1 | MODULE 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 | !! p4z_mort : Compute the mortality terms for phytoplankton |
---|
10 | !! p4z_mort_init : Initialize the mortality params for phytoplankton |
---|
11 | !!---------------------------------------------------------------------- |
---|
12 | USE oce_trc ! shared variables between ocean and passive tracers |
---|
13 | USE trc ! passive tracers common variables |
---|
14 | USE sms_pisces ! PISCES Source Minus Sink variables |
---|
15 | USE p4zprod ! Primary productivity |
---|
16 | USE p4zlim ! Phytoplankton limitation terms |
---|
17 | USE prtctl_trc ! print control for debugging |
---|
18 | |
---|
19 | IMPLICIT NONE |
---|
20 | PRIVATE |
---|
21 | |
---|
22 | PUBLIC p4z_mort |
---|
23 | PUBLIC p4z_mort_init |
---|
24 | |
---|
25 | REAL(wp), PUBLIC :: wchl !: |
---|
26 | REAL(wp), PUBLIC :: wchld !: |
---|
27 | REAL(wp), PUBLIC :: wchldm !: |
---|
28 | REAL(wp), PUBLIC :: mprat !: |
---|
29 | REAL(wp), PUBLIC :: mprat2 !: |
---|
30 | |
---|
31 | !! * Substitutions |
---|
32 | # include "do_loop_substitute.h90" |
---|
33 | !!---------------------------------------------------------------------- |
---|
34 | !! NEMO/TOP 4.0 , NEMO Consortium (2018) |
---|
35 | !! $Id$ |
---|
36 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
37 | !!---------------------------------------------------------------------- |
---|
38 | CONTAINS |
---|
39 | |
---|
40 | SUBROUTINE p4z_mort( kt, Kbb, Krhs ) |
---|
41 | !!--------------------------------------------------------------------- |
---|
42 | !! *** ROUTINE p4z_mort *** |
---|
43 | !! |
---|
44 | !! ** Purpose : Calls the different subroutine to initialize and compute |
---|
45 | !! the different phytoplankton mortality terms |
---|
46 | !! |
---|
47 | !! ** Method : - ??? |
---|
48 | !!--------------------------------------------------------------------- |
---|
49 | INTEGER, INTENT(in) :: kt ! ocean time step |
---|
50 | INTEGER, INTENT(in) :: Kbb, Krhs ! time level indices |
---|
51 | !!--------------------------------------------------------------------- |
---|
52 | ! |
---|
53 | CALL p4z_nano( Kbb, Krhs ) ! nanophytoplankton |
---|
54 | ! |
---|
55 | CALL p4z_diat( Kbb, Krhs ) ! diatoms |
---|
56 | ! |
---|
57 | END SUBROUTINE p4z_mort |
---|
58 | |
---|
59 | |
---|
60 | SUBROUTINE p4z_nano( Kbb, Krhs ) |
---|
61 | !!--------------------------------------------------------------------- |
---|
62 | !! *** ROUTINE p4z_nano *** |
---|
63 | !! |
---|
64 | !! ** Purpose : Compute the mortality terms for nanophytoplankton |
---|
65 | !! |
---|
66 | !! ** Method : - ??? |
---|
67 | !!--------------------------------------------------------------------- |
---|
68 | INTEGER, INTENT(in) :: Kbb, Krhs ! time level indices |
---|
69 | INTEGER :: ji, jj, jk |
---|
70 | REAL(wp) :: zsizerat, zcompaph |
---|
71 | REAL(wp) :: zfactfe, zfactch, zprcaca, zfracal |
---|
72 | REAL(wp) :: ztortp , zrespp , zmortp |
---|
73 | CHARACTER (len=25) :: charout |
---|
74 | !!--------------------------------------------------------------------- |
---|
75 | ! |
---|
76 | IF( ln_timing ) CALL timing_start('p4z_nano') |
---|
77 | ! |
---|
78 | prodcal(:,:,:) = 0._wp ! calcite production variable set to zero |
---|
79 | DO_3D_11_11( 1, jpkm1 ) |
---|
80 | zcompaph = MAX( ( tr(ji,jj,jk,jpphy,Kbb) - 1e-8 ), 0.e0 ) |
---|
81 | ! When highly limited by macronutrients, very small cells |
---|
82 | ! dominate the community. As a consequence, aggregation |
---|
83 | ! due to turbulence is negligible. Mortality is also set |
---|
84 | ! to 0 |
---|
85 | zsizerat = MIN(1., MAX( 0., (quotan(ji,jj,jk) - 0.2) / 0.3) ) * tr(ji,jj,jk,jpphy,Kbb) |
---|
86 | ! Squared mortality of Phyto similar to a sedimentation term during |
---|
87 | ! blooms (Doney et al. 1996) |
---|
88 | zrespp = wchl * 1.e6 * xstep * xdiss(ji,jj,jk) * zcompaph * zsizerat |
---|
89 | |
---|
90 | ! Phytoplankton mortality. This mortality loss is slightly |
---|
91 | ! increased when nutrients are limiting phytoplankton growth |
---|
92 | ! as observed for instance in case of iron limitation. |
---|
93 | ztortp = mprat * xstep * zcompaph / ( xkmort + tr(ji,jj,jk,jpphy,Kbb) ) * zsizerat |
---|
94 | |
---|
95 | zmortp = zrespp + ztortp |
---|
96 | |
---|
97 | ! Update the arrays TRA which contains the biological sources and sinks |
---|
98 | |
---|
99 | zfactfe = tr(ji,jj,jk,jpnfe,Kbb)/(tr(ji,jj,jk,jpphy,Kbb)+rtrn) |
---|
100 | zfactch = tr(ji,jj,jk,jpnch,Kbb)/(tr(ji,jj,jk,jpphy,Kbb)+rtrn) |
---|
101 | tr(ji,jj,jk,jpphy,Krhs) = tr(ji,jj,jk,jpphy,Krhs) - zmortp |
---|
102 | tr(ji,jj,jk,jpnch,Krhs) = tr(ji,jj,jk,jpnch,Krhs) - zmortp * zfactch |
---|
103 | tr(ji,jj,jk,jpnfe,Krhs) = tr(ji,jj,jk,jpnfe,Krhs) - zmortp * zfactfe |
---|
104 | zprcaca = xfracal(ji,jj,jk) * zmortp |
---|
105 | ! |
---|
106 | prodcal(ji,jj,jk) = prodcal(ji,jj,jk) + zprcaca ! prodcal=prodcal(nanophy)+prodcal(microzoo)+prodcal(mesozoo) |
---|
107 | ! |
---|
108 | zfracal = 0.5 * xfracal(ji,jj,jk) |
---|
109 | tr(ji,jj,jk,jpdic,Krhs) = tr(ji,jj,jk,jpdic,Krhs) - zprcaca |
---|
110 | tr(ji,jj,jk,jptal,Krhs) = tr(ji,jj,jk,jptal,Krhs) - 2. * zprcaca |
---|
111 | tr(ji,jj,jk,jpcal,Krhs) = tr(ji,jj,jk,jpcal,Krhs) + zprcaca |
---|
112 | tr(ji,jj,jk,jpgoc,Krhs) = tr(ji,jj,jk,jpgoc,Krhs) + zfracal * zmortp |
---|
113 | tr(ji,jj,jk,jppoc,Krhs) = tr(ji,jj,jk,jppoc,Krhs) + ( 1. - zfracal ) * zmortp |
---|
114 | prodpoc(ji,jj,jk) = prodpoc(ji,jj,jk) + ( 1. - zfracal ) * zmortp |
---|
115 | prodgoc(ji,jj,jk) = prodgoc(ji,jj,jk) + zfracal * zmortp |
---|
116 | tr(ji,jj,jk,jpsfe,Krhs) = tr(ji,jj,jk,jpsfe,Krhs) + ( 1. - zfracal ) * zmortp * zfactfe |
---|
117 | tr(ji,jj,jk,jpbfe,Krhs) = tr(ji,jj,jk,jpbfe,Krhs) + zfracal * zmortp * zfactfe |
---|
118 | END_3D |
---|
119 | ! |
---|
120 | IF(sn_cfctl%l_prttrc) THEN ! print mean trends (used for debugging) |
---|
121 | WRITE(charout, FMT="('nano')") |
---|
122 | CALL prt_ctl_trc_info(charout) |
---|
123 | CALL prt_ctl_trc(tab4d=tr(:,:,:,:,Krhs), mask=tmask, clinfo=ctrcnm) |
---|
124 | ENDIF |
---|
125 | ! |
---|
126 | IF( ln_timing ) CALL timing_stop('p4z_nano') |
---|
127 | ! |
---|
128 | END SUBROUTINE p4z_nano |
---|
129 | |
---|
130 | |
---|
131 | SUBROUTINE p4z_diat( Kbb, Krhs ) |
---|
132 | !!--------------------------------------------------------------------- |
---|
133 | !! *** ROUTINE p4z_diat *** |
---|
134 | !! |
---|
135 | !! ** Purpose : Compute the mortality terms for diatoms |
---|
136 | !! |
---|
137 | !! ** Method : - ??? |
---|
138 | !!--------------------------------------------------------------------- |
---|
139 | INTEGER, INTENT(in) :: Kbb, Krhs ! time level indices |
---|
140 | INTEGER :: ji, jj, jk |
---|
141 | REAL(wp) :: zfactfe,zfactsi,zfactch, zcompadi |
---|
142 | REAL(wp) :: zrespp2, ztortp2, zmortp2 |
---|
143 | REAL(wp) :: zlim2, zlim1 |
---|
144 | CHARACTER (len=25) :: charout |
---|
145 | !!--------------------------------------------------------------------- |
---|
146 | ! |
---|
147 | IF( ln_timing ) CALL timing_start('p4z_diat') |
---|
148 | ! |
---|
149 | ! Aggregation term for diatoms is increased in case of nutrient |
---|
150 | ! stress as observed in reality. The stressed cells become more |
---|
151 | ! sticky and coagulate to sink quickly out of the euphotic zone |
---|
152 | ! ------------------------------------------------------------ |
---|
153 | |
---|
154 | DO_3D_11_11( 1, jpkm1 ) |
---|
155 | |
---|
156 | zcompadi = MAX( ( tr(ji,jj,jk,jpdia,Kbb) - 1e-9), 0. ) |
---|
157 | |
---|
158 | ! Aggregation term for diatoms is increased in case of nutrient |
---|
159 | ! stress as observed in reality. The stressed cells become more |
---|
160 | ! sticky and coagulate to sink quickly out of the euphotic zone |
---|
161 | ! ------------------------------------------------------------ |
---|
162 | ! Phytoplankton respiration |
---|
163 | ! ------------------------ |
---|
164 | zlim2 = xlimdia(ji,jj,jk) * xlimdia(ji,jj,jk) |
---|
165 | zlim1 = 0.25 * ( 1. - zlim2 ) / ( 0.25 + zlim2 ) |
---|
166 | zrespp2 = 1.e6 * xstep * ( wchld + wchldm * zlim1 ) * xdiss(ji,jj,jk) * zcompadi * tr(ji,jj,jk,jpdia,Kbb) |
---|
167 | |
---|
168 | ! Phytoplankton mortality. |
---|
169 | ! ------------------------ |
---|
170 | ztortp2 = mprat2 * xstep * tr(ji,jj,jk,jpdia,Kbb) / ( xkmort + tr(ji,jj,jk,jpdia,Kbb) ) * zcompadi |
---|
171 | |
---|
172 | zmortp2 = zrespp2 + ztortp2 |
---|
173 | |
---|
174 | ! Update the arrays tr(:,:,:,:,Krhs) which contains the biological sources and sinks |
---|
175 | ! --------------------------------------------------------------------- |
---|
176 | zfactch = tr(ji,jj,jk,jpdch,Kbb) / ( tr(ji,jj,jk,jpdia,Kbb) + rtrn ) |
---|
177 | zfactfe = tr(ji,jj,jk,jpdfe,Kbb) / ( tr(ji,jj,jk,jpdia,Kbb) + rtrn ) |
---|
178 | zfactsi = tr(ji,jj,jk,jpdsi,Kbb) / ( tr(ji,jj,jk,jpdia,Kbb) + rtrn ) |
---|
179 | tr(ji,jj,jk,jpdia,Krhs) = tr(ji,jj,jk,jpdia,Krhs) - zmortp2 |
---|
180 | tr(ji,jj,jk,jpdch,Krhs) = tr(ji,jj,jk,jpdch,Krhs) - zmortp2 * zfactch |
---|
181 | tr(ji,jj,jk,jpdfe,Krhs) = tr(ji,jj,jk,jpdfe,Krhs) - zmortp2 * zfactfe |
---|
182 | tr(ji,jj,jk,jpdsi,Krhs) = tr(ji,jj,jk,jpdsi,Krhs) - zmortp2 * zfactsi |
---|
183 | tr(ji,jj,jk,jpgsi,Krhs) = tr(ji,jj,jk,jpgsi,Krhs) + zmortp2 * zfactsi |
---|
184 | tr(ji,jj,jk,jpgoc,Krhs) = tr(ji,jj,jk,jpgoc,Krhs) + zrespp2 + 0.5 * ztortp2 |
---|
185 | tr(ji,jj,jk,jppoc,Krhs) = tr(ji,jj,jk,jppoc,Krhs) + 0.5 * ztortp2 |
---|
186 | prodpoc(ji,jj,jk) = prodpoc(ji,jj,jk) + 0.5 * ztortp2 |
---|
187 | prodgoc(ji,jj,jk) = prodgoc(ji,jj,jk) + zrespp2 + 0.5 * ztortp2 |
---|
188 | tr(ji,jj,jk,jpsfe,Krhs) = tr(ji,jj,jk,jpsfe,Krhs) + 0.5 * ztortp2 * zfactfe |
---|
189 | tr(ji,jj,jk,jpbfe,Krhs) = tr(ji,jj,jk,jpbfe,Krhs) + ( zrespp2 + 0.5 * ztortp2 ) * zfactfe |
---|
190 | END_3D |
---|
191 | ! |
---|
192 | IF(sn_cfctl%l_prttrc) THEN ! print mean trends (used for debugging) |
---|
193 | WRITE(charout, FMT="('diat')") |
---|
194 | CALL prt_ctl_trc_info(charout) |
---|
195 | CALL prt_ctl_trc(tab4d=tr(:,:,:,:,Krhs), mask=tmask, clinfo=ctrcnm) |
---|
196 | ENDIF |
---|
197 | ! |
---|
198 | IF( ln_timing ) CALL timing_stop('p4z_diat') |
---|
199 | ! |
---|
200 | END SUBROUTINE p4z_diat |
---|
201 | |
---|
202 | |
---|
203 | SUBROUTINE p4z_mort_init |
---|
204 | !!---------------------------------------------------------------------- |
---|
205 | !! *** ROUTINE p4z_mort_init *** |
---|
206 | !! |
---|
207 | !! ** Purpose : Initialization of phytoplankton parameters |
---|
208 | !! |
---|
209 | !! ** Method : Read the nampismort namelist and check the parameters |
---|
210 | !! called at the first timestep |
---|
211 | !! |
---|
212 | !! ** input : Namelist nampismort |
---|
213 | !! |
---|
214 | !!---------------------------------------------------------------------- |
---|
215 | INTEGER :: ios ! Local integer |
---|
216 | ! |
---|
217 | NAMELIST/namp4zmort/ wchl, wchld, wchldm, mprat, mprat2 |
---|
218 | !!---------------------------------------------------------------------- |
---|
219 | ! |
---|
220 | IF(lwp) THEN |
---|
221 | WRITE(numout,*) |
---|
222 | WRITE(numout,*) 'p4z_mort_init : Initialization of phytoplankton mortality parameters' |
---|
223 | WRITE(numout,*) '~~~~~~~~~~~~~' |
---|
224 | ENDIF |
---|
225 | ! |
---|
226 | READ ( numnatp_ref, namp4zmort, IOSTAT = ios, ERR = 901) |
---|
227 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namp4zmort in reference namelist' ) |
---|
228 | READ ( numnatp_cfg, namp4zmort, IOSTAT = ios, ERR = 902 ) |
---|
229 | 902 IF( ios > 0 ) CALL ctl_nam ( ios , 'namp4zmort in configuration namelist' ) |
---|
230 | IF(lwm) WRITE( numonp, namp4zmort ) |
---|
231 | ! |
---|
232 | IF(lwp) THEN ! control print |
---|
233 | WRITE(numout,*) ' Namelist : namp4zmort' |
---|
234 | WRITE(numout,*) ' quadratic mortality of phytoplankton wchl =', wchl |
---|
235 | WRITE(numout,*) ' maximum quadratic mortality of diatoms wchld =', wchld |
---|
236 | WRITE(numout,*) ' maximum quadratic mortality of diatoms wchldm =', wchldm |
---|
237 | WRITE(numout,*) ' phytoplankton mortality rate mprat =', mprat |
---|
238 | WRITE(numout,*) ' Diatoms mortality rate mprat2 =', mprat2 |
---|
239 | ENDIF |
---|
240 | ! |
---|
241 | END SUBROUTINE p4z_mort_init |
---|
242 | |
---|
243 | !!====================================================================== |
---|
244 | END MODULE p4zmort |
---|