-
Notifications
You must be signed in to change notification settings - Fork 9
/
ThermalPhysics.hh
365 lines (316 loc) · 11.1 KB
/
ThermalPhysics.hh
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
/* Copyright (c) 2016 - 2024, the adamantine authors.
*
* This file is subject to the Modified BSD License and may not be distributed
* without copyright and license information. Please refer to the file LICENSE
* for the text and further information on this license.
*/
#ifndef THERMAL_PHYSICS_HH
#define THERMAL_PHYSICS_HH
#include <Geometry.hh>
#include <HeatSource.hh>
#include <ImplicitOperator.hh>
#include <ThermalOperatorBase.hh>
#include <ThermalPhysicsInterface.hh>
#include <deal.II/base/time_stepping.h>
#include <deal.II/base/time_stepping.templates.h>
#include <deal.II/distributed/cell_weights.h>
#include <deal.II/hp/fe_collection.h>
#include <boost/property_tree/ptree.hpp>
#include <memory>
namespace adamantine
{
/**
* This class takes care of building the linear operator and the
* right-hand-side. Also used to evolve the system in time.
*/
template <int dim, int fe_degree, typename MemorySpaceType,
typename QuadratureType>
class ThermalPhysics : public ThermalPhysicsInterface<dim, MemorySpaceType>
{
public:
/**
* Constructor.
*/
ThermalPhysics(MPI_Comm const &communicator,
boost::property_tree::ptree const &database,
Geometry<dim> &geometry,
MaterialProperty<dim, MemorySpaceType> &material_properties);
void setup() override;
void setup_dofs() override;
void compute_inverse_mass_matrix() override;
void add_material(
std::vector<std::vector<
typename dealii::DoFHandler<dim>::active_cell_iterator>> const
&elements_to_activate,
std::vector<double> const &new_deposition_cos,
std::vector<double> const &new_deposition_sin,
std::vector<bool> &new_has_melted, unsigned int const activation_start,
unsigned int const activation_end, double const initial_temperature,
dealii::LA::distributed::Vector<double, MemorySpaceType> &solution)
override;
/**
* For ThermalPhysics, update_physics_parameters is used to modify the heat
* sources in the middle of a simulation, e.g. for data assimilation with an
* augmented ensemble involving heat source parameters.
*/
void update_physics_parameters(
boost::property_tree::ptree const &heat_source_database) override;
double evolve_one_time_step(
double t, double delta_t,
dealii::LA::distributed::Vector<double, MemorySpaceType> &solution,
std::vector<Timer> &timers) override;
double get_delta_t_guess() const override;
void
initialize_dof_vector(double const value,
dealii::LA::distributed::Vector<double, MemorySpaceType>
&vector) const override;
void get_state_from_material_properties() override;
void set_state_to_material_properties() override;
void load_checkpoint(std::string const &filename,
dealii::LA::distributed::Vector<double, MemorySpaceType>
&temperature) override;
void save_checkpoint(std::string const &filename,
dealii::LA::distributed::Vector<double, MemorySpaceType>
&temperature) override;
void set_material_deposition_orientation(
std::vector<double> const &deposition_cos,
std::vector<double> const &deposition_sin) override;
double get_deposition_cos(unsigned int const i) const override;
double get_deposition_sin(unsigned int const i) const override;
void mark_has_melted(double const threshold_temperature,
dealii::LA::distributed::Vector<double, MemorySpaceType>
&temperature) override;
std::vector<bool> get_has_melted_vector() const override;
void set_has_melted_vector(std::vector<bool> const &has_melted) override;
bool get_has_melted(unsigned int const i) const override;
dealii::DoFHandler<dim> &get_dof_handler() override;
dealii::AffineConstraints<double> &get_affine_constraints() override;
std::vector<std::shared_ptr<HeatSource<dim>>> &get_heat_sources() override;
unsigned int get_fe_degree() const override;
/**
* Return the current height of the heat source.
*/
double get_current_source_height() const;
private:
using LA_Vector =
typename dealii::LA::distributed::Vector<double, MemorySpaceType>;
/**
* Update the depostion cosine and sine from the Physics object to the
* operator object.
*/
void update_material_deposition_orientation();
/**
* Compute the right-hand side and apply the TermalOperator.
*/
LA_Vector evaluate_thermal_physics(double const t, LA_Vector const &y,
std::vector<Timer> &timers) const;
/**
* Compute the inverse of the ImplicitOperator.
*/
LA_Vector id_minus_tau_J_inverse(double const t, double const tau,
LA_Vector const &y,
std::vector<Timer> &timers) const;
/**
* This flag is true if the time stepping method is embedded.
*/
bool _embedded_method = false;
/**
* This flag is true if the time stepping method is implicit.
*/
bool _implicit_method = false;
/**
* This flag is true if right preconditioning is used to invert the
* ImplicitOperator.
*/
bool _right_preconditioning;
/**
* Maximum number of iterations to invert the ImplicitOperator.
*/
unsigned int _max_iter;
/**
* Maximum number of temporary vectors when inverting the ImplicitOperator.
*/
unsigned int _max_n_tmp_vectors;
/**
* Guess of the next time step.
*/
double _delta_t_guess;
/**
* Tolerance to inverte the ImplicitOperator.
*/
double _tolerance;
/**
* Current height of the object.
*/
double _current_source_height = 0.;
/**
* Type of boundary.
*/
BoundaryType _boundary_type;
/**
* Associated geometry.
*/
Geometry<dim> &_geometry;
/**
* Associated Lagrange finite elements.
*/
dealii::hp::FECollection<dim> _fe_collection;
/**
* Associated DoFHandler.
*/
dealii::DoFHandler<dim> _dof_handler;
/**
* Associated AffineConstraints<double>.
*/
dealii::AffineConstraints<double> _affine_constraints;
/**
* Associated quadature, either Gauss or Gauss-Lobatto.
*/
dealii::hp::QCollection<1> _q_collection;
/**
* Object used to attach to each cell, a weight (used for load balancing)
* equal to the number of degrees of freedom associated with the cell.
*/
dealii::parallel::CellWeights<dim> _cell_weights;
/**
* Cosine of the material deposition angles.
*/
std::vector<double> _deposition_cos;
/**
* Sine of the material deposition angles.
*/
std::vector<double> _deposition_sin;
/**
* Indicator variable for whether a point has ever been above the solidus. The
* value is false for material that has not yet melted and true for material
* that has melted.
*/
std::vector<bool> _has_melted;
/**
* Associated material properties.
*/
MaterialProperty<dim, MemorySpaceType> &_material_properties;
/**
* Vector of heat sources.
*/
std::vector<std::shared_ptr<HeatSource<dim>>> _heat_sources;
/**
* Shared pointer to the underlying ThermalOperator.
*/
std::shared_ptr<ThermalOperatorBase<dim, MemorySpaceType>> _thermal_operator;
/**
* Unique pointer to the underlying ImplicitOperator.
*/
std::unique_ptr<ImplicitOperator<MemorySpaceType>> _implicit_operator;
/**
* Shared pointer to the underlying time stepping scheme.
*/
std::unique_ptr<dealii::TimeStepping::RungeKutta<LA_Vector>> _time_stepping;
};
template <int dim, int fe_degree, typename MemorySpaceType,
typename QuadratureType>
inline double ThermalPhysics<dim, fe_degree, MemorySpaceType,
QuadratureType>::get_delta_t_guess() const
{
return _delta_t_guess;
}
template <int dim, int fe_degree, typename MemorySpaceType,
typename QuadratureType>
inline void
ThermalPhysics<dim, fe_degree, MemorySpaceType,
QuadratureType>::update_material_deposition_orientation()
{
_thermal_operator->set_material_deposition_orientation(_deposition_cos,
_deposition_sin);
}
template <int dim, int fe_degree, typename MemorySpaceType,
typename QuadratureType>
inline void ThermalPhysics<dim, fe_degree, MemorySpaceType, QuadratureType>::
set_material_deposition_orientation(
std::vector<double> const &deposition_cos,
std::vector<double> const &deposition_sin)
{
_deposition_cos = deposition_cos;
_deposition_sin = deposition_sin;
update_material_deposition_orientation();
}
template <int dim, int fe_degree, typename MemorySpaceType,
typename QuadratureType>
inline double
ThermalPhysics<dim, fe_degree, MemorySpaceType,
QuadratureType>::get_deposition_cos(unsigned int const i) const
{
return _deposition_cos[i];
}
template <int dim, int fe_degree, typename MemorySpaceType,
typename QuadratureType>
inline double
ThermalPhysics<dim, fe_degree, MemorySpaceType,
QuadratureType>::get_deposition_sin(unsigned int const i) const
{
return _deposition_sin[i];
}
template <int dim, int fe_degree, typename MemorySpaceType,
typename QuadratureType>
inline std::vector<bool>
ThermalPhysics<dim, fe_degree, MemorySpaceType,
QuadratureType>::get_has_melted_vector() const
{
return _has_melted;
}
template <int dim, int fe_degree, typename MemorySpaceType,
typename QuadratureType>
inline void ThermalPhysics<dim, fe_degree, MemorySpaceType, QuadratureType>::
set_has_melted_vector(std::vector<bool> const &has_melted)
{
_has_melted = has_melted;
}
template <int dim, int fe_degree, typename MemorySpaceType,
typename QuadratureType>
inline bool
ThermalPhysics<dim, fe_degree, MemorySpaceType, QuadratureType>::get_has_melted(
unsigned int const i) const
{
return _has_melted[i];
}
template <int dim, int fe_degree, typename MemorySpaceType,
typename QuadratureType>
inline dealii::DoFHandler<dim> &
ThermalPhysics<dim, fe_degree, MemorySpaceType,
QuadratureType>::get_dof_handler()
{
return _dof_handler;
}
template <int dim, int fe_degree, typename MemorySpaceType,
typename QuadratureType>
inline dealii::AffineConstraints<double> &
ThermalPhysics<dim, fe_degree, MemorySpaceType,
QuadratureType>::get_affine_constraints()
{
return _affine_constraints;
}
template <int dim, int fe_degree, typename MemorySpaceType,
typename QuadratureType>
inline std::vector<std::shared_ptr<HeatSource<dim>>> &
ThermalPhysics<dim, fe_degree, MemorySpaceType,
QuadratureType>::get_heat_sources()
{
return _heat_sources;
}
template <int dim, int fe_degree, typename MemorySpaceType,
typename QuadratureType>
inline unsigned int
ThermalPhysics<dim, fe_degree, MemorySpaceType, QuadratureType>::get_fe_degree()
const
{
return fe_degree;
}
template <int dim, int fe_degree, typename MemorySpaceType,
typename QuadratureType>
inline double ThermalPhysics<dim, fe_degree, MemorySpaceType,
QuadratureType>::get_current_source_height() const
{
return _current_source_height;
}
} // namespace adamantine
#endif