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 NEMO/branches/2021/dev_r14116_HPC-10_mcastril_Mixed_Precision_implementation/src/TOP/PISCES/P4Z – NEMO

source: NEMO/branches/2021/dev_r14116_HPC-10_mcastril_Mixed_Precision_implementation/src/TOP/PISCES/P4Z/p5zmort.F90 @ 15814

Last change on this file since 15814 was 15648, checked in by sparonuz, 2 years ago

Updated name preprocessor function CASTWP to CASTDP

  • Property svn:keywords set to Id
File size: 14.5 KB
RevLine 
[7162]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   !!   p5z_mort       :   Compute the mortality terms for phytoplankton
11   !!   p5z_mort_init  :   Initialize the mortality params for phytoplankton
12   !!----------------------------------------------------------------------
13   USE oce_trc         !  shared variables between ocean and passive tracers
14   USE trc             !  passive tracers common variables
15   USE sms_pisces      !  PISCES Source Minus Sink variables
[10362]16   USE p4zlim
[10227]17   USE p5zlim          !  Phytoplankton limitation terms
[13286]18   USE prtctl          !  print control for debugging
[7162]19
20   IMPLICIT NONE
21   PRIVATE
22
23   PUBLIC   p5z_mort   
24   PUBLIC   p5z_mort_init   
25
26   !! * Shared module variables
27   REAL(wp), PUBLIC :: wchln    !:
28   REAL(wp), PUBLIC :: wchlp   !:
29   REAL(wp), PUBLIC :: wchld   !:
30   REAL(wp), PUBLIC :: wchldm  !:
31   REAL(wp), PUBLIC :: mpratn   !:
32   REAL(wp), PUBLIC :: mpratp  !:
33   REAL(wp), PUBLIC :: mpratd  !:
34
[12377]35   !! * Substitutions
36#  include "do_loop_substitute.h90"
[14219]37#  include "single_precision_substitute.h90"
[7162]38   !!----------------------------------------------------------------------
[10067]39   !! NEMO/TOP 4.0 , NEMO Consortium (2018)
[10068]40   !! $Id$
41   !! Software governed by the CeCILL license (see ./LICENSE)
[7162]42   !!----------------------------------------------------------------------
43
44CONTAINS
45
[12377]46   SUBROUTINE p5z_mort( kt, Kbb, Krhs )
[7162]47      !!---------------------------------------------------------------------
48      !!                     ***  ROUTINE p5z_mort  ***
49      !!
50      !! ** Purpose :   Calls the different subroutine to initialize and compute
51      !!                the different phytoplankton mortality terms
52      !!
53      !! ** Method  : - ???
54      !!---------------------------------------------------------------------
55      INTEGER, INTENT(in) ::   kt ! ocean time step
[12377]56      INTEGER, INTENT(in) ::   Kbb, Krhs  ! time level indices
[7162]57      !!---------------------------------------------------------------------
58
[12377]59      CALL p5z_nano( Kbb, Krhs )            ! nanophytoplankton
60      CALL p5z_pico( Kbb, Krhs )            ! picophytoplankton
61      CALL p5z_diat( Kbb, Krhs )            ! diatoms
[7162]62
63   END SUBROUTINE p5z_mort
64
65
[12377]66   SUBROUTINE p5z_nano( Kbb, Krhs )
[7162]67      !!---------------------------------------------------------------------
68      !!                     ***  ROUTINE p5z_nano  ***
69      !!
70      !! ** Purpose :   Compute the mortality terms for nanophytoplankton
71      !!
72      !! ** Method  : - ???
73      !!---------------------------------------------------------------------
[12377]74      INTEGER, INTENT(in) ::   Kbb, Krhs  ! time level indices
[7162]75      INTEGER  :: ji, jj, jk
76      REAL(wp) :: zcompaph
77      REAL(wp) :: zfactfe, zfactch, zfactn, zfactp, zprcaca
78      REAL(wp) :: ztortp , zrespp , zmortp
79      CHARACTER (len=25) :: charout
80      !!---------------------------------------------------------------------
81      !
[9124]82      IF( ln_timing )   CALL timing_start('p5z_nano')
[7162]83      !
84      prodcal(:,:,:) = 0.  !: calcite production variable set to zero
[13295]85      DO_3D( 1, 1, 1, 1, 1, jpkm1 )
[12377]86         zcompaph = MAX( ( tr(ji,jj,jk,jpphy,Kbb) - 1e-9 ), 0.e0 )
87         !   Squared mortality of Phyto similar to a sedimentation term during
88         !   blooms (Doney et al. 1996)
89         !   -----------------------------------------------------------------
90         zrespp = wchln * 1.e6 * xstep * xdiss(ji,jj,jk) * zcompaph * tr(ji,jj,jk,jpphy,Kbb)
[7162]91
[12377]92         !   Phytoplankton linear mortality
93         !   ------------------------------
94         ztortp = mpratn * xstep  * zcompaph
95         zmortp = zrespp + ztortp
[7162]96
[12377]97         !   Update the arrays TRA which contains the biological sources and sinks
[7162]98
[12377]99         zfactn  = tr(ji,jj,jk,jpnph,Kbb)/(tr(ji,jj,jk,jpphy,Kbb)+rtrn)
100         zfactp  = tr(ji,jj,jk,jppph,Kbb)/(tr(ji,jj,jk,jpphy,Kbb)+rtrn)
101         zfactfe = tr(ji,jj,jk,jpnfe,Kbb)/(tr(ji,jj,jk,jpphy,Kbb)+rtrn)
102         zfactch = tr(ji,jj,jk,jpnch,Kbb)/(tr(ji,jj,jk,jpphy,Kbb)+rtrn)
103         tr(ji,jj,jk,jpphy,Krhs) = tr(ji,jj,jk,jpphy,Krhs) - zmortp
104         tr(ji,jj,jk,jpnph,Krhs) = tr(ji,jj,jk,jpnph,Krhs) - zmortp * zfactn
105         tr(ji,jj,jk,jppph,Krhs) = tr(ji,jj,jk,jppph,Krhs) - zmortp * zfactp
106         tr(ji,jj,jk,jpnch,Krhs) = tr(ji,jj,jk,jpnch,Krhs) - zmortp * zfactch
107         tr(ji,jj,jk,jpnfe,Krhs) = tr(ji,jj,jk,jpnfe,Krhs) - zmortp * zfactfe
108         zprcaca = xfracal(ji,jj,jk) * zmortp
109         !
110         prodcal(ji,jj,jk) = prodcal(ji,jj,jk) + zprcaca  ! prodcal=prodcal(nanophy)+prodcal(microzoo)+prodcal(mesozoo)
111         !
112         tr(ji,jj,jk,jpdic,Krhs) = tr(ji,jj,jk,jpdic,Krhs) - zprcaca
113         tr(ji,jj,jk,jptal,Krhs) = tr(ji,jj,jk,jptal,Krhs) - 2. * zprcaca
114         tr(ji,jj,jk,jpcal,Krhs) = tr(ji,jj,jk,jpcal,Krhs) + zprcaca
115         tr(ji,jj,jk,jppoc,Krhs) = tr(ji,jj,jk,jppoc,Krhs) + zmortp
116         tr(ji,jj,jk,jppon,Krhs) = tr(ji,jj,jk,jppon,Krhs) + zmortp * zfactn
117         tr(ji,jj,jk,jppop,Krhs) = tr(ji,jj,jk,jppop,Krhs) + zmortp * zfactp
118         prodpoc(ji,jj,jk) = prodpoc(ji,jj,jk) + zmortp
119         tr(ji,jj,jk,jpsfe,Krhs) = tr(ji,jj,jk,jpsfe,Krhs) + zmortp * zfactfe
120      END_3D
[7162]121      !
[12377]122       IF(sn_cfctl%l_prttrc)   THEN  ! print mean trends (used for debugging)
[7162]123         WRITE(charout, FMT="('nano')")
[13286]124         CALL prt_ctl_info( charout, cdcomp = 'top' )
[15648]125         CALL prt_ctl(tab4d_1=CASTDP(tr(:,:,:,:,Krhs)), mask1=tmask, clinfo=ctrcnm)
[7162]126       ENDIF
127      !
[9124]128      IF( ln_timing )   CALL timing_stop('p5z_nano')
[7162]129      !
130   END SUBROUTINE p5z_nano
131
[9124]132
[12377]133   SUBROUTINE p5z_pico( Kbb, Krhs )
[7162]134      !!---------------------------------------------------------------------
135      !!                     ***  ROUTINE p5z_pico  ***
136      !!
137      !! ** Purpose :   Compute the mortality terms for picophytoplankton
138      !!
139      !! ** Method  : - ???
140      !!---------------------------------------------------------------------
[12377]141      INTEGER, INTENT(in) ::   Kbb, Krhs  ! time level indices
[7162]142      INTEGER  :: ji, jj, jk
143      REAL(wp) :: zcompaph
144      REAL(wp) :: zfactfe, zfactch, zfactn, zfactp
145      REAL(wp) :: ztortp , zrespp , zmortp 
146      CHARACTER (len=25) :: charout
147      !!---------------------------------------------------------------------
148      !
[9124]149      IF( ln_timing )   CALL timing_start('p5z_pico')
[7162]150      !
[13295]151      DO_3D( 1, 1, 1, 1, 1, jpkm1 )
[12377]152         zcompaph = MAX( ( tr(ji,jj,jk,jppic,Kbb) - 1e-9 ), 0.e0 )
153         !  Squared mortality of Phyto similar to a sedimentation term during
154         !  blooms (Doney et al. 1996)
155         !  -----------------------------------------------------------------
156         zrespp = wchlp * 1.e6 * xstep * xdiss(ji,jj,jk) * zcompaph * tr(ji,jj,jk,jppic,Kbb)
[7162]157
[12377]158         !     Phytoplankton mortality
159         ztortp = mpratp * xstep  * zcompaph
160         zmortp = zrespp + ztortp
[7162]161
[12377]162         !   Update the arrays TRA which contains the biological sources and sinks
[7162]163
[12377]164         zfactn = tr(ji,jj,jk,jpnpi,Kbb)/(tr(ji,jj,jk,jppic,Kbb)+rtrn)
165         zfactp = tr(ji,jj,jk,jpppi,Kbb)/(tr(ji,jj,jk,jppic,Kbb)+rtrn)
166         zfactfe = tr(ji,jj,jk,jppfe,Kbb)/(tr(ji,jj,jk,jppic,Kbb)+rtrn)
167         zfactch = tr(ji,jj,jk,jppch,Kbb)/(tr(ji,jj,jk,jppic,Kbb)+rtrn)
168         tr(ji,jj,jk,jppic,Krhs) = tr(ji,jj,jk,jppic,Krhs) - zmortp
169         tr(ji,jj,jk,jpnpi,Krhs) = tr(ji,jj,jk,jpnpi,Krhs) - zmortp * zfactn
170         tr(ji,jj,jk,jpppi,Krhs) = tr(ji,jj,jk,jpppi,Krhs) - zmortp * zfactp
171         tr(ji,jj,jk,jppch,Krhs) = tr(ji,jj,jk,jppch,Krhs) - zmortp * zfactch
172         tr(ji,jj,jk,jppfe,Krhs) = tr(ji,jj,jk,jppfe,Krhs) - zmortp * zfactfe
173         tr(ji,jj,jk,jppoc,Krhs) = tr(ji,jj,jk,jppoc,Krhs) + zmortp
174         tr(ji,jj,jk,jppon,Krhs) = tr(ji,jj,jk,jppon,Krhs) + zmortp * zfactn
175         tr(ji,jj,jk,jppop,Krhs) = tr(ji,jj,jk,jppop,Krhs) + zmortp * zfactp
176         tr(ji,jj,jk,jpsfe,Krhs) = tr(ji,jj,jk,jpsfe,Krhs) + zmortp * zfactfe
177         prodpoc(ji,jj,jk) = prodpoc(ji,jj,jk) + zmortp
178      END_3D
[7162]179      !
[12377]180       IF(sn_cfctl%l_prttrc)   THEN  ! print mean trends (used for debugging)
[7162]181         WRITE(charout, FMT="('pico')")
[13286]182         CALL prt_ctl_info( charout, cdcomp = 'top' )
[15648]183         CALL prt_ctl(tab4d_1=CASTDP(tr(:,:,:,:,Krhs)), mask1=tmask, clinfo=ctrcnm)
[7162]184       ENDIF
185      !
[9124]186      IF( ln_timing )   CALL timing_stop('p5z_pico')
[7162]187      !
188   END SUBROUTINE p5z_pico
189
190
[12377]191   SUBROUTINE p5z_diat( Kbb, Krhs )
[7162]192      !!---------------------------------------------------------------------
193      !!                     ***  ROUTINE p5z_diat  ***
194      !!
195      !! ** Purpose :   Compute the mortality terms for diatoms
196      !!
197      !! ** Method  : - ???
198      !!---------------------------------------------------------------------
[12377]199      INTEGER, INTENT(in) ::   Kbb, Krhs  ! time level indices
[7162]200      INTEGER  ::  ji, jj, jk
201      REAL(wp) ::  zfactfe,zfactsi,zfactch, zfactn, zfactp, zcompadi
202      REAL(wp) ::  zrespp2, ztortp2, zmortp2
203      REAL(wp) ::  zlim2, zlim1
204      CHARACTER (len=25) :: charout
205      !!---------------------------------------------------------------------
206      !
[9124]207      IF( ln_timing )   CALL timing_start('p5z_diat')
[7162]208      !
209
[13295]210      DO_3D( 1, 1, 1, 1, 1, jpkm1 )
[7162]211
[12377]212         zcompadi = MAX( ( tr(ji,jj,jk,jpdia,Kbb) - 1E-9), 0. )
[7162]213
[12377]214         !   Aggregation term for diatoms is increased in case of nutrient
215         !   stress as observed in reality. The stressed cells become more
216         !   sticky and coagulate to sink quickly out of the euphotic zone
217         !   -------------------------------------------------------------
218         !  Phytoplankton squared mortality
219         !  -------------------------------
220         zlim2   = xlimdia(ji,jj,jk) * xlimdia(ji,jj,jk)
221         zlim1   = 0.25 * ( 1. - zlim2 ) / ( 0.25 + zlim2 ) 
222         zrespp2 = 1.e6 * xstep * (  wchld + wchldm * zlim1 ) * xdiss(ji,jj,jk) * zcompadi * tr(ji,jj,jk,jpdia,Kbb)
[7162]223
[12377]224         !  Phytoplankton linear mortality
225         !  ------------------------------
226         ztortp2 = mpratd * xstep  * zcompadi
227         zmortp2 = zrespp2 + ztortp2
[7162]228
[12377]229         !   Update the arrays tr(:,:,:,:,Krhs) which contains the biological sources and sinks
230         !   ---------------------------------------------------------------------
231         zfactn  = tr(ji,jj,jk,jpndi,Kbb) / ( tr(ji,jj,jk,jpdia,Kbb) + rtrn )
232         zfactp  = tr(ji,jj,jk,jppdi,Kbb) / ( tr(ji,jj,jk,jpdia,Kbb) + rtrn )
233         zfactch = tr(ji,jj,jk,jpdch,Kbb) / ( tr(ji,jj,jk,jpdia,Kbb) + rtrn )
234         zfactfe = tr(ji,jj,jk,jpdfe,Kbb) / ( tr(ji,jj,jk,jpdia,Kbb) + rtrn )
235         zfactsi = tr(ji,jj,jk,jpdsi,Kbb) / ( tr(ji,jj,jk,jpdia,Kbb) + rtrn )
236         tr(ji,jj,jk,jpdia,Krhs) = tr(ji,jj,jk,jpdia,Krhs) - zmortp2 
237         tr(ji,jj,jk,jpndi,Krhs) = tr(ji,jj,jk,jpndi,Krhs) - zmortp2 * zfactn
238         tr(ji,jj,jk,jppdi,Krhs) = tr(ji,jj,jk,jppdi,Krhs) - zmortp2 * zfactp
239         tr(ji,jj,jk,jpdch,Krhs) = tr(ji,jj,jk,jpdch,Krhs) - zmortp2 * zfactch
240         tr(ji,jj,jk,jpdfe,Krhs) = tr(ji,jj,jk,jpdfe,Krhs) - zmortp2 * zfactfe
241         tr(ji,jj,jk,jpdsi,Krhs) = tr(ji,jj,jk,jpdsi,Krhs) - zmortp2 * zfactsi
242         tr(ji,jj,jk,jpgsi,Krhs) = tr(ji,jj,jk,jpgsi,Krhs) + zmortp2 * zfactsi
243         tr(ji,jj,jk,jpgoc,Krhs) = tr(ji,jj,jk,jpgoc,Krhs) + zrespp2 
244         tr(ji,jj,jk,jpgon,Krhs) = tr(ji,jj,jk,jpgon,Krhs) + zrespp2 * zfactn
245         tr(ji,jj,jk,jpgop,Krhs) = tr(ji,jj,jk,jpgop,Krhs) + zrespp2 * zfactp
246         tr(ji,jj,jk,jpbfe,Krhs) = tr(ji,jj,jk,jpbfe,Krhs) + zrespp2 * zfactfe
247         tr(ji,jj,jk,jppoc,Krhs) = tr(ji,jj,jk,jppoc,Krhs) + ztortp2
248         tr(ji,jj,jk,jppon,Krhs) = tr(ji,jj,jk,jppon,Krhs) + ztortp2 * zfactn
249         tr(ji,jj,jk,jppop,Krhs) = tr(ji,jj,jk,jppop,Krhs) + ztortp2 * zfactp
250         tr(ji,jj,jk,jpsfe,Krhs) = tr(ji,jj,jk,jpsfe,Krhs) + ztortp2 * zfactfe
251         prodpoc(ji,jj,jk)   = prodpoc(ji,jj,jk) + ztortp2
252         prodgoc(ji,jj,jk)   = prodgoc(ji,jj,jk) + zrespp2
253      END_3D
[7162]254      !
[12377]255      IF(sn_cfctl%l_prttrc)   THEN  ! print mean trends (used for debugging)
[7162]256         WRITE(charout, FMT="('diat')")
[13286]257         CALL prt_ctl_info( charout, cdcomp = 'top' )
[15648]258         CALL prt_ctl(tab4d_1=CASTDP(tr(:,:,:,:,Krhs)), mask1=tmask, clinfo=ctrcnm)
[7162]259      ENDIF
260      !
[9124]261      IF( ln_timing )   CALL timing_stop('p5z_diat')
[7162]262      !
263   END SUBROUTINE p5z_diat
264
[9124]265
[7162]266   SUBROUTINE p5z_mort_init
267      !!----------------------------------------------------------------------
268      !!                  ***  ROUTINE p5z_mort_init  ***
269      !!
270      !! ** Purpose :   Initialization of phytoplankton parameters
271      !!
272      !! ** Method  :   Read the nampismort namelist and check the parameters
273      !!      called at the first timestep
274      !!
275      !! ** input   :   Namelist nampismort
276      !!
277      !!----------------------------------------------------------------------
[9124]278      INTEGER :: ios                 ! Local integer output status for namelist read
279      !!
[7162]280      NAMELIST/namp5zmort/ wchln, wchlp, wchld, wchldm, mpratn, mpratp, mpratd
[9124]281      !!----------------------------------------------------------------------
[7162]282
283      READ  ( numnatp_ref, namp5zmort, IOSTAT = ios, ERR = 901)
[11536]284901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namp5zmort in reference namelist' )
[7162]285
286      READ  ( numnatp_cfg, namp5zmort, IOSTAT = ios, ERR = 902 )
[11536]287902   IF( ios >  0 ) CALL ctl_nam ( ios , 'namp5zmort in configuration namelist' )
[7162]288      IF(lwm) WRITE ( numonp, namp5zmort )
289
290      IF(lwp) THEN                         ! control print
291         WRITE(numout,*) ' '
292         WRITE(numout,*) ' Namelist parameters for phytoplankton mortality, namp5zmort'
293         WRITE(numout,*) ' ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~'
294         WRITE(numout,*) '    quadratic mortality of phytoplankton      wchln     =', wchln
295         WRITE(numout,*) '    quadratic mortality of picophyto.         wchlp     =', wchlp
296         WRITE(numout,*) '    quadratic mortality of diatoms            wchld     =', wchld
297         WRITE(numout,*) '    Additional quadratic mortality of diatoms wchldm    =', wchldm
298         WRITE(numout,*) '    nanophyto. mortality rate                 mpratn    =', mpratn
299         WRITE(numout,*) '    picophyto. mortality rate                 mpratp    =', mpratp
300         WRITE(numout,*) '    Diatoms mortality rate                    mpratd    =', mpratd
301      ENDIF
302
303   END SUBROUTINE p5z_mort_init
304
305   !!======================================================================
[9788]306END MODULE p5zmort
Note: See TracBrowser for help on using the repository browser.