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/icebergs_restart_single_file/NEMOGCM/TOOLS/REBUILD_NEMO/src – NEMO

source: branches/UKMO/icebergs_restart_single_file/NEMOGCM/TOOLS/REBUILD_NEMO/src/rebuild_nemo.f90 @ 6019

Last change on this file since 6019 was 6019, checked in by timgraham, 8 years ago

Reinstated svn keywords before upgrading to head of trunk

  • Property svn:keywords set to Id
File size: 53.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
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
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                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_i1 ) )
387               CASE( NF90_SHORT )
388                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_i2 ) )
389               CASE( NF90_INT )
390                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_i4 ) )
391               CASE( NF90_FLOAT )
392                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_sp ) )
393               CASE( NF90_DOUBLE )
394                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_dp ) )
395               CASE DEFAULT
396                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
397                  STOP
398            END SELECT
399
400         ELSEIF( ndims == 1 ) THEN
401
402            SELECT CASE( xtype )
403               CASE( NF90_BYTE )
404                  ALLOCATE(globaldata_1d_i1(indimlens(dimids(1))))
405                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_i1 ) )
406               CASE( NF90_SHORT )
407                  ALLOCATE(globaldata_1d_i2(indimlens(dimids(1))))
408                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_i2 ) )
409               CASE( NF90_INT )
410                  ALLOCATE(globaldata_1d_i4(indimlens(dimids(1))))
411                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_i4 ) )
412               CASE( NF90_FLOAT )
413                  ALLOCATE(globaldata_1d_sp(indimlens(dimids(1))))
414                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_sp ) )
415               CASE( NF90_DOUBLE )
416                  ALLOCATE(globaldata_1d_dp(indimlens(dimids(1))))
417                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_dp ) )
418               CASE DEFAULT
419                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
420                  STOP
421            END SELECT
422
423         ELSEIF( ndims == 2 ) THEN
424
425            SELECT CASE( xtype )
426               CASE( NF90_BYTE )
427                  ALLOCATE(globaldata_2d_i1(indimlens(dimids(1)),indimlens(dimids(2))))
428                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_i1 ) )
429               CASE( NF90_SHORT )
430                  ALLOCATE(globaldata_2d_i2(indimlens(dimids(1)),indimlens(dimids(2))))
431                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_i2 ) )
432               CASE( NF90_INT )
433                  ALLOCATE(globaldata_2d_i4(indimlens(dimids(1)),indimlens(dimids(2))))
434                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_i4 ) )
435               CASE( NF90_FLOAT )
436                  ALLOCATE(globaldata_2d_sp(indimlens(dimids(1)),indimlens(dimids(2))))
437                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_sp ) )
438               CASE( NF90_DOUBLE )
439                  ALLOCATE(globaldata_2d_dp(indimlens(dimids(1)),indimlens(dimids(2))))
440                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_dp ) )
441               CASE DEFAULT
442                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
443                  STOP
444            END SELECT
445
446         ELSEIF( ndims == 3 ) THEN
447
448            SELECT CASE( xtype )
449               CASE( NF90_BYTE )
450                  ALLOCATE(globaldata_3d_i1(indimlens(dimids(1)),indimlens(dimids(2)),       &
451                     &                      indimlens(dimids(3))))
452                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_i1 ) )
453               CASE( NF90_SHORT )
454                  ALLOCATE(globaldata_3d_i2(indimlens(dimids(1)),indimlens(dimids(2)),       &
455                     &                      indimlens(dimids(3))))
456                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_i2 ) )
457               CASE( NF90_INT )
458                  ALLOCATE(globaldata_3d_i4(indimlens(dimids(1)),indimlens(dimids(2)),       &
459                     &                      indimlens(dimids(3))))
460                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_i4 ) )
461               CASE( NF90_FLOAT )
462                  ALLOCATE(globaldata_3d_sp(indimlens(dimids(1)),indimlens(dimids(2)),       &
463                     &                      indimlens(dimids(3))))
464                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_sp ) )
465               CASE( NF90_DOUBLE )
466                  ALLOCATE(globaldata_3d_dp(indimlens(dimids(1)),indimlens(dimids(2)),       &
467                     &                      indimlens(dimids(3))))
468                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_dp ) )
469               CASE DEFAULT
470                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
471                  STOP
472            END SELECT
473
474         ELSEIF( ndims == 4 ) THEN
475
476            SELECT CASE( xtype )
477               CASE( NF90_BYTE )
478                  ALLOCATE(globaldata_4d_i1(indimlens(dimids(1)),indimlens(dimids(2)),       &
479                     &                      indimlens(dimids(3)),ntchunk))
480                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_i1, start=(/1,1,1,nt/) ) )
481               CASE( NF90_SHORT )
482                  ALLOCATE(globaldata_4d_i2(indimlens(dimids(1)),indimlens(dimids(2)),       &
483                     &                      indimlens(dimids(3)),ntchunk))
484                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_i2, start=(/1,1,1,nt/) ) )
485               CASE( NF90_INT )
486                  ALLOCATE(globaldata_4d_i4(indimlens(dimids(1)),indimlens(dimids(2)),       &
487                     &                      indimlens(dimids(3)),ntchunk))
488                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_i4, start=(/1,1,1,nt/) ) )
489               CASE( NF90_FLOAT )
490                  ALLOCATE(globaldata_4d_sp(indimlens(dimids(1)),indimlens(dimids(2)),       &
491                     &                      indimlens(dimids(3)),ntchunk))
492                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_sp, start=(/1,1,1,nt/) ) )
493               CASE( NF90_DOUBLE )
494                  ALLOCATE(globaldata_4d_dp(indimlens(dimids(1)),indimlens(dimids(2)),       &
495                     &                      indimlens(dimids(3)),ntchunk))
496                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_dp, start=(/1,1,1,nt/) ) )
497               CASE DEFAULT
498                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
499                  STOP
500            END SELECT
501
502         ENDIF
503
504      ELSE  ! l_noRebuild = .false.
505
506!3.2.2 For variables that require rebuilding we need to read in from all ndomain files
507!      Here we allocate global variables ahead of looping over files
508         IF( nchunksize == nmax_unlimited .OR. ndims <= 3 ) THEN
509            IF (l_verbose) WRITE(numout,*) 'Rebuilding data from variable '//TRIM(varname)//'...'
510         ELSE
511            IF (l_verbose) WRITE(numout,'(A,I3,A,I3,A)') ' Rebuilding data from variable '  &
512            &                 //TRIM(varname)//' for chunks ',nt,' to ',nt+ntchunk-1,' ...'
513         ENDIF
514         IF( ndims == 1 ) THEN
515
516            SELECT CASE( xtype )
517               CASE( NF90_BYTE )
518                  ALLOCATE(globaldata_1d_i1(outdimlens(dimids(1))))
519               CASE( NF90_SHORT )
520                  ALLOCATE(globaldata_1d_i2(outdimlens(dimids(1))))
521               CASE( NF90_INT )
522                  ALLOCATE(globaldata_1d_i4(outdimlens(dimids(1))))
523               CASE( NF90_FLOAT )
524                  ALLOCATE(globaldata_1d_sp(outdimlens(dimids(1))))
525               CASE( NF90_DOUBLE )
526                  ALLOCATE(globaldata_1d_dp(outdimlens(dimids(1))))
527               CASE DEFAULT
528                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
529                  STOP
530            END SELECT
531
532         ELSEIF( ndims == 2 ) THEN
533
534            SELECT CASE( xtype )
535               CASE( NF90_BYTE )
536                  ALLOCATE(globaldata_2d_i1(outdimlens(dimids(1)),outdimlens(dimids(2))))
537               CASE( NF90_SHORT )
538                  ALLOCATE(globaldata_2d_i2(outdimlens(dimids(1)),outdimlens(dimids(2))))
539               CASE( NF90_INT )
540                  ALLOCATE(globaldata_2d_i4(outdimlens(dimids(1)),outdimlens(dimids(2))))
541               CASE( NF90_FLOAT )
542                  ALLOCATE(globaldata_2d_sp(outdimlens(dimids(1)),outdimlens(dimids(2))))
543               CASE( NF90_DOUBLE )
544                  ALLOCATE(globaldata_2d_dp(outdimlens(dimids(1)),outdimlens(dimids(2))))
545               CASE DEFAULT
546                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
547                  STOP
548            END SELECT
549
550         ELSEIF( ndims == 3 ) THEN
551
552            SELECT CASE( xtype )
553               CASE( NF90_BYTE )
554                  ALLOCATE(globaldata_3d_i1(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
555                     &                      outdimlens(dimids(3))))
556               CASE( NF90_SHORT )
557                  ALLOCATE(globaldata_3d_i2(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
558                     &                      outdimlens(dimids(3))))
559               CASE( NF90_INT )
560                  ALLOCATE(globaldata_3d_i4(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
561                     &                      outdimlens(dimids(3))))
562               CASE( NF90_FLOAT )
563                  ALLOCATE(globaldata_3d_sp(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
564                     &                      outdimlens(dimids(3))))
565               CASE( NF90_DOUBLE )
566                  ALLOCATE(globaldata_3d_dp(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
567                     &                      outdimlens(dimids(3))))
568               CASE DEFAULT
569                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
570                  STOP
571            END SELECT
572
573         ELSEIF( ndims == 4 ) THEN
574
575            SELECT CASE( xtype )
576               CASE( NF90_BYTE )
577                  ALLOCATE(globaldata_4d_i1(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
578                     &                      outdimlens(dimids(3)),ntchunk))
579               CASE( NF90_SHORT )
580                  ALLOCATE(globaldata_4d_i2(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
581                     &                      outdimlens(dimids(3)),ntchunk))
582               CASE( NF90_INT )
583                  ALLOCATE(globaldata_4d_i4(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
584                     &                      outdimlens(dimids(3)),ntchunk))
585               CASE( NF90_FLOAT )
586                  ALLOCATE(globaldata_4d_sp(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
587                     &                      outdimlens(dimids(3)),ntchunk))
588               CASE( NF90_DOUBLE )
589                  ALLOCATE(globaldata_4d_dp(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
590                     &                      outdimlens(dimids(3)),ntchunk))
591               CASE DEFAULT
592                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
593                  STOP
594            END SELECT
595         ELSE
596            WRITE(numerr,*) 'ERROR! : A netcdf variable has more than 4 dimensions which is not taken into account'
597            STOP
598         ENDIF
599
600!$OMP  PARALLEL DO DEFAULT(NONE)                                                          &
601!$OMP& PRIVATE(ifile,ncid,xtype,start_pos,local_sizes,InMin,InMax,natts,                  &
602!$OMP&         ndims,attid,attname,dimids,idim,dimname,dimlen,unlimitedDimId,             &
603!$OMP&         halo_start,halo_end,idomain,jdomain,rdomain,di,dj,dr,                      &
604!$OMP&         localdata_1d_i2,localdata_1d_i4,localdata_1d_sp,localdata_1d_dp,           &
605!$OMP&         localdata_2d_i2,localdata_2d_i4,localdata_2d_sp,localdata_2d_dp,           &
606!$OMP&         localdata_3d_i2,localdata_3d_i4,localdata_3d_sp,localdata_3d_dp,           &
607!$OMP&         localdata_4d_i2,localdata_4d_i4,localdata_4d_sp,localdata_4d_dp,           &
608!$OMP&         localdata_1d_i1,localdata_2d_i1,localdata_3d_i1,localdata_4d_i1)           &
609!$OMP& SHARED(jv,nvars,varname,filenames,ValMin,ValMax,indimlens,outdimlens,rbdims,       &
610!$OMP&        ndomain,outid,chunksize,istop,l_valid,nthreads,inncids,rebuild_dims,        &
611!$OMP&        globaldata_1d_i2,globaldata_1d_i4,globaldata_1d_sp,globaldata_1d_dp,        &
612!$OMP&        globaldata_2d_i2,globaldata_2d_i4,globaldata_2d_sp,globaldata_2d_dp,        &
613!$OMP&        globaldata_3d_i2,globaldata_3d_i4,globaldata_3d_sp,globaldata_3d_dp,        &
614!$OMP&        globaldata_4d_i2,globaldata_4d_i4,globaldata_4d_sp,globaldata_4d_dp,        &
615!$OMP&        globaldata_1d_i1,globaldata_2d_i1,globaldata_3d_i1,globaldata_4d_i1,        &
616!$OMP&        ntchunk,nt,nmax_unlimited)
617
618         DO ifile = 1, ndomain
619
620            ncid = inncids(ifile)
621!$OMP CRITICAL
622            CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_size_local', local_sizes ), istop )
623            CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_position_first', start_pos ), istop )
624            CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_halo_size_start', halo_start ), istop )
625            CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_halo_size_end', halo_end ), istop )
626            CALL check_nf90( nf90_inquire_variable( ncid, jv, varname, xtype, ndims, dimids, natts ), istop )
627            CALL check_nf90( nf90_inquire( ncid, unlimitedDimId = unlimitedDimId ), istop )
628!$OMP END CRITICAL
629
630            ! set defaults for rebuilding so that i is 1st, j 2nd
631            di=1
632            dj=2
633
634            IF( rbdims == 1 ) THEN
635               ! override defaults above and set other variables
636               start_pos(2) = 1
637               local_sizes(2) = outdimlens(3-dimids(2))
638               halo_end(2) = 0
639               halo_start(2) = 0
640               di=rebuild_dims(1)
641               dj=3-di
642            ENDIF
643
644!3.3.1 Generate local domain interior sizes from local_sizes and halo sizes
645!      idomain defines the 1st and last interior points in the i direction and
646!      jdomain defines the 1st and last interior points in the j direction
647
648            idomain(1) = 1 + halo_start(di)
649            idomain(2) = local_sizes(di) - halo_end(di)
650            jdomain(1) = 1 + halo_start(dj)
651            jdomain(2) = local_sizes(dj) - halo_end(dj)
652
653!3.3.2 For rbdims or more dimensions put the data array from this input file into the correct
654!      part of the output data array. Assume the first dimensions are those to be rebuilt.
655
656            IF( ndims == 1 ) THEN
657
658               IF( rebuild_dims(1) == 1 ) THEN
659                  dr = di
660                  rdomain = idomain
661               ELSE
662                  dr = dj
663                  rdomain = jdomain
664               ENDIF
665
666              SELECT CASE( xtype )
667                  CASE( NF90_BYTE )
668                     ALLOCATE(localdata_1d_i1(local_sizes(dr)))
669                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_i1 ), istop )
670                     DO jr = rdomain(1), rdomain(2)
671                        globaldata_1d_i1(start_pos(dr) + jr - 1) = localdata_1d_i1(jr)
672                     END DO
673                     DEALLOCATE(localdata_1d_i1)
674                  CASE( NF90_SHORT )
675                     ALLOCATE(localdata_1d_i2(local_sizes(dr)))
676                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_i2 ), istop )
677                     DO jr = rdomain(1), rdomain(2)
678                        globaldata_1d_i2(start_pos(dr) + jr - 1) = localdata_1d_i2(jr)
679                     END DO
680                     DEALLOCATE(localdata_1d_i2)
681                  CASE( NF90_INT )
682                     ALLOCATE(localdata_1d_i4(local_sizes(dr)))
683                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_i4 ), istop )
684                     DO jr = rdomain(1), rdomain(2)
685                        globaldata_1d_i4(start_pos(dr) + jr - 1) = localdata_1d_i4(jr)
686                     END DO
687                     DEALLOCATE(localdata_1d_i4)
688                  CASE( NF90_FLOAT )
689                     ALLOCATE(localdata_1d_sp(local_sizes(dr)))
690                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_sp ), istop )
691                     DO jr = rdomain(1), rdomain(2)
692                        globaldata_1d_sp(start_pos(dr) + jr - 1) = localdata_1d_sp(jr)
693                     END DO
694                     DEALLOCATE(localdata_1d_sp)
695                  CASE( NF90_DOUBLE )
696                     ALLOCATE(localdata_1d_dp(local_sizes(dr)))
697                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_dp ), istop )
698                     DO jr = rdomain(1), rdomain(2)
699                        globaldata_1d_dp(start_pos(dr) + jr - 1) = localdata_1d_dp(jr)
700                     END DO
701                     DEALLOCATE(localdata_1d_dp)
702                  CASE DEFAULT
703                     WRITE(numerr,*) 'Unknown nf90 type: ', xtype
704                     istop = istop + 1
705               END SELECT
706
707            ELSEIF( ndims == 2 ) THEN
708       
709               SELECT CASE( xtype )
710                  CASE( NF90_BYTE )
711                     ALLOCATE(localdata_2d_i1(local_sizes(di),local_sizes(dj)))
712                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_i1 ), istop )
713                        DO jj = jdomain(1), jdomain(2)
714                           DO ji = idomain(1), idomain(2)
715                              globaldata_2d_i1(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_i1(ji,jj)
716                        END DO
717                     END DO
718                     DEALLOCATE(localdata_2d_i1)
719                  CASE( NF90_SHORT )
720                     ALLOCATE(localdata_2d_i2(local_sizes(di),local_sizes(dj)))
721                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_i2 ), istop )
722                     DO jj = jdomain(1), jdomain(2)
723                        DO ji = idomain(1), idomain(2)
724                           globaldata_2d_i2(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_i2(ji,jj)
725                        END DO
726                     END DO
727                     DEALLOCATE(localdata_2d_i2)
728                  CASE( NF90_INT )
729                     ALLOCATE(localdata_2d_i4(local_sizes(di),local_sizes(dj)))
730                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_i4 ), istop )
731                     DO jj = jdomain(1), jdomain(2)
732                        DO ji = idomain(1), idomain(2)
733                           globaldata_2d_i4(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_i4(ji,jj)
734                        END DO
735                     END DO
736                     DEALLOCATE(localdata_2d_i4)
737                  CASE( NF90_FLOAT )
738                     ALLOCATE(localdata_2d_sp(local_sizes(di),local_sizes(dj)))
739                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_sp ), istop )
740                     DO jj = jdomain(1), jdomain(2)
741                        DO ji = idomain(1), idomain(2)
742                           globaldata_2d_sp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_sp(ji,jj)
743                        END DO
744                     END DO
745                     DEALLOCATE(localdata_2d_sp)
746                  CASE( NF90_DOUBLE )
747                     ALLOCATE(localdata_2d_dp(local_sizes(di),local_sizes(dj)))
748                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_dp ), istop )
749                     DO jj = jdomain(1), jdomain(2)
750                        DO ji = idomain(1), idomain(2)
751                           globaldata_2d_dp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_dp(ji,jj)
752                        END DO
753                     END DO
754                     DEALLOCATE(localdata_2d_dp) 
755                  CASE DEFAULT
756                     WRITE(numerr,*) 'Unknown nf90 type: ', xtype
757                     istop = istop + 1
758               END SELECT
759
760            ELSEIF( ndims == 3 ) THEN
761
762               SELECT CASE( xtype )
763                  CASE( NF90_BYTE )
764                     ALLOCATE(localdata_3d_i1(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
765                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_i1 ), istop )
766!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
767!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_i1,localdata_3d_i1,di,dj)
768                     DO jk = 1, indimlens(dimids(3))
769                        DO jj = jdomain(1), jdomain(2)
770                           DO ji = idomain(1), idomain(2)
771                              globaldata_3d_i1(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_i1(ji,jj,jk)
772                           END DO
773                        END DO
774                     END DO
775!$OMP END PARALLEL DO
776                     DEALLOCATE(localdata_3d_i1)
777                  CASE( NF90_SHORT )
778                     ALLOCATE(localdata_3d_i2(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
779                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_i2 ), istop )
780!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
781!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_i2,localdata_3d_i2,di,dj)
782                     DO jk = 1, indimlens(dimids(3))
783                        DO jj = jdomain(1), jdomain(2)
784                           DO ji = idomain(1), idomain(2)
785                              globaldata_3d_i2(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_i2(ji,jj,jk)
786                           END DO
787                        END DO
788                     END DO
789!$OMP END PARALLEL DO
790                     DEALLOCATE(localdata_3d_i2)
791                  CASE( NF90_INT )
792                     ALLOCATE(localdata_3d_i4(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
793                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_i4 ), istop )
794!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
795!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_i4,localdata_3d_i4,di,dj)
796                     DO jk = 1, indimlens(dimids(3))
797                        DO jj = jdomain(1), jdomain(2)
798                           DO ji = idomain(1), idomain(2)
799                              globaldata_3d_i4(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_i4(ji,jj,jk)
800                           END DO
801                        END DO
802                     END DO
803!$OMP END PARALLEL DO
804                     DEALLOCATE(localdata_3d_i4)
805                  CASE( NF90_FLOAT )
806                     ALLOCATE(localdata_3d_sp(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
807                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_sp ), istop )
808!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
809!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_sp,localdata_3d_sp,di,dj)
810                     DO jk = 1, indimlens(dimids(3))
811                        DO jj = jdomain(1), jdomain(2)
812                           DO ji = idomain(1), idomain(2)
813                              globaldata_3d_sp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_sp(ji,jj,jk)
814                           END DO
815                        END DO
816                     END DO
817!$OMP END PARALLEL DO
818                     DEALLOCATE(localdata_3d_sp)
819                  CASE( NF90_DOUBLE )
820                     ALLOCATE(localdata_3d_dp(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
821                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_dp ), istop )
822!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
823!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_dp,localdata_3d_dp,di,dj)
824                     DO jk = 1, indimlens(dimids(3))
825                        DO jj = jdomain(1), jdomain(2)
826                           DO ji = idomain(1), idomain(2)
827                              globaldata_3d_dp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_dp(ji,jj,jk)
828                           END DO
829                        END DO
830                     END DO
831!$OMP END PARALLEL DO
832                     DEALLOCATE(localdata_3d_dp)
833                  CASE DEFAULT
834                     WRITE(numerr,*) 'Unknown nf90 type: ', xtype
835                     istop = istop + 1
836               END SELECT
837     
838            ELSEIF (ndims == 4) THEN
839       
840               SELECT CASE( xtype )
841                  CASE( NF90_BYTE )
842                     ALLOCATE(localdata_4d_i1(local_sizes(di),local_sizes(dj),               &
843                         &                     indimlens(dimids(3)),ntchunk))
844                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_i1, start=(/1,1,1,nt/) ), istop )
845!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
846!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_i1,localdata_4d_i1,di,dj,nt,ntchunk)
847                     DO jl = 1, ntchunk
848!$OMP DO
849                        DO jk = 1, indimlens(dimids(3))
850                           DO jj = jdomain(1), jdomain(2)
851                              DO ji = idomain(1), idomain(2)
852                                 globaldata_4d_i1(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_i1(ji,jj,jk,jl)
853                              END DO
854                           END DO
855                        END DO
856!$OMP END DO nowait
857                     END DO
858!$OMP END PARALLEL
859                     DEALLOCATE(localdata_4d_i1)
860                  CASE( NF90_SHORT )
861                     ALLOCATE(localdata_4d_i2(local_sizes(di),local_sizes(dj),               &
862                        &                     indimlens(dimids(3)),ntchunk))
863                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_i2, start=(/1,1,1,nt/) ), istop )
864!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
865!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_i2,localdata_4d_i2,di,dj,nt,ntchunk)
866                     DO jl = 1, ntchunk
867!$OMP DO
868                        DO jk = 1, indimlens(dimids(3))
869                           DO jj = jdomain(1), jdomain(2) 
870                              DO ji = idomain(1), idomain(2)
871                                 globaldata_4d_i2(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_i2(ji,jj,jk,jl)
872                              END DO
873                           END DO
874                        END DO
875!$OMP END DO nowait
876                     END DO
877!$OMP END PARALLEL
878                     DEALLOCATE(localdata_4d_i2)
879                  CASE( NF90_INT )
880                     ALLOCATE(localdata_4d_i4(local_sizes(di),local_sizes(dj),               &
881                        &                     indimlens(dimids(3)),ntchunk))
882                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_i4, start=(/1,1,1,nt/) ), istop )
883!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
884!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_i4,localdata_4d_i4,di,dj,nt,ntchunk)
885                     DO jl = 1, ntchunk
886!$OMP DO
887                        DO jk = 1, indimlens(dimids(3))
888                           DO jj = jdomain(1), jdomain(2) 
889                              DO ji = idomain(1), idomain(2)
890                                 globaldata_4d_i4(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_i4(ji,jj,jk,jl)
891                              END DO
892                           END DO
893                        END DO
894!$OMP END DO nowait
895                     END DO
896!$OMP END PARALLEL
897                     DEALLOCATE(localdata_4d_i4)
898                  CASE( NF90_FLOAT )
899                     ALLOCATE(localdata_4d_sp(local_sizes(di),local_sizes(dj),               &
900                        &                     indimlens(dimids(3)),ntchunk))
901                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_sp, start=(/1,1,1,nt/) ), istop )
902!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
903!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_sp,localdata_4d_sp,di,dj,nt,ntchunk)
904                     DO jl = 1, ntchunk
905!$OMP DO
906                        DO jk = 1, indimlens(dimids(3))
907                           DO jj = jdomain(1), jdomain(2) 
908                              DO ji = idomain(1), idomain(2)
909                                 globaldata_4d_sp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_sp(ji,jj,jk,jl)
910                              END DO
911                           END DO
912                        END DO
913!$OMP END DO nowait
914                     END DO
915!$OMP END PARALLEL
916                     DEALLOCATE(localdata_4d_sp)
917                  CASE( NF90_DOUBLE )
918                     ALLOCATE(localdata_4d_dp(local_sizes(di),local_sizes(dj),               &
919                        &                     indimlens(dimids(3)),ntchunk))
920                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_dp, start=(/1,1,1,nt/) ), istop )
921!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
922!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_dp,localdata_4d_dp,di,dj,nt,ntchunk)
923                     DO jl = 1, ntchunk
924!$OMP DO
925                        DO jk = 1, indimlens(dimids(3))
926                           DO jj = jdomain(1), jdomain(2) 
927                              DO ji = idomain(1), idomain(2)
928                                 globaldata_4d_dp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_dp(ji,jj,jk,jl)
929                              END DO
930                           END DO
931                        END DO
932!$OMP END DO nowait
933                     END DO
934!$OMP END PARALLEL
935                     DEALLOCATE(localdata_4d_dp)
936                  CASE DEFAULT
937                     WRITE(numerr,*) 'Unknown nf90 type: ', xtype
938                     istop = istop + 1
939               END SELECT
940
941            ENDIF ! l_noRebuild false
942
943!3.4 Work out if the valid_min and valid_max attributes exist for this variable.
944!    If they do then calculate the extrema over all input files.
945
946            DO attid = 1, natts
947               CALL check_nf90( nf90_inq_attname( ncid, jv, attid, attname ), istop )
948               IF( INDEX( attname, "valid_min" ) == 1 ) THEN
949                  CALL check_nf90( nf90_get_att( ncid, jv, attname, InMin), istop )
950                  l_valid = .true.
951               ENDIF
952               IF( INDEX( attname, "valid_max" ) == 1 ) THEN
953                  CALL check_nf90( nf90_get_att( ncid, jv, attname, InMax ), istop )
954                  l_valid = .true.
955               ENDIF
956            END DO
957
958            IF (l_valid) THEN
959!$OMP CRITICAL
960               IF( InMin < ValMin ) ValMin = InMin
961               IF( InMax > ValMax ) ValMax = InMax
962!$OMP END CRITICAL
963            ENDIF
964
965!3.5 Abort if failure and only 1 thread
966
967            IF( nthreads == 1 .AND. istop /= nf90_noerr )  THEN
968               WRITE(numerr,*) '*** NEMO rebuild failed! ***'
969               STOP
970            ENDIF
971       
972         END DO  ! loop over files
973!$OMP END PARALLEL DO
974
975!3.6 Abort if any of the OMP threads failed
976         IF( istop /= nf90_noerr )  THEN
977            WRITE(numerr,*) '*** NEMO rebuild failed! ***'
978            STOP
979         ENDIF
980
981      ENDIF ! ndims > 2
982
983!---------------------------------------------------------------------------
984!4. Write data to output file
985
986      IF (l_verbose) WRITE(numout,*) 'Writing variable '//TRIM(varname)//'...'
987
988!4.1 If the valid min and max attributes exist then update them in the file
989
990      IF( l_valid ) THEN
991         CALL check_nf90( nf90_put_att( outid, jv, "valid_min", ValMin ) )
992         CALL check_nf90( nf90_put_att( outid, jv, "valid_max", ValMax ) )
993      ENDIF
994
995!4.2 Write the data to the output file depending on how many dimensions it has
996
997      IF( ndims == 0 ) THEN
998
999         SELECT CASE( xtype )
1000            CASE( NF90_BYTE )
1001               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_i1 ) )
1002            CASE( NF90_SHORT )
1003               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_i2 ) )
1004            CASE( NF90_INT )
1005               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_i4 ) )
1006            CASE( NF90_FLOAT )
1007               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_sp ) )
1008            CASE( NF90_DOUBLE )
1009               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_dp ) )
1010         END SELECT
1011
1012      ELSEIF( ndims == 1 ) THEN
1013
1014         SELECT CASE( xtype )
1015            CASE( NF90_BYTE )
1016               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_i1 ) )
1017               DEALLOCATE(globaldata_1d_i1)
1018            CASE( NF90_SHORT )
1019               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_i2 ) )
1020               DEALLOCATE(globaldata_1d_i2)
1021            CASE( NF90_INT )
1022               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_i4 ) )
1023               DEALLOCATE(globaldata_1d_i4)
1024            CASE( NF90_FLOAT )
1025               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_sp ) )
1026               DEALLOCATE(globaldata_1d_sp)
1027            CASE( NF90_DOUBLE )
1028               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_dp ) )
1029               DEALLOCATE(globaldata_1d_dp)
1030         END SELECT
1031
1032      ELSEIF( ndims == 2 ) THEN 
1033     
1034         SELECT CASE( xtype )   
1035            CASE( NF90_BYTE )                   
1036               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_i1 ) )
1037               DEALLOCATE(globaldata_2d_i1)
1038            CASE( NF90_SHORT )                   
1039               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_i2 ) )
1040               DEALLOCATE(globaldata_2d_i2)
1041            CASE( NF90_INT )                             
1042               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_i4 ) )
1043               DEALLOCATE(globaldata_2d_i4)
1044            CASE( NF90_FLOAT )                             
1045               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_sp ) )
1046               DEALLOCATE(globaldata_2d_sp)
1047            CASE( NF90_DOUBLE )                                         
1048               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_dp ) )
1049               DEALLOCATE(globaldata_2d_dp)
1050            CASE DEFAULT   
1051               WRITE(numerr,*) 'Unknown nf90 type: ', xtype
1052               STOP
1053         END SELECT     
1054                     
1055      ELSEIF( ndims == 3 ) THEN
1056     
1057         SELECT CASE( xtype ) 
1058            CASE( NF90_BYTE )                   
1059               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_i1 ) )
1060               DEALLOCATE(globaldata_3d_i1)
1061            CASE( NF90_SHORT )                   
1062               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_i2 ) )
1063               DEALLOCATE(globaldata_3d_i2)
1064            CASE( NF90_INT )                             
1065               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_i4 ) )
1066               DEALLOCATE(globaldata_3d_i4)
1067            CASE( NF90_FLOAT )                             
1068               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_sp ) )
1069               DEALLOCATE(globaldata_3d_sp)
1070            CASE( NF90_DOUBLE )                                         
1071               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_dp ) )
1072               DEALLOCATE(globaldata_3d_dp)
1073            CASE DEFAULT   
1074               WRITE(numerr,*) 'Unknown nf90 type: ', xtype
1075               STOP
1076         END SELECT     
1077   
1078      ELSEIF( ndims == 4 ) THEN
1079
1080         SELECT CASE( xtype )   
1081            CASE( NF90_BYTE )                   
1082               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_i1, start=(/1,1,1,nt/) ) )
1083               DEALLOCATE(globaldata_4d_i1)
1084            CASE( NF90_SHORT )                   
1085               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_i2, start=(/1,1,1,nt/) ) )
1086               DEALLOCATE(globaldata_4d_i2)
1087            CASE( NF90_INT )                             
1088               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_i4, start=(/1,1,1,nt/) ) )
1089               DEALLOCATE(globaldata_4d_i4)
1090            CASE( NF90_FLOAT )                             
1091               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_sp, start=(/1,1,1,nt/) ) )
1092               DEALLOCATE(globaldata_4d_sp)
1093            CASE( NF90_DOUBLE )                                         
1094               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_dp, start=(/1,1,1,nt/) ) )
1095               DEALLOCATE(globaldata_4d_dp)
1096            CASE DEFAULT   
1097               WRITE(numerr,*) 'Unknown nf90 type: ', xtype
1098               STOP
1099         END SELECT     
1100   
1101      ENDIF
1102
1103         nt = nt + ntchunk
1104
1105      END DO ! WHILE loop
1106   
1107   END DO  ! loop over variables
1108
1109!---------------------------------------------------------------------------
1110!5. Close files
1111
1112!5.1 Close input files
1113
1114   IF (l_verbose) WRITE(numout,*) 'Closing input files...'
1115   DO ifile = 1, ndomain
1116      ncid = inncids(ifile)
1117      CALL check_nf90( nf90_close( ncid ) )
1118   END DO
1119
1120!5.2 Close output file
1121
1122   IF (l_verbose) WRITE(numout,*) 'Closing output file...'
1123   CALL check_nf90( nf90_close( outid ) )
1124
1125   IF (l_verbose) WRITE(numout,*) 'NEMO rebuild completed successfully'
1126   IF (l_verbose) WRITE(numout,*)
1127
1128CONTAINS
1129
1130
1131   SUBROUTINE check_nf90(status, errorFlag)
1132   !---------------------------------------------------------------------
1133   !  Checks return code from nf90 library calls and warns if needed
1134   !  If errorFlag is present then it just increments this flag (OMP use)
1135   !
1136   !---------------------------------------------------------------------
1137      INTEGER, INTENT(IN   ) :: status
1138      INTEGER, INTENT(INOUT), OPTIONAL :: errorFlag
1139   !---------------------------------------------------------------------
1140
1141      IF( status /= nf90_noerr ) THEN
1142         WRITE(numerr,*) 'ERROR! : '//TRIM(nf90_strerror(status))
1143         IF( PRESENT( errorFlag ) ) THEN
1144            errorFlag = errorFlag + status
1145         ELSE
1146            WRITE(numerr,*) "*** NEMO rebuild failed ***"
1147            WRITE(numerr,*)
1148            STOP
1149         ENDIF
1150      ENDIF
1151
1152   END SUBROUTINE check_nf90
1153
1154
1155END PROGRAM rebuild_nemo
Note: See TracBrowser for help on using the repository browser.