historyautoaddweight.f90

Path: gtool/gtool_historyauto/historyautoaddweight.f90
Last Update: Sun May 10 21:19:18 +0900 2009

座標の重みデータの設定

Settings of weights of axes.

Authors:Yasuhiro MORIKAWA
Version:$Id: historyautoaddweight.f90,v 1.1 2009-05-10 12:19:18 morikawa Exp $
Tag Name:$Name: gtool5-20101228-1 $
Copyright:Copyright (C) GFD Dennou Club, 2008-2009. All rights reserved.
License:See COPYRIGHT

Required files

Methods

Included Modules

gtool_historyauto_internal gtool_history dc_trace dc_error dc_types

Public Instance methods

Subroutine :
dim :character(*), intent(in)
: 座標重みを設定する座標の名称.

ただし, ここで指定するもの は, HistoryAutoCreate の dims 既に指定されていなければなりません.

Name of axis to which "weight" are set.

Note that this value must be set as "dims" of "HistoryAutoCreate".

weight(:) :real(DP), intent(in)
: 座標重みデータ.

データ型は整数, 単精度実数型, 倍精度実数型のどれでもかまいません. ただし, ファイルへ出力される際には, xtype もしくは座標データの型へと 変換されます.

Weight of axis.

Integer, single or double precision are acceptable as data type. Note that when this is output to a file, data type is converted into "xtype" or type of the axis.

units :character(*), intent(in), optional
: 座標重みの単位. 省略した場合には, 座標の単位が 使用されます.

Units of axis weight. If this argument is omitted, unit of the dimension is used.

xtype :character(*), intent(in), optional
: 座標重みのデータ型. 省略した場合には, 座標のデータ型が 使用されます.

Data type of weight of the dimension. If this argument is omitted, data type of the dimension is used.

座標の重みデータを設定します.

Set weights of axes.

[Source]

  subroutine HistoryAutoAddWeightDouble( dim, weight, units, xtype )
    !
    ! 座標の重みデータを設定します. 
    !
    ! Set weights of axes. 
    !

    use gtool_historyauto_internal, only: initialized, numdims, numwgts, wgtsuf, gthst_axes, data_weights, gthst_weights
    use gtool_history, only: HistoryAxisInquire, HistoryAxisAddAttr, HistoryVarinfoCreate
    use dc_trace, only: BeginSub, EndSub
    use dc_error, only: StoreError, DC_NOERR, GT_EARGSIZEMISMATCH, HST_ENOAXISNAME, DC_ENOTINIT
    use dc_types, only: DP, STRING, TOKEN

    implicit none
    character(*), intent(in):: dim
                              ! 座標重みを設定する座標の名称. 
                              !
                              ! ただし, ここで指定するもの
                              ! は, HistoryAutoCreate の *dims*
                              ! 既に指定されていなければなりません.
                              !
                              ! Name of axis to which "weight" are set. 
                              !
                              ! Note that this value must be set 
                              ! as "dims" of "HistoryAutoCreate". 
                              !
    real(DP), intent(in):: weight(:)
                              ! 座標重みデータ. 
                              !
                              ! データ型は整数, 単精度実数型, 
                              ! 倍精度実数型のどれでもかまいません. 
                              ! ただし, ファイルへ出力される際には, 
                              ! xtype もしくは座標データの型へと
                              ! 変換されます. 
                              ! 
                              ! Weight of axis. 
                              ! 
                              ! Integer, single or double precision are 
                              ! acceptable as data type. 
                              ! Note that when this is output to a file, 
                              ! data type is converted into "xtype" or 
                              ! type of the axis. 
                              ! 
    character(*), intent(in), optional:: units
                              ! 座標重みの単位. 
                              ! 省略した場合には, 座標の単位が
                              ! 使用されます. 
                              !
                              ! Units of axis weight.
                              ! If this argument is omitted, 
                              ! unit of the dimension is used. 
                              !
    character(*), intent(in),  optional:: xtype
                              ! 座標重みのデータ型. 
                              ! 省略した場合には, 座標のデータ型が
                              ! 使用されます. 
                              ! 
                              ! Data type of weight of the dimension. 
                              ! If this argument is omitted, 
                              ! data type of the dimension is used. 
                              ! 

    character(STRING):: name, longname
    character(TOKEN):: dim_units, dim_xtype
    integer:: dim_size
    integer:: stat, i
    character(STRING):: cause_c
    character(*), parameter:: subname = "HistoryAutoAddWeightDouble"
  continue
    call BeginSub(subname, 'dim=<%c>', c1=trim(dim) )
    stat = DC_NOERR
    cause_c = ""

    ! 初期設定チェック
    ! Check initialization
    !
    if ( .not. initialized ) then
      stat = DC_ENOTINIT
      cause_c = 'gtool_historyauto'
      goto 999
    end if

    do i = 1, numdims
      call HistoryAxisInquire( axis = gthst_axes(i), name = name, size = dim_size, longname = longname, units = dim_units, xtype = dim_xtype )      ! (out)
      if ( trim(dim) == trim(name) ) then
        if ( dim_size /= size(weight) ) then
          stat = GT_EARGSIZEMISMATCH
          cause_c = 'weight'
        end if
        if ( present(units) ) dim_units = units
        if ( present(xtype) ) dim_xtype = xtype
        call HistoryVarinfoCreate( varinfo = gthst_weights(numwgts + 1), name = trim(dim) // wgtsuf, dims = (/ dim /), longname = 'weight for integration or average in ' // trim(longname), units = dim_units, xtype = dim_xtype )           ! (in)

        call HistoryAxisAddAttr( axis = gthst_axes(i), attrname = 'gt_calc_weight', value = trim(dim) // wgtsuf )  ! (in)

        allocate( data_weights(numwgts + 1) % a_axis( dim_size ) )
        data_weights(numwgts + 1) % a_axis = weight

        numwgts = numwgts + 1
        goto 999
      end if
    end do

    stat = HST_ENOAXISNAME
    cause_c = dim
    
999 continue
    call StoreError(stat, subname, cause_c = cause_c)
    call EndSub(subname)
  end subroutine HistoryAutoAddWeightDouble
Subroutine :
dim :character(*), intent(in)
: 座標重みを設定する座標の名称.

ただし, ここで指定するもの は, HistoryAutoCreate の dims 既に指定されていなければなりません.

Name of axis to which "weight" are set.

Note that this value must be set as "dims" of "HistoryAutoCreate".

weight(:) :integer, intent(in)
: 座標重みデータ.

データ型は整数, 単精度実数型, 倍精度実数型のどれでもかまいません. ただし, ファイルへ出力される際には, xtype もしくは座標データの型へと 変換されます.

Weight of axis.

Integer, single or double precision are acceptable as data type. Note that when this is output to a file, data type is converted into "xtype" or type of the axis.

units :character(*), intent(in), optional
: 座標重みの単位. 省略した場合には, 座標の単位が 使用されます.

Units of axis weight. If this argument is omitted, unit of the dimension is used.

xtype :character(*), intent(in), optional
: 座標重みのデータ型. 省略した場合には, 座標のデータ型が 使用されます.

Data type of weight of the dimension. If this argument is omitted, data type of the dimension is used.

座標の重みデータを設定します.

Set weights of axes.

[Source]

  subroutine HistoryAutoAddWeightInt( dim, weight, units, xtype )
    !
    ! 座標の重みデータを設定します. 
    !
    ! Set weights of axes. 
    !

    use gtool_historyauto_internal, only: initialized, numdims, numwgts, wgtsuf, gthst_axes, data_weights, gthst_weights
    use gtool_history, only: HistoryAxisInquire, HistoryAxisAddAttr, HistoryVarinfoCreate
    use dc_trace, only: BeginSub, EndSub
    use dc_error, only: StoreError, DC_NOERR, GT_EARGSIZEMISMATCH, HST_ENOAXISNAME, DC_ENOTINIT
    use dc_types, only: DP, STRING, TOKEN

    implicit none
    character(*), intent(in):: dim
                              ! 座標重みを設定する座標の名称. 
                              !
                              ! ただし, ここで指定するもの
                              ! は, HistoryAutoCreate の *dims*
                              ! 既に指定されていなければなりません.
                              !
                              ! Name of axis to which "weight" are set. 
                              !
                              ! Note that this value must be set 
                              ! as "dims" of "HistoryAutoCreate". 
                              !
    integer, intent(in):: weight(:)
                              ! 座標重みデータ. 
                              !
                              ! データ型は整数, 単精度実数型, 
                              ! 倍精度実数型のどれでもかまいません. 
                              ! ただし, ファイルへ出力される際には, 
                              ! xtype もしくは座標データの型へと
                              ! 変換されます. 
                              ! 
                              ! Weight of axis. 
                              ! 
                              ! Integer, single or double precision are 
                              ! acceptable as data type. 
                              ! Note that when this is output to a file, 
                              ! data type is converted into "xtype" or 
                              ! type of the axis. 
                              ! 
    character(*), intent(in), optional:: units
                              ! 座標重みの単位. 
                              ! 省略した場合には, 座標の単位が
                              ! 使用されます. 
                              !
                              ! Units of axis weight.
                              ! If this argument is omitted, 
                              ! unit of the dimension is used. 
                              !
    character(*), intent(in),  optional:: xtype
                              ! 座標重みのデータ型. 
                              ! 省略した場合には, 座標のデータ型が
                              ! 使用されます. 
                              ! 
                              ! Data type of weight of the dimension. 
                              ! If this argument is omitted, 
                              ! data type of the dimension is used. 
                              ! 

    character(STRING):: name, longname
    character(TOKEN):: dim_units, dim_xtype
    integer:: dim_size
    integer:: stat, i
    character(STRING):: cause_c
    character(*), parameter:: subname = "HistoryAutoAddWeightInt"
  continue
    call BeginSub(subname, 'dim=<%c>', c1=trim(dim) )
    stat = DC_NOERR
    cause_c = ""

    ! 初期設定チェック
    ! Check initialization
    !
    if ( .not. initialized ) then
      stat = DC_ENOTINIT
      cause_c = 'gtool_historyauto'
      goto 999
    end if

    do i = 1, numdims
      call HistoryAxisInquire( axis = gthst_axes(i), name = name, size = dim_size, longname = longname, units = dim_units, xtype = dim_xtype )      ! (out)
      if ( trim(dim) == trim(name) ) then
        if ( dim_size /= size(weight) ) then
          stat = GT_EARGSIZEMISMATCH
          cause_c = 'weight'
        end if
        if ( present(units) ) dim_units = units
        if ( present(xtype) ) dim_xtype = xtype
        call HistoryVarinfoCreate( varinfo = gthst_weights(numwgts + 1), name = trim(dim) // wgtsuf, dims = (/ dim /), longname = 'weight for integration or average in ' // trim(longname), units = dim_units, xtype = dim_xtype )           ! (in)

        call HistoryAxisAddAttr( axis = gthst_axes(i), attrname = 'gt_calc_weight', value = trim(dim) // wgtsuf )  ! (in)

        allocate( data_weights(numwgts + 1) % a_axis( dim_size ) )
        data_weights(numwgts + 1) % a_axis = weight

        numwgts = numwgts + 1
        goto 999
      end if
    end do

    stat = HST_ENOAXISNAME
    cause_c = dim
    
999 continue
    call StoreError(stat, subname, cause_c = cause_c)
    call EndSub(subname)
  end subroutine HistoryAutoAddWeightInt
Subroutine :
dim :character(*), intent(in)
: 座標重みを設定する座標の名称.

ただし, ここで指定するもの は, HistoryAutoCreate の dims 既に指定されていなければなりません.

Name of axis to which "weight" are set.

Note that this value must be set as "dims" of "HistoryAutoCreate".

weight(:) :real, intent(in)
: 座標重みデータ.

データ型は整数, 単精度実数型, 倍精度実数型のどれでもかまいません. ただし, ファイルへ出力される際には, xtype もしくは座標データの型へと 変換されます.

Weight of axis.

Integer, single or double precision are acceptable as data type. Note that when this is output to a file, data type is converted into "xtype" or type of the axis.

units :character(*), intent(in), optional
: 座標重みの単位. 省略した場合には, 座標の単位が 使用されます.

Units of axis weight. If this argument is omitted, unit of the dimension is used.

xtype :character(*), intent(in), optional
: 座標重みのデータ型. 省略した場合には, 座標のデータ型が 使用されます.

Data type of weight of the dimension. If this argument is omitted, data type of the dimension is used.

座標の重みデータを設定します.

Set weights of axes.

[Source]

  subroutine HistoryAutoAddWeightReal( dim, weight, units, xtype )
    !
    ! 座標の重みデータを設定します. 
    !
    ! Set weights of axes. 
    !

    use gtool_historyauto_internal, only: initialized, numdims, numwgts, wgtsuf, gthst_axes, data_weights, gthst_weights
    use gtool_history, only: HistoryAxisInquire, HistoryAxisAddAttr, HistoryVarinfoCreate
    use dc_trace, only: BeginSub, EndSub
    use dc_error, only: StoreError, DC_NOERR, GT_EARGSIZEMISMATCH, HST_ENOAXISNAME, DC_ENOTINIT
    use dc_types, only: DP, STRING, TOKEN

    implicit none
    character(*), intent(in):: dim
                              ! 座標重みを設定する座標の名称. 
                              !
                              ! ただし, ここで指定するもの
                              ! は, HistoryAutoCreate の *dims*
                              ! 既に指定されていなければなりません.
                              !
                              ! Name of axis to which "weight" are set. 
                              !
                              ! Note that this value must be set 
                              ! as "dims" of "HistoryAutoCreate". 
                              !
    real, intent(in):: weight(:)
                              ! 座標重みデータ. 
                              !
                              ! データ型は整数, 単精度実数型, 
                              ! 倍精度実数型のどれでもかまいません. 
                              ! ただし, ファイルへ出力される際には, 
                              ! xtype もしくは座標データの型へと
                              ! 変換されます. 
                              ! 
                              ! Weight of axis. 
                              ! 
                              ! Integer, single or double precision are 
                              ! acceptable as data type. 
                              ! Note that when this is output to a file, 
                              ! data type is converted into "xtype" or 
                              ! type of the axis. 
                              ! 
    character(*), intent(in), optional:: units
                              ! 座標重みの単位. 
                              ! 省略した場合には, 座標の単位が
                              ! 使用されます. 
                              !
                              ! Units of axis weight.
                              ! If this argument is omitted, 
                              ! unit of the dimension is used. 
                              !
    character(*), intent(in),  optional:: xtype
                              ! 座標重みのデータ型. 
                              ! 省略した場合には, 座標のデータ型が
                              ! 使用されます. 
                              ! 
                              ! Data type of weight of the dimension. 
                              ! If this argument is omitted, 
                              ! data type of the dimension is used. 
                              ! 

    character(STRING):: name, longname
    character(TOKEN):: dim_units, dim_xtype
    integer:: dim_size
    integer:: stat, i
    character(STRING):: cause_c
    character(*), parameter:: subname = "HistoryAutoAddWeightReal"
  continue
    call BeginSub(subname, 'dim=<%c>', c1=trim(dim) )
    stat = DC_NOERR
    cause_c = ""

    ! 初期設定チェック
    ! Check initialization
    !
    if ( .not. initialized ) then
      stat = DC_ENOTINIT
      cause_c = 'gtool_historyauto'
      goto 999
    end if

    do i = 1, numdims
      call HistoryAxisInquire( axis = gthst_axes(i), name = name, size = dim_size, longname = longname, units = dim_units, xtype = dim_xtype )      ! (out)
      if ( trim(dim) == trim(name) ) then
        if ( dim_size /= size(weight) ) then
          stat = GT_EARGSIZEMISMATCH
          cause_c = 'weight'
        end if
        if ( present(units) ) dim_units = units
        if ( present(xtype) ) dim_xtype = xtype
        call HistoryVarinfoCreate( varinfo = gthst_weights(numwgts + 1), name = trim(dim) // wgtsuf, dims = (/ dim /), longname = 'weight for integration or average in ' // trim(longname), units = dim_units, xtype = dim_xtype )           ! (in)

        call HistoryAxisAddAttr( axis = gthst_axes(i), attrname = 'gt_calc_weight', value = trim(dim) // wgtsuf )  ! (in)

        allocate( data_weights(numwgts + 1) % a_axis( dim_size ) )
        data_weights(numwgts + 1) % a_axis = weight

        numwgts = numwgts + 1
        goto 999
      end if
    end do

    stat = HST_ENOAXISNAME
    cause_c = dim
    
999 continue
    call StoreError(stat, subname, cause_c = cause_c)
    call EndSub(subname)
  end subroutine HistoryAutoAddWeightReal