!-------------------------------------- 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 hzd_main - applies horizontal diffusion on a given set of fields
*
#include "model_macros_f.h"
*
subroutine hzd_main 3,10
*
implicit none
*
*author
* Joseph-Pierre Toviessi ( after version v1_03 of dif )
*
*revision
* v2_00 - Desgagne M. - initial MPI version
* v2_10 - Qaddouri&Desgagne - higher order diffusion operator
* v2_21 - Desgagne M. - new call to horwavg
* v2_30 - Edouard S. - adapt for vertical hybrid coordinate
* v3_00 - Desgagne & Lee - Lam configuration
* v3_01 - Toviessi J. P. - add call hzd_ho_parite
* v3_02 - Desgagne M. - correction for non-hydrostatic version
* v3_10 - Corbeil & Desgagne & Lee - AIXport+Opti+OpenMP
* v3_20 - Tanguay M. - Introduce Hzd_hzdmain_n_L
* v3_21 - Desgagne M. - added explicit horiz diff.
* v3_30 - Tanguay M. - activate Hzd_type_S='EXPLICIT'
*
*object
*
*arguments
* none
*
*implicits
#include "glb_ld.cdk"
#include "cstv.cdk"
#include "schm.cdk"
#include "lun.cdk"
#include "geomg.cdk"
#include "fft.cdk"
#include "hzd.cdk"
#include "eigv.cdk"
#include "vt0.cdk"
#include "vt1.cdk"
#include "vtx.cdk"
#include "p_geof.cdk"
#include "grd.cdk"
#include "lctl.cdk"
*
*modules
integer vmmlod,vmmuld,vmmget
external vmmlod,vmmuld,vmmget
*
integer err, nlod, key1(20), i, j, k
real*8 ONE_8, bbb_8, eta_8
parameter( ONE_8 = 1. )
**
* _________________________________________________________________
*
if (.not.G_lam.and.Hzd_type_S.eq.'EXPLICIT')
% call gem_stop
('ABORT .not.G_lam.and.Hzd_type=EXPLICIT',-1)
*
if ((.not.Hzd_fact_L).and.(.not.Hzd_ho_L).or..not.Hzd_hzdmain_0_L) then
if (Lun_debug_L) write(Lun_out,1001)
return
endif
*
key1(1) = VMM_KEY(ut1)
key1(2) = VMM_KEY(vt1)
key1(3) = VMM_KEY(psdt1)
key1(4) = VMM_KEY(tpt1)
key1(5) = VMM_KEY(fipt1)
key1(6) = VMM_KEY(tdt1)
key1(7) = VMM_KEY(tt1)
key1(8) = VMM_KEY(fit1)
key1(9) = VMM_KEY(topo)
key1(10) = VMM_KEY(tplt1)
key1(11) = VMM_KEY(qt1)
key1(12) = VMM_KEY(pipt0)
key1(13) = VMM_KEY(pipt1)
key1(14) = VMM_KEY(st1)
nlod = 14
if (.not. Schm_hydro_L) then
key1(nlod+1) = VMM_KEY(wt1)
key1(nlod+2) = VMM_KEY(qpt1)
key1(nlod+3) = VMM_KEY(mut1)
key1(nlod+4) = VMM_KEY(multx)
nlod = nlod+4
endif
* - - - - - - - - - - - - - - -
err = vmmlod(key1,nlod)
* - - - - - - - - - - - - - - -
err = VMM_GET_VAR(ut1)
err = VMM_GET_VAR(vt1)
err = VMM_GET_VAR(psdt1)
err = VMM_GET_VAR(tpt1)
err = VMM_GET_VAR(fipt1)
err = VMM_GET_VAR(tdt1)
err = VMM_GET_VAR(tt1)
err = VMM_GET_VAR(fit1)
err = VMM_GET_VAR(topo)
err = VMM_GET_VAR(tplt1)
err = VMM_GET_VAR(qt1)
err = VMM_GET_VAR(pipt0)
err = VMM_GET_VAR(pipt1)
err = VMM_GET_VAR(st1)
wt1_ = 0
qpt1_ = 0
mut1_ = 0
multx_= 0
if (.not. Schm_hydro_L) then
err = VMM_GET_VAR(wt1)
err = VMM_GET_VAR(qpt1)
err = VMM_GET_VAR(mut1)
err = VMM_GET_VAR(multx)
endif
*************************************
* 1. Implicit horizontal diffusion *
*************************************
if (Hzd_rwnd_L) call iw2rwnd
( ut1, vt1, LDIST_DIM, G_nk, 1 )
if (Hzd_fact_L) then
call hzd_fact
( ut1 , vt1 , psdt1, tpt1, fipt1, tdt1,
% tt1 , fit1, tplt1, qt1 , pipt1, wt1 ,
% qpt1, mut1, multx, st1 , topo ,
% LDIST_DIM, G_nk)
endif
*
if (Hzd_ho_L) then
if (Hzd_type_S.eq.'EXPLICIT') then
call hzd_exhrdif
( ut1 , vt1 , psdt1, tdt1, LDIST_DIM, G_nk)
else
if (Fft_fast_L) then
* use FFT in diffusion-solver
if (G_lam) then
call hzd_hof_lam
( ut1 , vt1 , psdt1, tpt1, fipt1, tdt1,
% tt1 , fit1, tplt1, pipt1, wt1 ,
% qpt1, mut1, multx, topo ,
% LDIST_DIM, G_nk)
else
call hzd_hof
( ut1 , vt1 , psdt1, tpt1, fipt1, tdt1,
% tt1 , fit1, tplt1, pipt1, wt1 ,
% qpt1, mut1, multx, topo ,
% LDIST_DIM, G_nk)
endif
else
* use MXMA in diffusion-solver
if ( .not. Eigv_parity_L) then
call hzd_ho
( ut1 , vt1 , psdt1, tpt1, fipt1, tdt1,
% tt1 , fit1, tplt1, pipt1, wt1 ,
% qpt1, mut1, multx, topo ,
% LDIST_DIM, G_nk)
else
call hzd_ho_parite
( ut1 , vt1 , psdt1, tpt1, fipt1, tdt1,
% tt1 , fit1, tplt1, pipt1, wt1 ,
% qpt1, mut1, multx, topo ,
% LDIST_DIM, G_nk)
endif
endif
endif
endif
if ((Hzd_fact_L.and.Cstv_phidf_8.gt.0.).or.(Hzd_ho_L)) then
*
* Adjust surface pressure and hence pi' (FOR CLIMATE APPLICATIONS)
* ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
!$omp parallel
if (Schm_psadj_L) then
*
!$omp single
call horwavg
( bbb_8, pipt0(l_minx,l_miny,G_nk),
$ pipt1(l_minx,l_miny,G_nk),LDIST_DIM )
!$omp end single
*
* Redistribute the average mass loss at the surface, ...
*
!$omp do
do j= 1+pil_s, l_nj-pil_n
do i= 1+pil_w, l_ni-pil_e
pipt1(i,j,G_nk) = pipt1(i,j,G_nk) + bbb_8
end do
! end do
!!$omp enddo
*
* ... correct s immediately and ...
*
!!$omp do
! do j= 1+pil_s, l_nj-pil_n
do i= 1+pil_w, l_ni-pil_e
st1(i,j)= log(ONE_8+dble(pipt1(i,j,G_nk)/geomg_pib(G_nk)))
end do
end do
!$omp enddo
*
* ... refine pi' everywhere else accordingly.
*
!$omp do
do k=1,G_nk-1
do j= 1+pil_s, l_nj-pil_n
do i= 1+pil_w, l_ni-pil_e
pipt1(i,j,k) = geomg_pib(k) *( exp(st1(i,j)) - 1 )
end do
end do
end do
!$omp enddo
*
else
*
*******************************************
* $. Indirect diffusion of s through pi' *
*******************************************
*
!$omp do
do j= 1+pil_s, l_nj-pil_n
do i= 1+pil_w, l_ni-pil_e
st1(i,j)= log(ONE_8+dble(pipt1(i,j,G_nk)/geomg_pib(G_nk)))
end do
end do
!$omp enddo
*
endif
*******************************************
* $. Indirect diffusion of q through pi' *
*******************************************
*
!$omp do
do k=1,G_nk
do j= 1+pil_s, l_nj-pil_n
do i= 1+pil_w, l_ni-pil_e
qt1(i,j,k) = log( Geomg_z_8(k) + dble(pipt1(i,j,k)) )
end do
end do
end do
!$omp enddo
*
if ( .not. Schm_hydro_L ) then
!$omp do
do k=1,G_nk
do j= 1+pil_s, l_nj-pil_n
do i= 1+pil_w, l_ni-pil_e
qt1(i,j,k) = qt1(i,j,k) + qpt1(i,j,k)
end do
end do
end do
!$omp enddo
endif
!$omp end parallel
*
endif
*
if (Hzd_rwnd_L) call iw2rwnd
( ut1, vt1, LDIST_DIM, G_nk, 2 )
*
1001 format(/,3X,'NO HORIZONTAL DIFFUSION REQUIRED')
*
return
end