forked from ESCOMP/CTSM
-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathCNGapMortalityMod.F90
485 lines (423 loc) · 33.8 KB
/
CNGapMortalityMod.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
module CNGapMortalityMod
#include "shr_assert.h"
!-----------------------------------------------------------------------
! !DESCRIPTION:
! Module holding routines used in gap mortality for coupled carbon
! nitrogen code.
!
! !USES:
use shr_kind_mod , only : r8 => shr_kind_r8
use decompMod , only : bounds_type
use abortutils , only : endrun
use shr_log_mod , only : errMsg => shr_log_errMsg
use pftconMod , only : pftcon
use CNDVType , only : dgvs_type
use CNVegCarbonStateType , only : cnveg_carbonstate_type, spinup_factor_deadwood
use CNVegCarbonFluxType , only : cnveg_carbonflux_type
use CNVegNitrogenStateType , only : cnveg_nitrogenstate_type
use CNVegNitrogenFluxType , only : cnveg_nitrogenflux_type
use SoilBiogeochemNitrogenFluxType , only : soilbiogeochem_nitrogenflux_type
use CanopyStateType , only : canopystate_type
use ColumnType , only : col
use PatchType , only : patch
use GridcellType , only : grc
use CNSharedParamsMod , only : use_matrixcn
implicit none
private
!
! !PUBLIC MEMBER FUNCTIONS:
public :: readParams
public :: CNGapMortality
type, private :: params_type
real(r8):: am ! mortality rate based on annual rate, fractional mortality (1/yr)
real(r8):: k_mort ! coeff. of growth efficiency in mortality equation
end type params_type
!
type(params_type), private :: params_inst
!
! !PRIVATE MEMBER FUNCTIONS:
private :: CNGap_PatchToColumn
character(len=*), parameter, private :: sourcefile = &
__FILE__
!-----------------------------------------------------------------------
contains
!-----------------------------------------------------------------------
subroutine readParams ( ncid )
!
! !DESCRIPTION:
! Read in parameters
!
! !USES:
use ncdio_pio , only : file_desc_t,ncd_io
!
! !ARGUMENTS:
implicit none
type(file_desc_t),intent(inout) :: ncid ! pio netCDF file id
!
! !LOCAL VARIABLES:
character(len=32) :: subname = 'CNGapMortParamsType'
character(len=100) :: errCode = '-Error reading in parameters file:'
logical :: readv ! has variable been read in or not
real(r8) :: tempr ! temporary to read in constant
character(len=100) :: tString ! temp. var for reading
!-----------------------------------------------------------------------
tString='r_mort'
call ncd_io(varname=trim(tString),data=tempr, flag='read', ncid=ncid, readvar=readv)
if ( .not. readv ) call endrun(msg=trim(errCode)//trim(tString)//errMsg(sourcefile, __LINE__))
params_inst%am=tempr
tString='k_mort'
call ncd_io(varname=trim(tString),data=tempr, flag='read', ncid=ncid, readvar=readv)
if ( .not. readv ) call endrun(msg=trim(errCode)//trim(tString)//errMsg(sourcefile, __LINE__))
params_inst%k_mort=tempr
end subroutine readParams
!-----------------------------------------------------------------------
subroutine CNGapMortality (bounds, num_soilp, filter_soilp, &
dgvs_inst, cnveg_carbonstate_inst, cnveg_nitrogenstate_inst, soilbiogeochem_nitrogenflux_inst,&
cnveg_carbonflux_inst, cnveg_nitrogenflux_inst, canopystate_inst, &
leaf_prof_patch, froot_prof_patch, croot_prof_patch, stem_prof_patch)
!
! !DESCRIPTION:
! Gap-phase mortality routine for coupled carbon-nitrogen code (CN)
!
! !USES:
use clm_time_manager , only: get_average_days_per_year, get_step_size
use clm_varpar , only: nlevdecomp_full
use clm_varcon , only: secspday
use clm_varctl , only: use_cndv, spinup_state
use pftconMod , only: npcropmin
!
! !ARGUMENTS:
type(bounds_type) , intent(in) :: bounds
integer , intent(in) :: num_soilp ! number of soil patches in filter
integer , intent(in) :: filter_soilp(:) ! patch filter for soil points
type(dgvs_type) , intent(inout) :: dgvs_inst
type(cnveg_carbonstate_type) , intent(in) :: cnveg_carbonstate_inst
type(cnveg_nitrogenstate_type) , intent(in) :: cnveg_nitrogenstate_inst
type(cnveg_carbonflux_type) , intent(inout) :: cnveg_carbonflux_inst
type(cnveg_nitrogenflux_type) , intent(inout) :: cnveg_nitrogenflux_inst
type(canopystate_type) , intent(in) :: canopystate_inst
type(soilbiogeochem_nitrogenflux_type) , intent(inout) :: soilbiogeochem_nitrogenflux_inst
real(r8) , intent(in) :: leaf_prof_patch(bounds%begp:,1:)
real(r8) , intent(in) :: froot_prof_patch(bounds%begp:,1:)
real(r8) , intent(in) :: croot_prof_patch(bounds%begp:,1:)
real(r8) , intent(in) :: stem_prof_patch(bounds%begp:,1:)
!
! !LOCAL VARIABLES:
integer :: p ! patch index
integer :: fp ! patch filter index
real(r8):: dt ! time step (sec)
real(r8):: am ! rate for fractional mortality (1/yr)
real(r8):: m ! rate for fractional mortality (1/s)
real(r8):: mort_max ! asymptotic max mortality rate (/yr)
real(r8):: k_mort = 0.3_r8 ! coeff of growth efficiency in mortality equation
logical,parameter :: matrixcheck_gm = .False. ! If matrix check should be done
!-----------------------------------------------------------------------
SHR_ASSERT_ALL_FL((ubound(leaf_prof_patch) == (/bounds%endp,nlevdecomp_full/)), sourcefile, __LINE__)
SHR_ASSERT_ALL_FL((ubound(froot_prof_patch) == (/bounds%endp,nlevdecomp_full/)), sourcefile, __LINE__)
SHR_ASSERT_ALL_FL((ubound(croot_prof_patch) == (/bounds%endp,nlevdecomp_full/)), sourcefile, __LINE__)
SHR_ASSERT_ALL_FL((ubound(stem_prof_patch) == (/bounds%endp,nlevdecomp_full/)), sourcefile, __LINE__)
associate( &
ivt => patch%itype , & ! Input: [integer (:) ] patch vegetation type
woody => pftcon%woody , & ! Input: binary flag for woody lifeform
greffic => dgvs_inst%greffic_patch , & ! Input: [real(r8) (:) ]
heatstress => dgvs_inst%heatstress_patch , & ! Input: [real(r8) (:) ]
leafcn => pftcon%leafcn , & ! Input: [real(r8) (:)] leaf C:N (gC/gN)
livewdcn => pftcon%livewdcn , & ! Input: [real(r8) (:)] live wood (phloem and ray parenchyma) C:N (gC/gN)
laisun => canopystate_inst%laisun_patch , & ! Input: [real(r8) (:) ] sunlit projected leaf area index
laisha => canopystate_inst%laisha_patch , & ! Input: [real(r8) (:) ] shaded projected leaf area index
nind => dgvs_inst%nind_patch & ! Output:[real(r8)(:)] number of individuals (#/m2) added by F. Li and S. Levis
)
dt = real( get_step_size(), r8 )
! set the mortality rate based on annual rate
am = params_inst%am
! set coeff of growth efficiency in mortality equation
k_mort = params_inst%k_mort
! patch loop
do fp = 1,num_soilp
p = filter_soilp(fp)
if (use_cndv) then
! Stress mortality from lpj's subr Mortality.
if (woody(ivt(p)) == 1._r8) then
if (ivt(p) == 8) then
mort_max = 0.03_r8 ! BDT boreal
else
mort_max = 0.01_r8 ! original value for all patches
end if
! heatstress and greffic calculated in Establishment once/yr
! Mortality rate inversely related to growth efficiency
! (Prentice et al 1993)
am = mort_max / (1._r8 + k_mort * greffic(p))
! Mortality rate inversely related to growth efficiency
! (Prentice et al 1993)
am = mort_max / (1._r8 + k_mort * greffic(p))
am = min(1._r8, am + heatstress(p))
else ! lpj didn't set this for grasses; cn does
! set the mortality rate based on annual rate
am = params_inst%am
end if
end if
m = am/(get_average_days_per_year() * secspday)
!------------------------------------------------------
! patch-level gap mortality carbon fluxes
!------------------------------------------------------
if(.not. use_matrixcn)then
! displayed pools
cnveg_carbonflux_inst%m_leafc_to_litter_patch(p) = cnveg_carbonstate_inst%leafc_patch(p) * m
cnveg_carbonflux_inst%m_frootc_to_litter_patch(p) = cnveg_carbonstate_inst%frootc_patch(p) * m
cnveg_carbonflux_inst%m_livestemc_to_litter_patch(p) = cnveg_carbonstate_inst%livestemc_patch(p) * m
cnveg_carbonflux_inst%m_livecrootc_to_litter_patch(p) = cnveg_carbonstate_inst%livecrootc_patch(p) * m
cnveg_carbonflux_inst%m_deadstemc_to_litter_patch(p) = cnveg_carbonstate_inst%deadstemc_patch(p) * m * spinup_factor_deadwood
cnveg_carbonflux_inst%m_deadcrootc_to_litter_patch(p) = cnveg_carbonstate_inst%deadcrootc_patch(p) * m * spinup_factor_deadwood
! storage pools
cnveg_carbonflux_inst%m_leafc_storage_to_litter_patch(p) = cnveg_carbonstate_inst%leafc_storage_patch(p) * m
cnveg_carbonflux_inst%m_frootc_storage_to_litter_patch(p) = cnveg_carbonstate_inst%frootc_storage_patch(p) * m
cnveg_carbonflux_inst%m_livestemc_storage_to_litter_patch(p) = cnveg_carbonstate_inst%livestemc_storage_patch(p) * m
cnveg_carbonflux_inst%m_deadstemc_storage_to_litter_patch(p) = cnveg_carbonstate_inst%deadstemc_storage_patch(p) * m
cnveg_carbonflux_inst%m_livecrootc_storage_to_litter_patch(p) = cnveg_carbonstate_inst%livecrootc_storage_patch(p) * m
cnveg_carbonflux_inst%m_deadcrootc_storage_to_litter_patch(p) = cnveg_carbonstate_inst%deadcrootc_storage_patch(p) * m
cnveg_carbonflux_inst%m_gresp_storage_to_litter_patch(p) = cnveg_carbonstate_inst%gresp_storage_patch(p) * m
! transfer pools
cnveg_carbonflux_inst%m_leafc_xfer_to_litter_patch(p) = cnveg_carbonstate_inst%leafc_xfer_patch(p) * m
cnveg_carbonflux_inst%m_frootc_xfer_to_litter_patch(p) = cnveg_carbonstate_inst%frootc_xfer_patch(p) * m
cnveg_carbonflux_inst%m_livestemc_xfer_to_litter_patch(p) = cnveg_carbonstate_inst%livestemc_xfer_patch(p) * m
cnveg_carbonflux_inst%m_deadstemc_xfer_to_litter_patch(p) = cnveg_carbonstate_inst%deadstemc_xfer_patch(p) * m
cnveg_carbonflux_inst%m_livecrootc_xfer_to_litter_patch(p) = cnveg_carbonstate_inst%livecrootc_xfer_patch(p) * m
cnveg_carbonflux_inst%m_deadcrootc_xfer_to_litter_patch(p) = cnveg_carbonstate_inst%deadcrootc_xfer_patch(p) * m
cnveg_carbonflux_inst%m_gresp_xfer_to_litter_patch(p) = cnveg_carbonstate_inst%gresp_xfer_patch(p) * m
else
! For the matrix solution the same mortality gets applied, but it may be limited by the matrix solution
! This could be unified, by not limiting matrix_update_gmc when use_matrixcn is true
! displayed pools
! storage pools
! transfer pools
end if !use_matrixcn
!------------------------------------------------------
! patch-level gap mortality nitrogen fluxes
!------------------------------------------------------
! displayed pools
if(.not. use_matrixcn)then
cnveg_nitrogenflux_inst%m_leafn_to_litter_patch(p) = cnveg_nitrogenstate_inst%leafn_patch(p) * m
cnveg_nitrogenflux_inst%m_frootn_to_litter_patch(p) = cnveg_nitrogenstate_inst%frootn_patch(p) * m
cnveg_nitrogenflux_inst%m_livestemn_to_litter_patch(p) = cnveg_nitrogenstate_inst%livestemn_patch(p) * m
cnveg_nitrogenflux_inst%m_livecrootn_to_litter_patch(p) = cnveg_nitrogenstate_inst%livecrootn_patch(p) * m
else
! For the matrix solution the same mortality gets applied, but it may be limited by the matrix solution
! This could be unified, by not limiting matrix_update_gmn when use_matrixcn is true
end if
if (spinup_state == 2 .and. .not. use_cndv) then !accelerate mortality of dead woody pools
if(.not. use_matrixcn)then
cnveg_nitrogenflux_inst%m_deadstemn_to_litter_patch(p) = cnveg_nitrogenstate_inst%deadstemn_patch(p) * m * spinup_factor_deadwood
cnveg_nitrogenflux_inst%m_deadcrootn_to_litter_patch(p) = cnveg_nitrogenstate_inst%deadcrootn_patch(p) * m * spinup_factor_deadwood
else
! For the matrix solution the same mortality gets applied, but it may be limited by the matrix solution
! This could be unified, by not limiting matrix_update_gmn when use_matrixcn is true
end if !.not. use_matrixcn
else
if (.not. use_matrixcn) then
cnveg_nitrogenflux_inst%m_deadstemn_to_litter_patch(p) = cnveg_nitrogenstate_inst%deadstemn_patch(p) * m
cnveg_nitrogenflux_inst%m_deadcrootn_to_litter_patch(p) = cnveg_nitrogenstate_inst%deadcrootn_patch(p) * m
else
! For the matrix solution the same mortality gets applied, but it may be limited by the matrix solution
! This could be unified, by not limiting matrix_update_gmn when use_matrixcn is true
end if !use_matrixcn
end if
if (ivt(p) < npcropmin) then
if(.not. use_matrixcn)then
cnveg_nitrogenflux_inst%m_retransn_to_litter_patch(p) = cnveg_nitrogenstate_inst%retransn_patch(p) * m
else
! For the matrix solution the same mortality gets applied, but it may be limited by the matrix solution
! This could be unified, by not limiting matrix_update_gmn when use_matrixcn is true
end if
end if
if(.not. use_matrixcn)then
! storage pools
cnveg_nitrogenflux_inst%m_leafn_storage_to_litter_patch(p) = cnveg_nitrogenstate_inst%leafn_storage_patch(p) * m
cnveg_nitrogenflux_inst%m_frootn_storage_to_litter_patch(p) = cnveg_nitrogenstate_inst%frootn_storage_patch(p) * m
cnveg_nitrogenflux_inst%m_livestemn_storage_to_litter_patch(p) = cnveg_nitrogenstate_inst%livestemn_storage_patch(p) * m
cnveg_nitrogenflux_inst%m_deadstemn_storage_to_litter_patch(p) = cnveg_nitrogenstate_inst%deadstemn_storage_patch(p) * m
cnveg_nitrogenflux_inst%m_livecrootn_storage_to_litter_patch(p) = cnveg_nitrogenstate_inst%livecrootn_storage_patch(p) * m
cnveg_nitrogenflux_inst%m_deadcrootn_storage_to_litter_patch(p) = cnveg_nitrogenstate_inst%deadcrootn_storage_patch(p) * m
! transfer pools
cnveg_nitrogenflux_inst%m_leafn_xfer_to_litter_patch(p) = cnveg_nitrogenstate_inst%leafn_xfer_patch(p) * m
cnveg_nitrogenflux_inst%m_frootn_xfer_to_litter_patch(p) = cnveg_nitrogenstate_inst%frootn_xfer_patch(p) * m
cnveg_nitrogenflux_inst%m_livestemn_xfer_to_litter_patch(p) = cnveg_nitrogenstate_inst%livestemn_xfer_patch(p) * m
cnveg_nitrogenflux_inst%m_deadstemn_xfer_to_litter_patch(p) = cnveg_nitrogenstate_inst%deadstemn_xfer_patch(p) * m
cnveg_nitrogenflux_inst%m_livecrootn_xfer_to_litter_patch(p) = cnveg_nitrogenstate_inst%livecrootn_xfer_patch(p) * m
cnveg_nitrogenflux_inst%m_deadcrootn_xfer_to_litter_patch(p) = cnveg_nitrogenstate_inst%deadcrootn_xfer_patch(p) * m
else
! For the matrix solution the same mortality gets applied, but it may be limited by the matrix solution
! This could be unified, by not limiting matrix_update_gmn when use_matrixcn is true
! storage pools
! transfer pools
end if !use_matrixcn
! added by F. Li and S. Levis
if (use_cndv) then
if (woody(ivt(p)) == 1._r8)then
if (cnveg_carbonstate_inst%livestemc_patch(p) + cnveg_carbonstate_inst%deadstemc_patch(p)> 0._r8)then
nind(p)=nind(p)*(1._r8-m)
else
nind(p) = 0._r8
end if
end if
end if
end do ! end of patch loop
! gather all patch-level litterfall fluxes to the column
! for litter C and N inputs
call CNGap_PatchToColumn(bounds, num_soilp, filter_soilp, &
cnveg_carbonflux_inst, cnveg_nitrogenflux_inst, &
leaf_prof_patch(bounds%begp:bounds%endp, 1:nlevdecomp_full), &
froot_prof_patch(bounds%begp:bounds%endp, 1:nlevdecomp_full), &
croot_prof_patch(bounds%begp:bounds%endp, 1:nlevdecomp_full), &
stem_prof_patch(bounds%begp:bounds%endp, 1:nlevdecomp_full))
end associate
end subroutine CNGapMortality
!-----------------------------------------------------------------------
subroutine CNGap_PatchToColumn (bounds, num_soilp, filter_soilp, &
cnveg_carbonflux_inst, cnveg_nitrogenflux_inst, &
leaf_prof_patch, froot_prof_patch, croot_prof_patch, stem_prof_patch)
!
! !DESCRIPTION:
! gathers all patch-level gap mortality fluxes to the column level and
! assigns them to the three litter pools
!
! !USES:
use clm_varpar , only : maxsoil_patches, nlevdecomp, nlevdecomp_full, i_litr_min, i_litr_max, i_met_lit
!
! !ARGUMENTS:
type(bounds_type) , intent(in) :: bounds
integer , intent(in) :: num_soilp ! number of soil patches in filter
integer , intent(in) :: filter_soilp(:) ! soil patch filter
type(cnveg_carbonflux_type) , intent(inout) :: cnveg_carbonflux_inst
type(cnveg_nitrogenflux_type) , intent(inout) :: cnveg_nitrogenflux_inst
real(r8) , intent(in) :: leaf_prof_patch(bounds%begp:,1:)
real(r8) , intent(in) :: froot_prof_patch(bounds%begp:,1:)
real(r8) , intent(in) :: croot_prof_patch(bounds%begp:,1:)
real(r8) , intent(in) :: stem_prof_patch(bounds%begp:,1:)
!
! !LOCAL VARIABLES:
integer :: fp,c,p,j,i ! indices
!-----------------------------------------------------------------------
SHR_ASSERT_ALL_FL((ubound(leaf_prof_patch) == (/bounds%endp,nlevdecomp_full/)), sourcefile, __LINE__)
SHR_ASSERT_ALL_FL((ubound(froot_prof_patch) == (/bounds%endp,nlevdecomp_full/)), sourcefile, __LINE__)
SHR_ASSERT_ALL_FL((ubound(croot_prof_patch) == (/bounds%endp,nlevdecomp_full/)), sourcefile, __LINE__)
SHR_ASSERT_ALL_FL((ubound(stem_prof_patch) == (/bounds%endp,nlevdecomp_full/)), sourcefile, __LINE__)
associate( &
leaf_prof => leaf_prof_patch , & ! Input: [real(r8) (:,:) ] (1/m) profile of leaves
froot_prof => froot_prof_patch , & ! Input: [real(r8) (:,:) ] (1/m) profile of fine roots
croot_prof => croot_prof_patch , & ! Input: [real(r8) (:,:) ] (1/m) profile of coarse roots
stem_prof => stem_prof_patch , & ! Input: [real(r8) (:,:) ] (1/m) profile of stems
ivt => patch%itype , & ! Input: [integer (:) ] patch vegetation type
wtcol => patch%wtcol , & ! Input: [real(r8) (:) ] patch weight relative to column (0-1)
lf_f => pftcon%lf_f , & ! Input: [real(r8) (:,:) ] leaf litter fractions
fr_f => pftcon%fr_f , & ! Input: [real(r8) (:,:) ] fine root litter fractions
m_leafc_to_litter => cnveg_carbonflux_inst%m_leafc_to_litter_patch , & ! Input: [real(r8) (:) ]
m_frootc_to_litter => cnveg_carbonflux_inst%m_frootc_to_litter_patch , & ! Input: [real(r8) (:) ]
m_livestemc_to_litter => cnveg_carbonflux_inst%m_livestemc_to_litter_patch , & ! Input: [real(r8) (:) ]
m_deadstemc_to_litter => cnveg_carbonflux_inst%m_deadstemc_to_litter_patch , & ! Input: [real(r8) (:) ]
m_livecrootc_to_litter => cnveg_carbonflux_inst%m_livecrootc_to_litter_patch , & ! Input: [real(r8) (:) ]
m_deadcrootc_to_litter => cnveg_carbonflux_inst%m_deadcrootc_to_litter_patch , & ! Input: [real(r8) (:) ]
m_leafc_storage_to_litter => cnveg_carbonflux_inst%m_leafc_storage_to_litter_patch , & ! Input: [real(r8) (:) ]
m_frootc_storage_to_litter => cnveg_carbonflux_inst%m_frootc_storage_to_litter_patch , & ! Input: [real(r8) (:) ]
m_livestemc_storage_to_litter => cnveg_carbonflux_inst%m_livestemc_storage_to_litter_patch , & ! Input: [real(r8) (:) ]
m_deadstemc_storage_to_litter => cnveg_carbonflux_inst%m_deadstemc_storage_to_litter_patch , & ! Input: [real(r8) (:) ]
m_livecrootc_storage_to_litter => cnveg_carbonflux_inst%m_livecrootc_storage_to_litter_patch , & ! Input: [real(r8) (:) ]
m_deadcrootc_storage_to_litter => cnveg_carbonflux_inst%m_deadcrootc_storage_to_litter_patch , & ! Input: [real(r8) (:) ]
m_gresp_storage_to_litter => cnveg_carbonflux_inst%m_gresp_storage_to_litter_patch , & ! Input: [real(r8) (:) ]
m_leafc_xfer_to_litter => cnveg_carbonflux_inst%m_leafc_xfer_to_litter_patch , & ! Input: [real(r8) (:) ]
m_frootc_xfer_to_litter => cnveg_carbonflux_inst%m_frootc_xfer_to_litter_patch , & ! Input: [real(r8) (:) ]
m_livestemc_xfer_to_litter => cnveg_carbonflux_inst%m_livestemc_xfer_to_litter_patch , & ! Input: [real(r8) (:) ]
m_deadstemc_xfer_to_litter => cnveg_carbonflux_inst%m_deadstemc_xfer_to_litter_patch , & ! Input: [real(r8) (:) ]
m_livecrootc_xfer_to_litter => cnveg_carbonflux_inst%m_livecrootc_xfer_to_litter_patch , & ! Input: [real(r8) (:) ]
m_deadcrootc_xfer_to_litter => cnveg_carbonflux_inst%m_deadcrootc_xfer_to_litter_patch , & ! Input: [real(r8) (:) ]
m_gresp_xfer_to_litter => cnveg_carbonflux_inst%m_gresp_xfer_to_litter_patch , & ! Input: [real(r8) (:) ]
gap_mortality_c_to_litr_c => cnveg_carbonflux_inst%gap_mortality_c_to_litr_c_col , & ! Output: [real(r8) (:,:,:) ] C fluxes associated with gap mortality to litter pools (gC/m3/s)
gap_mortality_c_to_cwdc => cnveg_carbonflux_inst%gap_mortality_c_to_cwdc_col , & ! Output: [real(r8) (:,:) ] C fluxes associated with gap mortality to CWD pool (gC/m3/s)
m_leafn_to_litter => cnveg_nitrogenflux_inst%m_leafn_to_litter_patch , & ! Input: [real(r8) (:) ]
m_frootn_to_litter => cnveg_nitrogenflux_inst%m_frootn_to_litter_patch , & ! Input: [real(r8) (:) ]
m_livestemn_to_litter => cnveg_nitrogenflux_inst%m_livestemn_to_litter_patch , & ! Input: [real(r8) (:) ]
m_deadstemn_to_litter => cnveg_nitrogenflux_inst%m_deadstemn_to_litter_patch , & ! Input: [real(r8) (:) ]
m_livecrootn_to_litter => cnveg_nitrogenflux_inst%m_livecrootn_to_litter_patch , & ! Input: [real(r8) (:) ]
m_deadcrootn_to_litter => cnveg_nitrogenflux_inst%m_deadcrootn_to_litter_patch , & ! Input: [real(r8) (:) ]
m_retransn_to_litter => cnveg_nitrogenflux_inst%m_retransn_to_litter_patch , & ! Input: [real(r8) (:) ]
m_leafn_storage_to_litter => cnveg_nitrogenflux_inst%m_leafn_storage_to_litter_patch , & ! Input: [real(r8) (:) ]
m_frootn_storage_to_litter => cnveg_nitrogenflux_inst%m_frootn_storage_to_litter_patch , & ! Input: [real(r8) (:) ]
m_livestemn_storage_to_litter => cnveg_nitrogenflux_inst%m_livestemn_storage_to_litter_patch , & ! Input: [real(r8) (:) ]
m_deadstemn_storage_to_litter => cnveg_nitrogenflux_inst%m_deadstemn_storage_to_litter_patch , & ! Input: [real(r8) (:) ]
m_livecrootn_storage_to_litter => cnveg_nitrogenflux_inst%m_livecrootn_storage_to_litter_patch , & ! Input: [real(r8) (:) ]
m_deadcrootn_storage_to_litter => cnveg_nitrogenflux_inst%m_deadcrootn_storage_to_litter_patch , & ! Input: [real(r8) (:) ]
m_leafn_xfer_to_litter => cnveg_nitrogenflux_inst%m_leafn_xfer_to_litter_patch , & ! Input: [real(r8) (:) ]
m_frootn_xfer_to_litter => cnveg_nitrogenflux_inst%m_frootn_xfer_to_litter_patch , & ! Input: [real(r8) (:) ]
m_livestemn_xfer_to_litter => cnveg_nitrogenflux_inst%m_livestemn_xfer_to_litter_patch , & ! Input: [real(r8) (:) ]
m_deadstemn_xfer_to_litter => cnveg_nitrogenflux_inst%m_deadstemn_xfer_to_litter_patch , & ! Input: [real(r8) (:) ]
m_livecrootn_xfer_to_litter => cnveg_nitrogenflux_inst%m_livecrootn_xfer_to_litter_patch , & ! Input: [real(r8) (:) ]
m_deadcrootn_xfer_to_litter => cnveg_nitrogenflux_inst%m_deadcrootn_xfer_to_litter_patch , & ! Input: [real(r8) (:) ]
gap_mortality_n_to_litr_n => cnveg_nitrogenflux_inst%gap_mortality_n_to_litr_n_col , & ! Output: [real(r8) (:,:,:)] N fluxes associated with gap mortality to litter pools (gN/m3/s)
gap_mortality_n_to_cwdn => cnveg_nitrogenflux_inst%gap_mortality_n_to_cwdn_col & ! Output: [real(r8) (:,:) ] N fluxes associated with gap mortality to CWD pool (gN/m3/s)
)
do j = 1,nlevdecomp
do fp = 1,num_soilp
p = filter_soilp(fp)
c = patch%column(p)
do i = i_litr_min, i_litr_max
gap_mortality_c_to_litr_c(c,j,i) = &
gap_mortality_c_to_litr_c(c,j,i) + &
! leaf gap mortality carbon fluxes
m_leafc_to_litter(p) * lf_f(ivt(p),i) * wtcol(p) * leaf_prof(p,j) + &
! fine root gap mortality carbon fluxes
m_frootc_to_litter(p) * fr_f(ivt(p),i) * wtcol(p) * froot_prof(p,j)
end do
! wood gap mortality carbon fluxes
gap_mortality_c_to_cwdc(c,j) = gap_mortality_c_to_cwdc(c,j) + &
(m_livestemc_to_litter(p) + m_deadstemc_to_litter(p)) * wtcol(p) * stem_prof(p,j)
gap_mortality_c_to_cwdc(c,j) = gap_mortality_c_to_cwdc(c,j) + &
(m_livecrootc_to_litter(p) + m_deadcrootc_to_litter(p)) * wtcol(p) * croot_prof(p,j)
! storage gap mortality carbon fluxes
! Metabolic litter is treated differently than other types
! of litter, so it gets this additional line after the
! most recent loop over all litter types
gap_mortality_c_to_litr_c(c,j,i_met_lit) = &
gap_mortality_c_to_litr_c(c,j,i_met_lit) + &
(m_leafc_storage_to_litter(p) + m_gresp_storage_to_litter(p)) * wtcol(p) * leaf_prof(p,j) + &
m_frootc_storage_to_litter(p) * wtcol(p) * froot_prof(p,j) + &
(m_livestemc_storage_to_litter(p) + m_deadstemc_storage_to_litter(p)) * wtcol(p) * stem_prof(p,j) + &
(m_livecrootc_storage_to_litter(p) + m_deadcrootc_storage_to_litter(p)) * wtcol(p) * croot_prof(p,j) + &
! transfer gap mortality carbon fluxes
(m_leafc_xfer_to_litter(p) + m_gresp_xfer_to_litter(p)) * wtcol(p) * leaf_prof(p,j) + &
m_frootc_xfer_to_litter(p) * wtcol(p) * froot_prof(p,j) + &
(m_livestemc_xfer_to_litter(p) + m_deadstemc_xfer_to_litter(p)) * wtcol(p) * stem_prof(p,j) + &
(m_livecrootc_xfer_to_litter(p) + m_deadcrootc_xfer_to_litter(p)) * wtcol(p) * croot_prof(p,j)
do i = i_litr_min, i_litr_max
gap_mortality_n_to_litr_n(c,j,i) = &
gap_mortality_n_to_litr_n(c,j,i) + &
! leaf gap mortality nitrogen fluxes
m_leafn_to_litter(p) * lf_f(ivt(p),i) * wtcol(p) * leaf_prof(p,j) + &
! fine root litter nitrogen fluxes
m_frootn_to_litter(p) * fr_f(ivt(p),i) * wtcol(p) * froot_prof(p,j)
end do
! wood gap mortality nitrogen fluxes
gap_mortality_n_to_cwdn(c,j) = gap_mortality_n_to_cwdn(c,j) + &
(m_livestemn_to_litter(p) + m_deadstemn_to_litter(p)) * wtcol(p) * stem_prof(p,j)
gap_mortality_n_to_cwdn(c,j) = gap_mortality_n_to_cwdn(c,j) + &
(m_livecrootn_to_litter(p) + m_deadcrootn_to_litter(p)) * wtcol(p) * croot_prof(p,j)
! Metabolic litter is treated differently than other types
! of litter, so it gets this additional line after the
! most recent loop over all litter types
gap_mortality_n_to_litr_n(c,j,i_met_lit) = &
gap_mortality_n_to_litr_n(c,j,i_met_lit) + &
! retranslocated N pool gap mortality fluxes
m_retransn_to_litter(p) * wtcol(p) * leaf_prof(p,j) + &
! storage gap mortality nitrogen fluxes
m_leafn_storage_to_litter(p) * wtcol(p) * leaf_prof(p,j) + &
m_frootn_storage_to_litter(p) * wtcol(p) * froot_prof(p,j) + &
(m_livestemn_storage_to_litter(p) + m_deadstemn_storage_to_litter(p)) * wtcol(p) * stem_prof(p,j) + &
(m_livecrootn_storage_to_litter(p) + m_deadcrootn_storage_to_litter(p)) * wtcol(p) * croot_prof(p,j) + &
! transfer gap mortality nitrogen fluxes
m_leafn_xfer_to_litter(p) * wtcol(p) * leaf_prof(p,j) + &
m_frootn_xfer_to_litter(p) * wtcol(p) * froot_prof(p,j) + &
(m_livestemn_xfer_to_litter(p) + m_deadstemn_xfer_to_litter(p)) * wtcol(p) * stem_prof(p,j) + &
(m_livecrootn_xfer_to_litter(p) + m_deadcrootn_xfer_to_litter(p)) * wtcol(p) * croot_prof(p,j)
end do
end do
end associate
end subroutine CNGap_PatchToColumn
end module CNGapMortalityMod