 |
ATLAS Offline Software
|
|
__DEVICE__ long long | getDDE (GeoGpu *geo, int sampling, float eta, float phi) |
|
__DEVICE__ int | find_index_f (float *array, int size, float value) |
|
__DEVICE__ int | find_index_uint32 (uint32_t *array, int size, uint32_t value) |
|
__DEVICE__ void | rnd_to_fct2d (float &valuex, float &valuey, float rnd0, float rnd1, FH2D *hf2d) |
|
__DEVICE__ float | rnd_to_fct1d (float rnd, uint32_t *contents, float *borders, int nbins, uint32_t s_MaxValue) |
|
__DEVICE__ void | CenterPositionCalculation_d (Hit &hit, const Chain0_Args args) |
|
__DEVICE__ void | HistoLateralShapeParametrization_d (Hit &hit, unsigned long t, Chain0_Args args, bool reweight) |
|
__DEVICE__ void | HitCellMapping_d (Hit &hit, unsigned long, Chain0_Args args) |
|
__DEVICE__ void | HitCellMappingWiggle_d (Hit &hit, Chain0_Args args, unsigned long t) |
|
◆ CenterPositionCalculation_d()
◆ find_index_f()
__DEVICE__ int CaloGpuGeneral_fnc::find_index_f |
( |
float * |
array, |
|
|
int |
size, |
|
|
float |
value |
|
) |
| |
Definition at line 76 of file CaloGpuGeneral_fnc.cxx.
80 int m_index = ( high - low ) / 2;
81 while ( high != low ) {
86 m_index = ( high + low ) / 2;
◆ find_index_uint32()
__DEVICE__ int CaloGpuGeneral_fnc::find_index_uint32 |
( |
uint32_t * |
array, |
|
|
int |
size, |
|
|
uint32_t |
value |
|
) |
| |
Definition at line 93 of file CaloGpuGeneral_fnc.cxx.
97 int m_index = ( high - low ) / 2;
98 while ( high != low ) {
105 m_index = ( high - low ) / 2 + low;
◆ getDDE()
__DEVICE__ long long CaloGpuGeneral_fnc::getDDE |
( |
GeoGpu * |
geo, |
|
|
int |
sampling, |
|
|
float |
eta, |
|
|
float |
phi |
|
) |
| |
Definition at line 16 of file CaloGpuGeneral_fnc.cxx.
21 int MAX_SAMPLING =
geo->max_sample;
25 if ( sampling < 0 )
return -1;
26 if ( sampling >= MAX_SAMPLING )
return -1;
28 int sample_size = SampleIdx[sampling].
size;
29 unsigned int sample_index = SampleIdx[sampling].
index;
32 if ( sample_size == 0 )
return -1;
34 long long bestDDE = -1;
40 beststeps = ( *steps );
44 if ( sampling < 21 ) {
45 for (
int skip_range_check = 0; skip_range_check <= 1; ++skip_range_check ) {
46 for (
unsigned int j = sample_index; j < sample_index + sample_size; ++j ) {
47 if ( !skip_range_check ) {
48 if (
eta <
gr[j].mineta() )
continue;
49 if (
eta >
gr[j].maxeta() )
continue;
52 intsteps = ( *steps );
56 long long newDDE =
gr[j].getDDE(
eta,
phi, &newdist, &intsteps );
60 if (
steps ) beststeps = intsteps;
61 if ( newdist < -0.1 )
break;
64 if ( bestDDE >= 0 )
break;
◆ HistoLateralShapeParametrization_d()
__DEVICE__ void CaloGpuGeneral_fnc::HistoLateralShapeParametrization_d |
( |
Hit & |
hit, |
|
|
unsigned long |
t, |
|
|
Chain0_Args |
args, |
|
|
bool |
reweight |
|
) |
| |
Definition at line 170 of file CaloGpuGeneral_fnc.cxx.
181 float alpha,
r, rnd1, rnd2;
185 if (
args.is_phi_symmetric ) {
204 if ( center_eta < 0. ) delta_eta_mm = -delta_eta_mm;
207 if (
charge < 0. ) delta_phi_mm = -delta_phi_mm;
209 float dist000 = sqrt( center_r * center_r + center_z * center_z );
210 float eta_jakobi = abs( 2.0 *
exp( -center_eta ) / ( 1.0 +
exp( -2 * center_eta ) ) );
212 float delta_eta = delta_eta_mm / eta_jakobi / dist000;
213 float delta_phi = delta_phi_mm / center_r;
216 float* histocontents = (
args.fh1d)->h_contents;
217 float* histoborders = (
args.fh1d)->h_borders;
218 float* histoerrors= (
args.fh1d)->h_errors;
220 float delta_r_mm = sqrt( delta_eta_mm * delta_eta_mm + delta_phi_mm * delta_phi_mm );
222 if( ibin < 1 ) ibin = 1 ;
224 float weight = histocontents[ibin];
225 float rms = histoerrors[ibin];
227 float rnd1 =
args.rand[
t + 2 *
args.nhits -100 ];
228 float rnd2 =
args.rand[
t + 2 *
args.nhits -10 ];
◆ HitCellMapping_d()
Definition at line 240 of file CaloGpuGeneral_fnc.cxx.
244 if ( cellele < 0 ) printf(
"cellele not found %lld \n", cellele );
246 atomicAdd( &
args.cells_energy[cellele], hit.
E() );
◆ HitCellMappingWiggle_d()
Definition at line 251 of file CaloGpuGeneral_fnc.cxx.
253 int nhist = ( *(
args.fhs ) ).nhist;
254 float* bin_low_edge = ( *(
args.fhs ) ).low_edge;
256 float eta = fabs( hit.
eta() );
260 for (
int i = 0;
i < nhist + 1; ++
i ) {
261 if ( bin_low_edge[
i] >
eta ) {
269 unsigned int mxsz =
args.fhs->mxsz;
271 float* borders = &(
args.fhs->d_borders1D[
bin * mxsz] );
272 int h_size = ( *(
args.fhs ) ).h_szs[
bin];
275 float rnd =
args.rand[
t + 2 *
args.nhits];
279 float hit_phi_shifted = hit.
phi() + wiggle;
◆ rnd_to_fct1d()
__DEVICE__ float CaloGpuGeneral_fnc::rnd_to_fct1d |
( |
float |
rnd, |
|
|
uint32_t * |
contents, |
|
|
float * |
borders, |
|
|
int |
nbins, |
|
|
uint32_t |
s_MaxValue |
|
) |
| |
Definition at line 139 of file CaloGpuGeneral_fnc.cxx.
141 uint32_t int_rnd = s_MaxValue * rnd;
148 if ( ibin > 0 ) basecont =
contents[ibin - 1];
152 return borders[binx] + ( ( borders[binx + 1] - borders[binx] ) * ( int_rnd - basecont ) ) / dcont;
154 return borders[binx] + ( borders[binx + 1] - borders[binx] ) / 2;
◆ rnd_to_fct2d()
__DEVICE__ void CaloGpuGeneral_fnc::rnd_to_fct2d |
( |
float & |
valuex, |
|
|
float & |
valuey, |
|
|
float |
rnd0, |
|
|
float |
rnd1, |
|
|
FH2D * |
hf2d |
|
) |
| |
Definition at line 112 of file CaloGpuGeneral_fnc.cxx.
114 int nbinsx = ( *hf2d ).nbinsx;
115 int nbinsy = ( *hf2d ).nbinsy;
116 float* HistoContents = ( *hf2d ).h_contents;
117 float* HistoBorders = ( *hf2d ).h_bordersx;
118 float* HistoBordersy = ( *hf2d ).h_bordersy;
120 int ibin =
find_index_f( HistoContents, nbinsx * nbinsy, rnd0 );
122 int biny = ibin / nbinsx;
123 int binx = ibin - nbinsx * biny;
126 if ( ibin > 0 ) basecont = HistoContents[ibin - 1];
128 float dcont = HistoContents[ibin] - basecont;
130 valuex = HistoBorders[binx] + ( HistoBorders[binx + 1] - HistoBorders[binx] ) * ( rnd0 - basecont ) / dcont;
132 valuex = HistoBorders[binx] + ( HistoBorders[binx + 1] - HistoBorders[binx] ) / 2;
134 valuey = HistoBordersy[biny] + ( HistoBordersy[biny + 1] - HistoBordersy[biny] ) * rnd1;
Scalar phi() const
phi method
Scalar eta() const
pseudorapidity method
__DEVICE__ int find_index_uint32(uint32_t *array, int size, uint32_t value)
__DEVICE__ void HitCellMapping_d(Hit &hit, unsigned long, Chain0_Args args)
CUDA_HOSTDEV float & center_eta()
__HOSTDEV__ double Phi_mpi_pi(double)
CUDA_HOSTDEV float & center_z()
CUDA_HOSTDEV float & center_phi()
__DEVICE__ long long getDDE(GeoGpu *geo, int sampling, float eta, float phi)
void contents(std::vector< std::string > &keys, TDirectory *td, const std::string &directory, const std::string &pattern, const std::string &path)
CUDA_HOSTDEV float & phi()
double charge(const T &p)
__DEVICE__ void rnd_to_fct2d(float &valuex, float &valuey, float rnd0, float rnd1, FH2D *hf2d)
CUDA_HOSTDEV void setCenter_phi(float phi)
__DEVICE__ int find_index_f(float *array, int size, float value)
CUDA_HOSTDEV void setEtaPhiZE(float eta, float phi, float z, float E)
def delta_phi(phi1, phi2)
CUDA_HOSTDEV void setCenter_eta(float eta)
CUDA_HOSTDEV float & eta()
CUDA_HOSTDEV void setCenter_r(float r)
__DEVICE__ float rnd_to_fct1d(float rnd, uint32_t *contents, float *borders, int nbins, uint32_t s_MaxValue)
CUDA_HOSTDEV void setCenter_z(float z)
float distance(const Amg::Vector3D &p1, const Amg::Vector3D &p2)
calculates the distance between two point in 3D space
CUDA_HOSTDEV float & center_r()