CGEM BOSS 6.6.5.h
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 "GaudiKernel/PropertyMgr.h"
#include "CLHEP/Units/PhysicalConstants.h"
#include "G4DigiManager.hh"
#include "Randomize.hh"
#include "G4ios.hh"
#include <iomanip>
#include <iostream>
#include <fstream>
#include <cmath>
#include "TTree.h"
#include "TRandom.h"
#include "TSystem.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 rectangle2_fast (double t, double *x, const double *y, int Npx)
 calc estimate using linear interpolation of curve {x,y} , at point t; mimic rectangle2() but faster version assuming x: from 0 to 1000 and 2000 points.
 
double rectangle2 (double t, std::vector< double > &x, const double *y, int Npx)
 
double T_branch2 (double t)
 
double E_branch2 (double t)
 
TH1F Convolution_Tbranch (const double *Input_Curr_plot_001, double xmin, double xmax, int Npx=2000)
 
TH1F Convolution_Tbranch (const double *Input_Curr_plot_001)
 
TH1F Convolution_Ebranch (const double *Input_Curr_plot_001, double xmin, double xmax, int Npx=2000)
 
TH1F Convolution_Ebranch (const double *Input_Curr_plot_001)
 

Typedef Documentation

◆ Vint

typedef std::vector<int> Vint

Definition at line 38 of file InductionGar.cxx.

Function Documentation

◆ Convolution_Ebranch() [1/2]

TH1F Convolution_Ebranch ( const double * Input_Curr_plot_001)

Definition at line 242 of file InductionGar.cxx.

243{
244 double xmin(0), xmax(1000);
245 const int Npx(2000);
246
247 TH1F h_signalE("h_signalE","",Npx,xmin,xmax);
248 double Input_Time_plot_001[Npx];
249 for(int i=0; i<Npx; i++) Input_Time_plot_001[i]=h_signalE.GetBinCenter(i+1);
250
251 double xmin_f1(0.), xmax_f1(1000);
252 int nStep_f1(200);
253 double step_f1=(xmax_f1-xmin_f1)/nStep_f1;
254 double x_f1(0.0);
255
256 double x[Npx];
257 double y_001[Npx];
258
259 for(int i=0; i<Npx; i++)
260 {
261 x[i]=xmin+(xmax-xmin)/Npx*(i+0.5);
262 double fun_001=0.0;
263
264 for(int j=0; j<nStep_f1; j++)
265 {
266 x_f1=xmin_f1+step_f1*(j+0.5);
267 fun_001+=rectangle2_fast(x_f1,Input_Time_plot_001,Input_Curr_plot_001,Npx)*E_branch2(x[i]-x_f1)*step_f1;
268
269 }
270 y_001[i]=fun_001;
271 }
272
273 for(int init=0; init<Npx; init++){
274 h_signalE.SetBinContent(init+1,-y_001[init]);
275 }
276
277 return h_signalE;
278}
Double_t x[10]
double rectangle2_fast(double t, double *x, const double *y, int Npx)
calc estimate using linear interpolation of curve {x,y} , at point t; mimic rectangle2() but faster v...
double E_branch2(double t)

◆ Convolution_Ebranch() [2/2]

TH1F Convolution_Ebranch ( const double * Input_Curr_plot_001,
double xmin,
double xmax,
int Npx = 2000 )

Definition at line 207 of file InductionGar.cxx.

208{
209 TH1F h_signalE("h_signalE","",Npx,xmin,xmax);
210 std::vector<double> Input_Time_plot_001(Npx);
211 for(int i=0; i<Npx; i++) Input_Time_plot_001[i]=h_signalE.GetBinCenter(i+1);
212
213 double xmin_f1(0.), xmax_f1(1000);
214 int nStep_f1(200);
215 double step_f1=(xmax_f1-xmin_f1)/nStep_f1;
216 double x_f1(0.0);
217
218 std::vector<double> x(Npx);
219 std::vector<double> y_001(Npx);
220
221 for(int i=0; i<Npx; i++)
222 {
223 x[i]=xmin+(xmax-xmin)/Npx*(i+0.5);
224 double fun_001=0.0;
225
226 for(int j=0; j<nStep_f1; j++)
227 {
228 x_f1=xmin_f1+step_f1*(j+0.5);
229 fun_001+=rectangle2(x_f1,Input_Time_plot_001,Input_Curr_plot_001,Npx)*E_branch2(x[i]-x_f1)*step_f1;
230
231 }
232 y_001[i]=fun_001;
233 }
234
235 for(int init=0; init<Npx; init++){
236 h_signalE.SetBinContent(init+1,-y_001[init]);
237 }
238
239 return h_signalE;
240}
double rectangle2(double t, double *x, double *y, int Npx)

◆ Convolution_Tbranch() [1/2]

TH1F Convolution_Tbranch ( const double * Input_Curr_plot_001)

Definition at line 169 of file InductionGar.cxx.

170{
171 double xmin(0), xmax(1000);
172 const int Npx(2000);
173
174 TH1F h_signalT("h_signalT","",Npx,xmin,xmax);
175 double Input_Time_plot_001[Npx];
176 for(int i=0; i<Npx; i++) Input_Time_plot_001[i]=h_signalT.GetBinCenter(i+1);
177
178 double xmin_f1(0.), xmax_f1(1000);
179 int nStep_f1(200);
180 double step_f1=(xmax_f1-xmin_f1)/nStep_f1;
181 double x_f1(0.0);
182
183 double x[Npx];
184 double y_001[Npx];
185
186 for(int i=0; i<Npx; i++)
187 {
188 x[i]=xmin+(xmax-xmin)/Npx*(i+0.5);
189 double fun_001=0.0;
190
191 for(int j=0; j<nStep_f1; j++)
192 {
193 x_f1=xmin_f1+step_f1*(j+0.5);
194 fun_001+=rectangle2_fast(x_f1,Input_Time_plot_001,Input_Curr_plot_001,Npx)*T_branch2(x[i]-x_f1)*step_f1;
195
196 }
197 y_001[i]=fun_001;
198 }
199
200 for(int init=0; init<Npx; init++){
201 h_signalT.SetBinContent(init+1,-y_001[init]);
202 }
203
204 return h_signalT;
205}
double T_branch2(double t)

◆ Convolution_Tbranch() [2/2]

TH1F Convolution_Tbranch ( const double * Input_Curr_plot_001,
double xmin,
double xmax,
int Npx = 2000 )

Definition at line 133 of file InductionGar.cxx.

134{
135
136 TH1F h_signalT("h_signalT","",Npx,xmin,xmax);
137 std::vector<double> Input_Time_plot_001(Npx);
138 for(int i=0; i<Npx; i++) Input_Time_plot_001[i]=h_signalT.GetBinCenter(i+1);
139
140 double xmin_f1(0.), xmax_f1(1000);
141 int nStep_f1(200);
142 double step_f1=(xmax_f1-xmin_f1)/nStep_f1;
143 double x_f1(0.0);
144
145 std::vector<double> x(Npx);
146 std::vector<double> y_001(Npx);
147
148 for(int i=0; i<Npx; i++)
149 {
150 x[i]=xmin+(xmax-xmin)/Npx*(i+0.5);
151 double fun_001=0.0;
152
153 for(int j=0; j<nStep_f1; j++)
154 {
155 x_f1=xmin_f1+step_f1*(j+0.5);
156 fun_001+=rectangle2(x_f1,Input_Time_plot_001,Input_Curr_plot_001,Npx)*T_branch2(x[i]-x_f1)*step_f1;
157
158 }
159 y_001[i]=fun_001;
160 }
161
162 for(int init=0; init<Npx; init++){
163 h_signalT.SetBinContent(init+1,-y_001[init]);
164 }
165
166 return h_signalT;
167}

◆ E_branch2()

double E_branch2 ( double t)

Definition at line 123 of file InductionGar.cxx.

124{
125 double result = 0.;
126 t=t-5;
127 if(t>0) {
128 result = 0.000627357*(1358.7*exp(-t*0.0385647)*t+1358.7*exp(-t*0.0114353)*t+100164.*exp(-t*0.0385647)-100164.*exp(-0.0114353*t));
129 }
130 return result;
131}
EvtComplex exp(const EvtComplex &c)
int t()
Definition t.c:1

Referenced by CEF::CEF(), CEF2::CEF2(), Convolution_Ebranch(), and Convolution_Ebranch().

◆ rectangle2() [1/2]

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

Definition at line 44 of file InductionGar.cxx.

45{
46
47 double I = y[0];
48 int id = -1;
49 for(int i = 0; i < Npx-1; i++){
50 if(t>=x[i]&&t<x[i+1]){id = i; break;}
51 }
52 if(id != -1){I=y[id]+(y[id+1]-y[id])*(t-x[id])/(x[id+1]-x[id]);}
53
54 return I;
55}
const DifComplex I

Referenced by Convolution_Ebranch(), and Convolution_Tbranch().

◆ rectangle2() [2/2]

double rectangle2 ( double t,
std::vector< double > & x,
const double * y,
int Npx )

Definition at line 100 of file InductionGar.cxx.

101{
102
103 double I = y[0];
104 int id = -1;
105 for(int i = 0; i < Npx-1; i++){
106 if(t>=x[i]&&t<x[i+1]){id = i; break;}
107 }
108 if(id != -1){I=y[id]+(y[id+1]-y[id])*(t-x[id])/(x[id+1]-x[id]);}
109
110 return I;
111}

◆ rectangle2_fast()

double rectangle2_fast ( double t,
double * x,
const double * y,
int Npx )

calc estimate using linear interpolation of curve {x,y} , at point t; mimic rectangle2() but faster version assuming x: from 0 to 1000 and 2000 points.

Parameters
t
x
y
Npx
Returns
double

Definition at line 66 of file InductionGar.cxx.

66 {
67
68 const int nbinsx=2000;
69 const double xmin=0;
70 const double xmax=1000;
71 const double binWidth=(xmax-xmin)/nbinsx;
72 double I = y[0];
73 int id = -1;
74
75 if (t > xmin+binWidth/2 && t< xmax-binWidth/2){
76 id=(int)floor((t-binWidth/2-xmin)/binWidth);
77 }
78#ifdef INDUCTIONGAR_DEBUG_RECTANGLE2_FAST_TEST // checking results of new and old version. for speed it should not run
79 int id2 = -1;
80 for (int i = 0; i < Npx - 1; i++) {
81 if (t >= x[i] && t < x[i + 1]) {
82 id2 = i;
83 break;
84 } // find the closist id, to make t=x[id]; find root.
85 }
86 if (id2 != id) {
87 cout<<
88 "CgemDigitizerSvc::InductionGar::rectangle2_fast: 2 method calced id doesnot correspond, "<<
89 "id new : "<< id2 <<" id old : "<<id<<
90 " t : "<<t<<" x at id old "<< x[id2]<< endl;
91 }
92#endif
93 if (id != -1) {
94 I = y[id] + (y[id + 1] - y[id]) * (t - x[id]) / (x[id + 1] - x[id]);
95 } // the tangent of y(i),x(i) at i=id, I is on tangent and at x=t;
96
97 return I;
98}

Referenced by Convolution_Ebranch(), and Convolution_Tbranch().

◆ T_branch2()

double T_branch2 ( double t)

Definition at line 113 of file InductionGar.cxx.

114{
115 double result = 0.;
116 t=t-10.;
117 if(t>0) {
118 result = 2000.*(0.00181928*exp(-t/3.)-0.0147059*exp(-t/20.)+0.0128866*exp(-t/100.));
119 }
120 return result;
121}

Referenced by Convolution_Tbranch(), Convolution_Tbranch(), CTF::CTF(), and CTF2::CTF2().