xref: /phasta/phSolver/common/genbkbSyncIO.f (revision 9d714148f3212e057272549b3ee56361da59830b)
1        subroutine genbkb (ibksz)
2c
3c----------------------------------------------------------------------
4c
5c  This routine reads the boundary elements, reorders them and
6c  generates traces for the gather/scatter operations.
7c
8c Zdenek Johan, Fall 1991.
9c----------------------------------------------------------------------
10c
11        use dtnmod
12        use pointer_data
13c
14        include "common.h"
15        include "mpif.h" !Required to determine the max for itpblk
16c
17
18        integer, allocatable :: ientp(:,:),iBCBtp(:,:)
19        real*8, allocatable :: BCBtp(:,:)
20        integer materb(ibksz)
21        integer intfromfile(50) ! integers read from headers
22        character*255 fname1
23
24cccccccccccccc New Phasta IO starts here cccccccccccccccccccccccccccccc
25
26        integer :: descriptor, descriptorG, GPID, color, nfiles, nfields
27        integer :: numparts, nppf, nppp, nprocs, writeLock
28        integer :: ierr_io, numprocs, itmp, itmp2
29        integer :: itpblktot,ierr
30
31        character*255 fnamer, fname2, temp2
32        character*64 temp1, temp3
33
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))
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        ieight=8
51        ieleven=11
52        itmp = int(log10(float(myrank+1)))+1
53
54
55cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
56
57        iel=1
58        itpblk=nelblb
59
60
61!MR CHANGE
62        ! Get the total number of different interior topologies in the whole domain.
63        ! Try to read from a field. If the field does not exist, scan the geombc file.
64        itpblktot=-1
65        write(temp1,
66     &   "('(''total number of boundary tpblocks@'',i',i1,',A1)')") itmp
67        write (fname2,temp1) (myrank+1),'?'
68        call readheader(igeom,fname2 // char(0) ,itpblktot,ione,
69     &  'integer' // char(0),iotype)
70
71!        write (*,*) 'Rank: ',myrank,' boundary itpblktot intermediate:',
72!     &               itpblktot
73
74        if (itpblktot == -1) then
75          ! The field 'total number of different boundary 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 boundary@xxx is not found, then
83            ! readheader will return intfromfile unchanged
84
85            intfromfile(:)=-1
86            iblk = iblk+1
87            write (temp1,"('connectivity boundary',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            !write(*,*) 'rank, fname2',myrank, trim(adjustl(fname2))
93            call readheader(igeom,fname2 // char(0),intfromfile,
94     &      ieight,'integer' // char(0),iotype)
95            neltp = intfromfile(1) ! -1 if fname2 was not found, >=0 otherwise
96          end do
97          itpblktot = iblk-1
98        end if
99
100        if (myrank == 0) then
101          write(*,*) 'Number of boundary topologies: ',itpblktot
102        endif
103!        write (*,*) 'Rank: ',myrank,' boundary itpblktot final:',
104!     &               itpblktot
105
106        nelblb=0
107        mattyp=0
108        ndofl = ndof
109        do iblk = 1, itpblktot
110           writeLock=0;
111
112           fname1='connectivity boundary?'
113
114!           print *, "Loop ",iblk, myrank, itpblk, trim(fnamer)
115
116ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
117
118           write (temp1,"('connectivity boundary',i1)") iblk
119           temp1 = trim(temp1)
120           write (temp3,"('(''@'',i',i1,',A1)')") itmp
121           write (fname2, temp3) (myrank+1), '?'
122           fname2 = trim(temp1)//trim(fname2)
123           fname2 = trim(fname2)
124
125ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
126
127           ! Synchronization for performance monitoring, as some parts do not include some topologies
128           call MPI_Barrier(MPI_COMM_WORLD,ierr)
129           call readheader(igeom,fname2 // char(0),intfromfile,ieight,
130     &                     'integer' // char(0),iotype)
131           neltp =intfromfile(1)
132           nenl  =intfromfile(2)
133           ipordl=intfromfile(3)
134           nshl  =intfromfile(4)
135           nshlb =intfromfile(5)
136           nenbl =intfromfile(6)
137           lcsyst=intfromfile(7)
138           numnbc=intfromfile(8)
139
140           allocate (ientp(neltp,nshl))
141           allocate (iBCBtp(neltp,ndiBCB))
142           allocate (BCBtp(neltp,ndBCB))
143           iientpsiz=neltp*nshl
144
145           if (neltp==0) then
146              writeLock=1;
147           endif
148
149!           print *, "neltp is ", neltp
150
151           call readdatablock(igeom,fname2 // char(0),ientp,iientpsiz,
152     &                     'integer' // char(0),iotype)
153
154
155c
156c.... Read the boundary flux codes
157c
158
159
160
161           fname1='nbc codes?'
162
163ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
164
165           write (temp1,"('nbc codes',i1)") iblk
166           temp1=trim(temp1)
167           write (temp3,"('(''@'',i',i1,',A1)')") itmp
168           write (fname2, temp3) (myrank+1), '?'
169           fname2 = trim(temp1)//trim(fname2)
170           call MPI_BARRIER(MPI_COMM_WORLD, ierr)
171
172ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
173
174           call readheader(igeom,fname2 // char(0) ,intfromfile,
175     &      ieight,'integer' // char(0),iotype)
176           iiBCBtpsiz=neltp*ndiBCB
177           call readdatablock(igeom,fname2 // char(0) ,iBCBtp,
178     &      iiBCBtpsiz,'integer' // char(0),iotype)
179
180c
181c.... read the boundary condition data
182c
183           fname1='nbc values?'
184
185ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
186
187           write (temp1,"('nbc values',i1)") iblk
188           temp1=trim(temp1)
189           write (temp3,"('(''@'',i',i1,',A1)')") itmp
190           write (fname2, temp3) (myrank+1), '?'
191           fname2 = trim(temp1)//trim(fname2)
192           call MPI_BARRIER(MPI_COMM_WORLD, ierr)
193ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
194
195           call readheader(igeom,fname2 // char(0) ,intfromfile,
196     &      ieight,'integer' // char(0) ,iotype)
197           BCBtp    = zero
198           iBCBtpsiz=neltp*ndBCB
199           call readdatablock(igeom,fname2 // char(0),BCBtp,iBCBtpsiz,
200     &                     'double' // char(0) ,iotype)
201
202
203c
204c This is a temporary fix until NSpre properly zeros this array where it
205c is not set.  DEC has indigestion with these arrays though the
206c result is never used (never effects solution).
207c
208
209
210           if(writeLock==0) then
211              where(.not.btest(iBCBtp(:,1),0)) BCBtp(:,1)=zero
212              where(.not.btest(iBCBtp(:,1),1)) BCBtp(:,2)=zero
213              where(.not.btest(iBCBtp(:,1),3)) BCBtp(:,6)=zero
214              if(ndBCB.gt.6) then
215                 do i=6,ndof
216                    where(.not.btest(iBCBtp(:,1),i-1)) BCBtp(:,i+1)=zero
217                 enddo
218              endif
219              where(.not.btest(iBCBtp(:,1),2))
220                 BCBtp(:,3)=zero
221                 BCBtp(:,4)=zero
222                 BCBtp(:,5)=zero
223              endwhere
224
225              do n=1,neltp,ibksz
226                 nelblb=nelblb+1
227                 npro= min(IBKSZ, neltp - n + 1)
228c
229                 lcblkb(1,nelblb)  = iel
230                 lcblkb(3,nelblb)  = lcsyst
231                 lcblkb(4,nelblb)  = ipordl
232                 lcblkb(5,nelblb)  = nenl
233                 lcblkb(6,nelblb)  = nenbl
234                 lcblkb(7,nelblb)  = mattyp
235                 lcblkb(8,nelblb)  = ndofl
236                 lcblkb(9,nelblb)  = nshl
237                 lcblkb(10,nelblb) = nshlb ! # of shape functions per elt
238c
239c.... save the element block
240c
241                 n1=n
242                 n2=n+npro-1
243                 materb=1       ! all one material for now
244c
245c.... allocate memory for stack arrays
246c
247
248                 allocate (mienb(nelblb)%p(npro,nshl))
249c
250                 allocate (miBCB(nelblb)%p(npro,ndiBCB))
251c
252                 allocate (mBCB(nelblb)%p(npro,nshlb,ndBCB))
253c
254                 allocate (mmatb(nelblb)%p(npro))
255c
256c.... save the boundary element block
257c
258                 call gensvb (ientp(n1:n2,1:nshl),
259     &                iBCBtp(n1:n2,:),      BCBtp(n1:n2,:),
260     &                materb,        mienb(nelblb)%p,
261     &                miBCB(nelblb)%p,        mBCB(nelblb)%p,
262     &                mmatb(nelblb)%p)
263c
264                 iel=iel+npro
265              enddo
266
267           endif
268           deallocate(ientp)
269           deallocate(iBCBtp)
270           deallocate(BCBtp)
271
272        enddo
273        lcblkb(1,nelblb+1) = iel
274
275c
276c.... return
277c
278        return
279c
280c.... end of file error handling
281c
282 911    call error ('genbcb  ','end file',igeomBAK)
283c
2841000    format(a80,//,
285     &  ' B o u n d a r y   E l e m e n t   C o n n e c t i v i t y',//,
286     &  '   Elem   BC codes',/,
287     &  '  Number  C P V H ',5x,27('Node',i1,:,2x))
2881100    format(2x,i5,2x,4i2,3x,27i7)
289c$$$2000    format(a80,//,
290c$$$     &  ' B o u n d a r y   E l e m e n t   B C   D a t a ',//,
291c$$$     &  '   Node   ',3x,'mass',/,
292c$$$     &  '  Number  ',3x,'flux',6x,'Pressure',6x,'Heat',6x,
293c$$$     &  3('Viscous',i1,:,4x))
2942100    format(2x,i5,1p,1x,6e12.4)
295c
296        end
297
298