22#ifndef EXADG_POISSON_THROUGHPUT_H_
23#define EXADG_POISSON_THROUGHPUT_H_
26#ifdef EXADG_WITH_LIKWID
31#include <deal.II/base/exceptions.h>
32#include <deal.II/base/parameter_handler.h>
35#include <exadg/operators/finite_element.h>
36#include <exadg/operators/hypercube_resolution_parameters.h>
37#include <exadg/operators/throughput_parameters.h>
38#include <exadg/poisson/driver.h>
39#include <exadg/poisson/user_interface/declare_get_application.h>
40#include <exadg/utilities/enum_patterns.h>
41#include <exadg/utilities/general_parameters.h>
46create_input_file(std::string
const & input_file)
48 dealii::ParameterHandler prm;
51 general.add_parameters(prm);
54 resolution.add_parameters(prm);
57 throughput.add_parameters(prm);
63 unsigned int const Dim = 2;
64 typedef double Number;
65 Poisson::get_application<Dim, 1, Number>(input_file, MPI_COMM_WORLD)->add_parameters(prm);
71 prm.print_parameters(input_file,
72 dealii::ParameterHandler::Short |
73 dealii::ParameterHandler::KeepDeclarationOrder);
76template<
int dim,
typename Number>
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,
86 std::shared_ptr<Poisson::ApplicationBase<dim, 1, Number>> application =
87 Poisson::get_application<dim, 1, Number>(input_file, mpi_comm);
89 application->set_parameters_throughput_study(degree, refine_space, n_cells_1d);
91 std::shared_ptr<Poisson::Driver<dim, Number>> driver =
92 std::make_shared<Poisson::Driver<dim, Number>>(mpi_comm, application, is_test,
true);
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);
101 throughput.wall_times.push_back(wall_time);
106main(
int argc,
char ** argv)
108#ifdef EXADG_WITH_LIKWID
112 dealii::Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
114 MPI_Comm mpi_comm(MPI_COMM_WORLD);
116 std::string input_file;
120 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
122 std::cout <<
"To run the program, use: ./throughput input_file" << std::endl
123 <<
"To setup the input file, use: ./throughput input_file --help" << std::endl;
130 input_file = std::string(argv[1]);
132 if(argc == 3 and std::string(argv[2]) ==
"--help")
134 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
135 ExaDG::create_input_file(input_file);
146 ExaDG::Poisson::SpatialDiscretization spatial_discretization =
147 ExaDG::Poisson::SpatialDiscretization::Undefined;
149 dealii::ParameterHandler prm;
150 prm.enter_subsection(
"Throughput");
152 prm.add_parameter(
"SpatialDiscretization",
153 spatial_discretization,
154 "Spatial discretization (CG vs. DG).",
155 ExaDG::Patterns::Enum<ExaDG::Poisson::SpatialDiscretization>(),
158 prm.leave_subsection();
160 prm.parse_input(input_file,
"",
true,
true);
162 auto const lambda_get_dofs_per_element =
163 [&](
unsigned int const dim,
unsigned int const degree, ExaDG::ElementType
const element_type) {
164 return ExaDG::get_dofs_per_element(element_type,
165 spatial_discretization ==
166 ExaDG::Poisson::SpatialDiscretization::DG,
173 resolution.fill_resolution_vector(lambda_get_dofs_per_element);
176 for(
auto iter = resolution.resolutions.begin(); iter != resolution.resolutions.end(); ++iter)
178 unsigned int const degree = std::get<0>(*iter);
179 unsigned int const refine_space = std::get<1>(*iter);
180 unsigned int const n_cells_1d = std::get<2>(*iter);
182 if(general.dim == 2 and general.precision ==
"float")
184 ExaDG::run<2, float>(
185 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
187 else if(general.dim == 2 and general.precision ==
"double")
189 ExaDG::run<2, double>(
190 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
192 else if(general.dim == 3 and general.precision ==
"float")
194 ExaDG::run<3, float>(
195 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
197 else if(general.dim == 3 and general.precision ==
"double")
199 ExaDG::run<3, double>(
200 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
205 dealii::ExcMessage(
"Only dim = 2|3 and precision=float|double implemented."));
209 if(not(general.is_test))
210 throughput.print_results(mpi_comm);
212#ifdef EXADG_WITH_LIKWID
Definition general_parameters.h:34
Definition hypercube_resolution_parameters.h:183
Definition throughput_parameters.h:91