22#ifndef INCLUDE_EXADG_COMPRESSIBLE_NAVIER_STOKES_SOLVER_H_
23#define INCLUDE_EXADG_COMPRESSIBLE_NAVIER_STOKES_SOLVER_H_
26#include <exadg/compressible_navier_stokes/driver.h>
29#include <exadg/operators/resolution_parameters.h>
30#include <exadg/time_integration/resolution_parameters.h>
31#include <exadg/utilities/general_parameters.h>
34#include <exadg/compressible_navier_stokes/user_interface/declare_get_application.h>
39create_input_file(std::string
const & input_file)
41 dealii::ParameterHandler prm;
44 general.add_parameters(prm);
47 spatial.add_parameters(prm);
50 temporal.add_parameters(prm);
54 unsigned int const Dim = 2;
55 typedef double Number;
56 CompNS::get_application<Dim, Number>(input_file, MPI_COMM_WORLD)->add_parameters(prm);
58 prm.print_parameters(input_file,
59 dealii::ParameterHandler::Short |
60 dealii::ParameterHandler::KeepDeclarationOrder);
63template<
int dim,
typename Number>
65run(std::string
const & input_file,
66 unsigned int const degree,
67 unsigned int const refine_space,
68 unsigned int const refine_time,
69 MPI_Comm
const & mpi_comm,
75 std::shared_ptr<CompNS::ApplicationBase<dim, Number>> application =
76 CompNS::get_application<dim, Number>(input_file, mpi_comm);
78 application->set_parameters_convergence_study(degree, refine_space, refine_time);
80 std::shared_ptr<CompNS::Driver<dim, Number>> driver =
81 std::make_shared<CompNS::Driver<dim, Number>>(mpi_comm, application, is_test,
false);
88 driver->print_performance_results(timer.wall_time());
93main(
int argc,
char ** argv)
95 dealii::Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
97 MPI_Comm mpi_comm(MPI_COMM_WORLD);
99 std::string input_file;
103 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
106 std::cout <<
"To run the program, use: ./solver input_file" << std::endl
107 <<
"To setup the input file, use: ./solver input_file --help" << std::endl;
115 input_file = std::string(argv[1]);
117 if(argc == 3 and std::string(argv[2]) ==
"--help")
119 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
120 ExaDG::create_input_file(input_file);
131 for(
unsigned int degree = spatial.degree_min; degree <= spatial.degree_max; ++degree)
134 for(
unsigned int refine_space = spatial.refine_space_min;
135 refine_space <= spatial.refine_space_max;
139 for(
unsigned int refine_time = temporal.refine_time_min;
140 refine_time <= temporal.refine_time_max;
144 if(general.dim == 2 and general.precision ==
"float")
146 ExaDG::run<2, float>(
147 input_file, degree, refine_space, refine_time, mpi_comm, general.is_test);
149 else if(general.dim == 2 and general.precision ==
"double")
151 ExaDG::run<2, double>(
152 input_file, degree, refine_space, refine_time, mpi_comm, general.is_test);
154 else if(general.dim == 3 and general.precision ==
"float")
156 ExaDG::run<3, float>(
157 input_file, degree, refine_space, refine_time, mpi_comm, general.is_test);
159 else if(general.dim == 3 and general.precision ==
"double")
161 ExaDG::run<3, double>(
162 input_file, degree, refine_space, refine_time, mpi_comm, general.is_test);
167 dealii::ExcMessage(
"Only dim = 2|3 and precision=float|double implemented."));
Definition general_parameters.h:32
Definition resolution_parameters.h:32
Definition resolution_parameters.h:32