22#ifndef EXADG_INCOMPRESSIBLE_NAVIER_STOKES_TIME_INTEGRATION_CREATE_TIME_INTEGRATOR_H_
23#define EXADG_INCOMPRESSIBLE_NAVIER_STOKES_TIME_INTEGRATION_CREATE_TIME_INTEGRATOR_H_
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>
40template<
int dim,
typename Number>
41std::shared_ptr<TimeIntBDF<dim, Number>>
43 std::shared_ptr<HelpersALE<dim, Number>
const> helpers_ale,
46 MPI_Comm
const & mpi_comm,
49 std::shared_ptr<TimeIntBDF<dim, Number>> time_integrator;
51 if(parameters.temporal_discretization == TemporalDiscretization::BDFCoupledSolution)
53 std::shared_ptr<OperatorCoupled<dim, Number>> operator_coupled =
54 std::dynamic_pointer_cast<OperatorCoupled<dim, Number>>(pde_operator);
56 time_integrator = std::make_shared<IncNS::TimeIntBDFCoupled<dim, Number>>(
57 operator_coupled, helpers_ale, postprocessor, parameters, mpi_comm, is_test);
59 else if(parameters.temporal_discretization == TemporalDiscretization::BDFDualSplitting)
61 std::shared_ptr<OperatorDualSplitting<dim, Number>> operator_dual_splitting =
62 std::dynamic_pointer_cast<OperatorDualSplitting<dim, Number>>(pde_operator);
64 time_integrator = std::make_shared<IncNS::TimeIntBDFDualSplitting<dim, Number>>(
65 operator_dual_splitting, helpers_ale, postprocessor, parameters, mpi_comm, is_test);
67 else if(parameters.temporal_discretization == TemporalDiscretization::BDFConsistentSplitting)
69 std::shared_ptr<OperatorConsistentSplitting<dim, Number>> operator_consistent_splitting =
70 std::dynamic_pointer_cast<OperatorConsistentSplitting<dim, Number>>(pde_operator);
72 time_integrator = std::make_shared<IncNS::TimeIntBDFConsistentSplitting<dim, Number>>(
73 operator_consistent_splitting, helpers_ale, postprocessor, parameters, mpi_comm, is_test);
75 else if(parameters.temporal_discretization == TemporalDiscretization::BDFPressureCorrection)
77 std::shared_ptr<OperatorPressureCorrection<dim, Number>> operator_pressure_correction =
78 std::dynamic_pointer_cast<OperatorPressureCorrection<dim, Number>>(pde_operator);
80 time_integrator = std::make_shared<IncNS::TimeIntBDFPressureCorrection<dim, Number>>(
81 operator_pressure_correction, helpers_ale, postprocessor, parameters, mpi_comm, is_test);
83 else if(parameters.temporal_discretization ==
84 TemporalDiscretization::InterpolateAnalyticalSolution)
86 time_integrator = std::make_shared<IncNS::TimeIntInterpolateAnalyticalSolution<dim, Number>>(
87 pde_operator, helpers_ale, postprocessor, parameters, mpi_comm, is_test);
91 AssertThrow(
false, dealii::ExcMessage(
"Not implemented."));
94 return time_integrator;
Definition parameters.h:46
Definition postprocessor_interface.h:37
Definition spatial_operator_base.h:68