Esempio n. 1
0
TEST(testEndToEnd,gracefulBadFit){
    pulsar _psr;
    pulsar *psr = &_psr;
    MAX_PSR=1;
    char timFile[MAX_PSR][MAX_FILELEN],parFile[MAX_PSR][MAX_FILELEN];
    char** argv=NULL;
    int npsr=1;
    initialise(psr,0); /* Initialise all */

    strcpy(parFile[0],DATDIR "/test3c.par");
    strcpy(timFile[0],DATDIR "/test3.tim");

    readParfile(psr,parFile,timFile,npsr);
    readTimfile(psr,timFile,npsr);
    preProcessSimple(psr);

    psr->noWarnings=2;
    formBatsAll(psr,npsr);
    formResiduals(psr,npsr,0);
    TK_errorCount=0;
    doFitAll(psr,npsr,"NULL");
//    ASSERT_GT(TK_errorCount,1u);
ASSERT_EQ(1,1);
    TK_errorCount=0;
    
}
Esempio n. 2
0
TEST(testEndToEnd,checkIdeal){
    pulsar _psr;
    pulsar *psr = &_psr;
    MAX_PSR=1;
    char timFile[MAX_PSR][MAX_FILELEN],parFile[MAX_PSR][MAX_FILELEN];
    char** argv=NULL;
    int npsr=1;
    initialise(psr,0); /* Initialise all */

    strcpy(parFile[0],DATDIR "/test3.par");
    strcpy(timFile[0],DATDIR "/test3.tim");

    readParfile(psr,parFile,timFile,npsr);
    readTimfile(psr,timFile,npsr);
    preProcessSimple(psr);

    psr->noWarnings=2;
    formBatsAll(psr,npsr);
    formResiduals(psr,npsr,0);
    for(int iobs = 0; iobs < psr->nobs; iobs++){
        ASSERT_LT(static_cast<double>(fabsl(psr->obsn[iobs].residual)),TEST_DELTA) << DATDIR "/test3.par test3.tim do not give idealised ToAs";
    }

    doFitAll(psr,npsr,"NULL");

    formBatsAll(psr,npsr);
    formResiduals(psr,npsr,0);
    for(int iobs = 0; iobs < psr->nobs; iobs++){
        ASSERT_LT(static_cast<double>(fabsl(psr->obsn[iobs].residual)),TEST_DELTA) << "Fitting has caused error in ideal";
    }

}
Esempio n. 3
0
/* The main function has to be called this, even if no
 * graphics are involved ;)
 */
extern "C" int graphicalInterface(int argc, char *argv[], 
        pulsar *psr, int *npsr) 
{
    char parFile[MAX_PSR][MAX_FILELEN];
    char timFile[MAX_PSR][MAX_FILELEN];
    *npsr = 1; // Necessary, I guess?

    printf("Chi2 grid plugin for tempo2\n");
    printf("Author:  Paul Demorest\n");
    printf("Version:           0.0\n");

    parFile[0][0]='\0';
    timFile[0][0]='\0';

    // Parse command line.. there has to be a better way to do this..
    double m2_min = 0.0, m2_max = 2.0;
    double sini_min = 0.5, sini_max = 1.0;
    int m2_steps = 128, sini_steps = 128;
    for (unsigned i=1; i<argc; i++) {
        if (strcmp(argv[i], "-f")==0) {
            strcpy(parFile[0], argv[i+1]);
            strcpy(timFile[0], argv[i+2]);
        } else if (strcmp(argv[i], "-sini_min")==0) {
            sini_min = atof(argv[i+1]);
        } else if (strcmp(argv[i], "-sini_max")==0) {
            sini_max = atof(argv[i+1]);
        } else if (strcmp(argv[i], "-sini_steps")==0) {
            sini_steps = atoi(argv[i+1]);
        } else if (strcmp(argv[i], "-m2_min")==0) {
            m2_min = atof(argv[i+1]);
        } else if (strcmp(argv[i], "-m2_max")==0) {
            m2_max = atof(argv[i+1]);
        } else if (strcmp(argv[i], "-m2_steps")==0) {
            m2_steps = atoi(argv[i+1]);
        }
    }

    if (parFile[0][0]=='\0' || timFile[0][0]=='\0') {
        fprintf(stderr, "Please provide .par and .tim files using -f\n");
        fprintf(stderr, "Optional arguments:\n");
        fprintf(stderr, "  -sini_min N    Minimum sin(i) value (0.5)\n");
        fprintf(stderr, "  -sini_max N    Maximum sin(i) value (1.0)\n");
        fprintf(stderr, "  -sini_steps N  Number of sin(i) steps (128)\n");
        fprintf(stderr, "  -m2_min N      Minimum m2 value (0.0)\n");
        fprintf(stderr, "  -m2_max N      Maximum m2 value (2.0)\n");
        fprintf(stderr, "  -m2_steps N    Number of m2 steps (128)\n");
        exit(1);
    }

    // Read the stuff
    readParfile(psr, parFile, timFile, *npsr);
    readTimfile(psr, timFile, *npsr);
    preProcess(psr, *npsr, argc, argv);

    // Do an initial fit to get best-fit params to use as a starting
    // place.
    formBatsAll(psr, *npsr);
    formResiduals(psr, *npsr, 0); // Final arg is "removeMean"
    doFit(psr, *npsr, 0);         // Final arg is "writeModel"
    formBatsAll(psr, *npsr);      // Run these again to get post-fit
    formResiduals(psr, *npsr, 0); // residuals.

    // Get initial chi2, etc.
    double init_chi2 = psr_chi2(psr);

    // step size, etc
    double delta_m2 = (m2_max - m2_min)/(double)m2_steps;
    double delta_sini = (sini_max - sini_min)/(double)sini_steps;

    // Store current parameters in psr[1], trial params in psr[0]
    store_params(psr);

    // Print header lines
    printf("# chi2_0=%.5e\n", init_chi2);
    printf("# SINI M2 chi2-chi2_0\n");

    // Main iteration loop
    signal(SIGINT, cc);
    double cur_m2, cur_sini, cur_chi2;
    for (cur_sini=sini_min; cur_sini<sini_max; cur_sini+=delta_sini) {
        for (cur_m2=m2_min; cur_m2<m2_max; cur_m2+=delta_m2) {

            reset_params(psr);

            psr[0].param[param_sini].val[0] = cur_sini;
            psr[0].param[param_sini].paramSet[0] = 1;
            psr[0].param[param_m2].val[0] = cur_m2;
            psr[0].param[param_m2].paramSet[0] = 1;

            // updateBats or FormBats?
            int ok=1;
            updateBatsAll(psr, *npsr);
            try {
                formResiduals(psr, *npsr, 0);
                doFit(psr, *npsr, 0);
            } catch (int e) {
                // Note, this is a hack to ignore spots in 
                // M2-sini parameter space where the fit does not
                // converge.  Plain tempo2 exit()'s here, killing
                // the whole process.  I've modified my copy to
                // throw an error instead so we can continue on
                // with the gridding.
                if (e==5) { ok=0; }
            }
            if (ok) {

                updateBatsAll(psr, *npsr);
                formResiduals(psr, *npsr, 0);
                cur_chi2 = psr_chi2(psr);

            } else {
                cur_chi2 = 1e6;
            }

            printf("%.9f %.5f %.10e\n", cur_sini, cur_m2, cur_chi2-init_chi2);

            if (run==0) break;
        }
        if (run==0) break;
    }

}
Esempio n. 4
0
extern "C" int graphicalInterface(int argc,char *argv[],pulsar *psr,int *npsr) 
{
  char parFile[MAX_PSR][MAX_FILELEN];
  char timFile[MAX_PSR][MAX_FILELEN];
  double epochs[MAX_OBSN];
  int i,nit,j,p;
  char fname[MAX_FILELEN];
  double globalParameter;
  long double result;
  long seed = TKsetSeed();
  
  //
  // For the output file
  //
  toasim_header_t* header;
  toasim_header_t* read_header;
  FILE* file;
  double offsets[MAX_OBSN]; // Will change to doubles - should use malloc
  double mean=0.0;
  // Create a set of corrections.
  toasim_corrections_t* corr = (toasim_corrections_t*)malloc(sizeof(toasim_corrections_t));
  
  double f1_1;
  double f1_2;
  double t0 = 0;
  double t;
  double phaseStart=0;
  double phase0=0;
  double f0;
  double phase;
  double nudot;
  double nu;
  int nStep;
  double stepTime[MAX_STEP];
  int sw;
  int nSwitch;
  double nu0,initialT,nudot_now;
  double nuSwitch[MAX_STEP+1];
  double fa,fb;
  double deltaT;
  int ii;

  nStep = 0;


  corr->offsets=offsets;
  corr->params=""; // Normally leave as NULL. Can store this along with each realisation. 
  // Same length string in every iteration - defined in r_param_length see below
  corr->a0=0; // constant
  corr->a1=0; // a1*x
  corr->a2=0; // a2*x*X
  
  *npsr = 0;
  nit = 1;

  printf("Graphical Interface: add2state\n");
  printf("Author:              G. Hobbs\n");
  printf("Version:             1.0\n");
  
  /* Obtain all parameters from the command line */
  for (i=2;i<argc;i++)
    {
      
      if (strcmp(argv[i],"-nreal")==0){
	nit=atoi(argv[++i]);
      }
      if (strcmp(argv[i],"-step")==0){
	sscanf(argv[++i],"%lf",&stepTime[nStep]);
	nStep++;
      }
      if (strcmp(argv[i],"-fa")==0){
	sscanf(argv[++i],"%lf",&f1_1);
      }
      if (strcmp(argv[i],"-fb")==0){
	sscanf(argv[++i],"%lf",&f1_2);
      }
      if (strcmp(argv[i],"-f")==0)
	{
	  strcpy(parFile[*npsr],argv[++i]); 
	  strcpy(timFile[*npsr],argv[++i]);
	  (*npsr)++;
	}
      if (strcmp(argv[i],"-seed")==0){
	sscanf(argv[++i],"%d",&seed);
      }
    }
  
  if (seed > 0)seed=-seed;
  
  readParfile(psr,parFile,timFile,*npsr); /* Load the parameters       */
  // Now read in all the .tim files
  readTimfile(psr,timFile,*npsr); /* Load the arrival times    */
  
  preProcess(psr,*npsr,argc,argv);
  
  for (p=0;p<*npsr;p++)
    {
      f0 = (double)psr[p].param[param_f].val[0];


      printf("NTOA = %d\n",psr[p].nobs);
      header = toasim_init_header();
      strcpy(header->short_desc,"add2state");
      strcpy(header->invocation,argv[0]);
      strcpy(header->timfile_name,timFile[p]);
      strcpy(header->parfile_name,"Unknown");
      header->idealised_toas="NotSet"; // What should this be
      header->orig_parfile="NA";
      header->gparam_desc=""; // Global parameters
      header->gparam_vals="";
      header->rparam_desc=""; // Desciprtion of the parameters
      header->rparam_len=0; // Size of the string
      header->seed = seed;
		
      header->ntoa = psr[p].nobs;
      header->nrealisations = nit;
      
      // First we write the header...
      sprintf(fname,"%s.add2state",timFile[p]);
      file = toasim_write_header(header,fname);
      
      
      for (i=0;i<nit;i++)
	{
	  mean=0;
	  if(i%10 == 0){
	    printf("\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b");
	    printf("\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b");
	    printf("Iteration %d/%d",i+1,nit);
	    fflush(stdout);
	  }
	  
	  t0 = (double)psr[p].obsn[0].sat;

	  for (j=0;j<psr[p].nobs;j++)
	    {
	      t= (psr[p].obsn[j].sat-t0)*86400.0;
	      // Calculate nudot at time t
	      sw=1;
	      for (ii=0;ii<nStep;ii++)
		{
		  if (psr[p].obsn[j].sat > stepTime[ii])
		    sw*=-1;
		}
	      if (sw == 1)
		nudot = f1_1;
	      else 
		nudot = f1_2;
	      
	      // Calculate nu at time t
	      // Integrate from time t0 to the current time
	      // Find necessary number of switches
	      nSwitch=0;
	      for (ii=0;ii<nStep;ii++)
		{
		  if (stepTime[ii] < psr[p].obsn[j].sat) {nSwitch++;}
		}
	      nu0 = f0;
	      initialT = t0;
	      nudot_now = f1_1;
	      sw=1;
	      nuSwitch[0] = nu0;
	      for (ii=0;ii<nSwitch;ii++)
		{
		  nu0 = nu0 + nudot_now*(stepTime[ii]-initialT)*86400.0;
		  nuSwitch[ii+1] = nu0;
		  sw*=-1;
		  if (sw==1)
		    nudot_now = f1_1;
		  else
		    nudot_now = f1_2;
		  initialT = stepTime[ii];
		}
	      // Now do the last bit
	      nu = nu0 + nudot_now*(psr[p].obsn[j].sat - initialT)*86400.0;
	      
	      // Now calculate phase at time t
	      // Dealing with all the switches
	      phase = phaseStart;
	      fb = f0;
	      initialT = t0;
	      nudot_now = f1_1;
	      sw=1;
	      
	      for (ii=0;ii<nSwitch;ii++)
		{
		  fa = nuSwitch[ii];
		  fb = nuSwitch[ii+1];
		  if (ii>0)
		    deltaT = stepTime[ii]-stepTime[ii-1];
		  else
		    deltaT = stepTime[ii] - t0;
		  sw*=-1;
		  if (sw==1)
		    nudot_now = f1_1;
		  else
		    nudot_now = f1_2;
		  
		  phase += (deltaT*86400.0*fb)+0.5*deltaT*86400.0*fabs(fa-fb);
		  initialT = stepTime[ii];
		}
	      // Now do last bit
	      deltaT = (psr[p].obsn[j].sat-initialT)*86400.0;
	      fa = fb;
	      fb = fa+deltaT*nudot_now;
	      phase += fb*deltaT+0.5*deltaT*fabs(fa-fb);
	      //	      phase = phase - floor(phase);
	      
	      printf("initial %g %g %d %.15f %.15f %g\n",(double)psr[p].obsn[j].sat,nudot,nSwitch,nu,phase,f0);
	      epochs[j] = (double)(psr[p].obsn[j].sat-psr[p].param[param_pepoch].val[0]);
	      offsets[j] = (double)(phase);
	      mean+=offsets[j];
	    }
	  for (j=0;j<psr[p].nobs;j++)
	    offsets[j]-=mean/(double)psr[p].nobs;
	  TKremovePoly_d(epochs,offsets,psr[p].nobs,3);

	  for (j=0;j<psr[p].nobs;j++)
	    printf("value = %g %.15f\n",epochs[j],offsets[j]);



  //			  //			  printf("Error = %g\n",(double)psr[p].obsn[j].toaErr);
  //			}
	    toasim_write_corrections(corr,header,file);

	  printf("\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b");
	  printf("\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b");
	  printf("Iteration %d/%d\n",i,nit);
	}	  
      printf("Close file\n");
      fclose(file);
    }
  return 0;
}
Esempio n. 5
0
extern "C" int graphicalInterface(int argc,char *argv[],pulsar *psr,int *npsr) 
{
	char parFile[MAX_PSR][MAX_FILELEN];
	char timFile[MAX_PSR][MAX_FILELEN];
	int i,nit,j,p;
	char fname[MAX_FILELEN];
	double globalParameter;
	long double result;
	long seed = TKsetSeed();

	//
	// For the output file
	//
	toasim_header_t* header;
	toasim_header_t* read_header;
	FILE* file;
	double offsets[MAX_OBSN]; // Will change to doubles - should use malloc
	// Create a set of corrections.
	toasim_corrections_t* corr = (toasim_corrections_t*)malloc(sizeof(toasim_corrections_t));

	corr->offsets=offsets;
	corr->params=""; // Normally leave as NULL. Can store this along with each realisation. 
	// Same length string in every iteration - defined in r_param_length see below
	corr->a0=0; // constant
	corr->a1=0; // a1*x
	corr->a2=0; // a2*x*X

	*npsr = 0;
	nit = 1;

	printf("Graphical Interface: addGaussian\n");
	printf("Author:              G. Hobbs, M. Keith\n");
	printf("Version:             1.0\n");

	/* Obtain all parameters from the command line */
	for (i=2;i<argc;i++)
	{

		if (strcmp(argv[i],"-nreal")==0){
			nit=atoi(argv[++i]);
		}
		if (strcmp(argv[i],"-f")==0)
		{
			strcpy(parFile[*npsr],argv[++i]); 
			strcpy(timFile[*npsr],argv[++i]);
			(*npsr)++;
		}
		if (strcmp(argv[i],"-seed")==0){
			sscanf(argv[++i],"%d",&seed);
		}
	}

	if (seed > 0)seed=-seed;

	readParfile(psr,parFile,timFile,*npsr); /* Load the parameters       */
	// Now read in all the .tim files
	readTimfile(psr,timFile,*npsr); /* Load the arrival times    */

	preProcess(psr,*npsr,argc,argv);

	for (p=0;p<*npsr;p++)
	{
		printf("NTOA = %d\n",psr[p].nobs);
		header = toasim_init_header();
		strcpy(header->short_desc,"addGaussian");
		strcpy(header->invocation,argv[0]);
		strcpy(header->timfile_name,timFile[p]);
		strcpy(header->parfile_name,"Unknown");
		header->idealised_toas="NotSet"; // What should this be
		header->orig_parfile="NA";
		header->gparam_desc=""; // Global parameters
		header->gparam_vals="";
		header->rparam_desc=""; // Desciprtion of the parameters
		header->rparam_len=0; // Size of the string
		header->seed = seed;

		header->ntoa = psr[p].nobs;
		header->nrealisations = nit;

		// First we write the header...
		sprintf(fname,"%s.addGauss",timFile[p]);
		file = toasim_write_header(header,fname);


		for (i=0;i<nit;i++)
		{
			if(i%10 == 0){
				printf("\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b");
				printf("\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b");
				printf("Iteration %d/%d",i+1,nit);
				fflush(stdout);
			}
			for (j=0;j<psr[p].nobs;j++){
					offsets[j] = (double)(psr[p].obsn[j].toaErr*1.0e-6*TKgaussDev(&seed));
			}
			toasim_write_corrections(corr,header,file);
		}
		printf("\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b");
		printf("\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b\b");
		printf("Iteration %d/%d\n",i,nit);

		printf("Close file\n");
		fclose(file);
	}
	return 0;
}
Esempio n. 6
0
extern "C" int graphicalInterface(int argc,char *argv[],pulsar *psr,int *npsr) 
{
  char parFile[MAX_PSR][MAX_FILELEN];
  char oparFile[MAX_PSR][MAX_FILELEN]; // Overlay files
  char timFile[MAX_PSR][MAX_FILELEN];
  char otimFile[MAX_PSR][MAX_FILELEN];
  char tpar[MAX_PSR][MAX_FILELEN];
  char ttim[MAX_PSR][MAX_FILELEN];
  int i;
  int plotUs=0;
  float scale[1000];
  int nScale=0;
  char grDev[100]="/vps";
  //  char grDev[100]="/xs";
  FILE *pin;
  FILE *fin;
  char str[1000];
  float fontsize = 1.2;
  float centreMJD=54500;
  int ptStyle=16;
  float ptSize = 1;
  int reverse=0;
  int error=1;
  float miny=-1;
  float maxy=-1;
  float minx=-1;
  float maxx=-1;
  int nOverlay=0;
  float labelsize=-1;
  float fracX = 0.85;
  *npsr = 0; 
  strcpy(covarFuncFile2,"NULL");

  printf("Graphical Interface: plotMany\n");
  printf("Author:              George Hobbs\n");
  printf("Version:             1.0\n");

  /* Obtain the .par and the .tim file from the command line */
  for (i=1;i<argc;i++)
    {
      if (strcmp(argv[i],"-f")==0)
	{	  
	  strcpy(parFile[*npsr],argv[i+1]); 
	  strcpy(timFile[*npsr],argv[i+2]);	  
	  (*npsr)++;
	  if (*npsr > MAX_PSR)
	    {
	      printf("Error, npsr > MAX_PSR.  Must increase MAX_PSR\n");
	      exit(1);
	    } 
	}
      else if (strcmp(argv[i],"-overlay")==0)
	sscanf(argv[i+1],"%d",&nOverlay);
      else if (strcmp(argv[i],"-dcf")==0){
        strcpy(covarFuncFile2,argv[++i]);
	  }
      else if (strcmp(argv[i],"-g")==0)
	strcpy(grDev,argv[++i]);
      else if (strcmp(argv[i],"-h")==0)
	help();
      else if (strcmp(argv[i],"-scale")==0)
	{
	  sscanf(argv[i+1],"%f",&scale[nScale]);
	  nScale++;
	}
      else if (strcmp(argv[i],"-fontsize")==0)
	sscanf(argv[i+1],"%f",&fontsize);
      else if (strcasecmp(argv[i],"-fracX")==0)
	sscanf(argv[i+1],"%f",&fracX);
      else if (strcasecmp(argv[i],"-labelSize")==0)
	sscanf(argv[i+1],"%f",&labelsize);
      else if (strcmp(argv[i],"-ptstyle")==0)
	  sscanf(argv[i+1],"%d",&ptStyle);
      else if (strcmp(argv[i],"-ptsize")==0)
	  sscanf(argv[i+1],"%f",&ptSize);
      else if (strcmp(argv[i],"-miny")==0)
	  sscanf(argv[i+1],"%f",&miny);
      else if (strcmp(argv[i],"-maxy")==0)
	  sscanf(argv[i+1],"%f",&maxy);
      else if (strcmp(argv[i],"-minx")==0)
	  sscanf(argv[i+1],"%f",&minx);
      else if (strcmp(argv[i],"-maxx")==0)
	  sscanf(argv[i+1],"%f",&maxx);
      else if (strcmp(argv[i],"-reverse")==0)
	reverse=1;
      else if (strcmp(argv[i],"-centremjd")==0)
	sscanf(argv[i+1],"%f",&centreMJD);
      else if (strcmp(argv[i],"-plotus")==0)
	plotUs=1;
      else if (strcmp(argv[i],"-noerror")==0)
	error=0;
    }

  if (*npsr==0) /* Select all files */
    {
      printf("Using all available .par and .tim files\n");
      sprintf(str,"ls `ls *.par | sed s/par/tim/` | sed s/.tim/\"\"/");
      pin = popen(str,"r");
      while (!feof(pin))
	{
	  if (fscanf(pin,"%s",str)==1)
	    {
	      sprintf(parFile[*npsr],"%s.par",str);
	      sprintf(timFile[*npsr],"%s.tim",str);
	      (*npsr)++;
	    }
	}
      pclose(pin);
      printf("Obtained files for %d pulsars\n",*npsr);
    }
 if (reverse==1) // Reverse order of .par and .tim files
   {
     for (i=0;i<*npsr;i++)
       {
     	 strcpy(tpar[i],parFile[i]);
     	 strcpy(ttim[i],timFile[i]);	 
       }
     for (i=0;i<*npsr;i++)
       {
	 //	 printf("Have %d %d %d >%s< >%s<\n",i,*npsr,(*npsr)-i-1,tpar[(*npsr)-i-1],ttim[(*npsr)-i-1]);
	 strcpy(parFile[i],tpar[(*npsr)-i-1]);
	 strcpy(timFile[i],ttim[(*npsr)-i-1]);
	 //	 strcpy(timFile[i],"J0900-3144.new2.tim");
       }
   }

 printf("Starting plugin1\n");
 //  initialise(psr,0);              /* Initialise the structures */
 readParfile(psr,parFile,timFile,*npsr); /* Load the parameters       */
 printf("Starting plugin2\n");
 readTimfile(psr,timFile,*npsr); /* Load the arrival times    */
 printf("Starting plugin3\n");
  preProcess(psr,*npsr,argc,argv);
 printf("Starting plugin4\n");
  callFit(psr,*npsr);             /* Do all the fitting routines */
 printf("Starting plugin5\n");
 doPlot(psr,*npsr,scale,nScale,grDev,plotUs,fontsize,centreMJD,ptStyle,ptSize,error,miny,maxy,minx,maxx,nOverlay,labelsize,fracX);              /* Do plot */
 printf("Starting plugin6\n");
}
Esempio n. 7
0
extern "C" int graphicalInterface(int argc,char *argv[],pulsar *psr,int *npsr) 
{
  longdouble imjd=-1.0,fmjd=-1.0,ra,grms=0.0;
  longdouble toa,ha0,almst,amjd,hnobs,mjd,trmjd,tstep=0.0;
  longdouble times[MAX_OBSN];
  longdouble out[MAX_OBSN];
  longdouble lowFreq,highFreq,alpha=-3.0,ampPL=1.0e-16;
  int setref=0;
  int nshots,npts;
  int j,count=0,i,k,ii,jj,kk;
  longdouble solsid  = 1.002737909;
  longdouble obslong = 149.0;  /* Longitude of Parkes */
  longdouble freq    = 1440.0;
  long iseed;
  int    site    = 7,p;
  int endit=0;
  long idum = 0;
  int nday = -1;
  float ngap=-1.0;
  char parFile[MAX_PSR][MAX_FILELEN];
  char timFile[MAX_PSR][MAX_FILELEN];
  char str[MAX_FILELEN],str2[MAX_FILELEN];
  double hamax =-1.0;
  char random[100];
  char read_f=0;
  FILE *fout,*fin;
  char addCubic[100]="n",temp[100];
  char formstr[50]="tempo2";
  char smooth[100]="n";
  int  psrNum,giveRMS=-1,setrand=-1,setred=-1,Npsr=0;
  int timesFile=0;
  char fake_fname[100];
  char have_outfile=0;
  char outfile[100];
  char timesfname[100];
  char telID[128]="7"; // Hardcode to Parkes
  int bunching=0; // flag on whether or not observations occur in groups
  // size of gap between observing runs, and length of observing runs.
  // These defaults give 7 observations every 28 days.
  long double gapsize=21,hillsize=7,gapstartmjd;
  // Flag whether or not to ask for red noise variables.
  *npsr = 1;

  strcpy(fake_fname,"fake.rf");

  for (i=0;i<argc;i++){
    if(strcmp(argv[i],"-ndobs")==0){
      sscanf(argv[i+1],"%f",&ngap);
      printf("Have >>%f<< days between observations\n",ngap);
    }
    if(strcmp(argv[i],"-nobsd")==0){
      sscanf(argv[i+1],"%d",&nday);
      printf("Have >>%d<< observations per day\n",nday);
    }
    if (strcmp(argv[i],"-idum")==0){
	sscanf(argv[i+1],"%d",&idum);
	printf("Have idum >>%d<<\n",idum);
    }
    if(strcmp(argv[i],"-ha")==0){
      sscanf(argv[i+1],"%lf",&hamax);
      printf("Have maximum absolute HA >>%lf<<\n",hamax);
    }
    if(strcmp(argv[i],"-randha")==0){
      strcpy(&random[0],argv[i+1]);
      printf("Have random >>%s<<\n",random);
      setrand=1;
    }
    if(strcmp(argv[i],"-start")==0){
      sscanf(argv[i+1],"%Lf",&imjd);
      printf("Have initial MJD >>%lf<<\n",(double)imjd);
    }
    if(strcmp(argv[i],"-tel")==0){
      sscanf(argv[i+1],"%s",telID);
    }
    if(strcmp(argv[i],"-end")==0){
      sscanf(argv[i+1],"%Lf",&fmjd);
      printf("Have final MJD >>%lf<<\n",(double)fmjd);
    }
    if(strcmp(argv[i],"-rms")==0){
      sscanf(argv[i+1],"%Lf",&grms);
      giveRMS=1;
      printf("Have Gaussian noise rms >>%lf<<\n",(double)grms);
    }
    if(strcmp(argv[i],"-format")==0){
	sscanf(argv[i+1],"%s",&formstr);
	printf("Have output format >>%s<<\n",formstr);
    }
    if (strcmp(argv[i],"-times")==0){
      timesFile=1;
      sscanf(argv[i+1],"%s",&timesfname);
      printf("Timesfile = %s\n",timesfname);
    }
    if (strcmp(argv[i],"-setref")==0)
      setref=1;

    if (strcmp(argv[i],"-o")==0){
      have_outfile=1;
      sscanf(argv[i+1],"%s",&outfile);
      printf("outfile = %s\n",outfile);
    }

    if (strcmp(argv[i],"-readtim")==0){
      read_f=1;
      printf("Read name,freq,mjd\n");
    }


    if(strcmp(argv[i],"-group")==0){
	bunching = 1;
	sscanf(argv[i+1],"%Lf",&hillsize);
	sscanf(argv[i+2],"%Lf",&gapsize);
	printf("Will simulate runs of %lg observing days long and leave a gap of %lg days between runs.\n",
	       (double)hillsize,(double)gapsize);
    }
    if(strcmp(argv[i],"-h")==0){
      printf("==========================================================================================\n");
      printf(" fake Tempo2 plugin - usage instructions.\n");
      printf(" tempo2 -gr fake -f file.par: The program will prompt you for parameters.\n");
      printf("\n Command-line arguments:\n");
      printf(" \t -f J0437-4715.par J1909-3744.par J1713+0747.par: specify a number of parfiles.\n");
      printf(" \t -ndobs xxx: specify number of days between observations.\n");
      printf(" \t -nobsd xxx: specify number of observations per day.\n");
      printf(" \t -ha xxx: specify maximal absolute Hour Angle.\n");
      printf(" \t -randha y: specify whether to use random HA coverage or not (y/n).\n");
      printf(" \t -start xxxxx: specify start MJD.\n");
      printf(" \t -end xxxxx: specify final MJD.\n");
      printf(" \t -rms xxx: specify Gaussian noise rms (in ms).\n");
      printf(" \t -times xxx: read observation times from specified file\n");
      printf(" \t           suppresses questions for red noise characteristics.\n");
      printf(" \t -format parkes : sets the tim-file format to tempo. (Default: tempo2).\n");
      printf(" \t -group 7 21 : simulate observations in groups of 7 days, with 21 days between groups.\n");
      printf(" \t               There will be nobsd observations every ndobs days during these 7 days.\n");
      printf("\n\n \t -idum xxx: specify random number seed (default = set from clock)\n");
      printf("\n The program will prompt you for the parameters not defined in the command line.\n");
      printf("\n\n Have a nice day!\n\n");
      printf("==========================================================================================\n");
      exit(0);
    }
    if(strcmp(argv[i],"-f")==0){
      Npsr=0;
      while(argv[i+Npsr+1][0]!='-'){
	strcpy(parFile[Npsr],argv[i+Npsr+1]);
	printf("Have parfile %d: >>%s<<.\n",Npsr+1,parFile[Npsr]);
	Npsr++;
	if(i+Npsr+1>=argc) break;
      }
      printf("Have %d parameter files.\n",Npsr);
    }
  }
  if (timesFile==1)
    {
      nday = 1;
      ngap = 1;
      hamax = 8;
      setrand = 1;
      imjd = 1;
      fmjd = 1;
    }
  printf("Simulate arrival times for parameter file %s\n",argv[3]);
  printf("----------------------------------------------------\n\n");
  if(ngap<0){    printf("Enter number of days between observations . "); scanf("%f",&ngap);}
  if(nday<0) {    printf("Enter number of observations/day .......... "); scanf("%d",&nday);}
  if(hamax<0 && nday!=1)  {  
      printf("Enter max absolute HA...................... "); scanf("%lf", &hamax);}
  if(setrand!=1){ printf("Random HA coverage? (y/n) ................. "); scanf("%s",&random);}
  if(imjd<0)   {  printf("Enter initial MJD ......................... "); scanf("%Lf",&imjd);}
  if(fmjd<0)   {  printf("Enter final MJD ........................... "); scanf("%Lf",&fmjd);}
  if(giveRMS<0){
    printf("Enter Gaussian noise rms  (ms/auto)........ "); scanf("%s",temp);
    giveRMS = sscanf(temp,"%Lf",&grms);
  }

  printf("GIVE RMS = %d\n",giveRMS);
  if (idum==0)
    {
      printf("Setting random number seed from the clock\n");      
      idum = TKsetSeed();
    }


  iseed = idum;
  psrNum = Npsr;
  if (timesFile==1)
    fin = fopen(timesfname,"r");

  hnobs = nday/2.0; 
  for(ii=0;ii<Npsr;ii++){
      count = 0;
      strcpy(parFile[0],parFile[ii]);
      printf("SET: %d\n",psr[0].param[param_pb].val[0]);
      psr[0].nJumps=0;
      psr[0].fitMode=0;
      psr[0].eclCoord=0;
  psr[0].nits=1;
  psr[0].clockFromOverride[0] = '\0';
  psr[0].nCompanion = 0;
  psr[0].bootStrap = 0;
  psr[0].units = SI_UNITS;
  psr[0].ne_sw  = NE_SW_DEFAULT; 
  psr[0].nWhite = 0;  /* No whitening by default */
  psr[0].timeEphemeris = IF99_TIMEEPH;
  psr[0].dilateFreq = 1;
  psr[0].planetShapiro = 1;
  psr[0].correctTroposphere = 1;
  psr[0].t2cMethod = T2C_IAU2000B;
  psr[0].fixedFormat=0;
  psr[0].nStorePrecision=0;
  strcpy(psr[0].deleteFileName,"NONE");
  strcpy(psr[0].tzrsite,"NULL");
  psr[0].calcShapiro=1;
  psr[0].ipm = 1;
  psr[0].swm = 0;
  psr[0].nPhaseJump=0;

      for (i=0;i<MAX_PARAMS;i++)
	{
	  for (j=0;j<psr[0].param[i].aSize;j++)
	    {
	      psr[0].param[i].fitFlag[j] = 0;
	      psr[0].param[i].paramSet[j] = 0;
	      psr[0].param[i].err[j] = 0;
	      psr[0].param[i].val[j] = 0;
	    }
	}
      //      initialise(psr,0);              /* Initialise the structures */      
      //      printf("SET AFTER 1: %d\n",psr[0].param[param_pb].val[0]);
      readParfile(psr,parFile,timFile,*npsr); /* Load the parameters       */
      //      printf("SET AFTER 2: %d %s\n",psr[0].param[param_pb].val[0],psr[0].name);
      ra = (double)psr[0].param[param_raj].val[0]/2.0/M_PI;
      
      /* Code based on fake1.f to calculate the TOA at transit for each of these observations */
      trmjd = imjd;
      /* 47892.0 = 1990??? */
      almst = fortran_mod((trmjd-47892.0)*solsid+0.276105324+obslong/360.0,(longdouble)1.0);
      
      /* Hour angle at 00h UT */
      ha0 = almst - ra;
      /* Approximate transit time */
      if (ha0 < 0.0) ha0+=1.0;
      trmjd += 1.0-ha0;
      
      if (nday > 1)tstep = hamax/12.0/nday;  /* Was 0.4/nday */
      amjd =  trmjd;
      gapstartmjd = imjd+(hillsize)/solsid;

      do {
	  for (j=0;j<nday;j++)
	      {
		if (timesFile==1)
		  {
		    if (read_f){
		      if (fscanf(fin,"%s %Lf %Lf\n",fake_fname,&freq,&mjd)==3)
			{
			  printf("Read %g %f\n",(double)mjd, (double)freq);
			  endit=0;
			}
		      else
			endit=1;
		    } else{
		      if (fscanf(fin,"%Lf",&mjd)==1)
			{
			  printf("Read %g\n",(double)mjd);
			  endit=0;
			}
		      else
			endit=1;
		    }
		  }
		else
		  {
		    if (random[0]=='y'||random[0]=='Y')
		      amjd=trmjd + (rand()/(longdouble)RAND_MAX - 0.5)*hamax/12.0;
		    else if (nday==1)
		      amjd=trmjd;
		    else
		      amjd=trmjd + ((j+1)-hnobs)*tstep;
		    
		    mjd=amjd;
		  }
		if (endit==0)
		  {
		    if (count==0 && setref==1)
		      {
			psr[0].obsn[count].sat    = psr[0].param[param_tzrmjd].val[0];
			strcpy(psr[0].obsn[count].fname,"reference");
			psr[0].obsn[count].freq   = psr[0].param[param_tzrfrq].val[0];
			if (giveRMS!=1) grms = psr[0].param[param_tres].val[0]/1e3;
			//	    else grms=0.0;
			psr[0].obsn[count].toaErr = grms*1000.0;
			psr[0].obsn[count].origErr = grms*1000.0;
			psr[0].obsn[count].phaseOffset = 0.0;
			strcpy(psr[0].obsn[count].telID, psr[0].tzrsite);
			psr[0].obsn[count].deleted = 0;
			psr[0].obsn[count].clockCorr=1;
			psr[0].obsn[count].delayCorr=1;
			psr[0].obsn[count].efac=1;
			count++;
		      }
		    psr[0].obsn[count].sat    = mjd;
		    strcpy(psr[0].obsn[count].fname,fake_fname);
		    psr[0].obsn[count].freq   = freq;
		    if (giveRMS!=1) grms = psr[0].param[param_tres].val[0]/1e3;
		    //	    else grms=0.0;
		    psr[0].obsn[count].toaErr = grms*1000.0;
		    psr[0].obsn[count].origErr = grms*1000.0;
		    psr[0].obsn[count].phaseOffset = 0.0;
		    strcpy(psr[0].obsn[count].telID, telID);
		    psr[0].obsn[count].deleted = 0;
		    psr[0].obsn[count].clockCorr=1;
		    psr[0].obsn[count].delayCorr=1;
		    psr[0].obsn[count].efac=1;
		    count++;
		    if (count>MAX_OBSN)
		      {
			printf("Number of TOAs > MAX_OBSN.\n");
			count--;
		      }
		  }
	      }
	  if((bunching == 1) && (trmjd >= (gapstartmjd-1))){
	    trmjd += gapsize/solsid;
	    gapstartmjd += (gapsize+hillsize)/solsid;
	  }
	  else{
	    trmjd+=ngap/solsid;
	  }
      }while ((timesFile == 0 && amjd<fmjd) || (timesFile == 1 && endit==0));
      if (timesFile==1)
    fclose(fin);


      psr[0].nobs=count;
      
      if (have_outfile){
	      strcpy(str,outfile);
	      strcpy(timFile[0],outfile);
      }else{
	      strcpy(str,parFile[0]);
	      str[strlen(str)-4]='\0';
	      strcat(str,".simulate");
	      strcpy(timFile[0],str);
      }
      
      /* Now run the tempo2 code */
      preProcess(psr,*npsr,argc,argv);
      callFit(psr,*npsr);             /* Do all the fitting routines */
      for (j=0;j<9;j++)
	{
	  /* Now update the site arrival times depending upon the residuals */
	  
	  for (i=0;i<psr[0].nobs;i++)  
	    {
	      psr[0].obsn[i].sat -= psr[0].obsn[i].prefitResidual/SECDAY; 
	      psr->obsn[i].nFlags = 0;
	    } 
	  writeTim(str,psr,"tempo2");
	  //	  initialise(&psr[ii],0);
	  // Reset the jumps
	  psr[ii].nJumps = 0;
	  for(kk=0;kk<MAX_JUMPS;kk++){
	      psr[ii].jumpVal[kk] = 0.0;
	      psr[ii].jumpValErr[kk] = 0.0;
	  }
	  for(jj=0;jj<MAX_PARAMS;jj++){
	      psr[ii].param[jj].nLinkTo = 0;
	      psr[ii].param[jj].nLinkFrom = 0;
	  }
	  readParfile(psr,parFile,timFile,*npsr); /* Load the parameters       */
	  readTimfile(psr,timFile,*npsr); 
	  preProcess(psr,1,argc,argv);
	  /* Now run the superTEMPO code again */
	  callFit(psr,*npsr);             /* Do all the fitting routines */
	}

      printf("Complete 10 iterations\n");
      for (i=0;i<psr[0].nobs;i++)  
	{
	  psr[0].obsn[i].sat -= psr[0].obsn[i].prefitResidual/SECDAY;  
	  psr->obsn[i].nFlags = 0;
	}
      
      
      for (i=0;i<psr[0].nobs;i++)
	{ 
	  times[i] = (psr[0].obsn[i].sat-psr[0].param[param_posepoch].val[0]);
	}
      npts = psr[0].nobs;
      
      /* Add Gaussian noise */
      if (giveRMS!=1) grms = psr[0].param[param_tres].val[0]/1e3;
      if (grms>0.0)
	{
	  printf("Adding Gaussian noise with rms = %f\n",(float)grms);
	  for (i=0;i<psr[0].nobs;i++)
	    psr[0].obsn[i].sat += TKgaussDev(&idum)*grms/1000.0/SECDAY;
	}
      
      printf("Output TOA file written to %s\n",str);
      writeTim(str,psr,formstr);      
    }
}
Esempio n. 8
0
extern "C" int graphicalInterface(int argc,char *argv[],pulsar *psr,int *npsr) 
{
  char parFile[MAX_PSR][MAX_FILELEN];
  char timFile[MAX_PSR][MAX_FILELEN];
  int i;
  int overlay=0;
  FILE *pin;
  char str[1000];
  *npsr = 0;  /* This graphical interface will only show results for one pulsar */


  strcpy(dcmFile,"NULL");
  strcpy(covarFuncFile,"NULL");

  printf("Graphical Interface: plk emulator\n");
  printf("Author:              George Hobbs (21 Dec 2003)\n");
  printf("CVS Version:         $Revision $\n");

  /* Obtain the .par and the .tim file from the command line */
  for (i=1;i<argc;i++)
    {
      if (strcmp(argv[i],"-f")==0)
	{	  
	  strcpy(parFile[*npsr],argv[i+1]); 
	  strcpy(timFile[*npsr],argv[i+2]);	  
	  (*npsr)++;
	  if (*npsr > MAX_PSR)
	    {
	      printf("Error, npsr > MAX_PSR.  Must increase MAX_PSR\n");
	      exit(1);
	    } 
	}
      else if (strcmp(argv[i],"-dcm")==0)
	strcpy(dcmFile,argv[++i]);
      else if (strcmp(argv[i],"-dcf")==0)
	strcpy(covarFuncFile,argv[++i]);
      else if (strcmp(argv[i],"-overlay")==0)
	overlay=1;
    }

if (*npsr==0) /* Select all files */
    {
      printf("Using all available .par and .tim files\n");
      sprintf(str,"ls `ls *.par | sed s/par/tim/` | sed s/.tim/\"\"/");
      pin = popen(str,"r");
      while (!feof(pin))
	{
	  if (fscanf(pin,"%s",str)==1)
	    {
	      sprintf(parFile[*npsr],"%s.par",str);
	      sprintf(timFile[*npsr],"%s.tim",str);
	      (*npsr)++;
	    }
	}
      pclose(pin);
      printf("Obtained files for %d pulsars\n",*npsr);
    }

//  initialise(psr,0);              /* Initialise the structures */
  readParfile(psr,parFile,timFile,*npsr); /* Load the parameters       */
  readTimfile(psr,timFile,*npsr); /* Load the arrival times    */
  preProcess(psr,*npsr,argc,argv);
  callFit(psr,*npsr);             /* Do all the fitting routines */
  doPlot(psr,*npsr,overlay);              /* Do plot */
  return 0;
}