]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
access dot shape gradients eliminated
authorguido <guido@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 30 May 2002 18:20:35 +0000 (18:20 +0000)
committerguido <guido@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 30 May 2002 18:20:35 +0000 (18:20 +0000)
git-svn-id: https://svn.dealii.org/trunk@5946 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/deal.II/source/numerics/matrices.cc
deal.II/deal.II/source/numerics/vectors.cc

index 3301ada8514f4c558c38bc889786e621ddfa9ce9..4f57642b4e174ecc7169549b091e0d3aa4dc9427 100644 (file)
@@ -980,9 +980,6 @@ void MatrixCreator::create_laplace_matrix_1 (const Mapping<dim>       &mapping,
       cell_matrix.clear ();
       cell->get_dof_indices (dof_indices);
       
-      const std::vector<std::vector<Tensor<1,dim> > >
-       &grads   = fe_values.get_shape_grads ();
-      
       if (coefficient != 0)
        {
          if (coefficient->n_components==1)
@@ -994,14 +991,14 @@ void MatrixCreator::create_laplace_matrix_1 (const Mapping<dim>       &mapping,
                  const double weight = fe_values.JxW(point);
                  for (unsigned int i=0; i<dofs_per_cell; ++i)
                    {
+                     const Tensor<1,dim>& Dv = fe_values.shape_grad(i,point);
                      for (unsigned int j=0; j<dofs_per_cell; ++j)
                        {
+                         const Tensor<1,dim>& Du = fe_values.shape_grad(j,point);
                          if ((n_components==1) ||
                              (fe.system_to_component_index(i).first ==
                               fe.system_to_component_index(j).first))
-                           cell_matrix(i,j) += (grads[i][point] *
-                                                grads[j][point] *
-                                                weight *
+                           cell_matrix(i,j) += (Du * Dv * weight *
                                                 coefficient_values[point]);
                        }
                    }
@@ -1016,15 +1013,15 @@ void MatrixCreator::create_laplace_matrix_1 (const Mapping<dim>       &mapping,
                  const double weight = fe_values.JxW(point);
                  for (unsigned int i=0; i<dofs_per_cell; ++i)
                    {
+                     const Tensor<1,dim>& Dv = fe_values.shape_grad(i,point);
                      const unsigned int component_i=
                        fe.system_to_component_index(i).first;
                      for (unsigned int j=0; j<dofs_per_cell; ++j)
                        {
+                         const Tensor<1,dim>& Du = fe_values.shape_grad(j,point);
                          if ((n_components==1) ||
                              (fe.system_to_component_index(j).first == component_i))
-                           cell_matrix(i,j) += (grads[i][point] *
-                                                grads[j][point] *
-                                                weight *
+                           cell_matrix(i,j) += (Du * Dv * weight *
                                                 coefficient_vector_values[point](component_i));
                          
                        }
@@ -1038,14 +1035,14 @@ void MatrixCreator::create_laplace_matrix_1 (const Mapping<dim>       &mapping,
            const double weight = fe_values.JxW(point);
            for (unsigned int i=0; i<dofs_per_cell; ++i)
              {
+               const Tensor<1,dim>& Dv = fe_values.shape_grad(i,point);
                for (unsigned int j=0; j<dofs_per_cell; ++j)
                  {
+                   const Tensor<1,dim>& Du = fe_values.shape_grad(j,point);
                    if ((n_components==1) ||
                        (fe.system_to_component_index(i).first ==
                         fe.system_to_component_index(j).first))
-                     cell_matrix(i,j) += (grads[i][point] *
-                                          grads[j][point] *
-                                          weight);
+                     cell_matrix(i,j) += (Du * Dv * weight);
                  }
              }
          }
@@ -1177,8 +1174,6 @@ MatrixCreator::create_laplace_matrix_2 (const Mapping<dim>       &mapping,
       local_rhs.clear ();
       cell->get_dof_indices (dof_indices);
       
-      const std::vector<std::vector<Tensor<1,dim> > >
-       &grads   = fe_values.get_shape_grads ();
       rhs.value_list (fe_values.get_quadrature_points(), rhs_values);
       
       if (coefficient != 0)
@@ -1193,14 +1188,14 @@ MatrixCreator::create_laplace_matrix_2 (const Mapping<dim>       &mapping,
                    for (unsigned int i=0; i<dofs_per_cell; ++i) 
                      {
                        const double v = fe_values.shape_value(i,point);
+                       const Tensor<1,dim>& Dv = fe_values.shape_grad(i,point);
                        for (unsigned int j=0; j<dofs_per_cell; ++j)
                          {
+                           const Tensor<1,dim>& Du = fe_values.shape_grad(j,point);
                            if ((n_components==1) ||
                                (fe.system_to_component_index(i).first ==
                                 fe.system_to_component_index(j).first))
-                             cell_matrix(i,j) += (grads[i][point] *
-                                                  grads[j][point] *
-                                                  weight *
+                             cell_matrix(i,j) += (Du * Dv * weight *
                                                   coefficient_values[point]);
                            local_rhs(i) += v * rhs_values[point] * weight;
                          }
@@ -1217,15 +1212,15 @@ MatrixCreator::create_laplace_matrix_2 (const Mapping<dim>       &mapping,
                    for (unsigned int i=0; i<dofs_per_cell; ++i) 
                      {
                        const double v = fe_values.shape_value(i,point);
+                       const Tensor<1,dim>& Dv = fe_values.shape_grad(i,point);
                        const unsigned int component_i=
                          fe.system_to_component_index(i).first;                    
                        for (unsigned int j=0; j<dofs_per_cell; ++j)
                          {
+                           const Tensor<1,dim>& Du = fe_values.shape_grad(j,point);
                            if ((n_components==1) ||
                                (fe.system_to_component_index(j).first == component_i))
-                             cell_matrix(i,j) += (grads[i][point] *
-                                                  grads[j][point] *
-                                                  weight *
+                             cell_matrix(i,j) += (Du * Dv * weight *
                                                   coefficient_vector_values[point](component_i));
                            local_rhs(i) += v * rhs_values[point] * weight;
                          }
@@ -1240,14 +1235,14 @@ MatrixCreator::create_laplace_matrix_2 (const Mapping<dim>       &mapping,
            for (unsigned int i=0; i<dofs_per_cell; ++i) 
              {
                const double v = fe_values.shape_value(i,point);
+               const Tensor<1,dim>& Dv = fe_values.shape_grad(i,point);
                for (unsigned int j=0; j<dofs_per_cell; ++j)
                  {
+                   const Tensor<1,dim>& Du = fe_values.shape_grad(j,point);
                    if ((n_components==1) ||
                        (fe.system_to_component_index(i).first ==
                         fe.system_to_component_index(j).first))
-                     cell_matrix(i,j) += (grads[i][point] *
-                                          grads[j][point] *
-                                          weight);
+                     cell_matrix(i,j) += (Du * Dv * weight);
                    local_rhs(i) += v * rhs_values[point] * weight;
                  }
              }
index 0d7a2ae07d2b9ddf39bd64bbbab82c8b58d81984..f1086644acc4c9ef7de820cb840b9e3b271c3d81 100644 (file)
@@ -484,7 +484,6 @@ void VectorTools::create_right_hand_side (const Mapping<dim>    &mapping,
        {
          fe_values.reinit(cell);
          
-         const FullMatrix<double>  &values    = fe_values.get_shape_values ();
          const std::vector<double> &weights   = fe_values.get_JxW_values ();
          rhs_function.value_list (fe_values.get_quadrature_points(), rhs_values);
          
@@ -492,7 +491,7 @@ void VectorTools::create_right_hand_side (const Mapping<dim>    &mapping,
          for (unsigned int point=0; point<n_q_points; ++point)
            for (unsigned int i=0; i<dofs_per_cell; ++i) 
              cell_vector(i) += rhs_values[point] *
-                               values(i,point) *
+                               fe_values.shape_value(i,point) *
                                weights[point];
        
          cell->get_dof_indices (dofs);
@@ -510,7 +509,6 @@ void VectorTools::create_right_hand_side (const Mapping<dim>    &mapping,
        {
          fe_values.reinit(cell);
          
-         const FullMatrix<double>  &values    = fe_values.get_shape_values ();
          const std::vector<double> &weights   = fe_values.get_JxW_values ();
          rhs_function.vector_value_list (fe_values.get_quadrature_points(), rhs_values);
          
@@ -522,7 +520,7 @@ void VectorTools::create_right_hand_side (const Mapping<dim>    &mapping,
                  = fe.system_to_component_index(i).first;
                
                cell_vector(i) += rhs_values[point](component) *
-                                 values(i,point) *
+                                 fe_values.shape_value(i,point) *
                                  weights[point];
              }
          
@@ -595,7 +593,6 @@ VectorTools::create_boundary_right_hand_side (const Mapping<dim>      &mapping,
            {
              fe_values.reinit(cell, face);
          
-             const FullMatrix<double>  &values    = fe_values.get_shape_values ();
              const std::vector<double> &weights   = fe_values.get_JxW_values ();
              rhs_function.value_list (fe_values.get_quadrature_points(), rhs_values);
              
@@ -603,7 +600,7 @@ VectorTools::create_boundary_right_hand_side (const Mapping<dim>      &mapping,
              for (unsigned int point=0; point<n_q_points; ++point)
                for (unsigned int i=0; i<dofs_per_cell; ++i) 
                  cell_vector(i) += rhs_values[point] *
-                                   values(i,point) *
+                                   fe_values.shape_value(i,point) *
                                    weights[point];
        
              cell->get_dof_indices (dofs);
@@ -625,7 +622,6 @@ VectorTools::create_boundary_right_hand_side (const Mapping<dim>      &mapping,
            {
              fe_values.reinit(cell, face);
          
-             const FullMatrix<double>  &values    = fe_values.get_shape_values ();
              const std::vector<double> &weights   = fe_values.get_JxW_values ();
              rhs_function.vector_value_list (fe_values.get_quadrature_points(), rhs_values);
          
@@ -637,7 +633,7 @@ VectorTools::create_boundary_right_hand_side (const Mapping<dim>      &mapping,
                      = fe.system_to_component_index(i).first;
                    
                    cell_vector(i) += rhs_values[point](component) *
-                                     values(i,point) *
+                                     fe_values.shape_value(i,point) *
                                      weights[point];
                  }
              

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.