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.
dtadyn.F90 in branches/dev_r2586_dynamic_mem/NEMOGCM/NEMO/OFF_SRC – NEMO

source: branches/dev_r2586_dynamic_mem/NEMOGCM/NEMO/OFF_SRC/dtadyn.F90 @ 2648

Last change on this file since 2648 was 2648, checked in by cetlod, 13 years ago

Changed OFF_SRC component to use dynamic memory

  • Property svn:keywords set to Id
File size: 36.6 KB
Line 
1MODULE dtadyn
2   !!======================================================================
3   !!                       ***  MODULE  dtadyn  ***
4   !! Off-line : interpolation of the physical fields
5   !!======================================================================
6   !! History :   OPA  ! 1992-01 (M. Imbard) Original code
7   !!             8.0  ! 1998-04 (L.Bopp MA Foujols) slopes for isopyc.
8   !!              -   ! 1998-05 (L. Bopp) read output of coupled run
9   !!             8.2  ! 2001-01 (M. Levy et M. Benjelloul) add netcdf FORMAT
10   !!   NEMO      1.0  ! 2005-03 (O. Aumont and A. El Moussaoui) F90
11   !!              -   ! 2005-12 (C. Ethe) Adapted for DEGINT
12   !!             3.0  ! 2007-06 (C. Ethe) use of iom module
13   !!              -   ! 2007-09  (C. Ethe)  add swap_dyn_data
14   !!             3.3  ! 2010-11 (C. Ethe) Full reorganization of the off-line: phasing with the on-line
15   !!----------------------------------------------------------------------
16
17   !!----------------------------------------------------------------------
18   !!   dta_dyn_init : initialization, namelist read, and parameters control
19   !!   dta_dyn      : Interpolation of the fields
20   !!----------------------------------------------------------------------
21   USE oce             ! ocean dynamics and tracers variables
22   USE c1d             ! 1D configuration: lk_c1d
23   USE dom_oce         ! ocean domain: variables
24   USE zdf_oce         ! ocean vertical physics: variables
25   USE sbc_oce         ! surface module: variables
26   USE phycst          ! physical constants
27   USE trabbl          ! active tracer: bottom boundary layer
28   USE ldfslp          ! lateral diffusion: iso-neutral slopes
29   USE ldfeiv          ! eddy induced velocity coef.
30   USE ldftra_oce      ! ocean tracer   lateral physics
31   USE zdfmxl          ! vertical physics: mixed layer depth
32   USE eosbn2          ! equation of state - Brunt Vaisala frequency
33   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
34   USE zpshde          ! z-coord. with partial steps: horizontal derivatives
35   USE in_out_manager  ! I/O manager
36   USE iom             ! I/O library
37   USE lib_mpp         ! distributed memory computing library
38   USE prtctl          !  print control
39
40   IMPLICIT NONE
41   PRIVATE
42
43   PUBLIC   dta_dyn_init   ! called by opa.F90
44   PUBLIC   dta_dyn        ! called by step.F90
45
46   LOGICAL, PUBLIC ::   lperdyn = .TRUE.   !: boolean for periodic fields or not
47   LOGICAL, PUBLIC ::   lfirdyn = .TRUE.   !: boolean for the first call or not
48
49   INTEGER, PUBLIC ::   ndtadyn = 73       !: Number of dat in one year
50   INTEGER, PUBLIC ::   ndtatot = 73       !: Number of data in the input field
51   INTEGER, PUBLIC ::   nsptint = 1        !: type of spatial interpolation
52
53   CHARACTER(len=45) ::   cfile_grid_T = 'dyna_grid_T.nc'   ! name of the grid_T file
54   CHARACTER(len=45) ::   cfile_grid_U = 'dyna_grid_U.nc'   ! name of the grid_U file
55   CHARACTER(len=45) ::   cfile_grid_V = 'dyna_grid_V.nc'   ! name of the grid_V file
56   CHARACTER(len=45) ::   cfile_grid_W = 'dyna_grid_W.nc'   ! name of the grid_W file
57   
58   REAL(wp) ::   rnspdta    ! number of time step per 2 consecutives data
59   REAL(wp) ::   rnspdta2   ! rnspdta * 0.5
60
61   INTEGER ::   ndyn1, ndyn2    !
62   INTEGER ::   nlecoff = 0     ! switch for the first read
63   INTEGER ::   numfl_t, numfl_u, numfl_v, numfl_w
64
65   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: tdta       ! temperature at two consecutive times
66   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: sdta       ! salinity at two consecutive times
67   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: udta       ! zonal velocity at two consecutive times
68   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: vdta       ! meridional velocity at two consecutive times
69   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: wdta       ! vertical velocity at two consecutive times
70   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: avtdta     ! vertical diffusivity coefficient
71
72   REAL(wp), ALLOCATABLE, SAVE,   DIMENSION(:,:,:) :: hmlddta    ! mixed layer depth at two consecutive times
73   REAL(wp), ALLOCATABLE, SAVE,   DIMENSION(:,:,:) :: wspddta    ! wind speed at two consecutive times
74   REAL(wp), ALLOCATABLE, SAVE,   DIMENSION(:,:,:) :: frlddta    ! sea-ice fraction at two consecutive times
75   REAL(wp), ALLOCATABLE, SAVE,   DIMENSION(:,:,:) :: empdta     ! E-P at two consecutive times
76   REAL(wp), ALLOCATABLE, SAVE,   DIMENSION(:,:,:) :: qsrdta     ! short wave heat flux at two consecutive times
77   REAL(wp), ALLOCATABLE, SAVE,   DIMENSION(:,:,:) :: bblxdta    ! frequency of bbl in the x direction at 2 consecutive times
78   REAL(wp), ALLOCATABLE, SAVE,   DIMENSION(:,:,:) :: bblydta    ! frequency of bbl in the y direction at 2 consecutive times
79   LOGICAL :: l_offbbl
80#if defined key_ldfslp
81   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: uslpdta    ! zonal isopycnal slopes
82   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: vslpdta    ! meridional isopycnal slopes
83   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: wslpidta   ! zonal diapycnal slopes
84   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: wslpjdta   ! meridional diapycnal slopes
85#endif
86#if ! defined key_degrad &&  defined key_traldf_c2d && defined key_traldf_eiv
87   REAL(wp), ALLOCATABLE, SAVE,   DIMENSION(:,:,:) :: aeiwdta    ! G&M coefficient
88#endif
89#if defined key_degrad
90   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: ahtudta, ahtvdta, ahtwdta   ! Lateral diffusivity
91# if defined key_traldf_eiv
92   REAL(wp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: aeiudta, aeivdta, aeiwdta   ! G&M coefficient
93# endif
94#endif
95
96   !! * Substitutions
97#  include "domzgr_substitute.h90"
98#  include "vectopt_loop_substitute.h90"
99   !!----------------------------------------------------------------------
100   !! NEMO/OFF 3.3 , NEMO Consortium (2010)
101   !! $Id$
102   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
103   !!----------------------------------------------------------------------
104CONTAINS
105
106   SUBROUTINE dta_dyn( kt )
107      !!----------------------------------------------------------------------
108      !!                  ***  ROUTINE dta_dyn  ***
109      !!
110      !! ** Purpose :   Prepares dynamics and physics fields from an NEMO run
111      !!              for an off-line simulation of passive tracers
112      !!
113      !! ** Method : calculates the position of DATA to read READ DATA
114      !!             (example month changement) computes slopes IF needed
115      !!             interpolates DATA IF needed
116      !!----------------------------------------------------------------------
117      INTEGER, INTENT(in) ::   kt   ! ocean time-step index
118      !!
119      INTEGER  ::   iper, iperm1, iswap, izt   ! local integers
120      REAL(wp) ::   zt, zweigh                 ! local scalars
121      !!----------------------------------------------------------------------
122
123      zt     = ( REAL(kt,wp) + rnspdta2 ) / rnspdta
124      izt    = INT( zt )
125      zweigh = zt - REAL( INT(zt), wp )
126
127      IF( lperdyn ) THEN   ;   iperm1 = MOD( izt, ndtadyn )
128      ELSE                 ;   iperm1 = MOD( izt, ndtatot - 1 ) + 1
129      ENDIF
130
131      iper = iperm1 + 1
132      IF( iperm1 == 0 ) THEN
133          IF( lperdyn ) THEN
134              iperm1 = ndtadyn
135          ELSE
136              IF( lfirdyn ) THEN
137                  IF(lwp) WRITE (numout,*) 'dta_dyn:  dynamic file is not periodic with or without interpolation    &
138                     &                                we take the first value for the previous period iperm1 = 0  '
139              END IF
140          END IF
141      END IF
142
143      iswap  = 0
144
145      ! 1. First call lfirdyn = true
146      ! ----------------------------
147
148      IF( lfirdyn ) THEN
149         ndyn1 = iperm1         ! store the information of the period read
150         ndyn2 = iper
151         
152         IF(lwp) THEN
153            WRITE (numout,*) ' dynamics data read for the period ndyn1 =', ndyn1,   &
154               &             ' and for the period ndyn2 = ', ndyn2
155            WRITE (numout,*) ' time step is : ', kt
156            WRITE (numout,*) ' we have ndtadyn = ', ndtadyn, ' records in the dynamic file for one year'
157         END IF
158         !
159         CALL dynrea( kt, MAX( 1, iperm1) )           ! data read for the iperm1 period
160         
161         IF( lk_ldfslp .AND. .NOT. lk_c1d ) THEN      ! Computes slopes (here tsn and avt are used as workspace)
162            tsn (:,:,:,jp_tem) = tdta  (:,:,:,2)
163            tsn (:,:,:,jp_sal) = sdta  (:,:,:,2)
164            avt(:,:,:)         = avtdta(:,:,:,2)
165         
166            CALL eos( tsn, rhd, rhop )   ! Time-filtered in situ density
167            CALL bn2( tsn, rn2 )         ! before Brunt-Vaisala frequency
168            IF( ln_zps )   &
169               &   CALL zps_hde( kt, jpts, tsn, gtsu, gtsv,  &  ! Partial steps: before Horizontal DErivative
170               &                           rhd, gru , grv   )    ! of t, s, rd at the bottom ocean level
171            CALL zdf_mxl( kt )           ! mixed layer depth
172            CALL ldf_slp( kt, rhd, rn2 )
173         
174            uslpdta (:,:,:,2) = uslp (:,:,:)
175            vslpdta (:,:,:,2) = vslp (:,:,:)
176            wslpidta(:,:,:,2) = wslpi(:,:,:)
177            wslpjdta(:,:,:,2) = wslpj(:,:,:)
178         END IF
179         !
180         CALL swap_dyn_data            ! swap from record 2 to 1
181         !
182         iswap = 1        !  indicates swap
183         !
184         CALL dynrea( kt, iper )       ! data read for the iper period
185         !
186         IF( lk_ldfslp .AND. .NOT. lk_c1d ) THEN      ! Computes slopes (here tsn and avt are used as workspace)
187            tsn (:,:,:,jp_tem) = tdta  (:,:,:,2)
188            tsn (:,:,:,jp_sal) = sdta  (:,:,:,2)
189            avt(:,:,:)         = avtdta(:,:,:,2)
190            !
191                           CALL eos( tsn, rhd, rhop )                   ! now in situ density
192                           CALL bn2( tsn, rn2 )                         ! now Brunt-Vaisala frequency
193            IF( ln_zps )   CALL zps_hde( kt, jpts, tsn, gtsu, gtsv,  &  ! Partial steps: before Horizontal DErivative
194               &                                   rhd, gru , grv   )   ! of t, s, rd at the bottom ocean level
195                           CALL zdf_mxl( kt )                           ! mixed layer depth
196                           CALL ldf_slp( kt, rhd, rn2 )                 ! slope of iso-neutral surfaces
197            !
198            uslpdta (:,:,:,2) = uslp (:,:,:)
199            vslpdta (:,:,:,2) = vslp (:,:,:)
200            wslpidta(:,:,:,2) = wslpi(:,:,:)
201            wslpjdta(:,:,:,2) = wslpj(:,:,:)
202         END IF
203         !
204         lfirdyn = .FALSE.    ! trace the first call
205      ENDIF
206      !
207      ! And now what we have to do at every time step
208      ! check the validity of the period in memory
209      !
210      IF( iperm1 /= ndyn1 ) THEN 
211         !
212         IF( iperm1 == 0 ) THEN
213            IF(lwp) THEN
214               WRITE (numout,*) ' dynamic file is not periodic with periodic interpolation'
215               WRITE (numout,*) ' we take the last value for the last period '
216               WRITE (numout,*) ' iperm1 = 12,   iper = 13  '
217            ENDIF
218            iperm1 = 12
219            iper   = 13
220         ENDIF
221         !
222         CALL swap_dyn_data         ! We have to prepare a new read of data : swap from record 2 to 1
223         !
224         iswap = 1                  !  indicates swap
225         !
226         CALL dynrea( kt, iper )    ! data read for the iper period
227         !
228         IF( lk_ldfslp .AND. .NOT. lk_c1d ) THEN
229            ! Computes slopes. Caution : here tsn and avt are used as workspace
230            tsn(:,:,:,jp_tem) = tdta  (:,:,:,2)
231            tsn(:,:,:,jp_sal) = sdta  (:,:,:,2)
232            avt(:,:,:)        = avtdta(:,:,:,2)
233            !
234                           CALL eos( tsn, rhd, rhop )                   ! now in situ density
235                           CALL bn2( tsn, rn2 )                         ! now Brunt-Vaisala frequency
236            IF( ln_zps )   CALL zps_hde( kt, jpts, tsn, gtsu, gtsv,  &  ! Partial steps: before Horizontal DErivative
237               &                                   rhd, gru , grv   )   ! of t, s, rd at the bottom ocean level
238            CALL zdf_mxl( kt )                                          ! mixed layer depth
239            CALL ldf_slp( kt, rhd, rn2 )                                ! slope of iso-neutral surfaces
240            !
241            uslpdta (:,:,:,2) = uslp (:,:,:)
242            vslpdta (:,:,:,2) = vslp (:,:,:)
243            wslpidta(:,:,:,2) = wslpi(:,:,:)
244            wslpjdta(:,:,:,2) = wslpj(:,:,:)
245         END IF
246         !
247         ndyn1 = ndyn2         ! store the information of the period read
248         ndyn2 = iper
249         !
250         IF(lwp) THEN
251            WRITE (numout,*) ' dynamics data read for the period ndyn1 =', ndyn1,   &
252               &             ' and for the period ndyn2 = ', ndyn2
253            WRITE (numout,*) ' time step is : ', kt
254         END IF
255         !
256      END IF
257      !
258      ! Compute the data at the given time step
259      !----------------------------------------     
260
261      IF( nsptint == 0 ) THEN          ! No space interpolation, data are probably correct
262         !                             ! We have to initialize data if we have changed the period         
263         CALL assign_dyn_data
264      ELSEIF( nsptint == 1 ) THEN      ! linear interpolation
265         CALL linear_interp_dyn_data( zweigh )
266      ELSE                             ! other interpolation
267         WRITE (numout,*) ' this kind of interpolation do not exist at the moment : we stop'
268         STOP 'dtadyn'         
269      END IF
270      !
271      CALL eos( tsn, rhd, rhop )       ! In any case, we need rhop
272      !
273#if ! defined key_degrad && defined key_traldf_c2d
274      !                                ! In case of 2D varying coefficients, we need aeiv and aeiu
275      IF( lk_traldf_eiv )   CALL dta_eiv( kt )      ! eddy induced velocity coefficient
276#endif
277      !
278      IF( .NOT. l_offbbl ) THEN       ! Compute bbl coefficients if needed
279         tsb(:,:,:,:) = tsn(:,:,:,:)
280         CALL bbl( kt, 'TRC')
281      END IF
282      IF(ln_ctl) THEN
283         CALL prt_ctl(tab3d_1=tsn(:,:,:,jp_tem), clinfo1=' tn      - : ', mask1=tmask, ovlap=1, kdim=jpk   )
284         CALL prt_ctl(tab3d_1=tsn(:,:,:,jp_sal), clinfo1=' sn      - : ', mask1=tmask, ovlap=1, kdim=jpk   )
285         CALL prt_ctl(tab3d_1=un               , clinfo1=' un      - : ', mask1=tmask, ovlap=1, kdim=jpk   )
286         CALL prt_ctl(tab3d_1=vn               , clinfo1=' vn      - : ', mask1=tmask, ovlap=1, kdim=jpk   )
287         CALL prt_ctl(tab3d_1=wn               , clinfo1=' wn      - : ', mask1=tmask, ovlap=1, kdim=jpk   )
288         CALL prt_ctl(tab3d_1=avt              , clinfo1=' kz      - : ', mask1=tmask, ovlap=1, kdim=jpk   )
289         CALL prt_ctl(tab2d_1=fr_i             , clinfo1=' fr_i    - : ', mask1=tmask, ovlap=1 )
290         CALL prt_ctl(tab2d_1=hmld             , clinfo1=' hmld    - : ', mask1=tmask, ovlap=1 )
291         CALL prt_ctl(tab2d_1=emps             , clinfo1=' emps    - : ', mask1=tmask, ovlap=1 )
292         CALL prt_ctl(tab2d_1=wndm             , clinfo1=' wspd    - : ', mask1=tmask, ovlap=1 )
293         CALL prt_ctl(tab2d_1=qsr              , clinfo1=' qsr     - : ', mask1=tmask, ovlap=1 )
294      ENDIF
295      !
296   END SUBROUTINE dta_dyn
297
298   INTEGER FUNCTION dta_dyn_alloc()
299      !!---------------------------------------------------------------------
300      !!                 ***  ROUTINE dta_dyn_alloc  ***
301      !!---------------------------------------------------------------------
302
303      ALLOCATE( tdta    (jpi,jpj,jpk,2), sdta    (jpi,jpj,jpk,2),    &
304       &        udta    (jpi,jpj,jpk,2), vdta    (jpi,jpj,jpk,2),    &
305       &        wdta    (jpi,jpj,jpk,2), avtdta  (jpi,jpj,jpk,2),    &
306#if defined key_ldfslp 
307       &        uslpdta (jpi,jpj,jpk,2), vslpdta (jpi,jpj,jpk,2),    &
308       &        wslpidta(jpi,jpj,jpk,2), wslpjdta(jpi,jpj,jpk,2),    &
309#endif
310#if defined key_degrad
311       &        ahtudta (jpi,jpj,jpk,2), ahtvdta (jpi,jpj,jpk,2),    &
312       &        ahtwdta (jpi,jpj,jpk,2),                             &
313# if defined key_traldf_eiv
314       &        aeiudta (jpi,jpj,jpk,2), aeivdta (jpi,jpj,jpk,2),    &
315       &        aeiwdta (jpi,jpj,jpk,2),                             &
316# endif
317#endif
318#if ! defined key_degrad &&  defined key_traldf_c2d && defined key_traldf_eiv
319       &        aeiwdta (jpi,jpj,    2),                             &
320#endif
321
322       &        hmlddta (jpi,jpj,    2), wspddta (jpi,jpj,    2),    &
323       &        frlddta (jpi,jpj,    2), qsrdta  (jpi,jpj,    2),    &
324       &        empdta  (jpi,jpj,    2),                         STAT=dta_dyn_alloc ) 
325
326      IF( dta_dyn_alloc /= 0 ) CALL ctl_warn('dta_dyn_alloc: failed to allocate facvol array.')
327
328   END FUNCTION dta_dyn_alloc
329
330   SUBROUTINE dynrea( kt, kenr )
331      !!----------------------------------------------------------------------
332      !!                  ***  ROUTINE dynrea  ***
333      !!
334      !! ** Purpose : READ dynamics fiels from OPA9 netcdf output
335      !!
336      !! ** Method : READ the kenr records of DATA and store in udta(...,2), .... 
337      !!----------------------------------------------------------------------
338      USE wrk_nemo, ONLY: wrk_in_use, wrk_not_released
339      USE wrk_nemo, ONLY: zu    => wrk_3d_1 , zv    => wrk_3d_2 , zw    => wrk_3d_3
340      USE wrk_nemo, ONLY: zt    => wrk_3d_4 , zs    => wrk_3d_5
341      USE wrk_nemo, ONLY: zavt  => wrk_3d_6 , zhdiv => wrk_3d_7
342      USE wrk_nemo, ONLY: zahtu => wrk_3d_8 , zahtv => wrk_3d_9 , zahtw => wrk_3d_10
343      USE wrk_nemo, ONLY: zaeiu => wrk_3d_11, zaeiv => wrk_3d_12, zaeiw => wrk_3d_13
344      !
345      USE wrk_nemo, ONLY: zemp  => wrk_2d_1 , zqsr  => wrk_2d_2 , zmld  => wrk_2d_3
346      USE wrk_nemo, ONLY: zice  => wrk_2d_4 , zwspd => wrk_2d_5 
347      USE wrk_nemo, ONLY: ztaux => wrk_2d_6 , ztauy => wrk_2d_7
348      USE wrk_nemo, ONLY: zbblx => wrk_2d_8 , zbbly => wrk_2d_9
349      USE wrk_nemo, ONLY: zaeiw2d => wrk_2d_10
350      !
351      INTEGER, INTENT(in) ::   kt, kenr   ! time index
352      !!
353      INTEGER ::  jkenr
354      !!----------------------------------------------------------------------
355
356      ! 0. Memory allocation
357      IF( wrk_in_use(3, 1,2,3,4,5,6,7,8,9,10,11,12,13) .OR. &
358          wrk_in_use(2, 1,2,3,4,5,6,7,8,9,10)               ) THEN
359         CALL ctl_stop('domrea/dta_dyn: requested workspace arrays unavailable')  ;  RETURN
360      END IF
361     
362      ! cas d'un fichier non periodique : on utilise deux fois le premier et
363      ! le dernier champ temporel
364
365      jkenr = kenr
366
367      IF(lwp) THEN
368         WRITE(numout,*)
369         WRITE(numout,*) 'Dynrea : read dynamical fields, kenr = ', jkenr
370         WRITE(numout,*) '~~~~~~~'
371#if defined key_degrad
372         WRITE(numout,*) ' Degraded fields'
373#endif
374         WRITE(numout,*)
375      ENDIF
376
377
378      IF( kt == nit000 .AND. nlecoff == 0 ) THEN
379         nlecoff = 1
380         CALL  iom_open ( cfile_grid_T, numfl_t )
381         CALL  iom_open ( cfile_grid_U, numfl_u )
382         CALL  iom_open ( cfile_grid_V, numfl_v )
383         CALL  iom_open ( cfile_grid_W, numfl_w )
384      ENDIF
385
386      ! file grid-T
387      !---------------
388      CALL iom_get( numfl_t, jpdom_data, 'votemper', zt   (:,:,:), jkenr )
389      CALL iom_get( numfl_t, jpdom_data, 'vosaline', zs   (:,:,:), jkenr )
390      CALL iom_get( numfl_t, jpdom_data, 'somixhgt', zmld (:,:  ), jkenr )
391      CALL iom_get( numfl_t, jpdom_data, 'sowaflcd', zemp (:,:  ), jkenr )
392      CALL iom_get( numfl_t, jpdom_data, 'soshfldo', zqsr (:,:  ), jkenr )
393      CALL iom_get( numfl_t, jpdom_data, 'soicecov', zice (:,:  ), jkenr )
394      IF( iom_varid( numfl_t, 'sowindsp', ldstop = .FALSE. ) > 0 ) THEN
395         CALL iom_get( numfl_t, jpdom_data, 'sowindsp', zwspd(:,:), jkenr ) 
396      ELSE
397         CALL iom_get( numfl_u, jpdom_data, 'sozotaux', ztaux(:,:), jkenr )
398         CALL iom_get( numfl_v, jpdom_data, 'sometauy', ztauy(:,:), jkenr )
399         CALL tau2wnd( ztaux, ztauy, zwspd )
400      ENDIF
401      ! files grid-U / grid_V
402      CALL iom_get( numfl_u, jpdom_data, 'vozocrtx', zu   (:,:,:), jkenr )
403      CALL iom_get( numfl_v, jpdom_data, 'vomecrty', zv   (:,:,:), jkenr )
404#if defined key_trabbl
405      IF( .NOT. lk_c1d .AND. nn_bbl_ldf == 1 ) THEN
406         IF( iom_varid( numfl_u, 'sobblcox', ldstop = .FALSE. ) > 0  .AND. &
407         &   iom_varid( numfl_v, 'sobblcoy', ldstop = .FALSE. ) > 0 ) THEN
408             CALL iom_get( numfl_u, jpdom_data, 'sobblcox', zbblx(:,:), jkenr )
409             CALL iom_get( numfl_v, jpdom_data, 'sobblcoy', zbbly(:,:), jkenr )
410             l_offbbl = .TRUE.
411         ENDIF
412      ENDIF
413#endif 
414
415      ! file grid-W
416!!      CALL iom_get ( numfl_w, jpdom_data, 'vovecrtz', zw   (:,:,:), jkenr )
417      ! Computation of vertical velocity using horizontal divergence
418      CALL wzv( zu, zv, zw, zhdiv )
419
420      IF( iom_varid( numfl_w, 'voddmavs', ldstop = .FALSE. ) > 0 ) THEN          ! avs exist: it is used
421         CALL iom_get( numfl_w, jpdom_data, 'voddmavs', zavt (:,:,:), jkenr )
422      ELSE                                                                       ! no avs: use avt
423         CALL iom_get( numfl_w, jpdom_data, 'votkeavt', zavt (:,:,:), jkenr )
424      ENDIF
425
426#if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv
427      CALL iom_get( numfl_w, jpdom_data, 'soleaeiw', zaeiw2d(:,: ), jkenr )
428#endif
429
430#if defined key_degrad
431      CALL iom_get( numfl_u, jpdom_data, 'vozoahtu', zahtu(:,:,:), jkenr )
432      CALL iom_get( numfl_v, jpdom_data, 'vomeahtv', zahtv(:,:,:), jkenr )
433      CALL iom_get( numfl_w, jpdom_data, 'voveahtw', zahtw(:,:,:), jkenr )
434#  if defined key_traldf_eiv
435      CALL iom_get( numfl_u, jpdom_data, 'vozoaeiu', zaeiu(:,:,:), jkenr )
436      CALL iom_get( numfl_v, jpdom_data, 'vomeaeiv', zaeiv(:,:,:), jkenr )
437      CALL iom_get( numfl_w, jpdom_data, 'voveaeiw', zaeiw(:,:,:), jkenr )
438#  endif
439#endif
440
441      udta(:,:,:,2) = zu(:,:,:) * umask(:,:,:)
442      vdta(:,:,:,2) = zv(:,:,:) * vmask(:,:,:) 
443      wdta(:,:,:,2) = zw(:,:,:) * tmask(:,:,:)
444
445      tdta(:,:,:,2)   = zt  (:,:,:) * tmask(:,:,:)
446      sdta(:,:,:,2)   = zs  (:,:,:) * tmask(:,:,:)
447      avtdta(:,:,:,2) = zavt(:,:,:) * tmask(:,:,:)
448
449#if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv
450      aeiwdta(:,:,2)  = zaeiw2d(:,:) * tmask(:,:,1)
451#endif
452
453#if defined key_degrad
454        ahtudta(:,:,:,2) = zahtu(:,:,:) * umask(:,:,:)
455        ahtvdta(:,:,:,2) = zahtv(:,:,:) * vmask(:,:,:)
456        ahtwdta(:,:,:,2) = zahtw(:,:,:) * tmask(:,:,:)
457#  if defined key_traldf_eiv
458        aeiudta(:,:,:,2) = zaeiu(:,:,:) * umask(:,:,:)
459        aeivdta(:,:,:,2) = zaeiv(:,:,:) * vmask(:,:,:)
460        aeiwdta(:,:,:,2) = zaeiw(:,:,:) * tmask(:,:,:)
461#  endif
462#endif
463
464      ! fluxes
465      !
466      wspddta(:,:,2)  = zwspd(:,:) * tmask(:,:,1)
467      frlddta(:,:,2)  = MIN( 1., zice(:,:) ) * tmask(:,:,1)
468      empdta (:,:,2)  = zemp(:,:) * tmask(:,:,1)
469      qsrdta (:,:,2)  = zqsr(:,:) * tmask(:,:,1)
470      hmlddta(:,:,2)  = zmld(:,:) * tmask(:,:,1)
471
472#if defined key_trabbl
473      IF( l_offbbl ) THEN
474         bblxdta(:,:,2) = MAX( 0., zbblx(:,:) )
475         bblydta(:,:,2) = MAX( 0., zbbly(:,:) )
476         WHERE( bblxdta(:,:,2) > 2. ) bblxdta(:,:,2) = 0.
477         WHERE( bblydta(:,:,2) > 2. ) bblydta(:,:,2) = 0.
478      ENDIF
479#endif
480     
481      IF( kt == nitend ) THEN
482         CALL iom_close ( numfl_t )
483         CALL iom_close ( numfl_u )
484         CALL iom_close ( numfl_v )
485         CALL iom_close ( numfl_w )
486      ENDIF
487      !     
488      IF( wrk_not_released(3, 1,2,3,4,5,6,7,8,9,10,11,12,13) .OR. &
489          wrk_not_released(2, 1,2,3,4,5,6,7,8,9,10)               ) THEN
490         CALL ctl_stop('domrea/dta_dyn: failed to release workspace arrays.')
491      END IF
492      !
493   END SUBROUTINE dynrea
494
495
496   SUBROUTINE dta_dyn_init
497      !!----------------------------------------------------------------------
498      !!                  ***  ROUTINE dta_dyn_init  ***
499      !!
500      !! ** Purpose :   initializations of parameters for the interpolation
501      !!
502      !! ** Method :
503      !!----------------------------------------------------------------------
504      REAL(wp) :: znspyr   !: number of time step per year
505      INTEGER  :: ierr
506      !!
507      NAMELIST/namdyn/ ndtadyn, ndtatot, nsptint, lperdyn,  &
508      &                cfile_grid_T, cfile_grid_U, cfile_grid_V, cfile_grid_W
509      !!----------------------------------------------------------------------
510
511      ierr = dta_dyn_alloc()
512      IF( ierr /= 0 )  CALL ctl_stop( 'STOP', 'dta_dyn_alloc : unable to allocate standard ocean arrays' )
513
514      !  Define the dynamical input parameters
515      ! ======================================
516
517      REWIND( numnam )              ! Read Namelist namdyn : Lateral physics on tracers
518      READ  ( numnam, namdyn )
519
520      IF(lwp) THEN                  ! control print
521         WRITE(numout,*)
522         WRITE(numout,*) 'namdyn : offline dynamical selection'
523         WRITE(numout,*) '~~~~~~~'
524         WRITE(numout,*) '  Namelist namdyn : set parameters for the lecture of the dynamical fields'
525         WRITE(numout,*) 
526         WRITE(numout,*) ' number of elements in the FILE for a year  ndtadyn = ' , ndtadyn
527         WRITE(numout,*) ' total number of elements in the FILE       ndtatot = ' , ndtatot
528         WRITE(numout,*) ' type of interpolation                      nsptint = ' , nsptint
529         WRITE(numout,*) ' loop on the same FILE                      lperdyn = ' , lperdyn
530         WRITE(numout,*) '  '
531         WRITE(numout,*) ' name of grid_T file                   cfile_grid_T = ', TRIM(cfile_grid_T)   
532         WRITE(numout,*) ' name of grid_U file                   cfile_grid_U = ', TRIM(cfile_grid_U) 
533         WRITE(numout,*) ' name of grid_V file                   cfile_grid_V = ', TRIM(cfile_grid_V) 
534         WRITE(numout,*) ' name of grid_W file                   cfile_grid_W = ', TRIM(cfile_grid_W)     
535         WRITE(numout,*) ' '
536      ENDIF
537      !
538      znspyr   = nyear_len(1) * rday / rdt 
539      rnspdta  = znspyr / FLOAT( ndtadyn )
540      rnspdta2 = rnspdta * 0.5 
541      !
542      CALL dta_dyn( nit000 )
543      !
544   END SUBROUTINE dta_dyn_init
545
546
547   SUBROUTINE wzv( pu, pv, pw, phdiv )
548      !!----------------------------------------------------------------------
549      !!                    ***  ROUTINE wzv  ***
550      !!
551      !! ** Purpose :   Compute the now vertical velocity after the array swap
552      !!
553      !! ** Method  : - compute the now divergence given by :
554      !!         * z-coordinate ONLY !!!!
555      !!         hdiv = 1/(e1t*e2t) [ di(e2u  u) + dj(e1v  v) ]
556      !!     - Using the incompressibility hypothesis, the vertical
557      !!      velocity is computed by integrating the horizontal divergence
558      !!      from the bottom to the surface.
559      !!        The boundary conditions are w=0 at the bottom (no flux).
560      !!----------------------------------------------------------------------
561      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(in   ) :: pu, pv    !:  horizontal velocities
562      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(  out) :: pw        !:  verticla velocity
563      REAL(wp), DIMENSION(jpi,jpj,jpk), INTENT(inout) :: phdiv     !:  horizontal divergence
564      !!
565      INTEGER  ::  ji, jj, jk
566      REAL(wp) ::  zu, zu1, zv, zv1, zet
567      !!----------------------------------------------------------------------
568      !
569      ! Computation of vertical velocity using horizontal divergence
570      phdiv(:,:,:) = 0.
571      DO jk = 1, jpkm1
572         DO jj = 2, jpjm1
573            DO ji = fs_2, fs_jpim1   ! vector opt.
574               zu  = pu(ji  ,jj  ,jk) * umask(ji  ,jj  ,jk) * e2u(ji  ,jj  ) * fse3u(ji  ,jj  ,jk)
575               zu1 = pu(ji-1,jj  ,jk) * umask(ji-1,jj  ,jk) * e2u(ji-1,jj  ) * fse3u(ji-1,jj  ,jk)
576               zv  = pv(ji  ,jj  ,jk) * vmask(ji  ,jj  ,jk) * e1v(ji  ,jj  ) * fse3v(ji  ,jj  ,jk)
577               zv1 = pv(ji  ,jj-1,jk) * vmask(ji  ,jj-1,jk) * e1v(ji  ,jj-1) * fse3v(ji  ,jj-1,jk)
578               zet = 1. / ( e1t(ji,jj) * e2t(ji,jj) * fse3t(ji,jj,jk) )
579               phdiv(ji,jj,jk) = ( zu - zu1 + zv - zv1 ) * zet 
580            END DO
581         END DO
582      END DO
583      CALL lbc_lnk( phdiv, 'T', 1. )      ! Lateral boundary conditions on phdiv
584      !
585      ! computation of vertical velocity from the bottom
586      pw(:,:,jpk) = 0._wp
587      DO jk = jpkm1, 1, -1
588         pw(:,:,jk) = pw(:,:,jk+1) - fse3t(:,:,jk) * phdiv(:,:,jk)
589      END DO
590      !
591   END SUBROUTINE wzv
592
593
594   SUBROUTINE dta_eiv( kt )
595      !!----------------------------------------------------------------------
596      !!                  ***  ROUTINE dta_eiv  ***
597      !!
598      !! ** Purpose :   Compute the eddy induced velocity coefficient from the
599      !!      growth rate of baroclinic instability.
600      !!
601      !! ** Method : Specific to the offline model. Computes the horizontal
602      !!             values from the vertical value
603      !!----------------------------------------------------------------------
604      INTEGER, INTENT( in ) ::   kt     ! ocean time-step inedx
605      !!
606      INTEGER ::   ji, jj           ! dummy loop indices
607      !!----------------------------------------------------------------------
608      !
609      IF( kt == nit000 ) THEN
610         IF(lwp) WRITE(numout,*)
611         IF(lwp) WRITE(numout,*) 'dta_eiv : eddy induced velocity coefficients'
612         IF(lwp) WRITE(numout,*) '~~~~~~~'
613      ENDIF
614      !
615#if defined key_ldfeiv
616      ! Average the diffusive coefficient at u- v- points
617      DO jj = 2, jpjm1
618         DO ji = fs_2, fs_jpim1   ! vector opt.
619            aeiu(ji,jj) = .5 * ( aeiw(ji,jj) + aeiw(ji+1,jj  ) )
620            aeiv(ji,jj) = .5 * ( aeiw(ji,jj) + aeiw(ji  ,jj+1) )
621         END DO
622      END DO
623      CALL lbc_lnk( aeiu, 'U', 1. )   ;   CALL lbc_lnk( aeiv, 'V', 1. )    ! lateral boundary condition
624#endif
625      !
626   END SUBROUTINE dta_eiv
627
628
629   SUBROUTINE tau2wnd( ptaux, ptauy, pwspd )
630      !!---------------------------------------------------------------------
631      !!                    ***  ROUTINE sbc_tau2wnd  ***
632      !!
633      !! ** Purpose : Estimation of wind speed as a function of wind stress
634      !!
635      !! ** Method  : |tau|=rhoa*Cd*|U|^2
636      !!---------------------------------------------------------------------
637      REAL(wp), DIMENSION(jpi,jpj), INTENT(in   ) ::   ptaux, ptauy   ! wind stress in i-j direction resp.
638      REAL(wp), DIMENSION(jpi,jpj), INTENT(  out) ::   pwspd          ! wind speed
639      !!
640      REAL(wp) ::   zrhoa  = 1.22_wp       ! Air density kg/m3
641      REAL(wp) ::   zcdrag = 1.5e-3_wp     ! drag coefficient
642      REAL(wp) ::   ztx, zty, ztau, zcoef  ! temporary variables
643      INTEGER  ::   ji, jj                 ! dummy indices
644      !!---------------------------------------------------------------------
645      zcoef = 1. / ( zrhoa * zcdrag )
646!CDIR NOVERRCHK
647      DO jj = 2, jpjm1
648!CDIR NOVERRCHK
649         DO ji = fs_2, fs_jpim1   ! vector opt.
650            ztx = ptaux(ji,jj) * umask(ji,jj,1) + ptaux(ji-1,jj  ) * umask(ji-1,jj  ,1)
651            zty = ptauy(ji,jj) * vmask(ji,jj,1) + ptauy(ji  ,jj-1) * vmask(ji  ,jj-1,1)
652            ztau = 0.5 * SQRT( ztx * ztx + zty * zty )
653            pwspd(ji,jj) = SQRT ( ztau * zcoef ) * tmask(ji,jj,1)
654         END DO
655      END DO
656      CALL lbc_lnk( pwspd(:,:), 'T', 1. )
657      !
658   END SUBROUTINE tau2wnd
659
660
661   SUBROUTINE swap_dyn_data
662      !!----------------------------------------------------------------------
663      !!                    ***  ROUTINE swap_dyn_data  ***
664      !!
665      !! ** Purpose :   swap array data
666      !!----------------------------------------------------------------------
667      !
668      ! swap from record 2 to 1
669      tdta   (:,:,:,1) = tdta   (:,:,:,2)
670      sdta   (:,:,:,1) = sdta   (:,:,:,2)
671      avtdta (:,:,:,1) = avtdta (:,:,:,2)
672      udta   (:,:,:,1) = udta   (:,:,:,2)
673      vdta   (:,:,:,1) = vdta   (:,:,:,2)
674      wdta   (:,:,:,1) = wdta   (:,:,:,2)
675#if defined key_ldfslp && ! defined key_c1d
676      uslpdta (:,:,:,1) = uslpdta (:,:,:,2)
677      vslpdta (:,:,:,1) = vslpdta (:,:,:,2)
678      wslpidta(:,:,:,1) = wslpidta(:,:,:,2)
679      wslpjdta(:,:,:,1) = wslpjdta(:,:,:,2)
680#endif
681      hmlddta(:,:,1) = hmlddta(:,:,2) 
682      wspddta(:,:,1) = wspddta(:,:,2) 
683      frlddta(:,:,1) = frlddta(:,:,2) 
684      empdta (:,:,1) = empdta (:,:,2) 
685      qsrdta (:,:,1) = qsrdta (:,:,2) 
686      IF( l_offbbl ) THEN
687         bblxdta(:,:,1) = bblxdta(:,:,2)
688         bblydta(:,:,1) = bblydta(:,:,2) 
689      ENDIF
690
691#if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv
692      aeiwdta(:,:,1) = aeiwdta(:,:,2)
693#endif
694
695#if defined key_degrad
696      ahtudta(:,:,:,1) = ahtudta(:,:,:,2)
697      ahtvdta(:,:,:,1) = ahtvdta(:,:,:,2)
698      ahtwdta(:,:,:,1) = ahtwdta(:,:,:,2)
699#  if defined key_traldf_eiv
700      aeiudta(:,:,:,1) = aeiudta(:,:,:,2)
701      aeivdta(:,:,:,1) = aeivdta(:,:,:,2)
702      aeiwdta(:,:,:,1) = aeiwdta(:,:,:,2)
703#  endif
704#endif
705      !
706   END SUBROUTINE swap_dyn_data
707
708
709   SUBROUTINE assign_dyn_data
710      !!----------------------------------------------------------------------
711      !!                    ***  ROUTINE assign_dyn_data  ***
712      !!
713      !! ** Purpose :   Assign dynamical data to the data that have been read
714      !!                without time interpolation
715      !!
716      !!----------------------------------------------------------------------
717     
718      tsn(:,:,:,jp_tem) = tdta  (:,:,:,2)
719      tsn(:,:,:,jp_sal) = sdta  (:,:,:,2)
720      avt(:,:,:)        = avtdta(:,:,:,2)
721     
722      un (:,:,:) = udta  (:,:,:,2) 
723      vn (:,:,:) = vdta  (:,:,:,2)
724      wn (:,:,:) = wdta  (:,:,:,2)
725     
726#if defined key_ldfslp && ! defined key_c1d
727      uslp (:,:,:) = uslpdta (:,:,:,2) 
728      vslp (:,:,:) = vslpdta (:,:,:,2) 
729      wslpi(:,:,:) = wslpidta(:,:,:,2) 
730      wslpj(:,:,:) = wslpjdta(:,:,:,2) 
731#endif
732
733      hmld(:,:) = hmlddta(:,:,2) 
734      wndm(:,:) = wspddta(:,:,2) 
735      fr_i(:,:) = frlddta(:,:,2) 
736      emp (:,:) = empdta (:,:,2) 
737      emps(:,:) = emp(:,:) 
738      qsr (:,:) = qsrdta (:,:,2) 
739#if defined key_trabbl
740      IF( l_offbbl ) THEN
741         ahu_bbl(:,:) = bblxdta(:,:,2)
742         ahv_bbl(:,:) = bblydta(:,:,2) 
743      ENDIF
744#endif
745#if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv
746      aeiw(:,:) = aeiwdta(:,:,2)
747#endif
748     
749#if defined key_degrad
750      ahtu(:,:,:) = ahtudta(:,:,:,2)
751      ahtv(:,:,:) = ahtvdta(:,:,:,2)
752      ahtw(:,:,:) = ahtwdta(:,:,:,2)
753#  if defined key_traldf_eiv
754      aeiu(:,:,:) = aeiudta(:,:,:,2)
755      aeiv(:,:,:) = aeivdta(:,:,:,2)
756      aeiw(:,:,:) = aeiwdta(:,:,:,2)
757#  endif
758#endif
759      !
760   END SUBROUTINE assign_dyn_data
761
762
763   SUBROUTINE linear_interp_dyn_data( pweigh )
764      !!----------------------------------------------------------------------
765      !!               ***  ROUTINE linear_interp_dyn_data  ***
766      !!
767      !! ** Purpose :   linear interpolation of data
768      !!----------------------------------------------------------------------
769      REAL(wp), INTENT(in) ::   pweigh   ! weigh
770      !!
771      REAL(wp) :: zweighm1
772      !!----------------------------------------------------------------------
773
774      zweighm1 = 1. - pweigh
775     
776      tsn(:,:,:,jp_tem) = zweighm1 * tdta  (:,:,:,1) + pweigh * tdta  (:,:,:,2)
777      tsn(:,:,:,jp_sal) = zweighm1 * sdta  (:,:,:,1) + pweigh * sdta  (:,:,:,2)
778      avt(:,:,:)        = zweighm1 * avtdta(:,:,:,1) + pweigh * avtdta(:,:,:,2)
779     
780      un (:,:,:) = zweighm1 * udta  (:,:,:,1) + pweigh * udta  (:,:,:,2) 
781      vn (:,:,:) = zweighm1 * vdta  (:,:,:,1) + pweigh * vdta  (:,:,:,2)
782      wn (:,:,:) = zweighm1 * wdta  (:,:,:,1) + pweigh * wdta  (:,:,:,2)
783     
784#if defined key_ldfslp && ! defined key_c1d
785      uslp (:,:,:) = zweighm1 * uslpdta (:,:,:,1) + pweigh * uslpdta (:,:,:,2) 
786      vslp (:,:,:) = zweighm1 * vslpdta (:,:,:,1) + pweigh * vslpdta (:,:,:,2) 
787      wslpi(:,:,:) = zweighm1 * wslpidta(:,:,:,1) + pweigh * wslpidta(:,:,:,2) 
788      wslpj(:,:,:) = zweighm1 * wslpjdta(:,:,:,1) + pweigh * wslpjdta(:,:,:,2) 
789#endif
790
791      hmld(:,:) = zweighm1 * hmlddta(:,:,1) + pweigh  * hmlddta(:,:,2) 
792      wndm(:,:) = zweighm1 * wspddta(:,:,1) + pweigh  * wspddta(:,:,2) 
793      fr_i(:,:) = zweighm1 * frlddta(:,:,1) + pweigh  * frlddta(:,:,2) 
794      emp (:,:) = zweighm1 * empdta (:,:,1) + pweigh  * empdta (:,:,2) 
795      emps(:,:) = emp(:,:) 
796      qsr (:,:) = zweighm1 * qsrdta (:,:,1) + pweigh  * qsrdta (:,:,2) 
797#if defined key_trabbl
798      IF( l_offbbl ) THEN
799         ahu_bbl(:,:) = zweighm1 * bblxdta(:,:,1) +  pweigh  * bblxdta(:,:,2)
800         ahv_bbl(:,:) = zweighm1 * bblydta(:,:,1) +  pweigh  * bblydta(:,:,2)
801      ENDIF
802#endif
803
804#if ! defined key_degrad && defined key_traldf_c2d && defined key_traldf_eiv
805      aeiw(:,:) = zweighm1 * aeiwdta(:,:,1) + pweigh * aeiwdta(:,:,2)
806#endif
807     
808#if defined key_degrad
809      ahtu(:,:,:) = zweighm1 * ahtudta(:,:,:,1) + pweigh * ahtudta(:,:,:,2)
810      ahtv(:,:,:) = zweighm1 * ahtvdta(:,:,:,1) + pweigh * ahtvdta(:,:,:,2)
811      ahtw(:,:,:) = zweighm1 * ahtwdta(:,:,:,1) + pweigh * ahtwdta(:,:,:,2)
812#  if defined key_traldf_eiv
813      aeiu(:,:,:) = zweighm1 * aeiudta(:,:,:,1) + pweigh * aeiudta(:,:,:,2)
814      aeiv(:,:,:) = zweighm1 * aeivdta(:,:,:,1) + pweigh * aeivdta(:,:,:,2)
815      aeiw(:,:,:) = zweighm1 * aeiwdta(:,:,:,1) + pweigh * aeiwdta(:,:,:,2)
816#  endif
817#endif
818      !     
819   END SUBROUTINE linear_interp_dyn_data
820
821   !!======================================================================
822END MODULE dtadyn
Note: See TracBrowser for help on using the repository browser.