22#ifndef EXADG_INCOMPRESSIBLE_NAVIER_STOKES_THROUGHPUT_H_
23#define EXADG_INCOMPRESSIBLE_NAVIER_STOKES_THROUGHPUT_H_
26#ifdef EXADG_WITH_LIKWID
31#include <exadg/incompressible_navier_stokes/driver.h>
32#include <exadg/incompressible_navier_stokes/user_interface/declare_get_application.h>
33#include <exadg/operators/hypercube_resolution_parameters.h>
34#include <exadg/operators/throughput_parameters.h>
35#include <exadg/utilities/general_parameters.h>
40create_input_file(std::string
const & input_file)
42 dealii::ParameterHandler prm;
45 general.add_parameters(prm);
48 resolution.add_parameters(prm);
51 throughput.add_parameters(prm);
57 unsigned int const Dim = 2;
58 typedef double Number;
59 IncNS::get_application<Dim, Number>(input_file, MPI_COMM_WORLD)->add_parameters(prm);
65 prm.print_parameters(input_file,
66 dealii::ParameterHandler::Short |
67 dealii::ParameterHandler::KeepDeclarationOrder);
70template<
int dim,
typename Number>
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,
80 std::shared_ptr<IncNS::ApplicationBase<dim, Number>> application =
81 IncNS::get_application<dim, Number>(input_file, mpi_comm);
83 application->set_parameters_throughput_study(degree, refine_space, n_cells_1d);
85 std::shared_ptr<IncNS::Driver<dim, Number>> driver =
86 std::make_shared<IncNS::Driver<dim, Number>>(mpi_comm, application, is_test,
true);
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);
95 throughput.wall_times.push_back(wall_time);
100main(
int argc,
char ** argv)
102#ifdef EXADG_WITH_LIKWID
106 dealii::Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
108 MPI_Comm mpi_comm(MPI_COMM_WORLD);
110 std::string input_file;
114 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
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;
126 input_file = std::string(argv[1]);
128 if(argc == 3 and std::string(argv[2]) ==
"--help")
130 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
131 ExaDG::create_input_file(input_file);
141 ExaDG::IncNS::PressureDegree pressure_degree = ExaDG::IncNS::PressureDegree::MixedOrder;
143 dealii::ParameterHandler prm;
145 prm.enter_subsection(
"Throughput");
147 prm.add_parameter(
"PressureDegree",
149 "Degree of pressure shape functions.",
150 ExaDG::Patterns::Enum<ExaDG::IncNS::PressureDegree>(),
153 prm.leave_subsection();
155 prm.parse_input(input_file,
"",
true,
true);
157 auto const lambda_get_dofs_per_element =
158 [&](
unsigned int const dim,
unsigned int const degree, ExaDG::ElementType
const element_type) {
159 return ExaDG::IncNS::get_dofs_per_element(
160 throughput.operator_type, pressure_degree, dim, degree, element_type);
164 resolution.fill_resolution_vector(lambda_get_dofs_per_element);
167 for(
auto iter = resolution.resolutions.begin(); iter != resolution.resolutions.end(); ++iter)
169 unsigned int const degree = std::get<0>(*iter);
170 unsigned int const refine_space = std::get<1>(*iter);
171 unsigned int const n_cells_1d = std::get<2>(*iter);
173 if(general.dim == 2 and general.precision ==
"float")
175 ExaDG::run<2, float>(
176 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
178 else if(general.dim == 2 and general.precision ==
"double")
180 ExaDG::run<2, double>(
181 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
183 else if(general.dim == 3 and general.precision ==
"float")
185 ExaDG::run<3, float>(
186 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
188 else if(general.dim == 3 and general.precision ==
"double")
190 ExaDG::run<3, double>(
191 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
196 dealii::ExcMessage(
"Only dim = 2|3 and precision=float|double implemented."));
200 if(not(general.is_test))
201 throughput.print_results(mpi_comm);
203#ifdef EXADG_WITH_LIKWID
Definition general_parameters.h:34
Definition hypercube_resolution_parameters.h:183
Definition throughput_parameters.h:91