CGEM BOSS 6.6.5.g
BESIII Offline Software System
Loading...
Searching...
No Matches
CTF Class Reference

#include <InductionGar.h>

Public Member Functions

 CTF ()
 
 ~CTF ()
 
TH1D Convolution_Tbranch_fft (const double *Input_Curr_plot_001)
 return convolve of Input_Curr_plot_001 and T_Branch2, a function; tries to mimic Convolution_Tbranch has a option to test correctness.
 

Public Attributes

ConvolveWithConstconv
 
std::vector< double > TBranch
 

Detailed Description

Definition at line 16 of file InductionGar.h.

Constructor & Destructor Documentation

◆ CTF()

CTF::CTF ( )

Definition at line 360 of file InductionGar.cxx.

360 {
361 double xmin(0), xmax(1000);
362 const int Npx(2000);
363 TBranch.resize(Npx);// TBranch at x<0 is 0
364 double *p_TBranch=&(TBranch.front());
365 for (int i = 0; i < Npx; i++) {
366 double x=xmin+(xmax-xmin)*i/Npx;
367 TBranch[i]=T_branch2(x);
368 }// set initial TBranch values.
369 // basically, a N and a (2N-1) needs (3N-2) bins to avoid overlap.
370 conv=new ConvolveWithConst( p_TBranch,Npx,Npx);
371 //TBranch_fft.resize(Npx*3-2);
372
373 //conv.convolve
374}
Double_t x[10]
double T_branch2(double t)
ConvolveWithConst * conv
Definition: InductionGar.h:18
std::vector< double > TBranch
Definition: InductionGar.h:19
basically, convolve with fft. convolves 1-d. if size of 2 inputs are fixed, it may yield optimal spee...

◆ ~CTF()

CTF::~CTF ( )

Definition at line 375 of file InductionGar.cxx.

375 {
376 delete conv;
377}

Member Function Documentation

◆ Convolution_Tbranch_fft()

TH1D CTF::Convolution_Tbranch_fft ( const double *  Input_Curr_plot_001)

return convolve of Input_Curr_plot_001 and T_Branch2, a function; tries to mimic Convolution_Tbranch has a option to test correctness.

Parameters
Input_Curr_plot_001
Returns
TH1D

Definition at line 385 of file InductionGar.cxx.

386{
387 double xmin(0), xmax(1000);
388 const int Npx(2000);
389
390 TH1D h_signal("h_signal", "", Npx, xmin, xmax);
391 double *p_signal=h_signal.GetArray()+1;//considering underflow.
392 this->conv->convolve(p_signal,Input_Curr_plot_001,0,Npx,-0.5);// the '-' before 0.5: from original function.
393
394#ifdef INDUCTIONGAR_DEBUG_CONVOLUTION_TBRANCH_FFT_TEST
395 const double dx_abs_threhold=1e-6;// think about the accurancy of float; warn if |x| > threhold
396 const double dx_ref_threhold=1e-6;//think about the accurancy of float; warn if |dx/x| > threhold
397 TH1F h_signal_ref=Convolution_Tbranch_noskip(Input_Curr_plot_001);
398 bool found_mismatch=false;
399 for (int i=1; i<=Npx; i++){
400 if (abs((h_signal[i]-h_signal_ref[i])/h_signal_ref[i]) > dx_ref_threhold
401 and abs(h_signal[i]-h_signal_ref[i]) > dx_abs_threhold ) {
402 if (not found_mismatch){
403 std::cout<<"CgemDigitizerSvc::InductionGar::Convolution_Tbranch_fft found results not match: ";
404 found_mismatch=true;
405 }
406 std::cout<<"ref: "<<h_signal_ref[i]<<"; this: "<<h_signal[i]<<"; at i="<<i<<";threhold abs="<<dx_abs_threhold<<"; ref="<<dx_ref_threhold<<std::endl;
407 }
408 }
409#endif
410
411 return h_signal;
412}
double abs(const EvtComplex &c)
Definition: EvtComplex.hh:212
void convolve(double *output, const double *B, const int leftIndex, const int sizeofOut, double factor=1) const
do a convolve of stored const A and B, and put results to output.

Referenced by InductionGar::setMultiElectrons().

Member Data Documentation

◆ conv

ConvolveWithConst* CTF::conv

Definition at line 18 of file InductionGar.h.

Referenced by Convolution_Tbranch_fft(), CTF(), and ~CTF().

◆ TBranch

std::vector<double> CTF::TBranch

Definition at line 19 of file InductionGar.h.

Referenced by CTF().


The documentation for this class was generated from the following files: