Reference documentation for deal.II version GIT relicensing-233-g802318d791 2024-03-28 20:20:02+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
patches.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2011 - 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_integrators_patches_h
16#define dealii_integrators_patches_h
17
18
19#include <deal.II/base/config.h>
20
23
25#include <deal.II/fe/mapping.h>
26
28
30
31namespace LocalIntegrators
32{
36 namespace Patches
37 {
38 template <int dim>
39 inline void
41 const FEValuesBase<dim> &fe,
42 const ArrayView<const std::vector<double>> &input)
43 {
44 const unsigned int n_comp = fe.get_fe().n_components();
46 AssertDimension(result.n_rows(), fe.n_quadrature_points);
47 AssertDimension(result.n_cols(), n_comp + dim);
48
49 for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
50 {
51 for (unsigned int d = 0; d < dim; ++d)
52 result(k, d) = fe.quadrature_point(k)[d];
53 for (unsigned int i = 0; i < n_comp; ++i)
54 result(k, dim + i) = input[i][k];
55 }
56 }
57 } // namespace Patches
58} // namespace LocalIntegrators
59
61
62#endif
const Point< spacedim > & quadrature_point(const unsigned int q_point) const
const unsigned int n_quadrature_points
const FiniteElement< dim, spacedim > & get_fe() const
unsigned int n_components() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
#define AssertVectorVectorDimension(VEC, DIM1, DIM2)
#define AssertDimension(dim1, dim2)
void points_and_values(Table< 2, double > &result, const FEValuesBase< dim > &fe, const ArrayView< const std::vector< double > > &input)
Definition patches.h:40
Library of integrals over cells and faces.
Definition advection.h:34