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 utils/tools/REBUILD_NEMO/src – NEMO

source: utils/tools/REBUILD_NEMO/src/rebuild_nemo.F90 @ 10338

Last change on this file since 10338 was 10338, checked in by mathiot, 5 years ago

Add deflation/chunking to rebuild_nemo tool. Fix #2165

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