CGEM BOSS 6.6.5.f
BESIII Offline Software System
Loading...
Searching...
No Matches
InductionGar.cxx File Reference
#include "CgemDigitizerSvc/InductionGar.h"
#include "GaudiKernel/ISvcLocator.h"
#include "GaudiKernel/Bootstrap.h"
#include "GaudiKernel/IDataProviderSvc.h"
#include "GaudiKernel/SmartDataPtr.h"
#include "GaudiKernel/DataSvc.h"
#include "CLHEP/Units/PhysicalConstants.h"
#include "G4DigiManager.hh"
#include "Randomize.hh"
#include "G4ios.hh"
#include <iomanip>
#include <iostream>
#include <fstream>
#include <cmath>
#include "TRandom.h"

Go to the source code of this file.

Typedefs

typedef std::vector< int > Vint
 

Functions

double rectangle2 (double t, double *x, double *y, int Npx)
 
double T_branch2 (double t)
 
TH1F Convolution_Tbranch (double *Input_Curr_plot_001)
 

Typedef Documentation

◆ Vint

typedef std::vector<int> Vint

Definition at line 24 of file InductionGar.cxx.

Function Documentation

◆ Convolution_Tbranch()

TH1F Convolution_Tbranch ( double *  Input_Curr_plot_001)

Definition at line 52 of file InductionGar.cxx.

53{
54 double xmin(0), xmax(1000);
55 const int Npx(2000);
56
57 TH1F h_signal("h_signal","",Npx,xmin,xmax);
58 double Input_Time_plot_001[Npx];
59 for(int i=0; i<Npx; i++) Input_Time_plot_001[i]=h_signal.GetBinCenter(i+1);
60
61 double xmin_f1(0.), xmax_f1(1000);
62 int nStep_f1(200);
63 double step_f1=(xmax_f1-xmin_f1)/nStep_f1;
64 double x_f1(0.0);
65
66 double x[Npx];
67 double y_001[Npx];
68
69 for(int i=0; i<Npx; i++)
70 {
71 x[i]=xmin+(xmax-xmin)/Npx*(i+0.5);
72 double fun_001=0.0;
73
74 for(int j=0; j<nStep_f1; j++)
75 {
76 x_f1=xmin_f1+step_f1*(j+0.5);
77 fun_001+=rectangle2(x_f1,Input_Time_plot_001,Input_Curr_plot_001,Npx)*T_branch2(x[i]-x_f1)*step_f1;
78
79 }
80 y_001[i]=fun_001;
81 }
82
83 for(int init=0; init<Npx; init++){
84 h_signal.SetBinContent(init+1,y_001[init]);
85 }
86
87 return h_signal;
88}
Double_t x[10]
double rectangle2(double t, double *x, double *y, int Npx)
double T_branch2(double t)

Referenced by InductionGar::setMultiElectrons().

◆ rectangle2()

double rectangle2 ( double  t,
double *  x,
double *  y,
int  Npx 
)

Definition at line 29 of file InductionGar.cxx.

30{
31
32 double I = y[0];
33 int id = -1;
34 for(int i = 0; i < Npx-1; i++){
35 if(t>=x[i]&&t<x[i+1]){id = i; break;}
36 }
37 if(id != -1){I=y[id]+(y[id+1]-y[id])*(t-x[id])/(x[id+1]-x[id]);}
38
39 return I;
40}
const DifComplex I
int t()
Definition: t.c:1

Referenced by Convolution_Tbranch().

◆ T_branch2()

double T_branch2 ( double  t)

Definition at line 42 of file InductionGar.cxx.

43{
44 double result = 0.;
45 t=t-10.;
46 if(t>0) {
47 result = 2000.*(0.00181928*exp(-t/3.)-0.0147059*exp(-t/20.)+0.0128866*exp(-t/100.));
48 }
49 return result;
50}
EvtComplex exp(const EvtComplex &c)
Definition: EvtComplex.hh:252

Referenced by Convolution_Tbranch().