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

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

cleaning and bug correction in CRS branch

  • Property svn:keywords set to Id
File size: 18.3 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
96      IF( kt /= nit000 )THEN
97         zfse3t(:,:,:) = e3t_b(:,:,:)
98         zfse3u(:,:,:) = e3u_b(:,:,:)
99         zfse3v(:,:,:) = e3v_b(:,:,:)
100         zfse3w(:,:,:) = e3w_b(:,:,:)
101
102         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)
103         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)
104         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)
105         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)
106
107         DO jk = 1, jpk
108            DO ji = 1, jpi_crs
109               DO jj = 1, jpj_crs
110                  IF( e3t_b_crs(ji,jj,jk) == 0._wp ) e3t_b_crs(ji,jj,jk) = e3t_1d(jk)
111                  IF( e3w_b_crs(ji,jj,jk) == 0._wp ) e3w_b_crs(ji,jj,jk) = e3w_1d(jk)
112                  IF( e3u_b_crs(ji,jj,jk) == 0._wp ) e3u_b_crs(ji,jj,jk) = e3t_1d(jk)
113                  IF( e3v_b_crs(ji,jj,jk) == 0._wp ) e3v_b_crs(ji,jj,jk) = e3t_1d(jk)
114               ENDDO
115            ENDDO
116         ENDDO
117
118         e3t_n_crs(:,:,:) = e3t_a_crs(:,:,:)
119         e3u_n_crs(:,:,:) = e3u_a_crs(:,:,:)
120         e3v_n_crs(:,:,:) = e3v_a_crs(:,:,:)
121         e3w_n_crs(:,:,:) = e3w_a_crs(:,:,:)
122
123      ENDIF
124#endif
125      !---------------------------------------------------------------------------------------------------
126      !variables domaine au temps before : swap
127      !---------------------------------------------------------------------------------------------------
128#if defined key_vvl
129      zfse3t(:,:,:) = e3t_b(:,:,:)
130      zfse3u(:,:,:) = e3u_b(:,:,:)
131      zfse3v(:,:,:) = e3v_b(:,:,:)
132      zfse3w(:,:,:) = e3w_b(:,:,:)
133#else
134      zfse3t(:,:,:) = e3t_0(:,:,:)
135      zfse3u(:,:,:) = e3u_0(:,:,:)
136      zfse3v(:,:,:) = e3v_0(:,:,:)
137      zfse3w(:,:,:) = e3w_0(:,:,:)
138#endif
139
140      IF( kt /= nit000 )THEN
141         emp_b_crs(:,:)        = emp_crs(:,:)
142         rnf_b_crs(:,:)        = rnf_crs(:,:)
143         hdivb_crs(:,:,:)      = hdivn_crs(:,:,:)
144      ELSE
145         emp_b_crs(:,:    ) = 0._wp
146         rnf_b_crs(:,:    ) = 0._wp
147         hdivb_crs(:,:,:  ) = 0._wp
148      ENDIF
149
150      !  Temperature
151      zt(:,:,:) = tsb(:,:,:,jp_tem)  ;      zt_crs(:,:,:) = 0._wp
152      CALL crs_dom_ope( zt, 'VOL', 'T', tmask, zt_crs, p_e12=e1e2t, p_e3=zfse3t, psgn=1.0 )
153      tsb_crs(:,:,:,jp_tem) = zt_crs(:,:,:)
154
155      !  Salinity
156      zs(:,:,:) = tsb(:,:,:,jp_sal)  ;      zs_crs(:,:,:) = 0._wp
157      CALL crs_dom_ope( zs, 'VOL', 'T', tmask, zs_crs, p_e12=e1e2t, p_e3=zfse3t, psgn=1.0 )
158      tsb_crs(:,:,:,jp_sal) = zs_crs(:,:,:)
159
160      !  U-velocity
161      CALL crs_dom_ope( ub, 'SUM', 'U', umask, ub_crs, p_e12=e2u, p_e3=zfse3u, p_surf_crs=e2e3u_msk, psgn=-1.0 )
162
163      !  V-velocity
164      CALL crs_dom_ope( vb, 'SUM', 'V', vmask, vb_crs, p_e12=e1v, p_e3=zfse3v, p_surf_crs=e1e3v_msk, psgn=-1.0 )
165
166      ! n2
167      CALL crs_dom_ope( rn2b, 'VOL', 'W', tmask, rb2_crs, p_e12=e1e2t, p_e3=zfse3t, psgn=1.0 )
168
169      !ssh
170      zfse3t(:,:,:) = 1._wp
171      CALL crs_dom_ope( sshb , 'VOL', 'T', tmask, sshb_crs , p_e12=e1e2t, p_e3=zfse3t         , psgn=1.0 )
172
173      !---------------------------------------------------------------------------------------------------
174      !variables at now time :
175      !---------------------------------------------------------------------------------------------------
176#if defined key_vvl
177      zfse3t(:,:,:) = e3t_n(:,:,:)
178      zfse3u(:,:,:) = e3u_n(:,:,:)
179      zfse3v(:,:,:) = e3v_n(:,:,:)
180      zfse3w(:,:,:) = e3w_n(:,:,:)
181      CALL iom_put("e3t",e3t_n_crs)
182      CALL iom_put("e3u",e3u_n_crs)
183      CALL iom_put("e3v",e3v_n_crs)
184      CALL iom_put("e3w",e3w_n_crs)
185#else
186      zfse3t(:,:,:) = e3t_0(:,:,:)
187      zfse3u(:,:,:) = e3u_0(:,:,:)
188      zfse3v(:,:,:) = e3v_0(:,:,:)
189      zfse3w(:,:,:) = e3w_0(:,:,:)
190#endif
191
192#if defined key_vvl
193
194      ! surfaces
195      CALL crs_dom_sfc( umask, 'U', e2e3u_crs, e2e3u_msk, p_e2=e2u, p_e3=zfse3u )
196      CALL crs_dom_sfc( vmask, 'V', e1e3v_crs, e1e3v_msk, p_e1=e1v, p_e3=zfse3v )
197      !cbr CALL iom_put("e2e3u_crs",e2e3u_crs)
198      !CALL iom_put("e2e3u_msk",e2e3u_msk)
199      !CALL iom_put("e1e3v_crs",e1e3v_crs)
200      !CALL iom_put("e1e3v_msk",e1e3v_msk)
201
202      ! vertical scale factors                                                                                 
203      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)
204      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)
205      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)
206      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)
207
208      DO jk = 1, jpk
209         DO ji = 1, jpi_crs
210            DO jj = 1, jpj_crs
211               IF( e3t_max_n_crs(ji,jj,jk) == 0._wp ) e3t_max_n_crs(ji,jj,jk) = e3t_1d(jk)
212               IF( e3w_max_n_crs(ji,jj,jk) == 0._wp ) e3w_max_n_crs(ji,jj,jk) = e3w_1d(jk)
213               IF( e3u_max_n_crs(ji,jj,jk) == 0._wp ) e3u_max_n_crs(ji,jj,jk) = e3t_1d(jk)
214               IF( e3v_max_n_crs(ji,jj,jk) == 0._wp ) e3v_max_n_crs(ji,jj,jk) = e3t_1d(jk)
215            ENDDO
216         ENDDO
217      ENDDO
218
219      ! depth
220      CALL crs_dom_ope( gdepw_n, 'MAX', 'T', tmask, gdept_n_crs, p_e3=zfse3t, psgn=1.0 )
221      CALL crs_dom_ope( gdepw_n, 'MAX', 'W', tmask, gdepw_n_crs, p_e3=zfse3w, psgn=1.0 )
222      DO jk = 1, jpk
223         DO ji = 1, jpi_crs
224            DO jj = 1, jpj_crs
225               IF( gdept_n_crs(ji,jj,jk) .LE. 0._wp ) gdept_n_crs(ji,jj,jk) = gdept_1d(jk)
226               IF( gdepw_n_crs(ji,jj,jk) .LE. 0._wp ) gdepw_n_crs(ji,jj,jk) = gdepw_1d(jk)
227            ENDDO
228         ENDDO
229      ENDDO
230
231      ! volume and facvol
232      CALL crs_dom_facvol( tmask, 'T', e1t, e2t, zfse3t, ocean_volume_crs_t, facvol_t )
233      !cbr CALL iom_put("cvol_crs_t",ocean_volume_crs_t)
234      !
235      bt_crs(:,:,:) = ocean_volume_crs_t(:,:,:) * facvol_t(:,:,:)*tmask_crs(:,:,:)
236      !
237      r1_bt_crs(:,:,:) = 0._wp
238      WHERE( bt_crs /= 0._wp ) r1_bt_crs(:,:,:) = 1._wp / bt_crs(:,:,:)
239
240      CALL crs_dom_facvol( tmask, 'W', e1t, e2t, zfse3w, ocean_volume_crs_w, facvol_w )
241
242#endif
243
244      !  Temperature
245      zt(:,:,:) = tsn(:,:,:,jp_tem)  ;      zt_crs(:,:,:) = 0._wp
246      CALL crs_dom_ope( zt, 'VOL', 'T', tmask, zt_crs, p_e12=e1e2t, p_e3=zfse3t, psgn=1.0 )
247      tsn_crs(:,:,:,jp_tem) = zt_crs(:,:,:)
248
249      CALL iom_put( "toce", tsn_crs(:,:,:,jp_tem) )    ! temp
250      CALL iom_put( "sst" , tsn_crs(:,:,1,jp_tem) )    ! sst
251
252      !  Salinity
253      zs(:,:,:) = tsn(:,:,:,jp_sal)  ;      zs_crs(:,:,:) = 0._wp
254      CALL crs_dom_ope( zs, 'VOL', 'T', tmask, zs_crs, p_e12=e1e2t, p_e3=zfse3t, psgn=1.0 )
255      tsn_crs(:,:,:,jp_sal) = zs_crs(:,:,:)
256
257      CALL iom_put( "soce" , tsn_crs(:,:,:,jp_sal) )    ! sal
258      CALL iom_put( "sss"  , tsn_crs(:,:,1,jp_sal) )    ! sss
259
260      !  U-velocity
261      CALL crs_dom_ope( un, 'SUM', 'U', umask, un_crs, p_e12=e2u, p_e3=zfse3u, p_surf_crs=e2e3u_msk, psgn=-1.0 )
262      CALL iom_put( "uoce"  , un_crs )   ! i-current
263
264      !  V-velocity
265      CALL crs_dom_ope( vn, 'SUM', 'V', vmask, vn_crs, p_e12=e1v, p_e3=zfse3v, p_surf_crs=e1e3v_msk, psgn=-1.0 )
266      CALL iom_put( "voce"  , vn_crs )   ! i-current
267
268      !n2
269      CALL crs_dom_ope( rn2 , 'VOL', 'W', tmask, rn2_crs, p_e12=e1e2t, p_e3=zfse3t, psgn=1.0 )
270     
271      !  Horizontal divergence ( following OPA_SRC/DYN/divcur.F90 )
272      hdivn_crs(:,:,:)=0._wp
273      DO jk = 1, jpkm1
274         DO jj = 2,jpj_crs
275            DO ji = 2,jpi_crs
276               z2dcrsu =  ( un_crs(ji  ,jj  ,jk) * e2e3u_msk(ji  ,jj  ,jk) ) &
277                 &      - ( un_crs(ji-1,jj  ,jk) * e2e3u_msk(ji-1,jj  ,jk) )
278               z2dcrsv =  ( vn_crs(ji  ,jj  ,jk) * e1e3v_msk(ji  ,jj  ,jk) ) &
279                 &      - ( vn_crs(ji  ,jj-1,jk) * e1e3v_msk(ji  ,jj-1,jk) )
280
281               hdivn_crs(ji,jj,jk) = ( z2dcrsu + z2dcrsv )
282            ENDDO
283         ENDDO
284      ENDDO
285      CALL crs_lbc_lnk( hdivn_crs, 'T', 1.0 )
286      !
287      CALL iom_put( "hdiv", hdivn_crs ) 
288
289
290      !  avt, avs
291      SELECT CASE ( nn_crs_kz )
292         CASE ( 0 )
293            CALL crs_dom_ope( avt, 'VOL', 'W', tmask, avt_crs, p_e12=e1e2t, p_e3=zfse3w, psgn=1.0 )
294         CASE ( 1 )
295            CALL crs_dom_ope( avt, 'MAX', 'W', tmask, avt_crs, p_e12=e1e2t, p_e3=zfse3w, psgn=1.0 )
296         CASE ( 2 )
297            CALL crs_dom_ope( avt, 'MIN', 'W', tmask, avt_crs, p_e12=e1e2t, p_e3=zfse3w, psgn=1.0 )
298         CASE ( 3 )
299            CALL crs_dom_ope( avt, 'LOGVOL', 'W', tmask, avt_crs, p_e12=e1e2t, p_e3=zfse3w, p_mask_crs=tmask_crs, psgn=1.0 )
300         CASE ( 4 )
301            CALL crs_dom_ope( avt, 'MED', 'W', tmask, avt_crs, p_e12=e1e2t, p_e3=zfse3w, psgn=1.0 )
302      END SELECT
303      !
304      CALL iom_put( "avt", avt_crs )   !  Kz
305     
306      !2D fields
307      CALL crs_dom_ope( utau , 'SUM', 'U', umask, utau_crs , p_e12=e2u  , p_surf_crs=e2u_crs  , psgn=1.0 )
308      CALL crs_dom_ope( vtau , 'SUM', 'V', vmask, vtau_crs , p_e12=e1v  , p_surf_crs=e1v_crs  , psgn=1.0 )
309      CALL crs_dom_ope( wndm , 'SUM', 'T', tmask, wndm_crs , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
310      CALL crs_dom_ope( rnf  , 'MAX', 'T', tmask, rnf_crs                                     , psgn=1.0 )
311      CALL crs_dom_ope( h_rnf, 'MAX', 'T', tmask, h_rnf_crs                                   , psgn=1.0 )
312
313      z2d=REAL(nk_rnf,wp)
314      CALL crs_dom_ope( z2d  , 'MAX', 'T', tmask, z2d_crs  , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
315      nk_rnf_crs=INT(z2d_crs)
316
317      CALL crs_dom_ope( qsr  , 'SUM', 'T', tmask, qsr_crs  , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
318
319      CALL crs_dom_ope( emp   ,'SUM', 'T', tmask, emp_crs   , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
320      CALL crs_dom_ope( fmmflx,'SUM', 'T', tmask, fmmflx_crs, p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
321      CALL crs_dom_ope( sfx   ,'SUM', 'T', tmask, sfx_crs   , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
322
323      CALL crs_dom_ope( fr_i  ,'SUM', 'T', tmask, fr_i_crs  , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
324      fr_i_crs=MAX( 0._wp, MIN( fr_i_crs , 1._wp ) )
325
326      z2d=REAL(nmln,wp)
327      CALL crs_dom_ope( z2d , 'MAX', 'T', tmask, z2d_crs , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
328      nmln_crs=INT(z2d_crs) 
329      nmln_crs=MAX(nlb10,nmln_crs)   
330
331      CALL iom_put( "utau"     , utau_crs )   ! i-tau output
332      CALL iom_put( "vtau"     , vtau_crs )   ! j-tau output
333      CALL iom_put( "wspd"     , wndm_crs )   ! wind speed output
334      CALL iom_put( "runoffs"  , rnf_crs  )   ! runoff output
335      CALL iom_put( "qsr"      , qsr_crs  )   ! qsr output
336      CALL iom_put( "empmr"    , emp_crs  )   ! water flux output
337      CALL iom_put( "saltflx"  , sfx_crs  )   ! salt flux output
338      CALL iom_put( "ice_cover", fr_i_crs )   ! ice cover output
339
340      zfse3t(:,:,:) = 1._wp
341      CALL crs_dom_ope( sshn , 'VOL', 'T', tmask, sshn_crs , p_e12=e1e2t, p_e3=zfse3t         , psgn=1.0 ) 
342      CALL iom_put( "ssh"      , sshn_crs )   ! ssh output
343
344
345#if defined key_vvl
346      !---------------------------------------------------------------------------------------------------
347      !variables au temps after
348      !---------------------------------------------------------------------------------------------------
349
350      !ssha
351      zfse3t(:,:,:) = 1._wp
352      zt(:,:,:) = tmask(:,:,:)
353      CALL crs_dom_ope( ssha , 'VOL', 'T', zt, ssha_crs , p_e12=e1e2t,  p_e3=zfse3t , psgn=1.0 )
354
355      !vertical scale factors
356      zfse3t(:,:,:) = e3t_a(:,:,:)
357      zfse3u(:,:,:) = e3u_a(:,:,:)
358      zfse3v(:,:,:) = e3v_a(:,:,:)
359      CALL dom_vvl_interpol( zfse3t(:,:,:), zfse3w(:,:,:), 'W'   )
360
361      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)
362      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)
363      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)
364      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)
365
366      DO jk = 1, jpk
367         DO ji = 1, jpi_crs
368            DO jj = 1, jpj_crs
369               IF( e3t_a_crs(ji,jj,jk) == 0._wp ) e3t_a_crs(ji,jj,jk) = e3t_1d(jk)
370               IF( e3w_a_crs(ji,jj,jk) == 0._wp ) e3w_a_crs(ji,jj,jk) = e3w_1d(jk)
371               IF( e3u_a_crs(ji,jj,jk) == 0._wp ) e3u_a_crs(ji,jj,jk) = e3t_1d(jk)
372               IF( e3v_a_crs(ji,jj,jk) == 0._wp ) e3v_a_crs(ji,jj,jk) = e3t_1d(jk)
373            ENDDO
374        ENDDO
375      ENDDO
376
377#endif
378
379#if defined key_vvl
380
381      z1_2dt = 1._wp / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
382      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1._wp / rdt
383
384      wn_crs(:,:,jpk) = 0._wp
385      DO jk = jpkm1, 1, -1
386         wn_crs(:,:,jk) = wn_crs(:,:,jk+1)*e1e2w_msk(:,:,jk+1) - (  hdivn_crs(:,:,jk)                                   &
387               &                          + z1_2dt * e1e2w_crs(:,:,jk) * ( e3t_a_crs(:,:,jk) - e3t_b_crs(:,:,jk) ) ) * tmask_crs(:,:,jk)
388         WHERE( e1e2w_msk(:,:,jk) .NE. 0._wp )  wn_crs(:,:,jk) =  wn_crs(:,:,jk) /e1e2w_msk(:,:,jk)
389
390
391      ENDDO
392#else
393      IF( ln_crs_wn ) THEN
394         CALL crs_dom_ope( wn, 'SUM', 'W', tmask, wn_crs, p_e12=e1e2t, p_surf_crs=e1e2w_msk, psgn=1.0 )
395      ELSE
396         wn_crs(:,:,jpk) = 0._wp
397         DO jk = jpkm1, 1, -1
398            wn_crs(:,:,jk) = e1e2w_msk(:,:,jk+1)*wn_crs(:,:,jk+1) - hdivn_crs(:,:,jk)
399            WHERE( e1e2w_msk(:,:,jk) .NE. 0._wp )  wn_crs(:,:,jk) =  wn_crs(:,:,jk) /e1e2w_msk(:,:,jk)
400         ENDDO
401       ENDIF
402
403#endif
404      CALL iom_put( "woce", wn_crs  )   ! vertical velocity
405
406      !---------------------------------------------------------------------------------------------------
407      !  free memory
408      CALL wrk_dealloc( jpi, jpj, jpk, zfse3t, zfse3w )
409      CALL wrk_dealloc( jpi, jpj, jpk, zfse3u, zfse3v )
410      CALL wrk_dealloc( jpi, jpj, jpk, zt, zs, ztmp   )
411      CALL wrk_dealloc( jpi, jpj, z2d                 )
412      CALL wrk_dealloc( jpi_crs, jpj_crs, jpk, zt_crs, zs_crs )
413      CALL wrk_dealloc( jpi_crs, jpj_crs, z2d_crs     )
414      !
415      CALL iom_swap( "nemo" )     ! return back on high-resolution grid
416      !
417      IF( nn_timing == 1 )   CALL timing_stop('crs_fld')
418      !
419   END SUBROUTINE crs_fld
420
421   !!======================================================================
422END MODULE crsfld
Note: See TracBrowser for help on using the repository browser.