Interoperable Physics Driver for NGGPS
gfs_physics_run_mod.f
Go to the documentation of this file.
1 !
2 ! !module: gfs_physics_run_mod --- run module of the grided
3 ! component of the gfs physics.
4 !
5 ! !description: gfs run module.
6 !
7 ! !revision history:
8 !
9 ! november 2004 weiyu yang initial code.
10 ! may 2005 weiyu yang, updated to the new version of gfs.
11 ! janusry 2007 hann-ming henry juang for gfs dynamics only
12 ! july 2007 shrinivas moorthi for gfs physics only
13 ! november 2007 hann-ming henry juang continue for gfs physics
14 ! october 2009 jun wang add nsout option
15 ! oct 11 2009 sarah lu, grid_gr replaced by grid_fld
16 ! oct 17 2009 sarah lu, q is replaced by tracers(1)
17 ! dec 08 2009 sarah lu, add g3d_fld to do_physics_one_step
18 ! calling argument
19 ! July 2010 Shrinivas Moorthi - Updated for new physics and added nst
20 ! eliminated calls to common_vars
21 ! jul 21 2010 sarah lu, add g2d_fld to do_physics_one_step
22 ! calling argument
23 ! Aug 03 2010 jun wang, fix lsout for dfi
24 ! Aug 25 2010 Jun Wang, add zhour_dfi for filtered dfi fields output
25 ! Oct 18 2010 s. moorthi added fscav to do tstep
26 ! Dec 23 2010 Sarah Lu, setup fscav from gfs_phy_tracer
27 ! Nov 27 2011 Sarah Lu, zerout fcld, dqdt, and wet1
28 ! Apr 06 2012 Henry Juang, add idea
29 ! Apr 09 2012 Jun Wang save phys state at 3hr and set back to
30 ! 3hr phys state aft dfi
31 ! Mar 09 2013 Jun Wang add restart step for idea
32 ! Mar 25 2014 Xingren Wu add aoi_fld for A/O/I coupling
33 ! Mar 31 2014 S Moorthi Add sstForGSM to do_physics_onestep argument
34 ! Jul -- 2014 S Moorthi update for new physics and semilag
35 ! Sep 18 2014 S Moorthi simplified argments for dfi_fixwr
36 ! Sep 30 2014 Sarah Lu, Remove fscav array
37 ! May -- 2015 S Moorthi - added SHOC option
38 ! Jun 09 2015 G Theurich, Generalize importData handling
39 ! Jul 29 2015 S Moorthi - added high frequency output capability
40 ! Jan 2016 P Tripp - NUOPC/GSM merge - importData
41 ! feb 2016 S Moorthi - grid-point digital filter fix to filter at initial
42 ! time and during restart
43 ! mar 05 2016 S Moorthi - add logic to change phour back to reset time after
44 ! clock was reset after digital filter initialization
45 ! Jun 10 2016 S Moorthi - add stratospheric h2o parameterization
46 ! Jun 2016 X.Li - modify to read dfi_fixwr (with flx_fld saved)
47 !
48 !
49 ! !interface:
50 !
52 !
53 !!uses:
54 !
55  use gfs_physics_internal_state_mod, ONLY: gfs_physics_internal_state
56  USE date_def, ONLY: fhour, zhour
57  USE namelist_physics_def, ONLY: nsout,ldfi,ndfi,nsout_hf, fhmax_hf
58  use gfs_phy_tracer_config, ONLY: gfs_phy_tracer
59  use layout1, only: me
60  use resol_def, only: kdt_start
61 
62  implicit none
63 
64  contains
65 
66  subroutine gfs_physics_run(gis_phy, rc)
67 
68  type(gfs_physics_internal_state), pointer, intent(inout) :: gis_phy
69  integer, optional, intent(out) :: rc
70 
71  real , save :: timestep=0.0
72  integer rc1, k , i1, i2, i3, kdt_dif
73  logical lsout1
74  real tem
75 
76 !***********************************************************************
77 !
78 ! lsfwd logical true during a first forward step
79 ! lssav logical true during a step for which
80 ! diagnostics are accumulated
81 ! lscca logical true during a step for which convective clouds
82 ! are calculated from convective precipitation rates
83 ! phour real forecast hour at the end of the previous timestep
84 !
85 ! lsout controls freq. of output
86 ! nsout controls freq. of output in time steps
87 ! fhout controls freq. of output in hours
88 ! nszer time steps between zeroing fluxes
89 !
90 !***********************************************************************
91 
92 ! print *,' enter gfs_physics_run '
93 !
94  rc1 = 0
95 !
96 ! print *,' uug=',gis_phy%grid_gr(1,gis_phy%g_u:gis_phy%g_u+gis_phy%levs-1)
97 ! print *,' pg=',gis_phy%grid_gr(1,gis_phy%g_p:gis_phy%g_p+gis_phy%levs-1)
98 ! print *,' dpg=',gis_phy%grid_gr(1,gis_phy%g_dp:gis_phy%g_dp+gis_phy%levs-1)
99 !
100 ! ---------------------------------------------------------------------
101 ! ======================================================================
102 ! do one physics time step
103 ! ---------------------------------------------------------------------
104 ! write(0,*)' in gfs_physics_run kdt=',gis_phy%kdt,' nsout=',nsout &
105 ! ,' kdt_start=',kdt_start
106 
107  if (fhmax_hf > 0 .and. nsout_hf > 0 .and. gis_phy%phour <= fhmax_hf &
108  .and. gis_phy%kdt*gis_phy%deltim <= fhmax_hf*3600.) then
109  lsout1 = mod(gis_phy%kdt ,nsout_hf) == 0
110  else
111  lsout1 = mod(gis_phy%kdt ,nsout) == 0
112  endif
113 
114  kdt_dif = gis_phy%kdt - kdt_start
115  if (.not. ldfi) then
116  gis_phy%lsout = lsout1 .or. gis_phy%kdt == 1
117  else
118  gis_phy%lsout = lsout1 .and. &
119  (kdt_dif <= ndfi/2 .or. kdt_dif > ndfi) .or. gis_phy%kdt == 1
120  endif
121 
122 ! if (me == 0) &
123 ! write(0,*)' in phy nsout_hf=',nsout_hf,' fhmax_hf=',fhmax_hf
124 ! write(0,*)' in phy ' ,&
125 ! 'gis_phy%lsout=',gis_phy%lsout,' kdt=',gis_phy%kdt
126 !
127 ! print *,' end of common_to_physics_vars,kdt=',gis_phy%kdt, &
128 ! 'nsout=',nsout,'lsout=',gis_phy%LSOUT,'zhour=',gis_phy%ZHOUR, &
129 ! 'ldfi=',ldfi,'ndfi=',ndfi,gis_phy%kdt<=ndfi/2,gis_phy%kdt>ndfi, &
130 ! gis_phy%kdt<=ndfi/2.or.gis_phy%kdt>ndfi
131 ! if(gis_phy%kdt==12.and.gis_phy%kdt<=13.or.gis_phy%kdt>=24.and.gis_phy%kdt<=25) then
132 ! print *,'be phys one,kdt=',gis_phy%kdt,'ps=',maxval(gis_phy%grid_fld%ps), &
133 ! minval(gis_phy%grid_fld%ps),'t=',maxval(gis_phy%grid_fld%t), &
134 ! minval(gis_phy%grid_fld%t),'spfh=',maxval(gis_phy%grid_fld%tracers(1)%flds), &
135 ! minval(gis_phy%grid_fld%tracers(1)%flds),'tsea=',maxval(gis_phy%sfc_fld%tsea),&
136 ! minval(gis_phy%sfc_fld%tsea),maxloc(gis_phy%sfc_fld%tsea),maxloc(gis_phy%grid_fld%ps)
137 ! print *,' ps1lp(',gis_phy%kdt,')= ',gis_phy%grid_fld%ps(154,58)
138 ! endif
139 
140  if( ndfi > 0 .and. kdt_dif == ndfi/2+1 .and. .not. ldfi ) then
141  call dfi_fixwr(2, gis_phy%sfc_fld, gis_phy%flx_fld, gis_phy%nst_fld)
142  endif
143 
144 
145  if ( gis_phy%lgocart ) then
146 ! i1 = gis_phy%gfs_phy_tracer%ntrac_met+1 ! 1st chemical tracer (excluding o3)
147 ! i2 = gis_phy%gfs_phy_tracer%ntrac ! last chemical tracer (excluding o3)
148  i1 = gfs_phy_tracer%ntrac_met+1 ! 1st chemical tracer (excluding o3)
149  i2 = gfs_phy_tracer%ntrac ! last chemical tracer (excluding o3)
150 !
151  i1 = gis_phy%lonr
152  i2 = gis_phy%lats_node_r_max
153  i3 = gis_phy%levs
154  gis_phy%flx_fld%wet1(1:i1,1:i2) = 0.
155  gis_phy%g3d_fld%fcld(1:i1,1:i2,1:i3) = 0.
156  gis_phy%g3d_fld%dqdt(1:i1,1:i2,1:i3) = 0.
157  gis_phy%g3d_fld%cnv_mfc(1:i1,1:i2,1:i3) = 0.
158  gis_phy%g3d_fld%cnv_mfd(1:i1,1:i2,1:i3) = 0.
159  gis_phy%g3d_fld%cnv_qc(1:i1,1:i2,1:i3) = 0.
160  endif
161 
162  tem = (gis_phy%kdt-1) * gis_phy%deltim / 3600.0
163  if (gis_phy%phour > tem) gis_phy%phour = tem
164 
165 ! gis_phy%phour = gis_phy%kdt * gis_phy%deltim / 3600.
166 ! if (me ==0) write(0,*)'gis_phy%phour=',gis_phy%phour,'gis_phy%kdt=',&
167 ! gis_phy%kdt,' gis_phy%deltim',gis_phy%deltim
168 !
169 ! ======================================================================
170  call do_physics_one_step( &
171  gis_phy%deltim, gis_phy%kdt, gis_phy%phour, &
172  gis_phy%grid_fld, gis_phy%sfc_fld, gis_phy%flx_fld, &
173  gis_phy%nst_fld, gis_phy%g3d_fld, gis_phy%g2d_fld, &
174  gis_phy%aoi_fld, gis_phy%importData, &
175  gis_phy%lats_nodes_r, gis_phy%global_lats_r, &
176  gis_phy%lonsperlar, &
177  gis_phy%XLON, gis_phy%XLAT, gis_phy%COSZDG, &
178  gis_phy%HPRIME, gis_phy%SWH, gis_phy%swhc, &
179  gis_phy%HLW, gis_phy%hlwc, &
180 ! idea add by hmhj - commented by moorthi since unused
181 ! gis_phy%HTRSWB, gis_phy%HTRLWB, &
182  gis_phy%FLUXR, gis_phy%SFALB, &
183  gis_phy%SLAG, gis_phy%SDEC, gis_phy%CDEC, &
184  gis_phy%OZPLIN, gis_phy%JINDX1, gis_phy%JINDX2, &
185  gis_phy%DDY, gis_phy%h2oplin, gis_phy%jindx1_h, &
186  gis_phy%jindx2_h, gis_phy%ddy_h, &
187  gis_phy%phy_f3d, gis_phy%phy_f2d, gis_phy%phy_fctd, &
188  gis_phy%num_ctp, &
189  gis_phy%NBLCK, gis_phy%ZHOUR_DFI, &
190 ! gis_phy%NBLCK, gis_phy%ZHOUR, gis_phy%ZHOUR_DFI, &
191  gis_phy%N3, gis_phy%N4, &
192  gis_phy%LSOUT, gis_phy%COLAT1, gis_phy%CFHOUR1, &
193  gis_phy%restart_step, gis_phy%mdl_parm)
194 
195 !
196 ! =======================================================================
197 !
198 ! save phys fields for digital filter
199 !
200  if( ldfi .and. kdt_dif == ndfi/2 ) then
201 ! write(0,*)'save phys state, at gis_phy%kdt=',gis_phy%kdt,'ldfi=',ldfi
202  call dfi_fixwr(1, gis_phy%sfc_fld, gis_phy%flx_fld, gis_phy%nst_fld)
203  endif
204 !
205  gis_phy%phour = fhour ! update hour
206  gis_phy%zhour = zhour ! update hour
207 !
208  if(present(rc)) then
209  rc = rc1
210  end if
211 
212  end subroutine gfs_physics_run
213 
214  end module gfs_physics_run_mod
subroutine do_physics_one_step(deltim, kdt, PHOUR,
subroutine gfs_physics_run(gis_phy, rc)