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.21 2000/01/17 19:41:17 fca
21 Revision 1.20 2000/01/12 11:29:27 fca
24 Revision 1.19 1999/12/17 09:03:12 fca
25 Introduce a names array
27 Revision 1.18 1999/11/26 16:55:39 fca
28 Reimplement CurrentVolName() to avoid memory leaks
30 Revision 1.17 1999/11/03 16:58:28 fca
31 Correct source of address violation in creating character string
33 Revision 1.16 1999/11/03 13:17:08 fca
34 Have ProdProcess return const char*
36 Revision 1.15 1999/10/26 06:04:50 fca
37 Introduce TLorentzVector in AliMC::GetSecondary. Thanks to I.Hrivnacova
39 Revision 1.14 1999/09/29 09:24:30 fca
40 Introduction of the Copyright and cvs Log
44 ///////////////////////////////////////////////////////////////////////////////
46 // Interface Class to the Geant3.21 MonteCarlo //
50 <img src="picts/TGeant3Class.gif">
55 ///////////////////////////////////////////////////////////////////////////////
61 #include <TDatabasePDG.h>
62 #include "AliCallf77.h"
65 # define gzebra gzebra_
66 # define grfile grfile_
67 # define gpcxyz gpcxyz_
68 # define ggclos ggclos_
71 # define gcinit gcinit_
74 # define gtrigc gtrigc_
75 # define gtrigi gtrigi_
77 # define gzinit gzinit_
78 # define gfmate gfmate_
79 # define gfpart gfpart_
80 # define gftmed gftmed_
81 # define gftmat gftmat_
85 # define gsmate gsmate_
86 # define gsmixt gsmixt_
87 # define gspart gspart_
88 # define gstmed gstmed_
89 # define gsckov gsckov_
90 # define gstpar gstpar_
91 # define gfkine gfkine_
92 # define gfvert gfvert_
93 # define gskine gskine_
94 # define gsvert gsvert_
95 # define gphysi gphysi_
96 # define gdebug gdebug_
97 # define gekbin gekbin_
98 # define gfinds gfinds_
99 # define gsking gsking_
100 # define gskpho gskpho_
101 # define gsstak gsstak_
102 # define gsxyz gsxyz_
103 # define gtrack gtrack_
104 # define gtreve gtreve_
105 # define gtreve_root gtreve_root_
106 # define grndm grndm_
107 # define grndmq grndmq_
108 # define gdtom gdtom_
109 # define glmoth glmoth_
110 # define gmedia gmedia_
111 # define gmtod gmtod_
112 # define gsdvn gsdvn_
113 # define gsdvn2 gsdvn2_
114 # define gsdvs gsdvs_
115 # define gsdvs2 gsdvs2_
116 # define gsdvt gsdvt_
117 # define gsdvt2 gsdvt2_
118 # define gsord gsord_
119 # define gspos gspos_
120 # define gsposp gsposp_
121 # define gsrotm gsrotm_
122 # define gprotm gprotm_
123 # define gsvolu gsvolu_
124 # define gprint gprint_
125 # define gdinit gdinit_
126 # define gdopt gdopt_
127 # define gdraw gdraw_
128 # define gdrayt gdrayt_
129 # define gdrawc gdrawc_
130 # define gdrawx gdrawx_
131 # define gdhead gdhead_
132 # define gdwmn1 gdwmn1_
133 # define gdwmn2 gdwmn2_
134 # define gdwmn3 gdwmn3_
135 # define gdxyz gdxyz_
136 # define gdcxyz gdcxyz_
137 # define gdman gdman_
138 # define gdspec gdspec_
139 # define gdtree gdtree_
140 # define gdelet gdelet_
141 # define gdclos gdclos_
142 # define gdshow gdshow_
143 # define gdopen gdopen_
144 # define dzshow dzshow_
145 # define gsatt gsatt_
146 # define gfpara gfpara_
147 # define gckpar gckpar_
148 # define gckmat gckmat_
149 # define geditv geditv_
150 # define mzdrop mzdrop_
152 # define ertrak ertrak_
153 # define ertrgo ertrgo_
155 # define setbomb setbomb_
156 # define setclip setclip_
157 # define gcomad gcomad_
159 # define gbrelm gbrelm_
160 # define gprelm gprelm_
162 # define gzebra GZEBRA
163 # define grfile GRFILE
164 # define gpcxyz GPCXYZ
165 # define ggclos GGCLOS
168 # define gcinit GCINIT
171 # define gtrigc GTRIGC
172 # define gtrigi GTRIGI
174 # define gzinit GZINIT
175 # define gfmate GFMATE
176 # define gfpart GFPART
177 # define gftmed GFTMED
178 # define gftmat GFTMAT
182 # define gsmate GSMATE
183 # define gsmixt GSMIXT
184 # define gspart GSPART
185 # define gstmed GSTMED
186 # define gsckov GSCKOV
187 # define gstpar GSTPAR
188 # define gfkine GFKINE
189 # define gfvert GFVERT
190 # define gskine GSKINE
191 # define gsvert GSVERT
192 # define gphysi GPHYSI
193 # define gdebug GDEBUG
194 # define gekbin GEKBIN
195 # define gfinds GFINDS
196 # define gsking GSKING
197 # define gskpho GSKPHO
198 # define gsstak GSSTAK
200 # define gtrack GTRACK
201 # define gtreve GTREVE
202 # define gtreve_root GTREVE_ROOT
204 # define grndmq GRNDMQ
206 # define glmoth GLMOTH
207 # define gmedia GMEDIA
210 # define gsdvn2 GSDVN2
212 # define gsdvs2 GSDVS2
214 # define gsdvt2 GSDVT2
217 # define gsposp GSPOSP
218 # define gsrotm GSROTM
219 # define gprotm GPROTM
220 # define gsvolu GSVOLU
221 # define gprint GPRINT
222 # define gdinit GDINIT
225 # define gdrayt GDRAYT
226 # define gdrawc GDRAWC
227 # define gdrawx GDRAWX
228 # define gdhead GDHEAD
229 # define gdwmn1 GDWMN1
230 # define gdwmn2 GDWMN2
231 # define gdwmn3 GDWMN3
233 # define gdcxyz GDCXYZ
235 # define gdfspc GDFSPC
236 # define gdspec GDSPEC
237 # define gdtree GDTREE
238 # define gdelet GDELET
239 # define gdclos GDCLOS
240 # define gdshow GDSHOW
241 # define gdopen GDOPEN
242 # define dzshow DZSHOW
244 # define gfpara GFPARA
245 # define gckpar GCKPAR
246 # define gckmat GCKMAT
247 # define geditv GEDITV
248 # define mzdrop MZDROP
250 # define ertrak ERTRAK
251 # define ertrgo ERTRGO
253 # define setbomb SETBOMB
254 # define setclip SETCLIP
255 # define gcomad GCOMAD
257 # define gbrelm GBRELM
258 # define gprelm GPRELM
262 //____________________________________________________________________________
266 // Prototypes for GEANT functions
268 void type_of_call gzebra(const int&);
270 void type_of_call gpcxyz();
272 void type_of_call ggclos();
274 void type_of_call glast();
276 void type_of_call ginit();
278 void type_of_call gcinit();
280 void type_of_call grun();
282 void type_of_call gtrig();
284 void type_of_call gtrigc();
286 void type_of_call gtrigi();
288 void type_of_call gwork(const int&);
290 void type_of_call gzinit();
292 void type_of_call gmate();
294 void type_of_call gpart();
296 void type_of_call gsdk(Int_t &, Float_t *, Int_t *);
298 void type_of_call gfkine(Int_t &, Float_t *, Float_t *, Int_t &,
299 Int_t &, Float_t *, Int_t &);
301 void type_of_call gfvert(Int_t &, Float_t *, Int_t &, Int_t &,
302 Float_t &, Float_t *, Int_t &);
304 void type_of_call gskine(Float_t *,Int_t &, Int_t &, Float_t *,
307 void type_of_call gsvert(Float_t *,Int_t &, Int_t &, Float_t *,
310 void type_of_call gphysi();
312 void type_of_call gdebug();
314 void type_of_call gekbin();
316 void type_of_call gfinds();
318 void type_of_call gsking(Int_t &);
320 void type_of_call gskpho(Int_t &);
322 void type_of_call gsstak(Int_t &);
324 void type_of_call gsxyz();
326 void type_of_call gtrack();
328 void type_of_call gtreve();
330 void type_of_call gtreve_root();
332 void type_of_call grndm(Float_t *, const Int_t &);
334 void type_of_call grndmq(Int_t &, Int_t &, const Int_t &,
337 void type_of_call gdtom(Float_t *, Float_t *, Int_t &);
339 void type_of_call glmoth(DEFCHARD, Int_t &, Int_t &, Int_t *,
340 Int_t *, Int_t * DEFCHARL);
342 void type_of_call gmedia(Float_t *, Int_t &);
344 void type_of_call gmtod(Float_t *, Float_t *, Int_t &);
346 void type_of_call gsrotm(const Int_t &, const Float_t &, const Float_t &,
347 const Float_t &, const Float_t &, const Float_t &,
350 void type_of_call gprotm(const Int_t &);
352 void type_of_call grfile(const Int_t&, DEFCHARD,
353 DEFCHARD DEFCHARL DEFCHARL);
355 void type_of_call gfmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
356 Float_t &, Float_t &, Float_t &, Float_t *,
359 void type_of_call gfpart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
360 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
362 void type_of_call gftmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
363 Float_t &, Float_t &, Float_t &, Float_t &,
364 Float_t &, Float_t &, Float_t *, Int_t * DEFCHARL);
366 void type_of_call gftmat(const Int_t&, const Int_t&, DEFCHARD, const Int_t&,
368 ,Float_t *, Int_t & DEFCHARL);
370 void type_of_call gsmate(const Int_t&, DEFCHARD, Float_t &, Float_t &,
371 Float_t &, Float_t &, Float_t &, Float_t *,
374 void type_of_call gsmixt(const Int_t&, DEFCHARD, Float_t *, Float_t *,
375 Float_t &, Int_t &, Float_t * DEFCHARL);
377 void type_of_call gspart(const Int_t&, DEFCHARD, Int_t &, Float_t &,
378 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
381 void type_of_call gstmed(const Int_t&, DEFCHARD, Int_t &, Int_t &, Int_t &,
382 Float_t &, Float_t &, Float_t &, Float_t &,
383 Float_t &, Float_t &, Float_t *, Int_t & DEFCHARL);
385 void type_of_call gsckov(Int_t &itmed, Int_t &npckov, Float_t *ppckov,
386 Float_t *absco, Float_t *effic, Float_t *rindex);
387 void type_of_call gstpar(const Int_t&, DEFCHARD, Float_t & DEFCHARL);
389 void type_of_call gsdvn(DEFCHARD,DEFCHARD, Int_t &, Int_t &
392 void type_of_call gsdvn2(DEFCHARD,DEFCHARD, Int_t &, Int_t &, Float_t &,
393 Int_t & DEFCHARL DEFCHARL);
395 void type_of_call gsdvs(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &
398 void type_of_call gsdvs2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t &,
399 Int_t & DEFCHARL DEFCHARL);
401 void type_of_call gsdvt(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Int_t &,
402 Int_t & DEFCHARL DEFCHARL);
404 void type_of_call gsdvt2(DEFCHARD,DEFCHARD, Float_t &, Int_t &, Float_t&,
405 Int_t &, Int_t & DEFCHARL DEFCHARL);
407 void type_of_call gsord(DEFCHARD, Int_t & DEFCHARL);
409 void type_of_call gspos(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
410 Float_t &, Int_t &, DEFCHARD DEFCHARL DEFCHARL
413 void type_of_call gsposp(DEFCHARD, Int_t &, DEFCHARD, Float_t &, Float_t &,
414 Float_t &, Int_t &, DEFCHARD,
415 Float_t *, Int_t & DEFCHARL DEFCHARL DEFCHARL);
417 void type_of_call gsvolu(DEFCHARD, DEFCHARD, Int_t &, Float_t *, Int_t &,
418 Int_t & DEFCHARL DEFCHARL);
420 void type_of_call gsatt(DEFCHARD, DEFCHARD, Int_t & DEFCHARL DEFCHARL);
422 void type_of_call gfpara(DEFCHARD , Int_t&, Int_t&, Int_t&, Int_t&, Float_t*,
425 void type_of_call gckpar(Int_t&, Int_t&, Float_t*);
427 void type_of_call gckmat(Int_t&, DEFCHARD DEFCHARL);
429 void type_of_call gprint(DEFCHARD,const int& DEFCHARL);
431 void type_of_call gdinit();
433 void type_of_call gdopt(DEFCHARD,DEFCHARD DEFCHARL DEFCHARL);
435 void type_of_call gdraw(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
436 Float_t &, Float_t &, Float_t & DEFCHARL);
437 void type_of_call gdrayt(DEFCHARD,Float_t &,Float_t &, Float_t &,Float_t &,
438 Float_t &, Float_t &, Float_t & DEFCHARL);
439 void type_of_call gdrawc(DEFCHARD,Int_t &, Float_t &, Float_t &, Float_t &,
440 Float_t &, Float_t & DEFCHARL);
441 void type_of_call gdrawx(DEFCHARD,Float_t &, Float_t &, Float_t &, Float_t &,
442 Float_t &, Float_t &, Float_t &, Float_t &,
444 void type_of_call gdhead(Int_t &,DEFCHARD, Float_t & DEFCHARL);
445 void type_of_call gdxyz(Int_t &);
446 void type_of_call gdcxyz();
447 void type_of_call gdman(Float_t &, Float_t &);
448 void type_of_call gdwmn1(Float_t &, Float_t &);
449 void type_of_call gdwmn2(Float_t &, Float_t &);
450 void type_of_call gdwmn3(Float_t &, Float_t &);
451 void type_of_call gdspec(DEFCHARD DEFCHARL);
452 void type_of_call gdfspc(DEFCHARD, Int_t &, Int_t & DEFCHARL) {;}
453 void type_of_call gdtree(DEFCHARD, Int_t &, Int_t & DEFCHARL);
455 void type_of_call gdopen(Int_t &);
456 void type_of_call gdclos();
457 void type_of_call gdelet(Int_t &);
458 void type_of_call gdshow(Int_t &);
459 void type_of_call geditv(Int_t &) {;}
462 void type_of_call dzshow(DEFCHARD,const int&,const int&,DEFCHARD,const int&,
463 const int&, const int&, const int& DEFCHARL
466 void type_of_call mzdrop(Int_t&, Int_t&, DEFCHARD DEFCHARL);
468 void type_of_call setbomb(Float_t &);
469 void type_of_call setclip(DEFCHARD, Float_t &,Float_t &,Float_t &,Float_t &,
470 Float_t &, Float_t & DEFCHARL);
471 void type_of_call gcomad(DEFCHARD, Int_t*& DEFCHARL);
473 void type_of_call ertrak(const Float_t *const x1, const Float_t *const p1,
474 const Float_t *x2, const Float_t *p2,
475 const Int_t &ipa, DEFCHARD DEFCHARL);
477 void type_of_call ertrgo();
479 float type_of_call gbrelm(const Float_t &z, const Float_t& t, const Float_t& cut);
480 float type_of_call gprelm(const Float_t &z, const Float_t& t, const Float_t& cut);
484 // Geant3 global pointer
486 static Int_t defSize = 600;
490 //____________________________________________________________________________
494 // Default constructor
498 //____________________________________________________________________________
499 TGeant3::TGeant3(const char *title, Int_t nwgeant)
500 :AliMC("TGeant3",title)
503 // Standard constructor for TGeant3 with ZEBRA initialisation
514 // Load Address of Geant3 commons
517 // Zero number of particles
521 //____________________________________________________________________________
522 Int_t TGeant3::CurrentMaterial(Float_t &a, Float_t &z, Float_t &dens,
523 Float_t &radl, Float_t &absl) const
526 // Return the parameters of the current material during transport
530 dens = fGcmate->dens;
531 radl = fGcmate->radl;
532 absl = fGcmate->absl;
533 return 1; //this could be the number of elements in mixture
536 //____________________________________________________________________________
537 void TGeant3::DefaultRange()
540 // Set range of current drawing pad to 20x20 cm
546 higz->Range(0,0,20,20);
549 //____________________________________________________________________________
550 void TGeant3::InitHIGZ()
561 //____________________________________________________________________________
562 void TGeant3::LoadAddress()
565 // Assigns the address of the GEANT common blocks to the structures
566 // that allow their access from C++
569 gcomad(PASSCHARD("QUEST"), (int*&) fQuest PASSCHARL("QUEST"));
570 gcomad(PASSCHARD("GCBANK"),(int*&) fGcbank PASSCHARL("GCBANK"));
571 gcomad(PASSCHARD("GCLINK"),(int*&) fGclink PASSCHARL("GCLINK"));
572 gcomad(PASSCHARD("GCCUTS"),(int*&) fGccuts PASSCHARL("GCCUTS"));
573 gcomad(PASSCHARD("GCMULO"),(int*&) fGcmulo PASSCHARL("GCMULO"));
574 gcomad(PASSCHARD("GCFLAG"),(int*&) fGcflag PASSCHARL("GCFLAG"));
575 gcomad(PASSCHARD("GCKINE"),(int*&) fGckine PASSCHARL("GCKINE"));
576 gcomad(PASSCHARD("GCKING"),(int*&) fGcking PASSCHARL("GCKING"));
577 gcomad(PASSCHARD("GCKIN2"),(int*&) fGckin2 PASSCHARL("GCKIN2"));
578 gcomad(PASSCHARD("GCKIN3"),(int*&) fGckin3 PASSCHARL("GCKIN3"));
579 gcomad(PASSCHARD("GCMATE"),(int*&) fGcmate PASSCHARL("GCMATE"));
580 gcomad(PASSCHARD("GCTMED"),(int*&) fGctmed PASSCHARL("GCTMED"));
581 gcomad(PASSCHARD("GCTRAK"),(int*&) fGctrak PASSCHARL("GCTRAK"));
582 gcomad(PASSCHARD("GCTPOL"),(int*&) fGctpol PASSCHARL("GCTPOL"));
583 gcomad(PASSCHARD("GCVOLU"),(int*&) fGcvolu PASSCHARL("GCVOLU"));
584 gcomad(PASSCHARD("GCNUM"), (int*&) fGcnum PASSCHARL("GCNUM"));
585 gcomad(PASSCHARD("GCSETS"),(int*&) fGcsets PASSCHARL("GCSETS"));
586 gcomad(PASSCHARD("GCPHYS"),(int*&) fGcphys PASSCHARL("GCPHYS"));
587 gcomad(PASSCHARD("GCOPTI"),(int*&) fGcopti PASSCHARL("GCOPTI"));
588 gcomad(PASSCHARD("GCTLIT"),(int*&) fGctlit PASSCHARL("GCTLIT"));
589 gcomad(PASSCHARD("GCVDMA"),(int*&) fGcvdma PASSCHARL("GCVDMA"));
592 gcomad(PASSCHARD("ERTRIO"),(int*&) fErtrio PASSCHARL("ERTRIO"));
593 gcomad(PASSCHARD("EROPTS"),(int*&) fEropts PASSCHARL("EROPTS"));
594 gcomad(PASSCHARD("EROPTC"),(int*&) fEroptc PASSCHARL("EROPTC"));
595 gcomad(PASSCHARD("ERWORK"),(int*&) fErwork PASSCHARL("ERWORK"));
597 // Variables for ZEBRA store
598 gcomad(PASSCHARD("IQ"), addr PASSCHARL("IQ"));
600 gcomad(PASSCHARD("LQ"), addr PASSCHARL("LQ"));
605 //_____________________________________________________________________________
606 void TGeant3::GeomIter()
609 // Geometry iterator for moving upward in the geometry tree
610 // Initialise the iterator
612 fNextVol=fGcvolu->nlevel;
615 //____________________________________________________________________________
616 Int_t TGeant3::NextVolUp(Text_t *name, Int_t ©)
619 // Geometry iterator for moving upward in the geometry tree
620 // Return next volume up
625 gname=fGcvolu->names[fNextVol];
626 copy=fGcvolu->number[fNextVol];
627 i=fGcvolu->lvolum[fNextVol];
628 name = fVolNames[i-1];
629 if(gname == fZiq[fGclink->jvolum+i]) return i;
630 else printf("GeomTree: Volume %s not found in bank\n",name);
635 //_____________________________________________________________________________
636 Int_t TGeant3::CurrentVolID(Int_t ©) const
639 // Returns the current volume ID and copy number
642 if( (i=fGcvolu->nlevel-1) < 0 ) {
643 Warning("CurrentVolID","Stack depth only %d\n",fGcvolu->nlevel);
645 gname=fGcvolu->names[i];
646 copy=fGcvolu->number[i];
647 i=fGcvolu->lvolum[i];
648 if(gname == fZiq[fGclink->jvolum+i]) return i;
649 else Warning("CurrentVolID","Volume %4s not found\n",(char*)&gname);
654 //_____________________________________________________________________________
655 Int_t TGeant3::CurrentVolOffID(Int_t off, Int_t ©) const
658 // Return the current volume "off" upward in the geometrical tree
659 // ID and copy number
662 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
663 Warning("CurrentVolOffID","Offset requested %d but stack depth %d\n",
664 off,fGcvolu->nlevel);
666 gname=fGcvolu->names[i];
667 copy=fGcvolu->number[i];
668 i=fGcvolu->lvolum[i];
669 if(gname == fZiq[fGclink->jvolum+i]) return i;
670 else Warning("CurrentVolOffID","Volume %4s not found\n",(char*)&gname);
675 //_____________________________________________________________________________
676 const char* TGeant3::CurrentVolName() const
679 // Returns the current volume name
682 if( (i=fGcvolu->nlevel-1) < 0 ) {
683 Warning("CurrentVolName","Stack depth %d\n",fGcvolu->nlevel);
685 gname=fGcvolu->names[i];
686 i=fGcvolu->lvolum[i];
687 if(gname == fZiq[fGclink->jvolum+i]) return fVolNames[i-1];
688 else Warning("CurrentVolName","Volume %4s not found\n",(char*) &gname);
693 //_____________________________________________________________________________
694 const char* TGeant3::CurrentVolOffName(Int_t off) const
697 // Return the current volume "off" upward in the geometrical tree
698 // ID, name and copy number
699 // if name=0 no name is returned
702 if( (i=fGcvolu->nlevel-off-1) < 0 ) {
703 Warning("CurrentVolOffName",
704 "Offset requested %d but stack depth %d\n",off,fGcvolu->nlevel);
706 gname=fGcvolu->names[i];
707 i=fGcvolu->lvolum[i];
708 if(gname == fZiq[fGclink->jvolum+i]) return fVolNames[i-1];
709 else Warning("CurrentVolOffName","Volume %4s not found\n",(char*)&gname);
714 //_____________________________________________________________________________
715 Int_t TGeant3::IdFromPDG(Int_t pdg) const
718 // Return Geant3 code from PDG and pseudo ENDF code
720 for(Int_t i=0;i<fNPDGCodes;++i)
721 if(pdg==fPDGCode[i]) return i;
725 //_____________________________________________________________________________
726 Int_t TGeant3::PDGFromId(Int_t id) const
728 if(id>0 && id<fNPDGCodes) return fPDGCode[id];
732 //_____________________________________________________________________________
733 void TGeant3::DefineParticles()
736 // Define standard Geant 3 particles
739 // Load standard numbers for GEANT particles and PDG conversion
740 fPDGCode[fNPDGCodes++]=-99; // 0 = unused location
741 fPDGCode[fNPDGCodes++]=22; // 1 = photon
742 fPDGCode[fNPDGCodes++]=-11; // 2 = positron
743 fPDGCode[fNPDGCodes++]=11; // 3 = electron
744 fPDGCode[fNPDGCodes++]=12; // 4 = neutrino e
745 fPDGCode[fNPDGCodes++]=-13; // 5 = muon +
746 fPDGCode[fNPDGCodes++]=13; // 6 = muon -
747 fPDGCode[fNPDGCodes++]=111; // 7 = pi0
748 fPDGCode[fNPDGCodes++]=211; // 8 = pi+
749 fPDGCode[fNPDGCodes++]=-211; // 9 = pi-
750 fPDGCode[fNPDGCodes++]=130; // 10 = Kaon Long
751 fPDGCode[fNPDGCodes++]=321; // 11 = Kaon +
752 fPDGCode[fNPDGCodes++]=-321; // 12 = Kaon -
753 fPDGCode[fNPDGCodes++]=2112; // 13 = Neutron
754 fPDGCode[fNPDGCodes++]=2212; // 14 = Proton
755 fPDGCode[fNPDGCodes++]=-2212; // 15 = Anti Proton
756 fPDGCode[fNPDGCodes++]=310; // 16 = Kaon Short
757 fPDGCode[fNPDGCodes++]=221; // 17 = Eta
758 fPDGCode[fNPDGCodes++]=3122; // 18 = Lambda
759 fPDGCode[fNPDGCodes++]=3222; // 19 = Sigma +
760 fPDGCode[fNPDGCodes++]=3212; // 20 = Sigma 0
761 fPDGCode[fNPDGCodes++]=3112; // 21 = Sigma -
762 fPDGCode[fNPDGCodes++]=3322; // 22 = Xi0
763 fPDGCode[fNPDGCodes++]=3312; // 23 = Xi-
764 fPDGCode[fNPDGCodes++]=3334; // 24 = Omega-
765 fPDGCode[fNPDGCodes++]=-2112; // 25 = Anti Proton
766 fPDGCode[fNPDGCodes++]=-3122; // 26 = Anti Proton
767 fPDGCode[fNPDGCodes++]=-3222; // 27 = Anti Sigma -
768 fPDGCode[fNPDGCodes++]=-3212; // 28 = Anti Sigma 0
769 fPDGCode[fNPDGCodes++]=-3112; // 29 = Anti Sigma 0
770 fPDGCode[fNPDGCodes++]=-3322; // 30 = Anti Xi 0
771 fPDGCode[fNPDGCodes++]=-3312; // 31 = Anti Xi +
772 fPDGCode[fNPDGCodes++]=-3334; // 32 = Anti Omega +
779 /* --- Define additional particles */
780 Gspart(33, "OMEGA(782)", 3, 0.782, 0., 7.836e-23);
781 fPDGCode[fNPDGCodes++]=223; // 33 = Omega(782)
783 Gspart(34, "PHI(1020)", 3, 1.019, 0., 1.486e-22);
784 fPDGCode[fNPDGCodes++]=333; // 34 = PHI (1020)
786 Gspart(35, "D +", 4, 1.87, 1., 1.066e-12);
787 fPDGCode[fNPDGCodes++]=411; // 35 = D+
789 Gspart(36, "D -", 4, 1.87, -1., 1.066e-12);
790 fPDGCode[fNPDGCodes++]=-411; // 36 = D-
792 Gspart(37, "D 0", 3, 1.865, 0., 4.2e-13);
793 fPDGCode[fNPDGCodes++]=421; // 37 = D0
795 Gspart(38, "ANTI D 0", 3, 1.865, 0., 4.2e-13);
796 fPDGCode[fNPDGCodes++]=-421; // 38 = D0 bar
798 fPDGCode[fNPDGCodes++]=-99; // 39 = unassigned
800 fPDGCode[fNPDGCodes++]=-99; // 40 = unassigned
802 fPDGCode[fNPDGCodes++]=-99; // 41 = unassigned
804 Gspart(42, "RHO +", 4, 0.768, 1., 4.353e-24);
805 fPDGCode[fNPDGCodes++]=213; // 42 = RHO+
807 Gspart(43, "RHO -", 4, 0.768, -1., 4.353e-24);
808 fPDGCode[fNPDGCodes++]=-213; // 40 = RHO-
810 Gspart(44, "RHO 0", 3, 0.768, 0., 4.353e-24);
811 fPDGCode[fNPDGCodes++]=113; // 37 = D0
814 // Use ENDF-6 mapping for ions = 10000*z+10*a+iso
816 // and numbers above 5 000 000 for special applications
819 const Int_t kion=10000000;
821 const Int_t kspe=50000000;
823 TDatabasePDG *pdgDB = TDatabasePDG::Instance();
825 const Double_t autogev=0.9314943228;
826 const Double_t hslash = 1.0545726663e-27;
827 const Double_t erggev = 1/1.6021773349e-3;
828 const Double_t hshgev = hslash*erggev;
829 const Double_t yearstosec = 3600*24*365.25;
832 pdgDB->AddParticle("Deuteron","Deuteron",2*autogev+8.071e-3,kTRUE,
833 0,1,"Ion",kion+10020);
834 fPDGCode[fNPDGCodes++]=kion+10020; // 45 = Deuteron
836 pdgDB->AddParticle("Triton","Triton",3*autogev+14.931e-3,kFALSE,
837 hshgev/(12.33*yearstosec),1,"Ion",kion+10030);
838 fPDGCode[fNPDGCodes++]=kion+10030; // 46 = Triton
840 pdgDB->AddParticle("Alpha","Alpha",4*autogev+2.424e-3,kTRUE,
841 hshgev/(12.33*yearstosec),2,"Ion",kion+20040);
842 fPDGCode[fNPDGCodes++]=kion+20040; // 47 = Alpha
844 fPDGCode[fNPDGCodes++]=0; // 48 = geantino mapped to rootino
846 pdgDB->AddParticle("HE3","HE3",3*autogev+14.931e-3,kFALSE,
847 0,2,"Ion",kion+20030);
848 fPDGCode[fNPDGCodes++]=kion+20030; // 49 = HE3
850 pdgDB->AddParticle("Cherenkov","Cherenkov",0,kFALSE,
851 0,0,"Special",kspe+50);
852 fPDGCode[fNPDGCodes++]=kspe+50; // 50 = Cherenkov
854 Gspart(51, "FeedbackPhoton", 7, 0., 0.,1.e20 );
855 pdgDB->AddParticle("FeedbackPhoton","FeedbackPhoton",0,kFALSE,
856 0,0,"Special",kspe+51);
857 fPDGCode[fNPDGCodes++]=kspe+51; // 51 = FeedbackPhoton
859 /* --- Define additional decay modes --- */
860 /* --- omega(783) --- */
861 for (kz = 0; kz < 6; ++kz) {
872 Gsdk(ipa, bratio, mode);
873 /* --- phi(1020) --- */
874 for (kz = 0; kz < 6; ++kz) {
889 Gsdk(ipa, bratio, mode);
891 for (kz = 0; kz < 6; ++kz) {
904 Gsdk(ipa, bratio, mode);
906 for (kz = 0; kz < 6; ++kz) {
919 Gsdk(ipa, bratio, mode);
921 for (kz = 0; kz < 6; ++kz) {
932 Gsdk(ipa, bratio, mode);
933 /* --- Anti D0 --- */
934 for (kz = 0; kz < 6; ++kz) {
945 Gsdk(ipa, bratio, mode);
947 for (kz = 0; kz < 6; ++kz) {
954 Gsdk(ipa, bratio, mode);
956 for (kz = 0; kz < 6; ++kz) {
963 Gsdk(ipa, bratio, mode);
965 for (kz = 0; kz < 6; ++kz) {
972 Gsdk(ipa, bratio, mode);
975 for (kz = 0; kz < 6; ++kz) {
984 Gsdk(ipa, bratio, mode);
987 Gsdk(ipa, bratio, mode);
990 Gsdk(ipa, bratio, mode);
995 //_____________________________________________________________________________
996 Int_t TGeant3::VolId(Text_t *name) const
999 // Return the unique numeric identifier for volume name
1002 strncpy((char *) &gname, name, 4);
1003 for(i=1; i<=fGcnum->nvolum; i++)
1004 if(gname == fZiq[fGclink->jvolum+i]) return i;
1005 printf("VolId: Volume %s not found\n",name);
1009 //_____________________________________________________________________________
1010 Int_t TGeant3::NofVolumes() const
1013 // Return total number of volumes in the geometry
1015 return fGcnum->nvolum;
1018 //_____________________________________________________________________________
1019 const char* TGeant3::VolName(Int_t id) const
1022 // Return the volume name given the volume identifier
1024 const char name[5]="NULL";
1025 if(id<1 || id > fGcnum->nvolum || fGclink->jvolum<=0)
1028 return fVolNames[id-1];
1031 //_____________________________________________________________________________
1032 Float_t TGeant3::Xsec(char* reac, Float_t energy, Int_t part, Int_t mate)
1034 Int_t gpart = IdFromPDG(part);
1035 if(!strcmp(reac,"PHOT"))
1038 Error("Xsec","Can calculate photoelectric only for photons\n");
1044 //_____________________________________________________________________________
1045 void TGeant3::TrackPosition(TLorentzVector &xyz) const
1048 // Return the current position in the master reference frame of the
1049 // track being transported
1051 xyz[0]=fGctrak->vect[0];
1052 xyz[1]=fGctrak->vect[1];
1053 xyz[2]=fGctrak->vect[2];
1054 xyz[3]=fGctrak->tofg;
1057 //_____________________________________________________________________________
1058 Float_t TGeant3::TrackTime() const
1061 // Return the current time of flight of the track being transported
1063 return fGctrak->tofg;
1066 //_____________________________________________________________________________
1067 void TGeant3::TrackMomentum(TLorentzVector &xyz) const
1070 // Return the direction and the momentum (GeV/c) of the track
1071 // currently being transported
1073 Double_t ptot=fGctrak->vect[6];
1074 xyz[0]=fGctrak->vect[3]*ptot;
1075 xyz[1]=fGctrak->vect[4]*ptot;
1076 xyz[2]=fGctrak->vect[5]*ptot;
1077 xyz[3]=fGctrak->getot;
1080 //_____________________________________________________________________________
1081 Float_t TGeant3::TrackCharge() const
1084 // Return charge of the track currently transported
1086 return fGckine->charge;
1089 //_____________________________________________________________________________
1090 Float_t TGeant3::TrackMass() const
1093 // Return the mass of the track currently transported
1095 return fGckine->amass;
1098 //_____________________________________________________________________________
1099 Int_t TGeant3::TrackPid() const
1102 // Return the id of the particle transported
1104 return PDGFromId(fGckine->ipart);
1107 //_____________________________________________________________________________
1108 Float_t TGeant3::TrackStep() const
1111 // Return the length in centimeters of the current step
1113 return fGctrak->step;
1116 //_____________________________________________________________________________
1117 Float_t TGeant3::TrackLength() const
1120 // Return the length of the current track from its origin
1122 return fGctrak->sleng;
1125 //_____________________________________________________________________________
1126 Bool_t TGeant3::IsTrackInside() const
1129 // True if the track is not at the boundary of the current volume
1131 return (fGctrak->inwvol==0);
1134 //_____________________________________________________________________________
1135 Bool_t TGeant3::IsTrackEntering() const
1138 // True if this is the first step of the track in the current volume
1140 return (fGctrak->inwvol==1);
1143 //_____________________________________________________________________________
1144 Bool_t TGeant3::IsTrackExiting() const
1147 // True if this is the last step of the track in the current volume
1149 return (fGctrak->inwvol==2);
1152 //_____________________________________________________________________________
1153 Bool_t TGeant3::IsTrackOut() const
1156 // True if the track is out of the setup
1158 return (fGctrak->inwvol==3);
1161 //_____________________________________________________________________________
1162 Bool_t TGeant3::IsTrackStop() const
1165 // True if the track energy has fallen below the threshold
1167 return (fGctrak->istop==2);
1170 //_____________________________________________________________________________
1171 Int_t TGeant3::NSecondaries() const
1174 // Number of secondary particles generated in the current step
1176 return fGcking->ngkine;
1179 //_____________________________________________________________________________
1180 Int_t TGeant3::CurrentEvent() const
1183 // Number of the current event
1185 return fGcflag->idevt;
1188 //_____________________________________________________________________________
1189 const char* TGeant3::ProdProcess() const
1192 // Name of the process that has produced the secondary particles
1193 // in the current step
1195 static char proc[5];
1196 const Int_t ipmec[13] = { 5,6,7,8,9,10,11,12,21,23,25,105,108 };
1199 if(fGcking->ngkine>0) {
1200 for (km = 0; km < fGctrak->nmec; ++km) {
1201 for (im = 0; im < 13; ++im) {
1202 if (fGctrak->lmec[km] == ipmec[im]) {
1203 mec = fGctrak->lmec[km];
1204 if (0 < mec && mec < 31) {
1205 strncpy(proc,(char *)&fGctrak->namec[mec - 1],4);
1206 } else if (mec - 100 <= 30 && mec - 100 > 0) {
1207 strncpy(proc,(char *)&fGctpol->namec1[mec - 101],4);
1214 strcpy(proc,"UNKN");
1215 } else strcpy(proc,"NONE");
1219 //_____________________________________________________________________________
1220 void TGeant3::GetSecondary(Int_t isec, Int_t& ipart,
1221 TLorentzVector &x, TLorentzVector &p)
1224 // Get the parameters of the secondary track number isec produced
1225 // in the current step
1228 if(-1<isec && isec<fGcking->ngkine) {
1229 ipart=Int_t (fGcking->gkin[isec][4] +0.5);
1231 x[i]=fGckin3->gpos[isec][i];
1232 p[i]=fGcking->gkin[isec][i];
1234 x[3]=fGcking->tofd[isec];
1235 p[3]=fGcking->gkin[isec][3];
1237 printf(" * TGeant3::GetSecondary * Secondary %d does not exist\n",isec);
1238 x[0]=x[1]=x[2]=x[3]=p[0]=p[1]=p[2]=p[3]=0;
1243 //_____________________________________________________________________________
1244 void TGeant3::InitLego()
1247 SetDEBU(0,0,0); //do not print a message
1250 //_____________________________________________________________________________
1251 Bool_t TGeant3::IsTrackDisappeared() const
1254 // True if the current particle has disappered
1255 // either because it decayed or because it underwent
1256 // an inelastic collision
1258 return (fGctrak->istop==1);
1261 //_____________________________________________________________________________
1262 Bool_t TGeant3::IsTrackAlive() const
1265 // True if the current particle is alive and will continue to be
1268 return (fGctrak->istop==0);
1271 //_____________________________________________________________________________
1272 void TGeant3::StopTrack()
1275 // Stop the transport of the current particle and skip to the next
1280 //_____________________________________________________________________________
1281 void TGeant3::StopEvent()
1284 // Stop simulation of the current event and skip to the next
1289 //_____________________________________________________________________________
1290 Float_t TGeant3::MaxStep() const
1293 // Return the maximum step length in the current medium
1295 return fGctmed->stemax;
1298 //_____________________________________________________________________________
1299 void TGeant3::SetColors()
1302 // Set the colors for all the volumes
1303 // this is done sequentially for all volumes
1304 // based on the number of their medium
1307 Int_t jvolum=fGclink->jvolum;
1308 //Int_t jtmed=fGclink->jtmed;
1309 //Int_t jmate=fGclink->jmate;
1310 Int_t nvolum=fGcnum->nvolum;
1313 // Now for all the volumes
1314 for(kv=1;kv<=nvolum;kv++) {
1315 // Get the tracking medium
1316 Int_t itm=Int_t (fZq[fZlq[jvolum-kv]+4]);
1318 //Int_t ima=Int_t (fZq[fZlq[jtmed-itm]+6]);
1320 //Float_t z=fZq[fZlq[jmate-ima]+7];
1321 // Find color number
1322 //icol = Int_t(z)%6+2;
1323 //icol = 17+Int_t(z*150./92.);
1326 strncpy(name,(char*)&fZiq[jvolum+kv],4);
1328 Gsatt(name,"COLO",icol);
1332 //_____________________________________________________________________________
1333 void TGeant3::SetMaxStep(Float_t maxstep)
1336 // Set the maximum step allowed till the particle is in the current medium
1338 fGctmed->stemax=maxstep;
1341 //_____________________________________________________________________________
1342 void TGeant3::SetMaxNStep(Int_t maxnstp)
1345 // Set the maximum number of steps till the particle is in the current medium
1347 fGctrak->maxnst=maxnstp;
1350 //_____________________________________________________________________________
1351 Int_t TGeant3::GetMaxNStep() const
1354 // Maximum number of steps allowed in current medium
1356 return fGctrak->maxnst;
1359 //_____________________________________________________________________________
1360 void TGeant3::Material(Int_t& kmat, const char* name, Float_t a, Float_t z,
1361 Float_t dens, Float_t radl, Float_t absl, Float_t* buf,
1365 // Defines a Material
1367 // kmat number assigned to the material
1368 // name material name
1369 // a atomic mass in au
1371 // dens density in g/cm3
1372 // absl absorbtion length in cm
1373 // if >=0 it is ignored and the program
1374 // calculates it, if <0. -absl is taken
1375 // radl radiation length in cm
1376 // if >=0 it is ignored and the program
1377 // calculates it, if <0. -radl is taken
1378 // buf pointer to an array of user words
1379 // nbuf number of user words
1381 Int_t jmate=fGclink->jmate;
1387 for(i=1; i<=ns; i++) {
1388 if(fZlq[jmate-i]==0) {
1394 gsmate(kmat,PASSCHARD(name), a, z, dens, radl, absl, buf,
1395 nwbuf PASSCHARL(name));
1398 //_____________________________________________________________________________
1399 void TGeant3::Mixture(Int_t& kmat, const char* name, Float_t* a, Float_t* z,
1400 Float_t dens, Int_t nlmat, Float_t* wmat)
1403 // Defines mixture OR COMPOUND IMAT as composed by
1404 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
1406 // If NLMAT > 0 then wmat contains the proportion by
1407 // weights of each basic material in the mixture.
1409 // If nlmat < 0 then WMAT contains the number of atoms
1410 // of a given kind into the molecule of the COMPOUND
1411 // In this case, WMAT in output is changed to relative
1414 Int_t jmate=fGclink->jmate;
1420 for(i=1; i<=ns; i++) {
1421 if(fZlq[jmate-i]==0) {
1427 gsmixt(kmat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
1430 //_____________________________________________________________________________
1431 void TGeant3::Medium(Int_t& kmed, const char* name, Int_t nmat, Int_t isvol,
1432 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
1433 Float_t stemax, Float_t deemax, Float_t epsil,
1434 Float_t stmin, Float_t* ubuf, Int_t nbuf)
1437 // kmed tracking medium number assigned
1438 // name tracking medium name
1439 // nmat material number
1440 // isvol sensitive volume flag
1441 // ifield magnetic field
1442 // fieldm max. field value (kilogauss)
1443 // tmaxfd max. angle due to field (deg/step)
1444 // stemax max. step allowed
1445 // deemax max. fraction of energy lost in a step
1446 // epsil tracking precision (cm)
1447 // stmin min. step due to continuos processes (cm)
1449 // ifield = 0 if no magnetic field; ifield = -1 if user decision in guswim;
1450 // ifield = 1 if tracking performed with grkuta; ifield = 2 if tracking
1451 // performed with ghelix; ifield = 3 if tracking performed with ghelx3.
1453 Int_t jtmed=fGclink->jtmed;
1459 for(i=1; i<=ns; i++) {
1460 if(fZlq[jtmed-i]==0) {
1466 gstmed(kmed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1467 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1470 //_____________________________________________________________________________
1471 void TGeant3::Matrix(Int_t& krot, Float_t thex, Float_t phix, Float_t they,
1472 Float_t phiy, Float_t thez, Float_t phiz)
1475 // krot rotation matrix number assigned
1476 // theta1 polar angle for axis i
1477 // phi1 azimuthal angle for axis i
1478 // theta2 polar angle for axis ii
1479 // phi2 azimuthal angle for axis ii
1480 // theta3 polar angle for axis iii
1481 // phi3 azimuthal angle for axis iii
1483 // it defines the rotation matrix number irot.
1485 Int_t jrotm=fGclink->jrotm;
1491 for(i=1; i<=ns; i++) {
1492 if(fZlq[jrotm-i]==0) {
1498 gsrotm(krot, thex, phix, they, phiy, thez, phiz);
1501 //_____________________________________________________________________________
1502 Int_t TGeant3::GetMedium() const
1505 // Return the number of the current medium
1507 return fGctmed->numed;
1510 //_____________________________________________________________________________
1511 Float_t TGeant3::Edep() const
1514 // Return the energy lost in the current step
1516 return fGctrak->destep;
1519 //_____________________________________________________________________________
1520 Float_t TGeant3::Etot() const
1523 // Return the total energy of the current track
1525 return fGctrak->getot;
1528 //_____________________________________________________________________________
1529 void TGeant3::Rndm(Float_t* r, const Int_t n) const
1532 // Return an array of n random numbers uniformly distributed
1533 // between 0 and 1 not included
1538 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1540 // Functions from GBASE
1542 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1544 //____________________________________________________________________________
1545 void TGeant3::Gfile(const char *filename, const char *option)
1548 // Routine to open a GEANT/RZ data base.
1550 // LUN logical unit number associated to the file
1552 // CHFILE RZ file name
1554 // CHOPT is a character string which may be
1555 // N To create a new file
1556 // U to open an existing file for update
1557 // " " to open an existing file for read only
1558 // Q The initial allocation (default 1000 records)
1559 // is given in IQUEST(10)
1560 // X Open the file in exchange format
1561 // I Read all data structures from file to memory
1562 // O Write all data structures from memory to file
1565 // If options "I" or "O" all data structures are read or
1566 // written from/to file and the file is closed.
1567 // See routine GRMDIR to create subdirectories
1568 // See routines GROUT,GRIN to write,read objects
1570 grfile(21, PASSCHARD(filename), PASSCHARD(option) PASSCHARL(filename)
1574 //____________________________________________________________________________
1575 void TGeant3::Gpcxyz()
1578 // Print track and volume parameters at current point
1583 //_____________________________________________________________________________
1584 void TGeant3::Ggclos()
1587 // Closes off the geometry setting.
1588 // Initializes the search list for the contents of each
1589 // volume following the order they have been positioned, and
1590 // inserting the content '0' when a call to GSNEXT (-1) has
1591 // been required by the user.
1592 // Performs the development of the JVOLUM structure for all
1593 // volumes with variable parameters, by calling GGDVLP.
1594 // Interprets the user calls to GSORD, through GGORD.
1595 // Computes and stores in a bank (next to JVOLUM mother bank)
1596 // the number of levels in the geometrical tree and the
1597 // maximum number of contents per level, by calling GGNLEV.
1598 // Sets status bit for CONCAVE volumes, through GGCAVE.
1599 // Completes the JSET structure with the list of volume names
1600 // which identify uniquely a given physical detector, the
1601 // list of bit numbers to pack the corresponding volume copy
1602 // numbers, and the generic path(s) in the JVOLUM tree,
1603 // through the routine GHCLOS.
1606 // Create internal list of volumes
1607 fVolNames = new char[fGcnum->nvolum][5];
1609 for(i=0; i<fGcnum->nvolum; ++i) {
1610 strncpy(fVolNames[i], (char *) &fZiq[fGclink->jvolum+i+1], 4);
1611 fVolNames[i][4]='\0';
1615 //_____________________________________________________________________________
1616 void TGeant3::Glast()
1619 // Finish a Geant run
1624 //_____________________________________________________________________________
1625 void TGeant3::Gprint(const char *name)
1628 // Routine to print data structures
1629 // CHNAME name of a data structure
1633 gprint(PASSCHARD(vname),0 PASSCHARL(vname));
1636 //_____________________________________________________________________________
1637 void TGeant3::Grun()
1640 // Steering function to process one run
1645 //_____________________________________________________________________________
1646 void TGeant3::Gtrig()
1649 // Steering function to process one event
1654 //_____________________________________________________________________________
1655 void TGeant3::Gtrigc()
1658 // Clear event partition
1663 //_____________________________________________________________________________
1664 void TGeant3::Gtrigi()
1667 // Initialises event partition
1672 //_____________________________________________________________________________
1673 void TGeant3::Gwork(Int_t nwork)
1676 // Allocates workspace in ZEBRA memory
1681 //_____________________________________________________________________________
1682 void TGeant3::Gzinit()
1685 // To initialise GEANT/ZEBRA data structures
1690 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1692 // Functions from GCONS
1694 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1696 //_____________________________________________________________________________
1697 void TGeant3::Gfmate(Int_t imat, char *name, Float_t &a, Float_t &z,
1698 Float_t &dens, Float_t &radl, Float_t &absl,
1699 Float_t* ubuf, Int_t& nbuf)
1702 // Return parameters for material IMAT
1704 gfmate(imat, PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
1708 //_____________________________________________________________________________
1709 void TGeant3::Gfpart(Int_t ipart, char *name, Int_t &itrtyp,
1710 Float_t &amass, Float_t &charge, Float_t &tlife)
1713 // Return parameters for particle of type IPART
1717 Int_t igpart = IdFromPDG(ipart);
1718 gfpart(igpart, PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
1722 //_____________________________________________________________________________
1723 void TGeant3::Gftmed(Int_t numed, char *name, Int_t &nmat, Int_t &isvol,
1724 Int_t &ifield, Float_t &fieldm, Float_t &tmaxfd,
1725 Float_t &stemax, Float_t &deemax, Float_t &epsil,
1726 Float_t &stmin, Float_t *ubuf, Int_t *nbuf)
1729 // Return parameters for tracking medium NUMED
1731 gftmed(numed, PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1732 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1736 void TGeant3::Gftmat(Int_t imate, Int_t ipart, char *chmeca, Int_t kdim,
1737 Float_t* tkin, Float_t* value, Float_t* pcut,
1741 // Return parameters for tracking medium NUMED
1743 gftmat(imate, ipart, PASSCHARD(chmeca), kdim,
1744 tkin, value, pcut, ixst PASSCHARL(chmeca));
1748 Float_t TGeant3::Gbrelm(Float_t z, Float_t t, Float_t bcut)
1750 return gbrelm(z,t,bcut);
1753 Float_t TGeant3::Gprelm(Float_t z, Float_t t, Float_t bcut)
1755 return gprelm(z,t,bcut);
1758 //_____________________________________________________________________________
1759 void TGeant3::Gmate()
1762 // Define standard GEANT materials
1767 //_____________________________________________________________________________
1768 void TGeant3::Gpart()
1771 // Define standard GEANT particles plus selected decay modes
1772 // and branching ratios.
1777 //_____________________________________________________________________________
1778 void TGeant3::Gsdk(Int_t ipart, Float_t *bratio, Int_t *mode)
1780 // Defines branching ratios and decay modes for standard
1782 gsdk(ipart,bratio,mode);
1785 //_____________________________________________________________________________
1786 void TGeant3::Gsmate(Int_t imat, const char *name, Float_t a, Float_t z,
1787 Float_t dens, Float_t radl, Float_t absl)
1790 // Defines a Material
1792 // kmat number assigned to the material
1793 // name material name
1794 // a atomic mass in au
1796 // dens density in g/cm3
1797 // absl absorbtion length in cm
1798 // if >=0 it is ignored and the program
1799 // calculates it, if <0. -absl is taken
1800 // radl radiation length in cm
1801 // if >=0 it is ignored and the program
1802 // calculates it, if <0. -radl is taken
1803 // buf pointer to an array of user words
1804 // nbuf number of user words
1808 gsmate(imat,PASSCHARD(name), a, z, dens, radl, absl, ubuf, nbuf
1812 //_____________________________________________________________________________
1813 void TGeant3::Gsmixt(Int_t imat, const char *name, Float_t *a, Float_t *z,
1814 Float_t dens, Int_t nlmat, Float_t *wmat)
1817 // Defines mixture OR COMPOUND IMAT as composed by
1818 // THE BASIC NLMAT materials defined by arrays A,Z and WMAT
1820 // If NLMAT.GT.0 then WMAT contains the PROPORTION BY
1821 // WEIGTHS OF EACH BASIC MATERIAL IN THE MIXTURE.
1823 // If NLMAT.LT.0 then WMAT contains the number of atoms
1824 // of a given kind into the molecule of the COMPOUND
1825 // In this case, WMAT in output is changed to relative
1828 gsmixt(imat,PASSCHARD(name), a, z,dens, nlmat,wmat PASSCHARL(name));
1831 //_____________________________________________________________________________
1832 void TGeant3::Gspart(Int_t ipart, const char *name, Int_t itrtyp,
1833 Float_t amass, Float_t charge, Float_t tlife)
1836 // Store particle parameters
1838 // ipart particle code
1839 // name particle name
1840 // itrtyp transport method (see GEANT manual)
1841 // amass mass in GeV/c2
1842 // charge charge in electron units
1843 // tlife lifetime in seconds
1847 gspart(ipart,PASSCHARD(name), itrtyp, amass, charge, tlife, ubuf, nbuf
1851 //_____________________________________________________________________________
1852 void TGeant3::Gstmed(Int_t numed, const char *name, Int_t nmat, Int_t isvol,
1853 Int_t ifield, Float_t fieldm, Float_t tmaxfd,
1854 Float_t stemax, Float_t deemax, Float_t epsil,
1858 // NTMED Tracking medium number
1859 // NAME Tracking medium name
1860 // NMAT Material number
1861 // ISVOL Sensitive volume flag
1862 // IFIELD Magnetic field
1863 // FIELDM Max. field value (Kilogauss)
1864 // TMAXFD Max. angle due to field (deg/step)
1865 // STEMAX Max. step allowed
1866 // DEEMAX Max. fraction of energy lost in a step
1867 // EPSIL Tracking precision (cm)
1868 // STMIN Min. step due to continuos processes (cm)
1870 // IFIELD = 0 if no magnetic field; IFIELD = -1 if user decision in GUSWIM;
1871 // IFIELD = 1 if tracking performed with GRKUTA; IFIELD = 2 if tracking
1872 // performed with GHELIX; IFIELD = 3 if tracking performed with GHELX3.
1876 gstmed(numed,PASSCHARD(name), nmat, isvol, ifield, fieldm, tmaxfd, stemax,
1877 deemax, epsil, stmin, ubuf, nbuf PASSCHARL(name));
1880 //_____________________________________________________________________________
1881 void TGeant3::Gsckov(Int_t itmed, Int_t npckov, Float_t *ppckov,
1882 Float_t *absco, Float_t *effic, Float_t *rindex)
1885 // Stores the tables for UV photon tracking in medium ITMED
1886 // Please note that it is the user's responsability to
1887 // provide all the coefficients:
1890 // ITMED Tracking medium number
1891 // NPCKOV Number of bins of each table
1892 // PPCKOV Value of photon momentum (in GeV)
1893 // ABSCO Absorbtion coefficients
1894 // dielectric: absorbtion length in cm
1895 // metals : absorbtion fraction (0<=x<=1)
1896 // EFFIC Detection efficiency for UV photons
1897 // RINDEX Refraction index (if=0 metal)
1899 gsckov(itmed,npckov,ppckov,absco,effic,rindex);
1902 //_____________________________________________________________________________
1903 void TGeant3::Gstpar(Int_t itmed, const char *param, Float_t parval)
1906 // To change the value of cut or mechanism "CHPAR"
1907 // to a new value PARVAL for tracking medium ITMED
1908 // The data structure JTMED contains the standard tracking
1909 // parameters (CUTS and flags to control the physics processes) which
1910 // are used by default for all tracking media. It is possible to
1911 // redefine individually with GSTPAR any of these parameters for a
1912 // given tracking medium.
1913 // ITMED tracking medium number
1914 // CHPAR is a character string (variable name)
1915 // PARVAL must be given as a floating point.
1917 gstpar(itmed,PASSCHARD(param), parval PASSCHARL(param));
1920 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1922 // Functions from GCONS
1924 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
1926 //_____________________________________________________________________________
1927 void TGeant3::Gfkine(Int_t itra, Float_t *vert, Float_t *pvert, Int_t &ipart,
1930 // Storing/Retrieving Vertex and Track parameters
1931 // ----------------------------------------------
1933 // Stores vertex parameters.
1934 // VERT array of (x,y,z) position of the vertex
1935 // NTBEAM beam track number origin of the vertex
1936 // =0 if none exists
1937 // NTTARG target track number origin of the vertex
1938 // UBUF user array of NUBUF floating point numbers
1940 // NVTX new vertex number (=0 in case of error).
1941 // Prints vertex parameters.
1942 // IVTX for vertex IVTX.
1943 // (For all vertices if IVTX=0)
1944 // Stores long life track parameters.
1945 // PLAB components of momentum
1946 // IPART type of particle (see GSPART)
1947 // NV vertex number origin of track
1948 // UBUF array of NUBUF floating point user parameters
1950 // NT track number (if=0 error).
1951 // Retrieves long life track parameters.
1952 // ITRA track number for which parameters are requested
1953 // VERT vector origin of the track
1954 // PVERT 4 momentum components at the track origin
1955 // IPART particle type (=0 if track ITRA does not exist)
1956 // NVERT vertex number origin of the track
1957 // UBUF user words stored in GSKINE.
1958 // Prints initial track parameters.
1959 // ITRA for track ITRA
1960 // (For all tracks if ITRA=0)
1964 gfkine(itra,vert,pvert,ipart,nvert,ubuf,nbuf);
1967 //_____________________________________________________________________________
1968 void TGeant3::Gfvert(Int_t nvtx, Float_t *v, Int_t &ntbeam, Int_t &nttarg,
1972 // Retrieves the parameter of a vertex bank
1973 // Vertex is generated from tracks NTBEAM NTTARG
1974 // NVTX is the new vertex number
1978 gfvert(nvtx,v,ntbeam,nttarg,tofg,ubuf,nbuf);
1981 //_____________________________________________________________________________
1982 Int_t TGeant3::Gskine(Float_t *plab, Int_t ipart, Int_t nv, Float_t *buf,
1986 // Store kinematics of track NT into data structure
1987 // Track is coming from vertex NV
1990 gskine(plab, ipart, nv, buf, nwbuf, nt);
1994 //_____________________________________________________________________________
1995 Int_t TGeant3::Gsvert(Float_t *v, Int_t ntbeam, Int_t nttarg, Float_t *ubuf,
1999 // Creates a new vertex bank
2000 // Vertex is generated from tracks NTBEAM NTTARG
2001 // NVTX is the new vertex number
2004 gsvert(v, ntbeam, nttarg, ubuf, nwbuf, nwtx);
2008 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2010 // Functions from GPHYS
2012 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2014 //_____________________________________________________________________________
2015 void TGeant3::Gphysi()
2018 // Initialise material constants for all the physics
2019 // mechanisms used by GEANT
2024 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2026 // Functions from GTRAK
2028 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2030 //_____________________________________________________________________________
2031 void TGeant3::Gdebug()
2034 // Debug the current step
2039 //_____________________________________________________________________________
2040 void TGeant3::Gekbin()
2043 // To find bin number in kinetic energy table
2044 // stored in ELOW(NEKBIN)
2049 //_____________________________________________________________________________
2050 void TGeant3::Gfinds()
2053 // Returns the set/volume parameters corresponding to
2054 // the current space point in /GCTRAK/
2055 // and fill common /GCSETS/
2057 // IHSET user set identifier
2058 // IHDET user detector identifier
2059 // ISET set number in JSET
2060 // IDET detector number in JS=LQ(JSET-ISET)
2061 // IDTYPE detector type (1,2)
2062 // NUMBV detector volume numbers (array of length NVNAME)
2063 // NVNAME number of volume levels
2068 //_____________________________________________________________________________
2069 void TGeant3::Gsking(Int_t igk)
2072 // Stores in stack JSTAK either the IGKth track of /GCKING/,
2073 // or the NGKINE tracks when IGK is 0.
2078 //_____________________________________________________________________________
2079 void TGeant3::Gskpho(Int_t igk)
2082 // Stores in stack JSTAK either the IGKth Cherenkov photon of
2083 // /GCKIN2/, or the NPHOT tracks when IGK is 0.
2088 //_____________________________________________________________________________
2089 void TGeant3::Gsstak(Int_t iflag)
2092 // Stores in auxiliary stack JSTAK the particle currently
2093 // described in common /GCKINE/.
2095 // On request, creates also an entry in structure JKINE :
2097 // 0 : No entry in JKINE structure required (user)
2098 // 1 : New entry in JVERTX / JKINE structures required (user)
2099 // <0 : New entry in JKINE structure at vertex -IFLAG (user)
2100 // 2 : Entry in JKINE structure exists already (from GTREVE)
2105 //_____________________________________________________________________________
2106 void TGeant3::Gsxyz()
2109 // Store space point VECT in banks JXYZ
2114 //_____________________________________________________________________________
2115 void TGeant3::Gtrack()
2118 // Controls tracking of current particle
2123 //_____________________________________________________________________________
2124 void TGeant3::Gtreve()
2127 // Controls tracking of all particles belonging to the current event
2132 //_____________________________________________________________________________
2133 void TGeant3::Gtreve_root()
2136 // Controls tracking of all particles belonging to the current event
2141 //_____________________________________________________________________________
2142 void TGeant3::Grndm(Float_t *rvec, const Int_t len) const
2145 // To generate a vector RVECV of LEN random numbers
2146 // Copy of the CERN Library routine RANECU
2150 //_____________________________________________________________________________
2151 void TGeant3::Grndmq(Int_t &is1, Int_t &is2, const Int_t iseq,
2152 const Text_t *chopt)
2155 // To set/retrieve the seed of the random number generator
2157 grndmq(is1,is2,iseq,PASSCHARD(chopt) PASSCHARL(chopt));
2160 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2162 // Functions from GDRAW
2164 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2166 //_____________________________________________________________________________
2167 void TGeant3::Gdxyz(Int_t it)
2170 // Draw the points stored with Gsxyz relative to track it
2175 //_____________________________________________________________________________
2176 void TGeant3::Gdcxyz()
2179 // Draw the position of the current track
2184 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2186 // Functions from GGEOM
2188 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2190 //_____________________________________________________________________________
2191 void TGeant3::Gdtom(Float_t *xd, Float_t *xm, Int_t iflag)
2194 // Computes coordinates XM (Master Reference System
2195 // knowing the coordinates XD (Detector Ref System)
2196 // The local reference system can be initialized by
2197 // - the tracking routines and GDTOM used in GUSTEP
2198 // - a call to GSCMED(NLEVEL,NAMES,NUMBER)
2199 // (inverse routine is GMTOD)
2201 // If IFLAG=1 convert coordinates
2202 // IFLAG=2 convert direction cosinus
2204 gdtom(xd, xm, iflag);
2207 //_____________________________________________________________________________
2208 void TGeant3::Glmoth(const char* iudet, Int_t iunum, Int_t &nlev, Int_t *lvols,
2212 // Loads the top part of the Volume tree in LVOLS (IVO's),
2213 // LINDX (IN indices) for a given volume defined through
2214 // its name IUDET and number IUNUM.
2216 // The routine stores only upto the last level where JVOLUM
2217 // data structure is developed. If there is no development
2218 // above the current level, it returns NLEV zero.
2220 glmoth(PASSCHARD(iudet), iunum, nlev, lvols, lindx, idum PASSCHARL(iudet));
2223 //_____________________________________________________________________________
2224 void TGeant3::Gmedia(Float_t *x, Int_t &numed)
2227 // Finds in which volume/medium the point X is, and updates the
2228 // common /GCVOLU/ and the structure JGPAR accordingly.
2230 // NUMED returns the tracking medium number, or 0 if point is
2231 // outside the experimental setup.
2236 //_____________________________________________________________________________
2237 void TGeant3::Gmtod(Float_t *xm, Float_t *xd, Int_t iflag)
2240 // Computes coordinates XD (in DRS)
2241 // from known coordinates XM in MRS
2242 // The local reference system can be initialized by
2243 // - the tracking routines and GMTOD used in GUSTEP
2244 // - a call to GMEDIA(XM,NUMED)
2245 // - a call to GLVOLU(NLEVEL,NAMES,NUMBER,IER)
2246 // (inverse routine is GDTOM)
2248 // If IFLAG=1 convert coordinates
2249 // IFLAG=2 convert direction cosinus
2251 gmtod(xm, xd, iflag);
2254 //_____________________________________________________________________________
2255 void TGeant3::Gsdvn(const char *name, const char *mother, Int_t ndiv,
2259 // Create a new volume by dividing an existing one
2262 // MOTHER Mother volume name
2263 // NDIV Number of divisions
2266 // X,Y,Z of CAXIS will be translated to 1,2,3 for IAXIS.
2267 // It divides a previously defined volume.
2272 Vname(mother,vmother);
2273 gsdvn(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis PASSCHARL(vname)
2274 PASSCHARL(vmother));
2277 //_____________________________________________________________________________
2278 void TGeant3::Gsdvn2(const char *name, const char *mother, Int_t ndiv,
2279 Int_t iaxis, Float_t c0i, Int_t numed)
2282 // Create a new volume by dividing an existing one
2284 // Divides mother into ndiv divisions called name
2285 // along axis iaxis starting at coordinate value c0.
2286 // the new volume created will be medium number numed.
2291 Vname(mother,vmother);
2292 gsdvn2(PASSCHARD(vname), PASSCHARD(vmother), ndiv, iaxis, c0i, numed
2293 PASSCHARL(vname) PASSCHARL(vmother));
2296 //_____________________________________________________________________________
2297 void TGeant3::Gsdvs(const char *name, const char *mother, Float_t step,
2298 Int_t iaxis, Int_t numed)
2301 // Create a new volume by dividing an existing one
2306 Vname(mother,vmother);
2307 gsdvs(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed
2308 PASSCHARL(vname) PASSCHARL(vmother));
2311 //_____________________________________________________________________________
2312 void TGeant3::Gsdvs2(const char *name, const char *mother, Float_t step,
2313 Int_t iaxis, Float_t c0, Int_t numed)
2316 // Create a new volume by dividing an existing one
2321 Vname(mother,vmother);
2322 gsdvs2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0, numed
2323 PASSCHARL(vname) PASSCHARL(vmother));
2326 //_____________________________________________________________________________
2327 void TGeant3::Gsdvt(const char *name, const char *mother, Float_t step,
2328 Int_t iaxis, Int_t numed, Int_t ndvmx)
2331 // Create a new volume by dividing an existing one
2333 // Divides MOTHER into divisions called NAME along
2334 // axis IAXIS in steps of STEP. If not exactly divisible
2335 // will make as many as possible and will centre them
2336 // with respect to the mother. Divisions will have medium
2337 // number NUMED. If NUMED is 0, NUMED of MOTHER is taken.
2338 // NDVMX is the expected maximum number of divisions
2339 // (If 0, no protection tests are performed)
2344 Vname(mother,vmother);
2345 gsdvt(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, numed, ndvmx
2346 PASSCHARL(vname) PASSCHARL(vmother));
2349 //_____________________________________________________________________________
2350 void TGeant3::Gsdvt2(const char *name, const char *mother, Float_t step,
2351 Int_t iaxis, Float_t c0, Int_t numed, Int_t ndvmx)
2354 // Create a new volume by dividing an existing one
2356 // Divides MOTHER into divisions called NAME along
2357 // axis IAXIS starting at coordinate value C0 with step
2359 // The new volume created will have medium number NUMED.
2360 // If NUMED is 0, NUMED of mother is taken.
2361 // NDVMX is the expected maximum number of divisions
2362 // (If 0, no protection tests are performed)
2367 Vname(mother,vmother);
2368 gsdvt2(PASSCHARD(vname), PASSCHARD(vmother), step, iaxis, c0,
2369 numed, ndvmx PASSCHARL(vname) PASSCHARL(vmother));
2372 //_____________________________________________________________________________
2373 void TGeant3::Gsord(const char *name, Int_t iax)
2376 // Flags volume CHNAME whose contents will have to be ordered
2377 // along axis IAX, by setting the search flag to -IAX
2381 // IAX = 4 Rxy (static ordering only -> GTMEDI)
2382 // IAX = 14 Rxy (also dynamic ordering -> GTNEXT)
2383 // IAX = 5 Rxyz (static ordering only -> GTMEDI)
2384 // IAX = 15 Rxyz (also dynamic ordering -> GTNEXT)
2385 // IAX = 6 PHI (PHI=0 => X axis)
2386 // IAX = 7 THETA (THETA=0 => Z axis)
2390 gsord(PASSCHARD(vname), iax PASSCHARL(vname));
2393 //_____________________________________________________________________________
2394 void TGeant3::Gspos(const char *name, Int_t nr, const char *mother, Float_t x,
2395 Float_t y, Float_t z, Int_t irot, const char *konly)
2398 // Position a volume into an existing one
2401 // NUMBER Copy number of the volume
2402 // MOTHER Mother volume name
2403 // X X coord. of the volume in mother ref. sys.
2404 // Y Y coord. of the volume in mother ref. sys.
2405 // Z Z coord. of the volume in mother ref. sys.
2406 // IROT Rotation matrix number w.r.t. mother ref. sys.
2407 // ONLY ONLY/MANY flag
2409 // It positions a previously defined volume in the mother.
2414 Vname(mother,vmother);
2415 gspos(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2416 PASSCHARD(konly) PASSCHARL(vname) PASSCHARL(vmother)
2420 //_____________________________________________________________________________
2421 void TGeant3::Gsposp(const char *name, Int_t nr, const char *mother,
2422 Float_t x, Float_t y, Float_t z, Int_t irot,
2423 const char *konly, Float_t *upar, Int_t np )
2426 // Place a copy of generic volume NAME with user number
2427 // NR inside MOTHER, with its parameters UPAR(1..NP)
2432 Vname(mother,vmother);
2433 gsposp(PASSCHARD(vname), nr, PASSCHARD(vmother), x, y, z, irot,
2434 PASSCHARD(konly), upar, np PASSCHARL(vname) PASSCHARL(vmother)
2438 //_____________________________________________________________________________
2439 void TGeant3::Gsrotm(Int_t nmat, Float_t theta1, Float_t phi1, Float_t theta2,
2440 Float_t phi2, Float_t theta3, Float_t phi3)
2443 // nmat Rotation matrix number
2444 // THETA1 Polar angle for axis I
2445 // PHI1 Azimuthal angle for axis I
2446 // THETA2 Polar angle for axis II
2447 // PHI2 Azimuthal angle for axis II
2448 // THETA3 Polar angle for axis III
2449 // PHI3 Azimuthal angle for axis III
2451 // It defines the rotation matrix number IROT.
2453 gsrotm(nmat, theta1, phi1, theta2, phi2, theta3, phi3);
2456 //_____________________________________________________________________________
2457 void TGeant3::Gprotm(Int_t nmat)
2460 // To print rotation matrices structure JROTM
2461 // nmat Rotation matrix number
2466 //_____________________________________________________________________________
2467 Int_t TGeant3::Gsvolu(const char *name, const char *shape, Int_t nmed,
2468 Float_t *upar, Int_t npar)
2472 // SHAPE Volume type
2473 // NUMED Tracking medium number
2474 // NPAR Number of shape parameters
2475 // UPAR Vector containing shape parameters
2477 // It creates a new volume in the JVOLUM data structure.
2483 Vname(shape,vshape);
2484 gsvolu(PASSCHARD(vname), PASSCHARD(vshape), nmed, upar, npar, ivolu
2485 PASSCHARL(vname) PASSCHARL(vshape));
2489 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2491 // T H E D R A W I N G P A C K A G E
2492 // ======================================
2493 // Drawing functions. These functions allow the visualization in several ways
2494 // of the volumes defined in the geometrical data structure. It is possible
2495 // to draw the logical tree of volumes belonging to the detector (DTREE),
2496 // to show their geometrical specification (DSPEC,DFSPC), to draw them
2497 // and their cut views (DRAW, DCUT). Moreover, it is possible to execute
2498 // these commands when the hidden line removal option is activated; in
2499 // this case, the volumes can be also either translated in the space
2500 // (SHIFT), or clipped by boolean operation (CVOL). In addition, it is
2501 // possible to fill the surfaces of the volumes
2502 // with solid colours when the shading option (SHAD) is activated.
2503 // Several tools (ZOOM, LENS) have been developed to zoom detailed parts
2504 // of the detectors or to scan physical events as well.
2505 // Finally, the command MOVE will allow the rotation, translation and zooming
2506 // on real time parts of the detectors or tracks and hits of a simulated event.
2507 // Ray-tracing commands. In case the command (DOPT RAYT ON) is executed,
2508 // the drawing is performed by the Geant ray-tracing;
2509 // automatically, the color is assigned according to the tracking medium of each
2510 // volume and the volumes with a density lower/equal than the air are considered
2511 // transparent; if the option (USER) is set (ON) (again via the command (DOPT)),
2512 // the user can set color and visibility for the desired volumes via the command
2513 // (SATT), as usual, relatively to the attributes (COLO) and (SEEN).
2514 // The resolution can be set via the command (SATT * FILL VALUE), where (VALUE)
2515 // is the ratio between the number of pixels drawn and 20 (user coordinates).
2516 // Parallel view and perspective view are possible (DOPT PROJ PARA/PERS); in the
2517 // first case, we assume that the first mother volume of the tree is a box with
2518 // dimensions 10000 X 10000 X 10000 cm and the view point (infinetely far) is
2519 // 5000 cm far from the origin along the Z axis of the user coordinates; in the
2520 // second case, the distance between the observer and the origin of the world
2521 // reference system is set in cm by the command (PERSP NAME VALUE); grand-angle
2522 // or telescopic effects can be achieved changing the scale factors in the command
2523 // (DRAW). When the final picture does not occupy the full window,
2524 // mapping the space before tracing can speed up the drawing, but can also
2525 // produce less precise results; values from 1 to 4 are allowed in the command
2526 // (DOPT MAPP VALUE), the mapping being more precise for increasing (VALUE); for
2527 // (VALUE = 0) no mapping is performed (therefore max precision and lowest speed).
2528 // The command (VALCUT) allows the cutting of the detector by three planes
2529 // ortogonal to the x,y,z axis. The attribute (LSTY) can be set by the command
2530 // SATT for any desired volume and can assume values from 0 to 7; it determines
2531 // the different light processing to be performed for different materials:
2532 // 0 = dark-matt, 1 = bright-matt, 2 = plastic, 3 = ceramic, 4 = rough-metals,
2533 // 5 = shiny-metals, 6 = glass, 7 = mirror. The detector is assumed to be in the
2534 // dark, the ambient light luminosity is 0.2 for each basic hue (the saturation
2535 // is 0.9) and the observer is assumed to have a light source (therefore he will
2536 // produce parallel light in the case of parallel view and point-like-source
2537 // light in the case of perspective view).
2539 //*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
2541 //_____________________________________________________________________________
2542 void TGeant3::Gsatt(const char *name, const char *att, Int_t val)
2546 // IOPT Name of the attribute to be set
2547 // IVAL Value to which the attribute is to be set
2549 // name= "*" stands for all the volumes.
2550 // iopt can be chosen among the following :
2552 // WORK 0=volume name is inactive for the tracking
2553 // 1=volume name is active for the tracking (default)
2555 // SEEN 0=volume name is invisible
2556 // 1=volume name is visible (default)
2557 // -1=volume invisible with all its descendants in the tree
2558 // -2=volume visible but not its descendants in the tree
2560 // LSTY line style 1,2,3,... (default=1)
2561 // LSTY=7 will produce a very precise approximation for
2562 // revolution bodies.
2564 // LWID line width -7,...,1,2,3,..7 (default=1)
2565 // LWID<0 will act as abs(LWID) was set for the volume
2566 // and for all the levels below it. When SHAD is 'ON', LWID
2567 // represent the linewidth of the scan lines filling the surfaces
2568 // (whereas the FILL value represent their number). Therefore
2569 // tuning this parameter will help to obtain the desired
2570 // quality/performance ratio.
2572 // COLO colour code -166,...,1,2,..166 (default=1)
2574 // n=2=red; n=17+m, m=0,25, increasing luminosity according to 'm';
2575 // n=3=green; n=67+m, m=0,25, increasing luminosity according to 'm';
2576 // n=4=blue; n=117+m, m=0,25, increasing luminosity according to 'm';
2577 // n=5=yellow; n=42+m, m=0,25, increasing luminosity according to 'm';
2578 // n=6=violet; n=142+m, m=0,25, increasing luminosity according to 'm';
2579 // n=7=lightblue; n=92+m, m=0,25, increasing luminosity according to 'm';
2580 // colour=n*10+m, m=1,2,...9, will produce the same colour
2581 // as 'n', but with increasing luminosity according to 'm';
2582 // COLO<0 will act as if abs(COLO) was set for the volume
2583 // and for all the levels below it.
2584 // When for a volume the attribute FILL is > 1 (and the
2585 // option SHAD is on), the ABS of its colour code must be < 8
2586 // because an automatic shading of its faces will be
2589 // FILL (1992) fill area -7,...,0,1,...7 (default=0)
2590 // when option SHAD is "on" the FILL attribute of any
2591 // volume can be set different from 0 (normal drawing);
2592 // if it is set to 1, the faces of such volume will be filled
2593 // with solid colours; if ABS(FILL) is > 1, then a light
2594 // source is placed along the observer line, and the faces of
2595 // such volumes will be painted by colours whose luminosity
2596 // will depend on the amount of light reflected;
2597 // if ABS(FILL) = 1, then it is possible to use all the 166
2598 // colours of the colour table, becouse the automatic shading
2599 // is not performed;
2600 // for increasing values of FILL the drawing will be performed
2601 // with higher and higher resolution improving the quality (the
2602 // number of scan lines used to fill the faces increases with FILL);
2603 // it is possible to set different values of FILL
2604 // for different volumes, in order to optimize at the same time
2605 // the performance and the quality of the picture;
2606 // FILL<0 will act as if abs(FILL) was set for the volume
2607 // and for all the levels below it.
2608 // This kind of drawing can be saved in 'picture files'
2609 // or in view banks.
2610 // 0=drawing without fill area
2611 // 1=faces filled with solid colours and resolution = 6
2612 // 2=lowest resolution (very fast)
2613 // 3=default resolution
2614 // 4=.................
2615 // 5=.................
2616 // 6=.................
2618 // Finally, if a coloured background is desired, the FILL
2619 // attribute for the first volume of the tree must be set
2620 // equal to -abs(colo), colo being >0 and <166.
2622 // SET set number associated to volume name
2623 // DET detector number associated to volume name
2624 // DTYP detector type (1,2)
2631 gsatt(PASSCHARD(vname), PASSCHARD(vatt), val PASSCHARL(vname)
2635 //_____________________________________________________________________________
2636 void TGeant3::Gfpara(const char *name, Int_t number, Int_t intext, Int_t& npar,
2637 Int_t& natt, Float_t* par, Float_t* att)
2640 // Find the parameters of a volume
2642 gfpara(PASSCHARD(name), number, intext, npar, natt, par, att
2646 //_____________________________________________________________________________
2647 void TGeant3::Gckpar(Int_t ish, Int_t npar, Float_t* par)
2650 // Check the parameters of a shape
2652 gckpar(ish,npar,par);
2655 //_____________________________________________________________________________
2656 void TGeant3::Gckmat(Int_t itmed, char* natmed)
2659 // Check the parameters of a tracking medium
2661 gckmat(itmed, PASSCHARD(natmed) PASSCHARL(natmed));
2664 //_____________________________________________________________________________
2665 void TGeant3::Gdelete(Int_t iview)
2668 // IVIEW View number
2670 // It deletes a view bank from memory.
2675 //_____________________________________________________________________________
2676 void TGeant3::Gdopen(Int_t iview)
2679 // IVIEW View number
2681 // When a drawing is very complex and requires a long time to be
2682 // executed, it can be useful to store it in a view bank: after a
2683 // call to DOPEN and the execution of the drawing (nothing will
2684 // appear on the screen), and after a necessary call to DCLOSE,
2685 // the contents of the bank can be displayed in a very fast way
2686 // through a call to DSHOW; therefore, the detector can be easily
2687 // zoomed many times in different ways. Please note that the pictures
2688 // with solid colours can now be stored in a view bank or in 'PICTURE FILES'
2695 //_____________________________________________________________________________
2696 void TGeant3::Gdclose()
2699 // It closes the currently open view bank; it must be called after the
2700 // end of the drawing to be stored.
2705 //_____________________________________________________________________________
2706 void TGeant3::Gdshow(Int_t iview)
2709 // IVIEW View number
2711 // It shows on the screen the contents of a view bank. It
2712 // can be called after a view bank has been closed.
2717 //_____________________________________________________________________________
2718 void TGeant3::Gdopt(const char *name,const char *value)
2722 // VALUE Option value
2724 // To set/modify the drawing options.
2727 // THRZ ON Draw tracks in R vs Z
2728 // OFF (D) Draw tracks in X,Y,Z
2731 // PROJ PARA (D) Parallel projection
2733 // TRAK LINE (D) Trajectory drawn with lines
2734 // POIN " " with markers
2735 // HIDE ON Hidden line removal using the CG package
2736 // OFF (D) No hidden line removal
2737 // SHAD ON Fill area and shading of surfaces.
2738 // OFF (D) Normal hidden line removal.
2739 // RAYT ON Ray-tracing on.
2740 // OFF (D) Ray-tracing off.
2741 // EDGE OFF Does not draw contours when shad is on.
2742 // ON (D) Normal shading.
2743 // MAPP 1,2,3,4 Mapping before ray-tracing.
2744 // 0 (D) No mapping.
2745 // USER ON User graphics options in the raytracing.
2746 // OFF (D) Automatic graphics options.
2752 Vname(value,vvalue);
2753 gdopt(PASSCHARD(vname), PASSCHARD(vvalue) PASSCHARL(vname)
2757 //_____________________________________________________________________________
2758 void TGeant3::Gdraw(const char *name,Float_t theta, Float_t phi, Float_t psi,
2759 Float_t u0,Float_t v0,Float_t ul,Float_t vl)
2764 // THETA Viewing angle theta (for 3D projection)
2765 // PHI Viewing angle phi (for 3D projection)
2766 // PSI Viewing angle psi (for 2D rotation)
2767 // U0 U-coord. (horizontal) of volume origin
2768 // V0 V-coord. (vertical) of volume origin
2769 // SU Scale factor for U-coord.
2770 // SV Scale factor for V-coord.
2772 // This function will draw the volumes,
2773 // selected with their graphical attributes, set by the Gsatt
2774 // facility. The drawing may be performed with hidden line removal
2775 // and with shading effects according to the value of the options HIDE
2776 // and SHAD; if the option SHAD is ON, the contour's edges can be
2777 // drawn or not. If the option HIDE is ON, the detector can be
2778 // exploded (BOMB), clipped with different shapes (CVOL), and some
2779 // of its parts can be shifted from their original
2780 // position (SHIFT). When HIDE is ON, if
2781 // the drawing requires more than the available memory, the program
2782 // will evaluate and display the number of missing words
2783 // (so that the user can increase the
2784 // size of its ZEBRA store). Finally, at the end of each drawing (with HIDE on),
2785 // the program will print messages about the memory used and
2786 // statistics on the volumes' visibility.
2787 // The following commands will produce the drawing of a green
2788 // volume, specified by NAME, without using the hidden line removal
2789 // technique, using the hidden line removal technique,
2790 // with different linewidth and colour (red), with
2791 // solid colour, with shading of surfaces, and without edges.
2792 // Finally, some examples are given for the ray-tracing. (A possible
2793 // string for the NAME of the volume can be found using the command DTREE).
2799 if (fGcvdma->raytra != 1) {
2800 gdraw(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
2802 gdrayt(PASSCHARD(vname), theta,phi,psi,u0,v0,ul,vl PASSCHARL(vname));
2806 //_____________________________________________________________________________
2807 void TGeant3::Gdrawc(const char *name,Int_t axis, Float_t cut,Float_t u0,
2808 Float_t v0,Float_t ul,Float_t vl)
2813 // CUTVAL Cut plane distance from the origin along the axis
2815 // U0 U-coord. (horizontal) of volume origin
2816 // V0 V-coord. (vertical) of volume origin
2817 // SU Scale factor for U-coord.
2818 // SV Scale factor for V-coord.
2820 // The cut plane is normal to caxis (X,Y,Z), corresponding to iaxis (1,2,3),
2821 // and placed at the distance cutval from the origin.
2822 // The resulting picture is seen from the the same axis.
2823 // When HIDE Mode is ON, it is possible to get the same effect with
2824 // the CVOL/BOX function.
2830 gdrawc(PASSCHARD(vname), axis,cut,u0,v0,ul,vl PASSCHARL(vname));
2833 //_____________________________________________________________________________
2834 void TGeant3::Gdrawx(const char *name,Float_t cutthe, Float_t cutphi,
2835 Float_t cutval, Float_t theta, Float_t phi, Float_t u0,
2836 Float_t v0,Float_t ul,Float_t vl)
2840 // CUTTHE Theta angle of the line normal to cut plane
2841 // CUTPHI Phi angle of the line normal to cut plane
2842 // CUTVAL Cut plane distance from the origin along the axis
2844 // THETA Viewing angle theta (for 3D projection)
2845 // PHI Viewing angle phi (for 3D projection)
2846 // U0 U-coord. (horizontal) of volume origin
2847 // V0 V-coord. (vertical) of volume origin
2848 // SU Scale factor for U-coord.
2849 // SV Scale factor for V-coord.
2851 // The cut plane is normal to the line given by the cut angles
2852 // cutthe and cutphi and placed at the distance cutval from the origin.
2853 // The resulting picture is seen from the viewing angles theta,phi.
2859 gdrawx(PASSCHARD(vname), cutthe,cutphi,cutval,theta,phi,u0,v0,ul,vl
2863 //_____________________________________________________________________________
2864 void TGeant3::Gdhead(Int_t isel, const char *name, Float_t chrsiz)
2869 // ISEL Option flag D=111110
2871 // CHRSIZ Character size (cm) of title NAME D=0.6
2874 // 0 to have only the header lines
2875 // xxxxx1 to add the text name centered on top of header
2876 // xxxx1x to add global detector name (first volume) on left
2877 // xxx1xx to add date on right
2878 // xx1xxx to select thick characters for text on top of header
2879 // x1xxxx to add the text 'EVENT NR x' on top of header
2880 // 1xxxxx to add the text 'RUN NR x' on top of header
2881 // NOTE that ISEL=x1xxx1 or ISEL=1xxxx1 are illegal choices,
2882 // i.e. they generate overwritten text.
2884 gdhead(isel,PASSCHARD(name),chrsiz PASSCHARL(name));
2887 //_____________________________________________________________________________
2888 void TGeant3::Gdman(Float_t u, Float_t v, const char *type)
2891 // Draw a 2D-man at position (U0,V0)
2893 // U U-coord. (horizontal) of the centre of man' R
2894 // V V-coord. (vertical) of the centre of man' R
2895 // TYPE D='MAN' possible values: 'MAN,WM1,WM2,WM3'
2897 // CALL GDMAN(u,v),CALL GDWMN1(u,v),CALL GDWMN2(u,v),CALL GDWMN2(u,v)
2898 // It superimposes the picure of a man or of a woman, chosen among
2899 // three different ones, with the same scale factors as the detector
2900 // in the current drawing.
2903 if (opt.Contains("WM1")) {
2905 } else if (opt.Contains("WM3")) {
2907 } else if (opt.Contains("WM2")) {
2914 //_____________________________________________________________________________
2915 void TGeant3::Gdspec(const char *name)
2920 // Shows 3 views of the volume (two cut-views and a 3D view), together with
2921 // its geometrical specifications. The 3D drawing will
2922 // be performed according the current values of the options HIDE and
2923 // SHAD and according the current SetClipBox clipping parameters for that
2930 gdspec(PASSCHARD(vname) PASSCHARL(vname));
2933 //_____________________________________________________________________________
2934 void TGeant3::DrawOneSpec(const char *name)
2937 // Function called when one double-clicks on a volume name
2938 // in a TPavelabel drawn by Gdtree.
2940 THIGZ *higzSave = higz;
2941 higzSave->SetName("higzSave");
2942 THIGZ *higzSpec = (THIGZ*)gROOT->FindObject("higzSpec");
2943 //printf("DrawOneSpec, higz=%x, higzSpec=%x\n",higz,higzSpec);
2944 if (higzSpec) higz = higzSpec;
2945 else higzSpec = new THIGZ(defSize);
2946 higzSpec->SetName("higzSpec");
2951 gdspec(PASSCHARD(vname) PASSCHARL(vname));
2954 higzSave->SetName("higz");
2958 //_____________________________________________________________________________
2959 void TGeant3::Gdtree(const char *name,Int_t levmax, Int_t isel)
2963 // LEVMAX Depth level
2966 // This function draws the logical tree,
2967 // Each volume in the tree is represented by a TPaveTree object.
2968 // Double-clicking on a TPaveTree draws the specs of the corresponding volume.
2969 // Use TPaveTree pop-up menu to select:
2972 // - drawing tree of parent
2978 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
2982 //_____________________________________________________________________________
2983 void TGeant3::GdtreeParent(const char *name,Int_t levmax, Int_t isel)
2987 // LEVMAX Depth level
2990 // This function draws the logical tree of the parent of name.
2994 // Scan list of volumes in JVOLUM
2996 Int_t gname, i, jvo, in, nin, jin, num;
2997 strncpy((char *) &gname, name, 4);
2998 for(i=1; i<=fGcnum->nvolum; i++) {
2999 jvo = fZlq[fGclink->jvolum-i];
3000 nin = Int_t(fZq[jvo+3]);
3001 if (nin == -1) nin = 1;
3002 for (in=1;in<=nin;in++) {
3004 num = Int_t(fZq[jin+2]);
3005 if(gname == fZiq[fGclink->jvolum+num]) {
3006 strncpy(vname,(char*)&fZiq[fGclink->jvolum+i],4);
3008 gdtree(PASSCHARD(vname), levmax, isel PASSCHARL(vname));
3016 //_____________________________________________________________________________
3017 void TGeant3::SetABAN(Int_t par)
3020 // par = 1 particles will be stopped according to their residual
3021 // range if they are not in a sensitive material and are
3022 // far enough from the boundary
3023 // 0 particles are transported normally
3025 fGcphys->dphys1 = par;
3029 //_____________________________________________________________________________
3030 void TGeant3::SetANNI(Int_t par)
3033 // To control positron annihilation.
3034 // par =0 no annihilation
3035 // =1 annihilation. Decays processed.
3036 // =2 annihilation. No decay products stored.
3038 fGcphys->ianni = par;
3042 //_____________________________________________________________________________
3043 void TGeant3::SetAUTO(Int_t par)
3046 // To control automatic calculation of tracking medium parameters:
3047 // par =0 no automatic calculation;
3048 // =1 automati calculation.
3050 fGctrak->igauto = par;
3054 //_____________________________________________________________________________
3055 void TGeant3::SetBOMB(Float_t boom)
3058 // BOOM : Exploding factor for volumes position
3060 // To 'explode' the detector. If BOOM is positive (values smaller
3061 // than 1. are suggested, but any value is possible)
3062 // all the volumes are shifted by a distance
3063 // proportional to BOOM along the direction between their centre
3064 // and the origin of the MARS; the volumes which are symmetric
3065 // with respect to this origin are simply not shown.
3066 // BOOM equal to 0 resets the normal mode.
3067 // A negative (greater than -1.) value of
3068 // BOOM will cause an 'implosion'; for even lower values of BOOM
3069 // the volumes' positions will be reflected respect to the origin.
3070 // This command can be useful to improve the 3D effect for very
3071 // complex detectors. The following commands will make explode the
3078 //_____________________________________________________________________________
3079 void TGeant3::SetBREM(Int_t par)
3082 // To control bremstrahlung.
3083 // par =0 no bremstrahlung
3084 // =1 bremstrahlung. Photon processed.
3085 // =2 bremstrahlung. No photon stored.
3087 fGcphys->ibrem = par;
3091 //_____________________________________________________________________________
3092 void TGeant3::SetCKOV(Int_t par)
3095 // To control Cerenkov production
3096 // par =0 no Cerenkov;
3098 // =2 Cerenkov with primary stopped at each step.
3100 fGctlit->itckov = par;
3104 //_____________________________________________________________________________
3105 void TGeant3::SetClipBox(const char *name,Float_t xmin,Float_t xmax,
3106 Float_t ymin,Float_t ymax,Float_t zmin,Float_t zmax)
3109 // The hidden line removal technique is necessary to visualize properly
3110 // very complex detectors. At the same time, it can be useful to visualize
3111 // the inner elements of a detector in detail. This function allows
3112 // subtractions (via boolean operation) of BOX shape from any part of
3113 // the detector, therefore showing its inner contents.
3114 // If "*" is given as the name of the
3115 // volume to be clipped, all volumes are clipped by the given box.
3116 // A volume can be clipped at most twice.
3117 // if a volume is explicitely clipped twice,
3118 // the "*" will not act on it anymore. Giving "." as the name
3119 // of the volume to be clipped will reset the clipping.
3121 // NAME Name of volume to be clipped
3123 // XMIN Lower limit of the Shape X coordinate
3124 // XMAX Upper limit of the Shape X coordinate
3125 // YMIN Lower limit of the Shape Y coordinate
3126 // YMAX Upper limit of the Shape Y coordinate
3127 // ZMIN Lower limit of the Shape Z coordinate
3128 // ZMAX Upper limit of the Shape Z coordinate
3130 // This function performs a boolean subtraction between the volume
3131 // NAME and a box placed in the MARS according the values of the given
3137 setclip(PASSCHARD(vname),xmin,xmax,ymin,ymax,zmin,zmax PASSCHARL(vname));
3140 //_____________________________________________________________________________
3141 void TGeant3::SetCOMP(Int_t par)
3144 // To control Compton scattering
3145 // par =0 no Compton
3146 // =1 Compton. Electron processed.
3147 // =2 Compton. No electron stored.
3150 fGcphys->icomp = par;
3153 //_____________________________________________________________________________
3154 void TGeant3::SetCUTS(Float_t cutgam,Float_t cutele,Float_t cutneu,
3155 Float_t cuthad,Float_t cutmuo ,Float_t bcute ,
3156 Float_t bcutm ,Float_t dcute ,Float_t dcutm ,
3157 Float_t ppcutm, Float_t tofmax)
3160 // CUTGAM Cut for gammas D=0.001
3161 // CUTELE Cut for electrons D=0.001
3162 // CUTHAD Cut for charged hadrons D=0.01
3163 // CUTNEU Cut for neutral hadrons D=0.01
3164 // CUTMUO Cut for muons D=0.01
3165 // BCUTE Cut for electron brems. D=-1.
3166 // BCUTM Cut for muon brems. D=-1.
3167 // DCUTE Cut for electron delta-rays D=-1.
3168 // DCUTM Cut for muon delta-rays D=-1.
3169 // PPCUTM Cut for e+e- pairs by muons D=0.01
3170 // TOFMAX Time of flight cut D=1.E+10
3172 // If the default values (-1.) for BCUTE ,BCUTM ,DCUTE ,DCUTM
3173 // are not modified, they will be set to CUTGAM,CUTGAM,CUTELE,CUTELE
3175 // If one of the parameters from CUTGAM to PPCUTM included
3176 // is modified, cross-sections and energy loss tables must be
3177 // recomputed via the function Gphysi.
3179 fGccuts->cutgam = cutgam;
3180 fGccuts->cutele = cutele;
3181 fGccuts->cutneu = cutneu;
3182 fGccuts->cuthad = cuthad;
3183 fGccuts->cutmuo = cutmuo;
3184 fGccuts->bcute = bcute;
3185 fGccuts->bcutm = bcutm;
3186 fGccuts->dcute = dcute;
3187 fGccuts->dcutm = dcutm;
3188 fGccuts->ppcutm = ppcutm;
3189 fGccuts->tofmax = tofmax;
3192 //_____________________________________________________________________________
3193 void TGeant3::SetDCAY(Int_t par)
3196 // To control Decay mechanism.
3197 // par =0 no decays.
3198 // =1 Decays. secondaries processed.
3199 // =2 Decays. No secondaries stored.
3201 fGcphys->idcay = par;
3205 //_____________________________________________________________________________
3206 void TGeant3::SetDEBU(Int_t emin, Int_t emax, Int_t emod)
3209 // Set the debug flag and frequency
3210 // Selected debug output will be printed from
3211 // event emin to even emax each emod event
3213 fGcflag->idemin = emin;
3214 fGcflag->idemax = emax;
3215 fGcflag->itest = emod;
3219 //_____________________________________________________________________________
3220 void TGeant3::SetDRAY(Int_t par)
3223 // To control delta rays mechanism.
3224 // par =0 no delta rays.
3225 // =1 Delta rays. secondaries processed.
3226 // =2 Delta rays. No secondaries stored.
3228 fGcphys->idray = par;
3231 //_____________________________________________________________________________
3232 void TGeant3::SetERAN(Float_t ekmin, Float_t ekmax, Int_t nekbin)
3235 // To control cross section tabulations
3236 // ekmin = minimum kinetic energy in GeV
3237 // ekmax = maximum kinetic energy in GeV
3238 // nekbin = number of logatithmic bins (<200)
3240 fGcmulo->ekmin = ekmin;
3241 fGcmulo->ekmax = ekmax;
3242 fGcmulo->nekbin = nekbin;
3245 //_____________________________________________________________________________
3246 void TGeant3::SetHADR(Int_t par)
3249 // To control hadronic interactions.
3250 // par =0 no hadronic interactions.
3251 // =1 Hadronic interactions. secondaries processed.
3252 // =2 Hadronic interactions. No secondaries stored.
3254 fGcphys->ihadr = par;
3257 //_____________________________________________________________________________
3258 void TGeant3::SetKINE(Int_t kine, Float_t xk1, Float_t xk2, Float_t xk3,
3259 Float_t xk4, Float_t xk5, Float_t xk6, Float_t xk7,
3260 Float_t xk8, Float_t xk9, Float_t xk10)
3263 // Set the variables in /GCFLAG/ IKINE, PKINE(10)
3264 // Their meaning is user defined
3266 fGckine->ikine = kine;
3267 fGckine->pkine[0] = xk1;
3268 fGckine->pkine[1] = xk2;
3269 fGckine->pkine[2] = xk3;
3270 fGckine->pkine[3] = xk4;
3271 fGckine->pkine[4] = xk5;
3272 fGckine->pkine[5] = xk6;
3273 fGckine->pkine[6] = xk7;
3274 fGckine->pkine[7] = xk8;
3275 fGckine->pkine[8] = xk9;
3276 fGckine->pkine[9] = xk10;
3279 //_____________________________________________________________________________
3280 void TGeant3::SetLOSS(Int_t par)
3283 // To control energy loss.
3284 // par =0 no energy loss;
3285 // =1 restricted energy loss fluctuations;
3286 // =2 complete energy loss fluctuations;
3288 // =4 no energy loss fluctuations.
3289 // If the value ILOSS is changed, then cross-sections and energy loss
3290 // tables must be recomputed via the command 'PHYSI'.
3292 fGcphys->iloss = par;
3296 //_____________________________________________________________________________
3297 void TGeant3::SetMULS(Int_t par)
3300 // To control multiple scattering.
3301 // par =0 no multiple scattering.
3302 // =1 Moliere or Coulomb scattering.
3303 // =2 Moliere or Coulomb scattering.
3304 // =3 Gaussian scattering.
3306 fGcphys->imuls = par;
3310 //_____________________________________________________________________________
3311 void TGeant3::SetMUNU(Int_t par)
3314 // To control muon nuclear interactions.
3315 // par =0 no muon-nuclear interactions.
3316 // =1 Nuclear interactions. Secondaries processed.
3317 // =2 Nuclear interactions. Secondaries not processed.
3319 fGcphys->imunu = par;
3322 //_____________________________________________________________________________
3323 void TGeant3::SetOPTI(Int_t par)
3326 // This flag controls the tracking optimisation performed via the
3328 // 1 no optimisation at all; GSORD calls disabled;
3329 // 0 no optimisation; only user calls to GSORD kept;
3330 // 1 all non-GSORDered volumes are ordered along the best axis;
3331 // 2 all volumes are ordered along the best axis.
3333 fGcopti->ioptim = par;
3336 //_____________________________________________________________________________
3337 void TGeant3::SetPAIR(Int_t par)
3340 // To control pair production mechanism.
3341 // par =0 no pair production.
3342 // =1 Pair production. secondaries processed.
3343 // =2 Pair production. No secondaries stored.
3345 fGcphys->ipair = par;
3349 //_____________________________________________________________________________
3350 void TGeant3::SetPFIS(Int_t par)
3353 // To control photo fission mechanism.
3354 // par =0 no photo fission.
3355 // =1 Photo fission. secondaries processed.
3356 // =2 Photo fission. No secondaries stored.
3358 fGcphys->ipfis = par;
3361 //_____________________________________________________________________________
3362 void TGeant3::SetPHOT(Int_t par)
3365 // To control Photo effect.
3366 // par =0 no photo electric effect.
3367 // =1 Photo effect. Electron processed.
3368 // =2 Photo effect. No electron stored.
3370 fGcphys->iphot = par;
3373 //_____________________________________________________________________________
3374 void TGeant3::SetRAYL(Int_t par)
3377 // To control Rayleigh scattering.
3378 // par =0 no Rayleigh scattering.
3381 fGcphys->irayl = par;
3384 //_____________________________________________________________________________
3385 void TGeant3::SetSWIT(Int_t sw, Int_t val)
3389 // val New switch value
3391 // Change one element of array ISWIT(10) in /GCFLAG/
3393 if (sw <= 0 || sw > 10) return;
3394 fGcflag->iswit[sw-1] = val;
3398 //_____________________________________________________________________________
3399 void TGeant3::SetTRIG(Int_t nevents)
3402 // Set number of events to be run
3404 fGcflag->nevent = nevents;
3407 //_____________________________________________________________________________
3408 void TGeant3::SetUserDecay(Int_t pdg)
3411 // Force the decays of particles to be done with Pythia
3412 // and not with the Geant routines.
3413 // just kill pointers doing mzdrop
3415 Int_t ipart = IdFromPDG(pdg);
3417 printf("Particle %d not in geant\n",pdg);
3420 Int_t jpart=fGclink->jpart;
3421 Int_t jpa=fZlq[jpart-ipart];
3424 Int_t jpa1=fZlq[jpa-1];
3426 mzdrop(fGcbank->ixcons,jpa1,PASSCHARD(" ") PASSCHARL(" "));
3427 Int_t jpa2=fZlq[jpa-2];
3429 mzdrop(fGcbank->ixcons,jpa2,PASSCHARD(" ") PASSCHARL(" "));
3433 //______________________________________________________________________________
3434 void TGeant3::Vname(const char *name, char *vname)
3437 // convert name to upper case. Make vname at least 4 chars
3439 Int_t l = strlen(name);
3442 for (i=0;i<l;i++) vname[i] = toupper(name[i]);
3443 for (i=l;i<4;i++) vname[i] = ' ';
3447 //______________________________________________________________________________
3448 void TGeant3::Ertrgo()
3453 //______________________________________________________________________________
3454 void TGeant3::Ertrak(const Float_t *const x1, const Float_t *const p1,
3455 const Float_t *x2, const Float_t *p2,
3456 Int_t ipa, Option_t *chopt)
3458 ertrak(x1,p1,x2,p2,ipa,PASSCHARD(chopt) PASSCHARL(chopt));
3461 //_____________________________________________________________________________
3462 void TGeant3::WriteEuclid(const char* filnam, const char* topvol,
3463 Int_t number, Int_t nlevel)
3467 // ******************************************************************
3469 // * Write out the geometry of the detector in EUCLID file format *
3471 // * filnam : will be with the extension .euc *
3472 // * topvol : volume name of the starting node *
3473 // * number : copy number of topvol (relevant for gsposp) *
3474 // * nlevel : number of levels in the tree structure *
3475 // * to be written out, starting from topvol *
3477 // * Author : M. Maire *
3479 // ******************************************************************
3481 // File filnam.tme is written out with the definitions of tracking
3482 // medias and materials.
3483 // As to restore original numbers for materials and medias, program
3484 // searches in the file euc_medi.dat and comparing main parameters of
3485 // the mat. defined inside geant and the one in file recognizes them
3486 // and is able to take number from file. If for any material or medium,
3487 // this procedure fails, ordering starts from 1.
3488 // Arrays IOTMED and IOMATE are used for this procedure
3490 const char shape[][5]={"BOX ","TRD1","TRD2","TRAP","TUBE","TUBS","CONE",
3491 "CONS","SPHE","PARA","PGON","PCON","ELTU","HYPE",
3493 Int_t i, end, itm, irm, jrm, k, nmed;
3497 char *filext, *filetme;
3498 char natmed[21], namate[21];
3499 char natmedc[21], namatec[21];
3500 char key[5], name[5], mother[5], konly[5];
3502 Int_t iadvol, iadtmd, iadrot, nwtot, iret;
3503 Int_t mlevel, numbr, natt, numed, nin, ndata;
3504 Int_t iname, ivo, ish, jvo, nvstak, ivstak;
3505 Int_t jdiv, ivin, in, jin, jvin, irot;
3506 Int_t jtm, imat, jma, flag=0, imatc;
3507 Float_t az, dens, radl, absl, a, step, x, y, z;
3508 Int_t npar, ndvmx, left;
3509 Float_t zc, densc, radlc, abslc, c0, tmaxfd;
3511 Int_t iomate[100], iotmed[100];
3512 Float_t par[50], att[20], ubuf[50];
3515 Int_t level, ndiv, iaxe;
3516 Int_t itmedc, nmatc, isvolc, ifieldc, nwbufc, isvol, nmat, ifield, nwbuf;
3517 Float_t fieldmc, tmaxfdc, stemaxc, deemaxc, epsilc, stminc, fieldm;
3518 Float_t tmaxf, stemax, deemax, epsil, stmin;
3519 const char *f10000="!\n%s\n!\n";
3520 //Open the input file
3522 for(i=0;i<end;i++) if(filnam[i]=='.') {
3526 filext=new char[end+5];
3527 filetme=new char[end+5];
3528 strncpy(filext,filnam,end);
3529 strncpy(filetme,filnam,end);
3531 // *** The output filnam name will be with extension '.euc'
3532 strcpy(&filext[end],".euc");
3533 strcpy(&filetme[end],".tme");
3534 lun=fopen(filext,"w");
3536 // *** Initialisation of the working space
3537 iadvol=fGcnum->nvolum;
3538 iadtmd=iadvol+fGcnum->nvolum;
3539 iadrot=iadtmd+fGcnum->ntmed;
3540 if(fGclink->jrotm) {
3541 fGcnum->nrotm=fZiq[fGclink->jrotm-2];
3545 nwtot=iadrot+fGcnum->nrotm;
3546 qws = new float[nwtot+1];
3547 for (i=0;i<nwtot+1;i++) qws[i]=0;
3550 if(nlevel==0) mlevel=20;
3552 // *** find the top volume and put it in the stak
3553 numbr = number>0 ? number : 1;
3554 Gfpara(topvol,numbr,1,npar,natt,par,att);
3556 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3561 // *** authorized shape ?
3562 strncpy((char *)&iname, topvol, 4);
3564 for(i=1; i<=fGcnum->nvolum; i++) if(fZiq[fGclink->jvolum+i]==iname) {
3568 jvo = fZlq[fGclink->jvolum-ivo];
3569 ish = Int_t (fZq[jvo+2]);
3571 printf(" *** GWEUCL *** top volume : %s number : %3d can not be a valid root\n",
3578 iws[iadvol+ivo] = level;
3581 //*** flag all volumes and fill the stak
3585 // pick the next volume in stak
3587 ivo = TMath::Abs(iws[ivstak]);
3588 jvo = fZlq[fGclink->jvolum - ivo];
3590 // flag the tracking medium
3591 numed = Int_t (fZq[jvo + 4]);
3592 iws[iadtmd + numed] = 1;
3594 // get the daughters ...
3595 level = iws[iadvol+ivo];
3596 if (level < mlevel) {
3598 nin = Int_t (fZq[jvo + 3]);
3600 // from division ...
3602 jdiv = fZlq[jvo - 1];
3603 ivin = Int_t (fZq[jdiv + 2]);
3605 iws[nvstak] = -ivin;
3606 iws[iadvol+ivin] = level;
3608 // from position ...
3609 } else if (nin > 0) {
3610 for(in=1; in<=nin; in++) {
3611 jin = fZlq[jvo - in];
3612 ivin = Int_t (fZq[jin + 2 ]);
3613 jvin = fZlq[fGclink->jvolum - ivin];
3614 ish = Int_t (fZq[jvin + 2]);
3615 // authorized shape ?
3617 // not yet flagged ?
3618 if (iws[iadvol+ivin]==0) {
3621 iws[iadvol+ivin] = level;
3623 // flag the rotation matrix
3624 irot = Int_t ( fZq[jin + 4 ]);
3625 if (irot > 0) iws[iadrot+irot] = 1;
3631 // next volume in stak ?
3632 if (ivstak < nvstak) goto L10;
3634 // *** restore original material and media numbers
3635 // file euc_medi.dat is needed to compare materials and medias
3637 FILE* luncor=fopen("euc_medi.dat","r");
3640 for(itm=1; itm<=fGcnum->ntmed; itm++) {
3641 if (iws[iadtmd+itm] > 0) {
3642 jtm = fZlq[fGclink->jtmed-itm];
3643 strncpy(natmed,(char *)&fZiq[jtm+1],20);
3644 imat = Int_t (fZq[jtm+6]);
3645 jma = fZlq[fGclink->jmate-imat];
3647 printf(" *** GWEUCL *** material not defined for tracking medium %5i %s\n",itm,natmed);
3650 strncpy(namate,(char *)&fZiq[jma+1],20);
3653 //** find the material original number
3656 iret=fscanf(luncor,"%4s,%130s",key,card);
3657 if(iret<=0) goto L26;
3659 if(!strcmp(key,"MATE")) {
3660 sscanf(card,"%d %s %f %f %f %f %f %d",&imatc,namatec,&az,&zc,&densc,&radlc,&abslc,&nparc);
3661 Gfmate(imat,namate,a,z,dens,radl,absl,par,npar);
3662 if(!strcmp(namatec,namate)) {
3663 if(az==a && zc==z && densc==dens && radlc==radl
3664 && abslc==absl && nparc==nparc) {
3667 printf("*** GWEUCL *** material : %3d '%s' restored as %3d\n",imat,namate,imatc);
3669 printf("*** GWEUCL *** different definitions for material: %s\n",namate);
3673 if(strcmp(key,"END") && !flag) goto L23;
3675 printf("*** GWEUCL *** cannot restore original number for material: %s\n",namate);
3679 //*** restore original tracking medium number
3682 iret=fscanf(luncor,"%4s,%130s",key,card);
3683 if(iret<=0) goto L26;
3685 if (!strcmp(key,"TMED")) {
3686 sscanf(card,"%d %s %d %d %d %f %f %f %f %f %f %d\n",
3687 &itmedc,natmedc,&nmatc,&isvolc,&ifieldc,&fieldmc,
3688 &tmaxfdc,&stemaxc,&deemaxc,&epsilc,&stminc,&nwbufc);
3689 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxf,stemax,deemax,
3690 epsil,stmin,ubuf,&nwbuf);
3691 if(!strcmp(natmedc,natmed)) {
3692 if (iomate[nmat]==nmatc && nwbuf==nwbufc) {
3695 printf("*** GWEUCL *** medium : %3d '%20s' restored as %3d\n",
3698 printf("*** GWEUCL *** different definitions for tracking medium: %s\n",natmed);
3702 if(strcmp(key,"END") && !flag) goto L24;
3704 printf("cannot restore original number for medium : %s\n",natmed);
3712 L26: printf("*** GWEUCL *** cannot read the data file\n");
3714 L29: if(luncor) fclose (luncor);
3717 // *** write down the tracking medium definition
3719 strcpy(card,"! Tracking medium");
3720 fprintf(lun,f10000,card);
3722 for(itm=1;itm<=fGcnum->ntmed;itm++) {
3723 if (iws[iadtmd+itm]>0) {
3724 jtm = fZlq[fGclink->jtmed-itm];
3725 strncpy(natmed,(char *)&fZiq[jtm+1],20);
3727 imat = Int_t (fZq[jtm+6]);
3728 jma = fZlq[fGclink->jmate-imat];
3729 //* order media from one, if comparing with database failed
3731 iotmed[itm]=++imxtmed;
3732 iomate[imat]=++imxmate;
3737 printf(" *** GWEUCL *** material not defined for tracking medium %5d %s\n",
3740 strncpy(namate,(char *)&fZiq[jma+1],20);
3743 fprintf(lun,"TMED %3d '%20s' %3d '%20s'\n",iotmed[itm],natmed,iomate[imat],namate);
3747 //* *** write down the rotation matrix
3749 strcpy(card,"! Reperes");
3750 fprintf(lun,f10000,card);
3752 for(irm=1;irm<=fGcnum->nrotm;irm++) {
3753 if (iws[iadrot+irm]>0) {
3754 jrm = fZlq[fGclink->jrotm-irm];
3755 fprintf(lun,"ROTM %3d",irm);
3756 for(k=11;k<=16;k++) fprintf(lun," %8.3f",fZq[jrm+k]);
3761 //* *** write down the volume definition
3763 strcpy(card,"! Volumes");
3764 fprintf(lun,f10000,card);
3766 for(ivstak=1;ivstak<=nvstak;ivstak++) {
3769 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivo],4);
3771 jvo = fZlq[fGclink->jvolum-ivo];
3772 ish = Int_t (fZq[jvo+2]);
3773 nmed = Int_t (fZq[jvo+4]);
3774 npar = Int_t (fZq[jvo+5]);
3776 if (ivstak>1) for(i=0;i<npar;i++) par[i]=fZq[jvo+7+i];
3777 Gckpar (ish,npar,par);
3778 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,shape[ish-1],iotmed[nmed],npar);
3779 for(i=0;i<(npar-1)/6+1;i++) {
3782 for(k=0;k<(left<6?left:6);k++) fprintf(lun," %11.5f",par[i*6+k]);
3786 fprintf(lun,"VOLU '%4s' '%4s' %3d %3d\n",name,shape[ish-1],iotmed[nmed],npar);
3791 //* *** write down the division of volumes
3793 fprintf(lun,f10000,"! Divisions");
3794 for(ivstak=1;ivstak<=nvstak;ivstak++) {
3795 ivo = TMath::Abs(iws[ivstak]);
3796 jvo = fZlq[fGclink->jvolum-ivo];
3797 ish = Int_t (fZq[jvo+2]);
3798 nin = Int_t (fZq[jvo+3]);
3799 //* this volume is divided ...
3802 iaxe = Int_t ( fZq[jdiv+1]);
3803 ivin = Int_t ( fZq[jdiv+2]);
3804 ndiv = Int_t ( fZq[jdiv+3]);
3807 jvin = fZlq[fGclink->jvolum-ivin];
3808 nmed = Int_t ( fZq[jvin+4]);
3809 strncpy(mother,(char *)&fZiq[fGclink->jvolum+ivo ],4);
3811 strncpy(name,(char *)&fZiq[fGclink->jvolum+ivin],4);
3813 if ((step<=0.)||(ish>=11)) {
3814 //* volume with negative parameter or gsposp or pgon ...
3815 fprintf(lun,"DIVN '%4s' '%4s' %3d %3d\n",name,mother,ndiv,iaxe);
3816 } else if ((ndiv<=0)||(ish==10)) {
3817 //* volume with negative parameter or gsposp or para ...
3818 ndvmx = TMath::Abs(ndiv);
3819 fprintf(lun,"DIVT '%4s' '%4s' %11.5f %3d %3d %3d\n",
3820 name,mother,step,iaxe,iotmed[nmed],ndvmx);
3822 //* normal volume : all kind of division are equivalent
3823 fprintf(lun,"DVT2 '%4s' '%4s' %11.5f %3d %11.5f %3d %3d\n",
3824 name,mother,step,iaxe,c0,iotmed[nmed],ndiv);
3829 //* *** write down the the positionnement of volumes
3831 fprintf(lun,f10000,"! Positionnements\n");
3833 for(ivstak = 1;ivstak<=nvstak;ivstak++) {
3834 ivo = TMath::Abs(iws[ivstak]);
3835 strncpy(mother,(char*)&fZiq[fGclink->jvolum+ivo ],4);
3837 jvo = fZlq[fGclink->jvolum-ivo];
3838 nin = Int_t( fZq[jvo+3]);
3839 //* this volume has daughters ...
3841 for (in=1;in<=nin;in++) {
3843 ivin = Int_t (fZq[jin +2]);
3844 numb = Int_t (fZq[jin +3]);
3845 irot = Int_t (fZq[jin +4]);
3849 strcpy(konly,"ONLY");
3850 if (fZq[jin+8]!=1.) strcpy(konly,"MANY");
3851 strncpy(name,(char*)&fZiq[fGclink->jvolum+ivin],4);
3853 jvin = fZlq[fGclink->jvolum-ivin];
3854 ish = Int_t (fZq[jvin+2]);
3855 //* gspos or gsposp ?
3856 ndata = fZiq[jin-1];
3858 fprintf(lun,"POSI '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s'\n",
3859 name,numb,mother,x,y,z,irot,konly);
3861 npar = Int_t (fZq[jin+9]);
3862 for(i=0;i<npar;i++) par[i]=fZq[jin+10+i];
3863 Gckpar (ish,npar,par);
3864 fprintf(lun,"POSP '%4s' %4d '%4s' %11.5f %11.5f %11.5f %3d '%4s' %3d\n",
3865 name,numb,mother,x,y,z,irot,konly,npar);
3867 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
3874 fprintf(lun,"END\n");
3877 //****** write down the materials and medias *****
3879 lun=fopen(filetme,"w");
3881 for(itm=1;itm<=fGcnum->ntmed;itm++) {
3882 if (iws[iadtmd+itm]>0) {
3883 jtm = fZlq[fGclink->jtmed-itm];
3884 strncpy(natmed,(char*)&fZiq[jtm+1],4);
3885 imat = Int_t (fZq[jtm+6]);
3886 jma = Int_t (fZlq[fGclink->jmate-imat]);
3888 Gfmate (imat,namate,a,z,dens,radl,absl,par,npar);
3889 fprintf(lun,"MATE %4d '%20s'%11.5E %11.5E %11.5E %11.5E %11.5E %3d\n",
3890 iomate[imat],namate,a,z,dens,radl,absl,npar);
3894 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
3898 Gftmed(itm,natmed,nmat,isvol,ifield,fieldm,tmaxfd,stemax,deemax,epsil,stmin,par,&npar);
3899 fprintf(lun,"TMED %4d '%20s' %3d %1d %3d %11.5f %11.5f %11.5f %11.5f %11.5f %11.5f %3d\n",
3900 iotmed[itm],natmed,iomate[nmat],isvol,ifield,
3901 fieldm,tmaxfd,stemax,deemax,epsil,stmin,npar);
3905 for(i=0;i<npar;i++) fprintf(lun," %11.5f",par[i]);
3911 fprintf(lun,"END\n");
3913 printf(" *** GWEUCL *** file: %s is now written out\n",filext);
3914 printf(" *** GWEUCL *** file: %s is now written out\n",filetme);
3923 //_____________________________________________________________________________
3924 void TGeant3::Streamer(TBuffer &R__b)
3927 // Stream an object of class TGeant3.
3929 if (R__b.IsReading()) {
3930 Version_t R__v = R__b.ReadVersion(); if (R__v) { }
3931 AliMC::Streamer(R__b);
3934 R__b.ReadStaticArray(fPDGCode);
3936 R__b.WriteVersion(TGeant3::IsA());
3937 AliMC::Streamer(R__b);
3940 R__b.WriteArray(fPDGCode, fNPDGCodes);