ExaDG
Loading...
Searching...
No Matches
throughput.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_THROUGHPUT_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_THROUGHPUT_H_
24
25// likwid
26#ifdef EXADG_WITH_LIKWID
27# include <likwid.h>
28#endif
29
30// ExaDG
31
32// driver
33#include <exadg/incompressible_navier_stokes/driver.h>
34
35// utilities
36#include <exadg/operators/hypercube_resolution_parameters.h>
37#include <exadg/operators/throughput_parameters.h>
38#include <exadg/utilities/general_parameters.h>
39
40// application
41#include <exadg/incompressible_navier_stokes/user_interface/declare_get_application.h>
42
43namespace ExaDG
44{
45void
46create_input_file(std::string const & input_file)
47{
48 dealii::ParameterHandler prm;
49
50 GeneralParameters general;
51 general.add_parameters(prm);
52
53 HypercubeResolutionParameters resolution;
54 resolution.add_parameters(prm);
55
56 ThroughputParameters<IncNS::OperatorType> throughput;
57 throughput.add_parameters(prm);
58
59 try
60 {
61 // we have to assume a default dimension and default Number type
62 // for the automatic generation of a default input file
63 unsigned int const Dim = 2;
64 typedef double Number;
65 IncNS::get_application<Dim, Number>(input_file, MPI_COMM_WORLD)->add_parameters(prm);
66 }
67 catch(...)
68 {
69 }
70
71 prm.print_parameters(input_file,
72 dealii::ParameterHandler::Short |
73 dealii::ParameterHandler::KeepDeclarationOrder);
74}
75
76template<int dim, typename Number>
77void
78run(ThroughputParameters<IncNS::OperatorType> const & throughput,
79 std::string const & input_file,
80 unsigned int const degree,
81 unsigned int const refine_space,
82 unsigned int const n_cells_1d,
83 MPI_Comm const & mpi_comm,
84 bool const is_test)
85{
86 std::shared_ptr<IncNS::ApplicationBase<dim, Number>> application =
87 IncNS::get_application<dim, Number>(input_file, mpi_comm);
88
89 application->set_parameters_throughput_study(degree, refine_space, n_cells_1d);
90
91 std::shared_ptr<IncNS::Driver<dim, Number>> driver =
92 std::make_shared<IncNS::Driver<dim, Number>>(mpi_comm, application, is_test, true);
93
94 driver->setup();
95
96 std::tuple<unsigned int, dealii::types::global_dof_index, double> wall_time =
97 driver->apply_operator(throughput.operator_type,
98 throughput.n_repetitions_inner,
99 throughput.n_repetitions_outer);
100
101 throughput.wall_times.push_back(wall_time);
102}
103} // namespace ExaDG
104
105int
106main(int argc, char ** argv)
107{
108#ifdef EXADG_WITH_LIKWID
109 LIKWID_MARKER_INIT;
110#endif
111
112 dealii::Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
113
114 MPI_Comm mpi_comm(MPI_COMM_WORLD);
115
116 std::string input_file;
117
118 if(argc == 1)
119 {
120 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
121 {
122 // clang-format off
123 std::cout << "To run the program, use: ./throughput input_file" << std::endl
124 << "To setup the input file, use: ./throughput input_file --help" << std::endl;
125 // clang-format on
126 }
127
128 return 0;
129 }
130 else if(argc >= 2)
131 {
132 input_file = std::string(argv[1]);
133
134 if(argc == 3 and std::string(argv[2]) == "--help")
135 {
136 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
137 ExaDG::create_input_file(input_file);
138
139 return 0;
140 }
141 }
142
143 ExaDG::GeneralParameters general(input_file);
144 ExaDG::HypercubeResolutionParameters resolution(input_file, general.dim);
146
147 ExaDG::IncNS::PressureDegree pressure_degree = ExaDG::IncNS::PressureDegree::MixedOrder;
148
149 dealii::ParameterHandler prm;
150
151 prm.enter_subsection("Throughput");
152 {
153 prm.add_parameter("PressureDegree",
154 pressure_degree,
155 "Degree of pressure shape functions.",
156 ExaDG::Patterns::Enum<ExaDG::IncNS::PressureDegree>(),
157 true);
158 }
159 prm.leave_subsection();
160
161 prm.parse_input(input_file, "", true, true);
162
163 auto const lambda_get_dofs_per_element =
164 [&](unsigned int const dim, unsigned int const degree, ExaDG::ElementType const element_type) {
165 return ExaDG::IncNS::get_dofs_per_element(
166 throughput.operator_type, pressure_degree, dim, degree, element_type);
167 };
168
169 // fill resolution vector depending on the operator_type
170 resolution.fill_resolution_vector(lambda_get_dofs_per_element);
171
172 // loop over resolutions vector and run simulations
173 for(auto iter = resolution.resolutions.begin(); iter != resolution.resolutions.end(); ++iter)
174 {
175 unsigned int const degree = std::get<0>(*iter);
176 unsigned int const refine_space = std::get<1>(*iter);
177 unsigned int const n_cells_1d = std::get<2>(*iter);
178
179 if(general.dim == 2 and general.precision == "float")
180 {
181 ExaDG::run<2, float>(
182 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
183 }
184 else if(general.dim == 2 and general.precision == "double")
185 {
186 ExaDG::run<2, double>(
187 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
188 }
189 else if(general.dim == 3 and general.precision == "float")
190 {
191 ExaDG::run<3, float>(
192 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
193 }
194 else if(general.dim == 3 and general.precision == "double")
195 {
196 ExaDG::run<3, double>(
197 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
198 }
199 else
200 {
201 AssertThrow(false,
202 dealii::ExcMessage("Only dim = 2|3 and precision=float|double implemented."));
203 }
204 }
205
206 if(not(general.is_test))
207 throughput.print_results(mpi_comm);
208
209#ifdef EXADG_WITH_LIKWID
210 LIKWID_MARKER_CLOSE;
211#endif
212
213 return 0;
214}
215
216#endif /* INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_THROUGHPUT_H_ */
Definition driver.cpp:33
Definition general_parameters.h:32
Definition hypercube_resolution_parameters.h:182
Definition throughput_parameters.h:92