22#ifndef INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_TIME_INTEGRATION_CREATE_TIME_INTEGRATOR_H_
23#define INCLUDE_EXADG_INCOMPRESSIBLE_NAVIER_STOKES_TIME_INTEGRATION_CREATE_TIME_INTEGRATOR_H_
25#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf.h>
26#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf_coupled_solver.h>
27#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf_dual_splitting.h>
28#include <exadg/incompressible_navier_stokes/time_integration/time_int_bdf_pressure_correction.h>
29#include <exadg/incompressible_navier_stokes/time_integration/time_int_interpolate_analytical_solution.h>
38template<
int dim,
typename Number>
39std::shared_ptr<TimeIntBDF<dim, Number>>
41 std::shared_ptr<HelpersALE<dim, Number>
const> helpers_ale,
44 MPI_Comm
const & mpi_comm,
47 std::shared_ptr<TimeIntBDF<dim, Number>> time_integrator;
49 if(parameters.temporal_discretization == TemporalDiscretization::BDFCoupledSolution)
51 std::shared_ptr<OperatorCoupled<dim, Number>> operator_coupled =
52 std::dynamic_pointer_cast<OperatorCoupled<dim, Number>>(pde_operator);
54 time_integrator = std::make_shared<IncNS::TimeIntBDFCoupled<dim, Number>>(
55 operator_coupled, helpers_ale, postprocessor, parameters, mpi_comm, is_test);
57 else if(parameters.temporal_discretization == TemporalDiscretization::BDFDualSplittingScheme)
59 std::shared_ptr<OperatorDualSplitting<dim, Number>> operator_dual_splitting =
60 std::dynamic_pointer_cast<OperatorDualSplitting<dim, Number>>(pde_operator);
62 time_integrator = std::make_shared<IncNS::TimeIntBDFDualSplitting<dim, Number>>(
63 operator_dual_splitting, helpers_ale, postprocessor, parameters, mpi_comm, is_test);
65 else if(parameters.temporal_discretization == TemporalDiscretization::BDFPressureCorrection)
67 std::shared_ptr<OperatorPressureCorrection<dim, Number>> operator_pressure_correction =
68 std::dynamic_pointer_cast<OperatorPressureCorrection<dim, Number>>(pde_operator);
70 time_integrator = std::make_shared<IncNS::TimeIntBDFPressureCorrection<dim, Number>>(
71 operator_pressure_correction, helpers_ale, postprocessor, parameters, mpi_comm, is_test);
73 else if(parameters.temporal_discretization ==
74 TemporalDiscretization::InterpolateAnalyticalSolution)
76 time_integrator = std::make_shared<IncNS::TimeIntInterpolateAnalyticalSolution<dim, Number>>(
77 pde_operator, helpers_ale, postprocessor, parameters, mpi_comm, is_test);
81 AssertThrow(
false, dealii::ExcMessage(
"Not implemented."));
84 return time_integrator;
Definition parameters.h:46
Definition postprocessor_interface.h:34
Definition spatial_operator_base.h:68