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

source: trunk/NEMO/OPA_SRC/TRA/traqsr.F90 @ 1423

Last change on this file since 1423 was 1423, checked in by ctlod, 15 years ago

add light penetration following 3 wavebands model (RGB) and the use of ocean color (chlorophyll), see ticket: #428

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