int main(int argc, char *argv[]) { FILE *fftfile, *candfile = NULL, *psfile = NULL; char filenm[100], candnm[100], psfilenm[120]; float locpow, norm, powargr, powargi; float *powr, *spreadpow, *minizoompow, *freqs; fcomplex *data, *minifft, *minizoom, *spread; fcomplex *resp, *kernel; double T, dr, ftobinp; int ii, nbins, ncands, candnum, lofreq = 0, nzoom, numsumpow = 1; int numbetween, numkern, kern_half_width; binaryprops binprops; infodata idata; if (argc < 3 || argc > 6) { usage(); exit(1); } printf("\n\n"); printf(" Binary Candidate Display Routine\n"); printf(" by Scott M. Ransom\n\n"); /* Initialize the filenames: */ sprintf(filenm, "%s.fft", argv[1]); sprintf(candnm, "%s_bin.cand", argv[1]); /* Read the info file */ readinf(&idata, argv[1]); if (idata.object) { printf("Plotting a %s candidate from '%s'.\n", idata.object, filenm); } else { printf("Plotting a candidate from '%s'.\n", filenm); } T = idata.N * idata.dt; /* Open the FFT file and get its length */ fftfile = chkfopen(filenm, "rb"); nbins = chkfilelen(fftfile, sizeof(fcomplex)); /* Open the candidate file and get its length */ candfile = chkfopen(candnm, "rb"); ncands = chkfilelen(candfile, sizeof(binaryprops)); /* The candidate number to examine */ candnum = atoi(argv[2]); /* Check that candnum is in range */ if ((candnum < 1) || (candnum > ncands)) { printf("\nThe candidate number is out of range.\n\n"); exit(1); } /* The lowest freq present in the FFT file */ if (argc >= 4) { lofreq = atoi(argv[3]); if ((lofreq < 0) || (lofreq > nbins - 1)) { printf("\n'lofreq' is out of range.\n\n"); exit(1); } } /* Is the original FFT a sum of other FFTs with the amplitudes added */ /* in quadrature? (i.e. an incoherent sum) */ if (argc >= 5) { numsumpow = atoi(argv[4]); if (numsumpow < 1) { printf("\nNumber of summed powers must be at least one.\n\n"); exit(1); } } /* Initialize PGPLOT using Postscript if requested */ if ((argc == 6) && (!strcmp(argv[5], "ps"))) { sprintf(psfilenm, "%s_bin_cand_%d.ps", argv[1], candnum); cpgstart_ps(psfilenm, "landscape"); } else { cpgstart_x("landscape"); } /* Read the binary candidate */ chkfileseek(candfile, (long) (candnum - 1), sizeof(binaryprops), SEEK_SET); chkfread(&binprops, sizeof(binaryprops), 1, candfile); fclose(candfile); /* Output the binary candidate */ print_bin_candidate(&binprops, 2); /* Allocate some memory */ powr = gen_fvect(binprops.nfftbins); minifft = gen_cvect(binprops.nfftbins / 2); spread = gen_cvect(binprops.nfftbins); spreadpow = gen_fvect(binprops.nfftbins); nzoom = 2 * ZOOMFACT * ZOOMNEIGHBORS; minizoom = gen_cvect(nzoom); minizoompow = gen_fvect(nzoom); /* Allocate and initialize our interpolation kernel */ numbetween = 2; kern_half_width = r_resp_halfwidth(LOWACC); numkern = 2 * numbetween * kern_half_width; resp = gen_r_response(0.0, numbetween, numkern); kernel = gen_cvect(binprops.nfftbins); place_complex_kernel(resp, numkern, kernel, binprops.nfftbins); COMPLEXFFT(kernel, binprops.nfftbins, -1); /* Read the data from the FFT file */ data = read_fcomplex_file(fftfile, binprops.lowbin - lofreq, binprops.nfftbins); /* Turn the Fourier amplitudes into powers */ for (ii = 0; ii < binprops.nfftbins; ii++) powr[ii] = POWER(data[ii].r, data[ii].i); /* Chop the powers that are way above the median level */ prune_powers(powr, binprops.nfftbins, numsumpow); /* Perform the minifft */ memcpy((float *) minifft, powr, sizeof(float) * binprops.nfftbins); realfft((float *) minifft, binprops.nfftbins, -1); /* Calculate the normalization constant */ norm = sqrt((double) binprops.nfftbins * (double) numsumpow) / minifft[0].r; locpow = minifft[0].r / binprops.nfftbins; /* Divide the original power spectrum by the local power level */ for (ii = 0; ii < binprops.nfftbins; ii++) powr[ii] /= locpow; /* Now normalize the miniFFT */ minifft[0].r = 1.0; minifft[0].i = 1.0; for (ii = 1; ii < binprops.nfftbins / 2; ii++) { minifft[ii].r *= norm; minifft[ii].i *= norm; } /* Interpolate the minifft and convert to power spectrum */ corr_complex(minifft, binprops.nfftbins / 2, RAW, kernel, binprops.nfftbins, FFT, spread, binprops.nfftbins, kern_half_width, numbetween, kern_half_width, CORR); for (ii = 0; ii < binprops.nfftbins; ii++) spreadpow[ii] = POWER(spread[ii].r, spread[ii].i); /* Plot the initial data set */ freqs = gen_freqs(binprops.nfftbins, binprops.lowbin / T, 1.0 / T); xyline(binprops.nfftbins, freqs, powr, "Pulsar Frequency (hz)", "Power / Local Power", 1); vect_free(freqs); printf("The initial data set (with high power outliers removed):\n\n"); /* Plot the miniFFT */ freqs = gen_freqs(binprops.nfftbins, 0.0, T / (2 * binprops.nfftbins)); xyline(binprops.nfftbins, freqs, spreadpow, "Binary Period (sec)", "Normalized Power", 1); vect_free(freqs); printf("The miniFFT:\n\n"); /* Interpolate and plot the actual candidate peak */ ftobinp = T / binprops.nfftbins; freqs = gen_freqs(nzoom, (binprops.rdetect - ZOOMNEIGHBORS) * ftobinp, ftobinp / (double) ZOOMFACT); for (ii = 0; ii < nzoom; ii++) { dr = -ZOOMNEIGHBORS + (double) ii / ZOOMFACT; rz_interp(minifft, binprops.nfftbins / 2, binprops.rdetect + dr, 0.0, kern_half_width, &minizoom[ii]); minizoompow[ii] = POWER(minizoom[ii].r, minizoom[ii].i); } xyline(nzoom, freqs, minizoompow, "Binary Period (sec)", "Normalized Power", 1); vect_free(freqs); printf("The candidate itself:\n\n"); printf("Done.\n\n"); /* Cleanup */ cpgend(); vect_free(data); vect_free(powr); vect_free(resp); vect_free(kernel); vect_free(minifft); vect_free(spread); vect_free(spreadpow); vect_free(minizoom); vect_free(minizoompow); fclose(fftfile); if ((argc == 6) && (!strcmp(argv[5], "ps"))) { fclose(psfile); } return (0); }
ffdotpows *subharm_ffdot_plane(int numharm, int harmnum, double fullrlo, double fullrhi, subharminfo * shi, accelobs * obs) { int ii, lobin, hibin, numdata, nice_numdata, nrs, fftlen, binoffset; static int numrs_full = 0, numzs_full = 0; float powargr, powargi; double drlo, drhi, harm_fract; ffdotpows *ffdot; fcomplex *data, **result; presto_datainf datainf; if (numrs_full == 0) { if (numharm == 1 && harmnum == 1) { numrs_full = ACCEL_USELEN; numzs_full = shi->numkern; } else { printf("You must call subharm_ffdot_plane() with numharm=1 and\n"); printf("harnum=1 before you use other values! Exiting.\n\n"); exit(0); } } ffdot = (ffdotpows *) malloc(sizeof(ffdotpows)); /* Calculate and get the required amplitudes */ harm_fract = (double) harmnum / (double) numharm; drlo = calc_required_r(harm_fract, fullrlo); drhi = calc_required_r(harm_fract, fullrhi); ffdot->rlo = (int) floor(drlo); ffdot->zlo = calc_required_z(harm_fract, obs->zlo); /* Initialize the lookup indices */ if (numharm > 1) { double rr, subr; for (ii = 0; ii < numrs_full; ii++) { rr = fullrlo + ii * ACCEL_DR; subr = calc_required_r(harm_fract, rr); shi->rinds[ii] = index_from_r(subr, ffdot->rlo); } } ffdot->rinds = shi->rinds; ffdot->numrs = (int) ((ceil(drhi) - floor(drlo)) * ACCEL_RDR + DBLCORRECT) + 1; if (numharm == 1 && harmnum == 1) { ffdot->numrs = ACCEL_USELEN; } else { if (ffdot->numrs % ACCEL_RDR) { ffdot->numrs = (ffdot->numrs / ACCEL_RDR + 1) * ACCEL_RDR; } } ffdot->numzs = shi->numkern; binoffset = shi->kern[0].kern_half_width; fftlen = shi->kern[0].fftlen; lobin = ffdot->rlo - binoffset; hibin = (int) ceil(drhi) + binoffset; numdata = hibin - lobin + 1; nice_numdata = next2_to_n(numdata); // for FFTs data = get_fourier_amplitudes(lobin, nice_numdata, obs); if (!obs->mmap_file && !obs->dat_input && 0) printf("This is newly malloc'd!\n"); // Normalize the Fourier amplitudes if (obs->nph > 0.0) { // Use freq 0 normalization if requested (i.e. photons) double norm = 1.0 / sqrt(obs->nph); for (ii = 0; ii < numdata; ii++) { data[ii].r *= norm; data[ii].i *= norm; } } else if (obs->norm_type == 0) { // old-style block median normalization float *powers; double norm; powers = gen_fvect(numdata); for (ii = 0; ii < numdata; ii++) powers[ii] = POWER(data[ii].r, data[ii].i); norm = 1.0 / sqrt(median(powers, numdata)/log(2.0)); free(powers); for (ii = 0; ii < numdata; ii++) { data[ii].r *= norm; data[ii].i *= norm; } } else { // new-style running double-tophat local-power normalization float *powers, *loc_powers; powers = gen_fvect(nice_numdata); for (ii = 0; ii < nice_numdata; ii++) { powers[ii] = POWER(data[ii].r, data[ii].i); } loc_powers = corr_loc_pow(powers, nice_numdata); for (ii = 0; ii < numdata; ii++) { float norm = invsqrt(loc_powers[ii]); data[ii].r *= norm; data[ii].i *= norm; } free(powers); free(loc_powers); } /* Perform the correlations */ result = gen_cmatrix(ffdot->numzs, ffdot->numrs); datainf = RAW; for (ii = 0; ii < ffdot->numzs; ii++) { nrs = corr_complex(data, numdata, datainf, shi->kern[ii].data, fftlen, FFT, result[ii], ffdot->numrs, binoffset, ACCEL_NUMBETWEEN, binoffset, CORR); datainf = SAME; } // Always free data free(data); /* Convert the amplitudes to normalized powers */ ffdot->powers = gen_fmatrix(ffdot->numzs, ffdot->numrs); for (ii = 0; ii < (ffdot->numzs * ffdot->numrs); ii++) ffdot->powers[0][ii] = POWER(result[0][ii].r, result[0][ii].i); free(result[0]); free(result); return ffdot; }
static void process_bird(double basebin, int harm, double *lofreq, double *hifreq) { int ii, plotnumpts = 1000, not_done_yet = 1, plotoffset; int lodatabin, firstcorrbin, numgoodpts, replot = 1; char inchar; float med, xx[2], yy[2], inx, iny; float powargr, powargi, pwr, maxpow = 0.0, maxbin = 0.0; double truebin, pred_freq, average; double firstbin = 0.0, lastbin = 0.0, numbins = 0.0; fcomplex *data, *result; /* 'bin' means normal resolution FFT amplitude */ /* 'pt' means an interpolated FFT amplitude */ /* 'pts'=='bins' only if NUMBETWEEN==1 */ *lofreq = *hifreq = 0.0; truebin = basebin * harm; pred_freq = truebin / T; xx[0] = xx[1] = pred_freq; data = get_rawbins(fftfile, truebin, BINSTOGET, &med, &lodatabin); if (lodatabin <= 0) { data[abs(lodatabin)].r = 1.0; data[abs(lodatabin)].i = 1.0; } firstcorrbin = (int) truebin - MAXBINSTOSHOW / 2; average = med / -log(0.5); result = gen_cvect(FFTLEN); numgoodpts = corr_complex(data, BINSTOGET, RAW, kernel, FFTLEN, FFT, result, MAXPTSTOSHOW, firstcorrbin - lodatabin, NUMBETWEEN, khw, CORR); for (ii = 0; ii < numgoodpts; ii++) { pwr = POWER(result[ii].r, result[ii].i) / average; if (pwr > maxpow) { maxpow = pwr; maxbin = firstcorrbin + dr * ii; } } printf("\nHarmonic %d of %.15g Hz (%.15g Hz, bin = %.15g)\n", harm, basebin / T, pred_freq, truebin); printf(" Max power = %.2f at %.15g Hz (bin = %.15g)\n", maxpow, maxbin / T, maxbin); do { cpgsci(1); if (replot) { plotoffset = MAXPTSTOSHOW / 2 - plotnumpts / 2; firstbin = firstcorrbin + dr * plotoffset; numbins = dr * plotnumpts; lastbin = firstbin + numbins; plot_spectrum(result + plotoffset, plotnumpts, firstbin, dr, T, average); cpgswin(0.0, 1.0, 0.0, 1.0); xx[0] = xx[1] = (truebin - firstbin) / numbins; yy[0] = 0.0; yy[1] = 1.0; cpgsci(2); /* red */ cpgline(2, xx, yy); /* Predicted freq */ cpgsci(7); /* yellow */ if (*lofreq) { xx[0] = xx[1] = ((*lofreq * T) - firstbin) / numbins; cpgline(2, xx, yy); /* Boundary */ } if (*hifreq) { xx[0] = xx[1] = ((*hifreq * T) - firstbin) / numbins; cpgline(2, xx, yy); /* Boundary */ } } replot = 1; cpgsci(7); /* yellow */ cpgcurs(&inx, &iny, &inchar); switch (inchar) { case ' ': case 'A': case 'a': xx[0] = xx[1] = inx; cpgline(2, xx, yy); /* New boundary */ if (*lofreq == 0.0) { *lofreq = (inx * numbins + firstbin) / T; printf(" Added 1st boundary at %.12g Hz\n", *lofreq); } else { *hifreq = (inx * numbins + firstbin) / T; printf(" Added 2nd boundary at %.12g Hz\n", *hifreq); } replot = 0; break; case 'I': /* Zoom in */ case 'i': plotnumpts /= 2; if (plotnumpts <= 8) plotnumpts = 8; printf(" Zooming in...\n"); break; case 'O': /* Zoom out */ case 'o': plotnumpts *= 2; if (plotnumpts > MAXPTSTOSHOW) plotnumpts = MAXPTSTOSHOW; printf(" Zooming out...\n"); break; case 'C': /* Clear/Delete the points */ case 'c': case 'D': case 'd': case 'X': case 'x': *lofreq = *hifreq = 0.0; printf(" Clearing boundaries.\n"); break; case 'Q': /* Quit/Next birdie */ case 'q': case 'N': case 'n': *lofreq = *hifreq = 0.0; free(data); vect_free(result); printf(" Skipping to next harmonic.\n"); return; default: printf(" Unrecognized option '%c'.\n", inchar); break; } if (*lofreq && *hifreq) not_done_yet = 0; } while (not_done_yet); if (*hifreq < *lofreq) { double tmpfreq; tmpfreq = *lofreq; *lofreq = *hifreq; *hifreq = tmpfreq; } free(data); vect_free(result); }