ExaDG
Loading...
Searching...
No Matches
throughput.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2023 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_ACOUSTIC_CONSERVATION_EQUATIONS_THROUGHPUT_H_
23#define EXADG_ACOUSTIC_CONSERVATION_EQUATIONS_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/acoustic_conservation_equations/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/acoustic_conservation_equations/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<Acoustics::OperatorType> throughput;
57 throughput.add_parameters(prm);
58
59 // we have to assume a default dimension and default Number type
60 // for the automatic generation of a default input file
61 unsigned int const Dim = 2;
62 using Number = double;
63 Acoustics::get_application<Dim, Number>(input_file, MPI_COMM_WORLD)->add_parameters(prm);
64
65 prm.print_parameters(input_file,
66 dealii::ParameterHandler::Short |
67 dealii::ParameterHandler::KeepDeclarationOrder);
68}
69
70template<int dim, typename Number>
71void
72run(ThroughputParameters<Acoustics::OperatorType> const & throughput,
73 std::string const & input_file,
74 unsigned int const degree,
75 unsigned int const refine_space,
76 unsigned int const n_cells_1d,
77 MPI_Comm const & mpi_comm,
78 bool const is_test)
79{
80 std::shared_ptr<Acoustics::ApplicationBase<dim, Number>> application =
81 Acoustics::get_application<dim, Number>(input_file, mpi_comm);
82
83 application->set_parameters_throughput_study(degree, refine_space, n_cells_1d);
84
85 std::shared_ptr<Acoustics::Driver<dim, Number>> driver =
86 std::make_shared<Acoustics::Driver<dim, Number>>(mpi_comm, application, is_test, true);
87
88 driver->setup();
89
90 std::tuple<unsigned int, dealii::types::global_dof_index, double> wall_time =
91 driver->apply_operator(throughput.operator_type,
92 throughput.n_repetitions_inner,
93 throughput.n_repetitions_outer);
94
95 throughput.wall_times.push_back(wall_time);
96}
97} // namespace ExaDG
98
99int
100main(int argc, char ** argv)
101{
102#ifdef EXADG_WITH_LIKWID
103 LIKWID_MARKER_INIT;
104#endif
105
106 dealii::Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
107
108 MPI_Comm mpi_comm(MPI_COMM_WORLD);
109
110 std::string input_file;
111
112 if(argc == 1)
113 {
114 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
115 {
116 // clang-format off
117 std::cout << "To run the program, use: ./throughput input_file" << std::endl
118 << "To setup the input file, use: ./throughput input_file --help" << std::endl;
119 // clang-format on
120 }
121
122 return 0;
123 }
124 else if(argc >= 2)
125 {
126 input_file = std::string(argv[1]);
127
128 if(argc == 3 and std::string(argv[2]) == "--help")
129 {
130 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
131 ExaDG::create_input_file(input_file);
132
133 return 0;
134 }
135 }
136
137 ExaDG::GeneralParameters general(input_file);
138 ExaDG::HypercubeResolutionParameters resolution(input_file, general.dim);
140
141 auto const lambda_get_dofs_per_element =
142 [&](unsigned int const dim, unsigned int const degree, ExaDG::ElementType const element_type) {
143 return ExaDG::Acoustics::get_dofs_per_element(dim, degree, element_type);
144 };
145
146 // fill resolution vector depending on the operator_type
147 resolution.fill_resolution_vector(lambda_get_dofs_per_element);
148
149 // loop over resolutions vector and run simulations
150 for(auto iter = resolution.resolutions.begin(); iter != resolution.resolutions.end(); ++iter)
151 {
152 unsigned int const degree = std::get<0>(*iter);
153 unsigned int const refine_space = std::get<1>(*iter);
154 unsigned int const n_cells_1d = std::get<2>(*iter);
155
156 if(general.dim == 2 and general.precision == "float")
157 {
158 ExaDG::run<2, float>(
159 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
160 }
161 else if(general.dim == 2 and general.precision == "double")
162 {
163 ExaDG::run<2, double>(
164 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
165 }
166 else if(general.dim == 3 and general.precision == "float")
167 {
168 ExaDG::run<3, float>(
169 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
170 }
171 else if(general.dim == 3 and general.precision == "double")
172 {
173 ExaDG::run<3, double>(
174 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
175 }
176 else
177 {
178 AssertThrow(false,
179 dealii::ExcMessage("Only dim = 2|3 and precision=float|double implemented."));
180 }
181 }
182
183 if(not(general.is_test))
184 throughput.print_results(mpi_comm);
185
186#ifdef EXADG_WITH_LIKWID
187 LIKWID_MARKER_CLOSE;
188#endif
189
190 return 0;
191}
192
193#endif /* EXADG_ACOUSTIC_CONSERVATION_EQUATIONS_THROUGHPUT_H_ */
Definition driver.cpp:33
Definition general_parameters.h:32
Definition hypercube_resolution_parameters.h:182
Definition throughput_parameters.h:92