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 @ 10131

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

add option to manage compression in rebuild_nemo, add compression and chunking in rebuild_nemo.F90

File size: 62.0 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
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           if( TRIM(indimnames(dimids(idim))) == 'x' )             &
371    &                             chunksizes(idim) = min(outdimlens(dimids(idim)), max(nc4_xchunk,1))
372           if( TRIM(indimnames(dimids(idim))) == 'y' )             &
373    &                             chunksizes(idim) = min(outdimlens(dimids(idim)), max(nc4_ychunk,1))
374! trick to find var in a list of suggestion (var0 and var1 : INDEX(|var0|var1|,|var|)
375           if( INDEX('|depth|z|nav_lev|','|'//TRIM(indimnames(dimids(idim)))//'|') > 0 ) &
376    &                             chunksizes(idim) = min(outdimlens(dimids(idim)), max(nc4_zchunk,1))
377           if( INDEX('|time|time_counter|','|'//TRIM(indimnames(dimids(idim)))//'|') > 0 ) &
378    &                             chunksizes(idim) = min(outdimlens(dimids(idim)), max(nc4_tchunk,1))
379        END DO
380#if defined key_netcdf4
381        CALL check_nf90( nf90_def_var( outid, varname, xtype, outdimids, varid, &
382                                       deflate_level=deflate_level ) )
383        IF (l_verbose) WRITE(numout,*) 'Dims    : ',ndims, outdimids(1:ndims)
384        IF (l_verbose) WRITE(numout,*) 'names   : ',(TRIM(indimnames(dimids(idim)))//' ',idim=1,ndims)
385        IF (l_verbose) WRITE(numout,*) 'lens   : ',(outdimlens(dimids(idim)),idim=1,ndims)
386        IF (l_verbose) WRITE(numout,*) 'Chunking: ',chunkalg, chunksizes
387        CALL check_nf90( nf90_def_var_chunking( outid, varid, chunkalg, &
388   &                                 chunksizes ) )
389      ELSE
390        CALL check_nf90( nf90_def_var( outid, varname, xtype, outdimids, varid ) )
391#else
392      CALL check_nf90( nf90_def_var( outid, varname, xtype, outdimids, varid ) )
393#endif
394      ENDIF
395      DEALLOCATE(outdimids)
396      DEALLOCATE(chunksizes)
397      IF (l_verbose) WRITE(numout,*) 'Defining variable '//TRIM(varname)//'...' 
398      IF( natts > 0 ) THEN
399         DO attid = 1, natts
400            CALL check_nf90( nf90_inq_attname( ncid, varid, attid, attname ) )
401            IF ( attname == "_FillValue" ) THEN
402               CALL check_nf90( nf90_get_att( ncid, varid, attname, rmdi ) )
403               mdiVals(jv)=rmdi
404            ENDIF
405            CALL check_nf90( nf90_copy_att( ncid, varid, attname, outid, varid ) )
406         END DO
407      ENDIF
408   END DO
409
410!2.3 End definitions in output file and copy 1st file ncid to the inncids array
411
412   CALL check_nf90( nf90_enddef( outid ) )
413   inncids(1) = ncid
414   IF (l_verbose) WRITE(numout,*) 'Finished defining output file.'
415 
416!---------------------------------------------------------------------------
417!3. Read in data from each file for each variable
418
419!3.1 Open each file and store the ncid in inncids array
420
421   IF (l_verbose) WRITE(numout,*) 'Opening input files...'
422
423   ! Set a file chunk cache size for the processor-domain files that scales with the number of processors
424   patchchunk = max(8192, fchunksize/ndomain)
425
426   ! open files
427   DO ifile = 2, ndomain
428      CALL check_nf90( nf90_open( TRIM(filenames(ifile)), nf90_share, ncid, chunksize=patchchunk ) )
429      inncids(ifile) = ncid
430   END DO
431   IF (l_verbose) WRITE(numout,*) 'All input files open.'
432
433   DO jv = 1, nvars
434
435      ValMin = 1.e10
436      ValMax = -1.e10
437      l_valid = .false.
438      istop = nf90_noerr
439      nt = 1
440      ntslice = nmax_unlimited
441      IF( nslicesize == 0 ) nslicesize = nmax_unlimited
442
443!3.2 Inquire variable to find out name and how many dimensions it has
444!    and importantly whether it contains the dimensions in rebuild_dims()
445
446      ncid = inncids(1)
447      CALL check_nf90( nf90_inquire_variable( ncid, jv, varname, xtype, ndims, dimids, natts ) )
448
449      l_noRebuild = .true.
450      IF( ANY( dimids(1:ndims) == rebuild_dims(1) )) l_noRebuild = .false.
451      IF( rbdims > 1 ) THEN
452         IF( ANY( dimids(1:ndims) == rebuild_dims(2) )) l_noRebuild = .false.
453      ENDIF
454
455!3.2.0 start while loop for time slicing
456
457      DO WHILE( nt <= nmax_unlimited )
458
459         IF( ndims > 3 ) THEN
460            ntslice = MIN( nslicesize, nmax_unlimited + 1 - nt )
461         ENDIF
462
463      IF (l_noRebuild) THEN
464
465         IF( nslicesize == nmax_unlimited .OR. ndims <= 3 ) THEN
466            IF (l_verbose) WRITE(numout,*) 'Copying data from variable '//TRIM(varname)//'...'
467         ELSE
468            IF (l_verbose) WRITE(numout,'(A,I3,A,I3,A)') ' Copying data from variable '  &
469            &                 //TRIM(varname)//' for slices ',nt,' to ',nt+ntslice-1,' ...'
470         ENDIF
471
472!3.2.1 If rebuilding not required then just need to read in variable
473!      for copying direct into output file after the OMP (files) loop.
474         IF( ndims == 0 ) THEN
475
476            SELECT CASE( xtype )
477               CASE( NF90_BYTE )
478                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_i1 ) )
479               CASE( NF90_SHORT )
480                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_i2 ) )
481               CASE( NF90_INT )
482                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_i4 ) )
483               CASE( NF90_FLOAT )
484                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_sp ) )
485               CASE( NF90_DOUBLE )
486                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_dp ) )
487               CASE DEFAULT
488                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
489                  STOP 4
490            END SELECT
491
492         ELSEIF( ndims == 1 ) THEN
493
494            SELECT CASE( xtype )
495               CASE( NF90_BYTE )
496                  ALLOCATE(globaldata_1d_i1(indimlens(dimids(1))))
497                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_i1 ) )
498               CASE( NF90_SHORT )
499                  ALLOCATE(globaldata_1d_i2(indimlens(dimids(1))))
500                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_i2 ) )
501               CASE( NF90_INT )
502                  ALLOCATE(globaldata_1d_i4(indimlens(dimids(1))))
503                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_i4 ) )
504               CASE( NF90_FLOAT )
505                  ALLOCATE(globaldata_1d_sp(indimlens(dimids(1))))
506                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_sp ) )
507               CASE( NF90_DOUBLE )
508                  ALLOCATE(globaldata_1d_dp(indimlens(dimids(1))))
509                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_dp ) )
510               CASE DEFAULT
511                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
512                  STOP 4
513            END SELECT
514
515         ELSEIF( ndims == 2 ) THEN
516
517            SELECT CASE( xtype )
518               CASE( NF90_BYTE )
519                  ALLOCATE(globaldata_2d_i1(indimlens(dimids(1)),indimlens(dimids(2))))
520                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_i1 ) )
521               CASE( NF90_SHORT )
522                  ALLOCATE(globaldata_2d_i2(indimlens(dimids(1)),indimlens(dimids(2))))
523                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_i2 ) )
524               CASE( NF90_INT )
525                  ALLOCATE(globaldata_2d_i4(indimlens(dimids(1)),indimlens(dimids(2))))
526                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_i4 ) )
527               CASE( NF90_FLOAT )
528                  ALLOCATE(globaldata_2d_sp(indimlens(dimids(1)),indimlens(dimids(2))))
529                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_sp ) )
530               CASE( NF90_DOUBLE )
531                  ALLOCATE(globaldata_2d_dp(indimlens(dimids(1)),indimlens(dimids(2))))
532                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_dp ) )
533               CASE DEFAULT
534                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
535                  STOP 4
536            END SELECT
537
538         ELSEIF( ndims == 3 ) THEN
539
540            SELECT CASE( xtype )
541               CASE( NF90_BYTE )
542                  ALLOCATE(globaldata_3d_i1(indimlens(dimids(1)),indimlens(dimids(2)),       &
543                     &                      indimlens(dimids(3))))
544                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_i1 ) )
545               CASE( NF90_SHORT )
546                  ALLOCATE(globaldata_3d_i2(indimlens(dimids(1)),indimlens(dimids(2)),       &
547                     &                      indimlens(dimids(3))))
548                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_i2 ) )
549               CASE( NF90_INT )
550                  ALLOCATE(globaldata_3d_i4(indimlens(dimids(1)),indimlens(dimids(2)),       &
551                     &                      indimlens(dimids(3))))
552                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_i4 ) )
553               CASE( NF90_FLOAT )
554                  ALLOCATE(globaldata_3d_sp(indimlens(dimids(1)),indimlens(dimids(2)),       &
555                     &                      indimlens(dimids(3))))
556                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_sp ) )
557               CASE( NF90_DOUBLE )
558                  ALLOCATE(globaldata_3d_dp(indimlens(dimids(1)),indimlens(dimids(2)),       &
559                     &                      indimlens(dimids(3))))
560                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_dp ) )
561               CASE DEFAULT
562                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
563                  STOP 4
564            END SELECT
565
566         ELSEIF( ndims == 4 ) THEN
567
568            SELECT CASE( xtype )
569               CASE( NF90_BYTE )
570                  ALLOCATE(globaldata_4d_i1(indimlens(dimids(1)),indimlens(dimids(2)),       &
571                     &                      indimlens(dimids(3)),ntslice))
572                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_i1, start=(/1,1,1,nt/) ) )
573               CASE( NF90_SHORT )
574                  ALLOCATE(globaldata_4d_i2(indimlens(dimids(1)),indimlens(dimids(2)),       &
575                     &                      indimlens(dimids(3)),ntslice))
576                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_i2, start=(/1,1,1,nt/) ) )
577               CASE( NF90_INT )
578                  ALLOCATE(globaldata_4d_i4(indimlens(dimids(1)),indimlens(dimids(2)),       &
579                     &                      indimlens(dimids(3)),ntslice))
580                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_i4, start=(/1,1,1,nt/) ) )
581               CASE( NF90_FLOAT )
582                  ALLOCATE(globaldata_4d_sp(indimlens(dimids(1)),indimlens(dimids(2)),       &
583                     &                      indimlens(dimids(3)),ntslice))
584                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_sp, start=(/1,1,1,nt/) ) )
585               CASE( NF90_DOUBLE )
586                  ALLOCATE(globaldata_4d_dp(indimlens(dimids(1)),indimlens(dimids(2)),       &
587                     &                      indimlens(dimids(3)),ntslice))
588                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_dp, start=(/1,1,1,nt/) ) )
589               CASE DEFAULT
590                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
591                  STOP 4
592            END SELECT
593
594         ENDIF
595
596      ELSE  ! l_noRebuild = .false.
597
598!3.2.2 For variables that require rebuilding we need to read in from all ndomain files
599!      Here we allocate global variables ahead of looping over files
600         IF( nslicesize == nmax_unlimited .OR. ndims <= 3 ) THEN
601            IF (l_verbose) WRITE(numout,*) 'Rebuilding data from variable '//TRIM(varname)//'...'
602         ELSE
603            IF (l_verbose) WRITE(numout,'(A,I3,A,I3,A)') ' Rebuilding data from variable '  &
604            &                 //TRIM(varname)//' for slices ',nt,' to ',nt+ntslice-1,' ...'
605         ENDIF
606         IF( ndims == 1 ) THEN
607
608            SELECT CASE( xtype )
609               CASE( NF90_BYTE )
610                  ALLOCATE(globaldata_1d_i1(outdimlens(dimids(1))))
611                  IF (l_maskout) globaldata_1d_i1(:)=mdiVals(jv)
612               CASE( NF90_SHORT )
613                  ALLOCATE(globaldata_1d_i2(outdimlens(dimids(1))))
614                  IF (l_maskout) globaldata_1d_i2(:)=mdiVals(jv)
615               CASE( NF90_INT )
616                  ALLOCATE(globaldata_1d_i4(outdimlens(dimids(1))))
617                  IF (l_maskout) globaldata_1d_i4(:)=mdiVals(jv)
618               CASE( NF90_FLOAT )
619                  ALLOCATE(globaldata_1d_sp(outdimlens(dimids(1))))
620                  IF (l_maskout) globaldata_1d_sp(:)=mdiVals(jv)
621               CASE( NF90_DOUBLE )
622                  ALLOCATE(globaldata_1d_dp(outdimlens(dimids(1))))
623                  IF (l_maskout) globaldata_1d_dp(:)=mdiVals(jv)
624               CASE DEFAULT
625                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
626                  STOP 4
627            END SELECT
628
629         ELSEIF( ndims == 2 ) THEN
630
631            SELECT CASE( xtype )
632               CASE( NF90_BYTE )
633                  ALLOCATE(globaldata_2d_i1(outdimlens(dimids(1)),outdimlens(dimids(2))))
634                  IF (l_maskout) globaldata_2d_i1(:,:)=mdiVals(jv)
635               CASE( NF90_SHORT )
636                  ALLOCATE(globaldata_2d_i2(outdimlens(dimids(1)),outdimlens(dimids(2))))
637                  IF (l_maskout) globaldata_2d_i2(:,:)=mdiVals(jv)
638               CASE( NF90_INT )
639                  ALLOCATE(globaldata_2d_i4(outdimlens(dimids(1)),outdimlens(dimids(2))))
640                  IF (l_maskout) globaldata_2d_i4(:,:)=mdiVals(jv)
641               CASE( NF90_FLOAT )
642                  ALLOCATE(globaldata_2d_sp(outdimlens(dimids(1)),outdimlens(dimids(2))))
643                  IF (l_maskout) globaldata_2d_sp(:,:)=mdiVals(jv)
644               CASE( NF90_DOUBLE )
645                  ALLOCATE(globaldata_2d_dp(outdimlens(dimids(1)),outdimlens(dimids(2))))
646                  IF (l_maskout) globaldata_2d_dp(:,:)=mdiVals(jv)
647               CASE DEFAULT
648                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
649                  STOP 4
650            END SELECT
651
652         ELSEIF( ndims == 3 ) THEN
653
654            SELECT CASE( xtype )
655               CASE( NF90_BYTE )
656                  ALLOCATE(globaldata_3d_i1(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
657                     &                      outdimlens(dimids(3))))
658                  IF (l_maskout) globaldata_3d_i1(:,:,:)=mdiVals(jv)
659               CASE( NF90_SHORT )
660                  ALLOCATE(globaldata_3d_i2(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
661                     &                      outdimlens(dimids(3))))
662                  IF (l_maskout) globaldata_3d_i2(:,:,:)=mdiVals(jv)
663               CASE( NF90_INT )
664                  ALLOCATE(globaldata_3d_i4(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
665                     &                      outdimlens(dimids(3))))
666                  IF (l_maskout) globaldata_3d_i4(:,:,:)=mdiVals(jv)
667               CASE( NF90_FLOAT )
668                  ALLOCATE(globaldata_3d_sp(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
669                     &                      outdimlens(dimids(3))))
670                  IF (l_maskout) globaldata_3d_sp(:,:,:)=mdiVals(jv)
671               CASE( NF90_DOUBLE )
672                  ALLOCATE(globaldata_3d_dp(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
673                     &                      outdimlens(dimids(3))))
674                  IF (l_maskout) globaldata_3d_dp(:,:,:)=mdiVals(jv)
675               CASE DEFAULT
676                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
677                  STOP 4
678            END SELECT
679
680         ELSEIF( ndims == 4 ) THEN
681
682            SELECT CASE( xtype )
683               CASE( NF90_BYTE )
684                  ALLOCATE(globaldata_4d_i1(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
685                     &                      outdimlens(dimids(3)),ntslice))
686                  IF (l_maskout) globaldata_4d_i1(:,:,:,:)=mdiVals(jv)
687               CASE( NF90_SHORT )
688                  ALLOCATE(globaldata_4d_i2(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
689                     &                      outdimlens(dimids(3)),ntslice))
690                  IF (l_maskout) globaldata_4d_i2(:,:,:,:)=mdiVals(jv)
691               CASE( NF90_INT )
692                  ALLOCATE(globaldata_4d_i4(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
693                     &                      outdimlens(dimids(3)),ntslice))
694                  IF (l_maskout) globaldata_4d_i4(:,:,:,:)=mdiVals(jv)
695               CASE( NF90_FLOAT )
696                  ALLOCATE(globaldata_4d_sp(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
697                     &                      outdimlens(dimids(3)),ntslice))
698                  IF (l_maskout) globaldata_4d_sp(:,:,:,:)=mdiVals(jv)
699               CASE( NF90_DOUBLE )
700                  ALLOCATE(globaldata_4d_dp(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
701                     &                      outdimlens(dimids(3)),ntslice))
702                  IF (l_maskout) globaldata_4d_dp(:,:,:,:)=mdiVals(jv)
703               CASE DEFAULT
704                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
705                  STOP 4
706            END SELECT
707         ELSE
708            WRITE(numerr,*) 'ERROR! : A netcdf variable has more than 4 dimensions which is not taken into account'
709            STOP 4
710         ENDIF
711
712!$OMP  PARALLEL DO DEFAULT(NONE)                                                          &
713!$OMP& PRIVATE(ifile,ncid,xtype,start_pos,local_sizes,InMin,InMax,natts,                  &
714!$OMP&         ndims,attid,attname,dimids,idim,dimname,dimlen,unlimitedDimId,             &
715!$OMP&         halo_start,halo_end,idomain,jdomain,rdomain,di,dj,dr,                      &
716!$OMP&         localdata_1d_i2,localdata_1d_i4,localdata_1d_sp,localdata_1d_dp,           &
717!$OMP&         localdata_2d_i2,localdata_2d_i4,localdata_2d_sp,localdata_2d_dp,           &
718!$OMP&         localdata_3d_i2,localdata_3d_i4,localdata_3d_sp,localdata_3d_dp,           &
719!$OMP&         localdata_4d_i2,localdata_4d_i4,localdata_4d_sp,localdata_4d_dp,           &
720!$OMP&         localdata_1d_i1,localdata_2d_i1,localdata_3d_i1,localdata_4d_i1)           &
721!$OMP& SHARED(jv,nvars,varname,filenames,ValMin,ValMax,indimlens,outdimlens,rbdims,       &
722!$OMP&        ndomain,outid,fchunksize,istop,l_valid,nthreads,inncids,rebuild_dims,       &
723!$OMP&        globaldata_1d_i2,globaldata_1d_i4,globaldata_1d_sp,globaldata_1d_dp,        &
724!$OMP&        globaldata_2d_i2,globaldata_2d_i4,globaldata_2d_sp,globaldata_2d_dp,        &
725!$OMP&        globaldata_3d_i2,globaldata_3d_i4,globaldata_3d_sp,globaldata_3d_dp,        &
726!$OMP&        globaldata_4d_i2,globaldata_4d_i4,globaldata_4d_sp,globaldata_4d_dp,        &
727!$OMP&        globaldata_1d_i1,globaldata_2d_i1,globaldata_3d_i1,globaldata_4d_i1,        &
728!$OMP&        ntslice,nt,nmax_unlimited,indimnames,dims,patchchunk)
729
730         DO ifile = 1, ndomain
731
732            ncid = inncids(ifile)
733!$OMP CRITICAL
734            CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_size_local', local_sizes ), istop )
735            CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_position_first', start_pos ), istop )
736            CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_halo_size_start', halo_start ), istop )
737            CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_halo_size_end', halo_end ), istop )
738            CALL check_nf90( nf90_inquire_variable( ncid, jv, varname, xtype, ndims, dimids, natts ), istop )
739            CALL check_nf90( nf90_inquire( ncid, unlimitedDimId = unlimitedDimId ), istop )
740!$OMP END CRITICAL
741
742            ! set defaults for rebuilding so that i is 1st, j 2nd
743            di=1
744            dj=2
745
746            IF( rbdims == 1 ) THEN
747               ! override defaults above and set other variables
748               start_pos(2) = 1
749               local_sizes(2) = outdimlens(3-dimids(2))
750               halo_end(2) = 0
751               halo_start(2) = 0
752               di=rebuild_dims(1)
753               dj=3-di
754            ENDIF
755
756!3.3.1 Generate local domain interior sizes from local_sizes and halo sizes
757!      idomain defines the 1st and last interior points in the i direction and
758!      jdomain defines the 1st and last interior points in the j direction
759
760            idomain(1) = 1 + halo_start(di)
761            idomain(2) = local_sizes(di) - halo_end(di)
762            jdomain(1) = 1 + halo_start(dj)
763            jdomain(2) = local_sizes(dj) - halo_end(dj)
764
765!3.3.2 For rbdims or more dimensions put the data array from this input file into the correct
766!      part of the output data array. Assume the first dimensions are those to be rebuilt.
767
768            IF( ndims == 1 ) THEN
769
770               IF( rebuild_dims(1) == 1 ) THEN
771                  dr = di
772                  rdomain = idomain
773               ELSE
774                  dr = dj
775                  rdomain = jdomain
776               ENDIF
777
778              SELECT CASE( xtype )
779                  CASE( NF90_BYTE )
780                     ALLOCATE(localdata_1d_i1(local_sizes(dr)))
781                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_i1 ), istop )
782                     DO jr = rdomain(1), rdomain(2)
783                        globaldata_1d_i1(start_pos(dr) + jr - 1) = localdata_1d_i1(jr)
784                     END DO
785                     DEALLOCATE(localdata_1d_i1)
786                  CASE( NF90_SHORT )
787                     ALLOCATE(localdata_1d_i2(local_sizes(dr)))
788                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_i2 ), istop )
789                     DO jr = rdomain(1), rdomain(2)
790                        globaldata_1d_i2(start_pos(dr) + jr - 1) = localdata_1d_i2(jr)
791                     END DO
792                     DEALLOCATE(localdata_1d_i2)
793                  CASE( NF90_INT )
794                     ALLOCATE(localdata_1d_i4(local_sizes(dr)))
795                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_i4 ), istop )
796                     DO jr = rdomain(1), rdomain(2)
797                        globaldata_1d_i4(start_pos(dr) + jr - 1) = localdata_1d_i4(jr)
798                     END DO
799                     DEALLOCATE(localdata_1d_i4)
800                  CASE( NF90_FLOAT )
801                     ALLOCATE(localdata_1d_sp(local_sizes(dr)))
802                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_sp ), istop )
803                     DO jr = rdomain(1), rdomain(2)
804                        globaldata_1d_sp(start_pos(dr) + jr - 1) = localdata_1d_sp(jr)
805                     END DO
806                     DEALLOCATE(localdata_1d_sp)
807                  CASE( NF90_DOUBLE )
808                     ALLOCATE(localdata_1d_dp(local_sizes(dr)))
809                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_dp ), istop )
810                     DO jr = rdomain(1), rdomain(2)
811                        globaldata_1d_dp(start_pos(dr) + jr - 1) = localdata_1d_dp(jr)
812                     END DO
813                     DEALLOCATE(localdata_1d_dp)
814                  CASE DEFAULT
815                     WRITE(numerr,*) 'Unknown nf90 type: ', xtype
816                     istop = istop + 1
817               END SELECT
818
819            ELSEIF( ndims == 2 ) THEN
820       
821               SELECT CASE( xtype )
822                  CASE( NF90_BYTE )
823                     ALLOCATE(localdata_2d_i1(local_sizes(di),local_sizes(dj)))
824                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_i1 ), istop )
825                        DO jj = jdomain(1), jdomain(2)
826                           DO ji = idomain(1), idomain(2)
827                              globaldata_2d_i1(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_i1(ji,jj)
828                        END DO
829                     END DO
830                     DEALLOCATE(localdata_2d_i1)
831                  CASE( NF90_SHORT )
832                     ALLOCATE(localdata_2d_i2(local_sizes(di),local_sizes(dj)))
833                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_i2 ), istop )
834                     DO jj = jdomain(1), jdomain(2)
835                        DO ji = idomain(1), idomain(2)
836                           globaldata_2d_i2(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_i2(ji,jj)
837                        END DO
838                     END DO
839                     DEALLOCATE(localdata_2d_i2)
840                  CASE( NF90_INT )
841                     ALLOCATE(localdata_2d_i4(local_sizes(di),local_sizes(dj)))
842                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_i4 ), istop )
843                     DO jj = jdomain(1), jdomain(2)
844                        DO ji = idomain(1), idomain(2)
845                           globaldata_2d_i4(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_i4(ji,jj)
846                        END DO
847                     END DO
848                     DEALLOCATE(localdata_2d_i4)
849                  CASE( NF90_FLOAT )
850                     ALLOCATE(localdata_2d_sp(local_sizes(di),local_sizes(dj)))
851                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_sp ), istop )
852                     DO jj = jdomain(1), jdomain(2)
853                        DO ji = idomain(1), idomain(2)
854                           globaldata_2d_sp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_sp(ji,jj)
855                        END DO
856                     END DO
857                     DEALLOCATE(localdata_2d_sp)
858                  CASE( NF90_DOUBLE )
859                     ALLOCATE(localdata_2d_dp(local_sizes(di),local_sizes(dj)))
860                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_dp ), istop )
861                     DO jj = jdomain(1), jdomain(2)
862                        DO ji = idomain(1), idomain(2)
863                           globaldata_2d_dp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_dp(ji,jj)
864                        END DO
865                     END DO
866                     DEALLOCATE(localdata_2d_dp) 
867                  CASE DEFAULT
868                     WRITE(numerr,*) 'Unknown nf90 type: ', xtype
869                     istop = istop + 1
870               END SELECT
871
872            ELSEIF( ndims == 3 ) THEN
873
874               SELECT CASE( xtype )
875                  CASE( NF90_BYTE )
876                     ALLOCATE(localdata_3d_i1(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
877                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_i1 ), istop )
878!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
879!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_i1,localdata_3d_i1,di,dj)
880                     DO jk = 1, indimlens(dimids(3))
881                        DO jj = jdomain(1), jdomain(2)
882                           DO ji = idomain(1), idomain(2)
883                              globaldata_3d_i1(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_i1(ji,jj,jk)
884                           END DO
885                        END DO
886                     END DO
887!$OMP END PARALLEL DO
888                     DEALLOCATE(localdata_3d_i1)
889                  CASE( NF90_SHORT )
890                     ALLOCATE(localdata_3d_i2(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
891                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_i2 ), istop )
892!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
893!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_i2,localdata_3d_i2,di,dj)
894                     DO jk = 1, indimlens(dimids(3))
895                        DO jj = jdomain(1), jdomain(2)
896                           DO ji = idomain(1), idomain(2)
897                              globaldata_3d_i2(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_i2(ji,jj,jk)
898                           END DO
899                        END DO
900                     END DO
901!$OMP END PARALLEL DO
902                     DEALLOCATE(localdata_3d_i2)
903                  CASE( NF90_INT )
904                     ALLOCATE(localdata_3d_i4(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
905                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_i4 ), istop )
906!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
907!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_i4,localdata_3d_i4,di,dj)
908                     DO jk = 1, indimlens(dimids(3))
909                        DO jj = jdomain(1), jdomain(2)
910                           DO ji = idomain(1), idomain(2)
911                              globaldata_3d_i4(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_i4(ji,jj,jk)
912                           END DO
913                        END DO
914                     END DO
915!$OMP END PARALLEL DO
916                     DEALLOCATE(localdata_3d_i4)
917                  CASE( NF90_FLOAT )
918                     ! TG: This if statement is added to check if the 1st dimension is the corners (for lon_bounds) variables
919                     ! TG: Had to add the unsatisfactory check for 'lon' as it failed for diaptr files
920                     ! TG: Would like to find a better assumption for this.
921                     IF ( trim(indimnames(dimids(1))) /= dims(1) .AND. indimnames(dimids(1)) .NE. 'lon' ) THEN
922                        ALLOCATE(localdata_3d_sp(indimlens(dimids(1)),local_sizes(di),local_sizes(dj)))
923                        WRITE(*,*) 'test', ifile, jv, indimlens(dimids(1)),local_sizes(di),local_sizes(dj)
924                        CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_sp ), istop )
925                        WRITE(*,*) 'test2'
926                        DO jj = jdomain(1), jdomain(2)
927                           DO ji = idomain(1), idomain(2)
928                              DO jk = 1, indimlens(dimids(1))
929                                 globaldata_3d_sp(jk, start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_3d_sp(jk,ji,jj)
930                              END DO
931                           END DO
932                        END DO
933                     ELSE
934                        ALLOCATE(localdata_3d_sp(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
935                        CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_sp ), istop ) 
936!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
937!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_sp,localdata_3d_sp,di,dj)
938                        DO jk = 1, indimlens(dimids(3))
939                           DO jj = jdomain(1), jdomain(2)
940                              DO ji = idomain(1), idomain(2)
941                                 globaldata_3d_sp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_sp(ji,jj,jk)
942                              END DO
943                           END DO
944                        END DO
945!$OMP END PARALLEL DO
946                     ENDIF
947                     DEALLOCATE(localdata_3d_sp)
948                  CASE( NF90_DOUBLE )
949                     IF ( trim(indimnames(dimids(1))) /= dims(1) ) THEN
950                        ALLOCATE(localdata_3d_dp(indimlens(dimids(1)),local_sizes(di),local_sizes(dj)))
951                        CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_dp ), istop ) 
952                        DO jj = jdomain(1), jdomain(2)
953                           DO ji = idomain(1), idomain(2)
954                              DO jk = 1, indimlens(dimids(1))
955                                 globaldata_3d_dp(jk, start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_3d_dp(jk,ji,jj)
956                              END DO
957                           END DO
958                        END DO
959                     ELSE
960                        ALLOCATE(localdata_3d_dp(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
961                        CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_dp ), istop ) 
962!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
963!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_dp,localdata_3d_dp,di,dj)
964                        DO jk = 1, indimlens(dimids(3))
965                           DO jj = jdomain(1), jdomain(2)
966                              DO ji = idomain(1), idomain(2)
967                                 globaldata_3d_dp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_dp(ji,jj,jk)
968                              END DO
969                           END DO
970                        END DO
971!$OMP END PARALLEL DO
972                     ENDIF
973                     DEALLOCATE(localdata_3d_dp)
974                  CASE DEFAULT
975                     WRITE(numerr,*) 'Unknown nf90 type: ', xtype
976                     istop = istop + 1
977               END SELECT
978     
979            ELSEIF (ndims == 4) THEN
980       
981               SELECT CASE( xtype )
982                  CASE( NF90_BYTE )
983                     ALLOCATE(localdata_4d_i1(local_sizes(di),local_sizes(dj),               &
984                         &                     indimlens(dimids(3)),ntslice))
985                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_i1, start=(/1,1,1,nt/) ), istop )
986!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
987!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_i1,localdata_4d_i1,di,dj,nt,ntslice)
988                     DO jl = 1, ntslice
989!$OMP DO
990                        DO jk = 1, indimlens(dimids(3))
991                           DO jj = jdomain(1), jdomain(2)
992                              DO ji = idomain(1), idomain(2)
993                                 globaldata_4d_i1(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_i1(ji,jj,jk,jl)
994                              END DO
995                           END DO
996                        END DO
997!$OMP END DO nowait
998                     END DO
999!$OMP END PARALLEL
1000                     DEALLOCATE(localdata_4d_i1)
1001                  CASE( NF90_SHORT )
1002                     ALLOCATE(localdata_4d_i2(local_sizes(di),local_sizes(dj),               &
1003                        &                     indimlens(dimids(3)),ntslice))
1004                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_i2, start=(/1,1,1,nt/) ), istop )
1005!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
1006!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_i2,localdata_4d_i2,di,dj,nt,ntslice)
1007                     DO jl = 1, ntslice
1008!$OMP DO
1009                        DO jk = 1, indimlens(dimids(3))
1010                           DO jj = jdomain(1), jdomain(2) 
1011                              DO ji = idomain(1), idomain(2)
1012                                 globaldata_4d_i2(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_i2(ji,jj,jk,jl)
1013                              END DO
1014                           END DO
1015                        END DO
1016!$OMP END DO nowait
1017                     END DO
1018!$OMP END PARALLEL
1019                     DEALLOCATE(localdata_4d_i2)
1020                  CASE( NF90_INT )
1021                     ALLOCATE(localdata_4d_i4(local_sizes(di),local_sizes(dj),               &
1022                        &                     indimlens(dimids(3)),ntslice))
1023                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_i4, start=(/1,1,1,nt/) ), istop )
1024!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
1025!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_i4,localdata_4d_i4,di,dj,nt,ntslice)
1026                     DO jl = 1, ntslice
1027!$OMP DO
1028                        DO jk = 1, indimlens(dimids(3))
1029                           DO jj = jdomain(1), jdomain(2) 
1030                              DO ji = idomain(1), idomain(2)
1031                                 globaldata_4d_i4(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_i4(ji,jj,jk,jl)
1032                              END DO
1033                           END DO
1034                        END DO
1035!$OMP END DO nowait
1036                     END DO
1037!$OMP END PARALLEL
1038                     DEALLOCATE(localdata_4d_i4)
1039                  CASE( NF90_FLOAT )
1040                     ALLOCATE(localdata_4d_sp(local_sizes(di),local_sizes(dj),               &
1041                        &                     indimlens(dimids(3)),ntslice))
1042                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_sp, start=(/1,1,1,nt/) ), istop )
1043!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
1044!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_sp,localdata_4d_sp,di,dj,nt,ntslice)
1045                     DO jl = 1, ntslice
1046!$OMP DO
1047                        DO jk = 1, indimlens(dimids(3))
1048                           DO jj = jdomain(1), jdomain(2) 
1049                              DO ji = idomain(1), idomain(2)
1050                                 globaldata_4d_sp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_sp(ji,jj,jk,jl)
1051                              END DO
1052                           END DO
1053                        END DO
1054!$OMP END DO nowait
1055                     END DO
1056!$OMP END PARALLEL
1057                     DEALLOCATE(localdata_4d_sp)
1058                  CASE( NF90_DOUBLE )
1059                     ALLOCATE(localdata_4d_dp(local_sizes(di),local_sizes(dj),               &
1060                        &                     indimlens(dimids(3)),ntslice))
1061                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_dp, start=(/1,1,1,nt/) ), istop )
1062!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
1063!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_dp,localdata_4d_dp,di,dj,nt,ntslice)
1064                     DO jl = 1, ntslice
1065!$OMP DO
1066                        DO jk = 1, indimlens(dimids(3))
1067                           DO jj = jdomain(1), jdomain(2) 
1068                              DO ji = idomain(1), idomain(2)
1069                                 globaldata_4d_dp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_dp(ji,jj,jk,jl)
1070                              END DO
1071                           END DO
1072                        END DO
1073!$OMP END DO nowait
1074                     END DO
1075!$OMP END PARALLEL
1076                     DEALLOCATE(localdata_4d_dp)
1077                  CASE DEFAULT
1078                     WRITE(numerr,*) 'Unknown nf90 type: ', xtype
1079                     istop = istop + 1
1080               END SELECT
1081
1082            ENDIF ! l_noRebuild false
1083
1084!3.4 Work out if the valid_min and valid_max attributes exist for this variable.
1085!    If they do then calculate the extrema over all input files.
1086
1087            DO attid = 1, natts
1088               CALL check_nf90( nf90_inq_attname( ncid, jv, attid, attname ), istop )
1089               IF( INDEX( attname, "valid_min" ) == 1 ) THEN
1090                  CALL check_nf90( nf90_get_att( ncid, jv, attname, InMin), istop )
1091                  l_valid = .true.
1092               ENDIF
1093               IF( INDEX( attname, "valid_max" ) == 1 ) THEN
1094                  CALL check_nf90( nf90_get_att( ncid, jv, attname, InMax ), istop )
1095                  l_valid = .true.
1096               ENDIF
1097            END DO
1098
1099            IF (l_valid) THEN
1100!$OMP CRITICAL
1101               IF( InMin < ValMin ) ValMin = InMin
1102               IF( InMax > ValMax ) ValMax = InMax
1103!$OMP END CRITICAL
1104            ENDIF
1105
1106!3.5 Abort if failure and only 1 thread
1107
1108            IF( nthreads == 1 .AND. istop /= nf90_noerr )  THEN
1109               WRITE(numerr,*) '*** NEMO rebuild failed! ***'
1110               STOP 5
1111            ENDIF
1112       
1113         END DO  ! loop over files
1114!$OMP END PARALLEL DO
1115
1116!3.6 Abort if any of the OMP threads failed
1117         IF( istop /= nf90_noerr )  THEN
1118            WRITE(numerr,*) '*** NEMO rebuild failed! ***'
1119            STOP 5
1120         ENDIF
1121
1122      ENDIF ! ndims > 2
1123
1124!---------------------------------------------------------------------------
1125!4. Write data to output file
1126
1127      IF (l_verbose) WRITE(numout,*) 'Writing variable '//TRIM(varname)//'...'
1128
1129!4.1 If the valid min and max attributes exist then update them in the file
1130
1131      IF( l_valid ) THEN
1132         CALL check_nf90( nf90_put_att( outid, jv, "valid_min", ValMin ) )
1133         CALL check_nf90( nf90_put_att( outid, jv, "valid_max", ValMax ) )
1134      ENDIF
1135
1136!4.2 Write the data to the output file depending on how many dimensions it has
1137
1138      IF( ndims == 0 ) THEN
1139
1140         SELECT CASE( xtype )
1141            CASE( NF90_BYTE )
1142               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_i1 ) )
1143            CASE( NF90_SHORT )
1144               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_i2 ) )
1145            CASE( NF90_INT )
1146               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_i4 ) )
1147            CASE( NF90_FLOAT )
1148               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_sp ) )
1149            CASE( NF90_DOUBLE )
1150               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_dp ) )
1151            CASE DEFAULT   
1152               WRITE(numerr,*) '0d Unknown nf90 type: ', xtype
1153               STOP 4
1154         END SELECT
1155
1156      ELSEIF( ndims == 1 ) THEN
1157
1158         SELECT CASE( xtype )
1159            CASE( NF90_BYTE )
1160               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_i1 ) )
1161               DEALLOCATE(globaldata_1d_i1)
1162            CASE( NF90_SHORT )
1163               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_i2 ) )
1164               DEALLOCATE(globaldata_1d_i2)
1165            CASE( NF90_INT )
1166               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_i4 ) )
1167               DEALLOCATE(globaldata_1d_i4)
1168            CASE( NF90_FLOAT )
1169               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_sp ) )
1170               DEALLOCATE(globaldata_1d_sp)
1171            CASE( NF90_DOUBLE )
1172               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_dp ) )
1173               DEALLOCATE(globaldata_1d_dp)
1174            CASE DEFAULT   
1175               WRITE(numerr,*) '1d Unknown nf90 type: ', xtype
1176               STOP 4
1177         END SELECT
1178
1179      ELSEIF( ndims == 2 ) THEN 
1180     
1181         SELECT CASE( xtype )   
1182            CASE( NF90_BYTE )                   
1183               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_i1 ) )
1184               DEALLOCATE(globaldata_2d_i1)
1185            CASE( NF90_SHORT )                   
1186               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_i2 ) )
1187               DEALLOCATE(globaldata_2d_i2)
1188            CASE( NF90_INT )                             
1189               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_i4 ) )
1190               DEALLOCATE(globaldata_2d_i4)
1191            CASE( NF90_FLOAT )                             
1192               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_sp ) )
1193               DEALLOCATE(globaldata_2d_sp)
1194            CASE( NF90_DOUBLE )                                         
1195               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_dp ) )
1196               DEALLOCATE(globaldata_2d_dp)
1197            CASE DEFAULT   
1198               WRITE(numerr,*) '2d Unknown nf90 type: ', xtype
1199               STOP 4
1200         END SELECT     
1201                     
1202      ELSEIF( ndims == 3 ) THEN
1203     
1204         SELECT CASE( xtype ) 
1205            CASE( NF90_BYTE )                   
1206               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_i1 ) )
1207               DEALLOCATE(globaldata_3d_i1)
1208            CASE( NF90_SHORT )                   
1209               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_i2 ) )
1210               DEALLOCATE(globaldata_3d_i2)
1211            CASE( NF90_INT )                             
1212               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_i4 ) )
1213               DEALLOCATE(globaldata_3d_i4)
1214            CASE( NF90_FLOAT )                             
1215               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_sp ) )
1216               DEALLOCATE(globaldata_3d_sp)
1217            CASE( NF90_DOUBLE )                                         
1218               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_dp ) )
1219               DEALLOCATE(globaldata_3d_dp)
1220            CASE DEFAULT   
1221               WRITE(numerr,*) '3d Unknown nf90 type: ', xtype
1222               STOP 4
1223         END SELECT     
1224   
1225      ELSEIF( ndims == 4 ) THEN
1226
1227         SELECT CASE( xtype )   
1228            CASE( NF90_BYTE )                   
1229               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_i1, start=(/1,1,1,nt/) ) )
1230               DEALLOCATE(globaldata_4d_i1)
1231            CASE( NF90_SHORT )                   
1232               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_i2, start=(/1,1,1,nt/) ) )
1233               DEALLOCATE(globaldata_4d_i2)
1234            CASE( NF90_INT )                             
1235               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_i4, start=(/1,1,1,nt/) ) )
1236               DEALLOCATE(globaldata_4d_i4)
1237            CASE( NF90_FLOAT )                             
1238               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_sp, start=(/1,1,1,nt/) ) )
1239               DEALLOCATE(globaldata_4d_sp)
1240            CASE( NF90_DOUBLE )                                         
1241               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_dp, start=(/1,1,1,nt/) ) )
1242               DEALLOCATE(globaldata_4d_dp)
1243            CASE DEFAULT   
1244               WRITE(numerr,*) '4d Unknown nf90 type: ', xtype
1245               STOP 4
1246         END SELECT     
1247         ! why only for big data set, test the cost.
1248         CALL check_nf90( nf90_sync( outid ) )    ! flush buffers to disk after writing big 4D datasets
1249   
1250      ENDIF
1251
1252         nt = nt + ntslice
1253
1254      END DO ! WHILE loop
1255   
1256   END DO  ! loop over variables
1257
1258!---------------------------------------------------------------------------
1259!5. Close files
1260
1261!5.1 Close input files
1262
1263   IF (l_verbose) WRITE(numout,*) 'Closing input files...'
1264   DO ifile = 1, ndomain
1265      ncid = inncids(ifile)
1266      CALL check_nf90( nf90_close( ncid ) )
1267   END DO
1268
1269!5.2 Close output file
1270
1271   IF (l_verbose) WRITE(numout,*) 'Closing output file...'
1272   CALL check_nf90( nf90_close( outid ) )
1273
1274   IF (l_verbose) WRITE(numout,*) 'NEMO rebuild completed successfully'
1275   IF (l_verbose) WRITE(numout,*)
1276
1277CONTAINS
1278
1279
1280   SUBROUTINE check_nf90(status, errorFlag)
1281   !---------------------------------------------------------------------
1282   !  Checks return code from nf90 library calls and warns if needed
1283   !  If errorFlag is present then it just increments this flag (OMP use)
1284   !
1285   !---------------------------------------------------------------------
1286      INTEGER, INTENT(IN   ) :: status
1287      INTEGER, INTENT(INOUT), OPTIONAL :: errorFlag
1288   !---------------------------------------------------------------------
1289
1290      IF( status /= nf90_noerr ) THEN
1291         WRITE(numerr,*) 'ERROR! : '//TRIM(nf90_strerror(status))
1292         IF( PRESENT( errorFlag ) ) THEN
1293            errorFlag = errorFlag + status
1294         ELSE
1295            WRITE(numerr,*) "*** NEMO rebuild failed ***"
1296            WRITE(numerr,*)
1297            STOP 5
1298         ENDIF
1299      ENDIF
1300
1301   END SUBROUTINE check_nf90
1302
1303
1304END PROGRAM rebuild_nemo
Note: See TracBrowser for help on using the repository browser.