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 **************************************************************************/
20 ///////////////////////////////////////////////////////////////////////////////
22 // Interface Class to the Geant3.21 MonteCarlo //
26 <img src="picts/TGeant3Class.gif">
31 ///////////////////////////////////////////////////////////////////////////////
37 #include <TDatabasePDG.h>
38 #include "AliCallf77.h"
41 # define gzebra gzebra_
42 # define grfile grfile_
43 # define gpcxyz gpcxyz_
44 # define ggclos ggclos_
47 # define gcinit gcinit_
50 # define gtrigc gtrigc_
51 # define gtrigi gtrigi_
53 # define gzinit gzinit_
54 # define gfmate gfmate_
55 # define gfpart gfpart_
56 # define gftmed gftmed_
60 # define gsmate gsmate_
61 # define gsmixt gsmixt_
62 # define gspart gspart_
63 # define gstmed gstmed_
64 # define gsckov gsckov_
65 # define gstpar gstpar_
66 # define gfkine gfkine_
67 # define gfvert gfvert_
68 # define gskine gskine_
69 # define gsvert gsvert_
70 # define gphysi gphysi_
71 # define gdebug gdebug_
72 # define gekbin gekbin_
73 # define gfinds gfinds_
74 # define gsking gsking_
75 # define gskpho gskpho_
76 # define gsstak gsstak_
78 # define gtrack gtrack_
79 # define gtreve gtreve_
80 # define gtreve_root gtreve_root_
82 # define grndmq grndmq_
84 # define glmoth glmoth_
85 # define gmedia gmedia_
88 # define gsdvn2 gsdvn2_
90 # define gsdvs2 gsdvs2_
92 # define gsdvt2 gsdvt2_
95 # define gsposp gsposp_
96 # define gsrotm gsrotm_
97 # define gprotm gprotm_
98 # define gsvolu gsvolu_
99 # define gprint gprint_
100 # define gdinit gdinit_
101 # define gdopt gdopt_
102 # define gdraw gdraw_
103 # define gdrayt gdrayt_
104 # define gdrawc gdrawc_
105 # define gdrawx gdrawx_
106 # define gdhead gdhead_
107 # define gdwmn1 gdwmn1_
108 # define gdwmn2 gdwmn2_
109 # define gdwmn3 gdwmn3_
110 # define gdxyz gdxyz_
111 # define gdcxyz gdcxyz_
112 # define gdman gdman_
113 # define gdspec gdspec_
114 # define gdtree gdtree_
115 # define gdelet gdelet_
116 # define gdclos gdclos_
117 # define gdshow gdshow_
118 # define gdopen gdopen_
119 # define dzshow dzshow_
120 # define gsatt gsatt_
121 # define gfpara gfpara_
122 # define gckpar gckpar_
123 # define gckmat gckmat_
124 # define geditv geditv_
125 # define mzdrop mzdrop_
127 # define ertrak ertrak_
128 # define ertrgo ertrgo_
130 # define setbomb setbomb_
131 # define setclip setclip_
132 # define gcomad gcomad_
135 # define gzebra GZEBRA
136 # define grfile GRFILE
137 # define gpcxyz GPCXYZ
138 # define ggclos GGCLOS
141 # define gcinit GCINIT
144 # define gtrigc GTRIGC
145 # define gtrigi GTRIGI
147 # define gzinit GZINIT
148 # define gfmate GFMATE
149 # define gfpart GFPART
150 # define gftmed GFTMED
154 # define gsmate GSMATE
155 # define gsmixt GSMIXT
156 # define gspart GSPART
157 # define gstmed GSTMED
158 # define gsckov GSCKOV
159 # define gstpar GSTPAR
160 # define gfkine GFKINE
161 # define gfvert GFVERT
162 # define gskine GSKINE
163 # define gsvert GSVERT
164 # define gphysi GPHYSI
165 # define gdebug GDEBUG
166 # define gekbin GEKBIN
167 # define gfinds GFINDS
168 # define gsking GSKING
169 # define gskpho GSKPHO
170 # define gsstak GSSTAK
172 # define gtrack GTRACK
173 # define gtreve GTREVE
174 # define gtreve_root GTREVE_ROOT
176 # define grndmq GRNDMQ
178 # define glmoth GLMOTH
179 # define gmedia GMEDIA
182 # define gsdvn2 GSDVN2
184 # define gsdvs2 GSDVS2
186 # define gsdvt2 GSDVT2
189 # define gsposp GSPOSP
190 # define gsrotm GSROTM
191 # define gprotm GPROTM
192 # define gsvolu GSVOLU
193 # define gprint GPRINT
194 # define gdinit GDINIT
197 # define gdrayt GDRAYT
198 # define gdrawc GDRAWC
199 # define gdrawx GDRAWX
200 # define gdhead GDHEAD
201 # define gdwmn1 GDWMN1
202 # define gdwmn2 GDWMN2
203 # define gdwmn3 GDWMN3
205 # define gdcxyz GDCXYZ
207 # define gdfspc GDFSPC
208 # define gdspec GDSPEC
209 # define gdtree GDTREE
210 # define gdelet GDELET
211 # define gdclos GDCLOS
212 # define gdshow GDSHOW
213 # define gdopen GDOPEN
214 # define dzshow DZSHOW
216 # define gfpara GFPARA
217 # define gckpar GCKPAR
218 # define gckmat GCKMAT
219 # define geditv GEDITV
220 # define mzdrop MZDROP
222 # define ertrak ERTRAK
223 # define ertrgo ERTRGO
225 # define setbomb SETBOMB
226 # define setclip SETCLIP
227 # define gcomad GCOMAD
231 //____________________________________________________________________________
235 // Prototypes for GEANT functions
237 void type_of_call gzebra(const int&);
239 void type_of_call gpcxyz();
241 void type_of_call ggclos();
243 void type_of_call glast();
245 void type_of_call ginit();
247 void type_of_call gcinit();
249 void type_of_call grun();
251 void type_of_call gtrig();
253 void type_of_call gtrigc();
255 void type_of_call gtrigi();
257 void type_of_call gwork(const int&);
259 void type_of_call gzinit();
261 void type_of_call gmate();
263 void type_of_call gpart();
265 void type_of_call gsdk(Int_t &, Float_t *, Int_t *);
267 void type_of_call gfkine(Int_t &, Float_t *, Float_t *, Int_t &,
268 Int_t &, Float_t *, Int_t &);
270 void type_of_call gfvert(Int_t &, Float_t *, Int_t &, Int_t &,
271 Float_t &, Float_t *, Int_t &);
273 void type_of_call gskine(Float_t *,Int_t &, Int_t &, Float_t *,
276 void type_of_call gsvert(Float_t *,Int_t &, Int_t &, Float_t *,
279 void type_of_call gphysi();
281 void type_of_call gdebug();
283 void type_of_call gekbin();
285 void type_of_call gfinds();
287 void type_of_call gsking(Int_t &);
289 void type_of_call gskpho(Int_t &);
291 void type_of_call gsstak(Int_t &);
293 void type_of_call gsxyz();
295 void type_of_call gtrack();
297 void type_of_call gtreve();
299 void type_of_call gtreve_root();
301 void type_of_call grndm(Float_t *, const Int_t &);
303 void type_of_call grndmq(Int_t &, Int_t &, const Int_t &,
306 void type_of_call gdtom(Float_t *, Float_t *, Int_t &);
308 void type_of_call glmoth(DEFCHARD, Int_t &, Int_t &, Int_t *,
309 Int_t *, Int_t * DEFCHARL);
311 void type_of_call gmedia(Float_t *, Int_t &);
313 void type_of_call gmtod(Float_t *, Float_t *, Int_t &);
315 void type_of_call gsrotm(const Int_t &, const Float_t &, const Float_t &,
316 const Float_t &, const Float_t &, const Float_t &,
319 void type_of_call gprotm(const Int_t &);
321 void type_of_call grfile(const Int_t&, DEFCHARD,
322 DEFCHARD DEFCHARL DEFCHARL);
324 void type_of_call gfmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
325 Float_t &, Float_t &, Float_t &, Float_t *,
328 void type_of_call gfpart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
329 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
331 void type_of_call gftmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
332 Float_t &, Float_t &, Float_t &, Float_t &,
333 Float_t &, Float_t &, Float_t *, Int_t * DEFCHARL);
335 void type_of_call gsmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
336 Float_t &, Float_t &, Float_t &, Float_t *,
339 void type_of_call gsmixt(const Int_t&, DEFCHARD, Float_t *, Float_t *,
340 Float_t &, Int_t &, Float_t * DEFCHARL);
342 void type_of_call gspart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
343 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
346 void type_of_call gstmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
347 Float_t &, Float_t &, Float_t &, Float_t &,
348 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
350 void type_of_call gsckov(Int_t &itmed, Int_t &npckov, Float_t *ppckov,
351 Float_t *absco, Float_t *effic, Float_t *rindex);
352 void type_of_call gstpar(const Int_t&, DEFCHARD, Float_t & DEFCHARL);
354 void type_of_call gsdvn(DEFCHARD,DEFCHARD, Int_t &, Int_t &
357 void type_of_call gsdvn2(DEFCHARD,DEFCHARD, Int_t &, Int_t &, Float_t &,
358 Int_t & DEFCHARL DEFCHARL);
360 void type_of_call gsdvs(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &
363 void type_of_call gsdvs2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t &,
364 Int_t & DEFCHARL DEFCHARL);
366 void type_of_call gsdvt(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &,
367 Int_t & DEFCHARL DEFCHARL);
369 void type_of_call gsdvt2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t&,
370 Int_t &, Int_t & DEFCHARL DEFCHARL);
372 void type_of_call gsord(DEFCHARD, Int_t & DEFCHARL);
374 void type_of_call gspos(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
375 Float_t &, Int_t &, DEFCHARD DEFCHARL DEFCHARL
378 void type_of_call gsposp(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
379 Float_t &, Int_t &, DEFCHARD,
380 Float_t *, Int_t & DEFCHARL DEFCHARL DEFCHARL);
382 void type_of_call gsvolu(DEFCHARD, DEFCHARD, Int_t &, Float_t *, Int_t &,
383 Int_t & DEFCHARL DEFCHARL);
385 void type_of_call gsatt(DEFCHARD, DEFCHARD, Int_t & DEFCHARL DEFCHARL);
387 void type_of_call gfpara(DEFCHARD , Int_t&, Int_t&, Int_t&, Int_t&, Float_t*,
390 void type_of_call gckpar(Int_t&, Int_t&, Float_t*);
392 void type_of_call gckmat(Int_t&, DEFCHARD DEFCHARL);
394 void type_of_call gprint(DEFCHARD,const int& DEFCHARL);
396 void type_of_call gdinit();
398 void type_of_call gdopt(DEFCHARD,DEFCHARD DEFCHARL DEFCHARL);
400 void type_of_call gdraw(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
401 Float_t &, Float_t &, Float_t & DEFCHARL);
402 void type_of_call gdrayt(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
403 Float_t &, Float_t &, Float_t & DEFCHARL);
404 void type_of_call gdrawc(DEFCHARD,Int_t &, Float_t &, Float_t &, Float_t &,
405 Float_t &, Float_t & DEFCHARL);
406 void type_of_call gdrawx(DEFCHARD,Float_t &, Float_t &, Float_t &, Float_t &,
407 Float_t &, Float_t &, Float_t &, Float_t &,
409 void type_of_call gdhead(Int_t &,DEFCHARD, Float_t & DEFCHARL);
410 void type_of_call gdxyz(Int_t &);
411 void type_of_call gdcxyz();
412 void type_of_call gdman(Float_t &, Float_t &);
413 void type_of_call gdwmn1(Float_t &, Float_t &);
414 void type_of_call gdwmn2(Float_t &, Float_t &);
415 void type_of_call gdwmn3(Float_t &, Float_t &);
416 void type_of_call gdspec(DEFCHARD DEFCHARL);
417 void type_of_call gdfspc(DEFCHARD, Int_t &, Int_t & DEFCHARL) {;}
418 void type_of_call gdtree(DEFCHARD, Int_t &, Int_t & DEFCHARL);
420 void type_of_call gdopen(Int_t &);
421 void type_of_call gdclos();
422 void type_of_call gdelet(Int_t &);
423 void type_of_call gdshow(Int_t &);
424 void type_of_call geditv(Int_t &) {;}
427 void type_of_call dzshow(DEFCHARD,const int&,const int&,DEFCHARD,const int&,
428 const int&, const int&, const int& DEFCHARL
431 void type_of_call mzdrop(Int_t&, Int_t&, DEFCHARD DEFCHARL);
433 void type_of_call setbomb(Float_t &);
434 void type_of_call setclip(DEFCHARD, Float_t &,Float_t &,Float_t &,Float_t &,
435 Float_t &, Float_t & DEFCHARL);
436 void type_of_call gcomad(DEFCHARD, Int_t*& DEFCHARL);
438 void type_of_call ertrak(const Float_t *const x1, const Float_t *const p1,
439 const Float_t *x2, const Float_t *p2,
440 const Int_t &ipa, DEFCHARD DEFCHARL);
442 void type_of_call ertrgo();
446 // Geant3 global pointer
448 static Int_t defSize = 600;
452 //____________________________________________________________________________
456 // Default constructor
460 //____________________________________________________________________________
461 TGeant3::TGeant3(const char *title, Int_t nwgeant)
462 :AliMC("TGeant3",title)
465 // Standard constructor for TGeant3 with ZEBRA initialisation
476 // Load Address of Geant3 commons
479 // Zero number of particles
483 //____________________________________________________________________________
484 Int_t TGeant3::CurrentMaterial(Float_t &a, Float_t &z, Float_t &dens,
485 Float_t &radl, Float_t &absl) const
488 // Return the parameters of the current material during transport
492 dens = fGcmate->dens;
493 radl = fGcmate->radl;
494 absl = fGcmate->absl;
495 return 1; //this could be the number of elements in mixture
498 //____________________________________________________________________________
499 void TGeant3::DefaultRange()
502 // Set range of current drawing pad to 20x20 cm
508 higz->Range(0,0,20,20);
511 //____________________________________________________________________________
512 void TGeant3::InitHIGZ()
523 //____________________________________________________________________________
524 void TGeant3::LoadAddress()
527 // Assigns the address of the GEANT common blocks to the structures
528 // that allow their access from C++
531 gcomad(PASSCHARD("QUEST"), (int*&) fQuest PASSCHARL("QUEST"));
532 gcomad(PASSCHARD("GCBANK"),(int*&) fGcbank PASSCHARL("GCBANK"));
533 gcomad(PASSCHARD("GCLINK"),(int*&) fGclink PASSCHARL("GCLINK"));
534 gcomad(PASSCHARD("GCCUTS"),(int*&) fGccuts PASSCHARL("GCCUTS"));
535 gcomad(PASSCHARD("GCFLAG"),(int*&) fGcflag PASSCHARL("GCFLAG"));
536 gcomad(PASSCHARD("GCKINE"),(int*&) fGckine PASSCHARL("GCKINE"));
537 gcomad(PASSCHARD("GCKING"),(int*&) fGcking PASSCHARL("GCKING"));
538 gcomad(PASSCHARD("GCKIN2"),(int*&) fGckin2 PASSCHARL("GCKIN2"));
539 gcomad(PASSCHARD("GCKIN3"),(int*&) fGckin3 PASSCHARL("GCKIN3"));
540 gcomad(PASSCHARD("GCMATE"),(int*&) fGcmate PASSCHARL("GCMATE"));
541 gcomad(PASSCHARD("GCTMED"),(int*&) fGctmed PASSCHARL("GCTMED"));
542 gcomad(PASSCHARD("GCTRAK"),(int*&) fGctrak PASSCHARL("GCTRAK"));
543 gcomad(PASSCHARD("GCTPOL"),(int*&) fGctpol PASSCHARL("GCTPOL"));
544 gcomad(PASSCHARD("GCVOLU"),(int*&) fGcvolu PASSCHARL("GCVOLU"));
545 gcomad(PASSCHARD("GCNUM"), (int*&) fGcnum PASSCHARL("GCNUM"));
546 gcomad(PASSCHARD("GCSETS"),(int*&) fGcsets PASSCHARL("GCSETS"));
547 gcomad(PASSCHARD("GCPHYS"),(int*&) fGcphys PASSCHARL("GCPHYS"));
548 gcomad(PASSCHARD("GCOPTI"),(int*&) fGcopti PASSCHARL("GCOPTI"));
549 gcomad(PASSCHARD("GCTLIT"),(int*&) fGctlit PASSCHARL("GCTLIT"));
550 gcomad(PASSCHARD("GCVDMA"),(int*&) fGcvdma PASSCHARL("GCVDMA"));
553 gcomad(PASSCHARD("ERTRIO"),(int*&) fErtrio PASSCHARL("ERTRIO"));
554 gcomad(PASSCHARD("EROPTS"),(int*&) fEropts PASSCHARL("EROPTS"));
555 gcomad(PASSCHARD("EROPTC"),(int*&) fEroptc PASSCHARL("EROPTC"));
556 gcomad(PASSCHARD("ERWORK"),(int*&) fErwork PASSCHARL("ERWORK"));
558 // Variables for ZEBRA store
559 gcomad(PASSCHARD("IQ"), addr PASSCHARL("IQ"));
561 gcomad(PASSCHARD("LQ"), addr PASSCHARL("LQ"));
566 //_____________________________________________________________________________
567 void TGeant3::GeomIter()
570 // Geometry iterator for moving upward in the geometry tree
571 // Initialise the iterator
573 fNextVol=fGcvolu->nlevel;
576 //____________________________________________________________________________
577 Int_t TGeant3::NextVolUp(Text_t *name, Int_t ©)
580 // Geometry iterator for moving upward in the geometry tree
581 // Return next volume up
586 gname=fGcvolu->names[fNextVol];
587 strncpy(name,(char *) &gname, 4);
589 copy=fGcvolu->number[fNextVol];
590 i=fGcvolu->lvolum[fNextVol];
591 if(gname == fZiq[fGclink->jvolum+i]) return i;
592 else printf("GeomTree: Volume %s not found in bank\n",name);
597 //_____________________________________________________________________________
598 Int_t TGeant3::CurrentVolID(Int_t ©) const
601 // Returns the current volume ID and copy number
604 if( (i=fGcvolu->nlevel-1) < 0 ) {
605 Warning("CurrentVolID","Stack depth only %d\n",fGcvolu->nlevel);
607 gname=fGcvolu->names[i];
608 copy=fGcvolu->number[i];
609 i=fGcvolu->lvolum[i];
610 if(gname == fZiq[fGclink->jvolum+i]) return i;
611 else Warning("CurrentVolID","Volume %4s not found\n",(char*)&gname);
616 //_____________________________________________________________________________
617 Int_t TGeant3::CurrentVolOffID(Int_t off, Int_t ©) const
620 // Return the current volume "off" upward in the geometrical tree
621 // ID and copy number
624 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
625 Warning("CurrentVolOffID","Offset requested %d but stack depth %d\n",
626 off,fGcvolu->nlevel);
628 gname=fGcvolu->names[i];
629 copy=fGcvolu->number[i];
630 i=fGcvolu->lvolum[i];
631 if(gname == fZiq[fGclink->jvolum+i]) return i;
632 else Warning("CurrentVolOffID","Volume %4s not found\n",(char*)&gname);
637 //_____________________________________________________________________________
638 const char* TGeant3::CurrentVolName() const
641 // Returns the current volume name
645 if( (i=fGcvolu->nlevel-1) < 0 ) {
646 Warning("CurrentVolName","Stack depth %d\n",fGcvolu->nlevel);
648 gname=fGcvolu->names[i];
650 strncpy(name,(char *) &gname, 4);
652 i=fGcvolu->lvolum[i];
653 if(gname == fZiq[fGclink->jvolum+i]) return name;
654 else Warning("CurrentVolName","Volume %4s not found\n",name);
659 //_____________________________________________________________________________
660 const char* TGeant3::CurrentVolOffName(Int_t off) const
663 // Return the current volume "off" upward in the geometrical tree
664 // ID, name and copy number
665 // if name=0 no name is returned
669 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
670 Warning("CurrentVolOffName",
671 "Offset requested %d but stack depth %d\n",off,fGcvolu->nlevel);
673 gname=fGcvolu->names[i];
675 strncpy(name,(char *) &gname, 4);
677 i=fGcvolu->lvolum[i];
678 if(gname == fZiq[fGclink->jvolum+i]) return name;
679 else Warning("CurrentVolOffName","Volume %4s not found\n",name);
684 //_____________________________________________________________________________
685 Int_t TGeant3::IdFromPDG(Int_t pdg) const
688 // Return Geant3 code from PDG and pseudo ENDF code
690 for(Int_t i=0;i<fNPDGCodes;++i)
691 if(pdg==fPDGCode[i]) return i;
695 //_____________________________________________________________________________
696 Int_t TGeant3::PDGFromId(Int_t id) const
698 if(id>0 && id<fNPDGCodes) return fPDGCode[id];
702 //_____________________________________________________________________________
703 void TGeant3::DefineParticles()
706 // Define standard Geant 3 particles
709 // Load standard numbers for GEANT particles and PDG conversion
710 fPDGCode[fNPDGCodes++]=-99; // 0 = unused location
711 fPDGCode[fNPDGCodes++]=22; // 1 = photon
712 fPDGCode[fNPDGCodes++]=-11; // 2 = positron
713 fPDGCode[fNPDGCodes++]=11; // 3 = electron
714 fPDGCode[fNPDGCodes++]=12; // 4 = neutrino e
715 fPDGCode[fNPDGCodes++]=-13; // 5 = muon +
716 fPDGCode[fNPDGCodes++]=13; // 6 = muon -
717 fPDGCode[fNPDGCodes++]=111; // 7 = pi0
718 fPDGCode[fNPDGCodes++]=211; // 8 = pi+
719 fPDGCode[fNPDGCodes++]=-211; // 9 = pi-
720 fPDGCode[fNPDGCodes++]=130; // 10 = Kaon Long
721 fPDGCode[fNPDGCodes++]=321; // 11 = Kaon +
722 fPDGCode[fNPDGCodes++]=-321; // 12 = Kaon -
723 fPDGCode[fNPDGCodes++]=2112; // 13 = Neutron
724 fPDGCode[fNPDGCodes++]=2212; // 14 = Proton
725 fPDGCode[fNPDGCodes++]=-2212; // 15 = Anti Proton
726 fPDGCode[fNPDGCodes++]=310; // 16 = Kaon Short
727 fPDGCode[fNPDGCodes++]=221; // 17 = Eta
728 fPDGCode[fNPDGCodes++]=3122; // 18 = Lambda
729 fPDGCode[fNPDGCodes++]=3222; // 19 = Sigma +
730 fPDGCode[fNPDGCodes++]=3212; // 20 = Sigma 0
731 fPDGCode[fNPDGCodes++]=3112; // 21 = Sigma -
732 fPDGCode[fNPDGCodes++]=3322; // 22 = Xi0
733 fPDGCode[fNPDGCodes++]=3312; // 23 = Xi-
734 fPDGCode[fNPDGCodes++]=3334; // 24 = Omega-
735 fPDGCode[fNPDGCodes++]=-2112; // 25 = Anti Proton
736 fPDGCode[fNPDGCodes++]=-3122; // 26 = Anti Proton
737 fPDGCode[fNPDGCodes++]=-3222; // 27 = Anti Sigma -
738 fPDGCode[fNPDGCodes++]=-3212; // 28 = Anti Sigma 0
739 fPDGCode[fNPDGCodes++]=-3112; // 29 = Anti Sigma 0
740 fPDGCode[fNPDGCodes++]=-3322; // 30 = Anti Xi 0
741 fPDGCode[fNPDGCodes++]=-3312; // 31 = Anti Xi +
742 fPDGCode[fNPDGCodes++]=-3334; // 32 = Anti Omega +
749 /* --- Define additional particles */
750 Gspart(33, "OMEGA(782)", 3, 0.782, 0., 7.836e-23);
751 fPDGCode[fNPDGCodes++]=223; // 33 = Omega(782)
753 Gspart(34, "PHI(1020)", 3, 1.019, 0., 1.486e-22);
754 fPDGCode[fNPDGCodes++]=333; // 34 = PHI (1020)
756 Gspart(35, "D +", 4, 1.87, 1., 1.066e-12);
757 fPDGCode[fNPDGCodes++]=411; // 35 = D+
759 Gspart(36, "D -", 4, 1.87, -1., 1.066e-12);
760 fPDGCode[fNPDGCodes++]=-411; // 36 = D-
762 Gspart(37, "D 0", 3, 1.865, 0., 4.2e-13);
763 fPDGCode[fNPDGCodes++]=421; // 37 = D0
765 Gspart(38, "ANTI D 0", 3, 1.865, 0., 4.2e-13);
766 fPDGCode[fNPDGCodes++]=-421; // 38 = D0 bar
768 fPDGCode[fNPDGCodes++]=-99; // 39 = unassigned
770 fPDGCode[fNPDGCodes++]=-99; // 40 = unassigned
772 fPDGCode[fNPDGCodes++]=-99; // 41 = unassigned
774 Gspart(42, "RHO +", 4, 0.768, 1., 4.353e-24);
775 fPDGCode[fNPDGCodes++]=213; // 42 = RHO+
777 Gspart(43, "RHO -", 4, 0.768, -1., 4.353e-24);
778 fPDGCode[fNPDGCodes++]=-213; // 40 = RHO-
780 Gspart(44, "RHO 0", 3, 0.768, 0., 4.353e-24);
781 fPDGCode[fNPDGCodes++]=113; // 37 = D0
784 // Use ENDF-6 mapping for ions = 10000*z+10*a+iso
786 // and numbers above 5 000 000 for special applications
789 const Int_t kion=10000000;
791 const Int_t kspe=50000000;
793 TDatabasePDG *pdgDB = TDatabasePDG::Instance();
795 const Double_t autogev=0.9314943228;
796 const Double_t hslash = 1.0545726663e-27;
797 const Double_t erggev = 1/1.6021773349e-3;
798 const Double_t hshgev = hslash*erggev;
799 const Double_t yearstosec = 3600*24*365.25;
802 pdgDB->AddParticle("Deuteron","Deuteron",2*autogev+8.071e-3,kTRUE,
803 0,1,"Ion",kion+10020);
804 fPDGCode[fNPDGCodes++]=kion+10020; // 45 = Deuteron
806 pdgDB->AddParticle("Triton","Triton",3*autogev+14.931e-3,kFALSE,
807 hshgev/(12.33*yearstosec),1,"Ion",kion+10030);
808 fPDGCode[fNPDGCodes++]=kion+10030; // 46 = Triton
810 pdgDB->AddParticle("Alpha","Alpha",4*autogev+2.424e-3,kTRUE,
811 hshgev/(12.33*yearstosec),2,"Ion",kion+20040);
812 fPDGCode[fNPDGCodes++]=kion+20040; // 47 = Alpha
814 fPDGCode[fNPDGCodes++]=0; // 48 = geantino mapped to rootino
816 pdgDB->AddParticle("HE3","HE3",3*autogev+14.931e-3,kFALSE,
817 0,2,"Ion",kion+20030);
818 fPDGCode[fNPDGCodes++]=kion+20030; // 49 = HE3
820 pdgDB->AddParticle("Cherenkov","Cherenkov",0,kFALSE,
821 0,0,"Special",kspe+50);
822 fPDGCode[fNPDGCodes++]=kspe+50; // 50 = Cherenkov
824 /* --- Define additional decay modes --- */
825 /* --- omega(783) --- */
826 for (kz = 0; kz < 6; ++kz) {
837 Gsdk(ipa, bratio, mode);
838 /* --- phi(1020) --- */
839 for (kz = 0; kz < 6; ++kz) {
854 Gsdk(ipa, bratio, mode);
856 for (kz = 0; kz < 6; ++kz) {
869 Gsdk(ipa, bratio, mode);
871 for (kz = 0; kz < 6; ++kz) {
884 Gsdk(ipa, bratio, mode);
886 for (kz = 0; kz < 6; ++kz) {
897 Gsdk(ipa, bratio, mode);
898 /* --- Anti D0 --- */
899 for (kz = 0; kz < 6; ++kz) {
910 Gsdk(ipa, bratio, mode);
912 for (kz = 0; kz < 6; ++kz) {
919 Gsdk(ipa, bratio, mode);
921 for (kz = 0; kz < 6; ++kz) {
928 Gsdk(ipa, bratio, mode);
930 for (kz = 0; kz < 6; ++kz) {
937 Gsdk(ipa, bratio, mode);
940 for (kz = 0; kz < 6; ++kz) {
949 Gsdk(ipa, bratio, mode);
952 Gsdk(ipa, bratio, mode);
955 Gsdk(ipa, bratio, mode);
960 //_____________________________________________________________________________
961 Int_t TGeant3::VolId(Text_t *name) const
964 // Return the unique numeric identifier for volume name
967 strncpy((char *) &gname, name, 4);
968 for(i=1; i<=fGcnum->nvolum; i++)
969 if(gname == fZiq[fGclink->jvolum+i]) return i;
970 printf("VolId: Volume %s not found\n",name);
974 //_____________________________________________________________________________
975 Int_t TGeant3::NofVolumes() const
978 // Return total number of volumes in the geometry
980 return fGcnum->nvolum;
983 //_____________________________________________________________________________
984 const char* TGeant3::VolName(Int_t id) const
987 // Return the volume name given the volume identifier
990 if(id<1 || id > fGcnum->nvolum || fGclink->jvolum<=0)
993 strncpy(name,(char *)&fZiq[fGclink->jvolum+id],4);
998 //_____________________________________________________________________________
999 Float_t TGeant3::Xsec(char* reac, Float_t energy, Int_t part, Int_t mate)
1001 Int_t gpart = IdFromPDG(part);
1002 if(!strcmp(reac,"PHOT"))
1005 Error("Xsec","Can calculate photoelectric only for photons\n");
1011 //_____________________________________________________________________________
1012 void TGeant3::TrackPosition(TLorentzVector &xyz) const
1015 // Return the current position in the master reference frame of the
1016 // track being transported
1018 xyz[0]=fGctrak->vect[0];
1019 xyz[1]=fGctrak->vect[1];
1020 xyz[2]=fGctrak->vect[2];
1021 xyz[3]=fGctrak->tofg;
1024 //_____________________________________________________________________________
1025 Float_t TGeant3::TrackTime() const
1028 // Return the current time of flight of the track being transported
1030 return fGctrak->tofg;
1033 //_____________________________________________________________________________
1034 void TGeant3::TrackMomentum(TLorentzVector &xyz) const
1037 // Return the direction and the momentum (GeV/c) of the track
1038 // currently being transported
1040 Double_t ptot=fGctrak->vect[6];
1041 xyz[0]=fGctrak->vect[3]*ptot;
1042 xyz[1]=fGctrak->vect[4]*ptot;
1043 xyz[2]=fGctrak->vect[5]*ptot;
1044 xyz[3]=fGctrak->getot;
1047 //_____________________________________________________________________________
1048 Float_t TGeant3::TrackCharge() const
1051 // Return charge of the track currently transported
1053 return fGckine->charge;
1056 //_____________________________________________________________________________
1057 Float_t TGeant3::TrackMass() const
1060 // Return the mass of the track currently transported
1062 return fGckine->amass;
1065 //_____________________________________________________________________________
1066 Int_t TGeant3::TrackPid() const
1069 // Return the id of the particle transported
1071 return PDGFromId(fGckine->ipart);
1074 //_____________________________________________________________________________
1075 Float_t TGeant3::TrackStep() const
1078 // Return the length in centimeters of the current step
1080 return fGctrak->step;
1083 //_____________________________________________________________________________
1084 Float_t TGeant3::TrackLength() const
1087 // Return the length of the current track from its origin
1089 return fGctrak->sleng;
1092 //_____________________________________________________________________________
1093 Bool_t TGeant3::IsTrackInside() const
1096 // True if the track is not at the boundary of the current volume
1098 return (fGctrak->inwvol==0);
1101 //_____________________________________________________________________________
1102 Bool_t TGeant3::IsTrackEntering() const
1105 // True if this is the first step of the track in the current volume
1107 return (fGctrak->inwvol==1);
1110 //_____________________________________________________________________________
1111 Bool_t TGeant3::IsTrackExiting() const
1114 // True if this is the last step of the track in the current volume
1116 return (fGctrak->inwvol==2);
1119 //_____________________________________________________________________________
1120 Bool_t TGeant3::IsTrackOut() const
1123 // True if the track is out of the setup
1125 return (fGctrak->inwvol==3);
1128 //_____________________________________________________________________________
1129 Bool_t TGeant3::IsTrackStop() const
1132 // True if the track energy has fallen below the threshold
1134 return (fGctrak->istop==2);
1137 //_____________________________________________________________________________
1138 Int_t TGeant3::NSecondaries() const
1141 // Number of secondary particles generated in the current step
1143 return fGcking->ngkine;
1146 //_____________________________________________________________________________
1147 Int_t TGeant3::CurrentEvent() const
1150 // Number of the current event
1152 return fGcflag->idevt;
1155 //_____________________________________________________________________________
1156 void TGeant3::ProdProcess(char* proc) const
1159 // Name of the process that has produced the secondary particles
1160 // in the current step
1162 const Int_t ipmec[13] = { 5,6,7,8,9,10,11,12,21,23,25,105,108 };
1165 if(fGcking->ngkine>0) {
1166 for (km = 0; km < fGctrak->nmec; ++km) {
1167 for (im = 0; im < 13; ++im) {
1168 if (fGctrak->lmec[km] == ipmec[im]) {
1169 mec = fGctrak->lmec[km];
1170 if (0 < mec && mec < 31) {
1171 strncpy(proc,(char *)&fGctrak->namec[mec - 1],4);
1172 } else if (mec - 100 <= 30 && mec - 100 > 0) {
1173 strncpy(proc,(char *)&fGctpol->namec1[mec - 101],4);
1180 strcpy(proc,"UNKN");
1181 } else strcpy(proc,"NONE");
1184 //_____________________________________________________________________________
1185 void TGeant3::GetSecondary(Int_t isec, Int_t& ipart, Float_t* x, Float_t* p)
1188 // Get the parameters of the secondary track number isec produced
1189 // in the current step
1192 if(-1<isec && isec<fGcking->ngkine) {
1193 ipart=Int_t (fGcking->gkin[isec][4] +0.5);
1195 x[i]=fGckin3->gpos[isec][i];
1196 p[i]=fGcking->gkin[isec][i];
1198 x[3]=fGcking->tofd[isec];
1199 p[3]=fGcking->gkin[isec][3];
1201 printf(" * TGeant3::GetSecondary * Secondary %d does not exist\n",isec);
1202 x[0]=x[1]=x[2]=x[3]=p[0]=p[1]=p[2]=p[3]=0;
1207 //_____________________________________________________________________________
1208 void TGeant3::InitLego()
1211 SetDEBU(0,0,0); //do not print a message
1214 //_____________________________________________________________________________
1215 Bool_t TGeant3::IsTrackDisappeared() const
1218 // True if the current particle has disappered
1219 // either because it decayed or because it underwent
1220 // an inelastic collision
1222 return (fGctrak->istop==1);
1225 //_____________________________________________________________________________
1226 Bool_t TGeant3::IsTrackAlive() const
1229 // True if the current particle is alive and will continue to be
1232 return (fGctrak->istop==0);
1235 //_____________________________________________________________________________
1236 void TGeant3::StopTrack()
1239 // Stop the transport of the current particle and skip to the next
1244 //_____________________________________________________________________________
1245 void TGeant3::StopEvent()
1248 // Stop simulation of the current event and skip to the next
1253 //_____________________________________________________________________________
1254 Float_t TGeant3::MaxStep() const
1257 // Return the maximum step length in the current medium
1259 return fGctmed->stemax;
1262 //_____________________________________________________________________________
1263 void TGeant3::SetColors()
1266 // Set the colors for all the volumes
1267 // this is done sequentially for all volumes
1268 // based on the number of their medium
1271 Int_t jvolum=fGclink->jvolum;
1272 //Int_t jtmed=fGclink->jtmed;
1273 //Int_t jmate=fGclink->jmate;
1274 Int_t nvolum=fGcnum->nvolum;
1277 // Now for all the volumes
1278 for(kv=1;kv<=nvolum;kv++) {
1279 // Get the tracking medium
1280 Int_t itm=Int_t (fZq[fZlq[jvolum-kv]+4]);
1282 //Int_t ima=Int_t (fZq[fZlq[jtmed-itm]+6]);
1284 //Float_t z=fZq[fZlq[jmate-ima]+7];
1285 // Find color number
1286 //icol = Int_t(z)%6+2;
1287 //icol = 17+Int_t(z*150./92.);
1290 strncpy(name,(char*)&fZiq[jvolum+kv],4);
1292 Gsatt(name,"COLO",icol);
1296 //_____________________________________________________________________________
1297 void TGeant3::SetMaxStep(Float_t maxstep)
1300 // Set the maximum step allowed till the particle is in the current medium
1302 fGctmed->stemax=maxstep;
1305 //_____________________________________________________________________________
1306 void TGeant3::SetMaxNStep(Int_t maxnstp)
1309 // Set the maximum number of steps till the particle is in the current medium
1311 fGctrak->maxnst=maxnstp;
1314 //_____________________________________________________________________________
1315 Int_t TGeant3::GetMaxNStep() const
1318 // Maximum number of steps allowed in current medium
1320 return fGctrak->maxnst;
1323 //_____________________________________________________________________________
1324 void TGeant3::Material(Int_t& kmat, const char* name, Float_t a, Float_t z,
1325 Float_t dens, Float_t radl, Float_t absl, Float_t* buf,
1329 // Defines a Material
1331 // kmat number assigned to the material
1332 // name material name
1333 // a atomic mass in au
1335 // dens density in g/cm3
1336 // absl absorbtion length in cm
1337 // if >=0 it is ignored and the program
1338 // calculates it, if <0. -absl is taken
1339 // radl radiation length in cm
1340 // if >=0 it is ignored and the program
1341 // calculates it, if <0. -radl is taken
1342 // buf pointer to an array of user words
1343 // nbuf number of user words
1345 Int_t jmate=fGclink->jmate;
1351 for(i=1; i<=ns; i++) {
1352 if(fZlq[jmate-i]==0) {
1358 gsmate(kmat,PASSCHARD(name), a, z, dens, radl, absl, buf,
1359 nwbuf PASSCHARL(name));
1362 //_____________________________________________________________________________
1363 void TGeant3::Mixture(Int_t& kmat, const char* name, Float_t* a, Float_t* z,
1364 Float_t dens, Int_t nlmat, Float_t* wmat)
1367 // Defines mixture OR COMPOUND IMAT as composed by
1368 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
1370 // If NLMAT > 0 then wmat contains the proportion by
1371 // weights of each basic material in the mixture.
1373 // If nlmat < 0 then WMAT contains the number of atoms
1374 // of a given kind into the molecule of the COMPOUND
1375 // In this case, WMAT in output is changed to relative
1378 Int_t jmate=fGclink->jmate;
1384 for(i=1; i<=ns; i++) {
1385 if(fZlq[jmate-i]==0) {
1391 gsmixt(kmat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
1394 //_____________________________________________________________________________
1395 void TGeant3::Medium(Int_t& kmed, const char* name, Int_t nmat, Int_t isvol,
1396 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
1397 Float_t stemax, Float_t deemax, Float_t epsil,
1398 Float_t stmin, Float_t* ubuf, Int_t nbuf)
1401 // kmed tracking medium number assigned
1402 // name tracking medium name
1403 // nmat material number
1404 // isvol sensitive volume flag
1405 // ifield magnetic field
1406 // fieldm max. field value (kilogauss)
1407 // tmaxfd max. angle due to field (deg/step)
1408 // stemax max. step allowed
1409 // deemax max. fraction of energy lost in a step
1410 // epsil tracking precision (cm)
1411 // stmin min. step due to continuos processes (cm)
1413 // ifield = 0 if no magnetic field; ifield = -1 if user decision in guswim;
1414 // ifield = 1 if tracking performed with grkuta; ifield = 2 if tracking
1415 // performed with ghelix; ifield = 3 if tracking performed with ghelx3.
1417 Int_t jtmed=fGclink->jtmed;
1423 for(i=1; i<=ns; i++) {
1424 if(fZlq[jtmed-i]==0) {
1430 gstmed(kmed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1431 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1434 //_____________________________________________________________________________
1435 void TGeant3::Matrix(Int_t& krot, Float_t thex, Float_t phix, Float_t they,
1436 Float_t phiy, Float_t thez, Float_t phiz)
1439 // krot rotation matrix number assigned
1440 // theta1 polar angle for axis i
1441 // phi1 azimuthal angle for axis i
1442 // theta2 polar angle for axis ii
1443 // phi2 azimuthal angle for axis ii
1444 // theta3 polar angle for axis iii
1445 // phi3 azimuthal angle for axis iii
1447 // it defines the rotation matrix number irot.
1449 Int_t jrotm=fGclink->jrotm;
1455 for(i=1; i<=ns; i++) {
1456 if(fZlq[jrotm-i]==0) {
1462 gsrotm(krot, thex, phix, they, phiy, thez, phiz);
1465 //_____________________________________________________________________________
1466 Int_t TGeant3::GetMedium() const
1469 // Return the number of the current medium
1471 return fGctmed->numed;
1474 //_____________________________________________________________________________
1475 Float_t TGeant3::Edep() const
1478 // Return the energy lost in the current step
1480 return fGctrak->destep;
1483 //_____________________________________________________________________________
1484 Float_t TGeant3::Etot() const
1487 // Return the total energy of the current track
1489 return fGctrak->getot;
1492 //_____________________________________________________________________________
1493 void TGeant3::Rndm(Float_t* r, const Int_t n) const
1496 // Return an array of n random numbers uniformly distributed
1497 // between 0 and 1 not included
1502 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1504 // Functions from GBASE
1506 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1508 //____________________________________________________________________________
1509 void TGeant3::Gfile(const char *filename, const char *option)
1512 // Routine to open a GEANT/RZ data base.
1514 // LUN logical unit number associated to the file
1516 // CHFILE RZ file name
1518 // CHOPT is a character string which may be
1519 // N To create a new file
1520 // U to open an existing file for update
1521 // " " to open an existing file for read only
1522 // Q The initial allocation (default 1000 records)
1523 // is given in IQUEST(10)
1524 // X Open the file in exchange format
1525 // I Read all data structures from file to memory
1526 // O Write all data structures from memory to file
1529 // If options "I" or "O" all data structures are read or
1530 // written from/to file and the file is closed.
1531 // See routine GRMDIR to create subdirectories
1532 // See routines GROUT,GRIN to write,read objects
1534 grfile(21, PASSCHARD(filename), PASSCHARD(option) PASSCHARL(filename)
1538 //____________________________________________________________________________
1539 void TGeant3::Gpcxyz()
1542 // Print track and volume parameters at current point
1547 //_____________________________________________________________________________
1548 void TGeant3::Ggclos()
1551 // Closes off the geometry setting.
1552 // Initializes the search list for the contents of each
1553 // volume following the order they have been positioned, and
1554 // inserting the content '0' when a call to GSNEXT (-1) has
1555 // been required by the user.
1556 // Performs the development of the JVOLUM structure for all
1557 // volumes with variable parameters, by calling GGDVLP.
1558 // Interprets the user calls to GSORD, through GGORD.
1559 // Computes and stores in a bank (next to JVOLUM mother bank)
1560 // the number of levels in the geometrical tree and the
1561 // maximum number of contents per level, by calling GGNLEV.
1562 // Sets status bit for CONCAVE volumes, through GGCAVE.
1563 // Completes the JSET structure with the list of volume names
1564 // which identify uniquely a given physical detector, the
1565 // list of bit numbers to pack the corresponding volume copy
1566 // numbers, and the generic path(s) in the JVOLUM tree,
1567 // through the routine GHCLOS.
1572 //_____________________________________________________________________________
1573 void TGeant3::Glast()
1576 // Finish a Geant run
1581 //_____________________________________________________________________________
1582 void TGeant3::Gprint(const char *name)
1585 // Routine to print data structures
1586 // CHNAME name of a data structure
1590 gprint(PASSCHARD(vname),0 PASSCHARL(vname));
1593 //_____________________________________________________________________________
1594 void TGeant3::Grun()
1597 // Steering function to process one run
1602 //_____________________________________________________________________________
1603 void TGeant3::Gtrig()
1606 // Steering function to process one event
1611 //_____________________________________________________________________________
1612 void TGeant3::Gtrigc()
1615 // Clear event partition
1620 //_____________________________________________________________________________
1621 void TGeant3::Gtrigi()
1624 // Initialises event partition
1629 //_____________________________________________________________________________
1630 void TGeant3::Gwork(Int_t nwork)
1633 // Allocates workspace in ZEBRA memory
1638 //_____________________________________________________________________________
1639 void TGeant3::Gzinit()
1642 // To initialise GEANT/ZEBRA data structures
1647 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1649 // Functions from GCONS
1651 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1653 //_____________________________________________________________________________
1654 void TGeant3::Gfmate(Int_t imat, char *name, Float_t &a, Float_t &z,
1655 Float_t &dens, Float_t &radl, Float_t &absl,
1656 Float_t* ubuf, Int_t& nbuf)
1659 // Return parameters for material IMAT
1661 gfmate(imat, PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
1665 //_____________________________________________________________________________
1666 void TGeant3::Gfpart(Int_t ipart, char *name, Int_t &itrtyp,
1667 Float_t &amass, Float_t &charge, Float_t &tlife)
1670 // Return parameters for particle of type IPART
1674 Int_t igpart = IdFromPDG(ipart);
1675 gfpart(igpart, PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
1679 //_____________________________________________________________________________
1680 void TGeant3::Gftmed(Int_t numed, char *name, Int_t &nmat, Int_t &isvol,
1681 Int_t &ifield, Float_t &fieldm, Float_t &tmaxfd,
1682 Float_t &stemax, Float_t &deemax, Float_t &epsil,
1683 Float_t &stmin, Float_t *ubuf, Int_t *nbuf)
1686 // Return parameters for tracking medium NUMED
1688 gftmed(numed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1689 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1692 //_____________________________________________________________________________
1693 void TGeant3::Gmate()
1696 // Define standard GEANT materials
1701 //_____________________________________________________________________________
1702 void TGeant3::Gpart()
1705 // Define standard GEANT particles plus selected decay modes
1706 // and branching ratios.
1711 //_____________________________________________________________________________
1712 void TGeant3::Gsdk(Int_t ipart, Float_t *bratio, Int_t *mode)
1714 // Defines branching ratios and decay modes for standard
1716 gsdk(ipart,bratio,mode);
1719 //_____________________________________________________________________________
1720 void TGeant3::Gsmate(Int_t imat, const char *name, Float_t a, Float_t z,
1721 Float_t dens, Float_t radl, Float_t absl)
1724 // Defines a Material
1726 // kmat number assigned to the material
1727 // name material name
1728 // a atomic mass in au
1730 // dens density in g/cm3
1731 // absl absorbtion length in cm
1732 // if >=0 it is ignored and the program
1733 // calculates it, if <0. -absl is taken
1734 // radl radiation length in cm
1735 // if >=0 it is ignored and the program
1736 // calculates it, if <0. -radl is taken
1737 // buf pointer to an array of user words
1738 // nbuf number of user words
1742 gsmate(imat,PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
1746 //_____________________________________________________________________________
1747 void TGeant3::Gsmixt(Int_t imat, const char *name, Float_t *a, Float_t *z,
1748 Float_t dens, Int_t nlmat, Float_t *wmat)
1751 // Defines mixture OR COMPOUND IMAT as composed by
1752 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
1754 // If NLMAT.GT.0 then WMAT contains the PROPORTION BY
1755 // WEIGTHS OF EACH BASIC MATERIAL IN THE MIXTURE.
1757 // If NLMAT.LT.0 then WMAT contains the number of atoms
1758 // of a given kind into the molecule of the COMPOUND
1759 // In this case, WMAT in output is changed to relative
1762 gsmixt(imat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
1765 //_____________________________________________________________________________
1766 void TGeant3::Gspart(Int_t ipart, const char *name, Int_t itrtyp,
1767 Float_t amass, Float_t charge, Float_t tlife)
1770 // Store particle parameters
1772 // ipart particle code
1773 // name particle name
1774 // itrtyp transport method (see GEANT manual)
1775 // amass mass in GeV/c2
1776 // charge charge in electron units
1777 // tlife lifetime in seconds
1781 gspart(ipart,PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
1785 //_____________________________________________________________________________
1786 void TGeant3::Gstmed(Int_t numed, const char *name, Int_t nmat, Int_t isvol,
1787 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
1788 Float_t stemax, Float_t deemax, Float_t epsil,
1792 // NTMED Tracking medium number
1793 // NAME Tracking medium name
1794 // NMAT Material number
1795 // ISVOL Sensitive volume flag
1796 // IFIELD Magnetic field
1797 // FIELDM Max. field value (Kilogauss)
1798 // TMAXFD Max. angle due to field (deg/step)
1799 // STEMAX Max. step allowed
1800 // DEEMAX Max. fraction of energy lost in a step
1801 // EPSIL Tracking precision (cm)
1802 // STMIN Min. step due to continuos processes (cm)
1804 // IFIELD = 0 if no magnetic field; IFIELD = -1 if user decision in GUSWIM;
1805 // IFIELD = 1 if tracking performed with GRKUTA; IFIELD = 2 if tracking
1806 // performed with GHELIX; IFIELD = 3 if tracking performed with GHELX3.
1810 gstmed(numed,PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1811 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1814 //_____________________________________________________________________________
1815 void TGeant3::Gsckov(Int_t itmed, Int_t npckov, Float_t *ppckov,
1816 Float_t *absco, Float_t *effic, Float_t *rindex)
1819 // Stores the tables for UV photon tracking in medium ITMED
1820 // Please note that it is the user's responsability to
1821 // provide all the coefficients:
1824 // ITMED Tracking medium number
1825 // NPCKOV Number of bins of each table
1826 // PPCKOV Value of photon momentum (in GeV)
1827 // ABSCO Absorbtion coefficients
1828 // dielectric: absorbtion length in cm
1829 // metals : absorbtion fraction (0<=x<=1)
1830 // EFFIC Detection efficiency for UV photons
1831 // RINDEX Refraction index (if=0 metal)
1833 gsckov(itmed,npckov,ppckov,absco,effic,rindex);
1836 //_____________________________________________________________________________
1837 void TGeant3::Gstpar(Int_t itmed, const char *param, Float_t parval)
1840 // To change the value of cut or mechanism "CHPAR"
1841 // to a new value PARVAL for tracking medium ITMED
1842 // The data structure JTMED contains the standard tracking
1843 // parameters (CUTS and flags to control the physics processes) which
1844 // are used by default for all tracking media. It is possible to
1845 // redefine individually with GSTPAR any of these parameters for a
1846 // given tracking medium.
1847 // ITMED tracking medium number
1848 // CHPAR is a character string (variable name)
1849 // PARVAL must be given as a floating point.
1851 gstpar(itmed,PASSCHARD(param), parval PASSCHARL(param));
1854 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1856 // Functions from GCONS
1858 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1860 //_____________________________________________________________________________
1861 void TGeant3::Gfkine(Int_t itra, Float_t *vert, Float_t *pvert, Int_t &ipart,
1864 // Storing/Retrieving Vertex and Track parameters
1865 // ----------------------------------------------
1867 // Stores vertex parameters.
1868 // VERT array of (x,y,z) position of the vertex
1869 // NTBEAM beam track number origin of the vertex
1870 // =0 if none exists
1871 // NTTARG target track number origin of the vertex
1872 // UBUF user array of NUBUF floating point numbers
1874 // NVTX new vertex number (=0 in case of error).
1875 // Prints vertex parameters.
1876 // IVTX for vertex IVTX.
1877 // (For all vertices if IVTX=0)
1878 // Stores long life track parameters.
1879 // PLAB components of momentum
1880 // IPART type of particle (see GSPART)
1881 // NV vertex number origin of track
1882 // UBUF array of NUBUF floating point user parameters
1884 // NT track number (if=0 error).
1885 // Retrieves long life track parameters.
1886 // ITRA track number for which parameters are requested
1887 // VERT vector origin of the track
1888 // PVERT 4 momentum components at the track origin
1889 // IPART particle type (=0 if track ITRA does not exist)
1890 // NVERT vertex number origin of the track
1891 // UBUF user words stored in GSKINE.
1892 // Prints initial track parameters.
1893 // ITRA for track ITRA
1894 // (For all tracks if ITRA=0)
1898 gfkine(itra,vert,pvert,ipart,nvert,ubuf,nbuf);
1901 //_____________________________________________________________________________
1902 void TGeant3::Gfvert(Int_t nvtx, Float_t *v, Int_t &ntbeam, Int_t &nttarg,
1906 // Retrieves the parameter of a vertex bank
1907 // Vertex is generated from tracks NTBEAM NTTARG
1908 // NVTX is the new vertex number
1912 gfvert(nvtx,v,ntbeam,nttarg,tofg,ubuf,nbuf);
1915 //_____________________________________________________________________________
1916 Int_t TGeant3::Gskine(Float_t *plab, Int_t ipart, Int_t nv, Float_t *buf,
1920 // Store kinematics of track NT into data structure
1921 // Track is coming from vertex NV
1924 gskine(plab, ipart, nv, buf, nwbuf, nt);
1928 //_____________________________________________________________________________
1929 Int_t TGeant3::Gsvert(Float_t *v, Int_t ntbeam, Int_t nttarg, Float_t *ubuf,
1933 // Creates a new vertex bank
1934 // Vertex is generated from tracks NTBEAM NTTARG
1935 // NVTX is the new vertex number
1938 gsvert(v, ntbeam, nttarg, ubuf, nwbuf, nwtx);
1942 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1944 // Functions from GPHYS
1946 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1948 //_____________________________________________________________________________
1949 void TGeant3::Gphysi()
1952 // Initialise material constants for all the physics
1953 // mechanisms used by GEANT
1958 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1960 // Functions from GTRAK
1962 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1964 //_____________________________________________________________________________
1965 void TGeant3::Gdebug()
1968 // Debug the current step
1973 //_____________________________________________________________________________
1974 void TGeant3::Gekbin()
1977 // To find bin number in kinetic energy table
1978 // stored in ELOW(NEKBIN)
1983 //_____________________________________________________________________________
1984 void TGeant3::Gfinds()
1987 // Returns the set/volume parameters corresponding to
1988 // the current space point in /GCTRAK/
1989 // and fill common /GCSETS/
1991 // IHSET user set identifier
1992 // IHDET user detector identifier
1993 // ISET set number in JSET
1994 // IDET detector number in JS=LQ(JSET-ISET)
1995 // IDTYPE detector type (1,2)
1996 // NUMBV detector volume numbers (array of length NVNAME)
1997 // NVNAME number of volume levels
2002 //_____________________________________________________________________________
2003 void TGeant3::Gsking(Int_t igk)
2006 // Stores in stack JSTAK either the IGKth track of /GCKING/,
2007 // or the NGKINE tracks when IGK is 0.
2012 //_____________________________________________________________________________
2013 void TGeant3::Gskpho(Int_t igk)
2016 // Stores in stack JSTAK either the IGKth Cherenkov photon of
2017 // /GCKIN2/, or the NPHOT tracks when IGK is 0.
2022 //_____________________________________________________________________________
2023 void TGeant3::Gsstak(Int_t iflag)
2026 // Stores in auxiliary stack JSTAK the particle currently
2027 // described in common /GCKINE/.
2029 // On request, creates also an entry in structure JKINE :
2031 // 0 : No entry in JKINE structure required (user)
2032 // 1 : New entry in JVERTX / JKINE structures required (user)
2033 // <0 : New entry in JKINE structure at vertex -IFLAG (user)
2034 // 2 : Entry in JKINE structure exists already (from GTREVE)
2039 //_____________________________________________________________________________
2040 void TGeant3::Gsxyz()
2043 // Store space point VECT in banks JXYZ
2048 //_____________________________________________________________________________
2049 void TGeant3::Gtrack()
2052 // Controls tracking of current particle
2057 //_____________________________________________________________________________
2058 void TGeant3::Gtreve()
2061 // Controls tracking of all particles belonging to the current event
2066 //_____________________________________________________________________________
2067 void TGeant3::Gtreve_root()
2070 // Controls tracking of all particles belonging to the current event
2075 //_____________________________________________________________________________
2076 void TGeant3::Grndm(Float_t *rvec, const Int_t len) const
2079 // To generate a vector RVECV of LEN random numbers
2080 // Copy of the CERN Library routine RANECU
2084 //_____________________________________________________________________________
2085 void TGeant3::Grndmq(Int_t &is1, Int_t &is2, const Int_t iseq,
2086 const Text_t *chopt)
2089 // To set/retrieve the seed of the random number generator
2091 grndmq(is1,is2,iseq,PASSCHARD(chopt) PASSCHARL(chopt));
2094 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2096 // Functions from GDRAW
2098 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2100 //_____________________________________________________________________________
2101 void TGeant3::Gdxyz(Int_t it)
2104 // Draw the points stored with Gsxyz relative to track it
2109 //_____________________________________________________________________________
2110 void TGeant3::Gdcxyz()
2113 // Draw the position of the current track
2118 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2120 // Functions from GGEOM
2122 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2124 //_____________________________________________________________________________
2125 void TGeant3::Gdtom(Float_t *xd, Float_t *xm, Int_t iflag)
2128 // Computes coordinates XM (Master Reference System
2129 // knowing the coordinates XD (Detector Ref System)
2130 // The local reference system can be initialized by
2131 // - the tracking routines and GDTOM used in GUSTEP
2132 // - a call to GSCMED(NLEVEL,NAMES,NUMBER)
2133 // (inverse routine is GMTOD)
2135 // If IFLAG=1 convert coordinates
2136 // IFLAG=2 convert direction cosinus
2138 gdtom(xd, xm, iflag);
2141 //_____________________________________________________________________________
2142 void TGeant3::Glmoth(const char* iudet, Int_t iunum, Int_t &nlev, Int_t *lvols,
2146 // Loads the top part of the Volume tree in LVOLS (IVO's),
2147 // LINDX (IN indices) for a given volume defined through
2148 // its name IUDET and number IUNUM.
2150 // The routine stores only upto the last level where JVOLUM
2151 // data structure is developed. If there is no development
2152 // above the current level, it returns NLEV zero.
2154 glmoth(PASSCHARD(iudet), iunum, nlev, lvols, lindx, idum PASSCHARL(iudet));
2157 //_____________________________________________________________________________
2158 void TGeant3::Gmedia(Float_t *x, Int_t &numed)
2161 // Finds in which volume/medium the point X is, and updates the
2162 // common /GCVOLU/ and the structure JGPAR accordingly.
2164 // NUMED returns the tracking medium number, or 0 if point is
2165 // outside the experimental setup.
2170 //_____________________________________________________________________________
2171 void TGeant3::Gmtod(Float_t *xm, Float_t *xd, Int_t iflag)
2174 // Computes coordinates XD (in DRS)
2175 // from known coordinates XM in MRS
2176 // The local reference system can be initialized by
2177 // - the tracking routines and GMTOD used in GUSTEP
2178 // - a call to GMEDIA(XM,NUMED)
2179 // - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER)
2180 // (inverse routine is GDTOM)
2182 // If IFLAG=1 convert coordinates
2183 // IFLAG=2 convert direction cosinus
2185 gmtod(xm, xd, iflag);
2188 //_____________________________________________________________________________
2189 void TGeant3::Gsdvn(const char *name, const char *mother, Int_t ndiv,
2193 // Create a new volume by dividing an existing one
2196 // MOTHER Mother volume name
2197 // NDIV Number of divisions
2200 // X,Y,Z of CAXIS will be translated to 1,2,3 for IAXIS.
2201 // It divides a previously defined volume.
2206 Vname(mother,vmother);
2207 gsdvn(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis PASSCHARL(vname)
2208 PASSCHARL(vmother));
2211 //_____________________________________________________________________________
2212 void TGeant3::Gsdvn2(const char *name, const char *mother, Int_t ndiv,
2213 Int_t iaxis, Float_t c0i, Int_t numed)
2216 // Create a new volume by dividing an existing one
2218 // Divides mother into ndiv divisions called name
2219 // along axis iaxis starting at coordinate value c0.
2220 // the new volume created will be medium number numed.
2225 Vname(mother,vmother);
2226 gsdvn2(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis, c0i, numed
2227 PASSCHARL(vname) PASSCHARL(vmother));
2230 //_____________________________________________________________________________
2231 void TGeant3::Gsdvs(const char *name, const char *mother, Float_t step,
2232 Int_t iaxis, Int_t numed)
2235 // Create a new volume by dividing an existing one
2240 Vname(mother,vmother);
2241 gsdvs(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed
2242 PASSCHARL(vname) PASSCHARL(vmother));
2245 //_____________________________________________________________________________
2246 void TGeant3::Gsdvs2(const char *name, const char *mother, Float_t step,
2247 Int_t iaxis, Float_t c0, Int_t numed)
2250 // Create a new volume by dividing an existing one
2255 Vname(mother,vmother);
2256 gsdvs2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0, numed
2257 PASSCHARL(vname) PASSCHARL(vmother));
2260 //_____________________________________________________________________________
2261 void TGeant3::Gsdvt(const char *name, const char *mother, Float_t step,
2262 Int_t iaxis, Int_t numed, Int_t ndvmx)
2265 // Create a new volume by dividing an existing one
2267 // Divides MOTHER into divisions called NAME along
2268 // axis IAXIS in steps of STEP. If not exactly divisible
2269 // will make as many as possible and will centre them
2270 // with respect to the mother. Divisions will have medium
2271 // number NUMED. If NUMED is 0, NUMED of MOTHER is taken.
2272 // NDVMX is the expected maximum number of divisions
2273 // (If 0, no protection tests are performed)
2278 Vname(mother,vmother);
2279 gsdvt(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed, ndvmx
2280 PASSCHARL(vname) PASSCHARL(vmother));
2283 //_____________________________________________________________________________
2284 void TGeant3::Gsdvt2(const char *name, const char *mother, Float_t step,
2285 Int_t iaxis, Float_t c0, Int_t numed, Int_t ndvmx)
2288 // Create a new volume by dividing an existing one
2290 // Divides MOTHER into divisions called NAME along
2291 // axis IAXIS starting at coordinate value C0 with step
2293 // The new volume created will have medium number NUMED.
2294 // If NUMED is 0, NUMED of mother is taken.
2295 // NDVMX is the expected maximum number of divisions
2296 // (If 0, no protection tests are performed)
2301 Vname(mother,vmother);
2302 gsdvt2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0,
2303 numed, ndvmx PASSCHARL(vname) PASSCHARL(vmother));
2306 //_____________________________________________________________________________
2307 void TGeant3::Gsord(const char *name, Int_t iax)
2310 // Flags volume CHNAME whose contents will have to be ordered
2311 // along axis IAX, by setting the search flag to -IAX
2315 // IAX = 4 Rxy (static ordering only -> GTMEDI)
2316 // IAX = 14 Rxy (also dynamic ordering -> GTNEXT)
2317 // IAX = 5 Rxyz (static ordering only -> GTMEDI)
2318 // IAX = 15 Rxyz (also dynamic ordering -> GTNEXT)
2319 // IAX = 6 PHI (PHI=0 => X axis)
2320 // IAX = 7 THETA (THETA=0 => Z axis)
2324 gsord(PASSCHARD(vname), iax PASSCHARL(vname));
2327 //_____________________________________________________________________________
2328 void TGeant3::Gspos(const char *name, Int_t nr, const char *mother, Float_t x,
2329 Float_t y, Float_t z, Int_t irot, const char *konly)
2332 // Position a volume into an existing one
2335 // NUMBER Copy number of the volume
2336 // MOTHER Mother volume name
2337 // X X coord. of the volume in mother ref. sys.
2338 // Y Y coord. of the volume in mother ref. sys.
2339 // Z Z coord. of the volume in mother ref. sys.
2340 // IROT Rotation matrix number w.r.t. mother ref. sys.
2341 // ONLY ONLY/MANY flag
2343 // It positions a previously defined volume in the mother.
2348 Vname(mother,vmother);
2349 gspos(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2350 PASSCHARD(konly) PASSCHARL(vname) PASSCHARL(vmother)
2354 //_____________________________________________________________________________
2355 void TGeant3::Gsposp(const char *name, Int_t nr, const char *mother,
2356 Float_t x, Float_t y, Float_t z, Int_t irot,
2357 const char *konly, Float_t *upar, Int_t np )
2360 // Place a copy of generic volume NAME with user number
2361 // NR inside MOTHER, with its parameters UPAR(1..NP)
2366 Vname(mother,vmother);
2367 gsposp(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2368 PASSCHARD(konly), upar, np PASSCHARL(vname) PASSCHARL(vmother)
2372 //_____________________________________________________________________________
2373 void TGeant3::Gsrotm(Int_t nmat, Float_t theta1, Float_t phi1, Float_t theta2,
2374 Float_t phi2, Float_t theta3, Float_t phi3)
2377 // nmat Rotation matrix number
2378 // THETA1 Polar angle for axis I
2379 // PHI1 Azimuthal angle for axis I
2380 // THETA2 Polar angle for axis II
2381 // PHI2 Azimuthal angle for axis II
2382 // THETA3 Polar angle for axis III
2383 // PHI3 Azimuthal angle for axis III
2385 // It defines the rotation matrix number IROT.
2387 gsrotm(nmat, theta1, phi1, theta2, phi2, theta3, phi3);
2390 //_____________________________________________________________________________
2391 void TGeant3::Gprotm(Int_t nmat)
2394 // To print rotation matrices structure JROTM
2395 // nmat Rotation matrix number
2400 //_____________________________________________________________________________
2401 Int_t TGeant3::Gsvolu(const char *name, const char *shape, Int_t nmed,
2402 Float_t *upar, Int_t npar)
2406 // SHAPE Volume type
2407 // NUMED Tracking medium number
2408 // NPAR Number of shape parameters
2409 // UPAR Vector containing shape parameters
2411 // It creates a new volume in the JVOLUM data structure.
2417 Vname(shape,vshape);
2418 gsvolu(PASSCHARD(vname), PASSCHARD(vshape), nmed, upar, npar, ivolu
2419 PASSCHARL(vname) PASSCHARL(vshape));
2423 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2425 // T H E D R A W I N G P A C K A G E
2426 // ======================================
2427 // Drawing functions. These functions allow the visualization in several ways
2428 // of the volumes defined in the geometrical data structure. It is possible
2429 // to draw the logical tree of volumes belonging to the detector (DTREE),
2430 // to show their geometrical specification (DSPEC,DFSPC), to draw them
2431 // and their cut views (DRAW, DCUT). Moreover, it is possible to execute
2432 // these commands when the hidden line removal option is activated; in
2433 // this case, the volumes can be also either translated in the space
2434 // (SHIFT), or clipped by boolean operation (CVOL). In addition, it is
2435 // possible to fill the surfaces of the volumes
2436 // with solid colours when the shading option (SHAD) is activated.
2437 // Several tools (ZOOM, LENS) have been developed to zoom detailed parts
2438 // of the detectors or to scan physical events as well.
2439 // Finally, the command MOVE will allow the rotation, translation and zooming
2440 // on real time parts of the detectors or tracks and hits of a simulated event.
2441 // Ray-tracing commands. In case the command (DOPT RAYT ON) is executed,
2442 // the drawing is performed by the Geant ray-tracing;
2443 // automatically, the color is assigned according to the tracking medium of each
2444 // volume and the volumes with a density lower/equal than the air are considered
2445 // transparent; if the option (USER) is set (ON) (again via the command (DOPT)),
2446 // the user can set color and visibility for the desired volumes via the command
2447 // (SATT), as usual, relatively to the attributes (COLO) and (SEEN).
2448 // The resolution can be set via the command (SATT * FILL VALUE), where (VALUE)
2449 // is the ratio between the number of pixels drawn and 20 (user coordinates).
2450 // Parallel view and perspective view are possible (DOPT PROJ PARA/PERS); in the
2451 // first case, we assume that the first mother volume of the tree is a box with
2452 // dimensions 10000 X 10000 X 10000 cm and the view point (infinetely far) is
2453 // 5000 cm far from the origin along the Z axis of the user coordinates; in the
2454 // second case, the distance between the observer and the origin of the world
2455 // reference system is set in cm by the command (PERSP NAME VALUE); grand-angle
2456 // or telescopic effects can be achieved changing the scale factors in the command
2457 // (DRAW). When the final picture does not occupy the full window,
2458 // mapping the space before tracing can speed up the drawing, but can also
2459 // produce less precise results; values from 1 to 4 are allowed in the command
2460 // (DOPT MAPP VALUE), the mapping being more precise for increasing (VALUE); for
2461 // (VALUE = 0) no mapping is performed (therefore max precision and lowest speed).
2462 // The command (VALCUT) allows the cutting of the detector by three planes
2463 // ortogonal to the x,y,z axis. The attribute (LSTY) can be set by the command
2464 // SATT for any desired volume and can assume values from 0 to 7; it determines
2465 // the different light processing to be performed for different materials:
2466 // 0 = dark-matt, 1 = bright-matt, 2 = plastic, 3 = ceramic, 4 = rough-metals,
2467 // 5 = shiny-metals, 6 = glass, 7 = mirror. The detector is assumed to be in the
2468 // dark, the ambient light luminosity is 0.2 for each basic hue (the saturation
2469 // is 0.9) and the observer is assumed to have a light source (therefore he will
2470 // produce parallel light in the case of parallel view and point-like-source
2471 // light in the case of perspective view).
2473 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2475 //_____________________________________________________________________________
2476 void TGeant3::Gsatt(const char *name, const char *att, Int_t val)
2480 // IOPT Name of the attribute to be set
2481 // IVAL Value to which the attribute is to be set
2483 // name= "*" stands for all the volumes.
2484 // iopt can be chosen among the following :
2486 // WORK 0=volume name is inactive for the tracking
2487 // 1=volume name is active for the tracking (default)
2489 // SEEN 0=volume name is invisible
2490 // 1=volume name is visible (default)
2491 // -1=volume invisible with all its descendants in the tree
2492 // -2=volume visible but not its descendants in the tree
2494 // LSTY line style 1,2,3,... (default=1)
2495 // LSTY=7 will produce a very precise approximation for
2496 // revolution bodies.
2498 // LWID line width -7,...,1,2,3,..7 (default=1)
2499 // LWID<0 will act as abs(LWID) was set for the volume
2500 // and for all the levels below it. When SHAD is 'ON', LWID
2501 // represent the linewidth of the scan lines filling the surfaces
2502 // (whereas the FILL value represent their number). Therefore
2503 // tuning this parameter will help to obtain the desired
2504 // quality/performance ratio.
2506 // COLO colour code -166,...,1,2,..166 (default=1)
2508 // n=2=red; n=17+m, m=0,25, increasing luminosity according to 'm';
2509 // n=3=green; n=67+m, m=0,25, increasing luminosity according to 'm';
2510 // n=4=blue; n=117+m, m=0,25, increasing luminosity according to 'm';
2511 // n=5=yellow; n=42+m, m=0,25, increasing luminosity according to 'm';
2512 // n=6=violet; n=142+m, m=0,25, increasing luminosity according to 'm';
2513 // n=7=lightblue; n=92+m, m=0,25, increasing luminosity according to 'm';
2514 // colour=n*10+m, m=1,2,...9, will produce the same colour
2515 // as 'n', but with increasing luminosity according to 'm';
2516 // COLO<0 will act as if abs(COLO) was set for the volume
2517 // and for all the levels below it.
2518 // When for a volume the attribute FILL is > 1 (and the
2519 // option SHAD is on), the ABS of its colour code must be < 8
2520 // because an automatic shading of its faces will be
2523 // FILL (1992) fill area -7,...,0,1,...7 (default=0)
2524 // when option SHAD is "on" the FILL attribute of any
2525 // volume can be set different from 0 (normal drawing);
2526 // if it is set to 1, the faces of such volume will be filled
2527 // with solid colours; if ABS(FILL) is > 1, then a light
2528 // source is placed along the observer line, and the faces of
2529 // such volumes will be painted by colours whose luminosity
2530 // will depend on the amount of light reflected;
2531 // if ABS(FILL) = 1, then it is possible to use all the 166
2532 // colours of the colour table, becouse the automatic shading
2533 // is not performed;
2534 // for increasing values of FILL the drawing will be performed
2535 // with higher and higher resolution improving the quality (the
2536 // number of scan lines used to fill the faces increases with FILL);
2537 // it is possible to set different values of FILL
2538 // for different volumes, in order to optimize at the same time
2539 // the performance and the quality of the picture;
2540 // FILL<0 will act as if abs(FILL) was set for the volume
2541 // and for all the levels below it.
2542 // This kind of drawing can be saved in 'picture files'
2543 // or in view banks.
2544 // 0=drawing without fill area
2545 // 1=faces filled with solid colours and resolution = 6
2546 // 2=lowest resolution (very fast)
2547 // 3=default resolution
2548 // 4=.................
2549 // 5=.................
2550 // 6=.................
2552 // Finally, if a coloured background is desired, the FILL
2553 // attribute for the first volume of the tree must be set
2554 // equal to -abs(colo), colo being >0 and <166.
2556 // SET set number associated to volume name
2557 // DET detector number associated to volume name
2558 // DTYP detector type (1,2)
2565 gsatt(PASSCHARD(vname), PASSCHARD(vatt), val PASSCHARL(vname)
2569 //_____________________________________________________________________________
2570 void TGeant3::Gfpara(const char *name, Int_t number, Int_t intext, Int_t& npar,
2571 Int_t& natt, Float_t* par, Float_t* att)
2574 // Find the parameters of a volume
2576 gfpara(PASSCHARD(name), number, intext, npar, natt, par, att
2580 //_____________________________________________________________________________
2581 void TGeant3::Gckpar(Int_t ish, Int_t npar, Float_t* par)
2584 // Check the parameters of a shape
2586 gckpar(ish,npar,par);
2589 //_____________________________________________________________________________
2590 void TGeant3::Gckmat(Int_t itmed, char* natmed)
2593 // Check the parameters of a tracking medium
2595 gckmat(itmed, PASSCHARD(natmed) PASSCHARL(natmed));
2598 //_____________________________________________________________________________
2599 void TGeant3::Gdelete(Int_t iview)
2602 // IVIEW View number
2604 // It deletes a view bank from memory.
2609 //_____________________________________________________________________________
2610 void TGeant3::Gdopen(Int_t iview)
2613 // IVIEW View number
2615 // When a drawing is very complex and requires a long time to be
2616 // executed, it can be useful to store it in a view bank: after a
2617 // call to DOPEN and the execution of the drawing (nothing will
2618 // appear on the screen), and after a necessary call to DCLOSE,
2619 // the contents of the bank can be displayed in a very fast way
2620 // through a call to DSHOW; therefore, the detector can be easily
2621 // zoomed many times in different ways. Please note that the pictures
2622 // with solid colours can now be stored in a view bank or in 'PICTURE FILES'
2629 //_____________________________________________________________________________
2630 void TGeant3::Gdclose()
2633 // It closes the currently open view bank; it must be called after the
2634 // end of the drawing to be stored.
2639 //_____________________________________________________________________________
2640 void TGeant3::Gdshow(Int_t iview)
2643 // IVIEW View number
2645 // It shows on the screen the contents of a view bank. It
2646 // can be called after a view bank has been closed.
2651 //_____________________________________________________________________________
2652 void TGeant3::Gdopt(const char *name,const char *value)
2656 // VALUE Option value
2658 // To set/modify the drawing options.
2661 // THRZ ON Draw tracks in R vs Z
2662 // OFF (D) Draw tracks in X,Y,Z
2665 // PROJ PARA (D) Parallel projection
2667 // TRAK LINE (D) Trajectory drawn with lines
2668 // POIN " " with markers
2669 // HIDE ON Hidden line removal using the CG package
2670 // OFF (D) No hidden line removal
2671 // SHAD ON Fill area and shading of surfaces.
2672 // OFF (D) Normal hidden line removal.
2673 // RAYT ON Ray-tracing on.
2674 // OFF (D) Ray-tracing off.
2675 // EDGE OFF Does not draw contours when shad is on.
2676 // ON (D) Normal shading.
2677 // MAPP 1,2,3,4 Mapping before ray-tracing.
2678 // 0 (D) No mapping.
2679 // USER ON User graphics options in the raytracing.
2680 // OFF (D) Automatic graphics options.
2686 Vname(value,vvalue);
2687 gdopt(PASSCHARD(vname), PASSCHARD(vvalue) PASSCHARL(vname)
2691 //_____________________________________________________________________________
2692 void TGeant3::Gdraw(const char *name,Float_t theta, Float_t phi, Float_t psi,
2693 Float_t u0,Float_t v0,Float_t ul,Float_t vl)
2698 // THETA Viewing angle theta (for 3D projection)
2699 // PHI Viewing angle phi (for 3D projection)
2700 // PSI Viewing angle psi (for 2D rotation)
2701 // U0 U-coord. (horizontal) of volume origin
2702 // V0 V-coord. (vertical) of volume origin
2703 // SU Scale factor for U-coord.
2704 // SV Scale factor for V-coord.
2706 // This function will draw the volumes,
2707 // selected with their graphical attributes, set by the Gsatt
2708 // facility. The drawing may be performed with hidden line removal
2709 // and with shading effects according to the value of the options HIDE
2710 // and SHAD; if the option SHAD is ON, the contour's edges can be
2711 // drawn or not. If the option HIDE is ON, the detector can be
2712 // exploded (BOMB), clipped with different shapes (CVOL), and some
2713 // of its parts can be shifted from their original
2714 // position (SHIFT). When HIDE is ON, if
2715 // the drawing requires more than the available memory, the program
2716 // will evaluate and display the number of missing words
2717 // (so that the user can increase the
2718 // size of its ZEBRA store). Finally, at the end of each drawing (with HIDE on),
2719 // the program will print messages about the memory used and
2720 // statistics on the volumes' visibility.
2721 // The following commands will produce the drawing of a green
2722 // volume, specified by NAME, without using the hidden line removal
2723 // technique, using the hidden line removal technique,
2724 // with different linewidth and colour (red), with
2725 // solid colour, with shading of surfaces, and without edges.
2726 // Finally, some examples are given for the ray-tracing. (A possible
2727 // string for the NAME of the volume can be found using the command DTREE).
2733 if (fGcvdma->raytra != 1) {
2734 gdraw(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
2736 gdrayt(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
2740 //_____________________________________________________________________________
2741 void TGeant3::Gdrawc(const char *name,Int_t axis, Float_t cut,Float_t u0,
2742 Float_t v0,Float_t ul,Float_t vl)
2747 // CUTVAL Cut plane distance from the origin along the axis
2749 // U0 U-coord. (horizontal) of volume origin
2750 // V0 V-coord. (vertical) of volume origin
2751 // SU Scale factor for U-coord.
2752 // SV Scale factor for V-coord.
2754 // The cut plane is normal to caxis (X,Y,Z), corresponding to iaxis (1,2,3),
2755 // and placed at the distance cutval from the origin.
2756 // The resulting picture is seen from the the same axis.
2757 // When HIDE Mode is ON, it is possible to get the same effect with
2758 // the CVOL/BOX function.
2764 gdrawc(PASSCHARD(vname), axis,cut,u0,v0,ul,vl PASSCHARL(vname));
2767 //_____________________________________________________________________________
2768 void TGeant3::Gdrawx(const char *name,Float_t cutthe, Float_t cutphi,
2769 Float_t cutval, Float_t theta, Float_t phi, Float_t u0,
2770 Float_t v0,Float_t ul,Float_t vl)
2774 // CUTTHE Theta angle of the line normal to cut plane
2775 // CUTPHI Phi angle of the line normal to cut plane
2776 // CUTVAL Cut plane distance from the origin along the axis
2778 // THETA Viewing angle theta (for 3D projection)
2779 // PHI Viewing angle phi (for 3D projection)
2780 // U0 U-coord. (horizontal) of volume origin
2781 // V0 V-coord. (vertical) of volume origin
2782 // SU Scale factor for U-coord.
2783 // SV Scale factor for V-coord.
2785 // The cut plane is normal to the line given by the cut angles
2786 // cutthe and cutphi and placed at the distance cutval from the origin.
2787 // The resulting picture is seen from the viewing angles theta,phi.
2793 gdrawx(PASSCHARD(vname), cutthe,cutphi,cutval,theta,phi,u0,v0,ul,vl
2797 //_____________________________________________________________________________
2798 void TGeant3::Gdhead(Int_t isel, const char *name, Float_t chrsiz)
2803 // ISEL Option flag D=111110
2805 // CHRSIZ Character size (cm) of title NAME D=0.6
2808 // 0 to have only the header lines
2809 // xxxxx1 to add the text name centered on top of header
2810 // xxxx1x to add global detector name (first volume) on left
2811 // xxx1xx to add date on right
2812 // xx1xxx to select thick characters for text on top of header
2813 // x1xxxx to add the text 'EVENT NR x' on top of header
2814 // 1xxxxx to add the text 'RUN NR x' on top of header
2815 // NOTE that ISEL=x1xxx1 or ISEL=1xxxx1 are illegal choices,
2816 // i.e. they generate overwritten text.
2818 gdhead(isel,PASSCHARD(name),chrsiz PASSCHARL(name));
2821 //_____________________________________________________________________________
2822 void TGeant3::Gdman(Float_t u, Float_t v, const char *type)
2825 // Draw a 2D-man at position (U0,V0)
2827 // U U-coord. (horizontal) of the centre of man' R
2828 // V V-coord. (vertical) of the centre of man' R
2829 // TYPE D='MAN' possible values: 'MAN,WM1,WM2,WM3'
2831 // CALL GDMAN(u,v),CALL GDWMN1(u,v),CALL GDWMN2(u,v),CALL GDWMN2(u,v)
2832 // It superimposes the picure of a man or of a woman, chosen among
2833 // three different ones, with the same scale factors as the detector
2834 // in the current drawing.
2837 if (opt.Contains("WM1")) {
2839 } else if (opt.Contains("WM3")) {
2841 } else if (opt.Contains("WM2")) {
2848 //_____________________________________________________________________________
2849 void TGeant3::Gdspec(const char *name)
2854 // Shows 3 views of the volume (two cut-views and a 3D view), together with
2855 // its geometrical specifications. The 3D drawing will
2856 // be performed according the current values of the options HIDE and
2857 // SHAD and according the current SetClipBox clipping parameters for that
2864 gdspec(PASSCHARD(vname) PASSCHARL(vname));
2867 //_____________________________________________________________________________
2868 void TGeant3::DrawOneSpec(const char *name)
2871 // Function called when one double-clicks on a volume name
2872 // in a TPavelabel drawn by Gdtree.
2874 THIGZ *higzSave = higz;
2875 higzSave->SetName("higzSave");
2876 THIGZ *higzSpec = (THIGZ*)gROOT->FindObject("higzSpec");
2877 //printf("DrawOneSpec, higz=%x, higzSpec=%x\n",higz,higzSpec);
2878 if (higzSpec) higz = higzSpec;
2879 else higzSpec = new THIGZ(defSize);
2880 higzSpec->SetName("higzSpec");
2885 gdspec(PASSCHARD(vname) PASSCHARL(vname));
2888 higzSave->SetName("higz");
2892 //_____________________________________________________________________________
2893 void TGeant3::Gdtree(const char *name,Int_t levmax, Int_t isel)
2897 // LEVMAX Depth level
2900 // This function draws the logical tree,
2901 // Each volume in the tree is represented by a TPaveTree object.
2902 // Double-clicking on a TPaveTree draws the specs of the corresponding volume.
2903 // Use TPaveTree pop-up menu to select:
2906 // - drawing tree of parent
2912 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
2916 //_____________________________________________________________________________
2917 void TGeant3::GdtreeParent(const char *name,Int_t levmax, Int_t isel)
2921 // LEVMAX Depth level
2924 // This function draws the logical tree of the parent of name.
2928 // Scan list of volumes in JVOLUM
2930 Int_t gname, i, jvo, in, nin, jin, num;
2931 strncpy((char *) &gname, name, 4);
2932 for(i=1; i<=fGcnum->nvolum; i++) {
2933 jvo = fZlq[fGclink->jvolum-i];
2934 nin = Int_t(fZq[jvo+3]);
2935 if (nin == -1) nin = 1;
2936 for (in=1;in<=nin;in++) {
2938 num = Int_t(fZq[jin+2]);
2939 if(gname == fZiq[fGclink->jvolum+num]) {
2940 strncpy(vname,(char*)&fZiq[fGclink->jvolum+i],4);
2942 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
2950 //_____________________________________________________________________________
2951 void TGeant3::SetABAN(Int_t par)
2954 // par = 1 particles will be stopped according to their residual
2955 // range if they are not in a sensitive material and are
2956 // far enough from the boundary
2957 // 0 particles are transported normally
2959 fGcphys->dphys1 = par;
2963 //_____________________________________________________________________________
2964 void TGeant3::SetANNI(Int_t par)
2967 // To control positron annihilation.
2968 // par =0 no annihilation
2969 // =1 annihilation. Decays processed.
2970 // =2 annihilation. No decay products stored.
2972 fGcphys->ianni = par;
2976 //_____________________________________________________________________________
2977 void TGeant3::SetAUTO(Int_t par)
2980 // To control automatic calculation of tracking medium parameters:
2981 // par =0 no automatic calculation;
2982 // =1 automati calculation.
2984 fGctrak->igauto = par;
2988 //_____________________________________________________________________________
2989 void TGeant3::SetBOMB(Float_t boom)
2992 // BOOM : Exploding factor for volumes position
2994 // To 'explode' the detector. If BOOM is positive (values smaller
2995 // than 1. are suggested, but any value is possible)
2996 // all the volumes are shifted by a distance
2997 // proportional to BOOM along the direction between their centre
2998 // and the origin of the MARS; the volumes which are symmetric
2999 // with respect to this origin are simply not shown.
3000 // BOOM equal to 0 resets the normal mode.
3001 // A negative (greater than -1.) value of
3002 // BOOM will cause an 'implosion'; for even lower values of BOOM
3003 // the volumes' positions will be reflected respect to the origin.
3004 // This command can be useful to improve the 3D effect for very
3005 // complex detectors. The following commands will make explode the
3012 //_____________________________________________________________________________
3013 void TGeant3::SetBREM(Int_t par)
3016 // To control bremstrahlung.
3017 // par =0 no bremstrahlung
3018 // =1 bremstrahlung. Photon processed.
3019 // =2 bremstrahlung. No photon stored.
3021 fGcphys->ibrem = par;
3025 //_____________________________________________________________________________
3026 void TGeant3::SetCKOV(Int_t par)
3029 // To control Cerenkov production
3030 // par =0 no Cerenkov;
3032 // =2 Cerenkov with primary stopped at each step.
3034 fGctlit->itckov = par;
3038 //_____________________________________________________________________________
3039 void TGeant3::SetClipBox(const char *name,Float_t xmin,Float_t xmax,
3040 Float_t ymin,Float_t ymax,Float_t zmin,Float_t zmax)
3043 // The hidden line removal technique is necessary to visualize properly
3044 // very complex detectors. At the same time, it can be useful to visualize
3045 // the inner elements of a detector in detail. This function allows
3046 // subtractions (via boolean operation) of BOX shape from any part of
3047 // the detector, therefore showing its inner contents.
3048 // If "*" is given as the name of the
3049 // volume to be clipped, all volumes are clipped by the given box.
3050 // A volume can be clipped at most twice.
3051 // if a volume is explicitely clipped twice,
3052 // the "*" will not act on it anymore. Giving "." as the name
3053 // of the volume to be clipped will reset the clipping.
3055 // NAME Name of volume to be clipped
3057 // XMIN Lower limit of the Shape X coordinate
3058 // XMAX Upper limit of the Shape X coordinate
3059 // YMIN Lower limit of the Shape Y coordinate
3060 // YMAX Upper limit of the Shape Y coordinate
3061 // ZMIN Lower limit of the Shape Z coordinate
3062 // ZMAX Upper limit of the Shape Z coordinate
3064 // This function performs a boolean subtraction between the volume
3065 // NAME and a box placed in the MARS according the values of the given
3071 setclip(PASSCHARD(vname),xmin,xmax,ymin,ymax,zmin,zmax PASSCHARL(vname));
3074 //_____________________________________________________________________________
3075 void TGeant3::SetCOMP(Int_t par)
3078 // To control Compton scattering
3079 // par =0 no Compton
3080 // =1 Compton. Electron processed.
3081 // =2 Compton. No electron stored.
3084 fGcphys->icomp = par;
3087 //_____________________________________________________________________________
3088 void TGeant3::SetCUTS(Float_t cutgam,Float_t cutele,Float_t cutneu,
3089 Float_t cuthad,Float_t cutmuo ,Float_t bcute ,
3090 Float_t bcutm ,Float_t dcute ,Float_t dcutm ,
3091 Float_t ppcutm, Float_t tofmax)
3094 // CUTGAM Cut for gammas D=0.001
3095 // CUTELE Cut for electrons D=0.001
3096 // CUTHAD Cut for charged hadrons D=0.01
3097 // CUTNEU Cut for neutral hadrons D=0.01
3098 // CUTMUO Cut for muons D=0.01
3099 // BCUTE Cut for electron brems. D=-1.
3100 // BCUTM Cut for muon brems. D=-1.
3101 // DCUTE Cut for electron delta-rays D=-1.
3102 // DCUTM Cut for muon delta-rays D=-1.
3103 // PPCUTM Cut for e+e- pairs by muons D=0.01
3104 // TOFMAX Time of flight cut D=1.E+10
3106 // If the default values (-1.) for BCUTE ,BCUTM ,DCUTE ,DCUTM
3107 // are not modified, they will be set to CUTGAM,CUTGAM,CUTELE,CUTELE
3109 // If one of the parameters from CUTGAM to PPCUTM included
3110 // is modified, cross-sections and energy loss tables must be
3111 // recomputed via the function Gphysi.
3113 fGccuts->cutgam = cutgam;
3114 fGccuts->cutele = cutele;
3115 fGccuts->cutneu = cutneu;
3116 fGccuts->cuthad = cuthad;
3117 fGccuts->cutmuo = cutmuo;
3118 fGccuts->bcute = bcute;
3119 fGccuts->bcutm = bcutm;
3120 fGccuts->dcute = dcute;
3121 fGccuts->dcutm = dcutm;
3122 fGccuts->ppcutm = ppcutm;
3123 fGccuts->tofmax = tofmax;
3126 //_____________________________________________________________________________
3127 void TGeant3::SetDCAY(Int_t par)
3130 // To control Decay mechanism.
3131 // par =0 no decays.
3132 // =1 Decays. secondaries processed.
3133 // =2 Decays. No secondaries stored.
3135 fGcphys->idcay = par;
3139 //_____________________________________________________________________________
3140 void TGeant3::SetDEBU(Int_t emin, Int_t emax, Int_t emod)
3143 // Set the debug flag and frequency
3144 // Selected debug output will be printed from
3145 // event emin to even emax each emod event
3147 fGcflag->idemin = emin;
3148 fGcflag->idemax = emax;
3149 fGcflag->itest = emod;
3153 //_____________________________________________________________________________
3154 void TGeant3::SetDRAY(Int_t par)
3157 // To control delta rays mechanism.
3158 // par =0 no delta rays.
3159 // =1 Delta rays. secondaries processed.
3160 // =2 Delta rays. No secondaries stored.
3162 fGcphys->idray = par;
3165 //_____________________________________________________________________________
3166 void TGeant3::SetHADR(Int_t par)
3169 // To control hadronic interactions.
3170 // par =0 no hadronic interactions.
3171 // =1 Hadronic interactions. secondaries processed.
3172 // =2 Hadronic interactions. No secondaries stored.
3174 fGcphys->ihadr = par;
3177 //_____________________________________________________________________________
3178 void TGeant3::SetKINE(Int_t kine, Float_t xk1, Float_t xk2, Float_t xk3,
3179 Float_t xk4, Float_t xk5, Float_t xk6, Float_t xk7,
3180 Float_t xk8, Float_t xk9, Float_t xk10)
3183 // Set the variables in /GCFLAG/ IKINE, PKINE(10)
3184 // Their meaning is user defined
3186 fGckine->ikine = kine;
3187 fGckine->pkine[0] = xk1;
3188 fGckine->pkine[1] = xk2;
3189 fGckine->pkine[2] = xk3;
3190 fGckine->pkine[3] = xk4;
3191 fGckine->pkine[4] = xk5;
3192 fGckine->pkine[5] = xk6;
3193 fGckine->pkine[6] = xk7;
3194 fGckine->pkine[7] = xk8;
3195 fGckine->pkine[8] = xk9;
3196 fGckine->pkine[9] = xk10;
3199 //_____________________________________________________________________________
3200 void TGeant3::SetLOSS(Int_t par)
3203 // To control energy loss.
3204 // par =0 no energy loss;
3205 // =1 restricted energy loss fluctuations;
3206 // =2 complete energy loss fluctuations;
3208 // =4 no energy loss fluctuations.
3209 // If the value ILOSS is changed, then cross-sections and energy loss
3210 // tables must be recomputed via the command 'PHYSI'.
3212 fGcphys->iloss = par;
3216 //_____________________________________________________________________________
3217 void TGeant3::SetMULS(Int_t par)
3220 // To control multiple scattering.
3221 // par =0 no multiple scattering.
3222 // =1 Moliere or Coulomb scattering.
3223 // =2 Moliere or Coulomb scattering.
3224 // =3 Gaussian scattering.
3226 fGcphys->imuls = par;
3230 //_____________________________________________________________________________
3231 void TGeant3::SetMUNU(Int_t par)
3234 // To control muon nuclear interactions.
3235 // par =0 no muon-nuclear interactions.
3236 // =1 Nuclear interactions. Secondaries processed.
3237 // =2 Nuclear interactions. Secondaries not processed.
3239 fGcphys->imunu = par;
3242 //_____________________________________________________________________________
3243 void TGeant3::SetOPTI(Int_t par)
3246 // This flag controls the tracking optimisation performed via the
3248 // 1 no optimisation at all; GSORD calls disabled;
3249 // 0 no optimisation; only user calls to GSORD kept;
3250 // 1 all non-GSORDered volumes are ordered along the best axis;
3251 // 2 all volumes are ordered along the best axis.
3253 fGcopti->ioptim = par;
3256 //_____________________________________________________________________________
3257 void TGeant3::SetPAIR(Int_t par)
3260 // To control pair production mechanism.
3261 // par =0 no pair production.
3262 // =1 Pair production. secondaries processed.
3263 // =2 Pair production. No secondaries stored.
3265 fGcphys->ipair = par;
3269 //_____________________________________________________________________________
3270 void TGeant3::SetPFIS(Int_t par)
3273 // To control photo fission mechanism.
3274 // par =0 no photo fission.
3275 // =1 Photo fission. secondaries processed.
3276 // =2 Photo fission. No secondaries stored.
3278 fGcphys->ipfis = par;
3281 //_____________________________________________________________________________
3282 void TGeant3::SetPHOT(Int_t par)
3285 // To control Photo effect.
3286 // par =0 no photo electric effect.
3287 // =1 Photo effect. Electron processed.
3288 // =2 Photo effect. No electron stored.
3290 fGcphys->iphot = par;
3293 //_____________________________________________________________________________
3294 void TGeant3::SetRAYL(Int_t par)
3297 // To control Rayleigh scattering.
3298 // par =0 no Rayleigh scattering.
3301 fGcphys->irayl = par;
3304 //_____________________________________________________________________________
3305 void TGeant3::SetSWIT(Int_t sw, Int_t val)
3309 // val New switch value
3311 // Change one element of array ISWIT(10) in /GCFLAG/
3313 if (sw <= 0 || sw > 10) return;
3314 fGcflag->iswit[sw-1] = val;
3318 //_____________________________________________________________________________
3319 void TGeant3::SetTRIG(Int_t nevents)
3322 // Set number of events to be run
3324 fGcflag->nevent = nevents;
3327 //_____________________________________________________________________________
3328 void TGeant3::SetUserDecay(Int_t pdg)
3331 // Force the decays of particles to be done with Pythia
3332 // and not with the Geant routines.
3333 // just kill pointers doing mzdrop
3335 Int_t ipart = IdFromPDG(pdg);
3337 printf("Particle %d not in geant\n",pdg);
3340 Int_t jpart=fGclink->jpart;
3341 Int_t jpa=fZlq[jpart-ipart];
3344 Int_t jpa1=fZlq[jpa-1];
3346 mzdrop(fGcbank->ixcons,jpa1,PASSCHARD(" ") PASSCHARL(" "));
3347 Int_t jpa2=fZlq[jpa-2];
3349 mzdrop(fGcbank->ixcons,jpa2,PASSCHARD(" ") PASSCHARL(" "));
3353 //______________________________________________________________________________
3354 void TGeant3::Vname(const char *name, char *vname)
3357 // convert name to upper case. Make vname at least 4 chars
3359 Int_t l = strlen(name);
3362 for (i=0;i<l;i++) vname[i] = toupper(name[i]);
3363 for (i=l;i<4;i++) vname[i] = ' ';
3367 //______________________________________________________________________________
3368 void TGeant3::Ertrgo()
3373 //______________________________________________________________________________
3374 void TGeant3::Ertrak(const Float_t *const x1, const Float_t *const p1,
3375 const Float_t *x2, const Float_t *p2,
3376 Int_t ipa, Option_t *chopt)
3378 ertrak(x1,p1,x2,p2,ipa,PASSCHARD(chopt) PASSCHARL(chopt));
3381 //_____________________________________________________________________________
3382 void TGeant3::WriteEuclid(const char* filnam, const char* topvol,
3383 Int_t number, Int_t nlevel)
3387 // ******************************************************************
3389 // * Write out the geometry of the detector in EUCLID file format *
3391 // * filnam : will be with the extension .euc *
3392 // * topvol : volume name of the starting node *
3393 // * number : copy number of topvol (relevant for gsposp) *
3394 // * nlevel : number of levels in the tree structure *
3395 // * to be written out, starting from topvol *
3397 // * Author : M. Maire *
3399 // ******************************************************************
3401 // File filnam.tme is written out with the definitions of tracking
3402 // medias and materials.
3403 // As to restore original numbers for materials and medias, program
3404 // searches in the file euc_medi.dat and comparing main parameters of
3405 // the mat. defined inside geant and the one in file recognizes them
3406 // and is able to take number from file. If for any material or medium,
3407 // this procedure fails, ordering starts from 1.
3408 // Arrays IOTMED and IOMATE are used for this procedure
3410 const char shape[][5]={"BOX ","TRD1","TRD2","TRAP","TUBE","TUBS","CONE",
3411 "CONS","SPHE","PARA","PGON","PCON","ELTU","HYPE",
3413 Int_t i, end, itm, irm, jrm, k, nmed;
3417 char *filext, *filetme;
3418 char natmed[21], namate[21];
3419 char natmedc[21], namatec[21];
3420 char key[5], name[5], mother[5], konly[5];
3422 Int_t iadvol, iadtmd, iadrot, nwtot, iret;
3423 Int_t mlevel, numbr, natt, numed, nin, ndata;
3424 Int_t iname, ivo, ish, jvo, nvstak, ivstak;
3425 Int_t jdiv, ivin, in, jin, jvin, irot;
3426 Int_t jtm, imat, jma, flag=0, imatc;
3427 Float_t az, dens, radl, absl, a, step, x, y, z;
3428 Int_t npar, ndvmx, left;
3429 Float_t zc, densc, radlc, abslc, c0, tmaxfd;
3431 Int_t iomate[100], iotmed[100];
3432 Float_t par[50], att[20], ubuf[50];
3435 Int_t level, ndiv, iaxe;
3436 Int_t itmedc, nmatc, isvolc, ifieldc, nwbufc, isvol, nmat, ifield, nwbuf;
3437 Float_t fieldmc, tmaxfdc, stemaxc, deemaxc, epsilc, stminc, fieldm;
3438 Float_t tmaxf, stemax, deemax, epsil, stmin;
3439 const char *f10000="!\n%s\n!\n";
3440 //Open the input file
3442 for(i=0;i<end;i++) if(filnam[i]=='.') {
3446 filext=new char[end+4];
3447 filetme=new char[end+4];
3448 strncpy(filext,filnam,end);
3449 strncpy(filetme,filnam,end);
3451 // *** The output filnam name will be with extension '.euc'
3452 strcpy(&filext[end],".euc");
3453 strcpy(&filetme[end],".tme");
3454 lun=fopen(filext,"w");
3456 // *** Initialisation of the working space
3457 iadvol=fGcnum->nvolum;
3458 iadtmd=iadvol+fGcnum->nvolum;
3459 iadrot=iadtmd+fGcnum->ntmed;
3460 if(fGclink->jrotm) {
3461 fGcnum->nrotm=fZiq[fGclink->jrotm-2];
3465 nwtot=iadrot+fGcnum->nrotm;
3466 qws = new float[nwtot+1];
3467 for (i=0;i<nwtot+1;i++) qws[i]=0;
3470 if(nlevel==0) mlevel=20;
3472 // *** find the top volume and put it in the stak
3473 numbr = number>0 ? number : 1;
3474 Gfpara(topvol,numbr,1,npar,natt,par,att);
3476 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3481 // *** authorized shape ?
3482 strncpy((char *)&iname, topvol, 4);
3484 for(i=1; i<=fGcnum->nvolum; i++) if(fZiq[fGclink->jvolum+i]==iname) {
3488 jvo = fZlq[fGclink->jvolum-ivo];
3489 ish = Int_t (fZq[jvo+2]);
3491 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3498 iws[iadvol+ivo] = level;
3501 //*** flag all volumes and fill the stak
3505 // pick the next volume in stak
3507 ivo = TMath::Abs(iws[ivstak]);
3508 jvo = fZlq[fGclink->jvolum - ivo];
3510 // flag the tracking medium
3511 numed = Int_t (fZq[jvo + 4]);
3512 iws[iadtmd + numed] = 1;
3514 // get the daughters ...
3515 level = iws[iadvol+ivo];
3516 if (level < mlevel) {
3518 nin = Int_t (fZq[jvo + 3]);
3520 // from division ...
3522 jdiv = fZlq[jvo - 1];
3523 ivin = Int_t (fZq[jdiv + 2]);
3525 iws[nvstak] = -ivin;
3526 iws[iadvol+ivin] = level;
3528 // from position ...
3529 } else if (nin > 0) {
3530 for(in=1; in<=nin; in++) {
3531 jin = fZlq[jvo - in];
3532 ivin = Int_t (fZq[jin + 2 ]);
3533 jvin = fZlq[fGclink->jvolum - ivin];
3534 ish = Int_t (fZq[jvin + 2]);
3535 // authorized shape ?
3537 // not yet flagged ?
3538 if (iws[iadvol+ivin]==0) {
3541 iws[iadvol+ivin] = level;
3543 // flag the rotation matrix
3544 irot = Int_t ( fZq[jin + 4 ]);
3545 if (irot > 0) iws[iadrot+irot] = 1;
3551 // next volume in stak ?
3552 if (ivstak < nvstak) goto L10;
3554 // *** restore original material and media numbers
3555 // file euc_medi.dat is needed to compare materials and medias
3557 FILE* luncor=fopen("euc_medi.dat","r");
3560 for(itm=1; itm<=fGcnum->ntmed; itm++) {
3561 if (iws[iadtmd+itm] > 0) {
3562 jtm = fZlq[fGclink->jtmed-itm];
3563 strncpy(natmed,(char *)&fZiq[jtm+1],20);
3564 imat = Int_t (fZq[jtm+6]);
3565 jma = fZlq[fGclink->jmate-imat];
3567 printf(" *** GWEUCL *** material not defined for tracking medium %5i %s\n",itm,natmed);
3570 strncpy(namate,(char *)&fZiq[jma+1],20);
3573 //** find the material original number
3576 iret=fscanf(luncor,"%4s,%130s",key,card);
3577 if(iret<=0) goto L26;
3579 if(!strcmp(key,"MATE")) {
3580 sscanf(card,"%d %s %f %f %f %f %f %d",&imatc,namatec,&az,&zc,&densc,&radlc,&abslc,&nparc);
3581 Gfmate(imat,namate,a,z,dens,radl,absl,par,npar);
3582 if(!strcmp(namatec,namate)) {
3583 if(az==a && zc==z && densc==dens && radlc==radl
3584 && abslc==absl && nparc==nparc) {
3587 printf("*** GWEUCL *** material : %3d '%s' restored as %3d\n",imat,namate,imatc);
3589 printf("*** GWEUCL *** different definitions for material: %s\n",namate);
3593 if(strcmp(key,"END") && !flag) goto L23;
3595 printf("*** GWEUCL *** cannot restore original number for material: %s\n",namate);
3599 //*** restore original tracking medium number
3602 iret=fscanf(luncor,"%4s,%130s",key,card);
3603 if(iret<=0) goto L26;
3605 if (!strcmp(key,"TMED")) {
3606 sscanf(card,"%d %s %d %d %d %f %f %f %f %f %f %d\n",
3607 &itmedc,natmedc,&nmatc,&isvolc,&ifieldc,&fieldmc,
3608 &tmaxfdc,&stemaxc,&deemaxc,&epsilc,&stminc,&nwbufc);
3609 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxf,stemax,deemax,
3610 epsil,stmin,ubuf,&nwbuf);
3611 if(!strcmp(natmedc,natmed)) {
3612 if (iomate[nmat]==nmatc && nwbuf==nwbufc) {
3615 printf("*** GWEUCL *** medium : %3d '%20s' restored as %3d\n",
3618 printf("*** GWEUCL *** different definitions for tracking medium: %s\n",natmed);
3622 if(strcmp(key,"END") && !flag) goto L24;
3624 printf("cannot restore original number for medium : %s\n",natmed);
3632 L26: printf("*** GWEUCL *** cannot read the data file\n");
3634 L29: if(luncor) fclose (luncor);
3637 // *** write down the tracking medium definition
3639 strcpy(card,"! Tracking medium");
3640 fprintf(lun,f10000,card);
3642 for(itm=1;itm<=fGcnum->ntmed;itm++) {
3643 if (iws[iadtmd+itm]>0) {
3644 jtm = fZlq[fGclink->jtmed-itm];
3645 strncpy(natmed,(char *)&fZiq[jtm+1],20);
3647 imat = Int_t (fZq[jtm+6]);
3648 jma = fZlq[fGclink->jmate-imat];
3649 //* order media from one, if comparing with database failed
3651 iotmed[itm]=++imxtmed;
3652 iomate[imat]=++imxmate;
3657 printf(" *** GWEUCL *** material not defined for tracking medium %5d %s\n",
3660 strncpy(namate,(char *)&fZiq[jma+1],20);
3663 fprintf(lun,"TMED %3d '%20s' %3d '%20s'\n",iotmed[itm],natmed,iomate[imat],namate);
3667 //* *** write down the rotation matrix
3669 strcpy(card,"! Reperes");
3670 fprintf(lun,f10000,card);
3672 for(irm=1;irm<=fGcnum->nrotm;irm++) {
3673 if (iws[iadrot+irm]>0) {
3674 jrm = fZlq[fGclink->jrotm-irm];
3675 fprintf(lun,"ROTM %3d",irm);
3676 for(k=11;k<=16;k++) fprintf(lun," %8.3f",fZq[jrm+k]);
3681 //* *** write down the volume definition
3683 strcpy(card,"! Volumes");
3684 fprintf(lun,f10000,card);
3686 for(ivstak=1;ivstak<=nvstak;ivstak++) {
3689 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivo],4);
3691 jvo = fZlq[fGclink->jvolum-ivo];
3692 ish = Int_t (fZq[jvo+2]);
3693 nmed = Int_t (fZq[jvo+4]);
3694 npar = Int_t (fZq[jvo+5]);
3696 if (ivstak>1) for(i=0;i<npar;i++) par[i]=fZq[jvo+7+i];
3697 Gckpar (ish,npar,par);
3698 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,shape[ish-1],iotmed[nmed],npar);
3699 for(i=0;i<(npar-1)/6+1;i++) {
3702 for(k=0;k<(left<6?left:6);k++) fprintf(lun," %11.5f",par[i*6+k]);
3706 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,shape[ish-1],iotmed[nmed],npar);
3711 //* *** write down the division of volumes
3713 fprintf(lun,f10000,"! Divisions");
3714 for(ivstak=1;ivstak<=nvstak;ivstak++) {
3715 ivo = TMath::Abs(iws[ivstak]);
3716 jvo = fZlq[fGclink->jvolum-ivo];
3717 ish = Int_t (fZq[jvo+2]);
3718 nin = Int_t (fZq[jvo+3]);
3719 //* this volume is divided ...
3722 iaxe = Int_t ( fZq[jdiv+1]);
3723 ivin = Int_t ( fZq[jdiv+2]);
3724 ndiv = Int_t ( fZq[jdiv+3]);
3727 jvin = fZlq[fGclink->jvolum-ivin];
3728 nmed = Int_t ( fZq[jvin+4]);
3729 strncpy(mother,(char *)&fZiq[fGclink->jvolum+ivo ],4);
3731 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivin],4);
3733 if ((step<=0.)||(ish>=11)) {
3734 //* volume with negative parameter or gsposp or pgon ...
3735 fprintf(lun,"DIVN '%4s' '%4s' %3d %3d\n",name,mother,ndiv,iaxe);
3736 } else if ((ndiv<=0)||(ish==10)) {
3737 //* volume with negative parameter or gsposp or para ...
3738 ndvmx = TMath::Abs(ndiv);
3739 fprintf(lun,"DIVT '%4s' '%4s' %11.5f %3d %3d %3d\n",
3740 name,mother,step,iaxe,iotmed[nmed],ndvmx);
3742 //* normal volume : all kind of division are equivalent
3743 fprintf(lun,"DVT2 '%4s' '%4s' %11.5f %3d %11.5f %3d %3d\n",
3744 name,mother,step,iaxe,c0,iotmed[nmed],ndiv);
3749 //* *** write down the the positionnement of volumes
3751 fprintf(lun,f10000,"! Positionnements\n");
3753 for(ivstak = 1;ivstak<=nvstak;ivstak++) {
3754 ivo = TMath::Abs(iws[ivstak]);
3755 strncpy(mother,(char*)&fZiq[fGclink->jvolum+ivo ],4);
3757 jvo = fZlq[fGclink->jvolum-ivo];
3758 nin = Int_t( fZq[jvo+3]);
3759 //* this volume has daughters ...
3761 for (in=1;in<=nin;in++) {
3763 ivin = Int_t (fZq[jin +2]);
3764 numb = Int_t (fZq[jin +3]);
3765 irot = Int_t (fZq[jin +4]);
3769 strcpy(konly,"ONLY");
3770 if (fZq[jin+8]!=1.) strcpy(konly,"MANY");
3771 strncpy(name,(char*)&fZiq[fGclink->jvolum+ivin],4);
3773 jvin = fZlq[fGclink->jvolum-ivin];
3774 ish = Int_t (fZq[jvin+2]);
3775 //* gspos or gsposp ?
3776 ndata = fZiq[jin-1];
3778 fprintf(lun,"POSI '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s'\n",
3779 name,numb,mother,x,y,z,irot,konly);
3781 npar = Int_t (fZq[jin+9]);
3782 for(i=0;i<npar;i++) par[i]=fZq[jin+10+i];
3783 Gckpar (ish,npar,par);
3784 fprintf(lun,"POSP '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s' %3d\n",
3785 name,numb,mother,x,y,z,irot,konly,npar);
3787 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
3794 fprintf(lun,"END\n");
3797 //****** write down the materials and medias *****
3799 lun=fopen(filetme,"w");
3801 for(itm=1;itm<=fGcnum->ntmed;itm++) {
3802 if (iws[iadtmd+itm]>0) {
3803 jtm = fZlq[fGclink->jtmed-itm];
3804 strncpy(natmed,(char*)&fZiq[jtm+1],4);
3805 imat = Int_t (fZq[jtm+6]);
3806 jma = Int_t (fZlq[fGclink->jmate-imat]);
3808 Gfmate (imat,namate,a,z,dens,radl,absl,par,npar);
3809 fprintf(lun,"MATE %4d '%20s'%11.5E %11.5E %11.5E %11.5E %11.5E %3d\n",
3810 iomate[imat],namate,a,z,dens,radl,absl,npar);
3814 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
3818 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxfd,stemax,deemax,epsil,stmin,par,&npar);
3819 fprintf(lun,"TMED %4d '%20s' %3d %1d %3d %11.5f %11.5f %11.5f %11.5f %11.5f %11.5f %3d\n",
3820 iotmed[itm],natmed,iomate[nmat],isvol,ifield,
3821 fieldm,tmaxfd,stemax,deemax,epsil,stmin,npar);
3825 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
3831 fprintf(lun,"END\n");
3832 printf(" *** GWEUCL *** file: %s is now written out\n",filext);
3833 printf(" *** GWEUCL *** file: %s is now written out\n",filetme);
3842 //_____________________________________________________________________________
3843 void TGeant3::Streamer(TBuffer &R__b)
3846 // Stream an object of class TGeant3.
3848 if (R__b.IsReading()) {
3849 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
3850 AliMC::Streamer(R__b);
3853 R__b.ReadStaticArray(fPDGCode);
3855 R__b.WriteVersion(TGeant3::IsA());
3856 AliMC::Streamer(R__b);
3859 R__b.WriteArray(fPDGCode, fNPDGCodes);