22#ifndef EXADG_CONVECTION_DIFFUSION_THROUGHPUT_H_
23#define EXADG_CONVECTION_DIFFUSION_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/convection_diffusion/driver.h>
36#include <exadg/convection_diffusion/user_interface/declare_get_application.h>
37#include <exadg/operators/finite_element.h>
38#include <exadg/operators/hypercube_resolution_parameters.h>
39#include <exadg/operators/throughput_parameters.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 ConvDiff::get_application<Dim, Number>(input_file, MPI_COMM_WORLD)->add_parameters(prm);
71 prm.print_parameters(input_file,
72 dealii::ParameterHandler::Short |
73 dealii::ParameterHandler::KeepDeclarationOrder);
77template<
int dim,
typename Number>
80 std::string
const & input_file,
81 unsigned int const degree,
82 unsigned int const refine_space,
83 unsigned int const n_cells_1d,
84 MPI_Comm
const & mpi_comm,
87 std::shared_ptr<ConvDiff::ApplicationBase<dim, Number>> application =
88 ConvDiff::get_application<dim, Number>(input_file, mpi_comm);
90 application->set_parameters_throughput_study(degree, refine_space, n_cells_1d);
92 std::shared_ptr<ConvDiff::Driver<dim, Number>> driver =
93 std::make_shared<ConvDiff::Driver<dim, Number>>(mpi_comm, application, is_test,
true);
97 std::tuple<unsigned int, dealii::types::global_dof_index, double> wall_time =
98 driver->apply_operator(throughput.operator_type,
99 throughput.n_repetitions_inner,
100 throughput.n_repetitions_outer);
102 throughput.wall_times.push_back(wall_time);
108main(
int argc,
char ** argv)
110#ifdef EXADG_WITH_LIKWID
114 dealii::Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
116 MPI_Comm mpi_comm(MPI_COMM_WORLD);
118 std::string input_file;
122 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
125 std::cout <<
"To run the program, use: ./throughput input_file" << std::endl
126 <<
"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);
149 auto const lambda_get_dofs_per_element =
150 [&](
unsigned int const dim,
unsigned int const degree, ExaDG::ElementType
const element_type) {
151 return ExaDG::get_dofs_per_element(
152 element_type,
true , 1 , degree, dim);
156 resolution.fill_resolution_vector(lambda_get_dofs_per_element);
159 for(
auto iter = resolution.resolutions.begin(); iter != resolution.resolutions.end(); ++iter)
161 unsigned int const degree = std::get<0>(*iter);
162 unsigned int const refine_space = std::get<1>(*iter);
163 unsigned int const n_cells_1d = std::get<2>(*iter);
165 if(general.dim == 2 and general.precision ==
"float")
167 ExaDG::run<2, float>(
168 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
170 else if(general.dim == 2 and general.precision ==
"double")
172 ExaDG::run<2, double>(
173 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
175 else if(general.dim == 3 and general.precision ==
"float")
177 ExaDG::run<3, float>(
178 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
180 else if(general.dim == 3 and general.precision ==
"double")
182 ExaDG::run<3, double>(
183 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
188 dealii::ExcMessage(
"Only dim = 2|3 and precision=float|double implemented."));
192 if(not(general.is_test))
193 throughput.print_results(mpi_comm);
195#ifdef EXADG_WITH_LIKWID
Definition general_parameters.h:34
Definition hypercube_resolution_parameters.h:183
Definition throughput_parameters.h:91