historycopy.f90 File Reference

Go to the source code of this file.

Functions/Subroutines

subroutine historycopy1 (hist_dest, file, hist_src, title, source, institution, origin, interval, conventions, gt_version)
 
subroutine historycopy2 (hist_dest, file, hist_src, title, source, institution, origin, interval, conventions, gt_version)
 

Function/Subroutine Documentation

◆ historycopy1()

subroutine historycopy1 ( type(gt_history), intent(out), target  hist_dest,
character(*), intent(in)  file,
type(gt_history), intent(in), optional, target  hist_src,
character(*), intent(in), optional  title,
character(*), intent(in), optional  source,
character(*), intent(in), optional  institution,
real, intent(in), optional  origin,
real, intent(in), optional  interval,
character(*), intent(in), optional  conventions,
character(*), intent(in), optional  gt_version 
)

Definition at line 14 of file historycopy.f90.

References dc_trace::beginsub(), dc_trace::dbgmessage(), dc_error::dc_noerr, gtool_history_internal::default, dc_types::dp, dc_trace::endsub(), dc_error::gt_eargsizemismatch, dc_error::storeerror(), dc_types::string, and dc_types::token.

14  !
15  ! 引数 *hist_src* の内容にコピーし, *hist_dest* へ返します. *hist_src*
16  ! が与えられない場合は, 引数 *history* を与えずに呼び出した
17  ! HistoryCreate の設定内容が参照されます.
18  ! HistoryCreate と同様に, 出力の初期設定を行います. *file*
19  ! は必ず与えなければならず, *hist_src* と同じファイルへ出力
20  ! しようとする場合はエラーを生じます.
21  ! HistoryAddVariable で設定される内容に関してはコピーされません.
22  !
23  ! それ以降の引数を与えることで, hist_src の設定を
24  ! 上書きすることが可能です.
25  !
30  use dc_present, only: present_select
31  use dc_types, only: string, token, dp
32  use dc_date, only: evalbyunit
33  use dc_trace, only: beginsub, endsub, dbgmessage
35  implicit none
36  type(gt_history), intent(out), target:: hist_dest
37  character(*), intent(in):: file
38  type(gt_history), intent(in), optional, target:: hist_src
39  character(*), intent(in), optional:: title, source, institution
40  real, intent(in), optional:: origin, interval
41  character(*), intent(in), optional:: conventions, gt_version
42 
43  ! Internal Work
44  type(gt_history), pointer:: src =>null()
45  character(STRING) :: title_src, source_src, institution_src
46  character(STRING) :: conventions_src, gt_version_src
47  character(STRING), pointer:: dims(:) => null()
48  integer , pointer:: dimsizes(:) => null()
49  character(STRING), pointer:: longnames(:) => null()
50  character(STRING), pointer:: units(:) => null()
51  character(STRING), pointer:: xtypes(:) => null()
52  real:: originw, intervalw
53  integer :: i, numdims
54  logical :: err
55  real(DP),pointer :: dimvalue(:) => null()
56  character(len = *),parameter:: subname = "HistoryCopy1"
57  continue
58  call beginsub(subname, 'file=<%c>', c1=trim(file))
59 
60  if (present(hist_src)) then
61  src => hist_src
62  else
63  src => default
64  endif
65 
66 
67 
68 
69 
70 
71  numdims = size(src % dimvars)
72 
73  call historyinquire(history=src, title=title_src, &
74  & source=source_src, institution=institution_src, &
75  & dims=dims, dimsizes=dimsizes, longnames=longnames, &
76  & units=units, xtypes=xtypes, &
77  & conventions=conventions_src, gt_version=gt_version_src)
78 
79  if ( present(origin) ) then
80  originw = origin
81  else
82  originw = src % origin
83 ! originw = EvalByUnit( src % origin, '', src % unlimited_units_symbol )
84  end if
85 
86  intervalw = src % interval
87 ! intervalw = EvalByUnit( src % interval, '', src % unlimited_units_symbol )
88  if ( present(interval) ) then
89  if ( interval /= 0.0 ) then
90  intervalw = interval
91  end if
92  end if
93 
94  call historycreate(file=trim(file), &
95  & title=trim(present_select('', title_src, title)), &
96  & source=trim(present_select('', source_src, source)), &
97  & institution=trim(present_select('', institution_src, institution)), &
98  & dims=dims, dimsizes=dimsizes, longnames=longnames, units=units, &
99  & origin=originw, interval=intervalw, &
100  & xtypes=xtypes, &
101  & history=hist_dest, &
102  & conventions=trim(present_select('', conventions_src, conventions)), &
103  & gt_version=trim(present_select('', gt_version_src, gt_version)) )
104 
105  !
106  ! 次元変数が属性を持っている場合のことも考え, 最後に直接
107  ! hist_dst % dimvars へ copy_attr (gtvarcopyattrall) をかける.
108  !
109  do i = 1, numdims
110  call copy_attr(hist_dest % dimvars(i), src % dimvars (i), global=.false.)
111  end do
112 
113  ! dimvars を Get してみて, 値を持っているようならデータを与えてしまう.
114  do i = 1, numdims
115  if (dimsizes(i) == 0) cycle
116  call get(src % dimvars(i), dimvalue, err)
117  if (err) cycle
118  call historyput(dims(i), dimvalue, hist_dest)
119  deallocate(dimvalue)
120  end do
121 
122  deallocate(dims, dimsizes, longnames, units, xtypes)
123 
124 
125 
126 
127 
128  call endsub(subname)
type(gt_history), target, save, public default
integer, parameter, public token
単語やキーワードを保持する文字型変数の種別型パラメタ
Definition: dc_types.f90:109
subroutine, public storeerror(number, where, err, cause_c, cause_i)
Definition: dc_error.f90:830
integer, parameter, public dc_noerr
Definition: dc_error.f90:509
integer, parameter, public dp
倍精度実数型変数
Definition: dc_types.f90:83
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
integer, parameter, public gt_eargsizemismatch
Definition: dc_error.f90:536
種別型パラメタを提供します。
Definition: dc_types.f90:49
subroutine, public endsub(name, fmt, i, r, d, L, n, c1, c2, c3, ca)
Definition: dc_trace.f90:446
integer, parameter, public string
文字列を保持する 文字型変数の種別型パラメタ
Definition: dc_types.f90:118
Here is the call graph for this function:

◆ historycopy2()

subroutine historycopy2 ( type(gt_history), intent(out), target  hist_dest,
character(*), intent(in)  file,
type(gt_history), intent(in), optional, target  hist_src,
character(*), intent(in), optional  title,
character(*), intent(in), optional  source,
character(*), intent(in), optional  institution,
real, intent(in), optional  origin,
real, intent(in), optional  interval,
character(*), intent(in), optional  conventions,
character(*), intent(in), optional  gt_version 
)

Definition at line 137 of file historycopy.f90.

References dc_trace::beginsub(), dc_trace::dbgmessage(), and dc_trace::endsub().

137  !
138  ! 使用方法は HistoryCopy と同様です.
139  !
140  ! Usage is same as "HistoryCopy".
141  !
142  !--
143  ! 総称名 Copy として提供するための関数です.
144  ! 機能は HistoryCopy1 と同じです.
145  !++
146  use dc_trace, only: beginsub, endsub, dbgmessage
149  implicit none
150  type(gt_history), intent(out), target:: hist_dest
151  character(*), intent(in):: file
152  type(gt_history), intent(in), optional, target:: hist_src
153  character(*), intent(in), optional:: title, source, institution
154  real, intent(in), optional:: origin, interval
155  character(*), intent(in), optional:: conventions, gt_version
156 
157  ! Internal Work
158  character(len = *),parameter:: subname = "HistoryCopy2"
159  continue
160  call beginsub(subname, 'file=<%c>', c1=trim(file))
161  call historycopy(hist_dest, file, hist_src, &
162  & title, source, institution, &
163  & origin, interval, &
164  & conventions, gt_version)
165  call endsub(subname)
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 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: