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\}}\)
evaluation_selector.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2017 - 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
17#ifndef dealii_matrix_free_evaluation_selector_h
18#define dealii_matrix_free_evaluation_selector_h
19
20#include <deal.II/base/config.h>
21
23
25
26
27
33template <int dim, int fe_degree, int n_q_points_1d, typename Number>
35{
43 static void
44 evaluate(const unsigned int n_components,
45 const EvaluationFlags::EvaluationFlags evaluation_flag,
47 Number *values_dofs_actual,
48 Number *values_quad,
49 Number *gradients_quad,
50 Number *hessians_quad,
51 Number *scratch_data);
52
60 static void
61 integrate(const unsigned int n_components,
62 const EvaluationFlags::EvaluationFlags integration_flag,
64 Number * values_dofs_actual,
65 Number * values_quad,
66 Number * gradients_quad,
67 Number * scratch_data,
68 const bool sum_into_values_array = false);
69};
70
71//----------------------Implementation for SelectEvaluator---------------------
72#ifndef DOXYGEN
73
74template <int dim, int fe_degree, int n_q_points_1d, typename Number>
75inline void
77 const unsigned int n_components,
78 const EvaluationFlags::EvaluationFlags evaluation_flag,
80 Number * values_dofs_actual,
81 Number * values_quad,
82 Number * gradients_quad,
83 Number * hessians_quad,
84 Number * scratch_data)
85{
86 Assert(fe_degree >= 0 && n_q_points_1d > 0, ExcInternalError());
87
89 template run<fe_degree, n_q_points_1d>(n_components,
90 evaluation_flag,
91 shape_info,
92 values_dofs_actual,
93 values_quad,
94 gradients_quad,
95 hessians_quad,
96 scratch_data);
97}
98
99
100
101template <int dim, int fe_degree, int n_q_points_1d, typename Number>
102inline void
104 const unsigned int n_components,
105 const EvaluationFlags::EvaluationFlags integration_flag,
107 Number * values_dofs_actual,
108 Number * values_quad,
109 Number * gradients_quad,
110 Number * scratch_data,
111 const bool sum_into_values_array)
112{
113 Assert(fe_degree >= 0 && n_q_points_1d > 0, ExcInternalError());
114
116 template run<fe_degree, n_q_points_1d>(n_components,
117 integration_flag,
118 shape_info,
119 values_dofs_actual,
120 values_quad,
121 gradients_quad,
122 scratch_data,
123 sum_into_values_array);
124}
125#endif // DOXYGEN
126
128
129#endif
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcInternalError()
EvaluationFlags
The EvaluationFlags enum.
static void evaluate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags evaluation_flag, const internal::MatrixFreeFunctions::ShapeInfo< Number > &shape_info, Number *values_dofs_actual, Number *values_quad, Number *gradients_quad, Number *hessians_quad, Number *scratch_data)
static void integrate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags integration_flag, const internal::MatrixFreeFunctions::ShapeInfo< Number > &shape_info, Number *values_dofs_actual, Number *values_quad, Number *gradients_quad, Number *scratch_data, const bool sum_into_values_array=false)