ROMS
Loading...
Searching...
No Matches
mod_netcdf::netcdf_get_fatt Interface Reference

Public Member Functions

subroutine netcdf_get_fatt_dp (ng, model, ncname, varid, attname, attvalue, foundit, ncid)
 
subroutine netcdf_get_fatt_r8 (ng, model, ncname, varid, attname, attvalue, foundit, ncid)
 

Detailed Description

Definition at line 38 of file mod_netcdf.F.

Member Function/Subroutine Documentation

◆ netcdf_get_fatt_dp()

subroutine mod_netcdf::netcdf_get_fatt::netcdf_get_fatt_dp ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
integer, intent(in) varid,
character (len=*), dimension(:), intent(in) attname,
real(dp), dimension(:), intent(out) attvalue,
logical, dimension(:), intent(out) foundit,
integer, intent(in), optional ncid )

Definition at line 1968 of file mod_netcdf.F.

1970!
1971!=======================================================================
1972! !
1973! This routine gets requested variable double-precision attribute(s). !
1974! !
1975! On Input: !
1976! !
1977! ng Nested grid number (integer) !
1978! model Calling model identifier (integer) !
1979! ncname NetCDF file name (string) !
1980! varid Variable ID for variable attribute or !
1981! NF90_GLOBAL for a global attribute (integer) !
1982! AttName Attribute name to read (string array) !
1983! ncid NetCDF file ID (integer, OPTIONAL) !
1984! !
1985! On Ouput: !
1986! !
1987! AttValue Attribute value (double precision array) !
1988! foundit Switch (T/F) activated when the requested !
1989! attribute is found (logical array) !
1990! !
1991!=======================================================================
1992!
1993! Imported variable declarations.
1994!
1995 integer, intent(in) :: ng, model, varid
1996
1997 integer, intent(in), optional :: ncid
1998!
1999 character (len=*), intent(in) :: ncname
2000 character (len=*), intent(in) :: AttName(:)
2001!
2002 logical, intent(out) :: foundit(:)
2003!
2004 real(dp), intent(out) :: AttValue(:)
2005!
2006! Local variable declarations.
2007!
2008 integer :: i, j, my_natts, my_ncid, natts, status
2009
2010# if !defined PARALLEL_IO && defined DISTRIBUTE
2011!
2012 real(dp), allocatable :: rbuffer(:)
2013# endif
2014!
2015 character (len=40) :: my_Aname
2016 character (len=40) :: my_Vname
2017
2018 character (len=*), parameter :: MyFile = &
2019 & __FILE__//", netcdf_get_fatt"
2020!
2021!-----------------------------------------------------------------------
2022! Inquire ID of requested variable.
2023!-----------------------------------------------------------------------
2024!
2025! Get number of variable attributes to process.
2026!
2027 natts=ubound(attname, dim=1)
2028 DO i=1,natts
2029 foundit(i)=.false.
2030 attvalue(i)=0.0_dp
2031 END DO
2032# if !defined PARALLEL_IO && defined DISTRIBUTE
2033 IF (.not.allocated(rbuffer)) THEN
2034 allocate ( rbuffer(2*natts+1) )
2035 END IF
2036# endif
2037!
2038! If appropriate, open file for reading.
2039!
2040 IF (.not.PRESENT(ncid)) THEN
2041 CALL netcdf_open (ng, model, trim(ncname), 0, my_ncid)
2042 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
2043 ELSE
2044 my_ncid=ncid
2045 END IF
2046!
2047! Inquire about requested attribute value.
2048!
2049 IF (inpthread) THEN
2050 IF (varid.eq.nf90_global) THEN
2051 status=nf90_inquire(my_ncid, &
2052 & nattributes = my_natts)
2053 ELSE
2054 status=nf90_inquire_variable(my_ncid, varid, &
2055 & name = my_vname, &
2056 & natts = my_natts)
2057 END IF
2058 IF (status.eq.nf90_noerr) THEN
2059 DO j=1,my_natts
2060 status=nf90_inq_attname(my_ncid, varid, j, my_aname)
2061 IF (status.eq.nf90_noerr) THEN
2062 DO i=1,natts
2063 IF (trim(my_aname).eq.trim(attname(i))) THEN
2064 status=nf90_get_att(my_ncid, varid, trim(attname(i)), &
2065 & attvalue(i))
2066 IF (founderror(status, nf90_noerr, &
2067 & __line__, myfile)) THEN
2068 IF (master) WRITE (stdout,10) trim(attname(i)), &
2069 & trim(my_vname), &
2070 & trim(ncname), &
2071 & trim(sourcefile), &
2072 & nf90_strerror(status)
2073 exit_flag=2
2074 ioerror=status
2075 END IF
2076 foundit(i)=.true.
2077 EXIT
2078 END IF
2079 END DO
2080 ELSE
2081 IF (master) WRITE (stdout,20) j, &
2082 & trim(my_vname), &
2083 & trim(ncname), &
2084 & trim(sourcefile), &
2085 & nf90_strerror(status)
2086 exit_flag=2
2087 ioerror=status
2088 EXIT
2089 END IF
2090 END DO
2091# if !defined PARALLEL_IO && defined DISTRIBUTE
2092!
2093 rbuffer=0.0_dp
2094 DO i=1,natts
2095 rbuffer(i)=attvalue(i)
2096 IF (foundit(i)) THEN
2097 rbuffer(i+natts)=1.0_dp
2098 END IF
2099 END DO
2100 rbuffer(2*natts+1)=real(ioerror, dp)
2101# endif
2102 ELSE
2103 IF (master) WRITE (stdout,30) trim(my_vname), &
2104 & trim(ncname), &
2105 & trim(sourcefile), &
2106 & nf90_strerror(status)
2107 exit_flag=2
2108 ioerror=status
2109 END IF
2110 END IF
2111
2112# if !defined PARALLEL_IO && defined DISTRIBUTE
2113!
2114! Broadcast information to all processors in the group. Pack all
2115! data into a real array for efficient communications.
2116!
2117 CALL mp_bcasti (ng, model, exit_flag)
2118 IF (exit_flag.eq.noerror) THEN
2119 CALL mp_bcastf (ng, model, rbuffer)
2120 DO i=1,natts
2121 attvalue(i)=rbuffer(i)
2122 IF (rbuffer(i+natts).gt.0.0_dp) THEN
2123 foundit(i)=.true.
2124 END IF
2125 END DO
2126 ioerror=int(rbuffer(2*natts+1))
2127 IF (allocated(rbuffer)) THEN
2128 deallocate (rbuffer)
2129 END IF
2130 END IF
2131# endif
2132!
2133! If applicable, close input NetCDF file.
2134!
2135 IF (.not.PRESENT(ncid)) THEN
2136 CALL netcdf_close (ng, model, my_ncid, ncname, .false.)
2137 END IF
2138!
2139 10 FORMAT (/,' NETCDF_GET_FATT_DP - error while reading attribute:', &
2140 & 1x,a,'for variable',1x,a,/,22x,'in input file:',2x,a,/, &
2141 & 22x,'call from:',2x,a,/,22x,a)
2142 20 FORMAT (/,' NETCDF_GET_FATT_DP - error while inquiring ', &
2143 & 'attribute: ',i2.2,'for variable',1x,a,/,22x, &
2144 & 'in input file:',2x,a,/,19x,'call from:',2x,a,/,19x,a)
2145 30 FORMAT (/,' NETCDF_GET_FATT_DP - error while inquiring number of',&
2146 & ' attributes for variable:',1x,a,/,22x,'in input file:', &
2147 & 2x,a,/,22x,'call from:',2x,a,/,22x,a)
2148!
2149 RETURN

References mod_scalars::exit_flag, mod_parallel::inpthread, mod_iounits::ioerror, mod_parallel::master, mod_netcdf::netcdf_close(), mod_netcdf::netcdf_open(), mod_scalars::noerror, mod_iounits::sourcefile, and mod_iounits::stdout.

Here is the call graph for this function:

◆ netcdf_get_fatt_r8()

subroutine mod_netcdf::netcdf_get_fatt::netcdf_get_fatt_r8 ( integer, intent(in) ng,
integer, intent(in) model,
character (len=*), intent(in) ncname,
integer, intent(in) varid,
character (len=*), dimension(:), intent(in) attname,
real(r8), dimension(:), intent(out) attvalue,
logical, dimension(:), intent(out) foundit,
integer, intent(in), optional ncid )

Definition at line 2153 of file mod_netcdf.F.

2155!
2156!=======================================================================
2157! !
2158! This routine gets requested variable floating-point attribute(s). !
2159! !
2160! On Input: !
2161! !
2162! ng Nested grid number (integer) !
2163! model Calling model identifier (integer) !
2164! ncname NetCDF file name (string) !
2165! varid Variable ID for variable attribute or !
2166! NF90_GLOBAL for a global attribute (integer) !
2167! AttName Attribute name to read (string array) !
2168! ncid NetCDF file ID (integer, OPTIONAL) !
2169! !
2170! On Ouput: !
2171! !
2172! AttValue Attribute value (real array) !
2173! foundit Switch (T/F) activated when the requested !
2174! attribute is found (logical array) !
2175! !
2176!=======================================================================
2177!
2178! Imported variable declarations.
2179!
2180 integer, intent(in) :: ng, model, varid
2181
2182 integer, intent(in), optional :: ncid
2183!
2184 character (len=*), intent(in) :: ncname
2185 character (len=*), intent(in) :: AttName(:)
2186!
2187 logical, intent(out) :: foundit(:)
2188!
2189 real(r8), intent(out) :: AttValue(:)
2190!
2191! Local variable declarations.
2192!
2193 integer :: i, j, my_natts, my_ncid, natts, status
2194
2195#if !defined PARALLEL_IO && defined DISTRIBUTE
2196!
2197 real(r8), allocatable :: rbuffer(:)
2198#endif
2199!
2200 character (len=40) :: my_Aname
2201 character (len=40) :: my_Vname
2202
2203 character (len=*), parameter :: MyFile = &
2204 & __FILE__//", netcdf_get_fatt_r8"
2205!
2206!-----------------------------------------------------------------------
2207! Inquire ID of requested variable.
2208!-----------------------------------------------------------------------
2209!
2210! Get number of variable attributes to process.
2211!
2212 natts=ubound(attname, dim=1)
2213 DO i=1,natts
2214 foundit(i)=.false.
2215 attvalue(i)=0.0_r8
2216 END DO
2217#if !defined PARALLEL_IO && defined DISTRIBUTE
2218 IF (.not.allocated(rbuffer)) THEN
2219 allocate ( rbuffer(2*natts+1) )
2220 END IF
2221#endif
2222!
2223! If appropriate, open file for reading.
2224!
2225 IF (.not.PRESENT(ncid)) THEN
2226 CALL netcdf_open (ng, model, trim(ncname), 0, my_ncid)
2227 IF (founderror(exit_flag, noerror, __line__, myfile)) RETURN
2228 ELSE
2229 my_ncid=ncid
2230 END IF
2231!
2232! Inquire about requested attribute value.
2233!
2234 IF (inpthread) THEN
2235 IF (varid.eq.nf90_global) THEN
2236 status=nf90_inquire(my_ncid, &
2237 & nattributes = my_natts)
2238 ELSE
2239 status=nf90_inquire_variable(my_ncid, varid, &
2240 & name = my_vname, &
2241 & natts = my_natts)
2242 END IF
2243 IF (status.eq.nf90_noerr) THEN
2244 DO j=1,my_natts
2245 status=nf90_inq_attname(my_ncid, varid, j, my_aname)
2246 IF (status.eq.nf90_noerr) THEN
2247 DO i=1,natts
2248 IF (trim(my_aname).eq.trim(attname(i))) THEN
2249 status=nf90_get_att(my_ncid, varid, trim(attname(i)), &
2250 & attvalue(i))
2251 IF (founderror(status, nf90_noerr, &
2252 & __line__, myfile)) THEN
2253 IF (master) WRITE (stdout,10) trim(attname(i)), &
2254 & trim(my_vname), &
2255 & trim(ncname), &
2256 & trim(sourcefile), &
2257 & nf90_strerror(status)
2258 exit_flag=2
2259 ioerror=status
2260 END IF
2261 foundit(i)=.true.
2262 EXIT
2263 END IF
2264 END DO
2265 ELSE
2266 IF (master) WRITE (stdout,20) j, &
2267 & trim(my_vname), &
2268 & trim(ncname), &
2269 & trim(sourcefile), &
2270 & nf90_strerror(status)
2271 exit_flag=2
2272 ioerror=status
2273 EXIT
2274 END IF
2275 END DO
2276#if !defined PARALLEL_IO && defined DISTRIBUTE
2277!
2278 rbuffer=0.0_r8
2279 DO i=1,natts
2280 rbuffer(i)=attvalue(i)
2281 IF (foundit(i)) THEN
2282 rbuffer(i+natts)=1.0_r8
2283 END IF
2284 END DO
2285 rbuffer(2*natts+1)=real(ioerror, r8)
2286#endif
2287 ELSE
2288 IF (master) WRITE (stdout,30) trim(my_vname), &
2289 & trim(ncname), &
2290 & trim(sourcefile), &
2291 & nf90_strerror(status)
2292 exit_flag=2
2293 ioerror=status
2294 END IF
2295 END IF
2296
2297#if !defined PARALLEL_IO && defined DISTRIBUTE
2298!
2299! Broadcast information to all processors in the group. Pack all
2300! data into a real array for efficient communications.
2301!
2302 CALL mp_bcasti (ng, model, exit_flag)
2303 IF (exit_flag.eq.noerror) THEN
2304 CALL mp_bcastf (ng, model, rbuffer)
2305 DO i=1,natts
2306 attvalue(i)=rbuffer(i)
2307 IF (rbuffer(i+natts).gt.0.0_r8) THEN
2308 foundit(i)=.true.
2309 END IF
2310 END DO
2311 ioerror=int(rbuffer(2*natts+1))
2312 IF (allocated(rbuffer)) THEN
2313 deallocate (rbuffer)
2314 END IF
2315 END IF
2316#endif
2317!
2318! If applicable, close input NetCDF file.
2319!
2320 IF (.not.PRESENT(ncid)) THEN
2321 CALL netcdf_close (ng, model, my_ncid, ncname, .false.)
2322 END IF
2323!
2324 10 FORMAT (/,' NETCDF_GET_FATT_R8 - error while reading attribute:', &
2325 & 1x,a,'for variable',1x,a,/,22x,'in input file:',2x,a,/, &
2326 & 22x,'call from:',2x,a,/,22x,a)
2327 20 FORMAT (/,' NETCDF_GET_FATT_R8 - error while inquiring ', &
2328 & 'attribute: ',i2.2,'for variable',1x,a,/,22x, &
2329 & 'in input file:',2x,a,/,19x,'call from:',2x,a,/,19x,a)
2330 30 FORMAT (/,' NETCDF_GET_FATT_R8 - error while inquiring number of',&
2331 & ' attributes for variable:',1x,a,/,22x,'in input file:', &
2332 & 2x,a,/,22x,'call from:',2x,a,/,22x,a)
2333!
2334 RETURN

References mod_scalars::exit_flag, mod_parallel::inpthread, mod_iounits::ioerror, mod_parallel::master, mod_netcdf::netcdf_close(), mod_netcdf::netcdf_open(), mod_scalars::noerror, mod_iounits::sourcefile, and mod_iounits::stdout.

Here is the call graph for this function:

The documentation for this interface was generated from the following file: