Logo ROOT  
Reference Guide
Loading...
Searching...
No Matches
g2root.f
Go to the documentation of this file.
1*CMZ : 2.22/00 08/04/99 18.41.30 by Rene Brun
2*CMZ : 2.21/05 08/02/99 12.19.15 by Rene Brun
3*CMZ : 2.00/11 22/08/98 15.25.26 by Rene Brun
4*CMZ : 2.00/10 28/07/98 13.32.42 by Rene Brun
5*-- Author :
6 Program g2root
7*
8**********************************************************************
9*
10* Program to convert an existing GEANT geometry/RZ file
11* into a ROOT macro (C++ file).
12*
13* Author: Rene Brun
14* modified by Nikolay I. Root <nroot@inp.nsk.su> to support map_names
15* modified by Mihaela Gheata for the new geometry format
16*
17* To use this conversion program (in $ROOTSYS/bin),
18* g2root [-f <map_names>] <geant_rzfile> <macro_name> [lrecl]
19* run g2root without parameters to see the usage info.
20*
21* for example
22* g2root brahms.geom brahms.C
23* will convert the GEANT RZ file brahms.geom into a ROOT macro brahms.C
24*
25* The default value for lrecl is 1024. The parameter lrecl must be specified
26* if the record length of the Zebra file is greater than 1024.
27*
28*
29* To generate the Geometry structure within Root, do:
30*
31* root[0]> .x brahms.C;
32* root[0]> new TBrowser();
33*
34* An interactive session
35* ------------------------
36*
37* Provided that a geometry was successfully built and closed , the manager class will
38* register itself to ROOT and the logical/physical structures will become immediately
39* browsable. The ROOT browser will display starting from the geometry folder : the list of
40* transformations and materials, the top volume and the top logical node. These last
41* two can be fully expanded, any intermediate volume/node in the browser being subject
42* of direct access context menu operations (right mouse button click). All user
43* utilities of classes TGeoManager, TGeoVolume and TGeoNode can be called via the
44* context menu.
45*
46* --- Drawing the geometry
47*
48* Any logical volume can be drawn via TGeoVolume::Draw() member function.
49* This can be direcly accessed from the context menu of the volume object
50* directly from the browser.
51* There are several drawing options that can be set with
52* TGeoManager::SetVisOption(Int_t opt) method :
53* opt=0 - only the content of the volume is drawn, N levels down (default N=3).
54* This is the default behavior. The number of levels to be drawn can be changed
55* via TGeoManager::SetVisLevel(Int_t level) method.
56*
57* opt=1 - the final leaves (e.g. daughters with no containment) of the branch
58* starting from volume are drawn down to the current number of levels.
59* WARNING : This mode is memory consuming
60* depending of the size of geometry, so drawing from top level within this mode
61* should be handled with care for expensive geometries. In future there will be
62* a limitation on the maximum number of nodes to be visualized.
63*
64* opt=2 - only the clicked volume is visualized. This is automatically set by
65* TGeoVolume::DrawOnly() method
66* opt=3 - only a given path is visualized. This is automatically set by
67* TGeoVolume::DrawPath(const char *path) method
68*
69* The current view can be exploded in cartesian, cylindrical or spherical
70* coordinates :
71* TGeoManager::SetExplodedView(Int_t opt). Options may be :
72* - 0 - default (no bombing)
73* - 1 - cartesian coordinates. The bomb factor on each axis can be set with
74* TGeoManager::SetBombX(Double_t bomb) and corresponding Y and Z.
75* - 2 - bomb in cylindrical coordinates. Only the bomb factors on Z and R
76* are considered
77*
78* - 3 - bomb in radial spherical coordinate : TGeoManager::SetBombR()
79*
80* Volumes themselves support different visualization settings :
81* - TGeoVolume::SetVisibility() : set volume visibility.
82* - TGeoVolume::VisibleDaughters() : set daughters visibility.
83* All these actions automatically updates the current view if any.
84*
85* --- Checking the geometry
86*
87* Several checking methods are accesible from the volume context menu. They
88* generally apply only to the visible parts of the drawn geometry in order to
89* ease geometry checking, and their implementation is in the TGeoChecker class
90* from the painting package.
91*
92* 1. Checking a given point.
93* Can be called from TGeoManager::CheckPoint(Double_t x, Double_t y, Double_t z).
94* This method is drawing the daughters of the volume containing the point one
95* level down, printing the path to the deepest physical node holding this point.
96* It also computes the closest distance to any boundary. The point will be drawn
97* in red.
98*
99* 2. Shooting random points.
100* Can be called from TGeoVolume::RandomPoints() (context menu function) and
101* it will draw this volume with current visualization settings. Random points
102* are generated in the bounding box of the top drawn volume. The points are
103* classified and drawn with the color of their deepest container. Only points
104* in visible nodes will be drawn.
105*
106*
107* 3. Raytracing.
108* Can be called from TGeoVolume::RandomRays() (context menu of volumes) and
109* will shoot rays from a given point in the local reference frame with random
110* directions. The intersections with displayed nodes will appear as segments
111* having the color of the touched node. Drawn geometry will be then made invisible
112* in order to enhance rays.
113*
114* IMPORTANT NOTE
115* To be compiled, this program requires a Fortran compiler supporting
116* recursive calls.
117*
118**********************************************************************
119* NOT YET MAPPED FROM OLD g2root
120**********************************************************************
121* The following lines starting with the 2 characters ** are an example
122* of (map_names> file.
123* To make a valid <map_names> file, remove these 2 characters.
124*
125**#
126**# formal definitions :
127**#
128**# 'comments' - the things that are ignored by the parser.
129**# The lines starting with '#' or ' ' - comments
130**# 'map' - first two words of 'non-comments'
131**# ==> the trailing part of this lines - also 'comments'
132**
133**# next lines are the 'map' examples :
134**
135**CMD2 Detector Names translation map for CMD2 detector
136**VALL Internals Internal part of detector.
137**VTBE TubeBe
138**VTAL TubeAL
139**VTML DCInner
140**DC0_ DCOuter
141**DRIF DCCell
142**
143** first part of map ('key') consists of exactly 4 characters,
144** including trailing blanks. This part of map is a key of
145** GEANT3 volumes/detectors.
146** The second part - 'name' - is just a sequence of non-blank chars.
147** 'name' is used as a replacement of 'key' in g2root output file.
148**
149** 'alias' - is a map with 'key' and 'name'
150** having one 'stab character' - '%'.
151**
152**DCC% DCSector% - example of alias.
153**
154** For above example - any 'keys' that have first 3 chars 'DCC'
155** would be matched by this alias.
156** The last char of 'key' - 'stab' - used as substitution for
157** '%' in 'name' (like gmake rules, but not so complicated).
158**
159** Keys are matched against aliases only if they do not
160** match any explicit 'map'. First found alias, that match
161** the key - is used for substitution.
162**
163** ==> The order of aliases is important for matching !
164**
165** The last alias may be like this :
166**
167**% Block% Match any key !
168**
169**% % <- assumed by default.
170*
171**********************************************************************
172
173 parameter(nwpaw=4000000)
174 common/pawc/paw(nwpaw)
175
176 character *80 gname
177 character *80 fname
178 character *8 chtop
179 character *8 crecl
180 integer npar, lrecl
181
182
183 call hlimit(nwpaw)
184
185 npar = iargc()
186 if (npar.eq.0) then
187 print *,
188 + 'Invoke g2root [-f map_name] geant_name macro_name [lrecl]'
189 go to 90
190 endif
191 narg = 1
192 call getarg(narg,gname)
193 narg = narg + 1
194 npar = npar - 1
195 if(gname .eq. "-f") then
196 if(npar .eq. 0) then
197 print *,'Invoke g2root [-f map_name] geant_name macro_name'
198 print *,' Parse error: you need specify <map_name>'
199 go to 90
200 endif
201 call getarg(narg,gname)
202 call create_map(gname)
203 narg = narg + 1
204 npar = npar - 1
205 if(npar .eq. 0) then
206 print *,'Invoke g2root [-f map_name] geant_name macro_name'
207 print *,' Parse error: you need specify <geant_name>'
208 go to 90
209 endif
210 call getarg(narg,gname)
211 narg = narg + 1
212 npar = npar - 1
213 endif
214 if (npar.ge.1) then
215 call getarg(narg,fname)
216 narg = narg + 1
217 else
218 idot=index(gname,'.')
219 fname = gname(1:idot-1)//'.C'
220 endif
221
222 lrecl = 1024
223 if (npar.ge.2) then
224 call getarg(narg,crecl)
225 read (crecl,'(I6)') lrecl
226 endif
227 call rzopen(1,chtop,gname,'W',lrecl,istat)
228 if (istat.ne.0) then
229 print *,'Cannot open file'
230 go to 90
231 endif
232 call rzfile(1,chtop,' ')
233* call rzldir(' ',' ')
234
235 call g2rin
236
237 call toroot(fname)
238 lg = lenocc(gname)
239 lf = lenocc(fname)
240 print 1000,gname(1:lg),fname(1:lf)
241 1000 format(' GEANT file: ',a, ' converted to ROOT macro: ',a)
242 90 continue
243 end
244 Subroutine toroot(fname)
245*
246**********************************************************************
247*
248* Rotation matrices (structure JROTM) are saved into a linked
249* list of class objects TRotMatrix.
250* The JVOLUM structure generates the following lists:
251* - the TMaterial list (material definition only).
252* - the TRotmatrix list (Rotation matrices definition only).
253* - the TShape list (volume definition only).
254* - the TNode list assembling all detector elements
255* the TNode list is a real tree.
256* The Node list contains two variants TNode and TNodeDiv
257* corresponding to the GSPOS and GSDIV mechanisms.
258*
259* Author: Rene Brun
260**********************************************************************
261*
262*KEEP,HCBOOK.
263 parameter(nwpaw=4000000)
264 common/pawc/paw(nwpaw)
265
266 INTEGER IQ(2), LQ(8000)
267 REAL Q(2)
268 equivalence(lq(1),paw(11)),(iq(1),paw(19)),(q(1),paw(19))
269
270 INTEGER JDIGI ,JDRAW ,JHEAD ,JHITS ,JKINE ,JMATE ,JPART
271 + ,JROTM ,JRUNG ,JSET ,JSTAK ,JGSTAT,JTMED ,JTRACK,JVERTX
272 + ,JVOLUM,JXYZ ,JGPAR ,JGPAR2,JSKLT
273C
274 common/gclink/jdigi ,jdraw ,jhead ,jhits ,jkine ,jmate ,jpart
275 + ,jrotm ,jrung ,jset ,jstak ,jgstat,jtmed ,jtrack,jvertx
276 + ,jvolum,jxyz ,jgpar ,jgpar2,jsklt
277C
278*KEEP,GCNUM.
279 common/gcnum/nmate ,nvolum,nrotm,ntmed,ntmult,ntrack,npart
280 + ,nstmax,nvertx,nhead,nbit
281*
282
283 common/cnodes/nnodes
284
285 parameter(maxpos=250000)
286* parameter (MAXPOS=50000)
287 common/cnpos/nodepos(maxpos),nodediv(maxpos),nvflags(maxpos),
288 +npflags(maxpos),nppflags(maxpos)
289
290 CHARACTER*4 KSHAP(30),klshap(30)
291 character*20 matname,medname
292 character*16 cname,mname,pname, rname
293 character*(*) fname
294 character*256 line
295 character*128 creals
296 character*16 astring,astring2,ptrname
297 dimension pmixt(3)
298 logical map_found
299
300 DATA kshap/'BRIK','TRD1','TRD2','TRAP','TUBE','TUBS','CONE',
301 + 'CONS','SPHE','PARA','PGON','PCON','ELTU','HYPE',
302 + 13*' ','GTRA','CTUB',' '/
303*________________________________________________________________________
304*
305
306 open(unit=51,file=fname,status='unknown')
307*
308* Create new Geometry object
309* ==========================
310 nch=lenocc(fname)
311 idot=index(fname,'.')
312 nct=idot-1
313 write(51,1111)fname(1:nct)
3141111 format('void ',a,'()',/,'{',/,
315 +'//',/,
316 +'// This file has been generated automatically via the root',/,
317 +'// utility g2root from an interactive version of GEANT',/,
318 +'// (see ROOT class TGeoManager for an example of use)',/,
319 +'//',/,
320 +'R__ASSERT(gSystem->Load("libGeom") >= 0);',/,
321 +'TGeoRotation *rot;',/,
322 +'TGeoNode *Node, *Node1;')
323
324 do 1 i=1,30
325 klshap(i) = kshap(i)
326 call cutol(klshap(i))
327 1 continue
328 do 2 i=1,maxpos
329 nodepos(i) = 0
330 nodediv(i) = 0
331 2 continue
332 nodepos(1) = 1
333 nodediv(1) = 1
334
335 write(51,490)fname(1:nct),fname(1:nct),fname(1:nch)
336 490 format(/,'TGeoManager *',a,' = new TGeoManager("',a,'","',a,'");'
337 +,/)
338 IF(jvolum.NE.0 ) nvolum = iq(jvolum-2)
339 IF(jmate.NE.0 ) nmate = iq(jmate-2)
340 IF(jtmed.NE.0 ) ntmed = iq(jtmed-2)
341 IF(jrotm.NE.0 ) nrotm = iq(jrotm-2)
342* Print Materials
343* =======================
344 write(51,3019)
345 3019 format(/,
346 +'//-----------List of Materials and Mixtures--------------',/)
347 do 300 imat=1,nmate
348 jma=lq(jmate-imat)
349 if(jma.eq.0)go to 300
350 nmixt=q(jma+11)
351 call uhtoc(iq(jma+1),4,matname,20)
352 ncn=lenocc(matname)
353 call toint(imat,astring,nc)
354 nm=abs(nmixt)
355*-* Case of a simple material
356 if (nm.le.1)then
357 call toreals(3,q(jma+6),creals,ncr)
358 if(q(jma+6).lt.1.and.q(jma+7).lt.1)then
359 creals=',0,0,0'
360 ncr=lenocc(creals)
361 endif
362 line=' '
363 write(line,3000)astring(1:nc),matname(1:ncn),creals(1:ncr)
364 3000 format('TGeoMaterial *mat',a,' = new TGeoMaterial("',a,
365 + '"',a,');')
366 nch = lenocc(line)
367 write(51,'(a)')line(1:nch)
368 write(line,3005) astring(1:nc),imat
369 3005 format(4x,'mat',a,'->SetUniqueID(',i4,');')
370 nch = lenocc(line)
371 write(51,'(a)')line(1:nch)
372*-* Case of a mixture
373 else
374 call toint(nm,creals,ncm)
375 jmixt=lq(jma-5)
376 if(nmixt.gt.0)then
377 mname=creals(1:ncm)
378 else
379 mname='-'//creals(1:ncm)
380 ncm=ncm+1
381 endif
382 line=' '
383 write(line,3010)astring(1:nc),matname(1:ncn),mname(1:ncm),
384 + q(jma+8)
385 3010 format('TGeoMixture *mat',a,' = new TGeoMixture("',a,'",',a,
386 + ',',g14.6,');')
387 nch = lenocc(line)
388 write(51,'(a)')line(1:nch)
389 write(line,3011) astring(1:nc),imat
390 3011 format(4x,'mat',a,'->SetUniqueID(',i4,');')
391 nch = lenocc(line)
392 write(51,'(a)')line(1:nch)
393 do 292 im=1,nm
394 call toint(im-1,astring2,nc2)
395 pmixt(1) = q(jmixt+im)
396 pmixt(2) = q(jmixt+nm+im)
397 pmixt(3) = q(jmixt+2*nm+im)
398 call toreals(3,pmixt,creals,ncr)
399 line=' '
400 write(line,3020)astring(1:nc),astring2(1:nc2),
401 + creals(1:ncr)
402 3020 format(4x,'mat',a,'->DefineElement(',a,a,');')
403 nch = lenocc(line)
404 write(51,'(a)')line(1:nch)
405 292 continue
406 endif
407 300 continue
408* Print Tracking Media
409* ====================
410 write(51,3069)
411 3069 format(/,
412 +'//-----------List of Tracking Media--------------',/)
413 do 350 itmed=1,ntmed
414 jtm=lq(jtmed-itmed)
415 if(jtm.eq.0)go to 350
416 imat=q(jtm+6)
417 call toint(imat,astring2,ncm)
418 call uhtoc(iq(jtm+1),4,medname,20)
419 ncn=lenocc(medname)
420 call toint(itmed,astring,nc)
421 call toreals(8,q(jtm+7),creals,ncr)
422 line=' '
423 write(line,3050)astring(1:nc),medname(1:ncn),astring(1:nc),
424 + astring2(1:ncm),creals(1:ncr)
425 3050 format('TGeoMedium *med',a,' = new TGeoMedium("',a,'",',a,
426 + ',',a,a,');')
427 nch = lenocc(line)
428 write(51,'(a)')line(1:nch)
429 350 continue
430* Print Rotation matrices
431* =======================
432 write(51,3021)
433 3021 format(/,'//-----------List of Rotation matrices--------------',/)
434 do 100 irot=1,nrotm
435 jr=lq(jrotm-irot)
436 if(jr.eq.0)go to 100
437 call toint(irot,astring,nc)
438 call toreals(6,q(jr+11),creals,ncr)
439 line=' '
440 ptrname = 'rot'//astring(1:nc)
441 nch = nc+3
442 write(line,1000)ptrname(1:nch),ptrname(1:nch),
443 + creals(1:ncr)
444 1000 format('TGeoRotation *',a,
445 + ' = new TGeoRotation("',a,'"',a,');')
446 nch = lenocc(line)
447 write(51,'(a)')line(1:nch)
448 100 continue
449
450* Print volume definition (ROOT shapes)
451* =======================
452 write(51,3022)
453 3022 format(/,'//-----------List of Volumes--------------',/)
454 print *,' nvolum= ',nvolum, ' jvolum=',jvolum
455C
456C ???? Convert GEANT3 keys to legal C++ names ?????
457C
458 do 50 ivo = 1,nvolum
459 if (lq(jvolum-ivo).eq.0)go to 50 ! <-- That's
460 cname=' ' ! 'real'
461* write(cname,'(a4)')iq(jvolum+ivo) ! trick !!!
462 call uhtoc(iq(jvolum+ivo),4,cname,4) ! <--
463 do 29 i=1,4
464 if (ichar(cname(i:i)).eq.0)cname(i:i) = ' '
465 29 continue
466 n1=lenocc(cname)
467 if(n1.lt.4)then
468 do 30 i=n1+1,4
469 cname(i:i)='_'
470 30 continue
471 endif
472 50 continue
473C----------------------------------------------
474 do 77 ivo = 1,nvolum
475 77 nvflags(ivo) = 0
476 nlevel = 0
477 call markdiv(1,1)
478
479 do 200 ivo = 1,nvolum
480 if (nvflags(ivo).eq.2) goto 200
481 jv=lq(jvolum-ivo)
482 if (jv.eq.0)go to 200
483 cname=' '
484 if(.not.map_found(iq(jvolum+ivo),cname)) then
485 write(cname,'(a4)')iq(jvolum+ivo)
486 endif
487 call volume(cname,q(jv+1),0,0)
488 200 continue
489
490 do 88 ivo = 1,nvolum
491 nvflags(ivo) = 0
492 npflags(ivo) = 0
493 nppflags(ivo) = 0
494 88 continue
495* Print volume positioning (ROOT nodes)
496* ========================
497
498 nnodes = 1
499 nlevel = 0
500 if (nvolum.gt.0) then
501 call node(1,1,1)
502 write(51,3023)
503 3023 format(/,'//-----------List of Nodes--------------',/)
504 do 89 ivo = 1,nvolum
505 89 nvflags(ivo) = 0
506 call node(1,1,0)
507 endif
508
509 write(51,2223)
510 write(51,2222)
5112222 format('}')
5122223 format(' gGeoManager->CloseGeometry();')
513 close(51)
514 end
515
516*_______________________________________________________________________
517 Subroutine markdiv2(ivo,nuserm)
518 call markdiv(ivo,nuserm)
519 End
520
521*_______________________________________________________________________
522 Subroutine markdiv(ivo,nuserm)
523*
524* Process one node (volume with contents)
525*KEEP,HCBOOK.
526 parameter(nwpaw=4000000)
527 parameter(maxpos=250000)
528 common/pawc/paw(nwpaw)
529
530 INTEGER IQ(2), LQ(8000)
531 REAL Q(2)
532 equivalence(lq(1),paw(11)),(iq(1),paw(19)),(q(1),paw(19))
533
534 INTEGER JDIGI ,JDRAW ,JHEAD ,JHITS ,JKINE ,JMATE ,JPART
535 + ,JROTM ,JRUNG ,JSET ,JSTAK ,JGSTAT,JTMED ,JTRACK,JVERTX
536 + ,JVOLUM,JXYZ ,JGPAR ,JGPAR2,JSKLT
537C
538 common/gclink/jdigi ,jdraw ,jhead ,jhits ,jkine ,jmate ,jpart
539 + ,jrotm ,jrung ,jset ,jstak ,jgstat,jtmed ,jtrack,jvertx
540 + ,jvolum,jxyz ,jgpar ,jgpar2,jsklt
541C
542 common/cnodes/nnodes
543 common/clevel/nodeold(20),nlevel
544
545 common/cnpos/nodepos(maxpos),nodediv(maxpos),nvflags(maxpos),
546 +npflags(maxpos),nppflags(maxpos)
547
548* dimension qjv(1000)
549 character*16 cname
550*---------------------------------------------------------------------
551 nlevel = nlevel + 1
552 nodeold(nlevel) = nnodes
553 jv=lq(jvolum-ivo)
554 ishape = q(jv+2)
555 nin = q(jv+3)
556*- Loop subnodes
557 if(nin.eq.0)then
558 nlevel=nlevel-1
559 return
560 endif
561 call cdnode(nodeold(nlevel))
562 if(nin.gt.0)then
563 if (nvflags(ivo).ne.0) then
564 goto 996
565 endif
566 nvflags(ivo)=1
567* Volume has positioned contents
568 do 300 in=1,nin
569 jin=lq(jv-in)
570 ivom=q(jin+2)
571 nuser = q(jin+3)
572 jinvom = lq(jvolum-ivom)
573 ninvom = q(jinvom+3)
574 cname=' '
575 write(cname,'(a4)')iq(jvolum+ivom)
576 n1=lenocc(cname)
577 if(ninvom.ge.0)then
578 nnodes = nnodes+1
579 if (nnodes.gt.maxpos) then
580 print *,'Too many nodes =',nnodes
581 go to 300
582 endif
583 if (nodepos(nnodes).eq.0) then
584 nodepos(nnodes) = 1
585 endif
586 endif
587 if(ninvom.ne.0) then
588 call markdiv2(ivom,nuser)
589 endif
590 300 continue
591 else
592 nnodes = nnodes+1
593 if (nodediv(nnodes).eq.0) then
594 nodediv(nnodes) = 1
595 endif
596 jin=lq(jv-1)
597 ivod=q(jin+2)
598 if (nvflags(ivod).gt.0) goto 996
599 cname=' '
600 write(cname,'(a4)')iq(jvolum+ivod)
601 n1=lenocc(cname)
602 print 200, cname(1:n1)
603 200 format('Division volume', a4)
604 call markdiv2(ivod,0)
605 nvflags(ivod) = 2
606 endif
607
608996 continue
609 nlevel = nlevel - 1
610 if (nlevel.gt.0)call cdnode(nodeold(nlevel))
611 end
612*_______________________________________________________________________
613 subroutine volume(cname,qjv,iposp,ifirst)
614*KEEP,HCBOOK.
615 parameter(nwpaw=4000000)
616 parameter(maxpos=250000)
617 common/pawc/paw(nwpaw)
618
619 INTEGER IQ(2), LQ(8000)
620 REAL Q(2)
621 equivalence(lq(1),paw(11)),(iq(1),paw(19)),(q(1),paw(19))
622
623 INTEGER JDIGI ,JDRAW ,JHEAD ,JHITS ,JKINE ,JMATE ,JPART
624 + ,JROTM ,JRUNG ,JSET ,JSTAK ,JGSTAT,JTMED ,JTRACK,JVERTX
625 + ,JVOLUM,JXYZ ,JGPAR ,JGPAR2,JSKLT
626C
627 common/gclink/jdigi ,jdraw ,jhead ,jhits ,jkine ,jmate ,jpart
628 + ,jrotm ,jrung ,jset ,jstak ,jgstat,jtmed ,jtrack,jvertx
629 + ,jvolum,jxyz ,jgpar ,jgpar2,jsklt
630C
631 common/cnpos/nodepos(maxpos),nodediv(maxpos),nvflags(maxpos),
632 +npflags(maxpos),nppflags(maxpos)
633 character *(*) cname
634 character*16 astring,cmater, pname, rname
635 character*128 creals
636 character*256 line
637 real qjv(100)
638 double precision RADDEG
639 dimension dummypars(100)
640 dimension npars(30)
641 CHARACTER*6 KSHAP(30)
642 data dummypars/100*0./
643
644 DATA kshap/'Box','Trd1','Trd2','Trap','Tube','Tubs','Cone',
645 + 'Cons','Sphere','Para','Pgon','Pcon','Eltu','Hype',
646 + 13*' ','Gtra','Ctub',' '/
647 DATA npars/3,4,5,11,3,5,5,7,6,6,10,9,3,4,13*0,12,11,0/
648*________________________________________________________________________
649*
650
651 raddeg = 57.2957795130823209
652 n1=lenocc(cname)
653
654** print *, 'VOLUME n1=',n1,' cname=',cname(1:n1)
655 ishape = qjv(2)
656 numed = qjv(4)
657 jtm = lq(jtmed-numed)
658 call toint(numed,astring,nc)
659 cmater='med'//astring(1:nc)
660 ncmat = lenocc(cmater)
661 nord = qjv(1)
662 nin = qjv(3)
663 npar = qjv(5)
664 npar0 = npars(ishape)
665**TRAP
666 jpar = 6
667 if (ishape.eq.4) then
668 ph = 0.
669 if (qjv(jpar+2).ne.0.)ph=atan2(qjv(jpar+3),qjv(jpar+2))*raddeg
670 tt = sqrt(qjv(jpar+2)**2+qjv(jpar+3)**2)
671 qjv(jpar+2) = atan(tt)*raddeg
672 if (ph.lt.0.0) ph = ph+360.0
673 qjv(jpar+3) = ph
674 qjv(jpar+7) = atan(qjv(jpar+7))*raddeg
675 if (qjv(jpar+7).gt.90.0) qjv(jpar+7) = qjv(jpar+7)-180.0
676 qjv(jpar+11)= atan(qjv(jpar+11))*raddeg
677 if (qjv(jpar+11).gt.90.0) qjv(jpar+11) = qjv(jpar+11)-180.0
678 endif
679**PARA
680 if (ishape.eq.10) then
681 ph = 0.
682 if (qjv(jpar+5).ne.0.)ph=atan2(qjv(jpar+6),qjv(jpar+5))*raddeg
683 tt = sqrt(qjv(jpar+5)**2+qjv(jpar+6)**2)
684 qjv(jpar+4) = atan(qjv(jpar+4))*raddeg
685 if (qjv(jpar+4).gt.90.0) qjv(jpar+4) = qjv(jpar+4)-180.0
686 qjv(jpar+5) = atan(tt)*raddeg
687 if (ph.lt.0.0) ph = ph+360.0
688 qjv(jpar+6) = ph
689 endif
690 if(ishape.eq.11)npar0=4
691 if(ishape.eq.12)npar0=3
692**HYPE
693 if (ishape.eq.14) then
694 hyrmin = qjv(jpar+1)
695 hyrmax = qjv(jpar+2)
696 hydz = qjv(jpar+3)
697 hyst = qjv(jpar+4)
698 dummypars(1) = hyrmin
699 dummypars(2) = hyst
700 dummypars(3) = hyrmax
701 dummypars(4) = hyst
702 dummypars(5) = hydz
703 npar0 = -5
704 endif
705
706* print 2351, cname(1:n1),kshap(ishape)
707* 2351 format('Volume:',a, ' shape=',a)
708 if (npar.le.0) then
709** print 2352, cname(1:n1),kshap(ishape)
710** 2352 format('Warning, volume with 0 parameters:',a, ' shape=',a)
711 return
712 endif
713 if(npar0.le.0)then
714 call toreals(-npar0,dummypars(1),creals,ncr)
715 else
716 call toreals(npar0,qjv(7),creals,ncr)
717 endif
718 line=' '
719 nshape = lenocc(kshap(ishape))
720 call ptname(cname, pname)
721 np = lenocc(pname)
722 call realname(cname, rname)
723 nrr = lenocc(rname)
724 if (iposp.eq.0) then
725 write(line,2000)pname(1:np),kshap(ishape)(1:nshape)
726 + ,rname(1:nrr),cmater(1:ncmat),creals(1:ncr)
727 else
728 if (ifirst.eq.1) then
729 write(line,2001)pname(1:np),rname(1:nrr),cmater(1:ncmat)
730 nch=lenocc(line)
731 write(51,'(a)')line(1:nch)
732 endif
733 line=' '
734 write(line,2002)pname(1:np),kshap(ishape)(1:nshape)
735 + ,rname(1:nrr),cmater(1:ncmat),creals(1:ncr)
736 nch=lenocc(line)
737 write(51,'(a)')line(1:nch)
738 endif
7392000 format(
740 + 'TGeoVolume',' *',a,' = gGeoManager->Make',a,'("',a,'",'
741 +,a,a,');')
7422001 format('TGeoVolumeMulti *',a,' = gGeoManager->MakeVolumeMulti("'
743 +,a,'", ',a,');')
7442002 format(' ',a,'->AddVolume(gGeoManager->Make',a,'("',a,'",',
745 +a,a,'));')
746 nch = lenocc(line)
747 if (iposp.eq.0) write(51,'(a)')line(1:nch)
748 if(ishape.eq.11)then
749 ndz=qjv(10)
750 do iz=1,ndz
751 call toreals(3,qjv(11+(iz-1)*3),creals,ncr)
752 line=' '
753 call toint(iz-1,astring,nci)
754 if (iposp.eq.0) then
755 write(line,2010)pname(1:np),astring(1:nci),creals(1:ncr)
756 else
757 write(line,2011)pname(1:np),astring(1:nci),creals(1:ncr)
758 endif
759 nch = lenocc(line)
760 write(51,'(a)')line(1:nch)
761 enddo
762 endif
763 if(ishape.eq.12)then
764 ndz=qjv(9)
765 do iz=1,ndz
766 call toreals(3,qjv(10+(iz-1)*3),creals,ncr)
767 line=' '
768 call toint(iz-1,astring,nci)
769 if (iposp.eq.0) then
770 write(line,2010)pname(1:np),astring(1:nci),creals(1:ncr)
771 else
772 write(line,2011)pname(1:np),astring(1:nci),creals(1:ncr)
773 endif
774 nch = lenocc(line)
775 write(51,'(a)')line(1:nch)
776 enddo
777 endif
7782010 format(2x,'((TGeoPcon*)',a,'->GetShape())->DefineSection(',
779 + a,a,');')
7802011 format(2x,'((TGeoPcon*)',a,'->GetLastShape())->DefineSection(',
781 + a,a,');')
782* Any attributes set ?
783 lseen = qjv(npar+8)
784 lstyle = qjv(npar+9)
785 lwidth = qjv(npar+10)
786 lcolor = qjv(npar+11)
787 lfill = qjv(npar+12)
788 if (lstyle.le.0) lstyle = 1
789 if (lwidth.le.0) lwidth = 1
790 if (lcolor.lt.0) lcolor = 1
791 if (lfill.lt.0) lfill = 0
792* if(ivo.eq.1)lseen=0
793* if(nord.lt.0)then
794* print *,'ordering : ',-nord
795* call toint(-nord,creals,ncr)
796* endif
797 if ((iposp.eq.0).or.(ifirst.eq.1)) then
798 if(lseen.ne.1)then
799 call toint(lseen,creals,ncr)
800 write(51,195)pname(1:np),creals(1:ncr)
801195 format(2x,a,'->SetVisibility(',a,');')
802 endif
803 if(lstyle.ne.1)then
804 call toint(lstyle,creals,ncr)
805 write(51,196)pname(1:np),creals(1:ncr)
806196 format(2x,a,'->SetLineStyle(',a,');')
807 endif
808 if(lwidth.ne.1)then
809 call toint(lwidth,creals,ncr)
810 write(51,197)pname(1:np),creals(1:ncr)
811197 format(2x,a,'->SetLineWidth(',a,');')
812 endif
813 if(lcolor.ne.1)then
814 call toint(lcolor,creals,ncr)
815 write(51,198)pname(1:np),creals(1:ncr)
816198 format(2x,a,'->SetLineColor(',a,');')
817 endif
818 if(lfill.ne.0)then
819 call toint(lfill,creals,ncr)
820 write(51,199)pname(1:np),creals(1:ncr)
821199 format(2x,a,'->SetFillStyle(',a,');')
822 endif
823 endif
824 end
825
826*_______________________________________________________________________
827 Subroutine node2(ivo,nuserm,iposp)
828 call node(ivo,nuserm,iposp)
829 End
830
831*_______________________________________________________________________
832 Subroutine node(ivo,nuserm,iposp)
833*
834* Process one node (volume with contents)
835*KEEP,HCBOOK.
836 parameter(nwpaw=4000000)
837 parameter(maxpos=250000)
838 common/pawc/paw(nwpaw)
839
840 INTEGER IQ(2), LQ(8000)
841 REAL Q(2)
842 equivalence(lq(1),paw(11)),(iq(1),paw(19)),(q(1),paw(19))
843
844 INTEGER JDIGI ,JDRAW ,JHEAD ,JHITS ,JKINE ,JMATE ,JPART
845 + ,JROTM ,JRUNG ,JSET ,JSTAK ,JGSTAT,JTMED ,JTRACK,JVERTX
846 + ,JVOLUM,JXYZ ,JGPAR ,JGPAR2,JSKLT
847C
848 common/gclink/jdigi ,jdraw ,jhead ,jhits ,jkine ,jmate ,jpart
849 + ,jrotm ,jrung ,jset ,jstak ,jgstat,jtmed ,jtrack,jvertx
850 + ,jvolum,jxyz ,jgpar ,jgpar2,jsklt
851C
852 common/cnodes/nnodes
853 common/clevel/nodeold(20),nlevel
854
855 common/cnpos/nodepos(maxpos),nodediv(maxpos),nvflags(maxpos),
856 +npflags(maxpos),nppflags(maxpos)
857
858 dimension qjv(1000)
859 character*16 cnode,cname,mname,anode,mother,pname, rname
860 integer nmother
861 character*256 line
862 character*128 creals
863 character*16 astring,astring1
864 character*16 matrix
865 character*256 matrixs
866 character *16 cblank
867 Logical map_found
868 data cblank/' '/
869*---------------------------------------------------------------------
870 cblank = ' '
871 nlevel = nlevel + 1
872 nodeold(nlevel) = nnodes
873 jv=lq(jvolum-ivo)
874 ishape = q(jv+2)
875 nin = q(jv+3)
876 mname=' '
877 if(.not.map_found(iq(jvolum+ivo),mname)) then
878 write(mname,'(a4)')iq(jvolum+ivo)
879 endif
880 n2=lenocc(mname)
881 call toint(nuserm,astring,nci)
882*- If top volume, create the top node
883 call ptname(mname, mother)
884 nmother = lenocc(mother)
885 if((ivo.eq.1).and.(iposp.eq.0))then
886 write(51,510)mother(1:nmother)
887* write(51,510)
888510 format('gGeoManager->SetTopVolume(',a,');')
889 endif
890*- Generate subnodes of this node (if any)
891** print 2346, iq(jvolum+ivo), nin
8922346 format('Processing node:',a4,' nin=',i4)
893 if(nin.eq.0)then
894 nlevel=nlevel-1
895 return
896 endif
897 call cdnode(nodeold(nlevel))
898* print 520, mother(1:nmother), ivo
899520 format('mother ',a,' index ',i9)
900 if(nin.gt.0)then
901 if (nvflags(ivo).ne.0) then
902 goto 996
903 endif
904 nvflags(ivo)=1
905* Volume has positioned contents
906 do 300 in=1,nin
907 ifirst = 0
908 icurrent = 0
909 imulti = 0
910 nci1 = 0
911 jin=lq(jv-in)
912 ivom=q(jin+2)
913 nuser = q(jin+3)
914 imany = q(jin+8)
915* print *,'in=',in,' nuser=',nuser
916 jinvom = lq(jvolum-ivom)
917 npar = q(jinvom+5)
918 ninvom = q(jinvom+3)
919 cname=' '
920 if(.not.map_found(iq(jvolum+ivom),cname)) then
921 write(cname,'(a4)')iq(jvolum+ivom)
922 endif
923 n1=lenocc(cname)
924 if (npar.gt.0) then
925 jpar = jinvom+6
926 else
927 jpar = jin+9
928 if (iposp.eq.1) then
929 if (npflags(ivom).eq.0) then
930 ifirst = 1
931 npflags(ivom) = 1
932 else
933 npflags(ivom) = npflags(ivom)+1
934 endif
935 else
936 icurrent = nppflags(ivom)
937 call toint(icurrent,astring1,nci1)
938 imulti = 1
939 nppflags(ivom) = nppflags(ivom)+1
940 endif
941 npar = q(jin+9)
942 call ucopy(q(jinvom+1),qjv(1),6)
943 qjv(5) = npar
944 call ucopy(q(jin+10),qjv(7),npar)
945 call ucopy(q(jinvom+7),qjv(7+npar),6)
946 call toint(in,astring,nci)
947 mname=cname(1:n1)//astring(1:nci)
948 if (iposp.eq.1) call volume(cname,qjv,iposp,ifirst)
949 endif
950 if(ninvom.ge.0)then
951 nnodes = nnodes+1
952 if (nnodes.gt.maxpos) then
953 print *,'Too many nodes =',nnodes
954 go to 300
955 endif
956 call toint(nnodes,anode,ncd)
957 cnode = 'Node'//anode(1:ncd)
958 if (nodepos(nnodes).eq.0) then
959* write(51,4444)cblank(1:nlevel),anode(1:ncd)
960 4444 format(a,'TNode *Node',a,';')
961 nodepos(nnodes) = 1
962 endif
963 else
964 cnode = 'Node'
965 endif
966 nd=lenocc(cnode)
967 call toreals(3,q(jin+5),creals,ncr)
968 itrans = 1
969 if ((abs(q(jin+5)).lt.1e-30).and.
970 + (abs(q(jin+6)).lt.1e-30).and.
971 + (abs(q(jin+7)).lt.1e-30)) then
972 itrans = 0
973 endif
974 irot=q(jin+4)
975 matrixs=' '
976 if(irot.eq.0)then
977 matrix='0'
978 ncmatrix=1
979 if (itrans.eq.0) then
980 matrixs='gGeoIdentity'
981 else
982 matrixs='new TGeoTranslation('//creals(2:ncr)//')'
983 endif
984 else
985 call toint(irot,astring,nci)
986 matrix='rot'//astring(1:nci)
987 ncmatrix=nci+3
988 if (itrans.eq.0) then
989 matrixs=matrix(1:ncmatrix)
990 else
991 matrixs='new TGeoCombiTrans('//creals(2:ncr)//','//
992 + matrix(1:ncmatrix)//')'
993 endif
994 endif
995 call toint(nuser,astring,nci)
996** print *,' cname=',cname(1:n1), ' astring=',astring(1:nci)
997 mname=cname(1:n1)//astring(1:nci)
998 n2=lenocc(mname)
999 ncmats=lenocc(matrixs)
1000 line=' '
1001
1002 call ptname(cname, pname)
1003 np = lenocc(pname)
1004 if (imany.eq.1) then
1005 if (imulti.eq.0) then
1006 write(line,3000)cblank(1:nlevel),mother(1:nmother),
1007 + pname(1:np), astring(1:nci), matrixs(1:ncmats)
1008 3000 format(a,a,'->AddNode(',a,',',a,',',a,');')
1009 else
1010 write(line,3002)cblank(1:nlevel),mother(1:nmother),
1011 + pname(1:np), astring1(1:nci1),astring(1:nci),
1012 + matrixs(1:ncmats)
1013 3002 format(a,a,'->AddNode(',a,'->GetVolume(',a,'),',a
1014 + ,',',a,');')
1015 endif
1016 else
1017 if (imulti.eq.0) then
1018 write(line,3001)cblank(1:nlevel),mother(1:nmother),
1019 + pname(1:np), astring(1:nci), matrixs(1:ncmats)
1020 3001 format(a,a,'->AddNodeOverlap(',a,',',a,',',a,');')
1021 else
1022 write(line,3003)cblank(1:nlevel),mother(1:nmother),
1023 + pname(1:np), astring1(1:nci1), astring(1:nci),
1024 + matrixs(1:ncmats)
1025 3003 format(a,a,'->AddNodeOverlap(',a,'->GetVolume(',a,'),',a
1026 + ,',',a,');')
1027 endif
1028 endif
1029 nch = lenocc(line)
1030 if (iposp.eq.0) write(51,'(a)')line(1:nch)
1031 npar=q(jv+5)
1032 if(ninvom.ne.0) then
1033 call node2(ivom,nuser,iposp)
1034 endif
1035 300 continue
1036 else
1037* Print *,'===== DIVISION ====='
1038* Print 4567,mother(1:nmother)
1039 nnodes = nnodes+1
1040 call toint(nnodes,anode,ncd)
1041 cnode = 'Nodiv'//anode(1:ncd)
1042 nd=lenocc(cnode)
1043 if (nodediv(nnodes).eq.0) then
1044 nodediv(nnodes) = 1
1045 endif
1046 jin=lq(jv-1)
1047 ivod=q(jin+2)
1048 if (nvflags(ivod).eq.1) goto 996
1049 cname=' '
1050* if(.not.map_found(iq(jvolum+ivod),cname)) then
1051 write(cname,'(a4)')iq(jvolum+ivod)
1052* endif
1053* Print 4445,iq(jvolum+ivod)
1054 4445 format('daughter division', a4)
1055 n1=lenocc(cname)
1056* if(cname(n1:n1).eq.'+')cname(n1:)='plus'
1057* if(cname(n1:n1).eq.'-')cname(n1:)='minus'
1058* n1=lenocc(cname)
1059* cname=cname(1:n1)//'_0'
1060* n2=lenocc(cname)
1061 iaxis=q(jin+1)
1062 call toint(iaxis,astring,nci)
1063 call toreals(3,q(jin+3),creals,ncr)
1064 line = ' '
1065 call ptname(cname, pname)
1066 np = lenocc(pname)
1067 call realname(cname, rname)
1068 nrr = lenocc(rname)
1069 write(line,995)cblank(1:nlevel),pname(1:np),mother(1:nmother),
1070 + rname(1:nrr),astring(1:nci), creals(1:ncr)
1071 995 format(a,'TGeoVolume *',a,' = ',a,'->Divide("',a,'",',a,a,');')
1072 nch = lenocc(line)
1073 if (iposp.eq.0) write(51,'(a)')line(1:nch)
1074
1075 call node2(ivod,0,iposp)
1076 nvflags(ivod) = 1
1077 endif
1078
1079996 continue
1080 nlevel = nlevel - 1
1081 if (nlevel.gt.0)call cdnode(nodeold(nlevel))
1082 end
1083
1084 subroutine realname(cname, pname)
1085 character *4 cname
1086 character *16 pname
1087 nind = 0
1088 pname = ' '
1089 write(pname,'(a4)')cname
1090 do i=1,4
1091 nind = nind+1
1092 pname(nind:nind)=cname(i:i)
1093 if(ichar(cname(i:i)).eq.0) then
1094 pname(nind:nind)=' '
1095 nind = nind+1
1096 pname(nind:nind)=' '
1097 endif
1098 if(ichar(cname(i:i)).eq.92) then
1099 pname(nind:nind)=char(92)
1100 nind = nind+1
1101 pname(nind:nind)=char(92)
1102 endif
1103 if(ichar(cname(i:i)).eq.34) then
1104 pname(nind:nind)=char(92)
1105 nind = nind+1
1106 pname(nind:nind)=char(34)
1107 endif
1108 enddo
1109*------ supress blanks
11102333 if (pname(lenocc(pname):lenocc(pname)).eq.' ') then
1111 pname = pname(1:lenocc(pname)-1)
1112 goto 2333
1113 endif
1114 end
1115
1116 subroutine ptname(cname, pname)
1117 character *4 cname
1118 character *16 pname
1119 pname = ' '
1120 write(pname,'(a4)')cname
1121 do i=1,4
1122 if(ichar(cname(i:i)).eq.0) then
1123 pname(i:i)='_'
1124 pname(5:5)='_'
1125 endif
1126 if(ichar(cname(i:i)).eq.92) then
1127 pname(i:i)='a'
1128 pname(5:5)='a'
1129 endif
1130 if(cname(i:i).eq.'?') then
1131 pname(i:i)='b'
1132 pname(5:5)='b'
1133 endif
1134 if(cname(i:i).eq.'`') then
1135 pname(i:i)='c'
1136 pname(5:5)='c'
1137 endif
1138 if(cname(i:i).eq.' ') then
1139 pname(i:i)='_'
1140 endif
1141 if(cname(i:i).eq.'+') then
1142 pname(i:i)='d'
1143 pname(5:5)='d'
1144 endif
1145 if(cname(i:i).eq.'-') then
1146 pname(i:i)='e'
1147 pname(5:5)='e'
1148 endif
1149 if(cname(i:i).eq.'*') then
1150 pname(i:i)='f'
1151 pname(5:5)='f'
1152 endif
1153 if(cname(i:i).eq.'/') then
1154 pname(i:i)='g'
1155 pname(5:5)='g'
1156 endif
1157 if(cname(i:i).eq.'.') then
1158 pname(i:i)='h'
1159 pname(5:5)='h'
1160 endif
1161 if(cname(i:i).eq.'''') then
1162 pname(i:i)='i'
1163 pname(5:5)='i'
1164 endif
1165 if(cname(i:i).eq.';') then
1166 pname(i:i)='j'
1167 pname(5:5)='j'
1168 endif
1169 if(cname(i:i).eq.':') then
1170 pname(i:i)='k'
1171 pname(5:5)='k'
1172 endif
1173 if(cname(i:i).eq.',') then
1174 pname(i:i)='l'
1175 pname(5:5)='l'
1176 endif
1177 if(cname(i:i).eq.'<') then
1178 pname(i:i)='m'
1179 pname(5:5)='m'
1180 endif
1181 if(cname(i:i).eq.'>') then
1182 pname(i:i)='n'
1183 pname(5:5)='n'
1184 endif
1185 if(cname(i:i).eq.'!') then
1186 pname(i:i)='o'
1187 pname(5:5)='o'
1188 endif
1189 if(cname(i:i).eq.'@') then
1190 pname(i:i)='p'
1191 pname(5:5)='p'
1192 endif
1193 if(cname(i:i).eq.'#') then
1194 pname(i:i)='q'
1195 pname(5:5)='q'
1196 endif
1197 if(cname(i:i).eq.'$') then
1198 pname(i:i)='r'
1199 pname(5:5)='r'
1200 endif
1201 if(cname(i:i).eq.'%') then
1202 pname(i:i)='s'
1203 pname(5:5)='s'
1204 endif
1205 if(cname(i:i).eq.'^') then
1206 pname(i:i)='t'
1207 pname(5:5)='t'
1208 endif
1209 if(cname(i:i).eq.'&') then
1210 pname(i:i)='u'
1211 pname(5:5)='u'
1212 endif
1213 if(cname(i:i).eq.'(') then
1214 pname(i:i)='v'
1215 pname(5:5)='v'
1216 endif
1217 if(cname(i:i).eq.')') then
1218 pname(i:i)='x'
1219 pname(5:5)='x'
1220 endif
1221 if(cname(i:i).eq.'[') then
1222 pname(i:i)='y'
1223 pname(5:5)='y'
1224 endif
1225 if(cname(i:i).eq.']') then
1226 pname(i:i)='z'
1227 pname(5:5)='z'
1228 endif
1229 if(cname(i:i).eq.'{') then
1230 pname(i:i)='c'
1231 pname(5:5)='a'
1232 endif
1233 if(cname(i:i).eq.'}') then
1234 pname(i:i)='c'
1235 pname(5:5)='b'
1236 endif
1237 if(cname(i:i).eq.'=') then
1238 pname(i:i)='c'
1239 pname(5:5)='d'
1240 endif
1241 if(cname(i:i).eq.'~') then
1242 pname(i:i)='c'
1243 pname(5:5)='e'
1244 endif
1245 if(cname(i:i).eq.'|') then
1246 pname(i:i)='c'
1247 pname(5:5)='f'
1248 endif
1249 enddo
1250 if ((ichar(pname(1:1)).ge.48).and.
1251 + (ichar(pname(1:1)).le.57)) then
1252 pname='Z'//pname(1:lenocc(pname))
1253 endif
1254 end
1255
1256 subroutine cdnode(node)
1257 common/clevel/nodeold(20),nlevel
1258 character*16 anode
1259 character*16 cblank
1260 data cblank/' '/
1261 call toint(node,anode,ncd)
1262 if (nlevel.gt.1)then
1263* write(51,1000)cblank(1:nlevel-1),anode(1:ncd)
12641000 format(a,'Node',a,'->cd();')
1265 else
1266* write(51,1001)anode(1:ncd)
12671001 format('Node',a,'->cd();')
1268 endif
1269 end
1270 subroutine toint(i,s,nc)
1271 character*16 s1,s
1272 s1=' '
1273 write(s1,'(i7)')i
1274 do j=1,16
1275 if(s1(j:j).ne.' ') then
1276 j1 = j
1277 go to 10
1278 endif
1279 enddo
1280 10 continue
1281 do j=16,1,-1
1282 if (s1(j:j).ne.' ') then
1283 s=s1(j1:j)
1284 nc=j-j1+1
1285 return
1286 endif
1287 enddo
1288 end
1289 subroutine toreals(n,r,s,nc)
1290 character*(*) s
1291 character*14 s1
1292 dimension r(200)
1293 s=' '
1294 nc=0
1295 do 10 i=1,n
1296 call toreal(r(i),s1,nch)
1297 if(nc+nch.gt.128)then
1298* print *,'n=',n,' nc=',nc,' nch=',nch
1299 return
1300 endif
1301 s(nc+1:)=','
1302 s(nc+2:)=s1
1303 nc=nc+nch+1
1304 10 continue
1305 end
1306 subroutine toreal(r,s,nc)
1307 character*14 s
1308 character*14 s1
1309 s=' '
1310 if(r.eq.0)then
1311 s1='0'
1312 jbeg=1
1313 jend=1
1314 else
1315 write(s1,'(g14.7)')r
1316 do k=1,14
1317 if(s1(k:k).ne.' ') then
1318 jbeg=k
1319 goto 10
1320 endif
1321 enddo
1322 jbeg=1
1323 10 continue
1324 do k=14,jbeg,-1
1325 if(s1(k:k).ne.' '.and.s1(k:k).ne.'0') then
1326 jend=k
1327 goto 20
1328 endif
1329 enddo
1330 jend=14
1331 20 continue
1332 if(s1(jend:jend).eq.'.') jend=jend-1
1333 endif
1334 nc=jend-jbeg+1
1335 if(nc.le.0) then
1336 print *, 'Should never happen'
1337 endif
1338 s(1:nc)=s1(jbeg:jend)
1339 read(s(1:nc),*)t
1340 if(abs(t-r).gt.5e-7*abs(t+r)) print *, s(1:nc), t, r
1341 end
1342 subroutine toreal_old(r,s,nc)
1343 character*16 s1,s
1344 if(r.eq.0)then
1345 s='0'
1346 nc=1
1347 return
1348 endif
1349 s1=' '
1350 write(s1,'(f14.7)')r
1351 j1=1
1352 do j=1,16
1353 if(s1(j:j).ne.' ') then
1354 j1=j
1355 go to 10
1356 endif
1357 enddo
1358 10 continue
1359 j2=j1+7
1360 if(j2.gt.16)j2=16
1361 do 20 j=j2,j1+1,-1
1362 if (s1(j:j).eq.' ')go to 20
1363 if (s1(j:j).ne.'0') then
1364 if(s1(j:j).eq.'.')then
1365 s=s1(j1:j-1)
1366 nc=j-j1
1367 go to 30
1368 endif
1369 s=s1(j1:j)
1370 nc=j-j1+1
1371 go to 30
1372 endif
1373 20 continue
1374 30 continue
1375 if(nc.eq.1.and.s(1:1).eq.'-')then
1376 s='0'
1377 endif
1378 if(nc.eq.0)then
1379 nc=1
1380 s='0'
1381 endif
1382 end
1383C--------------------------------------------------------------
1384 subroutine create_map(fname)
1385 character*80 fname
1386
1387 parameter(max_maps = 1000)
1388 integer nmap,nalias,id(max_MAPS),ialias(max_MAPS)
1389 character*16 names(max_MAPS)
1390 common /maps_pool/ nmap,nalias,id,ialias,names
1391
1392 character*80 line
1393 character*1 comment,blank
1394 character*4 chid
1395 comment='#'
1396 blank=' '
1397 nmap = 0
1398 open(52,file=fname,status='old',err=100)
1399 10 read(52,300,err=100,END=100) line
1400 300 format(a)
1401 if(line(1:1) .eq. comment) goto 10
1402 if(line(1:1) .eq. blank) goto 10
1403 nmap = nmap + 1
1404 if(nmap .gt. max_maps) then
1405 nmap = nmap - 1
1406 print *, 'Warning: Number of names exceed maximum:',nmap
1407 print *, 'Warning: Rest of file are ignored'
1408 goto 100
1409 endif
1410 chid=line(1:4)
1411 ialias(nmap) = index(chid,'%')
1412 if(ialias(nmap).ne.0) then
1413 nalias = nalias + 1
1414 endif
1415 call uctoh(chid,id(nmap),4,4)
1416*
1417* Find substitution word
1418*
1419 l1=5
1420 ll=lenocc(line)
1421 do while(line(l1:l1).le.blank .and. l1.le.ll)
1422 l1=l1+1
1423 enddo
1424 if(l1.gt.ll) then
1425 nmap = nmap - 1
1426 goto 10
1427 endif
1428 l2=l1
1429 do while(line(l2:l2).gt.blank .and. l2.le.ll)
1430 l2=l2+1
1431 enddo
1432 if(l2.gt.ll) l2=ll
1433*************************************************************
1434 names(nmap) = line(l1:l2)
1435 if(ialias(nmap).ne.0) then
1436* Check that 'name' also has '%'
1437 if(index(names(nmap),'%').eq.0) then
1438 nmap = nmap - 1
1439 nalias = nalias - 1
1440 print *, '* Error in line ',line(1:l2)
1441 print *, '* Both words should have stab symbol "%"'
1442 print *, '*==> Line ignored ***'
1443 endif
1444 endif
1445* print 200,ialias(nmap),id(nmap),names(nmap)
1446* 200 format('Stab index:',i4,' ID = ',a4,' ==> ',a)
1447 goto 10
1448 100 continue
1449 close(52)
1450 end
1451C-------------------------------------------------------------------
1452 logical function map_found(idv,name)
1453 integer idv
1454 character*16 name
1455
1456 parameter(max_maps = 1000)
1457 integer nmap,nalias,id(max_maps),ialias(max_maps)
1458 character*16 names(max_maps)
1459 common /maps_pool/ nmap,nalias,id,ialias,names
1460
1461 character*4 chid,chidv,stab
1462 integer i,j,l
1463 map_found = .false.
1464* First step: Search only non-alias
1465 do i = 1,nmap
1466 if(ialias(i).eq.0 .and. id(i).eq.idv) then
1467 name = names(i)
1468 map_found = .true.
1469 return
1470 endif
1471 enddo
1472 if(nalias.eq.0) return
1473* Second step: Search only alias
1474 call uhtoc(idv,4,chidv,4)
1475 do i = 1,nmap
1476 if(ialias(i).eq.0) goto 100
1477 call uhtoc(id(i),4,chid,4)
1478 l = lenocc(chid)
1479 if(ialias(i).eq.1) then
1480* case 1 : %aaa
1481 if(l.eq.1) then ! single '%' match any name
1482 stab = chidv
1483 goto 50 ! accept
1484 endif
1485 if(chid(2:l) .ne. chidv(4-l+2:4)) goto 100
1486 stab = chidv(1:4-l+1)
1487 elseif(ialias(i).eq.l) then
1488* case 2 : aaa%
1489 if(chid(1:l-1) .ne. chidv(1:l-1)) goto 100
1490 stab = chidv(l:4)
1491 else
1492* case 3 : a%aa
1493 j=ialias(i) ! index of '%'
1494 if(chid(j+1:l) .ne. chidv(4-(l-j-1):4)) goto 100
1495 if(chid(1:j-1) .ne. chidv(1:j-1)) goto 100
1496 stab = chidv(j:4-(l-j-1)-1)
1497 endif
1498 50 continue
1499* print *, chidv,' matched to ',chid,' with stab ',stab
1500 name = names(i)
1501 j = index(name,'%')
1502 l = lenocc(name)
1503 if(j.eq.1 .and. l.eq.1) then
1504 name = stab
1505 elseif(j.eq.1) then
1506 name = stab//name(2:l)
1507 elseif(j.eq.l) then
1508 name = name(1:l-1)//stab
1509 else
1510 name = name(1:j-1)//stab//name(j+1:l)
1511 endif
1512 map_found = .true.
1513* print *, names(i),' ==> ',name
1514 return
1515 100 continue
1516 enddo
1517 end
1518C----------------------------------------------------------------------
1519 SUBROUTINE g2rin
1520C.
1521C. ******************************************************************
1522C. * *
1523C. * Routine to read GEANT object(s) fromin the RZ file *
1524C. * at the Current Working Directory (See RZCDIR) *
1525C. * The data structures from disk are read in memory *
1526C. * (VOLU,ROTM,TMED,MATE,SETS,PART,SCAN) *
1527C. * *
1528C. * This routine is s short-cut of the GEANT routine GRIN *
1529C. * Author R.Brun ********* *
1530C. * *
1531C. ******************************************************************
1532C.
1533*KEEP,HCBOOK.
1534 parameter(nwpaw=4000000)
1535 common/pawc/paw(nwpaw)
1536
1537 INTEGER IQ(2), LQ(8000)
1538 REAL Q(2)
1539 equivalence(lq(1),paw(11)),(iq(1),paw(19)),(q(1),paw(19))
1540
1541 INTEGER JDIGI ,JDRAW ,JHEAD ,JHITS ,JKINE ,JMATE ,JPART
1542 + ,JROTM ,JRUNG ,JSET ,JSTAK ,JGSTAT,JTMED ,JTRACK,JVERTX
1543 + ,JVOLUM,JXYZ ,JGPAR ,JGPAR2,JSKLT
1544C
1545 common/gclink/jdigi ,jdraw ,jhead ,jhits ,jkine ,jmate ,jpart
1546 + ,jrotm ,jrung ,jset ,jstak ,jgstat,jtmed ,jtrack,jvertx
1547 + ,jvolum,jxyz ,jgpar ,jgpar2,jsklt
1548C
1549 common/quest/iquest(100)
1550 parameter(nlinit=9,nmkey=22)
1551 dimension jnames(20),keys(2),linit(nlinit)
1552 dimension link(nmkey)
1553 equivalence(jnames(1),jdigi)
1554 CHARACTER*4 CKEY,KNAMES(NMKEY)
1555
1556 DATA knames/'DIGI','DRAW','HEAD','HITS','KINE','MATE','PART',
1557 + 'ROTM','RUNG','SETS','STAK','STAT','TMED','NULL','VERT',
1558 + 'VOLU','JXYZ','NULL','NULL','NULL','SCAN','NULL'/
1559 DATA linit/2,6,7,8,9,10,13,16,21/
1560C.
1561C. ------------------------------------------------------------------
1562C.
1563 print *,' In g2rin, iquest(13)=',iquest(13)
1564 iquest(1)=0
1565 idiv=2
1566 kvol=jvolum
1567 CALL vzero(jnames,20)
1568C
1569C Create a permanent link area for master pointers
1570C
1571 CALL mzlink(0,'/GCLINK/',jdigi,jsklt,jdigi)
1572*
1573 DO 20 j=1, nlinit
1574 link(j)=linit(j)
1575 20 CONTINUE
1576 nlink=nlinit
1577*
1578 ikey=0
1579 130 CONTINUE
1580 ikey=ikey+1
1581 CALL rzink(ikey,0,'S')
1582 print *,' after rzink, ikey=',ikey,'iquest(1)=',iquest(1)
1583 IF(iquest(1).NE.0) THEN
1584 print *,' nkeys=',iquest(7),' nwkey=',iquest(8)
1585 if (iquest(1).ne.11) then
1586 iquest(1)=0
1587 GOTO 150
1588 endif
1589 ENDIF
1590 indkey = iquest(21)
1591 CALL uhtoc(indkey,4,ckey,4)
1592 print *,'trying to read:',ckey
1593 DO 140 i=1,nlink
1594 nkey=abs(link(i))
1595 IF(ckey.EQ.knames(nkey))THEN
1596 keys(1)=iquest(21)
1597 keys(2)=iquest(22)
1598 IF(nkey.LE.20)THEN
1599 IF(jnames(nkey).NE.0)THEN
1600 CALL mzdrop(idiv,jnames(nkey),'L')
1601 jnames(nkey)=0
1602 ENDIF
1603 CALL rzin(idiv,jnames(nkey),1,keys,0,' ')
1604 ENDIF
1605 ENDIF
1606 140 CONTINUE
1607 GOTO 130
1608*
1609 150 nin=0
1610*
1611 999 END
#define a(i)
Definition RSha256.hxx:99
#define e(i)
Definition RSha256.hxx:103
#define quest
#define rzink
#define hlimit
#define pawc
void file()
Definition file.C:11
subroutine markdiv(ivo, nuserm)
Definition g2root.f:523
subroutine ptname(cname, pname)
Definition g2root.f:1117
logical function map_found(idv, name)
Definition g2root.f:1453
subroutine toroot(fname)
Definition g2root.f:245
subroutine node(ivo, nuserm, iposp)
Definition g2root.f:833
subroutine markdiv2(ivo, nuserm)
Definition g2root.f:518
subroutine toreal(r, s, nc)
Definition g2root.f:1307
subroutine realname(cname, pname)
Definition g2root.f:1085
subroutine toint(i, s, nc)
Definition g2root.f:1271
subroutine cdnode(node)
Definition g2root.f:1257
subroutine volume(cname, qjv, iposp, ifirst)
Definition g2root.f:614
subroutine node2(ivo, nuserm, iposp)
Definition g2root.f:828
subroutine create_map(fname)
Definition g2root.f:1385
program g2root
Definition g2root.f:6
subroutine toreals(n, r, s, nc)
Definition g2root.f:1290
subroutine toreal_old(r, s, nc)
Definition g2root.f:1343
subroutine g2rin
Definition g2root.f:1520
#define uhtoc
Definition h2root.cxx:129
subroutine cutol(chv)
Definition kernlib.f:22
subroutine uctoh(ms, mt, npw, nch)
Definition kernlib.f:198
subroutine vzero(a, n)
Definition kernlib.f:106
subroutine ucopy(a, b, n)
Definition kernlib.f:317
function lenocc(chv)
Definition kernlib.f:46
subroutine rzopen(lunin, chdir, cfname, choptt, lrecl, istat)
Definition zebra.f:461
subroutine rzin(ixdiv, lsup, jbias, keyu, icycle, chopt)
Definition zebra.f:873
subroutine rzfile(lunin, chdir, chopt)
Definition zebra.f:2799
subroutine mzlink(ixstor, chname, larea, lref, lrefl)
Definition zebra.f:2103
subroutine mzdrop(ixstor, lheadp, chopt)
Definition zebra.f:4274