gtvarexchdim.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine gtvarexchdim (var, dimord1, dimord2, count_compact, err)
 

Function/Subroutine Documentation

◆ gtvarexchdim()

subroutine gtvarexchdim ( type(gt_variable), intent(in)  var,
integer, intent(in)  dimord1,
integer, intent(in)  dimord2,
logical, intent(in), optional  count_compact,
logical, intent(out)  err 
)

Definition at line 14 of file gtvarexchdim.f90.

References dc_trace::beginsub(), dc_trace::dbgmessage(), gtdata_internal_map::dimord_skip_compact(), dc_trace::endsub(), gtdata_internal_map::map_lookup(), gtdata_internal_map::map_set(), and gtdata_internal_map::map_set_ndims().

14  !
15  !== 次元順序番号の交換
16  !
17  ! 変数 *var* の次元順序番号 <b>dimord1</b>, <b>dimord2</b> のそれぞれに
18  ! 対応する次元を入れ替えます。
19  !
20  ! *count_compact* に .true. を渡すと、縮退した次元も含めて
21  ! 動作します。
22  !
23  ! エラーが生じた場合、メッセージを出力
24  ! してプログラムは強制終了します。*err* を与えてある場合には
25  ! の引数に .true. が返り、プログラムは終了しません。
26  !
27  use gtdata_types, only: gt_variable
30  use dc_trace, only: beginsub, endsub, dbgmessage
31  implicit none
32  type(gt_variable), intent(in):: var
33  integer, intent(in):: dimord1, dimord2
34  logical, intent(in), optional:: count_compact
35  logical, intent(out):: err
36  type(gt_dimmap), allocatable:: map(:)
37  type(gt_dimmap):: tmpmap
38  integer:: ndimsp, stat, idim1, idim2
39  logical:: direct_mode
40  character(*), parameter:: subname = 'GTVarExchDim'
41 continue
42  err = .true.
43  direct_mode = .false.
44  if (present(count_compact)) then
45  direct_mode = count_compact
46  endif
47  call beginsub(subname)
48  if (dimord1 < 1 .or. dimord2 < 1) then
49  call endsub(subname, "negative dimord=%d %d invalid", i=(/dimord1, dimord2/))
50  return
51  endif
52  call map_lookup(var, ndims=ndimsp)
53  if (ndimsp <= 0) then
54  call endsub(subname, "variable invalid")
55  return
56  else if (dimord1 > ndimsp .or. dimord2 > ndimsp) then
57  call endsub(subname, "dimord=%d %d not exist", i=(/dimord1, dimord2/))
58  return
59  endif
60 
61  allocate(map(ndimsp))
62  call map_lookup(var, map=map)
63 
64  if (.not. direct_mode) then
65  idim1 = dimord_skip_compact(dimord1, map)
66  idim2 = dimord_skip_compact(dimord2, map)
67  if (idim1 < 0 .or. idim2 < 0) then
68  call endsub(subname, "dimord=%d %d not found after compaction", &
69  & i=(/dimord1, dimord2/))
70  deallocate(map)
71  return
72  endif
73  else
74  idim1 = dimord1
75  idim2 = dimord2
76  endif
77 
78  tmpmap = map(idim1)
79  map(idim1) = map(idim2)
80  map(idim2) = tmpmap
81  call map_set(var, map, stat)
82  deallocate(map)
83 
84  err = stat /= 0
85  call endsub(subname)
integer function dimord_skip_compact(dimord, map)
subroutine map_set(var, map, stat)
subroutine map_set_ndims(var, ndims, stat)
subroutine, public dbgmessage(fmt, i, r, d, L, n, c1, c2, c3, ca)
Definition: dc_trace.f90:509
subroutine, public beginsub(name, fmt, i, r, d, L, n, c1, c2, c3, ca, version)
Definition: dc_trace.f90:351
subroutine, public map_lookup(var, vid, map, ndims)
subroutine, public endsub(name, fmt, i, r, d, L, n, c1, c2, c3, ca)
Definition: dc_trace.f90:446
Here is the call graph for this function: