Reference documentation for deal.II version 9.2.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\}}\)
patches.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2019 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 
25 #include <deal.II/fe/fe_values.h>
26 #include <deal.II/fe/mapping.h>
27 
29 
31 
32 namespace LocalIntegrators
33 {
40  namespace Patches
41  {
42  template <int dim>
43  inline void
45  const FEValuesBase<dim> & fe,
46  const ArrayView<const std::vector<double>> &input)
47  {
48  const unsigned int n_comp = fe.get_fe().n_components();
50  AssertDimension(result.n_rows(), fe.n_quadrature_points);
51  AssertDimension(result.n_cols(), n_comp + dim);
52 
53  for (unsigned int k = 0; k < fe.n_quadrature_points; ++k)
54  {
55  for (unsigned int d = 0; d < dim; ++d)
56  result(k, d) = fe.quadrature_point(k)[d];
57  for (unsigned int i = 0; i < n_comp; ++i)
58  result(k, dim + i) = input[i][k];
59  }
60  }
61  } // namespace Patches
62 } // namespace LocalIntegrators
63 
65 
66 #endif
fe_values.h
ArrayView
Definition: array_view.h:77
Physics::Elasticity::Kinematics::d
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
Table< 2, double >
AssertVectorVectorDimension
#define AssertVectorVectorDimension(VEC, DIM1, DIM2)
Definition: exceptions.h:1607
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
AssertDimension
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1579
LocalIntegrators
Library of integrals over cells and faces.
Definition: advection.h:34
dof_info.h
FEValuesBase
Definition: fe.h:36
exceptions.h
mapping.h
FEValuesBase::get_fe
const FiniteElement< dim, spacedim > & get_fe() const
quadrature.h
FEValuesBase::n_quadrature_points
const unsigned int n_quadrature_points
Definition: fe_values.h:2101
config.h
FEValuesBase::quadrature_point
const Point< spacedim > & quadrature_point(const unsigned int q) const
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
LocalIntegrators::Patches::points_and_values
void points_and_values(Table< 2, double > &result, const FEValuesBase< dim > &fe, const ArrayView< const std::vector< double >> &input)
Definition: patches.h:44