/[MITgcm]/MITgcm/pkg/bling/bling_dvm.F
ViewVC logotype

Contents of /MITgcm/pkg/bling/bling_dvm.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.7 - (show annotations) (download)
Thu Mar 16 17:03:26 2017 UTC (7 years, 2 months ago) by mmazloff
Branch: MAIN
CVS Tags: checkpoint66g, checkpoint66f, checkpoint66o, checkpoint66n, checkpoint66m, checkpoint66l, checkpoint66k, checkpoint66j, checkpoint66i, checkpoint66h, HEAD
Changes since 1.6: +1 -13 lines
spring cleaning. minor changes

1 C $Header: /u/gcmpack/MITgcm/pkg/bling/bling_dvm.F,v 1.6 2016/11/16 16:41:50 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_dvm
16 C | o Diel Vertical Migration
17 C =================================================================
18
19 implicit none
20
21 C === Global variables ===
22
23 #include "SIZE.h"
24 #include "DYNVARS.h"
25 #include "EEPARAMS.h"
26 #include "PARAMS.h"
27 #include "GRID.h"
28 #include "BLING_VARS.h"
29 #include "PTRACERS_SIZE.h"
30 #include "PTRACERS_PARAMS.h"
31 #ifdef ALLOW_AUTODIFF
32 # include "tamc.h"
33 #endif
34
35 C === Routine arguments ===
36 C bi,bj :: tile indices
37 C iMin,iMax :: computation domain: 1rst index range
38 C jMin,jMax :: computation domain: 2nd index range
39 C myTime :: current time
40 C myIter :: current timestep
41 C myThid :: thread Id. number
42 INTEGER bi, bj, imin, imax, jmin, jmax
43 _RL myTime
44 INTEGER myIter
45 INTEGER myThid
46 C === Input ===
47 C N_dvm :: vertical transport of nitrogen by DVM
48 C P_dvm :: vertical transport of phosphorus by DVM
49 C Fe_dvm :: vertical transport of iron by DVM
50 C PTR_O2 :: nitrate concentration
51 C mld :: mixed layer depth
52 _RL N_dvm (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
53 _RL P_dvm (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
54 _RL Fe_dvm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
55 _RL PTR_O2(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
56 _RL mld (1-OLx:sNx+OLx,1-OLy:sNy+OLy)
57
58 C === Output ===
59 C N_remindvm :: nitrogen remineralization due to diel vertical migration
60 C P_remindvm :: phosphorus remineralization due to diel vertical migration
61 C Fe_remindvm :: iron remineralization due to diel vertical migration
62 _RL N_remindvm (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
63 _RL P_remindvm (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
64 _RL Fe_remindvm (1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
65
66 #ifdef ALLOW_BLING
67 C === Local variables ===
68 C i,j,k :: loop indices
69 INTEGER i,j,k
70 INTEGER tmp
71 _RL depth_l
72 _RL o2_upper
73 _RL o2_lower
74 _RL dz_upper
75 _RL dz_lower
76 _RL temp_upper
77 _RL temp_lower
78 _RL z_dvm_regr
79 _RL frac_migr
80 _RL fdvm_migr
81 _RL fdvm_stat
82 _RL fdvmn_vint
83 _RL fdvmp_vint
84 _RL fdvmfe_vint
85 _RL z_dvm
86 _RL dvm(1-OLx:sNx+OLx,1-OLy:sNy+OLy,Nr)
87 _RL Fe_burial(1-OLx:sNx+OLx,1-OLy:sNy+OLy)
88 _RL x_erfcc,z_erfcc,t_erfcc,erfcc
89 CEOP
90
91 c ---------------------------------------------------------------------
92 c Initialize output and diagnostics
93 DO k=1,Nr
94 DO j=jmin,jmax
95 DO i=imin,imax
96 N_remindvm(i,j,k) = 0. _d 0
97 P_remindvm(i,j,k) = 0. _d 0
98 Fe_remindvm(i,j,k) = 0. _d 0
99 dvm(i,j,k) = 0. _d 0
100 ENDDO
101 Fe_burial(i,j) = 0. _d 0
102 ENDDO
103 ENDDO
104
105 C ---------------------------------------------------------------------
106 c DIEL VERTICAL MIGRATOR EXPORT
107 c The effect of vertically-migrating animals on the export flux of organic
108 c matter from the ocean surface is treated similarly to the scheme of
109 c Bianchi et al., Nature Geoscience 2013.
110 c This involves calculating the stationary depth of vertical migrators, using
111 c an empirical multivariate regression, and ensuring that this remains
112 c above the bottom as well as any suboxic waters.
113 c The total DVM export flux is partitioned between a swimming migratory
114 c component and the stationary component, and these are summed.
115
116 C$TAF LOOP = parallel
117 DO j=jmin,jmax
118 C$TAF LOOP = parallel
119 DO i=imin,imax
120
121 c Initialize
122 o2_upper = 0.
123 o2_lower = 0.
124 dz_upper = 0.
125 dz_lower = 0.
126 temp_upper = 0.
127 temp_lower = 0.
128 z_dvm_regr = 0.
129 z_dvm = 0.
130 frac_migr = 0.
131 fdvm_migr = 0.
132 fdvm_stat = 0.
133 fdvmn_vint = 0.
134 fdvmp_vint = 0.
135 fdvmfe_vint = 0.
136
137 DO k=1,Nr
138
139 IF ( hFacC(i,j,k,bi,bj).gt.0. _d 0 ) THEN
140
141 c Calculate the depth of migration based on linear regression.
142
143 depth_l=-rF(k+1)
144
145 c Average temperature and oxygen over upper 35 m, and 140-515m.
146 c Also convert O2 to mmol m-3.
147
148 if ( abs(depth_l) .lt. 35.) then
149 dz_upper = dz_upper + drf(k)
150 temp_upper = temp_upper + theta(i,j,k,bi,bj)*drf(k)
151 o2_upper = o2_upper + PTR_O2(i,j,k) * drf(k)*1.0 _d 3
152 endif
153 if ( (abs(depth_l) .gt. 140.0 _d 0) .and.
154 & (abs(depth_l) .lt. 515. _d 0)) then
155 dz_lower = dz_lower + drf(k)
156 temp_lower = temp_lower + theta(i,j,k,bi,bj)*drf(k)
157 o2_lower = o2_lower + PTR_O2(i,j,k) * drf(k)*1.0 _d 3
158 endif
159
160 ENDIF
161 ENDDO
162
163 o2_upper = o2_upper / (epsln + dz_upper)
164 temp_upper = temp_upper / (epsln + dz_upper)
165 o2_lower = o2_lower / (epsln + dz_lower)
166 temp_lower = temp_lower / (epsln + dz_lower)
167
168 c Calculate the regression, using the constants given in Bianchi et al. (2013).
169 c The variable values are bounded to lie within reasonable ranges:
170 c O2 gradient : [-10,300] mmol/m3
171 c Log10 Chl : [-1.8,0.85] log10(mg/m3)
172 c mld : [0,500] m
173 c T gradient : [-3,20] C
174
175 z_dvm_regr = 398. _d 0
176 & - 0.56 _d 0*min(300. _d 0,max(-10. _d 0,(o2_upper - o2_lower)))
177 & - 115. _d 0*min(0.85 _d 0,max(-1.80 _d 0,
178 & log10(max(chl(i,j,1,bi,bj),chl_min))))
179 & + 0.36 _d 0*min(500. _d 0,max(epsln,mld(i,j)))
180 & - 2.40 _d 0*min(20. _d 0,max(-3. _d 0,(temp_upper-temp_lower)))
181
182 c Limit the depth of migration in polar winter.
183 c Use irr_mem since this is averaged over multiple days, dampening the
184 c diurnal cycle.
185 c Tapers Z_DVM to the minimum when surface irradince is below a given
186 c threshold (here 10 W/m2).
187
188 if ( irr_mem(i,j,1,bi,bj) .lt. 10. ) then
189 z_dvm_regr = 150. _d 0 + (z_dvm_regr - 150. _d 0) *
190 & irr_mem(i,j,1,bi,bj) / 10. _d 0
191 endif
192
193 c Check for suboxic water within the column. If found, set dvm
194 c stationary depth to 2 layers above it. This is not meant to
195 c represent a cessation of downward migration, but rather the
196 c requirement for aerobic DVM respiration to occur above the suboxic
197 c water, where O2 is available.
198
199 tmp = 0
200 DO k=1,Nr-2
201
202 IF ( (hFacC(i,j,k,bi,bj).gt.0. _d 0) .and. (tmp.eq.0)) THEN
203
204 z_dvm = -rf(k+1)
205 if (PTR_O2(i,j,k+2) .lt. (5. _d 0*oxic_min)) tmp = 1
206
207 ENDIF
208
209 enddo
210
211 c The stationary depth is constrained between 150 and 700, above any
212 c anoxic waters found, and above the bottom.
213
214 z_dvm = min(700. _d 0,max(150. _d 0,z_dvm_regr),z_dvm,-rf(k+1))
215
216 c Calculate the fraction of migratory respiration that occurs
217 c during upwards and downwards swimming. The remainder is
218 c respired near the stationary depth.
219 c Constants for swimming speed and resting time are hard-coded
220 c after Bianchi et al, Nature Geoscience 2013.
221
222 frac_migr = max( 0.0 _d 0, min( 1.0 _d 0, (2.0 _d 0 * z_dvm) /
223 & (epsln + 0.05 _d 0 * 0.5 _d 0 * 86400. _d 0)))
224
225 c Calculate the vertical profile shapes of DVM fluxes.
226 c These are given as the downward organic flux due to migratory
227 c DVM remineralization, defined at the bottom of each layer k.
228
229 tmp = 0
230 DO k=1,Nr
231
232 IF ( (hFacC(i,j,k,bi,bj).gt.0. _d 0) .and. (tmp.eq.0)) THEN
233
234 ! First, calculate the part due to active migration above
235 ! the stationary depth.
236 if (-rf(k+1) .lt. z_dvm) then
237 fdvm_migr = frac_migr / (epsln + z_dvm - (-rf(2))) *
238 & (z_dvm - (-rf(k+1)) )
239 else
240 fdvm_migr = 0.0
241 endif
242
243 c Then, calculate the part at the stationary depth.
244
245 c Approximation of the complementary error function
246 c From Numerical Recipes (F90, Ch. 6, p. 216)
247 c Returns the complementary error function erfc(x)
248 c with fractional error everywhere less than 1.2e-7
249 x_erfcc = (-rf(k) - z_dvm) /
250 & ( (epsln + 2. _d 0 * sigma_dvm**2. _d 0)**0.5)
251
252 z_erfcc = abs(x_erfcc)
253
254 t_erfcc = 1. _d 0/(1. _d 0+0.5 _d 0*z_erfcc)
255
256 erfcc = t_erfcc*exp(-z_erfcc*z_erfcc-1.26551223+t_erfcc*
257 & (1.00002368+t_erfcc*(0.37409196+t_erfcc*
258 & (.09678418+t_erfcc*(-.18628806+t_erfcc*(.27886807+
259 & t_erfcc*(-1.13520398+t_erfcc*(1.48851587+
260 & t_erfcc*(-0.82215223+t_erfcc*0.17087277)))))))))
261
262 if (x_erfcc .lt. 0.0) then
263 erfcc = 2.0 - erfcc
264 endif
265
266 fdvm_stat = (1. _d 0 - frac_migr) / 2. _d 0 * erfcc
267
268 c Add the shapes, resulting in the 3-d DVM flux operator. If the
269 c current layer is the bottom layer, or the layer beneath the
270 c underlying layer is suboxic, all fluxes at and below the current
271 c layer remain at the initialized value of zero. This will cause all
272 c remaining DVM remineralization to occur in this layer.
273 IF (k.LT.NR-1) THEN
274 if (PTR_O2(i,j,k+2) .lt. (5. _d 0*oxic_min)) tmp = 1
275 ENDIF
276 c!! if (k .eq. grid_kmt(i,j)) exit
277 dvm(i,j,k) = fdvm_migr + fdvm_stat
278
279 ENDIF
280
281 enddo
282
283 c Sum up the total organic flux to be transported by DVM
284
285 do k = 1, nr
286 fdvmn_vint = fdvmn_vint + N_dvm(i,j,k) * drf(k)
287 fdvmp_vint = fdvmp_vint + P_dvm(i,j,k) * drf(k)
288 fdvmfe_vint = fdvmfe_vint + Fe_dvm(i,j,k) * drf(k)
289 enddo
290
291 c Calculate the remineralization terms as the divergence of the flux
292
293 N_remindvm(i,j,1) = fdvmn_vint * (1 - dvm(i,j,1)) /
294 & (epsln + drf(1))
295 P_remindvm(i,j,1) = fdvmp_vint * (1 - dvm(i,j,1)) /
296 & (epsln + drf(1))
297 Fe_remindvm(i,j,1) = fdvmfe_vint * (1 - dvm(i,j,1)) /
298 & (epsln + drf(1))
299
300 do k = 2, nr
301 N_remindvm(i,j,k) = fdvmn_vint *
302 & (dvm(i,j,k-1) - dvm(i,j,k)) / (epsln + drf(k))
303 P_remindvm(i,j,k) = fdvmp_vint *
304 & (dvm(i,j,k-1) - dvm(i,j,k)) / (epsln + drf(k))
305 Fe_remindvm(i,j,k) = fdvmfe_vint *
306 & (dvm(i,j,k-1) - dvm(i,j,k)) / (epsln + drf(k))
307 enddo
308
309 enddo
310 enddo
311
312 #endif /* ALLOW_BLING */
313
314 RETURN
315 END

  ViewVC Help
Powered by ViewVC 1.1.22