xref: /phasta/phSolver/common/genblkSyncIO.f (revision ab5b07a458969128644e941ce150fa049270fdf1)
1*ab5b07a4SKenneth E. Jansen        subroutine genblkSyncIO (IBKSZ)
29d714148SKenneth E. Jansenc
39d714148SKenneth E. Jansenc----------------------------------------------------------------------
49d714148SKenneth E. Jansenc
59d714148SKenneth E. Jansenc  This routine reads the interior elements and generates the
69d714148SKenneth E. Jansenc  appropriate blocks.
79d714148SKenneth E. Jansenc
89d714148SKenneth E. Jansenc Zdenek Johan, Fall 1991.
99d714148SKenneth E. Jansenc----------------------------------------------------------------------
109d714148SKenneth E. Jansenc
119d714148SKenneth E. Jansen        use pointer_data
12*ab5b07a4SKenneth E. Jansen        use phio
13*ab5b07a4SKenneth E. Jansen        use iso_c_binding
149d714148SKenneth E. Jansen        include "common.h"
159d714148SKenneth E. Jansen        include "mpif.h" !Required to determine the max for itpblk
16*ab5b07a4SKenneth E. Jansen
17*ab5b07a4SKenneth E. Jansen        integer, target, allocatable :: ientp(:,:)
189d714148SKenneth E. Jansen        integer mater(ibksz)
19*ab5b07a4SKenneth E. Jansen        integer, target :: intfromfile(50) ! integers read from headers
209d714148SKenneth E. Jansen        character*255 fname1
21*ab5b07a4SKenneth E. Jansen        integer :: descriptor, descriptorG, GPID, color
229d714148SKenneth E. Jansen        integer ::  numparts, writeLock
23*ab5b07a4SKenneth E. Jansen        integer :: ierr_io, numprocs
24*ab5b07a4SKenneth E. Jansen        integer, target :: itpblktot,ierr,iseven
25*ab5b07a4SKenneth E. Jansen        character*255 fname2
26*ab5b07a4SKenneth E. Jansen        character(len=30) :: dataInt
27*ab5b07a4SKenneth E. Jansen        dataInt = c_char_'integer'//c_null_char
289d714148SKenneth E. Jansen        numparts = numpe !This is the common settings. Beware if you try to compute several parts per process
299d714148SKenneth E. Jansen        ione=1
309d714148SKenneth E. Jansen        itwo=2
319d714148SKenneth E. Jansen        iseven=7
329d714148SKenneth E. Jansen        ieleven=11
339d714148SKenneth E. Jansen        iel=1
349d714148SKenneth E. Jansen        itpblk=nelblk
359d714148SKenneth E. Jansen
369d714148SKenneth E. Jansen        ! Get the total number of different interior topologies in the whole domain.
379d714148SKenneth E. Jansen        ! Try to read from a field. If the field does not exist, scan the geombc file.
389d714148SKenneth E. Jansen          itpblktot=-1
39*ab5b07a4SKenneth E. Jansen        call phio_readheader(fhandle,
40*ab5b07a4SKenneth E. Jansen     &   c_char_'total number of interior tpblocks' // char(0),
41*ab5b07a4SKenneth E. Jansen     &   c_loc(itpblktot), ione, dataInt, iotype)
429d714148SKenneth E. Jansen
439d714148SKenneth E. Jansen        if (itpblktot == -1) then
449d714148SKenneth E. Jansen          ! The field 'total number of different interior tpblocks' was not found in the geombc file.
459d714148SKenneth E. Jansen          ! Scan all the geombc file for the 'connectivity interior' fields to get this information.
469d714148SKenneth E. Jansen          iblk=0
479d714148SKenneth E. Jansen          neltp=0
489d714148SKenneth E. Jansen          do while(neltp .ne. -1)
499d714148SKenneth E. Jansen
509d714148SKenneth E. Jansen            ! intfromfile is reinitialized to -1 every time.
519d714148SKenneth E. Jansen            ! If connectivity interior@xxx is not found, then
529d714148SKenneth E. Jansen            ! readheader will return intfromfile unchanged
539d714148SKenneth E. Jansen
549d714148SKenneth E. Jansen            intfromfile(:)=-1
559d714148SKenneth E. Jansen            iblk = iblk+1
56*ab5b07a4SKenneth E. Jansen              write (fname2,"('connectivity interior',i1)") iblk
57*ab5b07a4SKenneth E. Jansen            call phio_readheader(fhandle, fname2 // char(0),
58*ab5b07a4SKenneth E. Jansen     &       c_loc(intfromfile), iseven, dataInt, iotype)
599d714148SKenneth E. Jansen            neltp = intfromfile(1) ! -1 if fname2 was not found, >=0 otherwise
609d714148SKenneth E. Jansen          end do
619d714148SKenneth E. Jansen          itpblktot = iblk-1
629d714148SKenneth E. Jansen        end if
639d714148SKenneth E. Jansen
649d714148SKenneth E. Jansen        if (myrank == 0) then
659d714148SKenneth E. Jansen          write(*,*) 'Number of interior topologies: ',itpblktot
669d714148SKenneth E. Jansen        endif
679d714148SKenneth E. Jansen!        write (*,*) 'Rank: ',myrank,' interior itpblktot final:',
689d714148SKenneth E. Jansen!     &               itpblktot
699d714148SKenneth E. Jansen
709d714148SKenneth E. Jansen
719d714148SKenneth E. Jansen        nelblk=0
729d714148SKenneth E. Jansen        mattyp = 0
739d714148SKenneth E. Jansen        ndofl = ndof
749d714148SKenneth E. Jansen        nsymdl = nsymdf
759d714148SKenneth E. Jansen
769d714148SKenneth E. Jansen        do iblk = 1, itpblktot
779d714148SKenneth E. Jansen           writeLock=0;
78*ab5b07a4SKenneth E. Jansen           write (fname2,"('connectivity interior',i1)") iblk
799d714148SKenneth E. Jansen
80*ab5b07a4SKenneth E. Jansen           ! Synchronization for performance monitoring, as some parts do not include some topologies
819d714148SKenneth E. Jansen           call MPI_Barrier(MPI_COMM_WORLD,ierr)
82*ab5b07a4SKenneth E. Jansen           call phio_readheader(fhandle, fname2 // char(0),
83*ab5b07a4SKenneth E. Jansen     &      c_loc(intfromfile), iseven, dataInt, iotype)
849d714148SKenneth E. Jansen           neltp  =intfromfile(1)
859d714148SKenneth E. Jansen           nenl   =intfromfile(2)
869d714148SKenneth E. Jansen           ipordl =intfromfile(3)
879d714148SKenneth E. Jansen           nshl   =intfromfile(4)
889d714148SKenneth E. Jansen           ijunk  =intfromfile(5)
899d714148SKenneth E. Jansen           ijunk  =intfromfile(6)
909d714148SKenneth E. Jansen           lcsyst =intfromfile(7)
919d714148SKenneth E. Jansen           allocate (ientp(neltp,nshl))
929d714148SKenneth E. Jansen           iientpsiz=neltp*nshl
939d714148SKenneth E. Jansen
949d714148SKenneth E. Jansen           if (neltp==0) then
959d714148SKenneth E. Jansen              writeLock=1;
969d714148SKenneth E. Jansen           endif
979d714148SKenneth E. Jansen
98*ab5b07a4SKenneth E. Jansen           call phio_readdatablock(fhandle,fname2 // char(0),
99*ab5b07a4SKenneth E. Jansen     &      c_loc(ientp), iientpsiz, dataInt, iotype)
1009d714148SKenneth E. Jansen
1019d714148SKenneth E. Jansen           if(writeLock==0) then
1029d714148SKenneth E. Jansen             do n=1,neltp,ibksz
1039d714148SKenneth E. Jansen                nelblk=nelblk+1
1049d714148SKenneth E. Jansen                npro= min(IBKSZ, neltp - n + 1)
1059d714148SKenneth E. Jansen                lcblk(1,nelblk)  = iel
1069d714148SKenneth E. Jansen                lcblk(3,nelblk)  = lcsyst
1079d714148SKenneth E. Jansen                lcblk(4,nelblk)  = ipordl
1089d714148SKenneth E. Jansen                lcblk(5,nelblk)  = nenl
1099d714148SKenneth E. Jansen                lcblk(6,nelblk)  = nfacel
1109d714148SKenneth E. Jansen                lcblk(7,nelblk)  = mattyp
1119d714148SKenneth E. Jansen                lcblk(8,nelblk)  = ndofl
1129d714148SKenneth E. Jansen                lcblk(9,nelblk)  = nsymdl
1139d714148SKenneth E. Jansen                lcblk(10,nelblk) = nshl ! # of shape functions per elt
1149d714148SKenneth E. Jansenc
1159d714148SKenneth E. Jansenc.... allocate memory for stack arrays
1169d714148SKenneth E. Jansenc
1179d714148SKenneth E. Jansen                allocate (mmat(nelblk)%p(npro))
1189d714148SKenneth E. Jansenc
1199d714148SKenneth E. Jansen                allocate (mien(nelblk)%p(npro,nshl))
1209d714148SKenneth E. Jansen                allocate (mxmudmi(nelblk)%p(npro,maxsh))
121*ab5b07a4SKenneth E. Jansen                if(usingpetsc.eq.0) then
122*ab5b07a4SKenneth E. Jansen                    allocate (mienG(nelblk)%p(1,1))
123*ab5b07a4SKenneth E. Jansen                else
124*ab5b07a4SKenneth E. Jansen                    allocate (mienG(nelblk)%p(npro,nshl))
125*ab5b07a4SKenneth E. Jansen                endif
126*ab5b07a4SKenneth E. Jansen                ! note mienG will be passed to gensav but nothing filled if not
127*ab5b07a4SKenneth E. Jansen                ! using PETSc so this is safe
1289d714148SKenneth E. Jansenc
1299d714148SKenneth E. Jansenc.... save the element block
1309d714148SKenneth E. Jansenc
1319d714148SKenneth E. Jansen                n1=n
1329d714148SKenneth E. Jansen                n2=n+npro-1
1339d714148SKenneth E. Jansen                mater=1   ! all one material for now
1349d714148SKenneth E. Jansen                call gensav (ientp(n1:n2,1:nshl),
1359d714148SKenneth E. Jansen     &                       mater,           mien(nelblk)%p,
136*ab5b07a4SKenneth E. Jansen     &                       mienG(nelblk)%p,
1379d714148SKenneth E. Jansen     &                       mmat(nelblk)%p)
1389d714148SKenneth E. Jansen                iel=iel+npro
1399d714148SKenneth E. Jansen             enddo
1409d714148SKenneth E. Jansen           endif
1419d714148SKenneth E. Jansen           deallocate(ientp)
1429d714148SKenneth E. Jansen        enddo
1439d714148SKenneth E. Jansen
1449d714148SKenneth E. Jansen        lcblk(1,nelblk+1) = iel
1459d714148SKenneth E. Jansenc
1469d714148SKenneth E. Jansenc.... return
1479d714148SKenneth E. Jansenc
1489d714148SKenneth E. Jansenc
1499d714148SKenneth E. Jansen        return
1509d714148SKenneth E. Jansenc
1519d714148SKenneth E. Jansen1000    format(a80,//,
1529d714148SKenneth E. Jansen     &  ' N o d a l   C o n n e c t i v i t y',//,
1539d714148SKenneth E. Jansen     &  '   Elem  ',/,
1549d714148SKenneth E. Jansen     &  '  Number  ',7x,27('Node',i2,:,2x))
1559d714148SKenneth E. Jansen1100    format(2x,i5,6x,27i8)
1569d714148SKenneth E. Jansen        end
157