]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add compatibility layer for 1D access 16776/head
authorMartin Kronbichler <martin.kronbichler@rub.de>
Mon, 25 Mar 2024 12:33:21 +0000 (13:33 +0100)
committerMartin Kronbichler <martin.kronbichler@rub.de>
Mon, 25 Mar 2024 12:33:21 +0000 (13:33 +0100)
include/deal.II/matrix_free/fe_evaluation.h

index afff5148bd8c1d1d672a8ab86702ef5edd592747..913a8418f0bc2d4f99f0fcb356089193afbe467d 100644 (file)
@@ -340,6 +340,19 @@ public:
   void
   submit_value(const value_type val_in, const unsigned int q_point);
 
+  /**
+   * In 1D, the value_type and gradient_type can be unintentionally be mixed
+   * up because FEEvaluationBase does not distinguish between scalar accessors
+   * and vector-valued accessors and the respective types, but solely in terms
+   * of the number of components and dimension. Thus, enable the use of
+   * submit_value in that case as well.
+   */
+  template <int dim_ = dim,
+            typename = std::enable_if_t<dim_ == 1 && n_components == dim_>>
+  void
+  submit_value(const Tensor<1, 1, VectorizedArrayType> val_in,
+               const unsigned int                      q_point);
+
   /**
    * Return the gradient of a finite element function at quadrature point
    * number @p q_point after a call to FEEvaluation::evaluate() with
@@ -373,6 +386,19 @@ public:
   void
   submit_gradient(const gradient_type grad_in, const unsigned int q_point);
 
+  /**
+   * In 1D, the value_type and gradient_type can be unintentionally be mixed
+   * up because FEEvaluationBase does not distinguish between scalar accessors
+   * and vector-valued accessors and the respective types, but solely in terms
+   * of the number of components and dimension. Thus, enable the use of
+   * submit_value in that case as well.
+   */
+  template <int dim_ = dim,
+            typename = std::enable_if_t<dim_ == 1 && n_components == dim_>>
+  void
+  submit_gradient(const Tensor<2, 1, VectorizedArrayType> val_in,
+                  const unsigned int                      q_point);
+
   /**
    * Write a contribution that is tested by the gradient to the field
    * containing the values on quadrature points with component @p
@@ -4832,6 +4858,22 @@ FEEvaluationBase<dim, n_components_, Number, is_face, VectorizedArrayType>::
 
 
 
+template <int dim,
+          int n_components_,
+          typename Number,
+          bool is_face,
+          typename VectorizedArrayType>
+template <int, typename>
+inline DEAL_II_ALWAYS_INLINE void
+FEEvaluationBase<dim, n_components_, Number, is_face, VectorizedArrayType>::
+  submit_value(const Tensor<1, 1, VectorizedArrayType> val_in,
+               const unsigned int                      q_point)
+{
+  submit_value(val_in[0], q_point);
+}
+
+
+
 template <int dim,
           int n_components_,
           typename Number,
@@ -5096,6 +5138,22 @@ FEEvaluationBase<dim, n_components_, Number, is_face, VectorizedArrayType>::
 
 
 
+template <int dim,
+          int n_components_,
+          typename Number,
+          bool is_face,
+          typename VectorizedArrayType>
+template <int, typename>
+inline DEAL_II_ALWAYS_INLINE void
+FEEvaluationBase<dim, n_components_, Number, is_face, VectorizedArrayType>::
+  submit_gradient(const Tensor<2, 1, VectorizedArrayType> grad_in,
+                  const unsigned int                      q_point)
+{
+  submit_value(grad_in[0], q_point);
+}
+
+
+
 template <int dim,
           int n_components_,
           typename Number,

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.