xref: /phasta/phSolver/common/genblkSyncIO.f (revision cc72a73fd2b79f4dd0a850fa4af718cd73554811)
1        subroutine genblkSyncIO (IBKSZ)
2c
3c----------------------------------------------------------------------
4c
5c  This routine reads the interior elements and generates the
6c  appropriate blocks.
7c
8c Zdenek Johan, Fall 1991.
9c----------------------------------------------------------------------
10c
11        use pointer_data
12        use phio
13        use iso_c_binding
14        include "common.h"
15        include "mpif.h" !Required to determine the max for itpblk
16
17        integer, target, allocatable :: ientp(:,:)
18        integer mater(ibksz)
19        integer, target :: intfromfile(50) ! integers read from headers
20        character*255 fname1
21        integer :: descriptor, descriptorG, GPID, color
22        integer ::  numparts, writeLock
23        integer :: ierr_io, numprocs
24        integer, target :: itpblktot,ierr,iseven
25        character*255 fname2
26        character(len=30) :: dataInt
27        dataInt = c_char_'integer'//c_null_char
28        numparts = numpe !This is the common settings. Beware if you try to compute several parts per process
29        ione=1
30        itwo=2
31        iseven=7
32        ieleven=11
33        iel=1
34        itpblk=nelblk
35
36        ! Get the total number of different interior topologies in the whole domain.
37        ! Try to read from a field. If the field does not exist, scan the geombc file.
38          itpblktot=-1
39        call phio_readheader(fhandle,
40     &   c_char_'total number of interior tpblocks' // char(0),
41     &   c_loc(itpblktot), ione, dataInt, iotype)
42
43        if (itpblktot == -1) then
44          ! The field 'total number of different interior tpblocks' was not found in the geombc file.
45          ! Scan all the geombc file for the 'connectivity interior' fields to get this information.
46          iblk=0
47          neltp=0
48          do while(neltp .ne. -1)
49
50            ! intfromfile is reinitialized to -1 every time.
51            ! If connectivity interior@xxx is not found, then
52            ! readheader will return intfromfile unchanged
53
54            intfromfile(:)=-1
55            iblk = iblk+1
56              write (fname2,"('connectivity interior',i1)") iblk
57            call phio_readheader(fhandle, fname2 // char(0),
58     &       c_loc(intfromfile), iseven, dataInt, iotype)
59            neltp = intfromfile(1) ! -1 if fname2 was not found, >=0 otherwise
60          end do
61          itpblktot = iblk-1
62        end if
63
64        if (myrank == 0) then
65          write(*,*) 'Number of interior topologies: ',itpblktot
66        endif
67!        write (*,*) 'Rank: ',myrank,' interior itpblktot final:',
68!     &               itpblktot
69
70
71        nelblk=0
72        mattyp = 0
73        ndofl = ndof
74        nsymdl = nsymdf
75
76        do iblk = 1, itpblktot
77           writeLock=0;
78           write (fname2,"('connectivity interior',i1)") iblk
79
80           ! Synchronization for performance monitoring, as some parts do not include some topologies
81           call MPI_Barrier(MPI_COMM_WORLD,ierr)
82           call phio_readheader(fhandle, fname2 // char(0),
83     &      c_loc(intfromfile), iseven, dataInt, iotype)
84           neltp  =intfromfile(1)
85           nenl   =intfromfile(2)
86           ipordl =intfromfile(3)
87           nshl   =intfromfile(4)
88           ijunk  =intfromfile(5)
89           ijunk  =intfromfile(6)
90           lcsyst =intfromfile(7)
91           allocate (ientp(neltp,nshl))
92           iientpsiz=neltp*nshl
93
94           if (neltp==0) then
95              writeLock=1;
96           endif
97
98           call phio_readdatablock(fhandle,fname2 // char(0),
99     &      c_loc(ientp), iientpsiz, dataInt, iotype)
100
101           if(writeLock==0) then
102             do n=1,neltp,ibksz
103                nelblk=nelblk+1
104                npro= min(IBKSZ, neltp - n + 1)
105                lcblk(1,nelblk)  = iel
106                lcblk(3,nelblk)  = lcsyst
107                lcblk(4,nelblk)  = ipordl
108                lcblk(5,nelblk)  = nenl
109                lcblk(6,nelblk)  = nfacel
110                lcblk(7,nelblk)  = mattyp
111                lcblk(8,nelblk)  = ndofl
112                lcblk(9,nelblk)  = nsymdl
113                lcblk(10,nelblk) = nshl ! # of shape functions per elt
114c
115c.... allocate memory for stack arrays
116c
117                allocate (mmat(nelblk)%p(npro))
118c
119                allocate (mien(nelblk)%p(npro,nshl))
120                allocate (mxmudmi(nelblk)%p(npro,maxsh))
121                if(usingpetsc.eq.0) then
122                    allocate (mienG(nelblk)%p(1,1))
123                else
124                    allocate (mienG(nelblk)%p(npro,nshl))
125                endif
126                ! note mienG will be passed to gensav but nothing filled if not
127                ! using PETSc so this is safe
128c
129c.... save the element block
130c
131                n1=n
132                n2=n+npro-1
133                mater=1   ! all one material for now
134                call gensav (ientp(n1:n2,1:nshl),
135     &                       mater,           mien(nelblk)%p,
136     &                       mienG(nelblk)%p,
137     &                       mmat(nelblk)%p)
138                iel=iel+npro
139             enddo
140           endif
141           deallocate(ientp)
142        enddo
143
144        lcblk(1,nelblk+1) = iel
145c
146c.... return
147c
148c
149        return
150c
1511000    format(a80,//,
152     &  ' N o d a l   C o n n e c t i v i t y',//,
153     &  '   Elem  ',/,
154     &  '  Number  ',7x,27('Node',i2,:,2x))
1551100    format(2x,i5,6x,27i8)
156        end
157