xref: /phasta/phSolver/common/genbkbSyncIO.f (revision ab5b07a458969128644e941ce150fa049270fdf1)
1        subroutine genbkbSyncIO (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
13        use phio
14        use iso_c_binding
15        include "common.h"
16        include "mpif.h" !Required to determine the max for itpblk
17
18        integer, target, allocatable :: ientp(:,:),iBCBtp(:,:)
19        real*8, target, allocatable :: BCBtp(:,:)
20        integer materb(ibksz)
21        integer, target :: intfromfile(50) ! integers read from headers
22        character*255 fname1
23        integer :: descriptor, descriptorG, GPID, color, nfields
24        integer :: numparts, nppp, nprocs, writeLock
25        integer :: ierr_io, numprocs, itmp, itmp2
26        integer, target :: itpblktot,ierr
27        character*255 fname2
28        character(len=30) :: dataInt, dataDbl
29        dataInt = c_char_'integer'//c_null_char
30        dataDbl = c_char_'double'//c_null_char
31
32        nfields = nsynciofieldsreadgeombc
33        numparts = numpe !This is the common settings. Beware if you try to compute several parts per process
34        nppp = numparts/numpe
35        ione=1
36        itwo=2
37        ieight=8
38        ieleven=11
39        itmp = int(log10(float(myrank+1)))+1
40        iel=1
41        itpblk=nelblb
42
43        ! Get the total number of different interior topologies in the whole domain.
44        ! Try to read from a field. If the field does not exist, scan the geombc file.
45        itpblktot=-1
46        call phio_readheader(fhandle,
47     &   c_char_'total number of boundary tpblocks' // char(0),
48     &   c_loc(itpblktot), ione, dataInt, iotype)
49
50        if (itpblktot == -1) then
51          ! The field 'total number of different boundary tpblocks' was not found in the geombc file.
52          ! Scan all the geombc file for the 'connectivity interior' fields to get this information.
53          iblk=0
54          neltp=0
55          do while(neltp .ne. -1)
56
57            ! intfromfile is reinitialized to -1 every time.
58            ! If connectivity boundary@xxx is not found, then
59            ! readheader will return intfromfile unchanged
60
61            intfromfile(:)=-1
62            iblk = iblk+1
63            write (fname2,"('connectivity boundary',i1)") iblk
64
65
66            call phio_readheader(fhandle, fname2 // char(0),
67     &       c_loc(intfromfile), ieight, dataInt, iotype)
68            neltp = intfromfile(1) ! -1 if fname2 was not found, >=0 otherwise
69          end do
70          itpblktot = iblk-1
71        end if
72
73        if (myrank == 0) then
74          write(*,*) 'Number of boundary topologies: ',itpblktot
75        endif
76!        write (*,*) 'Rank: ',myrank,' boundary itpblktot final:',
77!     &               itpblktot
78
79        nelblb=0
80        mattyp=0
81        ndofl = ndof
82
83        do iblk = 1, itpblktot
84           writeLock=0;
85           write (fname2,"('connectivity boundary',i1)") iblk
86           ! Synchronization for performance monitoring, as some parts do not include some topologies
87           call MPI_Barrier(MPI_COMM_WORLD,ierr)
88           call phio_readheader(fhandle, fname2 // char(0),
89     &      c_loc(intfromfile), ieight, dataInt, iotype)
90           neltp =intfromfile(1)
91           nenl  =intfromfile(2)
92           ipordl=intfromfile(3)
93           nshl  =intfromfile(4)
94           nshlb =intfromfile(5)
95           nenbl =intfromfile(6)
96           lcsyst=intfromfile(7)
97           numnbc=intfromfile(8)
98
99           allocate (ientp(neltp,nshl))
100           allocate (iBCBtp(neltp,ndiBCB))
101           allocate (BCBtp(neltp,ndBCB))
102           iientpsiz=neltp*nshl
103
104           if (neltp==0) then
105              writeLock=1;
106           endif
107
108           call phio_readdatablock(fhandle, fname2 // char(0),
109     &      c_loc(ientp),iientpsiz,dataInt,iotype)
110c
111c.... Read the boundary flux codes
112c
113           call MPI_BARRIER(MPI_COMM_WORLD, ierr)
114           write (fname2,"('nbc codes',i1)") iblk
115
116           call phio_readheader(fhandle, fname2 // char(0),
117     &      c_loc(intfromfile), ieight, dataInt, iotype)
118           iiBCBtpsiz=neltp*ndiBCB
119           call phio_readdatablock(fhandle, fname2 // char(0),
120     &      c_loc(iBCBtp),iiBCBtpsiz,dataInt,iotype)
121c
122c.... read the boundary condition data
123c
124           call MPI_BARRIER(MPI_COMM_WORLD, ierr)
125           write (fname2,"('nbc values',i1)") iblk
126
127           call phio_readheader(fhandle, fname2 // char(0),
128     &      c_loc(intfromfile), ieight, dataInt, iotype)
129           BCBtp    = zero
130           iBCBtpsiz=neltp*ndBCB
131           call phio_readdatablock(fhandle, fname2 // char(0),
132     &      c_loc(BCBtp),iBCBtpsiz,dataDbl,iotype)
133c
134c This is a temporary fix until NSpre properly zeros this array where it
135c is not set.  DEC has indigestion with these arrays though the
136c result is never used (never effects solution).
137c
138           if(writeLock==0) then
139              where(.not.btest(iBCBtp(:,1),0)) BCBtp(:,1)=zero
140              where(.not.btest(iBCBtp(:,1),1)) BCBtp(:,2)=zero
141              where(.not.btest(iBCBtp(:,1),3)) BCBtp(:,6)=zero
142              if(ndBCB.gt.6) then
143                 do i=6,ndof
144                    where(.not.btest(iBCBtp(:,1),i-1)) BCBtp(:,i+1)=zero
145                 enddo
146              endif
147              where(.not.btest(iBCBtp(:,1),2))
148                 BCBtp(:,3)=zero
149                 BCBtp(:,4)=zero
150                 BCBtp(:,5)=zero
151              endwhere
152
153              do n=1,neltp,ibksz
154                 nelblb=nelblb+1
155                 npro= min(IBKSZ, neltp - n + 1)
156                 lcblkb(1,nelblb)  = iel
157                 lcblkb(3,nelblb)  = lcsyst
158                 lcblkb(4,nelblb)  = ipordl
159                 lcblkb(5,nelblb)  = nenl
160                 lcblkb(6,nelblb)  = nenbl
161                 lcblkb(7,nelblb)  = mattyp
162                 lcblkb(8,nelblb)  = ndofl
163                 lcblkb(9,nelblb)  = nshl
164                 lcblkb(10,nelblb) = nshlb ! # of shape functions per elt
165c
166c.... save the element block
167c
168                 n1=n
169                 n2=n+npro-1
170                 materb=1       ! all one material for now
171c
172c.... allocate memory for stack arrays
173c
174                 allocate (mienb(nelblb)%p(npro,nshl))
175                 allocate (miBCB(nelblb)%p(npro,ndiBCB))
176                 allocate (mBCB(nelblb)%p(npro,nshlb,ndBCB))
177                 allocate (mmatb(nelblb)%p(npro))
178c
179c.... save the boundary element block
180c
181                 call gensvb (ientp(n1:n2,1:nshl),
182     &                iBCBtp(n1:n2,:),      BCBtp(n1:n2,:),
183     &                materb,        mienb(nelblb)%p,
184     &                miBCB(nelblb)%p,        mBCB(nelblb)%p,
185     &                mmatb(nelblb)%p)
186                 iel=iel+npro
187              enddo
188           endif
189           deallocate(ientp)
190           deallocate(iBCBtp)
191           deallocate(BCBtp)
192
193        enddo
194        lcblkb(1,nelblb+1) = iel
195
196c
197c.... return
198c
199        return
200c
201c.... end of file error handling
202c
203 911    call error ('genbcb  ','end file',igeomBAK)
204c
2051000    format(a80,//,
206     &  ' 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',//,
207     &  '   Elem   BC codes',/,
208     &  '  Number  C P V H ',5x,27('Node',i1,:,2x))
2091100    format(2x,i5,2x,4i2,3x,27i7)
210c$$$2000    format(a80,//,
211c$$$     &  ' B o u n d a r y   E l e m e n t   B C   D a t a ',//,
212c$$$     &  '   Node   ',3x,'mass',/,
213c$$$     &  '  Number  ',3x,'flux',6x,'Pressure',6x,'Heat',6x,
214c$$$     &  3('Viscous',i1,:,4x))
2152100    format(2x,i5,1p,1x,6e12.4)
216c
217        end
218
219