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

Last change on this file since 8894 was 8894, checked in by cbricaud, 6 years ago

crs improvments

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