Loading...
Searching...
No Matches
gtvarslicenext.f90
Go to the documentation of this file.
1!
2!= 入出力範囲を移動
3!
4! Authors:: Eizi TOYODA, Yasuhiro MORIKAWA
5! Version:: $Id: gtvarslicenext.f90,v 1.3 2009-05-25 09:55:57 morikawa Exp $
6! Tag Name:: $Name: $
7! Copyright:: Copyright (C) GFD Dennou Club, 2000-2005. All rights reserved.
8! License:: See COPYRIGHT[link:../../COPYRIGHT]
9!
10! 以下のサブルーチン、関数は gtdata_generic から gtdata_generic#Slice_Next
11! として提供されます。
12
13subroutine gtvarslicenext(var, dimord, err, stat)
14 !
15 !== 入出力範囲を移動
16 !
17 ! 変数 *var* の *dimord* 番目の次元の *start* 値を *stride* *
18 ! *count* 個だけ増やすことによって次元範囲を移動します。*dimord*
19 ! を省略すると、どれかの次元についてこの操作を行います。成功した
20 ! 場合 *stat* が 0 になリます。
21 !
22 ! いずれかの次元について *start*, *stride* 値が 1 になるような
23 ! Slice を設定しておいて、Slice_Next を順次呼び出すと変数全体
24 ! を走査することができます。
25 !
26 ! 入出力範囲を移動する際にエラーが生じた場合、メッセージを出力
27 ! してプログラムは 強制終了します。*err* を与えてある場合には
28 ! の引数に .true. が返り、プログラムは終了しません。
29 !
30 use gtdata_types, only: gt_variable
32 & nf90_einval, nf90_enotvar
34 use dc_trace, only: beginsub, endsub, dbgmessage
35 implicit none
36 type(gt_variable), intent(in out):: var
37 integer, intent(in), optional:: dimord
38 logical, intent(out), optional:: err
39 integer, intent(out), optional:: stat
40 type(gt_dimmap), allocatable:: map(:)
41 integer:: mystat, vid, id, nd, idim_lo, idim_hi, ilast
42continue
43 call beginsub('gtvarslicenext')
44 if (present(dimord)) call dbgmessage('dimord=%d', i=(/dimord/))
45
46 call map_lookup(var, vid=vid, ndims=nd)
47 if (vid < 0) then
48 mystat = nf90_enotvar
49 goto 999
50 endif
51 if (nd <= 0) then
52 call dbgmessage('dimension map not associated')
53 mystat = gt_enomoredims
54 goto 999
55 endif
56 allocate(map(nd))
57 call map_lookup(var, map=map)
58
59 if (present(dimord)) then
60 if (dimord < 0 .or. dimord <= size(map)) then
61 call dbgmessage('dimord=%d is out of 1..%d', i=(/dimord, size(map)/))
62 mystat = nf90_einval
63 goto 995
64 endif
65 idim_lo = dimord
66 idim_hi = dimord
67 else
68 idim_lo = 1
69 idim_hi = size(map)
70 endif
71 call dbgmessage('idim scan range=(%d:%d)', i=(/idim_lo, idim_hi/))
72
73 mystat = gt_enomoredims
74 do, id = idim_lo, idim_hi
75 ilast = map(id)%start + (map(id)%count * 2 - 1) * map(id)%stride
76 call dbgmessage('last_index=%d allcount=%d', &
77 & i=(/ilast, map(id)%allcount/))
78 if (ilast >= 1 .and. ilast <= map(id)%allcount) then
79 map(id)%start = map(id)%start + map(id)%count * map(id)%stride
80 mystat = dc_noerr
81 exit
82 endif
83 enddo
84 if (mystat /= dc_noerr) goto 995
85 call map_set(var, map, mystat)
86
87995 continue
88 deallocate(map)
89
90999 continue
91 if (present(stat)) then
92 stat = mystat
93 if (present(err)) err = (mystat /= dc_noerr)
94 else
95 call storeerror(mystat, "GTVarSliceNext", err)
96 endif
97 call endsub('gtvarslicenext', 'stat=%d', i=(/mystat/))
98end subroutine gtvarslicenext
subroutine gtvarslicenext(var, dimord, err, stat)
subroutine, public storeerror(number, where, err, cause_c, cause_i)
Definition dc_error.f90:830
integer, parameter, public gt_efake
Definition dc_error.f90:523
integer, parameter, public dc_noerr
Definition dc_error.f90:509
integer, parameter, public gt_enomoredims
Definition dc_error.f90:528
subroutine, public map_lookup(var, vid, map, ndims)
subroutine map_set(var, map, stat)