]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Avoid the use of determinant()
authorNiklas Wik <niiklaswiik@gmail.com>
Tue, 3 May 2022 09:24:42 +0000 (11:24 +0200)
committerNiklas Wik <niiklaswiik@gmail.com>
Thu, 12 May 2022 08:58:06 +0000 (10:58 +0200)
+ Bug fix for general submit_values

include/deal.II/matrix_free/fe_evaluation.h

index bad7578625cd64f8b282d6edfe472ea1af22a4fc..5b15ee71e2cf6bc9a15157e89344f1a3345e5b00 100644 (file)
@@ -5747,7 +5747,10 @@ FEEvaluationAccess<dim, dim, Number, is_face, VectorizedArrayType>::get_value(
           // Cartesian cell
           const Tensor<2, dim, dealii::VectorizedArray<Number>> jac =
             this->jacobian[1];
-          const VectorizedArrayType inv_det = determinant(this->jacobian[0]);
+          const VectorizedArrayType inv_det =
+            (dim == 2) ? this->jacobian[0][0][0] * this->jacobian[0][1][1] :
+                         this->jacobian[0][0][0] * this->jacobian[0][1][1] *
+                           this->jacobian[0][2][2];
 
           // J * u * det(J^-1)
           for (unsigned int comp = 0; comp < n_components; ++comp)
@@ -5818,7 +5821,10 @@ FEEvaluationAccess<dim, dim, Number, is_face, VectorizedArrayType>::
           const Tensor<2, dim, VectorizedArrayType> &inv_t_jac =
             this->jacobian[0];
           const Tensor<2, dim, VectorizedArrayType> &jac = this->jacobian[1];
-          const VectorizedArrayType inv_det = determinant(inv_t_jac);
+          const VectorizedArrayType                  inv_det =
+            (dim == 2) ? this->jacobian[0][0][0] * this->jacobian[0][1][1] :
+                                          this->jacobian[0][0][0] * this->jacobian[0][1][1] *
+                           this->jacobian[0][2][2];
 
           // J * grad_quad * J^-1 * det(J^-1)
           for (unsigned int d = 0; d < dim; ++d)
@@ -5891,7 +5897,22 @@ FEEvaluationAccess<dim, dim, Number, is_face, VectorizedArrayType>::
   if (this->data->element_type ==
       internal::MatrixFreeFunctions::ElementType::tensor_raviart_thomas)
     {
-      if (this->cell_type <= internal::MatrixFreeFunctions::affine)
+      if (!is_face &&
+          this->cell_type == internal::MatrixFreeFunctions::cartesian)
+        {
+          // Cartesian cell
+          const VectorizedArrayType inv_det =
+            (dim == 2) ? this->jacobian[0][0][0] * this->jacobian[0][1][1] :
+                         this->jacobian[0][0][0] * this->jacobian[0][1][1] *
+                           this->jacobian[0][2][2];
+
+          // div * det(J^-1)
+          divergence = this->gradients_quad[q_point] * inv_det;
+          for (unsigned int d = 1; d < dim; ++d)
+            divergence +=
+              this->gradients_quad[(dim * d + d) * nqp + q_point] * inv_det;
+        }
+      else if (this->cell_type <= internal::MatrixFreeFunctions::affine)
         {
           // Affine cell
           // Derivatives are reordered for faces. Need to take this into account
@@ -6075,7 +6096,7 @@ FEEvaluationAccess<dim, dim, Number, is_face, VectorizedArrayType>::
               this->jacobian[0];
           const Tensor<2, dim, VectorizedArrayType> &jac =
             (this->cell_type > internal::MatrixFreeFunctions::affine) ?
-              invert(inv_t_jac) :
+              transpose(invert(inv_t_jac)) :
               this->jacobian[1];
 
           // Derivatives are reordered for faces. Need to take this into account

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.