CCPP SciDoc v7.0.x  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
noahmpdrv.F90
1#ifndef CCPP
2#define CCPP
3#endif
4
6
12
14 module noahmpdrv
15
16 use module_sf_noahmplsm
17
18! These hold and apply Land IAU increments for soil temperature
19! (possibly will extend to soil moisture increments)
21 land_iau_mod_init, land_iau_mod_getiauforcing, land_iau_mod_finalize, calculate_landinc_mask
22
23 implicit none
24
25 integer, parameter :: psi_opt = 0 ! 0: MYNN or 1:GFS
26
27 private
28
29 public :: noahmpdrv_init, noahmpdrv_run, &
31
32 contains
33
40 subroutine noahmpdrv_init(lsm, lsm_noahmp, me, isot, ivegsrc, &
41 nlunit, pores, resid, &
42 do_mynnsfclay,do_mynnedmf, &
43 errmsg, errflg, &
44 Land_IAU_Control, Land_IAU_Data, Land_IAU_state)
45
46 use machine, only: kind_phys
50
51 implicit none
52 integer, intent(in) :: lsm
53 integer, intent(in) :: lsm_noahmp
54 integer, intent(in) :: me, isot, ivegsrc, nlunit
55
56 real (kind=kind_phys), dimension(:), intent(out) :: pores, resid
57
58 logical, intent(in) :: do_mynnsfclay
59 logical, intent(in) :: do_mynnedmf
60
61
62 character(len=*), intent(out) :: errmsg
63 integer, intent(out) :: errflg
64
65 ! Land iau mod DDTs ! made optional to allow NoahMP Component model call this function without having to deal with IAU
66
67 ! Land IAU Control holds settings' information, maily read from namelist
68 ! (e.g., block of global domain that belongs to current process,
69 ! whether to do IAU increment at this time step, time step informatoin, etc)
70 type(land_iau_control_type), intent(inout), optional :: Land_IAU_Control
71
72 ! land iau state holds increment data read from file (before interpolation)
73 type(land_iau_state_type), intent(inout), optional :: Land_IAU_state
74
75 ! Land IAU Data holds spatially and temporally interpolated increments per time step
76 type(land_iau_external_data_type), intent(inout), optional :: Land_IAU_Data ! arry of (number of blocks):each proc holds nblks
77
78 ! Initialize CCPP error handling variables
79 errmsg = ''
80 errflg = 0
81
82 ! Consistency checks
83 if (lsm/=lsm_noahmp) then
84 write(errmsg,'(*(a))') 'Logic error: namelist choice of ', &
85 & 'LSM is different from Noah'
86 errflg = 1
87 return
88 end if
89
90 if (ivegsrc /= 1) then
91 errmsg = 'The NOAHMP LSM expects that the ivegsrc physics '// &
92 'namelist parameter is 1. Exiting...'
93 errflg = 1
94 return
95 end if
96 if (isot /= 1) then
97 errmsg = 'The NOAHMP LSM expects that the isot physics '// &
98 'namelist parameter is 1. Exiting...'
99 errflg = 1
100 return
101 end if
102
103 if ( do_mynnsfclay .and. .not. do_mynnedmf) then
104 errmsg = 'Problem : do_mynnsfclay = .true.' // &
105 'but mynnpbl is .false.. Exiting ...'
106 errflg = 1
107 return
108 end if
109
110
111 !--- initialize soil vegetation
112 call set_soilveg(me, isot, ivegsrc, nlunit, errmsg, errflg)
113 if(errflg/=0) return
114
115 !--- read in noahmp table
116 call read_mp_table_parameters(errmsg, errflg)
117 if(errflg/=0) return
118
119 ! initialize psih and psim
120 if ( do_mynnsfclay ) then
121 call psi_init(psi_opt,errmsg,errflg)
122 if(errflg/=0) return
123 endif
124
125 pores(:) = maxsmc(:)
126 resid(:) = drysmc(:)
127
128 if (present(land_iau_control) .and. present(land_iau_data) .and. present(land_iau_state)) then
129
130 ! Initialize IAU for land--land_iau_control was set by host model
131 if (.not. land_iau_control%do_land_iau) return
132 call land_iau_mod_init (land_iau_control, land_iau_data, land_iau_state, errmsg, errflg)
133
134 endif
135
136 end subroutine noahmpdrv_init
137
144subroutine noahmpdrv_timestep_init (itime, fhour, delt, km, ncols, &
145 isot, ivegsrc, soiltyp, vegtype, weasd, &
146 land_iau_control, land_iau_data, land_iau_state, &
147 stc, slc, smc, errmsg, errflg, &
148 con_g, con_t0c, con_hfus)
149
150 use machine, only: kind_phys
152 ! use set_soilveg_snippet_mod, only: set_soilveg_noahmp
153 use noahmp_tables
154
155 implicit none
156
157 integer , intent(in) :: itime !current forecast iteration
158 real(kind=kind_phys) , intent(in) :: fhour !current forecast time (hr)
159 real(kind=kind_phys) , intent(in) :: delt ! time interval [s]
160 integer , intent(in) :: km !vertical soil layer dimension
161 integer, intent(in) :: ncols
162 integer, intent(in) :: isot
163 integer, intent(in) :: ivegsrc
164
165 integer , dimension(:) , intent(in) :: soiltyp ! soil type (integer index)
166 integer , dimension(:) , intent(in) :: vegtype ! vegetation type (integer index)
167 real(kind=kind_phys), dimension(:) , intent(inout) :: weasd ! water equivalent accumulated snow depth [mm]
168
169 type(land_iau_control_type) , intent(inout) :: land_iau_control
170 type(land_iau_external_data_type) , intent(inout) :: Land_IAU_Data
171 type(land_iau_state_type) , intent(inout) :: Land_IAU_State
172 real(kind=kind_phys), dimension(:,:) , intent(inout) :: stc ! soiltemp [K]
173 real(kind=kind_phys), dimension(:,:) , intent(inout) :: slc !liquid soil moisture [m3/m3]'
174 real(kind=kind_phys), dimension(:,:) , intent(inout) :: smc !
175 character(len=*), intent(out) :: errmsg
176 integer, intent(out) :: errflg
177 real(kind=kind_phys), intent(in) :: con_g ! grav
178 real(kind=kind_phys), intent(in) :: con_t0c ! tfreez
179 real(kind=kind_phys), intent(in) :: con_hfus ! hfus
180
181 ! IAU update
182 real(kind=kind_phys),allocatable, dimension(:,:) :: stc_inc_flat, slc_inc_flat
183 real(kind=kind_phys), dimension(km) :: dz ! layer thickness
184
185!TODO: This is hard-coded in noahmpdrv
186 real(kind=kind_phys) :: zsoil(4) = (/ -0.1, -0.4, -1.0, -2.0 /) !zsoil(km)
187
188 integer :: lsoil_incr
189 integer, allocatable :: mask_tile(:)
190 integer,allocatable :: stc_updated(:), slc_updated(:)
191 logical :: soil_freeze, soil_ice
192 integer :: soiltype, n_stc, n_slc
193 real(kind=kind_phys) :: slc_new
194
195 integer :: i, j, ij, l, k, ib
196 integer :: lensfc
197
198 real(kind=kind_phys) :: smp
199 real(kind=kind_phys) :: hc_incr
200
201 integer :: nother, nsnowupd
202 integer :: nstcupd, nslcupd, nfrozen, nfrozen_upd
203 logical :: print_update_stats = .false.
204
205 ! --- Initialize CCPP error handling variables
206 errmsg = ''
207 errflg = 0
208
209 if (.not. land_iau_control%do_land_iau) return
210
212 land_iau_control%fhour=fhour
213
215 call land_iau_mod_getiauforcing(land_iau_control, land_iau_data, land_iau_state, errmsg, errflg)
216 if (errflg .ne. 0) then
217 return
218 endif
219
221 if (.not. land_iau_data%in_interval) then
222 return
223 endif
224
225 if(land_iau_control%me == land_iau_control%mpi_root) then
226 print*, "adding land iau increments "
227 endif
228
229 if (land_iau_control%lsoil .ne. km) then
230 write(errmsg,*) 'noahmpdrv_timestep_init: Land_IAU_Data%lsoil ',land_iau_control%lsoil,' not equal to km ',km
231 errflg = 1
232 return
233 endif
234
235 ! local variable to copy blocked data Land_IAU_Data%stc_inc
236 allocate(stc_inc_flat(land_iau_control%nx * land_iau_control%ny, km)) !GFS_Control%ncols
237 allocate(slc_inc_flat(land_iau_control%nx * land_iau_control%ny, km)) !GFS_Control%ncols
238 allocate(stc_updated(land_iau_control%nx * land_iau_control%ny))
239 allocate(slc_updated(land_iau_control%nx * land_iau_control%ny))
240
241 !copy background stc
242 stc_updated = 0
243 slc_updated = 0
244 ib = 1
245 do j = 1, land_iau_control%ny
246 do k = 1, km
247 stc_inc_flat(ib:ib+land_iau_control%nx-1, k) = land_iau_data%stc_inc(:,j, k)
248 slc_inc_flat(ib:ib+land_iau_control%nx-1, k) = land_iau_data%slc_inc(:,j, k)
249 enddo
250 ib = ib + land_iau_control%nx
251 enddo
252
253 if ((land_iau_control%dtp - delt) > 0.0001) then
254 if(land_iau_control%me == land_iau_control%mpi_root) then
255 print*, "Warning! noahmpdrv_timestep_init delt ",delt," different from Land_IAU_Control%dtp ",land_iau_control%dtp
256 endif
257 endif
258
259 lsoil_incr = land_iau_control%lsoil_incr
260 lensfc = land_iau_control%nx * land_iau_control%ny
261
262 if(land_iau_control%me == land_iau_control%mpi_root) print*,' adjusting first ', lsoil_incr, ' surface layers only, delt ', delt
263 ! initialize variables for counts statitics to be zeros
264 nother = 0 ! grid cells not land
265 nsnowupd = 0 ! grid cells with snow (temperature not yet updated)
266 nstcupd = 0 ! grid cells that are updated stc
267 nslcupd = 0 ! grid cells that are updated slc
268 nfrozen = 0 ! not update as frozen soil
269 nfrozen_upd = 0 ! not update as frozen soil
270
271!TODO---if only fv3 increment files are used, this can be read from file
272 allocate(mask_tile(lensfc))
273 call calculate_landinc_mask(weasd, vegtype, soiltyp, lensfc, isice_table, mask_tile)
274
275 !IAU increments are in units of 1/sec !Land_IAU_Control%dtp
276 !* only updating soil temp for now
277 ij_loop : do ij = 1, lensfc
278 ! mask: 1 - soil, 2 - snow, 0 - land-ice, -1 - not land
279 if (mask_tile(ij) == 1) then
280
281 soil_freeze=.false.
282 soil_ice=.false.
283 do k = 1, lsoil_incr ! k = 1, km
284 if ( stc(ij,k) < con_t0c) soil_freeze=.true.
285 if ( smc(ij,k) - slc(ij,k) > 0.001 ) soil_ice=.true.
286
287 if (land_iau_control%upd_stc) then
288 stc(ij,k) = stc(ij,k) + stc_inc_flat(ij,k)*delt !Land_IAU_Control%dtp
289 if (k==1) then
290 stc_updated(ij) = 1
291 nstcupd = nstcupd + 1
292 endif
293 endif
294
295 if ( (stc(ij,k) < con_t0c) .and. (.not. soil_freeze) .and. (k==1) ) nfrozen_upd = nfrozen_upd + 1
296
297 ! do not do updates if this layer or any above is frozen
298 if ( (.not. soil_freeze ) .and. (.not. soil_ice ) ) then
299 if (land_iau_control%upd_slc) then
300 if (k==1) then
301 nslcupd = nslcupd + 1
302 slc_updated(ij) = 1
303 endif
304 ! apply zero limit here (higher, model-specific limits are later)
305 slc(ij,k) = max(slc(ij,k) + slc_inc_flat(ij,k)*delt, 0.0)
306 smc(ij,k) = max(smc(ij,k) + slc_inc_flat(ij,k)*delt, 0.0)
307 endif
308 else
309 if (k==1) nfrozen = nfrozen+1
310 endif
311 enddo
312 endif ! if soil/snow point
313 enddo ij_loop
314
315 !!do moisture/temperature adjustment for consistency after increment add
316 call read_mp_table_parameters(errmsg, errflg)
317 if (errflg .ne. 0) then
318 errmsg = 'FATAL ERROR in noahmpdrv_timestep_init: problem in set_soilveg_noahmp'
319 return
320 endif
321 n_stc = 0
322 n_slc = 0
323 if (land_iau_control%do_stcsmc_adjustment) then
324 if (land_iau_control%upd_stc) then
325 do i=1,lensfc
326 if (stc_updated(i) == 1 ) then ! soil-only location
327 n_stc = n_stc+1
328 soiltype = soiltyp(i)
329 do l = 1, lsoil_incr
330 if (abs(stc_inc_flat(i,l)) > land_iau_control%min_T_increment) then
331 !the following if case applies when updated stc > melting point, it handles both
332 !case 1: frz ==> frz, recalculate slc, smc remains
333 !case 2: unfrz ==> frz, recalculate slc, smc remains
334 if (stc(i,l) .LT. con_t0c )then
335 !recompute supercool liquid water,smc_anl remain unchanged
336 smp = con_hfus*(con_t0c-stc(i,l))/(con_g*stc(i,l)) !(m)
337 slc_new=maxsmc(soiltype)*(smp/satpsi(soiltype))**(-1./bb(soiltype))
338 slc(i,l) = max( min( slc_new, smc(i,l)), 0.0 )
339 endif
340 !case 3: frz ==> unfrz (or unfrz ==> unfrz), melt all soil ice (if any)
341 if (stc(i,l) .GT. con_t0c )then !do not rely on stc_bck
342 slc(i,l)=smc(i,l)
343 endif
344 endif
345 enddo
346 endif
347 enddo
348 endif
349
350 if (land_iau_control%upd_slc) then
351 dz(1) = -zsoil(1)
352 do l = 2, km
353 dz(l) = -zsoil(l) + zsoil(l-1)
354 enddo
355 do i=1,lensfc
356 if (slc_updated(i) == 1 ) then
357 n_slc = n_slc+1
358 ! apply SM bounds (later: add upper SMC limit)
359 do l = 1, lsoil_incr
360 if (abs(slc_inc_flat(i, l)) > land_iau_control%min_SLC_increment) then
361 ! noah-mp minimum is 1 mm per layer (in SMC)
362 ! no need to maintain frozen amount, would be v. small.
363 slc(i,l) = max( 0.001/dz(l), slc(i,l) )
364 smc(i,l) = max( 0.001/dz(l), smc(i,l) )
365 endif
366 enddo
367 endif
368 enddo
369 endif
370 endif
371
372 deallocate(stc_inc_flat, slc_inc_flat, stc_updated, slc_updated)
373 deallocate(mask_tile)
374
375 !Remove non-warning, non-error log write
376 !write(*,'(a,i4,a,i8)') 'noahmpdrv_timestep_init rank ', Land_IAU_Control%me, ' # of cells with stc update ', nstcupd
377
378
379end subroutine noahmpdrv_timestep_init
380
387 subroutine noahmpdrv_finalize (Land_IAU_Control, Land_IAU_Data, Land_IAU_State, errmsg, errflg)
388
389 use machine, only: kind_phys
390 implicit none
391 type(land_iau_control_type) , intent(in ) :: Land_IAU_Control
392 type(land_iau_external_data_type) , intent(inout) :: Land_IAU_Data
393 type(land_iau_state_type) , intent(inout) :: Land_IAU_State
394 character(len=*), intent(out) :: errmsg
395 integer, intent(out) :: errflg
396 integer :: j, k, ib
397 ! --- Initialize CCPP error handling variables
398 errmsg = ''
399 errflg = 0
400
401 if (.not. land_iau_control%do_land_iau) return
402 call land_iau_mod_finalize(land_iau_control, land_iau_data, land_iau_state, errmsg, errflg)
403
404 end subroutine noahmpdrv_finalize
405
427 subroutine noahmpdrv_run &
428!...................................
429! --- inputs:
430 ( im, km, lsnowl, itime, ps, u1, v1, t1, q1, soiltyp,soilcol,&
431 vegtype, sigmaf, dlwflx, dswsfc, snet, delt, tg3, cm, ch, &
432 prsl1, prslk1, prslki, prsik1, zf,pblh, dry, wind, slopetyp,&
433 shdmin, shdmax, snoalb, sfalb, flag_iter,con_g, &
434 idveg, iopt_crs, iopt_btr, iopt_run, iopt_sfc, iopt_frz, &
435 iopt_inf, iopt_rad, iopt_alb, iopt_snf, iopt_tbot,iopt_stc,&
436 iopt_trs,iopt_diag,xlatin, xcoszin, iyrlen, julian, garea, &
437 rainn_mp, rainc_mp, snow_mp, graupel_mp, ice_mp, rhonewsn1,&
438 con_hvap, con_cp, con_jcal, rhoh2o, con_eps, con_epsm1, &
439 con_fvirt, con_rd, con_hfus, thsfc_loc, cpllnd, cpllnd2atm,&
440
441! --- in/outs:
442 weasd, snwdph, tskin, tprcp, srflag, smc, stc, slc, &
443 canopy, trans, tsurf, zorl, &
444 rb1, fm1, fh1, ustar1, stress1, fm101, fh21, &
445 rmol1,flhc1,flqc1,do_mynnsfclay, &
446
447! --- Noah MP specific
448
449 snowxy, tvxy, tgxy, canicexy, canliqxy, eahxy, tahxy, cmxy,&
450 chxy, fwetxy, sneqvoxy, alboldxy, qsnowxy, wslakexy, zwtxy,&
451 waxy, wtxy, tsnoxy, zsnsoxy, snicexy, snliqxy, lfmassxy, &
452 rtmassxy, stmassxy, woodxy, stblcpxy, fastcpxy, xlaixy, &
453 xsaixy, taussxy, smoiseq, smcwtdxy, deeprechxy, rechxy, &
454 albdvis, albdnir, albivis, albinir,emiss, &
455
456! --- outputs:
457 sncovr1, qsurf, gflux, drain, evap, hflx, ep, runoff, &
458 cmm, chh, evbs, evcw, sbsno, pah, ecan, etran, edir, snowc,&
459 stm, snohf,smcwlt2, smcref2, wet1, t2mmp, q2mp,zvfun, &
460 ztmax, rca, errmsg, errflg, &
461 canopy_heat_storage_ccpp, &
462 rainfall_ccpp, &
463 sw_absorbed_total_ccpp, &
464 sw_reflected_total_ccpp, &
465 lw_absorbed_total_ccpp, &
466 temperature_bare_grd_ccpp, &
467 temperature_veg_grd_ccpp, &
468 temperature_veg_2m_ccpp, &
469 temperature_bare_2m_ccpp, &
470 spec_humidity_veg_2m_ccpp, &
471 spec_humidity_bare_2m_ccpp, &
472 sw_absorbed_veg_ccpp, &
473 sw_absorbed_ground_ccpp, &
474 snowmelt_out_ccpp, &
475 snowmelt_shallow_ccpp, &
476 albedo_direct_snow_ccpp, &
477 albedo_diffuse_snow_ccpp, &
478 ch_vegetated_ccpp, &
479 ch_bare_ground_ccpp, &
480 sensible_heat_grd_veg_ccpp, &
481 sensible_heat_leaf_ccpp, &
482 sensible_heat_grd_bar_ccpp, &
483 latent_heat_grd_veg_ccpp, &
484 latent_heat_grd_bare_ccpp, &
485 ground_heat_veg_ccpp, &
486 ground_heat_bare_ccpp, &
487 lw_absorbed_grd_veg_ccpp, &
488 lw_absorbed_leaf_ccpp, &
489 lw_absorbed_grd_bare_ccpp, &
490 latent_heat_trans_ccpp, &
491 latent_heat_leaf_ccpp, &
492 ch_leaf_ccpp, &
493 ch_below_canopy_ccpp, &
494 ch_vegetated_2m_ccpp, &
495 ch_bare_ground_2m_ccpp, &
496 precip_adv_heat_veg_ccpp, &
497 precip_adv_heat_grd_v_ccpp, &
498 precip_adv_heat_grd_b_ccpp, &
499 spec_humid_sfc_veg_ccpp, &
500 spec_humid_sfc_bare_ccpp &
501 )
502
503 use machine , only : kind_phys
504 use funcphys, only : fpvs
505
506 use module_sf_noahmplsm, only : gfs_stability
508 use noahmp_tables
509
510 implicit none
511
512 real(kind=kind_phys), parameter :: a2 = 17.2693882
513 real(kind=kind_phys), parameter :: a3 = 273.16
514 real(kind=kind_phys), parameter :: a4 = 35.86
515 real(kind=kind_phys), parameter :: a23m4 = a2*(a3-a4)
516 real(kind=kind_phys), intent(in) :: con_g
517
518 real, parameter :: undefined = 9.99e20_kind_phys
519
520 integer, parameter :: nsoil = 4 ! hardwired to Noah
521 integer, parameter :: nsnow = 3 ! max. snow layers
522
523 integer, parameter :: iz0tlnd = 0 ! z0t treatment option
524
525 real(kind=kind_phys), save :: zsoil(nsoil)
526 data zsoil / -0.1, -0.4, -1.0, -2.0 /
527
528!
529! --- CCPP interface fields (in call order)
530!
531
532 integer , intent(in) :: im ! horiz dimension and num of used pts
533 integer , intent(in) :: km ! vertical soil layer dimension
534 integer , intent(in) :: lsnowl ! lower bound for snow level arrays
535 integer , intent(in) :: itime ! NOT USED
536 real(kind=kind_phys), dimension(:) , intent(in) :: ps ! surface pressure [Pa]
537 real(kind=kind_phys), dimension(:) , intent(in) :: u1 ! u-component of wind [m/s]
538 real(kind=kind_phys), dimension(:) , intent(in) :: v1 ! u-component of wind [m/s]
539 real(kind=kind_phys), dimension(:) , intent(in) :: t1 ! layer 1 temperature [K]
540 real(kind=kind_phys), dimension(:) , intent(in) :: q1 ! layer 1 specific humidity [kg/kg]
541 integer , dimension(:) , intent(in) :: soiltyp ! soil type (integer index)
542 integer , dimension(:) , intent(in) :: soilcol ! soil color (integer index)
543 integer , dimension(:) , intent(in) :: vegtype ! vegetation type (integer index)
544 real(kind=kind_phys), dimension(:) , intent(in) :: sigmaf ! areal fractional cover of green vegetation
545 real(kind=kind_phys), dimension(:) , intent(in) :: dlwflx ! downward longwave radiation [W/m2]
546 real(kind=kind_phys), dimension(:) , intent(in) :: dswsfc ! downward shortwave radiation [W/m2]
547 real(kind=kind_phys), dimension(:) , intent(in) :: snet ! total sky sfc netsw flx into ground[W/m2]
548 real(kind=kind_phys) , intent(in) :: delt ! time interval [s]
549 real(kind=kind_phys), dimension(:) , intent(in) :: tg3 ! deep soil temperature [K]
550 real(kind=kind_phys), dimension(:) , intent(inout) :: cm ! surface exchange coeff for momentum [-]
551 real(kind=kind_phys), dimension(:) , intent(inout) :: ch ! surface exchange coeff heat & moisture[-]
552 real(kind=kind_phys), dimension(:) , intent(in) :: prsl1 ! sfc layer 1 mean pressure [Pa]
553 real(kind=kind_phys), dimension(:) , intent(in) :: prslk1 ! exner_function_at lowest model layer
554
555 real(kind=kind_phys), dimension(:) , intent(in) :: prslki ! Exner function bt midlayer and interface at 1st layer
556 real(kind=kind_phys), dimension(:) , intent(in) :: prsik1 ! Exner function at the ground surfac
557
558 real(kind=kind_phys), dimension(:) , intent(in) :: zf ! height of bottom layer [m]
559
560 logical , intent(in) :: do_mynnsfclay !flag for MYNN sfc layer scheme
561
562 real(kind=kind_phys), dimension(:) , intent(in) :: pblh ! height of pbl
563 real(kind=kind_phys), dimension(:) , intent(inout) :: rmol1 !
564 real(kind=kind_phys), dimension(:) , intent(inout) :: flhc1 !
565 real(kind=kind_phys), dimension(:) , intent(inout) :: flqc1 !
566
567
568 logical , dimension(:) , intent(in) :: dry ! = T if a point with any land
569 real(kind=kind_phys), dimension(:) , intent(in) :: wind ! wind speed [m/s]
570 integer , dimension(:) , intent(in) :: slopetyp ! surface slope classification
571 real(kind=kind_phys), dimension(:) , intent(in) :: shdmin ! min green vegetation coverage [fraction]
572 real(kind=kind_phys), dimension(:) , intent(in) :: shdmax ! max green vegetation coverage [fraction]
573 real(kind=kind_phys), dimension(:) , intent(in) :: snoalb ! upper bound on max albedo over deep snow
574 real(kind=kind_phys), dimension(:) , intent(inout) :: sfalb ! mean surface albedo [fraction]
575 logical , dimension(:) , intent(in) :: flag_iter !
576 integer , intent(in) :: idveg ! option for dynamic vegetation
577 integer , intent(in) :: iopt_crs ! option for canopy stomatal resistance
578 integer , intent(in) :: iopt_btr ! option for soil moisture factor for stomatal resistance
579 integer , intent(in) :: iopt_run ! option for runoff and groundwater
580 integer , intent(in) :: iopt_sfc ! option for surface layer drag coeff (ch & cm)
581 integer , intent(in) :: iopt_frz ! option for supercooled liquid water (or ice fraction)
582 integer , intent(in) :: iopt_inf ! option for frozen soil permeability
583 integer , intent(in) :: iopt_rad ! option for radiation transfer
584 integer , intent(in) :: iopt_alb ! option for ground snow surface albedo
585 integer , intent(in) :: iopt_snf ! option for partitioning precipitation into rainfall & snowfall
586 integer , intent(in) :: iopt_tbot ! option for lower boundary condition of soil temperature
587 integer , intent(in) :: iopt_stc ! option for snow/soil temperature time scheme (only layer 1)
588 integer , intent(in) :: iopt_trs ! option for thermal roughness scheme
589 integer , intent(in) :: iopt_diag ! option for surface diagnose approach
590 real(kind=kind_phys), dimension(:) , intent(in) :: xlatin ! latitude
591 real(kind=kind_phys), dimension(:) , intent(in) :: xcoszin ! cosine of zenith angle
592 integer , intent(in) :: iyrlen ! year length [days]
593 real(kind=kind_phys) , intent(in) :: julian ! julian day of year
594 real(kind=kind_phys), dimension(:) , intent(in) :: garea ! area of the grid cell
595 real(kind=kind_phys), dimension(:) , intent(in) :: rainn_mp ! microphysics non-convective precipitation [mm]
596 real(kind=kind_phys), dimension(:) , intent(in) :: rainc_mp ! microphysics convective precipitation [mm]
597 real(kind=kind_phys), dimension(:) , intent(in) :: snow_mp ! microphysics snow [mm]
598 real(kind=kind_phys), dimension(:) , intent(in) :: graupel_mp ! microphysics graupel [mm]
599 real(kind=kind_phys), dimension(:) , intent(in) :: ice_mp ! microphysics ice/hail [mm]
600 real(kind=kind_phys), dimension(:) , intent(in) :: rhonewsn1 ! precipitation ice density (kg/m^3)
601 real(kind=kind_phys) , intent(in) :: con_hvap ! latent heat condensation [J/kg]
602 real(kind=kind_phys) , intent(in) :: con_cp ! specific heat air [J/kg/K]
603 real(kind=kind_phys) , intent(in) :: con_jcal ! joules per calorie (not used)
604 real(kind=kind_phys) , intent(in) :: rhoh2o ! density of water [kg/m^3]
605 real(kind=kind_phys) , intent(in) :: con_eps ! Rd/Rv
606 real(kind=kind_phys) , intent(in) :: con_epsm1 ! Rd/Rv - 1
607 real(kind=kind_phys) , intent(in) :: con_fvirt ! Rv/Rd - 1
608 real(kind=kind_phys) , intent(in) :: con_rd ! gas constant air [J/kg/K]
609 real(kind=kind_phys) , intent(in) :: con_hfus ! lat heat H2O fusion [J/kg]
610
611 logical , intent(in) :: thsfc_loc ! Flag for reference pressure in theta calculation
612
613 logical , intent(in) :: cpllnd ! Flag for land coupling (atm->lnd)
614 logical , intent(in) :: cpllnd2atm ! Flag for land coupling (lnd->atm)
615
616 real(kind=kind_phys), dimension(:) , intent(inout) :: weasd ! water equivalent accumulated snow depth [mm]
617 real(kind=kind_phys), dimension(:) , intent(inout) :: snwdph ! snow depth [mm]
618 real(kind=kind_phys), dimension(:) , intent(inout) :: tskin ! ground surface skin temperature [K]
619 real(kind=kind_phys), dimension(:) , intent(inout) :: tprcp ! total precipitation [m]
620 real(kind=kind_phys), dimension(:) , intent(inout) :: srflag ! snow/rain flag for precipitation
621 real(kind=kind_phys), dimension(:,:) , intent(inout) :: smc ! total soil moisture content [m3/m3]
622 real(kind=kind_phys), dimension(:,:) , intent(inout) :: stc ! soil temp [K]
623 real(kind=kind_phys), dimension(:,:) , intent(inout) :: slc ! liquid soil moisture [m3/m3]
624 real(kind=kind_phys), dimension(:) , intent(inout) :: canopy ! canopy moisture content [mm]
625 real(kind=kind_phys), dimension(:) , intent(inout) :: trans ! total plant transpiration [m/s]
626 real(kind=kind_phys), dimension(:) , intent(inout) :: tsurf ! surface skin temperature [K]
627 real(kind=kind_phys), dimension(:) , intent(inout) :: zorl ! surface roughness [cm]
628
629 real(kind=kind_phys), dimension(:) , intent(inout) :: rb1 ! bulk richardson #
630 real(kind=kind_phys), dimension(:) , intent(inout) :: fm1 ! Monin_Obukhov_silarity_function for momentum
631 real(kind=kind_phys), dimension(:) , intent(inout) :: fh1 ! Monin_Obukhov_silarity_function for heat
632 real(kind=kind_phys), dimension(:) , intent(inout) :: ustar1 ! friction velocity m s-1
633 real(kind=kind_phys), dimension(:) , intent(inout) :: stress1 ! Wind stress m2 S-2
634 real(kind=kind_phys), dimension(:) , intent(inout) :: fm101 ! MOS function for momentum evaulated @ 10 m
635 real(kind=kind_phys), dimension(:) , intent(inout) :: fh21 ! MOS function for heat evaulated @ 2m
636
637 real(kind=kind_phys), dimension(:) , intent(inout) :: snowxy ! actual no. of snow layers
638 real(kind=kind_phys), dimension(:) , intent(inout) :: tvxy ! vegetation leaf temperature [K]
639 real(kind=kind_phys), dimension(:) , intent(inout) :: tgxy ! bulk ground surface temperature [K]
640 real(kind=kind_phys), dimension(:) , intent(inout) :: canicexy ! canopy-intercepted ice [mm]
641 real(kind=kind_phys), dimension(:) , intent(inout) :: canliqxy ! canopy-intercepted liquid water [mm]
642 real(kind=kind_phys), dimension(:) , intent(inout) :: eahxy ! canopy air vapor pressure [Pa]
643 real(kind=kind_phys), dimension(:) , intent(inout) :: tahxy ! canopy air temperature [K]
644 real(kind=kind_phys), dimension(:) , intent(inout) :: cmxy ! bulk momentum drag coefficient [m/s]
645 real(kind=kind_phys), dimension(:) , intent(inout) :: chxy ! bulk sensible heat exchange coefficient [m/s]
646 real(kind=kind_phys), dimension(:) , intent(inout) :: fwetxy ! wetted or snowed fraction of the canopy [-]
647 real(kind=kind_phys), dimension(:) , intent(inout) :: sneqvoxy ! snow mass at last time step[mm h2o]
648 real(kind=kind_phys), dimension(:) , intent(inout) :: alboldxy ! snow albedo at last time step [-]
649 real(kind=kind_phys), dimension(:) , intent(inout) :: qsnowxy ! snowfall on the ground [mm/s]
650 real(kind=kind_phys), dimension(:) , intent(inout) :: wslakexy ! lake water storage [mm]
651 real(kind=kind_phys), dimension(:) , intent(inout) :: zwtxy ! water table depth [m]
652 real(kind=kind_phys), dimension(:) , intent(inout) :: waxy ! water in the "aquifer" [mm]
653 real(kind=kind_phys), dimension(:) , intent(inout) :: wtxy ! groundwater storage [mm]
654 real(kind=kind_phys), dimension(:,lsnowl:), intent(inout) :: tsnoxy ! snow temperature [K]
655 real(kind=kind_phys), dimension(:,lsnowl:), intent(inout) :: zsnsoxy ! snow/soil layer depth [m]
656 real(kind=kind_phys), dimension(:,lsnowl:), intent(inout) :: snicexy ! snow layer ice [mm]
657 real(kind=kind_phys), dimension(:,lsnowl:), intent(inout) :: snliqxy ! snow layer liquid water [mm]
658 real(kind=kind_phys), dimension(:) , intent(inout) :: lfmassxy ! leaf mass [g/m2]
659 real(kind=kind_phys), dimension(:) , intent(inout) :: rtmassxy ! mass of fine roots [g/m2]
660 real(kind=kind_phys), dimension(:) , intent(inout) :: stmassxy ! stem mass [g/m2]
661 real(kind=kind_phys), dimension(:) , intent(inout) :: woodxy ! mass of wood (incl. woody roots) [g/m2]
662 real(kind=kind_phys), dimension(:) , intent(inout) :: stblcpxy ! stable carbon in deep soil [g/m2]
663 real(kind=kind_phys), dimension(:) , intent(inout) :: fastcpxy ! short-lived carbon, shallow soil [g/m2]
664 real(kind=kind_phys), dimension(:) , intent(inout) :: xlaixy ! leaf area index [m2/m2]
665 real(kind=kind_phys), dimension(:) , intent(inout) :: xsaixy ! stem area index [m2/m2]
666 real(kind=kind_phys), dimension(:) , intent(inout) :: taussxy ! snow age factor [-]
667 real(kind=kind_phys), dimension(:,:) , intent(inout) :: smoiseq ! eq volumetric soil moisture [m3/m3]
668 real(kind=kind_phys), dimension(:) , intent(inout) :: smcwtdxy ! soil moisture content in the layer to the water table when deep
669 real(kind=kind_phys), dimension(:) , intent(inout) :: deeprechxy ! recharge to the water table when deep
670 real(kind=kind_phys), dimension(:) , intent(inout) :: rechxy ! recharge to the water table
671 real(kind=kind_phys), dimension(:) , intent(out) :: albdvis ! albedo - direct visible [fraction]
672 real(kind=kind_phys), dimension(:) , intent(out) :: albdnir ! albedo - direct NIR [fraction]
673 real(kind=kind_phys), dimension(:) , intent(out) :: albivis ! albedo - diffuse visible [fraction]
674 real(kind=kind_phys), dimension(:) , intent(out) :: albinir ! albedo - diffuse NIR [fraction]
675 real(kind=kind_phys), dimension(:) , intent(out) :: emiss ! sfc lw emissivity [fraction]
676 real(kind=kind_phys), dimension(:) , intent(out) :: sncovr1 ! snow cover over land [fraction]
677 real(kind=kind_phys), dimension(:) , intent(out) :: qsurf ! specific humidity at sfc [kg/kg]
678 real(kind=kind_phys), dimension(:) , intent(out) :: gflux ! soil heat flux [W/m2]
679 real(kind=kind_phys), dimension(:) , intent(out) :: drain ! subsurface runoff [mm/s]
680 real(kind=kind_phys), dimension(:) , intent(out) :: evap ! total latent heat flux [W/m2]
681 real(kind=kind_phys), dimension(:) , intent(out) :: hflx ! sensible heat flux [W/m2]
682 real(kind=kind_phys), dimension(:) , intent(out) :: ep ! potential evaporation [mm/s?]
683 real(kind=kind_phys), dimension(:) , intent(out) :: runoff ! surface runoff [mm/s]
684 real(kind=kind_phys), dimension(:) , intent(out) :: cmm ! cm*U [m/s]
685 real(kind=kind_phys), dimension(:) , intent(out) :: chh ! ch*U*rho [kg/m2/s]
686 real(kind=kind_phys), dimension(:) , intent(out) :: evbs ! direct soil evaporation [m/s]
687 real(kind=kind_phys), dimension(:) , intent(out) :: evcw ! canopy water evaporation [m/s]
688 real(kind=kind_phys), dimension(:) , intent(out) :: sbsno ! sublimation/deposit from snopack [W/m2]
689 real(kind=kind_phys), dimension(:) , intent(out) :: pah ! precipitation advected heat - total (w/m2)
690 real(kind=kind_phys), dimension(:) , intent(out) :: ecan ! evaporation of intercepted water (mm/s)
691 real(kind=kind_phys), dimension(:) , intent(out) :: etran ! transpiration rate (mm/s)
692 real(kind=kind_phys), dimension(:) , intent(out) :: edir ! soil surface evaporation rate (mm/s)
693 real(kind=kind_phys), dimension(:) , intent(out) :: snowc ! fractional snow cover [-]
694 real(kind=kind_phys), dimension(:) , intent(out) :: stm ! total soil column moisture content [mm]
695 real(kind=kind_phys), dimension(:) , intent(out) :: snohf ! snow/freezing-rain latent heat flux [W/m2]
696 real(kind=kind_phys), dimension(:) , intent(out) :: smcwlt2 ! dry soil moisture threshold [m3/m3]
697 real(kind=kind_phys), dimension(:) , intent(out) :: smcref2 ! soil moisture threshold [m3/m3]
698 real(kind=kind_phys), dimension(:) , intent(out) :: wet1 ! normalized surface soil saturated fraction
699 real(kind=kind_phys), dimension(:) , intent(out) :: t2mmp ! combined T2m from tiles
700 real(kind=kind_phys), dimension(:) , intent(out) :: q2mp ! combined q2m from tiles
701 real(kind=kind_phys), dimension(:) , intent(out) :: zvfun !
702 real(kind=kind_phys), dimension(:) , intent(out) :: ztmax ! thermal roughness length
703 real(kind=kind_phys), dimension(:) , intent(out) :: rca ! total canopy/stomatal resistance (s/m)
704
705 character(len=*) , intent(out) :: errmsg
706 integer , intent(out) :: errflg
707
708 real(kind=kind_phys), dimension(:) , intent(out), optional :: canopy_heat_storage_ccpp ! within-canopy heat [W/m2]
709 real(kind=kind_phys), dimension(:) , intent(out), optional :: rainfall_ccpp
710 real(kind=kind_phys), dimension(:) , intent(out), optional :: sw_absorbed_total_ccpp
711 real(kind=kind_phys), dimension(:) , intent(out), optional :: sw_reflected_total_ccpp
712 real(kind=kind_phys), dimension(:) , intent(out), optional :: lw_absorbed_total_ccpp
713 real(kind=kind_phys), dimension(:) , intent(out), optional :: temperature_bare_grd_ccpp
714 real(kind=kind_phys), dimension(:) , intent(out), optional :: temperature_veg_grd_ccpp
715 real(kind=kind_phys), dimension(:) , intent(out), optional :: temperature_veg_2m_ccpp
716 real(kind=kind_phys), dimension(:) , intent(out), optional :: temperature_bare_2m_ccpp
717 real(kind=kind_phys), dimension(:) , intent(out), optional :: spec_humidity_veg_2m_ccpp
718 real(kind=kind_phys), dimension(:) , intent(out), optional :: spec_humidity_bare_2m_ccpp
719 real(kind=kind_phys), dimension(:) , intent(out), optional :: sw_absorbed_veg_ccpp
720 real(kind=kind_phys), dimension(:) , intent(out), optional :: sw_absorbed_ground_ccpp
721 real(kind=kind_phys), dimension(:) , intent(out), optional :: snowmelt_out_ccpp
722 real(kind=kind_phys), dimension(:) , intent(out), optional :: snowmelt_shallow_ccpp
723 real(kind=kind_phys), dimension(:,:), intent(out), optional :: albedo_direct_snow_ccpp
724 real(kind=kind_phys), dimension(:,:), intent(out), optional :: albedo_diffuse_snow_ccpp
725 real(kind=kind_phys), dimension(:) , intent(out), optional :: ch_vegetated_ccpp
726 real(kind=kind_phys), dimension(:) , intent(out), optional :: ch_bare_ground_ccpp
727 real(kind=kind_phys), dimension(:) , intent(out), optional :: sensible_heat_grd_veg_ccpp
728 real(kind=kind_phys), dimension(:) , intent(out), optional :: sensible_heat_leaf_ccpp
729 real(kind=kind_phys), dimension(:) , intent(out), optional :: sensible_heat_grd_bar_ccpp
730 real(kind=kind_phys), dimension(:) , intent(out), optional :: latent_heat_grd_veg_ccpp
731 real(kind=kind_phys), dimension(:) , intent(out), optional :: latent_heat_grd_bare_ccpp
732 real(kind=kind_phys), dimension(:) , intent(out), optional :: ground_heat_veg_ccpp
733 real(kind=kind_phys), dimension(:) , intent(out), optional :: ground_heat_bare_ccpp
734 real(kind=kind_phys), dimension(:) , intent(out), optional :: lw_absorbed_grd_veg_ccpp
735 real(kind=kind_phys), dimension(:) , intent(out), optional :: lw_absorbed_leaf_ccpp
736 real(kind=kind_phys), dimension(:) , intent(out), optional :: lw_absorbed_grd_bare_ccpp
737 real(kind=kind_phys), dimension(:) , intent(out), optional :: latent_heat_trans_ccpp
738 real(kind=kind_phys), dimension(:) , intent(out), optional :: latent_heat_leaf_ccpp
739 real(kind=kind_phys), dimension(:) , intent(out), optional :: ch_leaf_ccpp
740 real(kind=kind_phys), dimension(:) , intent(out), optional :: ch_below_canopy_ccpp
741 real(kind=kind_phys), dimension(:) , intent(out), optional :: ch_vegetated_2m_ccpp
742 real(kind=kind_phys), dimension(:) , intent(out), optional :: ch_bare_ground_2m_ccpp
743 real(kind=kind_phys), dimension(:) , intent(out), optional :: precip_adv_heat_veg_ccpp
744 real(kind=kind_phys), dimension(:) , intent(out), optional :: precip_adv_heat_grd_v_ccpp
745 real(kind=kind_phys), dimension(:) , intent(out), optional :: precip_adv_heat_grd_b_ccpp
746 real(kind=kind_phys), dimension(:) , intent(out), optional :: spec_humid_sfc_veg_ccpp
747 real(kind=kind_phys), dimension(:) , intent(out), optional :: spec_humid_sfc_bare_ccpp
748
749!
750! --- some new options, hard code for now
751!
752
753 integer :: iopt_rsf = 4 ! option for surface resistance
754 integer :: iopt_soil = 1 ! option for soil parameter treatment
755 integer :: iopt_pedo = 1 ! option for pedotransfer function
756 integer :: iopt_crop = 0 ! option for crop model
757 integer :: iopt_gla = 2 ! option for glacier treatment
758 integer :: iopt_z0m = 1 ! option for z0m treatment
759
760!
761! --- local inputs to noah-mp and glacier subroutines; listed in order in noah-mp call
762!
763 ! intent
764 integer :: i_location ! in | grid index
765 integer :: j_location ! in | grid index (not used in ccpp)
766 real (kind=kind_phys) :: latitude ! in | latitude [radians]
767 integer :: year_length ! in | number of days in the current year
768 real (kind=kind_phys) :: julian_day ! in | julian day of year [floating point]
769 real (kind=kind_phys) :: cosine_zenith ! in | cosine solar zenith angle [-1,1]
770 real (kind=kind_phys) :: timestep ! in | time step [sec]
771 real (kind=kind_phys) :: spatial_scale ! in | spatial scale [m] (not used in noah-mp)
772 real (kind=kind_phys) :: atmosphere_thickness ! in | thickness of lowest atmo layer [m] (not used in noah-mp)
773 integer :: soil_levels ! in | soil levels
774 real (kind=kind_phys), dimension( 1:nsoil) :: soil_interface_depth ! in | soil layer-bottom depth from surface [m]
775 integer :: max_snow_levels ! in | maximum number of snow levels
776 real (kind=kind_phys) :: vegetation_frac ! in | vegetation fraction [0.0-1.0]
777 real (kind=kind_phys) :: area_grid ! in |
778 real (kind=kind_phys) :: max_vegetation_frac ! in | annual maximum vegetation fraction [0.0-1.0]
779 integer :: vegetation_category ! in | vegetation category
780 integer :: ice_flag ! in | ice flag (1->ice)
781 integer :: surface_type ! in | surface type flag 1->soil; 2->lake
782 integer :: crop_type ! in | crop type category
783 real (kind=kind_phys), dimension( 1:nsoil) :: eq_soil_water_vol ! in | (opt_run=5) equilibrium soil water content [m3/m3]
784 real (kind=kind_phys) :: temperature_forcing ! in | forcing air temperature [K]
785 real (kind=kind_phys) :: air_pressure_surface ! in | surface air pressure [Pa]
786 real (kind=kind_phys) :: air_pressure_forcing ! in | forcing air pressure [Pa]
787 real (kind=kind_phys) :: uwind_forcing ! in | forcing u-wind [m/s]
788 real (kind=kind_phys) :: vwind_forcing ! in | forcing v-wind [m/s]
789 real (kind=kind_phys) :: spec_humidity_forcing ! in | forcing mixing ratio [kg/kg]
790 real (kind=kind_phys) :: cloud_water_forcing ! in | cloud water mixing ratio [kg/kg] (not used in noah-mp)
791 real (kind=kind_phys) :: sw_radiation_forcing ! in | forcing downward shortwave radiation [W/m2]
792 real (kind=kind_phys) :: radiation_lw_forcing ! in | forcing downward longwave radiation [W/m2]
793 real (kind=kind_phys) :: precipitation_forcing ! in | total precipitation [mm/s]
794 real (kind=kind_phys) :: precip_convective ! in | convective precipitation [mm/s]
795 real (kind=kind_phys) :: precip_non_convective ! in | non-convective precipitation [mm/s]
796 real (kind=kind_phys) :: precip_sh_convective ! in | shallow convective precipitation [mm/s]
797 real (kind=kind_phys) :: precip_snow ! in | snow precipitation [mm/s]
798 real (kind=kind_phys) :: precip_graupel ! in | graupel precipitation [mm/s]
799 real (kind=kind_phys) :: precip_hail ! in | hail precipitation [mm/s]
800 real (kind=kind_phys) :: temperature_soil_bot ! in | soil bottom boundary condition temperature [K]
801 real (kind=kind_phys) :: co2_air ! in | atmospheric co2 concentration [Pa]
802 real (kind=kind_phys) :: o2_air ! in | atmospheric o2 concentration [Pa]
803 real (kind=kind_phys) :: foliage_nitrogen ! in | foliage nitrogen [%] [1-saturated]
804 real (kind=kind_phys), dimension(-nsnow+1: 0) :: snow_ice_frac_old ! in | snow ice fraction at last timestep [-]
805 real (kind=kind_phys) :: forcing_height ! inout | forcing height [m]
806 real (kind=kind_phys) :: snow_albedo_old ! inout | snow albedo at last time step (class option) [-]
807 real (kind=kind_phys) :: snow_water_equiv_old ! inout | snow water equivalent at last time step [mm]
808 real (kind=kind_phys), dimension(-nsnow+1:nsoil) :: temperature_snow_soil ! inout | snow/soil temperature [K]
809 real (kind=kind_phys), dimension( 1:nsoil) :: soil_liquid_vol ! inout | volumetric liquid soil moisture [m3/m3]
810 real (kind=kind_phys), dimension( 1:nsoil) :: soil_moisture_vol ! inout | volumetric soil moisture (ice + liq.) [m3/m3]
811
812 real (kind=kind_phys) :: surface_temperature ! out | surface aerodynamic temp
813
814 real (kind=kind_phys) :: temperature_canopy_air! inout | canopy air tmeperature [K]
815 real (kind=kind_phys) :: vapor_pres_canopy_air ! inout | canopy air vapor pressure [Pa]
816 real (kind=kind_phys) :: canopy_wet_fraction ! inout | wetted or snowed fraction of canopy (-)
817 real (kind=kind_phys) :: canopy_liquid ! inout | canopy intercepted liquid [mm]
818 real (kind=kind_phys) :: canopy_ice ! inout | canopy intercepted ice [mm]
819 real (kind=kind_phys) :: temperature_leaf ! inout | leaf temperature [K]
820 real (kind=kind_phys) :: temperature_ground ! inout | grid ground surface temperature [K]
821 real (kind=kind_phys) :: spec_humidity_surface ! inout | surface specific humidty [kg/kg]
822 real (kind=kind_phys) :: snowfall ! inout | land model partitioned snowfall [mm/s]
823 real (kind=kind_phys) :: rainfall ! inout | land model partitioned rainfall [mm/s]
824 integer :: snow_levels ! inout | active snow levels [-]
825 real (kind=kind_phys), dimension(-nsnow+1:nsoil) :: interface_depth ! inout | layer-bottom depth from snow surf [m]
826 real (kind=kind_phys) :: snow_depth ! inout | snow depth [m]
827 real (kind=kind_phys) :: snow_water_equiv ! inout | snow water equivalent [mm]
828 real (kind=kind_phys), dimension(-nsnow+1: 0) :: snow_level_ice ! inout | snow level ice [mm]
829 real (kind=kind_phys), dimension(-nsnow+1: 0) :: snow_level_liquid ! inout | snow level liquid [mm]
830 real (kind=kind_phys) :: depth_water_table ! inout | depth to water table [m]
831 real (kind=kind_phys) :: aquifer_water ! inout | water storage in aquifer [mm]
832 real (kind=kind_phys) :: saturated_water ! inout | water in aquifer+saturated soil [mm]
833 real (kind=kind_phys) :: lake_water ! inout | lake water storage (can be neg.) [mm]
834 real (kind=kind_phys) :: leaf_carbon ! inout | leaf mass [g/m2]
835 real (kind=kind_phys) :: root_carbon ! inout | mass of fine roots [g/m2]
836 real (kind=kind_phys) :: stem_carbon ! inout | stem mass [g/m2]
837 real (kind=kind_phys) :: wood_carbon ! inout | mass of wood (incl. woody roots) [g/m2]
838 real (kind=kind_phys) :: soil_carbon_stable ! inout | stable soil carbon [g/m2]
839 real (kind=kind_phys) :: soil_carbon_fast ! inout | short-lived soil carbon [g/m2]
840 real (kind=kind_phys) :: leaf_area_index ! inout | leaf area index [-]
841 real (kind=kind_phys) :: stem_area_index ! inout | stem area index [-]
842 real (kind=kind_phys) :: cm_noahmp ! inout | grid momentum drag coefficient [m/s]
843 real (kind=kind_phys) :: ch_noahmp ! inout | grid heat exchange coefficient [m/s]
844 real (kind=kind_phys) :: snow_age ! inout | non-dimensional snow age [-]
845 real (kind=kind_phys) :: grain_carbon ! inout | grain mass [g/m2]
846 real (kind=kind_phys) :: growing_deg_days ! inout | growing degree days [-]
847 integer :: plant_growth_stage ! inout | plant growing stage [-]
848 real (kind=kind_phys) :: soil_moisture_wtd ! inout | (opt_run=5) soil water content between bottom of the soil and water table [m3/m3]
849 real (kind=kind_phys) :: deep_recharge ! inout | (opt_run=5) recharge to or from the water table when deep [m]
850 real (kind=kind_phys) :: recharge ! inout | (opt_run=5) recharge to or from the water table when shallow [m] (diagnostic)
851 real (kind=kind_phys) :: z0_total ! out | weighted z0 sent to coupled model [m]
852 real (kind=kind_phys) :: z0h_total ! out | weighted z0h sent to coupled model [m]
853
854 real (kind=kind_phys) :: sw_absorbed_total ! out | total absorbed solar radiation [W/m2]
855 real (kind=kind_phys) :: sw_reflected_total ! out | total reflected solar radiation [W/m2]
856 real (kind=kind_phys) :: lw_absorbed_total ! out | total net lw rad [W/m2] [+ to atm]
857 real (kind=kind_phys) :: sensible_heat_total ! out | total sensible heat [W/m2] [+ to atm]
858 real (kind=kind_phys) :: ground_heat_total ! out | ground heat flux [W/m2] [+ to soil]
859 real (kind=kind_phys) :: latent_heat_canopy ! out | canopy evaporation heat flux [W/m2] [+ to atm]
860 real (kind=kind_phys) :: latent_heat_ground ! out | ground evaporation heat flux [W/m2] [+ to atm]
861 real (kind=kind_phys) :: transpiration_heat ! out | transpiration heat flux [W/m2] [+ to atm]
862 real (kind=kind_phys) :: evaporation_canopy ! out | canopy evaporation [mm/s]
863 real (kind=kind_phys) :: transpiration ! out | transpiration [mm/s]
864 real (kind=kind_phys) :: evaporation_soil ! out | soil surface evaporation [mm/s]
865 real (kind=kind_phys) :: temperature_radiative ! out | surface radiative temperature [K]
866 real (kind=kind_phys) :: temperature_bare_grd ! out | bare ground surface temperature [K]
867 real (kind=kind_phys) :: temperature_veg_grd ! out | below_canopy ground surface temperature [K]
868 real (kind=kind_phys) :: temperature_veg_2m ! out | vegetated 2-m air temperature [K]
869 real (kind=kind_phys) :: temperature_bare_2m ! out | bare ground 2-m air temperature [K]
870 real (kind=kind_phys) :: spec_humidity_veg_2m ! out | vegetated 2-m air specific humidity [K]
871 real (kind=kind_phys) :: spec_humidity_bare_2m ! out | bare ground 2-m air specfic humidity [K]
872 real (kind=kind_phys) :: runoff_surface ! out | surface runoff [mm/s]
873 real (kind=kind_phys) :: runoff_baseflow ! out | baseflow runoff [mm/s]
874 real (kind=kind_phys) :: par_absorbed ! out | absorbed photosynthesis active radiation [W/m2]
875 real (kind=kind_phys) :: photosynthesis ! out | total photosynthesis [umol CO2/m2/s] [+ out]
876 real (kind=kind_phys) :: sw_absorbed_veg ! out | solar radiation absorbed by vegetation [W/m2]
877 real (kind=kind_phys) :: sw_absorbed_ground ! out | solar radiation absorbed by ground [W/m2]
878 real (kind=kind_phys) :: snow_cover_fraction ! out | snow cover fraction on the ground [-]
879 real (kind=kind_phys) :: net_eco_exchange ! out | net ecosystem exchange [g/m2/s CO2]
880 real (kind=kind_phys) :: global_prim_prod ! out | global primary production [g/m2/s C]
881 real (kind=kind_phys) :: net_prim_prod ! out | net primary productivity [g/m2/s C]
882 real (kind=kind_phys) :: vegetation_fraction ! out | vegetation fraction [0.0-1.0]
883 real (kind=kind_phys) :: albedo_total ! out | total surface albedo [-]
884 real (kind=kind_phys) :: snowmelt_out ! out | snowmelt out bottom of pack [mm/s]
885 real (kind=kind_phys) :: snowmelt_shallow ! out | shallow snow melt [mm]
886 real (kind=kind_phys) :: snowmelt_shallow_1 ! out | additional shallow snow melt [mm]
887 real (kind=kind_phys) :: snowmelt_shallow_2 ! out | additional shallow snow melt [mm]
888 real (kind=kind_phys) :: rs_sunlit ! out | sunlit leaf stomatal resistance [s/m]
889 real (kind=kind_phys) :: rs_shaded ! out | shaded leaf stomatal resistance [s/m]
890 real (kind=kind_phys), dimension(1:2) :: albedo_direct ! out | direct vis/nir albedo [-]
891 real (kind=kind_phys), dimension(1:2) :: albedo_diffuse ! out | diffuse vis/nir albedo [-]
892 real (kind=kind_phys), dimension(1:2) :: albedo_direct_snow ! out | direct vis/nir snow albedo [-]
893 real (kind=kind_phys), dimension(1:2) :: albedo_diffuse_snow ! out | diffuse vis/nir snow albedo [-]
894 real (kind=kind_phys) :: canopy_gap_fraction ! out | between canopy gap fraction [-]
895 real (kind=kind_phys) :: incanopy_gap_fraction ! out | within canopy gap fraction for beam [-]
896 real (kind=kind_phys) :: ch_vegetated ! out | vegetated heat exchange coefficient [m/s]
897 real (kind=kind_phys) :: ch_bare_ground ! out | bare-ground heat exchange coefficient [m/s]
898 real (kind=kind_phys) :: emissivity_total ! out | grid emissivity [-]
899 real (kind=kind_phys) :: sensible_heat_grd_veg ! out | below-canopy ground sensible heat flux [W/m2]
900 real (kind=kind_phys) :: sensible_heat_leaf ! out | leaf-to-canopy sensible heat flux [W/m2]
901 real (kind=kind_phys) :: sensible_heat_grd_bar ! out | bare ground sensible heat flux [W/m2]
902 real (kind=kind_phys) :: latent_heat_grd_veg ! out | below-canopy ground evaporation heat flux [W/m2]
903 real (kind=kind_phys) :: latent_heat_grd_bare ! out | bare ground evaporation heat flux [W/m2]
904 real (kind=kind_phys) :: ground_heat_veg ! out | below-canopy ground heat flux [W/m2]
905 real (kind=kind_phys) :: ground_heat_bare ! out | bare ground heat flux [W/m2]
906 real (kind=kind_phys) :: lw_absorbed_grd_veg ! out | below-canopy ground absorbed longwave radiation [W/m2]
907 real (kind=kind_phys) :: lw_absorbed_leaf ! out | leaf absorbed longwave radiation [W/m2]
908 real (kind=kind_phys) :: lw_absorbed_grd_bare ! out | bare ground net longwave radiation [W/m2]
909 real (kind=kind_phys) :: latent_heat_trans ! out | transpiration [W/m2]
910 real (kind=kind_phys) :: latent_heat_leaf ! out | leaf evaporation [W/m2]
911 real (kind=kind_phys) :: ch_leaf ! out | leaf exchange coefficient [m/s]
912 real (kind=kind_phys) :: ch_below_canopy ! out | below-canopy exchange coefficient [m/s]
913 real (kind=kind_phys) :: ch_vegetated_2m ! out | 2-m vegetated heat exchange coefficient [m/s]
914 real (kind=kind_phys) :: ch_bare_ground_2m ! out | 2-m bare-ground heat exchange coefficient [m/s]
915 real (kind=kind_phys) :: precip_frozen_frac ! out | precipitation snow fraction [-]
916 real (kind=kind_phys) :: precip_adv_heat_veg ! out | precipitation advected heat - vegetation net [W/m2]
917 real (kind=kind_phys) :: precip_adv_heat_grd_v ! out | precipitation advected heat - below-canopy net [W/m2]
918 real (kind=kind_phys) :: precip_adv_heat_grd_b ! out | precipitation advected heat - bare ground net [W/m2]
919 real (kind=kind_phys) :: precip_adv_heat_total ! out | precipitation advected heat - total [W/m2)
920 real (kind=kind_phys) :: snow_sublimation ! out | snow sublimation [W/m2]
921 real (kind=kind_phys) :: lai_sunlit ! out | sunlit leaf area index [m2/m2]
922 real (kind=kind_phys) :: lai_shaded ! out | shaded leaf area index [m2/m2]
923 real (kind=kind_phys) :: leaf_air_resistance ! out | leaf boundary layer resistance [s/m]
924
925 real (kind=kind_phys) :: canopy_heat_storage ! out | within-canopy heat [W/m2]
926 real (kind=kind_phys) :: spec_humid_sfc_veg ! out | surface specific humidty over vegetation [kg/kg]
927 real (kind=kind_phys) :: spec_humid_sfc_bare ! out | surface specific humidty over bare soil [kg/kg]
928
929 real (kind=kind_phys) :: ustarx ! inout |surface friction velocity
930 real (kind=kind_phys) :: prslkix ! in exner function
931 real (kind=kind_phys) :: prsik1x ! in exner function
932 real (kind=kind_phys) :: prslk1x ! in exner function
933
934 real (kind=kind_phys) :: ch2
935 real (kind=kind_phys) :: cq2
936 real (kind=kind_phys) :: qfx
937 real (kind=kind_phys) :: wspd1 ! wind speed with all components
938 real (kind=kind_phys) :: pblhx ! height of pbl
939 integer :: mnice
940
941 real (kind=kind_phys) :: rah_total !
942 real (kind=kind_phys) :: cah_total !
943
944
945!
946! --- local variable
947!
948
949 integer :: soil_category(nsoil)
950 integer :: slope_category
951 integer :: soil_color_category
952 character(len=256) :: dataset_identifier
953
954 real (kind=kind_phys) :: spec_humidity_sat ! saturation specific humidity
955 real (kind=kind_phys) :: vapor_pressure_sat ! saturation vapor pressure
956 real (kind=kind_phys) :: latent_heat_total ! total latent heat flux [W/m2]
957 real (kind=kind_phys) :: density ! air density
958 real (kind=kind_phys) :: virtual_temperature ! used for penman calculation and density
959 real (kind=kind_phys) :: potential_evaporation ! used for penman calculation
960 real (kind=kind_phys) :: potential_temperature ! used for penman calculation
961 real (kind=kind_phys) :: penman_radiation ! used for penman calculation
962 real (kind=kind_phys) :: dqsdt ! used for penman calculation
963 real (kind=kind_phys) :: precip_freeze_frac_in ! used for penman calculation
964
965 real (kind=kind_phys) :: virtfac1 ! virtual factor
966 real (kind=kind_phys) :: tflux ! surface flux temp
967 real (kind=kind_phys) :: tvs1 ! surface virtual temp
968 real (kind=kind_phys) :: vptemp ! virtual potential temp
969
970 real(kind=kind_phys) :: tem1,tem2,gdx
971 real(kind=kind_phys), parameter :: z0lo=0.1, z0up=1.0
972
973 logical :: is_snowing ! used for penman calculation
974 logical :: is_freeze_rain ! used for penman calculation
975 integer :: i, k
976
977!
978! --- local derived constants:
979!
980
981 type(noahmp_parameters) :: parameters
982
983!
984! --- end declaration
985!
986
987!
988! --- Initialize CCPP error handling variables
989!
990 errmsg = ''
991 errflg = 0
992
993!
994! --- Just return if external land component is activated for two-way interaction
995!
996 if (cpllnd .and. cpllnd2atm) return
997
998 do i = 1, im
999
1000 if (flag_iter(i) .and. dry(i)) then
1001
1002!
1003! --- variable checks and derived fields
1004!
1005
1006 if(vegtype(i) == isice_table ) then
1007 if(weasd(i) < 0.1) then
1008 weasd(i) = 0.1
1009 end if
1010 end if
1011
1012!
1013! --- noah-mp input variables (except snow_ice_frac_old done later)
1014!
1015
1016 dataset_identifier = "modified_igbp_modis_noah"
1017
1018 i_location = i
1019 j_location = -9999
1020 latitude = xlatin(i)
1021 year_length = iyrlen
1022 julian_day = julian
1023 cosine_zenith = xcoszin(i)
1024 timestep = delt
1025 spatial_scale = -9999.0
1026 atmosphere_thickness = -9999.0
1027 soil_levels = km
1028 soil_interface_depth = zsoil
1029 max_snow_levels = nsnow
1030 vegetation_frac = sigmaf(i)
1031 max_vegetation_frac = shdmax(i)
1032 vegetation_category = vegtype(i)
1033 surface_type = 1
1034 crop_type = 0
1035 eq_soil_water_vol = smoiseq(i,:) ! only need for run=5
1036 temperature_forcing = t1(i)
1037 air_pressure_surface = ps(i)
1038 air_pressure_forcing = prsl1(i)
1039 uwind_forcing = u1(i)
1040 vwind_forcing = v1(i)
1041 area_grid = garea(i)
1042
1043 pblhx = pblh(i)
1044
1045 prslkix = prslki(i)
1046 prsik1x = prsik1(i)
1047 prslk1x = prslk1(i)
1048
1049 spec_humidity_forcing = max(q1(i), 1.e-8) ! specific humidity at level 1 (kg/kg)
1050 virtual_temperature = temperature_forcing * &
1051 (1.0 + con_fvirt * spec_humidity_forcing) ! virtual temperature
1052 vapor_pressure_sat = fpvs( temperature_forcing ) ! sat. vapor pressure at level 1 (Pa)
1053 spec_humidity_sat = con_eps*vapor_pressure_sat / &
1054 (prsl1(i) + con_epsm1*vapor_pressure_sat) ! sat. specific humidity at level 1 (kg/kg)
1055 spec_humidity_sat = max(spec_humidity_sat, 1.e-8) ! lower limit sat. specific humidity (kg/kg)
1056 spec_humidity_forcing = min(spec_humidity_sat,spec_humidity_forcing) ! limit specific humidity at level 1 (kg/kg)
1057
1058 cloud_water_forcing = -9999.0
1059 sw_radiation_forcing = dswsfc(i)
1060 radiation_lw_forcing = dlwflx(i)
1061 precipitation_forcing = 1000.0 * tprcp(i) / delt
1062 precip_convective = rainc_mp(i)
1063 precip_non_convective = rainn_mp(i)
1064 precip_sh_convective = 0.
1065 precip_snow = snow_mp(i)
1066 precip_graupel = graupel_mp(i)
1067 precip_hail = ice_mp(i)
1068 temperature_soil_bot = tg3(i)
1069 co2_air = co2_table * air_pressure_forcing
1070 o2_air = o2_table * air_pressure_forcing
1071 foliage_nitrogen = 1.0
1072
1073!
1074! --- noah-mp inout variables
1075!
1076
1077 forcing_height = zf(i)
1078 snow_albedo_old = alboldxy(i)
1079 snow_water_equiv_old = sneqvoxy(i)
1080 temperature_snow_soil(-2: 0) = tsnoxy(i,:)
1081 temperature_snow_soil( 1:km) = stc(i,:)
1082 soil_liquid_vol = slc(i,:)
1083 soil_moisture_vol = smc(i,:)
1084 temperature_canopy_air = tahxy(i)
1085 vapor_pres_canopy_air = eahxy(i)
1086 canopy_wet_fraction = fwetxy(i)
1087 canopy_liquid = canliqxy(i)
1088 canopy_ice = canicexy(i)
1089 temperature_leaf = tvxy(i)
1090 temperature_ground = tgxy(i)
1091 spec_humidity_surface = undefined ! doesn't need inout; should be out
1092 snowfall = qsnowxy(i) ! doesn't need inout; should be out
1093 rainfall = -9999.0 ! doesn't need inout; should be out
1094 snow_levels = nint(snowxy(i))
1095 interface_depth = zsnsoxy(i,:)
1096 snow_depth = snwdph(i) * 0.001 ! convert from mm to m
1097 snow_water_equiv = weasd(i)
1098 if (snow_water_equiv /= 0.0 .and. snow_depth == 0.0) then
1099 snow_depth = 10.0 * snow_water_equiv /1000.0
1100 endif
1101 snow_level_ice = snicexy(i,:)
1102 snow_level_liquid = snliqxy(i,:)
1103 depth_water_table = zwtxy(i)
1104 aquifer_water = waxy(i)
1105 saturated_water = wtxy(i)
1106 lake_water = wslakexy(i)
1107 leaf_carbon = lfmassxy(i)
1108 root_carbon = rtmassxy(i)
1109 stem_carbon = stmassxy(i)
1110 wood_carbon = woodxy(i)
1111 soil_carbon_stable = stblcpxy(i)
1112 soil_carbon_fast = fastcpxy(i)
1113 leaf_area_index = xlaixy(i)
1114 stem_area_index = xsaixy(i)
1115 cm_noahmp = cmxy(i)
1116 ch_noahmp = chxy(i)
1117 snow_age = taussxy(i)
1118! grain_carbon ! new variable
1119! growing_deg_days ! new variable
1120! plant_growth_stage ! new variable
1121 soil_moisture_wtd = smcwtdxy(i)
1122 deep_recharge = deeprechxy(i)
1123 recharge = rechxy(i)
1124
1125 ustarx = ustar1(i)
1126
1127 snow_ice_frac_old = 0.0
1128 do k = snow_levels+1, 0
1129 if(snow_level_ice(k) > 0.0 ) &
1130 snow_ice_frac_old(k) = snow_level_ice(k) /(snow_level_ice(k)+snow_level_liquid(k))
1131 end do
1132
1133
1134 if (snow_depth .gt. 0.1 .or. vegetation_category == isice_table ) then
1135 mnice = 1
1136 else
1137 mnice = 0
1138 endif
1139
1140!
1141! --- some outputs for atm model?
1142!
1143 density = air_pressure_forcing / (con_rd * virtual_temperature)
1144 chh(i) = ch(i) * wind(i) * density
1145 cmm(i) = cm(i) * wind(i)
1146!
1147! --- noah-mp additional variables
1148!
1149
1150 soil_category = soiltyp(i)
1151 slope_category = slopetyp(i)
1152 soil_color_category = soilcol(i)
1153! soil_color_category = 4
1154
1155
1156 call transfer_mp_parameters(vegetation_category, soil_category, &
1157 slope_category, soil_color_category, crop_type,parameters)
1158 parameters%prcpiceden = rhonewsn1(i)
1159 call noahmp_options(idveg ,iopt_crs, iopt_btr , iopt_run, iopt_sfc, &
1160 iopt_frz, iopt_inf , iopt_rad, iopt_alb, &
1161 iopt_snf, iopt_tbot, iopt_stc, iopt_rsf, &
1162 iopt_soil,iopt_pedo, iopt_crop,iopt_trs, &
1163 iopt_diag,iopt_z0m)
1164
1165 if ( vegetation_category == isice_table ) then
1166
1167 if (precipitation_forcing > 0.0) then
1168 if (srflag(i) > 0.0) then
1169 snowfall = srflag(i) * precipitation_forcing ! need snowfall for glacier snow age
1170 endif
1171 endif
1172
1173 ice_flag = -1
1174 temperature_soil_bot = min(temperature_soil_bot,263.15)
1175
1176 call noahmp_options_glacier(iopt_alb, iopt_snf, iopt_tbot, iopt_stc, iopt_gla, &
1177 iopt_sfc ,iopt_trs)
1178 vegetation_frac = 0.0
1179 call noahmp_glacier ( &
1180 i_location ,1 ,cosine_zenith ,nsnow , &
1181 nsoil ,timestep , &
1182 temperature_forcing ,air_pressure_forcing ,uwind_forcing ,vwind_forcing , &
1183 spec_humidity_forcing,sw_radiation_forcing ,precipitation_forcing,radiation_lw_forcing , &
1184 temperature_soil_bot ,forcing_height ,snow_ice_frac_old ,zsoil , &
1185 thsfc_loc ,prslkix ,prsik1x ,prslk1x , &
1186 air_pressure_surface ,pblhx ,iz0tlnd ,itime , &
1187 vegetation_frac ,area_grid ,psi_opt , &
1188 con_fvirt ,con_eps ,con_epsm1 ,con_cp , &
1189 snowfall ,snow_water_equiv_old ,snow_albedo_old , &
1190 cm_noahmp ,ch_noahmp ,snow_levels ,snow_water_equiv , &
1191 soil_moisture_vol ,interface_depth ,snow_depth ,snow_level_ice , &
1192 snow_level_liquid ,temperature_ground ,temperature_snow_soil,soil_liquid_vol , &
1193 snow_age ,spec_humidity_surface,sw_absorbed_total ,sw_reflected_total , &
1194 lw_absorbed_total ,sensible_heat_total ,latent_heat_ground ,ground_heat_total , &
1195 temperature_radiative,evaporation_soil ,runoff_surface ,runoff_baseflow , &
1196 sw_absorbed_ground ,albedo_total ,snowmelt_out ,snowmelt_shallow , &
1197 snowmelt_shallow_1 ,snowmelt_shallow_2 ,temperature_bare_2m ,spec_humidity_bare_2m, &
1198 z0h_total , &
1199 emissivity_total ,precip_frozen_frac ,ch_bare_ground_2m ,snow_sublimation , &
1200#ifdef CCPP
1201 albedo_direct ,albedo_diffuse, errmsg ,errflg )
1202#else
1203 albedo_direct ,albedo_diffuse)
1204#endif
1205
1206#ifdef CCPP
1207 if (errflg /= 0) return
1208#endif
1209
1210!
1211! set some non-glacier fields over the glacier
1212!
1213
1214 snow_cover_fraction = 1.0
1215 temperature_leaf = temperature_radiative
1216 canopy_ice = 0.0
1217 canopy_liquid = 0.0
1218 vapor_pres_canopy_air = 2000.0
1219 temperature_canopy_air = temperature_radiative
1220 canopy_wet_fraction = 0.0
1221 lake_water = 0.0
1222 depth_water_table = 0.0
1223 aquifer_water = 0.0
1224 saturated_water = 0.0
1225 leaf_carbon = 0.0
1226 root_carbon = 0.0
1227 stem_carbon = 0.0
1228 wood_carbon = 0.0
1229 soil_carbon_stable = 0.0
1230 soil_carbon_fast = 0.0
1231 leaf_area_index = 0.0
1232 stem_area_index = 0.0
1233 evaporation_canopy = 0.0
1234 transpiration = 0.0
1235 aquifer_water = 0.0
1236 precip_adv_heat_total = 0.0
1237 soil_moisture_wtd = 0.0
1238 recharge = 0.0
1239 deep_recharge = 0.0
1240 eq_soil_water_vol = soil_moisture_vol
1241 transpiration_heat = 0.0
1242 latent_heat_canopy = 0.0
1243 z0_total = 0.002
1244 latent_heat_total = latent_heat_ground
1245 t2mmp(i) = temperature_bare_2m
1246 q2mp(i) = spec_humidity_bare_2m
1247
1248 tskin(i) = temperature_radiative
1249 tflux = temperature_ground
1250 surface_temperature = temperature_ground
1251 vegetation_fraction = vegetation_frac
1252 ch_vegetated = 0.0
1253 ch_bare_ground = ch_noahmp
1254 canopy_heat_storage = 0.0
1255 lai_sunlit = 0.0
1256 lai_shaded = 0.0
1257 rs_sunlit = 0.0
1258 rs_shaded = 0.0
1259
1260 else ! not glacier
1261
1262 ice_flag = 0
1263
1264 call noahmp_sflx (parameters , &
1265 i_location ,j_location ,latitude , &
1266 year_length ,julian_day ,cosine_zenith , &
1267 timestep ,spatial_scale ,atmosphere_thickness , &
1268 soil_levels ,soil_interface_depth ,max_snow_levels , &
1269 vegetation_frac ,max_vegetation_frac ,vegetation_category , &
1270 ice_flag ,surface_type ,crop_type , &
1271 eq_soil_water_vol ,temperature_forcing ,air_pressure_forcing , &
1272 air_pressure_surface ,uwind_forcing ,vwind_forcing , &
1273 spec_humidity_forcing ,area_grid ,cloud_water_forcing , &
1274 sw_radiation_forcing ,radiation_lw_forcing ,thsfc_loc , &
1275 prslkix ,prsik1x ,prslk1x , &
1276 pblhx ,iz0tlnd ,itime , &
1277 psi_opt , &
1278 precip_convective , &
1279 precip_non_convective ,precip_sh_convective ,precip_snow , &
1280 precip_graupel ,precip_hail ,temperature_soil_bot , &
1281 co2_air ,o2_air ,foliage_nitrogen , &
1282 snow_ice_frac_old ,forcing_height , &
1283 con_fvirt ,con_eps, con_epsm1 ,con_cp , &
1284 snow_albedo_old ,snow_water_equiv_old , &
1285 temperature_snow_soil ,soil_liquid_vol ,soil_moisture_vol , &
1286 temperature_canopy_air,vapor_pres_canopy_air ,canopy_wet_fraction , &
1287 canopy_liquid ,canopy_ice ,temperature_leaf , &
1288 temperature_ground ,spec_humidity_surface ,snowfall , &
1289 rainfall ,snow_levels ,interface_depth , &
1290 snow_depth ,snow_water_equiv ,snow_level_ice , &
1291 snow_level_liquid ,depth_water_table ,aquifer_water , &
1292 saturated_water , &
1293 lake_water ,leaf_carbon ,root_carbon , &
1294 stem_carbon ,wood_carbon ,soil_carbon_stable , &
1295 soil_carbon_fast ,leaf_area_index ,stem_area_index , &
1296 cm_noahmp ,ch_noahmp ,snow_age , &
1297 grain_carbon ,growing_deg_days ,plant_growth_stage , &
1298 soil_moisture_wtd ,deep_recharge ,recharge,ustarx , &
1299 z0_total ,z0h_total ,surface_temperature , &
1300 sw_absorbed_total ,sw_reflected_total , &
1301 lw_absorbed_total ,sensible_heat_total ,ground_heat_total , &
1302 latent_heat_canopy ,latent_heat_ground ,transpiration_heat , &
1303 evaporation_canopy ,transpiration ,evaporation_soil , &
1304 temperature_radiative ,temperature_bare_grd ,temperature_veg_grd , &
1305 temperature_veg_2m ,temperature_bare_2m ,spec_humidity_veg_2m , &
1306 spec_humidity_bare_2m ,runoff_surface ,runoff_baseflow , &
1307 par_absorbed ,photosynthesis ,sw_absorbed_veg , &
1308 sw_absorbed_ground ,snow_cover_fraction ,net_eco_exchange , &
1309 global_prim_prod ,net_prim_prod ,vegetation_fraction , &
1310 albedo_total ,snowmelt_out ,snowmelt_shallow , &
1311 snowmelt_shallow_1 ,snowmelt_shallow_2 ,rs_sunlit , &
1312 rs_shaded ,albedo_direct ,albedo_diffuse , &
1313 albedo_direct_snow ,albedo_diffuse_snow , &
1314 canopy_gap_fraction , &
1315 incanopy_gap_fraction ,ch_vegetated ,ch_bare_ground , &
1316 emissivity_total ,sensible_heat_grd_veg ,sensible_heat_leaf , &
1317 sensible_heat_grd_bar ,latent_heat_grd_veg ,latent_heat_grd_bare , &
1318 ground_heat_veg ,ground_heat_bare ,lw_absorbed_grd_veg , &
1319 lw_absorbed_leaf ,lw_absorbed_grd_bare ,latent_heat_trans , &
1320 latent_heat_leaf ,ch_leaf ,ch_below_canopy , &
1321 ch_vegetated_2m ,ch_bare_ground_2m ,precip_frozen_frac , &
1322 precip_adv_heat_veg ,precip_adv_heat_grd_v ,precip_adv_heat_grd_b , &
1323 precip_adv_heat_total ,snow_sublimation ,canopy_heat_storage , &
1324 lai_sunlit ,lai_shaded ,leaf_air_resistance , &
1325#ifdef CCPP
1326 spec_humid_sfc_veg ,spec_humid_sfc_bare , &
1327 errmsg ,errflg )
1328#else
1329 spec_humid_sfc_veg ,spec_humid_sfc_bare )
1330#endif
1331
1332#ifdef CCPP
1333 if (errflg /= 0) return
1334#endif
1335
1336 latent_heat_total = latent_heat_canopy + latent_heat_ground + transpiration_heat
1337
1338 t2mmp(i) = temperature_veg_2m * vegetation_fraction + &
1339 temperature_bare_2m * (1-vegetation_fraction)
1340 q2mp(i) = spec_humidity_veg_2m * vegetation_fraction + &
1341 spec_humidity_bare_2m * (1-vegetation_fraction)
1342
1343 tskin(i) = temperature_radiative
1344 tflux = surface_temperature
1345
1346 endif ! glacial split ends
1347
1348!
1349! --- noah-mp inout and out variables
1350!
1351
1352 tsnoxy(i,:) = temperature_snow_soil(-2: 0)
1353 stc(i,:) = temperature_snow_soil( 1:km)
1354 hflx(i) = sensible_heat_total !note unit change below
1355 evap(i) = latent_heat_total !note unit change below
1356 evbs(i) = latent_heat_ground
1357 evcw(i) = latent_heat_canopy
1358 trans(i) = transpiration_heat
1359 gflux(i) = -1.0*ground_heat_total ! opposite sign to be consistent with noah
1360 snohf(i) = snowmelt_out * con_hfus ! only snow that exits pack
1361 sbsno(i) = snow_sublimation
1362 pah(i) = precip_adv_heat_total
1363
1364 cmxy(i) = cm_noahmp
1365 chxy(i) = ch_noahmp
1366 zorl(i) = z0_total * 100.0 ! convert to cm
1367 ztmax(i) = z0h_total
1368
1369 !LAI-scale canopy resistance based on weighted sunlit shaded fraction
1370 if(rs_sunlit .le. 0.0 .or. rs_shaded .le. 0.0 .or. &
1371 lai_sunlit .eq. 0.0 .or. lai_shaded .eq. 0.0) then
1372 rca(i) = parameters%rsmax
1373 else !calculate LAI-scale canopy conductance (1/Rs)
1374 rca(i) = ((1.0/(rs_sunlit+leaf_air_resistance)*lai_sunlit) + &
1375 ((1.0/(rs_shaded+leaf_air_resistance))*lai_shaded))
1376 rca(i) = max((1.0/rca(i)),parameters%rsmin) !resistance
1377 end if
1378
1379 smc(i,:) = soil_moisture_vol
1380 slc(i,:) = soil_liquid_vol
1381 snowxy(i) = float(snow_levels)
1382 weasd(i) = snow_water_equiv
1383 snicexy(i,:) = snow_level_ice
1384 snliqxy(i,:) = snow_level_liquid
1385 snwdph(i) = snow_depth * 1000.0 ! convert from mm to m
1386 canopy(i) = canopy_ice + canopy_liquid
1387 canliqxy(i) = canopy_liquid
1388 canicexy(i) = canopy_ice
1389 zwtxy(i) = depth_water_table
1390 waxy(i) = aquifer_water
1391 wtxy(i) = saturated_water
1392 qsnowxy(i) = snowfall
1393 ecan(i) = evaporation_canopy
1394 etran(i) = transpiration
1395 edir(i) = evaporation_soil
1396 drain(i) = runoff_baseflow
1397 runoff(i) = runoff_surface
1398
1399 lfmassxy(i) = leaf_carbon
1400 rtmassxy(i) = root_carbon
1401 stmassxy(i) = stem_carbon
1402 woodxy(i) = wood_carbon
1403 stblcpxy(i) = soil_carbon_stable
1404 fastcpxy(i) = soil_carbon_fast
1405 xlaixy(i) = leaf_area_index
1406 xsaixy(i) = stem_area_index
1407
1408 snowc(i) = snow_cover_fraction
1409 sncovr1(i) = snow_cover_fraction
1410
1411 qsurf(i) = spec_humidity_surface
1412 tsurf(i) = tskin(i)
1413
1414 tvxy(i) = temperature_leaf
1415 tgxy(i) = temperature_ground
1416 tahxy(i) = temperature_canopy_air
1417 eahxy(i) = vapor_pres_canopy_air
1418 emiss(i) = emissivity_total
1419
1420 if(albedo_total > 0.0) then
1421 sfalb(i) = albedo_total
1422 albdvis(i) = albedo_direct(1)
1423 albdnir(i) = albedo_direct(2)
1424 albivis(i) = albedo_diffuse(1)
1425 albinir(i) = albedo_diffuse(2)
1426 end if
1427
1428 zsnsoxy(i,:) = interface_depth
1429
1430 if(present(canopy_heat_storage_ccpp )) canopy_heat_storage_ccpp(i) = canopy_heat_storage
1431 if(present(rainfall_ccpp )) rainfall_ccpp(i) = rainfall
1432 if(present(sw_absorbed_total_ccpp )) sw_absorbed_total_ccpp(i) = sw_absorbed_total
1433 if(present(sw_reflected_total_ccpp )) sw_reflected_total_ccpp(i) = sw_reflected_total
1434 if(present(lw_absorbed_total_ccpp )) lw_absorbed_total_ccpp(i) = lw_absorbed_total
1435 if(present(temperature_bare_grd_ccpp )) temperature_bare_grd_ccpp(i) = temperature_bare_grd
1436 if(present(temperature_veg_grd_ccpp )) temperature_veg_grd_ccpp(i) = temperature_veg_grd
1437 if(present(temperature_veg_2m_ccpp )) temperature_veg_2m_ccpp(i) = temperature_veg_2m
1438 if(present(temperature_bare_2m_ccpp )) temperature_bare_2m_ccpp(i) = temperature_bare_2m
1439 if(present(spec_humidity_veg_2m_ccpp )) spec_humidity_veg_2m_ccpp(i) = spec_humidity_veg_2m
1440 if(present(spec_humidity_bare_2m_ccpp)) spec_humidity_bare_2m_ccpp(i) = spec_humidity_bare_2m
1441 if(present(sw_absorbed_veg_ccpp )) sw_absorbed_veg_ccpp(i) = sw_absorbed_veg
1442 if(present(sw_absorbed_ground_ccpp )) sw_absorbed_ground_ccpp(i) = sw_absorbed_ground
1443 if(present(snowmelt_out_ccpp )) snowmelt_out_ccpp(i) = snowmelt_out
1444 if(present(snowmelt_shallow_ccpp )) snowmelt_shallow_ccpp(i) = snowmelt_shallow
1445 if(present(albedo_direct_snow_ccpp )) albedo_direct_snow_ccpp(i,:) = albedo_direct_snow
1446 if(present(albedo_diffuse_snow_ccpp )) albedo_diffuse_snow_ccpp(i,:) = albedo_diffuse_snow
1447 if(present(ch_vegetated_ccpp )) ch_vegetated_ccpp(i) = ch_vegetated
1448 if(present(ch_bare_ground_ccpp )) ch_bare_ground_ccpp(i) = ch_bare_ground
1449 if(present(sensible_heat_grd_veg_ccpp)) sensible_heat_grd_veg_ccpp(i) = sensible_heat_grd_veg
1450 if(present(sensible_heat_leaf_ccpp )) sensible_heat_leaf_ccpp(i) = sensible_heat_leaf
1451 if(present(sensible_heat_grd_bar_ccpp)) sensible_heat_grd_bar_ccpp(i) = sensible_heat_grd_bar
1452 if(present(latent_heat_grd_veg_ccpp )) latent_heat_grd_veg_ccpp(i) = latent_heat_grd_veg
1453 if(present(latent_heat_grd_bare_ccpp )) latent_heat_grd_bare_ccpp(i) = latent_heat_grd_bare
1454 if(present(ground_heat_veg_ccpp )) ground_heat_veg_ccpp(i) = ground_heat_veg
1455 if(present(ground_heat_bare_ccpp )) ground_heat_bare_ccpp(i) = ground_heat_bare
1456 if(present(lw_absorbed_grd_veg_ccpp )) lw_absorbed_grd_veg_ccpp(i) = lw_absorbed_grd_veg
1457 if(present(lw_absorbed_leaf_ccpp )) lw_absorbed_leaf_ccpp(i) = lw_absorbed_leaf
1458 if(present(lw_absorbed_grd_bare_ccpp )) lw_absorbed_grd_bare_ccpp(i) = lw_absorbed_grd_bare
1459 if(present(latent_heat_trans_ccpp )) latent_heat_trans_ccpp(i) = latent_heat_trans
1460 if(present(latent_heat_leaf_ccpp )) latent_heat_leaf_ccpp(i) = latent_heat_leaf
1461 if(present(ch_leaf_ccpp )) ch_leaf_ccpp(i) = ch_leaf
1462 if(present(ch_below_canopy_ccpp )) ch_below_canopy_ccpp(i) = ch_below_canopy
1463 if(present(ch_vegetated_2m_ccpp )) ch_vegetated_2m_ccpp(i) = ch_vegetated_2m
1464 if(present(ch_bare_ground_2m_ccpp )) ch_bare_ground_2m_ccpp(i) = ch_bare_ground_2m
1465 if(present(precip_adv_heat_veg_ccpp )) precip_adv_heat_veg_ccpp(i) = precip_adv_heat_veg
1466 if(present(precip_adv_heat_grd_v_ccpp)) precip_adv_heat_grd_v_ccpp(i) = precip_adv_heat_grd_v
1467 if(present(precip_adv_heat_grd_b_ccpp)) precip_adv_heat_grd_b_ccpp(i) = precip_adv_heat_grd_b
1468 if(present(spec_humid_sfc_veg_ccpp )) spec_humid_sfc_veg_ccpp(i) = spec_humid_sfc_veg
1469 if(present(spec_humid_sfc_bare_ccpp )) spec_humid_sfc_bare_ccpp(i) = spec_humid_sfc_bare
1470
1471 wslakexy(i) = lake_water ! not active
1472 fwetxy(i) = canopy_wet_fraction
1473 taussxy(i) = snow_age
1474 alboldxy(i) = snow_albedo_old
1475 sneqvoxy(i) = snow_water_equiv_old
1476
1477 smcwtdxy(i) = soil_moisture_wtd ! only need for run=5
1478 deeprechxy(i) = deep_recharge ! only need for run=5
1479 rechxy(i) = recharge ! only need for run=5
1480 smoiseq(i,:) = eq_soil_water_vol ! only need for run=5; listed as in
1481
1482 stm(i) = (0.1*soil_moisture_vol(1) + &
1483 0.3*soil_moisture_vol(2) + &
1484 0.6*soil_moisture_vol(3) + & ! clean up and use depths above
1485 1.0*soil_moisture_vol(4))*1000.0 ! unit conversion from m to kg m-2
1486
1487 wet1(i) = soil_moisture_vol(1) / smcmax_table(soil_category(1))
1488 smcwlt2(i) = smcdry_table(soil_category(1)) !!!change to wilt?
1489 smcref2(i) = smcref_table(soil_category(1))
1490
1491 virtfac1 = 1.0 + con_fvirt * max(q1(i), 1.e-8) !from forcing
1492
1493 if(thsfc_loc) then ! Use local potential temperature
1494 vptemp =temperature_forcing * prslki(i)*virtfac1 !virtual potential temperature @zlvl 1
1495 else ! Use potential temperature reference to 1000 hPa
1496 vptemp =temperature_forcing /prslk1(i) * virtfac1
1497 endif
1498
1499 if(thsfc_loc) then ! Use local potential temperature
1500 tvs1 = tflux * virtfac1
1501 else ! Use potential temperature referenced to 1000 hPa
1502 tvs1 = tflux/prsik1(i) * virtfac1
1503 endif
1504
1505 z0_total = max(min(z0_total,forcing_height),1.0e-6)
1506 z0h_total = max(z0h_total,1.0e-6)
1507
1508
1509 tem1 = (z0_total - z0lo) / (z0up - z0lo)
1510 tem1 = min(max(tem1, 0.0_kind_phys), 1.0_kind_phys)
1511 tem2 = max(vegetation_fraction, 0.1_kind_phys)
1512 zvfun(i) = sqrt(tem1 * tem2)
1513 gdx=sqrt(garea(i))
1514
1515! if ( .not. do_mynnsfclay) then !GFS sfcdiff
1516 if ( iopt_sfc .ne. 4 ) then !GFS sfcdiff
1517
1518 call gfs_stability &
1519 (zf(i), zvfun(i), gdx, virtual_temperature, vptemp,wind(i), z0_total, z0h_total, &
1520 tvs1, con_g, thsfc_loc, &
1521 rb1(i), fm1(i), fh1(i), fm101(i), fh21(i), cm(i), ch(i), stress1(i), ustar1(i))
1522
1523 rmol1(i) = undefined !not used in GFS sfcdif -> to satsify output
1524 flhc1(i) = undefined
1525 flqc1(i) = undefined
1526
1527 rah_total = max(1.0,1.0/( ch(i)*wind(i)) )
1528 cah_total = density * con_cp /rah_total
1529! tskin(i) = sensible_heat_total/cah_total + temperature_forcing ! test to use combined ch and SH to backout Ts
1530
1531! ch(i) = ch_vegetated * vegetation_frac + ch_bare_ground*(1.0-vegetation_frac)
1532
1533 else ! MYNN - note the GFS option is the same as sfcdif3; so removed.
1534
1535 qfx = evap(i) / con_hvap ! use flux from output
1536
1537 call sfcdif4(i_location ,j_location ,uwind_forcing ,vwind_forcing , &
1538 temperature_forcing, air_pressure_forcing ,air_pressure_surface , &
1539 pblhx,gdx,z0_total,con_fvirt,con_eps,con_cp,itime,snwdph(i),mnice, &
1540 psi_opt,surface_temperature, &
1541 spec_humidity_forcing,forcing_height,iz0tlnd,spec_humidity_surface,&
1542 sensible_heat_total,qfx,cm(i),ch(i),ch2,cq2,rmol1(i),ustar1(i), &
1543 rb1(i),fm1(i),fh1(i),stress1(i),fm101(i),fh21(i),wspd1,flhc1(i), &
1544 flqc1(i) )
1545
1546 ch(i)=ch(i)/wspd1
1547 cm(i)=cm(i)/wspd1
1548
1549 ch(i) = ch_vegetated * vegetation_fraction + ch_bare_ground*(1.0-vegetation_fraction)
1550
1551 rah_total = max(1.0,1.0/( ch(i)*wind(i)) )
1552 cah_total = density * con_cp /rah_total
1553
1554! tskin(i) = sensible_heat_total/cah_total + temperature_forcing !
1555
1556 endif
1557
1558
1559
1560 cmxy(i) = cm(i)
1561 chxy(i) = ch(i)
1562
1563 chh(i) = chxy(i) * wind(i) * density
1564 cmm(i) = cmxy(i) * wind(i)
1565
1566 snwdph(i) = snow_depth * 1000.0 ! convert from m to mm; wait after the stability call
1567! qsurf (i) = q1(i) + evap(i)/(con_hvap*density*ch(i)*wind(i))
1568
1569!
1570! --- change units for output
1571!
1572 hflx(i) = hflx(i) / density / con_cp
1573 evap(i) = evap(i) / density / con_hvap
1574
1575!
1576! --- calculate potential evaporation using noah code
1577!
1578 potential_temperature = temperature_forcing * prslki(i)
1579 virtual_temperature = temperature_forcing * (1.0 + 0.61*spec_humidity_forcing)
1580 penman_radiation = sw_absorbed_total + radiation_lw_forcing
1581 dqsdt = spec_humidity_sat * a23m4/(temperature_forcing-a4)**2
1582
1583 precip_freeze_frac_in = srflag(i)
1584 is_snowing = .false.
1585 is_freeze_rain = .false.
1586 if (precipitation_forcing > 0.0) then
1587 if (precip_freeze_frac_in > 0.0) then ! rain/snow flag, one condition is enough?
1588 is_snowing = .true.
1589 else
1590 if (temperature_forcing <= 275.15) is_freeze_rain = .true.
1591 end if
1592 end if
1593
1594!
1595! using new combined ch output to compute ep
1596!
1597 ch_noahmp = chxy(i) * wind(i)
1598
1599 call penman (temperature_forcing, air_pressure_forcing , ch_noahmp , &
1600 virtual_temperature, potential_temperature, precipitation_forcing, &
1601 penman_radiation , ground_heat_total , spec_humidity_forcing, &
1602 spec_humidity_sat , potential_evaporation, is_snowing , &
1603 is_freeze_rain , precip_freeze_frac_in, dqsdt , &
1604 emissivity_total , snow_cover_fraction )
1605
1606 ep(i) = potential_evaporation
1607
1608 end if ! flag_iter(i) .and. dry(i)
1609
1610 end do ! im loop
1611
1612 return
1613
1614 end subroutine noahmpdrv_run
1616!-----------------------------------
1617
1621 subroutine transfer_mp_parameters (vegtype,soiltype,slopetype, &
1622 soilcolor,croptype,parameters)
1623
1624 use noahmp_tables
1625 use module_sf_noahmplsm
1626
1627 implicit none
1628
1629 integer, intent(in) :: vegtype
1630 integer, intent(in) :: soiltype(4)
1631 integer, intent(in) :: slopetype
1632 integer, intent(in) :: soilcolor
1633 integer, intent(in) :: croptype
1634
1635 type (noahmp_parameters), intent(out) :: parameters
1636
1637 real :: refdk
1638 real :: refkdt
1639 real :: frzk
1640 real :: frzfact
1641 integer :: isoil
1642
1643 parameters%iswater = iswater_table
1644 parameters%isbarren = isbarren_table
1645 parameters%isice = isice_table
1646 parameters%iscrop = iscrop_table
1647 parameters%eblforest = eblforest_table
1648
1649!-----------------------------------------------------------------------&
1650 parameters%urban_flag = .false.
1651 if( vegtype == isurban_table .or. vegtype == 31 &
1652 & .or.vegtype == 32 .or. vegtype == 33) then
1653 parameters%urban_flag = .true.
1654 endif
1655
1656!------------------------------------------------------------------------------------------!
1657! transfer veg parameters
1658!------------------------------------------------------------------------------------------!
1659
1660 parameters%ch2op = ch2op_table(vegtype) !maximum intercepted h2o per unit lai+sai (mm)
1661 parameters%dleaf = dleaf_table(vegtype) !characteristic leaf dimension (m)
1662 parameters%z0mvt = z0mvt_table(vegtype) !momentum roughness length (m)
1663 parameters%hvt = hvt_table(vegtype) !top of canopy (m)
1664 parameters%hvb = hvb_table(vegtype) !bottom of canopy (m)
1665 parameters%z0mhvt = z0mhvt_table(vegtype) !momentum roughness length (m)
1666 parameters%den = den_table(vegtype) !tree density (no. of trunks per m2)
1667 parameters%rc = rc_table(vegtype) !tree crown radius (m)
1668 parameters%mfsno = mfsno_table(vegtype) !snowmelt m parameter ()
1669 parameters%scffac = scffac_table(vegtype) !snow cover factor
1670 parameters%cbiom = cbiom_table(vegtype) !canopy biomass heat capacity parameter (m)
1671 parameters%saim = saim_table(vegtype,:) !monthly stem area index, one-sided
1672 parameters%laim = laim_table(vegtype,:) !monthly leaf area index, one-sided
1673 parameters%sla = sla_table(vegtype) !single-side leaf area per kg [m2/kg]
1674 parameters%dilefc = dilefc_table(vegtype) !coeficient for leaf stress death [1/s]
1675 parameters%dilefw = dilefw_table(vegtype) !coeficient for leaf stress death [1/s]
1676 parameters%fragr = fragr_table(vegtype) !fraction of growth respiration !original was 0.3
1677 parameters%ltovrc = ltovrc_table(vegtype) !leaf turnover [1/s]
1678
1679 parameters%c3psn = c3psn_table(vegtype) !photosynthetic pathway: 0. = c4, 1. = c3
1680 parameters%kc25 = kc25_table(vegtype) !co2 michaelis-menten constant at 25c (pa)
1681 parameters%akc = akc_table(vegtype) !q10 for kc25
1682 parameters%ko25 = ko25_table(vegtype) !o2 michaelis-menten constant at 25c (pa)
1683 parameters%ako = ako_table(vegtype) !q10 for ko25
1684 parameters%vcmx25 = vcmx25_table(vegtype) !maximum rate of carboxylation at 25c (umol co2/m**2/s)
1685 parameters%avcmx = avcmx_table(vegtype) !q10 for vcmx25
1686 parameters%bp = bp_table(vegtype) !minimum leaf conductance (umol/m**2/s)
1687 parameters%mp = mp_table(vegtype) !slope of conductance-to-photosynthesis relationship
1688 parameters%qe25 = qe25_table(vegtype) !quantum efficiency at 25c (umol co2 / umol photon)
1689 parameters%aqe = aqe_table(vegtype) !q10 for qe25
1690 parameters%rmf25 = rmf25_table(vegtype) !leaf maintenance respiration at 25c (umol co2/m**2/s)
1691 parameters%rms25 = rms25_table(vegtype) !stem maintenance respiration at 25c (umol co2/kg bio/s)
1692 parameters%rmr25 = rmr25_table(vegtype) !root maintenance respiration at 25c (umol co2/kg bio/s)
1693 parameters%arm = arm_table(vegtype) !q10 for maintenance respiration
1694 parameters%folnmx = folnmx_table(vegtype) !foliage nitrogen concentration when f(n)=1 (%)
1695 parameters%tmin = tmin_table(vegtype) !minimum temperature for photosynthesis (k)
1696
1697 parameters%xl = xl_table(vegtype) !leaf/stem orientation index
1698 parameters%rhol = rhol_table(vegtype,:) !leaf reflectance: 1=vis, 2=nir
1699 parameters%rhos = rhos_table(vegtype,:) !stem reflectance: 1=vis, 2=nir
1700 parameters%taul = taul_table(vegtype,:) !leaf transmittance: 1=vis, 2=nir
1701 parameters%taus = taus_table(vegtype,:) !stem transmittance: 1=vis, 2=nir
1702
1703 parameters%mrp = mrp_table(vegtype) !microbial respiration parameter (umol co2 /kg c/ s)
1704 parameters%cwpvt = cwpvt_table(vegtype) !empirical canopy wind parameter
1705
1706 parameters%wrrat = wrrat_table(vegtype) !wood to non-wood ratio
1707 parameters%wdpool = wdpool_table(vegtype) !wood pool (switch 1 or 0) depending on woody or not [-]
1708 parameters%tdlef = tdlef_table(vegtype) !characteristic t for leaf freezing [k]
1709
1710 parameters%nroot = nroot_table(vegtype) !number of soil layers with root present
1711 parameters%rgl = rgl_table(vegtype) !parameter used in radiation stress function
1712 parameters%rsmin = rs_table(vegtype) !minimum stomatal resistance [s m-1]
1713 parameters%hs = hs_table(vegtype) !parameter used in vapor pressure deficit function
1714 parameters%topt = topt_table(vegtype) !optimum transpiration air temperature [k]
1715 parameters%rsmax = rsmax_table(vegtype) !maximal stomatal resistance [s m-1]
1716
1717!------------------------------------------------------------------------------------------!
1718! transfer rad parameters
1719!------------------------------------------------------------------------------------------!
1720
1721 parameters%albsat = albsat_table(soilcolor,:)
1722 parameters%albdry = albdry_table(soilcolor,:)
1723 parameters%albice = albice_table
1724 parameters%alblak = alblak_table
1725 parameters%omegas = omegas_table
1726 parameters%betads = betads_table
1727 parameters%betais = betais_table
1728 parameters%eg = eg_table
1729
1730!------------------------------------------------------------------------------------------!
1731! Transfer crop parameters
1732!------------------------------------------------------------------------------------------!
1733
1734 if(croptype > 0) then
1735 parameters%pltday = pltday_table(croptype) ! planting date
1736 parameters%hsday = hsday_table(croptype) ! harvest date
1737 parameters%plantpop = plantpop_table(croptype) ! plant density [per ha] - used?
1738 parameters%irri = irri_table(croptype) ! irrigation strategy 0= non-irrigation 1=irrigation (no water-stress)
1739 parameters%gddtbase = gddtbase_table(croptype) ! base temperature for gdd accumulation [c]
1740 parameters%gddtcut = gddtcut_table(croptype) ! upper temperature for gdd accumulation [c]
1741 parameters%gdds1 = gdds1_table(croptype) ! gdd from seeding to emergence
1742 parameters%gdds2 = gdds2_table(croptype) ! gdd from seeding to initial vegetative
1743 parameters%gdds3 = gdds3_table(croptype) ! gdd from seeding to post vegetative
1744 parameters%gdds4 = gdds4_table(croptype) ! gdd from seeding to intial reproductive
1745 parameters%gdds5 = gdds5_table(croptype) ! gdd from seeding to pysical maturity
1746 parameters%c3c4 = c3c4_table(croptype) ! photosynthetic pathway: 1. = c3 2. = c4
1747 parameters%aref = aref_table(croptype) ! reference maximum co2 assimulation rate
1748 parameters%psnrf = psnrf_table(croptype) ! co2 assimulation reduction factor(0-1) (e.g.pests, weeds)
1749 parameters%i2par = i2par_table(croptype) ! fraction of incoming solar radiation to photosynthetically active radiation
1750 parameters%tassim0 = tassim0_table(croptype) ! minimum temperature for co2 assimulation [c]
1751 parameters%tassim1 = tassim1_table(croptype) ! co2 assimulation linearly increasing until temperature reaches t1 [c]
1752 parameters%tassim2 = tassim2_table(croptype) ! co2 assmilation rate remain at aref until temperature reaches t2 [c]
1753 parameters%k = k_table(croptype) ! light extinction coefficient
1754 parameters%epsi = epsi_table(croptype) ! initial light use efficiency
1755 parameters%q10mr = q10mr_table(croptype) ! q10 for maintainance respiration
1756 parameters%foln_mx = foln_mx_table(croptype) ! foliage nitrogen concentration when f(n)=1 (%)
1757 parameters%lefreez = lefreez_table(croptype) ! characteristic t for leaf freezing [k]
1758 parameters%dile_fc = dile_fc_table(croptype,:) ! coeficient for temperature leaf stress death [1/s]
1759 parameters%dile_fw = dile_fw_table(croptype,:) ! coeficient for water leaf stress death [1/s]
1760 parameters%fra_gr = fra_gr_table(croptype) ! fraction of growth respiration
1761 parameters%lf_ovrc = lf_ovrc_table(croptype,:) ! fraction of leaf turnover [1/s]
1762 parameters%st_ovrc = st_ovrc_table(croptype,:) ! fraction of stem turnover [1/s]
1763 parameters%rt_ovrc = rt_ovrc_table(croptype,:) ! fraction of root tunrover [1/s]
1764 parameters%lfmr25 = lfmr25_table(croptype) ! leaf maintenance respiration at 25c [umol co2/m**2 /s]
1765 parameters%stmr25 = stmr25_table(croptype) ! stem maintenance respiration at 25c [umol co2/kg bio/s]
1766 parameters%rtmr25 = rtmr25_table(croptype) ! root maintenance respiration at 25c [umol co2/kg bio/s]
1767 parameters%grainmr25 = grainmr25_table(croptype) ! grain maintenance respiration at 25c [umol co2/kg bio/s]
1768 parameters%lfpt = lfpt_table(croptype,:) ! fraction of carbohydrate flux to leaf
1769 parameters%stpt = stpt_table(croptype,:) ! fraction of carbohydrate flux to stem
1770 parameters%rtpt = rtpt_table(croptype,:) ! fraction of carbohydrate flux to root
1771 parameters%grainpt = grainpt_table(croptype,:) ! fraction of carbohydrate flux to grain
1772 parameters%bio2lai = bio2lai_table(croptype) ! leaf are per living leaf biomass [m^2/kg]
1773 end if
1774
1775!------------------------------------------------------------------------------------------!
1776! transfer global parameters
1777!------------------------------------------------------------------------------------------!
1778
1779 parameters%co2 = co2_table
1780 parameters%o2 = o2_table
1781 parameters%timean = timean_table
1782 parameters%fsatmx = fsatmx_table
1783 parameters%z0sno = z0sno_table
1784 parameters%ssi = ssi_table
1785 parameters%snow_ret_fac = snow_ret_fac_table
1786 parameters%swemx = swemx_table
1787 parameters%tau0 = tau0_table
1788 parameters%grain_growth = grain_growth_table
1789 parameters%extra_growth = extra_growth_table
1790 parameters%dirt_soot = dirt_soot_table
1791 parameters%bats_cosz = bats_cosz_table
1792 parameters%bats_vis_new = bats_vis_new_table
1793 parameters%bats_nir_new = bats_nir_new_table
1794 parameters%bats_vis_age = bats_vis_age_table
1795 parameters%bats_nir_age = bats_nir_age_table
1796 parameters%bats_vis_dir = bats_vis_dir_table
1797 parameters%bats_nir_dir = bats_nir_dir_table
1798 parameters%rsurf_snow = rsurf_snow_table
1799 parameters%rsurf_exp = rsurf_exp_table
1800 parameters%snow_emis = snow_emis_table
1801
1802! ----------------------------------------------------------------------
1803! transfer soil parameters
1804! ----------------------------------------------------------------------
1805
1806 do isoil = 1, size(soiltype)
1807 parameters%bexp(isoil) = bexp_table(soiltype(isoil))
1808 parameters%dksat(isoil) = dksat_table(soiltype(isoil))
1809 parameters%dwsat(isoil) = dwsat_table(soiltype(isoil))
1810 parameters%psisat(isoil) = psisat_table(soiltype(isoil))
1811 parameters%quartz(isoil) = quartz_table(soiltype(isoil))
1812 parameters%smcdry(isoil) = smcdry_table(soiltype(isoil))
1813 parameters%smcmax(isoil) = smcmax_table(soiltype(isoil))
1814 parameters%smcref(isoil) = smcref_table(soiltype(isoil))
1815 parameters%smcwlt(isoil) = smcwlt_table(soiltype(isoil))
1816 end do
1817
1818 parameters%f1 = f1_table(soiltype(1))
1819 parameters%refdk = refdk_table
1820 parameters%refkdt = refkdt_table
1821
1822! ----------------------------------------------------------------------
1823! transfer genparm parameters
1824! ----------------------------------------------------------------------
1825 parameters%csoil = csoil_table
1826 parameters%zbot = zbot_table
1827 parameters%czil = czil_table
1828
1829 frzk = frzk_table
1830 parameters%kdt = parameters%refkdt * parameters%dksat(1) / parameters%refdk
1831 parameters%slope = slope_table(slopetype)
1832
1833 if(parameters%urban_flag)then ! hardcoding some urban parameters for soil
1834 parameters%smcmax = 0.45
1835 parameters%smcref = 0.42
1836 parameters%smcwlt = 0.40
1837 parameters%smcdry = 0.40
1838 parameters%csoil = 3.e6
1839 endif
1840
1841 ! adjust frzk parameter to actual soil type: frzk * frzfact
1842
1843!-----------------------------------------------------------------------&
1844 if(soiltype(1) /= 14) then
1845 frzfact = (parameters%smcmax(1) / parameters%smcref(1)) * (0.412 / 0.468)
1846 parameters%frzx = frzk * frzfact
1847 end if
1848
1849 end subroutine transfer_mp_parameters
1850
1853SUBROUTINE pedotransfer_sr2006(nsoil,sand,clay,orgm,parameters)
1854
1855 use module_sf_noahmplsm
1856 use noahmp_tables
1857
1858 implicit none
1859
1860 integer, intent(in ) :: nsoil ! number of soil layers
1861 real, dimension( 1:nsoil ), intent(inout) :: sand
1862 real, dimension( 1:nsoil ), intent(inout) :: clay
1863 real, dimension( 1:nsoil ), intent(inout) :: orgm
1864
1865 real, dimension( 1:nsoil ) :: theta_1500t
1866 real, dimension( 1:nsoil ) :: theta_1500
1867 real, dimension( 1:nsoil ) :: theta_33t
1868 real, dimension( 1:nsoil ) :: theta_33
1869 real, dimension( 1:nsoil ) :: theta_s33t
1870 real, dimension( 1:nsoil ) :: theta_s33
1871 real, dimension( 1:nsoil ) :: psi_et
1872 real, dimension( 1:nsoil ) :: psi_e
1873
1874 type(noahmp_parameters), intent(inout) :: parameters
1875 integer :: k
1876
1877 do k = 1,4
1878 if(sand(k) <= 0 .or. clay(k) <= 0) then
1879 sand(k) = 0.41
1880 clay(k) = 0.18
1881 end if
1882 if(orgm(k) <= 0 ) orgm(k) = 0.0
1883 end do
1884
1885 theta_1500t = sr2006_theta_1500t_a*sand &
1886 + sr2006_theta_1500t_b*clay &
1887 + sr2006_theta_1500t_c*orgm &
1888 + sr2006_theta_1500t_d*sand*orgm &
1889 + sr2006_theta_1500t_e*clay*orgm &
1890 + sr2006_theta_1500t_f*sand*clay &
1891 + sr2006_theta_1500t_g
1892
1893 theta_1500 = theta_1500t &
1894 + sr2006_theta_1500_a*theta_1500t &
1895 + sr2006_theta_1500_b
1896
1897 theta_33t = sr2006_theta_33t_a*sand &
1898 + sr2006_theta_33t_b*clay &
1899 + sr2006_theta_33t_c*orgm &
1900 + sr2006_theta_33t_d*sand*orgm &
1901 + sr2006_theta_33t_e*clay*orgm &
1902 + sr2006_theta_33t_f*sand*clay &
1903 + sr2006_theta_33t_g
1904
1905 theta_33 = theta_33t &
1906 + sr2006_theta_33_a*theta_33t*theta_33t &
1907 + sr2006_theta_33_b*theta_33t &
1908 + sr2006_theta_33_c
1909
1910 theta_s33t = sr2006_theta_s33t_a*sand &
1911 + sr2006_theta_s33t_b*clay &
1912 + sr2006_theta_s33t_c*orgm &
1913 + sr2006_theta_s33t_d*sand*orgm &
1914 + sr2006_theta_s33t_e*clay*orgm &
1915 + sr2006_theta_s33t_f*sand*clay &
1916 + sr2006_theta_s33t_g
1917
1918 theta_s33 = theta_s33t &
1919 + sr2006_theta_s33_a*theta_s33t &
1920 + sr2006_theta_s33_b
1921
1922 psi_et = sr2006_psi_et_a*sand &
1923 + sr2006_psi_et_b*clay &
1924 + sr2006_psi_et_c*theta_s33 &
1925 + sr2006_psi_et_d*sand*theta_s33 &
1926 + sr2006_psi_et_e*clay*theta_s33 &
1927 + sr2006_psi_et_f*sand*clay &
1928 + sr2006_psi_et_g
1929
1930 psi_e = psi_et &
1931 + sr2006_psi_e_a*psi_et*psi_et &
1932 + sr2006_psi_e_b*psi_et &
1933 + sr2006_psi_e_c
1934
1935 parameters%smcwlt = theta_1500
1936 parameters%smcref = theta_33
1937 parameters%smcmax = theta_33 &
1938 + theta_s33 &
1939 + sr2006_smcmax_a*sand &
1940 + sr2006_smcmax_b
1941
1942 parameters%bexp = 3.816712826 / (log(theta_33) - log(theta_1500) )
1943 parameters%psisat = psi_e
1944 parameters%dksat = 1930.0 * (parameters%smcmax - theta_33) ** (3.0 - 1.0/parameters%bexp)
1945 parameters%quartz = sand
1946
1947! Units conversion
1948
1949 parameters%psisat = max(0.1,parameters%psisat) ! arbitrarily impose a limit of 0.1kpa
1950 parameters%psisat = 0.101997 * parameters%psisat ! convert kpa to m
1951 parameters%dksat = parameters%dksat / 3600000.0 ! convert mm/h to m/s
1952 parameters%dwsat = parameters%dksat * parameters%psisat *parameters%bexp / parameters%smcmax ! units should be m*m/s
1953 parameters%smcdry = parameters%smcwlt
1954
1955! Introducing somewhat arbitrary limits (based on SOILPARM) to prevent bad things
1956
1957 parameters%smcmax = max(0.32 ,min(parameters%smcmax, 0.50 ))
1958 parameters%smcref = max(0.17 ,min(parameters%smcref,parameters%smcmax ))
1959 parameters%smcwlt = max(0.01 ,min(parameters%smcwlt,parameters%smcref ))
1960 parameters%smcdry = max(0.01 ,min(parameters%smcdry,parameters%smcref ))
1961 parameters%bexp = max(2.50 ,min(parameters%bexp, 12.0 ))
1962 parameters%psisat = max(0.03 ,min(parameters%psisat, 1.00 ))
1963 parameters%dksat = max(5.e-7,min(parameters%dksat, 1.e-5))
1964 parameters%dwsat = max(1.e-6,min(parameters%dwsat, 3.e-5))
1965 parameters%quartz = max(0.05 ,min(parameters%quartz, 0.95 ))
1966
1967 END SUBROUTINE pedotransfer_sr2006
1968
1969!-----------------------------------------------------------------------&
1970
1975 subroutine penman (sfctmp,sfcprs,ch,t2v,th2,prcp,fdown,ssoil, &
1976 & q2,q2sat,etp,snowng,frzgra,ffrozp, &
1977 & dqsdt2,emissi_in,sncovr)
1978
1979! etp is calcuated right after ssoil
1980
1981! ----------------------------------------------------------------------
1982! subroutine penman
1983! ----------------------------------------------------------------------
1984 use machine, only: kind_phys
1985 implicit none
1986 logical, intent(in) :: snowng, frzgra
1987 real(kind=kind_phys), intent(in) :: ch, dqsdt2,fdown,prcp,ffrozp, &
1988 & q2, q2sat,ssoil, sfcprs, sfctmp, &
1989 & t2v, th2,emissi_in,sncovr
1990 real(kind=kind_phys), intent(out) :: etp
1991 real(kind=kind_phys) :: epsca,flx2,rch,rr,t24
1992 real(kind=kind_phys) :: a, delta, fnet,rad,rho,emissi,elcp1,lvs
1993
1994 real(kind=kind_phys), parameter :: elcp = 2.4888e+3, lsubc = 2.501000e+6,cp = 1004.6
1995 real(kind=kind_phys), parameter :: lsubs = 2.83e+6, rd = 287.05, cph2o = 4.1855e+3
1996 real(kind=kind_phys), parameter :: cpice = 2.106e+3, lsubf = 3.335e5
1997 real(kind=kind_phys), parameter :: sigma = 5.6704e-8
1998
1999! ----------------------------------------------------------------------
2000! executable code begins here:
2001! ----------------------------------------------------------------------
2002! ----------------------------------------------------------------------
2003! prepare partial quantities for penman equation.
2004! ----------------------------------------------------------------------
2005 emissi=emissi_in
2006! elcp1 = (1.0-sncovr)*elcp + sncovr*elcp*lsubs/lsubc
2007 lvs = (1.0-sncovr)*lsubc + sncovr*lsubs
2008
2009 flx2 = 0.0
2010 delta = elcp * dqsdt2
2011! delta = elcp1 * dqsdt2
2012 t24 = sfctmp * sfctmp * sfctmp * sfctmp
2013 rr = t24 * 6.48e-8 / (sfcprs * ch) + 1.0
2014! rr = emissi*t24 * 6.48e-8 / (sfcprs * ch) + 1.0
2015 rho = sfcprs / (rd * t2v)
2016
2017! ----------------------------------------------------------------------
2018! adjust the partial sums / products with the latent heat
2019! effects caused by falling precipitation.
2020! ----------------------------------------------------------------------
2021 rch = rho * cp * ch
2022 if (.not. snowng) then
2023 if (prcp > 0.0) rr = rr + cph2o * prcp / rch
2024 else
2025! ---- ... fractional snowfall/rainfall
2026 rr = rr + (cpice*ffrozp+cph2o*(1.-ffrozp)) &
2027 & *prcp/rch
2028 end if
2029
2030! ----------------------------------------------------------------------
2031! include the latent heat effects of frzng rain converting to ice on
2032! impact in the calculation of flx2 and fnet.
2033! ----------------------------------------------------------------------
2034! fnet = fdown - sigma * t24- ssoil
2035 fnet = fdown - emissi*sigma * t24- ssoil
2036 if (frzgra) then
2037 flx2 = - lsubf * prcp
2038 fnet = fnet - flx2
2039! ----------------------------------------------------------------------
2040! finish penman equation calculations.
2041! ----------------------------------------------------------------------
2042 end if
2043 rad = fnet / rch + th2- sfctmp
2044 a = elcp * (q2sat - q2)
2045! a = elcp1 * (q2sat - q2)
2046 epsca = (a * rr + rad * delta) / (delta + rr)
2047 etp = epsca * rch / lsubc
2048! etp = epsca * rch / lvs
2049
2050! ----------------------------------------------------------------------
2051 end subroutine penman
2052
2053 end module noahmpdrv
subroutine penman
This subroutine calculates potential evaporation for the current point. various partial sums/products...
Definition sflx.f:1572
subroutine, public set_soilveg(me, isot, ivet, nlunit, errmsg, errflg)
This subroutine initializes soil and vegetation.
Definition set_soilveg.f:17
subroutine noahmpdrv_run(im, km, lsnowl, itime, ps, u1, v1, t1, q1, soiltyp, soilcol, vegtype, sigmaf, dlwflx, dswsfc, snet, delt, tg3, cm, ch, prsl1, prslk1, prslki, prsik1, zf, pblh, dry, wind, slopetyp, shdmin, shdmax, snoalb, sfalb, flag_iter, con_g, idveg, iopt_crs, iopt_btr, iopt_run, iopt_sfc, iopt_frz, iopt_inf, iopt_rad, iopt_alb, iopt_snf, iopt_tbot, iopt_stc, iopt_trs, iopt_diag, xlatin, xcoszin, iyrlen, julian, garea, rainn_mp, rainc_mp, snow_mp, graupel_mp, ice_mp, rhonewsn1, con_hvap, con_cp, con_jcal, rhoh2o, con_eps, con_epsm1, con_fvirt, con_rd, con_hfus, thsfc_loc, cpllnd, cpllnd2atm, weasd, snwdph, tskin, tprcp, srflag, smc, stc, slc, canopy, trans, tsurf, zorl, rb1, fm1, fh1, ustar1, stress1, fm101, fh21, rmol1, flhc1, flqc1, do_mynnsfclay, snowxy, tvxy, tgxy, canicexy, canliqxy, eahxy, tahxy, cmxy, chxy, fwetxy, sneqvoxy, alboldxy, qsnowxy, wslakexy, zwtxy, waxy, wtxy, tsnoxy, zsnsoxy, snicexy, snliqxy, lfmassxy, rtmassxy, stmassxy, woodxy, stblcpxy, fastcpxy, xlaixy, xsaixy, taussxy, smoiseq, smcwtdxy, deeprechxy, rechxy, albdvis, albdnir, albivis, albinir, emiss, sncovr1, qsurf, gflux, drain, evap, hflx, ep, runoff, cmm, chh, evbs, evcw, sbsno, pah, ecan, etran, edir, snowc, stm, snohf, smcwlt2, smcref2, wet1, t2mmp, q2mp, zvfun, ztmax, rca, errmsg, errflg, canopy_heat_storage_ccpp, rainfall_ccpp, sw_absorbed_total_ccpp, sw_reflected_total_ccpp, lw_absorbed_total_ccpp, temperature_bare_grd_ccpp, temperature_veg_grd_ccpp, temperature_veg_2m_ccpp, temperature_bare_2m_ccpp, spec_humidity_veg_2m_ccpp, spec_humidity_bare_2m_ccpp, sw_absorbed_veg_ccpp, sw_absorbed_ground_ccpp, snowmelt_out_ccpp, snowmelt_shallow_ccpp, albedo_direct_snow_ccpp, albedo_diffuse_snow_ccpp, ch_vegetated_ccpp, ch_bare_ground_ccpp, sensible_heat_grd_veg_ccpp, sensible_heat_leaf_ccpp, sensible_heat_grd_bar_ccpp, latent_heat_grd_veg_ccpp, latent_heat_grd_bare_ccpp, ground_heat_veg_ccpp, ground_heat_bare_ccpp, lw_absorbed_grd_veg_ccpp, lw_absorbed_leaf_ccpp, lw_absorbed_grd_bare_ccpp, latent_heat_trans_ccpp, latent_heat_leaf_ccpp, ch_leaf_ccpp, ch_below_canopy_ccpp, ch_vegetated_2m_ccpp, ch_bare_ground_2m_ccpp, precip_adv_heat_veg_ccpp, precip_adv_heat_grd_v_ccpp, precip_adv_heat_grd_b_ccpp, spec_humid_sfc_veg_ccpp, spec_humid_sfc_bare_ccpp)
This subroutine is the main CCPP entry point for the NoahMP LSM.
subroutine snowfall(parameters, nsoil, nsnow, dt, qsnow, snowhin, sfctmp, iloc, jloc, isnow, snowh, dzsnso, stc, snice, snliq, sneqv)
snow depth and density to account for the new snowfall. new values of snow depth & density returned.
subroutine noahmp_options(idveg, iopt_crs, iopt_btr, iopt_run, iopt_sfc, iopt_frz, iopt_inf, iopt_rad, iopt_alb, iopt_snf, iopt_tbot, iopt_stc, iopt_rsf, iopt_soil, iopt_pedo, iopt_crop, iopt_trs, iopt_diag, iopt_z0m)
subroutine psi_init(psi_opt, errmsg, errflg)
subroutine transfer_mp_parameters(vegtype, soiltype, slopetype, soilcolor, croptype, parameters)
This subroutine fills in a derived data type of type noahmp_parameters with data from the module noah...
subroutine sfcdif4(iloc, jloc, ux, vx, t1d, p1d, psfcpa, pblhx, dx, znt, ep_1, ep_2, cp, itime, snwh, isice, psi_opt, tsk, qx, zlvl, iz0tlnd, qsfc, hfx, qfx, cm, chs, chs2, cqs2, rmolx, ust, rbx, fmx, fhx, stressx, fm10x, fh2x, wspdx, flhcx, flqcx)
subroutine noahmpdrv_init(lsm, lsm_noahmp, me, isot, ivegsrc, nlunit, pores, resid, do_mynnsfclay, do_mynnedmf, errmsg, errflg, land_iau_control, land_iau_data, land_iau_state)
This subroutine is called during the CCPP initialization phase and calls set_soilveg() to initialize ...
Definition noahmpdrv.F90:45
subroutine pedotransfer_sr2006(nsoil, sand, clay, orgm, parameters)
This subroutine uses a pedotransfer method to calculate soil properties.
subroutine snow_age(parameters, dt, tg, sneqvo, sneqv, tauss, fage)
subroutine noahmpdrv_timestep_init(itime, fhour, delt, km, ncols, isot, ivegsrc, soiltyp, vegtype, weasd, land_iau_control, land_iau_data, land_iau_state, stc, slc, smc, errmsg, errflg, con_g, con_t0c, con_hfus)
This subroutine is called before noahmpdrv_run to update states with iau increments,...
subroutine noahmp_sflx(parameters, iloc, jloc, lat, yearlen, julian, cosz, dt, dx, dz8w, nsoil, zsoil, nsnow, shdfac, shdmax, vegtyp, ice, ist, croptype, smceq, sfctmp, sfcprs, psfc, uu, vv, q2, garea1, qc, soldn, lwdn, thsfc_loc, prslkix, prsik1x, prslk1x, pblhx, iz0tlnd, itime,psi_opt, prcpconv, prcpnonc, prcpshcv, prcpsnow, prcpgrpl, prcphail, tbot, co2air, o2air, foln, ficeold, zlvl, ep_1, ep_2, epsm1, cp, albold, sneqvo, stc, sh2o, smc, tah, eah, fwet, canliq, canice, tv, tg, qsfc, qsnow, qrain, isnow, zsnso, snowh, sneqv, snice, snliq, zwt, wa, wt, wslake, lfmass, rtmass, stmass, wood, stblcp, fastcp, lai, sai, cm, ch, tauss, grain, gdd, pgs, smcwtd,deeprech, rech, ustarx, z0wrf, z0hwrf, ts, fsa, fsr, fira, fsh, ssoil, fcev, fgev, fctr, ecan, etran, edir, trad, tgb, tgv, t2mv, t2mb, q2v, q2b, runsrf, runsub, apar, psn, sav, sag, fsno, nee, gpp, npp, fveg, albedo, qsnbot, ponding, ponding1, ponding2, rssun, rssha, albd, albi, albsnd, albsni, bgap, wgap, chv, chb, emissi, shg, shc, shb, evg, evb, ghv, ghb, irg, irc, irb, tr, evc, chleaf, chuc, chv2, chb2, fpice, pahv, pahg, pahb, pah, esnow, canhs, laisun, laisha, rb, qsfcveg, qsfcbare ifdef ccpp
subroutine noahmpdrv_finalize(land_iau_control, land_iau_data, land_iau_state, errmsg, errflg)
This subroutine mirrors noahmpdrv_init it calls land_iau_finalize which frees up allocated memory b...
TODO: replace with appropriate licence for CCPP.
This module contains the interface of noahmp_glacier_routines and noahmp_glacier_globals.
This module contains namelist options for Noah LSM.
brief Data from MPTABLE.TBL, SOILPARM.TBL, GENPARM.TBL for NoahMP
This module contains set_soilveg subroutine.
Definition set_soilveg.f:4