Class | netcdf_wrapper |
In: |
io/netcdf_wrapper.F90
|
Note that Japanese and English are described in parallel.
gtool で扱わない部分をカバーする NetCDF のラッパープログラム.
The wrapper routines in this module treat NetCDF input/output that are not covered by gtool.
NWInqDimLen : | 軸の長さを問い合わせる |
———— : | ———— |
NWInqDimLen : | Inquire length of a dimension |
NAMELIST#radiation_DennouAGCM_nml
Subroutine : | |||
ncfn : | character(*), intent(in )
| ||
varname : | character(*), intent(in )
| ||
attname : | character(*), intent(in )
| ||
att : | character(*), intent(out)
|
Alias for NWGetAttChar
Subroutine : | |||
ncfn : | character(*), intent(in )
| ||
varname : | character(*), intent(in )
| ||
attname : | character(*), intent(in )
| ||
att : | integer , intent(out)
|
Alias for NWGetAttInteger
Subroutine : | |||
ncfn : | character(*), intent(in )
| ||
dimname : | character(*), intent(in )
| ||
len : | integer , intent(out)
|
軸の長さを問い合わせます.
Inquire length of a dimension.
subroutine NWInqDimLen( ncfn, dimname, len ) ! ! 軸の長さを問い合わせます. ! ! Inquire length of a dimension. ! character(*), intent(in ) :: ncfn ! NetCDF filename character(*), intent(in ) :: dimname ! Dimension name integer , intent(out) :: len ! Length of a dimension ! ! Local variables ! character(STRING) :: mode integer :: ncid integer :: dimid integer :: status character(STRING) :: err_mes err_mes = "In NWInqDimLen" mode = "read" call NWOpen( NWMkNCFN(ncfn), mode, ncid ) status = NF90_inq_dimid( ncid, dimname, dimid ) call NWHandleErr( status, err_mes ) status = NF90_Inquire_Dimension( ncid, dimid, len = len ) call NWHandleErr( status, err_mes ) call NWClose( ncid ) end subroutine NWInqDimLen
Function : | |||
Flag : | logical | ||
ncfn : | character(*), intent(in )
| ||
varname : | character(*), intent(in )
|
function NWPresentAVarInFile( ncfn, varname ) result( Flag ) character(*), intent(in ) :: ncfn ! NetCDF filename character(*), intent(in ) :: varname ! Variable name logical :: Flag ! Local variables ! character(STRING) :: mode integer :: ncid integer :: varid integer :: status character(STRING) :: err_mes err_mes = "In NWChkAVarInFile" err_mes = trim( err_mes ) // ', Var: ' // trim( varname ) mode = "read" call NWOpen( NWMkNCFN(ncfn), mode, ncid ) status = nf90_inq_varid( ncid, varname, varid ) if( status == nf90_noerr ) then Flag = .true. else Flag = .false. call MessageNotify( 'M', module_name, 'Variable, %c, cannot be found in file, %c.', c1 = trim(varname), c2 = trim(NWMkNCFN(ncfn)) ) end if call NWClose( ncid ) end function NWPresentAVarInFile
Subroutine : | |
ncid : | integer, intent(in) |
subroutine NWClose( ncid ) integer, intent(in) :: ncid ! Local variables ! integer :: status character(STRING) :: err_mes err_mes = "In NWClose" status = nf90_close( ncid ) call NWHandleErr( status, err_mes ) end subroutine NWClose
Subroutine : | |
ncid : | integer, intent(in) |
subroutine NWEndDef( ncid ) integer, intent(in) :: ncid ! Local variables ! integer :: status character(STRING) :: err_mes err_mes = "In NWEndDef" status = nf90_enddef( ncid ) if( ( status .ne. nf90_noerr ) .and. ( status .ne. nf90_enotindefine ) ) call NWHandleErr( status, err_mes ) end subroutine NWEndDef
Subroutine : | |||
ncfn : | character(*), intent(in )
| ||
varname : | character(*), intent(in )
| ||
attname : | character(*), intent(in )
| ||
att : | character(*), intent(out)
|
subroutine NWGetAttChar( ncfn, varname, attname, att ) character(*), intent(in ) :: ncfn ! NetCDF filename character(*), intent(in ) :: varname ! Variable name character(*), intent(in ) :: attname ! Attribute name character(*), intent(out) :: att ! Attribute ! Local variables ! character(STRING) :: mode integer :: ncid integer :: varid integer :: status character(STRING) :: err_mes err_mes = "In NWGetAttChar" err_mes = trim( err_mes ) // ', Var: ' // trim( varname ) // ', Att: ' // trim( attname ) mode = "read" call NWOpen( NWMkNCFN(ncfn), mode, ncid ) if ( ( varname == 'global' ) .or. ( varname == 'GLOBAL' ) ) then varid = NF90_GLOBAL else status = nf90_inq_varid( ncid, varname, varid ) call NWHandleErr( status, err_mes ) end if status = nf90_get_att( ncid, varid, attname, att ) call NWHandleErr( status, err_mes ) call NWClose( ncid ) end subroutine NWGetAttChar
Subroutine : | |||
ncfn : | character(*), intent(in )
| ||
varname : | character(*), intent(in )
| ||
attname : | character(*), intent(in )
| ||
att : | integer , intent(out)
|
subroutine NWGetAttInteger( ncfn, varname, attname, att ) character(*), intent(in ) :: ncfn ! NetCDF filename character(*), intent(in ) :: varname ! Variable name character(*), intent(in ) :: attname ! Attribute name integer , intent(out) :: att ! Attribute ! Local variables ! character(STRING) :: mode integer :: ncid integer :: varid integer :: status character(STRING) :: err_mes err_mes = "In NWGetAttInteger" err_mes = trim( err_mes ) // ', Var: ' // trim( varname ) // ', Att: ' // trim( attname ) mode = "read" call NWOpen( NWMkNCFN(ncfn), mode, ncid ) if ( ( varname == 'global' ) .or. ( varname == 'GLOBAL' ) ) then varid = NF90_GLOBAL else status = nf90_inq_varid( ncid, varname, varid ) call NWHandleErr( status, err_mes ) end if status = nf90_get_att( ncid, varid, attname, att ) call NWHandleErr( status, err_mes ) call NWClose( ncid ) end subroutine NWGetAttInteger
Subroutine : | |
status : | integer , intent(in) |
err_mes : | character(*), intent(in), optional |
subroutine NWHandleErr( status, err_mes ) integer , intent(in) :: status character(*), intent(in), optional :: err_mes if( status .ne. nf90_noerr ) then if( present( err_mes ) ) print *, trim( err_mes ) print *, trim( nf90_strerror( status ) ) stop "STOP" end if end subroutine NWHandleErr
Function : | |||
ncfn_out : | character(STRING) | ||
ncfn_in : | character(*), intent(in )
|
character(STRING) function NWMkNCFN( ncfn_in ) result( ncfn_out ) ! MPI wrapper ! use mpi_wrapper, only : myrank character(*), intent(in ) :: ncfn_in ! NetCDF filename #ifdef LIB_MPI if ( ncfn_in(len_trim(ncfn_in)-2:len_trim(ncfn_in)) /= '.nc' ) then write( 6, * ) trim(ncfn_in), ' is inappropriate.' stop end if write( ncfn_out, '(a,a,i6.6,a)' ) ncfn_in(1:len_trim(ncfn_in)-3), '_rank', myrank, '.nc' #else ncfn_out = ncfn_in #endif end function NWMkNCFN
Subroutine : | |
path : | character(*), intent(in ) |
mode : | character(*), intent(in ) |
ncid : | integer , intent(out) |
subroutine NWOpen( path, mode, ncid ) character(*), intent(in ) :: path character(*), intent(in ) :: mode integer , intent(out) :: ncid ! Local variables ! integer :: cmode, omode integer :: status character(STRING) :: err_mes err_mes = "In NWOpen" err_mes = trim( err_mes ) // ', path: ' // trim( path ) // ', mode: ' // trim( mode ) if( ( mode .eq. "new" ) .or. ( mode .eq. "NEW" ) ) then cmode = NF90_CLOBBER status = nf90_create( path, cmode, ncid ) call NWHandleErr( status, err_mes ) else if( ( mode .eq. "read" ) .or. ( mode .eq. "READ" ) ) then omode = NF90_NOWRITE else if( ( mode .eq. "write" ) .or. ( mode .eq. "WRITE" ) ) then omode = NF90_WRITE else write( 6, * ) "Inproper argument, mode" write( 6, * ) "Argument mode should be either of words below." write( 6, * ) " 'new', 'NEW', 'read', 'READ', 'write', or 'WRITE'" stop end if status = nf90_open( path, omode, ncid ) call NWHandleErr( status, err_mes ) end if end subroutine NWOpen
Subroutine : | |
ncid : | integer, intent(in) |
subroutine NWReDef( ncid ) integer, intent(in) :: ncid ! Local variables ! integer :: status character(STRING) :: err_mes err_mes = "In NWReDef" status = nf90_redef( ncid ) if( ( status .ne. nf90_noerr ) .and. ( status .ne. nf90_eindefine ) ) call NWHandleErr( status, err_mes ) end subroutine NWReDef
Constant : | |||
module_name = ‘netcdf_wrapper‘ : | character(*), parameter
|
Constant : | |||
version = ’$Name: dcpam5-20130921 $’ // ’$Id: netcdf_wrapper.F90,v 1.3 2011-09-10 18:23:29 yot Exp $’ : | character(*), parameter
|