xref: /phasta/phSolver/common/genblk.f (revision 0bfd29162c763803e760e1c0b5b3d9fbb074d255)
1        subroutine genblk (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
12c
13        include "common.h"
14!MR CHANGE
15        include "mpif.h" !Required to determine the max for itpblk
16!MR CHANGE END
17c
18        integer, allocatable :: ientp(:,:)
19        integer mater(ibksz)
20        integer intfromfile(50) ! integers read from headers
21        character*255 fname1
22
23cccccccccccccc New Phasta IO starts here ccccccccccccccccccccccccc
24
25        integer :: descriptor, descriptorG, GPID, color, nfiles
26        integer ::  numparts, writeLock
27        integer :: ierr_io, numprocs, itmp, itmp2
28!MR CHANGE
29        integer :: itpblktot,ierr,iseven
30!MR CHANGE END
31        character*255 fnamer, fname2, temp2
32        character*64 temp1, temp3
33!THIS NEEDS TO BE CLEANED - MR
34        nfiles = nsynciofiles
35!        nfields = nsynciofieldsreadgeombc
36        numparts = numpe !This is the common settings. Beware if you try to compute several parts per process
37
38!        nppp = numparts/numpe
39!        nppf = numparts/nfiles
40
41        color = int(myrank/(numparts/nfiles)) !Should call the SyncIO routine here
42        itmp2 = int(log10(float(color+1)))+1
43        write (temp2,"('(''geombc-dat.'',i',i1,')')") itmp2
44        temp2=trim(temp2)
45        write (fnamer,temp2) (color+1)
46        fnamer=trim(fnamer)
47
48        ione=1
49        itwo=2
50        iseven=7
51        ieleven=11
52        itmp = int(log10(float(myrank+1)))+1
53
54cccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
55
56c
57        iel=1
58        itpblk=nelblk
59!MR CHANGE
60
61        ! Get the total number of different interior topologies in the whole domain.
62        ! Try to read from a field. If the field does not exist, scan the geombc file.
63        itpblktot=-1
64        call phio_readheader(igeom,
65     &   'total number of interior tpblocks' // char(0),
66     &   itpblktot,ione,'integer' // char(0),iotype)
67
68!        write (*,*) 'Rank: ',myrank,' interior itpblktot intermediate:',
69!     &               itpblktot
70
71        if (itpblktot == -1) then
72          ! The field 'total number of different interior tpblocks' was not found in the geombc file.
73          ! Scan all the geombc file for the 'connectivity interior' fields to get this information.
74          iblk=0
75          neltp=0
76          do while(neltp .ne. -1)
77
78            ! intfromfile is reinitialized to -1 every time.
79            ! If connectivity interior@xxx is not found, then
80            ! readheader will return intfromfile unchanged
81
82            intfromfile(:)=-1
83            iblk = iblk+1
84            write (fname2,"('connectivity interior',i1)") iblk
85
86            !write(*,*) 'rank, fname2',myrank, trim(adjustl(fname2))
87            call phio_readheader(igeom,fname2 // char(0),intfromfile,
88     &       iseven,'integer' // char(0),iotype)
89            neltp = intfromfile(1) ! -1 if fname2 was not found, >=0 otherwise
90          end do
91          itpblktot = iblk-1
92        end if
93
94        if (myrank == 0) then
95          write(*,*) 'Number of interior topologies: ',itpblktot
96        endif
97!        write (*,*) 'Rank: ',myrank,' interior itpblktot final:',
98!     &               itpblktot
99
100!MR CHANGE END
101
102        nelblk=0
103        mattyp = 0
104        ndofl = ndof
105        nsymdl = nsymdf
106
107!        call initphmpiio( nfields, nppf, nfiles, igeom )
108!        call openfile( fnamer, 'read', igeom )
109
110!         do iblk = 1, itpblk
111        do iblk = 1, itpblktot
112           writeLock=0;
113!MR CHANGE END
114c
115c           read(igeomBAK) neltp,nenl,ipordl,nshl, ijunk, ijunk, lcsyst
116c           call creadlist(igeomBAK,iseven,
117c     &          neltp,nenl,ipordl,nshl, ijunk, ijunk, lcsyst)
118
119ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
120
121           write (fname2,"('connectivity interior',i1)") iblk
122
123ccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccccc
124
125c           fname1='connectivity interior?'
126
127           ! Synchronization for performance monitoring, as some parts do not include some topologies
128           call MPI_Barrier(MPI_COMM_WORLD,ierr)
129           call phio_readheader(igeom,fname2 // char(0) ,intfromfile,
130     &     iseven,"integer" // char(0), iotype)
131           neltp  =intfromfile(1)
132           nenl   =intfromfile(2)
133           ipordl =intfromfile(3)
134           nshl   =intfromfile(4)
135           ijunk  =intfromfile(5)
136           ijunk  =intfromfile(6)
137           lcsyst =intfromfile(7)
138           allocate (ientp(neltp,nshl))
139c           read(igeomBAK) ientp
140           iientpsiz=neltp*nshl
141
142           if (neltp==0) then
143              writeLock=1;
144           endif
145
146           call phio_readdatablock(igeom,fname2 // char(0),ientp,iientpsiz,
147     &                     "integer" // char(0), iotype)
148
149!            call closefile( igeom, "read" // char(0) )
150!            call finalizephmpiio( igeom )
151
152!MR CHANGE
153           if(writeLock==0) then
154!MR CHANGE
155
156             do n=1,neltp,ibksz
157                nelblk=nelblk+1
158                npro= min(IBKSZ, neltp - n + 1)
159c
160                lcblk(1,nelblk)  = iel
161c                lcblk(2,nelblk)  = iopen ! available for later use
162                lcblk(3,nelblk)  = lcsyst
163                lcblk(4,nelblk)  = ipordl
164                lcblk(5,nelblk)  = nenl
165                lcblk(6,nelblk)  = nfacel
166                lcblk(7,nelblk)  = mattyp
167                lcblk(8,nelblk)  = ndofl
168                lcblk(9,nelblk)  = nsymdl
169                lcblk(10,nelblk) = nshl ! # of shape functions per elt
170c
171c.... allocate memory for stack arrays
172c
173                allocate (mmat(nelblk)%p(npro))
174c
175                allocate (mien(nelblk)%p(npro,nshl))
176                allocate (mxmudmi(nelblk)%p(npro,maxsh))
177c
178c.... save the element block
179c
180                n1=n
181                n2=n+npro-1
182                mater=1   ! all one material for now
183                call gensav (ientp(n1:n2,1:nshl),
184     &                       mater,           mien(nelblk)%p,
185     &                       mmat(nelblk)%p)
186                iel=iel+npro
187c
188             enddo
189!MR CHANGE
190           endif
191!MR CHANGE
192           deallocate(ientp)
193        enddo
194
195!        call closefile( igeom, "read" // char(0) )
196!        call finalizephmpiio( igeom )
197
198        lcblk(1,nelblk+1) = iel
199c
200c.... return
201c
202CAD        call timer ('Back    ')
203c
204        return
205c
2061000    format(a80,//,
207     &  ' N o d a l   C o n n e c t i v i t y',//,
208     &  '   Elem  ',/,
209     &  '  Number  ',7x,27('Node',i2,:,2x))
2101100    format(2x,i5,6x,27i8)
211        end
212