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 **************************************************************************/
18 Revision 1.48 2001/04/04 11:47:56 morsch
19 - muon and tau neutrinos added to g3 particle list (needed for D,B decays).
20 - some (up to now harmless) bugs in Gspart calls corrected.
22 Revision 1.47 2001/03/20 06:36:29 alibrary
23 100 parameters now allowed for geant shapes
25 Revision 1.46 2000/12/21 17:35:05 morsch
26 Last updates on the right version (1.44).
27 (1.45) does not compile.
29 Revision 1.45 2000/12/21 16:49:56 morsch
30 Adding particles to the PDG database delegated to AliPDG.
32 Revision 1.44 2000/12/20 09:46:51 alibrary
33 dlsym not supported on HP, reverting to gcomad
35 Revision 1.43 2000/12/20 08:39:39 fca
36 Support for Cerenkov and process list in Virtual MC
38 Revision 1.42 2000/12/19 08:37:48 alibrary
39 Using dlsym to retrieve address of commons
41 Revision 1.41 2000/12/18 11:33:50 alibrary
42 New call frequence histograms per module and volume
44 Revision 1.40 2000/12/06 10:06:58 morsch
45 Add all D and B baryons produced by HIJING to PDG DataBase.
47 Revision 1.39 2000/11/30 07:12:54 alibrary
48 Introducing new Rndm and QA classes
50 Revision 1.38 2000/10/30 15:19:06 morsch
51 Xi(b) (pdg code 5232) added to Pdg data base.
53 Revision 1.37 2000/10/02 21:28:16 fca
54 Removal of useless dependecies via forward declarations
56 Revision 1.36 2000/09/14 07:08:41 fca
57 Introducing glvolu in the interface
59 Revision 1.35 2000/09/12 14:27:10 morsch
60 No instance of AliDecayer created to initialize fDecayer.
62 Revision 1.34 2000/09/07 12:12:01 morsch
63 Comment inside comment removed.
65 Revision 1.33 2000/09/06 16:03:42 morsch
66 Set ExternalDecayer, Decayer and SetForceDecay methods added.
67 Gspart calls for charmed and bottom hadrons added.
68 Decay mode definitions for charmed and beauty hadrons have been taken out.
69 They will be handled by an external decayer.
71 Revision 1.32 2000/08/24 16:28:53 hristov
72 TGeant3::IsNewTrack corrected by F.Carminati
74 Revision 1.31 2000/07/13 16:19:10 fca
75 Mainly coding conventions + some small bug fixes
77 Revision 1.30 2000/07/12 08:56:30 fca
78 Coding convention correction and warning removal
80 Revision 1.29 2000/07/11 18:24:59 fca
81 Coding convention corrections + few minor bug fixes
83 Revision 1.28 2000/06/29 10:51:55 morsch
84 Add some charmed and bottom baryons to the particle list (TDatabasePDG). This
85 is needed by Hijing. Should be part of a future review of TDatabasePDG.
87 Revision 1.27 2000/06/21 17:40:15 fca
88 Adding possibility to set ISTRA, PAI model
90 Revision 1.26 2000/05/16 13:10:41 fca
91 New method IsNewTrack and fix for a problem in Father-Daughter relations
93 Revision 1.25 2000/04/07 11:12:35 fca
94 G4 compatibility changes
96 Revision 1.24 2000/02/28 21:03:57 fca
97 Some additions to improve the compatibility with G4
99 Revision 1.23 2000/02/23 16:25:25 fca
100 AliVMC and AliGeant3 classes introduced
101 ReadEuclid moved from AliRun to AliModule
103 Revision 1.22 2000/01/18 15:40:13 morsch
104 Interface to GEANT3 routines GFTMAT, GBRELM and GPRELM added
105 Define geant particle type 51: Feedback Photon with Cherenkov photon properties.
107 Revision 1.21 2000/01/17 19:41:17 fca
110 Revision 1.20 2000/01/12 11:29:27 fca
113 Revision 1.19 1999/12/17 09:03:12 fca
114 Introduce a names array
116 Revision 1.18 1999/11/26 16:55:39 fca
117 Reimplement CurrentVolName() to avoid memory leaks
119 Revision 1.17 1999/11/03 16:58:28 fca
120 Correct source of address violation in creating character string
122 Revision 1.16 1999/11/03 13:17:08 fca
123 Have ProdProcess return const char*
125 Revision 1.15 1999/10/26 06:04:50 fca
126 Introduce TLorentzVector in AliMC::GetSecondary. Thanks to I.Hrivnacova
128 Revision 1.14 1999/09/29 09:24:30 fca
129 Introduction of the Copyright and cvs Log
133 ///////////////////////////////////////////////////////////////////////////////
135 // Interface Class to the Geant3.21 MonteCarlo //
139 <img src="picts/TGeant3Class.gif">
144 ///////////////////////////////////////////////////////////////////////////////
149 #include "TDatabasePDG.h"
150 #include "TLorentzVector.h"
156 #include "AliCallf77.h"
157 #include "AliDecayer.h"
161 # define gzebra gzebra_
162 # define grfile grfile_
163 # define gpcxyz gpcxyz_
164 # define ggclos ggclos_
165 # define glast glast_
166 # define ginit ginit_
167 # define gcinit gcinit_
169 # define gtrig gtrig_
170 # define gtrigc gtrigc_
171 # define gtrigi gtrigi_
172 # define gwork gwork_
173 # define gzinit gzinit_
174 # define gfmate gfmate_
175 # define gfpart gfpart_
176 # define gftmed gftmed_
177 # define gftmat gftmat_
178 # define gmate gmate_
179 # define gpart gpart_
181 # define gsmate gsmate_
182 # define gsmixt gsmixt_
183 # define gspart gspart_
184 # define gstmed gstmed_
185 # define gsckov gsckov_
186 # define gstpar gstpar_
187 # define gfkine gfkine_
188 # define gfvert gfvert_
189 # define gskine gskine_
190 # define gsvert gsvert_
191 # define gphysi gphysi_
192 # define gdebug gdebug_
193 # define gekbin gekbin_
194 # define gfinds gfinds_
195 # define gsking gsking_
196 # define gskpho gskpho_
197 # define gsstak gsstak_
198 # define gsxyz gsxyz_
199 # define gtrack gtrack_
200 # define gtreve gtreve_
201 # define gtreveroot gtreveroot_
202 # define grndm grndm_
203 # define grndmq grndmq_
204 # define gdtom gdtom_
205 # define glmoth glmoth_
206 # define gmedia gmedia_
207 # define gmtod gmtod_
208 # define gsdvn gsdvn_
209 # define gsdvn2 gsdvn2_
210 # define gsdvs gsdvs_
211 # define gsdvs2 gsdvs2_
212 # define gsdvt gsdvt_
213 # define gsdvt2 gsdvt2_
214 # define gsord gsord_
215 # define gspos gspos_
216 # define gsposp gsposp_
217 # define gsrotm gsrotm_
218 # define gprotm gprotm_
219 # define gsvolu gsvolu_
220 # define gprint gprint_
221 # define gdinit gdinit_
222 # define gdopt gdopt_
223 # define gdraw gdraw_
224 # define gdrayt gdrayt_
225 # define gdrawc gdrawc_
226 # define gdrawx gdrawx_
227 # define gdhead gdhead_
228 # define gdwmn1 gdwmn1_
229 # define gdwmn2 gdwmn2_
230 # define gdwmn3 gdwmn3_
231 # define gdxyz gdxyz_
232 # define gdcxyz gdcxyz_
233 # define gdman gdman_
234 # define gdspec gdspec_
235 # define gdtree gdtree_
236 # define gdelet gdelet_
237 # define gdclos gdclos_
238 # define gdshow gdshow_
239 # define gdopen gdopen_
240 # define dzshow dzshow_
241 # define gsatt gsatt_
242 # define gfpara gfpara_
243 # define gckpar gckpar_
244 # define gckmat gckmat_
245 # define glvolu glvolu_
246 # define geditv geditv_
247 # define mzdrop mzdrop_
249 # define ertrak ertrak_
250 # define ertrgo ertrgo_
252 # define setbomb setbomb_
253 # define setclip setclip_
254 # define gcomad gcomad_
256 # define gbrelm gbrelm_
257 # define gprelm gprelm_
259 # define gzebra GZEBRA
260 # define grfile GRFILE
261 # define gpcxyz GPCXYZ
262 # define ggclos GGCLOS
265 # define gcinit GCINIT
268 # define gtrigc GTRIGC
269 # define gtrigi GTRIGI
271 # define gzinit GZINIT
272 # define gfmate GFMATE
273 # define gfpart GFPART
274 # define gftmed GFTMED
275 # define gftmat GFTMAT
279 # define gsmate GSMATE
280 # define gsmixt GSMIXT
281 # define gspart GSPART
282 # define gstmed GSTMED
283 # define gsckov GSCKOV
284 # define gstpar GSTPAR
285 # define gfkine GFKINE
286 # define gfvert GFVERT
287 # define gskine GSKINE
288 # define gsvert GSVERT
289 # define gphysi GPHYSI
290 # define gdebug GDEBUG
291 # define gekbin GEKBIN
292 # define gfinds GFINDS
293 # define gsking GSKING
294 # define gskpho GSKPHO
295 # define gsstak GSSTAK
297 # define gtrack GTRACK
298 # define gtreve GTREVE
299 # define gtreveroot GTREVEROOT
301 # define grndmq GRNDMQ
303 # define glmoth GLMOTH
304 # define gmedia GMEDIA
307 # define gsdvn2 GSDVN2
309 # define gsdvs2 GSDVS2
311 # define gsdvt2 GSDVT2
314 # define gsposp GSPOSP
315 # define gsrotm GSROTM
316 # define gprotm GPROTM
317 # define gsvolu GSVOLU
318 # define gprint GPRINT
319 # define gdinit GDINIT
322 # define gdrayt GDRAYT
323 # define gdrawc GDRAWC
324 # define gdrawx GDRAWX
325 # define gdhead GDHEAD
326 # define gdwmn1 GDWMN1
327 # define gdwmn2 GDWMN2
328 # define gdwmn3 GDWMN3
330 # define gdcxyz GDCXYZ
332 # define gdfspc GDFSPC
333 # define gdspec GDSPEC
334 # define gdtree GDTREE
335 # define gdelet GDELET
336 # define gdclos GDCLOS
337 # define gdshow GDSHOW
338 # define gdopen GDOPEN
339 # define dzshow DZSHOW
341 # define gfpara GFPARA
342 # define gckpar GCKPAR
343 # define gckmat GCKMAT
344 # define glvolu GLVOLU
345 # define geditv GEDITV
346 # define mzdrop MZDROP
348 # define ertrak ERTRAK
349 # define ertrgo ERTRGO
351 # define setbomb SETBOMB
352 # define setclip SETCLIP
353 # define gcomad GCOMAD
355 # define gbrelm GBRELM
356 # define gprelm GPRELM
360 //____________________________________________________________________________
364 // Prototypes for GEANT functions
366 void type_of_call gzebra(const int&);
368 void type_of_call gpcxyz();
370 void type_of_call ggclos();
372 void type_of_call glast();
374 void type_of_call ginit();
376 void type_of_call gcinit();
378 void type_of_call grun();
380 void type_of_call gtrig();
382 void type_of_call gtrigc();
384 void type_of_call gtrigi();
386 void type_of_call gwork(const int&);
388 void type_of_call gzinit();
390 void type_of_call gmate();
392 void type_of_call gpart();
394 void type_of_call gsdk(Int_t &, Float_t *, Int_t *);
396 void type_of_call gfkine(Int_t &, Float_t *, Float_t *, Int_t &,
397 Int_t &, Float_t *, Int_t &);
399 void type_of_call gfvert(Int_t &, Float_t *, Int_t &, Int_t &,
400 Float_t &, Float_t *, Int_t &);
402 void type_of_call gskine(Float_t *,Int_t &, Int_t &, Float_t *,
405 void type_of_call gsvert(Float_t *,Int_t &, Int_t &, Float_t *,
408 void type_of_call gphysi();
410 void type_of_call gdebug();
412 void type_of_call gekbin();
414 void type_of_call gfinds();
416 void type_of_call gsking(Int_t &);
418 void type_of_call gskpho(Int_t &);
420 void type_of_call gsstak(Int_t &);
422 void type_of_call gsxyz();
424 void type_of_call gtrack();
426 void type_of_call gtreve();
428 void type_of_call gtreveroot();
430 void type_of_call grndm(Float_t *r, const Int_t &n)
433 void type_of_call grndmq(Int_t &, Int_t &, const Int_t &,
435 {/*printf("Dummy grndmq called\n");*/}
437 void type_of_call gdtom(Float_t *, Float_t *, Int_t &);
439 void type_of_call glmoth(DEFCHARD, Int_t &, Int_t &, Int_t *,
440 Int_t *, Int_t * DEFCHARL);
442 void type_of_call gmedia(Float_t *, Int_t &);
444 void type_of_call gmtod(Float_t *, Float_t *, Int_t &);
446 void type_of_call gsrotm(const Int_t &, const Float_t &, const Float_t &,
447 const Float_t &, const Float_t &, const Float_t &,
450 void type_of_call gprotm(const Int_t &);
452 void type_of_call grfile(const Int_t&, DEFCHARD,
453 DEFCHARD DEFCHARL DEFCHARL);
455 void type_of_call gfmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
456 Float_t &, Float_t &, Float_t &, Float_t *,
459 void type_of_call gfpart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
460 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
462 void type_of_call gftmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
463 Float_t &, Float_t &, Float_t &, Float_t &,
464 Float_t &, Float_t &, Float_t *, Int_t * DEFCHARL);
466 void type_of_call gftmat(const Int_t&, const Int_t&, DEFCHARD, const Int_t&,
468 ,Float_t *, Int_t & DEFCHARL);
470 void type_of_call gsmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
471 Float_t &, Float_t &, Float_t &, Float_t *,
474 void type_of_call gsmixt(const Int_t&, DEFCHARD, Float_t *, Float_t *,
475 Float_t &, Int_t &, Float_t * DEFCHARL);
477 void type_of_call gspart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
478 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
481 void type_of_call gstmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
482 Float_t &, Float_t &, Float_t &, Float_t &,
483 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
485 void type_of_call gsckov(Int_t &itmed, Int_t &npckov, Float_t *ppckov,
486 Float_t *absco, Float_t *effic, Float_t *rindex);
487 void type_of_call gstpar(const Int_t&, DEFCHARD, Float_t & DEFCHARL);
489 void type_of_call gsdvn(DEFCHARD,DEFCHARD, Int_t &, Int_t &
492 void type_of_call gsdvn2(DEFCHARD,DEFCHARD, Int_t &, Int_t &, Float_t &,
493 Int_t & DEFCHARL DEFCHARL);
495 void type_of_call gsdvs(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &
498 void type_of_call gsdvs2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t &,
499 Int_t & DEFCHARL DEFCHARL);
501 void type_of_call gsdvt(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &,
502 Int_t & DEFCHARL DEFCHARL);
504 void type_of_call gsdvt2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t&,
505 Int_t &, Int_t & DEFCHARL DEFCHARL);
507 void type_of_call gsord(DEFCHARD, Int_t & DEFCHARL);
509 void type_of_call gspos(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
510 Float_t &, Int_t &, DEFCHARD DEFCHARL DEFCHARL
513 void type_of_call gsposp(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
514 Float_t &, Int_t &, DEFCHARD,
515 Float_t *, Int_t & DEFCHARL DEFCHARL DEFCHARL);
517 void type_of_call gsvolu(DEFCHARD, DEFCHARD, Int_t &, Float_t *, Int_t &,
518 Int_t & DEFCHARL DEFCHARL);
520 void type_of_call gsatt(DEFCHARD, DEFCHARD, Int_t & DEFCHARL DEFCHARL);
522 void type_of_call gfpara(DEFCHARD , Int_t&, Int_t&, Int_t&, Int_t&, Float_t*,
525 void type_of_call gckpar(Int_t&, Int_t&, Float_t*);
527 void type_of_call gckmat(Int_t&, DEFCHARD DEFCHARL);
529 void type_of_call glvolu(Int_t&, Int_t*, Int_t*, Int_t&);
531 void type_of_call gprint(DEFCHARD,const int& DEFCHARL);
533 void type_of_call gdinit();
535 void type_of_call gdopt(DEFCHARD,DEFCHARD DEFCHARL DEFCHARL);
537 void type_of_call gdraw(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
538 Float_t &, Float_t &, Float_t & DEFCHARL);
539 void type_of_call gdrayt(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
540 Float_t &, Float_t &, Float_t & DEFCHARL);
541 void type_of_call gdrawc(DEFCHARD,Int_t &, Float_t &, Float_t &, Float_t &,
542 Float_t &, Float_t & DEFCHARL);
543 void type_of_call gdrawx(DEFCHARD,Float_t &, Float_t &, Float_t &, Float_t &,
544 Float_t &, Float_t &, Float_t &, Float_t &,
546 void type_of_call gdhead(Int_t &,DEFCHARD, Float_t & DEFCHARL);
547 void type_of_call gdxyz(Int_t &);
548 void type_of_call gdcxyz();
549 void type_of_call gdman(Float_t &, Float_t &);
550 void type_of_call gdwmn1(Float_t &, Float_t &);
551 void type_of_call gdwmn2(Float_t &, Float_t &);
552 void type_of_call gdwmn3(Float_t &, Float_t &);
553 void type_of_call gdspec(DEFCHARD DEFCHARL);
554 void type_of_call gdfspc(DEFCHARD, Int_t &, Int_t & DEFCHARL) {;}
555 void type_of_call gdtree(DEFCHARD, Int_t &, Int_t & DEFCHARL);
557 void type_of_call gdopen(Int_t &);
558 void type_of_call gdclos();
559 void type_of_call gdelet(Int_t &);
560 void type_of_call gdshow(Int_t &);
561 void type_of_call geditv(Int_t &) {;}
564 void type_of_call dzshow(DEFCHARD,const int&,const int&,DEFCHARD,const int&,
565 const int&, const int&, const int& DEFCHARL
568 void type_of_call mzdrop(Int_t&, Int_t&, DEFCHARD DEFCHARL);
570 void type_of_call setbomb(Float_t &);
571 void type_of_call setclip(DEFCHARD, Float_t &,Float_t &,Float_t &,Float_t &,
572 Float_t &, Float_t & DEFCHARL);
573 void type_of_call gcomad(DEFCHARD, Int_t*& DEFCHARL);
575 void type_of_call ertrak(const Float_t *const x1, const Float_t *const p1,
576 const Float_t *x2, const Float_t *p2,
577 const Int_t &ipa, DEFCHARD DEFCHARL);
579 void type_of_call ertrgo();
581 float type_of_call gbrelm(const Float_t &z, const Float_t& t, const Float_t& cut);
582 float type_of_call gprelm(const Float_t &z, const Float_t& t, const Float_t& cut);
586 // Geant3 global pointer
588 static const Int_t kDefSize = 600;
592 //____________________________________________________________________________
596 // Default constructor
600 //____________________________________________________________________________
601 TGeant3::TGeant3(const char *title, Int_t nwgeant)
602 :AliMC("TGeant3",title)
605 // Standard constructor for TGeant3 with ZEBRA initialisation
616 // Load Address of Geant3 commons
619 // Zero number of particles
624 //____________________________________________________________________________
625 Int_t TGeant3::CurrentMaterial(Float_t &a, Float_t &z, Float_t &dens,
626 Float_t &radl, Float_t &absl) const
629 // Return the parameters of the current material during transport
633 dens = fGcmate->dens;
634 radl = fGcmate->radl;
635 absl = fGcmate->absl;
636 return 1; //this could be the number of elements in mixture
639 //____________________________________________________________________________
640 void TGeant3::DefaultRange()
643 // Set range of current drawing pad to 20x20 cm
649 gHigz->Range(0,0,20,20);
652 //____________________________________________________________________________
653 void TGeant3::InitHIGZ()
664 //____________________________________________________________________________
665 void TGeant3::LoadAddress()
668 // Assigns the address of the GEANT common blocks to the structures
669 // that allow their access from C++
672 gcomad(PASSCHARD("QUEST"), (int*&) fQuest PASSCHARL("QUEST"));
673 gcomad(PASSCHARD("GCBANK"),(int*&) fGcbank PASSCHARL("GCBANK"));
674 gcomad(PASSCHARD("GCLINK"),(int*&) fGclink PASSCHARL("GCLINK"));
675 gcomad(PASSCHARD("GCCUTS"),(int*&) fGccuts PASSCHARL("GCCUTS"));
676 gcomad(PASSCHARD("GCMULO"),(int*&) fGcmulo PASSCHARL("GCMULO"));
677 gcomad(PASSCHARD("GCFLAG"),(int*&) fGcflag PASSCHARL("GCFLAG"));
678 gcomad(PASSCHARD("GCKINE"),(int*&) fGckine PASSCHARL("GCKINE"));
679 gcomad(PASSCHARD("GCKING"),(int*&) fGcking PASSCHARL("GCKING"));
680 gcomad(PASSCHARD("GCKIN2"),(int*&) fGckin2 PASSCHARL("GCKIN2"));
681 gcomad(PASSCHARD("GCKIN3"),(int*&) fGckin3 PASSCHARL("GCKIN3"));
682 gcomad(PASSCHARD("GCMATE"),(int*&) fGcmate PASSCHARL("GCMATE"));
683 gcomad(PASSCHARD("GCTMED"),(int*&) fGctmed PASSCHARL("GCTMED"));
684 gcomad(PASSCHARD("GCTRAK"),(int*&) fGctrak PASSCHARL("GCTRAK"));
685 gcomad(PASSCHARD("GCTPOL"),(int*&) fGctpol PASSCHARL("GCTPOL"));
686 gcomad(PASSCHARD("GCVOLU"),(int*&) fGcvolu PASSCHARL("GCVOLU"));
687 gcomad(PASSCHARD("GCNUM"), (int*&) fGcnum PASSCHARL("GCNUM"));
688 gcomad(PASSCHARD("GCSETS"),(int*&) fGcsets PASSCHARL("GCSETS"));
689 gcomad(PASSCHARD("GCPHYS"),(int*&) fGcphys PASSCHARL("GCPHYS"));
690 gcomad(PASSCHARD("GCPHLT"),(int*&) fGcphlt PASSCHARL("GCPHLT"));
691 gcomad(PASSCHARD("GCOPTI"),(int*&) fGcopti PASSCHARL("GCOPTI"));
692 gcomad(PASSCHARD("GCTLIT"),(int*&) fGctlit PASSCHARL("GCTLIT"));
693 gcomad(PASSCHARD("GCVDMA"),(int*&) fGcvdma PASSCHARL("GCVDMA"));
696 gcomad(PASSCHARD("ERTRIO"),(int*&) fErtrio PASSCHARL("ERTRIO"));
697 gcomad(PASSCHARD("EROPTS"),(int*&) fEropts PASSCHARL("EROPTS"));
698 gcomad(PASSCHARD("EROPTC"),(int*&) fEroptc PASSCHARL("EROPTC"));
699 gcomad(PASSCHARD("ERWORK"),(int*&) fErwork PASSCHARL("ERWORK"));
701 // Variables for ZEBRA store
702 gcomad(PASSCHARD("IQ"), addr PASSCHARL("IQ"));
704 gcomad(PASSCHARD("LQ"), addr PASSCHARL("LQ"));
709 //_____________________________________________________________________________
710 void TGeant3::GeomIter()
713 // Geometry iterator for moving upward in the geometry tree
714 // Initialise the iterator
716 fNextVol=fGcvolu->nlevel;
719 //____________________________________________________________________________
720 void TGeant3::FinishGeometry()
722 //Close the geometry structure
726 //____________________________________________________________________________
727 Int_t TGeant3::NextVolUp(Text_t *name, Int_t ©)
730 // Geometry iterator for moving upward in the geometry tree
731 // Return next volume up
736 gname=fGcvolu->names[fNextVol];
737 copy=fGcvolu->number[fNextVol];
738 i=fGcvolu->lvolum[fNextVol];
739 name = fVolNames[i-1];
740 if(gname == fZiq[fGclink->jvolum+i]) return i;
741 else printf("GeomTree: Volume %s not found in bank\n",name);
746 //_____________________________________________________________________________
747 void TGeant3::BuildPhysics()
752 //_____________________________________________________________________________
753 Int_t TGeant3::CurrentVolID(Int_t ©) const
756 // Returns the current volume ID and copy number
759 if( (i=fGcvolu->nlevel-1) < 0 ) {
760 Warning("CurrentVolID","Stack depth only %d\n",fGcvolu->nlevel);
762 gname=fGcvolu->names[i];
763 copy=fGcvolu->number[i];
764 i=fGcvolu->lvolum[i];
765 if(gname == fZiq[fGclink->jvolum+i]) return i;
766 else Warning("CurrentVolID","Volume %4s not found\n",(char*)&gname);
771 //_____________________________________________________________________________
772 Int_t TGeant3::CurrentVolOffID(Int_t off, Int_t ©) const
775 // Return the current volume "off" upward in the geometrical tree
776 // ID and copy number
779 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
780 Warning("CurrentVolOffID","Offset requested %d but stack depth %d\n",
781 off,fGcvolu->nlevel);
783 gname=fGcvolu->names[i];
784 copy=fGcvolu->number[i];
785 i=fGcvolu->lvolum[i];
786 if(gname == fZiq[fGclink->jvolum+i]) return i;
787 else Warning("CurrentVolOffID","Volume %4s not found\n",(char*)&gname);
792 //_____________________________________________________________________________
793 const char* TGeant3::CurrentVolName() const
796 // Returns the current volume name
799 if( (i=fGcvolu->nlevel-1) < 0 ) {
800 Warning("CurrentVolName","Stack depth %d\n",fGcvolu->nlevel);
802 gname=fGcvolu->names[i];
803 i=fGcvolu->lvolum[i];
804 if(gname == fZiq[fGclink->jvolum+i]) return fVolNames[i-1];
805 else Warning("CurrentVolName","Volume %4s not found\n",(char*) &gname);
810 //_____________________________________________________________________________
811 const char* TGeant3::CurrentVolOffName(Int_t off) const
814 // Return the current volume "off" upward in the geometrical tree
815 // ID, name and copy number
816 // if name=0 no name is returned
819 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
820 Warning("CurrentVolOffName",
821 "Offset requested %d but stack depth %d\n",off,fGcvolu->nlevel);
823 gname=fGcvolu->names[i];
824 i=fGcvolu->lvolum[i];
825 if(gname == fZiq[fGclink->jvolum+i]) return fVolNames[i-1];
826 else Warning("CurrentVolOffName","Volume %4s not found\n",(char*)&gname);
831 //_____________________________________________________________________________
832 Int_t TGeant3::IdFromPDG(Int_t pdg) const
835 // Return Geant3 code from PDG and pseudo ENDF code
837 for(Int_t i=0;i<fNPDGCodes;++i)
838 if(pdg==fPDGCode[i]) return i;
842 //_____________________________________________________________________________
843 Int_t TGeant3::PDGFromId(Int_t id) const
846 // Return PDG code and pseudo ENDF code from Geant3 code
848 if(id>0 && id<fNPDGCodes) return fPDGCode[id];
852 //_____________________________________________________________________________
853 void TGeant3::DefineParticles()
856 // Define standard Geant 3 particles
859 // Load standard numbers for GEANT particles and PDG conversion
860 fPDGCode[fNPDGCodes++]=-99; // 0 = unused location
861 fPDGCode[fNPDGCodes++]=22; // 1 = photon
862 fPDGCode[fNPDGCodes++]=-11; // 2 = positron
863 fPDGCode[fNPDGCodes++]=11; // 3 = electron
864 fPDGCode[fNPDGCodes++]=12; // 4 = neutrino e
865 fPDGCode[fNPDGCodes++]=-13; // 5 = muon +
866 fPDGCode[fNPDGCodes++]=13; // 6 = muon -
867 fPDGCode[fNPDGCodes++]=111; // 7 = pi0
868 fPDGCode[fNPDGCodes++]=211; // 8 = pi+
869 fPDGCode[fNPDGCodes++]=-211; // 9 = pi-
870 fPDGCode[fNPDGCodes++]=130; // 10 = Kaon Long
871 fPDGCode[fNPDGCodes++]=321; // 11 = Kaon +
872 fPDGCode[fNPDGCodes++]=-321; // 12 = Kaon -
873 fPDGCode[fNPDGCodes++]=2112; // 13 = Neutron
874 fPDGCode[fNPDGCodes++]=2212; // 14 = Proton
875 fPDGCode[fNPDGCodes++]=-2212; // 15 = Anti Proton
876 fPDGCode[fNPDGCodes++]=310; // 16 = Kaon Short
877 fPDGCode[fNPDGCodes++]=221; // 17 = Eta
878 fPDGCode[fNPDGCodes++]=3122; // 18 = Lambda
879 fPDGCode[fNPDGCodes++]=3222; // 19 = Sigma +
880 fPDGCode[fNPDGCodes++]=3212; // 20 = Sigma 0
881 fPDGCode[fNPDGCodes++]=3112; // 21 = Sigma -
882 fPDGCode[fNPDGCodes++]=3322; // 22 = Xi0
883 fPDGCode[fNPDGCodes++]=3312; // 23 = Xi-
884 fPDGCode[fNPDGCodes++]=3334; // 24 = Omega-
885 fPDGCode[fNPDGCodes++]=-2112; // 25 = Anti Proton
886 fPDGCode[fNPDGCodes++]=-3122; // 26 = Anti Proton
887 fPDGCode[fNPDGCodes++]=-3222; // 27 = Anti Sigma -
888 fPDGCode[fNPDGCodes++]=-3212; // 28 = Anti Sigma 0
889 fPDGCode[fNPDGCodes++]=-3112; // 29 = Anti Sigma 0
890 fPDGCode[fNPDGCodes++]=-3322; // 30 = Anti Xi 0
891 fPDGCode[fNPDGCodes++]=-3312; // 31 = Anti Xi +
892 fPDGCode[fNPDGCodes++]=-3334; // 32 = Anti Omega +
899 /* --- Define additional particles */
900 Gspart(33, "OMEGA(782)", 3, 0.782, 0., 7.836e-23);
901 fPDGCode[fNPDGCodes++]=223; // 33 = Omega(782)
903 Gspart(34, "PHI(1020)", 3, 1.019, 0., 1.486e-22);
904 fPDGCode[fNPDGCodes++]=333; // 34 = PHI (1020)
906 Gspart(35, "D +", 4, 1.87, 1., 1.066e-12);
907 fPDGCode[fNPDGCodes++]=411; // 35 = D+
909 Gspart(36, "D -", 4, 1.87, -1., 1.066e-12);
910 fPDGCode[fNPDGCodes++]=-411; // 36 = D-
912 Gspart(37, "D 0", 3, 1.865, 0., 4.2e-13);
913 fPDGCode[fNPDGCodes++]=421; // 37 = D0
915 Gspart(38, "ANTI D 0", 3, 1.865, 0., 4.2e-13);
916 fPDGCode[fNPDGCodes++]=-421; // 38 = D0 bar
919 fPDGCode[fNPDGCodes++]=-99; // 39 = unassigned
921 fPDGCode[fNPDGCodes++]=-99; // 40 = unassigned
923 fPDGCode[fNPDGCodes++]=-99; // 41 = unassigned
925 Gspart(42, "RHO +", 4, 0.768, 1., 4.353e-24);
926 fPDGCode[fNPDGCodes++]=213; // 42 = RHO+
928 Gspart(43, "RHO -", 4, 0.768, -1., 4.353e-24);
929 fPDGCode[fNPDGCodes++]=-213; // 43 = RHO-
931 Gspart(44, "RHO 0", 3, 0.768, 0., 4.353e-24);
932 fPDGCode[fNPDGCodes++]=113; // 44 = RHO0
935 // Use ENDF-6 mapping for ions = 10000*z+10*a+iso
937 // and numbers above 5 000 000 for special applications
940 const Int_t kion=10000000;
942 const Int_t kspe=50000000;
947 fPDGCode[fNPDGCodes++]=kion+10020; // 45 = Deuteron
949 fPDGCode[fNPDGCodes++]=kion+10030; // 46 = Triton
951 fPDGCode[fNPDGCodes++]=kion+20040; // 47 = Alpha
953 fPDGCode[fNPDGCodes++]=0; // 48 = geantino mapped to rootino
955 fPDGCode[fNPDGCodes++]=kion+20030; // 49 = HE3
957 fPDGCode[fNPDGCodes++]=kspe+50; // 50 = Cherenkov
959 Gspart(51, "FeedbackPhoton", 7, 0., 0.,1.e20 );
960 fPDGCode[fNPDGCodes++]=kspe+51; // 51 = FeedbackPhoton
962 Gspart(52, "Lambda_c+", 4, 2.2849, +1., 2.06e-13);
963 fPDGCode[fNPDGCodes++]=4122; //52 = Lambda_c+
965 Gspart(53, "Lambda_c-", 4, 2.2849, -1., 2.06e-13);
966 fPDGCode[fNPDGCodes++]=-4122; //53 = Lambda_c-
968 Gspart(54, "D_s+", 4, 1.9685, +1., 4.67e-13);
969 fPDGCode[fNPDGCodes++]=431; //54 = D_s+
971 Gspart(55, "D_s-", 4, 1.9685, -1., 4.67e-13);
972 fPDGCode[fNPDGCodes++]=-431; //55 = D_s-
974 Gspart(56, "Tau+", 5, 1.77705, +1., 2.9e-13);
975 fPDGCode[fNPDGCodes++]=15; //56 = Tau+
977 Gspart(57, "Tau-", 5, 1.77705, -1., 2.9e-13);
978 fPDGCode[fNPDGCodes++]=-15; //57 = Tau-
980 Gspart(58, "B0", 3, 5.2792, +0., 1.56e-12);
981 fPDGCode[fNPDGCodes++]=511; //58 = B0
983 Gspart(59, "B0 bar", 3, 5.2792, -0., 1.56e-12);
984 fPDGCode[fNPDGCodes++]=-511; //58 = B0bar
986 Gspart(60, "B+", 4, 5.2789, +1., 1.65e-12);
987 fPDGCode[fNPDGCodes++]=521; //60 = B+
989 Gspart(61, "B-", 4, 5.2789, -1., 1.65e-12);
990 fPDGCode[fNPDGCodes++]=-521; //61 = B-
992 Gspart(62, "Bs", 3, 5.3693, +0., 1.54e-12);
993 fPDGCode[fNPDGCodes++]=531; //62 = B_s
995 Gspart(63, "Bs bar", 3, 5.3693, -0., 1.54e-12);
996 fPDGCode[fNPDGCodes++]=-531; //63 = B_s bar
998 Gspart(64, "Lambda_b", 3, 5.624, +0., 1.24e-12);
999 fPDGCode[fNPDGCodes++]=5122; //64 = Lambda_b
1001 Gspart(65, "Lambda_b bar", 3, 5.624, -0., 1.24e-12);
1002 fPDGCode[fNPDGCodes++]=-5122; //65 = Lambda_b bar
1004 Gspart(66, "J/Psi", 3, 3.09688, 0., 0.);
1005 fPDGCode[fNPDGCodes++]=443; // 66 = J/Psi
1007 Gspart(67, "Psi Prime", 3, 3.686, 0., 0.);
1008 fPDGCode[fNPDGCodes++]=20443; // 67 = Psi prime
1010 Gspart(68, "Upsilon(1S)", 3, 9.46037, 0., 0.);
1011 fPDGCode[fNPDGCodes++]=553; // 68 = Upsilon(1S)
1013 Gspart(69, "Upsilon(2S)", 3, 10.0233, 0., 0.);
1014 fPDGCode[fNPDGCodes++]=20553; // 69 = Upsilon(2S)
1016 Gspart(70, "Upsilon(3S)", 3, 10.3553, 0., 0.);
1017 fPDGCode[fNPDGCodes++]=30553; // 70 = Upsilon(3S)
1019 Gspart(71, "Anti Neutrino (e)", 3, 0., 0., 1.e20);
1020 fPDGCode[fNPDGCodes++]=-12; // 71 = anti electron neutrino
1022 Gspart(72, "Neutrino (mu)", 3, 0., 0., 1.e20);
1023 fPDGCode[fNPDGCodes++]=14; // 72 = muon neutrino
1025 Gspart(73, "Anti Neutrino (mu)", 3, 0., 0., 1.e20);
1026 fPDGCode[fNPDGCodes++]=-14; // 73 = anti muon neutrino
1028 Gspart(74, "Neutrino (tau)", 3, 0., 0., 1.e20);
1029 fPDGCode[fNPDGCodes++]=16; // 74 = tau neutrino
1031 Gspart(75, "Anti Neutrino (tau)",3, 0., 0., 1.e20);
1032 fPDGCode[fNPDGCodes++]=-16; // 75 = anti tau neutrino
1034 /* --- Define additional decay modes --- */
1035 /* --- omega(783) --- */
1036 for (kz = 0; kz < 6; ++kz) {
1047 Gsdk(ipa, bratio, mode);
1048 /* --- phi(1020) --- */
1049 for (kz = 0; kz < 6; ++kz) {
1064 Gsdk(ipa, bratio, mode);
1067 for (kz = 0; kz < 6; ++kz) {
1080 Gsdk(ipa, bratio, mode);
1084 for (kz = 0; kz < 6; ++kz) {
1097 Gsdk(ipa, bratio, mode);
1101 for (kz = 0; kz < 6; ++kz) {
1112 Gsdk(ipa, bratio, mode);
1114 /* --- Anti D0 --- */
1116 for (kz = 0; kz < 6; ++kz) {
1127 Gsdk(ipa, bratio, mode);
1130 for (kz = 0; kz < 6; ++kz) {
1137 Gsdk(ipa, bratio, mode);
1139 for (kz = 0; kz < 6; ++kz) {
1146 Gsdk(ipa, bratio, mode);
1148 for (kz = 0; kz < 6; ++kz) {
1155 Gsdk(ipa, bratio, mode);
1158 for (kz = 0; kz < 6; ++kz) {
1167 Gsdk(ipa, bratio, mode);
1170 Gsdk(ipa, bratio, mode);
1173 Gsdk(ipa, bratio, mode);
1176 AliPDG::AddParticlesToPdgDataBase();
1179 //_____________________________________________________________________________
1180 Int_t TGeant3::VolId(const Text_t *name) const
1183 // Return the unique numeric identifier for volume name
1186 strncpy((char *) &gname, name, 4);
1187 for(i=1; i<=fGcnum->nvolum; i++)
1188 if(gname == fZiq[fGclink->jvolum+i]) return i;
1189 printf("VolId: Volume %s not found\n",name);
1193 //_____________________________________________________________________________
1194 Int_t TGeant3::NofVolumes() const
1197 // Return total number of volumes in the geometry
1199 return fGcnum->nvolum;
1202 //_____________________________________________________________________________
1203 Int_t TGeant3::VolId2Mate(Int_t id) const
1206 // Return material number for a given volume id
1208 if(id<1 || id > fGcnum->nvolum || fGclink->jvolum<=0)
1211 Int_t jvo = fZlq[fGclink->jvolum-id];
1212 return Int_t(fZq[jvo+4]);
1216 //_____________________________________________________________________________
1217 const char* TGeant3::VolName(Int_t id) const
1220 // Return the volume name given the volume identifier
1222 if(id<1 || id > fGcnum->nvolum || fGclink->jvolum<=0)
1223 return fVolNames[fGcnum->nvolum];
1225 return fVolNames[id-1];
1228 //_____________________________________________________________________________
1229 void TGeant3::SetCut(const char* cutName, Float_t cutValue)
1232 // Set transport cuts for particles
1234 if(!strcmp(cutName,"CUTGAM"))
1235 fGccuts->cutgam=cutValue;
1236 else if(!strcmp(cutName,"CUTGAM"))
1237 fGccuts->cutele=cutValue;
1238 else if(!strcmp(cutName,"CUTELE"))
1239 fGccuts->cutneu=cutValue;
1240 else if(!strcmp(cutName,"CUTHAD"))
1241 fGccuts->cuthad=cutValue;
1242 else if(!strcmp(cutName,"CUTMUO"))
1243 fGccuts->cutmuo=cutValue;
1244 else if(!strcmp(cutName,"BCUTE"))
1245 fGccuts->bcute=cutValue;
1246 else if(!strcmp(cutName,"BCUTM"))
1247 fGccuts->bcutm=cutValue;
1248 else if(!strcmp(cutName,"DCUTE"))
1249 fGccuts->dcute=cutValue;
1250 else if(!strcmp(cutName,"DCUTM"))
1251 fGccuts->dcutm=cutValue;
1252 else if(!strcmp(cutName,"PPCUTM"))
1253 fGccuts->ppcutm=cutValue;
1254 else if(!strcmp(cutName,"TOFMAX"))
1255 fGccuts->tofmax=cutValue;
1256 else Warning("SetCut","Cut %s not implemented\n",cutName);
1259 //_____________________________________________________________________________
1260 void TGeant3::SetProcess(const char* flagName, Int_t flagValue)
1263 // Set thresholds for different processes
1265 if(!strcmp(flagName,"PAIR"))
1266 fGcphys->ipair=flagValue;
1267 else if(!strcmp(flagName,"COMP"))
1268 fGcphys->icomp=flagValue;
1269 else if(!strcmp(flagName,"PHOT"))
1270 fGcphys->iphot=flagValue;
1271 else if(!strcmp(flagName,"PFIS"))
1272 fGcphys->ipfis=flagValue;
1273 else if(!strcmp(flagName,"DRAY"))
1274 fGcphys->idray=flagValue;
1275 else if(!strcmp(flagName,"ANNI"))
1276 fGcphys->ianni=flagValue;
1277 else if(!strcmp(flagName,"BREM"))
1278 fGcphys->ibrem=flagValue;
1279 else if(!strcmp(flagName,"HADR"))
1280 fGcphys->ihadr=flagValue;
1281 else if(!strcmp(flagName,"MUNU"))
1282 fGcphys->imunu=flagValue;
1283 else if(!strcmp(flagName,"DCAY"))
1284 fGcphys->idcay=flagValue;
1285 else if(!strcmp(flagName,"LOSS"))
1286 fGcphys->iloss=flagValue;
1287 else if(!strcmp(flagName,"MULS"))
1288 fGcphys->imuls=flagValue;
1289 else if(!strcmp(flagName,"RAYL"))
1290 fGcphys->irayl=flagValue;
1291 else if(!strcmp(flagName,"STRA"))
1292 fGcphlt->istra=flagValue;
1293 else if(!strcmp(flagName,"SYNC"))
1294 fGcphlt->isync=flagValue;
1295 else Warning("SetFlag","Flag %s not implemented\n",flagName);
1298 //_____________________________________________________________________________
1299 Float_t TGeant3::Xsec(char* reac, Float_t /* energy */,
1300 Int_t part, Int_t /* mate */)
1303 // Calculate X-sections -- dummy for the moment
1305 if(!strcmp(reac,"PHOT"))
1308 Error("Xsec","Can calculate photoelectric only for photons\n");
1314 //_____________________________________________________________________________
1315 void TGeant3::TrackPosition(TLorentzVector &xyz) const
1318 // Return the current position in the master reference frame of the
1319 // track being transported
1321 xyz[0]=fGctrak->vect[0];
1322 xyz[1]=fGctrak->vect[1];
1323 xyz[2]=fGctrak->vect[2];
1324 xyz[3]=fGctrak->tofg;
1327 //_____________________________________________________________________________
1328 Float_t TGeant3::TrackTime() const
1331 // Return the current time of flight of the track being transported
1333 return fGctrak->tofg;
1336 //_____________________________________________________________________________
1337 void TGeant3::TrackMomentum(TLorentzVector &xyz) const
1340 // Return the direction and the momentum (GeV/c) of the track
1341 // currently being transported
1343 Double_t ptot=fGctrak->vect[6];
1344 xyz[0]=fGctrak->vect[3]*ptot;
1345 xyz[1]=fGctrak->vect[4]*ptot;
1346 xyz[2]=fGctrak->vect[5]*ptot;
1347 xyz[3]=fGctrak->getot;
1350 //_____________________________________________________________________________
1351 Float_t TGeant3::TrackCharge() const
1354 // Return charge of the track currently transported
1356 return fGckine->charge;
1359 //_____________________________________________________________________________
1360 Float_t TGeant3::TrackMass() const
1363 // Return the mass of the track currently transported
1365 return fGckine->amass;
1368 //_____________________________________________________________________________
1369 Int_t TGeant3::TrackPid() const
1372 // Return the id of the particle transported
1374 return PDGFromId(fGckine->ipart);
1377 //_____________________________________________________________________________
1378 Float_t TGeant3::TrackStep() const
1381 // Return the length in centimeters of the current step
1383 return fGctrak->step;
1386 //_____________________________________________________________________________
1387 Float_t TGeant3::TrackLength() const
1390 // Return the length of the current track from its origin
1392 return fGctrak->sleng;
1395 //_____________________________________________________________________________
1396 Bool_t TGeant3::IsNewTrack() const
1399 // True if the track is not at the boundary of the current volume
1401 return (fGctrak->sleng==0);
1404 //_____________________________________________________________________________
1405 Bool_t TGeant3::IsTrackInside() const
1408 // True if the track is not at the boundary of the current volume
1410 return (fGctrak->inwvol==0);
1413 //_____________________________________________________________________________
1414 Bool_t TGeant3::IsTrackEntering() const
1417 // True if this is the first step of the track in the current volume
1419 return (fGctrak->inwvol==1);
1422 //_____________________________________________________________________________
1423 Bool_t TGeant3::IsTrackExiting() const
1426 // True if this is the last step of the track in the current volume
1428 return (fGctrak->inwvol==2);
1431 //_____________________________________________________________________________
1432 Bool_t TGeant3::IsTrackOut() const
1435 // True if the track is out of the setup
1437 return (fGctrak->inwvol==3);
1440 //_____________________________________________________________________________
1441 Bool_t TGeant3::IsTrackStop() const
1444 // True if the track energy has fallen below the threshold
1446 return (fGctrak->istop==2);
1449 //_____________________________________________________________________________
1450 Int_t TGeant3::NSecondaries() const
1453 // Number of secondary particles generated in the current step
1455 return fGcking->ngkine;
1458 //_____________________________________________________________________________
1459 Int_t TGeant3::CurrentEvent() const
1462 // Number of the current event
1464 return fGcflag->idevt;
1467 //_____________________________________________________________________________
1468 AliMCProcess TGeant3::ProdProcess(Int_t ) const
1471 // Name of the process that has produced the secondary particles
1472 // in the current step
1474 const AliMCProcess kIpProc[13] = { kPDecay, kPPair, kPCompton,
1475 kPPhotoelectric, kPBrem, kPDeltaRay,
1476 kPAnnihilation, kPHadronic,
1477 kPMuonNuclear, kPPhotoFission,
1478 kPRayleigh, kPCerenkov, kPSynchrotron};
1481 if(fGcking->ngkine>0)
1482 for (km = 0; km < fGctrak->nmec; ++km)
1483 for (im = 0; im < 13; ++im)
1484 if (G3toVMC(fGctrak->lmec[km]) == kIpProc[im])
1490 //_____________________________________________________________________________
1491 Int_t TGeant3::StepProcesses(TArrayI &proc) const
1494 // Return processes active in the current step
1497 Int_t nproc=Gctrak()->nmec;
1502 for (i=0; i<nproc; ++i)
1503 if((proc[nvproc]=G3toVMC(Gctrak()->lmec[i]))!=kPNoProcess) nvproc++;
1510 //_____________________________________________________________________________
1511 AliMCProcess TGeant3::G3toVMC(Int_t iproc) const
1514 // Conversion between GEANT and AliMC processes
1517 const AliMCProcess kPG2MC1[30] = {kPNoProcess, kPMultipleScattering, kPEnergyLoss, kPMagneticFieldL, kPDecay,
1518 kPPair, kPCompton, kPPhotoelectric, kPBrem, kPDeltaRay,
1519 kPAnnihilation, kPHadronic, kPNoProcess, kPEvaporation, kPNuclearFission,
1520 kPNuclearAbsorption, kPPbarAnnihilation, kPNCapture, kPHElastic, kPHInhelastic,
1521 kPMuonNuclear, kPTOFlimit, kPPhotoFission, kPNoProcess, kPRayleigh,
1522 kPNoProcess, kPNoProcess, kPNoProcess, kPNull, kPStop};
1524 const AliMCProcess kPG2MC2[9] = {kPLightAbsorption, kPLightScattering, kStepMax, kPNoProcess, kPCerenkov,
1525 kPLightReflection, kPLightRefraction, kPSynchrotron, kPNoProcess};
1527 AliMCProcess proc=kPNoProcess;
1528 if(1<iproc && iproc<=30) proc= kPG2MC1[iproc-1];
1529 else if(101<=iproc && iproc<=109) proc= kPG2MC2[iproc-100-1];
1534 //_____________________________________________________________________________
1535 void TGeant3::GetSecondary(Int_t isec, Int_t& ipart,
1536 TLorentzVector &x, TLorentzVector &p)
1539 // Get the parameters of the secondary track number isec produced
1540 // in the current step
1543 if(-1<isec && isec<fGcking->ngkine) {
1544 ipart=Int_t (fGcking->gkin[isec][4] +0.5);
1546 x[i]=fGckin3->gpos[isec][i];
1547 p[i]=fGcking->gkin[isec][i];
1549 x[3]=fGcking->tofd[isec];
1550 p[3]=fGcking->gkin[isec][3];
1552 printf(" * TGeant3::GetSecondary * Secondary %d does not exist\n",isec);
1553 x[0]=x[1]=x[2]=x[3]=p[0]=p[1]=p[2]=p[3]=0;
1558 //_____________________________________________________________________________
1559 void TGeant3::InitLego()
1562 // Set switches for lego transport
1565 SetDEBU(0,0,0); //do not print a message
1568 //_____________________________________________________________________________
1569 Bool_t TGeant3::IsTrackDisappeared() const
1572 // True if the current particle has disappered
1573 // either because it decayed or because it underwent
1574 // an inelastic collision
1576 return (fGctrak->istop==1);
1579 //_____________________________________________________________________________
1580 Bool_t TGeant3::IsTrackAlive() const
1583 // True if the current particle is alive and will continue to be
1586 return (fGctrak->istop==0);
1589 //_____________________________________________________________________________
1590 void TGeant3::StopTrack()
1593 // Stop the transport of the current particle and skip to the next
1598 //_____________________________________________________________________________
1599 void TGeant3::StopEvent()
1602 // Stop simulation of the current event and skip to the next
1607 //_____________________________________________________________________________
1608 Float_t TGeant3::MaxStep() const
1611 // Return the maximum step length in the current medium
1613 return fGctmed->stemax;
1616 //_____________________________________________________________________________
1617 void TGeant3::SetMaxStep(Float_t maxstep)
1620 // Set the maximum step allowed till the particle is in the current medium
1622 fGctmed->stemax=maxstep;
1625 //_____________________________________________________________________________
1626 void TGeant3::SetMaxNStep(Int_t maxnstp)
1629 // Set the maximum number of steps till the particle is in the current medium
1631 fGctrak->maxnst=maxnstp;
1634 //_____________________________________________________________________________
1635 Int_t TGeant3::GetMaxNStep() const
1638 // Maximum number of steps allowed in current medium
1640 return fGctrak->maxnst;
1643 //_____________________________________________________________________________
1644 void TGeant3::Material(Int_t& kmat, const char* name, Float_t a, Float_t z,
1645 Float_t dens, Float_t radl, Float_t absl, Float_t* buf,
1649 // Defines a Material
1651 // kmat number assigned to the material
1652 // name material name
1653 // a atomic mass in au
1655 // dens density in g/cm3
1656 // absl absorbtion length in cm
1657 // if >=0 it is ignored and the program
1658 // calculates it, if <0. -absl is taken
1659 // radl radiation length in cm
1660 // if >=0 it is ignored and the program
1661 // calculates it, if <0. -radl is taken
1662 // buf pointer to an array of user words
1663 // nbuf number of user words
1665 Int_t jmate=fGclink->jmate;
1671 for(i=1; i<=ns; i++) {
1672 if(fZlq[jmate-i]==0) {
1678 gsmate(kmat,PASSCHARD(name), a, z, dens, radl, absl, buf,
1679 nwbuf PASSCHARL(name));
1682 //_____________________________________________________________________________
1683 void TGeant3::Mixture(Int_t& kmat, const char* name, Float_t* a, Float_t* z,
1684 Float_t dens, Int_t nlmat, Float_t* wmat)
1687 // Defines mixture OR COMPOUND IMAT as composed by
1688 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
1690 // If NLMAT > 0 then wmat contains the proportion by
1691 // weights of each basic material in the mixture.
1693 // If nlmat < 0 then WMAT contains the number of atoms
1694 // of a given kind into the molecule of the COMPOUND
1695 // In this case, WMAT in output is changed to relative
1698 Int_t jmate=fGclink->jmate;
1704 for(i=1; i<=ns; i++) {
1705 if(fZlq[jmate-i]==0) {
1711 gsmixt(kmat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
1714 //_____________________________________________________________________________
1715 void TGeant3::Medium(Int_t& kmed, const char* name, Int_t nmat, Int_t isvol,
1716 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
1717 Float_t stemax, Float_t deemax, Float_t epsil,
1718 Float_t stmin, Float_t* ubuf, Int_t nbuf)
1721 // kmed tracking medium number assigned
1722 // name tracking medium name
1723 // nmat material number
1724 // isvol sensitive volume flag
1725 // ifield magnetic field
1726 // fieldm max. field value (kilogauss)
1727 // tmaxfd max. angle due to field (deg/step)
1728 // stemax max. step allowed
1729 // deemax max. fraction of energy lost in a step
1730 // epsil tracking precision (cm)
1731 // stmin min. step due to continuos processes (cm)
1733 // ifield = 0 if no magnetic field; ifield = -1 if user decision in guswim;
1734 // ifield = 1 if tracking performed with grkuta; ifield = 2 if tracking
1735 // performed with ghelix; ifield = 3 if tracking performed with ghelx3.
1737 Int_t jtmed=fGclink->jtmed;
1743 for(i=1; i<=ns; i++) {
1744 if(fZlq[jtmed-i]==0) {
1750 gstmed(kmed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1751 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1754 //_____________________________________________________________________________
1755 void TGeant3::Matrix(Int_t& krot, Float_t thex, Float_t phix, Float_t they,
1756 Float_t phiy, Float_t thez, Float_t phiz)
1759 // krot rotation matrix number assigned
1760 // theta1 polar angle for axis i
1761 // phi1 azimuthal angle for axis i
1762 // theta2 polar angle for axis ii
1763 // phi2 azimuthal angle for axis ii
1764 // theta3 polar angle for axis iii
1765 // phi3 azimuthal angle for axis iii
1767 // it defines the rotation matrix number irot.
1769 Int_t jrotm=fGclink->jrotm;
1775 for(i=1; i<=ns; i++) {
1776 if(fZlq[jrotm-i]==0) {
1782 gsrotm(krot, thex, phix, they, phiy, thez, phiz);
1785 //_____________________________________________________________________________
1786 Int_t TGeant3::GetMedium() const
1789 // Return the number of the current medium
1791 return fGctmed->numed;
1794 //_____________________________________________________________________________
1795 Float_t TGeant3::Edep() const
1798 // Return the energy lost in the current step
1800 return fGctrak->destep;
1803 //_____________________________________________________________________________
1804 Float_t TGeant3::Etot() const
1807 // Return the total energy of the current track
1809 return fGctrak->getot;
1812 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1814 // Functions from GBASE
1816 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1818 //____________________________________________________________________________
1819 void TGeant3::Gfile(const char *filename, const char *option)
1822 // Routine to open a GEANT/RZ data base.
1824 // LUN logical unit number associated to the file
1826 // CHFILE RZ file name
1828 // CHOPT is a character string which may be
1829 // N To create a new file
1830 // U to open an existing file for update
1831 // " " to open an existing file for read only
1832 // Q The initial allocation (default 1000 records)
1833 // is given in IQUEST(10)
1834 // X Open the file in exchange format
1835 // I Read all data structures from file to memory
1836 // O Write all data structures from memory to file
1839 // If options "I" or "O" all data structures are read or
1840 // written from/to file and the file is closed.
1841 // See routine GRMDIR to create subdirectories
1842 // See routines GROUT,GRIN to write,read objects
1844 grfile(21, PASSCHARD(filename), PASSCHARD(option) PASSCHARL(filename)
1848 //____________________________________________________________________________
1849 void TGeant3::Gpcxyz()
1852 // Print track and volume parameters at current point
1857 //_____________________________________________________________________________
1858 void TGeant3::Ggclos()
1861 // Closes off the geometry setting.
1862 // Initializes the search list for the contents of each
1863 // volume following the order they have been positioned, and
1864 // inserting the content '0' when a call to GSNEXT (-1) has
1865 // been required by the user.
1866 // Performs the development of the JVOLUM structure for all
1867 // volumes with variable parameters, by calling GGDVLP.
1868 // Interprets the user calls to GSORD, through GGORD.
1869 // Computes and stores in a bank (next to JVOLUM mother bank)
1870 // the number of levels in the geometrical tree and the
1871 // maximum number of contents per level, by calling GGNLEV.
1872 // Sets status bit for CONCAVE volumes, through GGCAVE.
1873 // Completes the JSET structure with the list of volume names
1874 // which identify uniquely a given physical detector, the
1875 // list of bit numbers to pack the corresponding volume copy
1876 // numbers, and the generic path(s) in the JVOLUM tree,
1877 // through the routine GHCLOS.
1880 // Create internal list of volumes
1881 fVolNames = new char[fGcnum->nvolum+1][5];
1883 for(i=0; i<fGcnum->nvolum; ++i) {
1884 strncpy(fVolNames[i], (char *) &fZiq[fGclink->jvolum+i+1], 4);
1885 fVolNames[i][4]='\0';
1887 strcpy(fVolNames[fGcnum->nvolum],"NULL");
1890 //_____________________________________________________________________________
1891 void TGeant3::Glast()
1894 // Finish a Geant run
1899 //_____________________________________________________________________________
1900 void TGeant3::Gprint(const char *name)
1903 // Routine to print data structures
1904 // CHNAME name of a data structure
1908 gprint(PASSCHARD(vname),0 PASSCHARL(vname));
1911 //_____________________________________________________________________________
1912 void TGeant3::Grun()
1915 // Steering function to process one run
1920 //_____________________________________________________________________________
1921 void TGeant3::Gtrig()
1924 // Steering function to process one event
1929 //_____________________________________________________________________________
1930 void TGeant3::Gtrigc()
1933 // Clear event partition
1938 //_____________________________________________________________________________
1939 void TGeant3::Gtrigi()
1942 // Initialises event partition
1947 //_____________________________________________________________________________
1948 void TGeant3::Gwork(Int_t nwork)
1951 // Allocates workspace in ZEBRA memory
1956 //_____________________________________________________________________________
1957 void TGeant3::Gzinit()
1960 // To initialise GEANT/ZEBRA data structures
1965 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1967 // Functions from GCONS
1969 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1971 //_____________________________________________________________________________
1972 void TGeant3::Gfmate(Int_t imat, char *name, Float_t &a, Float_t &z,
1973 Float_t &dens, Float_t &radl, Float_t &absl,
1974 Float_t* ubuf, Int_t& nbuf)
1977 // Return parameters for material IMAT
1979 gfmate(imat, PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
1983 //_____________________________________________________________________________
1984 void TGeant3::Gfpart(Int_t ipart, char *name, Int_t &itrtyp,
1985 Float_t &amass, Float_t &charge, Float_t &tlife)
1988 // Return parameters for particle of type IPART
1992 Int_t igpart = IdFromPDG(ipart);
1993 gfpart(igpart, PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
1997 //_____________________________________________________________________________
1998 void TGeant3::Gftmed(Int_t numed, char *name, Int_t &nmat, Int_t &isvol,
1999 Int_t &ifield, Float_t &fieldm, Float_t &tmaxfd,
2000 Float_t &stemax, Float_t &deemax, Float_t &epsil,
2001 Float_t &stmin, Float_t *ubuf, Int_t *nbuf)
2004 // Return parameters for tracking medium NUMED
2006 gftmed(numed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
2007 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
2011 void TGeant3::Gftmat(Int_t imate, Int_t ipart, char *chmeca, Int_t kdim,
2012 Float_t* tkin, Float_t* value, Float_t* pcut,
2016 // Return parameters for tracking medium NUMED
2018 gftmat(imate, ipart, PASSCHARD(chmeca), kdim,
2019 tkin, value, pcut, ixst PASSCHARL(chmeca));
2023 //_____________________________________________________________________________
2024 Float_t TGeant3::Gbrelm(Float_t z, Float_t t, Float_t bcut)
2027 // To calculate energy loss due to soft muon BREMSSTRAHLUNG
2029 return gbrelm(z,t,bcut);
2032 //_____________________________________________________________________________
2033 Float_t TGeant3::Gprelm(Float_t z, Float_t t, Float_t bcut)
2036 // To calculate DE/DX in GeV*barn/atom for direct pair production by muons
2038 return gprelm(z,t,bcut);
2041 //_____________________________________________________________________________
2042 void TGeant3::Gmate()
2045 // Define standard GEANT materials
2050 //_____________________________________________________________________________
2051 void TGeant3::Gpart()
2054 // Define standard GEANT particles plus selected decay modes
2055 // and branching ratios.
2060 //_____________________________________________________________________________
2061 void TGeant3::Gsdk(Int_t ipart, Float_t *bratio, Int_t *mode)
2063 // Defines branching ratios and decay modes for standard
2065 gsdk(ipart,bratio,mode);
2068 //_____________________________________________________________________________
2069 void TGeant3::Gsmate(Int_t imat, const char *name, Float_t a, Float_t z,
2070 Float_t dens, Float_t radl, Float_t absl)
2073 // Defines a Material
2075 // kmat number assigned to the material
2076 // name material name
2077 // a atomic mass in au
2079 // dens density in g/cm3
2080 // absl absorbtion length in cm
2081 // if >=0 it is ignored and the program
2082 // calculates it, if <0. -absl is taken
2083 // radl radiation length in cm
2084 // if >=0 it is ignored and the program
2085 // calculates it, if <0. -radl is taken
2086 // buf pointer to an array of user words
2087 // nbuf number of user words
2091 gsmate(imat,PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
2095 //_____________________________________________________________________________
2096 void TGeant3::Gsmixt(Int_t imat, const char *name, Float_t *a, Float_t *z,
2097 Float_t dens, Int_t nlmat, Float_t *wmat)
2100 // Defines mixture OR COMPOUND IMAT as composed by
2101 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
2103 // If NLMAT.GT.0 then WMAT contains the PROPORTION BY
2104 // WEIGTHS OF EACH BASIC MATERIAL IN THE MIXTURE.
2106 // If NLMAT.LT.0 then WMAT contains the number of atoms
2107 // of a given kind into the molecule of the COMPOUND
2108 // In this case, WMAT in output is changed to relative
2111 gsmixt(imat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
2114 //_____________________________________________________________________________
2115 void TGeant3::Gspart(Int_t ipart, const char *name, Int_t itrtyp,
2116 Float_t amass, Float_t charge, Float_t tlife)
2119 // Store particle parameters
2121 // ipart particle code
2122 // name particle name
2123 // itrtyp transport method (see GEANT manual)
2124 // amass mass in GeV/c2
2125 // charge charge in electron units
2126 // tlife lifetime in seconds
2130 gspart(ipart,PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
2134 //_____________________________________________________________________________
2135 void TGeant3::Gstmed(Int_t numed, const char *name, Int_t nmat, Int_t isvol,
2136 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
2137 Float_t stemax, Float_t deemax, Float_t epsil,
2141 // NTMED Tracking medium number
2142 // NAME Tracking medium name
2143 // NMAT Material number
2144 // ISVOL Sensitive volume flag
2145 // IFIELD Magnetic field
2146 // FIELDM Max. field value (Kilogauss)
2147 // TMAXFD Max. angle due to field (deg/step)
2148 // STEMAX Max. step allowed
2149 // DEEMAX Max. fraction of energy lost in a step
2150 // EPSIL Tracking precision (cm)
2151 // STMIN Min. step due to continuos processes (cm)
2153 // IFIELD = 0 if no magnetic field; IFIELD = -1 if user decision in GUSWIM;
2154 // IFIELD = 1 if tracking performed with GRKUTA; IFIELD = 2 if tracking
2155 // performed with GHELIX; IFIELD = 3 if tracking performed with GHELX3.
2159 gstmed(numed,PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
2160 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
2163 //_____________________________________________________________________________
2164 void TGeant3::Gsckov(Int_t itmed, Int_t npckov, Float_t *ppckov,
2165 Float_t *absco, Float_t *effic, Float_t *rindex)
2168 // Stores the tables for UV photon tracking in medium ITMED
2169 // Please note that it is the user's responsability to
2170 // provide all the coefficients:
2173 // ITMED Tracking medium number
2174 // NPCKOV Number of bins of each table
2175 // PPCKOV Value of photon momentum (in GeV)
2176 // ABSCO Absorbtion coefficients
2177 // dielectric: absorbtion length in cm
2178 // metals : absorbtion fraction (0<=x<=1)
2179 // EFFIC Detection efficiency for UV photons
2180 // RINDEX Refraction index (if=0 metal)
2182 gsckov(itmed,npckov,ppckov,absco,effic,rindex);
2185 //_____________________________________________________________________________
2186 void TGeant3::SetCerenkov(Int_t itmed, Int_t npckov, Float_t *ppckov,
2187 Float_t *absco, Float_t *effic, Float_t *rindex)
2190 // Stores the tables for UV photon tracking in medium ITMED
2191 // Please note that it is the user's responsability to
2192 // provide all the coefficients:
2195 // ITMED Tracking medium number
2196 // NPCKOV Number of bins of each table
2197 // PPCKOV Value of photon momentum (in GeV)
2198 // ABSCO Absorbtion coefficients
2199 // dielectric: absorbtion length in cm
2200 // metals : absorbtion fraction (0<=x<=1)
2201 // EFFIC Detection efficiency for UV photons
2202 // RINDEX Refraction index (if=0 metal)
2204 gsckov(itmed,npckov,ppckov,absco,effic,rindex);
2207 //_____________________________________________________________________________
2208 void TGeant3::Gstpar(Int_t itmed, const char *param, Float_t parval)
2211 // To change the value of cut or mechanism "CHPAR"
2212 // to a new value PARVAL for tracking medium ITMED
2213 // The data structure JTMED contains the standard tracking
2214 // parameters (CUTS and flags to control the physics processes) which
2215 // are used by default for all tracking media. It is possible to
2216 // redefine individually with GSTPAR any of these parameters for a
2217 // given tracking medium.
2218 // ITMED tracking medium number
2219 // CHPAR is a character string (variable name)
2220 // PARVAL must be given as a floating point.
2222 gstpar(itmed,PASSCHARD(param), parval PASSCHARL(param));
2225 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2227 // Functions from GCONS
2229 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2231 //_____________________________________________________________________________
2232 void TGeant3::Gfkine(Int_t itra, Float_t *vert, Float_t *pvert, Int_t &ipart,
2235 // Storing/Retrieving Vertex and Track parameters
2236 // ----------------------------------------------
2238 // Stores vertex parameters.
2239 // VERT array of (x,y,z) position of the vertex
2240 // NTBEAM beam track number origin of the vertex
2241 // =0 if none exists
2242 // NTTARG target track number origin of the vertex
2243 // UBUF user array of NUBUF floating point numbers
2245 // NVTX new vertex number (=0 in case of error).
2246 // Prints vertex parameters.
2247 // IVTX for vertex IVTX.
2248 // (For all vertices if IVTX=0)
2249 // Stores long life track parameters.
2250 // PLAB components of momentum
2251 // IPART type of particle (see GSPART)
2252 // NV vertex number origin of track
2253 // UBUF array of NUBUF floating point user parameters
2255 // NT track number (if=0 error).
2256 // Retrieves long life track parameters.
2257 // ITRA track number for which parameters are requested
2258 // VERT vector origin of the track
2259 // PVERT 4 momentum components at the track origin
2260 // IPART particle type (=0 if track ITRA does not exist)
2261 // NVERT vertex number origin of the track
2262 // UBUF user words stored in GSKINE.
2263 // Prints initial track parameters.
2264 // ITRA for track ITRA
2265 // (For all tracks if ITRA=0)
2269 gfkine(itra,vert,pvert,ipart,nvert,ubuf,nbuf);
2272 //_____________________________________________________________________________
2273 void TGeant3::Gfvert(Int_t nvtx, Float_t *v, Int_t &ntbeam, Int_t &nttarg,
2277 // Retrieves the parameter of a vertex bank
2278 // Vertex is generated from tracks NTBEAM NTTARG
2279 // NVTX is the new vertex number
2283 gfvert(nvtx,v,ntbeam,nttarg,tofg,ubuf,nbuf);
2286 //_____________________________________________________________________________
2287 Int_t TGeant3::Gskine(Float_t *plab, Int_t ipart, Int_t nv, Float_t *buf,
2291 // Store kinematics of track NT into data structure
2292 // Track is coming from vertex NV
2295 gskine(plab, ipart, nv, buf, nwbuf, nt);
2299 //_____________________________________________________________________________
2300 Int_t TGeant3::Gsvert(Float_t *v, Int_t ntbeam, Int_t nttarg, Float_t *ubuf,
2304 // Creates a new vertex bank
2305 // Vertex is generated from tracks NTBEAM NTTARG
2306 // NVTX is the new vertex number
2309 gsvert(v, ntbeam, nttarg, ubuf, nwbuf, nwtx);
2313 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2315 // Functions from GPHYS
2317 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2319 //_____________________________________________________________________________
2320 void TGeant3::Gphysi()
2323 // Initialise material constants for all the physics
2324 // mechanisms used by GEANT
2329 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2331 // Functions from GTRAK
2333 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2335 //_____________________________________________________________________________
2336 void TGeant3::Gdebug()
2339 // Debug the current step
2344 //_____________________________________________________________________________
2345 void TGeant3::Gekbin()
2348 // To find bin number in kinetic energy table
2349 // stored in ELOW(NEKBIN)
2354 //_____________________________________________________________________________
2355 void TGeant3::Gfinds()
2358 // Returns the set/volume parameters corresponding to
2359 // the current space point in /GCTRAK/
2360 // and fill common /GCSETS/
2362 // IHSET user set identifier
2363 // IHDET user detector identifier
2364 // ISET set number in JSET
2365 // IDET detector number in JS=LQ(JSET-ISET)
2366 // IDTYPE detector type (1,2)
2367 // NUMBV detector volume numbers (array of length NVNAME)
2368 // NVNAME number of volume levels
2373 //_____________________________________________________________________________
2374 void TGeant3::Gsking(Int_t igk)
2377 // Stores in stack JSTAK either the IGKth track of /GCKING/,
2378 // or the NGKINE tracks when IGK is 0.
2383 //_____________________________________________________________________________
2384 void TGeant3::Gskpho(Int_t igk)
2387 // Stores in stack JSTAK either the IGKth Cherenkov photon of
2388 // /GCKIN2/, or the NPHOT tracks when IGK is 0.
2393 //_____________________________________________________________________________
2394 void TGeant3::Gsstak(Int_t iflag)
2397 // Stores in auxiliary stack JSTAK the particle currently
2398 // described in common /GCKINE/.
2400 // On request, creates also an entry in structure JKINE :
2402 // 0 : No entry in JKINE structure required (user)
2403 // 1 : New entry in JVERTX / JKINE structures required (user)
2404 // <0 : New entry in JKINE structure at vertex -IFLAG (user)
2405 // 2 : Entry in JKINE structure exists already (from GTREVE)
2410 //_____________________________________________________________________________
2411 void TGeant3::Gsxyz()
2414 // Store space point VECT in banks JXYZ
2419 //_____________________________________________________________________________
2420 void TGeant3::Gtrack()
2423 // Controls tracking of current particle
2428 //_____________________________________________________________________________
2429 void TGeant3::Gtreve()
2432 // Controls tracking of all particles belonging to the current event
2437 //_____________________________________________________________________________
2438 void TGeant3::GtreveRoot()
2441 // Controls tracking of all particles belonging to the current event
2446 //_____________________________________________________________________________
2447 void TGeant3::Grndm(Float_t *rvec, const Int_t len) const
2450 // To generate a vector RVECV of LEN random numbers
2451 // Copy of the CERN Library routine RANECU
2455 //_____________________________________________________________________________
2456 void TGeant3::Grndmq(Int_t &/*is1*/, Int_t &/*is2*/, const Int_t /*iseq*/,
2457 const Text_t */*chopt*/)
2460 // To set/retrieve the seed of the random number generator
2462 /*printf("Dummy grndmq called\n");*/
2465 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2467 // Functions from GDRAW
2469 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2471 //_____________________________________________________________________________
2472 void TGeant3::Gdxyz(Int_t it)
2475 // Draw the points stored with Gsxyz relative to track it
2480 //_____________________________________________________________________________
2481 void TGeant3::Gdcxyz()
2484 // Draw the position of the current track
2489 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2491 // Functions from GGEOM
2493 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2495 //_____________________________________________________________________________
2496 void TGeant3::Gdtom(Float_t *xd, Float_t *xm, Int_t iflag)
2499 // Computes coordinates XM (Master Reference System
2500 // knowing the coordinates XD (Detector Ref System)
2501 // The local reference system can be initialized by
2502 // - the tracking routines and GDTOM used in GUSTEP
2503 // - a call to GSCMED(NLEVEL,NAMES,NUMBER)
2504 // (inverse routine is GMTOD)
2506 // If IFLAG=1 convert coordinates
2507 // IFLAG=2 convert direction cosinus
2509 gdtom(xd, xm, iflag);
2512 //_____________________________________________________________________________
2513 void TGeant3::Glmoth(const char* iudet, Int_t iunum, Int_t &nlev, Int_t *lvols,
2517 // Loads the top part of the Volume tree in LVOLS (IVO's),
2518 // LINDX (IN indices) for a given volume defined through
2519 // its name IUDET and number IUNUM.
2521 // The routine stores only upto the last level where JVOLUM
2522 // data structure is developed. If there is no development
2523 // above the current level, it returns NLEV zero.
2525 glmoth(PASSCHARD(iudet), iunum, nlev, lvols, lindx, idum PASSCHARL(iudet));
2528 //_____________________________________________________________________________
2529 void TGeant3::Gmedia(Float_t *x, Int_t &numed)
2532 // Finds in which volume/medium the point X is, and updates the
2533 // common /GCVOLU/ and the structure JGPAR accordingly.
2535 // NUMED returns the tracking medium number, or 0 if point is
2536 // outside the experimental setup.
2541 //_____________________________________________________________________________
2542 void TGeant3::Gmtod(Float_t *xm, Float_t *xd, Int_t iflag)
2545 // Computes coordinates XD (in DRS)
2546 // from known coordinates XM in MRS
2547 // The local reference system can be initialized by
2548 // - the tracking routines and GMTOD used in GUSTEP
2549 // - a call to GMEDIA(XM,NUMED)
2550 // - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER)
2551 // (inverse routine is GDTOM)
2553 // If IFLAG=1 convert coordinates
2554 // IFLAG=2 convert direction cosinus
2556 gmtod(xm, xd, iflag);
2559 //_____________________________________________________________________________
2560 void TGeant3::Gsdvn(const char *name, const char *mother, Int_t ndiv,
2564 // Create a new volume by dividing an existing one
2567 // MOTHER Mother volume name
2568 // NDIV Number of divisions
2571 // X,Y,Z of CAXIS will be translated to 1,2,3 for IAXIS.
2572 // It divides a previously defined volume.
2577 Vname(mother,vmother);
2578 gsdvn(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis PASSCHARL(vname)
2579 PASSCHARL(vmother));
2582 //_____________________________________________________________________________
2583 void TGeant3::Gsdvn2(const char *name, const char *mother, Int_t ndiv,
2584 Int_t iaxis, Float_t c0i, Int_t numed)
2587 // Create a new volume by dividing an existing one
2589 // Divides mother into ndiv divisions called name
2590 // along axis iaxis starting at coordinate value c0.
2591 // the new volume created will be medium number numed.
2596 Vname(mother,vmother);
2597 gsdvn2(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis, c0i, numed
2598 PASSCHARL(vname) PASSCHARL(vmother));
2601 //_____________________________________________________________________________
2602 void TGeant3::Gsdvs(const char *name, const char *mother, Float_t step,
2603 Int_t iaxis, Int_t numed)
2606 // Create a new volume by dividing an existing one
2611 Vname(mother,vmother);
2612 gsdvs(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed
2613 PASSCHARL(vname) PASSCHARL(vmother));
2616 //_____________________________________________________________________________
2617 void TGeant3::Gsdvs2(const char *name, const char *mother, Float_t step,
2618 Int_t iaxis, Float_t c0, Int_t numed)
2621 // Create a new volume by dividing an existing one
2626 Vname(mother,vmother);
2627 gsdvs2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0, numed
2628 PASSCHARL(vname) PASSCHARL(vmother));
2631 //_____________________________________________________________________________
2632 void TGeant3::Gsdvt(const char *name, const char *mother, Float_t step,
2633 Int_t iaxis, Int_t numed, Int_t ndvmx)
2636 // Create a new volume by dividing an existing one
2638 // Divides MOTHER into divisions called NAME along
2639 // axis IAXIS in steps of STEP. If not exactly divisible
2640 // will make as many as possible and will centre them
2641 // with respect to the mother. Divisions will have medium
2642 // number NUMED. If NUMED is 0, NUMED of MOTHER is taken.
2643 // NDVMX is the expected maximum number of divisions
2644 // (If 0, no protection tests are performed)
2649 Vname(mother,vmother);
2650 gsdvt(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed, ndvmx
2651 PASSCHARL(vname) PASSCHARL(vmother));
2654 //_____________________________________________________________________________
2655 void TGeant3::Gsdvt2(const char *name, const char *mother, Float_t step,
2656 Int_t iaxis, Float_t c0, Int_t numed, Int_t ndvmx)
2659 // Create a new volume by dividing an existing one
2661 // Divides MOTHER into divisions called NAME along
2662 // axis IAXIS starting at coordinate value C0 with step
2664 // The new volume created will have medium number NUMED.
2665 // If NUMED is 0, NUMED of mother is taken.
2666 // NDVMX is the expected maximum number of divisions
2667 // (If 0, no protection tests are performed)
2672 Vname(mother,vmother);
2673 gsdvt2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0,
2674 numed, ndvmx PASSCHARL(vname) PASSCHARL(vmother));
2677 //_____________________________________________________________________________
2678 void TGeant3::Gsord(const char *name, Int_t iax)
2681 // Flags volume CHNAME whose contents will have to be ordered
2682 // along axis IAX, by setting the search flag to -IAX
2686 // IAX = 4 Rxy (static ordering only -> GTMEDI)
2687 // IAX = 14 Rxy (also dynamic ordering -> GTNEXT)
2688 // IAX = 5 Rxyz (static ordering only -> GTMEDI)
2689 // IAX = 15 Rxyz (also dynamic ordering -> GTNEXT)
2690 // IAX = 6 PHI (PHI=0 => X axis)
2691 // IAX = 7 THETA (THETA=0 => Z axis)
2695 gsord(PASSCHARD(vname), iax PASSCHARL(vname));
2698 //_____________________________________________________________________________
2699 void TGeant3::Gspos(const char *name, Int_t nr, const char *mother, Float_t x,
2700 Float_t y, Float_t z, Int_t irot, const char *konly)
2703 // Position a volume into an existing one
2706 // NUMBER Copy number of the volume
2707 // MOTHER Mother volume name
2708 // X X coord. of the volume in mother ref. sys.
2709 // Y Y coord. of the volume in mother ref. sys.
2710 // Z Z coord. of the volume in mother ref. sys.
2711 // IROT Rotation matrix number w.r.t. mother ref. sys.
2712 // ONLY ONLY/MANY flag
2714 // It positions a previously defined volume in the mother.
2720 Vname(mother,vmother);
2721 gspos(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2722 PASSCHARD(konly) PASSCHARL(vname) PASSCHARL(vmother)
2726 //_____________________________________________________________________________
2727 void TGeant3::Gsposp(const char *name, Int_t nr, const char *mother,
2728 Float_t x, Float_t y, Float_t z, Int_t irot,
2729 const char *konly, Float_t *upar, Int_t np )
2732 // Place a copy of generic volume NAME with user number
2733 // NR inside MOTHER, with its parameters UPAR(1..NP)
2738 Vname(mother,vmother);
2739 gsposp(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2740 PASSCHARD(konly), upar, np PASSCHARL(vname) PASSCHARL(vmother)
2744 //_____________________________________________________________________________
2745 void TGeant3::Gsrotm(Int_t nmat, Float_t theta1, Float_t phi1, Float_t theta2,
2746 Float_t phi2, Float_t theta3, Float_t phi3)
2749 // nmat Rotation matrix number
2750 // THETA1 Polar angle for axis I
2751 // PHI1 Azimuthal angle for axis I
2752 // THETA2 Polar angle for axis II
2753 // PHI2 Azimuthal angle for axis II
2754 // THETA3 Polar angle for axis III
2755 // PHI3 Azimuthal angle for axis III
2757 // It defines the rotation matrix number IROT.
2759 gsrotm(nmat, theta1, phi1, theta2, phi2, theta3, phi3);
2762 //_____________________________________________________________________________
2763 void TGeant3::Gprotm(Int_t nmat)
2766 // To print rotation matrices structure JROTM
2767 // nmat Rotation matrix number
2772 //_____________________________________________________________________________
2773 Int_t TGeant3::Gsvolu(const char *name, const char *shape, Int_t nmed,
2774 Float_t *upar, Int_t npar)
2778 // SHAPE Volume type
2779 // NUMED Tracking medium number
2780 // NPAR Number of shape parameters
2781 // UPAR Vector containing shape parameters
2783 // It creates a new volume in the JVOLUM data structure.
2789 Vname(shape,vshape);
2790 gsvolu(PASSCHARD(vname), PASSCHARD(vshape), nmed, upar, npar, ivolu
2791 PASSCHARL(vname) PASSCHARL(vshape));
2795 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2797 // T H E D R A W I N G P A C K A G E
2798 // ======================================
2799 // Drawing functions. These functions allow the visualization in several ways
2800 // of the volumes defined in the geometrical data structure. It is possible
2801 // to draw the logical tree of volumes belonging to the detector (DTREE),
2802 // to show their geometrical specification (DSPEC,DFSPC), to draw them
2803 // and their cut views (DRAW, DCUT). Moreover, it is possible to execute
2804 // these commands when the hidden line removal option is activated; in
2805 // this case, the volumes can be also either translated in the space
2806 // (SHIFT), or clipped by boolean operation (CVOL). In addition, it is
2807 // possible to fill the surfaces of the volumes
2808 // with solid colours when the shading option (SHAD) is activated.
2809 // Several tools (ZOOM, LENS) have been developed to zoom detailed parts
2810 // of the detectors or to scan physical events as well.
2811 // Finally, the command MOVE will allow the rotation, translation and zooming
2812 // on real time parts of the detectors or tracks and hits of a simulated event.
2813 // Ray-tracing commands. In case the command (DOPT RAYT ON) is executed,
2814 // the drawing is performed by the Geant ray-tracing;
2815 // automatically, the color is assigned according to the tracking medium of each
2816 // volume and the volumes with a density lower/equal than the air are considered
2817 // transparent; if the option (USER) is set (ON) (again via the command (DOPT)),
2818 // the user can set color and visibility for the desired volumes via the command
2819 // (SATT), as usual, relatively to the attributes (COLO) and (SEEN).
2820 // The resolution can be set via the command (SATT * FILL VALUE), where (VALUE)
2821 // is the ratio between the number of pixels drawn and 20 (user coordinates).
2822 // Parallel view and perspective view are possible (DOPT PROJ PARA/PERS); in the
2823 // first case, we assume that the first mother volume of the tree is a box with
2824 // dimensions 10000 X 10000 X 10000 cm and the view point (infinetely far) is
2825 // 5000 cm far from the origin along the Z axis of the user coordinates; in the
2826 // second case, the distance between the observer and the origin of the world
2827 // reference system is set in cm by the command (PERSP NAME VALUE); grand-angle
2828 // or telescopic effects can be achieved changing the scale factors in the command
2829 // (DRAW). When the final picture does not occupy the full window,
2830 // mapping the space before tracing can speed up the drawing, but can also
2831 // produce less precise results; values from 1 to 4 are allowed in the command
2832 // (DOPT MAPP VALUE), the mapping being more precise for increasing (VALUE); for
2833 // (VALUE = 0) no mapping is performed (therefore max precision and lowest speed).
2834 // The command (VALCUT) allows the cutting of the detector by three planes
2835 // ortogonal to the x,y,z axis. The attribute (LSTY) can be set by the command
2836 // SATT for any desired volume and can assume values from 0 to 7; it determines
2837 // the different light processing to be performed for different materials:
2838 // 0 = dark-matt, 1 = bright-matt, 2 = plastic, 3 = ceramic, 4 = rough-metals,
2839 // 5 = shiny-metals, 6 = glass, 7 = mirror. The detector is assumed to be in the
2840 // dark, the ambient light luminosity is 0.2 for each basic hue (the saturation
2841 // is 0.9) and the observer is assumed to have a light source (therefore he will
2842 // produce parallel light in the case of parallel view and point-like-source
2843 // light in the case of perspective view).
2845 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2847 //_____________________________________________________________________________
2848 void TGeant3::Gsatt(const char *name, const char *att, Int_t val)
2852 // IOPT Name of the attribute to be set
2853 // IVAL Value to which the attribute is to be set
2855 // name= "*" stands for all the volumes.
2856 // iopt can be chosen among the following :
2858 // WORK 0=volume name is inactive for the tracking
2859 // 1=volume name is active for the tracking (default)
2861 // SEEN 0=volume name is invisible
2862 // 1=volume name is visible (default)
2863 // -1=volume invisible with all its descendants in the tree
2864 // -2=volume visible but not its descendants in the tree
2866 // LSTY line style 1,2,3,... (default=1)
2867 // LSTY=7 will produce a very precise approximation for
2868 // revolution bodies.
2870 // LWID line width -7,...,1,2,3,..7 (default=1)
2871 // LWID<0 will act as abs(LWID) was set for the volume
2872 // and for all the levels below it. When SHAD is 'ON', LWID
2873 // represent the linewidth of the scan lines filling the surfaces
2874 // (whereas the FILL value represent their number). Therefore
2875 // tuning this parameter will help to obtain the desired
2876 // quality/performance ratio.
2878 // COLO colour code -166,...,1,2,..166 (default=1)
2880 // n=2=red; n=17+m, m=0,25, increasing luminosity according to 'm';
2881 // n=3=green; n=67+m, m=0,25, increasing luminosity according to 'm';
2882 // n=4=blue; n=117+m, m=0,25, increasing luminosity according to 'm';
2883 // n=5=yellow; n=42+m, m=0,25, increasing luminosity according to 'm';
2884 // n=6=violet; n=142+m, m=0,25, increasing luminosity according to 'm';
2885 // n=7=lightblue; n=92+m, m=0,25, increasing luminosity according to 'm';
2886 // colour=n*10+m, m=1,2,...9, will produce the same colour
2887 // as 'n', but with increasing luminosity according to 'm';
2888 // COLO<0 will act as if abs(COLO) was set for the volume
2889 // and for all the levels below it.
2890 // When for a volume the attribute FILL is > 1 (and the
2891 // option SHAD is on), the ABS of its colour code must be < 8
2892 // because an automatic shading of its faces will be
2895 // FILL (1992) fill area -7,...,0,1,...7 (default=0)
2896 // when option SHAD is "on" the FILL attribute of any
2897 // volume can be set different from 0 (normal drawing);
2898 // if it is set to 1, the faces of such volume will be filled
2899 // with solid colours; if ABS(FILL) is > 1, then a light
2900 // source is placed along the observer line, and the faces of
2901 // such volumes will be painted by colours whose luminosity
2902 // will depend on the amount of light reflected;
2903 // if ABS(FILL) = 1, then it is possible to use all the 166
2904 // colours of the colour table, becouse the automatic shading
2905 // is not performed;
2906 // for increasing values of FILL the drawing will be performed
2907 // with higher and higher resolution improving the quality (the
2908 // number of scan lines used to fill the faces increases with FILL);
2909 // it is possible to set different values of FILL
2910 // for different volumes, in order to optimize at the same time
2911 // the performance and the quality of the picture;
2912 // FILL<0 will act as if abs(FILL) was set for the volume
2913 // and for all the levels below it.
2914 // This kind of drawing can be saved in 'picture files'
2915 // or in view banks.
2916 // 0=drawing without fill area
2917 // 1=faces filled with solid colours and resolution = 6
2918 // 2=lowest resolution (very fast)
2919 // 3=default resolution
2920 // 4=.................
2921 // 5=.................
2922 // 6=.................
2924 // Finally, if a coloured background is desired, the FILL
2925 // attribute for the first volume of the tree must be set
2926 // equal to -abs(colo), colo being >0 and <166.
2928 // SET set number associated to volume name
2929 // DET detector number associated to volume name
2930 // DTYP detector type (1,2)
2937 gsatt(PASSCHARD(vname), PASSCHARD(vatt), val PASSCHARL(vname)
2941 //_____________________________________________________________________________
2942 void TGeant3::Gfpara(const char *name, Int_t number, Int_t intext, Int_t& npar,
2943 Int_t& natt, Float_t* par, Float_t* att)
2946 // Find the parameters of a volume
2948 gfpara(PASSCHARD(name), number, intext, npar, natt, par, att
2952 //_____________________________________________________________________________
2953 void TGeant3::Gckpar(Int_t ish, Int_t npar, Float_t* par)
2956 // Check the parameters of a shape
2958 gckpar(ish,npar,par);
2961 //_____________________________________________________________________________
2962 void TGeant3::Gckmat(Int_t itmed, char* natmed)
2965 // Check the parameters of a tracking medium
2967 gckmat(itmed, PASSCHARD(natmed) PASSCHARL(natmed));
2970 //_____________________________________________________________________________
2971 Int_t TGeant3::Glvolu(Int_t nlev, Int_t *lnam,Int_t *lnum)
2974 // nlev number of leveles deap into the volume tree
2975 // size of the arrays lnam and lnum
2976 // lnam an integer array whos 4 bytes contain the askii code for the
2978 // lnum an integer array containing the copy numbers for that specific
2981 // This routine fills the volulme paramters in common /gcvolu/ for a
2982 // physical tree, specified by the list lnam and lnum of volume names
2983 // and numbers, and for all its ascendants up to level 1. This routine
2984 // is optimsed and does not re-compute the part of the history already
2985 // available in GCVOLU. This means that if it is used in user programs
2986 // outside the usual framwork of the tracking, the user has to initilise
2987 // to zero NLEVEL in the common GCVOLU. It return 0 if there were no
2988 // problems in make the call.
2991 glvolu(nlev, lnam, lnum, ier);
2995 //_____________________________________________________________________________
2996 void TGeant3::Gdelete(Int_t iview)
2999 // IVIEW View number
3001 // It deletes a view bank from memory.
3006 //_____________________________________________________________________________
3007 void TGeant3::Gdopen(Int_t iview)
3010 // IVIEW View number
3012 // When a drawing is very complex and requires a long time to be
3013 // executed, it can be useful to store it in a view bank: after a
3014 // call to DOPEN and the execution of the drawing (nothing will
3015 // appear on the screen), and after a necessary call to DCLOSE,
3016 // the contents of the bank can be displayed in a very fast way
3017 // through a call to DSHOW; therefore, the detector can be easily
3018 // zoomed many times in different ways. Please note that the pictures
3019 // with solid colours can now be stored in a view bank or in 'PICTURE FILES'
3026 //_____________________________________________________________________________
3027 void TGeant3::Gdclose()
3030 // It closes the currently open view bank; it must be called after the
3031 // end of the drawing to be stored.
3036 //_____________________________________________________________________________
3037 void TGeant3::Gdshow(Int_t iview)
3040 // IVIEW View number
3042 // It shows on the screen the contents of a view bank. It
3043 // can be called after a view bank has been closed.
3048 //_____________________________________________________________________________
3049 void TGeant3::Gdopt(const char *name,const char *value)
3053 // VALUE Option value
3055 // To set/modify the drawing options.
3058 // THRZ ON Draw tracks in R vs Z
3059 // OFF (D) Draw tracks in X,Y,Z
3062 // PROJ PARA (D) Parallel projection
3064 // TRAK LINE (D) Trajectory drawn with lines
3065 // POIN " " with markers
3066 // HIDE ON Hidden line removal using the CG package
3067 // OFF (D) No hidden line removal
3068 // SHAD ON Fill area and shading of surfaces.
3069 // OFF (D) Normal hidden line removal.
3070 // RAYT ON Ray-tracing on.
3071 // OFF (D) Ray-tracing off.
3072 // EDGE OFF Does not draw contours when shad is on.
3073 // ON (D) Normal shading.
3074 // MAPP 1,2,3,4 Mapping before ray-tracing.
3075 // 0 (D) No mapping.
3076 // USER ON User graphics options in the raytracing.
3077 // OFF (D) Automatic graphics options.
3083 Vname(value,vvalue);
3084 gdopt(PASSCHARD(vname), PASSCHARD(vvalue) PASSCHARL(vname)
3088 //_____________________________________________________________________________
3089 void TGeant3::Gdraw(const char *name,Float_t theta, Float_t phi, Float_t psi,
3090 Float_t u0,Float_t v0,Float_t ul,Float_t vl)
3095 // THETA Viewing angle theta (for 3D projection)
3096 // PHI Viewing angle phi (for 3D projection)
3097 // PSI Viewing angle psi (for 2D rotation)
3098 // U0 U-coord. (horizontal) of volume origin
3099 // V0 V-coord. (vertical) of volume origin
3100 // SU Scale factor for U-coord.
3101 // SV Scale factor for V-coord.
3103 // This function will draw the volumes,
3104 // selected with their graphical attributes, set by the Gsatt
3105 // facility. The drawing may be performed with hidden line removal
3106 // and with shading effects according to the value of the options HIDE
3107 // and SHAD; if the option SHAD is ON, the contour's edges can be
3108 // drawn or not. If the option HIDE is ON, the detector can be
3109 // exploded (BOMB), clipped with different shapes (CVOL), and some
3110 // of its parts can be shifted from their original
3111 // position (SHIFT). When HIDE is ON, if
3112 // the drawing requires more than the available memory, the program
3113 // will evaluate and display the number of missing words
3114 // (so that the user can increase the
3115 // size of its ZEBRA store). Finally, at the end of each drawing (with HIDE on),
3116 // the program will print messages about the memory used and
3117 // statistics on the volumes' visibility.
3118 // The following commands will produce the drawing of a green
3119 // volume, specified by NAME, without using the hidden line removal
3120 // technique, using the hidden line removal technique,
3121 // with different linewidth and colour (red), with
3122 // solid colour, with shading of surfaces, and without edges.
3123 // Finally, some examples are given for the ray-tracing. (A possible
3124 // string for the NAME of the volume can be found using the command DTREE).
3130 if (fGcvdma->raytra != 1) {
3131 gdraw(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
3133 gdrayt(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
3137 //_____________________________________________________________________________
3138 void TGeant3::Gdrawc(const char *name,Int_t axis, Float_t cut,Float_t u0,
3139 Float_t v0,Float_t ul,Float_t vl)
3144 // CUTVAL Cut plane distance from the origin along the axis
3146 // U0 U-coord. (horizontal) of volume origin
3147 // V0 V-coord. (vertical) of volume origin
3148 // SU Scale factor for U-coord.
3149 // SV Scale factor for V-coord.
3151 // The cut plane is normal to caxis (X,Y,Z), corresponding to iaxis (1,2,3),
3152 // and placed at the distance cutval from the origin.
3153 // The resulting picture is seen from the the same axis.
3154 // When HIDE Mode is ON, it is possible to get the same effect with
3155 // the CVOL/BOX function.
3161 gdrawc(PASSCHARD(vname), axis,cut,u0,v0,ul,vl PASSCHARL(vname));
3164 //_____________________________________________________________________________
3165 void TGeant3::Gdrawx(const char *name,Float_t cutthe, Float_t cutphi,
3166 Float_t cutval, Float_t theta, Float_t phi, Float_t u0,
3167 Float_t v0,Float_t ul,Float_t vl)
3171 // CUTTHE Theta angle of the line normal to cut plane
3172 // CUTPHI Phi angle of the line normal to cut plane
3173 // CUTVAL Cut plane distance from the origin along the axis
3175 // THETA Viewing angle theta (for 3D projection)
3176 // PHI Viewing angle phi (for 3D projection)
3177 // U0 U-coord. (horizontal) of volume origin
3178 // V0 V-coord. (vertical) of volume origin
3179 // SU Scale factor for U-coord.
3180 // SV Scale factor for V-coord.
3182 // The cut plane is normal to the line given by the cut angles
3183 // cutthe and cutphi and placed at the distance cutval from the origin.
3184 // The resulting picture is seen from the viewing angles theta,phi.
3190 gdrawx(PASSCHARD(vname), cutthe,cutphi,cutval,theta,phi,u0,v0,ul,vl
3194 //_____________________________________________________________________________
3195 void TGeant3::Gdhead(Int_t isel, const char *name, Float_t chrsiz)
3200 // ISEL Option flag D=111110
3202 // CHRSIZ Character size (cm) of title NAME D=0.6
3205 // 0 to have only the header lines
3206 // xxxxx1 to add the text name centered on top of header
3207 // xxxx1x to add global detector name (first volume) on left
3208 // xxx1xx to add date on right
3209 // xx1xxx to select thick characters for text on top of header
3210 // x1xxxx to add the text 'EVENT NR x' on top of header
3211 // 1xxxxx to add the text 'RUN NR x' on top of header
3212 // NOTE that ISEL=x1xxx1 or ISEL=1xxxx1 are illegal choices,
3213 // i.e. they generate overwritten text.
3215 gdhead(isel,PASSCHARD(name),chrsiz PASSCHARL(name));
3218 //_____________________________________________________________________________
3219 void TGeant3::Gdman(Float_t u, Float_t v, const char *type)
3222 // Draw a 2D-man at position (U0,V0)
3224 // U U-coord. (horizontal) of the centre of man' R
3225 // V V-coord. (vertical) of the centre of man' R
3226 // TYPE D='MAN' possible values: 'MAN,WM1,WM2,WM3'
3228 // CALL GDMAN(u,v),CALL GDWMN1(u,v),CALL GDWMN2(u,v),CALL GDWMN2(u,v)
3229 // It superimposes the picure of a man or of a woman, chosen among
3230 // three different ones, with the same scale factors as the detector
3231 // in the current drawing.
3234 if (opt.Contains("WM1")) {
3236 } else if (opt.Contains("WM3")) {
3238 } else if (opt.Contains("WM2")) {
3245 //_____________________________________________________________________________
3246 void TGeant3::Gdspec(const char *name)
3251 // Shows 3 views of the volume (two cut-views and a 3D view), together with
3252 // its geometrical specifications. The 3D drawing will
3253 // be performed according the current values of the options HIDE and
3254 // SHAD and according the current SetClipBox clipping parameters for that
3261 gdspec(PASSCHARD(vname) PASSCHARL(vname));
3264 //_____________________________________________________________________________
3265 void TGeant3::DrawOneSpec(const char *name)
3268 // Function called when one double-clicks on a volume name
3269 // in a TPavelabel drawn by Gdtree.
3271 THIGZ *higzSave = gHigz;
3272 higzSave->SetName("higzSave");
3273 THIGZ *higzSpec = (THIGZ*)gROOT->FindObject("higzSpec");
3274 //printf("DrawOneSpec, gHigz=%x, higzSpec=%x\n",gHigz,higzSpec);
3275 if (higzSpec) gHigz = higzSpec;
3276 else higzSpec = new THIGZ(kDefSize);
3277 higzSpec->SetName("higzSpec");
3282 gdspec(PASSCHARD(vname) PASSCHARL(vname));
3285 higzSave->SetName("higz");
3289 //_____________________________________________________________________________
3290 void TGeant3::Gdtree(const char *name,Int_t levmax, Int_t isel)
3294 // LEVMAX Depth level
3297 // This function draws the logical tree,
3298 // Each volume in the tree is represented by a TPaveTree object.
3299 // Double-clicking on a TPaveTree draws the specs of the corresponding volume.
3300 // Use TPaveTree pop-up menu to select:
3303 // - drawing tree of parent
3309 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
3310 gHigz->SetPname("");
3313 //_____________________________________________________________________________
3314 void TGeant3::GdtreeParent(const char *name,Int_t levmax, Int_t isel)
3318 // LEVMAX Depth level
3321 // This function draws the logical tree of the parent of name.
3325 // Scan list of volumes in JVOLUM
3327 Int_t gname, i, jvo, in, nin, jin, num;
3328 strncpy((char *) &gname, name, 4);
3329 for(i=1; i<=fGcnum->nvolum; i++) {
3330 jvo = fZlq[fGclink->jvolum-i];
3331 nin = Int_t(fZq[jvo+3]);
3332 if (nin == -1) nin = 1;
3333 for (in=1;in<=nin;in++) {
3335 num = Int_t(fZq[jin+2]);
3336 if(gname == fZiq[fGclink->jvolum+num]) {
3337 strncpy(vname,(char*)&fZiq[fGclink->jvolum+i],4);
3339 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
3340 gHigz->SetPname("");
3347 //_____________________________________________________________________________
3348 void TGeant3::SetABAN(Int_t par)
3351 // par = 1 particles will be stopped according to their residual
3352 // range if they are not in a sensitive material and are
3353 // far enough from the boundary
3354 // 0 particles are transported normally
3356 fGcphys->dphys1 = par;
3360 //_____________________________________________________________________________
3361 void TGeant3::SetANNI(Int_t par)
3364 // To control positron annihilation.
3365 // par =0 no annihilation
3366 // =1 annihilation. Decays processed.
3367 // =2 annihilation. No decay products stored.
3369 fGcphys->ianni = par;
3373 //_____________________________________________________________________________
3374 void TGeant3::SetAUTO(Int_t par)
3377 // To control automatic calculation of tracking medium parameters:
3378 // par =0 no automatic calculation;
3379 // =1 automati calculation.
3381 fGctrak->igauto = par;
3385 //_____________________________________________________________________________
3386 void TGeant3::SetBOMB(Float_t boom)
3389 // BOOM : Exploding factor for volumes position
3391 // To 'explode' the detector. If BOOM is positive (values smaller
3392 // than 1. are suggested, but any value is possible)
3393 // all the volumes are shifted by a distance
3394 // proportional to BOOM along the direction between their centre
3395 // and the origin of the MARS; the volumes which are symmetric
3396 // with respect to this origin are simply not shown.
3397 // BOOM equal to 0 resets the normal mode.
3398 // A negative (greater than -1.) value of
3399 // BOOM will cause an 'implosion'; for even lower values of BOOM
3400 // the volumes' positions will be reflected respect to the origin.
3401 // This command can be useful to improve the 3D effect for very
3402 // complex detectors. The following commands will make explode the
3409 //_____________________________________________________________________________
3410 void TGeant3::SetBREM(Int_t par)
3413 // To control bremstrahlung.
3414 // par =0 no bremstrahlung
3415 // =1 bremstrahlung. Photon processed.
3416 // =2 bremstrahlung. No photon stored.
3418 fGcphys->ibrem = par;
3422 //_____________________________________________________________________________
3423 void TGeant3::SetCKOV(Int_t par)
3426 // To control Cerenkov production
3427 // par =0 no Cerenkov;
3429 // =2 Cerenkov with primary stopped at each step.
3431 fGctlit->itckov = par;
3435 //_____________________________________________________________________________
3436 void TGeant3::SetClipBox(const char *name,Float_t xmin,Float_t xmax,
3437 Float_t ymin,Float_t ymax,Float_t zmin,Float_t zmax)
3440 // The hidden line removal technique is necessary to visualize properly
3441 // very complex detectors. At the same time, it can be useful to visualize
3442 // the inner elements of a detector in detail. This function allows
3443 // subtractions (via boolean operation) of BOX shape from any part of
3444 // the detector, therefore showing its inner contents.
3445 // If "*" is given as the name of the
3446 // volume to be clipped, all volumes are clipped by the given box.
3447 // A volume can be clipped at most twice.
3448 // if a volume is explicitely clipped twice,
3449 // the "*" will not act on it anymore. Giving "." as the name
3450 // of the volume to be clipped will reset the clipping.
3452 // NAME Name of volume to be clipped
3454 // XMIN Lower limit of the Shape X coordinate
3455 // XMAX Upper limit of the Shape X coordinate
3456 // YMIN Lower limit of the Shape Y coordinate
3457 // YMAX Upper limit of the Shape Y coordinate
3458 // ZMIN Lower limit of the Shape Z coordinate
3459 // ZMAX Upper limit of the Shape Z coordinate
3461 // This function performs a boolean subtraction between the volume
3462 // NAME and a box placed in the MARS according the values of the given
3468 setclip(PASSCHARD(vname),xmin,xmax,ymin,ymax,zmin,zmax PASSCHARL(vname));
3471 //_____________________________________________________________________________
3472 void TGeant3::SetCOMP(Int_t par)
3475 // To control Compton scattering
3476 // par =0 no Compton
3477 // =1 Compton. Electron processed.
3478 // =2 Compton. No electron stored.
3481 fGcphys->icomp = par;
3484 //_____________________________________________________________________________
3485 void TGeant3::SetCUTS(Float_t cutgam,Float_t cutele,Float_t cutneu,
3486 Float_t cuthad,Float_t cutmuo ,Float_t bcute ,
3487 Float_t bcutm ,Float_t dcute ,Float_t dcutm ,
3488 Float_t ppcutm, Float_t tofmax)
3491 // CUTGAM Cut for gammas D=0.001
3492 // CUTELE Cut for electrons D=0.001
3493 // CUTHAD Cut for charged hadrons D=0.01
3494 // CUTNEU Cut for neutral hadrons D=0.01
3495 // CUTMUO Cut for muons D=0.01
3496 // BCUTE Cut for electron brems. D=-1.
3497 // BCUTM Cut for muon brems. D=-1.
3498 // DCUTE Cut for electron delta-rays D=-1.
3499 // DCUTM Cut for muon delta-rays D=-1.
3500 // PPCUTM Cut for e+e- pairs by muons D=0.01
3501 // TOFMAX Time of flight cut D=1.E+10
3503 // If the default values (-1.) for BCUTE ,BCUTM ,DCUTE ,DCUTM
3504 // are not modified, they will be set to CUTGAM,CUTGAM,CUTELE,CUTELE
3506 // If one of the parameters from CUTGAM to PPCUTM included
3507 // is modified, cross-sections and energy loss tables must be
3508 // recomputed via the function Gphysi.
3510 fGccuts->cutgam = cutgam;
3511 fGccuts->cutele = cutele;
3512 fGccuts->cutneu = cutneu;
3513 fGccuts->cuthad = cuthad;
3514 fGccuts->cutmuo = cutmuo;
3515 fGccuts->bcute = bcute;
3516 fGccuts->bcutm = bcutm;
3517 fGccuts->dcute = dcute;
3518 fGccuts->dcutm = dcutm;
3519 fGccuts->ppcutm = ppcutm;
3520 fGccuts->tofmax = tofmax;
3523 //_____________________________________________________________________________
3524 void TGeant3::SetDCAY(Int_t par)
3527 // To control Decay mechanism.
3528 // par =0 no decays.
3529 // =1 Decays. secondaries processed.
3530 // =2 Decays. No secondaries stored.
3532 fGcphys->idcay = par;
3536 //_____________________________________________________________________________
3537 void TGeant3::SetDEBU(Int_t emin, Int_t emax, Int_t emod)
3540 // Set the debug flag and frequency
3541 // Selected debug output will be printed from
3542 // event emin to even emax each emod event
3544 fGcflag->idemin = emin;
3545 fGcflag->idemax = emax;
3546 fGcflag->itest = emod;
3550 //_____________________________________________________________________________
3551 void TGeant3::SetDRAY(Int_t par)
3554 // To control delta rays mechanism.
3555 // par =0 no delta rays.
3556 // =1 Delta rays. secondaries processed.
3557 // =2 Delta rays. No secondaries stored.
3559 fGcphys->idray = par;
3562 //_____________________________________________________________________________
3563 void TGeant3::SetERAN(Float_t ekmin, Float_t ekmax, Int_t nekbin)
3566 // To control cross section tabulations
3567 // ekmin = minimum kinetic energy in GeV
3568 // ekmax = maximum kinetic energy in GeV
3569 // nekbin = number of logatithmic bins (<200)
3571 fGcmulo->ekmin = ekmin;
3572 fGcmulo->ekmax = ekmax;
3573 fGcmulo->nekbin = nekbin;
3576 //_____________________________________________________________________________
3577 void TGeant3::SetHADR(Int_t par)
3580 // To control hadronic interactions.
3581 // par =0 no hadronic interactions.
3582 // =1 Hadronic interactions. secondaries processed.
3583 // =2 Hadronic interactions. No secondaries stored.
3585 fGcphys->ihadr = par;
3588 //_____________________________________________________________________________
3589 void TGeant3::SetKINE(Int_t kine, Float_t xk1, Float_t xk2, Float_t xk3,
3590 Float_t xk4, Float_t xk5, Float_t xk6, Float_t xk7,
3591 Float_t xk8, Float_t xk9, Float_t xk10)
3594 // Set the variables in /GCFLAG/ IKINE, PKINE(10)
3595 // Their meaning is user defined
3597 fGckine->ikine = kine;
3598 fGckine->pkine[0] = xk1;
3599 fGckine->pkine[1] = xk2;
3600 fGckine->pkine[2] = xk3;
3601 fGckine->pkine[3] = xk4;
3602 fGckine->pkine[4] = xk5;
3603 fGckine->pkine[5] = xk6;
3604 fGckine->pkine[6] = xk7;
3605 fGckine->pkine[7] = xk8;
3606 fGckine->pkine[8] = xk9;
3607 fGckine->pkine[9] = xk10;
3610 //_____________________________________________________________________________
3611 void TGeant3::SetLOSS(Int_t par)
3614 // To control energy loss.
3615 // par =0 no energy loss;
3616 // =1 restricted energy loss fluctuations;
3617 // =2 complete energy loss fluctuations;
3619 // =4 no energy loss fluctuations.
3620 // If the value ILOSS is changed, then cross-sections and energy loss
3621 // tables must be recomputed via the command 'PHYSI'.
3623 fGcphys->iloss = par;
3627 //_____________________________________________________________________________
3628 void TGeant3::SetMULS(Int_t par)
3631 // To control multiple scattering.
3632 // par =0 no multiple scattering.
3633 // =1 Moliere or Coulomb scattering.
3634 // =2 Moliere or Coulomb scattering.
3635 // =3 Gaussian scattering.
3637 fGcphys->imuls = par;
3641 //_____________________________________________________________________________
3642 void TGeant3::SetMUNU(Int_t par)
3645 // To control muon nuclear interactions.
3646 // par =0 no muon-nuclear interactions.
3647 // =1 Nuclear interactions. Secondaries processed.
3648 // =2 Nuclear interactions. Secondaries not processed.
3650 fGcphys->imunu = par;
3653 //_____________________________________________________________________________
3654 void TGeant3::SetOPTI(Int_t par)
3657 // This flag controls the tracking optimisation performed via the
3659 // 1 no optimisation at all; GSORD calls disabled;
3660 // 0 no optimisation; only user calls to GSORD kept;
3661 // 1 all non-GSORDered volumes are ordered along the best axis;
3662 // 2 all volumes are ordered along the best axis.
3664 fGcopti->ioptim = par;
3667 //_____________________________________________________________________________
3668 void TGeant3::SetPAIR(Int_t par)
3671 // To control pair production mechanism.
3672 // par =0 no pair production.
3673 // =1 Pair production. secondaries processed.
3674 // =2 Pair production. No secondaries stored.
3676 fGcphys->ipair = par;
3680 //_____________________________________________________________________________
3681 void TGeant3::SetPFIS(Int_t par)
3684 // To control photo fission mechanism.
3685 // par =0 no photo fission.
3686 // =1 Photo fission. secondaries processed.
3687 // =2 Photo fission. No secondaries stored.
3689 fGcphys->ipfis = par;
3692 //_____________________________________________________________________________
3693 void TGeant3::SetPHOT(Int_t par)
3696 // To control Photo effect.
3697 // par =0 no photo electric effect.
3698 // =1 Photo effect. Electron processed.
3699 // =2 Photo effect. No electron stored.
3701 fGcphys->iphot = par;
3704 //_____________________________________________________________________________
3705 void TGeant3::SetRAYL(Int_t par)
3708 // To control Rayleigh scattering.
3709 // par =0 no Rayleigh scattering.
3712 fGcphys->irayl = par;
3715 //_____________________________________________________________________________
3716 void TGeant3::SetSTRA(Int_t par)
3719 // To control energy loss fluctuations
3720 // with the PhotoAbsorption Ionisation model.
3721 // par =0 no Straggling.
3722 // =1 Straggling yes => no Delta rays.
3724 fGcphlt->istra = par;
3727 //_____________________________________________________________________________
3728 void TGeant3::SetSWIT(Int_t sw, Int_t val)
3732 // val New switch value
3734 // Change one element of array ISWIT(10) in /GCFLAG/
3736 if (sw <= 0 || sw > 10) return;
3737 fGcflag->iswit[sw-1] = val;
3741 //_____________________________________________________________________________
3742 void TGeant3::SetTRIG(Int_t nevents)
3745 // Set number of events to be run
3747 fGcflag->nevent = nevents;
3750 //_____________________________________________________________________________
3751 void TGeant3::SetUserDecay(Int_t pdg)
3754 // Force the decays of particles to be done with Pythia
3755 // and not with the Geant routines.
3756 // just kill pointers doing mzdrop
3758 Int_t ipart = IdFromPDG(pdg);
3760 printf("Particle %d not in geant\n",pdg);
3763 Int_t jpart=fGclink->jpart;
3764 Int_t jpa=fZlq[jpart-ipart];
3767 Int_t jpa1=fZlq[jpa-1];
3769 mzdrop(fGcbank->ixcons,jpa1,PASSCHARD(" ") PASSCHARL(" "));
3770 Int_t jpa2=fZlq[jpa-2];
3772 mzdrop(fGcbank->ixcons,jpa2,PASSCHARD(" ") PASSCHARL(" "));
3776 //______________________________________________________________________________
3777 void TGeant3::Vname(const char *name, char *vname)
3780 // convert name to upper case. Make vname at least 4 chars
3782 Int_t l = strlen(name);
3785 for (i=0;i<l;i++) vname[i] = toupper(name[i]);
3786 for (i=l;i<4;i++) vname[i] = ' ';
3790 //______________________________________________________________________________
3791 void TGeant3::Ertrgo()
3794 // Perform the tracking of the track Track parameters are in VECT
3799 //______________________________________________________________________________
3800 void TGeant3::Ertrak(const Float_t *const x1, const Float_t *const p1,
3801 const Float_t *x2, const Float_t *p2,
3802 Int_t ipa, Option_t *chopt)
3804 //************************************************************************
3806 //* Perform the tracking of the track from point X1 to *
3808 //* (Before calling this routine the user should also provide *
3809 //* the input informations in /EROPTS/ and /ERTRIO/ *
3810 //* using subroutine EUFIL(L/P/V) *
3811 //* X1 - Starting coordinates (Cartesian) *
3812 //* P1 - Starting 3-momentum (Cartesian) *
3813 //* X2 - Final coordinates (Cartesian) *
3814 //* P2 - Final 3-momentum (Cartesian) *
3815 //* IPA - Particle code (a la GEANT) of the track *
3818 //* 'B' 'Backward tracking' - i.e. energy loss *
3819 //* added to the current energy *
3820 //* 'E' 'Exact' calculation of errors assuming *
3821 //* helix (i.e. pathlength not *
3822 //* assumed as infinitesimal) *
3823 //* 'L' Tracking upto prescribed Lengths reached *
3824 //* 'M' 'Mixed' prediction (not yet coded) *
3825 //* 'O' Tracking 'Only' without calculating errors *
3826 //* 'P' Tracking upto prescribed Planes reached *
3827 //* 'V' Tracking upto prescribed Volumes reached *
3828 //* 'X' Tracking upto prescribed Point approached *
3830 //* Interface with GEANT : *
3831 //* Track parameters are in /CGKINE/ and /GCTRAK/ *
3833 //* ==>Called by : USER *
3834 //* Authors M.Maire, E.Nagy ********//* *
3836 //************************************************************************
3837 ertrak(x1,p1,x2,p2,ipa,PASSCHARD(chopt) PASSCHARL(chopt));
3840 //_____________________________________________________________________________
3841 void TGeant3::WriteEuclid(const char* filnam, const char* topvol,
3842 Int_t number, Int_t nlevel)
3846 // ******************************************************************
3848 // * Write out the geometry of the detector in EUCLID file format *
3850 // * filnam : will be with the extension .euc *
3851 // * topvol : volume name of the starting node *
3852 // * number : copy number of topvol (relevant for gsposp) *
3853 // * nlevel : number of levels in the tree structure *
3854 // * to be written out, starting from topvol *
3856 // * Author : M. Maire *
3858 // ******************************************************************
3860 // File filnam.tme is written out with the definitions of tracking
3861 // medias and materials.
3862 // As to restore original numbers for materials and medias, program
3863 // searches in the file euc_medi.dat and comparing main parameters of
3864 // the mat. defined inside geant and the one in file recognizes them
3865 // and is able to take number from file. If for any material or medium,
3866 // this procedure fails, ordering starts from 1.
3867 // Arrays IOTMED and IOMATE are used for this procedure
3869 const char kShape[][5]={"BOX ","TRD1","TRD2","TRAP","TUBE","TUBS","CONE",
3870 "CONS","SPHE","PARA","PGON","PCON","ELTU","HYPE",
3872 Int_t i, end, itm, irm, jrm, k, nmed;
3876 char *filext, *filetme;
3877 char natmed[21], namate[21];
3878 char natmedc[21], namatec[21];
3879 char key[5], name[5], mother[5], konly[5];
3881 Int_t iadvol, iadtmd, iadrot, nwtot, iret;
3882 Int_t mlevel, numbr, natt, numed, nin, ndata;
3883 Int_t iname, ivo, ish, jvo, nvstak, ivstak;
3884 Int_t jdiv, ivin, in, jin, jvin, irot;
3885 Int_t jtm, imat, jma, flag=0, imatc;
3886 Float_t az, dens, radl, absl, a, step, x, y, z;
3887 Int_t npar, ndvmx, left;
3888 Float_t zc, densc, radlc, abslc, c0, tmaxfd;
3890 Int_t iomate[100], iotmed[100];
3891 Float_t par[100], att[20], ubuf[50];
3894 Int_t level, ndiv, iaxe;
3895 Int_t itmedc, nmatc, isvolc, ifieldc, nwbufc, isvol, nmat, ifield, nwbuf;
3896 Float_t fieldmc, tmaxfdc, stemaxc, deemaxc, epsilc, stminc, fieldm;
3897 Float_t tmaxf, stemax, deemax, epsil, stmin;
3898 const char *k10000="!\n%s\n!\n";
3899 //Open the input file
3901 for(i=0;i<end;i++) if(filnam[i]=='.') {
3905 filext=new char[end+5];
3906 filetme=new char[end+5];
3907 strncpy(filext,filnam,end);
3908 strncpy(filetme,filnam,end);
3910 // *** The output filnam name will be with extension '.euc'
3911 strcpy(&filext[end],".euc");
3912 strcpy(&filetme[end],".tme");
3913 lun=fopen(filext,"w");
3915 // *** Initialisation of the working space
3916 iadvol=fGcnum->nvolum;
3917 iadtmd=iadvol+fGcnum->nvolum;
3918 iadrot=iadtmd+fGcnum->ntmed;
3919 if(fGclink->jrotm) {
3920 fGcnum->nrotm=fZiq[fGclink->jrotm-2];
3924 nwtot=iadrot+fGcnum->nrotm;
3925 qws = new float[nwtot+1];
3926 for (i=0;i<nwtot+1;i++) qws[i]=0;
3929 if(nlevel==0) mlevel=20;
3931 // *** find the top volume and put it in the stak
3932 numbr = number>0 ? number : 1;
3933 Gfpara(topvol,numbr,1,npar,natt,par,att);
3935 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3940 // *** authorized shape ?
3941 strncpy((char *)&iname, topvol, 4);
3943 for(i=1; i<=fGcnum->nvolum; i++) if(fZiq[fGclink->jvolum+i]==iname) {
3947 jvo = fZlq[fGclink->jvolum-ivo];
3948 ish = Int_t (fZq[jvo+2]);
3950 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3957 iws[iadvol+ivo] = level;
3960 //*** flag all volumes and fill the stak
3964 // pick the next volume in stak
3966 ivo = TMath::Abs(iws[ivstak]);
3967 jvo = fZlq[fGclink->jvolum - ivo];
3969 // flag the tracking medium
3970 numed = Int_t (fZq[jvo + 4]);
3971 iws[iadtmd + numed] = 1;
3973 // get the daughters ...
3974 level = iws[iadvol+ivo];
3975 if (level < mlevel) {
3977 nin = Int_t (fZq[jvo + 3]);
3979 // from division ...
3981 jdiv = fZlq[jvo - 1];
3982 ivin = Int_t (fZq[jdiv + 2]);
3984 iws[nvstak] = -ivin;
3985 iws[iadvol+ivin] = level;
3987 // from position ...
3988 } else if (nin > 0) {
3989 for(in=1; in<=nin; in++) {
3990 jin = fZlq[jvo - in];
3991 ivin = Int_t (fZq[jin + 2 ]);
3992 jvin = fZlq[fGclink->jvolum - ivin];
3993 ish = Int_t (fZq[jvin + 2]);
3994 // authorized shape ?
3996 // not yet flagged ?
3997 if (iws[iadvol+ivin]==0) {
4000 iws[iadvol+ivin] = level;
4002 // flag the rotation matrix
4003 irot = Int_t ( fZq[jin + 4 ]);
4004 if (irot > 0) iws[iadrot+irot] = 1;
4010 // next volume in stak ?
4011 if (ivstak < nvstak) goto L10;
4013 // *** restore original material and media numbers
4014 // file euc_medi.dat is needed to compare materials and medias
4016 FILE* luncor=fopen("euc_medi.dat","r");
4019 for(itm=1; itm<=fGcnum->ntmed; itm++) {
4020 if (iws[iadtmd+itm] > 0) {
4021 jtm = fZlq[fGclink->jtmed-itm];
4022 strncpy(natmed,(char *)&fZiq[jtm+1],20);
4023 imat = Int_t (fZq[jtm+6]);
4024 jma = fZlq[fGclink->jmate-imat];
4026 printf(" *** GWEUCL *** material not defined for tracking medium %5i %s\n",itm,natmed);
4029 strncpy(namate,(char *)&fZiq[jma+1],20);
4032 //** find the material original number
4035 iret=fscanf(luncor,"%4s,%130s",key,card);
4036 if(iret<=0) goto L26;
4038 if(!strcmp(key,"MATE")) {
4039 sscanf(card,"%d %s %f %f %f %f %f %d",&imatc,namatec,&az,&zc,&densc,&radlc,&abslc,&nparc);
4040 Gfmate(imat,namate,a,z,dens,radl,absl,par,npar);
4041 if(!strcmp(namatec,namate)) {
4042 if(az==a && zc==z && densc==dens && radlc==radl
4043 && abslc==absl && nparc==nparc) {
4046 printf("*** GWEUCL *** material : %3d '%s' restored as %3d\n",imat,namate,imatc);
4048 printf("*** GWEUCL *** different definitions for material: %s\n",namate);
4052 if(strcmp(key,"END") && !flag) goto L23;
4054 printf("*** GWEUCL *** cannot restore original number for material: %s\n",namate);
4058 //*** restore original tracking medium number
4061 iret=fscanf(luncor,"%4s,%130s",key,card);
4062 if(iret<=0) goto L26;
4064 if (!strcmp(key,"TMED")) {
4065 sscanf(card,"%d %s %d %d %d %f %f %f %f %f %f %d\n",
4066 &itmedc,natmedc,&nmatc,&isvolc,&ifieldc,&fieldmc,
4067 &tmaxfdc,&stemaxc,&deemaxc,&epsilc,&stminc,&nwbufc);
4068 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxf,stemax,deemax,
4069 epsil,stmin,ubuf,&nwbuf);
4070 if(!strcmp(natmedc,natmed)) {
4071 if (iomate[nmat]==nmatc && nwbuf==nwbufc) {
4074 printf("*** GWEUCL *** medium : %3d '%20s' restored as %3d\n",
4077 printf("*** GWEUCL *** different definitions for tracking medium: %s\n",natmed);
4081 if(strcmp(key,"END") && !flag) goto L24;
4083 printf("cannot restore original number for medium : %s\n",natmed);
4091 L26: printf("*** GWEUCL *** cannot read the data file\n");
4093 L29: if(luncor) fclose (luncor);
4096 // *** write down the tracking medium definition
4098 strcpy(card,"! Tracking medium");
4099 fprintf(lun,k10000,card);
4101 for(itm=1;itm<=fGcnum->ntmed;itm++) {
4102 if (iws[iadtmd+itm]>0) {
4103 jtm = fZlq[fGclink->jtmed-itm];
4104 strncpy(natmed,(char *)&fZiq[jtm+1],20);
4106 imat = Int_t (fZq[jtm+6]);
4107 jma = fZlq[fGclink->jmate-imat];
4108 //* order media from one, if comparing with database failed
4110 iotmed[itm]=++imxtmed;
4111 iomate[imat]=++imxmate;
4116 printf(" *** GWEUCL *** material not defined for tracking medium %5d %s\n",
4119 strncpy(namate,(char *)&fZiq[jma+1],20);
4122 fprintf(lun,"TMED %3d '%20s' %3d '%20s'\n",iotmed[itm],natmed,iomate[imat],namate);
4126 //* *** write down the rotation matrix
4128 strcpy(card,"! Reperes");
4129 fprintf(lun,k10000,card);
4131 for(irm=1;irm<=fGcnum->nrotm;irm++) {
4132 if (iws[iadrot+irm]>0) {
4133 jrm = fZlq[fGclink->jrotm-irm];
4134 fprintf(lun,"ROTM %3d",irm);
4135 for(k=11;k<=16;k++) fprintf(lun," %8.3f",fZq[jrm+k]);
4140 //* *** write down the volume definition
4142 strcpy(card,"! Volumes");
4143 fprintf(lun,k10000,card);
4145 for(ivstak=1;ivstak<=nvstak;ivstak++) {
4148 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivo],4);
4150 jvo = fZlq[fGclink->jvolum-ivo];
4151 ish = Int_t (fZq[jvo+2]);
4152 nmed = Int_t (fZq[jvo+4]);
4153 npar = Int_t (fZq[jvo+5]);
4155 if (ivstak>1) for(i=0;i<npar;i++) par[i]=fZq[jvo+7+i];
4156 Gckpar (ish,npar,par);
4157 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,kShape[ish-1],iotmed[nmed],npar);
4158 for(i=0;i<(npar-1)/6+1;i++) {
4161 for(k=0;k<(left<6?left:6);k++) fprintf(lun," %11.5f",par[i*6+k]);
4165 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,kShape[ish-1],iotmed[nmed],npar);
4170 //* *** write down the division of volumes
4172 fprintf(lun,k10000,"! Divisions");
4173 for(ivstak=1;ivstak<=nvstak;ivstak++) {
4174 ivo = TMath::Abs(iws[ivstak]);
4175 jvo = fZlq[fGclink->jvolum-ivo];
4176 ish = Int_t (fZq[jvo+2]);
4177 nin = Int_t (fZq[jvo+3]);
4178 //* this volume is divided ...
4181 iaxe = Int_t ( fZq[jdiv+1]);
4182 ivin = Int_t ( fZq[jdiv+2]);
4183 ndiv = Int_t ( fZq[jdiv+3]);
4186 jvin = fZlq[fGclink->jvolum-ivin];
4187 nmed = Int_t ( fZq[jvin+4]);
4188 strncpy(mother,(char *)&fZiq[fGclink->jvolum+ivo ],4);
4190 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivin],4);
4192 if ((step<=0.)||(ish>=11)) {
4193 //* volume with negative parameter or gsposp or pgon ...
4194 fprintf(lun,"DIVN '%4s' '%4s' %3d %3d\n",name,mother,ndiv,iaxe);
4195 } else if ((ndiv<=0)||(ish==10)) {
4196 //* volume with negative parameter or gsposp or para ...
4197 ndvmx = TMath::Abs(ndiv);
4198 fprintf(lun,"DIVT '%4s' '%4s' %11.5f %3d %3d %3d\n",
4199 name,mother,step,iaxe,iotmed[nmed],ndvmx);
4201 //* normal volume : all kind of division are equivalent
4202 fprintf(lun,"DVT2 '%4s' '%4s' %11.5f %3d %11.5f %3d %3d\n",
4203 name,mother,step,iaxe,c0,iotmed[nmed],ndiv);
4208 //* *** write down the the positionnement of volumes
4210 fprintf(lun,k10000,"! Positionnements\n");
4212 for(ivstak = 1;ivstak<=nvstak;ivstak++) {
4213 ivo = TMath::Abs(iws[ivstak]);
4214 strncpy(mother,(char*)&fZiq[fGclink->jvolum+ivo ],4);
4216 jvo = fZlq[fGclink->jvolum-ivo];
4217 nin = Int_t( fZq[jvo+3]);
4218 //* this volume has daughters ...
4220 for (in=1;in<=nin;in++) {
4222 ivin = Int_t (fZq[jin +2]);
4223 numb = Int_t (fZq[jin +3]);
4224 irot = Int_t (fZq[jin +4]);
4228 strcpy(konly,"ONLY");
4229 if (fZq[jin+8]!=1.) strcpy(konly,"MANY");
4230 strncpy(name,(char*)&fZiq[fGclink->jvolum+ivin],4);
4232 jvin = fZlq[fGclink->jvolum-ivin];
4233 ish = Int_t (fZq[jvin+2]);
4234 //* gspos or gsposp ?
4235 ndata = fZiq[jin-1];
4237 fprintf(lun,"POSI '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s'\n",
4238 name,numb,mother,x,y,z,irot,konly);
4240 npar = Int_t (fZq[jin+9]);
4241 for(i=0;i<npar;i++) par[i]=fZq[jin+10+i];
4242 Gckpar (ish,npar,par);
4243 fprintf(lun,"POSP '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s' %3d\n",
4244 name,numb,mother,x,y,z,irot,konly,npar);
4246 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
4253 fprintf(lun,"END\n");
4256 //****** write down the materials and medias *****
4258 lun=fopen(filetme,"w");
4260 for(itm=1;itm<=fGcnum->ntmed;itm++) {
4261 if (iws[iadtmd+itm]>0) {
4262 jtm = fZlq[fGclink->jtmed-itm];
4263 strncpy(natmed,(char*)&fZiq[jtm+1],4);
4264 imat = Int_t (fZq[jtm+6]);
4265 jma = Int_t (fZlq[fGclink->jmate-imat]);
4267 Gfmate (imat,namate,a,z,dens,radl,absl,par,npar);
4268 fprintf(lun,"MATE %4d '%20s'%11.5E %11.5E %11.5E %11.5E %11.5E %3d\n",
4269 iomate[imat],namate,a,z,dens,radl,absl,npar);
4273 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
4277 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxfd,stemax,deemax,epsil,stmin,par,&npar);
4278 fprintf(lun,"TMED %4d '%20s' %3d %1d %3d %11.5f %11.5f %11.5f %11.5f %11.5f %11.5f %3d\n",
4279 iotmed[itm],natmed,iomate[nmat],isvol,ifield,
4280 fieldm,tmaxfd,stemax,deemax,epsil,stmin,npar);
4284 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
4290 fprintf(lun,"END\n");
4292 printf(" *** GWEUCL *** file: %s is now written out\n",filext);
4293 printf(" *** GWEUCL *** file: %s is now written out\n",filetme);
4302 //_____________________________________________________________________________
4303 void TGeant3::Streamer(TBuffer &R__b)
4306 // Stream an object of class TGeant3.
4308 if (R__b.IsReading()) {
4309 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
4310 AliMC::Streamer(R__b);
4313 R__b.ReadStaticArray(fPDGCode);
4315 R__b.WriteVersion(TGeant3::IsA());
4316 AliMC::Streamer(R__b);
4319 R__b.WriteArray(fPDGCode, fNPDGCodes);