ExaDG
Loading...
Searching...
No Matches
solver.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_CONVECTION_DIFFUSION_SOLVER_H_
23#define INCLUDE_EXADG_CONVECTION_DIFFUSION_SOLVER_H_
24
25// deal.II
26#include <deal.II/base/exceptions.h>
27#include <deal.II/base/parameter_handler.h>
28
29// driver
30#include <exadg/convection_diffusion/driver.h>
31
32// utilities
33#include <exadg/operators/resolution_parameters.h>
34#include <exadg/time_integration/resolution_parameters.h>
35#include <exadg/utilities/enum_patterns.h>
36#include <exadg/utilities/general_parameters.h>
37
38// application
39#include <exadg/convection_diffusion/user_interface/declare_get_application.h>
40
41namespace ExaDG
42{
43void
44create_input_file(std::string const & input_file)
45{
46 dealii::ParameterHandler prm;
47
48 GeneralParameters general;
49 general.add_parameters(prm);
50
51 SpatialResolutionParametersMinMax spatial;
52 spatial.add_parameters(prm);
53
54 TemporalResolutionParameters temporal;
55 temporal.add_parameters(prm);
56
57 // we have to assume a default dimension and default Number type
58 // for the automatic generation of a default input file
59 unsigned int const Dim = 2;
60 typedef double Number;
61 ConvDiff::get_application<Dim, Number>(input_file, MPI_COMM_WORLD)->add_parameters(prm);
62
63 prm.print_parameters(input_file,
64 dealii::ParameterHandler::Short |
65 dealii::ParameterHandler::KeepDeclarationOrder);
66}
67
68template<int dim, typename Number>
69void
70run(std::string const & input_file,
71 unsigned int const degree,
72 unsigned int const refine_space,
73 unsigned int const refine_time,
74 MPI_Comm const & mpi_comm,
75 bool const is_test)
76{
77 dealii::Timer timer;
78 timer.restart();
79
80 std::shared_ptr<ConvDiff::ApplicationBase<dim, Number>> application =
81 ConvDiff::get_application<dim, Number>(input_file, mpi_comm);
82
83 application->set_parameters_convergence_study(degree, refine_space, refine_time);
84
85 std::shared_ptr<ConvDiff::Driver<dim, Number>> driver =
86 std::make_shared<ConvDiff::Driver<dim, Number>>(mpi_comm, application, is_test, false);
87
88 driver->setup();
89
90 driver->solve();
91
92 if(not(is_test))
93 driver->print_performance_results(timer.wall_time());
94}
95
96} // namespace ExaDG
97
98int
99main(int argc, char ** argv)
100{
101 dealii::Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
102
103 MPI_Comm mpi_comm(MPI_COMM_WORLD);
104
105 std::string input_file;
106
107 if(argc == 1)
108 {
109 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
110 {
111 // clang-format off
112 std::cout << "To run the program, use: ./solver input_file" << std::endl
113 << "To setup the input file, use: ./solver input_file --help" << std::endl;
114 // clang-format on
115 }
116
117 return 0;
118 }
119 else if(argc >= 2)
120 {
121 input_file = std::string(argv[1]);
122
123 if(argc == 3 and std::string(argv[2]) == "--help")
124 {
125 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
126 ExaDG::create_input_file(input_file);
127
128 return 0;
129 }
130 }
131
132 ExaDG::GeneralParameters general(input_file);
134 ExaDG::TemporalResolutionParameters temporal(input_file);
135
136 // k-refinement
137 for(unsigned int degree = spatial.degree_min; degree <= spatial.degree_max; ++degree)
138 {
139 // h-refinement
140 for(unsigned int refine_space = spatial.refine_space_min;
141 refine_space <= spatial.refine_space_max;
142 ++refine_space)
143 {
144 // dt-refinement
145 for(unsigned int refine_time = temporal.refine_time_min;
146 refine_time <= temporal.refine_time_max;
147 ++refine_time)
148 {
149 // run the simulation
150 if(general.dim == 2 and general.precision == "float")
151 {
152 ExaDG::run<2, float>(
153 input_file, degree, refine_space, refine_time, mpi_comm, general.is_test);
154 }
155 else if(general.dim == 2 and general.precision == "double")
156 {
157 ExaDG::run<2, double>(
158 input_file, degree, refine_space, refine_time, mpi_comm, general.is_test);
159 }
160 else if(general.dim == 3 and general.precision == "float")
161 {
162 ExaDG::run<3, float>(
163 input_file, degree, refine_space, refine_time, mpi_comm, general.is_test);
164 }
165 else if(general.dim == 3 and general.precision == "double")
166 {
167 ExaDG::run<3, double>(
168 input_file, degree, refine_space, refine_time, mpi_comm, general.is_test);
169 }
170 else
171 {
172 AssertThrow(false,
173 dealii::ExcMessage("Only dim = 2|3 and precision=float|double implemented."));
174 }
175 }
176 }
177 }
178
179 return 0;
180}
181
182#endif /* INCLUDE_EXADG_CONVECTION_DIFFUSION_SOLVER_H_ */
Definition driver.cpp:33
Definition general_parameters.h:32
Definition resolution_parameters.h:32
Definition resolution_parameters.h:32