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.
p5zpoc.F90 in branches/CNRS/dev_r6270_PISCES_QUOTA/NEMOGCM/NEMO/TOP_SRC/PISCES/P5Z – NEMO

source: branches/CNRS/dev_r6270_PISCES_QUOTA/NEMOGCM/NEMO/TOP_SRC/PISCES/P5Z/p5zpoc.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

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