!-------------------------------------- LICENCE BEGIN ------------------------------------ !Environment Canada - Atmospheric Science and Technology License/Disclaimer, ! version 3; Last Modified: May 7, 2008. !This is free but copyrighted software; you can use/redistribute/modify it under the terms !of the Environment Canada - Atmospheric Science and Technology License/Disclaimer !version 3 or (at your option) any later version that should be found at: !http://collaboration.cmc.ec.gc.ca/science/rpn.comm/license.html ! !This software is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; !without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. !See the above mentioned License/Disclaimer for more details. !You should have received a copy of the License/Disclaimer along with this software; !if not, you can write to: EC-RPN COMM Group, 2121 TransCanada, suite 500, Dorval (Quebec), !CANADA, H9P 1J3; or send e-mail to service.rpn@ec.gc.ca !-------------------------------------- LICENCE END -------------------------------------- ***s/r adw_main_3_intlag - Interpolation of rhs * #include "model_macros_f.h"*
subroutine adw_main_3_intlag ( F_u, F_v, F_w ) 1,13 * implicit none real F_u(*), F_v(*), F_w(*) * *author * alain patoine * *revision * v2_31 - Desgagne & Tanguay - removed stkmemw, introduce tracers * v2_31 - tracers not monotone if V4dg_conf.ne.0 * v2_31 - R.Moffet - do precalc (zz1) for ZZ in blomega * v3_00 - Desgagne & Lee - Lam configuration * v3_02 - Tanguay - Restore tracers monotone if V4dg_conf.ne.0 * v3_02 - Lee V. - revert adw_exch_1 for GLB only, * v3_02 added adw_ckbd_lam,adw_cfl_lam for LAM only * v3_03 - Tanguay M. - stop if adw_exch_1 is activated when 4D-Var * v3_10 - Corbeil & Desgagne & Lee - AIXport+Opti+OpenMP * v3_11 - Gravel S. - introduce key Adw_mono_L * v3_20 - Gravel & Valin & Tanguay - Lagrange 3D * v3_20 - Tanguay M. - Improve alarm when points outside advection grid * v3_20 - Dugas B. - correct calculation for LAM when Glb_pil gt 7 * v3_30 - McTaggart-Cowan R. - add adw_comp_cub_L * v3_30 - Tanguay M. - Validation for LAM version * v3_31 - Desgagne M. - new scope for operator + adw_cliptraj (LAM) * v3_31 - Tanguay M. - Introduce time extrapolation * v3_31 - Tanguay M. - Cancel 'Increase Halo' when G_lam * v3_31 - Tanguay M. - Do adjoint of outsiders in advection * *language * fortran 77 * *object * see id section * *arguments *______________________________________________________________________ * | | | * NAME | DESCRIPTION | I/O | *--------|-------------------------------------------------------|-----| * | | | * | | | * F_u,F_v| 3 components of upstream positions at t1 at input | iw | * F_w | used as work field afterward | | *________|_______________________________________________________|_____| * *implicits #include "glb_ld.cdk"
#include "glb_pil.cdk"
#include "ptopo.cdk"
#include "geomg.cdk"
#include "orh.cdk"
#include "lctl.cdk"
#include "step.cdk"
#include "lun.cdk"
#include "schm.cdk"
#include "offc.cdk"
#include "cstv.cdk"
#include "rhsc.cdk"
#include "v4dg.cdk"
#include "adw.cdk"
#include "tr3d.cdk"
#include "vt1.cdk"
#include "vt0.cdk"
#include "adw_comp.cdk"
* *modules integer vmmlod, vmmget, vmmuld external vmmlod, vmmget, vmmuld integer i0,in,j0,jn integer pnerr, pnlkey1(30), key1(Tr3d_ntr), key0(Tr3d_ntr), $ key1_, key0_, pnlod, err, dim * integer outside,sum_outside,ier * integer n, nij, nijk, nijkag, cnt, unf, i,j,k integer*8 pnt_trt1(Tr3d_ntr),pnt_trt0(Tr3d_ntr) * real*8 aaa real tr,tr0,buf(l_ni,G_nk) pointer (patr, tr(LDIST_SHAPE,*)),(patr0, tr0(LDIST_SHAPE,*)) * logical done_once_L save done_once_L data done_once_L /.false./ * * __________________________________________________________________ * if (Lun_debug_L) write (Lun_out,1000) nij = l_ni *l_nj nijk = l_ni *l_nj *l_nk nijkag = Adw_nit*Adw_njt*l_nk call hpalloc(Adw_capx1_ ,nijk, err,1) call hpalloc(Adw_capy1_ ,nijk, err,1) call hpalloc(Adw_capz1_ ,nijk, err,1) call hpalloc(Adw_n1_ ,nijk, err,1) call hpalloc(Adw_xdd1_ ,nijk, err,1) call hpalloc(Adw_xgg1_ ,nijk, err,1) call hpalloc(Adw_ydd1_ ,nijk, err,1) call hpalloc(Adw_ygg1_ ,nijk, err,1) call hpalloc(Adw_cz1_ ,nijk, err,1) call hpalloc(Adw_c1_ ,nijk, err,1) call hpalloc(Adw_wrkb_ ,nijk, err,1) call hpalloc(Adw_wrkc_ ,nijk, err,1) * ************************************************************************ i0=1 in=l_ni j0=1 jn=l_nj if (G_lam) then if (l_west) i0=pil_w if (l_east) in=l_niu - pil_e + 2 if (l_south) j0=pil_s if (l_north) jn=l_njv - pil_n + 2 if ( ((.not.Orh_crank_L).or.(Orh_icn.eq.Schm_itcn)) .and. $ (mod(Lctl_step,Step_gstat).eq.0) ) $ call adw_cfl_lam
( F_u, F_v, F_w, i0, in, j0, jn ) call adw_cliptraj
( F_u, F_v, i0, in, j0, jn, 'INTERP') * endif * if (.not.G_lam) then call adw_exch_1
( Adw_n1, Adw_xgg1, Adw_xdd1, Adw_c1, $ F_u, F_v, F_w ) * if ( V4dg_conf.ne.0.0 ) then * outside = 0 if ( Adw_fro_a .gt. 0 ) outside = 1 * sum_outside = 0 call rpn_comm_Allreduce(outside,sum_outside,1,"MPI_INTEGER", $ "MPI_SUM","grid",ier) * if(sum_outside.ne.0.and..not.done_once_L.and.Ptopo_myproc.eq.0) then write(Lun_out,*) 'NUMBER OF PE WITH OUTSIDERS IN ADW_MAIN_3_LAG AT current TIME-CN = ',sum_outside call flush(Lun_out) endif * endif * dim = max(1,Adw_fro_a) call hpalloc(Adw_capx2_ ,dim, err,1) call hpalloc(Adw_capy2_ ,dim, err,1) call hpalloc(Adw_capz2_ ,dim, err,1) call hpalloc(Adw_n2_ ,dim, err,1) call hpalloc(Adw_xdd2_ ,dim, err,1) call hpalloc(Adw_xgg2_ ,dim, err,1) call hpalloc(Adw_ydd2_ ,dim, err,1) call hpalloc(Adw_ygg2_ ,dim, err,1) call hpalloc(Adw_cz2_ ,dim, err,1) call hpalloc(Adw_wrka_ ,dim, err,1) * call adw_exch_2
( Adw_capx2, Adw_capy2, Adw_capz2, % Adw_n1, Adw_xgg1, Adw_xdd1, % Adw_fro_n, Adw_fro_s, Adw_fro_a, % Adw_for_n, Adw_for_s, Adw_for_a, 3 ) * if ( Adw_fro_a .gt. 0 .and. Adw_ckbd_L ) $ call adw_ckbd
( Adw_capy2 ) * endif * pnlkey1(1) = VMM_KEY(ruw1) pnlkey1(2) = VMM_KEY(rvw1) pnlkey1(3) = VMM_KEY(ruw2) pnlkey1(4) = VMM_KEY(rvw2) pnlkey1(5) = VMM_KEY(rcn) pnlkey1(6) = VMM_KEY(rth) pnlkey1(7) = VMM_KEY(fit1) pnlkey1(8) = VMM_KEY(zz1) pnlod = 8 * if (.not. Schm_hydro_L) then pnlkey1(9) = VMM_KEY(rw) pnlkey1(10) = VMM_KEY(rvv) pnlod = 10 endif * pnerr = vmmlod(pnlkey1,pnlod) * pnerr = VMM_GET_VAR(ruw1) pnerr = VMM_GET_VAR(rvw1) pnerr = VMM_GET_VAR(ruw2) pnerr = VMM_GET_VAR(rvw2) pnerr = VMM_GET_VAR(rcn) pnerr = VMM_GET_VAR(rth) pnerr = VMM_GET_VAR(fit1) pnerr = VMM_GET_VAR(zz1) * if (.not. Schm_hydro_L) then pnerr = VMM_GET_VAR(rw) pnerr = VMM_GET_VAR(rvv) endif * if (Tr3d_ntr.gt.0) then key1_ = VMM_KEY (trt1) key0_ = VMM_KEY (trt0) do n=1,Tr3d_ntr key1(n) = key1_ + n key0(n) = key0_ + n end do err = vmmlod(key1,Tr3d_ntr) err = vmmlod(key0,Tr3d_ntr) do n=1,Tr3d_ntr err = vmmget(key1(n),patr,tr) pnt_trt1(n) = patr err = vmmget(key0(n),patr0,tr0) pnt_trt0(n) = patr0 end do endif * * ---------------------------- * Keep positions in CAP fields * ---------------------------- !$omp parallel private(n) !$omp do do k=1,l_nk do j=j0,jn do i=i0,in n = (k-1)*nij + ((j-1)*l_ni) + i Adw_capx1(n) = F_u(n) Adw_capy1(n) = F_v(n) Adw_capz1(n) = F_w(n) end do end do end do !$omp enddo * *********************************************************************** * Perform interpolation *********************************************************************** * adw_comp_cub_L = .true. !force recomputation of position * call adw_interp2
(ruw2, ruw1, F_u, F_v, % .true., .false., LDIST_DIM, l_nk,i0,in,j0,jn) * call adw_interp2
(rvw2, rvw1, F_u, F_v, % .true., .false., LDIST_DIM, l_nk,i0,in,j0,jn) * call adw_interp2
(rcn, rcn, F_u, F_v, % .false., .false., LDIST_DIM, l_nk,i0,in,j0,jn) * call adw_interp2
(rth, rth, F_u, F_v, % .false., .false., LDIST_DIM, l_nk,i0,in,j0,jn) * call adw_interp2
(zz1, fit1, F_u, F_v, % .false., .false., LDIST_DIM, l_nk,i0,in,j0,jn) * if (.not. Schm_hydro_L) then call adw_interp2
(rw, rw, F_u, F_v, % .false., .false., LDIST_DIM, l_nk,i0,in,j0,jn) * call adw_interp2
(rvv, rvv, F_u, F_v, % .false., .false., LDIST_DIM, l_nk,i0,in,j0,jn) endif * if ( Orh_icn .eq. Schm_itcn. or. .not.Orh_crank_L ) then * * tr3d advection * aaa = ( Offc_a1_8 / Offc_b0_8 )/ Cstv_dt_8 * do n=1,Tr3d_ntr patr = pnt_trt1(n) patr0= pnt_trt0(n) !$omp do do k=1,G_nk do j=1,l_nj do i=1,l_ni tr0(i,j,k) = - aaa*tr(i,j,k) end do end do end do !$omp enddo call adw_interp2
( tr0, tr0, F_u, F_v, % .false. , Adw_mono_L, LDIST_DIM, l_nk,i0,in,j0,jn) !$omp do do k=1,G_nk do j=1,l_nj do i=1,l_ni tr0(i,j,k) = Cstv_tau_8*tr0(i,j,k) end do end do end do !$omp enddo end do * endif !$omp end parallel *********************************************************************** * Deallocate *********************************************************************** call hpdeallc(Adw_capx1_ ,err,1) call hpdeallc(Adw_capy1_ ,err,1) call hpdeallc(Adw_capz1_ ,err,1) call hpdeallc(Adw_n1_ ,err,1) call hpdeallc(Adw_xdd1_ ,err,1) call hpdeallc(Adw_xgg1_ ,err,1) call hpdeallc(Adw_ydd1_ ,err,1) call hpdeallc(Adw_ygg1_ ,err,1) call hpdeallc(Adw_cz1_ ,err,1) call hpdeallc(Adw_c1_ ,err,1) call hpdeallc(Adw_wrkb_ ,err,1) call hpdeallc(Adw_wrkc_ ,err,1) if (.not.G_lam) then call hpdeallc(Adw_capx2_ ,err,1) call hpdeallc(Adw_capy2_ ,err,1) call hpdeallc(Adw_capz2_ ,err,1) call hpdeallc(Adw_n2_ ,err,1) call hpdeallc(Adw_xdd2_ ,err,1) call hpdeallc(Adw_xgg2_ ,err,1) call hpdeallc(Adw_ydd2_ ,err,1) call hpdeallc(Adw_ygg2_ ,err,1) call hpdeallc(Adw_cz2_ ,err,1) call hpdeallc(Adw_wrka_ ,err,1) endif * if (V4dg_conf.ne.0.and.Lctl_step.eq.Step_total.and.Orh_icn.eq.Schm_itcn) done_once_L = .true. * 1000 format(3X,'ADVECT THE RIGHT-HAND-SIDES: (S/R ADW_MAIN_3_INT)') * return end