Path: | gtvarattrsearch.f90 |
Last Update: | Thu Apr 27 17:08:18 JST 2006 |
Authors: | Eizi TOYODA, Yasuhiro MORIKAWA |
Version: | $Id: gtvarattrsearch.f90,v 1.4 2006/04/27 08:08:18 morikawa Exp $ |
Tag Name: | $Name: gt4f90io-20070616 $ |
Copyright: | Copyright (C) GFD Dennou Club, 2000-2005. All rights reserved. |
License: | See COPYRIGHT |
以下のサブルーチン, 関数は gtdata_generic から提供されます。
Subroutine : | |
var : | type(GT_VARIABLE), intent(inout), target |
name : | character(len = *), intent(out) |
end : | logical, intent(out), optional |
Attr_Rewind を参照してください。
subroutine GTVarAttrNext(var, name, end) ! !== 変数からの属性リスト取得 ! ! Attr_Rewind を参照してください。 ! use gtdata_types, only: GT_VARIABLE use gt_map, only: var_class, vtb_class_netcdf, vtb_class_memory use an_generic, only: attr_next, an_variable use gt_mem, only: attr_next, mem_variable implicit none type(GT_VARIABLE), intent(inout), target:: var character(len = *), intent(out):: name logical, intent(out), optional:: end integer:: class, cid continue call var_class(var, class, cid) select case(class) case(vtb_class_netcdf) call attr_next(an_variable(cid), name, end) case(vtb_class_memory) call attr_next(mem_variable(cid), name, end) end select end subroutine GTVarAttrNext
Subroutine : | |
var : | type(GT_VARIABLE), intent(inout), target |
var から属性名のリストを取得するために利用するサブルーチンです。 このサブルーチンと Attr_Next によって属性リスト一覧を取得できます。
ある変数 var について全ての属性を列挙するためには、まず Attr_Rewind を呼んだ後、Attr_Next を呼びます。最初の呼び出しで 最初の属性が、次の呼び出しで次の属性の名前が得られます。最後の 属性のあとでは end == .true. となります。
以下のサンプルソースコードを参照ください。
! 属性一覧の取得 use gt4f90io type(GT_VARIABLE):: var character(len = STRING):: attrname logical:: end call Attr_Rewind(var) do call Attr_Next(var, attrname, end) if (end) exit write(*,*) trim(attrname) enddo
subroutine GTVarAttrRewind(var) ! !== 変数からの属性リスト取得 (初期化用) ! ! *var* から属性名のリストを取得するために利用するサブルーチンです。 ! このサブルーチンと Attr_Next によって属性リスト一覧を取得できます。 ! ! ある変数 *var* について全ての属性を列挙するためには、まず ! Attr_Rewind を呼んだ後、Attr_Next を呼びます。最初の呼び出しで ! 最初の属性が、次の呼び出しで次の属性の名前が得られます。最後の ! 属性のあとでは end == .true. となります。 ! ! 以下のサンプルソースコードを参照ください。 ! ! ! ! 属性一覧の取得 ! use gt4f90io ! type(GT_VARIABLE):: var ! character(len = STRING):: attrname ! logical:: end ! ! call Attr_Rewind(var) ! do ! call Attr_Next(var, attrname, end) ! if (end) exit ! write(*,*) trim(attrname) ! enddo ! use gtdata_types, only: GT_VARIABLE use gt_map, only: var_class, vtb_class_netcdf, vtb_class_memory use an_generic, only: attr_rewind, an_variable use gt_mem, only: attr_rewind, mem_variable implicit none type(GT_VARIABLE), intent(inout), target:: var integer:: class, cid continue call var_class(var, class, cid) select case(class) case(vtb_class_netcdf) call attr_rewind(an_variable(cid)) case(vtb_class_memory) call attr_rewind(mem_variable(cid)) end select end subroutine GTVarAttrRewind