Loading [MathJax]/extensions/TeX/newcommand.js
 Reference documentation for deal.II version 9.3.3
\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\}}
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
patches.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2010 - 2020 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_integrators_patches_h
17#define dealii_integrators_patches_h
18
19
20#include <deal.II/base/config.h>
21
24
26#include <deal.II/fe/mapping.h>
27
29
31
32namespace LocalIntegrators
33{
37 namespace Patches
38 {
39 template <int dim>
40 inline void
42 const FEValuesBase<dim> & fe,
43 const ArrayView<const std::vector<double>> &input)
44 {
45 const unsigned int n_comp = fe.get_fe().n_components();
47 AssertDimension(result.n_rows(), fe.n_quadrature_points);
48 AssertDimension(result.n_cols(), n_comp + dim);
49
50 for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
51 {
52 for (unsigned int d = 0; d < dim; ++d)
53 result(k, d) = fe.quadrature_point(k)[d];
54 for (unsigned int i = 0; i < n_comp; ++i)
55 result(k, dim + i) = input[i][k];
56 }
57 }
58 } // namespace Patches
59} // namespace LocalIntegrators
60
62
63#endif
const unsigned int n_quadrature_points
Definition: fe_values.h:2432
const Point< spacedim > & quadrature_point(const unsigned int q) const
const FiniteElement< dim, spacedim > & get_fe() const
unsigned int n_components() const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define AssertVectorVectorDimension(VEC, DIM1, DIM2)
Definition: exceptions.h:1649
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1622
void points_and_values(Table< 2, double > &result, const FEValuesBase< dim > &fe, const ArrayView< const std::vector< double > > &input)
Definition: patches.h:41
Library of integrals over cells and faces.
Definition: advection.h:35
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)