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.
diawri.F90 in NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/tests/BENCH/MY_SRC – NEMO

source: NEMO/branches/2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps/tests/BENCH/MY_SRC/diawri.F90 @ 11758

Last change on this file since 11758 was 11758, checked in by acc, 5 years ago

Branch 2019/dev_r10721_KERNEL-02_Storkey_Coward_IMMERSE_first_steps. Removal of unnecessary cyclic dependancy between oce_trc and step and update to previously untested BENCH MY_SRC files

File size: 21.7 KB
Line 
1MODULE diawri
2   !!======================================================================
3   !!                     ***  MODULE  diawri  ***
4   !! Ocean diagnostics :  write ocean output files
5   !!=====================================================================
6   !! History :  OPA  ! 1991-03  (M.-A. Foujols)  Original code
7   !!            4.0  ! 1991-11  (G. Madec)
8   !!                 ! 1992-06  (M. Imbard)  correction restart file
9   !!                 ! 1992-07  (M. Imbard)  split into diawri and rstwri
10   !!                 ! 1993-03  (M. Imbard)  suppress writibm
11   !!                 ! 1998-01  (C. Levy)  NETCDF format using ioipsl INTERFACE
12   !!                 ! 1999-02  (E. Guilyardi)  name of netCDF files + variables
13   !!            8.2  ! 2000-06  (M. Imbard)  Original code (diabort.F)
14   !!   NEMO     1.0  ! 2002-06  (A.Bozec, E. Durand)  Original code (diainit.F)
15   !!             -   ! 2002-09  (G. Madec)  F90: Free form and module
16   !!             -   ! 2002-12  (G. Madec)  merge of diabort and diainit, F90
17   !!                 ! 2005-11  (V. Garnier) Surface pressure gradient organization
18   !!            3.2  ! 2008-11  (B. Lemaire) creation from old diawri
19   !!            3.7  ! 2014-01  (G. Madec) remove eddy induced velocity from no-IOM output
20   !!                 !                     change name of output variables in dia_wri_state
21   !!----------------------------------------------------------------------
22
23   !!----------------------------------------------------------------------
24   !!   dia_wri       : create the standart output files
25   !!   dia_wri_state : create an output NetCDF file for a single instantaeous ocean state and forcing fields
26   !!----------------------------------------------------------------------
27   USE oce            ! ocean dynamics and tracers
28   USE dom_oce        ! ocean space and time domain
29   USE phycst         ! physical constants
30   USE dianam         ! build name of file (routine)
31   USE diahth         ! thermocline diagnostics
32   USE dynadv   , ONLY: ln_dynadv_vec
33   USE icb_oce        ! Icebergs
34   USE icbdia         ! Iceberg budgets
35   USE ldftra         ! lateral physics: eddy diffusivity coef.
36   USE ldfdyn         ! lateral physics: eddy viscosity   coef.
37   USE sbc_oce        ! Surface boundary condition: ocean fields
38   USE sbc_ice        ! Surface boundary condition: ice fields
39   USE sbcssr         ! restoring term toward SST/SSS climatology
40   USE sbcwave        ! wave parameters
41   USE wet_dry        ! wetting and drying
42   USE zdf_oce        ! ocean vertical physics
43   USE zdfdrg         ! ocean vertical physics: top/bottom friction
44   USE zdfmxl         ! mixed layer
45   !
46   USE lbclnk         ! ocean lateral boundary conditions (or mpp link)
47   USE in_out_manager ! I/O manager
48   USE diatmb         ! Top,middle,bottom output
49   USE dia25h         ! 25h Mean output
50   USE iom            !
51   USE ioipsl         !
52
53#if defined key_si3
54   USE ice 
55   USE icewri 
56#endif
57   USE lib_mpp         ! MPP library
58   USE timing          ! preformance summary
59   USE diu_bulk        ! diurnal warm layer
60   USE diu_coolskin    ! Cool skin
61
62   IMPLICIT NONE
63   PRIVATE
64
65   PUBLIC   dia_wri                 ! routines called by step.F90
66   PUBLIC   dia_wri_state
67   PUBLIC   dia_wri_alloc           ! Called by nemogcm module
68
69   INTEGER ::   nid_T, nz_T, nh_T, ndim_T, ndim_hT   ! grid_T file
70   INTEGER ::          nb_T              , ndim_bT   ! grid_T file
71   INTEGER ::   nid_U, nz_U, nh_U, ndim_U, ndim_hU   ! grid_U file
72   INTEGER ::   nid_V, nz_V, nh_V, ndim_V, ndim_hV   ! grid_V file
73   INTEGER ::   nid_W, nz_W, nh_W                    ! grid_W file
74   INTEGER ::   ndex(1)                              ! ???
75   INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_hT, ndex_hU, ndex_hV
76   INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_T, ndex_U, ndex_V
77   INTEGER, SAVE, ALLOCATABLE, DIMENSION(:) :: ndex_bT
78
79   !! * Substitutions
80#  include "vectopt_loop_substitute.h90"
81   !!----------------------------------------------------------------------
82   !! NEMO/OCE 4.0 , NEMO Consortium (2018)
83   !! $Id: diawri.F90 9598 2018-05-15 22:47:16Z nicolasmartin $
84   !! Software governed by the CeCILL licence     (./LICENSE)
85   !!----------------------------------------------------------------------
86CONTAINS
87
88#if defined key_iomput
89   !!----------------------------------------------------------------------
90   !!   'key_iomput'                                        use IOM library
91   !!----------------------------------------------------------------------
92   INTEGER FUNCTION dia_wri_alloc()
93      !
94      dia_wri_alloc = 0
95      !
96   END FUNCTION dia_wri_alloc
97
98   
99   SUBROUTINE dia_wri( kt, Kmm )
100      !!---------------------------------------------------------------------
101      !!                  ***  ROUTINE dia_wri  ***
102      !!                   
103      !! ** Purpose :   Standard output of opa: dynamics and tracer fields
104      !!      NETCDF format is used by default
105      !!
106      !! ** Method  :  use iom_put
107      !!----------------------------------------------------------------------
108      INTEGER, INTENT( in ) ::   kt      ! ocean time-step index
109      INTEGER, INTENT( in ) ::   Kmm     ! ocean time level index
110      !!
111      INTEGER ::   ji, jj, jk       ! dummy loop indices
112      INTEGER ::   ikbot            ! local integer
113      REAL(wp)::   zztmp , zztmpx   ! local scalar
114      REAL(wp)::   zztmp2, zztmpy   !   -      -
115      REAL(wp), DIMENSION(jpi,jpj)     ::   z2d   ! 2D workspace
116      REAL(wp), DIMENSION(jpi,jpj,jpk) ::   z3d   ! 3D workspace
117      !!----------------------------------------------------------------------
118      !
119      IF( ln_timing )   CALL timing_start('dia_wri')
120      !
121      ! Output the initial state and forcings
122      IF( ninist == 1 ) THEN                       
123         CALL dia_wri_state( Kmm, 'output.init' )
124         ninist = 0
125      ENDIF
126
127      ! Output of initial vertical scale factor
128      CALL iom_put("e3t_0", e3t_0(:,:,:) )
129      CALL iom_put("e3u_0", e3u_0(:,:,:) )
130      CALL iom_put("e3v_0", e3v_0(:,:,:) )
131      !
132      CALL iom_put( "e3t" , e3t(:,:,:,Kmm) )
133      CALL iom_put( "e3u" , e3u(:,:,:,Kmm) )
134      CALL iom_put( "e3v" , e3v(:,:,:,Kmm) )
135      CALL iom_put( "e3w" , e3w(:,:,:,Kmm) )
136      IF( iom_use("e3tdef") )   &
137         CALL iom_put( "e3tdef"  , ( ( e3t(:,:,:,Kmm) - e3t_0(:,:,:) ) / e3t_0(:,:,:) * 100 * tmask(:,:,:) ) ** 2 )
138
139      IF( ll_wd ) THEN
140         CALL iom_put( "ssh" , (ssh(:,:,Kmm)+ssh_ref)*tmask(:,:,1) )   ! sea surface height (brought back to the reference used for wetting and drying)
141      ELSE
142         CALL iom_put( "ssh" , ssh(:,:,Kmm) )              ! sea surface height
143      ENDIF
144
145      IF( iom_use("wetdep") )   &                  ! wet depth
146         CALL iom_put( "wetdep" , ht_0(:,:) + ssh(:,:,Kmm) )
147     
148      CALL iom_put( "toce", ts(:,:,:,jp_tem,Kmm) )    ! 3D temperature
149      CALL iom_put(  "sst", ts(:,:,1,jp_tem,Kmm) )    ! surface temperature
150      IF ( iom_use("sbt") ) THEN
151         DO jj = 1, jpj
152            DO ji = 1, jpi
153               ikbot = mbkt(ji,jj)
154               z2d(ji,jj) = ts(ji,jj,ikbot,jp_tem,Kmm)
155            END DO
156         END DO
157         CALL iom_put( "sbt", z2d )                ! bottom temperature
158      ENDIF
159     
160      CALL iom_put( "soce", ts(:,:,:,jp_sal,Kmm) )    ! 3D salinity
161      CALL iom_put(  "sss", ts(:,:,1,jp_sal,Kmm) )    ! surface salinity
162      IF ( iom_use("sbs") ) THEN
163         DO jj = 1, jpj
164            DO ji = 1, jpi
165               ikbot = mbkt(ji,jj)
166               z2d(ji,jj) = ts(ji,jj,ikbot,jp_sal,Kmm)
167            END DO
168         END DO
169         CALL iom_put( "sbs", z2d )                ! bottom salinity
170      ENDIF
171
172      IF ( iom_use("taubot") ) THEN                ! bottom stress
173         zztmp = rau0 * 0.25
174         z2d(:,:) = 0._wp
175         DO jj = 2, jpjm1
176            DO ji = fs_2, fs_jpim1   ! vector opt.
177               zztmp2 = (  ( rCdU_bot(ji+1,jj)+rCdU_bot(ji  ,jj) ) * uu(ji  ,jj,mbku(ji  ,jj),Kmm)  )**2   &
178                  &   + (  ( rCdU_bot(ji  ,jj)+rCdU_bot(ji-1,jj) ) * uu(ji-1,jj,mbku(ji-1,jj),Kmm)  )**2   &
179                  &   + (  ( rCdU_bot(ji,jj+1)+rCdU_bot(ji,jj  ) ) * vv(ji,jj  ,mbkv(ji,jj  ),Kmm)  )**2   &
180                  &   + (  ( rCdU_bot(ji,jj  )+rCdU_bot(ji,jj-1) ) * vv(ji,jj-1,mbkv(ji,jj-1),Kmm)  )**2
181               z2d(ji,jj) = zztmp * SQRT( zztmp2 ) * tmask(ji,jj,1) 
182               !
183            END DO
184         END DO
185         CALL lbc_lnk( 'diawri', z2d, 'T', 1. )
186         CALL iom_put( "taubot", z2d )           
187      ENDIF
188         
189      CALL iom_put( "uoce", uu(:,:,:,Kmm) )            ! 3D i-current
190      CALL iom_put(  "ssu", uu(:,:,1,Kmm) )            ! surface i-current
191      IF ( iom_use("sbu") ) THEN
192         DO jj = 1, jpj
193            DO ji = 1, jpi
194               ikbot = mbku(ji,jj)
195               z2d(ji,jj) = uu(ji,jj,ikbot,Kmm)
196            END DO
197         END DO
198         CALL iom_put( "sbu", z2d )                ! bottom i-current
199      ENDIF
200     
201      CALL iom_put( "voce", vv(:,:,:,Kmm) )            ! 3D j-current
202      CALL iom_put(  "ssv", vv(:,:,1,Kmm) )            ! surface j-current
203      IF ( iom_use("sbv") ) THEN
204         DO jj = 1, jpj
205            DO ji = 1, jpi
206               ikbot = mbkv(ji,jj)
207               z2d(ji,jj) = vv(ji,jj,ikbot,Kmm)
208            END DO
209         END DO
210         CALL iom_put( "sbv", z2d )                ! bottom j-current
211      ENDIF
212
213      CALL iom_put( "woce", ww )                   ! vertical velocity
214      IF( iom_use('w_masstr') .OR. iom_use('w_masstr2') ) THEN   ! vertical mass transport & its square value
215         ! Caution: in the VVL case, it only correponds to the baroclinic mass transport.
216         z2d(:,:) = rau0 * e1e2t(:,:)
217         DO jk = 1, jpk
218            z3d(:,:,jk) = ww(:,:,jk) * z2d(:,:)
219         END DO
220         CALL iom_put( "w_masstr" , z3d ) 
221         IF( iom_use('w_masstr2') )   CALL iom_put( "w_masstr2", z3d(:,:,:) * z3d(:,:,:) )
222      ENDIF
223
224      CALL iom_put( "avt" , avt )                  ! T vert. eddy diff. coef.
225      CALL iom_put( "avs" , avs )                  ! S vert. eddy diff. coef.
226      CALL iom_put( "avm" , avm )                  ! T vert. eddy visc. coef.
227
228      IF( iom_use('logavt') )   CALL iom_put( "logavt", LOG( MAX( 1.e-20_wp, avt(:,:,:) ) ) )
229      IF( iom_use('logavs') )   CALL iom_put( "logavs", LOG( MAX( 1.e-20_wp, avs(:,:,:) ) ) )
230
231      IF ( iom_use("sstgrad") .OR. iom_use("sstgrad2") ) THEN
232         DO jj = 2, jpjm1                                    ! sst gradient
233            DO ji = fs_2, fs_jpim1   ! vector opt.
234               zztmp  = ts(ji,jj,1,jp_tem,Kmm)
235               zztmpx = ( ts(ji+1,jj,1,jp_tem,Kmm) - zztmp ) * r1_e1u(ji,jj) + ( zztmp - ts(ji-1,jj  ,1,jp_tem,Kmm) ) * r1_e1u(ji-1,jj)
236               zztmpy = ( ts(ji,jj+1,1,jp_tem,Kmm) - zztmp ) * r1_e2v(ji,jj) + ( zztmp - ts(ji  ,jj-1,1,jp_tem,Kmm) ) * r1_e2v(ji,jj-1)
237               z2d(ji,jj) = 0.25 * ( zztmpx * zztmpx + zztmpy * zztmpy )   &
238                  &              * umask(ji,jj,1) * umask(ji-1,jj,1) * vmask(ji,jj,1) * umask(ji,jj-1,1)
239            END DO
240         END DO
241         CALL lbc_lnk( 'diawri', z2d, 'T', 1. )
242         CALL iom_put( "sstgrad2",  z2d )          ! square of module of sst gradient
243         z2d(:,:) = SQRT( z2d(:,:) )
244         CALL iom_put( "sstgrad" ,  z2d )          ! module of sst gradient
245      ENDIF
246         
247      ! heat and salt contents
248      IF( iom_use("heatc") ) THEN
249         z2d(:,:)  = 0._wp 
250         DO jk = 1, jpkm1
251            DO jj = 1, jpj
252               DO ji = 1, jpi
253                  z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_tem,Kmm) * tmask(ji,jj,jk)
254               END DO
255            END DO
256         END DO
257         CALL iom_put( "heatc", rau0_rcp * z2d )   ! vertically integrated heat content (J/m2)
258      ENDIF
259
260      IF( iom_use("saltc") ) THEN
261         z2d(:,:)  = 0._wp 
262         DO jk = 1, jpkm1
263            DO jj = 1, jpj
264               DO ji = 1, jpi
265                  z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_sal,Kmm) * tmask(ji,jj,jk)
266               END DO
267            END DO
268         END DO
269         CALL iom_put( "saltc", rau0 * z2d )          ! vertically integrated salt content (PSU*kg/m2)
270      ENDIF
271      !
272      IF ( iom_use("eken") ) THEN
273         z3d(:,:,jpk) = 0._wp 
274         DO jk = 1, jpkm1
275            DO jj = 2, jpjm1
276               DO ji = fs_2, fs_jpim1   ! vector opt.
277                  zztmp  = 0.25_wp * r1_e1e2t(ji,jj) / e3t(ji,jj,jk,Kmm)
278                  z3d(ji,jj,jk) = zztmp * (  uu(ji-1,jj,jk,Kmm)**2 * e2u(ji-1,jj) * e3u(ji-1,jj,jk,Kmm)   &
279                     &                     + uu(ji  ,jj,jk,Kmm)**2 * e2u(ji  ,jj) * e3u(ji  ,jj,jk,Kmm)   &
280                     &                     + vv(ji,jj-1,jk,Kmm)**2 * e1v(ji,jj-1) * e3v(ji,jj-1,jk,Kmm)   &
281                     &                     + vv(ji,jj  ,jk,Kmm)**2 * e1v(ji,jj  ) * e3v(ji,jj  ,jk,Kmm)   )
282               END DO
283            END DO
284         END DO
285         CALL lbc_lnk( 'diawri', z3d, 'T', 1. )
286         CALL iom_put( "eken", z3d )                 ! kinetic energy
287      ENDIF
288      !
289      CALL iom_put( "hdiv", hdiv )                  ! Horizontal divergence
290      !
291      IF( iom_use("u_masstr") .OR. iom_use("u_masstr_vint") .OR. iom_use("u_heattr") .OR. iom_use("u_salttr") ) THEN
292         z3d(:,:,jpk) = 0.e0
293         z2d(:,:) = 0.e0
294         DO jk = 1, jpkm1
295            z3d(:,:,jk) = rau0 * uu(:,:,jk,Kmm) * e2u(:,:) * e3u(:,:,jk,Kmm) * umask(:,:,jk)
296            z2d(:,:) = z2d(:,:) + z3d(:,:,jk)
297         END DO
298         CALL iom_put( "u_masstr"     , z3d )         ! mass transport in i-direction
299         CALL iom_put( "u_masstr_vint", z2d )         ! mass transport in i-direction vertical sum
300      ENDIF
301     
302      IF( iom_use("u_heattr") ) THEN
303         z2d(:,:) = 0._wp 
304         DO jk = 1, jpkm1
305            DO jj = 2, jpjm1
306               DO ji = fs_2, fs_jpim1   ! vector opt.
307                  z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( ts(ji,jj,jk,jp_tem,Kmm) + ts(ji+1,jj,jk,jp_tem,Kmm) )
308               END DO
309            END DO
310         END DO
311         CALL lbc_lnk( 'diawri', z2d, 'U', -1. )
312         CALL iom_put( "u_heattr", 0.5*rcp * z2d )    ! heat transport in i-direction
313      ENDIF
314
315      IF( iom_use("u_salttr") ) THEN
316         z2d(:,:) = 0.e0 
317         DO jk = 1, jpkm1
318            DO jj = 2, jpjm1
319               DO ji = fs_2, fs_jpim1   ! vector opt.
320                  z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( ts(ji,jj,jk,jp_sal,Kmm) + ts(ji+1,jj,jk,jp_sal,Kmm) )
321               END DO
322            END DO
323         END DO
324         CALL lbc_lnk( 'diawri', z2d, 'U', -1. )
325         CALL iom_put( "u_salttr", 0.5 * z2d )        ! heat transport in i-direction
326      ENDIF
327
328     
329      IF( iom_use("v_masstr") .OR. iom_use("v_heattr") .OR. iom_use("v_salttr") ) THEN
330         z3d(:,:,jpk) = 0.e0
331         DO jk = 1, jpkm1
332            z3d(:,:,jk) = rau0 * vv(:,:,jk,Kmm) * e1v(:,:) * e3v(:,:,jk,Kmm) * vmask(:,:,jk)
333         END DO
334         CALL iom_put( "v_masstr", z3d )              ! mass transport in j-direction
335      ENDIF
336     
337      IF( iom_use("v_heattr") ) THEN
338         z2d(:,:) = 0.e0 
339         DO jk = 1, jpkm1
340            DO jj = 2, jpjm1
341               DO ji = fs_2, fs_jpim1   ! vector opt.
342                  z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( ts(ji,jj,jk,jp_tem,Kmm) + ts(ji,jj+1,jk,jp_tem,Kmm) )
343               END DO
344            END DO
345         END DO
346         CALL lbc_lnk( 'diawri', z2d, 'V', -1. )
347         CALL iom_put( "v_heattr", 0.5*rcp * z2d )    !  heat transport in j-direction
348      ENDIF
349
350      IF( iom_use("v_salttr") ) THEN
351         z2d(:,:) = 0._wp 
352         DO jk = 1, jpkm1
353            DO jj = 2, jpjm1
354               DO ji = fs_2, fs_jpim1   ! vector opt.
355                  z2d(ji,jj) = z2d(ji,jj) + z3d(ji,jj,jk) * ( ts(ji,jj,jk,jp_sal,Kmm) + ts(ji,jj+1,jk,jp_sal,Kmm) )
356               END DO
357            END DO
358         END DO
359         CALL lbc_lnk( 'diawri', z2d, 'V', -1. )
360         CALL iom_put( "v_salttr", 0.5 * z2d )        !  heat transport in j-direction
361      ENDIF
362
363      IF( iom_use("tosmint") ) THEN
364         z2d(:,:) = 0._wp
365         DO jk = 1, jpkm1
366            DO jj = 2, jpjm1
367               DO ji = fs_2, fs_jpim1   ! vector opt.
368                  z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) *  ts(ji,jj,jk,jp_tem,Kmm)
369               END DO
370            END DO
371         END DO
372         CALL lbc_lnk( 'diawri', z2d, 'T', -1. )
373         CALL iom_put( "tosmint", rau0 * z2d )        ! Vertical integral of temperature
374      ENDIF
375      IF( iom_use("somint") ) THEN
376         z2d(:,:)=0._wp
377         DO jk = 1, jpkm1
378            DO jj = 2, jpjm1
379               DO ji = fs_2, fs_jpim1   ! vector opt.
380                  z2d(ji,jj) = z2d(ji,jj) + e3t(ji,jj,jk,Kmm) * ts(ji,jj,jk,jp_sal,Kmm)
381               END DO
382            END DO
383         END DO
384         CALL lbc_lnk( 'diawri', z2d, 'T', -1. )
385         CALL iom_put( "somint", rau0 * z2d )         ! Vertical integral of salinity
386      ENDIF
387
388      CALL iom_put( "bn2", rn2 )                      ! Brunt-Vaisala buoyancy frequency (N^2)
389      !
390
391      IF (ln_diatmb)   CALL dia_tmb( Kmm )            ! tmb values
392         
393      IF (ln_dia25h)   CALL dia_25h( kt, Kmm )        ! 25h averaging
394
395      IF( ln_timing )   CALL timing_stop('dia_wri')
396      !
397   END SUBROUTINE dia_wri
398
399#else
400   !!----------------------------------------------------------------------
401   !!   Default option                                  use IOIPSL  library
402   !!----------------------------------------------------------------------
403
404   INTEGER FUNCTION dia_wri_alloc()
405      !
406      dia_wri_alloc = 0
407      !
408   END FUNCTION dia_wri_alloc
409
410   
411   SUBROUTINE dia_wri( kt, Kmm )
412
413      INTEGER, INTENT( in ) ::   kt   ! ocean time-step index
414      INTEGER, INTENT( in ) ::   Kmm  ! ocean time level index
415      IF( ninist == 1 ) THEN     !==  Output the initial state and forcings  ==!
416         CALL dia_wri_state( Kmm, 'output.init' )
417         ninist = 0
418      ENDIF
419      !
420      ! 0. Initialisation
421      ! -----------------
422
423   END SUBROUTINE dia_wri
424
425#endif
426
427   SUBROUTINE dia_wri_state( Kmm, cdfile_name )
428      !!---------------------------------------------------------------------
429      !!                 ***  ROUTINE dia_wri_state  ***
430      !!       
431      !! ** Purpose :   create a NetCDF file named cdfile_name which contains
432      !!      the instantaneous ocean state and forcing fields.
433      !!        Used to find errors in the initial state or save the last
434      !!      ocean state in case of abnormal end of a simulation
435      !!
436      !! ** Method  :   NetCDF files using ioipsl
437      !!      File 'output.init.nc'  is created if ninist = 1 (namelist)
438      !!      File 'output.abort.nc' is created in case of abnormal job end
439      !!----------------------------------------------------------------------
440      INTEGER           , INTENT( in ) ::   Kmm              ! time level index
441      CHARACTER (len=* ), INTENT( in ) ::   cdfile_name      ! name of the file created
442      !!
443      INTEGER :: inum
444      !!----------------------------------------------------------------------
445      !
446      IF(lwp) WRITE(numout,*)
447      IF(lwp) WRITE(numout,*) 'dia_wri_state : single instantaneous ocean state'
448      IF(lwp) WRITE(numout,*) '~~~~~~~~~~~~~   and forcing fields file created '
449      IF(lwp) WRITE(numout,*) '                and named :', cdfile_name, '...nc'
450
451#if defined key_si3
452     CALL iom_open( TRIM(cdfile_name), inum, ldwrt = .TRUE., kdlev = jpl )
453#else
454     CALL iom_open( TRIM(cdfile_name), inum, ldwrt = .TRUE. )
455#endif
456
457      CALL iom_rstput( 0, 0, inum, 'votemper', ts(:,:,:,jp_tem,Kmm) )    ! now temperature
458      CALL iom_rstput( 0, 0, inum, 'vosaline', ts(:,:,:,jp_sal,Kmm) )    ! now salinity
459      CALL iom_rstput( 0, 0, inum, 'sossheig', ssh(:,:,Kmm)              )    ! sea surface height
460      CALL iom_rstput( 0, 0, inum, 'vozocrtx', uu(:,:,:,Kmm)                )    ! now i-velocity
461      CALL iom_rstput( 0, 0, inum, 'vomecrty', vv(:,:,:,Kmm)                )    ! now j-velocity
462      CALL iom_rstput( 0, 0, inum, 'vovecrtz', ww                )    ! now k-velocity
463      IF( ALLOCATED(ahtu) ) THEN
464         CALL iom_rstput( 0, 0, inum,  'ahtu', ahtu              )    ! aht at u-point
465         CALL iom_rstput( 0, 0, inum,  'ahtv', ahtv              )    ! aht at v-point
466      ENDIF
467      IF( ALLOCATED(ahmt) ) THEN
468         CALL iom_rstput( 0, 0, inum,  'ahmt', ahmt              )    ! ahmt at u-point
469         CALL iom_rstput( 0, 0, inum,  'ahmf', ahmf              )    ! ahmf at v-point
470      ENDIF
471      CALL iom_rstput( 0, 0, inum, 'sowaflup', emp - rnf         )    ! freshwater budget
472      CALL iom_rstput( 0, 0, inum, 'sohefldo', qsr + qns         )    ! total heat flux
473      CALL iom_rstput( 0, 0, inum, 'soshfldo', qsr               )    ! solar heat flux
474      CALL iom_rstput( 0, 0, inum, 'soicecov', fr_i              )    ! ice fraction
475      CALL iom_rstput( 0, 0, inum, 'sozotaux', utau              )    ! i-wind stress
476      CALL iom_rstput( 0, 0, inum, 'sometauy', vtau              )    ! j-wind stress
477      IF(  .NOT.ln_linssh  ) THEN             
478         CALL iom_rstput( 0, 0, inum, 'vovvldep', gdept(:,:,:,Kmm)        )    !  T-cell depth
479         CALL iom_rstput( 0, 0, inum, 'vovvle3t', e3t(:,:,:,Kmm)          )    !  T-cell thickness 
480      END IF
481      IF( ln_wave .AND. ln_sdw ) THEN
482         CALL iom_rstput( 0, 0, inum, 'sdzocrtx', usd            )    ! now StokesDrift i-velocity
483         CALL iom_rstput( 0, 0, inum, 'sdmecrty', vsd            )    ! now StokesDrift j-velocity
484         CALL iom_rstput( 0, 0, inum, 'sdvecrtz', wsd            )    ! now StokesDrift k-velocity
485      ENDIF
486 
487#if defined key_si3
488      IF( nn_ice == 2 ) THEN   ! condition needed in case agrif + ice-model but no-ice in child grid
489         CALL ice_wri_state( inum )
490      ENDIF
491#endif
492      !
493      CALL iom_close( inum )
494      !
495   END SUBROUTINE dia_wri_state
496
497   !!======================================================================
498END MODULE diawri
Note: See TracBrowser for help on using the repository browser.