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.
p4zmort.F90 in branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/TOP_SRC/PISCES/P4Z – NEMO

source: branches/2016/dev_r6519_HPC_4/NEMOGCM/NEMO/TOP_SRC/PISCES/P4Z/p4zmort.F90 @ 7508

Last change on this file since 7508 was 7508, checked in by mocavero, 7 years ago

changes on code duplication and workshare construct

  • Property svn:keywords set to Id
File size: 13.1 KB
Line 
1MODULE p4zmort
2   !!======================================================================
3   !!                         ***  MODULE p4zmort  ***
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   !!----------------------------------------------------------------------
9#if defined key_pisces
10   !!----------------------------------------------------------------------
11   !!   'key_pisces'                                       PISCES bio-model
12   !!----------------------------------------------------------------------
13   !!   p4z_mort       :   Compute the mortality terms for phytoplankton
14   !!   p4z_mort_init  :   Initialize the mortality params for phytoplankton
15   !!----------------------------------------------------------------------
16   USE oce_trc         !  shared variables between ocean and passive tracers
17   USE trc             !  passive tracers common variables
18   USE sms_pisces      !  PISCES Source Minus Sink variables
19   USE p4zsink         !  vertical flux of particulate matter due to sinking
20   USE p4zprod         !  Primary productivity
21   USE prtctl_trc      !  print control for debugging
22
23   IMPLICIT NONE
24   PRIVATE
25
26   PUBLIC   p4z_mort   
27   PUBLIC   p4z_mort_init   
28
29   !! * Shared module variables
30   REAL(wp), PUBLIC :: wchl    !:
31   REAL(wp), PUBLIC :: wchld   !:
32   REAL(wp), PUBLIC :: wchldm  !:
33   REAL(wp), PUBLIC :: mprat   !:
34   REAL(wp), PUBLIC :: mprat2  !:
35
36
37   !!----------------------------------------------------------------------
38   !! NEMO/TOP 3.3 , NEMO Consortium (2010)
39   !! $Id$
40   !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt)
41   !!----------------------------------------------------------------------
42
43CONTAINS
44
45   SUBROUTINE p4z_mort( kt )
46      !!---------------------------------------------------------------------
47      !!                     ***  ROUTINE p4z_mort  ***
48      !!
49      !! ** Purpose :   Calls the different subroutine to initialize and compute
50      !!                the different phytoplankton mortality terms
51      !!
52      !! ** Method  : - ???
53      !!---------------------------------------------------------------------
54      INTEGER, INTENT(in) ::   kt ! ocean time step
55      !!---------------------------------------------------------------------
56
57      CALL p4z_nano            ! nanophytoplankton
58
59      CALL p4z_diat            ! diatoms
60
61   END SUBROUTINE p4z_mort
62
63
64   SUBROUTINE p4z_nano
65      !!---------------------------------------------------------------------
66      !!                     ***  ROUTINE p4z_nano  ***
67      !!
68      !! ** Purpose :   Compute the mortality terms for nanophytoplankton
69      !!
70      !! ** Method  : - ???
71      !!---------------------------------------------------------------------
72      INTEGER  :: ji, jj, jk
73      REAL(wp) :: zsizerat, zcompaph
74      REAL(wp) :: zfactfe, zfactch, zprcaca, zfracal
75      REAL(wp) :: ztortp , zrespp , zmortp , zstep
76      CHARACTER (len=25) :: charout
77      !!---------------------------------------------------------------------
78      !
79      IF( nn_timing == 1 )  CALL timing_start('p4z_nano')
80      !
81!$OMP PARALLEL
82!$OMP DO schedule(static) private(jk,jj,ji)
83      DO jk = 1, jpk
84         DO jj = 1, jpj
85            DO ji = 1, jpi
86               prodcal(ji,jj,jk) = 0.  !: calcite production variable set to zero
87            END DO
88         END DO
89      END DO
90!$OMP DO schedule(static) private(jk,jj,ji,zcompaph,zstep,zsizerat,zrespp,ztortp,zmortp,zfactfe,zfactch,zprcaca,zfracal)
91      DO jk = 1, jpkm1
92         DO jj = 1, jpj
93            DO ji = 1, jpi
94               zcompaph = MAX( ( trb(ji,jj,jk,jpphy) - 1e-8 ), 0.e0 )
95               zstep    = xstep
96# if defined key_degrad
97               zstep    = zstep * facvol(ji,jj,jk)
98# endif
99               !     When highly limited by macronutrients, very small cells
100               !     dominate the community. As a consequence, aggregation
101               !     due to turbulence is negligible. Mortality is also set
102               !     to 0
103               zsizerat = MIN(1., MAX( 0., (quotan(ji,jj,jk) - 0.2) / 0.3) ) * trb(ji,jj,jk,jpphy)
104               !     Squared mortality of Phyto similar to a sedimentation term during
105               !     blooms (Doney et al. 1996)
106               zrespp = wchl * 1.e6 * zstep * xdiss(ji,jj,jk) * zcompaph * zsizerat 
107
108               !     Phytoplankton mortality. This mortality loss is slightly
109               !     increased when nutrients are limiting phytoplankton growth
110               !     as observed for instance in case of iron limitation.
111               ztortp = mprat * xstep * zcompaph / ( xkmort + trb(ji,jj,jk,jpphy) ) * zsizerat
112
113               zmortp = zrespp + ztortp
114
115               !   Update the arrays TRA which contains the biological sources and sinks
116
117               zfactfe = trb(ji,jj,jk,jpnfe)/(trb(ji,jj,jk,jpphy)+rtrn)
118               zfactch = trb(ji,jj,jk,jpnch)/(trb(ji,jj,jk,jpphy)+rtrn)
119               tra(ji,jj,jk,jpphy) = tra(ji,jj,jk,jpphy) - zmortp
120               tra(ji,jj,jk,jpnch) = tra(ji,jj,jk,jpnch) - zmortp * zfactch
121               tra(ji,jj,jk,jpnfe) = tra(ji,jj,jk,jpnfe) - zmortp * zfactfe
122               zprcaca = xfracal(ji,jj,jk) * zmortp
123               !
124               prodcal(ji,jj,jk) = prodcal(ji,jj,jk) + zprcaca  ! prodcal=prodcal(nanophy)+prodcal(microzoo)+prodcal(mesozoo)
125               !
126               zfracal = 0.5 * xfracal(ji,jj,jk)
127               tra(ji,jj,jk,jpdic) = tra(ji,jj,jk,jpdic) - zprcaca
128               tra(ji,jj,jk,jptal) = tra(ji,jj,jk,jptal) - 2. * zprcaca
129               tra(ji,jj,jk,jpcal) = tra(ji,jj,jk,jpcal) + zprcaca
130#if defined key_kriest
131               tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + zmortp
132               tra(ji,jj,jk,jpnum) = tra(ji,jj,jk,jpnum) + ztortp * xkr_dnano + zrespp * xkr_ddiat
133               tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + zmortp * zfactfe
134#else
135               tra(ji,jj,jk,jpgoc) = tra(ji,jj,jk,jpgoc) + zfracal * zmortp
136               tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + ( 1. - zfracal ) * zmortp
137               tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + ( 1. - zfracal ) * zmortp * zfactfe
138               tra(ji,jj,jk,jpbfe) = tra(ji,jj,jk,jpbfe) + zfracal * zmortp * zfactfe
139#endif
140            END DO
141         END DO
142      END DO
143!$OMP END DO NOWAIT
144!$OMP END PARALLEL
145      !
146       IF(ln_ctl)   THEN  ! print mean trends (used for debugging)
147         WRITE(charout, FMT="('nano')")
148         CALL prt_ctl_trc_info(charout)
149         CALL prt_ctl_trc(tab4d=tra, mask=tmask, clinfo=ctrcnm)
150       ENDIF
151      !
152      IF( nn_timing == 1 )  CALL timing_stop('p4z_nano')
153      !
154   END SUBROUTINE p4z_nano
155
156   SUBROUTINE p4z_diat
157      !!---------------------------------------------------------------------
158      !!                     ***  ROUTINE p4z_diat  ***
159      !!
160      !! ** Purpose :   Compute the mortality terms for diatoms
161      !!
162      !! ** Method  : - ???
163      !!---------------------------------------------------------------------
164      INTEGER  ::  ji, jj, jk
165      REAL(wp) ::  zfactfe,zfactsi,zfactch, zcompadi
166      REAL(wp) ::  zrespp2, ztortp2, zmortp2, zstep
167      REAL(wp) ::  zlim2, zlim1
168      CHARACTER (len=25) :: charout
169      !!---------------------------------------------------------------------
170      !
171      IF( nn_timing == 1 )  CALL timing_start('p4z_diat')
172      !
173
174      !    Aggregation term for diatoms is increased in case of nutrient
175      !    stress as observed in reality. The stressed cells become more
176      !    sticky and coagulate to sink quickly out of the euphotic zone
177      !     ------------------------------------------------------------
178
179!$OMP PARALLEL DO schedule(static) private(jk,jj,ji,zcompadi,zstep,zlim2,zlim1,zrespp2,ztortp2,zmortp2,zfactfe,zfactch,zfactsi)
180      DO jk = 1, jpkm1
181         DO jj = 1, jpj
182            DO ji = 1, jpi
183
184               zcompadi = MAX( ( trb(ji,jj,jk,jpdia) - 1e-9), 0. )
185
186               !    Aggregation term for diatoms is increased in case of nutrient
187               !    stress as observed in reality. The stressed cells become more
188               !    sticky and coagulate to sink quickly out of the euphotic zone
189               !     ------------------------------------------------------------
190               zstep   = xstep
191# if defined key_degrad
192               zstep = zstep * facvol(ji,jj,jk)
193# endif
194               !  Phytoplankton respiration
195               !     ------------------------
196               zlim2   = xlimdia(ji,jj,jk) * xlimdia(ji,jj,jk)
197               zlim1   = 0.25 * ( 1. - zlim2 ) / ( 0.25 + zlim2 ) 
198               zrespp2 = 1.e6 * zstep * (  wchld + wchldm * zlim1 ) * xdiss(ji,jj,jk) * zcompadi * trb(ji,jj,jk,jpdia)
199
200               !     Phytoplankton mortality.
201               !     ------------------------
202               ztortp2 = mprat2 * zstep * trb(ji,jj,jk,jpdia)  / ( xkmort + trb(ji,jj,jk,jpdia) ) * zcompadi 
203
204               zmortp2 = zrespp2 + ztortp2
205
206               !   Update the arrays tra which contains the biological sources and sinks
207               !   ---------------------------------------------------------------------
208               zfactch = trb(ji,jj,jk,jpdch) / ( trb(ji,jj,jk,jpdia) + rtrn )
209               zfactfe = trb(ji,jj,jk,jpdfe) / ( trb(ji,jj,jk,jpdia) + rtrn )
210               zfactsi = trb(ji,jj,jk,jpdsi) / ( trb(ji,jj,jk,jpdia) + rtrn )
211               tra(ji,jj,jk,jpdia) = tra(ji,jj,jk,jpdia) - zmortp2 
212               tra(ji,jj,jk,jpdch) = tra(ji,jj,jk,jpdch) - zmortp2 * zfactch
213               tra(ji,jj,jk,jpdfe) = tra(ji,jj,jk,jpdfe) - zmortp2 * zfactfe
214               tra(ji,jj,jk,jpdsi) = tra(ji,jj,jk,jpdsi) - zmortp2 * zfactsi
215               tra(ji,jj,jk,jpgsi) = tra(ji,jj,jk,jpgsi) + zmortp2 * zfactsi
216#if defined key_kriest
217               tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + zmortp2 
218               tra(ji,jj,jk,jpnum) = tra(ji,jj,jk,jpnum) + ztortp2 * xkr_ddiat + zrespp2 * xkr_daggr
219               tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + zmortp2 * zfactfe
220#else
221               tra(ji,jj,jk,jpgoc) = tra(ji,jj,jk,jpgoc) + zrespp2 + 0.5 * ztortp2
222               tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + 0.5 * ztortp2
223               tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + 0.5 * ztortp2 * zfactfe
224               tra(ji,jj,jk,jpbfe) = tra(ji,jj,jk,jpbfe) + ( zrespp2 + 0.5 * ztortp2 ) * zfactfe
225#endif
226            END DO
227         END DO
228      END DO
229      !
230      IF(ln_ctl)   THEN  ! print mean trends (used for debugging)
231         WRITE(charout, FMT="('diat')")
232         CALL prt_ctl_trc_info(charout)
233         CALL prt_ctl_trc(tab4d=tra, mask=tmask, clinfo=ctrcnm)
234      ENDIF
235      !
236      IF( nn_timing == 1 )  CALL timing_stop('p4z_diat')
237      !
238   END SUBROUTINE p4z_diat
239
240   SUBROUTINE p4z_mort_init
241
242      !!----------------------------------------------------------------------
243      !!                  ***  ROUTINE p4z_mort_init  ***
244      !!
245      !! ** Purpose :   Initialization of phytoplankton parameters
246      !!
247      !! ** Method  :   Read the nampismort namelist and check the parameters
248      !!      called at the first timestep
249      !!
250      !! ** input   :   Namelist nampismort
251      !!
252      !!----------------------------------------------------------------------
253
254      NAMELIST/nampismort/ wchl, wchld, wchldm, mprat, mprat2
255      INTEGER :: ios                 ! Local integer output status for namelist read
256
257      REWIND( numnatp_ref )              ! Namelist nampismort in reference namelist : Pisces phytoplankton
258      READ  ( numnatp_ref, nampismort, IOSTAT = ios, ERR = 901)
259901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nampismort in reference namelist', lwp )
260
261      REWIND( numnatp_cfg )              ! Namelist nampismort in configuration namelist : Pisces phytoplankton
262      READ  ( numnatp_cfg, nampismort, IOSTAT = ios, ERR = 902 )
263902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'nampismort in configuration namelist', lwp )
264      IF(lwm) WRITE ( numonp, nampismort )
265
266      IF(lwp) THEN                         ! control print
267         WRITE(numout,*) ' '
268         WRITE(numout,*) ' Namelist parameters for phytoplankton mortality, nampismort'
269         WRITE(numout,*) ' ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~'
270         WRITE(numout,*) '    quadratic mortality of phytoplankton      wchl      =', wchl
271         WRITE(numout,*) '    maximum quadratic mortality of diatoms    wchld     =', wchld
272         WRITE(numout,*) '    maximum quadratic mortality of diatoms    wchldm    =', wchldm
273         WRITE(numout,*) '    phytoplankton mortality rate              mprat     =', mprat
274         WRITE(numout,*) '    Diatoms mortality rate                    mprat2    =', mprat2
275      ENDIF
276
277   END SUBROUTINE p4z_mort_init
278
279#else
280   !!======================================================================
281   !!  Dummy module :                                   No PISCES bio-model
282   !!======================================================================
283CONTAINS
284   SUBROUTINE p4z_mort                    ! Empty routine
285   END SUBROUTINE p4z_mort
286#endif 
287
288   !!======================================================================
289END MODULE p4zmort
Note: See TracBrowser for help on using the repository browser.