[8395] | 1 | MODULE detritus_mod |
---|
| 2 | !!====================================================================== |
---|
| 3 | !! *** MODULE detritus_mod *** |
---|
| 4 | !! Calculates detritus processes and fast-sinking detritus |
---|
| 5 | !!====================================================================== |
---|
| 6 | !! History : |
---|
| 7 | !! - ! 2017-04 (M. Stringer) Code taken from trcbio_medusa.F90 |
---|
| 8 | !!---------------------------------------------------------------------- |
---|
| 9 | #if defined key_medusa |
---|
| 10 | !!---------------------------------------------------------------------- |
---|
| 11 | !! MEDUSA bio-model |
---|
| 12 | !!---------------------------------------------------------------------- |
---|
| 13 | |
---|
| 14 | IMPLICIT NONE |
---|
| 15 | PRIVATE |
---|
| 16 | |
---|
| 17 | PUBLIC detritus ! Called in trcbio_medusa.F90 |
---|
| 18 | |
---|
| 19 | !!---------------------------------------------------------------------- |
---|
| 20 | !! NEMO/TOP 2.0 , LOCEAN-IPSL (2007) |
---|
| 21 | !! $Id$ |
---|
| 22 | !! Software governed by the CeCILL licence (modipsl/doc/NEMO_CeCILL.txt) |
---|
| 23 | !!---------------------------------------------------------------------- |
---|
| 24 | |
---|
| 25 | CONTAINS |
---|
| 26 | |
---|
| 27 | SUBROUTINE detritus( jk, iball ) |
---|
| 28 | !!------------------------------------------------------------------- |
---|
| 29 | !! *** ROUTINE detritus *** |
---|
| 30 | !! This called from TRC_BIO_MEDUSA and |
---|
| 31 | !! - Calculates detritus processes |
---|
| 32 | !! - Fast-sinking detritus |
---|
| 33 | !!------------------------------------------------------------------- |
---|
| 34 | USE bio_medusa_mod, ONLY: f_sbenin_c, f_sbenin_fe, & |
---|
| 35 | f_sbenin_n, fdd, & |
---|
| 36 | idf, idfval, & |
---|
| 37 | # if defined key_roam |
---|
| 38 | fddc, & |
---|
| 39 | # endif |
---|
| 40 | fun_T, fun_Q10, zdet, zdtc |
---|
| 41 | USE detritus_fast_sink_mod, ONLY: detritus_fast_sink |
---|
| 42 | USE dom_oce, ONLY: mbathy, tmask |
---|
| 43 | USE in_out_manager, ONLY: lwp, numout |
---|
| 44 | USE par_oce, ONLY: jpim1, jpjm1 |
---|
| 45 | USE sms_medusa, ONLY: jmd, jorgben, jsfd, vsed, & |
---|
| 46 | xrfn, xmd, xmdc, xthetad |
---|
| 47 | |
---|
| 48 | !! Level |
---|
| 49 | INTEGER, INTENT( in ) :: jk |
---|
| 50 | !! Fast detritus ballast scheme (0 = no; 1 = yes) |
---|
| 51 | INTEGER, INTENT( in ) :: iball |
---|
| 52 | |
---|
| 53 | INTEGER :: ji, jj |
---|
| 54 | |
---|
| 55 | !!------------------------------------------------------------------ |
---|
| 56 | !! Detritus remineralisation |
---|
| 57 | !! Constant or temperature-dependent |
---|
| 58 | !!------------------------------------------------------------------ |
---|
| 59 | DO jj = 2,jpjm1 |
---|
| 60 | DO ji = 2,jpim1 |
---|
| 61 | !! OPEN wet point IF..THEN loop |
---|
| 62 | if (tmask(ji,jj,jk) == 1) then |
---|
| 63 | !! |
---|
| 64 | if (jmd.eq.1) then |
---|
| 65 | !! temperature-dependent |
---|
| 66 | fdd(ji,jj) = xmd * fun_T(ji,jj) * zdet(ji,jj) |
---|
| 67 | # if defined key_roam |
---|
| 68 | fddc(ji,jj) = xmdc * fun_T(ji,jj) * zdtc(ji,jj) |
---|
| 69 | # endif |
---|
| 70 | elseif (jmd.eq.2) then |
---|
| 71 | !! AXY (16/05/13): add in Q10-based parameterisation |
---|
| 72 | !! (def in nmlst) |
---|
| 73 | !! temperature-dependent |
---|
| 74 | fdd(ji,jj) = xmd * fun_Q10(ji,jj) * zdet(ji,jj) |
---|
| 75 | #if defined key_roam |
---|
| 76 | fddc(ji,jj) = xmdc * fun_Q10(ji,jj) * zdtc(ji,jj) |
---|
| 77 | #endif |
---|
| 78 | else |
---|
| 79 | !! temperature-independent |
---|
| 80 | fdd(ji,jj) = xmd * zdet(ji,jj) |
---|
| 81 | # if defined key_roam |
---|
| 82 | fddc(ji,jj) = xmdc * zdtc(ji,jj) |
---|
| 83 | # endif |
---|
| 84 | endif |
---|
| 85 | !! |
---|
| 86 | !! AXY (22/07/09): accelerate detrital remineralisation |
---|
| 87 | !! in the bottom box |
---|
| 88 | if ((jk.eq.mbathy(ji,jj)) .and. jsfd.eq.1) then |
---|
| 89 | fdd(ji,jj) = 1.0 * zdet(ji,jj) |
---|
| 90 | # if defined key_roam |
---|
| 91 | fddc(ji,jj) = 1.0 * zdtc(ji,jj) |
---|
| 92 | # endif |
---|
| 93 | endif |
---|
| 94 | |
---|
| 95 | # if defined key_debug_medusa |
---|
| 96 | !! report plankton mortality and remineralisation |
---|
| 97 | if (idf.eq.1.AND.idfval.eq.1) then |
---|
| 98 | IF (lwp) write (numout,*) '------------------------------' |
---|
| 99 | ! I've removed the lines below, because the variables are not in this |
---|
| 100 | ! routine. If these debug prints need to stay, they should probably be |
---|
| 101 | ! moved - marc 27/4/17 |
---|
| 102 | ! IF (lwp) write (numout,*) 'fdpn2(',jk,') = ', fdpn2(ji,jj) |
---|
| 103 | ! IF (lwp) write (numout,*) 'fdpd2(',jk,') = ', fdpd2(ji,jj) |
---|
| 104 | ! IF (lwp) write (numout,*) 'fdpds2(',jk,')= ', fdpds2(ji,jj) |
---|
| 105 | ! IF (lwp) write (numout,*) 'fdzmi2(',jk,')= ', fdzmi2(ji,jj) |
---|
| 106 | ! IF (lwp) write (numout,*) 'fdzme2(',jk,')= ', fdzme2(ji,jj) |
---|
| 107 | ! IF (lwp) write (numout,*) 'fdpn(',jk,') = ', fdpn(ji,jj) |
---|
| 108 | ! IF (lwp) write (numout,*) 'fdpd(',jk,') = ', fdpd(ji,jj) |
---|
| 109 | ! IF (lwp) write (numout,*) 'fdpds(',jk,') = ', fdpds(ji,jj) |
---|
| 110 | ! IF (lwp) write (numout,*) 'fdzmi(',jk,') = ', fdzmi(ji,jj) |
---|
| 111 | ! IF (lwp) write (numout,*) 'fdzme(',jk,') = ', fdzme(ji,jj) |
---|
| 112 | IF (lwp) write (numout,*) 'fdd(',jk,') = ', fdd(ji,jj) |
---|
| 113 | # if defined key_roam |
---|
| 114 | IF (lwp) write (numout,*) 'fddc(',jk,') = ', fddc(ji,jj) |
---|
| 115 | # endif |
---|
| 116 | endif |
---|
| 117 | # endif |
---|
| 118 | ENDIF |
---|
| 119 | ENDDO |
---|
| 120 | ENDDO |
---|
| 121 | |
---|
| 122 | DO jj = 2,jpjm1 |
---|
| 123 | DO ji = 2,jpim1 |
---|
| 124 | if (tmask(ji,jj,jk) == 1) then |
---|
| 125 | !!--------------------------------------------------------- |
---|
| 126 | !! Detritus addition to benthos |
---|
| 127 | !! If activated, slow detritus in the bottom box will enter |
---|
| 128 | !! the benthic pool |
---|
| 129 | !!--------------------------------------------------------- |
---|
| 130 | !! |
---|
| 131 | if ((jk.eq.mbathy(ji,jj)) .and. jorgben.eq.1) then |
---|
| 132 | !! this is the BOTTOM OCEAN BOX -> into the benthic pool! |
---|
| 133 | !! |
---|
| 134 | f_sbenin_n(ji,jj) = (zdet(ji,jj) * vsed * 86400.) |
---|
| 135 | f_sbenin_fe(ji,jj) = (zdet(ji,jj) * vsed * 86400. * xrfn) |
---|
| 136 | # if defined key_roam |
---|
| 137 | f_sbenin_c(ji,jj) = (zdtc(ji,jj) * vsed * 86400.) |
---|
| 138 | # else |
---|
| 139 | f_sbenin_c(ji,jj) = (zdet(ji,jj) * vsed * 86400. * xthetad) |
---|
| 140 | # endif |
---|
| 141 | endif |
---|
| 142 | |
---|
| 143 | ENDIF |
---|
| 144 | ENDDO |
---|
| 145 | ENDDO |
---|
| 146 | |
---|
| 147 | !!------------------------------------------------------------------ |
---|
| 148 | !! Fast-sinking detritus |
---|
| 149 | !!------------------------------------------------------------------ |
---|
| 150 | CALL detritus_fast_sink( jk, iball ) |
---|
| 151 | |
---|
| 152 | END SUBROUTINE detritus |
---|
| 153 | |
---|
| 154 | #else |
---|
| 155 | !!====================================================================== |
---|
| 156 | !! Dummy module : No MEDUSA bio-model |
---|
| 157 | !!====================================================================== |
---|
| 158 | CONTAINS |
---|
| 159 | SUBROUTINE detritus( ) ! Empty routine |
---|
| 160 | WRITE(*,*) 'detritus: You should not have seen this print! error?' |
---|
| 161 | END SUBROUTINE detritus |
---|
| 162 | #endif |
---|
| 163 | |
---|
| 164 | !!====================================================================== |
---|
| 165 | END MODULE detritus_mod |
---|