]> git.uio.no Git - u/mrichter/AliRoot.git/blame - TFluka/TFlukaMCGeometry.cxx
Coding rule violations corrected (A. Gheata).
[u/mrichter/AliRoot.git] / TFluka / TFlukaMCGeometry.cxx
CommitLineData
e2e989f9 1/**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3 * *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
6 * *
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 **************************************************************************/
15
16// $Id$
8495a208 17// Author: Andrei Gheata 10/07/2003
18
19#include "TObjString.h"
89125dc7 20#include "TFluka.h"
8495a208 21#include "TFlukaMCGeometry.h"
22#include "TGeoManager.h"
23#include "TGeoVolume.h"
24
25#include "TCallf77.h"
26
27#ifndef WIN32
28# define idnrwr idnrwr_
29# define g1wr g1wr_
30# define g1rtwr g1rtwr_
31# define conhwr conhwr_
32# define inihwr inihwr_
33# define jomiwr jomiwr_
34# define lkdbwr lkdbwr_
35# define lkfxwr lkfxwr_
36# define lkmgwr lkmgwr_
37# define lkwr lkwr_
38# define magfld magfld_
39# define nrmlwr nrmlwr_
40# define rgrpwr rgrpwr_
41# define isvhwr isvhwr_
42
43#else
44
45# define idnrwr IDNRWR
46# define g1wr G1WR
47# define g1rtwr G1RTWR
48# define conhwr CONHWR
49# define inihwr INIHWR
50# define jomiwr JOMIWR
51# define lkdbwr LKDBWR
52# define lkfxwr LKFXWR
53# define lkmgwr LKMGWR
54# define lkwr LKWR
55# define magfld MAGFLD
56# define nrmlwr NRMLWR
57# define rgrpwr RGRPWR
58# define isvhwr ISVHWR
59
60#endif
61
62//____________________________________________________________________________
63extern "C"
64{
65 //
66 // Prototypes for FLUKA navigation methods
67 //
68 Int_t type_of_call idnrwr(const Int_t & /*nreg*/, const Int_t & /*mlat*/);
69 void type_of_call g1wr(Double_t & /*pSx*/, Double_t & /*pSy*/, Double_t & /*pSz*/,
70 Double_t * /*pV*/, Int_t & /*oldReg*/ , const Int_t & /*oldLttc*/, Double_t & /*propStep*/,
71 Int_t & /*nascFlag*/, Double_t & /*retStep*/, Int_t & /*newReg*/,
72 Double_t & /*saf*/, Int_t & /*newLttc*/, Int_t & /*LttcFlag*/,
73 Double_t *s /*Lt*/, Int_t * /*jrLt*/);
74
75 void type_of_call g1rtwr();
76 void type_of_call conhwr(Int_t & /*intHist*/, Int_t * /*incrCount*/);
77 void type_of_call inihwr(Int_t & /*intHist*/);
78 void type_of_call jomiwr(const Int_t & /*nge*/, const Int_t & /*lin*/, const Int_t & /*lou*/,
79 Int_t & /*flukaReg*/);
80 void type_of_call lkdbwr(Double_t & /*pSx*/, Double_t & /*pSy*/, Double_t & /*pSz*/,
81 Double_t * /*pV*/, const Int_t & /*oldReg*/, const Int_t & /*oldLttc*/,
82 Int_t & /*newReg*/, Int_t & /*flagErr*/, Int_t & /*newLttc*/);
83 void type_of_call lkfxwr(Double_t & /*pSx*/, Double_t & /*pSy*/, Double_t & /*pSz*/,
84 Double_t * /*pV*/, const Int_t & /*oldReg*/, const Int_t & /*oldLttc*/,
85 Int_t & /*newReg*/, Int_t & /*flagErr*/, Int_t & /*newLttc*/);
86 void type_of_call lkmgwr(Double_t & /*pSx*/, Double_t & /*pSy*/, Double_t & /*pSz*/,
87 Double_t * /*pV*/, const Int_t & /*oldReg*/, const Int_t & /*oldLttc*/,
88 Int_t & /*flagErr*/, Int_t & /*newReg*/, Int_t & /*newLttc*/);
89 void type_of_call lkwr(Double_t & /*pSx*/, Double_t & /*pSy*/, Double_t & /*pSz*/,
90 Double_t * /*pV*/, const Int_t & /*oldReg*/, const Int_t & /*oldLttc*/,
91 Int_t & /*newReg*/, Int_t & /*flagErr*/, Int_t & /*newLttc*/);
efde9b4d 92// void type_of_call magfld(const Double_t & /*pX*/, const Double_t & /*pY*/, const Double_t & /*pZ*/,
93// Double_t & /*cosBx*/, Double_t & /*cosBy*/, Double_t & /*cosBz*/,
94// Double_t & /*Bmag*/, Int_t & /*reg*/, Int_t & /*idiscflag*/);
8495a208 95 void type_of_call nrmlwr(Double_t & /*pSx*/, Double_t & /*pSy*/, Double_t & /*pSz*/,
96 Double_t & /*pVx*/, Double_t & /*pVy*/, Double_t & /*pVz*/,
97 Double_t * /*norml*/, const Int_t & /*oldReg*/,
98 const Int_t & /*newReg*/, Int_t & /*flagErr*/);
99 void type_of_call rgrpwr(const Int_t & /*flukaReg*/, const Int_t & /*ptrLttc*/, Int_t & /*g4Reg*/,
100 Int_t * /*indMother*/, Int_t * /*repMother*/, Int_t & /*depthFluka*/);
101 Int_t type_of_call isvhwr(const Int_t & /*fCheck*/, const Int_t & /*intHist*/);
102};
103
104// TFluka global pointer
89125dc7 105TFluka *gFluka = 0;
106TFlukaMCGeometry *gMCGeom = 0;
6f5667d1 107Int_t gNstep = 0;
8495a208 108
109ClassImp(TFlukaMCGeometry)
110
111TFlukaMCGeometry* TFlukaMCGeometry::fgInstance=0;
112
113//_____________________________________________________________________________
114TFlukaMCGeometry::TFlukaMCGeometry(const char *name, const char *title)
89125dc7 115 : TNamed(name, title)
8495a208 116{
117 //
118 // Standard constructor
119 //
2bc4c610 120 fDebug = kFALSE;
efde9b4d 121 fLastMaterial = 0;
b1536e91 122 fNextRegion = 0;
123 fNextLattice = 0;
124 fRegionList = 0;
89125dc7 125 gFluka = (TFluka*)gMC;
126 gMCGeom = this;
6f5667d1 127 gNstep = 0;
89125dc7 128 fMatList = new TObjArray(256);
129 fMatNames = new TObjArray(256);
8495a208 130}
131
132//_____________________________________________________________________________
133TFlukaMCGeometry::TFlukaMCGeometry()
8495a208 134{
135 //
136 // Default constructor
137 //
2bc4c610 138 fDebug = kFALSE;
efde9b4d 139 fLastMaterial = 0;
b1536e91 140 fNextRegion = 0;
141 fNextLattice = 0;
142 fRegionList = 0;
89125dc7 143 gFluka = (TFluka*)gMC;
144 gMCGeom = this;
6f5667d1 145 gNstep = 0;
89125dc7 146 fMatList = 0;
147 fMatNames = 0;
8495a208 148}
149
150//_____________________________________________________________________________
151TFlukaMCGeometry::~TFlukaMCGeometry()
152{
153 //
154 // Destructor
155 //
156 fgInstance=0;
b1536e91 157 if (fRegionList) delete [] fRegionList;
89125dc7 158 if (fMatList) delete fMatList;
159 if (fMatNames) {fMatNames->Delete(); delete fMatNames;}
8495a208 160 if (gGeoManager) delete gGeoManager;
161}
162
163//
164// private methods
165//
166//_____________________________________________________________________________
167TFlukaMCGeometry::TFlukaMCGeometry(const TFlukaMCGeometry &)
89125dc7 168 : TNamed()
8495a208 169{
170 //
171 // Copy constructor
172 //
173}
174
175//_____________________________________________________________________________
176Double_t* TFlukaMCGeometry::CreateDoubleArray(Float_t* array, Int_t size) const
177{
178// Converts Float_t* array to Double_t*,
179// !! The new array has to be deleted by user.
180// ---
181
182 Double_t* doubleArray;
183 if (size>0) {
184 doubleArray = new Double_t[size];
185 for (Int_t i=0; i<size; i++) doubleArray[i] = array[i];
186 }
187 else {
188 //doubleArray = 0;
189 doubleArray = new Double_t[1];
190 }
191 return doubleArray;
192}
193//
194// public methods
8495a208 195
8495a208 196
8495a208 197//_____________________________________________________________________________
198Int_t TFlukaMCGeometry::GetMedium() const
199{
200// Get current medium number
201 Int_t imed = 0;
202 TGeoNode *node = gGeoManager->GetCurrentNode();
203 if (!node) imed = gGeoManager->GetTopNode()->GetVolume()->GetMedium()->GetId();
204 else imed = node->GetVolume()->GetMedium()->GetId();
2bc4c610 205 if (fDebug) printf("GetMedium=%i\n", imed);
8495a208 206 return imed;
207}
208
0bc73b6c 209//_____________________________________________________________________________
210Int_t TFlukaMCGeometry::GetFlukaMaterial(Int_t imed) const
211{
212// Returns FLUKA material index for medium IMED
213 TGeoMedium *med = (TGeoMedium*)gGeoManager->GetListOfMedia()->At(imed-1);
214 if (!med) {
215 Error("GetFlukaMaterial", "MEDIUM %i nor found", imed);
216 return -1;
217 }
218 Int_t imatfl = med->GetMaterial()->GetIndex();
219 return imatfl;
220}
221
b1536e91 222//_____________________________________________________________________________
223Int_t *TFlukaMCGeometry::GetRegionList(Int_t imed, Int_t &nreg)
224{
225// Get an ordered list of regions matching a given medium number
226 nreg = 0;
227 if (!fRegionList) fRegionList = new Int_t[NofVolumes()+1];
228 TIter next(gGeoManager->GetListOfUVolumes());
229 TGeoVolume *vol;
230 Int_t imedium, ireg;
231 while ((vol = (TGeoVolume*)next())) {
232 imedium = vol->GetMedium()->GetId();
233 if (imedium == imed) {
234 ireg = vol->GetNumber();
235 fRegionList[nreg++] = ireg;
236 }
237 }
238 return fRegionList;
239}
240
241//_____________________________________________________________________________
242Int_t *TFlukaMCGeometry::GetMaterialList(Int_t imat, Int_t &nreg)
243{
244// Get an ordered list of regions matching a given medium number
245 nreg = 0;
246 if (!fRegionList) fRegionList = new Int_t[NofVolumes()+1];
247 TIter next(gGeoManager->GetListOfUVolumes());
248 TGeoVolume *vol;
249 Int_t imaterial, ireg;
250 while ((vol = (TGeoVolume*)next())) {
251 imaterial = vol->GetMedium()->GetMaterial()->GetIndex();
252 if (imaterial == imat) {
253 ireg = vol->GetNumber();
254 fRegionList[nreg++] = ireg;
255 }
256 }
257 return fRegionList;
258}
8495a208 259
8495a208 260//_____________________________________________________________________________
261Int_t TFlukaMCGeometry::NofVolumes() const
262{
263 //
264 // Return total number of volumes in the geometry
265 //
266
267 return gGeoManager->GetListOfUVolumes()->GetEntriesFast()-1;
268}
8495a208 269
270//_____________________________________________________________________________
89125dc7 271TGeoMaterial * TFlukaMCGeometry::GetMakeWrongMaterial(Double_t z)
272{
273// Try to replace a wrongly-defined material
274 static Double_t kz[23] = {7.3, 17.8184, 7.2167, 10.856, 8.875, 8.9, 7.177,
275 25.72, 6.2363, 7.1315, 47.7056, 10.6467, 7.8598, 2.10853, 10.6001, 9.1193,
276 15.3383, 4.55, 9.6502, 6.4561, 21.7963, 29.8246, 15.4021};
277
278 Int_t ind;
279 Double_t dz;
280 for (ind=0; ind<23; ind++) {
281 dz = TMath::Abs(z-kz[ind]);
282 if (dz<1E-4) break;
283 }
284 if (ind>22) {
285 printf("Cannot patch material with Z=%g\n", z);
286 return 0;
287 }
288 TGeoMixture *mix = 0;
289 TGeoElement *element;
290 TGeoElementTable *table = TGeoElementTable::Instance();
291 switch (ind) {
292 case 0: // AIR
293 mix = new TGeoMixture("AIR", 4, 0.001205);
294 element = table->GetElement(6); // C
295 mix->DefineElement(0, element, 0.000124);
296 element = table->GetElement(7); // N
297 mix->DefineElement(1, element, 0.755267);
298 element = table->GetElement(8); // O
299 mix->DefineElement(2, element, 0.231781);
300 element = table->GetElement(18); // AR
301 mix->DefineElement(3, element, 0.012827);
302 break;
303 case 1: //SDD SI CHIP
304 mix = new TGeoMixture("SDD_SI", 6, 2.4485);
305 element = table->GetElement(1);
306 mix->DefineElement(0, element, 0.004367771);
307 element = table->GetElement(6);
308 mix->DefineElement(1, element, 0.039730642);
309 element = table->GetElement(7);
310 mix->DefineElement(2, element, 0.001396798);
311 element = table->GetElement(8);
312 mix->DefineElement(3, element, 0.01169634);
313 element = table->GetElement(14);
314 mix->DefineElement(4, element, 0.844665);
315 element = table->GetElement(47);
316 mix->DefineElement(5, element, 0.09814344903);
317 break;
318 case 2: // WATER
319 mix = new TGeoMixture("WATER", 2, 1.0);
320 element = table->GetElement(1);
321 mix->DefineElement(0, element, 0.111898344);
322 element = table->GetElement(8);
323 mix->DefineElement(1, element, 0.888101656);
324 break;
325 case 3: // CERAMICS
326 mix = new TGeoMixture("CERAMICS", 5, 3.6);
327 element = table->GetElement(8);
328 mix->DefineElement(0, element, 0.59956);
329 element = table->GetElement(13);
330 mix->DefineElement(1, element, 0.3776);
331 element = table->GetElement(14);
332 mix->DefineElement(2, element, 0.00933);
333 element = table->GetElement(24);
334 mix->DefineElement(3, element, 0.002);
335 element = table->GetElement(25);
336 mix->DefineElement(4, element, 0.0115);
337 break;
338 case 4: // EPOXY
339 mix = new TGeoMixture("G10FR4", 4, 1.8);
340 element = table->GetElement(1);
341 mix->DefineElement(0, element, 0.19);
342 element = table->GetElement(6);
343 mix->DefineElement(1, element, 0.18);
344 element = table->GetElement(8);
345 mix->DefineElement(2, element, 0.35);
346 element = table->GetElement(14);
347 mix->DefineElement(3, element, 0.28);
348 break;
349 case 5: // EPOXY
350 mix = new TGeoMixture("G10FR4", 4, 1.8);
351 element = table->GetElement(1);
352 mix->DefineElement(0, element, 0.19);
353 element = table->GetElement(6);
354 mix->DefineElement(1, element, 0.18);
355 element = table->GetElement(8);
356 mix->DefineElement(2, element, 0.35);
357 element = table->GetElement(14);
358 mix->DefineElement(3, element, 0.28);
359 break;
360 case 6: // KAPTON
361 mix = new TGeoMixture("KAPTON", 4, 1.3);
362 element = table->GetElement(1);
363 mix->DefineElement(0, element, 0.026363415);
364 element = table->GetElement(6);
365 mix->DefineElement(1, element, 0.6911272);
366 element = table->GetElement(7);
367 mix->DefineElement(2, element, 0.073271325);
368 element = table->GetElement(8);
369 mix->DefineElement(3, element, 0.209238060);
370 break;
371 case 7: // INOX
372 mix = new TGeoMixture("INOX", 9, 7.9);
373 element = table->GetElement(6);
374 mix->DefineElement(0, element, 0.0003);
375 element = table->GetElement(14);
376 mix->DefineElement(1, element, 0.01);
377 element = table->GetElement(15);
378 mix->DefineElement(2, element, 0.00045);
379 element = table->GetElement(16);
380 mix->DefineElement(3, element, 0.0003);
381 element = table->GetElement(24);
382 mix->DefineElement(4, element, 0.17);
383 element = table->GetElement(25);
384 mix->DefineElement(5, element, 0.02);
385 element = table->GetElement(26);
386 mix->DefineElement(6, element, 0.654);
387 element = table->GetElement(28);
388 mix->DefineElement(7, element, 0.12);
389 element = table->GetElement(42);
390 mix->DefineElement(8, element, 0.025);
391 break;
392 case 8: // ROHACELL
393 mix = new TGeoMixture("ROHACELL", 4, 0.05);
394 element = table->GetElement(1);
395 mix->DefineElement(0, element, 0.07836617);
396 element = table->GetElement(6);
397 mix->DefineElement(1, element, 0.64648941);
398 element = table->GetElement(7);
399 mix->DefineElement(2, element, 0.08376983);
400 element = table->GetElement(8);
401 mix->DefineElement(3, element, 0.19137459);
402 break;
403 case 9: // SDD-C-AL
404 mix = new TGeoMixture("SDD-C-AL", 5, 1.9837);
405 element = table->GetElement(1);
406 mix->DefineElement(0, element, 0.022632);
407 element = table->GetElement(6);
408 mix->DefineElement(1, element, 0.8176579);
409 element = table->GetElement(7);
410 mix->DefineElement(2, element, 0.0093488);
411 element = table->GetElement(8);
412 mix->DefineElement(3, element, 0.0503618);
413 element = table->GetElement(13);
414 mix->DefineElement(4, element, 0.1);
415 break;
416 case 10: // X7R-CAP
417 mix = new TGeoMixture("X7R-CAP", 7, 6.72);
418 element = table->GetElement(8);
419 mix->DefineElement(0, element, 0.085975822);
420 element = table->GetElement(22);
421 mix->DefineElement(1, element, 0.084755042);
422 element = table->GetElement(28);
423 mix->DefineElement(2, element, 0.038244751);
424 element = table->GetElement(29);
425 mix->DefineElement(3, element, 0.009471271);
426 element = table->GetElement(50);
427 mix->DefineElement(4, element, 0.321736471);
428 element = table->GetElement(56);
429 mix->DefineElement(5, element, 0.251639432);
430 element = table->GetElement(82);
431 mix->DefineElement(6, element, 0.2081768);
432 break;
433 case 11: // SDD ruby sph. Al2O3
434 mix = new TGeoMixture("AL2O3", 2, 3.97);
435 element = table->GetElement(8);
436 mix->DefineElement(0, element, 0.5293);
437 element = table->GetElement(13);
438 mix->DefineElement(1, element, 0.4707);
439 break;
440 case 12: // SDD HV microcable
441 mix = new TGeoMixture("HV-CABLE", 5, 1.6087);
442 element = table->GetElement(1);
443 mix->DefineElement(0, element, 0.01983871336);
444 element = table->GetElement(6);
445 mix->DefineElement(1, element, 0.520088819984);
446 element = table->GetElement(7);
447 mix->DefineElement(2, element, 0.0551367996);
448 element = table->GetElement(8);
449 mix->DefineElement(3, element, 0.157399667056);
450 element = table->GetElement(13);
451 mix->DefineElement(4, element, 0.247536);
452 break;
453 case 13: //SDD LV+signal cable
454 mix = new TGeoMixture("LV-CABLE", 5, 2.1035);
455 element = table->GetElement(1);
456 mix->DefineElement(0, element, 0.0082859922);
457 element = table->GetElement(6);
458 mix->DefineElement(1, element, 0.21722436468);
459 element = table->GetElement(7);
460 mix->DefineElement(2, element, 0.023028867);
461 element = table->GetElement(8);
462 mix->DefineElement(3, element, 0.06574077612);
463 element = table->GetElement(13);
464 mix->DefineElement(4, element, 0.68572);
465 break;
466 case 14: //SDD hybrid microcab
467 mix = new TGeoMixture("HYB-CAB", 5, 2.0502);
468 element = table->GetElement(1);
469 mix->DefineElement(0, element, 0.00926228815);
470 element = table->GetElement(6);
471 mix->DefineElement(1, element, 0.24281879711);
472 element = table->GetElement(7);
473 mix->DefineElement(2, element, 0.02574224025);
474 element = table->GetElement(8);
475 mix->DefineElement(3, element, 0.07348667449);
476 element = table->GetElement(13);
477 mix->DefineElement(4, element, 0.64869);
478 break;
479 case 15: //SDD anode microcab
480 mix = new TGeoMixture("ANOD-CAB", 5, 1.7854);
481 element = table->GetElement(1);
482 mix->DefineElement(0, element, 0.0128595919215);
483 element = table->GetElement(6);
484 mix->DefineElement(1, element, 0.392653705471);
485 element = table->GetElement(7);
486 mix->DefineElement(2, element, 0.041626868025);
487 element = table->GetElement(8);
488 mix->DefineElement(3, element, 0.118832707289);
489 element = table->GetElement(13);
490 mix->DefineElement(4, element, 0.431909);
491 break;
492 case 16: // inox/alum
493 mix = new TGeoMixture("INOX-AL", 5, 3.0705);
494 element = table->GetElement(13);
495 mix->DefineElement(0, element, 0.816164);
496 element = table->GetElement(14);
497 mix->DefineElement(1, element, 0.000919182);
498 element = table->GetElement(24);
499 mix->DefineElement(2, element, 0.0330906);
500 element = table->GetElement(26);
501 mix->DefineElement(3, element, 0.131443);
502 element = table->GetElement(28);
503 mix->DefineElement(4, element, 0.0183836);
504 case 17: // MYLAR
505 mix = new TGeoMixture("MYLAR", 3, 1.39);
506 element = table->GetElement(1);
507 mix->DefineElement(0, element, 0.0416667);
508 element = table->GetElement(6);
509 mix->DefineElement(1, element, 0.625);
510 element = table->GetElement(8);
511 mix->DefineElement(2, element, 0.333333);
512 break;
513 case 18: // SPDBUS(AL+KPT+EPOX) - unknown composition
514 mix = new TGeoMixture("SPDBUS", 1, 1.906);
515 element = table->GetElement(9);
516 mix->DefineElement(0, element, 1.);
517 break;
518 case 19: // SDD/SSD rings - unknown composition
519 mix = new TGeoMixture("SDDRINGS", 1, 1.8097);
520 element = table->GetElement(6);
521 mix->DefineElement(0, element, 1.);
522 break;
523 case 20: // SPD end ladder - unknown composition
524 mix = new TGeoMixture("SPDEL", 1, 3.6374);
525 element = table->GetElement(22);
526 mix->DefineElement(0, element, 1.);
527 break;
528 case 21: // SDD end ladder - unknown composition
529 mix = new TGeoMixture("SDDEL", 1, 0.3824);
530 element = table->GetElement(30);
531 mix->DefineElement(0, element, 1.);
532 break;
533 case 22: // SSD end ladder - unknown composition
534 mix = new TGeoMixture("SSDEL", 1, 0.68);
535 element = table->GetElement(16);
536 mix->DefineElement(0, element, 1.);
537 break;
538 }
539 mix->SetZ(z);
540 printf("Patched with mixture %s\n", mix->GetName());
541 return mix;
8495a208 542}
543
8495a208 544//_____________________________________________________________________________
efde9b4d 545void TFlukaMCGeometry::CreateFlukaMatFile(const char *fname)
8495a208 546{
547 // ==== from FLUGG ====
548 // NAMES OF ELEMENTS AND COMPOUNDS: the names must be written in upper case,
549 // according to the fluka standard. In addition,. they must be equal to the
550 // names of the fluka materials - see fluka manual - in order that the
551 // program load the right cross sections, and equal to the names included in
552 // the .pemf. Otherwise the user must define the LOW-MAT CARDS, and make his
553 // own .pemf, in order to get the right cross sections loaded in memory.
89125dc7 554
555
8495a208 556 TString sname;
557 gGeoManager->Export("flgeom.root");
558 if (fname) sname = fname;
559 else sname = "flukaMat.inp";
560 ofstream out;
561 out.open(sname.Data(), ios::out);
562 if (!out.good()) {
563 Fatal("CreateFlukaMatFile", "could not open file %s for writing", sname.Data());
564 return;
565 }
566 PrintHeader(out, "MATERIALS AND COMPOUNDS");
567 PrintHeader(out, "MATERIALS");
89125dc7 568 Int_t i,j,idmat;
569 Int_t counttothree, nelem;
570 Double_t a,z,rho, w;
571 TGeoElementTable *table = TGeoElementTable::Instance();
572 TGeoElement *element;
573 element = table->GetElement(13);
574 element->SetTitle("ALUMINUM"); // this is how FLUKA likes it ...
575 element = table->GetElement(15);
576 element->SetTitle("PHOSPHO"); // same story ...
577// element = table->GetElement(10);
578// element->SetTitle("ARGON"); // NEON not in neutron xsec table
579 Int_t nelements = table->GetNelements();
8495a208 580 TList *matlist = gGeoManager->GetListOfMaterials();
89125dc7 581 TList *medlist = gGeoManager->GetListOfMedia();
582 Int_t nmed = medlist->GetSize();
8495a208 583 TIter next(matlist);
584 Int_t nmater = matlist->GetSize();
585 Int_t nfmater = 0;
89125dc7 586 TGeoMaterial *mat;
8495a208 587 TGeoMixture *mix = 0;
588 TString matname;
89125dc7 589 TObjString *objstr;
590 // Create all needed elements
591 for (Int_t i=1; i<nelements; i++) {
592 element = table->GetElement(i);
593 // skip elements which are not defined
594 if (!element->IsUsed() && !element->IsDefined()) continue;
595 matname = element->GetTitle();
596 ToFlukaString(matname);
597 rho = 0.999;
598 mat = new TGeoMaterial(matname, element->A(), element->Z(), rho);
599 mat->SetIndex(nfmater+3);
600 mat->SetUsed(kTRUE);
601 fMatList->Add(mat);
602 objstr = new TObjString(matname.Data());
603 fMatNames->Add(objstr);
604 nfmater++;
605 }
606 Int_t indmat = nfmater;
607 TGeoMedium *med;
608 // Adjust material names and add them to FLUKA list
8495a208 609 for (i=0; i<nmater; i++) {
610 mat = (TGeoMaterial*)matlist->At(i);
89125dc7 611 if (!mat->IsUsed()) continue;
612 z = mat->GetZ();
613 a = mat->GetA();
614 rho = mat->GetDensity();
615 if (mat->GetZ()<0.001) {
efde9b4d 616 mat->SetIndex(2); // vacuum, built-in inside FLUKA
617 continue;
89125dc7 618 }
619 matname = mat->GetName();
620 FlukaMatName(matname);
621 // material with one element: create it as mixture since it can be duplicated
622 if (!mat->IsMixture()) {
623 // normal materials
624 mix = new TGeoMixture(matname.Data(), 1, rho);
625 mix->DefineElement(0, mat->GetElement(), 1.);
8495a208 626 mat->SetIndex(nfmater+3);
89125dc7 627 for (j=0; j<nmed; j++) {
628 med = (TGeoMedium*)medlist->At(j);
629 if (med->GetMaterial() == mat) {
630 med->SetMaterial(mix);
631 if (mat->GetCerenkovProperties()) {
632 mix->SetCerenkovProperties(mat->GetCerenkovProperties());
633 mat->SetCerenkovProperties(0);
8495a208 634 }
89125dc7 635 break;
636 }
637 }
638 mat = (TGeoMaterial*)mix;
8495a208 639 }
89125dc7 640 mat->SetIndex(nfmater+3);
641 objstr = new TObjString(matname.Data());
642 fMatList->Add(mat);
643 fMatNames->Add(objstr);
644 nfmater++;
645 }
646
647 // Dump all elements with MATERIAL cards
8495a208 648 for (i=0; i<nfmater; i++) {
89125dc7 649 mat = (TGeoMaterial*)fMatList->At(i);
650// mat->SetUsed(kFALSE);
651 mix = 0;
652// out << "* " << mat->GetName() << endl;
8495a208 653 out << setw(10) << "MATERIAL ";
654 out.setf(static_cast<std::ios::fmtflags>(0),std::ios::floatfield);
89125dc7 655 objstr = (TObjString*)fMatNames->At(i);
b0853588 656 matname = objstr->GetString();
89125dc7 657 z = mat->GetZ();
658 a = mat->GetA();
659 rho = mat->GetDensity();
8495a208 660 if (mat->IsMixture()) {
661 out << setw(10) << " ";
662 out << setw(10) << " ";
663 mix = (TGeoMixture*)mat;
664 } else {
89125dc7 665 out << setw(10) << setiosflags(ios::fixed) << setprecision(1) << z;
666 out << setw(10) << setprecision(3) << a;
8495a208 667 }
668 out.setf(static_cast<std::ios::fmtflags>(0),std::ios::floatfield);
89125dc7 669 out << setw(10) << setiosflags(ios::scientific) << setprecision(3) << rho;
8495a208 670 out.setf(static_cast<std::ios::fmtflags>(0),std::ios::floatfield);
89125dc7 671 out << setw(10) << setiosflags(ios::fixed) << setprecision(1) << Double_t(mat->GetIndex());
8495a208 672 out << setw(10) << " ";
efde9b4d 673 out << setw(10) << " ";
674 out << setw(8) << matname.Data() << endl;
89125dc7 675 if (!mix) {
676 // add LOW-MAT card for NEON to associate with ARGON neutron xsec
677 if (z==10 && matname.Contains("NEON")) {
678 out << setw(10) << "LOW-MAT ";
679 out.setf(static_cast<std::ios::fmtflags>(0),std::ios::floatfield);
680 out << setw(10) << setiosflags(ios::fixed) << setprecision(1) << Double_t(mat->GetIndex());
681 out << setw(10) << setiosflags(ios::fixed) << setprecision(1) << 18.;
682 out << setw(10) << setiosflags(ios::fixed) << setprecision(1) << -2.;
683 out << setw(10) << setiosflags(ios::fixed) << setprecision(1) << 293.;
684 out << setw(10) << " ";
685 out << setw(10) << " ";
686// out << setw(8) << matname.Data() << endl;
687 out << setw(8) << " " << endl;
688 }
689 continue;
d11bf3ca 690 }
8495a208 691 counttothree = 0;
692 out << setw(10) << "COMPOUND ";
693 nelem = mix->GetNelements();
89125dc7 694 objstr = (TObjString*)fMatNames->At(i);
8495a208 695 matname = objstr->GetString();
8495a208 696 for (j=0; j<nelem; j++) {
89125dc7 697 w = (mix->GetWmixt())[j];
698 if (w<0.00001) w=0.00001;
699 z = (mix->GetZmixt())[j];
700 a = (mix->GetAmixt())[j];
701 idmat = GetElementIndex(Int_t(z));
702 if (!idmat) Error("CreateFlukaMatFile", "element with Z=%f not found", z);
8495a208 703 out.setf(static_cast<std::ios::fmtflags>(0),std::ios::floatfield);
89125dc7 704 out << setw(10) << setiosflags(ios::fixed) << setprecision(6) << -w;
8495a208 705 out.setf(static_cast<std::ios::fmtflags>(0),std::ios::floatfield);
706 out << setw(10) << setiosflags(ios::fixed) << setprecision(1) << Double_t(idmat);
707 counttothree++;
708 if (counttothree == 3) {
709 out << matname.Data();
710 out << endl;
711 if ( (j+1) != nelem) out << setw(10) << "COMPOUND ";
712 counttothree = 0;
89125dc7 713 }
714 }
8495a208 715 if (nelem%3) {
716 for (j=0; j<(3-(nelem%3)); j++)
717 out << setw(10) << " " << setw(10) << " ";
718 out << matname.Data();
719 out << endl;
89125dc7 720 }
721 }
8495a208 722 Int_t nvols = gGeoManager->GetListOfUVolumes()->GetEntriesFast()-1;
723 TGeoVolume *vol;
8495a208 724 // Now print the material assignments
d11bf3ca 725 Double_t flagfield = 0.;
726 printf("#############################################################\n");
89125dc7 727 if (gFluka->IsFieldEnabled()) {
d11bf3ca 728 flagfield = 1.;
729 printf("Magnetic field enabled\n");
730 } else printf("Magnetic field disabled\n");
731 printf("#############################################################\n");
732
8495a208 733 PrintHeader(out, "TGEO MATERIAL ASSIGNMENTS");
734 for (i=1; i<=nvols; i++) {
735 vol = gGeoManager->GetVolume(i);
736 mat = vol->GetMedium()->GetMaterial();
89125dc7 737// mat->SetUsed(kTRUE);
8495a208 738 idmat = mat->GetIndex();
89125dc7 739 for (Int_t j=0; j<nfmater; j++) {
740 mat = (TGeoMaterial*)fMatList->At(j);
741 if (mat->GetIndex() == idmat) mat->SetUsed(kTRUE);
742 }
8495a208 743 out << setw(10) << "ASSIGNMAT ";
744 out.setf(static_cast<std::ios::fmtflags>(0),std::ios::floatfield);
745 out << setw(10) << setiosflags(ios::fixed) << Double_t(idmat);
746 out << setw(10) << setiosflags(ios::fixed) << Double_t(i);
747 out << setw(10) << "0.0";
b0853588 748 out << setw(10) << "0.0";
8495a208 749 out << setw(10) << setiosflags(ios::fixed) << flagfield;
b0853588 750 out << setw(10) << "0.0";
8495a208 751 out << endl;
752 }
8495a208 753 out.close();
efde9b4d 754 fLastMaterial = nfmater+2;
89125dc7 755/*
756 TGeoMaterial *mat1 = 0;
757 for (i=1; i<=nvols; i++) {
758 vol = gGeoManager->GetVolume(i);
759 med = vol->GetMedium();
760 mat = med->GetMaterial();
761 printf("Region %d: %s\n", i, vol->GetName());
762 printf(" medium %d: %s\n", med->GetId(), med->GetName());
763 for (j=0; j<nfmater; j++) {
764 mat1 = (TGeoMaterial*)fMatList->At(j);
765 if (mat1 != mat) continue;
766 objstr = (TObjString*)fMatNames->At(j);
767 matname = objstr->GetString();
768 break;
769 }
770 if (mat1 != mat) printf(" (*) material not found in Fluka list\n");
771 printf(" material %s (at ind=%d): FlukaID=%d FlukaName=%s\n",
772 mat->GetName(), j, mat->GetIndex(), matname.Data());
773 if (mat->GetCerenkovProperties()) printf(" Cerenkov properties found\n");
774 }
775*/
776
777 if (!gFluka->IsGeneratePemf()) return;
778 // Write peg files
779 char number[20];
780 for (i=indmat; i<nfmater; i++) {
781 mat = (TGeoMaterial*)fMatList->At(i);
782 if (!mat->IsUsed()) continue;
783 sname = "mat";
784 sprintf(number, "%d", i);
785 sname.Append(number);
786 WritePegFile(i);
787 sname.Prepend("$FLUPRO/pemf/rpemf peg/");
788 gSystem->Exec(sname.Data());
789 }
790 sname = "cat peg/*.pemf > peg/alice.pemf";
791 gSystem->Exec(sname.Data());
792 sname = "mv peg/alice.pemf alice.pemf";
793 gSystem->Exec(sname.Data());
794}
795
796//_____________________________________________________________________________
797void TFlukaMCGeometry::WritePegFile(Int_t imat) const
798{
799// Write the .peg file for one material
800 TGeoMaterial *mat = (TGeoMaterial*)fMatList->At(imat);
801 TString name = ((TObjString*)fMatNames->At(imat))->GetString();
802 TString line;
803 char number[20];
804 TGeoElement *elem = mat->GetElement();
805 name = name.Strip();
806 TString sname = "mat";
807 sprintf(number, "%d", imat);
808 sname.Append(number);
809 sname.Append(".peg");
810 sname.Prepend("peg/");
811 ofstream out;
812 out.open(sname.Data(), ios::out);
813 if (!out.good()) return;
814 Double_t dens = mat->GetDensity();
815 TGeoMixture *mix = 0;
816 Int_t nel = 1;
817 Int_t i;
818 if (mat->IsMixture()) {
819 mix = (TGeoMixture*)mat;
820 nel = mix->GetNelements();
821 }
822 if (nel==1) {
823 out << "ELEM" << endl;
824 out << " &INP IRAYL=1, RHO=" << dens << ", " << endl;
825 if (dens<0.01) out << " GASP=1." << endl;
826 out << " &END" << endl;
827 out << name.Data() << endl;
828 out << elem->GetName() << endl;
829 } else {
830 out << "MIXT" << endl;
831 out << " &INP IRAYL=1, NE=" << nel << ", RHOZ=";
832 line = "";
833 for (i=0; i<nel; i++) {
834 sprintf(number, "%f", mix->GetWmixt()[i]);
835 line += number;
836 line += ", ";
837 if (line.Length() > 30) {
838 out << line.Data() << endl;
839 line = "";
840 }
841 }
842 if (line.Length()) out << " " << line.Data() << endl;
843 out << " RHO=" << dens;
844 if (dens<0.01) out << ", GASP=1." << endl;
845 out << " &END" << endl;
846 out << name.Data() << endl;
847 for (i=0; i<nel; i++) {
848 elem = mix->GetElement(i);
849 line = elem->GetName();
850 if (line.Length()==1) line.Append(" ");
851 out << line.Data() << " ";
852 }
853 out << endl;
854 }
855 out << "ENER" << endl;
856 out << " $INP AE=0.56099906, UE=3000000., AP=.03, UP=3000000. $END" << endl;
857 out << "PWLF" << endl;
858 out << " $INP NALE=300, NALG=400, NALR=100 $END" << endl;
859 out << "DECK" << endl;
860 out << " $INP $END" << endl;
861 out << "TEST" << endl;
862 out << " $INP $END" << endl;
863 out.close();
8495a208 864}
865
866//_____________________________________________________________________________
867void TFlukaMCGeometry::PrintHeader(ofstream &out, const char *text) const
868{
869// Print a FLUKA header.
870 out << "*\n" << "*\n" << "*\n";
871 out << "********************* " << text << " *********************\n"
872 << "*\n";
873 out << "*...+....1....+....2....+....3....+....4....+....5....+....6....+....7..."
874 << endl;
875 out << "*" << endl;
876}
877
878//_____________________________________________________________________________
879Int_t TFlukaMCGeometry::RegionId() const
880{
881// Returns current region id <-> TGeo node id
882 if (gGeoManager->IsOutside()) return 0;
883 return gGeoManager->GetCurrentNode()->GetUniqueID();
884}
89125dc7 885
886//_____________________________________________________________________________
887Int_t TFlukaMCGeometry::GetElementIndex(Int_t z) const
888{
6f5667d1 889// Get index of a material having a given Z element.
89125dc7 890 TIter next(fMatList);
891 TGeoMaterial *mat;
892 Int_t index = 0;
893 while ((mat=(TGeoMaterial*)next())) {
894 if (mat->IsMixture()) continue;
895 if (mat->GetElement()->Z() == z) return mat->GetIndex();
896 }
897 return index;
898}
899
05265ca9 900//_____________________________________________________________________________
901void TFlukaMCGeometry::SetMreg(Int_t mreg)
902{
903// Update if needed next history;
89125dc7 904 if (gFluka->GetDummyBoundary()==2) {
d11bf3ca 905 gGeoManager->CdNode(fNextLattice-1);
906 return;
907 }
89125dc7 908 Int_t curreg = (gGeoManager->IsOutside())?(gMCGeom->NofVolumes()+1):gGeoManager->GetCurrentVolume()->GetNumber();
05265ca9 909 if (mreg==curreg) return;
910 if (mreg==fNextRegion) {
911 if (fNextLattice!=999999999) gGeoManager->CdNode(fNextLattice-1);
912 return;
d11bf3ca 913 } else {
914 if (mreg == fCurrentRegion) {
915 if (fCurrentLattice!=999999999) gGeoManager->CdNode(fCurrentLattice-1);
916 return;
917 }
918 }
2bc4c610 919 if (fDebug) printf("ERROR: mreg=%i neither current nor next region\n", mreg);
05265ca9 920}
921
d11bf3ca 922//_____________________________________________________________________________
923void TFlukaMCGeometry::SetCurrentRegion(Int_t mreg, Int_t latt)
924{
925// Set index/history for next entered region
926 fCurrentRegion = mreg;
927 fCurrentLattice = latt;
928}
929
05265ca9 930//_____________________________________________________________________________
931void TFlukaMCGeometry::SetNextRegion(Int_t mreg, Int_t latt)
932{
933// Set index/history for next entered region
934 fNextRegion = mreg;
935 fNextLattice = latt;
936}
8495a208 937
938//_____________________________________________________________________________
939void TFlukaMCGeometry::ToFlukaString(TString &str) const
940{
941// ToFlukaString converts an string to something usefull in FLUKA:
942// * Capital letters
943// * Only 8 letters
944// * Replace ' ' by '_'
945 if (str.Length()<8) {
946 str += " ";
947 }
948 str.Remove(8);
949 Int_t ilast;
950 for (ilast=7; ilast>0; ilast--) if (str(ilast)!=' ') break;
951 str.ToUpper();
952 for (Int_t pos=0; pos<ilast; pos++)
953 if (str(pos)==' ') str.Replace(pos,1,"_",1);
954 return;
955}
89125dc7 956
957//_____________________________________________________________________________
958void TFlukaMCGeometry::FlukaMatName(TString &str) const
959{
6f5667d1 960// Convert a name to upper case 8 chars.
89125dc7 961 ToFlukaString(str);
962 Int_t ilast;
963 for (ilast=7; ilast>0; ilast--) if (str(ilast)!=' ') break;
964 if (ilast>5) ilast = 5;
965 char number[3];
966 TIter next(fMatNames);
967 TObjString *objstr;
968 TString matname;
969 Int_t index = 0;
970 while ((objstr=(TObjString*)next())) {
971 matname = objstr->GetString();
972 if (matname == str) {
973 index++;
974 if (index<10) {
975 number[0] = '0';
976 sprintf(&number[1], "%d", index);
977 } else if (index<100) {
978 sprintf(number, "%d", index);
979 } else {
980 Error("FlukaMatName", "Too many materials %s", str.Data());
981 return;
982 }
983 str.Replace(ilast+1, 2, number);
984 str.Remove(8);
985 }
986 }
987}
988
8495a208 989//______________________________________________________________________________
990void TFlukaMCGeometry::Vname(const char *name, char *vname) const
991{
992 //
993 // convert name to upper case. Make vname at least 4 chars
994 //
995 Int_t l = strlen(name);
996 Int_t i;
997 l = l < 4 ? l : 4;
998 for (i=0;i<l;i++) vname[i] = toupper(name[i]);
999 for (i=l;i<4;i++) vname[i] = ' ';
1000 vname[4] = 0;
1001}
1002
1003
1004// FLUKA GEOMETRY WRAPPERS - to replace FLUGG wrappers
1005
8495a208 1006//_____________________________________________________________________________
1007Int_t idnrwr(const Int_t & /*nreg*/, const Int_t & /*mlat*/)
1008{
1009// from FLUGG:
1010// Wrapper for setting DNEAR option on fluka side. Must return 0
1011// if user doesn't want Fluka to use DNEAR to compute the
1012// step (the same effect is obtained with the GLOBAL (WHAT(3)=-1)
1013// card in fluka input), returns 1 if user wants Fluka always to
1014// use DNEAR (in this case, be sure that GEANT4 DNEAR is unique,
1015// coming from all directions!!!)
89125dc7 1016 if (gMCGeom->IsDebugging()) printf("========== Dummy IDNRWR\n");
8495a208 1017 return 0;
1018}
1019
8495a208 1020//_____________________________________________________________________________
1021void g1wr(Double_t &pSx, Double_t &pSy, Double_t &pSz,
d11bf3ca 1022 Double_t *pV, Int_t &oldReg , const Int_t &oldLttc, Double_t &propStep,
1023 Int_t &nascFlag, Double_t &retStep, Int_t &newReg,
2573ac89 1024 Double_t &saf, Int_t &newLttc, Int_t &lttcFlag,
8495a208 1025 Double_t *sLt, Int_t *jrLt)
d11bf3ca 1026
8495a208 1027{
d11bf3ca 1028 // Initialize FLUKa point and direction;
6f5667d1 1029 gNstep++;
d11bf3ca 1030/*
1031 if (kNstep>0) {
89125dc7 1032 gMCGeom->SetDebugMode(kTRUE);
1033 gFluka->SetVerbosityLevel(3);
d11bf3ca 1034 }
1035 if (kNstep>6520) {
89125dc7 1036 gMCGeom->SetDebugMode(kFALSE);
1037 gFluka->SetVerbosityLevel(0);
d11bf3ca 1038 }
1039 if ((kNstep%10)==0) printf("step %i\n", kNstep);
1040*/
8495a208 1041
89125dc7 1042 if (gMCGeom->IsDebugging()) {
2bc4c610 1043 printf("========== Inside G1WR\n");
1044 printf(" point/dir:(%14.9f, %14.9f, %14.9f, %g, %g, %g)\n", pSx,pSy,pSz,pV[0],pV[1],pV[2]);
d11bf3ca 1045 printf(" oldReg=%i oldLttc=%i pstep=%f\n",oldReg, oldLttc, propStep);
2bc4c610 1046 }
8495a208 1047 gGeoManager->SetCurrentPoint(pSx, pSy, pSz);
1048 gGeoManager->SetCurrentDirection(pV);
89125dc7 1049 gMCGeom->SetCurrentRegion(oldReg, oldLttc);
d11bf3ca 1050 // Initialize default return values
1051 lttcFlag = 0;
1052 jrLt[lttcFlag] = oldLttc;
1053 sLt[lttcFlag] = propStep;
1054 jrLt[lttcFlag+1] = -1;
1055 sLt[lttcFlag+1] = 0.;
1056 newReg = oldReg;
1057 newLttc = oldLttc;
1058 // check if dummy boundary flag is set
1059 Int_t curLttc, curReg;
89125dc7 1060 if (gFluka->IsDummyBoundary()) {
d11bf3ca 1061 // printf("Dummy boundary intercepted. Point is: %f, %f, %f\n", pSx, pSy, pSz);
1062 Bool_t crossedDummy = (oldLttc == TFlukaMCGeometry::kLttcVirtual)?kTRUE:kFALSE;
1063 if (crossedDummy) {
1064 // FLUKA crossed the dummy boundary - update new region/history
1065 retStep = 0.;
1066 saf = 0.;
89125dc7 1067 gMCGeom->GetNextRegion(newReg, newLttc);
1068 gMCGeom->SetMreg(newReg);
1069 if (gMCGeom->IsDebugging()) printf(" virtual newReg=%i newLttc=%i\n", newReg, newLttc);
d11bf3ca 1070 sLt[lttcFlag] = 0.; // null step in current region
1071 lttcFlag++;
1072 jrLt[lttcFlag] = newLttc;
1073 sLt[lttcFlag] = 0.; // null step in next region
1074 jrLt[lttcFlag+1] = -1;
1075 sLt[lttcFlag+1] = 0.;
89125dc7 1076 gFluka->SetDummyBoundary(0);
d11bf3ca 1077 return;
1078 }
1079 }
1080
1081 // Reset outside flag
05265ca9 1082 if (gGeoManager->IsOutside()) {
1083 gGeoManager->SetOutside(kFALSE);
1084 gGeoManager->CdTop();
d11bf3ca 1085 }
1086
1087 // Reset dummy boundary flag
89125dc7 1088 gFluka->SetDummyBoundary(0);
d11bf3ca 1089
1090 curLttc = gGeoManager->GetCurrentNodeId()+1;
1091 curReg = gGeoManager->GetCurrentVolume()->GetNumber();
2573ac89 1092 if (oldLttc != curLttc) {
d11bf3ca 1093 // FLUKA crossed the boundary : we trust that the given point is really there,
1094 // so we just update TGeo state
2573ac89 1095 gGeoManager->CdNode(oldLttc-1);
1096 curLttc = gGeoManager->GetCurrentNodeId()+1;
d11bf3ca 1097 curReg = gGeoManager->GetCurrentVolume()->GetNumber();
89125dc7 1098 if (gMCGeom->IsDebugging()) printf(" re-initialized point: curReg=%i curLttc=%i\n", curReg, curLttc);
d11bf3ca 1099 }
1100 // Now the current TGeo state reflects the FLUKA state
89125dc7 1101 if (gMCGeom->IsDebugging()) printf(" current path: %s\n", gGeoManager->GetPath());
d11bf3ca 1102 Double_t extra = 1E-6;
1103 Double_t tmpStep = propStep + extra;
1104 gGeoManager->FindNextBoundary(-tmpStep);
1105 Double_t snext = gGeoManager->GetStep();
1106 // !!!!!
1107 if (snext<=0) {
1108 // FLUKA is in the wrong region, notify it
89125dc7 1109 if (gMCGeom->IsDebugging()) printf("ERROR: snext=%f\n", snext);
d11bf3ca 1110// newReg = -3;
1111// return;
1112 snext = extra;
1113 }
1114 saf = gGeoManager->GetSafeDistance();
1115 Bool_t cross = kFALSE;
1116 Bool_t onBound = kFALSE;
1117 if (snext<tmpStep) {
1118 // We have some boundary in the way
1119 Double_t dd = snext-propStep;
1120 if (dd < 0) {
1121 cross = kTRUE;
1122 dd = -dd;
1123 }
1124 if (dd < 1E-8) onBound = kTRUE;
1125 }
1126 snext += 1.E-8;
89125dc7 1127 if (gMCGeom->IsDebugging()) {
d11bf3ca 1128 if (!cross) printf(" physical step approved: %f\n", propStep);
1129 else printf(" boundary crossing at: %f\n", snext);
1130 if (onBound) printf(" step on boundary limit ! NASC=%i\n", nascFlag);
1131 }
1132 if (!cross) {
1133 // Next boundary further than proposed step, which is approved
1134 retStep = propStep;
1135 sLt[lttcFlag] = propStep;
1136 return;
1137 }
1138 // The next boundary is closer. We try to cross it.
2573ac89 1139 Double_t *point = gGeoManager->GetCurrentPoint();
1140 Double_t *dir = gGeoManager->GetCurrentDirection();
d11bf3ca 1141 Double_t pt[3];
1142 memcpy(pt, point, 3*sizeof(Double_t));
1143
2573ac89 1144 Int_t i;
d11bf3ca 1145 for (i=0;i<3;i++) point[i] += snext*dir[i];
1146 gGeoManager->FindNode();
1147 newLttc = (gGeoManager->IsOutside())?(TFlukaMCGeometry::kLttcOutside):gGeoManager->GetCurrentNodeId()+1;
1148 if (newLttc == oldLttc) {
1149 // brute force ...
1150 // Just try a fast extra small step
1151 snext += 1E-6;
1152 for (i=0;i<3;i++) point[i] = pt[i]+snext*dir[i];
1153 gGeoManager->FindNode();
1154 newLttc = (gGeoManager->IsOutside())?(TFlukaMCGeometry::kLttcOutside):gGeoManager->GetCurrentNodeId()+1;
1155 if (newLttc == oldLttc) {
1156 // check if state changes at the end of the proposed step
1157 for (i=0;i<3;i++) point[i] = pt[i]+propStep*dir[i];
2573ac89 1158 gGeoManager->FindNode();
d11bf3ca 1159 newLttc = (gGeoManager->IsOutside())?(TFlukaMCGeometry::kLttcOutside):gGeoManager->GetCurrentNodeId()+1;
1160 if (newLttc==oldLttc) {
1161 // approve step
1162 retStep = propStep;
1163 sLt[lttcFlag] = propStep;
1164 return;
05265ca9 1165 }
d11bf3ca 1166 // snext is underestimated - we will create a virtual one to overcome the error
1167// printf("some boundary in the way...\n");
1168 }
1169 }
1170 gGeoManager->SetCurrentPoint(pt);
1171// newLttc = (gGeoManager->IsOutside())?(TFlukaMCGeometry::kLttcOutside):gGeoManager->GetCurrentNodeId()+1;
89125dc7 1172 newReg = (gGeoManager->IsOutside())?(gMCGeom->NofVolumes()+1):gGeoManager->GetCurrentVolume()->GetNumber();
1173 if (gMCGeom->IsDebugging()) printf(" newReg=%i newLttc=%i\n", newReg, newLttc);
d11bf3ca 1174
1175 // We really crossed the boundary, but is it the same region ?
89125dc7 1176 gMCGeom->SetNextRegion(newReg, newLttc);
d11bf3ca 1177 if (newReg == oldReg) {
1178 // Virtual boundary between replicants
89125dc7 1179 if (gMCGeom->IsDebugging()) printf(" DUMMY boundary\n");
d11bf3ca 1180 newReg = 1; // cheat FLUKA telling it it crossed the TOP region
1181 newLttc = TFlukaMCGeometry::kLttcVirtual;
1182 // mark that next boundary is virtual
89125dc7 1183 gFluka->SetDummyBoundary(1);
d11bf3ca 1184 }
1185 retStep = snext;
1186 sLt[lttcFlag] = snext;
1187 lttcFlag++;
1188 jrLt[lttcFlag] = newLttc;
1189 sLt[lttcFlag] = snext;
1190 jrLt[lttcFlag+1] = -1;
1191 sLt[lttcFlag+1] = 0.;
1192
2573ac89 1193 if (newLttc!=oldLttc) {
1194 if (gGeoManager->IsOutside()) {
1195 gGeoManager->SetOutside(kFALSE);
1196 gGeoManager->CdTop();
d11bf3ca 1197 }
1198 gGeoManager->CdTop();
6f5667d1 1199 if (!gGeoManager->GetCurrentMatrix()->IsIdentity()) printf("ERROR at step %i\n", gNstep);
2573ac89 1200 gGeoManager->CdNode(oldLttc-1);
8495a208 1201 }
89125dc7 1202 if (gMCGeom->IsDebugging()) {
d11bf3ca 1203 printf("=> snext=%g safe=%g\n", snext, saf);
1204 for (Int_t i=0; i<lttcFlag+1; i++) printf(" jrLt[%i]=%i sLt[%i]=%g\n", i,jrLt[i],i,sLt[i]);
1205 }
89125dc7 1206 if (gMCGeom->IsDebugging()) printf("<= G1WR (in: %s)\n", gGeoManager->GetPath());
8495a208 1207}
1208
efde9b4d 1209//_____________________________________________________________________________
1210void g1rtwr()
1211{
89125dc7 1212 if (gMCGeom->IsDebugging()) printf("========== Dummy G1RTWR\n");
efde9b4d 1213}
1214
1215//_____________________________________________________________________________
1216void conhwr(Int_t & /*intHist*/, Int_t * /*incrCount*/)
1217{
89125dc7 1218 if (gMCGeom->IsDebugging()) printf("========== Dummy CONHWR\n");
efde9b4d 1219}
1220
1221//_____________________________________________________________________________
2573ac89 1222void inihwr(Int_t &intHist)
efde9b4d 1223{
89125dc7 1224 if (gMCGeom->IsDebugging()) printf("========== Inside INIHWR -> reinitializing history: %i\n", intHist);
2573ac89 1225 if (gGeoManager->IsOutside()) gGeoManager->CdTop();
1226 if (intHist<=0) {
1227// printf("=== wrong history number\n");
1228 return;
1229 }
1230 if (intHist==0) gGeoManager->CdTop();
1231 else gGeoManager->CdNode(intHist-1);
89125dc7 1232 if (gMCGeom->IsDebugging()) {
2bc4c610 1233 printf(" --- current path: %s\n", gGeoManager->GetPath());
1234 printf("<= INIHWR\n");
1235 }
efde9b4d 1236}
1237
1238//_____________________________________________________________________________
1239void jomiwr(const Int_t & /*nge*/, const Int_t & /*lin*/, const Int_t & /*lou*/,
1240 Int_t &flukaReg)
1241{
1242// Geometry initialization wrapper called by FLUKAM. Provides to FLUKA the
1243// number of regions (volumes in TGeo)
1244 // build application geometry
89125dc7 1245 if (gMCGeom->IsDebugging()) printf("========== Inside JOMIWR\n");
2573ac89 1246 flukaReg = gGeoManager->GetListOfUVolumes()->GetEntriesFast();
89125dc7 1247 if (gMCGeom->IsDebugging()) printf("<= JOMIWR: last region=%i\n", flukaReg);
efde9b4d 1248}
1249
1250//_____________________________________________________________________________
1251void lkdbwr(Double_t &pSx, Double_t &pSy, Double_t &pSz,
1252 Double_t * /*pV*/, const Int_t &oldReg, const Int_t &oldLttc,
1253 Int_t &newReg, Int_t &flagErr, Int_t &newLttc)
1254{
89125dc7 1255 if (gMCGeom->IsDebugging()) {
2bc4c610 1256 printf("========== Inside LKDBWR (%f, %f, %f)\n",pSx, pSy, pSz);
1257// printf(" in: pV=(%f, %f, %f)\n", pV[0], pV[1], pV[2]);
1258 printf(" in: oldReg=%i oldLttc=%i\n", oldReg, oldLttc);
1259 }
efde9b4d 1260 TGeoNode *node = gGeoManager->FindNode(pSx, pSy, pSz);
1261 if (gGeoManager->IsOutside()) {
89125dc7 1262 newReg = gMCGeom->NofVolumes()+1;
05265ca9 1263// newLttc = gGeoManager->GetCurrentNodeId();
1264 newLttc = 999999999;
89125dc7 1265 if (gMCGeom->IsDebugging()) {
2bc4c610 1266 printf("OUTSIDE\n");
1267 printf(" out: newReg=%i newLttc=%i\n", newReg, newLttc);
1268 printf("<= LKMGWR\n");
1269 }
05265ca9 1270 flagErr = newReg;
1271 return;
efde9b4d 1272 }
1273 newReg = node->GetVolume()->GetNumber();
1274 newLttc = gGeoManager->GetCurrentNodeId()+1;
89125dc7 1275 gMCGeom->SetNextRegion(newReg, newLttc);
efde9b4d 1276 flagErr = newReg;
89125dc7 1277 if (gMCGeom->IsDebugging()) {
2bc4c610 1278 printf(" out: newReg=%i newLttc=%i\n", newReg, newLttc);
1279 printf("<= LKDBWR\n");
1280 }
efde9b4d 1281}
1282
1283//_____________________________________________________________________________
1284void lkfxwr(Double_t &pSx, Double_t &pSy, Double_t &pSz,
1285 Double_t * /*pV*/, const Int_t &oldReg, const Int_t &oldLttc,
1286 Int_t &newReg, Int_t &flagErr, Int_t &newLttc)
1287{
89125dc7 1288 if (gMCGeom->IsDebugging()) {
2bc4c610 1289 printf("========== Inside LKFXWR (%f, %f, %f)\n",pSx, pSy, pSz);
1290// printf(" in: pV=(%f, %f, %f)\n", pV[0], pV[1], pV[2]);
1291 printf(" in: oldReg=%i oldLttc=%i\n", oldReg, oldLttc);
1292 }
efde9b4d 1293 TGeoNode *node = gGeoManager->FindNode(pSx, pSy, pSz);
1294 if (gGeoManager->IsOutside()) {
89125dc7 1295 newReg = gMCGeom->NofVolumes()+1;
05265ca9 1296// newLttc = gGeoManager->GetCurrentNodeId();
1297 newLttc = 999999999;
89125dc7 1298 if (gMCGeom->IsDebugging()) {
2bc4c610 1299 printf("OUTSIDE\n");
1300 printf(" out: newReg=%i newLttc=%i\n", newReg, newLttc);
1301 printf("<= LKMGWR\n");
1302 }
05265ca9 1303 flagErr = newReg;
1304 return;
efde9b4d 1305 }
1306 newReg = node->GetVolume()->GetNumber();
1307 newLttc = gGeoManager->GetCurrentNodeId()+1;
89125dc7 1308 gMCGeom->SetNextRegion(newReg, newLttc);
efde9b4d 1309 flagErr = newReg;
89125dc7 1310 if (gMCGeom->IsDebugging()) {
2bc4c610 1311 printf(" out: newReg=%i newLttc=%i\n", newReg, newLttc);
1312 printf("<= LKFXWR\n");
1313 }
efde9b4d 1314}
1315
1316//_____________________________________________________________________________
1317void lkmgwr(Double_t &pSx, Double_t &pSy, Double_t &pSz,
1318 Double_t * /*pV*/, const Int_t &oldReg, const Int_t &oldLttc,
1319 Int_t &flagErr, Int_t &newReg, Int_t &newLttc)
1320{
89125dc7 1321 if (gMCGeom->IsDebugging()) {
2bc4c610 1322 printf("========== Inside LKMGWR (%f, %f, %f)\n",pSx, pSy, pSz);
1323// printf(" in: pV=(%f, %f, %f)\n", pV[0], pV[1], pV[2]);
1324 printf(" in: oldReg=%i oldLttc=%i\n", oldReg, oldLttc);
1325 }
efde9b4d 1326 TGeoNode *node = gGeoManager->FindNode(pSx, pSy, pSz);
1327 if (gGeoManager->IsOutside()) {
89125dc7 1328 newReg = gMCGeom->NofVolumes()+1;
05265ca9 1329// newLttc = gGeoManager->GetCurrentNodeId();
1330 newLttc = 999999999;
89125dc7 1331 if (gMCGeom->IsDebugging()) {
2bc4c610 1332 printf("OUTSIDE\n");
1333 printf(" out: newReg=%i newLttc=%i\n", newReg, newLttc);
1334 printf("<= LKMGWR\n");
1335 }
05265ca9 1336 flagErr = newReg;
1337 return;
efde9b4d 1338 }
1339 newReg = node->GetVolume()->GetNumber();
1340 newLttc = gGeoManager->GetCurrentNodeId()+1;
89125dc7 1341 gMCGeom->SetNextRegion(newReg, newLttc);
efde9b4d 1342 flagErr = newReg;
89125dc7 1343 if (gMCGeom->IsDebugging()) {
2bc4c610 1344 printf(" out: newReg=%i newLttc=%i\n", newReg, newLttc);
1345 printf("<= LKMGWR\n");
1346 }
efde9b4d 1347}
1348
1349//_____________________________________________________________________________
1350void lkwr(Double_t &pSx, Double_t &pSy, Double_t &pSz,
1351 Double_t * /*pV*/, const Int_t &oldReg, const Int_t &oldLttc,
1352 Int_t &newReg, Int_t &flagErr, Int_t &newLttc)
1353{
89125dc7 1354 if (gMCGeom->IsDebugging()) {
2bc4c610 1355 printf("========== Inside LKWR (%f, %f, %f)\n",pSx, pSy, pSz);
1356// printf(" in: pV=(%f, %f, %f)\n", pV[0], pV[1], pV[2]);
1357 printf(" in: oldReg=%i oldLttc=%i\n", oldReg, oldLttc);
1358 }
efde9b4d 1359 TGeoNode *node = gGeoManager->FindNode(pSx, pSy, pSz);
1360 if (gGeoManager->IsOutside()) {
89125dc7 1361 newReg = gMCGeom->NofVolumes()+1;
05265ca9 1362// newLttc = gGeoManager->GetCurrentNodeId();
1363 newLttc = 999999999;
89125dc7 1364 if (gMCGeom->IsDebugging()) {
2bc4c610 1365 printf("OUTSIDE\n");
1366 printf(" out: newReg=%i newLttc=%i\n", newReg, newLttc);
1367 printf("<= LKMGWR\n");
1368 }
05265ca9 1369 flagErr = newReg;
1370 return;
efde9b4d 1371 }
1372 newReg = node->GetVolume()->GetNumber();
d11bf3ca 1373 newLttc = gGeoManager->GetCurrentNodeId()+1;
89125dc7 1374 gMCGeom->SetNextRegion(newReg, newLttc);
efde9b4d 1375 flagErr = newReg;
89125dc7 1376 if (gMCGeom->IsDebugging()) {
2bc4c610 1377 printf(" out: newReg=%i newLttc=%i in %s\n", newReg, newLttc, gGeoManager->GetPath());
1378 printf("<= LKWR\n");
1379 }
efde9b4d 1380}
1381
1382//_____________________________________________________________________________
2573ac89 1383void nrmlwr(Double_t &pSx, Double_t &pSy, Double_t &pSz,
1384 Double_t &pVx, Double_t &pVy, Double_t &pVz,
1385 Double_t *norml, const Int_t &oldReg,
1386 const Int_t &newReg, Int_t &flagErr)
efde9b4d 1387{
89125dc7 1388 if (gMCGeom->IsDebugging()) {
2bc4c610 1389 printf("========== Inside NRMLWR (%g, %g, %g, %g, %g, %g)\n", pSx,pSy,pSz,pVx,pVy,pVz);
1390 printf(" oldReg=%i, newReg=%i\n", oldReg,newReg);
1391 }
89125dc7 1392// Int_t curreg = (gGeoManager->IsOutside())?(gMCGeom->NofVolumes()+1):gGeoManager->GetCurrentVolume()->GetNumber();
d11bf3ca 1393// Int_t curLttc = gGeoManager->GetCurrentNodeId()+1;
89125dc7 1394// if (gMCGeom->IsDebugging()) printf(" curReg=%i, curLttc=%i in: %s\n", curreg, curLttc, gGeoManager->GetPath());
d11bf3ca 1395// Bool_t regsame = (curreg==oldReg)?kTRUE:kFALSE;
2573ac89 1396 gGeoManager->SetCurrentPoint(pSx, pSy, pSz);
1397 gGeoManager->SetCurrentDirection(pVx,pVy,pVz);
d11bf3ca 1398/*
2573ac89 1399 if (!regsame) {
89125dc7 1400 if (gMCGeom->IsDebugging()) printf(" REGIONS DOEN NOT MATCH\n");
2573ac89 1401 gGeoManager->FindNode();
89125dc7 1402 curreg = (gGeoManager->IsOutside())?(gMCGeom->NofVolumes()+1):gGeoManager->GetCurrentVolume()->GetNumber();
2573ac89 1403 curLttc = gGeoManager->GetCurrentNodeId()+1;
89125dc7 1404 if (gMCGeom->IsDebugging()) printf(" re-initialized point: curReg=%i curLttc=%i curPath=%s\n", curreg, curLttc, gGeoManager->GetPath());
2573ac89 1405 }
d11bf3ca 1406*/
2573ac89 1407 Double_t *dnorm = gGeoManager->FindNormalFast();
1408 flagErr = 0;
1409 if (!dnorm) {
1410 printf(" ERROR: Cannot compute fast normal\n");
1411 flagErr = 1;
1412 norml[0] = -pVx;
1413 norml[1] = -pVy;
1414 norml[2] = -pVz;
1415 }
1416 norml[0] = -dnorm[0];
1417 norml[1] = -dnorm[1];
1418 norml[2] = -dnorm[2];
89125dc7 1419 if (gMCGeom->IsDebugging()) printf(" normal to boundary: (%g, %g, %g)\n", norml[0], norml[1], norml[2]);
1420// curreg = (gGeoManager->IsOutside())?(gMCGeom->NofVolumes()+1):gGeoManager->GetCurrentVolume()->GetNumber();
d11bf3ca 1421// curLttc = gGeoManager->GetCurrentNodeId()+1;
89125dc7 1422 if (gMCGeom->IsDebugging()) {
d11bf3ca 1423// printf(" final location: curReg=%i, curLttc=%i in %s\n", curreg,curLttc,gGeoManager->GetPath());
2bc4c610 1424 printf("<= NRMLWR\n");
1425 }
efde9b4d 1426}
1427
1428//_____________________________________________________________________________
1429void rgrpwr(const Int_t & /*flukaReg*/, const Int_t & /*ptrLttc*/, Int_t & /*g4Reg*/,
1430 Int_t * /*indMother*/, Int_t * /*repMother*/, Int_t & /*depthFluka*/)
1431{
89125dc7 1432 if (gMCGeom->IsDebugging()) printf("=> Dummy RGRPWR\n");
efde9b4d 1433}
1434
1435//_____________________________________________________________________________
1436Int_t isvhwr(const Int_t &check, const Int_t & intHist)
1437{
1438// from FLUGG:
1439// Wrapper for saving current navigation history (fCheck=default)
1440// and returning its pointer. If fCheck=-1 copy of history pointed
1441// by intHist is made in NavHistWithCount object, and its pointer
1442// is returned. fCheck=1 and fCheck=2 cases are only in debugging
1443// version: an array is created by means of FGeometryInit functions
1444// (but could be a static int * ptrArray = new int[10000] with
1445// file scope as well) that stores a flag for deleted/undeleted
1446// histories and at the end of event is checked to verify that
1447// all saved history objects have been deleted.
1448
1449// For TGeo, just return the current node ID. No copy need to be made.
1450
89125dc7 1451 if (gMCGeom->IsDebugging()) printf("=> Inside ISVHWR\n");
efde9b4d 1452 if (check<0) return intHist;
2573ac89 1453 Int_t histInt = gGeoManager->GetCurrentNodeId()+1;
89125dc7 1454 if (gMCGeom->IsDebugging()) printf("<= ISVHWR: history is: %i in: %s\n", histInt, gGeoManager->GetPath());
efde9b4d 1455 return histInt;
1456}
1457
1458
8495a208 1459
1460