Interoperable Physics Driver for NGGPS
gloopb.f
Go to the documentation of this file.
1 
5 !
13 
14  subroutine gloopb(grid_fld, g3d_fld, sfc_fld,
15  & flx_fld, aoi_fld, nst_fld,
16  & lats_nodes_r, global_lats_r, lonsperlar,
17  & tstep, phour, sfalb, xlon,
18  & swh, swhc, hlw, hlwc,
19 ! & nbdsw, nbdlw, htrsw htrlwb,
20  & hprime, slag, sdec, cdec,
21  & ozplin, jindx1, jindx2, ddy,
22  & h2oplin, jindx1_h, jindx2_h, ddy_h,
23  & phy_f3d, phy_f2d, phy_fctd, nctp,
24  & xlat, nblck, kdt, restart_step,
25  & mdl_parm)
26 !!
27 !! Code Revision:
28 !! Sep 2009 Shrinivas Moorthi added nst_fld
29 !! Oct 11 2009 Sarah Lu, grid_gr replaced by grid_fld
30 !! Oct 16 2009 Sarah Lu, grid_fld%tracers used
31 !! Nov 18 2009 Sarah Lu, rain/rainc added to gbphys call arg
32 !! Dec 14 2009 Sarah Lu, add g3d_fld to calling argument,
33 !! update dqdt after gbphys returns dqdt_v
34 !! July 2010 Shrinivas Moorthi - Updated for new physics
35 !! Aug 2010 Shrinivas Moorthi - Recoded 3d diagnostic arrays so that
36 ! trap will not occur on call to gbphys
37 !! Oct 18 2010 Shrinivas Moorthi - Added fscav
38 !! Dec 23 2010 Sarah Lu, add lgocart to gbphys call arg
39 !! Apr 06 2012 Henry Juang, add idea related changes
40 !! Oct 18 2012 Shrinivas Moorthi - Added random number realted chages
41 !! Jul 26 2012 Jun Wang pass mpi info to idea_phys
42 !! Nov 30 2012 Jun Wang update idea_phys using whole band radiation
43 !! Dec 27 2012 Jun Wang move co2 init to gloopb
44 !! Mar 29 2013 Shrinivas Moorthi - Added dtphys option
45 !! Apr 08 2013 Jun Wang add idea init variables to restart file
46 !! Oct 21 2013 Henry Juang compute prsi from model top
47 !! Oct 31 2013 Xingren Wu add flx_fld%dusfci/dvsfci
48 !! Mar 26 2014 Xingren Wu add aoi_fld
49 !! May 05 2014 Jun Wang add cgwf,prslrd0 in gbphys argument list
50 !! jun 2014 y-t hou, revised sw sfc spectral component fluxes
51 !! and ocean albedoes (no ice contamination) for coupled mdl
52 !! Jun 2014 Xingren Wu update net SW fluxes over the ocean
53 !! (no ice contamination)
54 !! Jul 2014 Xingren Wu Add Sea/Land/Ice Mask - slimsk_cpl
55 !! Jul -- 2014 Shrinivas Moorthi updating to semi-lagrangian and latest gfs physics
56 !! and some cleanup - bugfix for random number
57 !! Aug 01 2014 Shrinivas Moorthi - add tracer fixer option
58 !! Sep 16 2014 Shrinivas Moorthi - cleanup and rearrange argumets
59 !! Sep 30 2014 Sarah Lu, remove fscav (the option to compute tracer
60 !! scavenging in GFS is disable)
61 !! Dec 21 2014 Jun Wang, add g3d_fld cnv_qc to gbphys argument,
62 !! update 3D fields after gbphys
63 !! --- -- 2014 Donald Dazlich - Added CS convection
64 !! Apr 09 2014 S Moorthi - ported CS convection to NEMS/GSM
65 !! May 28 2015 Xingren Wu add t/q/u/v/p/zbot for cice coupling
66 !! Jul/Aug2015 S Moorthi - added SHOC option
67 !! Aug 21 2015 Xu Li - change nst_fcst to be nstf_name
68 !! Sep 25 2015 Xingren Wu connect slimskin/dusfcin/dvsfcin/dtsfcin/dqsfcin/ulwsfcin
69 !! for coupling
70 !! Oct 15 2015 Xingren Wu add aoi_fld%snow for cice coupling
71 !! Oct 19 2015 Weiyu Yang- add the inputted f10.7 and kp data.
72 !! Jan 13 2016 S Moorthi - add facsw to account for fhswr in seconds
73 !! Jan 2016 P. Tripp - NUOPC/GSM merge
74 !! Feb 10 2016 Hang Lei - use physics driver
75 !! Mar 03 2016 J. Han - add ncnvcld3d integer
76 !! for convective cloudiness enhancement
77 !! Mar 04 2016 J. Han - change newsas & sashal to imfdeepcnv
78 !! & imfshalcnv, respectively
79 !! Apr 15 2016 S Moorthi - initialize sumtrc, thread random number call
80 !! May -- 2016 S Moorthi - optimize and cleanup stochastic physics part
81 !! Jun 10 2016 S Moorthi - add stratospheric h2o parameterization
82 !! Jun 15 2016 X. Li - Add to assign 5 more nst_fld (nst_z_c ...)
83 !! Jun 18 2016 Hang Lei - Add ozone and stratospheric h2o into the physics driver
84 !! Jul 20 2016 Hang Lei - Add stochastic physics into the physics driver
85 !
86 
87  use resol_def
88  use layout1
89  use gg_def
90  use vert_def
91  use date_def
92  use namelist_physics_def
93  use coordinate_def
94  use module_ras , only: ras_init
95  use physcons, grav => con_g , rerth => con_rerth, rk => con_rocp
96  use ozne_def , only : latsozp,levozp,
97  & pl_coeff,pl_pres,timeoz
98  use h2o_def , only : latsh2o, levh2o
99  &, h2o_coeff, h2o_pres, timeh2o
100  use d3d_def
101  use gfs_physics_sfc_flx_mod
102  use gfs_physics_nst_var_mod
103  use gfs_physics_aoi_var_mod
104  use gfs_physics_gridgr_mod, ONLY: grid_var_data
105  use gfs_physics_g3d_mod, ONLY: g3d_var_data
106 
109  use nuopc_physics, only: state_fields_in,
113  & model_parameters, nuopc_phys_run, nuopc_sppt_phys,
114  & tbd_ddt,
115 ! & dyn_param_setphys, state_fld_setphys_in,state_fld_setphys_out,
116 ! & diagnostics_setphys, interface_fld_setphys,rad_tend_set,
117 ! & sfc_prop_setphys, cld_prop_setphys, tbd_set,
118  & phys_run_savein, phys_run_readin,
119  & phys_run_saveout, phys_run_readout
120 
121  use mpi_def, only: mpi_r_io_r, mpi_comm_all
122  use wam_f107_kp_mod, ONLY: read_wam_f107_kp_txt,
123  & f107, kp, f107_kp_size
124  use idea_composition, only: prlog,pr_idea,amgm,amgms,prsilvl
125  &, nlev_co2,k43,nlevc_h2o,k71,gg
126  use efield, only: efield_init
127  use mersenne_twister
128 
129  implicit none
130  include 'mpif.h'
131 !
132  TYPE(grid_var_data) :: grid_fld
133  TYPE(sfc_var_data) :: sfc_fld
134  TYPE(flx_var_data) :: flx_fld
135  TYPE(nst_var_data) :: nst_fld
136  TYPE(g3d_var_data) :: g3d_fld
137  TYPE(aoi_var_data) :: aoi_fld
138 
139 !
140  integer nblck, kdt, nbdsw, nbdlw, nctp
141  logical restart_step
142 !!
143  real(kind=kind_phys) tstep, phour,slag, sdec, cdec
144  real(kind=kind_phys) plyr(levs)
145  real(kind=kind_phys), dimension(ngptc) :: dpshc, pgr
146  real(kind=kind_phys), dimension(ngptc,levs) :: prsl, prslk
147  &, phil, gu, gv, gt, adu, adv, adt
148  real(kind=kind_phys), dimension(ngptc,levs+1) :: prsi, prsik
149  &, phii
150 !!
151  real (kind=kind_rad), dimension(ngptc,levs,ntrac) :: gr, adr
152 !!
153  real (kind=kind_rad), dimension(lonr,lats_node_r) :: xlon, xlat
154  &, sfalb
155 !!
156  real (kind=kind_rad), dimension(ngptc,levs,nblck,lats_node_r) ::
157  & swh, swhc, hlw, hlwc
158 !!
159  real (kind=kind_rad) hprime(nmtvr,lonr,lats_node_r)
160 
161 !idea add by hmhj
162  real (kind=kind_rad) hlwd(ngptc,levs,6)
163 ! &, htrswb(ngptc,levs,nbdsw,nblck,lats_node_r)
164 ! &, htrlwb(ngptc,levs,nbdlw,nblck,lats_node_r)
165 !!
166  real (kind=kind_phys) phy_f3d(ngptc,levs,ntot3d,nblck,lats_node_r)
167  &, phy_f2d(lonr,lats_node_r,ntot2d)
168  &, phy_fctd(lonr,lats_node_r,nctp) ! for CS convection
169 !!
170  real (kind=kind_phys) dtp,dtf,solhr,clstp
171 !!
172  integer lats_nodes_r(nodes)
173  integer, dimension(latr) :: global_lats_r, lonsperlar
174 !
175  integer i,j,k,kk,n,l,lan,lat,ii,lonrbm,jj
176  &, lon_dim,lons_lat,nsphys
177 !
178 !timers______________________________________________________---
179 
180 ! real*8 rtc ,timer1,timer2
181 ! real(kind=kind_evod) global_times_b(latr,nodes)
182 
183 !timers______________________________________________________---
184 !flipv is declared in namelist_phys_def
185 ! logical, parameter :: flipv = .true.
186  real(kind=kind_phys), parameter :: pt01=0.01, pt00001=1.0e-5
187  &, thousnd=1000.0
188  &, cons_0=0.0, cons_24=24.0
189  &, cons_99=99.0, cons_1p0d9=1.0e9
190  &, qmin=1.0e-10
191 !
192 ! for nrl/nasa ozone production and distruction rates:(input through fixio)
193 ! ---------------------------------------------------
194  integer, dimension(lats_node_r) :: jindx1, jindx2 !for ozone interpolation
195  &, jindx1_h, jindx2_h
196  real(kind=kind_phys) ozplin(latsozp,levozp,pl_coeff,timeoz)
197  &, ddy(lats_node_r) !for ozone interpolation
198  &, ddy_h(lats_node_r) ! for h2o interpolation
199  &, ozplout(levozp,lats_node_r,pl_coeff)
200  &, h2oplin(latsh2o,levh2o,h2o_coeff,timeh2o)
201  &, h2oplout(levh2o,lats_node_r,h2o_coeff)
202  &, tmp(lonr,lats_node_r)
203 !!
204  real(kind=kind_phys), allocatable :: acv(:,:),acvb(:,:),acvt(:,:)
205  &, fscav(:), fswtr(:)
206  save acv,acvb,acvt,fscav,fswtr
207 !!
208 ! integer, parameter :: maxran=5000
209 ! integer, parameter :: maxran=3000
210 ! integer, parameter :: maxran=6000, maxsub=6, maxrs=maxran/maxsub
211 ! integer, parameter :: maxran=3000, maxsub=6, maxrs=maxran/maxsub
212  integer, parameter :: maxran=3000, maxsub=10, maxrs=maxran/maxsub
213  type(random_stat) :: stat(lats_node_r)
214 
215  real (kind=kind_phys), allocatable, save :: rannum_tank(:,:,:)
216  real (kind=kind_phys), allocatable :: rannum(:), wrkn(:)
217  integer, allocatable :: indxr(:)
218 ! integer, allocatable, save :: indxr(:)
219  integer iseed, nrc, seed0, kss, ksr, iseedl, latseed
220  &, nf0,nf1,ind,nt,indod,indev
221  real(kind=kind_evod) fd2, wrk(1), facsw
222 
223  logical first
224  data first/.true./
225  save first, seed0, facsw
226 !
227  integer nlons_v(ngptc)
228 
229  real(kind=kind_phys), dimension(ngptc) :: sinlat_v,coslat_v,rqtk
230 ! real(kind=kind_phys), dimension(njeff) :: sinlat_v,coslat_v,rqtk
231  real(kind=kind_phys), dimension(ngptc,lsoil) :: smc_v,stc_v,slc_v
232  real(kind=kind_phys), dimension(ngptc,levs) :: vvel,dtdt,dqdt_v,
233  & cnvqc_v
234 
235  real(kind=kind_phys) hprime_v(ngptc,nmtvr)
236  &, phy_f2dv(ngptc,ntot2d)
237  &, rannum_v(ngptc,nrcm)
238  &, ozplout_v(ngptc,levozp,pl_coeff)
239  &, h2oplout_v(ngptc,levh2o,h2o_coeff)
240  &, phy_fctdv(ngptc,nctp) ! for CS convection
241 
242  real(kind=kind_phys) trcg(latr,ntrac,2),trcj(lats_node_r,ntrac,2)
243  &, trcp(lonr,ntrac,2)
244 
245  real(kind=kind_phys), allocatable, save :: sumtrc(:,:), adjtrc(:)
246 
247 ! variables for stochastic physics
248 !-----------------------------------
249 
250  real (kind=kind_phys), dimension(ngptc,levs) :: gu0,gv0,gt0
251  real (kind=kind_phys), dimension(ngptc,levs) :: qphys_v
252  real (kind=kind_phys), dimension(ngptc,levs,ntrac) :: gr0
253  real (kind=kind_phys) :: sppt_wts, qphys
254  real (kind=kind_phys), dimension(ngptc) :: cplrain0,cplsnow0
255  &, totprcp0,cnvprcp0
256 !!!1
257 
258  real(kind=kind_phys) dt3dt_v(ngptc,levs,6), du3dt_v(ngptc,levs,4)
259  &, dv3dt_v(ngptc,levs,4)
260  &, dq3dt_v(ngptc,levs,5+pl_coeff)
261  real(kind=kind_phys) upd_mfv(ngptc,levs), dwn_mfv(ngptc,levs)
262  &, det_mfv(ngptc,levs)
263 ! &, det_mfv(ngptc,levs), dkh_v(ngptc,levs)
264 ! &, rnp_v(ngptc,levs)
265  real(kind=kind_phys), dimension(ngptc) :: nirbmdi, nirdfdi,
266  & visbmdi, visdfdi, nirbmui, nirdfui,
267  & visbmui, visdfui
268  &, aoi_du, aoi_dv, aoi_dt, aoi_dq
269  &, aoi_dlw, aoi_dsw, aoi_dnirbm
270  &, aoi_dnirdf, aoi_dvisbm, aoi_dvisdf
271  &, aoi_rain, aoi_nlw, aoi_nsw
272  &, aoi_nnirbm, aoi_nnirdf, aoi_nvisbm
273  &, aoi_nvisdf, aoi_snow
274 
275  &, aoi_dusfci, aoi_dvsfci, aoi_dtsfci
276  &, aoi_dqsfci, aoi_dlwsfci, aoi_dswsfci
277  &, aoi_dnirbmi, aoi_dnirdfi, aoi_dvisbmi
278  &, aoi_dvisdfi, aoi_nlwsfci, aoi_nswsfci
279  &, aoi_nnirbmi, aoi_nnirdfi, aoi_nvisbmi
280  &, aoi_nvisdfi, aoi_t2mi, aoi_q2mi
281  &, aoi_u10mi, aoi_v10mi, aoi_tseai
282  &, aoi_psurfi
283  &, aoi_slimskin,aoi_ulwsfcin,aoi_dusfcin
284  &, aoi_dvsfcin, aoi_dqsfcin, aoi_dtsfcin
285 
286  &, nst_xt, nst_xs, nst_xu, nst_xv, nst_xz
287  &, nst_zm, nst_xtts, nst_xzts, nst_d_conv
288  &, nst_ifd, nst_dt_cool, nst_qrain
289 
290  &, nst_tref, nst_z_c, nst_c_0, nst_c_d
291  &, nst_w_0, nst_w_d
292 
293  real(kind=kind_phys) work1, tem
294  integer njeff,lon,iblk,item, nn, nnr, nnrcm, dbgu
295 
296 ! idea local vars
297  real(kind=kind_phys) pmod(levs),gg1(levs),philco2(levs),
298  & qtrac(levs,ntrac),prsilvl1(levs+1)
299  real(kind=kind_phys),allocatable :: prpa(:)
300  integer info,pelat1,pelatall,lanlat1
301  real(kind=8) :: gb_ini_time=0, btime, timef
302 
303 ! NUOPC physics driver types - PT
304 !! Defined in nuopc_physics.F90
305  type(model_parameters), intent(in) :: mdl_parm
306 
309  type(state_fields_in) :: state_fldin
310  type(state_fields_out) :: state_fldout
311  type(sfc_properties) :: sfc_prop
312  type(diagnostics) :: diags
313  type(interface_fields) :: intrfc_fld
314  type(cloud_properties) :: cld_prop
315  type(radiation_tendencies) :: rad_tend
316  type(tbd_ddt) :: tbddata ! To be determined where
317  ! this data should go
318  type(dynamic_parameters) :: dyn_parm
319 
320  integer :: lonbnd ! upper lon dimension in lon/lan loop
321  logical :: savecon ! Save nuopc driver in/out states
322 !
323 !---------------------------------------------------------------------------
324 !
325 !---------------------------------------------------------------------------
326 !
327  if (first) then
328  allocate (acv(lonr,lats_node_r))
329  allocate (acvb(lonr,lats_node_r))
330  allocate (acvt(lonr,lats_node_r))
331  allocate (fscav(ntrac-ncld+2), fswtr(ntrac-ncld+2))
332 
333  allocate (sumtrc(ntrac,3), adjtrc(ntrac))
334  sumtrc(:,:) = 0.0
335  adjtrc = 1.0
336  fscav = 0.0
337  fswtr = 0.0
338  facsw = 1.0
339  if (fhswr > 3.0) facsw = 1.0 / 3600.0
340 !
341  if (imfdeepcnv <= 0 .or. cal_pre) then ! random number needed for RAS and old SAS
342  if (random_clds) then ! create random number tank
343 ! -------------------------
344  seed0 = idate(1) + idate(2) + idate(3) + idate(4)
345 
346  call random_setseed(seed0)
347  call random_number(wrk)
348  seed0 = seed0 + nint(wrk(1)*thousnd)
349  if (me == 0) print *,' seed0=',seed0,' idate=',idate,
350  & ' wrk=',wrk
351 !
352  if (.not. allocated(rannum_tank))
353  & allocate (rannum_tank(lonr,maxran,lats_node_r))
354 ! if (.not. allocated(rannum)) allocate (rannum(lonr*maxrs))
355  allocate (rannum(lonr*maxrs))
356  lonrbm = lonr / maxsub
357  if (me == 0) write(0,*)' maxran=',maxran,' maxrs=',maxrs,
358  & 'maxsub=',maxsub,' lonrbm=',lonrbm,
359  & ' lats_node_r=',lats_node_r
360 !$omp parallel do private(j,iseedl, nrc, nn, i, ii, k, kk, rannum)
361  do j=1,lats_node_r
362  iseedl = global_lats_r(ipt_lats_node_r-1+j) + seed0
363  call random_setseed(iseedl, stat(j))
364  call random_number(rannum, stat(j))
365 !!$omp parallel do shared(j,lonr,lonrbm,rannum,rannum_tank)
366 !!$omp+private(nrc,nn,i,ii,k,kk)
367 !!$omp parallel do private(nrc, nn, i, ii, k, kk)
368  do nrc=1,maxrs
369  nn = (nrc-1)*lonr
370  do k=1,maxsub
371  kk = k - 1
372  do i=1,lonr
373  ii = kk*lonrbm + i
374  if (ii > lonr) ii = ii - lonr
375  rannum_tank(i,nrc+kk*maxrs,j) = rannum(ii+nn)
376  enddo
377  enddo
378  enddo
379  enddo
380  if (allocated(rannum)) deallocate (rannum)
381  endif
382  endif
383 !
384  if (me == 0) then
385 ! write(0,*)' seed0=',seed0,' idate=',idate,' wrk=',wrk
386  if (num_p3d == 3)
387  & write(0,*)' USING Ferrier-MICROPHYSICS'
388  if (num_p3d == 4)
389  & write(0,*)' USING ZHAO-MICROPHYSICS'
390  endif
391  if (fhour == 0.0) then
392 !$omp parallel do private(i,j)
393  do j=1,lats_node_r
394  do i=1,lonr
395  phy_f2d(i,j,num_p2d) = 0.0
396  enddo
397  enddo
398  endif
399 
400  if (ras) call ras_init(levs, me)
401 !
402 !***************************************idea below*****************************************
403  if ( lsidea ) then
404 ! read the f10.7 and kp multi-time input data.
405 !---------------------------------------------
406  IF(.NOT. ASSOCIATED(f107)) ALLOCATE(f107(f107_kp_size))
407  IF(.NOT. ASSOCIATED(kp)) ALLOCATE(kp(f107_kp_size))
408  call read_wam_f107_kp_txt
409 !
410 !find PE with lat 1
411  pelat1 = -1
412  lanlat1 = -1
413  findlat1pe: do lan=1,lats_node_r
414  lat = global_lats_r(ipt_lats_node_r-1+lan)
415  if(lat == 1) then
416  pelat1 = me
417  lanlat1 = lan
418  print *,'pelat1=',pelat1
419  exit findlat1pe
420  endif
421  enddo findlat1pe
422  call mpi_reduce(pelat1,pelatall,1,mpi_integer,mpi_max,0,
423  & mpi_comm_all,info)
424  if(me == 0) then
425  print *,'pelatall=',pelatall
426  pelat1 = pelatall
427  endif
428  call mpi_bcast(pelat1,1,mpi_integer,0,mpi_comm_all,info)
429  print *,'pelat1=',pelat1,'lanlat1=',lanlat1,
430  & 'gen_coord_hybrid=',gen_coord_hybrid,
431  & 'thermodyn_id=',thermodyn_id
432 !
433 ! set plyr from lat1 pe
434  if(me == pelat1) then
435  do k=1,levs
436  plyr(k) = grid_fld%p(1,lanlat1,k)
437  enddo
438  print *,' plyr in gloopb ',(plyr(k),k=1,levs)
439  endif
440  call mpi_bcast(plyr,levs,mpi_r_io_r,pelat1,mpi_comm_all,info)
441  call idea_composition_init(levs,plyr)
442 !
443  call idea_solar_init(levs)
444  call idea_tracer_init(levs)
445 !h2ocin
446  print *,'in gloopb,nlevc_h2o=',nlevc_h2o,'k71=',k71,
447  & 'levs=',levs,levs-k71+1
448  allocate(prpa(nlevc_h2o))
449  prpa(1:nlevc_h2o) = 100.*pr_idea(k71:levs)
450  call h2ocin(prpa,nlevc_h2o,me,mpi_r_io_r,mpi_comm_all)
451  deallocate(prpa)
452 !ion
453  call efield_init()
454 !o3
455  call o3ini(levs)
456 !
457 !co2
458  if(.not. restart_step) then
459  if(me == pelat1) then
460 !compute phil
461  do n=1,ntrac
462  do k=1,levs
463  qtrac(k,n) = grid_fld%tracers(n)%flds(1,lanlat1,k)
464  enddo
465  enddo
466 ! print *,'in gloopb,ps=',grid_fld%ps(1,lanlat1)
467  call getphilvl(levs,ntrac, grid_fld%ps(1,lanlat1),
468  & grid_fld%t(1,lanlat1,1:levs),qtrac,
469  & grid_fld%dp(1,lanlat1,1:levs),gen_coord_hybrid,
470  & thermodyn_id,philco2,prsilvl1)
471 
472 !! change prsi from cb to pascal
473  prsilvl1 = prsilvl1*1000.
474 ! print *,'prsi=',prsi(1,103)*1000.
475 !!compute gravity
476  call gravco2(levs,philco2,sfc_fld%oro(1,lanlat1),gg1)
477  endif
478  call mpi_bcast(gg1,levs,mpi_r_io_r,pelat1,mpi_comm_all,info)
479 !
480  call mpi_bcast(prsilvl1,levs+1,mpi_r_io_r,pelat1,
481  & mpi_comm_all,info)
482 !
483  endif
484  pmod = pr_idea*100.
485  if(.not. allocated(prsilvl)) then
486  allocate(prsilvl(levs+1), gg(levs))
487  do k=1,levs
488  prsilvl(k) = prsilvl1(k)
489  gg(k) = gg1(k)
490  amgms(k) = amgm(k)
491  enddo
492  prsilvl(levs+1) = prsilvl1(levs+1)
493  endif
494 
495 ! print *,'bf co2cin,prlog=',prlog(k43:k43+3)
496  call co2cin(prlog(k43),pmod(k43),amgms(k43),gg(k43),
497  & nlev_co2,me,mpi_r_io_r,mpi_comm_all)
498 !
499  call ideaca_init(prsilvl,levs+1)
500 
501  endif
502 !***************************************idea above****************************************
503 
504  first = .false.
505 
506  gb_ini_time = gb_ini_time + (timef() - btime)
507 ! write(0,*)' gb_ini_time=',gb_ini_time*1.0e-3,' me=',me
508 
509  endif
510 !
511  if (semilag) then
512  nsphys = max(int(tstep/dtphys),1)
513  dtp = tstep / nsphys
514  dtf = dtp
515  else
516  nsphys = max(int((tstep+tstep)/dtphys+0.9999),1)
517  dtp = (tstep+tstep)/nsphys
518  dtf = 0.5*dtp
519  if(lsfwd) dtf = dtp
520  endif
521 
522  if (kdt == kdt_start .and. me == 0)
523  & write(0,*)' in gloopb nsphys=',nsphys
524  &, ' dtp=',dtp,' tstep=',tstep,' dtf=',dtf
525 
526 !
527  solhr = mod(phour+idate(1),cons_24)
528 
529 ! ************** Ken Campana / Yu-Tai Hou legacy Stuff ************
530 !... set switch for saving convective clouds
531  if(lscca .and. lsswr) then
532  clstp = 1100 + min(fhswr*facsw,fhour,cons_99) !initialize,accumulate,convert
533  elseif(lscca) then
534  clstp = 0100 + min(fhswr*facsw,fhour,cons_99) !accumulate,convert
535  elseif(lsswr) then
536  clstp = 1100 !initialize,accumulate
537  else
538  clstp = 0100 !accumulate
539  endif
540 ! ************** Ken Campana / Yu-Tai Hou legacy Stuff ************
541 !
542 !
543  if (imfdeepcnv <= 0 .or. cal_pre) then ! random number needed for RAS and old SAS
544  ! and when cal_pre=.true.
545 
546  if (random_clds) then
547  iseed = mod(100.0*sqrt(fhour*3600),cons_1p0d9) + 1 + seed0
548 
549 ! write(0,*)' After Initialization in gloopb iseed=',iseed
550 
551  nnrcm = nrcm*nsphys
552  allocate(wrkn(nnrcm), indxr(nnrcm))
553  call random_setseed(iseed)
554  call random_number(wrkn)
555  do nrc=1,nnrcm
556  indxr(nrc) = max(1, min(nint(wrkn(nrc)*maxran)+1,maxran))
557  enddo
558  endif
559  if (allocated(wrkn)) deallocate (wrkn)
560  endif
561 !
562 ! do ozone i/o and latitudinal interpolation to local gaussian lats
563 !
564  if (.not.use_nuopc) then
565  if (ntoz > 0) then
566  call ozinterpol(me,lats_node_r,lats_node_r,idate,fhour,
567  & jindx1,jindx2,ozplin,ozplout,ddy)
568  endif
569  if (h2o_phys) then
570  call h2ointerpol(me,lats_node_r,lats_node_r,idate,fhour,
571  & jindx1_h,jindx2_h,h2oplin,h2oplout,ddy_h)
572  endif
573  endif
574 !---------------------------main latitude loop starts---------------------------------
575 !
576 
577  do lan=1,lats_node_r
578  lat = global_lats_r(ipt_lats_node_r-1+lan)
579  lon_dim = lon_dims_r(lan)
580  lons_lat = lonsperlar(lat)
581 
582 !$omp parallel do private(i,j)
583  do n=1,ntrac
584  do i=1,lonr
585  trcp(i,n,1) = 0.0
586  trcp(i,n,2) = 0.0
587  enddo
588  enddo
589 
590 !$omp parallel do schedule(dynamic,1) private(lon)
591 !$omp+private(hprime_v,stc_v,smc_v,slc_v)
592 !$omp+private(nlons_v,sinlat_v,coslat_v,ozplout_v,rannum_v)
593 !$omp+private(h2oplout_v)
594 !$omp+private(prslk,prsl,prsik,prsi,phil,phii,dpshc)
595 !$omp+private(gu,gv,gt,gr,vvel)
596 !$omp+private(adt,adr,adu,adv,pgr,rqtk)
597 !$omp+private(qphys,sppt_wts)
598 !$omp+private(gt0,gr0,gu0,gv0,totprcp0,cnvprcp0,cplrain0,cplsnow0)
599 !$omp+private(phy_f2dv,dtdt,phy_fctdv)
600 !$omp+private(dt3dt_v,du3dt_v,dv3dt_v,dq3dt_v)
601 !$omp+private(upd_mfv,dwn_mfv,det_mfv)
602 !$omp+private(dqdt_v,cnvqc_v,hlwd)
603 !$omp+private(njeff,iblk,i,j,k,n,item,nn,nnr,tem)
604 !!$omp+private(njeff,iblk,i,j,k,n,item,nn,nnr,dbgu)
605 
606 !$omp+private(nirbmdi, nirdfdi,visbmdi, visdfdi, nirbmui, nirdfui)
607 !$omp+private(visbmui,visdfui,aoi_du,aoi_dv,aoi_dt,aoi_dq)
608 !$omp+private(aoi_dlw,aoi_dsw,aoi_dnirbm,aoi_dnirdf,aoi_dvisbm)
609 !$omp+private(aoi_dvisdf,aoi_rain,aoi_snow,aoi_nlw,aoi_nsw,aoi_nnirbm)
610 !$omp+private(aoi_nnirdf,aoi_nvisbm,aoi_nvisdf,aoi_dusfci)
611 !$omp+private(aoi_dvsfci,aoi_dtsfci,aoi_dqsfci,aoi_dlwsfci)
612 !$omp+private(aoi_dswsfci,aoi_dnirbmi,aoi_dnirdfi,aoi_dvisbmi)
613 !$omp+private(aoi_dvisdfi,aoi_nlwsfci,aoi_nswsfci,aoi_nnirbmi)
614 !$omp+private(aoi_nnirdfi,aoi_nvisbmi,aoi_nvisdfi,aoi_t2mi,aoi_q2mi)
615 !$omp+private(aoi_u10mi,aoi_v10mi,aoi_tseai,aoi_psurfi)
616 !$omp+private(aoi_slimskin,aoi_ulwsfcin)
617 !$omp+private(aoi_dusfcin,aoi_dvsfcin,aoi_dtsfcin,aoi_dqsfcin)
618 !$omp+private(nst_xt,nst_xs,nst_xu,nst_xv,nst_xz,nst_zm,nst_xtts)
619 !$omp+private(nst_xzts,nst_d_conv,nst_ifd,nst_dt_cool,nst_qrain)
620 !$omp+private(nst_tref,nst_z_c,nst_c_0,nst_c_d,nst_w_0,nst_w_d)
621 
622 !$omp+private(state_fldin,state_fldout)
623 !$omp+private(sfc_prop,diags,cld_prop,rad_tend)
624 !$omp+private(intrfc_fld,tbddata,dyn_parm)
625 !$omp+private(lonbnd)
626 
627 !!$omp+private(upd_mfv,dwn_mfv,det_mfv,dkh_v,rnp_v)
628 !!!$omp+private(njeff,iblk,i,j,k,n,item,nn,nnr,dbgu)
629 !!!$omp+private(njeff,iblk,ilan,i,j,k,n,item)
630 !!!!$omp+private(temlon,temlat,lprnt,ipt)
631 
632 
633 !---------------------------main longitude loop starts--------------------------------
634  do lon=1,lons_lat,ngptc
635 !!
636  njeff = min(ngptc,lons_lat-lon+1)
637  iblk = (lon-1)/ngptc + 1
638 !
639 ! dbgu = 1000 + lon
640 ! dbgu = 0
641 
642 ! write(dbgu,*)' GLOOPB : LON=',lon,' lons_lat=',lons_lat,
643 ! &' njeff=',njeff,' iblk=',iblk
644 ! &,' zorl=',sfc_fld%zorl(njeff-3:njeff,lan),' lan=',lan
645 ! write(dbgu,*)' dbgu=',dbgu,' lon=',lon,' lan=',lan
646 !!
647 ! write(dbgu,*)' lan=',lan,' pgr=',pgr(i),' i=',i,' njeff=',njeff
648 ! print *,' lan=',lan,' pgr=',pgr(i),' grid_gr=',grid_gr(ilan,g_ps)
649 ! &,' i=',i,' lan=',lan
650 !
651 
652  do k = 1, levs
653  do i = 1, njeff
654  item = lon+i-1
655  gu(i,k) = grid_fld%u(item,lan,k) ! real u
656  gv(i,k) = grid_fld%v(item,lan,k) ! real v
657  gt(i,k) = grid_fld%t(item,lan,k) ! temperature K
658  prsl(i,k) = grid_fld%p(item,lan,k) ! pascal
659  vvel(i,k) = grid_fld%dpdt(item,lan,k) ! pascal/sec
660  enddo
661  enddo
662 
663 ! write(0,*)' mintem=',minval(gt(1:njeff,:))
664 ! &,' maxtem=',maxval(gt(1:njeff,:)),' lan=',lan
665 
666 !hmhj prsi should be computed from model top for accuracy
667  do i=1,njeff
668  prsi(i,levs+1) = 0.0
669  prsik(i,levs+1) = 0.0
670  enddo
671  do k=levs,1,-1
672  do i = 1, njeff
673  item = lon+i-1
674  prsi(i,k) = prsi(i,k+1) + grid_fld%dp(item,lan,k) !pascal
675  enddo
676  enddo
677  do i=1,njeff
678  pgr(i) = prsi(i,1)
679  enddo
680  do n=1,ntrac
681  do k=1,levs
682  do i=1,njeff
683  gr(i,k,n)= grid_fld%tracers(n)%flds(lon+i-1,lan,k) ! kg/kg
684  enddo
685  enddo
686  enddo
687 
688  do i=1,njeff
689  phil(i,levs) = 0.0 ! will force calculation of geopotential in gbphys.
690  dpshc(i) = 0.3 * prsi(i,1)
691 !
692  nlons_v(i) = lons_lat
693  sinlat_v(i) = sinlat_r(lat)
694  coslat_v(i) = coslat_r(lat)
695  enddo
696 
697 ! if (gen_coord_hybrid .and. thermodyn_id == 3) then
698  do i=1,njeff
699  prslk(i,1) = 0.0 ! forces calculation of (p/p00)**kappa in gbphys
700  prsik(i,1) = 0.0 ! forces calculation of (p/p00)**kappa in gbphys
701  enddo
702 ! else
703 ! do k = 1, levs
704 ! do i = 1, njeff
705 ! prslk(i,k) = (prsl(i,k)*pt00001)**rk
706 ! prsik(i,k) = (prsi(i,k)*pt00001)**rk
707 ! enddo
708 ! enddo
709 ! endif
710 
711  if ( .not.use_nuopc) then
712  if (ntoz > 0) then
713  do j=1,pl_coeff
714  do k=1,levozp
715  do i=1,njeff
716  ozplout_v(i,k,j) = ozplout(k,lan,j)
717  enddo
718  enddo
719  enddo
720  endif
721  if (h2o_phys) then
722  do j=1,h2o_coeff
723  do k=1,levh2o
724  do i=1,ngptc
725  h2oplout_v(i,k,j) = h2oplout(k,lan,j)
726  enddo
727  enddo
728  enddo
729  endif
730  endif
731  do k=1,lsoil
732  do i=1,njeff
733  item = lon+i-1
734  smc_v(i,k) = sfc_fld%smc(k,item,lan)
735  stc_v(i,k) = sfc_fld%stc(k,item,lan)
736  slc_v(i,k) = sfc_fld%slc(k,item,lan)
737  enddo
738  enddo
739  do k=1,nmtvr
740  do i=1,njeff
741  hprime_v(i,k) = hprime(k,lon+i-1,lan)
742  enddo
743  enddo
744 !!
745  do j=1,ntot2d
746  do i=1,njeff
747  phy_f2dv(i,j) = phy_f2d(lon+i-1,lan,j)
748  enddo
749  enddo
750 
751  if(cscnv) then
752  do j=1,nctp
753  do i=1,njeff
754  phy_fctdv(i,j) = phy_fctd(lon+i-1,lan,j)
755  enddo
756  enddo
757  endif
758 
759  if (ldiag3d) then
760  do k=1,6
761  do j=1,levs
762  do i=1,njeff
763  dt3dt_v(i,j,k) = dt3dt(i,j,k,iblk,lan)
764  enddo
765  enddo
766  enddo
767  do k=1,4
768  do j=1,levs
769  do i=1,njeff
770  du3dt_v(i,j,k) = du3dt(i,j,k,iblk,lan)
771  dv3dt_v(i,j,k) = dv3dt(i,j,k,iblk,lan)
772  enddo
773  enddo
774  enddo
775  do k=1,5+pl_coeff
776  do j=1,levs
777  do i=1,njeff
778  dq3dt_v(i,j,k) = dq3dt(i,j,k,iblk,lan)
779  enddo
780  enddo
781  enddo
782  endif
783  if (ldiag3d .or. lgocart) then
784  do k=1,levs
785  do i=1,njeff
786  upd_mfv(i,k) = 0.
787  dwn_mfv(i,k) = 0.
788  det_mfv(i,k) = 0.
789  enddo
790  enddo
791  endif
792  if (lgocart) then
793  do k=1,levs
794  do i=1,njeff
795  dqdt_v(i,k) = 0.
796  cnvqc_v(i,k) = 0.
797  enddo
798  enddo
799  endif
800 !
801 !*************************************idea below**************************
802  if( lsidea ) then
803 !hmhj debug
804 ! do n=1,ntrac
805 ! call mymaxmin(gr(1,1,n),njeff,ngptc,levs,' gr b idea ')
806 ! enddo
807 !hmhj debug
808 ! call mymaxmin(gt,njeff,ngptc,levs,' gt before idea_phys ')
809  call idea_phys(njeff,ngptc,levs,prsi,prsl,
810  & gu,gv,gt,gr,ntrac,dtp,lat,
811  & solhr,slag,sdec,cdec,sinlat_v,coslat_v,
812  & xlon(lon,lan),xlat(lon,lan),
813  & sfc_fld%oro(lon,lan),flx_fld%coszen(lon,lan),
814  & swh(1,1,iblk,lan),hlw(1,1,iblk,lan),hlwd,
815  & thermodyn_id,sfcpress_id,gen_coord_hybrid,
816  & me,mpi_r_io_r,mpi_comm_all)
817 !hmhj debug
818 ! hlwd(:,:,6) = 0.0
819 ! do n=1,ntrac
820 ! print *,' =========== gr at n=',n,' ==============='
821 ! call mymaxmin(gr(1,1,n),njeff,ngptc,levs,' gr a idea ')
822 ! enddo
823 ! do n=1,6
824 ! print *,' =========== hlwd at n=',n,' ==============='
825 ! call mymaxmin(hlwd(1,1,n),njeff,ngptc,levs,' hlwd a idea ')
826 ! enddo
827 ! call mymaxmin(gu,njeff,ngptc,levs,' gu after idea_phys ')
828 ! call mymaxmin(gv,njeff,ngptc,levs,' gv after idea_phys ')
829 ! call mymaxmin(gt,njeff,ngptc,levs,' gt after idea_phys ')
830  endif
831 !*************************************idea above**************************
832 !
833 ! write(dbgu,*)' before gbphys:', njeff,ngptc,levs,lsoil,lsm, &
834 ! & ntrac,ncld,ntoz,ntcw, &
835 ! & nmtvr,nrcm,levozp,lonr,latr,jcap,num_p3d,num_p2d, &
836 ! & kdt,lat,me,pl_coeff,ncw,flgmin,crtrh,cdmbgwd
837 ! &,'sizer1=',size(phy_f3dv,1),'sizer2=',size(phy_f3dv,2)
838 ! &,'sizer3=',size(phy_f3dv,3)
839 ! &,' ccwf=',ccwf,' dlqf=',dlqf,' lsidea=',lsidea
840 ! &,' evpco=',evpco,' wminco=',wminco
841 ! write(dbgu,*)' tisfc=',sfc_fld%tisfc(1:20,lan),' lan=',lan,' lon='&
842 ! &, lon
843 !
844 !hmhj debug
845 ! call mymaxmin(gt,njeff,ngptc,levs,' gt before gbphys ')
846 !
847  do k=1,levs
848  do i=1,njeff
849  dtdt(i,k) = 0.0
850  enddo
851  enddo
852 
853 ! For tracer fixer
854 
855  do n=1,ntrac
856  if (fixtrc(n)) then
857  do k=1,levs
858  do i=1,njeff
859  trcp(lon+i-1,n,1) = trcp(lon+i-1,n,1)
860  & + gr(i,k,n) * (prsi(i,k) - prsi(i,k+1))
861  enddo
862  enddo
863  endif
864  enddo
865  do i=1,njeff
866  rqtk(i) = 0.0
867  enddo
868 
869  if (cplflx) then
870  do i=1,njeff
871  item = lon+i-1
872  nirbmdi(i) = aoi_fld%nirbmdi(item,lan)
873  nirdfdi(i) = aoi_fld%nirdfdi(item,lan)
874  visbmdi(i) = aoi_fld%visbmdi(item,lan)
875  visdfdi(i) = aoi_fld%visdfdi(item,lan)
876  nirbmui(i) = aoi_fld%nirbmui(item,lan)
877  nirdfui(i) = aoi_fld%nirdfui(item,lan)
878  visbmui(i) = aoi_fld%visbmui(item,lan)
879  visdfui(i) = aoi_fld%visdfui(item,lan)
880 !
881  aoi_du(i) = aoi_fld%dusfc(item,lan)
882  aoi_dv(i) = aoi_fld%dvsfc(item,lan)
883  aoi_dt(i) = aoi_fld%dtsfc(item,lan)
884  aoi_dq(i) = aoi_fld%dqsfc(item,lan)
885  aoi_dlw(i) = aoi_fld%dlwsfc(item,lan)
886  aoi_dsw(i) = aoi_fld%dswsfc(item,lan)
887  aoi_dnirbm(i) = aoi_fld%dnirbm(item,lan)
888  aoi_dnirdf(i) = aoi_fld%dnirdf(item,lan)
889  aoi_dvisbm(i) = aoi_fld%dvisbm(item,lan)
890  aoi_dvisdf(i) = aoi_fld%dvisdf(item,lan)
891  aoi_rain(i) = aoi_fld%rain(item,lan)
892  aoi_snow(i) = aoi_fld%snow(item,lan)
893  aoi_nlw(i) = aoi_fld%nlwsfc(item,lan)
894  aoi_nsw(i) = aoi_fld%nswsfc(item,lan)
895  aoi_nnirbm(i) = aoi_fld%nnirbm(item,lan)
896  aoi_nnirdf(i) = aoi_fld%nnirdf(item,lan)
897  aoi_nvisbm(i) = aoi_fld%nvisbm(item,lan)
898  aoi_nvisdf(i) = aoi_fld%nvisdf(item,lan)
899 !
900  aoi_slimskin(i) = aoi_fld%slimskin(item,lan)
901  aoi_dusfcin(i) = aoi_fld%dusfcin(item,lan)
902  aoi_dvsfcin(i) = aoi_fld%dvsfcin(item,lan)
903  aoi_dtsfcin(i) = aoi_fld%dtsfcin(item,lan)
904  aoi_dqsfcin(i) = aoi_fld%dqsfcin(item,lan)
905  aoi_ulwsfcin(i) = aoi_fld%ulwsfcin(item,lan)
906 
907  enddo
908  endif
909  if (nstf_name(1) > 0) then
910  do i=1,njeff
911  item = lon+i-1
912  nst_xt(i) = nst_fld%xt(item,lan)
913  nst_xs(i) = nst_fld%xs(item,lan)
914  nst_xu(i) = nst_fld%xu(item,lan)
915  nst_xv(i) = nst_fld%xv(item,lan)
916  nst_xz(i) = nst_fld%xz(item,lan)
917  nst_zm(i) = nst_fld%zm(item,lan)
918  nst_xtts(i) = nst_fld%xtts(item,lan)
919  nst_xzts(i) = nst_fld%xzts(item,lan)
920  nst_d_conv(i) = nst_fld%d_conv(item,lan)
921  nst_ifd(i) = nst_fld%ifd(item,lan)
922  nst_dt_cool(i) = nst_fld%dt_cool(item,lan)
923  nst_qrain(i) = nst_fld%qrain(item,lan)
924  nst_tref(i) = nst_fld%tref(item,lan)
925  nst_z_c(i) = nst_fld%z_c(item,lan)
926  nst_c_0(i) = nst_fld%c_0(item,lan)
927  nst_c_d(i) = nst_fld%c_d(item,lan)
928  nst_w_0(i) = nst_fld%w_0(item,lan)
929  nst_w_d(i) = nst_fld%w_d(item,lan)
930  enddo
931  endif
932 !
933 ! save a copy of the current state in order to calculate
934 ! physics tendencies for stochastic pertrubations (SPPT)
935  if (do_sppt)then
936  do k = 1,levs
937  do i = 1,njeff
938  gt0(i,k) = gt(i,k)
939  gu0(i,k) = gu(i,k)
940  gv0(i,k) = gv(i,k)
941  qphys_v(i,k) = 0.0
942  enddo
943  enddo
944  do n = 1,ntrac
945  do k = 1,levs
946  do i = 1,njeff
947  gr0(i,k,n) = gr(i,k,n)
948  enddo
949  enddo
950  enddo
951  do i = 1,njeff
952  cplrain0(i) = aoi_rain(i)
953  cplsnow0(i) = aoi_snow(i)
954  totprcp0(i) = flx_fld%totprcp(lon-1+i,lan)
955  cnvprcp0(i) = flx_fld%cnvprcp(lon-1+i,lan)
956  enddo
957  endif
958 
959  lonbnd=lon+njeff-1
960  if (use_nuopc) then
961  call rad_tend%setphys(
962  & swh(1:ngptc,1:levs,iblk,lan),
963  & sfalb(lon:lonbnd,lan),
964  & flx_fld%coszen(lon:lonbnd,lan),
965  & hlw(1:ngptc,1:levs,iblk,lan),
966  & flx_fld%tsflw(lon:lonbnd,lan),
967  & flx_fld%sfcemis(lon:lonbnd,lan),
968  & rqtk(1:njeff), hlwd(1:ngptc,1:levs,1:6),
969  & dtdt(1:ngptc,1:levs),
970  & swhc(1:ngptc,1:levs,iblk,lan),
971  & hlwc(1:ngptc,1:levs,iblk,lan)
972  & )
973  endif
974 
975  do nn=1,nsphys ! physics sub-steps
976 
977  if (imfdeepcnv <= 0 .or. cal_pre) then
978  if (random_clds) then
979  nnr = (nn-1)*nrcm
980  do j=1,nrcm
981  do i=1,njeff
982  rannum_v(i,j) = rannum_tank(lon+i-1,indxr(nnr+j),lan)
983  enddo
984  enddo
985  else
986  do j=1,nrcm
987  do i=1,njeff
988  rannum_v(i,j) = 0.6 ! This is useful for debugging
989  enddo
990  enddo
991  endif
992  endif
993 ! lonbnd=lon+njeff-1
994 
995  if (use_nuopc) then
1006 
1007  call dyn_parm%setphys(
1008  & xlon(lon:lonbnd,lan), xlat(lon:lonbnd,lan),
1009  & sinlat_v, coslat_v, solhr, ngptc, njeff, kdt,
1010  & lssav, lat, dtp, dtf, clstp, nn, nlons_v, fhour,
1011  & slag, sdec, cdec )
1012  call state_fldin%setphys( prsi, prsl, prslk, gt, gr,
1013  & vvel, pgr, gu, gv, prsik, phii, phil, adjtrc)
1014 
1015  call state_fldout%setphys(adt, adr, adu, adv)
1016 
1017  call diags%setphys(
1018  & flx_fld%srunoff (lon:lonbnd,lan),
1019  & flx_fld%evbsa (lon:lonbnd,lan),
1020  & flx_fld%evcwa (lon:lonbnd,lan),
1021  & flx_fld%snohfa (lon:lonbnd,lan),
1022  & flx_fld%transa (lon:lonbnd,lan),
1023  & flx_fld%sbsnoa (lon:lonbnd,lan),
1024  & flx_fld%snowca (lon:lonbnd,lan),
1025  & flx_fld%soilm (lon:lonbnd,lan),
1026  & flx_fld%tmpmin (lon:lonbnd,lan),
1027  & flx_fld%tmpmax (lon:lonbnd,lan),
1028 ! & flx_fld%slimsk (lon:lonbnd,lan),
1029  & flx_fld%dusfc (lon:lonbnd,lan),
1030  & flx_fld%dvsfc (lon:lonbnd,lan),
1031  & flx_fld%dtsfc (lon:lonbnd,lan),
1032  & flx_fld%dqsfc (lon:lonbnd,lan),
1033  & flx_fld%totprcp(lon:lonbnd,lan),
1034  & flx_fld%gflux (lon:lonbnd,lan),
1035  & flx_fld%dlwsfc (lon:lonbnd,lan),
1036  & flx_fld%ulwsfc (lon:lonbnd,lan),
1037  & flx_fld%suntim (lon:lonbnd,lan),
1038  & flx_fld%runoff (lon:lonbnd,lan),
1039  & flx_fld%ep (lon:lonbnd,lan),
1040  & flx_fld%cldwrk (lon:lonbnd,lan),
1041  & flx_fld%dugwd (lon:lonbnd,lan),
1042  & flx_fld%dvgwd (lon:lonbnd,lan),
1043  & flx_fld%psmean (lon:lonbnd,lan),
1044  & flx_fld%cnvprcp(lon:lonbnd,lan),
1045  & flx_fld%spfhmin(lon:lonbnd,lan),
1046  & flx_fld%spfhmax(lon:lonbnd,lan),
1047  & flx_fld%rain (lon:lonbnd,lan),
1048  & flx_fld%rainc (lon:lonbnd,lan),
1049 ! & flx_fld%snow (lon:lonbnd,lan),
1050  & dt3dt_v, dq3dt_v, du3dt_v, dv3dt_v, dqdt_v,
1051 ! Phys Outputs
1052  & flx_fld%u10m (lon:lonbnd,lan),
1053  & flx_fld%v10m (lon:lonbnd,lan),
1054  & flx_fld%zlvl (lon:lonbnd,lan),
1055  & flx_fld%psurf (lon:lonbnd,lan),
1056  & flx_fld%hpbl (lon:lonbnd,lan),
1057  & flx_fld%pwat (lon:lonbnd,lan),
1058  & flx_fld%t1 (lon:lonbnd,lan),
1059  & flx_fld%q1 (lon:lonbnd,lan),
1060  & flx_fld%u1 (lon:lonbnd,lan),
1061  & flx_fld%v1 (lon:lonbnd,lan),
1062  & flx_fld%chh(lon:lonbnd,lan),
1063  & flx_fld%cmm (lon:lonbnd,lan),
1064  & flx_fld%dlwsfci (lon:lonbnd,lan),
1065  & flx_fld%ulwsfci (lon:lonbnd,lan),
1066  & flx_fld%dswsfci (lon:lonbnd,lan),
1067  & flx_fld%uswsfci (lon:lonbnd,lan),
1068  & flx_fld%dusfci (lon:lonbnd,lan),
1069  & flx_fld%dvsfci (lon:lonbnd,lan),
1070  & flx_fld%dtsfci (lon:lonbnd,lan),
1071  & flx_fld%dqsfci (lon:lonbnd,lan),
1072  & flx_fld%gfluxi (lon:lonbnd,lan),
1073  & flx_fld%epi (lon:lonbnd,lan),
1074  & flx_fld%smcwlt2 (lon:lonbnd,lan),
1075  & flx_fld%smcref2 (lon:lonbnd,lan),
1076  & flx_fld%wet1 (lon:lonbnd,lan),
1077  & flx_fld%sr (lon:lonbnd,lan)
1078  & )
1079  call intrfc_fld%setphys(
1080  & flx_fld%sfcdsw (lon:lonbnd,lan),
1081  & flx_fld%sfcnsw (lon:lonbnd,lan),
1082  & flx_fld%sfcdlw (lon:lonbnd,lan),
1083  & nirbmui, nirdfui, visbmui, visdfui,
1084  & nirbmdi, nirdfdi, visbmdi, visdfdi, ! aoi
1085  & aoi_du, aoi_dv, aoi_dt, aoi_dq, aoi_dlw, aoi_dsw,
1086  & aoi_dnirbm, aoi_dnirdf, aoi_dvisbm, aoi_dvisdf,
1087  & aoi_rain, aoi_nlw, aoi_nsw, aoi_nnirbm, aoi_nnirdf,
1088  & aoi_nvisbm, aoi_nvisdf,
1089  & aoi_slimskin, aoi_ulwsfcin, aoi_dusfcin,
1090  & aoi_dvsfcin, aoi_dtsfcin, aoi_dqsfcin, aoi_snow,
1091  & nst_xt, nst_xs, nst_xu, nst_xv, nst_xz, nst_zm,
1092  & nst_xtts, nst_xzts, nst_d_conv,
1093  & nst_ifd, nst_dt_cool, nst_qrain,
1094  & aoi_dusfci, aoi_dvsfci, aoi_dtsfci,
1095  & aoi_dqsfci, aoi_dlwsfci, aoi_dswsfci,
1096  & aoi_dnirbmi, aoi_dnirdfi, aoi_dvisbmi,
1097  & aoi_dvisdfi, aoi_nlwsfci, aoi_nswsfci,
1098  & aoi_nnirbmi, aoi_nnirdfi, aoi_nvisbmi, aoi_nvisdfi,
1099  & aoi_t2mi, aoi_q2mi, aoi_u10mi, aoi_v10mi,
1100  & aoi_tseai, aoi_psurfi,
1101 ! & aoi_fld%oro (lon:lonbnd,lan),
1102 ! & aoi_fld%slimsk (lon:lonbnd,lan)
1103 ! & grid_fld%shum_wts(lon:lonbnd,lan,:),
1104  & grid_fld%sppt_wts(lon:lonbnd,lan,:),
1105 ! & grid_fld%skebu_wts(lon:lonbnd,lan,:),
1106 ! & grid_fld%skebv_wts(lon:lonbnd,lan,:),
1107 ! & grid_fld%vcu_wts (lon:lonbnd,lan,:),
1108 ! & grid_fld%vcv_wts (lon:lonbnd,lan,:),
1109 ! & uphys,vphys,tphys,
1110  & qphys_v,
1111 ! & tpphys,cpphys,
1112  & cplrain0,cplsnow0,
1113 ! & raincpl,snowcpl,
1114  & totprcp0,cnvprcp0,
1115  & gu0,gv0,gt0,gr0,
1116  & do_sppt, do_shum, do_skeb, do_vc
1117  & )
1118 
1119 ! call rad_tend%set(
1120 !! & swh(1,1,iblk,lan),
1121 ! & swh(1:ngptc,1:levs,iblk,lan),
1122 ! & sfalb(lon:lonbnd,lan),
1123 ! & flx_fld%coszen(lon:lonbnd,lan),
1124 !! & hlw(1,1,iblk,lan),
1125 ! & hlw(1:ngptc,1:levs,iblk,lan),
1126 ! & flx_fld%tsflw(lon:lonbnd,lan),
1127 ! & flx_fld%sfcemis(lon:lonbnd,lan),
1128 ! & rqtk=rqtk, hlwd=hlwd, dtdtr=dtdt,
1129 ! & swhc=swhc(1:ngptc,1:levs,iblk,lan),
1130 ! & hlwc=hlwc(1:ngptc,1:levs,iblk,lan)
1131 !! & swhc=swhc(1,1,iblk,lan),
1132 !! & hlwc=hlwc(1,1,iblk,lan)
1133 ! & )
1134 
1135  call sfc_prop%setphys(
1136  & hprime_v,
1137  & sfc_fld%slope (lon:lonbnd,lan),
1138  & sfc_fld%shdmin (lon:lonbnd,lan),
1139  & sfc_fld%shdmax (lon:lonbnd,lan),
1140  & sfc_fld%snoalb (lon:lonbnd,lan),
1141  & sfc_fld%tg3 (lon:lonbnd,lan),
1142  & sfc_fld%slmsk (lon:lonbnd,lan),
1143  & sfc_fld%vfrac (lon:lonbnd,lan),
1144  & sfc_fld%vtype (lon:lonbnd,lan),
1145  & sfc_fld%stype (lon:lonbnd,lan),
1146  & sfc_fld%uustar (lon:lonbnd,lan),
1147  & sfc_fld%oro (lon:lonbnd,lan),
1148  & sfc_fld%oro_uf (lon:lonbnd,lan),
1149  & sfc_fld%hice (lon:lonbnd,lan),
1150  & sfc_fld%fice (lon:lonbnd,lan),
1151  & sfc_fld%tisfc (lon:lonbnd,lan),
1152  & sfc_fld%tsea (lon:lonbnd,lan),
1153  & sfc_fld%snwdph (lon:lonbnd,lan),
1154  & sfc_fld%weasd (lon:lonbnd,lan), ! sheleg
1155  & sfc_fld%sncovr (lon:lonbnd,lan),
1156  & sfc_fld%zorl (lon:lonbnd,lan),
1157  & sfc_fld%canopy (lon:lonbnd,lan),
1158  & sfc_fld%ffmm (lon:lonbnd,lan),
1159  & sfc_fld%ffhh (lon:lonbnd,lan),
1160  & sfc_fld%f10m (lon:lonbnd,lan),
1161  & sfc_fld%t2m (lon:lonbnd,lan),
1162  & sfc_fld%q2m (lon:lonbnd,lan)
1163  & )
1164 
1165 
1166  call cld_prop%setphys(flgmin,sfc_fld%cv(lon:lonbnd,lan),
1167  & sfc_fld%cvt(lon:lonbnd,lan),
1168  & sfc_fld%cvb(lon:lonbnd,lan),
1169  & cnvqc_v, sup )
1170 
1171  call tbddata%set(
1172  & dpshc,
1173  & ozplout,
1174  & lan,
1175  & ozplout_v,
1176  & h2oplout,
1177  & h2oplout_v,
1178  & pl_pres,
1179  & rannum_v,
1180  & bkgd_vdif_m,
1181  & bkgd_vdif_h,
1182  & bkgd_vdif_s,
1183  & psautco, prautco, evpco,
1184  & wminco,
1185  & acv(lon:lonbnd,lan),
1186  & acvb(lon:lonbnd,lan),
1187  & acvt(lon:lonbnd,lan),
1188  & slc_v,
1189  & smc_v,
1190  & stc_v,
1191  & upd_mfv,
1192  & dwn_mfv,
1193  & det_mfv,
1194  & phy_f3d(1:ngptc,1:levs,1:ntot3d,iblk,lan),
1195 ! & phy_f3d(1,1,1,iblk,lan),
1196  & phy_f2dv,
1197  & sfc_fld%tprcp(lon:lonbnd,lan),
1198  & sfc_fld%srflag(lon:lonbnd,lan),
1199  & nst_tref,
1200  & nst_z_c,
1201  & nst_c_0,
1202  & nst_c_d,
1203  & nst_w_0,
1204  & nst_w_d,
1205  & fscav,
1206  & fswtr,
1207  & phy_fctdv
1208  & )
1209 
1210 ! if (savecon) then
1211 
1212 ! call phys_run_savein (state_fldin, sfc_prop,
1213 ! & diags, intrfc_fld, cld_prop, rad_tend,
1214 ! & mdl_parm, tbddata, dyn_parm)
1215 ! end if
1216 
1217  if (me == 0) then
1218  print *, 'NUOPC WRAPPER : Calling nuopc_phys_run...'
1219  end if
1220 
1223  call nuopc_phys_run(state_fldin, state_fldout, sfc_prop,
1224  & diags, intrfc_fld, cld_prop, rad_tend,
1225  & mdl_parm, tbddata, dyn_parm )
1226 
1227 ! if (savecon) then
1228 ! call phys_run_saveout (state_fldout, sfc_prop,
1229 ! & diags, intrfc_fld, cld_prop, rad_tend,
1230 ! & tbddata)
1231 ! end if
1232 
1233  else
1234 
1235 ! write(0,*)' calling gbphys for lan=',lan
1236 ! &,'dlw=', flx_fld%sfcdlw(lon,lan), flx_fld%tsflw (lon,lan)
1237 
1238  call gbphys &
1239 ! --- inputs
1240  & ( njeff,ngptc,levs,lsoil,lsm,ntrac,ncld,ntoz,ntcw,ntke, &
1241  & ntiw,ntlnc,ntinc, &
1242  & nmtvr,nrcm,levozp,lonr,latr,jcap, &
1243  & num_p3d,num_p2d,npdf3d,ncnvcld3d, &
1244  & kdt,lat,me,pl_coeff,nlons_v,ncw,flgmin,crtrh,cdmbgwd, &
1245  & ccwf,dlqf,ctei_rm,clstp,cgwf,prslrd0,ral_ts,dtp,dtf,fhour,&
1246  & solhr,slag,sdec,cdec,sinlat_v,coslat_v,pgr,gu,gv, &
1247  & gt,gr,vvel,prsi,prsl,prslk,prsik,phii,phil, &
1248  & rannum_v,ozplout_v,pl_pres,dpshc,fscav,fswtr, &
1249  & hprime_v, xlon(lon,lan),xlat(lon,lan), &
1250  & h2o_phys, levh2o, h2oplout_v, h2o_pres, h2o_coeff, &
1251  & isot, ivegsrc, &
1252  & sfc_fld%slope (lon,lan), sfc_fld%shdmin(lon,lan), &
1253  & sfc_fld%shdmax(lon,lan), sfc_fld%snoalb(lon,lan), &
1254  & sfc_fld%tg3 (lon,lan), sfc_fld%slmsk (lon,lan), &
1255  & sfc_fld%vfrac (lon,lan), sfc_fld%vtype (lon,lan), &
1256  & sfc_fld%stype (lon,lan), sfc_fld%uustar(lon,lan), &
1257  & sfc_fld%oro (lon,lan), sfc_fld%oro_uf(lon,lan), &
1258  & flx_fld%coszen(lon,lan), &
1259  & flx_fld%sfcdsw(lon,lan), flx_fld%sfcnsw(lon,lan), &
1260 
1261  & nirbmdi, nirdfdi, visbmdi, visdfdi, &
1262  & nirbmui, nirdfui, visbmui, visdfui, &
1263  & aoi_slimskin, aoi_ulwsfcin, &
1264  & aoi_dusfcin, aoi_dvsfcin, aoi_dtsfcin, aoi_dqsfcin, &
1265 
1266  & flx_fld%sfcdlw(lon,lan), flx_fld%tsflw (lon,lan), &
1267  & flx_fld%sfcemis(lon,lan), sfalb(lon,lan), &
1268  & swh(1,1,iblk,lan),swhc(1,1,iblk,lan), &
1269  & hlw(1,1,iblk,lan),hlwc(1,1,iblk,lan), hlwd, lsidea,
1270  & ras,pre_rad,ldiag3d,lgocart,lssav,cplflx, &
1271  & bkgd_vdif_m,bkgd_vdif_h,bkgd_vdif_s,psautco,prautco,evpco,&
1272  & wminco,pdfcld,shcnvcw,sup,redrag,hybedmf,dspheat, &
1273  & flipv,old_monin,cnvgwd,shal_cnv, &
1274  & imfshalcnv,imfdeepcnv,cal_pre,aero_in, &
1275  & mom4ice,mstrat,trans_trac,nstf_name,moist_adj, &
1276  & thermodyn_id,sfcpress_id,gen_coord_hybrid,levr,adjtrc,nn, &
1277  & cscnv,nctp,do_shoc,shocaftcnv,ntot3d,ntot2d, &
1278 ! --- input/outputs:
1279  & sfc_fld%hice (lon,lan), sfc_fld%fice (lon,lan), &
1280  & sfc_fld%tisfc (lon,lan), sfc_fld%tsea (lon,lan), &
1281  & sfc_fld%tprcp (lon,lan), sfc_fld%cv (lon,lan), &
1282  & sfc_fld%cvb (lon,lan), sfc_fld%cvt (lon,lan), &
1283  & sfc_fld%srflag(lon,lan), sfc_fld%snwdph(lon,lan), &
1284  & sfc_fld%weasd(lon,lan), sfc_fld%sncovr(lon,lan), &
1285  & sfc_fld%zorl (lon,lan), sfc_fld%canopy(lon,lan), &
1286  & sfc_fld%ffmm (lon,lan), sfc_fld%ffhh (lon,lan), &
1287  & sfc_fld%f10m (lon,lan), flx_fld%srunoff(lon,lan), &
1288  & flx_fld%evbsa (lon,lan), flx_fld%evcwa (lon,lan), &
1289  & flx_fld%snohfa(lon,lan), flx_fld%transa(lon,lan), &
1290  & flx_fld%sbsnoa(lon,lan), flx_fld%snowca(lon,lan), &
1291  & flx_fld%soilm (lon,lan), flx_fld%tmpmin(lon,lan), &
1292  & flx_fld%tmpmax(lon,lan), flx_fld%dusfc (lon,lan), &
1293  & flx_fld%dvsfc (lon,lan), flx_fld%dtsfc (lon,lan), &
1294  & flx_fld%dqsfc (lon,lan), flx_fld%totprcp(lon,lan), &
1295  & flx_fld%gflux (lon,lan), flx_fld%dlwsfc(lon,lan), &
1296  & flx_fld%ulwsfc(lon,lan), flx_fld%suntim(lon,lan), &
1297  & flx_fld%runoff(lon,lan), flx_fld%ep (lon,lan), &
1298  & flx_fld%cldwrk(lon,lan), flx_fld%dugwd (lon,lan), &
1299  & flx_fld%dvgwd (lon,lan), flx_fld%psmean(lon,lan), &
1300  & flx_fld%cnvprcp(lon,lan), flx_fld%spfhmin(lon,lan), &
1301  & flx_fld%spfhmax(lon,lan), &
1302  & flx_fld%rain(lon,lan), flx_fld%rainc(lon,lan), &
1303 
1304  & dt3dt_v, dq3dt_v, du3dt_v, dv3dt_v, dqdt_v,cnvqc_v, & ! added for GOCART
1305  & acv(lon,lan), acvb(lon,lan), acvt(lon,lan), &
1306  & slc_v, smc_v, stc_v, upd_mfv, dwn_mfv, det_mfv, &
1307  & phy_f3d(1,1,1,iblk,lan), phy_f2dv, &
1308  & aoi_du, aoi_dv, aoi_dt, aoi_dq, aoi_dlw, aoi_dsw, &
1309  & aoi_dnirbm, aoi_dnirdf, aoi_dvisbm, aoi_dvisdf, aoi_rain, &
1310  & aoi_nlw, aoi_nsw, aoi_nnirbm, aoi_nnirdf, &
1311  & aoi_nvisbm, aoi_nvisdf, aoi_snow, &
1312 
1313  & nst_xt, nst_xs, nst_xu, nst_xv, nst_xz, &
1314  & nst_zm, nst_xtts, nst_xzts, nst_d_conv, nst_ifd, &
1315  & nst_dt_cool, nst_qrain, &
1316  & nst_tref, nst_z_c, nst_c_0, nst_c_d, nst_w_0, nst_w_d, &
1317  & phy_fctdv, &
1318 
1319 ! --- outputs:
1320  & adt, adr, adu, adv, &
1321  & sfc_fld%t2m (lon,lan), sfc_fld%q2m (lon,lan), &
1322  & flx_fld%u10m (lon,lan), flx_fld%v10m (lon,lan), &
1323  & flx_fld%zlvl (lon,lan), flx_fld%psurf (lon,lan), &
1324  & flx_fld%hpbl (lon,lan), flx_fld%pwat (lon,lan), &
1325  & flx_fld%t1 (lon,lan), flx_fld%q1 (lon,lan), &
1326  & flx_fld%u1 (lon,lan), flx_fld%v1 (lon,lan), &
1327  & flx_fld%chh (lon,lan), flx_fld%cmm (lon,lan), &
1328  & flx_fld%dlwsfci(lon,lan), flx_fld%ulwsfci(lon,lan), &
1329  & flx_fld%dswsfci(lon,lan), flx_fld%uswsfci(lon,lan), &
1330  & flx_fld%dusfci(lon,lan), flx_fld%dvsfci(lon,lan), &
1331  & flx_fld%dtsfci(lon,lan), flx_fld%dqsfci(lon,lan), &
1332  & flx_fld%gfluxi(lon,lan), flx_fld%epi (lon,lan), &
1333  & flx_fld%smcwlt2(lon,lan), flx_fld%smcref2(lon,lan), &
1334  & flx_fld%wet1(lon,lan), flx_fld%sr(lon,lan), &
1335  & rqtk, &! rqtkD
1336  & dtdt, &
1337 
1338  & aoi_dusfci, aoi_dvsfci, aoi_dtsfci, aoi_dqsfci, &
1339  & aoi_dlwsfci, aoi_dswsfci, aoi_dnirbmi, aoi_dnirdfi, &
1340  & aoi_dvisbmi, aoi_dvisdfi, aoi_nlwsfci, aoi_nswsfci, &
1341  & aoi_nnirbmi, aoi_nnirdfi, aoi_nvisbmi, aoi_nvisdfi, &
1342  & aoi_t2mi, aoi_q2mi, aoi_u10mi, aoi_v10mi, &
1343  & aoi_tseai, aoi_psurfi &
1344  & )
1345 
1346  end if ! use_nuopc
1347 
1348  if (nn < nsphys) then
1349  do k=1,levs
1350  do i=1,njeff
1351  gt(i,k) = adt(i,k)
1352  gu(i,k) = adu(i,k)
1353  gv(i,k) = adv(i,k)
1354  enddo
1355  enddo
1356  do n=1,ntrac
1357  do k=1,levs
1358  do i=1,njeff
1359  gr(i,k,n) = adr(i,k,n)
1360  enddo
1361  enddo
1362  enddo
1363  else
1364  do n=1,ntrac
1365  if (fixtrc(n)) then
1366  do k=1,levs
1367  do i=1,njeff
1368  trcp(lon+i-1,n,2) = trcp(lon+i-1,n,2)
1369  & + adr(i,k,n) * (prsi(i,k) - prsi(i,k+1))
1370  enddo
1371  enddo
1372  endif
1373  enddo
1374  if (gg_tracers) then
1375  do i=1,njeff
1376  item = lon+i-1
1377  grid_fld%rqtk(item,lan) = rqtk(i) / pgr(i)
1378  enddo
1379 ! else
1380 ! do i=1,njeff
1381 ! item = lon+i-1
1382 ! grid_fld%rqtk(item,lan) = 0.0
1383 ! enddo
1384  endif
1385  endif
1386  enddo ! end of nsphys loop
1387 !-------------------------------------------------------------------------
1388  if (use_nuopc) then
1389 ! print*,'nuopc_sppt_phys'
1390  call nuopc_sppt_phys( state_fldout, diags, intrfc_fld, rad_tend,
1391  & mdl_parm, tbddata, dyn_parm )
1392 ! print*,'nuopc_sppt_phys done'
1393  else
1394 ! code section for SPPT stocastic perturbations
1395  if (do_sppt) then
1396  do k=1,levs
1397  do i=1,njeff
1398  sppt_wts = grid_fld%sppt_wts(lon+i-1,lan,k)
1399 ! perturb increments (adding radiation contribution back in)
1400 ! and put new grid back into adjusted variables
1401  adu(i,k) = gu0(i,k) + (adu(i,k)-gu0(i,k))*sppt_wts
1402  adv(i,k) = gv0(i,k) + (adv(i,k)-gv0(i,k))*sppt_wts
1403  qphys = gr0(i,k,1) + (adr(i,k,1)-gr0(i,k,1))*sppt_wts
1404 
1405  if (qphys > 0.0) then ! check for negative humidities
1406  tem = gt0(i,k) + dtdt(i,k)
1407  adt(i,k) = tem + (adt(i,k)-tem)*sppt_wts
1408  adr(i,k,1) = qphys
1409  endif
1410  enddo
1411  enddo
1412 ! precip perturbations
1413  do i=1,njeff
1414  item = lon+i-1
1415  sppt_wts = grid_fld%sppt_wts(item,lan,3)
1416  flx_fld%totprcp(item,lan) = totprcp0(i) + sppt_wts
1417  & * (flx_fld%totprcp(item,lan)-totprcp0(i))
1418  flx_fld%cnvprcp(item,lan) = cnvprcp0(i) + sppt_wts
1419  & * (flx_fld%cnvprcp(item,lan)-cnvprcp0(i))
1420  sfc_fld%tprcp(item,lan) = sfc_fld%tprcp(item,lan)
1421  & * sppt_wts
1422  aoi_rain(i) = cplrain0(i) + sppt_wts
1423  & * (aoi_rain(i)-cplrain0(i))
1424  aoi_snow(i) = cplsnow0(i) + sppt_wts
1425  & * (aoi_snow(i)-cplsnow0(i))
1426  nst_qrain(i) = nst_qrain(i)*sppt_wts
1427  enddo
1428 
1429  endif
1430 ! do_sppt
1431 ! code section for SHUM stocastic perturbations
1432  if (do_shum) then
1433  do k=1,levs
1434  do i=1,njeff
1435  adr(i,k,1) = adr(i,k,1)
1436  & * (1.0 + grid_fld%shum_wts(lon+i-1,lan,k))
1437  enddo
1438  enddo
1439  endif
1440 ! do_shum
1441 !!! code section for additive noise (SKEB) perturbation
1442  if (do_skeb) then
1443  do k=1,levs
1444  do i=1,njeff
1445  item = lon+i-1
1446  adu(i,k) = adu(i,k) + grid_fld%skebu_wts(item,lan,k)
1447  adv(i,k) = adv(i,k) + grid_fld%skebv_wts(item,lan,k)
1448  enddo
1449  enddo
1450  endif ! do_skeb
1451 !!! code section for vorticity confinement perturbation.
1452  if (do_vc) then
1453  do k=1,levs
1454  do i=1,njeff
1455  item = lon+i-1
1456  adu(i,k) = adu(i,k) + grid_fld%vcu_wts(item,lan,k)
1457  adv(i,k) = adv(i,k) + grid_fld%vcv_wts(item,lan,k)
1458  enddo
1459  enddo
1460  endif ! do_vc
1461 
1462 !! end of stochastic physics code
1463  endif
1464 !-------------------------------------------------------------------------
1465 
1466 ! if(kdt==100) then
1467 ! print *,'in gloopb,aft gbphys,kdt=',kdt,'lat=',lat,lon,'smcwlt=',
1468 ! & flx_fld%smcwlt2(lon:lon+3,lan),
1469 ! & 'loc=',minloc(flx_fld%smcwlt2(lon:lon+njeff-1,lan))
1470 ! endif
1471 !
1472 !!
1473 !hmhj debug
1474 ! do n=1,ntrac
1475 ! call mymaxmin(gr(1,1,n),njeff,ngptc,levs,' gr a gbphys ')
1476 ! call mymaxmin(adr(1,1,n),njeff,ngptc,levs,' adr a gbphys ')
1477 ! enddo
1478 ! write(0,*)' in gloopb dusfc=',flx_fld%dusfc(:,lan),' lan=',
1479 ! &lan,' kdt=',kdt
1480 
1481  do k=1,lsoil
1482  do i=1,njeff
1483  item = lon+i-1
1484  sfc_fld%smc(k,item,lan) = smc_v(i,k)
1485  sfc_fld%stc(k,item,lan) = stc_v(i,k)
1486  sfc_fld%slc(k,item,lan) = slc_v(i,k)
1487  enddo
1488  enddo
1489 !
1490  if (ldiag3d) then
1491  do k=1,6
1492  do j=1,levs
1493  do i=1,njeff
1494  dt3dt(i,j,k,iblk,lan) = dt3dt_v(i,j,k)
1495  enddo
1496  enddo
1497  enddo
1498  do k=1,4
1499  do j=1,levs
1500  do i=1,njeff
1501  du3dt(i,j,k,iblk,lan) = du3dt_v(i,j,k)
1502  dv3dt(i,j,k,iblk,lan) = dv3dt_v(i,j,k)
1503  enddo
1504  enddo
1505  enddo
1506  do k=1,5+pl_coeff
1507  do j=1,levs
1508  do i=1,njeff
1509  dq3dt(i,j,k,iblk,lan) = dq3dt_v(i,j,k)
1510  enddo
1511  enddo
1512  enddo
1513  do j=1,levs
1514  do i=1,njeff
1515  upd_mf(i,j,iblk,lan) = upd_mf(i,j,iblk,lan)+upd_mfv(i,j)
1516  dwn_mf(i,j,iblk,lan) = dwn_mf(i,j,iblk,lan)+dwn_mfv(i,j)
1517  det_mf(i,j,iblk,lan) = det_mf(i,j,iblk,lan)+det_mfv(i,j)
1518  enddo
1519  enddo
1520  endif
1521 !!
1522 !! total moist tendency (kg/kg/s): from local to global array
1523 !!
1524  if (lgocart) then
1525  tem = 1.0 /dtf
1526  do k=1,levs
1527  do i=1,njeff
1528  item = lon+i-1
1529  g3d_fld%dqdt(item,lan,k) = dqdt_v(i,k)
1530  g3d_fld%cnv_mfc(item,lan,k) = (upd_mfv(i,k)
1531  & + dwn_mfv(i,k)) * tem
1532  g3d_fld%cnv_mfd(item,lan,k) = det_mfv(i,k) * tem
1533  g3d_fld%cnv_qc(item,lan,k) = cnvqc_v(i,k)
1534  enddo
1535  enddo
1536  endif
1537 !
1538 
1539  if (cplflx) then
1540  do i=1,njeff
1541  item = lon+i-1
1542 !
1543  aoi_fld%dusfc(item,lan) = aoi_du(i)
1544  aoi_fld%dvsfc(item,lan) = aoi_dv(i)
1545  aoi_fld%dtsfc(item,lan) = aoi_dt(i)
1546  aoi_fld%dqsfc(item,lan) = aoi_dq(i)
1547  aoi_fld%dlwsfc(item,lan) = aoi_dlw(i)
1548  aoi_fld%dswsfc(item,lan) = aoi_dsw(i)
1549  aoi_fld%dnirbm(item,lan) = aoi_dnirbm(i)
1550  aoi_fld%dnirdf(item,lan) = aoi_dnirdf(i)
1551  aoi_fld%dvisbm(item,lan) = aoi_dvisbm(i)
1552  aoi_fld%dvisdf(item,lan) = aoi_dvisdf(i)
1553  aoi_fld%rain(item,lan) = aoi_rain(i)
1554  aoi_fld%snow(item,lan) = aoi_snow(i)
1555  aoi_fld%nlwsfc(item,lan) = aoi_nlw(i)
1556  aoi_fld%nswsfc(item,lan) = aoi_nsw(i)
1557  aoi_fld%nnirbm(item,lan) = aoi_nnirbm(i)
1558  aoi_fld%nnirdf(item,lan) = aoi_nnirdf(i)
1559  aoi_fld%nvisbm(item,lan) = aoi_nvisbm(i)
1560  aoi_fld%nvisdf(item,lan) = aoi_nvisdf(i)
1561 
1562  aoi_fld%dusfci(item,lan) = aoi_dusfci(i)
1563  aoi_fld%dvsfci(item,lan) = aoi_dvsfci(i)
1564  aoi_fld%dtsfci(item,lan) = aoi_dtsfci(i)
1565  aoi_fld%dqsfci(item,lan) = aoi_dqsfci(i)
1566  aoi_fld%dlwsfci(item,lan) = aoi_dlwsfci(i)
1567  aoi_fld%dswsfci(item,lan) = aoi_dswsfci(i)
1568  aoi_fld%dnirbmi(item,lan) = aoi_dnirbmi(i)
1569  aoi_fld%dnirdfi(item,lan) = aoi_dnirdfi(i)
1570  aoi_fld%dvisbmi(item,lan) = aoi_dvisbmi(i)
1571  aoi_fld%dvisdfi(item,lan) = aoi_dvisdfi(i)
1572 
1573  aoi_fld%nlwsfci(item,lan) = aoi_nlwsfci(i)
1574  aoi_fld%nswsfci(item,lan) = aoi_nswsfci(i)
1575  aoi_fld%nnirbmi(item,lan) = aoi_nnirbmi(i)
1576  aoi_fld%nnirdfi(item,lan) = aoi_nnirdfi(i)
1577  aoi_fld%nvisbmi(item,lan) = aoi_nvisbmi(i)
1578  aoi_fld%nvisdfi(item,lan) = aoi_nvisdfi(i)
1579 
1580  aoi_fld%t2mi(item,lan) = aoi_t2mi(i)
1581  aoi_fld%q2mi(item,lan) = aoi_q2mi(i)
1582  aoi_fld%u10mi(item,lan) = aoi_u10mi(i)
1583  aoi_fld%v10mi(item,lan) = aoi_v10mi(i)
1584  aoi_fld%tseai(item,lan) = aoi_tseai(i)
1585  aoi_fld%psurfi(item,lan) = aoi_psurfi(i)
1586 
1587  aoi_fld%tboti(item,lan) = adt(i,1)
1588  aoi_fld%qboti(item,lan) = adr(i,1,1)
1589  aoi_fld%uboti(item,lan) = adu(i,1)
1590  aoi_fld%vboti(item,lan) = adv(i,1)
1591  aoi_fld%pboti(item,lan) = prsl(i,1)
1592  aoi_fld%zboti(item,lan) = phil(i,1)/grav
1593 
1594  aoi_fld%oro(item,lan) = sfc_fld%oro(item,lan)
1595  aoi_fld%slimsk(item,lan) = sfc_fld%slmsk(item,lan)
1596  enddo
1597  endif
1598  if (nstf_name(1) > 0) then
1599  do i=1,njeff
1600  item = lon+i-1
1601  nst_fld%xt(item,lan) = nst_xt(i)
1602  nst_fld%xs(item,lan) = nst_xs(i)
1603  nst_fld%xu(item,lan) = nst_xu(i)
1604  nst_fld%xv(item,lan) = nst_xv(i)
1605  nst_fld%xz(item,lan) = nst_xz(i)
1606  nst_fld%zm(item,lan) = nst_zm(i)
1607  nst_fld%xtts(item,lan) = nst_xtts(i)
1608  nst_fld%xzts(item,lan) = nst_xzts(i)
1609  nst_fld%d_conv(item,lan) = nst_d_conv(i)
1610  nst_fld%ifd(item,lan) = nst_ifd(i)
1611  nst_fld%dt_cool(item,lan) = nst_dt_cool(i)
1612  nst_fld%qrain(item,lan) = nst_qrain(i)
1613 
1614  nst_fld%tref(item,lan) = nst_tref(i)
1615  nst_fld%z_c(item,lan) = nst_z_c(i)
1616  nst_fld%c_0(item,lan) = nst_c_0(i)
1617  nst_fld%c_d(item,lan) = nst_c_d(i)
1618  nst_fld%w_0(item,lan) = nst_w_0(i)
1619  nst_fld%w_d(item,lan) = nst_w_d(i)
1620  enddo
1621  endif
1622 !!
1623  do j=1,num_p2d+nshoc_2d
1624  do i=1,njeff
1625  phy_f2d(lon+i-1,lan,j) = phy_f2dv(i,j)
1626  enddo
1627  enddo
1628 
1629  if (cscnv) then
1630  do j=1,nctp
1631  do i=1,njeff
1632  phy_fctd(lon+i-1,lan,j) = phy_fctdv(i,j)
1633  enddo
1634  enddo
1635  endif
1636 
1637  do k=1,levs
1638  do i=1,njeff
1639  item = lon+i-1
1640  grid_fld%u(item,lan,k) = adu(i,k)
1641  grid_fld%v(item,lan,k) = adv(i,k)
1642  grid_fld%t(item,lan,k) = adt(i,k)
1643  enddo
1644  enddo
1645 
1646  do n=1,ntrac
1647  do k=1,levs
1648  do i=1,njeff
1649  grid_fld%tracers(n)%flds(lon+i-1,lan,k) = adr(i,k,n)
1650  enddo
1651  enddo
1652 ! if (n == ntrac .and. lan == lats_node_r)
1653 ! &write(6000+me,*) ' in gloob lan=',lan,' njeff=',njeff,
1654 ! &' adr=',adr(1:njeff,1,n),' kdt=',kdt
1655  enddo
1656 
1657 !hmhj debug
1658 ! do n=1,ntrac
1659 ! call mymaxmin(adr(1,1,n),njeff,ngptc,levs,' adr a gbphys ')
1660 ! enddo
1661 !!
1662  enddo !lon
1663 !---------------------------main longitude loop ends--------------------------------
1664 
1665  tem = 0.5 / lonsperlar(lat)
1666 !$omp parallel do private(j,n)
1667  do n=1,ntrac
1668  if (fixtrc(n)) then
1669  trcj(lan,n,1) = 0.0
1670  trcj(lan,n,2) = 0.0
1671  do j=1,lons_lat
1672  trcj(lan,n,1) = trcj(lan,n,1) + trcp(j,n,1)
1673  trcj(lan,n,2) = trcj(lan,n,2) + trcp(j,n,2)
1674  enddo
1675  trcj(lan,n,1) = trcj(lan,n,1) * tem
1676  trcj(lan,n,2) = trcj(lan,n,2) * tem
1677  endif
1678  enddo
1679 !
1680  enddo !lan
1681 !---------------------------main latitude loop ends---------------------------------
1682 !!
1683  call excht(lats_nodes_r, global_lats_r, trcj, trcg)
1684 
1685 
1686 !$omp parallel do private(n,lat,tem)
1687  do n=1,ntrac
1688  if (fixtrc(n)) then
1689  sumtrc(n,1) = 0.0
1690  sumtrc(n,2) = 0.0
1691  do lat=1,latr
1692  sumtrc(n,1) = sumtrc(n,1)
1693  & + wgt_r(min(lat,latr-lat+1))*trcg(lat,n,1)
1694  sumtrc(n,2) = sumtrc(n,2)
1695  & + wgt_r(min(lat,latr-lat+1))*trcg(lat,n,2)
1696 ! write(0,*)' kdt=',kdt,' lat=',lat,' sumtrc=',sumtrc(n,:)
1697 ! &,' trcg=',trcg(lat,n,1),trcg(lat,n,2),' n=',n,' me=',me
1698  enddo
1699 
1700 ! write(0,*)' kdt=',kdt,' n=',n,' sumtrc=',sumtrc(n,:)
1701 ! &,' kdt_start=',kdt_start,' me=',me
1702 
1703  if (abs(sumtrc(n,1)) > 1.0e-12 .and. kdt > kdt_start+1) then
1704 ! if (abs(sumtrc(n,1)) > 1.0e-12 .and. kdt > 1
1705 ! & .and. n == 2) then
1706  tem = (sumtrc(n,3)-sumtrc(n,1)) / sumtrc(n,1)
1707  if (abs(tem) < 0.1) then
1708  adjtrc(n) = 1 + tem
1709  else
1710  adjtrc(n) = 1.0
1711  endif
1712  else
1713  adjtrc(n) = 1.0
1714  endif
1715  sumtrc(n,3) = sumtrc(n,2)
1716 ! write(0,*)' kdt=',kdt,' n=',n,' adjtrc=',adjtrc(n)
1717 ! &, ' sumtrc=',sumtrc(n,:),' me=',me
1718  endif
1719  enddo
1720 
1721 !
1722  if (allocated(indxr)) deallocate (indxr)
1723  return
1724 !...................................
1725  end subroutine gloopb
1726 !-----------------------------------
DDT for fields typically only used for diagnostic output.
DDT for basic inputs of radiation and physics parameters.
DDT for radiation tendencies.
DDT for surface fields.
DDT for basic outputs from radiation and physics.
subroutine gbphys
Parameter descriptions include intent, name, description, and size.
Definition: gbphys.f:530
DDT for data that changes frequently (e.g. inner loops)
DDT for data used for coupling (e.g. land and ocean)
DDT for data that has not been pigeonholed and is left to be determined.
DDT for non-changing model parameters - set once in initialize.
DDT for cloud data and parameters. Used by grrad and gbphys with different intent.
subroutine gloopb(grid_fld, g3d_fld, sfc_fld,
Definition: gloopb.f:15
the interface between the dynamic core and the physics packages