10module noahmp_glacier_globals
12 use machine ,
only : kind_phys
14 use module_sf_noahmplsm,
only :
sfcdif4
23 real (kind=kind_phys),
parameter :: grav = 9.80616
24 real (kind=kind_phys),
parameter :: sb = 5.67e-08
25 real (kind=kind_phys),
parameter :: vkc = 0.40
26 real (kind=kind_phys),
parameter :: tfrz = 273.16
27 real (kind=kind_phys),
parameter :: hsub = 2.8440e06
28 real (kind=kind_phys),
parameter :: hvap = 2.5104e06
29 real (kind=kind_phys),
parameter :: hfus = 0.3336e06
30 real (kind=kind_phys),
parameter :: cwat = 4.188e06
31 real (kind=kind_phys),
parameter :: cice = 2.094e06
32 real (kind=kind_phys),
parameter :: cpair = 1004.64
33 real (kind=kind_phys),
parameter :: tkwat = 0.6
34 real (kind=kind_phys),
parameter :: tkice = 2.2
35 real (kind=kind_phys),
parameter :: tkair = 0.023
36 real (kind=kind_phys),
parameter :: rair = 287.04
37 real (kind=kind_phys),
parameter :: rw = 461.269
38 real (kind=kind_phys),
parameter :: denh2o = 1000.
39 real (kind=kind_phys),
parameter :: denice = 917.
74 REAL,
PARAMETER :: Z0SNO = 0.002
75 REAL,
PARAMETER :: SSI = 0.03
76 REAL,
PARAMETER :: SWEMX = 1.00
80end module noahmp_glacier_globals
87 use noahmp_glacier_globals
129 iloc ,jloc ,cosz ,nsnow ,nsoil ,dt , & ! in : time/space/model-related
130 sfctmp ,sfcprs ,uu ,vv ,q2 ,soldn , & ! in : forcing
131 prcp ,lwdn ,tbot ,zlvl ,ficeold ,zsoil , & ! in : forcing
132 thsfc_loc ,prslkix ,prsik1x ,prslk1x , &
133 psfc ,pblhx ,iz0tlnd ,itime , &
134 sigmaf1 ,garea1 ,psi_opt , & ! in :
135 ep_1 ,ep_2 ,epsm1 ,cp , &
136 qsnow ,sneqvo ,albold ,cm ,ch ,isnow , & ! in/out :
137 sneqv ,smc ,zsnso ,snowh ,snice ,snliq , & ! in/out :
138 tg ,stc ,sh2o ,tauss ,qsfc , & ! in/out :
139 fsa ,fsr ,fira ,fsh ,fgev ,ssoil , & ! out :
140 trad ,edir ,runsrf ,runsub ,sag ,albedo , & ! out :
141 qsnbot ,ponding ,ponding1 ,ponding2 ,t2m,q2e ,z0h_total , & ! out :
143 emissi ,fpice ,ch2b , esnow , albsnd , albsni , &
146 emissi ,fpice ,ch2b , esnow. , albsnd , albsni)
157 integer ,
intent(in) :: iloc
158 integer ,
intent(in) :: jloc
159 real (kind=kind_phys) ,
intent(in) :: cosz
160 integer ,
intent(in) :: nsnow
161 integer ,
intent(in) :: nsoil
162 integer ,
intent(in) :: psi_opt
164 real (kind=kind_phys) ,
intent(in) :: dt
165 real (kind=kind_phys) ,
intent(in) :: sfctmp
166 real (kind=kind_phys) ,
intent(in) :: sfcprs
167 real (kind=kind_phys) ,
intent(in) :: uu
168 real (kind=kind_phys) ,
intent(in) :: vv
169 real (kind=kind_phys) ,
intent(in) :: q2
170 real (kind=kind_phys) ,
intent(in) :: soldn
171 real (kind=kind_phys) ,
intent(in) :: prcp
172 real (kind=kind_phys) ,
intent(in) :: lwdn
173 real (kind=kind_phys) ,
intent(in) :: tbot
174 real (kind=kind_phys) ,
intent(in) :: zlvl
175 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(in) :: ficeold
176 real (kind=kind_phys),
dimension( 1:nsoil),
intent(in) :: zsoil
177 logical ,
intent(in) :: thsfc_loc
178 real (kind=kind_phys) ,
intent(in) :: prslkix
179 real (kind=kind_phys) ,
intent(in) :: prsik1x
180 real (kind=kind_phys) ,
intent(in) :: prslk1x
182 real (kind=kind_phys) ,
intent(in) :: psfc
183 real (kind=kind_phys) ,
intent(in) :: pblhx
184 real (kind=kind_phys) ,
intent(in) :: ep_1
185 real (kind=kind_phys) ,
intent(in) :: ep_2
186 real (kind=kind_phys) ,
intent(in) :: epsm1
187 real (kind=kind_phys) ,
intent(in) :: cp
188 integer ,
intent(in) :: iz0tlnd
189 integer ,
intent(in) :: itime
191 real (kind=kind_phys) ,
intent(in) :: sigmaf1
192 real (kind=kind_phys) ,
intent(in) :: garea1
197 real (kind=kind_phys) ,
intent(inout) :: qsnow
198 real (kind=kind_phys) ,
intent(inout) :: sneqvo
199 real (kind=kind_phys) ,
intent(inout) :: albold
200 real (kind=kind_phys) ,
intent(inout) :: cm
201 real (kind=kind_phys) ,
intent(inout) :: ch
204 integer ,
intent(inout) :: isnow
205 real (kind=kind_phys) ,
intent(inout) :: sneqv
206 real (kind=kind_phys),
dimension( 1:nsoil),
intent(inout) :: smc
207 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: zsnso
208 real (kind=kind_phys) ,
intent(inout) :: snowh
209 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(inout) :: snice
210 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(inout) :: snliq
211 real (kind=kind_phys) ,
intent(inout) :: tg
212 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: stc
213 real (kind=kind_phys),
dimension( 1:nsoil),
intent(inout) :: sh2o
214 real (kind=kind_phys) ,
intent(inout) :: tauss
215 real (kind=kind_phys) ,
intent(inout) :: qsfc
218 real (kind=kind_phys) ,
intent(out) :: fsa
219 real (kind=kind_phys) ,
intent(out) :: fsr
220 real (kind=kind_phys) ,
intent(out) :: fira
221 real (kind=kind_phys) ,
intent(out) :: fsh
222 real (kind=kind_phys) ,
intent(out) :: fgev
223 real (kind=kind_phys) ,
intent(out) :: ssoil
224 real (kind=kind_phys) ,
intent(out) :: trad
225 real (kind=kind_phys) ,
intent(out) :: edir
226 real (kind=kind_phys) ,
intent(out) :: runsrf
227 real (kind=kind_phys) ,
intent(out) :: runsub
228 real (kind=kind_phys) ,
intent(out) :: sag
229 real (kind=kind_phys) ,
intent(out) ::
albedo
230 real (kind=kind_phys) ,
intent(out) :: qsnbot
231 real (kind=kind_phys) ,
intent(out) :: ponding
232 real (kind=kind_phys) ,
intent(out) :: ponding1
233 real (kind=kind_phys) ,
intent(out) :: ponding2
234 real (kind=kind_phys) ,
intent(out) :: t2m
235 real (kind=kind_phys) ,
intent(out) :: q2e
236 real (kind=kind_phys) ,
intent(out) :: z0h_total
237 real (kind=kind_phys) ,
intent(out) :: emissi
238 real (kind=kind_phys) ,
intent(out) :: fpice
239 real (kind=kind_phys) ,
intent(out) :: ch2b
240 real (kind=kind_phys) ,
intent(out) :: esnow
241 real (kind=kind_phys),
dimension(1:2) ,
intent(out) :: albsnd
242 real (kind=kind_phys),
dimension(1:2) ,
intent(out) :: albsni
246 character(len=*),
intent(inout) :: errmsg
247 integer,
intent(inout) :: errflg
252 integer,
dimension(-nsnow+1:nsoil) :: imelt
253 real (kind=kind_phys) :: rhoair
254 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: dzsnso
255 real (kind=kind_phys) :: thair
256 real (kind=kind_phys) :: qair
257 real (kind=kind_phys) :: eair
258 real (kind=kind_phys),
dimension( 1: 2) :: solad
259 real (kind=kind_phys),
dimension( 1: 2) :: solai
260 real (kind=kind_phys),
dimension( 1:nsoil) :: sice
261 real (kind=kind_phys),
dimension(-nsnow+1: 0) :: snicev
262 real (kind=kind_phys),
dimension(-nsnow+1: 0) :: snliqv
263 real (kind=kind_phys),
dimension(-nsnow+1: 0) :: epore
264 real (kind=kind_phys) :: qdew
265 real (kind=kind_phys) :: qvap
266 real (kind=kind_phys) :: lathea
267 real (kind=kind_phys) :: qmelt
268 real (kind=kind_phys) :: swdown
269 real (kind=kind_phys) :: beg_wb
270 real (kind=kind_phys) :: zbot = -8.0
272 character*256 message
277 call atm_glacier (ep_2, epsm1,sfcprs ,sfctmp ,q2 ,soldn ,cosz ,thair , &
278 qair ,eair ,rhoair ,solad ,solai ,swdown )
284 do iz = isnow+1, nsoil
285 if(iz == isnow+1)
then
286 dzsnso(iz) = - zsnso(iz)
288 dzsnso(iz) = zsnso(iz-1) - zsnso(iz)
295 eair ,sfcprs ,qair ,sfctmp ,lwdn ,uu , &
296 vv ,solad ,solai ,cosz ,zlvl , &
297 tbot ,zbot ,zsnso ,dzsnso ,sigmaf1 ,garea1 , &
298 thsfc_loc ,prslkix ,prsik1x ,prslk1x , &
299 psfc ,pblhx ,iz0tlnd ,itime ,psi_opt , &
300 ep_1, ep_2, epsm1, cp, &
301 tg ,stc ,snowh ,sneqv ,sneqvo ,sh2o , &
302 smc ,snice ,snliq ,albold ,cm ,ch , &
304 tauss ,qsfc ,errmsg ,errflg , &
308 imelt ,snicev ,snliqv ,epore ,qmelt ,ponding , &
309 sag ,fsa ,fsr ,fira ,fsh ,fgev , &
310 trad ,t2m ,ssoil ,lathea ,q2e ,emissi , &
311 ch2b ,albsnd ,albsni ,z0h_total)
314 if (errflg /= 0)
return
317 sice = max(0.0, smc - sh2o)
320 qvap = max( fgev/lathea, 0.)
321 qdew = abs( min(fgev/lathea, 0.))
326 call water_glacier (nsnow ,nsoil ,imelt ,dt ,prcp ,sfctmp , &
327 qvap ,qdew ,ficeold ,zsoil , &
328 isnow ,snowh ,sneqv ,snice ,snliq ,stc , &
329 dzsnso ,sh2o ,sice ,ponding ,zsnso ,fsh , &
330 runsrf ,runsub ,qsnow ,ponding1 ,ponding2 ,qsnbot , &
333 if(opt_gla == 2)
then
346 fsh ,fgev ,ssoil ,sag ,prcp ,edir , &
348 runsrf ,runsub ,sneqv ,dt ,beg_wb ,errmsg , errflg )
350 runsrf ,runsub ,sneqv ,dt ,beg_wb )
354 if (errflg /= 0)
return
357 if(snowh <= 1.e-6 .or. sneqv <= 1.e-3)
then
362 if(swdown.ne.0.)
then
373 subroutine atm_glacier (ep_2, epsm1, sfcprs ,sfctmp ,q2 ,soldn ,cosz ,thair , &
374 qair ,eair ,rhoair ,solad ,solai , &
383 real (kind=kind_phys) ,
intent(in) :: ep_2
384 real (kind=kind_phys) ,
intent(in) :: epsm1
385 real (kind=kind_phys) ,
intent(in) :: sfcprs
386 real (kind=kind_phys) ,
intent(in) :: sfctmp
387 real (kind=kind_phys) ,
intent(in) :: q2
388 real (kind=kind_phys) ,
intent(in) :: soldn
389 real (kind=kind_phys) ,
intent(in) :: cosz
393 real (kind=kind_phys) ,
intent(out) :: thair
394 real (kind=kind_phys) ,
intent(out) :: qair
395 real (kind=kind_phys) ,
intent(out) :: eair
396 real (kind=kind_phys),
dimension( 1: 2),
intent(out) :: solad
397 real (kind=kind_phys),
dimension( 1: 2),
intent(out) :: solai
398 real (kind=kind_phys) ,
intent(out) :: rhoair
399 real (kind=kind_phys) ,
intent(out) :: swdown
403 real (kind=kind_phys) :: pair
407 thair = sfctmp * (sfcprs/pair)**(rair/cpair)
411 eair = qair*sfcprs / (ep_2-epsm1*qair)
412 rhoair = (sfcprs+epsm1*eair) / (rair*sfctmp)
420 solad(1) = swdown*0.7*0.5
421 solad(2) = swdown*0.7*0.5
422 solai(1) = swdown*0.3*0.5
423 solai(2) = swdown*0.3*0.5
431 eair ,sfcprs ,qair ,sfctmp ,lwdn ,uu , & !in
432 vv ,solad ,solai ,cosz ,zref , & !in
433 tbot ,zbot ,zsnso ,dzsnso ,sigmaf1 ,garea1 , & !in
434 thsfc_loc ,prslkix ,prsik1x ,prslk1x , & !in
435 psfc ,pblhx ,iz0tlnd ,itime ,psi_opt , &
436 ep_1, ep_2, epsm1, cp, &
437 tg ,stc ,snowh ,sneqv ,sneqvo ,sh2o , & !inout
438 smc ,snice ,snliq ,albold ,cm ,ch , & !inout
440 tauss ,qsfc ,errmsg ,errflg , & !inout
442 tauss ,qsfc , & !inout
444 imelt ,snicev ,snliqv ,epore ,qmelt ,ponding , & !out
445 sag ,fsa ,fsr ,fira ,fsh ,fgev , & !out
446 trad ,t2m ,ssoil ,lathea ,q2e ,emissi , & !out
447 ch2b ,albsnd ,albsni ,z0h_total)
457 integer ,
intent(in) :: nsnow
458 integer ,
intent(in) :: nsoil
459 integer ,
intent(in) :: psi_opt
461 integer ,
intent(in) :: isnow
462 real (kind=kind_phys) ,
intent(in) :: dt
463 real (kind=kind_phys) ,
intent(in) :: qsnow
464 real (kind=kind_phys) ,
intent(in) :: rhoair
465 real (kind=kind_phys) ,
intent(in) :: eair
466 real (kind=kind_phys) ,
intent(in) :: sfcprs
467 real (kind=kind_phys) ,
intent(in) :: qair
468 real (kind=kind_phys) ,
intent(in) :: sfctmp
469 real (kind=kind_phys) ,
intent(in) :: lwdn
470 real (kind=kind_phys) ,
intent(in) :: uu
471 real (kind=kind_phys) ,
intent(in) :: vv
472 real (kind=kind_phys) ,
dimension( 1: 2),
intent(in) :: solad
473 real (kind=kind_phys) ,
dimension( 1: 2),
intent(in) :: solai
474 real (kind=kind_phys) ,
intent(in) :: cosz
475 real (kind=kind_phys) ,
intent(in) :: zref
476 real (kind=kind_phys) ,
intent(in) :: tbot
477 real (kind=kind_phys) ,
intent(in) :: zbot
478 real (kind=kind_phys) ,
dimension(-nsnow+1:nsoil),
intent(in) :: zsnso
479 real (kind=kind_phys) ,
dimension(-nsnow+1:nsoil),
intent(in) :: dzsnso
481 logical ,
intent(in) :: thsfc_loc
482 real (kind=kind_phys) ,
intent(in) :: prslkix
483 real (kind=kind_phys) ,
intent(in) :: prsik1x
484 real (kind=kind_phys) ,
intent(in) :: prslk1x
486 real (kind=kind_phys) ,
intent(in) :: pblhx
487 real (kind=kind_phys) ,
intent(in) :: psfc
488 real (kind=kind_phys) ,
intent(in) :: ep_1
489 real (kind=kind_phys) ,
intent(in) :: ep_2
490 real (kind=kind_phys) ,
intent(in) :: epsm1
491 real (kind=kind_phys) ,
intent(in) :: cp
492 integer ,
intent(in) :: iz0tlnd
493 integer ,
intent(in) :: itime
495 real (kind=kind_phys) ,
intent(in) :: sigmaf1
496 real (kind=kind_phys) ,
intent(in) :: garea1
499 real (kind=kind_phys) ,
intent(inout) :: tg
500 real (kind=kind_phys) ,
dimension(-nsnow+1:nsoil),
intent(inout) :: stc
501 real (kind=kind_phys) ,
intent(inout) :: snowh
502 real (kind=kind_phys) ,
intent(inout) :: sneqv
503 real (kind=kind_phys) ,
intent(inout) :: sneqvo
504 real (kind=kind_phys) ,
dimension( 1:nsoil),
intent(inout) :: sh2o
505 real (kind=kind_phys) ,
dimension( 1:nsoil),
intent(inout) :: smc
506 real (kind=kind_phys) ,
dimension(-nsnow+1: 0),
intent(inout) :: snice
507 real (kind=kind_phys) ,
dimension(-nsnow+1: 0),
intent(inout) :: snliq
508 real (kind=kind_phys) ,
intent(inout) :: albold
509 real (kind=kind_phys) ,
intent(inout) :: cm
510 real (kind=kind_phys) ,
intent(inout) :: ch
511 real (kind=kind_phys) ,
intent(inout) :: tauss
512 real (kind=kind_phys) ,
intent(inout) :: qsfc
515 character(len=*) ,
intent(inout) :: errmsg
516 integer ,
intent(inout) :: errflg
520 integer,
dimension(-nsnow+1:nsoil) ,
intent(out) :: imelt
521 real (kind=kind_phys) ,
dimension(-nsnow+1: 0),
intent(out) :: snicev
522 real (kind=kind_phys) ,
dimension(-nsnow+1: 0),
intent(out) :: snliqv
523 real (kind=kind_phys) ,
dimension(-nsnow+1: 0),
intent(out) :: epore
524 real (kind=kind_phys) ,
intent(out) :: qmelt
525 real (kind=kind_phys) ,
intent(out) :: ponding
526 real (kind=kind_phys) ,
intent(out) :: sag
527 real (kind=kind_phys) ,
intent(out) :: fsa
528 real (kind=kind_phys) ,
intent(out) :: fsr
529 real (kind=kind_phys) ,
intent(out) :: fira
530 real (kind=kind_phys) ,
intent(out) :: fsh
531 real (kind=kind_phys) ,
intent(out) :: fgev
532 real (kind=kind_phys) ,
intent(out) :: trad
533 real (kind=kind_phys) ,
intent(out) :: t2m
534 real (kind=kind_phys) ,
intent(out) :: ssoil
535 real (kind=kind_phys) ,
intent(out) :: lathea
536 real (kind=kind_phys) ,
intent(out) :: q2e
537 real (kind=kind_phys) ,
intent(out) :: emissi
538 real (kind=kind_phys) ,
intent(out) :: ch2b
539 real (kind=kind_phys),
dimension(1:2) ,
intent(out) :: albsnd
540 real (kind=kind_phys),
dimension(1:2) ,
intent(out) :: albsni
541 real (kind=kind_phys) ,
intent(out) :: z0h_total
545 real (kind=kind_phys) :: ur
546 real (kind=kind_phys) :: zlvl
547 real (kind=kind_phys) :: rsurf
548 real (kind=kind_phys) :: zpd
549 real (kind=kind_phys) :: z0mg
550 real (kind=kind_phys) :: emg
551 real (kind=kind_phys) :: fire
552 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: fact
553 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: df
554 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: hcpct
555 real (kind=kind_phys) :: gamma
556 real (kind=kind_phys) :: rhsur
562 ur = max( sqrt(uu**2.+vv**2.), 1. )
574 dt ,snowh ,snice ,snliq , &
575 df ,hcpct ,snicev ,snliqv ,epore , &
581 qsnow ,solad ,solai , &
583 sag ,fsr ,fsa ,albsnd ,albsni)
597 gamma = cpair*sfcprs/(ep_2*lathea)
601 call glacier_flux (nsoil ,nsnow ,emg ,isnow ,df ,dzsnso ,z0mg , &
602 zlvl ,zpd ,qair ,sfctmp ,rhoair ,sfcprs , &
603 ur ,gamma ,rsurf ,lwdn ,rhsur ,smc , &
604 eair ,stc ,sag ,snowh ,lathea ,sh2o , &
605 thsfc_loc ,prslkix ,prsik1x ,prslk1x , &
606 psfc ,pblhx ,iz0tlnd ,itime ,uu ,vv , &
607 sigmaf1 ,garea1 ,psi_opt ,ep_1, ep_2, epsm1, cp, &
609 cm ,ch ,tg ,qsfc ,errmsg ,errflg , &
613 fira ,fsh ,fgev ,ssoil , &
614 t2m ,q2e ,ch2b ,z0h_total)
623 errmsg =
"stop in noah-mp: emitted longwave <0"
626 call wrf_error_fatal(
"stop in noah-mp: emitted longwave <0")
637 trad = ( ( fire - (1-emissi)*lwdn ) / (emissi*sb) ) ** 0.25
642 ssoil ,snowh ,zbot ,zsnso ,df , &
647 if(opt_stc == 2)
then
648 if (snowh > 0.05 .and. tg > tfrz) tg = tfrz
655 stc ,snice ,snliq ,sneqv ,snowh , &
657 qmelt ,imelt ,ponding )
665 dt ,snowh ,snice ,snliq , & !in
666 df ,hcpct ,snicev ,snliqv ,epore , & !out
673 integer ,
intent(in) :: nsoil
674 integer ,
intent(in) :: nsnow
675 integer ,
intent(in) :: isnow
676 real (kind=kind_phys) ,
intent(in) :: dt
677 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(in) :: snice
678 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(in) :: snliq
679 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in) :: dzsnso
680 real (kind=kind_phys) ,
intent(in) :: snowh
683 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(out) :: df
684 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(out) :: hcpct
685 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(out) :: snicev
686 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(out) :: snliqv
687 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(out) :: epore
688 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(out) :: fact
693 real (kind=kind_phys),
dimension(-nsnow+1: 0) :: cvsno
694 real (kind=kind_phys),
dimension(-nsnow+1: 0) :: tksno
695 real (kind=kind_phys) :: zmid
700 call csnow_glacier (isnow ,nsnow ,nsoil ,snice ,snliq ,dzsnso , &
701 tksno ,cvsno ,snicev ,snliqv ,epore )
705 hcpct(iz) = cvsno(iz)
711 zmid = 0.5 * (dzsnso(iz))
713 zmid = zmid + dzsnso(iz2)
715 hcpct(iz) = 1.e6 * ( 0.8194 + 0.1309*zmid )
716 df(iz) = 0.32333 + ( 0.10073 * zmid )
721 do iz = isnow+1,nsoil
722 fact(iz) = dt/(hcpct(iz)*dzsnso(iz))
728 df(1) = (df(1)*dzsnso(1)+0.35*snowh) / (snowh +dzsnso(1))
730 df(1) = (df(1)*dzsnso(1)+df(0)*dzsnso(0)) / (dzsnso(0)+dzsnso(1))
740 tksno ,cvsno ,snicev ,snliqv ,epore )
748 integer,
intent(in) :: isnow
749 integer ,
intent(in) :: nsnow
750 integer ,
intent(in) :: nsoil
751 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(in) :: snice
752 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(in) :: snliq
753 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in) :: dzsnso
757 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(out) :: cvsno
758 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(out) :: tksno
759 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(out) :: snicev
760 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(out) :: snliqv
761 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(out) :: epore
766 real (kind=kind_phys),
dimension(-nsnow+1: 0) :: bdsnoi
772 snicev(iz) = min(1., snice(iz)/(dzsnso(iz)*denice) )
773 epore(iz) = 1. - snicev(iz)
774 snliqv(iz) = min(epore(iz),snliq(iz)/(dzsnso(iz)*denh2o))
778 bdsnoi(iz) = (snice(iz)+snliq(iz))/dzsnso(iz)
779 cvsno(iz) = cice*snicev(iz)+cwat*snliqv(iz)
798 qsnow ,solad ,solai , & !in
799 albold ,tauss , & !inout
800 sag ,fsr ,fsa ,albsnd ,albsni)
805 real (kind=kind_phys),
intent(in) :: dt
806 real (kind=kind_phys),
intent(in) :: tg
807 real (kind=kind_phys),
intent(in) :: sneqvo
808 real (kind=kind_phys),
intent(in) :: sneqv
809 real (kind=kind_phys),
intent(in) :: cosz
810 real (kind=kind_phys),
intent(in) :: qsnow
811 real (kind=kind_phys),
dimension(1:2) ,
intent(in) :: solad
812 real (kind=kind_phys),
dimension(1:2) ,
intent(in) :: solai
815 real (kind=kind_phys),
intent(inout) :: albold
816 real (kind=kind_phys),
intent(inout) :: tauss
817 real (kind=kind_phys),
dimension(1:2) :: albsnd
818 real (kind=kind_phys),
dimension(1:2) :: albsni
821 real (kind=kind_phys),
intent(out) :: sag
822 real (kind=kind_phys),
intent(out) :: fsr
823 real (kind=kind_phys),
intent(out) :: fsa
828 real (kind=kind_phys) :: fage
829 real (kind=kind_phys) :: alb
830 real (kind=kind_phys) :: abs
831 real (kind=kind_phys) :: ref
832 real (kind=kind_phys) :: fsno
833 real (kind=kind_phys),
dimension(1:2) :: albice
835 real (kind=kind_phys),
parameter :: mpe = 1.e-6
854 if(opt_alb == 2)
then
867 if(sneqv > 0.0) fsno = 1.0
873 albsnd(ib) = albice(ib)*(1.-fsno) + albsnd(ib)*fsno
874 albsni(ib) = albice(ib)*(1.-fsno) + albsni(ib)*fsno
878 abs = solad(ib)*(1.-albsnd(ib)) + solai(ib)*(1.-albsni(ib))
882 ref = solad(ib)*albsnd(ib) + solai(ib)*albsni(ib)
897 real (kind=kind_phys),
intent(in) :: dt
898 real (kind=kind_phys),
intent(in) :: tg
899 real (kind=kind_phys),
intent(in) :: sneqvo
900 real (kind=kind_phys),
intent(in) :: sneqv
903 real (kind=kind_phys),
intent(inout) :: tauss
906 real (kind=kind_phys),
intent(out) :: fage
909 real (kind=kind_phys) :: tage
910 real (kind=kind_phys) :: age1
911 real (kind=kind_phys) :: age2
912 real (kind=kind_phys) :: age3
913 real (kind=kind_phys) :: dela
914 real (kind=kind_phys) :: sge
915 real (kind=kind_phys) :: dels
916 real (kind=kind_phys) :: dela0
917 real (kind=kind_phys) :: arg
921 if(sneqv.le.0.0)
then
923 else if (sneqv.gt.800.)
then
928 arg = 5.e3*(1./tfrz-1./tg)
930 age2 = exp(amin1(0.,10.*arg))
932 tage = age1+age2+age3
934 dels = amax1(0.0,sneqv-sneqvo) / swemx
935 sge = (tauss+dela)*(1.0-dels)
936 tauss = amax1(0.,sge)
939 fage= tauss/(tauss+1.)
951 integer,
intent(in) :: nband
953 real (kind=kind_phys),
intent(in) :: cosz
954 real (kind=kind_phys),
intent(in) :: fage
958 real (kind=kind_phys),
dimension(1:2),
intent(out) :: albsnd
959 real (kind=kind_phys),
dimension(1:2),
intent(out) :: albsni
962 real (kind=kind_phys) :: fzen
963 real (kind=kind_phys) :: cf1
964 real (kind=kind_phys) :: sl2
965 real (kind=kind_phys) :: sl1
966 real (kind=kind_phys) :: sl
967 real (kind=kind_phys),
parameter :: c1 = 0.2
968 real (kind=kind_phys),
parameter :: c2 = 0.5
974 albsnd(1: nband) = 0.
975 albsni(1: nband) = 0.
982 cf1=((1.+sl1)/(1.+sl2*cosz)-sl1)
988 albsnd(1)=albsni(1)+0.4*fzen*(1.-albsni(1))
989 albsnd(2)=albsni(2)+0.4*fzen*(1.-albsni(2))
1001 integer,
intent(in) :: nband
1003 real (kind=kind_phys),
intent(in) :: qsnow
1004 real (kind=kind_phys),
intent(in) :: dt
1005 real (kind=kind_phys),
intent(in) :: albold
1009 real (kind=kind_phys),
intent(inout) :: alb
1012 real (kind=kind_phys),
dimension(1:2),
intent(out) :: albsnd
1013 real (kind=kind_phys),
dimension(1:2),
intent(out) :: albsni
1019 albsnd(1: nband) = 0.
1020 albsni(1: nband) = 0.
1024 alb = 0.55 + (albold-0.55) * exp(-0.01*dt/3600.)
1029 if (qsnow > 0.)
then
1030 alb = alb + min(qsnow*dt,swemx) * (0.84-alb)/(swemx)
1044 zlvl ,zpd ,qair ,sfctmp ,rhoair ,sfcprs , & !in
1045 ur ,gamma ,rsurf ,lwdn ,rhsur ,smc , & !in
1046 eair ,stc ,sag ,snowh ,lathea ,sh2o , & !in
1047 thsfc_loc ,prslkix ,prsik1x ,prslk1x , &
1048 psfc ,pblhx ,iz0tlnd ,itime ,uu ,vv , &
1049 sigmaf1 ,garea1 ,psi_opt ,ep_1, ep_2, epsm1, cp, & !in
1051 cm ,ch ,tgb ,qsfc ,errmsg ,errflg , & !inout
1053 cm ,ch ,tgb ,qsfc , & !inout
1055 irb ,shb ,evb ,ghb , & !out
1056 t2mb ,q2b ,ehb2 ,z0h_total)
1070 integer,
intent(in) :: nsnow
1071 integer,
intent(in) :: nsoil
1072 integer,
intent(in) :: psi_opt
1074 real (kind=kind_phys),
intent(in) :: emg
1075 integer,
intent(in) :: isnow
1076 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in) :: df
1077 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in) :: dzsnso
1078 real (kind=kind_phys),
intent(in) :: z0m
1079 real (kind=kind_phys),
intent(in) :: zlvl
1080 real (kind=kind_phys),
intent(in) :: zpd
1081 real (kind=kind_phys),
intent(in) :: qair
1082 real (kind=kind_phys),
intent(in) :: sfctmp
1083 real (kind=kind_phys),
intent(in) :: rhoair
1084 real (kind=kind_phys),
intent(in) :: sfcprs
1085 real (kind=kind_phys),
intent(in) :: ur
1086 real (kind=kind_phys),
intent(in) :: gamma
1087 real (kind=kind_phys),
intent(in) :: rsurf
1088 real (kind=kind_phys),
intent(in) :: lwdn
1089 real (kind=kind_phys),
intent(in) :: rhsur
1090 real (kind=kind_phys),
intent(in) :: eair
1091 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in) :: stc
1092 real (kind=kind_phys),
dimension( 1:nsoil),
intent(in) :: smc
1093 real (kind=kind_phys),
dimension( 1:nsoil),
intent(in) :: sh2o
1094 real (kind=kind_phys),
intent(in) :: sag
1095 real (kind=kind_phys),
intent(in) :: snowh
1096 real (kind=kind_phys),
intent(in) :: lathea
1098 logical ,
intent(in) :: thsfc_loc
1099 real (kind=kind_phys),
intent(in) :: prslkix
1100 real (kind=kind_phys),
intent(in) :: prsik1x
1101 real (kind=kind_phys),
intent(in) :: prslk1x
1103 real (kind=kind_phys) ,
intent(in) :: pblhx
1104 real (kind=kind_phys) ,
intent(in) :: psfc
1105 real (kind=kind_phys) ,
intent(in) :: ep_1
1106 real (kind=kind_phys) ,
intent(in) :: ep_2
1107 real (kind=kind_phys) ,
intent(in) :: epsm1
1108 real (kind=kind_phys) ,
intent(in) :: cp
1109 integer ,
intent(in) :: iz0tlnd
1110 integer ,
intent(in) :: itime
1111 real (kind=kind_phys) ,
intent(in) :: uu
1112 real (kind=kind_phys) ,
intent(in) :: vv
1114 real (kind=kind_phys),
intent(in) :: sigmaf1
1115 real (kind=kind_phys),
intent(in) :: garea1
1118 real (kind=kind_phys),
intent(inout) :: cm
1119 real (kind=kind_phys),
intent(inout) :: ch
1120 real (kind=kind_phys),
intent(inout) :: tgb
1121 real (kind=kind_phys),
intent(inout) :: qsfc
1124 character(len=*),
intent(inout) :: errmsg
1125 integer,
intent(inout) :: errflg
1130 real (kind=kind_phys),
intent(out) :: irb
1131 real (kind=kind_phys),
intent(out) :: shb
1132 real (kind=kind_phys),
intent(out) :: evb
1133 real (kind=kind_phys),
intent(out) :: ghb
1134 real (kind=kind_phys),
intent(out) :: t2mb
1135 real (kind=kind_phys),
intent(out) :: q2b
1136 real (kind=kind_phys),
intent(out) :: ehb2
1137 real (kind=kind_phys),
intent(out) :: z0h_total
1144 real (kind=kind_phys) :: mpe
1145 real (kind=kind_phys) :: dtg
1147 real (kind=kind_phys) :: mozold
1148 real (kind=kind_phys) :: fm2
1149 real (kind=kind_phys) :: fh2
1150 real (kind=kind_phys) :: ch2
1151 real (kind=kind_phys) :: h
1152 real (kind=kind_phys) :: fv
1153 real (kind=kind_phys) :: cir
1154 real (kind=kind_phys) :: cgh
1155 real (kind=kind_phys) :: csh
1156 real (kind=kind_phys) :: cev
1157 real (kind=kind_phys) :: cq2b
1159 real (kind=kind_phys) :: z0h
1161 real (kind=kind_phys) :: qfx
1162 real (kind=kind_phys) :: cq2
1165 real(kind=kind_phys) :: rb1i
1166 real(kind=kind_phys) :: fm10i
1168 real(kind=kind_phys) :: stress1i
1170 real(kind=kind_phys) :: wspd1i
1171 real(kind=kind_phys) :: flhc1i
1172 real(kind=kind_phys) :: flqc1i
1174 real(kind=kind_phys) :: tv1i
1176 real(kind=kind_phys) :: thv1i
1177 real(kind=kind_phys) :: tvsi
1178 real(kind=kind_phys) :: zlvli
1180 real(kind=kind_phys) :: snwd
1182 real(kind=kind_phys) :: reyni
1183 real(kind=kind_phys) :: virtfaci
1185 real(kind=kind_phys) :: tem1,tem2,zvfun1,gdx
1186 real(kind=kind_phys),
parameter :: z0lo=0.1, z0up=1.0
1188 real (kind=kind_phys) :: moz
1189 real (kind=kind_phys) :: fm
1190 real (kind=kind_phys) :: fh
1191 real (kind=kind_phys) :: ramb
1192 real (kind=kind_phys) :: rahb
1193 real (kind=kind_phys) :: rawb
1194 real (kind=kind_phys) :: estg
1195 real (kind=kind_phys) :: destg
1196 real (kind=kind_phys) :: esatw
1197 real (kind=kind_phys) :: esati
1198 real (kind=kind_phys) :: dsatw
1199 real (kind=kind_phys) :: dsati
1200 real (kind=kind_phys) :: a
1201 real (kind=kind_phys) :: b
1202 real (kind=kind_phys) :: t, tdc
1203 real (kind=kind_phys),
dimension( 1:nsoil) :: sice
1204 real (kind=kind_phys) :: czil
1206 tdc(t) = min( 50., max(-50.,(t-tfrz)) )
1233 fv = ur*vkc/log(zlvli/z0m)
1234 reyni = fv*z0m/(1.5e-05)
1236 if (opt_trs == 1)
then
1238 elseif (opt_trs == 2)
then
1239 z0h = z0m*exp(-czil*0.4*258.2*sqrt(fv*z0m))
1240 elseif (opt_trs == 3)
then
1242 elseif (opt_trs == 4)
then
1243 if (reyni .gt. 2.0)
then
1244 z0h = z0m/exp(2.46*(reyni)**0.25 - log(7.4))
1246 z0h = z0m/exp(-log(0.397))
1252 virtfaci = 1.0 + 0.61 * max(qair, 1.e-8)
1253 tv1i = sfctmp * virtfaci
1256 thv1i = sfctmp * prslkix * virtfaci
1258 thv1i = sfctmp / prslk1x * virtfaci
1261 if ( ur < 2.0) niter = 2
1264 cgh = 2.*df(isnow+1)/dzsnso(isnow+1)
1267 tem1 = (z0m - z0lo) / (z0up - z0lo)
1268 tem1 = min(max(tem1, 0.0_kind_phys), 1.0_kind_phys)
1269 tem2 = max(sigmaf1, 0.1_kind_phys)
1270 zvfun1= sqrt(tem1 * tem2)
1273 if(opt_sfc == 1 .or. opt_sfc == 2 .or. opt_sfc == 4)
then
1274 loop3:
do iter = 1, niterb
1275 if(opt_sfc == 1 .or. opt_sfc == 2)
then
1280 qair ,sfctmp ,h ,rhoair ,mpe ,ur , &
1282 & moz ,mozsgn ,fm ,fh ,fm2 ,fh2,fv, errmsg, errflg, &
1284 & moz ,mozsgn ,fm ,fh ,fm2 ,fh2,fv , &
1289 if (errflg /= 0)
return
1293 if(opt_sfc == 4)
then
1295 call sfcdif4(1 ,1 ,uu ,vv ,sfctmp , &
1296 sfcprs ,psfc ,pblhx ,gdx ,z0m , &
1298 itime ,snwd ,1 ,psi_opt, &
1299 tgb ,qair ,zlvl ,iz0tlnd,qsfc , &
1300 h ,qfx ,cm ,ch ,ch2 , &
1301 cq2 ,moz ,fv ,rb1i, fm, fh, &
1302 stress1i,fm10i ,fh2 , wspd1i ,flhc1i ,flqc1i)
1322 if(opt_sfc == 1 .or. opt_sfc == 2 .or. opt_sfc == 3)
then
1323 ramb = max(1.,1./(cm*ur))
1324 rahb = max(1.,1./(ch*ur))
1325 elseif(opt_sfc == 4)
then
1326 ramb = max(1.,1./(cm*wspd1i) )
1327 rahb = max(1.,1./(ch*wspd1i) )
1335 call esat(t, esatw, esati, dsatw, dsati)
1344 csh = rhoair*cpair/rahb
1345 if(snowh > 0.0 .or. opt_gla == 1)
then
1346 cev = rhoair*cpair/gamma/(rsurf+rawb)
1353 irb = cir * tgb**4 - emg*lwdn
1354 shb = csh * (tgb - sfctmp )
1355 evb = cev * (estg*rhsur - eair )
1356 ghb = cgh * (tgb - stc(isnow+1))
1358 b = sag-irb-shb-evb-ghb
1359 a = 4.*cir*tgb**3 + csh + cev*destg + cgh
1362 irb = irb + 4.*cir*tgb**3*dtg
1364 evb = evb + cev*destg*dtg
1371 h = csh * (tgb - sfctmp)
1374 call esat(t, esatw, esati, dsatw, dsati)
1380 qsfc = ep_2*(estg*rhsur)/(sfcprs+epsm1*(estg*rhsur))
1381 qfx = (qsfc-qair)*cev*gamma/cpair
1386 if (opt_sfc == 3)
then
1391 tvsi = tgb * virtfaci
1393 tvsi = tgb/prsik1x * virtfaci
1397 (zlvli, zvfun1, gdx,tv1i,thv1i, ur, z0m, z0h, tvsi, grav,thsfc_loc, &
1398 rb1i, fm,fh,fm10i,fh2,cm,ch,stress1i,fv)
1402 ramb = max(1.,1./(cm*ur))
1403 rahb = max(1.,1./(ch*ur))
1409 call esat(t, esatw, esati, dsatw, dsati)
1418 csh = rhoair*cpair/rahb
1420 if(snowh > 0.0 .or. opt_gla == 1)
then
1421 cev = rhoair*cpair/gamma/(rsurf+rawb)
1428 irb = cir * tgb**4 - emg*lwdn
1429 shb = csh * (tgb - sfctmp )
1430 evb = cev * (estg*rhsur - eair )
1431 ghb = cgh * (tgb - stc(isnow+1))
1433 b = sag-irb-shb-evb-ghb
1434 a = 4.*cir*tgb**3 + csh + cev*destg + cgh
1437 irb = irb + 4.*cir*tgb**3*dtg
1439 evb = evb + cev*destg*dtg
1447 call esat(t, esatw, esati, dsatw, dsati)
1453 qsfc = ep_2*(estg*rhsur)/(sfcprs+epsm1*(estg*rhsur))
1463 if(opt_stc == 1 .or. opt_stc ==3)
then
1464 if ((maxval(sice) > 0.0 .or. snowh > 0.0) .and. tgb > tfrz .and. opt_gla == 1)
then
1467 call esat(t, esatw, esati, dsatw, dsati)
1469 qsfc = ep_2*(estg*rhsur)/(sfcprs+epsm1*(estg*rhsur))
1470 irb = cir * tgb**4 - emg*lwdn
1471 shb = csh * (tgb - sfctmp)
1472 evb = cev * (estg*rhsur - eair )
1473 ghb = sag - (irb+shb+evb)
1478 ehb2 = fv*vkc/(log((2.+z0h)/z0h)-fh2)
1481 if (opt_sfc ==3)
then
1486 if (opt_sfc == 4)
then
1491 if (ehb2.lt.1.e-5 )
then
1495 t2mb = tgb - shb/(rhoair*cpair) * 1./ehb2
1496 q2b = qsfc - evb/(lathea*rhoair)*(1./cq2b + rsurf)
1505 subroutine esat(t, esw, esi, desw, desi)
1513 real (kind=kind_phys),
intent(in) :: t
1517 real (kind=kind_phys),
intent(out) :: esw
1518 real (kind=kind_phys),
intent(out) :: esi
1519 real (kind=kind_phys),
intent(out) :: desw
1520 real (kind=kind_phys),
intent(out) :: desi
1524 real (kind=kind_phys) :: a0,a1,a2,a3,a4,a5,a6
1525 real (kind=kind_phys) :: b0,b1,b2,b3,b4,b5,b6
1526 real (kind=kind_phys) :: c0,c1,c2,c3,c4,c5,c6
1527 real (kind=kind_phys) :: d0,d1,d2,d3,d4,d5,d6
1529 parameter(a0=6.107799961 , a1=4.436518521e-01, &
1530 a2=1.428945805e-02, a3=2.650648471e-04, &
1531 a4=3.031240396e-06, a5=2.034080948e-08, &
1534 parameter(b0=6.109177956 , b1=5.034698970e-01, &
1535 b2=1.886013408e-02, b3=4.176223716e-04, &
1536 b4=5.824720280e-06, b5=4.838803174e-08, &
1539 parameter(c0= 4.438099984e-01, c1=2.857002636e-02, &
1540 c2= 7.938054040e-04, c3=1.215215065e-05, &
1541 c4= 1.036561403e-07, c5=3.532421810e-10, &
1542 c6=-7.090244804e-13)
1544 parameter(d0=5.030305237e-01, d1=3.773255020e-02, &
1545 d2=1.267995369e-03, d3=2.477563108e-05, &
1546 d4=3.005693132e-07, d5=2.158542548e-09, &
1549 esw = 100.*(a0+t*(a1+t*(a2+t*(a3+t*(a4+t*(a5+t*a6))))))
1550 esi = 100.*(b0+t*(b1+t*(b2+t*(b3+t*(b4+t*(b5+t*b6))))))
1551 desw = 100.*(c0+t*(c1+t*(c2+t*(c3+t*(c4+t*(c5+t*c6))))))
1552 desi = 100.*(d0+t*(d1+t*(d2+t*(d3+t*(d4+t*(d5+t*d6))))))
1559 qair ,sfctmp ,h ,rhoair ,mpe ,ur , & !in
1561 & moz ,mozsgn ,fm ,fh ,fm2 ,fh2 , fv, & !inout
1562 & errmsg ,errflg , & !inout
1564 & moz ,mozsgn ,fm ,fh ,fm2 ,fh2 , fv, & !inout
1573 integer,
intent(in) :: iter
1574 real (kind=kind_phys),
intent(in) :: zlvl
1575 real (kind=kind_phys),
intent(in) :: zpd
1576 real (kind=kind_phys),
intent(in) :: z0h
1577 real (kind=kind_phys),
intent(in) :: z0m
1578 real (kind=kind_phys),
intent(in) :: qair
1579 real (kind=kind_phys),
intent(in) :: sfctmp
1580 real (kind=kind_phys),
intent(in) :: h
1581 real (kind=kind_phys),
intent(in) :: rhoair
1582 real (kind=kind_phys),
intent(in) :: mpe
1583 real (kind=kind_phys),
intent(in) :: ur
1586 real (kind=kind_phys),
intent(inout) :: moz
1587 integer,
intent(inout) :: mozsgn
1588 real (kind=kind_phys),
intent(inout) :: fm
1589 real (kind=kind_phys),
intent(inout) :: fh
1590 real (kind=kind_phys),
intent(inout) :: fm2
1591 real (kind=kind_phys),
intent(inout) :: fh2
1592 real (kind=kind_phys),
intent(inout) :: fv
1595 character(len=*),
intent(inout) :: errmsg
1596 integer,
intent(inout) :: errflg
1600 real (kind=kind_phys),
intent(out) :: cm
1601 real (kind=kind_phys),
intent(out) :: ch
1602 real (kind=kind_phys),
intent(out) :: ch2
1605 real (kind=kind_phys) :: mozold
1606 real (kind=kind_phys) :: tmpcm
1607 real (kind=kind_phys) :: tmpch
1608 real (kind=kind_phys) :: mol
1609 real (kind=kind_phys) :: tvir
1610 real (kind=kind_phys) :: tmp1,tmp2,tmp3
1611 real (kind=kind_phys) :: fmnew
1612 real (kind=kind_phys) :: fhnew
1613 real (kind=kind_phys) :: moz2
1614 real (kind=kind_phys) :: tmpcm2
1615 real (kind=kind_phys) :: tmpch2
1616 real (kind=kind_phys) :: fm2new
1617 real (kind=kind_phys) :: fh2new
1618 real (kind=kind_phys) :: tmp12,tmp22,tmp32
1620 real (kind=kind_phys) :: cmfm, chfh, cm2fm2, ch2fh2
1628 if(zlvl <= zpd)
then
1629 write(*,*)
'critical glacier problem: zlvl <= zpd; model stops', zlvl, zpd
1632 errmsg =
"stop in noah-mp glacier"
1635 call wrf_error_fatal(
"stop in noah-mp glacier")
1639 tmpcm = log((zlvl-zpd) / z0m)
1640 tmpch = log((zlvl-zpd) / z0h)
1641 tmpcm2 = log((2.0 + z0m) / z0m)
1642 tmpch2 = log((2.0 + z0h) / z0h)
1650 tvir = (1. + 0.61*qair) * sfctmp
1651 tmp1 = vkc * (grav/tvir) * h/(rhoair*cpair)
1652 if (abs(tmp1) .le. mpe) tmp1 = mpe
1653 mol = -1. * fv**3 / tmp1
1654 moz = min( (zlvl-zpd)/mol, 1.)
1655 moz2 = min( (2.0 + z0h)/mol, 1.)
1660 if (mozold*moz .lt. 0.) mozsgn = mozsgn+1
1661 if (mozsgn .ge. 2)
then
1671 if (moz .lt. 0.)
then
1672 tmp1 = (1. - 16.*moz)**0.25
1673 tmp2 = log((1.+tmp1*tmp1)/2.)
1674 tmp3 = log((1.+tmp1)/2.)
1675 fmnew = 2.*tmp3 + tmp2 - 2.*atan(tmp1) + 1.5707963
1679 tmp12 = (1. - 16.*moz2)**0.25
1680 tmp22 = log((1.+tmp12*tmp12)/2.)
1681 tmp32 = log((1.+tmp12)/2.)
1682 fm2new = 2.*tmp32 + tmp22 - 2.*atan(tmp12) + 1.5707963
1700 fm = 0.5 * (fm+fmnew)
1701 fh = 0.5 * (fh+fhnew)
1702 fm2 = 0.5 * (fm2+fm2new)
1703 fh2 = 0.5 * (fh2+fh2new)
1708 fh = min(fh,0.9*tmpch)
1709 fm = min(fm,0.9*tmpcm)
1710 fh2 = min(fh2,0.9*tmpch2)
1711 fm2 = min(fm2,0.9*tmpcm2)
1717 if(abs(cmfm) <= mpe) cmfm = mpe
1718 if(abs(chfh) <= mpe) chfh = mpe
1719 if(abs(cm2fm2) <= mpe) cm2fm2 = mpe
1720 if(abs(ch2fh2) <= mpe) ch2fh2 = mpe
1721 cm = vkc*vkc/(cmfm*cmfm)
1722 ch = vkc*vkc/(cmfm*chfh)
1723 ch2 = vkc*vkc/(cm2fm2*ch2fh2)
1734 ssoil ,snowh ,zbot ,zsnso ,df , & !in
1746 integer,
intent(in) :: nsoil
1747 integer,
intent(in) :: nsnow
1748 integer,
intent(in) :: isnow
1750 real (kind=kind_phys),
intent(in) :: dt
1751 real (kind=kind_phys),
intent(in) :: tbot
1752 real (kind=kind_phys),
intent(in) :: ssoil
1753 real (kind=kind_phys),
intent(in) :: snowh
1754 real (kind=kind_phys),
intent(in) :: zbot
1755 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in) :: zsnso
1756 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in) :: df
1757 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in) :: hcpct
1761 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: stc
1766 real (kind=kind_phys) :: zbotsno
1767 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: ai, bi, ci, rhsts
1768 real (kind=kind_phys) :: eflxb
1769 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: phi
1775 phi(isnow+1:nsoil) = 0.
1779 zbotsno = zbot - snowh
1784 stc ,tbot ,zbotsno ,df , &
1785 hcpct ,ssoil ,phi , &
1786 ai ,bi ,ci ,rhsts , &
1790 ai ,bi ,ci ,rhsts , &
1798 stc ,tbot ,zbot ,df , & !in
1799 hcpct ,ssoil ,phi , & !in
1800 ai ,bi ,ci ,rhsts , & !out
1812 integer,
intent(in) :: nsoil
1813 integer,
intent(in) :: nsnow
1814 integer,
intent(in) :: isnow
1815 real (kind=kind_phys),
intent(in) :: tbot
1816 real (kind=kind_phys),
intent(in) :: zbot
1818 real (kind=kind_phys),
intent(in) :: ssoil
1819 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in) :: zsnso
1820 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in) :: stc
1821 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in) :: df
1822 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in) :: hcpct
1823 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in) :: phi
1827 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(out) :: rhsts
1828 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(out) :: ai
1829 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(out) :: bi
1830 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(out) :: ci
1831 real (kind=kind_phys),
intent(out) :: botflx
1836 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: ddz
1837 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: denom
1838 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: dtsdz
1839 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: eflux
1840 real (kind=kind_phys) :: temp1
1843 do k = isnow+1, nsoil
1844 if (k == isnow+1)
then
1845 denom(k) = - zsnso(k) * hcpct(k)
1846 temp1 = - zsnso(k+1)
1847 ddz(k) = 2.0 / temp1
1848 dtsdz(k) = 2.0 * (stc(k) - stc(k+1)) / temp1
1849 eflux(k) = df(k) * dtsdz(k) - ssoil - phi(k)
1850 else if (k < nsoil)
then
1851 denom(k) = (zsnso(k-1) - zsnso(k)) * hcpct(k)
1852 temp1 = zsnso(k-1) - zsnso(k+1)
1853 ddz(k) = 2.0 / temp1
1854 dtsdz(k) = 2.0 * (stc(k) - stc(k+1)) / temp1
1855 eflux(k) = (df(k)*dtsdz(k) - df(k-1)*dtsdz(k-1)) - phi(k)
1856 else if (k == nsoil)
then
1857 denom(k) = (zsnso(k-1) - zsnso(k)) * hcpct(k)
1858 temp1 = zsnso(k-1) - zsnso(k)
1859 if(opt_tbot == 1)
then
1862 if(opt_tbot == 2)
then
1863 dtsdz(k) = (stc(k) - tbot) / ( 0.5*(zsnso(k-1)+zsnso(k)) - zbot)
1864 botflx = -df(k) * dtsdz(k)
1866 eflux(k) = (-botflx - df(k-1)*dtsdz(k-1) ) - phi(k)
1870 do k = isnow+1, nsoil
1871 if (k == isnow+1)
then
1873 ci(k) = - df(k) * ddz(k) / denom(k)
1874 if (opt_stc == 1 .or. opt_stc == 3)
then
1877 if (opt_stc == 2)
then
1878 bi(k) = - ci(k) + df(k)/(0.5*zsnso(k)*zsnso(k)*hcpct(k))
1880 else if (k < nsoil)
then
1881 ai(k) = - df(k-1) * ddz(k-1) / denom(k)
1882 ci(k) = - df(k ) * ddz(k ) / denom(k)
1883 bi(k) = - (ai(k) + ci(k))
1884 else if (k == nsoil)
then
1885 ai(k) = - df(k-1) * ddz(k-1) / denom(k)
1887 bi(k) = - (ai(k) + ci(k))
1889 rhsts(k) = eflux(k)/ (-denom(k))
1897 ai ,bi ,ci ,rhsts , & !inout
1906 integer,
intent(in) :: nsoil
1907 integer,
intent(in) :: nsnow
1908 integer,
intent(in) :: isnow
1909 real (kind=kind_phys),
intent(in) :: dt
1912 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: ai
1913 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: bi
1914 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: ci
1915 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: stc
1916 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: rhsts
1920 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: rhstsin
1921 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: ciin
1924 do k = isnow+1,nsoil
1925 rhsts(k) = rhsts(k) * dt
1927 bi(k) = 1. + bi(k) * dt
1933 do k = isnow+1,nsoil
1934 rhstsin(k) = rhsts(k)
1940 call rosr12_glacier (ci,ai,bi,ciin,rhstsin,rhsts,isnow+1,nsoil,nsnow)
1944 do k = isnow+1,nsoil
1945 stc(k) = stc(k) + ci(k)
1972 integer,
intent(in) :: ntop
1973 integer,
intent(in) :: nsoil,nsnow
1976 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in):: a, b, d
1977 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout):: c,p,delta
1983 p(ntop) = - c(ntop) / b(ntop)
1987 delta(ntop) = d(ntop) / b(ntop)
1992 p(k) = - c(k) * ( 1.0 / (b(k) + a(k) * p(k -1)) )
1993 delta(k) = (d(k) - a(k)* delta(k -1))* (1.0/ (b(k) + a(k)&
1999 p(nsoil) = delta(nsoil)
2004 kk = nsoil - k + (ntop-1) + 1
2005 p(kk) = p(kk) * p(kk +1) + delta(kk)
2014 stc ,snice ,snliq ,sneqv ,snowh , & !inout
2015 smc ,sh2o , & !inout
2016 qmelt ,imelt ,ponding )
2024 integer,
intent(in) :: nsnow
2025 integer,
intent(in) :: nsoil
2026 integer,
intent(in) :: isnow
2027 real (kind=kind_phys),
intent(in) :: dt
2028 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in) :: fact
2029 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in) :: dzsnso
2033 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: stc
2034 real (kind=kind_phys),
dimension(-nsnow+1:0) ,
intent(inout) :: snice
2035 real (kind=kind_phys),
dimension(-nsnow+1:0) ,
intent(inout) :: snliq
2036 real (kind=kind_phys),
intent(inout) :: sneqv
2037 real (kind=kind_phys),
intent(inout) :: snowh
2038 real (kind=kind_phys),
dimension( 1:nsoil),
intent(inout) :: sh2o
2039 real (kind=kind_phys),
dimension( 1:nsoil),
intent(inout) :: smc
2042 real (kind=kind_phys),
intent(out) :: qmelt
2043 integer,
dimension(-nsnow+1:nsoil),
intent(out) :: imelt
2044 real (kind=kind_phys),
intent(out) :: ponding
2049 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: hm
2050 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: xm
2051 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: wmass0
2052 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: wice0
2053 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: wliq0
2054 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: mice
2055 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: mliq
2056 real (kind=kind_phys),
dimension(-nsnow+1:nsoil) :: heatr
2057 real (kind=kind_phys) :: temp1
2058 real (kind=kind_phys) :: propor
2059 real (kind=kind_phys) :: xmf
2079 wmass0(j) = mice(j) + mliq(j)
2083 if (mice(j) > 0. .and. stc(j) >= tfrz)
then
2086 if (mliq(j) > 0. .and. stc(j) < tfrz)
then
2095 if (imelt(j) > 0)
then
2096 hm(j) = (stc(j)-tfrz)/fact(j)
2100 if (imelt(j) == 1 .and. hm(j) < 0.)
then
2104 if (imelt(j) == 2 .and. hm(j) > 0.)
then
2108 xm(j) = hm(j)*dt/hfus
2113if (opt_gla == 2)
then
2115 if (isnow == 0 .and. sneqv > 0. .and. stc(1) >= tfrz)
then
2116 hm(1) = (stc(1)-tfrz)/fact(1)
2118 xm(1) = hm(1)*dt/hfus
2121 sneqv = max(0.,temp1-xm(1))
2122 propor = sneqv/temp1
2123 snowh = max(0.,propor * snowh)
2124 heatr(1) = hm(1) - hfus*(temp1-sneqv)/dt
2125 if (heatr(1) > 0.)
then
2126 xm(1) = heatr(1)*dt/hfus
2127 stc(1) = stc(1) + fact(1)*heatr(1)
2132 qmelt = max(0.,(temp1-sneqv))/dt
2134 ponding = temp1-sneqv
2142 if (imelt(j) > 0 .and. abs(hm(j)) > 0.)
then
2145 if (xm(j) > 0.)
then
2146 mice(j) = max(0., wice0(j)-xm(j))
2147 heatr(j) = hm(j) - hfus*(wice0(j)-mice(j))/dt
2148 else if (xm(j) < 0.)
then
2149 mice(j) = min(wmass0(j), wice0(j)-xm(j))
2150 heatr(j) = hm(j) - hfus*(wice0(j)-mice(j))/dt
2153 mliq(j) = max(0.,wmass0(j)-mice(j))
2155 if (abs(heatr(j)) > 0.)
then
2156 stc(j) = stc(j) + fact(j)*heatr(j)
2157 if (mliq(j)*mice(j)>0.) stc(j) = tfrz
2160 qmelt = qmelt + max(0.,(wice0(j)-mice(j)))/dt
2165if (opt_gla == 1)
then
2168 mliq(j) = sh2o(j) * dzsnso(j) * 1000.
2169 mice(j) = (smc(j) - sh2o(j)) * dzsnso(j) * 1000.
2178 wmass0(j) = mice(j) + mliq(j)
2182 if (mice(j) > 0. .and. stc(j) >= tfrz)
then
2185 if (mliq(j) > 0. .and. stc(j) < tfrz)
then
2190 if (isnow == 0 .and. sneqv > 0. .and. j == 1)
then
2191 if (stc(j) >= tfrz)
then
2200 if (imelt(j) > 0)
then
2201 hm(j) = (stc(j)-tfrz)/fact(j)
2205 if (imelt(j) == 1 .and. hm(j) < 0.)
then
2209 if (imelt(j) == 2 .and. hm(j) > 0.)
then
2213 xm(j) = hm(j)*dt/hfus
2218 if (isnow == 0 .and. sneqv > 0. .and. xm(1) > 0.)
then
2220 sneqv = max(0.,temp1-xm(1))
2221 propor = sneqv/temp1
2222 snowh = max(0.,propor * snowh)
2223 heatr(1) = hm(1) - hfus*(temp1-sneqv)/dt
2224 if (heatr(1) > 0.)
then
2225 xm(1) = heatr(1)*dt/hfus
2233 qmelt = max(0.,(temp1-sneqv))/dt
2235 ponding = temp1-sneqv
2241 if (imelt(j) > 0 .and. abs(hm(j)) > 0.)
then
2244 if (xm(j) > 0.)
then
2245 mice(j) = max(0., wice0(j)-xm(j))
2246 heatr(j) = hm(j) - hfus*(wice0(j)-mice(j))/dt
2247 else if (xm(j) < 0.)
then
2248 mice(j) = min(wmass0(j), wice0(j)-xm(j))
2249 heatr(j) = hm(j) - hfus*(wice0(j)-mice(j))/dt
2252 mliq(j) = max(0.,wmass0(j)-mice(j))
2254 if (abs(heatr(j)) > 0.)
then
2255 stc(j) = stc(j) + fact(j)*heatr(j)
2257 if (mliq(j)*mice(j)>0.) stc(j) = tfrz
2261 if (j > 0) xmf = xmf + hfus * (wice0(j)-mice(j))/dt
2264 qmelt = qmelt + max(0.,(wice0(j)-mice(j)))/dt
2275 if (any(stc(1:4) > tfrz) .and. any(stc(1:4) < tfrz))
then
2277 if ( stc(j) > tfrz )
then
2278 heatr(j) = (stc(j)-tfrz)/fact(j)
2280 if (j .ne. k .and. stc(k) < tfrz .and. heatr(j) > 0.1)
then
2281 heatr(k) = (stc(k)-tfrz)/fact(k)
2282 if (abs(heatr(k)) > heatr(j))
then
2283 heatr(k) = heatr(k) + heatr(j)
2284 stc(k) = tfrz + heatr(k)*fact(k)
2287 heatr(j) = heatr(j) + heatr(k)
2293 stc(j) = tfrz + heatr(j)*fact(j)
2300 if (any(stc(1:4) > tfrz) .and. any(stc(1:4) < tfrz))
then
2302 if ( stc(j) < tfrz )
then
2303 heatr(j) = (stc(j)-tfrz)/fact(j)
2305 if (j .ne. k .and. stc(k) > tfrz .and. heatr(j) < -0.1)
then
2306 heatr(k) = (stc(k)-tfrz)/fact(k)
2307 if (heatr(k) > abs(heatr(j)))
then
2308 heatr(k) = heatr(k) + heatr(j)
2309 stc(k) = tfrz + heatr(k)*fact(k)
2312 heatr(j) = heatr(j) + heatr(k)
2318 stc(j) = tfrz + heatr(j)*fact(j)
2325 if (any(stc(1:4) > tfrz) .and. any(mice(1:4) > 0.))
then
2327 if ( stc(j) > tfrz )
then
2328 heatr(j) = (stc(j)-tfrz)/fact(j)
2329 xm(j) = heatr(j)*dt/hfus
2331 if (j .ne. k .and. mice(k) > 0. .and. xm(j) > 0.1)
then
2332 if (mice(k) > xm(j))
then
2333 mice(k) = mice(k) - xm(j)
2334 xmf = xmf + hfus * xm(j)/dt
2338 xm(j) = xm(j) - mice(k)
2339 xmf = xmf + hfus * mice(k)/dt
2343 mliq(k) = max(0.,wmass0(k)-mice(k))
2346 heatr(j) = xm(j)*hfus/dt
2347 stc(j) = tfrz + heatr(j)*fact(j)
2354 if (any(stc(1:4) < tfrz) .and. any(mliq(1:4) > 0.))
then
2356 if ( stc(j) < tfrz )
then
2357 heatr(j) = (stc(j)-tfrz)/fact(j)
2358 xm(j) = heatr(j)*dt/hfus
2360 if (j .ne. k .and. mliq(k) > 0. .and. xm(j) < -0.1)
then
2361 if (mliq(k) > abs(xm(j)))
then
2362 mice(k) = mice(k) - xm(j)
2363 xmf = xmf + hfus * xm(j)/dt
2367 xm(j) = xm(j) + mliq(k)
2368 xmf = xmf - hfus * mliq(k)/dt
2372 mliq(k) = max(0.,wmass0(k)-mice(k))
2375 heatr(j) = xm(j)*hfus/dt
2376 stc(j) = tfrz + heatr(j)*fact(j)
2389 if(opt_gla == 1)
then
2390 sh2o(j) = mliq(j) / (1000. * dzsnso(j))
2391 sh2o(j) = max(0.0,min(1.0,sh2o(j)))
2393 elseif(opt_gla == 2)
then
2403 qvap ,qdew ,ficeold ,zsoil , & !in
2404 isnow ,snowh ,sneqv ,snice ,snliq ,stc , & !inout
2405 dzsnso ,sh2o ,sice ,ponding ,zsnso ,fsh , & !inout
2406 runsrf ,runsub ,qsnow ,ponding1 ,ponding2 ,qsnbot , & !out
2415 integer,
intent(in) :: nsnow
2416 integer,
intent(in) :: nsoil
2417 integer,
dimension(-nsnow+1:0) ,
intent(in) :: imelt
2418 real (kind=kind_phys),
intent(in) :: dt
2419 real (kind=kind_phys),
intent(in) :: prcp
2420 real (kind=kind_phys),
intent(in) :: sfctmp
2421 real (kind=kind_phys),
intent(inout) :: qvap
2422 real (kind=kind_phys),
intent(inout) :: qdew
2423 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(in) :: ficeold
2424 real (kind=kind_phys),
dimension( 1:nsoil),
intent(in) :: zsoil
2427 integer,
intent(inout) :: isnow
2428 real (kind=kind_phys),
intent(inout) :: snowh
2429 real (kind=kind_phys),
intent(inout) :: sneqv
2430 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(inout) :: snice
2431 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(inout) :: snliq
2432 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: stc
2433 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: dzsnso
2434 real (kind=kind_phys),
dimension( 1:nsoil),
intent(inout) :: sh2o
2435 real (kind=kind_phys),
dimension( 1:nsoil),
intent(inout) :: sice
2436 real (kind=kind_phys) ,
intent(inout) :: ponding
2437 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: zsnso
2438 real (kind=kind_phys) ,
intent(inout) :: fsh
2441 real (kind=kind_phys),
intent(out) :: runsrf
2442 real (kind=kind_phys),
intent(out) :: runsub
2443 real (kind=kind_phys),
intent(out) :: qsnow
2444 real (kind=kind_phys),
intent(out) :: ponding1
2445 real (kind=kind_phys),
intent(out) :: ponding2
2446 real (kind=kind_phys),
intent(out) :: qsnbot
2447 real (kind=kind_phys),
intent(out) :: fpice
2448 real (kind=kind_phys),
intent(out) :: esnow
2451 real (kind=kind_phys) :: qrain
2452 real (kind=kind_phys) :: qseva
2453 real (kind=kind_phys) :: qsdew
2454 real (kind=kind_phys) :: qsnfro
2455 real (kind=kind_phys) :: qsnsub
2456 real (kind=kind_phys) :: snowhin
2457 real (kind=kind_phys) :: snoflow
2458 real (kind=kind_phys) :: bdfall
2459 real (kind=kind_phys) :: replace
2460 real (kind=kind_phys),
dimension( 1:nsoil) :: sice_save
2461 real (kind=kind_phys),
dimension( 1:nsoil) :: sh2o_save
2479 if(opt_snf == 1 .or. opt_snf == 4)
then
2480 if(sfctmp > tfrz+2.5)
then
2483 if(sfctmp <= tfrz+0.5)
then
2485 else if(sfctmp <= tfrz+2.)
then
2486 fpice = 1.-(-54.632 + 0.2*sfctmp)
2493 if(opt_snf == 2)
then
2494 if(sfctmp >= tfrz+2.2)
then
2501 if(opt_snf == 3)
then
2502 if(sfctmp >= tfrz)
then
2513 bdfall = min(120.,67.92+51.25*exp((sfctmp-tfrz)/2.59))
2515 qrain = prcp * (1.-fpice)
2516 qsnow = prcp * fpice
2517 snowhin = qsnow/bdfall
2524 esnow = qsnsub*2.83e+6
2527 snowhin ,qsnow ,qsnfro ,qsnsub ,qrain , &
2529 isnow ,snowh ,sneqv ,snice ,snliq , &
2530 sh2o ,sice ,stc ,dzsnso ,zsnso , &
2532 qsnbot ,snoflow ,ponding1 ,ponding2)
2536 runsrf = (ponding+ponding1+ponding2)/dt
2539 runsrf = runsrf + qsnbot + qrain
2541 runsrf = runsrf + qsnbot
2545 if(opt_gla == 1)
then
2548 replace = replace + dzsnso(ilev)*(sice(ilev) - sice_save(ilev) + sh2o(ilev) - sh2o_save(ilev))
2550 replace = replace * 1000.0 / dt
2552 sice = min(1.0,sice_save)
2553 elseif(opt_gla == 2)
then
2561 if(opt_gla == 1)
then
2562 runsub = snoflow + replace
2563 elseif(opt_gla == 2)
then
2574 snowhin ,qsnow ,qsnfro ,qsnsub ,qrain , & !in
2575 ficeold ,zsoil , & !in
2576 isnow ,snowh ,sneqv ,snice ,snliq , & !inout
2577 sh2o ,sice ,stc ,dzsnso ,zsnso , & !inout
2579 qsnbot ,snoflow ,ponding1 ,ponding2)
2584 integer,
intent(in) :: nsnow
2585 integer,
intent(in) :: nsoil
2586 integer,
dimension(-nsnow+1:0) ,
intent(in) :: imelt
2587 real (kind=kind_phys),
intent(in) :: dt
2588 real (kind=kind_phys),
intent(in) :: sfctmp
2589 real (kind=kind_phys),
intent(in) :: snowhin
2590 real (kind=kind_phys),
intent(in) :: qsnow
2591 real (kind=kind_phys),
intent(inout) :: qsnfro
2592 real (kind=kind_phys),
intent(inout) :: qsnsub
2593 real (kind=kind_phys),
intent(in) :: qrain
2594 real (kind=kind_phys),
dimension(-nsnow+1:0) ,
intent(in) :: ficeold
2595 real (kind=kind_phys),
dimension( 1:nsoil),
intent(in) :: zsoil
2598 integer,
intent(inout) :: isnow
2599 real (kind=kind_phys),
intent(inout) :: snowh
2600 real (kind=kind_phys),
intent(inout) :: sneqv
2601 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(inout) :: snice
2602 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(inout) :: snliq
2603 real (kind=kind_phys),
dimension( 1:nsoil),
intent(inout) :: sh2o
2604 real (kind=kind_phys),
dimension( 1:nsoil),
intent(inout) :: sice
2605 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: stc
2606 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: dzsnso
2607 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: zsnso
2608 real (kind=kind_phys),
intent(inout) :: fsh
2611 real (kind=kind_phys),
intent(out) :: qsnbot
2612 real (kind=kind_phys),
intent(out) :: snoflow
2613 real (kind=kind_phys),
intent(out) :: ponding1
2614 real (kind=kind_phys),
intent(out) :: ponding2
2618 real (kind=kind_phys) :: bdsnow
2619 real (kind=kind_phys),
parameter :: mwd = 600.
2627 isnow ,snowh ,dzsnso ,stc ,snice , &
2632 snliq ,imelt ,ficeold, &
2637 isnow ,sh2o ,stc ,snice ,snliq , &
2638 dzsnso ,sice ,snowh ,sneqv , &
2643 isnow ,stc ,snice ,snliq ,dzsnso )
2648 isnow ,dzsnso ,snowh ,sneqv ,snice , &
2649 snliq ,sh2o ,sice ,stc , &
2650 ponding1 ,ponding2 ,fsh , &
2662 if(stc( 0) < 100.0) stc( 0) = stc( 1)
2663 if(stc(-1) < 100.0) stc(-1) = stc( 0)
2664 if(stc(-2) < 100.0) stc(-2) = stc(-1)
2672 do iz = -nsnow+1, isnow
2682 if(sneqv > mwd)
then
2683 bdsnow = snice(0) / dzsnso(0)
2684 snoflow = (sneqv - mwd)
2685 snice(0) = snice(0) - snoflow
2686 dzsnso(0) = dzsnso(0) - snoflow/bdsnow
2687 snoflow = snoflow / dt
2696 sneqv = sneqv + snice(iz) + snliq(iz)
2697 snowh = snowh + dzsnso(iz)
2704 dzsnso(iz) = -dzsnso(iz)
2707 dzsnso(1) = zsoil(1)
2709 dzsnso(iz) = (zsoil(iz) - zsoil(iz-1))
2712 zsnso(isnow+1) = dzsnso(isnow+1)
2713 do iz = isnow+2 ,nsoil
2714 zsnso(iz) = zsnso(iz-1) + dzsnso(iz)
2717 do iz = isnow+1 ,nsoil
2718 dzsnso(iz) = -dzsnso(iz)
2726 isnow ,snowh ,dzsnso ,stc ,snice , & !inout
2736 integer,
intent(in) :: nsoil
2737 integer,
intent(in) :: nsnow
2738 real (kind=kind_phys),
intent(in) :: dt
2739 real (kind=kind_phys),
intent(in) :: qsnow
2740 real (kind=kind_phys),
intent(in) :: snowhin
2741 real (kind=kind_phys),
intent(in) :: sfctmp
2745 integer,
intent(inout) :: isnow
2746 real (kind=kind_phys),
intent(inout) :: snowh
2747 real (kind=kind_phys),
intent(inout) :: sneqv
2748 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: dzsnso
2749 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: stc
2750 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(inout) :: snice
2751 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(inout) :: snliq
2761 if(isnow == 0 .and. qsnow > 0.)
then
2762 snowh = snowh + snowhin * dt
2763 sneqv = sneqv + qsnow * dt
2768 if(isnow == 0 .and. qsnow>0. .and. snowh >= 0.025)
then
2773 stc(0) = min(273.16, sfctmp)
2780 if(isnow < 0 .and. newnode == 0 .and. qsnow > 0.)
then
2781 snice(isnow+1) = snice(isnow+1) + qsnow * dt
2782 dzsnso(isnow+1) = dzsnso(isnow+1) + snowhin * dt
2791 snliq ,imelt ,ficeold, & !in
2798 integer,
intent(in) :: nsoil
2799 integer,
intent(in) :: nsnow
2800 integer,
dimension(-nsnow+1:0) ,
intent(in) :: imelt
2801 real (kind=kind_phys),
intent(in) :: dt
2802 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(in) :: stc
2803 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(in) :: snice
2804 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(in) :: snliq
2805 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(in) :: ficeold
2808 integer,
intent(inout) :: isnow
2809 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: dzsnso
2812 real (kind=kind_phys),
parameter :: c2 = 21.e-3
2813 real (kind=kind_phys),
parameter :: c3 = 2.5e-6
2814 real (kind=kind_phys),
parameter :: c4 = 0.04
2815 real (kind=kind_phys),
parameter :: c5 = 2.0
2816 real (kind=kind_phys),
parameter :: dm = 100.0
2817 real (kind=kind_phys),
parameter :: eta0 = 1.8e+6
2819 real (kind=kind_phys) :: burden
2820 real (kind=kind_phys) :: ddz1
2821 real (kind=kind_phys) :: ddz2
2822 real (kind=kind_phys) :: ddz3
2823 real (kind=kind_phys) :: dexpf
2824 real (kind=kind_phys) :: td
2825 real (kind=kind_phys) :: pdzdtc
2826 real (kind=kind_phys) :: void
2827 real (kind=kind_phys) :: wx
2828 real (kind=kind_phys) :: bi
2829 real (kind=kind_phys),
dimension(-nsnow+1:0) :: fice
2838 wx = snice(j) + snliq(j)
2839 fice(j) = snice(j) / wx
2840 void = 1. - (snice(j)/denice + snliq(j)/denh2o) / dzsnso(j)
2843 if (void > 0.001 .and. snice(j) > 0.1)
then
2844 bi = snice(j) / dzsnso(j)
2845 td = max(0.,tfrz-stc(j))
2852 if (bi > dm) ddz1 = ddz1*exp(-46.0e-3*(bi-dm))
2856 if (snliq(j) > 0.01*dzsnso(j)) ddz1=ddz1*c5
2860 ddz2 = -(burden+0.5*wx)*exp(-0.08*td-c2*bi)/eta0
2864 if (imelt(j) == 1)
then
2865 ddz3 = max(0.,(ficeold(j) - fice(j))/max(1.e-6,ficeold(j)))
2873 pdzdtc = (ddz1 + ddz2 + ddz3)*dt
2874 pdzdtc = max(-0.5,pdzdtc)
2878 dzsnso(j) = dzsnso(j)*(1.+pdzdtc)
2879 dzsnso(j) = min(max(dzsnso(j),(snliq(j)+snice(j))/500.0),(snliq(j)+snice(j))/50.0)
2884 burden = burden + wx
2892 isnow ,sh2o ,stc ,snice ,snliq , & !inout
2893 dzsnso ,sice ,snowh ,sneqv , & !inout
2900 integer,
intent(in) :: nsnow
2901 integer,
intent(in) :: nsoil
2905 integer,
intent(inout) :: isnow
2906 real (kind=kind_phys),
dimension( 1:nsoil),
intent(inout) :: sh2o
2907 real (kind=kind_phys),
dimension( 1:nsoil),
intent(inout) :: sice
2908 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: stc
2909 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(inout) :: snice
2910 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(inout) :: snliq
2911 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: dzsnso
2912 real (kind=kind_phys),
intent(inout) :: sneqv
2913 real (kind=kind_phys),
intent(inout) :: snowh
2914 real (kind=kind_phys),
intent(inout) :: ponding1
2915 real (kind=kind_phys),
intent(inout) :: ponding2
2920 integer :: isnow_old
2923 real (kind=kind_phys) :: zwice
2924 real (kind=kind_phys) :: zwliq
2925 real (kind=kind_phys) :: dzmin(3)
2927 data dzmin /0.025, 0.025, 0.1/
2932 do j = isnow_old+1,0
2933 if (snice(j) <= .1)
then
2935 snliq(j+1) = snliq(j+1) + snliq(j)
2936 snice(j+1) = snice(j+1) + snice(j)
2937 dzsnso(j+1) = dzsnso(j+1) + dzsnso(j)
2939 if (isnow_old < -1)
then
2940 snliq(j-1) = snliq(j-1) + snliq(j)
2941 snice(j-1) = snice(j-1) + snice(j)
2942 dzsnso(j-1) = dzsnso(j-1) + dzsnso(j)
2944 if(snice(j) >= 0.)
then
2949 ponding1 = snliq(j) + snice(j)
2950 if(ponding1 < 0.)
then
2951 sice(1) = max(0.0,sice(1)+ponding1/(dzsnso(1)*1000.))
2966 if (j > isnow+1 .and. isnow < -1)
then
2967 do i = j, isnow+2, -1
2969 snliq(i) = snliq(i-1)
2970 snice(i) = snice(i-1)
2971 dzsnso(i)= dzsnso(i-1)
2980 if(sice(1) < 0.)
then
2981 sh2o(1) = sh2o(1) + sice(1)
2985 if(isnow ==0)
return
2993 sneqv = sneqv + snice(j) + snliq(j)
2994 snowh = snowh + dzsnso(j)
2995 zwice = zwice + snice(j)
2996 zwliq = zwliq + snliq(j)
3002 if (snowh < 0.025 .and. isnow < 0 )
then
3006 ponding2 = ponding2 + zwliq
3007 if(sneqv <= 0.) snowh = 0.
3019 if (isnow < -1)
then
3024 do i = isnow_old+1,0
3025 if (dzsnso(i) < dzmin(mssi))
then
3027 if (i == isnow+1)
then
3029 else if (i == 0)
then
3033 if ((dzsnso(i-1)+dzsnso(i)) < (dzsnso(i+1)+dzsnso(i))) neibor = i-1
3037 if (neibor > i)
then
3046 stc(j), dzsnso(l), snliq(l), snice(l), stc(l) )
3049 if (j-1 > isnow+1)
then
3050 do k = j-1, isnow+2, -1
3052 snice(k) = snice(k-1)
3053 snliq(k) = snliq(k-1)
3054 dzsnso(k) = dzsnso(k-1)
3060 if (isnow >= -1)
exit
3084 real (kind=kind_phys),
intent(in) :: dz2
3085 real (kind=kind_phys),
intent(in) :: wliq2
3086 real (kind=kind_phys),
intent(in) :: wice2
3087 real (kind=kind_phys),
intent(in) :: t2
3088 real (kind=kind_phys),
intent(inout) :: dz
3089 real (kind=kind_phys),
intent(inout) :: wliq
3090 real (kind=kind_phys),
intent(inout) :: wice
3091 real (kind=kind_phys),
intent(inout) :: t
3095 real (kind=kind_phys) :: dzc
3096 real (kind=kind_phys) :: wliqc
3097 real (kind=kind_phys) :: wicec
3098 real (kind=kind_phys) :: tc
3099 real (kind=kind_phys) :: h
3100 real (kind=kind_phys) :: h2
3101 real (kind=kind_phys) :: hc
3106 wicec = (wice+wice2)
3107 wliqc = (wliq+wliq2)
3108 h = (cice*wice+cwat*wliq) * (t-tfrz)+hfus*wliq
3109 h2= (cice*wice2+cwat*wliq2) * (t2-tfrz)+hfus*wliq2
3113 tc = tfrz + hc/(cice*wicec + cwat*wliqc)
3114 else if (hc.le.hfus*wliqc)
then
3117 tc = tfrz + (hc - hfus*wliqc) / (cice*wicec + cwat*wliqc)
3129 isnow ,stc ,snice ,snliq ,dzsnso )
3135 integer,
intent(in) :: nsnow
3136 integer,
intent(in) :: nsoil
3140 integer ,
intent(inout) :: isnow
3141 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: stc
3142 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(inout) :: snice
3143 real (kind=kind_phys),
dimension(-nsnow+1: 0),
intent(inout) :: snliq
3144 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: dzsnso
3150 real (kind=kind_phys) :: drr
3151 real (kind=kind_phys),
dimension( 1:nsnow) :: dz
3152 real (kind=kind_phys),
dimension( 1:nsnow) :: swice
3153 real (kind=kind_phys),
dimension( 1:nsnow) :: swliq
3154 real (kind=kind_phys),
dimension( 1:nsnow) :: tsno
3155 real (kind=kind_phys) :: zwice
3156 real (kind=kind_phys) :: zwliq
3157 real (kind=kind_phys) :: propor
3158 real (kind=kind_phys) :: dtdz
3162 if (j <= abs(isnow))
then
3163 dz(j) = dzsnso(j+isnow)
3164 swice(j) = snice(j+isnow)
3165 swliq(j) = snliq(j+isnow)
3166 tsno(j) = stc(j+isnow)
3174 if (dz(1) > 0.05)
then
3177 swice(1) = swice(1)/2.
3178 swliq(1) = swliq(1)/2.
3187 if (dz(1) > 0.05)
then
3190 zwice = propor*swice(1)
3191 zwliq = propor*swliq(1)
3193 swice(1) = propor*swice(1)
3194 swliq(1) = propor*swliq(1)
3197 call combo_glacier (dz(2), swliq(2), swice(2), tsno(2), drr, &
3198 zwliq, zwice, tsno(1))
3201 if (msno <= 2 .and. dz(2) > 0.20)
then
3204 dtdz = (tsno(1) - tsno(2))/((dz(1)+dz(2))/2.)
3206 swice(2) = swice(2)/2.
3207 swliq(2) = swliq(2)/2.
3211 tsno(3) = tsno(2) - dtdz*dz(2)/2.
3212 if (tsno(3) >= tfrz)
then
3215 tsno(2) = tsno(2) + dtdz*dz(2)/2.
3223 if (dz(2) > 0.2)
then
3226 zwice = propor*swice(2)
3227 zwliq = propor*swliq(2)
3229 swice(2) = propor*swice(2)
3230 swliq(2) = propor*swliq(2)
3232 call combo_glacier (dz(3), swliq(3), swice(3), tsno(3), drr, &
3233 zwliq, zwice, tsno(2))
3240 dzsnso(j) = dz(j-isnow)
3241 snice(j) = swice(j-isnow)
3242 snliq(j) = swliq(j-isnow)
3243 stc(j) = tsno(j-isnow)
3256 isnow ,dzsnso ,snowh ,sneqv ,snice , & !inout
3257 snliq ,sh2o ,sice ,stc , & !inout
3258 ponding1 ,ponding2 ,fsh , & !inout
3268 integer,
intent(in) :: nsnow
3269 integer,
intent(in) :: nsoil
3270 real (kind=kind_phys),
intent(in) :: dt
3271 real (kind=kind_phys),
intent(inout) :: qsnfro
3272 real (kind=kind_phys),
intent(inout) :: qsnsub
3273 real (kind=kind_phys),
intent(in) :: qrain
3277 real (kind=kind_phys),
intent(out) :: qsnbot
3281 integer,
intent(inout) :: isnow
3282 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: dzsnso
3283 real (kind=kind_phys),
intent(inout) :: snowh
3284 real (kind=kind_phys),
intent(inout) :: sneqv
3285 real (kind=kind_phys),
dimension(-nsnow+1:0),
intent(inout) :: snice
3286 real (kind=kind_phys),
dimension(-nsnow+1:0),
intent(inout) :: snliq
3287 real (kind=kind_phys),
dimension( 1:nsoil),
intent(inout) :: sh2o
3288 real (kind=kind_phys),
dimension( 1:nsoil),
intent(inout) :: sice
3289 real (kind=kind_phys),
dimension(-nsnow+1:nsoil),
intent(inout) :: stc
3290 real (kind=kind_phys),
intent(inout) :: ponding1
3291 real (kind=kind_phys),
intent(inout) :: ponding2
3292 real (kind=kind_phys),
intent(inout) :: fsh
3297 real (kind=kind_phys) :: qin
3298 real (kind=kind_phys) :: qout
3299 real (kind=kind_phys) :: wgdif
3300 real (kind=kind_phys),
dimension(-nsnow+1:0) :: vol_liq
3301 real (kind=kind_phys),
dimension(-nsnow+1:0) :: vol_ice
3302 real (kind=kind_phys),
dimension(-nsnow+1:0) :: epore
3303 real (kind=kind_phys) :: propor, temp
3308 if(sneqv == 0.)
then
3309 if(opt_gla == 1)
then
3310 sice(1) = sice(1) + (qsnfro-qsnsub)*dt/(dzsnso(1)*1000.)
3311 elseif(opt_gla == 2)
then
3312 fsh = fsh - (qsnfro-qsnsub)*hsub
3323 if(isnow == 0 .and. sneqv > 0.)
then
3324 if(opt_gla == 1)
then
3326 sneqv = sneqv - qsnsub*dt + qsnfro*dt
3328 snowh = max(0.,propor * snowh)
3329 snowh = min(max(snowh,sneqv/500.0),sneqv/50.0)
3330 elseif(opt_gla == 2)
then
3331 fsh = fsh - (qsnfro-qsnsub)*hsub
3337 sice(1) = sice(1) + sneqv/(dzsnso(1)*1000.)
3341 if(sice(1) < 0.)
then
3342 sh2o(1) = sh2o(1) + sice(1)
3347 if(snowh <= 1.e-8 .or. sneqv <= 1.e-6)
then
3354 if ( isnow < 0 )
then
3356 wgdif = snice(isnow+1) - qsnsub*dt + qsnfro*dt
3357 snice(isnow+1) = wgdif
3358 if (wgdif < 1.e-6 .and. isnow <0)
then
3360 isnow ,sh2o ,stc ,snice ,snliq , &
3361 dzsnso ,sice ,snowh ,sneqv , &
3362 ponding1, ponding2 )
3365 if ( isnow < 0 )
then
3366 snliq(isnow+1) = snliq(isnow+1) + qrain * dt
3367 snliq(isnow+1) = max(0., snliq(isnow+1))
3377 if (j >= isnow+1)
then
3378 vol_ice(j) = min(1., snice(j)/(dzsnso(j)*denice))
3379 epore(j) = 1. - vol_ice(j)
3380 vol_liq(j) = min(epore(j),snliq(j)/(dzsnso(j)*denh2o))
3390 if (j >= isnow+1)
then
3391 snliq(j) = snliq(j) + qin
3393 if (epore(j) < 0.05 .or. epore(j+1) < 0.05)
then
3396 qout = max(0.,(vol_liq(j)-ssi*epore(j))*dzsnso(j))
3397 qout = min(qout,(1.-vol_ice(j+1)-vol_liq(j+1))*dzsnso(j+1))
3400 qout = max(0.,(vol_liq(j) - ssi*epore(j))*dzsnso(j))
3403 snliq(j) = snliq(j) - qout
3417 fsh ,fgev ,ssoil ,sag ,prcp ,edir , &
3419 runsrf ,runsub ,sneqv ,dt ,beg_wb, errmsg , errflg )
3421 runsrf ,runsub ,sneqv ,dt ,beg_wb )
3429 integer ,
intent(in) :: iloc
3430 integer ,
intent(in) :: jloc
3431 real (kind=kind_phys) ,
intent(in) :: swdown
3432 real (kind=kind_phys) ,
intent(in) :: fsa
3433 real (kind=kind_phys) ,
intent(in) :: fsr
3434 real (kind=kind_phys) ,
intent(in) :: fira
3435 real (kind=kind_phys) ,
intent(in) :: fsh
3436 real (kind=kind_phys) ,
intent(in) :: fgev
3437 real (kind=kind_phys) ,
intent(in) :: ssoil
3438 real (kind=kind_phys) ,
intent(in) :: sag
3440 real (kind=kind_phys) ,
intent(in) :: prcp
3441 real (kind=kind_phys) ,
intent(in) :: edir
3442 real (kind=kind_phys) ,
intent(in) :: runsrf
3443 real (kind=kind_phys) ,
intent(in) :: runsub
3444 real (kind=kind_phys) ,
intent(in) :: sneqv
3445 real (kind=kind_phys) ,
intent(in) :: dt
3446 real (kind=kind_phys) ,
intent(in) :: beg_wb
3449 character(len=*) ,
intent(inout) :: errmsg
3450 integer ,
intent(inout) :: errflg
3453 real (kind=kind_phys) :: end_wb
3454 real (kind=kind_phys) :: errwat
3455 real (kind=kind_phys) :: erreng
3456 real (kind=kind_phys) :: errsw
3457 character(len=256) :: message
3459 errsw = swdown - (fsa + fsr)
3460 if (errsw > 0.01)
then
3461 write(*,*)
"sag =",sag
3462 write(*,*)
"fsa =",fsa
3463 write(*,*)
"fsr =",fsr
3464 write(message,*)
'errsw =',errsw
3467 errmsg = trim(message)//new_line(
'A')//
"radiation budget problem in noahmp glacier"
3470 call wrf_message(trim(message))
3471 call wrf_error_fatal(
"radiation budget problem in noahmp glacier")
3475 erreng = sag-(fira+fsh+fgev+ssoil)
3476 if(erreng > 0.01)
then
3477 write(message,*)
'erreng =',erreng
3479 errmsg = trim(message)
3481 call wrf_message(trim(message))
3483 write(message,
'(i6,1x,i6,1x,5f10.4)')iloc,jloc,sag,fira,fsh,fgev,ssoil
3486 errmsg = trim(errmsg)//new_line(
'A')//
"energy budget problem in noahmp glacier"
3489 call wrf_message(trim(message))
3490 call wrf_error_fatal(
"energy budget problem in noahmp glacier")
3495 errwat = end_wb-beg_wb-(prcp-edir-runsrf-runsub)*dt
3508 integer,
intent(in) :: iopt_alb
3509 integer,
intent(in) :: iopt_snf
3510 integer,
intent(in) :: iopt_tbot
3511 integer,
intent(in) :: iopt_stc
3513 integer,
intent(in) :: iopt_gla
3514 integer,
intent(in) :: iopt_sfc
3515 integer,
intent(in) :: iopt_trs
3521 opt_tbot = iopt_tbot
3536 use noahmp_glacier_globals
subroutine, public stability(z1, zvfun, gdx, tv1, thv1, wind, z0max, ztmax, tvs, grav, thsfc_loc, rb, fm, fh, fm10, fh2, cm, ch, stress, ustar)
subroutine, private atm_glacier(ep_2, epsm1, sfcprs, sfctmp, q2, soldn, cosz, thair, qair, eair, rhoair, solad, solai, swdown)
re-process atmospheric forcing
subroutine esat(t, esw, esi, desw, desi)
subroutine, private thermoprop_glacier(nsoil, nsnow, isnow, dzsnso, dt, snowh, snice, snliq, df, hcpct, snicev, snliqv, epore, fact)
calculate thermal properties of soil, snow, lake, and frozen soil.
subroutine, private combo_glacier(dz, wliq, wice, t, dz2, wliq2, wice2, t2)
subroutine, private sfcdif1_glacier(iter,zlvl,zpd,z0h,z0m, qair,sfctmp,h,rhoair,mpe,ur, ifdef ccpp
compute surface drag coefficient cm for momentum and ch for heat
subroutine, private rosr12_glacier(p, a, b, c, d, delta, ntop, nsoil, nsnow)
subroutine, private snowh2o_glacier(nsnow, nsoil, dt, qsnfro, qsnsub, qrain, isnow, dzsnso, snowh, sneqv, snice, snliq, sh2o, sice, stc, ponding1, ponding2, fsh, qsnbot)
subroutine, private hrt_glacier(nsnow, nsoil, isnow, zsnso, stc, tbot, zbot, df, hcpct, ssoil, phi, ai, bi, ci, rhsts, botflx)
subroutine, private tsnosoi_glacier(nsoil, nsnow, isnow, dt, tbot, ssoil, snowh, zbot, zsnso, df, hcpct, stc)
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, private compact_glacier(nsnow, nsoil, dt, stc, snice, snliq, imelt, ficeold, isnow, dzsnso)
subroutine, private snowalb_bats_glacier(nband, cosz, fage, albsnd, albsni)
subroutine, private error_glacier(iloc,jloc,swdown,fsa,fsr,fira, fsh,fgev,ssoil,sag,prcp,edir, ifdef ccpp
subroutine, private combine_glacier(nsnow, nsoil, isnow, sh2o, stc, snice, snliq, dzsnso, sice, snowh, sneqv, ponding1, ponding2)
subroutine, private radiation_glacier(dt, tg, sneqvo, sneqv, cosz, qsnow, solad, solai, albold, tauss, sag, fsr, fsa, albsnd, albsni)
Compute solar radiation: absorbed & reflected by the ground.
subroutine, private divide_glacier(nsnow, nsoil, isnow, stc, snice, snliq, dzsnso)
subroutine, private snowwater_glacier(nsnow, nsoil, imelt, dt, sfctmp, snowhin, qsnow, qsnfro, qsnsub, qrain, ficeold, zsoil, isnow, snowh, sneqv, snice, snliq, sh2o, sice, stc, dzsnso, zsnso, fsh, qsnbot, snoflow, ponding1, ponding2)
subroutine, private water_glacier(nsnow, nsoil, imelt, dt, prcp, sfctmp, qvap, qdew, ficeold, zsoil, isnow, snowh, sneqv, snice, snliq, stc, dzsnso, sh2o, sice, ponding, zsnso, fsh, runsrf, runsub, qsnow, ponding1, ponding2, qsnbot, fpice, esnow)
subroutine, private snowfall_glacier(nsoil, nsnow, dt, qsnow, snowhin, sfctmp, isnow, snowh, dzsnso, stc, snice, snliq, sneqv)
subroutine, private phasechange_glacier(nsnow, nsoil, isnow, dt, fact, dzsnso, stc, snice, snliq, sneqv, snowh, smc, sh2o, qmelt, imelt, ponding)
subroutine, public noahmp_glacier( iloc,jloc,cosz,nsnow,nsoil,dt, sfctmp,sfcprs,uu,vv,q2,soldn, prcp,lwdn,tbot,zlvl,ficeold,zsoil, thsfc_loc,prslkix,prsik1x,prslk1x, psfc,pblhx,iz0tlnd,itime, sigmaf1,garea1,psi_opt, ep_1,ep_2,epsm1,cp, qsnow,sneqvo,albold,cm,ch,isnow, sneqv,smc,zsnso,snowh,snice,snliq, tg,stc,sh2o,tauss,qsfc, fsa,fsr,fira,fsh,fgev,ssoil, trad,edir,runsrf,runsub,sag,albedo, qsnbot,ponding,ponding1,ponding2,t2m, q2e,z0h_total, ifdef ccpp
subroutine, private csnow_glacier(isnow, nsnow, nsoil, snice, snliq, dzsnso, tksno, cvsno, snicev, snliqv, epore)
<> snow bulk density, volumetric capacity, and thermal conductivity
subroutine, private energy_glacier(nsnow,nsoil,isnow,dt,qsnow,rhoair, eair,sfcprs,qair,sfctmp,lwdn,uu, vv,solad,solai,cosz,zref, tbot,zbot,zsnso,dzsnso,sigmaf1,garea1, thsfc_loc,prslkix,prsik1x,prslk1x, psfc,pblhx,iz0tlnd,itime,psi_opt, ep_1, ep_2, epsm1, cp, tg,stc,snowh,sneqv,sneqvo,sh2o, smc,snice,snliq,albold,cm,ch, ifdef ccpp
Compute energy budget (momentum & energy fluxes and phase changes).
subroutine albedo(parameters, vegtyp, ist, ice, nsoil, dt, cosz, fage, elai, esai, tg, tv, snowh, fsno, fwet, smc, sneqvo, sneqv, qsnow, fveg, iloc, jloc, albold, tauss, albgrd, albgri, albd, albi, fabd, fabi, ftdd, ftid, ftii, fsun, frevi, frevd, fregd, fregi, bgap, wgap, albsnd, albsni)
surface albedos. also fluxes (per unit incoming direct and diffuse radiation) reflected,...
subroutine, private snowalb_class_glacier(nband, qsnow, dt, alb, albold, albsnd, albsni)
subroutine, private snow_age_glacier(dt, tg, sneqvo, sneqv, tauss, fage)
subroutine, private hstep_glacier(nsnow, nsoil, isnow, dt, ai, bi, ci, rhsts, stc)
subroutine, private glacier_flux(nsoil,nsnow,emg,isnow,df,dzsnso,z0m, zlvl,zpd,qair,sfctmp,rhoair,sfcprs, ur,gamma,rsurf,lwdn,rhsur,smc, eair,stc,sag,snowh,lathea,sh2o, thsfc_loc,prslkix,prsik1x,prslk1x, psfc,pblhx,iz0tlnd,itime,uu,vv, sigmaf1,garea1,psi_opt,ep_1, ep_2, epsm1, cp, ifdef ccpp
use newton-raphson iteration to solve ground (tg) temperature that balances the surface energy budget...
subroutine, public noahmp_options_glacier(iopt_alb, iopt_snf, iopt_tbot, iopt_stc, iopt_gla, iopt_sfc, iopt_trs)
This module contains the interface of noahmp_glacier_routines and noahmp_glacier_globals.
This module contains NoahMP glacier routines.
This module contains the CCPP-compliant GFS surface layer scheme.