ExaDG
Loading...
Searching...
No Matches
project_velocity.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_SPATIAL_DISCRETIZATION_PROJECT_VELOCITY_H_
23#define INCLUDE_EXADG_CONVECTION_DIFFUSION_SPATIAL_DISCRETIZATION_PROJECT_VELOCITY_H_
24
25#include <exadg/functions_and_boundary_conditions/evaluate_functions.h>
26#include <exadg/matrix_free/integrators.h>
27#include <exadg/operators/inverse_mass_operator.h>
28
29namespace ExaDG
30{
31template<int dim, typename Number>
33{
34private:
35 typedef dealii::LinearAlgebra::distributed::Vector<Number> VectorType;
36
37 typedef CellIntegrator<dim, dim, Number> IntegratorCell;
38
39 typedef std::pair<unsigned int, unsigned int> Range;
40
41public:
42 /*
43 * (v_h, u_h)_Omega^e = (v_h, f)_Omega^e -> M * U = RHS -> U = M^{-1} * RHS
44 */
45 void
46 apply(dealii::MatrixFree<dim, Number> const & matrix_free,
47 InverseMassOperatorData const inverse_mass_operator_data,
48 std::shared_ptr<dealii::Function<dim>> const function,
49 double const & time,
50 VectorType & vector)
51 {
52 this->dof_index = inverse_mass_operator_data.dof_index;
53 this->quad_index = inverse_mass_operator_data.quad_index;
54 this->function = function;
55 this->time = time;
56
57 // calculate RHS
58 VectorType src;
59 matrix_free.cell_loop(&VelocityProjection<dim, Number>::cell_loop, this, vector, src);
60
61 // apply M^{-1}
63 inverse_mass.initialize(matrix_free, inverse_mass_operator_data);
64 inverse_mass.apply(vector, vector);
65 }
66
67private:
68 void
69 cell_loop(dealii::MatrixFree<dim, Number> const & matrix_free,
70 VectorType & dst,
71 VectorType const & src,
72 Range const & cell_range) const
73 {
74 (void)src;
75
76 IntegratorCell integrator(matrix_free, dof_index, quad_index);
77
78 for(unsigned int cell = cell_range.first; cell < cell_range.second; ++cell)
79 {
80 integrator.reinit(cell);
81
82 for(unsigned int q = 0; q < integrator.n_q_points; ++q)
83 {
84 integrator.submit_value(
85 FunctionEvaluator<1, dim, Number>::value(*function, integrator.quadrature_point(q), time),
86 q);
87 }
88
89 integrator.integrate(dealii::EvaluationFlags::values);
90
91 integrator.distribute_local_to_global(dst);
92 }
93 }
94
95 unsigned int dof_index;
96 unsigned int quad_index;
97 std::shared_ptr<dealii::Function<dim>> function;
98 double time;
99};
100
101} // namespace ExaDG
102
103#endif /* INCLUDE_EXADG_CONVECTION_DIFFUSION_SPATIAL_DISCRETIZATION_PROJECT_VELOCITY_H_ */
Definition inverse_mass_operator.h:53
Definition project_velocity.h:33
Definition driver.cpp:33
Definition evaluate_functions.h:40
Definition inverse_mass_operator.h:39