22#ifndef INCLUDE_EXADG_COMPRESSIBLE_NAVIER_STOKES_THROUGHPUT_H_
23#define INCLUDE_EXADG_COMPRESSIBLE_NAVIER_STOKES_THROUGHPUT_H_
25#ifdef EXADG_WITH_LIKWID
30#include <exadg/compressible_navier_stokes/driver.h>
33#include <exadg/operators/finite_element.h>
34#include <exadg/operators/hypercube_resolution_parameters.h>
35#include <exadg/operators/throughput_parameters.h>
36#include <exadg/utilities/general_parameters.h>
39#include <exadg/compressible_navier_stokes/user_interface/declare_get_application.h>
44create_input_file(std::string
const & input_file)
46 dealii::ParameterHandler prm;
49 general.add_parameters(prm);
52 resolution.add_parameters(prm);
55 throughput.add_parameters(prm);
61 unsigned int const Dim = 2;
62 typedef double Number;
63 CompNS::get_application<Dim, Number>(input_file, MPI_COMM_WORLD)->add_parameters(prm);
69 prm.print_parameters(input_file,
70 dealii::ParameterHandler::Short |
71 dealii::ParameterHandler::KeepDeclarationOrder);
75template<
int dim,
typename Number>
78 std::string
const & input_file,
79 unsigned int const degree,
80 unsigned int const refine_space,
81 unsigned int const n_cells_1d,
82 MPI_Comm
const & mpi_comm,
85 std::shared_ptr<CompNS::ApplicationBase<dim, Number>> application =
86 CompNS::get_application<dim, Number>(input_file, mpi_comm);
88 application->set_parameters_throughput_study(degree, refine_space, n_cells_1d);
90 std::shared_ptr<CompNS::Driver<dim, Number>> driver =
91 std::make_shared<CompNS::Driver<dim, Number>>(mpi_comm, application, is_test,
true);
95 std::tuple<unsigned int, dealii::types::global_dof_index, double> wall_time =
96 driver->apply_operator(throughput.operator_type,
97 throughput.n_repetitions_inner,
98 throughput.n_repetitions_outer);
100 throughput.wall_times.push_back(wall_time);
105main(
int argc,
char ** argv)
107#ifdef EXADG_WITH_LIKWID
111 dealii::Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
113 MPI_Comm mpi_comm(MPI_COMM_WORLD);
115 std::string input_file;
119 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;
131 input_file = std::string(argv[1]);
133 if(argc == 3 and std::string(argv[2]) ==
"--help")
135 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
136 ExaDG::create_input_file(input_file);
146 auto const lambda_get_dofs_per_element =
147 [&](
unsigned int const dim,
unsigned int const degree, ExaDG::ElementType
const element_type) {
148 return ExaDG::get_dofs_per_element(
149 element_type,
true , dim + 2 , degree, dim);
153 resolution.fill_resolution_vector(lambda_get_dofs_per_element);
156 for(
auto iter = resolution.resolutions.begin(); iter != resolution.resolutions.end(); ++iter)
158 unsigned int const degree = std::get<0>(*iter);
159 unsigned int const refine_space = std::get<1>(*iter);
160 unsigned int const n_cells_1d = std::get<2>(*iter);
162 if(general.dim == 2 and general.precision ==
"float")
164 ExaDG::run<2, float>(
165 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
167 else if(general.dim == 2 and general.precision ==
"double")
169 ExaDG::run<2, double>(
170 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
172 else if(general.dim == 3 and general.precision ==
"float")
174 ExaDG::run<3, float>(
175 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
177 else if(general.dim == 3 and general.precision ==
"double")
179 ExaDG::run<3, double>(
180 throughput, input_file, degree, refine_space, n_cells_1d, mpi_comm, general.is_test);
185 dealii::ExcMessage(
"Only dim = 2|3 and precision=float|double implemented."));
189 if(not(general.is_test))
190 throughput.print_results(mpi_comm);
192#ifdef EXADG_WITH_LIKWID
Definition general_parameters.h:32
Definition hypercube_resolution_parameters.h:183
Definition throughput_parameters.h:92