Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
network.hh
Go to the documentation of this file.
1 //**************************************************************
2 // Wavelet Analysis Tool
3 // Sergey Klimenko, University of Florida
4 // class for coherent network analysis used with DMT and ROOT
5 //**************************************************************
6 
7 #ifndef NETWORK_HH
8 #define NETWORK_HH
9 
10 #include <iostream>
11 #include <stdlib.h>
12 #include <string.h>
13 #include <vector>
14 #include "TH2F.h"
15 #include "wavearray.hh"
16 #include "wseries.hh"
17 #include "detector.hh"
18 #include "skymap.hh"
19 #include "netcluster.hh"
20 #include "wat.hh"
21 #include "monster.hh"
22 #include "watplot.hh" // remove
23 #include "TMatrixDSym.h"
24 #include "TMatrixDSymEigen.h"
25 #include "TVectorD.h"
26 #ifndef __CINT__
27 #include "watsse.hh"
28 #include "watavx.hh"
29 #endif
30 
31 using namespace std;
32 
33 typedef std::vector<double> vectorD;
34 
35 struct waveSegment {
36  int index;
37  double start;
38  double stop;
39 };
40 
41 class network : public TNamed
42 {
43  public:
44 
45  // constructors
46 
47  //: Default constructor
48  network();
49 
50  //: Copy constructor
51  //!param: value - object to copy from
52  network(const network&);
53 
54  //: destructor
55  virtual ~network();
56 
57  // operators
58 
59  network& operator= (const network&);
60 
61  // accessors
62 
63  //: add detector to the network
64  //!param: detector structure
65  //!return number of detectors in the network
66  size_t add(detector*);
67 
68  //:do forward/inverse wavelet transformation by k steps
69  //!param: number of steps
70  // do not use with WDM transform !!!
71  void Forward(size_t k)
72  { for(size_t i=0; i<ifoList.size(); i++) ifoList[i]->TFmap.Forward(k); }
73  void Inverse(size_t k)
74  { for(size_t i=0; i<ifoList.size(); i++) ifoList[i]->TFmap.Inverse(k); }
75 
76  // set time shifts
77  //!param number of time lags
78  //!param time shift step in seconds
79  //!param first lag ID
80  //!param maximum lag ID
81  //!param file name for lag configurations
82  //!param r/w/s - read/write/string mode
83  int setTimeShifts(size_t=1, double=1., size_t=0, size_t=0,
84  const char* = NULL, const char* = "w", size_t* = NULL);
85 
86  // print wc_List
87  //param - time lag
88  void printwc(size_t);
89 
90  //: initialize cluster for selected TF area, put it in wc_List[0];
91  //!param: cluster start time relative to segment start
92  //!param: cluster duration
93  //!return cluster list size
94  // time shifts between detectors are controlled by sHIFt parameters in the detector objects
95  // cluster bandwidth is controlled by TFmap.low(), TFmap.high() parameters.
96  virtual size_t initwc(double, double);
97 
98  //: initialize network sky maps: fill netvork sensitivity and alignment factor
99  //!param: sky map granularity step, degrees
100  //!param: theta begin, degrees
101  //!param: theta end, degrees
102  //!param: phi begin, degrees
103  //!param: phi end, degrees
104  void setSkyMaps(double,double=0.,double=180.,double=0.,double=360.);
105 
106  //: initialize network sky maps: fill netvork sensitivity and alignment factor
107  //!param: healpix order
108  void setSkyMaps(int);
109 
110  //:calculate network data matrix (NDM)
111  //!param: cluster ID
112  //!param: lag index
113  //!param: statistic identificator
114  //!param: resolution idenificator
115  //!return: status
116  bool setndm(size_t, size_t, bool=true, int=1);
117  bool SETNDM(size_t, size_t, bool=true, int=1); // used with likelihoodI
118 
119  //:read NDM element defined by detectors i,j
120  //!param: first detector
121  //!param: second detector
122  inline double getNDM(size_t i, size_t j) { return NDM[i][j]; }
123 
124  //:set delay filters for a network
125  //:(-t,t)_i=0, (-t,t)_i=1, .....
126  //!param: detector
127  // return number of delays for each layer.
128  size_t setFilter(detector* = NULL); // from detector object
129 
130  //:set delay filters and index arrays for a network
131  void setDelayFilters(detector* = NULL); // from detector
132  void setDelayFilters(char*, char* = NULL); // from detector filter files
133  void setFilter(char*, char* = NULL); // from network filter files
134 
135  //: Dumps network filter to file *fname in binary format.
136  void writeFilter(const char *fname);
137 
138  //: reads network filter from file
139  void readFilter(const char*);
140 
141  // set veto array from the input list of DQ segments
142  //!param: time window around injections
143  double setVeto(double=5.);
144 
145  //: read injections
146  //!param: MDC log file
147  //!param: approximate gps time
148  // (injections in a window +-10000 sec around gps are selected)
149  //!param: position of time field
150  //!param: position of type field
151  // returns number of injections
152  size_t readMDClog(char*,double=0.,int=11,int=12);
153 
154  //: read segment list
155  //!param: segment list file
156  //!param: start time collumn number
157  // returns number of segments
158  size_t readSEGlist(char*,int=1);
159 
160  //: read DQ segments
161  //!param: MDC log file
162  // returns number of injections
163  //size_t readSegments(char*,int=11,int=12);
164 
165  // set optimal delay index
166  // new version which works with WDM delay filters
167  // the delay index is symmetric - can be negative.
168  // index=0 corresponds to zero delay
169  // index>0 corresponds to positive delay (shift right)
170  // index<0 corresponds to negative delay (shift left)
171  // rate: effective data rate (1/delay_step)
172  void setDelayIndex(double rate);
173 
174  // set optimal delay index
175  // work with dumb wavelet 1G delay filters
176  //!param: dummy
177  void setDelayIndex(int=0);
178 
179  //:set sky index array
180  //: mode: 0 - now downselection
181  //: 1 - exclude duplicate delay configurations
182  //: 2 - downselect by 2
183  //: 4 - downselect by 4
184  size_t setIndexMode(size_t=0);
185 
186  //:get sky index for given theta, phi
187  //!param: theta [deg]
188  //!param: phi [deg]
189  inline int getIndex(double theta, double phi) {
190  return getifo(0)->tau.getSkyIndex(theta,phi);
191  }
192 
193  //:set antenna pattern buffers
194  //!param: detector (use theta, phi index array)
195  void setAntenna(detector*);
196  // set antenna patterns for all detectors
197  void setAntenna();
198 
199  // 1G: delay detectors in the network with respect to reference
200  // to match sky location theta and phi
201  // index array should be setup
202  void delay(double theta, double phi);
203  // 1G delay detector in the network:
204  // time delay convention: + - shift TS right
205  // - - shift TS left
206  //!param: detector pointer
207  //!param: delay filter index
208  void delay(detector*, size_t);
209 
210  // 1G analysis algorithm for selection of significant network pixles
211  //:select TF pixels by setting a threshold on the pixel network likelihood.
212  //:integrated time shifts for estimation of the background
213  //:input parameters define threshold on TF pixel parameters
214  //!param: threshold on lognormal pixel energy (in units of noise rms)
215  //!param: threshold on total pixel energy (in units of noise rms)
216  //!param: threshold on likelihood (in units of amplitude SNR per detector)
217  //!return: number of selected samples.
218  long coherence(double, double=0., double=0.);
219 
220  // 2G analysis algorithm for selection of significant network pixles
221  // LAG - time shift lag defining how detectors are shifted wrt each other.
222  // Eo - pixel energy threshold
223  // DD - dummy
224  // hist- pointer to a diagnostic histogram showing network pixel energy.
225  long getNetworkPixels(int LAG, double Eo, double DD=1., TH1F* hist=NULL);
226 
227  //:selection of clusters based on:
228  // 'C' - network correlation coefficient
229  // 'l' - likelihood (biased)
230  // 'L' - likelihood (unbiased)
231  // 'A' - snr - null assymetry (double OR)
232  // 'E' - snr - null energy (double OR)
233  // 'a' - snr - null assymetry (strict)
234  // 'e' - snr - null energy (strict)
235  //!param: threshold
236  //!param: minimum cluster size processed by the corrcut
237  //!param: cluster type
238  //!return: number of rejected pixels.
239  size_t netcut(double, char='L', size_t=0, int=1);
240 
241  // apply subnetwork cut
242  // subnet: sub network threshold
243  // subcut: sub network threshold in the skyloop
244  // lag: lag index
245  // return number of processed pixels
246  // works only with TD amplitudes.
247  // Includes x regulator for background rejection
248  long subNetCut(int lag, float subnet=0.6, float subcut=0.33, TH2F* hist=NULL);
249 
250  // calculate network likelihood for reconstructed clusters when
251  // independent h+ and hx polarisations are reconstructed
252  // implementation for 2-5 detector network.
253  //!param: maximized statistic:
254  //!param: threshold to define core pixels (in units of noise rms)
255  // effective only if the second parameter is false
256  //!param: cluster ID, if ID=0 - fill likelihood field for all clusters
257  // otherwise, calculate likelihood only for specified cluster.
258  // fill nLikelihood skymap if ID<=0 - need for CED
259  //!param: lag index
260  //!param: sky index
261  //!param: true - for core pixels, false - for core & halo pixels
262  //!return number of processed pixels
263  // Options for sky statistics
264  // 'e'/'E' - power
265  // 'b'/'B' - un-modeled search with 0 phase data
266  long likelihoodB(char='E', double=sqrt(2.), int=0, size_t=0, int=-1, bool=false);
267 
268  // implementation for 2-5 detector network and elliptical constraint
269  // 'p'/'P' - un-modeled search with 0 and 90deg phase data
270  // 'i'/'I' - (inspiral) elliptical constraint
271  // 's'/'S' - (supernova) linear constraint
272  // 'g'/'G' - (short GRB) circular constraint
273  long likelihoodI(char='P', double=sqrt(2.), int=0, size_t=0, int=-1, bool=false);
274 
275  // 2G likelihood with multi-resolution cluster analysis
276  // mode: analysis mode
277  // lag: lag index
278  // ID: cluster ID, if ID=0 - process all clusters
279  // fill nLikelihood skymap if ID<=0 - need for CED
280  // hist: chirp histogram: If not needed, TGraphErrors* hist=NULL
281  // shold be used as input
282  //return number of processed pixels
283  long likelihood2G(char mode, int lag, int ID, TH2F* hist=NULL);
284  long likelihoodWP(char mode, int lag, int ID, TH2F* hist=NULL);
285 
286  // combines likelihood0() and likelihoodI() and likelihoodB()
287  long likelihood(char='E', double=sqrt(2.), int=0, size_t=0, int=-1, bool=false);
288 
289  // set pixel ranks in the network netcluster structure
290  // param: time window around pixel for rank calculation
291  // param: frequency window around pixel for rank calculation (not used)
292  size_t setRank(double, double=0.);
293 
294 
295  // read the wdm xtalk catalog
296  // param: catalog filename
297  inline void setMRAcatalog(char* fn){ wdmMRA.read(fn);}
298 
299 
300  //:reconstruct time clusters
301  //!param: time gap in pixels
302  //!return: number of reconstructed clusters
303  size_t cluster(int kt, int kf) {
304  if(!wc_List.size()) return 0;
305  size_t m = 0;
306  for(size_t n=0; n<nLag; n++) m += wc_List[n].cluster(kt,kf);
307  return m;
308  }
309 
310  //:return number of events
311  size_t events() {
312  if(!wc_List.size()) return 0;
313  size_t m = 0;
314  for(size_t n=0; n<nLag; n++) {
315  m += wc_List[n].esize();
316  }
317  return m;
318  }
319 
320  //:return number of events with specified type, lag and cID
321  // type : sCuts
322  // lag : lag<0 -> select all lags
323  size_t events(int type, int lag=-1) {
324  if(!wc_List.size()) return 0;
325  if(lag>=(int)nLag) lag=nLag-1;
326  size_t m = 0;
327  if(lag>=0) m += wc_List[lag].esize(type);
328  else for(size_t n=0; n<nLag; n++) m += wc_List[n].esize(type);
329  return m;
330  }
331 
332  // set noise rms fields for pixels in the network netcluster structure
333  inline void setRMS();
334 
335  // delink superclusters in wc_List
336  inline void delink(){
337  for(size_t i=0; i<nLag; i++) wc_List[i].delink();
338  }
339 
340  //: extract time series for detector responses
341  //!param: cluster ID
342  //!param: delay index
343  //!param: time series type
344  //!return: true if time series are extracted
345  bool getwave(size_t, size_t, char='W');
346 
347  // get MRA waveforms of type atype in time domain given lag nomber and cluster ID
348  // if tof = true, apply time-of-flight corrections
349  // fill in waveform arrays in the detector class
350  // atype = 'W' - get whitened detector Wavelet PC
351  // atype = 'w' - get detector Wavelet PC
352  // atype = 'S' - get whitened reconstructed response (Signal)
353  // atype = 's' - get reconstructed response (Signal)
354  // mode: -1/0/1 - return 90/mra/0 phase
355  // tof: false/true - time-of-flight correction
356  bool getMRAwave(size_t ID, size_t lag, char atype='S', int mode=0, bool tof=false);
357 
358  //: calculation of sky error regions
359  //!param: cluster id
360  //!param: lag
361  //!param: cluster time
362  void getSkyArea(size_t id, size_t lag, double T);
363  //!param: noise rms per DoF
364  void getSkyArea(size_t id, size_t lag, double T, double rms);
365 
366  // read earth skyMask coordinates from file (1G)
367  // parameter - fraction of skymap to be rejected
368  // parameter - file name
369  size_t setSkyMask(double f, char* fname);
370 
371  // read celestial/earth skyMask coordinates from file
372  // parameter - file name
373  // parameter - sky coordinates : 'e'=earth, 'c'=celestial
374  size_t setSkyMask(char* fname, char skycoord);
375 
376  // read celestial/earth skyMask coordinates from skymap
377  // parameter - skymap
378  // parameter - sky coordinates : 'e'=earth, 'c'=celestial
379  size_t setSkyMask(skymap sm, char skycoord);
380 
381  // set threshold on amplitude of core pixels
382  // set threshold on subnetwork energy
383  inline void setAcore(double a) {
384  this->acor = a;
385  a = 1-Gamma(ifoList.size(),a*a*ifoList.size()); // probability
386  this->e2or = iGamma(ifoList.size()-1,a); // subnetwork energy threshold
387  }
388 
389  //: get size of mdcList
390  inline size_t mdcListSize() { return mdcList.size(); }
391  //: get size of mdcType
392  inline size_t mdcTypeSize() { return mdcType.size(); }
393  //: get size of mdcTime
394  inline size_t mdcTimeSize() { return mdcTime.size(); }
395  //: get size of mdc__ID
396  inline size_t mdc__IDSize() { return mdc__ID.size(); }
397  //: get size of livTime
398  inline size_t livTimeSize() { return livTime.size(); }
399  //: get element of mdcList
400  inline string getmdcList(size_t n) { return mdcListSize()>n ? mdcList[n] : "\n"; }
401  //: get element of mdcType
402  inline string getmdcType(size_t n) { return mdcTypeSize()>n ? mdcType[n] : "\n"; }
403  //: get pointer to mdcTime
404  inline std::vector<double>* getmdcTime() { return &mdcTime; }
405  //: get element of mdcTime
406  inline double getmdcTime(size_t n) { return mdcTimeSize()>n ? mdcTime[n] : 0.; }
407  //: get element of mdc__ID
408  inline size_t getmdc__ID(size_t n) { return mdc__IDSize()>n ? mdc__ID[n] : 0; }
409  //: get element of livTime
410  inline double getliveTime(size_t n) { return livTimeSize()>n ? livTime[n] : 0.; }
411 
412  //: get size of ifoList
413  inline size_t ifoListSize() { return ifoList.size(); }
414  //: get size of wc_List
415  inline size_t wc_ListSize() { return wc_List.size(); }
416  //:return pointer to a detector n
417  //!param: detector index
418  inline detector* getifo(size_t n) { return ifoListSize()>n ? ifoList[n] : NULL; }
419  //:return pointer to a netcluster in wc_List
420  //!param: delay index
421  inline netcluster* getwc(size_t n) { return wc_ListSize()>n ? &wc_List[n] : NULL; }
422 
423  //: add wdm
424  //!param: pointer to wdm
425  //!return number of wdm tronsforms in the list
426  size_t add(WDM<double>* wdm) {if(wdm) wdmList.push_back(wdm); return wdmList.size();}
427  //: get size of wdmList
428  inline size_t wdmListSize() { return wdmList.size(); }
429  //!param: number of wdm layers
430  inline WDM<double>* getwdm(size_t M) {
431  for(size_t n=0;n<wdmListSize();n++) if(M==(wdmList[n]->maxLayer()+1)) return wdmList[n];
432  return NULL;
433  }
434 
435  // set run number
436  //!param: run
437  inline void setRunID(size_t n) { this->nRun=n; return; }
438 
439  // set wavelet boundary offset
440  //!param: run
441  inline void setOffset(double t) { this->Edge = t; return; }
442 
443  // set constraint parameter
444  //!param: constraint parameter, p=0 - no constraint
445  inline void constraint(double d=1., double g=0.0001) {
446  this->delta = d==0. ? 0.00001 : d;
447  this->gamma = g;
448  //this->gamma=g>0. ? g : 0.00001;
449  }
450 
451  // set threshold for double OR energy
452  //!param: threshold
453  inline void set2or(double p) { this->e2or = p; }
454 
455  // calculate WaveBurst threshold as a function of resolution
456  //!param: selected fraction of LTF pixels assuming Gaussian noise
457  //!param: maximum time delay between detectors
458  double threshold(double, double);
459 
460  // calculate WaveBurst threshold as a function of resolution
461  //!param: selected fraction of LTF pixels assuming Gaussian noise
462  double THRESHOLD(double bpp);
463 
464  // calculate WaveBurst threshold for patterns
465  //!param: selected fraction of LTF pixels assuming Gaussian noise
466  //!param: Gamma distribution shape parameter
467  double THRESHOLD(double bpp, double shape);
468 
469  // calculate maximum delay between detectors
470  // "max" - maximum delay
471  // "min" - minimum delay
472  // "MAX" - max(abs(max),abs(min)
473  // def - max-min
474  double getDelay(const char* c="");
475 
476  //:recalculate time delay skymaps from Earth center to
477  //:'L1' - Livingston
478  //:'H1' - Hanford
479  //:'V1' - Cashina
480  //:'G1' - Hannover
481  //:'T1' - Tama
482  //:'A1' - Australia
483  void setDelay(const char* = "L1");
484 
485  // extract accurate time delay amplitudes for a given sky location
486  // parameter 1 - sky location index
487  // parameter 2 - array for 0-phase delayed amplitudes
488  // parameter 3 - array for 0-phase delayed amplitudes
489  void updateTDamp(int, float**, float**);
490 
491  // return WDM used/unused
492  bool wdm() {return _WDM;}
493  // get likelihood type, X=likelihoodX, M=likelihoodM, ''=others
494  char like() {return _LIKE;}
495 
496  // print network parameters
497  void print(); // *MENU*
498  virtual void Browse(TBrowser*) {print();}
499 
500  static inline double sumx(double*);
501  static inline double dotx(double*, double*);
502  static inline double dotx(float*, float*);
503  static inline double dot4(double*, double*);
504  static inline double dotx(double*, double*, double*);
505  static inline double dotx(float*, float*, float*);
506  static inline double dot4(double*, double*, double*);
507  static inline double dotx(double*, double**, size_t);
508  static inline double dotx(double**, size_t, double*);
509  static inline double dotx(double**, size_t, double**, size_t);
510  static inline double dotx(double**, size_t, double**, size_t, double*);
511  static inline double dotx(double*, double*, double);
512  static inline double dotx(float*, float*, float);
513  static inline void addx(double*, double*, double*);
514  static inline void addx(double*, double**, size_t, double*);
515  static inline void addx(double**, size_t, double**, size_t, double*);
516  static inline double dotx(double*, double**, size_t, double*);
517  static inline double dot32(std::vector<float>*, double*, std::vector<short>*);
518  static inline double dot32(double*, double*, int*);
519  static inline double divx(double*, double*);
520  static inline double rotx(double*, double, double*, double, double*);
521  static inline double rotx(float*, float, float*, float, float*);
522  static inline double rot4(double*, double, double*, double, double*);
523  static inline float rots(float*, float, float*, float, float*);
524  static inline void mulx(double**, size_t, double**, size_t, double*);
525  static inline void mulx(double*, double, double*);
526  static inline void mulx(float*, float, float*);
527  static inline void mulx(double*, double);
528  static inline void mulx(float*, float);
529  static inline void inix(double**, size_t, double*);
530  static inline void inix(double*, double);
531  static inline void inix(float*, float);
532  static inline int netx(double*, double, double*, double, double);
533  static inline int netx(float*, float, float*, float, float);
534  static inline void pnt_(float**, float**, short**, int, int);
535  static inline void cpp_(float*&, float**);
536  static inline void cpf_(float*& a, double** p);
537  static inline void cpf_(float*& a, double** p, size_t);
538 
539  static inline void dpfx(float* fp, float* fx);
540  static inline void pnpx(float* fp, float* fx, float* am, float* AM, float* u, float* v);
541  static inline void dspx(float* u, float* v, float* am, float* AM);
542  static inline void dspx(float* fp, float* fx, float* am, float* AM, float* u, float* v);
543 
544  inline int _sse_MRA_ps(float*, float*, float, int);
545  inline int _sse_mra_ps(float*, float*, float, int);
546  inline wavearray<float> _avx_norm_ps(float**, float**, std::vector<float*> &, int);
547  inline wavearray<float> _avx_norm_ps(float**, float**, float*, int);
548  inline void _avx_saveGW_ps(float**, float**, int);
549 
550  void test_sse(int, int);
551 
552 // data members
553 
554  size_t nRun; // run number
555  size_t nLag; // number of time lags
556  long nSky; // number of pixels for sky probability area
557  size_t mIFO; // master IFO
558  double rTDF; // effective rate of time-delay filter
559  double Step; // time shift step
560  double Edge; // time offset at the boundaries
561  double gNET; // network sensitivity
562  double aNET; // network alignment
563  double iNET; // network index
564  double eCOR; // correlation energy
565  double norm; // norm factor
566  double e2or; // threshold on double OR energy
567  double acor; // threshold on coherent pixel energy
568  bool pOUT; // true/false printout flag
569  bool EFEC; // true/false - EFEC/selestial coordinate system
570  char tYPe; // likelihood type
571  bool local; // true/false - local/global normalization
572  bool optim; // true/false - process optimal/all resolutions
573  double delta; // weak constraint parameter:
574  double gamma; // hard constraint parameter:
575  double precision; // precision of energy calculation
576  double pSigma; // integration limit in sigmas for probability
577  double penalty; // penalty factor:
578  double netCC; // threshold on netcc:
579  double netRHO; // threshold on rho:
580  bool eDisbalance; // true/false - enable/disable energy disbalance ECED
581  bool MRA; // true/false - used/not-used likelihoodMRA
582  bool wfsave; // true/false - if false only simulated wf are saved (2G)
583  int pattern; // clustering pattern
584 
585  std::vector<vectorD> NDM; // network data matrix
586 
587  WSeries<double> whp; // + polarization
588  WSeries<double> whx; // x polarization
589 
590  std::vector<detector*> ifoList; // detectors
591  std::vector<char*> ifoName; // detector's names
592  std::vector<netcluster> wc_List; // netcluster structures for time shifts
593  std::vector<double> livTime; // live time for time shifts
594  std::vector<std::string> mdcList; // list of injections
595  std::vector<std::string> mdcType; // list of injection types
596  std::vector<double> mdcTime; // gps time of selected injections
597  std::vector<size_t> mdc__ID; // ID of selected injections
598  std::vector<waveSegment> segList; // DQ segment list
599  std::vector<WDM<double>*> wdmList; //! list of wdm tranformations
600 
601  skymap nSensitivity; // network sensitivity
602  skymap nAlignment; // network alignment factor
603  skymap nCorrelation; // network correlation coefficient
604  skymap nLikelihood; // network likelihood
605  skymap nNullEnergy; // network null energy
606  skymap nPenalty; // signal * noise penalty factor
607  skymap nCorrEnergy; // reduced correlated energy
608  skymap nNetIndex; // network index
609  skymap nDisbalance; // energy disbalance
610  skymap nSkyStat; // sky optimization statistic
611  skymap nEllipticity; // waveform ellipticity
612  skymap nPolarisation; // polarisation angle
613  skymap nProbability; // probability skymap
614  skymap nAntenaPrior; // network sensitivtiy used as a prior for skyloc
615 
616  WSeries<double> pixeLHood; // pixel likelihood statistic
617  WSeries<double> pixeLNull; // pixel null statistic
618 
619  std::vector<delayFilter> filter; // delay filter (1G)
620  std::vector<delayFilter> filter90; // phase shifted delay filter (1G)
621 
622  wavearray<int> index; // theta, phi mask index array
623  wavearray<short> skyMask; // index array for setting sky mask
624  wavearray<double> skyMaskCC; // index array for setting sky mask Celestial Coordinates
625  wavearray<double> skyHole; // static sky mask describing "holes"
626  wavearray<short> veto; // veto array for pixel selection
627  wavearray<double> skyProb; // sky probability
628  wavearray<double> skyENRG; // energy skymap
629 
630 // data arrays for MRA and likelihood analysis
631 
632  std::vector<netpixel*> pList; //! list of pixel pointers for MRA
633  monster wdmMRA; //! wdm multi-resolution analysis
634  wavearray<float> a_00; //! buffer for cluster sky 00 amplitude
635  wavearray<float> a_90; //! buffer for cluster sky 90 amplitudes
636  wavearray<float> rNRG; //! buffers for cluster residual energy
637  wavearray<float> pNRG; //! buffers for cluster MRA energy
638 
639 // data arrays for polar coordinates storage : [0,1] = [radius,angle]
640 
641  wavearray<double> p00_POL[2]; //! buffer for projection on network plane 00 ampl
642  wavearray<double> p90_POL[2]; //! buffer for projection on network plane 90 ampl
643  wavearray<double> r00_POL[2]; //! buffer for standard response 00 ampl
644  wavearray<double> r90_POL[2]; //! buffer for standard response 90 ampl
645 
646 private:
647 
648  void like(char _LIKE) {this->_LIKE=_LIKE;} // set likelihood type
649  void wdm(bool _WDM) {this->_WDM =_WDM; } // set wdm used/unused
650  bool _WDM; // true/false - used/not-used WDM
651  char _LIKE; // X=likelihoodX, M=likelihoodM, ''=others
652 
653  ClassDef(network,4)
654 
655 }; // class network
656 
657 //inline void network::setCatalogFile(char* fn)
658 //{ wdmOvlp.read(fn);}
659 
660 // set noise rms fields for pixels in the network netcluster structure
661 inline void network::setRMS() {
662  size_t n = ifoList.size();
663  if(!ifoList.size() || wc_List.size()!=nLag) return;
664  for(size_t i=0; i<n; i++) {
665  for(size_t j=0; j<nLag; j++) {
666  if(!getwc(j)->size()) continue;
667  if(!getifo(i)->setrms(getwc(j),i)) {
668  cout<<"network::setRMS() error\n";
669  exit(1);
670  }
671  }
672  }
673  return;
674 }
675 
676 // special functions
677 inline double network::sumx(double* a) {
678  double d=0.;
679  NETX(d+= a[0]; ,
680  d+= a[1]; ,
681  d+= a[2]; ,
682  d+= a[3]; ,
683  d+= a[4]; ,
684  d+= a[5]; ,
685  d+= a[6]; ,
686  d+= a[7]; )
687  return d;
688 }
689 
690 inline double network::dotx(double* a, double* b) {
691  double d=0.;
692  NETX(d+= a[0]*b[0]; ,
693  d+= a[1]*b[1]; ,
694  d+= a[2]*b[2]; ,
695  d+= a[3]*b[3]; ,
696  d+= a[4]*b[4]; ,
697  d+= a[5]*b[5]; ,
698  d+= a[6]*b[6]; ,
699  d+= a[7]*b[7]; )
700  return d;
701 }
702 
703 inline double network::dotx(float* a, float* b) {
704  float d=0.;
705  NETX(d+= a[0]*b[0]; ,
706  d+= a[1]*b[1]; ,
707  d+= a[2]*b[2]; ,
708  d+= a[3]*b[3]; ,
709  d+= a[4]*b[4]; ,
710  d+= a[5]*b[5]; ,
711  d+= a[6]*b[6]; ,
712  d+= a[7]*b[7]; )
713  return d;
714 }
715 
716 inline double network::dot4(double* a, double* b) {
717  double d=0.;
718  d+= a[0]*b[0];
719  d+= a[1]*b[1];
720  d+= a[2]*b[2];
721  d+= a[3]*b[3];
722  return d;
723 }
724 
725 inline double network::dotx(double* a, double* b, double* c) {
726  double d=0.;
727  NETX(c[0] = a[0]*b[0]; d+=c[0]; ,
728  c[1] = a[1]*b[1]; d+=c[1]; ,
729  c[2] = a[2]*b[2]; d+=c[2]; ,
730  c[3] = a[3]*b[3]; d+=c[3]; ,
731  c[4] = a[4]*b[4]; d+=c[4]; ,
732  c[5] = a[5]*b[5]; d+=c[5]; ,
733  c[6] = a[6]*b[6]; d+=c[6]; ,
734  c[7] = a[7]*b[7]; d+=c[7]; )
735  return d;
736 }
737 
738 inline double network::dotx(float* a, float* b, float* c) {
739  float d=0.;
740  NETX(c[0] = a[0]*b[0]; d+=c[0]; ,
741  c[1] = a[1]*b[1]; d+=c[1]; ,
742  c[2] = a[2]*b[2]; d+=c[2]; ,
743  c[3] = a[3]*b[3]; d+=c[3]; ,
744  c[4] = a[4]*b[4]; d+=c[4]; ,
745  c[5] = a[5]*b[5]; d+=c[5]; ,
746  c[6] = a[6]*b[6]; d+=c[6]; ,
747  c[7] = a[7]*b[7]; d+=c[7]; )
748  return d;
749 }
750 
751 inline double network::dot4(double* a, double* b, double* c) {
752  double d=0.;
753  c[0] = a[0]*b[0]; d+=c[0];
754  c[1] = a[1]*b[1]; d+=c[1];
755  c[2] = a[2]*b[2]; d+=c[2];
756  c[3] = a[3]*b[3]; d+=c[3];
757  return d;
758 }
759 
760 inline double network::dotx(double* a, double** b, size_t j) {
761  double d=0.;
762  NETX(d+= a[0]*b[0][j]; ,
763  d+= a[1]*b[1][j]; ,
764  d+= a[2]*b[2][j]; ,
765  d+= a[3]*b[3][j]; ,
766  d+= a[4]*b[4][j]; ,
767  d+= a[5]*b[5][j]; ,
768  d+= a[6]*b[6][j]; ,
769  d+= a[7]*b[7][j]; )
770  return d;
771 }
772 
773 inline double network::dotx(double** a, size_t i, double* b) {
774  double d=0.;
775  NETX(d+= a[0][i]*b[0]; ,
776  d+= a[1][i]*b[1]; ,
777  d+= a[2][i]*b[2]; ,
778  d+= a[3][i]*b[3]; ,
779  d+= a[4][i]*b[4]; ,
780  d+= a[5][i]*b[5]; ,
781  d+= a[6][i]*b[6]; ,
782  d+= a[7][i]*b[7]; )
783  return d;
784 }
785 
786 inline double network::dotx(double** a, size_t i, double** b, size_t j) {
787  double d=0.;
788  NETX(d+= a[0][i]*b[0][j]; ,
789  d+= a[1][i]*b[1][j]; ,
790  d+= a[2][i]*b[2][j]; ,
791  d+= a[3][i]*b[3][j]; ,
792  d+= a[4][i]*b[4][j]; ,
793  d+= a[5][i]*b[5][j]; ,
794  d+= a[6][i]*b[6][j]; ,
795  d+= a[7][i]*b[7][j]; )
796  return d;
797 }
798 
799 inline double network::divx(double* a, double* b) {
800  double d=0.;
801  NETX(d+= a[0]/b[0]; ,
802  d+= a[1]/b[1]; ,
803  d+= a[2]/b[2]; ,
804  d+= a[3]/b[3]; ,
805  d+= a[4]/b[4]; ,
806  d+= a[5]/b[5]; ,
807  d+= a[6]/b[6]; ,
808  d+= a[7]/b[7]; )
809  return d;
810 }
811 
812 inline int network::netx(double* u, double um, double* v, double vm, double g) {
813  double d=0.;
814  double q = (1.-g)*um;
815  NETX(d+= int(u[0]*u[0]>q) - int((u[0]*u[0]/um+v[0]*v[0]/vm)>g); ,
816  d+= int(u[1]*u[1]>q) - int((u[1]*u[1]/um+v[1]*v[1]/vm)>g); ,
817  d+= int(u[2]*u[2]>q) - int((u[2]*u[2]/um+v[2]*v[2]/vm)>g); ,
818  d+= int(u[3]*u[3]>q) - int((u[3]*u[3]/um+v[3]*v[3]/vm)>g); ,
819  d+= int(u[4]*u[4]>q) - int((u[4]*u[4]/um+v[4]*v[4]/vm)>g); ,
820  d+= int(u[5]*u[5]>q) - int((u[5]*u[5]/um+v[5]*v[5]/vm)>g); ,
821  d+= int(u[6]*u[6]>q) - int((u[6]*u[6]/um+v[6]*v[6]/vm)>g); ,
822  d+= int(u[7]*u[7]>q) - int((u[7]*u[7]/um+v[7]*v[7]/vm)>g); )
823  return d;
824 }
825 
826 inline int network::netx(float* u, float um, float* v, float vm, float g) {
827  float d=0.;
828  float q = (1.-g)*um;
829  NETX(d+= int(u[0]*u[0]>q) - int((u[0]*u[0]/um+v[0]*v[0]/vm)>g); ,
830  d+= int(u[1]*u[1]>q) - int((u[1]*u[1]/um+v[1]*v[1]/vm)>g); ,
831  d+= int(u[2]*u[2]>q) - int((u[2]*u[2]/um+v[2]*v[2]/vm)>g); ,
832  d+= int(u[3]*u[3]>q) - int((u[3]*u[3]/um+v[3]*v[3]/vm)>g); ,
833  d+= int(u[4]*u[4]>q) - int((u[4]*u[4]/um+v[4]*v[4]/vm)>g); ,
834  d+= int(u[5]*u[5]>q) - int((u[5]*u[5]/um+v[5]*v[5]/vm)>g); ,
835  d+= int(u[6]*u[6]>q) - int((u[6]*u[6]/um+v[6]*v[6]/vm)>g); ,
836  d+= int(u[7]*u[7]>q) - int((u[7]*u[7]/um+v[7]*v[7]/vm)>g); )
837  return d;
838 }
839 
840 inline double network::dotx(double** a, size_t i, double** b, size_t j, double* p) {
841  double d=0.;
842  NETX(p[0] = a[0][i]*b[0][j];d+=p[0]; ,
843  p[1] = a[1][i]*b[1][j];d+=p[1]; ,
844  p[2] = a[2][i]*b[2][j];d+=p[2]; ,
845  p[3] = a[3][i]*b[3][j];d+=p[3]; ,
846  p[4] = a[4][i]*b[4][j];d+=p[4]; ,
847  p[5] = a[5][i]*b[5][j];d+=p[5]; ,
848  p[6] = a[6][i]*b[6][j];d+=p[6]; ,
849  p[7] = a[7][i]*b[7][j];d+=p[7]; )
850  return d;
851 }
852 
853 inline double network::dotx(double* a, double* b, double c) {
854  double d=0.;
855  NETX(d+= a[0]*b[0]; ,
856  d+= a[1]*b[1]; ,
857  d+= a[2]*b[2]; ,
858  d+= a[3]*b[3]; ,
859  d+= a[4]*b[4]; ,
860  d+= a[5]*b[5]; ,
861  d+= a[6]*b[6]; ,
862  d+= a[7]*b[7]; )
863  return c*d;
864 }
865 
866 inline double network::dotx(float* a, float* b, float c) {
867  float d=0.;
868  NETX(d+= a[0]*b[0]; ,
869  d+= a[1]*b[1]; ,
870  d+= a[2]*b[2]; ,
871  d+= a[3]*b[3]; ,
872  d+= a[4]*b[4]; ,
873  d+= a[5]*b[5]; ,
874  d+= a[6]*b[6]; ,
875  d+= a[7]*b[7]; )
876  return c*d;
877 }
878 
879 inline double network::dotx(double* a, double** b, size_t j, double* p) {
880  double d=0.;
881  NETX(p[0] = a[0]*b[0][j];d+=p[0]; ,
882  p[1] = a[1]*b[1][j];d+=p[1]; ,
883  p[2] = a[2]*b[2][j];d+=p[2]; ,
884  p[3] = a[3]*b[3][j];d+=p[3]; ,
885  p[4] = a[4]*b[4][j];d+=p[4]; ,
886  p[5] = a[5]*b[5][j];d+=p[5]; ,
887  p[6] = a[6]*b[6][j];d+=p[6]; ,
888  p[7] = a[7]*b[7][j];d+=p[7]; )
889  return d;
890 }
891 
892 inline void network::addx(double* a, double* b, double* p) {
893  NETX(p[0] += a[0]*b[0]; ,
894  p[1] += a[1]*b[1]; ,
895  p[2] += a[2]*b[2]; ,
896  p[3] += a[3]*b[3]; ,
897  p[4] += a[4]*b[4]; ,
898  p[5] += a[5]*b[5]; ,
899  p[6] += a[6]*b[6]; ,
900  p[7] += a[7]*b[7]; )
901  return;
902 }
903 
904 inline void network::addx(double* a, double** b, size_t j, double* p) {
905  NETX(p[0] += a[0]*b[0][j]; ,
906  p[1] += a[1]*b[1][j]; ,
907  p[2] += a[2]*b[2][j]; ,
908  p[3] += a[3]*b[3][j]; ,
909  p[4] += a[4]*b[4][j]; ,
910  p[5] += a[5]*b[5][j]; ,
911  p[6] += a[6]*b[6][j]; ,
912  p[7] += a[7]*b[7][j]; )
913  return;
914 }
915 
916 inline void network::addx(double** a, size_t i, double** b, size_t j, double* p) {
917  NETX(p[0] += a[0][i]*b[0][j]; ,
918  p[1] += a[1][i]*b[1][j]; ,
919  p[2] += a[2][i]*b[2][j]; ,
920  p[3] += a[3][i]*b[3][j]; ,
921  p[4] += a[4][i]*b[4][j]; ,
922  p[5] += a[5][i]*b[5][j]; ,
923  p[6] += a[6][i]*b[6][j]; ,
924  p[7] += a[7][i]*b[7][j]; )
925  return;
926 }
927 
928 inline void network::mulx(double** a, size_t i, double** b, size_t j, double* p) {
929  NETX(p[0] = a[0][i]*b[0][j]; ,
930  p[1] = a[1][i]*b[1][j]; ,
931  p[2] = a[2][i]*b[2][j]; ,
932  p[3] = a[3][i]*b[3][j]; ,
933  p[4] = a[4][i]*b[4][j]; ,
934  p[5] = a[5][i]*b[5][j]; ,
935  p[6] = a[6][i]*b[6][j]; ,
936  p[7] = a[7][i]*b[7][j]; )
937  return;
938 }
939 
940 inline void network::mulx(double* a, double b, double* p) {
941  NETX(p[0] = a[0]*b; ,
942  p[1] = a[1]*b; ,
943  p[2] = a[2]*b; ,
944  p[3] = a[3]*b; ,
945  p[4] = a[4]*b; ,
946  p[5] = a[5]*b; ,
947  p[6] = a[6]*b; ,
948  p[7] = a[7]*b; )
949  return;
950 }
951 
952 inline void network::mulx(float* a, float b, float* p) {
953  NETX(p[0] = a[0]*b; ,
954  p[1] = a[1]*b; ,
955  p[2] = a[2]*b; ,
956  p[3] = a[3]*b; ,
957  p[4] = a[4]*b; ,
958  p[5] = a[5]*b; ,
959  p[6] = a[6]*b; ,
960  p[7] = a[7]*b; )
961  return;
962 }
963 
964 inline void network::mulx(double* a, double b) {
965  NETX(a[0] *= b; ,
966  a[1] *= b; ,
967  a[2] *= b; ,
968  a[3] *= b; ,
969  a[4] *= b; ,
970  a[5] *= b; ,
971  a[6] *= b; ,
972  a[7] *= b; )
973  return;
974 }
975 
976 inline void network::mulx(float* a, float b) {
977  NETX(a[0] *= b; ,
978  a[1] *= b; ,
979  a[2] *= b; ,
980  a[3] *= b; ,
981  a[4] *= b; ,
982  a[5] *= b; ,
983  a[6] *= b; ,
984  a[7] *= b; )
985  return;
986 }
987 
988 inline void network::inix(double** a, size_t j, double* p) {
989  NETX(p[0] = a[0][j]; ,
990  p[1] = a[1][j]; ,
991  p[2] = a[2][j]; ,
992  p[3] = a[3][j]; ,
993  p[4] = a[4][j]; ,
994  p[5] = a[5][j]; ,
995  p[6] = a[6][j]; ,
996  p[7] = a[7][j]; )
997  return;
998 }
999 
1000 inline void network::inix(double* p, double a) {
1001  NETX(p[0] = a; ,
1002  p[1] = a; ,
1003  p[2] = a; ,
1004  p[3] = a; ,
1005  p[4] = a; ,
1006  p[5] = a; ,
1007  p[6] = a; ,
1008  p[7] = a; )
1009  return;
1010 }
1011 
1012 inline void network::inix(float* p, float a) {
1013  NETX(p[0] = a; ,
1014  p[1] = a; ,
1015  p[2] = a; ,
1016  p[3] = a; ,
1017  p[4] = a; ,
1018  p[5] = a; ,
1019  p[6] = a; ,
1020  p[7] = a; )
1021  return;
1022 }
1023 
1024 inline double network::rotx(double* u, double c, double* v, double s, double* e) {
1025  double d=0.;
1026  NETX(e[0] = u[0]*c + v[0]*s;d+=e[0]*e[0]; ,
1027  e[1] = u[1]*c + v[1]*s;d+=e[1]*e[1]; ,
1028  e[2] = u[2]*c + v[2]*s;d+=e[2]*e[2]; ,
1029  e[3] = u[3]*c + v[3]*s;d+=e[3]*e[3]; ,
1030  e[4] = u[4]*c + v[4]*s;d+=e[4]*e[4]; ,
1031  e[5] = u[5]*c + v[5]*s;d+=e[5]*e[5]; ,
1032  e[6] = u[6]*c + v[6]*s;d+=e[6]*e[6]; ,
1033  e[7] = u[7]*c + v[7]*s;d+=e[7]*e[7]; )
1034  return d;
1035 }
1036 
1037 inline double network::rotx(float* u, float c, float* v, float s, float* e) {
1038  double d=0.;
1039  NETX(e[0] = u[0]*c + v[0]*s;d+=e[0]*e[0]; ,
1040  e[1] = u[1]*c + v[1]*s;d+=e[1]*e[1]; ,
1041  e[2] = u[2]*c + v[2]*s;d+=e[2]*e[2]; ,
1042  e[3] = u[3]*c + v[3]*s;d+=e[3]*e[3]; ,
1043  e[4] = u[4]*c + v[4]*s;d+=e[4]*e[4]; ,
1044  e[5] = u[5]*c + v[5]*s;d+=e[5]*e[5]; ,
1045  e[6] = u[6]*c + v[6]*s;d+=e[6]*e[6]; ,
1046  e[7] = u[7]*c + v[7]*s;d+=e[7]*e[7]; )
1047  return d;
1048 }
1049 
1050 inline double network::rot4(double* u, double c, double* v, double s, double* e) {
1051  double d=0.;
1052  e[0] = u[0]*c + v[0]*s;d+=e[0]*e[0];
1053  e[1] = u[1]*c + v[1]*s;d+=e[1]*e[1];
1054  e[2] = u[2]*c + v[2]*s;d+=e[2]*e[2];
1055  e[3] = u[3]*c + v[3]*s;d+=e[3]*e[3];
1056  return d;
1057 }
1058 
1059 inline float network::rots(float* u, float c, float* v, float s, float* e) {
1060  float d=0.;
1061  NETX(e[0] -= u[0]*c + v[0]*s; d+=e[0]*e[0]; ,
1062  e[1] -= u[1]*c + v[1]*s; d+=e[1]*e[1]; ,
1063  e[2] -= u[2]*c + v[2]*s; d+=e[2]*e[2]; ,
1064  e[3] -= u[3]*c + v[3]*s; d+=e[3]*e[3]; ,
1065  e[4] -= u[4]*c + v[4]*s; d+=e[4]*e[4]; ,
1066  e[5] -= u[5]*c + v[5]*s; d+=e[5]*e[5]; ,
1067  e[6] -= u[6]*c + v[6]*s; d+=e[6]*e[6]; ,
1068  e[7] -= u[7]*c + v[7]*s; d+=e[7]*e[7]; )
1069  return d/2.;
1070 }
1071 
1072 
1073 inline void network::pnt_(float** q, float** p, short** m, int l, int n) {
1074 // point 0-7 float pointers to first network pixel
1075  NETX(q[0] = (p[0] + m[0][l]*n);,
1076  q[1] = (p[1] + m[1][l]*n);,
1077  q[2] = (p[2] + m[2][l]*n);,
1078  q[3] = (p[3] + m[3][l]*n);,
1079  q[4] = (p[4] + m[4][l]*n);,
1080  q[5] = (p[5] + m[5][l]*n);,
1081  q[6] = (p[6] + m[6][l]*n);,
1082  q[7] = (p[7] + m[7][l]*n);)
1083  return;
1084 }
1085 
1086 inline void network::cpf_(float*& a, double** p) {
1087 // copy to a data defined by array of pointers p and increment target pointer
1088  for(int i=0;i<XIFO;i++) *(a++) = *p[i];
1089  a+=NIFO-XIFO;
1090  return;
1091 }
1092 
1093 inline void network::cpf_(float*& a, double** p, size_t i) {
1094 // copy to a data defined by array of pointers p and increment target pointer
1095  for(int k=0;k<XIFO;k++) *(a++) = p[k][i];
1096  a+=NIFO-XIFO;
1097  return;
1098 }
1099 
1100 inline void network::cpp_(float*& a, float** p) {
1101 // copy to a data defined by array of pointers p and increment pointers
1102  for(int i=0;i<XIFO;i++) *(a++) = *p[i]++;
1103  a+=NIFO-XIFO;
1104  return;
1105 }
1106 
1107 inline double network::dot32(std::vector<float>* F, double* p, std::vector<short>* J) {
1108  return (*F)[0] *p[(*J)[0]] + (*F)[1] *p[(*J)[1]]
1109  + (*F)[2] *p[(*J)[2]] + (*F)[3] *p[(*J)[3]]
1110  + (*F)[4] *p[(*J)[4]] + (*F)[5] *p[(*J)[5]]
1111  + (*F)[6] *p[(*J)[6]] + (*F)[7] *p[(*J)[7]]
1112  + (*F)[8] *p[(*J)[8]] + (*F)[9] *p[(*J)[9]]
1113  + (*F)[10]*p[(*J)[10]] + (*F)[11]*p[(*J)[11]]
1114  + (*F)[12]*p[(*J)[12]] + (*F)[13]*p[(*J)[13]]
1115  + (*F)[14]*p[(*J)[14]] + (*F)[15]*p[(*J)[15]]
1116  + (*F)[16]*p[(*J)[16]] + (*F)[17]*p[(*J)[17]]
1117  + (*F)[18]*p[(*J)[18]] + (*F)[19]*p[(*J)[19]]
1118  + (*F)[20]*p[(*J)[20]] + (*F)[21]*p[(*J)[21]]
1119  + (*F)[22]*p[(*J)[22]] + (*F)[23]*p[(*J)[23]]
1120  + (*F)[24]*p[(*J)[24]] + (*F)[25]*p[(*J)[25]]
1121  + (*F)[26]*p[(*J)[26]] + (*F)[27]*p[(*J)[27]]
1122  + (*F)[28]*p[(*J)[28]] + (*F)[29]*p[(*J)[29]]
1123  + (*F)[30]*p[(*J)[30]] + (*F)[31]*p[(*J)[31]];
1124 }
1125 
1126 inline double network::dot32(double* F, double* p, int* J) {
1127  return F[0] *p[J[0]] + F[1] *p[J[1]] + F[2] *p[J[2]] + F[3] *p[J[3]]
1128  + F[4] *p[J[4]] + F[5] *p[J[5]] + F[6] *p[J[6]] + F[7] *p[J[7]]
1129  + F[8] *p[J[8]] + F[9] *p[J[9]] + F[10]*p[J[10]] + F[11]*p[J[11]]
1130  + F[12]*p[J[12]] + F[13]*p[J[13]] + F[14]*p[J[14]] + F[15]*p[J[15]]
1131  + F[16]*p[J[16]] + F[17]*p[J[17]] + F[18]*p[J[18]] + F[19]*p[J[19]]
1132  + F[20]*p[J[20]] + F[21]*p[J[21]] + F[22]*p[J[22]] + F[23]*p[J[23]]
1133  + F[24]*p[J[24]] + F[25]*p[J[25]] + F[26]*p[J[26]] + F[27]*p[J[27]]
1134  + F[28]*p[J[28]] + F[29]*p[J[29]] + F[30]*p[J[30]] + F[31]*p[J[31]];
1135 }
1136 
1137 inline void network::dpfx(float* fp, float* fx) {
1138 // transformation to DPF for 4 consecutive pixels.
1139 // rotate vectors fp and fx into DPF
1140 
1141 #ifndef __CINT__
1142  float t[32],T[32];
1143 
1144  __m128* _t = (__m128*) t; __m128* _T = (__m128*) T;
1145  __m128* _fp = (__m128*) fp; __m128* _fx = (__m128*) fx;
1146 
1147  _sse_dpf4_ps(_fp,_fx,_t,_T);
1148  _sse_cpf4_ps(_fp,_t);
1149  _sse_cpf4_ps(_fx,_T);
1150 #endif
1151 }
1152 
1153 inline void network::pnpx(float* fp, float* fx, float* am, float* AM, float* u, float* v) {
1154 // projection to network plane (pnp)
1155 // project vectors am and AM on the network plane fp,fx
1156 
1157 #ifndef __CINT__
1158  __m128* _fp = (__m128*) fp; __m128* _fx = (__m128*) fx;
1159  __m128* _am = (__m128*) am; __m128* _AM = (__m128*) AM;
1160  __m128* _u = (__m128*) u; __m128* _v = (__m128*) v;
1161 
1162  _sse_pnp4_ps(_fp,_fx,_am,_AM,_u,_v);
1163 #endif
1164 }
1165 
1166 inline void network::dspx(float* u, float* v, float* am, float* AM) {
1167 // dual stream phase (dsp) transformation
1168 // take projection vectors u and v,
1169 // make them orthogonal,
1170 // apply phase transformation boths to data a,A and projections u,v
1171 
1172 #ifndef __CINT__
1173  float U[32],V[32];
1174 
1175  __m128* _am = (__m128*) am; __m128* _AM = (__m128*) AM;
1176  __m128* _u = (__m128*) u; __m128* _v = (__m128*) v;
1177  __m128* _U = (__m128*) U; __m128* _V = (__m128*) V;
1178 
1179  _sse_dsp4_ps(_u,_v,_am,_AM,_U,_V);
1180  _sse_cpf4_ps(_u,_U);
1181  _sse_cpf4_ps(_v,_V);
1182 #endif
1183 }
1184 
1185 inline void network::dspx(float* fp, float* fx, float* am, float* AM, float* u, float* v) {
1186 // DPF + DSP transformations
1187 // applied dpfx+pnpx+dspx
1188 
1189  dpfx(fp, fx);
1190  pnpx(fp, fx, am, AM, u, v);
1191  dspx(u, v, am, AM);
1192 }
1193 
1194 inline int network::_sse_MRA_ps(float* amp, float* AMP, float Eo, int K) {
1195 // fast multi-resolution analysis inside sky loop
1196 // select max E pixel and either scale or skip it based on the value of residual
1197 // pointer to 00 phase amplitude of monster pixels
1198 // pointer to 90 phase amplitude of monster pixels
1199 // Eo - energy threshold
1200 // K - number of principle components to extract
1201 // returns number of MRA pixels
1202 
1203 #ifndef __CINT__
1204  int j,n,mm;
1205  int k = 0;
1206  int m = 0;
1207  int f = NIFO/4;
1208  int V = (int)this->rNRG.size();
1209  float* ee = this->rNRG.data; // residual energy
1210  float* pp = this->pNRG.data; // residual energy
1211  float EE = 0.; // extracted energy
1212  float E;
1213  float mam[NIFO] _ALIGNED;
1214  float mAM[NIFO] _ALIGNED;
1215  this->pNRG=-1;
1216  for(j=0; j<V; ++j) if(ee[j]>Eo) pp[j]=0;
1217 
1218  __m128* _m00 = (__m128*) mam;
1219  __m128* _m90 = (__m128*) mAM;
1220  __m128* _amp = (__m128*) amp;
1221  __m128* _AMP = (__m128*) AMP;
1222  __m128* _a00 = (__m128*) a_00.data;
1223  __m128* _a90 = (__m128*) a_90.data;
1224 
1225  while(k<K){
1226 
1227  for(j=0; j<V; ++j) if(ee[j]>ee[m]) m=j; // find max pixel
1228  if(ee[m]<=Eo) break; mm = m*f;
1229 
1230  E = _sse_abs_ps(_a00+mm,_a90+mm); EE += E; // get PC energy
1231  int J = wdmMRA.size(m);
1232  float* c = wdmMRA.getXTalk(m); // c1*c2+c3*c4=c1*c3+c2*c4=0
1233 
1234  if(E/EE < 0.01) break; // ignore small PC
1235 
1236  _sse_cpf_ps(mam,_a00+mm); // store a00 for max pixel
1237  _sse_cpf_ps(mAM,_a90+mm); // store a90 for max pixel
1238  _sse_add_ps(_amp+mm,_m00); // update 00 PC
1239  _sse_add_ps(_AMP+mm,_m90); // update 90 PC
1240 
1241  for(j=0; j<J; j++) {
1242  n = int(c[0]+0.1);
1243  if(ee[n]>Eo) {
1244  ee[n] = _sse_rotsub_ps(_m00,c[4],_m90,c[5],_a00+n*f); // subtract PC from a00
1245  ee[n]+= _sse_rotsub_ps(_m00,c[6],_m90,c[7],_a90+n*f); // subtract PC from a90
1246  }
1247  c += 8;
1248  }
1249  pp[m] = _sse_abs_ps(_amp+mm,_AMP+mm); // store PC energy
1250  k++;
1251  }
1252  return k;
1253 #else
1254  return 0;
1255 #endif
1256 }
1257 
1258 inline int network::_sse_mra_ps(float* amp, float* AMP, float Eo, int K) {
1259 // fast multi-resolution analysis inside sky loop
1260 // select max E pixel and either scale or skip it based on the value of residual
1261 // pointer to 00 phase amplitude of monster pixels
1262 // pointer to 90 phase amplitude of monster pixels
1263 // Eo - energy threshold
1264 // K - max number of principle components to extract
1265 // returns number of MRA pixels
1266 
1267 #ifndef __CINT__
1268  int j,n,mm,J;
1269  int k = 0;
1270  int m = 0;
1271  int f = NIFO/4;
1272  int V = (int)this->rNRG.size();
1273  float* ee = this->rNRG.data; // residual energy
1274  float* pp = this->pNRG.data; // PC energy
1275  float* c = NULL;
1276  float E2 = Eo/2; // threshold
1277  float E;
1278  float mam[NIFO] _ALIGNED;
1279  float mAM[NIFO] _ALIGNED;
1280 
1281  __m128* _m00 = (__m128*) mam;
1282  __m128* _m90 = (__m128*) mAM;
1283  __m128* _amp = (__m128*) amp;
1284  __m128* _AMP = (__m128*) AMP;
1285  __m128* _a00 = (__m128*) a_00.data;
1286  __m128* _a90 = (__m128*) a_90.data;
1287 
1288  this->pNRG = 0;
1289 
1290  while(k<K){
1291 
1292  for(j=0; j<V; ++j) if(ee[j]>ee[m]) m=j; // find max pixel
1293  if(ee[m]<=Eo) break; mm = m*f;
1294 
1295  //cout<<k<<" "<<" V= "<<V<<" m="<<m<<" ee[m]="<<ee[m];
1296 
1297  float cc = 0.;
1298 
1299  _sse_zero_ps(_m00);
1300  _sse_zero_ps(_m90);
1301 
1302  J = wdmMRA.size(m);
1303  c = wdmMRA.getXTalk(m); // c1*c2+c3*c4=c1*c3+c2*c4=0
1304  for(j=0; j<J; j++) {
1305  n = int(c[0]+0.1);
1306  if(ee[n]>Eo) {
1307  _sse_rotadd_ps(_a00+n*f,c[4],_a90+n*f,c[6],_m00); // construct 00 vector
1308  _sse_rotadd_ps(_a00+n*f,c[5],_a90+n*f,c[7],_m90); // construct 90 vector
1309  }
1310  if(ee[n]>0) cc += c[1];
1311  c += 8;
1312  }
1313  _sse_mul_ps(_m00,cc>1?1./cc:0.7);
1314  _sse_mul_ps(_m90,cc>1?1./cc:0.7);
1315  E = _sse_abs_ps(_m00,_m90); // get PC energy
1316 
1317  if(E > ee[m]) { // correct overtuning PC
1318  _sse_cpf_ps(mam,_a00+mm); // store a00 for max pixel
1319  _sse_cpf_ps(mAM,_a90+mm); // store a90 for max pixel
1320  }
1321  _sse_add_ps(_amp+mm,_m00); // update 00 PC
1322  _sse_add_ps(_AMP+mm,_m90); // update 90 PC
1323 
1324  J = wdmMRA.size(m);
1325  c = wdmMRA.getXTalk(m); // c1*c2+c3*c4=c1*c3+c2*c4=0
1326  for(j=0; j<J; j++) {
1327  n = int(c[0]+0.1);
1328  if(E<E2 && n!=m) {c+=8; continue;}
1329  if(ee[n]>Eo) {
1330  ee[n] = _sse_rotsub_ps(_m00,c[4],_m90,c[5],_a00+n*f); // subtract PC from a00
1331  ee[n]+= _sse_rotsub_ps(_m00,c[6],_m90,c[7],_a90+n*f); // subtract PC from a90
1332  ee[n]+= 1.e-6;
1333  }
1334  c += 8;
1335  }
1336  pp[m] = _sse_abs_ps(_amp+mm,_AMP+mm)+E2/2; // store PC energy
1337  //cout<<" "<<ee[m]<<" "<<k<<" "<<E<<" "<<EE<<" "<<endl;
1338  //cout<<" "<<m<<" "<<cc<<" "<<ee[m]<<" "<<E<<" "<<pp[m]<<endl;
1339  k++;
1340  }
1341  return k;
1342 #else
1343  return 0;
1344 #endif
1345 }
1346 
1347 inline wavearray<float> network::_avx_norm_ps(float** p, float** q,
1348  std::vector<float*> &pAVX, int I) {
1349  wavearray<float> norm(NIFO+1); // output array for packet norms
1350  float* g = norm.data+1; norm=0.;
1351 
1352 #ifndef __CINT__
1353 
1354 // return packet norm for each detector
1355  int i,j,k,n,m;
1356  int M = abs(I);
1357  int II = abs(I*2);
1358  float o = 1.e-12;
1359  float* mk = pAVX[1]; // pixel energy mask
1360  float* rn = pAVX[22]; // halo noise
1361  wavearray<float> tmp(NIFO); // array to store data
1362  float* t = tmp.data; tmp=0.;
1363  float e,u,v;
1364 
1365  float am[4*8] _ALIGNED; __m128* _am = (__m128*)am;
1366  float x[4*8] _ALIGNED; __m128* _x = (__m128*)x;
1367  float h[4*8] _ALIGNED; __m128* _h = (__m128*)h; // halo
1368 
1369  float* an = pAVX[17]; // M*4*NIFO array
1370  NETX(__m128* _a0 = (__m128*)(an+4*M*0);, __m128* _a1 = (__m128*)(an+4*M*1);,
1371  __m128* _a2 = (__m128*)(an+4*M*2);, __m128* _a3 = (__m128*)(an+4*M*3);,
1372  __m128* _a4 = (__m128*)(an+4*M*4);, __m128* _a5 = (__m128*)(an+4*M*5);,
1373  __m128* _a6 = (__m128*)(an+4*M*6);, __m128* _a7 = (__m128*)(an+4*M*7);)
1374 
1375  for(m=0; m<M; m++) {
1376  if(I>0) rn[m] = 0.;
1377  NETX(_a0[m] = _mm_set_ps(q[0][m],q[0][m],p[0][m],p[0][m]); q[0][M+m]=0; ,
1378  _a1[m] = _mm_set_ps(q[1][m],q[1][m],p[1][m],p[1][m]); q[1][M+m]=0; ,
1379  _a2[m] = _mm_set_ps(q[2][m],q[2][m],p[2][m],p[2][m]); q[2][M+m]=0; ,
1380  _a3[m] = _mm_set_ps(q[3][m],q[3][m],p[3][m],p[3][m]); q[3][M+m]=0; ,
1381  _a4[m] = _mm_set_ps(q[4][m],q[4][m],p[4][m],p[4][m]); q[4][M+m]=0; ,
1382  _a5[m] = _mm_set_ps(q[5][m],q[5][m],p[5][m],p[5][m]); q[5][M+m]=0; ,
1383  _a6[m] = _mm_set_ps(q[6][m],q[6][m],p[6][m],p[6][m]); q[6][M+m]=0; ,
1384  _a7[m] = _mm_set_ps(q[7][m],q[7][m],p[7][m],p[7][m]); q[7][M+m]=0; )
1385  }
1386 
1387  for(m=0; m<M; m++) {
1388  if(mk[m]<=0.) continue;
1389 
1390  int J = wdmMRA.size(m)*2;
1391  float cc = 0;
1392  float* c = wdmMRA.getXTalk(m);
1393  __m128* _c = (__m128*)(c+4);
1394 
1395  NETX(u=p[0][m]; v=q[0][m]; _am[0]=_mm_set_ps(v,u,v,u); _x[0]=_mm_setzero_ps(); ,
1396  u=p[1][m]; v=q[1][m]; _am[1]=_mm_set_ps(v,u,v,u); _x[1]=_mm_setzero_ps(); ,
1397  u=p[2][m]; v=q[2][m]; _am[2]=_mm_set_ps(v,u,v,u); _x[2]=_mm_setzero_ps(); ,
1398  u=p[3][m]; v=q[3][m]; _am[3]=_mm_set_ps(v,u,v,u); _x[3]=_mm_setzero_ps(); ,
1399  u=p[4][m]; v=q[4][m]; _am[4]=_mm_set_ps(v,u,v,u); _x[4]=_mm_setzero_ps(); ,
1400  u=p[5][m]; v=q[5][m]; _am[5]=_mm_set_ps(v,u,v,u); _x[5]=_mm_setzero_ps(); ,
1401  u=p[6][m]; v=q[6][m]; _am[6]=_mm_set_ps(v,u,v,u); _x[6]=_mm_setzero_ps(); ,
1402  u=p[7][m]; v=q[7][m]; _am[7]=_mm_set_ps(v,u,v,u); _x[7]=_mm_setzero_ps(); )
1403 
1404  for(j=0; j<J; j+=2) {
1405  n = int(c[j*4]);
1406  NETX(_x[0]=_mm_add_ps(_x[0],_mm_mul_ps(_c[j],_a0[n]));,
1407  _x[1]=_mm_add_ps(_x[1],_mm_mul_ps(_c[j],_a1[n]));,
1408  _x[2]=_mm_add_ps(_x[2],_mm_mul_ps(_c[j],_a2[n]));,
1409  _x[3]=_mm_add_ps(_x[3],_mm_mul_ps(_c[j],_a3[n]));,
1410  _x[4]=_mm_add_ps(_x[4],_mm_mul_ps(_c[j],_a4[n]));,
1411  _x[5]=_mm_add_ps(_x[5],_mm_mul_ps(_c[j],_a5[n]));,
1412  _x[6]=_mm_add_ps(_x[6],_mm_mul_ps(_c[j],_a6[n]));,
1413  _x[7]=_mm_add_ps(_x[7],_mm_mul_ps(_c[j],_a7[n]));)
1414  }
1415 
1416  NETX(_h[0]=_mm_mul_ps(_x[0],_am[0]);,
1417  _h[1]=_mm_mul_ps(_x[1],_am[1]);,
1418  _h[2]=_mm_mul_ps(_x[2],_am[2]);,
1419  _h[3]=_mm_mul_ps(_x[3],_am[3]);,
1420  _h[4]=_mm_mul_ps(_x[4],_am[4]);,
1421  _h[5]=_mm_mul_ps(_x[5],_am[5]);,
1422  _h[6]=_mm_mul_ps(_x[6],_am[6]);,
1423  _h[7]=_mm_mul_ps(_x[7],_am[7]);)
1424 
1425  NETX(t[0]=h[ 0]+h[ 1]+h[ 2]+h[ 3]; t[0]=t[0]>0?t[0]:0; g[0]+=t[0]; ,
1426  t[1]=h[ 4]+h[ 5]+h[ 6]+h[ 7]; t[1]=t[1]>0?t[1]:0; g[1]+=t[1]; ,
1427  t[2]=h[ 8]+h[ 9]+h[10]+h[11]; t[2]=t[2]>0?t[2]:0; g[2]+=t[2]; ,
1428  t[3]=h[12]+h[13]+h[14]+h[15]; t[3]=t[3]>0?t[3]:0; g[3]+=t[3]; ,
1429  t[4]=h[16]+h[17]+h[18]+h[19]; t[4]=t[4]>0?t[4]:0; g[4]+=t[4]; ,
1430  t[5]=h[20]+h[21]+h[22]+h[23]; t[5]=t[5]>0?t[5]:0; g[5]+=t[5]; ,
1431  t[6]=h[24]+h[25]+h[26]+h[27]; t[6]=t[6]>0?t[6]:0; g[6]+=t[6]; ,
1432  t[7]=h[28]+h[29]+h[30]+h[31]; t[7]=t[7]>0?t[7]:0; g[7]+=t[7]; )
1433 
1434  if(I<0) continue;
1435 
1436  NETX(u=p[0][m]; v=q[0][m]; e=(u*u+v*v)/(t[0]+o); q[0][M+m]=(e>=1)?0:e; ,
1437  u=p[1][m]; v=q[1][m]; e=(u*u+v*v)/(t[1]+o); q[1][M+m]=(e>=1)?0:e; ,
1438  u=p[2][m]; v=q[2][m]; e=(u*u+v*v)/(t[2]+o); q[2][M+m]=(e>=1)?0:e; ,
1439  u=p[3][m]; v=q[3][m]; e=(u*u+v*v)/(t[3]+o); q[3][M+m]=(e>=1)?0:e; ,
1440  u=p[4][m]; v=q[4][m]; e=(u*u+v*v)/(t[4]+o); q[4][M+m]=(e>=1)?0:e; ,
1441  u=p[5][m]; v=q[5][m]; e=(u*u+v*v)/(t[5]+o); q[5][M+m]=(e>=1)?0:e; ,
1442  u=p[6][m]; v=q[6][m]; e=(u*u+v*v)/(t[6]+o); q[6][M+m]=(e>=1)?0:e; ,
1443  u=p[7][m]; v=q[7][m]; e=(u*u+v*v)/(t[7]+o); q[7][M+m]=(e>=1)?0:e; )
1444 
1445  NETX(u=x[ 0]+x[ 2]; v=x[ 1]+x[ 3]; rn[m]+=u*u+v*v; ,
1446  u=x[ 4]+x[ 6]; v=x[ 5]+x[ 7]; rn[m]+=u*u+v*v; ,
1447  u=x[ 8]+x[10]; v=x[ 9]+x[11]; rn[m]+=u*u+v*v; ,
1448  u=x[12]+x[14]; v=x[13]+x[15]; rn[m]+=u*u+v*v; ,
1449  u=x[16]+x[18]; v=x[17]+x[19]; rn[m]+=u*u+v*v; ,
1450  u=x[20]+x[22]; v=x[21]+x[23]; rn[m]+=u*u+v*v; ,
1451  u=x[24]+x[26]; v=x[25]+x[27]; rn[m]+=u*u+v*v; ,
1452  u=x[28]+x[30]; v=x[29]+x[31]; rn[m]+=u*u+v*v; )
1453 
1454  }
1455 
1456  for(n=1; n<=XIFO; n++) { // save norms
1457  if(I>0) {
1458  e = q[n-1][II+4]*q[n-1][II+4]; // TF-Domain SNR
1459  if(norm.data[n]<2.) norm.data[n]=2;
1460  q[n-1][II+5] = norm.data[n]; // save norms
1461  norm.data[n] = e/norm.data[n]; // detector {1:NIFO} SNR
1462  }
1463  norm.data[0] += norm.data[n]; // total SNR
1464  }
1465 #endif
1466  return norm;
1467 }
1468 
1469 inline wavearray<float> network::_avx_norm_ps(float** p, float** q, float* ec, int I) {
1470 // use GW norm from data packet
1471 // p - GW array
1472 // q - Data array
1473  float e;
1474  int II = abs(I*2);
1475  wavearray<float> norm(NIFO+1); // array for packet norms
1476  float* nn = norm.data; // array for packet norms
1477  norm = 0;
1478  for(int n=1; n<=XIFO; n++) { // save norms
1479  nn[n] = q[n-1][II+5]; // get data norms
1480  p[n-1][II+5] = nn[n]; // save norms
1481  e = p[n-1][II+4]*p[n-1][II+4]; // TF-Domain SNR
1482  nn[n] = e/nn[n]; // detector {1:NIFO} SNR
1483  nn[0] += nn[n]; // total SNR
1484  for(int i=0; i<I; i++) { // save signal norms
1485  p[n-1][I+i] = ec[i]>0 ? q[n-1][I+i] : 0.; // set signal norms
1486  }
1487  }
1488  return norm;
1489 }
1490 
1491 inline void network::_avx_saveGW_ps(float** p, float** q, int I) {
1492 // save GW strain amplitudes into a_00, a_90 arrays
1493 // p,q - input - GW warray
1494 // I - number of GW pixels
1495 // in likelihoodWP these arrays should be stored exactly in the same order.
1496 
1497  for(int i=0; i<I; i++) {
1498  for(int n=0; n<NIFO; n++) {
1499  a_00[i*NIFO+n]=p[n][i];
1500  a_90[i*NIFO+n]=q[n][i];
1501  }
1502  }
1503  return;
1504 }
1505 
1506 #endif // NETWORK_HH
wavearray< double > t(hp.size())
monster wdmMRA
list of pixel pointers for MRA
Definition: network.hh:633
std::vector< char * > ifoName
Definition: network.hh:591
size_t mdcTypeSize()
Definition: network.hh:392
wavearray< short > veto
Definition: network.hh:626
detector * getifo(size_t n)
param: detector index
Definition: network.hh:418
static float _sse_abs_ps(__m128 *_a)
Definition: watsse.hh:119
double netCC
Definition: network.hh:578
#define NIFO
Definition: wat.hh:56
size_t nLag
Definition: network.hh:555
double start
Definition: network.hh:37
static double g(double e)
Definition: GNGen.cc:98
wavearray< double > skyHole
Definition: network.hh:625
void like(char _LIKE)
buffer for standard response 90 ampl
Definition: network.hh:648
WSeries< double > pixeLHood
Definition: network.hh:616
std::vector< netcluster > wc_List
Definition: network.hh:592
tuple f
Definition: cwb_online.py:91
void set2or(double p)
param: threshold
Definition: network.hh:453
std::vector< double > * getmdcTime()
Definition: network.hh:404
static void cpp_(float *&, float **)
Definition: network.hh:1100
double delta
static float rots(float *, float, float *, float, float *)
Definition: network.hh:1059
std::vector< pixel > cluster
Definition: wavepath.hh:43
double THRESHOLD
#define _ALIGNED
Definition: wat.hh:53
wavearray< double > a(hp.size())
size_t setSkyMask(double, char *, bool *, int)
Definition: CreateSkyMask.C:70
WSeries< float > v[nIFO]
Definition: cwb_net.C:62
static void mulx(double **, size_t, double **, size_t, double *)
Definition: network.hh:928
static void inix(double **, size_t, double *)
Definition: network.hh:988
void wdm(bool _WDM)
Definition: network.hh:649
double Gamma(double r)
Definition: watfun.hh:175
static void pnpx(float *fp, float *fx, float *am, float *AM, float *u, float *v)
Definition: network.hh:1153
int n
Definition: cwb_net.C:10
static void _sse_zero_ps(__m128 *_p)
Definition: watsse.hh:26
skymap nSkyStat
Definition: network.hh:610
std::vector< vectorD > NDM
Definition: network.hh:585
wavearray< float > a_90
buffer for cluster sky 00 amplitude
Definition: network.hh:635
double iGamma(double r, double p)
Definition: watfun.hh:187
static void _sse_dsp4_ps(__m128 *u, __m128 *v, __m128 *_am, __m128 *_AM, __m128 *_u, __m128 *_v)
Definition: watsse.hh:644
size_t nRun
Definition: network.hh:554
double getliveTime(size_t n)
Definition: network.hh:410
int ID
Definition: TestMDC.C:70
double getmdcTime(size_t n)
Definition: network.hh:406
skymap nPenalty
Definition: network.hh:606
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 bpp
Definition: test_config1.C:22
static void _sse_cpf4_ps(__m128 *_aa, __m128 *_pp)
Definition: watsse.hh:299
int _sse_mra_ps(float *, float *, float, int)
Definition: network.hh:1258
float theta
WDM< double > wdm(nLAYERS, nLAYERS, 6, 10)
bool eDisbalance
Definition: network.hh:580
std::vector< netpixel * > pList
Definition: network.hh:632
string getmdcType(size_t n)
Definition: network.hh:402
skymap nAntenaPrior
Definition: network.hh:614
static float _sse_rotsub_ps(__m128 *_u, float c, __m128 *_v, float s, __m128 *_a)
Definition: watsse.hh:381
void _avx_saveGW_ps(float **, float **, int)
Definition: network.hh:1491
std::vector< std::string > mdcType
Definition: network.hh:595
size_t mIFO
Definition: network.hh:557
STL namespace.
skymap nPolarisation
Definition: network.hh:612
string getmdcList(size_t n)
Definition: network.hh:400
static void _sse_cpf_ps(float *a, __m128 *_p)
Definition: watsse.hh:289
#define M
Definition: UniqSLagsList.C:3
int m
Definition: cwb_net.C:10
std::vector< size_t > mdc__ID
Definition: network.hh:597
int j
Definition: cwb_net.C:10
static double rotx(double *, double, double *, double, double *)
Definition: network.hh:1024
i drho i
std::vector< double > mdcTime
Definition: network.hh:596
size_t cluster(int kt, int kf)
param: time gap in pixels return: number of reconstructed clusters
Definition: network.hh:303
void setRunID(size_t n)
param: run
Definition: network.hh:437
std::vector< double > vectorD
Definition: network.hh:33
std::vector< detector * > ifoList
Definition: network.hh:590
static void _sse_rotadd_ps(__m128 *_u, float c, __m128 *_v, float s, __m128 *_a)
Definition: watsse.hh:370
long subNetCut(network *net, int lag, float snc, TH2F *hist)
SUPERCLUSTER.
size_t mdcListSize()
Definition: network.hh:390
size_t mdcTimeSize()
Definition: network.hh:394
wavearray< double > skyENRG
Definition: network.hh:628
size_t wdmListSize()
Definition: network.hh:428
skymap nCorrEnergy
Definition: network.hh:607
double Edge
Definition: network.hh:560
size_t events(int type, int lag=-1)
Definition: network.hh:323
virtual void Browse(TBrowser *)
Definition: network.hh:498
size_t ifoListSize()
Definition: network.hh:413
void constraint(double d=1., double g=0.0001)
param: constraint parameter, p=0 - no constraint
Definition: network.hh:445
WDM< double > * getwdm(size_t M)
param: number of wdm layers
Definition: network.hh:430
static int netx(double *, double, double *, double, double)
Definition: network.hh:812
double aNET
Definition: network.hh:562
skymap nLikelihood
Definition: network.hh:604
bool optim
Definition: network.hh:572
size_t mode
wavearray< float > pNRG
buffers for cluster residual energy
Definition: network.hh:637
float phi
double precision
Definition: network.hh:575
std::vector< delayFilter > filter90
Definition: network.hh:620
static double divx(double *, double *)
Definition: network.hh:799
void Forward(size_t k)
param: number of steps
Definition: network.hh:71
int _sse_MRA_ps(network *net, float *amp, float *AMP, float Eo, int K)
static void dspx(float *u, float *v, float *am, float *AM)
Definition: network.hh:1166
std::vector< double > livTime
Definition: network.hh:593
wavearray< double > h
Definition: Regression_H1.C:25
static void addx(double *, double *, double *)
Definition: network.hh:892
bool local
Definition: network.hh:571
r setFilter(10)
i() int(T_cor *100))
std::vector< std::string > mdcList
Definition: network.hh:594
static void dpfx(float *fp, float *fx)
Definition: network.hh:1137
double * tmp
Definition: testWDM_5.C:31
double eCOR
Definition: network.hh:564
void setRMS()
Definition: network.hh:661
wavearray< float > a_00
wdm multi-resolution analysis
Definition: network.hh:634
std::vector< WDM< double > * > wdmList
Definition: network.hh:599
double gNET
Definition: network.hh:561
char fname[1024]
static void cpf_(float *&a, double **p)
Definition: network.hh:1086
size_t mdc__IDSize()
Definition: network.hh:396
double Step
Definition: network.hh:559
void delink()
Definition: network.hh:336
double acor
static void pnt_(float **, float **, short **, int, int)
Definition: network.hh:1073
int k
int pattern
Definition: network.hh:583
void setOffset(double t)
param: run
Definition: network.hh:441
void Inverse(size_t k)
Definition: network.hh:73
std::vector< waveSegment > segList
Definition: network.hh:598
int _sse_mra_ps(network *NET, float *amp, float *AMP, float Eo, int K)
r add(tfmap,"hchannel")
wavearray< int > index
Definition: network.hh:622
double acor
Definition: network.hh:567
static double dot32(std::vector< float > *, double *, std::vector< short > *)
Definition: network.hh:1107
double F
int getIndex(double theta, double phi)
param: theta [deg] param: phi [deg]
Definition: network.hh:189
Definition: skymap.hh:45
double e2or
Definition: network.hh:566
double e
wavearray< double > skyProb
Definition: network.hh:627
void setAcore(double a)
Definition: network.hh:383
WSeries< double > whp
Definition: network.hh:587
char like()
Definition: network.hh:494
size_t events()
Definition: network.hh:311
s s
Definition: cwb_net.C:137
static double rot4(double *, double, double *, double, double *)
Definition: network.hh:1050
int index
Definition: network.hh:36
bool _WDM
Definition: network.hh:650
std::vector< delayFilter > filter
Definition: network.hh:619
bool wfsave
Definition: network.hh:582
WSeries< double > pixeLNull
Definition: network.hh:617
double pSigma
Definition: network.hh:576
skymap nAlignment
Definition: network.hh:602
double T
Definition: testWDM_4.C:11
static void _sse_add_ps(__m128 *_a, __m128 *_b)
Definition: watsse.hh:230
WSeries< double > whx
Definition: network.hh:588
static double dot4(double *, double *)
Definition: network.hh:716
skymap nDisbalance
Definition: network.hh:609
double getNDM(size_t i, size_t j)
param: first detector param: second detector
Definition: network.hh:122
wavearray< float > _avx_norm_ps(float **, float **, std::vector< float * > &, int)
Definition: network.hh:1347
size_t add(WDM< double > *wdm)
param: pointer to wdm return number of wdm tronsforms in the list
Definition: network.hh:426
bool pOUT
Definition: network.hh:568
double norm
Definition: network.hh:565
iD print()
size_t livTimeSize()
Definition: network.hh:398
double gamma
Definition: network.hh:574
wavearray< double > skyMaskCC
Definition: network.hh:624
netcluster * getwc(size_t n)
param: delay index
Definition: network.hh:421
int l
Definition: cbc_plots.C:434
skymap nNetIndex
Definition: network.hh:608
skymap nSensitivity
list of wdm tranformations
Definition: network.hh:601
int _sse_MRA_ps(float *, float *, float, int)
Definition: network.hh:1194
static double dotx(double *, double *)
Definition: network.hh:690
static void _sse_mul_ps(__m128 *_a, float b)
Definition: watsse.hh:38
DataType_t * data
Definition: wavearray.hh:301
wavearray< float > rNRG
buffer for cluster sky 90 amplitudes
Definition: network.hh:636
long nSky
Definition: network.hh:556
size_t wc_ListSize()
Definition: network.hh:415
void setMRAcatalog(char *fn)
Definition: network.hh:297
double rTDF
Definition: network.hh:558
skymap nNullEnergy
Definition: network.hh:605
wavearray< short > skyMask
Definition: network.hh:623
skymap nCorrelation
Definition: network.hh:603
size_t getmdc__ID(size_t n)
Definition: network.hh:408
bool MRA
Definition: network.hh:581
double penalty
Definition: network.hh:577
double stop
Definition: network.hh:38
double delta
Definition: network.hh:573
bool wdm()
Definition: network.hh:492
char _LIKE
Definition: network.hh:651
char tYPe
Definition: network.hh:570
double netRHO
Definition: network.hh:579
static void _sse_dpf4_ps(__m128 *_Fp, __m128 *_Fx, __m128 *_fp, __m128 *_fx)
Definition: watsse.hh:613
skymap nProbability
Definition: network.hh:613
static double sumx(double *)
Definition: network.hh:677
static void _sse_pnp4_ps(__m128 *_fp, __m128 *_fx, __m128 *_am, __m128 *_AM, __m128 *_u, __m128 *_v)
Definition: watsse.hh:623
bool EFEC
Definition: network.hh:569
double iNET
Definition: network.hh:563
exit(0)
skymap nEllipticity
Definition: network.hh:611