-
Notifications
You must be signed in to change notification settings - Fork 26
/
Copy pathmarbl_ciso_interior_tendency_mod.F90
1295 lines (1027 loc) · 64.4 KB
/
marbl_ciso_interior_tendency_mod.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 marbl_ciso_interior_tendency_mod
use marbl_kinds_mod, only : r8
use marbl_kinds_mod, only : int_kind
use marbl_kinds_mod, only : char_len
use marbl_constants_mod, only : c0
use marbl_constants_mod, only : c1
use marbl_constants_mod, only : c2
use marbl_constants_mod, only : c1000
use marbl_constants_mod, only : mpercm
use marbl_settings_mod, only : autotroph_cnt
use marbl_settings_mod, only : autotroph_settings
use marbl_settings_mod, only : ciso_on
use marbl_logging, only : marbl_log_type
use marbl_interface_public_types, only : marbl_diagnostics_type
use marbl_interface_public_types, only : marbl_domain_type
use marbl_interface_private_types, only : autotroph_local_type
use marbl_interface_private_types, only : column_sinking_particle_type
use marbl_interface_private_types, only : marbl_interior_tendency_share_type
use marbl_interface_private_types, only : zooplankton_share_type
use marbl_interface_private_types, only : marbl_particulate_share_type
use marbl_interface_private_types, only : marbl_tracer_index_type
use marbl_interface_private_types, only : autotroph_derived_terms_type
implicit none
private
public :: marbl_ciso_interior_tendency_compute
public :: marbl_ciso_interior_tendency_autotroph_zero_consistency_enforce
!-----------------------------------------------------------------------
! scalar constants for 14C decay calculation
!-----------------------------------------------------------------------
real (r8), parameter :: c14_halflife_years = 5730.0_r8 !C14 half file
real (r8) :: c14_lambda_inv_sec ! Decay variable in 1/seconds (convert 5730 year decay to seconds)
contains
!***********************************************************************
subroutine marbl_ciso_interior_tendency_compute( &
marbl_domain, &
interior_tendency_share, &
zooplankton_share, &
marbl_particulate_share, &
tracer_local, &
autotroph_local, &
autotroph_derived_terms, &
temperature, &
marbl_tracer_indices, &
interior_tendencies, &
marbl_interior_diags, &
marbl_status_log)
! Compute time derivatives for 13C and 14C state variables.
! 13C code is based on code from X. Giraud, ETH Zürich, 2008, for pop1
! Also added biotic 14C
use marbl_settings_mod, only : ciso_fract_factors
use marbl_settings_mod, only : f_graze_CaCO3_REMIN
use marbl_constants_mod, only : R13C_std
use marbl_constants_mod, only : R14C_std
use marbl_constants_mod, only : spd
use marbl_constants_mod, only : spy
use marbl_ciso_diagnostics_mod, only : store_diagnostics_ciso_interior
type(marbl_domain_type), intent(in) :: marbl_domain
type(marbl_interior_tendency_share_type), intent(in) :: interior_tendency_share
type(zooplankton_share_type), intent(in) :: zooplankton_share
type(marbl_particulate_share_type), intent(in) :: marbl_particulate_share
real (r8), intent(in) :: tracer_local(:,:)
type(autotroph_local_type), intent(in) :: autotroph_local
type(autotroph_derived_terms_type), intent(in) :: autotroph_derived_terms
real (r8), intent(in) :: temperature(:)
type(marbl_tracer_index_type), intent(in) :: marbl_tracer_indices
real (r8), intent(inout) :: interior_tendencies(:,:) ! computed source/sink terms (inout because we don't touch non-ciso tracers)
type(marbl_diagnostics_type), intent(inout) :: marbl_interior_diags
type(marbl_log_type), intent(inout) :: marbl_status_log
!-----------------------------------------------------------------------
! local variables
!-----------------------------------------------------------------------
character(len=*), parameter :: subname = 'marbl_ciso_interior_tendency_mod:marbl_ciso_interior_tendency_compute'
real (r8) :: work1 ! temporaries
integer (int_kind) :: &
n, & ! tracer index
auto_ind, & ! autotroph functional group index
k ! index for looping over k levels
real (r8), dimension(autotroph_cnt) :: &
cell_active_C_uptake, & ! ratio of active carbon uptake to carbon fixation
cell_surf, & ! surface areas of cells ( m^2 )
cell_carb_cont, & ! cell carbon content ( mol C cell^-1 )
cell_radius, & ! cell radius ( um )
cell_permea, & ! cell wall permeability to CO2(aq) (m/s)
cell_eps_fix ! fractionation effect of carbon fixation
real(r8), parameter :: &
eps_carb = -2.0_r8 ! eps_carb = d13C(CaCO3) - d13C(DIC) Ziveri et al., 2003
type(column_sinking_particle_type) :: &
PO13C, & ! base units = nmol 13C
PO14C, & ! base units = nmol 14C
P_Ca13CO3, & ! base units = nmol CaCO3 13C
P_Ca14CO3 ! base units = nmol CaCO3 14C
real (r8), dimension(marbl_domain%km) :: &
mui_to_co2star_loc ! local carbon autotroph instanteous growth rate over [CO2*] (m^3 /mol C /s)
real (r8), dimension(marbl_domain%km) :: &
R13C_CaCO3_form, & ! 13C/12C in CaCO3 production of small phyto
R13C_CO2STAR, & ! 13C/12C in CO2* water
R13C_DIC, & ! 13C/12C in total DIC
R13C_DOCtot, & ! 13C/12C in total DOCtot
R13C_zoototC, & ! 13C/12C in total zooplankton
R14C_CaCO3_form, & ! 14C/12C in CaCO3 production of small phyto
R14C_CO2STAR, & ! 14C/12C in CO2* water
R14C_DIC, & ! 14C/12C in total DIC
R14C_DOCtot, & ! 14C/12C in total DOCtot
R14C_zoototC ! 14C/12C in total zooplankton
real (r8), dimension(autotroph_cnt, marbl_domain%km) :: &
Ca13CO3_PROD, & ! prod. of 13C CaCO3 by small phyto (mmol Ca13CO3/m^3/sec)
Ca14CO3_PROD, & ! prod. of 14C CaCO3 by small phyto (mmol Ca14CO3/m^3/sec)
eps_autotroph, & ! Permil fractionation (or discrimination factor) for Carbon autotroph types sp, diat, diaz
mui_to_co2star, & ! Carbon autotroph instanteous growth rate over [CO2*] (m^3 /mmol C /s)
R13C_photoC, & ! 13C/12C in Carbon autotroph C-fixation (mmol 13C/m^3/sec)
R13C_autotroph, & ! 13C/12C in total small phytoplankton (mmol 13C/m^3/sec)
photo13C, & ! Carbon autotroph 13C-fixation (mmol 13C/m^3/sec)
photo14C, & ! Carbon autotroph 14C-fixation (mmol 14C/m^3/sec)
R14C_photoC, & ! 14C/12C in Carbon autotroph C-fixation (mmol 14C/m^3/sec)
R14C_autotroph, & ! 14C/12C in total small phytoplankton (mmol 14C/m^3/sec)
autotrophCaCO3_d13C, & ! d13C of autotrophCaCO3
autotrophCaCO3_d14C, & ! d14C of autotrophCaCO3
autotroph_d13C, & ! d13C of autotroph C
autotroph_d14C, & ! d14C of autotroph C
R13C_autotrophCaCO3, & ! 13C/12C in total small phytoplankton carbonate
R14C_autotrophCaCO3 ! 14C/12C in total small phytoplankton carbonate
real (r8), dimension(marbl_domain%km) :: &
frac_co3, & ! carbonate fraction fCO3 = [CO3--]/DIC
CO2STAR_int, & ! [CO2*] water (mmol/m^3) in interior domain (not only surface)
DO13Ctot_prod, & ! production of 13C DOCtot (mmol 13C/m^3/sec)
DO13Ctot_remin, & ! remineralization of 13C DOCtot (mmol 13C/m^3/sec)
eps_aq_g, & ! equilibrium fractionation (CO2_gaseous <-> CO2_aq)
eps_dic_g, & ! equilibrium fractionation between total DIC and gaseous CO2
alpha_aq_g, & ! eps = ( alpha -1 ) * 1000
alpha_dic_g, & ! eps = ( alpha -1 ) * 1000
delta_C13_Corg, & ! deltaC13 of Net Primary Production
delta_C13_CO2STAR, & ! deltaC13 of CO2*
DO14Ctot_prod, & ! production of 13C DOCtot (mmol 14C/m^3/sec)
DO14Ctot_remin, & ! remineralization of 13C DOCtot (mmol 14C/m^3/sec)
alpha_aq_g_14c, & ! alpha for 14C, with fractionation twice as large as for 13C
alpha_dic_g_14c, & ! alpha for 14C, with fractionation twice as large as for 13C
delta_C14_CO2STAR, & ! deltaC14 of CO2*
DIC_d13C, & ! d13C of DIC
DOCtot_d13C, & ! d13C of DOCtot
zoototC_d13C, & ! d13C of zoototC
DIC_d14C, & ! d14C of DIC
DOCtot_d14C, & ! d14C of DOCtot
zoototC_d14C, & ! d14C of zoototC
decay_14Ctot ! 14C decay loss term
!-------------------------------------------------------------
! Return immediately if not running with carbon isotope tracer module
if (.not. ciso_on) return
associate( &
column_km => marbl_domain%km, &
column_kmt => marbl_domain%kmt, &
CO3 => interior_tendency_share%CO3_fields, & ! INPUT carbonate ion
HCO3 => interior_tendency_share%HCO3_fields, & ! INPUT bicarbonate ion
H2CO3 => interior_tendency_share%H2CO3_fields, & ! INPUT carbonic acid
DOCtot_remin => interior_tendency_share%DOCtot_remin_fields, & ! INPUT remineralization of DOCtot (mmol C/m^3/sec)
DOCtot_loc => interior_tendency_share%DOCtot_loc_fields, & ! INPUT local copy of model DOCtot
DO13Ctot_loc => tracer_local(marbl_tracer_indices%DO13Ctot_ind,:), & ! local copy of model DO13Ctot
DO14Ctot_loc => tracer_local(marbl_tracer_indices%DO14Ctot_ind,:), & ! local copy of model DO14Ctot
DIC_loc => tracer_local(marbl_tracer_indices%DIC_ind,:), & ! INPUT local copy of model DIC
DI13C_loc => tracer_local(marbl_tracer_indices%DI13C_ind,:), & ! local copy of model DI13C
DI14C_loc => tracer_local(marbl_tracer_indices%DI14C_ind,:), & ! local copy of model DI14C
zootot13C_loc => tracer_local(marbl_tracer_indices%zootot13C_ind,:), & ! local copy of model zootot13C
zootot14C_loc => tracer_local(marbl_tracer_indices%zootot14C_ind,:), & ! local copy of model zootot14C
QCaCO3 => autotroph_derived_terms%QCaCO3, & ! INPUT small phyto CaCO3/C ratio (mmol CaCO3/mmol C)
auto_graze => autotroph_derived_terms%auto_graze, & ! INPUT autotroph grazing rate (mmol C/m^3/sec)
auto_graze_zoo => autotroph_derived_terms%auto_graze_zootot, & ! INPUT auto_graze routed to zoo (mmol C/m^3/sec)
auto_graze_poc => autotroph_derived_terms%auto_graze_poc, & ! INPUT auto_graze routed to poc (mmol C/m^3/sec)
auto_graze_doc => autotroph_derived_terms%auto_graze_doc, & ! INPUT auto_graze routed to doc (mmol C/m^3/sec)
auto_graze_dic => autotroph_derived_terms%auto_graze_dic, & ! INPUT auto_graze routed to dic (mmol C/m^3/sec)
auto_loss => autotroph_derived_terms%auto_loss, & ! INPUT autotroph non-grazing mort (mmol C/m^3/sec)
auto_loss_poc => autotroph_derived_terms%auto_loss_poc, & ! INPUT auto_loss routed to poc (mmol C/m^3/sec)
auto_loss_doc => autotroph_derived_terms%auto_loss_doc, & ! INPUT auto_loss routed to doc (mmol C/m^3/sec)
auto_loss_dic => autotroph_derived_terms%auto_loss_dic, & ! INPUT auto_loss routed to dic (mmol C/m^3/sec)
auto_agg => autotroph_derived_terms%auto_agg, & ! INPUT autotroph aggregation (mmol C/m^3/sec)
photoC => autotroph_derived_terms%photoC, & ! INPUT C-fixation (mmol C/m^3/sec)
CaCO3_form => autotroph_derived_terms%CaCO3_form, & ! INPUT prod. of CaCO3 by small phyto (mmol CaCO3/m^3/sec)
PCphoto => autotroph_derived_terms%PCphoto, & ! INPUT C-specific rate of photosynth. (1/sec)
zoototC_loc => zooplankton_share%zoototC_loc_fields(:), & ! INPUT local copy of model zoototC
zootot_loss => zooplankton_share%zootot_loss_fields(:), & ! INPUT mortality & higher trophic grazing on zooplankton (mmol C/m^3/sec)
zootot_loss_poc => zooplankton_share%zootot_loss_poc_fields(:), & ! INPUT zootot_loss routed to large detrital pool (mmol C/m^3/sec)
zootot_loss_doc => zooplankton_share%zootot_loss_doc_fields(:), & ! INPUT zootot_loss routed to doc (mmol C/m^3/sec)
zootot_loss_dic => zooplankton_share%zootot_loss_dic_fields(:), & ! INPUT zootot_loss routed to dic (mmol C/m^3/sec)
zootot_graze => zooplankton_share%zootot_graze_fields(:), & ! INPUT zooplankton losses due to grazing (mmol C/m^3/sec)
zootot_graze_zoo => zooplankton_share%zootot_graze_zootot_fields(:), & ! INPUT grazing of zooplankton routed to total zoo (mmol C/m^3/sec)
zootot_graze_poc => zooplankton_share%zootot_graze_poc_fields(:), & ! INPUT grazing of zooplankton routed to poc (mmol C/m^3/sec)
zootot_graze_doc => zooplankton_share%zootot_graze_doc_fields(:), & ! INPUT grazing of zooplankton routed to doc (mmol C/m^3/sec)
zootot_graze_dic => zooplankton_share%zootot_graze_dic_fields(:), & ! INPUT grazing of zooplankton routed to dic (mmol C/m^3/sec)
POC => marbl_particulate_share%POC, & ! INPUT
P_CaCO3 => marbl_particulate_share%P_CaCO3, & ! INPUT
di13c_ind => marbl_tracer_indices%di13c_ind, &
do13ctot_ind => marbl_tracer_indices%do13ctot_ind, &
zootot13C_ind => marbl_tracer_indices%zootot13C_ind, &
di14c_ind => marbl_tracer_indices%di14c_ind, &
do14ctot_ind => marbl_tracer_indices%do14ctot_ind, &
zootot14C_ind => marbl_tracer_indices%zootot14C_ind &
)
!-----------------------------------------------------------------------
! Allocate memory for column_sinking_particle data types
!-----------------------------------------------------------------------
call PO13C%construct(num_levels=column_km)
call PO14C%construct(num_levels=column_km)
call P_Ca13CO3%construct(num_levels=column_km)
call P_Ca14CO3%construct(num_levels=column_km)
!-----------------------------------------------------------------------
! Set module variables
!-----------------------------------------------------------------------
! Define decay variable for DI14C, using earlier defined half-life of 14C
c14_lambda_inv_sec = log(c2) / (c14_halflife_years * spy)
!----------------------------------------------------------------------------------------
! Set cell attributes
!----------------------------------------------------------------------------------------
call setup_cell_attributes(ciso_fract_factors, &
cell_active_C_uptake, cell_surf, cell_carb_cont, &
cell_radius, cell_permea, cell_eps_fix, marbl_status_log)
if (marbl_status_log%labort_marbl) then
call marbl_status_log%log_error_trace("setup_cell_attributes", subname)
return
end if
!-----------------------------------------------------------------------
! Initialize Particulate terms for k=1
!-----------------------------------------------------------------------
call set_surface_particulate_terms(POC_ciso=PO13C, P_CaCO3_ciso=P_Ca13CO3)
call set_surface_particulate_terms(POC_ciso=PO14C, P_CaCO3_ciso=P_Ca14CO3)
!-----------------------------------------------------------------------
! Set ratios
!-----------------------------------------------------------------------
do k = 1, column_km
!-----------------------------------------------------------------------
! set local 13C/12C ratios, assuming ecosystem carries 12C (C=C12+C13+C14)
! If any Carbon boxes are zero, set corresponding 13C to zeros.
! Calculate fraction of CO3
!-----------------------------------------------------------------------
if (DOCtot_loc(k) > c0) then
R13C_DOCtot(k) = DO13Ctot_loc(k) / DOCtot_loc(k)
R14C_DOCtot(k) = DO14Ctot_loc(k) / DOCtot_loc(k)
else
R13C_DOCtot(k) = c0
R14C_DOCtot(k) = c0
end if
if (DIC_loc(k) > c0) then
R13C_DIC(k) = DI13C_loc(k) / DIC_loc(k)
R14C_DIC(k) = DI14C_loc(k) / DIC_loc(k)
frac_co3(k) = CO3(k) / DIC_loc(k)
else
R13C_DIC(k) = c0
R14C_DIC(k) = c0
frac_co3(k) = c0
end if
if (zoototC_loc(k) > c0) then
R13C_zoototC(k) = zootot13C_loc(k) / zoototC_loc(k)
R14C_zoototC(k) = zootot14C_loc(k) / zoototC_loc(k)
else
R13C_zoototC(k) = c0
R14C_zoototC(k) = c0
end if
do auto_ind = 1, autotroph_cnt
if (autotroph_local%C(auto_ind,k) > c0) then
R13C_autotroph(auto_ind,k) = autotroph_local%C13(auto_ind,k) / autotroph_local%C(auto_ind,k)
R14C_autotroph(auto_ind,k) = autotroph_local%C14(auto_ind,k) / autotroph_local%C(auto_ind,k)
else
R13C_autotroph(auto_ind,k) = c0
R14C_autotroph(auto_ind,k) = c0
end if
if (marbl_tracer_indices%auto_inds(auto_ind)%CaCO3_ind > 0) then
if (autotroph_local%CaCO3(auto_ind,k) > c0) then
R13C_autotrophCaCO3(auto_ind,k) = autotroph_local%Ca13CO3(auto_ind,k) / autotroph_local%CaCO3(auto_ind,k)
R14C_autotrophCaCO3(auto_ind,k) = autotroph_local%Ca14CO3(auto_ind,k) / autotroph_local%CaCO3(auto_ind,k)
else
R13C_autotrophCaCO3(auto_ind,k) = c0
R14C_autotrophCaCO3(auto_ind,k) = c0
end if
end if
end do
!-----------------------------------------------------------------------
! discrimination factors of carbone chemistry based on
! Zhang et al, 1995, Geochim. et Cosmochim. Acta, 59 (1), 107-114
!
! eps = permil fractionation and alpha is the fractionation factor
! with eps =(alpha - 1) *1000
!
! Fractionation is twice as large for 14C compared to 13C
!-----------------------------------------------------------------------
eps_aq_g(k) = -0.0049_r8 * temperature(k) - 1.31_r8
eps_dic_g(k) = 0.014_r8 * temperature(k) * frac_co3(k) - 0.105_r8 * temperature(k) + 10.53_r8
alpha_aq_g(k) = c1 + eps_aq_g(k) / c1000
alpha_dic_g(k) = c1 + eps_dic_g(k) / c1000
!fractionation is twice as large for 14C compared to 13C
alpha_aq_g_14c(k) = c1 + eps_aq_g(k) * 2.0_r8 / c1000
alpha_dic_g_14c(k) = c1 + eps_dic_g(k) * 2.0_r8 / c1000
!-----------------------------------------------------------------------
! 13C/12C ratios of CO2* (CO2STAR)
!-----------------------------------------------------------------------
R13C_CO2STAR(k) = R13C_DIC(k) * alpha_aq_g(k) / alpha_dic_g(k)
!-----------------------------------------------------------------------
! delta_13C of CO2* (CO2STAR)
!-----------------------------------------------------------------------
delta_C13_CO2STAR(k) = ( R13C_CO2STAR(k) / R13C_std - c1 ) * c1000
!-----------------------------------------------------------------------
! 14C/12C ratios of CO2* (CO2STAR)
!-----------------------------------------------------------------------
R14C_CO2STAR(k) = R14C_DIC(k) * alpha_aq_g_14c(k) / alpha_dic_g_14c(k)
!-----------------------------------------------------------------------
! delta_14C of CO2* (CO2STAR)
!-----------------------------------------------------------------------
delta_C14_CO2STAR(k) = ( R14C_CO2STAR(k) / R14C_std - c1 ) * c1000
!-----------------------------------------------------------------------
! [CO2STAR] = [CO2*] = [CO2(aq)] + [H2CO3]
! (this is eq 1.1.1 in Zeebe and Wolf-Gladrow, CO2 in seawater:
! equilibrium, kinetics, isotopes, Elseview Oceanography Series 65)
!
! DIC= [CO3] + [HCO3] + [CO2*] (eq 1.1.7)
!
! => CO2STAR_int = DIC_loc - HCO3 - CO3 !
!-----------------------------------------------------------------------
CO2STAR_int(k) = DIC_loc(k) - HCO3(k) - CO3(k)
!------------------------------------------------------------------------
! Loop over autotrophe types sp, diat, diaz and calculate fractionation
! for each type
!------------------------------------------------------------------------
do auto_ind = 1, autotroph_cnt
!------------------------------------------------------------------------
! mu(i) / [CO2*] of small phytoplankton ( m^3 / mmol C /s )
!-----------------------------------------------------------------------
if ( CO2STAR_int(k) /= c0 ) then
mui_to_co2star(auto_ind,k) = PCphoto(auto_ind,k) / CO2STAR_int(k)
else
mui_to_co2star(auto_ind,k) = c0
end if
!-----------------------------------------------------------------------
! fractionation factors for 13C fixation against CO2* in
! authotrophe types (sp, diaz, diat)
!-----------------------------------------------------------------------
select case (ciso_fract_factors)
case ('Rau')
!-----------------------------------------------------------------------
! Rau et al., 1989 ( see Gruber et al., 1998 )
! with restriction between -19 and -31 permil (see Rau et al., 1989)
!-----------------------------------------------------------------------
delta_C13_Corg(k) = -0.8_r8 * CO2STAR_int(k) - 12.6_r8
delta_C13_Corg(k) = min( delta_C13_Corg(k) , -18.0_r8 )
delta_C13_Corg(k) = max( delta_C13_Corg(k) , -32.0_r8 )
eps_autotroph(auto_ind,k) = c1000 * (delta_C13_CO2STAR(k) - delta_C13_Corg(k)) &
/(c1000 + delta_C13_Corg(k))
case ('Laws')
!-----------------------------------------------------------------------
! Laws et al, 1995
! with restriction between 10 and 26 for size effect (Tagliabue and Bopp, 2008)
! convert mui_to_co2star from m3/mmol/s to kg/mumol/d
!-----------------------------------------------------------------------
eps_autotroph(auto_ind,k) = (( mui_to_co2star(auto_ind,k) * spd ) - 0.371_r8 ) / (-0.015_r8)
!--------------------------------------------------------------------------
! uncomment the following two lines to restrict eps_sp between 10 and 26
!--------------------------------------------------------------------------
! eps_autotroph(auto_ind,k) = min( eps_autotroph(auto_ind,k), 26.0_r8 )
! eps_autotroph(auto_ind,k) = max( eps_autotroph(auto_ind,k), 10.0_r8 )
case ('KellerMorel')
!-----------------------------------------------------------------------
! Keller and morel, 1999
!-----------------------------------------------------------------------
! convert mui_to_co2start from m3/mmol/s to m3/mol/s
mui_to_co2star_loc(k) = mui_to_co2star(auto_ind,k) * c1000
call fract_keller_morel( &
mui_to_co2star_loc(k), &
cell_active_C_uptake(auto_ind), &
cell_surf(auto_ind), &
cell_carb_cont(auto_ind), &
cell_radius(auto_ind), &
cell_permea(auto_ind), &
cell_eps_fix(auto_ind), &
eps_autotroph(auto_ind,k) )
end select
!-----------------------------------------------------------------------
if (eps_autotroph(auto_ind,k) /= -c1000 ) then
R13C_photoC(auto_ind,k) = R13C_CO2STAR(k) *c1000 / (eps_autotroph(auto_ind,k) + c1000)
R14C_photoC(auto_ind,k) = R14C_CO2STAR(k) *c1000 / (2.0_r8* eps_autotroph(auto_ind,k) + c1000)
else
R13C_photoC(auto_ind,k) = c1
R14C_photoC(auto_ind,k) = c1
end if
!-----------------------------------------------------------------------
! Use R13/14C_photoC to determine small phytoplankton, Diatom, and
! Diaztroph 13C and 14C fixation
!-----------------------------------------------------------------------
photo13C(auto_ind,k) = photoC(auto_ind,k) * R13C_photoC(auto_ind,k)
photo14C(auto_ind,k) = photoC(auto_ind,k) * R14C_photoC(auto_ind,k)
!-----------------------------------------------------------------------
! C13 & C14 CaCO3 production
!-----------------------------------------------------------------------
if (autotroph_settings(auto_ind)%imp_calcifier) then
R13C_CaCO3_form(k) = R13C_DIC(k) + R13C_std * eps_carb / c1000
R14C_CaCO3_form(k) = R14C_DIC(k) + R14C_std * eps_carb * 2.0_r8 / c1000
Ca13CO3_PROD(auto_ind,k) = CaCO3_form(auto_ind,k) * R13C_CaCO3_form(k)
Ca14CO3_PROD(auto_ind,k) = CaCO3_form(auto_ind,k) * R14C_CaCO3_form(k)
end if
end do ! end loop over auto_ind
!-----------------------------------------------------------------------
! compute terms for DO13Ctot and DO14Ctot
!-----------------------------------------------------------------------
DO13Ctot_prod(k) = &
(zootot_loss_doc(k) + zootot_graze_doc(k))*R13C_zoototC(k) + &
sum((auto_loss_doc(:,k) + auto_graze_doc(:,k)) * R13C_autotroph(:,k),dim=1)
DO14Ctot_prod(k) = &
(zootot_loss_doc(k) + zootot_graze_doc(k))*R14C_zoototC(k) + &
sum((auto_loss_doc(:,k) + auto_graze_doc(:,k)) * R14C_autotroph(:,k),dim=1)
DO13Ctot_remin(k) = DOCtot_remin(k) * R13C_DOCtot(k)
DO14Ctot_remin(k) = DOCtot_remin(k) * R14C_DOCtot(k)
!-----------------------------------------------------------------------
! large detritus 13C and 14C
!-----------------------------------------------------------------------
PO13C%prod(k) = &
(zootot_loss_poc(k) + zootot_graze_poc(k))*R13C_zoototC(k) + &
sum((auto_graze_poc(:,k) + auto_agg(:,k) + auto_loss_poc(:,k)) * R13C_autotroph(:,k),dim=1)
PO14C%prod(k) = &
(zootot_loss_poc(k) + zootot_graze_poc(k))*R14C_zoototC(k) + &
sum((auto_graze_poc(:,k) + auto_agg(:,k) + auto_loss_poc(:,k)) * R14C_autotroph(:,k),dim=1)
!-----------------------------------------------------------------------
! large detrital Ca13CO3 and Ca14CCO3
!-----------------------------------------------------------------------
do auto_ind = 1, autotroph_cnt
if (marbl_tracer_indices%auto_inds(auto_ind)%CaCO3_ind > 0) then
P_Ca13CO3%prod(k) = P_CaCO3%prod(k) * R13C_autotrophCaCO3(auto_ind,k)
P_Ca14CO3%prod(k) = P_CaCO3%prod(k) * R14C_autotrophCaCO3(auto_ind,k)
endif
end do
!-----------------------------------------------------------------------
! Calculate oceanic D14C and D13C of carbon pools
!-----------------------------------------------------------------------
DIC_d13C(k) = ( R13C_DIC(k) / R13C_std - c1 ) * c1000
DIC_d14C(k) = ( R14C_DIC(k) / R14C_std - c1 ) * c1000
DOCtot_d13C(k) = ( R13C_DOCtot(k) / R13C_std - c1 ) * c1000
DOCtot_d14C(k) = ( R14C_DOCtot(k) / R14C_std - c1 ) * c1000
zoototC_d13C(k)= ( R13C_zoototC(k) / R13C_std - c1 ) * c1000
zoototC_d14C(k)= ( R14C_zoototC(k) / R14C_std - c1 ) * c1000
do auto_ind = 1, autotroph_cnt
if (marbl_tracer_indices%auto_inds(auto_ind)%CaCO3_ind > 0) then
autotrophCaCO3_d13C(auto_ind,k) = ( R13C_autotrophCaCO3(auto_ind,k) / R13C_std - c1 ) * c1000
autotrophCaCO3_d14C(auto_ind,k) = ( R14C_autotrophCaCO3(auto_ind,k) / R14C_std - c1 ) * c1000
else
autotrophCaCO3_d13C(auto_ind,k) = c0
autotrophCaCO3_d14C(auto_ind,k) = c0
endif
autotroph_d13C(auto_ind,k) = ( R13C_autotroph(auto_ind,k) / R13C_std - c1 ) * c1000
autotroph_d14C(auto_ind,k) = ( R14C_autotroph(auto_ind,k) / R14C_std - c1 ) * c1000
end do ! end loop over auto_ind
!-----------------------------------------------------------------------
! Compute carbon isotope particulate terms
!-----------------------------------------------------------------------
call compute_particulate_terms(k, marbl_domain, tracer_local(:,k), marbl_tracer_indices, &
interior_tendency_share, marbl_particulate_share, PO13C, P_Ca13CO3)
call compute_particulate_terms(k, marbl_domain, tracer_local(:,k), marbl_tracer_indices, &
interior_tendency_share, marbl_particulate_share, PO14C, P_Ca14CO3)
!-----------------------------------------------------------------------
! Update interior_tendencies for the 7 carbon pools for each Carbon isotope
!-----------------------------------------------------------------------
decay_14Ctot(k) = c0
!-----------------------------------------------------------------------
! dtracrs: autotroph Carbon (3 carbon pools), autotroph Ca13CO3 and Ca14CO3
!-----------------------------------------------------------------------
do auto_ind = 1, autotroph_cnt
work1 = auto_graze(auto_ind,k) + auto_loss(auto_ind,k) + auto_agg(auto_ind,k)
n = marbl_tracer_indices%auto_inds(auto_ind)%C13_ind
interior_tendencies(n,k) = photo13C(auto_ind,k) - work1 * R13C_autotroph(auto_ind,k)
n = marbl_tracer_indices%auto_inds(auto_ind)%C14_ind
interior_tendencies(n,k) = photo14C(auto_ind,k) - work1 * R14C_autotroph(auto_ind,k) - &
c14_lambda_inv_sec * autotroph_local%C14(auto_ind,k)
decay_14Ctot(k) = decay_14Ctot(k) + c14_lambda_inv_sec * autotroph_local%C14(auto_ind,k)
n = marbl_tracer_indices%auto_inds(auto_ind)%Ca13CO3_ind
if (n > 0) then
interior_tendencies(n,k) = Ca13CO3_PROD(auto_ind,k) - QCaCO3(auto_ind,k) &
* work1 * R13C_autotrophCaCO3(auto_ind,k)
endif
n = marbl_tracer_indices%auto_inds(auto_ind)%Ca14CO3_ind
if (n > 0) then
interior_tendencies(n,k) = Ca14CO3_PROD(auto_ind,k) - QCaCO3(auto_ind,k) &
* work1 * R14C_autotrophCaCO3(auto_ind,k) &
- c14_lambda_inv_sec * autotroph_local%Ca14CO3(auto_ind,k)
decay_14Ctot(k) = decay_14Ctot(k) + c14_lambda_inv_sec * autotroph_local%Ca14CO3(auto_ind,k)
endif
end do
!-----------------------------------------------------------------------
! interior_tendencies: zoo 13 and 14 Carbon
!-----------------------------------------------------------------------
interior_tendencies(zootot13C_ind,k) = &
sum(auto_graze_zoo(:,k) * R13C_autotroph(:,k),dim=1) &
+ (zootot_graze_zoo(k) - zootot_graze(k) - zootot_loss(k)) &
* R13C_zoototC(k)
interior_tendencies(zootot14C_ind,k) = &
sum(auto_graze_zoo(:,k) * R14C_autotroph(:,k),dim=1) &
+ (zootot_graze_zoo(k) - zootot_graze(k) - zootot_loss(k)) &
* R14C_zoototC(k) - c14_lambda_inv_sec * zootot14C_loc(k)
decay_14Ctot(k) = decay_14Ctot(k) + c14_lambda_inv_sec * zootot14C_loc(k)
!-----------------------------------------------------------------------
! interior_tendencies: dissolved organic Matter 13C and 14C
!-----------------------------------------------------------------------
interior_tendencies(do13ctot_ind,k) = DO13Ctot_prod(k) - DO13Ctot_remin(k)
interior_tendencies(do14ctot_ind,k) = DO14Ctot_prod(k) - DO14Ctot_remin(k) - c14_lambda_inv_sec * DO14Ctot_loc(k)
decay_14Ctot(k) = decay_14Ctot(k) + c14_lambda_inv_sec * DO14Ctot_loc(k)
!-----------------------------------------------------------------------
! interior_tendencies: dissolved inorganic Carbon 13 and 14
!-----------------------------------------------------------------------
interior_tendencies(di13c_ind,k) = &
sum((auto_loss_dic(:,k) + auto_graze_dic(:,k))*R13C_autotroph(:,k),dim=1) &
- sum(photo13C(:,k),dim=1) &
+ DO13Ctot_remin(k) + PO13C%remin(k) &
+ (zootot_loss_dic(k) + zootot_graze_dic(k)) * R13C_zoototC(k) &
+ P_Ca13CO3%remin(k)
interior_tendencies(di14c_ind,k) = &
sum((auto_loss_dic(:,k) + auto_graze_dic(:,k))*R14C_autotroph(:,k),dim=1) &
- sum(photo14C(:,k),dim=1) &
+ DO14Ctot_remin(k) + PO14C%remin(k) &
+ (zootot_loss_dic(k) + zootot_graze_dic(k)) * R14C_zoototC(k) &
+ P_Ca14CO3%remin(k) &
- c14_lambda_inv_sec * DI14C_loc(k)
decay_14Ctot(k) = decay_14Ctot(k) + c14_lambda_inv_sec * DI14C_loc(k)
do auto_ind = 1, autotroph_cnt
if (marbl_tracer_indices%auto_inds(auto_ind)%Ca13CO3_ind > 0) then
interior_tendencies(di13c_ind,k) = interior_tendencies(di13c_ind,k) &
+ f_graze_CaCO3_REMIN * auto_graze(auto_ind,k) &
* QCaCO3(auto_ind,k) * R13C_autotrophCaCO3(auto_ind,k) &
- Ca13CO3_PROD(auto_ind,k)
endif
if (marbl_tracer_indices%auto_inds(auto_ind)%Ca14CO3_ind > 0) then
interior_tendencies(di14c_ind,k) = interior_tendencies(di14c_ind,k) &
+ f_graze_CaCO3_REMIN * auto_graze(auto_ind,k) &
* QCaCO3(auto_ind,k) * R14C_autotrophCaCO3(auto_ind,k) &
- Ca14CO3_PROD(auto_ind,k)
endif
end do
!-----------------------------------------------------------------------
! Update particulate terms from prior level for next level
!-----------------------------------------------------------------------
if (k < column_km) then
call update_particulate_terms_from_prior_level(k+1, PO13C, P_Ca13CO3)
call update_particulate_terms_from_prior_level(k+1, PO14C, P_Ca14CO3)
endif
end do ! end of loop over k levels
end associate
! update carbon isotope diagnostics
! FIXME #18: the following arguments need to be group into a derived type
call store_diagnostics_ciso_interior(&
marbl_domain, &
autotroph_d13C, &
autotroph_d14C, &
autotrophCaCO3_d13C, &
autotrophCaCO3_d14C, &
photo13C, &
photo14C, &
eps_autotroph, &
mui_to_co2star, &
Ca13CO3_prod, &
Ca14CO3_prod, &
DIC_d13C, &
DIC_d14C, &
DOCtot_d13C, &
DOCtot_d14C, &
zoototC_d13C, &
zoototC_d14C, &
DO13Ctot_prod, &
DO14Ctot_prod, &
DO13Ctot_remin, &
DO14Ctot_remin, &
eps_aq_g, &
eps_dic_g, &
decay_14Ctot, &
PO13C, &
PO14C, &
P_Ca13CO3, &
P_Ca14CO3, &
interior_tendencies, &
marbl_tracer_indices,&
marbl_interior_diags,&
marbl_status_log)
if (marbl_status_log%labort_marbl) then
call marbl_status_log%log_error_trace("store_diagnostics_ciso_interior", subname)
return
end if
!-----------------------------------------------------------------------
! Deallocate memory for column_sinking_particle data types
!-----------------------------------------------------------------------
call PO13C%destruct()
call PO14C%destruct()
call P_Ca13CO3%destruct()
call P_Ca14CO3%destruct()
end subroutine marbl_ciso_interior_tendency_compute
!***********************************************************************
subroutine marbl_ciso_interior_tendency_autotroph_zero_consistency_enforce(auto_ind, column_kmt, zero_mask, &
autotroph_tracer_indices, autotroph_local)
use marbl_interface_private_types, only : marbl_living_tracer_index_type
integer, intent(in) :: auto_ind
integer, intent(in) :: column_kmt
logical, intent(in) :: zero_mask(column_kmt)
type(marbl_living_tracer_index_type), intent(in) :: autotroph_tracer_indices
type(autotroph_local_type), intent(inout) :: autotroph_local
if (.not. ciso_on) return
where (zero_mask)
autotroph_local%C13(auto_ind,1:column_kmt) = c0
autotroph_local%C14(auto_ind,1:column_kmt) = c0
end where
if (autotroph_tracer_indices%Ca13CO3_ind > 0) then
where (zero_mask)
autotroph_local%Ca13CO3(auto_ind,1:column_kmt) = c0
end where
end if
if (autotroph_tracer_indices%Ca14CO3_ind > 0) then
where (zero_mask)
autotroph_local%Ca14CO3(auto_ind,1:column_kmt) = c0
end where
end if
end subroutine marbl_ciso_interior_tendency_autotroph_zero_consistency_enforce
!***********************************************************************
subroutine setup_cell_attributes(ciso_fract_factors, &
cell_active_C_uptake, cell_surf, cell_carb_cont, &
cell_radius, cell_permea, cell_eps_fix, marbl_status_log)
!----------------------------------------------------------------------------------------
! For Keller and Morel, set cell attributes based on autotroph type (from observations)
!----------------------------------------------------------------------------------------
character(len=char_len), intent(in) :: ciso_fract_factors ! option for which biological fractionation calculation to use
real (r8), intent(out) :: cell_active_C_uptake(autotroph_cnt) ! ratio of active carbon uptake to carbon fixation
real (r8), intent(out) :: cell_surf(autotroph_cnt) ! surface areas of cells ( m2 )
real (r8), intent(out) :: cell_carb_cont(autotroph_cnt) ! cell carbon content ( mol C cell-1 )
real (r8), intent(out) :: cell_radius(autotroph_cnt) ! cell radius ( um )
real (r8), intent(out) :: cell_permea(autotroph_cnt) ! cell wall permeability to CO2(aq) (m/s)
real (r8), intent(out) :: cell_eps_fix(autotroph_cnt) ! fractionation effect of carbon fixation
type(marbl_log_type), intent(inout) :: marbl_status_log
!-----------------------------------------------------------------------
! local variables
!-----------------------------------------------------------------------
character(len=*), parameter :: subname = 'marbl_ciso_interior_tendency_mod:setup_cell_attributes'
character(len=char_len) :: log_message
integer(int_kind) :: auto_ind ! autotroph functional group index
!-----------------------------------------------------------------------
select case (ciso_fract_factors)
case ('KellerMorel')
do auto_ind = 1, autotroph_cnt
if (autotroph_settings(auto_ind)%silicifier) then
!----------------------------------------------------------------------------------------
! Diatom based on P. tricornumtum ( Keller and morel, 1999; Popp et al., 1998 )
!----------------------------------------------------------------------------------------
cell_active_C_uptake(auto_ind) = 2.3_r8 ! ratio of active carbon uptake to carbon fixation
cell_surf(auto_ind) = 100.6e-12_r8 ! surface areas of cells ( m2 )
cell_carb_cont(auto_ind) = 63.3e-14_r8 ! cell carbon content ( mol C cell-1 )
cell_radius(auto_ind) = 14.2_r8 ! cell radius ( um )
cell_permea(auto_ind) = 3.3e-5_r8 ! cell wall permeability to CO2(aq) (m/s)
cell_eps_fix(auto_ind) = 26.6_r8 ! fractionation effect of carbon fixation
else if (autotroph_settings(auto_ind)%Nfixer) then
!----------------------------------------------------------------------------------------
! Diazotroph based on Standard Phyto of Rau et al., (1996)
!----------------------------------------------------------------------------------------
!cell_active_C_uptake(auto_ind) = 0.0_r8 ! ratio of active carbon uptake to carbon fixation
!cell_surf(auto_ind) = -99.9_r8 ! surface areas of cells ( m2 ) - not used -
!cell_carb_cont(auto_ind) = -99.9_r8 ! cell carbon content ( mol C cell-1 ) - not used -
!cell_radius(auto_ind) = 10.0_r8 ! cell radius ( um )
!cell_permea(auto_ind) = 10.0e-5_r8 ! cell wall permeability to CO2(aq) (m/s)
!cell_eps_fix(auto_ind) = 25.0_r8 ! fractionation effect of carbon fixation
!----------------------------------------------------------------------------------------
! Diazotroph based on Synechococcus sp. ( Keller and morel, 1999; Popp et al., 1998 )
!----------------------------------------------------------------------------------------
cell_active_C_uptake(auto_ind) = 7.5_r8 ! ratio of active carbon uptake to carbon fixation
cell_surf(auto_ind) = 5.8e-12_r8 ! surface areas of cells ( m2 )
cell_carb_cont(auto_ind) = 3e-14_r8 ! cell carbon content ( mol C cell-1 )
cell_radius(auto_ind) = 0.68_r8 ! cell radius ( um )
cell_permea(auto_ind) = 3.0e-8_r8 ! cell wall permeability to CO2(aq) (m/s)
cell_eps_fix(auto_ind) = 30.0_r8 ! fractionation effect of carbon fixation
!else if (autotroph_settings(auto_ind)%exp_calcifier) then
!Currently not set up to separate exp_calcifiers, needs cell_radius value from data
!----------------------------------------------------------------------------------------
! Calcifier based on P. glacialis ( Keller and morel, 1999; Popp et al., 1998 )
! Popp et al express cell carbon content in ( pg C cell-1 ), here we convert in (mol C cell-1)
! convert pgC to molC : ! Mc = 12 g / mol ! Mc = 12 e12 pg / mol
!----------------------------------------------------------------------------------------
! cell_active_C_uptake(auto_ind) = 0.0_r9 ! ratio of active carbon uptake to carbon fixation
! cell_surf(auto_ind) = 3.886e-9_r8 ! surface areas of cells ( m2 )
! cell_carb_cont(auto_ind) = 1.68e-10_r8 ! cell carbon content ( mol C cell-1 )
! cell_radius(auto_ind) = ! cell radius ( um )
! cell_permea(auto_ind) = 1.1e-5_r8 ! cell wall permeability to CO2(aq) (m/s)
! cell_eps_fix(auto_ind) = 23.0_r8 ! fractionation effect of carbon fixation
else if (autotroph_settings(auto_ind)%Nfixer .and. &
autotroph_settings(auto_ind)%silicifier) then
log_message = "ciso: Currently Keller and Morel fractionation does not work for Diatoms-Diazotrophs"
call marbl_status_log%log_error(log_message, subname)
return
else
!----------------------------------------------------------------------------------------
! Small phytoplankton based on E. huxleyi ( Keller and morel, 1999; Popp et al., 1998 )
! Popp et al express cell carbon content in ( pg C cell-1 ), here we convert in (mol C cell-1)
! convert pgC to molC : ! Mc = 12 g / mol ! Mc = 12 e12 pg / mol
!----------------------------------------------------------------------------------------
cell_active_C_uptake(auto_ind) = 2.2_r8 ! ratio of active carbon uptake to carbon fixation
cell_surf(auto_ind) = 87.6e-12_r8 ! surface areas of cells ( m2 )
cell_carb_cont(auto_ind) = 69.2e-14_r8 ! cell carbon content ( mol C cell-1 )
cell_radius(auto_ind) = 2.6_r8 ! cell radius ( um )
cell_permea(auto_ind) = 1.8e-5_r8 ! cell wall permeability to CO2(aq) (m/s)
cell_eps_fix(auto_ind) = 25.3_r8 ! fractionation effect of carbon fixation
endif
end do
end select
end subroutine setup_cell_attributes
!***********************************************************************
subroutine fract_keller_morel( &
mui_to_co2star, &
cell_active_C_uptake, &
cell_surf, &
cell_carb_cont, &
cell_radius, &
cell_permea, &
cell_eps_fix, &
eps_p)
!---------------------------------------------------------------------------
! Calculate the carbon isotopic fractionation of phytoplankton photosynthesis : eps_p
!
! COMPUTATION : based on Keller and Morel, 1999
!
! eps_p = d13C(co2(aq)) - d13C(phyto)
!
! eps_p = eps_diff + (cell_active_C_uptake/(cell_active_C_uptake + 1/var)) * delta_d13C
! + ( (1 + (cell_active_C_uptake-1)*var )/(1+cell_active_C_uptake*var) )
! * ( cell_eps_fix - eps_diff )
!
! delta_d13C = d13C(CO2) - d13C(source)
! = 9 per mil
!
! var = mui_to_co2star * cell_carb_cont / ( cell_permea * cell_surf )
!
! mui_to_co2star = mu_i / [CO2*]
!
! --------
!
! Developed by X. Giraud, ETH Zürich, 21.07.2008
!---------------------------------------------------------------------------
use marbl_constants_mod, only : pi, c4, c3
real (r8), intent(in) :: &
mui_to_co2star, & ! mui_to_co2star = mu_i / [CO2*] (m3 / mol C / s)
cell_active_C_uptake, & ! ratio of active carbon uptake to carbon fixation
cell_eps_fix, & ! fractionation effect of carbon fixation
cell_surf, & ! surface areas of cells ( m2 )
cell_carb_cont, & ! cell carbon content ( mol C )
cell_radius, & ! cell radius ( um )
cell_permea ! cell wall permeability to CO2(aq) ( m /s )
real (r8), intent(out) :: eps_p ! = d13C(co2(aq)) - d13C(phyto)
!-----------------------------------------------------------------------
! local variables and parameters
!-----------------------------------------------------------------------
real (r8) :: &
var, theta, eps_up
real (r8) :: Vol, Qc, Surf, radius_m
real (r8) :: &
eps_diff = 0.7_r8, & ! fractionation by diffusion, O'Leary, 1984
delta_d13C = -9.0_r8 ! = d13C(CO2) - d13C(source), difference between the
! isotopic compositions of the external CO2 and the
! organic matter pools (Goericke et al. 1994).
! For active HCO3- uptake, the substrate for
! the carbon uptake mechanism has an isotopic
! composition which is around 9 permil higher
! than the external CO2, so D13CO2 -D13C_source
! is -9 permil ((Mook et al. 197)
!-----------------------------------------------------------------------
!---------------------------------------------------------------------
! cell surface in m^2
!---------------------------------------------------------------------
radius_m = cell_radius * 1e-6_r8 ! convert radius from um to m
if ( cell_surf > c0 ) then
Surf = cell_surf
else
Surf = c4 * pi * (radius_m ** 2)
endif
!---------------------------------------------------------------------
! cellular carbon content ( mol C )
! volume in um^3
!
! Qc based on Dieter A. Wolf-Gladrow, Ulf Riebeseel, Steffen Burkhardt and Jelle Buma (1999)
! Direct effects of CO2 concentration on growth and isotopic composition of marine plankton, based on
!
! Strathmann, 1967: Estimating the organic carbon content of phytoplankton from cell volume
! or plasma volume
!---------------------------------------------------------------------
if ( cell_carb_cont > c0 ) then
Qc = cell_carb_cont
else
Vol = c4 * pi * (cell_radius ** 3) / c3
Qc = 3.154e-14_r8 * (Vol ** (0.758_r8 ))
endif
!---------------------------------------------------------------------
! final expression of eps_p
!---------------------------------------------------------------------
if ( mui_to_co2star /= c0 ) then
var = mui_to_co2star * Qc / ( cell_permea * Surf )
theta = c1 + ( cell_active_C_uptake - c1 ) * var
theta = theta / ( c1 + cell_active_C_uptake * var )
eps_up = eps_diff + ( cell_active_C_uptake / &
( cell_active_C_uptake + c1 / var ) ) * delta_d13C
eps_p = eps_up + theta * ( cell_eps_fix - eps_diff )
else
eps_p = cell_eps_fix
end if
end subroutine fract_keller_morel
!***********************************************************************
subroutine set_surface_particulate_terms(POC_ciso, P_CaCO3_ciso)
!---------------------------------------------------------------------
! Set incoming fluxes (put into outgoing flux for first level usage).
! Set dissolution length, production fraction and mass terms.
!---------------------------------------------------------------------
type(column_sinking_particle_type), intent(inout) :: POC_ciso ! base units = nmol C_ciso
type(column_sinking_particle_type), intent(inout) :: P_CaCO3_ciso ! base units = nmol C_ciso