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/UKMO/dev_r5518_GO6_starthour_obsoper/NEMOGCM/NEMO/OPA_SRC/TRA – NEMO

source: branches/UKMO/dev_r5518_GO6_starthour_obsoper/NEMOGCM/NEMO/OPA_SRC/TRA/traqsr.F90 @ 15603

Last change on this file since 15603 was 15603, checked in by mattmartin, 3 years ago

Updated NEMO branch for coupled NWP at GO6 to include stochastic model perturbations.
For more info see ticket: https://code.metoffice.gov.uk/trac/nwpscience/ticket/1125.

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