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 branches/nemo_v3_3_beta/NEMOGCM/EXTERNAL/IOIPSL/src – NEMO

source: branches/nemo_v3_3_beta/NEMOGCM/EXTERNAL/IOIPSL/src/histcom.f90 @ 2364

Last change on this file since 2364 was 2364, checked in by acc, 13 years ago

Added basic NetCDF4 chunking and compression support (key_netcdf4). See ticket #754

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