xref: /phasta/phSolver/compressible/itrdrv.f (revision f0b806cbcd176bc4c068569bf06644e2ab71f1d7)
159599516SKenneth E. Jansen              subroutine itrdrv (y,         ac,   uold, x,
259599516SKenneth E. Jansen     &                   iBC,       BC,
359599516SKenneth E. Jansen     &                   iper,      ilwork,     shp,
459599516SKenneth E. Jansen     &                   shgl,      shpb,       shglb,
559599516SKenneth E. Jansen     &                   ifath,     velbar,     nsons )
659599516SKenneth E. Jansenc
759599516SKenneth E. Jansenc----------------------------------------------------------------------
859599516SKenneth E. Jansenc
959599516SKenneth E. Jansenc This iterative driver is the semi-discrete, predictor multi-corrector
1059599516SKenneth E. Jansenc algorithm. It contains the Hulbert Generalized Alpha method which
1159599516SKenneth E. Jansenc is 2nd order accurate for Rho_inf from 0 to 1.  The method can be
1259599516SKenneth E. Jansenc made  first-order accurate by setting Rho_inf=-1. It uses a
1359599516SKenneth E. Jansenc GMRES iterative solver.
1459599516SKenneth E. Jansenc
1559599516SKenneth E. Jansenc working arrays:
1659599516SKenneth E. Jansenc  y      (nshg,ndof)           : Y variables
1759599516SKenneth E. Jansenc  x      (nshg,nsd)            : node coordinates
1859599516SKenneth E. Jansenc  iBC    (nshg)                : BC codes
1959599516SKenneth E. Jansenc  BC     (nshg,ndofBC)         : BC constraint parameters
2059599516SKenneth E. Jansenc  iper   (nshg)                : periodicity table
2159599516SKenneth E. Jansenc
2259599516SKenneth E. Jansenc shape functions:
2359599516SKenneth E. Jansenc  shp    (nshape,ngauss)        : interior element shape functions
2459599516SKenneth E. Jansenc  shgl   (nsd,nshape,ngauss)    : local shape function gradients
2559599516SKenneth E. Jansenc  shpb   (nshapeb,ngaussb)      : boundary element shape functions
2659599516SKenneth E. Jansenc  shglb  (nsd,nshapeb,ngaussb)  : bdry. elt. shape gradients
2759599516SKenneth E. Jansenc
2859599516SKenneth E. Jansenc Zdenek Johan,  Winter 1991.  (Fortran 90)
2959599516SKenneth E. Jansenc----------------------------------------------------------------------
3059599516SKenneth E. Jansenc
3159599516SKenneth E. Jansen      use pvsQbi        !gives us splag (the spmass at the end of this run
3259599516SKenneth E. Jansen      use specialBC     !gives us itvn
330d32f9a8SKenneth E. Jansen      use timedataC      !allows collection of time series
34513954efSKenneth E. Jansen      use MachControl   !PID to control the inlet velocity.
35513954efSKenneth E. Jansen      use blowerControl !gives us BC_enable
3659599516SKenneth E. Jansen      use turbSA
37513954efSKenneth E. Jansen      use wallData
38513954efSKenneth E. Jansen      use fncorpmod
3959599516SKenneth E. Jansen
4059599516SKenneth E. Jansen        include "common.h"
4159599516SKenneth E. Jansen        include "mpif.h"
4259599516SKenneth E. Jansen        include "auxmpi.h"
4359599516SKenneth E. Jansen
4459599516SKenneth E. Jansenc
4559599516SKenneth E. Jansen        dimension y(nshg,ndof),            ac(nshg,ndof),
4659599516SKenneth E. Jansen     &            yold(nshg,ndof),         acold(nshg,ndof),
4759599516SKenneth E. Jansen     &            x(numnp,nsd),            iBC(nshg),
4859599516SKenneth E. Jansen     &            BC(nshg,ndofBC),         ilwork(nlwork),
4959599516SKenneth E. Jansen     &            iper(nshg),              uold(nshg,nsd)
50*f0b806cbSKenneth E. Jansen
51*f0b806cbSKenneth E. Jansen
5259599516SKenneth E. Jansenc
53513954efSKenneth E. Jansen        dimension res(nshg,nflow),
5459599516SKenneth E. Jansen     &            rest(nshg),              solinc(nshg,ndof)
5559599516SKenneth E. Jansenc
5659599516SKenneth E. Jansen        dimension shp(MAXTOP,maxsh,MAXQPT),
5759599516SKenneth E. Jansen     &            shgl(MAXTOP,nsd,maxsh,MAXQPT),
5859599516SKenneth E. Jansen     &            shpb(MAXTOP,maxsh,MAXQPT),
5959599516SKenneth E. Jansen     &            shglb(MAXTOP,nsd,maxsh,MAXQPT)
6059599516SKenneth E. Jansen        real*8   almit, alfit, gamit
61*f0b806cbSKenneth E. Jansen
62*f0b806cbSKenneth E. Jansen
63*f0b806cbSKenneth E. Jansen
6459599516SKenneth E. Jansen        dimension ifath(numnp),    velbar(nfath,ndof),  nsons(nfath)
6559599516SKenneth E. Jansen        real*8 rerr(nshg,10),ybar(nshg,ndof+8) ! 8 is for avg. of square as uu, vv, ww, pp, TT, uv, uw, and vw
66513954efSKenneth E. Jansen        real*8, allocatable, dimension(:,:) :: vortG
67513954efSKenneth E. Jansen        real*8, allocatable, dimension(:,:,:) :: BDiag
68513954efSKenneth E. Jansen
69513954efSKenneth E. Jansen!       integer, allocatable, dimension(:) :: iv_rank  !used with MRasquin's version of probe points
70513954efSKenneth E. Jansen!       integer :: iv_rankpernode, iv_totnodes, iv_totcores
71513954efSKenneth E. Jansen!       integer :: iv_node,        iv_core,     iv_thread
72513954efSKenneth E. Jansen
7359599516SKenneth E. Jansen! assuming three profiles to control (inlet, bottom FC and top FC)
7459599516SKenneth E. Jansen! store velocity profile set via BC
7559599516SKenneth E. Jansen        real*8 vbc_prof(nshg,3)
76513954efSKenneth E. Jansen        character(len=60) fvarts
7759599516SKenneth E. Jansen        integer ifuncs(6), iarray(10)
785124a526SKenneth E. Jansen        integer BCdtKW, tsBase
795124a526SKenneth E. Jansen
8059599516SKenneth E. Jansen        real*8 elDw(numel) ! element average of DES d variable
81e72fd12cSBen Matthews
82e72fd12cSBen Matthews        real*8, allocatable, dimension(:,:) :: HBrg
83e72fd12cSBen Matthews        real*8, allocatable, dimension(:) :: eBrg
84e72fd12cSBen Matthews        real*8, allocatable, dimension(:) :: yBrg
85e72fd12cSBen Matthews        real*8, allocatable, dimension(:) :: Rcos, Rsin
8659599516SKenneth E. Jansenc
8759599516SKenneth E. Jansenc  Here are the data structures for sparse matrix GMRES
8859599516SKenneth E. Jansenc
8959599516SKenneth E. Jansen        integer, allocatable, dimension(:,:) :: rowp
9059599516SKenneth E. Jansen        integer, allocatable, dimension(:) :: colm
9159599516SKenneth E. Jansen        real*8, allocatable, dimension(:,:) :: lhsK
9259599516SKenneth E. Jansen        real*8, allocatable, dimension(:,:) :: EGmass
9359599516SKenneth E. Jansen        real*8, allocatable, dimension(:,:) :: EGmasst
9459599516SKenneth E. Jansen
95513954efSKenneth E. Jansen        integer iTurbWall(nshg)
96513954efSKenneth E. Jansen        real*8 yInlet(3), yInletg(3)
97513954efSKenneth E. Jansen        integer imapped, imapInlet(nshg)  !for now, used for setting Blower conditions
98513954efSKenneth E. Jansen!        real*8 M_th, M_tc, M_tt
99513954efSKenneth E. Jansen!        logical  exMc
100513954efSKenneth E. Jansen!        real*8 vBC, vBCg
101513954efSKenneth E. Jansen        real*8 vortmax, vortmaxg
10259599516SKenneth E. Jansen
103513954efSKenneth E. Jansen       iprec=0 !PETSc - Disable PHASTA's BDiag. TODO: Preprocssor Switch
10459599516SKenneth E. Jansen
105513954efSKenneth E. Jansen       call findTurbWall(iTurbWall)
10659599516SKenneth E. Jansen
107513954efSKenneth E. Jansen!-------
108513954efSKenneth E. Jansen! SETUP
109513954efSKenneth E. Jansen!-------
11059599516SKenneth E. Jansen
111513954efSKenneth E. Jansen       !HACK for debugging suction
112513954efSKenneth E. Jansen!       call Write_Debug(myrank, 'wallNormal'//char(0),
113513954efSKenneth E. Jansen!     &                          'wnorm'//char(0), wnorm,
114513954efSKenneth E. Jansen!     &                          'd', nshg, 3, lstep)
115513954efSKenneth E. Jansen
116513954efSKenneth E. Jansen       !Probe Point Setup
117513954efSKenneth E. Jansen       call initProbePoints()
118513954efSKenneth E. Jansen       if(exts) then  !exts is set in initProbePoints
119513954efSKenneth E. Jansen         write(fvarts, "('./varts/varts.', I0, '.dat')") lstep
120513954efSKenneth E. Jansen         fvarts = trim(fvarts)
12159599516SKenneth E. Jansen
12259599516SKenneth E. Jansen         if(myrank .eq. master) then
123513954efSKenneth E. Jansen           call TD_writeHeader(fvarts)
124513954efSKenneth E. Jansen         endif
12559599516SKenneth E. Jansen       endif
12659599516SKenneth E. Jansen
127513954efSKenneth E. Jansen       !Mach Control Setup
128513954efSKenneth E. Jansen       call MC_init(Delt, lstep, BC)
129513954efSKenneth E. Jansen       exMC = exMC .and. exts   !If probe points aren't available, turn
130513954efSKenneth E. Jansen                                !the Mach Control off
131513954efSKenneth E. Jansen       if(exMC) then
132513954efSKenneth E. Jansen         call MC_applyBC(BC)
133513954efSKenneth E. Jansen         call MC_printState()
13459599516SKenneth E. Jansen       endif
13559599516SKenneth E. Jansen
136513954efSKenneth E. Jansen
13759599516SKenneth E. Jansenc
13859599516SKenneth E. Jansenc.... open history and aerodynamic forces files
13959599516SKenneth E. Jansenc
14059599516SKenneth E. Jansen        if (myrank .eq. master) then
14159599516SKenneth E. Jansen          open (unit=ihist,  file=fhist,  status='unknown')
14259599516SKenneth E. Jansen          open (unit=iforce, file=fforce, status='unknown')
14359599516SKenneth E. Jansen        endif
14459599516SKenneth E. Jansenc
14559599516SKenneth E. Jansenc
14659599516SKenneth E. Jansenc.... initialize
14759599516SKenneth E. Jansenc
14859599516SKenneth E. Jansen        ifuncs  = 0                      ! func. evaluation counter
14959599516SKenneth E. Jansen        istep  = 0
15059599516SKenneth E. Jansen        ntotGM = 0                      ! number of GMRES iterations
15159599516SKenneth E. Jansen        time   = 0
15259599516SKenneth E. Jansen        yold   = y
15359599516SKenneth E. Jansen        acold  = ac
154513954efSKenneth E. Jansen
155513954efSKenneth E. Jansen!Blower Setup
156513954efSKenneth E. Jansen       call BC_init(Delt, lstep, BC)  !Note: sets BC_enable
157513954efSKenneth E. Jansen! fix the yold values to the reset BC
158513954efSKenneth E. Jansen      if(BC_enable) call itrBC (yold,  ac,  iBC,  BC,  iper, ilwork)
159513954efSKenneth E. Jansen! without the above, second solve of first steps is fouled up
160513954efSKenneth E. Jansen!
161e72fd12cSBen Matthews        allocate(HBrg(Kspace+1,Kspace))
162e72fd12cSBen Matthews        allocate(eBrg(Kspace+1))
163e72fd12cSBen Matthews        allocate(yBrg(Kspace))
164e72fd12cSBen Matthews        allocate(Rcos(Kspace))
165e72fd12cSBen Matthews        allocate(Rsin(Kspace))
166e72fd12cSBen Matthews
16759599516SKenneth E. Jansen        if (mod(impl(1),100)/10 .eq. 1) then
16859599516SKenneth E. Jansenc
16959599516SKenneth E. Jansenc     generate the sparse data fill vectors
17059599516SKenneth E. Jansenc
17159599516SKenneth E. Jansen           allocate  (rowp(nshg,nnz))
17259599516SKenneth E. Jansen           allocate  (colm(nshg+1))
17359599516SKenneth E. Jansen           call genadj(colm, rowp, icnt ) ! preprocess the adjacency list
17459599516SKenneth E. Jansen
17559599516SKenneth E. Jansen           nnz_tot=icnt         ! this is exactly the number of non-zero
17659599516SKenneth E. Jansen                                ! blocks on this proc
177513954efSKenneth E. Jansen           if(usingpetsc.eq.1) then
178513954efSKenneth E. Jansen             allocate (BDiag(1,1,1))
179513954efSKenneth E. Jansen           else
180513954efSKenneth E. Jansen             allocate (BDiag(nshg,nflow,nflow))
18159599516SKenneth E. Jansen             allocate (lhsK(nflow*nflow,nnz_tot))
18259599516SKenneth E. Jansen           endif
183513954efSKenneth E. Jansen        endif
18459599516SKenneth E. Jansen        if (mod(impl(1),100)/10 .eq. 3) then
18559599516SKenneth E. Jansenc
18659599516SKenneth E. Jansenc     generate the ebe data fill vectors
18759599516SKenneth E. Jansenc
18859599516SKenneth E. Jansen           nedof=nflow*nshape
18959599516SKenneth E. Jansen           allocate  (EGmass(numel,nedof*nedof))
19059599516SKenneth E. Jansen        endif
19159599516SKenneth E. Jansen
19259599516SKenneth E. Jansenc..........................................
19359599516SKenneth E. Jansen        rerr = zero
19459599516SKenneth E. Jansen        ybar(:,1:ndof) = y(:,1:ndof)
19559599516SKenneth E. Jansen        do idof=1,5
19659599516SKenneth E. Jansen           ybar(:,ndof+idof) = y(:,idof)*y(:,idof)
19759599516SKenneth E. Jansen        enddo
19859599516SKenneth E. Jansen        ybar(:,ndof+6) = y(:,1)*y(:,2)
19959599516SKenneth E. Jansen        ybar(:,ndof+7) = y(:,1)*y(:,3)
20059599516SKenneth E. Jansen        ybar(:,ndof+8) = y(:,2)*y(:,3)
20159599516SKenneth E. Jansenc.........................................
20259599516SKenneth E. Jansen
203513954efSKenneth E. Jansen!  change the freestream and inflow eddy viscosity to reflect different
204513954efSKenneth E. Jansen!  levels of freestream turbulence
205513954efSKenneth E. Jansen!
206513954efSKenneth E. Jansen! First override boundary condition values
207513954efSKenneth E. Jansen!  USES imapInlet from Mach Control so if that gets conditionaled away
208513954efSKenneth E. Jansen!  it has to know if it is needed here
209513954efSKenneth E. Jansen!
210513954efSKenneth E. Jansen      if(isetEV_IC_BC.eq.1) then
211513954efSKenneth E. Jansen        allocate(vortG(nshg, 4))
212513954efSKenneth E. Jansen        call vortGLB(yold, x, shp, shgl, ilwork, vortG)
213513954efSKenneth E. Jansen        vortmax=maxval(abs(vortG(:,4)))  !column 4 is the magnitude of the shear tensor - should actually probably be calld shearmax instead of vortmax
214513954efSKenneth E. Jansen        write(*,*) "vortmax = ", vortmax
21559599516SKenneth E. Jansen
216513954efSKenneth E. Jansen        !Find the maximum shear in the simulation
217513954efSKenneth E. Jansen        if(numpe.gt.1) then
218513954efSKenneth E. Jansen           call MPI_ALLREDUCE(vortmax, vortmaxg, 1,
219513954efSKenneth E. Jansen     &          MPI_DOUBLE_PRECISION, MPI_MAX, MPI_COMM_WORLD, ierr )
220513954efSKenneth E. Jansen           vortmax = vortmaxg
221513954efSKenneth E. Jansen        endif
22259599516SKenneth E. Jansen
223513954efSKenneth E. Jansen        !Apply eddy viscosity at the inlet
224513954efSKenneth E. Jansen        do i=1,icountInlet ! for now only coded to catch primary inlet, not blower
225513954efSKenneth E. Jansen           BC(imapInlet(i),7)=evis_IC_BC
226513954efSKenneth E. Jansen        enddo
227513954efSKenneth E. Jansen
228513954efSKenneth E. Jansen        !Apply eddy viscosity through the quasi-inviscid portion of the domain
229513954efSKenneth E. Jansen        do i = 1,nshg
230513954efSKenneth E. Jansen          if(abs(vortG(i,4)).ge.vortmax*0.01) yold(i,6)=evis_IC_BC
231513954efSKenneth E. Jansen        enddo
232513954efSKenneth E. Jansen        isclr=1 ! fix scalar
233513954efSKenneth E. Jansen        call itrBCsclr(yold,ac,iBC,BC,iper,ilwork)
234513954efSKenneth E. Jansen
235513954efSKenneth E. Jansen        deallocate(vortG)
236513954efSKenneth E. Jansen      endif
23759599516SKenneth E. Jansenc
23859599516SKenneth E. Jansenc.... loop through the time sequences
23959599516SKenneth E. Jansenc
24059599516SKenneth E. Jansen        do 3000 itsq = 1, ntseq
24159599516SKenneth E. Jansen        itseq = itsq
24259599516SKenneth E. Jansenc
24359599516SKenneth E. Jansenc.... set up the current parameters
24459599516SKenneth E. Jansenc
24559599516SKenneth E. Jansen        nstp   = nstep(itseq)
24659599516SKenneth E. Jansen        nitr   = niter(itseq)
24759599516SKenneth E. Jansen        LCtime = loctim(itseq)
24859599516SKenneth E. Jansenc
24959599516SKenneth E. Jansen        call itrSetup ( y,  acold)
25059599516SKenneth E. Jansen        isclr=0
25159599516SKenneth E. Jansen
25259599516SKenneth E. Jansen        niter(itseq)=0          ! count number of flow solves in a step
25359599516SKenneth E. Jansen                                !  (# of iterations)
25459599516SKenneth E. Jansen        do i=1,seqsize
25559599516SKenneth E. Jansen           if(stepseq(i).eq.0) niter(itseq)=niter(itseq)+1
25659599516SKenneth E. Jansen        enddo
25759599516SKenneth E. Jansen        nitr = niter(itseq)
25859599516SKenneth E. Jansenc
25959599516SKenneth E. Jansenc.... determine how many scalar equations we are going to need to solve
26059599516SKenneth E. Jansenc
26159599516SKenneth E. Jansen        nsclrsol=nsclr          ! total number of scalars solved. At
26259599516SKenneth E. Jansen                                ! some point we probably want to create
26359599516SKenneth E. Jansen                                ! a map, considering stepseq(), to find
26459599516SKenneth E. Jansen                                ! what is actually solved and only
26559599516SKenneth E. Jansen                                ! dimension EGmasst to the appropriate
26659599516SKenneth E. Jansen                                ! size.
26759599516SKenneth E. Jansen
26859599516SKenneth E. Jansen        if(nsclrsol.gt.0)allocate  (EGmasst(numel*nshape*nshape
26959599516SKenneth E. Jansen     &                              ,nsclrsol))
27059599516SKenneth E. Jansen
27159599516SKenneth E. Jansenc
27259599516SKenneth E. Jansenc.... loop through the time steps
27359599516SKenneth E. Jansenc
27459599516SKenneth E. Jansen        ttim(1) = REAL(secs(0.0)) / 100.
27559599516SKenneth E. Jansen        ttim(2) = secs(0.0)
27659599516SKenneth E. Jansen
27759599516SKenneth E. Jansenc        tcorecp1 = REAL(secs(0.0)) / 100.
27859599516SKenneth E. Jansenc        tcorewc1 = secs(0.0)
27959599516SKenneth E. Jansen        if (numpe > 1) call MPI_BARRIER(MPI_COMM_WORLD, ierr)
280513954efSKenneth E. Jansen        if(myrank.eq.master)  then
28159599516SKenneth E. Jansen           tcorecp1 = TMRC()
28259599516SKenneth E. Jansen        endif
28359599516SKenneth E. Jansen
28459599516SKenneth E. Jansen        rmub=datmat(1,2,1)
28559599516SKenneth E. Jansen        if(rmutarget.gt.0) then
28659599516SKenneth E. Jansen           rmue=rmutarget
28759599516SKenneth E. Jansen           xmulfact=(rmue/rmub)**(1.0/nstp)
288513954efSKenneth E. Jansen           if(myrank.eq.master) then
28959599516SKenneth E. Jansen              write(*,*) 'viscosity will by multiplied by ', xmulfact
29059599516SKenneth E. Jansen              write(*,*) 'to bring it from ', rmub,' down to ', rmue
29159599516SKenneth E. Jansen           endif
29259599516SKenneth E. Jansen           datmat(1,2,1)=datmat(1,2,1)/xmulfact ! make first step right
29359599516SKenneth E. Jansen        else
29459599516SKenneth E. Jansen           rmue=datmat(1,2,1)   ! keep constant
29559599516SKenneth E. Jansen           xmulfact=one
29659599516SKenneth E. Jansen        endif
29759599516SKenneth E. Jansen        if(iramp.eq.1) then
29859599516SKenneth E. Jansen                call initBCprofileScale(vbc_prof,BC,yold,x)
29959599516SKenneth E. Jansen! fix the yold values to the reset BC
30059599516SKenneth E. Jansen                call itrBC (yold,  ac,  iBC,  BC,  iper, ilwork)
30159599516SKenneth E. Jansen                isclr=1 ! fix scalar
30259599516SKenneth E. Jansen                call itrBCsclr(yold,ac,iBC,BC,iper,ilwork)
30359599516SKenneth E. Jansen        endif
3045124a526SKenneth E. Jansenc  Time Varying BCs------------------------------------(Kyle W 6-6-13)
3055124a526SKenneth E. Jansenc        BCdtKW=0
3065124a526SKenneth E. Jansen        if(BCdtKW.gt.0) then
3075124a526SKenneth E. Jansen           call BCprofileInitKW(PresBase,VelBase,BC)
3085124a526SKenneth E. Jansen        endif
3095124a526SKenneth E. Jansenc  Time Varying BCs------------------------------------(Kyle W 6-6-13)
31059599516SKenneth E. Jansen
31159599516SKenneth E. Jansen867     continue
312513954efSKenneth E. Jansen
313513954efSKenneth E. Jansen
314513954efSKenneth E. Jansenc============ Start the loop of time steps============================c
315513954efSKenneth E. Jansen!        edamp2=.99
316513954efSKenneth E. Jansen!        edamp3=0.05
317513954efSKenneth E. Jansen        deltaInlInv=one/(0.125*0.0254)
31859599516SKenneth E. Jansen        do 2000 istp = 1, nstp
319*f0b806cbSKenneth E. Jansen           rerr=zero !extreme limit of 1 step window or error stats....later a variable
320513954efSKenneth E. Jansen
321*f0b806cbSKenneth E. Jansen!           if (myrank.eq.master) write(*,*) 'Time step of current run',
322*f0b806cbSKenneth E. Jansen!     &                                    istp
3235124a526SKenneth E. Jansen
3245124a526SKenneth E. Jansenc  Time Varying BCs------------------------------------(Kyle W 6-6-13)
3255124a526SKenneth E. Jansen           if(BCdtKW.gt.0) then
3265124a526SKenneth E. Jansen              call BCprofileScaleKW(PresBase,VelBase,BC,yold)
3275124a526SKenneth E. Jansen           endif
3285124a526SKenneth E. Jansenc  Time Varying BCs------------------------------------(Kyle W 6-6-13)
329513954efSKenneth E. Jansen
33059599516SKenneth E. Jansen           if(iramp.eq.1)
33159599516SKenneth E. Jansen     &        call BCprofileScale(vbc_prof,BC,yold)
33259599516SKenneth E. Jansen
33359599516SKenneth E. Jansenc           call rerun_check(stopjob)
33459599516SKenneth E. Jansencc          if(stopjob.ne.0) goto 2001
33559599516SKenneth E. Jansenc
33659599516SKenneth E. Jansenc Decay of scalars
33759599516SKenneth E. Jansenc
33859599516SKenneth E. Jansen           if(nsclr.gt.0 .and. tdecay.ne.1) then
33959599516SKenneth E. Jansen              yold(:,6:ndof)=y(:,6:ndof)*tdecay
34059599516SKenneth E. Jansen              BC(:,7:6+nsclr)= BC(:,7:6+nsclr)*tdecay
34159599516SKenneth E. Jansen           endif
34259599516SKenneth E. Jansen
34359599516SKenneth E. Jansen           if(nosource.eq.1) BC(:,7:6+nsclr)= BC(:,7:6+nsclr)*0.8
34459599516SKenneth E. Jansen
34559599516SKenneth E. Jansen
34659599516SKenneth E. Jansenc           xi=istp*one/nstp
34759599516SKenneth E. Jansenc           datmat(1,2,1)=rmub*(1.0-xi)+xi*rmue
34859599516SKenneth E. Jansen           datmat(1,2,1)=xmulfact*datmat(1,2,1)
34959599516SKenneth E. Jansen
35059599516SKenneth E. Jansenc.... if we have time varying boundary conditions update the values of BC.
35159599516SKenneth E. Jansenc     these will be for time step n+1 so use lstep+1
35259599516SKenneth E. Jansenc
35359599516SKenneth E. Jansen           if(itvn.gt.0) call BCint((lstep+1)*Delt(1), shp, shgl,
35459599516SKenneth E. Jansen     &                              shpb, shglb, x, BC, iBC)
35559599516SKenneth E. Jansen
35659599516SKenneth E. Jansen            if(iLES.gt.0) then
35759599516SKenneth E. Jansenc
35859599516SKenneth E. Jansenc.... get dynamic model coefficient
35959599516SKenneth E. Jansenc
36059599516SKenneth E. Jansen            ilesmod=iLES/10
36159599516SKenneth E. Jansenc
36259599516SKenneth E. Jansenc digit bit set filter rule, 10 bit set model
36359599516SKenneth E. Jansenc
36459599516SKenneth E. Jansen            if (ilesmod.eq.0) then ! 0 < iLES < 10 => dyn. model calculated
36559599516SKenneth E. Jansen                                   ! at nodes based on discrete filtering
36659599516SKenneth E. Jansen               call getdmc (yold,       shgl,      shp,
36759599516SKenneth E. Jansen     &                      iper,       ilwork,    nsons,
36859599516SKenneth E. Jansen     &                      ifath,      x)
36959599516SKenneth E. Jansen            endif
37059599516SKenneth E. Jansen            if (ilesmod .eq. 1) then ! 10 < iLES < 20 => dynamic-mixed
37159599516SKenneth E. Jansen                                     ! at nodes based on discrete filtering
37259599516SKenneth E. Jansen               call bardmc (yold,       shgl,      shp,
37359599516SKenneth E. Jansen     &                      iper,       ilwork,
37459599516SKenneth E. Jansen     &                      nsons,      ifath,     x)
37559599516SKenneth E. Jansen            endif
37659599516SKenneth E. Jansen            if (ilesmod .eq. 2) then ! 20 < iLES < 30 => dynamic at quad
37759599516SKenneth E. Jansen                                     ! pts based on lumped projection filt.
37859599516SKenneth E. Jansen               call projdmc (yold,       shgl,      shp,
37959599516SKenneth E. Jansen     &                       iper,       ilwork,    x)
38059599516SKenneth E. Jansen            endif
38159599516SKenneth E. Jansenc
382513954efSKenneth E. Jansen           endif ! endif of iLES
38359599516SKenneth E. Jansen
38459599516SKenneth E. Jansen
38559599516SKenneth E. Jansenc
38659599516SKenneth E. Jansenc.... set traction BCs for modeled walls
38759599516SKenneth E. Jansenc
38859599516SKenneth E. Jansen            if (itwmod.ne.0) then   !wallfn check
38959599516SKenneth E. Jansen               call asbwmod(yold,   acold,   x,      BC,     iBC,
39059599516SKenneth E. Jansen     &                      iper,   ilwork,  ifath,  velbar)
39159599516SKenneth E. Jansen            endif
39259599516SKenneth E. Jansenc
39359599516SKenneth E. Jansenc.... -----------------------> predictor phase <-----------------------
39459599516SKenneth E. Jansenc
39559599516SKenneth E. Jansen            call itrPredict(   yold,    acold,    y,   ac )
39659599516SKenneth E. Jansen            call itrBC (y,  ac,  iBC,  BC,  iper, ilwork)
39759599516SKenneth E. Jansen            isclr = zero
39859599516SKenneth E. Jansen            if (nsclr.gt.zero) then
39959599516SKenneth E. Jansen            do isclr=1,nsclr
40059599516SKenneth E. Jansen               call itrBCSclr (y, ac,  iBC, BC, iper, ilwork)
40159599516SKenneth E. Jansen            enddo
40259599516SKenneth E. Jansen            endif
40359599516SKenneth E. Jansenc
40459599516SKenneth E. Jansenc.... --------------------> multi-corrector phase <--------------------
40559599516SKenneth E. Jansenc
40659599516SKenneth E. Jansen           iter=0
40759599516SKenneth E. Jansen            ilss=0  ! this is a switch thrown on first solve of LS redistance
40859599516SKenneth E. Jansenc
40959599516SKenneth E. JansencHACK to make it keep solving RANS until tolerance is reached
41059599516SKenneth E. Jansenc
41159599516SKenneth E. Jansen       istop=0
412513954efSKenneth E. Jansen! blocking extra RANS steps for now       iMoreRANS=0
413513954efSKenneth E. Jansen       iMoreRANS=6
41459599516SKenneth E. Jansenc
41559599516SKenneth E. Jansenc find the the RANS portion of the sequence
41659599516SKenneth E. Jansenc
41759599516SKenneth E. Jansen       do istepc=1,seqsize
41859599516SKenneth E. Jansen          if(stepseq(istepc).eq.10) iLastRANS=istepc
41959599516SKenneth E. Jansen       enddo
42059599516SKenneth E. Jansen
42159599516SKenneth E. Jansen       iseqStart=1
42259599516SKenneth E. Jansen9876   continue
42359599516SKenneth E. Jansenc
42459599516SKenneth E. Jansen            do istepc=iseqStart,seqsize
42559599516SKenneth E. Jansen               icode=stepseq(istepc)
42659599516SKenneth E. Jansen               if(mod(icode,10).eq.0) then ! this is a solve
42759599516SKenneth E. Jansen                  isolve=icode/10
42859599516SKenneth E. Jansen                  if(isolve.eq.0) then   ! flow solve (encoded as 0)
42959599516SKenneth E. Jansenc
43059599516SKenneth E. Jansen                     etol=epstol(1)
43159599516SKenneth E. Jansen                     iter   = iter+1
43259599516SKenneth E. Jansen                     ifuncs(1)  = ifuncs(1) + 1
43359599516SKenneth E. Jansenc
43459599516SKenneth E. Jansenc.... reset the aerodynamic forces
43559599516SKenneth E. Jansenc
43659599516SKenneth E. Jansen                     Force(1) = zero
43759599516SKenneth E. Jansen                     Force(2) = zero
43859599516SKenneth E. Jansen                     Force(3) = zero
43959599516SKenneth E. Jansen                     HFlux    = zero
44059599516SKenneth E. Jansenc
44159599516SKenneth E. Jansenc.... form the element data and solve the matrix problem
44259599516SKenneth E. Jansenc
44359599516SKenneth E. Jansenc.... explicit solver
44459599516SKenneth E. Jansenc
44559599516SKenneth E. Jansen                     if (impl(itseq) .eq. 0) then
44659599516SKenneth E. Jansen                        if (myrank .eq. master)
44759599516SKenneth E. Jansen     &                       call error('itrdrv  ','impl ',impl(itseq))
44859599516SKenneth E. Jansen                     endif
44959599516SKenneth E. Jansen                     if (mod(impl(1),100)/10 .eq. 1) then  ! sparse solve
45059599516SKenneth E. Jansenc
45159599516SKenneth E. Jansenc.... preconditioned sparse matrix GMRES solver
45259599516SKenneth E. Jansenc
45359599516SKenneth E. Jansen                        lhs = 1 - min(1,mod(ifuncs(1)-1,LHSupd(1)))
45459599516SKenneth E. Jansen                        iprec=lhs
45559599516SKenneth E. Jansen                        nedof = nflow*nshape
45659599516SKenneth E. Jansenc                        write(*,*) 'lhs=',lhs
457513954efSKenneth E. Jansen                    if(usingpetsc.eq.1) then
45817860365SKenneth E. Jansen#if (HAVE_PETSC)
459513954efSKenneth E. Jansen               call SolGMRp (y,             ac,            yold,
460513954efSKenneth E. Jansen     &                       x,
461513954efSKenneth E. Jansen     &                       iBC,           BC,
462513954efSKenneth E. Jansen     &                       colm,          rowp,          lhsk,
463513954efSKenneth E. Jansen     &                       res,
464513954efSKenneth E. Jansen     &                       BDiag,
465513954efSKenneth E. Jansen     &                       iper,          ilwork,
466513954efSKenneth E. Jansen     &                       shp,           shgl,
467513954efSKenneth E. Jansen     &                       shpb,          shglb,         solinc,
468513954efSKenneth E. Jansen     &                       rerr,          fncorp )
46917860365SKenneth E. Jansen#else
47079f1763eSKenneth E. Jansen                     if(myrank.eq.0) write(*,*) 'exiting because run time input asked for PETSc, not linked in exec'
47117860365SKenneth E. Jansen                     call error('itrdrv  ','noPETSc',usingpetsc)
47217860365SKenneth E. Jansen#endif
473513954efSKenneth E. Jansen                     else
47459599516SKenneth E. Jansen                      call SolGMRs (y,             ac,            yold,
47559599516SKenneth E. Jansen     &                       acold,         x,
47659599516SKenneth E. Jansen     &                       iBC,           BC,
47759599516SKenneth E. Jansen     &                       colm,          rowp,          lhsk,
47859599516SKenneth E. Jansen     &                       res,
479e72fd12cSBen Matthews     &                       BDiag,         hBrg,          eBrg,
480e72fd12cSBen Matthews     &                       yBrg,          Rcos,          Rsin,
48159599516SKenneth E. Jansen     &                       iper,          ilwork,
48259599516SKenneth E. Jansen     &                       shp,           shgl,
48359599516SKenneth E. Jansen     &                       shpb,          shglb,         solinc,
48459599516SKenneth E. Jansen     &                       rerr)
485513954efSKenneth E. Jansen                    endif
48659599516SKenneth E. Jansen                      else if (mod(impl(1),100)/10 .eq. 2) then ! mfg solve
48759599516SKenneth E. Jansenc
48859599516SKenneth E. Jansenc.... preconditioned matrix-free GMRES solver
48959599516SKenneth E. Jansenc
49059599516SKenneth E. Jansen                        lhs=0
49159599516SKenneth E. Jansen                        iprec = 1 - min(1,mod(ifuncs(1)-1,LHSupd(1)))
49259599516SKenneth E. Jansen                        nedof = 0
49359599516SKenneth E. Jansen                        call SolMFG (y,             ac,            yold,
49459599516SKenneth E. Jansen     &                       acold,         x,
49559599516SKenneth E. Jansen     &                       iBC,           BC,
49659599516SKenneth E. Jansen     &                       res,
497e72fd12cSBen Matthews     &                       BDiag,         HBrg,          eBrg,
498e72fd12cSBen Matthews     &                       yBrg,          Rcos,          Rsin,
49959599516SKenneth E. Jansen     &                       iper,          ilwork,
50059599516SKenneth E. Jansen     &                       shp,           shgl,
50159599516SKenneth E. Jansen     &                       shpb,          shglb,         solinc,
50259599516SKenneth E. Jansen     &                       rerr)
50359599516SKenneth E. Jansenc
50459599516SKenneth E. Jansen                     else if (mod(impl(1),100)/10 .eq. 3) then ! ebe solve
50559599516SKenneth E. Jansenc.... preconditioned ebe matrix GMRES solver
50659599516SKenneth E. Jansenc
50759599516SKenneth E. Jansen
50859599516SKenneth E. Jansen                        lhs = 1 - min(1,mod(ifuncs(1)-1,LHSupd(1)))
50959599516SKenneth E. Jansen                        iprec = lhs
51059599516SKenneth E. Jansen                        nedof = nflow*nshape
51159599516SKenneth E. Jansenc                        write(*,*) 'lhs=',lhs
51259599516SKenneth E. Jansen                      call SolGMRe (y,             ac,            yold,
51359599516SKenneth E. Jansen     &                       acold,         x,
51459599516SKenneth E. Jansen     &                       iBC,           BC,
51559599516SKenneth E. Jansen     &                       EGmass,        res,
516e72fd12cSBen Matthews     &                       BDiag,         HBrg,          eBrg,
517e72fd12cSBen Matthews     &                       yBrg,          Rcos,          Rsin,
51859599516SKenneth E. Jansen     &                       iper,          ilwork,
51959599516SKenneth E. Jansen     &                       shp,           shgl,
52059599516SKenneth E. Jansen     &                       shpb,          shglb,         solinc,
52159599516SKenneth E. Jansen     &                       rerr)
52259599516SKenneth E. Jansen                     endif
52359599516SKenneth E. Jansenc
52459599516SKenneth E. Jansen                else          ! solve a scalar  (encoded at isclr*10)
52559599516SKenneth E. Jansen                     isclr=isolve
526513954efSKenneth E. Jansen                     etol=epstol(isclr+2) ! note that for both epstol and LHSupd 1 is flow 2 temp isclr+2 for scalars
527513954efSKenneth E. Jansen                     ifuncs(isclr+2)  = ifuncs(isclr+2) + 1
52859599516SKenneth E. Jansen                     if((iLSet.eq.2).and.(ilss.eq.0)
52959599516SKenneth E. Jansen     &                    .and.(isclr.eq.2)) then
53059599516SKenneth E. Jansen                        ilss=1  ! throw switch (once per step)
53159599516SKenneth E. Jansenc
53259599516SKenneth E. Jansenc... copy the first scalar at t=n+1 into the second scalar of the
53359599516SKenneth E. Jansenc... level sets
53459599516SKenneth E. Jansenc
53559599516SKenneth E. Jansen                     y(:,6)    = yold(:,6)  + (y(:,6)-yold(:,6))/alfi
53659599516SKenneth E. Jansen                     y(:,7)    =  y(:,6)
53759599516SKenneth E. Jansen                     yold(:,7) = y(:,7)
53859599516SKenneth E. Jansen                     ac(:,7)   = zero
53959599516SKenneth E. Jansenc
54059599516SKenneth E. Jansen                     call itrBCSclr (y, ac,  iBC, BC, iper, ilwork)
54159599516SKenneth E. Jansenc
54259599516SKenneth E. Jansenc....store the flow alpha, gamma parameter values and assigm them the
54359599516SKenneth E. Jansenc....Backward Euler parameters to solve the second levelset scalar
54459599516SKenneth E. Jansenc
54559599516SKenneth E. Jansen                        alfit=alfi
54659599516SKenneth E. Jansen                        gamit=gami
54759599516SKenneth E. Jansen                        almit=almi
54859599516SKenneth E. Jansen                        alfi = 1
54959599516SKenneth E. Jansen                        gami = 1
55059599516SKenneth E. Jansen                        almi = 1
55159599516SKenneth E. Jansen                     endif
55259599516SKenneth E. Jansenc
55359599516SKenneth E. Jansen                     lhs = 1 - min(1,mod(ifuncs(isclr+2)-1,
55459599516SKenneth E. Jansen     &                                       LHSupd(isclr+2)))
55559599516SKenneth E. Jansen                     iprec = lhs
55659599516SKenneth E. Jansen                     istop=0
557513954efSKenneth E. Jansen                 if(usingPETSc.eq.1) then
55817860365SKenneth E. Jansen#if (HAVE_PETSC)
559513954efSKenneth E. Jansen                     call SolGMRpSclr(y,             ac,
560513954efSKenneth E. Jansen     &                    x,             elDw,
561513954efSKenneth E. Jansen     &                    iBC,           BC,
562513954efSKenneth E. Jansen     &                    colm,           rowp,
563513954efSKenneth E. Jansen     &                    iper,          ilwork,
564513954efSKenneth E. Jansen     &                    shp,           shgl,
565513954efSKenneth E. Jansen     &                    shpb,          shglb,     rest,
566513954efSKenneth E. Jansen     &                    solinc(1,isclr+5),fncorp)
56717860365SKenneth E. Jansen#else
56817860365SKenneth E. Jansen                     write(*,*) 'exiting because run time input asked for PETSc, not linked in exec'
56917860365SKenneth E. Jansen                     call error('itrdrv  ','noPETSc',usingpetsc)
57017860365SKenneth E. Jansen#endif
571513954efSKenneth E. Jansen                 else
57259599516SKenneth E. Jansen                     call SolGMRSclr(y,             ac,         yold,
57359599516SKenneth E. Jansen     &                    acold,         EGmasst(1,isclr),
57459599516SKenneth E. Jansen     &                    x,             elDw,
57559599516SKenneth E. Jansen     &                    iBC,           BC,
57659599516SKenneth E. Jansen     &                    rest,
577e72fd12cSBen Matthews     &                    HBrg,          eBrg,
578e72fd12cSBen Matthews     &                    yBrg,          Rcos,      Rsin,
57959599516SKenneth E. Jansen     &                    iper,          ilwork,
58059599516SKenneth E. Jansen     &                    shp,           shgl,
58159599516SKenneth E. Jansen     &                    shpb,          shglb, solinc(1,isclr+5))
582513954efSKenneth E. Jansen                  endif
58359599516SKenneth E. Jansenc
58459599516SKenneth E. Jansen                  endif         ! end of scalar type solve
58559599516SKenneth E. Jansenc
58659599516SKenneth E. Jansenc
58759599516SKenneth E. Jansenc.... end of the multi-corrector loop
58859599516SKenneth E. Jansenc
58959599516SKenneth E. Jansen 1000             continue      !check this
59059599516SKenneth E. Jansen
59159599516SKenneth E. Jansen               else             ! this is an update  (mod did not equal zero)
59259599516SKenneth E. Jansen                  iupdate=icode/10 ! what to update
59359599516SKenneth E. Jansen                  if(iupdate.eq.0) then !update flow
59459599516SKenneth E. Jansen                     call itrCorrect ( y, ac, yold, acold, solinc)
59559599516SKenneth E. Jansen                     call itrBC (y,  ac,  iBC,  BC, iper, ilwork)
59659599516SKenneth E. Jansen                     call tnanq(y, 5, 'y_updbc')
59759599516SKenneth E. Jansenc Elaine-SPEBC
59859599516SKenneth E. Jansen                     if((irscale.ge.0).and.(myrank.eq.master)) then
59959599516SKenneth E. Jansen                        call genscale(y, x, iBC)
60059599516SKenneth E. Jansenc                       call itrBC (y,  ac,  iBC,  BC, iper, ilwork)
60159599516SKenneth E. Jansen                     endif
60259599516SKenneth E. Jansen                  else          ! update scalar
60359599516SKenneth E. Jansen                     isclr=iupdate !unless
60459599516SKenneth E. Jansen                     if(iupdate.eq.nsclr+1) isclr=0
60559599516SKenneth E. Jansen                     call itrCorrectSclr ( y, ac, yold, acold,
60659599516SKenneth E. Jansen     &                                     solinc(1,isclr+5))
60759599516SKenneth E. Jansen                     if (ilset.eq.2 .and. isclr.eq.2)  then
60859599516SKenneth E. Jansen                        fct2=one/almi
60959599516SKenneth E. Jansen                        fct3=one/alfi
61059599516SKenneth E. Jansen                        acold(:,7) = acold(:,7)
61159599516SKenneth E. Jansen     &                             + (ac(:,7)-acold(:,7))*fct2
61259599516SKenneth E. Jansen                        yold(:,7)  = yold(:,7)
61359599516SKenneth E. Jansen     &                             + (y(:,7)-yold(:,7))*fct3
61459599516SKenneth E. Jansen                        call itrBCSclr (  yold,  acold,  iBC,  BC,
61559599516SKenneth E. Jansen     &                                    iper,  ilwork)
61659599516SKenneth E. Jansen                        ac(:,7) = acold(:,7)*(one-almi/gami)
61759599516SKenneth E. Jansen                        y(:,7)  = yold(:,7)
61859599516SKenneth E. Jansen                        ac(:,7) = zero
61959599516SKenneth E. Jansen                        if (ivconstraint .eq. 1) then
62059599516SKenneth E. Jansen     &
62159599516SKenneth E. Jansenc ... applying the volume constraint
62259599516SKenneth E. Jansenc
62359599516SKenneth E. Jansen                           call solvecon (y,    x,      iBC,  BC,
62459599516SKenneth E. Jansen     &                                    iper, ilwork, shp,  shgl)
62559599516SKenneth E. Jansenc
62659599516SKenneth E. Jansen                        endif   ! end of volume constraint calculations
62759599516SKenneth E. Jansen                     endif
62859599516SKenneth E. Jansen                     call itrBCSclr (  y,  ac,  iBC,  BC, iper, ilwork)
62959599516SKenneth E. Jansen                  endif
63059599516SKenneth E. Jansen               endif            !end of switch between solve or update
63159599516SKenneth E. Jansen            enddo               ! loop over sequence in step
63259599516SKenneth E. Jansen        if((istop.lt.0).and.(iMoreRANS.lt.5)) then
63359599516SKenneth E. Jansen            iMoreRANS=iMoreRANS+1
634513954efSKenneth E. Jansen            if(myrank.eq.master) write(*,*) 'istop =', istop
63559599516SKenneth E. Jansen       iseqStart=iLastRANS
63659599516SKenneth E. Jansen       goto 9876
63759599516SKenneth E. Jansen       endif
63859599516SKenneth E. Jansenc
63959599516SKenneth E. Jansenc     Find the solution at the end of the timestep and move it to old
64059599516SKenneth E. Jansenc
64159599516SKenneth E. Jansenc.... First to reassign the parameters for the original time integrator scheme
64259599516SKenneth E. Jansenc
64359599516SKenneth E. Jansen            if((iLSet.eq.2).and.(ilss.eq.1)) then
64459599516SKenneth E. Jansen               alfi =alfit
64559599516SKenneth E. Jansen               gami =gamit
64659599516SKenneth E. Jansen               almi =almit
64759599516SKenneth E. Jansen            endif
64859599516SKenneth E. Jansen            call itrUpdate( yold,  acold,   y,    ac)
64959599516SKenneth E. Jansen            call itrBC (yold, acold,  iBC,  BC, iper,ilwork)
65059599516SKenneth E. Jansenc Elaine-SPEBC
65159599516SKenneth E. Jansen            if((irscale.ge.0).and.(myrank.eq.master)) then
65259599516SKenneth E. Jansen                call genscale(yold, x, iBC)
65359599516SKenneth E. Jansenc               call itrBC (y,  ac,  iBC,  BC, iper, ilwork)
65459599516SKenneth E. Jansen            endif
65559599516SKenneth E. Jansen            do isclr=1,nsclr
65659599516SKenneth E. Jansen               call itrBCSclr (yold, acold,  iBC, BC, iper, ilwork)
65759599516SKenneth E. Jansen            enddo
65859599516SKenneth E. Jansenc
65959599516SKenneth E. Jansen            istep = istep + 1
66059599516SKenneth E. Jansen            lstep = lstep + 1
66159599516SKenneth E. Jansen            ntoutv=max(ntout,100)
662513954efSKenneth E. Jansen            !boundary flux output moved after the error calculation so
663513954efSKenneth E. Jansen            !everything can be written out in a single chunk of code -
664513954efSKenneth E. Jansen            !Nicholas Mati
66559599516SKenneth E. Jansen
666513954efSKenneth E. Jansen            !dump TIME SERIES
66759599516SKenneth E. Jansen            if (exts) then
668513954efSKenneth E. Jansen              !Write the probe data to disc. Note that IO to disc only
669513954efSKenneth E. Jansen              !occurs when mod(lstep, nbuff) == 0. However, this
670513954efSKenneth E. Jansen              !function also does data buffering and must be called
671513954efSKenneth E. Jansen              !every time step.
672513954efSKenneth E. Jansen              call TD_bufferData()
673513954efSKenneth E. Jansen              call TD_writeData(fvarts, .false.)
67459599516SKenneth E. Jansen            endif
67559599516SKenneth E. Jansen
676513954efSKenneth E. Jansen            !Update the Mach Control
677513954efSKenneth E. Jansen            if(exts .and. exMC) then
678513954efSKenneth E. Jansen              !Note: the function MC_updateState must be called after
679513954efSKenneth E. Jansen              !the function TD_bufferData due to dependencies on
680513954efSKenneth E. Jansen              !vartsbuff(:,:,:).
681513954efSKenneth E. Jansen              call MC_updateState()
682513954efSKenneth E. Jansen              call MC_applyBC(BC)
683513954efSKenneth E. Jansen              call MC_printState()
68459599516SKenneth E. Jansen
685513954efSKenneth E. Jansen              !Write the state if a restart is also being written.
686513954efSKenneth E. Jansen              if(mod(lstep,ntout).eq.0 ) call MC_writeState(lstep)
687513954efSKenneth E. Jansen            endif
68859599516SKenneth E. Jansen
689513954efSKenneth E. Jansen            !update blower control
690513954efSKenneth E. Jansen            if(BC_enable) then
691513954efSKenneth E. Jansen              !Update the blower boundary conditions for the next
692513954efSKenneth E. Jansen              !iteration.
693513954efSKenneth E. Jansen              call BC_iter(BC)
69459599516SKenneth E. Jansen
695513954efSKenneth E. Jansen              !Also write the current phases of the blowers if a
696513954efSKenneth E. Jansen              !restart is also being written.
697513954efSKenneth E. Jansen              if(mod(lstep, ntout) == 0) call BC_writePhase(lstep)
698513954efSKenneth E. Jansen            endif
69959599516SKenneth E. Jansen
700513954efSKenneth E. Jansen            !.... Yi Chen Duct geometry8
701513954efSKenneth E. Jansen            if(isetBlowing_Duct.gt.0)then
702513954efSKenneth E. Jansen              if(ifixBlowingVel_Duct.eq.0)then
703513954efSKenneth E. Jansen                if(nstp.gt.nBlowingStepsDuct)then
704513954efSKenneth E. Jansen                  nBlowingStepsDuct = nstp-2
705513954efSKenneth E. Jansen                endif
706513954efSKenneth E. Jansen                call setBlowing_Duct2(x,BC,yold,iTurbWall,istp)
70759599516SKenneth E. Jansen              endif
70859599516SKenneth E. Jansen            endif
709513954efSKenneth E. Jansen          !... Yi Chen Duct geometry8
71059599516SKenneth E. Jansen
71159599516SKenneth E. Jansenc
71259599516SKenneth E. Jansenc.... -------------------> error calculation  <-----------------
71359599516SKenneth E. Jansen            if(ierrcalc.eq.1.or.ioybar.eq.1) then
71459599516SKenneth E. Jansen               tfact=one/istep
71559599516SKenneth E. Jansen               do idof=1,ndof
71659599516SKenneth E. Jansen                 ybar(:,idof) =tfact*yold(:,idof) +
71759599516SKenneth E. Jansen     &                         (one-tfact)*ybar(:,idof)
71859599516SKenneth E. Jansen               enddo
71959599516SKenneth E. Jansenc....compute average
72059599516SKenneth E. Jansenc...  ybar(:,ndof+1:ndof+8) is for avg. of square as uu, vv, ww, pp, TT, uv, uw, and vw
72159599516SKenneth E. Jansen               do idof=1,5 ! avg. of square for 5 flow variables
72259599516SKenneth E. Jansen                   ybar(:,ndof+idof) = tfact*yold(:,idof)**2 +
72359599516SKenneth E. Jansen     &                             (one-tfact)*ybar(:,ndof+idof)
72459599516SKenneth E. Jansen               enddo
72559599516SKenneth E. Jansen               ybar(:,ndof+6) = tfact*yold(:,1)*yold(:,2) + !uv
72659599516SKenneth E. Jansen     &                          (one-tfact)*ybar(:,ndof+6)
72759599516SKenneth E. Jansen               ybar(:,ndof+7) = tfact*yold(:,1)*yold(:,3) + !uw
72859599516SKenneth E. Jansen     &                          (one-tfact)*ybar(:,ndof+7)
72959599516SKenneth E. Jansen               ybar(:,ndof+8) = tfact*yold(:,2)*yold(:,3) + !vw
73059599516SKenneth E. Jansen     &                          (one-tfact)*ybar(:,ndof+8)
73159599516SKenneth E. Jansenc... compute err
732*f0b806cbSKenneth E. Jansenc hack ShockError
733*f0b806cbSKenneth E. Jansenc
734*f0b806cbSKenneth E. Jansen               errmax=maxval(rerr(:,6))
735*f0b806cbSKenneth E. Jansen               errswitch=0.1*errmax
736*f0b806cbSKenneth E. Jansen               where(rerr(:,6).gt.errswitch)
737*f0b806cbSKenneth E. Jansen                    rerr(:,6)=1.0
738*f0b806cbSKenneth E. Jansen               elsewhere
739*f0b806cbSKenneth E. Jansen                    rerr(:,6)=1e-10
740*f0b806cbSKenneth E. Jansen               endwhere
741513954efSKenneth E. Jansen               rerr(:, 7)=rerr(:, 7)+(yold(:,1)-ybar(:,1))**2
742513954efSKenneth E. Jansen               rerr(:, 8)=rerr(:, 8)+(yold(:,2)-ybar(:,2))**2
743513954efSKenneth E. Jansen               rerr(:, 9)=rerr(:, 9)+(yold(:,3)-ybar(:,3))**2
744513954efSKenneth E. Jansen               rerr(:,10)=rerr(:,10)+(yold(:,4)-ybar(:,4))**2
74559599516SKenneth E. Jansen            endif
746513954efSKenneth E. Jansen
747513954efSKenneth E. Jansenc.. writing ybar field if requested in each restart file
748513954efSKenneth E. Jansen
749513954efSKenneth E. Jansen            !here is where we save our averaged field.  In some cases we want to
750513954efSKenneth E. Jansen            !write it less frequently
751513954efSKenneth E. Jansen            if( (irs >= 1) .and. (
752513954efSKenneth E. Jansen     &        mod(lstep, ntout) == 0 .or. !Checkpoint
753513954efSKenneth E. Jansen     &        istep == nstp) )then        !End of simulation
754*f0b806cbSKenneth E. Jansen              if(output_mode .eq. -1 ) then ! this is an in-memory adapt case
755*f0b806cbSKenneth E. Jansen                if(istep == nstp) then ! go ahead and take care of it
756*f0b806cbSKenneth E. Jansen                  call checkpoint (nstp,yold, acold, ybar, rerr,  velbar,
757*f0b806cbSKenneth E. Jansen     &                       x, iper, ilwork, shp, shgl, iBC )
758513954efSKenneth E. Jansen                endif
759*f0b806cbSKenneth E. Jansen                if(ntout.le.lstep) then ! user also wants file output
760*f0b806cbSKenneth E. Jansen                  output_mode=0
761*f0b806cbSKenneth E. Jansen                  call checkpoint (nstp,yold, acold, ybar, rerr,  velbar,
762*f0b806cbSKenneth E. Jansen     &                       x, iper, ilwork, shp, shgl, iBC )
763*f0b806cbSKenneth E. Jansen                  output_mode=-1 ! reset to stream
764513954efSKenneth E. Jansen                endif
765*f0b806cbSKenneth E. Jansen              else
766*f0b806cbSKenneth E. Jansen                call checkpoint (nstp,yold, acold, ybar, rerr,  velbar,
767*f0b806cbSKenneth E. Jansen     &                       x, iper, ilwork, shp, shgl, iBC )
768513954efSKenneth E. Jansen              endif
769513954efSKenneth E. Jansen             endif
770513954efSKenneth E. Jansen
771513954efSKenneth E. Jansen 2000    continue  !end of NSTEP loop
77259599516SKenneth E. Jansen 2001    continue
77359599516SKenneth E. Jansen
77459599516SKenneth E. Jansen         ttim(1) = REAL(secs(0.0)) / 100. - ttim(1)
77559599516SKenneth E. Jansen         ttim(2) = secs(0.0)              - ttim(2)
77659599516SKenneth E. Jansen
77759599516SKenneth E. Jansenc         tcorecp2 = REAL(secs(0.0)) / 100.
77859599516SKenneth E. Jansenc         tcorewc2 = secs(0.0)
77959599516SKenneth E. Jansen         if (numpe > 1) call MPI_BARRIER(MPI_COMM_WORLD, ierr)
780513954efSKenneth E. Jansen         if(myrank.eq.master)  then
78159599516SKenneth E. Jansen            tcorecp2 = TMRC()
78259599516SKenneth E. Jansen            write(6,*) 'T(core) cpu = ',tcorecp2-tcorecp1
78359599516SKenneth E. Jansen         endif
78459599516SKenneth E. Jansen
78559599516SKenneth E. Jansenc     call wtime
78659599516SKenneth E. Jansen
78775f1c48cSCameron Smith      call destroyWallData
78892e15685SCameron Smith      call destroyfncorp
78975f1c48cSCameron Smith
790513954efSKenneth E. Jansen 3000 continue !end of NTSEQ loop
79159599516SKenneth E. Jansenc
79259599516SKenneth E. Jansenc.... ---------------------->  Post Processing  <----------------------
79359599516SKenneth E. Jansenc
79459599516SKenneth E. Jansenc.... print out the last step
79559599516SKenneth E. Jansenc
796513954efSKenneth E. Jansen!      if ( (irs .ge. 1) .and. ((mod(lstep, ntout) .ne. 0) .or.
797513954efSKenneth E. Jansen!     &    (nstp .eq. 0))) then
798513954efSKenneth E. Jansen!          if( (mod(lstep, ntoutv) .eq. 0) .and.
799513954efSKenneth E. Jansen!     &        ((irscale.ge.0).or.(itwmod.gt.0) .or.
800513954efSKenneth E. Jansen!     &        ((nsonmax.eq.1).and.(iLES.gt.0))))
801513954efSKenneth E. Jansen!     &        call rwvelb  ('out ',  velbar  ,ifail)
802513954efSKenneth E. Jansen!
803513954efSKenneth E. Jansen!          call Bflux  (yold,  acold,     x,
804513954efSKenneth E. Jansen!     &         shp,           shgl,      shpb,
805513954efSKenneth E. Jansen!     &         shglb,         nodflx,    ilwork)
806513954efSKenneth E. Jansen!      endif
80759599516SKenneth E. Jansen
808513954efSKenneth E. Jansen
809513954efSKenneth E. Jansen
810513954efSKenneth E. Jansenc      if(ioybar.eq.1) then
811513954efSKenneth E. Jansenc         call write_field(myrank,'a'//char(0),'ybar'//char(0),4,
812513954efSKenneth E. Jansenc     &                      ybar,'d'//char(0),nshg,ndof+8,lstep)
813513954efSKenneth E. Jansenc      endif
814513954efSKenneth E. Jansen
815513954efSKenneth E. Jansenc     if(iRANS.lt.0 .and. idistcalc.eq.1) then
816513954efSKenneth E. Jansenc        call write_field(myrank,'a'//char(0),'DESd'//char(0),4,
817513954efSKenneth E. Jansenc     &                      elDw,'d'//char(0),numel,1,lstep)
81859599516SKenneth E. Jansenc
819513954efSKenneth E. Jansenc         call write_field(myrank,'a'//char(0),'dwal'//char(0),4,
820513954efSKenneth E. Jansenc     &                    d2wall,'d'//char(0),nshg,1,lstep)
821513954efSKenneth E. Jansenc     endif
82259599516SKenneth E. Jansen
82359599516SKenneth E. Jansenc
82459599516SKenneth E. Jansenc.... close history and aerodynamic forces files
82559599516SKenneth E. Jansenc
82659599516SKenneth E. Jansen      if (myrank .eq. master) then
82759599516SKenneth E. Jansen         close (ihist)
82859599516SKenneth E. Jansen         close (iforce)
829513954efSKenneth E. Jansen
830513954efSKenneth E. Jansen         if(exMC) then
831513954efSKenneth E. Jansen           call MC_writeState(lstep)
832513954efSKenneth E. Jansen           call MC_finalize
833513954efSKenneth E. Jansen         endif
834513954efSKenneth E. Jansen
83559599516SKenneth E. Jansen         if(exts) then
836513954efSKenneth E. Jansen           call TD_writeData(fvarts, .true.)    !force the flush of the buffer.
837513954efSKenneth E. Jansen           call TD_finalize
83859599516SKenneth E. Jansen         endif
83959599516SKenneth E. Jansen      endif
840513954efSKenneth E. Jansen
841513954efSKenneth E. Jansen      if(BC_enable) then  !blower is allocated on all processes.
842513954efSKenneth E. Jansen        if(mod(lstep, ntout) /= 0) then !May have already written file.
843513954efSKenneth E. Jansen           call BC_writePhase(lstep)
844513954efSKenneth E. Jansen        endif
845513954efSKenneth E. Jansen
846513954efSKenneth E. Jansen        call BC_finalize
847513954efSKenneth E. Jansen      endif
848513954efSKenneth E. Jansen
84959599516SKenneth E. Jansen      close (iecho)
85059599516SKenneth E. Jansen      if(iabc==1) deallocate(acs)
85159599516SKenneth E. Jansenc
85259599516SKenneth E. Jansenc.... end
85359599516SKenneth E. Jansenc
85459599516SKenneth E. Jansen      return
85559599516SKenneth E. Jansen      end
856513954efSKenneth E. Jansen
857*f0b806cbSKenneth E. Jansen      subroutine checkpoint (nstp,yold, acold, ybar, rerr,  velbar,
858*f0b806cbSKenneth E. Jansen     &                       x, iper, ilwork, shp, shgl, iBC )
859*f0b806cbSKenneth E. Jansenc
860*f0b806cbSKenneth E. Jansen      use turbSA
861*f0b806cbSKenneth E. Jansen      include "common.h"
862*f0b806cbSKenneth E. Jansen      dimension shp(MAXTOP,maxsh,MAXQPT),
863*f0b806cbSKenneth E. Jansen     &            shgl(MAXTOP,nsd,maxsh,MAXQPT),
864*f0b806cbSKenneth E. Jansen     &            iper(nshg),              iBC(nshg),
865*f0b806cbSKenneth E. Jansen     &            x(nshg,nsd),         ilwork(nlwork)
866*f0b806cbSKenneth E. Jansen      real*8  velbar(nfath,ndof),
867*f0b806cbSKenneth E. Jansen     &        yold(nshg,ndof),      acold(nshg,ndof),
868*f0b806cbSKenneth E. Jansen     &        rerr(nshg,10),        ybar(nshg,ndof+8)
869*f0b806cbSKenneth E. Jansen! 8 is for avg. of square as uu, vv, ww, pp, TT, uv, uw, and vw
870*f0b806cbSKenneth E. Jansen
871*f0b806cbSKenneth E. Jansen      if( (mod(lstep, ntout) .eq. 0) .and.
872*f0b806cbSKenneth E. Jansen     &              ((irscale.ge.0).or.(itwmod.gt.0) .or.
873*f0b806cbSKenneth E. Jansen     &              ((nsonmax.eq.1).and.(iLES.gt.0))))
874*f0b806cbSKenneth E. Jansen     &              call rwvelb  ('out ',  velbar  ,ifail)
875*f0b806cbSKenneth E. Jansen
876*f0b806cbSKenneth E. Jansen!BUG: need to update new_interface to work with SyncIO.
877*f0b806cbSKenneth E. Jansen      !Bflux is presently completely crippled. Note that restar
878*f0b806cbSKenneth E. Jansen      !has also been moved here for readability.
879*f0b806cbSKenneth E. Jansen!              call Bflux  (yold,          acold,     x,  compute boundary fluxes and print out
880*f0b806cbSKenneth E. Jansen!    &              shp,           shgl,      shpb,
881*f0b806cbSKenneth E. Jansen!    &              shglb,         nodflx,    ilwork)
882*f0b806cbSKenneth E. Jansen
883*f0b806cbSKenneth E. Jansen      call timer ('Output  ')      !set up the timer
884*f0b806cbSKenneth E. Jansen
885*f0b806cbSKenneth E. Jansen      !write the solution and time derivative
886*f0b806cbSKenneth E. Jansen      call restar ('out ',  yold, acold)
887*f0b806cbSKenneth E. Jansen
888*f0b806cbSKenneth E. Jansen      !Write the distance to wall field in each restart
889*f0b806cbSKenneth E. Jansen      if((istep==nstp) .and. (irans < 0 )) then !d2wall is allocated
890*f0b806cbSKenneth E. Jansen                 call write_field(myrank,'a'//char(0),'dwal'//char(0),4,
891*f0b806cbSKenneth E. Jansen     &                            d2wall,'d'//char(0), nshg, 1, lstep)
892*f0b806cbSKenneth E. Jansen      endif
893*f0b806cbSKenneth E. Jansen
894*f0b806cbSKenneth E. Jansen      !Write the time average in each restart.
895*f0b806cbSKenneth E. Jansen      if(ioybar.eq.1)then
896*f0b806cbSKenneth E. Jansen                 call write_field(myrank,'a'//char(0),'ybar'//char(0),4,
897*f0b806cbSKenneth E. Jansen     &                              ybar,'d'//char(0),nshg,ndof+8,lstep)
898*f0b806cbSKenneth E. Jansen      endif
899*f0b806cbSKenneth E. Jansen
900*f0b806cbSKenneth E. Jansen      !Write the error feild at the end of each step sequence
901*f0b806cbSKenneth E. Jansen      if(ierrcalc.eq.1 .and. istep == nstp) then
902*f0b806cbSKenneth E. Jansen        !smooth the error indicators
903*f0b806cbSKenneth E. Jansen
904*f0b806cbSKenneth E. Jansen        do i=1,ierrsmooth
905*f0b806cbSKenneth E. Jansen         call errsmooth( rerr, x, iper, ilwork, shp, shgl, iBC )
906*f0b806cbSKenneth E. Jansen        enddo
907*f0b806cbSKenneth E. Jansen
908*f0b806cbSKenneth E. Jansen        call write_field( myrank, 'a'//char(0), 'errors'//char(0), 6,
909*f0b806cbSKenneth E. Jansen     &                        rerr, 'd'//char(0), nshg, 10, lstep)
910*f0b806cbSKenneth E. Jansen      endif
911*f0b806cbSKenneth E. Jansen
912*f0b806cbSKenneth E. Jansenc the following is a future way to have the number of steps in the header...done for posix but not yet for syncio
913*f0b806cbSKenneth E. Jansenc
914*f0b806cbSKenneth E. Jansenc              call write_field2(myrank,'a'//char(0),'ybar'//char(0),
915*f0b806cbSKenneth E. Jansenc     &                          4,ybar,'d'//char(0),nshg,ndof+8,
916*f0b806cbSKenneth E. Jansenc     &                         lstep,istep)
917*f0b806cbSKenneth E. Jansenc
918*f0b806cbSKenneth E. Jansenc.... end
919*f0b806cbSKenneth E. Jansenc
920*f0b806cbSKenneth E. Jansen      return
921*f0b806cbSKenneth E. Jansen      end
922513954efSKenneth E. Jansen
923