Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
deal.II
integrators
patches.h
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
22
#include <deal.II/base/exceptions.h>
23
#include <deal.II/base/quadrature.h>
24
25
#include <deal.II/fe/fe_values.h>
26
#include <deal.II/fe/mapping.h>
27
28
#include <deal.II/meshworker/dof_info.h>
29
30
DEAL_II_NAMESPACE_OPEN
31
32
namespace
LocalIntegrators
33
{
40
namespace
Patches
41
{
42
template
<
int
dim>
43
inline
void
44
points_and_values(
Table<2, double>
& result,
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();
49
AssertVectorVectorDimension
(input, n_comp, fe.
n_quadrature_points
);
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
64
DEAL_II_NAMESPACE_CLOSE
65
66
#endif
ArrayView
Definition:
array_view.h:77
Table< 2, double >
AssertVectorVectorDimension
#define AssertVectorVectorDimension(VEC, DIM1, DIM2)
Definition:
exceptions.h:1595
AssertDimension
#define AssertDimension(dim1, dim2)
Definition:
exceptions.h:1567
LocalIntegrators
Library of integrals over cells and faces.
Definition:
advection.h:34
FEValuesBase
Definition:
fe.h:36
FEValuesBase::get_fe
const FiniteElement< dim, spacedim > & get_fe() const
FEValuesBase::n_quadrature_points
const unsigned int n_quadrature_points
Definition:
fe_values.h:2099
FEValuesBase::quadrature_point
const Point< spacedim > & quadrature_point(const unsigned int q) const
Generated by
1.8.16