1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
16 ////////////////////////////////////////////////////////////
17 // Factory for muon chambers, segmentations and response //
18 ////////////////////////////////////////////////////////////
22 Revision 1.2 2001/05/16 14:57:17 alibrary
23 New files for folders and Stack
25 Revision 1.1 2001/04/06 11:24:43 morsch
26 Dependency on implementations of AliSegmentation and AliMUONResponse moved to AliMUONFactory class.
27 Static method Build() builds the MUON system out of chambers, segmentation and response.
30 #include "AliMUONFactory.h"
32 #include "AliMUONChamber.h"
33 #include "AliMUONResponseV0.h"
34 #include "AliMUONResponseTrigger.h"
35 #include "AliMUONSegmentationV0.h"
36 #include "AliMUONSegmentationV01.h"
37 #include "AliMUONSegmentationV02.h"
38 #include "AliMUONSegmentationV04.h"
39 #include "AliMUONSegmentationV05.h"
40 #include "AliMUONSegmentationSlat.h"
41 #include "AliMUONSegmentationSlatN.h"
42 #include "AliMUONSegmentationTrigger.h"
43 #include "AliMUONSegmentationTriggerX.h"
44 #include "AliMUONSegmentationTriggerY.h"
46 ClassImp(AliMUONFactory)
49 void AliMUONFactory::Build(AliMUON* where, const char* what)
52 // Construct MUON from chambers, segmentation and responses
55 AliMUON* pMUON = where;
58 if (strcmp(tmp, "default")==0) {
59 if(pMUON->GetDebug()) {
61 printf("\nAliMUONFactory: --------AliMUONFactory------------------------------");
62 printf("\nAliMUONFactory: Non default version of MUON selected ");
63 printf("\nAliMUONFactory: You have to construct yourself the MUON elements !!");
64 printf("\nAliMUONFactory: ----------------------------------------------------");
67 pMUON->SetMaxStepGas(0.1);
68 pMUON->SetMaxStepAlu(0.1);
72 // First define the number of planes that are segmented (1 or 2) by a call
74 // Then chose for each chamber (chamber plane) the segmentation
75 // and response model.
76 // They should be equal for the two chambers of each station. In a future
77 // version this will be enforced.
81 // Default response: 5 mm of gas
82 AliMUONResponseV0* response0 = new AliMUONResponseV0;
83 response0->SetSqrtKx3AndDeriveKx2Kx4(0.7131); // sqrt(0.5085)
84 response0->SetSqrtKy3AndDeriveKy2Ky4(0.7642); // sqrt(0.5840)
85 response0->SetPitch(0.25); // anode-cathode distance
86 response0->SetSigmaIntegration(10.);
87 response0->SetChargeSlope(50);
88 response0->SetChargeSpread(0.18, 0.18);
89 response0->SetMaxAdc(4096);
90 response0->SetZeroSuppression(6);
92 // Response for 4 mm of gas (station 1)
93 // automatic consistency with width of sensitive medium in CreateGeometry ????
94 AliMUONResponseV0* responseSt1 = new AliMUONResponseV0;
95 // Mathieson parameters from L.Kharmandarian's thesis, page 190
96 responseSt1->SetSqrtKx3AndDeriveKx2Kx4(0.7000); // sqrt(0.4900)
97 responseSt1->SetSqrtKy3AndDeriveKy2Ky4(0.7550); // sqrt(0.5700)
98 responseSt1->SetPitch(0.20); // anode-cathode distance
99 responseSt1->SetSigmaIntegration(10.);
100 // ChargeSlope larger to compensate for the smaller anode-cathode distance
101 // and keep the same most probable ADC channel for mip's
102 responseSt1->SetChargeSlope(62.5);
103 // assumed proportionality to anode-cathode distance for ChargeSpread
104 responseSt1->SetChargeSpread(0.144, 0.144);
105 responseSt1->SetMaxAdc(4096);
106 responseSt1->SetZeroSuppression(6);
108 //--------------------------------------------------------
109 // Configuration for Chamber TC1/2 (Station 1) ----------
110 //^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
111 Float_t rseg1[4]={17.5, 55.2, 71.3, 95.5};
112 Int_t nseg1[4]={4, 4, 2, 1};
116 pMUON->SetNsec(chamber-1,2);
118 AliMUONSegmentationV01 *seg11=new AliMUONSegmentationV01(4);
120 seg11->SetSegRadii(rseg1);
121 seg11->SetPadSize(2.4, 0.4); // smaller pad size
122 seg11->SetDAnod(0.20); // smaller distance between anode wires
123 seg11->SetPadDivision(nseg1);
125 pMUON->SetSegmentationModel(chamber-1, 1, seg11);
127 AliMUONSegmentationV02 *seg12=new AliMUONSegmentationV02(4);
128 seg12->SetSegRadii(rseg1);
129 seg12->SetPadSize(0.6, 1.6); // smaller pad size
130 seg12->SetDAnod(0.20); // smaller distance between anode wires
131 seg12->SetPadDivision(nseg1);
133 pMUON->SetSegmentationModel(chamber-1, 2, seg12);
135 pMUON->SetResponseModel(chamber-1, responseSt1); // special response
136 pMUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
141 pMUON->SetNsec(chamber-1,2);
143 AliMUONSegmentationV01 *seg21=new AliMUONSegmentationV01(4);
144 seg21->SetSegRadii(rseg1);
145 seg21->SetPadSize(2.4, 0.4); // smaller pad size
146 seg21->SetDAnod(0.20); // smaller distance between anode wires
147 seg21->SetPadDivision(nseg1);
148 pMUON->SetSegmentationModel(chamber-1, 1, seg21);
150 AliMUONSegmentationV02 *seg22=new AliMUONSegmentationV02(4);
151 seg22->SetSegRadii(rseg1);
152 seg22->SetPadSize(0.6, 1.6); // smaller pad size
153 seg22->SetDAnod(0.20); // smaller distance between anode wires
154 seg22->SetPadDivision(nseg1);
155 pMUON->SetSegmentationModel(chamber-1, 2, seg22);
157 pMUON->SetResponseModel(chamber-1, responseSt1); // special response
158 pMUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
161 //--------------------------------------------------------
162 // Configuration for Chamber TC3/4 (Station 2) -----------
163 ///^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
164 // Float_t rseg2[4]={23.5, 87.7, 122.4, 122.5};
165 Float_t rseg2[4]={23.5, 53.5, 90.5, 122.5};
166 Int_t nseg2[4]={4, 4, 2, 1};
170 pMUON->SetNsec(chamber-1,2);
172 AliMUONSegmentationV01 *seg31=new AliMUONSegmentationV01(4);
173 seg31->SetSegRadii(rseg2);
174 seg31->SetPadSize(3.0, 0.5);
175 seg31->SetDAnod(3.0/3./4);
176 seg31->SetPadDivision(nseg2);
177 pMUON->SetSegmentationModel(chamber-1, 1, seg31);
179 AliMUONSegmentationV02 *seg32=new AliMUONSegmentationV02(4);
180 seg32->SetSegRadii(rseg2);
181 seg32->SetPadSize(0.75, 2.0);
182 seg32->SetPadDivision(nseg2);
183 seg32->SetDAnod(3.0/3./4);
185 pMUON->SetSegmentationModel(chamber-1, 2, seg32);
187 pMUON->SetResponseModel(chamber-1, response0);
188 pMUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
193 pMUON->SetNsec(chamber-1,2);
195 AliMUONSegmentationV01 *seg41=new AliMUONSegmentationV01(4);
196 seg41->SetSegRadii(rseg2);
197 seg41->SetPadSize(3.0, 0.5);
198 seg41->SetDAnod(3.0/3./4);
199 seg41->SetPadDivision(nseg2);
200 pMUON->SetSegmentationModel(chamber-1, 1, seg41);
202 AliMUONSegmentationV02 *seg42=new AliMUONSegmentationV02(4);
203 seg42->SetSegRadii(rseg2);
204 seg42->SetPadSize(0.75, 2.0);
205 seg42->SetPadDivision(nseg2);
206 seg42->SetDAnod(3.0/3./4);
208 pMUON->SetSegmentationModel(chamber-1, 2, seg42);
210 pMUON->SetResponseModel(chamber-1, response0);
211 pMUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
214 //--------------------------------------------------------
215 // Configuration for Chamber TC5/6 (Station 3) ----------
216 //^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
217 Int_t nseg3[4]={4, 4, 2, 1};
218 Int_t npcb5[36] = {0,0,2,0,
228 Float_t shift = 1.5/2.;
229 Float_t xpos5[9] = {2., 2., 2., 2.,33., 2., 2., 2., 2.};
230 Float_t ypos5 = -(20.+4.*(40.-2.*shift));
233 pMUON->SetNsec(chamber-1,2);
234 AliMUONSegmentationSlat *seg51=new AliMUONSegmentationSlat(4);
236 seg51->SetShift(shift);
237 seg51->SetNPCBperSector(npcb5);
238 seg51->SetSlatXPositions(xpos5);
239 seg51->SetSlatYPosition(ypos5);
240 seg51->SetPadSize(10.,0.5);
241 seg51->SetDAnod(0.25);
242 seg51->SetPadDivision(nseg3);
243 pMUON->SetSegmentationModel(chamber-1, 1, seg51);
245 AliMUONSegmentationSlatN *seg52=new AliMUONSegmentationSlatN(4);
247 seg52->SetShift(shift);
248 seg52->SetNPCBperSector(npcb5);
249 seg52->SetSlatXPositions(xpos5);
250 seg52->SetSlatYPosition(ypos5);
251 seg52->SetPadSize(1., 10.); // DeltaX(non bending) = 2 * DeltaY(bending)
252 seg52->SetDAnod(0.25);
253 seg52->SetPadDivision(nseg3);
254 pMUON->SetSegmentationModel(chamber-1, 2, seg52);
255 pMUON->SetResponseModel(chamber-1, response0);
256 pMUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
259 pMUON->SetNsec(chamber-1,2);
260 AliMUONSegmentationSlat *seg61=new AliMUONSegmentationSlat(4);
262 seg61->SetShift(shift);
263 seg61->SetNPCBperSector(npcb5);
264 seg61->SetSlatXPositions(xpos5);
265 seg61->SetSlatYPosition(ypos5);
266 seg61->SetPadSize(10.,0.5);
267 seg61->SetDAnod(0.25);
268 seg61->SetPadDivision(nseg3);
269 pMUON->SetSegmentationModel(chamber-1, 1, seg61);
271 AliMUONSegmentationSlatN *seg62=new AliMUONSegmentationSlatN(4);
273 seg62->SetShift(shift);
274 seg62->SetNPCBperSector(npcb5);
275 seg62->SetSlatXPositions(xpos5);
276 seg62->SetSlatYPosition(ypos5);
277 seg62->SetPadSize(1., 10.); // DeltaX(non bending) = 2 * DeltaY(bending)
278 seg62->SetDAnod(0.25);
279 seg62->SetPadDivision(nseg3);
280 pMUON->SetSegmentationModel(chamber-1, 2, seg62);
281 pMUON->SetResponseModel(chamber-1, response0);
282 pMUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
284 //--------------------------------------------------------
285 // Configuration for Chamber TC7/8 (Station 4) ----------
286 //^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
288 Int_t nseg4[4]={4, 4, 2, 1};
293 pMUON->SetNsec(chamber-1,2);
295 AliMUONSegmentationSlat *seg71=new AliMUONSegmentationSlat(4);
296 Int_t npcb7[44] = {0,0,0,3,
307 Float_t xpos7[11] = {2., 2., 2., 2., 2., 40.5, 2., 2., 2., 2., 2.};
308 Float_t ypos7 = -(20.+5.*(40.-2.*shift));
310 seg71->SetNSlats(11);
311 seg71->SetShift(shift);
312 seg71->SetNPCBperSector(npcb7);
313 seg71->SetSlatXPositions(xpos7);
314 seg71->SetSlatYPosition(ypos7);
316 seg71->SetPadSize(10.,0.5);
317 seg71->SetDAnod(0.25);
318 seg71->SetPadDivision(nseg4);
319 pMUON->SetSegmentationModel(chamber-1, 1, seg71);
321 AliMUONSegmentationSlatN *seg72=new AliMUONSegmentationSlatN(4);
323 pMUON->SetSegmentationModel(chamber-1, 2, seg72);
324 seg72->SetNSlats(11);
325 seg72->SetShift(shift);
326 seg72->SetNPCBperSector(npcb7);
327 seg72->SetSlatXPositions(xpos7);
328 seg72->SetSlatYPosition(ypos7);
329 seg72->SetPadSize(1., 10.); // DeltaX(non bending) = 2 * DeltaY(bending)
330 seg72->SetDAnod(0.25);
331 seg72->SetPadDivision(nseg4);
333 pMUON->SetResponseModel(chamber-1, response0);
334 pMUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
338 pMUON->SetNsec(chamber-1,2);
340 AliMUONSegmentationSlat *seg81=new AliMUONSegmentationSlat(4);
342 seg81->SetNSlats(11);
343 seg81->SetShift(shift);
344 seg81->SetNPCBperSector(npcb7);
345 seg81->SetSlatXPositions(xpos7);
346 seg81->SetSlatYPosition(ypos7);
347 seg81->SetPadSize(10.,0.5);
348 seg81->SetDAnod(0.25);
349 seg81->SetPadDivision(nseg4);
350 pMUON->SetSegmentationModel(chamber-1, 1, seg81);
352 AliMUONSegmentationSlat *seg82=new AliMUONSegmentationSlatN(4);
354 pMUON->SetSegmentationModel(chamber-1, 2, seg82);
355 seg82->SetNSlats(11);
356 seg82->SetShift(shift);
357 seg82->SetNPCBperSector(npcb7);
358 seg82->SetSlatXPositions(xpos7);
359 seg82->SetSlatYPosition(ypos7);
360 seg82->SetPadSize(1., 10.); // DeltaX(non bending) = 2 * DeltaY(bending)
361 seg82->SetDAnod(0.25);
362 seg82->SetPadDivision(nseg4);
364 pMUON->SetResponseModel(chamber-1, response0);
365 pMUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
368 //--------------------------------------------------------
369 // Configuration for Chamber TC9/10 (Station 5) ---------
370 //^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
374 pMUON->SetNsec(chamber-1,2);
376 AliMUONSegmentationSlat *seg91=new AliMUONSegmentationSlat(4);
377 Int_t npcb9[52] = {0,0,0,3,
391 Float_t xpos9[13] = {2., 2., 2., 2., 2., 2., 40.5, 2., 2., 2., 2., 2., 2.};
392 Float_t ypos9 = -(20.+6.*(40.-2.*shift));
394 seg91->SetNSlats(13);
395 seg91->SetShift(shift);
396 seg91->SetNPCBperSector(npcb9);
397 seg91->SetSlatXPositions(xpos9);
398 seg91->SetSlatYPosition(ypos9);
399 seg91->SetPadSize(10.,0.5);
400 seg91->SetDAnod(0.25);
401 seg91->SetPadDivision(nseg4);
402 pMUON->SetSegmentationModel(chamber-1, 1, seg91);
404 AliMUONSegmentationSlatN *seg92=new AliMUONSegmentationSlatN(4);
406 pMUON->SetSegmentationModel(chamber-1, 2, seg92);
407 seg92->SetNSlats(13);
408 seg92->SetShift(shift);
409 seg92->SetNPCBperSector(npcb9);
410 seg92->SetSlatXPositions(xpos9);
411 seg92->SetSlatYPosition(ypos9);
412 seg92->SetPadSize(1., 10.); // DeltaX(non bending) = 2 * DeltaY(bending)
413 seg92->SetDAnod(0.25);
414 seg92->SetPadDivision(nseg4);
416 pMUON->SetResponseModel(chamber-1, response0);
417 pMUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
421 pMUON->SetNsec(chamber-1,2);
423 AliMUONSegmentationSlat *seg101=new AliMUONSegmentationSlat(4);
425 seg101->SetNSlats(13);
426 seg101->SetShift(shift);
427 seg101->SetNPCBperSector(npcb9);
428 seg101->SetSlatXPositions(xpos9);
429 seg101->SetSlatYPosition(ypos9);
430 seg101->SetPadSize(10.,0.5);
431 seg101->SetDAnod(0.25);
432 seg101->SetPadDivision(nseg4);
433 pMUON->SetSegmentationModel(chamber-1, 1, seg101);
435 AliMUONSegmentationSlatN *seg102=new AliMUONSegmentationSlatN(4);
437 pMUON->SetSegmentationModel(chamber-1, 2, seg102);
438 seg102->SetNSlats(13);
439 seg102->SetShift(shift);
440 seg102->SetNPCBperSector(npcb9);
441 seg102->SetSlatXPositions(xpos9);
442 seg102->SetSlatYPosition(ypos9);
443 seg102->SetPadSize(1., 10.); // DeltaX(non bending) = 2 * DeltaY(bending)
444 seg102->SetDAnod(0.25);
445 seg102->SetPadDivision(nseg4);
447 pMUON->SetResponseModel(chamber-1, response0);
448 pMUON->Chamber(chamber-1).SetChargeCorrel(0.11); // 11% charge spread
450 //--------------------------------------------------------
451 // Configuration for Trigger Stations --------------------
452 //^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
454 AliMUONResponseTrigger* responseTrigger0 = new AliMUONResponseTrigger;
456 // AliMUONResponseTriggerV1* responseTrigger0 = new AliMUONResponseTriggerV1;
459 pMUON->SetNsec(chamber-1,2);
460 AliMUONSegmentationTriggerX *seg111=new AliMUONSegmentationTriggerX;
461 pMUON->SetSegmentationModel(chamber-1, 1, seg111);
462 AliMUONSegmentationTriggerY *seg112=new AliMUONSegmentationTriggerY;
463 pMUON->SetSegmentationModel(chamber-1, 2, seg112);
465 pMUON->SetResponseModel(chamber-1, responseTrigger0);
466 pMUON->Chamber(chamber-1).SetChargeCorrel(0); // same charge on cathodes
470 pMUON->SetNsec(chamber-1,2);
471 AliMUONSegmentationTriggerX *seg121=new AliMUONSegmentationTriggerX;
472 pMUON->SetSegmentationModel(chamber-1, 1, seg121);
473 AliMUONSegmentationTriggerY *seg122=new AliMUONSegmentationTriggerY;
474 pMUON->SetSegmentationModel(chamber-1, 2, seg122);
476 pMUON->SetResponseModel(chamber-1, responseTrigger0);
477 pMUON->Chamber(chamber-1).SetChargeCorrel(0); // same charge on cathodes
480 pMUON->SetNsec(chamber-1,2);
481 AliMUONSegmentationTriggerX *seg131=new AliMUONSegmentationTriggerX;
482 pMUON->SetSegmentationModel(chamber-1, 1, seg131);
483 AliMUONSegmentationTriggerY *seg132=new AliMUONSegmentationTriggerY;
484 pMUON->SetSegmentationModel(chamber-1, 2, seg132);
485 pMUON->SetResponseModel(chamber-1, responseTrigger0);
486 pMUON->Chamber(chamber-1).SetChargeCorrel(0); // same charge on cathodes
489 pMUON->SetNsec(chamber-1,2);
490 AliMUONSegmentationTriggerX *seg141=new AliMUONSegmentationTriggerX;
491 pMUON->SetSegmentationModel(chamber-1, 1, seg141);
492 AliMUONSegmentationTriggerY *seg142=new AliMUONSegmentationTriggerY;
493 pMUON->SetSegmentationModel(chamber-1, 2, seg142);
495 pMUON->SetResponseModel(chamber-1, responseTrigger0);
496 pMUON->Chamber(chamber-1).SetChargeCorrel(0); // same charge on cathodes
498 if(pMUON->GetDebug()) {
499 printf("\nAliMUONFactory: --------AliMUONFactory------------------------------");
500 printf("\nAliMUONFactory: Non default version of MUON selected ");
501 printf("\nAliMUONFactory: You have to construct yourself the MUON elements !!");
502 printf("\nAliMUONFactory: ----------------------------------------------------");