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