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.
histcom.f90 in utils/tools/DOMAINcfg/src – NEMO

source: utils/tools/DOMAINcfg/src/histcom.f90 @ 14623

Last change on this file since 14623 was 14623, checked in by ldebreu, 3 years ago

AGFdomcfg: 1) Update DOMAINcfg to be compliant with the removal of halo cells 2) Update most of the LBC ... subroutines to a recent NEMO 4 version #2638

  • Property svn:executable set to *
File size: 80.8 KB
Line 
1MODULE histcom
2!$AGRIF_DO_NOT_TREAT
3!-
4!$Id: histcom.f90 2368 2010-11-09 15:38:45Z acc $
5!-
6! This software is governed by the CeCILL license
7! See IOIPSL/IOIPSL_License_CeCILL.txt
8!-
9  USE netcdf
10  USE nc4interface ! needed to allow compilation with netcdf3 libraries
11!-
12  USE stringop, ONLY : nocomma,cmpblank,findpos,find_str,strlowercase
13  USE mathelp,  ONLY : mathop,moycum,buildop
14  USE fliocom,  ONLY : flio_dom_file,flio_dom_att
15  USE calendar
16  USE errioipsl, ONLY : ipslerr,ipsldbg
17!-
18  IMPLICIT NONE
19!-
20  PRIVATE
21  PUBLIC :: histbeg,histdef,histhori,histvert,histend, &
22 &          histwrite,histclo,histsync,ioconf_modname
23!---------------------------------------------------------------------
24!- Some confusing vocabulary in this code !
25!- =========================================
26!-
27!- A REGULAR grid is a grid which is i,j indexes
28!- and thus it is stored in a 2D matrix.
29!- This is opposed to a IRREGULAR grid which is only in a vector
30!- and where we do not know which neighbors we have.
31!- As a consequence we need the bounds for each grid-cell.
32!-
33!- A RECTILINEAR grid is a special case of a regular grid
34!- in which all longitudes for i constant are equal
35!- and all latitudes for j constant.
36!- In other words we do not need the full 2D matrix
37!- to describe the grid, just two vectors.
38!---------------------------------------------------------------------
39!-
40  INTERFACE histbeg
41    MODULE PROCEDURE histb_reg1d,histb_reg2d,histb_irreg
42  END INTERFACE
43!-
44  INTERFACE histhori
45    MODULE PROCEDURE histh_reg1d,histh_reg2d,histh_irreg
46  END INTERFACE
47!-
48  INTERFACE histwrite
49!---------------------------------------------------------------------
50!- The "histwrite" routines will give the data to the I/O system.
51!- It will trigger the operations to be performed,
52!- and the writting to the file if needed
53!-
54!- We test for the work to be done at this time here so that at a
55!- later stage we can call different operation and write subroutine
56!- for the REAL and INTEGER interfaces
57!-
58!- INPUT
59!- idf      : The ID of the file on which this variable is to be,
60!-            written. The variable should have been defined in
61!-            this file before.
62!- pvarname : The short name of the variable
63!- pitau    : Current timestep
64!- pdata    : The variable, I mean the real data !
65!- nbindex  : The number of indexes provided. If it is equal to
66!-            the size of the full field as provided in histdef
67!-            then nothing is done.
68!- nindex   : The indices used to expand the variable (pdata)
69!-            onto the full field.
70!---------------------------------------------------------------------
71!- histwrite - we have to prepare different type of fields :
72!-             real and integer, 1,2 or 3D
73    MODULE PROCEDURE histwrite_r1d,histwrite_r2d,histwrite_r3d
74  END INTERFACE
75!-
76! Fixed parameter
77!-
78  INTEGER,PARAMETER :: nb_files_max=20,nb_var_max=400, &
79 &                     nb_hax_max=5,nb_zax_max=10,nbopp_max=10
80  REAL,PARAMETER :: missing_val=nf90_fill_real
81  INTEGER,PARAMETER,PUBLIC :: &
82 &  hist_r4=nf90_real4, hist_r8=nf90_real8
83!-
84! Variable derived type
85!-
86TYPE T_D_V
87  INTEGER :: ncvid
88  INTEGER :: nbopp
89  CHARACTER(LEN=20)  :: v_name,unit_name
90  CHARACTER(LEN=256) :: title,std_name
91  CHARACTER(LEN=80)  :: fullop
92  CHARACTER(LEN=7)   :: topp
93  CHARACTER(LEN=7),DIMENSION(nbopp_max) :: sopp
94  REAL,DIMENSION(nbopp_max) :: scal
95!-External type (for R4/R8)
96  INTEGER :: v_typ
97!-Sizes of the associated grid and zommed area
98  INTEGER,DIMENSION(3) :: scsize,zorig,zsize
99!-Sizes for the data as it goes through the various math operations
100  INTEGER,DIMENSION(3) :: datasz_in = -1
101  INTEGER :: datasz_max = -1
102!-
103  INTEGER :: h_axid,z_axid,t_axid
104!-
105  REAL,DIMENSION(2) :: hist_minmax
106  LOGICAL :: hist_calc_rng=.FALSE.,hist_wrt_rng=.FALSE.
107!-Book keeping of the axes
108  INTEGER :: tdimid,tbndid=-1,tax_last
109  LOGICAL :: l_bnd
110  CHARACTER(LEN=40) :: tax_name
111!-
112  REAL :: freq_opp,freq_wrt
113  INTEGER :: &
114 &  last_opp,last_wrt,last_opp_chk,last_wrt_chk,nb_opp,nb_wrt
115!- For future optimization
116  REAL,POINTER,DIMENSION(:) :: t_bf
117!#  REAL,ALLOCATABLE,DIMENSION(:) :: V_1_D
118!#  REAL,ALLOCATABLE,DIMENSION(:,:) :: V_2_D
119!#  REAL,ALLOCATABLE,DIMENSION(:,:,:) :: V_3_D
120END TYPE T_D_V
121!-
122! File derived type
123!-
124TYPE :: T_D_F
125!-NETCDF IDs for file
126  INTEGER :: ncfid=-1
127!-Time variables
128  INTEGER :: itau0=0
129  REAL :: date0,deltat
130!-Counter of elements (variables, time-horizontal-vertical axis
131  INTEGER :: n_var=0,n_tax=0,n_hax=0,n_zax=0
132!-NETCDF dimension IDs for time-[time_bounds]-longitude-latitude
133  INTEGER :: tid,bid,xid,yid
134!-General definitions in the NETCDF file
135  INTEGER,DIMENSION(2) :: full_size=0,slab_ori,slab_siz
136!-The horizontal axes
137  CHARACTER(LEN=25),DIMENSION(nb_hax_max,2) :: hax_name
138!-The vertical axes
139  INTEGER,DIMENSION(nb_zax_max) :: zax_size,zax_ids
140  CHARACTER(LEN=20),DIMENSION(nb_zax_max) :: zax_name
141!-
142  LOGICAL :: regular=.TRUE.
143!-DOMAIN ID
144  INTEGER :: dom_id_svg=-1
145!-
146  TYPE(T_D_V),DIMENSION(nb_var_max) :: W_V
147END TYPE T_D_F
148!-
149TYPE(T_D_F),DIMENSION(nb_files_max),SAVE :: W_F
150!-
151! A list of functions which require special action
152! (Needs to be updated when functions are added
153!  but they are well located here)
154!-
155  CHARACTER(LEN=30),SAVE :: fuchnbout = 'scatter, fill'
156!- Some configurable variables with locks
157  CHARACTER(LEN=80),SAVE :: model_name='An IPSL model'
158  LOGICAL,SAVE :: lock_modname=.FALSE.
159!-
160!===
161CONTAINS
162!===
163!-
164SUBROUTINE histb_reg1d                 &
165 & (pfilename,pim,plon,pjm,plat,       &
166 &  par_orix,par_szx,par_oriy,par_szy, &
167 &  pitau0,pdate0,pdeltat,phoriid,idf,domain_id,mode,snc4chunks)
168!---------------------------------------------------------------------
169!- histbeg for 1D regular horizontal coordinates (see histb_all)
170!---------------------------------------------------------------------
171  IMPLICIT NONE
172!-
173  CHARACTER(LEN=*) :: pfilename
174  INTEGER,INTENT(IN) :: pim,pjm
175  REAL,DIMENSION(pim),INTENT(IN) :: plon
176  REAL,DIMENSION(pjm),INTENT(IN) :: plat
177  INTEGER,INTENT(IN):: par_orix,par_szx,par_oriy,par_szy
178  INTEGER,INTENT(IN) :: pitau0
179  REAL,INTENT(IN) :: pdate0,pdeltat
180  INTEGER,INTENT(OUT) :: idf,phoriid
181  INTEGER,INTENT(IN),OPTIONAL :: domain_id
182  CHARACTER(LEN=*),OPTIONAL,INTENT(IN) :: mode
183  TYPE(snc4_ctl), OPTIONAL,INTENT(IN) :: snc4chunks
184!---------------------------------------------------------------------
185  CALL histb_all &
186 & (1,pfilename,pim,pjm,pitau0,pdate0,pdeltat,phoriid,idf, &
187 &  x_1d=plon,y_1d=plat, &
188 &  k_orx=par_orix,k_szx=par_szx,k_ory=par_oriy,k_szy=par_szy, &
189 &  domain_id=domain_id,mode=mode,snc4chunks=snc4chunks)
190!-------------------------
191END SUBROUTINE histb_reg1d
192!===
193SUBROUTINE histb_reg2d &
194 & (pfilename,pim,plon,pjm,plat,       &
195 &  par_orix,par_szx,par_oriy,par_szy, &
196 &  pitau0,pdate0,pdeltat,phoriid,idf,domain_id,mode,snc4chunks)
197!---------------------------------------------------------------------
198!- histbeg for 2D regular horizontal coordinates (see histb_all)
199!---------------------------------------------------------------------
200  IMPLICIT NONE
201!-
202  CHARACTER(LEN=*) :: pfilename
203  INTEGER,INTENT(IN) :: pim,pjm
204  REAL,DIMENSION(pim,pjm),INTENT(IN) :: plon,plat
205  INTEGER,INTENT(IN):: par_orix,par_szx,par_oriy,par_szy
206  INTEGER,INTENT(IN) :: pitau0
207  REAL,INTENT(IN) :: pdate0,pdeltat
208  INTEGER,INTENT(OUT) :: idf,phoriid
209  INTEGER,INTENT(IN),OPTIONAL :: domain_id
210  CHARACTER(LEN=*),OPTIONAL,INTENT(IN) :: mode
211  TYPE(snc4_ctl), OPTIONAL,INTENT(IN) :: snc4chunks
212!---------------------------------------------------------------------
213  CALL histb_all &
214 & (2,pfilename,pim,pjm,pitau0,pdate0,pdeltat,phoriid,idf,  &
215 &  x_2d=plon,y_2d=plat, &
216 &  k_orx=par_orix,k_szx=par_szx,k_ory=par_oriy,k_szy=par_szy,    &
217 &  domain_id=domain_id,mode=mode,snc4chunks=snc4chunks)
218!-------------------------
219END SUBROUTINE histb_reg2d
220!===
221SUBROUTINE histb_irreg &
222 &  (pfilename,pim,plon,plon_bounds,plat,plat_bounds, &
223 &   pitau0,pdate0,pdeltat,phoriid,idf,domain_id,mode,snc4chunks)
224!---------------------------------------------------------------------
225!- histbeg for irregular horizontal coordinates (see histb_all)
226!---------------------------------------------------------------------
227  IMPLICIT NONE
228!-
229  CHARACTER(LEN=*) :: pfilename
230  INTEGER,INTENT(IN) :: pim
231  REAL,DIMENSION(pim),INTENT(IN) :: plon,plat
232  REAL,DIMENSION(:,:),INTENT(IN) :: plon_bounds,plat_bounds
233  INTEGER,INTENT(IN) :: pitau0
234  REAL,INTENT(IN) :: pdate0,pdeltat
235  INTEGER,INTENT(OUT) :: idf,phoriid
236  INTEGER,INTENT(IN),OPTIONAL :: domain_id
237  CHARACTER(LEN=*),OPTIONAL,INTENT(IN) :: mode
238  TYPE(snc4_ctl), OPTIONAL,INTENT(IN) :: snc4chunks
239!---------------------------------------------------------------------
240  CALL histb_all &
241 & (3,pfilename,pim,pim,pitau0,pdate0,pdeltat,phoriid,idf,  &
242 &  x_1d=plon,y_1d=plat,x_bnds=plon_bounds,y_bnds=plat_bounds, &
243 &  domain_id=domain_id,mode=mode,snc4chunks=snc4chunks)
244!-------------------------
245END SUBROUTINE histb_irreg
246!===
247SUBROUTINE histb_all &
248 & (k_typ,nc_name,pim,pjm,pitau0,pdate0,pdeltat,phoriid,idf, &
249 &  x_1d,y_1d,x_2d,y_2d,k_orx,k_szx,k_ory,k_szy, &
250 &  x_bnds,y_bnds,domain_id,mode,snc4chunks)
251!---------------------------------------------------------------------
252!- General interface for horizontal grids.
253!- This subroutine initializes a netcdf file and returns the ID.
254!- It will set up the geographical space on which the data will be
255!- stored and offers the possibility of seting a zoom.
256!- In the case of irregular grids, all the data comes in as vectors
257!- and for the grid we have the coordinates of the 4 corners.
258!- It also gets the global parameters into the I/O subsystem.
259!-
260!- INPUT
261!-
262!- k_typ    : Type of the grid (1 rectilinear, 2 regular, 3 irregular)
263!- nc_name  : Name of the netcdf file to be created
264!- pim      : Size of arrays in longitude direction
265!- pjm      : Size of arrays in latitude direction (pjm=pim for type 3)
266!-
267!- pitau0   : time step at which the history tape starts
268!- pdate0   : The Julian date at which the itau was equal to 0
269!- pdeltat  : Time step, in seconds, of the counter itau
270!-            used in histwrite for instance
271!-
272!- OUTPUT
273!-
274!- phoriid  : Identifier of the horizontal grid
275!- idf      : Identifier of the file
276!-
277!- Optional INPUT arguments
278!-
279!- For rectilinear or irregular grid
280!- x_1d  : The longitudes
281!- y_1d  : The latitudes
282!- For regular grid
283!- x_2d  : The longitudes
284!- y_2d  : The latitudes
285!- One pair (x_1d,y_1d) or (x_2d,y_2d) must be supplied.
286!-
287!- For regular grid (reg1d or reg2d),
288!- the next 4 arguments allow to define a horizontal zoom
289!- for this file. It is assumed that all variables to come
290!- have the same index space. This can not be assumed for
291!- the z axis and thus we define the zoom in histdef.
292!- k_orx  : Origin of the slab of data within the X axis (pim)
293!- k_szx  : Size of the slab of data in X
294!- k_ory  : Origin of the slab of data within the Y axis (pjm)
295!- k_szy  : Size of the slab of data in Y
296!-
297!- For irregular grid.
298!- x_bnds : The boundaries of the grid in longitude
299!- y_bnds : The boundaries of the grid in latitude
300!-
301!- For all grids.
302!-
303!- domain_id  : Domain identifier
304!-
305!- mode       : String of (case insensitive) blank-separated words
306!-              defining the mode used to create the file.
307!-              Supported keywords : 32, 64
308!-              "32/64" defines the offset mode.
309!-              The default offset mode is 64 bits.
310!-              Keywords "NETCDF4" and "CLASSIC" are reserved
311!-              for future use.
312!-
313!- snc4chunks    : Structure containing chunk partitioning parameters
314!-              for 4-D variables and a logical switch to toggle
315!-              between netcdf3 o/p (false) and netcdf4 chunked
316!-              and compressed o/p (true)
317!---------------------------------------------------------------------
318  IMPLICIT NONE
319!-
320  INTEGER,INTENT(IN) :: k_typ
321  CHARACTER(LEN=*),INTENT(IN) :: nc_name
322  INTEGER,INTENT(IN) :: pim,pjm
323  INTEGER,INTENT(IN) :: pitau0
324  REAL,INTENT(IN) :: pdate0,pdeltat
325  INTEGER,INTENT(OUT) :: idf,phoriid
326  REAL,DIMENSION(:),INTENT(IN),OPTIONAL :: x_1d,y_1d
327  REAL,DIMENSION(:,:),INTENT(IN),OPTIONAL :: x_2d,y_2d
328  INTEGER,INTENT(IN),OPTIONAL :: k_orx,k_szx,k_ory,k_szy
329  REAL,DIMENSION(:,:),INTENT(IN),OPTIONAL :: x_bnds,y_bnds
330  INTEGER,INTENT(IN),OPTIONAL :: domain_id
331  CHARACTER(LEN=*),OPTIONAL,INTENT(IN) :: mode
332  TYPE(snc4_ctl), OPTIONAL,INTENT(IN) :: snc4chunks
333!-
334  INTEGER :: nfid,iret,m_c
335  CHARACTER(LEN=120) :: file
336  CHARACTER(LEN=30) :: timenow
337  CHARACTER(LEN=11) :: c_nam
338  LOGICAL :: l_dbg
339!---------------------------------------------------------------------
340  CALL ipsldbg (old_status=l_dbg)
341!-
342  IF     (k_typ == 1) THEN
343    c_nam = 'histb_reg1d'
344  ELSEIF (k_typ == 2) THEN
345    c_nam = 'histb_reg2d'
346  ELSEIF (k_typ == 3) THEN
347    c_nam = 'histb_irreg'
348  ELSE
349    CALL ipslerr (3,"histbeg", &
350 &    'Illegal value of k_typ argument','in internal interface','?')
351  ENDIF
352!-
353  IF (l_dbg) WRITE(*,*) c_nam//" 0.0"
354!-
355! Search for a free index
356!-
357  idf = -1
358  DO nfid=1,nb_files_max
359    IF (W_F(nfid)%ncfid < 0) THEN
360      idf = nfid; EXIT;
361    ENDIF
362  ENDDO
363  IF (idf < 0) THEN
364    CALL ipslerr (3,"histbeg", &
365   &  'Table of files too small. You should increase nb_files_max', &
366   &  'in histcom.f90 in order to accomodate all these files',' ')
367  ENDIF
368!-
369! 1.0 Transfering into the common for future use
370!-
371  IF (l_dbg) WRITE(*,*) c_nam//" 1.0"
372!-
373  W_F(idf)%itau0  = pitau0
374  W_F(idf)%date0  = pdate0
375  W_F(idf)%deltat = pdeltat
376!-
377! 2.0 Initializes all variables for this file
378!-
379  IF (l_dbg) WRITE(*,*) c_nam//" 2.0"
380!-
381  W_F(idf)%n_var = 0
382  W_F(idf)%n_tax = 0
383  W_F(idf)%n_hax = 0
384  W_F(idf)%n_zax = 0
385!-
386  IF ( (k_typ == 1).OR.(k_typ == 2) ) THEN
387    W_F(idf)%slab_ori(1:2) = (/ k_orx,k_ory /)
388    W_F(idf)%slab_siz(1:2)  = (/ k_szx,k_szy /)
389  ELSE
390    W_F(idf)%slab_ori(1:2) = (/ 1,1 /)
391    W_F(idf)%slab_siz(1:2) = (/ pim,1 /)
392  ENDIF
393!-
394! 3.0 Opening netcdf file and defining dimensions
395!-
396  IF (l_dbg) WRITE(*,*) c_nam//" 3.0"
397!-
398! Add DOMAIN number and ".nc" suffix in file name if needed
399!-
400  file = nc_name
401  CALL flio_dom_file (file,domain_id)
402!-
403! Check the mode
404!? See fliocom for HDF4 ????????????????????????????????????????????????
405!-
406  IF (PRESENT(mode)) THEN
407    SELECT CASE (TRIM(mode))
408    CASE('32')
409      m_c = NF90_CLOBBER
410    CASE('64')
411      m_c = IOR(NF90_CLOBBER,NF90_64BIT_OFFSET)
412    CASE DEFAULT
413      CALL ipslerr (3,"histbeg", &
414 &      'Invalid argument mode for file :',TRIM(file), &
415 &      'Supported values are 32 or 64')
416    END SELECT
417  ELSE
418    m_c = IOR(NF90_CLOBBER,NF90_64BIT_OFFSET)
419  ENDIF
420!-
421  IF (PRESENT(snc4chunks)) THEN
422    IF (snc4chunks%luse) CALL get_nf90_symbol("NF90_HDF5", m_c)
423  ENDIF
424!-
425! Create file
426!-
427  iret = NF90_CREATE(file,m_c,nfid)
428!-
429  IF     (k_typ == 1) THEN
430    iret = NF90_DEF_DIM(nfid,'lon',k_szx,W_F(idf)%xid)
431    iret = NF90_DEF_DIM(nfid,'lat',k_szy,W_F(idf)%yid)
432  ELSEIF (k_typ == 2) THEN
433    iret = NF90_DEF_DIM(nfid,'x',k_szx,W_F(idf)%xid)
434    iret = NF90_DEF_DIM(nfid,'y',k_szy,W_F(idf)%yid)
435  ELSEIF (k_typ == 3) THEN
436    iret = NF90_DEF_DIM(nfid,'x',pim,W_F(idf)%xid)
437    W_F(idf)%yid = W_F(idf)%xid
438  ENDIF
439!-
440! 4.0 Declaring the geographical coordinates and other attributes
441!-
442  IF (l_dbg) WRITE(*,*) c_nam//" 4.0"
443!-
444  iret = NF90_PUT_ATT(nfid,NF90_GLOBAL,'Conventions','CF-1.1')
445  iret = NF90_PUT_ATT(nfid,NF90_GLOBAL,'file_name',TRIM(file))
446  iret = NF90_PUT_ATT(nfid,NF90_GLOBAL,'production',TRIM(model_name))
447  lock_modname = .TRUE.
448  CALL ioget_timestamp (timenow)
449  iret = NF90_PUT_ATT(nfid,NF90_GLOBAL,'TimeStamp',TRIM(timenow))
450!-
451! 5.0 Saving some important information on this file in the common
452!-
453  IF (l_dbg) WRITE(*,*) c_nam//" 5.0"
454!-
455  IF (PRESENT(domain_id)) THEN
456    W_F(idf)%dom_id_svg = domain_id
457  ENDIF
458  W_F(idf)%ncfid = nfid
459  IF ( (k_typ == 1).OR.(k_typ == 2) ) THEN
460    W_F(idf)%full_size(1:2) = (/ pim,pjm /)
461    W_F(idf)%regular=.TRUE.
462  ELSEIF (k_typ == 3) THEN
463    W_F(idf)%full_size(1:2) = (/ pim,1 /)
464    W_F(idf)%regular=.FALSE.
465  ENDIF
466!-
467! 6.0 storing the geographical coordinates
468!-
469  IF     (k_typ == 1) THEN
470    CALL histh_all &
471 &   (k_typ,idf,pim,pjm,' ','Default grid',phoriid, &
472 &    x_1d=x_1d,y_1d=y_1d)
473  ELSEIF (k_typ == 2) THEN
474    CALL histh_all &
475 &   (k_typ,idf,pim,pjm,' ','Default grid',phoriid, &
476 &    x_2d=x_2d,y_2d=y_2d)
477  ELSEIF (k_typ == 3) THEN
478    CALL histh_all &
479 &   (k_typ,idf,pim,pim,' ','Default grid',phoriid, &
480 &    x_1d=x_1d,y_1d=y_1d,x_bnds=x_bnds,y_bnds=y_bnds)
481  ENDIF
482!-----------------------
483END SUBROUTINE histb_all
484!===
485SUBROUTINE histh_reg1d &
486 &  (idf,pim,plon,pjm,plat,phname,phtitle,phid)
487!---------------------------------------------------------------------
488!- histhori for 1d regular grid (see histh_all)
489!---------------------------------------------------------------------
490  IMPLICIT NONE
491!-
492  INTEGER,INTENT(IN) :: idf,pim,pjm
493  REAL,INTENT(IN),DIMENSION(:) :: plon,plat
494  CHARACTER(LEN=*),INTENT(IN) :: phname,phtitle
495  INTEGER,INTENT(OUT) :: phid
496!---------------------------------------------------------------------
497  CALL histh_all &
498 & (1,idf,pim,pjm,phname,phtitle,phid,x_1d=plon,y_1d=plat)
499!-------------------------
500END SUBROUTINE histh_reg1d
501!===
502SUBROUTINE histh_reg2d &
503 & (idf,pim,plon,pjm,plat,phname,phtitle,phid)
504!---------------------------------------------------------------------
505!- histhori for 2d regular grid (see histh_all)
506!---------------------------------------------------------------------
507  IMPLICIT NONE
508!-
509  INTEGER,INTENT(IN) :: idf,pim,pjm
510  REAL,INTENT(IN),DIMENSION(:,:) :: plon,plat
511  CHARACTER(LEN=*),INTENT(IN) :: phname,phtitle
512  INTEGER,INTENT(OUT) :: phid
513!---------------------------------------------------------------------
514  CALL histh_all &
515 & (2,idf,pim,pjm,phname,phtitle,phid,x_2d=plon,y_2d=plat)
516!-------------------------
517END SUBROUTINE histh_reg2d
518!===
519SUBROUTINE histh_irreg &
520 & (idf,pim,plon,plon_bounds,plat,plat_bounds,phname,phtitle,phid)
521!---------------------------------------------------------------------
522!- histhori for irregular grid (see histh_all)
523!---------------------------------------------------------------------
524  IMPLICIT NONE
525!-
526  INTEGER,INTENT(IN) :: idf,pim
527  REAL,DIMENSION(:),INTENT(IN) :: plon,plat
528  REAL,DIMENSION(:,:),INTENT(IN) :: plon_bounds,plat_bounds
529  CHARACTER(LEN=*),INTENT(IN) :: phname,phtitle
530  INTEGER,INTENT(OUT) :: phid
531!---------------------------------------------------------------------
532  CALL histh_all &
533 & (3,idf,pim,pim,phname,phtitle,phid, &
534 &  x_1d=plon,y_1d=plat,x_bnds=plon_bounds,y_bnds=plat_bounds)
535!-------------------------
536END SUBROUTINE histh_irreg
537!===
538SUBROUTINE histh_all &
539 & (k_typ,idf,pim,pjm,phname,phtitle,phid, &
540 &  x_1d,y_1d,x_2d,y_2d,x_bnds,y_bnds)
541!---------------------------------------------------------------------
542!- General interface for horizontal grids.
543!- This subroutine is made to declare a new horizontal grid.
544!- It has to have the same number of points as
545!- the original and thus in this routine we will only
546!- add two variable (longitude and latitude).
547!- Any variable in the file can thus point to this pair
548!- through an attribute. This routine is very usefull
549!- to allow staggered grids.
550!-
551!- INPUT
552!-
553!- k_typ   : Type of the grid (1 rectilinear, 2 regular, 3 irregular)
554!- idf     : The id of the file to which the grid should be added
555!- pim     : Size in the longitude direction
556!- pjm     : Size in the latitude direction (pjm=pim for type 3)
557!- phname  : The name of grid
558!- phtitle : The title of the grid
559!-
560!- OUTPUT
561!-
562!- phid    : Id of the created grid
563!-
564!- Optional INPUT arguments
565!-
566!- For rectilinear or irregular grid
567!- x_1d  : The longitudes
568!- y_1d  : The latitudes
569!- For regular grid
570!- x_2d  : The longitudes
571!- y_2d  : The latitudes
572!- One pair (x_1d,y_1d) or (x_2d,y_2d) must be supplied.
573!-
574!- For irregular grid.
575!- x_bnds : The boundaries of the grid in longitude
576!- y_bnds : The boundaries of the grid in latitude
577!---------------------------------------------------------------------
578  IMPLICIT NONE
579!-
580  INTEGER,INTENT(IN) :: k_typ
581  INTEGER,INTENT(IN) :: idf,pim,pjm
582  CHARACTER(LEN=*),INTENT(IN) :: phname,phtitle
583  INTEGER,INTENT(OUT) :: phid
584  REAL,DIMENSION(:),INTENT(IN),OPTIONAL :: x_1d,y_1d
585  REAL,DIMENSION(:,:),INTENT(IN),OPTIONAL :: x_2d,y_2d
586  REAL,DIMENSION(:,:),INTENT(IN),OPTIONAL :: x_bnds,y_bnds
587!-
588  CHARACTER(LEN=25) :: lon_name,lat_name
589  CHARACTER(LEN=30) :: lonbound_name,latbound_name
590  INTEGER :: i_s,i_e
591  INTEGER,DIMENSION(2) :: dims,dims_b
592  INTEGER :: nbbounds
593  INTEGER :: nlonidb,nlatidb,twoid
594  LOGICAL :: transp = .FALSE.
595  REAL,ALLOCATABLE,DIMENSION(:,:) :: bounds_trans
596  REAL :: wmn,wmx
597  INTEGER :: nlonid,nlatid
598  INTEGER :: o_x,o_y,s_x,s_y
599  INTEGER :: iret,nfid
600  CHARACTER(LEN=11) :: c_nam
601  LOGICAL :: l_dbg
602!---------------------------------------------------------------------
603  CALL ipsldbg (old_status=l_dbg)
604!-
605  IF     (k_typ == 1) THEN
606    c_nam = 'histh_reg1d'
607  ELSEIF (k_typ == 2) THEN
608    c_nam = 'histh_reg2d'
609  ELSEIF (k_typ == 3) THEN
610    c_nam = 'histh_irreg'
611  ELSE
612    CALL ipslerr (3,"histhori", &
613 &    'Illegal value of k_typ argument','in internal interface','?')
614  ENDIF
615!-
616! 1.0 Check that all fits in the buffers
617!-
618  IF (    (pim /= W_F(idf)%full_size(1)) &
619 &    .OR.(W_F(idf)%regular.AND.(pjm /= W_F(idf)%full_size(2)))  &
620 &    .OR.(.NOT.W_F(idf)%regular.AND.(W_F(idf)%full_size(2) /= 1)) ) THEN
621    CALL ipslerr (3,"histhori", &
622 &   'The new horizontal grid does not have the same size', &
623 &   'as the one provided to histbeg. This is not yet ', &
624 &   'possible in the hist package.')
625  ENDIF
626!-
627! 1.1 Create all the variables needed
628!-
629  IF (l_dbg) WRITE(*,*) c_nam//" 1.0"
630!-
631  nfid = W_F(idf)%ncfid
632!-
633  IF (k_typ == 3) THEN
634    IF     (SIZE(x_bnds,DIM=1) == pim) THEN
635      nbbounds = SIZE(x_bnds,DIM=2)
636      transp = .TRUE.
637    ELSEIF (SIZE(x_bnds,DIM=2) == pim) THEN
638      nbbounds = SIZE(x_bnds,DIM=1)
639      transp = .FALSE.
640    ELSE
641      CALL ipslerr (3,"histhori", &
642 &     'The boundary variable does not have any axis corresponding', &
643 &     'to the size of the longitude or latitude variable','.')
644    ENDIF
645    ALLOCATE(bounds_trans(nbbounds,pim))
646    iret = NF90_DEF_DIM(nfid,'nbnd',nbbounds,twoid)
647    dims_b(1:2) = (/ twoid,W_F(idf)%xid /)
648  ENDIF
649!-
650  dims(1:2) = (/ W_F(idf)%xid,W_F(idf)%yid /)
651!-
652  IF     (k_typ == 1) THEN
653    IF (W_F(idf)%n_hax == 0) THEN
654      lon_name = 'lon'
655      lat_name = 'lat'
656    ELSE
657      lon_name = 'lon_'//TRIM(phname)
658      lat_name = 'lat_'//TRIM(phname)
659    ENDIF
660  ELSEIF (k_typ == 2) THEN
661    IF (W_F(idf)%n_hax == 0) THEN
662      lon_name = 'nav_lon'
663      lat_name = 'nav_lat'
664    ELSE
665      lon_name = 'nav_lon_'//TRIM(phname)
666      lat_name = 'nav_lat_'//TRIM(phname)
667    ENDIF
668  ELSEIF (k_typ == 3) THEN
669    IF (W_F(idf)%n_hax == 0) THEN
670      lon_name = 'nav_lon'
671      lat_name = 'nav_lat'
672    ELSE
673      lon_name = 'nav_lon_'//TRIM(phname)
674      lat_name = 'nav_lat_'//TRIM(phname)
675    ENDIF
676    lonbound_name = TRIM(lon_name)//'_bounds'
677    latbound_name = TRIM(lat_name)//'_bounds'
678  ENDIF
679!-
680! 1.2 Save the informations
681!-
682  phid = W_F(idf)%n_hax+1
683  W_F(idf)%n_hax = phid
684  W_F(idf)%hax_name(phid,1:2) = (/ lon_name,lat_name /)
685!-
686! 2.0 Longitude
687!-
688  IF (l_dbg) WRITE(*,*) c_nam//" 2.0"
689!-
690  i_s = 1;
691  IF ( (k_typ == 1).OR.(k_typ == 3) ) THEN
692    i_e = 1; wmn = MINVAL(x_1d); wmx = MAXVAL(x_1d);
693  ELSEIF (k_typ == 2) THEN
694    i_e = 2; wmn = MINVAL(x_2d); wmx = MAXVAL(x_2d);
695  ENDIF
696  iret = NF90_DEF_VAR(nfid,lon_name,NF90_REAL4,dims(i_s:i_e),nlonid)
697  IF (k_typ == 1) THEN
698    iret = NF90_PUT_ATT(nfid,nlonid,'axis',"X")
699  ENDIF
700  iret = NF90_PUT_ATT(nfid,nlonid,'standard_name',"longitude")
701  iret = NF90_PUT_ATT(nfid,nlonid,'units',"degrees_east")
702  iret = NF90_PUT_ATT(nfid,nlonid,'valid_min',REAL(wmn,KIND=4))
703  iret = NF90_PUT_ATT(nfid,nlonid,'valid_max',REAL(wmx,KIND=4))
704  iret = NF90_PUT_ATT(nfid,nlonid,'long_name',"Longitude")
705  iret = NF90_PUT_ATT(nfid,nlonid,'nav_model',TRIM(phtitle))
706!-
707  IF (k_typ == 3) THEN
708!---
709!-- 2.1 Longitude bounds
710!---
711    iret = NF90_PUT_ATT(nfid,nlonid,'bounds',TRIM(lonbound_name))
712    iret = NF90_DEF_VAR(nfid,lonbound_name,NF90_REAL4,dims_b(1:2),nlonidb)
713    iret = NF90_PUT_ATT(nfid,nlonidb,'long_name', &
714 &          'Boundaries for coordinate variable '//TRIM(lon_name))
715  ENDIF
716!-
717! 3.0 Latitude
718!-
719  IF (l_dbg) WRITE(*,*) c_nam//" 3.0"
720!-
721  i_e = 2;
722  IF ( (k_typ == 1).OR.(k_typ == 3) ) THEN
723    i_s = 2; wmn = MINVAL(y_1d); wmx = MAXVAL(y_1d);
724  ELSEIF (k_typ == 2) THEN
725    i_s = 1; wmn = MINVAL(y_2d); wmx = MAXVAL(y_2d);
726  ENDIF
727  iret = NF90_DEF_VAR(nfid,lat_name,NF90_REAL4,dims(i_s:i_e),nlatid)
728  IF (k_typ == 1) THEN
729    iret = NF90_PUT_ATT(nfid,nlatid,'axis',"Y")
730  ENDIF
731!-
732  iret = NF90_PUT_ATT(nfid,nlatid,'standard_name',"latitude")
733  iret = NF90_PUT_ATT(nfid,nlatid,'units',"degrees_north")
734  iret = NF90_PUT_ATT(nfid,nlatid,'valid_min',REAL(wmn,KIND=4))
735  iret = NF90_PUT_ATT(nfid,nlatid,'valid_max',REAL(wmx,KIND=4))
736  iret = NF90_PUT_ATT(nfid,nlatid,'long_name',"Latitude")
737  iret = NF90_PUT_ATT(nfid,nlatid,'nav_model',TRIM(phtitle))
738!-
739  IF (k_typ == 3) THEN
740!---
741!-- 3.1 Latitude bounds
742!---
743    iret = NF90_PUT_ATT(nfid,nlatid,'bounds',TRIM(latbound_name))
744    iret = NF90_DEF_VAR(nfid,latbound_name,NF90_REAL4,dims_b(1:2),nlatidb)
745    iret = NF90_PUT_ATT(nfid,nlatidb,'long_name', &
746 &          'Boundaries for coordinate variable '//TRIM(lat_name))
747  ENDIF
748!-
749  iret = NF90_ENDDEF(nfid)
750!-
751! 4.0 storing the geographical coordinates
752!-
753  IF (l_dbg) WRITE(*,*) c_nam//" 4.0"
754!-
755  IF ( (k_typ == 1).OR.(k_typ == 2) ) THEN
756    o_x = W_F(idf)%slab_ori(1)
757    o_y = W_F(idf)%slab_ori(2)
758    s_x = W_F(idf)%slab_siz(1)
759    s_y = W_F(idf)%slab_siz(2)
760!---
761!-- Transfer the longitude and  the latitude
762!---
763    IF     (k_typ == 1) THEN
764      iret = NF90_PUT_VAR(nfid,nlonid,x_1d(o_x:o_x+s_x-1))
765      iret = NF90_PUT_VAR(nfid,nlatid,y_1d(o_y:o_y+s_y-1))
766    ELSEIF (k_typ == 2) THEN
767      iret = NF90_PUT_VAR(nfid,nlonid, &
768 &            x_2d(o_x:o_x+s_x-1,o_y:o_y+s_y-1))
769      iret = NF90_PUT_VAR(nfid,nlatid, &
770 &            y_2d(o_x:o_x+s_x-1,o_y:o_y+s_y-1))
771    ENDIF
772  ELSEIF (k_typ == 3) THEN
773!---
774!-- Transfer the longitude and the longitude bounds
775!---
776    iret = NF90_PUT_VAR(nfid,nlonid,x_1d(1:pim))
777!---
778    IF (transp) THEN
779      bounds_trans = TRANSPOSE(x_bnds)
780    ELSE
781      bounds_trans = x_bnds
782    ENDIF
783    iret = NF90_PUT_VAR(nfid,nlonidb,bounds_trans(1:nbbounds,1:pim))
784!---
785!-- Transfer the latitude and the latitude bounds
786!---
787    iret = NF90_PUT_VAR(nfid,nlatid,y_1d(1:pim))
788!---
789    IF (transp) THEN
790      bounds_trans = TRANSPOSE(y_bnds)
791    ELSE
792      bounds_trans = y_bnds
793    ENDIF
794    iret = NF90_PUT_VAR(nfid,nlatidb,bounds_trans(1:nbbounds,1:pim))
795!---
796    DEALLOCATE(bounds_trans)
797  ENDIF
798!-
799  iret = NF90_REDEF(nfid)
800!-----------------------
801END SUBROUTINE histh_all
802!===
803SUBROUTINE histvert (idf,pzaxname,pzaxtitle,pzaxunit, &
804 &                   pzsize,pzvalues,pzaxid,pdirect)
805!---------------------------------------------------------------------
806!- This subroutine defines a vertical axis and returns it s id.
807!- It gives the user the possibility to the user to define many
808!- different vertical axes. For each variable defined with histdef a
809!- vertical axis can be specified with by it s ID.
810!-
811!- INPUT
812!-
813!- idf      : ID of the file the variable should be archived in
814!- pzaxname : Name of the vertical axis
815!- pzaxtitle: title of the vertical axis
816!- pzaxunit : Units of the vertical axis (no units if blank string)
817!- pzsize   : size of the vertical axis
818!- pzvalues : Coordinate values of the vetical axis
819!-
820!- pdirect  : is an optional argument which allows to specify the
821!-            the positive direction of the axis : up or down.
822!- OUTPUT
823!-
824!- pzaxid   : Returns the ID of the axis.
825!-            Note that this is not the netCDF ID !
826!-
827!- VERSION
828!-
829!---------------------------------------------------------------------
830  IMPLICIT NONE
831!-
832  INTEGER,INTENT(IN) :: idf,pzsize
833  CHARACTER(LEN=*),INTENT(IN) :: pzaxname,pzaxunit,pzaxtitle
834  REAL,INTENT(IN) :: pzvalues(pzsize)
835  INTEGER,INTENT(OUT) :: pzaxid
836  CHARACTER(LEN=*),INTENT(IN),OPTIONAL :: pdirect
837!-
838  INTEGER :: pos,iv,zdimid,zaxid_tmp
839  CHARACTER(LEN=70) :: str71
840  CHARACTER(LEN=20) :: direction
841  INTEGER :: iret,leng,nfid
842  LOGICAL :: l_dbg
843!---------------------------------------------------------------------
844  CALL ipsldbg (old_status=l_dbg)
845!-
846! 1.0 Verifications :
847!    Do we have enough space for an extra axis ?
848!    Is the name already in use ?
849!-
850  IF (l_dbg) WRITE(*,*) "histvert : 1.0 Verifications", &
851 &                      pzaxname,'---',pzaxunit,'---',pzaxtitle
852!-
853! Direction of the vertical axis. Can we get if from the user.
854!-
855  IF (PRESENT(pdirect)) THEN
856    direction = TRIM(pdirect)
857    CALL strlowercase (direction)
858  ELSE
859    direction = 'unknown'
860  ENDIF
861!-
862! Check the consistency of the attribute
863!-
864  IF (     PRESENT(pdirect)    &
865 &    .AND.(direction /= 'up') &
866 &    .AND.(direction /= 'down') ) THEN
867    direction = 'unknown'
868    CALL ipslerr (2,"histvert",&
869 & "The specified positive direction for the vertical axis is invalid.",&
870 & "The value must be up or down.","The attribute will not be written.")
871  ENDIF
872!-
873  IF (W_F(idf)%n_zax+1 > nb_zax_max) THEN
874    CALL ipslerr (3,"histvert", &
875   &  'Table of vertical axes too small. You should increase ',&
876   &  'nb_zax_max in histcom.f90 in order to accomodate all ', &
877   &  'these variables ')
878  ENDIF
879!-
880  iv = W_F(idf)%n_zax
881  IF (iv > 1) THEN
882    CALL find_str (W_F(idf)%zax_name(1:iv-1),pzaxname,pos)
883  ELSE
884    pos = 0
885  ENDIF
886!-
887  IF (pos > 0) THEN
888    WRITE(str71,'("Check variable ",A," in file",I3)') &
889 &    TRIM(pzaxname),idf
890    CALL ipslerr (3,"histvert", &
891 &    "Vertical axis already exists",TRIM(str71), &
892 &    "Can also be a wrong file ID in another declaration")
893  ENDIF
894!-
895  iv = W_F(idf)%n_zax+1
896!-
897! 2.0 Add the information to the file
898!-
899  IF (l_dbg) &
900 &  WRITE(*,*) "histvert : 2.0 Add the information to the file"
901!-
902  nfid = W_F(idf)%ncfid
903!-
904  leng = MIN(LEN_TRIM(pzaxname),20)
905  iret = NF90_DEF_DIM (nfid,pzaxname(1:leng),pzsize,zaxid_tmp)
906  iret = NF90_DEF_VAR (nfid,pzaxname(1:leng),NF90_REAL4, &
907 &                     zaxid_tmp,zdimid)
908  iret = NF90_PUT_ATT (nfid,zdimid,'axis',"Z")
909  iret = NF90_PUT_ATT (nfid,zdimid,'standard_name',"model_level_number")
910  leng = MIN(LEN_TRIM(pzaxunit),20)
911  IF (leng > 0) THEN
912    iret = NF90_PUT_ATT (nfid,zdimid,'units',pzaxunit(1:leng))
913  ENDIF
914  IF (direction /= 'unknown') THEN
915    iret = NF90_PUT_ATT (nfid,zdimid,'positive',TRIM(direction))
916  ENDIF
917  iret = NF90_PUT_ATT (nfid,zdimid,'valid_min', &
918 &                     REAL(MINVAL(pzvalues(1:pzsize)),KIND=4))
919  iret = NF90_PUT_ATT (nfid,zdimid,'valid_max', &
920 &                     REAL(MAXVAL(pzvalues(1:pzsize)),KIND=4))
921  leng = MIN(LEN_TRIM(pzaxname),20)
922  iret = NF90_PUT_ATT (nfid,zdimid,'title',pzaxname(1:leng))
923  leng = MIN(LEN_TRIM(pzaxtitle),80)
924  iret = NF90_PUT_ATT (nfid,zdimid,'long_name',pzaxtitle(1:leng))
925!-
926  iret = NF90_ENDDEF (nfid)
927!-
928  iret = NF90_PUT_VAR (nfid,zdimid,pzvalues(1:pzsize))
929!-
930  iret = NF90_REDEF (nfid)
931!-
932!- 3.0 add the information to the common
933!-
934  IF (l_dbg) &
935  &  WRITE(*,*) "histvert : 3.0 add the information to the common"
936!-
937  W_F(idf)%n_zax = iv
938  W_F(idf)%zax_size(iv) = pzsize
939  W_F(idf)%zax_name(iv) = pzaxname
940  W_F(idf)%zax_ids(iv) = zaxid_tmp
941  pzaxid = iv
942!----------------------
943END SUBROUTINE histvert
944!===
945SUBROUTINE histdef &
946 &  (idf,pvarname,ptitle,punit, &
947 &   pxsize,pysize,phoriid,pzsize,par_oriz,par_szz,pzid, &
948 &   xtype,popp,pfreq_opp,pfreq_wrt,var_range,standard_name)
949!---------------------------------------------------------------------
950!- With this subroutine each variable to be archived on the history
951!- tape should be declared.
952!-
953!- It gives the user the choise of operation
954!- to be performed on the variables, the frequency of this operation
955!- and finaly the frequency of the archiving.
956!-
957!- INPUT
958!-
959!- idf      : ID of the file the variable should be archived in
960!- pvarname : Name of the variable, short and easy to remember
961!- ptitle   : Full name of the variable
962!- punit    : Units of the variable (no units if blank string)
963!-
964!- The next 3 arguments give the size of that data
965!- that will be passed to histwrite. The zoom will be
966!- done there with the horizontal information obtained
967!- in histbeg and the vertical information to follow.
968!-
969!- pxsize   : Size in X direction (size of the data that will be
970!-            given to histwrite)
971!- pysize   : Size in Y direction
972!- phoriid  : ID of the horizontal axis
973!-
974!- The next two arguments give the vertical zoom to use.
975!-
976!- pzsize   : Size in Z direction (If 1 then no axis is declared
977!-            for this variable and pzid is not used)
978!- par_oriz : Off set of the zoom
979!- par_szz  : Size of the zoom
980!-
981!- pzid     : ID of the vertical axis to use. It has to have
982!-            the size of the zoom.
983!- xtype    : External netCDF type (hist_r4/hist_r8)
984!- popp     : Operation to be performed. The following options
985!-            exist today :
986!- inst : keeps instantaneous values for writting
987!- ave  : Computes the average from call between writes
988!- pfreq_opp: Frequency of this operation (in seconds)
989!- pfreq_wrt: Frequency at which the variable should be
990!-            written (in seconds)
991!- var_range: Range of the variable.
992!-            If the minimum is greater than the maximum,
993!-            the values will be calculated.
994!-
995!- VERSION
996!---------------------------------------------------------------------
997  IMPLICIT NONE
998!-
999  INTEGER,INTENT(IN) :: idf,pxsize,pysize,pzsize,pzid
1000  INTEGER,INTENT(IN) :: par_oriz,par_szz,xtype,phoriid
1001  CHARACTER(LEN=*),INTENT(IN) :: pvarname,punit,popp,ptitle
1002  REAL,INTENT(IN) :: pfreq_opp,pfreq_wrt
1003  REAL,DIMENSION(2),OPTIONAL,INTENT(IN) :: var_range
1004  CHARACTER(LEN=*),OPTIONAL,INTENT(IN) :: standard_name
1005!-
1006  INTEGER :: iv
1007  CHARACTER(LEN=70) :: str70,str71,str72
1008  CHARACTER(LEN=20) :: tmp_name
1009  CHARACTER(LEN=40) :: str40
1010  CHARACTER(LEN=10) :: str10
1011  CHARACTER(LEN=120) :: ex_topps
1012  REAL :: un_an,un_jour,test_fopp,test_fwrt
1013  INTEGER :: pos,buff_sz
1014  LOGICAL :: l_dbg
1015!---------------------------------------------------------------------
1016  CALL ipsldbg (old_status=l_dbg)
1017!-
1018  ex_topps = 'ave, inst, t_min, t_max, t_sum, once, never, l_max, l_min'
1019!-
1020  W_F(idf)%n_var = W_F(idf)%n_var+1
1021  iv = W_F(idf)%n_var
1022!-
1023  IF (iv > nb_var_max) THEN
1024    CALL ipslerr (3,"histdef", &
1025   &  'Table of variables too small. You should increase nb_var_max',&
1026   &  'in histcom.f90 in order to accomodate all these variables', &
1027   &  ' ')
1028  ENDIF
1029!-
1030! 1.0 Transfer informations on the variable to the common
1031!     and verify that it does not already exist
1032!-
1033  IF (l_dbg) WRITE(*,*) "histdef : 1.0"
1034!-
1035  IF (iv > 1) THEN
1036    CALL find_str (W_F(idf)%W_V(1:iv-1)%v_name,pvarname,pos)
1037  ELSE
1038    pos = 0
1039  ENDIF
1040!-
1041  IF (pos > 0) THEN
1042    str70 = "Variable already exists"
1043    WRITE(str71,'("Check variable  ",a," in file",I3)') &
1044 &    TRIM(pvarname),idf
1045    str72 = "Can also be a wrong file ID in another declaration"
1046    CALL ipslerr (3,"histdef",str70,str71,str72)
1047  ENDIF
1048!-
1049  W_F(idf)%W_V(iv)%v_name = pvarname
1050  W_F(idf)%W_V(iv)%title = ptitle
1051  W_F(idf)%W_V(iv)%unit_name = punit
1052  IF (PRESENT(standard_name)) THEN
1053    W_F(idf)%W_V(iv)%std_name = standard_name
1054  ELSE
1055    W_F(idf)%W_V(iv)%std_name = ptitle
1056  ENDIF
1057  tmp_name = W_F(idf)%W_V(iv)%v_name
1058!-
1059! 1.1 decode the operations
1060!-
1061  W_F(idf)%W_V(iv)%fullop = popp
1062  CALL buildop &
1063 &  (TRIM(popp),ex_topps,W_F(idf)%W_V(iv)%topp,missing_val, &
1064 &   W_F(idf)%W_V(iv)%sopp,W_F(idf)%W_V(iv)%scal, &
1065 &   W_F(idf)%W_V(iv)%nbopp)
1066!-
1067! 1.2 If we have an even number of operations
1068!     then we need to add identity
1069!-
1070  IF ( MOD(W_F(idf)%W_V(iv)%nbopp,2) == 0) THEN
1071    W_F(idf)%W_V(iv)%nbopp = W_F(idf)%W_V(iv)%nbopp+1
1072    W_F(idf)%W_V(iv)%sopp(W_F(idf)%W_V(iv)%nbopp) = 'ident'
1073    W_F(idf)%W_V(iv)%scal(W_F(idf)%W_V(iv)%nbopp) = missing_val
1074  ENDIF
1075!-
1076! 1.3 External type of the variable
1077!-
1078  IF (xtype == hist_r8) THEN
1079    W_F(idf)%W_V(iv)%v_typ = hist_r8
1080  ELSE
1081    W_F(idf)%W_V(iv)%v_typ = hist_r4
1082  ENDIF
1083!-
1084! 2.0 Put the size of the variable in the common and check
1085!-
1086  IF (l_dbg) THEN
1087    WRITE(*,*) "histdef : 2.0",idf,iv,W_F(idf)%W_V(iv)%nbopp, &
1088 &    W_F(idf)%W_V(iv)%sopp(1:W_F(idf)%W_V(iv)%nbopp), &
1089 &    W_F(idf)%W_V(iv)%scal(1:W_F(idf)%W_V(iv)%nbopp)
1090  ENDIF
1091!-
1092  W_F(idf)%W_V(iv)%scsize(1:3) = (/ pxsize,pysize,pzsize /)
1093  W_F(idf)%W_V(iv)%zorig(1:3) = &
1094 &  (/ W_F(idf)%slab_ori(1),W_F(idf)%slab_ori(2),par_oriz /)
1095  W_F(idf)%W_V(iv)%zsize(1:3) = &
1096 &  (/ W_F(idf)%slab_siz(1),W_F(idf)%slab_siz(2),par_szz /)
1097!-
1098! Is the size of the full array the same as that of the coordinates ?
1099!-
1100  IF (    (pxsize > W_F(idf)%full_size(1)) &
1101 &    .OR.(pysize > W_F(idf)%full_size(2)) ) THEN
1102!-
1103    str70 = "The size of the variable is different "// &
1104 &          "from the one of the coordinates"
1105    WRITE(str71,'("Size of coordinates :",2I4)') &
1106 &   W_F(idf)%full_size(1),W_F(idf)%full_size(2)
1107    WRITE(str72,'("Size declared for variable ",a," :",2I4)') &
1108 &   TRIM(tmp_name),pxsize,pysize
1109    CALL ipslerr (3,"histdef",str70,str71,str72)
1110  ENDIF
1111!-
1112! Is the size of the zoom smaller than the coordinates ?
1113!-
1114  IF (    (W_F(idf)%full_size(1) < W_F(idf)%slab_siz(1)) &
1115 &    .OR.(W_F(idf)%full_size(2) < W_F(idf)%slab_siz(2)) ) THEN
1116    str70 = &
1117 &   "Size of variable should be greater or equal to those of the zoom"
1118    WRITE(str71,'("Size of XY zoom :",2I4)') &
1119 &   W_F(idf)%slab_siz(1),W_F(idf)%slab_siz(2)
1120    WRITE(str72,'("Size declared for variable ",A," :",2I4)') &
1121 &   TRIM(tmp_name),pxsize,pysize
1122    CALL ipslerr (3,"histdef",str70,str71,str72)
1123  ENDIF
1124!-
1125! 2.1 We store the horizontal grid information with minimal
1126!     and a fall back onto the default grid
1127!-
1128  IF ( (phoriid > 0).AND.(phoriid <= W_F(idf)%n_hax) ) THEN
1129    W_F(idf)%W_V(iv)%h_axid = phoriid
1130  ELSE
1131    W_F(idf)%W_V(iv)%h_axid = 1
1132    CALL ipslerr (2,"histdef", &
1133   &  'We use the default grid for variable as an invalide',&
1134   &  'ID was provided for variable : ',TRIM(pvarname))
1135  ENDIF
1136!-
1137! 2.2 Check the vertical coordinates if needed
1138!-
1139  IF (par_szz > 1) THEN
1140!-
1141!-- Does the vertical coordinate exist ?
1142!-
1143    IF (pzid > W_F(idf)%n_zax) THEN
1144      WRITE(str70, &
1145 &    '("The vertical coordinate chosen for variable ",A)') &
1146 &     TRIM(tmp_name)
1147      str71 = " Does not exist."
1148      CALL ipslerr (3,"histdef",str70,str71," ")
1149    ENDIF
1150!-
1151!-- Is the vertical size of the variable equal to that of the axis ?
1152!-
1153    IF (par_szz /= W_F(idf)%zax_size(pzid)) THEN
1154      str70 = "The size of the zoom does not correspond "// &
1155 &            "to the size of the chosen vertical axis"
1156      WRITE(str71,'("Size of zoom in z :",I4)') par_szz
1157      WRITE(str72,'("Size declared for axis ",A," :",I4)') &
1158 &     TRIM(W_F(idf)%zax_name(pzid)),W_F(idf)%zax_size(pzid)
1159      CALL ipslerr (3,"histdef",str70,str71,str72)
1160    ENDIF
1161!-
1162!-- Is the zoom smaller that the total size of the variable ?
1163!-
1164    IF (pzsize < par_szz) THEN
1165      str70 = "The vertical size of variable "// &
1166 &            "is smaller than that of the zoom."
1167      WRITE(str71,'("Declared vertical size of data :",I5)') pzsize
1168      WRITE(str72,'("Size of zoom for variable ",a," = ",I5)') &
1169 &     TRIM(tmp_name),par_szz
1170      CALL ipslerr (3,"histdef",str70,str71,str72)
1171    ENDIF
1172    W_F(idf)%W_V(iv)%z_axid = pzid
1173  ELSE
1174    W_F(idf)%W_V(iv)%z_axid = -99
1175  ENDIF
1176!-
1177! 3.0 We get the size of the arrays histwrite will get
1178!     and eventually allocate the time_buffer
1179!-
1180  IF (l_dbg) THEN
1181    WRITE(*,*) "histdef : 3.0"
1182  ENDIF
1183!-
1184  buff_sz = W_F(idf)%W_V(iv)%zsize(1) &
1185 &         *W_F(idf)%W_V(iv)%zsize(2) &
1186 &         *W_F(idf)%W_V(iv)%zsize(3)
1187!-
1188  IF (     (TRIM(W_F(idf)%W_V(iv)%topp) /= "inst") &
1189 &    .AND.(TRIM(W_F(idf)%W_V(iv)%topp) /= "once") &
1190 &    .AND.(TRIM(W_F(idf)%W_V(iv)%topp) /= "never") )THEN
1191    ALLOCATE(W_F(idf)%W_V(iv)%t_bf(buff_sz))
1192    W_F(idf)%W_V(iv)%t_bf(:) = 0.
1193    IF (l_dbg) THEN
1194      WRITE(*,*) "histdef : 3.0 allocating time_buffer for", &
1195 &      " idf = ",idf," iv = ",iv," size = ",buff_sz
1196    ENDIF
1197  ENDIF
1198!-
1199! 4.0 Transfer the frequency of the operations and check
1200!     for validity. We have to pay attention to negative values
1201!     of the frequency which indicate monthly time-steps.
1202!     The strategy is to bring it back to seconds for the tests
1203!-
1204  IF (l_dbg) WRITE(*,*) "histdef : 4.0"
1205!-
1206  W_F(idf)%W_V(iv)%freq_opp = pfreq_opp
1207  W_F(idf)%W_V(iv)%freq_wrt = pfreq_wrt
1208!-
1209  CALL ioget_calendar(un_an,un_jour)
1210  IF (pfreq_opp < 0) THEN
1211    CALL ioget_calendar(un_an)
1212    test_fopp = pfreq_opp*(-1.)*un_an/12.*un_jour
1213  ELSE
1214    test_fopp = pfreq_opp
1215  ENDIF
1216  IF (pfreq_wrt < 0) THEN
1217    CALL ioget_calendar(un_an)
1218    test_fwrt = pfreq_wrt*(-1.)*un_an/12.*un_jour
1219  ELSE
1220    test_fwrt = pfreq_wrt
1221  ENDIF
1222!-
1223! 4.1 Frequency of operations and output should be larger than deltat !
1224!-
1225  IF (test_fopp < W_F(idf)%deltat) THEN
1226    str70 = 'Frequency of operations should be larger than deltat'
1227    WRITE(str71,'("It is not the case for variable ",a," :",F10.4)') &
1228 &   TRIM(tmp_name),pfreq_opp
1229    str72 = "PATCH : frequency set to deltat"
1230!-
1231    CALL ipslerr (2,"histdef",str70,str71,str72)
1232!-
1233    W_F(idf)%W_V(iv)%freq_opp = W_F(idf)%deltat
1234  ENDIF
1235!-
1236  IF (test_fwrt < W_F(idf)%deltat) THEN
1237    str70 = 'Frequency of output should be larger than deltat'
1238    WRITE(str71,'("It is not the case for variable ",a," :",F10.4)') &
1239 &   TRIM(tmp_name),pfreq_wrt
1240    str72 = "PATCH : frequency set to deltat"
1241!-
1242    CALL ipslerr (2,"histdef",str70,str71,str72)
1243!-
1244    W_F(idf)%W_V(iv)%freq_wrt = W_F(idf)%deltat
1245  ENDIF
1246!-
1247! 4.2 First the existence of the operation is tested and then
1248!     its compaticility with the choice of frequencies
1249!-
1250  IF (TRIM(W_F(idf)%W_V(iv)%topp) == "inst") THEN
1251    IF (test_fopp /= test_fwrt) THEN
1252      str70 = 'For instantaneous output the frequency '// &
1253 &            'of operations and output'
1254      WRITE(str71, &
1255 &     '("should be the same, this was not case for variable ",a)') &
1256 &      TRIM(tmp_name)
1257      str72 = "PATCH : The smalest frequency of both is used"
1258      CALL ipslerr (2,"histdef",str70,str71,str72)
1259      IF (test_fopp < test_fwrt) THEN
1260        W_F(idf)%W_V(iv)%freq_opp = pfreq_opp
1261        W_F(idf)%W_V(iv)%freq_wrt = pfreq_opp
1262      ELSE
1263        W_F(idf)%W_V(iv)%freq_opp = pfreq_wrt
1264        W_F(idf)%W_V(iv)%freq_wrt = pfreq_wrt
1265      ENDIF
1266    ENDIF
1267  ELSE IF (INDEX(ex_topps,TRIM(W_F(idf)%W_V(iv)%topp)) > 0) THEN
1268    IF (test_fopp > test_fwrt) THEN
1269      str70 = 'For averages the frequency of operations '// &
1270 &            'should be smaller or equal'
1271      WRITE(str71, &
1272 &     '("to that of output. It is not the case for variable ",a)') &
1273 &     TRIM(tmp_name)
1274      str72 = 'PATCH : The output frequency is used for both'
1275      CALL ipslerr (2,"histdef",str70,str71,str72)
1276      W_F(idf)%W_V(iv)%freq_opp = pfreq_wrt
1277    ENDIF
1278  ELSE
1279    WRITE (str70,'("Operation on variable ",A," is unknown")') &
1280 &   TRIM(tmp_name)
1281    WRITE (str71,'("operation requested is :",A)') &
1282 &   W_F(idf)%W_V(iv)%topp
1283    WRITE (str72,'("File ID :",I3)') idf
1284    CALL ipslerr (3,"histdef",str70,str71,str72)
1285  ENDIF
1286!-
1287! 5.0 Initialize other variables of the common
1288!-
1289  IF (l_dbg) WRITE(*,*) "histdef : 5.0"
1290!-
1291  W_F(idf)%W_V(iv)%hist_wrt_rng = (PRESENT(var_range))
1292  IF (W_F(idf)%W_V(iv)%hist_wrt_rng) THEN
1293    W_F(idf)%W_V(iv)%hist_calc_rng = (var_range(1) > var_range(2))
1294    IF (W_F(idf)%W_V(iv)%hist_calc_rng) THEN
1295      W_F(idf)%W_V(iv)%hist_minmax(1:2) = &
1296 &      (/ ABS(missing_val),-ABS(missing_val) /)
1297    ELSE
1298      W_F(idf)%W_V(iv)%hist_minmax(1:2) = var_range(1:2)
1299    ENDIF
1300  ENDIF
1301!-
1302! - freq_opp(idf,iv)/2./deltat(idf)
1303  W_F(idf)%W_V(iv)%last_opp = W_F(idf)%itau0
1304! - freq_wrt(idf,iv)/2./deltat(idf)
1305  W_F(idf)%W_V(iv)%last_wrt = W_F(idf)%itau0
1306! - freq_opp(idf,iv)/2./deltat(idf)
1307  W_F(idf)%W_V(iv)%last_opp_chk = W_F(idf)%itau0
1308! - freq_wrt(idf,iv)/2./deltat(idf)
1309  W_F(idf)%W_V(iv)%last_wrt_chk = W_F(idf)%itau0
1310  W_F(idf)%W_V(iv)%nb_opp = 0
1311  W_F(idf)%W_V(iv)%nb_wrt = 0
1312!-
1313! 6.0 Get the time axis for this variable
1314!-
1315  IF (l_dbg) WRITE(*,*) "histdef : 6.0"
1316!-
1317! No time axis for once, l_max, l_min or never operation
1318!-
1319  IF (     (TRIM(W_F(idf)%W_V(iv)%topp) /= 'once')  &
1320 &    .AND.(TRIM(W_F(idf)%W_V(iv)%topp) /= 'never') &
1321 &    .AND.(TRIM(W_F(idf)%W_V(iv)%topp) /= 'l_max') &
1322 &    .AND.(TRIM(W_F(idf)%W_V(iv)%topp) /= 'l_min') ) THEN
1323    IF (TRIM(W_F(idf)%W_V(iv)%topp) == 'inst') THEN
1324      str10 = 't_inst_'
1325    ELSE
1326      str10 = 't_op_'
1327    ENDIF
1328    IF (W_F(idf)%W_V(iv)%freq_wrt > 0) THEN
1329      WRITE (UNIT=str40,FMT='(A,I8.8)') &
1330&      TRIM(str10),INT(W_F(idf)%W_V(iv)%freq_wrt)
1331    ELSE
1332      WRITE (UNIT=str40,FMT='(A,I2.2,"month")') &
1333&      TRIM(str10),ABS(INT(W_F(idf)%W_V(iv)%freq_wrt))
1334    ENDIF
1335    CALL find_str (W_F(idf)%W_V(1:W_F(idf)%n_tax)%tax_name,str40,pos)
1336    IF (pos < 0) THEN
1337      W_F(idf)%n_tax = W_F(idf)%n_tax+1
1338      W_F(idf)%W_V(iv)%l_bnd = &
1339 &      (TRIM(W_F(idf)%W_V(iv)%topp) /= 'inst')
1340      W_F(idf)%W_V(W_F(idf)%n_tax)%tax_name = str40
1341      W_F(idf)%W_V(W_F(idf)%n_tax)%tax_last = 0
1342      W_F(idf)%W_V(iv)%t_axid = W_F(idf)%n_tax
1343    ELSE
1344      W_F(idf)%W_V(iv)%t_axid = pos
1345    ENDIF
1346  ELSE
1347    IF (l_dbg) THEN
1348      WRITE(*,*) "histdef : 7.0 ",TRIM(W_F(idf)%W_V(iv)%topp),'----'
1349    ENDIF
1350    W_F(idf)%W_V(iv)%t_axid = -99
1351  ENDIF
1352!-
1353! 7.0 prepare frequence of writing and operation
1354!     for never or once operation
1355!-
1356  IF (    (TRIM(W_F(idf)%W_V(iv)%topp) == 'once')  &
1357 &    .OR.(TRIM(W_F(idf)%W_V(iv)%topp) == 'never') ) THEN
1358    W_F(idf)%W_V(iv)%freq_opp = 0.
1359    W_F(idf)%W_V(iv)%freq_wrt = 0.
1360  ENDIF
1361!---------------------
1362END SUBROUTINE histdef
1363!===
1364SUBROUTINE histend (idf, snc4chunks)
1365!---------------------------------------------------------------------
1366!- This subroutine end the decalaration of variables and sets the
1367!- time axes in the netcdf file and puts it into the write mode.
1368!-
1369!- INPUT
1370!-
1371!- idf : ID of the file to be worked on
1372!-
1373!- VERSION
1374!-
1375!---------------------------------------------------------------------
1376  IMPLICIT NONE
1377!-
1378  INTEGER,INTENT(IN) :: idf
1379  TYPE(snc4_ctl), OPTIONAL,INTENT(IN) :: snc4chunks
1380!-
1381  INTEGER :: nfid,nvid,iret,ndim,iv,itx,ziv,itax,dim_cnt
1382  INTEGER,DIMENSION(4) :: dims
1383  INTEGER :: year,month,day,hours,minutes
1384  REAL :: sec
1385  REAL :: rtime0
1386  CHARACTER(LEN=30) :: str30
1387  CHARACTER(LEN=35) :: str35
1388  CHARACTER(LEN=120) :: assoc
1389  CHARACTER(LEN=70) :: str70
1390  CHARACTER(LEN=3),DIMENSION(12) :: cal =   &
1391 &  (/ 'JAN','FEB','MAR','APR','MAY','JUN', &
1392 &     'JUL','AUG','SEP','OCT','NOV','DEC' /)
1393  CHARACTER(LEN=7) :: tmp_opp
1394  LOGICAL :: l_b
1395  LOGICAL :: l_dbg
1396  INTEGER, DIMENSION(4) :: ichunksz    ! NETCDF4 chunk sizes
1397  INTEGER :: ichunkalg, ishuffle,&
1398             ideflate, ideflate_level
1399  LOGICAL :: lchunk = .FALSE.   ! logical switch to activate chunking when appropriate
1400!-
1401  ! NetCDF4 chunking and compression parameters
1402  ichunkalg = 0
1403  ishuffle = 1
1404  ideflate = 1
1405  ideflate_level = 1
1406  !
1407!---------------------------------------------------------------------
1408  CALL ipsldbg (old_status=l_dbg)
1409!-
1410  nfid = W_F(idf)%ncfid
1411!-
1412! 1.0 Create the time axes
1413!-
1414  IF (l_dbg) WRITE(*,*) "histend : 1.0"
1415!-
1416! 1.1 Define the time dimensions needed for this file
1417!-
1418  iret = NF90_DEF_DIM (nfid,'time_counter', &
1419 &                     NF90_UNLIMITED,W_F(idf)%tid)
1420  DO iv=1,W_F(idf)%n_var
1421    IF (W_F(idf)%W_V(iv)%l_bnd) THEN
1422      iret = NF90_DEF_DIM (nfid,'tbnds',2,W_F(idf)%bid)
1423      EXIT
1424    ENDIF
1425  ENDDO
1426!-
1427! 1.2 Define all the time axes needed for this file
1428!-
1429  DO itx=1,W_F(idf)%n_tax
1430    dims(1) = W_F(idf)%tid
1431    l_b = (INDEX(W_F(idf)%W_V(itx)%tax_name,"t_op_") == 1)
1432    IF (itx > 1) THEN
1433      str30 = W_F(idf)%W_V(itx)%tax_name
1434    ELSE
1435      str30 = "time_counter"
1436    ENDIF
1437    IF (l_b) THEN
1438      str35 = TRIM(str30)//'_bnds'
1439    ENDIF
1440    iret = NF90_DEF_VAR (nfid,TRIM(str30),NF90_REAL8, &
1441 &          dims(1),W_F(idf)%W_V(itx)%tdimid)
1442    IF (itx <= 1) THEN
1443      iret = NF90_PUT_ATT (nfid,W_F(idf)%W_V(itx)%tdimid,'axis',"T")
1444    ENDIF
1445    iret = NF90_PUT_ATT (nfid,W_F(idf)%W_V(itx)%tdimid, &
1446 &          'standard_name',"time")
1447!---
1448!   To transform the current itau into a real date and take it
1449!   as the origin of the file requires the time counter to change.
1450!   Thus it is an operation the user has to ask for.
1451!   This function should thus only be re-instated
1452!   if there is a ioconf routine to control it.
1453!---
1454!-- rtime0 = itau2date(itau0(idf),date0(idf),deltat(idf))
1455    rtime0 = W_F(idf)%date0
1456!-
1457    CALL ju2ymds(rtime0,year,month,day,sec)
1458!---
1459!   Catch any error induced by a change in calendar !
1460!---
1461    IF (year < 0) THEN
1462      year = 2000+year
1463    ENDIF
1464!-
1465    hours = INT(sec/(60.*60.))
1466    minutes = INT((sec-hours*60.*60.)/60.)
1467    sec = sec-(hours*60.*60.+minutes*60.)
1468!-
1469    WRITE (UNIT=str70, &
1470 &   FMT='(A,I4.4,"-",I2.2,"-",I2.2," ",I2.2,":",I2.2,":",I2.2)') &
1471 &    'seconds since ',year,month,day,hours,minutes,INT(sec)
1472    iret = NF90_PUT_ATT (nfid,W_F(idf)%W_V(itx)%tdimid, &
1473 &           'units',TRIM(str70))
1474!-
1475    CALL ioget_calendar (str30)
1476    iret = NF90_PUT_ATT (nfid,W_F(idf)%W_V(itx)%tdimid, &
1477 &           'calendar',TRIM(str30))
1478!-
1479    iret = NF90_PUT_ATT (nfid,W_F(idf)%W_V(itx)%tdimid, &
1480 &           'title','Time')
1481!-
1482    iret = NF90_PUT_ATT (nfid,W_F(idf)%W_V(itx)%tdimid, &
1483 &           'long_name','Time axis')
1484!-
1485    WRITE (UNIT=str70, &
1486 &   FMT='(" ",I4.4,"-",A3,"-",I2.2," ",I2.2,":",I2.2,":",I2.2)') &
1487 &    year,cal(month),day,hours,minutes,INT(sec)
1488    iret = NF90_PUT_ATT (nfid,W_F(idf)%W_V(itx)%tdimid, &
1489 &           'time_origin',TRIM(str70))
1490!---
1491    IF (l_b) THEN
1492      iret = NF90_PUT_ATT (nfid,W_F(idf)%W_V(itx)%tdimid, &
1493 &             'bounds',TRIM(str35))
1494      dims(1:2) = (/ W_F(idf)%bid,W_F(idf)%tid /)
1495      iret = NF90_DEF_VAR (nfid,TRIM(str35),NF90_REAL8, &
1496 &            dims(1:2),W_F(idf)%W_V(itx)%tbndid)
1497    ENDIF
1498  ENDDO
1499!-
1500! 2.0 declare the variables
1501!-
1502  IF (l_dbg) WRITE(*,*) "histend : 2.0"
1503!-
1504  DO iv=1,W_F(idf)%n_var
1505!---
1506    itax = W_F(idf)%W_V(iv)%t_axid
1507!---
1508    IF (W_F(idf)%regular) THEN
1509      dims(1:2) = (/ W_F(idf)%xid,W_F(idf)%yid /)
1510      dim_cnt = 2
1511    ELSE
1512      dims(1) = W_F(idf)%xid
1513      dim_cnt = 1
1514    ENDIF
1515!---
1516    tmp_opp = W_F(idf)%W_V(iv)%topp
1517    ziv = W_F(idf)%W_V(iv)%z_axid
1518!---
1519!   2.1 dimension of field
1520!---
1521    IF ((TRIM(tmp_opp) /= 'never')) THEN
1522      IF (     (TRIM(tmp_opp) /= 'once')  &
1523     &    .AND.(TRIM(tmp_opp) /= 'l_max') &
1524     &    .AND.(TRIM(tmp_opp) /= 'l_min') ) THEN
1525        IF (ziv == -99) THEN
1526          ndim = dim_cnt+1
1527          dims(dim_cnt+1:dim_cnt+2) = (/ W_F(idf)%tid,0 /)
1528        ELSE
1529          ndim = dim_cnt+2
1530          dims(dim_cnt+1:dim_cnt+2) = &
1531 &         (/ W_F(idf)%zax_ids(ziv),W_F(idf)%tid /)
1532        ENDIF
1533      ELSE
1534        IF (ziv == -99) THEN
1535          ndim = dim_cnt
1536          dims(dim_cnt+1:dim_cnt+2) = (/ 0,0 /)
1537        ELSE
1538          ndim = dim_cnt+1
1539          dims(dim_cnt+1:dim_cnt+2) = (/ W_F(idf)%zax_ids(ziv),0 /)
1540        ENDIF
1541      ENDIF
1542!-
1543      iret = NF90_DEF_VAR (nfid,TRIM(W_F(idf)%W_V(iv)%v_name), &
1544 &             W_F(idf)%W_V(iv)%v_typ,dims(1:ABS(ndim)),nvid)
1545!-
1546      IF( ndim == 4 ) THEN
1547        IF( PRESENT( snc4chunks ) ) THEN
1548          IF( snc4chunks%luse ) THEN
1549           ichunksz = 1
1550           iret = NF90_INQUIRE_DIMENSION( nfid, W_F(idf)%xid, len = ichunksz(1) )
1551           iret = NF90_INQUIRE_DIMENSION( nfid, W_F(idf)%yid, len = ichunksz(2) )
1552           IF ( ziv .NE. -99 ) &
1553              iret = NF90_INQUIRE_DIMENSION( nfid, W_F(idf)%zax_ids(ziv), len = ichunksz(3) )
1554           ichunksz(1) = MIN(ichunksz(1), MAX((ichunksz(1)-1)/snc4chunks%ni + 1,16))
1555           ichunksz(2) = MIN(ichunksz(2), MAX((ichunksz(2)-1)/snc4chunks%nj + 1,16))
1556           ichunksz(3) = MIN(ichunksz(3), MAX((ichunksz(3)-1)/snc4chunks%nk + 1, 1))
1557           ! Always use a chunk size of 1 for the unlimited dimension
1558           iret = SET_NF90_DEF_VAR_CHUNKING(nfid, nvid, ichunkalg, ichunksz)
1559           iret = SET_NF90_DEF_VAR_DEFLATE(nfid, nvid, ishuffle, ideflate, ideflate_level)
1560          ENDIF
1561        ENDIF
1562      ENDIF
1563      W_F(idf)%W_V(iv)%ncvid = nvid
1564!-
1565      IF (LEN_TRIM(W_F(idf)%W_V(iv)%unit_name) > 0) THEN
1566        iret = NF90_PUT_ATT (nfid,nvid,'units', &
1567 &                           TRIM(W_F(idf)%W_V(iv)%unit_name))
1568      ENDIF
1569      iret = NF90_PUT_ATT (nfid,nvid,'standard_name', &
1570 &                         TRIM(W_F(idf)%W_V(iv)%std_name))
1571!-
1572      IF (W_F(idf)%W_V(iv)%v_typ == hist_r8) THEN
1573        iret = NF90_PUT_ATT (nfid,nvid,'_FillValue',NF90_FILL_REAL8)
1574      ELSE
1575        iret = NF90_PUT_ATT (nfid,nvid,'_FillValue',NF90_FILL_REAL4)
1576      ENDIF
1577      IF (W_F(idf)%W_V(iv)%hist_wrt_rng) THEN
1578        IF (W_F(idf)%W_V(iv)%v_typ == hist_r8) THEN
1579          iret = NF90_PUT_ATT (nfid,nvid,'valid_min', &
1580 &                 REAL(W_F(idf)%W_V(iv)%hist_minmax(1),KIND=8))
1581          iret = NF90_PUT_ATT (nfid,nvid,'valid_max', &
1582 &                 REAL(W_F(idf)%W_V(iv)%hist_minmax(2),KIND=8))
1583        ELSE
1584          iret = NF90_PUT_ATT (nfid,nvid,'valid_min', &
1585 &                 REAL(W_F(idf)%W_V(iv)%hist_minmax(1),KIND=4))
1586          iret = NF90_PUT_ATT (nfid,nvid,'valid_max', &
1587 &                 REAL(W_F(idf)%W_V(iv)%hist_minmax(2),KIND=4))
1588        ENDIF
1589      ENDIF
1590      iret = NF90_PUT_ATT (nfid,nvid,'long_name', &
1591 &                         TRIM(W_F(idf)%W_V(iv)%title))
1592      iret = NF90_PUT_ATT (nfid,nvid,'online_operation', &
1593 &                         TRIM(W_F(idf)%W_V(iv)%fullop))
1594!-
1595      SELECT CASE(ndim)
1596      CASE(-3,2:4)
1597      CASE DEFAULT
1598        CALL ipslerr (3,"histend", &
1599       &  'less than 2 or more than 4 dimensions are not', &
1600       &  'allowed at this stage',' ')
1601      END SELECT
1602!-
1603      assoc=TRIM(W_F(idf)%hax_name(W_F(idf)%W_V(iv)%h_axid,2)) &
1604 &   //' '//TRIM(W_F(idf)%hax_name(W_F(idf)%W_V(iv)%h_axid,1))
1605!-
1606      ziv = W_F(idf)%W_V(iv)%z_axid
1607      IF (ziv > 0) THEN
1608        str30 = W_F(idf)%zax_name(ziv)
1609        assoc = TRIM(str30)//' '//TRIM(assoc)
1610      ENDIF
1611!-
1612      IF (itax > 0) THEN
1613        IF (itax > 1) THEN
1614          str30 = W_F(idf)%W_V(itax)%tax_name
1615        ELSE
1616          str30 = "time_counter"
1617        ENDIF
1618        assoc = TRIM(str30)//' '//TRIM(assoc)
1619!-
1620        IF (l_dbg) THEN
1621          WRITE(*,*) "histend : 2.0.n, freq_opp, freq_wrt", &
1622 &          W_F(idf)%W_V(iv)%freq_opp,W_F(idf)%W_V(iv)%freq_wrt
1623        ENDIF
1624!-
1625        iret = NF90_PUT_ATT (nfid,nvid,'interval_operation', &
1626 &                           REAL(W_F(idf)%W_V(iv)%freq_opp,KIND=4))
1627        iret = NF90_PUT_ATT (nfid,nvid,'interval_write', &
1628 &                           REAL(W_F(idf)%W_V(iv)%freq_wrt,KIND=4))
1629      ENDIF
1630      iret = NF90_PUT_ATT (nfid,nvid,'coordinates',TRIM(assoc))
1631    ENDIF
1632  ENDDO
1633!-
1634! 2.2 Add DOMAIN attributes if needed
1635!-
1636  IF (W_F(idf)%dom_id_svg >= 0) THEN
1637    CALL flio_dom_att (nfid,W_F(idf)%dom_id_svg)
1638  ENDIF
1639!-
1640! 3.0 Put the netcdf file into write mode
1641!-
1642  IF (l_dbg) WRITE(*,*) "histend : 3.0"
1643!-
1644  iret = NF90_ENDDEF (nfid)
1645!-
1646! 4.0 Give some informations to the user
1647!-
1648  IF (l_dbg) WRITE(*,*) "histend : 4.0"
1649!-
1650!!$  WRITE(str70,'("All variables have been initialized on file :",I3)') idf
1651!!$  CALL ipslerr (1,'histend',str70,'',' ')
1652!---------------------
1653END SUBROUTINE histend
1654!===
1655SUBROUTINE histwrite_r1d (idf,pvarname,pitau,pdata,nbindex,nindex)
1656!---------------------------------------------------------------------
1657  IMPLICIT NONE
1658!-
1659  INTEGER,INTENT(IN) :: idf,pitau,nbindex
1660  INTEGER,DIMENSION(nbindex),INTENT(IN) :: nindex
1661  REAL,DIMENSION(:),INTENT(IN) :: pdata
1662  CHARACTER(LEN=*),INTENT(IN) :: pvarname
1663!---------------------------------------------------------------------
1664  CALL histw_rnd (idf,pvarname,pitau,nbindex,nindex,pdata_1d=pdata)
1665!---------------------------
1666END SUBROUTINE histwrite_r1d
1667!===
1668SUBROUTINE histwrite_r2d (idf,pvarname,pitau,pdata,nbindex,nindex)
1669!---------------------------------------------------------------------
1670  IMPLICIT NONE
1671!-
1672  INTEGER,INTENT(IN) :: idf,pitau,nbindex
1673  INTEGER,DIMENSION(nbindex),INTENT(IN) :: nindex
1674  REAL,DIMENSION(:,:),INTENT(IN) :: pdata
1675  CHARACTER(LEN=*),INTENT(IN) :: pvarname
1676!---------------------------------------------------------------------
1677  CALL histw_rnd (idf,pvarname,pitau,nbindex,nindex,pdata_2d=pdata)
1678!---------------------------
1679END SUBROUTINE histwrite_r2d
1680!===
1681SUBROUTINE histwrite_r3d (idf,pvarname,pitau,pdata,nbindex,nindex)
1682!---------------------------------------------------------------------
1683  IMPLICIT NONE
1684!-
1685  INTEGER,INTENT(IN) :: idf,pitau,nbindex
1686  INTEGER,DIMENSION(nbindex),INTENT(IN) :: nindex
1687  REAL,DIMENSION(:,:,:),INTENT(IN) :: pdata
1688  CHARACTER(LEN=*),INTENT(IN) :: pvarname
1689!---------------------------------------------------------------------
1690  CALL histw_rnd (idf,pvarname,pitau,nbindex,nindex,pdata_3d=pdata)
1691!---------------------------
1692END SUBROUTINE histwrite_r3d
1693!===
1694SUBROUTINE histw_rnd (idf,pvarname,pitau,nbindex,nindex, &
1695  &                   pdata_1d,pdata_2d,pdata_3d)
1696!---------------------------------------------------------------------
1697  IMPLICIT NONE
1698!-
1699  INTEGER,INTENT(IN) :: idf,pitau,nbindex
1700  INTEGER,DIMENSION(nbindex),INTENT(IN) :: nindex
1701  CHARACTER(LEN=*),INTENT(IN) :: pvarname
1702  REAL,DIMENSION(:),INTENT(IN),OPTIONAL     :: pdata_1d
1703  REAL,DIMENSION(:,:),INTENT(IN),OPTIONAL   :: pdata_2d
1704  REAL,DIMENSION(:,:,:),INTENT(IN),OPTIONAL :: pdata_3d
1705!-
1706  LOGICAL :: do_oper,do_write,largebuf,l1d,l2d,l3d
1707  INTEGER :: iv,io,nbpt_out
1708  INTEGER              :: nbpt_in1
1709  INTEGER,DIMENSION(2) :: nbpt_in2
1710  INTEGER,DIMENSION(3) :: nbpt_in3
1711  REAL,ALLOCATABLE,DIMENSION(:),SAVE :: tbf_1
1712  CHARACTER(LEN=7) :: tmp_opp
1713  CHARACTER(LEN=13) :: c_nam
1714  LOGICAL :: l_dbg
1715!---------------------------------------------------------------------
1716  CALL ipsldbg (old_status=l_dbg)
1717!-
1718  l1d=PRESENT(pdata_1d); l2d=PRESENT(pdata_2d); l3d=PRESENT(pdata_3d);
1719  IF      (l1d) THEN
1720    c_nam = 'histwrite_r1d'
1721  ELSE IF (l2d) THEN
1722    c_nam = 'histwrite_r2d'
1723  ELSE IF (l3d) THEN
1724    c_nam = 'histwrite_r3d'
1725  ENDIF
1726!-
1727  IF (l_dbg) THEN
1728    WRITE(*,*) "histwrite : ",c_nam
1729  ENDIF
1730!-
1731! 1.0 Try to catch errors like specifying the wrong file ID.
1732!     Thanks Marine for showing us what errors users can make !
1733!-
1734  IF ( (idf < 1).OR.(idf > nb_files_max) ) THEN
1735    CALL ipslerr (3,"histwrite", &
1736 &    'Illegal file ID in the histwrite of variable',pvarname,' ')
1737  ENDIF
1738!-
1739! 1.1 Find the id of the variable to be written and the real time
1740!-
1741  CALL histvar_seq (idf,pvarname,iv)
1742!-
1743! 2.0 do nothing for never operation
1744!-
1745  tmp_opp = W_F(idf)%W_V(iv)%topp
1746!-
1747  IF (TRIM(tmp_opp) == "never") THEN
1748    W_F(idf)%W_V(iv)%last_opp_chk = -99
1749    W_F(idf)%W_V(iv)%last_wrt_chk = -99
1750  ENDIF
1751!-
1752! 3.0 We check if we need to do an operation
1753!-
1754  IF (W_F(idf)%W_V(iv)%last_opp_chk == pitau) THEN
1755    CALL ipslerr (3,"histwrite", &
1756 &    'This variable has already been analysed at the present', &
1757 &    'time step',TRIM(pvarname))
1758  ENDIF
1759!-
1760  CALL isittime &
1761 &  (pitau,W_F(idf)%date0,W_F(idf)%deltat, &
1762 &   W_F(idf)%W_V(iv)%freq_opp, &
1763 &   W_F(idf)%W_V(iv)%last_opp, &
1764 &   W_F(idf)%W_V(iv)%last_opp_chk,do_oper)
1765!-
1766! 4.0 We check if we need to write the data
1767!-
1768  IF (W_F(idf)%W_V(iv)%last_wrt_chk == pitau) THEN
1769    CALL ipslerr (3,"histwrite", &
1770 &    'This variable as already been written for the present', &
1771 &    'time step',' ')
1772  ENDIF
1773!-
1774  CALL isittime &
1775 &  (pitau,W_F(idf)%date0,W_F(idf)%deltat, &
1776 &   W_F(idf)%W_V(iv)%freq_wrt, &
1777 &   W_F(idf)%W_V(iv)%last_wrt, &
1778 &   W_F(idf)%W_V(iv)%last_wrt_chk,do_write)
1779!-
1780! 5.0 histwrite called
1781!-
1782  IF (do_oper.OR.do_write) THEN
1783!-
1784!-- 5.1 Get the sizes of the data we will handle
1785!-
1786    IF (W_F(idf)%W_V(iv)%datasz_in(1) <= 0) THEN
1787!---- There is the risk here that the user has over-sized the array.
1788!---- But how can we catch this ?
1789!---- In the worst case we will do impossible operations
1790!---- on part of the data !
1791      W_F(idf)%W_V(iv)%datasz_in(1:3) = -1
1792      IF      (l1d) THEN
1793        W_F(idf)%W_V(iv)%datasz_in(1) = SIZE(pdata_1d)
1794      ELSE IF (l2d) THEN
1795        W_F(idf)%W_V(iv)%datasz_in(1) = SIZE(pdata_2d,DIM=1)
1796        W_F(idf)%W_V(iv)%datasz_in(2) = SIZE(pdata_2d,DIM=2)
1797      ELSE IF (l3d) THEN
1798        W_F(idf)%W_V(iv)%datasz_in(1) = SIZE(pdata_3d,DIM=1)
1799        W_F(idf)%W_V(iv)%datasz_in(2) = SIZE(pdata_3d,DIM=2)
1800        W_F(idf)%W_V(iv)%datasz_in(3) = SIZE(pdata_3d,DIM=3)
1801      ENDIF
1802    ENDIF
1803!-
1804!-- 5.2 The maximum size of the data will give the size of the buffer
1805!-
1806    IF (W_F(idf)%W_V(iv)%datasz_max <= 0) THEN
1807      largebuf = .FALSE.
1808      DO io=1,W_F(idf)%W_V(iv)%nbopp
1809        IF (INDEX(fuchnbout,W_F(idf)%W_V(iv)%sopp(io)) > 0) THEN
1810          largebuf = .TRUE.
1811        ENDIF
1812      ENDDO
1813      IF (largebuf) THEN
1814        W_F(idf)%W_V(iv)%datasz_max = &
1815 &        W_F(idf)%W_V(iv)%scsize(1) &
1816 &       *W_F(idf)%W_V(iv)%scsize(2) &
1817 &       *W_F(idf)%W_V(iv)%scsize(3)
1818      ELSE
1819        IF      (l1d) THEN
1820          W_F(idf)%W_V(iv)%datasz_max = &
1821 &          W_F(idf)%W_V(iv)%datasz_in(1)
1822        ELSE IF (l2d) THEN
1823          W_F(idf)%W_V(iv)%datasz_max = &
1824 &          W_F(idf)%W_V(iv)%datasz_in(1) &
1825 &         *W_F(idf)%W_V(iv)%datasz_in(2)
1826        ELSE IF (l3d) THEN
1827          W_F(idf)%W_V(iv)%datasz_max = &
1828 &          W_F(idf)%W_V(iv)%datasz_in(1) &
1829 &         *W_F(idf)%W_V(iv)%datasz_in(2) &
1830 &         *W_F(idf)%W_V(iv)%datasz_in(3)
1831        ENDIF
1832      ENDIF
1833    ENDIF
1834!-
1835    IF (.NOT.ALLOCATED(tbf_1)) THEN
1836      IF (l_dbg) THEN
1837        WRITE(*,*) &
1838 &       c_nam//" : allocate tbf_1 for size = ", &
1839 &       W_F(idf)%W_V(iv)%datasz_max
1840      ENDIF
1841      ALLOCATE(tbf_1(W_F(idf)%W_V(iv)%datasz_max))
1842    ELSE IF (W_F(idf)%W_V(iv)%datasz_max > SIZE(tbf_1)) THEN
1843      IF (l_dbg) THEN
1844        WRITE(*,*) &
1845 &       c_nam//" : re-allocate tbf_1 for size = ", &
1846 &       W_F(idf)%W_V(iv)%datasz_max
1847      ENDIF
1848      DEALLOCATE(tbf_1)
1849      ALLOCATE(tbf_1(W_F(idf)%W_V(iv)%datasz_max))
1850    ENDIF
1851!-
1852!-- We have to do the first operation anyway.
1853!-- Thus we do it here and change the ranke
1854!-- of the data at the same time. This should speed up things.
1855!-
1856    nbpt_out = W_F(idf)%W_V(iv)%datasz_max
1857    IF      (l1d) THEN
1858      nbpt_in1 = W_F(idf)%W_V(iv)%datasz_in(1)
1859      CALL mathop (W_F(idf)%W_V(iv)%sopp(1),nbpt_in1,pdata_1d, &
1860 &                 missing_val,nbindex,nindex, &
1861 &                 W_F(idf)%W_V(iv)%scal(1),nbpt_out,tbf_1)
1862    ELSE IF (l2d) THEN
1863      nbpt_in2(1:2) = W_F(idf)%W_V(iv)%datasz_in(1:2)
1864      CALL mathop (W_F(idf)%W_V(iv)%sopp(1),nbpt_in2,pdata_2d, &
1865 &                 missing_val,nbindex,nindex, &
1866 &                 W_F(idf)%W_V(iv)%scal(1),nbpt_out,tbf_1)
1867    ELSE IF (l3d) THEN
1868      nbpt_in3(1:3) = W_F(idf)%W_V(iv)%datasz_in(1:3)
1869      CALL mathop (W_F(idf)%W_V(iv)%sopp(1),nbpt_in3,pdata_3d, &
1870 &                 missing_val,nbindex,nindex, &
1871 &                 W_F(idf)%W_V(iv)%scal(1),nbpt_out,tbf_1)
1872    ENDIF
1873    CALL histwrite_real (idf,iv,pitau,nbpt_out, &
1874 &            tbf_1,nbindex,nindex,do_oper,do_write)
1875  ENDIF
1876!-
1877! 6.0 Manage time steps
1878!-
1879  IF ((TRIM(tmp_opp) /= "once").AND.(TRIM(tmp_opp) /= "never")) THEN
1880    W_F(idf)%W_V(iv)%last_opp_chk = pitau
1881    W_F(idf)%W_V(iv)%last_wrt_chk = pitau
1882  ELSE
1883    W_F(idf)%W_V(iv)%last_opp_chk = -99
1884    W_F(idf)%W_V(iv)%last_wrt_chk = -99
1885  ENDIF
1886!-----------------------
1887END SUBROUTINE histw_rnd
1888!===
1889SUBROUTINE histwrite_real &
1890 & (idf,iv,pitau,nbdpt,tbf_1,nbindex,nindex,do_oper,do_write)
1891!---------------------------------------------------------------------
1892!- This subroutine is internal and does the calculations and writing
1893!- if needed. At a later stage it should be split into an operation
1894!- and writing subroutines.
1895!---------------------------------------------------------------------
1896  IMPLICIT NONE
1897!-
1898  INTEGER,INTENT(IN) :: idf,pitau,iv, &
1899 &                      nbindex,nindex(nbindex),nbdpt
1900  REAL,DIMENSION(:)  :: tbf_1
1901  LOGICAL,INTENT(IN) :: do_oper,do_write
1902!-
1903  INTEGER :: tsz,nfid,nvid,iret,itax,io,nbin,nbout
1904  INTEGER :: nx,ny,nz,ky,kz,kt,kc
1905  INTEGER,DIMENSION(4) :: corner,edges
1906  INTEGER :: itime
1907!-
1908  REAL :: rtime
1909  REAL,DIMENSION(2) :: t_bnd
1910  CHARACTER(LEN=7) :: tmp_opp
1911  REAL,ALLOCATABLE,DIMENSION(:),SAVE :: tbf_2
1912  LOGICAL :: l_dbg
1913!---------------------------------------------------------------------
1914  CALL ipsldbg (old_status=l_dbg)
1915!-
1916  IF (l_dbg) THEN
1917    WRITE(*,*) "histwrite 0.0 : VAR : ",W_F(idf)%W_V(iv)%v_name
1918    WRITE(*,*) "histwrite 0.0 : nbindex :",nbindex
1919    WRITE(*,*) "histwrite 0.0 : nindex  :",nindex(1:MIN(3,nbindex)),'...'
1920  ENDIF
1921!-
1922! The sizes which can be encoutered
1923!-
1924  tsz =  W_F(idf)%W_V(iv)%zsize(1) &
1925 &      *W_F(idf)%W_V(iv)%zsize(2) &
1926 &      *W_F(idf)%W_V(iv)%zsize(3)
1927!-
1928! 1.0 We allocate and the temporary space needed for operations.
1929!     The buffers are only deallocated when more space is needed.
1930!     This reduces the umber of allocates but increases memory needs.
1931!-
1932  IF (.NOT.ALLOCATED(tbf_2)) THEN
1933    IF (l_dbg) THEN
1934      WRITE(*,*) "histwrite_real 1.1 allocate tbf_2 ",SIZE(tbf_1)
1935    ENDIF
1936    ALLOCATE(tbf_2(W_F(idf)%W_V(iv)%datasz_max))
1937  ELSE IF (W_F(idf)%W_V(iv)%datasz_max > SIZE(tbf_2)) THEN
1938    IF (l_dbg) THEN
1939      WRITE(*,*) "histwrite_real 1.2 re-allocate tbf_2 : ", &
1940     & SIZE(tbf_1)," instead of ",SIZE(tbf_2)
1941    ENDIF
1942    DEALLOCATE(tbf_2)
1943    ALLOCATE(tbf_2(W_F(idf)%W_V(iv)%datasz_max))
1944  ENDIF
1945!-
1946  rtime = pitau*W_F(idf)%deltat
1947  tmp_opp = W_F(idf)%W_V(iv)%topp
1948!-
1949! 3.0 Do the operations or transfer the slab of data into tbf_1
1950!-
1951  IF (l_dbg) THEN
1952    WRITE(*,*) "histwrite: 3.0",idf
1953  ENDIF
1954!-
1955! 3.1 DO the Operations only if needed
1956!-
1957  IF (do_oper) THEN
1958    nbout = nbdpt
1959!-
1960!-- 3.4 We continue the sequence of operations
1961!--     we started in the interface routine
1962!-
1963    DO io=2,W_F(idf)%W_V(iv)%nbopp,2
1964      nbin = nbout
1965      nbout = W_F(idf)%W_V(iv)%datasz_max
1966      CALL mathop(W_F(idf)%W_V(iv)%sopp(io),nbin,tbf_1, &
1967 &      missing_val,nbindex,nindex,W_F(idf)%W_V(iv)%scal(io), &
1968 &      nbout,tbf_2)
1969      IF (l_dbg) THEN
1970        WRITE(*,*) &
1971 &       "histwrite: 3.4a nbout : ",nbin,nbout,W_F(idf)%W_V(iv)%sopp(io)
1972      ENDIF
1973!-
1974      nbin = nbout
1975      nbout = W_F(idf)%W_V(iv)%datasz_max
1976      CALL mathop(W_F(idf)%W_V(iv)%sopp(io+1),nbin,tbf_2, &
1977 &      missing_val,nbindex,nindex,W_F(idf)%W_V(iv)%scal(io+1), &
1978 &      nbout,tbf_1)
1979      IF (l_dbg) THEN
1980        WRITE(*,*) &
1981 &     "histwrite: 3.4b nbout : ",nbin,nbout,W_F(idf)%W_V(iv)%sopp(io+1)
1982      ENDIF
1983    ENDDO
1984!-
1985!   3.5 Zoom into the data
1986!-
1987    IF (l_dbg) THEN
1988      WRITE(*,*) &
1989 &     "histwrite: 3.5 size(tbf_1) : ",SIZE(tbf_1)
1990      WRITE(*,*) &
1991 &     "histwrite: 3.5 slab in X :", &
1992 &     W_F(idf)%W_V(iv)%zorig(1),W_F(idf)%W_V(iv)%zsize(1)
1993      WRITE(*,*) &
1994 &     "histwrite: 3.5 slab in Y :", &
1995 &     W_F(idf)%W_V(iv)%zorig(2),W_F(idf)%W_V(iv)%zsize(2)
1996      WRITE(*,*) &
1997 &     "histwrite: 3.5 slab in Z :", &
1998 &     W_F(idf)%W_V(iv)%zorig(3),W_F(idf)%W_V(iv)%zsize(3)
1999      WRITE(*,*) &
2000 &     "histwrite: 3.5 slab of input:", &
2001 &     W_F(idf)%W_V(iv)%scsize(1), &
2002 &     W_F(idf)%W_V(iv)%scsize(2), &
2003 &     W_F(idf)%W_V(iv)%scsize(3)
2004    ENDIF
2005!---
2006!-- We have to consider blocks of contiguous data
2007!---
2008    nx=MAX(W_F(idf)%W_V(iv)%zsize(1),1)
2009    ny=MAX(W_F(idf)%W_V(iv)%zsize(2),1)
2010    nz=MAX(W_F(idf)%W_V(iv)%zsize(3),1)
2011    IF     (     (W_F(idf)%W_V(iv)%zorig(1) == 1) &
2012 &          .AND.(   W_F(idf)%W_V(iv)%zsize(1) &
2013 &                == W_F(idf)%W_V(iv)%scsize(1)) &
2014 &          .AND.(W_F(idf)%W_V(iv)%zorig(2) == 1) &
2015 &          .AND.(   W_F(idf)%W_V(iv)%zsize(2) &
2016 &                == W_F(idf)%W_V(iv)%scsize(2))) THEN
2017      kt = (W_F(idf)%W_V(iv)%zorig(3)-1)*nx*ny
2018      tbf_2(1:nx*ny*nz) = tbf_1(kt+1:kt+nx*ny*nz)
2019    ELSEIF (     (W_F(idf)%W_V(iv)%zorig(1) == 1) &
2020 &          .AND.(   W_F(idf)%W_V(iv)%zsize(1) &
2021 &                == W_F(idf)%W_V(iv)%scsize(1))) THEN
2022      kc = -nx*ny
2023      DO kz=W_F(idf)%W_V(iv)%zorig(3),W_F(idf)%W_V(iv)%zorig(3)+nz-1
2024        kc = kc+nx*ny
2025        kt = ( (kz-1)*W_F(idf)%W_V(iv)%scsize(2) &
2026 &            +W_F(idf)%W_V(iv)%zorig(2)-1)*nx
2027        tbf_2(kc+1:kc+nx*ny) = tbf_1(kt+1:kt+nx*ny)
2028      ENDDO
2029    ELSE
2030      kc = -nx
2031      DO kz=W_F(idf)%W_V(iv)%zorig(3),W_F(idf)%W_V(iv)%zorig(3)+nz-1
2032        DO ky=W_F(idf)%W_V(iv)%zorig(2),W_F(idf)%W_V(iv)%zorig(2)+ny-1
2033          kc = kc+nx
2034          kt = ((kz-1)*W_F(idf)%W_V(iv)%scsize(2)+ky-1) &
2035 &            *W_F(idf)%W_V(iv)%scsize(1) &
2036 &            +W_F(idf)%W_V(iv)%zorig(1)-1
2037          tbf_2(kc+1:kc+nx) = tbf_1(kt+1:kt+nx)
2038        ENDDO
2039      ENDDO
2040    ENDIF
2041!-
2042!-- 4.0 Get the min and max of the field
2043!-
2044    IF (l_dbg) THEN
2045      WRITE(*,*) "histwrite: 4.0 tbf_1",idf,iv, &
2046 &      TRIM(tmp_opp),' ---- ',LEN_TRIM(tmp_opp),nbindex
2047    ENDIF
2048!-
2049    IF (W_F(idf)%W_V(iv)%hist_calc_rng) THEN
2050      W_F(idf)%W_V(iv)%hist_minmax(1) = &
2051 &      MIN(W_F(idf)%W_V(iv)%hist_minmax(1), &
2052 &      MINVAL(tbf_2(1:tsz),MASK=tbf_2(1:tsz) /= missing_val))
2053      W_F(idf)%W_V(iv)%hist_minmax(2) = &
2054 &      MAX(W_F(idf)%W_V(iv)%hist_minmax(2), &
2055 &      MAXVAL(tbf_2(1:tsz),MASK=tbf_2(1:tsz) /= missing_val))
2056    ENDIF
2057!-
2058!-- 5.0 Do the operations if needed. In the case of instantaneous
2059!--     output we do not transfer to the time_buffer.
2060!-
2061    IF (l_dbg) THEN
2062      WRITE(*,*) "histwrite: 5.0 idf : ",idf," iv : ",iv," tsz : ",tsz
2063    ENDIF
2064!-
2065    IF (     (TRIM(tmp_opp) /= "inst") &
2066 &      .AND.(TRIM(tmp_opp) /= "once") ) THEN
2067      CALL moycum(tmp_opp,tsz,W_F(idf)%W_V(iv)%t_bf, &
2068 &           tbf_2,W_F(idf)%W_V(iv)%nb_opp)
2069    ENDIF
2070!-
2071    W_F(idf)%W_V(iv)%last_opp = pitau
2072    W_F(idf)%W_V(iv)%nb_opp = W_F(idf)%W_V(iv)%nb_opp+1
2073!-
2074  ENDIF
2075!-
2076! 6.0 Write to file if needed
2077!-
2078  IF (l_dbg) WRITE(*,*) "histwrite: 6.0",idf
2079!-
2080  IF (do_write) THEN
2081!-
2082    nfid = W_F(idf)%ncfid
2083    nvid = W_F(idf)%W_V(iv)%ncvid
2084!-
2085!-- 6.1 Do the operations that are needed before writting
2086!-
2087    IF (l_dbg) WRITE(*,*) "histwrite: 6.1",idf
2088!-
2089    IF (     (TRIM(tmp_opp) /= "inst") &
2090 &      .AND.(TRIM(tmp_opp) /= "once") ) THEN
2091      t_bnd(1:2) = (/ W_F(idf)%W_V(iv)%last_wrt*W_F(idf)%deltat,rtime /)
2092      rtime = (t_bnd(1)+t_bnd(2))/2.0
2093    ENDIF
2094!-
2095!-- 6.2 Add a value to the time axis of this variable if needed
2096!-
2097    IF (     (TRIM(tmp_opp) /= "l_max") &
2098 &      .AND.(TRIM(tmp_opp) /= "l_min") &
2099 &      .AND.(TRIM(tmp_opp) /= "once") ) THEN
2100!-
2101      IF (l_dbg) WRITE(*,*) "histwrite: 6.2",idf
2102!-
2103      itax  = W_F(idf)%W_V(iv)%t_axid
2104      itime = W_F(idf)%W_V(iv)%nb_wrt+1
2105!-
2106      IF (W_F(idf)%W_V(itax)%tax_last < itime) THEN
2107        iret = NF90_PUT_VAR (nfid,W_F(idf)%W_V(itax)%tdimid, &
2108 &               (/ rtime /),start=(/ itime /),count=(/ 1 /))
2109        IF (W_F(idf)%W_V(itax)%tbndid > 0) THEN
2110          iret = NF90_PUT_VAR (nfid,W_F(idf)%W_V(itax)%tbndid, &
2111 &                 t_bnd,start=(/ 1,itime /),count=(/ 2,1 /))
2112        ENDIF
2113        W_F(idf)%W_V(itax)%tax_last = itime
2114      ENDIF
2115    ELSE
2116      itime=1
2117    ENDIF
2118!-
2119!-- 6.3 Write the data. Only in the case of instantaneous output
2120!       we do not write the buffer.
2121!-
2122    IF (l_dbg) THEN
2123      WRITE(*,*) "histwrite: 6.3",idf,nfid,nvid,iv,itime
2124    ENDIF
2125!-
2126    IF (W_F(idf)%W_V(iv)%scsize(3) == 1) THEN
2127      IF (W_F(idf)%regular) THEN
2128        corner(1:4) = (/ 1,1,itime,0 /)
2129        edges(1:4) = (/ W_F(idf)%W_V(iv)%zsize(1), &
2130 &                      W_F(idf)%W_V(iv)%zsize(2),1,0 /)
2131      ELSE
2132        corner(1:4) = (/ 1,itime,0,0 /)
2133        edges(1:4) = (/ W_F(idf)%W_V(iv)%zsize(1),1,0,0 /)
2134      ENDIF
2135    ELSE
2136      IF (W_F(idf)%regular) THEN
2137        corner(1:4) = (/ 1,1,1,itime /)
2138        edges(1:4) = (/ W_F(idf)%W_V(iv)%zsize(1), &
2139 &                      W_F(idf)%W_V(iv)%zsize(2), &
2140 &                      W_F(idf)%W_V(iv)%zsize(3),1 /)
2141      ELSE
2142        corner(1:4) = (/ 1,1,itime,0 /)
2143        edges(1:4) = (/ W_F(idf)%W_V(iv)%zsize(1), &
2144 &                      W_F(idf)%W_V(iv)%zsize(3),1,0 /)
2145      ENDIF
2146    ENDIF
2147!-
2148    IF (     (TRIM(tmp_opp) /= "inst") &
2149 &      .AND.(TRIM(tmp_opp) /= "once") ) THEN
2150      iret = NF90_PUT_VAR (nfid,nvid,W_F(idf)%W_V(iv)%t_bf, &
2151 &                         start=corner(1:4),count=edges(1:4))
2152    ELSE
2153      iret = NF90_PUT_VAR (nfid,nvid,tbf_2, &
2154 &                         start=corner(1:4),count=edges(1:4))
2155    ENDIF
2156!-
2157    W_F(idf)%W_V(iv)%last_wrt = pitau
2158    W_F(idf)%W_V(iv)%nb_wrt = W_F(idf)%W_V(iv)%nb_wrt+1
2159    W_F(idf)%W_V(iv)%nb_opp = 0
2160!---
2161!   After the write the file can be synchronized so that no data is
2162!   lost in case of a crash. This feature gives up on the benefits of
2163!   buffering and should only be used in debuging mode. A flag is
2164!   needed here to switch to this mode.
2165!---
2166!   iret = NF90_SYNC (nfid)
2167!-
2168  ENDIF
2169!----------------------------
2170END SUBROUTINE histwrite_real
2171!===
2172SUBROUTINE histvar_seq (idf,pvarname,idv)
2173!---------------------------------------------------------------------
2174!- This subroutine optimize the search for the variable in the table.
2175!- In a first phase it will learn the succession of the variables
2176!- called and then it will use the table to guess what comes next.
2177!- It is the best solution to avoid lengthy searches through array
2178!- vectors.
2179!-
2180!- ARGUMENTS :
2181!-
2182!- idf      : id of the file on which we work
2183!- pvarname : The name of the variable we are looking for
2184!- idv      : The var id we found
2185!---------------------------------------------------------------------
2186  IMPLICIT NONE
2187!-
2188  INTEGER,INTENT(in)  :: idf
2189  CHARACTER(LEN=*),INTENT(IN) :: pvarname
2190  INTEGER,INTENT(out) :: idv
2191!-
2192  LOGICAL,SAVE :: learning(nb_files_max)=.TRUE.
2193  INTEGER,SAVE :: overlap(nb_files_max) = -1
2194  INTEGER,SAVE :: varseq(nb_files_max,nb_var_max*3)
2195  INTEGER,SAVE :: varseq_len(nb_files_max) = 0
2196  INTEGER,SAVE :: varseq_pos(nb_files_max)
2197  INTEGER,SAVE :: varseq_err(nb_files_max) = 0
2198  INTEGER      :: ib,sp,nn,pos
2199  CHARACTER(LEN=70) :: str70
2200  LOGICAL :: l_dbg
2201!---------------------------------------------------------------------
2202  CALL ipsldbg (old_status=l_dbg)
2203!-
2204  IF (l_dbg) THEN
2205    WRITE(*,*) 'histvar_seq, start of the subroutine :',learning(idf)
2206  ENDIF
2207!-
2208  IF (learning(idf)) THEN
2209!-
2210!-- 1.0 We compute the length over which we are going
2211!--     to check the overlap
2212!-
2213    IF (overlap(idf) <= 0) THEN
2214      IF (W_F(idf)%n_var > 6) THEN
2215        overlap(idf) = W_F(idf)%n_var/3*2
2216      ELSE
2217        overlap(idf) = W_F(idf)%n_var
2218      ENDIF
2219    ENDIF
2220!-
2221!-- 1.1 Find the position of this string
2222!-
2223    CALL find_str (W_F(idf)%W_V(1:W_F(idf)%n_var)%v_name,pvarname,pos)
2224    IF (pos > 0) THEN
2225      idv = pos
2226    ELSE
2227      CALL ipslerr (3,"histvar_seq", &
2228 &      'The name of the variable you gave has not been declared', &
2229 &      'You should use subroutine histdef for declaring variable', &
2230 &      TRIM(pvarname))
2231    ENDIF
2232!-
2233!-- 1.2 If we have not given up we store the position
2234!--     in the sequence of calls
2235!-
2236    IF (varseq_err(idf) >= 0) THEN
2237      sp = varseq_len(idf)+1
2238      IF (sp <= nb_var_max*3) THEN
2239        varseq(idf,sp) = idv
2240        varseq_len(idf) = sp
2241      ELSE
2242        CALL ipslerr (2,"histvar_seq",&
2243 &       'The learning process has failed and we give up. '// &
2244 &       'Either you sequence is',&
2245 &       'too complex or I am too dumb. '// &
2246 &       'This will only affect the efficiency',&
2247 &       'of your code. Thus if you wish to save time'// &
2248 &       ' contact the IOIPSL team. ')
2249        WRITE(*,*) 'The sequence we have found up to now :'
2250        WRITE(*,*) varseq(idf,1:sp-1)
2251        varseq_err(idf) = -1
2252      ENDIF
2253!-
2254!---- 1.3 Check if we have found the right overlap
2255!-
2256      IF (varseq_len(idf) >= overlap(idf)*2) THEN
2257!-
2258!------ We skip a few variables if needed as they could come
2259!------ from the initialisation of the model.
2260!-
2261        DO ib = 0,sp-overlap(idf)*2
2262          IF ( learning(idf) .AND.&
2263            & SUM(ABS(varseq(idf,ib+1:ib+overlap(idf)) -&
2264            & varseq(idf,sp-overlap(idf)+1:sp))) == 0 ) THEN
2265            learning(idf) = .FALSE.
2266            varseq_len(idf) = sp-overlap(idf)-ib
2267            varseq_pos(idf) = overlap(idf)+ib
2268            varseq(idf,1:varseq_len(idf)) = &
2269 &            varseq(idf,ib+1:ib+varseq_len(idf))
2270          ENDIF
2271        ENDDO
2272      ENDIF
2273    ENDIF
2274  ELSE
2275!-
2276!-- 2.0 Now we know how the calls to histwrite are sequenced
2277!--     and we can get a guess at the var ID
2278!-
2279    nn = varseq_pos(idf)+1
2280    IF (nn > varseq_len(idf)) nn = 1
2281!-
2282    idv = varseq(idf,nn)
2283!-
2284    IF (TRIM(W_F(idf)%W_V(idv)%v_name) /= TRIM(pvarname)) THEN
2285      CALL find_str (W_F(idf)%W_V(1:W_F(idf)%n_var)%v_name,pvarname,pos)
2286      IF (pos > 0) THEN
2287        idv = pos
2288      ELSE
2289        CALL ipslerr (3,"histvar_seq", &
2290 &  'The name of the variable you gave has not been declared',&
2291 &  'You should use subroutine histdef for declaring variable', &
2292 &  TRIM(pvarname))
2293      ENDIF
2294      varseq_err(idf) = varseq_err(idf)+1
2295    ELSE
2296!-
2297!---- We only keep the new position if we have found the variable
2298!---- this way. This way an out of sequence call to histwrite does
2299!---- not defeat the process.
2300!-
2301      varseq_pos(idf) = nn
2302    ENDIF
2303!-
2304!!$    IF (varseq_err(idf) >= 10) THEN
2305!!$      WRITE(str70,'("for file ",I3)') idf
2306!!$      CALL ipslerr (2,"histvar_seq", &
2307!!$ &  'There were 10 errors in the learned sequence of variables',&
2308!!$ &  str70,'This looks like a bug, please report it.')
2309!!$         varseq_err(idf) = 0
2310!!$    ENDIF
2311  ENDIF
2312!-
2313  IF (l_dbg) THEN
2314    WRITE(*,*) &
2315 &   'histvar_seq, end of the subroutine :',TRIM(pvarname),idv
2316  ENDIF
2317!-------------------------
2318END SUBROUTINE histvar_seq
2319!===
2320SUBROUTINE histsync (idf)
2321!---------------------------------------------------------------------
2322!- This subroutine will synchronise all
2323!- (or one if defined) opened files.
2324!-
2325!- VERSION
2326!-
2327!---------------------------------------------------------------------
2328  IMPLICIT NONE
2329!-
2330! idf  : optional argument for fileid
2331  INTEGER,INTENT(in),OPTIONAL :: idf
2332!-
2333  INTEGER :: ifile,iret,i_s,i_e
2334!-
2335  LOGICAL :: l_dbg
2336!---------------------------------------------------------------------
2337  CALL ipsldbg (old_status=l_dbg)
2338!-
2339  IF (l_dbg) THEN
2340    WRITE(*,*) "->histsync"
2341  ENDIF
2342!-
2343  IF (PRESENT(idf)) THEN
2344    IF ( (idf >= 1).AND.(idf <= nb_files_max) ) THEN
2345      IF (W_F(idf)%ncfid > 0) THEN
2346        i_s = idf
2347        i_e = idf
2348      ELSE
2349        i_s = 1
2350        i_e = 0
2351        CALL ipslerr (2,'histsync', &
2352 &       'Unable to synchronise the file :','probably','not opened')
2353      ENDIF
2354    ELSE
2355      CALL ipslerr (3,'histsync','Invalid file identifier',' ',' ')
2356    ENDIF
2357  ELSE
2358    i_s = 1
2359    i_e = nb_files_max
2360  ENDIF
2361!-
2362  DO ifile=i_s,i_e
2363    IF (W_F(ifile)%ncfid > 0) THEN
2364      IF (l_dbg) THEN
2365        WRITE(*,*) '  histsync - synchronising file number ',ifile
2366      ENDIF
2367      iret = NF90_SYNC(W_F(ifile)%ncfid)
2368    ENDIF
2369  ENDDO
2370!-
2371  IF (l_dbg) THEN
2372    WRITE(*,*) "<-histsync"
2373  ENDIF
2374!----------------------
2375END SUBROUTINE histsync
2376!===
2377SUBROUTINE histclo (idf)
2378!---------------------------------------------------------------------
2379!- This subroutine will close all (or one if defined) opened files
2380!-
2381!- VERSION
2382!-
2383!---------------------------------------------------------------------
2384  IMPLICIT NONE
2385!-
2386! idf  : optional argument for fileid
2387  INTEGER,INTENT(in),OPTIONAL :: idf
2388!-
2389  INTEGER :: ifile,nfid,nvid,iret,iv,i_s,i_e
2390  LOGICAL :: l_dbg
2391!---------------------------------------------------------------------
2392  CALL ipsldbg (old_status=l_dbg)
2393!-
2394  IF (l_dbg) THEN
2395    WRITE(*,*) "->histclo"
2396  ENDIF
2397!-
2398  IF (PRESENT(idf)) THEN
2399    IF ( (idf >= 1).AND.(idf <= nb_files_max) ) THEN
2400      IF (W_F(idf)%ncfid > 0) THEN
2401        i_s = idf
2402        i_e = idf
2403      ELSE
2404        i_s = 1
2405        i_e = 0
2406        CALL ipslerr (2,'histclo', &
2407 &       'Unable to close the file :','probably','not opened')
2408      ENDIF
2409    ELSE
2410      CALL ipslerr (3,'histclo','Invalid file identifier',' ',' ')
2411    ENDIF
2412  ELSE
2413    i_s = 1
2414    i_e = nb_files_max
2415  ENDIF
2416!-
2417  DO ifile=i_s,i_e
2418    IF (W_F(ifile)%ncfid > 0) THEN
2419      IF (l_dbg) THEN
2420        WRITE(*,*) '  histclo - closing specified file number :',ifile
2421      ENDIF
2422      nfid = W_F(ifile)%ncfid
2423      iret = NF90_REDEF(nfid)
2424!-----
2425!---- 1. Loop on the number of variables to add some final information
2426!-----
2427      IF (l_dbg) THEN
2428        WRITE(*,*) '  Entering loop on vars : ',W_F(ifile)%n_var
2429      ENDIF
2430      DO iv=1,W_F(ifile)%n_var
2431!------ Extrema
2432        IF (W_F(ifile)%W_V(iv)%hist_wrt_rng) THEN
2433          IF (l_dbg) THEN
2434            WRITE(*,*) 'min value for file :',ifile,' var n. :',iv, &
2435 &                     ' is : ',W_F(ifile)%W_V(iv)%hist_minmax(1)
2436            WRITE(*,*) 'max value for file :',ifile,' var n. :',iv, &
2437 &                     ' is : ',W_F(ifile)%W_V(iv)%hist_minmax(2)
2438          ENDIF
2439          IF (W_F(ifile)%W_V(iv)%hist_calc_rng) THEN
2440!---------- Put the min and max values on the file
2441            nvid = W_F(ifile)%W_V(iv)%ncvid
2442            IF (W_F(ifile)%W_V(iv)%v_typ == hist_r8) THEN
2443              iret = NF90_PUT_ATT(nfid,nvid,'valid_min', &
2444 &                     REAL(W_F(ifile)%W_V(iv)%hist_minmax(1),KIND=8))
2445              iret = NF90_PUT_ATT(nfid,nvid,'valid_max', &
2446 &                     REAL(W_F(ifile)%W_V(iv)%hist_minmax(2),KIND=8))
2447            ELSE
2448              iret = NF90_PUT_ATT(nfid,nvid,'valid_min', &
2449 &                     REAL(W_F(ifile)%W_V(iv)%hist_minmax(1),KIND=4))
2450              iret = NF90_PUT_ATT(nfid,nvid,'valid_max', &
2451 &                     REAL(W_F(ifile)%W_V(iv)%hist_minmax(2),KIND=4))
2452            ENDIF
2453          ENDIF
2454        ENDIF
2455!------ Time-Buffers
2456        IF (ASSOCIATED(W_F(ifile)%W_V(iv)%t_bf)) THEN
2457          DEALLOCATE(W_F(ifile)%W_V(iv)%t_bf)
2458        ENDIF
2459!------ Reinitialize the sizes
2460        W_F(ifile)%W_V(iv)%datasz_in(:) = -1
2461        W_F(ifile)%W_V(iv)%datasz_max = -1
2462      ENDDO
2463!-----
2464!---- 2. Close the file
2465!-----
2466      IF (l_dbg) WRITE(*,*) '  close file :',nfid
2467      iret = NF90_CLOSE(nfid)
2468      W_F(ifile)%ncfid = -1
2469      W_F(ifile)%dom_id_svg = -1
2470    ENDIF
2471  ENDDO
2472!-
2473  IF (l_dbg) THEN
2474    WRITE(*,*) "<-histclo"
2475  ENDIF
2476!---------------------
2477END SUBROUTINE histclo
2478!===
2479SUBROUTINE ioconf_modname (str)
2480!---------------------------------------------------------------------
2481!- This subroutine allows to configure the name
2482!- of the model written into the file
2483!---------------------------------------------------------------------
2484  IMPLICIT NONE
2485!-
2486  CHARACTER(LEN=*),INTENT(IN) :: str
2487!---------------------------------------------------------------------
2488  IF (.NOT.lock_modname) THEN
2489    model_name = str(1:MIN(LEN_TRIM(str),80))
2490    lock_modname = .TRUE.
2491  ELSE
2492    CALL ipslerr (2,"ioconf_modname", &
2493   &  'The model name can only be changed once and only', &
2494   &  'before it is used. It is now set to :',model_name)
2495  ENDIF
2496!----------------------------
2497END SUBROUTINE ioconf_modname
2498!-
2499!===
2500!-
2501!-----------------
2502!$AGRIF_END_DO_NOT_TREAT
2503END MODULE histcom
Note: See TracBrowser for help on using the repository browser.