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 @ 7180

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

various bug fixes on iron chemistry

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