CGEM BOSS 6.6.5.i
BESIII Offline Software System
Loading...
Searching...
No Matches
CEF2 Class Reference

#include <InductionGar2.h>

Public Member Functions

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

Public Attributes

ConvolveWithConstconv
 
std::vector< double > EBranch
 

Detailed Description

Definition at line 31 of file InductionGar2.h.

Constructor & Destructor Documentation

◆ CEF2()

CEF2::CEF2 ( )

Definition at line 483 of file InductionGar2.cxx.

483 {
484 double xmin(0), xmax(1000);
485 const int Npx(2000);
486 EBranch.resize(Npx); // EBranch at x<0 is 0
487 double *p_EBranch = &(EBranch.front());
488 for (int i = 0; i < Npx; i++) {
489 double x = xmin + (xmax - xmin) * i / Npx;
490 EBranch[i] = E_branch2(x);
491 } // set initial EBranch values.
492 // basically, a N and a (2N-1) needs (3N-2) bins to avoid overlap.
493 conv = new ConvolveWithConst(p_EBranch, Npx, Npx);
494 //EBranch_fft.resize(Npx*3-2);
495
496 //conv.convolve
497}
Double_t x[10]
double E_branch2(double t)
ConvolveWithConst * conv
std::vector< double > EBranch
basically, convolve with fft. convolves 1-d. if size of 2 inputs are fixed, it may yield optimal spee...

◆ ~CEF2()

CEF2::~CEF2 ( )

Definition at line 498 of file InductionGar2.cxx.

498 {
499 delete conv;
500}

Member Function Documentation

◆ Convolution_Ebranch_fft()

TH1D CEF2::Convolution_Ebranch_fft ( const double * Input_Curr_plot_001)

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

Parameters
Input_Curr_plot_001
Returns
TH1D

Definition at line 508 of file InductionGar2.cxx.

509{
510 double xmin(0), xmax(1000);
511 const int Npx(2000);
512
513 TH1D h_signal("h_signalE", "", Npx, xmin, xmax);
514 double *p_signal = h_signal.GetArray() + 1; //considering underflow.
515 this->conv->convolve(p_signal, Input_Curr_plot_001, 0, Npx, -0.5);
516
517#ifdef INDUCTIONGAR_DEBUG_CONVOLUTION_TBRANCH_FFT_TEST
518 const double dx_abs_threhold = 1e-6; // think about the accurancy of float; warn if |x| > threhold
519 const double dx_ref_threhold = 1e-6; //think about the accurancy of float; warn if |dx/x| > threhold
520 TH1F h_signal_ref = Convolution_Ebranch_noskip(Input_Curr_plot_001);
521 bool found_mismatch = false;
522 for (int i = 1; i <= Npx; i++) {
523 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) {
524 if (not found_mismatch) {
525 std::cout << "CgemDigitizerSvc::InductionGar2::Convolution_Ebranch_fft found results not match: ";
526 found_mismatch = true;
527 }
528 std::cout << "ref: " << h_signal_ref[i] << "; this: " << h_signal[i] << "; at i=" << i << ";threhold abs=" << dx_abs_threhold << "; ref=" << dx_ref_threhold << '\n';
529 }
530 }
531#endif
532
533 return h_signal;
534}
double abs(const EvtComplex &c)
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* CEF2::conv

Definition at line 34 of file InductionGar2.h.

Referenced by CEF2(), Convolution_Ebranch_fft(), and ~CEF2().

◆ EBranch

std::vector<double> CEF2::EBranch

Definition at line 35 of file InductionGar2.h.

Referenced by CEF2().


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