xref: /phasta/phSolver/common/proces.f (revision 6b966dd8220c388f9c18e3b5a7cb164734f542c6)
1        subroutine proces
2c
3c----------------------------------------------------------------------
4c
5c This subroutine generates the problem data and calls the solution
6c  driver.
7c
8c
9c Zdenek Johan, Winter 1991.  (Fortran 90)
10c----------------------------------------------------------------------
11c
12        use readarrays          ! used to access x, iper, ilwork
13        use turbsa          ! used to access d2wall
14        use dtnmod
15        include "common.h"
16        include "mpif.h"
17c
18c arrays in the following 2 lines are now dimensioned in readnblk
19c        dimension x(numnp,nsd)
20c        dimension iper(nshg), ilwork(nlwork)
21c
22        dimension y(nshg,ndof),
23     &            iBC(nshg),
24     &            BC(nshg,ndofBC),
25     &            ac(nshg,ndof)
26c
27c.... shape function declarations
28c
29        dimension shp(MAXTOP,maxsh,MAXQPT),
30     &            shgl(MAXTOP,nsd,maxsh,MAXQPT),
31     &            shpb(MAXTOP,maxsh,MAXQPT),
32     &            shglb(MAXTOP,nsd,maxsh,MAXQPT)
33c
34c  stuff for dynamic model s.w.avg and wall model
35c
36c        dimension ifath(numnp),    velbar(nfath,nflow),
37c     &            nsons(nfath)
38
39        dimension velbar(nfath,nflow)
40c
41c stuff to interpolate profiles at inlet
42c
43        real*8 bcinterp(100,ndof+1),interp_mask(ndof)
44        logical exlog
45
46c        if ((irscale .ge. 0) .and. (myrank.eq.master)) then
47c           call setSPEBC(numnp, point2nfath, nsonmax)
48c        endif
49c
50c.... generate the geometry and boundary conditions data
51c
52        call gendat (y,              ac,             point2x,
53     &               iBC,            BC,
54     &               point2iper,     point2ilwork,   shp,
55     &               shgl,           shpb,           shglb,
56     &               point2ifath,    velbar,         point2nsons )
57        call setper(nshg)
58        call perprep(iBC,point2iper,nshg)
59        if (iLES/10 .eq. 1) then
60        call keeplhsG ! Allocating space for the mass (Gram) matrix to be
61                      ! used for projection filtering and reconstruction
62                      ! of the strain rate tensor.
63
64        call setrls   ! Allocating space for the resolved Leonard stresses
65c                         See bardmc.f
66        endif
67c
68c.... time averaged statistics
69c
70        if (ioform .eq. 2) then
71           call initStats(point2x, iBC, point2iper, point2ilwork)
72        endif
73c
74c.... RANS turbulence model
75c
76        if (iRANS .lt. 0) then
77           call initTurb( point2x )
78        endif
79c
80c.... p vs. Q boundary
81c
82           call initNABI( point2x, shpb )
83c
84c.... check for execution mode
85c
86        if (iexec .eq. 0) then
87           lstep = 0
88           call restar ('out ',  y  ,ac)
89           return
90        endif
91c
92c.... initialize AutoSponge
93c
94        if(matflg(5,1).ge.4) then ! cool case (sponge)
95           call initSponge( y,point2x)
96        endif
97c
98c
99c.... adjust BC's to interpolate from file
100c
101
102        inquire(file="inlet.dat",exist=exlog)
103        if(exlog) then
104           open (unit=654,file="inlet.dat",status="old")
105           read(654,*) ninterp,ibcmatch,(interp_mask(j),j=1,ndof)
106           do i=1,ninterp
107              read(654,*) (bcinterp(i,j),j=1,ndof+1) ! distance to wall+
108                        ! ndof but note order of BC's
109                        ! p,t,u,v,w,scalars. Also note that must be in
110                        ! increasing distance from the wall order.
111
112           enddo
113           do i=1,nshg  ! only correct for linears at this time
114              if(mod(iBC(i),1024).eq.ibcmatch) then
115                 iupper=0
116                 do j=2,ninterp
117                    if(bcinterp(j,1).gt.d2wall(i)) then !bound found
118                       xi=(d2wall(i)-bcinterp(j-1,1))/
119     &                    (bcinterp(j,1)-bcinterp(j-1,1))
120                       iupper=j
121                       exit
122                    endif
123                 enddo
124                 if(iupper.eq.0) then
125                    write(*,*) 'failure in finterp, ynint=',d2wall(i)
126                    stop
127                 endif
128                 if(interp_mask(1).ne.zero) then
129                    BC(i,1)=(xi*bcinterp(iupper,2)
130     &                +(one-xi)*bcinterp(iupper-1,2))*interp_mask(1)
131                 endif
132                 if(interp_mask(2).ne.zero) then
133                    BC(i,2)=(xi*bcinterp(iupper,3)
134     &                +(one-xi)*bcinterp(iupper-1,3))*interp_mask(2)
135                 endif
136                 if(interp_mask(3).ne.zero) then
137                    BC(i,3)=(xi*bcinterp(iupper,4)
138     &                +(one-xi)*bcinterp(iupper-1,4))*interp_mask(3)
139                 endif
140                 if(interp_masK(4).ne.zero) then
141                    BC(i,4)=(xi*bcinterp(iupper,5)
142     &                +(one-xi)*bcinterp(iupper-1,5))*interp_mask(4)
143                 endif
144                 if(interp_mask(5).ne.zero) then
145                    BC(i,5)=(xi*bcinterp(iupper,6)
146     &                +(one-xi)*bcinterp(iupper-1,6))*interp_mask(5)
147                 endif
148                 if(interp_mask(6).ne.zero) then
149                    BC(i,7)=(xi*bcinterp(iupper,7)
150     &                +(one-xi)*bcinterp(iupper-1,7))*interp_mask(6)
151                 endif
152                 if(interp_mask(7).ne.zero) then
153                    BC(i,8)=(xi*bcinterp(iupper,8)
154     &                +(one-xi)*bcinterp(iupper-1,8))*interp_mask(7)
155                 endif
156              endif
157           enddo
158        endif
159c$$$$$$$$$$$$$$$$$$$$
160
161c
162c
163c.... call the semi-discrete predictor multi-corrector iterative driver
164c
165        call itrdrv (y,              ac,
166     &               uold,           point2x,
167     &               iBC,            BC,
168     &               point2iper,     point2ilwork,   shp,
169     &               shgl,           shpb,           shglb,
170     &               point2ifath,    velbar,         point2nsons )
171c
172c.... return
173c
174c
175c.... stop CPU-timer
176c
177CAD        call timer ('End     ')
178c
179c.... close echo file
180c
181
182!MR CHANGE
183      if (numpe > 1) call MPI_BARRIER(MPI_COMM_WORLD, ierr)
184      if(myrank.eq.0)  then
185          write(*,*) 'process - before closing iecho'
186      endif
187!MR CHANGE END
188
189        close (iecho)
190
191!MR CHANGE
192      if (numpe > 1) call MPI_BARRIER(MPI_COMM_WORLD, ierr)
193      if(myrank.eq.0)  then
194          write(*,*) 'process - after closing iecho'
195      endif
196!MR CHANGE END
197
198
199c
200c.... end of the program
201c
202CAD        write(6,*) 'Life: ', second(0) - ttim(100)
203        deallocate(point2iper)
204        if(numpe.gt.1) then
205          call Dctypes(point2ilwork(1))
206          deallocate(point2ilwork)
207        endif
208        deallocate(point2x)
209        deallocate(point2nsons)
210        deallocate(point2ifath)
211        deallocate(uold)
212        deallocate(wnrm)
213        deallocate(otwn)
214        call finalizeDtN
215
216        return
217        end
218
219
220