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 branches/UKMO/dev_r5518_nemo2cice_prints/NEMOGCM/TOOLS/REBUILD_NEMO/src – NEMO

source: branches/UKMO/dev_r5518_nemo2cice_prints/NEMOGCM/TOOLS/REBUILD_NEMO/src/rebuild_nemo.f90 @ 9817

Last change on this file since 9817 was 9817, checked in by dancopsey, 6 years ago

Merged in GO6 package branch up to revision 8356.

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