Logo Coherent WaveBurst  
Reference Guide
Logo
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
regression.cc
Go to the documentation of this file.
1 //---------------------------------------------------
2 // WAT class for regression analysis
3 // S. Klimenko, University of Florida
4 //---------------------------------------------------
5 
6 #include <time.h>
7 #include <iostream>
8 #include <stdexcept>
9 #include "regression.hh"
10 #include "Biorthogonal.hh"
11 #include "SymmArray.hh"
12 #include "WDM.hh"
13 
14 ClassImp(regression) // used by THtml doc
15 
16 regression::regression() : Edge(0.), pOUT(false)
17 {
18 //
19 // Default constructor
20 //
21 
22  this->clear();
23 }
24 
25 regression::regression(WSeries<double> &in, char* ch, double fL, double fH)
26 {
27 //
28 // Constructor from WSeries, add target channel
29 //
30 // in: TF transform of target
31 // ch: channel name
32 // fL, fH: low and high frequency
33 //
34 // It calls the add function
35 //
36 
37  this->add(in,ch,fL,fH);
38 }
39 
41  Edge(0.), pOUT(false)
42 {
43 //
44 // Copy constructor
45 //
46 
47  *this = value;
48 }
49 
51 {
52 //
53 // Copy constructor using operator
54 //
55 
56  std::cout<<"I am =\n";
57  this->chList = value.chList; // TF data: 0 - target, >0 - withess
58  this->FILTER = value.FILTER; // total Wiener filter
59  this->chName = value.chName; // channel names
60  this->chMask = value.chMask; // vectors of selected/rejected (1/0) layers
61  this->matrix = value.matrix; // symmetric matrix
62  this->vCROSS = value.vCROSS; // cross-correlation vector
63  this->vEIGEN = value.vEIGEN; // eigenvaue vector
64  this->target = value.target; // target time series
65  this->rnoise = value.rnoise; // regressed out noise
66  this->WNoise = value.WNoise; // Wavelet series for regressed out noise
67  this->Edge = value.Edge; // boundary buffer
68  this->pOUT = value.pOUT; // printout flag
69 
70  return *this;
71 }
72 
73 size_t regression::add(WSeries<double> &in, char* ch, double fL, double fH) {
74 //
75 // clear and add target channel to the regression list
76 //
77 // in: target channel, time-frequency series
78 // ch: name - channel name
79 // fL, fH: - fL <-- frequency band --> fH. If fL>=fH - set boundaries from input in
80 //
81 // return value: number of channels in the regression list
82 //
83 
84  this->clear();
85  size_t I = in.maxLayer();
86  wavearray<int> mask(I+1); mask=0;
88  ww = in;
89  ww.setlow(fL>=fH ? in.getlow() : fL);
90  ww.sethigh(fL>=fH ? in.gethigh() : fH);
91 
92  for(size_t i=1; i<I; i++) {
93  if(ww.frequency(i)<=ww.gethigh() && ww.frequency(i)>=ww.getlow())
94  mask.data[i] = 1;
95  }
96 
97  this->chList.push_back(ww);
98  this->chName.push_back(ch);
99  this->chMask.push_back(mask);
100  ww.Inverse();
101  this->target = ww;
102  return chList.size();
103 }
104 
105 size_t regression::add(wavearray<double> &in, char* ch, double fL, double fH) {
106 //
107 // add witness channel to the regression list
108 //
109 // in: witness time series
110 // ch: channel name, if ch name is the same as witness, build LPE filter
111 // fL, fH: fL <-- frequency band --> fH
112 //
113 // return value: number of channels in the regression list
114 //
115 
116  if(!chList.size()) {
117  std::cout<<"regression::add() ERROR: add target first\n";
118  exit(1);
119  }
120  if(chList.size()>1 && strcmp(this->chName[0],this->chName[1])==0) {
121  std::cout<<"regression::add() ERROR: second witness channel can not be added to LPEF\n";
122  exit(1);
123  }
124  in -= in.mean();
125  if(in.rms()==0.) {
126  std::cout<<"regression::add() warning: empty witness channel.\n";
127  return 0;
128  }
129  double rate = chList[0].rate();
130  size_t size = chList[0].size();
131  int n = in.rate()>rate ? int(in.rate())%int(rate) : int(rate)%int(in.rate());
132  int m = in.rate()>rate ? int(in.rate()/rate) : -int(rate/in.rate());
133  int l = int(log(fabs(m))/log(2.)+0.1);
134 
135  if(n || abs(m)!=int(pow(2.,l)+0.1)) {
136  std::cout<<"regression::add() ERROR: incompatible target and witness channels\n";
137  exit(1);
138  }
139 
140  Biorthogonal<double> Bio(512);
141  WSeries<double> ws(Bio);
143  WDM<double>* wdm = (WDM<double>*)this->chList[0].pWavelet->Clone();
144 
145  if(m==1) { // copy witness
146  ts = in;
147  }
148  else if(m>0) { // downsample witness
149  ws.Forward(in,l);
150  ws.getLayer(ts,0);
151  }
152  else { // upsample witness
153  ws.resize(abs(m)*in.size());
154  ws.rate(abs(m)*in.rate());
155  ws.Forward(l); ws=0.;
156  ws.putLayer(in,0);
157  ws.Inverse(); ts = ws;
158  }
159  ws.Forward(ts,*wdm);
160 
161  if(ws.rate()!=rate || ws.size() != size) {
162  std::cout<<"regression::add() ERROR: incompatible target and witness channels\n";
163  exit(1);
164  }
165 
166  ws.setlow(fL>=fH ? 0. : fL);
167  ws.sethigh(fL>=fH ? ws.rate()/2. : fH);
168 
169  size_t I = ws.maxLayer()+1;
170  wavearray<int> mask(I); mask=0;
171 
172  for(size_t i=0; i<I; i++) {
173  if(ws.frequency(i)<=ws.gethigh() && ws.frequency(i)>=ws.getlow())
174  mask.data[i] = 1;
175  //ws.getLayer(ts,i);
176  //if (ws.frequency(i)>500 && ws.frequency(i)<540) std::cout << ws.frequency(i) << " " << ts.mean() << " " << ts.rms() << endl;
177  //if (ts.rms()==0) mask.data[i]=0;
178  }
179 
180  this->chList.push_back(ws);
181  this->chName.push_back(ch);
182  this->chMask.push_back(mask);
183 
184  delete wdm;
185 
186  return chList.size();
187 }
188 
189 size_t regression::add(int n, int m, char* ch) {
190 //
191 // construct witness channel from 2 channels and add
192 //
193 // n, m channel indeces in chList
194 // ch: channel name
195 //
196 // return value: number of channels in the regression list
197 //
198 
199 
200  if(chList.size()>1 && strcmp(this->chName[0],this->chName[1])==0) {
201  std::cout<<"regression::add() ERROR: second witness channel can not be added to LPEF\n";
202  exit(1);
203  }
204  if(n>=(int)chList.size() || m>=(int)chList.size()) {
205  std::cout<<"regression::add() WARNING: can't construct witness channel.\n";
206  return chList.size();
207  }
208 
209  WSeries<double> wsn = this->chList[n];
210  WSeries<double> wsm = this->chList[m];
211  wavearray<int>* pMn = &(this->chMask[n]);
212  wavearray<int>* pMm = &(this->chMask[m]);
214  WDM<double>* wdm = (WDM<double>*)this->chList[0].pWavelet->Clone();
215 
216  size_t I = wsn.maxLayer()+1;
218  wavearray<int> mask(I); mask=0;
219 
220  for(size_t i=0; i<I; i++) {
221  if(!(pMn->data[i])) {
222  wsn.getLayer(wa,i); wa = 0;
223  wsn.putLayer(wa,i);
224  }
225  if(!(pMm->data[i])) {
226  wsm.getLayer(wa,i); wa = 0;
227  wsm.putLayer(wa,i);
228  }
229  if(!(pMn->data[i])) continue;
230  for(size_t j=0; j<I; j++) {
231  if(!(pMm->data[j])) continue;
232  if(i+j < I-1) mask.data[i+j] = 1;
233  if(i>j) mask.data[i-j] = 1;
234  if(i<j) mask.data[j-i] = 1;
235  }
236  }
237 
238  wsn.setlow(chList[0].getlow());
239  wsn.sethigh(chList[0].gethigh());
240  for(size_t i=0; i<I; i++) {
241  if(wsn.frequency(i)>wsn.gethigh() || wsn.frequency(i)<wsn.getlow())
242  mask.data[i] = 0;
243  }
244 
245  wsn.Inverse();
246  wsm.Inverse();
247  wsm *= wsn; ts = wsm; ts -= ts.mean();
248  wsn.Forward(ts, *wdm);
249  this->chList.push_back(wsn);
250  this->chName.push_back(ch);
251  this->chMask.push_back(mask);
252 
253  delete wdm;
254 
255  return chList.size();
256 }
257 
258 size_t regression::setFilter(size_t K) {
259 //
260 // set Wiener filter structure
261 //
262 // K: half-size length of a unit filter: each witness layer has filter 2(2*k+1) long
263 //
264 // return value: total number of witness layers
265 //
266 // FILTER.filter[0] is not used
267 //
268 
269  size_t i,j;
270  size_t I = this->chList[0].maxLayer();
271  size_t count = 0;
272  wavearray<int>* pT = &(this->chMask[0]);
273  wavearray<int>* pW;
274  vectorD f; f.resize(2*K+1);
275 
276  this->FILTER.clear(); std::vector<Wiener>().swap(FILTER);
277  this->matrix.clear(); std::vector<TMatrixDSym>().swap(matrix);
278  this->vCROSS.clear(); std::vector< wavearray<double> >().swap(vCROSS);
279  this->vEIGEN.clear(); std::vector< wavearray<double> >().swap(vEIGEN);
280 
281  this->kSIZE = K;
282  for(i=1; i<I; i++) {
283  if(!pT->data[i]) continue; // check target mask
284 
285  Wiener WF;
286 
287  for(j=0; j<chList.size(); j++) {
288  pW = &(this->chMask[j]);
289  if(!pW->data[i]) continue; // check witness mask
290  WF.channel.push_back(j);
291  WF.layer.push_back(i);
292  WF.norm.push_back(1.);
293  WF.filter00.push_back(f);
294  WF.filter90.push_back(f);
295  if(j) count++;
296  }
297  if(WF.channel.size()>1) this->FILTER.push_back(WF);
298  else pT->data[i] = 0; // mask target layer
299  }
300  return count;
301 }
302 
303 void regression::mask(int n, double flow, double fhigh) {
304 //
305 // mask (set to 0) Mask vector
306 //
307 // n: channel index
308 // flow, high: low and high frequency bounds
309 //
310 
311  int* pM = this->chMask[n].data;
312  WSeries<double>* pW = &(this->chList[n]);
313  size_t I = pW->maxLayer()+1;
314  for(size_t i=0; i<I; i++) {
315  if((pW->frequency(i)<fhigh && pW->frequency(i)>flow) || fhigh<=flow)
316  pM[i] = 0;
317  }
318 
319 }
320 
321 void regression::unmask(int n, double flow, double fhigh) {
322 //
323 // unmask (set to 1) Mask vector
324 //
325 // n: channel index
326 // flow, high: low and high frequency bounds
327 //
328 
329  int* pM = this->chMask[n].data;
330  WSeries<double>* pW = &(this->chList[n]);
331  size_t I = pW->maxLayer()+1;
332  for(size_t i=0; i<I; i++) {
333  if((pW->frequency(i)<fhigh && pW->frequency(i)>flow) || fhigh<=flow)
334  pM[i] = 1;
335  }
336 
337 }
338 
340 //
341 // extract eigenvalues
342 //
343 // n: channel index. If -1 return all eigenvalues
344 //
345 // return wavearray: eigen-values
346 //
347 
349  size_t I = this->vEIGEN.size();
350  if(n<(int)I && n>=0) return vEIGEN[n];
351  for(size_t i=0; i<I; i++) E.append(this->vEIGEN[i]);
352  return E;
353 }
354 
356 //
357 // extract filter for target index nT and witness nW
358 //
359 // c: 'f' or 'F' - extract 0-phase or 90-phase filter
360 // 'n' or 'N' - extract 0-phase or 90-phase norm
361 // 'a' or 'A' - extract 0-phase - 90-phase asymmetry
362 // nT: target index. If <0, return full vector
363 // nW: witness index
364 //
365 // return wavearray: filter values
366 //
367 
368  wavearray<double> OUT;
369  int N = int(this->FILTER.size());
370  int i,I,j,J;
371 
372  if(nT>=0 && nT<N) {i=nT; I=nT+1;}
373  else {i=0; I=N;}
374 
375  for(int n=i ;n<I;n++) {
376  int M = int(this->FILTER[n].filter00.size());
377 
378  if(nW>0 && nW<M) {j=nW; J=nW+1;}
379  else {j=1; J=M;}
380 
381  for(int m=j; m<J; m++) {
382  int K = int(this->FILTER[n].filter00[m].size());
383  wavearray<double> e(K);
384  wavearray<double> E(K);
385 
386  for(int k=0; k<K; k++) {
387  e.data[k] = this->FILTER[n].filter00[m][k];
388  E.data[k] = this->FILTER[n].filter90[m][k];
389  }
390 
391  if(c=='a') {
392  double re = e.rms();
393  double RE = E.rms();
394  OUT.append((re*re-RE*RE)/(re*re+RE*RE));
395  }
396  else if(c=='n') OUT.append(e.rms());
397  else if(c=='N') OUT.append(E.rms());
398  else {
399  if(c!='F') OUT.append(e);
400  if(c!='f') OUT.append(E);
401  }
402  }
403  }
404  return OUT;
405 }
406 
407 void regression::setMatrix(double edge, double f) {
408 //
409 // set system of linear equations: M * F = V
410 // M = matrix array, V = vector of free coefficients, F = filters
411 //
412 // edge: boundary interval in seconds excluded from training
413 // 1-f: tail fraction to be removed during training
414 //
415 
416  int i,j,n,m,nn,mm,N;
417  int I = this->FILTER.size(); // number of target layers
418  int K = (int)this->kSIZE; // unit filter half-size
419  int K2 = 2*K;
420  int K4 = 2*(2*K+1);
421  int k,ii,jj,kk;
422  int sIZe,sTEp,k00n,k90n,k00m,k90m;
423 
424  double FLTR = !strcmp(this->chName[0],this->chName[1]) ? 0. : 1.; // LPEF : Wiener
425  double fm = fabs(f);
426 
427  double fraction, rms;
428  WSeries<double>* pTFn;
429  WSeries<double>* pTFm;
433  std::slice s00n,s90n;
434  std::slice s00m,s90m;
435  Wiener* pW;
436  SymmArray<double> acf(K2);
437  SymmArray<double> ccf(K2);
438 
439  this->Edge = edge;
440  this->target.edge(edge);
441  for(i=0; i<(int)this->chList.size(); i++) { // loop on channels
442  this->chList[i].edge(edge); // set edge
443  }
444  if(!I) {std::cout<<"regression::nothing to clean.\n"; return;}
445 
446  this->matrix.clear(); std::vector<TMatrixDSym>().swap(matrix);
447  this->vCROSS.clear(); std::vector< wavearray<double> >().swap(vCROSS);
448  this->vEIGEN.clear(); std::vector< wavearray<double> >().swap(vEIGEN);
449 
450  for(i=0; i<I; i++) { // loop on target layers
451  pW = &FILTER[i];
452  N = pW->layer.size();
453  TMatrixDSym MS((N-1)*K4);
454  wavearray<double> VC((N-1)*K4);
455 
456  pTFm = &(this->chList[0]); // target
457  s00m = pTFm->pWavelet->getSlice(pW->layer[0]);
458  s90m = pTFm->pWavelet->getSlice(-pW->layer[0]);
459  k00m = s00m.start();
460  k90m = s90m.start();
461  sIZe = s00m.size(); // # of samples in a layer
462  sTEp = s00m.stride(); // layer stide
463 
464  for(n=0; n<N; n++) { // normalization & cross-vector
465  pTFn = &(this->chList[pW->channel[n]]);
466 
467  // normalize target & witness channels
468 
469  pTFn->getLayer(qq,pW->layer[n]);
470  pTFn->getLayer(ww,pW->layer[n]);
471  pTFn->getLayer(WW,-pW->layer[n]);
472 
473  for(j=0; j<(int)ww.size(); j++) {
474  qq.data[j] = ww.data[j]*ww.data[j]+WW.data[j]*WW.data[j];
475  }
476  rms = sqrt(qq.mean(fm));
477  pW->norm[n] = rms;
478 // std::cout<<i<<" "<<rms<<std::endl;
479  if(!n) continue; // skip target
480 
481  // calculate target-witness cross-vector
482 
483  s00n = pTFn->pWavelet->getSlice(pW->layer[n]);
484  s90n = pTFn->pWavelet->getSlice(-pW->layer[n]);
485  k00n = s00n.start();
486  k90n = s90n.start();
487 
488  for(k=-K; k<=K; k++) {
489 
490  if(k<0) { // n-channel is shifted
491  k00n = s00n.start();
492  k90n = s90n.start();
493  k00m = s00m.start()-k*sTEp;
494  k90m = s90m.start()-k*sTEp;
495  }
496  else { // m-channel is shifted
497  k00n = s00n.start()+k*sTEp;
498  k90n = s90n.start()+k*sTEp;
499  k00m = s00m.start();
500  k90m = s90m.start();
501  }
502  for(j=K; j<sIZe-K; j++) {
503  jj = j*sTEp;
504  ww.data[j] = pTFn->data[k00n+jj]*pTFm->data[k00m+jj];
505  ww.data[j]+= pTFn->data[k90n+jj]*pTFm->data[k90m+jj];
506  WW.data[j] = pTFm->data[k90m+jj]*pTFn->data[k00n+jj];
507  WW.data[j]-= pTFm->data[k00m+jj]*pTFn->data[k90n+jj];
508  }
509 
510  kk = (n-1)*K4+K+k;
511  VC.data[kk] = ww.mean(fm)/pW->norm[0]/pW->norm[n]; // Chx+Ch'x'
512  VC.data[kk+K4/2] = WW.mean(fm)/pW->norm[0]/pW->norm[n]; // Ch'x-Chx'
513  if(k==0) {VC.data[kk] *= FLTR; VC.data[kk+K4/2] *= FLTR;} // handle LPE filter
514 // printf("%f ",VC.data[kk+0*K4/2]);
515  }
516 // std::cout<<std::endl;
517  }
518  ww=0.; WW=0.;
519 
520  for(n=1; n<N; n++) { // loop on witness channels
521  nn = (n-1)*K4+K; // index of ZERO element of ii block
522  pTFn = &(this->chList[pW->channel[n]]); // witness channel i
523  s00n = pTFn->pWavelet->getSlice(pW->layer[n]);
524  s90n = pTFn->pWavelet->getSlice(-pW->layer[n]);
525  k00n = s00n.start();
526  k90n = s90n.start();
527  sIZe = s00n.size();
528  sTEp = s00n.stride();
529 
530  for(m=n; m<N; m++) {
531  mm = (m-1)*K4+K; // index of ZERO element of jj block
532  pTFm = &(this->chList[pW->channel[m]]); // witness channel j
533  s00m = pTFm->pWavelet->getSlice(pW->layer[m]);
534  s90m = pTFm->pWavelet->getSlice(-pW->layer[m]);
535  k00m = s00m.start();
536  k90m = s90m.start();
537 
538 // printf("n/m=%2d/%2d %d %d %d %d\n",n,m,sIZe,sTEp,(int)s00n.start(),(int)s90n.start());
539 
540  for(k=-K2; k<=K2; k++) {
541 
542  if(k<0) {
543  k00n = s00n.start();
544  k90n = s90n.start();
545  k00m = s00m.start()-k*sTEp;
546  k90m = s90m.start()-k*sTEp;
547  }
548  else {
549  k00n = s00n.start()+k*sTEp;
550  k90n = s90n.start()+k*sTEp;
551  k00m = s00m.start();
552  k90m = s90m.start();
553  }
554 
555  for(j=K2; j<sIZe-K2; j++) {
556  jj = j*sTEp;
557  ww.data[j] = pTFn->data[k00n+jj]*pTFm->data[k00m+jj];
558  ww.data[j]+= pTFn->data[k90n+jj]*pTFm->data[k90m+jj];
559  WW.data[j] = pTFm->data[k00m+jj]*pTFn->data[k90n+jj];
560  WW.data[j]-= pTFm->data[k90m+jj]*pTFn->data[k00n+jj];
561  }
562 
563  acf[k] = ww.mean(fm)/pW->norm[n]/pW->norm[m]; // "auto" function
564  ccf[k] = WW.mean(fm)/pW->norm[n]/pW->norm[m]; // "cros" function
565 
566 // printf("%f ",acf[k]);
567  }
568 // std::cout<<"\n";
569 
570  for(ii=-K; ii<=K; ii++) { // fill in matrix
571  for(jj=-K; jj<=K; jj++) {
572  MS[nn+ii][mm+jj] = (ii==0 || jj==0) ? acf[ii-jj]*FLTR : acf[ii-jj];
573  MS[mm+jj][nn+ii] = (ii==0 || jj==0) ? acf[ii-jj]*FLTR : acf[ii-jj];
574  MS[nn+ii][mm+jj+K4/2] = (ii==0 || jj==0) ? ccf[ii-jj]*FLTR : ccf[ii-jj];
575  MS[mm+jj+K4/2][nn+ii] = (ii==0 || jj==0) ? ccf[ii-jj]*FLTR : ccf[ii-jj];
576  MS[nn+ii+K4/2][mm+jj] = (ii==0 || jj==0) ? -ccf[ii-jj]*FLTR :-ccf[ii-jj];
577  MS[mm+jj][nn+ii+K4/2] = (ii==0 || jj==0) ? -ccf[ii-jj]*FLTR :-ccf[ii-jj];
578  MS[nn+ii+K4/2][mm+jj+K4/2] = (ii==0 || jj==0) ? acf[ii-jj]*FLTR : acf[ii-jj];
579  MS[mm+jj+K4/2][nn+ii+K4/2] = (ii==0 || jj==0) ? acf[ii-jj]*FLTR : acf[ii-jj];
580  }
581  }
582  }
583  }
584  // for(int i=0;i<VC.size();i++)
585  // std::cout<<VC[i]<<std::endl;
586  // MS.Print();
587  this->matrix.push_back(MS);
588  this->vCROSS.push_back(VC);
589  }
590 }
591 
592 void regression::solve(double th, int nE, char c) {
593 //
594 // solve for eigenvalues and calculate Wiener filters
595 //
596 // th: eigenvalue threshold. If <0, than in units of max eigenvalue
597 // nE: number of selected eigenvalues (minimum is 4)
598 // c: regulator h=mild, s=soft, h=hard
599 //
600 
601  int i, j, k;
602  int nA = this->FILTER.size();
603  int nR = this->matrix.size();
604  int nC = this->vCROSS.size();
605  int ne;
606 
607  if(!nA) {std::cout<<"regression::nothing to clean.\n"; return;}
608  if (nA!=nR || nA!=nC || nA==0) {
609  std::cout << "Error, filter not initialized\n"; exit(1);
610  }
611 
612  int K = (int) this->kSIZE;
613  int K4=2*(2*K+1);
614  double temp;
615 
616  this->vEIGEN.clear(); std::vector< wavearray<double> >().swap(vEIGEN);
617 
618  for (i=0; i<nA; i++) { //loop on wavelet layers
619  int J = FILTER[i].layer.size()-1; //number of witness channels
620  ne = nE<=0 ? K4*J : nE-1;
621  if(ne>=K4*J) ne = K4*J-1;
622  if(ne<1) ne = 1;
623 
624  wavearray<double> lambda(K4*J); //inverse of regulators
625  TMatrixDSym R(this->matrix[i]);
626  wavearray<double>* pC = &(this->vCROSS[i]);
627  TMatrixDSymEigen QP(R);
628  TVectorD eigenval(QP.GetEigenValues());
629  TMatrixD eigenvec(QP.GetEigenVectors());
630 
631  //select threshold
632  double last = 0.;
633  double TH = th<0. ? -th*eigenval[0] : th+1.e-12;
634  int nlast = -1;
635  for(j=0; j<K4*J; j++) {
636  lambda.data[j] = eigenval[j];
637  if(eigenval[j]>=TH) nlast++;
638  if(j<K4*J-1)
639  if(eigenval[j+1]>eigenval[j])
640  std::cout<<eigenval[j]<<" "<<eigenval[j+1]<<endl;
641  }
642  if(nlast<1) nlast = 1;
643 
644  this->vEIGEN.push_back(lambda);
645  if(nlast>ne) nlast = ne;
646  lambda = 0.;
647 
648  //std::cout<<i<<" "<<J<<" "<<TH<<" "<<nlast<<" "<<K<<" "<<eigenval[nlast]<<std::endl;
649 
650  switch(c) { //regulators
651  case 'h':
652  last = 0.;
653  break;
654  case 's':
655  last = eigenval[nlast]>0. ? 1./eigenval[nlast] : 0.;
656  break;
657  case 'm':
658  last = 1./eigenval[0];
659  break;
660  }
661  for(j=0;j<=nlast;j++) lambda[j] = eigenval[j]>0. ? 1./eigenval[j] : 0.;
662  for(j=nlast+1;j<K4*J;j++) lambda[j] = last;
663 
664  //calculte filters
665  wavearray<double> vv(K4*J); // lambda * eigenvector * cross vector
666  wavearray<double> aa(K4*J); // eigenvector^T * lambda * eigenvector * cross vector
667 
668  for(j=0;j<K4*J;j++) {
669  temp=0.0;
670  for(k=0;k<K4*J;k++) temp += eigenvec[k][j]*pC->data[k];
671  vv.data[j]=temp*lambda.data[j];
672  }
673 
674  for(j=0;j<K4*J;j++) {
675  temp=0.0;
676  for(k=0;k<K4*J;k++) temp += eigenvec[j][k]*vv.data[k];
677  aa.data[j]=temp;
678  }
679 
680  //save filters coefficients
681  for (j=0; j<J; j++) {
682  for(k=0;k<=2*K;k++) {
683  FILTER[i].filter00[j+1][k] = aa.data[j*K4+k];
684  FILTER[i].filter90[j+1][k] = aa[j*K4+k+K4/2];
685  }
686  }
687  }
688 }
689 
690 
691 void regression::apply(double threshold, char c) {
692 //
693 // apply filter to target channel and produce noise TS
694 //
695 // threshold: filter is not applied if regressed noise rms < threshold
696 // mask channels with rms<threshold, if c='m' or c='M'.
697 // mask changes configuration of witness channels
698 // c: 'n' or 'N' - either 0-phase or 90-phase noise
699 // 'a' - (0-phase + 90-phase)/2 noise, standard RMS threshold
700 // 'A' - (0-phase + 90-phase)/2 noise, differential RMS threshold
701 // 'm' - (0-phase + 90-phase)/2 noise, mask, standard RMS threshold
702 // 'M' - (0-phase + 90-phase)/2 noise, mask, differential RMS threshold
703 //
704 
705  int n,m;
706  int nA = this->FILTER.size();
707  int nR = this->matrix.size();
708  int nC = this->vCROSS.size();
709 
710  if(!nA) {std::cout<<"regression::nothing to clean.\n"; return;}
711  if (nA!=nR || nA!=nC || nA==0) {
712  std::cout << "Error, filter is not initialized\n"; exit(1);
713  }
714 
715  int K = (int) this->kSIZE;
716  int N = this->FILTER.size(); // number of target layers
717 
718  wavearray<double> tmp(N); tmp=0.; //RANK
719  int L = this->chList.size(); //RANK
720  for(n=0; n<L; n++) this->vrank.push_back(tmp); //RANK
721 
726  WSeries<double>* pT;
727  Wiener* pF;
728 
729  double sum, SUM;
730  WSeries<double> wnoise(this->chList[0]);
731  WSeries<double> WNOISE(this->chList[0]);
732  wnoise=0; WNOISE=0.;
733 
734  for(n=0; n<N; n++) { // loop on target layers
735  std::vector< wavearray<double> > wno;
736  std::vector< wavearray<double> > WNO;
737  _apply_(n,wno,WNO);
738 
739  pF = &FILTER[n]; // pointer to Weiner structure
740  pT = &(this->chList[0]); // pointer to target WS
741  pT->getLayer(ww,pF->layer[0]); // get target 00 data
742  pT->getLayer(WW,-pF->layer[0]); // get target 90 data
743  ww = 0.; WW = 0.; // zero total noise arrays
744 
745  double freq = pT->frequency(pF->layer[0]); //RANK
746  vfreq.append(freq); //RANK
747 
748  int KK = int(ww.rate()*this->Edge); // calculate # of edge samples
749  if(KK<K) KK = K; KK++;
750  std::slice S(KK,ww.size()-2*KK,1); // slice for RMS calculation
751 
752  double trms = wno[0].rms(S);
753  double TRMS = WNO[0].rms(S);
754  int M = pF->layer.size(); // number of witness channels + 1
755 
756  for(m=1; m<M; m++) { // loop over witnwss channels
757  int j = (int)pF->channel[m]; // witness channel index
758  int* pM = this->chMask[j].data; // poiner to j's mask
759 
760  if(c=='A' || c=='M') {
761  qq=wno[0]; qq-=wno[m]; // channel j subtracted
762  QQ=WNO[0]; QQ-=WNO[m]; // channel j subtracted
763  sum = qq.rms(S);
764  SUM = QQ.rms(S);
765  sum = trms*trms-sum*sum;
766  SUM = TRMS*TRMS-SUM*SUM;
767  }
768  else {
769  sum = pow(wno[m].rms(S),2); // rms 00
770  SUM = pow(WNO[m].rms(S),2); // rms 90
771  }
772  if(sum<0.) sum = 0.;
773  if(SUM<0.) SUM = 0.;
774 
775  this->vrank[j].data[n] = sum+SUM; //RANK
776  this->vrank[0].data[n] += sum+SUM; //RANK
777 
778  if(sum+SUM < threshold*threshold) {
779  if(c=='m' || c=='M') pM[pF->layer[m]] = 0;
780  continue;
781  }
782  ww += wno[m]; WW += WNO[m];
783  //std::cout << "RMS: M: " << m << " " << wno[m].rms(S) << " " << WNO[m].rms(S) << std::endl;
784  }
785  if(c=='A' || c=='M') this->vrank[0].data[n] = trms*trms+TRMS*TRMS; //RANK
786  //std::cout << "RMS: N: " << n << " " << ww.rms() << " " << WW.rms() << std::endl;
787  ww *= pF->norm[0];
788  WW *= pF->norm[0];
789  //std::cout << "RMS: N: " << n << " " << ww.rms() << " " << WW.rms() << " " << pF->norm[0] << std::endl;
790  wnoise.putLayer(ww,pF->layer[0]);
791  wnoise.putLayer(WW,-pF->layer[0]);
792  }
793  WNOISE = wnoise;
794  this->WNoise = wnoise;
795  wnoise.Inverse();
796  WNOISE.Inverse(-2);
797 
798  if(c=='n') this->rnoise = wnoise;
799  else if(c=='N') this->rnoise = WNOISE;
800  else{
801  this->rnoise = wnoise;
802  this->rnoise += WNOISE;
803  this->rnoise *= 0.5;
804  }
805 }
806 
807 wavearray<double> regression::rank(int nbins, double fL, double fH) {
808 //
809 // get RMS for predicted noise (per layer/bin) for all witness channels
810 // (including masked)
811 //
812 // nbins: number of loudest frequency layers to calculate rms/layer.
813 // nbins=0 - take all layers
814 // if positive, return standard RMS
815 // if negative, return differential RMS
816 // fL, fH: frequency range: low, high
817 //
818 // return wavearray: rank values
819 //
820 
821  int n,m;
822  int nA = this->FILTER.size();
823  int nR = this->matrix.size();
824  int nC = this->vCROSS.size();
825 
826  if(!nA) {
827  std::cout<<"regression::nothing to clean.\n";
828  return wavearray<double>();
829  }
830  if (nA!=nR || nA!=nC) {
831  std::cout << "Error, filter not initialized\n"; exit(1);
832  }
833 
834  int N = this->FILTER.size(); // number of target layers
835  int L = this->chList.size();
836 
837  wavearray<double> tmp(N); tmp=0.;
838  wavearray<double> rms(L); rms=0.;
839  std::vector< wavearray<double> > vrms;
840  for(n=0; n<L; n++) { // loop on channels
841  vrms.push_back(tmp);
842  int tsize=this->vfreq.size();
843  for (int i=0; i<tsize; i++) {
844  double freq = vfreq.data[i];
845  if(fL<fH && (freq>fH || freq<fL)) continue;
846  vrms[n].data[i]=vrank[n].data[i];
847  }
848  }
849 
850  for(n=0; n<L; n++) { // loop on channels
851  vrms[n].waveSort();
852  nA = 0;
853  for(m=N-1; m>=N-abs(nbins); m--) { // loop over loud layers
854  rms.data[n] += vrms[n].data[m];
855  if(vrms[n].data[m]>0.) nA++;
856  }
857  if(nA) rms.data[n] /= nA;
858  rms.data[n] = sqrt(rms.data[n]);
859  }
860  return rms;
861 }
862 
864  std::vector< wavearray<double> > &wno,
865  std::vector< wavearray<double> > &WNO)
866 {
867 //
868 // internal regression apply function
869 //
870 
871  int m,k,i;
876  WSeries<double>* pT;
877  WSeries<double>* pW;
878  Wiener* pF;
879  std::vector<double> *ff, *FF;
880  double val, VAL;
881 
882  pF = &FILTER[n]; // pointer to Weiner structure
883  pT = &(this->chList[0]); // pointer to target WS
884  pT->getLayer(ww,pF->layer[0]); // get target 00 data
885  pT->getLayer(WW,-pF->layer[0]); // get target 90 data
886 
887  ww = 0.; WW = 0.; // zero total noise arrays
888  wno.clear();
889  WNO.clear();
890  wno.push_back(ww);
891  WNO.push_back(WW);
892 
893  int K = (int) this->kSIZE;
894  int M = pF->layer.size(); // number of witness channels + 1
895  for(m=1; m<M; m++) { // loop over witnwss channels
896  int j = (int)pF->channel[m]; // witness channel index
897  pW = &(this->chList[j]); // pointer to witness WS
898  pW->getLayer(qq,pF->layer[m]); // get witness 00 data
899  pW->getLayer(QQ,-pF->layer[m]); // get witness 90 data
900  ff = &(pF->filter00[m]); // get 00 filter pointer
901  FF = &(pF->filter90[m]); // get 90 filter pointer
902  wavearray<double> nn(qq.size()); // array for witness prediction data
903  wavearray<double> NN(qq.size()); // array for witness prediction data
904  nn.rate(qq.rate());
905  NN.rate(QQ.rate());
906  qq *= 1./pF->norm[m];
907  QQ *= 1./pF->norm[m];
908 
909  for(i=0; i<(int)qq.size(); i++) {
910  val = VAL = 0.;
911  if(i>=K && i<int(qq.size())-K) {
912  for(k=-K; k<=K; k++) {
913  val += (*ff)[k+K]*qq.data[i+k] - (*FF)[k+K]*QQ.data[i+k];
914  VAL += (*FF)[k+K]*qq.data[i+k] + (*ff)[k+K]*QQ.data[i+k];
915  }
916  }
917  nn.data[i] = val;
918  NN.data[i] = VAL;
919  wno[0].data[i] += val;
920  WNO[0].data[i] += VAL;
921  }
922  wno.push_back(nn);
923  WNO.push_back(NN);
924  }
925 }
926 
927 
928 
929 
930 
void sethigh(double f)
Definition: wseries.hh:114
size_t append(const wavearray< DataType_t > &)
Definition: wavearray.cc:775
virtual void resize(unsigned int)
Definition: wseries.cc:883
virtual size_t size() const
Definition: wavearray.hh:127
char val[20]
void setMatrix(double edge=0., double f=1.)
Definition: regression.cc:407
par[0] value
void _apply_(int n, std::vector< wavearray< double > > &w, std::vector< wavearray< double > > &W)
Definition: regression.cc:863
tuple f
Definition: cwb_online.py:91
virtual void rate(double r)
Definition: wavearray.hh:123
wavearray< double > target
Definition: regression.hh:167
virtual void edge(double s)
Definition: wavearray.hh:125
int n
Definition: cwb_net.C:10
TMatrixTSym< double > TMatrixDSym
Definition: regression.hh:20
int count
Definition: compare_bkg.C:373
WDM< double > wdm(nLAYERS, nLAYERS, 6, 10)
void setlow(double f)
Definition: wseries.hh:107
std::vector< wavearray< double > > vCROSS
Definition: regression.hh:165
Long_t size
#define M
Definition: UniqSLagsList.C:3
int m
Definition: cwb_net.C:10
size_t add(WSeries< double > &target, char *name, double fL=0., double fH=0.)
Definition: regression.cc:73
int j
Definition: cwb_net.C:10
i drho i
std::vector< double > vectorD
Definition: network.hh:33
wc clear()
#define N
virtual double mean() const
Definition: wavearray.cc:1053
tuple ff
Definition: cwb_online.py:394
std::vector< wavearray< double > > vEIGEN
Definition: regression.hh:166
nT
Definition: cbc_plots.C:659
std::vector< Wiener > FILTER
Definition: regression.hh:163
std::vector< WSeries< double > > chList
Definition: regression.hh:160
virtual double rms()
Definition: wavearray.cc:1188
wavearray< double > rank(int nbins=0, double fL=0., double fH=0.)
Definition: regression.cc:807
double getlow() const
Definition: wseries.hh:111
std::vector< vectorD > filter00
Definition: regression.hh:27
size_t size() const
Definition: wslice.hh:71
i() int(T_cor *100))
bool log
Definition: WaveMDC.C:41
double fhigh
double * tmp
Definition: testWDM_5.C:31
void apply(double threshold=0., char c='a')
Definition: regression.cc:691
int getLayer(wavearray< DataType_t > &w, double n)
param: n - layer number
Definition: wseries.cc:175
std::vector< TMatrixDSym > matrix
Definition: regression.hh:164
std::vector< wavearray< int > > chMask
Definition: regression.hh:162
int k
void solve(double th, int nE=0, char c='s')
Definition: regression.cc:592
std::vector< int > channel
Definition: regression.hh:24
wavearray< double > vfreq
Definition: regression.hh:171
std::vector< vectorD > filter90
Definition: regression.hh:28
size_t setFilter(size_t)
Definition: regression.cc:258
wavearray< double > getVEIGEN(int n=-1)
Definition: regression.cc:339
double e
double Edge
Definition: regression.hh:157
WSeries< double > ww
Definition: Regression_H1.C:33
double gethigh() const
Definition: wseries.hh:118
size_t kSIZE
Definition: regression.hh:156
void clear()
Definition: regression.hh:142
std::vector< int > layer
Definition: regression.hh:25
double flow
std::vector< char * > chName
Definition: regression.hh:161
WSeries< double > WNoise
Definition: regression.hh:169
ifstream in
virtual std::slice getSlice(const double)
Definition: WaveDWT.cc:111
std::vector< double > norm
Definition: regression.hh:26
wavearray< double > ts(N)
#define SUM
Definition: FrDisplay.cc:112
double fabs(const Complex &x)
Definition: numpy.cc:37
wavearray< double > getFILTER(char c='a', int nT=-1, int nW=-1)
Definition: regression.cc:355
void Forward(int n=-1)
param: wavelet - n is number of steps (-1 means full decomposition)
Definition: wseries.cc:228
void unmask(int n, double flow=0., double fhigh=0.)
Definition: regression.cc:321
Meyer< double > S(1024, 2)
void mask(int n, double flow=0., double fhigh=0.)
Definition: regression.cc:303
int l
Definition: cbc_plots.C:434
DataType_t * data
Definition: wavearray.hh:301
std::vector< wavearray< double > > vrank
Definition: regression.hh:170
WaveDWT< DataType_t > * pWavelet
Definition: wseries.hh:438
regression & operator=(const regression &)
Definition: regression.cc:50
wavearray< double > rnoise
Definition: regression.hh:168
size_t stride() const
Definition: wslice.hh:75
void Inverse(int n=-1)
param: n - number of steps (-1 means full reconstruction)
Definition: wseries.cc:273
void putLayer(wavearray< DataType_t > &, double n)
param: n - layer number
Definition: wseries.cc:201
double frequency(int l)
Definition: wseries.cc:99
size_t start() const
Definition: wslice.hh:67
int maxLayer()
Definition: wseries.hh:121
exit(0)