-
Notifications
You must be signed in to change notification settings - Fork 51
/
modlsm.f90
2230 lines (1806 loc) · 76.9 KB
/
modlsm.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
!
! Copyright (c) 2020-2020 Wageningen University and Research (WUR)
!
! This file is part of DALES
!
! DALES is free software: you can redistribute it and/or modify
! it under the terms of the GNU General Public License as published by
! the Free Software Foundation, either version 3 of the License, or
! (at your option) any later version.
!
! DALES is distributed in the hope that it will be useful,
! but WITHOUT ANY WARRANTY; without even the implied warranty of
! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
! GNU General Public License for more details.
!
! You should have received a copy of the GNU General Public License
! along with DALES. If not, see <http://www.gnu.org/licenses/>.
!
module modlsm
use netcdf
implicit none
public :: initlsm, lsm, exitlsm, init_lsm_tiles, lags
logical :: llsm ! On/off switch LSM
logical :: lfreedrainage ! Free drainage bottom BC for soil moisture
logical :: lags ! Switch for A-Gs scheme
! Interpolation types soil from full to half level
integer :: iinterp_t, iinterp_theta
integer, parameter :: iinterp_amean = 1 ! val = 0.5*(v1+v2)
integer, parameter :: iinterp_gmean = 2 ! val = sqrt(v1*v2)
integer, parameter :: iinterp_hmean = 3 ! val = ((dz1+dz2)*v1*v2)/(dz1*v2+dz2*v1)
integer, parameter :: iinterp_max = 4 ! val = max(a,b)
! Soil grid
integer :: kmax_soil
real :: z_size_soil
real, allocatable :: z_soil(:), zh_soil(:)
real, allocatable :: dz_soil(:), dzh_soil(:)
real, allocatable :: dzi_soil(:), dzhi_soil(:)
! Soil index in `van_genuchten_parameters.nc` lookup table:
integer, allocatable :: soil_index(:,:,:)
! Source term in soil water budget
real, allocatable :: phiw_source(:,:,:)
! Precipitation, interception et al.
real, allocatable :: throughfall(:,:)
real, allocatable :: interception(:,:)
real, allocatable :: wl_max(:,:)
! Dependency factor canopy resistance on VPD (high veg only)
real, allocatable :: gD(:,:)
! Reduction functions canopy resistance
real, allocatable :: f1(:,:), f2_lv(:,:), f2_hv(:,:), f2b(:,:), f3(:,:)
! Random
real, allocatable :: du_tot(:,:), thv_1(:,:), land_frac(:,:)
! A-Gs
real, allocatable :: an_co2(:,:), resp_co2(:,:)
integer :: co2_index = -1
! Data structure for sub-grid tiles
type lsm_tile
! Static properties:
real, allocatable :: z0m(:,:), z0h(:,:)
! Base tile fraction (i.e. without liquid water)
real, allocatable :: base_frac(:,:)
! Dynamic tile fraction:
real, allocatable :: frac(:,:)
! Monin-obukhov / surface layer:
real, allocatable :: obuk(:,:), ustar(:,:), ra(:,:)
! Conductivity skin layer:
real, allocatable :: lambda_stable(:,:), lambda_unstable(:,:)
! Surface fluxes:
real, allocatable :: H(:,:), LE(:,:), G(:,:)
real, allocatable :: wthl(:,:), wqt(:,:)
! Surface (potential) temperature and humidity:
real, allocatable :: tskin(:,:), thlskin(:,:), qtskin(:,:)
! Buoyancy difference surface - atmosphere
real, allocatable :: db(:,:)
! Vegetation properties:
real, allocatable :: lai(:,:), rs_min(:,:), rs(:,:)
! Root fraction in soil
real, allocatable :: a_r(:,:), b_r(:,:)
real, allocatable :: root_frac(:,:,:)
! Root fraction weighted mean soil water content
real, allocatable :: phiw_mean(:,:)
end type lsm_tile
! Tiles for low veg (lv), high veg (hv), bare soil (bs), wet skin (ws), water (aq):
type(lsm_tile) tile_lv, tile_hv, tile_bs, tile_ws, tile_aq
! Land-surface / van Genuchten parameters from NetCDF input table.
real, allocatable :: &
theta_res(:), theta_wp(:), theta_fc(:), theta_sat(:), &
gamma_theta_sat(:), vg_a(:), vg_l(:), vg_n(:)
! Derived soil parameters
real, allocatable :: &
vg_m(:), &
lambda_theta_min(:), lambda_theta_max(:), &
gamma_theta_min(:), gamma_theta_max(:), &
gamma_t_dry(:), rho_C(:)
contains
subroutine lsm
implicit none
if (.not. llsm) return
! Calculate dynamic tile fractions,
! based on the amount of liquid water on vegetation.
call calc_tile_fractions
! Calculate root fraction weighted mean soil water content.
call calc_theta_mean(tile_lv)
call calc_theta_mean(tile_hv)
! Calculate canopy/soil resistances.
if (lags) then
call calc_canopy_resistance_ags
else
call calc_canopy_resistance_js
endif
! Calculate aerodynamic resistance (and u*, obuk).
call calc_stability
! Set grid point averaged boundary conditions (thls, qts, gradients, ..)
call calc_bulk_bcs
! Calculate soil tendencies
! Calc diffusivity heat:
call calc_thermal_properties
! Solve diffusion equation:
call integrate_t_soil
! Calc diffusivity and conductivity soil moisture:
call calc_hydraulic_properties
! Calculate tendency due to root water extraction
call calc_root_water_extraction
! Solve diffusion equation:
call integrate_theta_soil
! Update liquid water reservoir
call calc_liquid_reservoir
end subroutine lsm
!
! Calculate dynamic tile fractions, based on liquid water on vegetation
!
subroutine calc_tile_fractions
use modglobal, only : i1, j1
use modsurfdata, only : wl, wmax
implicit none
integer :: i, j
real :: c_liq
! Water tile fraction is not variable
tile_aq%frac(:,:) = tile_aq%base_frac(:,:)
do j=2, j1
do i=2, i1
c_liq = min(1., wl(i,j)/wl_max(i,j))
tile_ws%frac(i,j) = land_frac(i,j) * c_liq
tile_lv%frac(i,j) = (1.-c_liq) * tile_lv%base_frac(i,j)
tile_hv%frac(i,j) = (1.-c_liq) * tile_hv%base_frac(i,j)
tile_bs%frac(i,j) = (1.-c_liq) * (1.-tile_lv%base_frac(i,j)-tile_hv%base_frac(i,j)-tile_aq%base_frac(i,j))
end do
end do
end subroutine calc_tile_fractions
!
! Calculate changes in the liquid water reservoir
!
subroutine calc_liquid_reservoir
use modglobal, only : rk3step, rdt, i1, j1, rhow, rlv
use modsurfdata, only : wl, wlm
use modmicrodata, only : imicro, precep
use modfields, only : rhof
implicit none
integer :: i, j
real :: tend, rk3coef, rainrate, c_veg, wl_tend_max, wl_tend_min
real :: wl_tend_liq, wl_tend_dew, wl_tend_precip, wl_tend_sum, wl_tend_lim
real, parameter :: intercept_eff = 0.5
real, parameter :: to_ms = 1./(rhow*rlv)
rk3coef = rdt / (4. - dble(rk3step))
if(rk3step == 1) wlm(:,:) = wl(:,:)
do j=2, j1
do i=2, i1
c_veg = tile_lv%base_frac(i,j)+tile_hv%base_frac(i,j)
! Max and min possible tendencies
wl_tend_min = -wlm(i,j) / rk3coef
wl_tend_max = (wl_max(i,j) - wlm(i,j)) / rk3coef
! Tendency due to evaporation from liquid water reservoir/tile.
wl_tend_liq = -max(0., tile_ws%frac(i,j) * tile_ws%LE(i,j) * to_ms)
! Tendency due to dewfall into vegetation/soil/liquid water tiles
wl_tend_dew = &
-( min(0., tile_lv%frac(i,j) * tile_lv%LE(i,j) * to_ms) + &
min(0., tile_hv%frac(i,j) * tile_hv%LE(i,j) * to_ms) + &
min(0., tile_bs%frac(i,j) * tile_bs%LE(i,j) * to_ms) + &
min(0., tile_ws%frac(i,j) * tile_ws%LE(i,j) * to_ms) )
! Tendency due to interception of precipitation by vegetation
if (imicro == 0) then
rainrate = 0.
else
!rainrate = -sed_qr(i,j,1)/rhow
rainrate = -precep(i,j,1)*rhof(1)/rhow
end if
wl_tend_precip = intercept_eff * c_veg * rainrate
! Total and limited tendencies
wl_tend_sum = wl_tend_liq + wl_tend_dew + wl_tend_precip;
wl_tend_lim = min(wl_tend_max, max(wl_tend_min, wl_tend_sum));
! Diagnose interception and throughfall
throughfall(i,j) = &
-(1.-c_veg) * rainrate &
-(1.-intercept_eff) * c_veg * rainrate &
+ min(0., wl_tend_lim - wl_tend_sum)
interception(i,j) = max(0., wl_tend_lim)
! Integrate
wl(i,j) = wlm(i,j) + rk3coef * wl_tend_lim
end do
end do
end subroutine calc_liquid_reservoir
!
! Calculate root fraction weighted mean soil water content
!
subroutine calc_theta_mean(tile)
use modglobal, only : i1, j1
use modsurfdata, only : phiw
implicit none
type(lsm_tile), intent(inout) :: tile
integer :: i, j, k, si
real :: theta_lim
tile%phiw_mean(:,:) = 0.
do k=1, kmax_soil
do j=2,j1
do i=2,i1
si = soil_index(i,j,k)
theta_lim = max(phiw(i,j,k), theta_wp(si))
tile%phiw_mean(i,j) = tile%phiw_mean(i,j) + tile%root_frac(i,j,k) * &
(theta_lim - theta_wp(si)) / (theta_fc(si) - theta_wp(si))
end do
end do
end do
end subroutine calc_theta_mean
!
! Calculate canopy and soil resistances using Jarvis-Stewart method.
!
subroutine calc_canopy_resistance_js
use modglobal, only : i1, j1
use modfields, only : thl0, qt0, exnf, presf
use modsurface, only : ps
use modraddata, only : swd
use modsurfdata, only : phiw
implicit none
integer :: i, j, k, si
real :: swd_pos, T, esat, e, theta_min, theta_rel, c_veg
! Constants f1 calculation:
real, parameter :: a_f1 = 0.81
real, parameter :: b_f1 = 0.004
real, parameter :: c_f1 = 0.05
k = kmax_soil
do j=2,j1
do i=2,i1
si = soil_index(i,j,k)
! f1: reduction vegetation resistance as f(sw_in):
swd_pos = max(0., -swd(i,j,1))
f1(i,j) = 1./min(1., (b_f1*swd_pos + c_f1) / (a_f1 * (b_f1*swd_pos + 1.)))
! f2: reduction vegetation resistance as f(theta):
f2_lv(i,j) = 1./min(1., max(1.e-9, tile_lv%phiw_mean(i,j)))
f2_hv(i,j) = 1./min(1., max(1.e-9, tile_hv%phiw_mean(i,j)))
! f3: reduction vegetation resistance as f(VPD) (high veg only):
T = thl0(i,j,1) * exnf(1)
esat = 0.611e3 * exp(17.2694 * (T - 273.16) / (T - 35.86))
e = qt0(i,j,1) * presf(1) / 0.622
f3(i,j) = 1./exp(-gD(i,j) * (esat-e))
! f2b: reduction soil resistance as f(theta)
c_veg = tile_lv%base_frac(i,j) + tile_hv%base_frac(i,j)
theta_min = c_veg * theta_wp(si) + (1.-c_veg) * theta_res(si);
theta_rel = (phiw(i,j,k) - theta_min) / (theta_fc(si) - theta_min);
f2b(i,j) = 1./min(1., max(1.e-9, theta_rel))
! Calculate canopy and soil resistance
tile_lv%rs(i,j) = tile_lv%rs_min(i,j) / tile_lv%lai(i,j) * f1(i,j) * f2_lv(i,j)
tile_hv%rs(i,j) = tile_hv%rs_min(i,j) / tile_hv%lai(i,j) * f1(i,j) * f2_hv(i,j) * f3(i,j)
tile_bs%rs(i,j) = tile_bs%rs_min(i,j) / f2b(i,j)
tile_ws%rs(i,j) = 0.
end do
end do
end subroutine calc_canopy_resistance_js
!
! Calculate canopy and soil resistances using A-Gs (plant physiology).
! In addition, this calculates/sets the surface CO2 fluxes...
!
subroutine calc_canopy_resistance_ags
use modglobal, only : i1, j1
use modfields, only : rhof, exnh, qt0, presf, svm
use modsurface, only : E1, ps
use modsurfdata, only : svflux, tskin, phiw, tsoil
use modraddata, only : swd
use modemisdata, only : l_emission, svco2ags, svco2sum, svco2veg
implicit none
! NOTE: these should become a namelist switches...
logical, parameter :: lsplitleaf = .false.
logical, parameter :: lrelaxgc = .false.
real :: Ts, co2_comp, gm, fmin0, fmin, esatsurf, e, Ds, Dmax, cfrac, co2_abs, ci, to_ppb, from_ppb
real :: Ammax, fstr, Am, Rdark, PAR, alphac, AGSa1, Dstar, tempy, An, gc_inf, gcco2, fw, t_mean, th_mean, rs_co2
real :: cveg, cland, theta_min, theta_rel
integer :: i, j, k, si
! Fixed constants (** = same in DALES and IFS, !! = different in DALES and IFS)
real, parameter :: Q10gm = 2.0 ! (**) Parameter to calculate the mesophyll conductance
real, parameter :: Q10am = 2.0 ! (**) Parameter to calculate max primary productivity
real, parameter :: Q10lambda = 1.5 ! (!!) Parameter to calculate the CO2 compensation concentration. (2 in IFS, 1.5 in DALES)
! Reference temperatures calculation mesophyll conductance:
real, parameter :: T1gm = 278 ! (**)
real, parameter :: T2gm = 301 ! (!!: IFS=309, DALES=301 (default), 309 (C4))
! Reference temperatues calculation max primary productivity:
real, parameter :: T1Am = 286 ! (!!: IFS=281, DALES=286 (C4))
real, parameter :: T2Am = 311 ! (**)
real, parameter :: nuco2q = 1.6 ! Ratio molecular viscosity water to carbon dioxide
real, parameter :: gmin = 2.5e-4 ! Cuticular (minimum) conductance. NOTE: = g_cu in IFS?
real, parameter :: ad = 0.07 ! Regression coefficient to calculate Cfrac
real, parameter :: Kx = 0.7 ! Extinction coefficient PAR
!????
real, parameter :: alpha0 = 0.014 ! Initial low light conditions (?)
real, parameter :: Mair = 28.97
real, parameter :: Mco2 = 44
! Parameters respiration Jacobs (2006)
real, parameter :: Cw = 1.6e-3 ! Constant water stress correction
real, parameter :: wsmax = 0.55 ! Upper reference value soil water # NOTE: I guess these are soil dependant?
real, parameter :: wsmin = 0.005 ! Lower reference value soil water # NOTE: see line above :-)
real, parameter :: R10 = 0.23 ! Respiration at 10oC (Jacobs 2007)
real, parameter :: Eact0 = 53.3e3 ! Activation energy
! Vegetation specific constants
! TODO: link to lookup table IFS, for now hardcoded for single veg type.....
real, parameter :: gm298 = 7 !1.3 ! (mm s-1) NOTE: Much lower than DALES...
real, parameter :: Ammax298 = 1.7 ! CO2 maximal primary productivity
real, parameter :: f0 = 0.85 ! Maximum value Cfrac (constant or equation in IFS)
real, parameter :: co2_comp298 = 68.5 ! CO2 compensation concentration = Lambda(25) in IFS
! For sw_splitleaf
!nr_gauss = 3 ! Amount of bins to use for Gaussian integrations
!weight_g = np.array([0.2778,0.4444,0.2778]) ! Weights of the Gaussian bins (must add up to 1)
!angle_g = np.array([0.1127, 0.5,0.8873]) ! Sines of the leaf angles compared to the sun in the first Gaussian integration
!LAI_g = np.array([0.1127, 0.5,0.8873]) ! Ratio of integrated LAI at locations where shaded leaves are evaluated in the second Gaussian integration
!sigma = 0.2 ! Scattering coefficient
!kdfbl = 0.8 ! Diffuse radiation extinction coefficient for black leaves
! Conversion factors
! NOTE: this should use the surface density, not `rhof(1)`.
to_ppb = Mair/Mco2/rhof(1)*1000
from_ppb = 1./to_ppb
do j=2,j1
do i=2,i1
si = soil_index(i, j, kmax_soil)
Ts = tskin(i,j) * exnh(1)
! Calculate the CO2 compensation concentration (IFS eq. 8.92)
! "The compensation point Γ is defined as the CO2 concentration at which the net CO2 assimilation of a fully lit leaf becomes zero."
! NOTE: The old DALES LSM used the atmospheric `thl`, IFS uses the surface temperature.
co2_comp = rhof(1) * co2_comp298 * Q10lambda**(0.1 * (Ts - 298.0))
! Calculate the mesophyll conductance (IFS eq. 8.93)
! "The mesophyll conductance gm describes the transport of CO2 from the substomatal cavities to the mesophyll cells where the carbon is fixed."
! NOTE: The old DALES LSM used the atmospheric `thl`, IFS uses the surface temperature.
gm = gm298 * Q10gm**(0.1 * (Ts - 298.0)) / ((1. + exp(0.3 * (T1gm - Ts))) * (1. + exp(0.3 * (Ts - T2gm)))) / 1000.
! Calculate CO2 concentration inside the leaf (ci)
! NOTE: Differs from IFS
fmin0 = gmin/nuco2q - (1./9.)*gm
fmin = (-fmin0 + (fmin0**2 + 4*gmin/nuco2q*gm)**0.5) / (2.*gm)
! Calculate atmospheric moisture deficit
! NOTE: "Therefore Ci/Cs is specified as a function of atmospheric moisture deficit Ds at the leaf surface"
! In DALES, this uses a mix between esat(surface) and e(atmosphere)
! In IFS, this uses (in kg kg-1): qsat(Ts)-qs instead of qsat(Ts)-qa!
! NOTE: Old DALES LSM used the surface pressure in the calculation of `e`, not sure why...
esatsurf = 0.611e3 * exp(17.2694 * (Ts - 273.16) / (Ts - 35.86))
e = qt0(i,j,1) * presf(1) / 0.622
Ds = (esatsurf - e) / 1000.
! This seems to differ from IFS?
Dmax = (f0 - fmin) / ad
! Coupling factor (IFS eq. 8.101)
!cfrac = f0 * (1.0 - Ds/Dmax) + fmin * (Ds/Dmax)
cfrac = max(0.01, f0 * (1.0 - Ds/Dmax) + fmin * (Ds/Dmax))
! Absolute CO2 concentration.
if (l_emission) then
co2_abs = svm(i,j,1,svco2sum) * from_ppb
else
co2_abs = svm(i,j,1,co2_index) * from_ppb
endif
!if (lrelaxci) then
! if (ci_old_set) then
! ci_inf = cfrac * (co2_abs - co2_comp) + co2_comp
! ci = ci_old(i,j) + min(kci*rk3coef, 1.0) * (ci_inf - ci_old(i,j))
! if (rk3step == 3) then
! ci_old(i,j) = ci
! endif
! else
! ci = cfrac * (co2_abs - co2_comp) + co2_comp
! ci_old(i,j) = ci
! endif
!else
! ci = cfrac * (co2_abs - co2_comp) + co2_comp
!endif
! CO2 concentration in leaf (IFS eq. ???):
ci = cfrac * (co2_abs - co2_comp) + co2_comp
! Max gross primary production in high light conditions (Ag) (IFS eq. 8.94)
! NOTE: The old DALES LSM used the atmospheric `thl`, IFS uses the surface temperature.
Ammax = Ammax298 * Q10am**(0.1 * (Ts - 298.0)) / ((1.0 + exp(0.3 * (T1Am - Ts))) * (1. + exp(0.3 * (Ts - T2Am))))
! Effect of soil moisture stress on gross assimilation rate.
! NOTE: this seems to be different in IFS...
! NOTE: for now, this uses the relative soil moisture content from the low vegetaion only!
fstr = max(1.0e-3, min(1.0, tile_lv%phiw_mean(i,j)))
! Gross assimilation rate (Am, IFS eq. 8.97)
Am = Ammax * (1 - exp(-(gm * (ci - co2_comp) / Ammax)))
! Autotrophic dark respiration (IFS eq. 8.99)
Rdark = Am/9.
!PAR = 0.40 * max(0.1,-swdav * cveg(i,j))
PAR = 0.5 * max(0.1, -swd(i,j,1))
! Light use efficiency
alphac = alpha0 * (co2_abs - co2_comp) / (co2_abs + 2*co2_comp)
if (lsplitleaf) then
print*,'Splitleaf A-Gs not (yet) implemented!'
stop
! PARdir = 0.5 * max(0.1, sw_in[j,i])
! PARdif = 0.5 * max(0.1, sw_in[j,i])
!
! cos_sza = max(1.e-10, calc_zenith_angle(lat, lon, doy, time))
!
! kdrbl = 0.5 / cos_sza ! Direct radiation extinction coefficient for black leaves
! kdf = kdfbl * np.sqrt(1.0-sigma)
! kdr = kdrbl * np.sqrt(1.0-sigma)
! ref = (1.0 - np.sqrt(1.0-sigma)) / (1.0 + np.sqrt(1.0-sigma)) ! Reflection coefficient
! ref_dir = 2 * ref / (1.0 + 1.6 * cos_sza)
!
! Hleaf = np.zeros(nr_gauss+1)
! Fleaf = np.zeros(nr_gauss+1)
! Agl = np.zeros(nr_gauss+1)
!
! Fnet = np.zeros(nr_gauss)
! gnet = np.zeros(nr_gauss)
!
! ! Loop over different LAI locations
! for it in range(nr_gauss):
!
! iLAI = LAI * LAI_g[it] ! Integrated LAI between here and canopy top; Gaussian distributed
! fSL = np.exp(-kdrbl * iLAI) ! Fraction of sun-lit leaves
!
! PARdfD = PARdif * (1.0-ref) * np.exp(-kdf * iLAI) ! Total downward PAR due to diffuse radiation at canopy top
! PARdrD = PARdir * (1.0-ref_dir) * np.exp(-kdr * iLAI) ! Total downward PAR due to direct radiation at canopy top
! PARdfU = PARdif * (1.0-ref) * np.exp(-kdf * LAI) * \
! albedo * (1.0-ref) * np.exp(-kdf * (LAI-iLAI)) ! Total upward (reflected) PAR that originates as diffuse radiation
! PARdrU = PARdir * (1.0-ref_dir) * np.exp(-kdr * LAI) * \
! albedo * (1.0-ref) * np.exp(-kdf * (LAI-iLAI)) ! Total upward (reflected) PAR that originates as direct radiation
! PARdfT = PARdfD + PARdfU ! Total PAR due to diffuse radiation at canopy top
! PARdrT = PARdrD + PARdrU ! Total PAR due to direct radiation at canopy top
!
! dirPAR = (1.0-sigma) * PARdir * fSL ! Purely direct PAR (can only be downward)
! difPAR = PARdfT + PARdrT - dirPAR ! Total diffuse radiation
!
! HdfT = kdf * PARdfD + kdf * PARdfU
! HdrT = kdr * PARdrD + kdf * PARdrU
! dirH = kdrbl * dirPAR
! Hshad = HdfT + HdrT - dirH
!
! Hsun = Hshad + angle_g * (1.0-sigma) * kdrbl * PARdir / np.sum(angle_g * weight_g)
!
! Hleaf[0] = Hshad
! Hleaf[1:] = Hsun
!
! Agl = fstr * (Am + Rdark) * (1 - np.exp(-alphac*Hleaf/(Am + Rdark)))
! gleaf = gmin/nuco2q + Agl/(co2_abs-ci)
! Fleaf = Agl - Rdark
!
! Fshad = Fleaf[0]
! Fsun = np.sum(weight_g * Fleaf[1:])
! gshad = gleaf[0]
! gsun = np.sum(weight_g * gleaf[1:])
!
! Fnet[it] = Fsun * fSL + Fshad * (1 - fSL)
! gnet[it] = gsun * fSL + gshad * (1 - fSL)
!
! An = LAI * np.sum(weight_g * Fnet)
! gc_inf = LAI * np.sum(weight_g * gnet)
!
else
! Calculate upscaling from leaf to canopy: net flow CO2 into the plant (An)
! NOTE: this only uses LAI from low vegetation (for now..)
AGSa1 = 1.0 / (1 - f0)
Dstar = Dmax / (AGSa1 * (f0 - fmin))
tempy = alphac * Kx * PAR / (Am + Rdark)
An = (Am + Rdark) * (1 - 1.0 / (Kx * tile_lv%lai(i,j)) * (E1(tempy * exp(-Kx*tile_lv%lai(i,j))) - E1(tempy)))
gc_inf = tile_lv%lai(i,j) * (gmin/nuco2q + AGSa1 * fstr * An / ((co2_abs - co2_comp) * (1 + Ds / Dstar)))
endif
if (lrelaxgc) then
print*,'Relax GC A-Gs not (yet) implemented!'
stop
! if (gc_old_set) then
! gcco2 = gc_old(i,j) + min(kgc*rk3coef, 1.0) * (gc_inf - gc_old(i,j))
! if (rk3step ==3) then
! gc_old(i,j) = gcco2
! endif
! else
! gcco2 = gc_inf
! gc_old(i,j) = gcco2
! endif
else
gcco2 = gc_inf
endif
! Calculate mean t_soil and theta_soil for calculation soil respiration
t_mean = 0
th_mean = 0
k = kmax_soil
do k=1, kmax_soil
t_mean = t_mean + tsoil(i,j,k) * dz_soil(k)
th_mean = th_mean + phiw(i,j,k) * dz_soil(k)
enddo
t_mean = t_mean / (-zh_soil(1))
th_mean = th_mean / (-zh_soil(1))
! Sub-grid fractions of low+high veg, and low+high veg + bare soil
cveg = tile_lv%base_frac(i,j) + tile_hv%base_frac(i,j)
cland = cveg + tile_bs%base_frac(i,j)
! Water stress function:
fw = Cw * wsmax / (th_mean + wsmin)
! Soil resistance, using old/JS approach
theta_min = cveg * theta_wp(si) + (1.-cveg) * theta_res(si)
theta_rel = (phiw(i,j,kmax_soil) - theta_min) / (theta_fc(si) - theta_min);
f2b(i,j) = 1./min(1., max(1.e-9, theta_rel))
! Surface resistances for moisture and carbon dioxide
tile_lv%rs(i,j) = 1.0 / (1.6 * gcco2)
tile_hv%rs(i,j) = 1.0 / (1.6 * gcco2)
tile_bs%rs(i,j) = tile_bs%rs_min(i,j) / f2b(i,j)
tile_ws%rs(i,j) = 0.
rs_co2 = 1. / gcco2
! NOTE:Combined assimilation for low and high veg, using only low veg as vegetation type:
resp_co2(i,j) = R10 * (1.-fw) * exp(Eact0 / (283.15 * 8.314) * (1.0 - 283.15 / t_mean))
! Scale with vegetation fraction, and translate to `ppb m s-1`.
resp_co2(i,j) = cveg * resp_co2(i,j) * to_ppb
! Calculate net assimilation
! NOTE: this uses the aerodynamic resistance from low vegetation...
an_co2(i,j) = -(co2_abs - ci) / (tile_lv%ra(i,j) + rs_co2)
! Scale with vegetation + bare soil fraction, and translate to `ppb m s-1`.
an_co2(i,j) = cland * an_co2(i,j) * to_ppb
! Set CO2 fluxes:
if (l_emission) then
! Total flux into the CO2 field holding the sum of all CO2 concentrations:
svflux(i,j,svco2sum) = resp_co2(i,j) + an_co2(i,j)
! Respiration flux into the respiration specific field:
svflux(i,j,svco2ags) = resp_co2(i,j)
! Net update into the vegetation specific field:
svflux(i,j,svco2veg) = an_co2(i,j)
else
svflux(i,j,co2_index) = resp_co2(i,j) + an_co2(i,j)
endif
end do
end do
end subroutine calc_canopy_resistance_ags
!
! Calculate Obukhov length, ustar, and aerodynamic resistance, for all tiles
!
subroutine calc_stability
use modglobal, only : i1, j1, cu, cv, rv, rd
use modfields, only : u0, v0, thl0, qt0
implicit none
real, parameter :: du_min = 0.1
real :: du, dv
integer :: i, j
! Calculate properties shared by all tiles:
! Absolute wind speed difference, and virtual potential temperature atmosphere
do j=2,j1
do i=2,i1
du = 0.5*(u0(i,j,1) + u0(i+1,j,1)) + cu
dv = 0.5*(v0(i,j,1) + v0(i,j+1,1)) + cv
du_tot(i,j) = max(0.1, sqrt(du**2 + dv**2))
thv_1(i,j) = thl0(i,j,1) * (1.+(rv/rd-1.)*qt0(i,j,1))
end do
end do
call calc_obuk_ustar_ra(tile_lv)
call calc_obuk_ustar_ra(tile_hv)
call calc_obuk_ustar_ra(tile_bs)
call calc_obuk_ustar_ra(tile_ws)
call calc_obuk_ustar_ra(tile_aq)
end subroutine calc_stability
!
! Calculate Obukhov length and ustar, for single tile
!
subroutine calc_obuk_ustar_ra(tile)
use modglobal, only : i1, j1, rd, rv, grav, zf
use modfields, only : u0, v0
implicit none
type(lsm_tile), intent(inout) :: tile
integer :: i, j
real :: thvs
do j=2,j1
do i=2,i1
if (tile%frac(i,j) > 0) then
! Buoyancy difference surface - atmosphere
thvs = tile%thlskin(i,j) * (1.+(rv/rd-1.)*tile%qtskin(i,j))
!tile%db(i,j) = grav/thvs * (thvs - thv_1(i,j))
tile%db(i,j) = grav/thvs * (thv_1(i,j) - thvs)
! Iteratively find Obukhov length
tile%obuk(i,j) = calc_obuk_dirichlet( &
tile%obuk(i,j), du_tot(i,j), tile%db(i,j), zf(1), tile%z0m(i,j), tile%z0h(i,j))
end if
end do
end do
do j=2,j1
do i=2,i1
if (tile%frac(i,j) > 0) then
! Calculate friction velocity and aerodynamic resistance
tile%ustar(i,j) = du_tot(i,j) * fm(zf(1), tile%z0m(i,j), tile%obuk(i,j))
tile%ra(i,j) = 1./(tile%ustar(i,j) * fh(zf(1), tile%z0h(i,j), tile%obuk(i,j)))
end if
end do
end do
end subroutine calc_obuk_ustar_ra
!
! Calculate surface temperature for each tile, and calculate
! surface fluxes (H, LE, G0, wthl, wqt) and values (thlskin, qtskin)
!
subroutine calc_tile_bcs(tile)
use modglobal, only : i1, j1, cp, rlv, boltz
use modfields, only : exnh, exnf, presh, thl0, qt0, rhof
use modraddata, only : swd, swu, lwd, lwu
use modsurfdata, only : ps, tsoil
implicit none
type(lsm_tile), intent(inout) :: tile
integer :: i, j
real :: Ts, thvs, esats, qsats, desatdTs, dqsatdTs, &
rs_lim, fH, fLE, fG, num, denom, Ta, qsat_new, &
rhocp_i, rholv_i, Qnet
rhocp_i = 1. / (rhof(1) * cp)
rholv_i = 1. / (rhof(1) * rlv)
! tmp
!real :: HLEG, lwu_new, Qnet_new
do j=2, j1
do i=2, i1
if (tile%frac(i,j) > 0) then
! Disable canopy resistance in case of dew fall
Ts = tile%thlskin(i,j) * exnh(1)
esats = 0.611e3 * exp(17.2694 * (Ts - 273.16) / (Ts - 35.86))
qsats = 0.622 * esats / ps
if (qsats < qt0(i,j,1)) then
rs_lim = 0.
else
rs_lim = tile%rs(i,j)
end if
! Calculate recuring terms
! NOTE: this should use the surface density, not rhof(1).
! Not sure if rhoh is available somewhere...
fH = rhof(1) * cp / tile%ra(i,j)
fLE = rhof(1) * rlv / (tile%ra(i,j) + rs_lim)
if (tile%db(i,j) > 0) then
fG = tile%lambda_stable(i,j)
else
fG = tile%lambda_unstable(i,j)
end if
! Net radiation; negative sign = net input of radiation at surface
Qnet = swd(i,j,1) + swu(i,j,1) + lwd(i,j,1) + lwu(i,j,1)
! Solve new skin temperature from SEB
desatdTs = esats * (17.2694 / (Ts - 35.86) - 17.2694 * (Ts - 273.16) / (Ts - 35.86)**2.)
dqsatdTs = 0.622 * desatdTs / ps
Ta = thl0(i,j,1) * exnf(1)
num = -(Qnet - lwu(i,j,1) &
- fH * Ta + (qsats - dqsatdTs * Ts - qt0(i,j,1)) * fLE &
- fG * tsoil(i,j,kmax_soil) - 3.*boltz * Ts**4)
denom = (fH + fLE * dqsatdTs + fG + 4.*boltz * Ts**3)
tile%tskin(i,j) = num / denom
! Update qsat with linearised relation, to make sure that the SEB closes
qsat_new = qsats + dqsatdTs * (tile%tskin(i,j) - Ts)
! Calculate energetic surface fluxes
tile%H (i,j) = fH * (tile%tskin(i,j) - Ta)
tile%LE(i,j) = fLE * (qsat_new - qt0(i,j,1))
tile%G (i,j) = fG * (tsoil(i,j,kmax_soil) - tile%tskin(i,j))
! Calculate kinematic surface fluxes
tile%wthl(i,j) = tile%H (i,j) * rhocp_i
tile%wqt (i,j) = tile%LE(i,j) * rholv_i
! Calculate surface values
tile%thlskin(i,j) = thl0(i,j,1) + tile%wthl(i,j) * tile%ra(i,j)
tile%qtskin (i,j) = qt0 (i,j,1) + tile%wqt (i,j) * tile%ra(i,j)
end if
end do
end do
end subroutine calc_tile_bcs
!
! Calculate BCs and fluxes for the water tile
!
subroutine calc_water_bcs
use modglobal, only : i1, j1, cp, rlv
use modfields, only : exnh, thl0, qt0, rhof
use modsurfdata, only : ps
implicit none
integer :: i, j
real :: esats
do j=2, j1
do i=2, i1
if (tile_aq%frac(i,j) > 0) then
! Calculate BCs. `tile_aq%tskin` is fixed in time (for now...)
tile_aq%thlskin(i,j) = tile_aq%tskin(i,j) / exnh(1)
esats = 0.611e3 * exp(17.2694 * (tile_aq%tskin(i,j) - 273.16) / (tile_aq%tskin(i,j) - 35.86))
tile_aq%qtskin(i,j) = 0.622 * esats / ps
! Calculate kinematic fluxes
tile_aq%wthl(i,j) = 1./tile_aq%ra(i,j) * (tile_aq%thlskin(i,j) - thl0(i,j,1))
tile_aq%wqt (i,j) = 1./tile_aq%ra(i,j) * (tile_aq%qtskin (i,j) - qt0 (i,j,1))
! Calculate energetic fluxes
tile_aq%H (i,j) = tile_aq%wthl(i,j) * rhof(1) * cp
tile_aq%LE(i,j) = tile_aq%wqt (i,j) * rhof(1) * rlv
tile_aq%G (i,j) = 0.
end if
end do
end do
end subroutine calc_water_bcs
!
! Set grid point mean boundary conditions, as used by
! the diffusion scheme, thermodynamics, ...
!
subroutine calc_bulk_bcs
use modglobal, only : i1, j1, i2, j2, cp, rlv, fkar, zf, cu, cv, grav, rv, rd, lopenbc,lboundary,lperiodic
use modfields, only : rhof, thl0, qt0, u0, v0, thvh
use modsurface, only : phim, phih
use modmpi, only : excjs
use modopenboundary, only : openboundary_excjs
use modsurfdata, only : &
H, LE, G0, tskin, qskin, thlflux, qtflux, dthldz, dqtdz, &
dudz, dvdz, ustar, obl, cliq, ra, rsveg, rssoil
implicit none
integer :: i, j
real :: rhocp_i, rholv_i, ucu, vcv, cveg, bflux
real, pointer :: ustar_3D(:,:,:)
rhocp_i = 1. / (rhof(1) * cp)
rholv_i = 1. / (rhof(1) * rlv)
! Calculate surface temperature for each tile, and calculate
! surface fluxes (H, LE, G0, wthl, wqt) and values (thlskin, qtskin)
call calc_tile_bcs(tile_lv)
call calc_tile_bcs(tile_hv)
call calc_tile_bcs(tile_bs)
call calc_tile_bcs(tile_ws)
call calc_water_bcs
do j=2,j1
do i=2,i1
! Calc grid point averaged quantities
H(i,j) = tile_lv%frac(i,j) * tile_lv%H(i,j) + &
tile_hv%frac(i,j) * tile_hv%H(i,j) + &
tile_bs%frac(i,j) * tile_bs%H(i,j) + &
tile_ws%frac(i,j) * tile_ws%H(i,j) + &
tile_aq%frac(i,j) * tile_aq%H(i,j)
LE(i,j) = tile_lv%frac(i,j) * tile_lv%LE(i,j) + &
tile_hv%frac(i,j) * tile_hv%LE(i,j) + &
tile_bs%frac(i,j) * tile_bs%LE(i,j) + &
tile_ws%frac(i,j) * tile_ws%LE(i,j) + &
tile_aq%frac(i,j) * tile_aq%LE(i,j)
G0(i,j) = tile_lv%frac(i,j) * tile_lv%G(i,j) + &
tile_hv%frac(i,j) * tile_hv%G(i,j) + &
tile_bs%frac(i,j) * tile_bs%G(i,j) + &
tile_ws%frac(i,j) * tile_ws%G(i,j) + &
tile_aq%frac(i,j) * tile_aq%G(i,j)
ustar(i,j) = tile_lv%frac(i,j) * tile_lv%ustar(i,j) + &
tile_hv%frac(i,j) * tile_hv%ustar(i,j) + &
tile_bs%frac(i,j) * tile_bs%ustar(i,j) + &
tile_ws%frac(i,j) * tile_ws%ustar(i,j) + &
tile_aq%frac(i,j) * tile_aq%ustar(i,j)
tskin(i,j) = tile_lv%frac(i,j) * tile_lv%thlskin(i,j) + &
tile_hv%frac(i,j) * tile_hv%thlskin(i,j) + &
tile_bs%frac(i,j) * tile_bs%thlskin(i,j) + &
tile_ws%frac(i,j) * tile_ws%thlskin(i,j) + &
tile_aq%frac(i,j) * tile_aq%thlskin(i,j)
qskin(i,j) = tile_lv%frac(i,j) * tile_lv%qtskin(i,j) + &
tile_hv%frac(i,j) * tile_hv%qtskin(i,j) + &
tile_bs%frac(i,j) * tile_bs%qtskin(i,j) + &
tile_ws%frac(i,j) * tile_ws%qtskin(i,j) + &
tile_aq%frac(i,j) * tile_aq%qtskin(i,j)
! Kinematic surface fluxes
thlflux(i,j) = H(i,j) * rhocp_i
qtflux (i,j) = LE(i,j) * rholv_i
! Calculate mean Obukhov length from mean fluxes
bflux = grav/thvh(1) * (thlflux(i,j) * (1.-(1.-rv/rd)*qskin(i,j)) - &
(1.-rv/rd)*tskin(i,j)*qtflux(i,j))
obl(i,j) = -ustar(i,j)**3 / (fkar * bflux)
! MO gradients
dthldz(i,j) = -thlflux(i,j) / (fkar * zf(1) * ustar(i,j)) * phih(zf(1)/obl(i,j))
dqtdz (i,j) = -qtflux (i,j) / (fkar * zf(1) * ustar(i,j)) * phih(zf(1)/obl(i,j))
! NOTE: dudz, dvdz are at the grid center (full level), not the velocity locations.
ucu = 0.5*(u0(i,j,1) + u0(i+1,j,1))+cu
vcv = 0.5*(v0(i,j,1) + v0(i,j+1,1))+cv
dudz(i,j) = ustar(i,j) / (fkar * zf(1)) * phim(zf(1)/obl(i,j)) * (ucu/du_tot(i,j))
dvdz(i,j) = ustar(i,j) / (fkar * zf(1)) * phim(zf(1)/obl(i,j)) * (vcv/du_tot(i,j))
! Cyclic BCs where needed.
ustar_3D(1:i2,1:j2,1:1) => ustar
if(lopenbc) then ! Only use periodicity for non-domain boundaries when openboundaries are used
call openboundary_excjs(ustar_3D, 2,i1,2,j1,1,1,1,1, &
& (.not.lboundary(1:4)).or.lperiodic(1:4))
else
!call excjs(ustar,2,i1,2,j1,1,1,1,1)
call excjs(ustar_3D,2,i1,2,j1,1,1,1,1)
endif
! Just for diagnostics (modlsmcrosssection)
cliq(i,j) = tile_ws%frac(i,j) / land_frac(i,j)
! Calculate ra consistent with mean flux and temperature difference:
ra(i,j) = (tskin(i,j) - thl0(i,j,1)) / thlflux(i,j)
cveg = tile_lv%frac(i,j) + tile_hv%frac(i,j)
if (cveg > 0) then
rsveg(i,j) = (tile_lv%frac(i,j) * tile_lv%rs(i,j) + &
tile_hv%frac(i,j) * tile_hv%rs(i,j)) / cveg
else
rsveg(i,j) = 0.
end if
if (tile_bs%frac(i,j) > 0) then
rssoil(i,j) = tile_bs%rs(i,j)
else
rssoil(i,j) = 0.
end if
end do
end do
end subroutine calc_bulk_bcs
!
! Interpolation soil from full to half levels,
! using various interpolation methods
!
subroutine interpolate_soil(fieldh, field, iinterp)
use modglobal, only : i1, j1
implicit none
real, intent(inout) :: fieldh(:,:,:)
real, intent(in) :: field(:,:,:)
integer, intent(in) :: iinterp
integer i, j, k
if (iinterp == iinterp_amean) then
do k=2,kmax_soil
do j=2,j1
do i=2,i1
fieldh(i,j,k) = 0.5*(field(i,j,k-1) + field(i,j,k))
end do
end do
end do
else if (iinterp == iinterp_gmean) then
do k=2,kmax_soil
do j=2,j1
do i=2,i1
fieldh(i,j,k) = sqrt(field(i,j,k-1) * field(i,j,k))
end do
end do
end do
else if (iinterp == iinterp_hmean) then
do k=2,kmax_soil
do j=2,j1
do i=2,i1
fieldh(i,j,k) = ((dz_soil(k-1)+dz_soil(k))*field(i,j,k-1)*field(i,j,k)) / &
(dz_soil(k-1)*field(i,j,k) + dz_soil(k)*field(i,j,k-1))
end do
end do
end do
else if (iinterp == iinterp_max) then
do k=2,kmax_soil
do j=2,j1
do i=2,i1
fieldh(i,j,k) = max(field(i,j,k-1), field(i,j,k))
end do
end do
end do