Reference documentation for deal.II version 9.6.0
\(\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
derivative_approximation.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2000 - 2023 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#ifndef dealii_derivative_approximation_h
16#define dealii_derivative_approximation_h
17
18#include <deal.II/base/config.h>
19
22
24
26#include <deal.II/fe/mapping.h>
27
29
30#include <deal.II/lac/vector.h>
31#ifdef _MSC_VER
33#endif
34#include <utility>
35
37
161{
177 template <int dim, class InputVector, int spacedim>
178 void
180 const DoFHandler<dim, spacedim> &dof,
181 const InputVector &solution,
183 const unsigned int component = 0);
184
188 template <int dim, class InputVector, int spacedim>
189 void
191 const InputVector &solution,
193 const unsigned int component = 0);
194
212 template <int dim, class InputVector, int spacedim>
213 void
215 const DoFHandler<dim, spacedim> &dof,
216 const InputVector &solution,
218 const unsigned int component = 0);
219
223 template <int dim, class InputVector, int spacedim>
224 void
226 const InputVector &solution,
228 const unsigned int component = 0);
229
243 template <int dim, int spacedim, class InputVector, int order>
244 void
246 const Mapping<dim, spacedim> &mapping,
247 const DoFHandler<dim, spacedim> &dof,
248 const InputVector &solution,
249#ifndef _MSC_VER
251#else
253 &cell,
254#endif
255 Tensor<order, dim> &derivative,
256 const unsigned int component = 0);
257
261 template <int dim, int spacedim, class InputVector, int order>
262 void
264 const DoFHandler<dim, spacedim> &dof,
265 const InputVector &solution,
266#ifndef _MSC_VER
268#else
270 &cell,
271#endif
272 Tensor<order, dim> &derivative,
273 const unsigned int component = 0);
274
278 template <int dim, int order>
279 double
280 derivative_norm(const Tensor<order, dim> &derivative);
281
286 int,
287 int,
288 << "The output vector needs to have a size equal "
289 "to the number of active cells of your triangulation "
290 "but has length "
291 << arg1 << "There are " << arg2
292 << " active cells in your triangulation.");
297 "While computing a finite difference approximation to "
298 "derivatives, the algorithm encountered a cell on which "
299 "the number of linearly "
300 "independent directions that span the matrix Y (discussed "
301 "in the documentation of the DerivativeApproximation "
302 "class) is not equal to dim. The matrix Y then is "
303 "rank deficient and can not be inverted. A common reason "
304 "why this might be happening is if a cell has neither "
305 "left/right (or up/down, or front/back) neighbors, for "
306 "example because the mesh is too coarse.");
307} // namespace DerivativeApproximation
308
309
310
312
313#endif
Abstract base class for mapping classes.
Definition mapping.h:318
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:503
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:504
static ::ExceptionBase & ExcInsufficientDirections()
static ::ExceptionBase & ExcVectorLengthVsNActiveCells(int arg1, int arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition exceptions.h:539
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:494
typename ActiveSelector::active_cell_iterator active_cell_iterator
void approximate_gradient(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const InputVector &solution, Vector< float > &derivative_norm, const unsigned int component=0)
double derivative_norm(const Tensor< order, dim > &derivative)
void approximate_derivative_tensor(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const InputVector &solution, const typename DoFHandler< dim, spacedim >::active_cell_iterator &cell, Tensor< order, dim > &derivative, const unsigned int component=0)
void approximate_second_derivative(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const InputVector &solution, Vector< float > &derivative_norm, const unsigned int component=0)