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 branches/CNRS/dev_r6270_PISCES_QUOTA/NEMOGCM/NEMO/TOP_SRC/PISCES/P4Z – NEMO

source: branches/CNRS/dev_r6270_PISCES_QUOTA/NEMOGCM/NEMO/TOP_SRC/PISCES/P4Z/p4zpoc.F90 @ 6966

Last change on this file since 6966 was 6966, checked in by aumont, 8 years ago

various bug fixes and updates on carbon chemistry

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