xref: /phasta/phSolver/common/genbkb.f (revision 7acde132a6def0fe2daaec0d1a712dff0e5c6636)
159599516SKenneth E. Jansen        subroutine genbkb (ibksz)
259599516SKenneth E. Jansenc
359599516SKenneth E. Jansenc----------------------------------------------------------------------
459599516SKenneth E. Jansenc
559599516SKenneth E. Jansenc  This routine reads the boundary elements, reorders them and
659599516SKenneth E. Jansenc  generates traces for the gather/scatter operations.
759599516SKenneth E. Jansenc
859599516SKenneth E. Jansenc Zdenek Johan, Fall 1991.
959599516SKenneth E. Jansenc----------------------------------------------------------------------
1059599516SKenneth E. Jansenc
1159599516SKenneth E. Jansen        use dtnmod
1259599516SKenneth E. Jansen        use pointer_data
13e5afe575SCameron Smith        use phio
14e5afe575SCameron Smith        use iso_c_binding
1559599516SKenneth E. Jansen        include "common.h"
1659599516SKenneth E. Jansen        include "mpif.h" !Required to determine the max for itpblk
1759599516SKenneth E. Jansen
189a6935e5SKenneth E. Jansen        integer, target, allocatable :: ientp(:,:),iBCBtp(:,:)
199a6935e5SKenneth E. Jansen        real*8, target, allocatable :: BCBtp(:,:)
2059599516SKenneth E. Jansen        integer materb(ibksz)
219a6935e5SKenneth E. Jansen        integer, target :: intfromfile(50) ! integers read from headers
2259599516SKenneth E. Jansen        character*255 fname1
23fcf561c1SCameron Smith        integer :: descriptor, descriptorG, GPID, color, nfields
242efdc748SKenneth E. Jansen        integer :: numparts, nppp, nprocs, writeLock
2559599516SKenneth E. Jansen        integer :: ierr_io, numprocs, itmp, itmp2
269a6935e5SKenneth E. Jansen        integer, target :: itpblktot,ierr
272efdc748SKenneth E. Jansen        character*255 fname2
28bc62cfd4SCameron Smith        character(len=30) :: dataInt, dataDbl
29e5afe575SCameron Smith        dataInt = c_char_'integer'//c_null_char
30bc62cfd4SCameron Smith        dataDbl = c_char_'double'//c_null_char
31e5afe575SCameron Smith
3259599516SKenneth E. Jansen        nfields = nsynciofieldsreadgeombc
3359599516SKenneth E. Jansen        numparts = numpe !This is the common settings. Beware if you try to compute several parts per process
3459599516SKenneth E. Jansen        nppp = numparts/numpe
3559599516SKenneth E. Jansen        ione=1
3659599516SKenneth E. Jansen        itwo=2
3759599516SKenneth E. Jansen        ieight=8
3859599516SKenneth E. Jansen        ieleven=11
3959599516SKenneth E. Jansen        itmp = int(log10(float(myrank+1)))+1
4059599516SKenneth E. Jansen        iel=1
4159599516SKenneth E. Jansen        itpblk=nelblb
4259599516SKenneth E. Jansen
4359599516SKenneth E. Jansen        ! Get the total number of different interior topologies in the whole domain.
4459599516SKenneth E. Jansen        ! Try to read from a field. If the field does not exist, scan the geombc file.
452efdc748SKenneth E. Jansen        itpblktot=1  ! hardwired to monotopology for now
46d5d2f64dSCameron Smith        call phio_readheader(fhandle,
47e5afe575SCameron Smith     &   c_char_'total number of boundary tpblocks' // char(0),
48e5afe575SCameron Smith     &   c_loc(itpblktot), ione, dataInt, iotype)
4959599516SKenneth E. Jansen
5059599516SKenneth E. Jansen        if (itpblktot == -1) then
5159599516SKenneth E. Jansen          ! The field 'total number of different boundary tpblocks' was not found in the geombc file.
5259599516SKenneth E. Jansen          ! Scan all the geombc file for the 'connectivity interior' fields to get this information.
5359599516SKenneth E. Jansen          iblk=0
5459599516SKenneth E. Jansen          neltp=0
5559599516SKenneth E. Jansen          do while(neltp .ne. -1)
5659599516SKenneth E. Jansen
5759599516SKenneth E. Jansen            ! intfromfile is reinitialized to -1 every time.
5859599516SKenneth E. Jansen            ! If connectivity boundary@xxx is not found, then
5959599516SKenneth E. Jansen            ! readheader will return intfromfile unchanged
6059599516SKenneth E. Jansen
6159599516SKenneth E. Jansen            intfromfile(:)=-1
6259599516SKenneth E. Jansen            iblk = iblk+1
63*5b7f36ccSCameron Smith            if(input_mode.ge.1)then
642efdc748SKenneth E. Jansen               write (fname2,"('connectivity boundary',i1)") iblk
652efdc748SKenneth E. Jansen            else
662efdc748SKenneth E. Jansen               write (fname2,"('connectivity boundary linear tetrahedron')")
672efdc748SKenneth E. Jansen            endif
682efdc748SKenneth E. Jansen
692efdc748SKenneth E. Jansen            call phio_readheader(fhandle, fname2 // char(0),
70e5afe575SCameron Smith     &       c_loc(intfromfile), ieight, dataInt, iotype)
7159599516SKenneth E. Jansen            neltp = intfromfile(1) ! -1 if fname2 was not found, >=0 otherwise
7259599516SKenneth E. Jansen          end do
7359599516SKenneth E. Jansen          itpblktot = iblk-1
7459599516SKenneth E. Jansen        end if
7559599516SKenneth E. Jansen
7659599516SKenneth E. Jansen        if (myrank == 0) then
7759599516SKenneth E. Jansen          write(*,*) 'Number of boundary topologies: ',itpblktot
7859599516SKenneth E. Jansen        endif
7959599516SKenneth E. Jansen
8059599516SKenneth E. Jansen        nelblb=0
8159599516SKenneth E. Jansen        mattyp=0
8259599516SKenneth E. Jansen        ndofl = ndof
8359599516SKenneth E. Jansen
8459599516SKenneth E. Jansen        do iblk = 1, itpblktot
8559599516SKenneth E. Jansen           writeLock=0;
86*5b7f36ccSCameron Smith            if(input_mode.ge.1)then
872efdc748SKenneth E. Jansen               write (fname2,"('connectivity boundary',i1)") iblk
882efdc748SKenneth E. Jansen            else
892efdc748SKenneth E. Jansen               write (fname2,"('connectivity boundary linear tetrahedron')")
902efdc748SKenneth E. Jansen            endif
912efdc748SKenneth E. Jansen
9259599516SKenneth E. Jansen           ! Synchronization for performance monitoring, as some parts do not include some topologies
9359599516SKenneth E. Jansen           call MPI_Barrier(MPI_COMM_WORLD,ierr)
942efdc748SKenneth E. Jansen           call phio_readheader(fhandle, fname2 // char(0),
95e5afe575SCameron Smith     &      c_loc(intfromfile), ieight, dataInt, iotype)
9659599516SKenneth E. Jansen           neltp =intfromfile(1)
9759599516SKenneth E. Jansen           nenl  =intfromfile(2)
9859599516SKenneth E. Jansen           ipordl=intfromfile(3)
9959599516SKenneth E. Jansen           nshl  =intfromfile(4)
10059599516SKenneth E. Jansen           nshlb =intfromfile(5)
10159599516SKenneth E. Jansen           nenbl =intfromfile(6)
10259599516SKenneth E. Jansen           lcsyst=intfromfile(7)
10359599516SKenneth E. Jansen           numnbc=intfromfile(8)
10459599516SKenneth E. Jansen
10559599516SKenneth E. Jansen           allocate (ientp(neltp,nshl))
10659599516SKenneth E. Jansen           allocate (iBCBtp(neltp,ndiBCB))
10759599516SKenneth E. Jansen           allocate (BCBtp(neltp,ndBCB))
10859599516SKenneth E. Jansen           iientpsiz=neltp*nshl
10959599516SKenneth E. Jansen
11059599516SKenneth E. Jansen           if (neltp==0) then
11159599516SKenneth E. Jansen              writeLock=1;
11259599516SKenneth E. Jansen           endif
11359599516SKenneth E. Jansen
1142efdc748SKenneth E. Jansen           call phio_readdatablock(fhandle, fname2 // char(0),
115bc62cfd4SCameron Smith     &      c_loc(ientp),iientpsiz,dataInt,iotype)
11659599516SKenneth E. Jansenc
11759599516SKenneth E. Jansenc.... Read the boundary flux codes
11859599516SKenneth E. Jansenc
11959599516SKenneth E. Jansen           call MPI_BARRIER(MPI_COMM_WORLD, ierr)
120*5b7f36ccSCameron Smith            if(input_mode.ge.1)then
1212efdc748SKenneth E. Jansen               write (fname2,"('nbc codes',i1)") iblk
1222efdc748SKenneth E. Jansen            else
1232efdc748SKenneth E. Jansen               write (fname2,"('nbc codes linear tetrahedron')")
1242efdc748SKenneth E. Jansen            endif
12559599516SKenneth E. Jansen
1262efdc748SKenneth E. Jansen           call phio_readheader(fhandle, fname2 // char(0),
127e5afe575SCameron Smith     &      c_loc(intfromfile), ieight, dataInt, iotype)
12859599516SKenneth E. Jansen           iiBCBtpsiz=neltp*ndiBCB
1292efdc748SKenneth E. Jansen           call phio_readdatablock(fhandle, fname2 // char(0),
130bc62cfd4SCameron Smith     &      c_loc(iBCBtp),iiBCBtpsiz,dataInt,iotype)
13159599516SKenneth E. Jansenc
13259599516SKenneth E. Jansenc.... read the boundary condition data
13359599516SKenneth E. Jansenc
13459599516SKenneth E. Jansen           call MPI_BARRIER(MPI_COMM_WORLD, ierr)
135*5b7f36ccSCameron Smith            if(input_mode.ge.1)then
1362efdc748SKenneth E. Jansen               write (fname2,"('nbc values',i1)") iblk
1372efdc748SKenneth E. Jansen            else
1382efdc748SKenneth E. Jansen               write (fname2,"('nbc values linear tetrahedron')")
1392efdc748SKenneth E. Jansen            endif
14059599516SKenneth E. Jansen
1412efdc748SKenneth E. Jansen           call phio_readheader(fhandle, fname2 // char(0),
142e5afe575SCameron Smith     &      c_loc(intfromfile), ieight, dataInt, iotype)
14359599516SKenneth E. Jansen           BCBtp    = zero
14459599516SKenneth E. Jansen           iBCBtpsiz=neltp*ndBCB
1452efdc748SKenneth E. Jansen           call phio_readdatablock(fhandle, fname2 // char(0),
146bc62cfd4SCameron Smith     &      c_loc(BCBtp),iBCBtpsiz,dataDbl,iotype)
14759599516SKenneth E. Jansenc
14859599516SKenneth E. Jansenc This is a temporary fix until NSpre properly zeros this array where it
14959599516SKenneth E. Jansenc is not set.  DEC has indigestion with these arrays though the
15059599516SKenneth E. Jansenc result is never used (never effects solution).
15159599516SKenneth E. Jansenc
15259599516SKenneth E. Jansen           if(writeLock==0) then
15359599516SKenneth E. Jansen              where(.not.btest(iBCBtp(:,1),0)) BCBtp(:,1)=zero
15459599516SKenneth E. Jansen              where(.not.btest(iBCBtp(:,1),1)) BCBtp(:,2)=zero
15559599516SKenneth E. Jansen              where(.not.btest(iBCBtp(:,1),3)) BCBtp(:,6)=zero
15659599516SKenneth E. Jansen              if(ndBCB.gt.6) then
15759599516SKenneth E. Jansen                 do i=6,ndof
15859599516SKenneth E. Jansen                    where(.not.btest(iBCBtp(:,1),i-1)) BCBtp(:,i+1)=zero
15959599516SKenneth E. Jansen                 enddo
16059599516SKenneth E. Jansen              endif
16159599516SKenneth E. Jansen              where(.not.btest(iBCBtp(:,1),2))
16259599516SKenneth E. Jansen                 BCBtp(:,3)=zero
16359599516SKenneth E. Jansen                 BCBtp(:,4)=zero
16459599516SKenneth E. Jansen                 BCBtp(:,5)=zero
16559599516SKenneth E. Jansen              endwhere
16659599516SKenneth E. Jansen
16759599516SKenneth E. Jansen              do n=1,neltp,ibksz
16859599516SKenneth E. Jansen                 nelblb=nelblb+1
16959599516SKenneth E. Jansen                 npro= min(IBKSZ, neltp - n + 1)
17059599516SKenneth E. Jansen                 lcblkb(1,nelblb)  = iel
17159599516SKenneth E. Jansen                 lcblkb(3,nelblb)  = lcsyst
17259599516SKenneth E. Jansen                 lcblkb(4,nelblb)  = ipordl
17359599516SKenneth E. Jansen                 lcblkb(5,nelblb)  = nenl
17459599516SKenneth E. Jansen                 lcblkb(6,nelblb)  = nenbl
17559599516SKenneth E. Jansen                 lcblkb(7,nelblb)  = mattyp
17659599516SKenneth E. Jansen                 lcblkb(8,nelblb)  = ndofl
17759599516SKenneth E. Jansen                 lcblkb(9,nelblb)  = nshl
17859599516SKenneth E. Jansen                 lcblkb(10,nelblb) = nshlb ! # of shape functions per elt
17959599516SKenneth E. Jansenc
18059599516SKenneth E. Jansenc.... save the element block
18159599516SKenneth E. Jansenc
18259599516SKenneth E. Jansen                 n1=n
18359599516SKenneth E. Jansen                 n2=n+npro-1
18459599516SKenneth E. Jansen                 materb=1       ! all one material for now
18559599516SKenneth E. Jansenc
18659599516SKenneth E. Jansenc.... allocate memory for stack arrays
18759599516SKenneth E. Jansenc
18859599516SKenneth E. Jansen                 allocate (mienb(nelblb)%p(npro,nshl))
18959599516SKenneth E. Jansen                 allocate (miBCB(nelblb)%p(npro,ndiBCB))
19059599516SKenneth E. Jansen                 allocate (mBCB(nelblb)%p(npro,nshlb,ndBCB))
19159599516SKenneth E. Jansen                 allocate (mmatb(nelblb)%p(npro))
19259599516SKenneth E. Jansenc
19359599516SKenneth E. Jansenc.... save the boundary element block
19459599516SKenneth E. Jansenc
19559599516SKenneth E. Jansen                 call gensvb (ientp(n1:n2,1:nshl),
19659599516SKenneth E. Jansen     &                iBCBtp(n1:n2,:),      BCBtp(n1:n2,:),
19759599516SKenneth E. Jansen     &                materb,        mienb(nelblb)%p,
19859599516SKenneth E. Jansen     &                miBCB(nelblb)%p,        mBCB(nelblb)%p,
19959599516SKenneth E. Jansen     &                mmatb(nelblb)%p)
20059599516SKenneth E. Jansen                 iel=iel+npro
20159599516SKenneth E. Jansen              enddo
20259599516SKenneth E. Jansen           endif
20359599516SKenneth E. Jansen           deallocate(ientp)
20459599516SKenneth E. Jansen           deallocate(iBCBtp)
20559599516SKenneth E. Jansen           deallocate(BCBtp)
20659599516SKenneth E. Jansen
20759599516SKenneth E. Jansen        enddo
20859599516SKenneth E. Jansen        lcblkb(1,nelblb+1) = iel
20959599516SKenneth E. Jansen        return
21059599516SKenneth E. Jansenc
21159599516SKenneth E. Jansenc.... end of file error handling
21259599516SKenneth E. Jansenc
21359599516SKenneth E. Jansen 911    call error ('genbcb  ','end file',igeomBAK)
21459599516SKenneth E. Jansen1000    format(a80,//,
21559599516SKenneth 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',//,
21659599516SKenneth E. Jansen     &  '   Elem   BC codes',/,
21759599516SKenneth E. Jansen     &  '  Number  C P V H ',5x,27('Node',i1,:,2x))
21859599516SKenneth E. Jansen1100    format(2x,i5,2x,4i2,3x,27i7)
21959599516SKenneth E. Jansen2100    format(2x,i5,1p,1x,6e12.4)
22059599516SKenneth E. Jansen        end
221