xref: /phasta/phSolver/common/proces.f (revision bd5c8633bf620393147397a304ccd0eac0e12731)
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
15*bd5c8633SCameron Smith        use periodicity
16*bd5c8633SCameron 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
4859599516SKenneth E. Jansenc        if ((irscale .ge. 0) .and. (myrank.eq.master)) then
4959599516SKenneth E. Jansenc           call setSPEBC(numnp, point2nfath, nsonmax)
5059599516SKenneth E. Jansenc        endif
5159599516SKenneth E. Jansenc
5259599516SKenneth E. Jansenc.... generate the geometry and boundary conditions data
5359599516SKenneth E. Jansenc
5459599516SKenneth E. Jansen        call gendat (y,              ac,             point2x,
5559599516SKenneth E. Jansen     &               iBC,            BC,
5659599516SKenneth E. Jansen     &               point2iper,     point2ilwork,   shp,
5759599516SKenneth E. Jansen     &               shgl,           shpb,           shglb,
5859599516SKenneth E. Jansen     &               point2ifath,    velbar,         point2nsons )
5959599516SKenneth E. Jansen        call setper(nshg)
6059599516SKenneth E. Jansen        call perprep(iBC,point2iper,nshg)
6159599516SKenneth E. Jansen        if (iLES/10 .eq. 1) then
6259599516SKenneth E. Jansen        call keeplhsG ! Allocating space for the mass (Gram) matrix to be
6359599516SKenneth E. Jansen                      ! used for projection filtering and reconstruction
6459599516SKenneth E. Jansen                      ! of the strain rate tensor.
6559599516SKenneth E. Jansen
6659599516SKenneth E. Jansen        call setrls   ! Allocating space for the resolved Leonard stresses
6759599516SKenneth E. Jansenc                         See bardmc.f
6859599516SKenneth E. Jansen        endif
6959599516SKenneth E. Jansenc
7059599516SKenneth E. Jansenc.... time averaged statistics
7159599516SKenneth E. Jansenc
7259599516SKenneth E. Jansen        if (ioform .eq. 2) then
7359599516SKenneth E. Jansen           call initStats(point2x, iBC, point2iper, point2ilwork)
7459599516SKenneth E. Jansen        endif
7559599516SKenneth E. Jansenc
7659599516SKenneth E. Jansenc.... RANS turbulence model
7759599516SKenneth E. Jansenc
7859599516SKenneth E. Jansen        if (iRANS .lt. 0) then
7959599516SKenneth E. Jansen           call initTurb( point2x )
8059599516SKenneth E. Jansen        endif
8159599516SKenneth E. Jansenc
8259599516SKenneth E. Jansenc.... p vs. Q boundary
8359599516SKenneth E. Jansenc
8459599516SKenneth E. Jansen           call initNABI( point2x, shpb )
8559599516SKenneth E. Jansenc
8659599516SKenneth E. Jansenc.... check for execution mode
8759599516SKenneth E. Jansenc
8859599516SKenneth E. Jansen        if (iexec .eq. 0) then
8959599516SKenneth E. Jansen           lstep = 0
9059599516SKenneth E. Jansen           call restar ('out ',  y  ,ac)
9159599516SKenneth E. Jansen           return
9259599516SKenneth E. Jansen        endif
9359599516SKenneth E. Jansenc
9459599516SKenneth E. Jansenc.... initialize AutoSponge
9559599516SKenneth E. Jansenc
9659599516SKenneth E. Jansen        if(matflg(5,1).ge.4) then ! cool case (sponge)
9759599516SKenneth E. Jansen           call initSponge( y,point2x)
9859599516SKenneth E. Jansen        endif
9959599516SKenneth E. Jansenc
10059599516SKenneth E. Jansenc
10159599516SKenneth E. Jansenc.... adjust BC's to interpolate from file
10259599516SKenneth E. Jansenc
10359599516SKenneth E. Jansen
10459599516SKenneth E. Jansen        inquire(file="inlet.dat",exist=exlog)
10559599516SKenneth E. Jansen        if(exlog) then
10659599516SKenneth E. Jansen           open (unit=654,file="inlet.dat",status="old")
10759599516SKenneth E. Jansen           read(654,*) ninterp,ibcmatch,(interp_mask(j),j=1,ndof)
10859599516SKenneth E. Jansen           do i=1,ninterp
10959599516SKenneth E. Jansen              read(654,*) (bcinterp(i,j),j=1,ndof+1) ! distance to wall+
11059599516SKenneth E. Jansen                        ! ndof but note order of BC's
11159599516SKenneth E. Jansen                        ! p,t,u,v,w,scalars. Also note that must be in
11259599516SKenneth E. Jansen                        ! increasing distance from the wall order.
11359599516SKenneth E. Jansen
11459599516SKenneth E. Jansen           enddo
11559599516SKenneth E. Jansen           do i=1,nshg  ! only correct for linears at this time
11659599516SKenneth E. Jansen              if(mod(iBC(i),1024).eq.ibcmatch) then
11759599516SKenneth E. Jansen                 iupper=0
11859599516SKenneth E. Jansen                 do j=2,ninterp
11959599516SKenneth E. Jansen                    if(bcinterp(j,1).gt.d2wall(i)) then !bound found
12059599516SKenneth E. Jansen                       xi=(d2wall(i)-bcinterp(j-1,1))/
12159599516SKenneth E. Jansen     &                    (bcinterp(j,1)-bcinterp(j-1,1))
12259599516SKenneth E. Jansen                       iupper=j
12359599516SKenneth E. Jansen                       exit
12459599516SKenneth E. Jansen                    endif
12559599516SKenneth E. Jansen                 enddo
12659599516SKenneth E. Jansen                 if(iupper.eq.0) then
12759599516SKenneth E. Jansen                    write(*,*) 'failure in finterp, ynint=',d2wall(i)
12859599516SKenneth E. Jansen                    stop
12959599516SKenneth E. Jansen                 endif
13059599516SKenneth E. Jansen                 if(interp_mask(1).ne.zero) then
13159599516SKenneth E. Jansen                    BC(i,1)=(xi*bcinterp(iupper,2)
13259599516SKenneth E. Jansen     &                +(one-xi)*bcinterp(iupper-1,2))*interp_mask(1)
13359599516SKenneth E. Jansen                 endif
13459599516SKenneth E. Jansen                 if(interp_mask(2).ne.zero) then
13559599516SKenneth E. Jansen                    BC(i,2)=(xi*bcinterp(iupper,3)
13659599516SKenneth E. Jansen     &                +(one-xi)*bcinterp(iupper-1,3))*interp_mask(2)
13759599516SKenneth E. Jansen                 endif
13859599516SKenneth E. Jansen                 if(interp_mask(3).ne.zero) then
13959599516SKenneth E. Jansen                    BC(i,3)=(xi*bcinterp(iupper,4)
14059599516SKenneth E. Jansen     &                +(one-xi)*bcinterp(iupper-1,4))*interp_mask(3)
14159599516SKenneth E. Jansen                 endif
14259599516SKenneth E. Jansen                 if(interp_masK(4).ne.zero) then
14359599516SKenneth E. Jansen                    BC(i,4)=(xi*bcinterp(iupper,5)
14459599516SKenneth E. Jansen     &                +(one-xi)*bcinterp(iupper-1,5))*interp_mask(4)
14559599516SKenneth E. Jansen                 endif
14659599516SKenneth E. Jansen                 if(interp_mask(5).ne.zero) then
14759599516SKenneth E. Jansen                    BC(i,5)=(xi*bcinterp(iupper,6)
14859599516SKenneth E. Jansen     &                +(one-xi)*bcinterp(iupper-1,6))*interp_mask(5)
14959599516SKenneth E. Jansen                 endif
15059599516SKenneth E. Jansen                 if(interp_mask(6).ne.zero) then
15159599516SKenneth E. Jansen                    BC(i,7)=(xi*bcinterp(iupper,7)
15259599516SKenneth E. Jansen     &                +(one-xi)*bcinterp(iupper-1,7))*interp_mask(6)
15359599516SKenneth E. Jansen                 endif
15459599516SKenneth E. Jansen                 if(interp_mask(7).ne.zero) then
15559599516SKenneth E. Jansen                    BC(i,8)=(xi*bcinterp(iupper,8)
15659599516SKenneth E. Jansen     &                +(one-xi)*bcinterp(iupper-1,8))*interp_mask(7)
15759599516SKenneth E. Jansen                 endif
15859599516SKenneth E. Jansen              endif
15959599516SKenneth E. Jansen           enddo
16059599516SKenneth E. Jansen        endif
16159599516SKenneth E. Jansenc$$$$$$$$$$$$$$$$$$$$
16259599516SKenneth E. Jansen
16359599516SKenneth E. Jansenc
16459599516SKenneth E. Jansenc
16559599516SKenneth E. Jansenc.... call the semi-discrete predictor multi-corrector iterative driver
16659599516SKenneth E. Jansenc
16759599516SKenneth E. Jansen        call itrdrv (y,              ac,
16859599516SKenneth E. Jansen     &               uold,           point2x,
16959599516SKenneth E. Jansen     &               iBC,            BC,
17059599516SKenneth E. Jansen     &               point2iper,     point2ilwork,   shp,
17159599516SKenneth E. Jansen     &               shgl,           shpb,           shglb,
17259599516SKenneth E. Jansen     &               point2ifath,    velbar,         point2nsons )
17359599516SKenneth E. Jansenc
17459599516SKenneth E. Jansenc.... return
17559599516SKenneth E. Jansenc
17659599516SKenneth E. Jansenc
17759599516SKenneth E. Jansenc.... stop CPU-timer
17859599516SKenneth E. Jansenc
17959599516SKenneth E. JansenCAD        call timer ('End     ')
18059599516SKenneth E. Jansenc
18159599516SKenneth E. Jansenc.... close echo file
18259599516SKenneth E. Jansenc
18359599516SKenneth E. Jansen
18459599516SKenneth E. Jansen!MR CHANGE
18559599516SKenneth E. Jansen      if (numpe > 1) call MPI_BARRIER(MPI_COMM_WORLD, ierr)
18659599516SKenneth E. Jansen      if(myrank.eq.0)  then
18759599516SKenneth E. Jansen          write(*,*) 'process - before closing iecho'
18859599516SKenneth E. Jansen      endif
18959599516SKenneth E. Jansen!MR CHANGE END
19059599516SKenneth E. Jansen
19159599516SKenneth E. Jansen        close (iecho)
19259599516SKenneth E. Jansen
19359599516SKenneth E. Jansen!MR CHANGE
19459599516SKenneth E. Jansen      if (numpe > 1) call MPI_BARRIER(MPI_COMM_WORLD, ierr)
19559599516SKenneth E. Jansen      if(myrank.eq.0)  then
19659599516SKenneth E. Jansen          write(*,*) 'process - after closing iecho'
19759599516SKenneth E. Jansen      endif
19859599516SKenneth E. Jansen!MR CHANGE END
19959599516SKenneth E. Jansen
20059599516SKenneth E. Jansen
20159599516SKenneth E. Jansenc
20259599516SKenneth E. Jansenc.... end of the program
20359599516SKenneth E. Jansenc
20459599516SKenneth E. JansenCAD        write(6,*) 'Life: ', second(0) - ttim(100)
20559599516SKenneth E. Jansen        deallocate(point2iper)
206fcc77cc2SCameron Smith        if(numpe.gt.1) then
207fcc77cc2SCameron Smith          call Dctypes(point2ilwork(1))
208fcc77cc2SCameron Smith          deallocate(point2ilwork)
209fcc77cc2SCameron Smith        endif
21059599516SKenneth E. Jansen        deallocate(point2x)
21159599516SKenneth E. Jansen        deallocate(point2nsons)
21259599516SKenneth E. Jansen        deallocate(point2ifath)
2136b966dd8SCameron Smith        deallocate(uold)
2146b966dd8SCameron Smith        deallocate(wnrm)
2156b966dd8SCameron Smith        deallocate(otwn)
2166b966dd8SCameron Smith        call finalizeDtN
217*bd5c8633SCameron Smith        call clearper
218*bd5c8633SCameron Smith        call finalizeNABI
2196b966dd8SCameron Smith
22059599516SKenneth E. Jansen        return
22159599516SKenneth E. Jansen        end
22259599516SKenneth E. Jansen
22359599516SKenneth E. Jansen
224