[3443] | 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 | !!---------------------------------------------------------------------- |
---|
[9169] | 9 | !! p4z_mort : Compute the mortality terms for phytoplankton |
---|
| 10 | !! p4z_mort_init : Initialize the mortality params for phytoplankton |
---|
[3443] | 11 | !!---------------------------------------------------------------------- |
---|
[9169] | 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 |
---|
[10227] | 16 | USE p4zlim ! Phytoplankton limitation terms |
---|
[13286] | 17 | USE prtctl ! print control for debugging |
---|
[3443] | 18 | |
---|
| 19 | IMPLICIT NONE |
---|
| 20 | PRIVATE |
---|
| 21 | |
---|
[15459] | 22 | PUBLIC p4z_mort ! Called from p4zbio.F90 |
---|
| 23 | PUBLIC p4z_mort_init ! Called from trcini_pisces.F90 |
---|
[3443] | 24 | |
---|
[15459] | 25 | REAL(wp), PUBLIC :: wchln !: Quadratic mortality rate of nanophytoplankton |
---|
| 26 | REAL(wp), PUBLIC :: wchld !: Quadratic mortality rate of diatoms |
---|
| 27 | REAL(wp), PUBLIC :: mpratn !: Linear mortality rate of nanophytoplankton |
---|
| 28 | REAL(wp), PUBLIC :: mpratd !: Linear mortality rate of diatoms |
---|
[3443] | 29 | |
---|
[12377] | 30 | !! * Substitutions |
---|
| 31 | # include "do_loop_substitute.h90" |
---|
[3443] | 32 | !!---------------------------------------------------------------------- |
---|
[10067] | 33 | !! NEMO/TOP 4.0 , NEMO Consortium (2018) |
---|
[7753] | 34 | !! $Id$ |
---|
[10068] | 35 | !! Software governed by the CeCILL license (see ./LICENSE) |
---|
[3443] | 36 | !!---------------------------------------------------------------------- |
---|
| 37 | CONTAINS |
---|
| 38 | |
---|
[12377] | 39 | SUBROUTINE p4z_mort( kt, Kbb, Krhs ) |
---|
[3443] | 40 | !!--------------------------------------------------------------------- |
---|
| 41 | !! *** ROUTINE p4z_mort *** |
---|
| 42 | !! |
---|
[15459] | 43 | !! ** Purpose : Calls the different subroutine to compute |
---|
[3443] | 44 | !! the different phytoplankton mortality terms |
---|
| 45 | !! |
---|
| 46 | !! ** Method : - ??? |
---|
| 47 | !!--------------------------------------------------------------------- |
---|
| 48 | INTEGER, INTENT(in) :: kt ! ocean time step |
---|
[12377] | 49 | INTEGER, INTENT(in) :: Kbb, Krhs ! time level indices |
---|
[3443] | 50 | !!--------------------------------------------------------------------- |
---|
[9169] | 51 | ! |
---|
[15459] | 52 | CALL p4z_mort_nano( Kbb, Krhs ) ! nanophytoplankton |
---|
| 53 | CALL p4z_mort_diat( Kbb, Krhs ) ! diatoms |
---|
[9169] | 54 | ! |
---|
[3443] | 55 | END SUBROUTINE p4z_mort |
---|
| 56 | |
---|
| 57 | |
---|
[15459] | 58 | SUBROUTINE p4z_mort_nano( Kbb, Krhs ) |
---|
[3443] | 59 | !!--------------------------------------------------------------------- |
---|
[15459] | 60 | !! *** ROUTINE p4z_mort_nano *** |
---|
[3443] | 61 | !! |
---|
| 62 | !! ** Purpose : Compute the mortality terms for nanophytoplankton |
---|
| 63 | !! |
---|
[15459] | 64 | !! ** Method : Both quadratic and simili linear mortality terms |
---|
[3443] | 65 | !!--------------------------------------------------------------------- |
---|
[12377] | 66 | INTEGER, INTENT(in) :: Kbb, Krhs ! time level indices |
---|
[9169] | 67 | INTEGER :: ji, jj, jk |
---|
[15459] | 68 | REAL(wp) :: zcompaph |
---|
[9169] | 69 | REAL(wp) :: zfactfe, zfactch, zprcaca, zfracal |
---|
[15459] | 70 | REAL(wp) :: ztortp , zrespp , zmortp, zlim1, zlim2 |
---|
[9169] | 71 | CHARACTER (len=25) :: charout |
---|
[3443] | 72 | !!--------------------------------------------------------------------- |
---|
| 73 | ! |
---|
[15459] | 74 | IF( ln_timing ) CALL timing_start('p4z_mort_nano') |
---|
[3443] | 75 | ! |
---|
[9169] | 76 | prodcal(:,:,:) = 0._wp ! calcite production variable set to zero |
---|
[15090] | 77 | DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1) |
---|
[15459] | 78 | zcompaph = MAX( ( tr(ji,jj,jk,jpphy,Kbb) - 1e-9 ), 0.e0 ) |
---|
[3443] | 79 | |
---|
[15459] | 80 | ! Quadratic mortality of nano due to aggregation during |
---|
| 81 | ! blooms (Doney et al. 1996) |
---|
| 82 | ! ----------------------------------------------------- |
---|
| 83 | zlim2 = xlimphy(ji,jj,jk) * xlimphy(ji,jj,jk) |
---|
| 84 | zlim1 = 0.25 * ( 1. - zlim2 ) / ( 0.25 + zlim2 ) * tr(ji,jj,jk,jpphy,Kbb) |
---|
| 85 | zrespp = wchln * 1.e6 * xstep * zlim1 * xdiss(ji,jj,jk) * zcompaph |
---|
[3443] | 86 | |
---|
[15459] | 87 | ! Phytoplankton linear mortality |
---|
| 88 | ! A michaelis-menten like term is introduced to avoid |
---|
| 89 | ! extinction of nanophyto in highly limited areas |
---|
| 90 | ! ---------------------------------------------------- |
---|
| 91 | ztortp = mpratn * xstep * zcompaph / ( xkmort + tr(ji,jj,jk,jpphy,Kbb) ) * tr(ji,jj,jk,jpphy,Kbb) |
---|
| 92 | |
---|
[12377] | 93 | zmortp = zrespp + ztortp |
---|
[15459] | 94 | |
---|
[12377] | 95 | ! Update the arrays TRA which contains the biological sources and sinks |
---|
| 96 | zfactfe = tr(ji,jj,jk,jpnfe,Kbb)/(tr(ji,jj,jk,jpphy,Kbb)+rtrn) |
---|
| 97 | zfactch = tr(ji,jj,jk,jpnch,Kbb)/(tr(ji,jj,jk,jpphy,Kbb)+rtrn) |
---|
| 98 | tr(ji,jj,jk,jpphy,Krhs) = tr(ji,jj,jk,jpphy,Krhs) - zmortp |
---|
| 99 | tr(ji,jj,jk,jpnch,Krhs) = tr(ji,jj,jk,jpnch,Krhs) - zmortp * zfactch |
---|
| 100 | tr(ji,jj,jk,jpnfe,Krhs) = tr(ji,jj,jk,jpnfe,Krhs) - zmortp * zfactfe |
---|
[15459] | 101 | |
---|
| 102 | ! Production PIC particles due to mortality |
---|
[12377] | 103 | zprcaca = xfracal(ji,jj,jk) * zmortp |
---|
| 104 | prodcal(ji,jj,jk) = prodcal(ji,jj,jk) + zprcaca ! prodcal=prodcal(nanophy)+prodcal(microzoo)+prodcal(mesozoo) |
---|
[15459] | 105 | |
---|
| 106 | ! POC associated with the shell is supposed to be routed to |
---|
| 107 | ! big particles because of the ballasting effect |
---|
[12377] | 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 |
---|
[15459] | 112 | tr(ji,jj,jk,jppoc,Krhs) = tr(ji,jj,jk,jppoc,Krhs) + zmortp |
---|
| 113 | prodpoc(ji,jj,jk) = prodpoc(ji,jj,jk) + zmortp |
---|
| 114 | |
---|
| 115 | ! Update the arrays TRA which contains the biological sources and sinks |
---|
| 116 | tr(ji,jj,jk,jpsfe,Krhs) = tr(ji,jj,jk,jpsfe,Krhs) + zmortp * zfactfe |
---|
| 117 | ! |
---|
[12377] | 118 | END_3D |
---|
[3443] | 119 | ! |
---|
[12377] | 120 | IF(sn_cfctl%l_prttrc) THEN ! print mean trends (used for debugging) |
---|
[3443] | 121 | WRITE(charout, FMT="('nano')") |
---|
[13286] | 122 | CALL prt_ctl_info( charout, cdcomp = 'top' ) |
---|
| 123 | CALL prt_ctl(tab4d_1=tr(:,:,:,:,Krhs), mask1=tmask, clinfo=ctrcnm) |
---|
[3443] | 124 | ENDIF |
---|
| 125 | ! |
---|
[15459] | 126 | IF( ln_timing ) CALL timing_stop('p4z_mort_nano') |
---|
[3443] | 127 | ! |
---|
[15459] | 128 | END SUBROUTINE p4z_mort_nano |
---|
[3443] | 129 | |
---|
[9124] | 130 | |
---|
[15459] | 131 | SUBROUTINE p4z_mort_diat( Kbb, Krhs ) |
---|
[3443] | 132 | !!--------------------------------------------------------------------- |
---|
[15459] | 133 | !! *** ROUTINE p4z_mort_diat *** |
---|
[3443] | 134 | !! |
---|
| 135 | !! ** Purpose : Compute the mortality terms for diatoms |
---|
| 136 | !! |
---|
[15459] | 137 | !! ** Method : - Both quadratic and simili linear mortality terms |
---|
[3443] | 138 | !!--------------------------------------------------------------------- |
---|
[12377] | 139 | INTEGER, INTENT(in) :: Kbb, Krhs ! time level indices |
---|
[9169] | 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 |
---|
[3443] | 145 | !!--------------------------------------------------------------------- |
---|
| 146 | ! |
---|
[15459] | 147 | IF( ln_timing ) CALL timing_start('p4z_mort_diat') |
---|
[3443] | 148 | ! |
---|
[15459] | 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 | ! This is due to the production of EPS by stressed cells |
---|
| 153 | ! ------------------------------------------------------------- |
---|
[3443] | 154 | |
---|
[15090] | 155 | DO_3D( nn_hls, nn_hls, nn_hls, nn_hls, 1, jpkm1) |
---|
[3443] | 156 | |
---|
[12377] | 157 | zcompadi = MAX( ( tr(ji,jj,jk,jpdia,Kbb) - 1e-9), 0. ) |
---|
[3443] | 158 | |
---|
[15459] | 159 | ! Aggregation term for diatoms is increased in case of nutrient |
---|
| 160 | ! stress as observed in reality. The stressed cells become more |
---|
| 161 | ! sticky and coagulate to sink quickly out of the euphotic zone |
---|
| 162 | ! ------------------------------------------------------------ |
---|
[12377] | 163 | zlim2 = xlimdia(ji,jj,jk) * xlimdia(ji,jj,jk) |
---|
| 164 | zlim1 = 0.25 * ( 1. - zlim2 ) / ( 0.25 + zlim2 ) |
---|
[15459] | 165 | zrespp2 = 1.e6 * xstep * wchld * zlim1 * xdiss(ji,jj,jk) * zcompadi * tr(ji,jj,jk,jpdia,Kbb) |
---|
[3443] | 166 | |
---|
[15459] | 167 | ! Phytoplankton linear mortality |
---|
| 168 | ! A michaelis-menten like term is introduced to avoid |
---|
| 169 | ! extinction of diatoms in highly limited areas |
---|
| 170 | ! --------------------------------------------------- |
---|
| 171 | ztortp2 = mpratd * xstep * tr(ji,jj,jk,jpdia,Kbb) / ( xkmort + tr(ji,jj,jk,jpdia,Kbb) ) * zcompadi |
---|
[3443] | 172 | |
---|
[12377] | 173 | zmortp2 = zrespp2 + ztortp2 |
---|
[3443] | 174 | |
---|
[15459] | 175 | ! Update the arrays trends which contains the biological sources and sinks |
---|
[12377] | 176 | ! --------------------------------------------------------------------- |
---|
| 177 | zfactch = tr(ji,jj,jk,jpdch,Kbb) / ( tr(ji,jj,jk,jpdia,Kbb) + rtrn ) |
---|
| 178 | zfactfe = tr(ji,jj,jk,jpdfe,Kbb) / ( tr(ji,jj,jk,jpdia,Kbb) + rtrn ) |
---|
| 179 | zfactsi = tr(ji,jj,jk,jpdsi,Kbb) / ( tr(ji,jj,jk,jpdia,Kbb) + rtrn ) |
---|
| 180 | tr(ji,jj,jk,jpdia,Krhs) = tr(ji,jj,jk,jpdia,Krhs) - zmortp2 |
---|
| 181 | tr(ji,jj,jk,jpdch,Krhs) = tr(ji,jj,jk,jpdch,Krhs) - zmortp2 * zfactch |
---|
| 182 | tr(ji,jj,jk,jpdfe,Krhs) = tr(ji,jj,jk,jpdfe,Krhs) - zmortp2 * zfactfe |
---|
| 183 | tr(ji,jj,jk,jpdsi,Krhs) = tr(ji,jj,jk,jpdsi,Krhs) - zmortp2 * zfactsi |
---|
| 184 | tr(ji,jj,jk,jpgsi,Krhs) = tr(ji,jj,jk,jpgsi,Krhs) + zmortp2 * zfactsi |
---|
[15459] | 185 | |
---|
| 186 | ! Half of the linear mortality term is routed to big particles |
---|
| 187 | ! becaue of the ballasting effect |
---|
| 188 | tr(ji,jj,jk,jpgoc,Krhs) = tr(ji,jj,jk,jpgoc,Krhs) + zrespp2 |
---|
| 189 | tr(ji,jj,jk,jppoc,Krhs) = tr(ji,jj,jk,jppoc,Krhs) + ztortp2 |
---|
| 190 | prodpoc(ji,jj,jk) = prodpoc(ji,jj,jk) + ztortp2 |
---|
| 191 | prodgoc(ji,jj,jk) = prodgoc(ji,jj,jk) + zrespp2 |
---|
| 192 | tr(ji,jj,jk,jpsfe,Krhs) = tr(ji,jj,jk,jpsfe,Krhs) + ztortp2 * zfactfe |
---|
| 193 | tr(ji,jj,jk,jpbfe,Krhs) = tr(ji,jj,jk,jpbfe,Krhs) + zrespp2 * zfactfe |
---|
[12377] | 194 | END_3D |
---|
[3443] | 195 | ! |
---|
[12377] | 196 | IF(sn_cfctl%l_prttrc) THEN ! print mean trends (used for debugging) |
---|
[3443] | 197 | WRITE(charout, FMT="('diat')") |
---|
[13286] | 198 | CALL prt_ctl_info( charout, cdcomp = 'top' ) |
---|
| 199 | CALL prt_ctl(tab4d_1=tr(:,:,:,:,Krhs), mask1=tmask, clinfo=ctrcnm) |
---|
[3443] | 200 | ENDIF |
---|
| 201 | ! |
---|
[15459] | 202 | IF( ln_timing ) CALL timing_stop('p4z_mort_diat') |
---|
[3443] | 203 | ! |
---|
[15459] | 204 | END SUBROUTINE p4z_mort_diat |
---|
[3443] | 205 | |
---|
[9124] | 206 | |
---|
[3443] | 207 | SUBROUTINE p4z_mort_init |
---|
| 208 | !!---------------------------------------------------------------------- |
---|
| 209 | !! *** ROUTINE p4z_mort_init *** |
---|
| 210 | !! |
---|
| 211 | !! ** Purpose : Initialization of phytoplankton parameters |
---|
| 212 | !! |
---|
[15459] | 213 | !! ** Method : Read the namp4zmort namelist and check the parameters |
---|
[9169] | 214 | !! called at the first timestep |
---|
[3443] | 215 | !! |
---|
[15459] | 216 | !! ** input : Namelist namp4zmort |
---|
[3443] | 217 | !! |
---|
| 218 | !!---------------------------------------------------------------------- |
---|
[9124] | 219 | INTEGER :: ios ! Local integer |
---|
| 220 | ! |
---|
[15459] | 221 | NAMELIST/namp4zmort/ wchln, wchld, mpratn, mpratd |
---|
[9124] | 222 | !!---------------------------------------------------------------------- |
---|
| 223 | ! |
---|
[9169] | 224 | IF(lwp) THEN |
---|
| 225 | WRITE(numout,*) |
---|
| 226 | WRITE(numout,*) 'p4z_mort_init : Initialization of phytoplankton mortality parameters' |
---|
| 227 | WRITE(numout,*) '~~~~~~~~~~~~~' |
---|
| 228 | ENDIF |
---|
| 229 | ! |
---|
[7646] | 230 | READ ( numnatp_ref, namp4zmort, IOSTAT = ios, ERR = 901) |
---|
[11536] | 231 | 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namp4zmort in reference namelist' ) |
---|
[15459] | 232 | |
---|
[7646] | 233 | READ ( numnatp_cfg, namp4zmort, IOSTAT = ios, ERR = 902 ) |
---|
[11536] | 234 | 902 IF( ios > 0 ) CALL ctl_nam ( ios , 'namp4zmort in configuration namelist' ) |
---|
[9169] | 235 | IF(lwm) WRITE( numonp, namp4zmort ) |
---|
[9124] | 236 | ! |
---|
[3443] | 237 | IF(lwp) THEN ! control print |
---|
[9169] | 238 | WRITE(numout,*) ' Namelist : namp4zmort' |
---|
[15459] | 239 | WRITE(numout,*) ' quadratic mortality of phytoplankton wchln =', wchln |
---|
[9169] | 240 | WRITE(numout,*) ' maximum quadratic mortality of diatoms wchld =', wchld |
---|
[15459] | 241 | WRITE(numout,*) ' phytoplankton mortality rate mpratn =', mpratn |
---|
| 242 | WRITE(numout,*) ' Diatoms mortality rate mpratd =', mpratd |
---|
[3443] | 243 | ENDIF |
---|
[9124] | 244 | ! |
---|
[3443] | 245 | END SUBROUTINE p4z_mort_init |
---|
| 246 | |
---|
| 247 | !!====================================================================== |
---|
[5656] | 248 | END MODULE p4zmort |
---|