!*********************************************************************
subroutine run 1,42
! "Run" program for fim global model
! Alexander E. MacDonald 12/24/2005
! J. LEE 12/28/2005
!*********************************************************************
use module_constants
use module_control
,only: nts,numphr,nvl,nvlp1,nip,ntr,itsStart,nabl, &
CallPhysics,CallRadiation,CallOutqv,ipnDiag
use module_variables
,only: dpl,fu,fv,fdp,fdpl,flxp, &
nf,of,vof,us3d,vs3d,dp3d,pr3d,ph3d,ex3d,mp3d, &
tk3d,ssus,ssvs,ssdp,ssth,ue,ve,dpe,lpe, &
mke,adf,ab1,ab2,ab3,thtes, &
tr,tre,trpl,trdp,aftr,ftrl,ftr,ws3d,ssqv,ssqw,&
sdps,dpold,thold,rh3d,qs3d,rn2d,rc2d,pw2d, &
ts2d,us2d,hf2d,qf2d,sw2d,lw2d,st3d,sm3d,vor, &
physdiag3d,qmass,qmsta,qmstr,qmste,qmstd, &
qmsqv,qmsqw
use module_abstart
,only: abstart
use module_output
,only: output
use module_physics
,only: physics
use module_load_ls
,only: load_ls
use module_fctprs
,only: fctprs
use module_fctscl
,only: fctscl
use module_force
,only: force
use module_timedif
,only: timedif
use module_getscl
,only: getscl
use module_diag
,only: diag
use module_hybgen
,only: hybgen
use module_dp2pi
,only: dp2pi
use module_pi2dp
,only: pi2dp
use module_diagnoise
,only: diagnoise
use module_diagnoses
,only: diagnoses
use module_outqv
,only: outqv
use module_outqv_mn
,only: outqv_mn
use module_out4d_mn
,only: out4d_mn
use module_profout
,only: profout
implicit none
!SMS$insert INCLUDE "mpif.h"
!SMS$insert real*8 t0,t1
! Declare local variables:
integer :: its ! index time step
integer :: ipn
integer :: ivl
integer :: mype
real*8 :: sum
!SMS$insert call nnt_me(mype)
!SMS$insert t0 = mpi_wtime()
do its=itsStart,nts ! its=index time step, nts = num time steps
! print "('its=',i0)",its
!SMS$insert if(its==400) t0 = mpi_wtime()
!.............................................................
! Start the Adams Bash 3rd order time diff:
call abstart
( &
its, & ! index time step - in
nf,of,vof, & ! Adams Bash time dif indices-out
ab1,ab2,ab3 ) ! constants for Adams Bash - out
!
!.........................................................
! load_ls loads the "local stereographic" grid
call load_ls
( &
us3d,vs3d, & ! west, south wind on s level
dp3d,pr3d,mp3d, & ! delta p, pres, mont pot, kin energy
ue,ve,dpe,lpe, & ! u v dp and layer-pressure edge variables - out
mke,tr,tre ) ! qv and mont pot+kin ener, edge - out
!........................................................
! Calculate the "mass flux corrected transport":
call fctprs
( &
nf,of,vof, & ! Adams Bash time dif indices-in
ab1,ab2,ab3, & ! constants for Adams Bashforth - in
ue,ve,dpe, & ! edge west, so wind and vert dp-in
dp3d,dpl, & ! dp, dp low ord - inout
fdp,fdpl, & ! forcing for dp, low ord for dp-out
flxp,adf,ws3d ) ! anti diffusive flux - out, omega
!........................................................
! Calculate the "water vapor flux transport":
call fctscl
( &
nf,of,vof, & ! Adams Bash time dif indices-in
ab1,ab2,ab3, & ! constants for Adams Bashforth - in
ue,ve,tre, & ! edge west, so wind and vert dp - in
tr,trdp,trpl, & ! qv, qv*dp, qv*dp low ord - inout
ftr,ftrl, & ! forcing for qv, low ord for qv - out
flxp,aftr,its ) ! flxp - in, anti diffusive flux - out
!.........................................................
! Calculates forcing (tendency) terms:
!sms$compare_var(aftr, "run.F90 - aftr4 ")
call force
( &
nf, & ! Adams Bash time level index
us3d,vs3d,dp3d,ex3d, & ! west, south wind, delta p, exner funciton on layers
ue,ve,mke,lpe, & ! u v KE and layer-pressure edge variables - in
vor,ws3d , & ! relative vorticity, and omega=dp/dt
adf,aftr ,tre, & ! anti dif fluxes and edge-tracer variables
fu,fv,fdp ,ftr ) ! forcing for u, v, dp and tracers
!...........................................................
! Time differencing
!!!sms$compare_var(fu, "run.F90 - fu5 ")
!sms$compare_var(us3d, "run.F90 - us3d5 ")
!sms$compare_var(trdp, "run.F90 - trdp5 ")
!sms$compare_var(ftr, "run.F90 - ftr5 ")
call timedif
( &
nf,of,vof, & ! new field, old field, very o field
ab1,ab2,ab3, & ! adams bash constants - in
fu,fv,fdp, & ! forcing(tendendcy)for u,v,dp -in
dpl, & ! dp low order soln - in
us3d,vs3d, & ! west wind, south wind - inout
dp3d, & ! delta p, del p low ord - inout
trpl,ftr ,trdp,its ) ! delta qw, del qw low ord - inout
!sms$compare_var(us3d, "run.F90 - us3d6 ")
!sms$compare_var(trdp, "run.F90 - trdp6 ")
! Calculate tracer values by dividing dp
call getscl
(dp3d,trdp, tr,its) !dp3d,tp3d,tr*delp,theta*delp,theta
call dp2pi
( &
dp3d, & ! delta p between layers
pr3d, & ! pres at interface
ex3d ) ! exner at interface
!JLEE
call hybgen
(its, &
thetac, & ! target theta for hybgen- in
us3d,vs3d,tr, & ! zonal, meridional winds, specific humidity
ex3d ) ! interface exner - in/out
call pi2dp
( &
ex3d,pr3d,dp3d,tr,trdp) ! pot.temp, delta p, tr*dp3d
!JLEE
!
!...........................................................
! Do column calculations:
!...........................................................
! Condensational heating parameterizations
!
call physics
&
(its,us3d,vs3d,ws3d,dp3d, & ! west wind, south wind, vert wind, del p-in
pr3d, mp3d, & ! pres,theta,mont pot,temp(k)-in
tr ,ex3d,qs3d,rh3d,ph3d, & ! exner funciton-in
trdp, st3d,sm3d, & ! dp3d from SWE advec
ts2d,us2d,hf2d,qf2d, & ! surface temperature, friction velocity, sensible/latent heat fluxes
ssus,ssvs, & ! source/sink for u and v
ssdp,ssth,ssqv,ssqw, & ! source/sink for dp,th,qv
rn2d,rc2d,pw2d,sw2d,lw2d, &
physdiag3d, &
CallPhysics,CallRadiation,&
qmstr,qmstd ) ! Timestep interval to call physics,radiation
!...........................................................
! call hybgen
!
call hybgen
(its, &
thetac, & ! target theta for hybgen- in
us3d,vs3d,tr, & ! zonal, meridional winds, specific humidity
ex3d ) ! interface exner - in/out
!...........................................................
! Compute thickness (dp3d) from the hybgen adjusted pressure variables.
!
call pi2dp
( &
ex3d,pr3d,dp3d,tr,trdp) ! pot.temp, delta p, tr*dp3d
!
!...........................................................
! Diagnostic calculations such as hydrostatic calc
call diag
(its, &
ph3d, & ! phi (=g*z) geopotential
us3d,vs3d, & ! west wind and south wind on s - in
ex3d,mp3d,dp3d, & ! pres,exner,mont pot,kin energy - out
tr,trdp ) ! specific humidity
!
!...........................................................
! Evaluate noise parameter
!
#ifdef GLOBAL_DIAG
call diagnoise
( &
its, & ! model time step
dp3d ) ! layer thickness
call diagnoses
( &
its, & ! model time step
dp3d,tr,rn2d, &
pr3d,ex3d,qf2d,qmass, &
qmsqv,qmsqw,qmstr,qmste ) ! layer thickness
#endif
if(CallOutqv) then
! call outqv ( &
! tr(1,1,2) ,deg_lat,deg_lon,nip,nvl )
if (its.eq.1 .or. mod(its,6*numphr).eq.0) then
write (6,*)' Pressure'
call outqv_mn
(pr3d,nvlp1,nip,its,1.e-2)
write (6,*)' THETA'
call out4d_mn
(tr,nvl,nip,ntr,its,1.,1)
endif
endif
! if( mod(its,15).eq.0 ) call outqv_mn (tr(1,1,2),nvl,nip,its)
!sms$compare_var(mp3d, "run.F90 - mp3d10 ")
!
!............................................................
! Generate the output field of main variables:
!SMS$insert if(its==nts) t1 = mpi_wtime()
call output
( &
its,nts, & ! index time step, final timestep
us3d,vs3d,dp3d, & ! west wind, south wind, delta pres
pr3d, mp3d, & ! pressure, theta, mont pot, temp (k)
tr,rh3d,vor,ws3d, & ! specific and relative humidity
ph3d,tk3d,rn2d,rc2d,pw2d, &
ts2d,us2d,hf2d,qf2d, &
sw2d,lw2d,st3d,sm3d )
call profout
( &
its,ipnDiag, & ! index time step
us3d,vs3d,dp3d, & ! west wind, south wind, delta pres
pr3d,tr(:,:,1),mp3d,tk3d, & ! pressure, theta, mont pot, temp (k)
tr(:,:,2),rh3d, & ! specific and relative humidity
ph3d,rn2d,pw2d, & ! kin energy,vorticity,divergence, accumulated precipitation/rainfall
ts2d,us2d,hf2d,qf2d ) ! skin temperature, ustar, sensible heat hlux, water vapor flux
enddo ! time step loop
!SMS$insert print*,'time =',t1-t0
!#ifdef GLOBAL_DIAG
open (21,file='mst.dat',form="unformatted")
write(21) (its-1),qmass,qmsqv,qmsqw,qmstr ,qmstd, qmste
close(21)
!#endif
return
end subroutine run