]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Move jacobian_gradients into FEEvalData 13253/head
authorPeter Munch <peterrmuench@gmail.com>
Mon, 17 Jan 2022 20:27:18 +0000 (21:27 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Mon, 17 Jan 2022 22:11:28 +0000 (23:11 +0100)
include/deal.II/matrix_free/fe_evaluation.h
include/deal.II/matrix_free/fe_evaluation_data.h

index 86e54949e2e91ce574055b40c1162c3a7a95a86e..7b1d2a1482f6d6f1c6acc92928ad720dba99563e 100644 (file)
@@ -4811,10 +4811,7 @@ FEEvaluationBase<dim, n_components_, Number, is_face, VectorizedArrayType>::
   // cell with general Jacobian
   else
     {
-      const auto &jac_grad =
-        this->mapping_data->jacobian_gradients
-          [1 - this->is_interior_face]
-          [this->mapping_data->data_index_offsets[this->cell] + q_point];
+      const auto &jac_grad = this->jacobian_gradients[q_point];
       for (unsigned int comp = 0; comp < n_components; ++comp)
         {
           VectorizedArrayType tmp[dim][dim];
@@ -4915,10 +4912,7 @@ FEEvaluationBase<dim, n_components_, Number, is_face, VectorizedArrayType>::
   // cell with general Jacobian
   else
     {
-      const Tensor<1, dim *(dim + 1) / 2, Tensor<1, dim, VectorizedArrayType>>
-        &jac_grad =
-          this->mapping_data->jacobian_gradients
-            [0][this->mapping_data->data_index_offsets[this->cell] + q_point];
+      const auto &jac_grad = this->jacobian_gradients[q_point];
       for (unsigned int comp = 0; comp < n_components; ++comp)
         {
           // compute laplacian before the gradient because it needs to access
@@ -5244,10 +5238,7 @@ FEEvaluationBase<dim, n_components_, Number, is_face, VectorizedArrayType>::
     {
       const Tensor<2, dim, VectorizedArrayType> jac = this->jacobian[q_point];
       const VectorizedArrayType                 JxW = this->J_value[q_point];
-      const auto &                              jac_grad =
-        this->mapping_data->jacobian_gradients
-          [1 - this->is_interior_face]
-          [this->mapping_data->data_index_offsets[this->cell] + q_point];
+      const auto &jac_grad = this->jacobian_gradients[q_point];
       for (unsigned int comp = 0; comp < n_components; ++comp)
         {
           // 1. tmp = hessian_in(u) * J
@@ -6902,6 +6893,8 @@ FEEvaluation<dim,
     this->mapping_data->data_index_offsets[cell_index];
   this->jacobian = &this->mapping_data->jacobians[0][offsets];
   this->J_value  = &this->mapping_data->JxW_values[offsets];
+  this->jacobian_gradients =
+    this->mapping_data->jacobian_gradients[0].data() + offsets;
 
   for (unsigned int i = 0; i < VectorizedArrayType::size(); ++i)
     this->cell_ids[i] = cell_index * VectorizedArrayType::size() + i;
@@ -7638,6 +7631,9 @@ FEFaceEvaluation<dim,
   this->normal_x_jacobian =
     &this->mapping_data
        ->normals_times_jacobians[!this->is_interior_face][offsets];
+  this->jacobian_gradients =
+    this->mapping_data->jacobian_gradients[!this->is_interior_face].data() +
+    offsets;
 
 #  ifdef DEBUG
   this->dof_values_initialized     = false;
@@ -7774,6 +7770,9 @@ FEFaceEvaluation<dim,
     &this->matrix_free->get_mapping_info()
        .face_data_by_cells[this->quad_no]
        .normals_times_jacobians[!this->is_interior_face][offsets];
+  this->jacobian_gradients =
+    this->mapping_data->jacobian_gradients[!this->is_interior_face].data() +
+    offsets;
 
 #  ifdef DEBUG
   this->dof_values_initialized     = false;
index de3cb3753535a152fa04ddc3e7f504d6bd589c5e..6d6f90cc5593384dc37c83639f6795240ec4c208 100644 (file)
@@ -614,6 +614,13 @@ protected:
    */
   const Tensor<2, dim, Number> *jacobian;
 
+  /**
+   * A pointer to the gradients of the inverse Jacobian transformation of the
+   * present cell.
+   */
+  const Tensor<1, dim *(dim + 1) / 2, Tensor<1, dim, Number>>
+    *jacobian_gradients;
+
   /**
    * A pointer to the Jacobian determinant of the present cell. If on a
    * Cartesian cell or on a cell with constant Jacobian, this is just the

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.