deal.II version GIT relicensing-2890-gf173123fa3 2025-03-22 01:40:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
tensor_product_polynomials_const.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2013 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
18
19#include <memory>
20
22
23
24
25/* ------------------- TensorProductPolynomialsConst -------------- */
26
27
28
29template <int dim>
30void
32{
33 std::array<unsigned int, dim> ix;
34 for (unsigned int i = 0; i < tensor_polys.n(); ++i)
35 {
36 tensor_polys.compute_index(i, ix);
37 out << i << "\t";
38 for (unsigned int d = 0; d < dim; ++d)
39 out << ix[d] << " ";
40 out << std::endl;
41 }
42}
43
44
45
46template <int dim>
47void
49 const std::vector<unsigned int> &renumber)
50{
51 Assert(renumber.size() == index_map.size(),
52 ExcDimensionMismatch(renumber.size(), index_map.size()));
53
54 index_map = renumber;
55 for (unsigned int i = 0; i < index_map.size(); ++i)
56 index_map_inverse[index_map[i]] = i;
57
58 std::vector<unsigned int> renumber_base;
59 renumber_base.reserve(tensor_polys.n());
60 for (unsigned int i = 0; i < tensor_polys.n(); ++i)
61 renumber_base.push_back(renumber[i]);
62
63 tensor_polys.set_numbering(renumber_base);
64}
65
66
67template <int dim>
68double
70 const Point<dim> &p) const
71{
72 const unsigned int max_indices = tensor_polys.n();
73 Assert(i <= max_indices, ExcInternalError());
74
75 // treat the regular basis functions
76 if (i < max_indices)
77 return tensor_polys.compute_value(i, p);
78 else
79 // this is for the constant function
80 return 1.;
81}
82
83
84
85template <int dim>
88 const Point<dim> &p) const
89{
90 if constexpr (dim == 0)
91 {
92 (void)i;
93 (void)p;
95 return {};
96 }
97 else
98 {
99 const unsigned int max_indices = tensor_polys.n();
100 Assert(i <= max_indices, ExcInternalError());
101
102 // treat the regular basis functions
103 if (i < max_indices)
104 return tensor_polys.compute_grad(i, p);
105 else
106 // this is for the constant function
107 return Tensor<1, dim>();
108 }
109}
110
111template <int dim>
114 const Point<dim> &p) const
115{
116 const unsigned int max_indices = tensor_polys.n();
117 Assert(i <= max_indices, ExcInternalError());
118
119 // treat the regular basis functions
120 if (i < max_indices)
121 return tensor_polys.compute_grad_grad(i, p);
122 else
123 // this is for the constant function
124 return Tensor<2, dim>();
125}
126
127template <int dim>
128void
130 const Point<dim> &p,
131 std::vector<double> &values,
132 std::vector<Tensor<1, dim>> &grads,
133 std::vector<Tensor<2, dim>> &grad_grads,
134 std::vector<Tensor<3, dim>> &third_derivatives,
135 std::vector<Tensor<4, dim>> &fourth_derivatives) const
136{
137 Assert(values.size() == tensor_polys.n() + 1 || values.empty(),
138 ExcDimensionMismatch2(values.size(), tensor_polys.n() + 1, 0));
139 Assert(grads.size() == tensor_polys.n() + 1 || grads.empty(),
140 ExcDimensionMismatch2(grads.size(), tensor_polys.n() + 1, 0));
141 Assert(grad_grads.size() == tensor_polys.n() + 1 || grad_grads.empty(),
142 ExcDimensionMismatch2(grad_grads.size(), tensor_polys.n() + 1, 0));
143 Assert(third_derivatives.size() == tensor_polys.n() + 1 ||
144 third_derivatives.empty(),
145 ExcDimensionMismatch2(third_derivatives.size(),
146 tensor_polys.n() + 1,
147 0));
148 Assert(fourth_derivatives.size() == tensor_polys.n() + 1 ||
149 fourth_derivatives.empty(),
150 ExcDimensionMismatch2(fourth_derivatives.size(),
151 tensor_polys.n() + 1,
152 0));
153
154 // remove slot for const value, go into the base class compute method and
155 // finally append the const value again
156 bool do_values = false, do_grads = false, do_grad_grads = false;
157 bool do_3rd_derivatives = false, do_4th_derivatives = false;
158 if (values.empty() == false)
159 {
160 values.pop_back();
161 do_values = true;
162 }
163 if (grads.empty() == false)
164 {
165 grads.pop_back();
166 do_grads = true;
167 }
168 if (grad_grads.empty() == false)
169 {
170 grad_grads.pop_back();
171 do_grad_grads = true;
172 }
173 if (third_derivatives.empty() == false)
174 {
175 third_derivatives.resize(tensor_polys.n());
176 do_3rd_derivatives = true;
177 }
178 if (fourth_derivatives.empty() == false)
179 {
180 fourth_derivatives.resize(tensor_polys.n());
181 do_4th_derivatives = true;
182 }
183
184 tensor_polys.evaluate(
185 p, values, grads, grad_grads, third_derivatives, fourth_derivatives);
186
187 // for dgq node: values =1, grads=0, grads_grads=0, third_derivatives=0,
188 // fourth_derivatives=0
189 if (do_values)
190 values.push_back(1.);
191 if (do_grads)
192 grads.emplace_back();
193 if (do_grad_grads)
194 grad_grads.emplace_back();
195 if (do_3rd_derivatives)
196 third_derivatives.emplace_back();
197 if (do_4th_derivatives)
198 fourth_derivatives.emplace_back();
199}
200
201
202
203template <int dim>
204std::unique_ptr<ScalarPolynomialsBase<dim>>
206{
207 return std::make_unique<TensorProductPolynomialsConst<dim>>(*this);
208}
209
210
211/* ------------------- explicit instantiations -------------- */
215
Definition point.h:113
double compute_value(const unsigned int i, const Point< dim > &p) const override
Tensor< 1, dim > compute_grad(const unsigned int i, const Point< dim > &p) const override
virtual std::unique_ptr< ScalarPolynomialsBase< dim > > clone() const override
void set_numbering(const std::vector< unsigned int > &renumber)
Tensor< 2, dim > compute_grad_grad(const unsigned int i, const Point< dim > &p) const override
void evaluate(const Point< dim > &unit_point, std::vector< double > &values, std::vector< Tensor< 1, dim > > &grads, std::vector< Tensor< 2, dim > > &grad_grads, std::vector< Tensor< 3, dim > > &third_derivatives, std::vector< Tensor< 4, dim > > &fourth_derivatives) const override
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_NOT_IMPLEMENTED()
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch2(std::size_t arg1, std::size_t arg2, std::size_t arg3)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)