-
Notifications
You must be signed in to change notification settings - Fork 318
/
VOCEmissionMod.F90
1086 lines (931 loc) · 51.8 KB
/
VOCEmissionMod.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 VOCEmissionMod
!-----------------------------------------------------------------------
! !DESCRIPTION:
! Volatile organic compound emission
!
! !USES:
use shr_kind_mod , only : r8 => shr_kind_r8
use shr_log_mod , only : errMsg => shr_log_errMsg
use clm_varctl , only : iulog
use clm_varpar , only : maxveg, nlevcan
use pftconMod , only : ndllf_evr_tmp_tree, ndllf_evr_brl_tree
use pftconMod , only : ndllf_dcd_brl_tree, nbrdlf_evr_trp_tree
use pftconMod , only : nbrdlf_evr_tmp_tree, nbrdlf_dcd_brl_shrub
use pftconMod , only : nbrdlf_dcd_trp_tree, nbrdlf_dcd_tmp_tree
use pftconMod , only : nbrdlf_dcd_brl_tree, nbrdlf_evr_shrub
use pftconMod , only : nc3_arctic_grass , nc3crop
use pftconMod , only : nc4_grass, noveg
use shr_megan_mod , only : shr_megan_megcomps_n, shr_megan_megcomp_t, shr_megan_linkedlist
use shr_megan_mod , only : shr_megan_mechcomps_n, shr_megan_mechcomps, shr_megan_mapped_emisfctrs
use MEGANFactorsMod , only : Agro, Amat, Anew, Aold, betaT, ct1, ct2, LDF, Ceo
use decompMod , only : bounds_type
use abortutils , only : endrun
use fileutils , only : getfil
use clm_varcon , only : grlnd
use atm2lndType , only : atm2lnd_type
use CanopyStateType , only : canopystate_type
use PhotosynthesisMod , only : photosyns_type
use SoilStateType , only : soilstate_type
use SolarAbsorbedType , only : solarabs_type
use TemperatureType , only : temperature_type
use PatchType , only : patch
use EnergyFluxType , only : energyflux_type
!
implicit none
private
!
! !PUBLIC MEMBER FUNCTIONS:
public :: VOCEmission
!
! !PUBLIC TYPES:
type, public :: vocemis_type
real(r8) , pointer, private :: Eopt_out_patch (:) ! Eopt coefficient
real(r8) , pointer, private :: topt_out_patch (:) ! topt coefficient
real(r8) , pointer, private :: alpha_out_patch (:) ! alpha coefficient
real(r8) , pointer, private :: cp_out_patch (:) ! cp coefficient
real(r8) , pointer, private :: paru_out_patch (:) !
real(r8) , pointer, private :: par24u_out_patch (:) !
real(r8) , pointer, private :: par240u_out_patch (:) !
real(r8) , pointer, private :: para_out_patch (:) !
real(r8) , pointer, private :: par24a_out_patch (:) !
real(r8) , pointer, private :: par240a_out_patch (:) !
real(r8) , pointer, private :: gamma_out_patch (:) !
real(r8) , pointer, private :: gammaL_out_patch (:) !
real(r8) , pointer, private :: gammaT_out_patch (:) !
real(r8) , pointer, private :: gammaP_out_patch (:) !
real(r8) , pointer, private :: gammaA_out_patch (:) !
real(r8) , pointer, private :: gammaS_out_patch (:) !
real(r8) , pointer, private :: gammaC_out_patch (:) !
real(r8) , pointer, private :: vocflx_tot_patch (:) ! total VOC flux into atmosphere [moles/m2/sec]
real(r8) , pointer, PUBLIC :: vocflx_patch (:,:) ! (num_mech_comps) MEGAN flux [moles/m2/sec]
real(r8) , pointer, private :: efisop_grc (:,:) ! gridcell isoprene emission factors
contains
procedure, public :: Init
procedure, private :: InitAllocate
procedure, private :: InitHistory
procedure, private :: InitCold
end type vocemis_type
!
! !PRIVATE TYPES:
type :: megan_out_type
! VOC fluxes structure for CLM history output
real(r8), pointer, private :: flux_out(:) ! patch MEGAN flux [ug C m-2 h-1]
end type megan_out_type
type(megan_out_type), private, pointer :: meg_out(:) ! (n_megan_comps) points to output fluxes
!
logical, parameter :: debug = .false.
character(len=*), parameter, private :: sourcefile = &
__FILE__
!------------------------------------------------------------------------
contains
!------------------------------------------------------------------------
subroutine Init(this, bounds)
use clm_varctl , only : use_fates, use_fates_sp
class(vocemis_type) :: this
type(bounds_type), intent(in) :: bounds
if ( shr_megan_mechcomps_n > 0) then
if ( use_fates .and. (.not. use_fates_sp) ) then
call endrun( msg='ERROR: MEGAN currently does NOT work with FATES outside of FATES-SP mode (see github issue #115)'//&
errMsg(sourcefile, __LINE__))
end if
call this%InitAllocate(bounds)
call this%InitHistory(bounds)
call this%InitCold(bounds)
end if
end subroutine Init
!-----------------------------------------------------------------------
subroutine InitAllocate(this, bounds)
!
! Allocate memory for module datatypes
use shr_infnan_mod , only : nan => shr_infnan_nan, assignment(=)
use shr_megan_mod , only : shr_megan_factors_file
use MEGANFactorsMod , only : megan_factors_init, megan_factors_get
use clm_varpar , only : mxpft
!
! !ARGUMENTS:
class(vocemis_type) :: this
type(bounds_type) , intent(in) :: bounds
!
! !LOCAL VARIABLES:
integer :: i, imeg
integer :: class_num
real(r8) :: factors(mxpft+1)
real(r8) :: molec_wght
integer :: begg, endg
integer :: begp, endp
type(shr_megan_megcomp_t), pointer :: meg_cmp
!-----------------------------------------------------------------------
begg = bounds%begg; endg = bounds%endg
begp = bounds%begp; endp = bounds%endp
call megan_factors_init( shr_megan_factors_file )
meg_cmp => shr_megan_linkedlist
do while(associated(meg_cmp))
allocate(meg_cmp%emis_factors(maxveg))
call megan_factors_get( trim(meg_cmp%name), factors, class_num, molec_wght )
meg_cmp%emis_factors(1:maxveg) = factors(1:maxveg)
meg_cmp%class_number = class_num
meg_cmp%molec_weight = molec_wght
meg_cmp => meg_cmp%next_megcomp
enddo
allocate(this%Eopt_out_patch (begp:endp)) ; this%EOPT_out_patch (:) = nan
allocate(this%topt_out_patch (begp:endp)) ; this%topt_out_patch (:) = nan
allocate(this%topt_out_patch (begp:endp)) ; this%Eopt_out_patch (:) = nan
allocate(this%alpha_out_patch (begp:endp)) ; this%alpha_out_patch (:) = nan
allocate(this%cp_out_patch (begp:endp)) ; this%cp_out_patch (:) = nan
allocate(this%para_out_patch (begp:endp)) ; this%para_out_patch (:) = nan
allocate(this%par24a_out_patch (begp:endp)) ; this%par24a_out_patch (:) = nan
allocate(this%par240a_out_patch (begp:endp)) ; this%par240a_out_patch (:) = nan
allocate(this%paru_out_patch (begp:endp)) ; this%paru_out_patch (:) = nan
allocate(this%par24u_out_patch (begp:endp)) ; this%par24u_out_patch (:) = nan
allocate(this%par240u_out_patch (begp:endp)) ; this%par240u_out_patch (:) = nan
allocate(this%gamma_out_patch (begp:endp)) ; this%gamma_out_patch (:) = nan
allocate(this%gammaL_out_patch (begp:endp)) ; this%gammaL_out_patch (:) = nan
allocate(this%gammaT_out_patch (begp:endp)) ; this%gammaT_out_patch (:) = nan
allocate(this%gammaP_out_patch (begp:endp)) ; this%gammaP_out_patch (:) = nan
allocate(this%gammaA_out_patch (begp:endp)) ; this%gammaA_out_patch (:) = nan
allocate(this%gammaS_out_patch (begp:endp)) ; this%gammaS_out_patch (:) = nan
allocate(this%gammaC_out_patch (begp:endp)) ; this%gammaC_out_patch (:) = nan
allocate(this%vocflx_tot_patch (begp:endp)); this%vocflx_tot_patch (:) = nan
allocate(this%efisop_grc (6,begg:endg)); this%efisop_grc (:,:) = nan
allocate(meg_out(shr_megan_megcomps_n))
do i=1,shr_megan_megcomps_n
allocate(meg_out(i)%flux_out(begp:endp))
meg_out(i)%flux_out(:) = 0._r8
end do
allocate(this%vocflx_patch(begp:endp,1:shr_megan_mechcomps_n))
this%vocflx_patch(:,1:shr_megan_mechcomps_n)= nan
end subroutine InitAllocate
!-----------------------------------------------------------------------
subroutine InitHistory(this, bounds)
!
! !DESCRIPTION:
! Initialize history output fields for MEGAN emissions diagnositics
!
! !USES
use clm_varcon , only : spval
use histFileMod , only : hist_addfld1d
!
! !ARGUMENTS:
class(vocemis_type) :: this
type(bounds_type), intent(in) :: bounds
!
! !LOCAL VARIABLES
integer :: imeg, ii
integer :: begp, endp
type(shr_megan_megcomp_t), pointer :: meg_cmp
!---------------------------------------------------------------------
begp = bounds%begp; endp = bounds%endp
if (shr_megan_megcomps_n>0) then
! loop over megan compounds
meg_cmp => shr_megan_linkedlist
do while(associated(meg_cmp))
imeg = meg_cmp%index
call hist_addfld1d ( fname='MEG_'//trim(meg_cmp%name), units='kg/m2/sec', &
avgflag='A', long_name='MEGAN flux', &
ptr_patch=meg_out(imeg)%flux_out, set_lake=0._r8, set_urb=0._r8 )
meg_cmp => meg_cmp%next_megcomp
enddo
this%vocflx_tot_patch(begp:endp)= spval
call hist_addfld1d (fname='VOCFLXT', units='moles/m2/sec', &
avgflag='A', long_name='total VOC flux into atmosphere', &
ptr_patch=this%vocflx_tot_patch, set_lake=0._r8, set_urb=0._r8, default='inactive')
this%gamma_out_patch(begp:endp) = spval
call hist_addfld1d (fname='GAMMA', units='non', &
avgflag='A', long_name='total gamma for VOC calc', &
ptr_patch=this%gamma_out_patch, set_lake=0._r8, default='inactive')
this%gammaL_out_patch(begp:endp) = spval
call hist_addfld1d (fname='GAMMAL', units='non', &
avgflag='A', long_name='gamma L for VOC calc', &
ptr_patch=this%gammaL_out_patch, set_lake=0._r8, default='inactive')
this%gammaT_out_patch(begp:endp) = spval
call hist_addfld1d (fname='GAMMAT', units='non', &
avgflag='A', long_name='gamma T for VOC calc', &
ptr_patch=this%gammaT_out_patch, set_lake=0._r8, default='inactive')
this%gammaP_out_patch(begp:endp) = spval
call hist_addfld1d (fname='GAMMAP', units='non', &
avgflag='A', long_name='gamma P for VOC calc', &
ptr_patch=this%gammaP_out_patch, set_lake=0._r8, default='inactive')
this%gammaA_out_patch(begp:endp) = spval
call hist_addfld1d (fname='GAMMAA', units='non', &
avgflag='A', long_name='gamma A for VOC calc', &
ptr_patch=this%gammaA_out_patch, set_lake=0._r8, default='inactive')
this%gammaS_out_patch(begp:endp) = spval
call hist_addfld1d (fname='GAMMAS', units='non', &
avgflag='A', long_name='gamma S for VOC calc', &
ptr_patch=this%gammaS_out_patch, set_lake=0._r8, default='inactive')
this%gammaC_out_patch(begp:endp) = spval
call hist_addfld1d (fname='GAMMAC', units='non', &
avgflag='A', long_name='gamma C for VOC calc', &
ptr_patch=this%gammaC_out_patch, set_lake=0._r8, default='inactive')
this%EOPT_out_patch(begp:endp) = spval
call hist_addfld1d (fname='EOPT', units='non', &
avgflag='A', long_name='Eopt coefficient for VOC calc', &
ptr_patch=this%Eopt_out_patch, set_lake=0._r8, default='inactive')
this%topt_out_patch(begp:endp) = spval
call hist_addfld1d (fname='TOPT', units='non', &
avgflag='A', long_name='topt coefficient for VOC calc', &
ptr_patch=this%topt_out_patch, set_lake=0._r8, default='inactive')
this%alpha_out_patch(begp:endp) = spval
call hist_addfld1d (fname='ALPHA', units='non', &
avgflag='A', long_name='alpha coefficient for VOC calc', &
ptr_patch=this%alpha_out_patch, set_lake=0._r8, default='inactive')
this%cp_out_patch(begp:endp) = spval
call hist_addfld1d (fname='currentPatch', units='non', &
avgflag='A', long_name='currentPatch coefficient for VOC calc', &
ptr_patch=this%cp_out_patch, set_lake=0._r8, default='inactive')
this%paru_out_patch(begp:endp) = spval
call hist_addfld1d (fname='PAR_sun', units='umol/m2/s', &
avgflag='A', long_name='sunlit PAR', &
ptr_patch=this%paru_out_patch, set_lake=0._r8, default='inactive')
this%par24u_out_patch(begp:endp) = spval
call hist_addfld1d (fname='PAR24_sun', units='umol/m2/s', &
avgflag='A', long_name='sunlit PAR (24 hrs)', &
ptr_patch=this%par24u_out_patch, set_lake=0._r8, default='inactive')
this%par240u_out_patch(begp:endp) = spval
call hist_addfld1d (fname='PAR240_sun', units='umol/m2/s', &
avgflag='A', long_name='sunlit PAR (240 hrs)', &
ptr_patch=this%par240u_out_patch, set_lake=0._r8, default='inactive')
this%para_out_patch(begp:endp) = spval
call hist_addfld1d (fname='PAR_shade', units='umol/m2/s', &
avgflag='A', long_name='shade PAR', &
ptr_patch=this%para_out_patch, set_lake=0._r8, default='inactive')
this%par24a_out_patch(begp:endp) = spval
call hist_addfld1d (fname='PAR24_shade', units='umol/m2/s', &
avgflag='A', long_name='shade PAR (24 hrs)', &
ptr_patch=this%par24a_out_patch, set_lake=0._r8, default='inactive')
this%par240a_out_patch(begp:endp) = spval
call hist_addfld1d (fname='PAR240_shade', units='umol/m2/s', &
avgflag='A', long_name='shade PAR (240 hrs)', &
ptr_patch=this%par240a_out_patch, set_lake=0._r8, default='inactive')
end if
end subroutine InitHistory
!-----------------------------------------------------------------------
subroutine InitCold(this, bounds)
!
! !DESCRIPTION:
! Initialize cold start conditions for module variables
!
! !USES
use ncdio_pio
use clm_varctl, only : fsurdat
!
! !ARGUMENTS:
class(vocemis_type) :: this
type(bounds_type), intent(in) :: bounds
!
! !LOCAL VARIABLES:
logical :: readvar
integer :: begg, endg
type(file_desc_t) :: ncid ! netcdf id
character(len=256) :: locfn ! local filename
real(r8) ,pointer :: temp_ef(:) ! read in - temporary EFs
!-----------------------------------------------------------------------
begg = bounds%begg; endg = bounds%endg
! Time constant
allocate(temp_ef(begg:endg))
call getfil (fsurdat, locfn, 0)
call ncd_pio_openfile (ncid, locfn, 0)
call ncd_io(ncid=ncid, varname='EF1_BTR', flag='read', data=temp_ef, dim1name=grlnd, readvar=readvar)
if (.not. readvar) then
call endrun(msg='iniTimeConst: errror reading EF1_BTR'//errMsg(sourcefile, __LINE__))
end if
this%efisop_grc(1,begg:endg)=temp_ef(begg:endg)
call ncd_io(ncid=ncid, varname='EF1_FET', flag='read', data=temp_ef, dim1name=grlnd, readvar=readvar)
if (.not. readvar) then
call endrun(msg='iniTimeConst: errror reading EF1_FET'//errMsg(sourcefile, __LINE__))
end if
this%efisop_grc(2,begg:endg)=temp_ef(begg:endg)
call ncd_io(ncid=ncid, varname='EF1_FDT', flag='read', data=temp_ef, dim1name=grlnd, readvar=readvar)
if (.not. readvar) then
call endrun(msg='iniTimeConst: errror reading EF1_FDT'//errMsg(sourcefile, __LINE__))
end if
this%efisop_grc(3,begg:endg)=temp_ef(begg:endg)
call ncd_io(ncid=ncid, varname='EF1_SHR', flag='read', data=temp_ef, dim1name=grlnd, readvar=readvar)
if (.not. readvar) then
call endrun(msg='iniTimeConst: errror reading EF1_SHR'//errMsg(sourcefile, __LINE__))
end if
this%efisop_grc(4,begg:endg)=temp_ef(begg:endg)
call ncd_io(ncid=ncid, varname='EF1_GRS', flag='read', data=temp_ef, dim1name=grlnd, readvar=readvar)
if (.not. readvar) then
call endrun(msg='iniTimeConst: errror reading EF1_GRS'//errMsg(sourcefile, __LINE__))
end if
this%efisop_grc(5,begg:endg)=temp_ef(begg:endg)
call ncd_io(ncid=ncid, varname='EF1_CRP', flag='read', data=temp_ef, dim1name=grlnd, readvar=readvar)
if (.not. readvar) then
call endrun(msg='iniTimeConst: errror reading EF1_CRP'//errMsg(sourcefile, __LINE__))
end if
this%efisop_grc(6,begg:endg)=temp_ef(begg:endg)
deallocate(temp_ef)
call ncd_pio_closefile(ncid)
end subroutine InitCold
!-----------------------------------------------------------------------
subroutine VOCEmission (bounds, num_soilp, filter_soilp, &
atm2lnd_inst, canopystate_inst, photosyns_inst, temperature_inst, &
vocemis_inst, energyflux_inst)
!
! ! NEW DESCRIPTION
! Volatile organic compound emission
! This code simulates volatile organic compound emissions following
! MEGAN (Model of Emissions of Gases and Aerosols from Nature) v2.1
! for 20 compound classes. The original description of this
! algorithm (for isoprene only) can be found in Guenther et al., 2006
! (we follow equations 2-9, 16-17, 20 for explicit canopy).
! The model scheme came be described as:
! E= epsilon * gamma * rho
! VOC flux (E) [ug m-2 h-1] is calculated from baseline emission
! factors (epsilon) [ug m-2 h-1] which are specified for each of the 16
! CLM Patches (in input file) OR in the case of isoprene, from
! mapped EFs for each PATCH which reflect species divergence of emissions,
! particularly in North America.
! The emission activity factor (gamma) [unitless] for includes
! dependence on PPFT, temperature, LAI, leaf age and soil moisture.
! For isoprene only we also include the effect of CO2 inhibition as
! described by Heald et al., 2009.
! The canopy environment constant was calculated offline for CLM+CAM at
! standard conditions.
! We assume that the escape efficiency (rho) here is unity following
! Guenther et al., 2006.
! A manuscript describing MEGAN 2.1 and the implementation in CLM is
! in preparation: Guenther, Heald et al., 2012
! Subroutine written to operate at the patch level.
!
! Input: <filename> to be read in with EFs and some parameters.
! Currently these are set in procedure init_EF_params
! Output: vocflx(shr_megan_mechcomps_n) !VOC flux [moles/m2/sec]
!
! !USES:
use subgridAveMod , only : p2g
!
! !ARGUMENTS:
type(bounds_type) , intent(in) :: bounds
integer , intent(in) :: num_soilp ! number of columns in soil patch filter
integer , intent(in) :: filter_soilp(num_soilp) ! patch filter for soil
type(atm2lnd_type) , intent(in) :: atm2lnd_inst
type(canopystate_type) , intent(in) :: canopystate_inst
type(photosyns_type) , intent(in) :: photosyns_inst
type(temperature_type) , intent(in) :: temperature_inst
type(vocemis_type) , intent(inout) :: vocemis_inst
type(energyflux_type) , intent(in) :: energyflux_inst
!
! !REVISION HISTORY:
! 4/29/11: Colette L. Heald: expand MEGAN to 20 compound classes
! 7 Feb 2012: Francis Vitt: Implemented capability to specify MEGAN emissions in namelist
! and read in MEGAN factors from file.
!
! !LOCAL VARIABLES:
integer :: fp,p,g,c ! indices
real(r8) :: epsilon ! emission factor [ug m-2 h-1]
real(r8) :: gamma ! activity factor (accounting for light, T, age, LAI conditions)
real(r8) :: gamma_p ! activity factor for PPFD
real(r8) :: gamma_l ! activity factor for PPFD & LAI
real(r8) :: gamma_t ! activity factor for temperature
real(r8) :: gamma_a ! activity factor for leaf age
real(r8) :: gamma_sm ! activity factor for soil moisture
real(r8) :: gamma_c ! activity factor for CO2 (only isoprene)
real(r8) :: par_sun ! temporary
real(r8) :: par24_sun ! temporary
real(r8) :: par240_sun ! temporary
real(r8) :: par_sha ! temporary
real(r8) :: par24_sha ! temporary
real(r8) :: par240_sha ! temporary
integer :: class_num, n_meg_comps, imech, imeg, ii
character(len=16) :: mech_name
type(shr_megan_megcomp_t), pointer :: meg_cmp
real(r8) :: cp, alpha, Eopt, topt ! for history output
real(r8) :: co2_ppmv
real(r8) :: vocflx_meg(shr_megan_megcomps_n)
! factor used convert MEGAN units [micro-grams/m2/hr] to CAM srf emis units [g/m2/sec]
real(r8), parameter :: megemis_units_factor = 1._r8/3600._r8/1.e6_r8
! real(r8) :: root_depth(0:maxveg) ! Root depth [m]
character(len=32), parameter :: subname = "VOCEmission"
!-----------------------------------------------------------------------
! ! root depth (m) (defined based on Zeng et al., 2001, cf Guenther 2006)
! root_depth(noveg) = 0._r8 ! bare-soil
! root_depth(ndllf_evr_tmp_tree:ndllf_evr_brl_tree) = 1.8_r8 ! evergreen tree
! root_depth(ndllf_dcd_brl_tree) = 2.0_r8 ! needleleaf deciduous boreal tree
! root_depth(nbrdlf_evr_trp_tree:nbrdlf_evr_tmp_tree) = 3.0_r8 ! broadleaf evergreen tree
! root_depth(nbrdlf_dcd_trp_tree:nbrdlf_dcd_brl_tree) = 2.0_r8 ! broadleaf deciduous tree
! root_depth(nbrdlf_evr_shrub:nbrdlf_dcd_brl_shrub) = 2.5_r8 ! shrub
! root_depth(nc3_arctic_grass:maxveg) = 1.5_r8 ! grass/crop
!
if ( shr_megan_mechcomps_n < 1) return
if ( nlevcan /= 1 )then
call endrun( subname//' error: can NOT work without nlevcan == 1' )
end if
associate( &
btran => energyflux_inst%btran_patch , & ! Input: [real(r8) (:) ] transpiration wetness factor (0 to 1)
forc_solad => atm2lnd_inst%forc_solad_downscaled_col, & ! Input: [real(r8) (:,:) ] direct beam radiation (visible only)
forc_solai => atm2lnd_inst%forc_solai_grc , & ! Input: [real(r8) (:,:) ] diffuse radiation (visible only)
forc_pbot => atm2lnd_inst%forc_pbot_downscaled_col , & ! Input: [real(r8) (:) ] downscaled atmospheric pressure (Pa)
forc_pco2 => atm2lnd_inst%forc_pco2_grc , & ! Input: [real(r8) (:) ] partial pressure co2 (Pa)
forc_solad24 => atm2lnd_inst%fsd24_patch , & ! Input: [real(r8) (:) ] direct beam radiation last 24hrs (visible only)
forc_solad240 => atm2lnd_inst%fsd240_patch , & ! Input: [real(r8) (:) ] direct beam radiation last 240hrs (visible only)
forc_solai24 => atm2lnd_inst%fsi24_patch , & ! Input: [real(r8) (:) ] diffuse radiation last 24hrs (visible only)
forc_solai240 => atm2lnd_inst%fsi240_patch , & ! Input: [real(r8) (:) ] diffuse radiation last 240hrs (visible only)
fsun => canopystate_inst%fsun_patch , & ! Input: [real(r8) (:) ] sunlit fraction of canopy
fsun24 => canopystate_inst%fsun24_patch , & ! Input: [real(r8) (:) ] sunlit fraction of canopy last 24 hrs
fsun240 => canopystate_inst%fsun240_patch , & ! Input: [real(r8) (:) ] sunlit fraction of canopy last 240 hrs
elai => canopystate_inst%elai_patch , & ! Input: [real(r8) (:) ] one-sided leaf area index with burying by snow
elai240 => canopystate_inst%elai240_patch , & ! Input: [real(r8) (:) ] one-sided leaf area index with burying by snow last 240 hrs
cisun_z => photosyns_inst%cisun_z_patch , & ! Input: [real(r8) (:,:) ] sunlit intracellular CO2 (Pa)
cisha_z => photosyns_inst%cisha_z_patch , & ! Input: [real(r8) (:,:) ] shaded intracellular CO2 (Pa)
t_veg => temperature_inst%t_veg_patch , & ! Input: [real(r8) (:) ] patch vegetation temperature (Kelvin)
t_veg24 => temperature_inst%t_veg24_patch , & ! Input: [real(r8) (:) ] avg patch vegetation temperature for last 24 hrs
t_veg240 => temperature_inst%t_veg240_patch , & ! Input: [real(r8) (:) ] avg patch vegetation temperature for last 240 hrs
Eopt_out => vocemis_inst%Eopt_out_patch , & ! Output: [real(r8) (:) ]
topt_out => vocemis_inst%topt_out_patch , & ! Output: [real(r8) (:) ]
alpha_out => vocemis_inst%alpha_out_patch , & ! Output: [real(r8) (:) ]
cp_out => vocemis_inst%cp_out_patch , & ! Output: [real(r8) (:) ]
paru_out => vocemis_inst%paru_out_patch , & ! Output: [real(r8) (:) ]
par24u_out => vocemis_inst%par24u_out_patch , & ! Output: [real(r8) (:) ]
par240u_out => vocemis_inst%par240u_out_patch , & ! Output: [real(r8) (:) ]
para_out => vocemis_inst%para_out_patch , & ! Output: [real(r8) (:) ]
par24a_out => vocemis_inst%par24a_out_patch , & ! Output: [real(r8) (:) ]
par240a_out => vocemis_inst%par240a_out_patch , & ! Output: [real(r8) (:) ]
gammaL_out => vocemis_inst%gammaL_out_patch , & ! Output: [real(r8) (:) ]
gammaT_out => vocemis_inst%gammaT_out_patch , & ! Output: [real(r8) (:) ]
gammaP_out => vocemis_inst%gammaP_out_patch , & ! Output: [real(r8) (:) ]
gammaA_out => vocemis_inst%gammaA_out_patch , & ! Output: [real(r8) (:) ]
gammaS_out => vocemis_inst%gammaS_out_patch , & ! Output: [real(r8) (:) ]
gammaC_out => vocemis_inst%gammaC_out_patch , & ! Output: [real(r8) (:) ]
gamma_out => vocemis_inst%gamma_out_patch , & ! Output: [real(r8) (:) ]
vocflx => vocemis_inst%vocflx_patch , & ! Output: [real(r8) (:,:) ] VOC flux [moles/m2/sec]
vocflx_tot => vocemis_inst%vocflx_tot_patch & ! Output: [real(r8) (:) ] VOC flux [moles/m2/sec]
)
! initialize variables which get passed to the atmosphere
vocflx(bounds%begp:bounds%endp,:) = 0._r8
vocflx_tot(bounds%begp:bounds%endp) = 0._r8
do imeg=1,shr_megan_megcomps_n
meg_out(imeg)%flux_out(bounds%begp:bounds%endp) = 0._r8
enddo
! Begin loop over points
!_______________________________________________________________________________
do fp = 1,num_soilp
p = filter_soilp(fp)
g = patch%gridcell(p)
c = patch%column(p)
! initialize EF
epsilon=0._r8
! initalize to zero since this might not alway get set
! this needs to be within the fp loop ...
vocflx_meg(:) = 0._r8
! calculate VOC emissions for non-bare ground Patches
if (patch%itype(p) > 0) then
gamma=0._r8
! Calculate PAR: multiply w/m2 by 4.6 to get umol/m2/s for par (added 8/14/02)
!------------------------
! SUN:
par_sun = (forc_solad(c,1) + fsun(p) * forc_solai(g,1)) * 4.6_r8
par24_sun = (forc_solad24(p) + fsun24(p) * forc_solai24(p)) * 4.6_r8
par240_sun = (forc_solad240(p) + fsun240(p) * forc_solai240(p)) * 4.6_r8
! SHADE:
par_sha = ((1._r8 - fsun(p)) * forc_solai(g,1)) * 4.6_r8
par24_sha = ((1._r8 - fsun24(p)) * forc_solai24(p)) * 4.6_r8
par240_sha = ((1._r8 - fsun240(p)) * forc_solai240(p)) * 4.6_r8
! Activity factor for LAI (Guenther et al., 2006): all species
gamma_l = get_gamma_L(fsun240(p), elai(p))
! Impact of soil moisture on isoprene emission
gamma_sm = get_gamma_SM(btran(p))
! Loop through VOCs for light, temperature and leaf age activity factor & apply
! all final activity factors to baseline emission factors
!_______________________________________________________________________________
! loop over megan compounds
meg_cmp => shr_megan_linkedlist
meg_cmp_loop: do while(associated(meg_cmp))
imeg = meg_cmp%index
! set emis factor
! if specified, set EF for isoprene with mapped values
if ( trim(meg_cmp%name) == 'isoprene' .and. shr_megan_mapped_emisfctrs) then
epsilon = get_map_EF(patch%itype(p),g, vocemis_inst)
else
epsilon = meg_cmp%emis_factors(patch%itype(p))
end if
class_num = meg_cmp%class_number
! Activity factor for PPFD
gamma_p = get_gamma_P(par_sun, par24_sun, par240_sun, par_sha, par24_sha, par240_sha, &
fsun(p), fsun240(p), forc_solad240(p),forc_solai240(p), LDF(class_num), cp, alpha)
! Activity factor for T
gamma_t = get_gamma_T(t_veg240(p), t_veg24(p),t_veg(p), ct1(class_num), ct2(class_num),&
betaT(class_num),LDF(class_num), Ceo(class_num), Eopt, topt, patch%itype(p))
! Activity factor for Leaf Age
gamma_a = get_gamma_A(patch%itype(p), elai240(p),elai(p),class_num)
! Activity factor for CO2 (only for isoprene)
if (trim(meg_cmp%name) == 'isoprene') then
co2_ppmv = 1.e6_r8*forc_pco2(g)/forc_pbot(c)
gamma_c = get_gamma_C(cisun_z(p,1),cisha_z(p,1),forc_pbot(c),fsun(p), co2_ppmv)
else
gamma_c = 1._r8
end if
! Calculate total scaling factor
gamma = gamma_l * gamma_sm * gamma_a * gamma_p * gamma_T * gamma_c
if ( (gamma >=0.0_r8) .and. (gamma< 100._r8) ) then
vocflx_meg(imeg) = meg_cmp%coeff * epsilon * gamma * megemis_units_factor / meg_cmp%molec_weight ! moles/m2/sec
! assign to arrays for history file output (not weighted by landfrac)
meg_out(imeg)%flux_out(p) = meg_out(imeg)%flux_out(p) &
+ epsilon * gamma * megemis_units_factor*1.e-3_r8 ! Kg/m2/sec
if (imeg==1) then
!
gamma_out(p)=gamma
gammaP_out(p)=gamma_p
gammaT_out(p)=gamma_t
gammaA_out(p)=gamma_a
gammaS_out(p)=gamma_sm
gammaL_out(p)=gamma_l
gammaC_out(p)=gamma_c
paru_out(p)=par_sun
par24u_out(p)=par24_sun
par240u_out(p)=par240_sun
para_out(p)=par_sha
par24a_out(p)=par24_sha
par240a_out(p)=par240_sha
alpha_out(p)=alpha
cp_out(p)=cp
topt_out(p)=topt
Eopt_out(p)=Eopt
end if
endif
if (debug .and. gamma > 0.0_r8) then
write(iulog,*) 'MEGAN: n, megan name, epsilon, gamma, vocflx: ', &
imeg, meg_cmp%name, epsilon, gamma, vocflx_meg(imeg), gamma_p,gamma_t,gamma_a,gamma_sm,gamma_l
endif
meg_cmp => meg_cmp%next_megcomp
enddo meg_cmp_loop
! sum up the megan compound fluxes for the fluxes of chem mechanism compounds
do imech = 1,shr_megan_mechcomps_n
n_meg_comps = shr_megan_mechcomps(imech)%n_megan_comps
do imeg = 1,n_meg_comps ! loop over number of megan compounds that make up the nth mechanism compoud
ii = shr_megan_mechcomps(imech)%megan_comps(imeg)%ptr%index
vocflx(p,imech) = vocflx(p,imech) + vocflx_meg(ii)
enddo
vocflx_tot(p) = vocflx_tot(p) + vocflx(p,imech) ! moles/m2/sec
enddo
end if ! patch%itype(1:15 only)
enddo ! fp
end associate
end subroutine VOCEmission
!-----------------------------------------------------------------------
function get_map_EF(ivt_in, g_in, vocemis_inst)
!
! Get mapped EF for isoprene
! Use gridded values for 6 Patches specified by MEGAN following
! Guenther et al. (2006). Map the maxveg CLM Patches to these 6.
! Units: [ug m-2 h-1]
!
! !ARGUMENTS:
integer, intent(in) :: ivt_in
integer, intent(in) :: g_in
type(vocemis_type), intent(in) :: vocemis_inst
!
! !LOCAL VARIABLES:
real(r8) :: get_map_EF
!-----------------------------------------------------------------------
! vocemis_inst%efisop_patch ! Output: [real(r8) (:,:)] emission factors for isoprene for each patch [ug m-2 h-1]
get_map_EF = 0._r8
if ( ivt_in == ndllf_evr_tmp_tree &
.or. ivt_in == ndllf_evr_brl_tree) then !fineleaf evergreen
get_map_EF = vocemis_inst%efisop_grc(2,g_in)
else if (ivt_in == ndllf_dcd_brl_tree) then !fineleaf deciduous
get_map_EF = vocemis_inst%efisop_grc(3,g_in)
else if (ivt_in >= nbrdlf_evr_trp_tree &
.and. ivt_in <= nbrdlf_dcd_brl_tree) then !broadleaf trees
get_map_EF = vocemis_inst%efisop_grc(1,g_in)
else if (ivt_in >= nbrdlf_evr_shrub &
.and. ivt_in <= nbrdlf_dcd_brl_shrub) then !shrubs
get_map_EF = vocemis_inst%efisop_grc(4,g_in)
else if (ivt_in >= nc3_arctic_grass &
.and. ivt_in <= nc4_grass) then !grass
get_map_EF = vocemis_inst%efisop_grc(5,g_in)
else if (ivt_in >= nc3crop) then !crops
get_map_EF = vocemis_inst%efisop_grc(6,g_in)
end if
end function get_map_EF
!-----------------------------------------------------------------------
function get_gamma_P(par_sun_in, par24_sun_in, par240_sun_in, par_sha_in, par24_sha_in, par240_sha_in, &
fsun_in, fsun240_in, forc_solad240_in,forc_solai240_in, LDF_in, cp, alpha)
!
! Activity factor for PPFD (Guenther et al., 2006): all light dependent species
!-------------------------
! With distinction between sunlit and shaded leafs, weight scalings by
! fsun and fshade
! Scale total incident par by fraction of sunlit leaves (added on 1/2002)
! fvitt -- forc_solad240, forc_solai240 can be zero when CLM finidat is specified
! which will cause par240 to be zero and produce NaNs via log(par240)
! dml -- fsun240 can be equal to or greater than one before 10 day averages are
! set on startup or if a new patch comes online during land cover change.
! Avoid this problem by only doing calculations with fsun240 when fsun240 is
! between 0 and 1
!
! !ARGUMENTS:
implicit none
real(r8),intent(in) :: par_sun_in
real(r8),intent(in) :: par24_sun_in
real(r8),intent(in) :: par240_sun_in
real(r8),intent(in) :: par_sha_in
real(r8),intent(in) :: par24_sha_in
real(r8),intent(in) :: par240_sha_in
real(r8),intent(in) :: fsun_in
real(r8),intent(in) :: fsun240_in
real(r8),intent(in) :: forc_solad240_in
real(r8),intent(in) :: forc_solai240_in
real(r8),intent(in) :: LDF_in
real(r8),intent(out):: cp ! temporary
real(r8),intent(out):: alpha ! temporary
!
! !LOCAL VARIABLES:
real(r8) :: gamma_p_LDF ! activity factor for PPFD
real(r8) :: get_gamma_P ! return value
real(r8), parameter :: ca1 = 0.004_r8 ! empirical coefficent for alpha
real(r8), parameter :: ca2 = 0.0005_r8 ! empirical coefficent for alpha
real(r8), parameter :: ca3 = 0.0468_r8 ! empirical coefficent for cp
real(r8), parameter :: par0_sun = 200._r8 ! std conditions for past 24 hrs [umol/m2/s]
real(r8), parameter :: par0_shade = 50._r8 ! std conditions for past 24 hrs [umol/m2/s]
real(r8), parameter :: alpha_fix = 0.001_r8 ! empirical coefficient
real(r8), parameter :: cp_fix = 1.21_r8 ! empirical coefficient
!-----------------------------------------------------------------------
if ( (fsun240_in > 0._r8) .and. (fsun240_in < 1._r8) .and. (forc_solad240_in > 0._r8) &
.and. (forc_solai240_in > 0._r8)) then
! With alpha and cp calculated based on eq 6 and 7:
! Note indexing for accumulated variables is all at patch level
! SUN:
alpha = ca1 - ca2 * log(par240_sun_in)
cp = ca3 * exp(ca2 * (par24_sun_in-par0_sun))*par240_sun_in**(0.6_r8)
gamma_p_LDF = fsun_in * ( cp * alpha * par_sun_in * (1._r8 + alpha*alpha*par_sun_in*par_sun_in)**(-0.5_r8) )
! SHADE:
alpha = ca1 - ca2 * log(par240_sha_in)
cp = ca3 * exp(ca2 * (par_sha_in-par0_shade))*par240_sha_in**(0.6_r8)
gamma_p_LDF = gamma_p_LDF + (1._r8-fsun_in) * (cp*alpha*par_sha_in*(1._r8 + alpha*alpha*par_sha_in*par_sha_in)**(-0.5_r8))
else
! With fixed alpha and cp (from MEGAN User's Guide):
! SUN: direct + diffuse
alpha = alpha_fix
cp = cp_fix
gamma_p_LDF = fsun_in * ( cp * alpha*par_sun_in * (1._r8 + alpha*alpha*par_sun_in*par_sun_in)**(-0.5_r8) )
! SHADE: diffuse
gamma_p_LDF = gamma_p_LDF + (1._r8-fsun_in) * (cp*alpha*par_sha_in*(1._r8 + alpha*alpha*par_sha_in*par_sha_in)**(-0.5_r8))
end if
! Calculate total activity factor for PPFD accounting for light-dependent fraction
get_gamma_P = (1._r8 - LDF_in) + LDF_in * gamma_p_LDF
end function get_gamma_P
!-----------------------------------------------------------------------
function get_gamma_L(fsun240_in,elai_in)
!
! Activity factor for LAI (Guenther et al., 2006): all species
! Guenther et al., 2006 eq 3
!
! !USES:
use clm_varcon , only : denice
use clm_varpar , only : nlevsoi
!
! !ARGUMENTS:
implicit none
real(r8),intent(in) :: fsun240_in
real(r8),intent(in) :: elai_in
real(r8) :: get_gamma_L ! return value
!
! !LOCAL VARIABLES:
real(r8), parameter :: cce = 0.30_r8 ! factor to set emissions to unity @ std
real(r8), parameter :: cce1 = 0.24_r8 ! same as Cce but for non-accumulated vars
!-----------------------------------------------------------------------
if ( (fsun240_in > 0.0_r8) .and. (fsun240_in < 1.e30_r8) ) then
get_gamma_L = cce * elai_in
else
get_gamma_L = cce1 * elai_in
end if
end function get_gamma_L
!-----------------------------------------------------------------------
function get_gamma_SM(btran_in)
!---------------------------------------
! May 22, 2024
! Activity factor for soil moisture of Isoprene (Wang et al., 2022, JAMES)
! It is based on eq. (11) in the paper. Because the temperature response
! of isoprene has been explicitly included in CLM;
!ARGUMENTS:
implicit none
real(r8),intent(in) :: btran_in
!!!------- the drought algorithm--------
real(r8), parameter :: a1 = -7.4463_r8
real(r8), parameter :: b1 = 3.2552_r8
real(r8), parameter :: btran_threshold = 0.2_r8
real(r8) :: get_gamma_SM
!---------------------------------------
if (btran_in >= 1._r8) then
get_gamma_SM = 1._r8
else
get_gamma_SM = 1._r8 / (1._r8 + b1 * exp(a1 * (btran_in - btran_threshold)))
endif
end function get_gamma_SM
!-----------------------------------------------------------------------
function get_gamma_T(t_veg240_in, t_veg24_in,t_veg_in, ct1_in, ct2_in, betaT_in, LDF_in, Ceo_in, Eopt, topt, ivt_in)
! Activity factor for temperature
!--------------------------------
! May 24, 2024 Hui updated the temperature response curves of isoprene for
! Boreal Broadleaf Deciduous Shrub and Arctic C3 grass based on
! Wang et al., 2024 (GRL) and Wang et al., 2024 (Nature Communications)
!--------------------------------
! Calculate both a light-dependent fraction as in Guenther et al., 2006 for isoprene
! of a max saturation type form. Also caculate a light-independent fraction of the
! form of an exponential. Final activity factor depends on light dependent fraction
! of compound type.
!
! !USES:
use clm_varcon, only: tfrz
! !ARGUMENTS:
implicit none
integer,intent(in) :: ivt_in
real(r8),intent(in) :: t_veg240_in
real(r8),intent(in) :: t_veg24_in
real(r8),intent(in) :: t_veg_in
real(r8),intent(in) :: ct1_in
real(r8),intent(in) :: ct2_in
real(r8),intent(in) :: betaT_in
real(r8),intent(in) :: LDF_in
real(r8),intent(in) :: Ceo_in
real(r8),intent(out) :: Eopt ! temporary
real(r8),intent(out) :: topt ! temporary
!
! !LOCAL VARIABLES:
real(r8) :: get_gamma_T
real(r8) :: gamma_t_LDF ! activity factor for temperature
real(r8) :: gamma_t_LIF ! activity factor for temperature
real(r8) :: x ! temporary i
real(r8) :: bet_arc_c3 ! activity factor for temperature for arctic C3 grass
real(r8), parameter :: bet_arc_c3_max = 300._r8 ! max value, activity factor for temperature for arctic C3 graass
real(r8), parameter :: co1 = 313._r8 ! empirical coefficient
real(r8), parameter :: co2 = 0.6_r8 ! empirical coefficient
real(r8), parameter :: co4 = 0.05_r8 ! empirical coefficient
real(r8), parameter :: tstd0 = 297_r8 ! std temperature [K]
real(r8), parameter :: topt_fix = 317._r8 ! std temperature [K]
real(r8), parameter :: Eopt_fix = 2.26_r8 ! empirical coefficient
real(r8), parameter :: ct3 = 0.00831_r8 ! empirical coefficient (0.0083 in User's Guide)
real(r8), parameter :: tstd = 303.15_r8 ! std temperature [K]
real(r8), parameter :: bet = 0.09_r8 ! beta empirical coefficient [K-1]
real(r8), parameter :: std_act_energy_isopr = 95._r8 ! standard activation energy for isoprene
real(r8), parameter :: empirical_param_1 = 9.49_r8 ! empirical param for the activation energy in response to 10-day temperature change
real(r8), parameter :: empirical_param_2 = 0.53_r8 ! empirical param for the activation energy in response to 10-day temperature change
real(r8), parameter :: empirical_param_3 = 0.12_r8 ! empirical param for the emission factors of arctic C3 grass in response to 10-day temperature change
real(r8), parameter :: empirical_param_4 = 7.9_r8 ! empirical param for the emission factors of broadleaf deciduous boreal shrubs in response to 10-day temperature change
real(r8), parameter :: empirical_param_5 = 0.217_r8 ! empirical param for the emission factors of broadleaf deciduous boreal shrubs in response to 10-day temperature change
!-----------------------------------------------------------------------
! Light dependent fraction (Guenther et al., 2006)
if ( (t_veg240_in > 0.0_r8) .and. (t_veg240_in < 1.e30_r8) ) then
! topt and Eopt from eq 8 and 9:
topt = co1 + (co2 * (t_veg240_in-tstd0))
if ( (ivt_in == nbrdlf_dcd_brl_shrub) ) then ! boreal-deciduous-shrub
! coming from BEAR-oNS campaign willows results
Eopt = empirical_param_4 * exp (empirical_param_5 * (t_veg24_in - tfrz - 24.0_r8))
else if ( (ivt_in == nc3_arctic_grass ) ) then ! boreal-grass
Eopt = exp(empirical_param_3 * (t_veg240_in - tfrz - 15._r8))
else
Eopt = Ceo_in * exp (co4 * (t_veg24_in - tstd0)) * exp(co4 * (t_veg240_in - tstd0))
endif
else
topt = topt_fix
Eopt = Eopt_fix
endif
x = ( (1._r8/topt) - (1._r8/(t_veg_in)) ) / ct3
! for the boreal grass from BEAR-oNS campaign
if ( (ivt_in == nc3_arctic_grass ) ) then ! boreal-grass
bet_arc_c3 = std_act_energy_isopr + empirical_param_1 * exp(empirical_param_2 * (tfrz + 15._r8 - t_veg240_in))
bet_arc_c3 = min(bet_arc_c3, bet_arc_c3_max)
gamma_t_LDF = Eopt * exp(bet_arc_c3 * ((1._r8 / (tfrz + 30._r8) - 1._r8 / t_veg_in) / ct3))
else
gamma_t_LDF = Eopt * ( ct2_in * exp(ct1_in * x) / (ct2_in - ct1_in * (1._r8 - exp(ct2_in * x))) )
endif
! Light independent fraction (of exp(beta T) form)
gamma_t_LIF = exp(betaT_in * (t_veg_in - tstd))
! Calculate total activity factor for light as a function of light-dependent fraction
!--------------------------------
get_gamma_T = (1-LDF_in)*gamma_T_LIF + LDF_in*gamma_T_LDF
end function get_gamma_T
!-----------------------------------------------------------------------
function get_gamma_A(ivt_in, elai240_in, elai_in, nclass_in)
! Activity factor for leaf age (Guenther et al., 2006)
!-----------------------------
! If not CNDV elai is constant therefore gamma_a=1.0
! gamma_a set to unity for evergreens (Patches 1, 2, 4, 5)
! Note that we assume here that the time step is shorter than the number of
!days after budbreak required to induce isoprene emissions (ti=12 days) and
! the number of days after budbreak to reach peak emission (tm=28 days)
!
! !ARGUMENTS:
implicit none
integer,intent(in) :: ivt_in
integer,intent(in) :: nclass_in
real(r8),intent(in) :: elai240_in
real(r8),intent(in) :: elai_in
!
! !LOCAL VARIABLES:
real(r8) :: get_gamma_A
real(r8) :: elai_prev ! lai for previous timestep
real(r8) :: fnew, fgro, fmat, fold ! fractions of leaves at different phenological stages
!-----------------------------------------------------------------------
if ( (ivt_in == ndllf_dcd_brl_tree) .or. (ivt_in >= nbrdlf_dcd_trp_tree) ) then ! non-evergreen
if ( (elai240_in > 0.0_r8) .and. (elai240_in < 1.e30_r8) )then
elai_prev = 2._r8*elai240_in-elai_in ! have accumulated average lai over last 10 days
if (elai_prev == elai_in) then
fnew = 0.0_r8
fgro = 0.0_r8
fmat = 1.0_r8
fold = 0.0_r8
else if (elai_prev > elai_in) then
fnew = 0.0_r8
fgro = 0.0_r8
fmat = 1.0_r8 - (elai_prev - elai_in)/elai_prev
fold = (elai_prev - elai_in)/elai_prev
else if (elai_prev < elai_in) then
fnew = 1 - (elai_prev / elai_in)
fgro = 0.0_r8
fmat = (elai_prev / elai_in)
fold = 0.0_r8
end if
get_gamma_A = fnew*Anew(nclass_in) + fgro*Agro(nclass_in) + fmat*Amat(nclass_in) + fold*Aold(nclass_in)
else
get_gamma_A = 1.0_r8
end if
else
get_gamma_A = 1.0_r8
end if
end function get_gamma_A
!-----------------------------------------------------------------------
function get_gamma_C(cisun_in,cisha_in,forc_pbot_in,fsun_in, co2_ppmv)
! Activity factor for instantaneous CO2 changes (Heald et al., 2009)
!-------------------------
! With distinction between sunlit and shaded leaves, weight scalings by
! fsun and fshade
!
! !CALLED FROM: VOCEmission
!
! !REVISION HISTORY: