!-------------------------------------- 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 v4d_setscalp - Initialization of the inner product * #include "model_macros_f.h"*
subroutine v4d_setscalp 2,1 * implicit none * *author * M.Tanguay * *revision * v2_10 - Tanguay M. - initial MPI version * v2_31 - Tanguay M. - adapt for vertical hybrid coordinate * - preserve model fields content * - adapt for tracers in tr3d * v3_02 - Laroche S. - modified weight for hu (V4dj_trwt(n)) * v3_03 - Tanguay M. - Adjoint NoHyd configuration * *object * see id section * *arguments * none * *implicits #include "glb_ld.cdk"
#include "geomn.cdk"
#include "geomg.cdk"
#include "lun.cdk"
#include "vt1.cdk"
#include "dcst.cdk"
#include "cstv.cdk"
#include "tr3d.cdk"
#include "v4dj.cdk"
#include "v4dc.cdk"
#include "schm.cdk"
* *modules integer vmmlod,vmmget,vmmuld external vmmlod,vmmget,vmmuld integer pnerr, pnlkey1(5), key1(Tr3d_ntr), $ key1_, err, i, j, k, m, n real tr pointer (patr, tr(LDIST_SHAPE,*)) * real pr1 * real*8 ZERO_8 parameter (ZERO_8=0.0) * real wu (LDIST_SHAPE,l_nk),wv (LDIST_SHAPE,l_nk) real wtp (LDIST_SHAPE,l_nk),ws (LDIST_SHAPE) real wfip(LDIST_SHAPE,l_nk),wtr(LDIST_SHAPE,l_nk,Tr3d_ntr) * * ----------------------------------------------------------------- * if(Lun_out.gt.0) write(Lun_out,1000) * * Allocation memory for cost function coefficients * ------------------------------------------------ call hpalloc(V4dj_thickx_ , l_ni, pnerr, 1) call hpalloc(V4dj_thickxu_ , l_ni, pnerr, 1) call hpalloc(V4dj_thicky_ , l_nj, pnerr, 1) call hpalloc(V4dj_thickyv_ , l_njv, pnerr, 1) call hpalloc(V4dj_thickz_ , l_nk, pnerr, 1) * * ------------------------------------------- * Initialization of weights for each variable * ------------------------------------------- V4dj_uvwt = 1. V4dj_tpwt = Dcst_cpd_8/Cstv_tstr_8 if(.not. Schm_hydro_L) V4dj_fipwt = 1. * * ---------------------------------------------------- * NOTE: Weight for surf.pres. is given with respect to * the 3D area. * ---------------------------------------------------- V4dj_spwt = Geomg_z_8(l_nk)*Dcst_rgasv_8*Cstv_tstr_8 * if (Tr3d_ntr.gt.0) then do n=1,Tr3d_ntr V4dj_trwt(n) = ( Dcst_chlc_8**2 )/ ( Dcst_cpd_8 * 270. ) enddo endif * * ---------------------------------- * Integral with respect to longitude * ---------------------------------- * * Staggering scalar grid * ---------------------- do i=1,l_ni V4dj_thickx(i) = Geomg_hxu_8(i-1) end do * * Staggering u grid * ----------------- do i=1,l_ni V4dj_thickxu(i)= Geomg_hx_8(i) end do * * --------------------------------- * Integral with respect to latitude * --------------------------------- * * Staggering scalar grid * ---------------------- do j=1,l_nj V4dj_thicky(j) = Geomg_hsyv_8(j-1) end do * * Staggering v grid * ----------------- do j=1,l_njv V4dj_thickyv(j)= Geomg_hsy_8(j) end do * * --------------------------------- * Integral with respect to vertical * --------------------------------- V4dj_thickz(1) = .5*( Geomg_z_8(2 ) - Geomg_z_8(1 ) ) do k=2,l_nk-1 V4dj_thickz(k) = .5*( Geomg_z_8(k+1 ) - Geomg_z_8(k-1 ) ) end do V4dj_thickz(l_nk)= .5*( Geomg_z_8(l_nk) - Geomg_z_8(l_nk-1) ) * * ------------------------------------ * Inverse of total area without ray**2 * ------------------------------------ V4dj_invarea = 4.* Dcst_pi_8 **2 * (Geomg_z_8(l_nk) - Cstv_pitop_ %8) V4dj_invarea = 1./V4dj_invarea * * Get fields in memory * -------------------- pnlkey1(1) = vmmk_ut1 pnlkey1(2) = vmmk_vt1 pnlkey1(3) = vmmk_tpt1 pnlkey1(4) = vmmk_st1 if (.not. Schm_hydro_L) then pnlkey1(5) = vmmk_fipt1 pnerr = vmmlod(pnlkey1,5) else pnerr = vmmlod(pnlkey1,4) endif * pnerr = vmmget(vmmk_ut1 , ut1_ , ut1 ) pnerr = vmmget(vmmk_vt1 , vt1_ , vt1 ) pnerr = vmmget(vmmk_tpt1 , tpt1_ , tpt1 ) pnerr = vmmget(vmmk_st1 , st1_ , st1 ) if (.not. Schm_hydro_L) then pnerr = vmmget(vmmk_fipt1 , fipt1_, fipt1) endif * * ----------------------- * Keep original variables * ----------------------- do k=1,l_nk do j=l_miny,l_maxy do i=l_minx,l_maxx wu (i,j,k)= ut1 (i,j,k) wv (i,j,k)= vt1 (i,j,k) wtp(i,j,k)= tpt1(i,j,k) end do end do end do do j=l_miny,l_maxy do i=l_minx,l_maxx ws(i,j) = st1(i,j) end do end do * if (.not. Schm_hydro_L) then do k=1,l_nk do j=l_miny,l_maxy do i=l_minx,l_maxx wfip(i,j,k)= fipt1(i,j,k) end do end do end do endif * * ---------------------- * Reinitialize variables * ---------------------- do k=1,l_nk do j=l_miny,l_maxy do i=l_minx,l_maxx ut1 (i,j,k) = ZERO_8 vt1 (i,j,k) = ZERO_8 tpt1(i,j,k) = ZERO_8 end do end do end do do j=l_miny,l_maxy do i=l_minx,l_maxx st1(i,j) = ZERO_8 end do end do * if (.not. Schm_hydro_L) then do k=1,l_nk do j=l_miny,l_maxy do i=l_minx,l_maxx fipt1(i,j,k) = ZERO_8 end do end do end do endif * do k = 1,l_nk do j = 1,l_nj do i = 1,l_niu ut1(i,j,k) = V4dj_thickxu(i) * V4dj_thicky(j) * V4dj_thickz(k) * % V4dj_uvwt * V4dj_invarea end do end do end do * do k = 1,l_nk do j = 1,l_njv do i = 1,l_ni vt1(i,j,k) = V4dj_thickx(i) * V4dj_thickyv(j) * V4dj_thickz(k) * % V4dj_uvwt * V4dj_invarea end do end do end do * do k = 1,l_nk do j = 1,l_nj do i = 1,l_ni tpt1(i,j,k) = V4dj_thickx(i) * V4dj_thicky(j) * V4dj_thickz(k) * % V4dj_tpwt * V4dj_invarea end do end do end do * do j = 1,l_nj do i = 1,l_ni st1(i,j) = V4dj_thickx(i) * V4dj_thicky(j) * % V4dj_spwt * V4dj_invarea end do end do * if (.not. Schm_hydro_L) then do k = 1,l_nk do j = 1,l_nj do i = 1,l_ni fipt1(i,j,k) = V4dj_thickx(i) * V4dj_thicky(j) * V4dj_thickz(k) * % V4dj_tpwt * V4dj_invarea end do end do end do endif * pnerr = vmmuld(-1,0) * key1_ = vmmk_trt1 do n=1,Tr3d_ntr key1(n) = key1_ + n end do * if (Tr3d_ntr.gt.0) then err = vmmlod(key1,Tr3d_ntr) do n=1,Tr3d_ntr err = vmmget(key1(n),patr,tr) * * ----------------------- * Keep original variables * ----------------------- do k=1,l_nk do j=l_miny,l_maxy do i=l_minx,l_maxx wtr(i,j,k,n) = tr(i,j,k) end do end do end do * * ---------------------- * Reinitialize variables * ---------------------- do k=1,l_nk do j=l_miny,l_maxy do i=l_minx,l_maxx tr(i,j,k) = ZERO_8 end do end do end do do k = 1,l_nk do j = 1,l_nj do i = 1,l_ni tr(i,j,k) = V4dj_thickx(i) * V4dj_thicky(j) * V4dj_thickz %(k) * % V4dj_trwt(n) * V4dj_invarea end do end do end do * end do err = vmmuld(key1,Tr3d_ntr) endif * call v4d_cainin
(V4dc_ncv,V4dc_scalp) * * * Get fields in memory * -------------------- pnlkey1(1) = vmmk_ut1 pnlkey1(2) = vmmk_vt1 pnlkey1(3) = vmmk_tpt1 pnlkey1(4) = vmmk_st1 if (.not. Schm_hydro_L) then pnlkey1(5) = vmmk_fipt1 pnerr = vmmlod(pnlkey1,5) else pnerr = vmmlod(pnlkey1,4) endif * pnerr = vmmget(vmmk_ut1 , ut1_ , ut1 ) pnerr = vmmget(vmmk_vt1 , vt1_ , vt1 ) pnerr = vmmget(vmmk_tpt1 , tpt1_ , tpt1 ) pnerr = vmmget(vmmk_st1 , st1_ , st1 ) if (.not. Schm_hydro_L) then pnerr = vmmget(vmmk_fipt1 , fipt1_, fipt1) endif * * ------------------------ * Reset original variables * ------------------------ do k=1,l_nk do j=l_miny,l_maxy do i=l_minx,l_maxx ut1 (i,j,k) = wu (i,j,k) vt1 (i,j,k) = wv (i,j,k) tpt1(i,j,k) = wtp(i,j,k) end do end do end do do j=l_miny,l_maxy do i=l_minx,l_maxx st1(i,j) = ws(i,j) end do end do * if (.not. Schm_hydro_L) then * do k=1,l_nk do j=l_miny,l_maxy do i=l_minx,l_maxx fipt1 (i,j,k) = wfip (i,j,k) end do end do end do * endif * pnerr = vmmuld(-1,0) * key1_ = vmmk_trt1 do n=1,Tr3d_ntr key1(n) = key1_ + n end do * if (Tr3d_ntr.gt.0) then err = vmmlod(key1,Tr3d_ntr) do n=1,Tr3d_ntr err = vmmget(key1(n),patr,tr) do k=1,l_nk do j=l_miny,l_maxy do i=l_minx,l_maxx tr(i,j,k) = wtr(i,j,k,n) end do end do end do end do err = vmmuld(key1,Tr3d_ntr) endif * 1000 format( +//,'INITIALIZATION INNER PRODUCT (S/R V4D_SETSCALP)', + /,'===============================================', +//) * * --------------------------------------------------------------- * return end