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