xref: /phasta/phSolver/common/proces.f (revision 513954ef803c300cddba2bb96b4a5dac0b93489a)
159599516SKenneth E. Jansen        subroutine proces
259599516SKenneth E. Jansenc
359599516SKenneth E. Jansenc----------------------------------------------------------------------
459599516SKenneth E. Jansenc
559599516SKenneth E. Jansenc This subroutine generates the problem data and calls the solution
659599516SKenneth E. Jansenc  driver.
759599516SKenneth E. Jansenc
859599516SKenneth E. Jansenc
959599516SKenneth E. Jansenc Zdenek Johan, Winter 1991.  (Fortran 90)
1059599516SKenneth E. Jansenc----------------------------------------------------------------------
1159599516SKenneth E. Jansenc
1259599516SKenneth E. Jansen        use readarrays          ! used to access x, iper, ilwork
1359599516SKenneth E. Jansen        use turbsa          ! used to access d2wall
146b966dd8SCameron Smith        use dtnmod
15bd5c8633SCameron Smith        use periodicity
16bd5c8633SCameron Smith        use pvsQbi
1759599516SKenneth E. Jansen        include "common.h"
1859599516SKenneth E. Jansen        include "mpif.h"
1959599516SKenneth E. Jansenc
2059599516SKenneth E. Jansenc arrays in the following 2 lines are now dimensioned in readnblk
2159599516SKenneth E. Jansenc        dimension x(numnp,nsd)
2259599516SKenneth E. Jansenc        dimension iper(nshg), ilwork(nlwork)
2359599516SKenneth E. Jansenc
2459599516SKenneth E. Jansen        dimension y(nshg,ndof),
2559599516SKenneth E. Jansen     &            iBC(nshg),
2659599516SKenneth E. Jansen     &            BC(nshg,ndofBC),
2759599516SKenneth E. Jansen     &            ac(nshg,ndof)
2859599516SKenneth E. Jansenc
2959599516SKenneth E. Jansenc.... shape function declarations
3059599516SKenneth E. Jansenc
3159599516SKenneth E. Jansen        dimension shp(MAXTOP,maxsh,MAXQPT),
3259599516SKenneth E. Jansen     &            shgl(MAXTOP,nsd,maxsh,MAXQPT),
3359599516SKenneth E. Jansen     &            shpb(MAXTOP,maxsh,MAXQPT),
3459599516SKenneth E. Jansen     &            shglb(MAXTOP,nsd,maxsh,MAXQPT)
3559599516SKenneth E. Jansenc
3659599516SKenneth E. Jansenc  stuff for dynamic model s.w.avg and wall model
3759599516SKenneth E. Jansenc
3859599516SKenneth E. Jansenc        dimension ifath(numnp),    velbar(nfath,nflow),
3959599516SKenneth E. Jansenc     &            nsons(nfath)
4059599516SKenneth E. Jansen
4159599516SKenneth E. Jansen        dimension velbar(nfath,nflow)
4259599516SKenneth E. Jansenc
4359599516SKenneth E. Jansenc stuff to interpolate profiles at inlet
4459599516SKenneth E. Jansenc
4559599516SKenneth E. Jansen        real*8 bcinterp(100,ndof+1),interp_mask(ndof)
4659599516SKenneth E. Jansen        logical exlog
4759599516SKenneth E. Jansen
48*513954efSKenneth E. Jansen        !Duct
49*513954efSKenneth E. Jansen        real*8 c0, c1, c2, c3, x1, x2
50*513954efSKenneth E. Jansen        integer nn
51*513954efSKenneth E. Jansen
5259599516SKenneth E. Jansenc        if ((irscale .ge. 0) .and. (myrank.eq.master)) then
5359599516SKenneth E. Jansenc           call setSPEBC(numnp, point2nfath, nsonmax)
5459599516SKenneth E. Jansenc        endif
5559599516SKenneth E. Jansenc
5659599516SKenneth E. Jansenc.... generate the geometry and boundary conditions data
5759599516SKenneth E. Jansenc
5859599516SKenneth E. Jansen        call gendat (y,              ac,             point2x,
5959599516SKenneth E. Jansen     &               iBC,            BC,
6059599516SKenneth E. Jansen     &               point2iper,     point2ilwork,   shp,
6159599516SKenneth E. Jansen     &               shgl,           shpb,           shglb,
6259599516SKenneth E. Jansen     &               point2ifath,    velbar,         point2nsons )
6359599516SKenneth E. Jansen        call setper(nshg)
6459599516SKenneth E. Jansen        call perprep(iBC,point2iper,nshg)
6559599516SKenneth E. Jansen        if (iLES/10 .eq. 1) then
6659599516SKenneth E. Jansen        call keeplhsG ! Allocating space for the mass (Gram) matrix to be
6759599516SKenneth E. Jansen                      ! used for projection filtering and reconstruction
6859599516SKenneth E. Jansen                      ! of the strain rate tensor.
6959599516SKenneth E. Jansen
7059599516SKenneth E. Jansen        call setrls   ! Allocating space for the resolved Leonard stresses
7159599516SKenneth E. Jansenc                         See bardmc.f
7259599516SKenneth E. Jansen        endif
7359599516SKenneth E. Jansenc
7459599516SKenneth E. Jansenc.... time averaged statistics
7559599516SKenneth E. Jansenc
7659599516SKenneth E. Jansen        if (ioform .eq. 2) then
7759599516SKenneth E. Jansen           call initStats(point2x, iBC, point2iper, point2ilwork)
7859599516SKenneth E. Jansen        endif
7959599516SKenneth E. Jansenc
8059599516SKenneth E. Jansenc.... RANS turbulence model
8159599516SKenneth E. Jansenc
8259599516SKenneth E. Jansen        if (iRANS .lt. 0) then
8359599516SKenneth E. Jansen           call initTurb( point2x )
8459599516SKenneth E. Jansen        endif
8559599516SKenneth E. Jansenc
8659599516SKenneth E. Jansenc.... p vs. Q boundary
8759599516SKenneth E. Jansenc
8859599516SKenneth E. Jansen           call initNABI( point2x, shpb )
8959599516SKenneth E. Jansenc
9059599516SKenneth E. Jansenc.... check for execution mode
9159599516SKenneth E. Jansenc
9259599516SKenneth E. Jansen        if (iexec .eq. 0) then
9359599516SKenneth E. Jansen           lstep = 0
9459599516SKenneth E. Jansen           call restar ('out ',  y  ,ac)
9559599516SKenneth E. Jansen           return
9659599516SKenneth E. Jansen        endif
9759599516SKenneth E. Jansenc
9859599516SKenneth E. Jansenc.... initialize AutoSponge
9959599516SKenneth E. Jansenc
10059599516SKenneth E. Jansen        if(matflg(5,1).ge.4) then ! cool case (sponge)
10159599516SKenneth E. Jansen           call initSponge( y,point2x)
10259599516SKenneth E. Jansen        endif
10359599516SKenneth E. Jansenc
10459599516SKenneth E. Jansenc
10559599516SKenneth E. Jansenc.... adjust BC's to interpolate from file
10659599516SKenneth E. Jansenc
10759599516SKenneth E. Jansen
10859599516SKenneth E. Jansen        inquire(file="inlet.dat",exist=exlog)
10959599516SKenneth E. Jansen        if(exlog) then
11059599516SKenneth E. Jansen           open (unit=654,file="inlet.dat",status="old")
11159599516SKenneth E. Jansen           read(654,*) ninterp,ibcmatch,(interp_mask(j),j=1,ndof)
11259599516SKenneth E. Jansen           do i=1,ninterp
11359599516SKenneth E. Jansen              read(654,*) (bcinterp(i,j),j=1,ndof+1) ! distance to wall+
11459599516SKenneth E. Jansen                        ! ndof but note order of BC's
11559599516SKenneth E. Jansen                        ! p,t,u,v,w,scalars. Also note that must be in
11659599516SKenneth E. Jansen                        ! increasing distance from the wall order.
11759599516SKenneth E. Jansen
11859599516SKenneth E. Jansen           enddo
11959599516SKenneth E. Jansen           do i=1,nshg  ! only correct for linears at this time
12059599516SKenneth E. Jansen              if(mod(iBC(i),1024).eq.ibcmatch) then
12159599516SKenneth E. Jansen                 iupper=0
12259599516SKenneth E. Jansen                 do j=2,ninterp
12359599516SKenneth E. Jansen                    if(bcinterp(j,1).gt.d2wall(i)) then !bound found
12459599516SKenneth E. Jansen                       xi=(d2wall(i)-bcinterp(j-1,1))/
12559599516SKenneth E. Jansen     &                    (bcinterp(j,1)-bcinterp(j-1,1))
12659599516SKenneth E. Jansen                       iupper=j
12759599516SKenneth E. Jansen                       exit
12859599516SKenneth E. Jansen                    endif
12959599516SKenneth E. Jansen                 enddo
13059599516SKenneth E. Jansen                 if(iupper.eq.0) then
13159599516SKenneth E. Jansen                    write(*,*) 'failure in finterp, ynint=',d2wall(i)
13259599516SKenneth E. Jansen                    stop
13359599516SKenneth E. Jansen                 endif
13459599516SKenneth E. Jansen                 if(interp_mask(1).ne.zero) then
13559599516SKenneth E. Jansen                    BC(i,1)=(xi*bcinterp(iupper,2)
13659599516SKenneth E. Jansen     &                +(one-xi)*bcinterp(iupper-1,2))*interp_mask(1)
13759599516SKenneth E. Jansen                 endif
13859599516SKenneth E. Jansen                 if(interp_mask(2).ne.zero) then
13959599516SKenneth E. Jansen                    BC(i,2)=(xi*bcinterp(iupper,3)
14059599516SKenneth E. Jansen     &                +(one-xi)*bcinterp(iupper-1,3))*interp_mask(2)
14159599516SKenneth E. Jansen                 endif
14259599516SKenneth E. Jansen                 if(interp_mask(3).ne.zero) then
14359599516SKenneth E. Jansen                    BC(i,3)=(xi*bcinterp(iupper,4)
14459599516SKenneth E. Jansen     &                +(one-xi)*bcinterp(iupper-1,4))*interp_mask(3)
14559599516SKenneth E. Jansen                 endif
14659599516SKenneth E. Jansen                 if(interp_masK(4).ne.zero) then
14759599516SKenneth E. Jansen                    BC(i,4)=(xi*bcinterp(iupper,5)
14859599516SKenneth E. Jansen     &                +(one-xi)*bcinterp(iupper-1,5))*interp_mask(4)
14959599516SKenneth E. Jansen                 endif
15059599516SKenneth E. Jansen                 if(interp_mask(5).ne.zero) then
15159599516SKenneth E. Jansen                    BC(i,5)=(xi*bcinterp(iupper,6)
15259599516SKenneth E. Jansen     &                +(one-xi)*bcinterp(iupper-1,6))*interp_mask(5)
15359599516SKenneth E. Jansen                 endif
15459599516SKenneth E. Jansen                 if(interp_mask(6).ne.zero) then
15559599516SKenneth E. Jansen                    BC(i,7)=(xi*bcinterp(iupper,7)
15659599516SKenneth E. Jansen     &                +(one-xi)*bcinterp(iupper-1,7))*interp_mask(6)
15759599516SKenneth E. Jansen                 endif
15859599516SKenneth E. Jansen                 if(interp_mask(7).ne.zero) then
15959599516SKenneth E. Jansen                    BC(i,8)=(xi*bcinterp(iupper,8)
16059599516SKenneth E. Jansen     &                +(one-xi)*bcinterp(iupper-1,8))*interp_mask(7)
16159599516SKenneth E. Jansen                 endif
16259599516SKenneth E. Jansen              endif
16359599516SKenneth E. Jansen           enddo
16459599516SKenneth E. Jansen        endif
16559599516SKenneth E. Jansenc$$$$$$$$$$$$$$$$$$$$
16659599516SKenneth E. Jansen
167*513954efSKenneth E. Jansen!======================================================================
168*513954efSKenneth E. Jansen!Modifications for Duct. Need to encapsulate in a function call.
169*513954efSKenneth E. Jansen        !specify an initial eddy viscosity ramp
170*513954efSKenneth E. Jansen        if(isetEVramp .gt. 0) then
171*513954efSKenneth E. Jansen          if(myrank .eq. 0) then
172*513954efSKenneth E. Jansen            write(*,*) "Setting eddy viscosity ramp with:"
173*513954efSKenneth E. Jansen            write(*,*) "  - ramp X min = ", EVrampXmin
174*513954efSKenneth E. Jansen            write(*,*) "  - ramp X max = ", EVrampXmax
175*513954efSKenneth E. Jansen            write(*,*) "  - EV min = ", EVrampMin
176*513954efSKenneth E. Jansen            write(*,*) "  - EV max = ", EVrampMax
177*513954efSKenneth E. Jansen          endif
178*513954efSKenneth E. Jansen
179*513954efSKenneth E. Jansen          x1 = EVrampXmin  !stuff in a shorter variable name to
180*513954efSKenneth E. Jansen          x2 = EVrampXmax  !make the formulas cleaner
181*513954efSKenneth E. Jansen          !Newton Divided differences to generate a polynomial of
182*513954efSKenneth E. Jansen          !the form P(x) = c0 + x*(c1 + x*(c2 + (x - (x2 - x1))*c3))
183*513954efSKenneth E. Jansen          !satisfying P(x1) = EVrampMin, P(x2) = EVrampMax,
184*513954efSKenneth E. Jansen          ! P'(x1) = 0, and P'(x2) = 0
185*513954efSKenneth E. Jansen
186*513954efSKenneth E. Jansen          c0 = EVrampMin
187*513954efSKenneth E. Jansen          c1 = 0            !zero derivative
188*513954efSKenneth E. Jansen          c2 = (EVrampMax - EVrampMin)/(x2 - x1)
189*513954efSKenneth E. Jansen          c3 = 0            !zero derivative
190*513954efSKenneth E. Jansen          c3 = (c3 - c2)/(x2 - x1)
191*513954efSKenneth E. Jansen          c2 = (c2 - c1)/(x2 - x1)
192*513954efSKenneth E. Jansen          c3 = (c3 - c2)/(x2 - x1)
193*513954efSKenneth E. Jansen
194*513954efSKenneth E. Jansen          do nn = 1,nshg
195*513954efSKenneth E. Jansen            if(y(nn,6) .eq. 0) cycle  !don't change wall boundary conditions, should be iTurbWall == 1
196*513954efSKenneth E. Jansen
197*513954efSKenneth E. Jansen            if(point2x(nn,1) .gt. EVrampXmax) then !downstream of the ramp
198*513954efSKenneth E. Jansen              y(nn,6) = EVrampMax
199*513954efSKenneth E. Jansen            elseif(point2x(nn,1) .gt. EVrampXmin) then !and x(:,1) <= EVrampXmax
200*513954efSKenneth E. Jansen
201*513954efSKenneth E. Jansen              !P(x) = c0 + x*(c1 + x*(c2 + (x - (x2 - x1))*c3))
202*513954efSKenneth E. Jansen              !     = c0 + x*(c1 + x*(c2 - (x2 - x1)*c3 + x*c3
203*513954efSKenneth E. Jansen              y(nn,6) = c0                 + (point2x(nn,1) - x1)*(
204*513954efSKenneth E. Jansen     &                  c1                 + (point2x(nn,1) - x1)*(
205*513954efSKenneth E. Jansen     &                 (c2 - (x2 - x1)*c3) + (point2x(nn,1) - x1)*c3))
206*513954efSKenneth E. Jansen            else
207*513954efSKenneth E. Jansen              y(nn,6) = EVrampMin
208*513954efSKenneth E. Jansen            endif
209*513954efSKenneth E. Jansen          enddo
210*513954efSKenneth E. Jansen        endif
211*513954efSKenneth E. Jansen!End modifications for Duct
212*513954efSKenneth E. Jansen!======================================================================
21359599516SKenneth E. Jansenc
21459599516SKenneth E. Jansenc
21559599516SKenneth E. Jansenc.... call the semi-discrete predictor multi-corrector iterative driver
21659599516SKenneth E. Jansenc
21759599516SKenneth E. Jansen        call itrdrv (y,              ac,
21859599516SKenneth E. Jansen     &               uold,           point2x,
21959599516SKenneth E. Jansen     &               iBC,            BC,
22059599516SKenneth E. Jansen     &               point2iper,     point2ilwork,   shp,
22159599516SKenneth E. Jansen     &               shgl,           shpb,           shglb,
22259599516SKenneth E. Jansen     &               point2ifath,    velbar,         point2nsons )
22359599516SKenneth E. Jansenc
22459599516SKenneth E. Jansenc.... return
22559599516SKenneth E. Jansenc
22659599516SKenneth E. Jansenc
22759599516SKenneth E. Jansenc.... stop CPU-timer
22859599516SKenneth E. Jansenc
22959599516SKenneth E. JansenCAD        call timer ('End     ')
23059599516SKenneth E. Jansenc
23159599516SKenneth E. Jansenc.... close echo file
23259599516SKenneth E. Jansenc
23359599516SKenneth E. Jansen
23459599516SKenneth E. Jansen!MR CHANGE
23559599516SKenneth E. Jansen      if (numpe > 1) call MPI_BARRIER(MPI_COMM_WORLD, ierr)
23659599516SKenneth E. Jansen      if(myrank.eq.0)  then
23759599516SKenneth E. Jansen          write(*,*) 'process - before closing iecho'
23859599516SKenneth E. Jansen      endif
23959599516SKenneth E. Jansen!MR CHANGE END
24059599516SKenneth E. Jansen
24159599516SKenneth E. Jansen        close (iecho)
24259599516SKenneth E. Jansen
24359599516SKenneth E. Jansen!MR CHANGE
24459599516SKenneth E. Jansen      if (numpe > 1) call MPI_BARRIER(MPI_COMM_WORLD, ierr)
24559599516SKenneth E. Jansen      if(myrank.eq.0)  then
24659599516SKenneth E. Jansen          write(*,*) 'process - after closing iecho'
24759599516SKenneth E. Jansen      endif
24859599516SKenneth E. Jansen!MR CHANGE END
24959599516SKenneth E. Jansen
25059599516SKenneth E. Jansen
25159599516SKenneth E. Jansenc
25259599516SKenneth E. Jansenc.... end of the program
25359599516SKenneth E. Jansenc
25459599516SKenneth E. JansenCAD        write(6,*) 'Life: ', second(0) - ttim(100)
25559599516SKenneth E. Jansen        deallocate(point2iper)
256fcc77cc2SCameron Smith        if(numpe.gt.1) then
257fcc77cc2SCameron Smith          call Dctypes(point2ilwork(1))
258fcc77cc2SCameron Smith          deallocate(point2ilwork)
259fcc77cc2SCameron Smith        endif
26059599516SKenneth E. Jansen        deallocate(point2x)
26159599516SKenneth E. Jansen        deallocate(point2nsons)
26259599516SKenneth E. Jansen        deallocate(point2ifath)
2636b966dd8SCameron Smith        deallocate(uold)
2646b966dd8SCameron Smith        deallocate(wnrm)
2656b966dd8SCameron Smith        deallocate(otwn)
2666b966dd8SCameron Smith        call finalizeDtN
267bd5c8633SCameron Smith        call clearper
268bd5c8633SCameron Smith        call finalizeNABI
2696b966dd8SCameron Smith
27059599516SKenneth E. Jansen        return
27159599516SKenneth E. Jansen        end
27259599516SKenneth E. Jansen
27359599516SKenneth E. Jansen
274