1 |
mmazloff |
1.4 |
C $Header: /u/gcmpack/MITgcm/pkg/bling/bling_dvm.F,v 1.3 2016/05/23 21:47:56 jmc Exp $ |
2 |
mmazloff |
1.1 |
C $Name: $ |
3 |
|
|
|
4 |
|
|
#include "BLING_OPTIONS.h" |
5 |
|
|
|
6 |
|
|
CBOP |
7 |
|
|
subroutine BLING_DVM( |
8 |
jmc |
1.2 |
I N_dvm,P_dvm,Fe_dvm, |
9 |
mmazloff |
1.1 |
I PTR_O2, mld, |
10 |
|
|
O N_remindvm, P_remindvm, Fe_remindvm, |
11 |
|
|
I bi, bj, imin, imax, jmin, jmax, |
12 |
|
|
I myIter, myTime, myThid ) |
13 |
|
|
|
14 |
|
|
C ================================================================= |
15 |
|
|
C | subroutine bling_prod |
16 |
|
|
C | o Nutrient uptake and partitioning between organic pools. |
17 |
jmc |
1.2 |
C | - Phytoplankton biomass-specific growth rate is calculated |
18 |
|
|
C | as a function of light, nutrient limitation, and |
19 |
|
|
C | temperature. |
20 |
mmazloff |
1.1 |
C | - Biomass growth xxx |
21 |
|
|
C | o Organic matter export, remineralization, and recycling. |
22 |
jmc |
1.2 |
C | - Sinking particulate flux and diel migration contribute to |
23 |
|
|
C | export. |
24 |
mmazloff |
1.1 |
C | - Denitrification xxx |
25 |
|
|
C ================================================================= |
26 |
|
|
|
27 |
|
|
implicit none |
28 |
jmc |
1.2 |
|
29 |
mmazloff |
1.1 |
C === Global variables === |
30 |
|
|
C P_sm :: Small phytoplankton biomass |
31 |
|
|
C P_lg :: Large phytoplankton biomass |
32 |
|
|
C P_diaz :: Diazotroph phytoplankton biomass |
33 |
|
|
C irr_mem :: Phyto irradiance memory |
34 |
|
|
|
35 |
|
|
#include "SIZE.h" |
36 |
|
|
#include "DYNVARS.h" |
37 |
|
|
#include "EEPARAMS.h" |
38 |
|
|
#include "PARAMS.h" |
39 |
|
|
#include "GRID.h" |
40 |
|
|
#include "BLING_VARS.h" |
41 |
|
|
#include "PTRACERS_SIZE.h" |
42 |
|
|
#include "PTRACERS_PARAMS.h" |
43 |
|
|
#ifdef ALLOW_AUTODIFF |
44 |
|
|
# include "tamc.h" |
45 |
|
|
#endif |
46 |
|
|
|
47 |
|
|
C === Routine arguments === |
48 |
|
|
C bi,bj :: tile indices |
49 |
|
|
C iMin,iMax :: computation domain: 1rst index range |
50 |
|
|
C jMin,jMax :: computation domain: 2nd index range |
51 |
|
|
C myTime :: current time |
52 |
|
|
C myIter :: current timestep |
53 |
|
|
C myThid :: thread Id. number |
54 |
|
|
INTEGER bi, bj, imin, imax, jmin, jmax |
55 |
|
|
_RL myTime |
56 |
|
|
INTEGER myIter |
57 |
|
|
INTEGER myThid |
58 |
|
|
C === Input === |
59 |
|
|
C PTR_NO3 :: nitrate concentration |
60 |
|
|
C PTR_PO4 :: phosphate concentration |
61 |
|
|
C PTR_FE :: iron concentration |
62 |
|
|
C PTR_DON :: dissolved organic nitrogen concentration |
63 |
|
|
C PTR_DOP :: dissolved organic phosphorus concentration |
64 |
|
|
C PTR_O2 :: oxygen concentration |
65 |
|
|
_RL PTR_NO3(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
66 |
|
|
_RL PTR_PO4(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
67 |
|
|
_RL PTR_FE (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
68 |
|
|
_RL PTR_DON(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
69 |
|
|
_RL PTR_DOP(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
70 |
|
|
_RL PTR_O2 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
71 |
|
|
C === Output === |
72 |
|
|
C G_xxx :: Tendency of xxx |
73 |
|
|
_RL G_NO3 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
74 |
|
|
_RL G_PO4 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
75 |
|
|
_RL G_FE (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
76 |
|
|
_RL G_O2 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
77 |
|
|
_RL G_DON (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
78 |
|
|
_RL G_DOP (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
79 |
|
|
_RL G_CACO3 (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
80 |
jmc |
1.2 |
|
81 |
mmazloff |
1.1 |
#ifdef ALLOW_BLING |
82 |
|
|
C === Local variables === |
83 |
|
|
C i,j,k :: loop indices |
84 |
|
|
C irr_eff :: effective irradiance |
85 |
jmc |
1.2 |
C NO3_lim :: nitrate limitation |
86 |
|
|
C PO4_lim :: phosphate limitation |
87 |
|
|
C Fe_lim :: iron limitation for phytoplankton |
88 |
mmazloff |
1.1 |
C Fe_lim_diaz :: iron limitation for diazotrophs |
89 |
|
|
C alpha_Fe :: initial slope of the P-I curve |
90 |
|
|
C irrk :: nut-limited efficiency of algal photosystems |
91 |
|
|
C Pc_m :: light-saturated max photosynthesis rate for phyt |
92 |
|
|
C Pc_m_diaz :: light-saturated max photosynthesis rate for diaz |
93 |
|
|
C Pc_tot :: carbon-specific photosynthesis rate |
94 |
|
|
C expkT :: temperature function |
95 |
|
|
C mu :: net carbon-specific growth rate for phyt |
96 |
|
|
C mu_diaz :: net carbon-specific growth rate for diaz |
97 |
|
|
C biomass_sm :: nitrogen concentration in small phyto biomass |
98 |
jmc |
1.2 |
C biomass_lg :: nitrogen concentration in large phyto biomass |
99 |
|
|
C N_uptake :: nitrogen uptake |
100 |
mmazloff |
1.1 |
C N_fix :: nitrogen fixation |
101 |
|
|
C P_uptake :: phosphorus uptake |
102 |
jmc |
1.2 |
C POC_flux :: carbon export flux 3d field |
103 |
mmazloff |
1.1 |
C PtoN :: variable ratio of phosphorus to nitrogen |
104 |
|
|
C FetoN :: variable ratio of iron to nitrogen |
105 |
jmc |
1.2 |
C N_spm :: particulate sinking of nitrogen |
106 |
mmazloff |
1.1 |
C P_spm :: particulate sinking of phosphorus |
107 |
jmc |
1.2 |
C Fe_spm :: particulate sinking of iron |
108 |
|
|
C N_dvm :: vertical transport of nitrogen by DVM |
109 |
|
|
C P_dvm :: vertical transport of phosphorus by DVM |
110 |
|
|
C Fe_dvm :: vertical transport of iron by DVM |
111 |
mmazloff |
1.1 |
C N_recycle :: recycling of newly-produced organic nitrogen |
112 |
|
|
C P_recycle :: recycling of newly-produced organic phosphorus |
113 |
|
|
C Fe_recycle :: recycling of newly-produced organic iron |
114 |
|
|
c xxx to be completed |
115 |
jmc |
1.2 |
INTEGER i,j,k |
116 |
mmazloff |
1.1 |
INTEGER tmp |
117 |
|
|
_RL irr_eff(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
118 |
jmc |
1.2 |
_RL NO3_lim |
119 |
mmazloff |
1.1 |
_RL PO4_lim |
120 |
|
|
_RL Fe_lim |
121 |
jmc |
1.2 |
_RL Fe_lim_diaz |
122 |
mmazloff |
1.1 |
_RL expkT |
123 |
|
|
_RL Pc_m |
124 |
|
|
_RL Pc_m_diaz |
125 |
jmc |
1.2 |
_RL irrk |
126 |
mmazloff |
1.1 |
_RL mu(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
127 |
|
|
_RL mu_diaz(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
128 |
|
|
_RL PtoN(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
129 |
|
|
_RL FetoN(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
130 |
|
|
_RL N_uptake(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
131 |
|
|
_RL N_fix(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
132 |
|
|
_RL P_uptake(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
133 |
|
|
_RL Fe_uptake(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
134 |
|
|
_RL frac_exp |
135 |
|
|
_RL N_spm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
136 |
|
|
_RL P_spm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
137 |
|
|
_RL Fe_spm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
138 |
|
|
_RL N_dvm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
139 |
|
|
_RL P_dvm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
140 |
|
|
_RL Fe_dvm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
141 |
|
|
_RL DON_prod(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
142 |
|
|
_RL DOP_prod(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
143 |
|
|
_RL N_recycle(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
144 |
|
|
_RL P_recycle(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
145 |
|
|
_RL Fe_recycle(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
146 |
|
|
_RL POC_flux(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
147 |
|
|
_RL NPP(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
148 |
|
|
_RL NCP(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
149 |
|
|
_RL CaCO3_prod(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
150 |
|
|
_RL PONflux_u |
151 |
|
|
_RL POPflux_u |
152 |
|
|
_RL PFEflux_u |
153 |
|
|
_RL CaCO3flux_u |
154 |
|
|
_RL PONflux_l |
155 |
|
|
_RL POPflux_l |
156 |
|
|
_RL PFEflux_l |
157 |
|
|
_RL CaCO3flux_l |
158 |
|
|
_RL depth_l |
159 |
|
|
_RL zremin |
160 |
|
|
_RL zremin_caco3 |
161 |
|
|
_RL wsink |
162 |
|
|
_RL POC_sed |
163 |
|
|
_RL Fe_sed(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
164 |
|
|
_RL NO3_sed(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
165 |
|
|
_RL PO4_sed(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
166 |
|
|
_RL O2_sed(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
167 |
|
|
_RL lig_stability |
168 |
|
|
_RL FreeFe |
169 |
|
|
_RL Fe_ads_inorg(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
170 |
|
|
_RL Fe_ads_org(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
171 |
|
|
_RL N_den_pelag(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
172 |
|
|
_RL N_den_benthic(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
173 |
|
|
_RL log_btm_flx |
174 |
|
|
_RL N_reminp(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
175 |
|
|
_RL P_reminp(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
176 |
|
|
_RL Fe_reminp(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
177 |
|
|
_RL CacO3_diss(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
178 |
|
|
_RL o2_upper |
179 |
jmc |
1.2 |
_RL o2_lower |
180 |
|
|
_RL dz_upper |
181 |
|
|
_RL dz_lower |
182 |
|
|
_RL temp_upper |
183 |
|
|
_RL temp_lower |
184 |
|
|
_RL z_dvm_regr |
185 |
|
|
_RL frac_migr |
186 |
|
|
_RL fdvm_migr |
187 |
|
|
_RL fdvm_stat |
188 |
|
|
_RL fdvmn_vint |
189 |
|
|
_RL fdvmp_vint |
190 |
|
|
_RL fdvmfe_vint |
191 |
mmazloff |
1.1 |
_RL z_dvm |
192 |
|
|
_RL N_remindvm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
193 |
|
|
_RL P_remindvm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
194 |
|
|
_RL Fe_remindvm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
195 |
|
|
_RL dvm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr) |
196 |
|
|
_RL mld(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
197 |
|
|
_RL Fe_burial(1-OLx:sNx+OLx,1-OLy:sNy+OLy) |
198 |
|
|
_RL x_erfcc,z_erfcc,t_erfcc,erfcc |
199 |
|
|
cxx order |
200 |
|
|
CEOP |
201 |
jmc |
1.2 |
|
202 |
mmazloff |
1.1 |
c --------------------------------------------------------------------- |
203 |
|
|
c Initialize output and diagnostics |
204 |
|
|
DO k=1,Nr |
205 |
|
|
DO j=jmin,jmax |
206 |
|
|
DO i=imin,imax |
207 |
|
|
G_NO3(i,j,k) = 0. _d 0 |
208 |
|
|
G_PO4(i,j,k) = 0. _d 0 |
209 |
|
|
G_FE(i,j,k) = 0. _d 0 |
210 |
|
|
G_O2(i,j,k) = 0. _d 0 |
211 |
|
|
G_DON(i,j,k) = 0. _d 0 |
212 |
|
|
G_DOP(i,j,k) = 0. _d 0 |
213 |
|
|
G_CaCO3(i,j,k) = 0. _d 0 |
214 |
|
|
N_uptake(i,j,k) = 0. _d 0 |
215 |
|
|
P_uptake(i,j,k) = 0. _d 0 |
216 |
|
|
Fe_uptake(i,j,k) = 0. _d 0 |
217 |
|
|
Fe_ads_org(i,j,k) = 0. _d 0 |
218 |
|
|
Fe_ads_inorg(i,j,k) = 0. _d 0 |
219 |
|
|
mu(i,j,k) = 0. _d 0 |
220 |
|
|
mu_diaz(i,j,k) = 0. _d 0 |
221 |
|
|
irr_eff(i,j,k) = 0. _d 0 |
222 |
|
|
PtoN(i,j,k) = 0. _d 0 |
223 |
|
|
FetoN(i,j,k) = 0. _d 0 |
224 |
|
|
N_fix(i,j,k) = 0. _d 0 |
225 |
|
|
N_spm(i,j,k) = 0. _d 0 |
226 |
|
|
P_spm(i,j,k) = 0. _d 0 |
227 |
|
|
Fe_spm(i,j,k) = 0. _d 0 |
228 |
|
|
dvm(i,j,k) = 0. _d 0 |
229 |
|
|
N_reminp(i,j,k) = 0. _d 0 |
230 |
|
|
P_reminp(i,j,k) = 0. _d 0 |
231 |
|
|
Fe_reminp(i,j,k) = 0. _d 0 |
232 |
|
|
N_recycle(i,j,k) = 0. _d 0 |
233 |
|
|
P_recycle(i,j,k) = 0. _d 0 |
234 |
|
|
Fe_recycle(i,j,k) = 0. _d 0 |
235 |
|
|
N_remindvm(i,j,k) = 0. _d 0 |
236 |
|
|
P_remindvm(i,j,k) = 0. _d 0 |
237 |
|
|
Fe_remindvm(i,j,k) = 0. _d 0 |
238 |
|
|
N_den_benthic(i,j,k)= 0. _d 0 |
239 |
|
|
N_den_pelag(i,j,k) = 0. _d 0 |
240 |
|
|
DON_prod(i,j,k) = 0. _d 0 |
241 |
|
|
DOP_prod(i,j,k) = 0. _d 0 |
242 |
|
|
CaCO3_prod(i,j,k) = 0. _d 0 |
243 |
|
|
CaCO3_diss(i,j,k) = 0. _d 0 |
244 |
|
|
POC_flux(i,j,k) = 0. _d 0 |
245 |
|
|
NPP(i,j,k) = 0. _d 0 |
246 |
|
|
NCP(i,j,k) = 0. _d 0 |
247 |
|
|
ENDDO |
248 |
|
|
ENDDO |
249 |
|
|
ENDDO |
250 |
|
|
DO j=jmin,jmax |
251 |
|
|
DO i=imin,imax |
252 |
|
|
Fe_burial(i,j) = 0. _d 0 |
253 |
|
|
NO3_sed(i,j) = 0. _d 0 |
254 |
|
|
PO4_sed(i,j) = 0. _d 0 |
255 |
|
|
O2_sed(i,j) = 0. _d 0 |
256 |
|
|
ENDDO |
257 |
|
|
ENDDO |
258 |
|
|
cxx order |
259 |
|
|
|
260 |
|
|
C --------------------------------------------------------------------- |
261 |
|
|
c DIEL VERTICAL MIGRATOR EXPORT |
262 |
|
|
c The effect of vertically-migrating animals on the export flux of organic |
263 |
|
|
c matter from the ocean surface is treated similarly to the scheme of |
264 |
jmc |
1.2 |
c Bianchi et al., Nature Geoscience 2013. |
265 |
mmazloff |
1.1 |
c This involves calculating the stationary depth of vertical migrators, using |
266 |
jmc |
1.2 |
c an empirical multivariate regression, and ensuring that this remains |
267 |
mmazloff |
1.1 |
c above the bottom as well as any suboxic waters. |
268 |
jmc |
1.2 |
c The total DVM export flux is partitioned between a swimming migratory |
269 |
mmazloff |
1.1 |
c component and the stationary component, and these are summed. |
270 |
|
|
|
271 |
|
|
C$TAF LOOP = parallel |
272 |
|
|
DO j=jmin,jmax |
273 |
|
|
C$TAF LOOP = parallel |
274 |
|
|
DO i=imin,imax |
275 |
|
|
|
276 |
|
|
C ! Initialize the working variables to zero |
277 |
jmc |
1.2 |
o2_upper = 0. |
278 |
|
|
o2_lower = 0. |
279 |
|
|
dz_upper = 0. |
280 |
|
|
dz_lower = 0. |
281 |
|
|
temp_upper = 0. |
282 |
|
|
temp_lower = 0. |
283 |
|
|
z_dvm_regr = 0. |
284 |
jmc |
1.3 |
z_dvm = 0. |
285 |
mmazloff |
1.1 |
frac_migr = 0. |
286 |
jmc |
1.2 |
fdvm_migr = 0. |
287 |
|
|
fdvm_stat = 0. |
288 |
|
|
fdvmn_vint = 0. |
289 |
|
|
fdvmp_vint = 0. |
290 |
mmazloff |
1.1 |
fdvmfe_vint = 0. |
291 |
|
|
|
292 |
|
|
DO k=1,Nr |
293 |
|
|
|
294 |
|
|
IF ( hFacC(i,j,k,bi,bj).gt.0. _d 0 ) THEN |
295 |
|
|
|
296 |
|
|
! Calculate the depth of migration based on linear regression. |
297 |
|
|
|
298 |
|
|
depth_l=-rF(k+1) |
299 |
jmc |
1.2 |
|
300 |
jmc |
1.3 |
! Average temperature and oxygen over upper 35 m, and 140-515m. |
301 |
|
|
! Also convert O2 to mmol m-3. |
302 |
mmazloff |
1.1 |
|
303 |
|
|
if ( abs(depth_l) .lt. 35.) then |
304 |
|
|
dz_upper = dz_upper + drf(k) |
305 |
|
|
temp_upper = temp_upper + theta(i,j,k,bi,bj)*drf(k) |
306 |
|
|
o2_upper = o2_upper + PTR_O2(i,j,k) * drf(k)*1.0 _d 3 |
307 |
|
|
endif |
308 |
jmc |
1.2 |
if ( (abs(depth_l) .gt. 140.0 _d 0) .and. |
309 |
mmazloff |
1.1 |
& (abs(depth_l) .lt. 515. _d 0)) then |
310 |
|
|
dz_lower = dz_lower + drf(k) |
311 |
|
|
temp_lower = temp_lower + theta(i,j,k,bi,bj)*drf(k) |
312 |
|
|
o2_lower = o2_lower + PTR_O2(i,j,k) * drf(k)*1.0 _d 3 |
313 |
|
|
endif |
314 |
|
|
|
315 |
|
|
ENDIF |
316 |
|
|
ENDDO |
317 |
jmc |
1.2 |
|
318 |
mmazloff |
1.1 |
o2_upper = o2_upper / (epsln + dz_upper) |
319 |
|
|
temp_upper = temp_upper / (epsln + dz_upper) |
320 |
|
|
o2_lower = o2_lower / (epsln + dz_lower) |
321 |
|
|
temp_lower = temp_lower / (epsln + dz_lower) |
322 |
|
|
|
323 |
jmc |
1.3 |
! Calculate the regression, using the constants given |
324 |
|
|
! in Bianchi et al. (2013). |
325 |
|
|
! The variable values are bounded to lie within reasonable |
326 |
|
|
! ranges: O2 gradient : [-10,300] mmol/m3 |
327 |
|
|
! Log10 Chl : [-1.8,0.85] log10(mg/m3) |
328 |
|
|
! mld : [0,500] m |
329 |
|
|
! T gradient : [-3,20] C |
330 |
mmazloff |
1.1 |
|
331 |
jmc |
1.2 |
C!! I am replacing hblt_depth(i,j) with mld... not sure what hblt_depth is |
332 |
mmazloff |
1.1 |
|
333 |
|
|
#ifdef BLING_ADJOINT_SAFE |
334 |
|
|
z_dvm = 300. _d 0 |
335 |
|
|
|
336 |
|
|
#else |
337 |
|
|
|
338 |
jmc |
1.2 |
z_dvm_regr = 398. _d 0 |
339 |
mmazloff |
1.1 |
& - 0.56 _d 0*min(300. _d 0,max(-10. _d 0,(o2_upper - o2_lower))) |
340 |
|
|
& - 115. _d 0*min(0.85 _d 0,max(-1.80 _d 0, |
341 |
|
|
& log10(max(chl(i,j,1,bi,bj),chl_min)))) |
342 |
|
|
& + 0.36 _d 0*min(500. _d 0,max(epsln,mld(i,j))) |
343 |
|
|
& - 2.40 _d 0*min(20. _d 0,max(-3. _d 0,(temp_upper-temp_lower))) |
344 |
|
|
|
345 |
|
|
c ! Limit the depth of migration in polar winter. |
346 |
|
|
c ! Use irr_mem since this is averaged over multiple days, dampening the diurnal cycle. |
347 |
|
|
c ! Tapers Z_DVM to the minimum when surface irradince is below a given threshold (here 10 W/m2). |
348 |
|
|
|
349 |
jmc |
1.2 |
if ( irr_mem(i,j,1,bi,bj) .lt. 10. ) then |
350 |
|
|
z_dvm_regr = 150. _d 0 + (z_dvm_regr - 150. _d 0) * |
351 |
mmazloff |
1.1 |
& irr_mem(i,j,1,bi,bj) / 10. _d 0 |
352 |
jmc |
1.2 |
endif |
353 |
mmazloff |
1.1 |
|
354 |
jmc |
1.2 |
C Check for suboxic water within the column. If found, set dvm |
355 |
|
|
C stationary depth to 2 layers above it. This is not meant to |
356 |
|
|
C represent a cessation of downward migration, but rather the |
357 |
|
|
C requirement for aerobic DVM respiration to occur above the suboxic |
358 |
mmazloff |
1.1 |
C water, where O2 is available. |
359 |
jmc |
1.2 |
|
360 |
mmazloff |
1.1 |
tmp = 0 |
361 |
|
|
DO k=1,Nr-2 |
362 |
|
|
|
363 |
|
|
IF ( (hFacC(i,j,k,bi,bj).gt.0. _d 0) .and. (tmp.eq.0)) THEN |
364 |
|
|
|
365 |
|
|
z_dvm = -rf(k+1) |
366 |
|
|
if (PTR_O2(i,j,k+2) .lt. (5. _d 0*oxic_min)) tmp = 1 |
367 |
jmc |
1.2 |
|
368 |
mmazloff |
1.1 |
ENDIF |
369 |
|
|
|
370 |
jmc |
1.2 |
enddo |
371 |
|
|
|
372 |
|
|
C The stationary depth is constrained between 150 and 700, above any |
373 |
mmazloff |
1.1 |
C anoxic waters found, and above the bottom. |
374 |
|
|
|
375 |
|
|
z_dvm = min(700. _d 0,max(150. _d 0,z_dvm_regr),z_dvm,-rf(k+1)) |
376 |
|
|
c!! bling%zbot(i,j,grid_kmt(i,j))) * grid_tmask(i,j,1) |
377 |
|
|
c!! what is grid_mkt?? |
378 |
|
|
|
379 |
|
|
#endif |
380 |
|
|
|
381 |
jmc |
1.3 |
! Calculate the fraction of migratory respiration that occurs |
382 |
|
|
! during upwards and downwards swimming. The remainder is |
383 |
|
|
! respired near the stationary depth. |
384 |
|
|
! Constants for swimming speed and resting time are hard-coded |
385 |
|
|
! after Bianchi et al, Nature Geoscience 2013. |
386 |
mmazloff |
1.1 |
|
387 |
|
|
frac_migr = max( 0.0 _d 0, min( 1.0 _d 0, (2.0 _d 0 * z_dvm) / |
388 |
|
|
& (epsln + 0.05 _d 0 * 0.5 _d 0 * 86400. _d 0))) |
389 |
|
|
|
390 |
jmc |
1.3 |
! Calculate the vertical profile shapes of DVM fluxes. |
391 |
|
|
! These are given as the downward organic flux due to migratory |
392 |
|
|
! DVM remineralization, defined at the bottom of each layer k. |
393 |
mmazloff |
1.1 |
|
394 |
|
|
tmp = 0 |
395 |
|
|
DO k=1,Nr |
396 |
|
|
|
397 |
|
|
IF ( (hFacC(i,j,k,bi,bj).gt.0. _d 0) .and. (tmp.eq.0)) THEN |
398 |
|
|
|
399 |
jmc |
1.3 |
! First, calculate the part due to active migration above |
400 |
|
|
! the stationary depth. |
401 |
jmc |
1.2 |
if (-rf(k+1) .lt. z_dvm) then |
402 |
|
|
fdvm_migr = frac_migr / (epsln + z_dvm - (-rf(2))) * |
403 |
mmazloff |
1.1 |
& (z_dvm - (-rf(k+1)) ) |
404 |
jmc |
1.2 |
else |
405 |
mmazloff |
1.1 |
fdvm_migr = 0.0 |
406 |
jmc |
1.2 |
endif |
407 |
mmazloff |
1.1 |
|
408 |
|
|
! Then, calculate the part at the stationary depth. |
409 |
jmc |
1.2 |
c fdvm_stat = (1. - frac_migr) / 2. * erfcc((-rf(k) - z_dvm) / |
410 |
mmazloff |
1.1 |
c & ( (epsln + 2. * sigma_dvm**2.)**0.5)) |
411 |
|
|
|
412 |
|
|
c ! Approximation of the complementary error function |
413 |
|
|
c ! From Numerical Recipes (F90, Ch. 6, p. 216) |
414 |
|
|
c ! Returns the complementary error function erfc(x) |
415 |
jmc |
1.2 |
c with fractional error everywhere less than 1.2e-7 |
416 |
|
|
x_erfcc = (-rf(k) - z_dvm) / |
417 |
mmazloff |
1.1 |
& ( (epsln + 2. _d 0 * sigma_dvm**2. _d 0)**0.5) |
418 |
|
|
|
419 |
|
|
z_erfcc = abs(x_erfcc) |
420 |
|
|
|
421 |
|
|
t_erfcc = 1. _d 0/(1. _d 0+0.5 _d 0*z_erfcc) |
422 |
|
|
|
423 |
|
|
erfcc = t_erfcc*exp(-z_erfcc*z_erfcc-1.26551223+t_erfcc* |
424 |
|
|
& (1.00002368+t_erfcc*(0.37409196+t_erfcc* |
425 |
|
|
& (.09678418+t_erfcc*(-.18628806+t_erfcc*(.27886807+ |
426 |
|
|
& t_erfcc*(-1.13520398+t_erfcc*(1.48851587+ |
427 |
|
|
& t_erfcc*(-0.82215223+t_erfcc*0.17087277))))))))) |
428 |
|
|
|
429 |
|
|
if (x_erfcc .lt. 0.0) then |
430 |
|
|
erfcc = 2.0 - erfcc |
431 |
|
|
endif |
432 |
jmc |
1.2 |
|
433 |
|
|
fdvm_stat = (1. _d 0 - frac_migr) / 2. _d 0 * erfcc |
434 |
|
|
|
435 |
|
|
C Add the shapes, resulting in the 3-d DVM flux operator. If the |
436 |
|
|
C current layer is the bottom layer, or the layer beneath the |
437 |
|
|
C underlying layer is suboxic, all fluxes at and below the current |
438 |
|
|
C layer remain at the initialized value of zero. This will cause all |
439 |
mmazloff |
1.1 |
C remaining DVM remineralization to occur in this layer. |
440 |
|
|
IF (k.LT.NR-1) THEN |
441 |
|
|
if (PTR_O2(i,j,k+2) .lt. (5. _d 0*oxic_min)) tmp = 1 |
442 |
|
|
ENDIF |
443 |
|
|
c!! if (k .eq. grid_kmt(i,j)) exit |
444 |
|
|
dvm(i,j,k) = fdvm_migr + fdvm_stat |
445 |
jmc |
1.2 |
|
446 |
mmazloff |
1.1 |
ENDIF |
447 |
|
|
|
448 |
jmc |
1.2 |
enddo |
449 |
|
|
|
450 |
|
|
c Sum up the total organic flux to be transported by DVM |
451 |
mmazloff |
1.1 |
|
452 |
jmc |
1.2 |
do k = 1, nr |
453 |
mmazloff |
1.1 |
fdvmn_vint = fdvmn_vint + N_dvm(i,j,k) * drf(k) |
454 |
|
|
fdvmp_vint = fdvmp_vint + P_dvm(i,j,k) * drf(k) |
455 |
|
|
fdvmfe_vint = fdvmfe_vint + Fe_dvm(i,j,k) * drf(k) |
456 |
jmc |
1.2 |
enddo |
457 |
mmazloff |
1.1 |
|
458 |
|
|
c Calculate the remineralization terms as the divergence of the flux |
459 |
|
|
|
460 |
|
|
N_remindvm(i,j,1) = fdvmn_vint * (1 - dvm(i,j,1)) / |
461 |
|
|
& (epsln + drf(1)) |
462 |
|
|
P_remindvm(i,j,1) = fdvmp_vint * (1 - dvm(i,j,1)) / |
463 |
|
|
& (epsln + drf(1)) |
464 |
jmc |
1.2 |
Fe_remindvm(i,j,1) = fdvmfe_vint * (1 - dvm(i,j,1)) / |
465 |
mmazloff |
1.1 |
& (epsln + drf(1)) |
466 |
|
|
|
467 |
|
|
do k = 2, nr |
468 |
jmc |
1.2 |
N_remindvm(i,j,k) = fdvmn_vint * |
469 |
mmazloff |
1.1 |
& (dvm(i,j,k-1) - dvm(i,j,k)) / (epsln + drf(k)) |
470 |
jmc |
1.2 |
P_remindvm(i,j,k) = fdvmp_vint * |
471 |
mmazloff |
1.1 |
& (dvm(i,j,k-1) - dvm(i,j,k)) / (epsln + drf(k)) |
472 |
|
|
Fe_remindvm(i,j,k) = fdvmfe_vint * |
473 |
|
|
& (dvm(i,j,k-1) - dvm(i,j,k)) / (epsln + drf(k)) |
474 |
jmc |
1.2 |
enddo |
475 |
mmazloff |
1.1 |
|
476 |
jmc |
1.2 |
enddo |
477 |
mmazloff |
1.1 |
enddo |
478 |
|
|
|
479 |
|
|
c --------------------------------------------------------------------- |
480 |
|
|
|
481 |
|
|
#ifdef ALLOW_DIAGNOSTICS |
482 |
|
|
IF ( useDiagnostics ) THEN |
483 |
|
|
|
484 |
|
|
CALL DIAGNOSTICS_FILL(Fe_dvm,'BLGFEDVM',0,Nr,2,bi,bj,myThid) |
485 |
|
|
CALL DIAGNOSTICS_FILL(N_dvm,'BLGNDVM ',0,Nr,2,bi,bj,myThid) |
486 |
|
|
CALL DIAGNOSTICS_FILL(P_dvm,'BLGPDVM ',0,Nr,2,bi,bj,myThid) |
487 |
|
|
CALL DIAGNOSTICS_FILL(dvm,'BLGDVM ',0,Nr,2,bi,bj,myThid) |
488 |
|
|
|
489 |
|
|
ENDIF |
490 |
|
|
#endif /* ALLOW_DIAGNOSTICS */ |
491 |
|
|
|
492 |
|
|
#endif /* ALLOW_BLING */ |
493 |
|
|
|
494 |
|
|
RETURN |
495 |
|
|
END |