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.
rebuild_nemo.F90 in NEMO/branches/UKMO/dev_rebuild_nemo_compression/REBUILD_NEMO/src – NEMO

source: NEMO/branches/UKMO/dev_rebuild_nemo_compression/REBUILD_NEMO/src/rebuild_nemo.F90 @ 10220

Last change on this file since 10220 was 10220, checked in by mathiot, 6 years ago

changes suggested by Daley

File size: 62.2 KB
Line 
1PROGRAM rebuild_nemo
2#define key_netcdf4
3   !!=========================================================================
4   !!                        ***  rebuild_nemo  ***
5   !!=========================================================================
6   !!
7   !!  A routine to rebuild NEMO files from multiple processors into one file.
8   !!  This routine is designed to be much quicker than the old IOIPSL rebuild
9   !!  but at the cost of an increased memory usage.
10   !!
11   !!  NEMO rebuild has the following features:
12   !!     * dynamically works out what variables require rebuilding
13   !!     * does not copy subdomain halo regions
14   !!     * works for 1,2,3 and 4d arrays or types for all valid NetCDF types
15   !!     * utilises OMP shared memory parallelisation where applicable
16   !!     * time 'slicing' for lower memory use
17   !!       (only for 4D vars with unlimited dimension)
18   !!
19   !!  Ed Blockley - August 2011
20   !!  (based on original code by Matt Martin)
21   !!  Julien Palmieri and Andrew Coward - September 2018 (add compression and chunking)
22   !!
23   !!-------------------------------------------------------------------------
24   !!
25   !!  The code reads the filestem and number of subdomains from the namelist file nam_rebuild.
26   !!
27   !!  The 1st subdomain file is used to determine the dimensions and variables in all the input files.
28   !!  It is also used to find which dimensions (and hence which variables) require rebuilding
29   !!  as well as information about the global domain.
30   !!
31   !!  It then opens all the input files (unbuffered) and creates an array of netcdf identifiers
32   !!  before looping through all the variables and updating the rebuilt output file (either by direct
33   !!  copying or looping over the number of domains and rebuilding as appropriate).
34   !! 
35   !!  The code looks more complicated than it is because it has lots of case statements to deal with all
36   !!  the various NetCDF data types and with various data dimensions (up to 4d).
37   !!
38   !!  Diagnostic output is written to numout (default 6 - stdout)
39   !!  and errors are written to numerr (default 0 - stderr).
40   !!
41   !!  If time slicing is specified the code will use less memory but take a little longer.
42   !!  It does this by breaking down the 4D input variables over their 4th dimension
43   !!  (generally time) by way of a while loop.
44   !!
45   !!-------------------------------------------------------------------------------
46 
47   USE netcdf
48
49!$ USE omp_lib           ! Note OpenMP sentinel
50
51   IMPLICIT NONE
52
53   ! kind specifications
54   INTEGER,PARAMETER :: i1=SELECTED_INT_KIND(2)          ! NF90_BYTE
55   INTEGER,PARAMETER :: i2=SELECTED_INT_KIND(4)          ! NF90_SHORT
56   INTEGER,PARAMETER :: i4=SELECTED_INT_KIND(9)          ! NF90_INT
57   INTEGER,PARAMETER :: sp=SELECTED_REAL_KIND(6,37)      ! NF90_FLOAT
58   INTEGER,PARAMETER :: dp=SELECTED_REAL_KIND(12,307)    ! NF90_DOUBLE
59
60   INTEGER,PARAMETER :: numnam = 11
61   INTEGER,PARAMETER :: numout = 6
62   INTEGER,PARAMETER :: numerr = 0
63
64   LOGICAL, PARAMETER :: l_verbose = .true. 
65   
66   CHARACTER(LEN=nf90_max_name) :: filebase, suffix, attname, dimname, varname, time, date, zone, timestamp
67   CHARACTER(LEN=nf90_max_name), ALLOCATABLE :: filenames(:), indimnames(:)
68   CHARACTER(LEN=nf90_max_name), DIMENSION(2) :: dims
69   CHARACTER(LEN=256) :: cnampath, cdimlst, cdim
70   CHARACTER(LEN=50)  :: clibnc ! netcdf library version
71
72   INTEGER :: ndomain, ifile, ndomain_file, nslicesize, deflate_level
73   INTEGER :: ncid, outid, idim, istop
74   INTEGER :: natts, attid, xtype, varid, rbdims 
75   INTEGER :: jv, ndims, nvars, dimlen, dimids(4)
76   INTEGER :: dimid, unlimitedDimId, di, dj, dr
77   INTEGER :: nmax_unlimited, nt, ntslice 
78   INTEGER :: fchunksize = 32000000  ! NetCDF global file chunk cache size
79   INTEGER :: patchchunk             ! NetCDF processor-domain file chunk cache size
80   INTEGER :: nthreads = 1
81   INTEGER :: chunkalg = 0          ! NetCDF4 variable chunking algorithm
82                                    ! Default variable chunksizes (typical ORCA025
83                                    ! recommendations which can be adjusted via namelist
84                                    ! or will be bounded if too large for domain.)
85   INTEGER :: nc4_xchunk = 206      ! Default x (longitude) variable chunk size
86   INTEGER :: nc4_ychunk = 135      ! Default y (latitude)  variable chunk size
87   INTEGER :: nc4_zchunk = 1        ! Default z (depth) variable chunk size (almost always 1)
88   INTEGER :: nc4_tchunk = 1        ! Default t (time)  variable chunk size (almost always 1)
89   INTEGER, ALLOCATABLE  :: outdimids(:), outdimlens(:), indimlens(:), inncids(:)
90   INTEGER, ALLOCATABLE  :: chunksizes(:)
91   INTEGER, ALLOCATABLE  :: global_sizes(:), rebuild_dims(:)
92   INTEGER, DIMENSION(2) :: halo_start, halo_end, local_sizes
93   INTEGER, DIMENSION(2) :: idomain, jdomain, rdomain, start_pos
94   INTEGER :: ji, jj, jk, jl, jr
95   INTEGER :: nargs                 ! number of arguments
96   INTEGER, EXTERNAL :: iargc
97 
98   REAL(sp) :: ValMin, ValMax, InMin, InMax, rmdi
99   REAL(dp), ALLOCATABLE :: mdiVals(:)
100
101   ! NF90_BYTE local data arrays
102   INTEGER(i1), ALLOCATABLE, SAVE, DIMENSION(:) :: localdata_1d_i1
103   INTEGER(i1), ALLOCATABLE, SAVE, DIMENSION(:,:) :: localdata_2d_i1
104   INTEGER(i1), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: localdata_3d_i1
105   INTEGER(i1), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: localdata_4d_i1
106
107   ! NF90_SHORT local data arrays
108   INTEGER(i2), ALLOCATABLE, SAVE, DIMENSION(:) :: localdata_1d_i2
109   INTEGER(i2), ALLOCATABLE, SAVE, DIMENSION(:,:) :: localdata_2d_i2
110   INTEGER(i2), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: localdata_3d_i2
111   INTEGER(i2), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: localdata_4d_i2
112
113   ! NF90_INT local data arrays
114   INTEGER(i4), ALLOCATABLE, SAVE, DIMENSION(:) :: localdata_1d_i4
115   INTEGER(i4), ALLOCATABLE, SAVE, DIMENSION(:,:) :: localdata_2d_i4
116   INTEGER(i4), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: localdata_3d_i4
117   INTEGER(i4), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: localdata_4d_i4
118
119   ! NF90_FLOAT local data arrays
120   REAL(sp), ALLOCATABLE, SAVE, DIMENSION(:) :: localdata_1d_sp
121   REAL(sp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: localdata_2d_sp
122   REAL(sp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: localdata_3d_sp
123   REAL(sp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: localdata_4d_sp
124
125   ! NF90_DOUBLE local data arrays
126   REAL(dp), ALLOCATABLE, SAVE, DIMENSION(:) :: localdata_1d_dp
127   REAL(dp), ALLOCATABLE, SAVE, DIMENSION(:,:) :: localdata_2d_dp
128   REAL(dp), ALLOCATABLE, SAVE, DIMENSION(:,:,:) :: localdata_3d_dp
129   REAL(dp), ALLOCATABLE, SAVE, DIMENSION(:,:,:,:) :: localdata_4d_dp
130
131   ! NF90_BYTE global data arrays
132   INTEGER(i1) :: globaldata_0d_i1
133   INTEGER(i1), ALLOCATABLE, DIMENSION(:) :: globaldata_1d_i1
134   INTEGER(i1), ALLOCATABLE, DIMENSION(:,:) :: globaldata_2d_i1
135   INTEGER(i1), ALLOCATABLE, DIMENSION(:,:,:) :: globaldata_3d_i1
136   INTEGER(i1), ALLOCATABLE, DIMENSION(:,:,:,:) :: globaldata_4d_i1
137
138   ! NF90_SHORT global data arrays
139   INTEGER(i2) :: globaldata_0d_i2
140   INTEGER(i2), ALLOCATABLE, DIMENSION(:) :: globaldata_1d_i2
141   INTEGER(i2), ALLOCATABLE, DIMENSION(:,:) :: globaldata_2d_i2
142   INTEGER(i2), ALLOCATABLE, DIMENSION(:,:,:) :: globaldata_3d_i2
143   INTEGER(i2), ALLOCATABLE, DIMENSION(:,:,:,:) :: globaldata_4d_i2
144
145   ! NF90_INT global data arrays
146   INTEGER(i4) :: globaldata_0d_i4
147   INTEGER(i4), ALLOCATABLE, DIMENSION(:) :: globaldata_1d_i4
148   INTEGER(i4), ALLOCATABLE, DIMENSION(:,:) :: globaldata_2d_i4
149   INTEGER(i4), ALLOCATABLE, DIMENSION(:,:,:) :: globaldata_3d_i4
150   INTEGER(i4), ALLOCATABLE, DIMENSION(:,:,:,:) :: globaldata_4d_i4
151 
152   ! NF90_FLOAT global data arrays
153   REAL(sp) :: globaldata_0d_sp
154   REAL(sp), ALLOCATABLE, DIMENSION(:) :: globaldata_1d_sp
155   REAL(sp), ALLOCATABLE, DIMENSION(:,:) :: globaldata_2d_sp
156   REAL(sp), ALLOCATABLE, DIMENSION(:,:,:) :: globaldata_3d_sp
157   REAL(sp), ALLOCATABLE, DIMENSION(:,:,:,:) :: globaldata_4d_sp
158
159   ! NF90_DOUBLE global data arrays
160   REAL(dp) :: globaldata_0d_dp
161   REAL(dp), ALLOCATABLE, DIMENSION(:) :: globaldata_1d_dp
162   REAL(dp), ALLOCATABLE, DIMENSION(:,:) :: globaldata_2d_dp
163   REAL(dp), ALLOCATABLE, DIMENSION(:,:,:) :: globaldata_3d_dp
164   REAL(dp), ALLOCATABLE, DIMENSION(:,:,:,:) :: globaldata_4d_dp
165
166   LOGICAL :: l_valid = .false.
167   LOGICAL :: l_noRebuild = .false.
168   LOGICAL :: l_findDims = .true.
169   LOGICAL :: l_maskout = .false.
170
171   NAMELIST/nam_rebuild/ filebase, ndomain, dims, nslicesize, l_maskout, deflate_level, &
172                       & nc4_xchunk, nc4_ychunk, nc4_zchunk, nc4_tchunk, fchunksize         
173
174   external      :: getarg
175
176   !End of definitions
177
178!--------------------------------------------------------------------------------
179!0. OMP setup
180
181!$OMP PARALLEL DEFAULT(NONE) SHARED(nthreads)
182!$OMP MASTER
183!$      nthreads = omp_get_num_threads()
184!$      WRITE(numout,*) 'Running OMP with ',nthreads,' thread(s).'
185!$OMP END MASTER
186!$OMP END PARALLEL
187 
188!--------------------------------------------------------------------------------
189!1.0 Check netcdf version for warning
190   clibnc = TRIM(nf90_inq_libvers())
191   IF (ICHAR(clibnc(1:1)) <= 3) THEN
192      PRINT *, '=========================================================='
193      PRINT *, 'You are using old netcdf library (',TRIM(clibnc),').'
194      PRINT *, 'REBUILD_NEMO support of old netcdf library will end soon'
195      PRINT *, 'please consider moving to netcdf 4 or higher'
196      PRINT *, '=========================================================='
197   END IF
198
199!1.1 Get the namelist path
200   !Determine the number of arguments on the command line
201   nargs=iargc()
202   !Check that the required argument is present, if it is not then set it to the default value: nam_rebuild
203   IF (nargs == 0) THEN
204      WRITE(numout,*) 'Namelist path not supplied as command line argument. Using default, nam_rebuild.'
205      cnampath='nam_rebuild'
206   ELSE IF (nargs == 1) THEN
207      CALL getarg(1, cnampath)
208   ELSE
209      WRITE(numerr,*) 'ERROR! : Incorrect number of command line arguments. Please supply only'
210      WRITE(numerr,*) '         the path to the namelist file, or no arguments to use default value'
211      STOP 1
212   END IF
213
214!1.2 Read in the namelist
215
216   dims(:) = ""
217   nslicesize = 0
218   deflate_level = 0
219   OPEN( UNIT=numnam, FILE=TRIM(cnampath), FORM='FORMATTED', STATUS='OLD' )
220   READ( numnam, nam_rebuild )
221   CLOSE( numnam )
222   IF( .NOT. ALL(dims(:) == "") ) l_findDims = .false.
223
224!1.3 Set up the filenames and fileids
225
226   ALLOCATE(filenames(ndomain))
227   IF (l_verbose) WRITE(numout,*) 'Rebuilding the following files:'
228   DO ifile = 1, ndomain
229      WRITE(suffix,'(i4.4)') ifile-1
230      filenames(ifile) = TRIM(filebase)//'_'//TRIM(suffix)//'.nc'
231      IF (l_verbose) WRITE(numout,*) TRIM(filenames(ifile))
232   END DO
233   ALLOCATE(inncids(ndomain))
234 
235!---------------------------------------------------------------------------
236!2. Read in the global dimensions from the first input file and set up the output file
237 
238   CALL check_nf90( nf90_open( TRIM(filenames(1)), nf90_share, ncid ) )
239   CALL check_nf90( nf90_inquire( ncid, ndims, nvars, natts ) )
240   
241!2.0 Read in the total number of processors the file is expecting and check it's correct 
242 
243   CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_number_total', ndomain_file ) )
244   IF( ndomain /= ndomain_file ) THEN
245      WRITE(numerr,*) 'ERROR! : number of files to rebuild in file does not agree with namelist'
246      WRITE(numerr,*) 'Attribute DOMAIN_number_total is : ', ndomain_file
247      WRITE(numerr,*) 'Number of files specified in namelist is: ', ndomain
248      STOP 2
249   ENDIF
250 
251!2.1 Set up the output file
252#if defined key_netcdf4
253   CALL check_nf90( nf90_create( TRIM(filebase)//'.nc', nf90_netcdf4, outid, chunksize=fchunksize ) )
254#else
255   CALL check_nf90( nf90_create( TRIM(filebase)//'.nc', nf90_64bit_offset, outid, chunksize=fchunksize ) )
256#endif
257
258!2.2 Set up dimensions in output file
259
260!2.2.0 Find out how many dimensions are required to be rebuilt and which ones they are
261   CALL check_nf90( nf90_inquire_attribute( ncid, nf90_global, 'DOMAIN_dimensions_ids', xtype, rbdims, attid ) )
262
263   ALLOCATE(rebuild_dims(rbdims))
264   CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_dimensions_ids', rebuild_dims ) )
265
266   ALLOCATE(global_sizes(rbdims))
267   CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_size_global', global_sizes ) )
268   IF (l_verbose) WRITE(numout,*) 'Size of global arrays: ', global_sizes
269
270
271!2.2.1 Copy the dimensions into the output file apart from rebuild_dims() which are dimensioned globally
272   ALLOCATE(indimlens(ndims), indimnames(ndims), outdimlens(ndims))
273   CALL check_nf90( nf90_inquire( ncid, unlimitedDimId = unlimitedDimId ) )
274   istop = 0
275   DO idim = 1, ndims
276      CALL check_nf90( nf90_inquire_dimension( ncid, idim, dimname, dimlen ) )
277      CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_size_local', local_sizes ) )
278      indimlens(idim) = dimlen   
279      indimnames(idim) = dimname
280      IF (l_findDims) THEN
281         IF( idim == rebuild_dims(1) ) THEN
282            IF( dimlen == local_sizes(1) ) THEN
283               dimlen = global_sizes(1)
284               dims(1) = trim(dimname)
285            ELSE
286               istop = 1
287            ENDIF
288         ENDIF
289         IF( rbdims > 1 .AND. idim == rebuild_dims(2) ) THEN
290            IF( dimlen == local_sizes(2) ) THEN
291               dimlen = global_sizes(2)
292               dims(2) = trim(dimname)
293            ELSE
294               istop = 1
295            ENDIF
296         ENDIF
297      ELSE ! l_findDims = false
298         IF( TRIM(dimname) == TRIM(dims(1))) THEN
299            dimlen = global_sizes(1)
300            rebuild_dims(1) = idim
301         ENDIF
302         IF( rbdims > 1 .AND. TRIM(dimname) == TRIM(dims(2))) THEN
303            dimlen = global_sizes(2)
304            rebuild_dims(2) = idim
305         ENDIF
306      ENDIF
307
308      IF( idim == unlimitedDimId ) THEN
309         CALL check_nf90( nf90_def_dim( outid, dimname, nf90_unlimited, dimid) )
310         nmax_unlimited = dimlen
311      ELSE
312         CALL check_nf90( nf90_def_dim( outid, dimname, dimlen, dimid) )
313      ENDIF
314      outdimlens(idim) = dimlen
315   END DO
316   ! nmax_unlimited is only used for time-slicing so we set it to be at least 1 to
317   ! account for files with no record dimension or zero length record dimension(!)
318   nmax_unlimited = max(nmax_unlimited,1)
319
320   IF( istop == 1 ) THEN
321      WRITE(numerr,*) 'ERROR! : DOMAIN_local_sizes attribute does not match rebuild dimension lengths in the first file'
322      WRITE(numerr,*) 'Attribute DOMAIN_local_sizes is : ', local_sizes
323      WRITE(numerr,*) 'Dimensions to be rebuilt are of size : ', outdimlens(rebuild_dims(1)), outdimlens(rebuild_dims(2)) 
324      STOP 3
325   ENDIF
326
327   IF (l_findDims) THEN
328      IF (l_verbose) WRITE(numout,*) 'Finding rebuild dimensions from the first file...'
329   ELSE
330      IF (l_verbose) WRITE(numout,*) 'Using rebuild dimensions given in namelist...'
331   ENDIF
332
333   IF( rbdims > 1 ) THEN
334      IF (l_verbose) WRITE(numout,*) 'Rebuilding across dimensions '//TRIM(indimnames(rebuild_dims(1)))//  &
335         &                      ' and '//TRIM(indimnames(rebuild_dims(2)))
336   ELSE
337      IF (l_verbose) WRITE(numout,*) 'Rebuilding across dimension '//TRIM(indimnames(rebuild_dims(1)))
338   ENDIF
339     
340!2.2.2 Copy the global attributes into the output file, apart from those beginning with DOMAIN_ 
341!      Also need to change the file_name attribute and the TimeStamp attribute.
342   DO attid = 1, natts
343      CALL check_nf90( nf90_inq_attname( ncid, nf90_global, attid, attname ) )
344      IF( INDEX( attname, "DOMAIN_" ) == 1 ) CYCLE
345      IF( INDEX( attname, "file_name") == 1 ) CYCLE
346      IF( INDEX( attname, "associate_file") == 1 ) CYCLE
347      IF (l_verbose) WRITE(numout,*) 'Copying attribute '//TRIM(attname)//' into destination file...'
348      CALL check_nf90( nf90_copy_att( ncid, nf90_global, attname, outid, nf90_global ) )
349   END DO
350   CALL check_nf90( nf90_put_att( outid, nf90_global, "file_name", TRIM(filebase)//'.nc') )
351   IF (l_verbose) WRITE(numout,*) 'Writing new file_name attribute' 
352   CALL DATE_AND_TIME ( date=date, time=time, zone=zone )
353   timestamp = date(7:8) // "/" // date(5:6) // "/" // date(1:4) // " " // &
354               time(1:2) // ":" // time(3:4) // ":" // time(5:6) // " " // &
355               zone 
356   CALL check_nf90( nf90_put_att( outid, nf90_global, "TimeStamp", timestamp ) )
357   IF (l_verbose) WRITE(numout,*) 'Writing new TimeStamp attribute'
358 
359!2.2.3 Copy the variable definitions and attributes into the output file.
360   ALLOCATE(mdiVals(nvars))
361   mdiVals(:)=0
362   DO jv = 1, nvars
363      CALL check_nf90( nf90_inquire_variable( ncid, jv, varname, xtype, ndims, dimids, natts ) )
364      ALLOCATE(outdimids(ndims))
365      ALLOCATE(chunksizes(ndims))
366      IF( ndims > 0 ) then
367        DO idim = 1, ndims
368           outdimids(idim) = dimids(idim)
369           chunksizes(idim) = outdimlens(dimids(idim))
370           cdim='|'//TRIM(indimnames(dimids(idim)))//'|'
371
372! trick to find var in a list of suggestion (var0 and var1 : INDEX(|var0|var1|,|var|)
373           cdimlst='|x|x_grid_T|x_grid_U|x_grid_V|x_grid_W|'
374           if( INDEX(TRIM(cdimlst),TRIM(cdim)) > 0 ) &
375    &                             chunksizes(idim) = min(outdimlens(dimids(idim)), max(nc4_xchunk,1))
376
377           cdimlst='|y|y_grid_T|y_grid_U|y_grid_V|y_grid_W|'
378           if( INDEX(TRIM(cdimlst),TRIM(cdim)) > 0 ) &
379    &                             chunksizes(idim) = min(outdimlens(dimids(idim)), max(nc4_ychunk,1))
380
381           cdimlst='|z|deptht|depthu|depthv|depthw|depth|nav_lev|'
382           if( INDEX(TRIM(cdimlst),TRIM(cdim)) > 0 ) &
383    &                             chunksizes(idim) = min(outdimlens(dimids(idim)), max(nc4_zchunk,1))
384
385           cdimlst='|t|time|time_counter|'
386           if( INDEX(TRIM(cdimlst),TRIM(cdim)) > 0 ) &
387    &                             chunksizes(idim) = min(outdimlens(dimids(idim)), max(nc4_tchunk,1))
388
389        END DO
390#if defined key_netcdf4
391        CALL check_nf90( nf90_def_var( outid, varname, xtype, outdimids, varid, &
392                                       deflate_level=deflate_level ) )
393        IF (l_verbose) WRITE(numout,*) 'Dims    : ',ndims, outdimids(1:ndims)
394        IF (l_verbose) WRITE(numout,*) 'names   : ',(TRIM(indimnames(dimids(idim)))//' ',idim=1,ndims)
395        IF (l_verbose) WRITE(numout,*) 'lens   : ',(outdimlens(dimids(idim)),idim=1,ndims)
396        IF (l_verbose) WRITE(numout,*) 'Chunking: ',chunkalg, chunksizes
397        CALL check_nf90( nf90_def_var_chunking( outid, varid, chunkalg, &
398   &                                 chunksizes ) )
399      ELSE
400        CALL check_nf90( nf90_def_var( outid, varname, xtype, outdimids, varid ) )
401#else
402      CALL check_nf90( nf90_def_var( outid, varname, xtype, outdimids, varid ) )
403#endif
404      ENDIF
405      DEALLOCATE(outdimids)
406      DEALLOCATE(chunksizes)
407      IF (l_verbose) WRITE(numout,*) 'Defining variable '//TRIM(varname)//'...' 
408      IF( natts > 0 ) THEN
409         DO attid = 1, natts
410            CALL check_nf90( nf90_inq_attname( ncid, varid, attid, attname ) )
411            IF ( attname == "_FillValue" ) THEN
412               CALL check_nf90( nf90_get_att( ncid, varid, attname, rmdi ) )
413               mdiVals(jv)=rmdi
414            ENDIF
415            CALL check_nf90( nf90_copy_att( ncid, varid, attname, outid, varid ) )
416         END DO
417      ENDIF
418   END DO
419
420!2.3 End definitions in output file and copy 1st file ncid to the inncids array
421
422   CALL check_nf90( nf90_enddef( outid ) )
423   inncids(1) = ncid
424   IF (l_verbose) WRITE(numout,*) 'Finished defining output file.'
425 
426!---------------------------------------------------------------------------
427!3. Read in data from each file for each variable
428
429!3.1 Open each file and store the ncid in inncids array
430
431   IF (l_verbose) WRITE(numout,*) 'Opening input files...'
432
433   ! Set a file chunk cache size for the processor-domain files that scales with the number of processors
434   patchchunk = max(8192, fchunksize/ndomain)
435
436   ! open files
437   DO ifile = 2, ndomain
438      CALL check_nf90( nf90_open( TRIM(filenames(ifile)), nf90_share, ncid, chunksize=patchchunk ) )
439      inncids(ifile) = ncid
440   END DO
441   IF (l_verbose) WRITE(numout,*) 'All input files open.'
442
443   DO jv = 1, nvars
444
445      ValMin = 1.e10
446      ValMax = -1.e10
447      l_valid = .false.
448      istop = nf90_noerr
449      nt = 1
450      ntslice = nmax_unlimited
451      IF( nslicesize == 0 ) nslicesize = nmax_unlimited
452
453!3.2 Inquire variable to find out name and how many dimensions it has
454!    and importantly whether it contains the dimensions in rebuild_dims()
455
456      ncid = inncids(1)
457      CALL check_nf90( nf90_inquire_variable( ncid, jv, varname, xtype, ndims, dimids, natts ) )
458
459      l_noRebuild = .true.
460      IF( ANY( dimids(1:ndims) == rebuild_dims(1) )) l_noRebuild = .false.
461      IF( rbdims > 1 ) THEN
462         IF( ANY( dimids(1:ndims) == rebuild_dims(2) )) l_noRebuild = .false.
463      ENDIF
464
465!3.2.0 start while loop for time slicing
466
467      DO WHILE( nt <= nmax_unlimited )
468
469         IF( ndims > 3 ) THEN
470            ntslice = MIN( nslicesize, nmax_unlimited + 1 - nt )
471         ENDIF
472
473      IF (l_noRebuild) THEN
474
475         IF( nslicesize == nmax_unlimited .OR. ndims <= 3 ) THEN
476            IF (l_verbose) WRITE(numout,*) 'Copying data from variable '//TRIM(varname)//'...'
477         ELSE
478            IF (l_verbose) WRITE(numout,'(A,I3,A,I3,A)') ' Copying data from variable '  &
479            &                 //TRIM(varname)//' for slices ',nt,' to ',nt+ntslice-1,' ...'
480         ENDIF
481
482!3.2.1 If rebuilding not required then just need to read in variable
483!      for copying direct into output file after the OMP (files) loop.
484         IF( ndims == 0 ) THEN
485
486            SELECT CASE( xtype )
487               CASE( NF90_BYTE )
488                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_i1 ) )
489               CASE( NF90_SHORT )
490                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_i2 ) )
491               CASE( NF90_INT )
492                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_i4 ) )
493               CASE( NF90_FLOAT )
494                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_sp ) )
495               CASE( NF90_DOUBLE )
496                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_dp ) )
497               CASE DEFAULT
498                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
499                  STOP 4
500            END SELECT
501
502         ELSEIF( ndims == 1 ) THEN
503
504            SELECT CASE( xtype )
505               CASE( NF90_BYTE )
506                  ALLOCATE(globaldata_1d_i1(indimlens(dimids(1))))
507                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_i1 ) )
508               CASE( NF90_SHORT )
509                  ALLOCATE(globaldata_1d_i2(indimlens(dimids(1))))
510                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_i2 ) )
511               CASE( NF90_INT )
512                  ALLOCATE(globaldata_1d_i4(indimlens(dimids(1))))
513                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_i4 ) )
514               CASE( NF90_FLOAT )
515                  ALLOCATE(globaldata_1d_sp(indimlens(dimids(1))))
516                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_sp ) )
517               CASE( NF90_DOUBLE )
518                  ALLOCATE(globaldata_1d_dp(indimlens(dimids(1))))
519                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_dp ) )
520               CASE DEFAULT
521                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
522                  STOP 4
523            END SELECT
524
525         ELSEIF( ndims == 2 ) THEN
526
527            SELECT CASE( xtype )
528               CASE( NF90_BYTE )
529                  ALLOCATE(globaldata_2d_i1(indimlens(dimids(1)),indimlens(dimids(2))))
530                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_i1 ) )
531               CASE( NF90_SHORT )
532                  ALLOCATE(globaldata_2d_i2(indimlens(dimids(1)),indimlens(dimids(2))))
533                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_i2 ) )
534               CASE( NF90_INT )
535                  ALLOCATE(globaldata_2d_i4(indimlens(dimids(1)),indimlens(dimids(2))))
536                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_i4 ) )
537               CASE( NF90_FLOAT )
538                  ALLOCATE(globaldata_2d_sp(indimlens(dimids(1)),indimlens(dimids(2))))
539                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_sp ) )
540               CASE( NF90_DOUBLE )
541                  ALLOCATE(globaldata_2d_dp(indimlens(dimids(1)),indimlens(dimids(2))))
542                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_dp ) )
543               CASE DEFAULT
544                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
545                  STOP 4
546            END SELECT
547
548         ELSEIF( ndims == 3 ) THEN
549
550            SELECT CASE( xtype )
551               CASE( NF90_BYTE )
552                  ALLOCATE(globaldata_3d_i1(indimlens(dimids(1)),indimlens(dimids(2)),       &
553                     &                      indimlens(dimids(3))))
554                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_i1 ) )
555               CASE( NF90_SHORT )
556                  ALLOCATE(globaldata_3d_i2(indimlens(dimids(1)),indimlens(dimids(2)),       &
557                     &                      indimlens(dimids(3))))
558                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_i2 ) )
559               CASE( NF90_INT )
560                  ALLOCATE(globaldata_3d_i4(indimlens(dimids(1)),indimlens(dimids(2)),       &
561                     &                      indimlens(dimids(3))))
562                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_i4 ) )
563               CASE( NF90_FLOAT )
564                  ALLOCATE(globaldata_3d_sp(indimlens(dimids(1)),indimlens(dimids(2)),       &
565                     &                      indimlens(dimids(3))))
566                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_sp ) )
567               CASE( NF90_DOUBLE )
568                  ALLOCATE(globaldata_3d_dp(indimlens(dimids(1)),indimlens(dimids(2)),       &
569                     &                      indimlens(dimids(3))))
570                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_dp ) )
571               CASE DEFAULT
572                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
573                  STOP 4
574            END SELECT
575
576         ELSEIF( ndims == 4 ) THEN
577
578            SELECT CASE( xtype )
579               CASE( NF90_BYTE )
580                  ALLOCATE(globaldata_4d_i1(indimlens(dimids(1)),indimlens(dimids(2)),       &
581                     &                      indimlens(dimids(3)),ntslice))
582                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_i1, start=(/1,1,1,nt/) ) )
583               CASE( NF90_SHORT )
584                  ALLOCATE(globaldata_4d_i2(indimlens(dimids(1)),indimlens(dimids(2)),       &
585                     &                      indimlens(dimids(3)),ntslice))
586                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_i2, start=(/1,1,1,nt/) ) )
587               CASE( NF90_INT )
588                  ALLOCATE(globaldata_4d_i4(indimlens(dimids(1)),indimlens(dimids(2)),       &
589                     &                      indimlens(dimids(3)),ntslice))
590                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_i4, start=(/1,1,1,nt/) ) )
591               CASE( NF90_FLOAT )
592                  ALLOCATE(globaldata_4d_sp(indimlens(dimids(1)),indimlens(dimids(2)),       &
593                     &                      indimlens(dimids(3)),ntslice))
594                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_sp, start=(/1,1,1,nt/) ) )
595               CASE( NF90_DOUBLE )
596                  ALLOCATE(globaldata_4d_dp(indimlens(dimids(1)),indimlens(dimids(2)),       &
597                     &                      indimlens(dimids(3)),ntslice))
598                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_dp, start=(/1,1,1,nt/) ) )
599               CASE DEFAULT
600                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
601                  STOP 4
602            END SELECT
603
604         ENDIF
605
606      ELSE  ! l_noRebuild = .false.
607
608!3.2.2 For variables that require rebuilding we need to read in from all ndomain files
609!      Here we allocate global variables ahead of looping over files
610         IF( nslicesize == nmax_unlimited .OR. ndims <= 3 ) THEN
611            IF (l_verbose) WRITE(numout,*) 'Rebuilding data from variable '//TRIM(varname)//'...'
612         ELSE
613            IF (l_verbose) WRITE(numout,'(A,I3,A,I3,A)') ' Rebuilding data from variable '  &
614            &                 //TRIM(varname)//' for slices ',nt,' to ',nt+ntslice-1,' ...'
615         ENDIF
616         IF( ndims == 1 ) THEN
617
618            SELECT CASE( xtype )
619               CASE( NF90_BYTE )
620                  ALLOCATE(globaldata_1d_i1(outdimlens(dimids(1))))
621                  IF (l_maskout) globaldata_1d_i1(:)=mdiVals(jv)
622               CASE( NF90_SHORT )
623                  ALLOCATE(globaldata_1d_i2(outdimlens(dimids(1))))
624                  IF (l_maskout) globaldata_1d_i2(:)=mdiVals(jv)
625               CASE( NF90_INT )
626                  ALLOCATE(globaldata_1d_i4(outdimlens(dimids(1))))
627                  IF (l_maskout) globaldata_1d_i4(:)=mdiVals(jv)
628               CASE( NF90_FLOAT )
629                  ALLOCATE(globaldata_1d_sp(outdimlens(dimids(1))))
630                  IF (l_maskout) globaldata_1d_sp(:)=mdiVals(jv)
631               CASE( NF90_DOUBLE )
632                  ALLOCATE(globaldata_1d_dp(outdimlens(dimids(1))))
633                  IF (l_maskout) globaldata_1d_dp(:)=mdiVals(jv)
634               CASE DEFAULT
635                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
636                  STOP 4
637            END SELECT
638
639         ELSEIF( ndims == 2 ) THEN
640
641            SELECT CASE( xtype )
642               CASE( NF90_BYTE )
643                  ALLOCATE(globaldata_2d_i1(outdimlens(dimids(1)),outdimlens(dimids(2))))
644                  IF (l_maskout) globaldata_2d_i1(:,:)=mdiVals(jv)
645               CASE( NF90_SHORT )
646                  ALLOCATE(globaldata_2d_i2(outdimlens(dimids(1)),outdimlens(dimids(2))))
647                  IF (l_maskout) globaldata_2d_i2(:,:)=mdiVals(jv)
648               CASE( NF90_INT )
649                  ALLOCATE(globaldata_2d_i4(outdimlens(dimids(1)),outdimlens(dimids(2))))
650                  IF (l_maskout) globaldata_2d_i4(:,:)=mdiVals(jv)
651               CASE( NF90_FLOAT )
652                  ALLOCATE(globaldata_2d_sp(outdimlens(dimids(1)),outdimlens(dimids(2))))
653                  IF (l_maskout) globaldata_2d_sp(:,:)=mdiVals(jv)
654               CASE( NF90_DOUBLE )
655                  ALLOCATE(globaldata_2d_dp(outdimlens(dimids(1)),outdimlens(dimids(2))))
656                  IF (l_maskout) globaldata_2d_dp(:,:)=mdiVals(jv)
657               CASE DEFAULT
658                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
659                  STOP 4
660            END SELECT
661
662         ELSEIF( ndims == 3 ) THEN
663
664            SELECT CASE( xtype )
665               CASE( NF90_BYTE )
666                  ALLOCATE(globaldata_3d_i1(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
667                     &                      outdimlens(dimids(3))))
668                  IF (l_maskout) globaldata_3d_i1(:,:,:)=mdiVals(jv)
669               CASE( NF90_SHORT )
670                  ALLOCATE(globaldata_3d_i2(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
671                     &                      outdimlens(dimids(3))))
672                  IF (l_maskout) globaldata_3d_i2(:,:,:)=mdiVals(jv)
673               CASE( NF90_INT )
674                  ALLOCATE(globaldata_3d_i4(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
675                     &                      outdimlens(dimids(3))))
676                  IF (l_maskout) globaldata_3d_i4(:,:,:)=mdiVals(jv)
677               CASE( NF90_FLOAT )
678                  ALLOCATE(globaldata_3d_sp(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
679                     &                      outdimlens(dimids(3))))
680                  IF (l_maskout) globaldata_3d_sp(:,:,:)=mdiVals(jv)
681               CASE( NF90_DOUBLE )
682                  ALLOCATE(globaldata_3d_dp(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
683                     &                      outdimlens(dimids(3))))
684                  IF (l_maskout) globaldata_3d_dp(:,:,:)=mdiVals(jv)
685               CASE DEFAULT
686                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
687                  STOP 4
688            END SELECT
689
690         ELSEIF( ndims == 4 ) THEN
691
692            SELECT CASE( xtype )
693               CASE( NF90_BYTE )
694                  ALLOCATE(globaldata_4d_i1(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
695                     &                      outdimlens(dimids(3)),ntslice))
696                  IF (l_maskout) globaldata_4d_i1(:,:,:,:)=mdiVals(jv)
697               CASE( NF90_SHORT )
698                  ALLOCATE(globaldata_4d_i2(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
699                     &                      outdimlens(dimids(3)),ntslice))
700                  IF (l_maskout) globaldata_4d_i2(:,:,:,:)=mdiVals(jv)
701               CASE( NF90_INT )
702                  ALLOCATE(globaldata_4d_i4(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
703                     &                      outdimlens(dimids(3)),ntslice))
704                  IF (l_maskout) globaldata_4d_i4(:,:,:,:)=mdiVals(jv)
705               CASE( NF90_FLOAT )
706                  ALLOCATE(globaldata_4d_sp(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
707                     &                      outdimlens(dimids(3)),ntslice))
708                  IF (l_maskout) globaldata_4d_sp(:,:,:,:)=mdiVals(jv)
709               CASE( NF90_DOUBLE )
710                  ALLOCATE(globaldata_4d_dp(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
711                     &                      outdimlens(dimids(3)),ntslice))
712                  IF (l_maskout) globaldata_4d_dp(:,:,:,:)=mdiVals(jv)
713               CASE DEFAULT
714                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
715                  STOP 4
716            END SELECT
717         ELSE
718            WRITE(numerr,*) 'ERROR! : A netcdf variable has more than 4 dimensions which is not taken into account'
719            STOP 4
720         ENDIF
721
722!$OMP  PARALLEL DO DEFAULT(NONE)                                                          &
723!$OMP& PRIVATE(ifile,ncid,xtype,start_pos,local_sizes,InMin,InMax,natts,                  &
724!$OMP&         ndims,attid,attname,dimids,idim,dimname,dimlen,unlimitedDimId,             &
725!$OMP&         halo_start,halo_end,idomain,jdomain,rdomain,di,dj,dr,                      &
726!$OMP&         localdata_1d_i2,localdata_1d_i4,localdata_1d_sp,localdata_1d_dp,           &
727!$OMP&         localdata_2d_i2,localdata_2d_i4,localdata_2d_sp,localdata_2d_dp,           &
728!$OMP&         localdata_3d_i2,localdata_3d_i4,localdata_3d_sp,localdata_3d_dp,           &
729!$OMP&         localdata_4d_i2,localdata_4d_i4,localdata_4d_sp,localdata_4d_dp,           &
730!$OMP&         localdata_1d_i1,localdata_2d_i1,localdata_3d_i1,localdata_4d_i1)           &
731!$OMP& SHARED(jv,nvars,varname,filenames,ValMin,ValMax,indimlens,outdimlens,rbdims,       &
732!$OMP&        ndomain,outid,fchunksize,istop,l_valid,nthreads,inncids,rebuild_dims,       &
733!$OMP&        globaldata_1d_i2,globaldata_1d_i4,globaldata_1d_sp,globaldata_1d_dp,        &
734!$OMP&        globaldata_2d_i2,globaldata_2d_i4,globaldata_2d_sp,globaldata_2d_dp,        &
735!$OMP&        globaldata_3d_i2,globaldata_3d_i4,globaldata_3d_sp,globaldata_3d_dp,        &
736!$OMP&        globaldata_4d_i2,globaldata_4d_i4,globaldata_4d_sp,globaldata_4d_dp,        &
737!$OMP&        globaldata_1d_i1,globaldata_2d_i1,globaldata_3d_i1,globaldata_4d_i1,        &
738!$OMP&        ntslice,nt,nmax_unlimited,indimnames,dims,patchchunk)
739
740         DO ifile = 1, ndomain
741
742            ncid = inncids(ifile)
743!$OMP CRITICAL
744            CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_size_local', local_sizes ), istop )
745            CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_position_first', start_pos ), istop )
746            CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_halo_size_start', halo_start ), istop )
747            CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_halo_size_end', halo_end ), istop )
748            CALL check_nf90( nf90_inquire_variable( ncid, jv, varname, xtype, ndims, dimids, natts ), istop )
749            CALL check_nf90( nf90_inquire( ncid, unlimitedDimId = unlimitedDimId ), istop )
750!$OMP END CRITICAL
751
752            ! set defaults for rebuilding so that i is 1st, j 2nd
753            di=1
754            dj=2
755
756            IF( rbdims == 1 ) THEN
757               ! override defaults above and set other variables
758               start_pos(2) = 1
759               local_sizes(2) = outdimlens(3-dimids(2))
760               halo_end(2) = 0
761               halo_start(2) = 0
762               di=rebuild_dims(1)
763               dj=3-di
764            ENDIF
765
766!3.3.1 Generate local domain interior sizes from local_sizes and halo sizes
767!      idomain defines the 1st and last interior points in the i direction and
768!      jdomain defines the 1st and last interior points in the j direction
769
770            idomain(1) = 1 + halo_start(di)
771            idomain(2) = local_sizes(di) - halo_end(di)
772            jdomain(1) = 1 + halo_start(dj)
773            jdomain(2) = local_sizes(dj) - halo_end(dj)
774
775!3.3.2 For rbdims or more dimensions put the data array from this input file into the correct
776!      part of the output data array. Assume the first dimensions are those to be rebuilt.
777
778            IF( ndims == 1 ) THEN
779
780               IF( rebuild_dims(1) == 1 ) THEN
781                  dr = di
782                  rdomain = idomain
783               ELSE
784                  dr = dj
785                  rdomain = jdomain
786               ENDIF
787
788              SELECT CASE( xtype )
789                  CASE( NF90_BYTE )
790                     ALLOCATE(localdata_1d_i1(local_sizes(dr)))
791                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_i1 ), istop )
792                     DO jr = rdomain(1), rdomain(2)
793                        globaldata_1d_i1(start_pos(dr) + jr - 1) = localdata_1d_i1(jr)
794                     END DO
795                     DEALLOCATE(localdata_1d_i1)
796                  CASE( NF90_SHORT )
797                     ALLOCATE(localdata_1d_i2(local_sizes(dr)))
798                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_i2 ), istop )
799                     DO jr = rdomain(1), rdomain(2)
800                        globaldata_1d_i2(start_pos(dr) + jr - 1) = localdata_1d_i2(jr)
801                     END DO
802                     DEALLOCATE(localdata_1d_i2)
803                  CASE( NF90_INT )
804                     ALLOCATE(localdata_1d_i4(local_sizes(dr)))
805                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_i4 ), istop )
806                     DO jr = rdomain(1), rdomain(2)
807                        globaldata_1d_i4(start_pos(dr) + jr - 1) = localdata_1d_i4(jr)
808                     END DO
809                     DEALLOCATE(localdata_1d_i4)
810                  CASE( NF90_FLOAT )
811                     ALLOCATE(localdata_1d_sp(local_sizes(dr)))
812                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_sp ), istop )
813                     DO jr = rdomain(1), rdomain(2)
814                        globaldata_1d_sp(start_pos(dr) + jr - 1) = localdata_1d_sp(jr)
815                     END DO
816                     DEALLOCATE(localdata_1d_sp)
817                  CASE( NF90_DOUBLE )
818                     ALLOCATE(localdata_1d_dp(local_sizes(dr)))
819                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_dp ), istop )
820                     DO jr = rdomain(1), rdomain(2)
821                        globaldata_1d_dp(start_pos(dr) + jr - 1) = localdata_1d_dp(jr)
822                     END DO
823                     DEALLOCATE(localdata_1d_dp)
824                  CASE DEFAULT
825                     WRITE(numerr,*) 'Unknown nf90 type: ', xtype
826                     istop = istop + 1
827               END SELECT
828
829            ELSEIF( ndims == 2 ) THEN
830       
831               SELECT CASE( xtype )
832                  CASE( NF90_BYTE )
833                     ALLOCATE(localdata_2d_i1(local_sizes(di),local_sizes(dj)))
834                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_i1 ), istop )
835                        DO jj = jdomain(1), jdomain(2)
836                           DO ji = idomain(1), idomain(2)
837                              globaldata_2d_i1(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_i1(ji,jj)
838                        END DO
839                     END DO
840                     DEALLOCATE(localdata_2d_i1)
841                  CASE( NF90_SHORT )
842                     ALLOCATE(localdata_2d_i2(local_sizes(di),local_sizes(dj)))
843                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_i2 ), istop )
844                     DO jj = jdomain(1), jdomain(2)
845                        DO ji = idomain(1), idomain(2)
846                           globaldata_2d_i2(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_i2(ji,jj)
847                        END DO
848                     END DO
849                     DEALLOCATE(localdata_2d_i2)
850                  CASE( NF90_INT )
851                     ALLOCATE(localdata_2d_i4(local_sizes(di),local_sizes(dj)))
852                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_i4 ), istop )
853                     DO jj = jdomain(1), jdomain(2)
854                        DO ji = idomain(1), idomain(2)
855                           globaldata_2d_i4(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_i4(ji,jj)
856                        END DO
857                     END DO
858                     DEALLOCATE(localdata_2d_i4)
859                  CASE( NF90_FLOAT )
860                     ALLOCATE(localdata_2d_sp(local_sizes(di),local_sizes(dj)))
861                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_sp ), istop )
862                     DO jj = jdomain(1), jdomain(2)
863                        DO ji = idomain(1), idomain(2)
864                           globaldata_2d_sp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_sp(ji,jj)
865                        END DO
866                     END DO
867                     DEALLOCATE(localdata_2d_sp)
868                  CASE( NF90_DOUBLE )
869                     ALLOCATE(localdata_2d_dp(local_sizes(di),local_sizes(dj)))
870                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_dp ), istop )
871                     DO jj = jdomain(1), jdomain(2)
872                        DO ji = idomain(1), idomain(2)
873                           globaldata_2d_dp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_dp(ji,jj)
874                        END DO
875                     END DO
876                     DEALLOCATE(localdata_2d_dp) 
877                  CASE DEFAULT
878                     WRITE(numerr,*) 'Unknown nf90 type: ', xtype
879                     istop = istop + 1
880               END SELECT
881
882            ELSEIF( ndims == 3 ) THEN
883
884               SELECT CASE( xtype )
885                  CASE( NF90_BYTE )
886                     ALLOCATE(localdata_3d_i1(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
887                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_i1 ), istop )
888!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
889!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_i1,localdata_3d_i1,di,dj)
890                     DO jk = 1, indimlens(dimids(3))
891                        DO jj = jdomain(1), jdomain(2)
892                           DO ji = idomain(1), idomain(2)
893                              globaldata_3d_i1(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_i1(ji,jj,jk)
894                           END DO
895                        END DO
896                     END DO
897!$OMP END PARALLEL DO
898                     DEALLOCATE(localdata_3d_i1)
899                  CASE( NF90_SHORT )
900                     ALLOCATE(localdata_3d_i2(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
901                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_i2 ), istop )
902!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
903!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_i2,localdata_3d_i2,di,dj)
904                     DO jk = 1, indimlens(dimids(3))
905                        DO jj = jdomain(1), jdomain(2)
906                           DO ji = idomain(1), idomain(2)
907                              globaldata_3d_i2(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_i2(ji,jj,jk)
908                           END DO
909                        END DO
910                     END DO
911!$OMP END PARALLEL DO
912                     DEALLOCATE(localdata_3d_i2)
913                  CASE( NF90_INT )
914                     ALLOCATE(localdata_3d_i4(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
915                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_i4 ), istop )
916!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
917!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_i4,localdata_3d_i4,di,dj)
918                     DO jk = 1, indimlens(dimids(3))
919                        DO jj = jdomain(1), jdomain(2)
920                           DO ji = idomain(1), idomain(2)
921                              globaldata_3d_i4(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_i4(ji,jj,jk)
922                           END DO
923                        END DO
924                     END DO
925!$OMP END PARALLEL DO
926                     DEALLOCATE(localdata_3d_i4)
927                  CASE( NF90_FLOAT )
928                     ! TG: This if statement is added to check if the 1st dimension is the corners (for lon_bounds) variables
929                     ! TG: Had to add the unsatisfactory check for 'lon' as it failed for diaptr files
930                     ! TG: Would like to find a better assumption for this.
931                     IF ( trim(indimnames(dimids(1))) /= dims(1) .AND. indimnames(dimids(1)) .NE. 'lon' ) THEN
932                        ALLOCATE(localdata_3d_sp(indimlens(dimids(1)),local_sizes(di),local_sizes(dj)))
933                        WRITE(*,*) 'test', ifile, jv, indimlens(dimids(1)),local_sizes(di),local_sizes(dj)
934                        CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_sp ), istop )
935                        WRITE(*,*) 'test2'
936                        DO jj = jdomain(1), jdomain(2)
937                           DO ji = idomain(1), idomain(2)
938                              DO jk = 1, indimlens(dimids(1))
939                                 globaldata_3d_sp(jk, start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_3d_sp(jk,ji,jj)
940                              END DO
941                           END DO
942                        END DO
943                     ELSE
944                        ALLOCATE(localdata_3d_sp(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
945                        CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_sp ), istop ) 
946!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
947!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_sp,localdata_3d_sp,di,dj)
948                        DO jk = 1, indimlens(dimids(3))
949                           DO jj = jdomain(1), jdomain(2)
950                              DO ji = idomain(1), idomain(2)
951                                 globaldata_3d_sp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_sp(ji,jj,jk)
952                              END DO
953                           END DO
954                        END DO
955!$OMP END PARALLEL DO
956                     ENDIF
957                     DEALLOCATE(localdata_3d_sp)
958                  CASE( NF90_DOUBLE )
959                     IF ( trim(indimnames(dimids(1))) /= dims(1) ) THEN
960                        ALLOCATE(localdata_3d_dp(indimlens(dimids(1)),local_sizes(di),local_sizes(dj)))
961                        CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_dp ), istop ) 
962                        DO jj = jdomain(1), jdomain(2)
963                           DO ji = idomain(1), idomain(2)
964                              DO jk = 1, indimlens(dimids(1))
965                                 globaldata_3d_dp(jk, start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_3d_dp(jk,ji,jj)
966                              END DO
967                           END DO
968                        END DO
969                     ELSE
970                        ALLOCATE(localdata_3d_dp(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
971                        CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_dp ), istop ) 
972!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
973!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_dp,localdata_3d_dp,di,dj)
974                        DO jk = 1, indimlens(dimids(3))
975                           DO jj = jdomain(1), jdomain(2)
976                              DO ji = idomain(1), idomain(2)
977                                 globaldata_3d_dp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_dp(ji,jj,jk)
978                              END DO
979                           END DO
980                        END DO
981!$OMP END PARALLEL DO
982                     ENDIF
983                     DEALLOCATE(localdata_3d_dp)
984                  CASE DEFAULT
985                     WRITE(numerr,*) 'Unknown nf90 type: ', xtype
986                     istop = istop + 1
987               END SELECT
988     
989            ELSEIF (ndims == 4) THEN
990       
991               SELECT CASE( xtype )
992                  CASE( NF90_BYTE )
993                     ALLOCATE(localdata_4d_i1(local_sizes(di),local_sizes(dj),               &
994                         &                     indimlens(dimids(3)),ntslice))
995                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_i1, start=(/1,1,1,nt/) ), istop )
996!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
997!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_i1,localdata_4d_i1,di,dj,nt,ntslice)
998                     DO jl = 1, ntslice
999!$OMP DO
1000                        DO jk = 1, indimlens(dimids(3))
1001                           DO jj = jdomain(1), jdomain(2)
1002                              DO ji = idomain(1), idomain(2)
1003                                 globaldata_4d_i1(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_i1(ji,jj,jk,jl)
1004                              END DO
1005                           END DO
1006                        END DO
1007!$OMP END DO nowait
1008                     END DO
1009!$OMP END PARALLEL
1010                     DEALLOCATE(localdata_4d_i1)
1011                  CASE( NF90_SHORT )
1012                     ALLOCATE(localdata_4d_i2(local_sizes(di),local_sizes(dj),               &
1013                        &                     indimlens(dimids(3)),ntslice))
1014                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_i2, start=(/1,1,1,nt/) ), istop )
1015!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
1016!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_i2,localdata_4d_i2,di,dj,nt,ntslice)
1017                     DO jl = 1, ntslice
1018!$OMP DO
1019                        DO jk = 1, indimlens(dimids(3))
1020                           DO jj = jdomain(1), jdomain(2) 
1021                              DO ji = idomain(1), idomain(2)
1022                                 globaldata_4d_i2(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_i2(ji,jj,jk,jl)
1023                              END DO
1024                           END DO
1025                        END DO
1026!$OMP END DO nowait
1027                     END DO
1028!$OMP END PARALLEL
1029                     DEALLOCATE(localdata_4d_i2)
1030                  CASE( NF90_INT )
1031                     ALLOCATE(localdata_4d_i4(local_sizes(di),local_sizes(dj),               &
1032                        &                     indimlens(dimids(3)),ntslice))
1033                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_i4, start=(/1,1,1,nt/) ), istop )
1034!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
1035!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_i4,localdata_4d_i4,di,dj,nt,ntslice)
1036                     DO jl = 1, ntslice
1037!$OMP DO
1038                        DO jk = 1, indimlens(dimids(3))
1039                           DO jj = jdomain(1), jdomain(2) 
1040                              DO ji = idomain(1), idomain(2)
1041                                 globaldata_4d_i4(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_i4(ji,jj,jk,jl)
1042                              END DO
1043                           END DO
1044                        END DO
1045!$OMP END DO nowait
1046                     END DO
1047!$OMP END PARALLEL
1048                     DEALLOCATE(localdata_4d_i4)
1049                  CASE( NF90_FLOAT )
1050                     ALLOCATE(localdata_4d_sp(local_sizes(di),local_sizes(dj),               &
1051                        &                     indimlens(dimids(3)),ntslice))
1052                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_sp, start=(/1,1,1,nt/) ), istop )
1053!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
1054!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_sp,localdata_4d_sp,di,dj,nt,ntslice)
1055                     DO jl = 1, ntslice
1056!$OMP DO
1057                        DO jk = 1, indimlens(dimids(3))
1058                           DO jj = jdomain(1), jdomain(2) 
1059                              DO ji = idomain(1), idomain(2)
1060                                 globaldata_4d_sp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_sp(ji,jj,jk,jl)
1061                              END DO
1062                           END DO
1063                        END DO
1064!$OMP END DO nowait
1065                     END DO
1066!$OMP END PARALLEL
1067                     DEALLOCATE(localdata_4d_sp)
1068                  CASE( NF90_DOUBLE )
1069                     ALLOCATE(localdata_4d_dp(local_sizes(di),local_sizes(dj),               &
1070                        &                     indimlens(dimids(3)),ntslice))
1071                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_dp, start=(/1,1,1,nt/) ), istop )
1072!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
1073!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_dp,localdata_4d_dp,di,dj,nt,ntslice)
1074                     DO jl = 1, ntslice
1075!$OMP DO
1076                        DO jk = 1, indimlens(dimids(3))
1077                           DO jj = jdomain(1), jdomain(2) 
1078                              DO ji = idomain(1), idomain(2)
1079                                 globaldata_4d_dp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_dp(ji,jj,jk,jl)
1080                              END DO
1081                           END DO
1082                        END DO
1083!$OMP END DO nowait
1084                     END DO
1085!$OMP END PARALLEL
1086                     DEALLOCATE(localdata_4d_dp)
1087                  CASE DEFAULT
1088                     WRITE(numerr,*) 'Unknown nf90 type: ', xtype
1089                     istop = istop + 1
1090               END SELECT
1091
1092            ENDIF ! l_noRebuild false
1093
1094!3.4 Work out if the valid_min and valid_max attributes exist for this variable.
1095!    If they do then calculate the extrema over all input files.
1096
1097            DO attid = 1, natts
1098               CALL check_nf90( nf90_inq_attname( ncid, jv, attid, attname ), istop )
1099               IF( INDEX( attname, "valid_min" ) == 1 ) THEN
1100                  CALL check_nf90( nf90_get_att( ncid, jv, attname, InMin), istop )
1101                  l_valid = .true.
1102               ENDIF
1103               IF( INDEX( attname, "valid_max" ) == 1 ) THEN
1104                  CALL check_nf90( nf90_get_att( ncid, jv, attname, InMax ), istop )
1105                  l_valid = .true.
1106               ENDIF
1107            END DO
1108
1109            IF (l_valid) THEN
1110!$OMP CRITICAL
1111               IF( InMin < ValMin ) ValMin = InMin
1112               IF( InMax > ValMax ) ValMax = InMax
1113!$OMP END CRITICAL
1114            ENDIF
1115
1116!3.5 Abort if failure and only 1 thread
1117
1118            IF( nthreads == 1 .AND. istop /= nf90_noerr )  THEN
1119               WRITE(numerr,*) '*** NEMO rebuild failed! ***'
1120               STOP 5
1121            ENDIF
1122       
1123         END DO  ! loop over files
1124!$OMP END PARALLEL DO
1125
1126!3.6 Abort if any of the OMP threads failed
1127         IF( istop /= nf90_noerr )  THEN
1128            WRITE(numerr,*) '*** NEMO rebuild failed! ***'
1129            STOP 5
1130         ENDIF
1131
1132      ENDIF ! ndims > 2
1133
1134!---------------------------------------------------------------------------
1135!4. Write data to output file
1136
1137      IF (l_verbose) WRITE(numout,*) 'Writing variable '//TRIM(varname)//'...'
1138
1139!4.1 If the valid min and max attributes exist then update them in the file
1140
1141      IF( l_valid ) THEN
1142         CALL check_nf90( nf90_put_att( outid, jv, "valid_min", ValMin ) )
1143         CALL check_nf90( nf90_put_att( outid, jv, "valid_max", ValMax ) )
1144      ENDIF
1145
1146!4.2 Write the data to the output file depending on how many dimensions it has
1147
1148      IF( ndims == 0 ) THEN
1149
1150         SELECT CASE( xtype )
1151            CASE( NF90_BYTE )
1152               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_i1 ) )
1153            CASE( NF90_SHORT )
1154               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_i2 ) )
1155            CASE( NF90_INT )
1156               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_i4 ) )
1157            CASE( NF90_FLOAT )
1158               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_sp ) )
1159            CASE( NF90_DOUBLE )
1160               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_dp ) )
1161            CASE DEFAULT   
1162               WRITE(numerr,*) '0d Unknown nf90 type: ', xtype
1163               STOP 4
1164         END SELECT
1165
1166      ELSEIF( ndims == 1 ) THEN
1167
1168         SELECT CASE( xtype )
1169            CASE( NF90_BYTE )
1170               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_i1 ) )
1171               DEALLOCATE(globaldata_1d_i1)
1172            CASE( NF90_SHORT )
1173               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_i2 ) )
1174               DEALLOCATE(globaldata_1d_i2)
1175            CASE( NF90_INT )
1176               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_i4 ) )
1177               DEALLOCATE(globaldata_1d_i4)
1178            CASE( NF90_FLOAT )
1179               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_sp ) )
1180               DEALLOCATE(globaldata_1d_sp)
1181            CASE( NF90_DOUBLE )
1182               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_dp ) )
1183               DEALLOCATE(globaldata_1d_dp)
1184            CASE DEFAULT   
1185               WRITE(numerr,*) '1d Unknown nf90 type: ', xtype
1186               STOP 4
1187         END SELECT
1188
1189      ELSEIF( ndims == 2 ) THEN 
1190     
1191         SELECT CASE( xtype )   
1192            CASE( NF90_BYTE )                   
1193               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_i1 ) )
1194               DEALLOCATE(globaldata_2d_i1)
1195            CASE( NF90_SHORT )                   
1196               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_i2 ) )
1197               DEALLOCATE(globaldata_2d_i2)
1198            CASE( NF90_INT )                             
1199               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_i4 ) )
1200               DEALLOCATE(globaldata_2d_i4)
1201            CASE( NF90_FLOAT )                             
1202               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_sp ) )
1203               DEALLOCATE(globaldata_2d_sp)
1204            CASE( NF90_DOUBLE )                                         
1205               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_dp ) )
1206               DEALLOCATE(globaldata_2d_dp)
1207            CASE DEFAULT   
1208               WRITE(numerr,*) '2d Unknown nf90 type: ', xtype
1209               STOP 4
1210         END SELECT     
1211                     
1212      ELSEIF( ndims == 3 ) THEN
1213     
1214         SELECT CASE( xtype ) 
1215            CASE( NF90_BYTE )                   
1216               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_i1 ) )
1217               DEALLOCATE(globaldata_3d_i1)
1218            CASE( NF90_SHORT )                   
1219               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_i2 ) )
1220               DEALLOCATE(globaldata_3d_i2)
1221            CASE( NF90_INT )                             
1222               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_i4 ) )
1223               DEALLOCATE(globaldata_3d_i4)
1224            CASE( NF90_FLOAT )                             
1225               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_sp ) )
1226               DEALLOCATE(globaldata_3d_sp)
1227            CASE( NF90_DOUBLE )                                         
1228               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_dp ) )
1229               DEALLOCATE(globaldata_3d_dp)
1230            CASE DEFAULT   
1231               WRITE(numerr,*) '3d Unknown nf90 type: ', xtype
1232               STOP 4
1233         END SELECT     
1234   
1235      ELSEIF( ndims == 4 ) THEN
1236
1237         SELECT CASE( xtype )   
1238            CASE( NF90_BYTE )                   
1239               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_i1, start=(/1,1,1,nt/) ) )
1240               DEALLOCATE(globaldata_4d_i1)
1241            CASE( NF90_SHORT )                   
1242               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_i2, start=(/1,1,1,nt/) ) )
1243               DEALLOCATE(globaldata_4d_i2)
1244            CASE( NF90_INT )                             
1245               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_i4, start=(/1,1,1,nt/) ) )
1246               DEALLOCATE(globaldata_4d_i4)
1247            CASE( NF90_FLOAT )                             
1248               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_sp, start=(/1,1,1,nt/) ) )
1249               DEALLOCATE(globaldata_4d_sp)
1250            CASE( NF90_DOUBLE )                                         
1251               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_dp, start=(/1,1,1,nt/) ) )
1252               DEALLOCATE(globaldata_4d_dp)
1253            CASE DEFAULT   
1254               WRITE(numerr,*) '4d Unknown nf90 type: ', xtype
1255               STOP 4
1256         END SELECT     
1257         ! why only for big data set, test the cost.
1258         CALL check_nf90( nf90_sync( outid ) )    ! flush buffers to disk after writing big 4D datasets
1259   
1260      ENDIF
1261
1262         nt = nt + ntslice
1263
1264      END DO ! WHILE loop
1265   
1266   END DO  ! loop over variables
1267
1268!---------------------------------------------------------------------------
1269!5. Close files
1270
1271!5.1 Close input files
1272
1273   IF (l_verbose) WRITE(numout,*) 'Closing input files...'
1274   DO ifile = 1, ndomain
1275      ncid = inncids(ifile)
1276      CALL check_nf90( nf90_close( ncid ) )
1277   END DO
1278
1279!5.2 Close output file
1280
1281   IF (l_verbose) WRITE(numout,*) 'Closing output file...'
1282   CALL check_nf90( nf90_close( outid ) )
1283
1284   IF (l_verbose) WRITE(numout,*) 'NEMO rebuild completed successfully'
1285   IF (l_verbose) WRITE(numout,*)
1286
1287CONTAINS
1288
1289
1290   SUBROUTINE check_nf90(status, errorFlag)
1291   !---------------------------------------------------------------------
1292   !  Checks return code from nf90 library calls and warns if needed
1293   !  If errorFlag is present then it just increments this flag (OMP use)
1294   !
1295   !---------------------------------------------------------------------
1296      INTEGER, INTENT(IN   ) :: status
1297      INTEGER, INTENT(INOUT), OPTIONAL :: errorFlag
1298   !---------------------------------------------------------------------
1299
1300      IF( status /= nf90_noerr ) THEN
1301         WRITE(numerr,*) 'ERROR! : '//TRIM(nf90_strerror(status))
1302         IF( PRESENT( errorFlag ) ) THEN
1303            errorFlag = errorFlag + status
1304         ELSE
1305            WRITE(numerr,*) "*** NEMO rebuild failed ***"
1306            WRITE(numerr,*)
1307            STOP 5
1308         ENDIF
1309      ENDIF
1310
1311   END SUBROUTINE check_nf90
1312
1313
1314END PROGRAM rebuild_nemo
Note: See TracBrowser for help on using the repository browser.