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

#include <InductionGar2.h>

Public Member Functions

 CTF2 ()
 
 ~CTF2 ()
 
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 20 of file InductionGar2.h.

Constructor & Destructor Documentation

◆ CTF2()

CTF2::CTF2 ( )

Definition at line 430 of file InductionGar2.cxx.

430 {
431 double xmin(0), xmax(1000);
432 const int Npx(2000);
433 TBranch.resize(Npx); // TBranch at x<0 is 0
434 double *p_TBranch = &(TBranch.front());
435 for (int i = 0; i < Npx; i++) {
436 double x = xmin + (xmax - xmin) * i / Npx;
437 TBranch[i] = T_branch2(x);
438 } // set initial TBranch values.
439 // basically, a N and a (2N-1) needs (3N-2) bins to avoid overlap.
440 conv = new ConvolveWithConst(p_TBranch, Npx, Npx);
441 //TBranch_fft.resize(Npx*3-2);
442
443 //conv.convolve
444}
Double_t x[10]
double T_branch2(double t)
ConvolveWithConst * conv
Definition: InductionGar2.h:22
std::vector< double > TBranch
Definition: InductionGar2.h:23
basically, convolve with fft. convolves 1-d. if size of 2 inputs are fixed, it may yield optimal spee...

◆ ~CTF2()

CTF2::~CTF2 ( )

Definition at line 445 of file InductionGar2.cxx.

445 {
446 delete conv;
447}

Member Function Documentation

◆ Convolution_Tbranch_fft()

TH1D CTF2::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 455 of file InductionGar2.cxx.

456{
457 double xmin(0), xmax(1000);
458 const int Npx(2000);
459
460 TH1D h_signal("h_signal", "", Npx, xmin, xmax);
461 double *p_signal = h_signal.GetArray() + 1; //considering underflow.
462 this->conv->convolve(p_signal, Input_Curr_plot_001, 0, Npx, -0.5); // the '-' before 0.5: from original function.
463
464#ifdef INDUCTIONGAR_DEBUG_CONVOLUTION_TBRANCH_FFT_TEST
465 const double dx_abs_threhold = 1e-6; // think about the accurancy of float; warn if |x| > threhold
466 const double dx_ref_threhold = 1e-6; //think about the accurancy of float; warn if |dx/x| > threhold
467 TH1F h_signal_ref = Convolution_Tbranch_noskip(Input_Curr_plot_001);
468 bool found_mismatch = false;
469 for (int i = 1; i <= Npx; i++) {
470 if (abs((h_signal[i] - h_signal_ref[i]) / h_signal_ref[i]) > dx_ref_threhold and abs(h_signal[i] - h_signal_ref[i]) > dx_abs_threhold) {
471 if (not found_mismatch) {
472 std::cout << "CgemDigitizerSvc::InductionGar2::Convolution_Tbranch_fft found results not match: ";
473 found_mismatch = true;
474 }
475 std::cout << "ref: " << h_signal_ref[i] << "; this: " << h_signal[i] << "; at i=" << i << ";threhold abs=" << dx_abs_threhold << "; ref=" << dx_ref_threhold << '\n';
476 }
477 }
478#endif
479
480 return h_signal;
481}
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 InductionGar2::setMultiElectrons().

Member Data Documentation

◆ conv

ConvolveWithConst* CTF2::conv

Definition at line 22 of file InductionGar2.h.

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

◆ TBranch

std::vector<double> CTF2::TBranch

Definition at line 23 of file InductionGar2.h.

Referenced by CTF2().


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