Esempio n. 1
0
void Pulsar::SquareWave::get_transitions (const Profile* profile,
					  vector<unsigned>& up,
					  vector<unsigned>& down)
{
  unsigned nbin = profile->get_nbin();
  unsigned scrunch = 1;

  Reference::To<Profile> clone;
  if (use_nbin && nbin > use_nbin) {
    clone = profile->clone();
    scrunch = nbin/use_nbin;
    clone->bscrunch(scrunch);
    profile = clone;
    nbin = profile->get_nbin();
  }

  unsigned offset = (unsigned) (risetime * nbin);

cerr << "nbin=" << nbin << " offset=" << offset << endl;

  // differentiate the profile
  Reference::To<Profile> difference = differentiate (profile, offset);

  float* amps = difference->get_amps();

  // find the phase window in which the mean is closest to zero
  BaselineWindow window;
  window.set_find_mean (0.0);
  window.get_smooth()->set_turns (0.2);
  float zero = window.find_phase (nbin, amps);

  // get the noise statistics of the zero mean region
  double mean = 0;
  double variance = 0;
  difference->stats (zero, &mean, &variance);

cerr << "mean=" << mean << " rms=" << sqrt(variance) << endl;

  // check that the mean is actually zero
  double rms = sqrt(variance);
  if (mean > rms)
    throw Error (InvalidState, "Pulsar::SquareWave::get_transitions",
		 "mean=%lf > rms=%lf", mean, rms);

  float cutoff = threshold * rms;

  find_transitions (nbin, amps, up, cutoff);
  find_transitions (nbin, amps, down, -cutoff);
}
Esempio n. 2
0
void psrspa::matched_finder ( const Archive* arch )
{
  string name = arch->get_filename ();
  for ( unsigned isub = 0; isub < arch->get_nsubint(); isub++ )
  {
    Reference::To<Profile> profile = arch->get_Profile(isub,0,0)->clone();

    while (profile->get_nbin() > 256)
    {
      finder->set_Profile( profile );  // might finder optimize on &profile?

      PhaseWeight weight;
      finder->get_weight( &weight );

      matched_report (name, isub, weight, *profile);

      profile->bscrunch (2);
    }
  }
}
Esempio n. 3
0
int main (int argc, char *argv[]) try {
  
    bool verbose = false;
    char* metafile = 0;
 
    string ulpath;

    bool save = false;
    string ext;
  
    bool tscr = false;
    int tscr_fac = 0;

    bool fscr = false;
    int fscr_fac = 0;

    bool bscr = false;
    int bscr_fac = 0;

    bool newdm = false;
    double dm = 0.0;

    bool scattered_power_correction = false;

    bool defaraday = false;

    bool newrm = false;
    double rm = 0.0;

    bool reset_weights = false;
    float new_weight = 1.0;

    float smear_dc = 0.0;

    bool rotate = false;
    double rphase = 0.0;

    bool dedisperse = false;
    bool dededisperse = false;

    bool pscr = false;

    bool invint = false;

    bool stokesify = false;
    bool unstokesify = false;

    bool flipsb = false;
    bool flip_freq = false;
    double flip_freq_mhz = 0.0;

    Pulsar::Parameters* new_eph = 0;

    string command = "pam";

    char* archive_class = 0;

    int new_nchn = 0;
    int new_nsub = 0;
    int new_nbin = 0;

    float tsub = 0.0;

    bool circ = false;
    bool lin = false;

    unsigned ronsub = 0;
    bool cbppo = false;
    bool cbpao = false;
    bool cblpo = false;
    bool cblao = false;

    int subint_extract_start = -1;
    int subint_extract_end = -1;

    bool new_cfreq = false;
    double new_fr = 0.0;
    Signal::Source new_type = Signal::Unknown;
    string instrument;
    bool reverse_freqs = false;
    string site;
    string name;
    float mult = -1.0;
    double new_folding_period = -1.0;

    bool update_dm_from_eph = false;
    double aux_rm = 0.0;

    Reference::To<Pulsar::IntegrationOrder> myio;
    Reference::To<Pulsar::Receiver> install_receiver;

    Pulsar::ReflectStokes reflections;

    int c = 0;

    const int TYPE = 1208;
    const int INST = 1209;
    const int REVERSE_FREQS = 1210;
    const int SITE = 1211;
    const int NAME = 1212;
    const int DD   = 1213;
    const int RR   = 1214;
    const int SPC  = 1215;
    const int RM   = 1216;
    const int MULT = 1218;
    const int PERIOD=1219;
    const int SS   = 1220;
    const int FLIP = 1221;
    const int UPDATE_DM = 1222;
    const int AUX_RM = 1223;

    while (1) {

      int options_index = 0;

      static struct option long_options[] = {
	{"setnchn",    1, 0, 200},
	{"setnsub",    1, 0, 201},
	{"setnbin",    1, 0, 202},
	{"binphsperi", 1, 0, 203},
	{"binphsasc",  1, 0, 204},
	{"binlngperi", 1, 0, 205},
	{"binlngasc",  1, 0, 206},
	{"receiver",   1, 0, 207},
	{"settsub",    1, 0, 208},
	{"type",       1, 0, TYPE},
	{"inst",       1, 0, INST},
	{"reverse_freqs",no_argument,0,REVERSE_FREQS},
	{"flip",       1 ,0, FLIP},
	{"site",       1, 0, SITE},
	{"name",       1, 0, NAME},
	{"DD",         no_argument,      0,DD},
	{"RR",         no_argument,      0,RR},
	{"RM",         required_argument,0,RM},
	{"spc",        no_argument,      0,SPC},
	{"mult",       required_argument,0,MULT},
	{"period",     required_argument,0,PERIOD},
	{"SS",         no_argument,      0,SS},
	{"update_dm",   no_argument,      0,UPDATE_DM},
	{"aux_rm",    required_argument,0,AUX_RM},
	{0, 0, 0, 0}
      };

      c = getopt_long(argc, argv, "hqvViM:mn:a:e:E:TFpIt:f:b:d:o:s:r:u:w:DSBLCx:R:",
		      long_options, &options_index);

      if (c == -1)
	break;

      switch (c) {
      case 'h':
	usage();
	return (0);
	break;
      case 'q':
	Pulsar::Archive::set_verbosity(0);
	break;
      case 'v':
	verbose = true;
	Pulsar::Archive::set_verbosity(2);
	break;
      case 'V':
	verbose = true;
	Pulsar::Archive::set_verbosity(3);
	break;
      case 'i':
	cout << "$Id: pam.C,v 1.101 2010/10/05 23:59:50 jonathan_khoo Exp $" << endl;
	return 0;
      case 'm':
	save = true;
	break;
      case 'M':
        metafile = optarg;
        break;
      case 'L':
	lin = true;
	break;
      case 'C':
	circ = true;
	break;
      case 'a':
	archive_class = optarg;
	break;
      case 'e':
	ext = optarg;
	if( !ext.empty() )
	  save = true;
	break;
      case 'E':

	try {
	  new_eph = factory<Pulsar::Parameters> (optarg);
	}
	catch (Error& error) {
	  cerr << "Could not load new ephemeris from " << optarg << endl;
	  return -1;
	}

	command += " -E";
	break;
      case 'T':
	tscr = true;
	command += " -T";
	break;
      case 'F':
	fscr = true;
	command += " -F";
	break;
      case 'p':
	pscr = true;
	command += " -p";
	break;
      case 'I':
	invint = true;
	pscr = false;
	command += " -I";
	break;
      case 'f':
	fscr = true;
	if (sscanf(optarg, "%d", &fscr_fac) != 1) {
	  cout << "That is not a valid fscrunch factor" << endl;
	  return -1;
	}
	command += " -f ";
	command += optarg;
	break;
	
      case 'n':

	reflections.add_reflection( optarg[0] );

	command += " -n ";
	command += optarg;
	break;

      case 'o':
	new_cfreq = true;
	if (sscanf(optarg, "%lf", &new_fr) != 1) {
	  cout << "That is not a valid centre frequency" << endl;
	  return -1;
	}
	command += " -o ";
	command += optarg;
	break;
      case 't':
	tscr = true;
	if (sscanf(optarg, "%d", &tscr_fac) != 1) {
	  cout << "That is not a valid tscrunch factor" << endl;
	  return -1;
	}
	command += " -t ";
	command += optarg;
	break;
      case 'b':
	bscr = true;
	if (sscanf(optarg, "%d", &bscr_fac) != 1) {
	  cout << "That is not a valid bscrunch factor" << endl;
	  return -1;
	}
	if (bscr_fac <= 0) {
	  cout << "That is not a valid bscrunch factor" << endl;
	  return -1;
	}
	command += " -b ";
	command += optarg;
	break;
      case 'd':
	newdm = true;
	if (sscanf(optarg, "%lf", &dm) != 1) {
	  cout << "That is not a valid dispersion measure" << endl;
	  return -1;
	}
	command += " -d ";
	command += optarg;
	break;
      case 'D':
	dedisperse = true;
	command += " -D ";
	break;
      case 'R':
	if (sscanf(optarg, "%lf", &rm) != 1) {
	  cout << "That is not a valid rotation measure" << endl;
	  return -1;
	}
	newrm = true;
	defaraday = true;
	command += " -R ";
	command += optarg;
	break;
      case 's':
	if (sscanf(optarg, "%f", &smear_dc) != 1) {
	  cout << "That is not a valid smearing duty cycle" << endl;
	  return -1;
	}
	command += " -s ";
	command += optarg;
	break;
      case 'r':
	rotate = true;
	if (sscanf(optarg, "%lf", &rphase) != 1) {
	  cout << "That is not a valid rotation phase" << endl;
	  return -1;
	}
	if (rphase <= -1.0 || rphase >= 1.0) {
	  cout << "That is not a valid rotation phase" << endl;
	  return -1;
	}
	command += " -r ";
	command += optarg;
	break;
      case 'u':
	ulpath = optarg;
	if( !ulpath.empty() )
	{
	  save = true;
	  if (ulpath.substr(ulpath.length()-1,1) != "/")
	    ulpath += "/";
	}
	break;
      case 'w':
	reset_weights = true;
	if (sscanf(optarg, "%f", &new_weight) != 1) {
	  cout << "That is not a valid weight" << endl;
	  return -1;
	}
	command += " -w ";
	command += optarg;
	break;
      case 'S':
	stokesify = true;
	break;
      case SS:
        unstokesify = true;
        break;
      case 'B':
	flipsb = true;
	break;
      case 'x' :
	if (sscanf(optarg, "%d %d", 
		   &subint_extract_start, &subint_extract_end) !=2 ) {
	  cout << "That is not a valid subint range" << endl;
	  return -1;
	}
	subint_extract_end++;
	break;
      case 200:
	fscr = true;
	if (sscanf(optarg, "%d", &new_nchn) != 1) {
	  cout << "That is not a valid number of channels" << endl;
	  return -1;
	}
	if (new_nchn <= 0) {
	  cout << "That is not a valid number of channels" << endl;
	  return -1;
	}
	command += " --setnchn ";
	command += optarg;
	break;
      case 201:
	tscr = true;
	if (sscanf(optarg, "%d", &new_nsub) != 1) {
	  cout << "That is not a valid number of subints" << endl;
	  return -1;
	}
	if (new_nsub <= 0) {
	  cout << "That is not a valid number of subints" << endl;
	  return -1;
	}
	command += " --setnsub ";
	command += optarg;
	break;
      case 202:
	bscr = true;
	if (sscanf(optarg, "%d", &new_nbin) != 1) {
	  cout << "That is not a valid number of bins" << endl;
	  return -1;
	}
	if (new_nbin <= 0) {
	  cout << "That is not a valid number of bins" << endl;
	  return -1;
	}
	command += " --setnbin ";
	command += optarg;
	break;
      case 203: {
	if (cbpao || cblpo || cblao) {
	  cerr << "You can only specify one re-ordering scheme!"
	       << endl;
	  return -1;
	}
	if (sscanf(optarg, "%ud", &ronsub) != 1) {
	  cerr << "Invalid nsub given" << endl;
	  return -1;
	}
	cbppo = true;
	break;
      } 
      case 204: {
	if (cbppo || cblpo || cblao) {
	  cerr << "You can only specify one re-ordering scheme!"
	       << endl;
	  return -1;
	}
	if (sscanf(optarg, "%ud", &ronsub) != 1) {
	  cerr << "Invalid nsub given" << endl;
	  return -1;
	}
	cbpao = true;
	break;
      }      
      case 205: {
	if (cblao || cbppo || cbpao) {
	  cerr << "You can only specify one re-ordering scheme!"
	       << endl;
	  return -1;
	}
	if (sscanf(optarg, "%ud", &ronsub) != 1) {
	  cerr << "Invalid nsub given" << endl;
	  return -1;
	}
	cblpo = true;
	break;
      } 
      case 206: {
	if (cblpo || cbppo || cbpao) {
	  cerr << "You can only specify one re-ordering scheme!"
	       << endl;
	  return -1;
	}
	if (sscanf(optarg, "%ud", &ronsub) != 1) {
	  cerr << "Invalid nsub given" << endl;
	  return -1;
	}
	cblao = true;
	break;
      }

      case 207: try {
	install_receiver = Pulsar::Receiver::load (optarg);
	break;
      }
      catch (Error& error) {
	cerr << "pam: Error loading Receiver from " << optarg << endl
	     << error.get_message() << endl;
	return -1;
      }

      case 208: {
	if (sscanf(optarg, "%f", &tsub) != 1) {
	  cerr << "Invalid tsub given" << endl;
	  return -1;
	}
	tscr = true;
	break;
      }

      case TYPE:
	{
	  string s = optarg;
	  if(s=="Pulsar")     new_type = Signal::Pulsar;
	  else if(s=="PolnCal")    new_type = Signal::PolnCal;
	  else if(s=="FluxCalOn")  new_type = Signal::FluxCalOn;
	  else if(s=="FluxCalOff") new_type = Signal::FluxCalOff;
	  else if(s=="Calibrator") new_type = Signal::Calibrator;
	  else{
	    fprintf(stderr,"Unrecognised argument to --type: '%s'\n",optarg);
	    exit(-1);
	  }
	  command += " --type " + s;
	}
	break;

      case INST: instrument = optarg; break;

      case REVERSE_FREQS: reverse_freqs = true; break;

      case SITE: site = optarg; break;

      case NAME: name = optarg; break;

      case DD: dededisperse = true; break;

      case RM:
        aux_rm = fromstring<double>(optarg);
        newrm = true;
        command += " --RM ";
        command += optarg;
        break;

      case SPC: scattered_power_correction = true; break;
	
      case MULT: mult = atof(optarg); break;

      case PERIOD: new_folding_period = fromstring<double>(optarg); break;

      case FLIP: flip_freq = true; flip_freq_mhz = atof(optarg); break;

      case UPDATE_DM: update_dm_from_eph = true; break;

      case AUX_RM:
        aux_rm = fromstring<double>(optarg);
        command += " --aux_rm ";
        command += optarg;
        break;

      default:
	cout << "Unrecognised option" << endl;
      }
    }

   if (verbose)
     cerr << "pam: parsing filenames" << endl;
 
   vector <string> filenames;

    if (metafile)
      stringfload (&filenames, metafile);
    else
      for (int ai=optind; ai<argc; ai++)
        dirglob (&filenames, argv[ai]);
 
    if (filenames.empty())
    {
      cerr << "pam: no filenames were specified" << endl;
      exit(-1);
    } 
  
    Reference::To<Pulsar::Archive> arch;

    if (!save)
    {
      cout << "Changes will not be saved. Use -m, -u or -e to write results to disk"
	   << endl;
    }

    if (stokesify && unstokesify)
    {
      cerr << "pam: Both -S and --SS options were given.  Poln state will not be changed!" << endl;
      stokesify = false;
      unstokesify = false;
    }

    int flip_option_count=0;
    if (flipsb) flip_option_count++;
    if (flip_freq) flip_option_count++;
    if (reverse_freqs) flip_option_count++;
    if (flip_option_count > 1) {
      cerr << "pam: More than one band-flip option was given, exiting." << endl;
      exit(-1);
    }

    for (unsigned i = 0; i < filenames.size(); i++) try
    {
      if (verbose)
	cerr << "Loading " << filenames[i] << endl;
      
      arch = Pulsar::Archive::load(filenames[i]);

      if( mult > 0.0 ){
	for( unsigned isub=0; isub<arch->get_nsubint();isub++)
	  for( unsigned ichan=0; ichan<arch->get_nchan();ichan++)
	    for( unsigned ipol=0; ipol<arch->get_npol();ipol++)
	      arch->get_Profile(isub,ipol,ichan)->scale( mult );
      }

      if( new_folding_period > 0.0 ){
	Pulsar::counter_drift( arch, new_folding_period, 0.0);
	for( unsigned isub=0; isub<arch->get_nsubint();isub++)
	  arch->get_Integration(isub)->set_folding_period( new_folding_period );
      }

      if (install_receiver) {
	if (verbose)
	  cerr << "pam: Installing receiver: " << install_receiver->get_name()
	       << " in archive" << endl;

	arch->add_extension (install_receiver);
      }

      if (lin || circ) {
	Pulsar::Receiver* receiver = arch->get<Pulsar::Receiver>();

	if (!receiver)
	  cerr << "No Receiver Extension in " << filenames[i] << endl;
	else {
	  if (lin) {
	    receiver->set_basis (Signal::Linear);
	    cout << "Feed basis set to Linear" << endl;
	  }

	  if (circ) {
	    receiver->set_basis (Signal::Circular);
	    cout << "Feed basis set to Circular" << endl;
	  }
	}
      }

      reflections.transform( arch );

      if (new_cfreq)
      {
	double nc = arch->get_nchan();
	double bw = arch->get_bandwidth();
	double cw = bw / nc;

	double fr = new_fr - (bw / 2.0) + (cw / 2.0);
	
	for (unsigned i = 0; i < arch->get_nsubint(); i++) {
	  for (unsigned j = 0; j < arch->get_nchan(); j++) {
	    arch->get_Integration(i)->set_centre_frequency(j,(fr + (j*cw)));
	  }
	}
	
	arch->set_centre_frequency(new_fr);
      }

      if( new_type != Signal::Unknown )
	arch->set_type( new_type );

      if( instrument != string() ){
	Pulsar::Backend* b = arch->get<Pulsar::Backend>();
	if( !b )
	  fprintf(stderr,"Could not change instrument name- archive does not have Backend extension\n");
	else
	  b->set_name(instrument);
      }

      if( site != string() )
	arch->set_telescope( site );

      if( name != string() )
	arch->set_source( name );

      if (new_eph) try
      {
        arch->set_ephemeris(new_eph);

        if (update_dm_from_eph) {
          update_dm(arch);
        }
      }
      catch (Error& error)
      {
	cerr << "Error while installing new ephemeris: " 
	     << error.get_message() << endl;
        continue;
      }

      if (flipsb) {
	for (unsigned i = 0; i < arch->get_nsubint(); i++) {
	  vector<double> labels;
	  labels.resize(arch->get_nchan());
	  for (unsigned j = 0; j < arch->get_nchan(); j++) {
	    labels[j] = arch->get_Integration(i)->get_centre_frequency(j);
	  }
	  for (unsigned j = 0; j < arch->get_nchan(); j++) {
	    double new_frequency = labels[labels.size()-1-j];
	    arch->get_Integration(i)->set_centre_frequency(j,new_frequency);
	  }
	}
	arch->set_bandwidth(-1.0 * arch->get_bandwidth());
      }

      if (flip_freq) {
        for (unsigned isub = 0; isub < arch->get_nsubint(); isub++) {
          Reference::To<Pulsar::Integration> 
            subint = arch->get_Integration(isub);
          for (unsigned ichan = 0; ichan < arch->get_nchan(); ichan++) {
            double new_freq = flip_freq_mhz 
              - (subint->get_centre_frequency(ichan) - flip_freq_mhz);
            subint->set_centre_frequency(ichan, new_freq);
          }
        }
        arch->set_bandwidth(-1.0 * arch->get_bandwidth());
      }

      if( reverse_freqs ) {
	// Of course it would be nice to do this with pointers.... but oh well I guess copying will have to do HSK 27/8/04

	unsigned nchan = arch->get_nchan();

	for( unsigned isub=0; isub<arch->get_nsubint(); isub++){
	  for( unsigned ipol =0; ipol<arch->get_npol(); ipol++){
	    for( unsigned ichan=0; ichan<nchan/2; ichan++){
	      Reference::To<Pulsar::Profile> lo = arch->get_Profile(isub,ipol,ichan);	      
	      Reference::To<Pulsar::Profile> tmp = lo->clone();

	      Reference::To<Pulsar::Profile> hi = arch->get_Profile(isub,ipol,nchan-1-ichan);

	      lo->operator=(*hi);
	      hi->operator=(*tmp);
	    }
	  }
	}
	arch->set_bandwidth( -1.0 * arch->get_bandwidth() );
      }

      if (reset_weights) {
	arch->uniform_weight(new_weight);
	if (verbose)
	  cout << "All profile weights set to " << new_weight << endl;
      }
      
      if (rotate)
	arch->rotate_phase (rphase);

      if (scattered_power_correction) {

	Pulsar::ScatteredPowerCorrection spc;
	if (arch->get_state() == Signal::Stokes)
	  arch->convert_state(Signal::Coherence);

	spc.correct (arch);

      }

      if (newdm)
      {
	arch->set_dispersion_measure(dm);
	if (verbose)
	  cout << "Archive dispersion measure set to " << dm << endl;

	if (arch->get_dedispersed())
        {
	  arch->dedisperse();

	  if (verbose)
	    cout << "Archive re-dedipsersed" << endl;
        }
      }

      if (dedisperse)
      {
	arch->dedisperse();
	if (verbose)
	  cout << "Archive dedipsersed" << endl;
      }

      if (dededisperse)
      {
	Pulsar::Dispersion correction;
	correction.revert (arch);
      }

      if (stokesify) {
	if (arch->get_npol() != 4)
	  throw Error(InvalidState, "Convert to Stokes",
		      "Not enough polarisation information");
	arch->convert_state(Signal::Stokes);
	if (verbose)
	  cout << "Archive converted to Stokes parameters" << endl;
      }

      if (unstokesify) {
	if (arch->get_npol() != 4)
	  throw Error(InvalidState, "Convert to coherence",
		      "Not enough polarisation information");
	arch->convert_state(Signal::Coherence);
	if (verbose)
	  cout << "Archive converted to coherence parameters" << endl;
      }

      if (cbppo) {
	myio = new Pulsar::PeriastronOrder();
	arch->add_extension(myio);
	myio->organise(arch, ronsub);
      }
      
      if (cbpao) {
	myio = new Pulsar::BinaryPhaseOrder();
	arch->add_extension(myio);
	myio->organise(arch, ronsub);
      }
      
      if (cblpo) {
	myio = new Pulsar::BinLngPeriOrder();
	arch->add_extension(myio);
	myio->organise(arch, ronsub);
      }
      
      if (cblao) {
	myio = new Pulsar::BinLngAscOrder();
	arch->add_extension(myio);
	myio->organise(arch, ronsub);
      }
      
      if( subint_extract_start >= 0 && subint_extract_end >= 0 ) {
	vector<unsigned> subints;
	unsigned isub = subint_extract_start;

	while ( isub<arch->get_nsubint() && isub<unsigned(subint_extract_end) ) {
	  subints.push_back( isub );
	  isub++;
	}

	Reference::To<Pulsar::Archive> extracted( arch->extract(subints) );
	extracted->set_filename( arch->get_filename() );

	arch = extracted;
      }

      if (tscr) {
	if (tsub > 0.0) {
	  unsigned factor = 
	    unsigned (tsub / arch->get_Integration(0)->get_duration());
	  if (factor == 0) {
	    cerr << "Warning: subints already too long" << endl;
	  }
	  else {
	    arch->tscrunch(factor);
	  }
	  if (verbose)
	    cout << arch->get_filename() << " tscrunched by a factor of " 
		 << factor << endl;
	}
	else if (new_nsub > 0) {
	  arch->tscrunch_to_nsub(new_nsub);
	  if (verbose)
	    cout << arch->get_filename() << " tscrunched to " 
		 << new_nsub << " subints" << endl;
	}
	else if (tscr_fac > 0) {
	  arch->tscrunch(tscr_fac);
	  if (verbose)
	    cout << arch->get_filename() << " tscrunched by a factor of " 
		 << tscr_fac << endl;
	}
	else {
	  arch->tscrunch();
	  if (verbose)
	    cout << arch->get_filename() << " tscrunched" << endl;
	}
      }
      
      if (pscr) {
	arch->pscrunch();
	if (verbose)
	  cout << arch->get_filename() << " pscrunched" << endl;
      } 

      if (invint) {
	arch->invint();
	if (verbose)
	  cout << arch->get_filename() << " invinted" << endl;
      }

      if (newrm) {
	arch->set_rotation_measure (rm);
	if (verbose)
	  cout << arch->get_filename() << " RM set to " << rm << endl;
      }

      if (defaraday) {
	arch->defaraday();
	if (verbose)
	  cout << arch->get_filename() << " defaradayed" <<endl;
      }

      if (aux_rm)
      {
	if (verbose)
	  cout << "pam: correct auxiliary Faraday rotation; iono RM="
	       << aux_rm << endl;
	correct_auxiliary_rm (arch, aux_rm);
      }

      if (fscr) {
	if (new_nchn > 0) {
	  arch->fscrunch_to_nchan(new_nchn);
	  if (verbose)
	    cout << arch->get_filename() << " fscrunched to " 
		 << new_nchn << " channels" << endl;
	}
	else if (fscr_fac > 0) {
	  arch->fscrunch(fscr_fac);
	  if (verbose)
	    cout << arch->get_filename() << " fscrunched by a factor of " 
		 << fscr_fac << endl;
	}
	else {
	  arch->fscrunch();
	  if (verbose)
	    cout << arch->get_filename() << " fscrunched" << endl;
	}
      }

      if (bscr) {
	if (new_nbin > 0) {
	  arch->bscrunch_to_nbin(new_nbin);
	  if (verbose)
	    cout << arch->get_filename() << " bscrunched to " 
		 << new_nbin << " bins" << endl;
	}
	else {
	  arch->bscrunch(bscr_fac);
	  if (verbose)
	    cout << arch->get_filename() << " bscrunched by a factor of " 
		 << bscr_fac << endl;
	}
      }     

      if (smear_dc) {
	for (unsigned i = 0; i < arch->get_nsubint(); i++) {
	  for (unsigned j = 0; j < arch->get_npol(); j++) {
	    for (unsigned k = 0; k < arch->get_nchan(); k++) {
	      smear (arch->get_Profile(i,j,k), smear_dc);
	    }
	  }
	}
      }
      
      if (save) {

	if (archive_class)  {

	  // unload an archive of the specified class
	  Reference::To<Pulsar::Archive> output;
	  output = Pulsar::Archive::new_Archive (archive_class);
	  output -> copy (*arch);
	  output -> set_filename ( arch->get_filename() );

	  arch = output;

	}


	// See if the archive contains a history that should be updated:
	
	Pulsar::ProcHistory* fitsext = arch->get<Pulsar::ProcHistory>();
	
	if (fitsext) {
	  
	  if (command.length() > 80) {
	    cout << "WARNING: ProcHistory command string truncated to 80 chars" 
		 << endl;
	    fitsext->set_command_str(command.substr(0, 80));
	  }
	  else {
	    fitsext->set_command_str(command);
	  }
	  
	}
	
	string out_filename = arch->get_filename();
	
	if( !ext.empty() )
	  out_filename = replace_extension( out_filename, ext );
	
	if( !ulpath.empty() )
	  out_filename = ulpath + basename(out_filename);
	
	arch->unload( out_filename );
	cout << out_filename << " written to disk" << endl;
      }
    }  
    catch (Error& error) {
      cerr << error << endl;
    } 
  
    return 0;

}
catch(Error& er) 
{
  cerr << er << endl;
  return -1;
}
catch (string& error)
{
  cerr << "exception thrown: " << error << endl;
  return -1;
}
catch (bad_alloc& ba)
{
  cerr << "Caught a bad_alloc: '" << ba.what() << "'" << endl ;
  return -1;
}
catch (exception& e)
{
  cerr << "caught an exception of type '" 
				 << typeid(e).name() << "'" << endl; 
  return -1;
}
catch(...)
{
  fprintf(stderr,"Unknown exception caught\n");
  return -1;
}
Esempio n. 4
0
int main (int argc, char** argv) try
{
  // the multiple component model
  ComponentModel model;

  string model_filename_in;
  string model_filename_out = "paas.m";
  string details_filename   = "paas.txt";
  string std_filename       = "paas.std";

  bool fit = false;
  vector<string> new_components;
  string fit_flags;

  int bscrunch = 1;

  string pgdev;

  bool line_plot=false;
  bool interactive = false;
  bool centre_model = false, rotate_peak=false;
  float rotate_amount = 0.0;
  bool align = false;

  const char* args = "hb:r:w:c:fF:it:d:Dl:j:Ws:CpR:a";
  int c;

  while ((c = getopt(argc, argv, args)) != -1)
    switch (c) {

    case 'h':
      usage ();
      return 0;

    case 'r':
      model_filename_in = optarg;
      break;
      
    case 'w':
      model_filename_out = optarg;
      break;
      
    case 'f':
      fit = true;
      break;

    case 'W':
      model.set_fit_derivative (true);
      break;

    case 'c':
      new_components.push_back(optarg);
      break;

    case 'F':
      fit_flags = optarg;
      break;

    case 'b':
      bscrunch = atoi (optarg);
      break;

    case 't':
      model.set_threshold( atof(optarg) );
      break;

    case 'd':
      pgdev = optarg;
      break;

    case 'i':
      interactive = true;
    case 'D':
      pgdev = "/xs";
      break;
      
    case 'l':
      line_plot = true;
      break;

    case 's':
      std_filename = optarg;
      break;

    case 'C':
      centre_model = true;
      break;

    case 'p':
      rotate_peak = true;
      break;

    case 'R':
      rotate_amount = atof(optarg);
      break;

    case 'a':
      align = true;
      break;

    case 'j':
      details_filename = optarg;
      break;

   default:
      cerr << "invalid param '" << c << "'" << endl;
    }

  if (!pgdev.empty())
  {
    cpgopen(pgdev.c_str());
    cpgask(0);
    cpgsvp(0.1, 0.9, 0.1, 0.9);
  }
  
  Reference::To<Archive> archive = Archive::load (argv[optind]);

  // preprocess
  archive->fscrunch();
  archive->tscrunch();
  archive->pscrunch();

  // phase up as requested
  if (centre_model)
    archive->centre();
  else if (rotate_peak)
    archive->centre_max_bin(0.0);
  
  if (rotate_amount != 0.0)
    archive->rotate_phase(-rotate_amount);
  
  archive->remove_baseline();

  // load from file if specified
  if (!model_filename_in.empty())
  {
    cerr << "paas: loading model from " << model_filename_in << endl;
    model.load(model_filename_in.c_str());
    
    // align to profile first if asked
    if (align)
      model.align(archive->get_Integration(0)->get_Profile(0,0));
  }

  // add any new components specified
  unsigned ic, nc=new_components.size();
  double centre, concentration, height;
  int name_start;
  for (ic=0; ic < nc; ic++)
  {
    if (sscanf(new_components[ic].c_str(), 
	       "%lf %lf %lf %n", &centre, &concentration, &height,
	       &name_start)!=3)
    {
      cerr << "Could not parse component " << ic << endl;
      return -1;
    }

    model.add_component(centre, concentration, height, 
			new_components[ic].c_str()+name_start);
  }

  bool iterate = true;

  while (iterate)
  {
    iterate = false;

    // fit if specified
    if (fit)
    {
      // set fit flags
      model.set_infit(fit_flags.c_str());
      // fit
      model.fit(archive->get_Integration(0)->get_Profile(0,0));
    }
 
    // plot
    if (!pgdev.empty())
    {
      Reference::To<Pulsar::Archive> scrunched;
      scrunched = archive->clone();
      if (bscrunch > 1)
	scrunched->bscrunch(bscrunch);

      cpgpage();

      Profile *prof = scrunched->get_Integration(0)->get_Profile(0,0);
      float ymin = prof->min();
      float ymax = prof->max();
      float extra = 0.05*(ymax-ymin);

      ymin -= extra;
      ymax += extra;
      cpgswin(0.0, 1.0, ymin, ymax);
      cpgbox("bcnst", 0, 0, "bcnst", 0, 0);
      cpglab("Pulse phase", "Intensity", "");
      unsigned i, npts=prof->get_nbin();
      cpgsci(14);
      for (ic=0; ic < model.get_ncomponents(); ic++)
	plot(model, npts, true, ic);
      vector<float> xvals(npts);
      cpgsci(4);
      for (i=0; i < npts; i++)
	xvals[i] = i/((double)npts);
      if (line_plot)
	cpgline(npts, &xvals[0], prof->get_amps());
      else
	cpgbin(npts, &xvals[0], prof->get_amps(), 0);
      cpgsci(2);
      plot_difference(model, prof, line_plot);
      cpgsci(1);
      plot(model, npts, true);
    }
    
    
    while (interactive)
    {
      iterate = true;
      fit = false;

      cerr << "Enter command ('h' for help)" << endl;

      float curs_x=0, curs_y=0;
      char key = ' ';
      if (cpgband(6,0,0,0,&curs_x, &curs_y, &key) != 1 || key == 'q')
      {
	cerr << "Quitting ..." << endl;
	iterate = false;
	break;
      }

      if (key == 'h')
      {
	cerr << 
	  "\n"
	  "left click - add a component at cursor position\n"
	  "f key      - fit current set of components\n"
	  "q key      - quit\n"
	  "\n"
	  "After left click:\n"
	  "   1) any key    - select width of new component, then \n"
	  "   2) any key    - select height of new component \n"
	  "\n"
	  "   right click or <Escape> to abort addition of new component \n"
	     << endl;
	continue;
      }
      
      if (key == 'f')
      {
	cerr << "Fitting components" << endl;
	fit = true;
	break;
      }
      
      if (key == 'a' || key == 65)
      {
	cerr << "Adding component centred at pulse phase " << curs_x << endl;
	double centre = curs_x;

	cerr << "Select width and hit any key" << endl;
	if (cpgband(4, 1, curs_x, curs_y, &curs_x, &curs_y, &key) != 1) {
	  cerr << "paas: cpgband error" << endl;
	  return -1;
	}
	
	// right click or escape to abort
	if (key == 88 || key == 27)
	  continue;
	
	double width = fabs(centre - curs_x);
	
	cerr << "Select height and hit any key" << endl;
	if (cpgband(5,0,0,0,&curs_x, &curs_y, &key) != 1) {
	  cerr << "paas: cpgband error" << endl;
	  return -1;
	}
	
	// right click or escape to abort
	if (key == 88 || key == 27)
	  continue;
	
	double height = curs_y;
	
	model.add_component (centre, 0.25/(width*width), height, "");
	
	break;
	
      }
      
      cerr << "Unrecognized command" << endl;
      
    }
    
  }
  
  // write out if specified
  if (!model_filename_out.empty())
  {
    cerr << "paas: writing model to " << model_filename_out << endl;
    model.unload (model_filename_out.c_str());
  }
  
  Reference::To<Pulsar::Archive> copy = archive->clone();

  // write out standard
  model.evaluate (archive->get_Integration(0)->get_Profile(0,0)->get_amps(),
		  archive->get_nbin());

  write_details_to_file (model, copy, archive, details_filename);
  
  cerr << "paas: writing standard to " << std_filename << endl;
  archive->unload (std_filename);

  if (!pgdev.empty())
    cpgend();

  return 0;
}
catch (Error& error)
{
  cerr << error << endl;
  return -1;
}
Esempio n. 5
0
void psrspa::create_histograms ( Reference::To<Archive> archive )
{
  if ( verbose )
    cerr << "psrspa::create_histograms entered" << endl;
  // ensure Stokes parameters if creating polarisation histograms
  if ( create_polar_degree || create_polar_angle )
  {
    archive->convert_state ( Signal::Stokes );
    if ( verbose )
      cerr << "psrspa::create_histograms converted state of the archive to Stokes" << endl;
  }

  // auxillary vectors
  //vector< Estimate<float > > aux_vec_f;
  vector< Estimate<double > > aux_vec_d;

  // Full Stokes profile
  Reference::To<PolnProfile> profile;
  // Polarized flux
  Reference::To<Profile> P;
  P = new Profile;
  // Total flux
  float *T_amps = new float [ nbin ];
  float *P_amps = new float [ nbin ];

  unsigned bin_min, bin_max;

  if ( verbose && max_bscrunch > 1 )
    cerr << "psrspa::create_histograms entering the bscrunch loop for the " << archive->get_filename () << " archive" << endl;
  // TODO Uh, this should be fixed up. The first idea was to enable the phase resolved histograms to be bscrunch-aware as well, but I think it doesn't make much sense so I decided to keep only the max flux bscrunch aware. This can be done much more neatly probably in such a case
  for ( current_bscrunch = 1 ; current_bscrunch <= max_bscrunch ; current_bscrunch *= 2 )
  {
    // in each passage, we bscrunch by a factor of 2
    if ( current_bscrunch > 1 )
    {
      if ( verbose )
	cerr << "psrspa::create_histograms bscrunching the archive " << archive->get_filename () << " by a factor of 2" << endl;
      archive->bscrunch ( 2 );
    }
    if ( verbose )
      cerr << "psrspa::create_histograms entering the loop through subints of the " << archive->get_filename () << " archive" << endl;
    // loop through subints
    for ( unsigned isub = 0; isub < archive->get_nsubint (); isub++ )
    {
      if ( verbose )
	cerr << "psrspa::create_histograms creating necessary profiles for subint " << isub << " of " << archive->get_filename () << endl;
      if ( create_polar_angle || create_polar_degree && current_bscrunch == 1 )
      {
	profile = archive->get_Integration(isub)->new_PolnProfile(0);
	if ( verbose )
	  cerr << "psrspa::create_histograms retrieved PolnProfile for subint " << isub << " of " << archive->get_filename () << endl;
      }
      if ( create_polar_angle && current_bscrunch == 1 )
      {
	profile->get_orientation ( aux_vec_d, 0 );
	if ( verbose )
	  cerr << "psrspa::create_histograms retrieved polarisation angle for subint " << isub << " of " << archive->get_filename () << endl;
      }
      if ( create_polar_degree || create_flux || find_max_amp_in_range )
      {
	stats.set_profile ( archive->get_Profile ( isub, 0, 0 ) );
	b_sigma = sqrt ( stats.get_baseline_variance ().get_value () );
	T_amps = archive->get_Profile ( isub, 0, 0 )->get_amps (); 
	if ( verbose )
	  cerr << "psrspa::create_histograms retrieved total flux amps for subint " << isub << " of " << archive->get_filename () << endl;
	if ( create_polar_degree && current_bscrunch == 1 )
	{
	  profile->get_polarized ( P );
	  if ( verbose )
	    cerr << "psrspa::create_histograms retrieved polarized flux profile for subint "  << isub << " of " << archive->get_filename () << endl;
	  P_amps = P->get_amps ();
	  if ( verbose )
	    cerr << "psrspa::create_histograms retrieved polarized flux amps for subint " << isub << " of " << archive->get_filename () << endl;
	}
      }

      if ( verbose )
	cerr << "psrspa::create_histograms looping through the provided phase ranges for subint " << isub << " of " << archive->get_filename () << endl;
      unsigned curr_hist = 0;
      // loop through phase ranges
      for ( unsigned irange = 0; irange < phase_range.size () ; irange++ )
      {
	if ( irange%2 == 0)
	{
	  bin_min = unsigned( floor ( phase_range[irange] * float(nbin / current_bscrunch ) + 0.5 ) );
	  if ( verbose ) 
	    cerr << "psrspa::create_histograms set minimal bin to " << bin_min << endl;
	}
	else
	{
	  bin_max = unsigned( floor ( phase_range[irange] * float(nbin / current_bscrunch ) + 0.5 ) );
	  if ( bin_max == nbin )
	    bin_max = nbin - 1 ;
	  if ( verbose ) 
	    cerr << "psrspa::create_histograms set maximum bin to " << bin_max << endl;
	  // loop through bins in the given phase range
	  for ( unsigned ibin = bin_min; ibin <= bin_max; ibin++ )
	  {
	    if ( create_polar_angle && current_bscrunch == 1 )
	    {
	      int result = gsl_histogram_increment ( h_polar_angle_vec[curr_hist], aux_vec_d[ibin].get_value () / 180.0 * M_PI );
	      if ( result == GSL_EDOM )
	      {
		if ( dynamic_histogram )
		{
		  gsl_histogram *temp_hist = gsl_histogram_clone ( h_polar_angle_vec[curr_hist] );
		  h_flux_pr_vec[curr_hist] = update_histogram_range ( temp_hist, aux_vec_d[ibin].get_value () / 180.0 * M_PI );
		}
		else {
		  warn << "WARNING psrspa::create_histograms polarisation angle the histogram range for the bin " << ibin << " in the subint " << isub << " of archive " << archive->get_filename () << endl;
		}
	      }
	    }
	    if ( create_polar_degree && current_bscrunch == 1 )
	    {
	      // if P_amps[ibin] or T_amps[ibin] < 3 sigma, set polar degree to zero
	      int result = gsl_histogram_increment ( h_polar_degree_vec[curr_hist], ( fabs ( P_amps[ibin] ) < 3.0 * b_sigma || fabs ( T_amps[ibin] ) < 3.0 * b_sigma ) ? 0.0 : P_amps[ibin] / T_amps[ibin] );
	      if ( result == GSL_EDOM )
	      {
		if ( dynamic_histogram )
		{
		  gsl_histogram *temp_hist = gsl_histogram_clone ( h_polar_degree_vec[curr_hist] );
		  h_flux_pr_vec[curr_hist] = update_histogram_range ( temp_hist, ( fabs ( P_amps[ibin] ) < 3.0 * b_sigma || fabs ( T_amps[ibin] ) < 3.0 * b_sigma ) ? 0.0 : P_amps[ibin] / T_amps[ibin] );
		}
		else {
		warn << "WARNING psrspa::create_histograms polarisation degree outside the histogram range for the bin " << ibin << " in the subint " << isub << " of archive " << archive->get_filename () << endl;
		}
	      }
	    }
	    if ( create_flux && current_bscrunch == 1 )
	    {
	      int result = 0;
	      if ( log && T_amps[ibin] > 0.0 )
	      {
		result = gsl_histogram_increment ( h_flux_pr_vec[curr_hist], log10f ( T_amps[ibin] ) );
	      }
	      else if ( !log )
	      {
		result = gsl_histogram_increment ( h_flux_pr_vec[curr_hist], T_amps[ibin] );
	      }
	      if ( result == GSL_EDOM )
	      {
		if ( dynamic_histogram )
		{
		  gsl_histogram *temp_hist = gsl_histogram_clone ( h_flux_pr_vec[curr_hist] );
		  h_flux_pr_vec[curr_hist] = update_histogram_range ( temp_hist, log ? log10f ( T_amps[ibin] ) : T_amps[ibin] );
		}
		else {
		  warn << "WARNING psrspa::create_histograms phase resolved flux outside the histogram range for the bin " << ibin << " in the subint " << isub << " of archive " << archive->get_filename () << " flux = " << T_amps[ibin] << endl;
		}
	      }
	    }
	    // increment the histogram id
	    curr_hist ++;
	  } // loop through bins in the given phase range
	  // find the maximum amplitude in given phase range.
	  if ( find_max_amp_in_range )
	  {
	    if ( verbose )
	      cerr << "psrspa::create_histograms finding the maximum amplitude in the phase range " << irange << " of the subint " << isub << " (archive " << archive->get_filename () << ")" << endl;
	    int max_bin = archive->get_Profile ( isub, 0, 0 )->find_max_bin ( (int)bin_min, (int)bin_max );
	    identifier current_identifier = { archive->get_filename (), isub, bin_min, bin_max, current_bscrunch,
	      (unsigned)max_bin, T_amps[max_bin], b_sigma }; 
	    max_amp_info.push_back ( current_identifier );
	  } // find maximal amplitude in the given phase range
	} // handle given phase range
      } // loop through phase ranges
    } // loop through subints
  } // bscrunch loop
  if ( verbose )
    cerr << "psrspa::create_histograms finished" << endl;
}// create_histograms