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