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