void contract_simplex(double ** simplex, int dim, double ** fx, int ilo, int ihi, Search_settings *sett, Aux_arrays *aux, double *nS, double **sigaa, double **sigbb){ double sinalt, cosalt, sindelt, cosdelt; double * fx3; int i, j, k, o; for (i = 0; i < dim + 1; i++){ if (i != ilo){ for (j = 0; j < dim; j++) simplex[i][j] = (simplex[ilo][j]+simplex[i][j])*0.5; sinalt = sin(simplex[i][3]); cosalt = cos(simplex[i][3]); sindelt = sin(simplex[i][2]); cosdelt = cos(simplex[i][2]); nS[0] = cosalt*cosdelt; nS[1] = sinalt*cosdelt; nS[2] = sindelt; for (o = 0; o < sett->nifo; ++o){ modvir(sinalt, cosalt, sindelt, cosdelt, sett->N, &ifo[o], aux, sigaa[o], sigbb[o]); } fx3 = Fstatnet(sett, simplex[i], nS, sigaa, sigbb); for (k = 0; k < 11; k++) fx[i][k] = fx3[k]; } } }
int update_simplex(double ** simplex, int dim, double fmax, double ** fx, int ihi, double * midpoint, double * line, double scale, Search_settings *sett, Aux_arrays *aux, double *nS, double **sigaa, double **sigbb){ int i, o, j, update = 0; double * next = alloc_vector(dim); double * fx2; double sinalt, cosalt, sindelt, cosdelt; for (i = 0; i < dim; i++) next[i] = midpoint[i] + scale * line[i]; sinalt = sin(next[3]); cosalt = cos(next[3]); sindelt = sin(next[2]); cosdelt = cos(next[2]); nS[0] = cosalt*cosdelt; nS[1] = sinalt*cosdelt; nS[2] = sindelt; for (o = 0; o < sett->nifo; ++o){ modvir(sinalt, cosalt, sindelt, cosdelt, sett->N, &ifo[o], aux, sigaa[o], sigbb[o]); } fx2 = Fstatnet(sett, next, nS, sigaa, sigbb); if (fx2[5] < fmax){ for (i = 0; i < dim; i++) simplex[ihi][i] = next[i]; for (j = 0; j < 11; j++) fx[ihi][j] = fx2[j]; update = 1; } free_vector(next, dim); return update; }
void evaluate_simplex(double ** simplex, int dim, double ** fx, Search_settings *sett, Aux_arrays *aux, double *nS, double **sigaa, double **sigbb){ double sinalt, cosalt, sindelt, cosdelt; double *out; int i, o, j; for (i = 0; i < dim + 1; i++){ sinalt = sin(simplex[i][3]); cosalt = cos(simplex[i][3]); sindelt = sin(simplex[i][2]); cosdelt = cos(simplex[i][2]); nS[0] = cosalt*cosdelt; nS[1] = sinalt*cosdelt; nS[2] = sindelt; for (o = 0; o < sett->nifo; ++o){ modvir(sinalt, cosalt, sindelt, cosdelt, sett->N, &ifo[o], aux, sigaa[o], sigbb[o]); } out = Fstatnet(sett, simplex[i], nS, sigaa, sigbb); for (j = 0; j < 11; j++) fx[i][j] = out[j]; } }
int job_core(int pm, // Hemisphere int mm, // Grid 'sky position' int nn, // Second grid 'sky position' Search_settings *sett, // Search settings Command_line_opts *opts, // Search options Search_range *s_range, // Range for searching FFTW_plans *plans, // Plans for fftw FFTW_arrays *fftw_arr, // Arrays for fftw Aux_arrays *aux, // Auxiliary arrays int *sgnlc, // Candidate trigger parameters FLOAT_TYPE *sgnlv, // Candidate array int *FNum) { // Candidate signal number int i, j, n; int smin = s_range->sst, smax = s_range->spndr[1]; double al1, al2, sinalt, cosalt, sindelt, cosdelt, sgnlt[NPAR], nSource[3], het0, sgnl0, ft; double _tmp1[sett->nifo][sett->N]; #undef NORMTOMAX #ifdef NORMTOMAX double blkavg, threshold = 6.; int imax, imax0, iblk, blkstart, ihi; int blksize = 1024; int nfft = sett->nmax - sett->nmin; static int *Fmax; if (!Fmax) Fmax = (int *) malloc(nfft*sizeof(int)); #endif struct timespec tstart, tend; double spindown_timer = 0; int spindown_counter = 0; //tstart = get_current_time(CLOCK_REALTIME); /* Matrix M(.,.) (defined on page 22 of PolGrawCWAllSkyReview1.pdf file) defines the transformation form integers (bin, ss, nn, mm) determining a grid point to linear coordinates omega, omegadot, alpha_1, alpha_2), where bin is the frequency bin number and alpha_1 and alpha_2 are defined on p. 22 of PolGrawCWAllSkyReview1.pdf file. [omega] [bin] [omegadot] = M(.,.) \times [ss] [alpha_1/omega] [nn] [alpha_2/omega] [mm] Array M[.] is related to matrix M(.,.) in the following way; [ M[0] M[4] M[8] M[12] ] M(.,.) = [ M[1] M[5] M[9] M[13] ] [ M[2] M[6] M[10] M[14] ] [ M[3] M[7] M[11] M[15] ] and M[1] = M[2] = M[3] = M[6] = M[7] = 0 */ // Grid positions al1 = nn*sett->M[10] + mm*sett->M[14]; al2 = nn*sett->M[11] + mm*sett->M[15]; // check if the search is in an appropriate region of the grid // if not, returns NULL if ((sqr(al1)+sqr(al2))/sqr(sett->oms) > 1.) return 0; int ss; double shft1, phase, cp, sp; complex double exph; // Change linear (grid) coordinates to real coordinates lin2ast(al1/sett->oms, al2/sett->oms, pm, sett->sepsm, sett->cepsm, &sinalt, &cosalt, &sindelt, &cosdelt); // calculate declination and right ascention // written in file as candidate signal sky positions sgnlt[2] = asin(sindelt); sgnlt[3] = fmod(atan2(sinalt, cosalt) + 2.*M_PI, 2.*M_PI); het0 = fmod(nn*sett->M[8] + mm*sett->M[12], sett->M[0]); // Nyquist frequency int nyqst = (sett->nfft)/2 + 1; // Loop for each detector /* Amplitude modulation functions aa and bb * for each detector (in signal sub-struct * of _detector, ifo[n].sig.aa, ifo[n].sig.bb) */ for(n=0; n<sett->nifo; ++n) { modvir(sinalt, cosalt, sindelt, cosdelt, sett->N, &ifo[n], aux); // Calculate detector positions with respect to baricenter nSource[0] = cosalt*cosdelt; nSource[1] = sinalt*cosdelt; nSource[2] = sindelt; shft1 = nSource[0]*ifo[n].sig.DetSSB[0] + nSource[1]*ifo[n].sig.DetSSB[1] + nSource[2]*ifo[n].sig.DetSSB[2]; #define CHUNK 4 #pragma omp parallel default(shared) private(phase,cp,sp,exph) { #pragma omp for schedule(static,CHUNK) for(i=0; i<sett->N; ++i) { ifo[n].sig.shft[i] = nSource[0]*ifo[n].sig.DetSSB[i*3] + nSource[1]*ifo[n].sig.DetSSB[i*3+1] + nSource[2]*ifo[n].sig.DetSSB[i*3+2]; ifo[n].sig.shftf[i] = ifo[n].sig.shft[i] - shft1; _tmp1[n][i] = aux->t2[i] + (double)(2*i)*ifo[n].sig.shft[i]; } #pragma omp for schedule(static,CHUNK) for(i=0; i<sett->N; ++i) { // Phase modulation phase = het0*i + sett->oms*ifo[n].sig.shft[i]; #ifdef NOSINCOS cp = cos(phase); sp = sin(phase); #else sincos(phase, &sp, &cp); #endif exph = cp - I*sp; // Matched filter ifo[n].sig.xDatma[i] = ifo[n].sig.xDat[i]*ifo[n].sig.aa[i]*exph; ifo[n].sig.xDatmb[i] = ifo[n].sig.xDat[i]*ifo[n].sig.bb[i]*exph; } /* Resampling using spline interpolation: * This will double the sampling rate */ #pragma omp for schedule(static,CHUNK) for(i=0; i < sett->N; ++i) { fftw_arr->xa[i] = ifo[n].sig.xDatma[i]; fftw_arr->xb[i] = ifo[n].sig.xDatmb[i]; } // Zero-padding (filling with 0s up to sett->nfft, // the nearest power of 2) #pragma omp for schedule(static,CHUNK) for (i=sett->N; i<sett->nfft; ++i) { fftw_arr->xa[i] = 0.; fftw_arr->xb[i] = 0.; } } //omp parallel fftw_execute_dft(plans->pl_int,fftw_arr->xa,fftw_arr->xa); //forward fft (len nfft) fftw_execute_dft(plans->pl_int,fftw_arr->xb,fftw_arr->xb); //forward fft (len nfft) // move frequencies from second half of spectrum; // and zero frequencies higher than nyquist // loop length: nfft - nyqst = nfft - nfft/2 - 1 = nfft/2 - 1 for(i=nyqst + sett->Ninterp - sett->nfft, j=nyqst; i<sett->Ninterp; ++i, ++j) { fftw_arr->xa[i] = fftw_arr->xa[j]; fftw_arr->xa[j] = 0.; } for(i=nyqst + sett->Ninterp - sett->nfft, j=nyqst; i<sett->Ninterp; ++i, ++j) { fftw_arr->xb[i] = fftw_arr->xb[j]; fftw_arr->xb[j] = 0.; } // Backward fft (len Ninterp = nfft*interpftpad) fftw_execute_dft(plans->pl_inv,fftw_arr->xa,fftw_arr->xa); fftw_execute_dft(plans->pl_inv,fftw_arr->xb,fftw_arr->xb); ft = (double)sett->interpftpad / sett->Ninterp; //scale FFT for (i=0; i < sett->Ninterp; ++i) { fftw_arr->xa[i] *= ft; fftw_arr->xb[i] *= ft; } // struct timeval tstart = get_current_time(), tend; // Spline interpolation to xDatma, xDatmb arrays splintpad(fftw_arr->xa, ifo[n].sig.shftf, sett->N, sett->interpftpad, ifo[n].sig.xDatma); splintpad(fftw_arr->xb, ifo[n].sig.shftf, sett->N, sett->interpftpad, ifo[n].sig.xDatmb); } // end of detector loop // square sums of modulation factors double aa = 0., bb = 0.; for(n=0; n<sett->nifo; ++n) { double aatemp = 0., bbtemp = 0.; for(i=0; i<sett->N; ++i) { aatemp += sqr(ifo[n].sig.aa[i]); bbtemp += sqr(ifo[n].sig.bb[i]); } for(i=0; i<sett->N; ++i) { ifo[n].sig.xDatma[i] /= ifo[n].sig.sig2; ifo[n].sig.xDatmb[i] /= ifo[n].sig.sig2; } aa += aatemp/ifo[n].sig.sig2; bb += bbtemp/ifo[n].sig.sig2; } #ifdef YEPPP #define VLEN 1024 int bnd = (sett->N/VLEN)*VLEN; #endif // Check if the signal is added to the data or the range file is given: // if not, proceed with the wide range of spindowns // if yes, use smin = s_range->sst, smax = s_range->spndr[1] if(!strcmp(opts->addsig, "") && !strcmp(opts->range, "")) { // Spindown range defined using Smin and Smax (settings.c) smin = trunc((sett->Smin - nn*sett->M[9] - mm*sett->M[13])/sett->M[5]); smax = trunc(-(nn*sett->M[9] + mm*sett->M[13] + sett->Smax)/sett->M[5]); } if(opts->s0_flag) smin = smax; // if spindown parameter is taken into account, smin != smax printf ("\n>>%d\t%d\t%d\t[%d..%d]\n", *FNum, mm, nn, smin, smax); static fftw_complex *fxa, *fxb; static double *F; #pragma omp threadprivate(fxa,fxb, F) #pragma omp threadprivate(F) //private loop counter: ss //private (declared inside): ii,Fc,het1,k,veto_status,a,v,_p,_c,_s,status //shared default: nn,mm,sett,_tmp1,ifo,het0,bnd,plans,opts,aa,bb, // fftw_arr (zostawiamy i robimy nowe), FNum (atomic!) //we use shared plans and fftw_execute with 'new-array' interface #pragma omp parallel default(shared) \ private(i, j, n, sgnl0, exph, phase, cp, sp, tstart, tend) \ firstprivate(sgnlt) \ reduction(+ : spindown_timer, spindown_counter) { #ifdef YEPPP Yep64f _p[VLEN], _s[VLEN], _c[VLEN]; enum YepStatus status; #endif #ifdef SLEEF double _p[VECTLENDP], _c[VECTLENDP]; vdouble2 v; vdouble a; #endif if (!fxa) fxa = (fftw_complex *)fftw_malloc(fftw_arr->arr_len*sizeof(fftw_complex)); if (!fxb) fxb = (fftw_complex *)fftw_malloc(fftw_arr->arr_len*sizeof(fftw_complex)); if (!F) F = (double *)calloc(2*sett->nfft, sizeof(double)); /* Spindown loop */ #pragma omp for schedule(static,4) for(ss=smin; ss<=smax; ++ss) { #if TIMERS>2 tstart = get_current_time(CLOCK_PROCESS_CPUTIME_ID); #endif // Spindown parameter sgnlt[1] = ss*sett->M[5] + nn*sett->M[9] + mm*sett->M[13]; int ii; double Fc, het1; #ifdef VERBOSE //print a 'dot' every new spindown printf ("."); fflush (stdout); #endif het1 = fmod(ss*sett->M[4], sett->M[0]); if(het1<0) het1 += sett->M[0]; sgnl0 = het0 + het1; // phase modulation before fft #if defined(SLEEF) // use simd sincos from the SLEEF library; // VECTLENDP is a simd vector length defined in the SLEEF library // and it depends on selected instruction set e.g. -DENABLE_AVX for(i=0; i<sett->N; i+=VECTLENDP) { for(j=0; j<VECTLENDP; j++) _p[j] = het1*(i+j) + sgnlt[1]*_tmp1[0][i+j]; a = vloadu(_p); v = xsincos(a); vstoreu(_p, v.x); // reuse _p for sin vstoreu(_c, v.y); for(j=0; j<VECTLENDP; ++j){ exph = _c[j] - I*_p[j]; fxa[i+j] = ifo[0].sig.xDatma[i+j]*exph; //ifo[0].sig.sig2; fxb[i+j] = ifo[0].sig.xDatmb[i+j]*exph; //ifo[0].sig.sig2; } } #elif defined(YEPPP) // use yeppp! library; // VLEN is length of vector to be processed // for caches L1/L2 64/256kb optimal value is ~2048 for (j=0; j<bnd; j+=VLEN) { //double *_tmp2 = &_tmp1[0][j]; for (i=0; i<VLEN; ++i) //_p[i] = het1*(i+j) + sgnlt[1]*_tmp2[i]; _p[i] = het1*(i+j) + sgnlt[1]*_tmp1[0][i+j]; status = yepMath_Sin_V64f_V64f(_p, _s, VLEN); assert(status == YepStatusOk); status = yepMath_Cos_V64f_V64f(_p, _c, VLEN); assert(status == YepStatusOk); for (i=0; i<VLEN; ++i) { // exph = _c[i] - I*_s[i]; fxa[i+j] = ifo[0].sig.xDatma[i+j]*_c[i]-I*ifo[0].sig.xDatma[i+j]*_s[i]; fxb[i+j] = ifo[0].sig.xDatmb[i+j]*_c[i]-I*ifo[0].sig.xDatmb[i+j]*_s[i]; } } // remaining part is shorter than VLEN - no need to vectorize for (i=0; i<sett->N-bnd; ++i){ j = bnd + i; _p[i] = het1*j + sgnlt[1]*_tmp1[0][j]; } status = yepMath_Sin_V64f_V64f(_p, _s, sett->N-bnd); assert(status == YepStatusOk); status = yepMath_Cos_V64f_V64f(_p, _c, sett->N-bnd); assert(status == YepStatusOk); for (i=0; i<sett->N-bnd; ++i) { j = bnd + i; //exph = _c[i] - I*_s[i]; //fxa[j] = ifo[0].sig.xDatma[j]*exph; //fxb[j] = ifo[0].sig.xDatmb[j]*exph; fxa[j] = ifo[0].sig.xDatma[j]*_c[i]-I*ifo[0].sig.xDatma[j]*_s[i]; fxb[j] = ifo[0].sig.xDatmb[j]*_c[i]-I*ifo[0].sig.xDatmb[j]*_s[i]; } #elif defined(GNUSINCOS) for(i=sett->N-1; i!=-1; --i) { phase = het1*i + sgnlt[1]*_tmp1[0][i]; sincos(phase, &sp, &cp); exph = cp - I*sp; fxa[i] = ifo[0].sig.xDatma[i]*exph; //ifo[0].sig.sig2; fxb[i] = ifo[0].sig.xDatmb[i]*exph; //ifo[0].sig.sig2; } #else for(i=sett->N-1; i!=-1; --i) { phase = het1*i + sgnlt[1]*_tmp1[0][i]; cp = cos(phase); sp = sin(phase); exph = cp - I*sp; fxa[i] = ifo[0].sig.xDatma[i]*exph; //ifo[0].sig.sig2; fxb[i] = ifo[0].sig.xDatmb[i]*exph; //ifo[0].sig.sig2; } #endif for(n=1; n<sett->nifo; ++n) { #if defined(SLEEF) // use simd sincos from the SLEEF library; // VECTLENDP is a simd vector length defined in the SLEEF library // and it depends on selected instruction set e.g. -DENABLE_AVX for (i=0; i<sett->N; i+=VECTLENDP) { for(j=0; j<VECTLENDP; j++) _p[j] = het1*(i+j) + sgnlt[1]*_tmp1[n][i+j]; a = vloadu(_p); v = xsincos(a); vstoreu(_p, v.x); // reuse _p for sin vstoreu(_c, v.y); for(j=0; j<VECTLENDP; ++j){ exph = _c[j] - I*_p[j]; fxa[i+j] = ifo[n].sig.xDatma[i+j]*exph; fxb[i+j] = ifo[n].sig.xDatmb[i+j]*exph; } } #elif defined(YEPPP) // use yeppp! library; // VLEN is length of vector to be processed // for caches L1/L2 64/256kb optimal value is ~2048 for (j=0; j<bnd; j+=VLEN) { //double *_tmp2 = &_tmp1[n][j]; for (i=0; i<VLEN; ++i) // _p[i] = het1*(i+j) + sgnlt[1]*_tmp2[i]; _p[i] = het1*(j+i) + sgnlt[1]*_tmp1[n][j+i]; status = yepMath_Sin_V64f_V64f(_p, _s, VLEN); assert(status == YepStatusOk); status = yepMath_Cos_V64f_V64f(_p, _c, VLEN); assert(status == YepStatusOk); for (i=0; i<VLEN; ++i) { //exph = _c[i] - I*_s[i]; //fxa[j+i] += ifo[n].sig.xDatma[j+i]*exph; //fxb[j+i] += ifo[n].sig.xDatmb[j+i]*exph; fxa[i+j] += ifo[n].sig.xDatma[i+j]*_c[i]-I*ifo[n].sig.xDatma[i+j]*_s[i]; fxb[i+j] += ifo[n].sig.xDatmb[i+j]*_c[i]-I*ifo[n].sig.xDatmb[i+j]*_s[i]; } } // remaining part is shorter than VLEN - no need to vectorize for (i=0; i<sett->N-bnd; ++i){ j = bnd + i; _p[i] = het1*j + sgnlt[1]*_tmp1[n][j]; } status = yepMath_Sin_V64f_V64f(_p, _s, sett->N-bnd); assert(status == YepStatusOk); status = yepMath_Cos_V64f_V64f(_p, _c, sett->N-bnd); assert(status == YepStatusOk); for (i=0; i<sett->N-bnd; ++i) { j = bnd + i; //exph = _c[i] - I*_s[i]; //fxa[j] += ifo[n].sig.xDatma[j]*exph; //fxb[j] += ifo[n].sig.xDatmb[j]*exph; fxa[j] += ifo[n].sig.xDatma[j]*_c[i]-I*ifo[n].sig.xDatma[j]*_s[i]; fxb[j] += ifo[n].sig.xDatmb[j]*_c[i]-I*ifo[n].sig.xDatmb[j]*_s[i]; } #elif defined(GNUSINCOS) for(i=sett->N-1; i!=-1; --i) { phase = het1*i + sgnlt[1]*_tmp1[n][i]; sincos(phase, &sp, &cp); exph = cp - I*sp; fxa[i] += ifo[n].sig.xDatma[i]*exph; fxb[i] += ifo[n].sig.xDatmb[i]*exph; } #else for(i=sett->N-1; i!=-1; --i) { phase = het1*i + sgnlt[1]*_tmp1[n][i]; cp = cos(phase); sp = sin(phase); exph = cp - I*sp; fxa[i] += ifo[n].sig.xDatma[i]*exph; fxb[i] += ifo[n].sig.xDatmb[i]*exph; } #endif } // Zero-padding for(i = sett->fftpad*sett->nfft-1; i != sett->N-1; --i) fxa[i] = fxb[i] = 0.; fftw_execute_dft(plans->plan, fxa, fxa); fftw_execute_dft(plans->plan, fxb, fxb); // Computing F-statistic for (i=sett->nmin; i<sett->nmax; i++) { F[i] = (sqr(creal(fxa[i])) + sqr(cimag(fxa[i])))/aa + (sqr(creal(fxb[i])) + sqr(cimag(fxb[i])))/bb; } // for (i=sett->nmin; i<sett->nmax; i++) // F[i] += (sqr(creal(fxb[i])) + sqr(cimag(fxb[i])))/bb; #pragma omp atomic (*FNum)++; #if 0 FILE *f1 = fopen("fraw-1.dat", "w"); for(i=sett->nmin; i<sett->nmax; i++) fprintf(f1, "%d %lf %lf\n", i, F[i], 2.*M_PI*i/((double) sett->fftpad*sett->nfft) + sgnl0); fclose(f1); #endif #ifndef NORMTOMAX //#define NAVFSTAT 4096 // Normalize F-statistics if(!(opts->white_flag)) // if the noise is not white noise FStat(F + sett->nmin, sett->nmax - sett->nmin, NAVFSTAT, 0); // f1 = fopen("fnorm-4096-1.dat", "w"); //for(i=sett->nmin; i<sett->nmax; i++) //fprintf(f1, "%d %lf %lf\n", i, F[i], 2.*M_PI*i/((double) sett->fftpad*sett->nfft) + sgnl0); //fclose(f1); // exit(EXIT_SUCCESS); for(i=sett->nmin; i<sett->nmax; i++) { if ((Fc = F[i]) > opts->trl) { // if F-stat exceeds trl (critical value) // Find local maximum for neighboring signals ii = i; while (++i < sett->nmax && F[i] > opts->trl) { if(F[i] >= Fc) { ii = i; Fc = F[i]; } // if F[i] } // while i // Candidate signal frequency sgnlt[0] = 2.*M_PI*ii/((FLOAT_TYPE)sett->fftpad*sett->nfft) + sgnl0; // Signal-to-noise ratio sgnlt[4] = sqrt(2.*(Fc-sett->nd)); // Checking if signal is within a known instrumental line int k, veto_status = 0; for(k=0; k<sett->numlines_band; k++) if(sgnlt[0]>=sett->lines[k][0] && sgnlt[0]<=sett->lines[k][1]) { veto_status=1; break; } int _sgnlc; if(!veto_status) { /* #pragma omp critical { (*sgnlc)++; // increase found number // Add new parameters to output array for (j=0; j<NPAR; ++j) // save new parameters sgnlv[NPAR*(*sgnlc-1)+j] = (FLOAT_TYPE)sgnlt[j]; } */ #pragma omp atomic capture { (*sgnlc)++; // increase found number _sgnlc = *sgnlc; } // Add new parameters to output array for (j=0; j<NPAR; ++j) // save new parameters sgnlv[NPAR*(_sgnlc-1)+j] = (FLOAT_TYPE)sgnlt[j]; #ifdef VERBOSE printf ("\nSignal %d: %d %d %d %d %d snr=%.2f\n", *sgnlc, pm, mm, nn, ss, ii, sgnlt[4]); #endif } } // if Fc > trl } // for i #else // new version imax = -1; // find local maxima first //printf("nmin=%d nmax=%d nfft=%d nblocks=%d\n", sett->nmin, sett->nmax, nfft, nfft/blksize); for(iblk=0; iblk < nfft/blksize; ++iblk) { blkavg = 0.; blkstart = sett->nmin + iblk*blksize; // block start index in F // in case the last block is shorter than blksize, include its elements in the previous block if(iblk==(nfft/blksize-1)) {blksize = sett->nmax - blkstart;} imax0 = imax+1; // index of first maximum in current block //printf("\niblk=%d blkstart=%d blksize=%d imax0=%d\n", iblk, blkstart, blksize, imax0); for(i=1; i <= blksize; ++i) { // include first element of the next block ii = blkstart + i; if(ii < sett->nmax) {ihi=ii+1;} else {ihi = sett->nmax; /*printf("ihi=%d ii=%d\n", ihi, ii);*/}; if(F[ii] > F[ii-1] && F[ii] > F[ihi]) { blkavg += F[ii]; Fmax[++imax] = ii; ++i; // next element can't be maximum - skip it } } // i // now imax points to the last element of Fmax // normalize in blocks blkavg /= (double)(imax - imax0 + 1); for(i=imax0; i <= imax; ++i) F[Fmax[i]] /= blkavg; } // iblk //f1 = fopen("fmax.dat", "w"); //for(i=1; i < imax; i++) //fprintf(f1, "%d %lf \n", Fmax[i], F[Fmax[i]]); //fclose(f1); //exit(EXIT_SUCCESS); // apply threshold limit for(i=0; i <= imax; ++i){ //if(F[Fmax[i]] > opts->trl) { if(F[Fmax[i]] > threshold) { sgnlt[0] = 2.*M_PI*i/((FLOAT_TYPE)sett->fftpad*sett->nfft) + sgnl0; // Signal-to-noise ratio sgnlt[4] = sqrt(2.*(F[Fmax[i]] - sett->nd)); // Checking if signal is within a known instrumental line int k, veto_status=0; for(k=0; k<sett->numlines_band; k++) if(sgnlt[0]>=sett->lines[k][0] && sgnlt[0]<=sett->lines[k][1]) { veto_status=1; break; } if(!veto_status) { (*sgnlc)++; // increase number of found candidates // Add new parameters to buffer array for (j=0; j<NPAR; ++j) sgnlv[NPAR*(*sgnlc-1)+j] = (FLOAT_TYPE)sgnlt[j]; #ifdef VERBOSE printf ("\nSignal %d: %d %d %d %d %d snr=%.2f\n", *sgnlc, pm, mm, nn, ss, Fmax[i], sgnlt[4]); #endif } } } // i #endif // old/new version #if TIMERS>2 tend = get_current_time(CLOCK_PROCESS_CPUTIME_ID); spindown_timer += get_time_difference(tstart, tend); spindown_counter++; #endif } // for ss } // omp parallel #ifndef VERBOSE printf("Number of signals found: %d\n", *sgnlc); #endif // tend = get_current_time(CLOCK_REALTIME); //time_elapsed = get_time_difference(tstart, tend); //printf("Parallel part: %e ( per thread %e ) s\n", time_elapsed, time_elapsed/omp_get_max_threads()); #if TIMERS>2 printf("\nTotal spindown loop time: %e s, mean spindown cpu-time: %e s (%d runs)\n", spindown_timer, spindown_timer/spindown_counter, spindown_counter); #endif return 0; } // jobcore
void add_signal( Search_settings *sett, Command_line_opts *opts, Aux_arrays *aux_arr, Search_range *s_range) { int i, j, n, gsize, reffr, k; double snr, sum = 0., h0, cof, thsnr = 0; double sigma_noise = 1.0; double sinaadd, cosaadd, sindadd, cosdadd, phaseadd, shiftadd, signadd; double nSource[3], sgnlo[10], sgnlol[4]; double **sigaa, **sigbb; // aa[nifo][N] sigaa = (double **)malloc(sett->nifo*sizeof(double *)); for (k=0; k < sett->nifo; k++) sigaa[k] = (double *)calloc(sett->N, sizeof(double)); sigbb = (double **)malloc(sett->nifo*sizeof(double *)); for (k=0; k < sett->nifo; k++) sigbb[k] = (double *)calloc(sett->N, sizeof(double)); FILE *data; // Signal parameters are read if ((data=fopen (opts->addsig, "r")) != NULL) { // Fscanning for the GW amplitude, grid size, hemisphere // and the reference frame (for which the signal freq. is not spun-down/up) fscanf (data, "%le %d %d %d", &snr, &gsize, s_range->pmr, &reffr); // Fscanning signal parameters: f, fdot, delta, alpha (sgnlo[0], ..., sgnlo[3]) // four amplitudes sgnlo[4], ..., sgnlo[7] // (see sigen.c and Phys. Rev. D 82, 022005 2010, Eqs. 2.13a-d) // be1, be2 (sgnlo[8], sgnlo[9] to translate the sky position into grid position) for(i=0; i<10; i++) fscanf(data, "%le",i+sgnlo); fclose (data); } else { perror (opts->addsig); } // Search-specific parametrization of freq. // for the software injections // sgnlo[0]: frequency, sgnlo[1]: frequency. derivative //#mb For VSR1 reffr=67 sgnlo[0] += -2.*sgnlo[1]*(sett->N)*(reffr - opts->ident); // Check if the signal is in band if(sgnlo[0]<0) exit(171); // « else if (sgnlo[0]>M_PI) exit(187); // » cof = sett->oms + sgnlo[0]; for(i=0; i<2; i++) sgnlol[i] = sgnlo[i]; sgnlol[2] = sgnlo[8]*cof; sgnlol[3] = sgnlo[9]*cof; // solving a linear system in order to translate // sky position, frequency and spindown (sgnlo parameters) // into the position in the grid double *MM ; MM = (double *) calloc (16, sizeof (double)); for(i=0; i<16; i++) MM[i] = sett->M[i] ; gsl_vector *x = gsl_vector_alloc (4); int s; gsl_matrix_view m = gsl_matrix_view_array (MM, 4, 4); gsl_matrix_transpose (&m.matrix) ; gsl_vector_view b = gsl_vector_view_array (sgnlol, 4); gsl_permutation *p = gsl_permutation_alloc (4); gsl_linalg_LU_decomp (&m.matrix, p, &s); gsl_linalg_LU_solve (&m.matrix, p, &b.vector, x); s_range->spndr[0] = round(gsl_vector_get(x,1)); s_range->nr[0] = round(gsl_vector_get(x,2)); s_range->mr[0] = round(gsl_vector_get(x,3)); gsl_permutation_free (p); gsl_vector_free (x); free (MM); // Define the grid range in which the signal will be looked for s_range->spndr[1] = s_range->spndr[0] + gsize; s_range->spndr[0] -= gsize; s_range->nr[1] = s_range->nr[0] + gsize; s_range->nr[0] -= gsize; s_range->mr[1] = s_range->mr[0] + gsize; s_range->mr[0] -= gsize; s_range->pmr[1] = s_range->pmr[0]; printf("add_signal() - the following grid range is used\n"); printf("(spndr, nr, mr, pmr pairs): %d %d %d %d %d %d %d %d\n", \ s_range->spndr[0], s_range->spndr[1], s_range->nr[0], s_range->nr[1], s_range->mr[0], s_range->mr[1], s_range->pmr[0], s_range->pmr[1]); // sgnlo[2]: declination, sgnlo[3]: right ascension sindadd = sin(sgnlo[2]); cosdadd = cos(sgnlo[2]); sinaadd = sin(sgnlo[3]); cosaadd = cos(sgnlo[3]); // To keep coherent phase between time segments double phaseshift = sgnlo[0]*sett->N*(reffr - opts->ident) + sgnlo[1]*pow(sett->N*(reffr - opts->ident), 2); // Loop for each detector for(n=0; n<sett->nifo; n++) { modvir(sinaadd, cosaadd, sindadd, cosdadd, sett->N, &ifo[n], aux_arr, sigaa[n], sigbb[n]); nSource[0] = cosaadd*cosdadd; nSource[1] = sinaadd*cosdadd; nSource[2] = sindadd; // adding signal to data (point by point) for (i=0; i<sett->N; i++) { shiftadd = 0.; for (j=0; j<3; j++) shiftadd += nSource[j]*ifo[n].sig.DetSSB[i*3+j]; // Phase phaseadd = sgnlo[0]*i + sgnlo[1]*aux_arr->t2[i] + (cof + 2.*sgnlo[1]*i)*shiftadd - phaseshift; // The whole signal with 4 amplitudes and modulations signadd = sgnlo[4]*(sigaa[n][i])*cos(phaseadd) + sgnlo[6]*(sigaa[n][i])*sin(phaseadd) + sgnlo[5]*(sigbb[n][i])*cos(phaseadd) + sgnlo[7]*(sigbb[n][i])*sin(phaseadd); //#mb test printout //printf("%d %le\n", i + sett->N*(opts->ident-1), phaseadd); // Adding the signal to the data vector if(ifo[n].sig.xDat[i]) { ifo[n].sig.xDat[i] += signadd; sum += pow(signadd, 2.); // thsnr += pow(h0*signadd, 2.); } } // Write the data+signal to file /* FILE *dataout; char xxx[512]; sprintf(xxx, "%s/%03d/%s/xdatc_%03d_%04d%s.bin", opts->dtaprefix, opts->ident, ifo[n].name, opts->ident, opts->band, opts->label); printf("Flag 1: try to write xDat to fole %s\n", xxx); if((dataout = fopen(xxx, "wb")) != NULL) { fwrite(ifo[n].sig.xDat, sizeof(*ifo[n].sig.xDat), sett->N, dataout); } else{ printf("Problem with %s file!\n", xxx); } fclose(dataout); */ } h0 = (snr*sigma_noise)/(sqrt(sum)); for(n=0; n<sett->nifo; n++) { for (i=0; i<sett->N; i++) { ifo[n].sig.xDat[i] = ifo[n].sig.xDat[i]*h0; } } printf("%le %le\n", snr, h0); for (i = 0; i < sett->nifo; i++) free(sigaa[i]); free(sigaa); for (i = 0; i < sett->nifo; i++) free(sigbb[i]); free(sigbb); //exit(0); }
// Main programm int main (int argc, char *argv[]) { Search_settings sett; Command_line_opts opts; Search_range s_range; Aux_arrays aux_arr; double *F; // F-statistic array int i, j, r, c, a, b, g; int d, o, m, k; int bins = 2, ROW, dim = 4; // neighbourhood of point will be divide into defined number of bins double pc[4]; // % define neighbourhood around each parameter for initial grid double pc2[4]; // % define neighbourhood around each parameter for direct maximum search (MADS & Simplex) double tol = 1e-10; // double delta = 1e-5; // initial step in MADS function // double *results; // Vector with results from Fstatnet function // double *maximum; // True maximum of Fstat // double results_max[11]; double s1, s2, s3, s4; double sgnlo[4]; double **arr; // arr[ROW][COL], arrg[ROW][COL]; double nSource[3]; double sinalt, cosalt, sindelt, cosdelt; double F_min; char path[512]; double x, y; ROW = pow((bins+1),4); #ifdef YEPPP yepLibrary_Init(); Yep64f *results_max = (Yep64f*)malloc(sizeof(Yep64f)*11); Yep64f *results_first = (Yep64f*)malloc(sizeof(Yep64f)*11); Yep64f *results = (Yep64f*)malloc(sizeof(Yep64f)*11); Yep64f *maximum = (Yep64f*)malloc(sizeof(Yep64f)*11); // Yep64f *sgnlo = (Yep64f*)malloc(sizeof(Yep64f)*4); // Yep64f *nSource = (Yep64f*)malloc(sizeof(Yep64f)*3); Yep64f *mean = (Yep64f*)malloc(sizeof(Yep64f)*4); enum YepStatus status; #endif pc[0] = 0.015; pc[1] = 0.015; pc[2] = 0.015; pc[3] = 0.015; for (i = 0; i < 4; i++){ pc2[i] = 2*pc[i]/bins; } // Time tests double tdiff; clock_t tstart, tend; // Command line options handle_opts(&sett, &opts, argc, argv); // Output data handling /* struct stat buffer; if (stat(opts.prefix, &buffer) == -1) { if (errno == ENOENT) { // Output directory apparently does not exist, try to create one if(mkdir(opts.prefix, S_IRWXU | S_IRGRP | S_IXGRP | S_IROTH | S_IXOTH) == -1) { perror (opts.prefix); return 1; } } else { // can't access output directory perror (opts.prefix); return 1; } } */ sprintf(path, "%s/candidates.coi", opts.dtaprefix); //Glue function if(strlen(opts.glue)) { glue(&opts); sprintf(opts.dtaprefix, "./data_total"); sprintf(opts.dtaprefix, "%s/followup_total_data", opts.prefix); opts.ident = 000; } FILE *coi; int z; if ((coi = fopen(path, "r")) != NULL) { // while(!feof(coi)) { /* if(!fread(&w, sizeof(unsigned short int), 1, coi)) { break; } fread(&mean, sizeof(float), 5, coi); fread(&fra, sizeof(unsigned short int), w, coi); fread(&ops, sizeof(int), w, coi); if((fread(&mean, sizeof(float), 4, coi)) == 4){ */ while(fscanf(coi, "%le %le %le %le", &mean[0], &mean[1], &mean[2], &mean[3]) == 4){ //Time test // tstart = clock(); arr = matrix(ROW, 4); //Function neighbourhood - generating grid around point arr = neigh(mean, pc, bins); // Output data handling /* struct stat buffer; if (stat(opts.prefix, &buffer) == -1) { if (errno == ENOENT) { // Output directory apparently does not exist, try to create one if(mkdir(opts.prefix, S_IRWXU | S_IRGRP | S_IXGRP | S_IROTH | S_IXOTH) == -1) { perror (opts.prefix); return 1; } } else { // can't access output directory perror (opts.prefix); return 1; } } */ // Grid data if(strlen(opts.addsig)) { read_grid(&sett, &opts); } // Search settings search_settings(&sett); // Detector network settings detectors_settings(&sett, &opts); // Array initialization init_arrays(&sett, &opts, &aux_arr, &F); // Amplitude modulation functions for each detector for(i=0; i<sett.nifo; i++) rogcvir(&ifo[i]); // Adding signal from file if(strlen(opts.addsig)) { add_signal(&sett, &opts, &aux_arr, &s_range); } // Setting number of using threads (not required) omp_set_num_threads(1); results_max[5] = 0.; // ifo zostaje shared // ifo....shft i ifo....xdatm{a,b] prerobić na lokalne tablice w fstatnet // w regionie parallel wprowadzić tablice private aa i bb ; alokować i przekazywać je jako argumenty do fstatnet i amoeba // Main loop - over all parameters + parallelisation #pragma omp parallel default(shared) private(d, i, sgnlo, sinalt, cosalt, sindelt, cosdelt, nSource, results, maximum) { double **sigaa, **sigbb; // aa[nifo][N] sigaa = matrix(sett.nifo, sett.N); sigbb = matrix(sett.nifo, sett.N); #pragma omp for for (d = 0; d < ROW; ++d){ for (i = 0; i < 4; i++){ sgnlo[i] = arr[d][i]; // sgnlo[i] = mean[i]; } sinalt = sin(sgnlo[3]); cosalt = cos(sgnlo[3]); sindelt = sin(sgnlo[2]); cosdelt = cos(sgnlo[2]); nSource[0] = cosalt*cosdelt; nSource[1] = sinalt*cosdelt; nSource[2] = sindelt; for (i = 0; i < sett.nifo; ++i){ modvir(sinalt, cosalt, sindelt, cosdelt, sett.N, &ifo[i], &aux_arr, sigaa[i], sigbb[i]); } // F-statistic in given point results = Fstatnet(&sett, sgnlo, nSource, sigaa, sigbb); //printf("Fstatnet: %le %le %le %le %le %le\n", results[6], results[7], results[8], results[9], results[5], results[4]); #pragma omp critical if(results[5] < results_max[5]){ for (i = 0; i < 11; i++){ results_max[i] = results[i]; } } // Maximum search using simplex algorithm if(opts.simplex_flag){ // puts("Simplex"); maximum = amoeba(&sett, &aux_arr, sgnlo, nSource, results, dim, tol, pc2, sigaa, sigbb); printf("Amoeba: %le %le %le %le %le %le\n", maximum[6], maximum[7], maximum[8], maximum[9], maximum[5], maximum[4]); // Maximum value in points searching #pragma omp critical if(maximum[5] < results_max[5]){ for (i = 0; i < 11; i++){ results_max[i] = maximum[i]; } } } //simplex } // d - main outside loop free_matrix(sigaa, sett.nifo, sett.N); free_matrix(sigbb, sett.nifo, sett.N); } //pragma for(g = 0; g < 11; g++) results_first[g] = results_max[g]; // Maximum search using MADS algorithm if(opts.mads_flag) { // puts("MADS"); maximum = MADS(&sett, &aux_arr, results_max, mean, tol, pc2, bins); } //Time test // tend = clock(); // tdiff = (tend - tstart)/(double)CLOCKS_PER_SEC; printf("%le %le %le %le %le %le\n", results_max[6], results_max[7], results_max[8], results_max[9], results_max[5], results_max[4]); } // while fread coi // } } //if coi else { perror (path); return 1; } // Output information /* puts("**********************************************************************"); printf("*** Maximum value of F-statistic for grid is : (-)%.8le ***\n", -results_first[5]); printf("Sgnlo: %.8le %.8le %.8le %.8le\n", results_first[6], results_first[7], results_first[8], results_first[9]); printf("Amplitudes: %.8le %.8le %.8le %.8le\n", results_first[0], results_first[1], results_first[2], results_first[3]); printf("Signal-to-noise ratio: %.8le\n", results_first[4]); printf("Signal-to-noise ratio from estimated amplitudes (for h0 = 1): %.8le\n", results_first[10]); puts("**********************************************************************"); if((opts.mads_flag)||(opts.simplex_flag)){ printf("*** True maximum is : (-)%.8le ***\n", -maximum[5]); printf("Sgnlo for true maximum: %.8le %.8le %.8le %.8le\n", maximum[6], maximum[7], maximum[8], maximum[9]); printf("Amplitudes for true maximum: %.8le %.8le %.8le %.8le\n", maximum[0], maximum[1], maximum[2], maximum[3]); printf("Signal-to-noise ratio for true maximum: %.8le\n", maximum[4]); printf("Signal-to-noise ratio from estimated amplitudes (for h0 = 1) for true maximum: %.8le\n", maximum[10]); puts("**********************************************************************"); }*/ // Cleanup & memory free free(results_max); free(results_first); free(results); free(maximum); free(mean); free_matrix(arr, ROW, 4); cleanup_followup(&sett, &opts, &s_range, &aux_arr, F); return 0; }
//mesh adaptive direct search (MADS) maximum search declaration //double Yep64f* MADS(Search_settings *sett, Aux_arrays *aux, double* in, double *start, double delta, double *pc, int bins){ int i, j, k, l, m, n, o, r, a = 0; double sinalt, cosalt, sindelt, cosdelt; double param[4]; //initial size of mesh double smallest = 100; double **sigaa, **sigbb; // aa[nifo][N] sigaa = matrix(sett->nifo, sett->N); sigbb = matrix(sett->nifo, sett->N); for(r = 0; r < 4;r++){ param[r] = pc[r]/bins; } for(r = 0; r < 4; r++){ if(param[r] < smallest) smallest = param[r]; } #ifdef YEPPP yepLibrary_Init(); Yep64f *p = (Yep64f*)malloc(sizeof(Yep64f)*4); Yep64f *out = (Yep64f*)malloc(sizeof(Yep64f)*11); Yep64f *nSource = (Yep64f*)malloc(sizeof(Yep64f)*3); Yep64f *extr = (Yep64f*)malloc(sizeof(Yep64f)*11); Yep64f *res = (Yep64f*)malloc(sizeof(Yep64f)*11); enum YepStatus status; #endif // puts("MADS"); for(i = 0; i < 11; i++) extr[i] = in[i]; while(smallest >= delta){ //when to end computations k = 0; for(j = 0; j < 8; j++){ for(i = 0; i < 4; i++) p[i] = in[6+i]; if(j < 4){ p[k] = in[6+k] - start[k]*(param[k] - delta); k++; } else { k--; p[k] = in[6+k] + start[k]*(param[k] - delta); } sinalt = sin(p[3]); cosalt = cos(p[3]); sindelt = sin(p[2]); cosdelt = cos(p[2]); nSource[0] = cosalt*cosdelt; nSource[1] = sinalt*cosdelt; nSource[2] = sindelt; for (o = 0; o < sett->nifo; ++o){ modvir(sinalt, cosalt, sindelt, cosdelt, sett->N, &ifo[o], aux, sigaa[o], sigbb[o]); } res = Fstatnet(sett, p, nSource, sigaa, sigbb); //Fstat function for mesh points if (res[5] < extr[5]){ for (l = 0; l < 11; l++){ extr[l] = res[l]; } } } //j if (extr[5] == in[5]){ smallest = smallest - delta; for(r = 0; r < 4; r++){ param[r] = param[r] - delta; } } else{ for(m = 0; m < 4; m++) p[m] = extr[6+m]; for(m = 0; m < 11; m++) in[m] = extr[m]; } } // while for (n= 0; n < 11; n++){ out[n] = extr[n]; } free(p); free(nSource); free(extr); free(res); free_matrix(sigaa, sett->nifo, sett->N); free_matrix(sigbb, sett->nifo, sett->N); return out; }