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/DEV_r2106_LOCEAN2010/NEMO/OPA_SRC/TRA – NEMO

source: branches/DEV_r2106_LOCEAN2010/NEMO/OPA_SRC/TRA/traqsr.F90 @ 2236

Last change on this file since 2236 was 2236, checked in by cetlod, 13 years ago

First guess of NEMO_v3.3

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 24.4 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   !!----------------------------------------------------------------------
13
14   !!----------------------------------------------------------------------
15   !!   tra_qsr      : trend due to the solar radiation penetration
16   !!   tra_qsr_init : solar radiation penetration initialization
17   !!----------------------------------------------------------------------
18   USE oce             ! ocean dynamics and active tracers
19   USE dom_oce         ! ocean space and time domain
20   USE sbc_oce         ! surface boundary condition: ocean
21   USE trc_oce         ! share SMS/Ocean variables
22   USE trdmod_oce      ! ocean variables trends
23   USE trdtra          ! ocean active tracers trends
24   USE in_out_manager  ! I/O manager
25   USE phycst          ! physical constants
26   USE prtctl          ! Print control
27   USE iom             ! I/O manager
28   USE fldread         ! read input fields
29   USE restart         ! ocean restart
30
31   IMPLICIT NONE
32   PRIVATE
33
34   PUBLIC   tra_qsr       ! routine called by step.F90 (ln_traqsr=T)
35   PUBLIC   tra_qsr_init  ! routine called by opa.F90
36
37   !                                           !!* Namelist namtra_qsr: penetrative solar radiation
38   LOGICAL , PUBLIC ::   ln_traqsr  = .TRUE.    !: light absorption (qsr) flag
39   LOGICAL , PUBLIC ::   ln_qsr_rgb = .FALSE.   !: Red-Green-Blue light absorption flag 
40   LOGICAL , PUBLIC ::   ln_qsr_2bd = .TRUE.    !: 2 band         light absorption flag
41   LOGICAL , PUBLIC ::   ln_qsr_bio = .FALSE.   !: bio-model      light absorption flag
42   INTEGER , PUBLIC ::   nn_chldta  = 0         !: use Chlorophyll data (=1) or not (=0)
43   REAL(wp), PUBLIC ::   rn_abs     = 0.58_wp   !: fraction absorbed in the very near surface (RGB & 2 bands)
44   REAL(wp), PUBLIC ::   rn_si0     = 0.35_wp   !: very near surface depth of extinction      (RGB & 2 bands)
45   REAL(wp), PUBLIC ::   rn_si1     = 23.0_wp   !: deepest depth of extinction (water type I)       (2 bands)
46   REAL(wp), PUBLIC ::   rn_si2     = 61.8_wp   !: deepest depth of extinction (blue & 0.01 mg.m-3)     (RGB)
47   
48   ! Module variables
49!$AGRIF_DO_NOT_TREAT
50   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf_chl   ! structure of input Chl (file informations, fields read)
51   INTEGER, PUBLIC ::   nksr              ! levels below which the light cannot penetrate ( depth larger than 391 m)
52   REAL(wp), DIMENSION(3,61) ::   rkrgb   !: tabulated attenuation coefficients for RGB absorption
53!$AGRIF_END_DO_NOT_TREAT
54
55   !! * Substitutions
56#  include "domzgr_substitute.h90"
57#  include "vectopt_loop_substitute.h90"
58   !!----------------------------------------------------------------------
59   !! NEMO/OPA 3.2 , LOCEAN-IPSL (2009)
60   !! $Id$
61   !! Software governed by the CeCILL licence  (NEMOGCM/License_CeCILL.txt)
62   !!----------------------------------------------------------------------
63
64CONTAINS
65
66   SUBROUTINE tra_qsr( kt )
67      !!----------------------------------------------------------------------
68      !!                  ***  ROUTINE tra_qsr  ***
69      !!
70      !! ** Purpose :   Compute the temperature trend due to the solar radiation
71      !!      penetration and add it to the general temperature trend.
72      !!
73      !! ** Method  : The profile of the solar radiation within the ocean is defined
74      !!      through 2 wavebands (rn_si0,rn_si1) or 3 wavebands (RGB) and a ratio rn_abs
75      !!      Considering the 2 wavebands case:
76      !!         I(k) = Qsr*( rn_abs*EXP(z(k)/rn_si0) + (1.-rn_abs)*EXP(z(k)/rn_si1) )
77      !!         The temperature trend associated with the solar radiation penetration
78      !!         is given by : zta = 1/e3t dk[ I ] / (rau0*Cp)
79      !!         At the bottom, boudary condition for the radiation is no flux :
80      !!      all heat which has not been absorbed in the above levels is put
81      !!      in the last ocean level.
82      !!         In z-coordinate case, the computation is only done down to the
83      !!      level where I(k) < 1.e-15 W/m2. In addition, the coefficients
84      !!      used for the computation are calculated one for once as they
85      !!      depends on k only.
86      !!
87      !! ** Action  : - update ta with the penetrative solar radiation trend
88      !!              - save the trend in ttrd ('key_trdtra')
89      !!
90      !! Reference  : Jerlov, N. G., 1968 Optical Oceanography, Elsevier, 194pp.
91      !!              Lengaigne et al. 2007, Clim. Dyn., V28, 5, 503-516.
92      !!----------------------------------------------------------------------
93      !!
94      INTEGER, INTENT(in) ::   kt     ! ocean time-step
95      !!
96      INTEGER  ::   ji, jj, jk           ! dummy loop indices
97      INTEGER  ::   irgb                 ! temporary integers
98      REAL(wp) ::   zchl, zcoef, zsi0r   ! temporary scalars
99      REAL(wp) ::   zc0, zc1, zc2, zc3   !    -         -
100      REAL(wp) ::   z1_e3t, zfact        !    -         -
101      REAL(wp), DIMENSION(jpi,jpj)     ::   zekb, zekg, zekr            ! 2D workspace
102      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   ze0, ze1 , ze2, ze3, zea    ! 3D workspace
103      REAL(wp), DIMENSION(:,:,:), ALLOCATABLE ::  ztrdt,  ztrds
104      !!----------------------------------------------------------------------
105
106      IF( kt == nit000 ) THEN
107         IF(lwp) WRITE(numout,*)
108         IF(lwp) WRITE(numout,*) 'tra_qsr : penetration of the surface solar radiation'
109         IF(lwp) WRITE(numout,*) '~~~~~~~'
110         IF( .NOT.ln_traqsr )   RETURN
111      ENDIF
112
113      IF( l_trdtra ) THEN      ! Save ta and sa trends
114         ALLOCATE( ztrdt(jpi,jpj,jpk) )   ;    ztrdt(:,:,:) = tsa(:,:,:,jp_tem)
115         ALLOCATE( ztrds(jpi,jpj,jpk) )   ;    ztrds(:,:,:) = 0.
116      ENDIF
117
118      !                                        Set before qsr tracer content field
119      !                                        ***********************************
120      IF( kt == nit000 ) THEN                     ! Set the forcing field at nit000 - 1
121         !                                        ! -----------------------------------
122         IF( ln_rstart .AND.    &                    ! Restart: read in restart file
123              & iom_varid( numror, 'qsr_hc_b', ldstop = .FALSE. ) > 0 ) THEN
124            IF(lwp) WRITE(numout,*) '          nit000-1 qsr tracer content forcing field red in the restart file'
125            zfact = 0.5e0
126            CALL iom_get( numror, jpdom_autoglo, 'qsr_hc_b', qsr_hc_b )   ! before heat content trend due to Qsr flux
127         ELSE                                           ! No restart or restart not found: Euler forward time stepping
128            zfact = 1.e0
129            qsr_hc_b(:,:,:) = 0.e0
130         ENDIF
131      ELSE                                        ! Swap of forcing field
132         !                                        ! ---------------------
133         zfact = 0.5e0
134         qsr_hc_b(:,:,:) = qsr_hc_n(:,:,:)
135      ENDIF
136      !                                        Compute now qsr tracer content field
137      !                                        ************************************
138     
139      !                                           ! ============================================== !
140      IF( lk_qsr_bio .AND. ln_qsr_bio ) THEN      !  bio-model fluxes  : all vertical coordinates  !
141         !                                        ! ============================================== !
142         DO jk = 1, jpkm1
143            DO jj = 2, jpjm1
144               DO ji = fs_2, fs_jpim1   ! vector opt.
145                  qsr_hc_n(ji,jj,jk) =  ro0cpr * ( etot3(ji,jj,jk) - etot3(ji,jj,jk+1) ) / fse3t(ji,jj,jk) 
146               END DO
147            END DO
148         END DO
149         CALL iom_put( 'qsr3d', etot3 )   ! Shortwave Radiation 3D distribution
150         !                                        ! ============================================== !
151      ELSE                                        !  Ocean alone :
152         !                                        ! ============================================== !
153         !
154         !                                                ! ------------------------- !
155         IF( ln_qsr_rgb) THEN                             !  R-G-B  light penetration !
156            !                                             ! ------------------------- !
157            ! Set chlorophyl concentration
158            IF( nn_chldta ==1 ) THEN                             !*  Variable Chlorophyll
159               !
160               CALL fld_read( kt, 1, sf_chl )                         ! Read Chl data and provides it at the current time step
161               !         
162!CDIR COLLAPSE
163!CDIR NOVERRCHK
164               DO jj = 1, jpj                                         ! Separation in R-G-B depending of the surface Chl
165!CDIR NOVERRCHK
166                  DO ji = 1, jpi
167                     zchl = MIN( 10. , MAX( 0.03, sf_chl(1)%fnow(ji,jj,1) ) )
168                     irgb = NINT( 41 + 20.*LOG10(zchl) + 1.e-15 )
169                     zekb(ji,jj) = rkrgb(1,irgb)
170                     zekg(ji,jj) = rkrgb(2,irgb)
171                     zekr(ji,jj) = rkrgb(3,irgb)
172                  END DO
173               END DO
174               !
175               zsi0r = 1.e0 / rn_si0
176               zcoef  = ( 1. - rn_abs ) / 3.e0                        ! equi-partition in R-G-B
177               ze0(:,:,1) = rn_abs  * qsr(:,:)
178               ze1(:,:,1) = zcoef * qsr(:,:)
179               ze2(:,:,1) = zcoef * qsr(:,:)
180               ze3(:,:,1) = zcoef * qsr(:,:)
181               zea(:,:,1) =         qsr(:,:)
182               !
183               DO jk = 2, nksr+1
184!CDIR NOVERRCHK
185                  DO jj = 1, jpj
186!CDIR NOVERRCHK   
187                     DO ji = 1, jpi
188                        zc0 = ze0(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zsi0r     )
189                        zc1 = ze1(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zekb(ji,jj) )
190                        zc2 = ze2(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zekg(ji,jj) )
191                        zc3 = ze3(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zekr(ji,jj) )
192                        ze0(ji,jj,jk) = zc0
193                        ze1(ji,jj,jk) = zc1
194                        ze2(ji,jj,jk) = zc2
195                        ze3(ji,jj,jk) = zc3
196                        zea(ji,jj,jk) = ( zc0 + zc1 + zc2 + zc3 ) * tmask(ji,jj,jk)
197                     END DO
198                  END DO
199               END DO
200               !
201               DO jk = 1, nksr                                        ! compute and add qsr trend to ta
202                  qsr_hc_n(:,:,jk) = ro0cpr * ( zea(:,:,jk) - zea(:,:,jk+1) )
203               END DO
204               zea(:,:,nksr+1:jpk) = 0.e0     ! below 400m set to zero
205               CALL iom_put( 'qsr3d', zea )   ! Shortwave Radiation 3D distribution
206               !
207            ELSE                                                 !*  Constant Chlorophyll
208               DO jk = 1, nksr
209                  qsr_hc_n(:,:,jk) =  etot3(:,:,jk) * qsr(:,:)
210               END DO
211            ENDIF
212
213         ENDIF
214         !                                                ! ------------------------- !
215         IF( ln_qsr_2bd ) THEN                            !  2 band light penetration !
216            !                                             ! ------------------------- !
217            !
218            DO jk = 1, nksr
219               DO jj = 2, jpjm1
220                  DO ji = fs_2, fs_jpim1   ! vector opt.
221                     qsr_hc_n(ji,jj,jk) =  etot3(ji,jj,jk) * qsr(ji,jj)
222                  END DO
223               END DO
224            END DO
225            !
226         ENDIF
227         !
228      ENDIF
229      !                                        Add to the general trend
230      !                                        ************************
231      DO jk = 1, nksr
232         DO jj = 2, jpjm1 
233            DO ji = fs_2, fs_jpim1   ! vector opt.
234               z1_e3t = zfact / fse3t(ji,jj,jk)
235               tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem) + ( qsr_hc_b(ji,jj,jk) + qsr_hc_n(ji,jj,jk) ) * z1_e3t
236            END DO
237         END DO
238      END DO
239      !
240      IF( lrst_oce ) THEN   !                  Write in the ocean restart file
241         !                                     *******************************
242         IF(lwp) WRITE(numout,*)
243         IF(lwp) WRITE(numout,*) 'qsr tracer content forcing field written in ocean restart file ',   &
244            &                    'at it= ', kt,' date= ', ndastp
245         IF(lwp) WRITE(numout,*) '~~~~'
246         CALL iom_rstput( kt, nitrst, numrow, 'qsr_hc_b', qsr_hc_n )
247         !
248      ENDIF
249
250      IF( l_trdtra ) THEN     ! qsr tracers trends saved for diagnostics
251         ztrdt(:,:,:) = tsa(:,:,:,jp_tem) - ztrdt(:,:,:)
252         CALL trd_tra( kt, 'TRA', jp_tem, jptra_trd_qsr, ztrdt )
253         CALL trd_tra( kt, 'TRA', jp_sal, jptra_trd_qsr, ztrds )
254         DEALLOCATE( ztrdt )    ;        DEALLOCATE( ztrds )
255      ENDIF
256      !                       ! print mean trends (used for debugging)
257      IF(ln_ctl)   CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' qsr  - Ta: ', mask1=tmask, clinfo3='tra-ta' )
258      !
259   END SUBROUTINE tra_qsr
260
261
262   SUBROUTINE tra_qsr_init
263      !!----------------------------------------------------------------------
264      !!                  ***  ROUTINE tra_qsr_init  ***
265      !!
266      !! ** Purpose :   Initialization for the penetrative solar radiation
267      !!
268      !! ** Method  :   The profile of solar radiation within the ocean is set
269      !!      from two length scale of penetration (rn_si0,rn_si1) and a ratio
270      !!      (rn_abs). These parameters are read in the namtra_qsr namelist. The
271      !!      default values correspond to clear water (type I in Jerlov'
272      !!      (1968) classification.
273      !!         called by tra_qsr at the first timestep (nit000)
274      !!
275      !! ** Action  : - initialize rn_si0, rn_si1 and rn_abs
276      !!
277      !! Reference : Jerlov, N. G., 1968 Optical Oceanography, Elsevier, 194pp.
278      !!----------------------------------------------------------------------
279      INTEGER  ::   ji, jj, jk            ! dummy loop indices
280      INTEGER  ::   irgb, ierror          ! temporary integer
281      INTEGER  ::   ioptio, nqsr          ! temporary integer
282      REAL(wp) ::   zc0  , zc1            ! temporary scalars
283      REAL(wp) ::   zc2  , zc3  , zchl    !    -         -
284      REAL(wp) ::   zsi0r, zsi1r, zcoef   !    -         -
285      REAL(wp), DIMENSION(jpi,jpj)     ::   zekb, zekg, zekr              ! 2D workspace
286      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   ze0 , ze1 , ze2 , ze3 , zea   ! 3D workspace
287      !!
288      CHARACTER(len=100) ::   cn_dir   ! Root directory for location of ssr files
289      TYPE(FLD_N)        ::   sn_chl   ! informations about the chlorofyl field to be read
290      NAMELIST/namtra_qsr/  sn_chl, cn_dir, ln_traqsr, ln_qsr_rgb, ln_qsr_2bd, ln_qsr_bio,   &
291         &                  nn_chldta, rn_abs, rn_si0, rn_si1, rn_si2
292      !!----------------------------------------------------------------------
293
294      cn_dir = './'       ! directory in which the model is executed
295      ! ... default values (NB: frequency positive => hours, negative => months)
296      !            !     file       ! frequency !  variable  ! time interp !  clim   ! 'yearly' or ! weights  ! rotation   !
297      !            !     name       !  (hours)  !    name    !    (T/F)    !  (T/F)  ! 'monthly'   ! filename ! pairs      !
298      sn_chl = FLD_N( 'chlorophyll' ,    -1     ,  'CHLA'    ,  .true.     , .true.  ,   'yearly'  , ''       , ''         )
299      !
300      REWIND( numnam )            ! Read Namelist namtra_qsr : ratio and length of penetration
301      READ  ( numnam, namtra_qsr )
302      !
303      IF(lwp) THEN                ! control print
304         WRITE(numout,*)
305         WRITE(numout,*) 'tra_qsr_init : penetration of the surface solar radiation'
306         WRITE(numout,*) '~~~~~~~~~~~~'
307         WRITE(numout,*) '   Namelist namtra_qsr : set the parameter of penetration'
308         WRITE(numout,*) '      Light penetration (T) or not (F)         ln_traqsr  = ', ln_traqsr
309         WRITE(numout,*) '      RGB (Red-Green-Blue) light penetration   ln_qsr_rgb = ', ln_qsr_rgb
310         WRITE(numout,*) '      2 band               light penetration   ln_qsr_2bd = ', ln_qsr_2bd
311         WRITE(numout,*) '      bio-model            light penetration   ln_qsr_bio = ', ln_qsr_bio
312         WRITE(numout,*) '      RGB : Chl data (=1) or cst value (=0)    nn_chldta  = ', nn_chldta
313         WRITE(numout,*) '      RGB & 2 bands: fraction of light (rn_si1)    rn_abs = ', rn_abs
314         WRITE(numout,*) '      RGB & 2 bands: shortess depth of extinction  rn_si0 = ', rn_si0
315         WRITE(numout,*) '      2 bands: longest depth of extinction         rn_si1 = ', rn_si1
316         WRITE(numout,*) '      3 bands: longest depth of extinction         rn_si2 = ', rn_si2
317      ENDIF
318
319      IF( ln_traqsr ) THEN     ! control consistency
320         !                     
321         IF( .NOT.lk_qsr_bio .AND. ln_qsr_bio )   THEN
322            CALL ctl_warn( 'No bio model : force ln_qsr_bio = FALSE ' )
323            ln_qsr_bio = .FALSE.
324         ENDIF
325         !
326         ioptio = 0                      ! Parameter control
327         IF( ln_qsr_rgb  )   ioptio = ioptio + 1
328         IF( ln_qsr_2bd  )   ioptio = ioptio + 1
329         IF( ln_qsr_bio  )   ioptio = ioptio + 1
330         !
331         IF( ioptio /= 1 ) THEN
332            ln_qsr_rgb = .TRUE.
333            nn_chldta  = 0
334            ln_qsr_2bd = .FALSE.
335            ln_qsr_bio = .FALSE.
336            CALL ctl_warn( '          Choose ONE type of light penetration in namelist namtra_qsr',   &
337           &               ' otherwise, we force the model to run with RGB light penetration' )
338         ENDIF
339         !
340         IF( ln_qsr_rgb .AND. nn_chldta == 0 )   nqsr =  1 
341         IF( ln_qsr_rgb .AND. nn_chldta == 1 )   nqsr =  2
342         IF( ln_qsr_2bd                      )   nqsr =  3
343         IF( ln_qsr_bio                      )   nqsr =  4
344         !
345         IF(lwp) THEN                   ! Print the choice
346            WRITE(numout,*)
347            IF( nqsr ==  1 )   WRITE(numout,*) '         R-G-B  light penetration - Constant Chlorophyll'
348            IF( nqsr ==  2 )   WRITE(numout,*) '         R-G-B  light penetration - Chl data '
349            IF( nqsr ==  3 )   WRITE(numout,*) '         2 band light penetration'
350            IF( nqsr ==  4 )   WRITE(numout,*) '         bio-model light penetration'
351         ENDIF
352         !
353      ENDIF
354      !                          ! ===================================== !
355      IF( ln_traqsr  ) THEN      !  Initialisation of Light Penetration  ! 
356         !                       ! ===================================== !
357         !
358         zsi0r = 1.e0 / rn_si0
359         zsi1r = 1.e0 / rn_si1
360         !                                ! ---------------------------------- !
361         IF( ln_qsr_rgb ) THEN            !  Red-Green-Blue light penetration  !
362            !                             ! ---------------------------------- !
363            !
364            !                                ! level of light extinction
365            nksr = trc_oce_ext_lev( rn_si2, 0.33e2 )
366            IF(lwp) THEN
367               WRITE(numout,*)
368               WRITE(numout,*) '        level max of computation of qsr = ', nksr, ' ref depth = ', gdepw_0(nksr+1), ' m'
369            ENDIF
370            !
371            CALL trc_oce_rgb( rkrgb )           !* tabulated attenuation coef.
372!!gm            CALL trc_oce_rgb_read( rkrgb )           !* tabulated attenuation coef.
373            !
374            IF( nn_chldta == 1 ) THEN           !* Chl data : set sf_chl structure
375               IF(lwp) WRITE(numout,*)
376               IF(lwp) WRITE(numout,*) '        Chlorophyll read in a file'
377               ALLOCATE( sf_chl(1), STAT=ierror )
378               IF( ierror > 0 ) THEN
379                  CALL ctl_stop( 'tra_qsr_init: unable to allocate sf_chl structure' )   ;   RETURN
380               ENDIF
381               ALLOCATE( sf_chl(1)%fnow(jpi,jpj,1)   )
382               IF( sn_chl%ln_tint )ALLOCATE( sf_chl(1)%fdta(jpi,jpj,1,2) )
383               !                                        ! fill sf_chl with sn_chl and control print
384               CALL fld_fill( sf_chl, (/ sn_chl /), cn_dir, 'tra_qsr_init',   &
385                  &                                         'Solar penetration function of read chlorophyll', 'namtra_qsr' )
386               !
387            ELSE                                !* constant Chl : compute once for all the distribution of light (etot3)
388               IF(lwp) WRITE(numout,*)
389               IF(lwp) WRITE(numout,*) '        Constant Chlorophyll concentration = 0.05'
390               IF(lwp) WRITE(numout,*) '        light distribution computed once for all'
391               !
392               zchl = 0.05                                 ! constant chlorophyll
393               irgb = NINT( 41 + 20.*LOG10(zchl) + 1.e-15 )
394               zekb(:,:) = rkrgb(1,irgb)                   ! Separation in R-G-B depending of the chlorophyl concentration
395               zekg(:,:) = rkrgb(2,irgb)
396               zekr(:,:) = rkrgb(3,irgb)
397               !
398               zcoef = ( 1. - rn_abs ) / 3.e0              ! equi-partition in R-G-B
399               ze0(:,:,1) = rn_abs
400               ze1(:,:,1) = zcoef
401               ze2(:,:,1) = zcoef 
402               ze3(:,:,1) = zcoef
403               zea(:,:,1) = tmask(:,:,1)                   ! = ( ze0+ze1+z2+ze3 ) * tmask
404               
405               DO jk = 2, nksr+1
406!CDIR NOVERRCHK
407                  DO jj = 1, jpj
408!CDIR NOVERRCHK   
409                     DO ji = 1, jpi
410                        zc0 = ze0(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zsi0r     )
411                        zc1 = ze1(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zekb(ji,jj) )
412                        zc2 = ze2(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zekg(ji,jj) )
413                        zc3 = ze3(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zekr(ji,jj) )
414                        ze0(ji,jj,jk) = zc0                 
415                        ze1(ji,jj,jk) = zc1
416                        ze2(ji,jj,jk) = zc2     
417                        ze3(ji,jj,jk) = zc3
418                        zea(ji,jj,jk) = ( zc0 + zc1 + zc2 + zc3 ) * tmask(ji,jj,jk)
419                     END DO
420                  END DO
421               END DO 
422               !
423               DO jk = 1, nksr
424                  etot3(:,:,jk) = ro0cpr * ( zea(:,:,jk) - zea(:,:,jk+1) ) / fse3t(:,:,jk)
425               END DO
426               etot3(:,:,nksr+1:jpk) = 0.e0                ! below 400m set to zero
427            ENDIF
428            !
429         ENDIF
430            !                             ! ---------------------------------- !
431         IF( ln_qsr_2bd ) THEN            !    2 bands    light penetration    !
432            !                             ! ---------------------------------- !
433            !
434            !                                ! level of light extinction
435            nksr = trc_oce_ext_lev( rn_si1, 1.e2 )
436            IF(lwp) THEN
437               WRITE(numout,*)
438               WRITE(numout,*) '        level max of computation of qsr = ', nksr, ' ref depth = ', gdepw_0(nksr+1), ' m'
439            ENDIF
440            !
441            DO jk = 1, nksr                    !*  solar heat absorbed at T-point computed once for all
442               DO jj = 1, jpj                              ! top 400 meters
443                  DO ji = 1, jpi
444                     zc0 = rn_abs * EXP( -fsdepw(ji,jj,jk  )*zsi0r ) + (1.-rn_abs) * EXP( -fsdepw(ji,jj,jk  )*zsi1r )
445                     zc1 = rn_abs * EXP( -fsdepw(ji,jj,jk+1)*zsi0r ) + (1.-rn_abs) * EXP( -fsdepw(ji,jj,jk+1)*zsi1r )
446                     etot3(ji,jj,jk) = ro0cpr * (  zc0 * tmask(ji,jj,jk) - zc1 * tmask(ji,jj,jk+1)  ) / fse3t(ji,jj,jk)
447                  END DO
448               END DO
449            END DO
450            etot3(:,:,nksr+1:jpk) = 0.e0                   ! below 400m set to zero
451            !
452         ENDIF
453         !                       ! ===================================== !
454      ELSE                       !        No light penetration           !                   
455         !                       ! ===================================== !
456         IF(lwp) THEN
457            WRITE(numout,*)
458            WRITE(numout,*) 'tra_qsr_init : NO solar flux penetration'
459            WRITE(numout,*) '~~~~~~~~~~~~'
460         ENDIF
461      ENDIF
462      !
463   END SUBROUTINE tra_qsr_init
464
465   !!======================================================================
466END MODULE traqsr
Note: See TracBrowser for help on using the repository browser.