xref: /phasta/phSolver/common/genblk.f (revision 2efdc748c163fca218950001636e2694eb132b46)
159599516SKenneth E. Jansen        subroutine genblk (IBKSZ)
259599516SKenneth E. Jansenc
359599516SKenneth E. Jansenc----------------------------------------------------------------------
459599516SKenneth E. Jansenc
559599516SKenneth E. Jansenc  This routine reads the interior elements and generates the
659599516SKenneth E. Jansenc  appropriate blocks.
759599516SKenneth E. Jansenc
859599516SKenneth E. Jansenc Zdenek Johan, Fall 1991.
959599516SKenneth E. Jansenc----------------------------------------------------------------------
1059599516SKenneth E. Jansenc
1159599516SKenneth E. Jansen        use pointer_data
12e5afe575SCameron Smith        use phio
13e5afe575SCameron Smith        use iso_c_binding
1459599516SKenneth E. Jansenc
1559599516SKenneth E. Jansen        include "common.h"
1659599516SKenneth E. Jansen!MR CHANGE
1759599516SKenneth E. Jansen        include "mpif.h" !Required to determine the max for itpblk
1859599516SKenneth E. Jansen!MR CHANGE END
1959599516SKenneth E. Jansenc
209a6935e5SKenneth E. Jansen        integer, target, allocatable :: ientp(:,:)
2159599516SKenneth E. Jansen        integer mater(ibksz)
229a6935e5SKenneth E. Jansen        integer, target :: intfromfile(50) ! integers read from headers
2359599516SKenneth E. Jansen        character*255 fname1
2459599516SKenneth E. Jansen
2559599516SKenneth E. Jansencccccccccccccc New Phasta IO starts here ccccccccccccccccccccccccc
2659599516SKenneth E. Jansen
2759599516SKenneth E. Jansen        integer :: descriptor, descriptorG, GPID, color, nfiles
2859599516SKenneth E. Jansen        integer ::  numparts, writeLock
29*2efdc748SKenneth E. Jansen        integer :: ierr_io, numprocs
3059599516SKenneth E. Jansen!MR CHANGE
319a6935e5SKenneth E. Jansen        integer, target :: itpblktot,ierr,iseven
3259599516SKenneth E. Jansen!MR CHANGE END
33*2efdc748SKenneth E. Jansen        character*255 fname2
34e5afe575SCameron Smith
35e5afe575SCameron Smith        character(len=30) :: dataInt
36d5d2f64dSCameron Smith        dataInt = c_char_'integer'//c_null_char
37e5afe575SCameron Smith
3859599516SKenneth E. Jansen!THIS NEEDS TO BE CLEANED - MR
3959599516SKenneth E. Jansen        nfiles = nsynciofiles
4059599516SKenneth E. Jansen!        nfields = nsynciofieldsreadgeombc
4159599516SKenneth E. Jansen        numparts = numpe !This is the common settings. Beware if you try to compute several parts per process
4259599516SKenneth E. Jansen
4359599516SKenneth E. Jansen
4459599516SKenneth E. Jansen        ione=1
4559599516SKenneth E. Jansen        itwo=2
4659599516SKenneth E. Jansen        iseven=7
4759599516SKenneth E. Jansen        ieleven=11
4859599516SKenneth E. Jansen
4959599516SKenneth E. Jansencccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
5059599516SKenneth E. Jansen
5159599516SKenneth E. Jansenc
5259599516SKenneth E. Jansen        iel=1
5359599516SKenneth E. Jansen        itpblk=nelblk
5459599516SKenneth E. Jansen!MR CHANGE
5559599516SKenneth E. Jansen
5659599516SKenneth E. Jansen        ! Get the total number of different interior topologies in the whole domain.
5759599516SKenneth E. Jansen        ! Try to read from a field. If the field does not exist, scan the geombc file.
58*2efdc748SKenneth E. Jansen          itpblktot=1  ! hardwired to montopology for now
59d5d2f64dSCameron Smith        call phio_readheader(fhandle,
60e5afe575SCameron Smith     &   c_char_'total number of interior tpblocks' // char(0),
61e5afe575SCameron Smith     &   c_loc(itpblktot), ione, dataInt, iotype)
6259599516SKenneth E. Jansen
6359599516SKenneth E. Jansen!        write (*,*) 'Rank: ',myrank,' interior itpblktot intermediate:',
6459599516SKenneth E. Jansen!     &               itpblktot
6559599516SKenneth E. Jansen
6659599516SKenneth E. Jansen        if (itpblktot == -1) then
6759599516SKenneth E. Jansen          ! The field 'total number of different interior tpblocks' was not found in the geombc file.
6859599516SKenneth E. Jansen          ! Scan all the geombc file for the 'connectivity interior' fields to get this information.
6959599516SKenneth E. Jansen          iblk=0
7059599516SKenneth E. Jansen          neltp=0
7159599516SKenneth E. Jansen          do while(neltp .ne. -1)
7259599516SKenneth E. Jansen
7359599516SKenneth E. Jansen            ! intfromfile is reinitialized to -1 every time.
7459599516SKenneth E. Jansen            ! If connectivity interior@xxx is not found, then
7559599516SKenneth E. Jansen            ! readheader will return intfromfile unchanged
7659599516SKenneth E. Jansen
7759599516SKenneth E. Jansen            intfromfile(:)=-1
7859599516SKenneth E. Jansen            iblk = iblk+1
79*2efdc748SKenneth E. Jansen            if(nfiles.gt.0) then
80aa9d7345SCameron Smith              write (fname2,"('connectivity interior',i1)") iblk
81*2efdc748SKenneth E. Jansen            else
82*2efdc748SKenneth E. Jansen              write (fname2,"('connectivity interior linear tetrahedron')")
83*2efdc748SKenneth E. Jansen            endif
8459599516SKenneth E. Jansen
8559599516SKenneth E. Jansen            !write(*,*) 'rank, fname2',myrank, trim(adjustl(fname2))
86d5d2f64dSCameron Smith            call phio_readheader(fhandle, fname2 // char(0),
87e5afe575SCameron Smith     &       c_loc(intfromfile), iseven, dataInt, iotype)
8859599516SKenneth E. Jansen            neltp = intfromfile(1) ! -1 if fname2 was not found, >=0 otherwise
8959599516SKenneth E. Jansen          end do
9059599516SKenneth E. Jansen          itpblktot = iblk-1
9159599516SKenneth E. Jansen        end if
9259599516SKenneth E. Jansen
9359599516SKenneth E. Jansen        if (myrank == 0) then
9459599516SKenneth E. Jansen          write(*,*) 'Number of interior topologies: ',itpblktot
9559599516SKenneth E. Jansen        endif
9659599516SKenneth E. Jansen!        write (*,*) 'Rank: ',myrank,' interior itpblktot final:',
9759599516SKenneth E. Jansen!     &               itpblktot
9859599516SKenneth E. Jansen
9959599516SKenneth E. Jansen!MR CHANGE END
10059599516SKenneth E. Jansen
10159599516SKenneth E. Jansen        nelblk=0
10259599516SKenneth E. Jansen        mattyp = 0
10359599516SKenneth E. Jansen        ndofl = ndof
10459599516SKenneth E. Jansen        nsymdl = nsymdf
10559599516SKenneth E. Jansen
10659599516SKenneth E. Jansen!        call initphmpiio( nfields, nppf, nfiles, igeom )
10759599516SKenneth E. Jansen!        call openfile( fnamer, 'read', igeom )
10859599516SKenneth E. Jansen
10959599516SKenneth E. Jansen!         do iblk = 1, itpblk
11059599516SKenneth E. Jansen        do iblk = 1, itpblktot
11159599516SKenneth E. Jansen           writeLock=0;
11259599516SKenneth E. Jansen!MR CHANGE END
11359599516SKenneth E. Jansenc
11459599516SKenneth E. Jansenc           read(igeomBAK) neltp,nenl,ipordl,nshl, ijunk, ijunk, lcsyst
11559599516SKenneth E. Jansenc           call creadlist(igeomBAK,iseven,
11659599516SKenneth E. Jansenc     &          neltp,nenl,ipordl,nshl, ijunk, ijunk, lcsyst)
11759599516SKenneth E. Jansen
11859599516SKenneth E. Jansenccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
11959599516SKenneth E. Jansen
120*2efdc748SKenneth E. Jansen            if(nfiles.gt.0) then
121aa9d7345SCameron Smith              write (fname2,"('connectivity interior',i1)") iblk
122*2efdc748SKenneth E. Jansen            else
123*2efdc748SKenneth E. Jansen              write (fname2,"('connectivity interior linear tetrahedron')")
124*2efdc748SKenneth E. Jansen            endif
12559599516SKenneth E. Jansen
12659599516SKenneth E. Jansenccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
12759599516SKenneth E. Jansen
12859599516SKenneth E. Jansenc           fname1='connectivity interior?'
12959599516SKenneth E. Jansen
13059599516SKenneth E. Jansen           ! Synchronization for performance monitoring, as some parts do not include some topologies
13159599516SKenneth E. Jansen           call MPI_Barrier(MPI_COMM_WORLD,ierr)
132d5d2f64dSCameron Smith           call phio_readheader(fhandle, fname2 // char(0),
133e5afe575SCameron Smith     &      c_loc(intfromfile), iseven, dataInt, iotype)
13459599516SKenneth E. Jansen           neltp  =intfromfile(1)
13559599516SKenneth E. Jansen           nenl   =intfromfile(2)
13659599516SKenneth E. Jansen           ipordl =intfromfile(3)
13759599516SKenneth E. Jansen           nshl   =intfromfile(4)
13859599516SKenneth E. Jansen           ijunk  =intfromfile(5)
13959599516SKenneth E. Jansen           ijunk  =intfromfile(6)
14059599516SKenneth E. Jansen           lcsyst =intfromfile(7)
14159599516SKenneth E. Jansen           allocate (ientp(neltp,nshl))
14259599516SKenneth E. Jansenc           read(igeomBAK) ientp
14359599516SKenneth E. Jansen           iientpsiz=neltp*nshl
14459599516SKenneth E. Jansen
14559599516SKenneth E. Jansen           if (neltp==0) then
14659599516SKenneth E. Jansen              writeLock=1;
14759599516SKenneth E. Jansen           endif
14859599516SKenneth E. Jansen
149d5d2f64dSCameron Smith           call phio_readdatablock(fhandle,fname2 // char(0),
150bc62cfd4SCameron Smith     &      c_loc(ientp), iientpsiz, dataInt, iotype)
15159599516SKenneth E. Jansen
15259599516SKenneth E. Jansen!            call closefile( igeom, "read" // char(0) )
15359599516SKenneth E. Jansen!            call finalizephmpiio( igeom )
15459599516SKenneth E. Jansen
15559599516SKenneth E. Jansen!MR CHANGE
15659599516SKenneth E. Jansen           if(writeLock==0) then
15759599516SKenneth E. Jansen!MR CHANGE
15859599516SKenneth E. Jansen
15959599516SKenneth E. Jansen             do n=1,neltp,ibksz
16059599516SKenneth E. Jansen                nelblk=nelblk+1
16159599516SKenneth E. Jansen                npro= min(IBKSZ, neltp - n + 1)
16259599516SKenneth E. Jansenc
16359599516SKenneth E. Jansen                lcblk(1,nelblk)  = iel
16459599516SKenneth E. Jansenc                lcblk(2,nelblk)  = iopen ! available for later use
16559599516SKenneth E. Jansen                lcblk(3,nelblk)  = lcsyst
16659599516SKenneth E. Jansen                lcblk(4,nelblk)  = ipordl
16759599516SKenneth E. Jansen                lcblk(5,nelblk)  = nenl
16859599516SKenneth E. Jansen                lcblk(6,nelblk)  = nfacel
16959599516SKenneth E. Jansen                lcblk(7,nelblk)  = mattyp
17059599516SKenneth E. Jansen                lcblk(8,nelblk)  = ndofl
17159599516SKenneth E. Jansen                lcblk(9,nelblk)  = nsymdl
17259599516SKenneth E. Jansen                lcblk(10,nelblk) = nshl ! # of shape functions per elt
17359599516SKenneth E. Jansenc
17459599516SKenneth E. Jansenc.... allocate memory for stack arrays
17559599516SKenneth E. Jansenc
17659599516SKenneth E. Jansen                allocate (mmat(nelblk)%p(npro))
17759599516SKenneth E. Jansenc
17859599516SKenneth E. Jansen                allocate (mien(nelblk)%p(npro,nshl))
17959599516SKenneth E. Jansen                allocate (mxmudmi(nelblk)%p(npro,maxsh))
18059599516SKenneth E. Jansenc
18159599516SKenneth E. Jansenc.... save the element block
18259599516SKenneth E. Jansenc
18359599516SKenneth E. Jansen                n1=n
18459599516SKenneth E. Jansen                n2=n+npro-1
18559599516SKenneth E. Jansen                mater=1   ! all one material for now
18659599516SKenneth E. Jansen                call gensav (ientp(n1:n2,1:nshl),
18759599516SKenneth E. Jansen     &                       mater,           mien(nelblk)%p,
18859599516SKenneth E. Jansen     &                       mmat(nelblk)%p)
18959599516SKenneth E. Jansen                iel=iel+npro
19059599516SKenneth E. Jansenc
19159599516SKenneth E. Jansen             enddo
19259599516SKenneth E. Jansen!MR CHANGE
19359599516SKenneth E. Jansen           endif
19459599516SKenneth E. Jansen!MR CHANGE
19559599516SKenneth E. Jansen           deallocate(ientp)
19659599516SKenneth E. Jansen        enddo
19759599516SKenneth E. Jansen
19859599516SKenneth E. Jansen!        call closefile( igeom, "read" // char(0) )
19959599516SKenneth E. Jansen!        call finalizephmpiio( igeom )
20059599516SKenneth E. Jansen
20159599516SKenneth E. Jansen        lcblk(1,nelblk+1) = iel
20259599516SKenneth E. Jansenc
20359599516SKenneth E. Jansenc.... return
20459599516SKenneth E. Jansenc
20559599516SKenneth E. JansenCAD        call timer ('Back    ')
20659599516SKenneth E. Jansenc
20759599516SKenneth E. Jansen        return
20859599516SKenneth E. Jansenc
20959599516SKenneth E. Jansen1000    format(a80,//,
21059599516SKenneth E. Jansen     &  ' N o d a l   C o n n e c t i v i t y',//,
21159599516SKenneth E. Jansen     &  '   Elem  ',/,
21259599516SKenneth E. Jansen     &  '  Number  ',7x,27('Node',i2,:,2x))
21359599516SKenneth E. Jansen1100    format(2x,i5,6x,27i8)
21459599516SKenneth E. Jansen        end
215