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 @ 5866

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

#1613: vvl by default: add ln_linssh and remove key_vvl

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