xref: /phasta/phSolver/common/genbkb.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen        subroutine genbkb (ibksz)
2*59599516SKenneth E. Jansenc
3*59599516SKenneth E. Jansenc----------------------------------------------------------------------
4*59599516SKenneth E. Jansenc
5*59599516SKenneth E. Jansenc  This routine reads the boundary elements, reorders them and
6*59599516SKenneth E. Jansenc  generates traces for the gather/scatter operations.
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 dtnmod
12*59599516SKenneth E. Jansen        use pointer_data
13*59599516SKenneth E. Jansenc
14*59599516SKenneth E. Jansen        include "common.h"
15*59599516SKenneth E. Jansen!MR CHANGE
16*59599516SKenneth E. Jansen        include "mpif.h" !Required to determine the max for itpblk
17*59599516SKenneth E. Jansen!MR CHANGE END
18*59599516SKenneth E. Jansenc
19*59599516SKenneth E. Jansen
20*59599516SKenneth E. Jansen        integer, allocatable :: ientp(:,:),iBCBtp(:,:)
21*59599516SKenneth E. Jansen        real*8, allocatable :: BCBtp(:,:)
22*59599516SKenneth E. Jansen        integer materb(ibksz)
23*59599516SKenneth E. Jansen        integer intfromfile(50) ! integers read from headers
24*59599516SKenneth E. Jansen        character*255 fname1
25*59599516SKenneth E. Jansen
26*59599516SKenneth E. Jansencccccccccccccc New Phasta IO starts here cccccccccccccccccccccccccccccc
27*59599516SKenneth E. Jansen
28*59599516SKenneth E. Jansen        integer :: descriptor, descriptorG, GPID, color, nfiles, nfields
29*59599516SKenneth E. Jansen        integer :: numparts, nppf, nppp, nprocs, writeLock
30*59599516SKenneth E. Jansen        integer :: ierr_io, numprocs, itmp, itmp2
31*59599516SKenneth E. Jansen!        integer :: num_local_loop, num_global_loop
32*59599516SKenneth E. Jansen!MR CHANGE
33*59599516SKenneth E. Jansen        integer :: itpblktot,ierr
34*59599516SKenneth E. Jansen!MR CHANGE END
35*59599516SKenneth E. Jansen
36*59599516SKenneth E. Jansen        character*255 fnamer, fname2, temp2
37*59599516SKenneth E. Jansen        character*64 temp1, temp3
38*59599516SKenneth E. Jansen
39*59599516SKenneth E. Jansen        nfiles = nsynciofiles
40*59599516SKenneth E. Jansen        nfields = nsynciofieldsreadgeombc
41*59599516SKenneth E. Jansen        numparts = numpe !This is the common settings. Beware if you try to compute several parts per process
42*59599516SKenneth E. Jansen
43*59599516SKenneth E. Jansen        nppp = numparts/numpe
44*59599516SKenneth E. Jansen        nppf = numparts/nfiles
45*59599516SKenneth E. Jansen
46*59599516SKenneth E. Jansen        color = int(myrank/(numparts/nfiles))
47*59599516SKenneth E. Jansen        itmp2 = int(log10(float(color+1)))+1
48*59599516SKenneth E. Jansen        write (temp2,"('(''geombc-dat.'',i',i1,')')") itmp2
49*59599516SKenneth E. Jansen        temp2=trim(temp2)
50*59599516SKenneth E. Jansen        write (fnamer,temp2) (color+1)
51*59599516SKenneth E. Jansen        fnamer=trim(fnamer)
52*59599516SKenneth E. Jansen
53*59599516SKenneth E. Jansen        ione=1
54*59599516SKenneth E. Jansen        itwo=2
55*59599516SKenneth E. Jansen        ieight=8
56*59599516SKenneth E. Jansen        ieleven=11
57*59599516SKenneth E. Jansen        itmp = int(log10(float(myrank+1)))+1
58*59599516SKenneth E. Jansen
59*59599516SKenneth E. Jansen!        num_global_loop = 4
60*59599516SKenneth E. Jansen!        num_local_loop = nelblb
61*59599516SKenneth E. Jansen
62*59599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
63*59599516SKenneth E. Jansen
64*59599516SKenneth E. Jansen        iel=1
65*59599516SKenneth E. Jansen        itpblk=nelblb
66*59599516SKenneth E. Jansen
67*59599516SKenneth E. Jansen
68*59599516SKenneth E. Jansen!MR CHANGE
69*59599516SKenneth E. Jansen        ! Get the total number of different interior topologies in the whole domain.
70*59599516SKenneth E. Jansen        ! Try to read from a field. If the field does not exist, scan the geombc file.
71*59599516SKenneth E. Jansen        itpblktot=-1
72*59599516SKenneth E. Jansen        write(temp1,
73*59599516SKenneth E. Jansen     &   "('(''total number of boundary tpblocks@'',i',i1,',A1)')") itmp
74*59599516SKenneth E. Jansen        write (fname2,temp1) (myrank+1),'?'
75*59599516SKenneth E. Jansen        call readheader(igeom,fname2 // char(0) ,itpblktot,ione,
76*59599516SKenneth E. Jansen     &  'integer' // char(0),iotype)
77*59599516SKenneth E. Jansen
78*59599516SKenneth E. Jansen!        write (*,*) 'Rank: ',myrank,' boundary itpblktot intermediate:',
79*59599516SKenneth E. Jansen!     &               itpblktot
80*59599516SKenneth E. Jansen
81*59599516SKenneth E. Jansen        if (itpblktot == -1) then
82*59599516SKenneth E. Jansen          ! The field 'total number of different boundary tpblocks' was not found in the geombc file.
83*59599516SKenneth E. Jansen          ! Scan all the geombc file for the 'connectivity interior' fields to get this information.
84*59599516SKenneth E. Jansen          iblk=0
85*59599516SKenneth E. Jansen          neltp=0
86*59599516SKenneth E. Jansen          do while(neltp .ne. -1)
87*59599516SKenneth E. Jansen
88*59599516SKenneth E. Jansen            ! intfromfile is reinitialized to -1 every time.
89*59599516SKenneth E. Jansen            ! If connectivity boundary@xxx is not found, then
90*59599516SKenneth E. Jansen            ! readheader will return intfromfile unchanged
91*59599516SKenneth E. Jansen
92*59599516SKenneth E. Jansen            intfromfile(:)=-1
93*59599516SKenneth E. Jansen            iblk = iblk+1
94*59599516SKenneth E. Jansen            write (temp1,"('connectivity boundary',i1)") iblk
95*59599516SKenneth E. Jansen            temp1 = trim(temp1)
96*59599516SKenneth E. Jansen            write (temp3,"('(''@'',i',i1,',A1)')") itmp
97*59599516SKenneth E. Jansen            write (fname2, temp3) (myrank+1), '?'
98*59599516SKenneth E. Jansen            fname2 = trim(temp1)//trim(fname2)
99*59599516SKenneth E. Jansen            !write(*,*) 'rank, fname2',myrank, trim(adjustl(fname2))
100*59599516SKenneth E. Jansen            call readheader(igeom,fname2 // char(0),intfromfile,
101*59599516SKenneth E. Jansen     &      ieight,'integer' // char(0),iotype)
102*59599516SKenneth E. Jansen            neltp = intfromfile(1) ! -1 if fname2 was not found, >=0 otherwise
103*59599516SKenneth E. Jansen          end do
104*59599516SKenneth E. Jansen          itpblktot = iblk-1
105*59599516SKenneth E. Jansen        end if
106*59599516SKenneth E. Jansen
107*59599516SKenneth E. Jansen        if (myrank == 0) then
108*59599516SKenneth E. Jansen          write(*,*) 'Number of boundary topologies: ',itpblktot
109*59599516SKenneth E. Jansen        endif
110*59599516SKenneth E. Jansen!        write (*,*) 'Rank: ',myrank,' boundary itpblktot final:',
111*59599516SKenneth E. Jansen!     &               itpblktot
112*59599516SKenneth E. Jansen
113*59599516SKenneth E. Jansen!MR CHANGE END
114*59599516SKenneth E. Jansen        nelblb=0
115*59599516SKenneth E. Jansen        mattyp=0
116*59599516SKenneth E. Jansen        ndofl = ndof
117*59599516SKenneth E. Jansen!MR CHANGE
118*59599516SKenneth E. Jansen!        call initphmpiio( nfields, nppf, nfiles, igeom )
119*59599516SKenneth E. Jansen!        call openfile( fnamer, 'read', igeom )
120*59599516SKenneth E. Jansen
121*59599516SKenneth E. Jansen        do iblk = 1, itpblktot
122*59599516SKenneth E. Jansen           writeLock=0;
123*59599516SKenneth E. Jansen!MR CHANGE END
124*59599516SKenneth E. Jansen
125*59599516SKenneth E. Jansen           fname1='connectivity boundary?'
126*59599516SKenneth E. Jansen
127*59599516SKenneth E. Jansen!           print *, "Loop ",iblk, myrank, itpblk, trim(fnamer)
128*59599516SKenneth E. Jansen
129*59599516SKenneth E. Jansenccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
130*59599516SKenneth E. Jansen
131*59599516SKenneth E. Jansen           write (temp1,"('connectivity boundary',i1)") iblk
132*59599516SKenneth E. Jansen           temp1 = trim(temp1)
133*59599516SKenneth E. Jansen           write (temp3,"('(''@'',i',i1,',A1)')") itmp
134*59599516SKenneth E. Jansen           write (fname2, temp3) (myrank+1), '?'
135*59599516SKenneth E. Jansen           fname2 = trim(temp1)//trim(fname2)
136*59599516SKenneth E. Jansen           fname2 = trim(fname2)
137*59599516SKenneth E. Jansen!           write(*,*) 'rank, fname2',myrank, trim(adjustl(fname2))
138*59599516SKenneth E. Jansen
139*59599516SKenneth E. Jansenccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
140*59599516SKenneth E. Jansen
141*59599516SKenneth E. Jansen           ! Synchronization for performance monitoring, as some parts do not include some topologies
142*59599516SKenneth E. Jansen           call MPI_Barrier(MPI_COMM_WORLD,ierr)
143*59599516SKenneth E. Jansen           call readheader(igeom,fname2 // char(0),intfromfile,ieight,
144*59599516SKenneth E. Jansen     &                     'integer' // char(0),iotype)
145*59599516SKenneth E. Jansen           neltp =intfromfile(1)
146*59599516SKenneth E. Jansen           nenl  =intfromfile(2)
147*59599516SKenneth E. Jansen           ipordl=intfromfile(3)
148*59599516SKenneth E. Jansen           nshl  =intfromfile(4)
149*59599516SKenneth E. Jansen           nshlb =intfromfile(5)
150*59599516SKenneth E. Jansen           nenbl =intfromfile(6)
151*59599516SKenneth E. Jansen           lcsyst=intfromfile(7)
152*59599516SKenneth E. Jansen           numnbc=intfromfile(8)
153*59599516SKenneth E. Jansen
154*59599516SKenneth E. Jansen!MR CHANGE
155*59599516SKenneth E. Jansen!           write (temp1,"('connectivityBoundaryHeader_',i1,'_',i1)")
156*59599516SKenneth E. Jansen!     &                                              iblk,myrank
157*59599516SKenneth E. Jansen!           temp1=trim(temp1)
158*59599516SKenneth E. Jansen!           open(unit=14,file=temp1,status='unknown')
159*59599516SKenneth E. Jansen!           write(14,*) intfromfile(:)
160*59599516SKenneth E. Jansen!           close(14)
161*59599516SKenneth E. Jansen!MR CHANGE END
162*59599516SKenneth E. Jansen
163*59599516SKenneth E. Jansen           allocate (ientp(neltp,nshl))
164*59599516SKenneth E. Jansen           allocate (iBCBtp(neltp,ndiBCB))
165*59599516SKenneth E. Jansen           allocate (BCBtp(neltp,ndBCB))
166*59599516SKenneth E. Jansen           iientpsiz=neltp*nshl
167*59599516SKenneth E. Jansen
168*59599516SKenneth E. Jansen           if (neltp==0) then
169*59599516SKenneth E. Jansen              writeLock=1;
170*59599516SKenneth E. Jansen           endif
171*59599516SKenneth E. Jansen
172*59599516SKenneth E. Jansen!           print *, "neltp is ", neltp
173*59599516SKenneth E. Jansen
174*59599516SKenneth E. Jansen           call readdatablock(igeom,fname2 // char(0),ientp,iientpsiz,
175*59599516SKenneth E. Jansen     &                     'integer' // char(0),iotype)
176*59599516SKenneth E. Jansen
177*59599516SKenneth E. Jansen!MR CHANGE
178*59599516SKenneth E. Jansen!           write (temp1,"('connectivityBoundaryDatablock_',i1,'_',i1)")
179*59599516SKenneth E. Jansen!     &                                              iblk,myrank
180*59599516SKenneth E. Jansen!           temp1=trim(temp1)
181*59599516SKenneth E. Jansen
182*59599516SKenneth E. Jansen!           open(unit=14,file=temp1,status='unknown')
183*59599516SKenneth E. Jansen!           do i=1,neltp
184*59599516SKenneth E. Jansen!              do j=1,nshl
185*59599516SKenneth E. Jansen!                write(14,*) ientp(i,j)
186*59599516SKenneth E. Jansen!              enddo
187*59599516SKenneth E. Jansen!           enddo
188*59599516SKenneth E. Jansen!           write(14,*) iientpsiz
189*59599516SKenneth E. Jansen!           close(14)
190*59599516SKenneth E. Jansen!MR CHANGE END
191*59599516SKenneth E. Jansen
192*59599516SKenneth E. Jansen
193*59599516SKenneth E. Jansenc
194*59599516SKenneth E. Jansenc.... Read the boundary flux codes
195*59599516SKenneth E. Jansenc
196*59599516SKenneth E. Jansen
197*59599516SKenneth E. Jansen!           print *,"connectivity [] is ", trim(fname2),ientp(0,0)
198*59599516SKenneth E. Jansen
199*59599516SKenneth E. Jansen
200*59599516SKenneth E. Jansen           fname1='nbc codes?'
201*59599516SKenneth E. Jansen
202*59599516SKenneth E. Jansenccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
203*59599516SKenneth E. Jansen
204*59599516SKenneth E. Jansen           write (temp1,"('nbc codes',i1)") iblk
205*59599516SKenneth E. Jansen           temp1=trim(temp1)
206*59599516SKenneth E. Jansen           write (temp3,"('(''@'',i',i1,',A1)')") itmp
207*59599516SKenneth E. Jansen           write (fname2, temp3) (myrank+1), '?'
208*59599516SKenneth E. Jansen           fname2 = trim(temp1)//trim(fname2)
209*59599516SKenneth E. Jansen!           write(*,*) 'rank, fname2',myrank, trim(adjustl(fname2))
210*59599516SKenneth E. Jansen           call MPI_BARRIER(MPI_COMM_WORLD, ierr)
211*59599516SKenneth E. Jansen
212*59599516SKenneth E. Jansenccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
213*59599516SKenneth E. Jansen
214*59599516SKenneth E. Jansen           call readheader(igeom,fname2 // char(0) ,intfromfile,
215*59599516SKenneth E. Jansen     &      ieight,'integer' // char(0),iotype)
216*59599516SKenneth E. Jansen           iiBCBtpsiz=neltp*ndiBCB
217*59599516SKenneth E. Jansen           call readdatablock(igeom,fname2 // char(0) ,iBCBtp,
218*59599516SKenneth E. Jansen     &      iiBCBtpsiz,'integer' // char(0),iotype)
219*59599516SKenneth E. Jansen
220*59599516SKenneth E. Jansen!MR CHANGE
221*59599516SKenneth E. Jansen!           print *, "ndiBCB is ",ndiBCB
222*59599516SKenneth E. Jansen!           write (temp1,"('nbcCodesDatablock_',i1,'_',i1)")
223*59599516SKenneth E. Jansen!     &                                              iblk,myrank
224*59599516SKenneth E. Jansen!           temp1=trim(temp1)
225*59599516SKenneth E. Jansen!
226*59599516SKenneth E. Jansen!           open(unit=13,file=temp1,status='unknown')
227*59599516SKenneth E. Jansen!           do i=1,neltp
228*59599516SKenneth E. Jansen!              do j=1,ndiBCB
229*59599516SKenneth E. Jansen!                write(13,*) iBCBtp(i,j)
230*59599516SKenneth E. Jansen!              enddo
231*59599516SKenneth E. Jansen!           enddo
232*59599516SKenneth E. Jansen!           write(13,*) iiBCBtpsiz
233*59599516SKenneth E. Jansen!           close(13)
234*59599516SKenneth E. Jansen!!           iBCBtp(:,:) = 0 ! JUST FOR TEST
235*59599516SKenneth E. Jansen!MR CHANGE END
236*59599516SKenneth E. Jansen
237*59599516SKenneth E. Jansenc
238*59599516SKenneth E. Jansenc.... read the boundary condition data
239*59599516SKenneth E. Jansenc
240*59599516SKenneth E. Jansen           fname1='nbc values?'
241*59599516SKenneth E. Jansen
242*59599516SKenneth E. Jansenccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
243*59599516SKenneth E. Jansen
244*59599516SKenneth E. Jansen           write (temp1,"('nbc values',i1)") iblk
245*59599516SKenneth E. Jansen           temp1=trim(temp1)
246*59599516SKenneth E. Jansen           write (temp3,"('(''@'',i',i1,',A1)')") itmp
247*59599516SKenneth E. Jansen           write (fname2, temp3) (myrank+1), '?'
248*59599516SKenneth E. Jansen           fname2 = trim(temp1)//trim(fname2)
249*59599516SKenneth E. Jansen           call MPI_BARRIER(MPI_COMM_WORLD, ierr)
250*59599516SKenneth E. Jansenccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
251*59599516SKenneth E. Jansen
252*59599516SKenneth E. Jansen           call readheader(igeom,fname2 // char(0) ,intfromfile,
253*59599516SKenneth E. Jansen     &      ieight,'integer' // char(0) ,iotype)
254*59599516SKenneth E. Jansen           BCBtp    = zero
255*59599516SKenneth E. Jansen           iBCBtpsiz=neltp*ndBCB
256*59599516SKenneth E. Jansen           call readdatablock(igeom,fname2 // char(0),BCBtp,iBCBtpsiz,
257*59599516SKenneth E. Jansen     &                     'double' // char(0) ,iotype)
258*59599516SKenneth E. Jansen
259*59599516SKenneth E. Jansen!MR CHANGE
260*59599516SKenneth E. Jansen!           write (temp1,"('nbcValuesDatablock_',i1,'_',i1)")
261*59599516SKenneth E. Jansen!     &                                              iblk,myrank
262*59599516SKenneth E. Jansen!           temp1=trim(temp1)
263*59599516SKenneth E. Jansen!           open(unit=12,file=temp1,status='unknown')
264*59599516SKenneth E. Jansen!           do i=1,neltp
265*59599516SKenneth E. Jansen!              do j=1,ndBCB
266*59599516SKenneth E. Jansen!                write(12,*) BCBtp(i,j)
267*59599516SKenneth E. Jansen!              enddo
268*59599516SKenneth E. Jansen!           enddo
269*59599516SKenneth E. Jansen!           write(12,*) iBCBtpsiz
270*59599516SKenneth E. Jansen!           close(12)
271*59599516SKenneth E. Jansen!MR CHANGE END
272*59599516SKenneth E. Jansen
273*59599516SKenneth E. Jansenc
274*59599516SKenneth E. Jansenc This is a temporary fix until NSpre properly zeros this array where it
275*59599516SKenneth E. Jansenc is not set.  DEC has indigestion with these arrays though the
276*59599516SKenneth E. Jansenc result is never used (never effects solution).
277*59599516SKenneth E. Jansenc
278*59599516SKenneth E. Jansen
279*59599516SKenneth E. Jansen
280*59599516SKenneth E. Jansen!MR CHANGE
281*59599516SKenneth E. Jansen           if(writeLock==0) then
282*59599516SKenneth E. Jansen!MR CHANGE
283*59599516SKenneth E. Jansen!              print *,"In ASSIGN ASSIGN",myrank
284*59599516SKenneth E. Jansen              where(.not.btest(iBCBtp(:,1),0)) BCBtp(:,1)=zero
285*59599516SKenneth E. Jansen              where(.not.btest(iBCBtp(:,1),1)) BCBtp(:,2)=zero
286*59599516SKenneth E. Jansen              where(.not.btest(iBCBtp(:,1),3)) BCBtp(:,6)=zero
287*59599516SKenneth E. Jansen              if(ndBCB.gt.6) then
288*59599516SKenneth E. Jansen                 do i=6,ndof
289*59599516SKenneth E. Jansen                    where(.not.btest(iBCBtp(:,1),i-1)) BCBtp(:,i+1)=zero
290*59599516SKenneth E. Jansen                 enddo
291*59599516SKenneth E. Jansen              endif
292*59599516SKenneth E. Jansen              where(.not.btest(iBCBtp(:,1),2))
293*59599516SKenneth E. Jansen                 BCBtp(:,3)=zero
294*59599516SKenneth E. Jansen                 BCBtp(:,4)=zero
295*59599516SKenneth E. Jansen                 BCBtp(:,5)=zero
296*59599516SKenneth E. Jansen              endwhere
297*59599516SKenneth E. Jansen
298*59599516SKenneth E. Jansen              do n=1,neltp,ibksz
299*59599516SKenneth E. Jansen                 nelblb=nelblb+1
300*59599516SKenneth E. Jansen                 npro= min(IBKSZ, neltp - n + 1)
301*59599516SKenneth E. Jansenc
302*59599516SKenneth E. Jansen                 lcblkb(1,nelblb)  = iel
303*59599516SKenneth E. Jansenc     lcblkb(2,nelblb)  = iopen ! available for later use
304*59599516SKenneth E. Jansen                 lcblkb(3,nelblb)  = lcsyst
305*59599516SKenneth E. Jansen                 lcblkb(4,nelblb)  = ipordl
306*59599516SKenneth E. Jansen                 lcblkb(5,nelblb)  = nenl
307*59599516SKenneth E. Jansen                 lcblkb(6,nelblb)  = nenbl
308*59599516SKenneth E. Jansen                 lcblkb(7,nelblb)  = mattyp
309*59599516SKenneth E. Jansen                 lcblkb(8,nelblb)  = ndofl
310*59599516SKenneth E. Jansen                 lcblkb(9,nelblb)  = nshl
311*59599516SKenneth E. Jansen                 lcblkb(10,nelblb) = nshlb ! # of shape functions per elt
312*59599516SKenneth E. Jansenc
313*59599516SKenneth E. Jansenc.... save the element block
314*59599516SKenneth E. Jansenc
315*59599516SKenneth E. Jansen                 n1=n
316*59599516SKenneth E. Jansen                 n2=n+npro-1
317*59599516SKenneth E. Jansen                 materb=1       ! all one material for now
318*59599516SKenneth E. Jansenc
319*59599516SKenneth E. Jansenc.... allocate memory for stack arrays
320*59599516SKenneth E. Jansenc
321*59599516SKenneth E. Jansen
322*59599516SKenneth E. Jansen                 allocate (mienb(nelblb)%p(npro,nshl))
323*59599516SKenneth E. Jansenc
324*59599516SKenneth E. Jansen                 allocate (miBCB(nelblb)%p(npro,ndiBCB))
325*59599516SKenneth E. Jansenc
326*59599516SKenneth E. Jansen                 allocate (mBCB(nelblb)%p(npro,nshlb,ndBCB))
327*59599516SKenneth E. Jansenc
328*59599516SKenneth E. Jansen                 allocate (mmatb(nelblb)%p(npro))
329*59599516SKenneth E. Jansenc
330*59599516SKenneth E. Jansenc.... save the boundary element block
331*59599516SKenneth E. Jansenc
332*59599516SKenneth E. Jansen                 call gensvb (ientp(n1:n2,1:nshl),
333*59599516SKenneth E. Jansen     &                iBCBtp(n1:n2,:),      BCBtp(n1:n2,:),
334*59599516SKenneth E. Jansen     &                materb,        mienb(nelblb)%p,
335*59599516SKenneth E. Jansen     &                miBCB(nelblb)%p,        mBCB(nelblb)%p,
336*59599516SKenneth E. Jansen     &                mmatb(nelblb)%p)
337*59599516SKenneth E. Jansenc
338*59599516SKenneth E. Jansen                 iel=iel+npro
339*59599516SKenneth E. Jansen              enddo
340*59599516SKenneth E. Jansen
341*59599516SKenneth E. Jansen!MR CHANGE
342*59599516SKenneth E. Jansen           endif
343*59599516SKenneth E. Jansen!MR CHANGE
344*59599516SKenneth E. Jansen           deallocate(ientp)
345*59599516SKenneth E. Jansen           deallocate(iBCBtp)
346*59599516SKenneth E. Jansen           deallocate(BCBtp)
347*59599516SKenneth E. Jansen
348*59599516SKenneth E. Jansen        enddo
349*59599516SKenneth E. Jansen        lcblkb(1,nelblb+1) = iel
350*59599516SKenneth E. Jansen
351*59599516SKenneth E. Jansen!        call closefile( igeom, "read" )
352*59599516SKenneth E. Jansen!        call finalizephmpiio( igeom )
353*59599516SKenneth E. Jansen
354*59599516SKenneth E. Jansenc
355*59599516SKenneth E. Jansenc.... return
356*59599516SKenneth E. Jansenc
357*59599516SKenneth E. Jansen        return
358*59599516SKenneth E. Jansenc
359*59599516SKenneth E. Jansenc.... end of file error handling
360*59599516SKenneth E. Jansenc
361*59599516SKenneth E. Jansen 911    call error ('genbcb  ','end file',igeomBAK)
362*59599516SKenneth E. Jansenc
363*59599516SKenneth E. Jansen1000    format(a80,//,
364*59599516SKenneth E. Jansen     &  ' 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',//,
365*59599516SKenneth E. Jansen     &  '   Elem   BC codes',/,
366*59599516SKenneth E. Jansen     &  '  Number  C P V H ',5x,27('Node',i1,:,2x))
367*59599516SKenneth E. Jansen1100    format(2x,i5,2x,4i2,3x,27i7)
368*59599516SKenneth E. Jansenc$$$2000    format(a80,//,
369*59599516SKenneth E. Jansenc$$$     &  ' B o u n d a r y   E l e m e n t   B C   D a t a ',//,
370*59599516SKenneth E. Jansenc$$$     &  '   Node   ',3x,'mass',/,
371*59599516SKenneth E. Jansenc$$$     &  '  Number  ',3x,'flux',6x,'Pressure',6x,'Heat',6x,
372*59599516SKenneth E. Jansenc$$$     &  3('Viscous',i1,:,4x))
373*59599516SKenneth E. Jansen2100    format(2x,i5,1p,1x,6e12.4)
374*59599516SKenneth E. Jansenc
375*59599516SKenneth E. Jansen        end
376*59599516SKenneth E. Jansen
377