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