forked from NGEET/fates
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathTwoStreamMLPEMod.F90
1731 lines (1357 loc) · 71.3 KB
/
TwoStreamMLPEMod.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
Module TwoStreamMLPEMod
! This module holds the routines to calculate two-tream
! radiation scattering of vegetation in
! "M"ultiple "L"ayers with "P"arellel "E"lements "MLPE"
!
! In summary,
! there may numerous canopy layers. In each canopy layer,
! plant media for different functional types are grouped
! so that they inhabit their own exclusive footprint
! within the layer. Within these exclusive functional
! columns, there are further sub-layer discretizations,
! which are organized by top-down integrated vegetation
! area index.
!
! Note that there is a separate allocation and call
! sequence for each broad band. In other words, the
! two_stream_type is instantiated for each broad band.
!
!
!
! Assumptions: band index 1 = visible (vis)
! 2 = near infrared (nir)
! 3 = thermal (not used at the moment)
!
use shr_log_mod , only: errMsg => shr_log_errMsg
use shr_sys_mod , only: shr_sys_abort
use FatesConstantsMod, only : r8 => fates_r8
implicit none
private
real(r8),parameter :: nearzero = 1.e-20_r8
logical, parameter :: debug=.true.
real(r8), parameter :: unset_r8 = 1.e-36_r8
real(r8), parameter :: unset_int = -999
integer, parameter :: twostr_vis = 1 ! Named index of visible shortwave radiation
integer, parameter :: twostr_nir = 2 ! Named index for near infrared shortwave radiation
! Allowable error, as a fraction of total incident for total canopy
! radiation balance checks
real(r8), public, parameter :: rel_err_thresh = 1.e-6_r8
real(r8), public, parameter :: area_err_thresh = rel_err_thresh*0.1_r8
! These are the codes for how the upper boundary is specified, normalized or absolute
integer,public, parameter :: normalized_upper_boundary = 1
integer,public, parameter :: absolute_upper_boundary = 2
integer :: log_unit ! fortran output unit for logging
! These are parameter constants, ie things that are specific to the plant material
! and radiation band. Not all of these need to be used. 2-stream ultimately wants
! optical depth, scattering coefficient and backscatter fractions for diffuse and
! direct light. So there are various ways to get to these parameters, depending
! on the host model's available parameters. The rho,tau,xl and clumping parameters
! are standard elm/clm parameters, and provided as a convenience.
! Snow optical parameter constants for visible (index=1) and NIR (index=2)
real(r8), parameter :: betad_snow(1:2) = (/0.5, 0.5/) ! Diffuse backscatter fraction (CLM50 Tech Man)
real(r8), parameter :: betab_snow(1:2) = (/0.5, 0.5/) ! Beam backscatter fraction (CLM50 Tech Man)
real(r8), parameter :: om_snow(1:2) = (/0.8, 0.4/) ! Scattering coefficient for snow (CLM50 Tech Man)
!real(r8), parameter :: om_snow(1:2) = (/0.85, 0.75/) ! Tarboton 95
! Cap the maximum optical depth. After 30 or so, its
! so close to zero, if the values get too large, then
! it will blow up the exponents and cause math problems
real(r8), parameter :: kb_max = 30._r8
! For air, use a nominal values to prevent div0s
! the key is that vai = 0
real(r8), parameter :: k_air = 0.5_r8
real(r8), parameter :: om_air = 0.5_r8
real(r8), parameter :: beta_air = 0.5_r8
integer, public, parameter :: air_ft = 0
type, public :: rad_params_type
! From the parameter file
real(r8), allocatable :: rhol(:,:) ! leaf material reflectance: (band x pft)
real(r8), allocatable :: rhos(:,:) ! stem material reflectance: (band x pft)
real(r8), allocatable :: taul(:,:) ! leaf material transmittance: (band x pft)
real(r8), allocatable :: taus(:,:) ! stem material transmittance: (band x pft)
real(r8), allocatable :: xl(:) ! leaf/stem orientation (pft)
real(r8), allocatable :: clumping_index(:) ! clumping index 0-1, when
! leaves stick together (pft)
! Derived parameters
real(r8), allocatable :: phi1(:) ! intermediate term for kd and kb
real(r8), allocatable :: phi2(:) ! intermediate term for kd and kb
real(r8), allocatable :: avmu(:) ! average "av" inverse optical depth "mu" per unit leaf and stem area
real(r8), allocatable :: kd_leaf(:) ! Mean optical depth per unit area leaves in diffuse
real(r8), allocatable :: kd_stem(:) ! Mean optical depth per unit area stems in diffuse
real(r8), allocatable :: om_leaf(:,:) ! Leaf scattering coefficient (band x pft)
real(r8), allocatable :: om_stem(:,:) ! Stem scattering coefficient (band x pft)
end type rad_params_type
type(rad_params_type),public :: rad_params
! Information describing the scattering elements
! that is based on "g"eometry, and independent of wavelength
type scelg_type
integer :: pft ! pft index
real(r8) :: area ! m2 col/m2 ground
real(r8) :: lai ! m2 of leaf area / m2 col
real(r8) :: sai ! m2 of stem area / m2 col
real(r8) :: Kb ! Optical depth of beam radiation
real(r8) :: Kb_leaf ! Optical depth of just leaves in beam radiation
real(r8) :: Kd ! Optical depth of diffuse radiation
real(r8) :: area_squeeze ! This is the ratio of the element area to the
! the area of its constituents. Ideally this
! should be 1.0, but if the host model does not
! do a good job of filling up a canopy with 100% space,
! and instead is fractionally more than 100%, we must
! squeeze the area of 1 or more elements to get an exact
! space usage.
end type scelg_type
! Information describing the scattering elemnets that
! is dependent on wavelengths, ie "b"ands (this is allocated for each broad band)
type scelb_type
! Terms used in the final solution, also used for decomposing solution
real(r8) :: Au ! Compound intercept term
real(r8) :: Ad ! Compound intercept term
real(r8) :: B1 ! Compound term w/ lambdas (operates on e^{av})
real(r8) :: B2 ! Compound term w/ lambdas (operates on e^{-av})
real(r8) :: lambda1_diff ! Compount term w/ B for diffuse forcing
real(r8) :: lambda2_diff ! Compound term w/ B for diffuse forcing
real(r8) :: lambda1_beam ! Compount term w/ B for beam forcing
real(r8) :: lambda2_beam ! Compound term w/ B for beam forcing
real(r8) :: a ! Complex term operating on veg area index
real(r8) :: om ! scattering coefficient for media as a whole
real(r8) :: betad ! backscatter fraction of diffuse radiation for media as a whole
real(r8) :: betab ! backscatter fraction of beam radiation for media as a whole
real(r8) :: Rbeam0 ! Normalized downwelling beam radiation at
! top of the element (relative to downwelling atmospheric beam) [-]
end type scelb_type
type band_type
type(scelb_type), pointer :: scelb(:,:) ! array of scattering coefficients (layer, column)
! can be sparse, will only solve indices up to
integer :: ib ! band index, should be consistent with rad_params
real(r8) :: Rbeam_atm ! Downwelling beam radiation from atmosphere [W/m2 ground]
real(r8) :: Rdiff_atm ! Downwelling diffuse radiation from atmosphere [W/m2 ground]
real(r8) :: albedo_grnd_diff ! Ground albedo diffuse
real(r8) :: albedo_grnd_beam ! Ground albedo direct
end type band_type
! This type contains the pre-processed scattering coefficients
! and routines. This is the parent type that holds almost everything
! in the two-stream solver.
! The scelg structure describes the scattering elements, these are values
! that need to be defined by the ecosystem model, somewhat of
! an input to the solver. Since this is a Perfect Plasticity Approximation
! enabled system, we partition the scattering media into "columns" and "layers"
! Layers are canopy layers, think understory, mid-story and upper canopy. Columns
! are divisions of horizontal space, ie literal columns of space. The current
! implementation limits this space to media that has uniform scattering coefficients.
! So there could not be different PFTs in the same column, because they would undoubtedly
! have different joint scattering coefficients at different height levels in
! the column. Therefore, every column is connected with a PFT.
type, public :: twostream_type
type(scelg_type), pointer :: scelg(:,:) ! array of scattering elements (layer, column)
! can be sparse, will only solve indices up to
! n_lyr,n_col(n_lyr). This is for band (wavelength)
! independent information
type(band_type), pointer :: band(:) ! Holds scattering coefficients for each band
! vis,nir,etc (nothing that emits though, no thermal)
integer :: n_bands ! number of bands (allocation size of band(:))
integer :: n_lyr ! number of (vertical) scattering element layers
integer, allocatable :: n_col(:) ! number of (horizontal) scattering element columns per layer
integer :: n_scel ! total number of scattering elements
logical :: force_prep ! Some coefficients are only updated
! when the canopy composition changes, ie
! changes in leaf, stem or snow structure.
! If so, this sets to true, signalling that diffuse
! scattering coefficients should be updated.
! Otherwise, we only updated zenith dependent
! parameters on short sub-daily timesteps
real(r8) :: frac_snow ! Current mean snow-fraction of the canopy
real(r8) :: frac_snow_old ! Previous mean snow-fraction of the canopy
real(r8) :: cosz ! Current cosine of the zenith angle
contains
procedure :: ZenithPrep ! Update coefficients as zenith changes
procedure :: CanopyPrep ! Update coefficients as canopy changes
procedure :: Solve ! Perform the scattering solution
procedure :: Dump ! Dump out (print out) the site of interest
procedure :: GetNSCel
procedure :: AllocInitTwoStream
procedure :: DeallocTwoStream
procedure :: GetRdUp
procedure :: GetRdDn
procedure :: GetRb
procedure :: GetAbsRad
end type twostream_type
public :: ParamPrep
public :: AllocateRadParams
public :: TwoStreamLogInit
character(len=*), parameter, private :: sourcefile = &
__FILE__
contains
subroutine TwoStreamLogInit(log_unit_in)
integer,intent(in) :: log_unit_in
log_unit = log_unit_in
end subroutine TwoStreamLogInit
subroutine endrun(msg)
!-----------------------------------------------------------------------
! !DESCRIPTION:
! Abort the model for abnormal termination
! This subroutine was derived from CLM's
! endrun_vanilla() in abortutils.F90
!
!
! !ARGUMENTS:
implicit none
character(len=*), intent(in) :: msg ! string to be printed
!-----------------------------------------------------------------------
write(log_unit,*)'ENDRUN:', msg
call shr_sys_abort()
end subroutine endrun
! ===============================================================================================
subroutine AllocInitTwoStream(this,band_indices,ncan,ncol)
class(twostream_type) :: this
integer :: band_indices(:)
integer :: ncan
integer :: ncol
integer :: nbands
integer :: ib
nbands = ubound(band_indices,1)
allocate(this%n_col(ncan))
allocate(this%scelg(ncan,ncol))
allocate(this%band(nbands))
this%n_col(1:ncan) = unset_int
this%n_bands = nbands
this%n_lyr = ncan
this%frac_snow = unset_r8
this%frac_snow_old = unset_r8
do ib = 1,nbands
allocate(this%band(ib)%scelb(ncan,ncol))
this%band(ib)%albedo_grnd_diff = unset_r8
this%band(ib)%albedo_grnd_beam = unset_r8
this%band(ib)%ib = band_indices(ib)
end do
return
end subroutine AllocInitTwoStream
! ===============================================================================================
subroutine DeallocTwoStream(this)
class(twostream_type) :: this
integer :: nbands
integer :: ib
nbands = ubound(this%band,1)
deallocate(this%scelg)
deallocate(this%n_col)
do ib = 1,nbands
deallocate(this%band(ib)%scelb)
end do
deallocate(this%band)
return
end subroutine DeallocTwoStream
! ===============================================================================================
subroutine AllocateRadParams(n_pft,n_bands)
integer,intent(in) :: n_pft
integer,intent(in) :: n_bands
! Include the zeroth pft index for air
allocate(rad_params%rhol(n_bands,n_pft))
allocate(rad_params%rhos(n_bands,n_pft))
allocate(rad_params%taul(n_bands,n_pft))
allocate(rad_params%taus(n_bands,n_pft))
allocate(rad_params%xl(n_pft))
allocate(rad_params%clumping_index(n_pft))
allocate(rad_params%phi1(n_pft))
allocate(rad_params%phi2(n_pft))
allocate(rad_params%avmu(n_pft))
allocate(rad_params%kd_leaf(n_pft))
allocate(rad_params%kd_stem(n_pft))
allocate(rad_params%om_leaf(n_bands,n_pft))
allocate(rad_params%om_stem(n_bands,n_pft))
end subroutine AllocateRadParams
! ================================================================================================
function GetRdDn(this,ican,icol,ib,vai) result(r_diff_dn)
class(twostream_type) :: this
real(r8),intent(in) :: vai
integer,intent(in) :: ican
integer,intent(in) :: icol
integer,intent(in) :: ib
real(r8) :: r_diff_dn
! Rdn = Ad e−(Kbv) + Re + λ1 B2 e^(av) + λ2 B1 e^(−av)
associate(scelb => this%band(ib)%scelb(ican,icol), &
scelg => this%scelg(ican,icol) )
r_diff_dn = this%band(ib)%Rbeam_atm*( &
scelb%Ad*exp(-scelg%Kb*vai) + &
scelb%B2*scelb%lambda1_beam*exp(scelb%a*vai) + &
scelb%B1*scelb%lambda2_beam*exp(-scelb%a*vai)) + &
this%band(ib)%Rdiff_atm*( &
scelb%B2*scelb%lambda1_diff*exp(scelb%a*vai) + &
scelb%B1*scelb%lambda2_diff*exp(-scelb%a*vai))
if(debug)then
if(isnan(r_diff_dn))then
write(log_unit,*)"GETRDN"
write(log_unit,*)scelg%Kb
write(log_unit,*)scelb%a
write(log_unit,*)vai
write(log_unit,*)scelb%Ad
write(log_unit,*)scelb%B1,scelb%B2
write(log_unit,*)scelb%lambda1_beam,scelb%lambda2_beam
write(log_unit,*)scelb%lambda1_diff,scelb%lambda2_diff
write(log_unit,*)this%band(ib)%Rbeam_atm
write(log_unit,*)this%band(ib)%Rdiff_atm
write(log_unit,*)exp(-scelg%Kb*vai)
write(log_unit,*)exp(scelb%a*vai)
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
end if
end associate
end function GetRdDn
function GetRdUp(this,ican,icol,ib,vai) result(r_diff_up)
class(twostream_type) :: this
real(r8),intent(in) :: vai
integer,intent(in) :: ican
integer,intent(in) :: icol
integer,intent(in) :: ib
real(r8) :: r_diff_up
! Rup = Au e−(Kbv) + Re + λ1 B1 e^(av) + λ2 B2 e^(−av)
associate(scelb => this%band(ib)%scelb(ican,icol), &
scelg => this%scelg(ican,icol) )
r_diff_up = this%band(ib)%Rbeam_atm*( &
scelb%Au*exp(-scelg%Kb*vai) + &
scelb%B1*scelb%lambda1_beam*exp(scelb%a*vai) + &
scelb%B2*scelb%lambda2_beam*exp(-scelb%a*vai)) + &
this%band(ib)%Rdiff_atm*( &
scelb%B1*scelb%lambda1_diff*exp(scelb%a*vai) + &
scelb%B2*scelb%lambda2_diff*exp(-scelb%a*vai))
end associate
end function GetRdUp
function GetRb(this,ican,icol,ib,vai) result(r_beam_dn)
class(twostream_type) :: this
real(r8),intent(in) :: vai
integer,intent(in) :: ican
integer,intent(in) :: icol
integer,intent(in) :: ib
real(r8) :: r_beam_dn
r_beam_dn = this%band(ib)%Rbeam_atm * &
this%band(ib)%scelb(ican,icol)%Rbeam0*exp(-this%scelg(ican,icol)%Kb*vai)
end function GetRb
subroutine GetAbsRad(this,ican,icol,ib,vai_top,vai_bot, &
Rb_abs,Rd_abs,Rd_abs_leaf,Rb_abs_leaf,R_abs_stem,R_abs_snow,leaf_sun_frac)
! This routine is used to help decompose radiation scattering
! and return the amount of absorbed radiation. The canopy layer and column
! index identify the element of interest. The other arguments are the upper and
! lower bounds within the element over which to evaluate absorbed radiation.
! The assumption is that the vegetation area index is zero at the top of the
! element, and increases going downwards. As with all assumptions in this
! module, the scattering parameters are uniform within the element itself,
! which includes an assumption of the leaf/stem proportionality.
! ---------------------------------------------------------------------------
! Solution for radiative intensity of diffuse up and down at tai=v
! Rup = Au e−(Kbv) + Re + λ1 B1 e^(av) + λ2 B2 e^(−av)
! Rdn = Ad e−(Kbv) + Re + λ1 B2 e^(av) + λ2 B1 e^(−av)
! ---------------------------------------------------------------------------
! Arguments
class(twostream_type) :: this
integer,intent(in) :: ican
integer,intent(in) :: icol
integer, intent(in) :: ib ! broad band index
real(r8), intent(in) :: vai_top ! veg area index (from the top of element) to start
real(r8), intent(in) :: vai_bot ! veg area index (from the top of element) to finish
real(r8), intent(out) :: Rb_abs ! total absorbed beam radiation [W/m2 ground]
real(r8), intent(out) :: Rd_abs ! total absorbed diffuse radiation [W/m2 ground]
real(r8), intent(out) :: Rb_abs_leaf ! Absorbed beam radiation from leaves [W/m2 ground]
real(r8), intent(out) :: Rd_abs_leaf ! Absorbed diff radiation from leaves [W/m2 ground]
real(r8), intent(out) :: R_abs_stem ! Absorbed beam+diff radiation stems [W/m2 ground]
real(r8), intent(out) :: R_abs_snow ! Absorbed beam+diff radiation snow [W/m2 ground]
real(r8), intent(out) :: leaf_sun_frac ! Fraction of leaves in the interval exposed
! to sunlight
real(r8) :: dvai,dlai ! Amount of VAI and LAI in this interval [m2/m2]
real(r8) :: Rd_net ! Difference in diffuse radiation at upper and lower boundaries [W/m2]
real(r8) :: Rb_net ! Difference in beam radiation at upper and lower boundaries [W/m2]
real(r8) :: vai_max ! total integrated (leaf+stem) area index of the current element
real(r8) :: frac_abs_snow ! fraction of radiation absorbed by snow
real(r8) :: diff_wt_leaf ! diffuse absorption weighting for leaves
real(r8) :: diff_wt_stem ! diffuse absorption weighting for stems
real(r8) :: beam_wt_leaf ! beam absorption weighting for leaves
real(r8) :: beam_wt_stem ! beam absorption weighting for stems
real(r8) :: lai_bot,lai_top
associate(scelb => this%band(ib)%scelb(ican,icol), &
scelg => this%scelg(ican,icol), &
ft => this%scelg(ican,icol)%pft )
! If this is air, trivial solutions
if(ft==air_ft) then
Rb_abs = 0._r8
Rd_abs = 0._r8
Rb_abs_leaf = 0._r8
Rd_abs_leaf = 0._r8
R_abs_stem = 0._r8
R_abs_snow = 0._r8
leaf_sun_frac = 0._r8
return
end if
! The total vegetation area index of the element
vai_max = scelg%lai + scelg%sai
dvai = vai_bot - vai_top
lai_top = vai_top*scelg%lai/( scelg%lai+ scelg%sai)
lai_bot = vai_bot*scelg%lai/( scelg%lai+ scelg%sai)
dlai = dvai * scelg%lai/( scelg%lai+ scelg%sai)
if(dlai>nearzero)then
leaf_sun_frac = max(0.001_r8,min(0.999_r8,scelb%Rbeam0/(dlai*scelg%Kb_leaf/rad_params%clumping_index(ft)) &
*(exp(-scelg%Kb_leaf*lai_top) - exp(-scelg%Kb_leaf*lai_bot))))
else
leaf_sun_frac = 0001._r8
end if
!leaf_sun_frac = max(0.001_r8,min(0.999_r8,scelb%Rbeam0/(dvai*scelg%Kb/rad_params%clumping_index(ft)) &
! *(exp(-scelg%Kb*vai_top) - exp(-scelg%Kb*vai_bot))))
if(debug) then
if(leaf_sun_frac>1.0_r8 .or. leaf_sun_frac<0._r8) then
write(log_unit,*)"impossible leaf sun fraction"
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
end if
! We have to disentangle the absorption between leaves and stems, we give them both
! a weighting fraction of total absorption of area*K*(1-om)
frac_abs_snow = this%frac_snow*(1._r8-om_snow(ib)) / (1._r8-scelb%om)
diff_wt_leaf = scelg%lai*(1._r8-rad_params%om_leaf(ib,ft))*rad_params%Kd_leaf(ft)
diff_wt_stem = scelg%sai*(1._r8-rad_params%om_stem(ib,ft))*rad_params%Kd_stem(ft)
beam_wt_leaf = scelg%lai*(1._r8-rad_params%om_leaf(ib,ft))*scelg%Kb_leaf
beam_wt_stem = scelg%sai*(1._r8-rad_params%om_stem(ib,ft))*1._r8
! Mean element transmission coefficients adding snow scattering
if(debug) then
if( (vai_bot-vai_max)>rel_err_thresh)then
write(log_unit,*)"During decomposition of the 2-stream radiation solution"
write(log_unit,*)"A vegetation area index (VAI) was requested in GetAbsRad()"
write(log_unit,*)"that is larger than the total integrated VAI of the "
write(log_unit,*)"computation element of interest."
write(log_unit,*)"vai_max: ",vai_max
write(log_unit,*)"vai_bot: ",vai_bot
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
if( (vai_bot-vai_top)<-rel_err_thresh ) then
write(log_unit,*)"During decomposition of the 2-stream radiation solution"
write(log_unit,*)"the vegetation area index at the lower position was set"
write(log_unit,*)"as greater than the value at the upper position."
write(log_unit,*)"vai_max: ",vai_max
write(log_unit,*)"vai_bot: ",vai_bot
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
end if
! Amount of absorbed radiation is retrieved by doing an energy
! balance on this boundaries over the depth of interest (ie net)
! Result is Watts / m2 of the element's area footprint NOT
! per m2 of tissue (at least not in this step)
Rb_net = this%GetRb(ican,icol,ib,vai_top)-this%GetRb(ican,icol,ib,vai_bot)
Rd_net = (this%GetRdDn(ican,icol,ib,vai_top) - this%GetRdDn(ican,icol,ib,vai_bot)) + &
(this%GetRdUp(ican,icol,ib,vai_bot) - this%GetRdUp(ican,icol,ib,vai_top))
! The net beam radiation includes that which is absorbed, but also,
! that which is re-scattered, the re-scattered acts as a source
! to the net diffuse balance and adds to the absorbed, and a sink
! on the beam absorbed term.
Rb_abs = Rb_net * (1._r8-this%band(ib)%scelb(ican,icol)%om)
Rd_abs = Rd_net + Rb_net * this%band(ib)%scelb(ican,icol)%om
Rb_abs_leaf = (1._r8-frac_abs_snow)*Rb_abs * beam_wt_leaf / (beam_wt_leaf+beam_wt_stem)
Rd_abs_leaf = (1._r8-frac_abs_snow)*Rd_abs * diff_wt_leaf / (diff_wt_leaf+diff_wt_stem)
R_abs_snow = (Rb_abs+Rd_abs)*frac_abs_snow
R_abs_stem = (1._r8-frac_abs_snow)* &
(Rb_abs*beam_wt_stem / (beam_wt_leaf+beam_wt_stem) + &
Rd_abs*diff_wt_stem / (diff_wt_leaf+diff_wt_stem))
end associate
return
end subroutine GetAbsRad
! ================================================================================================
subroutine Dump(this,ib,lat,lon)
! Dump out everything we know about these two-stream elements
class(twostream_type) :: this
integer,intent(in) :: ib
real(r8),optional,intent(in) :: lat
real(r8),optional,intent(in) :: lon
integer :: ican
integer :: icol
write(log_unit,*) 'Dumping Two-stream elements for band ', ib
write(log_unit,*)
write(log_unit,*) 'rbeam atm: ',this%band(ib)%Rbeam_atm
write(log_unit,*) 'rdiff_atm: ',this%band(ib)%Rdiff_atm
write(log_unit,*) 'alb grnd diff: ',this%band(ib)%albedo_grnd_diff
write(log_unit,*) 'alb grnd beam: ',this%band(ib)%albedo_grnd_beam
write(log_unit,*) 'cosz: ',this%cosz
write(log_unit,*) 'snow fraction: ',this%frac_snow
if(present(lat)) write(log_unit,*) 'lat: ',lat
if(present(lon)) write(log_unit,*) 'lon: ',lon
do_can: do ican = 1,this%n_lyr
do_col: do icol = 1,this%n_col(ican)
associate(scelg => this%scelg(ican,icol), &
scelb => this%band(ib)%scelb(ican,icol))
write(log_unit,*) '--',ican,icol,'--'
write(log_unit,*) 'pft:',scelg%pft
write(log_unit,*) 'area: ',scelg%area
write(log_unit,*) 'lai,sai: ',scelg%lai,scelg%sai
write(log_unit,*) 'Kb: ',scelg%Kb
write(log_unit,*) 'Kb leaf: ',scelg%Kb_leaf
write(log_unit,*) 'Kd: ',scelg%Kd
write(log_unit,*) 'Rb0: ',scelb%Rbeam0
write(log_unit,*) 'om: ',scelb%om
write(log_unit,*) 'betad: ',scelb%betad
write(log_unit,*) 'betab:',scelb%betab
write(log_unit,*) 'a: ',scelb%a
this%band(ib)%Rbeam_atm = 1.0_r8
this%band(ib)%Rdiff_atm = 1.0_r8
write(log_unit,*)'RDiff Down @ bottom: ',this%GetRdDn(ican,icol,ib,scelg%lai+scelg%sai)
write(log_unit,*)'RDiff Up @ bottom: ',this%GetRdUp(ican,icol,ib,scelg%lai+scelg%sai)
write(log_unit,*)'Rbeam @ bottom: ',this%GetRb(ican,icol,ib,scelg%lai+scelg%sai)
end associate
end do do_col
end do do_can
end subroutine Dump
! ================================================================================================
subroutine ParamPrep()
integer :: ft
integer :: nbands
integer :: numpft
integer :: ib
numpft = ubound(rad_params%om_leaf,2)
nbands = ubound(rad_params%om_leaf,1)
do ft = 1,numpft
! The non-band specific parameters here will be re-derived for each
! band, which is inefficient, however this is an incredibly cheap
! routine to begin with, its only called during initialization, so
! just let it go, dont worry about it.
if(rad_params%xl(ft)<-0.4_r8 .or. rad_params%xl(ft)>0.6_r8) then
write(log_unit,*)"Leaf orientation factors (xl) should be between -0.4 and 0.6"
write(log_unit,*)"ft: ",ft,"xl: ",rad_params%xl(ft)
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
! There is a singularity of leaf orientation is exactly 0
! phi1 = 0.5
! phi2 = 0.0
! avmu = 1/0 (1 - 0.5/0 * ln(0.5/0.5) ) but the limit approaches 1
! a value of 0.0001 does not break numerics and generates an avmu of nearly 1
if( abs(rad_params%xl(ft)) <0.0001) rad_params%xl(ft)=0.0001_r8
! There must be protections on xl to prevent div0 and other weirdness
rad_params%phi1(ft) = 0.5_r8 - 0.633_r8*rad_params%xl(ft) - 0.330_r8*rad_params%xl(ft)*rad_params%xl(ft)
rad_params%phi2(ft) = 0.877_r8 * (1._r8 - 2._r8*rad_params%phi1(ft)) !0 = horiz leaves, 1 - vert leaves.
! Eq. 3.4 CLM50 Tech Man
rad_params%avmu(ft) = (1._r8/rad_params%phi2(ft))* &
(1._r8-(rad_params%phi1(ft)/rad_params%phi2(ft))* &
log((rad_params%phi2(ft)+rad_params%phi1(ft))/rad_params%phi1(ft)))
do ib = 1, nbands
rad_params%Kd_leaf(ft) = rad_params%clumping_index(ft)/rad_params%avmu(ft)
rad_params%Kd_stem(ft) = 1._r8
rad_params%om_leaf(ib,ft) = rad_params%rhol(ib,ft) + rad_params%taul(ib,ft)
rad_params%om_stem(ib,ft) = rad_params%rhos(ib,ft) + rad_params%taus(ib,ft)
if( rad_params%om_leaf(ib,ft) > 0.99_r8 ) then
write(log_unit,*) "In: TwoStreamMLPEMod.F90:ParamPrep()"
write(log_unit,*) "An extremely high leaf scattering coefficient was generated:"
write(log_unit,*) "om = tau + rho"
write(log_unit,*) "band = ",ib
write(log_unit,*) "pft = ",ft
write(log_unit,*) "om_leaf = ",rad_params%om_leaf(ib,ft)
write(log_unit,*) "rhol = ",rad_params%rhol(ib,ft)
write(log_unit,*) "taul = ",rad_params%taul(ib,ft)
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
if( rad_params%om_stem(ib,ft) > 0.99_r8 ) then
write(log_unit,*) "In: TwoStreamMLPEMod.F90:ParamPrep()"
write(log_unit,*) "An extremely high stem scattering coefficient was generated:"
write(log_unit,*) "om = tau + rho"
write(log_unit,*) "band = ",ib
write(log_unit,*) "pft = ",ft
write(log_unit,*) "om_stem = ",rad_params%om_stem(ib,ft)
write(log_unit,*) "rhos = ",rad_params%rhos(ib,ft)
write(log_unit,*) "taus = ",rad_params%taus(ib,ft)
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
end do
end do
return
end subroutine ParamPrep
! ================================================================================================
! ================================================================================================
subroutine CanopyPrep(this,frac_snow)
! Pre-process things that change with canopy-geometry or snow cover
! We try to only run this when necessary. For instance we only
! run this when the canopy vegetation composition changes, or
! when the amount of snow-cover changes.
class(twostream_type) :: this
real(r8) :: frac_snow ! The fraction (in terms of vegetation area index)
! of vegetation covered with snow
! But we check if the snow conditions
! change during the high frequency calls
! as well.
integer :: ib ! The band of interest
integer :: ican ! scattering element canopy layer index (top down)
integer :: icol ! scattering element column
real(r8) :: rho ! element mean material reflectance
real(r8) :: tau ! element mean material transmittance
real(r8) :: vai ! vegetation area index lai+sai
real(r8) :: om_veg ! scattering coefficient for vegetation (no snow)
real(r8) :: betad_veg ! diffuse backscatter for vegetation (no snow)
real(r8) :: betad_om ! multiplication of diffuse backscatter and reflectance
real(r8) :: area_check ! Checks to make sure each layer has 100% coverage
real(r8) :: a2 ! The "a" term squared
this%frac_snow = frac_snow
if(.not.this%force_prep) then
if(abs(this%frac_snow-this%frac_snow_old)<nearzero) then
this%frac_snow_old = this%frac_snow
return
end if
end if
this%frac_snow_old = this%frac_snow
do_can: do ican = 1,this%n_lyr
area_check = 0._r8
do_col: do icol = 1,this%n_col(ican)
associate(lai => this%scelg(ican,icol)%lai, &
sai => this%scelg(ican,icol)%sai, &
ft => this%scelg(ican,icol)%pft, &
scelg => this%scelg(ican,icol))
vai = lai + sai
! Mean element transmission coefficients w/o snow effects
if(ft==air_ft) then
scelg%Kd = k_air
else
if(debug)then
if(vai<nearzero)then
write(log_unit,*)"zero vai in non-air element?"
write(log_unit,*) lai,sai,ican,icol,this%n_lyr,this%n_col(ican)
write(log_unit,*)"TwoStreamMLPEMod.F90:CanopyPrep"
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
end if
scelg%Kd = (lai * rad_params%Kd_leaf(ft) + &
sai * rad_params%Kd_stem(ft))/vai
end if
area_check = area_check + scelg%area
do_bands: do ib = 1, this%n_bands
associate(scelb => this%band(ib)%scelb(ican,icol))
if (ft==air_ft) then
scelb%om = om_air
scelb%betad = beta_air
else
! Material reflectance (weighted average of leaf stem and snow)
! Eq. 3.11 and 3.12 ClM5.0 Tech Man
om_veg = (lai*rad_params%om_leaf(ib,ft) + &
sai*rad_params%om_stem(ib,ft))/vai
! Eq. 3.5 ClM5.0 Tech Man
scelb%om = this%frac_snow*om_snow(ib) + (1._r8-this%frac_snow)*om_veg
! Diffuse backscatter, taken from G. Bonan's code
rho = (lai * rad_params%rhol(ib,ft) + &
sai * rad_params%rhos(ib,ft))/vai
tau = (lai * rad_params%taul(ib,ft) + &
sai * rad_params%taus(ib,ft))/vai
! Eq 3.13 from CLM5.0 Tech Man
betad_veg = 0.5_r8 / scelb%om * &
( scelb%om + (rho-tau) * ((1._r8+rad_params%xl(ft))/2._r8)**2._r8 )
! Eq. 3.6 from CLM5.0 Tech Man
betad_om = betad_veg*om_veg*(1._r8-this%frac_snow) + &
om_snow(ib)*betad_snow(ib)*this%frac_snow
scelb%betad = betad_om / scelb%om
if(debug)then
if(isnan(scelb%betad))then
write(log_unit,*)"nans in canopy prep"
write(log_unit,*) ib,ican,icol,ft
write(log_unit,*) scelb%betad,scelb%om,lai,sai
write(log_unit,*) this%frac_snow,om_snow(ib),vai,om_veg
write(log_unit,*)"TwoStreamMLPEMod.F90:CanopyPrep"
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
end if
end if
a2 = scelg%Kd*scelg%Kd*(1._r8-scelb%om)*(1._r8-scelb%om+2._r8*scelb%om*scelb%betad)
if(a2<0._r8) then
write(log_unit,*)'a^2 is less than zero'
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
! We also have to avoid singularities, see Ad and Au below,
! where a^2-Kb^2 is in the denominator
scelb%a = sqrt(a2)
end associate
end do do_bands
end associate
end do do_col
! RE-ENABLE THIS CHECK WHEN FATES IS BETTER AT CONSERVING AREA!!
if(.false.)then
!if( abs(area_check-1._r8) > 10._r8*area_err_thresh )then
write(log_unit,*)"Only a partial canopy was specified"
write(log_unit,*)"Scattering elements must constitute 100% of the ground cover."
write(log_unit,*)"for open spaces, create an air element with the respective area."
write(log_unit,*)"total area (out of 1): ",area_check,ican
write(log_unit,*)"layer: ",ican," of: ",this%n_lyr
do icol = 1,this%n_col(ican)
write(log_unit,*)this%scelg(ican,icol)%area,this%scelg(ican,icol)%pft
end do
write(log_unit,*)"TwoStreamMLPEMod.F90:CanopyPrep"
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
end do do_can
return
end subroutine CanopyPrep
! ================================================================================================
subroutine ZenithPrep(this,cosz)
! Pre-process things that change with the zenith angle
! i.e. the beam optical properties
! Important !!!!
! This should always be called after CanopyPrep() has been
! called. This routine relies on the results of that routine
! notably the scattering coefficient "om".
class(twostream_type) :: this
integer :: ib ! band index, matches indexing of rad_params
real(r8) :: cosz ! Cosine of the zenith angle
integer :: ican ! scattering element canopy layer index (top down)
integer :: icol ! scattering element column
real(r8) :: asu ! single scattering albedo
real(r8) :: gdir
real(r8) :: tmp0,tmp1,tmp2
real(r8) :: betab_veg ! beam backscatter for vegetation (no snow)
real(r8) :: betab_om ! multiplication of beam backscatter and reflectance
real(r8) :: om_veg ! scattering coefficient for vegetation (no snow)
real(r8) :: Kb_sing ! the KB_leaf that would generate a singularity
! with the scelb%a parameter
real(r8) :: Kb_stem ! actual optical depth of stem with not planar geometry effects
! usually the base value
real(r8), parameter :: Kb_stem_base = 1.0_r8
real(r8), parameter :: sing_tol = 0.01_r8 ! allowable difference between
! the Kb_leaf that creates
! a singularity and the actual
if( (cosz-1.0) > nearzero ) then
write(log_unit,*)"The cosine of the zenith angle cannot exceed 1"
write(log_unit,*)"cosz: ",cosz
write(log_unit,*)"TwoStreamMLPEMod.F90:ZenithPrep"
call endrun(msg=errMsg(sourcefile, __LINE__))
elseif(cosz<0._r8)then
write(log_unit,*)"The cosine of the zenith angle should not be less than zero"
write(log_unit,*)"It can be exactly zero, but not less than"
write(log_unit,*)"cosz: ",cosz
write(log_unit,*)"TwoStreamMLPEMod.F90:ZenithPrep"
call endrun(msg=errMsg(sourcefile, __LINE__))
end if
cosz = max(0.001,cosz)
this%cosz = cosz
do_ican: do ican = 1,this%n_lyr
do_ical: do icol = 1,this%n_col(ican)
associate(ft => this%scelg(ican,icol)%pft, &
scelg => this%scelg(ican,icol))
if(ft==air_ft)then
! Simple provisions for a ghost element (air)
scelg%Kb_leaf = k_air
scelg%Kb = k_air
else
gdir = rad_params%phi1(ft) + rad_params%phi2(ft) * cosz
Kb_stem = Kb_stem_base
!how much direct light penetrates a singleunit of lai?
scelg%Kb_leaf = min(kb_max,rad_params%clumping_index(ft) * gdir / cosz)
! To avoid singularities, we need to make sure that Kb =/ a
! If they are too similar, it will create a very large
! term in the linear solution and generate solution errors
! Lets identify the Kb_leaf that gives a singularity.
! We don't need to include the min() function
! a will never be that large.
!
! kb = a = (lai*kb_leaf + sai*kb_stem)/(lai+sai)
! (a*(lai+sai) - sai*kb_stem)/lai = Kb_sing
! or.. adjust stem Kb?
! (a*(lai+sai) - lai*kb_leaf)/sai = kb_stem_sing
if(scelg%lai>nearzero) then
do ib = 1,this%n_bands
Kb_sing = (this%band(ib)%scelb(ican,icol)%a*(scelg%lai+scelg%sai) - scelg%sai*Kb_stem)/scelg%lai
if(abs(scelg%Kb_leaf - Kb_sing)<sing_tol)then
scelg%Kb_leaf = Kb_sing + sing_tol
end if
end do
else
do ib = 1,this%n_bands
Kb_sing = this%band(ib)%scelb(ican,icol)%a
if(abs(Kb_stem - Kb_sing)<sing_tol)then
Kb_stem = Kb_sing + sing_tol
end if
end do
end if
! RGK: My sense is that snow should be adding optical depth
! but we don't have any precedent for that in the FATES
! code or old CLM. Re-view this.
!!scelbp%Kb = this%frac_snow*k_snow + scelbp%Kb
scelg%Kb = min(kb_max,(scelg%lai*scelg%Kb_leaf + scelg%sai*Kb_stem)/(scelg%lai+scelg%sai))
! Component terms for asu (single scatering albedo)
tmp0 = gdir + rad_params%phi2(ft) * cosz
tmp1 = rad_params%phi1(ft) * cosz
tmp2 = 1._r8 - tmp1/tmp0 * log((tmp1+tmp0)/tmp1)
end if
do_ib: do ib = 1,this%n_bands
associate( scelb => this%band(ib)%scelb(ican,icol) )
if(ft==air_ft)then
! Simple provisions for a ghost element (air)
scelb%betab = beta_air
else
! betab - upscatter parameter for direct beam radiation, from G. Bonan
! Eq. 3.16 CLM50 Tech Man
! asu is the single scattering albedo per om_veg (material reflectance)