forked from NGEET/fates
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathSFMainMod.F90
1411 lines (1117 loc) · 60.5 KB
/
SFMainMod.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 SFMainMod
! ============================================================================
! All subroutines realted to the SPITFIRE fire routine.
! Code originally developed by Allan Spessa & Rosie Fisher as part of the NERC-QUEST project.
! ============================================================================
use FatesConstantsMod , only : r8 => fates_r8
use FatesConstantsMod , only : itrue, ifalse
use FatesInterfaceMod , only : hlm_masterproc ! 1= master process, 0=not master process
use EDTypesMod , only : numWaterMem
use FatesGlobals , only : fates_log
use FatesInterfaceMod , only : bc_in_type
use EDPftvarcon , only : EDPftvarcon_inst
use EDtypesMod , only : ed_site_type
use EDtypesMod , only : ed_patch_type
use EDtypesMod , only : ed_cohort_type
use EDtypesMod , only : AREA
use EDtypesMod , only : DL_SF
use EDtypesMod , only : FIRE_THRESHOLD
use EDTypesMod , only : TW_SF
use EDtypesMod , only : LB_SF
use EDtypesMod , only : LG_SF
use EDtypesMod , only : NCWD
use EDtypesMod , only : NFSC
use EDtypesMod , only : TR_SF
use PRTGenericMod, only : leaf_organ
use PRTGenericMod, only : all_carbon_elements
use PRTGenericMod, only : leaf_organ
use PRTGenericMod, only : fnrt_organ
use PRTGenericMod, only : sapw_organ
use PRTGenericMod, only : store_organ
use PRTGenericMod, only : repro_organ
use PRTGenericMod, only : struct_organ
use PRTGenericMod, only : SetState
use spmdMod, only : masterproc, mpicom, comp_id
use shr_kind_mod, only : CL => shr_kind_CL
use shr_strdata_mod, only : shr_strdata_type, shr_strdata_create
use shr_strdata_mod, only : shr_strdata_advance, shr_strdata_print
use shr_log_mod, only : errmsg => shr_log_errMsg
use abortutils, only : endrun
use clm_varctl, only : iulog
use fileutils, only : getavu, relavu
use decompMod, only : bounds_type, gsmap_lnd_gdc2glo
use domainMod, only : ldomain
implicit none
private
public :: fire_model
public :: fire_danger_index
public :: charecteristics_of_fuel
public :: rate_of_spread
public :: ground_fuel_consumption
public :: fire_intensity
public :: wind_effect
public :: area_burnt
public :: crown_scorching
public :: crown_damage
public :: cambial_damage_kill
public :: post_fire_mortality
public :: FATESFireInit ! Initialize fire-related inputs
public :: FATESFireInterp ! Interpolate fire-related inputs
type, public :: sfmain_type
real(r8), public, pointer :: lnfm24(:) ! Daily avg lightning by grid cell (#/km2/hr)
contains
procedure, public :: InitAccBuffer ! Initialize accumulation processes
procedure, public :: InitAccVars ! Initialize accumulation variables
procedure, public :: UpdateAccVars ! Update/extract accumulations vars
end type sfmain_type
! !PRIVATE MEMBER FUNCTIONS:
private :: lnfm_init ! position datasets for Lightning
private :: lnfm_interp ! interpolates between two years of Lightning file data
real(r8), pointer :: forc_lnfm(:) ! Lightning frequency from file (#/km2/hr)
type(shr_strdata_type) :: sdat_lnfm ! Lightning input data stream
character(len=CL) :: stream_fldFileName_lightng ! lightning stream filename to read
character(len=*), parameter, private :: sourcefile = &
__FILE__
integer :: write_SF = 0 ! for debugging
logical :: debug = .false. ! for debugging
! ============================================================================
! ============================================================================
contains
subroutine InitAccBuffer (this, bounds)
!
! !DESCRIPTION:
! Initialize accumulation buffer for all required module accumulated fields
! This routine set defaults values that are then overwritten by the
! restart file for restart or branch runs
!
! !USES
use clm_varcon, only : spval
use accumulMod, only : init_accum_field
!
! !ARGUMENTS:
class(sfmain_type) :: this
type(bounds_type), intent(in) :: bounds
! !LOCAL VARIABLES:
integer :: begg, endg
integer :: ier
!---------------------------------------------------------------------
begg = bounds%begg; endg = bounds%endg
allocate(this%lnfm24(begg:endg), stat=ier)
if (ier/=0) then
call endrun(msg="allocation error for lnfm24"//&
errMsg(sourcefile, __LINE__))
endif
this%lnfm24(:) = spval
call init_accum_field (name='lnfm24', units='strikes/km2/hr', &
desc='24hr average of lightning strikes', accum_type='runmean', &
accum_period=-1, subgrid_type='gridcell', numlev=1, init_value=0._r8)
end subroutine InitAccBuffer
!-----------------------------------------------------------------------
subroutine InitAccVars(this, bounds)
!
! !DESCRIPTION:
! Initialize module variables that are associated with
! time accumulated fields. This routine is called for both an initial run
! and a restart run (and must therefore must be called after the restart
! file is read in and the accumulation buffer is obtained)
!
! !USES
use accumulMod , only : extract_accum_field
use clm_time_manager , only : get_nstep
!
! !ARGUMENTS:
class(sfmain_type) :: this
type(bounds_type), intent(in) :: bounds
!
! !LOCAL VARIABLES:
integer :: begg, endg
integer :: nstep
integer :: ier
real(r8), pointer :: rbufslg(:) ! temporary
!---------------------------------------------------------------------
begg = bounds%begg; endg = bounds%endg
! Allocate needed dynamic memory for single level patch field
allocate(rbufslg(begg:endg), stat=ier)
if (ier/=0) then
call endrun(msg="allocation error for rbufslg"//&
errMsg(sourcefile, __LINE__))
endif
! Determine time step
nstep = get_nstep()
call extract_accum_field ('lnfm24', rbufslg, nstep)
this%lnfm24(begg:endg) = rbufslg(begg:endg)
deallocate(rbufslg)
end subroutine InitAccVars
!-----------------------------------------------------------------------
subroutine UpdateAccVars (this, bounds)
!
! USES
use clm_time_manager, only : get_nstep
use accumulMod , only : update_accum_field, extract_accum_field
use abortutils , only : endrun
!
! !ARGUMENTS:
class(sfmain_type) :: this
type(bounds_type), intent(in) :: bounds
!
! !LOCAL VARIABLES:
integer :: dtime ! timestep size [seconds]
integer :: nstep ! timestep number
integer :: ier ! error status
integer :: begg, endg
real(r8), pointer :: rbufslg(:) ! temporary single level - gridcell level
!---------------------------------------------------------------------
begg = bounds%begg; endg = bounds%endg
nstep = get_nstep()
! Allocate needed dynamic memory for single level gridcell field
allocate(rbufslg(begg:endg), stat=ier)
if (ier/=0) then
write(iulog,*)'UpdateAccVars allocation error for rbuf1dg'
call endrun(msg=errMsg(sourcefile, __LINE__))
endif
! Accumulate and extract lnfm24
rbufslg(begg:endg) = forc_lnfm(begg:endg)
call update_accum_field ('lnfm24', rbufslg, nstep)
call extract_accum_field ('lnfm24', this%lnfm24, nstep)
deallocate(rbufslg)
end subroutine UpdateAccVars
!-----------------------------------------------------------------------
subroutine FATESFireInit(bounds, NLFilename)
! !DESCRIPTION:
! Initialize FATES Fire module
! !USES:
use FatesInterfaceMod, only : hlm_use_spitfire
use shr_infnan_mod , only : nan => shr_infnan_nan, assignment(=)
!
! !ARGUMENTS:
type(bounds_type), intent(in) :: bounds
character(len=*), intent(in) :: NLFilename
! !LOCAL VARIABLES:
integer :: begg, endg
!-----------------------------------------------------------------------
begg = bounds%begg; endg = bounds%endg
! Allocate lightning forcing data
allocate( forc_lnfm(begg:endg) )
forc_lnfm(begg:) = nan
if( hlm_use_spitfire == itrue )then
call lnfm_init(bounds, NLFilename)
call lnfm_interp(bounds)
end if
end subroutine FATESFireInit
!*****************************************************************
subroutine FATESFireInterp(bounds)
! !DESCRIPTION:
! Interpolate FATES Fire datasets
! !USES:
use FatesInterfaceMod, only : hlm_use_spitfire
!
! !ARGUMENTS:
type(bounds_type), intent(in) :: bounds
!-----------------------------------------------------------------------
if( hlm_use_spitfire == itrue )then
call lnfm_interp(bounds)
end if
end subroutine FATESFireInterp
!-----------------------------------------------------------------------
subroutine lnfm_init( bounds, NLFilename )
!
! !DESCRIPTION:
!
! Initialize data stream information for Lightning.
!
! !USES:
use clm_varctl , only : inst_name
use clm_time_manager , only : get_calendar
use ncdio_pio , only : pio_subsystem
use shr_pio_mod , only : shr_pio_getiotype
use clm_nlUtilsMod , only : find_nlgroup_name
use ndepStreamMod , only : clm_domain_mct
use histFileMod , only : hist_addfld1d
use shr_mpi_mod, only : shr_mpi_bcast
use mct_mod, only: mct_ggrid
!
! !ARGUMENTS:
implicit none
type(bounds_type), intent(in) :: bounds
character(len=*), intent(in) :: NLFilename
!
! !LOCAL VARIABLES:
integer :: stream_year_first_lightng ! first year in Lightning stream to use
integer :: stream_year_last_lightng ! last year in Lightning stream to use
integer :: model_year_align_lightng ! align stream_year_first_lnfm with
integer :: nu_nml ! unit for namelist file
integer :: nml_error ! namelist i/o error flag
type(mct_ggrid) :: dom_clm ! domain information
character(len=CL) :: lightngmapalgo = 'bilinear'! Mapping alogrithm
character(*), parameter :: subName = "('lnfmdyn_init')"
character(*), parameter :: F00 = "('(lnfmdyn_init) ',4a)"
!-----------------------------------------------------------------------
namelist /light_streams/ &
stream_year_first_lightng, &
stream_year_last_lightng, &
model_year_align_lightng, &
lightngmapalgo, &
stream_fldFileName_lightng
! Default values for namelist
stream_year_first_lightng = 1 ! first year in stream to use
stream_year_last_lightng = 1 ! last year in stream to use
model_year_align_lightng = 1 ! align stream_year_first_lnfm with this model year
stream_fldFileName_lightng = ' '
! Read light_streams namelist
if (masterproc) then
nu_nml = getavu()
open( nu_nml, file=trim(NLFilename), status='old', iostat=nml_error )
call find_nlgroup_name(nu_nml, 'light_streams', status=nml_error)
if (nml_error == 0) then
read(nu_nml, nml=light_streams,iostat=nml_error)
if (nml_error /= 0) then
call endrun(msg='ERROR reading light_streams namelist'//errmsg(sourcefile, __LINE__))
end if
end if
close(nu_nml)
call relavu( nu_nml )
endif
call shr_mpi_bcast(stream_year_first_lightng, mpicom)
call shr_mpi_bcast(stream_year_last_lightng, mpicom)
call shr_mpi_bcast(model_year_align_lightng, mpicom)
call shr_mpi_bcast(stream_fldFileName_lightng, mpicom)
if (masterproc) then
write(iulog,*) ' '
write(iulog,*) 'light_stream settings:'
write(iulog,*) ' stream_year_first_lightng = ',stream_year_first_lightng
write(iulog,*) ' stream_year_last_lightng = ',stream_year_last_lightng
write(iulog,*) ' model_year_align_lightng = ',model_year_align_lightng
write(iulog,*) ' stream_fldFileName_lightng = ',stream_fldFileName_lightng
write(iulog,*) ' '
endif
if (index(stream_fldFileName_lightng, 'nofile') == 0) then
call clm_domain_mct (bounds, dom_clm)
call shr_strdata_create(sdat_lnfm,name="clmlnfm", &
pio_subsystem=pio_subsystem, &
pio_iotype=shr_pio_getiotype(inst_name), &
mpicom=mpicom, compid=comp_id, &
gsmap=gsmap_lnd_gdc2glo, ggrid=dom_clm, &
nxg=ldomain%ni, nyg=ldomain%nj, &
yearFirst=stream_year_first_lightng, &
yearLast=stream_year_last_lightng, &
yearAlign=model_year_align_lightng, &
offset=0, &
domFilePath='', &
domFileName=trim(stream_fldFileName_lightng), &
domTvarName='time', &
domXvarName='lon' , &
domYvarName='lat' , &
domAreaName='area', &
domMaskName='mask', &
filePath='', &
filename=(/trim(stream_fldFileName_lightng)/),&
fldListFile='lnfm', &
fldListModel='lnfm', &
fillalgo='none', &
mapalgo=lightngmapalgo, &
calendar=get_calendar(), &
taxmode='cycle' )
if (masterproc) then
call shr_strdata_print(sdat_lnfm,'Lightning data')
endif
! Add history fields
call hist_addfld1d (fname='LNFM', units='counts/km^2/hr', &
avgflag='A', long_name='Lightning frequency', &
ptr_lnd=forc_lnfm, default='inactive')
end if
end subroutine lnfm_init
!-----------------------------------------------------------------------
subroutine lnfm_interp( bounds )
!
! !DESCRIPTION:
! Interpolate data stream information for Lightning.
!
! !USES:
use clm_time_manager, only : get_curr_date
!
! !ARGUMENTS:
type(bounds_type), intent(in) :: bounds
!
! !LOCAL VARIABLES:
integer :: g, ig
integer :: year ! year (0, ...) for nstep+1
integer :: mon ! month (1, ..., 12) for nstep+1
integer :: day ! day of month (1, ..., 31) for nstep+1
integer :: sec ! seconds into current date for nstep+1
integer :: mcdate ! Current model date (yyyymmdd)
!-----------------------------------------------------------------------
call get_curr_date(year, mon, day, sec)
mcdate = year*10000 + mon*100 + day
if (index(stream_fldFileName_lightng, 'nofile') == 0) then
call shr_strdata_advance(sdat_lnfm, mcdate, sec, mpicom, 'lnfmdyn')
ig = 0
do g = bounds%begg,bounds%endg
ig = ig+1
forc_lnfm(g) = sdat_lnfm%avs(1)%rAttr(1,ig)
end do
end if
end subroutine lnfm_interp
! ============================================================================
! Area of site burned by fire
! ============================================================================
subroutine fire_model( currentSite, bc_in)
use FatesInterfaceMod, only : hlm_use_spitfire
type(ed_site_type) , intent(inout), target :: currentSite
type(bc_in_type) , intent(in) :: bc_in
type (ed_patch_type), pointer :: currentPatch
!zero fire things
currentPatch => currentSite%youngest_patch
do while(associated(currentPatch))
currentPatch%frac_burnt = 0.0_r8
currentPatch%AB = 0.0_r8
currentPatch%fire = 0
currentPatch => currentPatch%older
enddo
if(write_SF==1)then
write(fates_log(),*) 'use_spitfire',hlm_use_spitfire
endif
if( hlm_use_spitfire == itrue )then
call fire_danger_index(currentSite, bc_in)
call wind_effect(currentSite, bc_in)
call charecteristics_of_fuel(currentSite)
call rate_of_spread(currentSite)
call ground_fuel_consumption(currentSite)
call fire_intensity(currentSite)
call area_burnt(currentSite, bc_in)
call crown_scorching(currentSite)
call crown_damage(currentSite)
call cambial_damage_kill(currentSite)
call post_fire_mortality(currentSite)
end if
end subroutine fire_model
!*****************************************************************
subroutine fire_danger_index ( currentSite, bc_in)
!*****************************************************************
! currentSite%acc_NI is the accumulated Nesterov fire danger index
use SFParamsMod, only : SF_val_fdi_a, SF_val_fdi_b
use FatesConstantsMod , only : tfrz => t_water_freeze_k_1atm
use FatesConstantsMod , only : sec_per_day
type(ed_site_type) , intent(inout), target :: currentSite
type(bc_in_type) , intent(in) :: bc_in
real(r8) :: temp_in_C ! daily averaged temperature in celcius
real(r8) :: rainfall ! daily precip in mm/day
real(r8) :: rh ! daily rh
real yipsolon; !intermediate varable for dewpoint calculation
real dewpoint; !dewpoint in K
real d_NI; !daily change in Nesterov Index. C^2
integer :: iofp ! index of oldest the fates patch
! NOTE that the boundary conditions of temperature, precipitation and relative humidity
! are available at the patch level. We are currently using a simplification where the whole site
! is simply using the values associated with the first patch.
! which probably won't have much inpact, unless we decide to ever calculated the NI for each patch.
iofp = currentSite%oldest_patch%patchno
temp_in_C = bc_in%t_veg24_pa(iofp) - tfrz
rainfall = bc_in%precip24_pa(iofp)*sec_per_day
rh = bc_in%relhumid24_pa(iofp)
if (rainfall > 3.0_r8) then !rezero NI if it rains...
d_NI = 0.0_r8
currentSite%acc_NI = 0.0_r8
else
yipsolon = (SF_val_fdi_a* temp_in_C)/(SF_val_fdi_b+ temp_in_C)+log(rh/100.0_r8)
dewpoint = (SF_val_fdi_b*yipsolon)/(SF_val_fdi_a-yipsolon) !Standard met. formula
d_NI = ( temp_in_C-dewpoint)* temp_in_C !follows Nesterov 1968. Equation 5. Thonicke et al. 2010.
if (d_NI < 0.0_r8) then !Change in NI cannot be negative.
d_NI = 0.0_r8 !check
endif
endif
currentSite%acc_NI = currentSite%acc_NI + d_NI !Accumulate Nesterov index over the fire season.
end subroutine fire_danger_index
!*****************************************************************
subroutine charecteristics_of_fuel ( currentSite )
!*****************************************************************
use SFParamsMod, only : SF_val_drying_ratio, SF_val_SAV, SF_val_FBD
type(ed_site_type), intent(in), target :: currentSite
type(ed_patch_type), pointer :: currentPatch
type(ed_cohort_type), pointer :: currentCohort
real(r8) timeav_swc
real(r8) alpha_FMC(nfsc) ! Relative fuel moisture adjusted per drying ratio
real(r8) fuel_moisture(nfsc) ! Scaled moisture content of small litter fuels.
real(r8) MEF(nfsc) ! Moisture extinction factor of fuels integer n
fuel_moisture(:) = 0.0_r8
currentPatch => currentSite%oldest_patch;
do while(associated(currentPatch))
! How much live grass is there?
currentPatch%livegrass = 0.0_r8
currentCohort => currentPatch%tallest
do while(associated(currentCohort))
if(EDPftvarcon_inst%woody(currentCohort%pft) == 0)then
currentPatch%livegrass = currentPatch%livegrass + &
currentCohort%prt%GetState(leaf_organ, all_carbon_elements) * &
currentCohort%n/currentPatch%area
endif
currentCohort => currentCohort%shorter
enddo
! There are SIX fuel classes
! 1:4) four CWD_AG pools (twig, s branch, l branch, trunk), 5) dead leaves and 6) live grass
! NCWD =4 NFSC = 6
! tw_sf = 1, lb_sf = 3, tr_sf = 4, dl_sf = 5, lg_sf = 6,
! zero fire arrays.
currentPatch%fuel_eff_moist = 0.0_r8
currentPatch%fuel_bulkd = 0.0_r8 !this is kgBiomass/m2 for use in rate of spread equations
currentPatch%fuel_sav = 0.0_r8
currentPatch%fuel_frac(:) = 0.0_r8
currentPatch%fuel_mef = 0.0_r8
currentPatch%sum_fuel = 0.0_r8
currentPatch%fuel_frac = 0.0_r8
if(write_sf == itrue)then
if ( hlm_masterproc == itrue ) write(fates_log(),*) ' leaf_litter1 ',currentPatch%leaf_litter
if ( hlm_masterproc == itrue ) write(fates_log(),*) ' leaf_litter2 ',sum(currentPatch%CWD_AG)
if ( hlm_masterproc == itrue ) write(fates_log(),*) ' leaf_litter3 ',currentPatch%livegrass
if ( hlm_masterproc == itrue ) write(fates_log(),*) ' sum fuel', currentPatch%sum_fuel
endif
currentPatch%sum_fuel = sum(currentPatch%leaf_litter) + sum(currentPatch%CWD_AG) + currentPatch%livegrass
if(write_SF == itrue)then
if ( hlm_masterproc == itrue ) write(fates_log(),*) 'sum fuel', currentPatch%sum_fuel,currentPatch%area
endif
! ===============================================
! Average moisture, bulk density, surface area-volume and moisture extinction of fuel
! ================================================
if (currentPatch%sum_fuel > 0.0) then
! Fraction of fuel in litter classes
currentPatch%fuel_frac(dl_sf) = sum(currentPatch%leaf_litter)/ currentPatch%sum_fuel
currentPatch%fuel_frac(tw_sf:tr_sf) = currentPatch%CWD_AG / currentPatch%sum_fuel
if(write_sf == itrue)then
if ( hlm_masterproc == itrue ) write(fates_log(),*) 'ff1 ',currentPatch%fuel_frac
if ( hlm_masterproc == itrue ) write(fates_log(),*) 'ff2 ',currentPatch%fuel_frac
if ( hlm_masterproc == itrue ) write(fates_log(),*) 'ff2a ', &
lg_sf,currentPatch%livegrass,currentPatch%sum_fuel
endif
currentPatch%fuel_frac(lg_sf) = currentPatch%livegrass / currentPatch%sum_fuel
MEF(1:nfsc) = 0.524_r8 - 0.066_r8 * log10(SF_val_SAV(1:nfsc))
!--- weighted average of relative moisture content---
! Equation 6 in Thonicke et al. 2010. across twig, small branch, large branch, and dead leaves
! dead leaves and twigs included in 1hr pool per Thonicke (2010)
! Calculate fuel moisture for trunks to hold value for fuel consumption
alpha_FMC(tw_sf:dl_sf) = SF_val_SAV(tw_sf:dl_sf)/SF_val_drying_ratio
fuel_moisture(tw_sf:dl_sf) = exp(-1.0_r8 * alpha_FMC(tw_sf:dl_sf) * currentSite%acc_NI)
if(write_SF == itrue)then
if ( hlm_masterproc == itrue ) write(fates_log(),*) 'ff3 ',currentPatch%fuel_frac
if ( hlm_masterproc == itrue ) write(fates_log(),*) 'fm ',fuel_moisture
if ( hlm_masterproc == itrue ) write(fates_log(),*) 'csa ',currentSite%acc_NI
if ( hlm_masterproc == itrue ) write(fates_log(),*) 'sfv ',alpha_FMC
endif
! FIX(RF,032414): needs refactoring.
! average water content !is this the correct metric?
timeav_swc = sum(currentSite%water_memory(1:numWaterMem)) / real(numWaterMem,r8)
! Equation B2 in Thonicke et al. 2010
! live grass moisture content depends on upper soil layer
fuel_moisture(lg_sf) = max(0.0_r8, 10.0_r8/9._r8 * timeav_swc - 1.0_r8/9.0_r8)
! Average properties over the first three litter pools (twigs, s branches, l branches)
currentPatch%fuel_bulkd = sum(currentPatch%fuel_frac(tw_sf:lb_sf) * SF_val_FBD(tw_sf:lb_sf))
currentPatch%fuel_sav = sum(currentPatch%fuel_frac(tw_sf:lb_sf) * SF_val_SAV(tw_sf:lb_sf))
currentPatch%fuel_mef = sum(currentPatch%fuel_frac(tw_sf:lb_sf) * MEF(tw_sf:lb_sf))
currentPatch%fuel_eff_moist = sum(currentPatch%fuel_frac(tw_sf:lb_sf) * fuel_moisture(tw_sf:lb_sf))
if(write_sf == itrue)then
if ( hlm_masterproc == itrue ) write(fates_log(),*) 'ff4 ',currentPatch%fuel_eff_moist
endif
! Add on properties of dead leaves and live grass pools (5 & 6)
currentPatch%fuel_bulkd = currentPatch%fuel_bulkd + sum(currentPatch%fuel_frac(dl_sf:lg_sf) * SF_val_FBD(dl_sf:lg_sf))
currentPatch%fuel_sav = currentPatch%fuel_sav + sum(currentPatch%fuel_frac(dl_sf:lg_sf) * SF_val_SAV(dl_sf:lg_sf))
currentPatch%fuel_mef = currentPatch%fuel_mef + sum(currentPatch%fuel_frac(dl_sf:lg_sf) * MEF(dl_sf:lg_sf))
currentPatch%fuel_eff_moist = currentPatch%fuel_eff_moist+ sum(currentPatch%fuel_frac(dl_sf:lg_sf) * fuel_moisture(dl_sf:lg_sf))
! Correct averaging for the fact that we are not using the trunks pool for fire ROS and intensity (5)
! Consumption of fuel in trunk pool does not influence fire ROS or intensity (Pyne 1996)
currentPatch%fuel_bulkd = currentPatch%fuel_bulkd * (1.0_r8/(1.0_r8-currentPatch%fuel_frac(tr_sf)))
currentPatch%fuel_sav = currentPatch%fuel_sav * (1.0_r8/(1.0_r8-currentPatch%fuel_frac(tr_sf)))
currentPatch%fuel_mef = currentPatch%fuel_mef * (1.0_r8/(1.0_r8-currentPatch%fuel_frac(tr_sf)))
currentPatch%fuel_eff_moist = currentPatch%fuel_eff_moist * (1.0_r8/(1.0_r8-currentPatch%fuel_frac(tr_sf)))
! Pass litter moisture into the fuel burning routine (all fuels: twigs,s branch,l branch,trunk,dead leaves,live grass)
! (wo/me term in Thonicke et al. 2010)
currentPatch%litter_moisture(tw_sf:lb_sf) = fuel_moisture(tw_sf:lb_sf)/MEF(tw_sf:lb_sf)
currentPatch%litter_moisture(tr_sf) = fuel_moisture(tr_sf)/MEF(tr_sf)
currentPatch%litter_moisture(dl_sf) = fuel_moisture(dl_sf)/MEF(dl_sf)
currentPatch%litter_moisture(lg_sf) = fuel_moisture(lg_sf)/MEF(lg_sf)
else
if(write_SF == itrue)then
if ( hlm_masterproc == itrue ) write(fates_log(),*) 'no litter fuel at all',currentPatch%patchno, &
currentPatch%sum_fuel,sum(currentPatch%cwd_ag), &
sum(currentPatch%cwd_bg),sum(currentPatch%leaf_litter)
endif
currentPatch%fuel_sav = sum(SF_val_SAV(1:nfsc))/(nfsc) ! make average sav to avoid crashing code.
if ( hlm_masterproc == itrue ) write(fates_log(),*) 'problem with spitfire fuel averaging'
! FIX(SPM,032414) refactor...should not have 0 fuel unless everything is burnt
! off.
currentPatch%fuel_eff_moist = 0.0000000001_r8
currentPatch%fuel_bulkd = 0.0000000001_r8
currentPatch%fuel_frac(:) = 0.0000000001_r8
currentPatch%fuel_mef = 0.0000000001_r8
currentPatch%sum_fuel = 0.0000000001_r8
currentPatch%fuel_frac = 0.0000000001_r8
endif
! check values.
! FIX(SPM,032414) refactor...
if(write_SF == itrue.and.currentPatch%fuel_sav <= 0.0_r8.or.currentPatch%fuel_bulkd <= &
0.0_r8.or.currentPatch%fuel_mef <= 0.0_r8.or.currentPatch%fuel_eff_moist <= 0.0_r8)then
if ( hlm_masterproc == itrue ) write(fates_log(),*) 'problem with spitfire fuel averaging'
endif
currentPatch => currentPatch%younger
enddo !end patch loop
end subroutine charecteristics_of_fuel
!*****************************************************************
subroutine wind_effect ( currentSite, bc_in)
!*****************************************************************.
! Routine called daily from within ED within a site loop.
! Calculates the effective windspeed based on vegetation charecteristics.
! currentSite%wind is daily wind converted to m/min for Spitfire units
use FatesConstantsMod, only : sec_per_min
type(ed_site_type) , intent(inout), target :: currentSite
type(bc_in_type) , intent(in) :: bc_in
type(ed_patch_type) , pointer :: currentPatch
type(ed_cohort_type), pointer :: currentCohort
real(r8) :: total_grass_area ! per patch,in m2
real(r8) :: tree_fraction ! site level. no units
real(r8) :: grass_fraction ! site level. no units
real(r8) :: bare_fraction ! site level. no units
integer :: iofp ! index of oldest fates patch
! note - this is a patch level temperature, which probably won't have much inpact,
! unless we decide to ever calculated the NI for each patch.
iofp = currentSite%oldest_patch%patchno
currentSite%wind = bc_in%wind24_pa(iofp) * sec_per_min !Convert to m/min for SPITFIRE
if(write_SF == itrue)then
if ( hlm_masterproc == itrue ) write(fates_log(),*) 'wind24', currentSite%wind
endif
! --- influence of wind speed, corrected for surface roughness----
! --- averaged over the whole grid cell to prevent extreme divergence
! average_wspeed = 0.0_r8
tree_fraction = 0.0_r8
grass_fraction = 0.0_r8
currentPatch=>currentSite%oldest_patch;
do while(associated(currentPatch))
currentPatch%total_tree_area = 0.0_r8
total_grass_area = 0.0_r8
currentCohort => currentPatch%tallest
do while(associated(currentCohort))
if (debug) write(fates_log(),*) 'SF currentCohort%c_area ',currentCohort%c_area
if(EDPftvarcon_inst%woody(currentCohort%pft) == 1)then
currentPatch%total_tree_area = currentPatch%total_tree_area + currentCohort%c_area
else
total_grass_area = total_grass_area + currentCohort%c_area
endif
currentCohort => currentCohort%shorter
enddo
tree_fraction = tree_fraction + min(currentPatch%area,currentPatch%total_tree_area)/AREA
grass_fraction = grass_fraction + min(currentPatch%area,total_grass_area)/AREA
if(debug)then
write(fates_log(),*) 'SF currentPatch%area ',currentPatch%area
write(fates_log(),*) 'SF currentPatch%total_tree_area ',currentPatch%total_tree_area
write(fates_log(),*) 'SF total_grass_area ', total_grass_area
write(fates_log(),*) 'SF AREA ',AREA
endif
currentPatch => currentPatch%younger
enddo !currentPatch loop
!if there is a cover of more than one, then the grasses are under the trees
grass_fraction = min(grass_fraction,1.0_r8-tree_fraction)
bare_fraction = 1.0_r8 - tree_fraction - grass_fraction
if(write_sf == itrue)then
if ( hlm_masterproc == itrue ) write(fates_log(),*) 'grass, trees, bare', &
grass_fraction, tree_fraction, bare_fraction
endif
currentPatch=>currentSite%oldest_patch;
do while(associated(currentPatch))
currentPatch%total_tree_area = min(currentPatch%total_tree_area,currentPatch%area)
! effect_wspeed in units m/min
currentPatch%effect_wspeed = currentSite%wind * (tree_fraction*0.4_r8+(grass_fraction+bare_fraction)*0.6_r8)
currentPatch => currentPatch%younger
enddo !end patch loop
end subroutine wind_effect
!*****************************************************************
subroutine rate_of_spread ( currentSite )
!*****************************************************************.
!Routine called daily from within ED within a site loop.
!Returns the updated currentPatch%ROS_front value for each patch.
use SFParamsMod, only : SF_val_miner_total, &
SF_val_part_dens, &
SF_val_miner_damp, &
SF_val_fuel_energy
use FatesInterfaceMod, only : hlm_current_day, hlm_current_month
type(ed_site_type), intent(in), target :: currentSite
type(ed_patch_type), pointer :: currentPatch
! Rothermal fire spread model parameters.
real(r8) beta,beta_op ! weighted average of packing ratio (unitless)
real(r8) ir ! reaction intensity (kJ/m2/min)
real(r8) xi,eps,phi_wind ! all are unitless
real(r8) q_ig ! heat of pre-ignition (kJ/kg)
real(r8) reaction_v_opt,reaction_v_max !reaction velocity (per min)!optimum and maximum
real(r8) moist_damp,mw_weight ! moisture dampening coefficient and ratio fuel moisture to extinction
real(r8) beta_ratio ! ratio of beta/beta_op
real(r8) a_beta ! dummy variable for product of a* beta_ratio for react_v_opt equation
real(r8) a,b,c,e ! function of fuel sav
logical,parameter :: debug_windspeed = .false. !for debugging
currentPatch=>currentSite%oldest_patch;
do while(associated(currentPatch))
! ---initialise parameters to zero.---
beta_ratio = 0.0_r8; q_ig = 0.0_r8; eps = 0.0_r8; a = 0.0_r8; b = 0.0_r8; c = 0.0_r8; e = 0.0_r8
phi_wind = 0.0_r8; xi = 0.0_r8; reaction_v_max = 0.0_r8; reaction_v_opt = 0.0_r8; mw_weight = 0.0_r8
moist_damp = 0.0_r8; ir = 0.0_r8; a_beta = 0.0_r8;
currentPatch%ROS_front = 0.0_r8
! remove mineral content from net fuel load per Thonicke 2010 for ir calculation
currentPatch%sum_fuel = currentPatch%sum_fuel * (1.0_r8 - SF_val_miner_total) !net of minerals
! ----start spreading---
if ( hlm_masterproc == itrue .and.debug) write(fates_log(),*) &
'SF - currentPatch%fuel_bulkd ',currentPatch%fuel_bulkd
if ( hlm_masterproc == itrue .and.debug) write(fates_log(),*) &
'SF - SF_val_part_dens ',SF_val_part_dens
! beta = packing ratio (unitless)
! fraction of fuel array volume occupied by fuel or compactness of fuel bed
beta = currentPatch%fuel_bulkd / SF_val_part_dens
! Equation A6 in Thonicke et al. 2010
! packing ratio (unitless)
beta_op = 0.200395_r8 *(currentPatch%fuel_sav**(-0.8189_r8))
if ( hlm_masterproc == itrue .and.debug) write(fates_log(),*) 'SF - beta ',beta
if ( hlm_masterproc == itrue .and.debug) write(fates_log(),*) 'SF - beta_op ',beta_op
beta_ratio = beta/beta_op !unitless
if(write_sf == itrue)then
if ( hlm_masterproc == itrue ) write(fates_log(),*) 'esf ',currentPatch%fuel_eff_moist
endif
! ---heat of pre-ignition---
! Equation A4 in Thonicke et al. 2010
! conversion of Rohtermal (1972) equation 12 in BTU/lb to current kJ/kg
! q_ig in kJ/kg
q_ig = 581.0_r8 +2594.0_r8 * currentPatch%fuel_eff_moist
! ---effective heating number---
! Equation A3 in Thonicke et al. 2010.
eps = exp(-4.528_r8 / currentPatch%fuel_sav)
! Equation A7 in Thonicke et al. 2010
b = 0.15988_r8 * (currentPatch%fuel_sav**0.54_r8)
! Equation A8 in Thonicke et al. 2010
c = 7.47_r8 * (exp(-0.8711_r8 * (currentPatch%fuel_sav**0.55_r8)))
! Equation A9 in Thonicke et al. 2010.
e = 0.715_r8 * (exp(-0.01094_r8 * currentPatch%fuel_sav))
if (debug) then
if ( hlm_masterproc == itrue .and.debug) write(fates_log(),*) 'SF - c ',c
if ( hlm_masterproc == itrue .and.debug) write(fates_log(),*) 'SF - currentPatch%effect_wspeed ', &
currentPatch%effect_wspeed
if ( hlm_masterproc == itrue .and.debug) write(fates_log(),*) 'SF - b ',b
if ( hlm_masterproc == itrue .and.debug) write(fates_log(),*) 'SF - beta_ratio ',beta_ratio
if ( hlm_masterproc == itrue .and.debug) write(fates_log(),*) 'SF - e ',e
endif
! Equation A5 in Thonicke et al. 2010
! phi_wind (unitless)
! convert current_wspeed (wind at elev relevant to fire) from m/min to ft/min for Rothermel ROS eqn
phi_wind = c * ((3.281_r8*currentPatch%effect_wspeed)**b)*(beta_ratio**(-e))
! ---propagating flux----
! Equation A2 in Thonicke et al.2010 and Eq. 42 Rothermal 1972
! xi (unitless)
xi = (exp((0.792_r8 + 3.7597_r8 * (currentPatch%fuel_sav**0.5_r8)) * (beta+0.1_r8))) / &
(192_r8+7.9095_r8 * currentPatch%fuel_sav)
! ---reaction intensity----
! Equation in table A1 Thonicke et al. 2010.
a = 8.9033_r8 * (currentPatch%fuel_sav**(-0.7913_r8))
a_beta = exp(a*(1.0_r8-beta_ratio)) !dummy variable for reaction_v_opt equation
! Equation in table A1 Thonicke et al. 2010.
! reaction_v_max and reaction_v_opt = reaction velocity in units of per min
! reaction_v_max = Equation 36 in Rothermal 1972 and Fig 12
reaction_v_max = 1.0_r8 / (0.0591_r8 + 2.926_r8* (currentPatch%fuel_sav**(-1.5_r8)))
! reaction_v_opt = Equation 38 in Rothermal 1972 and Fig 11
reaction_v_opt = reaction_v_max*(beta_ratio**a)*a_beta
! mw_weight = relative fuel moisture/fuel moisture of extinction
! average values for litter pools (dead leaves, twigs, small and large branches) plus grass
mw_weight = currentPatch%fuel_eff_moist/currentPatch%fuel_mef
! Equation in table A1 Thonicke et al. 2010.
! moist_damp is unitless
moist_damp = max(0.0_r8,(1.0_r8 - (2.59_r8 * mw_weight) + (5.11_r8 * (mw_weight**2.0_r8)) - &
(3.52_r8*(mw_weight**3.0_r8))))
! FIX(SPM, 040114) ask RF if this should be an endrun
! if(write_SF == itrue)then
! write(fates_log(),*) 'moist_damp' ,moist_damp,mw_weight,currentPatch%fuel_eff_moist,currentPatch%fuel_mef
! endif
! ir = reaction intenisty in kJ/m2/min
! currentPatch%sum_fuel converted from kgC/m2 to kgBiomass/m2 for ir calculation
ir = reaction_v_opt*(currentPatch%sum_fuel/0.45_r8)*SF_val_fuel_energy*moist_damp*SF_val_miner_damp
! write(fates_log(),*) 'ir',gamma_aptr,moist_damp,SF_val_fuel_energy,SF_val_miner_damp
if (((currentPatch%fuel_bulkd) <= 0.0_r8).or.(eps <= 0.0_r8).or.(q_ig <= 0.0_r8)) then
currentPatch%ROS_front = 0.0_r8
else ! Equation 9. Thonicke et al. 2010.
! forward ROS in m/min
currentPatch%ROS_front = (ir*xi*(1.0_r8+phi_wind)) / (currentPatch%fuel_bulkd*eps*q_ig)
! write(fates_log(),*) 'ROS',currentPatch%ROS_front,phi_wind,currentPatch%effect_wspeed
! write(fates_log(),*) 'ros calcs',currentPatch%fuel_bulkd,ir,xi,eps,q_ig
endif
! Equation 10 in Thonicke et al. 2010
! backward ROS from Can FBP System (1992) in m/min
! backward ROS wind not changed by vegetation
currentPatch%ROS_back = currentPatch%ROS_front*exp(-0.012_r8*currentSite%wind)
currentPatch => currentPatch%younger
enddo !end patch loop
end subroutine rate_of_spread
!*****************************************************************
subroutine ground_fuel_consumption ( currentSite )
!*****************************************************************
!returns the the hypothetic fuel consumed by the fire
use SFParamsMod, only : SF_val_miner_total, SF_val_min_moisture, &
SF_val_mid_moisture, SF_val_low_moisture_Coeff, SF_val_low_moisture_Slope, &
SF_val_mid_moisture_Coeff, SF_val_mid_moisture_Slope
type(ed_site_type) , intent(in), target :: currentSite
type(ed_patch_type), pointer :: currentPatch
real(r8) :: moist ! effective fuel moisture
real(r8) :: tau_b(nfsc) ! lethal heating rates for each fuel class (min)
real(r8) :: fc_ground(nfsc) ! proportion of fuel consumed
integer :: c
currentPatch => currentSite%oldest_patch;
do while(associated(currentPatch))
currentPatch%burnt_frac_litter = 1.0_r8
! Calculate fraction of litter is burnt for all classes.
! Equation B1 in Thonicke et al. 2010---
do c = 1, nfsc !work out the burnt fraction for all pools, even if those pools dont exist.
moist = currentPatch%litter_moisture(c)
! 1. Very dry litter
if (moist <= SF_val_min_moisture(c)) then
currentPatch%burnt_frac_litter(c) = 1.0_r8
endif
! 2. Low to medium moistures
if (moist > SF_val_min_moisture(c).and.moist <= SF_val_mid_moisture(c)) then
currentPatch%burnt_frac_litter(c) = max(0.0_r8,min(1.0_r8,SF_val_low_moisture_Coeff(c)- &
SF_val_low_moisture_Slope(c)*moist))
else
! For medium to high moistures.
if (moist > SF_val_mid_moisture(c).and.moist <= 1.0_r8) then
currentPatch%burnt_frac_litter(c) = max(0.0_r8,min(1.0_r8,SF_val_mid_moisture_Coeff(c)- &
SF_val_mid_moisture_Slope(c)*moist))
endif
endif
! Very wet litter
if (moist >= 1.0_r8) then !this shouldn't happen?
currentPatch%burnt_frac_litter(c) = 0.0_r8
endif
enddo !c
! we can't ever kill -all- of the grass.
currentPatch%burnt_frac_litter(lg_sf) = min(0.8_r8,currentPatch%burnt_frac_litter(lg_sf ))
! reduce burnt amount for mineral content.
currentPatch%burnt_frac_litter = currentPatch%burnt_frac_litter * (1.0_r8-SF_val_miner_total)
!---Calculate amount of fuel burnt.---
FC_ground(tw_sf:tr_sf) = currentPatch%burnt_frac_litter(tw_sf:tr_sf) * currentPatch%CWD_AG
FC_ground(dl_sf) = currentPatch%burnt_frac_litter(dl_sf) * sum(currentPatch%leaf_litter)
FC_ground(lg_sf) = currentPatch%burnt_frac_litter(lg_sf) * currentPatch%livegrass
! Following used for determination of cambial kill follows from Peterson & Ryan (1986) scheme
! less empirical cf current scheme used in SPITFIRE which attempts to mesh Rothermel
! and P&R, and while solving potential inconsistencies, actually results in BIG values for
! fire residence time, thus lots of vegetation death!
! taul is the duration of the lethal heating.
! The /10 is to convert from kgC/m2 into gC/cm2, as in the Peterson and Ryan paper #Rosie,Jun 2013
do c = 1,nfsc
tau_b(c) = 39.4_r8 *(currentPatch%fuel_frac(c)*currentPatch%sum_fuel/0.45_r8/10._r8)* &
(1.0_r8-((1.0_r8-currentPatch%burnt_frac_litter(c))**0.5_r8))
enddo
tau_b(tr_sf) = 0.0_r8
! Cap the residence time to 8mins, as suggested by literature survey by P&R (1986).
currentPatch%tau_l = min(8.0_r8,sum(tau_b))
!---calculate overall fuel consumed by spreading fire ---
! ignore 1000hr fuels. Just interested in fuels affecting ROS
currentPatch%TFC_ROS = sum(FC_ground)-FC_ground(tr_sf)
currentPatch=>currentPatch%younger;
enddo !end patch loop
end subroutine ground_fuel_consumption