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.17 1999/11/03 16:58:28 fca
19 Correct source of address violation in creating character string
21 Revision 1.16 1999/11/03 13:17:08 fca
22 Have ProdProcess return const char*
24 Revision 1.15 1999/10/26 06:04:50 fca
25 Introduce TLorentzVector in AliMC::GetSecondary. Thanks to I.Hrivnacova
27 Revision 1.14 1999/09/29 09:24:30 fca
28 Introduction of the Copyright and cvs Log
32 ///////////////////////////////////////////////////////////////////////////////
34 // Interface Class to the Geant3.21 MonteCarlo //
38 <img src="picts/TGeant3Class.gif">
43 ///////////////////////////////////////////////////////////////////////////////
49 #include <TDatabasePDG.h>
50 #include "AliCallf77.h"
53 # define gzebra gzebra_
54 # define grfile grfile_
55 # define gpcxyz gpcxyz_
56 # define ggclos ggclos_
59 # define gcinit gcinit_
62 # define gtrigc gtrigc_
63 # define gtrigi gtrigi_
65 # define gzinit gzinit_
66 # define gfmate gfmate_
67 # define gfpart gfpart_
68 # define gftmed gftmed_
72 # define gsmate gsmate_
73 # define gsmixt gsmixt_
74 # define gspart gspart_
75 # define gstmed gstmed_
76 # define gsckov gsckov_
77 # define gstpar gstpar_
78 # define gfkine gfkine_
79 # define gfvert gfvert_
80 # define gskine gskine_
81 # define gsvert gsvert_
82 # define gphysi gphysi_
83 # define gdebug gdebug_
84 # define gekbin gekbin_
85 # define gfinds gfinds_
86 # define gsking gsking_
87 # define gskpho gskpho_
88 # define gsstak gsstak_
90 # define gtrack gtrack_
91 # define gtreve gtreve_
92 # define gtreve_root gtreve_root_
94 # define grndmq grndmq_
96 # define glmoth glmoth_
97 # define gmedia gmedia_
100 # define gsdvn2 gsdvn2_
101 # define gsdvs gsdvs_
102 # define gsdvs2 gsdvs2_
103 # define gsdvt gsdvt_
104 # define gsdvt2 gsdvt2_
105 # define gsord gsord_
106 # define gspos gspos_
107 # define gsposp gsposp_
108 # define gsrotm gsrotm_
109 # define gprotm gprotm_
110 # define gsvolu gsvolu_
111 # define gprint gprint_
112 # define gdinit gdinit_
113 # define gdopt gdopt_
114 # define gdraw gdraw_
115 # define gdrayt gdrayt_
116 # define gdrawc gdrawc_
117 # define gdrawx gdrawx_
118 # define gdhead gdhead_
119 # define gdwmn1 gdwmn1_
120 # define gdwmn2 gdwmn2_
121 # define gdwmn3 gdwmn3_
122 # define gdxyz gdxyz_
123 # define gdcxyz gdcxyz_
124 # define gdman gdman_
125 # define gdspec gdspec_
126 # define gdtree gdtree_
127 # define gdelet gdelet_
128 # define gdclos gdclos_
129 # define gdshow gdshow_
130 # define gdopen gdopen_
131 # define dzshow dzshow_
132 # define gsatt gsatt_
133 # define gfpara gfpara_
134 # define gckpar gckpar_
135 # define gckmat gckmat_
136 # define geditv geditv_
137 # define mzdrop mzdrop_
139 # define ertrak ertrak_
140 # define ertrgo ertrgo_
142 # define setbomb setbomb_
143 # define setclip setclip_
144 # define gcomad gcomad_
147 # define gzebra GZEBRA
148 # define grfile GRFILE
149 # define gpcxyz GPCXYZ
150 # define ggclos GGCLOS
153 # define gcinit GCINIT
156 # define gtrigc GTRIGC
157 # define gtrigi GTRIGI
159 # define gzinit GZINIT
160 # define gfmate GFMATE
161 # define gfpart GFPART
162 # define gftmed GFTMED
166 # define gsmate GSMATE
167 # define gsmixt GSMIXT
168 # define gspart GSPART
169 # define gstmed GSTMED
170 # define gsckov GSCKOV
171 # define gstpar GSTPAR
172 # define gfkine GFKINE
173 # define gfvert GFVERT
174 # define gskine GSKINE
175 # define gsvert GSVERT
176 # define gphysi GPHYSI
177 # define gdebug GDEBUG
178 # define gekbin GEKBIN
179 # define gfinds GFINDS
180 # define gsking GSKING
181 # define gskpho GSKPHO
182 # define gsstak GSSTAK
184 # define gtrack GTRACK
185 # define gtreve GTREVE
186 # define gtreve_root GTREVE_ROOT
188 # define grndmq GRNDMQ
190 # define glmoth GLMOTH
191 # define gmedia GMEDIA
194 # define gsdvn2 GSDVN2
196 # define gsdvs2 GSDVS2
198 # define gsdvt2 GSDVT2
201 # define gsposp GSPOSP
202 # define gsrotm GSROTM
203 # define gprotm GPROTM
204 # define gsvolu GSVOLU
205 # define gprint GPRINT
206 # define gdinit GDINIT
209 # define gdrayt GDRAYT
210 # define gdrawc GDRAWC
211 # define gdrawx GDRAWX
212 # define gdhead GDHEAD
213 # define gdwmn1 GDWMN1
214 # define gdwmn2 GDWMN2
215 # define gdwmn3 GDWMN3
217 # define gdcxyz GDCXYZ
219 # define gdfspc GDFSPC
220 # define gdspec GDSPEC
221 # define gdtree GDTREE
222 # define gdelet GDELET
223 # define gdclos GDCLOS
224 # define gdshow GDSHOW
225 # define gdopen GDOPEN
226 # define dzshow DZSHOW
228 # define gfpara GFPARA
229 # define gckpar GCKPAR
230 # define gckmat GCKMAT
231 # define geditv GEDITV
232 # define mzdrop MZDROP
234 # define ertrak ERTRAK
235 # define ertrgo ERTRGO
237 # define setbomb SETBOMB
238 # define setclip SETCLIP
239 # define gcomad GCOMAD
243 //____________________________________________________________________________
247 // Prototypes for GEANT functions
249 void type_of_call gzebra(const int&);
251 void type_of_call gpcxyz();
253 void type_of_call ggclos();
255 void type_of_call glast();
257 void type_of_call ginit();
259 void type_of_call gcinit();
261 void type_of_call grun();
263 void type_of_call gtrig();
265 void type_of_call gtrigc();
267 void type_of_call gtrigi();
269 void type_of_call gwork(const int&);
271 void type_of_call gzinit();
273 void type_of_call gmate();
275 void type_of_call gpart();
277 void type_of_call gsdk(Int_t &, Float_t *, Int_t *);
279 void type_of_call gfkine(Int_t &, Float_t *, Float_t *, Int_t &,
280 Int_t &, Float_t *, Int_t &);
282 void type_of_call gfvert(Int_t &, Float_t *, Int_t &, Int_t &,
283 Float_t &, Float_t *, Int_t &);
285 void type_of_call gskine(Float_t *,Int_t &, Int_t &, Float_t *,
288 void type_of_call gsvert(Float_t *,Int_t &, Int_t &, Float_t *,
291 void type_of_call gphysi();
293 void type_of_call gdebug();
295 void type_of_call gekbin();
297 void type_of_call gfinds();
299 void type_of_call gsking(Int_t &);
301 void type_of_call gskpho(Int_t &);
303 void type_of_call gsstak(Int_t &);
305 void type_of_call gsxyz();
307 void type_of_call gtrack();
309 void type_of_call gtreve();
311 void type_of_call gtreve_root();
313 void type_of_call grndm(Float_t *, const Int_t &);
315 void type_of_call grndmq(Int_t &, Int_t &, const Int_t &,
318 void type_of_call gdtom(Float_t *, Float_t *, Int_t &);
320 void type_of_call glmoth(DEFCHARD, Int_t &, Int_t &, Int_t *,
321 Int_t *, Int_t * DEFCHARL);
323 void type_of_call gmedia(Float_t *, Int_t &);
325 void type_of_call gmtod(Float_t *, Float_t *, Int_t &);
327 void type_of_call gsrotm(const Int_t &, const Float_t &, const Float_t &,
328 const Float_t &, const Float_t &, const Float_t &,
331 void type_of_call gprotm(const Int_t &);
333 void type_of_call grfile(const Int_t&, DEFCHARD,
334 DEFCHARD DEFCHARL DEFCHARL);
336 void type_of_call gfmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
337 Float_t &, Float_t &, Float_t &, Float_t *,
340 void type_of_call gfpart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
341 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
343 void type_of_call gftmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
344 Float_t &, Float_t &, Float_t &, Float_t &,
345 Float_t &, Float_t &, Float_t *, Int_t * DEFCHARL);
347 void type_of_call gsmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
348 Float_t &, Float_t &, Float_t &, Float_t *,
351 void type_of_call gsmixt(const Int_t&, DEFCHARD, Float_t *, Float_t *,
352 Float_t &, Int_t &, Float_t * DEFCHARL);
354 void type_of_call gspart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
355 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
358 void type_of_call gstmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
359 Float_t &, Float_t &, Float_t &, Float_t &,
360 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
362 void type_of_call gsckov(Int_t &itmed, Int_t &npckov, Float_t *ppckov,
363 Float_t *absco, Float_t *effic, Float_t *rindex);
364 void type_of_call gstpar(const Int_t&, DEFCHARD, Float_t & DEFCHARL);
366 void type_of_call gsdvn(DEFCHARD,DEFCHARD, Int_t &, Int_t &
369 void type_of_call gsdvn2(DEFCHARD,DEFCHARD, Int_t &, Int_t &, Float_t &,
370 Int_t & DEFCHARL DEFCHARL);
372 void type_of_call gsdvs(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &
375 void type_of_call gsdvs2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t &,
376 Int_t & DEFCHARL DEFCHARL);
378 void type_of_call gsdvt(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &,
379 Int_t & DEFCHARL DEFCHARL);
381 void type_of_call gsdvt2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t&,
382 Int_t &, Int_t & DEFCHARL DEFCHARL);
384 void type_of_call gsord(DEFCHARD, Int_t & DEFCHARL);
386 void type_of_call gspos(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
387 Float_t &, Int_t &, DEFCHARD DEFCHARL DEFCHARL
390 void type_of_call gsposp(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
391 Float_t &, Int_t &, DEFCHARD,
392 Float_t *, Int_t & DEFCHARL DEFCHARL DEFCHARL);
394 void type_of_call gsvolu(DEFCHARD, DEFCHARD, Int_t &, Float_t *, Int_t &,
395 Int_t & DEFCHARL DEFCHARL);
397 void type_of_call gsatt(DEFCHARD, DEFCHARD, Int_t & DEFCHARL DEFCHARL);
399 void type_of_call gfpara(DEFCHARD , Int_t&, Int_t&, Int_t&, Int_t&, Float_t*,
402 void type_of_call gckpar(Int_t&, Int_t&, Float_t*);
404 void type_of_call gckmat(Int_t&, DEFCHARD DEFCHARL);
406 void type_of_call gprint(DEFCHARD,const int& DEFCHARL);
408 void type_of_call gdinit();
410 void type_of_call gdopt(DEFCHARD,DEFCHARD DEFCHARL DEFCHARL);
412 void type_of_call gdraw(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
413 Float_t &, Float_t &, Float_t & DEFCHARL);
414 void type_of_call gdrayt(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
415 Float_t &, Float_t &, Float_t & DEFCHARL);
416 void type_of_call gdrawc(DEFCHARD,Int_t &, Float_t &, Float_t &, Float_t &,
417 Float_t &, Float_t & DEFCHARL);
418 void type_of_call gdrawx(DEFCHARD,Float_t &, Float_t &, Float_t &, Float_t &,
419 Float_t &, Float_t &, Float_t &, Float_t &,
421 void type_of_call gdhead(Int_t &,DEFCHARD, Float_t & DEFCHARL);
422 void type_of_call gdxyz(Int_t &);
423 void type_of_call gdcxyz();
424 void type_of_call gdman(Float_t &, Float_t &);
425 void type_of_call gdwmn1(Float_t &, Float_t &);
426 void type_of_call gdwmn2(Float_t &, Float_t &);
427 void type_of_call gdwmn3(Float_t &, Float_t &);
428 void type_of_call gdspec(DEFCHARD DEFCHARL);
429 void type_of_call gdfspc(DEFCHARD, Int_t &, Int_t & DEFCHARL) {;}
430 void type_of_call gdtree(DEFCHARD, Int_t &, Int_t & DEFCHARL);
432 void type_of_call gdopen(Int_t &);
433 void type_of_call gdclos();
434 void type_of_call gdelet(Int_t &);
435 void type_of_call gdshow(Int_t &);
436 void type_of_call geditv(Int_t &) {;}
439 void type_of_call dzshow(DEFCHARD,const int&,const int&,DEFCHARD,const int&,
440 const int&, const int&, const int& DEFCHARL
443 void type_of_call mzdrop(Int_t&, Int_t&, DEFCHARD DEFCHARL);
445 void type_of_call setbomb(Float_t &);
446 void type_of_call setclip(DEFCHARD, Float_t &,Float_t &,Float_t &,Float_t &,
447 Float_t &, Float_t & DEFCHARL);
448 void type_of_call gcomad(DEFCHARD, Int_t*& DEFCHARL);
450 void type_of_call ertrak(const Float_t *const x1, const Float_t *const p1,
451 const Float_t *x2, const Float_t *p2,
452 const Int_t &ipa, DEFCHARD DEFCHARL);
454 void type_of_call ertrgo();
458 // Geant3 global pointer
460 static Int_t defSize = 600;
464 //____________________________________________________________________________
468 // Default constructor
472 //____________________________________________________________________________
473 TGeant3::TGeant3(const char *title, Int_t nwgeant)
474 :AliMC("TGeant3",title)
477 // Standard constructor for TGeant3 with ZEBRA initialisation
488 // Load Address of Geant3 commons
491 // Zero number of particles
495 //____________________________________________________________________________
496 Int_t TGeant3::CurrentMaterial(Float_t &a, Float_t &z, Float_t &dens,
497 Float_t &radl, Float_t &absl) const
500 // Return the parameters of the current material during transport
504 dens = fGcmate->dens;
505 radl = fGcmate->radl;
506 absl = fGcmate->absl;
507 return 1; //this could be the number of elements in mixture
510 //____________________________________________________________________________
511 void TGeant3::DefaultRange()
514 // Set range of current drawing pad to 20x20 cm
520 higz->Range(0,0,20,20);
523 //____________________________________________________________________________
524 void TGeant3::InitHIGZ()
535 //____________________________________________________________________________
536 void TGeant3::LoadAddress()
539 // Assigns the address of the GEANT common blocks to the structures
540 // that allow their access from C++
543 gcomad(PASSCHARD("QUEST"), (int*&) fQuest PASSCHARL("QUEST"));
544 gcomad(PASSCHARD("GCBANK"),(int*&) fGcbank PASSCHARL("GCBANK"));
545 gcomad(PASSCHARD("GCLINK"),(int*&) fGclink PASSCHARL("GCLINK"));
546 gcomad(PASSCHARD("GCCUTS"),(int*&) fGccuts PASSCHARL("GCCUTS"));
547 gcomad(PASSCHARD("GCFLAG"),(int*&) fGcflag PASSCHARL("GCFLAG"));
548 gcomad(PASSCHARD("GCKINE"),(int*&) fGckine PASSCHARL("GCKINE"));
549 gcomad(PASSCHARD("GCKING"),(int*&) fGcking PASSCHARL("GCKING"));
550 gcomad(PASSCHARD("GCKIN2"),(int*&) fGckin2 PASSCHARL("GCKIN2"));
551 gcomad(PASSCHARD("GCKIN3"),(int*&) fGckin3 PASSCHARL("GCKIN3"));
552 gcomad(PASSCHARD("GCMATE"),(int*&) fGcmate PASSCHARL("GCMATE"));
553 gcomad(PASSCHARD("GCTMED"),(int*&) fGctmed PASSCHARL("GCTMED"));
554 gcomad(PASSCHARD("GCTRAK"),(int*&) fGctrak PASSCHARL("GCTRAK"));
555 gcomad(PASSCHARD("GCTPOL"),(int*&) fGctpol PASSCHARL("GCTPOL"));
556 gcomad(PASSCHARD("GCVOLU"),(int*&) fGcvolu PASSCHARL("GCVOLU"));
557 gcomad(PASSCHARD("GCNUM"), (int*&) fGcnum PASSCHARL("GCNUM"));
558 gcomad(PASSCHARD("GCSETS"),(int*&) fGcsets PASSCHARL("GCSETS"));
559 gcomad(PASSCHARD("GCPHYS"),(int*&) fGcphys PASSCHARL("GCPHYS"));
560 gcomad(PASSCHARD("GCOPTI"),(int*&) fGcopti PASSCHARL("GCOPTI"));
561 gcomad(PASSCHARD("GCTLIT"),(int*&) fGctlit PASSCHARL("GCTLIT"));
562 gcomad(PASSCHARD("GCVDMA"),(int*&) fGcvdma PASSCHARL("GCVDMA"));
565 gcomad(PASSCHARD("ERTRIO"),(int*&) fErtrio PASSCHARL("ERTRIO"));
566 gcomad(PASSCHARD("EROPTS"),(int*&) fEropts PASSCHARL("EROPTS"));
567 gcomad(PASSCHARD("EROPTC"),(int*&) fEroptc PASSCHARL("EROPTC"));
568 gcomad(PASSCHARD("ERWORK"),(int*&) fErwork PASSCHARL("ERWORK"));
570 // Variables for ZEBRA store
571 gcomad(PASSCHARD("IQ"), addr PASSCHARL("IQ"));
573 gcomad(PASSCHARD("LQ"), addr PASSCHARL("LQ"));
578 //_____________________________________________________________________________
579 void TGeant3::GeomIter()
582 // Geometry iterator for moving upward in the geometry tree
583 // Initialise the iterator
585 fNextVol=fGcvolu->nlevel;
588 //____________________________________________________________________________
589 Int_t TGeant3::NextVolUp(Text_t *name, Int_t ©)
592 // Geometry iterator for moving upward in the geometry tree
593 // Return next volume up
598 gname=fGcvolu->names[fNextVol];
599 strncpy(name,(char *) &gname, 4);
601 copy=fGcvolu->number[fNextVol];
602 i=fGcvolu->lvolum[fNextVol];
603 if(gname == fZiq[fGclink->jvolum+i]) return i;
604 else printf("GeomTree: Volume %s not found in bank\n",name);
609 //_____________________________________________________________________________
610 Int_t TGeant3::CurrentVolID(Int_t ©) const
613 // Returns the current volume ID and copy number
616 if( (i=fGcvolu->nlevel-1) < 0 ) {
617 Warning("CurrentVolID","Stack depth only %d\n",fGcvolu->nlevel);
619 gname=fGcvolu->names[i];
620 copy=fGcvolu->number[i];
621 i=fGcvolu->lvolum[i];
622 if(gname == fZiq[fGclink->jvolum+i]) return i;
623 else Warning("CurrentVolID","Volume %4s not found\n",(char*)&gname);
628 //_____________________________________________________________________________
629 Int_t TGeant3::CurrentVolOffID(Int_t off, Int_t ©) const
632 // Return the current volume "off" upward in the geometrical tree
633 // ID and copy number
636 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
637 Warning("CurrentVolOffID","Offset requested %d but stack depth %d\n",
638 off,fGcvolu->nlevel);
640 gname=fGcvolu->names[i];
641 copy=fGcvolu->number[i];
642 i=fGcvolu->lvolum[i];
643 if(gname == fZiq[fGclink->jvolum+i]) return i;
644 else Warning("CurrentVolOffID","Volume %4s not found\n",(char*)&gname);
649 //_____________________________________________________________________________
650 const char* TGeant3::CurrentVolName() const
653 // Returns the current volume name
657 if( (i=fGcvolu->nlevel-1) < 0 ) {
658 Warning("CurrentVolName","Stack depth %d\n",fGcvolu->nlevel);
660 gname=fGcvolu->names[i];
661 strncpy(name,(char *) &gname, 4);
663 i=fGcvolu->lvolum[i];
664 if(gname == fZiq[fGclink->jvolum+i]) return name;
665 else Warning("CurrentVolName","Volume %4s not found\n",name);
670 //_____________________________________________________________________________
671 const char* TGeant3::CurrentVolOffName(Int_t off) const
674 // Return the current volume "off" upward in the geometrical tree
675 // ID, name and copy number
676 // if name=0 no name is returned
680 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
681 Warning("CurrentVolOffName",
682 "Offset requested %d but stack depth %d\n",off,fGcvolu->nlevel);
684 gname=fGcvolu->names[i];
686 strncpy(name,(char *) &gname, 4);
688 i=fGcvolu->lvolum[i];
689 if(gname == fZiq[fGclink->jvolum+i]) return name;
690 else Warning("CurrentVolOffName","Volume %4s not found\n",name);
695 //_____________________________________________________________________________
696 Int_t TGeant3::IdFromPDG(Int_t pdg) const
699 // Return Geant3 code from PDG and pseudo ENDF code
701 for(Int_t i=0;i<fNPDGCodes;++i)
702 if(pdg==fPDGCode[i]) return i;
706 //_____________________________________________________________________________
707 Int_t TGeant3::PDGFromId(Int_t id) const
709 if(id>0 && id<fNPDGCodes) return fPDGCode[id];
713 //_____________________________________________________________________________
714 void TGeant3::DefineParticles()
717 // Define standard Geant 3 particles
720 // Load standard numbers for GEANT particles and PDG conversion
721 fPDGCode[fNPDGCodes++]=-99; // 0 = unused location
722 fPDGCode[fNPDGCodes++]=22; // 1 = photon
723 fPDGCode[fNPDGCodes++]=-11; // 2 = positron
724 fPDGCode[fNPDGCodes++]=11; // 3 = electron
725 fPDGCode[fNPDGCodes++]=12; // 4 = neutrino e
726 fPDGCode[fNPDGCodes++]=-13; // 5 = muon +
727 fPDGCode[fNPDGCodes++]=13; // 6 = muon -
728 fPDGCode[fNPDGCodes++]=111; // 7 = pi0
729 fPDGCode[fNPDGCodes++]=211; // 8 = pi+
730 fPDGCode[fNPDGCodes++]=-211; // 9 = pi-
731 fPDGCode[fNPDGCodes++]=130; // 10 = Kaon Long
732 fPDGCode[fNPDGCodes++]=321; // 11 = Kaon +
733 fPDGCode[fNPDGCodes++]=-321; // 12 = Kaon -
734 fPDGCode[fNPDGCodes++]=2112; // 13 = Neutron
735 fPDGCode[fNPDGCodes++]=2212; // 14 = Proton
736 fPDGCode[fNPDGCodes++]=-2212; // 15 = Anti Proton
737 fPDGCode[fNPDGCodes++]=310; // 16 = Kaon Short
738 fPDGCode[fNPDGCodes++]=221; // 17 = Eta
739 fPDGCode[fNPDGCodes++]=3122; // 18 = Lambda
740 fPDGCode[fNPDGCodes++]=3222; // 19 = Sigma +
741 fPDGCode[fNPDGCodes++]=3212; // 20 = Sigma 0
742 fPDGCode[fNPDGCodes++]=3112; // 21 = Sigma -
743 fPDGCode[fNPDGCodes++]=3322; // 22 = Xi0
744 fPDGCode[fNPDGCodes++]=3312; // 23 = Xi-
745 fPDGCode[fNPDGCodes++]=3334; // 24 = Omega-
746 fPDGCode[fNPDGCodes++]=-2112; // 25 = Anti Proton
747 fPDGCode[fNPDGCodes++]=-3122; // 26 = Anti Proton
748 fPDGCode[fNPDGCodes++]=-3222; // 27 = Anti Sigma -
749 fPDGCode[fNPDGCodes++]=-3212; // 28 = Anti Sigma 0
750 fPDGCode[fNPDGCodes++]=-3112; // 29 = Anti Sigma 0
751 fPDGCode[fNPDGCodes++]=-3322; // 30 = Anti Xi 0
752 fPDGCode[fNPDGCodes++]=-3312; // 31 = Anti Xi +
753 fPDGCode[fNPDGCodes++]=-3334; // 32 = Anti Omega +
760 /* --- Define additional particles */
761 Gspart(33, "OMEGA(782)", 3, 0.782, 0., 7.836e-23);
762 fPDGCode[fNPDGCodes++]=223; // 33 = Omega(782)
764 Gspart(34, "PHI(1020)", 3, 1.019, 0., 1.486e-22);
765 fPDGCode[fNPDGCodes++]=333; // 34 = PHI (1020)
767 Gspart(35, "D +", 4, 1.87, 1., 1.066e-12);
768 fPDGCode[fNPDGCodes++]=411; // 35 = D+
770 Gspart(36, "D -", 4, 1.87, -1., 1.066e-12);
771 fPDGCode[fNPDGCodes++]=-411; // 36 = D-
773 Gspart(37, "D 0", 3, 1.865, 0., 4.2e-13);
774 fPDGCode[fNPDGCodes++]=421; // 37 = D0
776 Gspart(38, "ANTI D 0", 3, 1.865, 0., 4.2e-13);
777 fPDGCode[fNPDGCodes++]=-421; // 38 = D0 bar
779 fPDGCode[fNPDGCodes++]=-99; // 39 = unassigned
781 fPDGCode[fNPDGCodes++]=-99; // 40 = unassigned
783 fPDGCode[fNPDGCodes++]=-99; // 41 = unassigned
785 Gspart(42, "RHO +", 4, 0.768, 1., 4.353e-24);
786 fPDGCode[fNPDGCodes++]=213; // 42 = RHO+
788 Gspart(43, "RHO -", 4, 0.768, -1., 4.353e-24);
789 fPDGCode[fNPDGCodes++]=-213; // 40 = RHO-
791 Gspart(44, "RHO 0", 3, 0.768, 0., 4.353e-24);
792 fPDGCode[fNPDGCodes++]=113; // 37 = D0
795 // Use ENDF-6 mapping for ions = 10000*z+10*a+iso
797 // and numbers above 5 000 000 for special applications
800 const Int_t kion=10000000;
802 const Int_t kspe=50000000;
804 TDatabasePDG *pdgDB = TDatabasePDG::Instance();
806 const Double_t autogev=0.9314943228;
807 const Double_t hslash = 1.0545726663e-27;
808 const Double_t erggev = 1/1.6021773349e-3;
809 const Double_t hshgev = hslash*erggev;
810 const Double_t yearstosec = 3600*24*365.25;
813 pdgDB->AddParticle("Deuteron","Deuteron",2*autogev+8.071e-3,kTRUE,
814 0,1,"Ion",kion+10020);
815 fPDGCode[fNPDGCodes++]=kion+10020; // 45 = Deuteron
817 pdgDB->AddParticle("Triton","Triton",3*autogev+14.931e-3,kFALSE,
818 hshgev/(12.33*yearstosec),1,"Ion",kion+10030);
819 fPDGCode[fNPDGCodes++]=kion+10030; // 46 = Triton
821 pdgDB->AddParticle("Alpha","Alpha",4*autogev+2.424e-3,kTRUE,
822 hshgev/(12.33*yearstosec),2,"Ion",kion+20040);
823 fPDGCode[fNPDGCodes++]=kion+20040; // 47 = Alpha
825 fPDGCode[fNPDGCodes++]=0; // 48 = geantino mapped to rootino
827 pdgDB->AddParticle("HE3","HE3",3*autogev+14.931e-3,kFALSE,
828 0,2,"Ion",kion+20030);
829 fPDGCode[fNPDGCodes++]=kion+20030; // 49 = HE3
831 pdgDB->AddParticle("Cherenkov","Cherenkov",0,kFALSE,
832 0,0,"Special",kspe+50);
833 fPDGCode[fNPDGCodes++]=kspe+50; // 50 = Cherenkov
835 /* --- Define additional decay modes --- */
836 /* --- omega(783) --- */
837 for (kz = 0; kz < 6; ++kz) {
848 Gsdk(ipa, bratio, mode);
849 /* --- phi(1020) --- */
850 for (kz = 0; kz < 6; ++kz) {
865 Gsdk(ipa, bratio, mode);
867 for (kz = 0; kz < 6; ++kz) {
880 Gsdk(ipa, bratio, mode);
882 for (kz = 0; kz < 6; ++kz) {
895 Gsdk(ipa, bratio, mode);
897 for (kz = 0; kz < 6; ++kz) {
908 Gsdk(ipa, bratio, mode);
909 /* --- Anti D0 --- */
910 for (kz = 0; kz < 6; ++kz) {
921 Gsdk(ipa, bratio, mode);
923 for (kz = 0; kz < 6; ++kz) {
930 Gsdk(ipa, bratio, mode);
932 for (kz = 0; kz < 6; ++kz) {
939 Gsdk(ipa, bratio, mode);
941 for (kz = 0; kz < 6; ++kz) {
948 Gsdk(ipa, bratio, mode);
951 for (kz = 0; kz < 6; ++kz) {
960 Gsdk(ipa, bratio, mode);
963 Gsdk(ipa, bratio, mode);
966 Gsdk(ipa, bratio, mode);
971 //_____________________________________________________________________________
972 Int_t TGeant3::VolId(Text_t *name) const
975 // Return the unique numeric identifier for volume name
978 strncpy((char *) &gname, name, 4);
979 for(i=1; i<=fGcnum->nvolum; i++)
980 if(gname == fZiq[fGclink->jvolum+i]) return i;
981 printf("VolId: Volume %s not found\n",name);
985 //_____________________________________________________________________________
986 Int_t TGeant3::NofVolumes() const
989 // Return total number of volumes in the geometry
991 return fGcnum->nvolum;
994 //_____________________________________________________________________________
995 const char* TGeant3::VolName(Int_t id) const
998 // Return the volume name given the volume identifier
1000 static char name[5];
1001 if(id<1 || id > fGcnum->nvolum || fGclink->jvolum<=0)
1002 strcpy(name,"NULL");
1004 strncpy(name,(char *)&fZiq[fGclink->jvolum+id],4);
1009 //_____________________________________________________________________________
1010 Float_t TGeant3::Xsec(char* reac, Float_t energy, Int_t part, Int_t mate)
1012 Int_t gpart = IdFromPDG(part);
1013 if(!strcmp(reac,"PHOT"))
1016 Error("Xsec","Can calculate photoelectric only for photons\n");
1022 //_____________________________________________________________________________
1023 void TGeant3::TrackPosition(TLorentzVector &xyz) const
1026 // Return the current position in the master reference frame of the
1027 // track being transported
1029 xyz[0]=fGctrak->vect[0];
1030 xyz[1]=fGctrak->vect[1];
1031 xyz[2]=fGctrak->vect[2];
1032 xyz[3]=fGctrak->tofg;
1035 //_____________________________________________________________________________
1036 Float_t TGeant3::TrackTime() const
1039 // Return the current time of flight of the track being transported
1041 return fGctrak->tofg;
1044 //_____________________________________________________________________________
1045 void TGeant3::TrackMomentum(TLorentzVector &xyz) const
1048 // Return the direction and the momentum (GeV/c) of the track
1049 // currently being transported
1051 Double_t ptot=fGctrak->vect[6];
1052 xyz[0]=fGctrak->vect[3]*ptot;
1053 xyz[1]=fGctrak->vect[4]*ptot;
1054 xyz[2]=fGctrak->vect[5]*ptot;
1055 xyz[3]=fGctrak->getot;
1058 //_____________________________________________________________________________
1059 Float_t TGeant3::TrackCharge() const
1062 // Return charge of the track currently transported
1064 return fGckine->charge;
1067 //_____________________________________________________________________________
1068 Float_t TGeant3::TrackMass() const
1071 // Return the mass of the track currently transported
1073 return fGckine->amass;
1076 //_____________________________________________________________________________
1077 Int_t TGeant3::TrackPid() const
1080 // Return the id of the particle transported
1082 return PDGFromId(fGckine->ipart);
1085 //_____________________________________________________________________________
1086 Float_t TGeant3::TrackStep() const
1089 // Return the length in centimeters of the current step
1091 return fGctrak->step;
1094 //_____________________________________________________________________________
1095 Float_t TGeant3::TrackLength() const
1098 // Return the length of the current track from its origin
1100 return fGctrak->sleng;
1103 //_____________________________________________________________________________
1104 Bool_t TGeant3::IsTrackInside() const
1107 // True if the track is not at the boundary of the current volume
1109 return (fGctrak->inwvol==0);
1112 //_____________________________________________________________________________
1113 Bool_t TGeant3::IsTrackEntering() const
1116 // True if this is the first step of the track in the current volume
1118 return (fGctrak->inwvol==1);
1121 //_____________________________________________________________________________
1122 Bool_t TGeant3::IsTrackExiting() const
1125 // True if this is the last step of the track in the current volume
1127 return (fGctrak->inwvol==2);
1130 //_____________________________________________________________________________
1131 Bool_t TGeant3::IsTrackOut() const
1134 // True if the track is out of the setup
1136 return (fGctrak->inwvol==3);
1139 //_____________________________________________________________________________
1140 Bool_t TGeant3::IsTrackStop() const
1143 // True if the track energy has fallen below the threshold
1145 return (fGctrak->istop==2);
1148 //_____________________________________________________________________________
1149 Int_t TGeant3::NSecondaries() const
1152 // Number of secondary particles generated in the current step
1154 return fGcking->ngkine;
1157 //_____________________________________________________________________________
1158 Int_t TGeant3::CurrentEvent() const
1161 // Number of the current event
1163 return fGcflag->idevt;
1166 //_____________________________________________________________________________
1167 const char* TGeant3::ProdProcess() const
1170 // Name of the process that has produced the secondary particles
1171 // in the current step
1173 static char proc[5];
1174 const Int_t ipmec[13] = { 5,6,7,8,9,10,11,12,21,23,25,105,108 };
1177 if(fGcking->ngkine>0) {
1178 for (km = 0; km < fGctrak->nmec; ++km) {
1179 for (im = 0; im < 13; ++im) {
1180 if (fGctrak->lmec[km] == ipmec[im]) {
1181 mec = fGctrak->lmec[km];
1182 if (0 < mec && mec < 31) {
1183 strncpy(proc,(char *)&fGctrak->namec[mec - 1],4);
1184 } else if (mec - 100 <= 30 && mec - 100 > 0) {
1185 strncpy(proc,(char *)&fGctpol->namec1[mec - 101],4);
1192 strcpy(proc,"UNKN");
1193 } else strcpy(proc,"NONE");
1197 //_____________________________________________________________________________
1198 void TGeant3::GetSecondary(Int_t isec, Int_t& ipart,
1199 TLorentzVector &x, TLorentzVector &p)
1202 // Get the parameters of the secondary track number isec produced
1203 // in the current step
1206 if(-1<isec && isec<fGcking->ngkine) {
1207 ipart=Int_t (fGcking->gkin[isec][4] +0.5);
1209 x[i]=fGckin3->gpos[isec][i];
1210 p[i]=fGcking->gkin[isec][i];
1212 x[3]=fGcking->tofd[isec];
1213 p[3]=fGcking->gkin[isec][3];
1215 printf(" * TGeant3::GetSecondary * Secondary %d does not exist\n",isec);
1216 x[0]=x[1]=x[2]=x[3]=p[0]=p[1]=p[2]=p[3]=0;
1221 //_____________________________________________________________________________
1222 void TGeant3::InitLego()
1225 SetDEBU(0,0,0); //do not print a message
1228 //_____________________________________________________________________________
1229 Bool_t TGeant3::IsTrackDisappeared() const
1232 // True if the current particle has disappered
1233 // either because it decayed or because it underwent
1234 // an inelastic collision
1236 return (fGctrak->istop==1);
1239 //_____________________________________________________________________________
1240 Bool_t TGeant3::IsTrackAlive() const
1243 // True if the current particle is alive and will continue to be
1246 return (fGctrak->istop==0);
1249 //_____________________________________________________________________________
1250 void TGeant3::StopTrack()
1253 // Stop the transport of the current particle and skip to the next
1258 //_____________________________________________________________________________
1259 void TGeant3::StopEvent()
1262 // Stop simulation of the current event and skip to the next
1267 //_____________________________________________________________________________
1268 Float_t TGeant3::MaxStep() const
1271 // Return the maximum step length in the current medium
1273 return fGctmed->stemax;
1276 //_____________________________________________________________________________
1277 void TGeant3::SetColors()
1280 // Set the colors for all the volumes
1281 // this is done sequentially for all volumes
1282 // based on the number of their medium
1285 Int_t jvolum=fGclink->jvolum;
1286 //Int_t jtmed=fGclink->jtmed;
1287 //Int_t jmate=fGclink->jmate;
1288 Int_t nvolum=fGcnum->nvolum;
1291 // Now for all the volumes
1292 for(kv=1;kv<=nvolum;kv++) {
1293 // Get the tracking medium
1294 Int_t itm=Int_t (fZq[fZlq[jvolum-kv]+4]);
1296 //Int_t ima=Int_t (fZq[fZlq[jtmed-itm]+6]);
1298 //Float_t z=fZq[fZlq[jmate-ima]+7];
1299 // Find color number
1300 //icol = Int_t(z)%6+2;
1301 //icol = 17+Int_t(z*150./92.);
1304 strncpy(name,(char*)&fZiq[jvolum+kv],4);
1306 Gsatt(name,"COLO",icol);
1310 //_____________________________________________________________________________
1311 void TGeant3::SetMaxStep(Float_t maxstep)
1314 // Set the maximum step allowed till the particle is in the current medium
1316 fGctmed->stemax=maxstep;
1319 //_____________________________________________________________________________
1320 void TGeant3::SetMaxNStep(Int_t maxnstp)
1323 // Set the maximum number of steps till the particle is in the current medium
1325 fGctrak->maxnst=maxnstp;
1328 //_____________________________________________________________________________
1329 Int_t TGeant3::GetMaxNStep() const
1332 // Maximum number of steps allowed in current medium
1334 return fGctrak->maxnst;
1337 //_____________________________________________________________________________
1338 void TGeant3::Material(Int_t& kmat, const char* name, Float_t a, Float_t z,
1339 Float_t dens, Float_t radl, Float_t absl, Float_t* buf,
1343 // Defines a Material
1345 // kmat number assigned to the material
1346 // name material name
1347 // a atomic mass in au
1349 // dens density in g/cm3
1350 // absl absorbtion length in cm
1351 // if >=0 it is ignored and the program
1352 // calculates it, if <0. -absl is taken
1353 // radl radiation length in cm
1354 // if >=0 it is ignored and the program
1355 // calculates it, if <0. -radl is taken
1356 // buf pointer to an array of user words
1357 // nbuf number of user words
1359 Int_t jmate=fGclink->jmate;
1365 for(i=1; i<=ns; i++) {
1366 if(fZlq[jmate-i]==0) {
1372 gsmate(kmat,PASSCHARD(name), a, z, dens, radl, absl, buf,
1373 nwbuf PASSCHARL(name));
1376 //_____________________________________________________________________________
1377 void TGeant3::Mixture(Int_t& kmat, const char* name, Float_t* a, Float_t* z,
1378 Float_t dens, Int_t nlmat, Float_t* wmat)
1381 // Defines mixture OR COMPOUND IMAT as composed by
1382 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
1384 // If NLMAT > 0 then wmat contains the proportion by
1385 // weights of each basic material in the mixture.
1387 // If nlmat < 0 then WMAT contains the number of atoms
1388 // of a given kind into the molecule of the COMPOUND
1389 // In this case, WMAT in output is changed to relative
1392 Int_t jmate=fGclink->jmate;
1398 for(i=1; i<=ns; i++) {
1399 if(fZlq[jmate-i]==0) {
1405 gsmixt(kmat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
1408 //_____________________________________________________________________________
1409 void TGeant3::Medium(Int_t& kmed, const char* name, Int_t nmat, Int_t isvol,
1410 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
1411 Float_t stemax, Float_t deemax, Float_t epsil,
1412 Float_t stmin, Float_t* ubuf, Int_t nbuf)
1415 // kmed tracking medium number assigned
1416 // name tracking medium name
1417 // nmat material number
1418 // isvol sensitive volume flag
1419 // ifield magnetic field
1420 // fieldm max. field value (kilogauss)
1421 // tmaxfd max. angle due to field (deg/step)
1422 // stemax max. step allowed
1423 // deemax max. fraction of energy lost in a step
1424 // epsil tracking precision (cm)
1425 // stmin min. step due to continuos processes (cm)
1427 // ifield = 0 if no magnetic field; ifield = -1 if user decision in guswim;
1428 // ifield = 1 if tracking performed with grkuta; ifield = 2 if tracking
1429 // performed with ghelix; ifield = 3 if tracking performed with ghelx3.
1431 Int_t jtmed=fGclink->jtmed;
1437 for(i=1; i<=ns; i++) {
1438 if(fZlq[jtmed-i]==0) {
1444 gstmed(kmed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1445 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1448 //_____________________________________________________________________________
1449 void TGeant3::Matrix(Int_t& krot, Float_t thex, Float_t phix, Float_t they,
1450 Float_t phiy, Float_t thez, Float_t phiz)
1453 // krot rotation matrix number assigned
1454 // theta1 polar angle for axis i
1455 // phi1 azimuthal angle for axis i
1456 // theta2 polar angle for axis ii
1457 // phi2 azimuthal angle for axis ii
1458 // theta3 polar angle for axis iii
1459 // phi3 azimuthal angle for axis iii
1461 // it defines the rotation matrix number irot.
1463 Int_t jrotm=fGclink->jrotm;
1469 for(i=1; i<=ns; i++) {
1470 if(fZlq[jrotm-i]==0) {
1476 gsrotm(krot, thex, phix, they, phiy, thez, phiz);
1479 //_____________________________________________________________________________
1480 Int_t TGeant3::GetMedium() const
1483 // Return the number of the current medium
1485 return fGctmed->numed;
1488 //_____________________________________________________________________________
1489 Float_t TGeant3::Edep() const
1492 // Return the energy lost in the current step
1494 return fGctrak->destep;
1497 //_____________________________________________________________________________
1498 Float_t TGeant3::Etot() const
1501 // Return the total energy of the current track
1503 return fGctrak->getot;
1506 //_____________________________________________________________________________
1507 void TGeant3::Rndm(Float_t* r, const Int_t n) const
1510 // Return an array of n random numbers uniformly distributed
1511 // between 0 and 1 not included
1516 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1518 // Functions from GBASE
1520 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1522 //____________________________________________________________________________
1523 void TGeant3::Gfile(const char *filename, const char *option)
1526 // Routine to open a GEANT/RZ data base.
1528 // LUN logical unit number associated to the file
1530 // CHFILE RZ file name
1532 // CHOPT is a character string which may be
1533 // N To create a new file
1534 // U to open an existing file for update
1535 // " " to open an existing file for read only
1536 // Q The initial allocation (default 1000 records)
1537 // is given in IQUEST(10)
1538 // X Open the file in exchange format
1539 // I Read all data structures from file to memory
1540 // O Write all data structures from memory to file
1543 // If options "I" or "O" all data structures are read or
1544 // written from/to file and the file is closed.
1545 // See routine GRMDIR to create subdirectories
1546 // See routines GROUT,GRIN to write,read objects
1548 grfile(21, PASSCHARD(filename), PASSCHARD(option) PASSCHARL(filename)
1552 //____________________________________________________________________________
1553 void TGeant3::Gpcxyz()
1556 // Print track and volume parameters at current point
1561 //_____________________________________________________________________________
1562 void TGeant3::Ggclos()
1565 // Closes off the geometry setting.
1566 // Initializes the search list for the contents of each
1567 // volume following the order they have been positioned, and
1568 // inserting the content '0' when a call to GSNEXT (-1) has
1569 // been required by the user.
1570 // Performs the development of the JVOLUM structure for all
1571 // volumes with variable parameters, by calling GGDVLP.
1572 // Interprets the user calls to GSORD, through GGORD.
1573 // Computes and stores in a bank (next to JVOLUM mother bank)
1574 // the number of levels in the geometrical tree and the
1575 // maximum number of contents per level, by calling GGNLEV.
1576 // Sets status bit for CONCAVE volumes, through GGCAVE.
1577 // Completes the JSET structure with the list of volume names
1578 // which identify uniquely a given physical detector, the
1579 // list of bit numbers to pack the corresponding volume copy
1580 // numbers, and the generic path(s) in the JVOLUM tree,
1581 // through the routine GHCLOS.
1586 //_____________________________________________________________________________
1587 void TGeant3::Glast()
1590 // Finish a Geant run
1595 //_____________________________________________________________________________
1596 void TGeant3::Gprint(const char *name)
1599 // Routine to print data structures
1600 // CHNAME name of a data structure
1604 gprint(PASSCHARD(vname),0 PASSCHARL(vname));
1607 //_____________________________________________________________________________
1608 void TGeant3::Grun()
1611 // Steering function to process one run
1616 //_____________________________________________________________________________
1617 void TGeant3::Gtrig()
1620 // Steering function to process one event
1625 //_____________________________________________________________________________
1626 void TGeant3::Gtrigc()
1629 // Clear event partition
1634 //_____________________________________________________________________________
1635 void TGeant3::Gtrigi()
1638 // Initialises event partition
1643 //_____________________________________________________________________________
1644 void TGeant3::Gwork(Int_t nwork)
1647 // Allocates workspace in ZEBRA memory
1652 //_____________________________________________________________________________
1653 void TGeant3::Gzinit()
1656 // To initialise GEANT/ZEBRA data structures
1661 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1663 // Functions from GCONS
1665 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1667 //_____________________________________________________________________________
1668 void TGeant3::Gfmate(Int_t imat, char *name, Float_t &a, Float_t &z,
1669 Float_t &dens, Float_t &radl, Float_t &absl,
1670 Float_t* ubuf, Int_t& nbuf)
1673 // Return parameters for material IMAT
1675 gfmate(imat, PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
1679 //_____________________________________________________________________________
1680 void TGeant3::Gfpart(Int_t ipart, char *name, Int_t &itrtyp,
1681 Float_t &amass, Float_t &charge, Float_t &tlife)
1684 // Return parameters for particle of type IPART
1688 Int_t igpart = IdFromPDG(ipart);
1689 gfpart(igpart, PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
1693 //_____________________________________________________________________________
1694 void TGeant3::Gftmed(Int_t numed, char *name, Int_t &nmat, Int_t &isvol,
1695 Int_t &ifield, Float_t &fieldm, Float_t &tmaxfd,
1696 Float_t &stemax, Float_t &deemax, Float_t &epsil,
1697 Float_t &stmin, Float_t *ubuf, Int_t *nbuf)
1700 // Return parameters for tracking medium NUMED
1702 gftmed(numed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1703 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1706 //_____________________________________________________________________________
1707 void TGeant3::Gmate()
1710 // Define standard GEANT materials
1715 //_____________________________________________________________________________
1716 void TGeant3::Gpart()
1719 // Define standard GEANT particles plus selected decay modes
1720 // and branching ratios.
1725 //_____________________________________________________________________________
1726 void TGeant3::Gsdk(Int_t ipart, Float_t *bratio, Int_t *mode)
1728 // Defines branching ratios and decay modes for standard
1730 gsdk(ipart,bratio,mode);
1733 //_____________________________________________________________________________
1734 void TGeant3::Gsmate(Int_t imat, const char *name, Float_t a, Float_t z,
1735 Float_t dens, Float_t radl, Float_t absl)
1738 // Defines a Material
1740 // kmat number assigned to the material
1741 // name material name
1742 // a atomic mass in au
1744 // dens density in g/cm3
1745 // absl absorbtion length in cm
1746 // if >=0 it is ignored and the program
1747 // calculates it, if <0. -absl is taken
1748 // radl radiation length in cm
1749 // if >=0 it is ignored and the program
1750 // calculates it, if <0. -radl is taken
1751 // buf pointer to an array of user words
1752 // nbuf number of user words
1756 gsmate(imat,PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
1760 //_____________________________________________________________________________
1761 void TGeant3::Gsmixt(Int_t imat, const char *name, Float_t *a, Float_t *z,
1762 Float_t dens, Int_t nlmat, Float_t *wmat)
1765 // Defines mixture OR COMPOUND IMAT as composed by
1766 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
1768 // If NLMAT.GT.0 then WMAT contains the PROPORTION BY
1769 // WEIGTHS OF EACH BASIC MATERIAL IN THE MIXTURE.
1771 // If NLMAT.LT.0 then WMAT contains the number of atoms
1772 // of a given kind into the molecule of the COMPOUND
1773 // In this case, WMAT in output is changed to relative
1776 gsmixt(imat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
1779 //_____________________________________________________________________________
1780 void TGeant3::Gspart(Int_t ipart, const char *name, Int_t itrtyp,
1781 Float_t amass, Float_t charge, Float_t tlife)
1784 // Store particle parameters
1786 // ipart particle code
1787 // name particle name
1788 // itrtyp transport method (see GEANT manual)
1789 // amass mass in GeV/c2
1790 // charge charge in electron units
1791 // tlife lifetime in seconds
1795 gspart(ipart,PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
1799 //_____________________________________________________________________________
1800 void TGeant3::Gstmed(Int_t numed, const char *name, Int_t nmat, Int_t isvol,
1801 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
1802 Float_t stemax, Float_t deemax, Float_t epsil,
1806 // NTMED Tracking medium number
1807 // NAME Tracking medium name
1808 // NMAT Material number
1809 // ISVOL Sensitive volume flag
1810 // IFIELD Magnetic field
1811 // FIELDM Max. field value (Kilogauss)
1812 // TMAXFD Max. angle due to field (deg/step)
1813 // STEMAX Max. step allowed
1814 // DEEMAX Max. fraction of energy lost in a step
1815 // EPSIL Tracking precision (cm)
1816 // STMIN Min. step due to continuos processes (cm)
1818 // IFIELD = 0 if no magnetic field; IFIELD = -1 if user decision in GUSWIM;
1819 // IFIELD = 1 if tracking performed with GRKUTA; IFIELD = 2 if tracking
1820 // performed with GHELIX; IFIELD = 3 if tracking performed with GHELX3.
1824 gstmed(numed,PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1825 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1828 //_____________________________________________________________________________
1829 void TGeant3::Gsckov(Int_t itmed, Int_t npckov, Float_t *ppckov,
1830 Float_t *absco, Float_t *effic, Float_t *rindex)
1833 // Stores the tables for UV photon tracking in medium ITMED
1834 // Please note that it is the user's responsability to
1835 // provide all the coefficients:
1838 // ITMED Tracking medium number
1839 // NPCKOV Number of bins of each table
1840 // PPCKOV Value of photon momentum (in GeV)
1841 // ABSCO Absorbtion coefficients
1842 // dielectric: absorbtion length in cm
1843 // metals : absorbtion fraction (0<=x<=1)
1844 // EFFIC Detection efficiency for UV photons
1845 // RINDEX Refraction index (if=0 metal)
1847 gsckov(itmed,npckov,ppckov,absco,effic,rindex);
1850 //_____________________________________________________________________________
1851 void TGeant3::Gstpar(Int_t itmed, const char *param, Float_t parval)
1854 // To change the value of cut or mechanism "CHPAR"
1855 // to a new value PARVAL for tracking medium ITMED
1856 // The data structure JTMED contains the standard tracking
1857 // parameters (CUTS and flags to control the physics processes) which
1858 // are used by default for all tracking media. It is possible to
1859 // redefine individually with GSTPAR any of these parameters for a
1860 // given tracking medium.
1861 // ITMED tracking medium number
1862 // CHPAR is a character string (variable name)
1863 // PARVAL must be given as a floating point.
1865 gstpar(itmed,PASSCHARD(param), parval PASSCHARL(param));
1868 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1870 // Functions from GCONS
1872 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1874 //_____________________________________________________________________________
1875 void TGeant3::Gfkine(Int_t itra, Float_t *vert, Float_t *pvert, Int_t &ipart,
1878 // Storing/Retrieving Vertex and Track parameters
1879 // ----------------------------------------------
1881 // Stores vertex parameters.
1882 // VERT array of (x,y,z) position of the vertex
1883 // NTBEAM beam track number origin of the vertex
1884 // =0 if none exists
1885 // NTTARG target track number origin of the vertex
1886 // UBUF user array of NUBUF floating point numbers
1888 // NVTX new vertex number (=0 in case of error).
1889 // Prints vertex parameters.
1890 // IVTX for vertex IVTX.
1891 // (For all vertices if IVTX=0)
1892 // Stores long life track parameters.
1893 // PLAB components of momentum
1894 // IPART type of particle (see GSPART)
1895 // NV vertex number origin of track
1896 // UBUF array of NUBUF floating point user parameters
1898 // NT track number (if=0 error).
1899 // Retrieves long life track parameters.
1900 // ITRA track number for which parameters are requested
1901 // VERT vector origin of the track
1902 // PVERT 4 momentum components at the track origin
1903 // IPART particle type (=0 if track ITRA does not exist)
1904 // NVERT vertex number origin of the track
1905 // UBUF user words stored in GSKINE.
1906 // Prints initial track parameters.
1907 // ITRA for track ITRA
1908 // (For all tracks if ITRA=0)
1912 gfkine(itra,vert,pvert,ipart,nvert,ubuf,nbuf);
1915 //_____________________________________________________________________________
1916 void TGeant3::Gfvert(Int_t nvtx, Float_t *v, Int_t &ntbeam, Int_t &nttarg,
1920 // Retrieves the parameter of a vertex bank
1921 // Vertex is generated from tracks NTBEAM NTTARG
1922 // NVTX is the new vertex number
1926 gfvert(nvtx,v,ntbeam,nttarg,tofg,ubuf,nbuf);
1929 //_____________________________________________________________________________
1930 Int_t TGeant3::Gskine(Float_t *plab, Int_t ipart, Int_t nv, Float_t *buf,
1934 // Store kinematics of track NT into data structure
1935 // Track is coming from vertex NV
1938 gskine(plab, ipart, nv, buf, nwbuf, nt);
1942 //_____________________________________________________________________________
1943 Int_t TGeant3::Gsvert(Float_t *v, Int_t ntbeam, Int_t nttarg, Float_t *ubuf,
1947 // Creates a new vertex bank
1948 // Vertex is generated from tracks NTBEAM NTTARG
1949 // NVTX is the new vertex number
1952 gsvert(v, ntbeam, nttarg, ubuf, nwbuf, nwtx);
1956 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1958 // Functions from GPHYS
1960 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1962 //_____________________________________________________________________________
1963 void TGeant3::Gphysi()
1966 // Initialise material constants for all the physics
1967 // mechanisms used by GEANT
1972 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1974 // Functions from GTRAK
1976 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1978 //_____________________________________________________________________________
1979 void TGeant3::Gdebug()
1982 // Debug the current step
1987 //_____________________________________________________________________________
1988 void TGeant3::Gekbin()
1991 // To find bin number in kinetic energy table
1992 // stored in ELOW(NEKBIN)
1997 //_____________________________________________________________________________
1998 void TGeant3::Gfinds()
2001 // Returns the set/volume parameters corresponding to
2002 // the current space point in /GCTRAK/
2003 // and fill common /GCSETS/
2005 // IHSET user set identifier
2006 // IHDET user detector identifier
2007 // ISET set number in JSET
2008 // IDET detector number in JS=LQ(JSET-ISET)
2009 // IDTYPE detector type (1,2)
2010 // NUMBV detector volume numbers (array of length NVNAME)
2011 // NVNAME number of volume levels
2016 //_____________________________________________________________________________
2017 void TGeant3::Gsking(Int_t igk)
2020 // Stores in stack JSTAK either the IGKth track of /GCKING/,
2021 // or the NGKINE tracks when IGK is 0.
2026 //_____________________________________________________________________________
2027 void TGeant3::Gskpho(Int_t igk)
2030 // Stores in stack JSTAK either the IGKth Cherenkov photon of
2031 // /GCKIN2/, or the NPHOT tracks when IGK is 0.
2036 //_____________________________________________________________________________
2037 void TGeant3::Gsstak(Int_t iflag)
2040 // Stores in auxiliary stack JSTAK the particle currently
2041 // described in common /GCKINE/.
2043 // On request, creates also an entry in structure JKINE :
2045 // 0 : No entry in JKINE structure required (user)
2046 // 1 : New entry in JVERTX / JKINE structures required (user)
2047 // <0 : New entry in JKINE structure at vertex -IFLAG (user)
2048 // 2 : Entry in JKINE structure exists already (from GTREVE)
2053 //_____________________________________________________________________________
2054 void TGeant3::Gsxyz()
2057 // Store space point VECT in banks JXYZ
2062 //_____________________________________________________________________________
2063 void TGeant3::Gtrack()
2066 // Controls tracking of current particle
2071 //_____________________________________________________________________________
2072 void TGeant3::Gtreve()
2075 // Controls tracking of all particles belonging to the current event
2080 //_____________________________________________________________________________
2081 void TGeant3::Gtreve_root()
2084 // Controls tracking of all particles belonging to the current event
2089 //_____________________________________________________________________________
2090 void TGeant3::Grndm(Float_t *rvec, const Int_t len) const
2093 // To generate a vector RVECV of LEN random numbers
2094 // Copy of the CERN Library routine RANECU
2098 //_____________________________________________________________________________
2099 void TGeant3::Grndmq(Int_t &is1, Int_t &is2, const Int_t iseq,
2100 const Text_t *chopt)
2103 // To set/retrieve the seed of the random number generator
2105 grndmq(is1,is2,iseq,PASSCHARD(chopt) PASSCHARL(chopt));
2108 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2110 // Functions from GDRAW
2112 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2114 //_____________________________________________________________________________
2115 void TGeant3::Gdxyz(Int_t it)
2118 // Draw the points stored with Gsxyz relative to track it
2123 //_____________________________________________________________________________
2124 void TGeant3::Gdcxyz()
2127 // Draw the position of the current track
2132 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2134 // Functions from GGEOM
2136 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2138 //_____________________________________________________________________________
2139 void TGeant3::Gdtom(Float_t *xd, Float_t *xm, Int_t iflag)
2142 // Computes coordinates XM (Master Reference System
2143 // knowing the coordinates XD (Detector Ref System)
2144 // The local reference system can be initialized by
2145 // - the tracking routines and GDTOM used in GUSTEP
2146 // - a call to GSCMED(NLEVEL,NAMES,NUMBER)
2147 // (inverse routine is GMTOD)
2149 // If IFLAG=1 convert coordinates
2150 // IFLAG=2 convert direction cosinus
2152 gdtom(xd, xm, iflag);
2155 //_____________________________________________________________________________
2156 void TGeant3::Glmoth(const char* iudet, Int_t iunum, Int_t &nlev, Int_t *lvols,
2160 // Loads the top part of the Volume tree in LVOLS (IVO's),
2161 // LINDX (IN indices) for a given volume defined through
2162 // its name IUDET and number IUNUM.
2164 // The routine stores only upto the last level where JVOLUM
2165 // data structure is developed. If there is no development
2166 // above the current level, it returns NLEV zero.
2168 glmoth(PASSCHARD(iudet), iunum, nlev, lvols, lindx, idum PASSCHARL(iudet));
2171 //_____________________________________________________________________________
2172 void TGeant3::Gmedia(Float_t *x, Int_t &numed)
2175 // Finds in which volume/medium the point X is, and updates the
2176 // common /GCVOLU/ and the structure JGPAR accordingly.
2178 // NUMED returns the tracking medium number, or 0 if point is
2179 // outside the experimental setup.
2184 //_____________________________________________________________________________
2185 void TGeant3::Gmtod(Float_t *xm, Float_t *xd, Int_t iflag)
2188 // Computes coordinates XD (in DRS)
2189 // from known coordinates XM in MRS
2190 // The local reference system can be initialized by
2191 // - the tracking routines and GMTOD used in GUSTEP
2192 // - a call to GMEDIA(XM,NUMED)
2193 // - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER)
2194 // (inverse routine is GDTOM)
2196 // If IFLAG=1 convert coordinates
2197 // IFLAG=2 convert direction cosinus
2199 gmtod(xm, xd, iflag);
2202 //_____________________________________________________________________________
2203 void TGeant3::Gsdvn(const char *name, const char *mother, Int_t ndiv,
2207 // Create a new volume by dividing an existing one
2210 // MOTHER Mother volume name
2211 // NDIV Number of divisions
2214 // X,Y,Z of CAXIS will be translated to 1,2,3 for IAXIS.
2215 // It divides a previously defined volume.
2220 Vname(mother,vmother);
2221 gsdvn(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis PASSCHARL(vname)
2222 PASSCHARL(vmother));
2225 //_____________________________________________________________________________
2226 void TGeant3::Gsdvn2(const char *name, const char *mother, Int_t ndiv,
2227 Int_t iaxis, Float_t c0i, Int_t numed)
2230 // Create a new volume by dividing an existing one
2232 // Divides mother into ndiv divisions called name
2233 // along axis iaxis starting at coordinate value c0.
2234 // the new volume created will be medium number numed.
2239 Vname(mother,vmother);
2240 gsdvn2(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis, c0i, numed
2241 PASSCHARL(vname) PASSCHARL(vmother));
2244 //_____________________________________________________________________________
2245 void TGeant3::Gsdvs(const char *name, const char *mother, Float_t step,
2246 Int_t iaxis, Int_t numed)
2249 // Create a new volume by dividing an existing one
2254 Vname(mother,vmother);
2255 gsdvs(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed
2256 PASSCHARL(vname) PASSCHARL(vmother));
2259 //_____________________________________________________________________________
2260 void TGeant3::Gsdvs2(const char *name, const char *mother, Float_t step,
2261 Int_t iaxis, Float_t c0, Int_t numed)
2264 // Create a new volume by dividing an existing one
2269 Vname(mother,vmother);
2270 gsdvs2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0, numed
2271 PASSCHARL(vname) PASSCHARL(vmother));
2274 //_____________________________________________________________________________
2275 void TGeant3::Gsdvt(const char *name, const char *mother, Float_t step,
2276 Int_t iaxis, Int_t numed, Int_t ndvmx)
2279 // Create a new volume by dividing an existing one
2281 // Divides MOTHER into divisions called NAME along
2282 // axis IAXIS in steps of STEP. If not exactly divisible
2283 // will make as many as possible and will centre them
2284 // with respect to the mother. Divisions will have medium
2285 // number NUMED. If NUMED is 0, NUMED of MOTHER is taken.
2286 // NDVMX is the expected maximum number of divisions
2287 // (If 0, no protection tests are performed)
2292 Vname(mother,vmother);
2293 gsdvt(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed, ndvmx
2294 PASSCHARL(vname) PASSCHARL(vmother));
2297 //_____________________________________________________________________________
2298 void TGeant3::Gsdvt2(const char *name, const char *mother, Float_t step,
2299 Int_t iaxis, Float_t c0, Int_t numed, Int_t ndvmx)
2302 // Create a new volume by dividing an existing one
2304 // Divides MOTHER into divisions called NAME along
2305 // axis IAXIS starting at coordinate value C0 with step
2307 // The new volume created will have medium number NUMED.
2308 // If NUMED is 0, NUMED of mother is taken.
2309 // NDVMX is the expected maximum number of divisions
2310 // (If 0, no protection tests are performed)
2315 Vname(mother,vmother);
2316 gsdvt2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0,
2317 numed, ndvmx PASSCHARL(vname) PASSCHARL(vmother));
2320 //_____________________________________________________________________________
2321 void TGeant3::Gsord(const char *name, Int_t iax)
2324 // Flags volume CHNAME whose contents will have to be ordered
2325 // along axis IAX, by setting the search flag to -IAX
2329 // IAX = 4 Rxy (static ordering only -> GTMEDI)
2330 // IAX = 14 Rxy (also dynamic ordering -> GTNEXT)
2331 // IAX = 5 Rxyz (static ordering only -> GTMEDI)
2332 // IAX = 15 Rxyz (also dynamic ordering -> GTNEXT)
2333 // IAX = 6 PHI (PHI=0 => X axis)
2334 // IAX = 7 THETA (THETA=0 => Z axis)
2338 gsord(PASSCHARD(vname), iax PASSCHARL(vname));
2341 //_____________________________________________________________________________
2342 void TGeant3::Gspos(const char *name, Int_t nr, const char *mother, Float_t x,
2343 Float_t y, Float_t z, Int_t irot, const char *konly)
2346 // Position a volume into an existing one
2349 // NUMBER Copy number of the volume
2350 // MOTHER Mother volume name
2351 // X X coord. of the volume in mother ref. sys.
2352 // Y Y coord. of the volume in mother ref. sys.
2353 // Z Z coord. of the volume in mother ref. sys.
2354 // IROT Rotation matrix number w.r.t. mother ref. sys.
2355 // ONLY ONLY/MANY flag
2357 // It positions a previously defined volume in the mother.
2362 Vname(mother,vmother);
2363 gspos(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2364 PASSCHARD(konly) PASSCHARL(vname) PASSCHARL(vmother)
2368 //_____________________________________________________________________________
2369 void TGeant3::Gsposp(const char *name, Int_t nr, const char *mother,
2370 Float_t x, Float_t y, Float_t z, Int_t irot,
2371 const char *konly, Float_t *upar, Int_t np )
2374 // Place a copy of generic volume NAME with user number
2375 // NR inside MOTHER, with its parameters UPAR(1..NP)
2380 Vname(mother,vmother);
2381 gsposp(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2382 PASSCHARD(konly), upar, np PASSCHARL(vname) PASSCHARL(vmother)
2386 //_____________________________________________________________________________
2387 void TGeant3::Gsrotm(Int_t nmat, Float_t theta1, Float_t phi1, Float_t theta2,
2388 Float_t phi2, Float_t theta3, Float_t phi3)
2391 // nmat Rotation matrix number
2392 // THETA1 Polar angle for axis I
2393 // PHI1 Azimuthal angle for axis I
2394 // THETA2 Polar angle for axis II
2395 // PHI2 Azimuthal angle for axis II
2396 // THETA3 Polar angle for axis III
2397 // PHI3 Azimuthal angle for axis III
2399 // It defines the rotation matrix number IROT.
2401 gsrotm(nmat, theta1, phi1, theta2, phi2, theta3, phi3);
2404 //_____________________________________________________________________________
2405 void TGeant3::Gprotm(Int_t nmat)
2408 // To print rotation matrices structure JROTM
2409 // nmat Rotation matrix number
2414 //_____________________________________________________________________________
2415 Int_t TGeant3::Gsvolu(const char *name, const char *shape, Int_t nmed,
2416 Float_t *upar, Int_t npar)
2420 // SHAPE Volume type
2421 // NUMED Tracking medium number
2422 // NPAR Number of shape parameters
2423 // UPAR Vector containing shape parameters
2425 // It creates a new volume in the JVOLUM data structure.
2431 Vname(shape,vshape);
2432 gsvolu(PASSCHARD(vname), PASSCHARD(vshape), nmed, upar, npar, ivolu
2433 PASSCHARL(vname) PASSCHARL(vshape));
2437 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2439 // T H E D R A W I N G P A C K A G E
2440 // ======================================
2441 // Drawing functions. These functions allow the visualization in several ways
2442 // of the volumes defined in the geometrical data structure. It is possible
2443 // to draw the logical tree of volumes belonging to the detector (DTREE),
2444 // to show their geometrical specification (DSPEC,DFSPC), to draw them
2445 // and their cut views (DRAW, DCUT). Moreover, it is possible to execute
2446 // these commands when the hidden line removal option is activated; in
2447 // this case, the volumes can be also either translated in the space
2448 // (SHIFT), or clipped by boolean operation (CVOL). In addition, it is
2449 // possible to fill the surfaces of the volumes
2450 // with solid colours when the shading option (SHAD) is activated.
2451 // Several tools (ZOOM, LENS) have been developed to zoom detailed parts
2452 // of the detectors or to scan physical events as well.
2453 // Finally, the command MOVE will allow the rotation, translation and zooming
2454 // on real time parts of the detectors or tracks and hits of a simulated event.
2455 // Ray-tracing commands. In case the command (DOPT RAYT ON) is executed,
2456 // the drawing is performed by the Geant ray-tracing;
2457 // automatically, the color is assigned according to the tracking medium of each
2458 // volume and the volumes with a density lower/equal than the air are considered
2459 // transparent; if the option (USER) is set (ON) (again via the command (DOPT)),
2460 // the user can set color and visibility for the desired volumes via the command
2461 // (SATT), as usual, relatively to the attributes (COLO) and (SEEN).
2462 // The resolution can be set via the command (SATT * FILL VALUE), where (VALUE)
2463 // is the ratio between the number of pixels drawn and 20 (user coordinates).
2464 // Parallel view and perspective view are possible (DOPT PROJ PARA/PERS); in the
2465 // first case, we assume that the first mother volume of the tree is a box with
2466 // dimensions 10000 X 10000 X 10000 cm and the view point (infinetely far) is
2467 // 5000 cm far from the origin along the Z axis of the user coordinates; in the
2468 // second case, the distance between the observer and the origin of the world
2469 // reference system is set in cm by the command (PERSP NAME VALUE); grand-angle
2470 // or telescopic effects can be achieved changing the scale factors in the command
2471 // (DRAW). When the final picture does not occupy the full window,
2472 // mapping the space before tracing can speed up the drawing, but can also
2473 // produce less precise results; values from 1 to 4 are allowed in the command
2474 // (DOPT MAPP VALUE), the mapping being more precise for increasing (VALUE); for
2475 // (VALUE = 0) no mapping is performed (therefore max precision and lowest speed).
2476 // The command (VALCUT) allows the cutting of the detector by three planes
2477 // ortogonal to the x,y,z axis. The attribute (LSTY) can be set by the command
2478 // SATT for any desired volume and can assume values from 0 to 7; it determines
2479 // the different light processing to be performed for different materials:
2480 // 0 = dark-matt, 1 = bright-matt, 2 = plastic, 3 = ceramic, 4 = rough-metals,
2481 // 5 = shiny-metals, 6 = glass, 7 = mirror. The detector is assumed to be in the
2482 // dark, the ambient light luminosity is 0.2 for each basic hue (the saturation
2483 // is 0.9) and the observer is assumed to have a light source (therefore he will
2484 // produce parallel light in the case of parallel view and point-like-source
2485 // light in the case of perspective view).
2487 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2489 //_____________________________________________________________________________
2490 void TGeant3::Gsatt(const char *name, const char *att, Int_t val)
2494 // IOPT Name of the attribute to be set
2495 // IVAL Value to which the attribute is to be set
2497 // name= "*" stands for all the volumes.
2498 // iopt can be chosen among the following :
2500 // WORK 0=volume name is inactive for the tracking
2501 // 1=volume name is active for the tracking (default)
2503 // SEEN 0=volume name is invisible
2504 // 1=volume name is visible (default)
2505 // -1=volume invisible with all its descendants in the tree
2506 // -2=volume visible but not its descendants in the tree
2508 // LSTY line style 1,2,3,... (default=1)
2509 // LSTY=7 will produce a very precise approximation for
2510 // revolution bodies.
2512 // LWID line width -7,...,1,2,3,..7 (default=1)
2513 // LWID<0 will act as abs(LWID) was set for the volume
2514 // and for all the levels below it. When SHAD is 'ON', LWID
2515 // represent the linewidth of the scan lines filling the surfaces
2516 // (whereas the FILL value represent their number). Therefore
2517 // tuning this parameter will help to obtain the desired
2518 // quality/performance ratio.
2520 // COLO colour code -166,...,1,2,..166 (default=1)
2522 // n=2=red; n=17+m, m=0,25, increasing luminosity according to 'm';
2523 // n=3=green; n=67+m, m=0,25, increasing luminosity according to 'm';
2524 // n=4=blue; n=117+m, m=0,25, increasing luminosity according to 'm';
2525 // n=5=yellow; n=42+m, m=0,25, increasing luminosity according to 'm';
2526 // n=6=violet; n=142+m, m=0,25, increasing luminosity according to 'm';
2527 // n=7=lightblue; n=92+m, m=0,25, increasing luminosity according to 'm';
2528 // colour=n*10+m, m=1,2,...9, will produce the same colour
2529 // as 'n', but with increasing luminosity according to 'm';
2530 // COLO<0 will act as if abs(COLO) was set for the volume
2531 // and for all the levels below it.
2532 // When for a volume the attribute FILL is > 1 (and the
2533 // option SHAD is on), the ABS of its colour code must be < 8
2534 // because an automatic shading of its faces will be
2537 // FILL (1992) fill area -7,...,0,1,...7 (default=0)
2538 // when option SHAD is "on" the FILL attribute of any
2539 // volume can be set different from 0 (normal drawing);
2540 // if it is set to 1, the faces of such volume will be filled
2541 // with solid colours; if ABS(FILL) is > 1, then a light
2542 // source is placed along the observer line, and the faces of
2543 // such volumes will be painted by colours whose luminosity
2544 // will depend on the amount of light reflected;
2545 // if ABS(FILL) = 1, then it is possible to use all the 166
2546 // colours of the colour table, becouse the automatic shading
2547 // is not performed;
2548 // for increasing values of FILL the drawing will be performed
2549 // with higher and higher resolution improving the quality (the
2550 // number of scan lines used to fill the faces increases with FILL);
2551 // it is possible to set different values of FILL
2552 // for different volumes, in order to optimize at the same time
2553 // the performance and the quality of the picture;
2554 // FILL<0 will act as if abs(FILL) was set for the volume
2555 // and for all the levels below it.
2556 // This kind of drawing can be saved in 'picture files'
2557 // or in view banks.
2558 // 0=drawing without fill area
2559 // 1=faces filled with solid colours and resolution = 6
2560 // 2=lowest resolution (very fast)
2561 // 3=default resolution
2562 // 4=.................
2563 // 5=.................
2564 // 6=.................
2566 // Finally, if a coloured background is desired, the FILL
2567 // attribute for the first volume of the tree must be set
2568 // equal to -abs(colo), colo being >0 and <166.
2570 // SET set number associated to volume name
2571 // DET detector number associated to volume name
2572 // DTYP detector type (1,2)
2579 gsatt(PASSCHARD(vname), PASSCHARD(vatt), val PASSCHARL(vname)
2583 //_____________________________________________________________________________
2584 void TGeant3::Gfpara(const char *name, Int_t number, Int_t intext, Int_t& npar,
2585 Int_t& natt, Float_t* par, Float_t* att)
2588 // Find the parameters of a volume
2590 gfpara(PASSCHARD(name), number, intext, npar, natt, par, att
2594 //_____________________________________________________________________________
2595 void TGeant3::Gckpar(Int_t ish, Int_t npar, Float_t* par)
2598 // Check the parameters of a shape
2600 gckpar(ish,npar,par);
2603 //_____________________________________________________________________________
2604 void TGeant3::Gckmat(Int_t itmed, char* natmed)
2607 // Check the parameters of a tracking medium
2609 gckmat(itmed, PASSCHARD(natmed) PASSCHARL(natmed));
2612 //_____________________________________________________________________________
2613 void TGeant3::Gdelete(Int_t iview)
2616 // IVIEW View number
2618 // It deletes a view bank from memory.
2623 //_____________________________________________________________________________
2624 void TGeant3::Gdopen(Int_t iview)
2627 // IVIEW View number
2629 // When a drawing is very complex and requires a long time to be
2630 // executed, it can be useful to store it in a view bank: after a
2631 // call to DOPEN and the execution of the drawing (nothing will
2632 // appear on the screen), and after a necessary call to DCLOSE,
2633 // the contents of the bank can be displayed in a very fast way
2634 // through a call to DSHOW; therefore, the detector can be easily
2635 // zoomed many times in different ways. Please note that the pictures
2636 // with solid colours can now be stored in a view bank or in 'PICTURE FILES'
2643 //_____________________________________________________________________________
2644 void TGeant3::Gdclose()
2647 // It closes the currently open view bank; it must be called after the
2648 // end of the drawing to be stored.
2653 //_____________________________________________________________________________
2654 void TGeant3::Gdshow(Int_t iview)
2657 // IVIEW View number
2659 // It shows on the screen the contents of a view bank. It
2660 // can be called after a view bank has been closed.
2665 //_____________________________________________________________________________
2666 void TGeant3::Gdopt(const char *name,const char *value)
2670 // VALUE Option value
2672 // To set/modify the drawing options.
2675 // THRZ ON Draw tracks in R vs Z
2676 // OFF (D) Draw tracks in X,Y,Z
2679 // PROJ PARA (D) Parallel projection
2681 // TRAK LINE (D) Trajectory drawn with lines
2682 // POIN " " with markers
2683 // HIDE ON Hidden line removal using the CG package
2684 // OFF (D) No hidden line removal
2685 // SHAD ON Fill area and shading of surfaces.
2686 // OFF (D) Normal hidden line removal.
2687 // RAYT ON Ray-tracing on.
2688 // OFF (D) Ray-tracing off.
2689 // EDGE OFF Does not draw contours when shad is on.
2690 // ON (D) Normal shading.
2691 // MAPP 1,2,3,4 Mapping before ray-tracing.
2692 // 0 (D) No mapping.
2693 // USER ON User graphics options in the raytracing.
2694 // OFF (D) Automatic graphics options.
2700 Vname(value,vvalue);
2701 gdopt(PASSCHARD(vname), PASSCHARD(vvalue) PASSCHARL(vname)
2705 //_____________________________________________________________________________
2706 void TGeant3::Gdraw(const char *name,Float_t theta, Float_t phi, Float_t psi,
2707 Float_t u0,Float_t v0,Float_t ul,Float_t vl)
2712 // THETA Viewing angle theta (for 3D projection)
2713 // PHI Viewing angle phi (for 3D projection)
2714 // PSI Viewing angle psi (for 2D rotation)
2715 // U0 U-coord. (horizontal) of volume origin
2716 // V0 V-coord. (vertical) of volume origin
2717 // SU Scale factor for U-coord.
2718 // SV Scale factor for V-coord.
2720 // This function will draw the volumes,
2721 // selected with their graphical attributes, set by the Gsatt
2722 // facility. The drawing may be performed with hidden line removal
2723 // and with shading effects according to the value of the options HIDE
2724 // and SHAD; if the option SHAD is ON, the contour's edges can be
2725 // drawn or not. If the option HIDE is ON, the detector can be
2726 // exploded (BOMB), clipped with different shapes (CVOL), and some
2727 // of its parts can be shifted from their original
2728 // position (SHIFT). When HIDE is ON, if
2729 // the drawing requires more than the available memory, the program
2730 // will evaluate and display the number of missing words
2731 // (so that the user can increase the
2732 // size of its ZEBRA store). Finally, at the end of each drawing (with HIDE on),
2733 // the program will print messages about the memory used and
2734 // statistics on the volumes' visibility.
2735 // The following commands will produce the drawing of a green
2736 // volume, specified by NAME, without using the hidden line removal
2737 // technique, using the hidden line removal technique,
2738 // with different linewidth and colour (red), with
2739 // solid colour, with shading of surfaces, and without edges.
2740 // Finally, some examples are given for the ray-tracing. (A possible
2741 // string for the NAME of the volume can be found using the command DTREE).
2747 if (fGcvdma->raytra != 1) {
2748 gdraw(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
2750 gdrayt(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
2754 //_____________________________________________________________________________
2755 void TGeant3::Gdrawc(const char *name,Int_t axis, Float_t cut,Float_t u0,
2756 Float_t v0,Float_t ul,Float_t vl)
2761 // CUTVAL Cut plane distance from the origin along the axis
2763 // U0 U-coord. (horizontal) of volume origin
2764 // V0 V-coord. (vertical) of volume origin
2765 // SU Scale factor for U-coord.
2766 // SV Scale factor for V-coord.
2768 // The cut plane is normal to caxis (X,Y,Z), corresponding to iaxis (1,2,3),
2769 // and placed at the distance cutval from the origin.
2770 // The resulting picture is seen from the the same axis.
2771 // When HIDE Mode is ON, it is possible to get the same effect with
2772 // the CVOL/BOX function.
2778 gdrawc(PASSCHARD(vname), axis,cut,u0,v0,ul,vl PASSCHARL(vname));
2781 //_____________________________________________________________________________
2782 void TGeant3::Gdrawx(const char *name,Float_t cutthe, Float_t cutphi,
2783 Float_t cutval, Float_t theta, Float_t phi, Float_t u0,
2784 Float_t v0,Float_t ul,Float_t vl)
2788 // CUTTHE Theta angle of the line normal to cut plane
2789 // CUTPHI Phi angle of the line normal to cut plane
2790 // CUTVAL Cut plane distance from the origin along the axis
2792 // THETA Viewing angle theta (for 3D projection)
2793 // PHI Viewing angle phi (for 3D projection)
2794 // U0 U-coord. (horizontal) of volume origin
2795 // V0 V-coord. (vertical) of volume origin
2796 // SU Scale factor for U-coord.
2797 // SV Scale factor for V-coord.
2799 // The cut plane is normal to the line given by the cut angles
2800 // cutthe and cutphi and placed at the distance cutval from the origin.
2801 // The resulting picture is seen from the viewing angles theta,phi.
2807 gdrawx(PASSCHARD(vname), cutthe,cutphi,cutval,theta,phi,u0,v0,ul,vl
2811 //_____________________________________________________________________________
2812 void TGeant3::Gdhead(Int_t isel, const char *name, Float_t chrsiz)
2817 // ISEL Option flag D=111110
2819 // CHRSIZ Character size (cm) of title NAME D=0.6
2822 // 0 to have only the header lines
2823 // xxxxx1 to add the text name centered on top of header
2824 // xxxx1x to add global detector name (first volume) on left
2825 // xxx1xx to add date on right
2826 // xx1xxx to select thick characters for text on top of header
2827 // x1xxxx to add the text 'EVENT NR x' on top of header
2828 // 1xxxxx to add the text 'RUN NR x' on top of header
2829 // NOTE that ISEL=x1xxx1 or ISEL=1xxxx1 are illegal choices,
2830 // i.e. they generate overwritten text.
2832 gdhead(isel,PASSCHARD(name),chrsiz PASSCHARL(name));
2835 //_____________________________________________________________________________
2836 void TGeant3::Gdman(Float_t u, Float_t v, const char *type)
2839 // Draw a 2D-man at position (U0,V0)
2841 // U U-coord. (horizontal) of the centre of man' R
2842 // V V-coord. (vertical) of the centre of man' R
2843 // TYPE D='MAN' possible values: 'MAN,WM1,WM2,WM3'
2845 // CALL GDMAN(u,v),CALL GDWMN1(u,v),CALL GDWMN2(u,v),CALL GDWMN2(u,v)
2846 // It superimposes the picure of a man or of a woman, chosen among
2847 // three different ones, with the same scale factors as the detector
2848 // in the current drawing.
2851 if (opt.Contains("WM1")) {
2853 } else if (opt.Contains("WM3")) {
2855 } else if (opt.Contains("WM2")) {
2862 //_____________________________________________________________________________
2863 void TGeant3::Gdspec(const char *name)
2868 // Shows 3 views of the volume (two cut-views and a 3D view), together with
2869 // its geometrical specifications. The 3D drawing will
2870 // be performed according the current values of the options HIDE and
2871 // SHAD and according the current SetClipBox clipping parameters for that
2878 gdspec(PASSCHARD(vname) PASSCHARL(vname));
2881 //_____________________________________________________________________________
2882 void TGeant3::DrawOneSpec(const char *name)
2885 // Function called when one double-clicks on a volume name
2886 // in a TPavelabel drawn by Gdtree.
2888 THIGZ *higzSave = higz;
2889 higzSave->SetName("higzSave");
2890 THIGZ *higzSpec = (THIGZ*)gROOT->FindObject("higzSpec");
2891 //printf("DrawOneSpec, higz=%x, higzSpec=%x\n",higz,higzSpec);
2892 if (higzSpec) higz = higzSpec;
2893 else higzSpec = new THIGZ(defSize);
2894 higzSpec->SetName("higzSpec");
2899 gdspec(PASSCHARD(vname) PASSCHARL(vname));
2902 higzSave->SetName("higz");
2906 //_____________________________________________________________________________
2907 void TGeant3::Gdtree(const char *name,Int_t levmax, Int_t isel)
2911 // LEVMAX Depth level
2914 // This function draws the logical tree,
2915 // Each volume in the tree is represented by a TPaveTree object.
2916 // Double-clicking on a TPaveTree draws the specs of the corresponding volume.
2917 // Use TPaveTree pop-up menu to select:
2920 // - drawing tree of parent
2926 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
2930 //_____________________________________________________________________________
2931 void TGeant3::GdtreeParent(const char *name,Int_t levmax, Int_t isel)
2935 // LEVMAX Depth level
2938 // This function draws the logical tree of the parent of name.
2942 // Scan list of volumes in JVOLUM
2944 Int_t gname, i, jvo, in, nin, jin, num;
2945 strncpy((char *) &gname, name, 4);
2946 for(i=1; i<=fGcnum->nvolum; i++) {
2947 jvo = fZlq[fGclink->jvolum-i];
2948 nin = Int_t(fZq[jvo+3]);
2949 if (nin == -1) nin = 1;
2950 for (in=1;in<=nin;in++) {
2952 num = Int_t(fZq[jin+2]);
2953 if(gname == fZiq[fGclink->jvolum+num]) {
2954 strncpy(vname,(char*)&fZiq[fGclink->jvolum+i],4);
2956 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
2964 //_____________________________________________________________________________
2965 void TGeant3::SetABAN(Int_t par)
2968 // par = 1 particles will be stopped according to their residual
2969 // range if they are not in a sensitive material and are
2970 // far enough from the boundary
2971 // 0 particles are transported normally
2973 fGcphys->dphys1 = par;
2977 //_____________________________________________________________________________
2978 void TGeant3::SetANNI(Int_t par)
2981 // To control positron annihilation.
2982 // par =0 no annihilation
2983 // =1 annihilation. Decays processed.
2984 // =2 annihilation. No decay products stored.
2986 fGcphys->ianni = par;
2990 //_____________________________________________________________________________
2991 void TGeant3::SetAUTO(Int_t par)
2994 // To control automatic calculation of tracking medium parameters:
2995 // par =0 no automatic calculation;
2996 // =1 automati calculation.
2998 fGctrak->igauto = par;
3002 //_____________________________________________________________________________
3003 void TGeant3::SetBOMB(Float_t boom)
3006 // BOOM : Exploding factor for volumes position
3008 // To 'explode' the detector. If BOOM is positive (values smaller
3009 // than 1. are suggested, but any value is possible)
3010 // all the volumes are shifted by a distance
3011 // proportional to BOOM along the direction between their centre
3012 // and the origin of the MARS; the volumes which are symmetric
3013 // with respect to this origin are simply not shown.
3014 // BOOM equal to 0 resets the normal mode.
3015 // A negative (greater than -1.) value of
3016 // BOOM will cause an 'implosion'; for even lower values of BOOM
3017 // the volumes' positions will be reflected respect to the origin.
3018 // This command can be useful to improve the 3D effect for very
3019 // complex detectors. The following commands will make explode the
3026 //_____________________________________________________________________________
3027 void TGeant3::SetBREM(Int_t par)
3030 // To control bremstrahlung.
3031 // par =0 no bremstrahlung
3032 // =1 bremstrahlung. Photon processed.
3033 // =2 bremstrahlung. No photon stored.
3035 fGcphys->ibrem = par;
3039 //_____________________________________________________________________________
3040 void TGeant3::SetCKOV(Int_t par)
3043 // To control Cerenkov production
3044 // par =0 no Cerenkov;
3046 // =2 Cerenkov with primary stopped at each step.
3048 fGctlit->itckov = par;
3052 //_____________________________________________________________________________
3053 void TGeant3::SetClipBox(const char *name,Float_t xmin,Float_t xmax,
3054 Float_t ymin,Float_t ymax,Float_t zmin,Float_t zmax)
3057 // The hidden line removal technique is necessary to visualize properly
3058 // very complex detectors. At the same time, it can be useful to visualize
3059 // the inner elements of a detector in detail. This function allows
3060 // subtractions (via boolean operation) of BOX shape from any part of
3061 // the detector, therefore showing its inner contents.
3062 // If "*" is given as the name of the
3063 // volume to be clipped, all volumes are clipped by the given box.
3064 // A volume can be clipped at most twice.
3065 // if a volume is explicitely clipped twice,
3066 // the "*" will not act on it anymore. Giving "." as the name
3067 // of the volume to be clipped will reset the clipping.
3069 // NAME Name of volume to be clipped
3071 // XMIN Lower limit of the Shape X coordinate
3072 // XMAX Upper limit of the Shape X coordinate
3073 // YMIN Lower limit of the Shape Y coordinate
3074 // YMAX Upper limit of the Shape Y coordinate
3075 // ZMIN Lower limit of the Shape Z coordinate
3076 // ZMAX Upper limit of the Shape Z coordinate
3078 // This function performs a boolean subtraction between the volume
3079 // NAME and a box placed in the MARS according the values of the given
3085 setclip(PASSCHARD(vname),xmin,xmax,ymin,ymax,zmin,zmax PASSCHARL(vname));
3088 //_____________________________________________________________________________
3089 void TGeant3::SetCOMP(Int_t par)
3092 // To control Compton scattering
3093 // par =0 no Compton
3094 // =1 Compton. Electron processed.
3095 // =2 Compton. No electron stored.
3098 fGcphys->icomp = par;
3101 //_____________________________________________________________________________
3102 void TGeant3::SetCUTS(Float_t cutgam,Float_t cutele,Float_t cutneu,
3103 Float_t cuthad,Float_t cutmuo ,Float_t bcute ,
3104 Float_t bcutm ,Float_t dcute ,Float_t dcutm ,
3105 Float_t ppcutm, Float_t tofmax)
3108 // CUTGAM Cut for gammas D=0.001
3109 // CUTELE Cut for electrons D=0.001
3110 // CUTHAD Cut for charged hadrons D=0.01
3111 // CUTNEU Cut for neutral hadrons D=0.01
3112 // CUTMUO Cut for muons D=0.01
3113 // BCUTE Cut for electron brems. D=-1.
3114 // BCUTM Cut for muon brems. D=-1.
3115 // DCUTE Cut for electron delta-rays D=-1.
3116 // DCUTM Cut for muon delta-rays D=-1.
3117 // PPCUTM Cut for e+e- pairs by muons D=0.01
3118 // TOFMAX Time of flight cut D=1.E+10
3120 // If the default values (-1.) for BCUTE ,BCUTM ,DCUTE ,DCUTM
3121 // are not modified, they will be set to CUTGAM,CUTGAM,CUTELE,CUTELE
3123 // If one of the parameters from CUTGAM to PPCUTM included
3124 // is modified, cross-sections and energy loss tables must be
3125 // recomputed via the function Gphysi.
3127 fGccuts->cutgam = cutgam;
3128 fGccuts->cutele = cutele;
3129 fGccuts->cutneu = cutneu;
3130 fGccuts->cuthad = cuthad;
3131 fGccuts->cutmuo = cutmuo;
3132 fGccuts->bcute = bcute;
3133 fGccuts->bcutm = bcutm;
3134 fGccuts->dcute = dcute;
3135 fGccuts->dcutm = dcutm;
3136 fGccuts->ppcutm = ppcutm;
3137 fGccuts->tofmax = tofmax;
3140 //_____________________________________________________________________________
3141 void TGeant3::SetDCAY(Int_t par)
3144 // To control Decay mechanism.
3145 // par =0 no decays.
3146 // =1 Decays. secondaries processed.
3147 // =2 Decays. No secondaries stored.
3149 fGcphys->idcay = par;
3153 //_____________________________________________________________________________
3154 void TGeant3::SetDEBU(Int_t emin, Int_t emax, Int_t emod)
3157 // Set the debug flag and frequency
3158 // Selected debug output will be printed from
3159 // event emin to even emax each emod event
3161 fGcflag->idemin = emin;
3162 fGcflag->idemax = emax;
3163 fGcflag->itest = emod;
3167 //_____________________________________________________________________________
3168 void TGeant3::SetDRAY(Int_t par)
3171 // To control delta rays mechanism.
3172 // par =0 no delta rays.
3173 // =1 Delta rays. secondaries processed.
3174 // =2 Delta rays. No secondaries stored.
3176 fGcphys->idray = par;
3179 //_____________________________________________________________________________
3180 void TGeant3::SetHADR(Int_t par)
3183 // To control hadronic interactions.
3184 // par =0 no hadronic interactions.
3185 // =1 Hadronic interactions. secondaries processed.
3186 // =2 Hadronic interactions. No secondaries stored.
3188 fGcphys->ihadr = par;
3191 //_____________________________________________________________________________
3192 void TGeant3::SetKINE(Int_t kine, Float_t xk1, Float_t xk2, Float_t xk3,
3193 Float_t xk4, Float_t xk5, Float_t xk6, Float_t xk7,
3194 Float_t xk8, Float_t xk9, Float_t xk10)
3197 // Set the variables in /GCFLAG/ IKINE, PKINE(10)
3198 // Their meaning is user defined
3200 fGckine->ikine = kine;
3201 fGckine->pkine[0] = xk1;
3202 fGckine->pkine[1] = xk2;
3203 fGckine->pkine[2] = xk3;
3204 fGckine->pkine[3] = xk4;
3205 fGckine->pkine[4] = xk5;
3206 fGckine->pkine[5] = xk6;
3207 fGckine->pkine[6] = xk7;
3208 fGckine->pkine[7] = xk8;
3209 fGckine->pkine[8] = xk9;
3210 fGckine->pkine[9] = xk10;
3213 //_____________________________________________________________________________
3214 void TGeant3::SetLOSS(Int_t par)
3217 // To control energy loss.
3218 // par =0 no energy loss;
3219 // =1 restricted energy loss fluctuations;
3220 // =2 complete energy loss fluctuations;
3222 // =4 no energy loss fluctuations.
3223 // If the value ILOSS is changed, then cross-sections and energy loss
3224 // tables must be recomputed via the command 'PHYSI'.
3226 fGcphys->iloss = par;
3230 //_____________________________________________________________________________
3231 void TGeant3::SetMULS(Int_t par)
3234 // To control multiple scattering.
3235 // par =0 no multiple scattering.
3236 // =1 Moliere or Coulomb scattering.
3237 // =2 Moliere or Coulomb scattering.
3238 // =3 Gaussian scattering.
3240 fGcphys->imuls = par;
3244 //_____________________________________________________________________________
3245 void TGeant3::SetMUNU(Int_t par)
3248 // To control muon nuclear interactions.
3249 // par =0 no muon-nuclear interactions.
3250 // =1 Nuclear interactions. Secondaries processed.
3251 // =2 Nuclear interactions. Secondaries not processed.
3253 fGcphys->imunu = par;
3256 //_____________________________________________________________________________
3257 void TGeant3::SetOPTI(Int_t par)
3260 // This flag controls the tracking optimisation performed via the
3262 // 1 no optimisation at all; GSORD calls disabled;
3263 // 0 no optimisation; only user calls to GSORD kept;
3264 // 1 all non-GSORDered volumes are ordered along the best axis;
3265 // 2 all volumes are ordered along the best axis.
3267 fGcopti->ioptim = par;
3270 //_____________________________________________________________________________
3271 void TGeant3::SetPAIR(Int_t par)
3274 // To control pair production mechanism.
3275 // par =0 no pair production.
3276 // =1 Pair production. secondaries processed.
3277 // =2 Pair production. No secondaries stored.
3279 fGcphys->ipair = par;
3283 //_____________________________________________________________________________
3284 void TGeant3::SetPFIS(Int_t par)
3287 // To control photo fission mechanism.
3288 // par =0 no photo fission.
3289 // =1 Photo fission. secondaries processed.
3290 // =2 Photo fission. No secondaries stored.
3292 fGcphys->ipfis = par;
3295 //_____________________________________________________________________________
3296 void TGeant3::SetPHOT(Int_t par)
3299 // To control Photo effect.
3300 // par =0 no photo electric effect.
3301 // =1 Photo effect. Electron processed.
3302 // =2 Photo effect. No electron stored.
3304 fGcphys->iphot = par;
3307 //_____________________________________________________________________________
3308 void TGeant3::SetRAYL(Int_t par)
3311 // To control Rayleigh scattering.
3312 // par =0 no Rayleigh scattering.
3315 fGcphys->irayl = par;
3318 //_____________________________________________________________________________
3319 void TGeant3::SetSWIT(Int_t sw, Int_t val)
3323 // val New switch value
3325 // Change one element of array ISWIT(10) in /GCFLAG/
3327 if (sw <= 0 || sw > 10) return;
3328 fGcflag->iswit[sw-1] = val;
3332 //_____________________________________________________________________________
3333 void TGeant3::SetTRIG(Int_t nevents)
3336 // Set number of events to be run
3338 fGcflag->nevent = nevents;
3341 //_____________________________________________________________________________
3342 void TGeant3::SetUserDecay(Int_t pdg)
3345 // Force the decays of particles to be done with Pythia
3346 // and not with the Geant routines.
3347 // just kill pointers doing mzdrop
3349 Int_t ipart = IdFromPDG(pdg);
3351 printf("Particle %d not in geant\n",pdg);
3354 Int_t jpart=fGclink->jpart;
3355 Int_t jpa=fZlq[jpart-ipart];
3358 Int_t jpa1=fZlq[jpa-1];
3360 mzdrop(fGcbank->ixcons,jpa1,PASSCHARD(" ") PASSCHARL(" "));
3361 Int_t jpa2=fZlq[jpa-2];
3363 mzdrop(fGcbank->ixcons,jpa2,PASSCHARD(" ") PASSCHARL(" "));
3367 //______________________________________________________________________________
3368 void TGeant3::Vname(const char *name, char *vname)
3371 // convert name to upper case. Make vname at least 4 chars
3373 Int_t l = strlen(name);
3376 for (i=0;i<l;i++) vname[i] = toupper(name[i]);
3377 for (i=l;i<4;i++) vname[i] = ' ';
3381 //______________________________________________________________________________
3382 void TGeant3::Ertrgo()
3387 //______________________________________________________________________________
3388 void TGeant3::Ertrak(const Float_t *const x1, const Float_t *const p1,
3389 const Float_t *x2, const Float_t *p2,
3390 Int_t ipa, Option_t *chopt)
3392 ertrak(x1,p1,x2,p2,ipa,PASSCHARD(chopt) PASSCHARL(chopt));
3395 //_____________________________________________________________________________
3396 void TGeant3::WriteEuclid(const char* filnam, const char* topvol,
3397 Int_t number, Int_t nlevel)
3401 // ******************************************************************
3403 // * Write out the geometry of the detector in EUCLID file format *
3405 // * filnam : will be with the extension .euc *
3406 // * topvol : volume name of the starting node *
3407 // * number : copy number of topvol (relevant for gsposp) *
3408 // * nlevel : number of levels in the tree structure *
3409 // * to be written out, starting from topvol *
3411 // * Author : M. Maire *
3413 // ******************************************************************
3415 // File filnam.tme is written out with the definitions of tracking
3416 // medias and materials.
3417 // As to restore original numbers for materials and medias, program
3418 // searches in the file euc_medi.dat and comparing main parameters of
3419 // the mat. defined inside geant and the one in file recognizes them
3420 // and is able to take number from file. If for any material or medium,
3421 // this procedure fails, ordering starts from 1.
3422 // Arrays IOTMED and IOMATE are used for this procedure
3424 const char shape[][5]={"BOX ","TRD1","TRD2","TRAP","TUBE","TUBS","CONE",
3425 "CONS","SPHE","PARA","PGON","PCON","ELTU","HYPE",
3427 Int_t i, end, itm, irm, jrm, k, nmed;
3431 char *filext, *filetme;
3432 char natmed[21], namate[21];
3433 char natmedc[21], namatec[21];
3434 char key[5], name[5], mother[5], konly[5];
3436 Int_t iadvol, iadtmd, iadrot, nwtot, iret;
3437 Int_t mlevel, numbr, natt, numed, nin, ndata;
3438 Int_t iname, ivo, ish, jvo, nvstak, ivstak;
3439 Int_t jdiv, ivin, in, jin, jvin, irot;
3440 Int_t jtm, imat, jma, flag=0, imatc;
3441 Float_t az, dens, radl, absl, a, step, x, y, z;
3442 Int_t npar, ndvmx, left;
3443 Float_t zc, densc, radlc, abslc, c0, tmaxfd;
3445 Int_t iomate[100], iotmed[100];
3446 Float_t par[50], att[20], ubuf[50];
3449 Int_t level, ndiv, iaxe;
3450 Int_t itmedc, nmatc, isvolc, ifieldc, nwbufc, isvol, nmat, ifield, nwbuf;
3451 Float_t fieldmc, tmaxfdc, stemaxc, deemaxc, epsilc, stminc, fieldm;
3452 Float_t tmaxf, stemax, deemax, epsil, stmin;
3453 const char *f10000="!\n%s\n!\n";
3454 //Open the input file
3456 for(i=0;i<end;i++) if(filnam[i]=='.') {
3460 filext=new char[end+5];
3461 filetme=new char[end+5];
3462 strncpy(filext,filnam,end);
3463 strncpy(filetme,filnam,end);
3465 // *** The output filnam name will be with extension '.euc'
3466 strcpy(&filext[end],".euc");
3467 strcpy(&filetme[end],".tme");
3468 lun=fopen(filext,"w");
3470 // *** Initialisation of the working space
3471 iadvol=fGcnum->nvolum;
3472 iadtmd=iadvol+fGcnum->nvolum;
3473 iadrot=iadtmd+fGcnum->ntmed;
3474 if(fGclink->jrotm) {
3475 fGcnum->nrotm=fZiq[fGclink->jrotm-2];
3479 nwtot=iadrot+fGcnum->nrotm;
3480 qws = new float[nwtot+1];
3481 for (i=0;i<nwtot+1;i++) qws[i]=0;
3484 if(nlevel==0) mlevel=20;
3486 // *** find the top volume and put it in the stak
3487 numbr = number>0 ? number : 1;
3488 Gfpara(topvol,numbr,1,npar,natt,par,att);
3490 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3495 // *** authorized shape ?
3496 strncpy((char *)&iname, topvol, 4);
3498 for(i=1; i<=fGcnum->nvolum; i++) if(fZiq[fGclink->jvolum+i]==iname) {
3502 jvo = fZlq[fGclink->jvolum-ivo];
3503 ish = Int_t (fZq[jvo+2]);
3505 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3512 iws[iadvol+ivo] = level;
3515 //*** flag all volumes and fill the stak
3519 // pick the next volume in stak
3521 ivo = TMath::Abs(iws[ivstak]);
3522 jvo = fZlq[fGclink->jvolum - ivo];
3524 // flag the tracking medium
3525 numed = Int_t (fZq[jvo + 4]);
3526 iws[iadtmd + numed] = 1;
3528 // get the daughters ...
3529 level = iws[iadvol+ivo];
3530 if (level < mlevel) {
3532 nin = Int_t (fZq[jvo + 3]);
3534 // from division ...
3536 jdiv = fZlq[jvo - 1];
3537 ivin = Int_t (fZq[jdiv + 2]);
3539 iws[nvstak] = -ivin;
3540 iws[iadvol+ivin] = level;
3542 // from position ...
3543 } else if (nin > 0) {
3544 for(in=1; in<=nin; in++) {
3545 jin = fZlq[jvo - in];
3546 ivin = Int_t (fZq[jin + 2 ]);
3547 jvin = fZlq[fGclink->jvolum - ivin];
3548 ish = Int_t (fZq[jvin + 2]);
3549 // authorized shape ?
3551 // not yet flagged ?
3552 if (iws[iadvol+ivin]==0) {
3555 iws[iadvol+ivin] = level;
3557 // flag the rotation matrix
3558 irot = Int_t ( fZq[jin + 4 ]);
3559 if (irot > 0) iws[iadrot+irot] = 1;
3565 // next volume in stak ?
3566 if (ivstak < nvstak) goto L10;
3568 // *** restore original material and media numbers
3569 // file euc_medi.dat is needed to compare materials and medias
3571 FILE* luncor=fopen("euc_medi.dat","r");
3574 for(itm=1; itm<=fGcnum->ntmed; itm++) {
3575 if (iws[iadtmd+itm] > 0) {
3576 jtm = fZlq[fGclink->jtmed-itm];
3577 strncpy(natmed,(char *)&fZiq[jtm+1],20);
3578 imat = Int_t (fZq[jtm+6]);
3579 jma = fZlq[fGclink->jmate-imat];
3581 printf(" *** GWEUCL *** material not defined for tracking medium %5i %s\n",itm,natmed);
3584 strncpy(namate,(char *)&fZiq[jma+1],20);
3587 //** find the material original number
3590 iret=fscanf(luncor,"%4s,%130s",key,card);
3591 if(iret<=0) goto L26;
3593 if(!strcmp(key,"MATE")) {
3594 sscanf(card,"%d %s %f %f %f %f %f %d",&imatc,namatec,&az,&zc,&densc,&radlc,&abslc,&nparc);
3595 Gfmate(imat,namate,a,z,dens,radl,absl,par,npar);
3596 if(!strcmp(namatec,namate)) {
3597 if(az==a && zc==z && densc==dens && radlc==radl
3598 && abslc==absl && nparc==nparc) {
3601 printf("*** GWEUCL *** material : %3d '%s' restored as %3d\n",imat,namate,imatc);
3603 printf("*** GWEUCL *** different definitions for material: %s\n",namate);
3607 if(strcmp(key,"END") && !flag) goto L23;
3609 printf("*** GWEUCL *** cannot restore original number for material: %s\n",namate);
3613 //*** restore original tracking medium number
3616 iret=fscanf(luncor,"%4s,%130s",key,card);
3617 if(iret<=0) goto L26;
3619 if (!strcmp(key,"TMED")) {
3620 sscanf(card,"%d %s %d %d %d %f %f %f %f %f %f %d\n",
3621 &itmedc,natmedc,&nmatc,&isvolc,&ifieldc,&fieldmc,
3622 &tmaxfdc,&stemaxc,&deemaxc,&epsilc,&stminc,&nwbufc);
3623 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxf,stemax,deemax,
3624 epsil,stmin,ubuf,&nwbuf);
3625 if(!strcmp(natmedc,natmed)) {
3626 if (iomate[nmat]==nmatc && nwbuf==nwbufc) {
3629 printf("*** GWEUCL *** medium : %3d '%20s' restored as %3d\n",
3632 printf("*** GWEUCL *** different definitions for tracking medium: %s\n",natmed);
3636 if(strcmp(key,"END") && !flag) goto L24;
3638 printf("cannot restore original number for medium : %s\n",natmed);
3646 L26: printf("*** GWEUCL *** cannot read the data file\n");
3648 L29: if(luncor) fclose (luncor);
3651 // *** write down the tracking medium definition
3653 strcpy(card,"! Tracking medium");
3654 fprintf(lun,f10000,card);
3656 for(itm=1;itm<=fGcnum->ntmed;itm++) {
3657 if (iws[iadtmd+itm]>0) {
3658 jtm = fZlq[fGclink->jtmed-itm];
3659 strncpy(natmed,(char *)&fZiq[jtm+1],20);
3661 imat = Int_t (fZq[jtm+6]);
3662 jma = fZlq[fGclink->jmate-imat];
3663 //* order media from one, if comparing with database failed
3665 iotmed[itm]=++imxtmed;
3666 iomate[imat]=++imxmate;
3671 printf(" *** GWEUCL *** material not defined for tracking medium %5d %s\n",
3674 strncpy(namate,(char *)&fZiq[jma+1],20);
3677 fprintf(lun,"TMED %3d '%20s' %3d '%20s'\n",iotmed[itm],natmed,iomate[imat],namate);
3681 //* *** write down the rotation matrix
3683 strcpy(card,"! Reperes");
3684 fprintf(lun,f10000,card);
3686 for(irm=1;irm<=fGcnum->nrotm;irm++) {
3687 if (iws[iadrot+irm]>0) {
3688 jrm = fZlq[fGclink->jrotm-irm];
3689 fprintf(lun,"ROTM %3d",irm);
3690 for(k=11;k<=16;k++) fprintf(lun," %8.3f",fZq[jrm+k]);
3695 //* *** write down the volume definition
3697 strcpy(card,"! Volumes");
3698 fprintf(lun,f10000,card);
3700 for(ivstak=1;ivstak<=nvstak;ivstak++) {
3703 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivo],4);
3705 jvo = fZlq[fGclink->jvolum-ivo];
3706 ish = Int_t (fZq[jvo+2]);
3707 nmed = Int_t (fZq[jvo+4]);
3708 npar = Int_t (fZq[jvo+5]);
3710 if (ivstak>1) for(i=0;i<npar;i++) par[i]=fZq[jvo+7+i];
3711 Gckpar (ish,npar,par);
3712 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,shape[ish-1],iotmed[nmed],npar);
3713 for(i=0;i<(npar-1)/6+1;i++) {
3716 for(k=0;k<(left<6?left:6);k++) fprintf(lun," %11.5f",par[i*6+k]);
3720 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,shape[ish-1],iotmed[nmed],npar);
3725 //* *** write down the division of volumes
3727 fprintf(lun,f10000,"! Divisions");
3728 for(ivstak=1;ivstak<=nvstak;ivstak++) {
3729 ivo = TMath::Abs(iws[ivstak]);
3730 jvo = fZlq[fGclink->jvolum-ivo];
3731 ish = Int_t (fZq[jvo+2]);
3732 nin = Int_t (fZq[jvo+3]);
3733 //* this volume is divided ...
3736 iaxe = Int_t ( fZq[jdiv+1]);
3737 ivin = Int_t ( fZq[jdiv+2]);
3738 ndiv = Int_t ( fZq[jdiv+3]);
3741 jvin = fZlq[fGclink->jvolum-ivin];
3742 nmed = Int_t ( fZq[jvin+4]);
3743 strncpy(mother,(char *)&fZiq[fGclink->jvolum+ivo ],4);
3745 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivin],4);
3747 if ((step<=0.)||(ish>=11)) {
3748 //* volume with negative parameter or gsposp or pgon ...
3749 fprintf(lun,"DIVN '%4s' '%4s' %3d %3d\n",name,mother,ndiv,iaxe);
3750 } else if ((ndiv<=0)||(ish==10)) {
3751 //* volume with negative parameter or gsposp or para ...
3752 ndvmx = TMath::Abs(ndiv);
3753 fprintf(lun,"DIVT '%4s' '%4s' %11.5f %3d %3d %3d\n",
3754 name,mother,step,iaxe,iotmed[nmed],ndvmx);
3756 //* normal volume : all kind of division are equivalent
3757 fprintf(lun,"DVT2 '%4s' '%4s' %11.5f %3d %11.5f %3d %3d\n",
3758 name,mother,step,iaxe,c0,iotmed[nmed],ndiv);
3763 //* *** write down the the positionnement of volumes
3765 fprintf(lun,f10000,"! Positionnements\n");
3767 for(ivstak = 1;ivstak<=nvstak;ivstak++) {
3768 ivo = TMath::Abs(iws[ivstak]);
3769 strncpy(mother,(char*)&fZiq[fGclink->jvolum+ivo ],4);
3771 jvo = fZlq[fGclink->jvolum-ivo];
3772 nin = Int_t( fZq[jvo+3]);
3773 //* this volume has daughters ...
3775 for (in=1;in<=nin;in++) {
3777 ivin = Int_t (fZq[jin +2]);
3778 numb = Int_t (fZq[jin +3]);
3779 irot = Int_t (fZq[jin +4]);
3783 strcpy(konly,"ONLY");
3784 if (fZq[jin+8]!=1.) strcpy(konly,"MANY");
3785 strncpy(name,(char*)&fZiq[fGclink->jvolum+ivin],4);
3787 jvin = fZlq[fGclink->jvolum-ivin];
3788 ish = Int_t (fZq[jvin+2]);
3789 //* gspos or gsposp ?
3790 ndata = fZiq[jin-1];
3792 fprintf(lun,"POSI '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s'\n",
3793 name,numb,mother,x,y,z,irot,konly);
3795 npar = Int_t (fZq[jin+9]);
3796 for(i=0;i<npar;i++) par[i]=fZq[jin+10+i];
3797 Gckpar (ish,npar,par);
3798 fprintf(lun,"POSP '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s' %3d\n",
3799 name,numb,mother,x,y,z,irot,konly,npar);
3801 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
3808 fprintf(lun,"END\n");
3811 //****** write down the materials and medias *****
3813 lun=fopen(filetme,"w");
3815 for(itm=1;itm<=fGcnum->ntmed;itm++) {
3816 if (iws[iadtmd+itm]>0) {
3817 jtm = fZlq[fGclink->jtmed-itm];
3818 strncpy(natmed,(char*)&fZiq[jtm+1],4);
3819 imat = Int_t (fZq[jtm+6]);
3820 jma = Int_t (fZlq[fGclink->jmate-imat]);
3822 Gfmate (imat,namate,a,z,dens,radl,absl,par,npar);
3823 fprintf(lun,"MATE %4d '%20s'%11.5E %11.5E %11.5E %11.5E %11.5E %3d\n",
3824 iomate[imat],namate,a,z,dens,radl,absl,npar);
3828 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
3832 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxfd,stemax,deemax,epsil,stmin,par,&npar);
3833 fprintf(lun,"TMED %4d '%20s' %3d %1d %3d %11.5f %11.5f %11.5f %11.5f %11.5f %11.5f %3d\n",
3834 iotmed[itm],natmed,iomate[nmat],isvol,ifield,
3835 fieldm,tmaxfd,stemax,deemax,epsil,stmin,npar);
3839 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
3845 fprintf(lun,"END\n");
3846 printf(" *** GWEUCL *** file: %s is now written out\n",filext);
3847 printf(" *** GWEUCL *** file: %s is now written out\n",filetme);
3856 //_____________________________________________________________________________
3857 void TGeant3::Streamer(TBuffer &R__b)
3860 // Stream an object of class TGeant3.
3862 if (R__b.IsReading()) {
3863 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
3864 AliMC::Streamer(R__b);
3867 R__b.ReadStaticArray(fPDGCode);
3869 R__b.WriteVersion(TGeant3::IsA());
3870 AliMC::Streamer(R__b);
3873 R__b.WriteArray(fPDGCode, fNPDGCodes);