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.
p5zmort.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/p5zmort.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: 15.2 KB
Line 
1MODULE p5zmort
2   !!======================================================================
3   !!                         ***  MODULE p5zmort  ***
4   !! TOP :   PISCES Compute the mortality terms for phytoplankton
5   !!======================================================================
6   !! History :   1.0  !  2002     (O. Aumont)  Original code
7   !!             2.0  !  2007-12  (C. Ethe, G. Madec)  F90
8   !!             3.6  !  2015-05  (O. Aumont) PISCES quota
9   !!----------------------------------------------------------------------
10#if defined key_pisces_quota
11   !!----------------------------------------------------------------------
12   !!   'key_pisces_quota'                                 PISCES bio-model
13   !!----------------------------------------------------------------------
14   !!   p5z_mort       :   Compute the mortality terms for phytoplankton
15   !!   p5z_mort_init  :   Initialize the mortality params for phytoplankton
16   !!----------------------------------------------------------------------
17   USE oce_trc         !  shared variables between ocean and passive tracers
18   USE trc             !  passive tracers common variables
19   USE sms_pisces      !  PISCES Source Minus Sink variables
20   USE p5zlim          !  Phytoplankton limitation terms
21   USE prtctl_trc      !  print control for debugging
22
23   IMPLICIT NONE
24   PRIVATE
25
26   PUBLIC   p5z_mort   
27   PUBLIC   p5z_mort_init   
28
29   !! * Shared module variables
30   REAL(wp), PUBLIC :: wchl    !:
31   REAL(wp), PUBLIC :: wchlp   !:
32   REAL(wp), PUBLIC :: wchld   !:
33   REAL(wp), PUBLIC :: wchldm  !:
34   REAL(wp), PUBLIC :: mprat   !:
35   REAL(wp), PUBLIC :: mpratp  !:
36   REAL(wp), PUBLIC :: mprat2  !:
37
38
39   !!* Substitution
40#  include "top_substitute.h90"
41   !!----------------------------------------------------------------------
42   !! NEMO/TOP 3.3 , NEMO Consortium (2010)
43   !! $Id: p4zmort.F90 3160 2011-11-20 14:27:18Z cetlod $
44   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
45   !!----------------------------------------------------------------------
46
47CONTAINS
48
49   SUBROUTINE p5z_mort( kt )
50      !!---------------------------------------------------------------------
51      !!                     ***  ROUTINE p5z_mort  ***
52      !!
53      !! ** Purpose :   Calls the different subroutine to initialize and compute
54      !!                the different phytoplankton mortality terms
55      !!
56      !! ** Method  : - ???
57      !!---------------------------------------------------------------------
58      INTEGER, INTENT(in) ::   kt ! ocean time step
59      !!---------------------------------------------------------------------
60
61      CALL p5z_nano            ! nanophytoplankton
62      CALL p5z_pico            ! picophytoplankton
63      CALL p5z_diat            ! diatoms
64
65   END SUBROUTINE p5z_mort
66
67
68   SUBROUTINE p5z_nano
69      !!---------------------------------------------------------------------
70      !!                     ***  ROUTINE p5z_nano  ***
71      !!
72      !! ** Purpose :   Compute the mortality terms for nanophytoplankton
73      !!
74      !! ** Method  : - ???
75      !!---------------------------------------------------------------------
76      INTEGER  :: ji, jj, jk
77      REAL(wp) :: zcompaph
78      REAL(wp) :: zfactfe, zfactch, zfactn, zfactp, zprcaca
79      REAL(wp) :: ztortp , zrespp , zmortp , zstep
80      CHARACTER (len=25) :: charout
81      !!---------------------------------------------------------------------
82      !
83      IF( nn_timing == 1 )  CALL timing_start('p5z_nano')
84      !
85      prodcal(:,:,:) = 0.  !: calcite production variable set to zero
86      DO jk = 1, jpkm1
87         DO jj = 1, jpj
88            DO ji = 1, jpi
89               zcompaph = MAX( ( trb(ji,jj,jk,jpphy) - 1e-9 ), 0.e0 )
90               zstep    = xstep
91               !   Squared mortality of Phyto similar to a sedimentation term during
92               !   blooms (Doney et al. 1996)
93               !   -----------------------------------------------------------------
94               zrespp = wchl * 1.e6 * zstep * xdiss(ji,jj,jk) * zcompaph * trb(ji,jj,jk,jpphy)
95
96               !   Phytoplankton linear mortality
97               !   ------------------------------
98               ztortp = mprat * xstep  * zcompaph
99               zmortp = zrespp + ztortp
100
101               !   Update the arrays TRA which contains the biological sources and sinks
102
103               zfactn  = trb(ji,jj,jk,jpnph)/(trb(ji,jj,jk,jpphy)+rtrn)
104               zfactp  = trb(ji,jj,jk,jppph)/(trb(ji,jj,jk,jpphy)+rtrn)
105               zfactfe = trb(ji,jj,jk,jpnfe)/(trb(ji,jj,jk,jpphy)+rtrn)
106               zfactch = trb(ji,jj,jk,jpnch)/(trb(ji,jj,jk,jpphy)+rtrn)
107               tra(ji,jj,jk,jpphy) = tra(ji,jj,jk,jpphy) - zmortp
108               tra(ji,jj,jk,jpnph) = tra(ji,jj,jk,jpnph) - zmortp * zfactn
109               tra(ji,jj,jk,jppph) = tra(ji,jj,jk,jppph) - zmortp * zfactp
110               tra(ji,jj,jk,jpnch) = tra(ji,jj,jk,jpnch) - zmortp * zfactch
111               tra(ji,jj,jk,jpnfe) = tra(ji,jj,jk,jpnfe) - zmortp * zfactfe
112               zprcaca = xfracal(ji,jj,jk) * zmortp
113               !
114               prodcal(ji,jj,jk) = prodcal(ji,jj,jk) + zprcaca  ! prodcal=prodcal(nanophy)+prodcal(microzoo)+prodcal(mesozoo)
115               !
116               tra(ji,jj,jk,jpdic) = tra(ji,jj,jk,jpdic) - zprcaca
117               tra(ji,jj,jk,jptal) = tra(ji,jj,jk,jptal) - 2. * zprcaca
118               tra(ji,jj,jk,jpcal) = tra(ji,jj,jk,jpcal) + zprcaca
119               tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + zmortp
120               tra(ji,jj,jk,jppon) = tra(ji,jj,jk,jppon) + zmortp * zfactn
121               tra(ji,jj,jk,jppop) = tra(ji,jj,jk,jppop) + zmortp * zfactp
122               prodpoc(ji,jj,jk) = prodpoc(ji,jj,jk) + zmortp
123               tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + zmortp * zfactfe
124            END DO
125         END DO
126      END DO
127      !
128       IF(ln_ctl)   THEN  ! print mean trends (used for debugging)
129         WRITE(charout, FMT="('nano')")
130         CALL prt_ctl_trc_info(charout)
131         CALL prt_ctl_trc(tab4d=tra, mask=tmask, clinfo=ctrcnm)
132       ENDIF
133      !
134      IF( nn_timing == 1 )  CALL timing_stop('p5z_nano')
135      !
136   END SUBROUTINE p5z_nano
137
138   SUBROUTINE p5z_pico
139      !!---------------------------------------------------------------------
140      !!                     ***  ROUTINE p5z_pico  ***
141      !!
142      !! ** Purpose :   Compute the mortality terms for picophytoplankton
143      !!
144      !! ** Method  : - ???
145      !!---------------------------------------------------------------------
146      INTEGER  :: ji, jj, jk
147      REAL(wp) :: zcompaph
148      REAL(wp) :: zfactfe, zfactch, zfactn, zfactp
149      REAL(wp) :: ztortp , zrespp , zmortp , zstep
150      CHARACTER (len=25) :: charout
151      !!---------------------------------------------------------------------
152      !
153      IF( nn_timing == 1 )  CALL timing_start('p5z_pico')
154      !
155      DO jk = 1, jpkm1
156         DO jj = 1, jpj
157            DO ji = 1, jpi
158               zcompaph = MAX( ( trb(ji,jj,jk,jppic) - 1e-9 ), 0.e0 )
159               zstep    = xstep
160               !  Squared mortality of Phyto similar to a sedimentation term during
161               !  blooms (Doney et al. 1996)
162               !  -----------------------------------------------------------------
163               zrespp = wchlp * 1.e6 * zstep * xdiss(ji,jj,jk) * zcompaph * trb(ji,jj,jk,jppic)
164
165               !     Phytoplankton mortality
166               ztortp = mpratp * xstep  * zcompaph
167               zmortp = zrespp + ztortp
168
169               !   Update the arrays TRA which contains the biological sources and sinks
170
171               zfactn = trb(ji,jj,jk,jpnpi)/(trb(ji,jj,jk,jppic)+rtrn)
172               zfactp = trb(ji,jj,jk,jpppi)/(trb(ji,jj,jk,jppic)+rtrn)
173               zfactfe = trb(ji,jj,jk,jppfe)/(trb(ji,jj,jk,jppic)+rtrn)
174               zfactch = trb(ji,jj,jk,jppch)/(trb(ji,jj,jk,jppic)+rtrn)
175               tra(ji,jj,jk,jppic) = tra(ji,jj,jk,jppic) - zmortp
176               tra(ji,jj,jk,jpnpi) = tra(ji,jj,jk,jpnpi) - zmortp * zfactn
177               tra(ji,jj,jk,jpppi) = tra(ji,jj,jk,jpppi) - zmortp * zfactp
178               tra(ji,jj,jk,jppch) = tra(ji,jj,jk,jppch) - zmortp * zfactch
179               tra(ji,jj,jk,jppfe) = tra(ji,jj,jk,jppfe) - zmortp * zfactfe
180               tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + zmortp
181               tra(ji,jj,jk,jppon) = tra(ji,jj,jk,jppon) + zmortp * zfactn
182               tra(ji,jj,jk,jppop) = tra(ji,jj,jk,jppop) + zmortp * zfactp
183               tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + zmortp * zfactfe
184               prodpoc(ji,jj,jk) = prodpoc(ji,jj,jk) + zmortp
185            END DO
186         END DO
187      END DO
188      !
189       IF(ln_ctl)   THEN  ! print mean trends (used for debugging)
190         WRITE(charout, FMT="('pico')")
191         CALL prt_ctl_trc_info(charout)
192         CALL prt_ctl_trc(tab4d=tra, mask=tmask, clinfo=ctrcnm)
193       ENDIF
194      !
195      IF( nn_timing == 1 )  CALL timing_stop('p5z_pico')
196      !
197   END SUBROUTINE p5z_pico
198
199
200   SUBROUTINE p5z_diat
201      !!---------------------------------------------------------------------
202      !!                     ***  ROUTINE p5z_diat  ***
203      !!
204      !! ** Purpose :   Compute the mortality terms for diatoms
205      !!
206      !! ** Method  : - ???
207      !!---------------------------------------------------------------------
208      INTEGER  ::  ji, jj, jk
209      REAL(wp) ::  zfactfe,zfactsi,zfactch, zfactn, zfactp, zcompadi
210      REAL(wp) ::  zrespp2, ztortp2, zmortp2, zstep
211      REAL(wp) ::  zlim2, zlim1
212      CHARACTER (len=25) :: charout
213      !!---------------------------------------------------------------------
214      !
215      IF( nn_timing == 1 )  CALL timing_start('p5z_diat')
216      !
217
218      DO jk = 1, jpkm1
219         DO jj = 1, jpj
220            DO ji = 1, jpi
221
222               zcompadi = MAX( ( trb(ji,jj,jk,jpdia) - 1E-9), 0. )
223
224               !   Aggregation term for diatoms is increased in case of nutrient
225               !   stress as observed in reality. The stressed cells become more
226               !   sticky and coagulate to sink quickly out of the euphotic zone
227               !   -------------------------------------------------------------
228               zstep   = xstep
229               !  Phytoplankton squared mortality
230               !  -------------------------------
231               zlim2   = xlimdia(ji,jj,jk) * xlimdia(ji,jj,jk)
232               zlim1   = 0.25 * ( 1. - zlim2 ) / ( 0.25 + zlim2 ) 
233               zrespp2 = 1.e6 * zstep * (  wchld + wchldm * zlim1 ) * xdiss(ji,jj,jk) * zcompadi * trb(ji,jj,jk,jpdia)
234
235               !  Phytoplankton linear mortality
236               !  ------------------------------
237               ztortp2 = mprat2 * xstep  * zcompadi
238               zmortp2 = zrespp2 + ztortp2
239
240               !   Update the arrays tra which contains the biological sources and sinks
241               !   ---------------------------------------------------------------------
242               zfactn  = trb(ji,jj,jk,jpndi) / ( trb(ji,jj,jk,jpdia) + rtrn )
243               zfactp  = trb(ji,jj,jk,jppdi) / ( trb(ji,jj,jk,jpdia) + rtrn )
244               zfactch = trb(ji,jj,jk,jpdch) / ( trb(ji,jj,jk,jpdia) + rtrn )
245               zfactfe = trb(ji,jj,jk,jpdfe) / ( trb(ji,jj,jk,jpdia) + rtrn )
246               zfactsi = trb(ji,jj,jk,jpdsi) / ( trb(ji,jj,jk,jpdia) + rtrn )
247               tra(ji,jj,jk,jpdia) = tra(ji,jj,jk,jpdia) - zmortp2 
248               tra(ji,jj,jk,jpndi) = tra(ji,jj,jk,jpndi) - zmortp2 * zfactn
249               tra(ji,jj,jk,jppdi) = tra(ji,jj,jk,jppdi) - zmortp2 * zfactp
250               tra(ji,jj,jk,jpdch) = tra(ji,jj,jk,jpdch) - zmortp2 * zfactch
251               tra(ji,jj,jk,jpdfe) = tra(ji,jj,jk,jpdfe) - zmortp2 * zfactfe
252               tra(ji,jj,jk,jpdsi) = tra(ji,jj,jk,jpdsi) - zmortp2 * zfactsi
253               tra(ji,jj,jk,jpgsi) = tra(ji,jj,jk,jpgsi) + zmortp2 * zfactsi
254               tra(ji,jj,jk,jpgoc) = tra(ji,jj,jk,jpgoc) + zrespp2 
255               tra(ji,jj,jk,jpgon) = tra(ji,jj,jk,jpgon) + zrespp2 * zfactn
256               tra(ji,jj,jk,jpgop) = tra(ji,jj,jk,jpgop) + zrespp2 * zfactp
257               tra(ji,jj,jk,jpbfe) = tra(ji,jj,jk,jpbfe) + zrespp2 * zfactfe
258               tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + ztortp2
259               tra(ji,jj,jk,jppon) = tra(ji,jj,jk,jppon) + ztortp2 * zfactn
260               tra(ji,jj,jk,jppop) = tra(ji,jj,jk,jppop) + ztortp2 * zfactp
261               tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + ztortp2 * zfactfe
262               prodpoc(ji,jj,jk) = prodpoc(ji,jj,jk) + ztortp2
263               prodgoc(ji,jj,jk) = prodgoc(ji,jj,jk) + zrespp2
264            END DO
265         END DO
266      END DO
267      !
268      IF(ln_ctl)   THEN  ! print mean trends (used for debugging)
269         WRITE(charout, FMT="('diat')")
270         CALL prt_ctl_trc_info(charout)
271         CALL prt_ctl_trc(tab4d=tra, mask=tmask, clinfo=ctrcnm)
272      ENDIF
273      !
274      IF( nn_timing == 1 )  CALL timing_stop('p5z_diat')
275      !
276   END SUBROUTINE p5z_diat
277
278   SUBROUTINE p5z_mort_init
279
280      !!----------------------------------------------------------------------
281      !!                  ***  ROUTINE p5z_mort_init  ***
282      !!
283      !! ** Purpose :   Initialization of phytoplankton parameters
284      !!
285      !! ** Method  :   Read the nampismort namelist and check the parameters
286      !!      called at the first timestep
287      !!
288      !! ** input   :   Namelist nampismort
289      !!
290      !!----------------------------------------------------------------------
291
292      NAMELIST/nampismort/ wchl, wchlp, wchld, wchldm, mprat, mpratp, mprat2
293      INTEGER :: ios                 ! Local integer output status for namelist read
294
295      REWIND( numnatp_ref )              ! Namelist nampismort in reference namelist : Pisces phytoplankton
296      READ  ( numnatp_ref, nampismort, IOSTAT = ios, ERR = 901)
297901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nampismort in reference namelist', lwp )
298
299      REWIND( numnatp_cfg )              ! Namelist nampismort in configuration namelist : Pisces phytoplankton
300      READ  ( numnatp_cfg, nampismort, IOSTAT = ios, ERR = 902 )
301902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nampismort in configuration namelist', lwp )
302      IF(lwm) WRITE ( numonp, nampismort )
303
304      IF(lwp) THEN                         ! control print
305         WRITE(numout,*) ' '
306         WRITE(numout,*) ' Namelist parameters for phytoplankton mortality, nampismort'
307         WRITE(numout,*) ' ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~'
308         WRITE(numout,*) '    quadratic mortality of phytoplankton      wchl      =', wchl
309         WRITE(numout,*) '    quadratic mortality of picophyto.         wchlp     =', wchlp
310         WRITE(numout,*) '    quadratic mortality of diatoms            wchld     =', wchld
311         WRITE(numout,*) '    Additional quadratic mortality of diatoms wchldm    =', wchldm
312         WRITE(numout,*) '    nanophyto. mortality rate                 mprat     =', mprat
313         WRITE(numout,*) '    picophyto. mortality rate                 mpratp    =', mpratp
314         WRITE(numout,*) '    Diatoms mortality rate                    mprat2    =', mprat2
315      ENDIF
316
317   END SUBROUTINE p5z_mort_init
318
319#else
320   !!======================================================================
321   !!  Dummy module :                                   No PISCES bio-model
322   !!======================================================================
323CONTAINS
324   SUBROUTINE p5z_mort                    ! Empty routine
325   END SUBROUTINE p5z_mort
326#endif 
327
328   !!======================================================================
329END MODULE  p5zmort
Note: See TracBrowser for help on using the repository browser.