New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
p4zpoc.F90 in NEMO/trunk/src/TOP/PISCES/P4Z – NEMO

source: NEMO/trunk/src/TOP/PISCES/P4Z/p4zpoc.F90

Last change on this file was 15459, checked in by cetlod, 2 years ago

Various bug fixes and more comments in PISCES routines ; sette test OK in debug mode, nn_hls=1/2, with tiling ; run.stat unchanged ; of course tracer.stat different

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