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

source: branches/CMIP5_IPSL/NEMO/OPA_SRC/TRA/traqsr.F90 @ 8335

Last change on this file since 8335 was 1829, checked in by cetlod, 14 years ago

New option to have R-B-G light penetration with 3D chlorophyll data

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