gtvarsearch.f90

Path: gtvarsearch.f90
Last Update: Sun Jan 15 16:10:30 JST 2006

変数探査

Authors:Yasuhiro MORIKAWA, Eizi TOYODA
Version:$Id: gtvarsearch.f90,v 1.3 2006/01/15 07:10:30 morikawa Exp $
Tag Name:$Name: gt4f90io-20070719 $
Copyright:Copyright (C) GFD Dennou Club, 2000-2005. All rights reserved.
License:See COPYRIGHT

以下のサブルーチン、関数は gtdata_generic から gtdata_generic#GTVarSearch として提供されます。

This file provides following module.

Required files

Methods

Included Modules

an_generic gtvarsearch dc_trace

Public Instance methods

Subroutine :
urlBase :character(len = *), intent(in)

変数探査初期化サブルーチン

上記の GTVarSearch を参照してください。

[Source]

subroutine GTVarSearchInit(urlBase)
  !
  !== 変数探査初期化サブルーチン
  !
  ! 上記の GTVarSearch を参照してください。
  !
  use an_generic, only: var_search
  use gtvarsearch
  use dc_trace, only: beginsub, endsub
  character(len = *), intent(in):: urlBase
  call beginsub('gtvarsearchinit', 'urlbase=<%c>', c1=trim(urlbase))
  call var_search(an, urlBase=urlBase)
  call endsub('gtvarsearchinit')
end subroutine
Subroutine :
url :character(len = *), intent(out)
end :logical, intent(out)

変数探査サブルーチン

あるファイル名 urlBase に依存する変数すべてを取得するには、 まず GTVarSearch(urlBase) (下記のサブルーチン) を呼び出し、 その後無限ループの中で GTVarSearch(url, end) を呼び出します。 そうすることで url に1つ1つの変数名が返ります。 end が真になったとき、すべての変数名を探索し終えたことになります。

   use gt4f90io
   character(len = STRING) :: filename, varname
   logical                 :: end

   write(*,*) "Enter file name: "
   read(*,*) filename

   call GTVarSearch(filename)
   do
     call GTVarSearch(varname, end)
     if (end) exit
     write(*, *) trim(varname)
   enddo

[Source]

subroutine GTVarSearchNext(url, end)
  !
  !== 変数探査サブルーチン
  !
  ! あるファイル名 urlBase に依存する変数すべてを取得するには、
  ! まず GTVarSearch(urlBase) (下記のサブルーチン) を呼び出し、
  ! その後無限ループの中で GTVarSearch(url, end) を呼び出します。
  ! そうすることで url に1つ1つの変数名が返ります。
  ! *end* が真になったとき、すべての変数名を探索し終えたことになります。
  !
  !=== 例
  !
  !    use gt4f90io
  !    character(len = STRING) :: filename, varname
  !    logical                 :: end
  !
  !    write(*,*) "Enter file name: "
  !    read(*,*) filename
  !
  !    call GTVarSearch(filename)
  !    do
  !      call GTVarSearch(varname, end)
  !      if (end) exit
  !      write(*, *) trim(varname)
  !    enddo
  !
  use an_generic, only: var_search
  use gtvarsearch
  use dc_trace, only: beginsub, endsub
  character(len = *), intent(out):: url
  logical, intent(out):: end
  call beginsub('gtvarsearchnext')
  call var_search(an, url=url, end=end)
  call endsub('gtvarsearchnext', 'url=%c end=%y', c1=trim(url), L=(/end/))
end subroutine

[Validate]