MODULE p4zmeso !!====================================================================== !! *** MODULE p4zmeso *** !! TOP : PISCES Compute the sources/sinks for mesozooplankton !!====================================================================== !! History : 1.0 ! 2002 (O. Aumont) Original code !! 2.0 ! 2007-12 (C. Ethe, G. Madec) F90 !! 3.4 ! 2011-06 (O. Aumont, C. Ethe) Quota model for iron !!---------------------------------------------------------------------- !! p4z_meso : Compute the sources/sinks for mesozooplankton !! p4z_meso_init : Initialization of the parameters for mesozooplankton !! p4z_meso_alloc : Allocate variables for mesozooplankton !!---------------------------------------------------------------------- USE oce_trc ! shared variables between ocean and passive tracers USE trc ! passive tracers common variables USE sms_pisces ! PISCES Source Minus Sink variables USE p4zprod ! production USE prtctl_trc ! print control for debugging USE iom ! I/O manager IMPLICIT NONE PRIVATE PUBLIC p4z_meso ! called in p4zbio.F90 PUBLIC p4z_meso_init ! called in trcsms_pisces.F90 PUBLIC p4z_meso_alloc ! called in trcini_pisces.F90 !! * Shared module variables REAL(wp), PUBLIC :: part2 !: part of calcite not dissolved in mesozoo guts REAL(wp), PUBLIC :: xpref2d !: mesozoo preference for diatoms REAL(wp), PUBLIC :: xpref2n !: mesozoo preference for nanophyto REAL(wp), PUBLIC :: xpref2z !: mesozoo preference for microzooplankton REAL(wp), PUBLIC :: xpref2c !: mesozoo preference for POC REAL(wp), PUBLIC :: xthresh2zoo !: zoo feeding threshold for mesozooplankton REAL(wp), PUBLIC :: xthresh2dia !: diatoms feeding threshold for mesozooplankton REAL(wp), PUBLIC :: xthresh2phy !: nanophyto feeding threshold for mesozooplankton REAL(wp), PUBLIC :: xthresh2poc !: poc feeding threshold for mesozooplankton REAL(wp), PUBLIC :: xthresh2 !: feeding threshold for mesozooplankton REAL(wp), PUBLIC :: resrat2 !: exsudation rate of mesozooplankton REAL(wp), PUBLIC :: mzrat2 !: microzooplankton mortality rate REAL(wp), PUBLIC :: grazrat2 !: maximal mesozoo grazing rate REAL(wp), PUBLIC :: xkgraz2 !: non assimilated fraction of P by mesozoo REAL(wp), PUBLIC :: unass2 !: Efficicency of mesozoo growth REAL(wp), PUBLIC :: sigma2 !: Fraction of mesozoo excretion as DOM REAL(wp), PUBLIC :: epsher2 !: growth efficiency REAL(wp), PUBLIC :: epsher2min !: minimum growth efficiency at high food for grazing 2 REAL(wp), PUBLIC :: xsigma2 !: Width of the predation window REAL(wp), PUBLIC :: xsigma2del !: Maximum width of the predation window at low food density REAL(wp), PUBLIC :: grazflux !: mesozoo flux feeding rate REAL(wp), PUBLIC :: xfracmig !: Fractional biomass of meso that performs DVM LOGICAL , PUBLIC :: ln_dvm_meso !: Boolean to activate DVM of mesozooplankton REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: depmig !: DVM of mesozooplankton : migration depth INTEGER , ALLOCATABLE, SAVE, DIMENSION(:,:) :: kmig !: Vertical indice of the the migration depth !!---------------------------------------------------------------------- !! NEMO/TOP 4.0 , NEMO Consortium (2018) !! $Id$ !! Software governed by the CeCILL license (see ./LICENSE) !!---------------------------------------------------------------------- CONTAINS SUBROUTINE p4z_meso( kt, knt ) !!--------------------------------------------------------------------- !! *** ROUTINE p4z_meso *** !! !! ** Purpose : Compute the sources/sinks for mesozooplankton !! This includes ingestion and assimilation, flux feeding !! and mortality. We use a passive prey switching !! parameterization. !! All living compartments smaller than mesozooplankton !! are potential preys of mesozooplankton as well as small !! sinking particles !! !! ** Method : - ??? !!--------------------------------------------------------------------- INTEGER, INTENT(in) :: kt, knt ! ocean time step and ??? ! INTEGER :: ji, jj, jk, jkt REAL(wp) :: zcompadi, zcompaph, zcompapoc, zcompaz, zcompam REAL(wp) :: zgraze2 , zdenom, zdenom2, zfact , zfood, zfoodlim, zproport, zbeta REAL(wp) :: zmortzgoc, zfrac, zfracfe, zratio, zratio2, zfracal, zgrazcal REAL(wp) :: zepsherf, zepshert, zepsherq, zepsherv, zgrarsig, zgraztotc, zgraztotn, zgraztotf REAL(wp) :: zmigreltime, zprcaca, zmortz, zgrasratf, zgrasratn REAL(wp) :: zrespz, ztortz, zgrazdc, zgrazz, zgrazpof, zgraznc, zgrazpoc, zgraznf, zgrazdf REAL(wp) :: zgrazfffp, zgrazfffg, zgrazffep, zgrazffeg, zrum, zcodel, zargu, zval REAL(wp) :: zsigma, zdiffdn, ztmp1, ztmp2, ztmp3, ztmp4, ztmptot CHARACTER (len=25) :: charout REAL(wp), DIMENSION(jpi,jpj,jpk) :: zgrazing, zfezoo2 REAL(wp), DIMENSION(jpi,jpj,jpk) :: zgrarem, zgraref, zgrapoc, zgrapof, zgrabsi REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zgramigrem, zgramigref, zgramigpoc, zgramigpof REAL(wp), ALLOCATABLE, DIMENSION(:,:) :: zstrn, zgramigbsi REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: zw3d, zz2ligprod !!--------------------------------------------------------------------- ! IF( ln_timing ) CALL timing_start('p4z_meso') ! zgrazing(:,:,:) = 0._wp ; zgrapoc(:,:,:) = 0._wp zfezoo2 (:,:,:) = 0._wp ; zgrarem(:,:,:) = 0._wp zgraref (:,:,:) = 0._wp ; zgrapof(:,:,:) = 0._wp zgrabsi (:,:,:) = 0._wp ! IF (ln_ligand) THEN ALLOCATE( zz2ligprod(jpi,jpj,jpk) ) zz2ligprod(:,:,:) = 0._wp ENDIF ! ! Diurnal vertical migration of mesozooplankton ! Computation of the migration depth ! --------------------------------------------- IF (ln_dvm_meso) CALL p4z_meso_depmig ! DO jk = 1, jpk DO jj = 1, jpj DO ji = 1, jpi zcompam = MAX( ( trb(ji,jj,jk,jpmes) - 1.e-9 ), 0.e0 ) zfact = xstep * tgfunc2(ji,jj,jk) * zcompam ! linear mortality of mesozooplankton ! A michaelis menten modulation term is used to avoid extinction of ! mesozooplankton at very low food concentration. Mortality is ! enhanced in low O2 waters ! ----------------------------------------------------------------- zrespz = resrat2 * zfact * ( trb(ji,jj,jk,jpmes) / ( xkmort + trb(ji,jj,jk,jpmes) ) & & + 3. * nitrfac(ji,jj,jk) ) ! Zooplankton quadratic mortality. A square function has been selected with ! to mimic predation and disease (density dependent mortality). It also tends ! to stabilise the model ! ------------------------------------------------------------------------- ztortz = mzrat2 * 1.e6 * zfact * trb(ji,jj,jk,jpmes) * (1. - nitrfac(ji,jj,jk) ) ! Computation of the abundance of the preys ! A threshold can be specified in the namelist ! -------------------------------------------- zcompadi = MAX( ( trb(ji,jj,jk,jpdia) - xthresh2dia ), 0.e0 ) zcompaz = MAX( ( trb(ji,jj,jk,jpzoo) - xthresh2zoo ), 0.e0 ) zcompapoc = MAX( ( trb(ji,jj,jk,jppoc) - xthresh2poc ), 0.e0 ) ! Size effect of nanophytoplankton on grazing : the smaller it is, the less prone ! it is to predation by mesozooplankton. We use a quota dependant parameterization ! as a low quota indicates oligotrophic conditions which are charatcerized by ! small cells ! ------------------------------------------------------------------------------- zcompaph = MAX( ( trb(ji,jj,jk,jpphy) - xthresh2phy ), 0.e0 ) & & * MIN(1., MAX( 0., ( quotan(ji,jj,jk) - 0.2) / 0.3 ) ) ! Mesozooplankton grazing ! The total amount of food is the sum of all preys accessible to mesozooplankton ! multiplied by their food preference ! A threshold can be specified in the namelist (xthresh2). However, when food ! concentration is close to this threshold, it is decreased to avoid the ! accumulation of food in the mesozoopelagic domain ! ------------------------------------------------------------------------------- zfood = xpref2d * zcompadi + xpref2z * zcompaz + xpref2n * zcompaph + xpref2c * zcompapoc zfoodlim = MAX( 0., zfood - MIN( 0.5 * zfood, xthresh2 ) ) zdenom = zfoodlim / ( xkgraz2 + zfoodlim ) zdenom2 = zdenom / ( zfood + rtrn ) zgraze2 = grazrat2 * xstep * tgfunc2(ji,jj,jk) * trb(ji,jj,jk,jpmes) * (1. - nitrfac(ji,jj,jk)) ! An active switching parameterization is used here. ! We don't use the KTW parameterization proposed by ! Vallina et al. because it tends to produce too steady biomass ! composition and the variance of Chl is too low as it grazes ! too strongly on winning organisms. We use a generalized ! switching parameterization proposed by Morozov and ! Petrovskii (2013) ! ------------------------------------------------------------ ! The width of the selection window is increased when preys ! have low abundance, .i.e. zooplankton become less specific ! to avoid starvation. ! ---------------------------------------------------------- zsigma = 1.0 - zdenom**2/(0.05**2+zdenom**2) zsigma = xsigma2 + xsigma2del * zsigma ! Nanophytoplankton and diatoms are the only preys considered ! to be close enough to have potential interference ! ----------------------------------------------------------- zdiffdn = exp( -ABS(log(1.67 * sizen(ji,jj,jk) / (5.0 * sized(ji,jj,jk) + rtrn )) )**2 / zsigma**2 ) ztmp1 = xpref2n * zcompaph * ( zcompaph + zdiffdn * zcompadi ) / ( 1.0 + zdiffdn ) ztmp2 = xpref2c * zcompapoc**2 ztmp3 = xpref2d * zcompadi * ( zdiffdn * zcompadi + zcompaph ) / ( 1.0 + zdiffdn ) ztmp4 = xpref2z * zcompaz**2 ztmptot = ztmp1 + ztmp2 + ztmp3 + ztmp4 + rtrn ztmp1 = ztmp1 / ztmptot ztmp2 = ztmp2 / ztmptot ztmp3 = ztmp3 / ztmptot ztmp4 = ztmp4 / ztmptot ! Mesozooplankton regular grazing on the different preys ! ------------------------------------------------------ zgrazdc = zgraze2 * ztmp3 * zdenom ! diatoms zgraznc = zgraze2 * ztmp1 * zdenom ! nanophytoplankton zgrazpoc = zgraze2 * ztmp2 * zdenom ! small POC zgrazz = zgraze2 * ztmp4 * zdenom ! microzooplankton ! Ingestion rates of the Fe content of the different preys zgraznf = zgraznc * trb(ji,jj,jk,jpnfe) / ( trb(ji,jj,jk,jpphy) + rtrn) zgrazdf = zgrazdc * trb(ji,jj,jk,jpdfe) / ( trb(ji,jj,jk,jpdia) + rtrn) zgrazpof = zgrazpoc * trb(ji,jj,jk,jpsfe) / ( trb(ji,jj,jk,jppoc) + rtrn) ! Mesozooplankton flux feeding on GOC and POC. The feeding pressure ! is proportional to the flux ! ------------------------------------------------------------------ zgrazffeg = grazflux * xstep * wsbio4(ji,jj,jk) & & * tgfunc2(ji,jj,jk) * trb(ji,jj,jk,jpgoc) * trb(ji,jj,jk,jpmes) & & * (1. - nitrfac(ji,jj,jk)) zgrazfffg = zgrazffeg * trb(ji,jj,jk,jpbfe) / (trb(ji,jj,jk,jpgoc) + rtrn) zgrazffep = grazflux * xstep * wsbio3(ji,jj,jk) & & * tgfunc2(ji,jj,jk) * trb(ji,jj,jk,jppoc) * trb(ji,jj,jk,jpmes) & & * (1. - nitrfac(ji,jj,jk)) zgrazfffp = zgrazffep * trb(ji,jj,jk,jpsfe) / (trb(ji,jj,jk,jppoc) + rtrn) zgraztotc = zgrazdc + zgrazz + zgraznc + zgrazpoc + zgrazffep + zgrazffeg ! Compute the proportion of filter feeders. It is assumed steady state. ! --------------------------------------------------------------------- zproport = (zgrazffep + zgrazffeg)/(rtrn + zgraztotc) ! Compute fractionation of aggregates. It is assumed that ! diatoms based aggregates are more prone to fractionation ! since they are more porous (marine snow instead of fecal pellets) ! ----------------------------------------------------------------- zratio = trb(ji,jj,jk,jpgsi) / ( trb(ji,jj,jk,jpgoc) + rtrn ) zratio2 = zratio * zratio zfrac = zproport * grazflux * xstep * wsbio4(ji,jj,jk) & & * trb(ji,jj,jk,jpgoc) * trb(ji,jj,jk,jpmes) & & * ( 0.2 + 3.8 * zratio2 / ( 1.**2 + zratio2 ) ) zfracfe = zfrac * trb(ji,jj,jk,jpbfe) / (trb(ji,jj,jk,jpgoc) + rtrn) ! Flux feeding is multiplied by the fractional biomass of flux feeders zgrazffep = zproport * zgrazffep zgrazffeg = zproport * zgrazffeg zgrazfffp = zproport * zgrazfffp zgrazfffg = zproport * zgrazfffg ! Total ingestion rates in C, N, Fe zgraztotc = zgrazdc + zgrazz + zgraznc + zgrazpoc + zgrazffep + zgrazffeg zgraztotn = zgrazdc * quotad(ji,jj,jk) + zgrazz + zgraznc * quotan(ji,jj,jk) & & + zgrazpoc + zgrazffep + zgrazffeg zgraztotf = zgrazdf + zgraznf + zgrazz * ferat3 + zgrazpof + zgrazfffp + zgrazfffg ! Total grazing ( grazing by microzoo is already computed in p4zmicro ) zgrazing(ji,jj,jk) = zgraztotc ! Mesozooplankton efficiency. ! We adopt a formulation proposed by Mitra et al. (2007) ! The gross growth efficiency is controled by the most limiting nutrient. ! Growth is also further decreased when the food quality is poor. This is currently ! hard coded : it can be decreased by up to 50% (zepsherq) ! GGE can also be decreased when food quantity is high, zepsherf (Montagnes and ! Fulton, 2012) ! ----------------------------------------------------------------------------------- zgrasratf = ( zgraztotf + rtrn )/ ( zgraztotc + rtrn ) zgrasratn = ( zgraztotn + rtrn )/ ( zgraztotc + rtrn ) zepshert = MIN( 1., zgrasratn, zgrasratf / ferat3) zbeta = MAX(0., (epsher2 - epsher2min) ) ! Food quantity deprivation of GGE zepsherf = epsher2min + zbeta / ( 1.0 + 0.04E6 * 12. * zfood * zbeta ) ! Food quality deprivation of GGE zepsherq = 0.5 + (1.0 - 0.5) * zepshert * ( 1.0 + 1.0 ) / ( zepshert + 1.0 ) ! Actual GGE zepsherv = zepsherf * zepshert * zepsherq ! ! Impact of grazing on the prognostic variables ! --------------------------------------------- zmortz = ztortz + zrespz ! Mortality induced by the upper trophic levels, ztortz, is allocated ! according to a infinite chain of predators (ANderson et al., 2013) zmortzgoc = unass2 / ( 1. - epsher2 ) * ztortz + zrespz ! Update of the trends tra(ji,jj,jk,jpmes) = tra(ji,jj,jk,jpmes) - zmortz + zepsherv * zgraztotc tra(ji,jj,jk,jpdia) = tra(ji,jj,jk,jpdia) - zgrazdc tra(ji,jj,jk,jpzoo) = tra(ji,jj,jk,jpzoo) - zgrazz tra(ji,jj,jk,jpphy) = tra(ji,jj,jk,jpphy) - zgraznc tra(ji,jj,jk,jpnch) = tra(ji,jj,jk,jpnch) - zgraznc * trb(ji,jj,jk,jpnch) / ( trb(ji,jj,jk,jpphy) + rtrn ) tra(ji,jj,jk,jpdch) = tra(ji,jj,jk,jpdch) - zgrazdc * trb(ji,jj,jk,jpdch) / ( trb(ji,jj,jk,jpdia) + rtrn ) tra(ji,jj,jk,jpdsi) = tra(ji,jj,jk,jpdsi) - zgrazdc * trb(ji,jj,jk,jpdsi) / ( trb(ji,jj,jk,jpdia) + rtrn ) zgrabsi(ji,jj,jk) = zgrazdc * trb(ji,jj,jk,jpdsi) / ( trb(ji,jj,jk,jpdia) + rtrn ) tra(ji,jj,jk,jpnfe) = tra(ji,jj,jk,jpnfe) - zgraznf tra(ji,jj,jk,jpdfe) = tra(ji,jj,jk,jpdfe) - zgrazdf tra(ji,jj,jk,jppoc) = tra(ji,jj,jk,jppoc) - zgrazpoc - zgrazffep + zfrac prodpoc(ji,jj,jk) = prodpoc(ji,jj,jk) + zfrac conspoc(ji,jj,jk) = conspoc(ji,jj,jk) - zgrazpoc - zgrazffep tra(ji,jj,jk,jpgoc) = tra(ji,jj,jk,jpgoc) - zgrazffeg - zfrac consgoc(ji,jj,jk) = consgoc(ji,jj,jk) - zgrazffeg - zfrac tra(ji,jj,jk,jpsfe) = tra(ji,jj,jk,jpsfe) - zgrazpof - zgrazfffp + zfracfe tra(ji,jj,jk,jpbfe) = tra(ji,jj,jk,jpbfe) - zgrazfffg - zfracfe ! Calcite remineralization due to zooplankton activity ! part2 of the ingested calcite is not dissolving in the ! acidic gut ! ------------------------------------------------------ zfracal = trb(ji,jj,jk,jpcal) / ( trb(ji,jj,jk,jpgoc) + rtrn ) zgrazcal = zgrazffeg * (1. - part2) * zfracal ! calcite production by zooplankton activity zprcaca = xfracal(ji,jj,jk) * zgraznc prodcal(ji,jj,jk) = prodcal(ji,jj,jk) + zprcaca ! prodcal=prodcal(nanophy)+prodcal(microzoo)+prodcal(mesozoo) ! zprcaca = part2 * zprcaca tra(ji,jj,jk,jpdic) = tra(ji,jj,jk,jpdic) + zgrazcal - zprcaca tra(ji,jj,jk,jptal) = tra(ji,jj,jk,jptal) - 2. * ( zgrazcal + zprcaca ) tra(ji,jj,jk,jpcal) = tra(ji,jj,jk,jpcal) - zgrazcal + zprcaca ! Computation of total excretion and egestion by mesozoo. ! --------------------------------------------------------- zgrarem(ji,jj,jk) = zgraztotc * ( 1. - zepsherv - unass2 ) & & + ( 1. - epsher2 - unass2 ) / ( 1. - epsher2 ) * ztortz zgraref(ji,jj,jk) = zgraztotc * MAX( 0. , ( 1. - unass2 ) * zgrasratf - ferat3 * zepsherv ) & & + ferat3 * ( ( 1. - epsher2 - unass2 ) /( 1. - epsher2 ) * ztortz ) zgrapoc(ji,jj,jk) = zgraztotc * unass2 + zmortzgoc zgrapof(ji,jj,jk) = zgraztotf * unass2 + ferat3 * zmortzgoc END DO END DO END DO ! Computation of the effect of DVM by mesozooplankton ! This part is only activated if ln_dvm_meso is set to true ! The parameterization has been published in Gorgues et al. (2019). ! ----------------------------------------------------------------- IF (ln_dvm_meso) THEN ALLOCATE( zgramigrem(jpi,jpj), zgramigref(jpi,jpj), zgramigpoc(jpi,jpj), zgramigpof(jpi,jpj) ) ALLOCATE( zgramigbsi(jpi,jpj) ) ALLOCATE( zstrn(jpi,jpj) ) zgramigrem(:,:) = 0.0 ; zgramigref(:,:) = 0.0 zgramigpoc(:,:) = 0.0 ; zgramigpof(:,:) = 0.0 zgramigbsi(:,:) = 0.0 ! compute the day length depending on latitude and the day zrum = REAL( nday_year - 80, wp ) / REAL( nyear_len(1), wp ) zcodel = ASIN( SIN( zrum * rpi * 2._wp ) * SIN( rad * 23.5_wp ) ) ! day length in hours zstrn(:,:) = 0. DO jj = 1, jpj DO ji = 1, jpi zargu = TAN( zcodel ) * TAN( gphit(ji,jj) * rad ) zargu = MAX( -1., MIN( 1., zargu ) ) zstrn(ji,jj) = MAX( 0.0, 24. - 2. * ACOS( zargu ) / rad / 15. ) zstrn(ji,jj) = MIN(0.75, MAX( 0.25, zstrn(ji,jj) / 24.) ) END DO END DO ! Compute the amount of materials that will go into vertical migration ! This fraction is sumed over the euphotic zone and is removed from ! the fluxes driven by mesozooplankton in the euphotic zone. ! -------------------------------------------------------------------- DO jk = 1, jpk DO jj = 1, jpj DO ji = 1, jpi zmigreltime = (1. - zstrn(ji,jj)) IF ( gdept_n(ji,jj,jk) <= heup(ji,jj) ) THEN zgramigrem(ji,jj) = zgramigrem(ji,jj) + xfracmig * zgrarem(ji,jj,jk) * (1. - zmigreltime ) & & * e3t_n(ji,jj,jk) * tmask(ji,jj,jk) zgramigref(ji,jj) = zgramigref(ji,jj) + xfracmig * zgraref(ji,jj,jk) * (1. - zmigreltime ) & & * e3t_n(ji,jj,jk) * tmask(ji,jj,jk) zgramigpoc(ji,jj) = zgramigpoc(ji,jj) + xfracmig * zgrapoc(ji,jj,jk) * (1. - zmigreltime ) & & * e3t_n(ji,jj,jk) * tmask(ji,jj,jk) zgramigpof(ji,jj) = zgramigpof(ji,jj) + xfracmig * zgrapof(ji,jj,jk) * (1. - zmigreltime ) & & * e3t_n(ji,jj,jk) * tmask(ji,jj,jk) zgramigbsi(ji,jj) = zgramigbsi(ji,jj) + xfracmig * zgrabsi(ji,jj,jk) * (1. - zmigreltime ) & & * e3t_n(ji,jj,jk) * tmask(ji,jj,jk) zgrarem(ji,jj,jk) = zgrarem(ji,jj,jk) * ( (1.0 - xfracmig) + xfracmig * zmigreltime ) zgraref(ji,jj,jk) = zgraref(ji,jj,jk) * ( (1.0 - xfracmig) + xfracmig * zmigreltime ) zgrapoc(ji,jj,jk) = zgrapoc(ji,jj,jk) * ( (1.0 - xfracmig) + xfracmig * zmigreltime ) zgrapof(ji,jj,jk) = zgrapof(ji,jj,jk) * ( (1.0 - xfracmig) + xfracmig * zmigreltime ) zgrabsi(ji,jj,jk) = zgrabsi(ji,jj,jk) * ( (1.0 - xfracmig) + xfracmig * zmigreltime ) ENDIF END DO END DO END DO ! The inorganic and organic fluxes induced by migrating organisms are added at the ! the migration depth (corresponding indice is set by kmig) ! -------------------------------------------------------------------------------- DO jj = 1, jpj DO ji = 1, jpi IF (tmask(ji,jj,1) == 1.) THEN jkt = kmig(ji,jj) zgrarem(ji,jj,jkt) = zgrarem(ji,jj,jkt) + zgramigrem(ji,jj) / e3t_n(ji,jj,jkt) zgraref(ji,jj,jkt) = zgraref(ji,jj,jkt) + zgramigref(ji,jj) / e3t_n(ji,jj,jkt) zgrapoc(ji,jj,jkt) = zgrapoc(ji,jj,jkt) + zgramigpoc(ji,jj) / e3t_n(ji,jj,jkt) zgrapof(ji,jj,jkt) = zgrapof(ji,jj,jkt) + zgramigpof(ji,jj) / e3t_n(ji,jj,jkt) zgrabsi(ji,jj,jkt) = zgrabsi(ji,jj,jkt) + zgramigbsi(ji,jj) / e3t_n(ji,jj,jkt) ENDIF END DO END DO ! ! Deallocate temporary variables ! ------------------------------ DEALLOCATE( zgramigrem, zgramigref, zgramigpoc, zgramigpof, zgramigbsi ) DEALLOCATE( zstrn ) ! End of the ln_dvm_meso part ENDIF DO jk = 1, jpk DO jj = 1, jpj DO ji = 1, jpi ! Update the arrays TRA which contain the biological sources and sinks ! This only concerns the variables which are affected by DVM (inorganic ! nutrients, DOC agands, and particulate organic carbon). zgrarsig = zgrarem(ji,jj,jk) * sigma2 tra(ji,jj,jk,jppo4) = tra(ji,jj,jk,jppo4) + zgrarsig tra(ji,jj,jk,jpnh4) = tra(ji,jj,jk,jpnh4) + zgrarsig tra(ji,jj,jk,jpdoc) = tra(ji,jj,jk,jpdoc) + zgrarem(ji,jj,jk) - zgrarsig ! IF( ln_ligand ) THEN tra(ji,jj,jk,jplgw) = tra(ji,jj,jk,jplgw) + (zgrarem(ji,jj,jk) - zgrarsig) * ldocz zz2ligprod(ji,jj,jk) = (zgrarem(ji,jj,jk) - zgrarsig) * ldocz ENDIF ! tra(ji,jj,jk,jpoxy) = tra(ji,jj,jk,jpoxy) - o2ut * zgrarsig tra(ji,jj,jk,jpfer) = tra(ji,jj,jk,jpfer) + zgraref(ji,jj,jk) zfezoo2(ji,jj,jk) = zgraref(ji,jj,jk) tra(ji,jj,jk,jpdic) = tra(ji,jj,jk,jpdic) + zgrarsig tra(ji,jj,jk,jptal) = tra(ji,jj,jk,jptal) + rno3 * zgrarsig tra(ji,jj,jk,jpgoc) = tra(ji,jj,jk,jpgoc) + zgrapoc(ji,jj,jk) prodgoc(ji,jj,jk) = prodgoc(ji,jj,jk) + zgrapoc(ji,jj,jk) tra(ji,jj,jk,jpbfe) = tra(ji,jj,jk,jpbfe) + zgrapof(ji,jj,jk) tra(ji,jj,jk,jpgsi) = tra(ji,jj,jk,jpgsi) + zgrabsi(ji,jj,jk) END DO END DO END DO ! ! Write the output IF( lk_iomput .AND. knt == nrdttrc ) THEN ALLOCATE( zw3d(jpi,jpj,jpk) ) IF( iom_use( "GRAZ2" ) ) THEN zw3d(:,:,:) = zgrazing(:,:,:) * 1.e+3 * rfact2r * tmask(:,:,:) ! Total grazing of phyto by zooplankton CALL iom_put( "GRAZ2", zw3d ) ENDIF IF( iom_use( "PCAL" ) ) THEN zw3d(:,:,:) = prodcal(:,:,:) * 1.e+3 * rfact2r * tmask(:,:,:) ! Calcite production CALL iom_put( "PCAL", zw3d ) ENDIF IF( iom_use( "FEZOO2" ) ) THEN zw3d(:,:,:) = zfezoo2(:,:,:) * 1e9 * 1.e+3 * rfact2r * tmask(:,:,:) ! CALL iom_put( "FEZOO2", zw3d ) ENDIF IF( iom_use( "LPRODZ2" ) .AND. ln_ligand ) THEN zw3d(:,:,:) = zz2ligprod(:,:,:) * 1e9 * 1.e+3 * rfact2r * tmask(:,:,:) CALL iom_put( "LPRODZ2" , zw3d ) ENDIF DEALLOCATE( zw3d ) ENDIF ! IF (ln_ligand) DEALLOCATE( zz2ligprod ) ! IF(ln_ctl) THEN ! print mean trends (used for debugging) WRITE(charout, FMT="('meso')") CALL prt_ctl_trc_info(charout) CALL prt_ctl_trc(tab4d=tra, mask=tmask, clinfo=ctrcnm) ENDIF ! IF( ln_timing ) CALL timing_stop('p4z_meso') ! END SUBROUTINE p4z_meso SUBROUTINE p4z_meso_init !!---------------------------------------------------------------------- !! *** ROUTINE p4z_meso_init *** !! !! ** Purpose : Initialization of mesozooplankton parameters !! !! ** Method : Read the namp4zmes namelist and check the parameters !! called at the first timestep (nittrc000) !! !! ** input : Namelist nampismes !!---------------------------------------------------------------------- INTEGER :: ios ! Local integer ! NAMELIST/namp4zmes/ part2, grazrat2, resrat2, mzrat2, xpref2n, xpref2d, xpref2z, & & xpref2c, xthresh2dia, xthresh2phy, xthresh2zoo, xthresh2poc, & & xthresh2, xkgraz2, epsher2, epsher2min, sigma2, unass2, grazflux, ln_dvm_meso, & & xsigma2, xsigma2del, xfracmig !!---------------------------------------------------------------------- ! IF(lwp) THEN WRITE(numout,*) WRITE(numout,*) 'p4z_meso_init : Initialization of mesozooplankton parameters' WRITE(numout,*) '~~~~~~~~~~~~~' ENDIF ! REWIND( numnatp_ref ) ! Namelist namp4zmes in reference namelist : Pisces mesozooplankton READ ( numnatp_ref, namp4zmes, IOSTAT = ios, ERR = 901) 901 IF( ios /= 0 ) CALL ctl_nam ( ios , 'namp4zmes in reference namelist' ) REWIND( numnatp_cfg ) ! Namelist namp4zmes in configuration namelist : Pisces mesozooplankton READ ( numnatp_cfg, namp4zmes, IOSTAT = ios, ERR = 902 ) 902 IF( ios > 0 ) CALL ctl_nam ( ios , 'namp4zmes in configuration namelist' ) IF(lwm) WRITE( numonp, namp4zmes ) ! IF(lwp) THEN ! control print WRITE(numout,*) ' Namelist : namp4zmes' WRITE(numout,*) ' part of calcite not dissolved in mesozoo guts part2 =', part2 WRITE(numout,*) ' mesozoo preference for phyto xpref2n =', xpref2n WRITE(numout,*) ' mesozoo preference for diatoms xpref2d =', xpref2d WRITE(numout,*) ' mesozoo preference for zoo xpref2z =', xpref2z WRITE(numout,*) ' mesozoo preference for poc xpref2c =', xpref2c WRITE(numout,*) ' microzoo feeding threshold for mesozoo xthresh2zoo =', xthresh2zoo WRITE(numout,*) ' diatoms feeding threshold for mesozoo xthresh2dia =', xthresh2dia WRITE(numout,*) ' nanophyto feeding threshold for mesozoo xthresh2phy =', xthresh2phy WRITE(numout,*) ' poc feeding threshold for mesozoo xthresh2poc =', xthresh2poc WRITE(numout,*) ' feeding threshold for mesozooplankton xthresh2 =', xthresh2 WRITE(numout,*) ' exsudation rate of mesozooplankton resrat2 =', resrat2 WRITE(numout,*) ' mesozooplankton mortality rate mzrat2 =', mzrat2 WRITE(numout,*) ' maximal mesozoo grazing rate grazrat2 =', grazrat2 WRITE(numout,*) ' mesozoo flux feeding rate grazflux =', grazflux WRITE(numout,*) ' non assimilated fraction of P by mesozoo unass2 =', unass2 WRITE(numout,*) ' Efficiency of Mesozoo growth epsher2 =', epsher2 WRITE(numout,*) ' Minimum Efficiency of Mesozoo growth epsher2min =', epsher2min WRITE(numout,*) ' Fraction of mesozoo excretion as DOM sigma2 =', sigma2 WRITE(numout,*) ' half sturation constant for grazing 2 xkgraz2 =', xkgraz2 WRITE(numout,*) ' Width of the grazing window xsigma2 =', xsigma2 WRITE(numout,*) ' Maximum additional width of the grazing window xsigma2del =', xsigma2del WRITE(numout,*) ' Diurnal vertical migration of mesozoo. ln_dvm_meso =', ln_dvm_meso WRITE(numout,*) ' Fractional biomass of meso that performs DVM xfracmig =', xfracmig ENDIF ! END SUBROUTINE p4z_meso_init SUBROUTINE p4z_meso_depmig !!---------------------------------------------------------------------- !! *** ROUTINE p4z_meso_depmig *** !! !! ** Purpose : Computation the migration depth of mesozooplankton !! !! ** Method : Computes the DVM depth of mesozooplankton from oxygen !! temperature and chlorophylle following the parameterization !! proposed by Bianchi et al. (2013) !!---------------------------------------------------------------------- INTEGER :: ji, jj, jk ! REAL(wp) :: totchl REAL(wp), DIMENSION(jpi,jpj) :: oxymoy, tempmoy, zdepmoy !!--------------------------------------------------------------------- ! IF( ln_timing == 1 ) CALL timing_start('p4z_meso_zdepmig') ! oxymoy(:,:) = 0. tempmoy(:,:) = 0. zdepmoy(:,:) = 0. depmig (:,:) = 5. kmig (:,:) = 1 ! ! Compute the averaged values of oxygen, temperature over the domain ! 150m to 500 m depth. ! ------------------------------------------------------------------ DO jk =1, jpk DO jj = 1, jpj DO ji = 1, jpi IF (tmask(ji,jj,jk) == 1.) THEN IF (gdept_n(ji,jj,jk) >= 150. .AND. gdept_n(ji,jj,jk) <= 500.) THEN oxymoy(ji,jj) = oxymoy(ji,jj) + trb(ji,jj,jk,jpoxy)*e3t_n(ji,jj,jk)*1E6 tempmoy(ji,jj) = tempmoy(ji,jj) + tsn(ji,jj,jk,jp_tem)*e3t_n(ji,jj,jk) zdepmoy(ji,jj) = zdepmoy(ji,jj) + e3t_n(ji,jj,jk) ENDIF ENDIF END DO END DO END DO ! Compute the difference between surface values and the mean values in the mesopelagic ! domain ! ------------------------------------------------------------------------------------ DO jj = 1, jpj DO ji = 1, jpi oxymoy(ji,jj) = trb(ji,jj,1,jpoxy)*1E6 - oxymoy(ji,jj) / (zdepmoy(ji,jj) + rtrn) tempmoy(ji,jj) = tsn(ji,jj,1,jp_tem)-tempmoy(ji,jj) / (zdepmoy(ji,jj) + rtrn) END DO END DO ! ! Computation of the migration depth based on the parameterization of ! Bianchi et al. (2013) ! ------------------------------------------------------------------- DO jj = 1, jpj DO ji = 1, jpi IF (tmask(ji,jj,1) == 1.) THEN totchl = (trb(ji,jj,1,jpnch)+trb(ji,jj,1,jpdch))*1E6 depmig(ji,jj) = 398. - 0.56 * oxymoy(ji,jj) -115. * log10(totchl) + 0.36 * hmld(ji,jj) -2.4 * tempmoy(ji,jj) ENDIF END DO END DO ! ! Computation of the corresponding jk indice ! ------------------------------------------ DO jk = 1, jpk-1 DO jj = 1, jpj DO ji = 1, jpi IF (depmig(ji,jj) .GE. gdepw_n(ji,jj,jk) .AND. depmig(ji,jj) .LT. gdepw_n(ji,jj,jk+1) ) THEN kmig(ji,jj) = jk ENDIF END DO END DO END DO ! ! Correction of the migration depth and indice based on O2 levels ! If O2 is too low, imposing a migration depth at this low O2 levels ! would lead to negative O2 concentrations (respiration while O2 is close ! to 0. Thus, to avoid that problem, the migration depth is adjusted so ! that it falls above the OMZ ! ----------------------------------------------------------------------- DO ji =1, jpi DO jj = 1, jpj IF (trb(ji,jj,kmig(ji,jj),jpoxy) < 5E-6) THEN DO jk = kmig(ji,jj),1,-1 IF (trb(ji,jj,jk,jpoxy) >= 5E-6 .AND. trb(ji,jj,jk+1,jpoxy) < 5E-6) THEN kmig(ji,jj) = jk depmig(ji,jj) = gdept_n(ji,jj,jk) ENDIF END DO ENDIF END DO END DO ! IF( ln_timing ) CALL timing_stop('p4z_meso_depmig') ! END SUBROUTINE p4z_meso_depmig INTEGER FUNCTION p4z_meso_alloc() !!---------------------------------------------------------------------- !! *** ROUTINE p4z_meso_alloc *** !!---------------------------------------------------------------------- ! ALLOCATE( depmig(jpi,jpj), kmig(jpi,jpj), STAT= p4z_meso_alloc ) ! IF( p4z_meso_alloc /= 0 ) CALL ctl_stop( 'STOP', 'p4z_meso_alloc : failed to allocate arrays.' ) ! END FUNCTION p4z_meso_alloc !!====================================================================== END MODULE p4zmeso