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 | #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 trc |
---|
17 | USE oce_trc ! |
---|
18 | USE trp_trc ! |
---|
19 | USE sms_pisces ! |
---|
20 | USE p4zsink |
---|
21 | USE prtctl_trc |
---|
22 | |
---|
23 | IMPLICIT NONE |
---|
24 | PRIVATE |
---|
25 | |
---|
26 | PUBLIC p4z_mort |
---|
27 | |
---|
28 | |
---|
29 | !! * Shared module variables |
---|
30 | REAL(wp), PUBLIC :: & |
---|
31 | wchl = 0.001_wp , & !: |
---|
32 | wchld = 0.02_wp , & !: |
---|
33 | mprat = 0.01_wp , & !: |
---|
34 | mprat2 = 0.01_wp , & !: |
---|
35 | mpratm = 0.01_wp !: |
---|
36 | |
---|
37 | !! * Module variables |
---|
38 | REAL(wp) :: zstep |
---|
39 | |
---|
40 | |
---|
41 | |
---|
42 | !!* Substitution |
---|
43 | # include "domzgr_substitute.h90" |
---|
44 | !!---------------------------------------------------------------------- |
---|
45 | !! NEMO/TOP 2.0 , LOCEAN-IPSL (2007) |
---|
46 | !! $Header:$ |
---|
47 | !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt) |
---|
48 | !!---------------------------------------------------------------------- |
---|
49 | |
---|
50 | CONTAINS |
---|
51 | |
---|
52 | SUBROUTINE p4z_mort( kt, jnt ) |
---|
53 | !!--------------------------------------------------------------------- |
---|
54 | !! *** ROUTINE p4z_mort *** |
---|
55 | !! |
---|
56 | !! ** Purpose : Calls the different subroutine to initialize and compute |
---|
57 | !! the different phytoplankton mortality terms |
---|
58 | !! |
---|
59 | !! ** Method : - ??? |
---|
60 | !!--------------------------------------------------------------------- |
---|
61 | INTEGER, INTENT(in) :: kt, jnt ! ocean time step |
---|
62 | !!--------------------------------------------------------------------- |
---|
63 | |
---|
64 | IF( ( kt * jnt ) == nittrc000 ) CALL p4z_mort_init ! Initialization (first time-step only) |
---|
65 | |
---|
66 | zstep = rfact2 / rjjss ! Time step duration for biology |
---|
67 | |
---|
68 | CALL p4z_nano ! nanophytoplankton |
---|
69 | |
---|
70 | CALL p4z_diat ! diatoms |
---|
71 | |
---|
72 | END SUBROUTINE p4z_mort |
---|
73 | |
---|
74 | |
---|
75 | SUBROUTINE p4z_nano |
---|
76 | !!--------------------------------------------------------------------- |
---|
77 | !! *** ROUTINE p4z_nano *** |
---|
78 | !! |
---|
79 | !! ** Purpose : Compute the mortality terms for nanophytoplankton |
---|
80 | !! |
---|
81 | !! ** Method : - ??? |
---|
82 | !!--------------------------------------------------------------------- |
---|
83 | INTEGER :: ji, jj, jk |
---|
84 | REAL(wp) :: zcompaph |
---|
85 | REAL(wp) :: zfactfe,zfactch,zprcaca,zfracal |
---|
86 | REAL(wp) :: ztortp,zrespp,zmortp |
---|
87 | CHARACTER (len=25) :: charout |
---|
88 | !!--------------------------------------------------------------------- |
---|
89 | |
---|
90 | |
---|
91 | DO jk = 1, jpkm1 |
---|
92 | DO jj = 1, jpj |
---|
93 | DO ji = 1, jpi |
---|
94 | |
---|
95 | zcompaph = MAX( ( trn(ji,jj,jk,jpphy) - 1e-8 ), 0.e0 ) |
---|
96 | |
---|
97 | ! Squared mortality of Phyto similar to a sedimentation term during |
---|
98 | ! blooms (Doney et al. 1996) |
---|
99 | ! ----------------------------------------------------------------- |
---|
100 | zrespp = wchl * 1.e6 * zstep * xdiss(ji,jj,jk) & |
---|
101 | # if defined key_off_degrad |
---|
102 | & * facvol(ji,jj,jk) & |
---|
103 | # endif |
---|
104 | & * zcompaph * trn(ji,jj,jk,jpphy) |
---|
105 | |
---|
106 | ! Phytoplankton mortality. This mortality loss is slightly |
---|
107 | ! increased when nutrients are limiting phytoplankton growth |
---|
108 | ! as observed for instance in case of iron limitation. |
---|
109 | ! ---------------------------------------------------------- |
---|
110 | ztortp = mprat * zstep * trn(ji,jj,jk,jpphy) & |
---|
111 | # if defined key_off_degrad |
---|
112 | & * facvol(ji,jj,jk) & |
---|
113 | # endif |
---|
114 | & / ( xkmort + trn(ji,jj,jk,jpphy) ) * zcompaph |
---|
115 | |
---|
116 | |
---|
117 | zmortp = zrespp + ztortp |
---|
118 | |
---|
119 | ! Update the arrays TRA which contains the biological sources and sinks |
---|
120 | |
---|
121 | zfactfe = trn(ji,jj,jk,jpnfe)/(trn(ji,jj,jk,jpphy)+rtrn) |
---|
122 | zfactch = trn(ji,jj,jk,jpnch)/(trn(ji,jj,jk,jpphy)+rtrn) |
---|
123 | |
---|
124 | tra(ji,jj,jk,jpphy) = tra(ji,jj,jk,jpphy) - zmortp |
---|
125 | tra(ji,jj,jk,jpnch) = tra(ji,jj,jk,jpnch) - zmortp * zfactch |
---|
126 | tra(ji,jj,jk,jpnfe) = tra(ji,jj,jk,jpnfe) - zmortp * zfactfe |
---|
127 | zprcaca = xfracal(ji,jj,jk) * zmortp |
---|
128 | zfracal = 0.5 * xfracal(ji,jj,jk) |
---|
129 | tra(ji,jj,jk,jpdic) = tra(ji,jj,jk,jpdic) - zprcaca |
---|
130 | tra(ji,jj,jk,jptal) = tra(ji,jj,jk,jptal) - 2. * zprcaca |
---|
131 | tra(ji,jj,jk,jpcal) = tra(ji,jj,jk,jpcal) + zprcaca |
---|
132 | #if defined key_kriest |
---|
133 | tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + zmortp |
---|
134 | tra(ji,jj,jk,jpnum) = tra(ji,jj,jk,jpnum) + ztortp * xkr_dnano + zrespp * xkr_ddiat |
---|
135 | tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + zmortp * zfactfe |
---|
136 | #else |
---|
137 | tra(ji,jj,jk,jpgoc) = tra(ji,jj,jk,jpgoc) + zfracal * zmortp |
---|
138 | tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + ( 1. - zfracal ) * zmortp |
---|
139 | tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + ( 1. - zfracal ) * zmortp * zfactfe |
---|
140 | tra(ji,jj,jk,jpbfe) = tra(ji,jj,jk,jpbfe) + zfracal * zmortp * zfactfe |
---|
141 | #endif |
---|
142 | END DO |
---|
143 | END DO |
---|
144 | END DO |
---|
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 | END SUBROUTINE p4z_nano |
---|
153 | |
---|
154 | SUBROUTINE p4z_diat |
---|
155 | !!--------------------------------------------------------------------- |
---|
156 | !! *** ROUTINE p4z_diat *** |
---|
157 | !! |
---|
158 | !! ** Purpose : Compute the mortality terms for diatoms |
---|
159 | !! |
---|
160 | !! ** Method : - ??? |
---|
161 | !!--------------------------------------------------------------------- |
---|
162 | INTEGER :: ji, jj, jk |
---|
163 | REAL(wp) :: zfactfe,zfactsi,zfactch, zcompadi |
---|
164 | REAL(wp) :: zrespp2, ztortp2, zmortp2 |
---|
165 | CHARACTER (len=25) :: charout |
---|
166 | |
---|
167 | !!--------------------------------------------------------------------- |
---|
168 | |
---|
169 | |
---|
170 | ! Aggregation term for diatoms is increased in case of nutrient |
---|
171 | ! stress as observed in reality. The stressed cells become more |
---|
172 | ! sticky and coagulate to sink quickly out of the euphotic zone |
---|
173 | ! ------------------------------------------------------------ |
---|
174 | |
---|
175 | DO jk = 1, jpkm1 |
---|
176 | DO jj = 1, jpj |
---|
177 | DO ji = 1, jpi |
---|
178 | |
---|
179 | zcompadi = MAX( ( trn(ji,jj,jk,jpdia) - 1e-8), 0. ) |
---|
180 | |
---|
181 | ! Aggregation term for diatoms is increased in case of nutrient |
---|
182 | ! stress as observed in reality. The stressed cells become more |
---|
183 | ! sticky and coagulate to sink quickly out of the euphotic zone |
---|
184 | ! ------------------------------------------------------------ |
---|
185 | |
---|
186 | zrespp2 = 1.e6 * zstep * ( wchl + wchld * ( 1.- xlimdia(ji,jj,jk) ) ) & |
---|
187 | # if defined key_off_degrad |
---|
188 | & * facvol(ji,jj,jk) & |
---|
189 | # endif |
---|
190 | & * xdiss(ji,jj,jk) * zcompadi * trn(ji,jj,jk,jpdia) |
---|
191 | |
---|
192 | |
---|
193 | ! Phytoplankton mortality. |
---|
194 | ! ------------------------ |
---|
195 | ztortp2 = mprat2 * zstep * trn(ji,jj,jk,jpdia) & |
---|
196 | # if defined key_off_degrad |
---|
197 | & * facvol(ji,jj,jk) & |
---|
198 | # endif |
---|
199 | & / ( xkmort + trn(ji,jj,jk,jpdia) ) * zcompadi |
---|
200 | |
---|
201 | zmortp2 = zrespp2 + ztortp2 |
---|
202 | |
---|
203 | ! Update the arrays tra which contains the biological sources and sinks |
---|
204 | ! --------------------------------------------------------------------- |
---|
205 | zfactch = trn(ji,jj,jk,jpdch) / ( trn(ji,jj,jk,jpdia) + rtrn ) |
---|
206 | zfactfe = trn(ji,jj,jk,jpdfe) / ( trn(ji,jj,jk,jpdia) + rtrn ) |
---|
207 | zfactsi = trn(ji,jj,jk,jpbsi) / ( trn(ji,jj,jk,jpdia) + rtrn ) |
---|
208 | |
---|
209 | tra(ji,jj,jk,jpdia) = tra(ji,jj,jk,jpdia) - zmortp2 |
---|
210 | tra(ji,jj,jk,jpdch) = tra(ji,jj,jk,jpdch) - zmortp2 * zfactch |
---|
211 | tra(ji,jj,jk,jpdfe) = tra(ji,jj,jk,jpdfe) - zmortp2 * zfactfe |
---|
212 | tra(ji,jj,jk,jpbsi) = tra(ji,jj,jk,jpbsi) - zmortp2 * zfactsi |
---|
213 | tra(ji,jj,jk,jpdsi) = tra(ji,jj,jk,jpdsi) + zmortp2 * zfactsi |
---|
214 | #if defined key_kriest |
---|
215 | tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + zmortp2 |
---|
216 | tra(ji,jj,jk,jpnum) = tra(ji,jj,jk,jpnum) + ztortp2 * xkr_ddiat + zrespp2 * xkr_daggr |
---|
217 | tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + zmortp2 * zfactfe |
---|
218 | #else |
---|
219 | tra(ji,jj,jk,jpgoc) = tra(ji,jj,jk,jpgoc) + zrespp2 + 0.5 * ztortp2 |
---|
220 | tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + 0.5 * ztortp2 |
---|
221 | tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + 0.5 * ztortp2 * zfactfe |
---|
222 | tra(ji,jj,jk,jpbfe) = tra(ji,jj,jk,jpbfe) + ( zrespp2 + 0.5 * ztortp2 ) * zfactfe |
---|
223 | #endif |
---|
224 | END DO |
---|
225 | END DO |
---|
226 | END DO |
---|
227 | ! |
---|
228 | IF(ln_ctl) THEN ! print mean trends (used for debugging) |
---|
229 | WRITE(charout, FMT="('diat')") |
---|
230 | CALL prt_ctl_trc_info(charout) |
---|
231 | CALL prt_ctl_trc(tab4d=tra, mask=tmask, clinfo=ctrcnm) |
---|
232 | ENDIF |
---|
233 | |
---|
234 | END SUBROUTINE p4z_diat |
---|
235 | |
---|
236 | SUBROUTINE p4z_mort_init |
---|
237 | |
---|
238 | !!---------------------------------------------------------------------- |
---|
239 | !! *** ROUTINE p4z_mort_init *** |
---|
240 | !! |
---|
241 | !! ** Purpose : Initialization of phytoplankton parameters |
---|
242 | !! |
---|
243 | !! ** Method : Read the natmort namelist and check the parameters |
---|
244 | !! called at the first timestep (nittrc000) |
---|
245 | !! |
---|
246 | !! ** input : Namelist natmort |
---|
247 | !! |
---|
248 | !!---------------------------------------------------------------------- |
---|
249 | |
---|
250 | NAMELIST/natmort/ wchl, wchld, mprat, mprat2, mpratm |
---|
251 | |
---|
252 | REWIND( numnat ) ! read numnat |
---|
253 | READ ( numnat, natmort ) |
---|
254 | |
---|
255 | IF(lwp) THEN ! control print |
---|
256 | WRITE(numout,*) ' ' |
---|
257 | WRITE(numout,*) ' Namelist parameters for phytoplankton mortality, natmort' |
---|
258 | WRITE(numout,*) ' ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~' |
---|
259 | WRITE(numout,*) ' quadratic mortality of phytoplankton wchl =', wchl |
---|
260 | WRITE(numout,*) ' maximum quadratic mortality of diatoms wchld =', wchld |
---|
261 | WRITE(numout,*) ' phytoplankton mortality rate mprat =', mprat |
---|
262 | WRITE(numout,*) ' Diatoms mortality rate mprat2 =', mprat2 |
---|
263 | WRITE(numout,*) ' Phytoplankton minimum mortality rate mpratm =', mpratm |
---|
264 | ENDIF |
---|
265 | |
---|
266 | END SUBROUTINE p4z_mort_init |
---|
267 | |
---|
268 | #else |
---|
269 | !!====================================================================== |
---|
270 | !! Dummy module : No PISCES bio-model |
---|
271 | !!====================================================================== |
---|
272 | CONTAINS |
---|
273 | SUBROUTINE p4z_mort ! Empty routine |
---|
274 | END SUBROUTINE p4z_mort |
---|
275 | #endif |
---|
276 | |
---|
277 | !!====================================================================== |
---|
278 | END MODULE p4zmort |
---|