[P]arallel [Hi]gh-order [Li]brary for [P]DEs  Latest
Parallel High-Order Library for PDEs through hp-adaptive Discontinuous Galerkin methods
main.cpp
1 #include <deal.II/base/utilities.h>
2 
3 #include <deal.II/base/logstream.h>
4 #include <deal.II/base/parameter_handler.h>
5 
6 #include <fenv.h> // catch nan
7 #include <iostream>
8 #include <fstream>
9 
10 #include "testing/tests.h"
11 #include "flow_solver/flow_solver_factory.h"
12 #include "parameters/all_parameters.h"
13 
14 #include "global_counter.hpp"
15 
16 int main (int argc, char *argv[])
17 {
18 // #if !defined(__APPLE__)
19 // feenableexcept(FE_INVALID | FE_OVERFLOW); // catch nan
20 // #endif
21 
22  n_vmult = 0;
23  dRdW_form = 0;
24  dRdW_mult = 0;
25  dRdX_mult = 0;
26  d2R_mult = 0;
27 
28  dealii::Utilities::MPI::MPI_InitFinalize mpi_initialization(argc, argv, 1);
29  const int n_mpi = dealii::Utilities::MPI::n_mpi_processes(MPI_COMM_WORLD);
30  const int mpi_rank = dealii::Utilities::MPI::this_mpi_process(MPI_COMM_WORLD);
31  if (n_mpi==1 || mpi_rank==0) {
32  dealii::deallog.depth_console(99);
33  } else {
34  dealii::deallog.depth_console(0);
35  }
36 
37  dealii::ConditionalOStream pcout(std::cout, mpi_rank==0);
38  pcout << "Starting program with " << n_mpi << " processors..." << std::endl;
39  if ((PHILIP_DIM==1) && !(n_mpi==1)) {
40  std::cout << "********************************************************" << std::endl;
41  std::cout << "Can't use mpirun -np X, where X>1, for 1D." << std::endl
42  << "Currently using " << n_mpi << " processors." << std::endl
43  << "Aborting..." << std::endl;
44  std::cout << "********************************************************" << std::endl;
45  std::abort();
46  }
47  int run_error = 1;
48  try
49  {
50  // Declare possible inputs
51  dealii::ParameterHandler parameter_handler;
53  PHiLiP::Parameters::parse_command_line (argc, argv, parameter_handler);
54 
55  // Read inputs from parameter file and set those values in AllParameters object
56  PHiLiP::Parameters::AllParameters all_parameters;
57  pcout << "Reading input..." << std::endl;
58  all_parameters.parse_parameters (parameter_handler);
59 
60  AssertDimension(all_parameters.dimension, PHILIP_DIM);
61 
62  const int max_dim = PHILIP_DIM;
63  const int max_nstate = 5;
64 
65  if(all_parameters.run_type == PHiLiP::Parameters::AllParameters::RunType::flow_simulation) {
66  std::unique_ptr<PHiLiP::FlowSolver::FlowSolverBase> flow_solver = PHiLiP::FlowSolver::FlowSolverFactory<max_dim,max_nstate>::create_flow_solver(&all_parameters,parameter_handler);
67  run_error = flow_solver->run();
68  pcout << "Flow simulation complete with run error code: " << run_error << std::endl;
69  }
70  else if(all_parameters.run_type == PHiLiP::Parameters::AllParameters::RunType::integration_test) {
71  std::unique_ptr<PHiLiP::Tests::TestsBase> test = PHiLiP::Tests::TestsFactory<max_dim,max_nstate>::create_test(&all_parameters,parameter_handler);
72  run_error = test->run_test();
73  pcout << "Finished integration test with run error code: " << run_error << std::endl;
74  }
75  }
76  catch (std::exception &exc)
77  {
78  std::cerr << std::endl << std::endl
79  << "----------------------------------------------------"
80  << std::endl
81  << "Exception on processing: " << std::endl
82  << exc.what() << std::endl
83  << "Aborting!" << std::endl
84  << "----------------------------------------------------"
85  << std::endl;
86  return 1;
87  }
88 
89  catch (...)
90  {
91  std::cerr << std::endl
92  << std::endl
93  << "----------------------------------------------------"
94  << std::endl
95  << "Unknown exception!" << std::endl
96  << "Aborting!" << std::endl
97  << "----------------------------------------------------"
98  << std::endl;
99  return 1;
100  }
101  //std::cout << "MPI process " << mpi_rank+1 << " out of " << n_mpi << "reached end of program." << std::endl;
102  pcout << "End of program" << std::endl;
103  return run_error;
104 }
static std::unique_ptr< TestsBase > create_test(const Parameters::AllParameters *const parameters_input, dealii::ParameterHandler &parameter_handler_input)
Recursive factory that will create TestBase<int dim, int nstate>
Definition: tests.cpp:351
static std::unique_ptr< FlowSolverBase > create_flow_solver(const Parameters::AllParameters *const parameters_input, const dealii::ParameterHandler &parameter_handler_input)
Recursive factory that will create FlowSolverBase (i.e. FlowSolver<dim,nstate>)
unsigned int dimension
Number of dimensions. Note that it has to match the executable PHiLiP_xD.
Main parameter class that contains the various other sub-parameter classes.
RunType run_type
Selected RunType from the input file.
void parse_parameters(dealii::ParameterHandler &prm)
Retrieve parameters from dealii::ParameterHandler.
static void declare_parameters(dealii::ParameterHandler &prm)
Declare parameters that can be set as inputs and set up the default options.