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.49 2001/04/06 14:04:25 morsch
19 Anti e-neutrino added to g3 particle list.
21 Revision 1.48 2001/04/04 11:47:56 morsch
22 - muon and tau neutrinos added to g3 particle list (needed for D,B decays).
23 - some (up to now harmless) bugs in Gspart calls corrected.
25 Revision 1.47 2001/03/20 06:36:29 alibrary
26 100 parameters now allowed for geant shapes
28 Revision 1.46 2000/12/21 17:35:05 morsch
29 Last updates on the right version (1.44).
30 (1.45) does not compile.
32 Revision 1.45 2000/12/21 16:49:56 morsch
33 Adding particles to the PDG database delegated to AliPDG.
35 Revision 1.44 2000/12/20 09:46:51 alibrary
36 dlsym not supported on HP, reverting to gcomad
38 Revision 1.43 2000/12/20 08:39:39 fca
39 Support for Cerenkov and process list in Virtual MC
41 Revision 1.42 2000/12/19 08:37:48 alibrary
42 Using dlsym to retrieve address of commons
44 Revision 1.41 2000/12/18 11:33:50 alibrary
45 New call frequence histograms per module and volume
47 Revision 1.40 2000/12/06 10:06:58 morsch
48 Add all D and B baryons produced by HIJING to PDG DataBase.
50 Revision 1.39 2000/11/30 07:12:54 alibrary
51 Introducing new Rndm and QA classes
53 Revision 1.38 2000/10/30 15:19:06 morsch
54 Xi(b) (pdg code 5232) added to Pdg data base.
56 Revision 1.37 2000/10/02 21:28:16 fca
57 Removal of useless dependecies via forward declarations
59 Revision 1.36 2000/09/14 07:08:41 fca
60 Introducing glvolu in the interface
62 Revision 1.35 2000/09/12 14:27:10 morsch
63 No instance of AliDecayer created to initialize fDecayer.
65 Revision 1.34 2000/09/07 12:12:01 morsch
66 Comment inside comment removed.
68 Revision 1.33 2000/09/06 16:03:42 morsch
69 Set ExternalDecayer, Decayer and SetForceDecay methods added.
70 Gspart calls for charmed and bottom hadrons added.
71 Decay mode definitions for charmed and beauty hadrons have been taken out.
72 They will be handled by an external decayer.
74 Revision 1.32 2000/08/24 16:28:53 hristov
75 TGeant3::IsNewTrack corrected by F.Carminati
77 Revision 1.31 2000/07/13 16:19:10 fca
78 Mainly coding conventions + some small bug fixes
80 Revision 1.30 2000/07/12 08:56:30 fca
81 Coding convention correction and warning removal
83 Revision 1.29 2000/07/11 18:24:59 fca
84 Coding convention corrections + few minor bug fixes
86 Revision 1.28 2000/06/29 10:51:55 morsch
87 Add some charmed and bottom baryons to the particle list (TDatabasePDG). This
88 is needed by Hijing. Should be part of a future review of TDatabasePDG.
90 Revision 1.27 2000/06/21 17:40:15 fca
91 Adding possibility to set ISTRA, PAI model
93 Revision 1.26 2000/05/16 13:10:41 fca
94 New method IsNewTrack and fix for a problem in Father-Daughter relations
96 Revision 1.25 2000/04/07 11:12:35 fca
97 G4 compatibility changes
99 Revision 1.24 2000/02/28 21:03:57 fca
100 Some additions to improve the compatibility with G4
102 Revision 1.23 2000/02/23 16:25:25 fca
103 AliVMC and AliGeant3 classes introduced
104 ReadEuclid moved from AliRun to AliModule
106 Revision 1.22 2000/01/18 15:40:13 morsch
107 Interface to GEANT3 routines GFTMAT, GBRELM and GPRELM added
108 Define geant particle type 51: Feedback Photon with Cherenkov photon properties.
110 Revision 1.21 2000/01/17 19:41:17 fca
113 Revision 1.20 2000/01/12 11:29:27 fca
116 Revision 1.19 1999/12/17 09:03:12 fca
117 Introduce a names array
119 Revision 1.18 1999/11/26 16:55:39 fca
120 Reimplement CurrentVolName() to avoid memory leaks
122 Revision 1.17 1999/11/03 16:58:28 fca
123 Correct source of address violation in creating character string
125 Revision 1.16 1999/11/03 13:17:08 fca
126 Have ProdProcess return const char*
128 Revision 1.15 1999/10/26 06:04:50 fca
129 Introduce TLorentzVector in AliMC::GetSecondary. Thanks to I.Hrivnacova
131 Revision 1.14 1999/09/29 09:24:30 fca
132 Introduction of the Copyright and cvs Log
136 ///////////////////////////////////////////////////////////////////////////////
138 // Interface Class to the Geant3.21 MonteCarlo //
142 <img src="picts/TGeant3Class.gif">
147 ///////////////////////////////////////////////////////////////////////////////
152 #include "TDatabasePDG.h"
153 #include "TLorentzVector.h"
159 #include "AliCallf77.h"
160 #include "AliDecayer.h"
164 # define gzebra gzebra_
165 # define grfile grfile_
166 # define gpcxyz gpcxyz_
167 # define ggclos ggclos_
168 # define glast glast_
169 # define ginit ginit_
170 # define gcinit gcinit_
172 # define gtrig gtrig_
173 # define gtrigc gtrigc_
174 # define gtrigi gtrigi_
175 # define gwork gwork_
176 # define gzinit gzinit_
177 # define gfmate gfmate_
178 # define gfpart gfpart_
179 # define gftmed gftmed_
180 # define gftmat gftmat_
181 # define gmate gmate_
182 # define gpart gpart_
184 # define gsmate gsmate_
185 # define gsmixt gsmixt_
186 # define gspart gspart_
187 # define gstmed gstmed_
188 # define gsckov gsckov_
189 # define gstpar gstpar_
190 # define gfkine gfkine_
191 # define gfvert gfvert_
192 # define gskine gskine_
193 # define gsvert gsvert_
194 # define gphysi gphysi_
195 # define gdebug gdebug_
196 # define gekbin gekbin_
197 # define gfinds gfinds_
198 # define gsking gsking_
199 # define gskpho gskpho_
200 # define gsstak gsstak_
201 # define gsxyz gsxyz_
202 # define gtrack gtrack_
203 # define gtreve gtreve_
204 # define gtreveroot gtreveroot_
205 # define grndm grndm_
206 # define grndmq grndmq_
207 # define gdtom gdtom_
208 # define glmoth glmoth_
209 # define gmedia gmedia_
210 # define gmtod gmtod_
211 # define gsdvn gsdvn_
212 # define gsdvn2 gsdvn2_
213 # define gsdvs gsdvs_
214 # define gsdvs2 gsdvs2_
215 # define gsdvt gsdvt_
216 # define gsdvt2 gsdvt2_
217 # define gsord gsord_
218 # define gspos gspos_
219 # define gsposp gsposp_
220 # define gsrotm gsrotm_
221 # define gprotm gprotm_
222 # define gsvolu gsvolu_
223 # define gprint gprint_
224 # define gdinit gdinit_
225 # define gdopt gdopt_
226 # define gdraw gdraw_
227 # define gdrayt gdrayt_
228 # define gdrawc gdrawc_
229 # define gdrawx gdrawx_
230 # define gdhead gdhead_
231 # define gdwmn1 gdwmn1_
232 # define gdwmn2 gdwmn2_
233 # define gdwmn3 gdwmn3_
234 # define gdxyz gdxyz_
235 # define gdcxyz gdcxyz_
236 # define gdman gdman_
237 # define gdspec gdspec_
238 # define gdtree gdtree_
239 # define gdelet gdelet_
240 # define gdclos gdclos_
241 # define gdshow gdshow_
242 # define gdopen gdopen_
243 # define dzshow dzshow_
244 # define gsatt gsatt_
245 # define gfpara gfpara_
246 # define gckpar gckpar_
247 # define gckmat gckmat_
248 # define glvolu glvolu_
249 # define geditv geditv_
250 # define mzdrop mzdrop_
252 # define ertrak ertrak_
253 # define ertrgo ertrgo_
255 # define setbomb setbomb_
256 # define setclip setclip_
257 # define gcomad gcomad_
259 # define gbrelm gbrelm_
260 # define gprelm gprelm_
262 # define gzebra GZEBRA
263 # define grfile GRFILE
264 # define gpcxyz GPCXYZ
265 # define ggclos GGCLOS
268 # define gcinit GCINIT
271 # define gtrigc GTRIGC
272 # define gtrigi GTRIGI
274 # define gzinit GZINIT
275 # define gfmate GFMATE
276 # define gfpart GFPART
277 # define gftmed GFTMED
278 # define gftmat GFTMAT
282 # define gsmate GSMATE
283 # define gsmixt GSMIXT
284 # define gspart GSPART
285 # define gstmed GSTMED
286 # define gsckov GSCKOV
287 # define gstpar GSTPAR
288 # define gfkine GFKINE
289 # define gfvert GFVERT
290 # define gskine GSKINE
291 # define gsvert GSVERT
292 # define gphysi GPHYSI
293 # define gdebug GDEBUG
294 # define gekbin GEKBIN
295 # define gfinds GFINDS
296 # define gsking GSKING
297 # define gskpho GSKPHO
298 # define gsstak GSSTAK
300 # define gtrack GTRACK
301 # define gtreve GTREVE
302 # define gtreveroot GTREVEROOT
304 # define grndmq GRNDMQ
306 # define glmoth GLMOTH
307 # define gmedia GMEDIA
310 # define gsdvn2 GSDVN2
312 # define gsdvs2 GSDVS2
314 # define gsdvt2 GSDVT2
317 # define gsposp GSPOSP
318 # define gsrotm GSROTM
319 # define gprotm GPROTM
320 # define gsvolu GSVOLU
321 # define gprint GPRINT
322 # define gdinit GDINIT
325 # define gdrayt GDRAYT
326 # define gdrawc GDRAWC
327 # define gdrawx GDRAWX
328 # define gdhead GDHEAD
329 # define gdwmn1 GDWMN1
330 # define gdwmn2 GDWMN2
331 # define gdwmn3 GDWMN3
333 # define gdcxyz GDCXYZ
335 # define gdfspc GDFSPC
336 # define gdspec GDSPEC
337 # define gdtree GDTREE
338 # define gdelet GDELET
339 # define gdclos GDCLOS
340 # define gdshow GDSHOW
341 # define gdopen GDOPEN
342 # define dzshow DZSHOW
344 # define gfpara GFPARA
345 # define gckpar GCKPAR
346 # define gckmat GCKMAT
347 # define glvolu GLVOLU
348 # define geditv GEDITV
349 # define mzdrop MZDROP
351 # define ertrak ERTRAK
352 # define ertrgo ERTRGO
354 # define setbomb SETBOMB
355 # define setclip SETCLIP
356 # define gcomad GCOMAD
358 # define gbrelm GBRELM
359 # define gprelm GPRELM
363 //____________________________________________________________________________
367 // Prototypes for GEANT functions
369 void type_of_call gzebra(const int&);
371 void type_of_call gpcxyz();
373 void type_of_call ggclos();
375 void type_of_call glast();
377 void type_of_call ginit();
379 void type_of_call gcinit();
381 void type_of_call grun();
383 void type_of_call gtrig();
385 void type_of_call gtrigc();
387 void type_of_call gtrigi();
389 void type_of_call gwork(const int&);
391 void type_of_call gzinit();
393 void type_of_call gmate();
395 void type_of_call gpart();
397 void type_of_call gsdk(Int_t &, Float_t *, Int_t *);
399 void type_of_call gfkine(Int_t &, Float_t *, Float_t *, Int_t &,
400 Int_t &, Float_t *, Int_t &);
402 void type_of_call gfvert(Int_t &, Float_t *, Int_t &, Int_t &,
403 Float_t &, Float_t *, Int_t &);
405 void type_of_call gskine(Float_t *,Int_t &, Int_t &, Float_t *,
408 void type_of_call gsvert(Float_t *,Int_t &, Int_t &, Float_t *,
411 void type_of_call gphysi();
413 void type_of_call gdebug();
415 void type_of_call gekbin();
417 void type_of_call gfinds();
419 void type_of_call gsking(Int_t &);
421 void type_of_call gskpho(Int_t &);
423 void type_of_call gsstak(Int_t &);
425 void type_of_call gsxyz();
427 void type_of_call gtrack();
429 void type_of_call gtreve();
431 void type_of_call gtreveroot();
433 void type_of_call grndm(Float_t *r, const Int_t &n)
436 void type_of_call grndmq(Int_t &, Int_t &, const Int_t &,
438 {/*printf("Dummy grndmq called\n");*/}
440 void type_of_call gdtom(Float_t *, Float_t *, Int_t &);
442 void type_of_call glmoth(DEFCHARD, Int_t &, Int_t &, Int_t *,
443 Int_t *, Int_t * DEFCHARL);
445 void type_of_call gmedia(Float_t *, Int_t &);
447 void type_of_call gmtod(Float_t *, Float_t *, Int_t &);
449 void type_of_call gsrotm(const Int_t &, const Float_t &, const Float_t &,
450 const Float_t &, const Float_t &, const Float_t &,
453 void type_of_call gprotm(const Int_t &);
455 void type_of_call grfile(const Int_t&, DEFCHARD,
456 DEFCHARD DEFCHARL DEFCHARL);
458 void type_of_call gfmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
459 Float_t &, Float_t &, Float_t &, Float_t *,
462 void type_of_call gfpart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
463 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
465 void type_of_call gftmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
466 Float_t &, Float_t &, Float_t &, Float_t &,
467 Float_t &, Float_t &, Float_t *, Int_t * DEFCHARL);
469 void type_of_call gftmat(const Int_t&, const Int_t&, DEFCHARD, const Int_t&,
471 ,Float_t *, Int_t & DEFCHARL);
473 void type_of_call gsmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
474 Float_t &, Float_t &, Float_t &, Float_t *,
477 void type_of_call gsmixt(const Int_t&, DEFCHARD, Float_t *, Float_t *,
478 Float_t &, Int_t &, Float_t * DEFCHARL);
480 void type_of_call gspart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
481 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
484 void type_of_call gstmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
485 Float_t &, Float_t &, Float_t &, Float_t &,
486 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
488 void type_of_call gsckov(Int_t &itmed, Int_t &npckov, Float_t *ppckov,
489 Float_t *absco, Float_t *effic, Float_t *rindex);
490 void type_of_call gstpar(const Int_t&, DEFCHARD, Float_t & DEFCHARL);
492 void type_of_call gsdvn(DEFCHARD,DEFCHARD, Int_t &, Int_t &
495 void type_of_call gsdvn2(DEFCHARD,DEFCHARD, Int_t &, Int_t &, Float_t &,
496 Int_t & DEFCHARL DEFCHARL);
498 void type_of_call gsdvs(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &
501 void type_of_call gsdvs2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t &,
502 Int_t & DEFCHARL DEFCHARL);
504 void type_of_call gsdvt(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &,
505 Int_t & DEFCHARL DEFCHARL);
507 void type_of_call gsdvt2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t&,
508 Int_t &, Int_t & DEFCHARL DEFCHARL);
510 void type_of_call gsord(DEFCHARD, Int_t & DEFCHARL);
512 void type_of_call gspos(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
513 Float_t &, Int_t &, DEFCHARD DEFCHARL DEFCHARL
516 void type_of_call gsposp(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
517 Float_t &, Int_t &, DEFCHARD,
518 Float_t *, Int_t & DEFCHARL DEFCHARL DEFCHARL);
520 void type_of_call gsvolu(DEFCHARD, DEFCHARD, Int_t &, Float_t *, Int_t &,
521 Int_t & DEFCHARL DEFCHARL);
523 void type_of_call gsatt(DEFCHARD, DEFCHARD, Int_t & DEFCHARL DEFCHARL);
525 void type_of_call gfpara(DEFCHARD , Int_t&, Int_t&, Int_t&, Int_t&, Float_t*,
528 void type_of_call gckpar(Int_t&, Int_t&, Float_t*);
530 void type_of_call gckmat(Int_t&, DEFCHARD DEFCHARL);
532 void type_of_call glvolu(Int_t&, Int_t*, Int_t*, Int_t&);
534 void type_of_call gprint(DEFCHARD,const int& DEFCHARL);
536 void type_of_call gdinit();
538 void type_of_call gdopt(DEFCHARD,DEFCHARD DEFCHARL DEFCHARL);
540 void type_of_call gdraw(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
541 Float_t &, Float_t &, Float_t & DEFCHARL);
542 void type_of_call gdrayt(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
543 Float_t &, Float_t &, Float_t & DEFCHARL);
544 void type_of_call gdrawc(DEFCHARD,Int_t &, Float_t &, Float_t &, Float_t &,
545 Float_t &, Float_t & DEFCHARL);
546 void type_of_call gdrawx(DEFCHARD,Float_t &, Float_t &, Float_t &, Float_t &,
547 Float_t &, Float_t &, Float_t &, Float_t &,
549 void type_of_call gdhead(Int_t &,DEFCHARD, Float_t & DEFCHARL);
550 void type_of_call gdxyz(Int_t &);
551 void type_of_call gdcxyz();
552 void type_of_call gdman(Float_t &, Float_t &);
553 void type_of_call gdwmn1(Float_t &, Float_t &);
554 void type_of_call gdwmn2(Float_t &, Float_t &);
555 void type_of_call gdwmn3(Float_t &, Float_t &);
556 void type_of_call gdspec(DEFCHARD DEFCHARL);
557 void type_of_call gdfspc(DEFCHARD, Int_t &, Int_t & DEFCHARL) {;}
558 void type_of_call gdtree(DEFCHARD, Int_t &, Int_t & DEFCHARL);
560 void type_of_call gdopen(Int_t &);
561 void type_of_call gdclos();
562 void type_of_call gdelet(Int_t &);
563 void type_of_call gdshow(Int_t &);
564 void type_of_call geditv(Int_t &) {;}
567 void type_of_call dzshow(DEFCHARD,const int&,const int&,DEFCHARD,const int&,
568 const int&, const int&, const int& DEFCHARL
571 void type_of_call mzdrop(Int_t&, Int_t&, DEFCHARD DEFCHARL);
573 void type_of_call setbomb(Float_t &);
574 void type_of_call setclip(DEFCHARD, Float_t &,Float_t &,Float_t &,Float_t &,
575 Float_t &, Float_t & DEFCHARL);
576 void type_of_call gcomad(DEFCHARD, Int_t*& DEFCHARL);
578 void type_of_call ertrak(const Float_t *const x1, const Float_t *const p1,
579 const Float_t *x2, const Float_t *p2,
580 const Int_t &ipa, DEFCHARD DEFCHARL);
582 void type_of_call ertrgo();
584 float type_of_call gbrelm(const Float_t &z, const Float_t& t, const Float_t& cut);
585 float type_of_call gprelm(const Float_t &z, const Float_t& t, const Float_t& cut);
589 // Geant3 global pointer
591 static const Int_t kDefSize = 600;
595 //____________________________________________________________________________
599 // Default constructor
603 //____________________________________________________________________________
604 TGeant3::TGeant3(const char *title, Int_t nwgeant)
605 :AliMC("TGeant3",title)
608 // Standard constructor for TGeant3 with ZEBRA initialisation
619 // Load Address of Geant3 commons
622 // Zero number of particles
627 //____________________________________________________________________________
628 Int_t TGeant3::CurrentMaterial(Float_t &a, Float_t &z, Float_t &dens,
629 Float_t &radl, Float_t &absl) const
632 // Return the parameters of the current material during transport
636 dens = fGcmate->dens;
637 radl = fGcmate->radl;
638 absl = fGcmate->absl;
639 return 1; //this could be the number of elements in mixture
642 //____________________________________________________________________________
643 void TGeant3::DefaultRange()
646 // Set range of current drawing pad to 20x20 cm
652 gHigz->Range(0,0,20,20);
655 //____________________________________________________________________________
656 void TGeant3::InitHIGZ()
667 //____________________________________________________________________________
668 void TGeant3::LoadAddress()
671 // Assigns the address of the GEANT common blocks to the structures
672 // that allow their access from C++
675 gcomad(PASSCHARD("QUEST"), (int*&) fQuest PASSCHARL("QUEST"));
676 gcomad(PASSCHARD("GCBANK"),(int*&) fGcbank PASSCHARL("GCBANK"));
677 gcomad(PASSCHARD("GCLINK"),(int*&) fGclink PASSCHARL("GCLINK"));
678 gcomad(PASSCHARD("GCCUTS"),(int*&) fGccuts PASSCHARL("GCCUTS"));
679 gcomad(PASSCHARD("GCMULO"),(int*&) fGcmulo PASSCHARL("GCMULO"));
680 gcomad(PASSCHARD("GCFLAG"),(int*&) fGcflag PASSCHARL("GCFLAG"));
681 gcomad(PASSCHARD("GCKINE"),(int*&) fGckine PASSCHARL("GCKINE"));
682 gcomad(PASSCHARD("GCKING"),(int*&) fGcking PASSCHARL("GCKING"));
683 gcomad(PASSCHARD("GCKIN2"),(int*&) fGckin2 PASSCHARL("GCKIN2"));
684 gcomad(PASSCHARD("GCKIN3"),(int*&) fGckin3 PASSCHARL("GCKIN3"));
685 gcomad(PASSCHARD("GCMATE"),(int*&) fGcmate PASSCHARL("GCMATE"));
686 gcomad(PASSCHARD("GCTMED"),(int*&) fGctmed PASSCHARL("GCTMED"));
687 gcomad(PASSCHARD("GCTRAK"),(int*&) fGctrak PASSCHARL("GCTRAK"));
688 gcomad(PASSCHARD("GCTPOL"),(int*&) fGctpol PASSCHARL("GCTPOL"));
689 gcomad(PASSCHARD("GCVOLU"),(int*&) fGcvolu PASSCHARL("GCVOLU"));
690 gcomad(PASSCHARD("GCNUM"), (int*&) fGcnum PASSCHARL("GCNUM"));
691 gcomad(PASSCHARD("GCSETS"),(int*&) fGcsets PASSCHARL("GCSETS"));
692 gcomad(PASSCHARD("GCPHYS"),(int*&) fGcphys PASSCHARL("GCPHYS"));
693 gcomad(PASSCHARD("GCPHLT"),(int*&) fGcphlt PASSCHARL("GCPHLT"));
694 gcomad(PASSCHARD("GCOPTI"),(int*&) fGcopti PASSCHARL("GCOPTI"));
695 gcomad(PASSCHARD("GCTLIT"),(int*&) fGctlit PASSCHARL("GCTLIT"));
696 gcomad(PASSCHARD("GCVDMA"),(int*&) fGcvdma PASSCHARL("GCVDMA"));
699 gcomad(PASSCHARD("ERTRIO"),(int*&) fErtrio PASSCHARL("ERTRIO"));
700 gcomad(PASSCHARD("EROPTS"),(int*&) fEropts PASSCHARL("EROPTS"));
701 gcomad(PASSCHARD("EROPTC"),(int*&) fEroptc PASSCHARL("EROPTC"));
702 gcomad(PASSCHARD("ERWORK"),(int*&) fErwork PASSCHARL("ERWORK"));
704 // Variables for ZEBRA store
705 gcomad(PASSCHARD("IQ"), addr PASSCHARL("IQ"));
707 gcomad(PASSCHARD("LQ"), addr PASSCHARL("LQ"));
712 //_____________________________________________________________________________
713 void TGeant3::GeomIter()
716 // Geometry iterator for moving upward in the geometry tree
717 // Initialise the iterator
719 fNextVol=fGcvolu->nlevel;
722 //____________________________________________________________________________
723 void TGeant3::FinishGeometry()
725 //Close the geometry structure
729 //____________________________________________________________________________
730 Int_t TGeant3::NextVolUp(Text_t *name, Int_t ©)
733 // Geometry iterator for moving upward in the geometry tree
734 // Return next volume up
739 gname=fGcvolu->names[fNextVol];
740 copy=fGcvolu->number[fNextVol];
741 i=fGcvolu->lvolum[fNextVol];
742 name = fVolNames[i-1];
743 if(gname == fZiq[fGclink->jvolum+i]) return i;
744 else printf("GeomTree: Volume %s not found in bank\n",name);
749 //_____________________________________________________________________________
750 void TGeant3::BuildPhysics()
755 //_____________________________________________________________________________
756 Int_t TGeant3::CurrentVolID(Int_t ©) const
759 // Returns the current volume ID and copy number
762 if( (i=fGcvolu->nlevel-1) < 0 ) {
763 Warning("CurrentVolID","Stack depth only %d\n",fGcvolu->nlevel);
765 gname=fGcvolu->names[i];
766 copy=fGcvolu->number[i];
767 i=fGcvolu->lvolum[i];
768 if(gname == fZiq[fGclink->jvolum+i]) return i;
769 else Warning("CurrentVolID","Volume %4s not found\n",(char*)&gname);
774 //_____________________________________________________________________________
775 Int_t TGeant3::CurrentVolOffID(Int_t off, Int_t ©) const
778 // Return the current volume "off" upward in the geometrical tree
779 // ID and copy number
782 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
783 Warning("CurrentVolOffID","Offset requested %d but stack depth %d\n",
784 off,fGcvolu->nlevel);
786 gname=fGcvolu->names[i];
787 copy=fGcvolu->number[i];
788 i=fGcvolu->lvolum[i];
789 if(gname == fZiq[fGclink->jvolum+i]) return i;
790 else Warning("CurrentVolOffID","Volume %4s not found\n",(char*)&gname);
795 //_____________________________________________________________________________
796 const char* TGeant3::CurrentVolName() const
799 // Returns the current volume name
802 if( (i=fGcvolu->nlevel-1) < 0 ) {
803 Warning("CurrentVolName","Stack depth %d\n",fGcvolu->nlevel);
805 gname=fGcvolu->names[i];
806 i=fGcvolu->lvolum[i];
807 if(gname == fZiq[fGclink->jvolum+i]) return fVolNames[i-1];
808 else Warning("CurrentVolName","Volume %4s not found\n",(char*) &gname);
813 //_____________________________________________________________________________
814 const char* TGeant3::CurrentVolOffName(Int_t off) const
817 // Return the current volume "off" upward in the geometrical tree
818 // ID, name and copy number
819 // if name=0 no name is returned
822 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
823 Warning("CurrentVolOffName",
824 "Offset requested %d but stack depth %d\n",off,fGcvolu->nlevel);
826 gname=fGcvolu->names[i];
827 i=fGcvolu->lvolum[i];
828 if(gname == fZiq[fGclink->jvolum+i]) return fVolNames[i-1];
829 else Warning("CurrentVolOffName","Volume %4s not found\n",(char*)&gname);
834 //_____________________________________________________________________________
835 Int_t TGeant3::IdFromPDG(Int_t pdg) const
838 // Return Geant3 code from PDG and pseudo ENDF code
840 for(Int_t i=0;i<fNPDGCodes;++i)
841 if(pdg==fPDGCode[i]) return i;
845 //_____________________________________________________________________________
846 Int_t TGeant3::PDGFromId(Int_t id) const
849 // Return PDG code and pseudo ENDF code from Geant3 code
851 if(id>0 && id<fNPDGCodes) return fPDGCode[id];
855 //_____________________________________________________________________________
856 void TGeant3::DefineParticles()
859 // Define standard Geant 3 particles
862 // Load standard numbers for GEANT particles and PDG conversion
863 fPDGCode[fNPDGCodes++]=-99; // 0 = unused location
864 fPDGCode[fNPDGCodes++]=22; // 1 = photon
865 fPDGCode[fNPDGCodes++]=-11; // 2 = positron
866 fPDGCode[fNPDGCodes++]=11; // 3 = electron
867 fPDGCode[fNPDGCodes++]=12; // 4 = neutrino e
868 fPDGCode[fNPDGCodes++]=-13; // 5 = muon +
869 fPDGCode[fNPDGCodes++]=13; // 6 = muon -
870 fPDGCode[fNPDGCodes++]=111; // 7 = pi0
871 fPDGCode[fNPDGCodes++]=211; // 8 = pi+
872 fPDGCode[fNPDGCodes++]=-211; // 9 = pi-
873 fPDGCode[fNPDGCodes++]=130; // 10 = Kaon Long
874 fPDGCode[fNPDGCodes++]=321; // 11 = Kaon +
875 fPDGCode[fNPDGCodes++]=-321; // 12 = Kaon -
876 fPDGCode[fNPDGCodes++]=2112; // 13 = Neutron
877 fPDGCode[fNPDGCodes++]=2212; // 14 = Proton
878 fPDGCode[fNPDGCodes++]=-2212; // 15 = Anti Proton
879 fPDGCode[fNPDGCodes++]=310; // 16 = Kaon Short
880 fPDGCode[fNPDGCodes++]=221; // 17 = Eta
881 fPDGCode[fNPDGCodes++]=3122; // 18 = Lambda
882 fPDGCode[fNPDGCodes++]=3222; // 19 = Sigma +
883 fPDGCode[fNPDGCodes++]=3212; // 20 = Sigma 0
884 fPDGCode[fNPDGCodes++]=3112; // 21 = Sigma -
885 fPDGCode[fNPDGCodes++]=3322; // 22 = Xi0
886 fPDGCode[fNPDGCodes++]=3312; // 23 = Xi-
887 fPDGCode[fNPDGCodes++]=3334; // 24 = Omega-
888 fPDGCode[fNPDGCodes++]=-2112; // 25 = Anti Proton
889 fPDGCode[fNPDGCodes++]=-3122; // 26 = Anti Proton
890 fPDGCode[fNPDGCodes++]=-3222; // 27 = Anti Sigma -
891 fPDGCode[fNPDGCodes++]=-3212; // 28 = Anti Sigma 0
892 fPDGCode[fNPDGCodes++]=-3112; // 29 = Anti Sigma 0
893 fPDGCode[fNPDGCodes++]=-3322; // 30 = Anti Xi 0
894 fPDGCode[fNPDGCodes++]=-3312; // 31 = Anti Xi +
895 fPDGCode[fNPDGCodes++]=-3334; // 32 = Anti Omega +
902 /* --- Define additional particles */
903 Gspart(33, "OMEGA(782)", 3, 0.782, 0., 7.836e-23);
904 fPDGCode[fNPDGCodes++]=223; // 33 = Omega(782)
906 Gspart(34, "PHI(1020)", 3, 1.019, 0., 1.486e-22);
907 fPDGCode[fNPDGCodes++]=333; // 34 = PHI (1020)
909 Gspart(35, "D +", 4, 1.87, 1., 1.066e-12);
910 fPDGCode[fNPDGCodes++]=411; // 35 = D+
912 Gspart(36, "D -", 4, 1.87, -1., 1.066e-12);
913 fPDGCode[fNPDGCodes++]=-411; // 36 = D-
915 Gspart(37, "D 0", 3, 1.865, 0., 4.2e-13);
916 fPDGCode[fNPDGCodes++]=421; // 37 = D0
918 Gspart(38, "ANTI D 0", 3, 1.865, 0., 4.2e-13);
919 fPDGCode[fNPDGCodes++]=-421; // 38 = D0 bar
922 fPDGCode[fNPDGCodes++]=-99; // 39 = unassigned
924 fPDGCode[fNPDGCodes++]=-99; // 40 = unassigned
926 fPDGCode[fNPDGCodes++]=-99; // 41 = unassigned
928 Gspart(42, "RHO +", 4, 0.768, 1., 4.353e-24);
929 fPDGCode[fNPDGCodes++]=213; // 42 = RHO+
931 Gspart(43, "RHO -", 4, 0.768, -1., 4.353e-24);
932 fPDGCode[fNPDGCodes++]=-213; // 43 = RHO-
934 Gspart(44, "RHO 0", 3, 0.768, 0., 4.353e-24);
935 fPDGCode[fNPDGCodes++]=113; // 44 = RHO0
938 // Use ENDF-6 mapping for ions = 10000*z+10*a+iso
940 // and numbers above 5 000 000 for special applications
943 const Int_t kion=10000000;
945 const Int_t kspe=50000000;
950 fPDGCode[fNPDGCodes++]=kion+10020; // 45 = Deuteron
952 fPDGCode[fNPDGCodes++]=kion+10030; // 46 = Triton
954 fPDGCode[fNPDGCodes++]=kion+20040; // 47 = Alpha
956 fPDGCode[fNPDGCodes++]=0; // 48 = geantino mapped to rootino
958 fPDGCode[fNPDGCodes++]=kion+20030; // 49 = HE3
960 fPDGCode[fNPDGCodes++]=kspe+50; // 50 = Cherenkov
962 Gspart(51, "FeedbackPhoton", 7, 0., 0.,1.e20 );
963 fPDGCode[fNPDGCodes++]=kspe+51; // 51 = FeedbackPhoton
965 Gspart(52, "Lambda_c+", 4, 2.2849, +1., 2.06e-13);
966 fPDGCode[fNPDGCodes++]=4122; //52 = Lambda_c+
968 Gspart(53, "Lambda_c-", 4, 2.2849, -1., 2.06e-13);
969 fPDGCode[fNPDGCodes++]=-4122; //53 = Lambda_c-
971 Gspart(54, "D_s+", 4, 1.9685, +1., 4.67e-13);
972 fPDGCode[fNPDGCodes++]=431; //54 = D_s+
974 Gspart(55, "D_s-", 4, 1.9685, -1., 4.67e-13);
975 fPDGCode[fNPDGCodes++]=-431; //55 = D_s-
977 Gspart(56, "Tau+", 5, 1.77705, +1., 2.9e-13);
978 fPDGCode[fNPDGCodes++]=15; //56 = Tau+
980 Gspart(57, "Tau-", 5, 1.77705, -1., 2.9e-13);
981 fPDGCode[fNPDGCodes++]=-15; //57 = Tau-
983 Gspart(58, "B0", 3, 5.2792, +0., 1.56e-12);
984 fPDGCode[fNPDGCodes++]=511; //58 = B0
986 Gspart(59, "B0 bar", 3, 5.2792, -0., 1.56e-12);
987 fPDGCode[fNPDGCodes++]=-511; //58 = B0bar
989 Gspart(60, "B+", 4, 5.2789, +1., 1.65e-12);
990 fPDGCode[fNPDGCodes++]=521; //60 = B+
992 Gspart(61, "B-", 4, 5.2789, -1., 1.65e-12);
993 fPDGCode[fNPDGCodes++]=-521; //61 = B-
995 Gspart(62, "Bs", 3, 5.3693, +0., 1.54e-12);
996 fPDGCode[fNPDGCodes++]=531; //62 = B_s
998 Gspart(63, "Bs bar", 3, 5.3693, -0., 1.54e-12);
999 fPDGCode[fNPDGCodes++]=-531; //63 = B_s bar
1001 Gspart(64, "Lambda_b", 3, 5.624, +0., 1.24e-12);
1002 fPDGCode[fNPDGCodes++]=5122; //64 = Lambda_b
1004 Gspart(65, "Lambda_b bar", 3, 5.624, -0., 1.24e-12);
1005 fPDGCode[fNPDGCodes++]=-5122; //65 = Lambda_b bar
1007 Gspart(66, "J/Psi", 3, 3.09688, 0., 0.);
1008 fPDGCode[fNPDGCodes++]=443; // 66 = J/Psi
1010 Gspart(67, "Psi Prime", 3, 3.686, 0., 0.);
1011 fPDGCode[fNPDGCodes++]=20443; // 67 = Psi prime
1013 Gspart(68, "Upsilon(1S)", 3, 9.46037, 0., 0.);
1014 fPDGCode[fNPDGCodes++]=553; // 68 = Upsilon(1S)
1016 Gspart(69, "Upsilon(2S)", 3, 10.0233, 0., 0.);
1017 fPDGCode[fNPDGCodes++]=20553; // 69 = Upsilon(2S)
1019 Gspart(70, "Upsilon(3S)", 3, 10.3553, 0., 0.);
1020 fPDGCode[fNPDGCodes++]=30553; // 70 = Upsilon(3S)
1022 Gspart(71, "Anti Neutrino (e)", 3, 0., 0., 1.e20);
1023 fPDGCode[fNPDGCodes++]=-12; // 71 = anti electron neutrino
1025 Gspart(72, "Neutrino (mu)", 3, 0., 0., 1.e20);
1026 fPDGCode[fNPDGCodes++]=14; // 72 = muon neutrino
1028 Gspart(73, "Anti Neutrino (mu)", 3, 0., 0., 1.e20);
1029 fPDGCode[fNPDGCodes++]=-14; // 73 = anti muon neutrino
1031 Gspart(74, "Neutrino (tau)", 3, 0., 0., 1.e20);
1032 fPDGCode[fNPDGCodes++]=16; // 74 = tau neutrino
1034 Gspart(75, "Anti Neutrino (tau)",3, 0., 0., 1.e20);
1035 fPDGCode[fNPDGCodes++]=-16; // 75 = anti tau neutrino
1037 /* --- Define additional decay modes --- */
1038 /* --- omega(783) --- */
1039 for (kz = 0; kz < 6; ++kz) {
1050 Gsdk(ipa, bratio, mode);
1051 /* --- phi(1020) --- */
1052 for (kz = 0; kz < 6; ++kz) {
1067 Gsdk(ipa, bratio, mode);
1070 for (kz = 0; kz < 6; ++kz) {
1083 Gsdk(ipa, bratio, mode);
1087 for (kz = 0; kz < 6; ++kz) {
1100 Gsdk(ipa, bratio, mode);
1104 for (kz = 0; kz < 6; ++kz) {
1115 Gsdk(ipa, bratio, mode);
1117 /* --- Anti D0 --- */
1119 for (kz = 0; kz < 6; ++kz) {
1130 Gsdk(ipa, bratio, mode);
1133 for (kz = 0; kz < 6; ++kz) {
1140 Gsdk(ipa, bratio, mode);
1142 for (kz = 0; kz < 6; ++kz) {
1149 Gsdk(ipa, bratio, mode);
1151 for (kz = 0; kz < 6; ++kz) {
1158 Gsdk(ipa, bratio, mode);
1161 for (kz = 0; kz < 6; ++kz) {
1170 Gsdk(ipa, bratio, mode);
1173 Gsdk(ipa, bratio, mode);
1176 Gsdk(ipa, bratio, mode);
1179 AliPDG::AddParticlesToPdgDataBase();
1182 //_____________________________________________________________________________
1183 Int_t TGeant3::VolId(const Text_t *name) const
1186 // Return the unique numeric identifier for volume name
1189 strncpy((char *) &gname, name, 4);
1190 for(i=1; i<=fGcnum->nvolum; i++)
1191 if(gname == fZiq[fGclink->jvolum+i]) return i;
1192 printf("VolId: Volume %s not found\n",name);
1196 //_____________________________________________________________________________
1197 Int_t TGeant3::NofVolumes() const
1200 // Return total number of volumes in the geometry
1202 return fGcnum->nvolum;
1205 //_____________________________________________________________________________
1206 Int_t TGeant3::VolId2Mate(Int_t id) const
1209 // Return material number for a given volume id
1211 if(id<1 || id > fGcnum->nvolum || fGclink->jvolum<=0)
1214 Int_t jvo = fZlq[fGclink->jvolum-id];
1215 return Int_t(fZq[jvo+4]);
1219 //_____________________________________________________________________________
1220 const char* TGeant3::VolName(Int_t id) const
1223 // Return the volume name given the volume identifier
1225 if(id<1 || id > fGcnum->nvolum || fGclink->jvolum<=0)
1226 return fVolNames[fGcnum->nvolum];
1228 return fVolNames[id-1];
1231 //_____________________________________________________________________________
1232 void TGeant3::SetCut(const char* cutName, Float_t cutValue)
1235 // Set transport cuts for particles
1237 if(!strcmp(cutName,"CUTGAM"))
1238 fGccuts->cutgam=cutValue;
1239 else if(!strcmp(cutName,"CUTELE"))
1240 fGccuts->cutele=cutValue;
1241 else if(!strcmp(cutName,"CUTNEU"))
1242 fGccuts->cutneu=cutValue;
1243 else if(!strcmp(cutName,"CUTHAD"))
1244 fGccuts->cuthad=cutValue;
1245 else if(!strcmp(cutName,"CUTMUO"))
1246 fGccuts->cutmuo=cutValue;
1247 else if(!strcmp(cutName,"BCUTE"))
1248 fGccuts->bcute=cutValue;
1249 else if(!strcmp(cutName,"BCUTM"))
1250 fGccuts->bcutm=cutValue;
1251 else if(!strcmp(cutName,"DCUTE"))
1252 fGccuts->dcute=cutValue;
1253 else if(!strcmp(cutName,"DCUTM"))
1254 fGccuts->dcutm=cutValue;
1255 else if(!strcmp(cutName,"PPCUTM"))
1256 fGccuts->ppcutm=cutValue;
1257 else if(!strcmp(cutName,"TOFMAX"))
1258 fGccuts->tofmax=cutValue;
1259 else Warning("SetCut","Cut %s not implemented\n",cutName);
1262 //_____________________________________________________________________________
1263 void TGeant3::SetProcess(const char* flagName, Int_t flagValue)
1266 // Set thresholds for different processes
1268 if(!strcmp(flagName,"PAIR"))
1269 fGcphys->ipair=flagValue;
1270 else if(!strcmp(flagName,"COMP"))
1271 fGcphys->icomp=flagValue;
1272 else if(!strcmp(flagName,"PHOT"))
1273 fGcphys->iphot=flagValue;
1274 else if(!strcmp(flagName,"PFIS"))
1275 fGcphys->ipfis=flagValue;
1276 else if(!strcmp(flagName,"DRAY"))
1277 fGcphys->idray=flagValue;
1278 else if(!strcmp(flagName,"ANNI"))
1279 fGcphys->ianni=flagValue;
1280 else if(!strcmp(flagName,"BREM"))
1281 fGcphys->ibrem=flagValue;
1282 else if(!strcmp(flagName,"HADR"))
1283 fGcphys->ihadr=flagValue;
1284 else if(!strcmp(flagName,"MUNU"))
1285 fGcphys->imunu=flagValue;
1286 else if(!strcmp(flagName,"DCAY"))
1287 fGcphys->idcay=flagValue;
1288 else if(!strcmp(flagName,"LOSS"))
1289 fGcphys->iloss=flagValue;
1290 else if(!strcmp(flagName,"MULS"))
1291 fGcphys->imuls=flagValue;
1292 else if(!strcmp(flagName,"RAYL"))
1293 fGcphys->irayl=flagValue;
1294 else if(!strcmp(flagName,"STRA"))
1295 fGcphlt->istra=flagValue;
1296 else if(!strcmp(flagName,"SYNC"))
1297 fGcphlt->isync=flagValue;
1298 else Warning("SetFlag","Flag %s not implemented\n",flagName);
1301 //_____________________________________________________________________________
1302 Float_t TGeant3::Xsec(char* reac, Float_t /* energy */,
1303 Int_t part, Int_t /* mate */)
1306 // Calculate X-sections -- dummy for the moment
1308 if(!strcmp(reac,"PHOT"))
1311 Error("Xsec","Can calculate photoelectric only for photons\n");
1317 //_____________________________________________________________________________
1318 void TGeant3::TrackPosition(TLorentzVector &xyz) const
1321 // Return the current position in the master reference frame of the
1322 // track being transported
1324 xyz[0]=fGctrak->vect[0];
1325 xyz[1]=fGctrak->vect[1];
1326 xyz[2]=fGctrak->vect[2];
1327 xyz[3]=fGctrak->tofg;
1330 //_____________________________________________________________________________
1331 Float_t TGeant3::TrackTime() const
1334 // Return the current time of flight of the track being transported
1336 return fGctrak->tofg;
1339 //_____________________________________________________________________________
1340 void TGeant3::TrackMomentum(TLorentzVector &xyz) const
1343 // Return the direction and the momentum (GeV/c) of the track
1344 // currently being transported
1346 Double_t ptot=fGctrak->vect[6];
1347 xyz[0]=fGctrak->vect[3]*ptot;
1348 xyz[1]=fGctrak->vect[4]*ptot;
1349 xyz[2]=fGctrak->vect[5]*ptot;
1350 xyz[3]=fGctrak->getot;
1353 //_____________________________________________________________________________
1354 Float_t TGeant3::TrackCharge() const
1357 // Return charge of the track currently transported
1359 return fGckine->charge;
1362 //_____________________________________________________________________________
1363 Float_t TGeant3::TrackMass() const
1366 // Return the mass of the track currently transported
1368 return fGckine->amass;
1371 //_____________________________________________________________________________
1372 Int_t TGeant3::TrackPid() const
1375 // Return the id of the particle transported
1377 return PDGFromId(fGckine->ipart);
1380 //_____________________________________________________________________________
1381 Float_t TGeant3::TrackStep() const
1384 // Return the length in centimeters of the current step
1386 return fGctrak->step;
1389 //_____________________________________________________________________________
1390 Float_t TGeant3::TrackLength() const
1393 // Return the length of the current track from its origin
1395 return fGctrak->sleng;
1398 //_____________________________________________________________________________
1399 Bool_t TGeant3::IsNewTrack() const
1402 // True if the track is not at the boundary of the current volume
1404 return (fGctrak->sleng==0);
1407 //_____________________________________________________________________________
1408 Bool_t TGeant3::IsTrackInside() const
1411 // True if the track is not at the boundary of the current volume
1413 return (fGctrak->inwvol==0);
1416 //_____________________________________________________________________________
1417 Bool_t TGeant3::IsTrackEntering() const
1420 // True if this is the first step of the track in the current volume
1422 return (fGctrak->inwvol==1);
1425 //_____________________________________________________________________________
1426 Bool_t TGeant3::IsTrackExiting() const
1429 // True if this is the last step of the track in the current volume
1431 return (fGctrak->inwvol==2);
1434 //_____________________________________________________________________________
1435 Bool_t TGeant3::IsTrackOut() const
1438 // True if the track is out of the setup
1440 return (fGctrak->inwvol==3);
1443 //_____________________________________________________________________________
1444 Bool_t TGeant3::IsTrackStop() const
1447 // True if the track energy has fallen below the threshold
1449 return (fGctrak->istop==2);
1452 //_____________________________________________________________________________
1453 Int_t TGeant3::NSecondaries() const
1456 // Number of secondary particles generated in the current step
1458 return fGcking->ngkine;
1461 //_____________________________________________________________________________
1462 Int_t TGeant3::CurrentEvent() const
1465 // Number of the current event
1467 return fGcflag->idevt;
1470 //_____________________________________________________________________________
1471 AliMCProcess TGeant3::ProdProcess(Int_t ) const
1474 // Name of the process that has produced the secondary particles
1475 // in the current step
1477 const AliMCProcess kIpProc[13] = { kPDecay, kPPair, kPCompton,
1478 kPPhotoelectric, kPBrem, kPDeltaRay,
1479 kPAnnihilation, kPHadronic,
1480 kPMuonNuclear, kPPhotoFission,
1481 kPRayleigh, kPCerenkov, kPSynchrotron};
1484 if(fGcking->ngkine>0)
1485 for (km = 0; km < fGctrak->nmec; ++km)
1486 for (im = 0; im < 13; ++im)
1487 if (G3toVMC(fGctrak->lmec[km]) == kIpProc[im])
1493 //_____________________________________________________________________________
1494 Int_t TGeant3::StepProcesses(TArrayI &proc) const
1497 // Return processes active in the current step
1500 Int_t nproc=Gctrak()->nmec;
1505 for (i=0; i<nproc; ++i)
1506 if((proc[nvproc]=G3toVMC(Gctrak()->lmec[i]))!=kPNoProcess) nvproc++;
1513 //_____________________________________________________________________________
1514 AliMCProcess TGeant3::G3toVMC(Int_t iproc) const
1517 // Conversion between GEANT and AliMC processes
1520 const AliMCProcess kPG2MC1[30] = {kPNoProcess, kPMultipleScattering, kPEnergyLoss, kPMagneticFieldL, kPDecay,
1521 kPPair, kPCompton, kPPhotoelectric, kPBrem, kPDeltaRay,
1522 kPAnnihilation, kPHadronic, kPNoProcess, kPEvaporation, kPNuclearFission,
1523 kPNuclearAbsorption, kPPbarAnnihilation, kPNCapture, kPHElastic, kPHInhelastic,
1524 kPMuonNuclear, kPTOFlimit, kPPhotoFission, kPNoProcess, kPRayleigh,
1525 kPNoProcess, kPNoProcess, kPNoProcess, kPNull, kPStop};
1527 const AliMCProcess kPG2MC2[9] = {kPLightAbsorption, kPLightScattering, kStepMax, kPNoProcess, kPCerenkov,
1528 kPLightReflection, kPLightRefraction, kPSynchrotron, kPNoProcess};
1530 AliMCProcess proc=kPNoProcess;
1531 if(1<iproc && iproc<=30) proc= kPG2MC1[iproc-1];
1532 else if(101<=iproc && iproc<=109) proc= kPG2MC2[iproc-100-1];
1537 //_____________________________________________________________________________
1538 void TGeant3::GetSecondary(Int_t isec, Int_t& ipart,
1539 TLorentzVector &x, TLorentzVector &p)
1542 // Get the parameters of the secondary track number isec produced
1543 // in the current step
1546 if(-1<isec && isec<fGcking->ngkine) {
1547 ipart=Int_t (fGcking->gkin[isec][4] +0.5);
1549 x[i]=fGckin3->gpos[isec][i];
1550 p[i]=fGcking->gkin[isec][i];
1552 x[3]=fGcking->tofd[isec];
1553 p[3]=fGcking->gkin[isec][3];
1555 printf(" * TGeant3::GetSecondary * Secondary %d does not exist\n",isec);
1556 x[0]=x[1]=x[2]=x[3]=p[0]=p[1]=p[2]=p[3]=0;
1561 //_____________________________________________________________________________
1562 void TGeant3::InitLego()
1565 // Set switches for lego transport
1568 SetDEBU(0,0,0); //do not print a message
1571 //_____________________________________________________________________________
1572 Bool_t TGeant3::IsTrackDisappeared() const
1575 // True if the current particle has disappered
1576 // either because it decayed or because it underwent
1577 // an inelastic collision
1579 return (fGctrak->istop==1);
1582 //_____________________________________________________________________________
1583 Bool_t TGeant3::IsTrackAlive() const
1586 // True if the current particle is alive and will continue to be
1589 return (fGctrak->istop==0);
1592 //_____________________________________________________________________________
1593 void TGeant3::StopTrack()
1596 // Stop the transport of the current particle and skip to the next
1601 //_____________________________________________________________________________
1602 void TGeant3::StopEvent()
1605 // Stop simulation of the current event and skip to the next
1610 //_____________________________________________________________________________
1611 Float_t TGeant3::MaxStep() const
1614 // Return the maximum step length in the current medium
1616 return fGctmed->stemax;
1619 //_____________________________________________________________________________
1620 void TGeant3::SetMaxStep(Float_t maxstep)
1623 // Set the maximum step allowed till the particle is in the current medium
1625 fGctmed->stemax=maxstep;
1628 //_____________________________________________________________________________
1629 void TGeant3::SetMaxNStep(Int_t maxnstp)
1632 // Set the maximum number of steps till the particle is in the current medium
1634 fGctrak->maxnst=maxnstp;
1637 //_____________________________________________________________________________
1638 Int_t TGeant3::GetMaxNStep() const
1641 // Maximum number of steps allowed in current medium
1643 return fGctrak->maxnst;
1646 //_____________________________________________________________________________
1647 void TGeant3::Material(Int_t& kmat, const char* name, Float_t a, Float_t z,
1648 Float_t dens, Float_t radl, Float_t absl, Float_t* buf,
1652 // Defines a Material
1654 // kmat number assigned to the material
1655 // name material name
1656 // a atomic mass in au
1658 // dens density in g/cm3
1659 // absl absorbtion length in cm
1660 // if >=0 it is ignored and the program
1661 // calculates it, if <0. -absl is taken
1662 // radl radiation length in cm
1663 // if >=0 it is ignored and the program
1664 // calculates it, if <0. -radl is taken
1665 // buf pointer to an array of user words
1666 // nbuf number of user words
1668 Int_t jmate=fGclink->jmate;
1674 for(i=1; i<=ns; i++) {
1675 if(fZlq[jmate-i]==0) {
1681 gsmate(kmat,PASSCHARD(name), a, z, dens, radl, absl, buf,
1682 nwbuf PASSCHARL(name));
1685 //_____________________________________________________________________________
1686 void TGeant3::Mixture(Int_t& kmat, const char* name, Float_t* a, Float_t* z,
1687 Float_t dens, Int_t nlmat, Float_t* wmat)
1690 // Defines mixture OR COMPOUND IMAT as composed by
1691 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
1693 // If NLMAT > 0 then wmat contains the proportion by
1694 // weights of each basic material in the mixture.
1696 // If nlmat < 0 then WMAT contains the number of atoms
1697 // of a given kind into the molecule of the COMPOUND
1698 // In this case, WMAT in output is changed to relative
1701 Int_t jmate=fGclink->jmate;
1707 for(i=1; i<=ns; i++) {
1708 if(fZlq[jmate-i]==0) {
1714 gsmixt(kmat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
1717 //_____________________________________________________________________________
1718 void TGeant3::Medium(Int_t& kmed, const char* name, Int_t nmat, Int_t isvol,
1719 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
1720 Float_t stemax, Float_t deemax, Float_t epsil,
1721 Float_t stmin, Float_t* ubuf, Int_t nbuf)
1724 // kmed tracking medium number assigned
1725 // name tracking medium name
1726 // nmat material number
1727 // isvol sensitive volume flag
1728 // ifield magnetic field
1729 // fieldm max. field value (kilogauss)
1730 // tmaxfd max. angle due to field (deg/step)
1731 // stemax max. step allowed
1732 // deemax max. fraction of energy lost in a step
1733 // epsil tracking precision (cm)
1734 // stmin min. step due to continuos processes (cm)
1736 // ifield = 0 if no magnetic field; ifield = -1 if user decision in guswim;
1737 // ifield = 1 if tracking performed with grkuta; ifield = 2 if tracking
1738 // performed with ghelix; ifield = 3 if tracking performed with ghelx3.
1740 Int_t jtmed=fGclink->jtmed;
1746 for(i=1; i<=ns; i++) {
1747 if(fZlq[jtmed-i]==0) {
1753 gstmed(kmed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1754 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1757 //_____________________________________________________________________________
1758 void TGeant3::Matrix(Int_t& krot, Float_t thex, Float_t phix, Float_t they,
1759 Float_t phiy, Float_t thez, Float_t phiz)
1762 // krot rotation matrix number assigned
1763 // theta1 polar angle for axis i
1764 // phi1 azimuthal angle for axis i
1765 // theta2 polar angle for axis ii
1766 // phi2 azimuthal angle for axis ii
1767 // theta3 polar angle for axis iii
1768 // phi3 azimuthal angle for axis iii
1770 // it defines the rotation matrix number irot.
1772 Int_t jrotm=fGclink->jrotm;
1778 for(i=1; i<=ns; i++) {
1779 if(fZlq[jrotm-i]==0) {
1785 gsrotm(krot, thex, phix, they, phiy, thez, phiz);
1788 //_____________________________________________________________________________
1789 Int_t TGeant3::GetMedium() const
1792 // Return the number of the current medium
1794 return fGctmed->numed;
1797 //_____________________________________________________________________________
1798 Float_t TGeant3::Edep() const
1801 // Return the energy lost in the current step
1803 return fGctrak->destep;
1806 //_____________________________________________________________________________
1807 Float_t TGeant3::Etot() const
1810 // Return the total energy of the current track
1812 return fGctrak->getot;
1815 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1817 // Functions from GBASE
1819 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1821 //____________________________________________________________________________
1822 void TGeant3::Gfile(const char *filename, const char *option)
1825 // Routine to open a GEANT/RZ data base.
1827 // LUN logical unit number associated to the file
1829 // CHFILE RZ file name
1831 // CHOPT is a character string which may be
1832 // N To create a new file
1833 // U to open an existing file for update
1834 // " " to open an existing file for read only
1835 // Q The initial allocation (default 1000 records)
1836 // is given in IQUEST(10)
1837 // X Open the file in exchange format
1838 // I Read all data structures from file to memory
1839 // O Write all data structures from memory to file
1842 // If options "I" or "O" all data structures are read or
1843 // written from/to file and the file is closed.
1844 // See routine GRMDIR to create subdirectories
1845 // See routines GROUT,GRIN to write,read objects
1847 grfile(21, PASSCHARD(filename), PASSCHARD(option) PASSCHARL(filename)
1851 //____________________________________________________________________________
1852 void TGeant3::Gpcxyz()
1855 // Print track and volume parameters at current point
1860 //_____________________________________________________________________________
1861 void TGeant3::Ggclos()
1864 // Closes off the geometry setting.
1865 // Initializes the search list for the contents of each
1866 // volume following the order they have been positioned, and
1867 // inserting the content '0' when a call to GSNEXT (-1) has
1868 // been required by the user.
1869 // Performs the development of the JVOLUM structure for all
1870 // volumes with variable parameters, by calling GGDVLP.
1871 // Interprets the user calls to GSORD, through GGORD.
1872 // Computes and stores in a bank (next to JVOLUM mother bank)
1873 // the number of levels in the geometrical tree and the
1874 // maximum number of contents per level, by calling GGNLEV.
1875 // Sets status bit for CONCAVE volumes, through GGCAVE.
1876 // Completes the JSET structure with the list of volume names
1877 // which identify uniquely a given physical detector, the
1878 // list of bit numbers to pack the corresponding volume copy
1879 // numbers, and the generic path(s) in the JVOLUM tree,
1880 // through the routine GHCLOS.
1883 // Create internal list of volumes
1884 fVolNames = new char[fGcnum->nvolum+1][5];
1886 for(i=0; i<fGcnum->nvolum; ++i) {
1887 strncpy(fVolNames[i], (char *) &fZiq[fGclink->jvolum+i+1], 4);
1888 fVolNames[i][4]='\0';
1890 strcpy(fVolNames[fGcnum->nvolum],"NULL");
1893 //_____________________________________________________________________________
1894 void TGeant3::Glast()
1897 // Finish a Geant run
1902 //_____________________________________________________________________________
1903 void TGeant3::Gprint(const char *name)
1906 // Routine to print data structures
1907 // CHNAME name of a data structure
1911 gprint(PASSCHARD(vname),0 PASSCHARL(vname));
1914 //_____________________________________________________________________________
1915 void TGeant3::Grun()
1918 // Steering function to process one run
1923 //_____________________________________________________________________________
1924 void TGeant3::Gtrig()
1927 // Steering function to process one event
1932 //_____________________________________________________________________________
1933 void TGeant3::Gtrigc()
1936 // Clear event partition
1941 //_____________________________________________________________________________
1942 void TGeant3::Gtrigi()
1945 // Initialises event partition
1950 //_____________________________________________________________________________
1951 void TGeant3::Gwork(Int_t nwork)
1954 // Allocates workspace in ZEBRA memory
1959 //_____________________________________________________________________________
1960 void TGeant3::Gzinit()
1963 // To initialise GEANT/ZEBRA data structures
1968 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1970 // Functions from GCONS
1972 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1974 //_____________________________________________________________________________
1975 void TGeant3::Gfmate(Int_t imat, char *name, Float_t &a, Float_t &z,
1976 Float_t &dens, Float_t &radl, Float_t &absl,
1977 Float_t* ubuf, Int_t& nbuf)
1980 // Return parameters for material IMAT
1982 gfmate(imat, PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
1986 //_____________________________________________________________________________
1987 void TGeant3::Gfpart(Int_t ipart, char *name, Int_t &itrtyp,
1988 Float_t &amass, Float_t &charge, Float_t &tlife)
1991 // Return parameters for particle of type IPART
1995 Int_t igpart = IdFromPDG(ipart);
1996 gfpart(igpart, PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
2000 //_____________________________________________________________________________
2001 void TGeant3::Gftmed(Int_t numed, char *name, Int_t &nmat, Int_t &isvol,
2002 Int_t &ifield, Float_t &fieldm, Float_t &tmaxfd,
2003 Float_t &stemax, Float_t &deemax, Float_t &epsil,
2004 Float_t &stmin, Float_t *ubuf, Int_t *nbuf)
2007 // Return parameters for tracking medium NUMED
2009 gftmed(numed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
2010 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
2014 void TGeant3::Gftmat(Int_t imate, Int_t ipart, char *chmeca, Int_t kdim,
2015 Float_t* tkin, Float_t* value, Float_t* pcut,
2019 // Return parameters for tracking medium NUMED
2021 gftmat(imate, ipart, PASSCHARD(chmeca), kdim,
2022 tkin, value, pcut, ixst PASSCHARL(chmeca));
2026 //_____________________________________________________________________________
2027 Float_t TGeant3::Gbrelm(Float_t z, Float_t t, Float_t bcut)
2030 // To calculate energy loss due to soft muon BREMSSTRAHLUNG
2032 return gbrelm(z,t,bcut);
2035 //_____________________________________________________________________________
2036 Float_t TGeant3::Gprelm(Float_t z, Float_t t, Float_t bcut)
2039 // To calculate DE/DX in GeV*barn/atom for direct pair production by muons
2041 return gprelm(z,t,bcut);
2044 //_____________________________________________________________________________
2045 void TGeant3::Gmate()
2048 // Define standard GEANT materials
2053 //_____________________________________________________________________________
2054 void TGeant3::Gpart()
2057 // Define standard GEANT particles plus selected decay modes
2058 // and branching ratios.
2063 //_____________________________________________________________________________
2064 void TGeant3::Gsdk(Int_t ipart, Float_t *bratio, Int_t *mode)
2066 // Defines branching ratios and decay modes for standard
2068 gsdk(ipart,bratio,mode);
2071 //_____________________________________________________________________________
2072 void TGeant3::Gsmate(Int_t imat, const char *name, Float_t a, Float_t z,
2073 Float_t dens, Float_t radl, Float_t absl)
2076 // Defines a Material
2078 // kmat number assigned to the material
2079 // name material name
2080 // a atomic mass in au
2082 // dens density in g/cm3
2083 // absl absorbtion length in cm
2084 // if >=0 it is ignored and the program
2085 // calculates it, if <0. -absl is taken
2086 // radl radiation length in cm
2087 // if >=0 it is ignored and the program
2088 // calculates it, if <0. -radl is taken
2089 // buf pointer to an array of user words
2090 // nbuf number of user words
2094 gsmate(imat,PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
2098 //_____________________________________________________________________________
2099 void TGeant3::Gsmixt(Int_t imat, const char *name, Float_t *a, Float_t *z,
2100 Float_t dens, Int_t nlmat, Float_t *wmat)
2103 // Defines mixture OR COMPOUND IMAT as composed by
2104 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
2106 // If NLMAT.GT.0 then WMAT contains the PROPORTION BY
2107 // WEIGTHS OF EACH BASIC MATERIAL IN THE MIXTURE.
2109 // If NLMAT.LT.0 then WMAT contains the number of atoms
2110 // of a given kind into the molecule of the COMPOUND
2111 // In this case, WMAT in output is changed to relative
2114 gsmixt(imat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
2117 //_____________________________________________________________________________
2118 void TGeant3::Gspart(Int_t ipart, const char *name, Int_t itrtyp,
2119 Float_t amass, Float_t charge, Float_t tlife)
2122 // Store particle parameters
2124 // ipart particle code
2125 // name particle name
2126 // itrtyp transport method (see GEANT manual)
2127 // amass mass in GeV/c2
2128 // charge charge in electron units
2129 // tlife lifetime in seconds
2133 gspart(ipart,PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
2137 //_____________________________________________________________________________
2138 void TGeant3::Gstmed(Int_t numed, const char *name, Int_t nmat, Int_t isvol,
2139 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
2140 Float_t stemax, Float_t deemax, Float_t epsil,
2144 // NTMED Tracking medium number
2145 // NAME Tracking medium name
2146 // NMAT Material number
2147 // ISVOL Sensitive volume flag
2148 // IFIELD Magnetic field
2149 // FIELDM Max. field value (Kilogauss)
2150 // TMAXFD Max. angle due to field (deg/step)
2151 // STEMAX Max. step allowed
2152 // DEEMAX Max. fraction of energy lost in a step
2153 // EPSIL Tracking precision (cm)
2154 // STMIN Min. step due to continuos processes (cm)
2156 // IFIELD = 0 if no magnetic field; IFIELD = -1 if user decision in GUSWIM;
2157 // IFIELD = 1 if tracking performed with GRKUTA; IFIELD = 2 if tracking
2158 // performed with GHELIX; IFIELD = 3 if tracking performed with GHELX3.
2162 gstmed(numed,PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
2163 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
2166 //_____________________________________________________________________________
2167 void TGeant3::Gsckov(Int_t itmed, Int_t npckov, Float_t *ppckov,
2168 Float_t *absco, Float_t *effic, Float_t *rindex)
2171 // Stores the tables for UV photon tracking in medium ITMED
2172 // Please note that it is the user's responsability to
2173 // provide all the coefficients:
2176 // ITMED Tracking medium number
2177 // NPCKOV Number of bins of each table
2178 // PPCKOV Value of photon momentum (in GeV)
2179 // ABSCO Absorbtion coefficients
2180 // dielectric: absorbtion length in cm
2181 // metals : absorbtion fraction (0<=x<=1)
2182 // EFFIC Detection efficiency for UV photons
2183 // RINDEX Refraction index (if=0 metal)
2185 gsckov(itmed,npckov,ppckov,absco,effic,rindex);
2188 //_____________________________________________________________________________
2189 void TGeant3::SetCerenkov(Int_t itmed, Int_t npckov, Float_t *ppckov,
2190 Float_t *absco, Float_t *effic, Float_t *rindex)
2193 // Stores the tables for UV photon tracking in medium ITMED
2194 // Please note that it is the user's responsability to
2195 // provide all the coefficients:
2198 // ITMED Tracking medium number
2199 // NPCKOV Number of bins of each table
2200 // PPCKOV Value of photon momentum (in GeV)
2201 // ABSCO Absorbtion coefficients
2202 // dielectric: absorbtion length in cm
2203 // metals : absorbtion fraction (0<=x<=1)
2204 // EFFIC Detection efficiency for UV photons
2205 // RINDEX Refraction index (if=0 metal)
2207 gsckov(itmed,npckov,ppckov,absco,effic,rindex);
2210 //_____________________________________________________________________________
2211 void TGeant3::Gstpar(Int_t itmed, const char *param, Float_t parval)
2214 // To change the value of cut or mechanism "CHPAR"
2215 // to a new value PARVAL for tracking medium ITMED
2216 // The data structure JTMED contains the standard tracking
2217 // parameters (CUTS and flags to control the physics processes) which
2218 // are used by default for all tracking media. It is possible to
2219 // redefine individually with GSTPAR any of these parameters for a
2220 // given tracking medium.
2221 // ITMED tracking medium number
2222 // CHPAR is a character string (variable name)
2223 // PARVAL must be given as a floating point.
2225 gstpar(itmed,PASSCHARD(param), parval PASSCHARL(param));
2228 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2230 // Functions from GCONS
2232 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2234 //_____________________________________________________________________________
2235 void TGeant3::Gfkine(Int_t itra, Float_t *vert, Float_t *pvert, Int_t &ipart,
2238 // Storing/Retrieving Vertex and Track parameters
2239 // ----------------------------------------------
2241 // Stores vertex parameters.
2242 // VERT array of (x,y,z) position of the vertex
2243 // NTBEAM beam track number origin of the vertex
2244 // =0 if none exists
2245 // NTTARG target track number origin of the vertex
2246 // UBUF user array of NUBUF floating point numbers
2248 // NVTX new vertex number (=0 in case of error).
2249 // Prints vertex parameters.
2250 // IVTX for vertex IVTX.
2251 // (For all vertices if IVTX=0)
2252 // Stores long life track parameters.
2253 // PLAB components of momentum
2254 // IPART type of particle (see GSPART)
2255 // NV vertex number origin of track
2256 // UBUF array of NUBUF floating point user parameters
2258 // NT track number (if=0 error).
2259 // Retrieves long life track parameters.
2260 // ITRA track number for which parameters are requested
2261 // VERT vector origin of the track
2262 // PVERT 4 momentum components at the track origin
2263 // IPART particle type (=0 if track ITRA does not exist)
2264 // NVERT vertex number origin of the track
2265 // UBUF user words stored in GSKINE.
2266 // Prints initial track parameters.
2267 // ITRA for track ITRA
2268 // (For all tracks if ITRA=0)
2272 gfkine(itra,vert,pvert,ipart,nvert,ubuf,nbuf);
2275 //_____________________________________________________________________________
2276 void TGeant3::Gfvert(Int_t nvtx, Float_t *v, Int_t &ntbeam, Int_t &nttarg,
2280 // Retrieves the parameter of a vertex bank
2281 // Vertex is generated from tracks NTBEAM NTTARG
2282 // NVTX is the new vertex number
2286 gfvert(nvtx,v,ntbeam,nttarg,tofg,ubuf,nbuf);
2289 //_____________________________________________________________________________
2290 Int_t TGeant3::Gskine(Float_t *plab, Int_t ipart, Int_t nv, Float_t *buf,
2294 // Store kinematics of track NT into data structure
2295 // Track is coming from vertex NV
2298 gskine(plab, ipart, nv, buf, nwbuf, nt);
2302 //_____________________________________________________________________________
2303 Int_t TGeant3::Gsvert(Float_t *v, Int_t ntbeam, Int_t nttarg, Float_t *ubuf,
2307 // Creates a new vertex bank
2308 // Vertex is generated from tracks NTBEAM NTTARG
2309 // NVTX is the new vertex number
2312 gsvert(v, ntbeam, nttarg, ubuf, nwbuf, nwtx);
2316 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2318 // Functions from GPHYS
2320 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2322 //_____________________________________________________________________________
2323 void TGeant3::Gphysi()
2326 // Initialise material constants for all the physics
2327 // mechanisms used by GEANT
2332 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2334 // Functions from GTRAK
2336 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2338 //_____________________________________________________________________________
2339 void TGeant3::Gdebug()
2342 // Debug the current step
2347 //_____________________________________________________________________________
2348 void TGeant3::Gekbin()
2351 // To find bin number in kinetic energy table
2352 // stored in ELOW(NEKBIN)
2357 //_____________________________________________________________________________
2358 void TGeant3::Gfinds()
2361 // Returns the set/volume parameters corresponding to
2362 // the current space point in /GCTRAK/
2363 // and fill common /GCSETS/
2365 // IHSET user set identifier
2366 // IHDET user detector identifier
2367 // ISET set number in JSET
2368 // IDET detector number in JS=LQ(JSET-ISET)
2369 // IDTYPE detector type (1,2)
2370 // NUMBV detector volume numbers (array of length NVNAME)
2371 // NVNAME number of volume levels
2376 //_____________________________________________________________________________
2377 void TGeant3::Gsking(Int_t igk)
2380 // Stores in stack JSTAK either the IGKth track of /GCKING/,
2381 // or the NGKINE tracks when IGK is 0.
2386 //_____________________________________________________________________________
2387 void TGeant3::Gskpho(Int_t igk)
2390 // Stores in stack JSTAK either the IGKth Cherenkov photon of
2391 // /GCKIN2/, or the NPHOT tracks when IGK is 0.
2396 //_____________________________________________________________________________
2397 void TGeant3::Gsstak(Int_t iflag)
2400 // Stores in auxiliary stack JSTAK the particle currently
2401 // described in common /GCKINE/.
2403 // On request, creates also an entry in structure JKINE :
2405 // 0 : No entry in JKINE structure required (user)
2406 // 1 : New entry in JVERTX / JKINE structures required (user)
2407 // <0 : New entry in JKINE structure at vertex -IFLAG (user)
2408 // 2 : Entry in JKINE structure exists already (from GTREVE)
2413 //_____________________________________________________________________________
2414 void TGeant3::Gsxyz()
2417 // Store space point VECT in banks JXYZ
2422 //_____________________________________________________________________________
2423 void TGeant3::Gtrack()
2426 // Controls tracking of current particle
2431 //_____________________________________________________________________________
2432 void TGeant3::Gtreve()
2435 // Controls tracking of all particles belonging to the current event
2440 //_____________________________________________________________________________
2441 void TGeant3::GtreveRoot()
2444 // Controls tracking of all particles belonging to the current event
2449 //_____________________________________________________________________________
2450 void TGeant3::Grndm(Float_t *rvec, const Int_t len) const
2453 // To generate a vector RVECV of LEN random numbers
2454 // Copy of the CERN Library routine RANECU
2458 //_____________________________________________________________________________
2459 void TGeant3::Grndmq(Int_t &/*is1*/, Int_t &/*is2*/, const Int_t /*iseq*/,
2460 const Text_t */*chopt*/)
2463 // To set/retrieve the seed of the random number generator
2465 /*printf("Dummy grndmq called\n");*/
2468 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2470 // Functions from GDRAW
2472 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2474 //_____________________________________________________________________________
2475 void TGeant3::Gdxyz(Int_t it)
2478 // Draw the points stored with Gsxyz relative to track it
2483 //_____________________________________________________________________________
2484 void TGeant3::Gdcxyz()
2487 // Draw the position of the current track
2492 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2494 // Functions from GGEOM
2496 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2498 //_____________________________________________________________________________
2499 void TGeant3::Gdtom(Float_t *xd, Float_t *xm, Int_t iflag)
2502 // Computes coordinates XM (Master Reference System
2503 // knowing the coordinates XD (Detector Ref System)
2504 // The local reference system can be initialized by
2505 // - the tracking routines and GDTOM used in GUSTEP
2506 // - a call to GSCMED(NLEVEL,NAMES,NUMBER)
2507 // (inverse routine is GMTOD)
2509 // If IFLAG=1 convert coordinates
2510 // IFLAG=2 convert direction cosinus
2512 gdtom(xd, xm, iflag);
2515 //_____________________________________________________________________________
2516 void TGeant3::Glmoth(const char* iudet, Int_t iunum, Int_t &nlev, Int_t *lvols,
2520 // Loads the top part of the Volume tree in LVOLS (IVO's),
2521 // LINDX (IN indices) for a given volume defined through
2522 // its name IUDET and number IUNUM.
2524 // The routine stores only upto the last level where JVOLUM
2525 // data structure is developed. If there is no development
2526 // above the current level, it returns NLEV zero.
2528 glmoth(PASSCHARD(iudet), iunum, nlev, lvols, lindx, idum PASSCHARL(iudet));
2531 //_____________________________________________________________________________
2532 void TGeant3::Gmedia(Float_t *x, Int_t &numed)
2535 // Finds in which volume/medium the point X is, and updates the
2536 // common /GCVOLU/ and the structure JGPAR accordingly.
2538 // NUMED returns the tracking medium number, or 0 if point is
2539 // outside the experimental setup.
2544 //_____________________________________________________________________________
2545 void TGeant3::Gmtod(Float_t *xm, Float_t *xd, Int_t iflag)
2548 // Computes coordinates XD (in DRS)
2549 // from known coordinates XM in MRS
2550 // The local reference system can be initialized by
2551 // - the tracking routines and GMTOD used in GUSTEP
2552 // - a call to GMEDIA(XM,NUMED)
2553 // - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER)
2554 // (inverse routine is GDTOM)
2556 // If IFLAG=1 convert coordinates
2557 // IFLAG=2 convert direction cosinus
2559 gmtod(xm, xd, iflag);
2562 //_____________________________________________________________________________
2563 void TGeant3::Gsdvn(const char *name, const char *mother, Int_t ndiv,
2567 // Create a new volume by dividing an existing one
2570 // MOTHER Mother volume name
2571 // NDIV Number of divisions
2574 // X,Y,Z of CAXIS will be translated to 1,2,3 for IAXIS.
2575 // It divides a previously defined volume.
2580 Vname(mother,vmother);
2581 gsdvn(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis PASSCHARL(vname)
2582 PASSCHARL(vmother));
2585 //_____________________________________________________________________________
2586 void TGeant3::Gsdvn2(const char *name, const char *mother, Int_t ndiv,
2587 Int_t iaxis, Float_t c0i, Int_t numed)
2590 // Create a new volume by dividing an existing one
2592 // Divides mother into ndiv divisions called name
2593 // along axis iaxis starting at coordinate value c0.
2594 // the new volume created will be medium number numed.
2599 Vname(mother,vmother);
2600 gsdvn2(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis, c0i, numed
2601 PASSCHARL(vname) PASSCHARL(vmother));
2604 //_____________________________________________________________________________
2605 void TGeant3::Gsdvs(const char *name, const char *mother, Float_t step,
2606 Int_t iaxis, Int_t numed)
2609 // Create a new volume by dividing an existing one
2614 Vname(mother,vmother);
2615 gsdvs(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed
2616 PASSCHARL(vname) PASSCHARL(vmother));
2619 //_____________________________________________________________________________
2620 void TGeant3::Gsdvs2(const char *name, const char *mother, Float_t step,
2621 Int_t iaxis, Float_t c0, Int_t numed)
2624 // Create a new volume by dividing an existing one
2629 Vname(mother,vmother);
2630 gsdvs2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0, numed
2631 PASSCHARL(vname) PASSCHARL(vmother));
2634 //_____________________________________________________________________________
2635 void TGeant3::Gsdvt(const char *name, const char *mother, Float_t step,
2636 Int_t iaxis, Int_t numed, Int_t ndvmx)
2639 // Create a new volume by dividing an existing one
2641 // Divides MOTHER into divisions called NAME along
2642 // axis IAXIS in steps of STEP. If not exactly divisible
2643 // will make as many as possible and will centre them
2644 // with respect to the mother. Divisions will have medium
2645 // number NUMED. If NUMED is 0, NUMED of MOTHER is taken.
2646 // NDVMX is the expected maximum number of divisions
2647 // (If 0, no protection tests are performed)
2652 Vname(mother,vmother);
2653 gsdvt(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed, ndvmx
2654 PASSCHARL(vname) PASSCHARL(vmother));
2657 //_____________________________________________________________________________
2658 void TGeant3::Gsdvt2(const char *name, const char *mother, Float_t step,
2659 Int_t iaxis, Float_t c0, Int_t numed, Int_t ndvmx)
2662 // Create a new volume by dividing an existing one
2664 // Divides MOTHER into divisions called NAME along
2665 // axis IAXIS starting at coordinate value C0 with step
2667 // The new volume created will have medium number NUMED.
2668 // If NUMED is 0, NUMED of mother is taken.
2669 // NDVMX is the expected maximum number of divisions
2670 // (If 0, no protection tests are performed)
2675 Vname(mother,vmother);
2676 gsdvt2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0,
2677 numed, ndvmx PASSCHARL(vname) PASSCHARL(vmother));
2680 //_____________________________________________________________________________
2681 void TGeant3::Gsord(const char *name, Int_t iax)
2684 // Flags volume CHNAME whose contents will have to be ordered
2685 // along axis IAX, by setting the search flag to -IAX
2689 // IAX = 4 Rxy (static ordering only -> GTMEDI)
2690 // IAX = 14 Rxy (also dynamic ordering -> GTNEXT)
2691 // IAX = 5 Rxyz (static ordering only -> GTMEDI)
2692 // IAX = 15 Rxyz (also dynamic ordering -> GTNEXT)
2693 // IAX = 6 PHI (PHI=0 => X axis)
2694 // IAX = 7 THETA (THETA=0 => Z axis)
2698 gsord(PASSCHARD(vname), iax PASSCHARL(vname));
2701 //_____________________________________________________________________________
2702 void TGeant3::Gspos(const char *name, Int_t nr, const char *mother, Float_t x,
2703 Float_t y, Float_t z, Int_t irot, const char *konly)
2706 // Position a volume into an existing one
2709 // NUMBER Copy number of the volume
2710 // MOTHER Mother volume name
2711 // X X coord. of the volume in mother ref. sys.
2712 // Y Y coord. of the volume in mother ref. sys.
2713 // Z Z coord. of the volume in mother ref. sys.
2714 // IROT Rotation matrix number w.r.t. mother ref. sys.
2715 // ONLY ONLY/MANY flag
2717 // It positions a previously defined volume in the mother.
2723 Vname(mother,vmother);
2724 gspos(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2725 PASSCHARD(konly) PASSCHARL(vname) PASSCHARL(vmother)
2729 //_____________________________________________________________________________
2730 void TGeant3::Gsposp(const char *name, Int_t nr, const char *mother,
2731 Float_t x, Float_t y, Float_t z, Int_t irot,
2732 const char *konly, Float_t *upar, Int_t np )
2735 // Place a copy of generic volume NAME with user number
2736 // NR inside MOTHER, with its parameters UPAR(1..NP)
2741 Vname(mother,vmother);
2742 gsposp(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2743 PASSCHARD(konly), upar, np PASSCHARL(vname) PASSCHARL(vmother)
2747 //_____________________________________________________________________________
2748 void TGeant3::Gsrotm(Int_t nmat, Float_t theta1, Float_t phi1, Float_t theta2,
2749 Float_t phi2, Float_t theta3, Float_t phi3)
2752 // nmat Rotation matrix number
2753 // THETA1 Polar angle for axis I
2754 // PHI1 Azimuthal angle for axis I
2755 // THETA2 Polar angle for axis II
2756 // PHI2 Azimuthal angle for axis II
2757 // THETA3 Polar angle for axis III
2758 // PHI3 Azimuthal angle for axis III
2760 // It defines the rotation matrix number IROT.
2762 gsrotm(nmat, theta1, phi1, theta2, phi2, theta3, phi3);
2765 //_____________________________________________________________________________
2766 void TGeant3::Gprotm(Int_t nmat)
2769 // To print rotation matrices structure JROTM
2770 // nmat Rotation matrix number
2775 //_____________________________________________________________________________
2776 Int_t TGeant3::Gsvolu(const char *name, const char *shape, Int_t nmed,
2777 Float_t *upar, Int_t npar)
2781 // SHAPE Volume type
2782 // NUMED Tracking medium number
2783 // NPAR Number of shape parameters
2784 // UPAR Vector containing shape parameters
2786 // It creates a new volume in the JVOLUM data structure.
2792 Vname(shape,vshape);
2793 gsvolu(PASSCHARD(vname), PASSCHARD(vshape), nmed, upar, npar, ivolu
2794 PASSCHARL(vname) PASSCHARL(vshape));
2798 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2800 // T H E D R A W I N G P A C K A G E
2801 // ======================================
2802 // Drawing functions. These functions allow the visualization in several ways
2803 // of the volumes defined in the geometrical data structure. It is possible
2804 // to draw the logical tree of volumes belonging to the detector (DTREE),
2805 // to show their geometrical specification (DSPEC,DFSPC), to draw them
2806 // and their cut views (DRAW, DCUT). Moreover, it is possible to execute
2807 // these commands when the hidden line removal option is activated; in
2808 // this case, the volumes can be also either translated in the space
2809 // (SHIFT), or clipped by boolean operation (CVOL). In addition, it is
2810 // possible to fill the surfaces of the volumes
2811 // with solid colours when the shading option (SHAD) is activated.
2812 // Several tools (ZOOM, LENS) have been developed to zoom detailed parts
2813 // of the detectors or to scan physical events as well.
2814 // Finally, the command MOVE will allow the rotation, translation and zooming
2815 // on real time parts of the detectors or tracks and hits of a simulated event.
2816 // Ray-tracing commands. In case the command (DOPT RAYT ON) is executed,
2817 // the drawing is performed by the Geant ray-tracing;
2818 // automatically, the color is assigned according to the tracking medium of each
2819 // volume and the volumes with a density lower/equal than the air are considered
2820 // transparent; if the option (USER) is set (ON) (again via the command (DOPT)),
2821 // the user can set color and visibility for the desired volumes via the command
2822 // (SATT), as usual, relatively to the attributes (COLO) and (SEEN).
2823 // The resolution can be set via the command (SATT * FILL VALUE), where (VALUE)
2824 // is the ratio between the number of pixels drawn and 20 (user coordinates).
2825 // Parallel view and perspective view are possible (DOPT PROJ PARA/PERS); in the
2826 // first case, we assume that the first mother volume of the tree is a box with
2827 // dimensions 10000 X 10000 X 10000 cm and the view point (infinetely far) is
2828 // 5000 cm far from the origin along the Z axis of the user coordinates; in the
2829 // second case, the distance between the observer and the origin of the world
2830 // reference system is set in cm by the command (PERSP NAME VALUE); grand-angle
2831 // or telescopic effects can be achieved changing the scale factors in the command
2832 // (DRAW). When the final picture does not occupy the full window,
2833 // mapping the space before tracing can speed up the drawing, but can also
2834 // produce less precise results; values from 1 to 4 are allowed in the command
2835 // (DOPT MAPP VALUE), the mapping being more precise for increasing (VALUE); for
2836 // (VALUE = 0) no mapping is performed (therefore max precision and lowest speed).
2837 // The command (VALCUT) allows the cutting of the detector by three planes
2838 // ortogonal to the x,y,z axis. The attribute (LSTY) can be set by the command
2839 // SATT for any desired volume and can assume values from 0 to 7; it determines
2840 // the different light processing to be performed for different materials:
2841 // 0 = dark-matt, 1 = bright-matt, 2 = plastic, 3 = ceramic, 4 = rough-metals,
2842 // 5 = shiny-metals, 6 = glass, 7 = mirror. The detector is assumed to be in the
2843 // dark, the ambient light luminosity is 0.2 for each basic hue (the saturation
2844 // is 0.9) and the observer is assumed to have a light source (therefore he will
2845 // produce parallel light in the case of parallel view and point-like-source
2846 // light in the case of perspective view).
2848 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2850 //_____________________________________________________________________________
2851 void TGeant3::Gsatt(const char *name, const char *att, Int_t val)
2855 // IOPT Name of the attribute to be set
2856 // IVAL Value to which the attribute is to be set
2858 // name= "*" stands for all the volumes.
2859 // iopt can be chosen among the following :
2861 // WORK 0=volume name is inactive for the tracking
2862 // 1=volume name is active for the tracking (default)
2864 // SEEN 0=volume name is invisible
2865 // 1=volume name is visible (default)
2866 // -1=volume invisible with all its descendants in the tree
2867 // -2=volume visible but not its descendants in the tree
2869 // LSTY line style 1,2,3,... (default=1)
2870 // LSTY=7 will produce a very precise approximation for
2871 // revolution bodies.
2873 // LWID line width -7,...,1,2,3,..7 (default=1)
2874 // LWID<0 will act as abs(LWID) was set for the volume
2875 // and for all the levels below it. When SHAD is 'ON', LWID
2876 // represent the linewidth of the scan lines filling the surfaces
2877 // (whereas the FILL value represent their number). Therefore
2878 // tuning this parameter will help to obtain the desired
2879 // quality/performance ratio.
2881 // COLO colour code -166,...,1,2,..166 (default=1)
2883 // n=2=red; n=17+m, m=0,25, increasing luminosity according to 'm';
2884 // n=3=green; n=67+m, m=0,25, increasing luminosity according to 'm';
2885 // n=4=blue; n=117+m, m=0,25, increasing luminosity according to 'm';
2886 // n=5=yellow; n=42+m, m=0,25, increasing luminosity according to 'm';
2887 // n=6=violet; n=142+m, m=0,25, increasing luminosity according to 'm';
2888 // n=7=lightblue; n=92+m, m=0,25, increasing luminosity according to 'm';
2889 // colour=n*10+m, m=1,2,...9, will produce the same colour
2890 // as 'n', but with increasing luminosity according to 'm';
2891 // COLO<0 will act as if abs(COLO) was set for the volume
2892 // and for all the levels below it.
2893 // When for a volume the attribute FILL is > 1 (and the
2894 // option SHAD is on), the ABS of its colour code must be < 8
2895 // because an automatic shading of its faces will be
2898 // FILL (1992) fill area -7,...,0,1,...7 (default=0)
2899 // when option SHAD is "on" the FILL attribute of any
2900 // volume can be set different from 0 (normal drawing);
2901 // if it is set to 1, the faces of such volume will be filled
2902 // with solid colours; if ABS(FILL) is > 1, then a light
2903 // source is placed along the observer line, and the faces of
2904 // such volumes will be painted by colours whose luminosity
2905 // will depend on the amount of light reflected;
2906 // if ABS(FILL) = 1, then it is possible to use all the 166
2907 // colours of the colour table, becouse the automatic shading
2908 // is not performed;
2909 // for increasing values of FILL the drawing will be performed
2910 // with higher and higher resolution improving the quality (the
2911 // number of scan lines used to fill the faces increases with FILL);
2912 // it is possible to set different values of FILL
2913 // for different volumes, in order to optimize at the same time
2914 // the performance and the quality of the picture;
2915 // FILL<0 will act as if abs(FILL) was set for the volume
2916 // and for all the levels below it.
2917 // This kind of drawing can be saved in 'picture files'
2918 // or in view banks.
2919 // 0=drawing without fill area
2920 // 1=faces filled with solid colours and resolution = 6
2921 // 2=lowest resolution (very fast)
2922 // 3=default resolution
2923 // 4=.................
2924 // 5=.................
2925 // 6=.................
2927 // Finally, if a coloured background is desired, the FILL
2928 // attribute for the first volume of the tree must be set
2929 // equal to -abs(colo), colo being >0 and <166.
2931 // SET set number associated to volume name
2932 // DET detector number associated to volume name
2933 // DTYP detector type (1,2)
2940 gsatt(PASSCHARD(vname), PASSCHARD(vatt), val PASSCHARL(vname)
2944 //_____________________________________________________________________________
2945 void TGeant3::Gfpara(const char *name, Int_t number, Int_t intext, Int_t& npar,
2946 Int_t& natt, Float_t* par, Float_t* att)
2949 // Find the parameters of a volume
2951 gfpara(PASSCHARD(name), number, intext, npar, natt, par, att
2955 //_____________________________________________________________________________
2956 void TGeant3::Gckpar(Int_t ish, Int_t npar, Float_t* par)
2959 // Check the parameters of a shape
2961 gckpar(ish,npar,par);
2964 //_____________________________________________________________________________
2965 void TGeant3::Gckmat(Int_t itmed, char* natmed)
2968 // Check the parameters of a tracking medium
2970 gckmat(itmed, PASSCHARD(natmed) PASSCHARL(natmed));
2973 //_____________________________________________________________________________
2974 Int_t TGeant3::Glvolu(Int_t nlev, Int_t *lnam,Int_t *lnum)
2977 // nlev number of leveles deap into the volume tree
2978 // size of the arrays lnam and lnum
2979 // lnam an integer array whos 4 bytes contain the askii code for the
2981 // lnum an integer array containing the copy numbers for that specific
2984 // This routine fills the volulme paramters in common /gcvolu/ for a
2985 // physical tree, specified by the list lnam and lnum of volume names
2986 // and numbers, and for all its ascendants up to level 1. This routine
2987 // is optimsed and does not re-compute the part of the history already
2988 // available in GCVOLU. This means that if it is used in user programs
2989 // outside the usual framwork of the tracking, the user has to initilise
2990 // to zero NLEVEL in the common GCVOLU. It return 0 if there were no
2991 // problems in make the call.
2994 glvolu(nlev, lnam, lnum, ier);
2998 //_____________________________________________________________________________
2999 void TGeant3::Gdelete(Int_t iview)
3002 // IVIEW View number
3004 // It deletes a view bank from memory.
3009 //_____________________________________________________________________________
3010 void TGeant3::Gdopen(Int_t iview)
3013 // IVIEW View number
3015 // When a drawing is very complex and requires a long time to be
3016 // executed, it can be useful to store it in a view bank: after a
3017 // call to DOPEN and the execution of the drawing (nothing will
3018 // appear on the screen), and after a necessary call to DCLOSE,
3019 // the contents of the bank can be displayed in a very fast way
3020 // through a call to DSHOW; therefore, the detector can be easily
3021 // zoomed many times in different ways. Please note that the pictures
3022 // with solid colours can now be stored in a view bank or in 'PICTURE FILES'
3029 //_____________________________________________________________________________
3030 void TGeant3::Gdclose()
3033 // It closes the currently open view bank; it must be called after the
3034 // end of the drawing to be stored.
3039 //_____________________________________________________________________________
3040 void TGeant3::Gdshow(Int_t iview)
3043 // IVIEW View number
3045 // It shows on the screen the contents of a view bank. It
3046 // can be called after a view bank has been closed.
3051 //_____________________________________________________________________________
3052 void TGeant3::Gdopt(const char *name,const char *value)
3056 // VALUE Option value
3058 // To set/modify the drawing options.
3061 // THRZ ON Draw tracks in R vs Z
3062 // OFF (D) Draw tracks in X,Y,Z
3065 // PROJ PARA (D) Parallel projection
3067 // TRAK LINE (D) Trajectory drawn with lines
3068 // POIN " " with markers
3069 // HIDE ON Hidden line removal using the CG package
3070 // OFF (D) No hidden line removal
3071 // SHAD ON Fill area and shading of surfaces.
3072 // OFF (D) Normal hidden line removal.
3073 // RAYT ON Ray-tracing on.
3074 // OFF (D) Ray-tracing off.
3075 // EDGE OFF Does not draw contours when shad is on.
3076 // ON (D) Normal shading.
3077 // MAPP 1,2,3,4 Mapping before ray-tracing.
3078 // 0 (D) No mapping.
3079 // USER ON User graphics options in the raytracing.
3080 // OFF (D) Automatic graphics options.
3086 Vname(value,vvalue);
3087 gdopt(PASSCHARD(vname), PASSCHARD(vvalue) PASSCHARL(vname)
3091 //_____________________________________________________________________________
3092 void TGeant3::Gdraw(const char *name,Float_t theta, Float_t phi, Float_t psi,
3093 Float_t u0,Float_t v0,Float_t ul,Float_t vl)
3098 // THETA Viewing angle theta (for 3D projection)
3099 // PHI Viewing angle phi (for 3D projection)
3100 // PSI Viewing angle psi (for 2D rotation)
3101 // U0 U-coord. (horizontal) of volume origin
3102 // V0 V-coord. (vertical) of volume origin
3103 // SU Scale factor for U-coord.
3104 // SV Scale factor for V-coord.
3106 // This function will draw the volumes,
3107 // selected with their graphical attributes, set by the Gsatt
3108 // facility. The drawing may be performed with hidden line removal
3109 // and with shading effects according to the value of the options HIDE
3110 // and SHAD; if the option SHAD is ON, the contour's edges can be
3111 // drawn or not. If the option HIDE is ON, the detector can be
3112 // exploded (BOMB), clipped with different shapes (CVOL), and some
3113 // of its parts can be shifted from their original
3114 // position (SHIFT). When HIDE is ON, if
3115 // the drawing requires more than the available memory, the program
3116 // will evaluate and display the number of missing words
3117 // (so that the user can increase the
3118 // size of its ZEBRA store). Finally, at the end of each drawing (with HIDE on),
3119 // the program will print messages about the memory used and
3120 // statistics on the volumes' visibility.
3121 // The following commands will produce the drawing of a green
3122 // volume, specified by NAME, without using the hidden line removal
3123 // technique, using the hidden line removal technique,
3124 // with different linewidth and colour (red), with
3125 // solid colour, with shading of surfaces, and without edges.
3126 // Finally, some examples are given for the ray-tracing. (A possible
3127 // string for the NAME of the volume can be found using the command DTREE).
3133 if (fGcvdma->raytra != 1) {
3134 gdraw(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
3136 gdrayt(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
3140 //_____________________________________________________________________________
3141 void TGeant3::Gdrawc(const char *name,Int_t axis, Float_t cut,Float_t u0,
3142 Float_t v0,Float_t ul,Float_t vl)
3147 // CUTVAL Cut plane distance from the origin along the axis
3149 // U0 U-coord. (horizontal) of volume origin
3150 // V0 V-coord. (vertical) of volume origin
3151 // SU Scale factor for U-coord.
3152 // SV Scale factor for V-coord.
3154 // The cut plane is normal to caxis (X,Y,Z), corresponding to iaxis (1,2,3),
3155 // and placed at the distance cutval from the origin.
3156 // The resulting picture is seen from the the same axis.
3157 // When HIDE Mode is ON, it is possible to get the same effect with
3158 // the CVOL/BOX function.
3164 gdrawc(PASSCHARD(vname), axis,cut,u0,v0,ul,vl PASSCHARL(vname));
3167 //_____________________________________________________________________________
3168 void TGeant3::Gdrawx(const char *name,Float_t cutthe, Float_t cutphi,
3169 Float_t cutval, Float_t theta, Float_t phi, Float_t u0,
3170 Float_t v0,Float_t ul,Float_t vl)
3174 // CUTTHE Theta angle of the line normal to cut plane
3175 // CUTPHI Phi angle of the line normal to cut plane
3176 // CUTVAL Cut plane distance from the origin along the axis
3178 // THETA Viewing angle theta (for 3D projection)
3179 // PHI Viewing angle phi (for 3D projection)
3180 // U0 U-coord. (horizontal) of volume origin
3181 // V0 V-coord. (vertical) of volume origin
3182 // SU Scale factor for U-coord.
3183 // SV Scale factor for V-coord.
3185 // The cut plane is normal to the line given by the cut angles
3186 // cutthe and cutphi and placed at the distance cutval from the origin.
3187 // The resulting picture is seen from the viewing angles theta,phi.
3193 gdrawx(PASSCHARD(vname), cutthe,cutphi,cutval,theta,phi,u0,v0,ul,vl
3197 //_____________________________________________________________________________
3198 void TGeant3::Gdhead(Int_t isel, const char *name, Float_t chrsiz)
3203 // ISEL Option flag D=111110
3205 // CHRSIZ Character size (cm) of title NAME D=0.6
3208 // 0 to have only the header lines
3209 // xxxxx1 to add the text name centered on top of header
3210 // xxxx1x to add global detector name (first volume) on left
3211 // xxx1xx to add date on right
3212 // xx1xxx to select thick characters for text on top of header
3213 // x1xxxx to add the text 'EVENT NR x' on top of header
3214 // 1xxxxx to add the text 'RUN NR x' on top of header
3215 // NOTE that ISEL=x1xxx1 or ISEL=1xxxx1 are illegal choices,
3216 // i.e. they generate overwritten text.
3218 gdhead(isel,PASSCHARD(name),chrsiz PASSCHARL(name));
3221 //_____________________________________________________________________________
3222 void TGeant3::Gdman(Float_t u, Float_t v, const char *type)
3225 // Draw a 2D-man at position (U0,V0)
3227 // U U-coord. (horizontal) of the centre of man' R
3228 // V V-coord. (vertical) of the centre of man' R
3229 // TYPE D='MAN' possible values: 'MAN,WM1,WM2,WM3'
3231 // CALL GDMAN(u,v),CALL GDWMN1(u,v),CALL GDWMN2(u,v),CALL GDWMN2(u,v)
3232 // It superimposes the picure of a man or of a woman, chosen among
3233 // three different ones, with the same scale factors as the detector
3234 // in the current drawing.
3237 if (opt.Contains("WM1")) {
3239 } else if (opt.Contains("WM3")) {
3241 } else if (opt.Contains("WM2")) {
3248 //_____________________________________________________________________________
3249 void TGeant3::Gdspec(const char *name)
3254 // Shows 3 views of the volume (two cut-views and a 3D view), together with
3255 // its geometrical specifications. The 3D drawing will
3256 // be performed according the current values of the options HIDE and
3257 // SHAD and according the current SetClipBox clipping parameters for that
3264 gdspec(PASSCHARD(vname) PASSCHARL(vname));
3267 //_____________________________________________________________________________
3268 void TGeant3::DrawOneSpec(const char *name)
3271 // Function called when one double-clicks on a volume name
3272 // in a TPavelabel drawn by Gdtree.
3274 THIGZ *higzSave = gHigz;
3275 higzSave->SetName("higzSave");
3276 THIGZ *higzSpec = (THIGZ*)gROOT->FindObject("higzSpec");
3277 //printf("DrawOneSpec, gHigz=%x, higzSpec=%x\n",gHigz,higzSpec);
3278 if (higzSpec) gHigz = higzSpec;
3279 else higzSpec = new THIGZ(kDefSize);
3280 higzSpec->SetName("higzSpec");
3285 gdspec(PASSCHARD(vname) PASSCHARL(vname));
3288 higzSave->SetName("higz");
3292 //_____________________________________________________________________________
3293 void TGeant3::Gdtree(const char *name,Int_t levmax, Int_t isel)
3297 // LEVMAX Depth level
3300 // This function draws the logical tree,
3301 // Each volume in the tree is represented by a TPaveTree object.
3302 // Double-clicking on a TPaveTree draws the specs of the corresponding volume.
3303 // Use TPaveTree pop-up menu to select:
3306 // - drawing tree of parent
3312 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
3313 gHigz->SetPname("");
3316 //_____________________________________________________________________________
3317 void TGeant3::GdtreeParent(const char *name,Int_t levmax, Int_t isel)
3321 // LEVMAX Depth level
3324 // This function draws the logical tree of the parent of name.
3328 // Scan list of volumes in JVOLUM
3330 Int_t gname, i, jvo, in, nin, jin, num;
3331 strncpy((char *) &gname, name, 4);
3332 for(i=1; i<=fGcnum->nvolum; i++) {
3333 jvo = fZlq[fGclink->jvolum-i];
3334 nin = Int_t(fZq[jvo+3]);
3335 if (nin == -1) nin = 1;
3336 for (in=1;in<=nin;in++) {
3338 num = Int_t(fZq[jin+2]);
3339 if(gname == fZiq[fGclink->jvolum+num]) {
3340 strncpy(vname,(char*)&fZiq[fGclink->jvolum+i],4);
3342 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
3343 gHigz->SetPname("");
3350 //_____________________________________________________________________________
3351 void TGeant3::SetABAN(Int_t par)
3354 // par = 1 particles will be stopped according to their residual
3355 // range if they are not in a sensitive material and are
3356 // far enough from the boundary
3357 // 0 particles are transported normally
3359 fGcphys->dphys1 = par;
3363 //_____________________________________________________________________________
3364 void TGeant3::SetANNI(Int_t par)
3367 // To control positron annihilation.
3368 // par =0 no annihilation
3369 // =1 annihilation. Decays processed.
3370 // =2 annihilation. No decay products stored.
3372 fGcphys->ianni = par;
3376 //_____________________________________________________________________________
3377 void TGeant3::SetAUTO(Int_t par)
3380 // To control automatic calculation of tracking medium parameters:
3381 // par =0 no automatic calculation;
3382 // =1 automati calculation.
3384 fGctrak->igauto = par;
3388 //_____________________________________________________________________________
3389 void TGeant3::SetBOMB(Float_t boom)
3392 // BOOM : Exploding factor for volumes position
3394 // To 'explode' the detector. If BOOM is positive (values smaller
3395 // than 1. are suggested, but any value is possible)
3396 // all the volumes are shifted by a distance
3397 // proportional to BOOM along the direction between their centre
3398 // and the origin of the MARS; the volumes which are symmetric
3399 // with respect to this origin are simply not shown.
3400 // BOOM equal to 0 resets the normal mode.
3401 // A negative (greater than -1.) value of
3402 // BOOM will cause an 'implosion'; for even lower values of BOOM
3403 // the volumes' positions will be reflected respect to the origin.
3404 // This command can be useful to improve the 3D effect for very
3405 // complex detectors. The following commands will make explode the
3412 //_____________________________________________________________________________
3413 void TGeant3::SetBREM(Int_t par)
3416 // To control bremstrahlung.
3417 // par =0 no bremstrahlung
3418 // =1 bremstrahlung. Photon processed.
3419 // =2 bremstrahlung. No photon stored.
3421 fGcphys->ibrem = par;
3425 //_____________________________________________________________________________
3426 void TGeant3::SetCKOV(Int_t par)
3429 // To control Cerenkov production
3430 // par =0 no Cerenkov;
3432 // =2 Cerenkov with primary stopped at each step.
3434 fGctlit->itckov = par;
3438 //_____________________________________________________________________________
3439 void TGeant3::SetClipBox(const char *name,Float_t xmin,Float_t xmax,
3440 Float_t ymin,Float_t ymax,Float_t zmin,Float_t zmax)
3443 // The hidden line removal technique is necessary to visualize properly
3444 // very complex detectors. At the same time, it can be useful to visualize
3445 // the inner elements of a detector in detail. This function allows
3446 // subtractions (via boolean operation) of BOX shape from any part of
3447 // the detector, therefore showing its inner contents.
3448 // If "*" is given as the name of the
3449 // volume to be clipped, all volumes are clipped by the given box.
3450 // A volume can be clipped at most twice.
3451 // if a volume is explicitely clipped twice,
3452 // the "*" will not act on it anymore. Giving "." as the name
3453 // of the volume to be clipped will reset the clipping.
3455 // NAME Name of volume to be clipped
3457 // XMIN Lower limit of the Shape X coordinate
3458 // XMAX Upper limit of the Shape X coordinate
3459 // YMIN Lower limit of the Shape Y coordinate
3460 // YMAX Upper limit of the Shape Y coordinate
3461 // ZMIN Lower limit of the Shape Z coordinate
3462 // ZMAX Upper limit of the Shape Z coordinate
3464 // This function performs a boolean subtraction between the volume
3465 // NAME and a box placed in the MARS according the values of the given
3471 setclip(PASSCHARD(vname),xmin,xmax,ymin,ymax,zmin,zmax PASSCHARL(vname));
3474 //_____________________________________________________________________________
3475 void TGeant3::SetCOMP(Int_t par)
3478 // To control Compton scattering
3479 // par =0 no Compton
3480 // =1 Compton. Electron processed.
3481 // =2 Compton. No electron stored.
3484 fGcphys->icomp = par;
3487 //_____________________________________________________________________________
3488 void TGeant3::SetCUTS(Float_t cutgam,Float_t cutele,Float_t cutneu,
3489 Float_t cuthad,Float_t cutmuo ,Float_t bcute ,
3490 Float_t bcutm ,Float_t dcute ,Float_t dcutm ,
3491 Float_t ppcutm, Float_t tofmax)
3494 // CUTGAM Cut for gammas D=0.001
3495 // CUTELE Cut for electrons D=0.001
3496 // CUTHAD Cut for charged hadrons D=0.01
3497 // CUTNEU Cut for neutral hadrons D=0.01
3498 // CUTMUO Cut for muons D=0.01
3499 // BCUTE Cut for electron brems. D=-1.
3500 // BCUTM Cut for muon brems. D=-1.
3501 // DCUTE Cut for electron delta-rays D=-1.
3502 // DCUTM Cut for muon delta-rays D=-1.
3503 // PPCUTM Cut for e+e- pairs by muons D=0.01
3504 // TOFMAX Time of flight cut D=1.E+10
3506 // If the default values (-1.) for BCUTE ,BCUTM ,DCUTE ,DCUTM
3507 // are not modified, they will be set to CUTGAM,CUTGAM,CUTELE,CUTELE
3509 // If one of the parameters from CUTGAM to PPCUTM included
3510 // is modified, cross-sections and energy loss tables must be
3511 // recomputed via the function Gphysi.
3513 fGccuts->cutgam = cutgam;
3514 fGccuts->cutele = cutele;
3515 fGccuts->cutneu = cutneu;
3516 fGccuts->cuthad = cuthad;
3517 fGccuts->cutmuo = cutmuo;
3518 fGccuts->bcute = bcute;
3519 fGccuts->bcutm = bcutm;
3520 fGccuts->dcute = dcute;
3521 fGccuts->dcutm = dcutm;
3522 fGccuts->ppcutm = ppcutm;
3523 fGccuts->tofmax = tofmax;
3526 //_____________________________________________________________________________
3527 void TGeant3::SetDCAY(Int_t par)
3530 // To control Decay mechanism.
3531 // par =0 no decays.
3532 // =1 Decays. secondaries processed.
3533 // =2 Decays. No secondaries stored.
3535 fGcphys->idcay = par;
3539 //_____________________________________________________________________________
3540 void TGeant3::SetDEBU(Int_t emin, Int_t emax, Int_t emod)
3543 // Set the debug flag and frequency
3544 // Selected debug output will be printed from
3545 // event emin to even emax each emod event
3547 fGcflag->idemin = emin;
3548 fGcflag->idemax = emax;
3549 fGcflag->itest = emod;
3553 //_____________________________________________________________________________
3554 void TGeant3::SetDRAY(Int_t par)
3557 // To control delta rays mechanism.
3558 // par =0 no delta rays.
3559 // =1 Delta rays. secondaries processed.
3560 // =2 Delta rays. No secondaries stored.
3562 fGcphys->idray = par;
3565 //_____________________________________________________________________________
3566 void TGeant3::SetERAN(Float_t ekmin, Float_t ekmax, Int_t nekbin)
3569 // To control cross section tabulations
3570 // ekmin = minimum kinetic energy in GeV
3571 // ekmax = maximum kinetic energy in GeV
3572 // nekbin = number of logatithmic bins (<200)
3574 fGcmulo->ekmin = ekmin;
3575 fGcmulo->ekmax = ekmax;
3576 fGcmulo->nekbin = nekbin;
3579 //_____________________________________________________________________________
3580 void TGeant3::SetHADR(Int_t par)
3583 // To control hadronic interactions.
3584 // par =0 no hadronic interactions.
3585 // =1 Hadronic interactions. secondaries processed.
3586 // =2 Hadronic interactions. No secondaries stored.
3588 fGcphys->ihadr = par;
3591 //_____________________________________________________________________________
3592 void TGeant3::SetKINE(Int_t kine, Float_t xk1, Float_t xk2, Float_t xk3,
3593 Float_t xk4, Float_t xk5, Float_t xk6, Float_t xk7,
3594 Float_t xk8, Float_t xk9, Float_t xk10)
3597 // Set the variables in /GCFLAG/ IKINE, PKINE(10)
3598 // Their meaning is user defined
3600 fGckine->ikine = kine;
3601 fGckine->pkine[0] = xk1;
3602 fGckine->pkine[1] = xk2;
3603 fGckine->pkine[2] = xk3;
3604 fGckine->pkine[3] = xk4;
3605 fGckine->pkine[4] = xk5;
3606 fGckine->pkine[5] = xk6;
3607 fGckine->pkine[6] = xk7;
3608 fGckine->pkine[7] = xk8;
3609 fGckine->pkine[8] = xk9;
3610 fGckine->pkine[9] = xk10;
3613 //_____________________________________________________________________________
3614 void TGeant3::SetLOSS(Int_t par)
3617 // To control energy loss.
3618 // par =0 no energy loss;
3619 // =1 restricted energy loss fluctuations;
3620 // =2 complete energy loss fluctuations;
3622 // =4 no energy loss fluctuations.
3623 // If the value ILOSS is changed, then cross-sections and energy loss
3624 // tables must be recomputed via the command 'PHYSI'.
3626 fGcphys->iloss = par;
3630 //_____________________________________________________________________________
3631 void TGeant3::SetMULS(Int_t par)
3634 // To control multiple scattering.
3635 // par =0 no multiple scattering.
3636 // =1 Moliere or Coulomb scattering.
3637 // =2 Moliere or Coulomb scattering.
3638 // =3 Gaussian scattering.
3640 fGcphys->imuls = par;
3644 //_____________________________________________________________________________
3645 void TGeant3::SetMUNU(Int_t par)
3648 // To control muon nuclear interactions.
3649 // par =0 no muon-nuclear interactions.
3650 // =1 Nuclear interactions. Secondaries processed.
3651 // =2 Nuclear interactions. Secondaries not processed.
3653 fGcphys->imunu = par;
3656 //_____________________________________________________________________________
3657 void TGeant3::SetOPTI(Int_t par)
3660 // This flag controls the tracking optimisation performed via the
3662 // 1 no optimisation at all; GSORD calls disabled;
3663 // 0 no optimisation; only user calls to GSORD kept;
3664 // 1 all non-GSORDered volumes are ordered along the best axis;
3665 // 2 all volumes are ordered along the best axis.
3667 fGcopti->ioptim = par;
3670 //_____________________________________________________________________________
3671 void TGeant3::SetPAIR(Int_t par)
3674 // To control pair production mechanism.
3675 // par =0 no pair production.
3676 // =1 Pair production. secondaries processed.
3677 // =2 Pair production. No secondaries stored.
3679 fGcphys->ipair = par;
3683 //_____________________________________________________________________________
3684 void TGeant3::SetPFIS(Int_t par)
3687 // To control photo fission mechanism.
3688 // par =0 no photo fission.
3689 // =1 Photo fission. secondaries processed.
3690 // =2 Photo fission. No secondaries stored.
3692 fGcphys->ipfis = par;
3695 //_____________________________________________________________________________
3696 void TGeant3::SetPHOT(Int_t par)
3699 // To control Photo effect.
3700 // par =0 no photo electric effect.
3701 // =1 Photo effect. Electron processed.
3702 // =2 Photo effect. No electron stored.
3704 fGcphys->iphot = par;
3707 //_____________________________________________________________________________
3708 void TGeant3::SetRAYL(Int_t par)
3711 // To control Rayleigh scattering.
3712 // par =0 no Rayleigh scattering.
3715 fGcphys->irayl = par;
3718 //_____________________________________________________________________________
3719 void TGeant3::SetSTRA(Int_t par)
3722 // To control energy loss fluctuations
3723 // with the PhotoAbsorption Ionisation model.
3724 // par =0 no Straggling.
3725 // =1 Straggling yes => no Delta rays.
3727 fGcphlt->istra = par;
3730 //_____________________________________________________________________________
3731 void TGeant3::SetSWIT(Int_t sw, Int_t val)
3735 // val New switch value
3737 // Change one element of array ISWIT(10) in /GCFLAG/
3739 if (sw <= 0 || sw > 10) return;
3740 fGcflag->iswit[sw-1] = val;
3744 //_____________________________________________________________________________
3745 void TGeant3::SetTRIG(Int_t nevents)
3748 // Set number of events to be run
3750 fGcflag->nevent = nevents;
3753 //_____________________________________________________________________________
3754 void TGeant3::SetUserDecay(Int_t pdg)
3757 // Force the decays of particles to be done with Pythia
3758 // and not with the Geant routines.
3759 // just kill pointers doing mzdrop
3761 Int_t ipart = IdFromPDG(pdg);
3763 printf("Particle %d not in geant\n",pdg);
3766 Int_t jpart=fGclink->jpart;
3767 Int_t jpa=fZlq[jpart-ipart];
3770 Int_t jpa1=fZlq[jpa-1];
3772 mzdrop(fGcbank->ixcons,jpa1,PASSCHARD(" ") PASSCHARL(" "));
3773 Int_t jpa2=fZlq[jpa-2];
3775 mzdrop(fGcbank->ixcons,jpa2,PASSCHARD(" ") PASSCHARL(" "));
3779 //______________________________________________________________________________
3780 void TGeant3::Vname(const char *name, char *vname)
3783 // convert name to upper case. Make vname at least 4 chars
3785 Int_t l = strlen(name);
3788 for (i=0;i<l;i++) vname[i] = toupper(name[i]);
3789 for (i=l;i<4;i++) vname[i] = ' ';
3793 //______________________________________________________________________________
3794 void TGeant3::Ertrgo()
3797 // Perform the tracking of the track Track parameters are in VECT
3802 //______________________________________________________________________________
3803 void TGeant3::Ertrak(const Float_t *const x1, const Float_t *const p1,
3804 const Float_t *x2, const Float_t *p2,
3805 Int_t ipa, Option_t *chopt)
3807 //************************************************************************
3809 //* Perform the tracking of the track from point X1 to *
3811 //* (Before calling this routine the user should also provide *
3812 //* the input informations in /EROPTS/ and /ERTRIO/ *
3813 //* using subroutine EUFIL(L/P/V) *
3814 //* X1 - Starting coordinates (Cartesian) *
3815 //* P1 - Starting 3-momentum (Cartesian) *
3816 //* X2 - Final coordinates (Cartesian) *
3817 //* P2 - Final 3-momentum (Cartesian) *
3818 //* IPA - Particle code (a la GEANT) of the track *
3821 //* 'B' 'Backward tracking' - i.e. energy loss *
3822 //* added to the current energy *
3823 //* 'E' 'Exact' calculation of errors assuming *
3824 //* helix (i.e. pathlength not *
3825 //* assumed as infinitesimal) *
3826 //* 'L' Tracking upto prescribed Lengths reached *
3827 //* 'M' 'Mixed' prediction (not yet coded) *
3828 //* 'O' Tracking 'Only' without calculating errors *
3829 //* 'P' Tracking upto prescribed Planes reached *
3830 //* 'V' Tracking upto prescribed Volumes reached *
3831 //* 'X' Tracking upto prescribed Point approached *
3833 //* Interface with GEANT : *
3834 //* Track parameters are in /CGKINE/ and /GCTRAK/ *
3836 //* ==>Called by : USER *
3837 //* Authors M.Maire, E.Nagy ********//* *
3839 //************************************************************************
3840 ertrak(x1,p1,x2,p2,ipa,PASSCHARD(chopt) PASSCHARL(chopt));
3843 //_____________________________________________________________________________
3844 void TGeant3::WriteEuclid(const char* filnam, const char* topvol,
3845 Int_t number, Int_t nlevel)
3849 // ******************************************************************
3851 // * Write out the geometry of the detector in EUCLID file format *
3853 // * filnam : will be with the extension .euc *
3854 // * topvol : volume name of the starting node *
3855 // * number : copy number of topvol (relevant for gsposp) *
3856 // * nlevel : number of levels in the tree structure *
3857 // * to be written out, starting from topvol *
3859 // * Author : M. Maire *
3861 // ******************************************************************
3863 // File filnam.tme is written out with the definitions of tracking
3864 // medias and materials.
3865 // As to restore original numbers for materials and medias, program
3866 // searches in the file euc_medi.dat and comparing main parameters of
3867 // the mat. defined inside geant and the one in file recognizes them
3868 // and is able to take number from file. If for any material or medium,
3869 // this procedure fails, ordering starts from 1.
3870 // Arrays IOTMED and IOMATE are used for this procedure
3872 const char kShape[][5]={"BOX ","TRD1","TRD2","TRAP","TUBE","TUBS","CONE",
3873 "CONS","SPHE","PARA","PGON","PCON","ELTU","HYPE",
3875 Int_t i, end, itm, irm, jrm, k, nmed;
3879 char *filext, *filetme;
3880 char natmed[21], namate[21];
3881 char natmedc[21], namatec[21];
3882 char key[5], name[5], mother[5], konly[5];
3884 Int_t iadvol, iadtmd, iadrot, nwtot, iret;
3885 Int_t mlevel, numbr, natt, numed, nin, ndata;
3886 Int_t iname, ivo, ish, jvo, nvstak, ivstak;
3887 Int_t jdiv, ivin, in, jin, jvin, irot;
3888 Int_t jtm, imat, jma, flag=0, imatc;
3889 Float_t az, dens, radl, absl, a, step, x, y, z;
3890 Int_t npar, ndvmx, left;
3891 Float_t zc, densc, radlc, abslc, c0, tmaxfd;
3893 Int_t iomate[100], iotmed[100];
3894 Float_t par[100], att[20], ubuf[50];
3897 Int_t level, ndiv, iaxe;
3898 Int_t itmedc, nmatc, isvolc, ifieldc, nwbufc, isvol, nmat, ifield, nwbuf;
3899 Float_t fieldmc, tmaxfdc, stemaxc, deemaxc, epsilc, stminc, fieldm;
3900 Float_t tmaxf, stemax, deemax, epsil, stmin;
3901 const char *k10000="!\n%s\n!\n";
3902 //Open the input file
3904 for(i=0;i<end;i++) if(filnam[i]=='.') {
3908 filext=new char[end+5];
3909 filetme=new char[end+5];
3910 strncpy(filext,filnam,end);
3911 strncpy(filetme,filnam,end);
3913 // *** The output filnam name will be with extension '.euc'
3914 strcpy(&filext[end],".euc");
3915 strcpy(&filetme[end],".tme");
3916 lun=fopen(filext,"w");
3918 // *** Initialisation of the working space
3919 iadvol=fGcnum->nvolum;
3920 iadtmd=iadvol+fGcnum->nvolum;
3921 iadrot=iadtmd+fGcnum->ntmed;
3922 if(fGclink->jrotm) {
3923 fGcnum->nrotm=fZiq[fGclink->jrotm-2];
3927 nwtot=iadrot+fGcnum->nrotm;
3928 qws = new float[nwtot+1];
3929 for (i=0;i<nwtot+1;i++) qws[i]=0;
3932 if(nlevel==0) mlevel=20;
3934 // *** find the top volume and put it in the stak
3935 numbr = number>0 ? number : 1;
3936 Gfpara(topvol,numbr,1,npar,natt,par,att);
3938 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3943 // *** authorized shape ?
3944 strncpy((char *)&iname, topvol, 4);
3946 for(i=1; i<=fGcnum->nvolum; i++) if(fZiq[fGclink->jvolum+i]==iname) {
3950 jvo = fZlq[fGclink->jvolum-ivo];
3951 ish = Int_t (fZq[jvo+2]);
3953 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3960 iws[iadvol+ivo] = level;
3963 //*** flag all volumes and fill the stak
3967 // pick the next volume in stak
3969 ivo = TMath::Abs(iws[ivstak]);
3970 jvo = fZlq[fGclink->jvolum - ivo];
3972 // flag the tracking medium
3973 numed = Int_t (fZq[jvo + 4]);
3974 iws[iadtmd + numed] = 1;
3976 // get the daughters ...
3977 level = iws[iadvol+ivo];
3978 if (level < mlevel) {
3980 nin = Int_t (fZq[jvo + 3]);
3982 // from division ...
3984 jdiv = fZlq[jvo - 1];
3985 ivin = Int_t (fZq[jdiv + 2]);
3987 iws[nvstak] = -ivin;
3988 iws[iadvol+ivin] = level;
3990 // from position ...
3991 } else if (nin > 0) {
3992 for(in=1; in<=nin; in++) {
3993 jin = fZlq[jvo - in];
3994 ivin = Int_t (fZq[jin + 2 ]);
3995 jvin = fZlq[fGclink->jvolum - ivin];
3996 ish = Int_t (fZq[jvin + 2]);
3997 // authorized shape ?
3999 // not yet flagged ?
4000 if (iws[iadvol+ivin]==0) {
4003 iws[iadvol+ivin] = level;
4005 // flag the rotation matrix
4006 irot = Int_t ( fZq[jin + 4 ]);
4007 if (irot > 0) iws[iadrot+irot] = 1;
4013 // next volume in stak ?
4014 if (ivstak < nvstak) goto L10;
4016 // *** restore original material and media numbers
4017 // file euc_medi.dat is needed to compare materials and medias
4019 FILE* luncor=fopen("euc_medi.dat","r");
4022 for(itm=1; itm<=fGcnum->ntmed; itm++) {
4023 if (iws[iadtmd+itm] > 0) {
4024 jtm = fZlq[fGclink->jtmed-itm];
4025 strncpy(natmed,(char *)&fZiq[jtm+1],20);
4026 imat = Int_t (fZq[jtm+6]);
4027 jma = fZlq[fGclink->jmate-imat];
4029 printf(" *** GWEUCL *** material not defined for tracking medium %5i %s\n",itm,natmed);
4032 strncpy(namate,(char *)&fZiq[jma+1],20);
4035 //** find the material original number
4038 iret=fscanf(luncor,"%4s,%130s",key,card);
4039 if(iret<=0) goto L26;
4041 if(!strcmp(key,"MATE")) {
4042 sscanf(card,"%d %s %f %f %f %f %f %d",&imatc,namatec,&az,&zc,&densc,&radlc,&abslc,&nparc);
4043 Gfmate(imat,namate,a,z,dens,radl,absl,par,npar);
4044 if(!strcmp(namatec,namate)) {
4045 if(az==a && zc==z && densc==dens && radlc==radl
4046 && abslc==absl && nparc==nparc) {
4049 printf("*** GWEUCL *** material : %3d '%s' restored as %3d\n",imat,namate,imatc);
4051 printf("*** GWEUCL *** different definitions for material: %s\n",namate);
4055 if(strcmp(key,"END") && !flag) goto L23;
4057 printf("*** GWEUCL *** cannot restore original number for material: %s\n",namate);
4061 //*** restore original tracking medium number
4064 iret=fscanf(luncor,"%4s,%130s",key,card);
4065 if(iret<=0) goto L26;
4067 if (!strcmp(key,"TMED")) {
4068 sscanf(card,"%d %s %d %d %d %f %f %f %f %f %f %d\n",
4069 &itmedc,natmedc,&nmatc,&isvolc,&ifieldc,&fieldmc,
4070 &tmaxfdc,&stemaxc,&deemaxc,&epsilc,&stminc,&nwbufc);
4071 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxf,stemax,deemax,
4072 epsil,stmin,ubuf,&nwbuf);
4073 if(!strcmp(natmedc,natmed)) {
4074 if (iomate[nmat]==nmatc && nwbuf==nwbufc) {
4077 printf("*** GWEUCL *** medium : %3d '%20s' restored as %3d\n",
4080 printf("*** GWEUCL *** different definitions for tracking medium: %s\n",natmed);
4084 if(strcmp(key,"END") && !flag) goto L24;
4086 printf("cannot restore original number for medium : %s\n",natmed);
4094 L26: printf("*** GWEUCL *** cannot read the data file\n");
4096 L29: if(luncor) fclose (luncor);
4099 // *** write down the tracking medium definition
4101 strcpy(card,"! Tracking medium");
4102 fprintf(lun,k10000,card);
4104 for(itm=1;itm<=fGcnum->ntmed;itm++) {
4105 if (iws[iadtmd+itm]>0) {
4106 jtm = fZlq[fGclink->jtmed-itm];
4107 strncpy(natmed,(char *)&fZiq[jtm+1],20);
4109 imat = Int_t (fZq[jtm+6]);
4110 jma = fZlq[fGclink->jmate-imat];
4111 //* order media from one, if comparing with database failed
4113 iotmed[itm]=++imxtmed;
4114 iomate[imat]=++imxmate;
4119 printf(" *** GWEUCL *** material not defined for tracking medium %5d %s\n",
4122 strncpy(namate,(char *)&fZiq[jma+1],20);
4125 fprintf(lun,"TMED %3d '%20s' %3d '%20s'\n",iotmed[itm],natmed,iomate[imat],namate);
4129 //* *** write down the rotation matrix
4131 strcpy(card,"! Reperes");
4132 fprintf(lun,k10000,card);
4134 for(irm=1;irm<=fGcnum->nrotm;irm++) {
4135 if (iws[iadrot+irm]>0) {
4136 jrm = fZlq[fGclink->jrotm-irm];
4137 fprintf(lun,"ROTM %3d",irm);
4138 for(k=11;k<=16;k++) fprintf(lun," %8.3f",fZq[jrm+k]);
4143 //* *** write down the volume definition
4145 strcpy(card,"! Volumes");
4146 fprintf(lun,k10000,card);
4148 for(ivstak=1;ivstak<=nvstak;ivstak++) {
4151 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivo],4);
4153 jvo = fZlq[fGclink->jvolum-ivo];
4154 ish = Int_t (fZq[jvo+2]);
4155 nmed = Int_t (fZq[jvo+4]);
4156 npar = Int_t (fZq[jvo+5]);
4158 if (ivstak>1) for(i=0;i<npar;i++) par[i]=fZq[jvo+7+i];
4159 Gckpar (ish,npar,par);
4160 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,kShape[ish-1],iotmed[nmed],npar);
4161 for(i=0;i<(npar-1)/6+1;i++) {
4164 for(k=0;k<(left<6?left:6);k++) fprintf(lun," %11.5f",par[i*6+k]);
4168 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,kShape[ish-1],iotmed[nmed],npar);
4173 //* *** write down the division of volumes
4175 fprintf(lun,k10000,"! Divisions");
4176 for(ivstak=1;ivstak<=nvstak;ivstak++) {
4177 ivo = TMath::Abs(iws[ivstak]);
4178 jvo = fZlq[fGclink->jvolum-ivo];
4179 ish = Int_t (fZq[jvo+2]);
4180 nin = Int_t (fZq[jvo+3]);
4181 //* this volume is divided ...
4184 iaxe = Int_t ( fZq[jdiv+1]);
4185 ivin = Int_t ( fZq[jdiv+2]);
4186 ndiv = Int_t ( fZq[jdiv+3]);
4189 jvin = fZlq[fGclink->jvolum-ivin];
4190 nmed = Int_t ( fZq[jvin+4]);
4191 strncpy(mother,(char *)&fZiq[fGclink->jvolum+ivo ],4);
4193 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivin],4);
4195 if ((step<=0.)||(ish>=11)) {
4196 //* volume with negative parameter or gsposp or pgon ...
4197 fprintf(lun,"DIVN '%4s' '%4s' %3d %3d\n",name,mother,ndiv,iaxe);
4198 } else if ((ndiv<=0)||(ish==10)) {
4199 //* volume with negative parameter or gsposp or para ...
4200 ndvmx = TMath::Abs(ndiv);
4201 fprintf(lun,"DIVT '%4s' '%4s' %11.5f %3d %3d %3d\n",
4202 name,mother,step,iaxe,iotmed[nmed],ndvmx);
4204 //* normal volume : all kind of division are equivalent
4205 fprintf(lun,"DVT2 '%4s' '%4s' %11.5f %3d %11.5f %3d %3d\n",
4206 name,mother,step,iaxe,c0,iotmed[nmed],ndiv);
4211 //* *** write down the the positionnement of volumes
4213 fprintf(lun,k10000,"! Positionnements\n");
4215 for(ivstak = 1;ivstak<=nvstak;ivstak++) {
4216 ivo = TMath::Abs(iws[ivstak]);
4217 strncpy(mother,(char*)&fZiq[fGclink->jvolum+ivo ],4);
4219 jvo = fZlq[fGclink->jvolum-ivo];
4220 nin = Int_t( fZq[jvo+3]);
4221 //* this volume has daughters ...
4223 for (in=1;in<=nin;in++) {
4225 ivin = Int_t (fZq[jin +2]);
4226 numb = Int_t (fZq[jin +3]);
4227 irot = Int_t (fZq[jin +4]);
4231 strcpy(konly,"ONLY");
4232 if (fZq[jin+8]!=1.) strcpy(konly,"MANY");
4233 strncpy(name,(char*)&fZiq[fGclink->jvolum+ivin],4);
4235 jvin = fZlq[fGclink->jvolum-ivin];
4236 ish = Int_t (fZq[jvin+2]);
4237 //* gspos or gsposp ?
4238 ndata = fZiq[jin-1];
4240 fprintf(lun,"POSI '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s'\n",
4241 name,numb,mother,x,y,z,irot,konly);
4243 npar = Int_t (fZq[jin+9]);
4244 for(i=0;i<npar;i++) par[i]=fZq[jin+10+i];
4245 Gckpar (ish,npar,par);
4246 fprintf(lun,"POSP '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s' %3d\n",
4247 name,numb,mother,x,y,z,irot,konly,npar);
4249 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
4256 fprintf(lun,"END\n");
4259 //****** write down the materials and medias *****
4261 lun=fopen(filetme,"w");
4263 for(itm=1;itm<=fGcnum->ntmed;itm++) {
4264 if (iws[iadtmd+itm]>0) {
4265 jtm = fZlq[fGclink->jtmed-itm];
4266 strncpy(natmed,(char*)&fZiq[jtm+1],4);
4267 imat = Int_t (fZq[jtm+6]);
4268 jma = Int_t (fZlq[fGclink->jmate-imat]);
4270 Gfmate (imat,namate,a,z,dens,radl,absl,par,npar);
4271 fprintf(lun,"MATE %4d '%20s'%11.5E %11.5E %11.5E %11.5E %11.5E %3d\n",
4272 iomate[imat],namate,a,z,dens,radl,absl,npar);
4276 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
4280 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxfd,stemax,deemax,epsil,stmin,par,&npar);
4281 fprintf(lun,"TMED %4d '%20s' %3d %1d %3d %11.5f %11.5f %11.5f %11.5f %11.5f %11.5f %3d\n",
4282 iotmed[itm],natmed,iomate[nmat],isvol,ifield,
4283 fieldm,tmaxfd,stemax,deemax,epsil,stmin,npar);
4287 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
4293 fprintf(lun,"END\n");
4295 printf(" *** GWEUCL *** file: %s is now written out\n",filext);
4296 printf(" *** GWEUCL *** file: %s is now written out\n",filetme);
4305 //_____________________________________________________________________________
4306 void TGeant3::Streamer(TBuffer &R__b)
4309 // Stream an object of class TGeant3.
4311 if (R__b.IsReading()) {
4312 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
4313 AliMC::Streamer(R__b);
4316 R__b.ReadStaticArray(fPDGCode);
4318 R__b.WriteVersion(TGeant3::IsA());
4319 AliMC::Streamer(R__b);
4322 R__b.WriteArray(fPDGCode, fNPDGCodes);