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

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

change computing of rnf_b_crs and emp_b_crs to ensure restartability of passive trcacer in CRS mode

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