-
Notifications
You must be signed in to change notification settings - Fork 0
/
recom_sms.F
executable file
·2631 lines (2504 loc) · 88 KB
/
recom_sms.F
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
C$Header: /csys/software/cvsroot/REcoM/pkg/recom/recom_sms.F,v 1.61 2007/10/03 06:15:18 mlosch Exp $
C$Name: $
#include "PTRACERS_OPTIONS.h"
#include "GCHEM_OPTIONS.h"
#include "RECOM_OPTIONS.h"
CBOP
C !ROUTINE: RECOM_SMS
C !INTERFACE:
subroutine REcoM_sms(
& iMin, iMax, jMin, jMax, bi, bj,
& Nz, kLowC, bgc_num, ndiags3d, ndiags2d, delta_t,
& maskC, hFacC, recip_hFacC, drF, recip_drF, recip_drC,
& temp, dicFlux_loc,
& surf_light,cobeta,
& state,
& benthos,
& sms, diags3d, diags2d
#ifdef RECOM_WAVEBANDS
& ,PARlocal
& ,PARwlocal
& ,C_phot_nl
& ,C_phot_nl_dia
& ,Ek_nl
& ,Ek_nl_dia
#endif
& )
C !DESCRIPTION:
C=====================================================================
C Carbon and Nitrogen Regulated Ecosystem Model (CN-REcoM)
C version 2 (2 Phytoplankton size classes)
C NOTE: Phytoplankton growth equations
C mainly derived from Geider et al., (1998) L&O
C
C Model code started by Markus Schartau (AWI Bremerhaven, October 2002)
C Maintained by Martin Losch, Martin.Losch@awi.de
C=====================================================================
C
C The field state(Nz,bgc_num) contains all biological tracers,
C with the following units and typical values in the open ocean
C (coastal regions may differ):
C
C state(:,1): DIN (dissolved nitrate NO3)
C units: mmol/m^3, typical range: 0 < DIN < 54
C
C state(:,2): DIC (dissolved inorganic carbon, sometimes TCO2)
C units: mmol/m^3, typical range 2000 < DIC < 2400
C
C state(:,3): alk (alkalinity)
C units: mmol/m^3, typical range 2000 < alk < 2500
C
C state(:,4): phyN (small phytoplankton nitrogen)
C units: mmol/m^3,
C typical range: 0 < phyN < 5 (open ocean)
C < 10 (some coastal regions)
C
C state(:,5): phyC (small phytoplankton carbon)
C units: mmol/m^3:
C typical range: 0 < phyC < 35 (open ocean)
C 0 < phyC < 70 (coastal regions)
C
C state(:,6): CHLA (small phytoplankton chlorophyll a)
C units: mg/m^3, typical range: 0 < CHLA < 4
C
C state(:,7): detN (nitrate-based detritus)
C units: mmol/m^3
C typical range: 0 < detN < (hmm, can't really tell... approx.) 10
C
C state(:,8): detC (carbon-based detritus)
C units: mmol/m^3
C typical range: 0 < detC < (even worse to tell... approx.) 70
C
C state(:,9): hetN (heterotrophic zooplankton on N-basis)
C units: mmol/m^3
C typical range: 0 < hetN < (depends on size classes)
C
C state(:,10): hetC (heterotrophic zooplankton on C-basis)
C units: mmol/m^3
C typical range: 0 < hetC < (depends on size classes)
C
C state(:,11): DON (extracellular dissolved organic nitrogen)
C units: mmol/m^3
C typical range: model considers only "fresh" DON, which can range from
C 5-20, background values are 10-20
C
C state(:,12): EOC (extracellular dissolved organic carbon)
C units: mmol/m^3
C typical range: 0 < EOC < 100
C (this refers to "fresh" labile DOC); however, in the model part
C of of EOC is particulate material (POC); we call the remaining EOC,
C which does not contain any polysaccharides (PCHO), "residual DOC".
C EOC=DOC+TEPC=PCHO+residual EOC+TEPC
C
C state(:,13): diaN (diatom nitrogen)
C units: mmol/m^3,
C typical range: 0 < phyN < 5 (open ocean)
C < 10 (some coastal regions)
C
C state(:,14): diaC (diatom carbon)
C units: mmol/m^3:
C typical range: 0 < phyC < 35 (open ocean)
C 0 < phyC < 70 (coastal regions)
C
C state(:,15): DiaCHL (diatom chlorophyll a)
C units: mg/m^3, typical range: 0 < CHLA < 4
C
C state(:,16): diaSi (diatom silica)
C units: mmol/m^3
C typical range: 0-??
C
C state(:,17): detSi (silica-based detritus)
C units: mmol/m^3
C typical range: 0-??
C
C state(:,18): Si (disolved silicate)
C units: mmol/m^3
C typical range: 0 < Si < 170
C
C state(:,19): Fe (disolved iron)
C units: mumol/m^3
C typical range: 0.1 < Fe < 0.6
C
C state(:,20): phyCalc (phytoplankton calcium carbonate)
C units: mumol/m^3
C typical range: 0.0 < phyCalc < 1.0
C
C state(:,21): detCalc (detritus calcium carbonate)
C units: mumol/m^3
C typical range: 0.1 < detCalc < 1.0
C
#ifdef RECOM_ZOO2
CMS state(:,24): zoo2N (heterotrophic zooplankton on N-basis)
C units: mmol/m^3
C typical range: 0 < zoo2N < (depends on size classes)
CMS state(:,25): zoo2C (heterotrophic zooplankton on C-basis)
C units: mmol/m^3
C typical range: 0 < zoo2C < (depends on size classes)
#endif
C !USES:
C --- definitions
#ifdef ALLOW_MODULES
use mod_REcoM_para_def
implicit none
#else /* not ALLOW_MODULES */
implicit none
#include "RECOM_PARAMS.h"
#ifdef ALLOW_MITGCM
#include "SIZE.h"
#endif
#ifdef RECOM_WAVEBANDS
#include "SPECTRAL_SIZE.h"
#include "SPECTRAL.h"
#include "WAVEBANDS_PARAMS.h"
#endif
C----------
CCV The next two lines are necessary to use the pointers to
CCV biogeochemical tracers, defined in RECOM.h
#include "EEPARAMS.h"
#include "RECOM.h"
C----------
#endif /* not ALLOW_MODULES */
C
C !INPUT/OUTPUT PARAMETERS:
C
C loop boundaries
integer iMin, iMax, jMin, jMax, bi, bj
C external timestep in seconds (physics !)
_RL delta_t
C number of layers
integer Nz
C number of tracers passed
integer bgc_num
C number of diagnostics
integer ndiags3d, ndiags2d, ndiags3d_used
C
integer kLowC (1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
C thickness of grid boxes
_RS drF(Nz), recip_drF(Nz), recip_drC(Nz)
_RS maskC (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz,nSx,nSy)
_RS hFacC (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz,nSx,nSy)
_RS recip_hFacC(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz,nSx,nSy)
C temperature profile (in degC)
_RL temp(1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz,nSx,nSy)
C atmospheric CO2 flux
_RL dicFlux_loc(1-Olx:sNx+Olx,1-Oly:sNy+Oly,nSx,nSy)
C radiation at surface
_RL surf_light(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
C-----
_RL cobeta(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
C-----
_RL state (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz,nSx,nSy,bgc_num)
_RL benthos(1-Olx:sNx+Olx,1-Oly:sNy+Oly,4)
C-----
C auxillary variables for the flux into the benthic layer
_RL wFluxDet(4)
_RL wFluxPhy(4)
_RL wFluxDia(4)
_RL decayBenthos(4)
C [mmol m^{-3}]
_RL sms (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz,nSx,nSy,bgc_num)
C-----
C generic array for diagnostics
_RL diags3d (1-Olx:sNx+Olx,1-Oly:sNy+Oly,Nz,ndiags3d)
_RL diags2d (1-Olx:sNx+Olx,1-Oly:sNy+Oly,ndiags2d)
C-----
#ifdef RECOM_WAVEBANDS
C spectral PAR is used to compute specral CPhot
C all the other variables are used in forcing and just pass to sms to be exported
_RL PARlocal(1-Olx:sNx+Olx,1-Oly:sNy+Oly, Nz) !average PAR at midpoint of gridcell
_RL PARwlocal(1-Olx:sNx+Olx,1-Oly:sNy+Oly, tlam, Nz) !PAR at midpoint of previous(in) and local(out) gridcell
_RL C_phot_nl(tlam)
_RL C_phot_nl_dia(tlam)
_RL Ek_nl(tlam)
_RL Ek_nl_dia(tlam)
#endif
C !LOCAL VARIABLES:
integer kLoc
C sinking velocities are positive downwards
_RL vSinkPhy(1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nz)
_RL vSinkDia(1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nz)
_RL vSinkDet(1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nz)
C tendency due to sinking
_RL sink (1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nz) ! [mmol m^{-3}]
#ifdef RECOM_EXPORT_DIAGNOSTICS
C diagnostic dummy variable for calculating vertical export fluxes
_RL export (1-Olx:sNx+Olx,1-Oly:sNy+Oly)
#endif
#ifdef RECOM_MAREMIP
_RL export3d(1-Olx:sNx+Olx,1-Oly:sNy+Oly,1:Nz)
#endif
_RL arrFunc, rTloc, rTref
C-----
C small phytoplankton process variables
C-----
C carbon specific rate of photosynthesis [d^{-1}]
_RL C_phot
#ifdef RECOM_WAVEBANDS
_RL Ek
_RL alpha_I
#endif
C maximum photosynthetic rate
_RL pMax
C carbon specific nitrogen utilization rate [mmol N mmol C^{-1} d^{-1}]
_RL N_assim
C CHL a synthesis regulation term [mg CHL mmol N^{-1}]
_RL chlSynth
C molar cell quota [mmol N mmol C^{-1}]
_RL quota, recipQuota
C cellular CHL a : carbon/nitrogen ratio [mg CHL mmol C/N^{-1}]
_RL CHL2C, CHL2N
#ifdef RECOM_PHOTODAMAGE
_RL CHL2C_plast
#endif
C phytoplankton respiration rate [d^{-1}]
_RL phyRespRate
C coefficient for damage to the photosynthetic apparatus
_RL KOchl
#ifdef RECOM_MARSHALL
C quantum yield of photosynthesis [mmolC J-1]
_RL QY
C antena-based non-photochemical quenching [rel]
_RL NPQ
C damage rate [d-1]
_RL damage
C repair rate [d-1]
_RL repair
C initial slope of the PE curve [m2 mmolC-1 mgChl-1 J-1]
_RL ALPHAmar
C#ifdef RECOM_WAVEBANDS
C _RL ALPHAmar_nl(tlam)
C#endif
#endif
C
C-----
C diatom process variables
C-----
C carbon specific rate of photosynthesis [d^{-1}]
_RL C_phot_dia
#ifdef RECOM_WAVEBANDS
_RL Ek_dia
_RL alpha_I_dia
#endif
C maximum photosynthetic rate
_RL pMax_dia
C carbon specific nitrogen utilization rate [mmol N mmol C^{-1} d^{-1}]
_RL N_assim_dia
C CHL a synthesis regulation term [mg CHL mmol N^{-1}]
_RL chlSynth_dia
C molar cell quota [mmol N mmol C^{-1}]
_RL quota_dia, recipQuota_dia
C cellular CHL a : carbon/nitrogen ratio [mg CHL mmol C/N^{-1}]
_RL CHL2C_dia, CHL2N_dia
#ifdef RECOM_PHOTODAMAGE
_RL CHL2C_plast_dia
#endif
C phytoplankton respiration rate [d^{-1}]
_RL phyRespRate_dia
C carbon specific silicate uptake rate [mmol Si mmol C^{-1} d^{-1}]
_RL Si_assim
C molar cell quotas [mmol Si mmol N^{-1} and mmol Si mmol C^{-1}]
_RL qSiN, qSiC
C coefficient for damage to the photosynthetic apparatus
_RL KOchl_dia
#ifdef RECOM_MARSHALL
C quantum yield of photosynthesis [mmolC J-1]
_RL QY_dia
C antena-based non-photochemical quenching [rel]
_RL NPQ_dia
C damage rate [d-1]
_RL damage_dia
C repair rate [d-1]
_RL repair_dia
C initial slope of the PE curve [m2 mmolC-1 mgChl-1 J-1]
_RL ALPHAmar_dia
C#ifdef RECOM_WAVEBANDS
C _RL ALPHAmar_dia_nl(tlam)
C#endif
#endif
C
C-----
C zooplankton process variables
C-----
C inverse molar cell quota of zooplankton [(mmol N mmol C^{-1})^{-1}]
_RL recipQZoo
C grazing flux [mmol N m^{-3} d^{-1}]
_RL fdiaN
#ifdef RECOM_GRAZING_VARIABLE_PREFERENCE
_RL varpzdia, diaNsq
#endif
_RL grazingFlux
_RL grazingFlux_dia, grazingFlux_phy, food, foodsq
#ifdef RECOM_ZOO2
_RL grazingFlux_zoo2
#endif
C zooplankton respiration
_RL hetRespFlux
C zooplankton mortality (quadratic loss)
_RL hetLossFlux
#ifdef RECOM_ZOO2
C-----
C Second zooplankton process variables
C-----
C inverse molar cell quota of zooplankton [(mmol N mmol C^{-1})^{-1}]
_RL recipQZoo2, recipQuota_z2
C grazing flux [mmol N m^{-3} d^{-1}]
_RL fdiaN2
#ifdef RECOM_GRAZING_VARIABLE_PREFERENCE
_RL varpzdia2, diaNsq2
#endif
#endif
#ifdef RECOM_ZOO2
_RL grazingFluxz2
_RL grazingzoo2Flux_dia, grazingzoo2Flux_phy, food2, foodsq2
C zooplankton respiration
_RL zoo2RespFlux
C zooplankton mortality (quadratic loss)
_RL zoo2LossFlux
#endif
C aggregationRate (of nitrogen) [d^{-1}]
_RL aggregationRate
C tendency due to carbon air-sea flux
_RL dicSurfTend
C averaged photosynthetic available radiation (PAR)
_RL PARave
C cumulative vertical integral of kappa (light attenuation)
_RL kdzUpper(1-Olx:sNx+Olx,1-Oly:sNy+Oly)
_RL kdzLower, upper_light, lower_light
_RL kappa, kappaStar
_RL din, dic, alk, phyN, phyC, CHL
_RL detN, detC, hetN, hetC, DON, eoC
_RL diaN, diaC, diaCHL, diaSi, detSi, Si, Fe, freeFe
_RL calc_diss
#ifdef RECOM_CDOM
_RL cdomC, cdom_photo_rate
#endif
#ifdef RECOM_ZOO2
_RL zoo2N, zoo2C
#endif
#ifdef RECOM_MARSHALL
_RL D1, diaD1
_RL PPC, diaPPC
C#ifdef RECOM_WAVEBANDS
C _RL PPC_nl(tlam), diaPPC_nl(tlam)
C#endif
#endif
C temperature dependent remineralisation rate of Si
_RL reminSiT
C some loop parameters
integer i, j, k, step, ii
_RL dt
C light limitation variables
_RL qLimitFac, qLimitFacTmp, feLimitFac
#ifdef RECOM_MARSHALL
_RL qLimitFac_dia, qLimitFacTmp_dia, feLimitFac_dia
#endif
C downregulation of metabolic processes
_RL limitFacN
_RL limitFacN_dia, limitFacSi
C Temperature dependent maximum of C-specific nitrogen uptake
C [mmol N (mmol C)^{-1} d^{-1}]
_RL V_cm
C Limiter function
_RL recom_limiter
external recom_limiter
C external function to compute free iron
_RL iron_chemistry
external iron_chemistry
_RL phyCalc, detCalc
_RL calcification, calc_loss_agg, calc_loss_gra
#ifdef RECOM_LENA
_RL P1, P2, Q0, Pslope, flena
#endif
#ifdef RECOM_WAVEBANDS
integer ilam
integer nl
_RL cu_area
_RL PARwdn(tlam) !light at bottom of local gridcell
_RL attenwl(tlam) !attenuation (m-1)
_RL sumaphy_nl(tlam) !total phyto absorption at each wavelength
#endif
integer bgc_num_expected
parameter ( bgc_num_expected = 21
#ifdef RECOM_O2
& + 1
#endif
#ifdef RECOM_ZOO2
& + 2
#endif
#ifdef RECOM_CDOM
& + 1
#endif
#ifdef RECOM_MARSHALL
& + 2
#endif
& )
CEOP
C check some parameters
if ( bgc_num .ne. bgc_num_expected ) then
write(*,*) 'S/R REcoM_sms: bgc_num does not have the '//
& 'expected value'
write(*,*) 'bgc_num = ', bgc_num, ', bgc_num_expected = ',
& bgc_num_expected
stop 'UNEXPECTED IN RECOM_SMS'
endif
C
C---- INITIAL set to zero ------------------------------------------
C
DO k = 1,4
wFluxDet(k) = 0.
wFluxPhy(k) = 0.
wFluxDia(k) = 0.
decayBenthos(k) = 0.
ENDDO
C calc_diss_rate = 20.0/(3500*86400)
CMS Probably not necessary, but to keep the structure:
quota = c0
CHL2C = c0
CHL2N = c0
recipQuota = c0
qSiN = c0
qSiC = c0
dicSurfTend = c0
N_assim = c0
C_phot = c0
chlSynth = c0
aggregationRate= c0
phyRespRate = c0
hetRespFlux = c0
grazingFlux = c0
#ifdef RECOM_ZOO2
zoo2RespFlux = c0
grazingFluxz2 = c0
#endif
PARave = c0
do ii=1,bgc_num
do k=1,Nz
do j=jMin,jMax
do i=iMin,iMax
sms(i,j,k,bi,bj,ii) = c0
enddo
enddo
enddo
enddo
c---
do k=1,Nz
do j=jMin,jMax
do i=iMin,iMax
sink(i,j,k) = c0
enddo
enddo
enddo
c---
do ii=1,ndiags3d
do k=1,Nz
do j=jMin,jMax
do i=iMin,iMax
diags3d(i,j,k,ii) = c0
enddo
enddo
enddo
enddo
c---
do ii=1,ndiags2d
do j=jMin,jMax
do i=iMin,iMax
diags2d(i,j,ii) = c0
enddo
enddo
enddo
c---
do j=jMin,jMax
do i=iMin,iMax
vSinkPhy(i,j,1) = c0
vSinkDia(i,j,1) = c0
vSinkDet(i,j,1) = c0
enddo
enddo
do k=2,Nz
do j=jMin,jMax
do i=iMin,iMax
vSinkPhy(i,j,k) = Vphy*maskC(i,j,k,bi,bj)
vSinkDia(i,j,k) = Vdia*maskC(i,j,k,bi,bj)
CTW increased sinking speed
#ifdef ALLOW_SINK_INCREASE
vSinkDet(i,j,k) = Vdetfast(k)*maskC(i,j,k,bi,bj)
#else
vSinkDet(i,j,k) = Vdet*maskC(i,j,k,bi,bj)
#endif /* ALLOW_SINK_INCREASE */
enddo
enddo
enddo
#ifdef RECOM_WAVEBANDS
do ilam=1,tlam
C_phot_nl(ilam)=0. _d 0
C_phot_nl_dia(ilam)=0. _d 0
Ek_nl(ilam)=0. _d 0
Ek_nl_dia(ilam)=0. _d 0
enddo
#endif
C
C BIOLOGICAL LOOP ------------------------------------------------------------
C
rTref = c1/recom_Tref
dt = one_day/delta_t
dt = c1/(dt*real(bio_step)) ! time increment unit 'day'
do step = 1, bio_step
C before starting the vertical loop, set light attenuation,
C integral(kappa*dz), to zero at the surface
do j = jMin, jMax
do i = iMin, iMax
kdzUpper(i,j) = 0.
end do
end do
C set tiny tendencies to zero to avoid numerical problems
C (is this really necessary?)
do ii = 1, bgc_num
do k = 1, Nz
do j = jMin, jMax
do i = iMin, iMax
if (abs(sms(i,j,k,bi,bj,ii)).le.tiny) sms(i,j,k,bi,bj,ii)=0.
end do
end do
end do
end do
do k = 1, Nz
do j = jMin, jMax
do i = iMin, iMax
C friendly attachment of state variables to names, avoid divisions
C by zero by setting a lower limit
C state: 1. DIN, 2. DIC, 3. ALK, 4. PhyN, 5. PhyC, 6. Chl,
C 7. DetN, 8.DetC, 9. HetN, 10. HetC, 11. DON, 12. EOC
C 13. DiaN, 14. DiaC, 15. DiaCHL, 16.DiaSi, 17. DetSi, 18. Si 19. Fe
din = max(tiny,state(i,j,k,bi,bj,idin) +
& sms(i,j,k,bi,bj, idin))
dic = max(tiny,state(i,j,k,bi,bj,idic) +
& sms(i,j,k,bi,bj, idic))
alk = max(tiny,state(i,j,k,bi,bj,ialk) +
& sms(i,j,k,bi,bj, ialk))
phyN = max(tiny,state(i,j,k,bi,bj,iphyn) +
& sms(i,j,k,bi,bj, iphyn))
phyC = max(tiny,state(i,j,k,bi,bj,iphyc) +
& sms(i,j,k,bi,bj, iphyc))
CHL = max(tiny,state(i,j,k,bi,bj,ipchl) +
& sms(i,j,k,bi,bj, ipchl))
detN = max(tiny,state(i,j,k,bi,bj,idetn) +
& sms(i,j,k,bi,bj, idetn))
detC = max(tiny,state(i,j,k,bi,bj,idetc) +
& sms(i,j,k,bi,bj, idetc))
hetN = max(tiny,state(i,j,k,bi,bj,ihetn) +
& sms(i,j,k,bi,bj, ihetn))
hetC = max(tiny,state(i,j,k,bi,bj,ihetc) +
& sms(i,j,k,bi,bj, ihetc))
DON = max(tiny,state(i,j,k,bi,bj,idon) +
& sms(i,j,k,bi,bj, idon))
eoC = max(tiny,state(i,j,k,bi,bj,idoc) +
& sms(i,j,k,bi,bj, idoc))
diaN = max(tiny,state(i,j,k,bi,bj,idian) +
& sms(i,j,k,bi,bj,idian))
diaC = max(tiny,state(i,j,k,bi,bj,idiac) +
& sms(i,j,k,bi,bj,idiac))
diaChl = max(tiny,state(i,j,k,bi,bj,idchl) +
& sms(i,j,k,bi,bj,idchl))
diaSi= max(tiny,state(i,j,k,bi,bj,idiasi) +
& sms(i,j,k,bi,bj,idiasi))
detSi= max(tiny,state(i,j,k,bi,bj,idetsi) +
& sms(i,j,k,bi,bj,idetsi))
Si = max(tiny,state(i,j,k,bi,bj,isi)+sms(i,j,k,bi,bj,isi))
Fe = max(tiny,state(i,j,k,bi,bj,ife)+sms(i,j,k,bi,bj,ife))
freeFe = 0.
phyCalc = max(tiny,state(i,j,k,bi,bj,iphycalc) +
& sms(i,j,k,bi,bj,iphycalc))
detCalc = max(tiny,state(i,j,k,bi,bj,idetcalc) +
& sms(i,j,k,bi,bj,idetcalc))
CJH calc_diss_rate was initilaized for constant sinking speed of 20m/d
CJH with increasing sinking speed we have to correct for that to obtain
CJH the same profile
calc_diss = calc_diss_rate * vSinkDet(i,j,k) /20
#ifdef RECOM_CDOM
cdomC = max(tiny,state(i,j,k,bi,bj,icdom) +
& sms(i,j,k,bi,bj,icdom))
#endif
#ifdef RECOM_ZOO2
zoo2N = max(tiny,state(i,j,k,bi,bj,izoo2n) +
& sms(i,j,k,bi,bj, izoo2n))
zoo2C = max(tiny,state(i,j,k,bi,bj,izoo2c) +
& sms(i,j,k,bi,bj, izoo2c))
#endif /* RECOM_ZOO2 */
#ifdef RECOM_MARSHALL
D1 = max(tiny,state(i,j,k,bi,bj,id1) +
& sms(i,j,k,bi,bj,id1))
diaD1 = max(tiny,state(i,j,k,bi,bj,id1d) +
& sms(i,j,k,bi,bj,id1d))
#endif
C-----
C small phytoplankton, diatom and zooplankton cell quotas
C-----
quota = phyN / phyC
#ifdef RECOM_WAVEBANDS
quota = max(NCmin, quota)
quota = min(quota, NCmax)
#endif
recipQuota = 1./quota
CHL2C = CHL / phyC
CHL2N = CHL / phyN
#ifdef RECOM_PHOTODAMAGE
if (quota .gt. NCmin) then
CHL2C_plast = CHL2C * (quota/(quota - NCmin))
else
CHL2C_plast = CHL2C
end if
#endif
quota_dia = diaN / diaC
#ifdef RECOM_WAVEBANDS
quota_dia = max(NCmin_d, quota_dia)
quota_dia = min(quota_dia, NCmax_d)
#endif
recipQuota_dia = 1./quota_dia
CHL2C_dia = diaCHL / diaC
CHL2N_dia = diaCHL / diaN
#ifdef RECOM_PHOTODAMAGE
if (quota_dia .gt. NCmin_d) then
CHL2C_plast_dia = CHL2C_dia*(quota_dia/(quota_dia-NCmin_d))
else
CHL2C_plast_dia = CHL2C_dia
end if
#endif
qSiC = diaSi / diaC
qSiC = max(SiCmin, qSiC)
qSiC = min(qSiC, SiCmax)
qSiN = diaSi / diaN
recipQZoo = hetC / hetN
#ifdef RECOM_ZOO2
recipQZoo2 = zoo2C / zoo2N
recipQuota_z2 = recipQZoo2
#endif
C-----
C temperature dependence of rates
C-----
rTloc = c1/(temp(i,j,k,bi,bj) + recom_celsius2K)
arrFunc = exp( -Ae * ( rTloc - rTref ) )*maskC(i,j,k,bi,bj)
C silicate remineralisation has a different temperature dependence
reminSiT = reminSi*arrFunc
if ( useReminSiT )
& reminSiT = min(1.32E16*exp(-11200.*rTloc),reminSi)
& *maskC(i,j,k,bi,bj)
C----- upper boundary condition for carbon: exchange with atmopshere
C----- mmol C m^{-2} d^{-1} --> mmol C m^{-3} d^{-1}
dicSurfTend = dicFlux_loc(i,j,bi,bj)
& *recip_drF(k)*recip_hfacc(i,j,k,bi,bj) * max( 2-k, 0 )
CML if (k.eq.1) then
CML dicSurfTend=dicFlux_loc*recip_drF(k)*recip_hfacc(i,j,k,bi,bj)
CML else
CML dicSurfTend=c0
CML end if
C-------------------------------------------------------------------------
C---- photosynthesis section, light parameters and rates
C----
C small phytoplankton pMax
C----
qLimitFac = recom_limiter(NMinSlope,NCmin,quota)
if (FeLimit) then
feLimitFac = Fe/(k_Fe + Fe)
qLimitFac = min(qLimitFac,feLimitFac)
end if
#ifdef RECOM_LENA
P2 = P_cm
P1 = 0.0
Q0 = 0.04
Pslope = 100.0
flena = (P1 + P2)*0.5 + Pslope*(quota - Q0)
flena = max(P1,flena)
flena = min(P2,flena)
pMax = flena*qLimitFac*arrFunc
#else
pMax = P_cm*qLimitFac*arrFunc
#endif
C----
C diatom pMax
C----
qLimitFac_dia = recom_limiter(NMinSlope,NCmin_d,quota_dia)
qLimitFacTmp_dia = recom_limiter(SiMinSlope,SiCmin,qSiC)
qLimitFac_dia = min(qLimitFac_dia,qLimitFacTmp_dia)
if (FeLimit) then
feLimitFac_dia = Fe/(k_Fe_d + Fe)
qLimitFac_dia = min(qLimitFac_dia,feLimitFac_dia)
end if
#ifdef RECOM_LENA
P2 = P_cm_d
P1 = 0.0
Q0 = 0.04
Pslope = 100.0
flena = (P1 + P2)*0.5 + Pslope*(quota_dia - Q0)
flena = max(P1,flena)
flena = min(P2,flena)
pMax_dia = flena * qLimitFac_dia * arrFunc
#else
pMax_dia = P_cm_d * qLimitFac_dia * arrFunc
#endif
C----
C light
C----
#ifdef RECOM_WAVEBANDS
PARave = 0.
do nl = 1,tlam
PARave = PARave + (PARwlocal(i,j,nl,k) / WtouEins(nl))
end do
PARave = max(tiny, PARave)
C PARave = max(tiny,(PARlocal(i,j,k)/4.596))
#else
C attenuation coefficient
kappa = k_w + a_chl*(CHL + diaCHL)
kappaStar = kappa/cobeta(i,j)
C vertical light profile, averaged over grid cell
C - first, integrate kappaStar*deltaZ for a another layer
kdzLower = kdzUpper(i,j) + kappaStar*drF(k)*hFacC(i,j,k,bi,bj)
C - compute available light at upper and lower face of current layer
upper_light = surf_light(i,j)*exp( -kdzUpper(i,j) )
lower_light = surf_light(i,j)*exp( -kdzLower )
C - vertical average over box
PARave = max(tiny,(upper_light-lower_light)/kappaStar
& * recip_drF(k)*recip_hfacc(i,j,k,bi,bj))
C - store kdzLower for next layer
kdzUpper(i,j) = kdzLower
#endif
#ifdef RECOM_WAVEBANDS
c Get value for alphachl_nl * PARwlocal (value will depend on matchup between spectra of
c alphachl_nl and PARwlocal) and integrate the product over wavebands to get alpha_I
alpha_I = 0.d0
alpha_I_dia = 0.d0
do nl = 1,tlam
alpha_I = alpha_I
& + alphachl_nl(nl) * PARwlocal(i,j,nl,k)
alpha_I_dia = alpha_I_dia
& + alphachl_nl_dia(nl) * PARwlocal(i,j,nl,k)
end do
c Units of alpha_nl and alpha_mean m2 mgChla-1 mmolC uE-1 (x86400 s d-1)
alpha_I = alpha_I*86400
alpha_I_dia = alpha_I_dia*86400
c#else
c Without spectral alpha, the single value of alpha (parameter) is used.
c It should be comparable to alpha_mean precomputed in recom_init_vari.F
C Check the units of PARave and PARlocal.
c alpha=alpha_mean
c alpha_d=alpha_mean_dia
#endif
#ifdef RECOM_MARSHALL
C----
C Photoinhibition (Marshall et al 2000)-------------------------
C----
C Phyto
C qLimitFac = recom_limiter(NMinSlope,NCmin,quota)
C if (FeLimit) then
C feLimitFac = Fe/(k_Fe + Fe)
C qLimitFac = min(qLimitFac,feLimitFac)
C end if
if (D1 .gt. c1) then
QY = QYmax
else
QY = min(((QYmax/Drel) * D1) , QYmax)
endif
NPQ = 1-(QY/QYmax)
#ifdef RECOM_WAVEBANDS
astar = 0.d0
do nl = 1,tlam
C ALPHAmar_nl(nl) = aphy_chl_ps(nl) * (c1-NPQ) * QY * 86400
astar = astar + (wb_width(nl) * aphy_chl_ps(nl))
end do
astar = astar / wb_totalWidth
#endif
ALPHAmar = astar * (c1-NPQ) * QY * 86400
if (PARave .lt. tiny) then
damage = tiny
else
damage = (PARave * k_deg) * D1 * (c1-NPQ) * 86400
endif
if (D1 .gt. c1) then
repair = tiny
else
repair = ( ( r_max * (c1 - D1) )
& / ( k_rep + (c1 - D1) )
& ) * qLimitFac * arrFunc
endif
C
C Diatoms
C qLimitFac_dia = recom_limiter(NMinSlope,NCmin_d,quota_dia)
C qLimitFacTmp_dia = recom_limiter(SiMinSlope,SiCmin,qSiC)
C qLimitFac_dia = min(qLimitFac_dia,qLimitFacTmp_dia)
C if (FeLimit) then
C feLimitFac_dia = Fe/(k_Fe_d + Fe)
C qLimitFac_dia = min(qLimitFac_dia,feLimitFac_dia)
C end if
if (diaD1 .gt. c1) then
QY_dia = QYmax_d
else
QY_dia = min( (QYmax_d/Drel)* diaD1 , QYmax_d)
endif
NPQ_dia = 1-(QY_dia/QYmax_d)
#ifdef RECOM_WAVEBANDS
astar_d = 0.d0
do nl = 1,tlam
C ALPHAmar_dia_nl(nl) = aphy_chl_ps_dia(nl) * (c1-NPQ_dia)
C & * QY_dia * 86400
astar_d = astar_d + (wb_width(nl) * aphy_chl_ps_dia(nl))
end do
astar_d = astar_d / wb_totalWidth
#endif
ALPHAmar_dia = astar_d * (c1-NPQ_dia) * QY_dia * 86400
if (PARave .lt. tiny) then
damage_dia = tiny
else
damage_dia = (PARave * k_deg_d) * diaD1
& * (c1-NPQ_dia) * 86400
endif
if (diaD1 .gt. c1) then
repair_dia = tiny
else
repair_dia = ( ( r_max_d * (c1 - diaD1) )
& / ( k_rep_d + (c1 - diaD1) )
& ) * qLimitFac_dia * arrFunc
endif
CEA Proxy for non-photosynthetic pigments
PPC = 1-(ALPHAmar / (astar * QYmax * 86400) )
diaPPC = 1-(ALPHAmar_dia / (astar_d * QYmax_d * 86400) )
#endif /* RECOM_MARSHALL */
C----
C small phytoplankton photosynthesis rate
C----
#ifdef RECOM_WAVEBANDS
if (pMax.ge.tiny .AND. alpha_I.ge.tiny
& .AND. CHL2C.ge.tiny) then
C_phot = pMax * ( c1 - exp(-alpha_I * CHL2C / pMax))
if (PARlocal(i,j,k) .ge. tiny
& .AND. alpha_mean .ge. tiny) then
Ek = c1 - exp((-alpha_mean * 86400 * CHL2C
& * PARlocal(i,j,k)) / pMax)
else
Ek = c0
end if
do nl=1, tlam
if (PARwlocal(i,j,nl,k) .ge. tiny
& .AND. alphachl_nl(nl) .ge. tiny) then
C_phot_nl(nl)= pMax * ( c1 - exp((-alphachl_nl(nl)
& * 86400 * CHL2C * PARwlocal(i,j,nl,k)) / pMax))
Ek_nl(nl)= c1 - exp((-alphachl_nl(nl) * 86400
& * CHL2C * PARwlocal(i,j,nl,k)) / pMax)
else
C_phot_nl(nl) = c0
Ek_nl(nl) = c0
end if
enddo
else
C_phot = c0
Ek = c0
end if
alpha = alpha_mean
#else
if (pMax .lt. tiny) then
C_phot = c0
else
C_phot = pMax * ( c1 - exp(-alpha * CHL2C * PARave / pMax) )
end if
if (C_phot .lt. tiny) C_phot = c0
#endif
C
C----
C diatom photosynthesis rate
C----
#ifdef RECOM_WAVEBANDS
if (pMax_dia.ge.tiny .AND. alpha_I_dia.ge.tiny
& .AND. CHL2C_dia.ge.tiny) then
C_phot_dia = pMax_dia *
& ( c1 - exp(-alpha_I_dia * CHL2C_dia / pMax_dia) )
if (PARlocal(i,j,k).ge.tiny
& .AND. alpha_mean_dia.ge.tiny) then
Ek_dia = c1 - exp((-alpha_mean_dia * 86400
& * CHL2C_dia * PARlocal(i,j,k)) / pMax_dia)
else
Ek_dia = c0
end if
do nl=1, tlam
if (PARwlocal(i,j,nl,k) .ge. tiny
& .AND. alphachl_nl_dia(nl) .ge. tiny) then
C_phot_nl_dia(nl)= pMax_dia*(c1-exp((-alphachl_nl_dia(nl)
& * 86400 * CHL2C_dia * PARwlocal(i,j,nl,k))/pMax_dia))
Ek_nl_dia(nl) = c1 - exp((-alphachl_nl_dia(nl) * 86400
& * CHL2C_dia * PARwlocal(i,j,nl,k)) / pMax_dia)
else
C_phot_nl_dia(nl) = c0
Ek_nl_dia(nl) = c0
end if
enddo
else
C_phot_dia = c0
Ek_dia = c0
end if
alpha_d = alpha_mean_dia
#else
if (pMax_dia .lt. tiny) then
C_phot_dia = c0
else
C_phot_dia = pMax_dia *
& ( c1 - exp(-alpha_d * CHL2C_dia * PARave / pMax_dia) )
end if
if (C_phot_dia .lt. tiny) C_phot_dia = c0
#endif
C
C----
C chlorophyll degradation
C----
KOchl = deg_CHL
KOchl_dia = deg_CHL_d