xref: /phasta/phSolver/common/genblk.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen        subroutine genblk (IBKSZ)
2*59599516SKenneth E. Jansenc
3*59599516SKenneth E. Jansenc----------------------------------------------------------------------
4*59599516SKenneth E. Jansenc
5*59599516SKenneth E. Jansenc  This routine reads the interior elements and generates the
6*59599516SKenneth E. Jansenc  appropriate blocks.
7*59599516SKenneth E. Jansenc
8*59599516SKenneth E. Jansenc Zdenek Johan, Fall 1991.
9*59599516SKenneth E. Jansenc----------------------------------------------------------------------
10*59599516SKenneth E. Jansenc
11*59599516SKenneth E. Jansen        use pointer_data
12*59599516SKenneth E. Jansenc
13*59599516SKenneth E. Jansen        include "common.h"
14*59599516SKenneth E. Jansen!MR CHANGE
15*59599516SKenneth E. Jansen        include "mpif.h" !Required to determine the max for itpblk
16*59599516SKenneth E. Jansen!MR CHANGE END
17*59599516SKenneth E. Jansenc
18*59599516SKenneth E. Jansen        integer, allocatable :: ientp(:,:)
19*59599516SKenneth E. Jansen        integer mater(ibksz)
20*59599516SKenneth E. Jansen        integer intfromfile(50) ! integers read from headers
21*59599516SKenneth E. Jansen        character*255 fname1
22*59599516SKenneth E. Jansen
23*59599516SKenneth E. Jansencccccccccccccc New Phasta IO starts here ccccccccccccccccccccccccc
24*59599516SKenneth E. Jansen
25*59599516SKenneth E. Jansen        integer :: descriptor, descriptorG, GPID, color, nfiles
26*59599516SKenneth E. Jansen        integer ::  numparts, writeLock
27*59599516SKenneth E. Jansen        integer :: ierr_io, numprocs, itmp, itmp2
28*59599516SKenneth E. Jansen!MR CHANGE
29*59599516SKenneth E. Jansen        integer :: itpblktot,ierr,iseven
30*59599516SKenneth E. Jansen!MR CHANGE END
31*59599516SKenneth E. Jansen        character*255 fnamer, fname2, temp2
32*59599516SKenneth E. Jansen        character*64 temp1, temp3
33*59599516SKenneth E. Jansen!THIS NEEDS TO BE CLEANED - MR
34*59599516SKenneth E. Jansen        nfiles = nsynciofiles
35*59599516SKenneth E. Jansen!        nfields = nsynciofieldsreadgeombc
36*59599516SKenneth E. Jansen        numparts = numpe !This is the common settings. Beware if you try to compute several parts per process
37*59599516SKenneth E. Jansen
38*59599516SKenneth E. Jansen!        nppp = numparts/numpe
39*59599516SKenneth E. Jansen!        nppf = numparts/nfiles
40*59599516SKenneth E. Jansen
41*59599516SKenneth E. Jansen        color = int(myrank/(numparts/nfiles)) !Should call the SyncIO routine here
42*59599516SKenneth E. Jansen        itmp2 = int(log10(float(color+1)))+1
43*59599516SKenneth E. Jansen        write (temp2,"('(''geombc-dat.'',i',i1,')')") itmp2
44*59599516SKenneth E. Jansen        temp2=trim(temp2)
45*59599516SKenneth E. Jansen        write (fnamer,temp2) (color+1)
46*59599516SKenneth E. Jansen        fnamer=trim(fnamer)
47*59599516SKenneth E. Jansen
48*59599516SKenneth E. Jansen        ione=1
49*59599516SKenneth E. Jansen        itwo=2
50*59599516SKenneth E. Jansen        iseven=7
51*59599516SKenneth E. Jansen        ieleven=11
52*59599516SKenneth E. Jansen        itmp = int(log10(float(myrank+1)))+1
53*59599516SKenneth E. Jansen
54*59599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
55*59599516SKenneth E. Jansen
56*59599516SKenneth E. Jansenc
57*59599516SKenneth E. Jansen        iel=1
58*59599516SKenneth E. Jansen        itpblk=nelblk
59*59599516SKenneth E. Jansen!MR CHANGE
60*59599516SKenneth E. Jansen
61*59599516SKenneth E. Jansen        ! Get the total number of different interior topologies in the whole domain.
62*59599516SKenneth E. Jansen        ! Try to read from a field. If the field does not exist, scan the geombc file.
63*59599516SKenneth E. Jansen        itpblktot=-1
64*59599516SKenneth E. Jansen        write(temp1,
65*59599516SKenneth E. Jansen     &   "('(''total number of interior tpblocks@'',i',i1,',A1)')") itmp
66*59599516SKenneth E. Jansen
67*59599516SKenneth E. Jansen        write (fname2,temp1) (myrank+1),'?'
68*59599516SKenneth E. Jansen        call readheader(igeom,fname2 // char(0) ,itpblktot,ione,
69*59599516SKenneth E. Jansen     &  'integer' // char(0),iotype)
70*59599516SKenneth E. Jansen
71*59599516SKenneth E. Jansen!        write (*,*) 'Rank: ',myrank,' interior itpblktot intermediate:',
72*59599516SKenneth E. Jansen!     &               itpblktot
73*59599516SKenneth E. Jansen
74*59599516SKenneth E. Jansen        if (itpblktot == -1) then
75*59599516SKenneth E. Jansen          ! The field 'total number of different interior tpblocks' was not found in the geombc file.
76*59599516SKenneth E. Jansen          ! Scan all the geombc file for the 'connectivity interior' fields to get this information.
77*59599516SKenneth E. Jansen          iblk=0
78*59599516SKenneth E. Jansen          neltp=0
79*59599516SKenneth E. Jansen          do while(neltp .ne. -1)
80*59599516SKenneth E. Jansen
81*59599516SKenneth E. Jansen            ! intfromfile is reinitialized to -1 every time.
82*59599516SKenneth E. Jansen            ! If connectivity interior@xxx is not found, then
83*59599516SKenneth E. Jansen            ! readheader will return intfromfile unchanged
84*59599516SKenneth E. Jansen
85*59599516SKenneth E. Jansen            intfromfile(:)=-1
86*59599516SKenneth E. Jansen            iblk = iblk+1
87*59599516SKenneth E. Jansen            write (temp1,"('connectivity interior',i1)") iblk
88*59599516SKenneth E. Jansen            temp1 = trim(temp1)
89*59599516SKenneth E. Jansen            write (temp3,"('(''@'',i',i1,',A1)')") itmp
90*59599516SKenneth E. Jansen            write (fname2, temp3) (myrank+1), '?'
91*59599516SKenneth E. Jansen            fname2 = trim(temp1)//trim(fname2)
92*59599516SKenneth E. Jansen
93*59599516SKenneth E. Jansen            !write(*,*) 'rank, fname2',myrank, trim(adjustl(fname2))
94*59599516SKenneth E. Jansen            call readheader(igeom,fname2 // char(0),intfromfile,
95*59599516SKenneth E. Jansen     &       iseven,'integer' // char(0),iotype)
96*59599516SKenneth E. Jansen            neltp = intfromfile(1) ! -1 if fname2 was not found, >=0 otherwise
97*59599516SKenneth E. Jansen          end do
98*59599516SKenneth E. Jansen          itpblktot = iblk-1
99*59599516SKenneth E. Jansen        end if
100*59599516SKenneth E. Jansen
101*59599516SKenneth E. Jansen        if (myrank == 0) then
102*59599516SKenneth E. Jansen          write(*,*) 'Number of interior topologies: ',itpblktot
103*59599516SKenneth E. Jansen        endif
104*59599516SKenneth E. Jansen!        write (*,*) 'Rank: ',myrank,' interior itpblktot final:',
105*59599516SKenneth E. Jansen!     &               itpblktot
106*59599516SKenneth E. Jansen
107*59599516SKenneth E. Jansen!MR CHANGE END
108*59599516SKenneth E. Jansen
109*59599516SKenneth E. Jansen        nelblk=0
110*59599516SKenneth E. Jansen        mattyp = 0
111*59599516SKenneth E. Jansen        ndofl = ndof
112*59599516SKenneth E. Jansen        nsymdl = nsymdf
113*59599516SKenneth E. Jansen
114*59599516SKenneth E. Jansen!        call initphmpiio( nfields, nppf, nfiles, igeom )
115*59599516SKenneth E. Jansen!        call openfile( fnamer, 'read', igeom )
116*59599516SKenneth E. Jansen
117*59599516SKenneth E. Jansen!         do iblk = 1, itpblk
118*59599516SKenneth E. Jansen        do iblk = 1, itpblktot
119*59599516SKenneth E. Jansen           writeLock=0;
120*59599516SKenneth E. Jansen!MR CHANGE END
121*59599516SKenneth E. Jansenc
122*59599516SKenneth E. Jansenc           read(igeomBAK) neltp,nenl,ipordl,nshl, ijunk, ijunk, lcsyst
123*59599516SKenneth E. Jansenc           call creadlist(igeomBAK,iseven,
124*59599516SKenneth E. Jansenc     &          neltp,nenl,ipordl,nshl, ijunk, ijunk, lcsyst)
125*59599516SKenneth E. Jansen
126*59599516SKenneth E. Jansenccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
127*59599516SKenneth E. Jansen
128*59599516SKenneth E. Jansen           write (temp1,"('connectivity interior',i1)") iblk
129*59599516SKenneth E. Jansen           temp1=trim(temp1)
130*59599516SKenneth E. Jansen           write (temp3,"('(''@'',i',i1,',A1)')") itmp
131*59599516SKenneth E. Jansen           write (fname2, temp3) (myrank+1), '?'
132*59599516SKenneth E. Jansen           fname2 = trim(temp1)//trim(fname2)
133*59599516SKenneth E. Jansen
134*59599516SKenneth E. Jansenccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
135*59599516SKenneth E. Jansen
136*59599516SKenneth E. Jansenc           fname1='connectivity interior?'
137*59599516SKenneth E. Jansen
138*59599516SKenneth E. Jansen           ! Synchronization for performance monitoring, as some parts do not include some topologies
139*59599516SKenneth E. Jansen           call MPI_Barrier(MPI_COMM_WORLD,ierr)
140*59599516SKenneth E. Jansen           call readheader(igeom,fname2 // char(0) ,intfromfile,
141*59599516SKenneth E. Jansen     &     iseven,"integer" // char(0), iotype)
142*59599516SKenneth E. Jansen           neltp  =intfromfile(1)
143*59599516SKenneth E. Jansen           nenl   =intfromfile(2)
144*59599516SKenneth E. Jansen           ipordl =intfromfile(3)
145*59599516SKenneth E. Jansen           nshl   =intfromfile(4)
146*59599516SKenneth E. Jansen           ijunk  =intfromfile(5)
147*59599516SKenneth E. Jansen           ijunk  =intfromfile(6)
148*59599516SKenneth E. Jansen           lcsyst =intfromfile(7)
149*59599516SKenneth E. Jansen           allocate (ientp(neltp,nshl))
150*59599516SKenneth E. Jansenc           read(igeomBAK) ientp
151*59599516SKenneth E. Jansen           iientpsiz=neltp*nshl
152*59599516SKenneth E. Jansen
153*59599516SKenneth E. Jansen           if (neltp==0) then
154*59599516SKenneth E. Jansen              writeLock=1;
155*59599516SKenneth E. Jansen           endif
156*59599516SKenneth E. Jansen
157*59599516SKenneth E. Jansen           call readdatablock(igeom,fname2 // char(0),ientp,iientpsiz,
158*59599516SKenneth E. Jansen     &                     "integer" // char(0), iotype)
159*59599516SKenneth E. Jansen
160*59599516SKenneth E. Jansen!            call closefile( igeom, "read" // char(0) )
161*59599516SKenneth E. Jansen!            call finalizephmpiio( igeom )
162*59599516SKenneth E. Jansen
163*59599516SKenneth E. Jansen!MR CHANGE
164*59599516SKenneth E. Jansen           if(writeLock==0) then
165*59599516SKenneth E. Jansen!MR CHANGE
166*59599516SKenneth E. Jansen
167*59599516SKenneth E. Jansen             do n=1,neltp,ibksz
168*59599516SKenneth E. Jansen                nelblk=nelblk+1
169*59599516SKenneth E. Jansen                npro= min(IBKSZ, neltp - n + 1)
170*59599516SKenneth E. Jansenc
171*59599516SKenneth E. Jansen                lcblk(1,nelblk)  = iel
172*59599516SKenneth E. Jansenc                lcblk(2,nelblk)  = iopen ! available for later use
173*59599516SKenneth E. Jansen                lcblk(3,nelblk)  = lcsyst
174*59599516SKenneth E. Jansen                lcblk(4,nelblk)  = ipordl
175*59599516SKenneth E. Jansen                lcblk(5,nelblk)  = nenl
176*59599516SKenneth E. Jansen                lcblk(6,nelblk)  = nfacel
177*59599516SKenneth E. Jansen                lcblk(7,nelblk)  = mattyp
178*59599516SKenneth E. Jansen                lcblk(8,nelblk)  = ndofl
179*59599516SKenneth E. Jansen                lcblk(9,nelblk)  = nsymdl
180*59599516SKenneth E. Jansen                lcblk(10,nelblk) = nshl ! # of shape functions per elt
181*59599516SKenneth E. Jansenc
182*59599516SKenneth E. Jansenc.... allocate memory for stack arrays
183*59599516SKenneth E. Jansenc
184*59599516SKenneth E. Jansen                allocate (mmat(nelblk)%p(npro))
185*59599516SKenneth E. Jansenc
186*59599516SKenneth E. Jansen                allocate (mien(nelblk)%p(npro,nshl))
187*59599516SKenneth E. Jansen                allocate (mxmudmi(nelblk)%p(npro,maxsh))
188*59599516SKenneth E. Jansenc
189*59599516SKenneth E. Jansenc.... save the element block
190*59599516SKenneth E. Jansenc
191*59599516SKenneth E. Jansen                n1=n
192*59599516SKenneth E. Jansen                n2=n+npro-1
193*59599516SKenneth E. Jansen                mater=1   ! all one material for now
194*59599516SKenneth E. Jansen                call gensav (ientp(n1:n2,1:nshl),
195*59599516SKenneth E. Jansen     &                       mater,           mien(nelblk)%p,
196*59599516SKenneth E. Jansen     &                       mmat(nelblk)%p)
197*59599516SKenneth E. Jansen                iel=iel+npro
198*59599516SKenneth E. Jansenc
199*59599516SKenneth E. Jansen             enddo
200*59599516SKenneth E. Jansen!MR CHANGE
201*59599516SKenneth E. Jansen           endif
202*59599516SKenneth E. Jansen!MR CHANGE
203*59599516SKenneth E. Jansen           deallocate(ientp)
204*59599516SKenneth E. Jansen        enddo
205*59599516SKenneth E. Jansen
206*59599516SKenneth E. Jansen!        call closefile( igeom, "read" // char(0) )
207*59599516SKenneth E. Jansen!        call finalizephmpiio( igeom )
208*59599516SKenneth E. Jansen
209*59599516SKenneth E. Jansen        lcblk(1,nelblk+1) = iel
210*59599516SKenneth E. Jansenc
211*59599516SKenneth E. Jansenc.... return
212*59599516SKenneth E. Jansenc
213*59599516SKenneth E. JansenCAD        call timer ('Back    ')
214*59599516SKenneth E. Jansenc
215*59599516SKenneth E. Jansen        return
216*59599516SKenneth E. Jansenc
217*59599516SKenneth E. Jansen1000    format(a80,//,
218*59599516SKenneth E. Jansen     &  ' N o d a l   C o n n e c t i v i t y',//,
219*59599516SKenneth E. Jansen     &  '   Elem  ',/,
220*59599516SKenneth E. Jansen     &  '  Number  ',7x,27('Node',i2,:,2x))
221*59599516SKenneth E. Jansen1100    format(2x,i5,6x,27i8)
222*59599516SKenneth E. Jansen        end
223