xref: /phasta/phSolver/common/turbsa.f (revision 513954ef803c300cddba2bb96b4a5dac0b93489a)
159599516SKenneth E. Jansenc-----------------------------------------------------------------------
259599516SKenneth E. Jansenc
359599516SKenneth E. Jansenc     Spalart-Allmaras turbulence model constants
459599516SKenneth E. Jansenc
559599516SKenneth E. Jansenc-----------------------------------------------------------------------
659599516SKenneth E. Jansen      module turbSA
759599516SKenneth E. Jansen
859599516SKenneth E. Jansen      real*8, allocatable ::  d2wall(:)
959599516SKenneth E. Jansen      real*8, allocatable ::  wnrm(:,:)
1059599516SKenneth E. Jansen      integer, allocatable :: otwn(:)
1159599516SKenneth E. Jansen      real*8  saCb1, saCb2, saCw1, saCw2, saCw3, saCv1, saSigma,
1259599516SKenneth E. Jansen     &        saKappa, saKappaP2Inv, saCv2P3, saCw3P6, saSigmaInv
1359599516SKenneth E. Jansen      integer, allocatable :: sidmapg(:) ! list of all surfID's, low to high
1459599516SKenneth E. Jansen
1559599516SKenneth E. Jansen      parameter (
1659599516SKenneth E. Jansen     &     saCb1        = 0.1355d0,
1759599516SKenneth E. Jansen     &     saCb2        = 0.622d0,
1859599516SKenneth E. Jansen     &     saCw1        = 3.239067817d0,
1959599516SKenneth E. Jansen     &     saCw2        = 0.3d0,
2059599516SKenneth E. Jansen     &     saCw3        = 2.0d0,
2159599516SKenneth E. Jansen     &     saKappa      = 0.41d0,
2259599516SKenneth E. Jansen     &     saSigma      = 0.666666666666666667d0,
2359599516SKenneth E. Jansen     &     saCv1        = 7.1d0,
2459599516SKenneth E. Jansen     &     saKappaP2Inv = 5.94883997620464d0,
2559599516SKenneth E. Jansen     &     saCv1P3      = 3.579109999999999d+02,
2659599516SKenneth E. Jansen     &     saCw3P6      = 64.0d0,
2759599516SKenneth E. Jansen     &     saSigmaInv   = 1.50d0
2859599516SKenneth E. Jansen     &     )
2959599516SKenneth E. Jansen
3059599516SKenneth E. Jansen
3159599516SKenneth E. Jansen      end module
3259599516SKenneth E. Jansen
3359599516SKenneth E. Jansenc-----------------------------------------------------------------------
3459599516SKenneth E. Jansenc
3559599516SKenneth E. Jansenc     Initialize: compute the distance to the nearest wall for
3659599516SKenneth E. Jansenc     each vertex in the mesh.
3759599516SKenneth E. Jansenc
3859599516SKenneth E. Jansenc     Michael Yaworski (fall 1998)
3959599516SKenneth E. Jansenc
4059599516SKenneth E. Jansenc-----------------------------------------------------------------------
4159599516SKenneth E. Jansen      subroutine initTurb( x )
4259599516SKenneth E. Jansen
4359599516SKenneth E. Jansen      use     pointer_data
4459599516SKenneth E. Jansen      use     turbSA
4559599516SKenneth E. Jansen      include "common.h"
4659599516SKenneth E. Jansen      include "mpif.h"
4759599516SKenneth E. Jansen
48*513954efSKenneth E. Jansen      character(len=20) fname1,  fmt1
4959599516SKenneth E. Jansen      real*8   x(numnp,nsd)
5059599516SKenneth E. Jansen      integer  nwall(numpe),      idisp(numpe)
51*513954efSKenneth E. Jansen      character(len=5)  cname
5259599516SKenneth E. Jansen      real*8, allocatable :: xwi(:,:,:), xw(:,:,:)
5359599516SKenneth E. Jansen
5459599516SKenneth E. Jansen!MR CHANGE
5559599516SKenneth E. Jansen      integer :: ierr
5659599516SKenneth E. Jansen      integer :: numparts, nfile, itmp2, id2wall, itwo, ifoundd2wall
5759599516SKenneth E. Jansen      integer iarray(10)
58*513954efSKenneth E. Jansen      character(len=255) :: fnamer, fname2, temp2
59*513954efSKenneth E. Jansen      character(len=64) :: temp1
6059599516SKenneth E. Jansen!MR CHANGE END
6159599516SKenneth E. Jansen
6259599516SKenneth E. Jansen      allocate ( d2wall(numnp) )
6359599516SKenneth E. Jansen
6459599516SKenneth E. Jansen!MR CHANGE
6559599516SKenneth E. Jansen      if(myrank.eq.master) then
6659599516SKenneth E. Jansen        write (*,*) 'entering initTurb'
6759599516SKenneth E. Jansen      endif
6859599516SKenneth E. Jansen      call MPI_BARRIER(MPI_COMM_WORLD,ierr)
6959599516SKenneth E. Jansen
7059599516SKenneth E. Jansen      ! First we try to read the d2wall field from either the restart files or the d2wall files
7159599516SKenneth E. Jansen      call read_d2wall(myrank,numnp,d2wall,ifoundd2wall)
7259599516SKenneth E. Jansen!MR CHANGE END
7359599516SKenneth E. Jansen
7459599516SKenneth E. Jansen      if(ifoundd2wall.eq.0) then   ! d2wall was not found so calculate the distance
7559599516SKenneth E. Jansen        if(myrank.eq.master) then
7659599516SKenneth E. Jansen          write (*,*) 'Computing the d2wall field'
7759599516SKenneth E. Jansen        endif
7859599516SKenneth E. Jansenc
7959599516SKenneth E. Jansenc   hard code trip point until we are sure it is worth doing
8059599516SKenneth E. Jansenc
8159599516SKenneth E. Jansen
8259599516SKenneth E. Jansenc
8359599516SKenneth E. Jansenc.... Count the welts (wall-elements)
8459599516SKenneth E. Jansenc
8559599516SKenneth E. Jansen         nwalli=0
8659599516SKenneth E. Jansen         do iblk = 1, nelblb    ! loop over boundary elt blocks
8759599516SKenneth E. Jansen            npro = lcblkb(1,iblk+1) - lcblkb(1,iblk)
8859599516SKenneth E. Jansen            do j = 1, npro
8959599516SKenneth E. Jansen               if(btest(miBCB(iblk)%p(j,1),4)) nwalli=nwalli+1
9059599516SKenneth E. Jansen            enddo
9159599516SKenneth E. Jansen         enddo
9259599516SKenneth E. Jansenc
9359599516SKenneth E. Jansenc.... Create wallnode-coord list for welts on processor
9459599516SKenneth E. Jansenc
9559599516SKenneth E. Jansen         if (nwalli.eq.0) nwalli=1 !  patch for mpi's lack of imagination
9659599516SKenneth E. Jansen         allocate (xwi(nwalli,nenb+1,nsd))
9759599516SKenneth E. Jansen         xwi = 1.0d32
9859599516SKenneth E. Jansen         xwi(:,nenb+1,:)=zero
9959599516SKenneth E. Jansen         nwalli = 0
10059599516SKenneth E. Jansen         do iblk = 1, nelblb    ! loop over boundary elt blocks
10159599516SKenneth E. Jansenc
10259599516SKenneth E. Jansen            iel    = lcblkb(1,iblk)
10359599516SKenneth E. Jansen            nenbl  = lcblkb(6,iblk) ! no. of vertices per bdry. face
10459599516SKenneth E. Jansen            npro   = lcblkb(1,iblk+1) - iel
10559599516SKenneth E. Jansenc
10659599516SKenneth E. Jansen            do j = 1, npro      ! loop over belts in this blk
10759599516SKenneth E. Jansen               if(btest(miBCB(iblk)%p(j,1),4)) then
10859599516SKenneth E. Jansen                  nwalli=nwalli+1
10959599516SKenneth E. Jansenc assemble local coord list
11059599516SKenneth E. Jansen                  do node = 1, nenbl
11159599516SKenneth E. Jansen                     xwi(nwalli,node,1:3)=x(mienb(iblk)%p(j,node),:)
11259599516SKenneth E. Jansen                  enddo
11359599516SKenneth E. Jansenc put the centroid coordinates in the last slot
11459599516SKenneth E. Jansen                  do node = 1, nenbl
11559599516SKenneth E. Jansen                     xwi(nwalli,nenb+1,:)=xwi(nwalli,nenb+1,:)
11659599516SKenneth E. Jansen     &                    +xwi(nwalli,node,:)
11759599516SKenneth E. Jansen                  enddo
11859599516SKenneth E. Jansen                  xwi(nwalli,nenb+1,:)=xwi(nwalli,nenb+1,:)/nenbl
11959599516SKenneth E. Jansenc
12059599516SKenneth E. Jansen               endif
12159599516SKenneth E. Jansen            enddo               ! loop over belts in this blk
12259599516SKenneth E. Jansenc
12359599516SKenneth E. Jansen         enddo                  ! loop over boundary elt blocks
12459599516SKenneth E. Jansenc
12559599516SKenneth E. Jansen         if (nwalli.eq.0) xwi=1.0e32 ! fix for mpi's lack of imagination
12659599516SKenneth E. Jansen         if (nwalli.eq.0) nwalli=1 !  patch for mpi's lack of imagination
12759599516SKenneth E. Jansenc
12859599516SKenneth E. Jansenc  Pool "number of welts" info from all processors
12959599516SKenneth E. Jansenc
13059599516SKenneth E. JansencMPI_ALLGATHER(sendbuf,sendcount,sendtype,recvbuf,recvcount,recvtype,comm)
13159599516SKenneth E. Jansenc[ IN sendbuf] starting address of send buffer (choice)
13259599516SKenneth E. Jansenc[ IN sendcount] number of elements in send buffer (integer)
13359599516SKenneth E. Jansenc[ IN sendtype] data type of send buffer elements (handle)
13459599516SKenneth E. Jansenc[ OUT recvbuf] address of receive buffer (choice)
13559599516SKenneth E. Jansenc[ IN recvcount] number of elements received from any process (integer)
13659599516SKenneth E. Jansenc[ IN recvtype] data type of receive buffer elements (handle)
13759599516SKenneth E. Jansenc[ IN comm] communicator (handle)
13859599516SKenneth E. Jansen         if (numpe.gt.1) then   ! multiple processors
13959599516SKenneth E. Jansenc write the number of wall elts on the jth processor to slot j of nwall
14059599516SKenneth E. Jansen            call MPI_ALLGATHER(nwalli,1,MPI_INTEGER,nwall,1,
14159599516SKenneth E. Jansen     &          MPI_INTEGER,MPI_COMM_WORLD,ierr)
14259599516SKenneth E. Jansenc count up the total number of wall elts among all processes
14359599516SKenneth E. Jansen            nwallt=0
14459599516SKenneth E. Jansen            do j=1,numpe
14559599516SKenneth E. Jansen               nwallt=nwallt+nwall(j)
14659599516SKenneth E. Jansen            enddo
14759599516SKenneth E. Jansen         else                   ! single processor
14859599516SKenneth E. Jansenc the local information is the global information for single-processor
14959599516SKenneth E. Jansen            nwall=nwalli
15059599516SKenneth E. Jansen            nwallt=nwalli
15159599516SKenneth E. Jansen         endif                  ! if-else for multiple processors
15259599516SKenneth E. Jansenc
15359599516SKenneth E. Jansenc  Make all-processor wallnode-coord collage
15459599516SKenneth E. Jansenc
15559599516SKenneth E. Jansen         allocate (xw(nwallt,nenb+1,nsd))
15659599516SKenneth E. Jansen         if (numpe.gt.1) then   ! multiple processors
15759599516SKenneth E. Jansenc we will gather coordinates from local on-proc sets to a global set
15859599516SKenneth E. Jansenc we will stack each processor's coordinate list atop that of the
15959599516SKenneth E. Jansenc previous processor.  If the coordinate list for processor i is
16059599516SKenneth E. Jansenc called xwi, then our global coordinate list xw will look like this:
16159599516SKenneth E. Jansenc ---------------------------------------------------------------
16259599516SKenneth E. Jansenc | xw1            | xw2                | xw3        |   ...    |
16359599516SKenneth E. Jansenc ---------------------------------------------------------------
16459599516SKenneth E. Jansenc  <---nwall(1)---> <-----nwall(2)-----> <-nwall(3)->
16559599516SKenneth E. Jansenc  <------------------------nwallt-----------------------...---->
16659599516SKenneth E. Jansenc To accomplish this with MPI, we use MPI_ALLGATHERV, summarized as:
16759599516SKenneth E. JansencMPI_ALLGATHERV(sendbuf,sendcount,sendtype,recvbuf,recvcount,disp,recvtype,comm)
16859599516SKenneth E. Jansenc[ IN sendbuf] starting address of send buffer (choice)
16959599516SKenneth E. Jansenc[ IN sendcount] number of elements in send buffer (integer)
17059599516SKenneth E. Jansenc[ IN sendtype] data type of send buffer elements (handle)
17159599516SKenneth E. Jansenc[ OUT recvbuf] address of receive buffer (choice)
17259599516SKenneth E. Jansenc[ IN recvcount] number of elements received from any process (integer)
17359599516SKenneth E. Jansenc[ IN disp] displacement array
17459599516SKenneth E. Jansenc[ IN recvtype] data type of receive buffer elements (handle)
17559599516SKenneth E. Jansenc[ IN comm] communicator (handle)
17659599516SKenneth E. Jansenc The displacement array disp is an array of integers in which the jth
17759599516SKenneth E. Jansenc entry indicates which slot of xw marks beginning of xwj
17859599516SKenneth E. Jansenc So, first we will build this displacement array
17959599516SKenneth E. Jansen            idisp(:)=0 ! starting with zero, since MPI likes C-numbering
18059599516SKenneth E. Jansen            do j=2,numpe
18159599516SKenneth E. Jansen               idisp(j)=idisp(j-1)+nwall(j-1) ! see diagram above
18259599516SKenneth E. Jansen            enddo
18359599516SKenneth E. Jansenc Now, we gather the data one slice at a time (1:nwalli)
18459599516SKenneth E. Jansen            do j=1,nenb+1
18559599516SKenneth E. Jansen               do k=1,nsd
18659599516SKenneth E. Jansen                  call MPI_ALLGATHERV(xwi(:,j,k),nwalli,
18759599516SKenneth E. Jansen     &                 MPI_DOUBLE_PRECISION,xw(:,j,k),nwall,idisp,
18859599516SKenneth E. Jansen     &                 MPI_DOUBLE_PRECISION,MPI_COMM_WORLD,ierr)
18959599516SKenneth E. Jansen               enddo
19059599516SKenneth E. Jansen            enddo
19159599516SKenneth E. Jansen         else                   ! single-processor
19259599516SKenneth E. Jansenc global data is local data in single processor case
19359599516SKenneth E. Jansen            xw=xwi
19459599516SKenneth E. Jansen         endif
19559599516SKenneth E. Jansen
19659599516SKenneth E. Jansenc
19759599516SKenneth E. Jansenc  For each point, loop over wall nodes and calculate distance;
19859599516SKenneth E. Jansenc  save the distance in this node's position of d2wall if it's
19959599516SKenneth E. Jansenc  shorter than currently stored distance
20059599516SKenneth E. Jansenc
20159599516SKenneth E. Jansen         d2wall=1.0e32
20259599516SKenneth E. Jansen         do i=1,numnp
20359599516SKenneth E. Jansen            do j=1, nwallt
20459599516SKenneth E. Jansen               do k=1,nenb+1
20559599516SKenneth E. Jansen                  distance =  ( x(i,1) - xw(j,k,1) )**2
20659599516SKenneth E. Jansen     &                 +( x(i,2) - xw(j,k,2) )**2
20759599516SKenneth E. Jansen     &                 +( x(i,3) - xw(j,k,3) )**2
20859599516SKenneth E. Jansen                  if ( d2wall(i).gt.distance ) d2wall(i) = distance
20959599516SKenneth E. Jansen               enddo
21059599516SKenneth E. Jansen            enddo
21159599516SKenneth E. Jansen         enddo
21259599516SKenneth E. Jansen         d2wall=sqrt(d2wall)
21359599516SKenneth E. Jansenc
21459599516SKenneth E. Jansen         deallocate(xwi)
21559599516SKenneth E. Jansen         deallocate(xw)
21659599516SKenneth E. Jansenc
21759599516SKenneth E. Jansenc.... write d2wall to a file so we don't have to do this again
21859599516SKenneth E. Jansenc
21959599516SKenneth E. Jansen
22059599516SKenneth E. Jansen!MR CHANGE #Fix this with SyncIO!!!!!
22159599516SKenneth E. Jansen!         write (fmt1,"('(''d2wall.'',i',i1,',1x)')") 1
22259599516SKenneth E. Jansen!         write (fname1,fmt1) 0
22359599516SKenneth E. Jansen!         fname1 = trim(fname1) // cname(myrank+1)
22459599516SKenneth E. Jansen!         if(myrank.eq.master) write (*,*) 'Writing dist file : ', fname1
22559599516SKenneth E. Jansen!         open (unit=72, file=fname1, status='unknown',
22659599516SKenneth E. Jansen!     &                                    form='unformatted', err=996)
22759599516SKenneth E. Jansen!
22859599516SKenneth E. Jansen!         write (72) d2wall
22959599516SKenneth E. Jansen!         close (72)
23059599516SKenneth E. Jansen
23159599516SKenneth E. Jansen!MR CHANGE END
23259599516SKenneth E. Jansen
23359599516SKenneth E. Jansenc         write(*,*) "make sure to: echo 0 > distcalc.dat"
23459599516SKenneth E. Jansenc         call MPI_BARRIER (MPI_COMM_WORLD,ierr)
23559599516SKenneth E. Jansenc         call error ('distcalc','complete',0)
23659599516SKenneth E. Jansen      endif
23759599516SKenneth E. Jansen
23859599516SKenneth E. Jansen!MR CHANGE
23959599516SKenneth E. Jansen      call MPI_BARRIER(MPI_COMM_WORLD,ierr)
24059599516SKenneth E. Jansen      if(myrank.eq.master) then
24159599516SKenneth E. Jansen        write (*,*) 'leaving initTurb'
24259599516SKenneth E. Jansen      endif
24359599516SKenneth E. Jansen!MR CHANGE
24459599516SKenneth E. Jansen
24559599516SKenneth E. Jansen      return
24659599516SKenneth E. Jansen995     call error ('d2wall  ','opening ', 72)
24759599516SKenneth E. Jansen996     call error ('d2wall  ','opening ', 72)
24859599516SKenneth E. Jansen
24959599516SKenneth E. Jansen      end subroutine
25059599516SKenneth E. Jansen
25159599516SKenneth E. Jansen
252