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.32 2000/08/24 16:28:53 hristov
19 TGeant3::IsNewTrack corrected by F.Carminati
21 Revision 1.31 2000/07/13 16:19:10 fca
22 Mainly coding conventions + some small bug fixes
24 Revision 1.30 2000/07/12 08:56:30 fca
25 Coding convention correction and warning removal
27 Revision 1.29 2000/07/11 18:24:59 fca
28 Coding convention corrections + few minor bug fixes
30 Revision 1.28 2000/06/29 10:51:55 morsch
31 Add some charmed and bottom baryons to the particle list (TDatabasePDG). This
32 is needed by Hijing. Should be part of a future review of TDatabasePDG.
34 Revision 1.27 2000/06/21 17:40:15 fca
35 Adding possibility to set ISTRA, PAI model
37 Revision 1.26 2000/05/16 13:10:41 fca
38 New method IsNewTrack and fix for a problem in Father-Daughter relations
40 Revision 1.25 2000/04/07 11:12:35 fca
41 G4 compatibility changes
43 Revision 1.24 2000/02/28 21:03:57 fca
44 Some additions to improve the compatibility with G4
46 Revision 1.23 2000/02/23 16:25:25 fca
47 AliVMC and AliGeant3 classes introduced
48 ReadEuclid moved from AliRun to AliModule
50 Revision 1.22 2000/01/18 15:40:13 morsch
51 Interface to GEANT3 routines GFTMAT, GBRELM and GPRELM added
52 Define geant particle type 51: Feedback Photon with Cherenkov photon properties.
54 Revision 1.21 2000/01/17 19:41:17 fca
57 Revision 1.20 2000/01/12 11:29:27 fca
60 Revision 1.19 1999/12/17 09:03:12 fca
61 Introduce a names array
63 Revision 1.18 1999/11/26 16:55:39 fca
64 Reimplement CurrentVolName() to avoid memory leaks
66 Revision 1.17 1999/11/03 16:58:28 fca
67 Correct source of address violation in creating character string
69 Revision 1.16 1999/11/03 13:17:08 fca
70 Have ProdProcess return const char*
72 Revision 1.15 1999/10/26 06:04:50 fca
73 Introduce TLorentzVector in AliMC::GetSecondary. Thanks to I.Hrivnacova
75 Revision 1.14 1999/09/29 09:24:30 fca
76 Introduction of the Copyright and cvs Log
80 ///////////////////////////////////////////////////////////////////////////////
82 // Interface Class to the Geant3.21 MonteCarlo //
86 <img src="picts/TGeant3Class.gif">
91 ///////////////////////////////////////////////////////////////////////////////
97 #include <TDatabasePDG.h>
98 #include "AliDecayer.h"
99 #include "AliDecayerPythia.h"
101 #include "AliCallf77.h"
104 # define gzebra gzebra_
105 # define grfile grfile_
106 # define gpcxyz gpcxyz_
107 # define ggclos ggclos_
108 # define glast glast_
109 # define ginit ginit_
110 # define gcinit gcinit_
112 # define gtrig gtrig_
113 # define gtrigc gtrigc_
114 # define gtrigi gtrigi_
115 # define gwork gwork_
116 # define gzinit gzinit_
117 # define gfmate gfmate_
118 # define gfpart gfpart_
119 # define gftmed gftmed_
120 # define gftmat gftmat_
121 # define gmate gmate_
122 # define gpart gpart_
124 # define gsmate gsmate_
125 # define gsmixt gsmixt_
126 # define gspart gspart_
127 # define gstmed gstmed_
128 # define gsckov gsckov_
129 # define gstpar gstpar_
130 # define gfkine gfkine_
131 # define gfvert gfvert_
132 # define gskine gskine_
133 # define gsvert gsvert_
134 # define gphysi gphysi_
135 # define gdebug gdebug_
136 # define gekbin gekbin_
137 # define gfinds gfinds_
138 # define gsking gsking_
139 # define gskpho gskpho_
140 # define gsstak gsstak_
141 # define gsxyz gsxyz_
142 # define gtrack gtrack_
143 # define gtreve gtreve_
144 # define gtreveroot gtreveroot_
145 # define grndm grndm_
146 # define grndmq grndmq_
147 # define gdtom gdtom_
148 # define glmoth glmoth_
149 # define gmedia gmedia_
150 # define gmtod gmtod_
151 # define gsdvn gsdvn_
152 # define gsdvn2 gsdvn2_
153 # define gsdvs gsdvs_
154 # define gsdvs2 gsdvs2_
155 # define gsdvt gsdvt_
156 # define gsdvt2 gsdvt2_
157 # define gsord gsord_
158 # define gspos gspos_
159 # define gsposp gsposp_
160 # define gsrotm gsrotm_
161 # define gprotm gprotm_
162 # define gsvolu gsvolu_
163 # define gprint gprint_
164 # define gdinit gdinit_
165 # define gdopt gdopt_
166 # define gdraw gdraw_
167 # define gdrayt gdrayt_
168 # define gdrawc gdrawc_
169 # define gdrawx gdrawx_
170 # define gdhead gdhead_
171 # define gdwmn1 gdwmn1_
172 # define gdwmn2 gdwmn2_
173 # define gdwmn3 gdwmn3_
174 # define gdxyz gdxyz_
175 # define gdcxyz gdcxyz_
176 # define gdman gdman_
177 # define gdspec gdspec_
178 # define gdtree gdtree_
179 # define gdelet gdelet_
180 # define gdclos gdclos_
181 # define gdshow gdshow_
182 # define gdopen gdopen_
183 # define dzshow dzshow_
184 # define gsatt gsatt_
185 # define gfpara gfpara_
186 # define gckpar gckpar_
187 # define gckmat gckmat_
188 # define geditv geditv_
189 # define mzdrop mzdrop_
191 # define ertrak ertrak_
192 # define ertrgo ertrgo_
194 # define setbomb setbomb_
195 # define setclip setclip_
196 # define gcomad gcomad_
198 # define gbrelm gbrelm_
199 # define gprelm gprelm_
201 # define gzebra GZEBRA
202 # define grfile GRFILE
203 # define gpcxyz GPCXYZ
204 # define ggclos GGCLOS
207 # define gcinit GCINIT
210 # define gtrigc GTRIGC
211 # define gtrigi GTRIGI
213 # define gzinit GZINIT
214 # define gfmate GFMATE
215 # define gfpart GFPART
216 # define gftmed GFTMED
217 # define gftmat GFTMAT
221 # define gsmate GSMATE
222 # define gsmixt GSMIXT
223 # define gspart GSPART
224 # define gstmed GSTMED
225 # define gsckov GSCKOV
226 # define gstpar GSTPAR
227 # define gfkine GFKINE
228 # define gfvert GFVERT
229 # define gskine GSKINE
230 # define gsvert GSVERT
231 # define gphysi GPHYSI
232 # define gdebug GDEBUG
233 # define gekbin GEKBIN
234 # define gfinds GFINDS
235 # define gsking GSKING
236 # define gskpho GSKPHO
237 # define gsstak GSSTAK
239 # define gtrack GTRACK
240 # define gtreve GTREVE
241 # define gtreveroot GTREVEROOT
243 # define grndmq GRNDMQ
245 # define glmoth GLMOTH
246 # define gmedia GMEDIA
249 # define gsdvn2 GSDVN2
251 # define gsdvs2 GSDVS2
253 # define gsdvt2 GSDVT2
256 # define gsposp GSPOSP
257 # define gsrotm GSROTM
258 # define gprotm GPROTM
259 # define gsvolu GSVOLU
260 # define gprint GPRINT
261 # define gdinit GDINIT
264 # define gdrayt GDRAYT
265 # define gdrawc GDRAWC
266 # define gdrawx GDRAWX
267 # define gdhead GDHEAD
268 # define gdwmn1 GDWMN1
269 # define gdwmn2 GDWMN2
270 # define gdwmn3 GDWMN3
272 # define gdcxyz GDCXYZ
274 # define gdfspc GDFSPC
275 # define gdspec GDSPEC
276 # define gdtree GDTREE
277 # define gdelet GDELET
278 # define gdclos GDCLOS
279 # define gdshow GDSHOW
280 # define gdopen GDOPEN
281 # define dzshow DZSHOW
283 # define gfpara GFPARA
284 # define gckpar GCKPAR
285 # define gckmat GCKMAT
286 # define geditv GEDITV
287 # define mzdrop MZDROP
289 # define ertrak ERTRAK
290 # define ertrgo ERTRGO
292 # define setbomb SETBOMB
293 # define setclip SETCLIP
294 # define gcomad GCOMAD
296 # define gbrelm GBRELM
297 # define gprelm GPRELM
301 //____________________________________________________________________________
305 // Prototypes for GEANT functions
307 void type_of_call gzebra(const int&);
309 void type_of_call gpcxyz();
311 void type_of_call ggclos();
313 void type_of_call glast();
315 void type_of_call ginit();
317 void type_of_call gcinit();
319 void type_of_call grun();
321 void type_of_call gtrig();
323 void type_of_call gtrigc();
325 void type_of_call gtrigi();
327 void type_of_call gwork(const int&);
329 void type_of_call gzinit();
331 void type_of_call gmate();
333 void type_of_call gpart();
335 void type_of_call gsdk(Int_t &, Float_t *, Int_t *);
337 void type_of_call gfkine(Int_t &, Float_t *, Float_t *, Int_t &,
338 Int_t &, Float_t *, Int_t &);
340 void type_of_call gfvert(Int_t &, Float_t *, Int_t &, Int_t &,
341 Float_t &, Float_t *, Int_t &);
343 void type_of_call gskine(Float_t *,Int_t &, Int_t &, Float_t *,
346 void type_of_call gsvert(Float_t *,Int_t &, Int_t &, Float_t *,
349 void type_of_call gphysi();
351 void type_of_call gdebug();
353 void type_of_call gekbin();
355 void type_of_call gfinds();
357 void type_of_call gsking(Int_t &);
359 void type_of_call gskpho(Int_t &);
361 void type_of_call gsstak(Int_t &);
363 void type_of_call gsxyz();
365 void type_of_call gtrack();
367 void type_of_call gtreve();
369 void type_of_call gtreveroot();
371 void type_of_call grndm(Float_t *, const Int_t &);
373 void type_of_call grndmq(Int_t &, Int_t &, const Int_t &,
376 void type_of_call gdtom(Float_t *, Float_t *, Int_t &);
378 void type_of_call glmoth(DEFCHARD, Int_t &, Int_t &, Int_t *,
379 Int_t *, Int_t * DEFCHARL);
381 void type_of_call gmedia(Float_t *, Int_t &);
383 void type_of_call gmtod(Float_t *, Float_t *, Int_t &);
385 void type_of_call gsrotm(const Int_t &, const Float_t &, const Float_t &,
386 const Float_t &, const Float_t &, const Float_t &,
389 void type_of_call gprotm(const Int_t &);
391 void type_of_call grfile(const Int_t&, DEFCHARD,
392 DEFCHARD DEFCHARL DEFCHARL);
394 void type_of_call gfmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
395 Float_t &, Float_t &, Float_t &, Float_t *,
398 void type_of_call gfpart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
399 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
401 void type_of_call gftmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
402 Float_t &, Float_t &, Float_t &, Float_t &,
403 Float_t &, Float_t &, Float_t *, Int_t * DEFCHARL);
405 void type_of_call gftmat(const Int_t&, const Int_t&, DEFCHARD, const Int_t&,
407 ,Float_t *, Int_t & DEFCHARL);
409 void type_of_call gsmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
410 Float_t &, Float_t &, Float_t &, Float_t *,
413 void type_of_call gsmixt(const Int_t&, DEFCHARD, Float_t *, Float_t *,
414 Float_t &, Int_t &, Float_t * DEFCHARL);
416 void type_of_call gspart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
417 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
420 void type_of_call gstmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
421 Float_t &, Float_t &, Float_t &, Float_t &,
422 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
424 void type_of_call gsckov(Int_t &itmed, Int_t &npckov, Float_t *ppckov,
425 Float_t *absco, Float_t *effic, Float_t *rindex);
426 void type_of_call gstpar(const Int_t&, DEFCHARD, Float_t & DEFCHARL);
428 void type_of_call gsdvn(DEFCHARD,DEFCHARD, Int_t &, Int_t &
431 void type_of_call gsdvn2(DEFCHARD,DEFCHARD, Int_t &, Int_t &, Float_t &,
432 Int_t & DEFCHARL DEFCHARL);
434 void type_of_call gsdvs(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &
437 void type_of_call gsdvs2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t &,
438 Int_t & DEFCHARL DEFCHARL);
440 void type_of_call gsdvt(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &,
441 Int_t & DEFCHARL DEFCHARL);
443 void type_of_call gsdvt2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t&,
444 Int_t &, Int_t & DEFCHARL DEFCHARL);
446 void type_of_call gsord(DEFCHARD, Int_t & DEFCHARL);
448 void type_of_call gspos(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
449 Float_t &, Int_t &, DEFCHARD DEFCHARL DEFCHARL
452 void type_of_call gsposp(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
453 Float_t &, Int_t &, DEFCHARD,
454 Float_t *, Int_t & DEFCHARL DEFCHARL DEFCHARL);
456 void type_of_call gsvolu(DEFCHARD, DEFCHARD, Int_t &, Float_t *, Int_t &,
457 Int_t & DEFCHARL DEFCHARL);
459 void type_of_call gsatt(DEFCHARD, DEFCHARD, Int_t & DEFCHARL DEFCHARL);
461 void type_of_call gfpara(DEFCHARD , Int_t&, Int_t&, Int_t&, Int_t&, Float_t*,
464 void type_of_call gckpar(Int_t&, Int_t&, Float_t*);
466 void type_of_call gckmat(Int_t&, DEFCHARD DEFCHARL);
468 void type_of_call gprint(DEFCHARD,const int& DEFCHARL);
470 void type_of_call gdinit();
472 void type_of_call gdopt(DEFCHARD,DEFCHARD DEFCHARL DEFCHARL);
474 void type_of_call gdraw(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
475 Float_t &, Float_t &, Float_t & DEFCHARL);
476 void type_of_call gdrayt(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
477 Float_t &, Float_t &, Float_t & DEFCHARL);
478 void type_of_call gdrawc(DEFCHARD,Int_t &, Float_t &, Float_t &, Float_t &,
479 Float_t &, Float_t & DEFCHARL);
480 void type_of_call gdrawx(DEFCHARD,Float_t &, Float_t &, Float_t &, Float_t &,
481 Float_t &, Float_t &, Float_t &, Float_t &,
483 void type_of_call gdhead(Int_t &,DEFCHARD, Float_t & DEFCHARL);
484 void type_of_call gdxyz(Int_t &);
485 void type_of_call gdcxyz();
486 void type_of_call gdman(Float_t &, Float_t &);
487 void type_of_call gdwmn1(Float_t &, Float_t &);
488 void type_of_call gdwmn2(Float_t &, Float_t &);
489 void type_of_call gdwmn3(Float_t &, Float_t &);
490 void type_of_call gdspec(DEFCHARD DEFCHARL);
491 void type_of_call gdfspc(DEFCHARD, Int_t &, Int_t & DEFCHARL) {;}
492 void type_of_call gdtree(DEFCHARD, Int_t &, Int_t & DEFCHARL);
494 void type_of_call gdopen(Int_t &);
495 void type_of_call gdclos();
496 void type_of_call gdelet(Int_t &);
497 void type_of_call gdshow(Int_t &);
498 void type_of_call geditv(Int_t &) {;}
501 void type_of_call dzshow(DEFCHARD,const int&,const int&,DEFCHARD,const int&,
502 const int&, const int&, const int& DEFCHARL
505 void type_of_call mzdrop(Int_t&, Int_t&, DEFCHARD DEFCHARL);
507 void type_of_call setbomb(Float_t &);
508 void type_of_call setclip(DEFCHARD, Float_t &,Float_t &,Float_t &,Float_t &,
509 Float_t &, Float_t & DEFCHARL);
510 void type_of_call gcomad(DEFCHARD, Int_t*& DEFCHARL);
512 void type_of_call ertrak(const Float_t *const x1, const Float_t *const p1,
513 const Float_t *x2, const Float_t *p2,
514 const Int_t &ipa, DEFCHARD DEFCHARL);
516 void type_of_call ertrgo();
518 float type_of_call gbrelm(const Float_t &z, const Float_t& t, const Float_t& cut);
519 float type_of_call gprelm(const Float_t &z, const Float_t& t, const Float_t& cut);
523 // Geant3 global pointer
525 static const Int_t kDefSize = 600;
529 //____________________________________________________________________________
533 // Default constructor
537 //____________________________________________________________________________
538 TGeant3::TGeant3(const char *title, Int_t nwgeant)
539 :AliMC("TGeant3",title)
542 // Standard constructor for TGeant3 with ZEBRA initialisation
553 // Load Address of Geant3 commons
556 // Zero number of particles
558 fDecayer=new AliDecayerPythia();
562 //____________________________________________________________________________
563 Int_t TGeant3::CurrentMaterial(Float_t &a, Float_t &z, Float_t &dens,
564 Float_t &radl, Float_t &absl) const
567 // Return the parameters of the current material during transport
571 dens = fGcmate->dens;
572 radl = fGcmate->radl;
573 absl = fGcmate->absl;
574 return 1; //this could be the number of elements in mixture
577 //____________________________________________________________________________
578 void TGeant3::DefaultRange()
581 // Set range of current drawing pad to 20x20 cm
587 gHigz->Range(0,0,20,20);
590 //____________________________________________________________________________
591 void TGeant3::InitHIGZ()
602 //____________________________________________________________________________
603 void TGeant3::LoadAddress()
606 // Assigns the address of the GEANT common blocks to the structures
607 // that allow their access from C++
610 gcomad(PASSCHARD("QUEST"), (int*&) fQuest PASSCHARL("QUEST"));
611 gcomad(PASSCHARD("GCBANK"),(int*&) fGcbank PASSCHARL("GCBANK"));
612 gcomad(PASSCHARD("GCLINK"),(int*&) fGclink PASSCHARL("GCLINK"));
613 gcomad(PASSCHARD("GCCUTS"),(int*&) fGccuts PASSCHARL("GCCUTS"));
614 gcomad(PASSCHARD("GCMULO"),(int*&) fGcmulo PASSCHARL("GCMULO"));
615 gcomad(PASSCHARD("GCFLAG"),(int*&) fGcflag PASSCHARL("GCFLAG"));
616 gcomad(PASSCHARD("GCKINE"),(int*&) fGckine PASSCHARL("GCKINE"));
617 gcomad(PASSCHARD("GCKING"),(int*&) fGcking PASSCHARL("GCKING"));
618 gcomad(PASSCHARD("GCKIN2"),(int*&) fGckin2 PASSCHARL("GCKIN2"));
619 gcomad(PASSCHARD("GCKIN3"),(int*&) fGckin3 PASSCHARL("GCKIN3"));
620 gcomad(PASSCHARD("GCMATE"),(int*&) fGcmate PASSCHARL("GCMATE"));
621 gcomad(PASSCHARD("GCTMED"),(int*&) fGctmed PASSCHARL("GCTMED"));
622 gcomad(PASSCHARD("GCTRAK"),(int*&) fGctrak PASSCHARL("GCTRAK"));
623 gcomad(PASSCHARD("GCTPOL"),(int*&) fGctpol PASSCHARL("GCTPOL"));
624 gcomad(PASSCHARD("GCVOLU"),(int*&) fGcvolu PASSCHARL("GCVOLU"));
625 gcomad(PASSCHARD("GCNUM"), (int*&) fGcnum PASSCHARL("GCNUM"));
626 gcomad(PASSCHARD("GCSETS"),(int*&) fGcsets PASSCHARL("GCSETS"));
627 gcomad(PASSCHARD("GCPHYS"),(int*&) fGcphys PASSCHARL("GCPHYS"));
628 gcomad(PASSCHARD("GCPHLT"),(int*&) fGcphlt PASSCHARL("GCPHLT"));
629 gcomad(PASSCHARD("GCOPTI"),(int*&) fGcopti PASSCHARL("GCOPTI"));
630 gcomad(PASSCHARD("GCTLIT"),(int*&) fGctlit PASSCHARL("GCTLIT"));
631 gcomad(PASSCHARD("GCVDMA"),(int*&) fGcvdma PASSCHARL("GCVDMA"));
634 gcomad(PASSCHARD("ERTRIO"),(int*&) fErtrio PASSCHARL("ERTRIO"));
635 gcomad(PASSCHARD("EROPTS"),(int*&) fEropts PASSCHARL("EROPTS"));
636 gcomad(PASSCHARD("EROPTC"),(int*&) fEroptc PASSCHARL("EROPTC"));
637 gcomad(PASSCHARD("ERWORK"),(int*&) fErwork PASSCHARL("ERWORK"));
639 // Variables for ZEBRA store
640 gcomad(PASSCHARD("IQ"), addr PASSCHARL("IQ"));
642 gcomad(PASSCHARD("LQ"), addr PASSCHARL("LQ"));
647 //_____________________________________________________________________________
648 void TGeant3::GeomIter()
651 // Geometry iterator for moving upward in the geometry tree
652 // Initialise the iterator
654 fNextVol=fGcvolu->nlevel;
657 //____________________________________________________________________________
658 void TGeant3::FinishGeometry()
660 //Close the geometry structure
664 //____________________________________________________________________________
665 Int_t TGeant3::NextVolUp(Text_t *name, Int_t ©)
668 // Geometry iterator for moving upward in the geometry tree
669 // Return next volume up
674 gname=fGcvolu->names[fNextVol];
675 copy=fGcvolu->number[fNextVol];
676 i=fGcvolu->lvolum[fNextVol];
677 name = fVolNames[i-1];
678 if(gname == fZiq[fGclink->jvolum+i]) return i;
679 else printf("GeomTree: Volume %s not found in bank\n",name);
684 //_____________________________________________________________________________
685 void TGeant3::BuildPhysics()
690 //_____________________________________________________________________________
691 Int_t TGeant3::CurrentVolID(Int_t ©) const
694 // Returns the current volume ID and copy number
697 if( (i=fGcvolu->nlevel-1) < 0 ) {
698 Warning("CurrentVolID","Stack depth only %d\n",fGcvolu->nlevel);
700 gname=fGcvolu->names[i];
701 copy=fGcvolu->number[i];
702 i=fGcvolu->lvolum[i];
703 if(gname == fZiq[fGclink->jvolum+i]) return i;
704 else Warning("CurrentVolID","Volume %4s not found\n",(char*)&gname);
709 //_____________________________________________________________________________
710 Int_t TGeant3::CurrentVolOffID(Int_t off, Int_t ©) const
713 // Return the current volume "off" upward in the geometrical tree
714 // ID and copy number
717 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
718 Warning("CurrentVolOffID","Offset requested %d but stack depth %d\n",
719 off,fGcvolu->nlevel);
721 gname=fGcvolu->names[i];
722 copy=fGcvolu->number[i];
723 i=fGcvolu->lvolum[i];
724 if(gname == fZiq[fGclink->jvolum+i]) return i;
725 else Warning("CurrentVolOffID","Volume %4s not found\n",(char*)&gname);
730 //_____________________________________________________________________________
731 const char* TGeant3::CurrentVolName() const
734 // Returns the current volume name
737 if( (i=fGcvolu->nlevel-1) < 0 ) {
738 Warning("CurrentVolName","Stack depth %d\n",fGcvolu->nlevel);
740 gname=fGcvolu->names[i];
741 i=fGcvolu->lvolum[i];
742 if(gname == fZiq[fGclink->jvolum+i]) return fVolNames[i-1];
743 else Warning("CurrentVolName","Volume %4s not found\n",(char*) &gname);
748 //_____________________________________________________________________________
749 const char* TGeant3::CurrentVolOffName(Int_t off) const
752 // Return the current volume "off" upward in the geometrical tree
753 // ID, name and copy number
754 // if name=0 no name is returned
757 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
758 Warning("CurrentVolOffName",
759 "Offset requested %d but stack depth %d\n",off,fGcvolu->nlevel);
761 gname=fGcvolu->names[i];
762 i=fGcvolu->lvolum[i];
763 if(gname == fZiq[fGclink->jvolum+i]) return fVolNames[i-1];
764 else Warning("CurrentVolOffName","Volume %4s not found\n",(char*)&gname);
769 //_____________________________________________________________________________
770 Int_t TGeant3::IdFromPDG(Int_t pdg) const
773 // Return Geant3 code from PDG and pseudo ENDF code
775 for(Int_t i=0;i<fNPDGCodes;++i)
776 if(pdg==fPDGCode[i]) return i;
780 //_____________________________________________________________________________
781 Int_t TGeant3::PDGFromId(Int_t id) const
784 // Return PDG code and pseudo ENDF code from Geant3 code
786 if(id>0 && id<fNPDGCodes) return fPDGCode[id];
790 //_____________________________________________________________________________
791 void TGeant3::DefineParticles()
794 // Define standard Geant 3 particles
797 // Load standard numbers for GEANT particles and PDG conversion
798 fPDGCode[fNPDGCodes++]=-99; // 0 = unused location
799 fPDGCode[fNPDGCodes++]=22; // 1 = photon
800 fPDGCode[fNPDGCodes++]=-11; // 2 = positron
801 fPDGCode[fNPDGCodes++]=11; // 3 = electron
802 fPDGCode[fNPDGCodes++]=12; // 4 = neutrino e
803 fPDGCode[fNPDGCodes++]=-13; // 5 = muon +
804 fPDGCode[fNPDGCodes++]=13; // 6 = muon -
805 fPDGCode[fNPDGCodes++]=111; // 7 = pi0
806 fPDGCode[fNPDGCodes++]=211; // 8 = pi+
807 fPDGCode[fNPDGCodes++]=-211; // 9 = pi-
808 fPDGCode[fNPDGCodes++]=130; // 10 = Kaon Long
809 fPDGCode[fNPDGCodes++]=321; // 11 = Kaon +
810 fPDGCode[fNPDGCodes++]=-321; // 12 = Kaon -
811 fPDGCode[fNPDGCodes++]=2112; // 13 = Neutron
812 fPDGCode[fNPDGCodes++]=2212; // 14 = Proton
813 fPDGCode[fNPDGCodes++]=-2212; // 15 = Anti Proton
814 fPDGCode[fNPDGCodes++]=310; // 16 = Kaon Short
815 fPDGCode[fNPDGCodes++]=221; // 17 = Eta
816 fPDGCode[fNPDGCodes++]=3122; // 18 = Lambda
817 fPDGCode[fNPDGCodes++]=3222; // 19 = Sigma +
818 fPDGCode[fNPDGCodes++]=3212; // 20 = Sigma 0
819 fPDGCode[fNPDGCodes++]=3112; // 21 = Sigma -
820 fPDGCode[fNPDGCodes++]=3322; // 22 = Xi0
821 fPDGCode[fNPDGCodes++]=3312; // 23 = Xi-
822 fPDGCode[fNPDGCodes++]=3334; // 24 = Omega-
823 fPDGCode[fNPDGCodes++]=-2112; // 25 = Anti Proton
824 fPDGCode[fNPDGCodes++]=-3122; // 26 = Anti Proton
825 fPDGCode[fNPDGCodes++]=-3222; // 27 = Anti Sigma -
826 fPDGCode[fNPDGCodes++]=-3212; // 28 = Anti Sigma 0
827 fPDGCode[fNPDGCodes++]=-3112; // 29 = Anti Sigma 0
828 fPDGCode[fNPDGCodes++]=-3322; // 30 = Anti Xi 0
829 fPDGCode[fNPDGCodes++]=-3312; // 31 = Anti Xi +
830 fPDGCode[fNPDGCodes++]=-3334; // 32 = Anti Omega +
837 /* --- Define additional particles */
838 Gspart(33, "OMEGA(782)", 3, 0.782, 0., 7.836e-23);
839 fPDGCode[fNPDGCodes++]=223; // 33 = Omega(782)
841 Gspart(34, "PHI(1020)", 3, 1.019, 0., 1.486e-22);
842 fPDGCode[fNPDGCodes++]=333; // 34 = PHI (1020)
844 Gspart(35, "D +", 4, 1.87, 1., 1.066e-12);
845 fPDGCode[fNPDGCodes++]=411; // 35 = D+
847 Gspart(36, "D -", 4, 1.87, -1., 1.066e-12);
848 fPDGCode[fNPDGCodes++]=-411; // 36 = D-
850 Gspart(37, "D 0", 3, 1.865, 0., 4.2e-13);
851 fPDGCode[fNPDGCodes++]=421; // 37 = D0
853 Gspart(38, "ANTI D 0", 3, 1.865, 0., 4.2e-13);
854 fPDGCode[fNPDGCodes++]=-421; // 38 = D0 bar
856 fPDGCode[fNPDGCodes++]=-99; // 39 = unassigned
858 fPDGCode[fNPDGCodes++]=-99; // 40 = unassigned
860 fPDGCode[fNPDGCodes++]=-99; // 41 = unassigned
862 Gspart(42, "RHO +", 4, 0.768, 1., 4.353e-24);
863 fPDGCode[fNPDGCodes++]=213; // 42 = RHO+
865 Gspart(43, "RHO -", 4, 0.768, -1., 4.353e-24);
866 fPDGCode[fNPDGCodes++]=-213; // 43 = RHO-
868 Gspart(44, "RHO 0", 3, 0.768, 0., 4.353e-24);
869 fPDGCode[fNPDGCodes++]=113; // 44 = RHO0
872 // Use ENDF-6 mapping for ions = 10000*z+10*a+iso
874 // and numbers above 5 000 000 for special applications
877 const Int_t kion=10000000;
879 const Int_t kspe=50000000;
881 TDatabasePDG *pdgDB = TDatabasePDG::Instance();
883 const Double_t kAu2Gev=0.9314943228;
884 const Double_t khSlash = 1.0545726663e-27;
885 const Double_t kErg2Gev = 1/1.6021773349e-3;
886 const Double_t khShGev = khSlash*kErg2Gev;
887 const Double_t kYear2Sec = 3600*24*365.25;
890 // mass and life-time from PDG
891 pdgDB->AddParticle("B(s)*0","B(s)*0",
892 5.4163, kTRUE, 0.047, +0.,"Meson", 533);
894 pdgDB->AddParticle("B(s)*0 bar","B(s)*0 bar",
895 5.4163, kTRUE, 0.047, -0.,"Meson", -533);
899 // value for mass used by Hijing
900 pdgDB->AddParticle("Sigma(c)*+","Sigma(c)*+",
901 2.4536, kTRUE, -1., +1.,"Baryon", 4214);
903 pdgDB->AddParticle("Sigma(c)*-","Sigma(c)*-",
904 2.4536, kTRUE, -1., -1.,"Baryon", -4214);
905 // equivalent to 4312 ? Hijing uses m=2.55
906 pdgDB->AddParticle("Xsi(c)0","Xsi(c)0",
907 2.4703, kTRUE, -1., +0.,"Baryon", 4132);
909 pdgDB->AddParticle("Xsi(c)0 bar","Xsi(c)0 bar",
910 2.4703, kTRUE, -1., -0.,"Baryon", -4132);
911 // equivalent to 4322 ? Hijing uses m=2.55
912 pdgDB->AddParticle("Xi(c)+","Xi(c)+",
913 2.4656, kFALSE, -1., +1.,"Baryon", 4232);
915 pdgDB->AddParticle("Xi(c)-","Xi(c)-",
916 2.4656, kFALSE, -1., -1.,"Baryon", -4232);
917 // mass values from Hijing
919 pdgDB->AddParticle("Xsi(c)*0","Xsi(c)*0",
920 2.63, kTRUE, -1., +0.,"Baryon", 4314);
922 pdgDB->AddParticle("Xsi(c)*0 bar","Xsi(c)*0 bar",
923 2.63, kTRUE, -1., -0.,"Baryon", -4314);
925 pdgDB->AddParticle("Xsi(c)*+","Xsi(c)*+",
926 2.63, kTRUE, -1., +1.,"Baryon", 4324);
928 pdgDB->AddParticle("Xsi(c)*-","Xsi(c)*-",
929 2.63, kTRUE, -1., -1.,"Baryon", -4324);
931 // pdg mass value, Hijing uses m=2.73.
932 pdgDB->AddParticle("Omega(c)0","Omega(c)0",
933 2.7040, kFALSE, khShGev/0.064e-12, +0.,"Baryon", 4332);
935 pdgDB->AddParticle("Omega(c)0 bar","Omega(c)0 bar",
936 2.7040, kFALSE, khShGev/0.064e-12, -0.,"Baryon", -4332);
937 // mass value from Hijing
938 pdgDB->AddParticle("Omega(c)*0","Omega(c)*0",
939 2.8000, kFALSE, -1., +0.,"Baryon", 4334);
941 pdgDB->AddParticle("Omega(c)*0 bar","Omega(c)*0",
942 2.8000, kFALSE, -1., -0.,"Baryon", -4334);
946 // mass value from Hijing
947 pdgDB->AddParticle("Sigma(b)*+","Sigma(b)*+",
948 5.8100, kFALSE, -1., +1.,"Baryon", 5224);
950 pdgDB->AddParticle("Sigma(b)*-","Sigma(b)*-",
951 5.8100, kFALSE, -1., -1.,"Baryon", -5224);
955 pdgDB->AddParticle("Deuteron","Deuteron",2*kAu2Gev+8.071e-3,kTRUE,
956 0,1,"Ion",kion+10020);
957 fPDGCode[fNPDGCodes++]=kion+10020; // 45 = Deuteron
959 pdgDB->AddParticle("Triton","Triton",3*kAu2Gev+14.931e-3,kFALSE,
960 khShGev/(12.33*kYear2Sec),1,"Ion",kion+10030);
961 fPDGCode[fNPDGCodes++]=kion+10030; // 46 = Triton
963 pdgDB->AddParticle("Alpha","Alpha",4*kAu2Gev+2.424e-3,kTRUE,
964 khShGev/(12.33*kYear2Sec),2,"Ion",kion+20040);
965 fPDGCode[fNPDGCodes++]=kion+20040; // 47 = Alpha
967 fPDGCode[fNPDGCodes++]=0; // 48 = geantino mapped to rootino
969 pdgDB->AddParticle("HE3","HE3",3*kAu2Gev+14.931e-3,kFALSE,
970 0,2,"Ion",kion+20030);
971 fPDGCode[fNPDGCodes++]=kion+20030; // 49 = HE3
973 pdgDB->AddParticle("Cherenkov","Cherenkov",0,kFALSE,
974 0,0,"Special",kspe+50);
975 fPDGCode[fNPDGCodes++]=kspe+50; // 50 = Cherenkov
977 Gspart(51, "FeedbackPhoton", 7, 0., 0.,1.e20 );
978 pdgDB->AddParticle("FeedbackPhoton","FeedbackPhoton",0,kFALSE,
979 0,0,"Special",kspe+51);
980 fPDGCode[fNPDGCodes++]=kspe+51; // 51 = FeedbackPhoton
981 Gspart(52, "Lambda_c+", 4, 2.2849, +1., 2.06e-13);
982 fPDGCode[fNPDGCodes++]=4122; //52 = Lambda_c+
984 Gspart(53, "Lambda_c-", 4, 2.2849, -1., 2.06e-13);
985 fPDGCode[fNPDGCodes++]=-4122; //53 = Lambda_c-
987 Gspart(54, "D_s+", 4, 1.9685, +1., 4.67e-13);
988 fPDGCode[fNPDGCodes++]=431; //54 = D_s+
990 Gspart(55, "D_s-", 4, 1.9685, -1., 4.67e-13);
991 fPDGCode[fNPDGCodes++]=-431; //55 = D_s-
993 Gspart(56, "Tau+", 5, 1.77705, +1., 2.9e-13);
994 fPDGCode[fNPDGCodes++]=15; //56 = Tau+
996 Gspart(57, "Tau-", 5, 1.77705, -1., 2.9e-13);
997 fPDGCode[fNPDGCodes++]=-15; //57 = Tau-
999 Gspart(58, "B0", 3, 5.2792, +0., 1.56e-12);
1000 fPDGCode[fNPDGCodes++]=511; //58 = B0
1002 Gspart(59, "B0 bar", 3, 5.2792, -0., 1.56e-12);
1003 fPDGCode[fNPDGCodes++]=-511; //58 = B0bar
1005 Gspart(60, "B+", 4, 5.2789, +1., 1.65e-12);
1006 fPDGCode[fNPDGCodes++]=521; //60 = B+
1008 Gspart(61, "B-", 4, 5.2789, -1., 1.65e-12);
1009 fPDGCode[fNPDGCodes++]=-521; //61 = B-
1011 Gspart(62, "Bs", 3, 5.3693, +0., 1.54e-12);
1012 fPDGCode[fNPDGCodes++]=521; //62 = B_s
1014 Gspart(63, "Bs bar", 3, 5.3693, -0., 1.54e-12);
1015 fPDGCode[fNPDGCodes++]=-521; //63 = B_s bar
1017 Gspart(64, "Lambda_b", 3, 5.624, +0., 1.24e-12);
1018 fPDGCode[fNPDGCodes++]=5122; //64 = Lambda_b
1020 Gspart(65, "Lambda_b bar", 3, 5.624, -0., 1.24e-12);
1021 fPDGCode[fNPDGCodes++]=-5122; //65 = Lambda_b bar
1023 Gspart(66, "J/Psi", 3.09688, 3, 0., 0.);
1024 fPDGCode[fNPDGCodes++]=443; // 66 = J/Psi
1026 Gspart(67, "Psi Prime", 3, 3.686, 0., 0.);
1027 fPDGCode[fNPDGCodes++]=20443; // 67 = Psi prime
1029 Gspart(68, "Upsilon(1S)", 9.46037, 3, 0., 0.);
1030 fPDGCode[fNPDGCodes++]=553; // 68 = Upsilon(1S)
1032 Gspart(69, "Upsilon(2S)", 10.0233, 3, 0., 0.);
1033 fPDGCode[fNPDGCodes++]=20553; // 69 = Upsilon(2S)
1035 Gspart(70, "Upsilon(3S)", 10.3553, 3, 0., 0.);
1036 fPDGCode[fNPDGCodes++]=30553; // 70 = Upsilon(3S)
1038 /* --- Define additional decay modes --- */
1039 /* --- omega(783) --- */
1040 for (kz = 0; kz < 6; ++kz) {
1051 Gsdk(ipa, bratio, mode);
1052 /* --- phi(1020) --- */
1053 for (kz = 0; kz < 6; ++kz) {
1068 Gsdk(ipa, bratio, mode);
1071 for (kz = 0; kz < 6; ++kz) {
1084 Gsdk(ipa, bratio, mode);
1087 for (kz = 0; kz < 6; ++kz) {
1100 Gsdk(ipa, bratio, mode);
1104 for (kz = 0; kz < 6; ++kz) {
1115 Gsdk(ipa, bratio, mode);
1117 /* --- Anti D0 --- */
1119 for (kz = 0; kz < 6; ++kz) {
1130 Gsdk(ipa, bratio, mode);
1133 for (kz = 0; kz < 6; ++kz) {
1140 Gsdk(ipa, bratio, mode);
1142 for (kz = 0; kz < 6; ++kz) {
1149 Gsdk(ipa, bratio, mode);
1151 for (kz = 0; kz < 6; ++kz) {
1158 Gsdk(ipa, bratio, mode);
1161 for (kz = 0; kz < 6; ++kz) {
1170 Gsdk(ipa, bratio, mode);
1173 Gsdk(ipa, bratio, mode);
1176 Gsdk(ipa, bratio, mode);
1181 //_____________________________________________________________________________
1182 Int_t TGeant3::VolId(const Text_t *name) const
1185 // Return the unique numeric identifier for volume name
1188 strncpy((char *) &gname, name, 4);
1189 for(i=1; i<=fGcnum->nvolum; i++)
1190 if(gname == fZiq[fGclink->jvolum+i]) return i;
1191 printf("VolId: Volume %s not found\n",name);
1195 //_____________________________________________________________________________
1196 Int_t TGeant3::NofVolumes() const
1199 // Return total number of volumes in the geometry
1201 return fGcnum->nvolum;
1204 //_____________________________________________________________________________
1205 const char* TGeant3::VolName(Int_t id) const
1208 // Return the volume name given the volume identifier
1210 if(id<1 || id > fGcnum->nvolum || fGclink->jvolum<=0)
1211 return fVolNames[fGcnum->nvolum];
1213 return fVolNames[id-1];
1216 //_____________________________________________________________________________
1217 void TGeant3::SetCut(const char* cutName, Float_t cutValue)
1220 // Set transport cuts for particles
1222 if(!strcmp(cutName,"CUTGAM"))
1223 fGccuts->cutgam=cutValue;
1224 else if(!strcmp(cutName,"CUTGAM"))
1225 fGccuts->cutele=cutValue;
1226 else if(!strcmp(cutName,"CUTELE"))
1227 fGccuts->cutneu=cutValue;
1228 else if(!strcmp(cutName,"CUTHAD"))
1229 fGccuts->cuthad=cutValue;
1230 else if(!strcmp(cutName,"CUTMUO"))
1231 fGccuts->cutmuo=cutValue;
1232 else if(!strcmp(cutName,"BCUTE"))
1233 fGccuts->bcute=cutValue;
1234 else if(!strcmp(cutName,"BCUTM"))
1235 fGccuts->bcutm=cutValue;
1236 else if(!strcmp(cutName,"DCUTE"))
1237 fGccuts->dcute=cutValue;
1238 else if(!strcmp(cutName,"DCUTM"))
1239 fGccuts->dcutm=cutValue;
1240 else if(!strcmp(cutName,"PPCUTM"))
1241 fGccuts->ppcutm=cutValue;
1242 else if(!strcmp(cutName,"TOFMAX"))
1243 fGccuts->tofmax=cutValue;
1244 else Warning("SetCut","Cut %s not implemented\n",cutName);
1247 //_____________________________________________________________________________
1248 void TGeant3::SetProcess(const char* flagName, Int_t flagValue)
1251 // Set thresholds for different processes
1253 if(!strcmp(flagName,"PAIR"))
1254 fGcphys->ipair=flagValue;
1255 else if(!strcmp(flagName,"COMP"))
1256 fGcphys->icomp=flagValue;
1257 else if(!strcmp(flagName,"PHOT"))
1258 fGcphys->iphot=flagValue;
1259 else if(!strcmp(flagName,"PFIS"))
1260 fGcphys->ipfis=flagValue;
1261 else if(!strcmp(flagName,"DRAY"))
1262 fGcphys->idray=flagValue;
1263 else if(!strcmp(flagName,"ANNI"))
1264 fGcphys->ianni=flagValue;
1265 else if(!strcmp(flagName,"BREM"))
1266 fGcphys->ibrem=flagValue;
1267 else if(!strcmp(flagName,"HADR"))
1268 fGcphys->ihadr=flagValue;
1269 else if(!strcmp(flagName,"MUNU"))
1270 fGcphys->imunu=flagValue;
1271 else if(!strcmp(flagName,"DCAY"))
1272 fGcphys->idcay=flagValue;
1273 else if(!strcmp(flagName,"LOSS"))
1274 fGcphys->iloss=flagValue;
1275 else if(!strcmp(flagName,"MULS"))
1276 fGcphys->imuls=flagValue;
1277 else if(!strcmp(flagName,"RAYL"))
1278 fGcphys->irayl=flagValue;
1279 else if(!strcmp(flagName,"STRA"))
1280 fGcphlt->istra=flagValue;
1281 else if(!strcmp(flagName,"SYNC"))
1282 fGcphlt->isync=flagValue;
1283 else Warning("SetFlag","Flag %s not implemented\n",flagName);
1286 //_____________________________________________________________________________
1287 Float_t TGeant3::Xsec(char* reac, Float_t /* energy */,
1288 Int_t part, Int_t /* mate */)
1291 // Calculate X-sections -- dummy for the moment
1293 if(!strcmp(reac,"PHOT"))
1296 Error("Xsec","Can calculate photoelectric only for photons\n");
1302 //_____________________________________________________________________________
1303 void TGeant3::TrackPosition(TLorentzVector &xyz) const
1306 // Return the current position in the master reference frame of the
1307 // track being transported
1309 xyz[0]=fGctrak->vect[0];
1310 xyz[1]=fGctrak->vect[1];
1311 xyz[2]=fGctrak->vect[2];
1312 xyz[3]=fGctrak->tofg;
1315 //_____________________________________________________________________________
1316 Float_t TGeant3::TrackTime() const
1319 // Return the current time of flight of the track being transported
1321 return fGctrak->tofg;
1324 //_____________________________________________________________________________
1325 void TGeant3::TrackMomentum(TLorentzVector &xyz) const
1328 // Return the direction and the momentum (GeV/c) of the track
1329 // currently being transported
1331 Double_t ptot=fGctrak->vect[6];
1332 xyz[0]=fGctrak->vect[3]*ptot;
1333 xyz[1]=fGctrak->vect[4]*ptot;
1334 xyz[2]=fGctrak->vect[5]*ptot;
1335 xyz[3]=fGctrak->getot;
1338 //_____________________________________________________________________________
1339 Float_t TGeant3::TrackCharge() const
1342 // Return charge of the track currently transported
1344 return fGckine->charge;
1347 //_____________________________________________________________________________
1348 Float_t TGeant3::TrackMass() const
1351 // Return the mass of the track currently transported
1353 return fGckine->amass;
1356 //_____________________________________________________________________________
1357 Int_t TGeant3::TrackPid() const
1360 // Return the id of the particle transported
1362 return PDGFromId(fGckine->ipart);
1365 //_____________________________________________________________________________
1366 Float_t TGeant3::TrackStep() const
1369 // Return the length in centimeters of the current step
1371 return fGctrak->step;
1374 //_____________________________________________________________________________
1375 Float_t TGeant3::TrackLength() const
1378 // Return the length of the current track from its origin
1380 return fGctrak->sleng;
1383 //_____________________________________________________________________________
1384 Bool_t TGeant3::IsNewTrack() const
1387 // True if the track is not at the boundary of the current volume
1389 return (fGctrak->sleng==0);
1392 //_____________________________________________________________________________
1393 Bool_t TGeant3::IsTrackInside() const
1396 // True if the track is not at the boundary of the current volume
1398 return (fGctrak->inwvol==0);
1401 //_____________________________________________________________________________
1402 Bool_t TGeant3::IsTrackEntering() const
1405 // True if this is the first step of the track in the current volume
1407 return (fGctrak->inwvol==1);
1410 //_____________________________________________________________________________
1411 Bool_t TGeant3::IsTrackExiting() const
1414 // True if this is the last step of the track in the current volume
1416 return (fGctrak->inwvol==2);
1419 //_____________________________________________________________________________
1420 Bool_t TGeant3::IsTrackOut() const
1423 // True if the track is out of the setup
1425 return (fGctrak->inwvol==3);
1428 //_____________________________________________________________________________
1429 Bool_t TGeant3::IsTrackStop() const
1432 // True if the track energy has fallen below the threshold
1434 return (fGctrak->istop==2);
1437 //_____________________________________________________________________________
1438 Int_t TGeant3::NSecondaries() const
1441 // Number of secondary particles generated in the current step
1443 return fGcking->ngkine;
1446 //_____________________________________________________________________________
1447 Int_t TGeant3::CurrentEvent() const
1450 // Number of the current event
1452 return fGcflag->idevt;
1455 //_____________________________________________________________________________
1456 const char* TGeant3::ProdProcess() const
1459 // Name of the process that has produced the secondary particles
1460 // in the current step
1462 static char proc[5];
1463 const Int_t kIpMec[13] = { 5,6,7,8,9,10,11,12,21,23,25,105,108 };
1466 if(fGcking->ngkine>0) {
1467 for (km = 0; km < fGctrak->nmec; ++km) {
1468 for (im = 0; im < 13; ++im) {
1469 if (fGctrak->lmec[km] == kIpMec[im]) {
1470 mec = fGctrak->lmec[km];
1471 if (0 < mec && mec < 31) {
1472 strncpy(proc,(char *)&fGctrak->namec[mec - 1],4);
1473 } else if (mec - 100 <= 30 && mec - 100 > 0) {
1474 strncpy(proc,(char *)&fGctpol->namec1[mec - 101],4);
1481 strcpy(proc,"UNKN");
1482 } else strcpy(proc,"NONE");
1486 //_____________________________________________________________________________
1487 void TGeant3::GetSecondary(Int_t isec, Int_t& ipart,
1488 TLorentzVector &x, TLorentzVector &p)
1491 // Get the parameters of the secondary track number isec produced
1492 // in the current step
1495 if(-1<isec && isec<fGcking->ngkine) {
1496 ipart=Int_t (fGcking->gkin[isec][4] +0.5);
1498 x[i]=fGckin3->gpos[isec][i];
1499 p[i]=fGcking->gkin[isec][i];
1501 x[3]=fGcking->tofd[isec];
1502 p[3]=fGcking->gkin[isec][3];
1504 printf(" * TGeant3::GetSecondary * Secondary %d does not exist\n",isec);
1505 x[0]=x[1]=x[2]=x[3]=p[0]=p[1]=p[2]=p[3]=0;
1510 //_____________________________________________________________________________
1511 void TGeant3::InitLego()
1514 // Set switches for lego transport
1517 SetDEBU(0,0,0); //do not print a message
1520 //_____________________________________________________________________________
1521 Bool_t TGeant3::IsTrackDisappeared() const
1524 // True if the current particle has disappered
1525 // either because it decayed or because it underwent
1526 // an inelastic collision
1528 return (fGctrak->istop==1);
1531 //_____________________________________________________________________________
1532 Bool_t TGeant3::IsTrackAlive() const
1535 // True if the current particle is alive and will continue to be
1538 return (fGctrak->istop==0);
1541 //_____________________________________________________________________________
1542 void TGeant3::StopTrack()
1545 // Stop the transport of the current particle and skip to the next
1550 //_____________________________________________________________________________
1551 void TGeant3::StopEvent()
1554 // Stop simulation of the current event and skip to the next
1559 //_____________________________________________________________________________
1560 Float_t TGeant3::MaxStep() const
1563 // Return the maximum step length in the current medium
1565 return fGctmed->stemax;
1568 //_____________________________________________________________________________
1569 void TGeant3::SetMaxStep(Float_t maxstep)
1572 // Set the maximum step allowed till the particle is in the current medium
1574 fGctmed->stemax=maxstep;
1577 //_____________________________________________________________________________
1578 void TGeant3::SetMaxNStep(Int_t maxnstp)
1581 // Set the maximum number of steps till the particle is in the current medium
1583 fGctrak->maxnst=maxnstp;
1586 //_____________________________________________________________________________
1587 Int_t TGeant3::GetMaxNStep() const
1590 // Maximum number of steps allowed in current medium
1592 return fGctrak->maxnst;
1595 //_____________________________________________________________________________
1596 void TGeant3::Material(Int_t& kmat, const char* name, Float_t a, Float_t z,
1597 Float_t dens, Float_t radl, Float_t absl, Float_t* buf,
1601 // Defines a Material
1603 // kmat number assigned to the material
1604 // name material name
1605 // a atomic mass in au
1607 // dens density in g/cm3
1608 // absl absorbtion length in cm
1609 // if >=0 it is ignored and the program
1610 // calculates it, if <0. -absl is taken
1611 // radl radiation length in cm
1612 // if >=0 it is ignored and the program
1613 // calculates it, if <0. -radl is taken
1614 // buf pointer to an array of user words
1615 // nbuf number of user words
1617 Int_t jmate=fGclink->jmate;
1623 for(i=1; i<=ns; i++) {
1624 if(fZlq[jmate-i]==0) {
1630 gsmate(kmat,PASSCHARD(name), a, z, dens, radl, absl, buf,
1631 nwbuf PASSCHARL(name));
1634 //_____________________________________________________________________________
1635 void TGeant3::Mixture(Int_t& kmat, const char* name, Float_t* a, Float_t* z,
1636 Float_t dens, Int_t nlmat, Float_t* wmat)
1639 // Defines mixture OR COMPOUND IMAT as composed by
1640 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
1642 // If NLMAT > 0 then wmat contains the proportion by
1643 // weights of each basic material in the mixture.
1645 // If nlmat < 0 then WMAT contains the number of atoms
1646 // of a given kind into the molecule of the COMPOUND
1647 // In this case, WMAT in output is changed to relative
1650 Int_t jmate=fGclink->jmate;
1656 for(i=1; i<=ns; i++) {
1657 if(fZlq[jmate-i]==0) {
1663 gsmixt(kmat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
1666 //_____________________________________________________________________________
1667 void TGeant3::Medium(Int_t& kmed, const char* name, Int_t nmat, Int_t isvol,
1668 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
1669 Float_t stemax, Float_t deemax, Float_t epsil,
1670 Float_t stmin, Float_t* ubuf, Int_t nbuf)
1673 // kmed tracking medium number assigned
1674 // name tracking medium name
1675 // nmat material number
1676 // isvol sensitive volume flag
1677 // ifield magnetic field
1678 // fieldm max. field value (kilogauss)
1679 // tmaxfd max. angle due to field (deg/step)
1680 // stemax max. step allowed
1681 // deemax max. fraction of energy lost in a step
1682 // epsil tracking precision (cm)
1683 // stmin min. step due to continuos processes (cm)
1685 // ifield = 0 if no magnetic field; ifield = -1 if user decision in guswim;
1686 // ifield = 1 if tracking performed with grkuta; ifield = 2 if tracking
1687 // performed with ghelix; ifield = 3 if tracking performed with ghelx3.
1689 Int_t jtmed=fGclink->jtmed;
1695 for(i=1; i<=ns; i++) {
1696 if(fZlq[jtmed-i]==0) {
1702 gstmed(kmed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1703 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1706 //_____________________________________________________________________________
1707 void TGeant3::Matrix(Int_t& krot, Float_t thex, Float_t phix, Float_t they,
1708 Float_t phiy, Float_t thez, Float_t phiz)
1711 // krot rotation matrix number assigned
1712 // theta1 polar angle for axis i
1713 // phi1 azimuthal angle for axis i
1714 // theta2 polar angle for axis ii
1715 // phi2 azimuthal angle for axis ii
1716 // theta3 polar angle for axis iii
1717 // phi3 azimuthal angle for axis iii
1719 // it defines the rotation matrix number irot.
1721 Int_t jrotm=fGclink->jrotm;
1727 for(i=1; i<=ns; i++) {
1728 if(fZlq[jrotm-i]==0) {
1734 gsrotm(krot, thex, phix, they, phiy, thez, phiz);
1737 //_____________________________________________________________________________
1738 Int_t TGeant3::GetMedium() const
1741 // Return the number of the current medium
1743 return fGctmed->numed;
1746 //_____________________________________________________________________________
1747 Float_t TGeant3::Edep() const
1750 // Return the energy lost in the current step
1752 return fGctrak->destep;
1755 //_____________________________________________________________________________
1756 Float_t TGeant3::Etot() const
1759 // Return the total energy of the current track
1761 return fGctrak->getot;
1764 //_____________________________________________________________________________
1765 void TGeant3::Rndm(Float_t* r, const Int_t n) const
1768 // Return an array of n random numbers uniformly distributed
1769 // between 0 and 1 not included
1774 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1776 // Functions from GBASE
1778 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1780 //____________________________________________________________________________
1781 void TGeant3::Gfile(const char *filename, const char *option)
1784 // Routine to open a GEANT/RZ data base.
1786 // LUN logical unit number associated to the file
1788 // CHFILE RZ file name
1790 // CHOPT is a character string which may be
1791 // N To create a new file
1792 // U to open an existing file for update
1793 // " " to open an existing file for read only
1794 // Q The initial allocation (default 1000 records)
1795 // is given in IQUEST(10)
1796 // X Open the file in exchange format
1797 // I Read all data structures from file to memory
1798 // O Write all data structures from memory to file
1801 // If options "I" or "O" all data structures are read or
1802 // written from/to file and the file is closed.
1803 // See routine GRMDIR to create subdirectories
1804 // See routines GROUT,GRIN to write,read objects
1806 grfile(21, PASSCHARD(filename), PASSCHARD(option) PASSCHARL(filename)
1810 //____________________________________________________________________________
1811 void TGeant3::Gpcxyz()
1814 // Print track and volume parameters at current point
1819 //_____________________________________________________________________________
1820 void TGeant3::Ggclos()
1823 // Closes off the geometry setting.
1824 // Initializes the search list for the contents of each
1825 // volume following the order they have been positioned, and
1826 // inserting the content '0' when a call to GSNEXT (-1) has
1827 // been required by the user.
1828 // Performs the development of the JVOLUM structure for all
1829 // volumes with variable parameters, by calling GGDVLP.
1830 // Interprets the user calls to GSORD, through GGORD.
1831 // Computes and stores in a bank (next to JVOLUM mother bank)
1832 // the number of levels in the geometrical tree and the
1833 // maximum number of contents per level, by calling GGNLEV.
1834 // Sets status bit for CONCAVE volumes, through GGCAVE.
1835 // Completes the JSET structure with the list of volume names
1836 // which identify uniquely a given physical detector, the
1837 // list of bit numbers to pack the corresponding volume copy
1838 // numbers, and the generic path(s) in the JVOLUM tree,
1839 // through the routine GHCLOS.
1842 // Create internal list of volumes
1843 fVolNames = new char[fGcnum->nvolum+1][5];
1845 for(i=0; i<fGcnum->nvolum; ++i) {
1846 strncpy(fVolNames[i], (char *) &fZiq[fGclink->jvolum+i+1], 4);
1847 fVolNames[i][4]='\0';
1849 strcpy(fVolNames[fGcnum->nvolum],"NULL");
1852 //_____________________________________________________________________________
1853 void TGeant3::Glast()
1856 // Finish a Geant run
1861 //_____________________________________________________________________________
1862 void TGeant3::Gprint(const char *name)
1865 // Routine to print data structures
1866 // CHNAME name of a data structure
1870 gprint(PASSCHARD(vname),0 PASSCHARL(vname));
1873 //_____________________________________________________________________________
1874 void TGeant3::Grun()
1877 // Steering function to process one run
1882 //_____________________________________________________________________________
1883 void TGeant3::Gtrig()
1886 // Steering function to process one event
1889 fDecayer->ForceDecay(fForceDecay);
1893 //_____________________________________________________________________________
1894 void TGeant3::Gtrigc()
1897 // Clear event partition
1902 //_____________________________________________________________________________
1903 void TGeant3::Gtrigi()
1906 // Initialises event partition
1911 //_____________________________________________________________________________
1912 void TGeant3::Gwork(Int_t nwork)
1915 // Allocates workspace in ZEBRA memory
1920 //_____________________________________________________________________________
1921 void TGeant3::Gzinit()
1924 // To initialise GEANT/ZEBRA data structures
1929 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1931 // Functions from GCONS
1933 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1935 //_____________________________________________________________________________
1936 void TGeant3::Gfmate(Int_t imat, char *name, Float_t &a, Float_t &z,
1937 Float_t &dens, Float_t &radl, Float_t &absl,
1938 Float_t* ubuf, Int_t& nbuf)
1941 // Return parameters for material IMAT
1943 gfmate(imat, PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
1947 //_____________________________________________________________________________
1948 void TGeant3::Gfpart(Int_t ipart, char *name, Int_t &itrtyp,
1949 Float_t &amass, Float_t &charge, Float_t &tlife)
1952 // Return parameters for particle of type IPART
1956 Int_t igpart = IdFromPDG(ipart);
1957 gfpart(igpart, PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
1961 //_____________________________________________________________________________
1962 void TGeant3::Gftmed(Int_t numed, char *name, Int_t &nmat, Int_t &isvol,
1963 Int_t &ifield, Float_t &fieldm, Float_t &tmaxfd,
1964 Float_t &stemax, Float_t &deemax, Float_t &epsil,
1965 Float_t &stmin, Float_t *ubuf, Int_t *nbuf)
1968 // Return parameters for tracking medium NUMED
1970 gftmed(numed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1971 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1975 void TGeant3::Gftmat(Int_t imate, Int_t ipart, char *chmeca, Int_t kdim,
1976 Float_t* tkin, Float_t* value, Float_t* pcut,
1980 // Return parameters for tracking medium NUMED
1982 gftmat(imate, ipart, PASSCHARD(chmeca), kdim,
1983 tkin, value, pcut, ixst PASSCHARL(chmeca));
1987 //_____________________________________________________________________________
1988 Float_t TGeant3::Gbrelm(Float_t z, Float_t t, Float_t bcut)
1991 // To calculate energy loss due to soft muon BREMSSTRAHLUNG
1993 return gbrelm(z,t,bcut);
1996 //_____________________________________________________________________________
1997 Float_t TGeant3::Gprelm(Float_t z, Float_t t, Float_t bcut)
2000 // To calculate DE/DX in GeV*barn/atom for direct pair production by muons
2002 return gprelm(z,t,bcut);
2005 //_____________________________________________________________________________
2006 void TGeant3::Gmate()
2009 // Define standard GEANT materials
2014 //_____________________________________________________________________________
2015 void TGeant3::Gpart()
2018 // Define standard GEANT particles plus selected decay modes
2019 // and branching ratios.
2024 //_____________________________________________________________________________
2025 void TGeant3::Gsdk(Int_t ipart, Float_t *bratio, Int_t *mode)
2027 // Defines branching ratios and decay modes for standard
2029 gsdk(ipart,bratio,mode);
2032 //_____________________________________________________________________________
2033 void TGeant3::Gsmate(Int_t imat, const char *name, Float_t a, Float_t z,
2034 Float_t dens, Float_t radl, Float_t absl)
2037 // Defines a Material
2039 // kmat number assigned to the material
2040 // name material name
2041 // a atomic mass in au
2043 // dens density in g/cm3
2044 // absl absorbtion length in cm
2045 // if >=0 it is ignored and the program
2046 // calculates it, if <0. -absl is taken
2047 // radl radiation length in cm
2048 // if >=0 it is ignored and the program
2049 // calculates it, if <0. -radl is taken
2050 // buf pointer to an array of user words
2051 // nbuf number of user words
2055 gsmate(imat,PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
2059 //_____________________________________________________________________________
2060 void TGeant3::Gsmixt(Int_t imat, const char *name, Float_t *a, Float_t *z,
2061 Float_t dens, Int_t nlmat, Float_t *wmat)
2064 // Defines mixture OR COMPOUND IMAT as composed by
2065 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
2067 // If NLMAT.GT.0 then WMAT contains the PROPORTION BY
2068 // WEIGTHS OF EACH BASIC MATERIAL IN THE MIXTURE.
2070 // If NLMAT.LT.0 then WMAT contains the number of atoms
2071 // of a given kind into the molecule of the COMPOUND
2072 // In this case, WMAT in output is changed to relative
2075 gsmixt(imat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
2078 //_____________________________________________________________________________
2079 void TGeant3::Gspart(Int_t ipart, const char *name, Int_t itrtyp,
2080 Float_t amass, Float_t charge, Float_t tlife)
2083 // Store particle parameters
2085 // ipart particle code
2086 // name particle name
2087 // itrtyp transport method (see GEANT manual)
2088 // amass mass in GeV/c2
2089 // charge charge in electron units
2090 // tlife lifetime in seconds
2094 gspart(ipart,PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
2098 //_____________________________________________________________________________
2099 void TGeant3::Gstmed(Int_t numed, const char *name, Int_t nmat, Int_t isvol,
2100 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
2101 Float_t stemax, Float_t deemax, Float_t epsil,
2105 // NTMED Tracking medium number
2106 // NAME Tracking medium name
2107 // NMAT Material number
2108 // ISVOL Sensitive volume flag
2109 // IFIELD Magnetic field
2110 // FIELDM Max. field value (Kilogauss)
2111 // TMAXFD Max. angle due to field (deg/step)
2112 // STEMAX Max. step allowed
2113 // DEEMAX Max. fraction of energy lost in a step
2114 // EPSIL Tracking precision (cm)
2115 // STMIN Min. step due to continuos processes (cm)
2117 // IFIELD = 0 if no magnetic field; IFIELD = -1 if user decision in GUSWIM;
2118 // IFIELD = 1 if tracking performed with GRKUTA; IFIELD = 2 if tracking
2119 // performed with GHELIX; IFIELD = 3 if tracking performed with GHELX3.
2123 gstmed(numed,PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
2124 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
2127 //_____________________________________________________________________________
2128 void TGeant3::Gsckov(Int_t itmed, Int_t npckov, Float_t *ppckov,
2129 Float_t *absco, Float_t *effic, Float_t *rindex)
2132 // Stores the tables for UV photon tracking in medium ITMED
2133 // Please note that it is the user's responsability to
2134 // provide all the coefficients:
2137 // ITMED Tracking medium number
2138 // NPCKOV Number of bins of each table
2139 // PPCKOV Value of photon momentum (in GeV)
2140 // ABSCO Absorbtion coefficients
2141 // dielectric: absorbtion length in cm
2142 // metals : absorbtion fraction (0<=x<=1)
2143 // EFFIC Detection efficiency for UV photons
2144 // RINDEX Refraction index (if=0 metal)
2146 gsckov(itmed,npckov,ppckov,absco,effic,rindex);
2149 //_____________________________________________________________________________
2150 void TGeant3::Gstpar(Int_t itmed, const char *param, Float_t parval)
2153 // To change the value of cut or mechanism "CHPAR"
2154 // to a new value PARVAL for tracking medium ITMED
2155 // The data structure JTMED contains the standard tracking
2156 // parameters (CUTS and flags to control the physics processes) which
2157 // are used by default for all tracking media. It is possible to
2158 // redefine individually with GSTPAR any of these parameters for a
2159 // given tracking medium.
2160 // ITMED tracking medium number
2161 // CHPAR is a character string (variable name)
2162 // PARVAL must be given as a floating point.
2164 gstpar(itmed,PASSCHARD(param), parval PASSCHARL(param));
2167 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2169 // Functions from GCONS
2171 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2173 //_____________________________________________________________________________
2174 void TGeant3::Gfkine(Int_t itra, Float_t *vert, Float_t *pvert, Int_t &ipart,
2177 // Storing/Retrieving Vertex and Track parameters
2178 // ----------------------------------------------
2180 // Stores vertex parameters.
2181 // VERT array of (x,y,z) position of the vertex
2182 // NTBEAM beam track number origin of the vertex
2183 // =0 if none exists
2184 // NTTARG target track number origin of the vertex
2185 // UBUF user array of NUBUF floating point numbers
2187 // NVTX new vertex number (=0 in case of error).
2188 // Prints vertex parameters.
2189 // IVTX for vertex IVTX.
2190 // (For all vertices if IVTX=0)
2191 // Stores long life track parameters.
2192 // PLAB components of momentum
2193 // IPART type of particle (see GSPART)
2194 // NV vertex number origin of track
2195 // UBUF array of NUBUF floating point user parameters
2197 // NT track number (if=0 error).
2198 // Retrieves long life track parameters.
2199 // ITRA track number for which parameters are requested
2200 // VERT vector origin of the track
2201 // PVERT 4 momentum components at the track origin
2202 // IPART particle type (=0 if track ITRA does not exist)
2203 // NVERT vertex number origin of the track
2204 // UBUF user words stored in GSKINE.
2205 // Prints initial track parameters.
2206 // ITRA for track ITRA
2207 // (For all tracks if ITRA=0)
2211 gfkine(itra,vert,pvert,ipart,nvert,ubuf,nbuf);
2214 //_____________________________________________________________________________
2215 void TGeant3::Gfvert(Int_t nvtx, Float_t *v, Int_t &ntbeam, Int_t &nttarg,
2219 // Retrieves the parameter of a vertex bank
2220 // Vertex is generated from tracks NTBEAM NTTARG
2221 // NVTX is the new vertex number
2225 gfvert(nvtx,v,ntbeam,nttarg,tofg,ubuf,nbuf);
2228 //_____________________________________________________________________________
2229 Int_t TGeant3::Gskine(Float_t *plab, Int_t ipart, Int_t nv, Float_t *buf,
2233 // Store kinematics of track NT into data structure
2234 // Track is coming from vertex NV
2237 gskine(plab, ipart, nv, buf, nwbuf, nt);
2241 //_____________________________________________________________________________
2242 Int_t TGeant3::Gsvert(Float_t *v, Int_t ntbeam, Int_t nttarg, Float_t *ubuf,
2246 // Creates a new vertex bank
2247 // Vertex is generated from tracks NTBEAM NTTARG
2248 // NVTX is the new vertex number
2251 gsvert(v, ntbeam, nttarg, ubuf, nwbuf, nwtx);
2255 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2257 // Functions from GPHYS
2259 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2261 //_____________________________________________________________________________
2262 void TGeant3::Gphysi()
2265 // Initialise material constants for all the physics
2266 // mechanisms used by GEANT
2271 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2273 // Functions from GTRAK
2275 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2277 //_____________________________________________________________________________
2278 void TGeant3::Gdebug()
2281 // Debug the current step
2286 //_____________________________________________________________________________
2287 void TGeant3::Gekbin()
2290 // To find bin number in kinetic energy table
2291 // stored in ELOW(NEKBIN)
2296 //_____________________________________________________________________________
2297 void TGeant3::Gfinds()
2300 // Returns the set/volume parameters corresponding to
2301 // the current space point in /GCTRAK/
2302 // and fill common /GCSETS/
2304 // IHSET user set identifier
2305 // IHDET user detector identifier
2306 // ISET set number in JSET
2307 // IDET detector number in JS=LQ(JSET-ISET)
2308 // IDTYPE detector type (1,2)
2309 // NUMBV detector volume numbers (array of length NVNAME)
2310 // NVNAME number of volume levels
2315 //_____________________________________________________________________________
2316 void TGeant3::Gsking(Int_t igk)
2319 // Stores in stack JSTAK either the IGKth track of /GCKING/,
2320 // or the NGKINE tracks when IGK is 0.
2325 //_____________________________________________________________________________
2326 void TGeant3::Gskpho(Int_t igk)
2329 // Stores in stack JSTAK either the IGKth Cherenkov photon of
2330 // /GCKIN2/, or the NPHOT tracks when IGK is 0.
2335 //_____________________________________________________________________________
2336 void TGeant3::Gsstak(Int_t iflag)
2339 // Stores in auxiliary stack JSTAK the particle currently
2340 // described in common /GCKINE/.
2342 // On request, creates also an entry in structure JKINE :
2344 // 0 : No entry in JKINE structure required (user)
2345 // 1 : New entry in JVERTX / JKINE structures required (user)
2346 // <0 : New entry in JKINE structure at vertex -IFLAG (user)
2347 // 2 : Entry in JKINE structure exists already (from GTREVE)
2352 //_____________________________________________________________________________
2353 void TGeant3::Gsxyz()
2356 // Store space point VECT in banks JXYZ
2361 //_____________________________________________________________________________
2362 void TGeant3::Gtrack()
2365 // Controls tracking of current particle
2370 //_____________________________________________________________________________
2371 void TGeant3::Gtreve()
2374 // Controls tracking of all particles belonging to the current event
2379 //_____________________________________________________________________________
2380 void TGeant3::GtreveRoot()
2383 // Controls tracking of all particles belonging to the current event
2388 //_____________________________________________________________________________
2389 void TGeant3::Grndm(Float_t *rvec, const Int_t len) const
2392 // To generate a vector RVECV of LEN random numbers
2393 // Copy of the CERN Library routine RANECU
2397 //_____________________________________________________________________________
2398 void TGeant3::Grndmq(Int_t &is1, Int_t &is2, const Int_t iseq,
2399 const Text_t *chopt)
2402 // To set/retrieve the seed of the random number generator
2404 grndmq(is1,is2,iseq,PASSCHARD(chopt) PASSCHARL(chopt));
2407 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2409 // Functions from GDRAW
2411 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2413 //_____________________________________________________________________________
2414 void TGeant3::Gdxyz(Int_t it)
2417 // Draw the points stored with Gsxyz relative to track it
2422 //_____________________________________________________________________________
2423 void TGeant3::Gdcxyz()
2426 // Draw the position of the current track
2431 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2433 // Functions from GGEOM
2435 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2437 //_____________________________________________________________________________
2438 void TGeant3::Gdtom(Float_t *xd, Float_t *xm, Int_t iflag)
2441 // Computes coordinates XM (Master Reference System
2442 // knowing the coordinates XD (Detector Ref System)
2443 // The local reference system can be initialized by
2444 // - the tracking routines and GDTOM used in GUSTEP
2445 // - a call to GSCMED(NLEVEL,NAMES,NUMBER)
2446 // (inverse routine is GMTOD)
2448 // If IFLAG=1 convert coordinates
2449 // IFLAG=2 convert direction cosinus
2451 gdtom(xd, xm, iflag);
2454 //_____________________________________________________________________________
2455 void TGeant3::Glmoth(const char* iudet, Int_t iunum, Int_t &nlev, Int_t *lvols,
2459 // Loads the top part of the Volume tree in LVOLS (IVO's),
2460 // LINDX (IN indices) for a given volume defined through
2461 // its name IUDET and number IUNUM.
2463 // The routine stores only upto the last level where JVOLUM
2464 // data structure is developed. If there is no development
2465 // above the current level, it returns NLEV zero.
2467 glmoth(PASSCHARD(iudet), iunum, nlev, lvols, lindx, idum PASSCHARL(iudet));
2470 //_____________________________________________________________________________
2471 void TGeant3::Gmedia(Float_t *x, Int_t &numed)
2474 // Finds in which volume/medium the point X is, and updates the
2475 // common /GCVOLU/ and the structure JGPAR accordingly.
2477 // NUMED returns the tracking medium number, or 0 if point is
2478 // outside the experimental setup.
2483 //_____________________________________________________________________________
2484 void TGeant3::Gmtod(Float_t *xm, Float_t *xd, Int_t iflag)
2487 // Computes coordinates XD (in DRS)
2488 // from known coordinates XM in MRS
2489 // The local reference system can be initialized by
2490 // - the tracking routines and GMTOD used in GUSTEP
2491 // - a call to GMEDIA(XM,NUMED)
2492 // - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER)
2493 // (inverse routine is GDTOM)
2495 // If IFLAG=1 convert coordinates
2496 // IFLAG=2 convert direction cosinus
2498 gmtod(xm, xd, iflag);
2501 //_____________________________________________________________________________
2502 void TGeant3::Gsdvn(const char *name, const char *mother, Int_t ndiv,
2506 // Create a new volume by dividing an existing one
2509 // MOTHER Mother volume name
2510 // NDIV Number of divisions
2513 // X,Y,Z of CAXIS will be translated to 1,2,3 for IAXIS.
2514 // It divides a previously defined volume.
2519 Vname(mother,vmother);
2520 gsdvn(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis PASSCHARL(vname)
2521 PASSCHARL(vmother));
2524 //_____________________________________________________________________________
2525 void TGeant3::Gsdvn2(const char *name, const char *mother, Int_t ndiv,
2526 Int_t iaxis, Float_t c0i, Int_t numed)
2529 // Create a new volume by dividing an existing one
2531 // Divides mother into ndiv divisions called name
2532 // along axis iaxis starting at coordinate value c0.
2533 // the new volume created will be medium number numed.
2538 Vname(mother,vmother);
2539 gsdvn2(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis, c0i, numed
2540 PASSCHARL(vname) PASSCHARL(vmother));
2543 //_____________________________________________________________________________
2544 void TGeant3::Gsdvs(const char *name, const char *mother, Float_t step,
2545 Int_t iaxis, Int_t numed)
2548 // Create a new volume by dividing an existing one
2553 Vname(mother,vmother);
2554 gsdvs(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed
2555 PASSCHARL(vname) PASSCHARL(vmother));
2558 //_____________________________________________________________________________
2559 void TGeant3::Gsdvs2(const char *name, const char *mother, Float_t step,
2560 Int_t iaxis, Float_t c0, Int_t numed)
2563 // Create a new volume by dividing an existing one
2568 Vname(mother,vmother);
2569 gsdvs2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0, numed
2570 PASSCHARL(vname) PASSCHARL(vmother));
2573 //_____________________________________________________________________________
2574 void TGeant3::Gsdvt(const char *name, const char *mother, Float_t step,
2575 Int_t iaxis, Int_t numed, Int_t ndvmx)
2578 // Create a new volume by dividing an existing one
2580 // Divides MOTHER into divisions called NAME along
2581 // axis IAXIS in steps of STEP. If not exactly divisible
2582 // will make as many as possible and will centre them
2583 // with respect to the mother. Divisions will have medium
2584 // number NUMED. If NUMED is 0, NUMED of MOTHER is taken.
2585 // NDVMX is the expected maximum number of divisions
2586 // (If 0, no protection tests are performed)
2591 Vname(mother,vmother);
2592 gsdvt(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed, ndvmx
2593 PASSCHARL(vname) PASSCHARL(vmother));
2596 //_____________________________________________________________________________
2597 void TGeant3::Gsdvt2(const char *name, const char *mother, Float_t step,
2598 Int_t iaxis, Float_t c0, Int_t numed, Int_t ndvmx)
2601 // Create a new volume by dividing an existing one
2603 // Divides MOTHER into divisions called NAME along
2604 // axis IAXIS starting at coordinate value C0 with step
2606 // The new volume created will have medium number NUMED.
2607 // If NUMED is 0, NUMED of mother is taken.
2608 // NDVMX is the expected maximum number of divisions
2609 // (If 0, no protection tests are performed)
2614 Vname(mother,vmother);
2615 gsdvt2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0,
2616 numed, ndvmx PASSCHARL(vname) PASSCHARL(vmother));
2619 //_____________________________________________________________________________
2620 void TGeant3::Gsord(const char *name, Int_t iax)
2623 // Flags volume CHNAME whose contents will have to be ordered
2624 // along axis IAX, by setting the search flag to -IAX
2628 // IAX = 4 Rxy (static ordering only -> GTMEDI)
2629 // IAX = 14 Rxy (also dynamic ordering -> GTNEXT)
2630 // IAX = 5 Rxyz (static ordering only -> GTMEDI)
2631 // IAX = 15 Rxyz (also dynamic ordering -> GTNEXT)
2632 // IAX = 6 PHI (PHI=0 => X axis)
2633 // IAX = 7 THETA (THETA=0 => Z axis)
2637 gsord(PASSCHARD(vname), iax PASSCHARL(vname));
2640 //_____________________________________________________________________________
2641 void TGeant3::Gspos(const char *name, Int_t nr, const char *mother, Float_t x,
2642 Float_t y, Float_t z, Int_t irot, const char *konly)
2645 // Position a volume into an existing one
2648 // NUMBER Copy number of the volume
2649 // MOTHER Mother volume name
2650 // X X coord. of the volume in mother ref. sys.
2651 // Y Y coord. of the volume in mother ref. sys.
2652 // Z Z coord. of the volume in mother ref. sys.
2653 // IROT Rotation matrix number w.r.t. mother ref. sys.
2654 // ONLY ONLY/MANY flag
2656 // It positions a previously defined volume in the mother.
2661 Vname(mother,vmother);
2662 gspos(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2663 PASSCHARD(konly) PASSCHARL(vname) PASSCHARL(vmother)
2667 //_____________________________________________________________________________
2668 void TGeant3::Gsposp(const char *name, Int_t nr, const char *mother,
2669 Float_t x, Float_t y, Float_t z, Int_t irot,
2670 const char *konly, Float_t *upar, Int_t np )
2673 // Place a copy of generic volume NAME with user number
2674 // NR inside MOTHER, with its parameters UPAR(1..NP)
2679 Vname(mother,vmother);
2680 gsposp(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2681 PASSCHARD(konly), upar, np PASSCHARL(vname) PASSCHARL(vmother)
2685 //_____________________________________________________________________________
2686 void TGeant3::Gsrotm(Int_t nmat, Float_t theta1, Float_t phi1, Float_t theta2,
2687 Float_t phi2, Float_t theta3, Float_t phi3)
2690 // nmat Rotation matrix number
2691 // THETA1 Polar angle for axis I
2692 // PHI1 Azimuthal angle for axis I
2693 // THETA2 Polar angle for axis II
2694 // PHI2 Azimuthal angle for axis II
2695 // THETA3 Polar angle for axis III
2696 // PHI3 Azimuthal angle for axis III
2698 // It defines the rotation matrix number IROT.
2700 gsrotm(nmat, theta1, phi1, theta2, phi2, theta3, phi3);
2703 //_____________________________________________________________________________
2704 void TGeant3::Gprotm(Int_t nmat)
2707 // To print rotation matrices structure JROTM
2708 // nmat Rotation matrix number
2713 //_____________________________________________________________________________
2714 Int_t TGeant3::Gsvolu(const char *name, const char *shape, Int_t nmed,
2715 Float_t *upar, Int_t npar)
2719 // SHAPE Volume type
2720 // NUMED Tracking medium number
2721 // NPAR Number of shape parameters
2722 // UPAR Vector containing shape parameters
2724 // It creates a new volume in the JVOLUM data structure.
2730 Vname(shape,vshape);
2731 gsvolu(PASSCHARD(vname), PASSCHARD(vshape), nmed, upar, npar, ivolu
2732 PASSCHARL(vname) PASSCHARL(vshape));
2736 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2738 // T H E D R A W I N G P A C K A G E
2739 // ======================================
2740 // Drawing functions. These functions allow the visualization in several ways
2741 // of the volumes defined in the geometrical data structure. It is possible
2742 // to draw the logical tree of volumes belonging to the detector (DTREE),
2743 // to show their geometrical specification (DSPEC,DFSPC), to draw them
2744 // and their cut views (DRAW, DCUT). Moreover, it is possible to execute
2745 // these commands when the hidden line removal option is activated; in
2746 // this case, the volumes can be also either translated in the space
2747 // (SHIFT), or clipped by boolean operation (CVOL). In addition, it is
2748 // possible to fill the surfaces of the volumes
2749 // with solid colours when the shading option (SHAD) is activated.
2750 // Several tools (ZOOM, LENS) have been developed to zoom detailed parts
2751 // of the detectors or to scan physical events as well.
2752 // Finally, the command MOVE will allow the rotation, translation and zooming
2753 // on real time parts of the detectors or tracks and hits of a simulated event.
2754 // Ray-tracing commands. In case the command (DOPT RAYT ON) is executed,
2755 // the drawing is performed by the Geant ray-tracing;
2756 // automatically, the color is assigned according to the tracking medium of each
2757 // volume and the volumes with a density lower/equal than the air are considered
2758 // transparent; if the option (USER) is set (ON) (again via the command (DOPT)),
2759 // the user can set color and visibility for the desired volumes via the command
2760 // (SATT), as usual, relatively to the attributes (COLO) and (SEEN).
2761 // The resolution can be set via the command (SATT * FILL VALUE), where (VALUE)
2762 // is the ratio between the number of pixels drawn and 20 (user coordinates).
2763 // Parallel view and perspective view are possible (DOPT PROJ PARA/PERS); in the
2764 // first case, we assume that the first mother volume of the tree is a box with
2765 // dimensions 10000 X 10000 X 10000 cm and the view point (infinetely far) is
2766 // 5000 cm far from the origin along the Z axis of the user coordinates; in the
2767 // second case, the distance between the observer and the origin of the world
2768 // reference system is set in cm by the command (PERSP NAME VALUE); grand-angle
2769 // or telescopic effects can be achieved changing the scale factors in the command
2770 // (DRAW). When the final picture does not occupy the full window,
2771 // mapping the space before tracing can speed up the drawing, but can also
2772 // produce less precise results; values from 1 to 4 are allowed in the command
2773 // (DOPT MAPP VALUE), the mapping being more precise for increasing (VALUE); for
2774 // (VALUE = 0) no mapping is performed (therefore max precision and lowest speed).
2775 // The command (VALCUT) allows the cutting of the detector by three planes
2776 // ortogonal to the x,y,z axis. The attribute (LSTY) can be set by the command
2777 // SATT for any desired volume and can assume values from 0 to 7; it determines
2778 // the different light processing to be performed for different materials:
2779 // 0 = dark-matt, 1 = bright-matt, 2 = plastic, 3 = ceramic, 4 = rough-metals,
2780 // 5 = shiny-metals, 6 = glass, 7 = mirror. The detector is assumed to be in the
2781 // dark, the ambient light luminosity is 0.2 for each basic hue (the saturation
2782 // is 0.9) and the observer is assumed to have a light source (therefore he will
2783 // produce parallel light in the case of parallel view and point-like-source
2784 // light in the case of perspective view).
2786 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2788 //_____________________________________________________________________________
2789 void TGeant3::Gsatt(const char *name, const char *att, Int_t val)
2793 // IOPT Name of the attribute to be set
2794 // IVAL Value to which the attribute is to be set
2796 // name= "*" stands for all the volumes.
2797 // iopt can be chosen among the following :
2799 // WORK 0=volume name is inactive for the tracking
2800 // 1=volume name is active for the tracking (default)
2802 // SEEN 0=volume name is invisible
2803 // 1=volume name is visible (default)
2804 // -1=volume invisible with all its descendants in the tree
2805 // -2=volume visible but not its descendants in the tree
2807 // LSTY line style 1,2,3,... (default=1)
2808 // LSTY=7 will produce a very precise approximation for
2809 // revolution bodies.
2811 // LWID line width -7,...,1,2,3,..7 (default=1)
2812 // LWID<0 will act as abs(LWID) was set for the volume
2813 // and for all the levels below it. When SHAD is 'ON', LWID
2814 // represent the linewidth of the scan lines filling the surfaces
2815 // (whereas the FILL value represent their number). Therefore
2816 // tuning this parameter will help to obtain the desired
2817 // quality/performance ratio.
2819 // COLO colour code -166,...,1,2,..166 (default=1)
2821 // n=2=red; n=17+m, m=0,25, increasing luminosity according to 'm';
2822 // n=3=green; n=67+m, m=0,25, increasing luminosity according to 'm';
2823 // n=4=blue; n=117+m, m=0,25, increasing luminosity according to 'm';
2824 // n=5=yellow; n=42+m, m=0,25, increasing luminosity according to 'm';
2825 // n=6=violet; n=142+m, m=0,25, increasing luminosity according to 'm';
2826 // n=7=lightblue; n=92+m, m=0,25, increasing luminosity according to 'm';
2827 // colour=n*10+m, m=1,2,...9, will produce the same colour
2828 // as 'n', but with increasing luminosity according to 'm';
2829 // COLO<0 will act as if abs(COLO) was set for the volume
2830 // and for all the levels below it.
2831 // When for a volume the attribute FILL is > 1 (and the
2832 // option SHAD is on), the ABS of its colour code must be < 8
2833 // because an automatic shading of its faces will be
2836 // FILL (1992) fill area -7,...,0,1,...7 (default=0)
2837 // when option SHAD is "on" the FILL attribute of any
2838 // volume can be set different from 0 (normal drawing);
2839 // if it is set to 1, the faces of such volume will be filled
2840 // with solid colours; if ABS(FILL) is > 1, then a light
2841 // source is placed along the observer line, and the faces of
2842 // such volumes will be painted by colours whose luminosity
2843 // will depend on the amount of light reflected;
2844 // if ABS(FILL) = 1, then it is possible to use all the 166
2845 // colours of the colour table, becouse the automatic shading
2846 // is not performed;
2847 // for increasing values of FILL the drawing will be performed
2848 // with higher and higher resolution improving the quality (the
2849 // number of scan lines used to fill the faces increases with FILL);
2850 // it is possible to set different values of FILL
2851 // for different volumes, in order to optimize at the same time
2852 // the performance and the quality of the picture;
2853 // FILL<0 will act as if abs(FILL) was set for the volume
2854 // and for all the levels below it.
2855 // This kind of drawing can be saved in 'picture files'
2856 // or in view banks.
2857 // 0=drawing without fill area
2858 // 1=faces filled with solid colours and resolution = 6
2859 // 2=lowest resolution (very fast)
2860 // 3=default resolution
2861 // 4=.................
2862 // 5=.................
2863 // 6=.................
2865 // Finally, if a coloured background is desired, the FILL
2866 // attribute for the first volume of the tree must be set
2867 // equal to -abs(colo), colo being >0 and <166.
2869 // SET set number associated to volume name
2870 // DET detector number associated to volume name
2871 // DTYP detector type (1,2)
2878 gsatt(PASSCHARD(vname), PASSCHARD(vatt), val PASSCHARL(vname)
2882 //_____________________________________________________________________________
2883 void TGeant3::Gfpara(const char *name, Int_t number, Int_t intext, Int_t& npar,
2884 Int_t& natt, Float_t* par, Float_t* att)
2887 // Find the parameters of a volume
2889 gfpara(PASSCHARD(name), number, intext, npar, natt, par, att
2893 //_____________________________________________________________________________
2894 void TGeant3::Gckpar(Int_t ish, Int_t npar, Float_t* par)
2897 // Check the parameters of a shape
2899 gckpar(ish,npar,par);
2902 //_____________________________________________________________________________
2903 void TGeant3::Gckmat(Int_t itmed, char* natmed)
2906 // Check the parameters of a tracking medium
2908 gckmat(itmed, PASSCHARD(natmed) PASSCHARL(natmed));
2911 //_____________________________________________________________________________
2912 void TGeant3::Gdelete(Int_t iview)
2915 // IVIEW View number
2917 // It deletes a view bank from memory.
2922 //_____________________________________________________________________________
2923 void TGeant3::Gdopen(Int_t iview)
2926 // IVIEW View number
2928 // When a drawing is very complex and requires a long time to be
2929 // executed, it can be useful to store it in a view bank: after a
2930 // call to DOPEN and the execution of the drawing (nothing will
2931 // appear on the screen), and after a necessary call to DCLOSE,
2932 // the contents of the bank can be displayed in a very fast way
2933 // through a call to DSHOW; therefore, the detector can be easily
2934 // zoomed many times in different ways. Please note that the pictures
2935 // with solid colours can now be stored in a view bank or in 'PICTURE FILES'
2942 //_____________________________________________________________________________
2943 void TGeant3::Gdclose()
2946 // It closes the currently open view bank; it must be called after the
2947 // end of the drawing to be stored.
2952 //_____________________________________________________________________________
2953 void TGeant3::Gdshow(Int_t iview)
2956 // IVIEW View number
2958 // It shows on the screen the contents of a view bank. It
2959 // can be called after a view bank has been closed.
2964 //_____________________________________________________________________________
2965 void TGeant3::Gdopt(const char *name,const char *value)
2969 // VALUE Option value
2971 // To set/modify the drawing options.
2974 // THRZ ON Draw tracks in R vs Z
2975 // OFF (D) Draw tracks in X,Y,Z
2978 // PROJ PARA (D) Parallel projection
2980 // TRAK LINE (D) Trajectory drawn with lines
2981 // POIN " " with markers
2982 // HIDE ON Hidden line removal using the CG package
2983 // OFF (D) No hidden line removal
2984 // SHAD ON Fill area and shading of surfaces.
2985 // OFF (D) Normal hidden line removal.
2986 // RAYT ON Ray-tracing on.
2987 // OFF (D) Ray-tracing off.
2988 // EDGE OFF Does not draw contours when shad is on.
2989 // ON (D) Normal shading.
2990 // MAPP 1,2,3,4 Mapping before ray-tracing.
2991 // 0 (D) No mapping.
2992 // USER ON User graphics options in the raytracing.
2993 // OFF (D) Automatic graphics options.
2999 Vname(value,vvalue);
3000 gdopt(PASSCHARD(vname), PASSCHARD(vvalue) PASSCHARL(vname)
3004 //_____________________________________________________________________________
3005 void TGeant3::Gdraw(const char *name,Float_t theta, Float_t phi, Float_t psi,
3006 Float_t u0,Float_t v0,Float_t ul,Float_t vl)
3011 // THETA Viewing angle theta (for 3D projection)
3012 // PHI Viewing angle phi (for 3D projection)
3013 // PSI Viewing angle psi (for 2D rotation)
3014 // U0 U-coord. (horizontal) of volume origin
3015 // V0 V-coord. (vertical) of volume origin
3016 // SU Scale factor for U-coord.
3017 // SV Scale factor for V-coord.
3019 // This function will draw the volumes,
3020 // selected with their graphical attributes, set by the Gsatt
3021 // facility. The drawing may be performed with hidden line removal
3022 // and with shading effects according to the value of the options HIDE
3023 // and SHAD; if the option SHAD is ON, the contour's edges can be
3024 // drawn or not. If the option HIDE is ON, the detector can be
3025 // exploded (BOMB), clipped with different shapes (CVOL), and some
3026 // of its parts can be shifted from their original
3027 // position (SHIFT). When HIDE is ON, if
3028 // the drawing requires more than the available memory, the program
3029 // will evaluate and display the number of missing words
3030 // (so that the user can increase the
3031 // size of its ZEBRA store). Finally, at the end of each drawing (with HIDE on),
3032 // the program will print messages about the memory used and
3033 // statistics on the volumes' visibility.
3034 // The following commands will produce the drawing of a green
3035 // volume, specified by NAME, without using the hidden line removal
3036 // technique, using the hidden line removal technique,
3037 // with different linewidth and colour (red), with
3038 // solid colour, with shading of surfaces, and without edges.
3039 // Finally, some examples are given for the ray-tracing. (A possible
3040 // string for the NAME of the volume can be found using the command DTREE).
3046 if (fGcvdma->raytra != 1) {
3047 gdraw(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
3049 gdrayt(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
3053 //_____________________________________________________________________________
3054 void TGeant3::Gdrawc(const char *name,Int_t axis, Float_t cut,Float_t u0,
3055 Float_t v0,Float_t ul,Float_t vl)
3060 // CUTVAL Cut plane distance from the origin along the axis
3062 // U0 U-coord. (horizontal) of volume origin
3063 // V0 V-coord. (vertical) of volume origin
3064 // SU Scale factor for U-coord.
3065 // SV Scale factor for V-coord.
3067 // The cut plane is normal to caxis (X,Y,Z), corresponding to iaxis (1,2,3),
3068 // and placed at the distance cutval from the origin.
3069 // The resulting picture is seen from the the same axis.
3070 // When HIDE Mode is ON, it is possible to get the same effect with
3071 // the CVOL/BOX function.
3077 gdrawc(PASSCHARD(vname), axis,cut,u0,v0,ul,vl PASSCHARL(vname));
3080 //_____________________________________________________________________________
3081 void TGeant3::Gdrawx(const char *name,Float_t cutthe, Float_t cutphi,
3082 Float_t cutval, Float_t theta, Float_t phi, Float_t u0,
3083 Float_t v0,Float_t ul,Float_t vl)
3087 // CUTTHE Theta angle of the line normal to cut plane
3088 // CUTPHI Phi angle of the line normal to cut plane
3089 // CUTVAL Cut plane distance from the origin along the axis
3091 // THETA Viewing angle theta (for 3D projection)
3092 // PHI Viewing angle phi (for 3D projection)
3093 // U0 U-coord. (horizontal) of volume origin
3094 // V0 V-coord. (vertical) of volume origin
3095 // SU Scale factor for U-coord.
3096 // SV Scale factor for V-coord.
3098 // The cut plane is normal to the line given by the cut angles
3099 // cutthe and cutphi and placed at the distance cutval from the origin.
3100 // The resulting picture is seen from the viewing angles theta,phi.
3106 gdrawx(PASSCHARD(vname), cutthe,cutphi,cutval,theta,phi,u0,v0,ul,vl
3110 //_____________________________________________________________________________
3111 void TGeant3::Gdhead(Int_t isel, const char *name, Float_t chrsiz)
3116 // ISEL Option flag D=111110
3118 // CHRSIZ Character size (cm) of title NAME D=0.6
3121 // 0 to have only the header lines
3122 // xxxxx1 to add the text name centered on top of header
3123 // xxxx1x to add global detector name (first volume) on left
3124 // xxx1xx to add date on right
3125 // xx1xxx to select thick characters for text on top of header
3126 // x1xxxx to add the text 'EVENT NR x' on top of header
3127 // 1xxxxx to add the text 'RUN NR x' on top of header
3128 // NOTE that ISEL=x1xxx1 or ISEL=1xxxx1 are illegal choices,
3129 // i.e. they generate overwritten text.
3131 gdhead(isel,PASSCHARD(name),chrsiz PASSCHARL(name));
3134 //_____________________________________________________________________________
3135 void TGeant3::Gdman(Float_t u, Float_t v, const char *type)
3138 // Draw a 2D-man at position (U0,V0)
3140 // U U-coord. (horizontal) of the centre of man' R
3141 // V V-coord. (vertical) of the centre of man' R
3142 // TYPE D='MAN' possible values: 'MAN,WM1,WM2,WM3'
3144 // CALL GDMAN(u,v),CALL GDWMN1(u,v),CALL GDWMN2(u,v),CALL GDWMN2(u,v)
3145 // It superimposes the picure of a man or of a woman, chosen among
3146 // three different ones, with the same scale factors as the detector
3147 // in the current drawing.
3150 if (opt.Contains("WM1")) {
3152 } else if (opt.Contains("WM3")) {
3154 } else if (opt.Contains("WM2")) {
3161 //_____________________________________________________________________________
3162 void TGeant3::Gdspec(const char *name)
3167 // Shows 3 views of the volume (two cut-views and a 3D view), together with
3168 // its geometrical specifications. The 3D drawing will
3169 // be performed according the current values of the options HIDE and
3170 // SHAD and according the current SetClipBox clipping parameters for that
3177 gdspec(PASSCHARD(vname) PASSCHARL(vname));
3180 //_____________________________________________________________________________
3181 void TGeant3::DrawOneSpec(const char *name)
3184 // Function called when one double-clicks on a volume name
3185 // in a TPavelabel drawn by Gdtree.
3187 THIGZ *higzSave = gHigz;
3188 higzSave->SetName("higzSave");
3189 THIGZ *higzSpec = (THIGZ*)gROOT->FindObject("higzSpec");
3190 //printf("DrawOneSpec, gHigz=%x, higzSpec=%x\n",gHigz,higzSpec);
3191 if (higzSpec) gHigz = higzSpec;
3192 else higzSpec = new THIGZ(kDefSize);
3193 higzSpec->SetName("higzSpec");
3198 gdspec(PASSCHARD(vname) PASSCHARL(vname));
3201 higzSave->SetName("higz");
3205 //_____________________________________________________________________________
3206 void TGeant3::Gdtree(const char *name,Int_t levmax, Int_t isel)
3210 // LEVMAX Depth level
3213 // This function draws the logical tree,
3214 // Each volume in the tree is represented by a TPaveTree object.
3215 // Double-clicking on a TPaveTree draws the specs of the corresponding volume.
3216 // Use TPaveTree pop-up menu to select:
3219 // - drawing tree of parent
3225 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
3226 gHigz->SetPname("");
3229 //_____________________________________________________________________________
3230 void TGeant3::GdtreeParent(const char *name,Int_t levmax, Int_t isel)
3234 // LEVMAX Depth level
3237 // This function draws the logical tree of the parent of name.
3241 // Scan list of volumes in JVOLUM
3243 Int_t gname, i, jvo, in, nin, jin, num;
3244 strncpy((char *) &gname, name, 4);
3245 for(i=1; i<=fGcnum->nvolum; i++) {
3246 jvo = fZlq[fGclink->jvolum-i];
3247 nin = Int_t(fZq[jvo+3]);
3248 if (nin == -1) nin = 1;
3249 for (in=1;in<=nin;in++) {
3251 num = Int_t(fZq[jin+2]);
3252 if(gname == fZiq[fGclink->jvolum+num]) {
3253 strncpy(vname,(char*)&fZiq[fGclink->jvolum+i],4);
3255 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
3256 gHigz->SetPname("");
3263 //_____________________________________________________________________________
3264 void TGeant3::SetABAN(Int_t par)
3267 // par = 1 particles will be stopped according to their residual
3268 // range if they are not in a sensitive material and are
3269 // far enough from the boundary
3270 // 0 particles are transported normally
3272 fGcphys->dphys1 = par;
3276 //_____________________________________________________________________________
3277 void TGeant3::SetANNI(Int_t par)
3280 // To control positron annihilation.
3281 // par =0 no annihilation
3282 // =1 annihilation. Decays processed.
3283 // =2 annihilation. No decay products stored.
3285 fGcphys->ianni = par;
3289 //_____________________________________________________________________________
3290 void TGeant3::SetAUTO(Int_t par)
3293 // To control automatic calculation of tracking medium parameters:
3294 // par =0 no automatic calculation;
3295 // =1 automati calculation.
3297 fGctrak->igauto = par;
3301 //_____________________________________________________________________________
3302 void TGeant3::SetBOMB(Float_t boom)
3305 // BOOM : Exploding factor for volumes position
3307 // To 'explode' the detector. If BOOM is positive (values smaller
3308 // than 1. are suggested, but any value is possible)
3309 // all the volumes are shifted by a distance
3310 // proportional to BOOM along the direction between their centre
3311 // and the origin of the MARS; the volumes which are symmetric
3312 // with respect to this origin are simply not shown.
3313 // BOOM equal to 0 resets the normal mode.
3314 // A negative (greater than -1.) value of
3315 // BOOM will cause an 'implosion'; for even lower values of BOOM
3316 // the volumes' positions will be reflected respect to the origin.
3317 // This command can be useful to improve the 3D effect for very
3318 // complex detectors. The following commands will make explode the
3325 //_____________________________________________________________________________
3326 void TGeant3::SetBREM(Int_t par)
3329 // To control bremstrahlung.
3330 // par =0 no bremstrahlung
3331 // =1 bremstrahlung. Photon processed.
3332 // =2 bremstrahlung. No photon stored.
3334 fGcphys->ibrem = par;
3338 //_____________________________________________________________________________
3339 void TGeant3::SetCKOV(Int_t par)
3342 // To control Cerenkov production
3343 // par =0 no Cerenkov;
3345 // =2 Cerenkov with primary stopped at each step.
3347 fGctlit->itckov = par;
3351 //_____________________________________________________________________________
3352 void TGeant3::SetClipBox(const char *name,Float_t xmin,Float_t xmax,
3353 Float_t ymin,Float_t ymax,Float_t zmin,Float_t zmax)
3356 // The hidden line removal technique is necessary to visualize properly
3357 // very complex detectors. At the same time, it can be useful to visualize
3358 // the inner elements of a detector in detail. This function allows
3359 // subtractions (via boolean operation) of BOX shape from any part of
3360 // the detector, therefore showing its inner contents.
3361 // If "*" is given as the name of the
3362 // volume to be clipped, all volumes are clipped by the given box.
3363 // A volume can be clipped at most twice.
3364 // if a volume is explicitely clipped twice,
3365 // the "*" will not act on it anymore. Giving "." as the name
3366 // of the volume to be clipped will reset the clipping.
3368 // NAME Name of volume to be clipped
3370 // XMIN Lower limit of the Shape X coordinate
3371 // XMAX Upper limit of the Shape X coordinate
3372 // YMIN Lower limit of the Shape Y coordinate
3373 // YMAX Upper limit of the Shape Y coordinate
3374 // ZMIN Lower limit of the Shape Z coordinate
3375 // ZMAX Upper limit of the Shape Z coordinate
3377 // This function performs a boolean subtraction between the volume
3378 // NAME and a box placed in the MARS according the values of the given
3384 setclip(PASSCHARD(vname),xmin,xmax,ymin,ymax,zmin,zmax PASSCHARL(vname));
3387 //_____________________________________________________________________________
3388 void TGeant3::SetCOMP(Int_t par)
3391 // To control Compton scattering
3392 // par =0 no Compton
3393 // =1 Compton. Electron processed.
3394 // =2 Compton. No electron stored.
3397 fGcphys->icomp = par;
3400 //_____________________________________________________________________________
3401 void TGeant3::SetCUTS(Float_t cutgam,Float_t cutele,Float_t cutneu,
3402 Float_t cuthad,Float_t cutmuo ,Float_t bcute ,
3403 Float_t bcutm ,Float_t dcute ,Float_t dcutm ,
3404 Float_t ppcutm, Float_t tofmax)
3407 // CUTGAM Cut for gammas D=0.001
3408 // CUTELE Cut for electrons D=0.001
3409 // CUTHAD Cut for charged hadrons D=0.01
3410 // CUTNEU Cut for neutral hadrons D=0.01
3411 // CUTMUO Cut for muons D=0.01
3412 // BCUTE Cut for electron brems. D=-1.
3413 // BCUTM Cut for muon brems. D=-1.
3414 // DCUTE Cut for electron delta-rays D=-1.
3415 // DCUTM Cut for muon delta-rays D=-1.
3416 // PPCUTM Cut for e+e- pairs by muons D=0.01
3417 // TOFMAX Time of flight cut D=1.E+10
3419 // If the default values (-1.) for BCUTE ,BCUTM ,DCUTE ,DCUTM
3420 // are not modified, they will be set to CUTGAM,CUTGAM,CUTELE,CUTELE
3422 // If one of the parameters from CUTGAM to PPCUTM included
3423 // is modified, cross-sections and energy loss tables must be
3424 // recomputed via the function Gphysi.
3426 fGccuts->cutgam = cutgam;
3427 fGccuts->cutele = cutele;
3428 fGccuts->cutneu = cutneu;
3429 fGccuts->cuthad = cuthad;
3430 fGccuts->cutmuo = cutmuo;
3431 fGccuts->bcute = bcute;
3432 fGccuts->bcutm = bcutm;
3433 fGccuts->dcute = dcute;
3434 fGccuts->dcutm = dcutm;
3435 fGccuts->ppcutm = ppcutm;
3436 fGccuts->tofmax = tofmax;
3439 //_____________________________________________________________________________
3440 void TGeant3::SetDCAY(Int_t par)
3443 // To control Decay mechanism.
3444 // par =0 no decays.
3445 // =1 Decays. secondaries processed.
3446 // =2 Decays. No secondaries stored.
3448 fGcphys->idcay = par;
3452 //_____________________________________________________________________________
3453 void TGeant3::SetDEBU(Int_t emin, Int_t emax, Int_t emod)
3456 // Set the debug flag and frequency
3457 // Selected debug output will be printed from
3458 // event emin to even emax each emod event
3460 fGcflag->idemin = emin;
3461 fGcflag->idemax = emax;
3462 fGcflag->itest = emod;
3466 //_____________________________________________________________________________
3467 void TGeant3::SetDRAY(Int_t par)
3470 // To control delta rays mechanism.
3471 // par =0 no delta rays.
3472 // =1 Delta rays. secondaries processed.
3473 // =2 Delta rays. No secondaries stored.
3475 fGcphys->idray = par;
3478 //_____________________________________________________________________________
3479 void TGeant3::SetERAN(Float_t ekmin, Float_t ekmax, Int_t nekbin)
3482 // To control cross section tabulations
3483 // ekmin = minimum kinetic energy in GeV
3484 // ekmax = maximum kinetic energy in GeV
3485 // nekbin = number of logatithmic bins (<200)
3487 fGcmulo->ekmin = ekmin;
3488 fGcmulo->ekmax = ekmax;
3489 fGcmulo->nekbin = nekbin;
3492 //_____________________________________________________________________________
3493 void TGeant3::SetHADR(Int_t par)
3496 // To control hadronic interactions.
3497 // par =0 no hadronic interactions.
3498 // =1 Hadronic interactions. secondaries processed.
3499 // =2 Hadronic interactions. No secondaries stored.
3501 fGcphys->ihadr = par;
3504 //_____________________________________________________________________________
3505 void TGeant3::SetKINE(Int_t kine, Float_t xk1, Float_t xk2, Float_t xk3,
3506 Float_t xk4, Float_t xk5, Float_t xk6, Float_t xk7,
3507 Float_t xk8, Float_t xk9, Float_t xk10)
3510 // Set the variables in /GCFLAG/ IKINE, PKINE(10)
3511 // Their meaning is user defined
3513 fGckine->ikine = kine;
3514 fGckine->pkine[0] = xk1;
3515 fGckine->pkine[1] = xk2;
3516 fGckine->pkine[2] = xk3;
3517 fGckine->pkine[3] = xk4;
3518 fGckine->pkine[4] = xk5;
3519 fGckine->pkine[5] = xk6;
3520 fGckine->pkine[6] = xk7;
3521 fGckine->pkine[7] = xk8;
3522 fGckine->pkine[8] = xk9;
3523 fGckine->pkine[9] = xk10;
3526 //_____________________________________________________________________________
3527 void TGeant3::SetLOSS(Int_t par)
3530 // To control energy loss.
3531 // par =0 no energy loss;
3532 // =1 restricted energy loss fluctuations;
3533 // =2 complete energy loss fluctuations;
3535 // =4 no energy loss fluctuations.
3536 // If the value ILOSS is changed, then cross-sections and energy loss
3537 // tables must be recomputed via the command 'PHYSI'.
3539 fGcphys->iloss = par;
3543 //_____________________________________________________________________________
3544 void TGeant3::SetMULS(Int_t par)
3547 // To control multiple scattering.
3548 // par =0 no multiple scattering.
3549 // =1 Moliere or Coulomb scattering.
3550 // =2 Moliere or Coulomb scattering.
3551 // =3 Gaussian scattering.
3553 fGcphys->imuls = par;
3557 //_____________________________________________________________________________
3558 void TGeant3::SetMUNU(Int_t par)
3561 // To control muon nuclear interactions.
3562 // par =0 no muon-nuclear interactions.
3563 // =1 Nuclear interactions. Secondaries processed.
3564 // =2 Nuclear interactions. Secondaries not processed.
3566 fGcphys->imunu = par;
3569 //_____________________________________________________________________________
3570 void TGeant3::SetOPTI(Int_t par)
3573 // This flag controls the tracking optimisation performed via the
3575 // 1 no optimisation at all; GSORD calls disabled;
3576 // 0 no optimisation; only user calls to GSORD kept;
3577 // 1 all non-GSORDered volumes are ordered along the best axis;
3578 // 2 all volumes are ordered along the best axis.
3580 fGcopti->ioptim = par;
3583 //_____________________________________________________________________________
3584 void TGeant3::SetPAIR(Int_t par)
3587 // To control pair production mechanism.
3588 // par =0 no pair production.
3589 // =1 Pair production. secondaries processed.
3590 // =2 Pair production. No secondaries stored.
3592 fGcphys->ipair = par;
3596 //_____________________________________________________________________________
3597 void TGeant3::SetPFIS(Int_t par)
3600 // To control photo fission mechanism.
3601 // par =0 no photo fission.
3602 // =1 Photo fission. secondaries processed.
3603 // =2 Photo fission. No secondaries stored.
3605 fGcphys->ipfis = par;
3608 //_____________________________________________________________________________
3609 void TGeant3::SetPHOT(Int_t par)
3612 // To control Photo effect.
3613 // par =0 no photo electric effect.
3614 // =1 Photo effect. Electron processed.
3615 // =2 Photo effect. No electron stored.
3617 fGcphys->iphot = par;
3620 //_____________________________________________________________________________
3621 void TGeant3::SetRAYL(Int_t par)
3624 // To control Rayleigh scattering.
3625 // par =0 no Rayleigh scattering.
3628 fGcphys->irayl = par;
3631 //_____________________________________________________________________________
3632 void TGeant3::SetSTRA(Int_t par)
3635 // To control energy loss fluctuations
3636 // with the PhotoAbsorption Ionisation model.
3637 // par =0 no Straggling.
3638 // =1 Straggling yes => no Delta rays.
3640 fGcphlt->istra = par;
3643 //_____________________________________________________________________________
3644 void TGeant3::SetSWIT(Int_t sw, Int_t val)
3648 // val New switch value
3650 // Change one element of array ISWIT(10) in /GCFLAG/
3652 if (sw <= 0 || sw > 10) return;
3653 fGcflag->iswit[sw-1] = val;
3657 //_____________________________________________________________________________
3658 void TGeant3::SetTRIG(Int_t nevents)
3661 // Set number of events to be run
3663 fGcflag->nevent = nevents;
3666 //_____________________________________________________________________________
3667 void TGeant3::SetUserDecay(Int_t pdg)
3670 // Force the decays of particles to be done with Pythia
3671 // and not with the Geant routines.
3672 // just kill pointers doing mzdrop
3674 Int_t ipart = IdFromPDG(pdg);
3676 printf("Particle %d not in geant\n",pdg);
3679 Int_t jpart=fGclink->jpart;
3680 Int_t jpa=fZlq[jpart-ipart];
3683 Int_t jpa1=fZlq[jpa-1];
3685 mzdrop(fGcbank->ixcons,jpa1,PASSCHARD(" ") PASSCHARL(" "));
3686 Int_t jpa2=fZlq[jpa-2];
3688 mzdrop(fGcbank->ixcons,jpa2,PASSCHARD(" ") PASSCHARL(" "));
3692 //______________________________________________________________________________
3693 void TGeant3::Vname(const char *name, char *vname)
3696 // convert name to upper case. Make vname at least 4 chars
3698 Int_t l = strlen(name);
3701 for (i=0;i<l;i++) vname[i] = toupper(name[i]);
3702 for (i=l;i<4;i++) vname[i] = ' ';
3706 //______________________________________________________________________________
3707 void TGeant3::Ertrgo()
3710 // Perform the tracking of the track Track parameters are in VECT
3715 //______________________________________________________________________________
3716 void TGeant3::Ertrak(const Float_t *const x1, const Float_t *const p1,
3717 const Float_t *x2, const Float_t *p2,
3718 Int_t ipa, Option_t *chopt)
3720 //************************************************************************
3722 //* Perform the tracking of the track from point X1 to *
3724 //* (Before calling this routine the user should also provide *
3725 //* the input informations in /EROPTS/ and /ERTRIO/ *
3726 //* using subroutine EUFIL(L/P/V) *
3727 //* X1 - Starting coordinates (Cartesian) *
3728 //* P1 - Starting 3-momentum (Cartesian) *
3729 //* X2 - Final coordinates (Cartesian) *
3730 //* P2 - Final 3-momentum (Cartesian) *
3731 //* IPA - Particle code (a la GEANT) of the track *
3734 //* 'B' 'Backward tracking' - i.e. energy loss *
3735 //* added to the current energy *
3736 //* 'E' 'Exact' calculation of errors assuming *
3737 //* helix (i.e. pathlength not *
3738 //* assumed as infinitesimal) *
3739 //* 'L' Tracking upto prescribed Lengths reached *
3740 //* 'M' 'Mixed' prediction (not yet coded) *
3741 //* 'O' Tracking 'Only' without calculating errors *
3742 //* 'P' Tracking upto prescribed Planes reached *
3743 //* 'V' Tracking upto prescribed Volumes reached *
3744 //* 'X' Tracking upto prescribed Point approached *
3746 //* Interface with GEANT : *
3747 //* Track parameters are in /CGKINE/ and /GCTRAK/ *
3749 //* ==>Called by : USER *
3750 //* Authors M.Maire, E.Nagy ********//* *
3752 //************************************************************************
3753 ertrak(x1,p1,x2,p2,ipa,PASSCHARD(chopt) PASSCHARL(chopt));
3756 //_____________________________________________________________________________
3757 void TGeant3::WriteEuclid(const char* filnam, const char* topvol,
3758 Int_t number, Int_t nlevel)
3762 // ******************************************************************
3764 // * Write out the geometry of the detector in EUCLID file format *
3766 // * filnam : will be with the extension .euc *
3767 // * topvol : volume name of the starting node *
3768 // * number : copy number of topvol (relevant for gsposp) *
3769 // * nlevel : number of levels in the tree structure *
3770 // * to be written out, starting from topvol *
3772 // * Author : M. Maire *
3774 // ******************************************************************
3776 // File filnam.tme is written out with the definitions of tracking
3777 // medias and materials.
3778 // As to restore original numbers for materials and medias, program
3779 // searches in the file euc_medi.dat and comparing main parameters of
3780 // the mat. defined inside geant and the one in file recognizes them
3781 // and is able to take number from file. If for any material or medium,
3782 // this procedure fails, ordering starts from 1.
3783 // Arrays IOTMED and IOMATE are used for this procedure
3785 const char kShape[][5]={"BOX ","TRD1","TRD2","TRAP","TUBE","TUBS","CONE",
3786 "CONS","SPHE","PARA","PGON","PCON","ELTU","HYPE",
3788 Int_t i, end, itm, irm, jrm, k, nmed;
3792 char *filext, *filetme;
3793 char natmed[21], namate[21];
3794 char natmedc[21], namatec[21];
3795 char key[5], name[5], mother[5], konly[5];
3797 Int_t iadvol, iadtmd, iadrot, nwtot, iret;
3798 Int_t mlevel, numbr, natt, numed, nin, ndata;
3799 Int_t iname, ivo, ish, jvo, nvstak, ivstak;
3800 Int_t jdiv, ivin, in, jin, jvin, irot;
3801 Int_t jtm, imat, jma, flag=0, imatc;
3802 Float_t az, dens, radl, absl, a, step, x, y, z;
3803 Int_t npar, ndvmx, left;
3804 Float_t zc, densc, radlc, abslc, c0, tmaxfd;
3806 Int_t iomate[100], iotmed[100];
3807 Float_t par[50], att[20], ubuf[50];
3810 Int_t level, ndiv, iaxe;
3811 Int_t itmedc, nmatc, isvolc, ifieldc, nwbufc, isvol, nmat, ifield, nwbuf;
3812 Float_t fieldmc, tmaxfdc, stemaxc, deemaxc, epsilc, stminc, fieldm;
3813 Float_t tmaxf, stemax, deemax, epsil, stmin;
3814 const char *k10000="!\n%s\n!\n";
3815 //Open the input file
3817 for(i=0;i<end;i++) if(filnam[i]=='.') {
3821 filext=new char[end+5];
3822 filetme=new char[end+5];
3823 strncpy(filext,filnam,end);
3824 strncpy(filetme,filnam,end);
3826 // *** The output filnam name will be with extension '.euc'
3827 strcpy(&filext[end],".euc");
3828 strcpy(&filetme[end],".tme");
3829 lun=fopen(filext,"w");
3831 // *** Initialisation of the working space
3832 iadvol=fGcnum->nvolum;
3833 iadtmd=iadvol+fGcnum->nvolum;
3834 iadrot=iadtmd+fGcnum->ntmed;
3835 if(fGclink->jrotm) {
3836 fGcnum->nrotm=fZiq[fGclink->jrotm-2];
3840 nwtot=iadrot+fGcnum->nrotm;
3841 qws = new float[nwtot+1];
3842 for (i=0;i<nwtot+1;i++) qws[i]=0;
3845 if(nlevel==0) mlevel=20;
3847 // *** find the top volume and put it in the stak
3848 numbr = number>0 ? number : 1;
3849 Gfpara(topvol,numbr,1,npar,natt,par,att);
3851 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3856 // *** authorized shape ?
3857 strncpy((char *)&iname, topvol, 4);
3859 for(i=1; i<=fGcnum->nvolum; i++) if(fZiq[fGclink->jvolum+i]==iname) {
3863 jvo = fZlq[fGclink->jvolum-ivo];
3864 ish = Int_t (fZq[jvo+2]);
3866 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3873 iws[iadvol+ivo] = level;
3876 //*** flag all volumes and fill the stak
3880 // pick the next volume in stak
3882 ivo = TMath::Abs(iws[ivstak]);
3883 jvo = fZlq[fGclink->jvolum - ivo];
3885 // flag the tracking medium
3886 numed = Int_t (fZq[jvo + 4]);
3887 iws[iadtmd + numed] = 1;
3889 // get the daughters ...
3890 level = iws[iadvol+ivo];
3891 if (level < mlevel) {
3893 nin = Int_t (fZq[jvo + 3]);
3895 // from division ...
3897 jdiv = fZlq[jvo - 1];
3898 ivin = Int_t (fZq[jdiv + 2]);
3900 iws[nvstak] = -ivin;
3901 iws[iadvol+ivin] = level;
3903 // from position ...
3904 } else if (nin > 0) {
3905 for(in=1; in<=nin; in++) {
3906 jin = fZlq[jvo - in];
3907 ivin = Int_t (fZq[jin + 2 ]);
3908 jvin = fZlq[fGclink->jvolum - ivin];
3909 ish = Int_t (fZq[jvin + 2]);
3910 // authorized shape ?
3912 // not yet flagged ?
3913 if (iws[iadvol+ivin]==0) {
3916 iws[iadvol+ivin] = level;
3918 // flag the rotation matrix
3919 irot = Int_t ( fZq[jin + 4 ]);
3920 if (irot > 0) iws[iadrot+irot] = 1;
3926 // next volume in stak ?
3927 if (ivstak < nvstak) goto L10;
3929 // *** restore original material and media numbers
3930 // file euc_medi.dat is needed to compare materials and medias
3932 FILE* luncor=fopen("euc_medi.dat","r");
3935 for(itm=1; itm<=fGcnum->ntmed; itm++) {
3936 if (iws[iadtmd+itm] > 0) {
3937 jtm = fZlq[fGclink->jtmed-itm];
3938 strncpy(natmed,(char *)&fZiq[jtm+1],20);
3939 imat = Int_t (fZq[jtm+6]);
3940 jma = fZlq[fGclink->jmate-imat];
3942 printf(" *** GWEUCL *** material not defined for tracking medium %5i %s\n",itm,natmed);
3945 strncpy(namate,(char *)&fZiq[jma+1],20);
3948 //** find the material original number
3951 iret=fscanf(luncor,"%4s,%130s",key,card);
3952 if(iret<=0) goto L26;
3954 if(!strcmp(key,"MATE")) {
3955 sscanf(card,"%d %s %f %f %f %f %f %d",&imatc,namatec,&az,&zc,&densc,&radlc,&abslc,&nparc);
3956 Gfmate(imat,namate,a,z,dens,radl,absl,par,npar);
3957 if(!strcmp(namatec,namate)) {
3958 if(az==a && zc==z && densc==dens && radlc==radl
3959 && abslc==absl && nparc==nparc) {
3962 printf("*** GWEUCL *** material : %3d '%s' restored as %3d\n",imat,namate,imatc);
3964 printf("*** GWEUCL *** different definitions for material: %s\n",namate);
3968 if(strcmp(key,"END") && !flag) goto L23;
3970 printf("*** GWEUCL *** cannot restore original number for material: %s\n",namate);
3974 //*** restore original tracking medium number
3977 iret=fscanf(luncor,"%4s,%130s",key,card);
3978 if(iret<=0) goto L26;
3980 if (!strcmp(key,"TMED")) {
3981 sscanf(card,"%d %s %d %d %d %f %f %f %f %f %f %d\n",
3982 &itmedc,natmedc,&nmatc,&isvolc,&ifieldc,&fieldmc,
3983 &tmaxfdc,&stemaxc,&deemaxc,&epsilc,&stminc,&nwbufc);
3984 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxf,stemax,deemax,
3985 epsil,stmin,ubuf,&nwbuf);
3986 if(!strcmp(natmedc,natmed)) {
3987 if (iomate[nmat]==nmatc && nwbuf==nwbufc) {
3990 printf("*** GWEUCL *** medium : %3d '%20s' restored as %3d\n",
3993 printf("*** GWEUCL *** different definitions for tracking medium: %s\n",natmed);
3997 if(strcmp(key,"END") && !flag) goto L24;
3999 printf("cannot restore original number for medium : %s\n",natmed);
4007 L26: printf("*** GWEUCL *** cannot read the data file\n");
4009 L29: if(luncor) fclose (luncor);
4012 // *** write down the tracking medium definition
4014 strcpy(card,"! Tracking medium");
4015 fprintf(lun,k10000,card);
4017 for(itm=1;itm<=fGcnum->ntmed;itm++) {
4018 if (iws[iadtmd+itm]>0) {
4019 jtm = fZlq[fGclink->jtmed-itm];
4020 strncpy(natmed,(char *)&fZiq[jtm+1],20);
4022 imat = Int_t (fZq[jtm+6]);
4023 jma = fZlq[fGclink->jmate-imat];
4024 //* order media from one, if comparing with database failed
4026 iotmed[itm]=++imxtmed;
4027 iomate[imat]=++imxmate;
4032 printf(" *** GWEUCL *** material not defined for tracking medium %5d %s\n",
4035 strncpy(namate,(char *)&fZiq[jma+1],20);
4038 fprintf(lun,"TMED %3d '%20s' %3d '%20s'\n",iotmed[itm],natmed,iomate[imat],namate);
4042 //* *** write down the rotation matrix
4044 strcpy(card,"! Reperes");
4045 fprintf(lun,k10000,card);
4047 for(irm=1;irm<=fGcnum->nrotm;irm++) {
4048 if (iws[iadrot+irm]>0) {
4049 jrm = fZlq[fGclink->jrotm-irm];
4050 fprintf(lun,"ROTM %3d",irm);
4051 for(k=11;k<=16;k++) fprintf(lun," %8.3f",fZq[jrm+k]);
4056 //* *** write down the volume definition
4058 strcpy(card,"! Volumes");
4059 fprintf(lun,k10000,card);
4061 for(ivstak=1;ivstak<=nvstak;ivstak++) {
4064 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivo],4);
4066 jvo = fZlq[fGclink->jvolum-ivo];
4067 ish = Int_t (fZq[jvo+2]);
4068 nmed = Int_t (fZq[jvo+4]);
4069 npar = Int_t (fZq[jvo+5]);
4071 if (ivstak>1) for(i=0;i<npar;i++) par[i]=fZq[jvo+7+i];
4072 Gckpar (ish,npar,par);
4073 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,kShape[ish-1],iotmed[nmed],npar);
4074 for(i=0;i<(npar-1)/6+1;i++) {
4077 for(k=0;k<(left<6?left:6);k++) fprintf(lun," %11.5f",par[i*6+k]);
4081 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,kShape[ish-1],iotmed[nmed],npar);
4086 //* *** write down the division of volumes
4088 fprintf(lun,k10000,"! Divisions");
4089 for(ivstak=1;ivstak<=nvstak;ivstak++) {
4090 ivo = TMath::Abs(iws[ivstak]);
4091 jvo = fZlq[fGclink->jvolum-ivo];
4092 ish = Int_t (fZq[jvo+2]);
4093 nin = Int_t (fZq[jvo+3]);
4094 //* this volume is divided ...
4097 iaxe = Int_t ( fZq[jdiv+1]);
4098 ivin = Int_t ( fZq[jdiv+2]);
4099 ndiv = Int_t ( fZq[jdiv+3]);
4102 jvin = fZlq[fGclink->jvolum-ivin];
4103 nmed = Int_t ( fZq[jvin+4]);
4104 strncpy(mother,(char *)&fZiq[fGclink->jvolum+ivo ],4);
4106 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivin],4);
4108 if ((step<=0.)||(ish>=11)) {
4109 //* volume with negative parameter or gsposp or pgon ...
4110 fprintf(lun,"DIVN '%4s' '%4s' %3d %3d\n",name,mother,ndiv,iaxe);
4111 } else if ((ndiv<=0)||(ish==10)) {
4112 //* volume with negative parameter or gsposp or para ...
4113 ndvmx = TMath::Abs(ndiv);
4114 fprintf(lun,"DIVT '%4s' '%4s' %11.5f %3d %3d %3d\n",
4115 name,mother,step,iaxe,iotmed[nmed],ndvmx);
4117 //* normal volume : all kind of division are equivalent
4118 fprintf(lun,"DVT2 '%4s' '%4s' %11.5f %3d %11.5f %3d %3d\n",
4119 name,mother,step,iaxe,c0,iotmed[nmed],ndiv);
4124 //* *** write down the the positionnement of volumes
4126 fprintf(lun,k10000,"! Positionnements\n");
4128 for(ivstak = 1;ivstak<=nvstak;ivstak++) {
4129 ivo = TMath::Abs(iws[ivstak]);
4130 strncpy(mother,(char*)&fZiq[fGclink->jvolum+ivo ],4);
4132 jvo = fZlq[fGclink->jvolum-ivo];
4133 nin = Int_t( fZq[jvo+3]);
4134 //* this volume has daughters ...
4136 for (in=1;in<=nin;in++) {
4138 ivin = Int_t (fZq[jin +2]);
4139 numb = Int_t (fZq[jin +3]);
4140 irot = Int_t (fZq[jin +4]);
4144 strcpy(konly,"ONLY");
4145 if (fZq[jin+8]!=1.) strcpy(konly,"MANY");
4146 strncpy(name,(char*)&fZiq[fGclink->jvolum+ivin],4);
4148 jvin = fZlq[fGclink->jvolum-ivin];
4149 ish = Int_t (fZq[jvin+2]);
4150 //* gspos or gsposp ?
4151 ndata = fZiq[jin-1];
4153 fprintf(lun,"POSI '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s'\n",
4154 name,numb,mother,x,y,z,irot,konly);
4156 npar = Int_t (fZq[jin+9]);
4157 for(i=0;i<npar;i++) par[i]=fZq[jin+10+i];
4158 Gckpar (ish,npar,par);
4159 fprintf(lun,"POSP '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s' %3d\n",
4160 name,numb,mother,x,y,z,irot,konly,npar);
4162 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
4169 fprintf(lun,"END\n");
4172 //****** write down the materials and medias *****
4174 lun=fopen(filetme,"w");
4176 for(itm=1;itm<=fGcnum->ntmed;itm++) {
4177 if (iws[iadtmd+itm]>0) {
4178 jtm = fZlq[fGclink->jtmed-itm];
4179 strncpy(natmed,(char*)&fZiq[jtm+1],4);
4180 imat = Int_t (fZq[jtm+6]);
4181 jma = Int_t (fZlq[fGclink->jmate-imat]);
4183 Gfmate (imat,namate,a,z,dens,radl,absl,par,npar);
4184 fprintf(lun,"MATE %4d '%20s'%11.5E %11.5E %11.5E %11.5E %11.5E %3d\n",
4185 iomate[imat],namate,a,z,dens,radl,absl,npar);
4189 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
4193 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxfd,stemax,deemax,epsil,stmin,par,&npar);
4194 fprintf(lun,"TMED %4d '%20s' %3d %1d %3d %11.5f %11.5f %11.5f %11.5f %11.5f %11.5f %3d\n",
4195 iotmed[itm],natmed,iomate[nmat],isvol,ifield,
4196 fieldm,tmaxfd,stemax,deemax,epsil,stmin,npar);
4200 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
4206 fprintf(lun,"END\n");
4208 printf(" *** GWEUCL *** file: %s is now written out\n",filext);
4209 printf(" *** GWEUCL *** file: %s is now written out\n",filetme);
4218 //_____________________________________________________________________________
4219 void TGeant3::Streamer(TBuffer &R__b)
4222 // Stream an object of class TGeant3.
4224 if (R__b.IsReading()) {
4225 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
4226 AliMC::Streamer(R__b);
4229 R__b.ReadStaticArray(fPDGCode);
4231 R__b.WriteVersion(TGeant3::IsA());
4232 AliMC::Streamer(R__b);
4235 R__b.WriteArray(fPDGCode, fNPDGCodes);