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