copyright (C) 2001 MSC-RPN COMM %%%MC2%%% *subroutine rhs3_p_1d ( rhs, ww, nx, ny ) 1 implicit none * integer nx, ny real*8 rhs (nx, ny, *), ww (nx, ny, *) * *AUTHORs C. Girard & M. Desgagne * *OBJECT * ******************************************************************* * * * CALCULATES * * modified right-hand-side terms * * Qw*, Qp* * * FROM * * right-hand-side terms: * * Qu, Qv, Qw, Qb, Qp * * stored in * * up, vp, wp, bp, pp * * STORES * * results in * * wp=Qw*, rhs=Qp* * * * * PERFORMS essentially the following calculations: * * * * * * wp = [ wp + dtp * bp /(1+alfa) ] * nu0 * * * * ~ __Z * * rhs = - pp + dtp * [ d(wp)/dZ - gc02 * wp ] * * * * * * ~ * * SEE rhs3_pm for definitions of d/dZ * * * * * ******************************************************************* * #include "grd.cdk"
#include "dynmem.cdk"
#include "physcom.cdk"
#include "yomdyn.cdk"
#include "yomdyn1.cdk"
#include "dynpar.cdk"
#include "dtmdtp.cdk"
#include "nbcpu.cdk"
#include "alfmem.cdk"
* integer i, j, k, km1, kp1, nkt real*8 pt5, t1 parameter (pt5 = 0.5d0) * *----------------------------------------------------------------------- * nkt = gnk-1 if(flextop) nkt = gnk * do k = 1, gnk km1 = max(k-1, 1 ) do j = 1, ny do i = 1, nx wwp(i,j,k) = (wwp(i,j,k) + dtp*bbp(i,j,k)*ap1r(i,j,k)) $ * nu0b(i,j,k) ww (i,j,k) = wwp(i,j,k) end do end do if (k.eq.1) then do j = 1, ny do i = 1, nx t1 = sbxy(i,j) * pt5 * $ (-gg1(i+1,j,1)*uup(i+1,j,1) - gg1(i,j,1)*uup(i,j,1) $ - gg2(i,j+1,1)*vvp(i,j+1,1) - gg2(i,j,1)*vvp(i,j,1) ) ww(i,j,1) = dhdt(i,j,1) - t1 end do end do endif if ((k.eq.gnk).and.(.not.flextop)) then do j = 1, ny do i = 1, nx ww(i,j,gnk) = 0. end do end do endif end do * do k = 1, nkt kp1 = min( k+1 , gnk ) do j = 1, ny do i = 1, nx rhs(i,j,k) = - ppp(i,j,k) $ + dtp*gg0r(i,j,k)*(ww(i,j,kp1)-ww(i,j,k)) $ - dtp*pt5*gc02(i,j,k)*(ww(i,j,kp1)+ww(i,j,k)) end do end do end do * *----------------------------------------------------------------------- return end