22#ifndef INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SOLVER_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_SOLVER_H_
26#include <exadg/incompressible_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/incompressible_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 IncNS::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<IncNS::ApplicationBase<dim, Number>> application =
76 IncNS::get_application<dim, Number>(input_file, mpi_comm);
78 application->set_parameters_convergence_study(degree, refine_space, refine_time);
80 std::shared_ptr<IncNS::Driver<dim, Number>> driver =
81 std::make_shared<IncNS::Driver<dim, Number>>(mpi_comm, application, is_test,
false);
88 driver->print_performance_results(timer.wall_time());
96main(
int argc,
char ** argv)
98 dealii::Utilities::MPI::MPI_InitFinalize mpi(argc, argv, 1);
100 MPI_Comm mpi_comm(MPI_COMM_WORLD);
105#ifdef USE_SUB_COMMUNICATOR
109 int const rank = dealii::Utilities::MPI::this_mpi_process(mpi_comm);
110 int const size = dealii::Utilities::MPI::n_mpi_processes(mpi_comm);
112 int const new_rank = rank + (rank % n) * size;
115 MPI_Comm_split(mpi_comm, flag, new_rank, &sub_comm);
118 std::cout << std::endl <<
"Created sub communicator with stride of " << n << std::endl;
123 std::string input_file;
127 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
130 std::cout <<
"To run the program, use: ./solver input_file" << std::endl
131 <<
"To setup the input file, use: ./solver input_file --help" << std::endl;
139 input_file = std::string(argv[1]);
141 if(argc == 3 and std::string(argv[2]) ==
"--help")
143 if(dealii::Utilities::MPI::this_mpi_process(mpi_comm) == 0)
144 ExaDG::create_input_file(input_file);
155 for(
unsigned int degree = spatial.degree_min; degree <= spatial.degree_max; ++degree)
158 for(
unsigned int refine_space = spatial.refine_space_min;
159 refine_space <= spatial.refine_space_max;
163 for(
unsigned int refine_time = temporal.refine_time_min;
164 refine_time <= temporal.refine_time_max;
168 if(general.dim == 2 and general.precision ==
"float")
170 ExaDG::run<2, float>(
171 input_file, degree, refine_space, refine_time, sub_comm, general.is_test);
173 else if(general.dim == 2 and general.precision ==
"double")
175 ExaDG::run<2, double>(
176 input_file, degree, refine_space, refine_time, sub_comm, general.is_test);
178 else if(general.dim == 3 and general.precision ==
"float")
180 ExaDG::run<3, float>(
181 input_file, degree, refine_space, refine_time, sub_comm, general.is_test);
183 else if(general.dim == 3 and general.precision ==
"double")
185 ExaDG::run<3, double>(
186 input_file, degree, refine_space, refine_time, sub_comm, general.is_test);
191 false, dealii::ExcMessage(
"Only dim = 2|3 and precision = float|double implemented."));
197#ifdef USE_SUB_COMMUNICATOR
199 MPI_Comm_free(&sub_comm);
Definition general_parameters.h:32
Definition resolution_parameters.h:32
Definition resolution_parameters.h:32