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

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

fix a small compilation error related to the implementation of RGB, see ticket: #428

  • Property svn:eol-style set to native
  • Property svn:keywords set to Id
File size: 19.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   !!----------------------------------------------------------------------
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                     irgb = NINT( 41 + 20.*LOG10(zchl) + 1.e-15 )
141                     zekb(ji,jj) = rkrgb(1,irgb)
142                     zekg(ji,jj) = rkrgb(2,irgb)
143                     zekr(ji,jj) = rkrgb(3,irgb)
144                  END DO
145               END DO
146               !
147               zsi0r = 1.e0 / rn_si0
148               zcoef  = ( 1. - rn_abs ) / 3.e0                        ! equi-partition in R-G-B
149               ze0(:,:,1) = rn_abs  * qsr(:,:)
150               ze1(:,:,1) = zcoef * qsr(:,:)
151               ze2(:,:,1) = zcoef * qsr(:,:)
152               ze3(:,:,1) = zcoef * qsr(:,:)
153               zea(:,:,1) =         qsr(:,:)
154               !
155               DO jk = 2, nksr+1
156!CDIR NOVERRCHK
157                  DO jj = 1, jpj
158!CDIR NOVERRCHK   
159                     DO ji = 1, jpi
160                        zc0 = ze0(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zsi0r     )
161                        zc1 = ze1(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zekb(ji,jj) )
162                        zc2 = ze2(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zekg(ji,jj) )
163                        zc3 = ze3(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zekr(ji,jj) )
164                        ze0(ji,jj,jk) = zc0
165                        ze1(ji,jj,jk) = zc1
166                        ze2(ji,jj,jk) = zc2
167                        ze3(ji,jj,jk) = zc3
168                        zea(ji,jj,jk) = ( zc0 + zc1 + zc2 + zc3 ) * tmask(ji,jj,jk)
169                     END DO
170                  END DO
171               END DO
172               !
173               DO jk = 1, nksr                                        ! compute and add qsr trend to ta
174                  ta(:,:,jk) = ta(:,:,jk) + ro0cpr * ( zea(:,:,jk) - zea(:,:,jk+1) ) / fse3t(:,:,jk)
175               END DO
176               !
177            ELSE                                                 !*  Constant Chlorophyll
178               DO jk = 1, nksr
179                  ta(:,:,jk) = ta(:,:,jk) + etot3(:,:,jk) * qsr(:,:)
180               END DO
181            ENDIF
182
183!!gm BUG  the case key_vvl is missing: etot3 should be recomputed at each time step !!!
184
185            !                                             ! ------------------------- !
186         ELSE                                             !  2 band light penetration !
187            !                                             ! ------------------------- !
188            !
189            DO jk = 1, nksr
190               DO jj = 2, jpjm1
191                  DO ji = fs_2, fs_jpim1   ! vector opt.
192                     ta(ji,jj,jk) = ta(ji,jj,jk) + etot3(ji,jj,jk) * qsr(ji,jj)
193                  END DO
194               END DO
195            END DO
196            !
197         ENDIF
198
199!!gm BUG  the case key_vvl is missing: etot3 should be recomputed at each time step !!!
200
201         !
202      ENDIF
203
204      IF( l_trdtra ) THEN     ! qsr tracers trends saved for diagnostics
205         ztrdt(:,:,:) = ta(:,:,:) - ztrdt(:,:,:)
206         CALL trd_mod( ztrdt, ztrds, jptra_trd_qsr, 'TRA', kt )
207      ENDIF
208      !                       ! print mean trends (used for debugging)
209      IF(ln_ctl)   CALL prt_ctl( tab3d_1=ta, clinfo1=' qsr  - Ta: ', mask1=tmask, clinfo3='tra-ta' )
210      !
211   END SUBROUTINE tra_qsr
212
213
214   SUBROUTINE tra_qsr_init
215      !!----------------------------------------------------------------------
216      !!                  ***  ROUTINE tra_qsr_init  ***
217      !!
218      !! ** Purpose :   Initialization for the penetrative solar radiation
219      !!
220      !! ** Method  :   The profile of solar radiation within the ocean is set
221      !!      from two length scale of penetration (rn_si0,rn_si1) and a ratio
222      !!      (rn_abs). These parameters are read in the namqsr namelist. The
223      !!      default values correspond to clear water (type I in Jerlov'
224      !!      (1968) classification.
225      !!         called by tra_qsr at the first timestep (nit000)
226      !!
227      !! ** Action  : - initialize rn_si0, rn_si1 and rn_abs
228      !!
229      !! Reference : Jerlov, N. G., 1968 Optical Oceanography, Elsevier, 194pp.
230      !!----------------------------------------------------------------------
231      INTEGER  ::   ji, jj, jk     ! dummy loop indices
232      INTEGER  ::   irgb, ierror   ! temporary integer
233      REAL(wp) ::   zc0  , zc1            ! temporary scalars
234      REAL(wp) ::   zc2  , zc3  , zchl    !    -         -
235      REAL(wp) ::   zsi0r, zsi1r, zcoef   !    -         -
236      REAL(wp), DIMENSION(jpi,jpj)     ::   zekb, zekg, zekr              ! 2D workspace
237      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   ze0 , ze1 , ze2 , ze3 , zea   ! 3D workspace
238      !!
239      CHARACTER(len=100) ::   cn_dir   ! Root directory for location of ssr files
240      TYPE(FLD_N)        ::   sn_chl   ! informations about the chlorofyl field to be read
241      NAMELIST/namqsr/  sn_chl, cn_dir, ln_traqsr, ln_qsr_rgb, ln_qsr_bio,   &
242         &              nn_chldta, rn_abs, rn_si0, rn_si1, rn_si2
243      !!----------------------------------------------------------------------
244
245      cn_dir = './'       ! directory in which the model is executed
246      ! ... default values (NB: frequency positive => hours, negative => months)
247      !            !     file       ! frequency !  variable  ! time interp !  clim   ! 'yearly' or ! weights  ! rotation   !
248      !            !     name       !  (hours)  !    name    !    (T/F)    !  (T/F)  ! 'monthly'   ! filename ! pairs      !
249      sn_chl = FLD_N( 'chlorophyll' ,    -1     ,  'CHLA'    ,  .true.     , .true.  ,   'yearly'  , ''       , ''         )
250      !
251      REWIND( numnam )            ! Read Namelist namqsr : ratio and length of penetration
252      READ  ( numnam, namqsr )
253      !
254      IF(lwp) THEN                ! control print
255         WRITE(numout,*)
256         WRITE(numout,*) 'tra_qsr_init : penetration of the surface solar radiation'
257         WRITE(numout,*) '~~~~~~~~~~~~'
258         WRITE(numout,*) '    Namelist namqsr : set the parameter of penetration'
259         WRITE(numout,*) '        Light penetration (T) or not (F)         ln_traqsr  = ', ln_traqsr
260         WRITE(numout,*) '        RGB (Red-Green-Blue) light penetration   ln_qsr_rgb = ', ln_qsr_rgb
261         WRITE(numout,*) '        bio-model            light penetration   ln_qsr_bio = ', ln_qsr_bio
262         WRITE(numout,*) '        RGB : Chl data (=1) or cst value (=0)    nn_chldta  = ', nn_chldta
263         WRITE(numout,*) '        RGB & 2 bands: fraction of light (rn_si1)    rn_abs = ', rn_abs
264         WRITE(numout,*) '        RGB & 2 bands: shortess depth of extinction  rn_si0 = ', rn_si0
265         WRITE(numout,*) '        2 bands: longest depth of extinction         rn_si1 = ', rn_si1
266         WRITE(numout,*) '        3 bands: longest depth of extinction         rn_si2 = ', rn_si2
267      ENDIF
268      !                         ! control consistency
269      IF( lk_qsr_bio .AND. .NOT.ln_qsr_bio ) THEN
270         ln_qsr_bio = .true.
271         CALL ctl_warn( 'Force bio-model light penetraton ln_qsr_bio  = TRUE ' )
272      ENDIF
273     
274      !                          ! ===================================== !
275      IF( ln_traqsr  ) THEN      !  Initialisation of Light Penetration  ! 
276         !                       ! ===================================== !
277         !
278         zsi0r = 1.e0 / rn_si0
279         zsi1r = 1.e0 / rn_si1
280         !                                ! ---------------------------------- !
281         IF( ln_qsr_rgb ) THEN            !  Red-Green-Blue light penetration  !
282            !                             ! ---------------------------------- !
283            !
284            !                                ! level of light extinction
285            nksr = trc_oce_ext_lev( rn_si2, 0.33e2 )
286            IF(lwp) THEN
287               WRITE(numout,*)
288               WRITE(numout,*) '        level max of computation of qsr = ', nksr, ' ref depth = ', gdepw_0(nksr+1), ' m'
289            ENDIF
290            !
291            CALL trc_oce_rgb( rkrgb )           !* tabulated attenuation coef.
292!!gm            CALL trc_oce_rgb_read( rkrgb )           !* tabulated attenuation coef.
293            !
294            IF( nn_chldta == 1 ) THEN           !* Chl data : set sf_chl structure
295               IF(lwp) WRITE(numout,*)
296               IF(lwp) WRITE(numout,*) '        Chlorophyll read in a file'
297               ALLOCATE( sf_chl(1), STAT=ierror )
298               IF( ierror > 0 ) THEN
299                  CALL ctl_stop( 'tra_qsr_init: unable to allocate sf_chl structure' )   ;   RETURN
300               ENDIF
301               ALLOCATE( sf_chl(1)%fnow(jpi,jpj)   )
302               ALLOCATE( sf_chl(1)%fdta(jpi,jpj,2) )
303               !                                        ! fill sf_chl with sn_chl and control print
304               CALL fld_fill( sf_chl, (/ sn_chl /), cn_dir, 'tra_qsr_init',   &
305                  &                                         'Solar penetration function of read chlorophyll', 'namqsr' )
306               !
307            ELSE                                !* constant Chl : compute once for all the distribution of light (etot3)
308               IF(lwp) WRITE(numout,*)
309               IF(lwp) WRITE(numout,*) '        Constant Chlorophyll concentration = 0.05'
310               IF(lwp) WRITE(numout,*) '        light distribution computed once for all'
311               !
312               zchl = 0.05                                 ! constant chlorophyll
313               irgb = NINT( 41 + 20.*LOG10(zchl) + 1.e-15 )
314               zekb(:,:) = rkrgb(1,irgb)                   ! Separation in R-G-B depending of the chlorophyl concentration
315               zekg(:,:) = rkrgb(2,irgb)
316               zekr(:,:) = rkrgb(3,irgb)
317               !
318               zcoef = ( 1. - rn_abs ) / 3.e0              ! equi-partition in R-G-B
319               ze0(:,:,1) = rn_abs
320               ze1(:,:,1) = zcoef
321               ze2(:,:,1) = zcoef 
322               ze3(:,:,1) = zcoef
323               zea(:,:,1) = tmask(:,:,1)                   ! = ( ze0+ze1+z2+ze3 ) * tmask
324               
325               DO jk = 2, nksr+1
326!CDIR NOVERRCHK
327                  DO jj = 1, jpj
328!CDIR NOVERRCHK   
329                     DO ji = 1, jpi
330                        zc0 = ze0(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zsi0r     )
331                        zc1 = ze1(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zekb(ji,jj) )
332                        zc2 = ze2(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zekg(ji,jj) )
333                        zc3 = ze3(ji,jj,jk-1) * EXP( - fse3t(ji,jj,jk-1) * zekr(ji,jj) )
334                        ze0(ji,jj,jk) = zc0                 
335                        ze1(ji,jj,jk) = zc1
336                        ze2(ji,jj,jk) = zc2     
337                        ze3(ji,jj,jk) = zc3
338                        zea(ji,jj,jk) = ( zc0 + zc1 + zc2 + zc3 ) * tmask(ji,jj,jk)
339                     END DO
340                  END DO
341               END DO 
342               !
343               DO jk = 1, nksr
344                  etot3(:,:,jk) = ro0cpr * ( zea(:,:,jk) - zea(:,:,jk+1) ) / fse3t(:,:,jk)
345               END DO
346               etot3(:,:,nksr+1:jpk) = 0.e0                ! below 400m set to zero
347            ENDIF
348            !
349            !                             ! ---------------------------------- !
350         ELSE                             !    2 bands    light penetration    !
351            !                             ! ---------------------------------- !
352            !
353            !                                ! level of light extinction
354            nksr = trc_oce_ext_lev( rn_si1, 1.e2 )
355            IF(lwp) THEN
356               WRITE(numout,*)
357               WRITE(numout,*) '        level max of computation of qsr = ', nksr, ' ref depth = ', gdepw_0(nksr+1), ' m'
358            ENDIF
359            !
360            DO jk = 1, nksr                    !*  solar heat absorbed at T-point computed once for all
361               DO jj = 1, jpj                              ! top 400 meters
362                  DO ji = 1, jpi
363                     zc0 = rn_abs * EXP( -fsdepw(ji,jj,jk  )*zsi0r ) + (1.-rn_abs) * EXP( -fsdepw(ji,jj,jk  )*zsi1r )
364                     zc1 = rn_abs * EXP( -fsdepw(ji,jj,jk+1)*zsi0r ) + (1.-rn_abs) * EXP( -fsdepw(ji,jj,jk+1)*zsi1r )
365                     etot3(ji,jj,jk) = ro0cpr * (  zc0 * tmask(ji,jj,jk) - zc1 * tmask(ji,jj,jk+1)  ) / fse3t(ji,jj,jk)
366                  END DO
367               END DO
368            END DO
369            etot3(:,:,nksr+1:jpk) = 0.e0                   ! below 400m set to zero
370            !
371         ENDIF
372         !                       ! ===================================== !
373      ELSE                       !        No light penetration           !                   
374         !                       ! ===================================== !
375         IF(lwp) THEN
376            WRITE(numout,*)
377            WRITE(numout,*) 'tra_qsr_init : NO solar flux penetration'
378            WRITE(numout,*) '~~~~~~~~~~~~'
379         ENDIF
380      ENDIF
381      !
382   END SUBROUTINE tra_qsr_init
383
384   !!======================================================================
385END MODULE traqsr
Note: See TracBrowser for help on using the repository browser.