Reference documentation for deal.II version GIT 472adc65c1 2022-06-27 15:55:01+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\}}\)
evaluation_selector.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2017 - 2021 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 
25 
27 
28 
29 
35 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
37 {
45  static void
46  evaluate(const unsigned int n_components,
47  const EvaluationFlags::EvaluationFlags evaluation_flag,
48  const Number * values_dofs,
50 
58  static void
59  integrate(const unsigned int n_components,
60  const EvaluationFlags::EvaluationFlags integration_flag,
61  Number * values_dofs,
63  const bool sum_into_values_array = false);
64 };
65 
66 //----------------------Implementation for SelectEvaluator---------------------
67 #ifndef DOXYGEN
68 
69 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
70 inline void
72  const unsigned int n_components,
73  const EvaluationFlags::EvaluationFlags evaluation_flag,
74  const Number * values_dofs,
76 {
77  Assert(fe_degree >= 0 && n_q_points_1d > 0, ExcInternalError());
78 
80  fe_degree,
81  n_q_points_1d>(n_components, evaluation_flag, values_dofs, eval);
82 }
83 
84 
85 
86 template <int dim, int fe_degree, int n_q_points_1d, typename Number>
87 inline void
89  const unsigned int n_components,
90  const EvaluationFlags::EvaluationFlags integration_flag,
91  Number * values_dofs,
93  const bool sum_into_values_array)
94 {
95  Assert(fe_degree >= 0 && n_q_points_1d > 0, ExcInternalError());
96 
98  template run<fe_degree, n_q_points_1d>(
99  n_components, integration_flag, values_dofs, eval, sum_into_values_array);
100 }
101 #endif // DOXYGEN
102 
104 
105 #endif
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
static ::ExceptionBase & ExcInternalError()
#define Assert(cond, exc)
Definition: exceptions.h:1473
EvaluationFlags
The EvaluationFlags enum.
void run(const Iterator &begin, const typename identity< Iterator >::type &end, Worker worker, Copier copier, const ScratchData &sample_scratch_data, const CopyData &sample_copy_data, const unsigned int queue_length, const unsigned int chunk_size)
Definition: work_stream.h:474
static void integrate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags integration_flag, Number *values_dofs, FEEvaluationData< dim, Number, false > &eval, const bool sum_into_values_array=false)
static void evaluate(const unsigned int n_components, const EvaluationFlags::EvaluationFlags evaluation_flag, const Number *values_dofs, FEEvaluationData< dim, Number, false > &eval)