ExaDG
Loading...
Searching...
No Matches
boundary_conditions.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#ifndef INCLUDE_EXADG_STRUCTURE_SPATIAL_DISCRETIZATION_OPERATORS_BOUNDARY_CONDITIONS_H_
22#define INCLUDE_EXADG_STRUCTURE_SPATIAL_DISCRETIZATION_OPERATORS_BOUNDARY_CONDITIONS_H_
23
24#include <exadg/functions_and_boundary_conditions/evaluate_functions.h>
25#include <exadg/matrix_free/integrators.h>
26#include <exadg/structure/user_interface/boundary_descriptor.h>
27
28namespace ExaDG
29{
30namespace Structure
31{
32/*
33 * This function calculates the Neumann boundary value.
34 */
35template<int dim, typename Number>
36inline DEAL_II_ALWAYS_INLINE //
37 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>>
38 calculate_neumann_value(unsigned int const q,
39 FaceIntegrator<dim, dim, Number> const & integrator,
40 BoundaryType const & boundary_type,
41 dealii::types::boundary_id const boundary_id,
42 std::shared_ptr<BoundaryDescriptor<dim> const> boundary_descriptor,
43 double const & time)
44{
45 dealii::Tensor<1, dim, dealii::VectorizedArray<Number>> traction;
46
47 if(boundary_type == BoundaryType::Neumann)
48 {
49 auto bc = boundary_descriptor->neumann_bc.find(boundary_id)->second;
50 auto q_points = integrator.quadrature_point(q);
51
52 traction = FunctionEvaluator<1, dim, Number>::value(*bc, q_points, time);
53 }
54 else if(boundary_type == BoundaryType::NeumannCached)
55 {
56 auto bc = boundary_descriptor->get_neumann_cached_data();
57
58 traction = FunctionEvaluator<1, dim, Number>::value(*bc,
59 integrator.get_current_cell_index(),
60 q,
61 integrator.get_quadrature_index());
62 }
63 else
64 {
65 // do nothing
66
67 AssertThrow(boundary_type == BoundaryType::Dirichlet or
68 boundary_type == BoundaryType::DirichletCached,
69 dealii::ExcMessage("Boundary type of face is invalid or not implemented."));
70 }
71
72 return traction;
73}
74
75} // namespace Structure
76} // namespace ExaDG
77
78
79#endif /* INCLUDE_EXADG_STRUCTURE_SPATIAL_DISCRETIZATION_OPERATORS_BOUNDARY_CONDITIONS_H_ */
Definition driver.cpp:33