Loading [MathJax]/extensions/tex2jax.js
Plugins
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
simple_integrator_impl.hpp
Go to the documentation of this file.
1/*
2 * SPDX-FileCopyrightText: Copyright (c) 2025: Goethe University Frankfurt
3 * SPDX-License-Identifier: LicenseRef-UG4-LGPL-3.0
4 *
5 * Author: Arne Naegel
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__SIMPLE_INTEGRATOR_IMPL_HPP__
36#define __H__LIMEX__SIMPLE_INTEGRATOR_IMPL_HPP__
37
38namespace ug {
39
40template<typename TDomain, typename TAlgebra>
42(
44 number t1,
46 number t0
47)
48{
50
51 // short-cuts
52 GridLevel const &gl = u0->grid_level();
53 time_disc_type &tdisc = *tdisc_dep_type::m_spTimeDisc;
54 typename base_type::solver_type &solver = *base_type::m_spSolver;
55
56 // create solution vector & right hand side
58
59 // init solution time series
60 SmartPtr<vector_time_series_type> m_spSolTimeSeries;
61 m_spSolTimeSeries=make_sp(new vector_time_series_type());
62 m_spSolTimeSeries->clear();
63 m_spSolTimeSeries->push(u0->clone(), t0);
64
65 // init solver (and matrix operator)
67 spAssOp = make_sp(new typename base_type::assembled_operator_type(tdisc_dep_type::m_spTimeDisc, gl));
68 solver.init(spAssOp);
69
70 // integrate
71 double t = t0;
72 number currdt = base_type::m_dt;
73 int step = 1;
74
75 double final_dt = base_type::m_dt;
76
77 if(!base_type::m_bNoLogOut)
78 UG_LOG("+++ Integrating: [\t"<< t0 <<"\t, \t"<< t1 <<"\t] with\t" << currdt <<"\n");
79
80 while(!hasTerminated(t, t0, t1)) {
81 if(!base_type::m_bNoLogOut)
82 UG_LOG("+++ Timestep +++" << step << "\n");
83
84 if (this->debug_writer_valid())
85 {
86 char debug_name_ext[16]; snprintf(debug_name_ext, 16, "%04d", step);
87 this->enter_debug_writer_section(std::string("SimpleTimeIntegrator_step") + debug_name_ext);
88 }
89
90 // determine step size
91 UG_COND_THROW(currdt < base_type::get_dt_min(), "Time step size below minimum. ABORTING!")
92 number dt = std::min(currdt, t1-t);
93 final_dt = dt;
94
95 // prepare step
96 tdisc.prepare_step(m_spSolTimeSeries, dt);
97 if (solver.prepare(*u1) == false)
98 {
99 if(!base_type::m_bNoLogOut)
100 UG_LOG("Initialization failed! RETRY");
101
102 currdt *= base_type::get_reduction_factor();
103 continue;
104 }
105 //UG_LOG("m_spSolTimeSeries.size="<< m_spSolTimeSeries->size());
106
107 // execute step
108 if (solver.apply(*u1))
109 {
110 //
111 // ACCEPT step
112 //
113
114 // post prcess (e.g. physics)
115 if(!base_type::m_bNoLogOut)
116 {
117 // m_spSolTimeSeries->oldest() actually holds a pointer to a grid function
118 // but as the time series does not know this, we have to cast ourselves
119 // SmartPtr<grid_function_type> tmpOld = m_spSolTimeSeries->oldest().template cast_static<grid_function_type>();
120 this->notify_finalize_step(u1, step, t+dt, dt);
121 }
122
123 // update time
124 t += dt;
125
126 // push updated solution into time series (and continue)
127 //SmartPtr<typename base_type::vector_type> utmp = m_spSolTimeSeries->oldest();
128 //VecAssign(*utmp, static_cast<typename base_type::vector_type> (*u1) );
129 uold = m_spSolTimeSeries->push_discard_oldest(u1->clone(), t).template cast_static<grid_function_type>();
130 }
131 else
132 {
133 //
134 // REJECT step
135 //
136 UG_LOG("Solution failed! RETRY");
137 currdt *= base_type::get_reduction_factor();
138 continue;
139 }
140
141 // consistency check
142 if (step == 1 && m_spDerivative.valid())
143 {
144 UG_LOG("Computing consistency error: "<< std::endl);
145 UG_ASSERT(static_cast<typename base_type::vector_type*> (&*u1) != &(*u0),
146 "Huhh: Different vectors required!");
147 *m_spDerivative = *u0;
148 m_initial_consistency_error = m_spBanachSpace->distance(*m_spDerivative, *u1);
149 }
150
151 if (this->debug_writer_valid())
152 {
153 this->leave_debug_writer_section();
154 char debug_name_ext[16]; snprintf(debug_name_ext, 16, "%04d", step);
155 this->write_debug(*u1, std::string("TimeSolution_step") + debug_name_ext);
156 }
157
158 step++;
159 // tdisc.finish_step_elem(m_spSolTimeSeries, dt);
160 }
161
162 if(base_type::m_bNoLogOut)
163 {
164 this->notify_finalize_step(u1, /*uold,*/ step, t, final_dt);
165 }
166
167
168 if (m_spDerivative.valid())
169 {
170
171 //
172 // approximate derivative (by forward difference)
173 //
174 UG_ASSERT(static_cast<typename base_type::vector_type*> (&*u1) != &(*uold),
175 "Huhh: Different vectors required!");
176
177 VecScaleAdd(static_cast<typename TAlgebra::vector_type&>(*m_spDerivative),
178 1.0/final_dt, static_cast<typename TAlgebra::vector_type&>(*u1),
179 -1.0/final_dt, static_cast<typename TAlgebra::vector_type&>(*uold));
180
181 }
182
183 m_spSolTimeSeries->clear();
184
185 return true;
186
187};
188
189template<typename TDomain, typename TAlgebra>
191(
193 number t1,
195 number t0
196)
197{
198
200
201 // short-cuts
202 GridLevel const &gl = u0->grid_level();
203 time_disc_type &tdisc = *tdisc_dep_type::m_spTimeDisc;
204 typename base_type::solver_type &solver = *base_type::m_spSolver;
205
206 //using TimeIntegratorSubject<TDomain,TAlgebra>::notify_step_postprocess;
207
208 // create solution vector & right hand side
209 SmartPtr<grid_function_type> uold = u0->clone();
210
211 // init solution time series
212 SmartPtr<vector_time_series_type> m_spSolTimeSeries;
213 m_spSolTimeSeries=make_sp(new vector_time_series_type());
214 m_spSolTimeSeries->clear();
215 m_spSolTimeSeries->push(u0->clone(), t0);
216
217 // init solver (and matrix operator)
218 SmartPtr<typename base_type::assembled_operator_type> spAssOp=make_sp(new typename base_type::assembled_operator_type(tdisc_dep_type::m_spTimeDisc, gl));
219 solver.init(spAssOp);
220
221 // integrate
222 double t = t0;
223 number currdt = base_type::m_dt;
224 int step = 1;
225
226 double final_dt = base_type::m_dt;
227
228 if(!base_type::m_bNoLogOut)
229 UG_LOG("+++ Integrating: ["<< t0 <<", "<< t1 <<"] with " << currdt <<"\n");
230
231 while(!hasTerminated(t, t0, t1))
232 {
233 if(!base_type::m_bNoLogOut)
234 UG_LOG("++++++ TIMESTEP " << step++ << " BEGIN (current time: " << t << ") ++++++\n");
235
236 // determine step size
237 UG_COND_THROW(currdt < base_type::get_dt_min(), "Time step size below minimum. ABORTING!")
238 number dt = std::min(currdt, t1-t);
239
240 final_dt = dt;
241
242 double told = t;
243
244 const int num_stages = tdisc.num_stages();
245 int s=1;
246 do // 'for-each-stage' loop
247 {
248
249 // for (int s=1; s<=num_stages; ++s)
250 if(!base_type::m_bNoLogOut)
251 UG_LOG("+++ STAGE "<< s << " BEGIN +++\n");
252
253 // set stage
254 tdisc.set_stage(s);
255
256 // prepare step
257 tdisc.prepare_step(m_spSolTimeSeries, dt);
258 if (solver.prepare(*u1) == false) break;
259
260 // execute step
261 if (!solver.apply(*u1)) break;
262
263 // stage was successful:
264 // a. update (intermediate) time
265 t = tdisc.future_time();
266
267 // b. push updated solution into time series (and continue)
268 SmartPtr<typename base_type::vector_type> oldest= m_spSolTimeSeries->oldest();
269 VecAssign(*oldest, static_cast<typename base_type::vector_type> (*u1));
270 m_spSolTimeSeries->push_discard_oldest(oldest, t);
271
272 // c. output
273 if(!base_type::m_bNoLogOut)
274 UG_LOG("+++ STAGE "<< s << " END +++\n");
275
276 } while ((++s) <=num_stages);
277
278
279 if (s<=num_stages)
280 {
281 // REJECT time step
282 if(!base_type::m_bNoLogOut)
283 UG_LOG("Solution failed! RETRY");
284
285 currdt *= this->get_reduction_factor();
286 t = told;
287 continue;
288 }
289 else
290 {
291 if(!base_type::m_bNoLogOut)
292 {
293 this->notify_finalize_step(u1, /*uold, */step, t, dt);
294 }
295
296 // ACCEPT time step
297 if (!hasTerminated(t, t0, t1))
298 *uold = *u1; // save solution (but not in last step)
299 // tdisc.finish_step_elem(m_spSolTimeSeries, dt);
300
301 if(!base_type::m_bNoLogOut)
302 UG_LOG("++++++ TIMESTEP " << step++ << " END (current time: " << t << ") ++++++\n");
303 }
304 }
305
306 if(base_type::m_bNoLogOut)
307 {
308 this->notify_finalize_step(u1, /*uold,*/ step, t, final_dt);
309 }
310 return true;
311};
312
313} // ug
314
315#endif // __H__LIMEX__SIMPLE_INTEGRATOR_IMPL_HPP__
parameterString s
Definition Biogas.lua:2
base_type::vector_type vector_type
Definition time_integrator.hpp:354
virtual bool apply(X &u)=0
virtual bool prepare(X &u)=0
virtual bool init(SmartPtr< IOperator< Y, X > > N)=0
virtual void set_stage(size_t stage)=0
virtual void prepare_step(SmartPtr< VectorTimeSeries< vector_type > > prevSol, number dt)=0
virtual number future_time() const=0
virtual size_t num_stages() const=0
bool apply_single_stage(SmartPtr< grid_function_type > u1, number t1, ConstSmartPtr< grid_function_type > u0, number t0)
Definition simple_integrator_impl.hpp:42
bool apply_multi_stage(SmartPtr< grid_function_type > u1, number t1, ConstSmartPtr< grid_function_type > u0, number t0)
Definition simple_integrator_impl.hpp:191
void VecAssign(GPUVector< T > &dest, const GPUVector< T > &v1)
#define UG_ASSERT(expr, msg)
#define UG_LOG(msg)
#define UG_COND_THROW(cond, msg)
double number
#define LIMEX_PROFILE_FUNC()
Definition limex_tools.h:41
void VecScaleAdd(double &dest, double alpha1, const double &v1, double alpha2, const double &v2)
SmartPtr< T, FreePolicy > make_sp(T *inst)