ChanTrigger.h
1 #ifndef CHANTRIGGER_HH_
2 #define CHANTRIGGER_HH_
3 
5 //class ChanTrigger:
7 
8 class RX;
9 class Anita;
10 using std::vector;
11 class Balloon;
12 class IceModel;
13 class Settings;
14 class Screen;
15 class GlobalTrigger;
16 
18 class ChanTrigger {
19 
20  private:
21 
34  double ADCCountstoPowerThreshold(Anita *anita1, int ipol, int iant);
35 
36  const static int NSURF=9;
37  const static int NSURFPLUSONE=10;
38  const static int NSURFMINUSONE=8;
39  const static int NCHANNELS=32;
40  const static int NPOINTS=4073;
41  static const unsigned NFOUR = 1024;
42  static const unsigned HALFNFOUR = 512;
43 
44  double thisrate;
45  double thispowerthresh;
46  double e_component;
47  double h_component;
48  double n_component;
49  double e_component_kvector;
50  double h_component_kvector;
51  double n_component_kvector;
52  double hitangle_e;
53  double hitangle_h;
54 
55  public:
56 
58  ChanTrigger();
59 
60 
62 
67  void InitializeEachBand(Anita *anita1);
68 
70 
86  void ApplyAntennaGain(Settings *settings1, Anita *anita1, Balloon *bn1, Screen *panel1, int ant, Vector &n_eplane, Vector &n_hplane, Vector &n_normal);
87 
89 
99  void TriggerPath(Settings *settings1, Anita *anita1, int ant, Balloon *bn1);
100 
102 
112  void DigitizerPath(Settings *settings1, Anita *anita1, int ant, Balloon *bn1);
113 
115 
126  void TimeShiftAndSignalFluct(Settings *settings1, Anita *anita1, int ilayer, int ifold, double volts_rx_rfcm_lab_e_all[48][512], double volts_rx_rfcm_lab_h_all[48][512]);
127 
129 
137  static void ConvertEHtoLREfield(double,double,double&,double&);
138 
140 
148  static void ConvertEHtoLREnergy(double,double,double&,double&);
149 
151 
161  void L1Trigger(Anita *anita1,double timedomain_output_1[5][Anita::NFOUR],double timedomain_output_2[5][Anita::NFOUR],double powerthreshold[2][5],int *channels_passing_e_forglob,int *channels_passing_h_forglob,int &npass);
162 
163 
165 
169  double getRate();
170 
172 
185  double rateToThreshold(double rate, int band);
186 
188 
198  static double GetNoise(Settings *settings1,double altitude_bn,double geoid,double theta,double bw,double temp);
199 
201 
214  void WhichBandsPass(Settings *settings1, Anita *anita1, GlobalTrigger *globaltrig1, Balloon *bn1, int ilayer, int ifold, double dangle, double emfrac, double hadfrac, double thresholds[2][5]);
215 
217 
227  void WhichBandsPassTrigger1(Settings *settings1, Anita *anita1, GlobalTrigger *globaltrig1, Balloon *bn1, int ilayer, int ifold, double thresholds[2][5]);
228 
230 
243  void WhichBandsPassTrigger2(Settings *settings1, Anita *anita1, GlobalTrigger *globaltrig1, Balloon *bn1, int ilayer, int ifold, double dangle, double emfrac, double hadfrac, double thresholds[2][5]);
244 
246 
250  static double FindPeak(double *waveform,int n);
251 
252 
254 
261  void GetThresholds(Settings *settings1,Anita *anita1,int ilayer,double thresholds[2][5]); // get thresholds for this layer
262 
264 
278  void DiodeConvolution(Settings *settings1, Anita *anita1, GlobalTrigger *globaltrig1, int ilayer, int ifold, double mindiodeconvl[5], double onediodeconvl[5], double psignal[5][Anita::NFOUR], double timedomain_output[5][Anita::NFOUR], double timedomain_output_justNoise[5][Anita::NFOUR], int ibinshift, int ipol, double thresholds[2][5]);
279 
280 
282 
289  void addToChannelSums(Settings *settings1,Anita *anita1,int ibw,int k);
290 
292 
302  static int IsItUnmasked(unsigned short surfTrigBandMask[9][2],int ibw,int ilayer, int ifold, int ipol);
303 
305 
316  void applyImpulseResponseDigitizer(Settings *settings1, Anita *anita1, int nPoints, int ant, double *x, double y[512], bool pol);
317 
319 
329  void applyImpulseResponseTrigger(Settings *settings1, Anita *anita1, int ant, double y[512], double *vhz, bool pol);
330 
332 
342  void getNoiseFromFlight(Settings* settings1, Anita* anita1, int ant, bool also_digi = true);
343 
345 
352  void injectImpulseAfterAntenna(Anita *anita1, int ant);
353 
355 
361  TGraph *getPulserAtAMPA(Anita *anita1, int ant);
362 
364 
369  void saveTriggerWaveforms(Anita *anita1, double sig0[48], double sig1[48], double noise0[48], double noise1[48]);
370 
372 
377  void saveDigitizerWaveforms(Anita *anita1, double sig0[48], double sig1[48], double noise0[48], double noise1[48]);
378 
380 
390  void injectImpulseAtSurf(Anita *anita1, double volts_triggerPath_e[Anita::HALFNFOUR], double volts_triggerPath_h[Anita::HALFNFOUR], int ant);
391 
393 
401  void calculateCW(Anita *anita1, double frequency, double phase, double amplitude);
403 
411  double applyButterworthFilter(double ff, double ampl, int notchStatus[3]);
412 
413 
414  double vhz_rx[2][5][Anita::NFREQ];
415  double volts_rx_forfft[2][5][Anita::HALFNFOUR];
416  vector<int> flag_e[5];
417  vector<int> flag_h[5];
418  double bwslice_volts_pol0[5];
419  double bwslice_volts_pol1[5];
424  double bwslice_volts_pole[5];
426  double bwslice_volts_polh[5];
428  double volts_rx_rfcm_lab[2][Anita::HALFNFOUR];
429  double volts_rx_rfcm_lab_all[2][48][Anita::HALFNFOUR];
430  double volts_rx_rfcm[2][Anita::HALFNFOUR];
431  double justNoise_digPath[2][Anita::HALFNFOUR];
432  double justNoise_trigPath[2][Anita::HALFNFOUR];
433  double cw_digPath[2][Anita::HALFNFOUR];
434  double justSig_trigPath[2][Anita::HALFNFOUR];
435  double justSig_digPath[2][Anita::HALFNFOUR];
436 
437  // these are filled for triggerscheme==0 and triggerscheme==1
438  // frequency domain voltage and energy based
443  vector<double> vsignal_eachband[2];
444  vector<double> vthreshold_eachband[2];
445  vector<double> vnoise_eachband[2];
446  vector<int> vpasses_eachband[2];
447 
448  double v_banding_rfcm[2][5][Anita::NFREQ];
449  double v_banding_rfcm_forfft[2][5][HALFNFOUR];
450  double vm_banding_rfcm_forfft[2][5][HALFNFOUR];
451  double vm_banding_rfcm_forfft_justNoise[2][5][HALFNFOUR];
452  double v_banding_rfcm_forfft_temp[2][5][HALFNFOUR];
453  double integral_vmmhz;
454  int unwarned;
455 }; //class ChanTrigger
456 
457 #endif
458 
void injectImpulseAtSurf(Anita *anita1, double volts_triggerPath_e[Anita::HALFNFOUR], double volts_triggerPath_h[Anita::HALFNFOUR], int ant)
Inject pulse at the surf (used for trigger efficiency scans)
double cw_digPath[2][Anita::HALFNFOUR]
For digitizer path, time domain cw.
Definition: ChanTrigger.h:433
void WhichBandsPassTrigger2(Settings *settings1, Anita *anita1, GlobalTrigger *globaltrig1, Balloon *bn1, int ilayer, int ifold, double dangle, double emfrac, double hadfrac, double thresholds[2][5])
Which bands passes the trigger (for trigger scheme larger than 2)
Definition: ChanTrigger.cc:285
double justNoise_digPath[2][Anita::HALFNFOUR]
For digitizer path, time domain noise from flight.
Definition: ChanTrigger.h:431
double bwslice_energy_polh[5]
Square the sum of voltage for each slice in bandwidth for h polarization. The 5th element is the full...
Definition: ChanTrigger.h:427
ChanTrigger()
Channel trigger constructur.
double justNoise_trigPath[2][Anita::HALFNFOUR]
For trigger path, time domain noise from flight.
Definition: ChanTrigger.h:432
static const int NBANDS_MAX
max number of bands
Definition: anita.hh:50
double getRate()
Returns the thisrate variable value (in MHz)
Global Trigger.
Definition: GlobalTrigger.h:12
void InitializeEachBand(Anita *anita1)
Initialize trigger bands.
Definition: ChanTrigger.cc:711
double bwslice_volts_pol0[5]
Sum voltage for each slice in bandwidth for the lcp polarization.
Definition: ChanTrigger.h:418
double bwslice_energy_pol0[5]
Square the sum of voltage for each slice in bandwidth for the 0th polarization.
Definition: ChanTrigger.h:420
void TriggerPath(Settings *settings1, Anita *anita1, int ant, Balloon *bn1)
Apply trigger path.
Definition: ChanTrigger.cc:871
void WhichBandsPassTrigger1(Settings *settings1, Anita *anita1, GlobalTrigger *globaltrig1, Balloon *bn1, int ilayer, int ifold, double thresholds[2][5])
Which bands passes the trigger (for trigger scheme 0 and 1)
Definition: ChanTrigger.cc:103
double rateToThreshold(double rate, int band)
Calculates the trigger threshold for an antenna via a fit from the "singles" rate and band identifier...
double bwslice_energy_pole[5]
Square the sum of voltage for each slice in bandwidth for e polarization. The 5th element is the full...
Definition: ChanTrigger.h:425
vector< double > vthreshold_eachband[2]
Threshold in each band.
Definition: ChanTrigger.h:444
void L1Trigger(Anita *anita1, double timedomain_output_1[5][Anita::NFOUR], double timedomain_output_2[5][Anita::NFOUR], double powerthreshold[2][5], int *channels_passing_e_forglob, int *channels_passing_h_forglob, int &npass)
The L1 trigger of the Anita trigger scheme.
vector< int > flag_e[5]
Which bands pass trigger e.
Definition: ChanTrigger.h:416
double bwslice_volts_pol0_em[5]
Component of the voltage that comes from the em shower for 0th polarization.
Definition: ChanTrigger.h:422
void applyImpulseResponseDigitizer(Settings *settings1, Anita *anita1, int nPoints, int ant, double *x, double y[512], bool pol)
Apply impulse response to digitizer path.
Definition: screen.hh:22
void ApplyAntennaGain(Settings *settings1, Anita *anita1, Balloon *bn1, Screen *panel1, int ant, Vector &n_eplane, Vector &n_hplane, Vector &n_normal)
Apply the antenna gain.
Definition: ChanTrigger.cc:743
double vhz_rx[2][5][Anita::NFREQ]
Array of amplitudes in the Fourier domain (V/Hz) after the antenna gain. Indeces stand for [ipol][iba...
Definition: ChanTrigger.h:414
void applyImpulseResponseTrigger(Settings *settings1, Anita *anita1, int ant, double y[512], double *vhz, bool pol)
Apply impulse response to trigger path.
void TimeShiftAndSignalFluct(Settings *settings1, Anita *anita1, int ilayer, int ifold, double volts_rx_rfcm_lab_e_all[48][512], double volts_rx_rfcm_lab_h_all[48][512])
Time shift and fluctuate signal.
double bwslice_volts_pol1[5]
Sum voltage for each slice in bandwidth for the rcp polarization.
Definition: ChanTrigger.h:419
static void ConvertEHtoLREnergy(double, double, double &, double &)
Convert E and H to left and right energy.
Definition: ChanTrigger.cc:45
double volts_rx_forfft[2][5][Anita::HALFNFOUR]
Array of time domain after the antenna gain. Indeces stand for [ipol][iband][itime].
Definition: ChanTrigger.h:415
vector< double > vnoise_eachband[2]
Noise in each band.
Definition: ChanTrigger.h:445
void saveTriggerWaveforms(Anita *anita1, double sig0[48], double sig1[48], double noise0[48], double noise1[48])
Save signal and noise waveforms at trigger.
double noise_eachband[2][Anita::NBANDS_MAX]
Noise in each band.
Definition: ChanTrigger.h:441
double volts_rx_rfcm_lab[2][Anita::HALFNFOUR]
For digitizer path, time domain voltage vs. time after rx, rfcm&#39;s and lab.
Definition: ChanTrigger.h:428
Class that handles the channel trigger.
Definition: ChanTrigger.h:18
double applyButterworthFilter(double ff, double ampl, int notchStatus[3])
Apply Butterworth Filter.
Reads in and stores input settings for the run.
Definition: Settings.h:35
static double GetNoise(Settings *settings1, double altitude_bn, double geoid, double theta, double bw, double temp)
Get Noise.
void saveDigitizerWaveforms(Anita *anita1, double sig0[48], double sig1[48], double noise0[48], double noise1[48])
Save signal and noise waveforms at digitizer.
void getNoiseFromFlight(Settings *settings1, Anita *anita1, int ant, bool also_digi=true)
Add noise from ANITA-3 flight to the time domain waveforms.
double v_banding_rfcm[2][5][Anita::NFREQ]
This is Volts/m as a function of frequency after rfcm&#39;s and banding.
Definition: ChanTrigger.h:448
double integral_vmmhz
Electric field integral.
Definition: ChanTrigger.h:453
static void ConvertEHtoLREfield(double, double, double &, double &)
Convert E and H to left and right e field.
Definition: ChanTrigger.cc:39
vector< int > vpasses_eachband[2]
Whether the signal passes or not each band.
Definition: ChanTrigger.h:446
void DigitizerPath(Settings *settings1, Anita *anita1, int ant, Balloon *bn1)
Apply digitizer path.
Contains everything about positions within payload and signals it sees for each event, in both the trigger and signal paths.
Definition: anita.hh:32
void addToChannelSums(Settings *settings1, Anita *anita1, int ibw, int k)
Increment the volts in each band.
double justSig_trigPath[2][Anita::HALFNFOUR]
Just signal in trigger path.
Definition: ChanTrigger.h:434
double threshold_eachband[2][Anita::NBANDS_MAX]
Threshold in each band.
Definition: ChanTrigger.h:440
double bwslice_volts_pol1_em[5]
Component of the voltage that comes from the em shower for 1st polarization.
Definition: ChanTrigger.h:423
static double FindPeak(double *waveform, int n)
Find peak voltage of a waveform.
double v_banding_rfcm_forfft[2][5][HALFNFOUR]
Starts out as V/s vs. freq after banding, rfcm, after fft it is V vs. t.
Definition: ChanTrigger.h:449
double justSig_digPath[2][Anita::HALFNFOUR]
Just signal in trigger path.
Definition: ChanTrigger.h:435
int passes_eachband[2][Anita::NBANDS_MAX]
Whether the signal passes or not each band.
Definition: ChanTrigger.h:442
double vm_banding_rfcm_forfft[2][5][HALFNFOUR]
Tunnel diode input (signal + noise)
Definition: ChanTrigger.h:450
int unwarned
Whether we have warned the user about resetting thresholds when they are beyond the measured bounds...
Definition: ChanTrigger.h:454
double bwslice_volts_polh[5]
Sum voltage for each slice in bandwidth for the h polarization.
Definition: ChanTrigger.h:426
TGraph * getPulserAtAMPA(Anita *anita1, int ant)
Get time domain graph of pulse at AMPA (used for trigger efficiency scans)
void WhichBandsPass(Settings *settings1, Anita *anita1, GlobalTrigger *globaltrig1, Balloon *bn1, int ilayer, int ifold, double dangle, double emfrac, double hadfrac, double thresholds[2][5])
Which bands passes the trigger.
Definition: ChanTrigger.cc:54
void DiodeConvolution(Settings *settings1, Anita *anita1, GlobalTrigger *globaltrig1, int ilayer, int ifold, double mindiodeconvl[5], double onediodeconvl[5], double psignal[5][Anita::NFOUR], double timedomain_output[5][Anita::NFOUR], double timedomain_output_justNoise[5][Anita::NFOUR], int ibinshift, int ipol, double thresholds[2][5])
Apply the diode convolution.
Definition: ChanTrigger.cc:546
static int IsItUnmasked(unsigned short surfTrigBandMask[9][2], int ibw, int ilayer, int ifold, int ipol)
Returns whether the indicated antenna and band are "masked".
vector< double > vsignal_eachband[2]
Signal in each band.
Definition: ChanTrigger.h:443
Handles everything related to balloon positions, payload orientation over the course of a flight...
Definition: balloon.hh:30
This class represents a three-vector. Operators are overloaded to provide for the familiar operations...
Definition: vector.hh:27
void GetThresholds(Settings *settings1, Anita *anita1, int ilayer, double thresholds[2][5])
Sets the threshold values based on which payload and where the antenna is located physically...
double bwslice_volts_pole[5]
Sum voltage for each slice in bandwidth for the e polarization.
Definition: ChanTrigger.h:424
double volts_rx_rfcm_lab_all[2][48][Anita::HALFNFOUR]
For digitizer path, time domain voltage vs. time after rx, rfcm&#39;s and lab.
Definition: ChanTrigger.h:429
void injectImpulseAfterAntenna(Anita *anita1, int ant)
Inject pulse after the antenna (used for trigger efficiency scans)
vector< int > flag_h[5]
Which bands pass trigger h.
Definition: ChanTrigger.h:417
Definition: rx.hpp:21
double bwslice_energy_pol1[5]
Square the sum of voltage for each slice in bandwidth for the 1st polarization.
Definition: ChanTrigger.h:421
double v_banding_rfcm_forfft_temp[2][5][HALFNFOUR]
Use for the averaging over 10 neighboring bins.
Definition: ChanTrigger.h:452
double signal_eachband[2][Anita::NBANDS_MAX]
Signal in each band.
Definition: ChanTrigger.h:439
double vm_banding_rfcm_forfft_justNoise[2][5][HALFNFOUR]
Tunnel diode input (just noise)
Definition: ChanTrigger.h:451
double volts_rx_rfcm[2][Anita::HALFNFOUR]
For digitizer path, time domain voltage vs. time after rx, rfcm&#39;s.
Definition: ChanTrigger.h:430
void calculateCW(Anita *anita1, double frequency, double phase, double amplitude)
Add CW.
Ice thicknesses and water depth.
Definition: icemodel.hh:88