OGS
TESReactionAdaptor.cpp
Go to the documentation of this file.
1
11#include "TESReactionAdaptor.h"
12
13#include <cassert>
14
21
22namespace ProcessLib
23{
24namespace TES
25{
26std::unique_ptr<TESFEMReactionAdaptor> TESFEMReactionAdaptor::newInstance(
27 TESLocalAssemblerData const& data)
28{
29 auto const* ads = data.ap.react_sys.get();
30 if (dynamic_cast<Adsorption::AdsorptionReaction const*>(ads) != nullptr)
31 {
32 return std::make_unique<TESFEMReactionAdaptorAdsorption>(data);
33 }
34 if (dynamic_cast<Adsorption::ReactionInert const*>(ads) != nullptr)
35 {
36 return std::make_unique<TESFEMReactionAdaptorInert>(data);
37 }
38 if (dynamic_cast<Adsorption::ReactionSinusoidal const*>(ads) != nullptr)
39 {
40 return std::make_unique<TESFEMReactionAdaptorSinusoidal>(data);
41 }
42 if (dynamic_cast<Adsorption::ReactionCaOH2 const*>(ads) != nullptr)
43 {
44 return std::make_unique<TESFEMReactionAdaptorCaOH2>(data);
45 }
46
47 OGS_FATAL("No suitable TESFEMReactionAdaptor found. Aborting.");
48 return nullptr;
49}
50
52 TESLocalAssemblerData const& data)
53 // caution fragile: this relies in this constructor b eing called __after__
54 // data.solid_density has been properly set up!
55 : _bounds_violation(data.solid_density.size(), false), _d(data)
56{
57 assert(dynamic_cast<Adsorption::AdsorptionReaction const*>(
58 data.ap.react_sys.get()) != nullptr &&
59 "Reactive system has wrong type.");
60 assert(!_bounds_violation.empty());
61}
62
65 const unsigned int_pt)
66{
67 assert(_d.ap.number_of_try_of_iteration <= 20);
68
69 const double loading = Adsorption::AdsorptionReaction::getLoading(
71
72 double react_rate_R =
73 _d.ap.react_sys->getReactionRate(_d.p_V, _d.T, _d.ap.M_react, loading) *
75
76 // set reaction rate based on current damping factor
77 react_rate_R = (_reaction_damping_factor > 1e-3)
78 ? _reaction_damping_factor * react_rate_R
79 : 0.0;
80
81 if (_d.p_V <
83 _d.T) &&
84 react_rate_R > 0.0)
85 {
86 react_rate_R = 0.0;
87 }
88 else if (_d.p_V < 100.0 ||
90 getEquilibriumVapourPressure(_d.T))
91 {
92 // use equilibrium reaction for dry regime
93
94 // in the case of zeroth try in zeroth iteration: _p_V and loading are
95 // the values
96 // at the end of the previous timestep
97
98 const double pV_eq = estimateAdsorptionEquilibrium(_d.p_V, loading);
99 // TODO [CL]: it would be more correct to subtract pV from the previous
100 // timestep here
101 const double delta_pV = pV_eq - _d.p_V;
102 const double delta_rhoV =
103 delta_pV * _d.ap.M_react /
105 const double delta_rhoSR = delta_rhoV / (_d.ap.poro - 1.0);
106 double react_rate_R2 = delta_rhoSR / _d.ap.delta_t;
107
108 if (_bounds_violation[int_pt])
109 {
110 react_rate_R2 *= 0.5;
111 }
112
113 // 0th try: make sure reaction is not slower than allowed by local
114 // estimation
115 // nth try: make sure reaction is not made faster by local estimation
116 if ((_d.ap.number_of_try_of_iteration == 1 &&
117 std::abs(react_rate_R2) > std::abs(react_rate_R)) ||
119 std::abs(react_rate_R2) < std::abs(react_rate_R)))
120 {
121 react_rate_R = react_rate_R2;
122 }
123 }
124
125 // smooth out readjustment of reaction rate
127 {
129 {
130 // update reaction rate for for five iterations
131 const auto N = _d.ap.iteration_in_current_timestep - 4;
132
133 // take average s.t. does not oscillate so much
134 react_rate_R = 1.0 / (1.0 + N) *
135 (N * _d.reaction_rate[int_pt] + 1.0 * react_rate_R);
136 }
137 else
138 {
139 // afterwards no update anymore
140 react_rate_R = _d.reaction_rate[int_pt];
141 }
142 }
143
145 {
146 // assert that within tries reaction does not get faster
147 // (e.g. due to switch equilibrium reaction <--> kinetic reaction)
148
149 // factor of 0.9*N: in fact, even slow down reaction over tries
150 const double r = std::pow(0.9, _d.ap.number_of_try_of_iteration - 1) *
151 _d.reaction_rate[int_pt];
152 if (std::abs(react_rate_R) > std::abs(r))
153 {
154 react_rate_R = r;
155 }
156 }
157
158 return {react_rate_R,
159 _d.solid_density_prev_ts[int_pt] + react_rate_R * _d.ap.delta_t};
160}
161
163 const double p_V0, const double C0) const
164{
165 auto f = [this, p_V0, C0](double pV) -> double
166 {
167 // pV0 := _p_V
168 const double C_eq =
169 _d.ap.react_sys->getEquilibriumLoading(pV, _d.T, _d.ap.M_react);
170 return (pV - p_V0) * _d.ap.M_react /
172 _d.ap.poro +
173 (1.0 - _d.ap.poro) * (C_eq - C0) * _d.ap.rho_SR_dry;
174 };
175
176 // range where to search for roots of f
177 const double C_eq0 =
178 _d.ap.react_sys->getEquilibriumLoading(p_V0, _d.T, _d.ap.M_react);
179 const double limit =
180 (C_eq0 > C0)
181 ? 1e-8
183 _d.T);
184
185 // search for roots
186 auto rf = MathLib::Nonlinear::makeRegulaFalsi<MathLib::Nonlinear::Pegasus>(
187 f, p_V0, limit);
188 rf.step(3);
189
190 // set vapour pressure
191 return rf.getResult();
192}
193
195 std::vector<double> const& local_x,
196 std::vector<double> const& local_x_prev_ts)
197{
198 double alpha = 1.0;
199
200 const double min_xmV = 1e-6;
201 const std::size_t nnodes = local_x.size() / NODAL_DOF;
202 const std::size_t xmV_offset = COMPONENT_ID_MASS_FRACTION * nnodes;
203
204 for (std::size_t i = 0; i < nnodes; ++i)
205 {
206 auto const xnew = local_x[xmV_offset + i];
207 auto const xold = local_x_prev_ts[xmV_offset + i];
208
209 if (xnew < min_xmV)
210 {
211 const auto a = xold / (xold - xnew);
212 alpha = std::min(alpha, a);
213 _bounds_violation[i] = true;
214 }
215 else if (xnew > 1.0)
216 {
217 const auto a = xold / (xnew - xold);
218 alpha = std::min(alpha, a);
219 _bounds_violation[i] = true;
220 }
221 else
222 {
223 _bounds_violation[i] = false;
224 }
225 }
226
227 assert(alpha > 0.0);
228
229 if (alpha != 1.0)
230 {
231 if (alpha > 0.5)
232 {
233 alpha = 0.5;
234 }
235 if (alpha < 0.05)
236 {
237 alpha = 0.05;
238 }
239
241 {
242 _reaction_damping_factor *= sqrt(alpha);
243 }
244 else
245 {
247 }
248 }
249
250 return alpha == 1.0;
251}
252
263
269
271{
272 return {0.0, _d.solid_density_prev_ts[int_pt]};
273 // _d._qR = 0.0;
274 // _d._reaction_rate[int_pt] = 0.0;
275}
276
278 TESLocalAssemblerData const& data)
279 : _d(data)
280{
281 assert(dynamic_cast<Adsorption::ReactionSinusoidal const*>(
282 data.ap.react_sys.get()) != nullptr &&
283 "Reactive system has wrong type.");
284}
285
287 const unsigned /*int_pt*/)
288{
289 const double t = _d.ap.current_time;
290
291 // Cf. OGS5
292 const double rhoSR0 = 1.0;
293 const double rhoTil = 0.1;
294 const double omega = 2.0 * 3.1416;
295 const double poro = _d.ap.poro;
296
297 return {rhoTil * omega * cos(omega * t) / (1.0 - poro),
298 rhoSR0 + rhoTil * std::sin(omega * t) / (1.0 - poro)};
299}
300
302 TESLocalAssemblerData const& data)
303 : _d(data),
304 _react(dynamic_cast<Adsorption::ReactionCaOH2&>(*data.ap.react_sys))
305{
306 _ode_solver = MathLib::ODE::createODESolver<1>(_react.getOdeSolverConfig());
307 // TODO invalidate config
308
309 _ode_solver->setTolerance(1e-10, 1e-10);
310
311 auto f = [this](const double /*t*/,
314 ydot) -> bool
315 {
316 ydot[0] = _react.getReactionRate(y[0]);
317 return true;
318 };
319
320 _ode_solver->setFunction(f, nullptr);
321}
322
324{
325 // TODO if the first holds, the second also has to hold
328 {
329 return {_d.reaction_rate[int_pt], _d.solid_density[int_pt]};
330 }
331
332 // TODO: double check!
333 // const double xv_NR = SolidProp->non_reactive_solid_volume_fraction;
334 // const double rho_NR = SolidProp->non_reactive_solid_density;
335 const double xv_NR = 0.0;
336 const double rho_NR = 0.0;
337
338 const double t0 = 0.0;
339 const double y0 =
340 (_d.solid_density_prev_ts[int_pt] - xv_NR * rho_NR) / (1.0 - xv_NR);
341
342 const double t_end = _d.ap.delta_t;
343
345 _d.solid_density_prev_ts[int_pt]);
346
347 _ode_solver->setIC(t0, {y0});
348 _ode_solver->preSolve();
349 _ode_solver->solve(t_end);
350
351 const double time_reached = _ode_solver->getTime();
352 (void)time_reached;
353 assert(std::abs(t_end - time_reached) <
354 std::numeric_limits<double>::epsilon());
355
356 auto const& y_new = _ode_solver->getSolution();
357 auto const& y_dot_new = _ode_solver->getYDot(t_end, y_new);
358
359 double rho_react;
360
361 // cut off when limits are reached
363 {
365 }
366 else if (y_new[0] >
368 {
370 }
371 else
372 {
373 rho_react = y_new[0];
374 }
375
376 return {y_dot_new[0] * (1.0 - xv_NR),
377 (1.0 - xv_NR) * rho_react + xv_NR * rho_NR};
378}
379
380} // namespace TES
381} // namespace ProcessLib
#define OGS_FATAL(...)
Definition Error.h:26
static double getEquilibriumVapourPressure(const double T_Ads)
static double getLoading(const double rho_curr, const double rho_dry)
static MATERIALLIB_EXPORT constexpr double rho_low
lower density limit
void updateParam(double T_solid, double p_gas, double x_react, double rho_s_initial)
static MATERIALLIB_EXPORT constexpr double rho_up
upper density limit
double getReactionRate(const double, const double, const double, const double) const override
const BaseLib::ConfigTree & getOdeSolverConfig() const
ReactionRate initReaction_slowDownUndershootStrategy(const unsigned int_pt)
TESFEMReactionAdaptorAdsorption(TESLocalAssemblerData const &data)
double estimateAdsorptionEquilibrium(const double p_V0, const double C0) const
bool checkBounds(std::vector< double > const &local_x, std::vector< double > const &local_x_prev_ts) override
ReactionRate initReaction(const unsigned) override
TESFEMReactionAdaptorCaOH2(TESLocalAssemblerData const &data)
std::unique_ptr< MathLib::ODE::ODESolver< 1 > > _ode_solver
ReactionRate initReaction(const unsigned int_pt) override
TESFEMReactionAdaptorInert(TESLocalAssemblerData const &)
ReactionRate initReaction(const unsigned) override
TESFEMReactionAdaptorSinusoidal(TESLocalAssemblerData const &data)
static std::unique_ptr< TESFEMReactionAdaptor > newInstance(TESLocalAssemblerData const &data)
MappedConstMatrix< N, 1 > MappedConstVector
MappedMatrix< N, 1 > MappedVector
const unsigned NODAL_DOF
const unsigned COMPONENT_ID_MASS_FRACTION
std::unique_ptr< Adsorption::Reaction > react_sys