ExaDG
Loading...
Searching...
No Matches
driver.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_NAVIER_STOKES_DRIVER_PRECURSOR_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_DRIVER_PRECURSOR_H_
24
25#include <exadg/functions_and_boundary_conditions/verify_boundary_conditions.h>
26#include <exadg/incompressible_navier_stokes/postprocessor/postprocessor_base.h>
27#include <exadg/incompressible_navier_stokes/precursor/user_interface/application_base.h>
28#include <exadg/incompressible_navier_stokes/spatial_discretization/create_operator.h>
29#include <exadg/incompressible_navier_stokes/spatial_discretization/operator_coupled.h>
30#include <exadg/incompressible_navier_stokes/spatial_discretization/operator_dual_splitting.h>
31#include <exadg/incompressible_navier_stokes/spatial_discretization/operator_pressure_correction.h>
32#include <exadg/incompressible_navier_stokes/time_integration/create_time_integrator.h>
33#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf_coupled_solver.h>
34#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf_dual_splitting.h>
35#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf_pressure_correction.h>
36#include <exadg/matrix_free/matrix_free_data.h>
37#include <exadg/utilities/print_general_infos.h>
38
39namespace ExaDG
40{
41namespace IncNS
42{
43namespace Precursor
44{
45template<int dim, typename Number>
46class Solver
47{
48public:
49 void
50 setup(std::shared_ptr<Domain<dim, Number>> & domain,
51 std::vector<std::string> const & subsection_names_parameters,
52 std::string const & field,
53 MPI_Comm const & mpi_comm,
54 bool const is_test)
55 {
56 // setup application
57 domain->setup(grid, mapping, multigrid_mappings, subsection_names_parameters);
58
59 // ALE is not used for this solver
60 std::shared_ptr<HelpersALE<dim, Number>> helpers_ale_dummy;
61
62 // initialize pde_operator
63 pde_operator = create_operator<dim, Number>(grid,
64 mapping,
65 multigrid_mappings,
66 domain->get_boundary_descriptor(),
67 domain->get_field_functions(),
68 domain->get_parameters(),
69 field,
70 mpi_comm);
71
72 // initialize matrix_free
73 matrix_free_data = std::make_shared<MatrixFreeData<dim, Number>>();
74 matrix_free_data->append(pde_operator);
75
76 matrix_free = std::make_shared<dealii::MatrixFree<dim, Number>>();
77 if(domain->get_parameters().use_cell_based_face_loops)
78 Categorization::do_cell_based_loops(*grid->triangulation, matrix_free_data->data);
79 matrix_free->reinit(*mapping,
80 matrix_free_data->get_dof_handler_vector(),
81 matrix_free_data->get_constraint_vector(),
82 matrix_free_data->get_quadrature_vector(),
83 matrix_free_data->data);
84
85 // setup Navier-Stokes operator
86 pde_operator->setup(matrix_free, matrix_free_data);
87
88 // setup postprocessor
89 postprocessor = domain->create_postprocessor();
90 postprocessor->setup(*pde_operator);
91
92 // Setup time integrator
93 time_integrator = create_time_integrator<dim, Number>(
94 pde_operator, helpers_ale_dummy, postprocessor, domain->get_parameters(), mpi_comm, is_test);
95
96 time_integrator->setup(domain->get_parameters().restarted_simulation);
97 }
98
99 /*
100 * Grid and mapping
101 */
102 std::shared_ptr<Grid<dim>> grid;
103 std::shared_ptr<dealii::Mapping<dim>> mapping;
104
105 std::shared_ptr<MultigridMappings<dim, Number>> multigrid_mappings;
106
107 /*
108 * Spatial discretization
109 */
110 std::shared_ptr<SpatialOperatorBase<dim, Number>> pde_operator;
111
112 /*
113 * Postprocessor
114 */
116
117 std::shared_ptr<Postprocessor> postprocessor;
118
119 /*
120 * Temporal discretization
121 */
122 std::shared_ptr<TimeIntBDF<dim, Number>> time_integrator;
123
124private:
125 /*
126 * MatrixFree
127 */
128 std::shared_ptr<MatrixFreeData<dim, Number>> matrix_free_data;
129 std::shared_ptr<dealii::MatrixFree<dim, Number>> matrix_free;
130};
131
132template<int dim, typename Number>
134{
135public:
136 Driver(MPI_Comm const & mpi_comm,
137 std::shared_ptr<ApplicationBase<dim, Number>> application,
138 bool const is_test);
139
140 void
141 setup();
142
143 void
144 solve() const;
145
146 void
147 print_performance_results(double const total_time) const;
148
149private:
150 void
151 set_start_time() const;
152
153 void
154 synchronize_time_step_size() const;
155
156 void
157 consistency_checks() const;
158
159 // MPI communicator
160 MPI_Comm const mpi_comm;
161
162 // output to std::cout
163 dealii::ConditionalOStream pcout;
164
165 // do not print wall times if is_test
166 bool const is_test;
167
168 // application
169 std::shared_ptr<ApplicationBase<dim, Number>> application;
170
171 Solver<dim, Number> solver_main, solver_precursor;
172
173 bool use_adaptive_time_stepping;
174
175 /*
176 * Computation time (wall clock time).
177 */
178 mutable TimerTree timer_tree;
179};
180
181} // namespace Precursor
182} // namespace IncNS
183} // namespace ExaDG
184
185#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_DRIVER_PRECURSOR_H_ */
Definition postprocessor_base.h:39
Definition application_base.h:175
Definition application_base.h:51
Definition driver.h:134
Definition driver.h:47
Definition timer_tree.h:36
Definition driver.cpp:33