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

Last change on this file since 8003 was 8003, checked in by aumont, 7 years ago

modification in the code to remove unnecessary parts such as kriest and non iomput options

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