deal.II version GIT relicensing-2167-g9622207b8f 2024-11-21 12: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 for (unsigned int i = 0; i < tensor_polys.n(); ++i)
60 renumber_base.push_back(renumber[i]);
61
62 tensor_polys.set_numbering(renumber_base);
63}
64
65
66template <int dim>
67double
69 const Point<dim> &p) const
70{
71 const unsigned int max_indices = tensor_polys.n();
72 Assert(i <= max_indices, ExcInternalError());
73
74 // treat the regular basis functions
75 if (i < max_indices)
76 return tensor_polys.compute_value(i, p);
77 else
78 // this is for the constant function
79 return 1.;
80}
81
82
83
84template <>
85double
87 const Point<0> &) const
88{
90 return 0.;
91}
92
93
94template <int dim>
97 const Point<dim> &p) const
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
110template <int dim>
113 const Point<dim> &p) const
114{
115 const unsigned int max_indices = tensor_polys.n();
116 Assert(i <= max_indices, ExcInternalError());
117
118 // treat the regular basis functions
119 if (i < max_indices)
120 return tensor_polys.compute_grad_grad(i, p);
121 else
122 // this is for the constant function
123 return Tensor<2, dim>();
124}
125
126template <int dim>
127void
129 const Point<dim> &p,
130 std::vector<double> &values,
131 std::vector<Tensor<1, dim>> &grads,
132 std::vector<Tensor<2, dim>> &grad_grads,
133 std::vector<Tensor<3, dim>> &third_derivatives,
134 std::vector<Tensor<4, dim>> &fourth_derivatives) const
135{
136 Assert(values.size() == tensor_polys.n() + 1 || values.empty(),
137 ExcDimensionMismatch2(values.size(), tensor_polys.n() + 1, 0));
138 Assert(grads.size() == tensor_polys.n() + 1 || grads.empty(),
139 ExcDimensionMismatch2(grads.size(), tensor_polys.n() + 1, 0));
140 Assert(grad_grads.size() == tensor_polys.n() + 1 || grad_grads.empty(),
141 ExcDimensionMismatch2(grad_grads.size(), tensor_polys.n() + 1, 0));
142 Assert(third_derivatives.size() == tensor_polys.n() + 1 ||
143 third_derivatives.empty(),
144 ExcDimensionMismatch2(third_derivatives.size(),
145 tensor_polys.n() + 1,
146 0));
147 Assert(fourth_derivatives.size() == tensor_polys.n() + 1 ||
148 fourth_derivatives.empty(),
149 ExcDimensionMismatch2(fourth_derivatives.size(),
150 tensor_polys.n() + 1,
151 0));
152
153 // remove slot for const value, go into the base class compute method and
154 // finally append the const value again
155 bool do_values = false, do_grads = false, do_grad_grads = false;
156 bool do_3rd_derivatives = false, do_4th_derivatives = false;
157 if (values.empty() == false)
158 {
159 values.pop_back();
160 do_values = true;
161 }
162 if (grads.empty() == false)
163 {
164 grads.pop_back();
165 do_grads = true;
166 }
167 if (grad_grads.empty() == false)
168 {
169 grad_grads.pop_back();
170 do_grad_grads = true;
171 }
172 if (third_derivatives.empty() == false)
173 {
174 third_derivatives.resize(tensor_polys.n());
175 do_3rd_derivatives = true;
176 }
177 if (fourth_derivatives.empty() == false)
178 {
179 fourth_derivatives.resize(tensor_polys.n());
180 do_4th_derivatives = true;
181 }
182
183 tensor_polys.evaluate(
184 p, values, grads, grad_grads, third_derivatives, fourth_derivatives);
185
186 // for dgq node: values =1, grads=0, grads_grads=0, third_derivatives=0,
187 // fourth_derivatives=0
188 if (do_values)
189 values.push_back(1.);
190 if (do_grads)
191 grads.emplace_back();
192 if (do_grad_grads)
193 grad_grads.emplace_back();
194 if (do_3rd_derivatives)
195 third_derivatives.emplace_back();
196 if (do_4th_derivatives)
197 fourth_derivatives.emplace_back();
198}
199
200
201
202template <int dim>
203std::unique_ptr<ScalarPolynomialsBase<dim>>
205{
206 return std::make_unique<TensorProductPolynomialsConst<dim>>(*this);
207}
208
209
210/* ------------------- explicit instantiations -------------- */
214
Definition point.h:111
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:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
#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)
#define DEAL_II_NOT_IMPLEMENTED()