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