copyright (C) 2001 MSC-RPN COMM %%%MC2%%% ***s/r exhrdif2 *subroutine exhrdif2 (t1,t2_8) 1,6 implicit none * real t1(*) real*8 t2_8(*) * *AUTHORs C. Girard & M. Desgagne * *OBJECT * * ******************************************************************* * * * * * perform EXPLICIT HORIZONTAL DIFFUSION * * * * * * * * * of the type * * * * * * * * * ---------------- * * * | n | * * * | -(-nu_DEL2) | * * * | | * * * ---------------- * * * * * * * * * calling SUBROUTINE nu_del2_n n times * * * * * * * * ******************************************************************* * ** *EXTERNALS #include "dynmem.cdk"
#include "yomdyn1.cdk"
#include "levels.cdk"
#include "partopo.cdk"
#include "topo.cdk"
#include "consdyn_8.cdk"
#include "dynpar.cdk"
* integer i,j,n,k,mm,nn,ng real*8 visco(5),zero,pt25,pt5,psi,con1,con2,v1,v2, $ recv1(minx:ldni+hx,miny:ldnj+hy) pointer (pav1, v1(minx:ldni+hx,miny:ldnj+hy,*)), $ (pav2, v2(minx:ldni+hx,miny:ldnj+hy,*)) * real smuup,smvvp,smwwp,smbbp,smhup,smtrp pointer (pau, smuup(*)),(pav, smvvp(*)),(paw, smwwp(*)), $ (pab, smbbp(*)),(pah, smhup(*)),(pat, smtrp(*)) * parameter(zero=0.d0,pt25=0.25d0,pt5=0.5d0) *---------------------------------------------------------------------- * pau = loc(t1( 1)) pav = loc(t1(1*dim3d+1)) paw = loc(t1(2*dim3d+1)) pab = loc(t1(3*dim3d+1)) pah = loc(t1(4*dim3d+1)) if (ntr.gt.0) pat = loc(t1(5*dim3d+1)) pav1 = loc(t2_8( 1)) pav2 = loc(t2_8(dim3d+1)) * if (kh(1).lt.1.0e-10) return * * * Smoothing type: -(-del2)**n; n=1,2,3..., hor_del=2,4,6,... * visco(1)=min( hord_fuv * dble(kh(1)) , pt25 ) visco(2)=min( hord_fww * dble(kh(1)) , pt25 ) visco(3)=min( hord_ftt * dble(kh(1)) , pt25 ) visco(4)=min( hord_fhu * dble(kh(1)) , pt25 ) visco(5)=min( hord_ftr * dble(kh(1)) , pt25 ) * if(hord_ftt.ne.0.) then * * diffusion of potential temperature * ng = (ldni+2*hx)*(ldnj+2*hy) con1= grav_8/grtstar con2= gama_star*pt5 !$omp do do k = 1, gnk do j = 1-hy,ldnj+hy do i = 1-hx,ldni+hx v1(i,j,k) = -gama_star*ht(i,j,k) end do end do call vexp (recv1, v1(minx,miny,k), ng) do j = 1-hy,ldnj+hy do i = 1-hx,ldni+hx v1(i,j,k) = con1*recv1(i,j) v2(i,j,k) = con2*(ppp(i,j,k)+ppp(i,j,k-1)) end do end do call vrec (recv1, v1(minx,miny,k), ng) do j = 1-hy,ldnj+hy do i = 1-hx,ldni+hx bbp(i,j,k) = (grav_8 + bbp(i,j,k) - v2(i,j,k))*recv1(i,j) end do end do end do !$omp end do * endif * !$omp single call rpn_comm_xch_halo (uup,minx,maxx,miny,maxy,ldni,ldnj, $ (ndynvar-1)*gnk,hx,hy,period_x,period_y,ldni,0) !$omp end single * nn=hord_del/2 * do mm=1,nn * if(visco(1).ne.zero) $ call nu_del2_n
(uup,smuup,minx,maxx,miny,maxy,gnk-1, $ visco(1),1,0,mm,nn) if(visco(1).ne.zero) $ call nu_del2_n
(vvp,smvvp,minx,maxx,miny,maxy,gnk-1, $ visco(1),0,1,mm,nn) if(visco(2).ne.zero) $ call nu_del2_n
(wwp,smwwp,minx,maxx,miny,maxy,gnk-1, $ visco(2),0,0,mm,nn) if(visco(3).ne.zero) $ call nu_del2_n
(bbp,smbbp,minx,maxx,miny,maxy,gnk-1, $ visco(3),0,0,mm,nn) if(visco(4).ne.zero) $ call nu_del2_n
(hup,smhup,minx,maxx,miny,maxy,gnk-1, $ visco(4),0,0,mm,nn) do n=1,ntr if(visco(5).ne.zero) $ call nu_del2_n
(trp(1-hx,1-hy,1,n),smtrp((n-1)*dim3d+1), $ minx,maxx,miny,maxy,gnk-1,visco(5),0,0,mm,nn) end do * if(mm.ne.nn) then !$omp single call rpn_comm_xch_halo (smuup,minx,maxx,miny,maxy,ldni,ldnj, $ (5+ntr)*gnk,hx,hy,period_x,period_y,ldni,0) !$omp end single endif * end do * if(hord_ftt.ne.0.) then !$omp do do k = 1, gnk do j = 1-hy,ldnj+hy do i = 1-hx,ldni+hx bbp(i,j,k) = bbp(i,j,k)*v1(i,j,k) - grav_8 + v2(i,j,k) end do end do end do !$omp end do endif * *---------------------------------------------------------------------- return end