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

Last change on this file since 1445 was 1445, checked in by cetlod, 15 years ago

add the use of bio-optical retroaction on dynamics when coupling with PISCES, see ticket:428

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