Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
time_integrator_impl.hpp
Go to the documentation of this file.
1/*
2 * SPDX-FileCopyrightText: Copyright (c) 2014: G-CSC, Goethe University Frankfurt
3 * SPDX-License-Identifier: LicenseRef-UG4-LGPL-3.0
4 *
5 * Author: Arne Naegel, Andreas Kreienbuehl
6 *
7 * This file is part of UG4.
8 *
9 * UG4 is free software: you can redistribute it and/or modify it under the
10 * terms of the GNU Lesser General Public License version 3 (as published by the
11 * Free Software Foundation) with the following additional attribution
12 * requirements (according to LGPL/GPL v3 §7):
13 *
14 * (1) The following notice must be displayed in the Appropriate Legal Notices
15 * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
16 *
17 * (2) The following notice must be displayed at a prominent place in the
18 * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
19 *
20 * (3) The following bibliography is recommended for citation and must be
21 * preserved in all covered files:
22 * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
23 * parallel geometric multigrid solver on hierarchically distributed grids.
24 * Computing and visualization in science 16, 4 (2013), 151-164"
25 * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
26 * flexible software system for simulating pde based models on high performance
27 * computers. Computing and visualization in science 16, 4 (2013), 165-179"
28 *
29 * This program is distributed in the hope that it will be useful,
30 * but WITHOUT ANY WARRANTY; without even the implied warranty of
31 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
32 * GNU Lesser General Public License for more details.
33 */
34
35#ifndef __H__LIMEX__TIME_INTEGRATOR_IMPL_HPP__
36#define __H__LIMEX__TIME_INTEGRATOR_IMPL_HPP__
37
38namespace ug {
39
40template<typename TDomain, typename TAlgebra>
42(
44 number t1,
46 number t0
47)
48{
49
51 UG_COND_THROW(!base_type::m_spLinearSolver.valid(), "Linear solver invalid");
52
53 // short-cuts
54 GridLevel const &gl = u0->grid_level();
55 time_disc_type &tdisc = *tdisc_dep_type::m_spTimeDisc;
56
57 // create solution vector & right hand side
59 SmartPtr<typename base_type::vector_type> b= u0->clone_without_values();
60
61 // solution time series
62 SmartPtr<vector_time_series_type> m_spSolTimeSeries;
63 m_spSolTimeSeries=make_sp(new vector_time_series_type());
64 m_spSolTimeSeries->clear();
65 m_spSolTimeSeries->push(uold, t0);
66
67 // create matrix operator
68 SmartPtr<typename base_type::assembled_operator_type> spAssOp=make_sp(new typename base_type::assembled_operator_type(tdisc_dep_type::m_spTimeDisc, gl));
69
70 // integrate
71 if(!base_type::m_bNoLogOut)
72 UG_LOG("+++ Integrating: ["<< t0 <<", "<< t1 <<"]\n");
73
74 double t = t0;
75 number dt_assembled = -1.0; // invalid
76 int step = 1;
77
78 number currdt = base_type::m_dt;
79
80 while((t < t1) && (t1-t > base_type::m_precisionBound))
81 {
82
83 if(!base_type::m_bNoLogOut)
84 { UG_LOG("+++ Timestep +++" << step << "\n"); }
85
86 // determine step size
87 number dt = std::min(currdt, t1-t);
88
89 // prepare step
90 tdisc.prepare_step(m_spSolTimeSeries, dt);
91 if (fabs(dt-dt_assembled) > base_type::m_precisionBound)
92 {
93 // re-assemble operator
94 if(!base_type::m_bNoLogOut)
95 UG_LOG("+++ Reassemble (t=" << t << ", dt=" << dt <<")\n");
96
97 tdisc.assemble_linear(*spAssOp, *b, gl);
98 (base_type::m_spLinearSolver)->init(spAssOp, *u1);
99 dt_assembled = dt;
100 }
101 else
102 {
103 // keep old operator
104 tdisc.assemble_rhs(*b, gl);
105 }
106
107 // execute step
108 if (base_type::m_spLinearSolver->apply(*u1, *b))
109 {
110 // ACCEPTING:
111 // push updated solution into time series
112 t += dt;
113 SmartPtr<typename base_type::vector_type> tmp = m_spSolTimeSeries->oldest();
114 VecAssign(*tmp, *u1.template cast_dynamic<typename base_type::vector_type>());
115 m_spSolTimeSeries->push_discard_oldest(tmp, t);
116
117 this->notify_finalize_step(u1, step++, t+dt, dt);
118 }
119 else
120 {
121 // DISCARDING
122 currdt *= 0.5;
123 }
124
125 }
126
127 this->notify_end(u1, step, t1, currdt);
128
129 return true;
130};
131
132
133template<typename TDomain, typename TAlgebra>
135(
137 number t1,
139 number t0
140)
141{
143 UG_COND_THROW(!base_type::m_spLinearSolver.valid(), "Linear solver invalid");
144
145 // short-cuts
146 GridLevel const &gl = u0->grid_level();
147 time_disc_type &tdisc = *tdisc_dep_type::m_spTimeDisc;
148
149 // create solution vector & right hand side
151 SmartPtr<typename base_type::vector_type> b= u0->clone_without_values();
152
153 // solution time series
154 SmartPtr<vector_time_series_type> m_spSolTimeSeries;
155 m_spSolTimeSeries=make_sp(new vector_time_series_type());
156 m_spSolTimeSeries->push(uold, t0);
157
159
160 // select number of steps
161 double t = t0;
162 int numSteps = round((t1-t0) / base_type::m_dt);
163 number currdt = (t1-t0) / numSteps;
164
165
166
167 //std::cerr << "+++ Integrating: ["<< t0 <<", "<< t1 <<"] with dt=" << currdt << "("<< numSteps<< " iters)\n";
168 if(!base_type::m_bNoLogOut)
169 {
170 UG_LOG("+++ Integrating: [\t"<< t0 <<"\t, \t"<< t1 <<"\t] with dt=\t" << currdt << "("<< numSteps<< " iters)" << std::endl);
171 }
172
173 // integrate
174 for(int step = 1; step<=numSteps; ++step)
175 {
176 // determine step size
177 // number dt = std::min(currdt, t1-t);
178 const number dt = currdt;
179
180 if(!base_type::m_bNoLogOut)
181 {
182 UG_LOG("+++ Const timestep +++" << step<< "(t=" << t << ", dt=" << dt << ")"<< std::endl);
183 }
184 this->notify_init_step(u1, step, t, dt);
185
186 // prepare step
187 tdisc.prepare_step(m_spSolTimeSeries, dt);
188 if (spAssOp==SPNULL)
189 {
190 // Assemble operator.
191 if(!base_type::m_bNoLogOut) UG_LOG("+++ Assemble (t=" << t << ", dt=" << dt <<")" << std::endl);
192
193 spAssOp=make_sp(new typename base_type::assembled_operator_type(tdisc_dep_type::m_spTimeDisc, gl));
194 tdisc.assemble_linear(*spAssOp, *b, gl);
195 (base_type::m_spLinearSolver)->init(spAssOp, *u1);
196 }
197 else
198 {
199 // Recycle existing operator.
200 // std::cerr << "Recycling timestep " << step << "\n";
201 tdisc.assemble_rhs(*b, gl);
202 }
203
204 // execute step
205 if (base_type::m_spLinearSolver->apply(*u1, *b))
206 {
207 // ACCEPTING: push updated solution into time series
208 t += dt;
209 SmartPtr<typename base_type::vector_type> tmp = m_spSolTimeSeries->oldest();
210 VecAssign(*tmp, *u1.template cast_dynamic<typename base_type::vector_type>());
211 m_spSolTimeSeries->push_discard_oldest(tmp, t);
212 this->notify_finalize_step(u1, step, t, dt); // t updated already!
213 }
214 else
215 {
216 UG_THROW("ConstStepLinearTimeIntegrator::apply failed!!!");
217 this->notify_rewind_step(u1, step, t+dt, dt);
218 }
219
220 }
221
222 this->notify_end(u1, numSteps, t1, currdt);
223
224 return true;
225};
226
227template<typename TDomain, typename TAlgebra>
229(
231 number t1,
233 number t0
234)
235{
236 // short-cuts
237 GridLevel const &gl = u0->grid_level();
238 time_disc_type &tdisc = *tdisc_dep_type::m_spTimeDisc;
239 time_disc_type &tdisc2 = *m_spTimeDisc2;
240
241 // create solution vector & right hand side
243 SmartPtr<typename base_type::vector_type> b= u0->clone_without_values();
244
245 // create matrix operator
246 SmartPtr<typename base_type::assembled_operator_type> spAssOp=make_sp(new typename base_type::assembled_operator_type(tdisc_dep_type::m_spTimeDisc, gl));
247
248 // create additional solutions
251
252 // solution time series
253 SmartPtr<vector_time_series_type> m_spSolTimeSeries;
254 m_spSolTimeSeries=make_sp(new vector_time_series_type());
255 m_spSolTimeSeries->push(uold, t0);
256
257 SmartPtr<vector_time_series_type> m_spSolTimeSeries2;
258 m_spSolTimeSeries2=make_sp(new vector_time_series_type());
259 m_spSolTimeSeries2->push(u2old, t0);
260
261 // automatic selection of min/max time step
262 if (m_dtmin <= 0.0) { m_dtmin = (t1 - t0)/1.0e+5; }
263 if (m_dtmax <= 0.0) { m_dtmax = (t1 - t0)/10.0; }
264
265 // Aitken Neville extrapolation
266 const size_t tsteps[2] = {1,2};
267 std::vector<size_t> vsteps (tsteps, tsteps+2);
269
270 // integrate
271 if(!base_type::m_bNoLogOut)
272 UG_LOG("+++ Integrating: ["<< t0 <<", "<< t1 <<"]\n");
273
274 double t = t0;
275 int step = 0;
276
277 number dt = base_type::m_dt;
278 while((t < t1) && (t1-t > base_type::m_precisionBound))
279 {
280 // step: t -> t+dt
281 bool bSuccess = false;
282 while (!bSuccess){
283
284 // determine step size
285 if (dt<m_dtmin){
286 if(!base_type::m_bNoLogOut)
287 UG_LOG("Step size below minimum")
288 }
289 dt = std::min(dt, t1-t);
290
291 // basic step
292 if(!base_type::m_bNoLogOut)
293 UG_LOG("+++ Timestep: " << ++step << "\n");
294
295 tdisc.prepare_step(m_spSolTimeSeries, dt);
296 tdisc.assemble_linear(*spAssOp, *b, gl);
297 base_type::m_spLinearSolver->init(spAssOp, *u1);
298 base_type::m_spLinearSolver->apply(*u1, *b);
299
300
301 // control 1/2
302 if(!base_type::m_bNoLogOut)
303 UG_LOG("+++ Control: " << step << "\n");
304
305 tdisc2.prepare_step(m_spSolTimeSeries, 0.5*dt);
306 tdisc2.assemble_linear(*spAssOp, *b, gl);
307 base_type::m_spLinearSolver->init(spAssOp, *u2);
308 base_type::m_spLinearSolver->apply(*u2, *b);
309
310 // push back solution
311 SmartPtr<typename base_type::vector_type> tmp2 = m_spSolTimeSeries2->oldest();
312 (*tmp2)=*(u2.template cast_static<typename base_type::vector_type>());
313 m_spSolTimeSeries2->push_discard_oldest(tmp2, t + 0.5*dt);
314
315 // control 2/2
316 tdisc2.prepare_step(m_spSolTimeSeries2, 0.5*dt);
317 tdisc2.assemble_linear(*spAssOp, *b, gl);
318 base_type::m_spLinearSolver->init(spAssOp, *u2);
319 base_type::m_spLinearSolver->apply(*u2, *b);
320
321 // obtain extrapolated solution
322 timex.set_solution(u1, 0);
323 timex.set_solution(u2, 1);
324 timex.apply();
325
326 // predict (subsequent) time step
327 number eps = timex.get_error_estimates()[0];
328 number lambda = std::pow(0.8* m_tol/eps, 0.5);
329
330 number dtEst= dt*lambda;
331 dtEst = std::min(dtEst, 1.5*dt);
332 dtEst = std::min(dtEst, m_dtmax);
333
334 dt = dtEst;
335 if (eps <= m_tol)
336 {
337 // ACCEPT STEP (and thus solution u2)
338 if(!base_type::m_bNoLogOut)
339 UG_LOG("ACCEPTING solution, dtnew=" << dt);
340
341 bSuccess = true;
342 }
343 else
344 {
345 // DISCARD step
346 if(!base_type::m_bNoLogOut)
347 UG_LOG("DISCARDING solution, dtnew=" << dt);
348
349 // => reset solutions
350 VecAssign(*u1.template cast_dynamic<typename base_type::vector_type>(), *uold);
351
352 // => timeSeries2 has been updated...
353 SmartPtr<typename base_type::vector_type> tmp = m_spSolTimeSeries2->oldest();
354 VecAssign(*tmp, *uold);
355 m_spSolTimeSeries2->push_discard_oldest(tmp, t);
356
357 }
358
359 }
360
361
362 // prepare next loop
363 t += dt;
364
365 // push updated solution into time series (and continue)
366 SmartPtr<typename base_type::vector_type> tmp = m_spSolTimeSeries->oldest();
367 VecAssign(*tmp, static_cast<typename base_type::vector_type> (*u2));
368 m_spSolTimeSeries->push_discard_oldest(tmp, t);
369
370 }
371
372 return true;
373};
374
375} // namespace ug
376
377#endif /* __H__LIMEX__TIME_INTEGRATOR_IMPL_HPP__ */
Definition time_extrapolation.h:986
const std::vector< number > & get_error_estimates() const
Definition time_extrapolation.h:1030
void apply(size_t nstages, bool with_error=true)
best error estimate
Definition time_extrapolation.h:1065
void set_solution(SmartPtr< vector_type > soli, int i)
set solution (for stage i)
Definition time_extrapolation.h:1016
bool apply(SmartPtr< grid_function_type > u1, number t1, ConstSmartPtr< grid_function_type > u0, number t0)
Definition time_integrator_impl.hpp:135
void assemble_linear(matrix_type &A, vector_type &b)
void assemble_rhs(vector_type &b)
base_type::vector_type vector_type
Definition time_integrator.hpp:174
virtual void prepare_step(SmartPtr< VectorTimeSeries< vector_type > > prevSol, number dt)=0
bool apply(SmartPtr< grid_function_type > u1, number t1, ConstSmartPtr< grid_function_type > u0, number t0)
Definition time_integrator_impl.hpp:42
bool apply(SmartPtr< grid_function_type > u1, number t1, ConstSmartPtr< grid_function_type > u0, number t0)
Definition time_integrator_impl.hpp:229
void VecAssign(GPUVector< T > &dest, const GPUVector< T > &v1)
const NullSmartPtr SPNULL
#define UG_THROW(msg)
#define UG_LOG(msg)
#define UG_COND_THROW(cond, msg)
double number
#define LIMEX_PROFILE_FUNC()
Definition limex_tools.h:41
SmartPtr< T, FreePolicy > make_sp(T *inst)