1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 Revision 1.24 2000/02/28 21:03:57 fca
19 Some additions to improve the compatibility with G4
21 Revision 1.23 2000/02/23 16:25:25 fca
22 AliVMC and AliGeant3 classes introduced
23 ReadEuclid moved from AliRun to AliModule
25 Revision 1.22 2000/01/18 15:40:13 morsch
26 Interface to GEANT3 routines GFTMAT, GBRELM and GPRELM added
27 Define geant particle type 51: Feedback Photon with Cherenkov photon properties.
29 Revision 1.21 2000/01/17 19:41:17 fca
32 Revision 1.20 2000/01/12 11:29:27 fca
35 Revision 1.19 1999/12/17 09:03:12 fca
36 Introduce a names array
38 Revision 1.18 1999/11/26 16:55:39 fca
39 Reimplement CurrentVolName() to avoid memory leaks
41 Revision 1.17 1999/11/03 16:58:28 fca
42 Correct source of address violation in creating character string
44 Revision 1.16 1999/11/03 13:17:08 fca
45 Have ProdProcess return const char*
47 Revision 1.15 1999/10/26 06:04:50 fca
48 Introduce TLorentzVector in AliMC::GetSecondary. Thanks to I.Hrivnacova
50 Revision 1.14 1999/09/29 09:24:30 fca
51 Introduction of the Copyright and cvs Log
55 ///////////////////////////////////////////////////////////////////////////////
57 // Interface Class to the Geant3.21 MonteCarlo //
61 <img src="picts/TGeant3Class.gif">
66 ///////////////////////////////////////////////////////////////////////////////
72 #include <TDatabasePDG.h>
73 #include "AliCallf77.h"
76 # define gzebra gzebra_
77 # define grfile grfile_
78 # define gpcxyz gpcxyz_
79 # define ggclos ggclos_
82 # define gcinit gcinit_
85 # define gtrigc gtrigc_
86 # define gtrigi gtrigi_
88 # define gzinit gzinit_
89 # define gfmate gfmate_
90 # define gfpart gfpart_
91 # define gftmed gftmed_
92 # define gftmat gftmat_
96 # define gsmate gsmate_
97 # define gsmixt gsmixt_
98 # define gspart gspart_
99 # define gstmed gstmed_
100 # define gsckov gsckov_
101 # define gstpar gstpar_
102 # define gfkine gfkine_
103 # define gfvert gfvert_
104 # define gskine gskine_
105 # define gsvert gsvert_
106 # define gphysi gphysi_
107 # define gdebug gdebug_
108 # define gekbin gekbin_
109 # define gfinds gfinds_
110 # define gsking gsking_
111 # define gskpho gskpho_
112 # define gsstak gsstak_
113 # define gsxyz gsxyz_
114 # define gtrack gtrack_
115 # define gtreve gtreve_
116 # define gtreve_root gtreve_root_
117 # define grndm grndm_
118 # define grndmq grndmq_
119 # define gdtom gdtom_
120 # define glmoth glmoth_
121 # define gmedia gmedia_
122 # define gmtod gmtod_
123 # define gsdvn gsdvn_
124 # define gsdvn2 gsdvn2_
125 # define gsdvs gsdvs_
126 # define gsdvs2 gsdvs2_
127 # define gsdvt gsdvt_
128 # define gsdvt2 gsdvt2_
129 # define gsord gsord_
130 # define gspos gspos_
131 # define gsposp gsposp_
132 # define gsrotm gsrotm_
133 # define gprotm gprotm_
134 # define gsvolu gsvolu_
135 # define gprint gprint_
136 # define gdinit gdinit_
137 # define gdopt gdopt_
138 # define gdraw gdraw_
139 # define gdrayt gdrayt_
140 # define gdrawc gdrawc_
141 # define gdrawx gdrawx_
142 # define gdhead gdhead_
143 # define gdwmn1 gdwmn1_
144 # define gdwmn2 gdwmn2_
145 # define gdwmn3 gdwmn3_
146 # define gdxyz gdxyz_
147 # define gdcxyz gdcxyz_
148 # define gdman gdman_
149 # define gdspec gdspec_
150 # define gdtree gdtree_
151 # define gdelet gdelet_
152 # define gdclos gdclos_
153 # define gdshow gdshow_
154 # define gdopen gdopen_
155 # define dzshow dzshow_
156 # define gsatt gsatt_
157 # define gfpara gfpara_
158 # define gckpar gckpar_
159 # define gckmat gckmat_
160 # define geditv geditv_
161 # define mzdrop mzdrop_
163 # define ertrak ertrak_
164 # define ertrgo ertrgo_
166 # define setbomb setbomb_
167 # define setclip setclip_
168 # define gcomad gcomad_
170 # define gbrelm gbrelm_
171 # define gprelm gprelm_
173 # define gzebra GZEBRA
174 # define grfile GRFILE
175 # define gpcxyz GPCXYZ
176 # define ggclos GGCLOS
179 # define gcinit GCINIT
182 # define gtrigc GTRIGC
183 # define gtrigi GTRIGI
185 # define gzinit GZINIT
186 # define gfmate GFMATE
187 # define gfpart GFPART
188 # define gftmed GFTMED
189 # define gftmat GFTMAT
193 # define gsmate GSMATE
194 # define gsmixt GSMIXT
195 # define gspart GSPART
196 # define gstmed GSTMED
197 # define gsckov GSCKOV
198 # define gstpar GSTPAR
199 # define gfkine GFKINE
200 # define gfvert GFVERT
201 # define gskine GSKINE
202 # define gsvert GSVERT
203 # define gphysi GPHYSI
204 # define gdebug GDEBUG
205 # define gekbin GEKBIN
206 # define gfinds GFINDS
207 # define gsking GSKING
208 # define gskpho GSKPHO
209 # define gsstak GSSTAK
211 # define gtrack GTRACK
212 # define gtreve GTREVE
213 # define gtreve_root GTREVE_ROOT
215 # define grndmq GRNDMQ
217 # define glmoth GLMOTH
218 # define gmedia GMEDIA
221 # define gsdvn2 GSDVN2
223 # define gsdvs2 GSDVS2
225 # define gsdvt2 GSDVT2
228 # define gsposp GSPOSP
229 # define gsrotm GSROTM
230 # define gprotm GPROTM
231 # define gsvolu GSVOLU
232 # define gprint GPRINT
233 # define gdinit GDINIT
236 # define gdrayt GDRAYT
237 # define gdrawc GDRAWC
238 # define gdrawx GDRAWX
239 # define gdhead GDHEAD
240 # define gdwmn1 GDWMN1
241 # define gdwmn2 GDWMN2
242 # define gdwmn3 GDWMN3
244 # define gdcxyz GDCXYZ
246 # define gdfspc GDFSPC
247 # define gdspec GDSPEC
248 # define gdtree GDTREE
249 # define gdelet GDELET
250 # define gdclos GDCLOS
251 # define gdshow GDSHOW
252 # define gdopen GDOPEN
253 # define dzshow DZSHOW
255 # define gfpara GFPARA
256 # define gckpar GCKPAR
257 # define gckmat GCKMAT
258 # define geditv GEDITV
259 # define mzdrop MZDROP
261 # define ertrak ERTRAK
262 # define ertrgo ERTRGO
264 # define setbomb SETBOMB
265 # define setclip SETCLIP
266 # define gcomad GCOMAD
268 # define gbrelm GBRELM
269 # define gprelm GPRELM
273 //____________________________________________________________________________
277 // Prototypes for GEANT functions
279 void type_of_call gzebra(const int&);
281 void type_of_call gpcxyz();
283 void type_of_call ggclos();
285 void type_of_call glast();
287 void type_of_call ginit();
289 void type_of_call gcinit();
291 void type_of_call grun();
293 void type_of_call gtrig();
295 void type_of_call gtrigc();
297 void type_of_call gtrigi();
299 void type_of_call gwork(const int&);
301 void type_of_call gzinit();
303 void type_of_call gmate();
305 void type_of_call gpart();
307 void type_of_call gsdk(Int_t &, Float_t *, Int_t *);
309 void type_of_call gfkine(Int_t &, Float_t *, Float_t *, Int_t &,
310 Int_t &, Float_t *, Int_t &);
312 void type_of_call gfvert(Int_t &, Float_t *, Int_t &, Int_t &,
313 Float_t &, Float_t *, Int_t &);
315 void type_of_call gskine(Float_t *,Int_t &, Int_t &, Float_t *,
318 void type_of_call gsvert(Float_t *,Int_t &, Int_t &, Float_t *,
321 void type_of_call gphysi();
323 void type_of_call gdebug();
325 void type_of_call gekbin();
327 void type_of_call gfinds();
329 void type_of_call gsking(Int_t &);
331 void type_of_call gskpho(Int_t &);
333 void type_of_call gsstak(Int_t &);
335 void type_of_call gsxyz();
337 void type_of_call gtrack();
339 void type_of_call gtreve();
341 void type_of_call gtreve_root();
343 void type_of_call grndm(Float_t *, const Int_t &);
345 void type_of_call grndmq(Int_t &, Int_t &, const Int_t &,
348 void type_of_call gdtom(Float_t *, Float_t *, Int_t &);
350 void type_of_call glmoth(DEFCHARD, Int_t &, Int_t &, Int_t *,
351 Int_t *, Int_t * DEFCHARL);
353 void type_of_call gmedia(Float_t *, Int_t &);
355 void type_of_call gmtod(Float_t *, Float_t *, Int_t &);
357 void type_of_call gsrotm(const Int_t &, const Float_t &, const Float_t &,
358 const Float_t &, const Float_t &, const Float_t &,
361 void type_of_call gprotm(const Int_t &);
363 void type_of_call grfile(const Int_t&, DEFCHARD,
364 DEFCHARD DEFCHARL DEFCHARL);
366 void type_of_call gfmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
367 Float_t &, Float_t &, Float_t &, Float_t *,
370 void type_of_call gfpart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
371 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
373 void type_of_call gftmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
374 Float_t &, Float_t &, Float_t &, Float_t &,
375 Float_t &, Float_t &, Float_t *, Int_t * DEFCHARL);
377 void type_of_call gftmat(const Int_t&, const Int_t&, DEFCHARD, const Int_t&,
379 ,Float_t *, Int_t & DEFCHARL);
381 void type_of_call gsmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
382 Float_t &, Float_t &, Float_t &, Float_t *,
385 void type_of_call gsmixt(const Int_t&, DEFCHARD, Float_t *, Float_t *,
386 Float_t &, Int_t &, Float_t * DEFCHARL);
388 void type_of_call gspart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
389 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
392 void type_of_call gstmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
393 Float_t &, Float_t &, Float_t &, Float_t &,
394 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
396 void type_of_call gsckov(Int_t &itmed, Int_t &npckov, Float_t *ppckov,
397 Float_t *absco, Float_t *effic, Float_t *rindex);
398 void type_of_call gstpar(const Int_t&, DEFCHARD, Float_t & DEFCHARL);
400 void type_of_call gsdvn(DEFCHARD,DEFCHARD, Int_t &, Int_t &
403 void type_of_call gsdvn2(DEFCHARD,DEFCHARD, Int_t &, Int_t &, Float_t &,
404 Int_t & DEFCHARL DEFCHARL);
406 void type_of_call gsdvs(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &
409 void type_of_call gsdvs2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t &,
410 Int_t & DEFCHARL DEFCHARL);
412 void type_of_call gsdvt(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &,
413 Int_t & DEFCHARL DEFCHARL);
415 void type_of_call gsdvt2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t&,
416 Int_t &, Int_t & DEFCHARL DEFCHARL);
418 void type_of_call gsord(DEFCHARD, Int_t & DEFCHARL);
420 void type_of_call gspos(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
421 Float_t &, Int_t &, DEFCHARD DEFCHARL DEFCHARL
424 void type_of_call gsposp(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
425 Float_t &, Int_t &, DEFCHARD,
426 Float_t *, Int_t & DEFCHARL DEFCHARL DEFCHARL);
428 void type_of_call gsvolu(DEFCHARD, DEFCHARD, Int_t &, Float_t *, Int_t &,
429 Int_t & DEFCHARL DEFCHARL);
431 void type_of_call gsatt(DEFCHARD, DEFCHARD, Int_t & DEFCHARL DEFCHARL);
433 void type_of_call gfpara(DEFCHARD , Int_t&, Int_t&, Int_t&, Int_t&, Float_t*,
436 void type_of_call gckpar(Int_t&, Int_t&, Float_t*);
438 void type_of_call gckmat(Int_t&, DEFCHARD DEFCHARL);
440 void type_of_call gprint(DEFCHARD,const int& DEFCHARL);
442 void type_of_call gdinit();
444 void type_of_call gdopt(DEFCHARD,DEFCHARD DEFCHARL DEFCHARL);
446 void type_of_call gdraw(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
447 Float_t &, Float_t &, Float_t & DEFCHARL);
448 void type_of_call gdrayt(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
449 Float_t &, Float_t &, Float_t & DEFCHARL);
450 void type_of_call gdrawc(DEFCHARD,Int_t &, Float_t &, Float_t &, Float_t &,
451 Float_t &, Float_t & DEFCHARL);
452 void type_of_call gdrawx(DEFCHARD,Float_t &, Float_t &, Float_t &, Float_t &,
453 Float_t &, Float_t &, Float_t &, Float_t &,
455 void type_of_call gdhead(Int_t &,DEFCHARD, Float_t & DEFCHARL);
456 void type_of_call gdxyz(Int_t &);
457 void type_of_call gdcxyz();
458 void type_of_call gdman(Float_t &, Float_t &);
459 void type_of_call gdwmn1(Float_t &, Float_t &);
460 void type_of_call gdwmn2(Float_t &, Float_t &);
461 void type_of_call gdwmn3(Float_t &, Float_t &);
462 void type_of_call gdspec(DEFCHARD DEFCHARL);
463 void type_of_call gdfspc(DEFCHARD, Int_t &, Int_t & DEFCHARL) {;}
464 void type_of_call gdtree(DEFCHARD, Int_t &, Int_t & DEFCHARL);
466 void type_of_call gdopen(Int_t &);
467 void type_of_call gdclos();
468 void type_of_call gdelet(Int_t &);
469 void type_of_call gdshow(Int_t &);
470 void type_of_call geditv(Int_t &) {;}
473 void type_of_call dzshow(DEFCHARD,const int&,const int&,DEFCHARD,const int&,
474 const int&, const int&, const int& DEFCHARL
477 void type_of_call mzdrop(Int_t&, Int_t&, DEFCHARD DEFCHARL);
479 void type_of_call setbomb(Float_t &);
480 void type_of_call setclip(DEFCHARD, Float_t &,Float_t &,Float_t &,Float_t &,
481 Float_t &, Float_t & DEFCHARL);
482 void type_of_call gcomad(DEFCHARD, Int_t*& DEFCHARL);
484 void type_of_call ertrak(const Float_t *const x1, const Float_t *const p1,
485 const Float_t *x2, const Float_t *p2,
486 const Int_t &ipa, DEFCHARD DEFCHARL);
488 void type_of_call ertrgo();
490 float type_of_call gbrelm(const Float_t &z, const Float_t& t, const Float_t& cut);
491 float type_of_call gprelm(const Float_t &z, const Float_t& t, const Float_t& cut);
495 // Geant3 global pointer
497 static Int_t defSize = 600;
501 //____________________________________________________________________________
505 // Default constructor
509 //____________________________________________________________________________
510 TGeant3::TGeant3(const char *title, Int_t nwgeant)
511 :AliMC("TGeant3",title)
514 // Standard constructor for TGeant3 with ZEBRA initialisation
525 // Load Address of Geant3 commons
528 // Zero number of particles
532 //____________________________________________________________________________
533 Int_t TGeant3::CurrentMaterial(Float_t &a, Float_t &z, Float_t &dens,
534 Float_t &radl, Float_t &absl) const
537 // Return the parameters of the current material during transport
541 dens = fGcmate->dens;
542 radl = fGcmate->radl;
543 absl = fGcmate->absl;
544 return 1; //this could be the number of elements in mixture
547 //____________________________________________________________________________
548 void TGeant3::DefaultRange()
551 // Set range of current drawing pad to 20x20 cm
557 higz->Range(0,0,20,20);
560 //____________________________________________________________________________
561 void TGeant3::InitHIGZ()
572 //____________________________________________________________________________
573 void TGeant3::LoadAddress()
576 // Assigns the address of the GEANT common blocks to the structures
577 // that allow their access from C++
580 gcomad(PASSCHARD("QUEST"), (int*&) fQuest PASSCHARL("QUEST"));
581 gcomad(PASSCHARD("GCBANK"),(int*&) fGcbank PASSCHARL("GCBANK"));
582 gcomad(PASSCHARD("GCLINK"),(int*&) fGclink PASSCHARL("GCLINK"));
583 gcomad(PASSCHARD("GCCUTS"),(int*&) fGccuts PASSCHARL("GCCUTS"));
584 gcomad(PASSCHARD("GCMULO"),(int*&) fGcmulo PASSCHARL("GCMULO"));
585 gcomad(PASSCHARD("GCFLAG"),(int*&) fGcflag PASSCHARL("GCFLAG"));
586 gcomad(PASSCHARD("GCKINE"),(int*&) fGckine PASSCHARL("GCKINE"));
587 gcomad(PASSCHARD("GCKING"),(int*&) fGcking PASSCHARL("GCKING"));
588 gcomad(PASSCHARD("GCKIN2"),(int*&) fGckin2 PASSCHARL("GCKIN2"));
589 gcomad(PASSCHARD("GCKIN3"),(int*&) fGckin3 PASSCHARL("GCKIN3"));
590 gcomad(PASSCHARD("GCMATE"),(int*&) fGcmate PASSCHARL("GCMATE"));
591 gcomad(PASSCHARD("GCTMED"),(int*&) fGctmed PASSCHARL("GCTMED"));
592 gcomad(PASSCHARD("GCTRAK"),(int*&) fGctrak PASSCHARL("GCTRAK"));
593 gcomad(PASSCHARD("GCTPOL"),(int*&) fGctpol PASSCHARL("GCTPOL"));
594 gcomad(PASSCHARD("GCVOLU"),(int*&) fGcvolu PASSCHARL("GCVOLU"));
595 gcomad(PASSCHARD("GCNUM"), (int*&) fGcnum PASSCHARL("GCNUM"));
596 gcomad(PASSCHARD("GCSETS"),(int*&) fGcsets PASSCHARL("GCSETS"));
597 gcomad(PASSCHARD("GCPHYS"),(int*&) fGcphys PASSCHARL("GCPHYS"));
598 gcomad(PASSCHARD("GCOPTI"),(int*&) fGcopti PASSCHARL("GCOPTI"));
599 gcomad(PASSCHARD("GCTLIT"),(int*&) fGctlit PASSCHARL("GCTLIT"));
600 gcomad(PASSCHARD("GCVDMA"),(int*&) fGcvdma PASSCHARL("GCVDMA"));
603 gcomad(PASSCHARD("ERTRIO"),(int*&) fErtrio PASSCHARL("ERTRIO"));
604 gcomad(PASSCHARD("EROPTS"),(int*&) fEropts PASSCHARL("EROPTS"));
605 gcomad(PASSCHARD("EROPTC"),(int*&) fEroptc PASSCHARL("EROPTC"));
606 gcomad(PASSCHARD("ERWORK"),(int*&) fErwork PASSCHARL("ERWORK"));
608 // Variables for ZEBRA store
609 gcomad(PASSCHARD("IQ"), addr PASSCHARL("IQ"));
611 gcomad(PASSCHARD("LQ"), addr PASSCHARL("LQ"));
616 //_____________________________________________________________________________
617 void TGeant3::GeomIter()
620 // Geometry iterator for moving upward in the geometry tree
621 // Initialise the iterator
623 fNextVol=fGcvolu->nlevel;
626 //____________________________________________________________________________
627 void TGeant3::FinishGeometry()
629 //Close the geometry structure
633 //____________________________________________________________________________
634 Int_t TGeant3::NextVolUp(Text_t *name, Int_t ©)
637 // Geometry iterator for moving upward in the geometry tree
638 // Return next volume up
643 gname=fGcvolu->names[fNextVol];
644 copy=fGcvolu->number[fNextVol];
645 i=fGcvolu->lvolum[fNextVol];
646 name = fVolNames[i-1];
647 if(gname == fZiq[fGclink->jvolum+i]) return i;
648 else printf("GeomTree: Volume %s not found in bank\n",name);
653 //_____________________________________________________________________________
654 void TGeant3::BuildPhysics()
659 //_____________________________________________________________________________
660 Int_t TGeant3::CurrentVolID(Int_t ©) const
663 // Returns the current volume ID and copy number
666 if( (i=fGcvolu->nlevel-1) < 0 ) {
667 Warning("CurrentVolID","Stack depth only %d\n",fGcvolu->nlevel);
669 gname=fGcvolu->names[i];
670 copy=fGcvolu->number[i];
671 i=fGcvolu->lvolum[i];
672 if(gname == fZiq[fGclink->jvolum+i]) return i;
673 else Warning("CurrentVolID","Volume %4s not found\n",(char*)&gname);
678 //_____________________________________________________________________________
679 Int_t TGeant3::CurrentVolOffID(Int_t off, Int_t ©) const
682 // Return the current volume "off" upward in the geometrical tree
683 // ID and copy number
686 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
687 Warning("CurrentVolOffID","Offset requested %d but stack depth %d\n",
688 off,fGcvolu->nlevel);
690 gname=fGcvolu->names[i];
691 copy=fGcvolu->number[i];
692 i=fGcvolu->lvolum[i];
693 if(gname == fZiq[fGclink->jvolum+i]) return i;
694 else Warning("CurrentVolOffID","Volume %4s not found\n",(char*)&gname);
699 //_____________________________________________________________________________
700 const char* TGeant3::CurrentVolName() const
703 // Returns the current volume name
706 if( (i=fGcvolu->nlevel-1) < 0 ) {
707 Warning("CurrentVolName","Stack depth %d\n",fGcvolu->nlevel);
709 gname=fGcvolu->names[i];
710 i=fGcvolu->lvolum[i];
711 if(gname == fZiq[fGclink->jvolum+i]) return fVolNames[i-1];
712 else Warning("CurrentVolName","Volume %4s not found\n",(char*) &gname);
717 //_____________________________________________________________________________
718 const char* TGeant3::CurrentVolOffName(Int_t off) const
721 // Return the current volume "off" upward in the geometrical tree
722 // ID, name and copy number
723 // if name=0 no name is returned
726 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
727 Warning("CurrentVolOffName",
728 "Offset requested %d but stack depth %d\n",off,fGcvolu->nlevel);
730 gname=fGcvolu->names[i];
731 i=fGcvolu->lvolum[i];
732 if(gname == fZiq[fGclink->jvolum+i]) return fVolNames[i-1];
733 else Warning("CurrentVolOffName","Volume %4s not found\n",(char*)&gname);
738 //_____________________________________________________________________________
739 Int_t TGeant3::IdFromPDG(Int_t pdg) const
742 // Return Geant3 code from PDG and pseudo ENDF code
744 for(Int_t i=0;i<fNPDGCodes;++i)
745 if(pdg==fPDGCode[i]) return i;
749 //_____________________________________________________________________________
750 Int_t TGeant3::PDGFromId(Int_t id) const
752 if(id>0 && id<fNPDGCodes) return fPDGCode[id];
756 //_____________________________________________________________________________
757 void TGeant3::DefineParticles()
760 // Define standard Geant 3 particles
763 // Load standard numbers for GEANT particles and PDG conversion
764 fPDGCode[fNPDGCodes++]=-99; // 0 = unused location
765 fPDGCode[fNPDGCodes++]=22; // 1 = photon
766 fPDGCode[fNPDGCodes++]=-11; // 2 = positron
767 fPDGCode[fNPDGCodes++]=11; // 3 = electron
768 fPDGCode[fNPDGCodes++]=12; // 4 = neutrino e
769 fPDGCode[fNPDGCodes++]=-13; // 5 = muon +
770 fPDGCode[fNPDGCodes++]=13; // 6 = muon -
771 fPDGCode[fNPDGCodes++]=111; // 7 = pi0
772 fPDGCode[fNPDGCodes++]=211; // 8 = pi+
773 fPDGCode[fNPDGCodes++]=-211; // 9 = pi-
774 fPDGCode[fNPDGCodes++]=130; // 10 = Kaon Long
775 fPDGCode[fNPDGCodes++]=321; // 11 = Kaon +
776 fPDGCode[fNPDGCodes++]=-321; // 12 = Kaon -
777 fPDGCode[fNPDGCodes++]=2112; // 13 = Neutron
778 fPDGCode[fNPDGCodes++]=2212; // 14 = Proton
779 fPDGCode[fNPDGCodes++]=-2212; // 15 = Anti Proton
780 fPDGCode[fNPDGCodes++]=310; // 16 = Kaon Short
781 fPDGCode[fNPDGCodes++]=221; // 17 = Eta
782 fPDGCode[fNPDGCodes++]=3122; // 18 = Lambda
783 fPDGCode[fNPDGCodes++]=3222; // 19 = Sigma +
784 fPDGCode[fNPDGCodes++]=3212; // 20 = Sigma 0
785 fPDGCode[fNPDGCodes++]=3112; // 21 = Sigma -
786 fPDGCode[fNPDGCodes++]=3322; // 22 = Xi0
787 fPDGCode[fNPDGCodes++]=3312; // 23 = Xi-
788 fPDGCode[fNPDGCodes++]=3334; // 24 = Omega-
789 fPDGCode[fNPDGCodes++]=-2112; // 25 = Anti Proton
790 fPDGCode[fNPDGCodes++]=-3122; // 26 = Anti Proton
791 fPDGCode[fNPDGCodes++]=-3222; // 27 = Anti Sigma -
792 fPDGCode[fNPDGCodes++]=-3212; // 28 = Anti Sigma 0
793 fPDGCode[fNPDGCodes++]=-3112; // 29 = Anti Sigma 0
794 fPDGCode[fNPDGCodes++]=-3322; // 30 = Anti Xi 0
795 fPDGCode[fNPDGCodes++]=-3312; // 31 = Anti Xi +
796 fPDGCode[fNPDGCodes++]=-3334; // 32 = Anti Omega +
803 /* --- Define additional particles */
804 Gspart(33, "OMEGA(782)", 3, 0.782, 0., 7.836e-23);
805 fPDGCode[fNPDGCodes++]=223; // 33 = Omega(782)
807 Gspart(34, "PHI(1020)", 3, 1.019, 0., 1.486e-22);
808 fPDGCode[fNPDGCodes++]=333; // 34 = PHI (1020)
810 Gspart(35, "D +", 4, 1.87, 1., 1.066e-12);
811 fPDGCode[fNPDGCodes++]=411; // 35 = D+
813 Gspart(36, "D -", 4, 1.87, -1., 1.066e-12);
814 fPDGCode[fNPDGCodes++]=-411; // 36 = D-
816 Gspart(37, "D 0", 3, 1.865, 0., 4.2e-13);
817 fPDGCode[fNPDGCodes++]=421; // 37 = D0
819 Gspart(38, "ANTI D 0", 3, 1.865, 0., 4.2e-13);
820 fPDGCode[fNPDGCodes++]=-421; // 38 = D0 bar
822 fPDGCode[fNPDGCodes++]=-99; // 39 = unassigned
824 fPDGCode[fNPDGCodes++]=-99; // 40 = unassigned
826 fPDGCode[fNPDGCodes++]=-99; // 41 = unassigned
828 Gspart(42, "RHO +", 4, 0.768, 1., 4.353e-24);
829 fPDGCode[fNPDGCodes++]=213; // 42 = RHO+
831 Gspart(43, "RHO -", 4, 0.768, -1., 4.353e-24);
832 fPDGCode[fNPDGCodes++]=-213; // 40 = RHO-
834 Gspart(44, "RHO 0", 3, 0.768, 0., 4.353e-24);
835 fPDGCode[fNPDGCodes++]=113; // 37 = D0
838 // Use ENDF-6 mapping for ions = 10000*z+10*a+iso
840 // and numbers above 5 000 000 for special applications
843 const Int_t kion=10000000;
845 const Int_t kspe=50000000;
847 TDatabasePDG *pdgDB = TDatabasePDG::Instance();
849 const Double_t autogev=0.9314943228;
850 const Double_t hslash = 1.0545726663e-27;
851 const Double_t erggev = 1/1.6021773349e-3;
852 const Double_t hshgev = hslash*erggev;
853 const Double_t yearstosec = 3600*24*365.25;
856 pdgDB->AddParticle("Deuteron","Deuteron",2*autogev+8.071e-3,kTRUE,
857 0,1,"Ion",kion+10020);
858 fPDGCode[fNPDGCodes++]=kion+10020; // 45 = Deuteron
860 pdgDB->AddParticle("Triton","Triton",3*autogev+14.931e-3,kFALSE,
861 hshgev/(12.33*yearstosec),1,"Ion",kion+10030);
862 fPDGCode[fNPDGCodes++]=kion+10030; // 46 = Triton
864 pdgDB->AddParticle("Alpha","Alpha",4*autogev+2.424e-3,kTRUE,
865 hshgev/(12.33*yearstosec),2,"Ion",kion+20040);
866 fPDGCode[fNPDGCodes++]=kion+20040; // 47 = Alpha
868 fPDGCode[fNPDGCodes++]=0; // 48 = geantino mapped to rootino
870 pdgDB->AddParticle("HE3","HE3",3*autogev+14.931e-3,kFALSE,
871 0,2,"Ion",kion+20030);
872 fPDGCode[fNPDGCodes++]=kion+20030; // 49 = HE3
874 pdgDB->AddParticle("Cherenkov","Cherenkov",0,kFALSE,
875 0,0,"Special",kspe+50);
876 fPDGCode[fNPDGCodes++]=kspe+50; // 50 = Cherenkov
878 Gspart(51, "FeedbackPhoton", 7, 0., 0.,1.e20 );
879 pdgDB->AddParticle("FeedbackPhoton","FeedbackPhoton",0,kFALSE,
880 0,0,"Special",kspe+51);
881 fPDGCode[fNPDGCodes++]=kspe+51; // 51 = FeedbackPhoton
883 /* --- Define additional decay modes --- */
884 /* --- omega(783) --- */
885 for (kz = 0; kz < 6; ++kz) {
896 Gsdk(ipa, bratio, mode);
897 /* --- phi(1020) --- */
898 for (kz = 0; kz < 6; ++kz) {
913 Gsdk(ipa, bratio, mode);
915 for (kz = 0; kz < 6; ++kz) {
928 Gsdk(ipa, bratio, mode);
930 for (kz = 0; kz < 6; ++kz) {
943 Gsdk(ipa, bratio, mode);
945 for (kz = 0; kz < 6; ++kz) {
956 Gsdk(ipa, bratio, mode);
957 /* --- Anti D0 --- */
958 for (kz = 0; kz < 6; ++kz) {
969 Gsdk(ipa, bratio, mode);
971 for (kz = 0; kz < 6; ++kz) {
978 Gsdk(ipa, bratio, mode);
980 for (kz = 0; kz < 6; ++kz) {
987 Gsdk(ipa, bratio, mode);
989 for (kz = 0; kz < 6; ++kz) {
996 Gsdk(ipa, bratio, mode);
999 for (kz = 0; kz < 6; ++kz) {
1008 Gsdk(ipa, bratio, mode);
1011 Gsdk(ipa, bratio, mode);
1014 Gsdk(ipa, bratio, mode);
1019 //_____________________________________________________________________________
1020 Int_t TGeant3::VolId(const Text_t *name) const
1023 // Return the unique numeric identifier for volume name
1026 strncpy((char *) &gname, name, 4);
1027 for(i=1; i<=fGcnum->nvolum; i++)
1028 if(gname == fZiq[fGclink->jvolum+i]) return i;
1029 printf("VolId: Volume %s not found\n",name);
1033 //_____________________________________________________________________________
1034 Int_t TGeant3::NofVolumes() const
1037 // Return total number of volumes in the geometry
1039 return fGcnum->nvolum;
1042 //_____________________________________________________________________________
1043 const char* TGeant3::VolName(Int_t id) const
1046 // Return the volume name given the volume identifier
1048 const char name[5]="NULL";
1049 if(id<1 || id > fGcnum->nvolum || fGclink->jvolum<=0)
1052 return fVolNames[id-1];
1055 //_____________________________________________________________________________
1056 void TGeant3::SetCut(const char* cutName, Float_t cutValue)
1058 if(!strcmp(cutName,"CUTGAM"))
1059 fGccuts->cutgam=cutValue;
1060 else if(!strcmp(cutName,"CUTGAM"))
1061 fGccuts->cutele=cutValue;
1062 else if(!strcmp(cutName,"CUTELE"))
1063 fGccuts->cutneu=cutValue;
1064 else if(!strcmp(cutName,"CUTHAD"))
1065 fGccuts->cuthad=cutValue;
1066 else if(!strcmp(cutName,"CUTMUO"))
1067 fGccuts->cutmuo=cutValue;
1068 else if(!strcmp(cutName,"BCUTE"))
1069 fGccuts->bcute=cutValue;
1070 else if(!strcmp(cutName,"BCUTM"))
1071 fGccuts->bcutm=cutValue;
1072 else if(!strcmp(cutName,"DCUTE"))
1073 fGccuts->dcute=cutValue;
1074 else if(!strcmp(cutName,"DCUTM"))
1075 fGccuts->dcutm=cutValue;
1076 else if(!strcmp(cutName,"PPCUTM"))
1077 fGccuts->ppcutm=cutValue;
1078 else if(!strcmp(cutName,"TOFMAX"))
1079 fGccuts->tofmax=cutValue;
1080 else Warning("SetCut","Cut %s not implemented\n",cutName);
1083 //_____________________________________________________________________________
1084 void TGeant3::SetProcess(const char* flagName, Int_t flagValue)
1086 if(!strcmp(flagName,"PAIR"))
1087 fGcphys->ipair=flagValue;
1088 else if(!strcmp(flagName,"COMP"))
1089 fGcphys->icomp=flagValue;
1090 else if(!strcmp(flagName,"PHOT"))
1091 fGcphys->iphot=flagValue;
1092 else if(!strcmp(flagName,"PFIS"))
1093 fGcphys->ipfis=flagValue;
1094 else if(!strcmp(flagName,"DRAY"))
1095 fGcphys->idray=flagValue;
1096 else if(!strcmp(flagName,"ANNI"))
1097 fGcphys->ianni=flagValue;
1098 else if(!strcmp(flagName,"BREM"))
1099 fGcphys->ibrem=flagValue;
1100 else if(!strcmp(flagName,"HADR"))
1101 fGcphys->ihadr=flagValue;
1102 else if(!strcmp(flagName,"MUNU"))
1103 fGcphys->imunu=flagValue;
1104 else if(!strcmp(flagName,"DCAY"))
1105 fGcphys->idcay=flagValue;
1106 else if(!strcmp(flagName,"LOSS"))
1107 fGcphys->iloss=flagValue;
1108 else if(!strcmp(flagName,"MULS"))
1109 fGcphys->imuls=flagValue;
1110 else if(!strcmp(flagName,"RAYL"))
1111 fGcphys->irayl=flagValue;
1112 else Warning("SetFlag","Flag %s not implemented\n",flagName);
1115 //_____________________________________________________________________________
1116 Float_t TGeant3::Xsec(char* reac, Float_t energy, Int_t part, Int_t mate)
1118 Int_t gpart = IdFromPDG(part);
1119 if(!strcmp(reac,"PHOT"))
1122 Error("Xsec","Can calculate photoelectric only for photons\n");
1128 //_____________________________________________________________________________
1129 void TGeant3::TrackPosition(TLorentzVector &xyz) const
1132 // Return the current position in the master reference frame of the
1133 // track being transported
1135 xyz[0]=fGctrak->vect[0];
1136 xyz[1]=fGctrak->vect[1];
1137 xyz[2]=fGctrak->vect[2];
1138 xyz[3]=fGctrak->tofg;
1141 //_____________________________________________________________________________
1142 Float_t TGeant3::TrackTime() const
1145 // Return the current time of flight of the track being transported
1147 return fGctrak->tofg;
1150 //_____________________________________________________________________________
1151 void TGeant3::TrackMomentum(TLorentzVector &xyz) const
1154 // Return the direction and the momentum (GeV/c) of the track
1155 // currently being transported
1157 Double_t ptot=fGctrak->vect[6];
1158 xyz[0]=fGctrak->vect[3]*ptot;
1159 xyz[1]=fGctrak->vect[4]*ptot;
1160 xyz[2]=fGctrak->vect[5]*ptot;
1161 xyz[3]=fGctrak->getot;
1164 //_____________________________________________________________________________
1165 Float_t TGeant3::TrackCharge() const
1168 // Return charge of the track currently transported
1170 return fGckine->charge;
1173 //_____________________________________________________________________________
1174 Float_t TGeant3::TrackMass() const
1177 // Return the mass of the track currently transported
1179 return fGckine->amass;
1182 //_____________________________________________________________________________
1183 Int_t TGeant3::TrackPid() const
1186 // Return the id of the particle transported
1188 return PDGFromId(fGckine->ipart);
1191 //_____________________________________________________________________________
1192 Float_t TGeant3::TrackStep() const
1195 // Return the length in centimeters of the current step
1197 return fGctrak->step;
1200 //_____________________________________________________________________________
1201 Float_t TGeant3::TrackLength() const
1204 // Return the length of the current track from its origin
1206 return fGctrak->sleng;
1209 //_____________________________________________________________________________
1210 Bool_t TGeant3::IsTrackInside() const
1213 // True if the track is not at the boundary of the current volume
1215 return (fGctrak->inwvol==0);
1218 //_____________________________________________________________________________
1219 Bool_t TGeant3::IsTrackEntering() const
1222 // True if this is the first step of the track in the current volume
1224 return (fGctrak->inwvol==1);
1227 //_____________________________________________________________________________
1228 Bool_t TGeant3::IsTrackExiting() const
1231 // True if this is the last step of the track in the current volume
1233 return (fGctrak->inwvol==2);
1236 //_____________________________________________________________________________
1237 Bool_t TGeant3::IsTrackOut() const
1240 // True if the track is out of the setup
1242 return (fGctrak->inwvol==3);
1245 //_____________________________________________________________________________
1246 Bool_t TGeant3::IsTrackStop() const
1249 // True if the track energy has fallen below the threshold
1251 return (fGctrak->istop==2);
1254 //_____________________________________________________________________________
1255 Int_t TGeant3::NSecondaries() const
1258 // Number of secondary particles generated in the current step
1260 return fGcking->ngkine;
1263 //_____________________________________________________________________________
1264 Int_t TGeant3::CurrentEvent() const
1267 // Number of the current event
1269 return fGcflag->idevt;
1272 //_____________________________________________________________________________
1273 const char* TGeant3::ProdProcess() const
1276 // Name of the process that has produced the secondary particles
1277 // in the current step
1279 static char proc[5];
1280 const Int_t ipmec[13] = { 5,6,7,8,9,10,11,12,21,23,25,105,108 };
1283 if(fGcking->ngkine>0) {
1284 for (km = 0; km < fGctrak->nmec; ++km) {
1285 for (im = 0; im < 13; ++im) {
1286 if (fGctrak->lmec[km] == ipmec[im]) {
1287 mec = fGctrak->lmec[km];
1288 if (0 < mec && mec < 31) {
1289 strncpy(proc,(char *)&fGctrak->namec[mec - 1],4);
1290 } else if (mec - 100 <= 30 && mec - 100 > 0) {
1291 strncpy(proc,(char *)&fGctpol->namec1[mec - 101],4);
1298 strcpy(proc,"UNKN");
1299 } else strcpy(proc,"NONE");
1303 //_____________________________________________________________________________
1304 void TGeant3::GetSecondary(Int_t isec, Int_t& ipart,
1305 TLorentzVector &x, TLorentzVector &p)
1308 // Get the parameters of the secondary track number isec produced
1309 // in the current step
1312 if(-1<isec && isec<fGcking->ngkine) {
1313 ipart=Int_t (fGcking->gkin[isec][4] +0.5);
1315 x[i]=fGckin3->gpos[isec][i];
1316 p[i]=fGcking->gkin[isec][i];
1318 x[3]=fGcking->tofd[isec];
1319 p[3]=fGcking->gkin[isec][3];
1321 printf(" * TGeant3::GetSecondary * Secondary %d does not exist\n",isec);
1322 x[0]=x[1]=x[2]=x[3]=p[0]=p[1]=p[2]=p[3]=0;
1327 //_____________________________________________________________________________
1328 void TGeant3::InitLego()
1331 SetDEBU(0,0,0); //do not print a message
1334 //_____________________________________________________________________________
1335 Bool_t TGeant3::IsTrackDisappeared() const
1338 // True if the current particle has disappered
1339 // either because it decayed or because it underwent
1340 // an inelastic collision
1342 return (fGctrak->istop==1);
1345 //_____________________________________________________________________________
1346 Bool_t TGeant3::IsTrackAlive() const
1349 // True if the current particle is alive and will continue to be
1352 return (fGctrak->istop==0);
1355 //_____________________________________________________________________________
1356 void TGeant3::StopTrack()
1359 // Stop the transport of the current particle and skip to the next
1364 //_____________________________________________________________________________
1365 void TGeant3::StopEvent()
1368 // Stop simulation of the current event and skip to the next
1373 //_____________________________________________________________________________
1374 Float_t TGeant3::MaxStep() const
1377 // Return the maximum step length in the current medium
1379 return fGctmed->stemax;
1382 //_____________________________________________________________________________
1383 void TGeant3::SetMaxStep(Float_t maxstep)
1386 // Set the maximum step allowed till the particle is in the current medium
1388 fGctmed->stemax=maxstep;
1391 //_____________________________________________________________________________
1392 void TGeant3::SetMaxNStep(Int_t maxnstp)
1395 // Set the maximum number of steps till the particle is in the current medium
1397 fGctrak->maxnst=maxnstp;
1400 //_____________________________________________________________________________
1401 Int_t TGeant3::GetMaxNStep() const
1404 // Maximum number of steps allowed in current medium
1406 return fGctrak->maxnst;
1409 //_____________________________________________________________________________
1410 void TGeant3::Material(Int_t& kmat, const char* name, Float_t a, Float_t z,
1411 Float_t dens, Float_t radl, Float_t absl, Float_t* buf,
1415 // Defines a Material
1417 // kmat number assigned to the material
1418 // name material name
1419 // a atomic mass in au
1421 // dens density in g/cm3
1422 // absl absorbtion length in cm
1423 // if >=0 it is ignored and the program
1424 // calculates it, if <0. -absl is taken
1425 // radl radiation length in cm
1426 // if >=0 it is ignored and the program
1427 // calculates it, if <0. -radl is taken
1428 // buf pointer to an array of user words
1429 // nbuf number of user words
1431 Int_t jmate=fGclink->jmate;
1437 for(i=1; i<=ns; i++) {
1438 if(fZlq[jmate-i]==0) {
1444 gsmate(kmat,PASSCHARD(name), a, z, dens, radl, absl, buf,
1445 nwbuf PASSCHARL(name));
1448 //_____________________________________________________________________________
1449 void TGeant3::Mixture(Int_t& kmat, const char* name, Float_t* a, Float_t* z,
1450 Float_t dens, Int_t nlmat, Float_t* wmat)
1453 // Defines mixture OR COMPOUND IMAT as composed by
1454 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
1456 // If NLMAT > 0 then wmat contains the proportion by
1457 // weights of each basic material in the mixture.
1459 // If nlmat < 0 then WMAT contains the number of atoms
1460 // of a given kind into the molecule of the COMPOUND
1461 // In this case, WMAT in output is changed to relative
1464 Int_t jmate=fGclink->jmate;
1470 for(i=1; i<=ns; i++) {
1471 if(fZlq[jmate-i]==0) {
1477 gsmixt(kmat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
1480 //_____________________________________________________________________________
1481 void TGeant3::Medium(Int_t& kmed, const char* name, Int_t nmat, Int_t isvol,
1482 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
1483 Float_t stemax, Float_t deemax, Float_t epsil,
1484 Float_t stmin, Float_t* ubuf, Int_t nbuf)
1487 // kmed tracking medium number assigned
1488 // name tracking medium name
1489 // nmat material number
1490 // isvol sensitive volume flag
1491 // ifield magnetic field
1492 // fieldm max. field value (kilogauss)
1493 // tmaxfd max. angle due to field (deg/step)
1494 // stemax max. step allowed
1495 // deemax max. fraction of energy lost in a step
1496 // epsil tracking precision (cm)
1497 // stmin min. step due to continuos processes (cm)
1499 // ifield = 0 if no magnetic field; ifield = -1 if user decision in guswim;
1500 // ifield = 1 if tracking performed with grkuta; ifield = 2 if tracking
1501 // performed with ghelix; ifield = 3 if tracking performed with ghelx3.
1503 Int_t jtmed=fGclink->jtmed;
1509 for(i=1; i<=ns; i++) {
1510 if(fZlq[jtmed-i]==0) {
1516 gstmed(kmed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1517 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1520 //_____________________________________________________________________________
1521 void TGeant3::Matrix(Int_t& krot, Float_t thex, Float_t phix, Float_t they,
1522 Float_t phiy, Float_t thez, Float_t phiz)
1525 // krot rotation matrix number assigned
1526 // theta1 polar angle for axis i
1527 // phi1 azimuthal angle for axis i
1528 // theta2 polar angle for axis ii
1529 // phi2 azimuthal angle for axis ii
1530 // theta3 polar angle for axis iii
1531 // phi3 azimuthal angle for axis iii
1533 // it defines the rotation matrix number irot.
1535 Int_t jrotm=fGclink->jrotm;
1541 for(i=1; i<=ns; i++) {
1542 if(fZlq[jrotm-i]==0) {
1548 gsrotm(krot, thex, phix, they, phiy, thez, phiz);
1551 //_____________________________________________________________________________
1552 Int_t TGeant3::GetMedium() const
1555 // Return the number of the current medium
1557 return fGctmed->numed;
1560 //_____________________________________________________________________________
1561 Float_t TGeant3::Edep() const
1564 // Return the energy lost in the current step
1566 return fGctrak->destep;
1569 //_____________________________________________________________________________
1570 Float_t TGeant3::Etot() const
1573 // Return the total energy of the current track
1575 return fGctrak->getot;
1578 //_____________________________________________________________________________
1579 void TGeant3::Rndm(Float_t* r, const Int_t n) const
1582 // Return an array of n random numbers uniformly distributed
1583 // between 0 and 1 not included
1588 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1590 // Functions from GBASE
1592 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1594 //____________________________________________________________________________
1595 void TGeant3::Gfile(const char *filename, const char *option)
1598 // Routine to open a GEANT/RZ data base.
1600 // LUN logical unit number associated to the file
1602 // CHFILE RZ file name
1604 // CHOPT is a character string which may be
1605 // N To create a new file
1606 // U to open an existing file for update
1607 // " " to open an existing file for read only
1608 // Q The initial allocation (default 1000 records)
1609 // is given in IQUEST(10)
1610 // X Open the file in exchange format
1611 // I Read all data structures from file to memory
1612 // O Write all data structures from memory to file
1615 // If options "I" or "O" all data structures are read or
1616 // written from/to file and the file is closed.
1617 // See routine GRMDIR to create subdirectories
1618 // See routines GROUT,GRIN to write,read objects
1620 grfile(21, PASSCHARD(filename), PASSCHARD(option) PASSCHARL(filename)
1624 //____________________________________________________________________________
1625 void TGeant3::Gpcxyz()
1628 // Print track and volume parameters at current point
1633 //_____________________________________________________________________________
1634 void TGeant3::Ggclos()
1637 // Closes off the geometry setting.
1638 // Initializes the search list for the contents of each
1639 // volume following the order they have been positioned, and
1640 // inserting the content '0' when a call to GSNEXT (-1) has
1641 // been required by the user.
1642 // Performs the development of the JVOLUM structure for all
1643 // volumes with variable parameters, by calling GGDVLP.
1644 // Interprets the user calls to GSORD, through GGORD.
1645 // Computes and stores in a bank (next to JVOLUM mother bank)
1646 // the number of levels in the geometrical tree and the
1647 // maximum number of contents per level, by calling GGNLEV.
1648 // Sets status bit for CONCAVE volumes, through GGCAVE.
1649 // Completes the JSET structure with the list of volume names
1650 // which identify uniquely a given physical detector, the
1651 // list of bit numbers to pack the corresponding volume copy
1652 // numbers, and the generic path(s) in the JVOLUM tree,
1653 // through the routine GHCLOS.
1656 // Create internal list of volumes
1657 fVolNames = new char[fGcnum->nvolum][5];
1659 for(i=0; i<fGcnum->nvolum; ++i) {
1660 strncpy(fVolNames[i], (char *) &fZiq[fGclink->jvolum+i+1], 4);
1661 fVolNames[i][4]='\0';
1665 //_____________________________________________________________________________
1666 void TGeant3::Glast()
1669 // Finish a Geant run
1674 //_____________________________________________________________________________
1675 void TGeant3::Gprint(const char *name)
1678 // Routine to print data structures
1679 // CHNAME name of a data structure
1683 gprint(PASSCHARD(vname),0 PASSCHARL(vname));
1686 //_____________________________________________________________________________
1687 void TGeant3::Grun()
1690 // Steering function to process one run
1695 //_____________________________________________________________________________
1696 void TGeant3::Gtrig()
1699 // Steering function to process one event
1704 //_____________________________________________________________________________
1705 void TGeant3::Gtrigc()
1708 // Clear event partition
1713 //_____________________________________________________________________________
1714 void TGeant3::Gtrigi()
1717 // Initialises event partition
1722 //_____________________________________________________________________________
1723 void TGeant3::Gwork(Int_t nwork)
1726 // Allocates workspace in ZEBRA memory
1731 //_____________________________________________________________________________
1732 void TGeant3::Gzinit()
1735 // To initialise GEANT/ZEBRA data structures
1740 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1742 // Functions from GCONS
1744 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1746 //_____________________________________________________________________________
1747 void TGeant3::Gfmate(Int_t imat, char *name, Float_t &a, Float_t &z,
1748 Float_t &dens, Float_t &radl, Float_t &absl,
1749 Float_t* ubuf, Int_t& nbuf)
1752 // Return parameters for material IMAT
1754 gfmate(imat, PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
1758 //_____________________________________________________________________________
1759 void TGeant3::Gfpart(Int_t ipart, char *name, Int_t &itrtyp,
1760 Float_t &amass, Float_t &charge, Float_t &tlife)
1763 // Return parameters for particle of type IPART
1767 Int_t igpart = IdFromPDG(ipart);
1768 gfpart(igpart, PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
1772 //_____________________________________________________________________________
1773 void TGeant3::Gftmed(Int_t numed, char *name, Int_t &nmat, Int_t &isvol,
1774 Int_t &ifield, Float_t &fieldm, Float_t &tmaxfd,
1775 Float_t &stemax, Float_t &deemax, Float_t &epsil,
1776 Float_t &stmin, Float_t *ubuf, Int_t *nbuf)
1779 // Return parameters for tracking medium NUMED
1781 gftmed(numed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1782 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1786 void TGeant3::Gftmat(Int_t imate, Int_t ipart, char *chmeca, Int_t kdim,
1787 Float_t* tkin, Float_t* value, Float_t* pcut,
1791 // Return parameters for tracking medium NUMED
1793 gftmat(imate, ipart, PASSCHARD(chmeca), kdim,
1794 tkin, value, pcut, ixst PASSCHARL(chmeca));
1798 Float_t TGeant3::Gbrelm(Float_t z, Float_t t, Float_t bcut)
1800 return gbrelm(z,t,bcut);
1803 Float_t TGeant3::Gprelm(Float_t z, Float_t t, Float_t bcut)
1805 return gprelm(z,t,bcut);
1808 //_____________________________________________________________________________
1809 void TGeant3::Gmate()
1812 // Define standard GEANT materials
1817 //_____________________________________________________________________________
1818 void TGeant3::Gpart()
1821 // Define standard GEANT particles plus selected decay modes
1822 // and branching ratios.
1827 //_____________________________________________________________________________
1828 void TGeant3::Gsdk(Int_t ipart, Float_t *bratio, Int_t *mode)
1830 // Defines branching ratios and decay modes for standard
1832 gsdk(ipart,bratio,mode);
1835 //_____________________________________________________________________________
1836 void TGeant3::Gsmate(Int_t imat, const char *name, Float_t a, Float_t z,
1837 Float_t dens, Float_t radl, Float_t absl)
1840 // Defines a Material
1842 // kmat number assigned to the material
1843 // name material name
1844 // a atomic mass in au
1846 // dens density in g/cm3
1847 // absl absorbtion length in cm
1848 // if >=0 it is ignored and the program
1849 // calculates it, if <0. -absl is taken
1850 // radl radiation length in cm
1851 // if >=0 it is ignored and the program
1852 // calculates it, if <0. -radl is taken
1853 // buf pointer to an array of user words
1854 // nbuf number of user words
1858 gsmate(imat,PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
1862 //_____________________________________________________________________________
1863 void TGeant3::Gsmixt(Int_t imat, const char *name, Float_t *a, Float_t *z,
1864 Float_t dens, Int_t nlmat, Float_t *wmat)
1867 // Defines mixture OR COMPOUND IMAT as composed by
1868 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
1870 // If NLMAT.GT.0 then WMAT contains the PROPORTION BY
1871 // WEIGTHS OF EACH BASIC MATERIAL IN THE MIXTURE.
1873 // If NLMAT.LT.0 then WMAT contains the number of atoms
1874 // of a given kind into the molecule of the COMPOUND
1875 // In this case, WMAT in output is changed to relative
1878 gsmixt(imat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
1881 //_____________________________________________________________________________
1882 void TGeant3::Gspart(Int_t ipart, const char *name, Int_t itrtyp,
1883 Float_t amass, Float_t charge, Float_t tlife)
1886 // Store particle parameters
1888 // ipart particle code
1889 // name particle name
1890 // itrtyp transport method (see GEANT manual)
1891 // amass mass in GeV/c2
1892 // charge charge in electron units
1893 // tlife lifetime in seconds
1897 gspart(ipart,PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
1901 //_____________________________________________________________________________
1902 void TGeant3::Gstmed(Int_t numed, const char *name, Int_t nmat, Int_t isvol,
1903 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
1904 Float_t stemax, Float_t deemax, Float_t epsil,
1908 // NTMED Tracking medium number
1909 // NAME Tracking medium name
1910 // NMAT Material number
1911 // ISVOL Sensitive volume flag
1912 // IFIELD Magnetic field
1913 // FIELDM Max. field value (Kilogauss)
1914 // TMAXFD Max. angle due to field (deg/step)
1915 // STEMAX Max. step allowed
1916 // DEEMAX Max. fraction of energy lost in a step
1917 // EPSIL Tracking precision (cm)
1918 // STMIN Min. step due to continuos processes (cm)
1920 // IFIELD = 0 if no magnetic field; IFIELD = -1 if user decision in GUSWIM;
1921 // IFIELD = 1 if tracking performed with GRKUTA; IFIELD = 2 if tracking
1922 // performed with GHELIX; IFIELD = 3 if tracking performed with GHELX3.
1926 gstmed(numed,PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1927 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1930 //_____________________________________________________________________________
1931 void TGeant3::Gsckov(Int_t itmed, Int_t npckov, Float_t *ppckov,
1932 Float_t *absco, Float_t *effic, Float_t *rindex)
1935 // Stores the tables for UV photon tracking in medium ITMED
1936 // Please note that it is the user's responsability to
1937 // provide all the coefficients:
1940 // ITMED Tracking medium number
1941 // NPCKOV Number of bins of each table
1942 // PPCKOV Value of photon momentum (in GeV)
1943 // ABSCO Absorbtion coefficients
1944 // dielectric: absorbtion length in cm
1945 // metals : absorbtion fraction (0<=x<=1)
1946 // EFFIC Detection efficiency for UV photons
1947 // RINDEX Refraction index (if=0 metal)
1949 gsckov(itmed,npckov,ppckov,absco,effic,rindex);
1952 //_____________________________________________________________________________
1953 void TGeant3::Gstpar(Int_t itmed, const char *param, Float_t parval)
1956 // To change the value of cut or mechanism "CHPAR"
1957 // to a new value PARVAL for tracking medium ITMED
1958 // The data structure JTMED contains the standard tracking
1959 // parameters (CUTS and flags to control the physics processes) which
1960 // are used by default for all tracking media. It is possible to
1961 // redefine individually with GSTPAR any of these parameters for a
1962 // given tracking medium.
1963 // ITMED tracking medium number
1964 // CHPAR is a character string (variable name)
1965 // PARVAL must be given as a floating point.
1967 gstpar(itmed,PASSCHARD(param), parval PASSCHARL(param));
1970 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1972 // Functions from GCONS
1974 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1976 //_____________________________________________________________________________
1977 void TGeant3::Gfkine(Int_t itra, Float_t *vert, Float_t *pvert, Int_t &ipart,
1980 // Storing/Retrieving Vertex and Track parameters
1981 // ----------------------------------------------
1983 // Stores vertex parameters.
1984 // VERT array of (x,y,z) position of the vertex
1985 // NTBEAM beam track number origin of the vertex
1986 // =0 if none exists
1987 // NTTARG target track number origin of the vertex
1988 // UBUF user array of NUBUF floating point numbers
1990 // NVTX new vertex number (=0 in case of error).
1991 // Prints vertex parameters.
1992 // IVTX for vertex IVTX.
1993 // (For all vertices if IVTX=0)
1994 // Stores long life track parameters.
1995 // PLAB components of momentum
1996 // IPART type of particle (see GSPART)
1997 // NV vertex number origin of track
1998 // UBUF array of NUBUF floating point user parameters
2000 // NT track number (if=0 error).
2001 // Retrieves long life track parameters.
2002 // ITRA track number for which parameters are requested
2003 // VERT vector origin of the track
2004 // PVERT 4 momentum components at the track origin
2005 // IPART particle type (=0 if track ITRA does not exist)
2006 // NVERT vertex number origin of the track
2007 // UBUF user words stored in GSKINE.
2008 // Prints initial track parameters.
2009 // ITRA for track ITRA
2010 // (For all tracks if ITRA=0)
2014 gfkine(itra,vert,pvert,ipart,nvert,ubuf,nbuf);
2017 //_____________________________________________________________________________
2018 void TGeant3::Gfvert(Int_t nvtx, Float_t *v, Int_t &ntbeam, Int_t &nttarg,
2022 // Retrieves the parameter of a vertex bank
2023 // Vertex is generated from tracks NTBEAM NTTARG
2024 // NVTX is the new vertex number
2028 gfvert(nvtx,v,ntbeam,nttarg,tofg,ubuf,nbuf);
2031 //_____________________________________________________________________________
2032 Int_t TGeant3::Gskine(Float_t *plab, Int_t ipart, Int_t nv, Float_t *buf,
2036 // Store kinematics of track NT into data structure
2037 // Track is coming from vertex NV
2040 gskine(plab, ipart, nv, buf, nwbuf, nt);
2044 //_____________________________________________________________________________
2045 Int_t TGeant3::Gsvert(Float_t *v, Int_t ntbeam, Int_t nttarg, Float_t *ubuf,
2049 // Creates a new vertex bank
2050 // Vertex is generated from tracks NTBEAM NTTARG
2051 // NVTX is the new vertex number
2054 gsvert(v, ntbeam, nttarg, ubuf, nwbuf, nwtx);
2058 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2060 // Functions from GPHYS
2062 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2064 //_____________________________________________________________________________
2065 void TGeant3::Gphysi()
2068 // Initialise material constants for all the physics
2069 // mechanisms used by GEANT
2074 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2076 // Functions from GTRAK
2078 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2080 //_____________________________________________________________________________
2081 void TGeant3::Gdebug()
2084 // Debug the current step
2089 //_____________________________________________________________________________
2090 void TGeant3::Gekbin()
2093 // To find bin number in kinetic energy table
2094 // stored in ELOW(NEKBIN)
2099 //_____________________________________________________________________________
2100 void TGeant3::Gfinds()
2103 // Returns the set/volume parameters corresponding to
2104 // the current space point in /GCTRAK/
2105 // and fill common /GCSETS/
2107 // IHSET user set identifier
2108 // IHDET user detector identifier
2109 // ISET set number in JSET
2110 // IDET detector number in JS=LQ(JSET-ISET)
2111 // IDTYPE detector type (1,2)
2112 // NUMBV detector volume numbers (array of length NVNAME)
2113 // NVNAME number of volume levels
2118 //_____________________________________________________________________________
2119 void TGeant3::Gsking(Int_t igk)
2122 // Stores in stack JSTAK either the IGKth track of /GCKING/,
2123 // or the NGKINE tracks when IGK is 0.
2128 //_____________________________________________________________________________
2129 void TGeant3::Gskpho(Int_t igk)
2132 // Stores in stack JSTAK either the IGKth Cherenkov photon of
2133 // /GCKIN2/, or the NPHOT tracks when IGK is 0.
2138 //_____________________________________________________________________________
2139 void TGeant3::Gsstak(Int_t iflag)
2142 // Stores in auxiliary stack JSTAK the particle currently
2143 // described in common /GCKINE/.
2145 // On request, creates also an entry in structure JKINE :
2147 // 0 : No entry in JKINE structure required (user)
2148 // 1 : New entry in JVERTX / JKINE structures required (user)
2149 // <0 : New entry in JKINE structure at vertex -IFLAG (user)
2150 // 2 : Entry in JKINE structure exists already (from GTREVE)
2155 //_____________________________________________________________________________
2156 void TGeant3::Gsxyz()
2159 // Store space point VECT in banks JXYZ
2164 //_____________________________________________________________________________
2165 void TGeant3::Gtrack()
2168 // Controls tracking of current particle
2173 //_____________________________________________________________________________
2174 void TGeant3::Gtreve()
2177 // Controls tracking of all particles belonging to the current event
2182 //_____________________________________________________________________________
2183 void TGeant3::Gtreve_root()
2186 // Controls tracking of all particles belonging to the current event
2191 //_____________________________________________________________________________
2192 void TGeant3::Grndm(Float_t *rvec, const Int_t len) const
2195 // To generate a vector RVECV of LEN random numbers
2196 // Copy of the CERN Library routine RANECU
2200 //_____________________________________________________________________________
2201 void TGeant3::Grndmq(Int_t &is1, Int_t &is2, const Int_t iseq,
2202 const Text_t *chopt)
2205 // To set/retrieve the seed of the random number generator
2207 grndmq(is1,is2,iseq,PASSCHARD(chopt) PASSCHARL(chopt));
2210 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2212 // Functions from GDRAW
2214 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2216 //_____________________________________________________________________________
2217 void TGeant3::Gdxyz(Int_t it)
2220 // Draw the points stored with Gsxyz relative to track it
2225 //_____________________________________________________________________________
2226 void TGeant3::Gdcxyz()
2229 // Draw the position of the current track
2234 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2236 // Functions from GGEOM
2238 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2240 //_____________________________________________________________________________
2241 void TGeant3::Gdtom(Float_t *xd, Float_t *xm, Int_t iflag)
2244 // Computes coordinates XM (Master Reference System
2245 // knowing the coordinates XD (Detector Ref System)
2246 // The local reference system can be initialized by
2247 // - the tracking routines and GDTOM used in GUSTEP
2248 // - a call to GSCMED(NLEVEL,NAMES,NUMBER)
2249 // (inverse routine is GMTOD)
2251 // If IFLAG=1 convert coordinates
2252 // IFLAG=2 convert direction cosinus
2254 gdtom(xd, xm, iflag);
2257 //_____________________________________________________________________________
2258 void TGeant3::Glmoth(const char* iudet, Int_t iunum, Int_t &nlev, Int_t *lvols,
2262 // Loads the top part of the Volume tree in LVOLS (IVO's),
2263 // LINDX (IN indices) for a given volume defined through
2264 // its name IUDET and number IUNUM.
2266 // The routine stores only upto the last level where JVOLUM
2267 // data structure is developed. If there is no development
2268 // above the current level, it returns NLEV zero.
2270 glmoth(PASSCHARD(iudet), iunum, nlev, lvols, lindx, idum PASSCHARL(iudet));
2273 //_____________________________________________________________________________
2274 void TGeant3::Gmedia(Float_t *x, Int_t &numed)
2277 // Finds in which volume/medium the point X is, and updates the
2278 // common /GCVOLU/ and the structure JGPAR accordingly.
2280 // NUMED returns the tracking medium number, or 0 if point is
2281 // outside the experimental setup.
2286 //_____________________________________________________________________________
2287 void TGeant3::Gmtod(Float_t *xm, Float_t *xd, Int_t iflag)
2290 // Computes coordinates XD (in DRS)
2291 // from known coordinates XM in MRS
2292 // The local reference system can be initialized by
2293 // - the tracking routines and GMTOD used in GUSTEP
2294 // - a call to GMEDIA(XM,NUMED)
2295 // - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER)
2296 // (inverse routine is GDTOM)
2298 // If IFLAG=1 convert coordinates
2299 // IFLAG=2 convert direction cosinus
2301 gmtod(xm, xd, iflag);
2304 //_____________________________________________________________________________
2305 void TGeant3::Gsdvn(const char *name, const char *mother, Int_t ndiv,
2309 // Create a new volume by dividing an existing one
2312 // MOTHER Mother volume name
2313 // NDIV Number of divisions
2316 // X,Y,Z of CAXIS will be translated to 1,2,3 for IAXIS.
2317 // It divides a previously defined volume.
2322 Vname(mother,vmother);
2323 gsdvn(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis PASSCHARL(vname)
2324 PASSCHARL(vmother));
2327 //_____________________________________________________________________________
2328 void TGeant3::Gsdvn2(const char *name, const char *mother, Int_t ndiv,
2329 Int_t iaxis, Float_t c0i, Int_t numed)
2332 // Create a new volume by dividing an existing one
2334 // Divides mother into ndiv divisions called name
2335 // along axis iaxis starting at coordinate value c0.
2336 // the new volume created will be medium number numed.
2341 Vname(mother,vmother);
2342 gsdvn2(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis, c0i, numed
2343 PASSCHARL(vname) PASSCHARL(vmother));
2346 //_____________________________________________________________________________
2347 void TGeant3::Gsdvs(const char *name, const char *mother, Float_t step,
2348 Int_t iaxis, Int_t numed)
2351 // Create a new volume by dividing an existing one
2356 Vname(mother,vmother);
2357 gsdvs(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed
2358 PASSCHARL(vname) PASSCHARL(vmother));
2361 //_____________________________________________________________________________
2362 void TGeant3::Gsdvs2(const char *name, const char *mother, Float_t step,
2363 Int_t iaxis, Float_t c0, Int_t numed)
2366 // Create a new volume by dividing an existing one
2371 Vname(mother,vmother);
2372 gsdvs2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0, numed
2373 PASSCHARL(vname) PASSCHARL(vmother));
2376 //_____________________________________________________________________________
2377 void TGeant3::Gsdvt(const char *name, const char *mother, Float_t step,
2378 Int_t iaxis, Int_t numed, Int_t ndvmx)
2381 // Create a new volume by dividing an existing one
2383 // Divides MOTHER into divisions called NAME along
2384 // axis IAXIS in steps of STEP. If not exactly divisible
2385 // will make as many as possible and will centre them
2386 // with respect to the mother. Divisions will have medium
2387 // number NUMED. If NUMED is 0, NUMED of MOTHER is taken.
2388 // NDVMX is the expected maximum number of divisions
2389 // (If 0, no protection tests are performed)
2394 Vname(mother,vmother);
2395 gsdvt(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed, ndvmx
2396 PASSCHARL(vname) PASSCHARL(vmother));
2399 //_____________________________________________________________________________
2400 void TGeant3::Gsdvt2(const char *name, const char *mother, Float_t step,
2401 Int_t iaxis, Float_t c0, Int_t numed, Int_t ndvmx)
2404 // Create a new volume by dividing an existing one
2406 // Divides MOTHER into divisions called NAME along
2407 // axis IAXIS starting at coordinate value C0 with step
2409 // The new volume created will have medium number NUMED.
2410 // If NUMED is 0, NUMED of mother is taken.
2411 // NDVMX is the expected maximum number of divisions
2412 // (If 0, no protection tests are performed)
2417 Vname(mother,vmother);
2418 gsdvt2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0,
2419 numed, ndvmx PASSCHARL(vname) PASSCHARL(vmother));
2422 //_____________________________________________________________________________
2423 void TGeant3::Gsord(const char *name, Int_t iax)
2426 // Flags volume CHNAME whose contents will have to be ordered
2427 // along axis IAX, by setting the search flag to -IAX
2431 // IAX = 4 Rxy (static ordering only -> GTMEDI)
2432 // IAX = 14 Rxy (also dynamic ordering -> GTNEXT)
2433 // IAX = 5 Rxyz (static ordering only -> GTMEDI)
2434 // IAX = 15 Rxyz (also dynamic ordering -> GTNEXT)
2435 // IAX = 6 PHI (PHI=0 => X axis)
2436 // IAX = 7 THETA (THETA=0 => Z axis)
2440 gsord(PASSCHARD(vname), iax PASSCHARL(vname));
2443 //_____________________________________________________________________________
2444 void TGeant3::Gspos(const char *name, Int_t nr, const char *mother, Float_t x,
2445 Float_t y, Float_t z, Int_t irot, const char *konly)
2448 // Position a volume into an existing one
2451 // NUMBER Copy number of the volume
2452 // MOTHER Mother volume name
2453 // X X coord. of the volume in mother ref. sys.
2454 // Y Y coord. of the volume in mother ref. sys.
2455 // Z Z coord. of the volume in mother ref. sys.
2456 // IROT Rotation matrix number w.r.t. mother ref. sys.
2457 // ONLY ONLY/MANY flag
2459 // It positions a previously defined volume in the mother.
2464 Vname(mother,vmother);
2465 gspos(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2466 PASSCHARD(konly) PASSCHARL(vname) PASSCHARL(vmother)
2470 //_____________________________________________________________________________
2471 void TGeant3::Gsposp(const char *name, Int_t nr, const char *mother,
2472 Float_t x, Float_t y, Float_t z, Int_t irot,
2473 const char *konly, Float_t *upar, Int_t np )
2476 // Place a copy of generic volume NAME with user number
2477 // NR inside MOTHER, with its parameters UPAR(1..NP)
2482 Vname(mother,vmother);
2483 gsposp(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2484 PASSCHARD(konly), upar, np PASSCHARL(vname) PASSCHARL(vmother)
2488 //_____________________________________________________________________________
2489 void TGeant3::Gsrotm(Int_t nmat, Float_t theta1, Float_t phi1, Float_t theta2,
2490 Float_t phi2, Float_t theta3, Float_t phi3)
2493 // nmat Rotation matrix number
2494 // THETA1 Polar angle for axis I
2495 // PHI1 Azimuthal angle for axis I
2496 // THETA2 Polar angle for axis II
2497 // PHI2 Azimuthal angle for axis II
2498 // THETA3 Polar angle for axis III
2499 // PHI3 Azimuthal angle for axis III
2501 // It defines the rotation matrix number IROT.
2503 gsrotm(nmat, theta1, phi1, theta2, phi2, theta3, phi3);
2506 //_____________________________________________________________________________
2507 void TGeant3::Gprotm(Int_t nmat)
2510 // To print rotation matrices structure JROTM
2511 // nmat Rotation matrix number
2516 //_____________________________________________________________________________
2517 Int_t TGeant3::Gsvolu(const char *name, const char *shape, Int_t nmed,
2518 Float_t *upar, Int_t npar)
2522 // SHAPE Volume type
2523 // NUMED Tracking medium number
2524 // NPAR Number of shape parameters
2525 // UPAR Vector containing shape parameters
2527 // It creates a new volume in the JVOLUM data structure.
2533 Vname(shape,vshape);
2534 gsvolu(PASSCHARD(vname), PASSCHARD(vshape), nmed, upar, npar, ivolu
2535 PASSCHARL(vname) PASSCHARL(vshape));
2539 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2541 // T H E D R A W I N G P A C K A G E
2542 // ======================================
2543 // Drawing functions. These functions allow the visualization in several ways
2544 // of the volumes defined in the geometrical data structure. It is possible
2545 // to draw the logical tree of volumes belonging to the detector (DTREE),
2546 // to show their geometrical specification (DSPEC,DFSPC), to draw them
2547 // and their cut views (DRAW, DCUT). Moreover, it is possible to execute
2548 // these commands when the hidden line removal option is activated; in
2549 // this case, the volumes can be also either translated in the space
2550 // (SHIFT), or clipped by boolean operation (CVOL). In addition, it is
2551 // possible to fill the surfaces of the volumes
2552 // with solid colours when the shading option (SHAD) is activated.
2553 // Several tools (ZOOM, LENS) have been developed to zoom detailed parts
2554 // of the detectors or to scan physical events as well.
2555 // Finally, the command MOVE will allow the rotation, translation and zooming
2556 // on real time parts of the detectors or tracks and hits of a simulated event.
2557 // Ray-tracing commands. In case the command (DOPT RAYT ON) is executed,
2558 // the drawing is performed by the Geant ray-tracing;
2559 // automatically, the color is assigned according to the tracking medium of each
2560 // volume and the volumes with a density lower/equal than the air are considered
2561 // transparent; if the option (USER) is set (ON) (again via the command (DOPT)),
2562 // the user can set color and visibility for the desired volumes via the command
2563 // (SATT), as usual, relatively to the attributes (COLO) and (SEEN).
2564 // The resolution can be set via the command (SATT * FILL VALUE), where (VALUE)
2565 // is the ratio between the number of pixels drawn and 20 (user coordinates).
2566 // Parallel view and perspective view are possible (DOPT PROJ PARA/PERS); in the
2567 // first case, we assume that the first mother volume of the tree is a box with
2568 // dimensions 10000 X 10000 X 10000 cm and the view point (infinetely far) is
2569 // 5000 cm far from the origin along the Z axis of the user coordinates; in the
2570 // second case, the distance between the observer and the origin of the world
2571 // reference system is set in cm by the command (PERSP NAME VALUE); grand-angle
2572 // or telescopic effects can be achieved changing the scale factors in the command
2573 // (DRAW). When the final picture does not occupy the full window,
2574 // mapping the space before tracing can speed up the drawing, but can also
2575 // produce less precise results; values from 1 to 4 are allowed in the command
2576 // (DOPT MAPP VALUE), the mapping being more precise for increasing (VALUE); for
2577 // (VALUE = 0) no mapping is performed (therefore max precision and lowest speed).
2578 // The command (VALCUT) allows the cutting of the detector by three planes
2579 // ortogonal to the x,y,z axis. The attribute (LSTY) can be set by the command
2580 // SATT for any desired volume and can assume values from 0 to 7; it determines
2581 // the different light processing to be performed for different materials:
2582 // 0 = dark-matt, 1 = bright-matt, 2 = plastic, 3 = ceramic, 4 = rough-metals,
2583 // 5 = shiny-metals, 6 = glass, 7 = mirror. The detector is assumed to be in the
2584 // dark, the ambient light luminosity is 0.2 for each basic hue (the saturation
2585 // is 0.9) and the observer is assumed to have a light source (therefore he will
2586 // produce parallel light in the case of parallel view and point-like-source
2587 // light in the case of perspective view).
2589 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2591 //_____________________________________________________________________________
2592 void TGeant3::Gsatt(const char *name, const char *att, Int_t val)
2596 // IOPT Name of the attribute to be set
2597 // IVAL Value to which the attribute is to be set
2599 // name= "*" stands for all the volumes.
2600 // iopt can be chosen among the following :
2602 // WORK 0=volume name is inactive for the tracking
2603 // 1=volume name is active for the tracking (default)
2605 // SEEN 0=volume name is invisible
2606 // 1=volume name is visible (default)
2607 // -1=volume invisible with all its descendants in the tree
2608 // -2=volume visible but not its descendants in the tree
2610 // LSTY line style 1,2,3,... (default=1)
2611 // LSTY=7 will produce a very precise approximation for
2612 // revolution bodies.
2614 // LWID line width -7,...,1,2,3,..7 (default=1)
2615 // LWID<0 will act as abs(LWID) was set for the volume
2616 // and for all the levels below it. When SHAD is 'ON', LWID
2617 // represent the linewidth of the scan lines filling the surfaces
2618 // (whereas the FILL value represent their number). Therefore
2619 // tuning this parameter will help to obtain the desired
2620 // quality/performance ratio.
2622 // COLO colour code -166,...,1,2,..166 (default=1)
2624 // n=2=red; n=17+m, m=0,25, increasing luminosity according to 'm';
2625 // n=3=green; n=67+m, m=0,25, increasing luminosity according to 'm';
2626 // n=4=blue; n=117+m, m=0,25, increasing luminosity according to 'm';
2627 // n=5=yellow; n=42+m, m=0,25, increasing luminosity according to 'm';
2628 // n=6=violet; n=142+m, m=0,25, increasing luminosity according to 'm';
2629 // n=7=lightblue; n=92+m, m=0,25, increasing luminosity according to 'm';
2630 // colour=n*10+m, m=1,2,...9, will produce the same colour
2631 // as 'n', but with increasing luminosity according to 'm';
2632 // COLO<0 will act as if abs(COLO) was set for the volume
2633 // and for all the levels below it.
2634 // When for a volume the attribute FILL is > 1 (and the
2635 // option SHAD is on), the ABS of its colour code must be < 8
2636 // because an automatic shading of its faces will be
2639 // FILL (1992) fill area -7,...,0,1,...7 (default=0)
2640 // when option SHAD is "on" the FILL attribute of any
2641 // volume can be set different from 0 (normal drawing);
2642 // if it is set to 1, the faces of such volume will be filled
2643 // with solid colours; if ABS(FILL) is > 1, then a light
2644 // source is placed along the observer line, and the faces of
2645 // such volumes will be painted by colours whose luminosity
2646 // will depend on the amount of light reflected;
2647 // if ABS(FILL) = 1, then it is possible to use all the 166
2648 // colours of the colour table, becouse the automatic shading
2649 // is not performed;
2650 // for increasing values of FILL the drawing will be performed
2651 // with higher and higher resolution improving the quality (the
2652 // number of scan lines used to fill the faces increases with FILL);
2653 // it is possible to set different values of FILL
2654 // for different volumes, in order to optimize at the same time
2655 // the performance and the quality of the picture;
2656 // FILL<0 will act as if abs(FILL) was set for the volume
2657 // and for all the levels below it.
2658 // This kind of drawing can be saved in 'picture files'
2659 // or in view banks.
2660 // 0=drawing without fill area
2661 // 1=faces filled with solid colours and resolution = 6
2662 // 2=lowest resolution (very fast)
2663 // 3=default resolution
2664 // 4=.................
2665 // 5=.................
2666 // 6=.................
2668 // Finally, if a coloured background is desired, the FILL
2669 // attribute for the first volume of the tree must be set
2670 // equal to -abs(colo), colo being >0 and <166.
2672 // SET set number associated to volume name
2673 // DET detector number associated to volume name
2674 // DTYP detector type (1,2)
2681 gsatt(PASSCHARD(vname), PASSCHARD(vatt), val PASSCHARL(vname)
2685 //_____________________________________________________________________________
2686 void TGeant3::Gfpara(const char *name, Int_t number, Int_t intext, Int_t& npar,
2687 Int_t& natt, Float_t* par, Float_t* att)
2690 // Find the parameters of a volume
2692 gfpara(PASSCHARD(name), number, intext, npar, natt, par, att
2696 //_____________________________________________________________________________
2697 void TGeant3::Gckpar(Int_t ish, Int_t npar, Float_t* par)
2700 // Check the parameters of a shape
2702 gckpar(ish,npar,par);
2705 //_____________________________________________________________________________
2706 void TGeant3::Gckmat(Int_t itmed, char* natmed)
2709 // Check the parameters of a tracking medium
2711 gckmat(itmed, PASSCHARD(natmed) PASSCHARL(natmed));
2714 //_____________________________________________________________________________
2715 void TGeant3::Gdelete(Int_t iview)
2718 // IVIEW View number
2720 // It deletes a view bank from memory.
2725 //_____________________________________________________________________________
2726 void TGeant3::Gdopen(Int_t iview)
2729 // IVIEW View number
2731 // When a drawing is very complex and requires a long time to be
2732 // executed, it can be useful to store it in a view bank: after a
2733 // call to DOPEN and the execution of the drawing (nothing will
2734 // appear on the screen), and after a necessary call to DCLOSE,
2735 // the contents of the bank can be displayed in a very fast way
2736 // through a call to DSHOW; therefore, the detector can be easily
2737 // zoomed many times in different ways. Please note that the pictures
2738 // with solid colours can now be stored in a view bank or in 'PICTURE FILES'
2745 //_____________________________________________________________________________
2746 void TGeant3::Gdclose()
2749 // It closes the currently open view bank; it must be called after the
2750 // end of the drawing to be stored.
2755 //_____________________________________________________________________________
2756 void TGeant3::Gdshow(Int_t iview)
2759 // IVIEW View number
2761 // It shows on the screen the contents of a view bank. It
2762 // can be called after a view bank has been closed.
2767 //_____________________________________________________________________________
2768 void TGeant3::Gdopt(const char *name,const char *value)
2772 // VALUE Option value
2774 // To set/modify the drawing options.
2777 // THRZ ON Draw tracks in R vs Z
2778 // OFF (D) Draw tracks in X,Y,Z
2781 // PROJ PARA (D) Parallel projection
2783 // TRAK LINE (D) Trajectory drawn with lines
2784 // POIN " " with markers
2785 // HIDE ON Hidden line removal using the CG package
2786 // OFF (D) No hidden line removal
2787 // SHAD ON Fill area and shading of surfaces.
2788 // OFF (D) Normal hidden line removal.
2789 // RAYT ON Ray-tracing on.
2790 // OFF (D) Ray-tracing off.
2791 // EDGE OFF Does not draw contours when shad is on.
2792 // ON (D) Normal shading.
2793 // MAPP 1,2,3,4 Mapping before ray-tracing.
2794 // 0 (D) No mapping.
2795 // USER ON User graphics options in the raytracing.
2796 // OFF (D) Automatic graphics options.
2802 Vname(value,vvalue);
2803 gdopt(PASSCHARD(vname), PASSCHARD(vvalue) PASSCHARL(vname)
2807 //_____________________________________________________________________________
2808 void TGeant3::Gdraw(const char *name,Float_t theta, Float_t phi, Float_t psi,
2809 Float_t u0,Float_t v0,Float_t ul,Float_t vl)
2814 // THETA Viewing angle theta (for 3D projection)
2815 // PHI Viewing angle phi (for 3D projection)
2816 // PSI Viewing angle psi (for 2D rotation)
2817 // U0 U-coord. (horizontal) of volume origin
2818 // V0 V-coord. (vertical) of volume origin
2819 // SU Scale factor for U-coord.
2820 // SV Scale factor for V-coord.
2822 // This function will draw the volumes,
2823 // selected with their graphical attributes, set by the Gsatt
2824 // facility. The drawing may be performed with hidden line removal
2825 // and with shading effects according to the value of the options HIDE
2826 // and SHAD; if the option SHAD is ON, the contour's edges can be
2827 // drawn or not. If the option HIDE is ON, the detector can be
2828 // exploded (BOMB), clipped with different shapes (CVOL), and some
2829 // of its parts can be shifted from their original
2830 // position (SHIFT). When HIDE is ON, if
2831 // the drawing requires more than the available memory, the program
2832 // will evaluate and display the number of missing words
2833 // (so that the user can increase the
2834 // size of its ZEBRA store). Finally, at the end of each drawing (with HIDE on),
2835 // the program will print messages about the memory used and
2836 // statistics on the volumes' visibility.
2837 // The following commands will produce the drawing of a green
2838 // volume, specified by NAME, without using the hidden line removal
2839 // technique, using the hidden line removal technique,
2840 // with different linewidth and colour (red), with
2841 // solid colour, with shading of surfaces, and without edges.
2842 // Finally, some examples are given for the ray-tracing. (A possible
2843 // string for the NAME of the volume can be found using the command DTREE).
2849 if (fGcvdma->raytra != 1) {
2850 gdraw(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
2852 gdrayt(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
2856 //_____________________________________________________________________________
2857 void TGeant3::Gdrawc(const char *name,Int_t axis, Float_t cut,Float_t u0,
2858 Float_t v0,Float_t ul,Float_t vl)
2863 // CUTVAL Cut plane distance from the origin along the axis
2865 // U0 U-coord. (horizontal) of volume origin
2866 // V0 V-coord. (vertical) of volume origin
2867 // SU Scale factor for U-coord.
2868 // SV Scale factor for V-coord.
2870 // The cut plane is normal to caxis (X,Y,Z), corresponding to iaxis (1,2,3),
2871 // and placed at the distance cutval from the origin.
2872 // The resulting picture is seen from the the same axis.
2873 // When HIDE Mode is ON, it is possible to get the same effect with
2874 // the CVOL/BOX function.
2880 gdrawc(PASSCHARD(vname), axis,cut,u0,v0,ul,vl PASSCHARL(vname));
2883 //_____________________________________________________________________________
2884 void TGeant3::Gdrawx(const char *name,Float_t cutthe, Float_t cutphi,
2885 Float_t cutval, Float_t theta, Float_t phi, Float_t u0,
2886 Float_t v0,Float_t ul,Float_t vl)
2890 // CUTTHE Theta angle of the line normal to cut plane
2891 // CUTPHI Phi angle of the line normal to cut plane
2892 // CUTVAL Cut plane distance from the origin along the axis
2894 // THETA Viewing angle theta (for 3D projection)
2895 // PHI Viewing angle phi (for 3D projection)
2896 // U0 U-coord. (horizontal) of volume origin
2897 // V0 V-coord. (vertical) of volume origin
2898 // SU Scale factor for U-coord.
2899 // SV Scale factor for V-coord.
2901 // The cut plane is normal to the line given by the cut angles
2902 // cutthe and cutphi and placed at the distance cutval from the origin.
2903 // The resulting picture is seen from the viewing angles theta,phi.
2909 gdrawx(PASSCHARD(vname), cutthe,cutphi,cutval,theta,phi,u0,v0,ul,vl
2913 //_____________________________________________________________________________
2914 void TGeant3::Gdhead(Int_t isel, const char *name, Float_t chrsiz)
2919 // ISEL Option flag D=111110
2921 // CHRSIZ Character size (cm) of title NAME D=0.6
2924 // 0 to have only the header lines
2925 // xxxxx1 to add the text name centered on top of header
2926 // xxxx1x to add global detector name (first volume) on left
2927 // xxx1xx to add date on right
2928 // xx1xxx to select thick characters for text on top of header
2929 // x1xxxx to add the text 'EVENT NR x' on top of header
2930 // 1xxxxx to add the text 'RUN NR x' on top of header
2931 // NOTE that ISEL=x1xxx1 or ISEL=1xxxx1 are illegal choices,
2932 // i.e. they generate overwritten text.
2934 gdhead(isel,PASSCHARD(name),chrsiz PASSCHARL(name));
2937 //_____________________________________________________________________________
2938 void TGeant3::Gdman(Float_t u, Float_t v, const char *type)
2941 // Draw a 2D-man at position (U0,V0)
2943 // U U-coord. (horizontal) of the centre of man' R
2944 // V V-coord. (vertical) of the centre of man' R
2945 // TYPE D='MAN' possible values: 'MAN,WM1,WM2,WM3'
2947 // CALL GDMAN(u,v),CALL GDWMN1(u,v),CALL GDWMN2(u,v),CALL GDWMN2(u,v)
2948 // It superimposes the picure of a man or of a woman, chosen among
2949 // three different ones, with the same scale factors as the detector
2950 // in the current drawing.
2953 if (opt.Contains("WM1")) {
2955 } else if (opt.Contains("WM3")) {
2957 } else if (opt.Contains("WM2")) {
2964 //_____________________________________________________________________________
2965 void TGeant3::Gdspec(const char *name)
2970 // Shows 3 views of the volume (two cut-views and a 3D view), together with
2971 // its geometrical specifications. The 3D drawing will
2972 // be performed according the current values of the options HIDE and
2973 // SHAD and according the current SetClipBox clipping parameters for that
2980 gdspec(PASSCHARD(vname) PASSCHARL(vname));
2983 //_____________________________________________________________________________
2984 void TGeant3::DrawOneSpec(const char *name)
2987 // Function called when one double-clicks on a volume name
2988 // in a TPavelabel drawn by Gdtree.
2990 THIGZ *higzSave = higz;
2991 higzSave->SetName("higzSave");
2992 THIGZ *higzSpec = (THIGZ*)gROOT->FindObject("higzSpec");
2993 //printf("DrawOneSpec, higz=%x, higzSpec=%x\n",higz,higzSpec);
2994 if (higzSpec) higz = higzSpec;
2995 else higzSpec = new THIGZ(defSize);
2996 higzSpec->SetName("higzSpec");
3001 gdspec(PASSCHARD(vname) PASSCHARL(vname));
3004 higzSave->SetName("higz");
3008 //_____________________________________________________________________________
3009 void TGeant3::Gdtree(const char *name,Int_t levmax, Int_t isel)
3013 // LEVMAX Depth level
3016 // This function draws the logical tree,
3017 // Each volume in the tree is represented by a TPaveTree object.
3018 // Double-clicking on a TPaveTree draws the specs of the corresponding volume.
3019 // Use TPaveTree pop-up menu to select:
3022 // - drawing tree of parent
3028 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
3032 //_____________________________________________________________________________
3033 void TGeant3::GdtreeParent(const char *name,Int_t levmax, Int_t isel)
3037 // LEVMAX Depth level
3040 // This function draws the logical tree of the parent of name.
3044 // Scan list of volumes in JVOLUM
3046 Int_t gname, i, jvo, in, nin, jin, num;
3047 strncpy((char *) &gname, name, 4);
3048 for(i=1; i<=fGcnum->nvolum; i++) {
3049 jvo = fZlq[fGclink->jvolum-i];
3050 nin = Int_t(fZq[jvo+3]);
3051 if (nin == -1) nin = 1;
3052 for (in=1;in<=nin;in++) {
3054 num = Int_t(fZq[jin+2]);
3055 if(gname == fZiq[fGclink->jvolum+num]) {
3056 strncpy(vname,(char*)&fZiq[fGclink->jvolum+i],4);
3058 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
3066 //_____________________________________________________________________________
3067 void TGeant3::SetABAN(Int_t par)
3070 // par = 1 particles will be stopped according to their residual
3071 // range if they are not in a sensitive material and are
3072 // far enough from the boundary
3073 // 0 particles are transported normally
3075 fGcphys->dphys1 = par;
3079 //_____________________________________________________________________________
3080 void TGeant3::SetANNI(Int_t par)
3083 // To control positron annihilation.
3084 // par =0 no annihilation
3085 // =1 annihilation. Decays processed.
3086 // =2 annihilation. No decay products stored.
3088 fGcphys->ianni = par;
3092 //_____________________________________________________________________________
3093 void TGeant3::SetAUTO(Int_t par)
3096 // To control automatic calculation of tracking medium parameters:
3097 // par =0 no automatic calculation;
3098 // =1 automati calculation.
3100 fGctrak->igauto = par;
3104 //_____________________________________________________________________________
3105 void TGeant3::SetBOMB(Float_t boom)
3108 // BOOM : Exploding factor for volumes position
3110 // To 'explode' the detector. If BOOM is positive (values smaller
3111 // than 1. are suggested, but any value is possible)
3112 // all the volumes are shifted by a distance
3113 // proportional to BOOM along the direction between their centre
3114 // and the origin of the MARS; the volumes which are symmetric
3115 // with respect to this origin are simply not shown.
3116 // BOOM equal to 0 resets the normal mode.
3117 // A negative (greater than -1.) value of
3118 // BOOM will cause an 'implosion'; for even lower values of BOOM
3119 // the volumes' positions will be reflected respect to the origin.
3120 // This command can be useful to improve the 3D effect for very
3121 // complex detectors. The following commands will make explode the
3128 //_____________________________________________________________________________
3129 void TGeant3::SetBREM(Int_t par)
3132 // To control bremstrahlung.
3133 // par =0 no bremstrahlung
3134 // =1 bremstrahlung. Photon processed.
3135 // =2 bremstrahlung. No photon stored.
3137 fGcphys->ibrem = par;
3141 //_____________________________________________________________________________
3142 void TGeant3::SetCKOV(Int_t par)
3145 // To control Cerenkov production
3146 // par =0 no Cerenkov;
3148 // =2 Cerenkov with primary stopped at each step.
3150 fGctlit->itckov = par;
3154 //_____________________________________________________________________________
3155 void TGeant3::SetClipBox(const char *name,Float_t xmin,Float_t xmax,
3156 Float_t ymin,Float_t ymax,Float_t zmin,Float_t zmax)
3159 // The hidden line removal technique is necessary to visualize properly
3160 // very complex detectors. At the same time, it can be useful to visualize
3161 // the inner elements of a detector in detail. This function allows
3162 // subtractions (via boolean operation) of BOX shape from any part of
3163 // the detector, therefore showing its inner contents.
3164 // If "*" is given as the name of the
3165 // volume to be clipped, all volumes are clipped by the given box.
3166 // A volume can be clipped at most twice.
3167 // if a volume is explicitely clipped twice,
3168 // the "*" will not act on it anymore. Giving "." as the name
3169 // of the volume to be clipped will reset the clipping.
3171 // NAME Name of volume to be clipped
3173 // XMIN Lower limit of the Shape X coordinate
3174 // XMAX Upper limit of the Shape X coordinate
3175 // YMIN Lower limit of the Shape Y coordinate
3176 // YMAX Upper limit of the Shape Y coordinate
3177 // ZMIN Lower limit of the Shape Z coordinate
3178 // ZMAX Upper limit of the Shape Z coordinate
3180 // This function performs a boolean subtraction between the volume
3181 // NAME and a box placed in the MARS according the values of the given
3187 setclip(PASSCHARD(vname),xmin,xmax,ymin,ymax,zmin,zmax PASSCHARL(vname));
3190 //_____________________________________________________________________________
3191 void TGeant3::SetCOMP(Int_t par)
3194 // To control Compton scattering
3195 // par =0 no Compton
3196 // =1 Compton. Electron processed.
3197 // =2 Compton. No electron stored.
3200 fGcphys->icomp = par;
3203 //_____________________________________________________________________________
3204 void TGeant3::SetCUTS(Float_t cutgam,Float_t cutele,Float_t cutneu,
3205 Float_t cuthad,Float_t cutmuo ,Float_t bcute ,
3206 Float_t bcutm ,Float_t dcute ,Float_t dcutm ,
3207 Float_t ppcutm, Float_t tofmax)
3210 // CUTGAM Cut for gammas D=0.001
3211 // CUTELE Cut for electrons D=0.001
3212 // CUTHAD Cut for charged hadrons D=0.01
3213 // CUTNEU Cut for neutral hadrons D=0.01
3214 // CUTMUO Cut for muons D=0.01
3215 // BCUTE Cut for electron brems. D=-1.
3216 // BCUTM Cut for muon brems. D=-1.
3217 // DCUTE Cut for electron delta-rays D=-1.
3218 // DCUTM Cut for muon delta-rays D=-1.
3219 // PPCUTM Cut for e+e- pairs by muons D=0.01
3220 // TOFMAX Time of flight cut D=1.E+10
3222 // If the default values (-1.) for BCUTE ,BCUTM ,DCUTE ,DCUTM
3223 // are not modified, they will be set to CUTGAM,CUTGAM,CUTELE,CUTELE
3225 // If one of the parameters from CUTGAM to PPCUTM included
3226 // is modified, cross-sections and energy loss tables must be
3227 // recomputed via the function Gphysi.
3229 fGccuts->cutgam = cutgam;
3230 fGccuts->cutele = cutele;
3231 fGccuts->cutneu = cutneu;
3232 fGccuts->cuthad = cuthad;
3233 fGccuts->cutmuo = cutmuo;
3234 fGccuts->bcute = bcute;
3235 fGccuts->bcutm = bcutm;
3236 fGccuts->dcute = dcute;
3237 fGccuts->dcutm = dcutm;
3238 fGccuts->ppcutm = ppcutm;
3239 fGccuts->tofmax = tofmax;
3242 //_____________________________________________________________________________
3243 void TGeant3::SetDCAY(Int_t par)
3246 // To control Decay mechanism.
3247 // par =0 no decays.
3248 // =1 Decays. secondaries processed.
3249 // =2 Decays. No secondaries stored.
3251 fGcphys->idcay = par;
3255 //_____________________________________________________________________________
3256 void TGeant3::SetDEBU(Int_t emin, Int_t emax, Int_t emod)
3259 // Set the debug flag and frequency
3260 // Selected debug output will be printed from
3261 // event emin to even emax each emod event
3263 fGcflag->idemin = emin;
3264 fGcflag->idemax = emax;
3265 fGcflag->itest = emod;
3269 //_____________________________________________________________________________
3270 void TGeant3::SetDRAY(Int_t par)
3273 // To control delta rays mechanism.
3274 // par =0 no delta rays.
3275 // =1 Delta rays. secondaries processed.
3276 // =2 Delta rays. No secondaries stored.
3278 fGcphys->idray = par;
3281 //_____________________________________________________________________________
3282 void TGeant3::SetERAN(Float_t ekmin, Float_t ekmax, Int_t nekbin)
3285 // To control cross section tabulations
3286 // ekmin = minimum kinetic energy in GeV
3287 // ekmax = maximum kinetic energy in GeV
3288 // nekbin = number of logatithmic bins (<200)
3290 fGcmulo->ekmin = ekmin;
3291 fGcmulo->ekmax = ekmax;
3292 fGcmulo->nekbin = nekbin;
3295 //_____________________________________________________________________________
3296 void TGeant3::SetHADR(Int_t par)
3299 // To control hadronic interactions.
3300 // par =0 no hadronic interactions.
3301 // =1 Hadronic interactions. secondaries processed.
3302 // =2 Hadronic interactions. No secondaries stored.
3304 fGcphys->ihadr = par;
3307 //_____________________________________________________________________________
3308 void TGeant3::SetKINE(Int_t kine, Float_t xk1, Float_t xk2, Float_t xk3,
3309 Float_t xk4, Float_t xk5, Float_t xk6, Float_t xk7,
3310 Float_t xk8, Float_t xk9, Float_t xk10)
3313 // Set the variables in /GCFLAG/ IKINE, PKINE(10)
3314 // Their meaning is user defined
3316 fGckine->ikine = kine;
3317 fGckine->pkine[0] = xk1;
3318 fGckine->pkine[1] = xk2;
3319 fGckine->pkine[2] = xk3;
3320 fGckine->pkine[3] = xk4;
3321 fGckine->pkine[4] = xk5;
3322 fGckine->pkine[5] = xk6;
3323 fGckine->pkine[6] = xk7;
3324 fGckine->pkine[7] = xk8;
3325 fGckine->pkine[8] = xk9;
3326 fGckine->pkine[9] = xk10;
3329 //_____________________________________________________________________________
3330 void TGeant3::SetLOSS(Int_t par)
3333 // To control energy loss.
3334 // par =0 no energy loss;
3335 // =1 restricted energy loss fluctuations;
3336 // =2 complete energy loss fluctuations;
3338 // =4 no energy loss fluctuations.
3339 // If the value ILOSS is changed, then cross-sections and energy loss
3340 // tables must be recomputed via the command 'PHYSI'.
3342 fGcphys->iloss = par;
3346 //_____________________________________________________________________________
3347 void TGeant3::SetMULS(Int_t par)
3350 // To control multiple scattering.
3351 // par =0 no multiple scattering.
3352 // =1 Moliere or Coulomb scattering.
3353 // =2 Moliere or Coulomb scattering.
3354 // =3 Gaussian scattering.
3356 fGcphys->imuls = par;
3360 //_____________________________________________________________________________
3361 void TGeant3::SetMUNU(Int_t par)
3364 // To control muon nuclear interactions.
3365 // par =0 no muon-nuclear interactions.
3366 // =1 Nuclear interactions. Secondaries processed.
3367 // =2 Nuclear interactions. Secondaries not processed.
3369 fGcphys->imunu = par;
3372 //_____________________________________________________________________________
3373 void TGeant3::SetOPTI(Int_t par)
3376 // This flag controls the tracking optimisation performed via the
3378 // 1 no optimisation at all; GSORD calls disabled;
3379 // 0 no optimisation; only user calls to GSORD kept;
3380 // 1 all non-GSORDered volumes are ordered along the best axis;
3381 // 2 all volumes are ordered along the best axis.
3383 fGcopti->ioptim = par;
3386 //_____________________________________________________________________________
3387 void TGeant3::SetPAIR(Int_t par)
3390 // To control pair production mechanism.
3391 // par =0 no pair production.
3392 // =1 Pair production. secondaries processed.
3393 // =2 Pair production. No secondaries stored.
3395 fGcphys->ipair = par;
3399 //_____________________________________________________________________________
3400 void TGeant3::SetPFIS(Int_t par)
3403 // To control photo fission mechanism.
3404 // par =0 no photo fission.
3405 // =1 Photo fission. secondaries processed.
3406 // =2 Photo fission. No secondaries stored.
3408 fGcphys->ipfis = par;
3411 //_____________________________________________________________________________
3412 void TGeant3::SetPHOT(Int_t par)
3415 // To control Photo effect.
3416 // par =0 no photo electric effect.
3417 // =1 Photo effect. Electron processed.
3418 // =2 Photo effect. No electron stored.
3420 fGcphys->iphot = par;
3423 //_____________________________________________________________________________
3424 void TGeant3::SetRAYL(Int_t par)
3427 // To control Rayleigh scattering.
3428 // par =0 no Rayleigh scattering.
3431 fGcphys->irayl = par;
3434 //_____________________________________________________________________________
3435 void TGeant3::SetSWIT(Int_t sw, Int_t val)
3439 // val New switch value
3441 // Change one element of array ISWIT(10) in /GCFLAG/
3443 if (sw <= 0 || sw > 10) return;
3444 fGcflag->iswit[sw-1] = val;
3448 //_____________________________________________________________________________
3449 void TGeant3::SetTRIG(Int_t nevents)
3452 // Set number of events to be run
3454 fGcflag->nevent = nevents;
3457 //_____________________________________________________________________________
3458 void TGeant3::SetUserDecay(Int_t pdg)
3461 // Force the decays of particles to be done with Pythia
3462 // and not with the Geant routines.
3463 // just kill pointers doing mzdrop
3465 Int_t ipart = IdFromPDG(pdg);
3467 printf("Particle %d not in geant\n",pdg);
3470 Int_t jpart=fGclink->jpart;
3471 Int_t jpa=fZlq[jpart-ipart];
3474 Int_t jpa1=fZlq[jpa-1];
3476 mzdrop(fGcbank->ixcons,jpa1,PASSCHARD(" ") PASSCHARL(" "));
3477 Int_t jpa2=fZlq[jpa-2];
3479 mzdrop(fGcbank->ixcons,jpa2,PASSCHARD(" ") PASSCHARL(" "));
3483 //______________________________________________________________________________
3484 void TGeant3::Vname(const char *name, char *vname)
3487 // convert name to upper case. Make vname at least 4 chars
3489 Int_t l = strlen(name);
3492 for (i=0;i<l;i++) vname[i] = toupper(name[i]);
3493 for (i=l;i<4;i++) vname[i] = ' ';
3497 //______________________________________________________________________________
3498 void TGeant3::Ertrgo()
3503 //______________________________________________________________________________
3504 void TGeant3::Ertrak(const Float_t *const x1, const Float_t *const p1,
3505 const Float_t *x2, const Float_t *p2,
3506 Int_t ipa, Option_t *chopt)
3508 ertrak(x1,p1,x2,p2,ipa,PASSCHARD(chopt) PASSCHARL(chopt));
3511 //_____________________________________________________________________________
3512 void TGeant3::WriteEuclid(const char* filnam, const char* topvol,
3513 Int_t number, Int_t nlevel)
3517 // ******************************************************************
3519 // * Write out the geometry of the detector in EUCLID file format *
3521 // * filnam : will be with the extension .euc *
3522 // * topvol : volume name of the starting node *
3523 // * number : copy number of topvol (relevant for gsposp) *
3524 // * nlevel : number of levels in the tree structure *
3525 // * to be written out, starting from topvol *
3527 // * Author : M. Maire *
3529 // ******************************************************************
3531 // File filnam.tme is written out with the definitions of tracking
3532 // medias and materials.
3533 // As to restore original numbers for materials and medias, program
3534 // searches in the file euc_medi.dat and comparing main parameters of
3535 // the mat. defined inside geant and the one in file recognizes them
3536 // and is able to take number from file. If for any material or medium,
3537 // this procedure fails, ordering starts from 1.
3538 // Arrays IOTMED and IOMATE are used for this procedure
3540 const char shape[][5]={"BOX ","TRD1","TRD2","TRAP","TUBE","TUBS","CONE",
3541 "CONS","SPHE","PARA","PGON","PCON","ELTU","HYPE",
3543 Int_t i, end, itm, irm, jrm, k, nmed;
3547 char *filext, *filetme;
3548 char natmed[21], namate[21];
3549 char natmedc[21], namatec[21];
3550 char key[5], name[5], mother[5], konly[5];
3552 Int_t iadvol, iadtmd, iadrot, nwtot, iret;
3553 Int_t mlevel, numbr, natt, numed, nin, ndata;
3554 Int_t iname, ivo, ish, jvo, nvstak, ivstak;
3555 Int_t jdiv, ivin, in, jin, jvin, irot;
3556 Int_t jtm, imat, jma, flag=0, imatc;
3557 Float_t az, dens, radl, absl, a, step, x, y, z;
3558 Int_t npar, ndvmx, left;
3559 Float_t zc, densc, radlc, abslc, c0, tmaxfd;
3561 Int_t iomate[100], iotmed[100];
3562 Float_t par[50], att[20], ubuf[50];
3565 Int_t level, ndiv, iaxe;
3566 Int_t itmedc, nmatc, isvolc, ifieldc, nwbufc, isvol, nmat, ifield, nwbuf;
3567 Float_t fieldmc, tmaxfdc, stemaxc, deemaxc, epsilc, stminc, fieldm;
3568 Float_t tmaxf, stemax, deemax, epsil, stmin;
3569 const char *f10000="!\n%s\n!\n";
3570 //Open the input file
3572 for(i=0;i<end;i++) if(filnam[i]=='.') {
3576 filext=new char[end+5];
3577 filetme=new char[end+5];
3578 strncpy(filext,filnam,end);
3579 strncpy(filetme,filnam,end);
3581 // *** The output filnam name will be with extension '.euc'
3582 strcpy(&filext[end],".euc");
3583 strcpy(&filetme[end],".tme");
3584 lun=fopen(filext,"w");
3586 // *** Initialisation of the working space
3587 iadvol=fGcnum->nvolum;
3588 iadtmd=iadvol+fGcnum->nvolum;
3589 iadrot=iadtmd+fGcnum->ntmed;
3590 if(fGclink->jrotm) {
3591 fGcnum->nrotm=fZiq[fGclink->jrotm-2];
3595 nwtot=iadrot+fGcnum->nrotm;
3596 qws = new float[nwtot+1];
3597 for (i=0;i<nwtot+1;i++) qws[i]=0;
3600 if(nlevel==0) mlevel=20;
3602 // *** find the top volume and put it in the stak
3603 numbr = number>0 ? number : 1;
3604 Gfpara(topvol,numbr,1,npar,natt,par,att);
3606 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3611 // *** authorized shape ?
3612 strncpy((char *)&iname, topvol, 4);
3614 for(i=1; i<=fGcnum->nvolum; i++) if(fZiq[fGclink->jvolum+i]==iname) {
3618 jvo = fZlq[fGclink->jvolum-ivo];
3619 ish = Int_t (fZq[jvo+2]);
3621 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3628 iws[iadvol+ivo] = level;
3631 //*** flag all volumes and fill the stak
3635 // pick the next volume in stak
3637 ivo = TMath::Abs(iws[ivstak]);
3638 jvo = fZlq[fGclink->jvolum - ivo];
3640 // flag the tracking medium
3641 numed = Int_t (fZq[jvo + 4]);
3642 iws[iadtmd + numed] = 1;
3644 // get the daughters ...
3645 level = iws[iadvol+ivo];
3646 if (level < mlevel) {
3648 nin = Int_t (fZq[jvo + 3]);
3650 // from division ...
3652 jdiv = fZlq[jvo - 1];
3653 ivin = Int_t (fZq[jdiv + 2]);
3655 iws[nvstak] = -ivin;
3656 iws[iadvol+ivin] = level;
3658 // from position ...
3659 } else if (nin > 0) {
3660 for(in=1; in<=nin; in++) {
3661 jin = fZlq[jvo - in];
3662 ivin = Int_t (fZq[jin + 2 ]);
3663 jvin = fZlq[fGclink->jvolum - ivin];
3664 ish = Int_t (fZq[jvin + 2]);
3665 // authorized shape ?
3667 // not yet flagged ?
3668 if (iws[iadvol+ivin]==0) {
3671 iws[iadvol+ivin] = level;
3673 // flag the rotation matrix
3674 irot = Int_t ( fZq[jin + 4 ]);
3675 if (irot > 0) iws[iadrot+irot] = 1;
3681 // next volume in stak ?
3682 if (ivstak < nvstak) goto L10;
3684 // *** restore original material and media numbers
3685 // file euc_medi.dat is needed to compare materials and medias
3687 FILE* luncor=fopen("euc_medi.dat","r");
3690 for(itm=1; itm<=fGcnum->ntmed; itm++) {
3691 if (iws[iadtmd+itm] > 0) {
3692 jtm = fZlq[fGclink->jtmed-itm];
3693 strncpy(natmed,(char *)&fZiq[jtm+1],20);
3694 imat = Int_t (fZq[jtm+6]);
3695 jma = fZlq[fGclink->jmate-imat];
3697 printf(" *** GWEUCL *** material not defined for tracking medium %5i %s\n",itm,natmed);
3700 strncpy(namate,(char *)&fZiq[jma+1],20);
3703 //** find the material original number
3706 iret=fscanf(luncor,"%4s,%130s",key,card);
3707 if(iret<=0) goto L26;
3709 if(!strcmp(key,"MATE")) {
3710 sscanf(card,"%d %s %f %f %f %f %f %d",&imatc,namatec,&az,&zc,&densc,&radlc,&abslc,&nparc);
3711 Gfmate(imat,namate,a,z,dens,radl,absl,par,npar);
3712 if(!strcmp(namatec,namate)) {
3713 if(az==a && zc==z && densc==dens && radlc==radl
3714 && abslc==absl && nparc==nparc) {
3717 printf("*** GWEUCL *** material : %3d '%s' restored as %3d\n",imat,namate,imatc);
3719 printf("*** GWEUCL *** different definitions for material: %s\n",namate);
3723 if(strcmp(key,"END") && !flag) goto L23;
3725 printf("*** GWEUCL *** cannot restore original number for material: %s\n",namate);
3729 //*** restore original tracking medium number
3732 iret=fscanf(luncor,"%4s,%130s",key,card);
3733 if(iret<=0) goto L26;
3735 if (!strcmp(key,"TMED")) {
3736 sscanf(card,"%d %s %d %d %d %f %f %f %f %f %f %d\n",
3737 &itmedc,natmedc,&nmatc,&isvolc,&ifieldc,&fieldmc,
3738 &tmaxfdc,&stemaxc,&deemaxc,&epsilc,&stminc,&nwbufc);
3739 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxf,stemax,deemax,
3740 epsil,stmin,ubuf,&nwbuf);
3741 if(!strcmp(natmedc,natmed)) {
3742 if (iomate[nmat]==nmatc && nwbuf==nwbufc) {
3745 printf("*** GWEUCL *** medium : %3d '%20s' restored as %3d\n",
3748 printf("*** GWEUCL *** different definitions for tracking medium: %s\n",natmed);
3752 if(strcmp(key,"END") && !flag) goto L24;
3754 printf("cannot restore original number for medium : %s\n",natmed);
3762 L26: printf("*** GWEUCL *** cannot read the data file\n");
3764 L29: if(luncor) fclose (luncor);
3767 // *** write down the tracking medium definition
3769 strcpy(card,"! Tracking medium");
3770 fprintf(lun,f10000,card);
3772 for(itm=1;itm<=fGcnum->ntmed;itm++) {
3773 if (iws[iadtmd+itm]>0) {
3774 jtm = fZlq[fGclink->jtmed-itm];
3775 strncpy(natmed,(char *)&fZiq[jtm+1],20);
3777 imat = Int_t (fZq[jtm+6]);
3778 jma = fZlq[fGclink->jmate-imat];
3779 //* order media from one, if comparing with database failed
3781 iotmed[itm]=++imxtmed;
3782 iomate[imat]=++imxmate;
3787 printf(" *** GWEUCL *** material not defined for tracking medium %5d %s\n",
3790 strncpy(namate,(char *)&fZiq[jma+1],20);
3793 fprintf(lun,"TMED %3d '%20s' %3d '%20s'\n",iotmed[itm],natmed,iomate[imat],namate);
3797 //* *** write down the rotation matrix
3799 strcpy(card,"! Reperes");
3800 fprintf(lun,f10000,card);
3802 for(irm=1;irm<=fGcnum->nrotm;irm++) {
3803 if (iws[iadrot+irm]>0) {
3804 jrm = fZlq[fGclink->jrotm-irm];
3805 fprintf(lun,"ROTM %3d",irm);
3806 for(k=11;k<=16;k++) fprintf(lun," %8.3f",fZq[jrm+k]);
3811 //* *** write down the volume definition
3813 strcpy(card,"! Volumes");
3814 fprintf(lun,f10000,card);
3816 for(ivstak=1;ivstak<=nvstak;ivstak++) {
3819 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivo],4);
3821 jvo = fZlq[fGclink->jvolum-ivo];
3822 ish = Int_t (fZq[jvo+2]);
3823 nmed = Int_t (fZq[jvo+4]);
3824 npar = Int_t (fZq[jvo+5]);
3826 if (ivstak>1) for(i=0;i<npar;i++) par[i]=fZq[jvo+7+i];
3827 Gckpar (ish,npar,par);
3828 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,shape[ish-1],iotmed[nmed],npar);
3829 for(i=0;i<(npar-1)/6+1;i++) {
3832 for(k=0;k<(left<6?left:6);k++) fprintf(lun," %11.5f",par[i*6+k]);
3836 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,shape[ish-1],iotmed[nmed],npar);
3841 //* *** write down the division of volumes
3843 fprintf(lun,f10000,"! Divisions");
3844 for(ivstak=1;ivstak<=nvstak;ivstak++) {
3845 ivo = TMath::Abs(iws[ivstak]);
3846 jvo = fZlq[fGclink->jvolum-ivo];
3847 ish = Int_t (fZq[jvo+2]);
3848 nin = Int_t (fZq[jvo+3]);
3849 //* this volume is divided ...
3852 iaxe = Int_t ( fZq[jdiv+1]);
3853 ivin = Int_t ( fZq[jdiv+2]);
3854 ndiv = Int_t ( fZq[jdiv+3]);
3857 jvin = fZlq[fGclink->jvolum-ivin];
3858 nmed = Int_t ( fZq[jvin+4]);
3859 strncpy(mother,(char *)&fZiq[fGclink->jvolum+ivo ],4);
3861 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivin],4);
3863 if ((step<=0.)||(ish>=11)) {
3864 //* volume with negative parameter or gsposp or pgon ...
3865 fprintf(lun,"DIVN '%4s' '%4s' %3d %3d\n",name,mother,ndiv,iaxe);
3866 } else if ((ndiv<=0)||(ish==10)) {
3867 //* volume with negative parameter or gsposp or para ...
3868 ndvmx = TMath::Abs(ndiv);
3869 fprintf(lun,"DIVT '%4s' '%4s' %11.5f %3d %3d %3d\n",
3870 name,mother,step,iaxe,iotmed[nmed],ndvmx);
3872 //* normal volume : all kind of division are equivalent
3873 fprintf(lun,"DVT2 '%4s' '%4s' %11.5f %3d %11.5f %3d %3d\n",
3874 name,mother,step,iaxe,c0,iotmed[nmed],ndiv);
3879 //* *** write down the the positionnement of volumes
3881 fprintf(lun,f10000,"! Positionnements\n");
3883 for(ivstak = 1;ivstak<=nvstak;ivstak++) {
3884 ivo = TMath::Abs(iws[ivstak]);
3885 strncpy(mother,(char*)&fZiq[fGclink->jvolum+ivo ],4);
3887 jvo = fZlq[fGclink->jvolum-ivo];
3888 nin = Int_t( fZq[jvo+3]);
3889 //* this volume has daughters ...
3891 for (in=1;in<=nin;in++) {
3893 ivin = Int_t (fZq[jin +2]);
3894 numb = Int_t (fZq[jin +3]);
3895 irot = Int_t (fZq[jin +4]);
3899 strcpy(konly,"ONLY");
3900 if (fZq[jin+8]!=1.) strcpy(konly,"MANY");
3901 strncpy(name,(char*)&fZiq[fGclink->jvolum+ivin],4);
3903 jvin = fZlq[fGclink->jvolum-ivin];
3904 ish = Int_t (fZq[jvin+2]);
3905 //* gspos or gsposp ?
3906 ndata = fZiq[jin-1];
3908 fprintf(lun,"POSI '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s'\n",
3909 name,numb,mother,x,y,z,irot,konly);
3911 npar = Int_t (fZq[jin+9]);
3912 for(i=0;i<npar;i++) par[i]=fZq[jin+10+i];
3913 Gckpar (ish,npar,par);
3914 fprintf(lun,"POSP '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s' %3d\n",
3915 name,numb,mother,x,y,z,irot,konly,npar);
3917 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
3924 fprintf(lun,"END\n");
3927 //****** write down the materials and medias *****
3929 lun=fopen(filetme,"w");
3931 for(itm=1;itm<=fGcnum->ntmed;itm++) {
3932 if (iws[iadtmd+itm]>0) {
3933 jtm = fZlq[fGclink->jtmed-itm];
3934 strncpy(natmed,(char*)&fZiq[jtm+1],4);
3935 imat = Int_t (fZq[jtm+6]);
3936 jma = Int_t (fZlq[fGclink->jmate-imat]);
3938 Gfmate (imat,namate,a,z,dens,radl,absl,par,npar);
3939 fprintf(lun,"MATE %4d '%20s'%11.5E %11.5E %11.5E %11.5E %11.5E %3d\n",
3940 iomate[imat],namate,a,z,dens,radl,absl,npar);
3944 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
3948 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxfd,stemax,deemax,epsil,stmin,par,&npar);
3949 fprintf(lun,"TMED %4d '%20s' %3d %1d %3d %11.5f %11.5f %11.5f %11.5f %11.5f %11.5f %3d\n",
3950 iotmed[itm],natmed,iomate[nmat],isvol,ifield,
3951 fieldm,tmaxfd,stemax,deemax,epsil,stmin,npar);
3955 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
3961 fprintf(lun,"END\n");
3963 printf(" *** GWEUCL *** file: %s is now written out\n",filext);
3964 printf(" *** GWEUCL *** file: %s is now written out\n",filetme);
3973 //_____________________________________________________________________________
3974 void TGeant3::Streamer(TBuffer &R__b)
3977 // Stream an object of class TGeant3.
3979 if (R__b.IsReading()) {
3980 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
3981 AliMC::Streamer(R__b);
3984 R__b.ReadStaticArray(fPDGCode);
3986 R__b.WriteVersion(TGeant3::IsA());
3987 AliMC::Streamer(R__b);
3990 R__b.WriteArray(fPDGCode, fNPDGCodes);