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

Last change on this file since 6772 was 6772, checked in by cbricaud, 8 years ago

clean in coarsening branch

  • Property svn:keywords set to Id
File size: 22.2 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, zerr_crs,zmax_crs
78      REAL(wp), POINTER, DIMENSION(:,:,:) :: ztmp_crs
79      REAL(wp), POINTER, DIMENSION(:,:,:,:) :: avte_crs
80      REAL(wp)       :: z2dcrsu, z2dcrsv
81      REAL(wp)       :: z1_2dt
82      REAL(wp)       :: icnt1,icnt2
83      INTEGER :: i,j,ijis,ijie,ijjs,ijje
84      REAL(wp)       :: zw,zwp1,zum1,zu,zvm1,zv,zsnm,zsm,z
85      REAL(wp)       :: zerr, zerr0, zerr1, zmean
86      INTEGER,DIMENSION(4,3) :: ind
87      REAL(wp),DIMENSION(4) :: zwgt
88      INTEGER ::  iji,ijj
89      INTEGER :: jl,jm,jn
90      REAL(wp)       :: zmin,zmax,zsuma0,zsuma1,zsuma2,zsuma3,zsumb0,zsumb1,zsumb2,zsumb3,zsumb4
91      !!
92      !!----------------------------------------------------------------------
93
94      IF( nn_timing == 1 )   CALL timing_start('crs_fld')
95
96      !  Initialize arrays
97      CALL wrk_alloc( jpi, jpj, jpk, zfse3t, zfse3w )
98      CALL wrk_alloc( jpi, jpj, jpk, zfse3u, zfse3v )
99      CALL wrk_alloc( jpi, jpj, jpk, zt, zs , ztmp        )
100      CALL wrk_alloc( jpi, jpj,      z2d            )
101      !
102      CALL wrk_alloc( jpi_crs, jpj_crs, jpk, zt_crs, zs_crs, zerr_crs,zmax_crs )
103      CALL wrk_alloc( jpi_crs, jpj_crs, jpk, ztmp_crs )
104      CALL wrk_alloc( jpi_crs, jpj_crs, jpk, 4, avte_crs )
105      CALL wrk_alloc( jpi_crs, jpj_crs, z2d_crs     )
106
107
108      IF( kt == nit000  ) THEN
109         tsn_crs  (:,:,:,:) = 0._wp    ! temp/sal  array, now
110         un_crs   (:,:,:  ) = 0._wp    ! u-velocity
111         vn_crs   (:,:,:  ) = 0._wp    ! v-velocity
112         wn_crs   (:,:,:  ) = 0._wp    ! w
113         avt_crs  (:,:,:  ) = 0._wp    ! avt
114         hdivb_crs(:,:,:  ) = 0._wp    ! hdiv
115         hdivn_crs(:,:,:  ) = 0._wp    ! hdiv
116         rke_crs  (:,:,:  ) = 0._wp    ! rke
117         sshn_crs (:,:    ) = 0._wp    ! ssh
118         utau_crs (:,:    ) = 0._wp    ! taux
119         vtau_crs (:,:    ) = 0._wp    ! tauy
120         wndm_crs (:,:    ) = 0._wp    ! wind speed
121         qsr_crs  (:,:    ) = 0._wp    ! qsr
122         emp_crs  (:,:    ) = 0._wp    ! emp
123         emp_b_crs(:,:    ) = 0._wp    ! emp
124         rnf_crs  (:,:    ) = 0._wp    ! runoff
125         rnf_b_crs(:,:    ) = 0._wp    ! runoff
126         fr_i_crs (:,:    ) = 0._wp    ! ice cover
127      ENDIF
128
129      CALL iom_swap( "nemo_crs" )    ! swap on the coarse grid
130
131      !---------------------------------------------------------------------------------------------------
132      !variables domaine au temps before : swap
133      !---------------------------------------------------------------------------------------------------
134#if defined key_vvl
135      e3t_b_crs(:,:,:) = e3t_n_crs(:,:,:)
136      e3u_b_crs(:,:,:) = e3u_n_crs(:,:,:)
137      e3v_b_crs(:,:,:) = e3v_n_crs(:,:,:)
138      e3w_b_crs(:,:,:) = e3w_n_crs(:,:,:)
139      e3t_n_crs(:,:,:) = e3t_a_crs(:,:,:)
140      e3u_n_crs(:,:,:) = e3u_a_crs(:,:,:)
141      e3v_n_crs(:,:,:) = e3v_a_crs(:,:,:)
142      e3w_n_crs(:,:,:) = e3w_a_crs(:,:,:)
143#endif
144
145      IF( kt /= nit000 )THEN
146         tsb_crs(:,:,:,jp_tem) = tsn_crs(:,:,:,jp_tem) 
147         tsb_crs(:,:,:,jp_sal) = tsn_crs(:,:,:,jp_sal) 
148         ub_crs(:,:,:)         = un_crs(:,:,:) 
149         vb_crs(:,:,:)         = vn_crs(:,:,:) 
150         sshb_crs(:,:)         = sshb_crs(:,:)
151         emp_b_crs(:,:)        = emp_crs(:,:)
152         rnf_b_crs(:,:)        = rnf_crs(:,:)
153         rb2_crs(:,:,:)        = rn2_crs(:,:,:)
154      ENDIF
155
156      !---------------------------------------------------------------------------------------------------
157      !variables domaine au temps now :
158      !---------------------------------------------------------------------------------------------------
159#if defined key_vvl
160      zfse3t(:,:,:) = e3t_n(:,:,:)
161      zfse3u(:,:,:) = e3u_n(:,:,:)
162      zfse3v(:,:,:) = e3v_n(:,:,:)
163      zfse3w(:,:,:) = e3w_n(:,:,:)
164
165      CALL crs_dom_sfc( umask, 'U', e2e3u_crs, e2e3u_msk, p_e2=e2u, p_e3=zfse3u )
166      CALL crs_dom_sfc( vmask, 'V', e1e3v_crs, e1e3v_msk, p_e1=e1v, p_e3=zfse3v )
167      !                                                                                 
168      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_0_crs)
169      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_0_crs)
170      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_0_crs)
171      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_0_crs)
172
173      CALL crs_dom_ope( gdepw_n, 'MAX', 'T', tmask, gdept_n_crs, p_e3=zfse3t, psgn=1.0 )
174      CALL crs_dom_ope( gdepw_n, 'MAX', 'W', tmask, gdepw_n_crs, p_e3=zfse3w, psgn=1.0 )
175
176      CALL crs_dom_facvol( tmask, 'T', e1t, e2t, zfse3t, ocean_volume_crs_t, facvol_t )
177      CALL iom_put("ocean_volume_crs_t",ocean_volume_crs_t)
178      !
179      bt_crs(:,:,:) = ocean_volume_crs_t(:,:,:) * facvol_t(:,:,:)*tmask_crs(:,:,:)
180      !
181      r1_bt_crs(:,:,:) = 0._wp
182      WHERE( bt_crs /= 0._wp ) r1_bt_crs(:,:,:) = 1._wp / bt_crs(:,:,:)
183
184      CALL crs_dom_facvol( tmask, 'W', e1t, e2t, zfse3w, ocean_volume_crs_w, facvol_w )
185
186#endif
187
188#if defined key_vvl
189      zfse3t(:,:,:) = e3t_n(:,:,:)
190      zfse3u(:,:,:) = e3u_n(:,:,:)
191      zfse3v(:,:,:) = e3v_n(:,:,:)
192      zfse3w(:,:,:) = e3w_n(:,:,:)
193      CALL iom_put("e3t",e3t_n_crs)
194      CALL iom_put("e3u",e3u_n_crs)
195      CALL iom_put("e3v",e3v_n_crs)
196      CALL iom_put("e3w",e3w_n_crs)
197#else
198      zfse3t(:,:,:) = e3t_0(:,:,:)
199      zfse3u(:,:,:) = e3u_0(:,:,:)
200      zfse3v(:,:,:) = e3v_0(:,:,:)
201      zfse3w(:,:,:) = e3w_0(:,:,:)
202#endif
203
204      !  Temperature
205      zt(:,:,:) = tsn(:,:,:,jp_tem)  ;      zt_crs(:,:,:) = 0._wp
206      CALL crs_dom_ope( zt, 'VOL', 'T', tmask, zt_crs, p_e12=e1e2t, p_e3=zfse3t, psgn=1.0 )
207      tsn_crs(:,:,:,jp_tem) = zt_crs(:,:,:)
208
209      CALL iom_put( "toce", tsn_crs(:,:,:,jp_tem) )    ! temp
210      CALL iom_put( "sst" , tsn_crs(:,:,1,jp_tem) )    ! sst
211
212      !  Salinity
213      zs(:,:,:) = tsn(:,:,:,jp_sal)  ;      zs_crs(:,:,:) = 0._wp
214      CALL crs_dom_ope( zs, 'VOL', 'T', tmask, zs_crs, p_e12=e1e2t, p_e3=zfse3t, psgn=1.0 )
215      tsn_crs(:,:,:,jp_sal) = zs_crs(:,:,:)
216
217      CALL iom_put( "soce" , tsn_crs(:,:,:,jp_sal) )    ! sal
218      CALL iom_put( "sss"  , tsn_crs(:,:,1,jp_sal) )    ! sss
219
220      !  U-velocity
221      CALL crs_dom_ope( un, 'SUM', 'U', umask, un_crs, p_e12=e2u, p_e3=zfse3u, p_surf_crs=e2e3u_msk, psgn=-1.0 )
222      un_crs(:,:,:) = un_crs(:,:,:)*umask_crs(:,:,:)
223      CALL iom_put( "uoce"  , un_crs )   ! i-current
224
225      !  V-velocity
226      CALL crs_dom_ope( vn, 'SUM', 'V', vmask, vn_crs, p_e12=e1v, p_e3=zfse3v, p_surf_crs=e1e3v_msk, psgn=-1.0 )
227      vn_crs(:,:,:) = vn_crs(:,:,:)*vmask_crs(:,:,:)
228      CALL iom_put( "voce"  , vn_crs )   ! i-current
229     
230      !  Horizontal divergence ( following OPA_SRC/DYN/divcur.F90 )
231      hdivn_crs(:,:,:)=0._wp
232
233      DO jk = 1, jpkm1
234         DO jj = 2,jpj_crs
235            DO ji = 2,jpi_crs
236               z2dcrsu =  ( un_crs(ji  ,jj  ,jk) * e2e3u_msk(ji  ,jj  ,jk) ) &
237                 &      - ( un_crs(ji-1,jj  ,jk) * e2e3u_msk(ji-1,jj  ,jk) )
238               z2dcrsv =  ( vn_crs(ji  ,jj  ,jk) * e1e3v_msk(ji  ,jj  ,jk) ) &
239                 &      - ( vn_crs(ji  ,jj-1,jk) * e1e3v_msk(ji  ,jj-1,jk) )
240
241               hdivn_crs(ji,jj,jk) = ( z2dcrsu + z2dcrsv )
242            ENDDO
243         ENDDO
244      ENDDO
245      CALL crs_lbc_lnk( hdivn_crs, 'T', 1.0 )
246      !
247      CALL iom_put( "hdiv", hdivn_crs ) 
248
249
250      !  avt, avs
251      SELECT CASE ( nn_crs_kz )
252         CASE ( 0 )
253            CALL crs_dom_ope( avt, 'VOL', 'W', tmask, avt_crs, p_e12=e1e2t, p_e3=zfse3w, psgn=1.0 )
254         CASE ( 1 )
255            CALL crs_dom_ope( avt, 'MAX', 'W', tmask, avt_crs, p_e12=e1e2t, p_e3=zfse3w, psgn=1.0 )
256         CASE ( 2 )
257            CALL crs_dom_ope( avt, 'MIN', 'W', tmask, avt_crs, p_e12=e1e2t, p_e3=zfse3w, psgn=1.0 )
258         CASE ( 3 )
259            CALL crs_dom_ope( avt, 'LOGVOL', 'W', tmask, avt_crs, p_e12=e1e2t, p_e3=zfse3w, p_mask_crs=tmask_crs, psgn=1.0 )
260         CASE ( 4 )
261            CALL crs_dom_ope( avt, 'MED', 'W', tmask, avt_crs, p_e12=e1e2t, p_e3=zfse3w, psgn=1.0 )
262         CASE ( 5 )
263#if defined key_zdftke
264            CALL crs_dom_ope( en , 'VOL', 'W', tmask, en_crs, p_e12=e1e2t, p_e3=zfse3w, psgn=1.0 )
265            CALL crs_dom_ope( taum , 'SUM', 'T', tmask, taum_crs , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
266            CALL crs_dom_ope( rn2(:,:,:), 'VOL', 'W', tmask, rn2_crs, p_e12=e1e2t, p_e3=zfse3t, psgn=1.0 )
267            IF( kt==nit000 )CALL tke_avn_ini_crs
268            CALL tke_avn_crs
269            CALL zdf_evd_crs(kt)
270#endif
271         CASE ( 6 )
272
273            avte_crs(:,:,:,:) = 0._wp
274            ztmp(:,:,:)=1.
275
276            zt(:,:,:) = 0._wp
277            zs(:,:,:) = 0._wp
278            DO jk=2,jpk 
279               WHERE( fse3w(:,:,jk) .NE. 0._wp) zs(:,:,jk) = ( tsn(:,:,jk-1,jp_tem) - tsn(:,:,jk,jp_tem) ) / fse3w(:,:,jk)
280               zt(:,:,jk)=  avt(:,:,jk) *  zs(:,:,jk)
281            ENDDO
282            CALL crs_dom_ope( zt, 'VOL', 'W', tmask, zt_crs, p_e12=e1e2t, p_e3=ztmp, psgn=1.0 )
283            CALL crs_dom_ope( zs, 'VOL', 'W', tmask, zs_crs, p_e12=e1e2t, p_e3=ztmp, psgn=1.0 )
284            zt_crs=tmask_crs*zt_crs
285            WHERE( zs_crs .NE. 0._wp ) avte_crs(:,:,:,1) = zt_crs / zs_crs
286
287            zt(:,:,:) = 0._wp
288            zs(:,:,:) = 0._wp
289            DO jk=2,jpk 
290               WHERE( fse3w(:,:,jk) .NE. 0._wp) zs(:,:,jk) = ( tsn(:,:,jk-1,jp_sal) - tsn(:,:,jk,jp_sal) ) / fse3w(:,:,jk)
291               zt(:,:,jk)=  avt(:,:,jk) *  zs(:,:,jk)
292            ENDDO
293            CALL crs_dom_ope( zt, 'VOL', 'W', tmask, zt_crs, p_e12=e1e2t, p_e3=ztmp, psgn=1.0 )
294            CALL crs_dom_ope( zs, 'VOL', 'W', tmask, zs_crs, p_e12=e1e2t, p_e3=ztmp, psgn=1.0 )
295            zt_crs=tmask_crs*zt_crs
296            WHERE( zs_crs .NE. 0._wp ) avte_crs(:,:,:,2) = zt_crs / zs_crs
297
298            zt(:,:,:) = 0._wp
299            zs(:,:,:) = 0._wp
300            DO jk=2,jpk
301                WHERE( fse3w(:,:,jk) .NE. 0._wp ) zs(:,:,jk)= ( rn_a0 * ( tsn(:,:,jk-1,jp_tem) - tsn(:,:,jk,jp_tem) ) +  &
302                                                                  &   rn_b0 * ( tsn(:,:,jk-1,jp_sal) - tsn(:,:,jk,jp_sal) ) )/ fse3w(:,:,jk) 
303                zt(:,:,jk)=  avt(:,:,jk) *  zs(:,:,jk)                                                           
304            ENDDO
305            CALL crs_dom_ope( zt, 'VOL', 'W', tmask, zt_crs, p_e12=e1e2t, p_e3=ztmp, psgn=1.0 )
306            CALL crs_dom_ope( zs, 'VOL', 'W', tmask, zs_crs, p_e12=e1e2t, p_e3=ztmp, psgn=1.0 )
307            zt_crs=tmask_crs*zt_crs
308            WHERE( zs_crs .NE. 0._wp ) avte_crs(:,:,:,3) = zt_crs / zs_crs
309
310            zt(:,:,:) = 0._wp
311            zs(:,:,:) = 0._wp
312            DO jk=2,jpk
313                WHERE( fse3w(:,:,jk) .NE. 0._wp ) zs(:,:,jk)= ( rn_a0 * ( tsn(:,:,jk-1,jp_tem) - tsn(:,:,jk,jp_tem) ) -  &
314                                                                  &   rn_b0 * ( tsn(:,:,jk-1,jp_sal) - tsn(:,:,jk,jp_sal) ) )/ fse3w(:,:,jk) 
315                zt(:,:,jk)=  avt(:,:,jk) *  zs(:,:,jk)
316            ENDDO
317            CALL crs_dom_ope( zt, 'VOL', 'W', tmask, zt_crs, p_e12=e1e2t, p_e3=ztmp, psgn=1.0 )
318            CALL crs_dom_ope( zs, 'VOL', 'W', tmask, zs_crs, p_e12=e1e2t, p_e3=ztmp, psgn=1.0 )
319            zt_crs=tmask_crs*zt_crs
320            WHERE( zs_crs .NE. 0._wp ) avte_crs(:,:,:,4) = zt_crs / zs_crs
321
322            CALL iom_put( "avte_crs1", avte_crs(:,:,:,1) )   !  Kz
323            CALL iom_put( "avte_crs2", avte_crs(:,:,:,2) )   !  Kz
324            CALL iom_put( "avte_crs3", avte_crs(:,:,:,3) )   !  Kz
325            CALL iom_put( "avte_crs4", avte_crs(:,:,:,4) )   !  Kz
326
327            CALL crs_dom_ope( avt, 'MED', 'W', tmask, zs_crs, p_e12=e1e2t, p_e3=zfse3w, p_mask_crs=tmask_crs, psgn=1.0 )
328            CALL iom_put( "zs_crs", zs_crs )   !  Kzlogvol
329
330            CALL crs_dom_ope( avt, 'VOL', 'W', tmask, zmax_crs, p_e12=e1e2t, p_e3=zfse3w,  psgn=1.0 )
331            CALL iom_put( "zmax_crs", zmax_crs )   !  Kzlogvol
332            avt_crs=zs_crs
333
334
335            zerr0=0.01
336
337            icnt1=0
338            icnt2=0
339
340            zt_crs(:,:,:)=0._wp
341            zerr_crs(:,:,:)=0._wp
342            DO ji=1,jpi_crs 
343            DO jj=1,jpj_crs 
344            DO jk=1,jpk
345
346 
347               zwgt(1:4)=0._wp
348               DO jm=1,4 ; IF( avte_crs(ji,jj,jk,jm)  .GE. 0._wp .AND.  avte_crs(ji,jj,jk,jm)  .LE. zmax_crs(ji,jj,jk)  ) zwgt(jm) = 1._wp ; ENDDO
349               IF( SUM(zwgt(1:4)) .NE. 0._wp ) THEN 
350                  zmean = SUM( zwgt(1:4)*avte_crs(ji,jj,jk,1:4)) / SUM(zwgt(1:4) )
351                  zerr  = SQRT(SUM( zwgt(1:4)*(avte_crs(ji,jj,jk,1:4)-zmean)*(avte_crs(ji,jj,jk,1:4)-zmean) ) / SUM(zwgt(1:4) ) )
352               ELSE
353                  zmean=0._wp
354                 zerr=1.e10
355               ENDIF
356
357               zerr_crs(ji,jj,jk)=zerr
358
359               IF( zerr .LE. zerr0 .AND. zmean .GT. 0._wp )zt_crs(ji,jj,jk)=zmean
360               IF( zerr .LE. zerr0 .AND. zmean .GT. 0._wp )avt_crs(ji,jj,jk)=zmean
361
362               IF( tmask_crs(ji,jj,jk) == 1 ) icnt1=icnt1+1
363               IF( tmask_crs(ji,jj,jk) == 1 .AND.  zerr .LE. zerr0 .AND. zmean .GT. 0._wp ) icnt2=icnt2+1
364
365            ENDDO
366            ENDDO
367            ENDDO
368
369            CALL iom_put( "zt_crs", zt_crs )   !  Kz
370            CALL iom_put( "zerr_crs", zerr_crs )   !  Kz
371
372      END SELECT
373      !
374      CALL iom_put( "avt", avt_crs )   !  Kz
375     
376      CALL crs_dom_ope( sshn , 'VOL', 'T', tmask, sshn_crs , p_e12=e1e2t, p_e3=zfse3t         , psgn=1.0 ) 
377      CALL crs_dom_ope( utau , 'SUM', 'U', umask, utau_crs , p_e12=e2u  , p_surf_crs=e2u_crs  , psgn=1.0 )
378      CALL crs_dom_ope( vtau , 'SUM', 'V', vmask, vtau_crs , p_e12=e1v  , p_surf_crs=e1v_crs  , psgn=1.0 )
379      CALL crs_dom_ope( wndm , 'SUM', 'T', tmask, wndm_crs , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
380      CALL crs_dom_ope( rnf  , 'MAX', 'T', tmask, rnf_crs                                     , psgn=1.0 )
381      CALL crs_dom_ope( qsr  , 'SUM', 'T', tmask, qsr_crs  , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
382#if defined key_vvl
383      CALL crs_dom_ope( gdepw_n, 'MAX', 'W', tmask, gdepw_n_crs, p_e3=zfse3w, psgn=1.0 )
384#else
385      CALL crs_dom_ope( gdepw_0, 'MAX', 'W', tmask, gdepw_0_crs, p_e3=zfse3w, psgn=1.0 )
386#endif
387
388      CALL crs_dom_ope( emp  , 'SUM', 'T', tmask, emp_crs  , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
389
390      CALL crs_dom_ope( fmmflx,'SUM', 'T', tmask, fmmflx_crs, p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
391      CALL crs_dom_ope( sfx  , 'SUM', 'T', tmask, sfx_crs  , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
392      CALL crs_dom_ope( fr_i , 'SUM', 'T', tmask, fr_i_crs , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
393
394      z2d=REAL(nmln,wp)
395      CALL crs_dom_ope( z2d , 'MAX', 'T', tmask, z2d_crs , p_e12=e1e2t, p_surf_crs=e1e2t_crs, psgn=1.0 )
396      nmln_crs=INT(z2d_crs) 
397      nmln_crs=MAX(nlb10,nmln_crs)   
398
399      CALL iom_put( "ssh"      , sshn_crs )   ! ssh output
400      CALL iom_put( "utau"     , utau_crs )   ! i-tau output
401      CALL iom_put( "vtau"     , vtau_crs )   ! j-tau output
402      CALL iom_put( "wspd"     , wndm_crs )   ! wind speed output
403      CALL iom_put( "runoffs"  , rnf_crs  )   ! runoff output
404      CALL iom_put( "qsr"      , qsr_crs  )   ! qsr output
405      CALL iom_put( "empmr"    , emp_crs  )   ! water flux output
406      CALL iom_put( "saltflx"  , sfx_crs  )   ! salt flux output
407      CALL iom_put( "ice_cover", fr_i_crs )   ! ice cover output
408
409#if defined key_vvl
410     !---------------------------------------------------------------------------------------------------
411     !variables au temps after
412     !---------------------------------------------------------------------------------------------------
413
414     zfse3t(:,:,:) = 1._wp
415     zt(:,:,:) = tmask(:,:,:)
416     ssha(:,:) = ssha(:,:) * tmask(:,:,1)
417     CALL crs_dom_ope( ssha , 'VOL', 'T', zt, ssha_crs , p_e12=e1e2t,  p_e3=zfse3t , psgn=1.0 )
418     CALL crs_lbc_lnk( ssha_crs, 'T', 1.0 ) !!!!!!!!!!!!!!!!!!! pas utile !!!!!!!!!!!!!!!!!!!!!!!!!
419
420     zfse3t(:,:,:) = e3t_a(:,:,:)
421     zfse3u(:,:,:) = e3u_a(:,:,:)
422     zfse3v(:,:,:) = e3v_a(:,:,:)
423     CALL dom_vvl_interpol( zfse3t(:,:,:), zfse3w(:,:,:), 'W'   )
424
425     CALL crs_dom_sfc( umask, 'U', zt_crs, zs_crs, p_e2=e2u, p_e3=zfse3u ) ! zt_crs=e2e3u_crs,zs_crs=e2e3u_msk
426     CALL crs_dom_sfc( vmask, 'V', zt_crs, zs_crs, p_e1=e2v, p_e3=zfse3v ) ! zt_crs=e1e3v_crs,zs_crs=e1e3v_msk
427     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)
428     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)
429     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)
430     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)
431
432
433     DO jk = 1, jpk
434        DO ji = 1, jpi_crs
435           DO jj = 1, jpj_crs
436              IF( e3t_a_crs(ji,jj,jk) == 0._wp ) e3t_a_crs(ji,jj,jk) = e3t_1d(jk)
437              IF( e3w_a_crs(ji,jj,jk) == 0._wp ) e3w_a_crs(ji,jj,jk) = e3w_1d(jk)
438              IF( e3u_a_crs(ji,jj,jk) == 0._wp ) e3u_a_crs(ji,jj,jk) = e3t_1d(jk)
439              IF( e3v_a_crs(ji,jj,jk) == 0._wp ) e3v_a_crs(ji,jj,jk) = e3t_1d(jk)
440           ENDDO
441       ENDDO
442     ENDDO
443
444     !zt_crs=ocean_volume_crs_t ; zs_crs=facvol_t after time !!!
445     CALL crs_dom_facvol( tmask, 'T', e1t, e2t, zfse3t, zt_crs, zs_crs )
446
447#endif
448
449#if defined key_vvl
450      z1_2dt = 1._wp / ( 2. * rdt )                         ! set time step size (Euler/Leapfrog)
451      IF( neuler == 0 .AND. kt == nit000 )   z1_2dt = 1._wp / rdt
452      wn_crs(:,:,jpk) = 0._wp
453      DO jk = jpkm1, 1, -1
454         wn_crs(:,:,jk) = wn_crs(:,:,jk+1)*e1e2w_msk(:,:,jk+1) - (  hdivn_crs(:,:,jk)                                   &
455               &                          + z1_2dt * e1e2w_crs(:,:,jk) * ( e3t_a_crs(:,:,jk) - e3t_b_crs(:,:,jk) ) ) * tmask_crs(:,:,jk)
456         WHERE( e1e2w_msk(:,:,jk) .NE. 0._wp )  wn_crs(:,:,jk) =  wn_crs(:,:,jk) /e1e2w_msk(:,:,jk)
457      ENDDO
458#else
459      IF( ln_crs_wn ) THEN
460         CALL crs_dom_ope( wn, 'SUM', 'W', tmask, wn_crs, p_e12=e1e2t, p_surf_crs=e1e2w_msk, psgn=1.0 )
461      ELSE
462         wn_crs(:,:,jpk) = 0._wp
463         DO jk = jpkm1, 1, -1
464            wn_crs(:,:,jk) = e1e2w_msk(:,:,jk+1)*wn_crs(:,:,jk+1) - hdivn_crs(:,:,jk)
465            WHERE( e1e2w_msk(:,:,jk) .NE. 0._wp )  wn_crs(:,:,jk) =  wn_crs(:,:,jk) /e1e2w_msk(:,:,jk)
466         ENDDO
467       ENDIF
468
469#endif
470      CALL crs_lbc_lnk( wn_crs, 'W', 1.0 )   !!!!!!!pas utile, nan ??????????????????????
471      wn_crs(:,:,:) = wn_crs(:,:,:) * tmask_crs(:,:,:)
472      CALL iom_put( "woce", wn_crs  )   ! vertical velocity
473
474      !  free memory
475      CALL wrk_dealloc( jpi, jpj, jpk, zfse3t, zfse3w )
476      CALL wrk_dealloc( jpi, jpj, jpk, zfse3u, zfse3v )
477      CALL wrk_dealloc( jpi, jpj, jpk, zt, zs, ztmp   )
478      CALL wrk_dealloc( jpi, jpj, z2d                 )
479      CALL wrk_dealloc( jpi_crs, jpj_crs, jpk, zt_crs, zs_crs, zerr_crs,zmax_crs )
480      CALL wrk_dealloc( jpi_crs, jpj_crs, jpk, ztmp_crs )
481      CALL wrk_dealloc( jpi_crs, jpj_crs, jpk, 4, avte_crs )
482      CALL wrk_dealloc( jpi_crs, jpj_crs, z2d_crs     )
483      !
484      CALL iom_swap( "nemo" )     ! return back on high-resolution grid
485      !
486      IF( nn_timing == 1 )   CALL timing_stop('crs_fld')
487      !
488   END SUBROUTINE crs_fld
489
490   !!======================================================================
491END MODULE crsfld
Note: See TracBrowser for help on using the repository browser.