xref: /phasta/phSolver/common/settauw.f (revision 595995161822a203c8467e0e4a253d7bd7d6df32)
1*59599516SKenneth E. Jansen        subroutine settauw (y,  x,  BC,
2*59599516SKenneth E. Jansen     &                      ifath,   velbar)
3*59599516SKenneth E. Jansenc
4*59599516SKenneth E. Jansenc----------------------------------------------------------------------
5*59599516SKenneth E. Jansenc
6*59599516SKenneth E. Jansenc This routine computes the time varying viscous flux for turbulence wall
7*59599516SKenneth E. Jansenc  boundary elements.
8*59599516SKenneth E. Jansenc
9*59599516SKenneth E. Jansenc Zdenek Johan, Winter 1991.  (Fortran 90)
10*59599516SKenneth E. Jansenc----------------------------------------------------------------------
11*59599516SKenneth E. Jansenc
12*59599516SKenneth E. Jansen      use pointer_data
13*59599516SKenneth E. Jansen      use turbSA
14*59599516SKenneth E. Jansen      include "common.h"
15*59599516SKenneth E. Jansenc
16*59599516SKenneth E. Jansen      dimension y(nshg,ndofl),            x(numnp,nsd),
17*59599516SKenneth E. Jansen     &          BC(nshg,ndofBC),
18*59599516SKenneth E. Jansen     &          ifath(numnp),             velbar(nfath,nflow),
19*59599516SKenneth E. Jansen     &          ull(nsd),                 trx(numnp,nsd),
20*59599516SKenneth E. Jansen     &          ullb(nsd),                dull(nsd),
21*59599516SKenneth E. Jansen     &          evisc(numnp)
22*59599516SKenneth E. Jansen      real*8 outvec(nshg,nsd+nsd+nsd)
23*59599516SKenneth E. Jansenc
24*59599516SKenneth E. Jansenc calculate the traction magnitude for all nodes on the wall
25*59599516SKenneth E. Jansenc
26*59599516SKenneth E. Jansen      trx=zero
27*59599516SKenneth E. Jansen      evisc=zero
28*59599516SKenneth E. Jansen      rm=datmat(1,2,1)
29*59599516SKenneth E. Jansen      do nodw = 1, numnp ! loop over nodes
30*59599516SKenneth E. Jansen      if ( otwn(nodw).ne.nodw ) then ! wall node check
31*59599516SKenneth E. Jansen         if (iLES.gt.0) then
32*59599516SKenneth E. Jansen            ull(:)=velbar(ifath(otwn(nodw)),2:4) ! recall old #ing
33*59599516SKenneth E. Jansen         endif
34*59599516SKenneth E. Jansen         if (iRANS.lt.0) then
35*59599516SKenneth E. Jansen            ull(:)=y(otwn(nodw),1:3)
36*59599516SKenneth E. Jansen         endif
37*59599516SKenneth E. Jansen         ub=ull(1)*wnrm(nodw,1) !
38*59599516SKenneth E. Jansen     &     +ull(2)*wnrm(nodw,2) ! store u.n here for now
39*59599516SKenneth E. Jansen     &     +ull(3)*wnrm(nodw,3) !
40*59599516SKenneth E. Jansenc     u_parallel_to_boundary=u-(u.n)n  :
41*59599516SKenneth E. Jansen         ull(:)=ull(:)-ub*wnrm(nodw,:) ! ull in flow
42*59599516SKenneth E. Jansenc     u_b is || u_ll ||
43*59599516SKenneth E. Jansen         ub=sqrt(ull(1)**2+ull(2)**2+ull(3)**2)
44*59599516SKenneth E. Jansenc     perp d2wall is (x2-x1).n   :
45*59599516SKenneth E. Jansen         dw=abs(
46*59599516SKenneth E. Jansen     &         (x(otwn(nodw),1)-x(nodw,1))*wnrm(nodw,1)
47*59599516SKenneth E. Jansen     &        +(x(otwn(nodw),2)-x(nodw,2))*wnrm(nodw,2)
48*59599516SKenneth E. Jansen     &        +(x(otwn(nodw),3)-x(nodw,3))*wnrm(nodw,3)
49*59599516SKenneth E. Jansen     &         )
50*59599516SKenneth E. Jansen         if(ub.eq.0) then
51*59599516SKenneth E. Jansen            ut=zero
52*59599516SKenneth E. Jansen            twoub=zero
53*59599516SKenneth E. Jansen         else
54*59599516SKenneth E. Jansen            ut=utau(ub,dw,rm,nodw,x(nodw,:))   ! find utau
55*59599516SKenneth E. Jansen            twoub=-ut*ut/ub     ! find -tau_w/ub
56*59599516SKenneth E. Jansen         endif
57*59599516SKenneth E. Jansenc
58*59599516SKenneth E. Jansenc for LES ull has the mean-parallel velocity vector.  We want the
59*59599516SKenneth E. Jansenc instantaneous-parallel velocity vector
60*59599516SKenneth E. Jansenc
61*59599516SKenneth E. Jansen         if (iLES.gt.0) then
62*59599516SKenneth E. Jansen            ullb=ull
63*59599516SKenneth E. Jansen            ull(:)=y(otwn(nodw),1:3)
64*59599516SKenneth E. Jansen            ubn=ull(1)*wnrm(nodw,1) !
65*59599516SKenneth E. Jansen     &        +ull(2)*wnrm(nodw,2) ! store u.n here for now
66*59599516SKenneth E. Jansen     &        +ull(3)*wnrm(nodw,3) !
67*59599516SKenneth E. Jansenc
68*59599516SKenneth E. Jansenc     u_parallel_to_boundary=u-(u.n)n  :
69*59599516SKenneth E. Jansenc
70*59599516SKenneth E. Jansen            ull(:)=ull(:)-ubn*wnrm(nodw,:) ! ull in flow
71*59599516SKenneth E. Jansenc
72*59599516SKenneth E. Jansenc hack a limiter into this fluctuating vector. Early transients have
73*59599516SKenneth E. Jansenc huge differences from mean values
74*59599516SKenneth E. Jansenc
75*59599516SKenneth E. Jansen            dull=ull-ullb ! the current vector difference
76*59599516SKenneth E. Jansen            dullm=sqrt(dull(1)*dull(1)+dull(2)*dull(2)+dull(3)*dull(3))
77*59599516SKenneth E. Jansen     &           + 1.0e-9
78*59599516SKenneth E. Jansen            ullbm=ub ! mag of ullb still there
79*59599516SKenneth E. Jansenc
80*59599516SKenneth E. Jansenc limit the magnitude of the difference to a 40% change from the mean.
81*59599516SKenneth E. Jansenc if less than that already we will take the whole difference, otherwise
82*59599516SKenneth E. Jansenc only take a 40% change.
83*59599516SKenneth E. Jansenc
84*59599516SKenneth E. Jansen
85*59599516SKenneth E. Jansen            dullmod=min(one,0.4*ullbm/dullm)
86*59599516SKenneth E. Jansen            ull=ullb+dullmod*dull
87*59599516SKenneth E. Jansen         endif
88*59599516SKenneth E. Jansen
89*59599516SKenneth E. Jansen         trx(nodw,:)=twoub*ull(:)
90*59599516SKenneth E. Jansen         if(itwmod.eq.-2) then ! effective-viscosity
91*59599516SKenneth E. Jansen            tauw=ut*ut
92*59599516SKenneth E. Jansen            BC(nodw,7)=tauw*dw/ub-rm
93*59599516SKenneth E. Jansen         endif
94*59599516SKenneth E. Jansen         if(itwmod.eq.2) then ! effective-viscosity
95*59599516SKenneth E. Jansenc
96*59599516SKenneth E. Jansenc  mag of u instantaneous
97*59599516SKenneth E. Jansenc
98*59599516SKenneth E. Jansen            ullm=sqrt(ull(1)*ull(1)+ull(2)*ull(2)+ull(3)*ull(3))
99*59599516SKenneth E. Jansen            tauw=ut*ut*ullm/ub
100*59599516SKenneth E. Jansen            evisc(nodw)=tauw*dw/ub-rm
101*59599516SKenneth E. Jansen         endif
102*59599516SKenneth E. Jansen         if((itwmod.eq.-1)) then ! slip-velocity RANS
103*59599516SKenneth E. Jansen            up=sqrt(
104*59599516SKenneth E. Jansen     &           y(nodw,1)**2   !
105*59599516SKenneth E. Jansen     &           +y(nodw,2)**2  ! flow is auto-|| at boundary
106*59599516SKenneth E. Jansen     &           +y(nodw,3)**2  !
107*59599516SKenneth E. Jansen     &           )/ut
108*59599516SKenneth E. Jansen            BC(nodw,7)=savarw(up,rm,saCv1,nodw,x(nodw,:))
109*59599516SKenneth E. Jansen         endif
110*59599516SKenneth E. Jansen      endif ! wallnode check
111*59599516SKenneth E. Jansen      enddo ! loop over nodes
112*59599516SKenneth E. Jansenc
113*59599516SKenneth E. Jansenc Write the traction vectors to a file restart.4077.n
114*59599516SKenneth E. Jansenc
115*59599516SKenneth E. Jansenc$$$      ilstep=4077
116*59599516SKenneth E. Jansenc$$$      outvec(:,1:3)=trx(:,1:3)
117*59599516SKenneth E. Jansenc$$$      outvec(:,4:6)=0
118*59599516SKenneth E. Jansenc$$$      do i=1,numnp
119*59599516SKenneth E. Jansenc$$$         if(otwn(i).ne. i ) outvec(i,4:6)=y(otwn(i),1:3)
120*59599516SKenneth E. Jansenc$$$      enddo
121*59599516SKenneth E. Jansenc$$$      outvec(:,7:9)=wnrm(:,1:3)
122*59599516SKenneth E. Jansenc$$$      call write_restart(myrank,ilstep,numnp,nsd*3,outvec,outvec)
123*59599516SKenneth E. Jansenc$$$      write(*,*) 'Traction dumped to restart.4077.*'
124*59599516SKenneth E. Jansenc
125*59599516SKenneth E. Jansenc Put traction calculations into BCB
126*59599516SKenneth E. Jansenc
127*59599516SKenneth E. Jansen      do iblk = 1, nelblb
128*59599516SKenneth E. Jansen         iel    = lcblkb(1,iblk)
129*59599516SKenneth E. Jansen         nenl   = lcblkb(5,iblk) ! no. of vertices per element
130*59599516SKenneth E. Jansen         nenbl  = lcblkb(6,iblk) ! no. of vertices per bdry. face
131*59599516SKenneth E. Jansen         ndofl  = lcblkb(8,iblk)
132*59599516SKenneth E. Jansen         npro   = lcblkb(1,iblk+1) - iel
133*59599516SKenneth E. Jansenc For all elements in this block that lie on a wall, assign the traction
134*59599516SKenneth E. Jansen         do i=1,npro
135*59599516SKenneth E. Jansen            if(btest(miBCB(iblk)%p(i,1),4)) then ! wall elt
136*59599516SKenneth E. Jansen               do j = 1, nenbl
137*59599516SKenneth E. Jansen                  if(itwmod.eq.-2) then ! effective-viscosity
138*59599516SKenneth E. Jansen                     mBCB(iblk)%p(i,j,3:5)=0.0
139*59599516SKenneth E. Jansen                  endif
140*59599516SKenneth E. Jansen                  if((itwmod.eq.-1).or.(itwmod.eq.1)) then ! slip-velocity
141*59599516SKenneth E. Jansen                     mBCB(iblk)%p(i,j,3:5)=trx(mienb(iblk)%p(i,j),:)
142*59599516SKenneth E. Jansen                  endif
143*59599516SKenneth E. Jansen               enddo
144*59599516SKenneth E. Jansen            endif
145*59599516SKenneth E. Jansen         enddo
146*59599516SKenneth E. Jansen      enddo
147*59599516SKenneth E. Jansen
148*59599516SKenneth E. Jansen
149*59599516SKenneth E. Jansen      if(itwmod.eq.2) then   !effective viscosity
150*59599516SKenneth E. Jansenc
151*59599516SKenneth E. Jansenc For the elements which touch a modeled wall,
152*59599516SKenneth E. Jansenc modify the eddy viscosity at all quadrature points to be the average
153*59599516SKenneth E. Jansenc of the "optimal" nodal values.  This is an element-wise constant.
154*59599516SKenneth E. Jansenc
155*59599516SKenneth E. Jansen         do iblk = 1,nelblk
156*59599516SKenneth E. Jansen            nenl   = lcblk(5,iblk) ! no. of vertices per element
157*59599516SKenneth E. Jansen            iel    = lcblk(1,iblk)
158*59599516SKenneth E. Jansen            npro   = lcblk(1,iblk+1) - iel
159*59599516SKenneth E. Jansen            lcsyst = lcblk(3,iblk)
160*59599516SKenneth E. Jansen            ngauss = nint(lcsyst)
161*59599516SKenneth E. Jansen            do i=1,npro
162*59599516SKenneth E. Jansen               xnave=zero
163*59599516SKenneth E. Jansen               avevisc=zero
164*59599516SKenneth E. Jansen               do j=1,nenl
165*59599516SKenneth E. Jansen                  ev = evisc(mien(iblk)%p(i,j))
166*59599516SKenneth E. Jansen                  avevisc=avevisc + ev
167*59599516SKenneth E. Jansen                  if(ev.ne.zero) xnave=xnave+1
168*59599516SKenneth E. Jansen               enddo
169*59599516SKenneth E. Jansen               if(xnave.ne.0) mxmudmi(iblk)%p(i,1:ngauss)=avevisc/xnave
170*59599516SKenneth E. Jansen            enddo
171*59599516SKenneth E. Jansen         enddo
172*59599516SKenneth E. Jansen      endif
173*59599516SKenneth E. Jansenc
174*59599516SKenneth E. Jansenc.... end
175*59599516SKenneth E. Jansenc
176*59599516SKenneth E. Jansen        return
177*59599516SKenneth E. Jansen        end
178*59599516SKenneth E. Jansen
179*59599516SKenneth E. Jansen        function utaul(u,y,rm)
180*59599516SKenneth E. Jansen        real*8 utaul
181*59599516SKenneth E. Jansen        utaul=.2
182*59599516SKenneth E. Jansen        return
183*59599516SKenneth E. Jansen        end
184*59599516SKenneth E. Jansen
185*59599516SKenneth E. Jansen        function utau(u,y,rm,nodw,x)
186*59599516SKenneth E. Jansen        implicit none
187*59599516SKenneth E. Jansen        real*8 u,err,utau,yrmi,y,rm,yp,up,kup,f,dfds,rat,efac
188*59599516SKenneth E. Jansen        real*8 kup2,kup3,pt5,sxth,kappa
189*59599516SKenneth E. Jansen        integer iter
190*59599516SKenneth E. Jansen        integer nodw
191*59599516SKenneth E. Jansen        real*8 x(3)
192*59599516SKenneth E. Jansen        pt5=0.5
193*59599516SKenneth E. Jansen        sxth=0.1666666666666667
194*59599516SKenneth E. Jansen        err=1.0d-6
195*59599516SKenneth E. Jansen        utau=0.04
196*59599516SKenneth E. Jansen        yrmi=y/rm
197*59599516SKenneth E. Jansen        kappa=0.4
198*59599516SKenneth E. Jansenc$$$        B=5.5
199*59599516SKenneth E. Jansen        efac=0.1108 ! exp(-kappa*B)
200*59599516SKenneth E. Jansen        do iter=1,500
201*59599516SKenneth E. Jansen           yp=yrmi*utau
202*59599516SKenneth E. Jansen           up=u/utau
203*59599516SKenneth E. Jansen           kup=kappa*up
204*59599516SKenneth E. Jansen           kup2=kup*kup
205*59599516SKenneth E. Jansen           kup3=kup*kup2
206*59599516SKenneth E. Jansen           f=   up-yp+efac*(exp(kup)-1.0-kup - kup2*pt5 - kup3*sxth)
207*59599516SKenneth E. Jansen           dfds=up+yp+efac*(exp(kup)*kup-kup - kup2 - kup3*pt5) !/-utau in rat
208*59599516SKenneth E. Jansen           rat=f*utau/dfds ! this is -f/dfds but we did not do 1/(-utau) yet.
209*59599516SKenneth E. Jansen           utau=utau+rat
210*59599516SKenneth E. Jansen           if(abs(rat).le.err) goto 20
211*59599516SKenneth E. Jansen        enddo
212*59599516SKenneth E. Jansen        write(*,*)'utau failed to converge at ',nodw,x(1),x(2),x(3)
213*59599516SKenneth E. Jansen        write(*,*) 'u,          dwallperp,      mu'
214*59599516SKenneth E. Jansen        write(*,*) u, y, rm
215*59599516SKenneth E. Jansen        write(*,*) 'dfds,         rat,         utau'
216*59599516SKenneth E. Jansen        write(*,*) dfds,rat,utau
217*59599516SKenneth E. Jansenc        stop
218*59599516SKenneth E. Jansenc        utau=0.0
219*59599516SKenneth E. Jansenc
220*59599516SKenneth E. Jansenc  if the above fails then try a simple no-slip traction
221*59599516SKenneth E. Jansenc
222*59599516SKenneth E. Jansen        utau=sqrt(u/y*rm)
223*59599516SKenneth E. Jansen
224*59599516SKenneth E. Jansen 20     continue
225*59599516SKenneth E. Jansen        return
226*59599516SKenneth E. Jansen        end
227*59599516SKenneth E. Jansen
228*59599516SKenneth E. Jansenc$$$        function utau_log_layer_only(u,y,rm)
229*59599516SKenneth E. Jansenc$$$        real*8 u,err,utau,yrmi,y,rm,lnyp,f,dfds,rat
230*59599516SKenneth E. Jansenc$$$        err=1.0d-6
231*59599516SKenneth E. Jansenc$$$        utau=0.04
232*59599516SKenneth E. Jansenc$$$        yrmi=y/rm
233*59599516SKenneth E. Jansenc$$$        do iter=1,50
234*59599516SKenneth E. Jansenc$$$           lnyp=log(yrmi*utau)
235*59599516SKenneth E. Jansenc$$$           f=u-utau*(2.5*lnyp+5.2)
236*59599516SKenneth E. Jansenc$$$           dfds=-2.5*(lnyp+6.2)
237*59599516SKenneth E. Jansenc$$$           rat=-f/dfds
238*59599516SKenneth E. Jansenc$$$           utau=utau+rat
239*59599516SKenneth E. Jansenc$$$           if(utau.gt.0.5) then !flow is laminar for now give the parabolic
240*59599516SKenneth E. Jansenc$$$              utau=sqrt(2.0*rm*u)/y
241*59599516SKenneth E. Jansenc$$$              goto 20
242*59599516SKenneth E. Jansenc$$$           endif
243*59599516SKenneth E. Jansenc$$$
244*59599516SKenneth E. Jansenc$$$           if(abs(rat).le.err) goto 20
245*59599516SKenneth E. Jansenc$$$        enddo
246*59599516SKenneth E. Jansenc$$$        write(*,*)'utau failed to converge',r,dfds,rat,utau,y,u,rm
247*59599516SKenneth E. Jansenc$$$        stop
248*59599516SKenneth E. Jansenc$$$ 20     continue
249*59599516SKenneth E. Jansenc$$$        return
250*59599516SKenneth E. Jansenc$$$        end
251*59599516SKenneth E. Jansen
252*59599516SKenneth E. Jansen        function savarw(up,rm,cv1,nodw,x)
253*59599516SKenneth E. Jansen        implicit none
254*59599516SKenneth E. Jansen        real*8 err,savarw,rm,up,cv1,f,dfds,rat,efac
255*59599516SKenneth E. Jansen        real*8 pt5,kappa,B,xmut,chi3,denom,cv1_3
256*59599516SKenneth E. Jansen        integer iter
257*59599516SKenneth E. Jansen        integer nodw
258*59599516SKenneth E. Jansen        real*8 x(3)
259*59599516SKenneth E. Jansen        pt5=0.5
260*59599516SKenneth E. Jansen        err=1.0d-6
261*59599516SKenneth E. Jansen        savarw=rm*cv1*1.2599 ! inflection point chi=cv1*cuberoot(2)
262*59599516SKenneth E. Jansen        kappa=0.4
263*59599516SKenneth E. Jansenc$$$        B=5.5
264*59599516SKenneth E. Jansen        efac=0.1108 ! exp(-kappa*B)
265*59599516SKenneth E. Jansen        xmut=rm*kappa*efac*(exp(kappa*up)-1.0-kappa*up
266*59599516SKenneth E. Jansen     &       -pt5*(kappa**2)*(up**2))
267*59599516SKenneth E. Jansen        do iter=1,50
268*59599516SKenneth E. Jansen           chi3=savarw/rm
269*59599516SKenneth E. Jansen           chi3=chi3*chi3*chi3
270*59599516SKenneth E. Jansen           cv1_3=cv1**3
271*59599516SKenneth E. Jansen           denom=chi3+cv1_3
272*59599516SKenneth E. Jansen
273*59599516SKenneth E. Jansen           f=savarw*chi3/denom - xmut
274*59599516SKenneth E. Jansen           dfds=chi3*(chi3+4.0*cv1_3)/(denom**2)
275*59599516SKenneth E. Jansen           rat=-f/dfds
276*59599516SKenneth E. Jansen           savarw=savarw+rat
277*59599516SKenneth E. Jansen           if(abs(rat).le.err) goto 20
278*59599516SKenneth E. Jansen        enddo
279*59599516SKenneth E. Jansen        write(*,*)'savarw failed to converge at ',nodw,x(1),x(2),x(3)
280*59599516SKenneth E. Jansen        write(*,*) 'uplus,           mu'
281*59599516SKenneth E. Jansen        write(*,*) up,rm
282*59599516SKenneth E. Jansen        write(*,*) 'dfds,        rat,        savarw'
283*59599516SKenneth E. Jansen        write(*,*) dfds,rat,savarw
284*59599516SKenneth E. Jansen        savarw=10000.0*rm
285*59599516SKenneth E. Jansen 20     continue
286*59599516SKenneth E. Jansen        return
287*59599516SKenneth E. Jansen        end
288