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 branches/2015/dev_r5836_NOC3_vvl_by_default/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/2015/dev_r5836_NOC3_vvl_by_default/NEMOGCM/NEMO/OPA_SRC/TRA/traqsr.F90 @ 5883

Last change on this file since 5883 was 5883, checked in by gm, 8 years ago

#1613: vvl by default: TRA/TRC remove optimization associated with linear free surface

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