Geant4_10
tog4.F
Go to the documentation of this file.
1 *
2 * ********************************************************************
3 * * License and Disclaimer *
4 * * *
5 * * The Geant4 software is copyright of the Copyright Holders of *
6 * * the Geant4 Collaboration. It is provided under the terms and *
7 * * conditions of the Geant4 Software License, included in the file *
8 * * LICENSE and available at http://cern.ch/geant4/license . These *
9 * * include a list of copyright holders. *
10 * * *
11 * * Neither the authors of this software system, nor their employing *
12 * * institutes,nor the agencies providing financial support for this *
13 * * work make any representation or warranty, express or implied, *
14 * * regarding this software system or assume any liability for its *
15 * * use. Please see the license in the file LICENSE and URL above *
16 * * for the full disclaimer and the limitation of liability. *
17 * * *
18 * * This code implementation is the result of the scientific and *
19 * * technical work of the GEANT4 collaboration. *
20 * * By using, copying, modifying or distributing the software (or *
21 * * any work based on the software) you agree to acknowledge its *
22 * * use in resulting scientific publications, and indicate your *
23 * * acceptance of all terms of the Geant4 Software license. *
24 * ********************************************************************
25 *
26 *
27 * $Id: tog4.F 67982 2013-03-13 10:36:03Z gcosmo $
28 *
29  subroutine tog4
30 ************************************************************************
31 *
32 * tog4
33 *
34 * Perform the translation to G4
35 *
36 ************************************************************************
37  implicit none
38 #include "gcbank.inc"
39  integer maxdivols
40  parameter(maxdivols=20000)
41  integer nvol, nrotm, nmate, ntmed, nset, i, jma, nmixt, k, nin,
42  > jdiv, jd, iaxis, ivo, ndiv, numed, npar, natt, ivol, jin,
43  > nparv, npard, nr, irot, konly, nwbuf, isvol, nmat, ifield,
44  > nbits(5000), idtyp, nwhi, nwdi, iset, idet, j, in, jmx,
45  > jdh, jdd, jdu, ndet, nn, nupar, npos, ndvol, ndivols, ii,
46  > npositioned, iia(10000), imate, smixt
47 
48  real c0, step, x, y, z, a, dens, radl, absl, fact(5000),
49  > fieldm, tmaxfd, stemax, deemax, epsil, stmin, orig(5000),
50  > upar(5000)
51 
52  character shape*4, name*4, dname*4, chonly*4, chmat*20, chtmed*20,
53  > chset*4, chdet*4, chnms(5000)*4, divols(maxdivols)*4
54 *
55  npositioned = 0
56 *
57 *** count materials and convert
58  call bankcnt(jmate,iia, nmate)
59  print *,'Materials: ',nmate
60  do imate=1,nmate
61  ii=iia(imate)
62  jma = lq(jmate-ii)
63  call uhtoc(iq(jma+1),4,chmat,20)
64  a = q(jma+6)
65  z = q(jma+7)
66  dens = q(jma+8)
67  radl = q(jma+9)
68  absl = q(jma+10)
69  nwbuf = iq(jma-1)-11
70  if (jma.gt.0) then
71  smixt=q(jma+11)
72  nmixt=abs(smixt)
73  if (nmixt.le.1) then
74  write(6,101) imate, chmat, a, z, dens, radl, absl
75  call ksmate(ii, chmat, a, z, dens, radl, absl,
76  > q(jma+12), nwbuf)
77  else
78  jmx = lq(jma-5)
79  write(6,102) imate, chmat, a, z, dens, radl, absl,
80  > (j,q(jmx+j),q(jmx+nmixt+j),q(jmx+2*nmixt+j),
81  > j=1,nmixt)
82  call ksmixt(ii, chmat, q(jmx+1), q(jmx+nmixt+1),
83  > dens, smixt, q(jmx+2*nmixt+1))
84  end if
85  end if
86  enddo
87  101 format(1x,i5,1x,a12,f6.2,f5.1,f8.2,2f9.2)
88  102 format(1x,i5,1x,a12,f6.2,f5.1,f8.2,2f9.2,1x,i2, f6.2, f5.1,
89  > f6.2/(57x, i2, f6.2, f5.1, f6.2))
90 *
91 *** count tracking media and convert
92  call bankcnt(jtmed,iia, ntmed)
93  print *,'Media: ',ntmed
94  do i=1,ntmed
95  ii=iia(i)
96  j = lq(jtmed-ii)
97  call uhtoc(iq(j+1),4,chtmed,20)
98  nmat = q(j+6)
99  isvol = q(j+7)
100  ifield = q(j+8)
101  fieldm = q(j+9)
102  tmaxfd = q(j+10)
103  stemax = q(j+11)
104  deemax = q(j+12)
105  epsil = q(j+13)
106  stmin = q(j+14)
107  nwbuf = iq(j-1) -14
108  call kstmed(ii,chtmed,nmat,isvol,ifield,fieldm,tmaxfd,stemax,
109  + deemax,epsil,stmin,q(j+15),nwbuf)
110  enddo
111 *
112 *** count rotation matrices and convert
113  call bankcnt(jrotm,iia, nrotm)
114  print *,'Rotms: ',nrotm
115  do i=1,nrotm
116  ii=iia(i)
117  j = lq(jrotm-ii)
118  call ksrotm(ii,q(j+11),q(j+12),q(j+13),q(j+14),q(j+15),q(j+16))
119  enddo
120 *
121 *** count volumes
122  npos = 0
123  call bankcnt(jvolum,iia, nvol)
124  print *,'Volumes: ',nvol
125 *** pull out the names of the volumes which are subvolumes of
126 *** divided volumes (gsvolu should not be called on these)
127  ndivols = 0
128  do i=1, nvol
129  ii=iia(i)
130  j = lq(jvolum-ii)
131  nin = q(j+3)
132  if (nin.lt.0) then
133  jdiv = lq(j-1)
134  ivo = q(jdiv+2)
135  call uhtoc(iq(jvolum+ivo),4,dname,4)
136  ndivols = ndivols +1
137  if (ndivols.gt.maxdivols) then
138  ndivols = maxdivols
139  print *,
140  + '!!!ERROR!!! ndivols array exhausted. ',
141  + 'Too many divisions.'
142  endif
143  divols(ndivols) = dname
144  endif
145  enddo
146 *** create the logical volumes (gsvolu)
147  ndvol = 0
148  do i=1, nvol
149  ii=iia(i)
150  j = lq(jvolum-ii)
151  call uhtoc(iq(jvolum+ii),4,name,4)
152  call jshape(q(j+2),shape)
153  nin = q(j+3)
154  numed = q(j+4)
155  npar = q(j+5)
156  natt = q(j+6)
157  do k=1, ndivols
158  if (divols(k).eq.name) then
159  ndvol = ndvol +1
160 c print *,'Division volume ',name,'; no gsvolu call.'
161  goto 11
162  endif
163  enddo
164  call ksvolu(name, shape, numed, q(j+7), npar, ivol)
165  11 continue
166  enddo
167  print *,'Divided volumes: ',ndvol
168 
169 *** properties of the mother volume
170  call uhtoc(iq(jvolum+1),4,name,4)
171  j=lq(jvolum-1)
172  call jshape(q(j+2),shape)
173  print *,'mother volume: ',name,' shape: ',shape
174 *** convert physical volumes
175  do i=1, nvol
176  ii=iia(i)
177  j = lq(jvolum-ii)
178  call uhtoc(iq(jvolum+ii),4,name,4)
179  nin = q(j+3)
180  numed = q(j+4)
181  npar = q(j+5)
182  if (nin.lt.0) then
183 * ! divided volume
184  jdiv = lq(j-1)
185  iaxis = q(jdiv+1)
186  ivo = q(jdiv+2)
187  call uhtoc(iq(jvolum+ivo),4,dname,4)
188  jd = lq(jvolum-ivo)
189  numed = q(jd+4)
190  ndiv = q(jdiv+3)
191  c0 = q(jdiv+4)
192  step = q(jdiv+5)
193  call ksdvn2(dname, name, ndiv, iaxis, c0, numed)
194  else if (nin.gt.0) then
195 * ! volume not divided. Handle positioning of daughter vols
196  do in=1, nin
197  jin = lq(j-in)
198  ivo = q(jin+2)
199  call uhtoc(iq(jvolum+ivo),4,dname,4)
200  jd = lq(jvolum-ivo)
201  nparv = q(jd+5) ! NPAR declared in the GSVOLU call
202  nr = q(jin+3)
203  irot = q(jin+4)
204  x = q(jin+5)
205  y = q(jin+6)
206  z = q(jin+7)
207  konly = q(jin+8)
208  if (konly.ne.0) then
209  chonly = 'ONLY'
210  else
211  chonly = 'MANY'
212  endif
213  npard = iq(jin-1) -9
214  npositioned = npositioned +1
215  if (nparv.eq.0) then
216 * ! use GSPOSP
217  call ksposp(dname, nr, name, x, y, z, irot, chonly,
218  + q(jin+10), npard)
219  else
220 * ! GSPOS
221  call kspos(dname, nr, name, x, y, z, irot, chonly)
222  endif
223  enddo
224  endif
225  enddo
226 *
227 *** count sensitive detectors
228  call bankcnt(jset,iia, nset)
229  print *,'Sets: ',nset
230  do i=1,nset
231  ii=iia(i)
232  j = lq(jset-ii)
233  call uhtoc(iq(jset+i),4,chset,4)
234  ndet = iq(j-1)
235  do k=1,ndet
236  jd = lq(j-k)
237  call uhtoc(iq(j+k),4,chdet,4)
238  call gfdet(chset, chdet, nn, chnms, nbits, idtyp,
239  + nwhi, nwdi, iset, idet)
240  call ksdet(chset, chdet, nn, chnms, nbits, idtyp,
241  + nwhi, nwdi, iset, idet)
242  jdh = lq(jd-1)
243  if (jdh.ne.0) then
244  call gfdeth(chset,chdet,nn,chnms,nbits,orig,fact)
245  call ksdeth(chset,chdet,nn,chnms,nbits,orig,fact)
246  endif
247  jdd = lq(jd-2)
248  if (jdd.ne.0) then
249  call gfdetd(chset,chdet,nn,chnms,nbits)
250  call ksdetd(chset,chdet,nn,chnms,nbits)
251  endif
252  jdu = lq(jd-3)
253  if (jdu.ne.0) then
254  call gfdetu(chset,chdet,100,nupar,upar)
255  call ksdetu(chset,chdet,nupar,upar)
256  endif
257  enddo
258  enddo
259  print *,'Positioned volumes (gspos, gsposp):',npositioned
260 *
261  call kgclos
262 *
263  end
264 
265  subroutine bankcnt(link,iia,nbanks)
266 ************************************************************************
267 ************************************************************************
268  implicit none
269 #include "gcbank.inc"
270  integer i, link, nbanks, iia(*)
271 *
272  nbanks=0
273  if (link.eq.0) return
274 C* do i=1,9999999
275  do i=1,iq(link-2)
276 C* if(lq(link-nbanks-1).eq.0.or.iq(link-2).eq.nbanks) goto 10
277  if(lq(link-i).ne.0)then
278  nbanks = nbanks +1
279  iia(nbanks)=i
280  endif
281  enddo
282  10 continue
283  end
Double_t z
Definition: plot.C:279
subroutine ksposp(name, num, moth, x, y, z, irot, only, par, npar)
Definition: g3routines.F:113
subroutine ksmixt(imate, name, a, z, dens, nlmat, wmat)
Definition: g3routines.F:464
ifstream in
Definition: comparison.C:7
subroutine ksmate(imate, name, a, z, dens, radl, absl, ubf, nwbf)
Definition: g3routines.F:421
subroutine ksdetu(chset, chdet, nupar, upar)
Definition: g3routines.F:886
const XML_Char * name
Definition: expat.h:151
G4double a
Definition: TRTMaterials.hh:39
Double_t y
Definition: plot.C:279
subroutine ksvolu(name, shape, nmed, par, npar, ivol)
Definition: g3routines.F:34
subroutine bankcnt(link, iia, nbanks)
Definition: tog4.F:265
subroutine ksrotm(irot, theta1, phi1, theta2, phi2, theta3, phi3)
Definition: g3routines.F:202
subroutine kspos(name, num, moth, x, y, z, irot, only)
Definition: g3routines.F:76
Double_t x
Definition: plot.C:279
subroutine ksdet(chset, chdet, nv, chnam, nbits, idtyp, nwhi, nwdi, iset, idet)
Definition: g3routines.F:690
subroutine ksdetd(chset, chdet, nd, chnam, nbits)
Definition: g3routines.F:845
subroutine ksdvn2(name, moth, ndiv, iaxis, c0, numed)
Definition: g3routines.F:350
subroutine jshape(RSHAPE, SHAPE)
Definition: jshape.F:32
subroutine kstmed(itmed, name, nmat, isvol, ifield, fieldm, tmaxfd, stemax, deemax, epsil, stmin, ubuf, nwbuf)
Definition: g3routines.F:514
subroutine kgclos
Definition: g3routines.F:924
subroutine tog4
Definition: tog4.F:29
void print(const std::vector< T > &data)
Definition: DicomRun.hh:111
subroutine ksdeth(chset, chdet, nh, chnam, nbits, orig, fact)
Definition: g3routines.F:798