source: trunk/source/g3tog4/src/tog4.F @ 1249

Last change on this file since 1249 was 965, checked in by garnier, 15 years ago

update g3tog4

File size: 9.8 KB
Line 
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,v 1.5 2006/06/29 18:15:21 gunter Exp $
28*     GEANT4 tag $Name: geant4-09-02-ref-02 $
29*   
30      subroutine tog4
31************************************************************************
32*     
33*     tog4
34*     
35*     Perform the translation to G4
36*     
37************************************************************************
38      implicit none
39#include "gcbank.inc"
40      integer maxdivols
41      parameter (maxdivols=20000)
42      integer nvol, nrotm, nmate, ntmed, nset, i, jma, nmixt, k, nin,
43     >     jdiv, jd, iaxis, ivo, ndiv, numed, npar, natt, ivol, jin,
44     >     nparv, npard, nr, irot, konly, nwbuf, isvol, nmat, ifield,
45     >     nbits(5000), idtyp, nwhi, nwdi, iset, idet, j, in, jmx,
46     >     jdh, jdd, jdu, ndet, nn, nupar, npos, ndvol, ndivols, ii,
47     >     npositioned, iia(10000), imate, smixt
48
49      real c0, step, x, y, z, a, dens, radl, absl, fact(5000),
50     >     fieldm, tmaxfd, stemax, deemax, epsil, stmin, orig(5000),
51     >     upar(5000)
52
53      character shape*4, name*4, dname*4, chonly*4, chmat*20, chtmed*20,
54     >     chset*4, chdet*4, chnms(5000)*4, divols(maxdivols)*4
55*     
56      npositioned = 0
57*     
58***   count materials and convert
59      call bankcnt(jmate,iia, nmate)
60      print *,'Materials: ',nmate
61      do imate=1,nmate
62         ii=iia(imate)
63         jma = lq(jmate-ii)
64         call uhtoc(iq(jma+1),4,chmat,20)
65         a = q(jma+6)
66         z = q(jma+7)
67         dens = q(jma+8)
68         radl = q(jma+9)
69         absl = q(jma+10)
70         nwbuf = iq(jma-1)-11
71         if (jma.gt.0) then
72            smixt=q(jma+11)
73            nmixt=abs(smixt)
74            if (nmixt.le.1) then
75               write(6,101) imate, chmat, a, z, dens, radl, absl
76               call ksmate(ii, chmat, a, z, dens, radl, absl,
77     >              q(jma+12), nwbuf)
78            else
79               jmx = lq(jma-5)
80               write(6,102) imate, chmat, a, z, dens, radl, absl,
81     >              (j,q(jmx+j),q(jmx+nmixt+j),q(jmx+2*nmixt+j),
82     >              j=1,nmixt)
83               call ksmixt(ii, chmat, q(jmx+1), q(jmx+nmixt+1),
84     >              dens, smixt, q(jmx+2*nmixt+1))
85            end if
86         end if
87      enddo
88 101  format(1x,i5,1x,A12,f6.2,f5.1,f8.2,2f9.2)
89 102  format(1x,i5,1x,A12,f6.2,f5.1,f8.2,2f9.2,1x,i2, f6.2, f5.1,
90     >     f6.2/(57x, i2, f6.2, f5.1, f6.2))
91*     
92***   count tracking media and convert
93      call bankcnt(jtmed,iia, ntmed)
94      print *,'Media: ',ntmed
95      do i=1,ntmed
96         ii=iia(i)
97         j = lq(jtmed-ii)
98         call uhtoc(iq(j+1),4,chtmed,20)
99         nmat = q(j+6)
100         isvol = q(j+7)
101         ifield = q(j+8)
102         fieldm = q(j+9)
103         tmaxfd = q(j+10)
104         stemax = q(j+11)
105         deemax = q(j+12)
106         epsil = q(j+13)
107         stmin = q(j+14)
108         nwbuf = iq(j-1) -14
109         call kstmed(ii,chtmed,nmat,isvol,ifield,fieldm,tmaxfd,stemax,
110     +        deemax,epsil,stmin,q(j+15),nwbuf)
111      enddo
112*     
113***   count rotation matrices and convert
114      call bankcnt(jrotm,iia, nrotm)
115      print *,'Rotms: ',nrotm
116      do i=1,nrotm
117         ii=iia(i)
118         j = lq(jrotm-ii)
119         call ksrotm(ii,q(j+11),q(j+12),q(j+13),q(j+14),q(j+15),q(j+16))
120      enddo
121*     
122***   count volumes
123      npos = 0
124      call bankcnt(jvolum,iia, nvol)
125      print *,'Volumes: ',nvol
126***   pull out the names of the volumes which are subvolumes of
127***   divided volumes (gsvolu should not be called on these)
128      ndivols = 0
129      do i=1, nvol
130         ii=iia(i)
131         j = lq(jvolum-ii)
132         nin = q(j+3)
133         if (nin.lt.0) then
134            jdiv = lq(j-1)
135            ivo = q(jdiv+2)
136            call uhtoc(iq(jvolum+ivo),4,dname,4)
137            ndivols = ndivols +1
138            if (ndivols.gt.maxdivols) then
139               ndivols = maxdivols
140               print *,
141     +              '!!!ERROR!!! ndivols array exhausted. ',
142     +              'Too many divisions.'
143            endif
144            divols(ndivols) = dname
145         endif
146      enddo
147***   create the logical volumes (gsvolu)
148      ndvol = 0
149      do i=1, nvol
150         ii=iia(i)
151         j = lq(jvolum-ii)
152         call uhtoc(iq(jvolum+ii),4,name,4)
153         call jshape(q(j+2),shape)
154         nin = q(j+3)
155         numed = q(j+4)
156         npar = q(j+5)
157         natt = q(j+6)
158         do k=1, ndivols
159            if (divols(k).eq.name) then
160               ndvol = ndvol +1
161c     print *,'Division volume ',name,'; no gsvolu call.'
162               goto 11
163            endif
164         enddo
165         call ksvolu(name, shape, numed, q(j+7), npar, ivol)
166 11      continue
167      enddo
168      print *,'Divided volumes: ',ndvol
169
170***   properties of the mother volume
171      call uhtoc(iq(jvolum+1),4,name,4)
172      j=lq(jvolum-1)
173      call jshape(q(j+2),shape)
174      print *,'mother volume: ',name,' shape: ',shape
175***   convert physical volumes
176      do i=1, nvol
177         ii=iia(i)
178         j = lq(jvolum-ii)
179         call uhtoc(iq(jvolum+ii),4,name,4)
180         nin = q(j+3)
181         numed = q(j+4)
182         npar = q(j+5)
183         if (nin.lt.0) then
184*     ! divided volume
185            jdiv = lq(j-1)
186            iaxis = q(jdiv+1)
187            ivo = q(jdiv+2)
188            call uhtoc(iq(jvolum+ivo),4,dname,4)
189            jd = lq(jvolum-ivo)
190            numed = q(jd+4)
191            ndiv = q(jdiv+3)
192            c0 = q(jdiv+4)
193            step = q(jdiv+5)
194            call ksdvn2(dname, name, ndiv, iaxis, c0, numed)
195         else if (nin.gt.0) then
196*     ! volume not divided. Handle positioning of daughter vols
197            do in=1, nin
198               jin = lq(j-in)
199               ivo = q(jin+2)
200               call uhtoc(iq(jvolum+ivo),4,dname,4)
201               jd = lq(jvolum-ivo)
202               nparv = q(jd+5)  ! NPAR declared in the GSVOLU call
203               nr = q(jin+3)
204               irot = q(jin+4)
205               x = q(jin+5)
206               y = q(jin+6)
207               z = q(jin+7)
208               konly = q(jin+8)
209               if (konly.ne.0) then
210                  chonly = 'ONLY'
211               else
212                  chonly = 'MANY'
213               endif
214               npard = iq(jin-1) -9
215               npositioned = npositioned +1
216               if (nparv.eq.0) then
217*     ! use GSPOSP
218                  call ksposp(dname, nr, name, x, y, z, irot, chonly,
219     +                 q(jin+10), npard)
220               else
221*     ! GSPOS
222                  call kspos(dname, nr, name, x, y, z, irot, chonly)
223               endif
224            enddo
225         endif
226      enddo
227*     
228***   count sensitive detectors
229      call bankcnt(jset,iia, nset)
230      print *,'Sets: ',nset
231      do i=1,nset
232         ii=iia(i)
233         j = lq(jset-ii)
234         call uhtoc(iq(jset+i),4,chset,4)
235         ndet = iq(j-1)
236         do k=1,ndet
237            jd = lq(j-k)
238            call uhtoc(iq(j+k),4,chdet,4)
239            call gfdet(chset, chdet, nn, chnms, nbits, idtyp,
240     +           nwhi, nwdi, iset, idet)
241            call ksdet(chset, chdet, nn, chnms, nbits, idtyp,
242     +           nwhi, nwdi, iset, idet)
243            jdh = lq(jd-1)
244            if (jdh.ne.0) then
245               call gfdeth(chset,chdet,nn,chnms,nbits,orig,fact)
246               call ksdeth(chset,chdet,nn,chnms,nbits,orig,fact)
247            endif
248            jdd = lq(jd-2)
249            if (jdd.ne.0) then
250               call gfdetd(chset,chdet,nn,chnms,nbits)
251               call ksdetd(chset,chdet,nn,chnms,nbits)
252            endif
253            jdu = lq(jd-3)
254            if (jdu.ne.0) then
255               call gfdetu(chset,chdet,100,nupar,upar)
256               call ksdetu(chset,chdet,nupar,upar)
257            endif
258         enddo
259      enddo
260      print *,'Positioned volumes (gspos, gsposp):',npositioned
261*     
262      call kgclos
263*     
264      end
265
266      subroutine bankcnt(link,iia,nbanks)
267************************************************************************
268************************************************************************
269      implicit none
270#include "gcbank.inc"
271      integer i, link, nbanks, iia(*)
272*     
273      nbanks=0
274      if (link.eq.0) return
275C*      do i=1,9999999
276      do i=1,iq(link-2)
277C*         if(lq(link-nbanks-1).eq.0.or.iq(link-2).eq.nbanks) goto 10
278         if(lq(link-i).ne.0)then
279           nbanks = nbanks +1
280           iia(nbanks)=i
281         endif
282      enddo
283 10   continue
284      end
Note: See TracBrowser for help on using the repository browser.