[6453] | 1 | MODULE 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 | #if defined key_pisces_quota |
---|
| 11 | !!---------------------------------------------------------------------- |
---|
| 12 | !! 'key_pisces_quota' PISCES bio-model |
---|
| 13 | !!---------------------------------------------------------------------- |
---|
| 14 | !! p5z_mort : Compute the mortality terms for phytoplankton |
---|
| 15 | !! p5z_mort_init : Initialize the mortality params for phytoplankton |
---|
| 16 | !!---------------------------------------------------------------------- |
---|
| 17 | USE oce_trc ! shared variables between ocean and passive tracers |
---|
| 18 | USE trc ! passive tracers common variables |
---|
| 19 | USE sms_pisces ! PISCES Source Minus Sink variables |
---|
| 20 | USE p5zlim ! Phytoplankton limitation terms |
---|
| 21 | USE prtctl_trc ! print control for debugging |
---|
| 22 | |
---|
| 23 | IMPLICIT NONE |
---|
| 24 | PRIVATE |
---|
| 25 | |
---|
| 26 | PUBLIC p5z_mort |
---|
| 27 | PUBLIC p5z_mort_init |
---|
| 28 | |
---|
| 29 | !! * Shared module variables |
---|
| 30 | REAL(wp), PUBLIC :: wchl !: |
---|
| 31 | REAL(wp), PUBLIC :: wchlp !: |
---|
| 32 | REAL(wp), PUBLIC :: wchld !: |
---|
| 33 | REAL(wp), PUBLIC :: wchldm !: |
---|
| 34 | REAL(wp), PUBLIC :: mprat !: |
---|
| 35 | REAL(wp), PUBLIC :: mpratp !: |
---|
| 36 | REAL(wp), PUBLIC :: mprat2 !: |
---|
| 37 | |
---|
| 38 | |
---|
| 39 | !!* Substitution |
---|
| 40 | # include "top_substitute.h90" |
---|
| 41 | !!---------------------------------------------------------------------- |
---|
| 42 | !! NEMO/TOP 3.3 , NEMO Consortium (2010) |
---|
| 43 | !! $Id: p4zmort.F90 3160 2011-11-20 14:27:18Z cetlod $ |
---|
| 44 | !! Software governed by the CeCILL licence (NEMOGCM/NEMO_CeCILL.txt) |
---|
| 45 | !!---------------------------------------------------------------------- |
---|
| 46 | |
---|
| 47 | CONTAINS |
---|
| 48 | |
---|
| 49 | SUBROUTINE p5z_mort( kt ) |
---|
| 50 | !!--------------------------------------------------------------------- |
---|
| 51 | !! *** ROUTINE p5z_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 p5z_nano ! nanophytoplankton |
---|
| 62 | CALL p5z_pico ! picophytoplankton |
---|
| 63 | CALL p5z_diat ! diatoms |
---|
| 64 | |
---|
| 65 | END SUBROUTINE p5z_mort |
---|
| 66 | |
---|
| 67 | |
---|
| 68 | SUBROUTINE p5z_nano |
---|
| 69 | !!--------------------------------------------------------------------- |
---|
| 70 | !! *** ROUTINE p5z_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, zfactn, zfactp, zprcaca |
---|
| 79 | REAL(wp) :: ztortp , zrespp , zmortp , zstep |
---|
| 80 | CHARACTER (len=25) :: charout |
---|
| 81 | !!--------------------------------------------------------------------- |
---|
| 82 | ! |
---|
| 83 | IF( nn_timing == 1 ) CALL timing_start('p5z_nano') |
---|
| 84 | ! |
---|
| 85 | prodcal(:,:,:) = 0. !: calcite production variable set to zero |
---|
| 86 | DO jk = 1, jpkm1 |
---|
| 87 | DO jj = 1, jpj |
---|
| 88 | DO ji = 1, jpi |
---|
| 89 | zcompaph = MAX( ( trb(ji,jj,jk,jpphy) - 1e-9 ), 0.e0 ) |
---|
| 90 | zstep = xstep |
---|
| 91 | # if defined key_degrad |
---|
| 92 | zstep = zstep * facvol(ji,jj,jk) |
---|
| 93 | # endif |
---|
| 94 | ! Squared mortality of Phyto similar to a sedimentation term during |
---|
| 95 | ! blooms (Doney et al. 1996) |
---|
| 96 | ! ----------------------------------------------------------------- |
---|
| 97 | zrespp = wchl * 1.e6 * zstep * xdiss(ji,jj,jk) * zcompaph * trb(ji,jj,jk,jpphy) |
---|
| 98 | |
---|
| 99 | ! Phytoplankton linear mortality |
---|
| 100 | ! ------------------------------ |
---|
| 101 | ztortp = mprat * xstep * zcompaph |
---|
| 102 | zmortp = zrespp + ztortp |
---|
| 103 | |
---|
| 104 | ! Update the arrays TRA which contains the biological sources and sinks |
---|
| 105 | |
---|
| 106 | zfactn = trb(ji,jj,jk,jpnph)/(trb(ji,jj,jk,jpphy)+rtrn) |
---|
| 107 | zfactp = trb(ji,jj,jk,jppph)/(trb(ji,jj,jk,jpphy)+rtrn) |
---|
| 108 | zfactfe = trb(ji,jj,jk,jpnfe)/(trb(ji,jj,jk,jpphy)+rtrn) |
---|
| 109 | zfactch = trb(ji,jj,jk,jpnch)/(trb(ji,jj,jk,jpphy)+rtrn) |
---|
| 110 | tra(ji,jj,jk,jpphy) = tra(ji,jj,jk,jpphy) - zmortp |
---|
| 111 | tra(ji,jj,jk,jpnph) = tra(ji,jj,jk,jpnph) - zmortp * zfactn |
---|
| 112 | tra(ji,jj,jk,jppph) = tra(ji,jj,jk,jppph) - zmortp * zfactp |
---|
| 113 | tra(ji,jj,jk,jpnch) = tra(ji,jj,jk,jpnch) - zmortp * zfactch |
---|
| 114 | tra(ji,jj,jk,jpnfe) = tra(ji,jj,jk,jpnfe) - zmortp * zfactfe |
---|
| 115 | zprcaca = xfracal(ji,jj,jk) * zmortp |
---|
| 116 | ! |
---|
| 117 | prodcal(ji,jj,jk) = prodcal(ji,jj,jk) + zprcaca ! prodcal=prodcal(nanophy)+prodcal(microzoo)+prodcal(mesozoo) |
---|
| 118 | ! |
---|
| 119 | tra(ji,jj,jk,jpdic) = tra(ji,jj,jk,jpdic) - zprcaca |
---|
| 120 | tra(ji,jj,jk,jptal) = tra(ji,jj,jk,jptal) - 2. * zprcaca |
---|
| 121 | tra(ji,jj,jk,jpcal) = tra(ji,jj,jk,jpcal) + zprcaca |
---|
| 122 | #if defined key_kriest |
---|
| 123 | tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + zmortp |
---|
| 124 | tra(ji,jj,jk,jppon) = tra(ji,jj,jk,jppon) + zmortp * zfactn |
---|
| 125 | tra(ji,jj,jk,jppop) = tra(ji,jj,jk,jppop) + zmortp * zfactp |
---|
| 126 | tra(ji,jj,jk,jpnum) = tra(ji,jj,jk,jpnum) + ztortp * xkr_dnano + zrespp * xkr_ddiat |
---|
| 127 | tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + zmortp * zfactfe |
---|
| 128 | #else |
---|
| 129 | tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + zmortp |
---|
| 130 | tra(ji,jj,jk,jppon) = tra(ji,jj,jk,jppon) + zmortp * zfactn |
---|
| 131 | tra(ji,jj,jk,jppop) = tra(ji,jj,jk,jppop) + zmortp * zfactp |
---|
| 132 | prodpoc(ji,jj,jk) = prodpoc(ji,jj,jk) + zmortp |
---|
| 133 | tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + zmortp * zfactfe |
---|
| 134 | #endif |
---|
| 135 | END DO |
---|
| 136 | END DO |
---|
| 137 | END DO |
---|
| 138 | ! |
---|
| 139 | IF(ln_ctl) THEN ! print mean trends (used for debugging) |
---|
| 140 | WRITE(charout, FMT="('nano')") |
---|
| 141 | CALL prt_ctl_trc_info(charout) |
---|
| 142 | CALL prt_ctl_trc(tab4d=tra, mask=tmask, clinfo=ctrcnm) |
---|
| 143 | ENDIF |
---|
| 144 | ! |
---|
| 145 | IF( nn_timing == 1 ) CALL timing_stop('p5z_nano') |
---|
| 146 | ! |
---|
| 147 | END SUBROUTINE p5z_nano |
---|
| 148 | |
---|
| 149 | SUBROUTINE p5z_pico |
---|
| 150 | !!--------------------------------------------------------------------- |
---|
| 151 | !! *** ROUTINE p5z_pico *** |
---|
| 152 | !! |
---|
| 153 | !! ** Purpose : Compute the mortality terms for picophytoplankton |
---|
| 154 | !! |
---|
| 155 | !! ** Method : - ??? |
---|
| 156 | !!--------------------------------------------------------------------- |
---|
| 157 | INTEGER :: ji, jj, jk |
---|
| 158 | REAL(wp) :: zcompaph |
---|
| 159 | REAL(wp) :: zfactfe, zfactch, zfactn, zfactp |
---|
| 160 | REAL(wp) :: ztortp , zrespp , zmortp , zstep |
---|
| 161 | CHARACTER (len=25) :: charout |
---|
| 162 | !!--------------------------------------------------------------------- |
---|
| 163 | ! |
---|
| 164 | IF( nn_timing == 1 ) CALL timing_start('p5z_pico') |
---|
| 165 | ! |
---|
| 166 | DO jk = 1, jpkm1 |
---|
| 167 | DO jj = 1, jpj |
---|
| 168 | DO ji = 1, jpi |
---|
| 169 | zcompaph = MAX( ( trb(ji,jj,jk,jppic) - 1e-9 ), 0.e0 ) |
---|
| 170 | zstep = xstep |
---|
| 171 | # if defined key_degrad |
---|
| 172 | zstep = zstep * facvol(ji,jj,jk) |
---|
| 173 | # endif |
---|
| 174 | ! Squared mortality of Phyto similar to a sedimentation term during |
---|
| 175 | ! blooms (Doney et al. 1996) |
---|
| 176 | ! ----------------------------------------------------------------- |
---|
| 177 | zrespp = wchlp * 1.e6 * zstep * xdiss(ji,jj,jk) * zcompaph * trb(ji,jj,jk,jppic) |
---|
| 178 | |
---|
| 179 | ! Phytoplankton mortality |
---|
| 180 | ztortp = mpratp * xstep * zcompaph |
---|
| 181 | zmortp = zrespp + ztortp |
---|
| 182 | |
---|
| 183 | ! Update the arrays TRA which contains the biological sources and sinks |
---|
| 184 | |
---|
| 185 | zfactn = trb(ji,jj,jk,jpnpi)/(trb(ji,jj,jk,jppic)+rtrn) |
---|
| 186 | zfactp = trb(ji,jj,jk,jpppi)/(trb(ji,jj,jk,jppic)+rtrn) |
---|
| 187 | zfactfe = trb(ji,jj,jk,jppfe)/(trb(ji,jj,jk,jppic)+rtrn) |
---|
| 188 | zfactch = trb(ji,jj,jk,jppch)/(trb(ji,jj,jk,jppic)+rtrn) |
---|
| 189 | tra(ji,jj,jk,jppic) = tra(ji,jj,jk,jppic) - zmortp |
---|
| 190 | tra(ji,jj,jk,jpnpi) = tra(ji,jj,jk,jpnpi) - zmortp * zfactn |
---|
| 191 | tra(ji,jj,jk,jpppi) = tra(ji,jj,jk,jpppi) - zmortp * zfactp |
---|
| 192 | tra(ji,jj,jk,jppch) = tra(ji,jj,jk,jppch) - zmortp * zfactch |
---|
| 193 | tra(ji,jj,jk,jppfe) = tra(ji,jj,jk,jppfe) - zmortp * zfactfe |
---|
| 194 | #if defined key_kriest |
---|
| 195 | tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + zmortp |
---|
| 196 | tra(ji,jj,jk,jppon) = tra(ji,jj,jk,jppon) + zmortp * zfactn |
---|
| 197 | tra(ji,jj,jk,jppop) = tra(ji,jj,jk,jppop) + zmortp * zfactp |
---|
| 198 | tra(ji,jj,jk,jpnum) = tra(ji,jj,jk,jpnum) + ztortp * xkr_dnano + zrespp * xkr_ddiat |
---|
| 199 | tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + zmortp * zfactfe |
---|
| 200 | #else |
---|
| 201 | tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + zmortp |
---|
| 202 | tra(ji,jj,jk,jppon) = tra(ji,jj,jk,jppon) + zmortp * zfactn |
---|
| 203 | tra(ji,jj,jk,jppop) = tra(ji,jj,jk,jppop) + zmortp * zfactp |
---|
| 204 | tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + zmortp * zfactfe |
---|
| 205 | prodpoc(ji,jj,jk) = prodpoc(ji,jj,jk) + zmortp |
---|
| 206 | prodgoc(ji,jj,jk) = prodgoc(ji,jj,jk) + zmortp |
---|
| 207 | #endif |
---|
| 208 | END DO |
---|
| 209 | END DO |
---|
| 210 | END DO |
---|
| 211 | ! |
---|
| 212 | IF(ln_ctl) THEN ! print mean trends (used for debugging) |
---|
| 213 | WRITE(charout, FMT="('pico')") |
---|
| 214 | CALL prt_ctl_trc_info(charout) |
---|
| 215 | CALL prt_ctl_trc(tab4d=tra, mask=tmask, clinfo=ctrcnm) |
---|
| 216 | ENDIF |
---|
| 217 | ! |
---|
| 218 | IF( nn_timing == 1 ) CALL timing_stop('p5z_pico') |
---|
| 219 | ! |
---|
| 220 | END SUBROUTINE p5z_pico |
---|
| 221 | |
---|
| 222 | |
---|
| 223 | SUBROUTINE p5z_diat |
---|
| 224 | !!--------------------------------------------------------------------- |
---|
| 225 | !! *** ROUTINE p5z_diat *** |
---|
| 226 | !! |
---|
| 227 | !! ** Purpose : Compute the mortality terms for diatoms |
---|
| 228 | !! |
---|
| 229 | !! ** Method : - ??? |
---|
| 230 | !!--------------------------------------------------------------------- |
---|
| 231 | INTEGER :: ji, jj, jk |
---|
| 232 | REAL(wp) :: zfactfe,zfactsi,zfactch, zfactn, zfactp, zcompadi |
---|
| 233 | REAL(wp) :: zrespp2, ztortp2, zmortp2, zstep |
---|
| 234 | REAL(wp) :: zlim2, zlim1 |
---|
| 235 | CHARACTER (len=25) :: charout |
---|
| 236 | !!--------------------------------------------------------------------- |
---|
| 237 | ! |
---|
| 238 | IF( nn_timing == 1 ) CALL timing_start('p5z_diat') |
---|
| 239 | ! |
---|
| 240 | |
---|
| 241 | DO jk = 1, jpkm1 |
---|
| 242 | DO jj = 1, jpj |
---|
| 243 | DO ji = 1, jpi |
---|
| 244 | |
---|
| 245 | zcompadi = MAX( ( trb(ji,jj,jk,jpdia) - 1E-9), 0. ) |
---|
| 246 | |
---|
| 247 | ! Aggregation term for diatoms is increased in case of nutrient |
---|
| 248 | ! stress as observed in reality. The stressed cells become more |
---|
| 249 | ! sticky and coagulate to sink quickly out of the euphotic zone |
---|
| 250 | ! ------------------------------------------------------------- |
---|
| 251 | zstep = xstep |
---|
| 252 | # if defined key_degrad |
---|
| 253 | zstep = zstep * facvol(ji,jj,jk) |
---|
| 254 | # endif |
---|
| 255 | ! Phytoplankton squared mortality |
---|
| 256 | ! ------------------------------- |
---|
| 257 | zlim2 = xlimdia(ji,jj,jk) * xlimdia(ji,jj,jk) |
---|
| 258 | zlim1 = 0.25 * ( 1. - zlim2 ) / ( 0.25 + zlim2 ) |
---|
| 259 | zrespp2 = 1.e6 * zstep * ( wchld + wchldm * zlim1 ) * xdiss(ji,jj,jk) * zcompadi * trb(ji,jj,jk,jpdia) |
---|
| 260 | |
---|
| 261 | ! Phytoplankton linear mortality |
---|
| 262 | ! ------------------------------ |
---|
| 263 | ztortp2 = mprat2 * xstep * zcompadi |
---|
| 264 | zmortp2 = zrespp2 + ztortp2 |
---|
| 265 | |
---|
| 266 | ! Update the arrays tra which contains the biological sources and sinks |
---|
| 267 | ! --------------------------------------------------------------------- |
---|
| 268 | zfactn = trb(ji,jj,jk,jpndi) / ( trb(ji,jj,jk,jpdia) + rtrn ) |
---|
| 269 | zfactp = trb(ji,jj,jk,jppdi) / ( trb(ji,jj,jk,jpdia) + rtrn ) |
---|
| 270 | zfactch = trb(ji,jj,jk,jpdch) / ( trb(ji,jj,jk,jpdia) + rtrn ) |
---|
| 271 | zfactfe = trb(ji,jj,jk,jpdfe) / ( trb(ji,jj,jk,jpdia) + rtrn ) |
---|
| 272 | zfactsi = trb(ji,jj,jk,jpdsi) / ( trb(ji,jj,jk,jpdia) + rtrn ) |
---|
| 273 | tra(ji,jj,jk,jpdia) = tra(ji,jj,jk,jpdia) - zmortp2 |
---|
| 274 | tra(ji,jj,jk,jpndi) = tra(ji,jj,jk,jpndi) - zmortp2 * zfactn |
---|
| 275 | tra(ji,jj,jk,jppdi) = tra(ji,jj,jk,jppdi) - zmortp2 * zfactp |
---|
| 276 | tra(ji,jj,jk,jpdch) = tra(ji,jj,jk,jpdch) - zmortp2 * zfactch |
---|
| 277 | tra(ji,jj,jk,jpdfe) = tra(ji,jj,jk,jpdfe) - zmortp2 * zfactfe |
---|
| 278 | tra(ji,jj,jk,jpdsi) = tra(ji,jj,jk,jpdsi) - zmortp2 * zfactsi |
---|
| 279 | tra(ji,jj,jk,jpgsi) = tra(ji,jj,jk,jpgsi) + zmortp2 * zfactsi |
---|
| 280 | #if defined key_kriest |
---|
| 281 | tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + zmortp2 |
---|
| 282 | tra(ji,jj,jk,jppon) = tra(ji,jj,jk,jppon) + zmortp2 * zfactn |
---|
| 283 | tra(ji,jj,jk,jppop) = tra(ji,jj,jk,jppop) + zmortp2 * zfactp |
---|
| 284 | tra(ji,jj,jk,jpnum) = tra(ji,jj,jk,jpnum) + ztortp2 * xkr_ddiat + zrespp2 * xkr_daggr |
---|
| 285 | tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + zmortp2 * zfactfe |
---|
| 286 | #else |
---|
| 287 | tra(ji,jj,jk,jpgoc) = tra(ji,jj,jk,jpgoc) + zrespp2 + 0.5 * ztortp2 |
---|
| 288 | tra(ji,jj,jk,jpgon) = tra(ji,jj,jk,jpgon) + (zrespp2 + 0.5 * ztortp2) * zfactn |
---|
| 289 | tra(ji,jj,jk,jpgop) = tra(ji,jj,jk,jpgop) + (zrespp2 + 0.5 * ztortp2) * zfactp |
---|
| 290 | tra(ji,jj,jk,jpbfe) = tra(ji,jj,jk,jpbfe) + (zrespp2 + 0.5 * ztortp2) * zfactfe |
---|
| 291 | tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) + 0.5 * ztortp2 |
---|
| 292 | tra(ji,jj,jk,jppon) = tra(ji,jj,jk,jppon) + 0.5 * ztortp2 * zfactn |
---|
| 293 | tra(ji,jj,jk,jppop) = tra(ji,jj,jk,jppop) + 0.5 * ztortp2 * zfactp |
---|
| 294 | tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) + 0.5 * ztortp2 * zfactfe |
---|
| 295 | prodpoc(ji,jj,jk) = prodpoc(ji,jj,jk) + 0.5 * ztortp2 |
---|
| 296 | prodgoc(ji,jj,jk) = prodgoc(ji,jj,jk) + zrespp2 + 0.5 * ztortp2 |
---|
| 297 | #endif |
---|
| 298 | END DO |
---|
| 299 | END DO |
---|
| 300 | END DO |
---|
| 301 | ! |
---|
| 302 | IF(ln_ctl) THEN ! print mean trends (used for debugging) |
---|
| 303 | WRITE(charout, FMT="('diat')") |
---|
| 304 | CALL prt_ctl_trc_info(charout) |
---|
| 305 | CALL prt_ctl_trc(tab4d=tra, mask=tmask, clinfo=ctrcnm) |
---|
| 306 | ENDIF |
---|
| 307 | ! |
---|
| 308 | IF( nn_timing == 1 ) CALL timing_stop('p5z_diat') |
---|
| 309 | ! |
---|
| 310 | END SUBROUTINE p5z_diat |
---|
| 311 | |
---|
| 312 | SUBROUTINE p5z_mort_init |
---|
| 313 | |
---|
| 314 | !!---------------------------------------------------------------------- |
---|
| 315 | !! *** ROUTINE p5z_mort_init *** |
---|
| 316 | !! |
---|
| 317 | !! ** Purpose : Initialization of phytoplankton parameters |
---|
| 318 | !! |
---|
| 319 | !! ** Method : Read the nampismort namelist and check the parameters |
---|
| 320 | !! called at the first timestep |
---|
| 321 | !! |
---|
| 322 | !! ** input : Namelist nampismort |
---|
| 323 | !! |
---|
| 324 | !!---------------------------------------------------------------------- |
---|
| 325 | |
---|
| 326 | NAMELIST/nampismort/ wchl, wchlp, wchld, wchldm, mprat, mpratp, mprat2 |
---|
| 327 | INTEGER :: ios ! Local integer output status for namelist read |
---|
| 328 | |
---|
| 329 | REWIND( numnatp_ref ) ! Namelist nampismort in reference namelist : Pisces phytoplankton |
---|
| 330 | READ ( numnatp_ref, nampismort, IOSTAT = ios, ERR = 901) |
---|
| 331 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'nampismort in reference namelist', lwp ) |
---|
| 332 | |
---|
| 333 | REWIND( numnatp_cfg ) ! Namelist nampismort in configuration namelist : Pisces phytoplankton |
---|
| 334 | READ ( numnatp_cfg, nampismort, IOSTAT = ios, ERR = 902 ) |
---|
| 335 | 902 IF( ios /= 0 ) CALL ctl_nam ( ios , 'nampismort in configuration namelist', lwp ) |
---|
| 336 | IF(lwm) WRITE ( numonp, nampismort ) |
---|
| 337 | |
---|
| 338 | IF(lwp) THEN ! control print |
---|
| 339 | WRITE(numout,*) ' ' |
---|
| 340 | WRITE(numout,*) ' Namelist parameters for phytoplankton mortality, nampismort' |
---|
| 341 | WRITE(numout,*) ' ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~' |
---|
| 342 | WRITE(numout,*) ' quadratic mortality of phytoplankton wchl =', wchl |
---|
| 343 | WRITE(numout,*) ' quadratic mortality of picophyto. wchlp =', wchlp |
---|
| 344 | WRITE(numout,*) ' quadratic mortality of diatoms wchld =', wchld |
---|
| 345 | WRITE(numout,*) ' Additional quadratic mortality of diatoms wchldm =', wchldm |
---|
| 346 | WRITE(numout,*) ' nanophyto. mortality rate mprat =', mprat |
---|
| 347 | WRITE(numout,*) ' picophyto. mortality rate mpratp =', mpratp |
---|
| 348 | WRITE(numout,*) ' Diatoms mortality rate mprat2 =', mprat2 |
---|
| 349 | ENDIF |
---|
| 350 | |
---|
| 351 | END SUBROUTINE p5z_mort_init |
---|
| 352 | |
---|
| 353 | #else |
---|
| 354 | !!====================================================================== |
---|
| 355 | !! Dummy module : No PISCES bio-model |
---|
| 356 | !!====================================================================== |
---|
| 357 | CONTAINS |
---|
| 358 | SUBROUTINE p5z_mort ! Empty routine |
---|
| 359 | END SUBROUTINE p5z_mort |
---|
| 360 | #endif |
---|
| 361 | |
---|
| 362 | !!====================================================================== |
---|
| 363 | END MODULE p5zmort |
---|