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.
crsfld.F90 in branches/2015/dev_r5003_MERCATOR6_CRS/NEMOGCM/NEMO/OPA_SRC/CRS – NEMO

source: branches/2015/dev_r5003_MERCATOR6_CRS/NEMOGCM/NEMO/OPA_SRC/CRS/crsfld.F90 @ 7217

Last change on this file since 7217 was 7217, checked in by cbricaud, 7 years ago

CRS branch: code cleaning

  • Property svn:keywords set to Id
File size: 18.4 KB
Line 
1MODULE crsfld
2   !!======================================================================
3   !!                     ***  MODULE  crsdfld  ***
4   !!  Ocean coarsening :  coarse ocean fields
5   !!=====================================================================
6   !!   2012-07  (J. Simeon, C. Calone, G. Madec, C. Ethe)
7   !!----------------------------------------------------------------------
8
9   !!----------------------------------------------------------------------
10   !!   crs_fld       : create the standard output files for coarse grid and prep
11   !!                       other variables needed to be passed to TOP
12   !!----------------------------------------------------------------------
13   USE oce             ! ocean dynamics and tracers
14   USE dom_oce         ! ocean space and time domain
15   USE ldftra_oce      ! ocean active tracers: lateral physics
16   USE sbc_oce         ! Surface boundary condition: ocean fields
17   USE zdf_oce         ! vertical  physics: ocean fields
18   USE zdfddm          ! vertical  physics: double diffusion
19   USe zdfmxl
20   USE lbclnk          ! ocean lateral boundary conditions (or mpp link)
21   USE in_out_manager  ! I/O manager
22   USE timing          ! preformance summary
23   USE wrk_nemo        ! working array
24   USE crs
25   USE crsdom
26   USE domvvl
27   USE domvvl_crs
28   USE crslbclnk
29   USE iom
30   USE zdfmxl_crs
31   USE eosbn2
32   USE zdfevd_crs
33   USE zdftke
34   USE zdftke_crs
35
36   USE ieee_arithmetic
37
38   IMPLICIT NONE
39   PRIVATE
40
41   PUBLIC   crs_fld                 ! routines called by step.F90
42
43
44   !! * Substitutions
45#  include "zdfddm_substitute.h90"
46#  include "domzgr_substitute.h90"
47#  include "vectopt_loop_substitute.h90"
48   !!----------------------------------------------------------------------
49   !! NEMO/OPA 3.3 , NEMO Consortium (2010)
50   !! $Id $
51   !! Software governed by the CeCILL licence     (NEMOGCM/NEMO_CeCILL.txt)
52   !!----------------------------------------------------------------------
53CONTAINS
54
55   SUBROUTINE crs_fld( kt )
56      !!---------------------------------------------------------------------
57      !!                  ***  ROUTINE crs_fld  ***
58      !!                   
59      !! ** Purpose :   Basic output of coarsened dynamics and tracer fields
60      !!      NETCDF format is used by default
61      !!      1. Accumulate in time the dimensionally-weighted fields
62      !!      2. At time of output, rescale [1] by dimension and time
63      !!         to yield the spatial and temporal average.
64      !!  See. diawri_dimg.h90, sbcmod.F90
65      !!
66      !! ** Method  : 
67      !!----------------------------------------------------------------------
68      !!
69     
70      INTEGER, INTENT( in ) ::   kt      ! ocean time-step index
71      !!
72      INTEGER               ::   ji, jj, jk              ! dummy loop indices
73      !!
74      REAL(wp), POINTER, DIMENSION(:,:,:) :: zfse3t, zfse3u, zfse3v, zfse3w ! 3D workspace for e3
75      REAL(wp), POINTER, DIMENSION(:,:,:) :: zt, zs , ztmp
76      REAL(wp), POINTER, DIMENSION(:,:)   :: z2d,z2d_crs
77      REAL(wp), POINTER, DIMENSION(:,:,:) :: zt_crs, zs_crs
78      REAL(wp):: z2dcrsu, z2dcrsv
79      REAL(wp):: z1_2dt
80      REAL(wp):: zmin,zmax
81      INTEGER :: i,j,ijis,ijie,ijjs,ijje
82      INTEGER ::  iji,ijj
83      INTEGER :: jl,jm,jn
84      !!
85      !!----------------------------------------------------------------------
86
87      IF( nn_timing == 1 )   CALL timing_start('crs_fld')
88
89      !  Initialize arrays
90      CALL wrk_alloc( jpi, jpj, jpk, zfse3t, zfse3w )
91      CALL wrk_alloc( jpi, jpj, jpk, zfse3u, zfse3v )
92      CALL wrk_alloc( jpi, jpj, jpk, zt, zs , ztmp        )
93      CALL wrk_alloc( jpi, jpj,      z2d            )
94      !
95      CALL wrk_alloc( jpi_crs, jpj_crs, jpk, zt_crs, zs_crs )
96      CALL wrk_alloc( jpi_crs, jpj_crs, z2d_crs     )
97
98      CALL iom_swap( "nemo_crs" )    ! swap on the coarse grid
99
100      !---------------------------------------------------------------------------------------------------
101      !scale factors: before and now
102      !---------------------------------------------------------------------------------------------------
103#if defined key_vvl
104 
105
106      IF( kt /= nit000 )THEN
107         zfse3t(:,:,:) = e3t_b(:,:,:)
108         zfse3u(:,:,:) = e3u_b(:,:,:)
109         zfse3v(:,:,:) = e3v_b(:,:,:)
110         zfse3w(:,:,:) = e3w_b(:,:,:)
111
112         CALL crs_dom_e3( e1t, e2t, zfse3t, p_sfc_3d_crs=e1e2w_crs, cd_type='T', p_mask=tmask, p_e3_crs=e3t_b_crs, p_e3_max_crs=zs_crs)
113         CALL crs_dom_e3( e1t, e2t, zfse3w, p_sfc_3d_crs=e1e2w_crs, cd_type='W', p_mask=tmask, p_e3_crs=e3w_b_crs, p_e3_max_crs=zs_crs)
114         CALL crs_dom_e3( e1u, e2u, zfse3u, p_sfc_2d_crs=e2u_crs  , cd_type='U', p_mask=umask, p_e3_crs=e3u_b_crs, p_e3_max_crs=zs_crs)
115         CALL crs_dom_e3( e1v, e2v, zfse3v, p_sfc_2d_crs=e1v_crs  , cd_type='V', p_mask=vmask, p_e3_crs=e3v_b_crs, p_e3_max_crs=zs_crs)
116
117         DO jk = 1, jpk
118            DO ji = 1, jpi_crs
119               DO jj = 1, jpj_crs
120                  IF( e3t_b_crs(ji,jj,jk) == 0._wp ) e3t_b_crs(ji,jj,jk) = e3t_1d(jk)
121                  IF( e3w_b_crs(ji,jj,jk) == 0._wp ) e3w_b_crs(ji,jj,jk) = e3w_1d(jk)
122                  IF( e3u_b_crs(ji,jj,jk) == 0._wp ) e3u_b_crs(ji,jj,jk) = e3t_1d(jk)
123                  IF( e3v_b_crs(ji,jj,jk) == 0._wp ) e3v_b_crs(ji,jj,jk) = e3t_1d(jk)
124               ENDDO
125            ENDDO
126         ENDDO
127
128         e3t_n_crs(:,:,:) = e3t_a_crs(:,:,:)
129         e3u_n_crs(:,:,:) = e3u_a_crs(:,:,:)
130         e3v_n_crs(:,:,:) = e3v_a_crs(:,:,:)
131         e3w_n_crs(:,:,:) = e3w_a_crs(:,:,:)
132
133      ENDIF
134#endif
135      !---------------------------------------------------------------------------------------------------
136      !variables domaine au temps before : swap
137      !---------------------------------------------------------------------------------------------------
138#if defined key_vvl
139      zfse3t(:,:,:) = e3t_b(:,:,:)
140      zfse3u(:,:,:) = e3u_b(:,:,:)
141      zfse3v(:,:,:) = e3v_b(:,:,:)
142      zfse3w(:,:,:) = e3w_b(:,:,:)
143#else
144      zfse3t(:,:,:) = e3t_0(:,:,:)
145      zfse3u(:,:,:) = e3u_0(:,:,:)
146      zfse3v(:,:,:) = e3v_0(:,:,:)
147      zfse3w(:,:,:) = e3w_0(:,:,:)
148#endif
149
150      IF( kt /= nit000 )THEN
151         emp_b_crs(:,:)        = emp_crs(:,:)
152         rnf_b_crs(:,:)        = rnf_crs(:,:)
153         hdivb_crs(:,:,:)      = hdivn_crs(:,:,:)
154      ELSE
155         emp_b_crs(:,:    ) = 0._wp
156         rnf_b_crs(:,:    ) = 0._wp
157         hdivb_crs(:,:,:  ) = 0._wp
158      ENDIF
159
160      !  Temperature
161      zt(:,:,:) = tsb(:,:,:,jp_tem)  ;      zt_crs(:,:,:) = 0._wp
162      CALL crs_dom_ope( zt, 'VOL', 'T', tmask, zt_crs, p_e12=e1e2t, p_e3=zfse3t, psgn=1.0 )
163      tsb_crs(:,:,:,jp_tem) = zt_crs(:,:,:)
164
165      !  Salinity
166      zs(:,:,:) = tsb(:,:,:,jp_sal)  ;      zs_crs(:,:,:) = 0._wp
167      CALL crs_dom_ope( zs, 'VOL', 'T', tmask, zs_crs, p_e12=e1e2t, p_e3=zfse3t, psgn=1.0 )
168      tsb_crs(:,:,:,jp_sal) = zs_crs(:,:,:)
169
170      !  U-velocity
171      CALL crs_dom_ope( ub, 'SUM', 'U', umask, ub_crs, p_e12=e2u, p_e3=zfse3u, p_surf_crs=e2e3u_msk, psgn=-1.0 )
172
173      !  V-velocity
174      CALL crs_dom_ope( vb, 'SUM', 'V', vmask, vb_crs, p_e12=e1v, p_e3=zfse3v, p_surf_crs=e1e3v_msk, psgn=-1.0 )
175
176      ! n2
177      CALL crs_dom_ope( rn2b, 'VOL', 'W', tmask, rb2_crs, p_e12=e1e2t, p_e3=zfse3t, psgn=1.0 )
178
179      !ssh
180      zfse3t(:,:,:) = 1._wp
181      CALL crs_dom_ope( sshb , 'VOL', 'T', tmask, sshb_crs , p_e12=e1e2t, p_e3=zfse3t         , psgn=1.0 )
182
183      !---------------------------------------------------------------------------------------------------
184      !variables at now time :
185      !---------------------------------------------------------------------------------------------------
186#if defined key_vvl
187      zfse3t(:,:,:) = e3t_n(:,:,:)
188      zfse3u(:,:,:) = e3u_n(:,:,:)
189      zfse3v(:,:,:) = e3v_n(:,:,:)
190      zfse3w(:,:,:) = e3w_n(:,:,:)
191      CALL iom_put("e3t",e3t_n_crs)
192      CALL iom_put("e3u",e3u_n_crs)
193      CALL iom_put("e3v",e3v_n_crs)
194      CALL iom_put("e3w",e3w_n_crs)
195#else
196      zfse3t(:,:,:) = e3t_0(:,:,:)
197      zfse3u(:,:,:) = e3u_0(:,:,:)
198      zfse3v(:,:,:) = e3v_0(:,:,:)
199      zfse3w(:,:,:) = e3w_0(:,:,:)
200#endif
201
202#if defined key_vvl
203
204      ! surfaces
205      CALL crs_dom_sfc( umask, 'U', e2e3u_crs, e2e3u_msk, p_e2=e2u, p_e3=zfse3u )
206      CALL crs_dom_sfc( vmask, 'V', e1e3v_crs, e1e3v_msk, p_e1=e1v, p_e3=zfse3v )
207      CALL iom_put("e2e3u_crs",e2e3u_crs)
208      CALL iom_put("e2e3u_msk",e2e3u_msk)
209      CALL iom_put("e1e3v_crs",e1e3v_crs)
210      CALL iom_put("e1e3v_msk",e1e3v_msk)
211
212      ! vertical scale factors                                                                                 
213      CALL crs_dom_e3( e1t, e2t, zfse3t, p_sfc_3d_crs=e1e2w_crs, cd_type='T', p_mask=tmask, p_e3_crs=zs_crs, p_e3_max_crs=e3t_max_n_crs)
214      CALL crs_dom_e3( e1t, e2t, zfse3w, p_sfc_3d_crs=e1e2w_crs, cd_type='W', p_mask=tmask, p_e3_crs=zs_crs, p_e3_max_crs=e3w_max_n_crs)
215      CALL crs_dom_e3( e1u, e2u, zfse3u, p_sfc_2d_crs=e2u_crs  , cd_type='U', p_mask=umask, p_e3_crs=zs_crs, p_e3_max_crs=e3u_max_n_crs)
216      CALL crs_dom_e3( e1v, e2v, zfse3v, p_sfc_2d_crs=e1v_crs  , cd_type='V', p_mask=vmask, p_e3_crs=zs_crs, p_e3_max_crs=e3v_max_n_crs)
217
218      !cbr ??????????????????? faut pas mettre le profile 1d plutot ???????????
219      WHERE(e3t_max_n_crs == 0._wp) e3t_max_n_crs=r_inf
220      WHERE(e3u_max_n_crs == 0._wp) e3u_max_n_crs=r_inf
221      WHERE(e3v_max_n_crs == 0._wp) e3v_max_n_crs=r_inf
222      WHERE(e3w_max_n_crs == 0._wp) e3w_max_n_crs=r_inf
223
224      ! depth
225      CALL crs_dom_ope( gdepw_n, 'MAX', 'T', tmask, gdept_n_crs, p_e3=zfse3t, psgn=1.0 )
226      CALL crs_dom_ope( gdepw_n, 'MAX', 'W', tmask, gdepw_n_crs, p_e3=zfse3w, psgn=1.0 )
227
228      ! volume and facvol
229      CALL crs_dom_facvol( tmask, 'T', e1t, e2t, zfse3t, ocean_volume_crs_t, facvol_t )
230      CALL iom_put("cvol_crs_t",ocean_volume_crs_t)
231      !
232      bt_crs(:,:,:) = ocean_volume_crs_t(:,:,:) * facvol_t(:,:,:)*tmask_crs(:,:,:)
233      !
234      r1_bt_crs(:,:,:) = 0._wp
235      WHERE( bt_crs /= 0._wp ) r1_bt_crs(:,:,:) = 1._wp / bt_crs(:,:,:)
236
237      CALL crs_dom_facvol( tmask, 'W', e1t, e2t, zfse3w, ocean_volume_crs_w, facvol_w )
238
239#endif
240
241      !  Temperature
242      zt(:,:,:) = tsn(:,:,:,jp_tem)  ;      zt_crs(:,:,:) = 0._wp
243      CALL crs_dom_ope( zt, 'VOL', 'T', tmask, zt_crs, p_e12=e1e2t, p_e3=zfse3t, psgn=1.0 )
244      tsn_crs(:,:,:,jp_tem) = zt_crs(:,:,:)
245
246      CALL iom_put( "toce", tsn_crs(:,:,:,jp_tem) )    ! temp
247      CALL iom_put( "sst" , tsn_crs(:,:,1,jp_tem) )    ! sst
248
249      !  Salinity
250      zs(:,:,:) = tsn(:,:,:,jp_sal)  ;      zs_crs(:,:,:) = 0._wp
251      CALL crs_dom_ope( zs, 'VOL', 'T', tmask, zs_crs, p_e12=e1e2t, p_e3=zfse3t, psgn=1.0 )
252      tsn_crs(:,:,:,jp_sal) = zs_crs(:,:,:)
253
254      CALL iom_put( "soce" , tsn_crs(:,:,:,jp_sal) )    ! sal
255      CALL iom_put( "sss"  , tsn_crs(:,:,1,jp_sal) )    ! sss
256
257      !  U-velocity
258      CALL crs_dom_ope( un, 'SUM', 'U', umask, un_crs, p_e12=e2u, p_e3=zfse3u, p_surf_crs=e2e3u_msk, psgn=-1.0 )
259      un_crs(:,:,:) = un_crs(:,:,:)*umask_crs(:,:,:) !cbr utile ??????????????????
260      CALL iom_put( "uoce"  , un_crs )   ! i-current
261
262      !  V-velocity
263      CALL crs_dom_ope( vn, 'SUM', 'V', vmask, vn_crs, p_e12=e1v, p_e3=zfse3v, p_surf_crs=e1e3v_msk, psgn=-1.0 )
264      vn_crs(:,:,:) = vn_crs(:,:,:)*vmask_crs(:,:,:) !cbr utile ??????????????????
265      CALL iom_put( "voce"  , vn_crs )   ! i-current
266
267      !n2
268      CALL crs_dom_ope( rn2 , 'VOL', 'W', tmask, rn2_crs, p_e12=e1e2t, p_e3=zfse3t, psgn=1.0 )
269     
270      !  Horizontal divergence ( following OPA_SRC/DYN/divcur.F90 )
271      hdivn_crs(:,:,:)=0._wp
272      DO jk = 1, jpkm1
273         DO jj = 2,jpj_crs
274            DO ji = 2,jpi_crs
275               z2dcrsu =  ( un_crs(ji  ,jj  ,jk) * e2e3u_msk(ji  ,jj  ,jk) ) &
276                 &      - ( un_crs(ji-1,jj  ,jk) * e2e3u_msk(ji-1,jj  ,jk) )
277               z2dcrsv =  ( vn_crs(ji  ,jj  ,jk) * e1e3v_msk(ji  ,jj  ,jk) ) &
278                 &      - ( vn_crs(ji  ,jj-1,jk) * e1e3v_msk(ji  ,jj-1,jk) )
279
280               hdivn_crs(ji,jj,jk) = ( z2dcrsu + z2dcrsv )
281            ENDDO
282         ENDDO
283      ENDDO
284      CALL crs_lbc_lnk( hdivn_crs, 'T', 1.0 )
285      !
286      CALL iom_put( "hdiv", hdivn_crs ) 
287
288
289      !  avt, avs
290      SELECT CASE ( nn_crs_kz )
291         CASE ( 0 )
292            CALL crs_dom_ope( avt, 'VOL', 'W', tmask, avt_crs, p_e12=e1e2t, p_e3=zfse3w, psgn=1.0 )
293         CASE ( 1 )
294            CALL crs_dom_ope( avt, 'MAX', 'W', tmask, avt_crs, p_e12=e1e2t, p_e3=zfse3w, psgn=1.0 )
295         CASE ( 2 )
296            CALL crs_dom_ope( avt, 'MIN', 'W', tmask, avt_crs, p_e12=e1e2t, p_e3=zfse3w, psgn=1.0 )
297         CASE ( 3 )
298            CALL crs_dom_ope( avt, 'LOGVOL', 'W', tmask, avt_crs, p_e12=e1e2t, p_e3=zfse3w, p_mask_crs=tmask_crs, psgn=1.0 )
299         CASE ( 4 )
300            CALL crs_dom_ope( avt, 'MED', 'W', tmask, avt_crs, p_e12=e1e2t, p_e3=zfse3w, psgn=1.0 )
301      END SELECT
302      !
303      CALL iom_put( "avt", avt_crs )   !  Kz
304     
305      !2D fields
306      CALL crs_dom_ope( utau , 'SUM', 'U', umask, utau_crs , p_e12=e2u  , p_surf_crs=e2u_crs  , psgn=1.0 )
307      CALL crs_dom_ope( vtau , 'SUM', 'V', vmask, vtau_crs , p_e12=e1v  , p_surf_crs=e1v_crs  , psgn=1.0 )
308      CALL crs_dom_ope( wndm , 'SUM', 'T', tmask, wndm_crs , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
309      CALL crs_dom_ope( rnf  , 'MAX', 'T', tmask, rnf_crs                                     , psgn=1.0 )
310      CALL crs_dom_ope( qsr  , 'SUM', 'T', tmask, qsr_crs  , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
311#if defined key_vvl
312      CALL crs_dom_ope( gdepw_n, 'MAX', 'W', tmask, gdepw_n_crs, p_e3=zfse3w, psgn=1.0 )
313#else
314      CALL crs_dom_ope( gdepw_0, 'MAX', 'W', tmask, gdepw_0_crs, p_e3=zfse3w, psgn=1.0 )
315#endif
316      CALL crs_dom_ope( emp   ,'SUM', 'T', tmask, emp_crs   , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
317      CALL crs_dom_ope( fmmflx,'SUM', 'T', tmask, fmmflx_crs, p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
318      CALL crs_dom_ope( sfx   ,'SUM', 'T', tmask, sfx_crs   , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
319      CALL crs_dom_ope( fr_i  ,'SUM', 'T', tmask, fr_i_crs  , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
320
321      z2d=REAL(nmln,wp)
322      CALL crs_dom_ope( z2d , 'MAX', 'T', tmask, z2d_crs , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
323      nmln_crs=INT(z2d_crs) 
324      nmln_crs=MAX(nlb10,nmln_crs)   
325
326      CALL iom_put( "utau"     , utau_crs )   ! i-tau output
327      CALL iom_put( "vtau"     , vtau_crs )   ! j-tau output
328      CALL iom_put( "wspd"     , wndm_crs )   ! wind speed output
329      CALL iom_put( "runoffs"  , rnf_crs  )   ! runoff output
330      CALL iom_put( "qsr"      , qsr_crs  )   ! qsr output
331      CALL iom_put( "empmr"    , emp_crs  )   ! water flux output
332      CALL iom_put( "saltflx"  , sfx_crs  )   ! salt flux output
333      CALL iom_put( "ice_cover", fr_i_crs )   ! ice cover output
334
335      zfse3t(:,:,:) = 1._wp
336      CALL crs_dom_ope( sshn , 'VOL', 'T', tmask, sshn_crs , p_e12=e1e2t, p_e3=zfse3t         , psgn=1.0 ) 
337      CALL iom_put( "ssh"      , sshn_crs )   ! ssh output
338
339
340#if defined key_vvl
341      !---------------------------------------------------------------------------------------------------
342      !variables au temps after
343      !---------------------------------------------------------------------------------------------------
344
345      !ssha
346      zfse3t(:,:,:) = 1._wp
347      zt(:,:,:) = tmask(:,:,:)
348      ssha(:,:) = ssha(:,:) * tmask(:,:,1)  !cbr utile ??????????????????
349      CALL crs_dom_ope( ssha , 'VOL', 'T', zt, ssha_crs , p_e12=e1e2t,  p_e3=zfse3t , psgn=1.0 )
350
351      !vertical scale factors
352      zfse3t(:,:,:) = e3t_a(:,:,:)
353      zfse3u(:,:,:) = e3u_a(:,:,:)
354      zfse3v(:,:,:) = e3v_a(:,:,:)
355      CALL dom_vvl_interpol( zfse3t(:,:,:), zfse3w(:,:,:), 'W'   )
356
357      CALL crs_dom_e3( e1t, e2t, zfse3t, p_sfc_3d_crs=e1e2w_crs, cd_type='T', p_mask=tmask, p_e3_crs=e3t_a_crs, p_e3_max_crs=zs_crs)
358      CALL crs_dom_e3( e1t, e2t, zfse3w, p_sfc_3d_crs=e1e2w_crs, cd_type='W', p_mask=tmask, p_e3_crs=e3w_a_crs, p_e3_max_crs=zs_crs)
359      CALL crs_dom_e3( e1u, e2u, zfse3u, p_sfc_2d_crs=e2u_crs  , cd_type='U', p_mask=umask, p_e3_crs=e3u_a_crs, p_e3_max_crs=zs_crs)
360      CALL crs_dom_e3( e1v, e2v, zfse3v, p_sfc_2d_crs=e1v_crs  , cd_type='V', p_mask=vmask, p_e3_crs=e3v_a_crs, p_e3_max_crs=zs_crs)
361
362      DO jk = 1, jpk
363         DO ji = 1, jpi_crs
364            DO jj = 1, jpj_crs
365               IF( e3t_a_crs(ji,jj,jk) == 0._wp ) e3t_a_crs(ji,jj,jk) = e3t_1d(jk)
366               IF( e3w_a_crs(ji,jj,jk) == 0._wp ) e3w_a_crs(ji,jj,jk) = e3w_1d(jk)
367               IF( e3u_a_crs(ji,jj,jk) == 0._wp ) e3u_a_crs(ji,jj,jk) = e3t_1d(jk)
368               IF( e3v_a_crs(ji,jj,jk) == 0._wp ) e3v_a_crs(ji,jj,jk) = e3t_1d(jk)
369            ENDDO
370        ENDDO
371      ENDDO
372
373      !zt_crs=ocean_volume_crs_t ; zs_crs=facvol_t after time !!! ça sert à quoi ???????????????????????????????????????????
374      CALL crs_dom_facvol( tmask, 'T', e1t, e2t, zfse3t, zt_crs, zs_crs )
375
376#endif
377
378#if defined key_vvl
379      z1_2dt = 1._wp / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
380      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1._wp / rdt
381      wn_crs(:,:,jpk) = 0._wp
382      DO jk = jpkm1, 1, -1
383         wn_crs(:,:,jk) = wn_crs(:,:,jk+1)*e1e2w_msk(:,:,jk+1) - (  hdivn_crs(:,:,jk)                                   &
384               &                          + z1_2dt * e1e2w_crs(:,:,jk) * ( e3t_a_crs(:,:,jk) - e3t_b_crs(:,:,jk) ) ) * tmask_crs(:,:,jk)
385         WHERE( e1e2w_msk(:,:,jk) .NE. 0._wp )  wn_crs(:,:,jk) =  wn_crs(:,:,jk) /e1e2w_msk(:,:,jk)
386
387
388      ENDDO
389#else
390      IF( ln_crs_wn ) THEN
391         CALL crs_dom_ope( wn, 'SUM', 'W', tmask, wn_crs, p_e12=e1e2t, p_surf_crs=e1e2w_msk, psgn=1.0 )
392      ELSE
393         wn_crs(:,:,jpk) = 0._wp
394         DO jk = jpkm1, 1, -1
395            wn_crs(:,:,jk) = e1e2w_msk(:,:,jk+1)*wn_crs(:,:,jk+1) - hdivn_crs(:,:,jk)
396            WHERE( e1e2w_msk(:,:,jk) .NE. 0._wp )  wn_crs(:,:,jk) =  wn_crs(:,:,jk) /e1e2w_msk(:,:,jk)
397         ENDDO
398       ENDIF
399
400#endif
401      CALL crs_lbc_lnk( wn_crs, 'W', 1.0 )   !!!!!!!pas utile, nan ??????????????????????
402      wn_crs(:,:,:) = wn_crs(:,:,:) * tmask_crs(:,:,:)   !!!!!!!pas utile, nan ??????????????????????
403
404      CALL iom_put( "woce", wn_crs  )   ! vertical velocity
405
406      !  free memory
407      CALL wrk_dealloc( jpi, jpj, jpk, zfse3t, zfse3w )
408      CALL wrk_dealloc( jpi, jpj, jpk, zfse3u, zfse3v )
409      CALL wrk_dealloc( jpi, jpj, jpk, zt, zs, ztmp   )
410      CALL wrk_dealloc( jpi, jpj, z2d                 )
411      CALL wrk_dealloc( jpi_crs, jpj_crs, jpk, zt_crs, zs_crs )
412      CALL wrk_dealloc( jpi_crs, jpj_crs, z2d_crs     )
413      !
414      CALL iom_swap( "nemo" )     ! return back on high-resolution grid
415      !
416      IF( nn_timing == 1 )   CALL timing_stop('crs_fld')
417      !
418   END SUBROUTINE crs_fld
419
420   !!======================================================================
421END MODULE crsfld
Note: See TracBrowser for help on using the repository browser.