1 | MODULE p4zpoc |
---|
2 | !!====================================================================== |
---|
3 | !! *** MODULE p4zpoc *** |
---|
4 | !! TOP : PISCES Compute remineralization of organic particles |
---|
5 | !!========================================================================= |
---|
6 | !! History : 1.0 ! 2004 (O. Aumont) Original code |
---|
7 | !! 2.0 ! 2007-12 (C. Ethe, G. Madec) F90 |
---|
8 | !! 3.4 ! 2011-06 (O. Aumont, C. Ethe) Quota model for iron |
---|
9 | !! 3.6 ! 2016-03 (O. Aumont) Quota model and diverse |
---|
10 | !!---------------------------------------------------------------------- |
---|
11 | #if defined key_pisces |
---|
12 | !!---------------------------------------------------------------------- |
---|
13 | !! 'key_top' and TOP models |
---|
14 | !! 'key_pisces' PISCES bio-model |
---|
15 | !!---------------------------------------------------------------------- |
---|
16 | !! p4z_poc : Compute remineralization/dissolution of organic compounds |
---|
17 | !! p4z_poc_init : Initialisation of parameters for remineralisation |
---|
18 | !!---------------------------------------------------------------------- |
---|
19 | USE oce_trc ! shared variables between ocean and passive tracers |
---|
20 | USE trc ! passive tracers common variables |
---|
21 | USE sms_pisces ! PISCES Source Minus Sink variables |
---|
22 | USE p4zsink ! Sedimentation of organic particles |
---|
23 | USE prtctl_trc ! print control for debugging |
---|
24 | USE iom ! I/O manager |
---|
25 | |
---|
26 | |
---|
27 | IMPLICIT NONE |
---|
28 | PRIVATE |
---|
29 | |
---|
30 | PUBLIC p4z_poc ! called in p4zbio.F90 |
---|
31 | PUBLIC p4z_poc_init ! called in trcsms_pisces.F90 |
---|
32 | PUBLIC alngam |
---|
33 | PUBLIC gamain |
---|
34 | |
---|
35 | !! * Shared module variables |
---|
36 | REAL(wp), PUBLIC :: xremip !: remineralisation rate of POC |
---|
37 | INTEGER , PUBLIC :: jcpoc !: number of lability classes |
---|
38 | REAL(wp), PUBLIC :: rshape !: shape factor of the gamma distribution |
---|
39 | |
---|
40 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:) :: alphan, reminp |
---|
41 | REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: alphap |
---|
42 | |
---|
43 | |
---|
44 | !!* Substitution |
---|
45 | # include "top_substitute.h90" |
---|
46 | !!---------------------------------------------------------------------- |
---|
47 | !! NEMO/TOP 3.3 , NEMO Consortium (2010) |
---|
48 | !! $Id: p4zrem.F90 3160 2011-11-20 14:27:18Z cetlod $ |
---|
49 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
50 | !!---------------------------------------------------------------------- |
---|
51 | CONTAINS |
---|
52 | |
---|
53 | SUBROUTINE p4z_poc( kt, jnt ) |
---|
54 | !!--------------------------------------------------------------------- |
---|
55 | !! *** ROUTINE p4z_poc *** |
---|
56 | !! |
---|
57 | !! ** Purpose : Compute remineralization of organic particles |
---|
58 | !! |
---|
59 | !! ** Method : - ??? |
---|
60 | !!--------------------------------------------------------------------- |
---|
61 | ! |
---|
62 | INTEGER, INTENT(in) :: kt, jnt ! ocean time step |
---|
63 | ! |
---|
64 | INTEGER :: ji, jj, jk, jn |
---|
65 | REAL(wp) :: zremip, zremig, zdep, zorem, zorem2, zofer |
---|
66 | REAL(wp) :: zsizek, zsizek1, alphat, remint, solgoc, zpoc |
---|
67 | #if ! defined key_kriest |
---|
68 | REAL(wp) :: zofer2, zofer3 |
---|
69 | #endif |
---|
70 | REAL(wp) :: zstep, zrfact2 |
---|
71 | CHARACTER (len=25) :: charout |
---|
72 | REAL(wp), POINTER, DIMENSION(:,: ) :: totprod, totthick, totcons |
---|
73 | REAL(wp), POINTER, DIMENSION(:,:,:) :: zremipoc, zremigoc, zorem3 |
---|
74 | REAL(wp), POINTER, DIMENSION(:,:,:,:) :: alphag |
---|
75 | REAL(wp), POINTER, DIMENSION(:,:,:) :: zfolimi |
---|
76 | |
---|
77 | !!--------------------------------------------------------------------- |
---|
78 | ! |
---|
79 | IF( nn_timing == 1 ) CALL timing_start('p4z_poc') |
---|
80 | ! |
---|
81 | ! Allocate temporary workspace |
---|
82 | CALL wrk_alloc( jpi, jpj, totprod, totthick, totcons ) |
---|
83 | CALL wrk_alloc( jpi, jpj, jpk, zremipoc, zremigoc, zorem3 ) |
---|
84 | CALL wrk_alloc( jpi, jpj, jpk, zfolimi ) |
---|
85 | ALLOCATE( alphag(jpi,jpj,jpk,jcpoc) ) |
---|
86 | |
---|
87 | ! Initialization of local variables |
---|
88 | ! --------------------------------- |
---|
89 | |
---|
90 | ! Here we compute the GOC -> POC rate due to the shrinking |
---|
91 | ! of the fecal pellets/aggregates as a result of bacterial |
---|
92 | ! solubilization |
---|
93 | ! This is based on a fractal dimension of 2.56 and a spectral |
---|
94 | ! slope of -3.6 (identical to what is used in p4zsink to compute |
---|
95 | ! aggregation |
---|
96 | solgoc = 0.04/ 2.56 * 1./ ( 1.-50**(-0.04) ) |
---|
97 | |
---|
98 | ! Initialisation of temprary arrys |
---|
99 | zremipoc(:,:,:) = xremip |
---|
100 | zremigoc(:,:,:) = xremip |
---|
101 | zorem3(:,:,:) = 0. |
---|
102 | orem (:,:,:) = 0. |
---|
103 | zfolimi(:,:,:) = 0. |
---|
104 | |
---|
105 | DO jn = 1, jcpoc |
---|
106 | alphag(:,:,:,jn) = alphan(jn) |
---|
107 | alphap(:,:,:,jn) = alphan(jn) |
---|
108 | END DO |
---|
109 | |
---|
110 | #if ! defined key_kriest |
---|
111 | ! ----------------------------------------------------------------------- |
---|
112 | ! Lability parameterization. This is the big particles part (GOC) |
---|
113 | ! This lability parameterization can be activated only with the standard |
---|
114 | ! particle scheme. Does not work with Kriest parameterization. |
---|
115 | ! ----------------------------------------------------------------------- |
---|
116 | DO jk = 2, jpkm1 |
---|
117 | DO jj = 1, jpj |
---|
118 | DO ji = 1, jpi |
---|
119 | IF (tmask(ji,jj,jk) == 1.) THEN |
---|
120 | zdep = hmld(ji,jj) |
---|
121 | ! |
---|
122 | ! In the case of GOC, lability is constant in the mixed layer |
---|
123 | ! It is computed only below the mixed layer depth |
---|
124 | ! ------------------------------------------------------------ |
---|
125 | ! |
---|
126 | IF( fsdept(ji,jj,jk) > zdep ) THEN |
---|
127 | alphat = 0. |
---|
128 | remint = 0. |
---|
129 | ! |
---|
130 | zsizek1 = fse3t(ji,jj,jk-1) / 2. / (wsbio4(ji,jj,jk-1) + rtrn) * tgfunc(ji,jj,jk-1) |
---|
131 | zsizek = fse3t(ji,jj,jk) / 2. / (wsbio4(ji,jj,jk) + rtrn) * tgfunc(ji,jj,jk) |
---|
132 | ! |
---|
133 | IF ( fsdept(ji,jj,jk-1) <= zdep ) THEN |
---|
134 | ! |
---|
135 | ! The first level just below the mixed layer needs a |
---|
136 | ! specific treatment because lability is supposed constant |
---|
137 | ! everywhere within the mixed layer. This means that |
---|
138 | ! change in lability in the bottom part of the previous cell |
---|
139 | ! should not be computed |
---|
140 | ! ---------------------------------------------------------- |
---|
141 | ! |
---|
142 | ! POC concentration is computed using the lagrangian |
---|
143 | ! framework. It is only used for the lability param |
---|
144 | zpoc = trb(ji,jj,jk-1,jpgoc) + consgoc(ji,jj,jk) * rday / rfact2 & |
---|
145 | & * fse3t(ji,jj,jk) / 2. / (wsbio4(ji,jj,jk) + rtrn) |
---|
146 | zpoc = max(0., zpoc) |
---|
147 | ! |
---|
148 | DO jn = 1, jcpoc |
---|
149 | ! |
---|
150 | ! Lagrangian based algorithm. The fraction of each |
---|
151 | ! lability class is computed starting from the previous |
---|
152 | ! level |
---|
153 | ! ----------------------------------------------------- |
---|
154 | ! |
---|
155 | ! the concentration of each lability class is calculated |
---|
156 | ! as the sum of the different sources and sinks |
---|
157 | ! Please note that production of new GOC experiences |
---|
158 | ! degradation |
---|
159 | alphag(ji,jj,jk,jn) = alphag(ji,jj,jk-1,jn) * exp( -reminp(jn) * zsizek ) * zpoc & |
---|
160 | & + prodgoc(ji,jj,jk) * alphan(jn) / tgfunc(ji,jj,jk) / reminp(jn) & |
---|
161 | & * ( 1. - exp( -reminp(jn) * zsizek ) ) * rday / rfact2 |
---|
162 | alphat = alphat + alphag(ji,jj,jk,jn) |
---|
163 | remint = remint + alphag(ji,jj,jk,jn) * reminp(jn) |
---|
164 | END DO |
---|
165 | ELSE |
---|
166 | ! |
---|
167 | ! standard algorithm in the rest of the water column |
---|
168 | ! See the comments in the previous block. |
---|
169 | ! --------------------------------------------------- |
---|
170 | ! |
---|
171 | zpoc = trb(ji,jj,jk-1,jpgoc) + consgoc(ji,jj,jk-1) * rday / rfact2 & |
---|
172 | & * fse3t(ji,jj,jk-1) / 2. / (wsbio4(ji,jj,jk-1) + rtrn) + consgoc(ji,jj,jk) & |
---|
173 | & * rday / rfact2 * fse3t(ji,jj,jk) / 2. / (wsbio4(ji,jj,jk) + rtrn) |
---|
174 | zpoc = max(0., zpoc) |
---|
175 | ! |
---|
176 | DO jn = 1, jcpoc |
---|
177 | alphag(ji,jj,jk,jn) = alphag(ji,jj,jk-1,jn) * exp( -reminp(jn) * ( zsizek & |
---|
178 | & + zsizek1 ) ) * zpoc + ( prodgoc(ji,jj,jk-1) / tgfunc(ji,jj,jk-1) * ( 1. & |
---|
179 | & - exp( -reminp(jn) * zsizek1 ) ) * exp( -reminp(jn) * zsizek ) + prodgoc(ji,jj,jk) & |
---|
180 | & / tgfunc(ji,jj,jk) * ( 1. - exp( -reminp(jn) * zsizek ) ) ) * rday / rfact2 / reminp(jn) |
---|
181 | alphat = alphat + alphag(ji,jj,jk,jn) |
---|
182 | remint = remint + alphag(ji,jj,jk,jn) * reminp(jn) |
---|
183 | END DO |
---|
184 | ENDIF |
---|
185 | DO jn = 1, jcpoc |
---|
186 | ! The contribution of each lability class at the current |
---|
187 | ! level is computed |
---|
188 | alphag(ji,jj,jk,jn) = alphag(ji,jj,jk,jn) / ( alphat + rtrn) |
---|
189 | END DO |
---|
190 | ! Computation of the mean remineralisation rate |
---|
191 | zremigoc(ji,jj,jk) = MIN(xremip, MAX(0., remint / ( alphat + rtrn) )) |
---|
192 | ENDIF |
---|
193 | ENDIF |
---|
194 | END DO |
---|
195 | END DO |
---|
196 | END DO |
---|
197 | |
---|
198 | DO jk = 1, jpkm1 |
---|
199 | DO jj = 1, jpj |
---|
200 | DO ji = 1, jpi |
---|
201 | zstep = xstep |
---|
202 | # if defined key_degrad |
---|
203 | zstep = zstep * facvol(ji,jj,jk) |
---|
204 | # endif |
---|
205 | ! POC disaggregation by turbulence and bacterial activity. |
---|
206 | ! -------------------------------------------------------- |
---|
207 | zremig = zremigoc(ji,jj,jk) * zstep * tgfunc(ji,jj,jk) |
---|
208 | zorem2 = zremig * trb(ji,jj,jk,jpgoc) |
---|
209 | orem(ji,jj,jk) = zorem2 |
---|
210 | zorem3(ji,jj,jk) = zremig * solgoc * trb(ji,jj,jk,jpgoc) |
---|
211 | zofer2 = zremig * trb(ji,jj,jk,jpbfe) |
---|
212 | zofer3 = zremig * solgoc * trb(ji,jj,jk,jpbfe) |
---|
213 | |
---|
214 | ! Update the appropriate tracers trends |
---|
215 | ! ------------------------------------- |
---|
216 | |
---|
217 | tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + zorem3(ji,jj,jk) |
---|
218 | tra(ji,jj,jk,jpgoc) = tra(ji,jj,jk,jpgoc) - zorem2 - zorem3(ji,jj,jk) |
---|
219 | tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + zofer3 |
---|
220 | tra(ji,jj,jk,jpbfe) = tra(ji,jj,jk,jpbfe) - zofer2 - zofer3 |
---|
221 | tra(ji,jj,jk,jpdoc) = tra(ji,jj,jk,jpdoc) + zorem2 |
---|
222 | tra(ji,jj,jk,jpfer) = tra(ji,jj,jk,jpfer) + zofer2 |
---|
223 | zfolimi(ji,jj,jk) = zofer2 |
---|
224 | END DO |
---|
225 | END DO |
---|
226 | END DO |
---|
227 | |
---|
228 | IF(ln_ctl) THEN ! print mean trends (used for debugging) |
---|
229 | WRITE(charout, FMT="('poc1')") |
---|
230 | CALL prt_ctl_trc_info(charout) |
---|
231 | CALL prt_ctl_trc(tab4d=tra, mask=tmask, clinfo=ctrcnm) |
---|
232 | ENDIF |
---|
233 | |
---|
234 | ! ------------------------------------------------------------------ |
---|
235 | ! Lability parameterization for the small OM particles. This param |
---|
236 | ! is based on the same theoretical background as the big particles. |
---|
237 | ! However, because of its low sinking speed, lability is not supposed |
---|
238 | ! to be equal to its initial value (the value of the freshly produced |
---|
239 | ! organic matter). It is however uniform in the mixed layer. |
---|
240 | ! ------------------------------------------------------------------- |
---|
241 | ! |
---|
242 | totprod(:,:) = 0. |
---|
243 | totthick(:,:) = 0. |
---|
244 | totcons(:,:) = 0. |
---|
245 | ! intregrated production and consumption of POC in the mixed layer |
---|
246 | ! ---------------------------------------------------------------- |
---|
247 | ! |
---|
248 | DO jk = 1, jpkm1 |
---|
249 | DO jj = 1, jpj |
---|
250 | DO ji = 1, jpi |
---|
251 | zdep = hmld(ji,jj) |
---|
252 | IF (tmask(ji,jj,jk) == 1. .AND. fsdept(ji,jj,jk) <= zdep ) THEN |
---|
253 | totprod(ji,jj) = totprod(ji,jj) + prodpoc(ji,jj,jk) * fse3t(ji,jj,jk) * rday/ rfact2 |
---|
254 | ! The temperature effect is included here |
---|
255 | totthick(ji,jj) = totthick(ji,jj) + fse3t(ji,jj,jk)* tgfunc(ji,jj,jk) |
---|
256 | totcons(ji,jj) = totcons(ji,jj) - conspoc(ji,jj,jk) * fse3t(ji,jj,jk) * rday/ rfact2 & |
---|
257 | & / ( trb(ji,jj,jk,jppoc) + rtrn ) |
---|
258 | ENDIF |
---|
259 | END DO |
---|
260 | END DO |
---|
261 | END DO |
---|
262 | |
---|
263 | ! Computation of the lability spectrum in the mixed layer. In the mixed |
---|
264 | ! layer, this spectrum is supposed to be uniform. |
---|
265 | ! --------------------------------------------------------------------- |
---|
266 | DO jk = 1, jpkm1 |
---|
267 | DO jj = 1, jpj |
---|
268 | DO ji = 1, jpi |
---|
269 | IF (tmask(ji,jj,jk) == 1.) THEN |
---|
270 | zdep = hmld(ji,jj) |
---|
271 | alphat = 0.0 |
---|
272 | remint = 0.0 |
---|
273 | IF( fsdept(ji,jj,jk) <= zdep ) THEN |
---|
274 | DO jn = 1, jcpoc |
---|
275 | ! For each lability class, the system is supposed to be |
---|
276 | ! at equilibrium: Prod - Sink - w alphap = 0. |
---|
277 | alphap(ji,jj,jk,jn) = totprod(ji,jj) * alphan(jn) / ( reminp(jn) & |
---|
278 | & * totthick(ji,jj) + totcons(ji,jj) + wsbio + rtrn ) |
---|
279 | alphat = alphat + alphap(ji,jj,jk,jn) |
---|
280 | END DO |
---|
281 | DO jn = 1, jcpoc |
---|
282 | alphap(ji,jj,jk,jn) = alphap(ji,jj,jk,jn) / ( alphat + rtrn) |
---|
283 | remint = remint + alphap(ji,jj,jk,jn) * reminp(jn) |
---|
284 | END DO |
---|
285 | ! Mean remineralization rate in the mixed layer |
---|
286 | zremipoc(ji,jj,jk) = MIN(xremip, MAX(0., remint )) |
---|
287 | ENDIF |
---|
288 | ENDIF |
---|
289 | END DO |
---|
290 | END DO |
---|
291 | END DO |
---|
292 | ! |
---|
293 | ! ----------------------------------------------------------------------- |
---|
294 | ! The lability parameterization is used here. The code is here |
---|
295 | ! almost identical to what is done for big particles. The only difference |
---|
296 | ! is that an additional source from GOC to POC is included. This means |
---|
297 | ! that since we need the lability spectrum of GOC, GOC spectrum |
---|
298 | ! should be determined before. |
---|
299 | ! ----------------------------------------------------------------------- |
---|
300 | ! |
---|
301 | DO jk = 2, jpkm1 |
---|
302 | DO jj = 1, jpj |
---|
303 | DO ji = 1, jpi |
---|
304 | IF (tmask(ji,jj,jk) == 1.) THEN |
---|
305 | zdep = hmld(ji,jj) |
---|
306 | IF( fsdept(ji,jj,jk) > zdep ) THEN |
---|
307 | alphat = 0. |
---|
308 | remint = 0. |
---|
309 | ! |
---|
310 | ! the scale factors are corrected with temperature |
---|
311 | zsizek1 = fse3t(ji,jj,jk-1) / 2. / (wsbio3(ji,jj,jk-1) + rtrn) * tgfunc(ji,jj,jk-1) |
---|
312 | zsizek = fse3t(ji,jj,jk) / 2. / (wsbio3(ji,jj,jk) + rtrn) * tgfunc(ji,jj,jk) |
---|
313 | ! |
---|
314 | ! Special treatment of the level just below the MXL |
---|
315 | ! See the comments in the GOC section |
---|
316 | ! --------------------------------------------------- |
---|
317 | ! |
---|
318 | IF ( fsdept(ji,jj,jk-1) <= zdep ) THEN |
---|
319 | ! |
---|
320 | ! Computation of the POC concentration using the |
---|
321 | ! lagrangian algorithm |
---|
322 | zpoc = trb(ji,jj,jk-1,jppoc) + conspoc(ji,jj,jk) * rday / rfact2 & |
---|
323 | & * fse3t(ji,jj,jk) / 2. / (wsbio3(ji,jj,jk) + rtrn) |
---|
324 | zpoc = max(0., zpoc) |
---|
325 | ! |
---|
326 | DO jn = 1, jcpoc |
---|
327 | ! computation of the lability spectrum applying the |
---|
328 | ! different sources and sinks |
---|
329 | alphap(ji,jj,jk,jn) = alphap(ji,jj,jk-1,jn) * exp( -reminp(jn) * zsizek ) * zpoc & |
---|
330 | & + ( prodpoc(ji,jj,jk) * alphan(jn) + zorem3(ji,jj,jk) * alphag(ji,jj,jk,jn) ) & |
---|
331 | & / tgfunc(ji,jj,jk) / reminp(jn) * rday / rfact2 * ( 1. - exp( -reminp(jn) & |
---|
332 | & * zsizek ) ) |
---|
333 | alphap(ji,jj,jk,jn) = MAX( 0., alphap(ji,jj,jk,jn) ) |
---|
334 | alphat = alphat + alphap(ji,jj,jk,jn) |
---|
335 | END DO |
---|
336 | ELSE |
---|
337 | ! |
---|
338 | ! Lability parameterization for the interior of the ocean |
---|
339 | ! This is very similar to what is done in the previous |
---|
340 | ! block |
---|
341 | ! -------------------------------------------------------- |
---|
342 | ! |
---|
343 | zpoc = trb(ji,jj,jk-1,jppoc) + conspoc(ji,jj,jk-1) * rday / rfact2 & |
---|
344 | & * fse3t(ji,jj,jk-1) / 2. / (wsbio3(ji,jj,jk-1) + rtrn) + conspoc(ji,jj,jk) & |
---|
345 | & * rday / rfact2 * fse3t(ji,jj,jk) / 2. / (wsbio3(ji,jj,jk) + rtrn) |
---|
346 | zpoc = max(0., zpoc) |
---|
347 | ! |
---|
348 | DO jn = 1, jcpoc |
---|
349 | alphap(ji,jj,jk,jn) = alphap(ji,jj,jk-1,jn) * exp( -reminp(jn) & |
---|
350 | & * ( zsizek + zsizek1 ) ) * zpoc + ( prodpoc(ji,jj,jk-1) * alphan(jn) & |
---|
351 | & + zorem3(ji,jj,jk-1) * alphag(ji,jj,jk-1,jn) ) * rday / rfact2 / reminp(jn) & |
---|
352 | & / tgfunc(ji,jj,jk-1) * ( 1. - exp( -reminp(jn) * zsizek1 ) ) * exp( -reminp(jn) & |
---|
353 | & * zsizek ) + ( prodpoc(ji,jj,jk) * alphan(jn) + zorem3(ji,jj,jk) & |
---|
354 | & * alphag(ji,jj,jk,jn) ) * rday / rfact2 / reminp(jn) / tgfunc(ji,jj,jk) * ( 1. & |
---|
355 | & - exp( -reminp(jn) * zsizek ) ) |
---|
356 | alphap(ji,jj,jk,jn) = max(0., alphap(ji,jj,jk,jn) ) |
---|
357 | alphat = alphat + alphap(ji,jj,jk,jn) |
---|
358 | END DO |
---|
359 | ENDIF |
---|
360 | ! Normalization of the lability spectrum so that the |
---|
361 | ! integral is equal to 1 |
---|
362 | DO jn = 1, jcpoc |
---|
363 | alphap(ji,jj,jk,jn) = alphap(ji,jj,jk,jn) / ( alphat + rtrn) |
---|
364 | remint = remint + alphap(ji,jj,jk,jn) * reminp(jn) |
---|
365 | END DO |
---|
366 | ! Mean remineralization rate in the water column |
---|
367 | zremipoc(ji,jj,jk) = MIN(xremip, MAX(0., remint )) |
---|
368 | ENDIF |
---|
369 | ENDIF |
---|
370 | END DO |
---|
371 | END DO |
---|
372 | END DO |
---|
373 | #endif |
---|
374 | |
---|
375 | |
---|
376 | DO jk = 1, jpkm1 |
---|
377 | DO jj = 1, jpj |
---|
378 | DO ji = 1, jpi |
---|
379 | IF (tmask(ji,jj,jk) == 1.) THEN |
---|
380 | zstep = xstep |
---|
381 | # if defined key_degrad |
---|
382 | zstep = zstep * facvol(ji,jj,jk) |
---|
383 | # endif |
---|
384 | ! POC disaggregation by turbulence and bacterial activity. |
---|
385 | ! -------------------------------------------------------- |
---|
386 | zremip = zremipoc(ji,jj,jk) * zstep * tgfunc(ji,jj,jk) |
---|
387 | zorem = zremip * trb(ji,jj,jk,jppoc) |
---|
388 | zofer = zremip * trb(ji,jj,jk,jpsfe) |
---|
389 | #if defined key_kriest |
---|
390 | zorem2 = zremip * trb(ji,jj,jk,jpnum) |
---|
391 | #endif |
---|
392 | |
---|
393 | ! Update the appropriate tracers trends |
---|
394 | ! ------------------------------------- |
---|
395 | |
---|
396 | tra(ji,jj,jk,jpdoc) = tra(ji,jj,jk,jpdoc) + zorem |
---|
397 | orem(ji,jj,jk) = orem(ji,jj,jk) + zorem |
---|
398 | tra(ji,jj,jk,jpfer) = tra(ji,jj,jk,jpfer) + zofer |
---|
399 | zfolimi(ji,jj,jk) = zfolimi(ji,jj,jk) + zofer |
---|
400 | tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) - zorem |
---|
401 | #if defined key_kriest |
---|
402 | tra(ji,jj,jk,jpnum) = tra(ji,jj,jk,jpnum) - zorem2 |
---|
403 | #endif |
---|
404 | tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) - zofer |
---|
405 | |
---|
406 | ENDIF |
---|
407 | END DO |
---|
408 | END DO |
---|
409 | END DO |
---|
410 | |
---|
411 | IF( ln_diatrc .AND. lk_iomput .AND. jnt == nrdttrc ) THEN |
---|
412 | zrfact2 = 1.e3 * rfact2r |
---|
413 | CALL iom_put( "REMINP" , zremipoc(:,:,:) * tmask(:,:,:) ) ! Remineralisation rate |
---|
414 | CALL iom_put( "REMING" , zremigoc(:,:,:) * tmask(:,:,:) ) ! Remineralisation rate |
---|
415 | CALL iom_put( "REMINF" , zfolimi(:,:,:) * tmask(:,:,:) * 1.e+3 * rfact2r ) ! Remineralisation rate |
---|
416 | ENDIF |
---|
417 | |
---|
418 | IF(ln_ctl) THEN ! print mean trends (used for debugging) |
---|
419 | WRITE(charout, FMT="('poc2')") |
---|
420 | CALL prt_ctl_trc_info(charout) |
---|
421 | CALL prt_ctl_trc(tab4d=tra, mask=tmask, clinfo=ctrcnm) |
---|
422 | ENDIF |
---|
423 | ! |
---|
424 | CALL wrk_dealloc( jpi, jpj, totprod, totthick, totcons ) |
---|
425 | CALL wrk_dealloc( jpi, jpj, jpk, zremipoc, zremigoc, zorem3, zfolimi ) |
---|
426 | DEALLOCATE( alphag ) |
---|
427 | ! |
---|
428 | IF( nn_timing == 1 ) CALL timing_stop('p4z_poc') |
---|
429 | ! |
---|
430 | END SUBROUTINE p4z_poc |
---|
431 | |
---|
432 | |
---|
433 | SUBROUTINE p4z_poc_init |
---|
434 | !!---------------------------------------------------------------------- |
---|
435 | !! *** ROUTINE p4z_poc_init *** |
---|
436 | !! |
---|
437 | !! ** Purpose : Initialization of remineralization parameters |
---|
438 | !! |
---|
439 | !! ** Method : Read the nampispoc namelist and check the parameters |
---|
440 | !! called at the first timestep |
---|
441 | !! |
---|
442 | !! ** input : Namelist nampispoc |
---|
443 | !! |
---|
444 | !!---------------------------------------------------------------------- |
---|
445 | INTEGER :: jn |
---|
446 | REAL(wp) :: remindelta, reminup, remindown |
---|
447 | INTEGER :: ifault |
---|
448 | |
---|
449 | NAMELIST/nampispoc/ xremip, jcpoc, rshape |
---|
450 | INTEGER :: ios ! Local integer output status for namelist read |
---|
451 | |
---|
452 | REWIND( numnatp_ref ) ! Namelist nampisrem in reference namelist : Pisces remineralization |
---|
453 | READ ( numnatp_ref, nampispoc, IOSTAT = ios, ERR = 901) |
---|
454 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'nampispoc in reference namelist', lwp ) |
---|
455 | |
---|
456 | REWIND( numnatp_cfg ) ! Namelist nampisrem in configuration namelist : Pisces remineralization |
---|
457 | READ ( numnatp_cfg, nampispoc, IOSTAT = ios, ERR = 902 ) |
---|
458 | 902 IF( ios /= 0 ) CALL ctl_nam ( ios , 'nampispoc in configuration namelist', lwp ) |
---|
459 | IF(lwm) WRITE ( numonp, nampispoc ) |
---|
460 | |
---|
461 | IF(lwp) THEN ! control print |
---|
462 | WRITE(numout,*) ' ' |
---|
463 | WRITE(numout,*) ' Namelist parameters for remineralization, nampispoc' |
---|
464 | WRITE(numout,*) ' ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~' |
---|
465 | WRITE(numout,*) ' remineralisation rate of POC xremip =', xremip |
---|
466 | WRITE(numout,*) ' Number of lability classes for POC jcpoc =', jcpoc |
---|
467 | WRITE(numout,*) ' Shape factor of the gamma distribution rshape =', rshape |
---|
468 | ENDIF |
---|
469 | ! |
---|
470 | ! Discretization along the lability space |
---|
471 | ! --------------------------------------- |
---|
472 | ! |
---|
473 | ALLOCATE( alphan(jcpoc), reminp(jcpoc) ) |
---|
474 | ALLOCATE( alphap(jpi,jpj,jpk,jcpoc) ) |
---|
475 | ! |
---|
476 | IF (jcpoc > 1) THEN |
---|
477 | ! |
---|
478 | remindelta = log(4. * 1000. ) / float(jcpoc-1) |
---|
479 | reminup = 1./ 400. * exp(remindelta) |
---|
480 | ! |
---|
481 | ! Discretization based on incomplete gamma functions |
---|
482 | ! As incomplete gamma functions are not available in standard |
---|
483 | ! fortran 95, they have been coded as functions in this module (gamain) |
---|
484 | ! --------------------------------------------------------------------- |
---|
485 | ! |
---|
486 | alphan(1) = gamain(reminup, rshape, ifault) |
---|
487 | reminp(1) = gamain(reminup, rshape+1.0, ifault) * xremip / alphan(1) |
---|
488 | DO jn = 2, jcpoc-1 |
---|
489 | reminup = 1./ 400. * exp(float(jn) * remindelta) |
---|
490 | remindown = 1. / 400. * exp(float(jn-1) * remindelta) |
---|
491 | alphan(jn) = gamain(reminup, rshape, ifault) - gamain(remindown, rshape, ifault) |
---|
492 | reminp(jn) = gamain(reminup, rshape+1.0, ifault) - gamain(remindown, rshape+1.0, ifault) |
---|
493 | reminp(jn) = reminp(jn) * xremip / alphan(jn) |
---|
494 | END DO |
---|
495 | remindown = 1. / 400. * exp(float(jcpoc-1) * remindelta) |
---|
496 | alphan(jcpoc) = 1.0 - gamain(remindown, rshape, ifault) |
---|
497 | reminp(jcpoc) = 1.0 - gamain(remindown, rshape+1.0, ifault) |
---|
498 | reminp(jcpoc) = reminp(jcpoc) * xremip / alphan(jcpoc) |
---|
499 | |
---|
500 | ELSE |
---|
501 | alphan(jcpoc) = 1. |
---|
502 | reminp(jcpoc) = xremip |
---|
503 | ENDIF |
---|
504 | |
---|
505 | DO jn = 1, jcpoc |
---|
506 | alphap(:,:,:,jn) = alphan(jn) |
---|
507 | END DO |
---|
508 | |
---|
509 | END SUBROUTINE p4z_poc_init |
---|
510 | |
---|
511 | REAL FUNCTION alngam( xvalue, ifault ) |
---|
512 | |
---|
513 | !*****************************************************************************80 |
---|
514 | ! |
---|
515 | !! ALNGAM computes the logarithm of the gamma function. |
---|
516 | ! |
---|
517 | ! Modified: |
---|
518 | ! |
---|
519 | ! 13 January 2008 |
---|
520 | ! |
---|
521 | ! Author: |
---|
522 | ! |
---|
523 | ! Allan Macleod |
---|
524 | ! FORTRAN90 version by John Burkardt |
---|
525 | ! |
---|
526 | ! Reference: |
---|
527 | ! |
---|
528 | ! Allan Macleod, |
---|
529 | ! Algorithm AS 245, |
---|
530 | ! A Robust and Reliable Algorithm for the Logarithm of the Gamma Function, |
---|
531 | ! Applied Statistics, |
---|
532 | ! Volume 38, Number 2, 1989, pages 397-402. |
---|
533 | ! |
---|
534 | ! Parameters: |
---|
535 | ! |
---|
536 | ! Input, real ( kind = 8 ) XVALUE, the argument of the Gamma function. |
---|
537 | ! |
---|
538 | ! Output, integer ( kind = 4 ) IFAULT, error flag. |
---|
539 | ! 0, no error occurred. |
---|
540 | ! 1, XVALUE is less than or equal to 0. |
---|
541 | ! 2, XVALUE is too big. |
---|
542 | ! |
---|
543 | ! Output, real ( kind = 8 ) ALNGAM, the logarithm of the gamma function of X. |
---|
544 | ! |
---|
545 | implicit none |
---|
546 | |
---|
547 | real(wp), parameter :: alr2pi = 0.918938533204673E+00 |
---|
548 | integer:: ifault |
---|
549 | real(wp), dimension ( 9 ) :: r1 = (/ & |
---|
550 | -2.66685511495E+00, & |
---|
551 | -24.4387534237E+00, & |
---|
552 | -21.9698958928E+00, & |
---|
553 | 11.1667541262E+00, & |
---|
554 | 3.13060547623E+00, & |
---|
555 | 0.607771387771E+00, & |
---|
556 | 11.9400905721E+00, & |
---|
557 | 31.4690115749E+00, & |
---|
558 | 15.2346874070E+00 /) |
---|
559 | real(wp), dimension ( 9 ) :: r2 = (/ & |
---|
560 | -78.3359299449E+00, & |
---|
561 | -142.046296688E+00, & |
---|
562 | 137.519416416E+00, & |
---|
563 | 78.6994924154E+00, & |
---|
564 | 4.16438922228E+00, & |
---|
565 | 47.0668766060E+00, & |
---|
566 | 313.399215894E+00, & |
---|
567 | 263.505074721E+00, & |
---|
568 | 43.3400022514E+00 /) |
---|
569 | real(wp), dimension ( 9 ) :: r3 = (/ & |
---|
570 | -2.12159572323E+05, & |
---|
571 | 2.30661510616E+05, & |
---|
572 | 2.74647644705E+04, & |
---|
573 | -4.02621119975E+04, & |
---|
574 | -2.29660729780E+03, & |
---|
575 | -1.16328495004E+05, & |
---|
576 | -1.46025937511E+05, & |
---|
577 | -2.42357409629E+04, & |
---|
578 | -5.70691009324E+02 /) |
---|
579 | real(wp), dimension ( 5 ) :: r4 = (/ & |
---|
580 | 0.279195317918525E+00, & |
---|
581 | 0.4917317610505968E+00, & |
---|
582 | 0.0692910599291889E+00, & |
---|
583 | 3.350343815022304E+00, & |
---|
584 | 6.012459259764103E+00 /) |
---|
585 | real (wp) :: x |
---|
586 | real (wp) :: x1 |
---|
587 | real (wp) :: x2 |
---|
588 | real (wp), parameter :: xlge = 5.10E+05 |
---|
589 | real (wp), parameter :: xlgst = 1.0E+30 |
---|
590 | real (wp) :: xvalue |
---|
591 | real (wp) :: y |
---|
592 | |
---|
593 | x = xvalue |
---|
594 | alngam = 0.0E+00 |
---|
595 | ! |
---|
596 | ! Check the input. |
---|
597 | ! |
---|
598 | if ( xlgst <= x ) then |
---|
599 | ifault = 2 |
---|
600 | return |
---|
601 | end if |
---|
602 | if ( x <= 0.0E+00 ) then |
---|
603 | ifault = 1 |
---|
604 | return |
---|
605 | end if |
---|
606 | |
---|
607 | ifault = 0 |
---|
608 | ! |
---|
609 | ! Calculation for 0 < X < 0.5 and 0.5 <= X < 1.5 combined. |
---|
610 | ! |
---|
611 | if ( x < 1.5E+00 ) then |
---|
612 | |
---|
613 | if ( x < 0.5E+00 ) then |
---|
614 | alngam = - log ( x ) |
---|
615 | y = x + 1.0E+00 |
---|
616 | ! |
---|
617 | ! Test whether X < machine epsilon. |
---|
618 | ! |
---|
619 | if ( y == 1.0E+00 ) then |
---|
620 | return |
---|
621 | end if |
---|
622 | |
---|
623 | else |
---|
624 | |
---|
625 | alngam = 0.0E+00 |
---|
626 | y = x |
---|
627 | x = ( x - 0.5E+00 ) - 0.5E+00 |
---|
628 | |
---|
629 | end if |
---|
630 | |
---|
631 | alngam = alngam + x * (((( & |
---|
632 | r1(5) * y & |
---|
633 | + r1(4) ) * y & |
---|
634 | + r1(3) ) * y & |
---|
635 | + r1(2) ) * y & |
---|
636 | + r1(1) ) / (((( & |
---|
637 | y & |
---|
638 | + r1(9) ) * y & |
---|
639 | + r1(8) ) * y & |
---|
640 | + r1(7) ) * y & |
---|
641 | + r1(6) ) |
---|
642 | |
---|
643 | return |
---|
644 | |
---|
645 | end if |
---|
646 | ! |
---|
647 | ! Calculation for 1.5 <= X < 4.0. |
---|
648 | ! |
---|
649 | if ( x < 4.0E+00 ) then |
---|
650 | |
---|
651 | y = ( x - 1.0E+00 ) - 1.0E+00 |
---|
652 | |
---|
653 | alngam = y * (((( & |
---|
654 | r2(5) * x & |
---|
655 | + r2(4) ) * x & |
---|
656 | + r2(3) ) * x & |
---|
657 | + r2(2) ) * x & |
---|
658 | + r2(1) ) / (((( & |
---|
659 | x & |
---|
660 | + r2(9) ) * x & |
---|
661 | + r2(8) ) * x & |
---|
662 | + r2(7) ) * x & |
---|
663 | + r2(6) ) |
---|
664 | ! |
---|
665 | ! Calculation for 4.0 <= X < 12.0. |
---|
666 | ! |
---|
667 | else if ( x < 12.0E+00 ) then |
---|
668 | |
---|
669 | alngam = (((( & |
---|
670 | r3(5) * x & |
---|
671 | + r3(4) ) * x & |
---|
672 | + r3(3) ) * x & |
---|
673 | + r3(2) ) * x & |
---|
674 | + r3(1) ) / (((( & |
---|
675 | x & |
---|
676 | + r3(9) ) * x & |
---|
677 | + r3(8) ) * x & |
---|
678 | + r3(7) ) * x & |
---|
679 | + r3(6) ) |
---|
680 | ! |
---|
681 | ! Calculation for 12.0 <= X. |
---|
682 | ! |
---|
683 | else |
---|
684 | |
---|
685 | y = log ( x ) |
---|
686 | alngam = x * ( y - 1.0E+00 ) - 0.5E+00 * y + alr2pi |
---|
687 | |
---|
688 | if ( x <= xlge ) then |
---|
689 | |
---|
690 | x1 = 1.0E+00 / x |
---|
691 | x2 = x1 * x1 |
---|
692 | |
---|
693 | alngam = alngam + x1 * ( ( & |
---|
694 | r4(3) * & |
---|
695 | x2 + r4(2) ) * & |
---|
696 | x2 + r4(1) ) / ( ( & |
---|
697 | x2 + r4(5) ) * & |
---|
698 | x2 + r4(4) ) |
---|
699 | |
---|
700 | end if |
---|
701 | |
---|
702 | end if |
---|
703 | |
---|
704 | END FUNCTION alngam |
---|
705 | |
---|
706 | REAL FUNCTION gamain( x, p, ifault ) |
---|
707 | |
---|
708 | !*****************************************************************************80 |
---|
709 | ! |
---|
710 | !! GAMAIN computes the incomplete gamma ratio. |
---|
711 | ! |
---|
712 | ! Discussion: |
---|
713 | ! |
---|
714 | ! A series expansion is used if P > X or X <= 1. Otherwise, a |
---|
715 | ! continued fraction approximation is used. |
---|
716 | ! |
---|
717 | ! Modified: |
---|
718 | ! |
---|
719 | ! 17 January 2008 |
---|
720 | ! |
---|
721 | ! Author: |
---|
722 | ! |
---|
723 | ! G Bhattacharjee |
---|
724 | ! FORTRAN90 version by John Burkardt |
---|
725 | ! |
---|
726 | ! Reference: |
---|
727 | ! |
---|
728 | ! G Bhattacharjee, |
---|
729 | ! Algorithm AS 32: |
---|
730 | ! The Incomplete Gamma Integral, |
---|
731 | ! Applied Statistics, |
---|
732 | ! Volume 19, Number 3, 1970, pages 285-287. |
---|
733 | ! |
---|
734 | ! Parameters: |
---|
735 | ! |
---|
736 | ! Input, real ( kind = 8 ) X, P, the parameters of the incomplete |
---|
737 | ! gamma ratio. 0 <= X, and 0 < P. |
---|
738 | ! |
---|
739 | ! Output, integer ( kind = 4 ) IFAULT, error flag. |
---|
740 | ! 0, no errors. |
---|
741 | ! 1, P <= 0. |
---|
742 | ! 2, X < 0. |
---|
743 | ! 3, underflow. |
---|
744 | ! 4, error return from the Log Gamma routine. |
---|
745 | ! |
---|
746 | ! Output, real ( kind = 8 ) GAMAIN, the value of the incomplete |
---|
747 | ! gamma ratio. |
---|
748 | ! |
---|
749 | implicit none |
---|
750 | |
---|
751 | real (wp) a |
---|
752 | real (wp), parameter :: acu = 1.0E-08 |
---|
753 | real (wp) an |
---|
754 | real (wp) arg |
---|
755 | real (wp) b |
---|
756 | real (wp) dif |
---|
757 | real (wp) factor |
---|
758 | real (wp) g |
---|
759 | real (wp) gin |
---|
760 | integer i |
---|
761 | integer ifault |
---|
762 | real (wp), parameter :: oflo = 1.0E+37 |
---|
763 | real (wp) p |
---|
764 | real (wp) pn(6) |
---|
765 | real (wp) rn |
---|
766 | real (wp) term |
---|
767 | real (wp), parameter :: uflo = 1.0E-37 |
---|
768 | real (wp) x |
---|
769 | ! |
---|
770 | ! Check the input. |
---|
771 | ! |
---|
772 | if ( p <= 0.0E+00 ) then |
---|
773 | ifault = 1 |
---|
774 | gamain = 0.0E+00 |
---|
775 | return |
---|
776 | end if |
---|
777 | |
---|
778 | if ( x < 0.0E+00 ) then |
---|
779 | ifault = 2 |
---|
780 | gamain = 0.0E+00 |
---|
781 | return |
---|
782 | end if |
---|
783 | |
---|
784 | if ( x == 0.0E+00 ) then |
---|
785 | ifault = 0 |
---|
786 | gamain = 0.0E+00 |
---|
787 | return |
---|
788 | end if |
---|
789 | |
---|
790 | g = alngam ( p, ifault ) |
---|
791 | |
---|
792 | if ( ifault /= 0 ) then |
---|
793 | ifault = 4 |
---|
794 | gamain = 0.0E+00 |
---|
795 | return |
---|
796 | end if |
---|
797 | |
---|
798 | arg = p * log ( x ) - x - g |
---|
799 | |
---|
800 | if ( arg < log ( uflo ) ) then |
---|
801 | ifault = 3 |
---|
802 | gamain = 0.0E+00 |
---|
803 | return |
---|
804 | end if |
---|
805 | |
---|
806 | ifault = 0 |
---|
807 | factor = exp ( arg ) |
---|
808 | ! |
---|
809 | ! Calculation by series expansion. |
---|
810 | ! |
---|
811 | if ( x <= 1.0E+00 .or. x < p ) then |
---|
812 | |
---|
813 | gin = 1.0E+00 |
---|
814 | term = 1.0E+00 |
---|
815 | rn = p |
---|
816 | |
---|
817 | do |
---|
818 | |
---|
819 | rn = rn + 1.0E+00 |
---|
820 | term = term * x / rn |
---|
821 | gin = gin + term |
---|
822 | |
---|
823 | if ( term <= acu ) then |
---|
824 | exit |
---|
825 | end if |
---|
826 | |
---|
827 | end do |
---|
828 | |
---|
829 | gamain = gin * factor / p |
---|
830 | return |
---|
831 | |
---|
832 | end if |
---|
833 | ! |
---|
834 | ! Calculation by continued fraction. |
---|
835 | ! |
---|
836 | a = 1.0E+00 - p |
---|
837 | b = a + x + 1.0E+00 |
---|
838 | term = 0.0E+00 |
---|
839 | |
---|
840 | pn(1) = 1.0E+00 |
---|
841 | pn(2) = x |
---|
842 | pn(3) = x + 1.0E+00 |
---|
843 | pn(4) = x * b |
---|
844 | |
---|
845 | gin = pn(3) / pn(4) |
---|
846 | |
---|
847 | do |
---|
848 | |
---|
849 | a = a + 1.0E+00 |
---|
850 | b = b + 2.0E+00 |
---|
851 | term = term + 1.0E+00 |
---|
852 | an = a * term |
---|
853 | do i = 1, 2 |
---|
854 | pn(i+4) = b * pn(i+2) - an * pn(i) |
---|
855 | end do |
---|
856 | |
---|
857 | if ( pn(6) /= 0.0E+00 ) then |
---|
858 | |
---|
859 | rn = pn(5) / pn(6) |
---|
860 | dif = abs ( gin - rn ) |
---|
861 | ! |
---|
862 | ! Absolute error tolerance satisfied? |
---|
863 | ! |
---|
864 | if ( dif <= acu ) then |
---|
865 | ! |
---|
866 | ! Relative error tolerance satisfied? |
---|
867 | ! |
---|
868 | if ( dif <= acu * rn ) then |
---|
869 | gamain = 1.0E+00 - factor * gin |
---|
870 | exit |
---|
871 | end if |
---|
872 | |
---|
873 | end if |
---|
874 | |
---|
875 | gin = rn |
---|
876 | |
---|
877 | end if |
---|
878 | |
---|
879 | do i = 1, 4 |
---|
880 | pn(i) = pn(i+2) |
---|
881 | end do |
---|
882 | if ( oflo <= abs ( pn(5) ) ) then |
---|
883 | |
---|
884 | do i = 1, 4 |
---|
885 | pn(i) = pn(i) / oflo |
---|
886 | end do |
---|
887 | |
---|
888 | end if |
---|
889 | |
---|
890 | end do |
---|
891 | |
---|
892 | END FUNCTION gamain |
---|
893 | |
---|
894 | #else |
---|
895 | !!====================================================================== |
---|
896 | !! Dummy module : No PISCES bio-model |
---|
897 | !!====================================================================== |
---|
898 | CONTAINS |
---|
899 | SUBROUTINE p4z_poc ! Empty routine |
---|
900 | END SUBROUTINE p4z_poc |
---|
901 | #endif |
---|
902 | |
---|
903 | !!====================================================================== |
---|
904 | END MODULE p4zpoc |
---|