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.14 1999/09/29 09:24:30 fca
19 Introduction of the Copyright and cvs Log
23 ///////////////////////////////////////////////////////////////////////////////
25 // Interface Class to the Geant3.21 MonteCarlo //
29 <img src="picts/TGeant3Class.gif">
34 ///////////////////////////////////////////////////////////////////////////////
40 #include <TDatabasePDG.h>
41 #include "AliCallf77.h"
44 # define gzebra gzebra_
45 # define grfile grfile_
46 # define gpcxyz gpcxyz_
47 # define ggclos ggclos_
50 # define gcinit gcinit_
53 # define gtrigc gtrigc_
54 # define gtrigi gtrigi_
56 # define gzinit gzinit_
57 # define gfmate gfmate_
58 # define gfpart gfpart_
59 # define gftmed gftmed_
63 # define gsmate gsmate_
64 # define gsmixt gsmixt_
65 # define gspart gspart_
66 # define gstmed gstmed_
67 # define gsckov gsckov_
68 # define gstpar gstpar_
69 # define gfkine gfkine_
70 # define gfvert gfvert_
71 # define gskine gskine_
72 # define gsvert gsvert_
73 # define gphysi gphysi_
74 # define gdebug gdebug_
75 # define gekbin gekbin_
76 # define gfinds gfinds_
77 # define gsking gsking_
78 # define gskpho gskpho_
79 # define gsstak gsstak_
81 # define gtrack gtrack_
82 # define gtreve gtreve_
83 # define gtreve_root gtreve_root_
85 # define grndmq grndmq_
87 # define glmoth glmoth_
88 # define gmedia gmedia_
91 # define gsdvn2 gsdvn2_
93 # define gsdvs2 gsdvs2_
95 # define gsdvt2 gsdvt2_
98 # define gsposp gsposp_
99 # define gsrotm gsrotm_
100 # define gprotm gprotm_
101 # define gsvolu gsvolu_
102 # define gprint gprint_
103 # define gdinit gdinit_
104 # define gdopt gdopt_
105 # define gdraw gdraw_
106 # define gdrayt gdrayt_
107 # define gdrawc gdrawc_
108 # define gdrawx gdrawx_
109 # define gdhead gdhead_
110 # define gdwmn1 gdwmn1_
111 # define gdwmn2 gdwmn2_
112 # define gdwmn3 gdwmn3_
113 # define gdxyz gdxyz_
114 # define gdcxyz gdcxyz_
115 # define gdman gdman_
116 # define gdspec gdspec_
117 # define gdtree gdtree_
118 # define gdelet gdelet_
119 # define gdclos gdclos_
120 # define gdshow gdshow_
121 # define gdopen gdopen_
122 # define dzshow dzshow_
123 # define gsatt gsatt_
124 # define gfpara gfpara_
125 # define gckpar gckpar_
126 # define gckmat gckmat_
127 # define geditv geditv_
128 # define mzdrop mzdrop_
130 # define ertrak ertrak_
131 # define ertrgo ertrgo_
133 # define setbomb setbomb_
134 # define setclip setclip_
135 # define gcomad gcomad_
138 # define gzebra GZEBRA
139 # define grfile GRFILE
140 # define gpcxyz GPCXYZ
141 # define ggclos GGCLOS
144 # define gcinit GCINIT
147 # define gtrigc GTRIGC
148 # define gtrigi GTRIGI
150 # define gzinit GZINIT
151 # define gfmate GFMATE
152 # define gfpart GFPART
153 # define gftmed GFTMED
157 # define gsmate GSMATE
158 # define gsmixt GSMIXT
159 # define gspart GSPART
160 # define gstmed GSTMED
161 # define gsckov GSCKOV
162 # define gstpar GSTPAR
163 # define gfkine GFKINE
164 # define gfvert GFVERT
165 # define gskine GSKINE
166 # define gsvert GSVERT
167 # define gphysi GPHYSI
168 # define gdebug GDEBUG
169 # define gekbin GEKBIN
170 # define gfinds GFINDS
171 # define gsking GSKING
172 # define gskpho GSKPHO
173 # define gsstak GSSTAK
175 # define gtrack GTRACK
176 # define gtreve GTREVE
177 # define gtreve_root GTREVE_ROOT
179 # define grndmq GRNDMQ
181 # define glmoth GLMOTH
182 # define gmedia GMEDIA
185 # define gsdvn2 GSDVN2
187 # define gsdvs2 GSDVS2
189 # define gsdvt2 GSDVT2
192 # define gsposp GSPOSP
193 # define gsrotm GSROTM
194 # define gprotm GPROTM
195 # define gsvolu GSVOLU
196 # define gprint GPRINT
197 # define gdinit GDINIT
200 # define gdrayt GDRAYT
201 # define gdrawc GDRAWC
202 # define gdrawx GDRAWX
203 # define gdhead GDHEAD
204 # define gdwmn1 GDWMN1
205 # define gdwmn2 GDWMN2
206 # define gdwmn3 GDWMN3
208 # define gdcxyz GDCXYZ
210 # define gdfspc GDFSPC
211 # define gdspec GDSPEC
212 # define gdtree GDTREE
213 # define gdelet GDELET
214 # define gdclos GDCLOS
215 # define gdshow GDSHOW
216 # define gdopen GDOPEN
217 # define dzshow DZSHOW
219 # define gfpara GFPARA
220 # define gckpar GCKPAR
221 # define gckmat GCKMAT
222 # define geditv GEDITV
223 # define mzdrop MZDROP
225 # define ertrak ERTRAK
226 # define ertrgo ERTRGO
228 # define setbomb SETBOMB
229 # define setclip SETCLIP
230 # define gcomad GCOMAD
234 //____________________________________________________________________________
238 // Prototypes for GEANT functions
240 void type_of_call gzebra(const int&);
242 void type_of_call gpcxyz();
244 void type_of_call ggclos();
246 void type_of_call glast();
248 void type_of_call ginit();
250 void type_of_call gcinit();
252 void type_of_call grun();
254 void type_of_call gtrig();
256 void type_of_call gtrigc();
258 void type_of_call gtrigi();
260 void type_of_call gwork(const int&);
262 void type_of_call gzinit();
264 void type_of_call gmate();
266 void type_of_call gpart();
268 void type_of_call gsdk(Int_t &, Float_t *, Int_t *);
270 void type_of_call gfkine(Int_t &, Float_t *, Float_t *, Int_t &,
271 Int_t &, Float_t *, Int_t &);
273 void type_of_call gfvert(Int_t &, Float_t *, Int_t &, Int_t &,
274 Float_t &, Float_t *, Int_t &);
276 void type_of_call gskine(Float_t *,Int_t &, Int_t &, Float_t *,
279 void type_of_call gsvert(Float_t *,Int_t &, Int_t &, Float_t *,
282 void type_of_call gphysi();
284 void type_of_call gdebug();
286 void type_of_call gekbin();
288 void type_of_call gfinds();
290 void type_of_call gsking(Int_t &);
292 void type_of_call gskpho(Int_t &);
294 void type_of_call gsstak(Int_t &);
296 void type_of_call gsxyz();
298 void type_of_call gtrack();
300 void type_of_call gtreve();
302 void type_of_call gtreve_root();
304 void type_of_call grndm(Float_t *, const Int_t &);
306 void type_of_call grndmq(Int_t &, Int_t &, const Int_t &,
309 void type_of_call gdtom(Float_t *, Float_t *, Int_t &);
311 void type_of_call glmoth(DEFCHARD, Int_t &, Int_t &, Int_t *,
312 Int_t *, Int_t * DEFCHARL);
314 void type_of_call gmedia(Float_t *, Int_t &);
316 void type_of_call gmtod(Float_t *, Float_t *, Int_t &);
318 void type_of_call gsrotm(const Int_t &, const Float_t &, const Float_t &,
319 const Float_t &, const Float_t &, const Float_t &,
322 void type_of_call gprotm(const Int_t &);
324 void type_of_call grfile(const Int_t&, DEFCHARD,
325 DEFCHARD DEFCHARL DEFCHARL);
327 void type_of_call gfmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
328 Float_t &, Float_t &, Float_t &, Float_t *,
331 void type_of_call gfpart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
332 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
334 void type_of_call gftmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
335 Float_t &, Float_t &, Float_t &, Float_t &,
336 Float_t &, Float_t &, Float_t *, Int_t * DEFCHARL);
338 void type_of_call gsmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
339 Float_t &, Float_t &, Float_t &, Float_t *,
342 void type_of_call gsmixt(const Int_t&, DEFCHARD, Float_t *, Float_t *,
343 Float_t &, Int_t &, Float_t * DEFCHARL);
345 void type_of_call gspart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
346 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
349 void type_of_call gstmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
350 Float_t &, Float_t &, Float_t &, Float_t &,
351 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
353 void type_of_call gsckov(Int_t &itmed, Int_t &npckov, Float_t *ppckov,
354 Float_t *absco, Float_t *effic, Float_t *rindex);
355 void type_of_call gstpar(const Int_t&, DEFCHARD, Float_t & DEFCHARL);
357 void type_of_call gsdvn(DEFCHARD,DEFCHARD, Int_t &, Int_t &
360 void type_of_call gsdvn2(DEFCHARD,DEFCHARD, Int_t &, Int_t &, Float_t &,
361 Int_t & DEFCHARL DEFCHARL);
363 void type_of_call gsdvs(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &
366 void type_of_call gsdvs2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t &,
367 Int_t & DEFCHARL DEFCHARL);
369 void type_of_call gsdvt(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &,
370 Int_t & DEFCHARL DEFCHARL);
372 void type_of_call gsdvt2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t&,
373 Int_t &, Int_t & DEFCHARL DEFCHARL);
375 void type_of_call gsord(DEFCHARD, Int_t & DEFCHARL);
377 void type_of_call gspos(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
378 Float_t &, Int_t &, DEFCHARD DEFCHARL DEFCHARL
381 void type_of_call gsposp(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
382 Float_t &, Int_t &, DEFCHARD,
383 Float_t *, Int_t & DEFCHARL DEFCHARL DEFCHARL);
385 void type_of_call gsvolu(DEFCHARD, DEFCHARD, Int_t &, Float_t *, Int_t &,
386 Int_t & DEFCHARL DEFCHARL);
388 void type_of_call gsatt(DEFCHARD, DEFCHARD, Int_t & DEFCHARL DEFCHARL);
390 void type_of_call gfpara(DEFCHARD , Int_t&, Int_t&, Int_t&, Int_t&, Float_t*,
393 void type_of_call gckpar(Int_t&, Int_t&, Float_t*);
395 void type_of_call gckmat(Int_t&, DEFCHARD DEFCHARL);
397 void type_of_call gprint(DEFCHARD,const int& DEFCHARL);
399 void type_of_call gdinit();
401 void type_of_call gdopt(DEFCHARD,DEFCHARD DEFCHARL DEFCHARL);
403 void type_of_call gdraw(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
404 Float_t &, Float_t &, Float_t & DEFCHARL);
405 void type_of_call gdrayt(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
406 Float_t &, Float_t &, Float_t & DEFCHARL);
407 void type_of_call gdrawc(DEFCHARD,Int_t &, Float_t &, Float_t &, Float_t &,
408 Float_t &, Float_t & DEFCHARL);
409 void type_of_call gdrawx(DEFCHARD,Float_t &, Float_t &, Float_t &, Float_t &,
410 Float_t &, Float_t &, Float_t &, Float_t &,
412 void type_of_call gdhead(Int_t &,DEFCHARD, Float_t & DEFCHARL);
413 void type_of_call gdxyz(Int_t &);
414 void type_of_call gdcxyz();
415 void type_of_call gdman(Float_t &, Float_t &);
416 void type_of_call gdwmn1(Float_t &, Float_t &);
417 void type_of_call gdwmn2(Float_t &, Float_t &);
418 void type_of_call gdwmn3(Float_t &, Float_t &);
419 void type_of_call gdspec(DEFCHARD DEFCHARL);
420 void type_of_call gdfspc(DEFCHARD, Int_t &, Int_t & DEFCHARL) {;}
421 void type_of_call gdtree(DEFCHARD, Int_t &, Int_t & DEFCHARL);
423 void type_of_call gdopen(Int_t &);
424 void type_of_call gdclos();
425 void type_of_call gdelet(Int_t &);
426 void type_of_call gdshow(Int_t &);
427 void type_of_call geditv(Int_t &) {;}
430 void type_of_call dzshow(DEFCHARD,const int&,const int&,DEFCHARD,const int&,
431 const int&, const int&, const int& DEFCHARL
434 void type_of_call mzdrop(Int_t&, Int_t&, DEFCHARD DEFCHARL);
436 void type_of_call setbomb(Float_t &);
437 void type_of_call setclip(DEFCHARD, Float_t &,Float_t &,Float_t &,Float_t &,
438 Float_t &, Float_t & DEFCHARL);
439 void type_of_call gcomad(DEFCHARD, Int_t*& DEFCHARL);
441 void type_of_call ertrak(const Float_t *const x1, const Float_t *const p1,
442 const Float_t *x2, const Float_t *p2,
443 const Int_t &ipa, DEFCHARD DEFCHARL);
445 void type_of_call ertrgo();
449 // Geant3 global pointer
451 static Int_t defSize = 600;
455 //____________________________________________________________________________
459 // Default constructor
463 //____________________________________________________________________________
464 TGeant3::TGeant3(const char *title, Int_t nwgeant)
465 :AliMC("TGeant3",title)
468 // Standard constructor for TGeant3 with ZEBRA initialisation
479 // Load Address of Geant3 commons
482 // Zero number of particles
486 //____________________________________________________________________________
487 Int_t TGeant3::CurrentMaterial(Float_t &a, Float_t &z, Float_t &dens,
488 Float_t &radl, Float_t &absl) const
491 // Return the parameters of the current material during transport
495 dens = fGcmate->dens;
496 radl = fGcmate->radl;
497 absl = fGcmate->absl;
498 return 1; //this could be the number of elements in mixture
501 //____________________________________________________________________________
502 void TGeant3::DefaultRange()
505 // Set range of current drawing pad to 20x20 cm
511 higz->Range(0,0,20,20);
514 //____________________________________________________________________________
515 void TGeant3::InitHIGZ()
526 //____________________________________________________________________________
527 void TGeant3::LoadAddress()
530 // Assigns the address of the GEANT common blocks to the structures
531 // that allow their access from C++
534 gcomad(PASSCHARD("QUEST"), (int*&) fQuest PASSCHARL("QUEST"));
535 gcomad(PASSCHARD("GCBANK"),(int*&) fGcbank PASSCHARL("GCBANK"));
536 gcomad(PASSCHARD("GCLINK"),(int*&) fGclink PASSCHARL("GCLINK"));
537 gcomad(PASSCHARD("GCCUTS"),(int*&) fGccuts PASSCHARL("GCCUTS"));
538 gcomad(PASSCHARD("GCFLAG"),(int*&) fGcflag PASSCHARL("GCFLAG"));
539 gcomad(PASSCHARD("GCKINE"),(int*&) fGckine PASSCHARL("GCKINE"));
540 gcomad(PASSCHARD("GCKING"),(int*&) fGcking PASSCHARL("GCKING"));
541 gcomad(PASSCHARD("GCKIN2"),(int*&) fGckin2 PASSCHARL("GCKIN2"));
542 gcomad(PASSCHARD("GCKIN3"),(int*&) fGckin3 PASSCHARL("GCKIN3"));
543 gcomad(PASSCHARD("GCMATE"),(int*&) fGcmate PASSCHARL("GCMATE"));
544 gcomad(PASSCHARD("GCTMED"),(int*&) fGctmed PASSCHARL("GCTMED"));
545 gcomad(PASSCHARD("GCTRAK"),(int*&) fGctrak PASSCHARL("GCTRAK"));
546 gcomad(PASSCHARD("GCTPOL"),(int*&) fGctpol PASSCHARL("GCTPOL"));
547 gcomad(PASSCHARD("GCVOLU"),(int*&) fGcvolu PASSCHARL("GCVOLU"));
548 gcomad(PASSCHARD("GCNUM"), (int*&) fGcnum PASSCHARL("GCNUM"));
549 gcomad(PASSCHARD("GCSETS"),(int*&) fGcsets PASSCHARL("GCSETS"));
550 gcomad(PASSCHARD("GCPHYS"),(int*&) fGcphys PASSCHARL("GCPHYS"));
551 gcomad(PASSCHARD("GCOPTI"),(int*&) fGcopti PASSCHARL("GCOPTI"));
552 gcomad(PASSCHARD("GCTLIT"),(int*&) fGctlit PASSCHARL("GCTLIT"));
553 gcomad(PASSCHARD("GCVDMA"),(int*&) fGcvdma PASSCHARL("GCVDMA"));
556 gcomad(PASSCHARD("ERTRIO"),(int*&) fErtrio PASSCHARL("ERTRIO"));
557 gcomad(PASSCHARD("EROPTS"),(int*&) fEropts PASSCHARL("EROPTS"));
558 gcomad(PASSCHARD("EROPTC"),(int*&) fEroptc PASSCHARL("EROPTC"));
559 gcomad(PASSCHARD("ERWORK"),(int*&) fErwork PASSCHARL("ERWORK"));
561 // Variables for ZEBRA store
562 gcomad(PASSCHARD("IQ"), addr PASSCHARL("IQ"));
564 gcomad(PASSCHARD("LQ"), addr PASSCHARL("LQ"));
569 //_____________________________________________________________________________
570 void TGeant3::GeomIter()
573 // Geometry iterator for moving upward in the geometry tree
574 // Initialise the iterator
576 fNextVol=fGcvolu->nlevel;
579 //____________________________________________________________________________
580 Int_t TGeant3::NextVolUp(Text_t *name, Int_t ©)
583 // Geometry iterator for moving upward in the geometry tree
584 // Return next volume up
589 gname=fGcvolu->names[fNextVol];
590 strncpy(name,(char *) &gname, 4);
592 copy=fGcvolu->number[fNextVol];
593 i=fGcvolu->lvolum[fNextVol];
594 if(gname == fZiq[fGclink->jvolum+i]) return i;
595 else printf("GeomTree: Volume %s not found in bank\n",name);
600 //_____________________________________________________________________________
601 Int_t TGeant3::CurrentVolID(Int_t ©) const
604 // Returns the current volume ID and copy number
607 if( (i=fGcvolu->nlevel-1) < 0 ) {
608 Warning("CurrentVolID","Stack depth only %d\n",fGcvolu->nlevel);
610 gname=fGcvolu->names[i];
611 copy=fGcvolu->number[i];
612 i=fGcvolu->lvolum[i];
613 if(gname == fZiq[fGclink->jvolum+i]) return i;
614 else Warning("CurrentVolID","Volume %4s not found\n",(char*)&gname);
619 //_____________________________________________________________________________
620 Int_t TGeant3::CurrentVolOffID(Int_t off, Int_t ©) const
623 // Return the current volume "off" upward in the geometrical tree
624 // ID and copy number
627 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
628 Warning("CurrentVolOffID","Offset requested %d but stack depth %d\n",
629 off,fGcvolu->nlevel);
631 gname=fGcvolu->names[i];
632 copy=fGcvolu->number[i];
633 i=fGcvolu->lvolum[i];
634 if(gname == fZiq[fGclink->jvolum+i]) return i;
635 else Warning("CurrentVolOffID","Volume %4s not found\n",(char*)&gname);
640 //_____________________________________________________________________________
641 const char* TGeant3::CurrentVolName() const
644 // Returns the current volume name
648 if( (i=fGcvolu->nlevel-1) < 0 ) {
649 Warning("CurrentVolName","Stack depth %d\n",fGcvolu->nlevel);
651 gname=fGcvolu->names[i];
653 strncpy(name,(char *) &gname, 4);
655 i=fGcvolu->lvolum[i];
656 if(gname == fZiq[fGclink->jvolum+i]) return name;
657 else Warning("CurrentVolName","Volume %4s not found\n",name);
662 //_____________________________________________________________________________
663 const char* TGeant3::CurrentVolOffName(Int_t off) const
666 // Return the current volume "off" upward in the geometrical tree
667 // ID, name and copy number
668 // if name=0 no name is returned
672 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
673 Warning("CurrentVolOffName",
674 "Offset requested %d but stack depth %d\n",off,fGcvolu->nlevel);
676 gname=fGcvolu->names[i];
678 strncpy(name,(char *) &gname, 4);
680 i=fGcvolu->lvolum[i];
681 if(gname == fZiq[fGclink->jvolum+i]) return name;
682 else Warning("CurrentVolOffName","Volume %4s not found\n",name);
687 //_____________________________________________________________________________
688 Int_t TGeant3::IdFromPDG(Int_t pdg) const
691 // Return Geant3 code from PDG and pseudo ENDF code
693 for(Int_t i=0;i<fNPDGCodes;++i)
694 if(pdg==fPDGCode[i]) return i;
698 //_____________________________________________________________________________
699 Int_t TGeant3::PDGFromId(Int_t id) const
701 if(id>0 && id<fNPDGCodes) return fPDGCode[id];
705 //_____________________________________________________________________________
706 void TGeant3::DefineParticles()
709 // Define standard Geant 3 particles
712 // Load standard numbers for GEANT particles and PDG conversion
713 fPDGCode[fNPDGCodes++]=-99; // 0 = unused location
714 fPDGCode[fNPDGCodes++]=22; // 1 = photon
715 fPDGCode[fNPDGCodes++]=-11; // 2 = positron
716 fPDGCode[fNPDGCodes++]=11; // 3 = electron
717 fPDGCode[fNPDGCodes++]=12; // 4 = neutrino e
718 fPDGCode[fNPDGCodes++]=-13; // 5 = muon +
719 fPDGCode[fNPDGCodes++]=13; // 6 = muon -
720 fPDGCode[fNPDGCodes++]=111; // 7 = pi0
721 fPDGCode[fNPDGCodes++]=211; // 8 = pi+
722 fPDGCode[fNPDGCodes++]=-211; // 9 = pi-
723 fPDGCode[fNPDGCodes++]=130; // 10 = Kaon Long
724 fPDGCode[fNPDGCodes++]=321; // 11 = Kaon +
725 fPDGCode[fNPDGCodes++]=-321; // 12 = Kaon -
726 fPDGCode[fNPDGCodes++]=2112; // 13 = Neutron
727 fPDGCode[fNPDGCodes++]=2212; // 14 = Proton
728 fPDGCode[fNPDGCodes++]=-2212; // 15 = Anti Proton
729 fPDGCode[fNPDGCodes++]=310; // 16 = Kaon Short
730 fPDGCode[fNPDGCodes++]=221; // 17 = Eta
731 fPDGCode[fNPDGCodes++]=3122; // 18 = Lambda
732 fPDGCode[fNPDGCodes++]=3222; // 19 = Sigma +
733 fPDGCode[fNPDGCodes++]=3212; // 20 = Sigma 0
734 fPDGCode[fNPDGCodes++]=3112; // 21 = Sigma -
735 fPDGCode[fNPDGCodes++]=3322; // 22 = Xi0
736 fPDGCode[fNPDGCodes++]=3312; // 23 = Xi-
737 fPDGCode[fNPDGCodes++]=3334; // 24 = Omega-
738 fPDGCode[fNPDGCodes++]=-2112; // 25 = Anti Proton
739 fPDGCode[fNPDGCodes++]=-3122; // 26 = Anti Proton
740 fPDGCode[fNPDGCodes++]=-3222; // 27 = Anti Sigma -
741 fPDGCode[fNPDGCodes++]=-3212; // 28 = Anti Sigma 0
742 fPDGCode[fNPDGCodes++]=-3112; // 29 = Anti Sigma 0
743 fPDGCode[fNPDGCodes++]=-3322; // 30 = Anti Xi 0
744 fPDGCode[fNPDGCodes++]=-3312; // 31 = Anti Xi +
745 fPDGCode[fNPDGCodes++]=-3334; // 32 = Anti Omega +
752 /* --- Define additional particles */
753 Gspart(33, "OMEGA(782)", 3, 0.782, 0., 7.836e-23);
754 fPDGCode[fNPDGCodes++]=223; // 33 = Omega(782)
756 Gspart(34, "PHI(1020)", 3, 1.019, 0., 1.486e-22);
757 fPDGCode[fNPDGCodes++]=333; // 34 = PHI (1020)
759 Gspart(35, "D +", 4, 1.87, 1., 1.066e-12);
760 fPDGCode[fNPDGCodes++]=411; // 35 = D+
762 Gspart(36, "D -", 4, 1.87, -1., 1.066e-12);
763 fPDGCode[fNPDGCodes++]=-411; // 36 = D-
765 Gspart(37, "D 0", 3, 1.865, 0., 4.2e-13);
766 fPDGCode[fNPDGCodes++]=421; // 37 = D0
768 Gspart(38, "ANTI D 0", 3, 1.865, 0., 4.2e-13);
769 fPDGCode[fNPDGCodes++]=-421; // 38 = D0 bar
771 fPDGCode[fNPDGCodes++]=-99; // 39 = unassigned
773 fPDGCode[fNPDGCodes++]=-99; // 40 = unassigned
775 fPDGCode[fNPDGCodes++]=-99; // 41 = unassigned
777 Gspart(42, "RHO +", 4, 0.768, 1., 4.353e-24);
778 fPDGCode[fNPDGCodes++]=213; // 42 = RHO+
780 Gspart(43, "RHO -", 4, 0.768, -1., 4.353e-24);
781 fPDGCode[fNPDGCodes++]=-213; // 40 = RHO-
783 Gspart(44, "RHO 0", 3, 0.768, 0., 4.353e-24);
784 fPDGCode[fNPDGCodes++]=113; // 37 = D0
787 // Use ENDF-6 mapping for ions = 10000*z+10*a+iso
789 // and numbers above 5 000 000 for special applications
792 const Int_t kion=10000000;
794 const Int_t kspe=50000000;
796 TDatabasePDG *pdgDB = TDatabasePDG::Instance();
798 const Double_t autogev=0.9314943228;
799 const Double_t hslash = 1.0545726663e-27;
800 const Double_t erggev = 1/1.6021773349e-3;
801 const Double_t hshgev = hslash*erggev;
802 const Double_t yearstosec = 3600*24*365.25;
805 pdgDB->AddParticle("Deuteron","Deuteron",2*autogev+8.071e-3,kTRUE,
806 0,1,"Ion",kion+10020);
807 fPDGCode[fNPDGCodes++]=kion+10020; // 45 = Deuteron
809 pdgDB->AddParticle("Triton","Triton",3*autogev+14.931e-3,kFALSE,
810 hshgev/(12.33*yearstosec),1,"Ion",kion+10030);
811 fPDGCode[fNPDGCodes++]=kion+10030; // 46 = Triton
813 pdgDB->AddParticle("Alpha","Alpha",4*autogev+2.424e-3,kTRUE,
814 hshgev/(12.33*yearstosec),2,"Ion",kion+20040);
815 fPDGCode[fNPDGCodes++]=kion+20040; // 47 = Alpha
817 fPDGCode[fNPDGCodes++]=0; // 48 = geantino mapped to rootino
819 pdgDB->AddParticle("HE3","HE3",3*autogev+14.931e-3,kFALSE,
820 0,2,"Ion",kion+20030);
821 fPDGCode[fNPDGCodes++]=kion+20030; // 49 = HE3
823 pdgDB->AddParticle("Cherenkov","Cherenkov",0,kFALSE,
824 0,0,"Special",kspe+50);
825 fPDGCode[fNPDGCodes++]=kspe+50; // 50 = Cherenkov
827 /* --- Define additional decay modes --- */
828 /* --- omega(783) --- */
829 for (kz = 0; kz < 6; ++kz) {
840 Gsdk(ipa, bratio, mode);
841 /* --- phi(1020) --- */
842 for (kz = 0; kz < 6; ++kz) {
857 Gsdk(ipa, bratio, mode);
859 for (kz = 0; kz < 6; ++kz) {
872 Gsdk(ipa, bratio, mode);
874 for (kz = 0; kz < 6; ++kz) {
887 Gsdk(ipa, bratio, mode);
889 for (kz = 0; kz < 6; ++kz) {
900 Gsdk(ipa, bratio, mode);
901 /* --- Anti D0 --- */
902 for (kz = 0; kz < 6; ++kz) {
913 Gsdk(ipa, bratio, mode);
915 for (kz = 0; kz < 6; ++kz) {
922 Gsdk(ipa, bratio, mode);
924 for (kz = 0; kz < 6; ++kz) {
931 Gsdk(ipa, bratio, mode);
933 for (kz = 0; kz < 6; ++kz) {
940 Gsdk(ipa, bratio, mode);
943 for (kz = 0; kz < 6; ++kz) {
952 Gsdk(ipa, bratio, mode);
955 Gsdk(ipa, bratio, mode);
958 Gsdk(ipa, bratio, mode);
963 //_____________________________________________________________________________
964 Int_t TGeant3::VolId(Text_t *name) const
967 // Return the unique numeric identifier for volume name
970 strncpy((char *) &gname, name, 4);
971 for(i=1; i<=fGcnum->nvolum; i++)
972 if(gname == fZiq[fGclink->jvolum+i]) return i;
973 printf("VolId: Volume %s not found\n",name);
977 //_____________________________________________________________________________
978 Int_t TGeant3::NofVolumes() const
981 // Return total number of volumes in the geometry
983 return fGcnum->nvolum;
986 //_____________________________________________________________________________
987 const char* TGeant3::VolName(Int_t id) const
990 // Return the volume name given the volume identifier
993 if(id<1 || id > fGcnum->nvolum || fGclink->jvolum<=0)
996 strncpy(name,(char *)&fZiq[fGclink->jvolum+id],4);
1001 //_____________________________________________________________________________
1002 Float_t TGeant3::Xsec(char* reac, Float_t energy, Int_t part, Int_t mate)
1004 Int_t gpart = IdFromPDG(part);
1005 if(!strcmp(reac,"PHOT"))
1008 Error("Xsec","Can calculate photoelectric only for photons\n");
1014 //_____________________________________________________________________________
1015 void TGeant3::TrackPosition(TLorentzVector &xyz) const
1018 // Return the current position in the master reference frame of the
1019 // track being transported
1021 xyz[0]=fGctrak->vect[0];
1022 xyz[1]=fGctrak->vect[1];
1023 xyz[2]=fGctrak->vect[2];
1024 xyz[3]=fGctrak->tofg;
1027 //_____________________________________________________________________________
1028 Float_t TGeant3::TrackTime() const
1031 // Return the current time of flight of the track being transported
1033 return fGctrak->tofg;
1036 //_____________________________________________________________________________
1037 void TGeant3::TrackMomentum(TLorentzVector &xyz) const
1040 // Return the direction and the momentum (GeV/c) of the track
1041 // currently being transported
1043 Double_t ptot=fGctrak->vect[6];
1044 xyz[0]=fGctrak->vect[3]*ptot;
1045 xyz[1]=fGctrak->vect[4]*ptot;
1046 xyz[2]=fGctrak->vect[5]*ptot;
1047 xyz[3]=fGctrak->getot;
1050 //_____________________________________________________________________________
1051 Float_t TGeant3::TrackCharge() const
1054 // Return charge of the track currently transported
1056 return fGckine->charge;
1059 //_____________________________________________________________________________
1060 Float_t TGeant3::TrackMass() const
1063 // Return the mass of the track currently transported
1065 return fGckine->amass;
1068 //_____________________________________________________________________________
1069 Int_t TGeant3::TrackPid() const
1072 // Return the id of the particle transported
1074 return PDGFromId(fGckine->ipart);
1077 //_____________________________________________________________________________
1078 Float_t TGeant3::TrackStep() const
1081 // Return the length in centimeters of the current step
1083 return fGctrak->step;
1086 //_____________________________________________________________________________
1087 Float_t TGeant3::TrackLength() const
1090 // Return the length of the current track from its origin
1092 return fGctrak->sleng;
1095 //_____________________________________________________________________________
1096 Bool_t TGeant3::IsTrackInside() const
1099 // True if the track is not at the boundary of the current volume
1101 return (fGctrak->inwvol==0);
1104 //_____________________________________________________________________________
1105 Bool_t TGeant3::IsTrackEntering() const
1108 // True if this is the first step of the track in the current volume
1110 return (fGctrak->inwvol==1);
1113 //_____________________________________________________________________________
1114 Bool_t TGeant3::IsTrackExiting() const
1117 // True if this is the last step of the track in the current volume
1119 return (fGctrak->inwvol==2);
1122 //_____________________________________________________________________________
1123 Bool_t TGeant3::IsTrackOut() const
1126 // True if the track is out of the setup
1128 return (fGctrak->inwvol==3);
1131 //_____________________________________________________________________________
1132 Bool_t TGeant3::IsTrackStop() const
1135 // True if the track energy has fallen below the threshold
1137 return (fGctrak->istop==2);
1140 //_____________________________________________________________________________
1141 Int_t TGeant3::NSecondaries() const
1144 // Number of secondary particles generated in the current step
1146 return fGcking->ngkine;
1149 //_____________________________________________________________________________
1150 Int_t TGeant3::CurrentEvent() const
1153 // Number of the current event
1155 return fGcflag->idevt;
1158 //_____________________________________________________________________________
1159 void TGeant3::ProdProcess(char* proc) const
1162 // Name of the process that has produced the secondary particles
1163 // in the current step
1165 const Int_t ipmec[13] = { 5,6,7,8,9,10,11,12,21,23,25,105,108 };
1168 if(fGcking->ngkine>0) {
1169 for (km = 0; km < fGctrak->nmec; ++km) {
1170 for (im = 0; im < 13; ++im) {
1171 if (fGctrak->lmec[km] == ipmec[im]) {
1172 mec = fGctrak->lmec[km];
1173 if (0 < mec && mec < 31) {
1174 strncpy(proc,(char *)&fGctrak->namec[mec - 1],4);
1175 } else if (mec - 100 <= 30 && mec - 100 > 0) {
1176 strncpy(proc,(char *)&fGctpol->namec1[mec - 101],4);
1183 strcpy(proc,"UNKN");
1184 } else strcpy(proc,"NONE");
1187 //_____________________________________________________________________________
1188 void TGeant3::GetSecondary(Int_t isec, Int_t& ipart,
1189 TLorentzVector &x, TLorentzVector &p)
1192 // Get the parameters of the secondary track number isec produced
1193 // in the current step
1196 if(-1<isec && isec<fGcking->ngkine) {
1197 ipart=Int_t (fGcking->gkin[isec][4] +0.5);
1199 x[i]=fGckin3->gpos[isec][i];
1200 p[i]=fGcking->gkin[isec][i];
1202 x[3]=fGcking->tofd[isec];
1203 p[3]=fGcking->gkin[isec][3];
1205 printf(" * TGeant3::GetSecondary * Secondary %d does not exist\n",isec);
1206 x[0]=x[1]=x[2]=x[3]=p[0]=p[1]=p[2]=p[3]=0;
1211 //_____________________________________________________________________________
1212 void TGeant3::InitLego()
1215 SetDEBU(0,0,0); //do not print a message
1218 //_____________________________________________________________________________
1219 Bool_t TGeant3::IsTrackDisappeared() const
1222 // True if the current particle has disappered
1223 // either because it decayed or because it underwent
1224 // an inelastic collision
1226 return (fGctrak->istop==1);
1229 //_____________________________________________________________________________
1230 Bool_t TGeant3::IsTrackAlive() const
1233 // True if the current particle is alive and will continue to be
1236 return (fGctrak->istop==0);
1239 //_____________________________________________________________________________
1240 void TGeant3::StopTrack()
1243 // Stop the transport of the current particle and skip to the next
1248 //_____________________________________________________________________________
1249 void TGeant3::StopEvent()
1252 // Stop simulation of the current event and skip to the next
1257 //_____________________________________________________________________________
1258 Float_t TGeant3::MaxStep() const
1261 // Return the maximum step length in the current medium
1263 return fGctmed->stemax;
1266 //_____________________________________________________________________________
1267 void TGeant3::SetColors()
1270 // Set the colors for all the volumes
1271 // this is done sequentially for all volumes
1272 // based on the number of their medium
1275 Int_t jvolum=fGclink->jvolum;
1276 //Int_t jtmed=fGclink->jtmed;
1277 //Int_t jmate=fGclink->jmate;
1278 Int_t nvolum=fGcnum->nvolum;
1281 // Now for all the volumes
1282 for(kv=1;kv<=nvolum;kv++) {
1283 // Get the tracking medium
1284 Int_t itm=Int_t (fZq[fZlq[jvolum-kv]+4]);
1286 //Int_t ima=Int_t (fZq[fZlq[jtmed-itm]+6]);
1288 //Float_t z=fZq[fZlq[jmate-ima]+7];
1289 // Find color number
1290 //icol = Int_t(z)%6+2;
1291 //icol = 17+Int_t(z*150./92.);
1294 strncpy(name,(char*)&fZiq[jvolum+kv],4);
1296 Gsatt(name,"COLO",icol);
1300 //_____________________________________________________________________________
1301 void TGeant3::SetMaxStep(Float_t maxstep)
1304 // Set the maximum step allowed till the particle is in the current medium
1306 fGctmed->stemax=maxstep;
1309 //_____________________________________________________________________________
1310 void TGeant3::SetMaxNStep(Int_t maxnstp)
1313 // Set the maximum number of steps till the particle is in the current medium
1315 fGctrak->maxnst=maxnstp;
1318 //_____________________________________________________________________________
1319 Int_t TGeant3::GetMaxNStep() const
1322 // Maximum number of steps allowed in current medium
1324 return fGctrak->maxnst;
1327 //_____________________________________________________________________________
1328 void TGeant3::Material(Int_t& kmat, const char* name, Float_t a, Float_t z,
1329 Float_t dens, Float_t radl, Float_t absl, Float_t* buf,
1333 // Defines a Material
1335 // kmat number assigned to the material
1336 // name material name
1337 // a atomic mass in au
1339 // dens density in g/cm3
1340 // absl absorbtion length in cm
1341 // if >=0 it is ignored and the program
1342 // calculates it, if <0. -absl is taken
1343 // radl radiation length in cm
1344 // if >=0 it is ignored and the program
1345 // calculates it, if <0. -radl is taken
1346 // buf pointer to an array of user words
1347 // nbuf number of user words
1349 Int_t jmate=fGclink->jmate;
1355 for(i=1; i<=ns; i++) {
1356 if(fZlq[jmate-i]==0) {
1362 gsmate(kmat,PASSCHARD(name), a, z, dens, radl, absl, buf,
1363 nwbuf PASSCHARL(name));
1366 //_____________________________________________________________________________
1367 void TGeant3::Mixture(Int_t& kmat, const char* name, Float_t* a, Float_t* z,
1368 Float_t dens, Int_t nlmat, Float_t* wmat)
1371 // Defines mixture OR COMPOUND IMAT as composed by
1372 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
1374 // If NLMAT > 0 then wmat contains the proportion by
1375 // weights of each basic material in the mixture.
1377 // If nlmat < 0 then WMAT contains the number of atoms
1378 // of a given kind into the molecule of the COMPOUND
1379 // In this case, WMAT in output is changed to relative
1382 Int_t jmate=fGclink->jmate;
1388 for(i=1; i<=ns; i++) {
1389 if(fZlq[jmate-i]==0) {
1395 gsmixt(kmat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
1398 //_____________________________________________________________________________
1399 void TGeant3::Medium(Int_t& kmed, const char* name, Int_t nmat, Int_t isvol,
1400 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
1401 Float_t stemax, Float_t deemax, Float_t epsil,
1402 Float_t stmin, Float_t* ubuf, Int_t nbuf)
1405 // kmed tracking medium number assigned
1406 // name tracking medium name
1407 // nmat material number
1408 // isvol sensitive volume flag
1409 // ifield magnetic field
1410 // fieldm max. field value (kilogauss)
1411 // tmaxfd max. angle due to field (deg/step)
1412 // stemax max. step allowed
1413 // deemax max. fraction of energy lost in a step
1414 // epsil tracking precision (cm)
1415 // stmin min. step due to continuos processes (cm)
1417 // ifield = 0 if no magnetic field; ifield = -1 if user decision in guswim;
1418 // ifield = 1 if tracking performed with grkuta; ifield = 2 if tracking
1419 // performed with ghelix; ifield = 3 if tracking performed with ghelx3.
1421 Int_t jtmed=fGclink->jtmed;
1427 for(i=1; i<=ns; i++) {
1428 if(fZlq[jtmed-i]==0) {
1434 gstmed(kmed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1435 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1438 //_____________________________________________________________________________
1439 void TGeant3::Matrix(Int_t& krot, Float_t thex, Float_t phix, Float_t they,
1440 Float_t phiy, Float_t thez, Float_t phiz)
1443 // krot rotation matrix number assigned
1444 // theta1 polar angle for axis i
1445 // phi1 azimuthal angle for axis i
1446 // theta2 polar angle for axis ii
1447 // phi2 azimuthal angle for axis ii
1448 // theta3 polar angle for axis iii
1449 // phi3 azimuthal angle for axis iii
1451 // it defines the rotation matrix number irot.
1453 Int_t jrotm=fGclink->jrotm;
1459 for(i=1; i<=ns; i++) {
1460 if(fZlq[jrotm-i]==0) {
1466 gsrotm(krot, thex, phix, they, phiy, thez, phiz);
1469 //_____________________________________________________________________________
1470 Int_t TGeant3::GetMedium() const
1473 // Return the number of the current medium
1475 return fGctmed->numed;
1478 //_____________________________________________________________________________
1479 Float_t TGeant3::Edep() const
1482 // Return the energy lost in the current step
1484 return fGctrak->destep;
1487 //_____________________________________________________________________________
1488 Float_t TGeant3::Etot() const
1491 // Return the total energy of the current track
1493 return fGctrak->getot;
1496 //_____________________________________________________________________________
1497 void TGeant3::Rndm(Float_t* r, const Int_t n) const
1500 // Return an array of n random numbers uniformly distributed
1501 // between 0 and 1 not included
1506 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1508 // Functions from GBASE
1510 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1512 //____________________________________________________________________________
1513 void TGeant3::Gfile(const char *filename, const char *option)
1516 // Routine to open a GEANT/RZ data base.
1518 // LUN logical unit number associated to the file
1520 // CHFILE RZ file name
1522 // CHOPT is a character string which may be
1523 // N To create a new file
1524 // U to open an existing file for update
1525 // " " to open an existing file for read only
1526 // Q The initial allocation (default 1000 records)
1527 // is given in IQUEST(10)
1528 // X Open the file in exchange format
1529 // I Read all data structures from file to memory
1530 // O Write all data structures from memory to file
1533 // If options "I" or "O" all data structures are read or
1534 // written from/to file and the file is closed.
1535 // See routine GRMDIR to create subdirectories
1536 // See routines GROUT,GRIN to write,read objects
1538 grfile(21, PASSCHARD(filename), PASSCHARD(option) PASSCHARL(filename)
1542 //____________________________________________________________________________
1543 void TGeant3::Gpcxyz()
1546 // Print track and volume parameters at current point
1551 //_____________________________________________________________________________
1552 void TGeant3::Ggclos()
1555 // Closes off the geometry setting.
1556 // Initializes the search list for the contents of each
1557 // volume following the order they have been positioned, and
1558 // inserting the content '0' when a call to GSNEXT (-1) has
1559 // been required by the user.
1560 // Performs the development of the JVOLUM structure for all
1561 // volumes with variable parameters, by calling GGDVLP.
1562 // Interprets the user calls to GSORD, through GGORD.
1563 // Computes and stores in a bank (next to JVOLUM mother bank)
1564 // the number of levels in the geometrical tree and the
1565 // maximum number of contents per level, by calling GGNLEV.
1566 // Sets status bit for CONCAVE volumes, through GGCAVE.
1567 // Completes the JSET structure with the list of volume names
1568 // which identify uniquely a given physical detector, the
1569 // list of bit numbers to pack the corresponding volume copy
1570 // numbers, and the generic path(s) in the JVOLUM tree,
1571 // through the routine GHCLOS.
1576 //_____________________________________________________________________________
1577 void TGeant3::Glast()
1580 // Finish a Geant run
1585 //_____________________________________________________________________________
1586 void TGeant3::Gprint(const char *name)
1589 // Routine to print data structures
1590 // CHNAME name of a data structure
1594 gprint(PASSCHARD(vname),0 PASSCHARL(vname));
1597 //_____________________________________________________________________________
1598 void TGeant3::Grun()
1601 // Steering function to process one run
1606 //_____________________________________________________________________________
1607 void TGeant3::Gtrig()
1610 // Steering function to process one event
1615 //_____________________________________________________________________________
1616 void TGeant3::Gtrigc()
1619 // Clear event partition
1624 //_____________________________________________________________________________
1625 void TGeant3::Gtrigi()
1628 // Initialises event partition
1633 //_____________________________________________________________________________
1634 void TGeant3::Gwork(Int_t nwork)
1637 // Allocates workspace in ZEBRA memory
1642 //_____________________________________________________________________________
1643 void TGeant3::Gzinit()
1646 // To initialise GEANT/ZEBRA data structures
1651 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1653 // Functions from GCONS
1655 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1657 //_____________________________________________________________________________
1658 void TGeant3::Gfmate(Int_t imat, char *name, Float_t &a, Float_t &z,
1659 Float_t &dens, Float_t &radl, Float_t &absl,
1660 Float_t* ubuf, Int_t& nbuf)
1663 // Return parameters for material IMAT
1665 gfmate(imat, PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
1669 //_____________________________________________________________________________
1670 void TGeant3::Gfpart(Int_t ipart, char *name, Int_t &itrtyp,
1671 Float_t &amass, Float_t &charge, Float_t &tlife)
1674 // Return parameters for particle of type IPART
1678 Int_t igpart = IdFromPDG(ipart);
1679 gfpart(igpart, PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
1683 //_____________________________________________________________________________
1684 void TGeant3::Gftmed(Int_t numed, char *name, Int_t &nmat, Int_t &isvol,
1685 Int_t &ifield, Float_t &fieldm, Float_t &tmaxfd,
1686 Float_t &stemax, Float_t &deemax, Float_t &epsil,
1687 Float_t &stmin, Float_t *ubuf, Int_t *nbuf)
1690 // Return parameters for tracking medium NUMED
1692 gftmed(numed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1693 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1696 //_____________________________________________________________________________
1697 void TGeant3::Gmate()
1700 // Define standard GEANT materials
1705 //_____________________________________________________________________________
1706 void TGeant3::Gpart()
1709 // Define standard GEANT particles plus selected decay modes
1710 // and branching ratios.
1715 //_____________________________________________________________________________
1716 void TGeant3::Gsdk(Int_t ipart, Float_t *bratio, Int_t *mode)
1718 // Defines branching ratios and decay modes for standard
1720 gsdk(ipart,bratio,mode);
1723 //_____________________________________________________________________________
1724 void TGeant3::Gsmate(Int_t imat, const char *name, Float_t a, Float_t z,
1725 Float_t dens, Float_t radl, Float_t absl)
1728 // Defines a Material
1730 // kmat number assigned to the material
1731 // name material name
1732 // a atomic mass in au
1734 // dens density in g/cm3
1735 // absl absorbtion length in cm
1736 // if >=0 it is ignored and the program
1737 // calculates it, if <0. -absl is taken
1738 // radl radiation length in cm
1739 // if >=0 it is ignored and the program
1740 // calculates it, if <0. -radl is taken
1741 // buf pointer to an array of user words
1742 // nbuf number of user words
1746 gsmate(imat,PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
1750 //_____________________________________________________________________________
1751 void TGeant3::Gsmixt(Int_t imat, const char *name, Float_t *a, Float_t *z,
1752 Float_t dens, Int_t nlmat, Float_t *wmat)
1755 // Defines mixture OR COMPOUND IMAT as composed by
1756 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
1758 // If NLMAT.GT.0 then WMAT contains the PROPORTION BY
1759 // WEIGTHS OF EACH BASIC MATERIAL IN THE MIXTURE.
1761 // If NLMAT.LT.0 then WMAT contains the number of atoms
1762 // of a given kind into the molecule of the COMPOUND
1763 // In this case, WMAT in output is changed to relative
1766 gsmixt(imat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
1769 //_____________________________________________________________________________
1770 void TGeant3::Gspart(Int_t ipart, const char *name, Int_t itrtyp,
1771 Float_t amass, Float_t charge, Float_t tlife)
1774 // Store particle parameters
1776 // ipart particle code
1777 // name particle name
1778 // itrtyp transport method (see GEANT manual)
1779 // amass mass in GeV/c2
1780 // charge charge in electron units
1781 // tlife lifetime in seconds
1785 gspart(ipart,PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
1789 //_____________________________________________________________________________
1790 void TGeant3::Gstmed(Int_t numed, const char *name, Int_t nmat, Int_t isvol,
1791 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
1792 Float_t stemax, Float_t deemax, Float_t epsil,
1796 // NTMED Tracking medium number
1797 // NAME Tracking medium name
1798 // NMAT Material number
1799 // ISVOL Sensitive volume flag
1800 // IFIELD Magnetic field
1801 // FIELDM Max. field value (Kilogauss)
1802 // TMAXFD Max. angle due to field (deg/step)
1803 // STEMAX Max. step allowed
1804 // DEEMAX Max. fraction of energy lost in a step
1805 // EPSIL Tracking precision (cm)
1806 // STMIN Min. step due to continuos processes (cm)
1808 // IFIELD = 0 if no magnetic field; IFIELD = -1 if user decision in GUSWIM;
1809 // IFIELD = 1 if tracking performed with GRKUTA; IFIELD = 2 if tracking
1810 // performed with GHELIX; IFIELD = 3 if tracking performed with GHELX3.
1814 gstmed(numed,PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1815 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1818 //_____________________________________________________________________________
1819 void TGeant3::Gsckov(Int_t itmed, Int_t npckov, Float_t *ppckov,
1820 Float_t *absco, Float_t *effic, Float_t *rindex)
1823 // Stores the tables for UV photon tracking in medium ITMED
1824 // Please note that it is the user's responsability to
1825 // provide all the coefficients:
1828 // ITMED Tracking medium number
1829 // NPCKOV Number of bins of each table
1830 // PPCKOV Value of photon momentum (in GeV)
1831 // ABSCO Absorbtion coefficients
1832 // dielectric: absorbtion length in cm
1833 // metals : absorbtion fraction (0<=x<=1)
1834 // EFFIC Detection efficiency for UV photons
1835 // RINDEX Refraction index (if=0 metal)
1837 gsckov(itmed,npckov,ppckov,absco,effic,rindex);
1840 //_____________________________________________________________________________
1841 void TGeant3::Gstpar(Int_t itmed, const char *param, Float_t parval)
1844 // To change the value of cut or mechanism "CHPAR"
1845 // to a new value PARVAL for tracking medium ITMED
1846 // The data structure JTMED contains the standard tracking
1847 // parameters (CUTS and flags to control the physics processes) which
1848 // are used by default for all tracking media. It is possible to
1849 // redefine individually with GSTPAR any of these parameters for a
1850 // given tracking medium.
1851 // ITMED tracking medium number
1852 // CHPAR is a character string (variable name)
1853 // PARVAL must be given as a floating point.
1855 gstpar(itmed,PASSCHARD(param), parval PASSCHARL(param));
1858 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1860 // Functions from GCONS
1862 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1864 //_____________________________________________________________________________
1865 void TGeant3::Gfkine(Int_t itra, Float_t *vert, Float_t *pvert, Int_t &ipart,
1868 // Storing/Retrieving Vertex and Track parameters
1869 // ----------------------------------------------
1871 // Stores vertex parameters.
1872 // VERT array of (x,y,z) position of the vertex
1873 // NTBEAM beam track number origin of the vertex
1874 // =0 if none exists
1875 // NTTARG target track number origin of the vertex
1876 // UBUF user array of NUBUF floating point numbers
1878 // NVTX new vertex number (=0 in case of error).
1879 // Prints vertex parameters.
1880 // IVTX for vertex IVTX.
1881 // (For all vertices if IVTX=0)
1882 // Stores long life track parameters.
1883 // PLAB components of momentum
1884 // IPART type of particle (see GSPART)
1885 // NV vertex number origin of track
1886 // UBUF array of NUBUF floating point user parameters
1888 // NT track number (if=0 error).
1889 // Retrieves long life track parameters.
1890 // ITRA track number for which parameters are requested
1891 // VERT vector origin of the track
1892 // PVERT 4 momentum components at the track origin
1893 // IPART particle type (=0 if track ITRA does not exist)
1894 // NVERT vertex number origin of the track
1895 // UBUF user words stored in GSKINE.
1896 // Prints initial track parameters.
1897 // ITRA for track ITRA
1898 // (For all tracks if ITRA=0)
1902 gfkine(itra,vert,pvert,ipart,nvert,ubuf,nbuf);
1905 //_____________________________________________________________________________
1906 void TGeant3::Gfvert(Int_t nvtx, Float_t *v, Int_t &ntbeam, Int_t &nttarg,
1910 // Retrieves the parameter of a vertex bank
1911 // Vertex is generated from tracks NTBEAM NTTARG
1912 // NVTX is the new vertex number
1916 gfvert(nvtx,v,ntbeam,nttarg,tofg,ubuf,nbuf);
1919 //_____________________________________________________________________________
1920 Int_t TGeant3::Gskine(Float_t *plab, Int_t ipart, Int_t nv, Float_t *buf,
1924 // Store kinematics of track NT into data structure
1925 // Track is coming from vertex NV
1928 gskine(plab, ipart, nv, buf, nwbuf, nt);
1932 //_____________________________________________________________________________
1933 Int_t TGeant3::Gsvert(Float_t *v, Int_t ntbeam, Int_t nttarg, Float_t *ubuf,
1937 // Creates a new vertex bank
1938 // Vertex is generated from tracks NTBEAM NTTARG
1939 // NVTX is the new vertex number
1942 gsvert(v, ntbeam, nttarg, ubuf, nwbuf, nwtx);
1946 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1948 // Functions from GPHYS
1950 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1952 //_____________________________________________________________________________
1953 void TGeant3::Gphysi()
1956 // Initialise material constants for all the physics
1957 // mechanisms used by GEANT
1962 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1964 // Functions from GTRAK
1966 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1968 //_____________________________________________________________________________
1969 void TGeant3::Gdebug()
1972 // Debug the current step
1977 //_____________________________________________________________________________
1978 void TGeant3::Gekbin()
1981 // To find bin number in kinetic energy table
1982 // stored in ELOW(NEKBIN)
1987 //_____________________________________________________________________________
1988 void TGeant3::Gfinds()
1991 // Returns the set/volume parameters corresponding to
1992 // the current space point in /GCTRAK/
1993 // and fill common /GCSETS/
1995 // IHSET user set identifier
1996 // IHDET user detector identifier
1997 // ISET set number in JSET
1998 // IDET detector number in JS=LQ(JSET-ISET)
1999 // IDTYPE detector type (1,2)
2000 // NUMBV detector volume numbers (array of length NVNAME)
2001 // NVNAME number of volume levels
2006 //_____________________________________________________________________________
2007 void TGeant3::Gsking(Int_t igk)
2010 // Stores in stack JSTAK either the IGKth track of /GCKING/,
2011 // or the NGKINE tracks when IGK is 0.
2016 //_____________________________________________________________________________
2017 void TGeant3::Gskpho(Int_t igk)
2020 // Stores in stack JSTAK either the IGKth Cherenkov photon of
2021 // /GCKIN2/, or the NPHOT tracks when IGK is 0.
2026 //_____________________________________________________________________________
2027 void TGeant3::Gsstak(Int_t iflag)
2030 // Stores in auxiliary stack JSTAK the particle currently
2031 // described in common /GCKINE/.
2033 // On request, creates also an entry in structure JKINE :
2035 // 0 : No entry in JKINE structure required (user)
2036 // 1 : New entry in JVERTX / JKINE structures required (user)
2037 // <0 : New entry in JKINE structure at vertex -IFLAG (user)
2038 // 2 : Entry in JKINE structure exists already (from GTREVE)
2043 //_____________________________________________________________________________
2044 void TGeant3::Gsxyz()
2047 // Store space point VECT in banks JXYZ
2052 //_____________________________________________________________________________
2053 void TGeant3::Gtrack()
2056 // Controls tracking of current particle
2061 //_____________________________________________________________________________
2062 void TGeant3::Gtreve()
2065 // Controls tracking of all particles belonging to the current event
2070 //_____________________________________________________________________________
2071 void TGeant3::Gtreve_root()
2074 // Controls tracking of all particles belonging to the current event
2079 //_____________________________________________________________________________
2080 void TGeant3::Grndm(Float_t *rvec, const Int_t len) const
2083 // To generate a vector RVECV of LEN random numbers
2084 // Copy of the CERN Library routine RANECU
2088 //_____________________________________________________________________________
2089 void TGeant3::Grndmq(Int_t &is1, Int_t &is2, const Int_t iseq,
2090 const Text_t *chopt)
2093 // To set/retrieve the seed of the random number generator
2095 grndmq(is1,is2,iseq,PASSCHARD(chopt) PASSCHARL(chopt));
2098 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2100 // Functions from GDRAW
2102 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2104 //_____________________________________________________________________________
2105 void TGeant3::Gdxyz(Int_t it)
2108 // Draw the points stored with Gsxyz relative to track it
2113 //_____________________________________________________________________________
2114 void TGeant3::Gdcxyz()
2117 // Draw the position of the current track
2122 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2124 // Functions from GGEOM
2126 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2128 //_____________________________________________________________________________
2129 void TGeant3::Gdtom(Float_t *xd, Float_t *xm, Int_t iflag)
2132 // Computes coordinates XM (Master Reference System
2133 // knowing the coordinates XD (Detector Ref System)
2134 // The local reference system can be initialized by
2135 // - the tracking routines and GDTOM used in GUSTEP
2136 // - a call to GSCMED(NLEVEL,NAMES,NUMBER)
2137 // (inverse routine is GMTOD)
2139 // If IFLAG=1 convert coordinates
2140 // IFLAG=2 convert direction cosinus
2142 gdtom(xd, xm, iflag);
2145 //_____________________________________________________________________________
2146 void TGeant3::Glmoth(const char* iudet, Int_t iunum, Int_t &nlev, Int_t *lvols,
2150 // Loads the top part of the Volume tree in LVOLS (IVO's),
2151 // LINDX (IN indices) for a given volume defined through
2152 // its name IUDET and number IUNUM.
2154 // The routine stores only upto the last level where JVOLUM
2155 // data structure is developed. If there is no development
2156 // above the current level, it returns NLEV zero.
2158 glmoth(PASSCHARD(iudet), iunum, nlev, lvols, lindx, idum PASSCHARL(iudet));
2161 //_____________________________________________________________________________
2162 void TGeant3::Gmedia(Float_t *x, Int_t &numed)
2165 // Finds in which volume/medium the point X is, and updates the
2166 // common /GCVOLU/ and the structure JGPAR accordingly.
2168 // NUMED returns the tracking medium number, or 0 if point is
2169 // outside the experimental setup.
2174 //_____________________________________________________________________________
2175 void TGeant3::Gmtod(Float_t *xm, Float_t *xd, Int_t iflag)
2178 // Computes coordinates XD (in DRS)
2179 // from known coordinates XM in MRS
2180 // The local reference system can be initialized by
2181 // - the tracking routines and GMTOD used in GUSTEP
2182 // - a call to GMEDIA(XM,NUMED)
2183 // - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER)
2184 // (inverse routine is GDTOM)
2186 // If IFLAG=1 convert coordinates
2187 // IFLAG=2 convert direction cosinus
2189 gmtod(xm, xd, iflag);
2192 //_____________________________________________________________________________
2193 void TGeant3::Gsdvn(const char *name, const char *mother, Int_t ndiv,
2197 // Create a new volume by dividing an existing one
2200 // MOTHER Mother volume name
2201 // NDIV Number of divisions
2204 // X,Y,Z of CAXIS will be translated to 1,2,3 for IAXIS.
2205 // It divides a previously defined volume.
2210 Vname(mother,vmother);
2211 gsdvn(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis PASSCHARL(vname)
2212 PASSCHARL(vmother));
2215 //_____________________________________________________________________________
2216 void TGeant3::Gsdvn2(const char *name, const char *mother, Int_t ndiv,
2217 Int_t iaxis, Float_t c0i, Int_t numed)
2220 // Create a new volume by dividing an existing one
2222 // Divides mother into ndiv divisions called name
2223 // along axis iaxis starting at coordinate value c0.
2224 // the new volume created will be medium number numed.
2229 Vname(mother,vmother);
2230 gsdvn2(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis, c0i, numed
2231 PASSCHARL(vname) PASSCHARL(vmother));
2234 //_____________________________________________________________________________
2235 void TGeant3::Gsdvs(const char *name, const char *mother, Float_t step,
2236 Int_t iaxis, Int_t numed)
2239 // Create a new volume by dividing an existing one
2244 Vname(mother,vmother);
2245 gsdvs(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed
2246 PASSCHARL(vname) PASSCHARL(vmother));
2249 //_____________________________________________________________________________
2250 void TGeant3::Gsdvs2(const char *name, const char *mother, Float_t step,
2251 Int_t iaxis, Float_t c0, Int_t numed)
2254 // Create a new volume by dividing an existing one
2259 Vname(mother,vmother);
2260 gsdvs2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0, numed
2261 PASSCHARL(vname) PASSCHARL(vmother));
2264 //_____________________________________________________________________________
2265 void TGeant3::Gsdvt(const char *name, const char *mother, Float_t step,
2266 Int_t iaxis, Int_t numed, Int_t ndvmx)
2269 // Create a new volume by dividing an existing one
2271 // Divides MOTHER into divisions called NAME along
2272 // axis IAXIS in steps of STEP. If not exactly divisible
2273 // will make as many as possible and will centre them
2274 // with respect to the mother. Divisions will have medium
2275 // number NUMED. If NUMED is 0, NUMED of MOTHER is taken.
2276 // NDVMX is the expected maximum number of divisions
2277 // (If 0, no protection tests are performed)
2282 Vname(mother,vmother);
2283 gsdvt(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed, ndvmx
2284 PASSCHARL(vname) PASSCHARL(vmother));
2287 //_____________________________________________________________________________
2288 void TGeant3::Gsdvt2(const char *name, const char *mother, Float_t step,
2289 Int_t iaxis, Float_t c0, Int_t numed, Int_t ndvmx)
2292 // Create a new volume by dividing an existing one
2294 // Divides MOTHER into divisions called NAME along
2295 // axis IAXIS starting at coordinate value C0 with step
2297 // The new volume created will have medium number NUMED.
2298 // If NUMED is 0, NUMED of mother is taken.
2299 // NDVMX is the expected maximum number of divisions
2300 // (If 0, no protection tests are performed)
2305 Vname(mother,vmother);
2306 gsdvt2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0,
2307 numed, ndvmx PASSCHARL(vname) PASSCHARL(vmother));
2310 //_____________________________________________________________________________
2311 void TGeant3::Gsord(const char *name, Int_t iax)
2314 // Flags volume CHNAME whose contents will have to be ordered
2315 // along axis IAX, by setting the search flag to -IAX
2319 // IAX = 4 Rxy (static ordering only -> GTMEDI)
2320 // IAX = 14 Rxy (also dynamic ordering -> GTNEXT)
2321 // IAX = 5 Rxyz (static ordering only -> GTMEDI)
2322 // IAX = 15 Rxyz (also dynamic ordering -> GTNEXT)
2323 // IAX = 6 PHI (PHI=0 => X axis)
2324 // IAX = 7 THETA (THETA=0 => Z axis)
2328 gsord(PASSCHARD(vname), iax PASSCHARL(vname));
2331 //_____________________________________________________________________________
2332 void TGeant3::Gspos(const char *name, Int_t nr, const char *mother, Float_t x,
2333 Float_t y, Float_t z, Int_t irot, const char *konly)
2336 // Position a volume into an existing one
2339 // NUMBER Copy number of the volume
2340 // MOTHER Mother volume name
2341 // X X coord. of the volume in mother ref. sys.
2342 // Y Y coord. of the volume in mother ref. sys.
2343 // Z Z coord. of the volume in mother ref. sys.
2344 // IROT Rotation matrix number w.r.t. mother ref. sys.
2345 // ONLY ONLY/MANY flag
2347 // It positions a previously defined volume in the mother.
2352 Vname(mother,vmother);
2353 gspos(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2354 PASSCHARD(konly) PASSCHARL(vname) PASSCHARL(vmother)
2358 //_____________________________________________________________________________
2359 void TGeant3::Gsposp(const char *name, Int_t nr, const char *mother,
2360 Float_t x, Float_t y, Float_t z, Int_t irot,
2361 const char *konly, Float_t *upar, Int_t np )
2364 // Place a copy of generic volume NAME with user number
2365 // NR inside MOTHER, with its parameters UPAR(1..NP)
2370 Vname(mother,vmother);
2371 gsposp(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2372 PASSCHARD(konly), upar, np PASSCHARL(vname) PASSCHARL(vmother)
2376 //_____________________________________________________________________________
2377 void TGeant3::Gsrotm(Int_t nmat, Float_t theta1, Float_t phi1, Float_t theta2,
2378 Float_t phi2, Float_t theta3, Float_t phi3)
2381 // nmat Rotation matrix number
2382 // THETA1 Polar angle for axis I
2383 // PHI1 Azimuthal angle for axis I
2384 // THETA2 Polar angle for axis II
2385 // PHI2 Azimuthal angle for axis II
2386 // THETA3 Polar angle for axis III
2387 // PHI3 Azimuthal angle for axis III
2389 // It defines the rotation matrix number IROT.
2391 gsrotm(nmat, theta1, phi1, theta2, phi2, theta3, phi3);
2394 //_____________________________________________________________________________
2395 void TGeant3::Gprotm(Int_t nmat)
2398 // To print rotation matrices structure JROTM
2399 // nmat Rotation matrix number
2404 //_____________________________________________________________________________
2405 Int_t TGeant3::Gsvolu(const char *name, const char *shape, Int_t nmed,
2406 Float_t *upar, Int_t npar)
2410 // SHAPE Volume type
2411 // NUMED Tracking medium number
2412 // NPAR Number of shape parameters
2413 // UPAR Vector containing shape parameters
2415 // It creates a new volume in the JVOLUM data structure.
2421 Vname(shape,vshape);
2422 gsvolu(PASSCHARD(vname), PASSCHARD(vshape), nmed, upar, npar, ivolu
2423 PASSCHARL(vname) PASSCHARL(vshape));
2427 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2429 // T H E D R A W I N G P A C K A G E
2430 // ======================================
2431 // Drawing functions. These functions allow the visualization in several ways
2432 // of the volumes defined in the geometrical data structure. It is possible
2433 // to draw the logical tree of volumes belonging to the detector (DTREE),
2434 // to show their geometrical specification (DSPEC,DFSPC), to draw them
2435 // and their cut views (DRAW, DCUT). Moreover, it is possible to execute
2436 // these commands when the hidden line removal option is activated; in
2437 // this case, the volumes can be also either translated in the space
2438 // (SHIFT), or clipped by boolean operation (CVOL). In addition, it is
2439 // possible to fill the surfaces of the volumes
2440 // with solid colours when the shading option (SHAD) is activated.
2441 // Several tools (ZOOM, LENS) have been developed to zoom detailed parts
2442 // of the detectors or to scan physical events as well.
2443 // Finally, the command MOVE will allow the rotation, translation and zooming
2444 // on real time parts of the detectors or tracks and hits of a simulated event.
2445 // Ray-tracing commands. In case the command (DOPT RAYT ON) is executed,
2446 // the drawing is performed by the Geant ray-tracing;
2447 // automatically, the color is assigned according to the tracking medium of each
2448 // volume and the volumes with a density lower/equal than the air are considered
2449 // transparent; if the option (USER) is set (ON) (again via the command (DOPT)),
2450 // the user can set color and visibility for the desired volumes via the command
2451 // (SATT), as usual, relatively to the attributes (COLO) and (SEEN).
2452 // The resolution can be set via the command (SATT * FILL VALUE), where (VALUE)
2453 // is the ratio between the number of pixels drawn and 20 (user coordinates).
2454 // Parallel view and perspective view are possible (DOPT PROJ PARA/PERS); in the
2455 // first case, we assume that the first mother volume of the tree is a box with
2456 // dimensions 10000 X 10000 X 10000 cm and the view point (infinetely far) is
2457 // 5000 cm far from the origin along the Z axis of the user coordinates; in the
2458 // second case, the distance between the observer and the origin of the world
2459 // reference system is set in cm by the command (PERSP NAME VALUE); grand-angle
2460 // or telescopic effects can be achieved changing the scale factors in the command
2461 // (DRAW). When the final picture does not occupy the full window,
2462 // mapping the space before tracing can speed up the drawing, but can also
2463 // produce less precise results; values from 1 to 4 are allowed in the command
2464 // (DOPT MAPP VALUE), the mapping being more precise for increasing (VALUE); for
2465 // (VALUE = 0) no mapping is performed (therefore max precision and lowest speed).
2466 // The command (VALCUT) allows the cutting of the detector by three planes
2467 // ortogonal to the x,y,z axis. The attribute (LSTY) can be set by the command
2468 // SATT for any desired volume and can assume values from 0 to 7; it determines
2469 // the different light processing to be performed for different materials:
2470 // 0 = dark-matt, 1 = bright-matt, 2 = plastic, 3 = ceramic, 4 = rough-metals,
2471 // 5 = shiny-metals, 6 = glass, 7 = mirror. The detector is assumed to be in the
2472 // dark, the ambient light luminosity is 0.2 for each basic hue (the saturation
2473 // is 0.9) and the observer is assumed to have a light source (therefore he will
2474 // produce parallel light in the case of parallel view and point-like-source
2475 // light in the case of perspective view).
2477 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2479 //_____________________________________________________________________________
2480 void TGeant3::Gsatt(const char *name, const char *att, Int_t val)
2484 // IOPT Name of the attribute to be set
2485 // IVAL Value to which the attribute is to be set
2487 // name= "*" stands for all the volumes.
2488 // iopt can be chosen among the following :
2490 // WORK 0=volume name is inactive for the tracking
2491 // 1=volume name is active for the tracking (default)
2493 // SEEN 0=volume name is invisible
2494 // 1=volume name is visible (default)
2495 // -1=volume invisible with all its descendants in the tree
2496 // -2=volume visible but not its descendants in the tree
2498 // LSTY line style 1,2,3,... (default=1)
2499 // LSTY=7 will produce a very precise approximation for
2500 // revolution bodies.
2502 // LWID line width -7,...,1,2,3,..7 (default=1)
2503 // LWID<0 will act as abs(LWID) was set for the volume
2504 // and for all the levels below it. When SHAD is 'ON', LWID
2505 // represent the linewidth of the scan lines filling the surfaces
2506 // (whereas the FILL value represent their number). Therefore
2507 // tuning this parameter will help to obtain the desired
2508 // quality/performance ratio.
2510 // COLO colour code -166,...,1,2,..166 (default=1)
2512 // n=2=red; n=17+m, m=0,25, increasing luminosity according to 'm';
2513 // n=3=green; n=67+m, m=0,25, increasing luminosity according to 'm';
2514 // n=4=blue; n=117+m, m=0,25, increasing luminosity according to 'm';
2515 // n=5=yellow; n=42+m, m=0,25, increasing luminosity according to 'm';
2516 // n=6=violet; n=142+m, m=0,25, increasing luminosity according to 'm';
2517 // n=7=lightblue; n=92+m, m=0,25, increasing luminosity according to 'm';
2518 // colour=n*10+m, m=1,2,...9, will produce the same colour
2519 // as 'n', but with increasing luminosity according to 'm';
2520 // COLO<0 will act as if abs(COLO) was set for the volume
2521 // and for all the levels below it.
2522 // When for a volume the attribute FILL is > 1 (and the
2523 // option SHAD is on), the ABS of its colour code must be < 8
2524 // because an automatic shading of its faces will be
2527 // FILL (1992) fill area -7,...,0,1,...7 (default=0)
2528 // when option SHAD is "on" the FILL attribute of any
2529 // volume can be set different from 0 (normal drawing);
2530 // if it is set to 1, the faces of such volume will be filled
2531 // with solid colours; if ABS(FILL) is > 1, then a light
2532 // source is placed along the observer line, and the faces of
2533 // such volumes will be painted by colours whose luminosity
2534 // will depend on the amount of light reflected;
2535 // if ABS(FILL) = 1, then it is possible to use all the 166
2536 // colours of the colour table, becouse the automatic shading
2537 // is not performed;
2538 // for increasing values of FILL the drawing will be performed
2539 // with higher and higher resolution improving the quality (the
2540 // number of scan lines used to fill the faces increases with FILL);
2541 // it is possible to set different values of FILL
2542 // for different volumes, in order to optimize at the same time
2543 // the performance and the quality of the picture;
2544 // FILL<0 will act as if abs(FILL) was set for the volume
2545 // and for all the levels below it.
2546 // This kind of drawing can be saved in 'picture files'
2547 // or in view banks.
2548 // 0=drawing without fill area
2549 // 1=faces filled with solid colours and resolution = 6
2550 // 2=lowest resolution (very fast)
2551 // 3=default resolution
2552 // 4=.................
2553 // 5=.................
2554 // 6=.................
2556 // Finally, if a coloured background is desired, the FILL
2557 // attribute for the first volume of the tree must be set
2558 // equal to -abs(colo), colo being >0 and <166.
2560 // SET set number associated to volume name
2561 // DET detector number associated to volume name
2562 // DTYP detector type (1,2)
2569 gsatt(PASSCHARD(vname), PASSCHARD(vatt), val PASSCHARL(vname)
2573 //_____________________________________________________________________________
2574 void TGeant3::Gfpara(const char *name, Int_t number, Int_t intext, Int_t& npar,
2575 Int_t& natt, Float_t* par, Float_t* att)
2578 // Find the parameters of a volume
2580 gfpara(PASSCHARD(name), number, intext, npar, natt, par, att
2584 //_____________________________________________________________________________
2585 void TGeant3::Gckpar(Int_t ish, Int_t npar, Float_t* par)
2588 // Check the parameters of a shape
2590 gckpar(ish,npar,par);
2593 //_____________________________________________________________________________
2594 void TGeant3::Gckmat(Int_t itmed, char* natmed)
2597 // Check the parameters of a tracking medium
2599 gckmat(itmed, PASSCHARD(natmed) PASSCHARL(natmed));
2602 //_____________________________________________________________________________
2603 void TGeant3::Gdelete(Int_t iview)
2606 // IVIEW View number
2608 // It deletes a view bank from memory.
2613 //_____________________________________________________________________________
2614 void TGeant3::Gdopen(Int_t iview)
2617 // IVIEW View number
2619 // When a drawing is very complex and requires a long time to be
2620 // executed, it can be useful to store it in a view bank: after a
2621 // call to DOPEN and the execution of the drawing (nothing will
2622 // appear on the screen), and after a necessary call to DCLOSE,
2623 // the contents of the bank can be displayed in a very fast way
2624 // through a call to DSHOW; therefore, the detector can be easily
2625 // zoomed many times in different ways. Please note that the pictures
2626 // with solid colours can now be stored in a view bank or in 'PICTURE FILES'
2633 //_____________________________________________________________________________
2634 void TGeant3::Gdclose()
2637 // It closes the currently open view bank; it must be called after the
2638 // end of the drawing to be stored.
2643 //_____________________________________________________________________________
2644 void TGeant3::Gdshow(Int_t iview)
2647 // IVIEW View number
2649 // It shows on the screen the contents of a view bank. It
2650 // can be called after a view bank has been closed.
2655 //_____________________________________________________________________________
2656 void TGeant3::Gdopt(const char *name,const char *value)
2660 // VALUE Option value
2662 // To set/modify the drawing options.
2665 // THRZ ON Draw tracks in R vs Z
2666 // OFF (D) Draw tracks in X,Y,Z
2669 // PROJ PARA (D) Parallel projection
2671 // TRAK LINE (D) Trajectory drawn with lines
2672 // POIN " " with markers
2673 // HIDE ON Hidden line removal using the CG package
2674 // OFF (D) No hidden line removal
2675 // SHAD ON Fill area and shading of surfaces.
2676 // OFF (D) Normal hidden line removal.
2677 // RAYT ON Ray-tracing on.
2678 // OFF (D) Ray-tracing off.
2679 // EDGE OFF Does not draw contours when shad is on.
2680 // ON (D) Normal shading.
2681 // MAPP 1,2,3,4 Mapping before ray-tracing.
2682 // 0 (D) No mapping.
2683 // USER ON User graphics options in the raytracing.
2684 // OFF (D) Automatic graphics options.
2690 Vname(value,vvalue);
2691 gdopt(PASSCHARD(vname), PASSCHARD(vvalue) PASSCHARL(vname)
2695 //_____________________________________________________________________________
2696 void TGeant3::Gdraw(const char *name,Float_t theta, Float_t phi, Float_t psi,
2697 Float_t u0,Float_t v0,Float_t ul,Float_t vl)
2702 // THETA Viewing angle theta (for 3D projection)
2703 // PHI Viewing angle phi (for 3D projection)
2704 // PSI Viewing angle psi (for 2D rotation)
2705 // U0 U-coord. (horizontal) of volume origin
2706 // V0 V-coord. (vertical) of volume origin
2707 // SU Scale factor for U-coord.
2708 // SV Scale factor for V-coord.
2710 // This function will draw the volumes,
2711 // selected with their graphical attributes, set by the Gsatt
2712 // facility. The drawing may be performed with hidden line removal
2713 // and with shading effects according to the value of the options HIDE
2714 // and SHAD; if the option SHAD is ON, the contour's edges can be
2715 // drawn or not. If the option HIDE is ON, the detector can be
2716 // exploded (BOMB), clipped with different shapes (CVOL), and some
2717 // of its parts can be shifted from their original
2718 // position (SHIFT). When HIDE is ON, if
2719 // the drawing requires more than the available memory, the program
2720 // will evaluate and display the number of missing words
2721 // (so that the user can increase the
2722 // size of its ZEBRA store). Finally, at the end of each drawing (with HIDE on),
2723 // the program will print messages about the memory used and
2724 // statistics on the volumes' visibility.
2725 // The following commands will produce the drawing of a green
2726 // volume, specified by NAME, without using the hidden line removal
2727 // technique, using the hidden line removal technique,
2728 // with different linewidth and colour (red), with
2729 // solid colour, with shading of surfaces, and without edges.
2730 // Finally, some examples are given for the ray-tracing. (A possible
2731 // string for the NAME of the volume can be found using the command DTREE).
2737 if (fGcvdma->raytra != 1) {
2738 gdraw(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
2740 gdrayt(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
2744 //_____________________________________________________________________________
2745 void TGeant3::Gdrawc(const char *name,Int_t axis, Float_t cut,Float_t u0,
2746 Float_t v0,Float_t ul,Float_t vl)
2751 // CUTVAL Cut plane distance from the origin along the axis
2753 // U0 U-coord. (horizontal) of volume origin
2754 // V0 V-coord. (vertical) of volume origin
2755 // SU Scale factor for U-coord.
2756 // SV Scale factor for V-coord.
2758 // The cut plane is normal to caxis (X,Y,Z), corresponding to iaxis (1,2,3),
2759 // and placed at the distance cutval from the origin.
2760 // The resulting picture is seen from the the same axis.
2761 // When HIDE Mode is ON, it is possible to get the same effect with
2762 // the CVOL/BOX function.
2768 gdrawc(PASSCHARD(vname), axis,cut,u0,v0,ul,vl PASSCHARL(vname));
2771 //_____________________________________________________________________________
2772 void TGeant3::Gdrawx(const char *name,Float_t cutthe, Float_t cutphi,
2773 Float_t cutval, Float_t theta, Float_t phi, Float_t u0,
2774 Float_t v0,Float_t ul,Float_t vl)
2778 // CUTTHE Theta angle of the line normal to cut plane
2779 // CUTPHI Phi angle of the line normal to cut plane
2780 // CUTVAL Cut plane distance from the origin along the axis
2782 // THETA Viewing angle theta (for 3D projection)
2783 // PHI Viewing angle phi (for 3D projection)
2784 // U0 U-coord. (horizontal) of volume origin
2785 // V0 V-coord. (vertical) of volume origin
2786 // SU Scale factor for U-coord.
2787 // SV Scale factor for V-coord.
2789 // The cut plane is normal to the line given by the cut angles
2790 // cutthe and cutphi and placed at the distance cutval from the origin.
2791 // The resulting picture is seen from the viewing angles theta,phi.
2797 gdrawx(PASSCHARD(vname), cutthe,cutphi,cutval,theta,phi,u0,v0,ul,vl
2801 //_____________________________________________________________________________
2802 void TGeant3::Gdhead(Int_t isel, const char *name, Float_t chrsiz)
2807 // ISEL Option flag D=111110
2809 // CHRSIZ Character size (cm) of title NAME D=0.6
2812 // 0 to have only the header lines
2813 // xxxxx1 to add the text name centered on top of header
2814 // xxxx1x to add global detector name (first volume) on left
2815 // xxx1xx to add date on right
2816 // xx1xxx to select thick characters for text on top of header
2817 // x1xxxx to add the text 'EVENT NR x' on top of header
2818 // 1xxxxx to add the text 'RUN NR x' on top of header
2819 // NOTE that ISEL=x1xxx1 or ISEL=1xxxx1 are illegal choices,
2820 // i.e. they generate overwritten text.
2822 gdhead(isel,PASSCHARD(name),chrsiz PASSCHARL(name));
2825 //_____________________________________________________________________________
2826 void TGeant3::Gdman(Float_t u, Float_t v, const char *type)
2829 // Draw a 2D-man at position (U0,V0)
2831 // U U-coord. (horizontal) of the centre of man' R
2832 // V V-coord. (vertical) of the centre of man' R
2833 // TYPE D='MAN' possible values: 'MAN,WM1,WM2,WM3'
2835 // CALL GDMAN(u,v),CALL GDWMN1(u,v),CALL GDWMN2(u,v),CALL GDWMN2(u,v)
2836 // It superimposes the picure of a man or of a woman, chosen among
2837 // three different ones, with the same scale factors as the detector
2838 // in the current drawing.
2841 if (opt.Contains("WM1")) {
2843 } else if (opt.Contains("WM3")) {
2845 } else if (opt.Contains("WM2")) {
2852 //_____________________________________________________________________________
2853 void TGeant3::Gdspec(const char *name)
2858 // Shows 3 views of the volume (two cut-views and a 3D view), together with
2859 // its geometrical specifications. The 3D drawing will
2860 // be performed according the current values of the options HIDE and
2861 // SHAD and according the current SetClipBox clipping parameters for that
2868 gdspec(PASSCHARD(vname) PASSCHARL(vname));
2871 //_____________________________________________________________________________
2872 void TGeant3::DrawOneSpec(const char *name)
2875 // Function called when one double-clicks on a volume name
2876 // in a TPavelabel drawn by Gdtree.
2878 THIGZ *higzSave = higz;
2879 higzSave->SetName("higzSave");
2880 THIGZ *higzSpec = (THIGZ*)gROOT->FindObject("higzSpec");
2881 //printf("DrawOneSpec, higz=%x, higzSpec=%x\n",higz,higzSpec);
2882 if (higzSpec) higz = higzSpec;
2883 else higzSpec = new THIGZ(defSize);
2884 higzSpec->SetName("higzSpec");
2889 gdspec(PASSCHARD(vname) PASSCHARL(vname));
2892 higzSave->SetName("higz");
2896 //_____________________________________________________________________________
2897 void TGeant3::Gdtree(const char *name,Int_t levmax, Int_t isel)
2901 // LEVMAX Depth level
2904 // This function draws the logical tree,
2905 // Each volume in the tree is represented by a TPaveTree object.
2906 // Double-clicking on a TPaveTree draws the specs of the corresponding volume.
2907 // Use TPaveTree pop-up menu to select:
2910 // - drawing tree of parent
2916 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
2920 //_____________________________________________________________________________
2921 void TGeant3::GdtreeParent(const char *name,Int_t levmax, Int_t isel)
2925 // LEVMAX Depth level
2928 // This function draws the logical tree of the parent of name.
2932 // Scan list of volumes in JVOLUM
2934 Int_t gname, i, jvo, in, nin, jin, num;
2935 strncpy((char *) &gname, name, 4);
2936 for(i=1; i<=fGcnum->nvolum; i++) {
2937 jvo = fZlq[fGclink->jvolum-i];
2938 nin = Int_t(fZq[jvo+3]);
2939 if (nin == -1) nin = 1;
2940 for (in=1;in<=nin;in++) {
2942 num = Int_t(fZq[jin+2]);
2943 if(gname == fZiq[fGclink->jvolum+num]) {
2944 strncpy(vname,(char*)&fZiq[fGclink->jvolum+i],4);
2946 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
2954 //_____________________________________________________________________________
2955 void TGeant3::SetABAN(Int_t par)
2958 // par = 1 particles will be stopped according to their residual
2959 // range if they are not in a sensitive material and are
2960 // far enough from the boundary
2961 // 0 particles are transported normally
2963 fGcphys->dphys1 = par;
2967 //_____________________________________________________________________________
2968 void TGeant3::SetANNI(Int_t par)
2971 // To control positron annihilation.
2972 // par =0 no annihilation
2973 // =1 annihilation. Decays processed.
2974 // =2 annihilation. No decay products stored.
2976 fGcphys->ianni = par;
2980 //_____________________________________________________________________________
2981 void TGeant3::SetAUTO(Int_t par)
2984 // To control automatic calculation of tracking medium parameters:
2985 // par =0 no automatic calculation;
2986 // =1 automati calculation.
2988 fGctrak->igauto = par;
2992 //_____________________________________________________________________________
2993 void TGeant3::SetBOMB(Float_t boom)
2996 // BOOM : Exploding factor for volumes position
2998 // To 'explode' the detector. If BOOM is positive (values smaller
2999 // than 1. are suggested, but any value is possible)
3000 // all the volumes are shifted by a distance
3001 // proportional to BOOM along the direction between their centre
3002 // and the origin of the MARS; the volumes which are symmetric
3003 // with respect to this origin are simply not shown.
3004 // BOOM equal to 0 resets the normal mode.
3005 // A negative (greater than -1.) value of
3006 // BOOM will cause an 'implosion'; for even lower values of BOOM
3007 // the volumes' positions will be reflected respect to the origin.
3008 // This command can be useful to improve the 3D effect for very
3009 // complex detectors. The following commands will make explode the
3016 //_____________________________________________________________________________
3017 void TGeant3::SetBREM(Int_t par)
3020 // To control bremstrahlung.
3021 // par =0 no bremstrahlung
3022 // =1 bremstrahlung. Photon processed.
3023 // =2 bremstrahlung. No photon stored.
3025 fGcphys->ibrem = par;
3029 //_____________________________________________________________________________
3030 void TGeant3::SetCKOV(Int_t par)
3033 // To control Cerenkov production
3034 // par =0 no Cerenkov;
3036 // =2 Cerenkov with primary stopped at each step.
3038 fGctlit->itckov = par;
3042 //_____________________________________________________________________________
3043 void TGeant3::SetClipBox(const char *name,Float_t xmin,Float_t xmax,
3044 Float_t ymin,Float_t ymax,Float_t zmin,Float_t zmax)
3047 // The hidden line removal technique is necessary to visualize properly
3048 // very complex detectors. At the same time, it can be useful to visualize
3049 // the inner elements of a detector in detail. This function allows
3050 // subtractions (via boolean operation) of BOX shape from any part of
3051 // the detector, therefore showing its inner contents.
3052 // If "*" is given as the name of the
3053 // volume to be clipped, all volumes are clipped by the given box.
3054 // A volume can be clipped at most twice.
3055 // if a volume is explicitely clipped twice,
3056 // the "*" will not act on it anymore. Giving "." as the name
3057 // of the volume to be clipped will reset the clipping.
3059 // NAME Name of volume to be clipped
3061 // XMIN Lower limit of the Shape X coordinate
3062 // XMAX Upper limit of the Shape X coordinate
3063 // YMIN Lower limit of the Shape Y coordinate
3064 // YMAX Upper limit of the Shape Y coordinate
3065 // ZMIN Lower limit of the Shape Z coordinate
3066 // ZMAX Upper limit of the Shape Z coordinate
3068 // This function performs a boolean subtraction between the volume
3069 // NAME and a box placed in the MARS according the values of the given
3075 setclip(PASSCHARD(vname),xmin,xmax,ymin,ymax,zmin,zmax PASSCHARL(vname));
3078 //_____________________________________________________________________________
3079 void TGeant3::SetCOMP(Int_t par)
3082 // To control Compton scattering
3083 // par =0 no Compton
3084 // =1 Compton. Electron processed.
3085 // =2 Compton. No electron stored.
3088 fGcphys->icomp = par;
3091 //_____________________________________________________________________________
3092 void TGeant3::SetCUTS(Float_t cutgam,Float_t cutele,Float_t cutneu,
3093 Float_t cuthad,Float_t cutmuo ,Float_t bcute ,
3094 Float_t bcutm ,Float_t dcute ,Float_t dcutm ,
3095 Float_t ppcutm, Float_t tofmax)
3098 // CUTGAM Cut for gammas D=0.001
3099 // CUTELE Cut for electrons D=0.001
3100 // CUTHAD Cut for charged hadrons D=0.01
3101 // CUTNEU Cut for neutral hadrons D=0.01
3102 // CUTMUO Cut for muons D=0.01
3103 // BCUTE Cut for electron brems. D=-1.
3104 // BCUTM Cut for muon brems. D=-1.
3105 // DCUTE Cut for electron delta-rays D=-1.
3106 // DCUTM Cut for muon delta-rays D=-1.
3107 // PPCUTM Cut for e+e- pairs by muons D=0.01
3108 // TOFMAX Time of flight cut D=1.E+10
3110 // If the default values (-1.) for BCUTE ,BCUTM ,DCUTE ,DCUTM
3111 // are not modified, they will be set to CUTGAM,CUTGAM,CUTELE,CUTELE
3113 // If one of the parameters from CUTGAM to PPCUTM included
3114 // is modified, cross-sections and energy loss tables must be
3115 // recomputed via the function Gphysi.
3117 fGccuts->cutgam = cutgam;
3118 fGccuts->cutele = cutele;
3119 fGccuts->cutneu = cutneu;
3120 fGccuts->cuthad = cuthad;
3121 fGccuts->cutmuo = cutmuo;
3122 fGccuts->bcute = bcute;
3123 fGccuts->bcutm = bcutm;
3124 fGccuts->dcute = dcute;
3125 fGccuts->dcutm = dcutm;
3126 fGccuts->ppcutm = ppcutm;
3127 fGccuts->tofmax = tofmax;
3130 //_____________________________________________________________________________
3131 void TGeant3::SetDCAY(Int_t par)
3134 // To control Decay mechanism.
3135 // par =0 no decays.
3136 // =1 Decays. secondaries processed.
3137 // =2 Decays. No secondaries stored.
3139 fGcphys->idcay = par;
3143 //_____________________________________________________________________________
3144 void TGeant3::SetDEBU(Int_t emin, Int_t emax, Int_t emod)
3147 // Set the debug flag and frequency
3148 // Selected debug output will be printed from
3149 // event emin to even emax each emod event
3151 fGcflag->idemin = emin;
3152 fGcflag->idemax = emax;
3153 fGcflag->itest = emod;
3157 //_____________________________________________________________________________
3158 void TGeant3::SetDRAY(Int_t par)
3161 // To control delta rays mechanism.
3162 // par =0 no delta rays.
3163 // =1 Delta rays. secondaries processed.
3164 // =2 Delta rays. No secondaries stored.
3166 fGcphys->idray = par;
3169 //_____________________________________________________________________________
3170 void TGeant3::SetHADR(Int_t par)
3173 // To control hadronic interactions.
3174 // par =0 no hadronic interactions.
3175 // =1 Hadronic interactions. secondaries processed.
3176 // =2 Hadronic interactions. No secondaries stored.
3178 fGcphys->ihadr = par;
3181 //_____________________________________________________________________________
3182 void TGeant3::SetKINE(Int_t kine, Float_t xk1, Float_t xk2, Float_t xk3,
3183 Float_t xk4, Float_t xk5, Float_t xk6, Float_t xk7,
3184 Float_t xk8, Float_t xk9, Float_t xk10)
3187 // Set the variables in /GCFLAG/ IKINE, PKINE(10)
3188 // Their meaning is user defined
3190 fGckine->ikine = kine;
3191 fGckine->pkine[0] = xk1;
3192 fGckine->pkine[1] = xk2;
3193 fGckine->pkine[2] = xk3;
3194 fGckine->pkine[3] = xk4;
3195 fGckine->pkine[4] = xk5;
3196 fGckine->pkine[5] = xk6;
3197 fGckine->pkine[6] = xk7;
3198 fGckine->pkine[7] = xk8;
3199 fGckine->pkine[8] = xk9;
3200 fGckine->pkine[9] = xk10;
3203 //_____________________________________________________________________________
3204 void TGeant3::SetLOSS(Int_t par)
3207 // To control energy loss.
3208 // par =0 no energy loss;
3209 // =1 restricted energy loss fluctuations;
3210 // =2 complete energy loss fluctuations;
3212 // =4 no energy loss fluctuations.
3213 // If the value ILOSS is changed, then cross-sections and energy loss
3214 // tables must be recomputed via the command 'PHYSI'.
3216 fGcphys->iloss = par;
3220 //_____________________________________________________________________________
3221 void TGeant3::SetMULS(Int_t par)
3224 // To control multiple scattering.
3225 // par =0 no multiple scattering.
3226 // =1 Moliere or Coulomb scattering.
3227 // =2 Moliere or Coulomb scattering.
3228 // =3 Gaussian scattering.
3230 fGcphys->imuls = par;
3234 //_____________________________________________________________________________
3235 void TGeant3::SetMUNU(Int_t par)
3238 // To control muon nuclear interactions.
3239 // par =0 no muon-nuclear interactions.
3240 // =1 Nuclear interactions. Secondaries processed.
3241 // =2 Nuclear interactions. Secondaries not processed.
3243 fGcphys->imunu = par;
3246 //_____________________________________________________________________________
3247 void TGeant3::SetOPTI(Int_t par)
3250 // This flag controls the tracking optimisation performed via the
3252 // 1 no optimisation at all; GSORD calls disabled;
3253 // 0 no optimisation; only user calls to GSORD kept;
3254 // 1 all non-GSORDered volumes are ordered along the best axis;
3255 // 2 all volumes are ordered along the best axis.
3257 fGcopti->ioptim = par;
3260 //_____________________________________________________________________________
3261 void TGeant3::SetPAIR(Int_t par)
3264 // To control pair production mechanism.
3265 // par =0 no pair production.
3266 // =1 Pair production. secondaries processed.
3267 // =2 Pair production. No secondaries stored.
3269 fGcphys->ipair = par;
3273 //_____________________________________________________________________________
3274 void TGeant3::SetPFIS(Int_t par)
3277 // To control photo fission mechanism.
3278 // par =0 no photo fission.
3279 // =1 Photo fission. secondaries processed.
3280 // =2 Photo fission. No secondaries stored.
3282 fGcphys->ipfis = par;
3285 //_____________________________________________________________________________
3286 void TGeant3::SetPHOT(Int_t par)
3289 // To control Photo effect.
3290 // par =0 no photo electric effect.
3291 // =1 Photo effect. Electron processed.
3292 // =2 Photo effect. No electron stored.
3294 fGcphys->iphot = par;
3297 //_____________________________________________________________________________
3298 void TGeant3::SetRAYL(Int_t par)
3301 // To control Rayleigh scattering.
3302 // par =0 no Rayleigh scattering.
3305 fGcphys->irayl = par;
3308 //_____________________________________________________________________________
3309 void TGeant3::SetSWIT(Int_t sw, Int_t val)
3313 // val New switch value
3315 // Change one element of array ISWIT(10) in /GCFLAG/
3317 if (sw <= 0 || sw > 10) return;
3318 fGcflag->iswit[sw-1] = val;
3322 //_____________________________________________________________________________
3323 void TGeant3::SetTRIG(Int_t nevents)
3326 // Set number of events to be run
3328 fGcflag->nevent = nevents;
3331 //_____________________________________________________________________________
3332 void TGeant3::SetUserDecay(Int_t pdg)
3335 // Force the decays of particles to be done with Pythia
3336 // and not with the Geant routines.
3337 // just kill pointers doing mzdrop
3339 Int_t ipart = IdFromPDG(pdg);
3341 printf("Particle %d not in geant\n",pdg);
3344 Int_t jpart=fGclink->jpart;
3345 Int_t jpa=fZlq[jpart-ipart];
3348 Int_t jpa1=fZlq[jpa-1];
3350 mzdrop(fGcbank->ixcons,jpa1,PASSCHARD(" ") PASSCHARL(" "));
3351 Int_t jpa2=fZlq[jpa-2];
3353 mzdrop(fGcbank->ixcons,jpa2,PASSCHARD(" ") PASSCHARL(" "));
3357 //______________________________________________________________________________
3358 void TGeant3::Vname(const char *name, char *vname)
3361 // convert name to upper case. Make vname at least 4 chars
3363 Int_t l = strlen(name);
3366 for (i=0;i<l;i++) vname[i] = toupper(name[i]);
3367 for (i=l;i<4;i++) vname[i] = ' ';
3371 //______________________________________________________________________________
3372 void TGeant3::Ertrgo()
3377 //______________________________________________________________________________
3378 void TGeant3::Ertrak(const Float_t *const x1, const Float_t *const p1,
3379 const Float_t *x2, const Float_t *p2,
3380 Int_t ipa, Option_t *chopt)
3382 ertrak(x1,p1,x2,p2,ipa,PASSCHARD(chopt) PASSCHARL(chopt));
3385 //_____________________________________________________________________________
3386 void TGeant3::WriteEuclid(const char* filnam, const char* topvol,
3387 Int_t number, Int_t nlevel)
3391 // ******************************************************************
3393 // * Write out the geometry of the detector in EUCLID file format *
3395 // * filnam : will be with the extension .euc *
3396 // * topvol : volume name of the starting node *
3397 // * number : copy number of topvol (relevant for gsposp) *
3398 // * nlevel : number of levels in the tree structure *
3399 // * to be written out, starting from topvol *
3401 // * Author : M. Maire *
3403 // ******************************************************************
3405 // File filnam.tme is written out with the definitions of tracking
3406 // medias and materials.
3407 // As to restore original numbers for materials and medias, program
3408 // searches in the file euc_medi.dat and comparing main parameters of
3409 // the mat. defined inside geant and the one in file recognizes them
3410 // and is able to take number from file. If for any material or medium,
3411 // this procedure fails, ordering starts from 1.
3412 // Arrays IOTMED and IOMATE are used for this procedure
3414 const char shape[][5]={"BOX ","TRD1","TRD2","TRAP","TUBE","TUBS","CONE",
3415 "CONS","SPHE","PARA","PGON","PCON","ELTU","HYPE",
3417 Int_t i, end, itm, irm, jrm, k, nmed;
3421 char *filext, *filetme;
3422 char natmed[21], namate[21];
3423 char natmedc[21], namatec[21];
3424 char key[5], name[5], mother[5], konly[5];
3426 Int_t iadvol, iadtmd, iadrot, nwtot, iret;
3427 Int_t mlevel, numbr, natt, numed, nin, ndata;
3428 Int_t iname, ivo, ish, jvo, nvstak, ivstak;
3429 Int_t jdiv, ivin, in, jin, jvin, irot;
3430 Int_t jtm, imat, jma, flag=0, imatc;
3431 Float_t az, dens, radl, absl, a, step, x, y, z;
3432 Int_t npar, ndvmx, left;
3433 Float_t zc, densc, radlc, abslc, c0, tmaxfd;
3435 Int_t iomate[100], iotmed[100];
3436 Float_t par[50], att[20], ubuf[50];
3439 Int_t level, ndiv, iaxe;
3440 Int_t itmedc, nmatc, isvolc, ifieldc, nwbufc, isvol, nmat, ifield, nwbuf;
3441 Float_t fieldmc, tmaxfdc, stemaxc, deemaxc, epsilc, stminc, fieldm;
3442 Float_t tmaxf, stemax, deemax, epsil, stmin;
3443 const char *f10000="!\n%s\n!\n";
3444 //Open the input file
3446 for(i=0;i<end;i++) if(filnam[i]=='.') {
3450 filext=new char[end+4];
3451 filetme=new char[end+4];
3452 strncpy(filext,filnam,end);
3453 strncpy(filetme,filnam,end);
3455 // *** The output filnam name will be with extension '.euc'
3456 strcpy(&filext[end],".euc");
3457 strcpy(&filetme[end],".tme");
3458 lun=fopen(filext,"w");
3460 // *** Initialisation of the working space
3461 iadvol=fGcnum->nvolum;
3462 iadtmd=iadvol+fGcnum->nvolum;
3463 iadrot=iadtmd+fGcnum->ntmed;
3464 if(fGclink->jrotm) {
3465 fGcnum->nrotm=fZiq[fGclink->jrotm-2];
3469 nwtot=iadrot+fGcnum->nrotm;
3470 qws = new float[nwtot+1];
3471 for (i=0;i<nwtot+1;i++) qws[i]=0;
3474 if(nlevel==0) mlevel=20;
3476 // *** find the top volume and put it in the stak
3477 numbr = number>0 ? number : 1;
3478 Gfpara(topvol,numbr,1,npar,natt,par,att);
3480 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3485 // *** authorized shape ?
3486 strncpy((char *)&iname, topvol, 4);
3488 for(i=1; i<=fGcnum->nvolum; i++) if(fZiq[fGclink->jvolum+i]==iname) {
3492 jvo = fZlq[fGclink->jvolum-ivo];
3493 ish = Int_t (fZq[jvo+2]);
3495 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3502 iws[iadvol+ivo] = level;
3505 //*** flag all volumes and fill the stak
3509 // pick the next volume in stak
3511 ivo = TMath::Abs(iws[ivstak]);
3512 jvo = fZlq[fGclink->jvolum - ivo];
3514 // flag the tracking medium
3515 numed = Int_t (fZq[jvo + 4]);
3516 iws[iadtmd + numed] = 1;
3518 // get the daughters ...
3519 level = iws[iadvol+ivo];
3520 if (level < mlevel) {
3522 nin = Int_t (fZq[jvo + 3]);
3524 // from division ...
3526 jdiv = fZlq[jvo - 1];
3527 ivin = Int_t (fZq[jdiv + 2]);
3529 iws[nvstak] = -ivin;
3530 iws[iadvol+ivin] = level;
3532 // from position ...
3533 } else if (nin > 0) {
3534 for(in=1; in<=nin; in++) {
3535 jin = fZlq[jvo - in];
3536 ivin = Int_t (fZq[jin + 2 ]);
3537 jvin = fZlq[fGclink->jvolum - ivin];
3538 ish = Int_t (fZq[jvin + 2]);
3539 // authorized shape ?
3541 // not yet flagged ?
3542 if (iws[iadvol+ivin]==0) {
3545 iws[iadvol+ivin] = level;
3547 // flag the rotation matrix
3548 irot = Int_t ( fZq[jin + 4 ]);
3549 if (irot > 0) iws[iadrot+irot] = 1;
3555 // next volume in stak ?
3556 if (ivstak < nvstak) goto L10;
3558 // *** restore original material and media numbers
3559 // file euc_medi.dat is needed to compare materials and medias
3561 FILE* luncor=fopen("euc_medi.dat","r");
3564 for(itm=1; itm<=fGcnum->ntmed; itm++) {
3565 if (iws[iadtmd+itm] > 0) {
3566 jtm = fZlq[fGclink->jtmed-itm];
3567 strncpy(natmed,(char *)&fZiq[jtm+1],20);
3568 imat = Int_t (fZq[jtm+6]);
3569 jma = fZlq[fGclink->jmate-imat];
3571 printf(" *** GWEUCL *** material not defined for tracking medium %5i %s\n",itm,natmed);
3574 strncpy(namate,(char *)&fZiq[jma+1],20);
3577 //** find the material original number
3580 iret=fscanf(luncor,"%4s,%130s",key,card);
3581 if(iret<=0) goto L26;
3583 if(!strcmp(key,"MATE")) {
3584 sscanf(card,"%d %s %f %f %f %f %f %d",&imatc,namatec,&az,&zc,&densc,&radlc,&abslc,&nparc);
3585 Gfmate(imat,namate,a,z,dens,radl,absl,par,npar);
3586 if(!strcmp(namatec,namate)) {
3587 if(az==a && zc==z && densc==dens && radlc==radl
3588 && abslc==absl && nparc==nparc) {
3591 printf("*** GWEUCL *** material : %3d '%s' restored as %3d\n",imat,namate,imatc);
3593 printf("*** GWEUCL *** different definitions for material: %s\n",namate);
3597 if(strcmp(key,"END") && !flag) goto L23;
3599 printf("*** GWEUCL *** cannot restore original number for material: %s\n",namate);
3603 //*** restore original tracking medium number
3606 iret=fscanf(luncor,"%4s,%130s",key,card);
3607 if(iret<=0) goto L26;
3609 if (!strcmp(key,"TMED")) {
3610 sscanf(card,"%d %s %d %d %d %f %f %f %f %f %f %d\n",
3611 &itmedc,natmedc,&nmatc,&isvolc,&ifieldc,&fieldmc,
3612 &tmaxfdc,&stemaxc,&deemaxc,&epsilc,&stminc,&nwbufc);
3613 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxf,stemax,deemax,
3614 epsil,stmin,ubuf,&nwbuf);
3615 if(!strcmp(natmedc,natmed)) {
3616 if (iomate[nmat]==nmatc && nwbuf==nwbufc) {
3619 printf("*** GWEUCL *** medium : %3d '%20s' restored as %3d\n",
3622 printf("*** GWEUCL *** different definitions for tracking medium: %s\n",natmed);
3626 if(strcmp(key,"END") && !flag) goto L24;
3628 printf("cannot restore original number for medium : %s\n",natmed);
3636 L26: printf("*** GWEUCL *** cannot read the data file\n");
3638 L29: if(luncor) fclose (luncor);
3641 // *** write down the tracking medium definition
3643 strcpy(card,"! Tracking medium");
3644 fprintf(lun,f10000,card);
3646 for(itm=1;itm<=fGcnum->ntmed;itm++) {
3647 if (iws[iadtmd+itm]>0) {
3648 jtm = fZlq[fGclink->jtmed-itm];
3649 strncpy(natmed,(char *)&fZiq[jtm+1],20);
3651 imat = Int_t (fZq[jtm+6]);
3652 jma = fZlq[fGclink->jmate-imat];
3653 //* order media from one, if comparing with database failed
3655 iotmed[itm]=++imxtmed;
3656 iomate[imat]=++imxmate;
3661 printf(" *** GWEUCL *** material not defined for tracking medium %5d %s\n",
3664 strncpy(namate,(char *)&fZiq[jma+1],20);
3667 fprintf(lun,"TMED %3d '%20s' %3d '%20s'\n",iotmed[itm],natmed,iomate[imat],namate);
3671 //* *** write down the rotation matrix
3673 strcpy(card,"! Reperes");
3674 fprintf(lun,f10000,card);
3676 for(irm=1;irm<=fGcnum->nrotm;irm++) {
3677 if (iws[iadrot+irm]>0) {
3678 jrm = fZlq[fGclink->jrotm-irm];
3679 fprintf(lun,"ROTM %3d",irm);
3680 for(k=11;k<=16;k++) fprintf(lun," %8.3f",fZq[jrm+k]);
3685 //* *** write down the volume definition
3687 strcpy(card,"! Volumes");
3688 fprintf(lun,f10000,card);
3690 for(ivstak=1;ivstak<=nvstak;ivstak++) {
3693 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivo],4);
3695 jvo = fZlq[fGclink->jvolum-ivo];
3696 ish = Int_t (fZq[jvo+2]);
3697 nmed = Int_t (fZq[jvo+4]);
3698 npar = Int_t (fZq[jvo+5]);
3700 if (ivstak>1) for(i=0;i<npar;i++) par[i]=fZq[jvo+7+i];
3701 Gckpar (ish,npar,par);
3702 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,shape[ish-1],iotmed[nmed],npar);
3703 for(i=0;i<(npar-1)/6+1;i++) {
3706 for(k=0;k<(left<6?left:6);k++) fprintf(lun," %11.5f",par[i*6+k]);
3710 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,shape[ish-1],iotmed[nmed],npar);
3715 //* *** write down the division of volumes
3717 fprintf(lun,f10000,"! Divisions");
3718 for(ivstak=1;ivstak<=nvstak;ivstak++) {
3719 ivo = TMath::Abs(iws[ivstak]);
3720 jvo = fZlq[fGclink->jvolum-ivo];
3721 ish = Int_t (fZq[jvo+2]);
3722 nin = Int_t (fZq[jvo+3]);
3723 //* this volume is divided ...
3726 iaxe = Int_t ( fZq[jdiv+1]);
3727 ivin = Int_t ( fZq[jdiv+2]);
3728 ndiv = Int_t ( fZq[jdiv+3]);
3731 jvin = fZlq[fGclink->jvolum-ivin];
3732 nmed = Int_t ( fZq[jvin+4]);
3733 strncpy(mother,(char *)&fZiq[fGclink->jvolum+ivo ],4);
3735 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivin],4);
3737 if ((step<=0.)||(ish>=11)) {
3738 //* volume with negative parameter or gsposp or pgon ...
3739 fprintf(lun,"DIVN '%4s' '%4s' %3d %3d\n",name,mother,ndiv,iaxe);
3740 } else if ((ndiv<=0)||(ish==10)) {
3741 //* volume with negative parameter or gsposp or para ...
3742 ndvmx = TMath::Abs(ndiv);
3743 fprintf(lun,"DIVT '%4s' '%4s' %11.5f %3d %3d %3d\n",
3744 name,mother,step,iaxe,iotmed[nmed],ndvmx);
3746 //* normal volume : all kind of division are equivalent
3747 fprintf(lun,"DVT2 '%4s' '%4s' %11.5f %3d %11.5f %3d %3d\n",
3748 name,mother,step,iaxe,c0,iotmed[nmed],ndiv);
3753 //* *** write down the the positionnement of volumes
3755 fprintf(lun,f10000,"! Positionnements\n");
3757 for(ivstak = 1;ivstak<=nvstak;ivstak++) {
3758 ivo = TMath::Abs(iws[ivstak]);
3759 strncpy(mother,(char*)&fZiq[fGclink->jvolum+ivo ],4);
3761 jvo = fZlq[fGclink->jvolum-ivo];
3762 nin = Int_t( fZq[jvo+3]);
3763 //* this volume has daughters ...
3765 for (in=1;in<=nin;in++) {
3767 ivin = Int_t (fZq[jin +2]);
3768 numb = Int_t (fZq[jin +3]);
3769 irot = Int_t (fZq[jin +4]);
3773 strcpy(konly,"ONLY");
3774 if (fZq[jin+8]!=1.) strcpy(konly,"MANY");
3775 strncpy(name,(char*)&fZiq[fGclink->jvolum+ivin],4);
3777 jvin = fZlq[fGclink->jvolum-ivin];
3778 ish = Int_t (fZq[jvin+2]);
3779 //* gspos or gsposp ?
3780 ndata = fZiq[jin-1];
3782 fprintf(lun,"POSI '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s'\n",
3783 name,numb,mother,x,y,z,irot,konly);
3785 npar = Int_t (fZq[jin+9]);
3786 for(i=0;i<npar;i++) par[i]=fZq[jin+10+i];
3787 Gckpar (ish,npar,par);
3788 fprintf(lun,"POSP '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s' %3d\n",
3789 name,numb,mother,x,y,z,irot,konly,npar);
3791 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
3798 fprintf(lun,"END\n");
3801 //****** write down the materials and medias *****
3803 lun=fopen(filetme,"w");
3805 for(itm=1;itm<=fGcnum->ntmed;itm++) {
3806 if (iws[iadtmd+itm]>0) {
3807 jtm = fZlq[fGclink->jtmed-itm];
3808 strncpy(natmed,(char*)&fZiq[jtm+1],4);
3809 imat = Int_t (fZq[jtm+6]);
3810 jma = Int_t (fZlq[fGclink->jmate-imat]);
3812 Gfmate (imat,namate,a,z,dens,radl,absl,par,npar);
3813 fprintf(lun,"MATE %4d '%20s'%11.5E %11.5E %11.5E %11.5E %11.5E %3d\n",
3814 iomate[imat],namate,a,z,dens,radl,absl,npar);
3818 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
3822 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxfd,stemax,deemax,epsil,stmin,par,&npar);
3823 fprintf(lun,"TMED %4d '%20s' %3d %1d %3d %11.5f %11.5f %11.5f %11.5f %11.5f %11.5f %3d\n",
3824 iotmed[itm],natmed,iomate[nmat],isvol,ifield,
3825 fieldm,tmaxfd,stemax,deemax,epsil,stmin,npar);
3829 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
3835 fprintf(lun,"END\n");
3836 printf(" *** GWEUCL *** file: %s is now written out\n",filext);
3837 printf(" *** GWEUCL *** file: %s is now written out\n",filetme);
3846 //_____________________________________________________________________________
3847 void TGeant3::Streamer(TBuffer &R__b)
3850 // Stream an object of class TGeant3.
3852 if (R__b.IsReading()) {
3853 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
3854 AliMC::Streamer(R__b);
3857 R__b.ReadStaticArray(fPDGCode);
3859 R__b.WriteVersion(TGeant3::IsA());
3860 AliMC::Streamer(R__b);
3863 R__b.WriteArray(fPDGCode, fNPDGCodes);