subroutine cnuity(m,n) use mod_xc ! HYCOM communication interface use mod_pipe ! HYCOM debugging interface use mod_floats ! HYCOM synthetic floats, drifters and moorings #if defined(STOKES) use mod_stokes ! Stokes Drift Velocity Module #endif c c --- hycom version 1.0 implicit none c include 'common_blocks.h' c integer m,n c c --- ------------------------------------------------------ c --- continuity equation (flux-corrected transport version) c --- ------------------------------------------------------ c --- on entry: c --- dp( :,:,:,n) = time step t-1 c --- dp( :,:,:,m) = time step t c c --- dpv( :,:,:,n) = time step t-1 c --- dpv( :,:,:,m) = time step t c --- dpu( :,:,:,n) = time step t-1 c --- dpu( :,:,:,m) = time step t c c --- on exit: c --- dpo(:,:,:,n) = time step t-1 c --- dpo(:,:,:,m) = time step t c --- dp( :,:,:,m) = time step t with RA time smoothing c --- dp( :,:,:,n) = time step t+1 c --- ------------------------------------------------------ c real dpfatal parameter (dpfatal=-10.0) !fatal negative dp in meters c logical lpipe_cnuity parameter (lpipe_cnuity=.false.) c integer, save, dimension (1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy) :: & masku,maskv real, save, dimension (1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy) :: & pold real, save, dimension (1-nbdy:jdm+nbdy) :: & dpmn c integer i,iflip,iprint,isave,j,jsave,k,l,ia,ib,ja,jb,mbdy real q,dpmin,clip,flxhi,flxlo,dtinv,dpup,dpdn,thkdfu,thkdfv real dpkmin(2*kdm) c character*12 text,textu,textv c mbdy = 6 c call xctilr(dpmixl( 1-nbdy,1-nbdy, n),1, 1, 6,6, halo_ps) call xctilr(dp( 1-nbdy,1-nbdy,1,n),1,kk, 6,6, halo_ps) call xctilr(dp( 1-nbdy,1-nbdy,1,m),1,kk, 6,6, halo_ps) call xctilr(dpu( 1-nbdy,1-nbdy,1,m),1,kk, 6,6, halo_us) call xctilr(dpv( 1-nbdy,1-nbdy,1,m),1,kk, 6,6, halo_vs) call xctilr(u( 1-nbdy,1-nbdy,1,m),1,kk, 6,6, halo_uv) call xctilr(v( 1-nbdy,1-nbdy,1,m),1,kk, 6,6, halo_vv) call xctilr(ubavg( 1-nbdy,1-nbdy, m),1, 1, 6,6, halo_uv) call xctilr(vbavg( 1-nbdy,1-nbdy, m),1, 1, 6,6, halo_vv) c c --- rhs: dpmixl.n c --- lhs: util3, dpold, utotn, vtotn c margin = mbdy c !$OMP PARALLEL DO PRIVATE(j,l,i) !$OMP& SCHEDULE(STATIC,jblk) do j=1-margin,jj+margin do l=1,isp(j) do i=max(1-margin,ifp(j,l)),min(ii+margin,ilp(j,l)) util3(i,j)=0. dpmold( i,j)=dpmixl(i,j,n) ! save for Asselin filter enddo enddo do l=1,isu(j) do i=max(1-margin,ifu(j,l)),min(ii+margin,ilu(j,l)) utotn(i,j)=0. enddo enddo do l=1,isv(j) do i=max(1-margin,ifv(j,l)),min(ii+margin,ilv(j,l)) vtotn(i,j)=0. enddo enddo enddo !$OMP END PARALLEL DO do 76 k=1,kk c c --- uflux/vflux = low-order (diffusive) mass fluxes at old time level. c --- uflux2/vflux2 = 'antidiffusive' fluxes, defined as high-order minus low- c --- order fluxes. high-order fluxes are second-order in space, time-centered. c c --- rhs: depthu+, util3, dp.n, ubavg.m c --- lhs: uflux c margin = mbdy - 1 c do j=1-margin,jj+margin do l=1,isu(j) i=ifu(j,l)-1 if (i.ge.1-margin) then if (iuopn(i,j).ne.0) then q=min(dp(i ,j,k,n),max(0.,depthu(i+1,j)-util3(i ,j))) utotm(i,j)=(u(i+1,j,k,m)+ubavg(i,j,m))*scuy(i,j) #if defined(STOKES) utotm(i,j)=utotm(i,j)+usd(i+1,j,k)*scuy(i,j) #endif uflux(i,j)=utotm(i,j)*q endif endif i=ilu(j,l)+1 if (i.le.ii+margin) then if (iuopn(i,j).ne.0) then q=min(dp(i-1,j,k,n),max(0.,depthu(i-1,j)-util3(i-1,j))) utotm(i,j)=(u(i-1,j,k,m)+ubavg(i,j,m))*scuy(i,j) #if defined(STOKES) utotm(i,j)=utotm(i,j)+usd(i-1,j,k)*scuy(i,j) #endif uflux(i,j)=utotm(i,j)*q endif endif enddo enddo c c --- rhs: depthv+, util3, dp.n, vbavg.m c --- lhs: vflux c margin = mbdy - 1 c do i=1-margin,ii+margin do l=1,jsv(i) j=jfv(i,l)-1 if (j.ge.1-margin) then if (ivopn(i,j).ne.0) then q=min(dp(i,j ,k,n),max(0.,depthv(i,j+1)-util3(i,j ))) vtotm(i,j)=(v(i,j+1,k,m)+vbavg(i,j,m))*scvx(i,j) #if defined(STOKES) vtotm(i,j)=vtotm(i,j)+vsd(i,j+1,k)*scvx(i,j) #endif vflux(i,j)=vtotm(i,j)*q endif endif j=jlv(i,l)+1 if (j.le.jj+margin) then if (ivopn(i,j).ne.0) then q=min(dp(i,j-1,k,n),max(0.,depthv(i,j-1)-util3(i,j-1))) vtotm(i,j)=(v(i,j-1,k,m)+vbavg(i,j,m))*scvx(i,j) #if defined(STOKES) vtotm(i,j)=vtotm(i,j)+vsd(i,j-1,k)*scvx(i,j) #endif vflux(i,j)=vtotm(i,j)*q endif endif enddo enddo c c --- rhs: u.m, ubavg.m, depthu, dp.n+, util3+, dpu.m, uflux c --- rhs: v.m, vbavg.m, depthv, dp.n+, util3+, dpv.m, vflux c --- lhs: utotm,uflux,uflux2,uflx c --- lhs: vtotm,vflux,vflux2,vflx c margin = mbdy - 1 c !$OMP PARALLEL DO PRIVATE(j,l,i,q) !$OMP& SCHEDULE(STATIC,jblk) do j=1-margin,jj+margin c do l=1,isu(j) do i=max(1-margin,ifu(j,l)),min(ii+margin,ilu(j,l)) utotm(i,j)=(u(i,j,k,m)+ubavg(i,j,m))*scuy(i,j) #if defined(STOKES) utotm(i,j)=utotm(i,j)+usd(i,j,k)*scuy(i,j) #endif if (utotm(i,j).ge.0.) then q=min(dp(i-1,j,k,n),max(0.,depthu(i,j)-util3(i-1,j))) else q=min(dp(i ,j,k,n),max(0.,depthu(i,j)-util3(i ,j))) endif uflux(i,j)=utotm(i,j)*q uflux2(i,j)=utotm(i,j)*dpu(i,j,k,m)-uflux(i,j) uflx(i,j,k)=uflux(i,j) enddo enddo c do l=1,isv(j) do i=max(1-margin,ifv(j,l)),min(ii+margin,ilv(j,l)) c vtotm(i,j)=(v(i,j,k,m)+vbavg(i,j,m))*scvx(i,j) #if defined(STOKES) vtotm(i,j)=vtotm(i,j)+vsd(i,j,k)*scvx(i,j) #endif if (vtotm(i,j).ge.0.) then q=min(dp(i,j-1,k,n),max(0.,depthv(i,j)-util3(i,j-1))) else q=min(dp(i,j ,k,n),max(0.,depthv(i,j)-util3(i,j ))) endif vflux(i,j)=vtotm(i,j)*q vflux2(i,j)=vtotm(i,j)*dpv(i,j,k,m)-vflux(i,j) vflx(i,j,k)=vflux(i,j) enddo enddo enddo !$OMP END PARALLEL DO c c --- advance -dp- field using low-order (diffusive) flux values c --- rhs: dp.n, dp.m, util3, uflux+, vflux+ c --- lhs: dpold,dpoldm,util3,dp.n c margin = mbdy - 2 c !$OMP PARALLEL DO PRIVATE(j,l,i,dpmin) !$OMP& SCHEDULE(STATIC,jblk) do j=1-margin,jj+margin dpmin=999. do l=1,isp(j) do i=max(1-margin,ifp(j,l)),min(ii+margin,ilp(j,l)) dpold( i,j,k)=dp(i,j,k,n) util3(i,j)=util3(i,j)+dp(i,j,k,n) dp(i,j,k,n)=dp(i,j,k,n)- & ((uflux(i+1,j)-uflux(i,j))+ & (vflux(i,j+1)-vflux(i,j)))*delt1*scp2i(i,j) dpoldm(i,j,k)=dp(i,j,k,n) ! save for loop 19 test dpmin=min(dpmin,dp(i,j,k,n)) enddo enddo dpmn(j)=dpmin ! minimizes false sharing enddo ! loop 19 !$OMP END PARALLEL DO c dpmin=999. do j=1,jj dpmin=min(dpmin,dpmn(j)) enddo dpkmin(k)=dpmin c if (lpipe .and. lpipe_cnuity) then c --- compare two model runs. write (text,'(a9,i3)') 'dp.low k=',k call pipe_compare_sym1(dp(1-nbdy,1-nbdy,k,n),ip,text) endif c do j=1-margin,jj+margin do l=1,isu(j) i=ifu(j,l)-1 if (i.ge.1-margin) then if (iuopn(i,j).ne.0) then uflux(i,j)=0.0 endif endif i=ilu(j,l)+1 if (i.le.ii+margin) then if (iuopn(i,j).ne.0) then uflux(i,j)=0.0 endif endif enddo enddo c do i=1-margin,ii+margin do l=1,jsv(i) j=jfv(i,l)-1 if (j.ge.1-margin) then if (ivopn(i,j).ne.0) then vflux(i,j)=0.0 endif endif j=jlv(i,l)+1 if (j.le.jj+margin) then if (ivopn(i,j).ne.0) then vflux(i,j)=0.0 endif endif enddo enddo c if (lpipe .and. lpipe_cnuity) then c --- compare two model runs. do j=1-margin,jj+margin do i=1-margin,ii+margin masku(i,j)=iu(i,j) if (i.gt. 1) masku(i,j)=masku(i,j)+iu(i-1,j) if (i.lt.ii) masku(i,j)=masku(i,j)+iu(i+1,j) maskv(i,j)=iv(i,j) if (j.gt. 1) maskv(i,j)=maskv(i,j)+iv(i,j-1) if (j.lt.jj) maskv(i,j)=maskv(i,j)+iv(i,j+1) enddo enddo write (textu,'(a9,i3)') 'uflux k=',k write (textv,'(a9,i3)') 'vflux k=',k call pipe_compare_sym2(uflux,masku,textu, & vflux,maskv,textv) write (textu,'(a9,i3)') 'uflux2 k=',k write (textv,'(a9,i3)') 'vflux2 k=',k call pipe_compare_sym2(uflux2,masku,textu, & vflux2,maskv,textv) endif c cdiag if (mod(k,15).eq.1) then cdiag do i=itest-1,itest+1 cdiag do j=jtest-1,jtest+1 cdiag write (lp,101) nstep,i+i0,j+j0,k,dpold(i-1,j,k),uflux(i,j), cdiag. 'old dp''s, fluxes:',dpold(i,j-1,k),dpold(i,j,k),dpold(i,j+1,k) cdiag. ,vflux(i,j),dp(i,j,k,n),vflux(i,j+1),dpold(i+1,j,k),uflux(i+1,j) cdiag enddo cdiag enddo cdiag endif 101 format (i9,2i5,i3,1p,e15.2,e30.2/a17,6e10.2/e37.2,e30.2) c c --- at each grid point, determine the ratio of the largest permissible c --- pos. (neg.) change in -dp- to the sum of all incoming (outgoing) fluxes c c --- rhs: dp.n+, uflux2+, vflux2+ c --- lhs: util1,util2 c margin = mbdy - 2 c !$OMP PARALLEL DO PRIVATE(j,l,i,ia,ib,ja,jb) !$OMP& SCHEDULE(STATIC,jblk) do j=1-margin,jj+margin do l=1,isp(j) do i=max(1-margin,ifp(j,l)),min(ii+margin,ilp(j,l)) c --- assume margin