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 INCLUDE_EXADG_CONVECTION_DIFFUSION_TIME_INTEGRATION_CREATE_TIME_INTEGRATOR_H_
23#define INCLUDE_EXADG_CONVECTION_DIFFUSION_TIME_INTEGRATION_CREATE_TIME_INTEGRATOR_H_
24
25#include <exadg/convection_diffusion/time_integration/time_int_bdf.h>
26#include <exadg/convection_diffusion/time_integration/time_int_explicit_runge_kutta.h>
27#include <exadg/convection_diffusion/user_interface/parameters.h>
28
29namespace ExaDG
30{
31namespace ConvDiff
32{
36template<int dim, typename Number>
37std::shared_ptr<TimeIntBase>
38create_time_integrator(std::shared_ptr<Operator<dim, Number>> pde_operator,
39 std::shared_ptr<HelpersALE<dim, Number> const> helpers_ale,
40 std::shared_ptr<PostProcessorInterface<Number>> postprocessor,
41 Parameters const & parameters,
42 MPI_Comm const & mpi_comm,
43 bool const is_test)
44{
45 std::shared_ptr<TimeIntBase> time_integrator;
46
47 if(parameters.temporal_discretization == TemporalDiscretization::ExplRK)
48 {
49 time_integrator = std::make_shared<TimeIntExplRK<Number>>(
50 pde_operator, postprocessor, parameters, mpi_comm, is_test);
51 }
52 else if(parameters.temporal_discretization == TemporalDiscretization::BDF)
53 {
54 time_integrator = std::make_shared<TimeIntBDF<dim, Number>>(
55 pde_operator, helpers_ale, postprocessor, parameters, mpi_comm, is_test);
56 }
57 else
58 {
59 AssertThrow(parameters.temporal_discretization == TemporalDiscretization::ExplRK or
60 parameters.temporal_discretization == TemporalDiscretization::BDF,
61 dealii::ExcMessage("Specified time integration scheme is not implemented!"));
62 }
63
64 return time_integrator;
65}
66
67} // namespace ConvDiff
68} // namespace ExaDG
69
70
71
72#endif /* INCLUDE_EXADG_CONVECTION_DIFFUSION_TIME_INTEGRATION_CREATE_TIME_INTEGRATOR_H_ */
Definition driver.cpp:33