FairShip
Loading...
Searching...
No Matches
splitcal Class Reference

#include <splitcal.h>

Inheritance diagram for splitcal:
Collaboration diagram for splitcal:

Public Member Functions

 splitcal (const char *Name, Bool_t Active)
 
 splitcal ()
 
Bool_t ProcessHits (FairVolume *v=0) override
 
void SetZStart (Double_t ZStart)
 
void SetEmpty (Double_t Empty, Double_t BigGap, Double_t ActiveECAL_gas_gap, Double_t first_precision_layer, Double_t second_precision_layer, Double_t third_precision_layer, Double_t num_precision_layers)
 
void SetThickness (Double_t ActiveECALThickness, Double_t ActiveHCALThickness, Double_t FilterECALThickness, Double_t FilterECALThickness_first, Double_t FilterHCALThickness, Double_t ActiveECAL_gas_Thickness)
 
void SetMaterial (Double_t ActiveECALMaterial, Double_t ActiveHCALMaterial, Double_t FilterECALMaterial, Double_t FilterHCALMaterial)
 
void SetNSamplings (Int_t nECALSamplings, Int_t nHCALSamplings, Double_t ActiveHCAL)
 
void SetNModules (Int_t nModulesInX, Int_t nModulesInY)
 
void SetNStrips (Int_t nStrips)
 
void SetStripSize (Double_t stripHalfWidth, Double_t stripHalfLength)
 
void SetXMax (Double_t xMax)
 
void SetYMax (Double_t yMax)
 
void ConstructGeometry () override
 
- Public Member Functions inherited from SHiP::Detector< splitcalPoint >
 Detector ()=default
 
 Detector (const char *Name, Bool_t Active, Int_t detID)
 
 Detector (const char *Name, Bool_t Active)
 
 ~Detector () override=default
 
splitcalPointAddHit (Args &&... args)
 
void ConstructGeometry () override=0
 
void Initialize () override
 
void Reset () override
 
void EndOfEvent () override
 
void Register () override
 
TClonesArray * GetCollection (Int_t iColl) const override
 
void UpdatePointTrackIndices (const std::map< Int_t, Int_t > &indexMap) override
 Update track indices in point collection after track filtering.
 
void SetSpecialPhysicsCuts () override
 
void FinishPrimary () override
 
void FinishRun () override
 
void BeginPrimary () override
 
void PostTrack () override
 
void PreTrack () override
 
void BeginEvent () override
 
void CopyClones (TClonesArray *cl1, TClonesArray *cl2, Int_t offset) override
 
- Public Member Functions inherited from ISTLPointContainer
virtual void UpdatePointTrackIndices (const std::map< Int_t, Int_t > &indexMap)=0
 Update track indices in point collection after track filtering.
 
virtual ~ISTLPointContainer ()=default
 

Private Member Functions

 splitcal (const splitcal &)=delete
 
splitcaloperator= (const splitcal &)=delete
 

Private Attributes

Double_t fActiveECALThickness
 
Double_t fActiveHCALThickness
 
Double_t fFilterECALThickness
 
Double_t xfFilterECALThickness
 
Double_t fFilterECALThickness_first
 
Double_t fFilterHCALThickness
 
Double_t fActiveECALMaterial
 
Double_t fActiveHCALMaterial
 
Double_t fFilterECALMaterial
 
Double_t fFilterHCALMaterial
 
Double_t fActiveECAL_gas_gap
 
Double_t fActiveECAL_gas_Thickness
 
Int_t fnECALSamplings
 
Int_t fnHCALSamplings
 
Double_t fActiveHCAL
 
Double_t fZStart
 
Double_t fEmpty
 
Double_t fBigGap
 
Double_t fXMax
 
Double_t fYMax
 
Double_t ffirst_precision_layer
 
Double_t fsecond_precision_layer
 
Double_t fthird_precision_layer
 
Double_t fnum_precision_layers
 
Int_t fNModulesInX
 
Int_t fNModulesInY
 
Int_t fNStripsPerModule
 
Double_t fStripHalfWidth
 
Double_t fStripHalfLength
 

Additional Inherited Members

- Protected Attributes inherited from SHiP::Detector< splitcalPoint >
Int_t fEventID
 
Int_t fTrackID
 event index
 
Int_t fVolumeID
 track index
 
TLorentzVector fPos
 volume id
 
TLorentzVector fMom
 position at entrance
 
Double_t fTime
 momentum at entrance
 
Double_t fLength
 time
 
Double_t fELoss
 length
 
std::vector< splitcalPoint > * fDetPoints
 energy loss
 
TGeoVolume * fDetector
 

Detailed Description

Definition at line 19 of file splitcal.h.

Constructor & Destructor Documentation

◆ splitcal() [1/3]

splitcal::splitcal ( const char *  Name,
Bool_t  Active 
)

Name : Detector Name Active: kTRUE for active detectors (ProcessHits() will be called) kFALSE for inactive detectors

Definition at line 41 of file splitcal.cxx.

42 : Detector(name, active, kSplitCal) {}
@ kSplitCal

◆ splitcal() [2/3]

splitcal::splitcal ( )

default constructor

Definition at line 39 of file splitcal.cxx.

39: Detector("splitcal", kTRUE, kSplitCal) {}

◆ splitcal() [3/3]

splitcal::splitcal ( const splitcal )
privatedelete

container for data points

Member Function Documentation

◆ ConstructGeometry()

void splitcal::ConstructGeometry ( )
overridevirtual

Create the detector geometry

If you are using the standard ASCII input for the geometry just copy this and use it for your detector, otherwise you can implement here you own way of constructing the geometry.

Implements SHiP::Detector< splitcalPoint >.

Definition at line 145 of file splitcal.cxx.

145 {
150 TGeoVolume* top = gGeoManager->GetTopVolume();
151 TGeoVolume* tSplitCal = new TGeoVolumeAssembly("SplitCalDetector");
152
153 ShipGeo::InitMedium("iron");
154 ShipGeo::InitMedium("lead");
155 ShipGeo::InitMedium("Scintillator");
156 ShipGeo::InitMedium("argon");
157 ShipGeo::InitMedium("GEMmixture");
158
159 TGeoMedium* A2 = gGeoManager->GetMedium("iron");
160 TGeoMedium* A3 = gGeoManager->GetMedium("lead");
161 TGeoMedium* A4 = gGeoManager->GetMedium("GEMmixture");
162 TGeoMedium* A1 = gGeoManager->GetMedium("Scintillator");
163
164 Double_t zStartSplitCal = fZStart;
165
166 TGeoVolume* newECALfilter_first; // first layer can have different thickeness
167 TGeoVolume* newECALfilter;
168 TGeoVolume* newECALdet_gas;
169 TGeoVolume* stripGivingX;
170 TGeoVolume* stripGivingY;
171
172 TGeoVolume* newHCALfilter[100];
173 TGeoVolume* newHCALdet[100];
174 const char* char_labelHCALfilter[100];
175 TString labelHCALfilter;
176 const char* char_labelHCALdet[100];
177 TString labelHCALdet;
178
179 Double_t z_splitcal = 0;
180
181 // logical volume for the absorbing layers
182 // first absorbing layer can have different thinkens from the others
183 newECALfilter_first = gGeoManager->MakeBox(
184 "ECALfilter_first", A3, fXMax, fYMax, fFilterECALThickness_first / 2);
185 newECALfilter_first->SetLineColor(kGray);
186 newECALfilter = gGeoManager->MakeBox("ECALfilter", A3, fXMax, fYMax,
188 newECALfilter->SetLineColor(kGray);
189
190 stripGivingX =
191 gGeoManager->MakeBox("stripGivingX", A1, fStripHalfWidth,
193 stripGivingX->SetVisibility(kTRUE);
194 AddSensitiveVolume(stripGivingX);
195 stripGivingX->SetLineColor(kGreen);
196
197 stripGivingY =
198 gGeoManager->MakeBox("stripGivingY", A1, fStripHalfLength,
200 stripGivingY->SetVisibility(kTRUE);
201 AddSensitiveVolume(stripGivingY);
202 stripGivingY->SetLineColor(kGreen);
203
204 // logical volume for the high precision sensitive layers
205 newECALdet_gas = gGeoManager->MakeBox("ECALdet_gas", A4, fXMax, fYMax,
207 AddSensitiveVolume(newECALdet_gas);
208 newECALdet_gas->SetLineColor(kRed);
209
210 // now position layer volumes in tSplitCal mother volume
211 for (Int_t i_nlayECAL = 0; i_nlayECAL < fnECALSamplings; i_nlayECAL++) {
212 // position absorber layers
213 // thinkness of first layer can be different from others
214 if (i_nlayECAL == 0) {
215 z_splitcal += fFilterECALThickness_first / 2;
216 tSplitCal->AddNode(newECALfilter_first, i_nlayECAL * 1e5,
217 new TGeoTranslation(0, 0, z_splitcal));
218 z_splitcal += fFilterECALThickness_first / 2;
219 } else {
220 z_splitcal += fFilterECALThickness / 2;
221 tSplitCal->AddNode(newECALfilter, i_nlayECAL * 1e5,
222 new TGeoTranslation(0, 0, z_splitcal));
223 z_splitcal += fFilterECALThickness / 2;
224 }
225
226 if (i_nlayECAL == 0)
227 z_splitcal += fEmpty; // space after first layer? set to 0 in the
228 // config file? for whar is it for?
229 if (i_nlayECAL == 7) z_splitcal += fBigGap;
230
231 // position high precision sensitive layers
232 if (i_nlayECAL == ffirst_precision_layer ||
233 i_nlayECAL == fsecond_precision_layer ||
234 i_nlayECAL == fthird_precision_layer) {
235 z_splitcal += fActiveECAL_gas_Thickness / 2;
236
237 tSplitCal->AddNode(newECALdet_gas, 1e8 + (i_nlayECAL + 1) * 1e5,
238 new TGeoTranslation(0, 0, z_splitcal));
239 z_splitcal += fActiveECAL_gas_Thickness / 2;
240 if (fnum_precision_layers == 2) {
241 z_splitcal += fActiveECAL_gas_gap;
242 z_splitcal += fActiveECAL_gas_Thickness / 2;
243 tSplitCal->AddNode(newECALdet_gas, 1e8 + (i_nlayECAL + 1) * 1e5,
244 new TGeoTranslation(0, 0, z_splitcal));
245 z_splitcal += fActiveECAL_gas_Thickness / 2;
246 }
247 } else {
248 // position sensitive layers
249 z_splitcal += fActiveECALThickness / 2;
250 if (i_nlayECAL % 2 == 0) {
251 // strips giving x information
252 for (int mx = 0; mx < fNModulesInX; mx++) {
253 for (int my = 0; my < fNModulesInY; my++) {
254 for (int j = 0; j < fNStripsPerModule; j++) {
255 int index = (i_nlayECAL + 1) * 1e5 + (mx + 1) * 1e4 +
256 (my + 1) * 1e3 + j + 1;
257 double xCoordinate =
258 -fXMax + (fNStripsPerModule * mx + j + 0.5) *
260 2; // the times 2 is to get the total width
261 // from the half-width
262 double yCoordinate =
263 -fYMax + (my + 0.5) * fStripHalfLength *
264 2; // the times 2 is to get the total length
265 // from the half-length
266 tSplitCal->AddNode(
267 stripGivingX, index,
268 new TGeoTranslation(xCoordinate, yCoordinate, z_splitcal));
269
270 } // end loop on strips
271 } // end loop on modules in y
272 } // end loop on modules in x
273 } // end layer stripped in X
274 else {
275 // strips giving y information
276 for (int mx = 0; mx < fNModulesInX; mx++) {
277 for (int my = 0; my < fNModulesInY; my++) {
278 for (int j = 0; j < fNStripsPerModule; j++) {
279 int index = (i_nlayECAL + 1) * 1e5 + (mx + 1) * 1e4 +
280 (my + 1) * 1e3 + j + 1;
281 double xCoordinate =
282 -fXMax + (mx + 0.5) * fStripHalfLength *
283 2; // the times 2 is to get the total length
284 // from the half-length
285 double yCoordinate =
286 -fYMax + (fNStripsPerModule * my + j + 0.5) *
288 2; // the times 2 is to get the total width
289 // from the half-width
290 tSplitCal->AddNode(
291 stripGivingY, index,
292 new TGeoTranslation(xCoordinate, yCoordinate, z_splitcal));
293
294 } // end loop on strips
295 } // end loop on modules in y
296 } // end loop on modules in x
297 } // end layer stripped in Y
298
299 z_splitcal += fActiveECALThickness / 2;
300 } // end loop on sensitive scintillator layers
301
302 } // end loop in ecal layers
303
304 for (Int_t i_nlayHCAL = 0; i_nlayHCAL < 1; i_nlayHCAL++) {
305 labelHCALfilter = "HCALfilter_";
306 labelHCALfilter += i_nlayHCAL;
307 char_labelHCALfilter[i_nlayHCAL] = labelHCALfilter;
308 labelHCALdet = "HCAL_det";
309 labelHCALdet += i_nlayHCAL;
310 char_labelHCALdet[i_nlayHCAL] = labelHCALdet;
311 newHCALfilter[i_nlayHCAL] =
312 gGeoManager->MakeBox(char_labelHCALfilter[i_nlayHCAL], A2, fXMax, fYMax,
314 if (fActiveHCAL) {
315 newHCALdet[i_nlayHCAL] =
316 gGeoManager->MakeBox(char_labelHCALdet[i_nlayHCAL], A4, fXMax, fYMax,
318
319 AddSensitiveVolume(newHCALdet[i_nlayHCAL]);
320
321 newHCALdet[i_nlayHCAL]->SetLineColor(kRed);
322 }
323 newHCALfilter[i_nlayHCAL]->SetLineColor(kBlue);
324 }
325 for (Int_t i_nlayHCAL = 0; i_nlayHCAL < fnHCALSamplings; i_nlayHCAL++) {
326 z_splitcal += fFilterHCALThickness / 2;
327 tSplitCal->AddNode(newHCALfilter[i_nlayHCAL], 1,
328 new TGeoTranslation(0, 0, z_splitcal));
329 z_splitcal += fFilterHCALThickness / 2;
330
331 z_splitcal += fActiveHCALThickness / 2;
332 if (fActiveHCAL)
333 tSplitCal->AddNode(newHCALdet[i_nlayHCAL], 1,
334 new TGeoTranslation(0, 0, z_splitcal));
335 z_splitcal += fActiveHCALThickness / 2;
336 }
337
338 // finish assembly and position
339 TGeoShapeAssembly* asmb =
340 dynamic_cast<TGeoShapeAssembly*>(tSplitCal->GetShape());
341 Double_t totLength = asmb->GetDZ();
342 top->AddNode(tSplitCal, 1,
343 new TGeoTranslation(0, 0, zStartSplitCal + totLength));
344}
Double_t fFilterHCALThickness
Definition: splitcal.h:70
Double_t fnum_precision_layers
Definition: splitcal.h:81
Double_t fActiveHCALThickness
Definition: splitcal.h:69
Int_t fNModulesInY
Definition: splitcal.h:82
Double_t fFilterECALThickness_first
Definition: splitcal.h:70
Double_t fActiveECALThickness
Definition: splitcal.h:69
Double_t fthird_precision_layer
Definition: splitcal.h:81
Double_t fActiveHCAL
Definition: splitcal.h:75
Int_t fNStripsPerModule
Definition: splitcal.h:83
Int_t fNModulesInX
Definition: splitcal.h:82
Double_t fActiveECAL_gas_gap
Definition: splitcal.h:73
Double_t fBigGap
Definition: splitcal.h:77
Double_t fActiveECAL_gas_Thickness
Definition: splitcal.h:73
Double_t fStripHalfLength
Definition: splitcal.h:84
Double_t fEmpty
Definition: splitcal.h:77
Double_t fZStart
Definition: splitcal.h:76
Double_t fYMax
Definition: splitcal.h:79
Int_t fnHCALSamplings
Definition: splitcal.h:74
Double_t fXMax
Definition: splitcal.h:78
Double_t fsecond_precision_layer
Definition: splitcal.h:80
Double_t ffirst_precision_layer
Definition: splitcal.h:80
Double_t fStripHalfWidth
Definition: splitcal.h:84
Double_t fFilterECALThickness
Definition: splitcal.h:69
Int_t fnECALSamplings
Definition: splitcal.h:74
Int_t InitMedium(const char *name)
Definition: ShipGeoUtil.h:20

◆ operator=()

splitcal & splitcal::operator= ( const splitcal )
privatedelete

◆ ProcessHits()

Bool_t splitcal::ProcessHits ( FairVolume *  v = 0)
override

this method is called for each step during simulation (see FairMCApplication::Stepping())

This method is called from the MC stepping

Definition at line 44 of file splitcal.cxx.

44 {
46 // Set parameters at entrance of volume. Reset ELoss.
47 if (gMC->IsTrackEntering()) {
48 fELoss = 0.;
49 fTime = gMC->TrackTime() * 1.0e09;
50 fLength = gMC->TrackLength();
51 gMC->TrackPosition(fPos);
52 gMC->TrackMomentum(fMom);
53 }
54
55 // Sum energy loss for all steps in the active volume
56 fELoss += gMC->Edep();
57
58 // Create splitcalPoint at exit of active volume
59 if (gMC->IsTrackExiting() || gMC->IsTrackStop() ||
60 gMC->IsTrackDisappeared()) {
61 fEventID = gMC->CurrentEvent();
62 fTrackID = gMC->GetStack()->GetCurrentTrackNumber();
63 Int_t detID = 0;
64 gMC->CurrentVolID(detID);
65
67 if (fELoss == 0.) {
68 return kFALSE;
69 }
70 TParticle* p = gMC->GetStack()->GetCurrentTrack();
71 Int_t pdgCode = p->GetPdgCode();
73 TVector3(fPos.X(), fPos.Y(), fPos.Z()),
74 TVector3(fMom.Px(), fMom.Py(), fMom.Pz()), fTime, fLength, fELoss,
75 pdgCode);
76
77 // Increment number of splitcal det points in TParticle
78 ShipStack* stack = dynamic_cast<ShipStack*>(gMC->GetStack());
79 stack->AddPoint(kSplitCal);
80 }
81
82 return kTRUE;
83}
Int_t fTrackID
event index
Definition: Detector.h:85
Int_t fVolumeID
track index
Definition: Detector.h:86
TLorentzVector fPos
volume id
Definition: Detector.h:87
Double_t fTime
momentum at entrance
Definition: Detector.h:89
splitcalPoint * AddHit(Args &&... args)
Definition: Detector.h:38
TLorentzVector fMom
position at entrance
Definition: Detector.h:88
ROOT p
Definition: makeDecay.py:76

◆ SetEmpty()

void splitcal::SetEmpty ( Double_t  Empty,
Double_t  BigGap,
Double_t  ActiveECAL_gas_gap,
Double_t  first_precision_layer,
Double_t  second_precision_layer,
Double_t  third_precision_layer,
Double_t  num_precision_layers 
)

Definition at line 86 of file splitcal.cxx.

91 {
92 fEmpty = Empty;
93 fBigGap = BigGap;
94 fActiveECAL_gas_gap = ActiveECAL_gas_gap;
95 ffirst_precision_layer = first_precision_layer;
96 fsecond_precision_layer = second_precision_layer;
97 fthird_precision_layer = third_precision_layer;
98 fnum_precision_layers = num_precision_layers;
99}

◆ SetMaterial()

void splitcal::SetMaterial ( Double_t  ActiveECALMaterial,
Double_t  ActiveHCALMaterial,
Double_t  FilterECALMaterial,
Double_t  FilterHCALMaterial 
)

Definition at line 114 of file splitcal.cxx.

117 {
118 fActiveECALMaterial = ActiveECALMaterial;
119 fActiveHCALMaterial = ActiveHCALMaterial;
120 fFilterECALMaterial = FilterECALMaterial;
121 fFilterHCALMaterial = FilterHCALMaterial;
122}
Double_t fFilterHCALMaterial
Definition: splitcal.h:72
Double_t fActiveHCALMaterial
Definition: splitcal.h:71
Double_t fFilterECALMaterial
Definition: splitcal.h:71
Double_t fActiveECALMaterial
Definition: splitcal.h:71

◆ SetNModules()

void splitcal::SetNModules ( Int_t  nModulesInX,
Int_t  nModulesInY 
)

Definition at line 131 of file splitcal.cxx.

131 {
132 fNModulesInX = nModulesInX;
133 fNModulesInY = nModulesInY;
134}

◆ SetNSamplings()

void splitcal::SetNSamplings ( Int_t  nECALSamplings,
Int_t  nHCALSamplings,
Double_t  ActiveHCAL 
)

Definition at line 124 of file splitcal.cxx.

125 {
126 fnHCALSamplings = nHCALSamplings;
127 fnECALSamplings = nECALSamplings;
128 fActiveHCAL = ActiveHCAL;
129}

◆ SetNStrips()

void splitcal::SetNStrips ( Int_t  nStrips)

Definition at line 136 of file splitcal.cxx.

136{ fNStripsPerModule = nStrips; }

◆ SetStripSize()

void splitcal::SetStripSize ( Double_t  stripHalfWidth,
Double_t  stripHalfLength 
)

Definition at line 138 of file splitcal.cxx.

138 {
139 fStripHalfWidth = stripHalfWidth;
140 fStripHalfLength = stripHalfLength;
141}

◆ SetThickness()

void splitcal::SetThickness ( Double_t  ActiveECALThickness,
Double_t  ActiveHCALThickness,
Double_t  FilterECALThickness,
Double_t  FilterECALThickness_first,
Double_t  FilterHCALThickness,
Double_t  ActiveECAL_gas_Thickness 
)

Definition at line 101 of file splitcal.cxx.

106 {
107 fActiveECALThickness = ActiveECALThickness;
108 fActiveHCALThickness = ActiveHCALThickness;
109 fFilterECALThickness = FilterECALThickness;
110 fFilterECALThickness_first = FilterECALThickness_first;
111 fFilterHCALThickness = FilterHCALThickness;
112 fActiveECAL_gas_Thickness = ActiveECAL_gas_Thickness;
113}

◆ SetXMax()

void splitcal::SetXMax ( Double_t  xMax)

Definition at line 143 of file splitcal.cxx.

143{ fXMax = xMax; }

◆ SetYMax()

void splitcal::SetYMax ( Double_t  yMax)

Definition at line 144 of file splitcal.cxx.

144{ fYMax = yMax; }

◆ SetZStart()

void splitcal::SetZStart ( Double_t  ZStart)

Definition at line 85 of file splitcal.cxx.

85{ fZStart = ZStart; }

Member Data Documentation

◆ fActiveECAL_gas_gap

Double_t splitcal::fActiveECAL_gas_gap
private

Definition at line 73 of file splitcal.h.

◆ fActiveECAL_gas_Thickness

Double_t splitcal::fActiveECAL_gas_Thickness
private

Definition at line 73 of file splitcal.h.

◆ fActiveECALMaterial

Double_t splitcal::fActiveECALMaterial
private

Definition at line 71 of file splitcal.h.

◆ fActiveECALThickness

Double_t splitcal::fActiveECALThickness
private

Track information to be stored until the track leaves the active volume.

Definition at line 69 of file splitcal.h.

◆ fActiveHCAL

Double_t splitcal::fActiveHCAL
private

Definition at line 75 of file splitcal.h.

◆ fActiveHCALMaterial

Double_t splitcal::fActiveHCALMaterial
private

Definition at line 71 of file splitcal.h.

◆ fActiveHCALThickness

Double_t splitcal::fActiveHCALThickness
private

Definition at line 69 of file splitcal.h.

◆ fBigGap

Double_t splitcal::fBigGap
private

Definition at line 77 of file splitcal.h.

◆ fEmpty

Double_t splitcal::fEmpty
private

Definition at line 77 of file splitcal.h.

◆ fFilterECALMaterial

Double_t splitcal::fFilterECALMaterial
private

Definition at line 71 of file splitcal.h.

◆ fFilterECALThickness

Double_t splitcal::fFilterECALThickness
private

Definition at line 69 of file splitcal.h.

◆ fFilterECALThickness_first

Double_t splitcal::fFilterECALThickness_first
private

Definition at line 70 of file splitcal.h.

◆ fFilterHCALMaterial

Double_t splitcal::fFilterHCALMaterial
private

Definition at line 72 of file splitcal.h.

◆ fFilterHCALThickness

Double_t splitcal::fFilterHCALThickness
private

Definition at line 70 of file splitcal.h.

◆ ffirst_precision_layer

Double_t splitcal::ffirst_precision_layer
private

Definition at line 80 of file splitcal.h.

◆ fnECALSamplings

Int_t splitcal::fnECALSamplings
private

Definition at line 74 of file splitcal.h.

◆ fnHCALSamplings

Int_t splitcal::fnHCALSamplings
private

Definition at line 74 of file splitcal.h.

◆ fNModulesInX

Int_t splitcal::fNModulesInX
private

Definition at line 82 of file splitcal.h.

◆ fNModulesInY

Int_t splitcal::fNModulesInY
private

Definition at line 82 of file splitcal.h.

◆ fNStripsPerModule

Int_t splitcal::fNStripsPerModule
private

Definition at line 83 of file splitcal.h.

◆ fnum_precision_layers

Double_t splitcal::fnum_precision_layers
private

Definition at line 81 of file splitcal.h.

◆ fsecond_precision_layer

Double_t splitcal::fsecond_precision_layer
private

Definition at line 80 of file splitcal.h.

◆ fStripHalfLength

Double_t splitcal::fStripHalfLength
private

Definition at line 84 of file splitcal.h.

◆ fStripHalfWidth

Double_t splitcal::fStripHalfWidth
private

Definition at line 84 of file splitcal.h.

◆ fthird_precision_layer

Double_t splitcal::fthird_precision_layer
private

Definition at line 81 of file splitcal.h.

◆ fXMax

Double_t splitcal::fXMax
private

Definition at line 78 of file splitcal.h.

◆ fYMax

Double_t splitcal::fYMax
private

Definition at line 79 of file splitcal.h.

◆ fZStart

Double_t splitcal::fZStart
private

Definition at line 76 of file splitcal.h.

◆ xfFilterECALThickness

Double_t splitcal::xfFilterECALThickness
private

Definition at line 70 of file splitcal.h.


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