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

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

add extra depth dimensions names to deals with restart and output

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