[3182] | 1 | !> @file src/get_netcdf_variable.inc |
---|
| 2 | !------------------------------------------------------------------------------! |
---|
| 3 | ! This file is part of the PALM model system. |
---|
| 4 | ! |
---|
| 5 | ! PALM is free software: you can redistribute it and/or modify it under the |
---|
| 6 | ! terms of the GNU General Public License as published by the Free Software |
---|
| 7 | ! Foundation, either version 3 of the License, or (at your option) any later |
---|
| 8 | ! version. |
---|
| 9 | ! |
---|
| 10 | ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY |
---|
| 11 | ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR |
---|
| 12 | ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. |
---|
| 13 | ! |
---|
| 14 | ! You should have received a copy of the GNU General Public License along with |
---|
| 15 | ! PALM. If not, see <http://www.gnu.org/licenses/>. |
---|
| 16 | ! |
---|
| 17 | ! Copyright 2017-2018 Leibniz Universitaet Hannover |
---|
| 18 | ! Copyright 2017-2018 Deutscher Wetterdienst Offenbach |
---|
| 19 | !------------------------------------------------------------------------------! |
---|
| 20 | ! |
---|
| 21 | ! Current revisions: |
---|
| 22 | ! ----------------- |
---|
| 23 | ! |
---|
| 24 | ! |
---|
| 25 | ! Former revisions: |
---|
| 26 | ! ----------------- |
---|
| 27 | ! $Id$ |
---|
[3183] | 28 | ! Initial revision |
---|
| 29 | ! |
---|
| 30 | ! |
---|
| 31 | ! |
---|
[3182] | 32 | ! |
---|
| 33 | ! |
---|
| 34 | ! |
---|
| 35 | ! Authors: |
---|
| 36 | ! -------- |
---|
| 37 | ! @author Eckhard Kadasch |
---|
| 38 | ! |
---|
| 39 | ! Description: |
---|
| 40 | ! ------------ |
---|
| 41 | !> This file provides the implementation of the io.f90 subroutines |
---|
| 42 | !> get_netcdf_variable_int() and get_netcdf_variable_real(). Both routines |
---|
| 43 | !> differ only in the type of the netCDF variable to be read and otherwise share |
---|
| 44 | !> their functionality. In order to avoid duplicating code, both routines source |
---|
| 45 | !> the implementation below using Fortran's INCLUDE statement. |
---|
| 46 | !> |
---|
| 47 | !------------------------------------------------------------------------------! |
---|
| 48 | INTEGER :: i, ncid, start(3) |
---|
| 49 | |
---|
| 50 | ! Read in_var NetCDF attributes |
---|
| 51 | IF ( nf90_open( TRIM(in_file), NF90_NOWRITE, ncid ) .EQ. NF90_NOERR .AND. & |
---|
| 52 | nf90_inq_varid( ncid, in_var % name, in_var % varid ) .EQ. NF90_NOERR ) THEN |
---|
| 53 | |
---|
| 54 | CALL check(nf90_get_att(ncid, in_var % varid, "long_name", in_var % long_name)) |
---|
| 55 | CALL check(nf90_get_att(ncid, in_var % varid, "units", in_var % units)) |
---|
| 56 | |
---|
| 57 | ! Read in_var NetCDF dimensions |
---|
| 58 | CALL check(nf90_inquire_variable( ncid, in_var % varid, & |
---|
| 59 | ndims = in_var % ndim, & |
---|
| 60 | dimids = in_var % dimids )) |
---|
| 61 | |
---|
| 62 | DO i = 1, in_var % ndim |
---|
| 63 | CALL check(nf90_inquire_dimension( ncid, in_var % dimids(i), & |
---|
| 64 | name = in_var % dimname(i), & |
---|
| 65 | len = in_var % dimlen(i) )) |
---|
| 66 | END DO |
---|
| 67 | |
---|
| 68 | start = (/ 1, 1, 1 /) |
---|
| 69 | IF ( TRIM(in_var % name) .EQ. 'T_SO' ) THEN |
---|
| 70 | ! Skip depth = 0.0 for T_SO and reduce number of depths from 9 to 8 |
---|
| 71 | in_var % dimlen(3) = in_var % dimlen(3) - 1 |
---|
| 72 | |
---|
| 73 | ! Start reading from second level, e.g. depth = 0.005 instead of 0.0 |
---|
| 74 | start(3) = 2 |
---|
| 75 | END IF |
---|
| 76 | |
---|
| 77 | SELECT CASE(in_var % ndim) |
---|
| 78 | |
---|
| 79 | CASE (2) |
---|
| 80 | |
---|
| 81 | ALLOCATE( buffer( in_var % dimlen(1), & |
---|
| 82 | in_var % dimlen(2), & |
---|
| 83 | 1 ) ) |
---|
| 84 | |
---|
| 85 | CASE (3, 4) |
---|
| 86 | |
---|
| 87 | ALLOCATE( buffer( in_var % dimlen(1), & |
---|
| 88 | in_var % dimlen(2), & |
---|
| 89 | in_var % dimlen(3) ) ) |
---|
| 90 | |
---|
| 91 | CASE DEFAULT |
---|
| 92 | |
---|
| 93 | message = "Failed reading NetCDF variable " // & |
---|
| 94 | TRIM(in_var % name) // " with " // TRIM(str(in_var%ndim)) // & |
---|
| 95 | " dimensions because only two- and and three-dimensional" // & |
---|
| 96 | " variables are supported." |
---|
| 97 | CALL abort('get_netcdf_variable', message) |
---|
| 98 | |
---|
| 99 | END SELECT |
---|
| 100 | CALL run_control('time', 'alloc') |
---|
| 101 | |
---|
| 102 | ! TODO: Check for matching dimensions of buffer and var |
---|
| 103 | CALL check(nf90_get_var( ncid, in_var % varid, buffer, & |
---|
| 104 | start = start, & |
---|
| 105 | count = in_var % dimlen(1:3) ) ) |
---|
| 106 | |
---|
| 107 | CALL run_control('time', 'read') |
---|
| 108 | ELSE |
---|
| 109 | |
---|
| 110 | message = "Failed to read '" // TRIM(in_var % name) // & |
---|
| 111 | "' from file '" // TRIM(in_file) // "'." |
---|
| 112 | CALL report('get_netcdf_variable', message) |
---|
| 113 | |
---|
| 114 | END IF |
---|
| 115 | |
---|
| 116 | CALL check(nf90_close(ncid)) |
---|
| 117 | |
---|
| 118 | CALL run_control('time', 'read') |
---|