gtvarcreatecopy.f90

Path: src/gtvarcreatecopy.f90
Last Update: Wed Jul 20 18:22:24 JST 2005

Copyright (C) GFD Dennou Club, 2000. All rights reserved.

Methods

Included Modules

gtdata_types dc_types gtdata_generic dc_url dc_trace dc_error dc_string

Public Instance methods

to :type(GT_VARIABLE), intent(inout)
from :type(GT_VARIABLE), intent(inout)

すでに存在する変数について、値をコピーする。

[Source]

    subroutine GTVarCopyValue(to, from)

        type(GT_VARIABLE), intent(inout):: to
        type(GT_VARIABLE), intent(inout):: from
        real, allocatable:: rbuffer(:)
        logical:: err
        integer:: siz, stat
        !
        call beginsub('gtvarcopyvalue')
        ! 値のコピー
        call Slice(from)
        call Slice(to, compatible=from)
        call Inquire(from, size=siz)
        allocate (rbuffer(siz))
        do
            call GTVarGetReal(from, rbuffer, siz, err)
            if (err) call DumpError()
            call GtVarPutReal(to, rbuffer, siz, err)
            if (err) call DumpError()
            call Slice_Next(from, stat=stat)
            if (stat /= 0) exit
            call Slice_Next(to, stat=stat)
        enddo
        deallocate (rbuffer)
        call endsub('gtvarcopyvalue')
    end subroutine
var :type(GT_VARIABLE), intent(out)
url :character(len = *), intent(in)
url :character(len = string)
copyfrom :type(GT_VARIABLE), intent(inout)
copyvalue :logical, intent(in), optional
overwrite :logical, intent(in), optional
err :logical, intent(out), optional

なお、次元変数の複製は copyfrom と url が異なるファイルに 載っている場合に行なわれる。これは netCDF/an を想定したものだが ほかのファイル形式が追加されたときには変更を要するかもしれない。

[Source]


subroutine GTVarCreateCopyC(var, url, copyfrom, copyvalue,  overwrite, err)

    implicit none
    intrinsic trim
    type(GT_VARIABLE),    intent(out)   :: var
    character(len = *),   intent(in)    :: url
    type(GT_VARIABLE),    intent(inout) :: copyfrom
    logical, intent(in),  optional      :: copyvalue
    logical, intent(in),  optional      :: overwrite
    logical, intent(out), optional      :: err
    type(GT_VARIABLE),    allocatable   :: vDimSource(:)
    type(GT_VARIABLE),    allocatable   :: vDimDest(:)
    integer                             :: i, nd, stat
    logical                             :: myerr
    character(STRING)                   :: vpart, upart, desturl
    character(TOKEN)                    :: xtype
continue
    call beginsub('gtvarcreatecopy', 'url=%c copyfrom=%d',  c1=trim(url), i=(/copyfrom%mapid/))
    stat = 0
    myerr = .FALSE.
    !-----------------------------------------------------------------
    !  コピーする変数の次元をコピー先のファイルに作成
    !-----------------------------------------------------------------
    !----- コピー元 copyfrom の次元変数の取得 -----
    call Inquire(copyfrom, alldims=nd)
    allocate(vDimSource(nd), vDimDest(nd), stat=stat)
    if (stat /= 0) goto 999
    desturl = url
    !----- コピー元 copyfrom の各次元情報を vDimSource に取り出し, -----
    !----- それをコピー先 desturl へコピーしてその次元 ID を       -----
    !----- vDimDest に返してもらう.                                -----
    do, i = 1, nd
        call Open(vDimSource(i), copyfrom, dimord=i,  count_compact=.TRUE., err=myerr)
        call GTVarCopyDim(to=vDimDest(i), from=vDimSource(i),  target=desturl)
    end do
    !-----------------------------------------------------------------
    !  変数作成
    !-----------------------------------------------------------------
    !----- url に変数名が無い場合, コピー元の変数名を使用 -----
    call UrlSplit(url, var=vpart)
    if (vpart == "") then
        call Inquire(copyfrom, url=upart)
        call UrlSplit(upart, var=vpart)
        desturl = trim(desturl) // GT_ATMARK // trim(vpart)
    end if
    !----- 実際に変数作成 -----
    call Inquire(copyfrom, xtype=xtype)
    call Create(var, trim(desturl), dims=vDimDest, xtype=xtype,       overwrite=overwrite, err=myerr)
    if (myerr) goto 990
    call copy_attr(to=var, from=copyfrom, err=myerr)
    if (myerr) goto 990
    if (present(copyvalue)) then
        if (copyvalue) then
            call GTVarCopyValue(to=var, from=copyfrom)
        endif
    endif
    do, i = 1, nd
        call Close(vDimSource(i))
        call Close(vDimDest(i))
    end do
990 continue
    deallocate(vDimSource, vDimDest, stat=stat)
999 continue
    if (stat /= 0) then
        call StoreError(GT_ENOMEM, "GTVarCreateCopy", err)
    else if (present(err)) then
        err = myerr
    else if (myerr) then
        call DumpError
    end if
    call endsub('gtvarcreatecopy', 'result=%d', i=(/var%mapid/))
contains

    ! from と同じ内容の次元変数を URL target で示される変数の作成時に
    ! 次元として使えるように to に複写。
    ! なるべく再オープンで済まそうとする。
    ! 複写する場合もなるべく次元名を合わせようとする。
    !
    subroutine GTVarCopyDim(to, from, target)

        type(GT_VARIABLE), intent(out):: to
        type(GT_VARIABLE), intent(inout):: from
        character(len = *), intent(in):: target
        character(len = string):: url, file, dimname
        character(len = token):: xtype
        logical:: growable, myerr
        integer:: length
        call beginsub('gtvarcopydim', 'from=%d target=<%c>',  i=(/from%mapid/), c1=trim(target))
        !----- 同じファイル上にコピーする場合は参照カウンタを1つ回すだけ -----
        call Inquire(var=from, url=url)
        if (trim(url) .onthesamefile. trim(target)) then
            call Open(to, from, dimord=0)
            call endsub('gtvarcopydim', 'dup-handle')
            return
        endif
        !----- 異なるファイル上にコピーする場合, 既に次元変数 from が -----
        !----- target の次元変数として含まれるかチェック              -----
        call UrlSplit(target, file=file)
        if (LookupEquivalent(to, from, file)) then
           !----- 含まれる場合はそれで終了 -----
           call endsub('gtvarcopydim', 'equivalent-exists')
           return
        else
           !----- 含まれない場合次元変数 from を target 上に作成 -----
           ! 次元変数 from が無制限次元である場合には長さを 0 に
           call Inquire(var=from, growable=growable, allcount=length)
           if (growable) length = 0
           call Inquire(var=from, xtype=xtype, name=dimname)
           !
           url = urlmerge(file, dimname)
           call Create(to, trim(url), length, xtype, err=myerr)
           if (myerr) then
              ! 指定名称でうまくいかない場合は自動生成名にする
              call Create(to, trim(file), length, xtype)
           endif
           call copy_attr(to, from, myerr)
           call GTVarCopyValue(to, from)
           call endsub('gtvarcopydim', 'created')
           return
        endif
    end subroutine
result :logical
: —— 次元変数 from のサイズと単位, 無制限次元かどうかを探査 —— —— コピー先 file の変数情報を探査 —— とりあえずは次元だけでなく全ての変数について開く 次元変数のサイズと, 無制限次元かどうかを取得
  (次元変数でないもののサイズは, 依存する次元変数のサイズを
   掛け合わせたものとなるので, もしかすると誤動作するかも).

次元変数 from が無制限次元で, 且つ file 内の次元変数も 無制限次元の場合は, 同じ次元変数と考える. 次元変数 from のサイズと file 内の次元変数のサイズが 異なる場合はスキップ 本当は dc_units で比較すべきだがとりあえず文字列比較

to :type(GT_VARIABLE), intent(out)
from :type(GT_VARIABLE), intent(in)
file :character(len = *), intent(in)
     すること.
     ※ もしかすると条件が足りないかも知れない.

[Source]

    logical function LookupEquivalent(to, from, file) result(result)

        type(GT_VARIABLE), intent(out):: to
        type(GT_VARIABLE), intent(in):: from
        character(len = *), intent(in):: file
        character(len = string):: url, units1, units2, reason
        logical:: end, growable1, growable2
        integer:: len1, len2
        character(len = *), parameter:: subnam = "lookupequivalent"
        call beginsub(subnam, 'from=%d file=<%c>',  i=(/from%mapid/), c1=trim(file))
        result = .FALSE.
        !----- 次元変数 from のサイズと単位, 無制限次元かどうかを探査 -----
        call Inquire(from, allcount=len1, growable=growable1)
        call get_attr(from, 'units', units1, default='')
        !----- コピー先 file の変数情報を探査 -----
        ! とりあえずは次元だけでなく全ての変数について開く
        call GTVarSearch(file)
        do
            call GTVarSearch(url, end)
            if (end) exit
            call Open(to, url, writable=.TRUE., err=end)
            if (end) exit
            ! 次元変数のサイズと, 無制限次元かどうかを取得
            !   (次元変数でないもののサイズは, 依存する次元変数のサイズを
            !    掛け合わせたものとなるので, もしかすると誤動作するかも).
            call Inquire(to, allcount=len2, growable=growable2)
            ! 次元変数 from が無制限次元で, 且つ file 内の次元変数も
            ! 無制限次元の場合は, 同じ次元変数と考える.
            if (.not. growable1 .or. .not. growable2) then
               ! 次元変数 from のサイズと file 内の次元変数のサイズが
               ! 異なる場合はスキップ
               if (len1 /= len2) then
                  call Close(to)
                  cycle
               endif
               call get_attr(to, 'units', units2, default='')
               ! 本当は dc_units で比較すべきだがとりあえず文字列比較
               if (units1 /= units2) then
                  call Close(to)
                  cycle
               else
                  reason = 'length of from is ' // trim(toChar(len1)) //    '. units of from is ' // "[" //                  trim(units1) // "]" //                           '. And file has same length and units.'
               endif
            else
               reason = 'from is UNLIMITED dimension, and file has it'
            endif
            result = .TRUE.
            call endsub(subnam, 'found (%c)', c1=trim(reason))
            return
        enddo
        call endsub(subnam, 'not found')
    end function

[Validate]