Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
detector.cc
Go to the documentation of this file.
1 // ++++++++++++++++++++++++++++++++++++++++++++++
2 // S. Klimenko, University of Florida
3 // WAT detector class
4 // ++++++++++++++++++++++++++++++++++++++++++++++
5 
6 #define DETECTOR_CC
7 #include "detector.hh"
8 #include "Meyer.hh"
9 #include "Symlet.hh"
10 
11 #include "TVector3.h"
12 #include "TRotation.h"
13 #include "Math/Rotation3D.h"
14 #include "Math/Vector3Dfwd.h"
15 
16 using namespace ROOT::Math;
17 using namespace std;
18 
19 // All LIGO coordinates are verbatim from LIGO-P000006-D-E:
20 // Rev. Sci. Instrum., Vol. 72, No. 7, July 2001
21 
22 // LIGO Hanford 4k interferometer
23 extern const double _rH1[3] = {-2.161414928e6, -3.834695183e6, 4.600350224e6};
24 extern const double _eH1x[3] = {-0.223891216, 0.799830697, 0.556905359};
25 extern const double _eH1y[3] = {-0.913978490, 0.026095321, -0.404922650};
26 
27 // LIGO Hanford 2k interferometer
28 extern const double _rH2[3] = {-2.161414928e6, -3.834695183e6, 4.600350224e6};
29 extern const double _eH2x[3] = {-0.223891216, 0.799830697, 0.556905359};
30 extern const double _eH2y[3] = {-0.913978490, 0.026095321, -0.404922650};
31 
32 // LIGO Livingston interferometer
33 extern const double _rL1[3] = {-7.427604192e4, -5.496283721e6, 3.224257016e6};
34 extern const double _eL1x[3] = {-0.954574615, -0.141579994, -0.262187738};
35 extern const double _eL1y[3] = { 0.297740169, -0.487910627, -0.820544948};
36 
37 // GEO coordinates on
38 // http://www.geo600.uni-hannover.de/geo600/project/location.html
39 // x any y arms swaped compare to
40 // Anderson, Brady, Creighton, and Flanagan, PRD 63 042003, 2001,
41 
42 // GEO-600 interferometer
43 extern const double _rG1[3] = {3.8563112e6, 6.665978e5, 5.0196406e6};
44 extern const double _eG1x[3] = {-0.445184239, 0.866534205, 0.225675575};
45 extern const double _eG1y[3] = {-0.626000687, -0.552167273, 0.550667271};
46 
47 // Virgo interferometer
48 extern const double _rV1[3] = {4.5463741e6, 8.429897e5, 4.378577e6};
49 extern const double _eV1x[3] = {-0.700458309, 0.208487795, 0.682562083};
50 extern const double _eV1y[3] = {-0.053791331, -0.969082169, 0.240803326};
51 
52 // TAMA interferometer
53 extern const double _rT1[3] = {-3.946409e6, 3.366259e6, 3.6991507e6};
54 extern const double _eT1x[3] = {0.648969405, 0.760814505, 0};
55 extern const double _eT1y[3] = {-0.443713769, 0.378484715, -0.812322234};
56 
57 // LCGT interferometer
58 extern const double _rJ1[3] = {-3.7685777e6, 3.49218552e6, 3.76722999e6};
59 extern const double _eJ1x[3] = {-0.5014606643, -0.824380550, 0.262552681483};
60 extern const double _eJ1y[3] = {0.6313793503,-0.1412130545258,0.762508353526};
61 
62 // AIGO interferometer
63 extern const double _rA1[3] = {-2.3567784e6, 4.8970238e6, -3.3173147e6};
64 extern const double _eA1x[3] = {-9.01077021322091554e-01, -4.33659084587544319e-01, 0};
65 extern const double _eA1y[3] = {-2.25940560005277846e-01, 4.69469807139026196e-01, 8.53550797275327455e-01};
66 
67 // AURIGA bar
68 extern const double _rO1[3] = {4.392467e6, 0.9295086e6, 4.515029e6};
69 extern const double _eO1x[3] = {-0.644504130, 0.573655377, 0.50550364};
70 extern const double _eO1y[3] = {0., 0., 0.};
71 
72 // NAUTILUS bar
73 extern const double _rN1[3] = {4.64410999868e6, 1.04425342477e6, 4.23104713307e6};
74 extern const double _eN1x[3] = {-0.62792641437, 0.56480832712, 0.53544371484};
75 extern const double _eN1y[3] = {0., 0., 0.};
76 
77 // EXPLORER bar
78 extern const double _rE1[3] = {4.37645395452e6, 4.75435044067e5, 4.59985274450e6};
79 extern const double _eE1x[3] = {-0.62792641437, 0.56480832712, 0.53544371484};
80 extern const double _eE1y[3] = {0., 0., 0.};
81 
82 
83 // see also detector coordinates used by GravEn
84 // https://gravity.psu.edu/~s4/sims/BurstMDC/Validation/getifo.m
85 
86 ClassImp(detector)
87 
88 // constructors
89 
91 {
92  for(int i=0; i<3; i++){
93  this->Rv[i] = _rL1[i];
94  this->Ex[i] = _eL1x[i];
95  this->Ey[i] = _eL1y[i];
96  }
97  init();
98  this->sHIFt = 0.;
99  this->null = 0.;
100  this->sSNR = 0.;
101  this->xSNR = 0.;
102  this->ekXk = 0.;
103  this->ifoID = 0;
104  this->rate = 16384.;
105  this->TFmap.rate(4096);
106  detectorParams null_dP = {"",0.0,0,0,0,0,0,0};
107  this->dP = null_dP;
108  this->polarization = TENSOR;
109  this->wfSAVE = 0;
110 }
111 
112 detector::detector(char* name, double t)
113 {
114  const double* pRv;
115  const double* pEx;
116  const double* pEy;
117 
118  bool xifo=false;
119 
120  if(strstr(name,"L1")) {pRv=_rL1; pEx=_eL1x; pEy=_eL1y;xifo=true;}
121  if(strstr(name,"H1")) {pRv=_rH1; pEx=_eH1x; pEy=_eH1y;xifo=true;}
122  if(strstr(name,"H2")) {pRv=_rH2; pEx=_eH2x; pEy=_eH2y;xifo=true;}
123  if(strstr(name,"G1")) {pRv=_rG1; pEx=_eG1x; pEy=_eG1y;xifo=true;}
124  if(strstr(name,"T1")) {pRv=_rT1; pEx=_eT1x; pEy=_eT1y;xifo=true;}
125  if(strstr(name,"V1")) {pRv=_rV1; pEx=_eV1x; pEy=_eV1y;xifo=true;}
126  if(strstr(name,"A1")) {pRv=_rA1; pEx=_eA1x; pEy=_eA1y;xifo=true;}
127  if(strstr(name,"A2")) {pRv=_rA1; pEx=_eA1x; pEy=_eA1y;xifo=true;}
128  if(strstr(name,"Virgo")) {pRv=_rV1; pEx=_eV1x; pEy=_eV1y;xifo=true;}
129  if(strstr(name,"VIRGO")) {pRv=_rV1; pEx=_eV1x; pEy=_eV1y;xifo=true;}
130  if(strstr(name,"GEO")) {pRv=_rG1; pEx=_eG1x; pEy=_eG1y;xifo=true;}
131  if(strstr(name,"TAMA")) {pRv=_rT1; pEx=_eT1x; pEy=_eT1y;xifo=true;}
132  if(strstr(name,"LCGT")) {pRv=_rJ1; pEx=_eJ1x; pEy=_eJ1y;xifo=true;}
133  if(strstr(name,"J1")) {pRv=_rJ1; pEx=_eJ1x; pEy=_eJ1y;xifo=true;}
134  if(strstr(name,"O1")) {pRv=_rO1; pEx=_eO1x; pEy=_eO1y;xifo=true;}
135  if(strstr(name,"N1")) {pRv=_rN1; pEx=_eN1x; pEy=_eN1y;xifo=true;}
136  if(strstr(name,"E1")) {pRv=_rE1; pEx=_eE1x; pEy=_eE1y;xifo=true;}
137 
138  if(!xifo) {
139  cout << "detector::detector - Error : detector " << name
140  << " is not in present in the builtin list" << endl;
141  exit(1);
142  }
143 
144  sprintf(this->Name,"%s",name);
145  SetName(name);
146  for(int i=0; i<3; i++){
147  this->Rv[i] = pRv[i];
148  this->Ex[i] = pEx[i];
149  this->Ey[i] = pEy[i];
150  }
151  if(strstr(name,"A2")) this->rotate(45); // rotate arms for A2
152  init(); // fill detector tensor.
153  this->sHIFt = t;
154  this->null = 0.;
155  this->sSNR = 0.;
156  this->xSNR = 0.;
157  this->ekXk = 0.;
158  this->ifoID = 0;
159  this->rate = 16384.;
160  this->TFmap.rate(4096);
161  detectorParams null_dP = {"",0.0,0,0,0,0,0,0};
162  this->dP = null_dP;
163  this->polarization = TENSOR;
164  this->wfSAVE = 0;
165 }
166 
168 {
169  double rad2deg = 180./TMath::Pi();
170  double deg2rad = TMath::Pi()/180.;
171 
172  double uRv[3];
173  double uEx[3];
174  double uEy[3];
175 
176  GeodeticToGeocentric(dP.latitude*deg2rad,dP.longitude*deg2rad,dP.elevation,uRv[0],uRv[1],uRv[2]);
177  GetCartesianComponents(uEx,dP.AltX*deg2rad,dP.AzX*deg2rad, dP.latitude*deg2rad,dP.longitude*deg2rad);
178  GetCartesianComponents(uEy,dP.AltY*deg2rad,dP.AzY*deg2rad, dP.latitude*deg2rad,dP.longitude*deg2rad);
179 
180  sprintf(this->Name,"%s",dP.name);
181  SetName(dP.name);
182  for(int i=0; i<3; i++){
183  this->Rv[i] = uRv[i];
184  this->Ex[i] = uEx[i];
185  this->Ey[i] = uEy[i];
186  }
187  init(); // fill detector tensor.
188  this->sHIFt = t;
189  this->null = 0.;
190  this->sSNR = 0.;
191  this->xSNR = 0.;
192  this->ekXk = 0.;
193  this->ifoID = 0;
194  this->rate = 16384.;
195  this->TFmap.rate(4096);
196  this->dP = dP;
197  this->polarization = TENSOR;
198  this->wfSAVE = 0;
199 }
200 
202 {
203  // -----------------------------------------
204  // user define detector
205  // -----------------------------------------
206  if(strlen(this->dP.name)>0) return this->dP;
207 
208  // -----------------------------------------
209  // builtin detector
210  // -----------------------------------------
211  double rad2deg = 180./TMath::Pi();
212  double deg2rad = TMath::Pi()/180.;
213 
214  XYZVector iRv(this->Rv[0],this->Rv[1],this->Rv[2]);
215  TVector3 vRv(iRv.X(),iRv.Y(),iRv.Z());
216 
217  // Ex angle respect to east direction
218  XYZVector iEx(this->Ex[0],this->Ex[1],this->Ex[2]);
219  XYZVector iEy(this->Ey[0],this->Ey[1],this->Ey[2]);
220  XYZVector iEZ(0.,0.,1.); // Zeta cartesian axis
221 
222  TVector3 vEx(iEx.X(),iEx.Y(),iEx.Z());
223  TVector3 vEy(iEy.X(),iEy.Y(),iEy.Z());
224  TVector3 vEZ(iEZ.X(),iEZ.Y(),iEZ.Z());
225 
226  vEx*=1./vEx.Mag();
227  vEy*=1./vEy.Mag();
228  vEZ*=1./vEZ.Mag();
229  vRv*=1./vRv.Mag();
230 
231  TVector3 vEe=vEZ.Cross(vRv); // vEe point to east in the local detector frame
232  vEe*=1./vEe.Mag();
233 
234  double cosExAngle=vEe.Dot(vEx); // ExAngle is the angle of Ex respect to East counter-clockwise
235  if(cosExAngle>1.) cosExAngle=1.;
236  if(cosExAngle<-1.) cosExAngle=-1.;
237  double cosEyAngle=vEe.Dot(vEy); // EyAngle is the angle of Ey respect to Eas counter-clockwiset
238  if(cosEyAngle>1.) cosEyAngle=1.;
239  if(cosEyAngle<-1.) cosEyAngle=-1.;
240 
241  double ExAngle=acos(cosExAngle)*rad2deg;
242  double EyAngle=acos(cosEyAngle)*rad2deg;
243 
244  // fix sign of ExAngle,EyAngle
245  TVector3 vEn;
246  vEn=vEe.Cross(vEx);
247  if(vEn.Dot(vRv)<0) ExAngle*=-1;
248  vEn=vEe.Cross(vEy);
249  if(vEn.Dot(vRv)<0) EyAngle*=-1;
250 
251  // Convert ExAngle to the angle of Ex respect to North clockwise
252  ExAngle=fmod(90-ExAngle,360.);
253  // Convert EyAngle to the angle of Ey respect to North clockwise
254  EyAngle=fmod(90-EyAngle,360.);
255 
256  double latitude,longitude,elevation;
257  GeocentricToGeodetic(iRv.X(),iRv.Y(),iRv.Z(),latitude,longitude,elevation);
258 
259  detectorParams idP;
260 
261  sprintf(idP.name,"%s",this->Name);
262  idP.latitude = latitude*rad2deg;
263  idP.longitude = longitude*rad2deg;
264  idP.elevation = elevation;
265  idP.AltX = 0.;
266  idP.AzX = ExAngle;
267  idP.AltY = 0.;
268  idP.AzY = EyAngle;
269 
270  return idP;
271 }
272 
273 // rotate arms in the detector plane by angle a in degrees counter-clockwise
274 void detector::rotate(double a)
275 {
276  double ax[3];
277  double ay[3];
278  double si = sin(a*PI/180.);
279  double co = cos(a*PI/180.);
280 
281  for(int i=0; i<3; i++) {
282  ax[i] = this->Ex[i];
283  ay[i] = this->Ey[i];
284  }
285 /*
286  for(int i=0; i<3; i++) {
287  this->Ex[i] = ax[i]*co+ay[i]*si;
288  this->Ey[i] = ay[i]*co-ax[i]*si;
289  }
290 */
291  double aww=0.;
292  double aw[3];
293  double axy=0.;
294 
295  for(int i=0; i<3; i++) axy+=ax[i]*ay[i];
296 
297  // compute vector aw in the plane Ex,Ey ortogonal to Ex and rotate Ex
298  aww=0.;
299  for(int i=0; i<3; i++) {aw[i]=-ax[i]*axy+ay[i]; aww+=aw[i]*aw[i];}
300  for(int i=0; i<3; i++) aw[i]/=sqrt(aww);
301  for(int i=0; i<3; i++) this->Ex[i]=ax[i]*co+aw[i]*si; // rotate Ex
302 
303  // compute vector aw in the plane Ex,Ey ortogonal to Ey and rotate Ey
304  aww=0.;
305  for(int i=0; i<3; i++) {aw[i]=-ax[i]+ay[i]*axy; aww+=aw[i]*aw[i];}
306  for(int i=0; i<3; i++) aw[i]/=sqrt(aww);
307  for(int i=0; i<3; i++) this->Ey[i]=ay[i]*co+aw[i]*si; // rotate Ey
308 
309  init(); // fill detector tensor.
310 
311  // update user define detector
312  if(strlen(this->dP.name)>0) {
313  this->dP.AzX = fmod(this->dP.AzX-a,360.);
314  this->dP.AzY = fmod(this->dP.AzY-a,360.);
315  }
316 }
317 
318 
320 {
321  *this = value;
322 }
323 
324 // destructor
325 
327 
328  int n;
329 
330  n = IWFP.size();
331  for (int i=0;i<n;i++) {
333  delete wf;
334  }
335  IWFP.clear();
336  IWFID.clear();
337 
338  n = RWFP.size();
339  for (int i=0;i<n;i++) {
341  delete wf;
342  }
343  RWFP.clear();
344  RWFID.clear();
345 
346 }
347 
348 //: operator =
349 //: !!! not fully implemented
350 
352 {
353  sprintf(Name,"%s",value.Name);
354  SetName(Name);
355  for(int i=0; i<3; i++){
356  this->Rv[i] = value.Rv[i];
357  this->Ex[i] = value.Ex[i];
358  this->Ey[i] = value.Ey[i];
359  }
360  init();
361 
362  tau = value.tau;
363  mFp = value.mFp;
364  mFx = value.mFx;
365 
366  TFmap = value.TFmap;
367  waveForm = value.waveForm;
368  sHIFt = value.sHIFt;
369  null = value.null;
370  nRMS = value.nRMS;
371  nVAR = value.nVAR;
372 
373  dP = value.dP;
374 
375  polarization = value.polarization;
376 
377  wfSAVE = value.wfSAVE;
378 
379  return *this;
380 }
381 
382 
384 {
385  double tsRate =TFmap.wavearray<double>::rate();
386  TFmap = value;
387  waveForm.resize(size_t(tsRate));
388  waveForm.rate(tsRate);
389  waveForm.setWavelet(*(TFmap.pWavelet));
390  return *this;
391 }
392 
393 // copy 'from' injection stuff
395 {
396  wfSAVE = value.wfSAVE;
397  HRSS = value.HRSS;
398  ISNR = value.ISNR;
399  FREQ = value.FREQ;
400  BAND = value.BAND;
401  TIME = value.TIME;
402  TDUR = value.TDUR;
403 
404  IWFID = value.IWFID;
405 
406  for (int i=0;i<IWFP.size();i++) {
408  delete wf;
409  }
410  IWFP.clear();
411  for(int i=0;i<value.IWFP.size();i++) {
413  *wf = *value.IWFP[i];
414  IWFP.push_back(wf);
415  }
416 
417  return *this;
418 }
419 
420 // copy 'to' injection stuff
422 {
423  value.wfSAVE = wfSAVE;
424  value.HRSS = HRSS;
425  value.ISNR = ISNR;
426  value.FREQ = FREQ;
427  value.BAND = BAND;
428  value.TIME = TIME;
429  value.TDUR = TDUR;
430 
431  value.IWFID = IWFID;
432 
433  for (int i=0;i<value.IWFP.size();i++) {
435  delete wf;
436  }
437  value.IWFP.clear();
438  for(int i=0;i<IWFP.size();i++) {
440  *wf = *IWFP[i];
441  value.IWFP.push_back(wf);
442  }
443 
444  return *this;
445 }
446 
447 
448 //**************************************************************************
449 // initialize detector tenzor
450 //**************************************************************************
452 {
453  DT[0] = Ex[0]*Ex[0]-Ey[0]*Ey[0];
454  DT[1] = Ex[0]*Ex[1]-Ey[0]*Ey[1];
455  DT[2] = Ex[0]*Ex[2]-Ey[0]*Ey[2];
456 
457  DT[3] = DT[1];
458  DT[4] = Ex[1]*Ex[1]-Ey[1]*Ey[1];
459  DT[5] = Ex[1]*Ex[2]-Ey[1]*Ey[2];
460 
461  DT[6] = DT[2];
462  DT[7] = DT[5];
463  DT[8] = Ex[2]*Ex[2]-Ey[2]*Ey[2];
464 
465  if (strcmp(Name,"O1")==0) for (int i=0;i<9;i++) DT[i]*=2;
466  if (strcmp(Name,"N1")==0) for (int i=0;i<9;i++) DT[i]*=2;
467  if (strcmp(Name,"E1")==0) for (int i=0;i<9;i++) DT[i]*=2;
468 }
469 
470 //**************************************************************************
471 // return antenna pattern
472 //**************************************************************************
473 wavecomplex detector::antenna(double theta, double phi, double psi)
474 {
475  double a,b;
476 
477  theta *= PI/180.; phi *= PI/180.; psi *= PI/180.;
478 
479  double cT = cos(theta);
480  double sT = sin(theta);
481  double cP = cos(phi);
482  double sP = sin(phi);
483 
484  double d11 = DT[0];
485  double d12 = DT[1];
486  double d13 = DT[2];
487 
488  double d21 = DT[3];
489  double d22 = DT[4];
490  double d23 = DT[5];
491 
492  double d31 = DT[6];
493  double d32 = DT[7];
494  double d33 = DT[8];
495 
496 
497  double fp = 0.;
498  double fx = 0.;
499 
500  if(polarization==TENSOR) {
501 
502  fp = (cT*cP*d11 + cT*sP*d21 - sT*d31)*cT*cP
503  + (cT*cP*d12 + cT*sP*d22 - sT*d32)*cT*sP
504  - (cT*cP*d13 + cT*sP*d23 - sT*d33)*sT
505  + (cP*d21-sP*d11)*sP
506  - (cP*d22-sP*d12)*cP;
507 
508  fx = -(cT*cP*d11 + cT*sP*d21 - sT*d31)*sP
509  + (cT*cP*d12 + cT*sP*d22 - sT*d32)*cP
510  + (cP*d21-sP*d11)*cT*cP
511  + (cP*d22-sP*d12)*cT*sP
512  - (cP*d23-sP*d13)*sT;
513 
514  fp = -fp; // to follow convention in LIGO-T010110 and Anderson et al.
515 
516  if(fabs(psi)>0.) { // rotate in the waveframe A'=exp(-2*i*psi)A
517  a = fp*cos(2*psi)+fx*sin(2*psi);
518  b = -fp*sin(2*psi)+fx*cos(2*psi);
519  fp = a; fx = b;
520  }
521  }
522 
523  if(polarization==SCALAR) {
524 
525  fp = -(cT*cP*d11 + cT*sP*d21 - sT*d31)*cT*cP
526  - (cT*cP*d12 + cT*sP*d22 - sT*d32)*cT*sP
527  + (cT*cP*d13 + cT*sP*d23 - sT*d33)*sT
528  + (cP*d21-sP*d11)*sP
529  - (cP*d22-sP*d12)*cP;
530 
531  fp = 2*fp;
532  fx = 0;
533  }
534 
535  wavecomplex z(fp/2.,fx/2.);
536  return z;
537 }
538 
539 
540 
541 //**************************************************************************
542 // reconstruct wavelet series for a cluster, put it in waveForm
543 //**************************************************************************
544 double detector::getwave(int ID, netcluster& wc, char atype, size_t index)
545 {
546  int i,j,n,m,k,l;
547  double a,b,rms,fl,fh;
548  double R = this->TFmap.rate();
549  int L = int(this->TFmap.maxLayer()+1);
550 
551  waveForm.setWavelet(*(TFmap.pWavelet));
552  waveForm.rate(R);
553  rms = wc.getwave(ID,waveForm,atype,index);
554 
555  if(rms==0.) return rms;
556 
557 // create bandlimited detector output
558 
559  waveBand = waveForm; waveBand = 0.;
560  fh = waveBand.gethigh();
561  fl = waveBand.getlow();
562  l = int(this->TFmap.getLevel()) - int(waveBand.getLevel());
563 
564 // adjust waveBand resolution to match TFmap
565 
566  if(l<0) { waveBand.Inverse(-l); }
567  else { waveBand.Forward(l); }
568 
569  a = waveBand.start()*R;
570  i = int(a + ((a>0) ? 0.1 : -0.1));
571  k = waveBand.size(); j = 0;
572  if(i<0) { k += i; j = -i; i = 0; }
573  if((i/L)*L != i) cout<<"detector::getwave() time mismatch: "<<L<<" "<<i<<"\n";
574 
575  waveBand.cpf(this->TFmap,k,i,j);
576  waveNull = waveBand;
577 
578  n = int(2*L*fl/R+0.1)-1; // first layer
579  m = int(2*L*fh/R+0.1)+1; // last layer
580  if(n<=0) n=0;
581  if(m>=int(L)) m=L;
582  if(m<=n) { n=0; m=L; }
583 
584  WSeries<double> w = waveBand;
586 
587 // cout<<i<<" L="<<L<<" Band start="<<waveBand.start()<<endl;
588 
589  waveBand = 0;
590  for(k=n; k<m; k++) { w.getLayer(x,k); waveBand.putLayer(x,k); }
591 
592  waveBand.Inverse();
593  waveForm.Inverse();
594  waveNull.Inverse();
595  waveForm *= atype=='w' ? 1./sqrt(this->rate/R) : 1.; // rescale waveform
596  w = waveForm;
597 
598 // window the waveforms
599 
600  double sTARt= waveForm.start();
601  size_t I = waveForm.size();
602  size_t M = I/2;
603  double sum = waveForm.data[M]*waveForm.data[M];
604 
605  a = waveForm.rms();
606  double hrss = a*a*I;
607 
608  for(i=1; i<int(M); i++) {
609  a = waveForm.data[M-i];
610  b = waveForm.data[M+i];
611  sum += a*a+b*b;
612  if(sum/hrss > 0.999 && i/R>0.05) break;
613  }
614  n = i+int(0.05*R);
615  if(n < int(M-1)) i = size_t(n);
616  i = M - ((M-i)/L)*L; // sink with wavelet resolution.
617 
618 // cout<<"M="<<M<<" 2i="<<2*i<<" i/R"<<i/R<<endl;
619 
620  waveForm.cpf(w,2*i,M-i);
621  waveForm.resize(2*i);
622  waveForm.start(sTARt+(M-i)/R);
623 
624  w = waveBand;
625  waveBand.cpf(w,2*i,M-i);
626  waveBand.resize(2*i);
627  waveBand.start(sTARt+(M-i)/R);
628 
629  return rms;
630 }
631 
632 
633 //**************************************************************************
634 // set time delays
635 // time delay convention: t_detector-tau - arrival time at the center of Earth
636 //**************************************************************************
637 void detector::setTau(double sms,double t1,double t2,double p1,double p2)
638 {
639  size_t i;
640  skymap SM(sms,t1,t2,p1,p2);
641  size_t n = SM.size();
642  double x,y,z;
643 
644  for(i=0; i<n; i++) {
645  x = SM.getTheta(i)*PI/180.;
646  y = SM.getPhi(i)*PI/180.;
647  z = Rv[0]*sin(x)*cos(y) + Rv[1]*sin(x)*sin(y) + Rv[2]*cos(x);
648  SM.set(i,-z/speedlight);
649  }
650 
651  tau = SM;
652  return;
653 }
654 
655 
656 //**************************************************************************
657 // set time delays
658 // time delay convention: t_detector-tau - arrival time at the center of Earth
659 //**************************************************************************
660 void detector::setTau(int order)
661 {
662  size_t i;
663  skymap SM(order);
664  size_t n = SM.size();
665  double x,y,z;
666 
667  for(i=0; i<n; i++) {
668  x = SM.getTheta(i)*PI/180.;
669  y = SM.getPhi(i)*PI/180.;
670  z = Rv[0]*sin(x)*cos(y) + Rv[1]*sin(x)*sin(y) + Rv[2]*cos(x);
671  SM.set(i,-z/speedlight);
672  }
673 
674  tau = SM;
675  return;
676 }
677 
678 //**************************************************************************
679 // return detector time delay for specified source location
680 //**************************************************************************
681 double detector::getTau(double theta, double phi)
682 {
683  double x = theta*PI/180.;
684  double y = phi*PI/180.;
685  double z = Rv[0]*sin(x)*cos(y) + Rv[1]*sin(x)*sin(y) + Rv[2]*cos(x);
686  return -z/speedlight;
687 }
688 
689 //**************************************************************************
690 // set antenna patterns
691 //**************************************************************************
692 void detector::setFpFx(double sms,double t1,double t2,double p1,double p2)
693 {
694  size_t i;
695  skymap Sp(sms,t1,t2,p1,p2);
696  skymap Sx(sms,t1,t2,p1,p2);
697  size_t n = Sp.size();
698  double x,y;
699  wavecomplex a;
700 
701  for(i=0; i<n; i++) {
702  x = Sp.getTheta(i);
703  y = Sp.getPhi(i);
704  a = antenna(x,y);
705  Sp.set(i,a.real());
706  Sx.set(i,a.imag());
707  }
708 
709  mFp = Sp;
710  mFx = Sx;
711  return;
712 }
713 
714 //**************************************************************************
715 // set antenna patterns
716 //**************************************************************************
717 void detector::setFpFx(int order)
718 {
719  size_t i;
720  skymap Sp(order);
721  skymap Sx(order);
722  size_t n = Sp.size();
723  double x,y;
724  wavecomplex a;
725 
726  for(i=0; i<n; i++) {
727  x = Sp.getTheta(i);
728  y = Sp.getPhi(i);
729  a = antenna(x,y);
730  Sp.set(i,a.real());
731  Sx.set(i,a.imag());
732  }
733 
734  mFp = Sp;
735  mFx = Sx;
736  return;
737 }
738 
739 //: initialize delay filter
740 // delay index n: 0 1 2 3 4 5 6 ... M-3 M-2 M-1 M
741 // sample delay: 0 -1 -2 -3 -4 -5 -6 3 2 1 0
742 size_t detector::setFilter(size_t K, double phase, size_t upL)
743 {
744  if(TFmap.isWDM()) {
745  cout<<"wseries::setFilter(): not applicable to WDM TFmaps\n";
746  return 0;
747  }
748  size_t i,j,k,n,ii,jj;
749  size_t M = TFmap.maxLayer()+1; // number of wavelet layers
750  size_t L = TFmap.getLevel(); // wavelet decomposition depth
751  size_t m = M*TFmap.pWavelet->m_H; // buffer length
752  size_t N = M*(1<<upL); // number of time delays
753 
754 // K - total length of the delay filter
755  std::vector<delayFilter> F; // delay filter buffer
756  delayFilter v; v.index.resize(K); v.value.resize(K);
757  for(i=0; i<K; i++) v.value[i] = 0.;
758  for(i=0; i<M; i++) F.push_back(v);
759  slice S;
760  size_t s;
761  short inDex;
762  float vaLue;
763  int offst;
764 
765 
766 // set wavelet buffer
767 
768  Wavelet* pW = TFmap.pWavelet->Clone(); // wavelet used to calculate delay filter
769  Meyer<double> wM(1024,1); // up-sample wavelet
770  WSeries<double> w(*pW);
771  WSeries<double> W(wM); // up-sampled wavelet series
772 
773  cout<<"w.pWavelet->m_H "<<w.pWavelet->m_H<<" level "<<w.pWavelet->m_Level<<" ";
774 
775  w.resize(1024);
776  while(int(L)>w.getMaxLevel() || w.size()<m) w.resize(w.size()*2);
777  w.setLevel(L);
778  W.resize(w.size()*N/M);
779  W.setLevel(upL);
780 
781  cout<<"wsize: "<<w.size()<<endl;
782 
783  S = w.getSlice(0);
784  j = M*S.size()/2;
785 
786  double* pb = (double * )malloc(m*sizeof(double));
787  double** pp = (double **)malloc(m*sizeof(double*));
788  for(i=0; i<m; i++) pp[i] = w.data + i + (int(j) - int(m/2));
789  double* p0 = pp[0];
790  double* p;
791  double sum;
792  wavecomplex Z(cos(phase*PI/180.),sin(phase*PI/180.)); // phase shift
793  wavecomplex z;
794 
795  filter.clear();
796  this->nDFS = N; // store number of Delay Filter Samples
797  this->nDFL = M; // store number of Delay Filter Layers
798 
799  for(n=0; n<N; n++) { // loop over delays
800 
801  for(i=0; i<M; i++) { // loop over wavelet layers
802 
803  w = 0.;
804  S = w.getSlice(i);
805  p = w.data+S.start()+j;
806  s = S.start();
807  *p = 1.;
808  w.Inverse();
809 
810 // up-sample
811 
812  W = 0.;
813  W.putLayer(w,0);
814  W.Inverse();
815 
816 // phase shift
817  if(phase != 0.) {
818  W.FFTW(1);
819  for(k=2; k<W.size(); k+=2) {
820  z.set(W.data[k],W.data[k+1]); // complex amplitude
821  z *= Z;
822  W.data[k] = z.real();
823  W.data[k+1] = z.imag();
824  }
825  W.FFTW(-1);
826  }
827 
828 // time shift by integer number of samples
829 
830  W.cpf(W,W.size()-n,n);
831 
832 // down-sample
833 
834  W.Forward(upL);
835  W.getLayer(w,0);
836 
837 // get filter coefficients
838 
839  w.Forward(L);
840  if(n >= N/2) p -= M; // positive shift
841 
842  for(k=0; k<m; k++) {
843  *(pb+k) = *(p0+k); // save data in the buffer
844  *(p0+k) *= *(p0+k); // square
845  }
846 
847  w.waveSort(pp,0,m-1);
848 
849  for(k=m-1; k>=0; k--) {
850  offst = pp[k]-p;
851  if(abs(offst) >32767) continue;
852  inDex = short(offst);
853  vaLue = float(pb[pp[k]-p0]);
854  if(fabs(vaLue)<1.e-4) break;
855  offst += s;
856  offst -= (offst/M)*M;
857  if(offst<0) offst += M; // calculate offset
858  ii = pW->convertO2F(L,offst); // convert offset into frequency index
859 
860  if(offst != pW->getOffset(L,pW->convertF2L(L,ii))) cout<<"setFilter error 1\n";
861  offst -= inDex;
862  offst -= (offst/M)*M;
863  if(offst<0) offst += M; // calculate offset
864  if(offst != int(s)) cout<<"setFilter error 2: "<<offst<<" "<<s<<endl;
865 
866  jj = minDFindex(F[ii])-1; // index of least significant element in F[ii]
867 
868  for(size_t kk=0; kk<K; kk++)
869  if(fabs(F[ii].value[jj])>fabs(F[ii].value[kk])) cout<<"setFilter error 3:\n";
870 
871  if(jj>=K) {cout<<jj<<endl; continue;}
872  if(fabs(F[ii].value[jj]) < fabs(vaLue)){
873  F[ii].value[jj] = vaLue;
874  F[ii].index[jj] = -inDex;
875  }
876  }
877 
878  }
879 
880  for(i=0; i<M; i++) {
881  filter.push_back(F[i]);
882 
883 // if(n==3) {
884 // S = w.getSlice(i);
885 // printf("%3d %3d %1d %7.5f %7.5f %7.5f %7.5f %7.5f %7.5f %7.5f %7.5f %7.5f %7.5f\n",
886 // S.start(),i,n,F[i].value[0],F[i].value[1],F[i].value[2],F[i].value[3],F[i].value[4],
887 // F[i].value[5],F[i].value[6],F[i].value[7],F[i].value[8],F[i].value[9]);
888 // printf("%3d %3d %1d %7d %7d %7d %7d %7d %7d %7d %7d %7d %7d\n",
889 // S.start(),i,n,F[i].index[0],F[i].index[1],F[i].index[2],F[i].index[3],F[i].index[4],
890 // F[i].index[5],F[i].index[6],F[i].index[7],F[i].index[8],F[i].index[9]);
891 // }
892 
893  sum = 0;
894  v = filter[n*M+i];
895  for(k=0; k<K; k++) {
896  sum += F[i].value[k]*F[i].value[k];
897  if(n && F[i].value[k] == 0.) printf("%4d %4d %d4 \n",int(n),int(i),int(k));
898  if(v.value[k] != F[i].value[k] ||
899  v.index[k] != F[i].index[k]) cout<<"setFilter error 4\n";
900  F[i].value[k] = 0.;
901  }
902  if(sum<0.97) printf("%4d %4d %8.5f \n",int(n),int(i),sum);
903 
904  }
905 
906  }
907 
908  delete pW;
909  free(pp);
910  free(pb);
911  return filter.size();
912 }
913 
914 
915 //: initialize delay filter from another detector
917  size_t K = d.filter.size();
918  filter.clear();
919  std::vector<delayFilter>().swap(filter);
920  filter.reserve(K);
921 
922  for(size_t k=0; k<K; k++) {
923  filter.push_back(d.filter[k]);
924  }
925  return filter.size();
926 }
927 
928 
929 //: Dumps filters to file *fname in binary format.
930 void detector::writeFilter(const char *fname)
931 {
932  size_t i,j,k;
933  FILE *fp;
934 
935  if ( (fp=fopen(fname, "wb")) == NULL ) {
936  cout << " DumpBinary() error : cannot open file " << fname <<". \n";
937  exit(1);
938  }
939 
940  size_t M = size_t(TFmap.maxLayer()+1); // number of wavelet layers
941  size_t K = size_t(filter[0].index.size()); // delay filter length
942  size_t N = this->nDFS; // number of delays
943  size_t n = K * sizeof(float);
944  size_t m = K * sizeof(short);
945 
948 
949  fwrite(&K, sizeof(size_t), 1, fp); // write filter length
950  fwrite(&M, sizeof(size_t), 1, fp); // number of layers
951  fwrite(&N, sizeof(size_t), 1, fp); // number of delays
952 
953  for(i=0; i<N; i++) { // loop over delays
954  for(j=0; j<M; j++) { // loop over wavelet layers
955  for(k=0; k<K; k++) { // loop over filter coefficients
956  value.data[k] = filter[i*M+j].value[k];
957  index.data[k] = filter[i*M+j].index[k];
958  }
959  fwrite(value.data, n, 1, fp);
960  fwrite(index.data, m, 1, fp);
961  }
962  }
963  fclose(fp);
964 }
965 
966 //+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
967 //: Read filters from file *fname.
968 //+++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
969 void detector::readFilter(const char *fname)
970 {
971  size_t i,j,k;
972  FILE *fp;
973 
974  if ( (fp=fopen(fname, "rb")) == NULL ) {
975  cout << " DumpBinary() error : cannot open file " << fname <<". \n";
976  exit(1);
977  }
978 
979  size_t M; // number of wavelet layers
980  size_t K; // delay filter length
981  size_t N; // number of delays
982 
983  fread(&K, sizeof(size_t), 1, fp); // read filter length
984  fread(&M, sizeof(size_t), 1, fp); // read number of layers
985  fread(&N, sizeof(size_t), 1, fp); // read number of delays
986 
987  size_t n = K * sizeof(float);
988  size_t m = K * sizeof(short);
991  delayFilter v;
992 
993  v.value.clear(); v.value.reserve(K);
994  v.index.clear(); v.index.reserve(K);
995  this->clearFilter(); filter.reserve(N*M);
996  this->nDFS = N; // set number of delay samples
997  this->nDFL = M; // set number of delay layers
998 
999  for(k=0; k<K; k++) { // loop over filter coefficients
1000  v.value.push_back(0.);
1001  v.index.push_back(0);
1002  }
1003 
1004  for(i=0; i<N; i++) { // loop over delays
1005  for(j=0; j<M; j++) { // loop over wavelet layers
1006  fread(value.data, n, 1, fp);
1007  fread(index.data, m, 1, fp);
1008  for(k=0; k<K; k++) { // loop over filter coefficients
1009  v.value[k] = value.data[k];
1010  v.index[k] = index.data[k];
1011  }
1012 
1013 // if(i<3){
1014 // printf("%6.5f %6.5f %6.5f %6.5f %6.5f %6.5f %6.5f %6.5f %6.5f %6.5f \n",
1015 // v.value[0],v.value[1],v.value[2],v.value[3],v.value[4],
1016 // v.value[5],v.value[6],v.value[7],v.value[8],v.value[9]);
1017 // printf("%4d %4d %4d %4d %4d %4d %4d %4d %4d %4d \n",
1018 // v.index[0],v.index[1],v.index[2],v.index[3],v.index[4],
1019 // v.index[5],v.index[6],v.index[7],v.index[8],v.index[9]);
1020 // }
1021 
1022  filter.push_back(v);
1023  }
1024  }
1025  fclose(fp);
1026 }
1027 
1028 
1029 //: apply delay filter to input WSeries and put result in TFmap
1031 {
1032  int i,j,jb,je;
1033 
1034  register int k;
1035  register double* p;
1036  register double* q;
1037 
1038  if(TFmap.isWDM()) {
1039  cout<<"wseries::delay(): not applicable to WDM TFmaps\n";
1040  return;
1041  }
1042 
1043  slice S;
1044  delayFilter v = filter[0]; // delay filter
1045 
1046  int M = this->nDFL; // number of wavelet layers
1047  int N = this->nDFS; // number of delay samples per pixel
1048  int K = int(v.index.size()); // delay filter length
1049  int n = int(t*TFmap.wavearray<double>::rate());
1050  int m = n>0 ? (n+N/2-1)/N : (n-N/2)/N; // delay in wavelet pixels
1051  int l = TFmap.pWavelet->m_H/4+2;
1052 
1053  n = n - m*N; // n - delay in samples
1054  if(n <= 0) n = -n; // filter index for negative delays
1055  else n = N-n; // filter index for positive delays
1056 
1057  cout<<"delay="<<t<<" m="<<m<<" n="<<n<<" M="<<M<<" K="<<K<<" N="<<N<<endl;
1058 
1059  int mM = m*M;
1060  double* A = (double*)malloc(K*sizeof(double));
1061  int* I = (int*)malloc(K*sizeof(int));
1062 
1063  for(i=0; i<M; i++) {
1064 
1065  S = w.getSlice(i);
1066  v = filter[n*M+i];
1067  jb = m>0 ? S.start()+(l+m)*M : S.start()+l*M;
1068  je = m<0 ? w.size() +(m-l)*M : w.size() -l*M;
1069 
1070  for(k=0; k<K; k++) { A[k]=double(v.value[k]); I[k]=int(v.index[k]); }
1071 
1072  if(K==16){
1073  for(j=jb; j<je; j+=M) {
1074  p = w.data+j-mM;
1075  TFmap.data[j] = A[0]*p[I[0]] + A[1]*p[I[1]] + A[2]*p[I[2]] + A[3]*p[I[3]]
1076  + A[4]*p[I[4]] + A[5]*p[I[5]] + A[6]*p[I[6]] + A[7]*p[I[7]]
1077  + A[8]*p[I[8]] + A[9]*p[I[9]] + A[10]*p[I[10]] + A[11]*p[I[11]]
1078  + A[12]*p[I[12]] + A[13]*p[I[13]] + A[14]*p[I[14]] + A[15]*p[I[15]];
1079  }
1080  }
1081 
1082  else if(K==32){
1083  for(j=jb; j<je; j+=M) {
1084  p = w.data+j-mM;
1085  TFmap.data[j] = A[0]*p[I[0]] + A[1]*p[I[1]] + A[2]*p[I[2]] + A[3]*p[I[3]]
1086  + A[4]*p[I[4]] + A[5]*p[I[5]] + A[6]*p[I[6]] + A[7]*p[I[7]]
1087  + A[8]*p[I[8]] + A[9]*p[I[9]] + A[10]*p[I[10]] + A[11]*p[I[11]]
1088  + A[12]*p[I[12]] + A[13]*p[I[13]] + A[14]*p[I[14]] + A[15]*p[I[15]]
1089  + A[16]*p[I[16]] + A[17]*p[I[17]] + A[18]*p[I[18]] + A[19]*p[I[19]]
1090  + A[20]*p[I[20]] + A[21]*p[I[21]] + A[22]*p[I[22]] + A[23]*p[I[23]]
1091  + A[24]*p[I[24]] + A[25]*p[I[25]] + A[26]*p[I[26]] + A[27]*p[I[27]]
1092  + A[28]*p[I[28]] + A[29]*p[I[29]] + A[30]*p[I[30]] + A[31]*p[I[31]];
1093  }
1094  }
1095 
1096  else {
1097  for(j=jb; j<je; j+=M) {
1098  p = w.data+j-mM;
1099  q = TFmap.data+j;
1100  k = K; *q = 0.;
1101  while(k-- > 0) { *q += *(A++) * p[*(I++)]; }
1102  A -= K; I -= K;
1103  }
1104  }
1105 
1106  }
1107  free(A);
1108  free(I);
1109 }
1110 
1111 //: return noise variance for selected wavelet layer or pixel
1112 double detector::getNoise(size_t I, int J)
1113 {
1114  if(!nRMS.size()) return 0.;
1115  if(int(I)>TFmap.maxLayer()) return 0.;
1116 
1117  size_t l,j;
1118  double x, g;
1119  size_t N = nRMS.maxLayer()+1; // number of layers in nRMS
1120  size_t M = TFmap.maxLayer()+1; // number of layers in TFmap
1121  size_t n = I*N/M; // layer low index in nRMS
1122  size_t m = (I+1)*(N/M); // layer high index in nRMS
1123  double t = TFmap.start();
1124  double T = nRMS.start();
1125  slice S = TFmap.getSlice(I);
1126  double rATe = TFmap.wrate(); // map rate
1127 
1128  double rESn = nRMS.frequency(1)-nRMS.frequency(0); // noise F resolution
1129 
1130  // printf("%4d %4d %16.3f %16.3f\n",M,N,t,T);
1131 
1132  if(M>N || t>T) {
1133  cout<<"detector::getNoise(): invalid noise rms array nRMS\n";
1134  return 0.;
1135  }
1136 
1137  size_t L = size_t(TFmap.getlow()/rESn+0.1); // low F boundary
1138  size_t H = size_t(TFmap.gethigh()/rESn+0.1); // high F boundary
1139  if(n>=H || m<=L) return 0.; // out of boundaries
1140 
1141  // printf("%4d %4d\n",n,m);
1142 
1143  if(n >= m) {
1144  cout<<"detector::getNoise():: invalid noise rms array nRMS\n";
1145  return 0.;
1146  }
1147 
1148  S = nRMS.getSlice(n);
1149  size_t K = S.size(); // size of noise vector in each layer
1150 
1151  if(!K) return 0.;
1152 
1153  if(J < 0) { // get noise rms for specified layer
1154  wavearray<double> rms(K);
1155  rms = 0.;
1156 
1157  for(l=n; l<m; l++) {
1158  S = nRMS.getSlice(l);
1159  g = (n<L || m>H) ? 1.e10 : 1.; // supress layers below HP cut-off frequency
1160  for(j=0; j<rms.size(); j++) {
1161  x = nRMS.data[S.start()+j*S.stride()]*g;
1162  rms.data[j] += 1./x/x;
1163  }
1164  }
1165 
1166  for(j=0; j<rms.size(); j++) {
1167  rms.data[j] = sqrt((double(m)-double(n))/rms.data[j]);
1168  }
1169 
1170  return rms.mean();
1171  }
1172 
1173  else { // get noise rms for specified pixel
1174 
1175  double RMS = 0;
1176 
1177  t += J/rATe; // pixel GPS time
1178 
1179  int inRMS = int((t-nRMS.start())*nRMS.rate());
1180  int inVAR = nVAR.size() ? int((t-nVAR.start())*nVAR.rate()) : 0;
1181 
1182  if(inRMS < 0) inRMS = 0;
1183  else if(size_t(inRMS)>=K && K) inRMS = K-1;
1184 
1185  if(inVAR <= 0) inVAR = 0;
1186  else if(size_t(inVAR)>=nVAR.size()) inVAR = nVAR.size()-1;
1187 
1188  for(l=n; l<m; l++) { // get noise vector for specified fl-fh
1189  S = nRMS.getSlice(l);
1190  g = (n<L || m>H) ? 1.e10 : 1.; // supress layers below HP cut-off
1191  x = nRMS.data[S.start()+inRMS*S.stride()]*g;
1192  RMS += 1./x/x;
1193  }
1194 
1195  RMS /= double(m)-double(n);
1196  RMS = sqrt(1./RMS);
1197 
1198  if(!nVAR.size() || rESn*n<nVAR.getlow() || rESn*m>nVAR.gethigh()) return RMS;
1199 
1200  return RMS*double(nVAR.data[inVAR]);
1201  }
1202 }
1203 
1204 
1205 //**************************************************************************
1206 // set noise rms in pixel data array in cluster structure
1207 //**************************************************************************
1209 {
1210  size_t i,j,n,m;
1211  size_t M = wc->size();
1212 
1213  if(!M) return false;
1214 
1215  if(!nRMS.size()) return false;
1216 
1217  size_t max_layer = nRMS.maxLayer();
1218  netpixel* p = NULL; // pointer to pixel structure
1219  slice S;
1220 
1221  int k;
1222  int K = nRMS.size()/(max_layer+1); // number of RMS measurements per layer
1223  double To = nRMS.start();
1224  double Ro = nRMS.wrate();
1225  double Fo = nRMS.frequency(0); // central frequency of zero layer
1226  double dF = nRMS.frequency(1)-Fo; // nRMS frequency resolution
1227  double fl = wc->getlow()-0.1;
1228  double fh = wc->gethigh()+0.1;
1229 
1230  double x,f,t,r,F,g;
1231  Fo = Fo==0. ? 0.5 : 0.; // WDM : wavelet frequency correction
1232 
1233  for(i=0; i<M; i++){
1234  p = wc->getPixel(0,i);
1235 
1236  if(p->frequency > max_layer ||
1237  int(p->rate/Ro+0.01) < 1 ||
1238  p->frequency == 0) { // never use zero layer
1239  cout<<"detector::setrms() - illegal pixel from zero level\n";
1240  exit(0);
1241  }
1242 
1243  x = p->frequency-Fo; // fractional frequency index for wavelet and WDM
1244  f = x*p->rate/2.;
1245  n = size_t(f/dF+0.6); // first layer in nRMS
1246  F = (x+1)*p->rate/2.;
1247  m = size_t(F/dF+0.6); // last layer in nRMS
1248  if(m>max_layer) m=max_layer+1;
1249  t = p->getdata('I',I)/p->rate/p->layers; // takes into account time lag
1250  t+= wc->start; // gps time
1251  k = int((t-To)*Ro); // time index in the noise array
1252 
1253  if(k>=K) k -= k ? 1 : 0;
1254  if(k<0 || n>=m || k>=K) {
1255  cout<<"detector::setrms() - invalid input: ";
1256  cout<<k<<" "<<n<<" "<<m<<" "<<f<<" "<<F<<" "<<t<<endl;
1257  cout<<p->frequency<<" "<<p->rate/2.<<" "<<dF<<" "<<fl<<endl;
1258  exit(0);
1259  }
1260 
1261 // get noise rms for specified pixel
1262 
1263  r = 0.;
1264  for(j=n; j<m; j++) {
1265  S = nRMS.getSlice(j);
1266  g = (f<fl || F>fh) ? 1.e10 : 1.; // supress layers below HP cut-off
1267  x = nRMS.data[S.start()+k*S.stride()];
1268  r += 1./x/x;
1269  }
1270  r = sqrt((m-n)/r);
1271  if(nVAR.size() && f<nVAR.gethigh() && F>nVAR.getlow()) r *= nVAR.get(t);
1272  p->setdata(r,'N',I);
1273 
1274  }
1275  return true;
1276 }
1277 
1278 //**************************************************************************
1279 // apply band pass filter with cut-offs specified by parameters (used by 1G)
1280 //**************************************************************************
1281 void detector::bandPass1G(double f1, double f2)
1282 {
1283  int i;
1284  double dF = TFmap.frequency(1)-TFmap.frequency(0); // frequency resolution
1285  double fl = fabs(f1)>0. ? fabs(f1) : this->TFmap.getlow();
1286  double fh = fabs(f2)>0. ? fabs(f2) : this->TFmap.gethigh();
1287  size_t n = TFmap.pWavelet->m_WaveType==WDMT ? size_t((fl+dF/2.)/dF+0.1) : size_t(fl/dF+0.1);
1288  size_t m = TFmap.pWavelet->m_WaveType==WDMT ? size_t((fh+dF/2.)/dF+0.1)-1 : size_t(fh/dF+0.1)-1;
1289  size_t M = this->TFmap.maxLayer()+1;
1291 
1292  if(n>m) return;
1293 
1294  for(i=0; i<int(M); i++) { // ......f1......f2......
1295 
1296  if((f1>=0 && i>=n) && (f2>=0 && i<=m)) continue; // zzzzzz..........zzzzzz band pass
1297  if((f1<0 && i<n) || (f2<0 && i>m)) continue; // ......zzzzzzzzzz...... band cut
1298  if((f1<0 && f2>=0 && i<n)) continue; // ......zzzzzzzzzzzzzzzz low pass
1299  if((f1>=0 && f2<0 && i>=m)) continue; // zzzzzzzzzzzzzzzz...... high pass
1300 
1301  this->TFmap.getLayer(w,i+0.01); w=0.; this->TFmap.putLayer(w,i+0.01);
1302  this->TFmap.getLayer(w,-i-0.01); w=0.; this->TFmap.putLayer(w,-i-0.01);
1303  }
1304  return;
1305 }
1306 
1307 //**************************************************************************
1308 // apply band pass filter with cut-offs specified by parameters
1309 //**************************************************************************
1310 /*
1311 void detector::bandPass(double f1, double f2, double a)
1312 {
1313 // assign constant value a to wseries layer coefficients
1314 // 0utside of the band defined by frequencies f1 and f2
1315 // f1>0, f2>0 - zzzzzz..........zzzzzz band pass
1316 // f1<0, f2<0 - ......zzzzzzzzzz...... band cut
1317 // f1<0, f2>0 - ......zzzzzzzzzzzzzzzz low pass
1318 // f1>0, f2<0 - zzzzzzzzzzzzzzzz...... high pass
1319  int i;
1320  double dF = TFmap.frequency(1)-TFmap.frequency(0); // frequency resolution
1321  double fl = fabs(f1)>0. ? fabs(f1) : this->TFmap.getlow();
1322  double fh = fabs(f2)>0. ? fabs(f2) : this->TFmap.gethigh();
1323  size_t n = TFmap.pWavelet->m_WaveType==WDMT ? size_t((fl+dF/2.)/dF+0.1) : size_t(fl/dF+0.1);
1324  size_t m = TFmap.pWavelet->m_WaveType==WDMT ? size_t((fh+dF/2.)/dF+0.1)-1 : size_t(fh/dF+0.1)-1;
1325  size_t M = this->TFmap.maxLayer()+1;
1326  wavearray<double> w;
1327 
1328  if(n>m) return;
1329 
1330  for(i=0; i<int(M); i++) { // ......f1......f2......
1331 
1332  if((f1>=0 && i>=n) && (f2>=0 && i<=m)) continue; // zzzzzz..........zzzzzz band pass
1333  if((f1<0 && i<n) || (f2<0 && i>m)) continue; // ......zzzzzzzzzz...... band cut
1334  if((f1<0 && f2>=0 && i<n)) continue; // ......zzzzzzzzzzzzzzzz low pass
1335  if((f1>=0 && f2<0 && i>=m)) continue; // zzzzzzzzzzzzzzzz...... high pass
1336 
1337  this->TFmap.getLayer(w,i+0.01); w=a; this->TFmap.putLayer(w,i+0.01);
1338  this->TFmap.getLayer(w,-i-0.01); w=a; this->TFmap.putLayer(w,-i-0.01);
1339  }
1340  return;
1341 }
1342 */
1343 //**************************************************************************
1344 // calculate hrss of injected responses
1345 // returns number of eligible injections
1346 // update MDC series with whitened time series
1347 //**************************************************************************
1348 size_t detector::setsim(WSeries<double> &wi, std::vector<double>* pT, double dT, double offset, bool saveWF)
1349 {
1350  int j,nstop,nstrt,n,m,J;
1351  size_t i,k;
1352  double a,b,T,E,D,H,f;
1353  double R = this->rate; // original data rate
1354  size_t K = pT->size();
1355  size_t I = wi.maxLayer()+1;
1356  size_t M = 0;
1357  bool pOUT = dT>0. ? false : true; // printout flag
1358 
1359  dT = fabs(dT);
1360 
1361  if(wi.size() != TFmap.size()) {
1362  cout<<"setsim(): mismatch between MDC size "
1363  <<wi.size()<<" and data size "<<TFmap.size()<<"\n";
1364  exit(1);
1365  }
1366 
1367  if(!K) return K;
1368  if(!this->nRMS.size() || !this->TFmap.size()) return 0;
1369 
1370  if(this->HRSS.size() != K) {
1371  this->HRSS.resize(K);
1372  this->ISNR.resize(K);
1373  this->FREQ.resize(K);
1374  this->BAND.resize(K);
1375  this->TIME.resize(K);
1376  this->TDUR.resize(K);
1377  }
1378  this->HRSS = 0.;
1379  this->ISNR = 0.;
1380  this->FREQ = 0.;
1381  this->BAND = 0.;
1382  this->TIME = 0.;
1383  this->TDUR = 0.;
1384 
1385  WSeries<double> W; W=wi;
1387  WSeries<double> hot; hot=wi; hot.Inverse();
1389  std::slice S;
1390 
1391 // whiten injection if noise array is filled in
1392 
1393  if(W.pWavelet->m_WaveType==WDMT) { // WDM type
1394  if(!W.white(this->nRMS,1)) { // whiten 0 phase WSeries
1395  cout<<"detector::setsim error: invalid noise array\n"; exit(1);
1396  }
1397  if(!W.white(this->nRMS,-1)) { // whiten 90 phase WSeries
1398  cout<<"detector::setsim error: invalid noise array\n"; exit(1);
1399  }
1400  w = W; // whitened WS
1401  WSeries<double> wtmp = w;
1402  w.Inverse();
1403  wtmp.Inverse(-2);
1404  w += wtmp;
1405  w *= 0.5;
1406  } else { // wavelet type
1407  if(!W.white(this->nRMS)) {
1408  cout<<"detector::setsim error: invalid noise array\n"; exit(1);
1409  }
1410  w = W; // whitened WS
1411  w.Inverse(); // whitened TS
1412  }
1413 
1414 
1415 // isolate injections in time series w
1416 
1417  size_t N = w.size(); // MDC size
1418  double rate = w.wavearray<double>::rate(); // simulation rate
1419  double bgps = w.start()+offset+1.; // analysed start time
1420  double sgps = w.start()-offset+N/rate-1.; // analysed stop time
1421 
1422  for(k=0; k<K; k++) {
1423 
1424  T = (*pT)[k] - w.start();
1425 
1426  nstrt = int((T - dT)*rate);
1427  nstop = int((T + dT)*rate);
1428  if(nstrt<=0) nstrt = 0;
1429  if(nstop>=int(N)) nstop = N;
1430  if(nstop<=0) continue; // outside of the segment
1431  if(nstrt>=int(N)) continue; // outside of the segment
1432 
1433  E = T = 0.;
1434  for(j=nstrt; j<nstop; j++) {
1435  a = w.data[j];
1436  T += a*a*j/rate; // central time
1437  E += a*a; // SNR
1438  }
1439  T /= E; // central time for k-th injection
1440 
1441 // zoom in
1442 
1443  nstrt = int((T - dT)*rate);
1444  nstop = int((T + dT)*rate);
1445  if(nstrt<=0) nstrt = 0;
1446  if(nstop>=int(N)) nstop = N;
1447  if(nstop<=0) continue; // outside of the segment
1448  if(nstrt>=int(N)) continue; // outside of the segment
1449 
1450  E = T = D = H = 0.;
1451  for(j=nstrt; j<nstop; j++) {
1452  a = w.data[j];
1453  T += a*a*j/rate; // central time
1454  E += a*a; // SNR
1455  }
1456  T /= E;
1457 
1458  m = int((T - dT)*W.wrate());
1459  n = int((T + dT)*W.wrate());
1460  S = W.getSlice(0); // zero layer
1461  if(m<=0) m = 0; // check left
1462  if(n>=int(S.size())) n = S.size()-1; // check right
1463 
1464  for(j=nstrt; j<nstop; j++) {
1465  a = w.data[j]*(j/rate-T);
1466  D += a*a; // duration
1467  a = hot.data[j];
1468  H += a*a; // hrss
1469  }
1470 
1471  D = sqrt(D/E);
1472  H = sqrt(H/R);
1473  T += w.start();
1474  if(T<bgps || T>sgps) continue; // outside of the segment
1475 
1476  this->ISNR.data[k] = E;
1477  this->TIME.data[k] = T;
1478  this->TDUR.data[k] = D;
1479  this->HRSS.data[k] = H;
1480 
1481  E = 0.;
1482  for(i=0; i<I; i++) {
1483  S = W.getSlice(i+0.001);
1484  J = S.stride();
1485  a = W.rms(slice(S.start()+m*J,n-m,J)); // get rms
1486  f = W.frequency(i); // layer frequency
1487  this->FREQ.data[k] += f*a*a; // central frequency
1488  this->BAND.data[k] += f*f*a*a; // bandwidth
1489  E += a*a;
1490  }
1491 
1492  this->FREQ.data[k] /= E;
1493  a = this->FREQ.data[k];
1494  b = this->BAND.data[k]/E - a*a;
1495  this->BAND.data[k] = b>0. ? sqrt(b) : 0.;
1496 
1497 // save waveform
1498 
1499  double time = this->TIME.data[k];
1500  double tdur = this->TDUR.data[k];
1501  double tDur = 6*tdur;
1502  if (tDur > dT) tDur = dT;
1503  int nStrt = int((time-tDur-w.start())*rate);
1504  int nStop = int((time+tDur-w.start())*rate);
1505  if(nStrt<0) nStrt = 0;
1506  if(nStop>int(N)) nStop = N;
1507 
1508 // window the injected waveforms
1509 
1510  size_t I = int(2.*dT*rate);
1511  int OS = int((time - dT - w.start())*rate);
1512  double ms = 0;
1513 
1514  for (j=0;j<I;j++) ms += ((OS+j>=0)&&(OS+j<N)) ? w.data[OS+j]*w.data[OS+j] : 0.;
1515 
1516  OS += I/2;
1517  double a,b;
1518  double sum = ((OS>=0)&&(OS<N)) ? w.data[OS]*w.data[OS] : 0.;
1519  for(j=1; j<int(I/2); j++) {
1520  a = ((OS-j>=0)&&(OS-j<N)) ? w.data[OS-j] : 0.;
1521  b = ((OS+j>=0)&&(OS+j<N)) ? w.data[OS+j] : 0.;
1522  sum += a*a+b*b;
1523  if(sum/ms > 0.999) break;
1524  }
1525 
1526  nStrt = int((time-w.start())*rate)-j;
1527  nStop = int((time-w.start())*rate)+j;
1528  if(nStrt<0) nStrt = 0;
1529  if(nStop>int(N)) nStop = N;
1530 
1532  wf->rate(rate);
1533  wf->start(w.start()+nStrt/rate);
1534  wf->resize(nStop-nStrt);
1535  for(j=nStrt; j<nStop; j++) {
1536  wf->data[j-nStrt] = w.data[j];
1537  }
1538 
1539 // apply freq cuts
1540 
1541  double f_low = this->TFmap.getlow();
1542  double f_high = this->TFmap.gethigh();
1543  //cout << "f_low : " << f_low << " f_high : " << f_high << endl;
1544  wf->FFTW(1);
1545  double Fs=((double)wf->rate()/(double)wf->size())/2.;
1546  for (j=0;j<wf->size()/2;j+=2) {
1547  double f=j*Fs;
1548  if((f<f_low)||(f>f_high)) {wf->data[j]=0.;wf->data[j+1]=0.;}
1549  }
1550  wf->FFTW(-1);
1551 
1552 // compute SNR,TIME,DUR within the search frequency band
1553 
1554  E = T = D = 0.;
1555  for(j=0;j<wf->size();j++) {
1556  a = wf->data[j];
1557  T += a*a*j/rate; // central time
1558  E += a*a; // SNR
1559  }
1560  T /= E;
1561 
1562  for(j=0;j<wf->size();j++) {
1563  a = wf->data[j]*(j/rate-T);
1564  D += a*a; // duration
1565  }
1566 
1567  D = sqrt(D/E);
1568  T += wf->start();
1569 
1570  this->ISNR.data[k] = E;
1571  this->TIME.data[k] = T;
1572  this->TDUR.data[k] = D;
1573 
1574  if (saveWF) {
1575  wf->resetFFTW();
1576  IWFID.push_back(k);
1577  IWFP.push_back(wf);
1578  } else {
1579  delete wf;
1580  }
1581 
1582  // save strain waveform
1583  if (saveWF) {
1584 
1585  // compute central time
1586  T = E = 0.;
1587  for(j=nstrt; j<nstop; j++) {
1588  a = hot.data[j];
1589  T += a*a*j/rate; // central time
1590  E += a*a; // energy
1591  }
1592  T /= E;
1593 
1594  // compute the time range containing the 0.999 of the total energy
1595  I = int(2.*dT*rate);
1596  OS = int((T - dT)*rate);
1597  ms = 0;
1598 
1599  for (j=0;j<I;j++) ms += ((OS+j>=0)&&(OS+j<N)) ? hot.data[OS+j]*hot.data[OS+j] : 0.;
1600 
1601  OS += I/2;
1602  sum = ((OS>=0)&&(OS<N)) ? hot.data[OS]*hot.data[OS] : 0.;
1603  for(j=1; j<int(I/2); j++) {
1604  a = ((OS-j>=0)&&(OS-j<N)) ? hot.data[OS-j] : 0.;
1605  b = ((OS+j>=0)&&(OS+j<N)) ? hot.data[OS+j] : 0.;
1606  sum += a*a+b*b;
1607  if(sum/ms > 0.999) break;
1608  }
1609 
1610  nStrt = int(T*rate)-j;
1611  nStop = int(T*rate)+j;
1612  if(nStrt<0) nStrt = 0;
1613  if(nStop>int(N)) nStop = N;
1614 
1615  // select strain mdc data
1617  WF->rate(rate);
1618  WF->start(hot.start()+nStrt/rate);
1619  WF->resize(nStop-nStrt);
1620  for(j=nStrt; j<nStop; j++) WF->data[j-nStrt] = hot.data[j];
1621 
1622  // store strain data (use ID=-(k+1))
1623  IWFID.push_back(-(k+1));
1624  IWFP.push_back(WF);
1625  }
1626 
1627  if(pOUT)
1628  printf("%3d T+-dT: %8.3f +-%5.3f, F+-dF: %4.0f +-%4.0f, SNR: %6.1e, hrss: %6.1e\n",
1629  int(M), T-bgps, D, FREQ.data[k], BAND.data[k], E, H);
1630 
1631 // (*pT)[k] = T;
1632  M++;
1633  }
1634  wi = w;
1635  return M;
1636 }
1637 
1638 //**************************************************************************
1639 // modify input signals (wi) at times pT according the factor pF
1640 //**************************************************************************
1641 size_t detector::setsnr(wavearray<double> &wi, std::vector<double>* pT, std::vector<double>* pF, double dT, double offset)
1642 {
1643  int j,nstop,nstrt;
1644  size_t k;
1645  double F,T;
1646  size_t K = pT->size();
1647  size_t N = wi.size();
1648  size_t M = 0;
1649  double rate = wi.rate(); // simulation rate
1650 
1651  dT = fabs(dT);
1652 
1653  wavearray<double> w; w=wi;
1654 // isolate injections
1655 
1656  for(k=0; k<K; k++) {
1657 
1658  F = (*pF)[k];
1659  T = (*pT)[k] - w.start();
1660 
1661  nstrt = int((T - dT)*rate);
1662  nstop = int((T + dT)*rate);
1663  if(nstrt<=0) nstrt = 0;
1664  if(nstop>=int(N)) nstop = N;
1665  if(nstop<=0) continue; // outside of the segment
1666  if(nstrt>=int(N)) continue; // outside of the segment
1667 
1668  for(j=nstrt; j<nstop; j++) {
1669  w.data[j]*=F;
1670  }
1671 
1672  M++;
1673  }
1674  wi = w;
1675  return M;
1676 }
1677 
1678 
1679 //**************************************************************************
1680 // apply sample shift to time series in TFmap
1681 //**************************************************************************
1682 void detector::delay(double theta, double phi)
1683 {
1684  if(!TFmap.size()) return;
1685  double R = this->TFmap.wavearray<double>::rate();
1686  double T = this->getTau(theta,phi); // time delay: +/- increase/decrease gps
1687  size_t n = size_t(fabs(T)*R); // shift in samples
1688  size_t m = this->TFmap.size();
1690  w = this->TFmap;
1691  TFmap = 0.;
1692 
1693  if(T<0) TFmap.cpf(w,m-n,0,n); // shift forward
1694  else TFmap.cpf(w,m-n,n,0); // shift backward
1695  return;
1696 }
1697 
1698 //**************************************************************************
1699 // apply sample shift to input time series
1700 //**************************************************************************
1702 {
1703  if(!x.size()) return;
1704  double R = x.rate();
1705  double T = this->getTau(theta,phi); // time delay: +/- increase/decrease gps
1706  size_t n = size_t(fabs(T)*R); // shift in samples
1707  size_t m = x.size();
1709  w = x;
1710  x = 0.;
1711 
1712  if(T<0) x.cpf(w,m-n,0,n); // shift forward
1713  else x.cpf(w,m-n,n,0); // shift backward
1714  return;
1715 }
1716 
1717 //**************************************************************************
1718 // apply time shift T to input time series
1719 //**************************************************************************
1721 {
1722  if(!x.size()) return;
1723  double R = x.rate();
1724  size_t n = size_t(fabs(T)*R+0.5); // shift in samples
1725  size_t m = x.size();
1727  w = x;
1728  x = 0.;
1729 
1730  if(T<0) x.cpf(w,m-n,0,n); // shift forward
1731  else x.cpf(w,m-n,n,0); // shift backward
1732  return;
1733 }
1734 
1735 double detector::getWFtime(char atype)
1736 {
1737 // returns central time of reconstructed waveform
1738  double e;
1739  double T = 0.;
1740  double E = 0.;
1741  WSeries<double>* wf = atype=='S' ? &this->waveForm : &this->waveBand;
1742  for(size_t i=0; i<wf->size(); i++) {
1743  e = wf->data[i]*wf->data[i];
1744  T += e*i;
1745  E += e;
1746  }
1747  return E>0. ? wf->start()+T/E/wf->rate() : 0.;
1748 }
1749 
1750 double detector::getWFfreq(char atype)
1751 {
1752 // returns central frequency of reconstructed waveform
1753  double e;
1754  double F = 0.;
1755  double E = 0.;
1756  WSeries<double>* wf = atype=='S' ? &this->waveForm : &this->waveBand;
1757  wf->FFTW(1);
1758  for(size_t i=0; i<wf->size(); i+=2) {
1759  e = wf->data[i]*wf->data[i];
1760  e += wf->data[i+1]*wf->data[i+1];
1761  F += e*i;
1762  E += e;
1763  }
1764  return E>0. ? 0.5*F*wf->rate()/E/wf->size() : 0.;
1765 }
1766 
1767 //______________________________________________________________________________
1769 {
1770  detectorParams xdP = getDetectorParams();
1771 
1772  char LAT;
1773  double theta_t=xdP.latitude;
1774  if(theta_t>0) LAT='N'; else {LAT='S';theta_t=-theta_t;}
1775  int theta_d = int(theta_t);
1776  int theta_m = int((theta_t-theta_d)*60);
1777  float theta_s = (theta_t-theta_d-theta_m/60.)*3600.;
1778 
1779  char LON;
1780  double phi_t=xdP.longitude;
1781  if(phi_t>0) LON='E'; else {LON='W';phi_t=-phi_t;}
1782  int phi_d = int(phi_t);
1783  int phi_m = int((phi_t-phi_d)*60);
1784  float phi_s = (phi_t-phi_d-phi_m/60.)*3600.;
1785 
1786  cout << endl;
1787  cout << "----------------------------------------------" << endl;
1788  cout << "IFO : " << xdP.name << " (Geographic Coordinates) " << endl;
1789  cout << "----------------------------------------------" << endl;
1790  cout << endl;
1791  cout << "latitude : " << xdP.latitude << " longitude : " << xdP.longitude << endl;
1792  cout << endl;
1793  cout << "LAT : " << LAT << " " << theta_d << ", " << theta_m << ", " << theta_s << endl;
1794  cout << "LON : " << LON << " " << phi_d << ", " << phi_m << ", " << phi_s << endl;
1795  cout << endl;
1796  int precision=cout.precision(12);
1797  // radius vector to beam splitter
1798  cout << "Rv : " << Rv[0] << " " << Rv[1] << " " << Rv[2] << " " << endl;
1799  // vector along x-arm
1800  cout << "Ex : " << Ex[0] << " " << Ex[1] << " " << Ex[2] << " " << endl;
1801  // vector along y-arm
1802  cout << "Ey : " << Ey[0] << " " << Ey[1] << " " << Ey[2] << " " << endl;
1803  cout << endl;
1804  cout.precision(precision);
1805  cout << "Ex-North Angle Clockwise (deg) : " << xdP.AzX << endl;
1806  cout << "Ey-North Angle Clockwise (deg) : " << xdP.AzY << endl;
1807  cout << endl;
1808  if(polarization==TENSOR) cout << "GW Polarization = TENSOR" << endl;
1809  if(polarization==SCALAR) cout << "GW Polarization = SCALAR" << endl;
1810  cout << endl;
1811 
1812  return;
1813 }
1814 
1815 //______________________________________________________________________________
1816 void detector::Streamer(TBuffer &R__b)
1817 {
1818  // Stream an object of class detector.
1819 
1820  UInt_t R__s, R__c;
1821  if (R__b.IsReading()) {
1822  Version_t R__v = R__b.ReadVersion(&R__s, &R__c); if (R__v) { }
1823  TNamed::Streamer(R__b);
1824  R__b.ReadStaticArray((char*)Name);
1825  R__b.StreamObject(&(dP),typeid(detectorParams));
1826  R__b.ReadStaticArray((double*)Rv);
1827  R__b.ReadStaticArray((double*)Ex);
1828  R__b.ReadStaticArray((double*)Ey);
1829  R__b.ReadStaticArray((double*)DT);
1830  R__b.ReadStaticArray((double*)ED);
1831  if(R__v > 2) R__b >> ifoID;
1832  R__b >> sHIFt;
1833  if(R__v > 1) {
1834  void *ptr_polarization = (void*)&polarization;
1835  R__b >> *reinterpret_cast<Int_t*>(ptr_polarization);
1836  }
1837 
1838  if(R__v > 3) {
1839  R__b >> wfSAVE;
1840  if(wfSAVE) {
1841  HRSS.Streamer(R__b);
1842  ISNR.Streamer(R__b);
1843  FREQ.Streamer(R__b);
1844  BAND.Streamer(R__b);
1845  TIME.Streamer(R__b);
1846  TDUR.Streamer(R__b);
1847  if(wfSAVE==1||wfSAVE==3) {
1848  {
1849  vector<int> &R__stl = IWFID;
1850  R__stl.clear();
1851  int R__i, R__n;
1852  R__b >> R__n;
1853  R__stl.reserve(R__n);
1854  for (R__i = 0; R__i < R__n; R__i++) {
1855  int R__t;
1856  R__b >> R__t;
1857  R__stl.push_back(R__t);
1858  }
1859  }
1860  {
1861  vector<wavearray<double>*> &R__stl = IWFP;
1862  R__stl.clear();
1863  int R__i, R__n;
1864  R__b >> R__n;
1865  R__stl.reserve(R__n);
1866  for (R__i = 0; R__i < R__n; R__i++) {
1868  R__t->Streamer(R__b);
1869  R__stl.push_back(R__t);
1870  }
1871  }
1872  }
1873  if(wfSAVE==2||wfSAVE==3) {
1874  {
1875  vector<int> &R__stl = RWFID;
1876  R__stl.clear();
1877  int R__i, R__n;
1878  R__b >> R__n;
1879  R__stl.reserve(R__n);
1880  for (R__i = 0; R__i < R__n; R__i++) {
1881  int R__t;
1882  R__b >> R__t;
1883  R__stl.push_back(R__t);
1884  }
1885  }
1886  {
1887  vector<wavearray<double>*> &R__stl = RWFP;
1888  R__stl.clear();
1889  int R__i, R__n;
1890  R__b >> R__n;
1891  R__stl.reserve(R__n);
1892  for (R__i = 0; R__i < R__n; R__i++) {
1894  R__t->Streamer(R__b);
1895  R__stl.push_back(R__t);
1896  }
1897  }
1898  }
1899  }
1900  }
1901 
1902  R__b.CheckByteCount(R__s, R__c, detector::IsA());
1903  } else {
1904  R__c = R__b.WriteVersion(detector::IsA(), kTRUE);
1905  TNamed::Streamer(R__b);
1906  R__b.WriteArray(Name, 16);
1907  R__b.StreamObject(&(dP),typeid(detectorParams));
1908  R__b.WriteArray(Rv, 3);
1909  R__b.WriteArray(Ex, 3);
1910  R__b.WriteArray(Ey, 3);
1911  R__b.WriteArray(DT, 9);
1912  R__b.WriteArray(ED, 5);
1913  R__b << ifoID;
1914  R__b << sHIFt;
1915  R__b << (Int_t)polarization;
1916 
1917  R__b << wfSAVE;
1918  if(wfSAVE) {
1919  HRSS.Streamer(R__b);
1920  ISNR.Streamer(R__b);
1921  FREQ.Streamer(R__b);
1922  BAND.Streamer(R__b);
1923  TIME.Streamer(R__b);
1924  TDUR.Streamer(R__b);
1925  if(wfSAVE==1||wfSAVE==3) {
1926  {
1927  vector<int> &R__stl = IWFID;
1928  int R__n=(&R__stl) ? int(R__stl.size()) : 0;
1929  R__b << R__n;
1930  if(R__n) {
1931  vector<int>::iterator R__k;
1932  for (R__k = R__stl.begin(); R__k != R__stl.end(); ++R__k) {
1933  R__b << (*R__k);
1934  }
1935  }
1936  }
1937  {
1938  vector<wavearray<double> > IWF(IWFP.size());
1939  for(int i=0;i<IWFP.size();i++) IWF[i] = *IWFP[i];
1940  vector<wavearray<double> > &R__stl = IWF;
1941  int R__n=(&R__stl) ? int(R__stl.size()) : 0;
1942  R__b << R__n;
1943  if(R__n) {
1944  vector<wavearray<double> >::iterator R__k;
1945  for (R__k = R__stl.begin(); R__k != R__stl.end(); ++R__k) {
1946  ((wavearray<double>&)(*R__k)).Streamer(R__b);
1947  }
1948  }
1949  }
1950  }
1951  if(wfSAVE==2||wfSAVE==3) {
1952  {
1953  vector<int> &R__stl = RWFID;
1954  int R__n=(&R__stl) ? int(R__stl.size()) : 0;
1955  R__b << R__n;
1956  if(R__n) {
1957  vector<int>::iterator R__k;
1958  for (R__k = R__stl.begin(); R__k != R__stl.end(); ++R__k) {
1959  R__b << (*R__k);
1960  }
1961  }
1962  }
1963  {
1964  vector<wavearray<double> > RWF(RWFP.size());
1965  for(int i=0;i<RWFP.size();i++) RWF[i] = *RWFP[i];
1966  vector<wavearray<double> > &R__stl = RWF;
1967  int R__n=(&R__stl) ? int(R__stl.size()) : 0;
1968  R__b << R__n;
1969  if(R__n) {
1970  vector<wavearray<double> >::iterator R__k;
1971  for (R__k = R__stl.begin(); R__k != R__stl.end(); ++R__k) {
1972  ((wavearray<double>&)(*R__k)).Streamer(R__b);
1973  }
1974  }
1975  }
1976  }
1977  }
1978 
1979  R__b.SetByteCount(R__c, kTRUE);
1980  }
1981 }
1982 
wavearray< double > t(hp.size())
detectorParams dP
Definition: detector.hh:311
const double _eT1x[3]
const double _eH1x[3]
void writeFilter(const char *)
param: file name
Definition: detector.cc:930
double sHIFt
Definition: detector.hh:317
const double _rE1[3]
detectorParams getDetectorParams()
Definition: detector.cc:201
virtual void resize(unsigned int)
Definition: wseries.cc:883
double sTARt
virtual size_t size() const
Definition: wavearray.hh:127
int wfSAVE
Definition: detector.hh:326
double imag() const
Definition: wavecomplex.hh:52
double getWFfreq(char atype='S')
Definition: detector.cc:1750
const double _eA1y[3]
static double g(double e)
Definition: GNGen.cc:98
par[0] value
detector D
Definition: TestBandPass.C:15
const double _eL1x[3]
int getMaxLevel()
Definition: wseries.cc:85
tuple f
Definition: cwb_online.py:91
void setFpFx(double, double=0., double=180., double=0., double=360.)
param - step on phi and theta param - theta begin param - theta end param - phi begin param - phi end...
Definition: detector.cc:692
printf("total live time: non-zero lags = %10.1f \n", liveTot)
size_t setsnr(wavearray< double > &, std::vector< double > *, std::vector< double > *, double=5., double=8.)
Definition: detector.cc:1641
int offset
Definition: TestSTFT_2.C:19
virtual void rate(double r)
Definition: wavearray.hh:123
char name[32]
Definition: detector.hh:32
plot hist2D SetName("WSeries-1")
wavearray< double > a(hp.size())
size_t frequency
Definition: netpixel.hh:93
double Ey[3]
Definition: detector.hh:314
double AzX
Definition: detector.hh:37
WSeries< float > v[nIFO]
Definition: cwb_net.C:62
double getWFtime(char atype='S')
Definition: detector.cc:1735
par[0] name
int n
Definition: cwb_net.C:10
const double _eE1x[3]
wavearray< double > z
Definition: Test10.C:32
int ID
Definition: TestMDC.C:70
cout<< endl;cout<< "ts size = "<< ts.size()<< " ts rate = "<< ts.rate()<< endl;tf.Forward(ts, wdm);int levels=tf.getLevel();cout<< "tf size = "<< tf.size()<< endl;double dF=tf.resolution();double dT=1./(2 *dF);cout<< "rate(hz) : "<< RATE<< "\t layers : "<< nLAYERS<< "\t dF(hz) : "<< dF<< "\t dT(ms) : "<< dT *1000.<< endl;int itime=TIME_PIXEL_INDEX;int ifreq=FREQ_PIXEL_INDEX;int index=(levels+1)*itime+ifreq;double time=itime *dT;double freq=(ifreq >0)?ifreq *dF:dF/4;cout<< endl;cout<< "PIXEL TIME = "<< time<< " sec "<< endl;cout<< "PIXEL FREQ = "<< freq<< " Hz "<< endl;cout<< endl;wavearray< double > x
double phase
std::vector< delayFilter > filter
Definition: detector.hh:343
wavearray< double > HRSS
Definition: detector.hh:353
float theta
double getlow()
Definition: netcluster.hh:329
double real() const
Definition: wavecomplex.hh:51
STL namespace.
double getTheta(size_t i)
Definition: skymap.hh:206
std::slice getSlice(double n)
Definition: wseries.hh:134
const double _rT1[3]
size_t setFilter(size_t, double=0., size_t=0)
param: filter length param: filter phase delay in degrees param: number of up-sample layers return fi...
Definition: detector.cc:742
double AltX
Definition: detector.hh:36
double AzY
Definition: detector.hh:39
waveform wf
int polarization
bool setrms(netcluster *, size_t=0)
Definition: detector.cc:1208
const double _rH1[3]
watplot p1(const_cast< char * >("TFMap1"))
#define M
Definition: UniqSLagsList.C:3
size_t layers
Definition: netpixel.hh:94
int m
Definition: cwb_net.C:10
virtual void start(double s)
Definition: wavearray.hh:119
int j
Definition: cwb_net.C:10
size_t minDFindex(delayFilter &F)
Definition: detector.hh:382
void delay(double, double)
param: theta param: phi
Definition: detector.cc:1682
i drho i
double longitude
Definition: detector.hh:34
skymap tau
Definition: detector.hh:328
#define N
virtual double mean() const
Definition: wavearray.cc:1053
void init()
param: Rx,Ry,Rz in ECEF frame
Definition: detector.cc:451
double getNoise(size_t, int=-1)
param: wavelet layer index (frequency) param: wavelet time index if param 2 is specified - return noi...
Definition: detector.cc:1112
#define PI
Definition: watfun.hh:14
const double _rN1[3]
skymap mFx
Definition: detector.hh:330
const double _eO1x[3]
virtual double rms()
Definition: wavearray.cc:1188
wavearray< double > w
Definition: Test1.C:27
void wrate(double r)
Definition: wseries.hh:102
void readFilter(const char *)
Definition: detector.cc:969
float phi
const double _eN1y[3]
const double _eV1y[3]
wavecomplex antenna(double, double, double=0.)
param: source theta,phi, polarization angle psi in degrees
Definition: detector.cc:473
double hrss
Definition: TestMDC.C:70
float psi
double time[6]
Definition: cbc_plots.C:435
const double _eJ1x[3]
void set(double x, double y)
Definition: wavecomplex.hh:56
const double _eE1y[3]
detector & operator=(const detector &)
Definition: detector.cc:351
void bandPass1G(double f1=0., double f2=0.)
Definition: detector.cc:1281
const double _eT1y[3]
wavearray< double > hot[2]
double deg2rad
double getwave(int, WSeries< double > &, char='W', size_t=0)
param: cluster ID param: WSeries where to put the waveform return: noise rms
Definition: netcluster.cc:2743
size_t size() const
Definition: wslice.hh:71
WSeries< double > nRMS
Definition: detector.hh:340
const double _rG1[3]
const double _rL1[3]
static const double SM
Definition: GNGen.cc:8
wavearray< double > TIME
Definition: detector.hh:357
void GeocentricToGeodetic(double X, double Y, double Z, double &latitude, double &longitude, double &elevation)
Definition: skycoord.hh:289
i() int(T_cor *100))
double Pi
TF1 * f2
Definition: cbc_plots.C:1710
virtual int convertF2L(int, int)
Definition: Wavelet.cc:75
const double _eG1x[3]
const double _eO1y[3]
double start
Definition: netcluster.hh:361
int getLayer(wavearray< DataType_t > &w, double n)
param: n - layer number
Definition: wseries.cc:175
double elevation
Definition: detector.hh:35
static double tau
Definition: geodesics.cc:8
const double _rV1[3]
const double _eJ1y[3]
char fname[1024]
POLARIZATION polarization
Definition: detector.hh:366
Definition: Wavelet.hh:31
const double _eH2x[3]
double precision
#define speedlight
Definition: watfun.hh:15
detector & operator>>(detector &)
Definition: detector.cc:421
int k
virtual void resetFFTW()
Definition: wavearray.cc:959
const double _eN1x[3]
std::vector< wavearray< double > * > IWFP
Definition: detector.hh:361
detector()
Definition: detector.cc:90
const double _rO1[3]
static double A
Definition: geodesics.cc:8
virtual int convertO2F(int, int)
Definition: Wavelet.cc:87
int m_Level
Definition: Wavelet.hh:97
double F
Definition: skymap.hh:45
double latitude
Definition: detector.hh:33
wavearray< double > TDUR
Definition: detector.hh:358
std::vector< int > IWFID
Definition: detector.hh:360
double e
void setLevel(size_t n)
Definition: wseries.hh:94
WSeries< double > wM
Definition: cwb_job_obj.C:23
size_t size()
Definition: netcluster.hh:129
double gethigh()
Definition: netcluster.hh:330
wavearray< double > BAND
Definition: detector.hh:356
netpixel * getPixel(size_t n, size_t i)
Definition: netcluster.hh:395
double Ex[3]
Definition: detector.hh:313
virtual void FFTW(int=1)
Definition: wavearray.cc:878
regression r
Definition: Regression_H1.C:44
double Rv[3]
Definition: detector.hh:312
s s
Definition: cwb_net.C:137
double rad2deg
char filter[1024]
double getPhi(size_t i)
Definition: skymap.hh:146
WSeries< double > TFmap
Definition: detector.hh:336
double T
Definition: testWDM_4.C:11
void GeodeticToGeocentric(double latitude, double longitude, double elevation, double &X, double &Y, double &Z)
Definition: skycoord.hh:197
virtual Wavelet * Clone() const
return: Wavelet* - duplicate of *this, allocated on heap
Definition: Wavelet.cc:42
const double _eV1x[3]
virtual void waveSort(DataType_t **pp, size_t l=0, size_t r=0) const
Definition: wavearray.cc:1403
WSeries< double > waveForm
Definition: detector.hh:337
double AltY
Definition: detector.hh:38
const double _eL1y[3]
watplot p2(const_cast< char * >("TFMap2"))
const double _rH2[3]
const double _rJ1[3]
wavearray< int > index
Definition: Meyer.hh:18
char Name[16]
Definition: detector.hh:309
skymap mFp
Definition: detector.hh:329
double fabs(const Complex &x)
Definition: numpy.cc:37
void rotate(double)
Definition: detector.cc:274
virtual ~detector()
Definition: detector.cc:326
const double _eG1y[3]
void print()
Definition: detector.cc:1768
std::vector< short > index
Definition: detector.hh:27
virtual int getOffset(int, int)
Definition: Wavelet.cc:51
void Forward(int n=-1)
param: wavelet - n is number of steps (-1 means full decomposition)
Definition: wseries.cc:228
TString OS
Definition: cwb_rootlogon.C:25
Meyer< double > S(1024, 2)
int l
Definition: cbc_plots.C:434
sprintf(tfres,"(1/%g)x(%g) (sec)x(Hz)", 2 *df, df)
const double _eA1x[3]
double getTau(double, double)
param: source theta,phi angles in degrees
Definition: detector.cc:681
const double _eH2y[3]
void set(size_t i, double a)
param: sky index param: value to set
Definition: skymap.hh:104
double getwave(int, netcluster &, char, size_t)
param: no parameters
Definition: detector.cc:544
DataType_t * data
Definition: wavearray.hh:301
double null
Definition: detector.hh:318
netcluster wc
wavearray< double > FREQ
Definition: detector.hh:355
enum WAVETYPE m_WaveType
Definition: Wavelet.hh:88
TH1 * t1
WaveDWT< DataType_t > * pWavelet
Definition: wseries.hh:438
double dT
Definition: testWDM_5.C:12
fclose(ftrig)
void GetCartesianComponents(double u[3], double Alt, double Az, double Lat, double Lon)
Definition: skycoord.hh:369
float rate
Definition: netpixel.hh:95
size_t size()
Definition: skymap.hh:118
double getdata(char type='R', size_t n=0)
Definition: netpixel.hh:56
size_t stride() const
Definition: wslice.hh:75
WSeries< float > nVAR
Definition: detector.hh:341
std::vector< float > value
Definition: detector.hh:28
virtual void resize(unsigned int)
Definition: wavearray.cc:445
void Inverse(int n=-1)
param: n - number of steps (-1 means full reconstruction)
Definition: wseries.cc:273
wavearray< double > y
Definition: Test10.C:31
void putLayer(wavearray< DataType_t > &, double n)
param: n - layer number
Definition: wseries.cc:201
size_t setsim(WSeries< double > &, std::vector< double > *, double=5., double=8., bool saveWF=false)
Definition: detector.cc:1348
const double _eH1y[3]
double frequency(int l)
Definition: wseries.cc:99
const double _rA1[3]
wavearray< double > ISNR
Definition: detector.hh:354
size_t start() const
Definition: wslice.hh:67
int maxLayer()
Definition: wseries.hh:121
int m_H
Definition: Wavelet.hh:103
detector & operator<<(detector &)
Definition: detector.cc:394
void setTau(double, double=0., double=180., double=0., double=360.)
param - step on phi and theta param - theta begin param - theta end param - phi begin param - phi end...
Definition: detector.cc:637
init()
Definition: revMonster.cc:12
exit(0)
void cpf(const wavearray< DataType_t > &, int=0, int=0, int=0)
Definition: wavearray.cc:699
virtual WSeries< double > white(double, int, double=0., double=0.)
what it does: each wavelet layer is devided into k intervals.
Definition: wseries.cc:1128
bool setdata(double a, char type='R', size_t n=0)
Definition: netpixel.hh:40