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.37 2000/10/02 21:28:16 fca
19 Removal of useless dependecies via forward declarations
21 Revision 1.36 2000/09/14 07:08:41 fca
22 Introducing glvolu in the interface
24 Revision 1.35 2000/09/12 14:27:10 morsch
25 No instance of AliDecayer created to initialize fDecayer.
27 Revision 1.34 2000/09/07 12:12:01 morsch
28 Comment inside comment removed.
30 Revision 1.33 2000/09/06 16:03:42 morsch
31 Set ExternalDecayer, Decayer and SetForceDecay methods added.
32 Gspart calls for charmed and bottom hadrons added.
33 Decay mode definitions for charmed and beauty hadrons have been taken out.
34 They will be handled by an external decayer.
36 Revision 1.32 2000/08/24 16:28:53 hristov
37 TGeant3::IsNewTrack corrected by F.Carminati
39 Revision 1.31 2000/07/13 16:19:10 fca
40 Mainly coding conventions + some small bug fixes
42 Revision 1.30 2000/07/12 08:56:30 fca
43 Coding convention correction and warning removal
45 Revision 1.29 2000/07/11 18:24:59 fca
46 Coding convention corrections + few minor bug fixes
48 Revision 1.28 2000/06/29 10:51:55 morsch
49 Add some charmed and bottom baryons to the particle list (TDatabasePDG). This
50 is needed by Hijing. Should be part of a future review of TDatabasePDG.
52 Revision 1.27 2000/06/21 17:40:15 fca
53 Adding possibility to set ISTRA, PAI model
55 Revision 1.26 2000/05/16 13:10:41 fca
56 New method IsNewTrack and fix for a problem in Father-Daughter relations
58 Revision 1.25 2000/04/07 11:12:35 fca
59 G4 compatibility changes
61 Revision 1.24 2000/02/28 21:03:57 fca
62 Some additions to improve the compatibility with G4
64 Revision 1.23 2000/02/23 16:25:25 fca
65 AliVMC and AliGeant3 classes introduced
66 ReadEuclid moved from AliRun to AliModule
68 Revision 1.22 2000/01/18 15:40:13 morsch
69 Interface to GEANT3 routines GFTMAT, GBRELM and GPRELM added
70 Define geant particle type 51: Feedback Photon with Cherenkov photon properties.
72 Revision 1.21 2000/01/17 19:41:17 fca
75 Revision 1.20 2000/01/12 11:29:27 fca
78 Revision 1.19 1999/12/17 09:03:12 fca
79 Introduce a names array
81 Revision 1.18 1999/11/26 16:55:39 fca
82 Reimplement CurrentVolName() to avoid memory leaks
84 Revision 1.17 1999/11/03 16:58:28 fca
85 Correct source of address violation in creating character string
87 Revision 1.16 1999/11/03 13:17:08 fca
88 Have ProdProcess return const char*
90 Revision 1.15 1999/10/26 06:04:50 fca
91 Introduce TLorentzVector in AliMC::GetSecondary. Thanks to I.Hrivnacova
93 Revision 1.14 1999/09/29 09:24:30 fca
94 Introduction of the Copyright and cvs Log
98 ///////////////////////////////////////////////////////////////////////////////
100 // Interface Class to the Geant3.21 MonteCarlo //
104 <img src="picts/TGeant3Class.gif">
109 ///////////////////////////////////////////////////////////////////////////////
115 #include <TDatabasePDG.h>
116 #include "AliCallf77.h"
117 #include "AliDecayer.h"
118 #include "TLorentzVector.h"
121 # define gzebra gzebra_
122 # define grfile grfile_
123 # define gpcxyz gpcxyz_
124 # define ggclos ggclos_
125 # define glast glast_
126 # define ginit ginit_
127 # define gcinit gcinit_
129 # define gtrig gtrig_
130 # define gtrigc gtrigc_
131 # define gtrigi gtrigi_
132 # define gwork gwork_
133 # define gzinit gzinit_
134 # define gfmate gfmate_
135 # define gfpart gfpart_
136 # define gftmed gftmed_
137 # define gftmat gftmat_
138 # define gmate gmate_
139 # define gpart gpart_
141 # define gsmate gsmate_
142 # define gsmixt gsmixt_
143 # define gspart gspart_
144 # define gstmed gstmed_
145 # define gsckov gsckov_
146 # define gstpar gstpar_
147 # define gfkine gfkine_
148 # define gfvert gfvert_
149 # define gskine gskine_
150 # define gsvert gsvert_
151 # define gphysi gphysi_
152 # define gdebug gdebug_
153 # define gekbin gekbin_
154 # define gfinds gfinds_
155 # define gsking gsking_
156 # define gskpho gskpho_
157 # define gsstak gsstak_
158 # define gsxyz gsxyz_
159 # define gtrack gtrack_
160 # define gtreve gtreve_
161 # define gtreveroot gtreveroot_
162 # define grndm grndm_
163 # define grndmq grndmq_
164 # define gdtom gdtom_
165 # define glmoth glmoth_
166 # define gmedia gmedia_
167 # define gmtod gmtod_
168 # define gsdvn gsdvn_
169 # define gsdvn2 gsdvn2_
170 # define gsdvs gsdvs_
171 # define gsdvs2 gsdvs2_
172 # define gsdvt gsdvt_
173 # define gsdvt2 gsdvt2_
174 # define gsord gsord_
175 # define gspos gspos_
176 # define gsposp gsposp_
177 # define gsrotm gsrotm_
178 # define gprotm gprotm_
179 # define gsvolu gsvolu_
180 # define gprint gprint_
181 # define gdinit gdinit_
182 # define gdopt gdopt_
183 # define gdraw gdraw_
184 # define gdrayt gdrayt_
185 # define gdrawc gdrawc_
186 # define gdrawx gdrawx_
187 # define gdhead gdhead_
188 # define gdwmn1 gdwmn1_
189 # define gdwmn2 gdwmn2_
190 # define gdwmn3 gdwmn3_
191 # define gdxyz gdxyz_
192 # define gdcxyz gdcxyz_
193 # define gdman gdman_
194 # define gdspec gdspec_
195 # define gdtree gdtree_
196 # define gdelet gdelet_
197 # define gdclos gdclos_
198 # define gdshow gdshow_
199 # define gdopen gdopen_
200 # define dzshow dzshow_
201 # define gsatt gsatt_
202 # define gfpara gfpara_
203 # define gckpar gckpar_
204 # define gckmat gckmat_
205 # define glvolu glvolu_
206 # define geditv geditv_
207 # define mzdrop mzdrop_
209 # define ertrak ertrak_
210 # define ertrgo ertrgo_
212 # define setbomb setbomb_
213 # define setclip setclip_
214 # define gcomad gcomad_
216 # define gbrelm gbrelm_
217 # define gprelm gprelm_
219 # define gzebra GZEBRA
220 # define grfile GRFILE
221 # define gpcxyz GPCXYZ
222 # define ggclos GGCLOS
225 # define gcinit GCINIT
228 # define gtrigc GTRIGC
229 # define gtrigi GTRIGI
231 # define gzinit GZINIT
232 # define gfmate GFMATE
233 # define gfpart GFPART
234 # define gftmed GFTMED
235 # define gftmat GFTMAT
239 # define gsmate GSMATE
240 # define gsmixt GSMIXT
241 # define gspart GSPART
242 # define gstmed GSTMED
243 # define gsckov GSCKOV
244 # define gstpar GSTPAR
245 # define gfkine GFKINE
246 # define gfvert GFVERT
247 # define gskine GSKINE
248 # define gsvert GSVERT
249 # define gphysi GPHYSI
250 # define gdebug GDEBUG
251 # define gekbin GEKBIN
252 # define gfinds GFINDS
253 # define gsking GSKING
254 # define gskpho GSKPHO
255 # define gsstak GSSTAK
257 # define gtrack GTRACK
258 # define gtreve GTREVE
259 # define gtreveroot GTREVEROOT
261 # define grndmq GRNDMQ
263 # define glmoth GLMOTH
264 # define gmedia GMEDIA
267 # define gsdvn2 GSDVN2
269 # define gsdvs2 GSDVS2
271 # define gsdvt2 GSDVT2
274 # define gsposp GSPOSP
275 # define gsrotm GSROTM
276 # define gprotm GPROTM
277 # define gsvolu GSVOLU
278 # define gprint GPRINT
279 # define gdinit GDINIT
282 # define gdrayt GDRAYT
283 # define gdrawc GDRAWC
284 # define gdrawx GDRAWX
285 # define gdhead GDHEAD
286 # define gdwmn1 GDWMN1
287 # define gdwmn2 GDWMN2
288 # define gdwmn3 GDWMN3
290 # define gdcxyz GDCXYZ
292 # define gdfspc GDFSPC
293 # define gdspec GDSPEC
294 # define gdtree GDTREE
295 # define gdelet GDELET
296 # define gdclos GDCLOS
297 # define gdshow GDSHOW
298 # define gdopen GDOPEN
299 # define dzshow DZSHOW
301 # define gfpara GFPARA
302 # define gckpar GCKPAR
303 # define gckmat GCKMAT
304 # define glvolu GLVOLU
305 # define geditv GEDITV
306 # define mzdrop MZDROP
308 # define ertrak ERTRAK
309 # define ertrgo ERTRGO
311 # define setbomb SETBOMB
312 # define setclip SETCLIP
313 # define gcomad GCOMAD
315 # define gbrelm GBRELM
316 # define gprelm GPRELM
320 //____________________________________________________________________________
324 // Prototypes for GEANT functions
326 void type_of_call gzebra(const int&);
328 void type_of_call gpcxyz();
330 void type_of_call ggclos();
332 void type_of_call glast();
334 void type_of_call ginit();
336 void type_of_call gcinit();
338 void type_of_call grun();
340 void type_of_call gtrig();
342 void type_of_call gtrigc();
344 void type_of_call gtrigi();
346 void type_of_call gwork(const int&);
348 void type_of_call gzinit();
350 void type_of_call gmate();
352 void type_of_call gpart();
354 void type_of_call gsdk(Int_t &, Float_t *, Int_t *);
356 void type_of_call gfkine(Int_t &, Float_t *, Float_t *, Int_t &,
357 Int_t &, Float_t *, Int_t &);
359 void type_of_call gfvert(Int_t &, Float_t *, Int_t &, Int_t &,
360 Float_t &, Float_t *, Int_t &);
362 void type_of_call gskine(Float_t *,Int_t &, Int_t &, Float_t *,
365 void type_of_call gsvert(Float_t *,Int_t &, Int_t &, Float_t *,
368 void type_of_call gphysi();
370 void type_of_call gdebug();
372 void type_of_call gekbin();
374 void type_of_call gfinds();
376 void type_of_call gsking(Int_t &);
378 void type_of_call gskpho(Int_t &);
380 void type_of_call gsstak(Int_t &);
382 void type_of_call gsxyz();
384 void type_of_call gtrack();
386 void type_of_call gtreve();
388 void type_of_call gtreveroot();
390 void type_of_call grndm(Float_t *, const Int_t &);
392 void type_of_call grndmq(Int_t &, Int_t &, const Int_t &,
395 void type_of_call gdtom(Float_t *, Float_t *, Int_t &);
397 void type_of_call glmoth(DEFCHARD, Int_t &, Int_t &, Int_t *,
398 Int_t *, Int_t * DEFCHARL);
400 void type_of_call gmedia(Float_t *, Int_t &);
402 void type_of_call gmtod(Float_t *, Float_t *, Int_t &);
404 void type_of_call gsrotm(const Int_t &, const Float_t &, const Float_t &,
405 const Float_t &, const Float_t &, const Float_t &,
408 void type_of_call gprotm(const Int_t &);
410 void type_of_call grfile(const Int_t&, DEFCHARD,
411 DEFCHARD DEFCHARL DEFCHARL);
413 void type_of_call gfmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
414 Float_t &, Float_t &, Float_t &, Float_t *,
417 void type_of_call gfpart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
418 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
420 void type_of_call gftmed(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 gftmat(const Int_t&, const Int_t&, DEFCHARD, const Int_t&,
426 ,Float_t *, Int_t & DEFCHARL);
428 void type_of_call gsmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
429 Float_t &, Float_t &, Float_t &, Float_t *,
432 void type_of_call gsmixt(const Int_t&, DEFCHARD, Float_t *, Float_t *,
433 Float_t &, Int_t &, Float_t * DEFCHARL);
435 void type_of_call gspart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
436 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
439 void type_of_call gstmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
440 Float_t &, Float_t &, Float_t &, Float_t &,
441 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
443 void type_of_call gsckov(Int_t &itmed, Int_t &npckov, Float_t *ppckov,
444 Float_t *absco, Float_t *effic, Float_t *rindex);
445 void type_of_call gstpar(const Int_t&, DEFCHARD, Float_t & DEFCHARL);
447 void type_of_call gsdvn(DEFCHARD,DEFCHARD, Int_t &, Int_t &
450 void type_of_call gsdvn2(DEFCHARD,DEFCHARD, Int_t &, Int_t &, Float_t &,
451 Int_t & DEFCHARL DEFCHARL);
453 void type_of_call gsdvs(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &
456 void type_of_call gsdvs2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t &,
457 Int_t & DEFCHARL DEFCHARL);
459 void type_of_call gsdvt(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &,
460 Int_t & DEFCHARL DEFCHARL);
462 void type_of_call gsdvt2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t&,
463 Int_t &, Int_t & DEFCHARL DEFCHARL);
465 void type_of_call gsord(DEFCHARD, Int_t & DEFCHARL);
467 void type_of_call gspos(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
468 Float_t &, Int_t &, DEFCHARD DEFCHARL DEFCHARL
471 void type_of_call gsposp(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
472 Float_t &, Int_t &, DEFCHARD,
473 Float_t *, Int_t & DEFCHARL DEFCHARL DEFCHARL);
475 void type_of_call gsvolu(DEFCHARD, DEFCHARD, Int_t &, Float_t *, Int_t &,
476 Int_t & DEFCHARL DEFCHARL);
478 void type_of_call gsatt(DEFCHARD, DEFCHARD, Int_t & DEFCHARL DEFCHARL);
480 void type_of_call gfpara(DEFCHARD , Int_t&, Int_t&, Int_t&, Int_t&, Float_t*,
483 void type_of_call gckpar(Int_t&, Int_t&, Float_t*);
485 void type_of_call gckmat(Int_t&, DEFCHARD DEFCHARL);
487 void type_of_call glvolu(Int_t&, Int_t*, Int_t*, Int_t&);
489 void type_of_call gprint(DEFCHARD,const int& DEFCHARL);
491 void type_of_call gdinit();
493 void type_of_call gdopt(DEFCHARD,DEFCHARD DEFCHARL DEFCHARL);
495 void type_of_call gdraw(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
496 Float_t &, Float_t &, Float_t & DEFCHARL);
497 void type_of_call gdrayt(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
498 Float_t &, Float_t &, Float_t & DEFCHARL);
499 void type_of_call gdrawc(DEFCHARD,Int_t &, Float_t &, Float_t &, Float_t &,
500 Float_t &, Float_t & DEFCHARL);
501 void type_of_call gdrawx(DEFCHARD,Float_t &, Float_t &, Float_t &, Float_t &,
502 Float_t &, Float_t &, Float_t &, Float_t &,
504 void type_of_call gdhead(Int_t &,DEFCHARD, Float_t & DEFCHARL);
505 void type_of_call gdxyz(Int_t &);
506 void type_of_call gdcxyz();
507 void type_of_call gdman(Float_t &, Float_t &);
508 void type_of_call gdwmn1(Float_t &, Float_t &);
509 void type_of_call gdwmn2(Float_t &, Float_t &);
510 void type_of_call gdwmn3(Float_t &, Float_t &);
511 void type_of_call gdspec(DEFCHARD DEFCHARL);
512 void type_of_call gdfspc(DEFCHARD, Int_t &, Int_t & DEFCHARL) {;}
513 void type_of_call gdtree(DEFCHARD, Int_t &, Int_t & DEFCHARL);
515 void type_of_call gdopen(Int_t &);
516 void type_of_call gdclos();
517 void type_of_call gdelet(Int_t &);
518 void type_of_call gdshow(Int_t &);
519 void type_of_call geditv(Int_t &) {;}
522 void type_of_call dzshow(DEFCHARD,const int&,const int&,DEFCHARD,const int&,
523 const int&, const int&, const int& DEFCHARL
526 void type_of_call mzdrop(Int_t&, Int_t&, DEFCHARD DEFCHARL);
528 void type_of_call setbomb(Float_t &);
529 void type_of_call setclip(DEFCHARD, Float_t &,Float_t &,Float_t &,Float_t &,
530 Float_t &, Float_t & DEFCHARL);
531 void type_of_call gcomad(DEFCHARD, Int_t*& DEFCHARL);
533 void type_of_call ertrak(const Float_t *const x1, const Float_t *const p1,
534 const Float_t *x2, const Float_t *p2,
535 const Int_t &ipa, DEFCHARD DEFCHARL);
537 void type_of_call ertrgo();
539 float type_of_call gbrelm(const Float_t &z, const Float_t& t, const Float_t& cut);
540 float type_of_call gprelm(const Float_t &z, const Float_t& t, const Float_t& cut);
544 // Geant3 global pointer
546 static const Int_t kDefSize = 600;
550 //____________________________________________________________________________
554 // Default constructor
558 //____________________________________________________________________________
559 TGeant3::TGeant3(const char *title, Int_t nwgeant)
560 :AliMC("TGeant3",title)
563 // Standard constructor for TGeant3 with ZEBRA initialisation
574 // Load Address of Geant3 commons
577 // Zero number of particles
582 //____________________________________________________________________________
583 Int_t TGeant3::CurrentMaterial(Float_t &a, Float_t &z, Float_t &dens,
584 Float_t &radl, Float_t &absl) const
587 // Return the parameters of the current material during transport
591 dens = fGcmate->dens;
592 radl = fGcmate->radl;
593 absl = fGcmate->absl;
594 return 1; //this could be the number of elements in mixture
597 //____________________________________________________________________________
598 void TGeant3::DefaultRange()
601 // Set range of current drawing pad to 20x20 cm
607 gHigz->Range(0,0,20,20);
610 //____________________________________________________________________________
611 void TGeant3::InitHIGZ()
622 //____________________________________________________________________________
623 void TGeant3::LoadAddress()
626 // Assigns the address of the GEANT common blocks to the structures
627 // that allow their access from C++
630 gcomad(PASSCHARD("QUEST"), (int*&) fQuest PASSCHARL("QUEST"));
631 gcomad(PASSCHARD("GCBANK"),(int*&) fGcbank PASSCHARL("GCBANK"));
632 gcomad(PASSCHARD("GCLINK"),(int*&) fGclink PASSCHARL("GCLINK"));
633 gcomad(PASSCHARD("GCCUTS"),(int*&) fGccuts PASSCHARL("GCCUTS"));
634 gcomad(PASSCHARD("GCMULO"),(int*&) fGcmulo PASSCHARL("GCMULO"));
635 gcomad(PASSCHARD("GCFLAG"),(int*&) fGcflag PASSCHARL("GCFLAG"));
636 gcomad(PASSCHARD("GCKINE"),(int*&) fGckine PASSCHARL("GCKINE"));
637 gcomad(PASSCHARD("GCKING"),(int*&) fGcking PASSCHARL("GCKING"));
638 gcomad(PASSCHARD("GCKIN2"),(int*&) fGckin2 PASSCHARL("GCKIN2"));
639 gcomad(PASSCHARD("GCKIN3"),(int*&) fGckin3 PASSCHARL("GCKIN3"));
640 gcomad(PASSCHARD("GCMATE"),(int*&) fGcmate PASSCHARL("GCMATE"));
641 gcomad(PASSCHARD("GCTMED"),(int*&) fGctmed PASSCHARL("GCTMED"));
642 gcomad(PASSCHARD("GCTRAK"),(int*&) fGctrak PASSCHARL("GCTRAK"));
643 gcomad(PASSCHARD("GCTPOL"),(int*&) fGctpol PASSCHARL("GCTPOL"));
644 gcomad(PASSCHARD("GCVOLU"),(int*&) fGcvolu PASSCHARL("GCVOLU"));
645 gcomad(PASSCHARD("GCNUM"), (int*&) fGcnum PASSCHARL("GCNUM"));
646 gcomad(PASSCHARD("GCSETS"),(int*&) fGcsets PASSCHARL("GCSETS"));
647 gcomad(PASSCHARD("GCPHYS"),(int*&) fGcphys PASSCHARL("GCPHYS"));
648 gcomad(PASSCHARD("GCPHLT"),(int*&) fGcphlt PASSCHARL("GCPHLT"));
649 gcomad(PASSCHARD("GCOPTI"),(int*&) fGcopti PASSCHARL("GCOPTI"));
650 gcomad(PASSCHARD("GCTLIT"),(int*&) fGctlit PASSCHARL("GCTLIT"));
651 gcomad(PASSCHARD("GCVDMA"),(int*&) fGcvdma PASSCHARL("GCVDMA"));
654 gcomad(PASSCHARD("ERTRIO"),(int*&) fErtrio PASSCHARL("ERTRIO"));
655 gcomad(PASSCHARD("EROPTS"),(int*&) fEropts PASSCHARL("EROPTS"));
656 gcomad(PASSCHARD("EROPTC"),(int*&) fEroptc PASSCHARL("EROPTC"));
657 gcomad(PASSCHARD("ERWORK"),(int*&) fErwork PASSCHARL("ERWORK"));
659 // Variables for ZEBRA store
660 gcomad(PASSCHARD("IQ"), addr PASSCHARL("IQ"));
662 gcomad(PASSCHARD("LQ"), addr PASSCHARL("LQ"));
667 //_____________________________________________________________________________
668 void TGeant3::GeomIter()
671 // Geometry iterator for moving upward in the geometry tree
672 // Initialise the iterator
674 fNextVol=fGcvolu->nlevel;
677 //____________________________________________________________________________
678 void TGeant3::FinishGeometry()
680 //Close the geometry structure
684 //____________________________________________________________________________
685 Int_t TGeant3::NextVolUp(Text_t *name, Int_t ©)
688 // Geometry iterator for moving upward in the geometry tree
689 // Return next volume up
694 gname=fGcvolu->names[fNextVol];
695 copy=fGcvolu->number[fNextVol];
696 i=fGcvolu->lvolum[fNextVol];
697 name = fVolNames[i-1];
698 if(gname == fZiq[fGclink->jvolum+i]) return i;
699 else printf("GeomTree: Volume %s not found in bank\n",name);
704 //_____________________________________________________________________________
705 void TGeant3::BuildPhysics()
710 //_____________________________________________________________________________
711 Int_t TGeant3::CurrentVolID(Int_t ©) const
714 // Returns the current volume ID and copy number
717 if( (i=fGcvolu->nlevel-1) < 0 ) {
718 Warning("CurrentVolID","Stack depth only %d\n",fGcvolu->nlevel);
720 gname=fGcvolu->names[i];
721 copy=fGcvolu->number[i];
722 i=fGcvolu->lvolum[i];
723 if(gname == fZiq[fGclink->jvolum+i]) return i;
724 else Warning("CurrentVolID","Volume %4s not found\n",(char*)&gname);
729 //_____________________________________________________________________________
730 Int_t TGeant3::CurrentVolOffID(Int_t off, Int_t ©) const
733 // Return the current volume "off" upward in the geometrical tree
734 // ID and copy number
737 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
738 Warning("CurrentVolOffID","Offset requested %d but stack depth %d\n",
739 off,fGcvolu->nlevel);
741 gname=fGcvolu->names[i];
742 copy=fGcvolu->number[i];
743 i=fGcvolu->lvolum[i];
744 if(gname == fZiq[fGclink->jvolum+i]) return i;
745 else Warning("CurrentVolOffID","Volume %4s not found\n",(char*)&gname);
750 //_____________________________________________________________________________
751 const char* TGeant3::CurrentVolName() const
754 // Returns the current volume name
757 if( (i=fGcvolu->nlevel-1) < 0 ) {
758 Warning("CurrentVolName","Stack depth %d\n",fGcvolu->nlevel);
760 gname=fGcvolu->names[i];
761 i=fGcvolu->lvolum[i];
762 if(gname == fZiq[fGclink->jvolum+i]) return fVolNames[i-1];
763 else Warning("CurrentVolName","Volume %4s not found\n",(char*) &gname);
768 //_____________________________________________________________________________
769 const char* TGeant3::CurrentVolOffName(Int_t off) const
772 // Return the current volume "off" upward in the geometrical tree
773 // ID, name and copy number
774 // if name=0 no name is returned
777 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
778 Warning("CurrentVolOffName",
779 "Offset requested %d but stack depth %d\n",off,fGcvolu->nlevel);
781 gname=fGcvolu->names[i];
782 i=fGcvolu->lvolum[i];
783 if(gname == fZiq[fGclink->jvolum+i]) return fVolNames[i-1];
784 else Warning("CurrentVolOffName","Volume %4s not found\n",(char*)&gname);
789 //_____________________________________________________________________________
790 Int_t TGeant3::IdFromPDG(Int_t pdg) const
793 // Return Geant3 code from PDG and pseudo ENDF code
795 for(Int_t i=0;i<fNPDGCodes;++i)
796 if(pdg==fPDGCode[i]) return i;
800 //_____________________________________________________________________________
801 Int_t TGeant3::PDGFromId(Int_t id) const
804 // Return PDG code and pseudo ENDF code from Geant3 code
806 if(id>0 && id<fNPDGCodes) return fPDGCode[id];
810 //_____________________________________________________________________________
811 void TGeant3::DefineParticles()
814 // Define standard Geant 3 particles
817 // Load standard numbers for GEANT particles and PDG conversion
818 fPDGCode[fNPDGCodes++]=-99; // 0 = unused location
819 fPDGCode[fNPDGCodes++]=22; // 1 = photon
820 fPDGCode[fNPDGCodes++]=-11; // 2 = positron
821 fPDGCode[fNPDGCodes++]=11; // 3 = electron
822 fPDGCode[fNPDGCodes++]=12; // 4 = neutrino e
823 fPDGCode[fNPDGCodes++]=-13; // 5 = muon +
824 fPDGCode[fNPDGCodes++]=13; // 6 = muon -
825 fPDGCode[fNPDGCodes++]=111; // 7 = pi0
826 fPDGCode[fNPDGCodes++]=211; // 8 = pi+
827 fPDGCode[fNPDGCodes++]=-211; // 9 = pi-
828 fPDGCode[fNPDGCodes++]=130; // 10 = Kaon Long
829 fPDGCode[fNPDGCodes++]=321; // 11 = Kaon +
830 fPDGCode[fNPDGCodes++]=-321; // 12 = Kaon -
831 fPDGCode[fNPDGCodes++]=2112; // 13 = Neutron
832 fPDGCode[fNPDGCodes++]=2212; // 14 = Proton
833 fPDGCode[fNPDGCodes++]=-2212; // 15 = Anti Proton
834 fPDGCode[fNPDGCodes++]=310; // 16 = Kaon Short
835 fPDGCode[fNPDGCodes++]=221; // 17 = Eta
836 fPDGCode[fNPDGCodes++]=3122; // 18 = Lambda
837 fPDGCode[fNPDGCodes++]=3222; // 19 = Sigma +
838 fPDGCode[fNPDGCodes++]=3212; // 20 = Sigma 0
839 fPDGCode[fNPDGCodes++]=3112; // 21 = Sigma -
840 fPDGCode[fNPDGCodes++]=3322; // 22 = Xi0
841 fPDGCode[fNPDGCodes++]=3312; // 23 = Xi-
842 fPDGCode[fNPDGCodes++]=3334; // 24 = Omega-
843 fPDGCode[fNPDGCodes++]=-2112; // 25 = Anti Proton
844 fPDGCode[fNPDGCodes++]=-3122; // 26 = Anti Proton
845 fPDGCode[fNPDGCodes++]=-3222; // 27 = Anti Sigma -
846 fPDGCode[fNPDGCodes++]=-3212; // 28 = Anti Sigma 0
847 fPDGCode[fNPDGCodes++]=-3112; // 29 = Anti Sigma 0
848 fPDGCode[fNPDGCodes++]=-3322; // 30 = Anti Xi 0
849 fPDGCode[fNPDGCodes++]=-3312; // 31 = Anti Xi +
850 fPDGCode[fNPDGCodes++]=-3334; // 32 = Anti Omega +
857 /* --- Define additional particles */
858 Gspart(33, "OMEGA(782)", 3, 0.782, 0., 7.836e-23);
859 fPDGCode[fNPDGCodes++]=223; // 33 = Omega(782)
861 Gspart(34, "PHI(1020)", 3, 1.019, 0., 1.486e-22);
862 fPDGCode[fNPDGCodes++]=333; // 34 = PHI (1020)
864 Gspart(35, "D +", 4, 1.87, 1., 1.066e-12);
865 fPDGCode[fNPDGCodes++]=411; // 35 = D+
867 Gspart(36, "D -", 4, 1.87, -1., 1.066e-12);
868 fPDGCode[fNPDGCodes++]=-411; // 36 = D-
870 Gspart(37, "D 0", 3, 1.865, 0., 4.2e-13);
871 fPDGCode[fNPDGCodes++]=421; // 37 = D0
873 Gspart(38, "ANTI D 0", 3, 1.865, 0., 4.2e-13);
874 fPDGCode[fNPDGCodes++]=-421; // 38 = D0 bar
876 fPDGCode[fNPDGCodes++]=-99; // 39 = unassigned
878 fPDGCode[fNPDGCodes++]=-99; // 40 = unassigned
880 fPDGCode[fNPDGCodes++]=-99; // 41 = unassigned
882 Gspart(42, "RHO +", 4, 0.768, 1., 4.353e-24);
883 fPDGCode[fNPDGCodes++]=213; // 42 = RHO+
885 Gspart(43, "RHO -", 4, 0.768, -1., 4.353e-24);
886 fPDGCode[fNPDGCodes++]=-213; // 43 = RHO-
888 Gspart(44, "RHO 0", 3, 0.768, 0., 4.353e-24);
889 fPDGCode[fNPDGCodes++]=113; // 44 = RHO0
892 // Use ENDF-6 mapping for ions = 10000*z+10*a+iso
894 // and numbers above 5 000 000 for special applications
897 const Int_t kion=10000000;
899 const Int_t kspe=50000000;
901 TDatabasePDG *pdgDB = TDatabasePDG::Instance();
903 const Double_t kAu2Gev=0.9314943228;
904 const Double_t khSlash = 1.0545726663e-27;
905 const Double_t kErg2Gev = 1/1.6021773349e-3;
906 const Double_t khShGev = khSlash*kErg2Gev;
907 const Double_t kYear2Sec = 3600*24*365.25;
910 // mass and life-time from PDG
911 pdgDB->AddParticle("B(s)*0","B(s)*0",
912 5.4163, kTRUE, 0.047, +0.,"Meson", 533);
914 pdgDB->AddParticle("B(s)*0 bar","B(s)*0 bar",
915 5.4163, kTRUE, 0.047, -0.,"Meson", -533);
919 // value for mass used by Hijing
920 pdgDB->AddParticle("Sigma(c)*+","Sigma(c)*+",
921 2.4536, kTRUE, -1., +1.,"Baryon", 4214);
923 pdgDB->AddParticle("Sigma(c)*-","Sigma(c)*-",
924 2.4536, kTRUE, -1., -1.,"Baryon", -4214);
925 // equivalent to 4312 ? Hijing uses m=2.55
926 pdgDB->AddParticle("Xsi(c)0","Xsi(c)0",
927 2.4703, kTRUE, -1., +0.,"Baryon", 4132);
929 pdgDB->AddParticle("Xsi(c)0 bar","Xsi(c)0 bar",
930 2.4703, kTRUE, -1., -0.,"Baryon", -4132);
931 // equivalent to 4322 ? Hijing uses m=2.55
932 pdgDB->AddParticle("Xi(c)+","Xi(c)+",
933 2.4656, kFALSE, -1., +1.,"Baryon", 4232);
935 pdgDB->AddParticle("Xi(c)-","Xi(c)-",
936 2.4656, kFALSE, -1., -1.,"Baryon", -4232);
937 // mass values from Hijing
939 pdgDB->AddParticle("Xsi(c)*0","Xsi(c)*0",
940 2.63, kTRUE, -1., +0.,"Baryon", 4314);
942 pdgDB->AddParticle("Xsi(c)*0 bar","Xsi(c)*0 bar",
943 2.63, kTRUE, -1., -0.,"Baryon", -4314);
945 pdgDB->AddParticle("Xsi(c)*+","Xsi(c)*+",
946 2.63, kTRUE, -1., +1.,"Baryon", 4324);
948 pdgDB->AddParticle("Xsi(c)*-","Xsi(c)*-",
949 2.63, kTRUE, -1., -1.,"Baryon", -4324);
951 // pdg mass value, Hijing uses m=2.73.
952 pdgDB->AddParticle("Omega(c)0","Omega(c)0",
953 2.7040, kFALSE, khShGev/0.064e-12, +0.,"Baryon", 4332);
955 pdgDB->AddParticle("Omega(c)0 bar","Omega(c)0 bar",
956 2.7040, kFALSE, khShGev/0.064e-12, -0.,"Baryon", -4332);
957 // mass value from Hijing
958 pdgDB->AddParticle("Omega(c)*0","Omega(c)*0",
959 2.8000, kFALSE, -1., +0.,"Baryon", 4334);
961 pdgDB->AddParticle("Omega(c)*0 bar","Omega(c)*0",
962 2.8000, kFALSE, -1., -0.,"Baryon", -4334);
966 // mass value from Hijing
967 pdgDB->AddParticle("Sigma(b)*+","Sigma(b)*+",
968 5.8100, kFALSE, -1., +1.,"Baryon", 5224);
970 pdgDB->AddParticle("Sigma(b)*-","Sigma(b)*-",
971 5.8100, kFALSE, -1., -1.,"Baryon", -5224);
974 pdgDB->AddParticle("Xi(b)0","Xi(b)0",
975 5.8400, kFALSE, -1., +0.,"Baryon", 5232);
977 pdgDB->AddParticle("Xi(b)0 bar","Xi(b)0 bar",
978 5.8100, kFALSE, -1., -0.,"Baryon", -5232);
982 pdgDB->AddParticle("Deuteron","Deuteron",2*kAu2Gev+8.071e-3,kTRUE,
983 0,1,"Ion",kion+10020);
984 fPDGCode[fNPDGCodes++]=kion+10020; // 45 = Deuteron
986 pdgDB->AddParticle("Triton","Triton",3*kAu2Gev+14.931e-3,kFALSE,
987 khShGev/(12.33*kYear2Sec),1,"Ion",kion+10030);
988 fPDGCode[fNPDGCodes++]=kion+10030; // 46 = Triton
990 pdgDB->AddParticle("Alpha","Alpha",4*kAu2Gev+2.424e-3,kTRUE,
991 khShGev/(12.33*kYear2Sec),2,"Ion",kion+20040);
992 fPDGCode[fNPDGCodes++]=kion+20040; // 47 = Alpha
994 fPDGCode[fNPDGCodes++]=0; // 48 = geantino mapped to rootino
996 pdgDB->AddParticle("HE3","HE3",3*kAu2Gev+14.931e-3,kFALSE,
997 0,2,"Ion",kion+20030);
998 fPDGCode[fNPDGCodes++]=kion+20030; // 49 = HE3
1000 pdgDB->AddParticle("Cherenkov","Cherenkov",0,kFALSE,
1001 0,0,"Special",kspe+50);
1002 fPDGCode[fNPDGCodes++]=kspe+50; // 50 = Cherenkov
1004 Gspart(51, "FeedbackPhoton", 7, 0., 0.,1.e20 );
1005 pdgDB->AddParticle("FeedbackPhoton","FeedbackPhoton",0,kFALSE,
1006 0,0,"Special",kspe+51);
1007 fPDGCode[fNPDGCodes++]=kspe+51; // 51 = FeedbackPhoton
1008 Gspart(52, "Lambda_c+", 4, 2.2849, +1., 2.06e-13);
1009 fPDGCode[fNPDGCodes++]=4122; //52 = Lambda_c+
1011 Gspart(53, "Lambda_c-", 4, 2.2849, -1., 2.06e-13);
1012 fPDGCode[fNPDGCodes++]=-4122; //53 = Lambda_c-
1014 Gspart(54, "D_s+", 4, 1.9685, +1., 4.67e-13);
1015 fPDGCode[fNPDGCodes++]=431; //54 = D_s+
1017 Gspart(55, "D_s-", 4, 1.9685, -1., 4.67e-13);
1018 fPDGCode[fNPDGCodes++]=-431; //55 = D_s-
1020 Gspart(56, "Tau+", 5, 1.77705, +1., 2.9e-13);
1021 fPDGCode[fNPDGCodes++]=15; //56 = Tau+
1023 Gspart(57, "Tau-", 5, 1.77705, -1., 2.9e-13);
1024 fPDGCode[fNPDGCodes++]=-15; //57 = Tau-
1026 Gspart(58, "B0", 3, 5.2792, +0., 1.56e-12);
1027 fPDGCode[fNPDGCodes++]=511; //58 = B0
1029 Gspart(59, "B0 bar", 3, 5.2792, -0., 1.56e-12);
1030 fPDGCode[fNPDGCodes++]=-511; //58 = B0bar
1032 Gspart(60, "B+", 4, 5.2789, +1., 1.65e-12);
1033 fPDGCode[fNPDGCodes++]=521; //60 = B+
1035 Gspart(61, "B-", 4, 5.2789, -1., 1.65e-12);
1036 fPDGCode[fNPDGCodes++]=-521; //61 = B-
1038 Gspart(62, "Bs", 3, 5.3693, +0., 1.54e-12);
1039 fPDGCode[fNPDGCodes++]=521; //62 = B_s
1041 Gspart(63, "Bs bar", 3, 5.3693, -0., 1.54e-12);
1042 fPDGCode[fNPDGCodes++]=-521; //63 = B_s bar
1044 Gspart(64, "Lambda_b", 3, 5.624, +0., 1.24e-12);
1045 fPDGCode[fNPDGCodes++]=5122; //64 = Lambda_b
1047 Gspart(65, "Lambda_b bar", 3, 5.624, -0., 1.24e-12);
1048 fPDGCode[fNPDGCodes++]=-5122; //65 = Lambda_b bar
1050 Gspart(66, "J/Psi", 3.09688, 3, 0., 0.);
1051 fPDGCode[fNPDGCodes++]=443; // 66 = J/Psi
1053 Gspart(67, "Psi Prime", 3, 3.686, 0., 0.);
1054 fPDGCode[fNPDGCodes++]=20443; // 67 = Psi prime
1056 Gspart(68, "Upsilon(1S)", 9.46037, 3, 0., 0.);
1057 fPDGCode[fNPDGCodes++]=553; // 68 = Upsilon(1S)
1059 Gspart(69, "Upsilon(2S)", 10.0233, 3, 0., 0.);
1060 fPDGCode[fNPDGCodes++]=20553; // 69 = Upsilon(2S)
1062 Gspart(70, "Upsilon(3S)", 10.3553, 3, 0., 0.);
1063 fPDGCode[fNPDGCodes++]=30553; // 70 = Upsilon(3S)
1065 /* --- Define additional decay modes --- */
1066 /* --- omega(783) --- */
1067 for (kz = 0; kz < 6; ++kz) {
1078 Gsdk(ipa, bratio, mode);
1079 /* --- phi(1020) --- */
1080 for (kz = 0; kz < 6; ++kz) {
1095 Gsdk(ipa, bratio, mode);
1098 for (kz = 0; kz < 6; ++kz) {
1111 Gsdk(ipa, bratio, mode);
1115 for (kz = 0; kz < 6; ++kz) {
1128 Gsdk(ipa, bratio, mode);
1132 for (kz = 0; kz < 6; ++kz) {
1143 Gsdk(ipa, bratio, mode);
1145 /* --- Anti D0 --- */
1147 for (kz = 0; kz < 6; ++kz) {
1158 Gsdk(ipa, bratio, mode);
1161 for (kz = 0; kz < 6; ++kz) {
1168 Gsdk(ipa, bratio, mode);
1170 for (kz = 0; kz < 6; ++kz) {
1177 Gsdk(ipa, bratio, mode);
1179 for (kz = 0; kz < 6; ++kz) {
1186 Gsdk(ipa, bratio, mode);
1189 for (kz = 0; kz < 6; ++kz) {
1198 Gsdk(ipa, bratio, mode);
1201 Gsdk(ipa, bratio, mode);
1204 Gsdk(ipa, bratio, mode);
1209 //_____________________________________________________________________________
1210 Int_t TGeant3::VolId(const Text_t *name) const
1213 // Return the unique numeric identifier for volume name
1216 strncpy((char *) &gname, name, 4);
1217 for(i=1; i<=fGcnum->nvolum; i++)
1218 if(gname == fZiq[fGclink->jvolum+i]) return i;
1219 printf("VolId: Volume %s not found\n",name);
1223 //_____________________________________________________________________________
1224 Int_t TGeant3::NofVolumes() const
1227 // Return total number of volumes in the geometry
1229 return fGcnum->nvolum;
1232 //_____________________________________________________________________________
1233 const char* TGeant3::VolName(Int_t id) const
1236 // Return the volume name given the volume identifier
1238 if(id<1 || id > fGcnum->nvolum || fGclink->jvolum<=0)
1239 return fVolNames[fGcnum->nvolum];
1241 return fVolNames[id-1];
1244 //_____________________________________________________________________________
1245 void TGeant3::SetCut(const char* cutName, Float_t cutValue)
1248 // Set transport cuts for particles
1250 if(!strcmp(cutName,"CUTGAM"))
1251 fGccuts->cutgam=cutValue;
1252 else if(!strcmp(cutName,"CUTGAM"))
1253 fGccuts->cutele=cutValue;
1254 else if(!strcmp(cutName,"CUTELE"))
1255 fGccuts->cutneu=cutValue;
1256 else if(!strcmp(cutName,"CUTHAD"))
1257 fGccuts->cuthad=cutValue;
1258 else if(!strcmp(cutName,"CUTMUO"))
1259 fGccuts->cutmuo=cutValue;
1260 else if(!strcmp(cutName,"BCUTE"))
1261 fGccuts->bcute=cutValue;
1262 else if(!strcmp(cutName,"BCUTM"))
1263 fGccuts->bcutm=cutValue;
1264 else if(!strcmp(cutName,"DCUTE"))
1265 fGccuts->dcute=cutValue;
1266 else if(!strcmp(cutName,"DCUTM"))
1267 fGccuts->dcutm=cutValue;
1268 else if(!strcmp(cutName,"PPCUTM"))
1269 fGccuts->ppcutm=cutValue;
1270 else if(!strcmp(cutName,"TOFMAX"))
1271 fGccuts->tofmax=cutValue;
1272 else Warning("SetCut","Cut %s not implemented\n",cutName);
1275 //_____________________________________________________________________________
1276 void TGeant3::SetProcess(const char* flagName, Int_t flagValue)
1279 // Set thresholds for different processes
1281 if(!strcmp(flagName,"PAIR"))
1282 fGcphys->ipair=flagValue;
1283 else if(!strcmp(flagName,"COMP"))
1284 fGcphys->icomp=flagValue;
1285 else if(!strcmp(flagName,"PHOT"))
1286 fGcphys->iphot=flagValue;
1287 else if(!strcmp(flagName,"PFIS"))
1288 fGcphys->ipfis=flagValue;
1289 else if(!strcmp(flagName,"DRAY"))
1290 fGcphys->idray=flagValue;
1291 else if(!strcmp(flagName,"ANNI"))
1292 fGcphys->ianni=flagValue;
1293 else if(!strcmp(flagName,"BREM"))
1294 fGcphys->ibrem=flagValue;
1295 else if(!strcmp(flagName,"HADR"))
1296 fGcphys->ihadr=flagValue;
1297 else if(!strcmp(flagName,"MUNU"))
1298 fGcphys->imunu=flagValue;
1299 else if(!strcmp(flagName,"DCAY"))
1300 fGcphys->idcay=flagValue;
1301 else if(!strcmp(flagName,"LOSS"))
1302 fGcphys->iloss=flagValue;
1303 else if(!strcmp(flagName,"MULS"))
1304 fGcphys->imuls=flagValue;
1305 else if(!strcmp(flagName,"RAYL"))
1306 fGcphys->irayl=flagValue;
1307 else if(!strcmp(flagName,"STRA"))
1308 fGcphlt->istra=flagValue;
1309 else if(!strcmp(flagName,"SYNC"))
1310 fGcphlt->isync=flagValue;
1311 else Warning("SetFlag","Flag %s not implemented\n",flagName);
1314 //_____________________________________________________________________________
1315 Float_t TGeant3::Xsec(char* reac, Float_t /* energy */,
1316 Int_t part, Int_t /* mate */)
1319 // Calculate X-sections -- dummy for the moment
1321 if(!strcmp(reac,"PHOT"))
1324 Error("Xsec","Can calculate photoelectric only for photons\n");
1330 //_____________________________________________________________________________
1331 void TGeant3::TrackPosition(TLorentzVector &xyz) const
1334 // Return the current position in the master reference frame of the
1335 // track being transported
1337 xyz[0]=fGctrak->vect[0];
1338 xyz[1]=fGctrak->vect[1];
1339 xyz[2]=fGctrak->vect[2];
1340 xyz[3]=fGctrak->tofg;
1343 //_____________________________________________________________________________
1344 Float_t TGeant3::TrackTime() const
1347 // Return the current time of flight of the track being transported
1349 return fGctrak->tofg;
1352 //_____________________________________________________________________________
1353 void TGeant3::TrackMomentum(TLorentzVector &xyz) const
1356 // Return the direction and the momentum (GeV/c) of the track
1357 // currently being transported
1359 Double_t ptot=fGctrak->vect[6];
1360 xyz[0]=fGctrak->vect[3]*ptot;
1361 xyz[1]=fGctrak->vect[4]*ptot;
1362 xyz[2]=fGctrak->vect[5]*ptot;
1363 xyz[3]=fGctrak->getot;
1366 //_____________________________________________________________________________
1367 Float_t TGeant3::TrackCharge() const
1370 // Return charge of the track currently transported
1372 return fGckine->charge;
1375 //_____________________________________________________________________________
1376 Float_t TGeant3::TrackMass() const
1379 // Return the mass of the track currently transported
1381 return fGckine->amass;
1384 //_____________________________________________________________________________
1385 Int_t TGeant3::TrackPid() const
1388 // Return the id of the particle transported
1390 return PDGFromId(fGckine->ipart);
1393 //_____________________________________________________________________________
1394 Float_t TGeant3::TrackStep() const
1397 // Return the length in centimeters of the current step
1399 return fGctrak->step;
1402 //_____________________________________________________________________________
1403 Float_t TGeant3::TrackLength() const
1406 // Return the length of the current track from its origin
1408 return fGctrak->sleng;
1411 //_____________________________________________________________________________
1412 Bool_t TGeant3::IsNewTrack() const
1415 // True if the track is not at the boundary of the current volume
1417 return (fGctrak->sleng==0);
1420 //_____________________________________________________________________________
1421 Bool_t TGeant3::IsTrackInside() const
1424 // True if the track is not at the boundary of the current volume
1426 return (fGctrak->inwvol==0);
1429 //_____________________________________________________________________________
1430 Bool_t TGeant3::IsTrackEntering() const
1433 // True if this is the first step of the track in the current volume
1435 return (fGctrak->inwvol==1);
1438 //_____________________________________________________________________________
1439 Bool_t TGeant3::IsTrackExiting() const
1442 // True if this is the last step of the track in the current volume
1444 return (fGctrak->inwvol==2);
1447 //_____________________________________________________________________________
1448 Bool_t TGeant3::IsTrackOut() const
1451 // True if the track is out of the setup
1453 return (fGctrak->inwvol==3);
1456 //_____________________________________________________________________________
1457 Bool_t TGeant3::IsTrackStop() const
1460 // True if the track energy has fallen below the threshold
1462 return (fGctrak->istop==2);
1465 //_____________________________________________________________________________
1466 Int_t TGeant3::NSecondaries() const
1469 // Number of secondary particles generated in the current step
1471 return fGcking->ngkine;
1474 //_____________________________________________________________________________
1475 Int_t TGeant3::CurrentEvent() const
1478 // Number of the current event
1480 return fGcflag->idevt;
1483 //_____________________________________________________________________________
1484 const char* TGeant3::ProdProcess() const
1487 // Name of the process that has produced the secondary particles
1488 // in the current step
1490 static char proc[5];
1491 const Int_t kIpMec[13] = { 5,6,7,8,9,10,11,12,21,23,25,105,108 };
1494 if(fGcking->ngkine>0) {
1495 for (km = 0; km < fGctrak->nmec; ++km) {
1496 for (im = 0; im < 13; ++im) {
1497 if (fGctrak->lmec[km] == kIpMec[im]) {
1498 mec = fGctrak->lmec[km];
1499 if (0 < mec && mec < 31) {
1500 strncpy(proc,(char *)&fGctrak->namec[mec - 1],4);
1501 } else if (mec - 100 <= 30 && mec - 100 > 0) {
1502 strncpy(proc,(char *)&fGctpol->namec1[mec - 101],4);
1509 strcpy(proc,"UNKN");
1510 } else strcpy(proc,"NONE");
1514 //_____________________________________________________________________________
1515 void TGeant3::GetSecondary(Int_t isec, Int_t& ipart,
1516 TLorentzVector &x, TLorentzVector &p)
1519 // Get the parameters of the secondary track number isec produced
1520 // in the current step
1523 if(-1<isec && isec<fGcking->ngkine) {
1524 ipart=Int_t (fGcking->gkin[isec][4] +0.5);
1526 x[i]=fGckin3->gpos[isec][i];
1527 p[i]=fGcking->gkin[isec][i];
1529 x[3]=fGcking->tofd[isec];
1530 p[3]=fGcking->gkin[isec][3];
1532 printf(" * TGeant3::GetSecondary * Secondary %d does not exist\n",isec);
1533 x[0]=x[1]=x[2]=x[3]=p[0]=p[1]=p[2]=p[3]=0;
1538 //_____________________________________________________________________________
1539 void TGeant3::InitLego()
1542 // Set switches for lego transport
1545 SetDEBU(0,0,0); //do not print a message
1548 //_____________________________________________________________________________
1549 Bool_t TGeant3::IsTrackDisappeared() const
1552 // True if the current particle has disappered
1553 // either because it decayed or because it underwent
1554 // an inelastic collision
1556 return (fGctrak->istop==1);
1559 //_____________________________________________________________________________
1560 Bool_t TGeant3::IsTrackAlive() const
1563 // True if the current particle is alive and will continue to be
1566 return (fGctrak->istop==0);
1569 //_____________________________________________________________________________
1570 void TGeant3::StopTrack()
1573 // Stop the transport of the current particle and skip to the next
1578 //_____________________________________________________________________________
1579 void TGeant3::StopEvent()
1582 // Stop simulation of the current event and skip to the next
1587 //_____________________________________________________________________________
1588 Float_t TGeant3::MaxStep() const
1591 // Return the maximum step length in the current medium
1593 return fGctmed->stemax;
1596 //_____________________________________________________________________________
1597 void TGeant3::SetMaxStep(Float_t maxstep)
1600 // Set the maximum step allowed till the particle is in the current medium
1602 fGctmed->stemax=maxstep;
1605 //_____________________________________________________________________________
1606 void TGeant3::SetMaxNStep(Int_t maxnstp)
1609 // Set the maximum number of steps till the particle is in the current medium
1611 fGctrak->maxnst=maxnstp;
1614 //_____________________________________________________________________________
1615 Int_t TGeant3::GetMaxNStep() const
1618 // Maximum number of steps allowed in current medium
1620 return fGctrak->maxnst;
1623 //_____________________________________________________________________________
1624 void TGeant3::Material(Int_t& kmat, const char* name, Float_t a, Float_t z,
1625 Float_t dens, Float_t radl, Float_t absl, Float_t* buf,
1629 // Defines a Material
1631 // kmat number assigned to the material
1632 // name material name
1633 // a atomic mass in au
1635 // dens density in g/cm3
1636 // absl absorbtion length in cm
1637 // if >=0 it is ignored and the program
1638 // calculates it, if <0. -absl is taken
1639 // radl radiation length in cm
1640 // if >=0 it is ignored and the program
1641 // calculates it, if <0. -radl is taken
1642 // buf pointer to an array of user words
1643 // nbuf number of user words
1645 Int_t jmate=fGclink->jmate;
1651 for(i=1; i<=ns; i++) {
1652 if(fZlq[jmate-i]==0) {
1658 gsmate(kmat,PASSCHARD(name), a, z, dens, radl, absl, buf,
1659 nwbuf PASSCHARL(name));
1662 //_____________________________________________________________________________
1663 void TGeant3::Mixture(Int_t& kmat, const char* name, Float_t* a, Float_t* z,
1664 Float_t dens, Int_t nlmat, Float_t* wmat)
1667 // Defines mixture OR COMPOUND IMAT as composed by
1668 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
1670 // If NLMAT > 0 then wmat contains the proportion by
1671 // weights of each basic material in the mixture.
1673 // If nlmat < 0 then WMAT contains the number of atoms
1674 // of a given kind into the molecule of the COMPOUND
1675 // In this case, WMAT in output is changed to relative
1678 Int_t jmate=fGclink->jmate;
1684 for(i=1; i<=ns; i++) {
1685 if(fZlq[jmate-i]==0) {
1691 gsmixt(kmat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
1694 //_____________________________________________________________________________
1695 void TGeant3::Medium(Int_t& kmed, const char* name, Int_t nmat, Int_t isvol,
1696 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
1697 Float_t stemax, Float_t deemax, Float_t epsil,
1698 Float_t stmin, Float_t* ubuf, Int_t nbuf)
1701 // kmed tracking medium number assigned
1702 // name tracking medium name
1703 // nmat material number
1704 // isvol sensitive volume flag
1705 // ifield magnetic field
1706 // fieldm max. field value (kilogauss)
1707 // tmaxfd max. angle due to field (deg/step)
1708 // stemax max. step allowed
1709 // deemax max. fraction of energy lost in a step
1710 // epsil tracking precision (cm)
1711 // stmin min. step due to continuos processes (cm)
1713 // ifield = 0 if no magnetic field; ifield = -1 if user decision in guswim;
1714 // ifield = 1 if tracking performed with grkuta; ifield = 2 if tracking
1715 // performed with ghelix; ifield = 3 if tracking performed with ghelx3.
1717 Int_t jtmed=fGclink->jtmed;
1723 for(i=1; i<=ns; i++) {
1724 if(fZlq[jtmed-i]==0) {
1730 gstmed(kmed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1731 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1734 //_____________________________________________________________________________
1735 void TGeant3::Matrix(Int_t& krot, Float_t thex, Float_t phix, Float_t they,
1736 Float_t phiy, Float_t thez, Float_t phiz)
1739 // krot rotation matrix number assigned
1740 // theta1 polar angle for axis i
1741 // phi1 azimuthal angle for axis i
1742 // theta2 polar angle for axis ii
1743 // phi2 azimuthal angle for axis ii
1744 // theta3 polar angle for axis iii
1745 // phi3 azimuthal angle for axis iii
1747 // it defines the rotation matrix number irot.
1749 Int_t jrotm=fGclink->jrotm;
1755 for(i=1; i<=ns; i++) {
1756 if(fZlq[jrotm-i]==0) {
1762 gsrotm(krot, thex, phix, they, phiy, thez, phiz);
1765 //_____________________________________________________________________________
1766 Int_t TGeant3::GetMedium() const
1769 // Return the number of the current medium
1771 return fGctmed->numed;
1774 //_____________________________________________________________________________
1775 Float_t TGeant3::Edep() const
1778 // Return the energy lost in the current step
1780 return fGctrak->destep;
1783 //_____________________________________________________________________________
1784 Float_t TGeant3::Etot() const
1787 // Return the total energy of the current track
1789 return fGctrak->getot;
1792 //_____________________________________________________________________________
1793 void TGeant3::Rndm(Float_t* r, const Int_t n) const
1796 // Return an array of n random numbers uniformly distributed
1797 // between 0 and 1 not included
1802 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1804 // Functions from GBASE
1806 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1808 //____________________________________________________________________________
1809 void TGeant3::Gfile(const char *filename, const char *option)
1812 // Routine to open a GEANT/RZ data base.
1814 // LUN logical unit number associated to the file
1816 // CHFILE RZ file name
1818 // CHOPT is a character string which may be
1819 // N To create a new file
1820 // U to open an existing file for update
1821 // " " to open an existing file for read only
1822 // Q The initial allocation (default 1000 records)
1823 // is given in IQUEST(10)
1824 // X Open the file in exchange format
1825 // I Read all data structures from file to memory
1826 // O Write all data structures from memory to file
1829 // If options "I" or "O" all data structures are read or
1830 // written from/to file and the file is closed.
1831 // See routine GRMDIR to create subdirectories
1832 // See routines GROUT,GRIN to write,read objects
1834 grfile(21, PASSCHARD(filename), PASSCHARD(option) PASSCHARL(filename)
1838 //____________________________________________________________________________
1839 void TGeant3::Gpcxyz()
1842 // Print track and volume parameters at current point
1847 //_____________________________________________________________________________
1848 void TGeant3::Ggclos()
1851 // Closes off the geometry setting.
1852 // Initializes the search list for the contents of each
1853 // volume following the order they have been positioned, and
1854 // inserting the content '0' when a call to GSNEXT (-1) has
1855 // been required by the user.
1856 // Performs the development of the JVOLUM structure for all
1857 // volumes with variable parameters, by calling GGDVLP.
1858 // Interprets the user calls to GSORD, through GGORD.
1859 // Computes and stores in a bank (next to JVOLUM mother bank)
1860 // the number of levels in the geometrical tree and the
1861 // maximum number of contents per level, by calling GGNLEV.
1862 // Sets status bit for CONCAVE volumes, through GGCAVE.
1863 // Completes the JSET structure with the list of volume names
1864 // which identify uniquely a given physical detector, the
1865 // list of bit numbers to pack the corresponding volume copy
1866 // numbers, and the generic path(s) in the JVOLUM tree,
1867 // through the routine GHCLOS.
1870 // Create internal list of volumes
1871 fVolNames = new char[fGcnum->nvolum+1][5];
1873 for(i=0; i<fGcnum->nvolum; ++i) {
1874 strncpy(fVolNames[i], (char *) &fZiq[fGclink->jvolum+i+1], 4);
1875 fVolNames[i][4]='\0';
1877 strcpy(fVolNames[fGcnum->nvolum],"NULL");
1880 //_____________________________________________________________________________
1881 void TGeant3::Glast()
1884 // Finish a Geant run
1889 //_____________________________________________________________________________
1890 void TGeant3::Gprint(const char *name)
1893 // Routine to print data structures
1894 // CHNAME name of a data structure
1898 gprint(PASSCHARD(vname),0 PASSCHARL(vname));
1901 //_____________________________________________________________________________
1902 void TGeant3::Grun()
1905 // Steering function to process one run
1910 //_____________________________________________________________________________
1911 void TGeant3::Gtrig()
1914 // Steering function to process one event
1919 //_____________________________________________________________________________
1920 void TGeant3::Gtrigc()
1923 // Clear event partition
1928 //_____________________________________________________________________________
1929 void TGeant3::Gtrigi()
1932 // Initialises event partition
1937 //_____________________________________________________________________________
1938 void TGeant3::Gwork(Int_t nwork)
1941 // Allocates workspace in ZEBRA memory
1946 //_____________________________________________________________________________
1947 void TGeant3::Gzinit()
1950 // To initialise GEANT/ZEBRA data structures
1955 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1957 // Functions from GCONS
1959 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1961 //_____________________________________________________________________________
1962 void TGeant3::Gfmate(Int_t imat, char *name, Float_t &a, Float_t &z,
1963 Float_t &dens, Float_t &radl, Float_t &absl,
1964 Float_t* ubuf, Int_t& nbuf)
1967 // Return parameters for material IMAT
1969 gfmate(imat, PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
1973 //_____________________________________________________________________________
1974 void TGeant3::Gfpart(Int_t ipart, char *name, Int_t &itrtyp,
1975 Float_t &amass, Float_t &charge, Float_t &tlife)
1978 // Return parameters for particle of type IPART
1982 Int_t igpart = IdFromPDG(ipart);
1983 gfpart(igpart, PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
1987 //_____________________________________________________________________________
1988 void TGeant3::Gftmed(Int_t numed, char *name, Int_t &nmat, Int_t &isvol,
1989 Int_t &ifield, Float_t &fieldm, Float_t &tmaxfd,
1990 Float_t &stemax, Float_t &deemax, Float_t &epsil,
1991 Float_t &stmin, Float_t *ubuf, Int_t *nbuf)
1994 // Return parameters for tracking medium NUMED
1996 gftmed(numed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1997 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
2001 void TGeant3::Gftmat(Int_t imate, Int_t ipart, char *chmeca, Int_t kdim,
2002 Float_t* tkin, Float_t* value, Float_t* pcut,
2006 // Return parameters for tracking medium NUMED
2008 gftmat(imate, ipart, PASSCHARD(chmeca), kdim,
2009 tkin, value, pcut, ixst PASSCHARL(chmeca));
2013 //_____________________________________________________________________________
2014 Float_t TGeant3::Gbrelm(Float_t z, Float_t t, Float_t bcut)
2017 // To calculate energy loss due to soft muon BREMSSTRAHLUNG
2019 return gbrelm(z,t,bcut);
2022 //_____________________________________________________________________________
2023 Float_t TGeant3::Gprelm(Float_t z, Float_t t, Float_t bcut)
2026 // To calculate DE/DX in GeV*barn/atom for direct pair production by muons
2028 return gprelm(z,t,bcut);
2031 //_____________________________________________________________________________
2032 void TGeant3::Gmate()
2035 // Define standard GEANT materials
2040 //_____________________________________________________________________________
2041 void TGeant3::Gpart()
2044 // Define standard GEANT particles plus selected decay modes
2045 // and branching ratios.
2050 //_____________________________________________________________________________
2051 void TGeant3::Gsdk(Int_t ipart, Float_t *bratio, Int_t *mode)
2053 // Defines branching ratios and decay modes for standard
2055 gsdk(ipart,bratio,mode);
2058 //_____________________________________________________________________________
2059 void TGeant3::Gsmate(Int_t imat, const char *name, Float_t a, Float_t z,
2060 Float_t dens, Float_t radl, Float_t absl)
2063 // Defines a Material
2065 // kmat number assigned to the material
2066 // name material name
2067 // a atomic mass in au
2069 // dens density in g/cm3
2070 // absl absorbtion length in cm
2071 // if >=0 it is ignored and the program
2072 // calculates it, if <0. -absl is taken
2073 // radl radiation length in cm
2074 // if >=0 it is ignored and the program
2075 // calculates it, if <0. -radl is taken
2076 // buf pointer to an array of user words
2077 // nbuf number of user words
2081 gsmate(imat,PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
2085 //_____________________________________________________________________________
2086 void TGeant3::Gsmixt(Int_t imat, const char *name, Float_t *a, Float_t *z,
2087 Float_t dens, Int_t nlmat, Float_t *wmat)
2090 // Defines mixture OR COMPOUND IMAT as composed by
2091 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
2093 // If NLMAT.GT.0 then WMAT contains the PROPORTION BY
2094 // WEIGTHS OF EACH BASIC MATERIAL IN THE MIXTURE.
2096 // If NLMAT.LT.0 then WMAT contains the number of atoms
2097 // of a given kind into the molecule of the COMPOUND
2098 // In this case, WMAT in output is changed to relative
2101 gsmixt(imat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
2104 //_____________________________________________________________________________
2105 void TGeant3::Gspart(Int_t ipart, const char *name, Int_t itrtyp,
2106 Float_t amass, Float_t charge, Float_t tlife)
2109 // Store particle parameters
2111 // ipart particle code
2112 // name particle name
2113 // itrtyp transport method (see GEANT manual)
2114 // amass mass in GeV/c2
2115 // charge charge in electron units
2116 // tlife lifetime in seconds
2120 gspart(ipart,PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
2124 //_____________________________________________________________________________
2125 void TGeant3::Gstmed(Int_t numed, const char *name, Int_t nmat, Int_t isvol,
2126 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
2127 Float_t stemax, Float_t deemax, Float_t epsil,
2131 // NTMED Tracking medium number
2132 // NAME Tracking medium name
2133 // NMAT Material number
2134 // ISVOL Sensitive volume flag
2135 // IFIELD Magnetic field
2136 // FIELDM Max. field value (Kilogauss)
2137 // TMAXFD Max. angle due to field (deg/step)
2138 // STEMAX Max. step allowed
2139 // DEEMAX Max. fraction of energy lost in a step
2140 // EPSIL Tracking precision (cm)
2141 // STMIN Min. step due to continuos processes (cm)
2143 // IFIELD = 0 if no magnetic field; IFIELD = -1 if user decision in GUSWIM;
2144 // IFIELD = 1 if tracking performed with GRKUTA; IFIELD = 2 if tracking
2145 // performed with GHELIX; IFIELD = 3 if tracking performed with GHELX3.
2149 gstmed(numed,PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
2150 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
2153 //_____________________________________________________________________________
2154 void TGeant3::Gsckov(Int_t itmed, Int_t npckov, Float_t *ppckov,
2155 Float_t *absco, Float_t *effic, Float_t *rindex)
2158 // Stores the tables for UV photon tracking in medium ITMED
2159 // Please note that it is the user's responsability to
2160 // provide all the coefficients:
2163 // ITMED Tracking medium number
2164 // NPCKOV Number of bins of each table
2165 // PPCKOV Value of photon momentum (in GeV)
2166 // ABSCO Absorbtion coefficients
2167 // dielectric: absorbtion length in cm
2168 // metals : absorbtion fraction (0<=x<=1)
2169 // EFFIC Detection efficiency for UV photons
2170 // RINDEX Refraction index (if=0 metal)
2172 gsckov(itmed,npckov,ppckov,absco,effic,rindex);
2175 //_____________________________________________________________________________
2176 void TGeant3::Gstpar(Int_t itmed, const char *param, Float_t parval)
2179 // To change the value of cut or mechanism "CHPAR"
2180 // to a new value PARVAL for tracking medium ITMED
2181 // The data structure JTMED contains the standard tracking
2182 // parameters (CUTS and flags to control the physics processes) which
2183 // are used by default for all tracking media. It is possible to
2184 // redefine individually with GSTPAR any of these parameters for a
2185 // given tracking medium.
2186 // ITMED tracking medium number
2187 // CHPAR is a character string (variable name)
2188 // PARVAL must be given as a floating point.
2190 gstpar(itmed,PASSCHARD(param), parval PASSCHARL(param));
2193 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2195 // Functions from GCONS
2197 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2199 //_____________________________________________________________________________
2200 void TGeant3::Gfkine(Int_t itra, Float_t *vert, Float_t *pvert, Int_t &ipart,
2203 // Storing/Retrieving Vertex and Track parameters
2204 // ----------------------------------------------
2206 // Stores vertex parameters.
2207 // VERT array of (x,y,z) position of the vertex
2208 // NTBEAM beam track number origin of the vertex
2209 // =0 if none exists
2210 // NTTARG target track number origin of the vertex
2211 // UBUF user array of NUBUF floating point numbers
2213 // NVTX new vertex number (=0 in case of error).
2214 // Prints vertex parameters.
2215 // IVTX for vertex IVTX.
2216 // (For all vertices if IVTX=0)
2217 // Stores long life track parameters.
2218 // PLAB components of momentum
2219 // IPART type of particle (see GSPART)
2220 // NV vertex number origin of track
2221 // UBUF array of NUBUF floating point user parameters
2223 // NT track number (if=0 error).
2224 // Retrieves long life track parameters.
2225 // ITRA track number for which parameters are requested
2226 // VERT vector origin of the track
2227 // PVERT 4 momentum components at the track origin
2228 // IPART particle type (=0 if track ITRA does not exist)
2229 // NVERT vertex number origin of the track
2230 // UBUF user words stored in GSKINE.
2231 // Prints initial track parameters.
2232 // ITRA for track ITRA
2233 // (For all tracks if ITRA=0)
2237 gfkine(itra,vert,pvert,ipart,nvert,ubuf,nbuf);
2240 //_____________________________________________________________________________
2241 void TGeant3::Gfvert(Int_t nvtx, Float_t *v, Int_t &ntbeam, Int_t &nttarg,
2245 // Retrieves the parameter of a vertex bank
2246 // Vertex is generated from tracks NTBEAM NTTARG
2247 // NVTX is the new vertex number
2251 gfvert(nvtx,v,ntbeam,nttarg,tofg,ubuf,nbuf);
2254 //_____________________________________________________________________________
2255 Int_t TGeant3::Gskine(Float_t *plab, Int_t ipart, Int_t nv, Float_t *buf,
2259 // Store kinematics of track NT into data structure
2260 // Track is coming from vertex NV
2263 gskine(plab, ipart, nv, buf, nwbuf, nt);
2267 //_____________________________________________________________________________
2268 Int_t TGeant3::Gsvert(Float_t *v, Int_t ntbeam, Int_t nttarg, Float_t *ubuf,
2272 // Creates a new vertex bank
2273 // Vertex is generated from tracks NTBEAM NTTARG
2274 // NVTX is the new vertex number
2277 gsvert(v, ntbeam, nttarg, ubuf, nwbuf, nwtx);
2281 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2283 // Functions from GPHYS
2285 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2287 //_____________________________________________________________________________
2288 void TGeant3::Gphysi()
2291 // Initialise material constants for all the physics
2292 // mechanisms used by GEANT
2297 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2299 // Functions from GTRAK
2301 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2303 //_____________________________________________________________________________
2304 void TGeant3::Gdebug()
2307 // Debug the current step
2312 //_____________________________________________________________________________
2313 void TGeant3::Gekbin()
2316 // To find bin number in kinetic energy table
2317 // stored in ELOW(NEKBIN)
2322 //_____________________________________________________________________________
2323 void TGeant3::Gfinds()
2326 // Returns the set/volume parameters corresponding to
2327 // the current space point in /GCTRAK/
2328 // and fill common /GCSETS/
2330 // IHSET user set identifier
2331 // IHDET user detector identifier
2332 // ISET set number in JSET
2333 // IDET detector number in JS=LQ(JSET-ISET)
2334 // IDTYPE detector type (1,2)
2335 // NUMBV detector volume numbers (array of length NVNAME)
2336 // NVNAME number of volume levels
2341 //_____________________________________________________________________________
2342 void TGeant3::Gsking(Int_t igk)
2345 // Stores in stack JSTAK either the IGKth track of /GCKING/,
2346 // or the NGKINE tracks when IGK is 0.
2351 //_____________________________________________________________________________
2352 void TGeant3::Gskpho(Int_t igk)
2355 // Stores in stack JSTAK either the IGKth Cherenkov photon of
2356 // /GCKIN2/, or the NPHOT tracks when IGK is 0.
2361 //_____________________________________________________________________________
2362 void TGeant3::Gsstak(Int_t iflag)
2365 // Stores in auxiliary stack JSTAK the particle currently
2366 // described in common /GCKINE/.
2368 // On request, creates also an entry in structure JKINE :
2370 // 0 : No entry in JKINE structure required (user)
2371 // 1 : New entry in JVERTX / JKINE structures required (user)
2372 // <0 : New entry in JKINE structure at vertex -IFLAG (user)
2373 // 2 : Entry in JKINE structure exists already (from GTREVE)
2378 //_____________________________________________________________________________
2379 void TGeant3::Gsxyz()
2382 // Store space point VECT in banks JXYZ
2387 //_____________________________________________________________________________
2388 void TGeant3::Gtrack()
2391 // Controls tracking of current particle
2396 //_____________________________________________________________________________
2397 void TGeant3::Gtreve()
2400 // Controls tracking of all particles belonging to the current event
2405 //_____________________________________________________________________________
2406 void TGeant3::GtreveRoot()
2409 // Controls tracking of all particles belonging to the current event
2414 //_____________________________________________________________________________
2415 void TGeant3::Grndm(Float_t *rvec, const Int_t len) const
2418 // To generate a vector RVECV of LEN random numbers
2419 // Copy of the CERN Library routine RANECU
2423 //_____________________________________________________________________________
2424 void TGeant3::Grndmq(Int_t &is1, Int_t &is2, const Int_t iseq,
2425 const Text_t *chopt)
2428 // To set/retrieve the seed of the random number generator
2430 grndmq(is1,is2,iseq,PASSCHARD(chopt) PASSCHARL(chopt));
2433 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2435 // Functions from GDRAW
2437 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2439 //_____________________________________________________________________________
2440 void TGeant3::Gdxyz(Int_t it)
2443 // Draw the points stored with Gsxyz relative to track it
2448 //_____________________________________________________________________________
2449 void TGeant3::Gdcxyz()
2452 // Draw the position of the current track
2457 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2459 // Functions from GGEOM
2461 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2463 //_____________________________________________________________________________
2464 void TGeant3::Gdtom(Float_t *xd, Float_t *xm, Int_t iflag)
2467 // Computes coordinates XM (Master Reference System
2468 // knowing the coordinates XD (Detector Ref System)
2469 // The local reference system can be initialized by
2470 // - the tracking routines and GDTOM used in GUSTEP
2471 // - a call to GSCMED(NLEVEL,NAMES,NUMBER)
2472 // (inverse routine is GMTOD)
2474 // If IFLAG=1 convert coordinates
2475 // IFLAG=2 convert direction cosinus
2477 gdtom(xd, xm, iflag);
2480 //_____________________________________________________________________________
2481 void TGeant3::Glmoth(const char* iudet, Int_t iunum, Int_t &nlev, Int_t *lvols,
2485 // Loads the top part of the Volume tree in LVOLS (IVO's),
2486 // LINDX (IN indices) for a given volume defined through
2487 // its name IUDET and number IUNUM.
2489 // The routine stores only upto the last level where JVOLUM
2490 // data structure is developed. If there is no development
2491 // above the current level, it returns NLEV zero.
2493 glmoth(PASSCHARD(iudet), iunum, nlev, lvols, lindx, idum PASSCHARL(iudet));
2496 //_____________________________________________________________________________
2497 void TGeant3::Gmedia(Float_t *x, Int_t &numed)
2500 // Finds in which volume/medium the point X is, and updates the
2501 // common /GCVOLU/ and the structure JGPAR accordingly.
2503 // NUMED returns the tracking medium number, or 0 if point is
2504 // outside the experimental setup.
2509 //_____________________________________________________________________________
2510 void TGeant3::Gmtod(Float_t *xm, Float_t *xd, Int_t iflag)
2513 // Computes coordinates XD (in DRS)
2514 // from known coordinates XM in MRS
2515 // The local reference system can be initialized by
2516 // - the tracking routines and GMTOD used in GUSTEP
2517 // - a call to GMEDIA(XM,NUMED)
2518 // - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER)
2519 // (inverse routine is GDTOM)
2521 // If IFLAG=1 convert coordinates
2522 // IFLAG=2 convert direction cosinus
2524 gmtod(xm, xd, iflag);
2527 //_____________________________________________________________________________
2528 void TGeant3::Gsdvn(const char *name, const char *mother, Int_t ndiv,
2532 // Create a new volume by dividing an existing one
2535 // MOTHER Mother volume name
2536 // NDIV Number of divisions
2539 // X,Y,Z of CAXIS will be translated to 1,2,3 for IAXIS.
2540 // It divides a previously defined volume.
2545 Vname(mother,vmother);
2546 gsdvn(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis PASSCHARL(vname)
2547 PASSCHARL(vmother));
2550 //_____________________________________________________________________________
2551 void TGeant3::Gsdvn2(const char *name, const char *mother, Int_t ndiv,
2552 Int_t iaxis, Float_t c0i, Int_t numed)
2555 // Create a new volume by dividing an existing one
2557 // Divides mother into ndiv divisions called name
2558 // along axis iaxis starting at coordinate value c0.
2559 // the new volume created will be medium number numed.
2564 Vname(mother,vmother);
2565 gsdvn2(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis, c0i, numed
2566 PASSCHARL(vname) PASSCHARL(vmother));
2569 //_____________________________________________________________________________
2570 void TGeant3::Gsdvs(const char *name, const char *mother, Float_t step,
2571 Int_t iaxis, Int_t numed)
2574 // Create a new volume by dividing an existing one
2579 Vname(mother,vmother);
2580 gsdvs(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed
2581 PASSCHARL(vname) PASSCHARL(vmother));
2584 //_____________________________________________________________________________
2585 void TGeant3::Gsdvs2(const char *name, const char *mother, Float_t step,
2586 Int_t iaxis, Float_t c0, Int_t numed)
2589 // Create a new volume by dividing an existing one
2594 Vname(mother,vmother);
2595 gsdvs2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0, numed
2596 PASSCHARL(vname) PASSCHARL(vmother));
2599 //_____________________________________________________________________________
2600 void TGeant3::Gsdvt(const char *name, const char *mother, Float_t step,
2601 Int_t iaxis, Int_t numed, Int_t ndvmx)
2604 // Create a new volume by dividing an existing one
2606 // Divides MOTHER into divisions called NAME along
2607 // axis IAXIS in steps of STEP. If not exactly divisible
2608 // will make as many as possible and will centre them
2609 // with respect to the mother. Divisions will have medium
2610 // number NUMED. If NUMED is 0, NUMED of MOTHER is taken.
2611 // NDVMX is the expected maximum number of divisions
2612 // (If 0, no protection tests are performed)
2617 Vname(mother,vmother);
2618 gsdvt(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed, ndvmx
2619 PASSCHARL(vname) PASSCHARL(vmother));
2622 //_____________________________________________________________________________
2623 void TGeant3::Gsdvt2(const char *name, const char *mother, Float_t step,
2624 Int_t iaxis, Float_t c0, Int_t numed, Int_t ndvmx)
2627 // Create a new volume by dividing an existing one
2629 // Divides MOTHER into divisions called NAME along
2630 // axis IAXIS starting at coordinate value C0 with step
2632 // The new volume created will have medium number NUMED.
2633 // If NUMED is 0, NUMED of mother is taken.
2634 // NDVMX is the expected maximum number of divisions
2635 // (If 0, no protection tests are performed)
2640 Vname(mother,vmother);
2641 gsdvt2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0,
2642 numed, ndvmx PASSCHARL(vname) PASSCHARL(vmother));
2645 //_____________________________________________________________________________
2646 void TGeant3::Gsord(const char *name, Int_t iax)
2649 // Flags volume CHNAME whose contents will have to be ordered
2650 // along axis IAX, by setting the search flag to -IAX
2654 // IAX = 4 Rxy (static ordering only -> GTMEDI)
2655 // IAX = 14 Rxy (also dynamic ordering -> GTNEXT)
2656 // IAX = 5 Rxyz (static ordering only -> GTMEDI)
2657 // IAX = 15 Rxyz (also dynamic ordering -> GTNEXT)
2658 // IAX = 6 PHI (PHI=0 => X axis)
2659 // IAX = 7 THETA (THETA=0 => Z axis)
2663 gsord(PASSCHARD(vname), iax PASSCHARL(vname));
2666 //_____________________________________________________________________________
2667 void TGeant3::Gspos(const char *name, Int_t nr, const char *mother, Float_t x,
2668 Float_t y, Float_t z, Int_t irot, const char *konly)
2671 // Position a volume into an existing one
2674 // NUMBER Copy number of the volume
2675 // MOTHER Mother volume name
2676 // X X coord. of the volume in mother ref. sys.
2677 // Y Y coord. of the volume in mother ref. sys.
2678 // Z Z coord. of the volume in mother ref. sys.
2679 // IROT Rotation matrix number w.r.t. mother ref. sys.
2680 // ONLY ONLY/MANY flag
2682 // It positions a previously defined volume in the mother.
2687 Vname(mother,vmother);
2688 gspos(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2689 PASSCHARD(konly) PASSCHARL(vname) PASSCHARL(vmother)
2693 //_____________________________________________________________________________
2694 void TGeant3::Gsposp(const char *name, Int_t nr, const char *mother,
2695 Float_t x, Float_t y, Float_t z, Int_t irot,
2696 const char *konly, Float_t *upar, Int_t np )
2699 // Place a copy of generic volume NAME with user number
2700 // NR inside MOTHER, with its parameters UPAR(1..NP)
2705 Vname(mother,vmother);
2706 gsposp(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2707 PASSCHARD(konly), upar, np PASSCHARL(vname) PASSCHARL(vmother)
2711 //_____________________________________________________________________________
2712 void TGeant3::Gsrotm(Int_t nmat, Float_t theta1, Float_t phi1, Float_t theta2,
2713 Float_t phi2, Float_t theta3, Float_t phi3)
2716 // nmat Rotation matrix number
2717 // THETA1 Polar angle for axis I
2718 // PHI1 Azimuthal angle for axis I
2719 // THETA2 Polar angle for axis II
2720 // PHI2 Azimuthal angle for axis II
2721 // THETA3 Polar angle for axis III
2722 // PHI3 Azimuthal angle for axis III
2724 // It defines the rotation matrix number IROT.
2726 gsrotm(nmat, theta1, phi1, theta2, phi2, theta3, phi3);
2729 //_____________________________________________________________________________
2730 void TGeant3::Gprotm(Int_t nmat)
2733 // To print rotation matrices structure JROTM
2734 // nmat Rotation matrix number
2739 //_____________________________________________________________________________
2740 Int_t TGeant3::Gsvolu(const char *name, const char *shape, Int_t nmed,
2741 Float_t *upar, Int_t npar)
2745 // SHAPE Volume type
2746 // NUMED Tracking medium number
2747 // NPAR Number of shape parameters
2748 // UPAR Vector containing shape parameters
2750 // It creates a new volume in the JVOLUM data structure.
2756 Vname(shape,vshape);
2757 gsvolu(PASSCHARD(vname), PASSCHARD(vshape), nmed, upar, npar, ivolu
2758 PASSCHARL(vname) PASSCHARL(vshape));
2762 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2764 // T H E D R A W I N G P A C K A G E
2765 // ======================================
2766 // Drawing functions. These functions allow the visualization in several ways
2767 // of the volumes defined in the geometrical data structure. It is possible
2768 // to draw the logical tree of volumes belonging to the detector (DTREE),
2769 // to show their geometrical specification (DSPEC,DFSPC), to draw them
2770 // and their cut views (DRAW, DCUT). Moreover, it is possible to execute
2771 // these commands when the hidden line removal option is activated; in
2772 // this case, the volumes can be also either translated in the space
2773 // (SHIFT), or clipped by boolean operation (CVOL). In addition, it is
2774 // possible to fill the surfaces of the volumes
2775 // with solid colours when the shading option (SHAD) is activated.
2776 // Several tools (ZOOM, LENS) have been developed to zoom detailed parts
2777 // of the detectors or to scan physical events as well.
2778 // Finally, the command MOVE will allow the rotation, translation and zooming
2779 // on real time parts of the detectors or tracks and hits of a simulated event.
2780 // Ray-tracing commands. In case the command (DOPT RAYT ON) is executed,
2781 // the drawing is performed by the Geant ray-tracing;
2782 // automatically, the color is assigned according to the tracking medium of each
2783 // volume and the volumes with a density lower/equal than the air are considered
2784 // transparent; if the option (USER) is set (ON) (again via the command (DOPT)),
2785 // the user can set color and visibility for the desired volumes via the command
2786 // (SATT), as usual, relatively to the attributes (COLO) and (SEEN).
2787 // The resolution can be set via the command (SATT * FILL VALUE), where (VALUE)
2788 // is the ratio between the number of pixels drawn and 20 (user coordinates).
2789 // Parallel view and perspective view are possible (DOPT PROJ PARA/PERS); in the
2790 // first case, we assume that the first mother volume of the tree is a box with
2791 // dimensions 10000 X 10000 X 10000 cm and the view point (infinetely far) is
2792 // 5000 cm far from the origin along the Z axis of the user coordinates; in the
2793 // second case, the distance between the observer and the origin of the world
2794 // reference system is set in cm by the command (PERSP NAME VALUE); grand-angle
2795 // or telescopic effects can be achieved changing the scale factors in the command
2796 // (DRAW). When the final picture does not occupy the full window,
2797 // mapping the space before tracing can speed up the drawing, but can also
2798 // produce less precise results; values from 1 to 4 are allowed in the command
2799 // (DOPT MAPP VALUE), the mapping being more precise for increasing (VALUE); for
2800 // (VALUE = 0) no mapping is performed (therefore max precision and lowest speed).
2801 // The command (VALCUT) allows the cutting of the detector by three planes
2802 // ortogonal to the x,y,z axis. The attribute (LSTY) can be set by the command
2803 // SATT for any desired volume and can assume values from 0 to 7; it determines
2804 // the different light processing to be performed for different materials:
2805 // 0 = dark-matt, 1 = bright-matt, 2 = plastic, 3 = ceramic, 4 = rough-metals,
2806 // 5 = shiny-metals, 6 = glass, 7 = mirror. The detector is assumed to be in the
2807 // dark, the ambient light luminosity is 0.2 for each basic hue (the saturation
2808 // is 0.9) and the observer is assumed to have a light source (therefore he will
2809 // produce parallel light in the case of parallel view and point-like-source
2810 // light in the case of perspective view).
2812 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2814 //_____________________________________________________________________________
2815 void TGeant3::Gsatt(const char *name, const char *att, Int_t val)
2819 // IOPT Name of the attribute to be set
2820 // IVAL Value to which the attribute is to be set
2822 // name= "*" stands for all the volumes.
2823 // iopt can be chosen among the following :
2825 // WORK 0=volume name is inactive for the tracking
2826 // 1=volume name is active for the tracking (default)
2828 // SEEN 0=volume name is invisible
2829 // 1=volume name is visible (default)
2830 // -1=volume invisible with all its descendants in the tree
2831 // -2=volume visible but not its descendants in the tree
2833 // LSTY line style 1,2,3,... (default=1)
2834 // LSTY=7 will produce a very precise approximation for
2835 // revolution bodies.
2837 // LWID line width -7,...,1,2,3,..7 (default=1)
2838 // LWID<0 will act as abs(LWID) was set for the volume
2839 // and for all the levels below it. When SHAD is 'ON', LWID
2840 // represent the linewidth of the scan lines filling the surfaces
2841 // (whereas the FILL value represent their number). Therefore
2842 // tuning this parameter will help to obtain the desired
2843 // quality/performance ratio.
2845 // COLO colour code -166,...,1,2,..166 (default=1)
2847 // n=2=red; n=17+m, m=0,25, increasing luminosity according to 'm';
2848 // n=3=green; n=67+m, m=0,25, increasing luminosity according to 'm';
2849 // n=4=blue; n=117+m, m=0,25, increasing luminosity according to 'm';
2850 // n=5=yellow; n=42+m, m=0,25, increasing luminosity according to 'm';
2851 // n=6=violet; n=142+m, m=0,25, increasing luminosity according to 'm';
2852 // n=7=lightblue; n=92+m, m=0,25, increasing luminosity according to 'm';
2853 // colour=n*10+m, m=1,2,...9, will produce the same colour
2854 // as 'n', but with increasing luminosity according to 'm';
2855 // COLO<0 will act as if abs(COLO) was set for the volume
2856 // and for all the levels below it.
2857 // When for a volume the attribute FILL is > 1 (and the
2858 // option SHAD is on), the ABS of its colour code must be < 8
2859 // because an automatic shading of its faces will be
2862 // FILL (1992) fill area -7,...,0,1,...7 (default=0)
2863 // when option SHAD is "on" the FILL attribute of any
2864 // volume can be set different from 0 (normal drawing);
2865 // if it is set to 1, the faces of such volume will be filled
2866 // with solid colours; if ABS(FILL) is > 1, then a light
2867 // source is placed along the observer line, and the faces of
2868 // such volumes will be painted by colours whose luminosity
2869 // will depend on the amount of light reflected;
2870 // if ABS(FILL) = 1, then it is possible to use all the 166
2871 // colours of the colour table, becouse the automatic shading
2872 // is not performed;
2873 // for increasing values of FILL the drawing will be performed
2874 // with higher and higher resolution improving the quality (the
2875 // number of scan lines used to fill the faces increases with FILL);
2876 // it is possible to set different values of FILL
2877 // for different volumes, in order to optimize at the same time
2878 // the performance and the quality of the picture;
2879 // FILL<0 will act as if abs(FILL) was set for the volume
2880 // and for all the levels below it.
2881 // This kind of drawing can be saved in 'picture files'
2882 // or in view banks.
2883 // 0=drawing without fill area
2884 // 1=faces filled with solid colours and resolution = 6
2885 // 2=lowest resolution (very fast)
2886 // 3=default resolution
2887 // 4=.................
2888 // 5=.................
2889 // 6=.................
2891 // Finally, if a coloured background is desired, the FILL
2892 // attribute for the first volume of the tree must be set
2893 // equal to -abs(colo), colo being >0 and <166.
2895 // SET set number associated to volume name
2896 // DET detector number associated to volume name
2897 // DTYP detector type (1,2)
2904 gsatt(PASSCHARD(vname), PASSCHARD(vatt), val PASSCHARL(vname)
2908 //_____________________________________________________________________________
2909 void TGeant3::Gfpara(const char *name, Int_t number, Int_t intext, Int_t& npar,
2910 Int_t& natt, Float_t* par, Float_t* att)
2913 // Find the parameters of a volume
2915 gfpara(PASSCHARD(name), number, intext, npar, natt, par, att
2919 //_____________________________________________________________________________
2920 void TGeant3::Gckpar(Int_t ish, Int_t npar, Float_t* par)
2923 // Check the parameters of a shape
2925 gckpar(ish,npar,par);
2928 //_____________________________________________________________________________
2929 void TGeant3::Gckmat(Int_t itmed, char* natmed)
2932 // Check the parameters of a tracking medium
2934 gckmat(itmed, PASSCHARD(natmed) PASSCHARL(natmed));
2937 //_____________________________________________________________________________
2938 Int_t TGeant3::Glvolu(Int_t nlev, Int_t *lnam,Int_t *lnum)
2941 // nlev number of leveles deap into the volume tree
2942 // size of the arrays lnam and lnum
2943 // lnam an integer array whos 4 bytes contain the askii code for the
2945 // lnum an integer array containing the copy numbers for that specific
2948 // This routine fills the volulme paramters in common /gcvolu/ for a
2949 // physical tree, specified by the list lnam and lnum of volume names
2950 // and numbers, and for all its ascendants up to level 1. This routine
2951 // is optimsed and does not re-compute the part of the history already
2952 // available in GCVOLU. This means that if it is used in user programs
2953 // outside the usual framwork of the tracking, the user has to initilise
2954 // to zero NLEVEL in the common GCVOLU. It return 0 if there were no
2955 // problems in make the call.
2958 glvolu(nlev, lnam, lnum, ier);
2962 //_____________________________________________________________________________
2963 void TGeant3::Gdelete(Int_t iview)
2966 // IVIEW View number
2968 // It deletes a view bank from memory.
2973 //_____________________________________________________________________________
2974 void TGeant3::Gdopen(Int_t iview)
2977 // IVIEW View number
2979 // When a drawing is very complex and requires a long time to be
2980 // executed, it can be useful to store it in a view bank: after a
2981 // call to DOPEN and the execution of the drawing (nothing will
2982 // appear on the screen), and after a necessary call to DCLOSE,
2983 // the contents of the bank can be displayed in a very fast way
2984 // through a call to DSHOW; therefore, the detector can be easily
2985 // zoomed many times in different ways. Please note that the pictures
2986 // with solid colours can now be stored in a view bank or in 'PICTURE FILES'
2993 //_____________________________________________________________________________
2994 void TGeant3::Gdclose()
2997 // It closes the currently open view bank; it must be called after the
2998 // end of the drawing to be stored.
3003 //_____________________________________________________________________________
3004 void TGeant3::Gdshow(Int_t iview)
3007 // IVIEW View number
3009 // It shows on the screen the contents of a view bank. It
3010 // can be called after a view bank has been closed.
3015 //_____________________________________________________________________________
3016 void TGeant3::Gdopt(const char *name,const char *value)
3020 // VALUE Option value
3022 // To set/modify the drawing options.
3025 // THRZ ON Draw tracks in R vs Z
3026 // OFF (D) Draw tracks in X,Y,Z
3029 // PROJ PARA (D) Parallel projection
3031 // TRAK LINE (D) Trajectory drawn with lines
3032 // POIN " " with markers
3033 // HIDE ON Hidden line removal using the CG package
3034 // OFF (D) No hidden line removal
3035 // SHAD ON Fill area and shading of surfaces.
3036 // OFF (D) Normal hidden line removal.
3037 // RAYT ON Ray-tracing on.
3038 // OFF (D) Ray-tracing off.
3039 // EDGE OFF Does not draw contours when shad is on.
3040 // ON (D) Normal shading.
3041 // MAPP 1,2,3,4 Mapping before ray-tracing.
3042 // 0 (D) No mapping.
3043 // USER ON User graphics options in the raytracing.
3044 // OFF (D) Automatic graphics options.
3050 Vname(value,vvalue);
3051 gdopt(PASSCHARD(vname), PASSCHARD(vvalue) PASSCHARL(vname)
3055 //_____________________________________________________________________________
3056 void TGeant3::Gdraw(const char *name,Float_t theta, Float_t phi, Float_t psi,
3057 Float_t u0,Float_t v0,Float_t ul,Float_t vl)
3062 // THETA Viewing angle theta (for 3D projection)
3063 // PHI Viewing angle phi (for 3D projection)
3064 // PSI Viewing angle psi (for 2D rotation)
3065 // U0 U-coord. (horizontal) of volume origin
3066 // V0 V-coord. (vertical) of volume origin
3067 // SU Scale factor for U-coord.
3068 // SV Scale factor for V-coord.
3070 // This function will draw the volumes,
3071 // selected with their graphical attributes, set by the Gsatt
3072 // facility. The drawing may be performed with hidden line removal
3073 // and with shading effects according to the value of the options HIDE
3074 // and SHAD; if the option SHAD is ON, the contour's edges can be
3075 // drawn or not. If the option HIDE is ON, the detector can be
3076 // exploded (BOMB), clipped with different shapes (CVOL), and some
3077 // of its parts can be shifted from their original
3078 // position (SHIFT). When HIDE is ON, if
3079 // the drawing requires more than the available memory, the program
3080 // will evaluate and display the number of missing words
3081 // (so that the user can increase the
3082 // size of its ZEBRA store). Finally, at the end of each drawing (with HIDE on),
3083 // the program will print messages about the memory used and
3084 // statistics on the volumes' visibility.
3085 // The following commands will produce the drawing of a green
3086 // volume, specified by NAME, without using the hidden line removal
3087 // technique, using the hidden line removal technique,
3088 // with different linewidth and colour (red), with
3089 // solid colour, with shading of surfaces, and without edges.
3090 // Finally, some examples are given for the ray-tracing. (A possible
3091 // string for the NAME of the volume can be found using the command DTREE).
3097 if (fGcvdma->raytra != 1) {
3098 gdraw(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
3100 gdrayt(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
3104 //_____________________________________________________________________________
3105 void TGeant3::Gdrawc(const char *name,Int_t axis, Float_t cut,Float_t u0,
3106 Float_t v0,Float_t ul,Float_t vl)
3111 // CUTVAL Cut plane distance from the origin along the axis
3113 // U0 U-coord. (horizontal) of volume origin
3114 // V0 V-coord. (vertical) of volume origin
3115 // SU Scale factor for U-coord.
3116 // SV Scale factor for V-coord.
3118 // The cut plane is normal to caxis (X,Y,Z), corresponding to iaxis (1,2,3),
3119 // and placed at the distance cutval from the origin.
3120 // The resulting picture is seen from the the same axis.
3121 // When HIDE Mode is ON, it is possible to get the same effect with
3122 // the CVOL/BOX function.
3128 gdrawc(PASSCHARD(vname), axis,cut,u0,v0,ul,vl PASSCHARL(vname));
3131 //_____________________________________________________________________________
3132 void TGeant3::Gdrawx(const char *name,Float_t cutthe, Float_t cutphi,
3133 Float_t cutval, Float_t theta, Float_t phi, Float_t u0,
3134 Float_t v0,Float_t ul,Float_t vl)
3138 // CUTTHE Theta angle of the line normal to cut plane
3139 // CUTPHI Phi angle of the line normal to cut plane
3140 // CUTVAL Cut plane distance from the origin along the axis
3142 // THETA Viewing angle theta (for 3D projection)
3143 // PHI Viewing angle phi (for 3D projection)
3144 // U0 U-coord. (horizontal) of volume origin
3145 // V0 V-coord. (vertical) of volume origin
3146 // SU Scale factor for U-coord.
3147 // SV Scale factor for V-coord.
3149 // The cut plane is normal to the line given by the cut angles
3150 // cutthe and cutphi and placed at the distance cutval from the origin.
3151 // The resulting picture is seen from the viewing angles theta,phi.
3157 gdrawx(PASSCHARD(vname), cutthe,cutphi,cutval,theta,phi,u0,v0,ul,vl
3161 //_____________________________________________________________________________
3162 void TGeant3::Gdhead(Int_t isel, const char *name, Float_t chrsiz)
3167 // ISEL Option flag D=111110
3169 // CHRSIZ Character size (cm) of title NAME D=0.6
3172 // 0 to have only the header lines
3173 // xxxxx1 to add the text name centered on top of header
3174 // xxxx1x to add global detector name (first volume) on left
3175 // xxx1xx to add date on right
3176 // xx1xxx to select thick characters for text on top of header
3177 // x1xxxx to add the text 'EVENT NR x' on top of header
3178 // 1xxxxx to add the text 'RUN NR x' on top of header
3179 // NOTE that ISEL=x1xxx1 or ISEL=1xxxx1 are illegal choices,
3180 // i.e. they generate overwritten text.
3182 gdhead(isel,PASSCHARD(name),chrsiz PASSCHARL(name));
3185 //_____________________________________________________________________________
3186 void TGeant3::Gdman(Float_t u, Float_t v, const char *type)
3189 // Draw a 2D-man at position (U0,V0)
3191 // U U-coord. (horizontal) of the centre of man' R
3192 // V V-coord. (vertical) of the centre of man' R
3193 // TYPE D='MAN' possible values: 'MAN,WM1,WM2,WM3'
3195 // CALL GDMAN(u,v),CALL GDWMN1(u,v),CALL GDWMN2(u,v),CALL GDWMN2(u,v)
3196 // It superimposes the picure of a man or of a woman, chosen among
3197 // three different ones, with the same scale factors as the detector
3198 // in the current drawing.
3201 if (opt.Contains("WM1")) {
3203 } else if (opt.Contains("WM3")) {
3205 } else if (opt.Contains("WM2")) {
3212 //_____________________________________________________________________________
3213 void TGeant3::Gdspec(const char *name)
3218 // Shows 3 views of the volume (two cut-views and a 3D view), together with
3219 // its geometrical specifications. The 3D drawing will
3220 // be performed according the current values of the options HIDE and
3221 // SHAD and according the current SetClipBox clipping parameters for that
3228 gdspec(PASSCHARD(vname) PASSCHARL(vname));
3231 //_____________________________________________________________________________
3232 void TGeant3::DrawOneSpec(const char *name)
3235 // Function called when one double-clicks on a volume name
3236 // in a TPavelabel drawn by Gdtree.
3238 THIGZ *higzSave = gHigz;
3239 higzSave->SetName("higzSave");
3240 THIGZ *higzSpec = (THIGZ*)gROOT->FindObject("higzSpec");
3241 //printf("DrawOneSpec, gHigz=%x, higzSpec=%x\n",gHigz,higzSpec);
3242 if (higzSpec) gHigz = higzSpec;
3243 else higzSpec = new THIGZ(kDefSize);
3244 higzSpec->SetName("higzSpec");
3249 gdspec(PASSCHARD(vname) PASSCHARL(vname));
3252 higzSave->SetName("higz");
3256 //_____________________________________________________________________________
3257 void TGeant3::Gdtree(const char *name,Int_t levmax, Int_t isel)
3261 // LEVMAX Depth level
3264 // This function draws the logical tree,
3265 // Each volume in the tree is represented by a TPaveTree object.
3266 // Double-clicking on a TPaveTree draws the specs of the corresponding volume.
3267 // Use TPaveTree pop-up menu to select:
3270 // - drawing tree of parent
3276 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
3277 gHigz->SetPname("");
3280 //_____________________________________________________________________________
3281 void TGeant3::GdtreeParent(const char *name,Int_t levmax, Int_t isel)
3285 // LEVMAX Depth level
3288 // This function draws the logical tree of the parent of name.
3292 // Scan list of volumes in JVOLUM
3294 Int_t gname, i, jvo, in, nin, jin, num;
3295 strncpy((char *) &gname, name, 4);
3296 for(i=1; i<=fGcnum->nvolum; i++) {
3297 jvo = fZlq[fGclink->jvolum-i];
3298 nin = Int_t(fZq[jvo+3]);
3299 if (nin == -1) nin = 1;
3300 for (in=1;in<=nin;in++) {
3302 num = Int_t(fZq[jin+2]);
3303 if(gname == fZiq[fGclink->jvolum+num]) {
3304 strncpy(vname,(char*)&fZiq[fGclink->jvolum+i],4);
3306 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
3307 gHigz->SetPname("");
3314 //_____________________________________________________________________________
3315 void TGeant3::SetABAN(Int_t par)
3318 // par = 1 particles will be stopped according to their residual
3319 // range if they are not in a sensitive material and are
3320 // far enough from the boundary
3321 // 0 particles are transported normally
3323 fGcphys->dphys1 = par;
3327 //_____________________________________________________________________________
3328 void TGeant3::SetANNI(Int_t par)
3331 // To control positron annihilation.
3332 // par =0 no annihilation
3333 // =1 annihilation. Decays processed.
3334 // =2 annihilation. No decay products stored.
3336 fGcphys->ianni = par;
3340 //_____________________________________________________________________________
3341 void TGeant3::SetAUTO(Int_t par)
3344 // To control automatic calculation of tracking medium parameters:
3345 // par =0 no automatic calculation;
3346 // =1 automati calculation.
3348 fGctrak->igauto = par;
3352 //_____________________________________________________________________________
3353 void TGeant3::SetBOMB(Float_t boom)
3356 // BOOM : Exploding factor for volumes position
3358 // To 'explode' the detector. If BOOM is positive (values smaller
3359 // than 1. are suggested, but any value is possible)
3360 // all the volumes are shifted by a distance
3361 // proportional to BOOM along the direction between their centre
3362 // and the origin of the MARS; the volumes which are symmetric
3363 // with respect to this origin are simply not shown.
3364 // BOOM equal to 0 resets the normal mode.
3365 // A negative (greater than -1.) value of
3366 // BOOM will cause an 'implosion'; for even lower values of BOOM
3367 // the volumes' positions will be reflected respect to the origin.
3368 // This command can be useful to improve the 3D effect for very
3369 // complex detectors. The following commands will make explode the
3376 //_____________________________________________________________________________
3377 void TGeant3::SetBREM(Int_t par)
3380 // To control bremstrahlung.
3381 // par =0 no bremstrahlung
3382 // =1 bremstrahlung. Photon processed.
3383 // =2 bremstrahlung. No photon stored.
3385 fGcphys->ibrem = par;
3389 //_____________________________________________________________________________
3390 void TGeant3::SetCKOV(Int_t par)
3393 // To control Cerenkov production
3394 // par =0 no Cerenkov;
3396 // =2 Cerenkov with primary stopped at each step.
3398 fGctlit->itckov = par;
3402 //_____________________________________________________________________________
3403 void TGeant3::SetClipBox(const char *name,Float_t xmin,Float_t xmax,
3404 Float_t ymin,Float_t ymax,Float_t zmin,Float_t zmax)
3407 // The hidden line removal technique is necessary to visualize properly
3408 // very complex detectors. At the same time, it can be useful to visualize
3409 // the inner elements of a detector in detail. This function allows
3410 // subtractions (via boolean operation) of BOX shape from any part of
3411 // the detector, therefore showing its inner contents.
3412 // If "*" is given as the name of the
3413 // volume to be clipped, all volumes are clipped by the given box.
3414 // A volume can be clipped at most twice.
3415 // if a volume is explicitely clipped twice,
3416 // the "*" will not act on it anymore. Giving "." as the name
3417 // of the volume to be clipped will reset the clipping.
3419 // NAME Name of volume to be clipped
3421 // XMIN Lower limit of the Shape X coordinate
3422 // XMAX Upper limit of the Shape X coordinate
3423 // YMIN Lower limit of the Shape Y coordinate
3424 // YMAX Upper limit of the Shape Y coordinate
3425 // ZMIN Lower limit of the Shape Z coordinate
3426 // ZMAX Upper limit of the Shape Z coordinate
3428 // This function performs a boolean subtraction between the volume
3429 // NAME and a box placed in the MARS according the values of the given
3435 setclip(PASSCHARD(vname),xmin,xmax,ymin,ymax,zmin,zmax PASSCHARL(vname));
3438 //_____________________________________________________________________________
3439 void TGeant3::SetCOMP(Int_t par)
3442 // To control Compton scattering
3443 // par =0 no Compton
3444 // =1 Compton. Electron processed.
3445 // =2 Compton. No electron stored.
3448 fGcphys->icomp = par;
3451 //_____________________________________________________________________________
3452 void TGeant3::SetCUTS(Float_t cutgam,Float_t cutele,Float_t cutneu,
3453 Float_t cuthad,Float_t cutmuo ,Float_t bcute ,
3454 Float_t bcutm ,Float_t dcute ,Float_t dcutm ,
3455 Float_t ppcutm, Float_t tofmax)
3458 // CUTGAM Cut for gammas D=0.001
3459 // CUTELE Cut for electrons D=0.001
3460 // CUTHAD Cut for charged hadrons D=0.01
3461 // CUTNEU Cut for neutral hadrons D=0.01
3462 // CUTMUO Cut for muons D=0.01
3463 // BCUTE Cut for electron brems. D=-1.
3464 // BCUTM Cut for muon brems. D=-1.
3465 // DCUTE Cut for electron delta-rays D=-1.
3466 // DCUTM Cut for muon delta-rays D=-1.
3467 // PPCUTM Cut for e+e- pairs by muons D=0.01
3468 // TOFMAX Time of flight cut D=1.E+10
3470 // If the default values (-1.) for BCUTE ,BCUTM ,DCUTE ,DCUTM
3471 // are not modified, they will be set to CUTGAM,CUTGAM,CUTELE,CUTELE
3473 // If one of the parameters from CUTGAM to PPCUTM included
3474 // is modified, cross-sections and energy loss tables must be
3475 // recomputed via the function Gphysi.
3477 fGccuts->cutgam = cutgam;
3478 fGccuts->cutele = cutele;
3479 fGccuts->cutneu = cutneu;
3480 fGccuts->cuthad = cuthad;
3481 fGccuts->cutmuo = cutmuo;
3482 fGccuts->bcute = bcute;
3483 fGccuts->bcutm = bcutm;
3484 fGccuts->dcute = dcute;
3485 fGccuts->dcutm = dcutm;
3486 fGccuts->ppcutm = ppcutm;
3487 fGccuts->tofmax = tofmax;
3490 //_____________________________________________________________________________
3491 void TGeant3::SetDCAY(Int_t par)
3494 // To control Decay mechanism.
3495 // par =0 no decays.
3496 // =1 Decays. secondaries processed.
3497 // =2 Decays. No secondaries stored.
3499 fGcphys->idcay = par;
3503 //_____________________________________________________________________________
3504 void TGeant3::SetDEBU(Int_t emin, Int_t emax, Int_t emod)
3507 // Set the debug flag and frequency
3508 // Selected debug output will be printed from
3509 // event emin to even emax each emod event
3511 fGcflag->idemin = emin;
3512 fGcflag->idemax = emax;
3513 fGcflag->itest = emod;
3517 //_____________________________________________________________________________
3518 void TGeant3::SetDRAY(Int_t par)
3521 // To control delta rays mechanism.
3522 // par =0 no delta rays.
3523 // =1 Delta rays. secondaries processed.
3524 // =2 Delta rays. No secondaries stored.
3526 fGcphys->idray = par;
3529 //_____________________________________________________________________________
3530 void TGeant3::SetERAN(Float_t ekmin, Float_t ekmax, Int_t nekbin)
3533 // To control cross section tabulations
3534 // ekmin = minimum kinetic energy in GeV
3535 // ekmax = maximum kinetic energy in GeV
3536 // nekbin = number of logatithmic bins (<200)
3538 fGcmulo->ekmin = ekmin;
3539 fGcmulo->ekmax = ekmax;
3540 fGcmulo->nekbin = nekbin;
3543 //_____________________________________________________________________________
3544 void TGeant3::SetHADR(Int_t par)
3547 // To control hadronic interactions.
3548 // par =0 no hadronic interactions.
3549 // =1 Hadronic interactions. secondaries processed.
3550 // =2 Hadronic interactions. No secondaries stored.
3552 fGcphys->ihadr = par;
3555 //_____________________________________________________________________________
3556 void TGeant3::SetKINE(Int_t kine, Float_t xk1, Float_t xk2, Float_t xk3,
3557 Float_t xk4, Float_t xk5, Float_t xk6, Float_t xk7,
3558 Float_t xk8, Float_t xk9, Float_t xk10)
3561 // Set the variables in /GCFLAG/ IKINE, PKINE(10)
3562 // Their meaning is user defined
3564 fGckine->ikine = kine;
3565 fGckine->pkine[0] = xk1;
3566 fGckine->pkine[1] = xk2;
3567 fGckine->pkine[2] = xk3;
3568 fGckine->pkine[3] = xk4;
3569 fGckine->pkine[4] = xk5;
3570 fGckine->pkine[5] = xk6;
3571 fGckine->pkine[6] = xk7;
3572 fGckine->pkine[7] = xk8;
3573 fGckine->pkine[8] = xk9;
3574 fGckine->pkine[9] = xk10;
3577 //_____________________________________________________________________________
3578 void TGeant3::SetLOSS(Int_t par)
3581 // To control energy loss.
3582 // par =0 no energy loss;
3583 // =1 restricted energy loss fluctuations;
3584 // =2 complete energy loss fluctuations;
3586 // =4 no energy loss fluctuations.
3587 // If the value ILOSS is changed, then cross-sections and energy loss
3588 // tables must be recomputed via the command 'PHYSI'.
3590 fGcphys->iloss = par;
3594 //_____________________________________________________________________________
3595 void TGeant3::SetMULS(Int_t par)
3598 // To control multiple scattering.
3599 // par =0 no multiple scattering.
3600 // =1 Moliere or Coulomb scattering.
3601 // =2 Moliere or Coulomb scattering.
3602 // =3 Gaussian scattering.
3604 fGcphys->imuls = par;
3608 //_____________________________________________________________________________
3609 void TGeant3::SetMUNU(Int_t par)
3612 // To control muon nuclear interactions.
3613 // par =0 no muon-nuclear interactions.
3614 // =1 Nuclear interactions. Secondaries processed.
3615 // =2 Nuclear interactions. Secondaries not processed.
3617 fGcphys->imunu = par;
3620 //_____________________________________________________________________________
3621 void TGeant3::SetOPTI(Int_t par)
3624 // This flag controls the tracking optimisation performed via the
3626 // 1 no optimisation at all; GSORD calls disabled;
3627 // 0 no optimisation; only user calls to GSORD kept;
3628 // 1 all non-GSORDered volumes are ordered along the best axis;
3629 // 2 all volumes are ordered along the best axis.
3631 fGcopti->ioptim = par;
3634 //_____________________________________________________________________________
3635 void TGeant3::SetPAIR(Int_t par)
3638 // To control pair production mechanism.
3639 // par =0 no pair production.
3640 // =1 Pair production. secondaries processed.
3641 // =2 Pair production. No secondaries stored.
3643 fGcphys->ipair = par;
3647 //_____________________________________________________________________________
3648 void TGeant3::SetPFIS(Int_t par)
3651 // To control photo fission mechanism.
3652 // par =0 no photo fission.
3653 // =1 Photo fission. secondaries processed.
3654 // =2 Photo fission. No secondaries stored.
3656 fGcphys->ipfis = par;
3659 //_____________________________________________________________________________
3660 void TGeant3::SetPHOT(Int_t par)
3663 // To control Photo effect.
3664 // par =0 no photo electric effect.
3665 // =1 Photo effect. Electron processed.
3666 // =2 Photo effect. No electron stored.
3668 fGcphys->iphot = par;
3671 //_____________________________________________________________________________
3672 void TGeant3::SetRAYL(Int_t par)
3675 // To control Rayleigh scattering.
3676 // par =0 no Rayleigh scattering.
3679 fGcphys->irayl = par;
3682 //_____________________________________________________________________________
3683 void TGeant3::SetSTRA(Int_t par)
3686 // To control energy loss fluctuations
3687 // with the PhotoAbsorption Ionisation model.
3688 // par =0 no Straggling.
3689 // =1 Straggling yes => no Delta rays.
3691 fGcphlt->istra = par;
3694 //_____________________________________________________________________________
3695 void TGeant3::SetSWIT(Int_t sw, Int_t val)
3699 // val New switch value
3701 // Change one element of array ISWIT(10) in /GCFLAG/
3703 if (sw <= 0 || sw > 10) return;
3704 fGcflag->iswit[sw-1] = val;
3708 //_____________________________________________________________________________
3709 void TGeant3::SetTRIG(Int_t nevents)
3712 // Set number of events to be run
3714 fGcflag->nevent = nevents;
3717 //_____________________________________________________________________________
3718 void TGeant3::SetUserDecay(Int_t pdg)
3721 // Force the decays of particles to be done with Pythia
3722 // and not with the Geant routines.
3723 // just kill pointers doing mzdrop
3725 Int_t ipart = IdFromPDG(pdg);
3727 printf("Particle %d not in geant\n",pdg);
3730 Int_t jpart=fGclink->jpart;
3731 Int_t jpa=fZlq[jpart-ipart];
3734 Int_t jpa1=fZlq[jpa-1];
3736 mzdrop(fGcbank->ixcons,jpa1,PASSCHARD(" ") PASSCHARL(" "));
3737 Int_t jpa2=fZlq[jpa-2];
3739 mzdrop(fGcbank->ixcons,jpa2,PASSCHARD(" ") PASSCHARL(" "));
3743 //______________________________________________________________________________
3744 void TGeant3::Vname(const char *name, char *vname)
3747 // convert name to upper case. Make vname at least 4 chars
3749 Int_t l = strlen(name);
3752 for (i=0;i<l;i++) vname[i] = toupper(name[i]);
3753 for (i=l;i<4;i++) vname[i] = ' ';
3757 //______________________________________________________________________________
3758 void TGeant3::Ertrgo()
3761 // Perform the tracking of the track Track parameters are in VECT
3766 //______________________________________________________________________________
3767 void TGeant3::Ertrak(const Float_t *const x1, const Float_t *const p1,
3768 const Float_t *x2, const Float_t *p2,
3769 Int_t ipa, Option_t *chopt)
3771 //************************************************************************
3773 //* Perform the tracking of the track from point X1 to *
3775 //* (Before calling this routine the user should also provide *
3776 //* the input informations in /EROPTS/ and /ERTRIO/ *
3777 //* using subroutine EUFIL(L/P/V) *
3778 //* X1 - Starting coordinates (Cartesian) *
3779 //* P1 - Starting 3-momentum (Cartesian) *
3780 //* X2 - Final coordinates (Cartesian) *
3781 //* P2 - Final 3-momentum (Cartesian) *
3782 //* IPA - Particle code (a la GEANT) of the track *
3785 //* 'B' 'Backward tracking' - i.e. energy loss *
3786 //* added to the current energy *
3787 //* 'E' 'Exact' calculation of errors assuming *
3788 //* helix (i.e. pathlength not *
3789 //* assumed as infinitesimal) *
3790 //* 'L' Tracking upto prescribed Lengths reached *
3791 //* 'M' 'Mixed' prediction (not yet coded) *
3792 //* 'O' Tracking 'Only' without calculating errors *
3793 //* 'P' Tracking upto prescribed Planes reached *
3794 //* 'V' Tracking upto prescribed Volumes reached *
3795 //* 'X' Tracking upto prescribed Point approached *
3797 //* Interface with GEANT : *
3798 //* Track parameters are in /CGKINE/ and /GCTRAK/ *
3800 //* ==>Called by : USER *
3801 //* Authors M.Maire, E.Nagy ********//* *
3803 //************************************************************************
3804 ertrak(x1,p1,x2,p2,ipa,PASSCHARD(chopt) PASSCHARL(chopt));
3807 //_____________________________________________________________________________
3808 void TGeant3::WriteEuclid(const char* filnam, const char* topvol,
3809 Int_t number, Int_t nlevel)
3813 // ******************************************************************
3815 // * Write out the geometry of the detector in EUCLID file format *
3817 // * filnam : will be with the extension .euc *
3818 // * topvol : volume name of the starting node *
3819 // * number : copy number of topvol (relevant for gsposp) *
3820 // * nlevel : number of levels in the tree structure *
3821 // * to be written out, starting from topvol *
3823 // * Author : M. Maire *
3825 // ******************************************************************
3827 // File filnam.tme is written out with the definitions of tracking
3828 // medias and materials.
3829 // As to restore original numbers for materials and medias, program
3830 // searches in the file euc_medi.dat and comparing main parameters of
3831 // the mat. defined inside geant and the one in file recognizes them
3832 // and is able to take number from file. If for any material or medium,
3833 // this procedure fails, ordering starts from 1.
3834 // Arrays IOTMED and IOMATE are used for this procedure
3836 const char kShape[][5]={"BOX ","TRD1","TRD2","TRAP","TUBE","TUBS","CONE",
3837 "CONS","SPHE","PARA","PGON","PCON","ELTU","HYPE",
3839 Int_t i, end, itm, irm, jrm, k, nmed;
3843 char *filext, *filetme;
3844 char natmed[21], namate[21];
3845 char natmedc[21], namatec[21];
3846 char key[5], name[5], mother[5], konly[5];
3848 Int_t iadvol, iadtmd, iadrot, nwtot, iret;
3849 Int_t mlevel, numbr, natt, numed, nin, ndata;
3850 Int_t iname, ivo, ish, jvo, nvstak, ivstak;
3851 Int_t jdiv, ivin, in, jin, jvin, irot;
3852 Int_t jtm, imat, jma, flag=0, imatc;
3853 Float_t az, dens, radl, absl, a, step, x, y, z;
3854 Int_t npar, ndvmx, left;
3855 Float_t zc, densc, radlc, abslc, c0, tmaxfd;
3857 Int_t iomate[100], iotmed[100];
3858 Float_t par[50], att[20], ubuf[50];
3861 Int_t level, ndiv, iaxe;
3862 Int_t itmedc, nmatc, isvolc, ifieldc, nwbufc, isvol, nmat, ifield, nwbuf;
3863 Float_t fieldmc, tmaxfdc, stemaxc, deemaxc, epsilc, stminc, fieldm;
3864 Float_t tmaxf, stemax, deemax, epsil, stmin;
3865 const char *k10000="!\n%s\n!\n";
3866 //Open the input file
3868 for(i=0;i<end;i++) if(filnam[i]=='.') {
3872 filext=new char[end+5];
3873 filetme=new char[end+5];
3874 strncpy(filext,filnam,end);
3875 strncpy(filetme,filnam,end);
3877 // *** The output filnam name will be with extension '.euc'
3878 strcpy(&filext[end],".euc");
3879 strcpy(&filetme[end],".tme");
3880 lun=fopen(filext,"w");
3882 // *** Initialisation of the working space
3883 iadvol=fGcnum->nvolum;
3884 iadtmd=iadvol+fGcnum->nvolum;
3885 iadrot=iadtmd+fGcnum->ntmed;
3886 if(fGclink->jrotm) {
3887 fGcnum->nrotm=fZiq[fGclink->jrotm-2];
3891 nwtot=iadrot+fGcnum->nrotm;
3892 qws = new float[nwtot+1];
3893 for (i=0;i<nwtot+1;i++) qws[i]=0;
3896 if(nlevel==0) mlevel=20;
3898 // *** find the top volume and put it in the stak
3899 numbr = number>0 ? number : 1;
3900 Gfpara(topvol,numbr,1,npar,natt,par,att);
3902 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3907 // *** authorized shape ?
3908 strncpy((char *)&iname, topvol, 4);
3910 for(i=1; i<=fGcnum->nvolum; i++) if(fZiq[fGclink->jvolum+i]==iname) {
3914 jvo = fZlq[fGclink->jvolum-ivo];
3915 ish = Int_t (fZq[jvo+2]);
3917 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3924 iws[iadvol+ivo] = level;
3927 //*** flag all volumes and fill the stak
3931 // pick the next volume in stak
3933 ivo = TMath::Abs(iws[ivstak]);
3934 jvo = fZlq[fGclink->jvolum - ivo];
3936 // flag the tracking medium
3937 numed = Int_t (fZq[jvo + 4]);
3938 iws[iadtmd + numed] = 1;
3940 // get the daughters ...
3941 level = iws[iadvol+ivo];
3942 if (level < mlevel) {
3944 nin = Int_t (fZq[jvo + 3]);
3946 // from division ...
3948 jdiv = fZlq[jvo - 1];
3949 ivin = Int_t (fZq[jdiv + 2]);
3951 iws[nvstak] = -ivin;
3952 iws[iadvol+ivin] = level;
3954 // from position ...
3955 } else if (nin > 0) {
3956 for(in=1; in<=nin; in++) {
3957 jin = fZlq[jvo - in];
3958 ivin = Int_t (fZq[jin + 2 ]);
3959 jvin = fZlq[fGclink->jvolum - ivin];
3960 ish = Int_t (fZq[jvin + 2]);
3961 // authorized shape ?
3963 // not yet flagged ?
3964 if (iws[iadvol+ivin]==0) {
3967 iws[iadvol+ivin] = level;
3969 // flag the rotation matrix
3970 irot = Int_t ( fZq[jin + 4 ]);
3971 if (irot > 0) iws[iadrot+irot] = 1;
3977 // next volume in stak ?
3978 if (ivstak < nvstak) goto L10;
3980 // *** restore original material and media numbers
3981 // file euc_medi.dat is needed to compare materials and medias
3983 FILE* luncor=fopen("euc_medi.dat","r");
3986 for(itm=1; itm<=fGcnum->ntmed; itm++) {
3987 if (iws[iadtmd+itm] > 0) {
3988 jtm = fZlq[fGclink->jtmed-itm];
3989 strncpy(natmed,(char *)&fZiq[jtm+1],20);
3990 imat = Int_t (fZq[jtm+6]);
3991 jma = fZlq[fGclink->jmate-imat];
3993 printf(" *** GWEUCL *** material not defined for tracking medium %5i %s\n",itm,natmed);
3996 strncpy(namate,(char *)&fZiq[jma+1],20);
3999 //** find the material original number
4002 iret=fscanf(luncor,"%4s,%130s",key,card);
4003 if(iret<=0) goto L26;
4005 if(!strcmp(key,"MATE")) {
4006 sscanf(card,"%d %s %f %f %f %f %f %d",&imatc,namatec,&az,&zc,&densc,&radlc,&abslc,&nparc);
4007 Gfmate(imat,namate,a,z,dens,radl,absl,par,npar);
4008 if(!strcmp(namatec,namate)) {
4009 if(az==a && zc==z && densc==dens && radlc==radl
4010 && abslc==absl && nparc==nparc) {
4013 printf("*** GWEUCL *** material : %3d '%s' restored as %3d\n",imat,namate,imatc);
4015 printf("*** GWEUCL *** different definitions for material: %s\n",namate);
4019 if(strcmp(key,"END") && !flag) goto L23;
4021 printf("*** GWEUCL *** cannot restore original number for material: %s\n",namate);
4025 //*** restore original tracking medium number
4028 iret=fscanf(luncor,"%4s,%130s",key,card);
4029 if(iret<=0) goto L26;
4031 if (!strcmp(key,"TMED")) {
4032 sscanf(card,"%d %s %d %d %d %f %f %f %f %f %f %d\n",
4033 &itmedc,natmedc,&nmatc,&isvolc,&ifieldc,&fieldmc,
4034 &tmaxfdc,&stemaxc,&deemaxc,&epsilc,&stminc,&nwbufc);
4035 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxf,stemax,deemax,
4036 epsil,stmin,ubuf,&nwbuf);
4037 if(!strcmp(natmedc,natmed)) {
4038 if (iomate[nmat]==nmatc && nwbuf==nwbufc) {
4041 printf("*** GWEUCL *** medium : %3d '%20s' restored as %3d\n",
4044 printf("*** GWEUCL *** different definitions for tracking medium: %s\n",natmed);
4048 if(strcmp(key,"END") && !flag) goto L24;
4050 printf("cannot restore original number for medium : %s\n",natmed);
4058 L26: printf("*** GWEUCL *** cannot read the data file\n");
4060 L29: if(luncor) fclose (luncor);
4063 // *** write down the tracking medium definition
4065 strcpy(card,"! Tracking medium");
4066 fprintf(lun,k10000,card);
4068 for(itm=1;itm<=fGcnum->ntmed;itm++) {
4069 if (iws[iadtmd+itm]>0) {
4070 jtm = fZlq[fGclink->jtmed-itm];
4071 strncpy(natmed,(char *)&fZiq[jtm+1],20);
4073 imat = Int_t (fZq[jtm+6]);
4074 jma = fZlq[fGclink->jmate-imat];
4075 //* order media from one, if comparing with database failed
4077 iotmed[itm]=++imxtmed;
4078 iomate[imat]=++imxmate;
4083 printf(" *** GWEUCL *** material not defined for tracking medium %5d %s\n",
4086 strncpy(namate,(char *)&fZiq[jma+1],20);
4089 fprintf(lun,"TMED %3d '%20s' %3d '%20s'\n",iotmed[itm],natmed,iomate[imat],namate);
4093 //* *** write down the rotation matrix
4095 strcpy(card,"! Reperes");
4096 fprintf(lun,k10000,card);
4098 for(irm=1;irm<=fGcnum->nrotm;irm++) {
4099 if (iws[iadrot+irm]>0) {
4100 jrm = fZlq[fGclink->jrotm-irm];
4101 fprintf(lun,"ROTM %3d",irm);
4102 for(k=11;k<=16;k++) fprintf(lun," %8.3f",fZq[jrm+k]);
4107 //* *** write down the volume definition
4109 strcpy(card,"! Volumes");
4110 fprintf(lun,k10000,card);
4112 for(ivstak=1;ivstak<=nvstak;ivstak++) {
4115 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivo],4);
4117 jvo = fZlq[fGclink->jvolum-ivo];
4118 ish = Int_t (fZq[jvo+2]);
4119 nmed = Int_t (fZq[jvo+4]);
4120 npar = Int_t (fZq[jvo+5]);
4122 if (ivstak>1) for(i=0;i<npar;i++) par[i]=fZq[jvo+7+i];
4123 Gckpar (ish,npar,par);
4124 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,kShape[ish-1],iotmed[nmed],npar);
4125 for(i=0;i<(npar-1)/6+1;i++) {
4128 for(k=0;k<(left<6?left:6);k++) fprintf(lun," %11.5f",par[i*6+k]);
4132 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,kShape[ish-1],iotmed[nmed],npar);
4137 //* *** write down the division of volumes
4139 fprintf(lun,k10000,"! Divisions");
4140 for(ivstak=1;ivstak<=nvstak;ivstak++) {
4141 ivo = TMath::Abs(iws[ivstak]);
4142 jvo = fZlq[fGclink->jvolum-ivo];
4143 ish = Int_t (fZq[jvo+2]);
4144 nin = Int_t (fZq[jvo+3]);
4145 //* this volume is divided ...
4148 iaxe = Int_t ( fZq[jdiv+1]);
4149 ivin = Int_t ( fZq[jdiv+2]);
4150 ndiv = Int_t ( fZq[jdiv+3]);
4153 jvin = fZlq[fGclink->jvolum-ivin];
4154 nmed = Int_t ( fZq[jvin+4]);
4155 strncpy(mother,(char *)&fZiq[fGclink->jvolum+ivo ],4);
4157 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivin],4);
4159 if ((step<=0.)||(ish>=11)) {
4160 //* volume with negative parameter or gsposp or pgon ...
4161 fprintf(lun,"DIVN '%4s' '%4s' %3d %3d\n",name,mother,ndiv,iaxe);
4162 } else if ((ndiv<=0)||(ish==10)) {
4163 //* volume with negative parameter or gsposp or para ...
4164 ndvmx = TMath::Abs(ndiv);
4165 fprintf(lun,"DIVT '%4s' '%4s' %11.5f %3d %3d %3d\n",
4166 name,mother,step,iaxe,iotmed[nmed],ndvmx);
4168 //* normal volume : all kind of division are equivalent
4169 fprintf(lun,"DVT2 '%4s' '%4s' %11.5f %3d %11.5f %3d %3d\n",
4170 name,mother,step,iaxe,c0,iotmed[nmed],ndiv);
4175 //* *** write down the the positionnement of volumes
4177 fprintf(lun,k10000,"! Positionnements\n");
4179 for(ivstak = 1;ivstak<=nvstak;ivstak++) {
4180 ivo = TMath::Abs(iws[ivstak]);
4181 strncpy(mother,(char*)&fZiq[fGclink->jvolum+ivo ],4);
4183 jvo = fZlq[fGclink->jvolum-ivo];
4184 nin = Int_t( fZq[jvo+3]);
4185 //* this volume has daughters ...
4187 for (in=1;in<=nin;in++) {
4189 ivin = Int_t (fZq[jin +2]);
4190 numb = Int_t (fZq[jin +3]);
4191 irot = Int_t (fZq[jin +4]);
4195 strcpy(konly,"ONLY");
4196 if (fZq[jin+8]!=1.) strcpy(konly,"MANY");
4197 strncpy(name,(char*)&fZiq[fGclink->jvolum+ivin],4);
4199 jvin = fZlq[fGclink->jvolum-ivin];
4200 ish = Int_t (fZq[jvin+2]);
4201 //* gspos or gsposp ?
4202 ndata = fZiq[jin-1];
4204 fprintf(lun,"POSI '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s'\n",
4205 name,numb,mother,x,y,z,irot,konly);
4207 npar = Int_t (fZq[jin+9]);
4208 for(i=0;i<npar;i++) par[i]=fZq[jin+10+i];
4209 Gckpar (ish,npar,par);
4210 fprintf(lun,"POSP '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s' %3d\n",
4211 name,numb,mother,x,y,z,irot,konly,npar);
4213 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
4220 fprintf(lun,"END\n");
4223 //****** write down the materials and medias *****
4225 lun=fopen(filetme,"w");
4227 for(itm=1;itm<=fGcnum->ntmed;itm++) {
4228 if (iws[iadtmd+itm]>0) {
4229 jtm = fZlq[fGclink->jtmed-itm];
4230 strncpy(natmed,(char*)&fZiq[jtm+1],4);
4231 imat = Int_t (fZq[jtm+6]);
4232 jma = Int_t (fZlq[fGclink->jmate-imat]);
4234 Gfmate (imat,namate,a,z,dens,radl,absl,par,npar);
4235 fprintf(lun,"MATE %4d '%20s'%11.5E %11.5E %11.5E %11.5E %11.5E %3d\n",
4236 iomate[imat],namate,a,z,dens,radl,absl,npar);
4240 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
4244 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxfd,stemax,deemax,epsil,stmin,par,&npar);
4245 fprintf(lun,"TMED %4d '%20s' %3d %1d %3d %11.5f %11.5f %11.5f %11.5f %11.5f %11.5f %3d\n",
4246 iotmed[itm],natmed,iomate[nmat],isvol,ifield,
4247 fieldm,tmaxfd,stemax,deemax,epsil,stmin,npar);
4251 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
4257 fprintf(lun,"END\n");
4259 printf(" *** GWEUCL *** file: %s is now written out\n",filext);
4260 printf(" *** GWEUCL *** file: %s is now written out\n",filetme);
4269 //_____________________________________________________________________________
4270 void TGeant3::Streamer(TBuffer &R__b)
4273 // Stream an object of class TGeant3.
4275 if (R__b.IsReading()) {
4276 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
4277 AliMC::Streamer(R__b);
4280 R__b.ReadStaticArray(fPDGCode);
4282 R__b.WriteVersion(TGeant3::IsA());
4283 AliMC::Streamer(R__b);
4286 R__b.WriteArray(fPDGCode, fNPDGCodes);