ExaDG
Loading...
Searching...
No Matches
application_base.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef INCLUDE_EXADG_INCOMPRESSIBLE_FLOW_WITH_TRANSPORT_USER_INTERFACE_APPLICATION_BASE_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_FLOW_WITH_TRANSPORT_USER_INTERFACE_APPLICATION_BASE_H_
24
25// deal.II
26#include <deal.II/grid/grid_generator.h>
27#include <deal.II/grid/manifold_lib.h>
28#include <deal.II/grid/tria_description.h>
29
30// ExaDG
31#include <exadg/grid/grid.h>
32#include <exadg/grid/grid_utilities.h>
33
34#include <exadg/incompressible_navier_stokes/postprocessor/postprocessor.h>
35#include <exadg/incompressible_navier_stokes/user_interface/boundary_descriptor.h>
36#include <exadg/incompressible_navier_stokes/user_interface/field_functions.h>
37#include <exadg/incompressible_navier_stokes/user_interface/parameters.h>
38
39#include <exadg/convection_diffusion/postprocessor/postprocessor.h>
40#include <exadg/convection_diffusion/user_interface/boundary_descriptor.h>
41#include <exadg/convection_diffusion/user_interface/field_functions.h>
42#include <exadg/convection_diffusion/user_interface/parameters.h>
43
44#include <exadg/operators/resolution_parameters.h>
45#include <exadg/postprocessor/output_parameters.h>
46
47namespace ExaDG
48{
49namespace FTI
50{
51template<int dim, typename Number>
53{
54public:
55 virtual void
56 add_parameters(dealii::ParameterHandler & prm, std::vector<std::string> const & subsection_names)
57 {
58 for(auto & name : subsection_names)
59 {
60 prm.enter_subsection(name);
61 }
62
63 resolution.add_parameters(prm);
64
65 output_parameters.add_parameters(prm);
66
67 for(auto & name : subsection_names)
68 {
69 (void)name;
70 prm.leave_subsection();
71 }
72 }
73
74 FluidBase(std::string parameter_file, MPI_Comm const & comm)
75 : mpi_comm(comm),
76 pcout(std::cout, dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0),
77 parameter_file(parameter_file)
78 {
79 }
80
81 virtual ~FluidBase()
82 {
83 }
84
85 virtual void
86 setup(std::shared_ptr<Grid<dim>> & grid,
87 std::shared_ptr<dealii::Mapping<dim>> & mapping,
88 std::shared_ptr<MultigridMappings<dim, Number>> & multigrid_mappings,
89 std::vector<std::string> const & subsection_names)
90 {
91 parse_parameters(subsection_names);
92
93 // set resolution parameters
94 param.degree_u = resolution.degree;
95 param.grid.n_refine_global = resolution.refine_space;
96
97 set_parameters();
98 param.check(pcout);
99 param.print(pcout, "List of parameters:");
100
101 // grid
102 grid = std::make_shared<Grid<dim>>();
103 create_grid(*grid, mapping, multigrid_mappings);
104 print_grid_info(pcout, *grid);
105
106 // boundary conditions
107 boundary_descriptor = std::make_shared<IncNS::BoundaryDescriptor<dim>>();
108 set_boundary_descriptor();
109 verify_boundary_conditions<dim>(*boundary_descriptor, *grid);
110
111 // field functions
112 field_functions = std::make_shared<IncNS::FieldFunctions<dim>>();
113 set_field_functions();
114 }
115
116 virtual std::shared_ptr<IncNS::PostProcessorBase<dim, Number>>
117 create_postprocessor() = 0;
118
119 IncNS::Parameters const &
120 get_parameters() const
121 {
122 return param;
123 }
124
125 std::shared_ptr<IncNS::BoundaryDescriptor<dim> const>
126 get_boundary_descriptor() const
127 {
128 return boundary_descriptor;
129 }
130
131 std::shared_ptr<IncNS::FieldFunctions<dim> const>
132 get_field_functions() const
133 {
134 return field_functions;
135 }
136
137 // Analytical mesh motion
138 virtual std::shared_ptr<dealii::Function<dim>>
139 create_mesh_movement_function()
140 {
141 std::shared_ptr<dealii::Function<dim>> mesh_motion =
142 std::make_shared<dealii::Functions::ZeroFunction<dim>>(dim);
143
144 return mesh_motion;
145 }
146
147protected:
148 virtual void
149 parse_parameters(std::vector<std::string> const & subsection_names)
150 {
151 dealii::ParameterHandler prm;
152 this->add_parameters(prm, subsection_names);
153 prm.parse_input(parameter_file, "", true, true);
154 }
155
156 MPI_Comm const mpi_comm;
157
158 dealii::ConditionalOStream pcout;
159
160 IncNS::Parameters param;
161
162 std::shared_ptr<IncNS::FieldFunctions<dim>> field_functions;
163 std::shared_ptr<IncNS::BoundaryDescriptor<dim>> boundary_descriptor;
164
165 std::string parameter_file;
166
167 OutputParameters output_parameters;
168
169private:
170 virtual void
171 set_parameters() = 0;
172
173 virtual void
174 create_grid(Grid<dim> & grid,
175 std::shared_ptr<dealii::Mapping<dim>> & mapping,
176 std::shared_ptr<MultigridMappings<dim, Number>> & multigrid_mappings) = 0;
177
178 virtual void
179 set_boundary_descriptor() = 0;
180
181 virtual void
182 set_field_functions() = 0;
183
185};
186
187template<int dim, typename Number>
189{
190public:
191 virtual void
192 add_parameters(dealii::ParameterHandler & prm, std::vector<std::string> const & subsection_names)
193 {
194 for(auto & name : subsection_names)
195 {
196 prm.enter_subsection(name);
197 }
198
199 prm.enter_subsection("SpatialResolution");
200 {
201 prm.add_parameter("Degree",
202 degree,
203 "Polynomial degree of shape functions.",
204 dealii::Patterns::Integer(1),
205 true);
206 }
207 prm.leave_subsection();
208
209 output_parameters.add_parameters(prm);
210
211 for(auto & name : subsection_names)
212 {
213 (void)name;
214 prm.leave_subsection();
215 }
216 }
217
218 ScalarBase(std::string parameter_file, MPI_Comm const & comm)
219 : mpi_comm(comm),
220 pcout(std::cout, dealii::Utilities::MPI::this_mpi_process(comm) == 0),
221 parameter_file(parameter_file),
222 degree(1)
223 {
224 }
225
226 virtual ~ScalarBase()
227 {
228 }
229
230 void
231 setup(std::vector<std::string> const & subsection_names, Grid<dim> const & grid)
232 {
233 parse_parameters(subsection_names);
234
235 // set degree of shape functions (note that the refinement level is identical to fluid domain)
236 param.degree = degree;
237
238 set_parameters();
239 param.check();
240
241 param.print(this->pcout, "List of parameters for quantity " + subsection_names.back() + ":");
242
243 // boundary conditions
244 boundary_descriptor = std::make_shared<ConvDiff::BoundaryDescriptor<dim>>();
245 set_boundary_descriptor();
246 verify_boundary_conditions(*boundary_descriptor, grid);
247
248 // field functions
249 field_functions = std::make_shared<ConvDiff::FieldFunctions<dim>>();
250 set_field_functions();
251 }
252
253 virtual std::shared_ptr<ConvDiff::PostProcessorBase<dim, Number>>
254 create_postprocessor() = 0;
255
257 get_parameters() const
258 {
259 return param;
260 }
261
262 std::shared_ptr<ConvDiff::BoundaryDescriptor<dim> const>
263 get_boundary_descriptor() const
264 {
265 return boundary_descriptor;
266 }
267
268 std::shared_ptr<ConvDiff::FieldFunctions<dim> const>
269 get_field_functions() const
270 {
271 return field_functions;
272 }
273
274protected:
275 virtual void
276 parse_parameters(std::vector<std::string> const & subsection_names)
277 {
278 dealii::ParameterHandler prm;
279 this->add_parameters(prm, subsection_names);
280 prm.parse_input(parameter_file, "", true, true);
281 }
282
283 MPI_Comm const mpi_comm;
284
285 dealii::ConditionalOStream pcout;
286
288 std::shared_ptr<ConvDiff::FieldFunctions<dim>> field_functions;
289 std::shared_ptr<ConvDiff::BoundaryDescriptor<dim>> boundary_descriptor;
290
291 std::string parameter_file;
292
293 OutputParameters output_parameters;
294
295private:
296 virtual void
297 set_parameters() = 0;
298
299 virtual void
300 set_boundary_descriptor() = 0;
301
302 virtual void
303 set_field_functions() = 0;
304
305 unsigned int degree;
306};
307
308template<int dim, typename Number>
310{
311public:
312 void
313 add_parameters(dealii::ParameterHandler & prm)
314 {
315 AssertThrow(fluid.get(), dealii::ExcMessage("fluid has not been initialized."));
316
317 fluid->add_parameters(prm, {"Fluid"});
318
319 for(unsigned int i = 0; i < scalars.size(); ++i)
320 {
321 AssertThrow(scalars[i].get(),
322 dealii::ExcMessage("scalar[" + std::to_string(i) +
323 "] has not been initialized."));
324
325 scalars[i]->add_parameters(prm, {"Scalar" + std::to_string(i)});
326 }
327 }
328
329 ApplicationBase(std::string parameter_file, MPI_Comm const & comm)
330 : mpi_comm(comm),
331 pcout(std::cout, dealii::Utilities::MPI::this_mpi_process(comm) == 0),
332 parameter_file(parameter_file)
333 {
334 }
335
336 virtual ~ApplicationBase()
337 {
338 }
339
340 void
341 setup(std::shared_ptr<Grid<dim>> & grid,
342 std::shared_ptr<dealii::Mapping<dim>> & mapping,
343 std::shared_ptr<MultigridMappings<dim, Number>> & multigrid_mappings)
344 {
345 AssertThrow(fluid.get(), dealii::ExcMessage("fluid has not been initialized."));
346
347 // The fluid field is defined as the field that creates the grid and the mapping, while all
348 // scalar fields use the same grid/mapping
349 fluid->setup(grid, mapping, multigrid_mappings, {"Fluid"});
350
351 for(unsigned int i = 0; i < scalars.size(); ++i)
352 {
353 AssertThrow(scalars[i].get(),
354 dealii::ExcMessage("scalar[" + std::to_string(i) +
355 "] has not been initialized."));
356
357 scalars[i]->setup({"Scalar" + std::to_string(i)}, *grid);
358
359 // do additional parameter checks
360 AssertThrow(scalars[i]->get_parameters().ale_formulation ==
361 fluid->get_parameters().ale_formulation,
362 dealii::ExcMessage(
363 "Parameter ale_formulation is different for fluid field and scalar field"));
364
365 AssertThrow(
366 scalars[i]->get_parameters().adaptive_time_stepping ==
367 fluid->get_parameters().adaptive_time_stepping,
368 dealii::ExcMessage(
369 "The option adaptive_time_stepping has to be consistent for fluid and scalar transport solvers."));
370 }
371 }
372
373 std::shared_ptr<FluidBase<dim, Number>> fluid;
374 std::vector<std::shared_ptr<ScalarBase<dim, Number>>> scalars;
375
376protected:
377 MPI_Comm const mpi_comm;
378
379 dealii::ConditionalOStream pcout;
380
381private:
382 std::string parameter_file;
383};
384
385} // namespace FTI
386} // namespace ExaDG
387
388
389#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_FLOW_WITH_TRANSPORT_USER_INTERFACE_APPLICATION_BASE_H_ */
Definition parameters.h:46
Definition application_base.h:310
Definition application_base.h:53
Definition application_base.h:189
Definition grid.h:40
Definition parameters.h:46
Definition grid.h:84
Definition driver.cpp:33
Definition output_parameters.h:32
Definition resolution_parameters.h:81