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