]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
First round of eliminating calls to FEValues that expose internal data structures.
authorwolf <wolf@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 30 May 2002 12:32:27 +0000 (12:32 +0000)
committerwolf <wolf@0785d39b-7218-0410-832d-ea1e28bc413d>
Thu, 30 May 2002 12:32:27 +0000 (12:32 +0000)
git-svn-id: https://svn.dealii.org/trunk@5934 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/examples/step-14/step-14.cc
deal.II/examples/step-7/step-7.cc
deal.II/examples/step-8/step-8.cc
deal.II/examples/step-9/step-9.cc

index be09ab869faf0ea1ec0d6dcd76b4443bcfdea261..fc24ffbdf89d6443284ec3a3dede17fa412ddce8 100644 (file)
@@ -699,17 +699,13 @@ namespace LaplaceSolver
        cell_matrix.clear ();
 
        fe_values.reinit (cell);
-       const std::vector<std::vector<Tensor<1,dim> > >
-         & shape_grads  = fe_values.get_shape_grads();
-       const std::vector<double>
-         & JxW_values   = fe_values.get_JxW_values();
 
        for (unsigned int q_point=0; q_point<n_q_points; ++q_point)
          for (unsigned int i=0; i<dofs_per_cell; ++i)
            for (unsigned int j=0; j<dofs_per_cell; ++j)
-             cell_matrix(i,j) += (shape_grads[i][q_point] *
-                                  shape_grads[j][q_point] *
-                                  JxW_values[q_point]);
+             cell_matrix(i,j) += (fe_values.shape_grad(i,q_point) *
+                                  fe_values.shape_grad(j,q_point) *
+                                  fe_values.JxW(q_point));
 
 
        cell->get_dof_indices (local_dof_indices);
@@ -959,20 +955,15 @@ namespace LaplaceSolver
        cell_rhs.clear ();
 
        fe_values.reinit (cell);
-       const FullMatrix<double> 
-         & shape_values = fe_values.get_shape_values();
-       const std::vector<double>
-         & JxW_values   = fe_values.get_JxW_values();
-       const std::vector<Point<dim> >
-         & q_points     = fe_values.get_quadrature_points();
-
-       rhs_function->value_list (q_points, rhs_values);
+
+       rhs_function->value_list (fe_values.get_quadrature_points(),
+                                 rhs_values);
       
        for (unsigned int q_point=0; q_point<n_q_points; ++q_point)
          for (unsigned int i=0; i<dofs_per_cell; ++i)
-           cell_rhs(i) += (shape_values (i,q_point) *
+           cell_rhs(i) += (fe_values.shape_value(i,q_point) *
                            rhs_values[q_point] *
-                           JxW_values[q_point]);
+                           fe_values.JxW(q_point));
 
        cell->get_dof_indices (local_dof_indices);
        for (unsigned int i=0; i<dofs_per_cell; ++i)
index a27f10c6bd59ead3da994e29985817a0db8cb5cb..cf6cc13131bd6a7d45c3be48a18757560c8d5a60 100644 (file)
@@ -856,16 +856,8 @@ void LaplaceProblem<dim>::assemble_system ()
       cell_rhs.clear ();
 
       fe_values.reinit (cell);
-      const FullMatrix<double> 
-       & shape_values = fe_values.get_shape_values();
-      const std::vector<std::vector<Tensor<1,dim> > >
-       & shape_grads  = fe_values.get_shape_grads();
-      const std::vector<double>
-       & JxW_values   = fe_values.get_JxW_values();
-      const std::vector<Point<dim> >
-       & q_points     = fe_values.get_quadrature_points();
-
-      right_hand_side.value_list (q_points, rhs_values);
+
+      right_hand_side.value_list (fe_values.get_quadrature_points(), rhs_values);
       
       for (unsigned int q_point=0; q_point<n_q_points; ++q_point)
        for (unsigned int i=0; i<dofs_per_cell; ++i)
@@ -884,17 +876,17 @@ void LaplaceProblem<dim>::assemble_system ()
                                               // their gradients,
                                               // which is the second
                                               // term below:
-             cell_matrix(i,j) += ((shape_grads[i][q_point] *
-                                   shape_grads[j][q_point] *
-                                   JxW_values[q_point])
+             cell_matrix(i,j) += ((fe_values.shape_grad(i,q_point) *
+                                   fe_values.shape_grad(j,q_point) *
+                                   fe_values.JxW(q_point))
                                   +
-                                  (shape_values[i][q_point] *
-                                   shape_values[j][q_point] *
-                                   JxW_values[q_point]));
+                                  (fe_values.shape_value(i,q_point) *
+                                   fe_values.shape_value(j,q_point) *
+                                   fe_values.JxW(q_point)));
 
-           cell_rhs(i) += (shape_values (i,q_point) *
+           cell_rhs(i) += (fe_values.shape_value(i,q_point) *
                            rhs_values [q_point] *
-                           JxW_values[q_point]);
+                           fe_values.JxW(q_point));
          };
 
                                       // Then there is that second
@@ -945,19 +937,6 @@ void LaplaceProblem<dim>::assemble_system ()
                                             // class:
            fe_face_values.reinit (cell, face);
 
-                                            // Then, for simpler
-                                            // access, we alias the
-                                            // various quantities to
-                                            // local variables:
-           const FullMatrix<double> 
-             & face_shape_values   = fe_face_values.get_shape_values();
-           const std::vector<double>
-             & face_JxW_values     = fe_face_values.get_JxW_values();
-           const std::vector<Point<dim> >
-             & face_q_points       = fe_face_values.get_quadrature_points();
-           const std::vector<Point<dim> >
-             & face_normal_vectors = fe_face_values.get_normal_vectors ();
-
                                             // And we can then
                                             // perform the
                                             // integration by using a
@@ -978,8 +957,8 @@ void LaplaceProblem<dim>::assemble_system ()
                                                 // present quadrature
                                                 // point:
                const double neumann_value
-                 = (exact_solution.gradient (face_q_points[q_point]) *
-                    face_normal_vectors[q_point]);
+                 = (exact_solution.gradient (fe_face_values.quadrature_point(q_point)) *
+                    fe_face_values.normal_vector(q_point));
 
                                                 // Using this, we can
                                                 // compute the
@@ -988,8 +967,8 @@ void LaplaceProblem<dim>::assemble_system ()
                                                 // shape function:
                for (unsigned int i=0; i<dofs_per_cell; ++i)
                  cell_rhs(i) += (neumann_value *
-                                 face_shape_values[i][q_point] *
-                                 face_JxW_values[q_point]);
+                                 fe_face_values.shape_value(i,q_point) *
+                                 fe_face_values.JxW(q_point));
              };
          };
 
index fe24ce0cd1290c99af91c13833e3e7246b215451..1a84c60b9254a37ac8adbef0a8acc612e806e32d 100644 (file)
@@ -521,26 +521,12 @@ void ElasticProblem<dim>::assemble_system ()
       cell_rhs.clear ();
 
       fe_values.reinit (cell);
-
-                                      // As in previous examples, we
-                                      // define some abbreviations
-                                      // for the various data that
-                                      // the ``FEValues'' class
-                                      // offers:
-      const FullMatrix<double> 
-       & shape_values = fe_values.get_shape_values();
-      const std::vector<std::vector<Tensor<1,dim> > >
-       & shape_grads  = fe_values.get_shape_grads();
-      const std::vector<double>
-       & JxW_values   = fe_values.get_JxW_values();
-      const std::vector<Point<dim> >
-       & q_points     = fe_values.get_quadrature_points();
       
                                       // Next we get the values of
                                       // the coefficients at the
                                       // quadrature points:
-      lambda.value_list (q_points, lambda_values);
-      mu.value_list     (q_points, mu_values);
+      lambda.value_list (fe_values.get_quadrature_points(), lambda_values);
+      mu.value_list     (fe_values.get_quadrature_points(), mu_values);
 
                                       // Then assemble the entries of
                                       // the local stiffness matrix
@@ -600,7 +586,7 @@ void ElasticProblem<dim>::assemble_system ()
                                                     // This first term is
                                                     // ((lambda+mu) d_i u_i, d_j v_j).
                                                     // Note that
-                                                    // ``shape_grads[i][q_point]''
+                                                    // ``shape_grad(i,q_point)''
                                                     // returns the
                                                     // gradient of
                                                     // the i-th shape
@@ -621,8 +607,8 @@ void ElasticProblem<dim>::assemble_system ()
                                                     // the appended
                                                     // brackets.
                    (
-                     (shape_grads[i][q_point][component_i] *
-                      shape_grads[j][q_point][component_j] *
+                     (fe_values.shape_grad(i,q_point)[component_i] *
+                      fe_values.shape_grad(j,q_point)[component_j] *
                       (lambda_values[q_point] +
                        mu_values[q_point]))
                      +
@@ -666,13 +652,13 @@ void ElasticProblem<dim>::assemble_system ()
                                                       // away by the
                                                       // compiler).
                      ((component_i == component_j) ?
-                      (shape_grads[i][q_point] *
-                       shape_grads[j][q_point] *
+                      (fe_values.shape_grad(i,q_point) *
+                       fe_values.shape_grad(j,q_point) *
                        mu_values[q_point])  :
                       0)
                    )
                    *
-                   JxW_values[q_point];
+                   fe_values.JxW(q_point);
                };
            };
        };
@@ -683,16 +669,17 @@ void ElasticProblem<dim>::assemble_system ()
                                       // introduction. We will
                                       // therefore not discuss it
                                       // further.
-      right_hand_side.vector_value_list (q_points, rhs_values);
+      right_hand_side.vector_value_list (fe_values.get_quadrature_points(),
+                                        rhs_values);
       for (unsigned int i=0; i<dofs_per_cell; ++i)
        {
          const unsigned int 
            component_i = fe.system_to_component_index(i).first;
          
          for (unsigned int q_point=0; q_point<n_q_points; ++q_point)
-           cell_rhs(i) += shape_values[i][q_point] *
+           cell_rhs(i) += fe_values.shape_value(i,q_point) *
                           rhs_values[q_point](component_i) *
-                          JxW_values[q_point];
+                          fe_values.JxW(q_point);
        };
 
                                       // The transfer from local
index 50ee66aec3c435d3fa1bffc62b4dc034d3d34faa..1bf590056ff3943cd1ee274e49dc2317f62e3e50 100644 (file)
@@ -1009,25 +1009,17 @@ assemble_system_interval (const typename DoFHandler<dim>::active_cell_iterator &
       cell_rhs.clear ();
 
                                       // ... then initialize
-                                      // ``FEValues'' object and
-                                      // define aliases to the data
-                                      // it provides...
+                                      // the ``FEValues'' object...
       fe_values.reinit (cell);
-      const FullMatrix<double> 
-       & shape_values = fe_values.get_shape_values();
-      const std::vector<std::vector<Tensor<1,dim> > >
-       & shape_grads  = fe_values.get_shape_grads();
-      const std::vector<double>
-       & JxW_values   = fe_values.get_JxW_values();
-      const std::vector<Point<dim> >
-       & q_points     = fe_values.get_quadrature_points();
 
                                       // ... obtain the values of
                                       // right hand side and
                                       // advection directions at the
                                       // quadrature points...
-      advection_field.value_list (q_points, advection_directions);
-      right_hand_side.value_list (q_points, rhs_values);
+      advection_field.value_list (fe_values.get_quadrature_points(),
+                                 advection_directions);
+      right_hand_side.value_list (fe_values.get_quadrature_points(),
+                                 rhs_values);
 
                                       // ... set the value of the
                                       // streamline diffusion
@@ -1044,17 +1036,17 @@ assemble_system_interval (const typename DoFHandler<dim>::active_cell_iterator &
          {
            for (unsigned int j=0; j<dofs_per_cell; ++j)
              cell_matrix(i,j) += ((advection_directions[q_point] *
-                                   shape_grads[j][q_point]    *
-                                   (shape_values[i][q_point] +
+                                   fe_values.shape_grad(j,q_point)   *
+                                   (fe_values.shape_value(i,q_point) +
                                     delta *
                                     (advection_directions[q_point] *
-                                     shape_grads[i][q_point]))) *
-                                  JxW_values[q_point]);
+                                     fe_values.shape_grad(i,q_point)))) *
+                                  fe_values.JxW(q_point));
 
-           cell_rhs(i) += ((shape_values[i][q_point] +
+           cell_rhs(i) += ((fe_values.shape_value(i,q_point) +
                             delta *
                             (advection_directions[q_point] *
-                             shape_grads[i][q_point])        ) *
+                             fe_values.shape_grad(i,q_point))        ) *
                            rhs_values[i] *
                            fe_values.JxW (q_point));
          };
@@ -1102,32 +1094,19 @@ assemble_system_interval (const typename DoFHandler<dim>::active_cell_iterator &
                                             // above, we have to
                                             // reinitialize the
                                             // FEFaceValues object
-                                            // for the present face,
-                                            // and we also define the
-                                            // usual aliases to the
-                                            // fields holding values
-                                            // of shape functions,
-                                            // normal vectors, or
-                                            // quadrature points.
+                                            // for the present face:
            fe_face_values.reinit (cell, face);
            
-           const FullMatrix<double> 
-             & face_shape_values = fe_face_values.get_shape_values();
-           const std::vector<double>
-             & face_JxW_values   = fe_face_values.get_JxW_values();
-           const std::vector<Point<dim> >
-             & face_q_points     = fe_face_values.get_quadrature_points();
-           const std::vector<Point<dim> >
-             & normal_vectors    = fe_face_values.get_normal_vectors();
-           
                                             // For the quadrature
                                             // points at hand, we ask
                                             // for the values of the
                                             // inflow function and
                                             // for the direction of
                                             // flow:
-           boundary_values.value_list (face_q_points, face_boundary_values);
-           advection_field.value_list (face_q_points, face_advection_directions);
+           boundary_values.value_list (fe_face_values.get_quadrature_points(),
+                                       face_boundary_values);
+           advection_field.value_list (fe_face_values.get_quadrature_points(),
+                                       face_advection_directions);
            
                                             // Now loop over all
                                             // quadrature points and
@@ -1154,7 +1133,9 @@ assemble_system_interval (const typename DoFHandler<dim>::active_cell_iterator &
                                             // normal vector must be
                                             // negative):
            for (unsigned int q_point=0; q_point<n_face_q_points; ++q_point)
-             if (normal_vectors[q_point] * face_advection_directions[q_point] < 0)
+             if (fe_face_values.normal_vector(q_point) *
+                 face_advection_directions[q_point]
+                 < 0)
                                                 // If the is part of
                                                 // the inflow
                                                 // boundary, then
@@ -1174,16 +1155,16 @@ assemble_system_interval (const typename DoFHandler<dim>::active_cell_iterator &
                  {
                    for (unsigned int j=0; j<dofs_per_cell; ++j)
                      cell_matrix(i,j) -= (face_advection_directions[q_point] *
-                                          normal_vectors[q_point] *
-                                          face_shape_values[i][q_point] *
-                                          face_shape_values[j][q_point] *
-                                          face_JxW_values[q_point]);
+                                          fe_face_values.normal_vector(q_point) *
+                                          fe_face_values.shape_value(i,q_point) *
+                                          fe_face_values.shape_value(j,q_point) *
+                                          fe_face_values.JxW(q_point));
                    
                    cell_rhs(i) -= (face_advection_directions[q_point] *
-                                   normal_vectors[q_point] *
-                                   face_boundary_values[q_point] *
-                                   face_shape_values[i][q_point] *
-                                   face_JxW_values[q_point]);
+                                   fe_face_values.normal_vector(q_point) *
+                                   face_boundary_values[q_point]         *
+                                   fe_face_values.shape_value(i,q_point) *
+                                   fe_face_values.JxW(q_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.