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

Last change on this file since 3019 was 3019, checked in by edblockley, 12 years ago

5th commit for rebuild branch; Fixing a small bug in the rebuild to ensure that the rebuild dimension no.s are properly set when the names are specified in the namelist. Also changing indimlens to SHARED rather than PRIVATE for OMP loops to satisfy certain compilers. see ticket:#871

File size: 51.8 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))) THEN
240            dimlen = global_sizes(1)
241            rebuild_dims(1) = idim
242         ENDIF
243         IF( rbdims > 1 .AND. TRIM(dimname) == TRIM(dims(2))) THEN
244            dimlen = global_sizes(2)
245            rebuild_dims(2) = idim
246         ENDIF
247      ENDIF
248
249      IF( idim == unlimitedDimId ) THEN
250         CALL check_nf90( nf90_def_dim( outid, dimname, nf90_unlimited, dimid) )
251      ELSE
252         CALL check_nf90( nf90_def_dim( outid, dimname, dimlen, dimid) )
253      ENDIF
254      outdimlens(idim) = dimlen
255   END DO
256
257   IF( istop == 1 ) THEN
258      WRITE(numerr,*) 'ERROR! : DOMAIN_local_sizes attribute does not match rebuild dimension lengths in the first file'
259      WRITE(numerr,*) 'Attribute DOMAIN_local_sizes is : ', local_sizes
260      WRITE(numerr,*) 'Dimensions to be rebuilt are of size : ', outdimlens(rebuild_dims(1)), outdimlens(rebuild_dims(2)) 
261      STOP
262   ENDIF
263
264   IF (l_findDims) THEN
265      IF (l_verbose) WRITE(numout,*) 'Finding rebuild dimensions from the first file...'
266   ELSE
267      IF (l_verbose) WRITE(numout,*) 'Using rebuild dimensions given in namelist...'
268   ENDIF
269
270   IF( rbdims > 1 ) THEN
271      IF (l_verbose) WRITE(numout,*) 'Rebuilding across dimensions '//TRIM(indimnames(rebuild_dims(1)))//  &
272         &                      ' and '//TRIM(indimnames(rebuild_dims(2)))
273   ELSE
274      IF (l_verbose) WRITE(numout,*) 'Rebuilding across dimension '//TRIM(indimnames(rebuild_dims(1)))
275   ENDIF
276     
277!2.2.2 Copy the global attributes into the output file, apart from those beginning with DOMAIN_ 
278!      Also need to change the file_name attribute and the TimeStamp attribute.
279   DO attid = 1, natts
280      CALL check_nf90( nf90_inq_attname( ncid, nf90_global, attid, attname ) )
281      IF( INDEX( attname, "DOMAIN_" ) == 1 ) CYCLE
282      IF( INDEX( attname, "file_name") == 1 ) CYCLE
283      IF( INDEX( attname, "associate_file") == 1 ) CYCLE
284      IF (l_verbose) WRITE(numout,*) 'Copying attribute '//TRIM(attname)//' into destination file...'
285      CALL check_nf90( nf90_copy_att( ncid, nf90_global, attname, outid, nf90_global ) )
286   END DO
287   CALL check_nf90( nf90_put_att( outid, nf90_global, "file_name", TRIM(filebase)//'.nc') )
288   IF (l_verbose) WRITE(numout,*) 'Writing new file_name attribute' 
289   CALL DATE_AND_TIME ( date=date, time=time, zone=zone )
290   timestamp = date(7:8) // "/" // date(5:6) // "/" // date(1:4) // " " // &
291               time(1:2) // ":" // time(3:4) // ":" // time(5:6) // " " // &
292               zone 
293   CALL check_nf90( nf90_put_att( outid, nf90_global, "TimeStamp", timestamp ) )
294   IF (l_verbose) WRITE(numout,*) 'Writing new TimeStamp attribute'
295 
296!2.2.3 Copy the variable definitions and attributes into the output file.
297   DO jv = 1, nvars
298      CALL check_nf90( nf90_inquire_variable( ncid, jv, varname, xtype, ndims, dimids, natts ) )
299      ALLOCATE(outdimids(ndims))
300      DO idim = 1, ndims
301         outdimids(idim) = dimids(idim)
302      END DO
303      CALL check_nf90( nf90_def_var( outid, varname, xtype, outdimids, varid ) )
304      DEALLOCATE(outdimids)
305      IF (l_verbose) WRITE(numout,*) 'Defining variable '//TRIM(varname)//'...' 
306      IF( natts > 0 ) THEN
307         DO attid = 1, natts
308            CALL check_nf90( nf90_inq_attname( ncid, varid, attid, attname ) )
309            CALL check_nf90( nf90_copy_att( ncid, varid, attname, outid, varid ) )
310         END DO
311      ENDIF
312   END DO
313 
314!2.3 End definitions in output file and copy 1st file ncid to the inncids array
315
316   CALL check_nf90( nf90_enddef( outid ) )
317   inncids(1) = ncid
318   IF (l_verbose) WRITE(numout,*) 'Finished defining output file.'
319 
320!---------------------------------------------------------------------------
321!3. Read in data from each file for each variable
322
323!3.1 Open each file and store the ncid in inncids array
324
325   IF (l_verbose) WRITE(numout,*) 'Opening input files...'
326   DO ifile = 2, ndomain
327      CALL check_nf90( nf90_open( TRIM(filenames(ifile)), nf90_share, ncid, chunksize=chunksize ) )
328      inncids(ifile) = ncid
329   END DO
330   IF (l_verbose) WRITE(numout,*) 'All input files open.'
331
332   DO jv = 1, nvars
333
334      ValMin = 1.e10
335      ValMax = -1.e10
336      l_valid = .false.
337      istop = nf90_noerr
338
339!3.2 Inquire variable to find out name and how many dimensions it has
340!    and importantly whether it contains the dimensions in rebuild_dims()
341
342      ncid = inncids(1)
343      CALL check_nf90( nf90_inquire_variable( ncid, jv, varname, xtype, ndims, dimids, natts ) )
344
345      l_noRebuild = .true.
346      IF( ANY( dimids(1:ndims) == rebuild_dims(1) )) l_noRebuild = .false.
347      IF( rbdims > 1 ) THEN
348         IF( ANY( dimids(1:ndims) == rebuild_dims(2) )) l_noRebuild = .false.
349      ENDIF
350
351      IF (l_noRebuild) THEN
352
353         IF (l_verbose) WRITE(numout,*) 'Copying data from variable '//TRIM(varname)//'...'
354
355!3.2.1 If rebuilding not required then just need to read in variable
356!      for copying direct into output file after the OMP (files) loop.
357         IF( ndims == 0 ) THEN
358
359            SELECT CASE( xtype )
360               CASE( NF90_BYTE )
361                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_i1 ) )
362               CASE( NF90_SHORT )
363                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_i2 ) )
364               CASE( NF90_INT )
365                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_i4 ) )
366               CASE( NF90_FLOAT )
367                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_sp ) )
368               CASE( NF90_DOUBLE )
369                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_0d_dp ) )
370               CASE DEFAULT
371                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
372                  STOP
373            END SELECT
374
375         ELSEIF( ndims == 1 ) THEN
376
377            SELECT CASE( xtype )
378               CASE( NF90_BYTE )
379                  ALLOCATE(globaldata_1d_i1(indimlens(dimids(1))))
380                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_i1 ) )
381               CASE( NF90_SHORT )
382                  ALLOCATE(globaldata_1d_i2(indimlens(dimids(1))))
383                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_i2 ) )
384               CASE( NF90_INT )
385                  ALLOCATE(globaldata_1d_i4(indimlens(dimids(1))))
386                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_i4 ) )
387               CASE( NF90_FLOAT )
388                  ALLOCATE(globaldata_1d_sp(indimlens(dimids(1))))
389                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_sp ) )
390               CASE( NF90_DOUBLE )
391                  ALLOCATE(globaldata_1d_dp(indimlens(dimids(1))))
392                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_1d_dp ) )
393               CASE DEFAULT
394                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
395                  STOP
396            END SELECT
397
398         ELSEIF( ndims == 2 ) THEN
399
400            SELECT CASE( xtype )
401               CASE( NF90_BYTE )
402                  ALLOCATE(globaldata_2d_i1(indimlens(dimids(1)),indimlens(dimids(2))))
403                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_i1 ) )
404               CASE( NF90_SHORT )
405                  ALLOCATE(globaldata_2d_i2(indimlens(dimids(1)),indimlens(dimids(2))))
406                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_i2 ) )
407               CASE( NF90_INT )
408                  ALLOCATE(globaldata_2d_i4(indimlens(dimids(1)),indimlens(dimids(2))))
409                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_i4 ) )
410               CASE( NF90_FLOAT )
411                  ALLOCATE(globaldata_2d_sp(indimlens(dimids(1)),indimlens(dimids(2))))
412                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_sp ) )
413               CASE( NF90_DOUBLE )
414                  ALLOCATE(globaldata_2d_dp(indimlens(dimids(1)),indimlens(dimids(2))))
415                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_2d_dp ) )
416               CASE DEFAULT
417                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
418                  STOP
419            END SELECT
420
421         ELSEIF( ndims == 3 ) THEN
422
423            SELECT CASE( xtype )
424               CASE( NF90_BYTE )
425                  ALLOCATE(globaldata_3d_i1(indimlens(dimids(1)),indimlens(dimids(2)),       &
426                     &                      indimlens(dimids(3))))
427                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_i1 ) )
428               CASE( NF90_SHORT )
429                  ALLOCATE(globaldata_3d_i2(indimlens(dimids(1)),indimlens(dimids(2)),       &
430                     &                      indimlens(dimids(3))))
431                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_i2 ) )
432               CASE( NF90_INT )
433                  ALLOCATE(globaldata_3d_i4(indimlens(dimids(1)),indimlens(dimids(2)),       &
434                     &                      indimlens(dimids(3))))
435                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_i4 ) )
436               CASE( NF90_FLOAT )
437                  ALLOCATE(globaldata_3d_sp(indimlens(dimids(1)),indimlens(dimids(2)),       &
438                     &                      indimlens(dimids(3))))
439                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_sp ) )
440               CASE( NF90_DOUBLE )
441                  ALLOCATE(globaldata_3d_dp(indimlens(dimids(1)),indimlens(dimids(2)),       &
442                     &                      indimlens(dimids(3))))
443                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_3d_dp ) )
444               CASE DEFAULT
445                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
446                  STOP
447            END SELECT
448
449         ELSEIF( ndims == 4 ) THEN
450
451            SELECT CASE( xtype )
452               CASE( NF90_BYTE )
453                  ALLOCATE(globaldata_4d_i1(indimlens(dimids(1)),indimlens(dimids(2)),       &
454                     &                      indimlens(dimids(3)),indimlens(dimids(4))))
455                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_i1 ) )
456               CASE( NF90_SHORT )
457                  ALLOCATE(globaldata_4d_i2(indimlens(dimids(1)),indimlens(dimids(2)),       &
458                     &                      indimlens(dimids(3)),indimlens(dimids(4))))
459                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_i2 ) )
460               CASE( NF90_INT )
461                  ALLOCATE(globaldata_4d_i4(indimlens(dimids(1)),indimlens(dimids(2)),       &
462                     &                      indimlens(dimids(3)),indimlens(dimids(4))))
463                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_i4 ) )
464               CASE( NF90_FLOAT )
465                  ALLOCATE(globaldata_4d_sp(indimlens(dimids(1)),indimlens(dimids(2)),       &
466                     &                      indimlens(dimids(3)),indimlens(dimids(4))))
467                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_sp ) )
468               CASE( NF90_DOUBLE )
469                  ALLOCATE(globaldata_4d_dp(indimlens(dimids(1)),indimlens(dimids(2)),       &
470                     &                      indimlens(dimids(3)),indimlens(dimids(4))))
471                  CALL check_nf90( nf90_get_var( ncid, jv, globaldata_4d_dp ) )
472               CASE DEFAULT
473                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
474                  STOP
475            END SELECT
476
477         ENDIF
478
479      ELSE  ! l_noRebuild = .false.
480
481!3.2.2 For variables that require rebuilding we need to read in from all ndomain files
482!      Here we allocate global variables ahead of looping over files
483         IF (l_verbose) WRITE(numout,*) 'Rebuilding data from variable '//TRIM(varname)//'...'
484
485         IF( ndims == 1 ) THEN
486
487            SELECT CASE( xtype )
488               CASE( NF90_BYTE )
489                  ALLOCATE(globaldata_1d_i1(outdimlens(dimids(1))))
490               CASE( NF90_SHORT )
491                  ALLOCATE(globaldata_1d_i2(outdimlens(dimids(1))))
492               CASE( NF90_INT )
493                  ALLOCATE(globaldata_1d_i4(outdimlens(dimids(1))))
494               CASE( NF90_FLOAT )
495                  ALLOCATE(globaldata_1d_sp(outdimlens(dimids(1))))
496               CASE( NF90_DOUBLE )
497                  ALLOCATE(globaldata_1d_dp(outdimlens(dimids(1))))
498               CASE DEFAULT
499                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
500                  STOP
501            END SELECT
502
503         ELSEIF( ndims == 2 ) THEN
504
505            SELECT CASE( xtype )
506               CASE( NF90_BYTE )
507                  ALLOCATE(globaldata_2d_i1(outdimlens(dimids(1)),outdimlens(dimids(2))))
508               CASE( NF90_SHORT )
509                  ALLOCATE(globaldata_2d_i2(outdimlens(dimids(1)),outdimlens(dimids(2))))
510               CASE( NF90_INT )
511                  ALLOCATE(globaldata_2d_i4(outdimlens(dimids(1)),outdimlens(dimids(2))))
512               CASE( NF90_FLOAT )
513                  ALLOCATE(globaldata_2d_sp(outdimlens(dimids(1)),outdimlens(dimids(2))))
514               CASE( NF90_DOUBLE )
515                  ALLOCATE(globaldata_2d_dp(outdimlens(dimids(1)),outdimlens(dimids(2))))
516               CASE DEFAULT
517                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
518                  STOP
519            END SELECT
520
521         ELSEIF( ndims == 3 ) THEN
522
523            SELECT CASE( xtype )
524               CASE( NF90_BYTE )
525                  ALLOCATE(globaldata_3d_i1(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
526                     &                      outdimlens(dimids(3))))
527               CASE( NF90_SHORT )
528                  ALLOCATE(globaldata_3d_i2(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
529                     &                      outdimlens(dimids(3))))
530               CASE( NF90_INT )
531                  ALLOCATE(globaldata_3d_i4(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
532                     &                      outdimlens(dimids(3))))
533               CASE( NF90_FLOAT )
534                  ALLOCATE(globaldata_3d_sp(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
535                     &                      outdimlens(dimids(3))))
536               CASE( NF90_DOUBLE )
537                  ALLOCATE(globaldata_3d_dp(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
538                     &                      outdimlens(dimids(3))))
539               CASE DEFAULT
540                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
541                  STOP
542            END SELECT
543
544         ELSEIF( ndims == 4 ) THEN
545
546            SELECT CASE( xtype )
547               CASE( NF90_BYTE )
548                  ALLOCATE(globaldata_4d_i1(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
549                     &                      outdimlens(dimids(3)),outdimlens(dimids(4))))
550               CASE( NF90_SHORT )
551                  ALLOCATE(globaldata_4d_i2(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
552                     &                      outdimlens(dimids(3)),outdimlens(dimids(4))))
553               CASE( NF90_INT )
554                  ALLOCATE(globaldata_4d_i4(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
555                     &                      outdimlens(dimids(3)),outdimlens(dimids(4))))
556               CASE( NF90_FLOAT )
557                  ALLOCATE(globaldata_4d_sp(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
558                     &                      outdimlens(dimids(3)),outdimlens(dimids(4))))
559               CASE( NF90_DOUBLE )
560                  ALLOCATE(globaldata_4d_dp(outdimlens(dimids(1)),outdimlens(dimids(2)),     &
561                     &                      outdimlens(dimids(3)),outdimlens(dimids(4))))
562               CASE DEFAULT
563                  WRITE(numerr,*) 'Unknown nf90 type: ', xtype
564                  STOP
565            END SELECT
566         ELSE
567            WRITE(numerr,*) 'ERROR! : A netcdf variable has more than 4 dimensions which is not taken into account'
568            STOP
569         ENDIF
570
571!$OMP  PARALLEL DO DEFAULT(NONE)                                                          &
572!$OMP& PRIVATE(ifile,ncid,xtype,start_pos,local_sizes,InMin,InMax,natts,                  &
573!$OMP&         ndims,attid,attname,dimids,idim,dimname,dimlen,unlimitedDimId,             &
574!$OMP&         halo_start,halo_end,idomain,jdomain,rdomain,di,dj,dr,                      &
575!$OMP&         localdata_1d_i2,localdata_1d_i4,localdata_1d_sp,localdata_1d_dp,           &
576!$OMP&         localdata_2d_i2,localdata_2d_i4,localdata_2d_sp,localdata_2d_dp,           &
577!$OMP&         localdata_3d_i2,localdata_3d_i4,localdata_3d_sp,localdata_3d_dp,           &
578!$OMP&         localdata_4d_i2,localdata_4d_i4,localdata_4d_sp,localdata_4d_dp,           &
579!$OMP&         localdata_1d_i1,localdata_2d_i1,localdata_3d_i1,localdata_4d_i1)           &
580!$OMP& SHARED(jv,nvars,varname,filenames,ValMin,ValMax,indimlens,outdimlens,rbdims,       &
581!$OMP&        ndomain,outid,chunksize,istop,l_valid,nthreads,inncids,rebuild_dims,        &
582!$OMP&        globaldata_1d_i2,globaldata_1d_i4,globaldata_1d_sp,globaldata_1d_dp,        &
583!$OMP&        globaldata_2d_i2,globaldata_2d_i4,globaldata_2d_sp,globaldata_2d_dp,        &
584!$OMP&        globaldata_3d_i2,globaldata_3d_i4,globaldata_3d_sp,globaldata_3d_dp,        &
585!$OMP&        globaldata_4d_i2,globaldata_4d_i4,globaldata_4d_sp,globaldata_4d_dp,        &
586!$OMP&        globaldata_1d_i1,globaldata_2d_i1,globaldata_3d_i1,globaldata_4d_i1)
587
588         DO ifile = 1, ndomain
589
590            ncid = inncids(ifile)
591!$OMP CRITICAL
592            CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_size_local', local_sizes ), istop )
593            CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_position_first', start_pos ), istop )
594            CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_halo_size_start', halo_start ), istop )
595            CALL check_nf90( nf90_get_att( ncid, nf90_global, 'DOMAIN_halo_size_end', halo_end ), istop )
596            CALL check_nf90( nf90_inquire_variable( ncid, jv, varname, xtype, ndims, dimids, natts ), istop )
597            CALL check_nf90( nf90_inquire( ncid, unlimitedDimId = unlimitedDimId ), istop )
598!$OMP END CRITICAL
599
600            ! set defaults for rebuilding so that i is 1st, j 2nd
601            di=1
602            dj=2
603
604            IF( rbdims == 1 ) THEN
605               ! override defaults above and set other variables
606               start_pos(2) = 1
607               local_sizes(2) = outdimlens(3-dimids(2))
608               halo_end(2) = 0
609               halo_start(2) = 0
610               di=rebuild_dims(1)
611               dj=3-di
612            ENDIF
613
614!3.3.1 Generate local domain interior sizes from local_sizes and halo sizes
615!      idomain defines the 1st and last interior points in the i direction and
616!      jdomain defines the 1st and last interior points in the j direction
617
618            idomain(1) = 1 + halo_start(di)
619            idomain(2) = local_sizes(di) - halo_end(di)
620            jdomain(1) = 1 + halo_start(dj)
621            jdomain(2) = local_sizes(dj) - halo_end(dj)
622
623!3.3.2 For rbdims or more dimensions put the data array from this input file into the correct
624!      part of the output data array. Assume the first dimensions are those to be rebuilt.
625
626            IF( ndims == 1 ) THEN
627
628               IF( rebuild_dims(1) == 1 ) THEN
629                  dr = di
630                  rdomain = idomain
631               ELSE
632                  dr = dj
633                  rdomain = jdomain
634               ENDIF
635
636              SELECT CASE( xtype )
637                  CASE( NF90_BYTE )
638                     ALLOCATE(localdata_1d_i1(local_sizes(dr)))
639                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_i1 ), istop )
640                     DO jr = rdomain(1), rdomain(2)
641                        globaldata_1d_i1(start_pos(dr) + jr - 1) = localdata_1d_i1(jr)
642                     END DO
643                     DEALLOCATE(localdata_1d_i1)
644                  CASE( NF90_SHORT )
645                     ALLOCATE(localdata_1d_i2(local_sizes(dr)))
646                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_i2 ), istop )
647                     DO jr = rdomain(1), rdomain(2)
648                        globaldata_1d_i2(start_pos(dr) + jr - 1) = localdata_1d_i2(jr)
649                     END DO
650                     DEALLOCATE(localdata_1d_i2)
651                  CASE( NF90_INT )
652                     ALLOCATE(localdata_1d_i4(local_sizes(dr)))
653                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_i4 ), istop )
654                     DO jr = rdomain(1), rdomain(2)
655                        globaldata_1d_i4(start_pos(dr) + jr - 1) = localdata_1d_i4(jr)
656                     END DO
657                     DEALLOCATE(localdata_1d_i4)
658                  CASE( NF90_FLOAT )
659                     ALLOCATE(localdata_1d_sp(local_sizes(dr)))
660                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_sp ), istop )
661                     DO jr = rdomain(1), rdomain(2)
662                        globaldata_1d_sp(start_pos(dr) + jr - 1) = localdata_1d_sp(jr)
663                     END DO
664                     DEALLOCATE(localdata_1d_sp)
665                  CASE( NF90_DOUBLE )
666                     ALLOCATE(localdata_1d_dp(local_sizes(dr)))
667                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_1d_dp ), istop )
668                     DO jr = rdomain(1), rdomain(2)
669                        globaldata_1d_dp(start_pos(dr) + jr - 1) = localdata_1d_dp(jr)
670                     END DO
671                     DEALLOCATE(localdata_1d_dp)
672                  CASE DEFAULT
673                     WRITE(numerr,*) 'Unknown nf90 type: ', xtype
674                     istop = istop + 1
675               END SELECT
676
677            ELSEIF( ndims == 2 ) THEN
678       
679               SELECT CASE( xtype )
680                  CASE( NF90_BYTE )
681                     ALLOCATE(localdata_2d_i1(local_sizes(di),local_sizes(dj)))
682                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_i1 ), istop )
683                        DO jj = jdomain(1), jdomain(2)
684                           DO ji = idomain(1), idomain(2)
685                              globaldata_2d_i1(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_i1(ji,jj)
686                        END DO
687                     END DO
688                     DEALLOCATE(localdata_2d_i1)
689                  CASE( NF90_SHORT )
690                     ALLOCATE(localdata_2d_i2(local_sizes(di),local_sizes(dj)))
691                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_i2 ), istop )
692                     DO jj = jdomain(1), jdomain(2)
693                        DO ji = idomain(1), idomain(2)
694                           globaldata_2d_i2(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_i2(ji,jj)
695                        END DO
696                     END DO
697                     DEALLOCATE(localdata_2d_i2)
698                  CASE( NF90_INT )
699                     ALLOCATE(localdata_2d_i4(local_sizes(di),local_sizes(dj)))
700                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_i4 ), istop )
701                     DO jj = jdomain(1), jdomain(2)
702                        DO ji = idomain(1), idomain(2)
703                           globaldata_2d_i4(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_i4(ji,jj)
704                        END DO
705                     END DO
706                     DEALLOCATE(localdata_2d_i4)
707                  CASE( NF90_FLOAT )
708                     ALLOCATE(localdata_2d_sp(local_sizes(di),local_sizes(dj)))
709                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_sp ), istop )
710                     DO jj = jdomain(1), jdomain(2)
711                        DO ji = idomain(1), idomain(2)
712                           globaldata_2d_sp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_sp(ji,jj)
713                        END DO
714                     END DO
715                     DEALLOCATE(localdata_2d_sp)
716                  CASE( NF90_DOUBLE )
717                     ALLOCATE(localdata_2d_dp(local_sizes(di),local_sizes(dj)))
718                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_2d_dp ), istop )
719                     DO jj = jdomain(1), jdomain(2)
720                        DO ji = idomain(1), idomain(2)
721                           globaldata_2d_dp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1) = localdata_2d_dp(ji,jj)
722                        END DO
723                     END DO
724                     DEALLOCATE(localdata_2d_dp) 
725                  CASE DEFAULT
726                     WRITE(numerr,*) 'Unknown nf90 type: ', xtype
727                     istop = istop + 1
728               END SELECT
729
730            ELSEIF( ndims == 3 ) THEN
731
732               SELECT CASE( xtype )
733                  CASE( NF90_BYTE )
734                     ALLOCATE(localdata_3d_i1(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
735                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_i1 ), istop )
736!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
737!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_i1,localdata_3d_i1,di,dj)
738                     DO jk = 1, indimlens(dimids(3))
739                        DO jj = jdomain(1), jdomain(2)
740                           DO ji = idomain(1), idomain(2)
741                              globaldata_3d_i1(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_i1(ji,jj,jk)
742                           END DO
743                        END DO
744                     END DO
745!$OMP END PARALLEL DO
746                     DEALLOCATE(localdata_3d_i1)
747                  CASE( NF90_SHORT )
748                     ALLOCATE(localdata_3d_i2(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
749                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_i2 ), istop )
750!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
751!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_i2,localdata_3d_i2,di,dj)
752                     DO jk = 1, indimlens(dimids(3))
753                        DO jj = jdomain(1), jdomain(2)
754                           DO ji = idomain(1), idomain(2)
755                              globaldata_3d_i2(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_i2(ji,jj,jk)
756                           END DO
757                        END DO
758                     END DO
759!$OMP END PARALLEL DO
760                     DEALLOCATE(localdata_3d_i2)
761                  CASE( NF90_INT )
762                     ALLOCATE(localdata_3d_i4(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
763                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_i4 ), istop )
764!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
765!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_i4,localdata_3d_i4,di,dj)
766                     DO jk = 1, indimlens(dimids(3))
767                        DO jj = jdomain(1), jdomain(2)
768                           DO ji = idomain(1), idomain(2)
769                              globaldata_3d_i4(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_i4(ji,jj,jk)
770                           END DO
771                        END DO
772                     END DO
773!$OMP END PARALLEL DO
774                     DEALLOCATE(localdata_3d_i4)
775                  CASE( NF90_FLOAT )
776                     ALLOCATE(localdata_3d_sp(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
777                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_sp ), istop )
778!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
779!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_sp,localdata_3d_sp,di,dj)
780                     DO jk = 1, indimlens(dimids(3))
781                        DO jj = jdomain(1), jdomain(2)
782                           DO ji = idomain(1), idomain(2)
783                              globaldata_3d_sp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_sp(ji,jj,jk)
784                           END DO
785                        END DO
786                     END DO
787!$OMP END PARALLEL DO
788                     DEALLOCATE(localdata_3d_sp)
789                  CASE( NF90_DOUBLE )
790                     ALLOCATE(localdata_3d_dp(local_sizes(di),local_sizes(dj),indimlens(dimids(3))))
791                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_3d_dp ), istop )
792!$OMP  PARALLEL DO DEFAULT(NONE) PRIVATE(ji,jj,jk)   &
793!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_3d_dp,localdata_3d_dp,di,dj)
794                     DO jk = 1, indimlens(dimids(3))
795                        DO jj = jdomain(1), jdomain(2)
796                           DO ji = idomain(1), idomain(2)
797                              globaldata_3d_dp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk) = localdata_3d_dp(ji,jj,jk)
798                           END DO
799                        END DO
800                     END DO
801!$OMP END PARALLEL DO
802                     DEALLOCATE(localdata_3d_dp)
803                  CASE DEFAULT
804                     WRITE(numerr,*) 'Unknown nf90 type: ', xtype
805                     istop = istop + 1
806               END SELECT
807     
808            ELSEIF (ndims == 4) THEN
809       
810               SELECT CASE( xtype )
811                  CASE( NF90_BYTE )
812                     ALLOCATE(localdata_4d_i1(local_sizes(di),local_sizes(dj),               &
813                         &                     indimlens(dimids(3)),indimlens(dimids(4))))
814                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_i1 ), istop )
815!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
816!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_i1,localdata_4d_i1,di,dj)
817                     DO jl = 1, indimlens(dimids(4))
818!$OMP DO
819                        DO jk = 1, indimlens(dimids(3))
820                           DO jj = jdomain(1), jdomain(2)
821                              DO ji = idomain(1), idomain(2)
822                                 globaldata_4d_i1(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_i1(ji,jj,jk,jl)
823                              END DO
824                           END DO
825                        END DO
826!$OMP END DO nowait
827                     END DO
828!$OMP END PARALLEL
829                     DEALLOCATE(localdata_4d_i1)
830                  CASE( NF90_SHORT )
831                     ALLOCATE(localdata_4d_i2(local_sizes(di),local_sizes(dj),               &
832                        &                     indimlens(dimids(3)),indimlens(dimids(4))))
833                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_i2 ), istop )
834!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
835!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_i2,localdata_4d_i2,di,dj)
836                     DO jl = 1, indimlens(dimids(4))
837!$OMP DO
838                        DO jk = 1, indimlens(dimids(3))
839                           DO jj = jdomain(1), jdomain(2) 
840                              DO ji = idomain(1), idomain(2)
841                                 globaldata_4d_i2(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_i2(ji,jj,jk,jl)
842                              END DO
843                           END DO
844                        END DO
845!$OMP END DO nowait
846                     END DO
847!$OMP END PARALLEL
848                     DEALLOCATE(localdata_4d_i2)
849                  CASE( NF90_INT )
850                     ALLOCATE(localdata_4d_i4(local_sizes(di),local_sizes(dj),               &
851                        &                     indimlens(dimids(3)),indimlens(dimids(4))))
852                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_i4 ), istop )
853!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
854!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_i4,localdata_4d_i4,di,dj)
855                     DO jl = 1, indimlens(dimids(4))
856!$OMP DO
857                        DO jk = 1, indimlens(dimids(3))
858                           DO jj = jdomain(1), jdomain(2) 
859                              DO ji = idomain(1), idomain(2)
860                                 globaldata_4d_i4(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_i4(ji,jj,jk,jl)
861                              END DO
862                           END DO
863                        END DO
864!$OMP END DO nowait
865                     END DO
866!$OMP END PARALLEL
867                     DEALLOCATE(localdata_4d_i4)
868                  CASE( NF90_FLOAT )
869                     ALLOCATE(localdata_4d_sp(local_sizes(di),local_sizes(dj),               &
870                        &                     indimlens(dimids(3)),indimlens(dimids(4))))
871                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_sp ), istop )
872!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
873!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_sp,localdata_4d_sp,di,dj)
874                     DO jl = 1, indimlens(dimids(4))
875!$OMP DO
876                        DO jk = 1, indimlens(dimids(3))
877                           DO jj = jdomain(1), jdomain(2) 
878                              DO ji = idomain(1), idomain(2)
879                                 globaldata_4d_sp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_sp(ji,jj,jk,jl)
880                              END DO
881                           END DO
882                        END DO
883!$OMP END DO nowait
884                     END DO
885!$OMP END PARALLEL
886                     DEALLOCATE(localdata_4d_sp)
887                  CASE( NF90_DOUBLE )
888                     ALLOCATE(localdata_4d_dp(local_sizes(di),local_sizes(dj),               &
889                        &                     indimlens(dimids(3)),indimlens(dimids(4))))
890                     CALL check_nf90( nf90_get_var( ncid, jv, localdata_4d_dp ), istop )
891!$OMP  PARALLEL DEFAULT(NONE) PRIVATE(ji,jj,jk,jl)   &
892!$OMP& SHARED(idomain,jdomain,indimlens,dimids,start_pos,globaldata_4d_dp,localdata_4d_dp,di,dj)
893                     DO jl = 1, indimlens(dimids(4))
894!$OMP DO
895                        DO jk = 1, indimlens(dimids(3))
896                           DO jj = jdomain(1), jdomain(2) 
897                              DO ji = idomain(1), idomain(2)
898                                 globaldata_4d_dp(start_pos(di) + ji - 1, start_pos(dj) + jj - 1, jk, jl) = localdata_4d_dp(ji,jj,jk,jl)
899                              END DO
900                           END DO
901                        END DO
902!$OMP END DO nowait
903                     END DO
904!$OMP END PARALLEL
905                     DEALLOCATE(localdata_4d_dp)
906                  CASE DEFAULT
907                     WRITE(numerr,*) 'Unknown nf90 type: ', xtype
908                     istop = istop + 1
909               END SELECT
910
911            ENDIF ! l_noRebuild false
912
913!3.4 Work out if the valid_min and valid_max attributes exist for this variable.
914!    If they do then calculate the extrema over all input files.
915
916            DO attid = 1, natts
917               CALL check_nf90( nf90_inq_attname( ncid, jv, attid, attname ), istop )
918               IF( INDEX( attname, "valid_min" ) == 1 ) THEN
919                  CALL check_nf90( nf90_get_att( ncid, jv, attname, InMin), istop )
920                  l_valid = .true.
921               ENDIF
922               IF( INDEX( attname, "valid_max" ) == 1 ) THEN
923                  CALL check_nf90( nf90_get_att( ncid, jv, attname, InMax ), istop )
924                  l_valid = .true.
925               ENDIF
926            END DO
927
928            IF (l_valid) THEN
929!$OMP CRITICAL
930               IF( InMin < ValMin ) ValMin = InMin
931               IF( InMax > ValMax ) ValMax = InMax
932!$OMP END CRITICAL
933            ENDIF
934
935!3.5 Abort if failure and only 1 thread
936
937            IF( nthreads == 1 .AND. istop /= nf90_noerr )  THEN
938               WRITE(numerr,*) '*** NEMO rebuild failed! ***'
939               STOP
940            ENDIF
941       
942         END DO  ! loop over files
943!$OMP END PARALLEL DO
944
945!3.6 Abort if any of the OMP threads failed
946         IF( istop /= nf90_noerr )  THEN
947            WRITE(numerr,*) '*** NEMO rebuild failed! ***'
948            STOP
949         ENDIF
950
951      ENDIF ! ndims > 2
952
953!---------------------------------------------------------------------------
954!4. Write data to output file
955
956      IF (l_verbose) WRITE(numout,*) 'Writing variable '//TRIM(varname)//'...'
957
958!4.1 If the valid min and max attributes exist then update them in the file
959
960      IF( l_valid ) THEN
961         CALL check_nf90( nf90_put_att( outid, jv, "valid_min", ValMin ) )
962         CALL check_nf90( nf90_put_att( outid, jv, "valid_max", ValMax ) )
963      ENDIF
964
965!4.2 Write the data to the output file depending on how many dimensions it has
966
967      IF( ndims == 0 ) THEN
968
969         SELECT CASE( xtype )
970            CASE( NF90_BYTE )
971               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_i1 ) )
972            CASE( NF90_SHORT )
973               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_i2 ) )
974            CASE( NF90_INT )
975               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_i4 ) )
976            CASE( NF90_FLOAT )
977               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_sp ) )
978            CASE( NF90_DOUBLE )
979               CALL check_nf90( nf90_put_var( outid, jv, globaldata_0d_dp ) )
980         END SELECT
981
982      ELSEIF( ndims == 1 ) THEN
983
984         SELECT CASE( xtype )
985            CASE( NF90_BYTE )
986               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_i1 ) )
987               DEALLOCATE(globaldata_1d_i1)
988            CASE( NF90_SHORT )
989               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_i2 ) )
990               DEALLOCATE(globaldata_1d_i2)
991            CASE( NF90_INT )
992               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_i4 ) )
993               DEALLOCATE(globaldata_1d_i4)
994            CASE( NF90_FLOAT )
995               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_sp ) )
996               DEALLOCATE(globaldata_1d_sp)
997            CASE( NF90_DOUBLE )
998               CALL check_nf90( nf90_put_var( outid, jv, globaldata_1d_dp ) )
999               DEALLOCATE(globaldata_1d_dp)
1000         END SELECT
1001
1002      ELSEIF( ndims == 2 ) THEN 
1003     
1004         SELECT CASE( xtype )   
1005            CASE( NF90_BYTE )                   
1006               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_i1 ) )
1007               DEALLOCATE(globaldata_2d_i1)
1008            CASE( NF90_SHORT )                   
1009               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_i2 ) )
1010               DEALLOCATE(globaldata_2d_i2)
1011            CASE( NF90_INT )                             
1012               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_i4 ) )
1013               DEALLOCATE(globaldata_2d_i4)
1014            CASE( NF90_FLOAT )                             
1015               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_sp ) )
1016               DEALLOCATE(globaldata_2d_sp)
1017            CASE( NF90_DOUBLE )                                         
1018               CALL check_nf90( nf90_put_var( outid, jv, globaldata_2d_dp ) )
1019               DEALLOCATE(globaldata_2d_dp)
1020            CASE DEFAULT   
1021               WRITE(numerr,*) 'Unknown nf90 type: ', xtype
1022               STOP
1023         END SELECT     
1024                     
1025      ELSEIF( ndims == 3 ) THEN
1026     
1027         SELECT CASE( xtype ) 
1028            CASE( NF90_BYTE )                   
1029               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_i1 ) )
1030               DEALLOCATE(globaldata_3d_i1)
1031            CASE( NF90_SHORT )                   
1032               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_i2 ) )
1033               DEALLOCATE(globaldata_3d_i2)
1034            CASE( NF90_INT )                             
1035               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_i4 ) )
1036               DEALLOCATE(globaldata_3d_i4)
1037            CASE( NF90_FLOAT )                             
1038               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_sp ) )
1039               DEALLOCATE(globaldata_3d_sp)
1040            CASE( NF90_DOUBLE )                                         
1041               CALL check_nf90( nf90_put_var( outid, jv, globaldata_3d_dp ) )
1042               DEALLOCATE(globaldata_3d_dp)
1043            CASE DEFAULT   
1044               WRITE(numerr,*) 'Unknown nf90 type: ', xtype
1045               STOP
1046         END SELECT     
1047   
1048      ELSEIF( ndims == 4 ) THEN
1049     
1050         SELECT CASE( xtype )   
1051            CASE( NF90_BYTE )                   
1052               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_i1 ) )
1053               DEALLOCATE(globaldata_4d_i1)
1054            CASE( NF90_SHORT )                   
1055               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_i2 ) )
1056               DEALLOCATE(globaldata_4d_i2)
1057            CASE( NF90_INT )                             
1058               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_i4 ) )
1059               DEALLOCATE(globaldata_4d_i4)
1060            CASE( NF90_FLOAT )                             
1061               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_sp ) )
1062               DEALLOCATE(globaldata_4d_sp)
1063            CASE( NF90_DOUBLE )                                         
1064               CALL check_nf90( nf90_put_var( outid, jv, globaldata_4d_dp ) )
1065               DEALLOCATE(globaldata_4d_dp)
1066            CASE DEFAULT   
1067               WRITE(numerr,*) 'Unknown nf90 type: ', xtype
1068               STOP
1069         END SELECT     
1070   
1071      ENDIF
1072   
1073   END DO  ! loop over variables
1074
1075!---------------------------------------------------------------------------
1076!5. Close files
1077
1078!5.1 Close input files
1079
1080   IF (l_verbose) WRITE(numout,*) 'Closing input files...'
1081   DO ifile = 1, ndomain
1082      ncid = inncids(ifile)
1083      CALL check_nf90( nf90_close( ncid ) )
1084   END DO
1085
1086!5.2 Close output file
1087
1088   IF (l_verbose) WRITE(numout,*) 'Closing output file...'
1089   CALL check_nf90( nf90_close( outid ) )
1090
1091   IF (l_verbose) WRITE(numout,*) 'NEMO rebuild completed successfully'
1092   IF (l_verbose) WRITE(numout,*)
1093
1094CONTAINS
1095
1096
1097   SUBROUTINE check_nf90(status, errorFlag)
1098   !---------------------------------------------------------------------
1099   !  Checks return code from nf90 library calls and warns if needed
1100   !  If errorFlag is present then it just increments this flag (OMP use)
1101   !
1102   !---------------------------------------------------------------------
1103      INTEGER, INTENT(IN   ) :: status
1104      INTEGER, INTENT(INOUT), OPTIONAL :: errorFlag
1105   !---------------------------------------------------------------------
1106
1107      IF( status /= nf90_noerr ) THEN
1108         WRITE(numerr,*) 'ERROR! : '//TRIM(nf90_strerror(status))
1109         IF( PRESENT( errorFlag ) ) THEN
1110            errorFlag = errorFlag + status
1111         ELSE
1112            WRITE(numerr,*) "*** NEMO rebuild failed ***"
1113            WRITE(numerr,*)
1114            STOP
1115         ENDIF
1116      ENDIF
1117
1118   END SUBROUTINE check_nf90
1119
1120
1121END PROGRAM rebuild_nemo
Note: See TracBrowser for help on using the repository browser.