22#ifndef INCLUDE_EXADG_POISSON_THROUGHPUT_H_
23#define INCLUDE_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/poisson/driver.h>
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>
45#include <exadg/poisson/user_interface/declare_get_application.h>
50create_input_file(std::string
const & input_file)
52 dealii::ParameterHandler prm;
55 general.add_parameters(prm);
58 resolution.add_parameters(prm);
61 throughput.add_parameters(prm);
67 unsigned int const Dim = 2;
68 typedef double Number;
69 Poisson::get_application<Dim, 1, Number>(input_file, MPI_COMM_WORLD)->add_parameters(prm);
75 prm.print_parameters(input_file,
76 dealii::ParameterHandler::Short |
77 dealii::ParameterHandler::KeepDeclarationOrder);
80template<
int dim,
typename Number>
83 std::string
const & input_file,
84 unsigned int const degree,
85 unsigned int const refine_space,
86 unsigned int const n_cells_1d,
87 MPI_Comm
const & mpi_comm,
90 std::shared_ptr<Poisson::ApplicationBase<dim, 1, Number>> application =
91 Poisson::get_application<dim, 1, Number>(input_file, mpi_comm);
93 application->set_parameters_throughput_study(degree, refine_space, n_cells_1d);
95 std::shared_ptr<Poisson::Driver<dim, Number>> driver =
96 std::make_shared<Poisson::Driver<dim, Number>>(mpi_comm, application, is_test,
true);
100 std::tuple<unsigned int, dealii::types::global_dof_index, double> wall_time =
101 driver->apply_operator(throughput.operator_type,
102 throughput.n_repetitions_inner,
103 throughput.n_repetitions_outer);
105 throughput.wall_times.push_back(wall_time);
110main(
int argc,
char ** argv)
112#ifdef EXADG_WITH_LIKWID
116 dealii::Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
118 MPI_Comm mpi_comm(MPI_COMM_WORLD);
120 std::string input_file;
124 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
126 std::cout <<
"To run the program, use: ./throughput input_file" << std::endl
127 <<
"To setup the input file, use: ./throughput input_file --help" << std::endl;
134 input_file = std::string(argv[1]);
136 if(argc == 3 and std::string(argv[2]) ==
"--help")
138 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
139 ExaDG::create_input_file(input_file);
150 ExaDG::Poisson::SpatialDiscretization spatial_discretization =
151 ExaDG::Poisson::SpatialDiscretization::Undefined;
153 dealii::ParameterHandler prm;
154 prm.enter_subsection(
"Throughput");
156 prm.add_parameter(
"SpatialDiscretization",
157 spatial_discretization,
158 "Spatial discretization (CG vs. DG).",
159 ExaDG::Patterns::Enum<ExaDG::Poisson::SpatialDiscretization>(),
162 prm.leave_subsection();
164 prm.parse_input(input_file,
"",
true,
true);
166 auto const lambda_get_dofs_per_element =
167 [&](
unsigned int const dim,
unsigned int const degree, ExaDG::ElementType
const element_type) {
168 return ExaDG::get_dofs_per_element(element_type,
169 spatial_discretization ==
170 ExaDG::Poisson::SpatialDiscretization::DG,
177 resolution.fill_resolution_vector(lambda_get_dofs_per_element);
180 for(
auto iter = resolution.resolutions.begin(); iter != resolution.resolutions.end(); ++iter)
182 unsigned int const degree = std::get<0>(*iter);
183 unsigned int const refine_space = std::get<1>(*iter);
184 unsigned int const n_cells_1d = std::get<2>(*iter);
186 if(general.dim == 2 and general.precision ==
"float")
188 ExaDG::run<2, float>(
189 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
191 else if(general.dim == 2 and general.precision ==
"double")
193 ExaDG::run<2, double>(
194 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
196 else if(general.dim == 3 and general.precision ==
"float")
198 ExaDG::run<3, float>(
199 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
201 else if(general.dim == 3 and general.precision ==
"double")
203 ExaDG::run<3, double>(
204 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
209 dealii::ExcMessage(
"Only dim = 2|3 and precision=float|double implemented."));
213 if(not(general.is_test))
214 throughput.print_results(mpi_comm);
216#ifdef EXADG_WITH_LIKWID
Definition general_parameters.h:32
Definition hypercube_resolution_parameters.h:183
Definition throughput_parameters.h:92