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/2014/dev_CNRS1_10_TEOS10_Ediag/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/2014/dev_CNRS1_10_TEOS10_Ediag/NEMOGCM/NEMO/OPA_SRC/TRA/traqsr.F90 @ 4915

Last change on this file since 4915 was 4619, checked in by gm, 10 years ago

#1294 : TEOS-10 and Ediag

  • Property svn:keywords set to Id
File size: 29.9 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   !!----------------------------------------------------------------------
14
15   !!----------------------------------------------------------------------
16   !!   tra_qsr      : trend due to the solar radiation penetration
17   !!   tra_qsr_init : solar radiation penetration initialization
18   !!----------------------------------------------------------------------
19   USE oce            ! ocean dynamics and active tracers
20   USE dom_oce        ! ocean space and time domain
21   USE sbc_oce        ! surface boundary condition: ocean
22   USE trc_oce        ! share SMS/Ocean variables
23   USE trd_oce        ! trends: ocean variables
24   USE trdtra         ! trends manager: tracers
25   USE phycst         ! physical constants
26   USE sbc_ice,  ONLY : lk_lim3
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 lib_mpp        ! MPP library
33   USE wrk_nemo       ! Memory Allocation
34   USE timing         ! Timing
35
36   IMPLICIT NONE
37   PRIVATE
38
39   PUBLIC   tra_qsr       ! routine called by step.F90 (ln_traqsr=T)
40   PUBLIC   tra_qsr_init  ! routine called by opa.F90
41
42   !                                 !!* Namelist namtra_qsr: penetrative solar radiation
43   LOGICAL , PUBLIC ::   ln_traqsr    !: light absorption (qsr) flag
44   LOGICAL , PUBLIC ::   ln_qsr_rgb   !: Red-Green-Blue light absorption flag 
45   LOGICAL , PUBLIC ::   ln_qsr_2bd   !: 2 band         light absorption flag
46   LOGICAL , PUBLIC ::   ln_qsr_bio   !: bio-model      light absorption flag
47   LOGICAL , PUBLIC ::   ln_qsr_ice   !: light penetration for ice-model LIM3 (clem)
48   INTEGER , PUBLIC ::   nn_chldta    !: use Chlorophyll data (=1) or not (=0)
49   REAL(wp), PUBLIC ::   rn_abs       !: fraction absorbed in the very near surface (RGB & 2 bands)
50   REAL(wp), PUBLIC ::   rn_si0       !: very near surface depth of extinction      (RGB & 2 bands)
51   REAL(wp), PUBLIC ::   rn_si1       !: deepest depth of extinction (water type I)       (2 bands)
52   
53   INTEGER , PUBLIC ::   nksr   !: levels below which the light cannot penetrate ( depth larger than 391 m)
54
55   REAL(wp)                  ::   xsi0r, xsi1r        ! inverse of rn_si0 and rn_si1, resp.
56   REAL(wp), DIMENSION(3,61) ::   rkrgb               ! tabulated attenuation coefficients for RGB absorption
57   TYPE(FLD), ALLOCATABLE, DIMENSION(:) ::   sf_chl   ! structure of input Chl (file informations, fields read)
58
59   !! * Substitutions
60#  include "domzgr_substitute.h90"
61#  include "vectopt_loop_substitute.h90"
62   !!----------------------------------------------------------------------
63   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
64   !! $Id$
65   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
66   !!----------------------------------------------------------------------
67CONTAINS
68
69   SUBROUTINE tra_qsr( kt )
70      !!----------------------------------------------------------------------
71      !!                  ***  ROUTINE tra_qsr  ***
72      !!
73      !! ** Purpose :   Compute the temperature trend due to the solar radiation
74      !!      penetration and add it to the general temperature trend.
75      !!
76      !! ** Method  : The profile of the solar radiation within the ocean is defined
77      !!      through 2 wavebands (rn_si0,rn_si1) or 3 wavebands (RGB) and a ratio rn_abs
78      !!      Considering the 2 wavebands case:
79      !!         I(k) = Qsr*( rn_abs*EXP(z(k)/rn_si0) + (1.-rn_abs)*EXP(z(k)/rn_si1) )
80      !!         The temperature trend associated with the solar radiation penetration
81      !!         is given by : zta = 1/e3t dk[ I ] / (rau0*Cp)
82      !!         At the bottom, boudary condition for the radiation is no flux :
83      !!      all heat which has not been absorbed in the above levels is put
84      !!      in the last ocean level.
85      !!         In z-coordinate case, the computation is only done down to the
86      !!      level where I(k) < 1.e-15 W/m2. In addition, the coefficients
87      !!      used for the computation are calculated one for once as they
88      !!      depends on k only.
89      !!
90      !! ** Action  : - update ta with the penetrative solar radiation trend
91      !!              - send the trend to trdtra (l_trdtra=T)
92      !!
93      !! Reference  : Jerlov, N. G., 1968 Optical Oceanography, Elsevier, 194pp.
94      !!              Lengaigne et al. 2007, Clim. Dyn., V28, 5, 503-516.
95      !!----------------------------------------------------------------------
96      INTEGER, INTENT(in) ::   kt     ! ocean time-step
97      !
98      INTEGER  ::   ji, jj, jk           ! dummy loop indices
99      INTEGER  ::   irgb                 ! local integers
100      REAL(wp) ::   zchl, zcoef, zfact   ! local scalars
101      REAL(wp) ::   zc0, zc1, zc2, zc3   !    -         -
102      REAL(wp) ::   zzc0, zzc1, zzc2, zzc3   !    -         -
103      REAL(wp) ::   zz0, zz1, z1_e3t     !    -         -
104      REAL(wp), POINTER, DIMENSION(:,:  ) :: zekb, zekg, zekr
105      REAL(wp), POINTER, DIMENSION(:,:,:) :: ze0, ze1, ze2, ze3, zea, ztrdt
106      !!----------------------------------------------------------------------
107      !
108      IF( nn_timing == 1 )  CALL timing_start('tra_qsr')
109      !
110      CALL wrk_alloc( jpi, jpj,      zekb, zekg, zekr        ) 
111      CALL wrk_alloc( jpi, jpj, jpk, ze0, ze1, ze2, ze3, zea ) 
112      !
113      IF( kt == nit000 ) THEN
114         IF(lwp) WRITE(numout,*)
115         IF(lwp) WRITE(numout,*) 'tra_qsr : penetration of the surface solar radiation'
116         IF(lwp) WRITE(numout,*) '~~~~~~~'
117         IF( .NOT.ln_traqsr )   RETURN
118      ENDIF
119
120      IF( l_trdtra ) THEN      ! Save temperature trends
121         CALL wrk_alloc( jpi, jpj, jpk, ztrdt ) 
122         ztrdt(:,:,:) = tsa(:,:,:,jp_tem)
123      ENDIF
124
125      !                                        Set before qsr tracer content field
126      !                                        ***********************************
127      IF( kt == nit000 ) THEN                     ! Set the forcing field at nit000 - 1
128         !                                        ! -----------------------------------
129         IF( ln_rstart .AND.    &                    ! Restart: read in restart file
130              & iom_varid( numror, 'qsr_hc_b', ldstop = .FALSE. ) > 0 ) THEN
131            IF(lwp) WRITE(numout,*) '          nit000-1 qsr tracer content forcing field red in the restart file'
132            zfact = 0.5e0
133            CALL iom_get( numror, jpdom_autoglo, 'qsr_hc_b', qsr_hc_b )   ! before heat content trend due to Qsr flux
134         ELSE                                           ! No restart or restart not found: Euler forward time stepping
135            zfact = 1.e0
136            qsr_hc_b(:,:,:) = 0.e0
137         ENDIF
138      ELSE                                        ! Swap of forcing field
139         !                                        ! ---------------------
140         zfact = 0.5e0
141         qsr_hc_b(:,:,:) = qsr_hc(:,:,:)
142      ENDIF
143      !                                        Compute now qsr tracer content field
144      !                                        ************************************
145      !                                           ! ============================================== !
146      IF( lk_qsr_bio .AND. ln_qsr_bio ) THEN      !  bio-model fluxes  : all vertical coordinates  !
147         !                                        ! ============================================== !
148         DO jk = 1, jpkm1
149            qsr_hc(:,:,jk) = r1_rau0_rcp * ( etot3(:,:,jk) - etot3(:,:,jk+1) )
150         END DO
151         !                                        Add to the general trend
152         DO jk = 1, jpkm1
153            DO jj = 2, jpjm1 
154               DO ji = fs_2, fs_jpim1   ! vector opt.
155                  z1_e3t = zfact / fse3t(ji,jj,jk)
156                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem) + ( qsr_hc_b(ji,jj,jk) + qsr_hc(ji,jj,jk) ) * z1_e3t
157               END DO
158            END DO
159         END DO
160         CALL iom_put( 'qsr3d', etot3 )   ! Shortwave Radiation 3D distribution
161         ! clem: store attenuation coefficient of the first ocean level
162         IF ( lk_lim3 .AND. ln_qsr_ice ) THEN
163            DO jj = 1, jpj
164               DO ji = 1, jpi
165                  IF ( qsr(ji,jj) /= 0._wp ) THEN
166                     oatte(ji,jj) = ( qsr_hc(ji,jj,1) / ( r1_rau0_rcp * qsr(ji,jj) ) )
167                     iatte(ji,jj) = oatte(ji,jj)
168                  ENDIF
169               END DO
170            END DO
171         ENDIF
172         !                                        ! ============================================== !
173      ELSE                                        !  Ocean alone :
174         !                                        ! ============================================== !
175         !
176         !                                                ! ------------------------- !
177         IF( ln_qsr_rgb) THEN                             !  R-G-B  light penetration !
178            !                                             ! ------------------------- !
179            ! Set chlorophyl concentration
180            IF( nn_chldta == 1 .OR. lk_vvl ) THEN            !*  Variable Chlorophyll or ocean volume
181               !
182               IF( nn_chldta == 1 ) THEN                             !- Variable Chlorophyll
183                  !
184                  CALL fld_read( kt, 1, sf_chl )                         ! Read Chl data and provides it at the current time step
185                  !         
186!CDIR COLLAPSE
187!CDIR NOVERRCHK
188                  DO jj = 1, jpj                                         ! Separation in R-G-B depending of the surface Chl
189!CDIR NOVERRCHK
190                     DO ji = 1, jpi
191                        zchl = MIN( 10. , MAX( 0.03, sf_chl(1)%fnow(ji,jj,1) ) )
192                        irgb = NINT( 41 + 20.*LOG10(zchl) + 1.e-15 )
193                        zekb(ji,jj) = rkrgb(1,irgb)
194                        zekg(ji,jj) = rkrgb(2,irgb)
195                        zekr(ji,jj) = rkrgb(3,irgb)
196                     END DO
197                  END DO
198               ELSE                                                  !- Variable ocean volume but constant chrlorophyll
199                  zchl = 0.05                                           ! constant chlorophyll
200                  irgb = NINT( 41 + 20.*LOG10( zchl ) + 1.e-15 )
201                  zekb(:,:) = rkrgb(1,irgb)                             ! Separation in R-G-B depending of the chlorophyll
202                  zekg(:,:) = rkrgb(2,irgb)
203                  zekr(:,:) = rkrgb(3,irgb)
204               ENDIF
205               !
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) * xsi0r     )
219                        zc1 = ze1(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zekb(ji,jj) )
220                        zc2 = ze2(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zekg(ji,jj) )
221                        zc3 = ze3(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zekr(ji,jj) )
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               ! clem: store attenuation coefficient of the first ocean level
231               IF ( lk_lim3 .AND. ln_qsr_ice ) THEN
232                  DO jj = 1, jpj
233                     DO ji = 1, jpi
234                        zzc0 = rn_abs * EXP( - fse3t(ji,jj,1) * xsi0r     )
235                        zzc1 = zcoef  * EXP( - fse3t(ji,jj,1) * zekb(ji,jj) )
236                        zzc2 = zcoef  * EXP( - fse3t(ji,jj,1) * zekg(ji,jj) )
237                        zzc3 = zcoef  * EXP( - fse3t(ji,jj,1) * zekr(ji,jj) )
238                        oatte(ji,jj) = 1.0 - ( zzc0 + zzc1 + zzc2  + zzc3  ) * tmask(ji,jj,2) 
239                        iatte(ji,jj) = 1.0 - ( zzc0 + zzc1 + zcoef + zcoef ) * tmask(ji,jj,2)
240                     END DO
241                  END DO
242               ENDIF
243               !
244               DO jk = 1, nksr                                        ! compute and add qsr trend to ta
245                  qsr_hc(:,:,jk) = r1_rau0_rcp * ( zea(:,:,jk) - zea(:,:,jk+1) )
246               END DO
247               zea(:,:,nksr+1:jpk) = 0.e0     ! below 400m set to zero
248               CALL iom_put( 'qsr3d', zea )   ! Shortwave Radiation 3D distribution
249               !
250            ELSE                                                 !*  Constant Chlorophyll
251               DO jk = 1, nksr
252                  qsr_hc(:,:,jk) =  etot3(:,:,jk) * qsr(:,:)
253               END DO
254               ! clem: store attenuation coefficient of the first ocean level
255               IF ( lk_lim3 .AND. ln_qsr_ice ) THEN
256               
257!!gm  BUG ??????   ? ?  ?
258                  oatte(:,:) = etot3(:,:,1) / r1_rau0_rcp
259                  iatte(:,:) = oatte(:,:)
260               ENDIF
261           ENDIF
262
263         ENDIF
264         !                                                ! ------------------------- !
265         IF( ln_qsr_2bd ) THEN                            !  2 band light penetration !
266            !                                             ! ------------------------- !
267            !
268            IF( lk_vvl ) THEN                                  !* variable volume
269               zz0   =        rn_abs   * r1_rau0_rcp
270               zz1   = ( 1. - rn_abs ) * r1_rau0_rcp
271               DO jk = 1, nksr                    ! solar heat absorbed at T-point in the top 400m
272                  DO jj = 1, jpj
273                     DO ji = 1, jpi
274                        zc0 = zz0 * EXP( -fsdepw(ji,jj,jk  )*xsi0r ) + zz1 * EXP( -fsdepw(ji,jj,jk  )*xsi1r )
275                        zc1 = zz0 * EXP( -fsdepw(ji,jj,jk+1)*xsi0r ) + zz1 * EXP( -fsdepw(ji,jj,jk+1)*xsi1r )
276                        qsr_hc(ji,jj,jk) = qsr(ji,jj) * ( zc0*tmask(ji,jj,jk) - zc1*tmask(ji,jj,jk+1) ) 
277                     END DO
278                  END DO
279               END DO
280               ! clem: store attenuation coefficient of the first ocean level
281               IF ( lk_lim3 .AND. ln_qsr_ice ) THEN
282                  DO jj = 1, jpj
283                     DO ji = 1, jpi
284                        zc0 = zz0 * EXP( -fsdepw(ji,jj,1)*xsi0r ) + zz1 * EXP( -fsdepw(ji,jj,1)*xsi1r )
285                        zc1 = zz0 * EXP( -fsdepw(ji,jj,2)*xsi0r ) + zz1 * EXP( -fsdepw(ji,jj,2)*xsi1r )
286                        oatte(ji,jj) = ( zc0*tmask(ji,jj,1) - zc1*tmask(ji,jj,2) ) / r1_rau0_rcp
287                        iatte(ji,jj) = oatte(ji,jj)
288                     END DO
289                  END DO
290               ENDIF
291            ELSE                                               !* constant volume: coef. computed one for all
292               DO jk = 1, nksr
293                  DO jj = 2, jpjm1
294                     DO ji = fs_2, fs_jpim1   ! vector opt.
295                        qsr_hc(ji,jj,jk) =  etot3(ji,jj,jk) * qsr(ji,jj)
296                     END DO
297                  END DO
298               END DO
299               ! clem: store attenuation coefficient of the first ocean level
300               IF ( lk_lim3 .AND. ln_qsr_ice ) THEN
301                  oatte(:,:) = etot3(:,:,1) / r1_rau0_rcp
302                  iatte(:,:) = oatte(:,:)
303               ENDIF
304               !
305            ENDIF
306            !
307         ENDIF
308         !
309         !                                        Add to the general trend
310         DO jk = 1, nksr
311            DO jj = 2, jpjm1 
312               DO ji = fs_2, fs_jpim1   ! vector opt.
313                  z1_e3t = zfact / fse3t(ji,jj,jk)
314                  tsa(ji,jj,jk,jp_tem) = tsa(ji,jj,jk,jp_tem) + ( qsr_hc_b(ji,jj,jk) + qsr_hc(ji,jj,jk) ) * z1_e3t
315               END DO
316            END DO
317         END DO
318         !
319      ENDIF
320      !
321      IF( lrst_oce ) THEN   !                  Write in the ocean restart file
322         !                                     *******************************
323         IF(lwp) WRITE(numout,*)
324         IF(lwp) WRITE(numout,*) 'qsr tracer content forcing field written in ocean restart file ',   &
325            &                    'at it= ', kt,' date= ', ndastp
326         IF(lwp) WRITE(numout,*) '~~~~'
327         CALL iom_rstput( kt, nitrst, numrow, 'qsr_hc_b', qsr_hc )
328         !
329      ENDIF
330
331      IF( l_trdtra ) THEN     ! qsr tracers trends saved for diagnostics
332         ztrdt(:,:,:) = tsa(:,:,:,jp_tem) - ztrdt(:,:,:)
333         CALL trd_tra( kt, 'TRA', jp_tem, jptra_qsr, ztrdt )
334         CALL wrk_dealloc( jpi, jpj, jpk, ztrdt ) 
335      ENDIF
336      !                       ! print mean trends (used for debugging)
337      IF(ln_ctl)   CALL prt_ctl( tab3d_1=tsa(:,:,:,jp_tem), clinfo1=' qsr  - Ta: ', mask1=tmask, clinfo3='tra-ta' )
338      !
339      CALL wrk_dealloc( jpi, jpj,      zekb, zekg, zekr        ) 
340      CALL wrk_dealloc( jpi, jpj, jpk, ze0, ze1, ze2, ze3, zea ) 
341      !
342      IF( nn_timing == 1 )  CALL timing_stop('tra_qsr')
343      !
344   END SUBROUTINE tra_qsr
345
346
347   SUBROUTINE tra_qsr_init
348      !!----------------------------------------------------------------------
349      !!                  ***  ROUTINE tra_qsr_init  ***
350      !!
351      !! ** Purpose :   Initialization for the penetrative solar radiation
352      !!
353      !! ** Method  :   The profile of solar radiation within the ocean is set
354      !!      from two length scale of penetration (rn_si0,rn_si1) and a ratio
355      !!      (rn_abs). These parameters are read in the namtra_qsr namelist. The
356      !!      default values correspond to clear water (type I in Jerlov'
357      !!      (1968) classification.
358      !!         called by tra_qsr at the first timestep (nit000)
359      !!
360      !! ** Action  : - initialize rn_si0, rn_si1 and rn_abs
361      !!
362      !! Reference : Jerlov, N. G., 1968 Optical Oceanography, Elsevier, 194pp.
363      !!----------------------------------------------------------------------
364      INTEGER  ::   ji, jj, jk                   ! dummy loop indices
365      INTEGER  ::   irgb, ierror, ioptio, nqsr   ! local integer
366      INTEGER  ::   ios                          ! Local integer output status for namelist read
367      REAL(wp) ::   zz0, zc0  , zc1, zcoef       ! local scalars
368      REAL(wp) ::   zz1, zc2  , zc3, zchl        !   -      -
369      REAL(wp), POINTER, DIMENSION(:,:  ) :: zekb, zekg, zekr
370      REAL(wp), POINTER, DIMENSION(:,:,:) :: ze0, ze1, ze2, ze3, zea
371      !
372      CHARACTER(len=100) ::   cn_dir   ! Root directory for location of ssr files
373      TYPE(FLD_N)        ::   sn_chl   ! informations about the chlorofyl field to be read
374      !!
375      NAMELIST/namtra_qsr/  sn_chl, cn_dir, ln_traqsr, ln_qsr_rgb, ln_qsr_2bd, ln_qsr_bio, ln_qsr_ice,  &
376         &                  nn_chldta, rn_abs, rn_si0, rn_si1
377      !!----------------------------------------------------------------------
378
379      !
380      IF( nn_timing == 1 )  CALL timing_start('tra_qsr_init')
381      !
382      ! clem init for oatte and iatte
383      IF( .NOT. ln_rstart ) THEN
384         oatte(:,:) = 1._wp
385         iatte(:,:) = 1._wp
386      ENDIF
387      !
388      CALL wrk_alloc( jpi, jpj,      zekb, zekg, zekr        ) 
389      CALL wrk_alloc( jpi, jpj, jpk, ze0, ze1, ze2, ze3, zea ) 
390      !
391
392      REWIND( numnam_ref )              ! Namelist namtra_qsr in reference namelist : Ratio and length of penetration
393      READ  ( numnam_ref, namtra_qsr, IOSTAT = ios, ERR = 901)
394901   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtra_qsr in reference namelist', lwp )
395
396      REWIND( numnam_cfg )              !  Namelist namtra_qsr in configuration namelist : Ratio and length of penetration
397      READ  ( numnam_cfg, namtra_qsr, IOSTAT = ios, ERR = 902 )
398902   IF( ios /= 0 ) CALL ctl_nam ( ios , 'namtra_qsr in configuration namelist', lwp )
399      WRITE ( numond, namtra_qsr )
400      !
401      IF(lwp) THEN                ! control print
402         WRITE(numout,*)
403         WRITE(numout,*) 'tra_qsr_init : penetration of the surface solar radiation'
404         WRITE(numout,*) '~~~~~~~~~~~~'
405         WRITE(numout,*) '   Namelist namtra_qsr : set the parameter of penetration'
406         WRITE(numout,*) '      Light penetration (T) or not (F)         ln_traqsr  = ', ln_traqsr
407         WRITE(numout,*) '      RGB (Red-Green-Blue) light penetration   ln_qsr_rgb = ', ln_qsr_rgb
408         WRITE(numout,*) '      2 band               light penetration   ln_qsr_2bd = ', ln_qsr_2bd
409         WRITE(numout,*) '      bio-model            light penetration   ln_qsr_bio = ', ln_qsr_bio
410         WRITE(numout,*) '      light penetration for ice-model LIM3     ln_qsr_ice = ', ln_qsr_ice
411         WRITE(numout,*) '      RGB : Chl data (=1) or cst value (=0)    nn_chldta  = ', nn_chldta
412         WRITE(numout,*) '      RGB & 2 bands: fraction of light (rn_si1)    rn_abs = ', rn_abs
413         WRITE(numout,*) '      RGB & 2 bands: shortess depth of extinction  rn_si0 = ', rn_si0
414         WRITE(numout,*) '      2 bands: longest depth of extinction         rn_si1 = ', rn_si1
415         WRITE(numout,*) '      light penetration for ice-model LIM3     ln_qsr_ice = ', ln_qsr_ice   
416      ENDIF
417
418      IF( ln_traqsr ) THEN     ! control consistency
419         !                     
420         IF( .NOT.lk_qsr_bio .AND. ln_qsr_bio )   THEN
421            CALL ctl_warn( 'No bio model : force ln_qsr_bio = FALSE ' )
422            ln_qsr_bio = .FALSE.
423         ENDIF
424         !
425         ioptio = 0                      ! Parameter control
426         IF( ln_qsr_rgb  )   ioptio = ioptio + 1
427         IF( ln_qsr_2bd  )   ioptio = ioptio + 1
428         IF( ln_qsr_bio  )   ioptio = ioptio + 1
429         !
430         IF( ioptio /= 1 ) &
431            CALL ctl_stop( '          Choose ONE type of light penetration in namelist namtra_qsr',  &
432            &              ' 2 bands, 3 RGB bands or bio-model light penetration' )
433         !
434         IF( ln_qsr_rgb .AND. nn_chldta == 0 )   nqsr =  1 
435         IF( ln_qsr_rgb .AND. nn_chldta == 1 )   nqsr =  2
436         IF( ln_qsr_2bd                      )   nqsr =  3
437         IF( ln_qsr_bio                      )   nqsr =  4
438         !
439         IF(lwp) THEN                   ! Print the choice
440            WRITE(numout,*)
441            IF( nqsr ==  1 )   WRITE(numout,*) '         R-G-B   light penetration - Constant Chlorophyll'
442            IF( nqsr ==  2 )   WRITE(numout,*) '         R-G-B   light penetration - Chl data '
443            IF( nqsr ==  3 )   WRITE(numout,*) '         2 bands light penetration'
444            IF( nqsr ==  4 )   WRITE(numout,*) '         bio-model light penetration'
445         ENDIF
446         !
447      ENDIF
448      !                          ! ===================================== !
449      IF( ln_traqsr  ) THEN      !  Initialisation of Light Penetration  ! 
450         !                       ! ===================================== !
451         !
452         xsi0r = 1.e0 / rn_si0
453         xsi1r = 1.e0 / rn_si1
454         !                                ! ---------------------------------- !
455         IF( ln_qsr_rgb ) THEN            !  Red-Green-Blue light penetration  !
456            !                             ! ---------------------------------- !
457            !
458            CALL trc_oce_rgb( rkrgb )           !* tabulated attenuation coef.
459            !
460            !                                   !* level of light extinction
461            IF(  ln_sco ) THEN   ;   nksr = jpkm1
462            ELSE                 ;   nksr = trc_oce_ext_lev( r_si2, 0.33e2 )
463            ENDIF
464
465            IF(lwp) WRITE(numout,*) '        level of light extinction = ', nksr, ' ref depth = ', gdepw_1d(nksr+1), ' m'
466            !
467            IF( nn_chldta == 1 ) THEN           !* Chl data : set sf_chl structure
468               IF(lwp) WRITE(numout,*)
469               IF(lwp) WRITE(numout,*) '        Chlorophyll read in a file'
470               ALLOCATE( sf_chl(1), STAT=ierror )
471               IF( ierror > 0 ) THEN
472                  CALL ctl_stop( 'tra_qsr_init: unable to allocate sf_chl structure' )   ;   RETURN
473               ENDIF
474               ALLOCATE( sf_chl(1)%fnow(jpi,jpj,1)   )
475               IF( sn_chl%ln_tint )ALLOCATE( sf_chl(1)%fdta(jpi,jpj,1,2) )
476               !                                        ! fill sf_chl with sn_chl and control print
477               CALL fld_fill( sf_chl, (/ sn_chl /), cn_dir, 'tra_qsr_init',   &
478                  &                                         'Solar penetration function of read chlorophyll', 'namtra_qsr' )
479               !
480            ELSE                                !* constant Chl : compute once for all the distribution of light (etot3)
481               IF(lwp) WRITE(numout,*)
482               IF(lwp) WRITE(numout,*) '        Constant Chlorophyll concentration = 0.05'
483               IF( lk_vvl ) THEN                   ! variable volume
484                  IF(lwp) WRITE(numout,*) '        key_vvl: light distribution will be computed at each time step'
485               ELSE                                ! constant volume: computes one for all
486                  IF(lwp) WRITE(numout,*) '        fixed volume: light distribution computed one for all'
487                  !
488                  zchl = 0.05                                 ! constant chlorophyll
489                  irgb = NINT( 41 + 20.*LOG10(zchl) + 1.e-15 )
490                  zekb(:,:) = rkrgb(1,irgb)                   ! Separation in R-G-B depending of the chlorophyll
491                  zekg(:,:) = rkrgb(2,irgb)
492                  zekr(:,:) = rkrgb(3,irgb)
493                  !
494                  zcoef = ( 1. - rn_abs ) / 3.e0              ! equi-partition in R-G-B
495                  ze0(:,:,1) = rn_abs
496                  ze1(:,:,1) = zcoef
497                  ze2(:,:,1) = zcoef 
498                  ze3(:,:,1) = zcoef
499                  zea(:,:,1) = tmask(:,:,1)                   ! = ( ze0+ze1+z2+ze3 ) * tmask
500               
501                  DO jk = 2, nksr+1
502!CDIR NOVERRCHK
503                     DO jj = 1, jpj
504!CDIR NOVERRCHK   
505                        DO ji = 1, jpi
506                           zc0 = ze0(ji,jj,jk-1) * EXP( - e3t_0(ji,jj,jk-1) * xsi0r     )
507                           zc1 = ze1(ji,jj,jk-1) * EXP( - e3t_0(ji,jj,jk-1) * zekb(ji,jj) )
508                           zc2 = ze2(ji,jj,jk-1) * EXP( - e3t_0(ji,jj,jk-1) * zekg(ji,jj) )
509                           zc3 = ze3(ji,jj,jk-1) * EXP( - e3t_0(ji,jj,jk-1) * zekr(ji,jj) )
510                           ze0(ji,jj,jk) = zc0
511                           ze1(ji,jj,jk) = zc1
512                           ze2(ji,jj,jk) = zc2
513                           ze3(ji,jj,jk) = zc3
514                           zea(ji,jj,jk) = ( zc0 + zc1 + zc2 + zc3 ) * tmask(ji,jj,jk)
515                        END DO
516                     END DO
517                  END DO 
518                  !
519                  DO jk = 1, nksr
520                     etot3(:,:,jk) = r1_rau0_rcp * ( zea(:,:,jk) - zea(:,:,jk+1) ) 
521                  END DO
522                  etot3(:,:,nksr+1:jpk) = 0.e0                ! below 400m set to zero
523               ENDIF
524            ENDIF
525            !
526         ENDIF
527            !                             ! ---------------------------------- !
528         IF( ln_qsr_2bd ) THEN            !    2 bands    light penetration    !
529            !                             ! ---------------------------------- !
530            !
531            !                                ! level of light extinction
532            nksr = trc_oce_ext_lev( rn_si1, 1.e2 )
533            IF(lwp) THEN
534               WRITE(numout,*)
535            IF(lwp) WRITE(numout,*) '        level of light extinction = ', nksr, ' ref depth = ', gdepw_1d(nksr+1), ' m'
536            ENDIF
537            !
538            IF( lk_vvl ) THEN                   ! variable volume
539               IF(lwp) WRITE(numout,*) '        key_vvl: light distribution will be computed at each time step'
540            ELSE                                ! constant volume: computes one for all
541               zz0 =        rn_abs   * r1_rau0_rcp
542               zz1 = ( 1. - rn_abs ) * r1_rau0_rcp
543               DO jk = 1, nksr                    !*  solar heat absorbed at T-point computed once for all
544                  DO jj = 1, jpj                              ! top 400 meters
545                     DO ji = 1, jpi
546                        zc0 = zz0 * EXP( -fsdepw(ji,jj,jk  )*xsi0r ) + zz1 * EXP( -fsdepw(ji,jj,jk  )*xsi1r )
547                        zc1 = zz0 * EXP( -fsdepw(ji,jj,jk+1)*xsi0r ) + zz1 * EXP( -fsdepw(ji,jj,jk+1)*xsi1r )
548                        etot3(ji,jj,jk) = (  zc0 * tmask(ji,jj,jk) - zc1 * tmask(ji,jj,jk+1)  ) 
549                     END DO
550                  END DO
551               END DO
552               etot3(:,:,nksr+1:jpk) = 0.e0                   ! below 400m set to zero
553               !
554            ENDIF
555         ENDIF
556         !                       ! ===================================== !
557      ELSE                       !        No light penetration           !                   
558         !                       ! ===================================== !
559         IF(lwp) THEN
560            WRITE(numout,*)
561            WRITE(numout,*) 'tra_qsr_init : NO solar flux penetration'
562            WRITE(numout,*) '~~~~~~~~~~~~'
563         ENDIF
564      ENDIF
565      !
566      CALL wrk_dealloc( jpi, jpj,      zekb, zekg, zekr        ) 
567      CALL wrk_dealloc( jpi, jpj, jpk, ze0, ze1, ze2, ze3, zea ) 
568      !
569      IF( nn_timing == 1 )  CALL timing_stop('tra_qsr_init')
570      !
571   END SUBROUTINE tra_qsr_init
572
573   !!======================================================================
574END MODULE traqsr
Note: See TracBrowser for help on using the repository browser.