New URL for NEMO forge!   http://forge.nemo-ocean.eu

Since March 2022 along with NEMO 4.2 release, the code development moved to a self-hosted GitLab.
This present forge is now archived and remained online for history.
p4zsbc.F90 in NEMO/releases/release-4.0/src/TOP/PISCES/P4Z – NEMO

source: NEMO/releases/release-4.0/src/TOP/PISCES/P4Z/p4zsbc.F90 @ 10868

Last change on this file since 10868 was 10868, checked in by cetlod, 5 years ago

v4.0:Change units of Nitrogen deposition data in PISCES, see ticket:2271

  • Property svn:keywords set to Id
File size: 25.8 KB
RevLine 
[3443]1MODULE p4zsbc
2   !!======================================================================
3   !!                         ***  MODULE p4sbc  ***
4   !! TOP :   PISCES surface boundary conditions of external inputs of nutrients
5   !!======================================================================
6   !! History :   3.5  !  2012-07 (O. Aumont, C. Ethe) Original code
7   !!----------------------------------------------------------------------
8   !!   p4z_sbc        :  Read and interpolate time-varying nutrients fluxes
9   !!   p4z_sbc_init   :  Initialization of p4z_sbc
10   !!----------------------------------------------------------------------
11   USE oce_trc         !  shared variables between ocean and passive tracers
12   USE trc             !  passive tracers common variables
13   USE sms_pisces      !  PISCES Source Minus Sink variables
14   USE iom             !  I/O manager
15   USE fldread         !  time interpolation
16
17   IMPLICIT NONE
18   PRIVATE
19
20   PUBLIC   p4z_sbc
21   PUBLIC   p4z_sbc_init   
22
[9169]23   LOGICAL , PUBLIC ::   ln_dust      !: boolean for dust input from the atmosphere
24   LOGICAL , PUBLIC ::   ln_solub     !: boolean for variable solubility of atmospheric iron
25   LOGICAL , PUBLIC ::   ln_river     !: boolean for river input of nutrients
26   LOGICAL , PUBLIC ::   ln_ndepo     !: boolean for atmospheric deposition of N
27   LOGICAL , PUBLIC ::   ln_ironsed   !: boolean for Fe input from sediments
28   LOGICAL , PUBLIC ::   ln_hydrofe   !: boolean for Fe input from hydrothermal vents
29   REAL(wp), PUBLIC ::   sedfeinput   !: Coastal release of Iron
30   REAL(wp), PUBLIC ::   dustsolub    !: Solubility of the dust
31   REAL(wp), PUBLIC ::   mfrac        !: Mineral Content of the dust
32   REAL(wp), PUBLIC ::   icefeinput   !: Iron concentration in sea ice
33   REAL(wp), PUBLIC ::   wdust        !: Sinking speed of the dust
34   REAL(wp), PUBLIC ::   nitrfix      !: Nitrogen fixation rate   
35   REAL(wp), PUBLIC ::   diazolight   !: Nitrogen fixation sensitivty to light
36   REAL(wp), PUBLIC ::   concfediaz   !: Fe half-saturation Cste for diazotrophs
37   REAL(wp)         ::   hratio       !: Fe:3He ratio assumed for vent iron supply
[10111]38   REAL(wp)         ::   distcoast    !: Distance off the coast for Iron from sediments
[9169]39   REAL(wp), PUBLIC ::   lgw_rath     !: Weak ligand ratio from hydro sources
[3443]40
[9169]41   LOGICAL , PUBLIC ::   ll_sbc
42   LOGICAL          ::   ll_solub
[7646]43
[3443]44   INTEGER , PARAMETER  :: jpriv  = 7   !: Maximum number of river input fields
45   INTEGER , PARAMETER  :: jr_dic = 1   !: index of dissolved inorganic carbon
46   INTEGER , PARAMETER  :: jr_doc = 2   !: index of dissolved organic carbon
47   INTEGER , PARAMETER  :: jr_din = 3   !: index of dissolved inorganic nitrogen
48   INTEGER , PARAMETER  :: jr_don = 4   !: index of dissolved organic nitrogen
49   INTEGER , PARAMETER  :: jr_dip = 5   !: index of dissolved inorganic phosporus
50   INTEGER , PARAMETER  :: jr_dop = 6   !: index of dissolved organic phosphorus
51   INTEGER , PARAMETER  :: jr_dsi = 7   !: index of dissolved silicate
52
53   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf_dust      ! structure of input dust
[9169]54   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf_solub     ! structure of input dust
55   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf_river     ! structure of input riverdic
[3443]56   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf_ndepo     ! structure of input nitrogen deposition
57   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf_ironsed   ! structure of input iron from sediment
58   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf_hydrofe   ! structure of input iron from hydrothermal vents
59
[9169]60   INTEGER , PARAMETER ::   nbtimes = 365                          ! maximum number of times record in a file
61   INTEGER             ::   ntimes_dust, ntimes_riv, ntimes_ndep   ! number of time steps in a file
62   INTEGER             ::   ntimes_solub, ntimes_hydro             ! number of time steps in a file
[3443]63
[9169]64   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   dust  , solub    !: dust fields
65   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   rivdic, rivalk   !: river input fields
66   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   rivdin, rivdip   !: river input fields
67   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   rivdon, rivdop   !: river input fields
68   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   rivdoc           !: river input fields
69   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   rivdsi           !: river input fields
70   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:)   ::   nitdep           !: atmospheric N deposition
71   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   ironsed          !: Coastal supply of iron
72   REAL(wp), PUBLIC, ALLOCATABLE, SAVE, DIMENSION(:,:,:) ::   hydrofe          !: Hydrothermal vent supply of iron
[3443]73
[10127]74   REAL(wp), PUBLIC :: sedsilfrac, sedcalfrac
75   REAL(wp), PUBLIC :: rivalkinput, rivdicinput
76   REAL(wp), PUBLIC :: rivdininput, rivdipinput, rivdsiinput
[3443]77
[5836]78   !! * Substitutions
79#  include "vectopt_loop_substitute.h90"
[3443]80   !!----------------------------------------------------------------------
[10067]81   !! NEMO/TOP 4.0 , NEMO Consortium (2018)
[5215]82   !! $Id$
[10068]83   !! Software governed by the CeCILL license (see ./LICENSE)
[3443]84   !!----------------------------------------------------------------------
85CONTAINS
86
87   SUBROUTINE p4z_sbc( kt )
88      !!----------------------------------------------------------------------
89      !!                  ***  routine p4z_sbc  ***
90      !!
91      !! ** purpose :   read and interpolate the external sources of nutrients
92      !!
93      !! ** method  :   read the files and interpolate the appropriate variables
94      !!
95      !! ** input   :   external netcdf files
96      !!
97      !!----------------------------------------------------------------------
[9169]98      INTEGER, INTENT(in) ::   kt   ! ocean time step
99      !
100      INTEGER  ::   ji, jj 
101      REAL(wp) ::   zcoef, zyyss
[3443]102      !!---------------------------------------------------------------------
103      !
[9169]104      IF( ln_timing )   CALL timing_start('p4z_sbc')
[3443]105      !
106      ! Compute dust at nit000 or only if there is more than 1 time record in dust file
107      IF( ln_dust ) THEN
108         IF( kt == nit000 .OR. ( kt /= nit000 .AND. ntimes_dust > 1 ) ) THEN
109            CALL fld_read( kt, 1, sf_dust )
[10791]110            dust(:,:) = MAX( rtrn, sf_dust(1)%fnow(:,:,1) ) * ( 1.0 - fr_i(:,:) )
[3443]111         ENDIF
112      ENDIF
[9169]113      !
[3443]114      IF( ll_solub ) THEN
115         IF( kt == nit000 .OR. ( kt /= nit000 .AND. ntimes_solub > 1 ) ) THEN
116            CALL fld_read( kt, 1, sf_solub )
[7753]117            solub(:,:) = sf_solub(1)%fnow(:,:,1)
[3443]118         ENDIF
119      ENDIF
120
121      ! N/P and Si releases due to coastal rivers
122      ! Compute river at nit000 or only if there is more than 1 time record in river file
123      ! -----------------------------------------
124      IF( ln_river ) THEN
125         IF( kt == nit000 .OR. ( kt /= nit000 .AND. ntimes_riv > 1 ) ) THEN
126            CALL fld_read( kt, 1, sf_river )
[7646]127            IF( ln_p4z ) THEN
128               DO jj = 1, jpj
129                  DO ji = 1, jpi
130                     zcoef = ryyss * e1e2t(ji,jj) * h_rnf(ji,jj) 
[10362]131                     rivalk(ji,jj) =   sf_river(jr_dic)%fnow(ji,jj,1)  &
[7646]132                        &              * 1.E3        / ( 12. * zcoef + rtrn )
[10362]133                     rivdic(ji,jj) =   sf_river(jr_dic)%fnow(ji,jj,1)  &
[7646]134                        &              * 1.E3         / ( 12. * zcoef + rtrn )
[10362]135                     rivdin(ji,jj) =   sf_river(jr_din)%fnow(ji,jj,1)  &
[7646]136                        &              * 1.E3 / rno3 / ( 14. * zcoef + rtrn )
[10362]137                     rivdip(ji,jj) =   sf_river(jr_dip)%fnow(ji,jj,1)  &
[7646]138                        &              * 1.E3 / po4r / ( 31. * zcoef + rtrn )
[10362]139                     rivdsi(ji,jj) =   sf_river(jr_dsi)%fnow(ji,jj,1)  &
[7646]140                        &              * 1.E3        / ( 28.1 * zcoef + rtrn )
[10362]141                     rivdoc(ji,jj) =   sf_river(jr_doc)%fnow(ji,jj,1)  &
142                        &              * 1.E3        / ( 12. * zcoef + rtrn ) 
[7646]143                  END DO
[3443]144               END DO
[7646]145            ELSE    !  ln_p5z
146               DO jj = 1, jpj
147                  DO ji = 1, jpi
148                     zcoef = ryyss * e1e2t(ji,jj) * h_rnf(ji,jj) 
149                     rivalk(ji,jj) =   sf_river(jr_dic)%fnow(ji,jj,1)                                    &
150                        &              * 1.E3        / ( 12. * zcoef + rtrn )
151                     rivdic(ji,jj) = ( sf_river(jr_dic)%fnow(ji,jj,1) ) &
152                        &              * 1.E3 / ( 12. * zcoef + rtrn ) * tmask(ji,jj,1)
153                     rivdin(ji,jj) = ( sf_river(jr_din)%fnow(ji,jj,1) ) &
154                        &              * 1.E3 / rno3 / ( 14. * zcoef + rtrn ) * tmask(ji,jj,1)
155                     rivdip(ji,jj) = ( sf_river(jr_dip)%fnow(ji,jj,1) ) &
156                        &              * 1.E3 / po4r / ( 31. * zcoef + rtrn ) * tmask(ji,jj,1)
157                     rivdon(ji,jj) = ( sf_river(jr_don)%fnow(ji,jj,1) ) &
158                        &              * 1.E3 / rno3 / ( 14. * zcoef + rtrn ) * tmask(ji,jj,1)
159                     rivdop(ji,jj) = ( sf_river(jr_dop)%fnow(ji,jj,1) ) &
160                        &              * 1.E3 / po4r / ( 31. * zcoef + rtrn ) * tmask(ji,jj,1)
[10362]161                     rivdsi(ji,jj) =   sf_river(jr_dsi)%fnow(ji,jj,1)  &
162                        &              * 1.E3        / ( 28.1 * zcoef + rtrn )
163                     rivdoc(ji,jj) =   sf_river(jr_doc)%fnow(ji,jj,1)  &
164                        &              * 1.E3        / ( 12. * zcoef + rtrn )
[7646]165                  END DO
166               END DO
167            ENDIF
[3443]168         ENDIF
169      ENDIF
170
171      ! Compute N deposition at nit000 or only if there is more than 1 time record in N deposition file
172      IF( ln_ndepo ) THEN
173         IF( kt == nit000 .OR. ( kt /= nit000 .AND. ntimes_ndep > 1 ) ) THEN
[10868]174             zcoef = 14. * rno3
[7753]175             CALL fld_read( kt, 1, sf_ndepo )
[10127]176             nitdep(:,:) = MAX( rtrn, sf_ndepo(1)%fnow(:,:,1) / zcoef / e3t_n(:,:,1) )
[3443]177         ENDIF
[6962]178         IF( .NOT.ln_linssh ) THEN
[10868]179            zcoef = 14. * rno3
180            nitdep(:,:) = MAX( rtrn, sf_ndepo(1)%fnow(:,:,1) / zcoef / e3t_n(:,:,1) )
[6945]181         ENDIF
[3443]182      ENDIF
183      !
[9124]184      IF( ln_timing )  CALL timing_stop('p4z_sbc')
[3443]185      !
186   END SUBROUTINE p4z_sbc
187
[9124]188
[3443]189   SUBROUTINE p4z_sbc_init
190      !!----------------------------------------------------------------------
191      !!                  ***  routine p4z_sbc_init  ***
192      !!
193      !! ** purpose :   initialization of the external sources of nutrients
194      !!
195      !! ** method  :   read the files and compute the budget
196      !!                called at the first timestep (nittrc000)
197      !!
198      !! ** input   :   external netcdf files
199      !!
200      !!----------------------------------------------------------------------
201      INTEGER  :: ji, jj, jk, jm, ifpr
[3446]202      INTEGER  :: ii0, ii1, ij0, ij1
[3443]203      INTEGER  :: numdust, numsolub, numriv, numiron, numdepo, numhydro
204      INTEGER  :: ierr, ierr1, ierr2, ierr3
[4147]205      INTEGER  :: ios                 ! Local integer output status for namelist read
[5385]206      INTEGER  :: ik50                !  last level where depth less than 50 m
207      INTEGER  :: isrow             ! index for ORCA1 starting row
[10111]208      REAL(wp) :: zexpide, zdenitide, zmaskt, zsurfc, zsurfp,ze3t, ze3t2, zcslp
[3443]209      REAL(wp) :: ztimes_dust, ztimes_riv, ztimes_ndep 
210      REAL(wp), DIMENSION(:), ALLOCATABLE :: rivinput
[10127]211      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE :: zriver, zcmask
[3443]212      !
213      CHARACTER(len=100) ::  cn_dir          ! Root directory for location of ssr files
214      TYPE(FLD_N), DIMENSION(jpriv) ::  slf_river    ! array of namelist informations on the fields to read
215      TYPE(FLD_N) ::   sn_dust, sn_solub, sn_ndepo, sn_ironsed, sn_hydrofe   ! informations about the fields to be read
216      TYPE(FLD_N) ::   sn_riverdoc, sn_riverdic, sn_riverdsi   ! informations about the fields to be read
217      TYPE(FLD_N) ::   sn_riverdin, sn_riverdon, sn_riverdip, sn_riverdop
[9169]218      !!
[3443]219      NAMELIST/nampissbc/cn_dir, sn_dust, sn_solub, sn_riverdic, sn_riverdoc, sn_riverdin, sn_riverdon,     &
220        &                sn_riverdip, sn_riverdop, sn_riverdsi, sn_ndepo, sn_ironsed, sn_hydrofe, &
221        &                ln_dust, ln_solub, ln_river, ln_ndepo, ln_ironsed, ln_ironice, ln_hydrofe,    &
[10111]222        &                sedfeinput, distcoast, dustsolub, icefeinput, wdust, mfrac, nitrfix, diazolight, concfediaz, &
[10416]223        &                hratio, lgw_rath
[3443]224      !!----------------------------------------------------------------------
225      !
[9169]226      IF(lwp) THEN
227         WRITE(numout,*)
228         WRITE(numout,*) 'p4z_sbc_init : initialization of the external sources of nutrients '
229         WRITE(numout,*) '~~~~~~~~~~~~ '
230      ENDIF
[3443]231      !                            !* set file information
[4147]232      REWIND( numnatp_ref )              ! Namelist nampissbc in reference namelist : Pisces external sources of nutrients
233      READ  ( numnatp_ref, nampissbc, IOSTAT = ios, ERR = 901)
[9169]234901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'nampissbc in reference namelist', lwp )
[4147]235      REWIND( numnatp_cfg )              ! Namelist nampissbc in configuration namelist : Pisces external sources of nutrients
236      READ  ( numnatp_cfg, nampissbc, IOSTAT = ios, ERR = 902 )
[9169]237902   IF( ios >  0 )   CALL ctl_nam ( ios , 'nampissbc in configuration namelist', lwp )
[4624]238      IF(lwm) WRITE ( numonp, nampissbc )
[4147]239
[3443]240      IF(lwp) THEN
[9169]241         WRITE(numout,*) '   Namelist : nampissbc '
242         WRITE(numout,*) '      dust input from the atmosphere           ln_dust     = ', ln_dust
243         WRITE(numout,*) '      Variable solubility of iron input        ln_solub    = ', ln_solub
244         WRITE(numout,*) '      river input of nutrients                 ln_river    = ', ln_river
245         WRITE(numout,*) '      atmospheric deposition of n              ln_ndepo    = ', ln_ndepo
246         WRITE(numout,*) '      Fe input from sediments                  ln_ironsed  = ', ln_ironsed
247         WRITE(numout,*) '      Fe input from seaice                     ln_ironice  = ', ln_ironice
248         WRITE(numout,*) '      fe input from hydrothermal vents         ln_hydrofe  = ', ln_hydrofe
249         WRITE(numout,*) '      coastal release of iron                  sedfeinput  = ', sedfeinput
[10111]250         WRITE(numout,*) '      distance off the coast                   distcoast   = ', distcoast
[9169]251         WRITE(numout,*) '      solubility of the dust                   dustsolub   = ', dustsolub
252         WRITE(numout,*) '      Mineral Fe content of the dust           mfrac       = ', mfrac
253         WRITE(numout,*) '      Iron concentration in sea ice            icefeinput  = ', icefeinput
254         WRITE(numout,*) '      sinking speed of the dust                wdust       = ', wdust
255         WRITE(numout,*) '      nitrogen fixation rate                   nitrfix     = ', nitrfix
256         WRITE(numout,*) '      nitrogen fixation sensitivty to light    diazolight  = ', diazolight
257         WRITE(numout,*) '      Fe half-saturation cste for diazotrophs  concfediaz  = ', concfediaz
258         WRITE(numout,*) '      Fe to 3He ratio assumed for vent iron supply hratio  = ', hratio
[7646]259         IF( ln_ligand ) THEN
[9169]260            WRITE(numout,*) '      Weak ligand ratio from sed hydro sources  lgw_rath   = ', lgw_rath
[7646]261         ENDIF
[3443]262      END IF
263
[9169]264      IF( ln_dust .OR. ln_river .OR. ln_ndepo ) THEN   ;   ll_sbc = .TRUE.
265      ELSE                                             ;   ll_sbc = .FALSE.
[3443]266      ENDIF
267
[9169]268      IF( ln_dust .AND. ln_solub ) THEN                ;   ll_solub = .TRUE.
269      ELSE                                             ;   ll_solub = .FALSE.
270      ENDIF
271
[5385]272      ! set the number of level over which river runoffs are applied
273      ! online configuration : computed in sbcrnf
[7646]274      IF( l_offline ) THEN
[7753]275        nk_rnf(:,:) = 1
276        h_rnf (:,:) = gdept_n(:,:,1)
[5385]277      ENDIF
278
[3443]279      ! dust input from the atmosphere
280      ! ------------------------------
281      IF( ln_dust ) THEN 
282         !
283         IF(lwp) WRITE(numout,*) '    initialize dust input from atmosphere '
284         IF(lwp) WRITE(numout,*) '    ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~ '
285         !
286         ALLOCATE( dust(jpi,jpj) )    ! allocation
287         !
288         ALLOCATE( sf_dust(1), STAT=ierr )           !* allocate and fill sf_sst (forcing structure) with sn_sst
289         IF( ierr > 0 )   CALL ctl_stop( 'STOP', 'p4z_sed_init: unable to allocate sf_dust structure' )
290         !
291         CALL fld_fill( sf_dust, (/ sn_dust /), cn_dir, 'p4z_sed_init', 'Atmospheric dust deposition', 'nampissed' )
292                                   ALLOCATE( sf_dust(1)%fnow(jpi,jpj,1)   )
293         IF( sn_dust%ln_tint )     ALLOCATE( sf_dust(1)%fdta(jpi,jpj,1,2) )
294         !
[3557]295         IF( Agrif_Root() ) THEN   !  Only on the master grid
296            ! Get total input dust ; need to compute total atmospheric supply of Si in a year
297            CALL iom_open (  TRIM( sn_dust%clname ) , numdust )
[10522]298            ntimes_dust = iom_getszuld( numdust )   ! get number of record in file
[10127]299         END IF
[3443]300      END IF
301
302      ! Solubility of dust deposition of iron
303      ! Only if ln_dust and ln_solubility set to true (ll_solub = .true.)
304      ! -----------------------------------------------------------------
305      IF( ll_solub ) THEN
306         !
[9169]307         IF(lwp) WRITE(numout,*)
308         IF(lwp) WRITE(numout,*) '   ==>>>   ll_solub=T , initialize variable solubility of Fe '
[3443]309         !
310         ALLOCATE( solub(jpi,jpj) )    ! allocation
311         !
312         ALLOCATE( sf_solub(1), STAT=ierr )           !* allocate and fill sf_sst (forcing structure) with sn_sst
313         IF( ierr > 0 )   CALL ctl_stop( 'STOP', 'p4z_sed_init: unable to allocate sf_solub structure' )
314         !
315         CALL fld_fill( sf_solub, (/ sn_solub /), cn_dir, 'p4z_sed_init', 'Solubility of atm. iron ', 'nampissed' )
316                                   ALLOCATE( sf_solub(1)%fnow(jpi,jpj,1)   )
317         IF( sn_solub%ln_tint )    ALLOCATE( sf_solub(1)%fdta(jpi,jpj,1,2) )
318         ! get number of record in file
319         CALL iom_open (  TRIM( sn_solub%clname ) , numsolub )
[10522]320         ntimes_solub = iom_getszuld( numsolub )   ! get number of record in file
[3443]321         CALL iom_close( numsolub )
322      ENDIF
323
324      ! nutrient input from rivers
325      ! --------------------------
326      IF( ln_river ) THEN
327         !
[9169]328         slf_river(jr_dic) = sn_riverdic   ;   slf_river(jr_doc) = sn_riverdoc   ;   slf_river(jr_din) = sn_riverdin 
329         slf_river(jr_don) = sn_riverdon   ;   slf_river(jr_dip) = sn_riverdip   ;   slf_river(jr_dop) = sn_riverdop
[3443]330         slf_river(jr_dsi) = sn_riverdsi 
331         !
[10362]332         ALLOCATE( rivdic(jpi,jpj), rivalk(jpi,jpj), rivdin(jpi,jpj), rivdip(jpi,jpj), rivdsi(jpi,jpj), rivdoc(jpi,jpj) ) 
333         IF( ln_p5z )  ALLOCATE( rivdon(jpi,jpj), rivdop(jpi,jpj) )
[3443]334         !
[9169]335         ALLOCATE( sf_river(jpriv), rivinput(jpriv), STAT=ierr1 )    !* allocate and fill sf_river (forcing structure) with sn_river_
336         rivinput(:) = 0._wp
[3443]337
338         IF( ierr1 > 0 )   CALL ctl_stop( 'STOP', 'p4z_sed_init: unable to allocate sf_irver structure' )
339         !
340         CALL fld_fill( sf_river, slf_river, cn_dir, 'p4z_sed_init', 'Input from river ', 'nampissed' )
341         DO ifpr = 1, jpriv
342                                          ALLOCATE( sf_river(ifpr)%fnow(jpi,jpj,1  ) )
343            IF( slf_river(ifpr)%ln_tint ) ALLOCATE( sf_river(ifpr)%fdta(jpi,jpj,1,2) )
344         END DO
[3557]345         IF( Agrif_Root() ) THEN   !  Only on the master grid
346            ! Get total input rivers ; need to compute total river supply in a year
347            DO ifpr = 1, jpriv
348               CALL iom_open ( TRIM( slf_river(ifpr)%clname ), numriv )
[10522]349               ntimes_riv = iom_getszuld( numriv )
[3557]350               ALLOCATE( zriver(jpi,jpj,ntimes_riv) )
351               DO jm = 1, ntimes_riv
352                  CALL iom_get( numriv, jpdom_data, TRIM( slf_river(ifpr)%clvar ), zriver(:,:,jm), jm )
353               END DO
354               CALL iom_close( numriv )
[7646]355               ztimes_riv = 1._wp / REAL(ntimes_riv, wp) 
[3557]356               DO jm = 1, ntimes_riv
[10425]357                  rivinput(ifpr) = rivinput(ifpr) + glob_sum( 'p4zsbc', zriver(:,:,jm) * tmask(:,:,1) * ztimes_riv ) 
[3557]358               END DO
359               DEALLOCATE( zriver)
[3443]360            END DO
[3557]361            ! N/P and Si releases due to coastal rivers
362            ! -----------------------------------------
363            rivdicinput = (rivinput(jr_dic) + rivinput(jr_doc) ) * 1E3 / 12._wp
364            rivdininput = (rivinput(jr_din) + rivinput(jr_don) ) * 1E3 / rno3 / 14._wp
365            rivdipinput = (rivinput(jr_dip) + rivinput(jr_dop) ) * 1E3 / po4r / 31._wp
366            rivdsiinput = rivinput(jr_dsi) * 1E3 / 28.1_wp
367            rivalkinput = rivinput(jr_dic) * 1E3 / 12._wp
368            !
369         ENDIF
[3443]370      ELSE
371         rivdicinput = 0._wp
372         rivdininput = 0._wp
373         rivdipinput = 0._wp
374         rivdsiinput = 0._wp
375         rivalkinput = 0._wp
376      END IF 
377      ! nutrient input from dust
378      ! ------------------------
379      IF( ln_ndepo ) THEN
380         !
[9169]381         IF(lwp) WRITE(numout,*)
382         IF(lwp) WRITE(numout,*) '   ==>>>   ln_ndepo=T , initialize the nutrient input by dust from NetCDF file'
[3443]383         !
384         ALLOCATE( nitdep(jpi,jpj) )    ! allocation
385         !
386         ALLOCATE( sf_ndepo(1), STAT=ierr3 )           !* allocate and fill sf_sst (forcing structure) with sn_sst
387         IF( ierr3 > 0 )   CALL ctl_stop( 'STOP', 'p4z_sed_init: unable to allocate sf_ndepo structure' )
388         !
389         CALL fld_fill( sf_ndepo, (/ sn_ndepo /), cn_dir, 'p4z_sed_init', 'Nutrient atmospheric depositon ', 'nampissed' )
390                                   ALLOCATE( sf_ndepo(1)%fnow(jpi,jpj,1)   )
391         IF( sn_ndepo%ln_tint )    ALLOCATE( sf_ndepo(1)%fdta(jpi,jpj,1,2) )
392         !
[3557]393         IF( Agrif_Root() ) THEN   !  Only on the master grid
394            ! Get total input dust ; need to compute total atmospheric supply of N in a year
395            CALL iom_open ( TRIM( sn_ndepo%clname ), numdepo )
[10522]396            ntimes_ndep = iom_getszuld( numdepo )
[3557]397         ENDIF
[3443]398      ENDIF
399
400      ! coastal and island masks
401      ! ------------------------
402      IF( ln_ironsed ) THEN     
403         !
[9169]404         IF(lwp) WRITE(numout,*)
405         IF(lwp) WRITE(numout,*) '   ==>>>   ln_ironsed=T , computation of an island mask to enhance coastal supply of iron'
[3443]406         !
407         ALLOCATE( ironsed(jpi,jpj,jpk) )    ! allocation
408         !
409         CALL iom_open ( TRIM( sn_ironsed%clname ), numiron )
410         ALLOCATE( zcmask(jpi,jpj,jpk) )
411         CALL iom_get  ( numiron, jpdom_data, TRIM( sn_ironsed%clvar ), zcmask(:,:,:), 1 )
412         CALL iom_close( numiron )
413         !
[5385]414         ik50 = 5        !  last level where depth less than 50 m
415         DO jk = jpkm1, 1, -1
[9169]416            IF( gdept_1d(jk) > 50. )   ik50 = jk - 1
[5385]417         END DO
[9169]418         IF(lwp) WRITE(numout,*)
419         IF(lwp) WRITE(numout,*) ' Level corresponding to 50m depth ',  ik50,' ', gdept_1d(ik50+1)
[5385]420         DO jk = 1, ik50
[3443]421            DO jj = 2, jpjm1
422               DO ji = fs_2, fs_jpim1
[10111]423                  ze3t   = e3t_0(ji,jj,jk)
424                  zsurfc =  e1u(ji,jj) * ( 1. - umask(ji  ,jj  ,jk) )   &
425                          + e1u(ji,jj) * ( 1. - umask(ji-1,jj  ,jk) )   &
426                          + e2v(ji,jj) * ( 1. - vmask(ji  ,jj  ,jk) )   &
427                          + e2v(ji,jj) * ( 1. - vmask(ji  ,jj-1,jk) )
428                  zsurfp = zsurfc * ze3t / e1e2t(ji,jj)
429                  ! estimation of the coastal slope : 5 km off the coast
430                  ze3t2 = ze3t * ze3t
431                  zcslp = SQRT( ( distcoast*distcoast + ze3t2 ) / ze3t2 )
432                  !
433                  zcmask(ji,jj,jk) = zcmask(ji,jj,jk) + zcslp * zsurfp
[3443]434               END DO
435            END DO
436         END DO
[5385]437         !
[10425]438         CALL lbc_lnk( 'p4zsbc', zcmask , 'T', 1. )      ! lateral boundary conditions on cmask   (sign unchanged)
[5385]439         !
[3443]440         DO jk = 1, jpk
441            DO jj = 1, jpj
442               DO ji = 1, jpi
[6140]443                  zexpide   = MIN( 8.,( gdept_n(ji,jj,jk) / 500. )**(-1.5) )
[3443]444                  zdenitide = -0.9543 + 0.7662 * LOG( zexpide ) - 0.235 * LOG( zexpide )**2
445                  zcmask(ji,jj,jk) = zcmask(ji,jj,jk) * MIN( 1., EXP( zdenitide ) / 0.5 )
446               END DO
447            END DO
448         END DO
449         ! Coastal supply of iron
450         ! -------------------------
[7753]451         ironsed(:,:,jpk) = 0._wp
[3443]452         DO jk = 1, jpkm1
[7753]453            ironsed(:,:,jk) = sedfeinput * zcmask(:,:,jk) / ( e3t_0(:,:,jk) * rday )
[3443]454         END DO
455         DEALLOCATE( zcmask)
456      ENDIF
457      !
458      ! Iron from Hydrothermal vents
459      ! ------------------------
460      IF( ln_hydrofe ) THEN
461         !
[9169]462         IF(lwp) WRITE(numout,*)
463         IF(lwp) WRITE(numout,*) '   ==>>>   ln_hydrofe=T , Input of iron from hydrothermal vents'
[3443]464         !
465         ALLOCATE( hydrofe(jpi,jpj,jpk) )    ! allocation
466         !
467         CALL iom_open ( TRIM( sn_hydrofe%clname ), numhydro )
468         CALL iom_get  ( numhydro, jpdom_data, TRIM( sn_hydrofe%clvar ), hydrofe(:,:,:), 1 )
469         CALL iom_close( numhydro )
470         !
[6945]471         DO jk = 1, jpk
472            hydrofe(:,:,jk) = ( hydrofe(:,:,jk) * hratio ) / ( e1e2t(:,:) * e3t_0(:,:,jk) * ryyss + rtrn ) / 1000._wp
473         ENDDO
[3451]474         !
[3443]475      ENDIF
476      !
477      IF( ll_sbc ) CALL p4z_sbc( nit000 ) 
478      !
479      IF(lwp) THEN
480         WRITE(numout,*)
481         WRITE(numout,*) '    Total input of elements from river supply'
482         WRITE(numout,*) '    ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~'
483         WRITE(numout,*) '    N Supply   : ', rivdininput*rno3*1E3/1E12*14.,' TgN/yr'
[9169]484         WRITE(numout,*) '    Si Supply  : ', rivdsiinput*1E3/1E12*28.1    ,' TgSi/yr'
[3443]485         WRITE(numout,*) '    P Supply   : ', rivdipinput*1E3*po4r/1E12*31.,' TgP/yr'
[9169]486         WRITE(numout,*) '    Alk Supply : ', rivalkinput*1E3/1E12         ,' Teq/yr'
487         WRITE(numout,*) '    DIC Supply : ', rivdicinput*1E3*12./1E12     ,' TgC/yr'
[3443]488         WRITE(numout,*) 
489      ENDIF
490      !
[10127]491      sedsilfrac = 0.03     ! percentage of silica loss in the sediments
492      sedcalfrac = 0.6      ! percentage of calcite loss in the sediments
493      !
[3443]494   END SUBROUTINE p4z_sbc_init
495
496   !!======================================================================
[5656]497END MODULE p4zsbc
Note: See TracBrowser for help on using the repository browser.