p4zmort.F90 13 KB

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