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

source: branches/2011/dev_r2802_UKMET3_rebuild/NEMOGCM/TOOLS/REBUILD_NEMO/src/rebuild_nemo.f90 @ 2898

Last change on this file since 2898 was 2898, checked in by edblockley, 13 years ago

4th commit for rebuild branch; Moving rebuild_nemo code from TOOLS/REBUILD to TOOLS/REBUILD_NEMO (NB. this could be renamed REBUILD at a later date if that directory were to be removed) . see ticket:#871

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