ExaDG
Loading...
Searching...
No Matches
create_time_integrator.h
1/* ______________________________________________________________________
2 *
3 * ExaDG - High-Order Discontinuous Galerkin for the Exa-Scale
4 *
5 * Copyright (C) 2021 by the ExaDG authors
6 *
7 * This program is free software: you can redistribute it and/or modify
8 * it under the terms of the GNU General Public License as published by
9 * the Free Software Foundation, either version 3 of the License, or
10 * (at your option) any later version.
11 *
12 * This program is distributed in the hope that it will be useful,
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 * GNU General Public License for more details.
16 *
17 * You should have received a copy of the GNU General Public License
18 * along with this program. If not, see <https://www.gnu.org/licenses/>.
19 * ______________________________________________________________________
20 */
21
22#ifndef EXADG_INCOMPRESSIBLE_NAVIER_STOKES_TIME_INTEGRATION_CREATE_TIME_INTEGRATOR_H_
23#define EXADG_INCOMPRESSIBLE_NAVIER_STOKES_TIME_INTEGRATION_CREATE_TIME_INTEGRATOR_H_
24
25// ExaDG
26#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf.h>
27#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf_consistent_splitting.h>
28#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf_coupled_solver.h>
29#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf_dual_splitting.h>
30#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf_pressure_correction.h>
31#include <exadg/incompressible_navier_stokes/time_integration/time_int_interpolate_analytical_solution.h>
32
33namespace ExaDG
34{
35namespace IncNS
36{
40template<int dim, typename Number>
41std::shared_ptr<TimeIntBDF<dim, Number>>
42create_time_integrator(std::shared_ptr<SpatialOperatorBase<dim, Number>> pde_operator,
43 std::shared_ptr<HelpersALE<dim, Number> const> helpers_ale,
44 std::shared_ptr<PostProcessorInterface<Number>> postprocessor,
45 Parameters const & parameters,
46 MPI_Comm const & mpi_comm,
47 bool const is_test)
48{
49 std::shared_ptr<TimeIntBDF<dim, Number>> time_integrator;
50
51 if(parameters.temporal_discretization == TemporalDiscretization::BDFCoupledSolution)
52 {
53 std::shared_ptr<OperatorCoupled<dim, Number>> operator_coupled =
54 std::dynamic_pointer_cast<OperatorCoupled<dim, Number>>(pde_operator);
55
56 time_integrator = std::make_shared<IncNS::TimeIntBDFCoupled<dim, Number>>(
57 operator_coupled, helpers_ale, postprocessor, parameters, mpi_comm, is_test);
58 }
59 else if(parameters.temporal_discretization == TemporalDiscretization::BDFDualSplitting)
60 {
61 std::shared_ptr<OperatorDualSplitting<dim, Number>> operator_dual_splitting =
62 std::dynamic_pointer_cast<OperatorDualSplitting<dim, Number>>(pde_operator);
63
64 time_integrator = std::make_shared<IncNS::TimeIntBDFDualSplitting<dim, Number>>(
65 operator_dual_splitting, helpers_ale, postprocessor, parameters, mpi_comm, is_test);
66 }
67 else if(parameters.temporal_discretization == TemporalDiscretization::BDFConsistentSplitting)
68 {
69 std::shared_ptr<OperatorConsistentSplitting<dim, Number>> operator_consistent_splitting =
70 std::dynamic_pointer_cast<OperatorConsistentSplitting<dim, Number>>(pde_operator);
71
72 time_integrator = std::make_shared<IncNS::TimeIntBDFConsistentSplitting<dim, Number>>(
73 operator_consistent_splitting, helpers_ale, postprocessor, parameters, mpi_comm, is_test);
74 }
75 else if(parameters.temporal_discretization == TemporalDiscretization::BDFPressureCorrection)
76 {
77 std::shared_ptr<OperatorPressureCorrection<dim, Number>> operator_pressure_correction =
78 std::dynamic_pointer_cast<OperatorPressureCorrection<dim, Number>>(pde_operator);
79
80 time_integrator = std::make_shared<IncNS::TimeIntBDFPressureCorrection<dim, Number>>(
81 operator_pressure_correction, helpers_ale, postprocessor, parameters, mpi_comm, is_test);
82 }
83 else if(parameters.temporal_discretization ==
84 TemporalDiscretization::InterpolateAnalyticalSolution)
85 {
86 time_integrator = std::make_shared<IncNS::TimeIntInterpolateAnalyticalSolution<dim, Number>>(
87 pde_operator, helpers_ale, postprocessor, parameters, mpi_comm, is_test);
88 }
89 else
90 {
91 AssertThrow(false, dealii::ExcMessage("Not implemented."));
92 }
93
94 return time_integrator;
95}
96
97} // namespace IncNS
98} // namespace ExaDG
99
100#endif /* EXADG_INCOMPRESSIBLE_NAVIER_STOKES_TIME_INTEGRATION_CREATE_TIME_INTEGRATOR_H_ */
Definition parameters.h:46
Definition postprocessor_interface.h:37
Definition spatial_operator_base.h:68
Definition driver.cpp:33