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.
traqsr.F90 in NEMO/branches/UKMO/NEMO_4.0.4_change_chlorophyll/src/OCE/TRA – NEMO

source: NEMO/branches/UKMO/NEMO_4.0.4_change_chlorophyll/src/OCE/TRA/traqsr.F90 @ 15696

Last change on this file since 15696 was 15685, checked in by dancopsey, 2 years ago

Change chlorophyll ancil to be in the same units

File size: 22.8 KB
RevLine 
[3]1MODULE traqsr
2   !!======================================================================
3   !!                       ***  MODULE  traqsr  ***
[6140]4   !! Ocean physics:   solar radiation penetration in the top ocean levels
[3]5   !!======================================================================
[1423]6   !! History :  OPA  !  1990-10  (B. Blanke)  Original code
7   !!            7.0  !  1991-11  (G. Madec)
8   !!                 !  1996-01  (G. Madec)  s-coordinates
9   !!   NEMO     1.0  !  2002-06  (G. Madec)  F90: Free form and module
10   !!             -   !  2005-11  (G. Madec) zco, zps, sco coordinate
11   !!            3.2  !  2009-04  (G. Madec & NEMO team)
[6140]12   !!            3.6  !  2012-05  (C. Rousset) store attenuation coef for use in ice model
[6403]13   !!            3.6  !  2015-12  (O. Aumont, J. Jouanno, C. Ethe) use vertical profile of chlorophyll
[6140]14   !!            3.7  !  2015-11  (G. Madec, A. Coward)  remove optimisation for fix volume
[3]15   !!----------------------------------------------------------------------
[503]16
17   !!----------------------------------------------------------------------
[6140]18   !!   tra_qsr       : temperature trend due to the penetration of solar radiation
19   !!   tra_qsr_init  : initialization of the qsr penetration
[3]20   !!----------------------------------------------------------------------
[6140]21   USE oce            ! ocean dynamics and active tracers
22   USE phycst         ! physical constants
23   USE dom_oce        ! ocean space and time domain
24   USE sbc_oce        ! surface boundary condition: ocean
25   USE trc_oce        ! share SMS/Ocean variables
[4990]26   USE trd_oce        ! trends: ocean variables
27   USE trdtra         ! trends manager: tracers
[6140]28   !
29   USE in_out_manager ! I/O manager
30   USE prtctl         ! Print control
[9019]31   USE iom            ! I/O library
[6140]32   USE fldread        ! read input fields
33   USE restart        ! ocean restart
34   USE lib_mpp        ! MPP library
35   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
[3294]36   USE timing         ! Timing
[3]37
38   IMPLICIT NONE
39   PRIVATE
40
[2528]41   PUBLIC   tra_qsr       ! routine called by step.F90 (ln_traqsr=T)
[5407]42   PUBLIC   tra_qsr_init  ! routine called by nemogcm.F90
[3]43
[4147]44   !                                 !!* Namelist namtra_qsr: penetrative solar radiation
45   LOGICAL , PUBLIC ::   ln_traqsr    !: light absorption (qsr) flag
46   LOGICAL , PUBLIC ::   ln_qsr_rgb   !: Red-Green-Blue light absorption flag 
47   LOGICAL , PUBLIC ::   ln_qsr_2bd   !: 2 band         light absorption flag
48   LOGICAL , PUBLIC ::   ln_qsr_bio   !: bio-model      light absorption flag
49   INTEGER , PUBLIC ::   nn_chldta    !: use Chlorophyll data (=1) or not (=0)
[15678]50   REAL(wp), PUBLIC ::   rn_chl_conc  !: Chlorophyll concentration (for nn_chldta=0)
[4147]51   REAL(wp), PUBLIC ::   rn_abs       !: fraction absorbed in the very near surface (RGB & 2 bands)
52   REAL(wp), PUBLIC ::   rn_si0       !: very near surface depth of extinction      (RGB & 2 bands)
53   REAL(wp), PUBLIC ::   rn_si1       !: deepest depth of extinction (water type I)       (2 bands)
[6140]54   !
55   INTEGER , PUBLIC ::   nksr         !: levels below which the light cannot penetrate (depth larger than 391 m)
[5407]56 
[6140]57   INTEGER, PARAMETER ::   np_RGB  = 1   ! R-G-B     light penetration with constant Chlorophyll
58   INTEGER, PARAMETER ::   np_RGBc = 2   ! R-G-B     light penetration with Chlorophyll data
59   INTEGER, PARAMETER ::   np_2BD  = 3   ! 2 bands   light penetration
60   INTEGER, PARAMETER ::   np_BIO  = 4   ! bio-model light penetration
61   !
62   INTEGER  ::   nqsr    ! user choice of the type of light penetration
63   REAL(wp) ::   xsi0r   ! inverse of rn_si0
64   REAL(wp) ::   xsi1r   ! inverse of rn_si1
65   !
[13331]66   REAL(wp) , PUBLIC, DIMENSION(3,61)   ::   rkrgb    ! tabulated attenuation coefficients for RGB absorption
[1423]67   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf_chl   ! structure of input Chl (file informations, fields read)
[3]68
69   !! * Substitutions
70#  include "vectopt_loop_substitute.h90"
71   !!----------------------------------------------------------------------
[9598]72   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
[888]73   !! $Id$
[10068]74   !! Software governed by the CeCILL license (see ./LICENSE)
[3]75   !!----------------------------------------------------------------------
76CONTAINS
77
78   SUBROUTINE tra_qsr( kt )
79      !!----------------------------------------------------------------------
80      !!                  ***  ROUTINE tra_qsr  ***
81      !!
82      !! ** Purpose :   Compute the temperature trend due to the solar radiation
[6140]83      !!              penetration and add it to the general temperature trend.
[3]84      !!
[1423]85      !! ** Method  : The profile of the solar radiation within the ocean is defined
86      !!      through 2 wavebands (rn_si0,rn_si1) or 3 wavebands (RGB) and a ratio rn_abs
87      !!      Considering the 2 wavebands case:
88      !!         I(k) = Qsr*( rn_abs*EXP(z(k)/rn_si0) + (1.-rn_abs)*EXP(z(k)/rn_si1) )
89      !!         The temperature trend associated with the solar radiation penetration
90      !!         is given by : zta = 1/e3t dk[ I ] / (rau0*Cp)
[3]91      !!         At the bottom, boudary condition for the radiation is no flux :
92      !!      all heat which has not been absorbed in the above levels is put
93      !!      in the last ocean level.
[6140]94      !!         The computation is only done down to the level where
95      !!      I(k) < 1.e-15 W/m2 (i.e. over the top nksr levels) .
[3]96      !!
97      !! ** Action  : - update ta with the penetrative solar radiation trend
[6140]98      !!              - send  trend for further diagnostics (l_trdtra=T)
[1423]99      !!
100      !! Reference  : Jerlov, N. G., 1968 Optical Oceanography, Elsevier, 194pp.
101      !!              Lengaigne et al. 2007, Clim. Dyn., V28, 5, 503-516.
[6403]102      !!              Morel, A. et Berthon, JF, 1989, Limnol Oceanogr 34(8), 1545-1562
[503]103      !!----------------------------------------------------------------------
104      INTEGER, INTENT(in) ::   kt     ! ocean time-step
[2715]105      !
[6140]106      INTEGER  ::   ji, jj, jk               ! dummy loop indices
107      INTEGER  ::   irgb                     ! local integers
108      REAL(wp) ::   zchl, zcoef, z1_2        ! local scalars
109      REAL(wp) ::   zc0 , zc1 , zc2 , zc3    !    -         -
[4161]110      REAL(wp) ::   zzc0, zzc1, zzc2, zzc3   !    -         -
[6140]111      REAL(wp) ::   zz0 , zz1                !    -         -
[6403]112      REAL(wp) ::   zCb, zCmax, zze, zpsi, zpsimax, zdelpsi, zCtot, zCze
113      REAL(wp) ::   zlogc, zlogc2, zlogc3 
[9019]114      REAL(wp), ALLOCATABLE, DIMENSION(:,:)   :: zekb, zekg, zekr
115      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: ze0, ze1, ze2, ze3, zea, ztrdt
116      REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: zetot, zchl3d
[3]117      !!----------------------------------------------------------------------
[3294]118      !
[9019]119      IF( ln_timing )   CALL timing_start('tra_qsr')
[3294]120      !
[3]121      IF( kt == nit000 ) THEN
[503]122         IF(lwp) WRITE(numout,*)
123         IF(lwp) WRITE(numout,*) 'tra_qsr : penetration of the surface solar radiation'
124         IF(lwp) WRITE(numout,*) '~~~~~~~'
[3]125      ENDIF
[6140]126      !
127      IF( l_trdtra ) THEN      ! trends diagnostic: save the input temperature trend
[9019]128         ALLOCATE( ztrdt(jpi,jpj,jpk) ) 
[7753]129         ztrdt(:,:,:) = tsa(:,:,:,jp_tem)
[216]130      ENDIF
[6140]131      !
132      !                         !-----------------------------------!
133      !                         !  before qsr induced heat content  !
134      !                         !-----------------------------------!
135      IF( kt == nit000 ) THEN          !==  1st time step  ==!
136!!gm case neuler  not taken into account....
137         IF( ln_rstart .AND. iom_varid( numror, 'qsr_hc_b', ldstop = .FALSE. ) > 0 ) THEN    ! read in restart
138            IF(lwp) WRITE(numout,*) '          nit000-1 qsr tracer content forcing field read in the restart file'
139            z1_2 = 0.5_wp
[9367]140            CALL iom_get( numror, jpdom_autoglo, 'qsr_hc_b', qsr_hc_b, ldxios = lrxios )   ! before heat content trend due to Qsr flux
[2528]141         ELSE                                           ! No restart or restart not found: Euler forward time stepping
[6140]142            z1_2 = 1._wp
[7753]143            qsr_hc_b(:,:,:) = 0._wp
[2528]144         ENDIF
[6140]145      ELSE                             !==  Swap of qsr heat content  ==!
146         z1_2 = 0.5_wp
[7753]147         qsr_hc_b(:,:,:) = qsr_hc(:,:,:)
[2528]148      ENDIF
[6140]149      !
150      !                         !--------------------------------!
151      SELECT CASE( nqsr )       !  now qsr induced heat content  !
152      !                         !--------------------------------!
153      !
154      CASE( np_BIO )                   !==  bio-model fluxes  ==!
155         !
156         DO jk = 1, nksr
[7753]157            qsr_hc(:,:,jk) = r1_rau0_rcp * ( etot3(:,:,jk) - etot3(:,:,jk+1) )
[2528]158         END DO
[6140]159         !
160      CASE( np_RGB , np_RGBc )         !==  R-G-B fluxes  ==!
161         !
[9019]162         ALLOCATE( zekb(jpi,jpj)     , zekg(jpi,jpj)     , zekr  (jpi,jpj)     , &
163            &      ze0 (jpi,jpj,jpk) , ze1 (jpi,jpj,jpk) , ze2   (jpi,jpj,jpk) , &
164            &      ze3 (jpi,jpj,jpk) , zea (jpi,jpj,jpk) , zchl3d(jpi,jpj,jpk)   ) 
[6140]165         !
166         IF( nqsr == np_RGBc ) THEN          !*  Variable Chlorophyll
167            CALL fld_read( kt, 1, sf_chl )         ! Read Chl data and provides it at the current time step
[6403]168            DO jk = 1, nksr + 1
169               DO jj = 2, jpjm1                       ! Separation in R-G-B depending of the surface Chl
170                  DO ji = fs_2, fs_jpim1
[15685]171
172                     ! Commenting out this bit of code as I want to load in chlorophyll concentrations
173                     ! with the same units as rn_chl_conc
174                     ! zchl    = MIN( 10. , MAX( 0.03, sf_chl(1)%fnow(ji,jj,1) ) )
175                     ! zCtot   = 40.6  * zchl**0.459
176                     ! zze     = 568.2 * zCtot**(-0.746)
177                     ! IF( zze > 102. ) zze = 200.0 * zCtot**(-0.293)
178                     ! zpsi    = gdepw_n(ji,jj,jk) / zze
[6403]179                     !
[15685]180                     ! zlogc   = LOG( zchl )
181                     ! zlogc2  = zlogc * zlogc
182                     ! zlogc3  = zlogc * zlogc * zlogc
183                     ! zCb     = 0.768 + 0.087 * zlogc - 0.179 * zlogc2 - 0.025 * zlogc3
184                     ! zCmax   = 0.299 - 0.289 * zlogc + 0.579 * zlogc2
185                     ! zpsimax = 0.6   - 0.640 * zlogc + 0.021 * zlogc2 + 0.115 * zlogc3
186                     ! zdelpsi = 0.710 + 0.159 * zlogc + 0.021 * zlogc2
187                     ! zCze    = 1.12  * (zchl)**0.803
[6403]188                     !
[15685]189                     ! zchl3d(ji,jj,jk) = zCze * ( zCb + zCmax * EXP( -( (zpsi - zpsimax) / zdelpsi )**2 ) )
190
191                     zchl3d(ji,jj,jk) = MIN( 50. , MAX( 0.01, sf_chl(1)%fnow(ji,jj,1) ) )
[6403]192                  END DO
193                  !
[3]194               END DO
195            END DO
[6140]196         ELSE                                !* constant chrlorophyll
[6403]197           DO jk = 1, nksr + 1
[15678]198              zchl3d(:,:,jk) = rn_chl_conc
[6403]199            ENDDO
[4161]200         ENDIF
[1423]201         !
[6140]202         zcoef  = ( 1. - rn_abs ) / 3._wp    !* surface equi-partition in R-G-B
203         DO jj = 2, jpjm1
204            DO ji = fs_2, fs_jpim1
205               ze0(ji,jj,1) = rn_abs * qsr(ji,jj)
206               ze1(ji,jj,1) = zcoef  * qsr(ji,jj)
207               ze2(ji,jj,1) = zcoef  * qsr(ji,jj)
208               ze3(ji,jj,1) = zcoef  * qsr(ji,jj)
209               zea(ji,jj,1) =          qsr(ji,jj)
210            END DO
211         END DO
212         !
[6403]213         DO jk = 2, nksr+1                   !* interior equi-partition in R-G-B depending of vertical profile of Chl
[6140]214            DO jj = 2, jpjm1
215               DO ji = fs_2, fs_jpim1
[6403]216                  zchl = MIN( 10. , MAX( 0.03, zchl3d(ji,jj,jk) ) )
217                  irgb = NINT( 41 + 20.*LOG10(zchl) + 1.e-15 )
218                  zekb(ji,jj) = rkrgb(1,irgb)
219                  zekg(ji,jj) = rkrgb(2,irgb)
220                  zekr(ji,jj) = rkrgb(3,irgb)
221               END DO
222            END DO
223
224            DO jj = 2, jpjm1
225               DO ji = fs_2, fs_jpim1
[6140]226                  zc0 = ze0(ji,jj,jk-1) * EXP( - e3t_n(ji,jj,jk-1) * xsi0r       )
227                  zc1 = ze1(ji,jj,jk-1) * EXP( - e3t_n(ji,jj,jk-1) * zekb(ji,jj) )
228                  zc2 = ze2(ji,jj,jk-1) * EXP( - e3t_n(ji,jj,jk-1) * zekg(ji,jj) )
229                  zc3 = ze3(ji,jj,jk-1) * EXP( - e3t_n(ji,jj,jk-1) * zekr(ji,jj) )
230                  ze0(ji,jj,jk) = zc0
231                  ze1(ji,jj,jk) = zc1
232                  ze2(ji,jj,jk) = zc2
233                  ze3(ji,jj,jk) = zc3
234                  zea(ji,jj,jk) = ( zc0 + zc1 + zc2 + zc3 ) * wmask(ji,jj,jk)
[1423]235               END DO
[6140]236            END DO
237         END DO
238         !
239         DO jk = 1, nksr                     !* now qsr induced heat content
240            DO jj = 2, jpjm1
241               DO ji = fs_2, fs_jpim1
242                  qsr_hc(ji,jj,jk) = r1_rau0_rcp * ( zea(ji,jj,jk) - zea(ji,jj,jk+1) )
[1423]243               END DO
[6140]244            END DO
245         END DO
[187]246         !
[9019]247         DEALLOCATE( zekb , zekg , zekr , ze0 , ze1 , ze2 , ze3 , zea , zchl3d ) 
[6140]248         !
249      CASE( np_2BD  )            !==  2-bands fluxes  ==!
250         !
251         zz0 =        rn_abs   * r1_rau0_rcp      ! surface equi-partition in 2-bands
252         zz1 = ( 1. - rn_abs ) * r1_rau0_rcp
253         DO jk = 1, nksr                          ! solar heat absorbed at T-point in the top 400m
254            DO jj = 2, jpjm1
255               DO ji = fs_2, fs_jpim1
256                  zc0 = zz0 * EXP( -gdepw_n(ji,jj,jk  )*xsi0r ) + zz1 * EXP( -gdepw_n(ji,jj,jk  )*xsi1r )
257                  zc1 = zz0 * EXP( -gdepw_n(ji,jj,jk+1)*xsi0r ) + zz1 * EXP( -gdepw_n(ji,jj,jk+1)*xsi1r )
258                  qsr_hc(ji,jj,jk) = qsr(ji,jj) * ( zc0 * wmask(ji,jj,jk) - zc1 * wmask(ji,jj,jk+1) ) 
[2528]259               END DO
260            END DO
261         END DO
262         !
[6140]263      END SELECT
264      !
265      !                          !-----------------------------!
266      DO jk = 1, nksr            !  update to the temp. trend  !
267         DO jj = 2, jpjm1        !-----------------------------!
268            DO ji = fs_2, fs_jpim1   ! vector opt.
269               tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem)   &
270                  &                 + z1_2 * ( qsr_hc_b(ji,jj,jk) + qsr_hc(ji,jj,jk) ) / e3t_n(ji,jj,jk)
271            END DO
272         END DO
273      END DO
274      !
[9019]275      ! sea-ice: store the 1st ocean level attenuation coefficient
276      DO jj = 2, jpjm1 
277         DO ji = fs_2, fs_jpim1   ! vector opt.
278            IF( qsr(ji,jj) /= 0._wp ) THEN   ;   fraqsr_1lev(ji,jj) = qsr_hc(ji,jj,1) / ( r1_rau0_rcp * qsr(ji,jj) )
279            ELSE                             ;   fraqsr_1lev(ji,jj) = 1._wp
280            ENDIF
[6140]281         END DO
[9019]282      END DO
[10425]283      CALL lbc_lnk( 'traqsr', fraqsr_1lev(:,:), 'T', 1._wp )
[2528]284      !
[6140]285      IF( iom_use('qsr3d') ) THEN      ! output the shortwave Radiation distribution
[9019]286         ALLOCATE( zetot(jpi,jpj,jpk) )
[7753]287         zetot(:,:,nksr+1:jpk) = 0._wp     ! below ~400m set to zero
[6140]288         DO jk = nksr, 1, -1
[9019]289            zetot(:,:,jk) = zetot(:,:,jk+1) + qsr_hc(:,:,jk) * rau0_rcp
[6140]290         END DO         
291         CALL iom_put( 'qsr3d', zetot )   ! 3D distribution of shortwave Radiation
[9019]292         DEALLOCATE( zetot ) 
[2528]293      ENDIF
[6140]294      !
295      IF( lrst_oce ) THEN     ! write in the ocean restart file
[9367]296         IF( lwxios ) CALL iom_swap(      cwxios_context          )
297         CALL iom_rstput( kt, nitrst, numrow, 'qsr_hc_b'   , qsr_hc     , ldxios = lwxios )
298         CALL iom_rstput( kt, nitrst, numrow, 'fraqsr_1lev', fraqsr_1lev, ldxios = lwxios ) 
299         IF( lwxios ) CALL iom_swap(      cxios_context          )
[6140]300      ENDIF
301      !
[503]302      IF( l_trdtra ) THEN     ! qsr tracers trends saved for diagnostics
[7753]303         ztrdt(:,:,:) = tsa(:,:,:,jp_tem) - ztrdt(:,:,:)
[4990]304         CALL trd_tra( kt, 'TRA', jp_tem, jptra_qsr, ztrdt )
[9019]305         DEALLOCATE( ztrdt ) 
[3]306      ENDIF
[457]307      !                       ! print mean trends (used for debugging)
[2528]308      IF(ln_ctl)   CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' qsr  - Ta: ', mask1=tmask, clinfo3='tra-ta' )
[503]309      !
[9019]310      IF( ln_timing )   CALL timing_stop('tra_qsr')
[3294]311      !
[3]312   END SUBROUTINE tra_qsr
313
314
315   SUBROUTINE tra_qsr_init
316      !!----------------------------------------------------------------------
317      !!                  ***  ROUTINE tra_qsr_init  ***
318      !!
319      !! ** Purpose :   Initialization for the penetrative solar radiation
320      !!
321      !! ** Method  :   The profile of solar radiation within the ocean is set
[1423]322      !!      from two length scale of penetration (rn_si0,rn_si1) and a ratio
[1601]323      !!      (rn_abs). These parameters are read in the namtra_qsr namelist. The
[3]324      !!      default values correspond to clear water (type I in Jerlov'
325      !!      (1968) classification.
326      !!         called by tra_qsr at the first timestep (nit000)
327      !!
[1423]328      !! ** Action  : - initialize rn_si0, rn_si1 and rn_abs
[3]329      !!
[503]330      !! Reference : Jerlov, N. G., 1968 Optical Oceanography, Elsevier, 194pp.
[3]331      !!----------------------------------------------------------------------
[6140]332      INTEGER  ::   ji, jj, jk                  ! dummy loop indices
333      INTEGER  ::   ios, irgb, ierror, ioptio   ! local integer
334      REAL(wp) ::   zz0, zc0 , zc1, zcoef      ! local scalars
335      REAL(wp) ::   zz1, zc2 , zc3, zchl       !   -      -
[2715]336      !
[1423]337      CHARACTER(len=100) ::   cn_dir   ! Root directory for location of ssr files
338      TYPE(FLD_N)        ::   sn_chl   ! informations about the chlorofyl field to be read
[2715]339      !!
[9019]340      NAMELIST/namtra_qsr/  sn_chl, cn_dir, ln_qsr_rgb, ln_qsr_2bd, ln_qsr_bio,  &
[15678]341         &                  nn_chldta, rn_abs, rn_chl_conc, rn_si0, rn_si1
[3]342      !!----------------------------------------------------------------------
[3294]343      !
[6140]344      REWIND( numnam_ref )              ! Namelist namtra_qsr in reference     namelist
345      READ  ( numnam_ref, namtra_qsr, IOSTAT = ios, ERR = 901)
[11536]346901   IF( ios /= 0 )   CALL ctl_nam ( ios , 'namtra_qsr in reference namelist' )
[3294]347      !
[6140]348      REWIND( numnam_cfg )              ! Namelist namtra_qsr in configuration namelist
[4147]349      READ  ( numnam_cfg, namtra_qsr, IOSTAT = ios, ERR = 902 )
[11536]350902   IF( ios >  0 )   CALL ctl_nam ( ios , 'namtra_qsr in configuration namelist' )
[4624]351      IF(lwm) WRITE ( numond, namtra_qsr )
[1423]352      !
353      IF(lwp) THEN                ! control print
354         WRITE(numout,*)
355         WRITE(numout,*) 'tra_qsr_init : penetration of the surface solar radiation'
356         WRITE(numout,*) '~~~~~~~~~~~~'
[1601]357         WRITE(numout,*) '   Namelist namtra_qsr : set the parameter of penetration'
[6140]358         WRITE(numout,*) '      RGB (Red-Green-Blue) light penetration       ln_qsr_rgb = ', ln_qsr_rgb
359         WRITE(numout,*) '      2 band               light penetration       ln_qsr_2bd = ', ln_qsr_2bd
360         WRITE(numout,*) '      bio-model            light penetration       ln_qsr_bio = ', ln_qsr_bio
361         WRITE(numout,*) '      RGB : Chl data (=1) or cst value (=0)        nn_chldta  = ', nn_chldta
[15678]362         WRITE(numout,*) '      Chlorophyll concentration (for nn_chldta=0)  rn_chl_conc = ', rn_chl_conc
[6140]363         WRITE(numout,*) '      RGB & 2 bands: fraction of light (rn_si1)    rn_abs     = ', rn_abs
364         WRITE(numout,*) '      RGB & 2 bands: shortess depth of extinction  rn_si0     = ', rn_si0
365         WRITE(numout,*) '      2 bands: longest depth of extinction         rn_si1     = ', rn_si1
366         WRITE(numout,*)
[1423]367      ENDIF
[6140]368      !
369      ioptio = 0                    ! Parameter control
370      IF( ln_qsr_rgb  )   ioptio = ioptio + 1
371      IF( ln_qsr_2bd  )   ioptio = ioptio + 1
372      IF( ln_qsr_bio  )   ioptio = ioptio + 1
373      !
374      IF( ioptio /= 1 )   CALL ctl_stop( 'Choose ONE type of light penetration in namelist namtra_qsr',  &
375         &                               ' 2 bands, 3 RGB bands or bio-model light penetration' )
376      !
377      IF( ln_qsr_rgb .AND. nn_chldta == 0 )   nqsr = np_RGB 
378      IF( ln_qsr_rgb .AND. nn_chldta == 1 )   nqsr = np_RGBc
379      IF( ln_qsr_2bd                      )   nqsr = np_2BD
380      IF( ln_qsr_bio                      )   nqsr = np_BIO
381      !
382      !                             ! Initialisation
383      xsi0r = 1._wp / rn_si0
384      xsi1r = 1._wp / rn_si1
385      !
386      SELECT CASE( nqsr )
387      !                               
388      CASE( np_RGB , np_RGBc )         !==  Red-Green-Blue light penetration  ==!
389         !                             
[9169]390         IF(lwp)   WRITE(numout,*) '   ==>>>   R-G-B   light penetration '
[6140]391         !
392         CALL trc_oce_rgb( rkrgb )                 ! tabulated attenuation coef.
393         !                                   
394         nksr = trc_oce_ext_lev( r_si2, 33._wp )   ! level of light extinction
395         !
396         IF(lwp) WRITE(numout,*) '        level of light extinction = ', nksr, ' ref depth = ', gdepw_1d(nksr+1), ' m'
397         !
398         IF( nqsr == np_RGBc ) THEN                ! Chl data : set sf_chl structure
[9169]399            IF(lwp) WRITE(numout,*) '   ==>>>   Chlorophyll read in a file'
[6140]400            ALLOCATE( sf_chl(1), STAT=ierror )
401            IF( ierror > 0 ) THEN
402               CALL ctl_stop( 'tra_qsr_init: unable to allocate sf_chl structure' )   ;   RETURN
403            ENDIF
404            ALLOCATE( sf_chl(1)%fnow(jpi,jpj,1)   )
405            IF( sn_chl%ln_tint )   ALLOCATE( sf_chl(1)%fdta(jpi,jpj,1,2) )
406            !                                        ! fill sf_chl with sn_chl and control print
407            CALL fld_fill( sf_chl, (/ sn_chl /), cn_dir, 'tra_qsr_init',   &
[7646]408               &           'Solar penetration function of read chlorophyll', 'namtra_qsr' , no_print )
[1448]409         ENDIF
[6140]410         IF( nqsr == np_RGB ) THEN                 ! constant Chl
[15685]411            IF(lwp) WRITE(numout,*) '   ==>>>   Constant Chlorophyll concentration = ',rn_chl_conc
[6140]412         ENDIF
[1448]413         !
[6140]414      CASE( np_2BD )                   !==  2 bands light penetration  ==!
[1448]415         !
[9169]416         IF(lwp)  WRITE(numout,*) '   ==>>>   2 bands light penetration'
[1448]417         !
[6140]418         nksr = trc_oce_ext_lev( rn_si1, 100._wp )    ! level of light extinction
419         IF(lwp) WRITE(numout,*) '        level of light extinction = ', nksr, ' ref depth = ', gdepw_1d(nksr+1), ' m'
[1455]420         !
[6140]421      CASE( np_BIO )                   !==  BIO light penetration  ==!
[1448]422         !
[9169]423         IF(lwp) WRITE(numout,*) '   ==>>>   bio-model light penetration'
[7646]424         IF( .NOT.lk_top )   CALL ctl_stop( 'No bio model : ln_qsr_bio = true impossible ' )
[13331]425         !                                   
426         CALL trc_oce_rgb( rkrgb )                 ! tabulated attenuation coef.
427         !                                   
428         nksr = trc_oce_ext_lev( r_si2, 33._wp )   ! level of light extinction
[1423]429         !
[13331]430         IF(lwp) WRITE(numout,*) '        level of light extinction = ', nksr, ' ref depth = ', gdepw_1d(nksr+1), ' m'
431         !
[6140]432      END SELECT
[503]433      !
[7753]434      qsr_hc(:,:,:) = 0._wp     ! now qsr heat content set to zero where it will not be computed
[6140]435      !
436      ! 1st ocean level attenuation coefficient (used in sbcssm)
[5407]437      IF( iom_varid( numror, 'fraqsr_1lev', ldstop = .FALSE. ) > 0 ) THEN
[9367]438         CALL iom_get( numror, jpdom_autoglo, 'fraqsr_1lev'  , fraqsr_1lev, ldxios = lrxios  )
[5407]439      ELSE
[7753]440         fraqsr_1lev(:,:) = 1._wp   ! default : no penetration
[5407]441      ENDIF
442      !
[9367]443      IF( lwxios ) THEN
444         CALL iom_set_rstw_var_active('qsr_hc_b')
445         CALL iom_set_rstw_var_active('fraqsr_1lev')
446      ENDIF
447      !
[3]448   END SUBROUTINE tra_qsr_init
449
450   !!======================================================================
451END MODULE traqsr
Note: See TracBrowser for help on using the repository browser.