]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix some null reference warnings. 3715/head
authorDavid Wells <wellsd2@rpi.edu>
Wed, 28 Dec 2016 00:37:22 +0000 (19:37 -0500)
committerDavid Wells <wellsd2@rpi.edu>
Wed, 28 Dec 2016 01:11:57 +0000 (20:11 -0500)
An equivalent way to do this is to use a SmartPointer instead of a
reference: this lets us initialize with something that is not undefined
behavior.

include/deal.II/fe/fe_values.h
source/fe/fe_values.cc

index 59a57a4405d09fc4b44e17a92889de4c67388b27..94f1c6cb9d7cce5e9570b78030c74eb84a2a8f8b 100644 (file)
@@ -380,9 +380,9 @@ namespace FEValuesViews
 
   private:
     /**
-     * A reference to the FEValuesBase object we operate on.
+     * A pointer to the FEValuesBase object we operate on.
      */
-    const FEValuesBase<dim,spacedim> &fe_values;
+    const SmartPointer<const FEValuesBase<dim,spacedim> > fe_values;
 
     /**
      * The single scalar component this view represents of the FEValuesBase
@@ -854,9 +854,9 @@ namespace FEValuesViews
 
   private:
     /**
-     * A reference to the FEValuesBase object we operate on.
+     * A pointer to the FEValuesBase object we operate on.
      */
-    const FEValuesBase<dim,spacedim> &fe_values;
+    const SmartPointer<const FEValuesBase<dim,spacedim> > fe_values;
 
     /**
      * The first component of the vector this view represents of the
@@ -1070,9 +1070,9 @@ namespace FEValuesViews
 
   private:
     /**
-     * A reference to the FEValuesBase object we operate on.
+     * A pointer to the FEValuesBase object we operate on.
      */
-    const FEValuesBase<dim, spacedim> &fe_values;
+    const SmartPointer<const FEValuesBase<dim, spacedim> > fe_values;
 
     /**
      * The first component of the vector this view represents of the
@@ -1275,9 +1275,9 @@ namespace FEValuesViews
 
   private:
     /**
-     * A reference to the FEValuesBase object we operate on.
+     * A pointer to the FEValuesBase object we operate on.
      */
-    const FEValuesBase<dim, spacedim> &fe_values;
+    const SmartPointer<const FEValuesBase<dim, spacedim> > fe_values;
 
     /**
      * The first component of the vector this view represents of the
@@ -3204,18 +3204,18 @@ namespace FEValuesViews
                                const unsigned int q_point) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (shape_function < fe_values.fe->dofs_per_cell,
-            ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
-    Assert (fe_values.update_flags & update_values,
+    Assert (shape_function < fe_values->fe->dofs_per_cell,
+            ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
+    Assert (fe_values->update_flags & update_values,
             typename FVB::ExcAccessToUninitializedField("update_values"));
 
     // an adaptation of the FEValuesBase::shape_value_component function
     // except that here we know the component as fixed and we have
     // pre-computed and cached a bunch of information. See the comments there.
     if (shape_function_data[shape_function].is_nonzero_shape_function_component)
-      return fe_values.finite_element_output.shape_values(shape_function_data[shape_function]
-                                                          .row_index,
-                                                          q_point);
+      return fe_values->finite_element_output.shape_values(shape_function_data[shape_function]
+                                                           .row_index,
+                                                           q_point);
     else
       return 0;
   }
@@ -3230,9 +3230,9 @@ namespace FEValuesViews
                                   const unsigned int q_point) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (shape_function < fe_values.fe->dofs_per_cell,
-            ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
-    Assert (fe_values.update_flags & update_gradients,
+    Assert (shape_function < fe_values->fe->dofs_per_cell,
+            ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
+    Assert (fe_values->update_flags & update_gradients,
             typename FVB::ExcAccessToUninitializedField("update_gradients"));
 
     // an adaptation of the
@@ -3242,8 +3242,8 @@ namespace FEValuesViews
     // pre-computed and cached a bunch of
     // information. See the comments there.
     if (shape_function_data[shape_function].is_nonzero_shape_function_component)
-      return fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function]
-                                                             .row_index][q_point];
+      return fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function]
+                                                              .row_index][q_point];
     else
       return gradient_type();
   }
@@ -3257,9 +3257,9 @@ namespace FEValuesViews
                                  const unsigned int q_point) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (shape_function < fe_values.fe->dofs_per_cell,
-            ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
-    Assert (fe_values.update_flags & update_hessians,
+    Assert (shape_function < fe_values->fe->dofs_per_cell,
+            ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
+    Assert (fe_values->update_flags & update_hessians,
             typename FVB::ExcAccessToUninitializedField("update_hessians"));
 
     // an adaptation of the
@@ -3269,7 +3269,7 @@ namespace FEValuesViews
     // pre-computed and cached a bunch of
     // information. See the comments there.
     if (shape_function_data[shape_function].is_nonzero_shape_function_component)
-      return fe_values.finite_element_output.shape_hessians[shape_function_data[shape_function].row_index][q_point];
+      return fe_values->finite_element_output.shape_hessians[shape_function_data[shape_function].row_index][q_point];
     else
       return hessian_type();
   }
@@ -3283,9 +3283,9 @@ namespace FEValuesViews
                                           const unsigned int q_point) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (shape_function < fe_values.fe->dofs_per_cell,
-            ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
-    Assert (fe_values.update_flags & update_3rd_derivatives,
+    Assert (shape_function < fe_values->fe->dofs_per_cell,
+            ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
+    Assert (fe_values->update_flags & update_3rd_derivatives,
             typename FVB::ExcAccessToUninitializedField("update_3rd_derivatives"));
 
     // an adaptation of the
@@ -3295,7 +3295,7 @@ namespace FEValuesViews
     // pre-computed and cached a bunch of
     // information. See the comments there.
     if (shape_function_data[shape_function].is_nonzero_shape_function_component)
-      return fe_values.finite_element_output.shape_3rd_derivatives[shape_function_data[shape_function].row_index][q_point];
+      return fe_values->finite_element_output.shape_3rd_derivatives[shape_function_data[shape_function].row_index][q_point];
     else
       return third_derivative_type();
   }
@@ -3309,9 +3309,9 @@ namespace FEValuesViews
                                const unsigned int q_point) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (shape_function < fe_values.fe->dofs_per_cell,
-            ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
-    Assert (fe_values.update_flags & update_values,
+    Assert (shape_function < fe_values->fe->dofs_per_cell,
+            ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
+    Assert (fe_values->update_flags & update_values,
             typename FVB::ExcAccessToUninitializedField("update_values"));
 
     // same as for the scalar case except
@@ -3323,7 +3323,7 @@ namespace FEValuesViews
       {
         value_type return_value;
         return_value[shape_function_data[shape_function].single_nonzero_component_index]
-          = fe_values.finite_element_output.shape_values(snc,q_point);
+          = fe_values->finite_element_output.shape_values(snc,q_point);
         return return_value;
       }
     else
@@ -3332,7 +3332,7 @@ namespace FEValuesViews
         for (unsigned int d=0; d<dim; ++d)
           if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
             return_value[d]
-              = fe_values.finite_element_output.shape_values(shape_function_data[shape_function].row_index[d],q_point);
+              = fe_values->finite_element_output.shape_values(shape_function_data[shape_function].row_index[d],q_point);
 
         return return_value;
       }
@@ -3347,9 +3347,9 @@ namespace FEValuesViews
                                   const unsigned int q_point) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (shape_function < fe_values.fe->dofs_per_cell,
-            ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
-    Assert (fe_values.update_flags & update_gradients,
+    Assert (shape_function < fe_values->fe->dofs_per_cell,
+            ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
+    Assert (fe_values->update_flags & update_gradients,
             typename FVB::ExcAccessToUninitializedField("update_gradients"));
 
     // same as for the scalar case except
@@ -3361,7 +3361,7 @@ namespace FEValuesViews
       {
         gradient_type return_value;
         return_value[shape_function_data[shape_function].single_nonzero_component_index]
-          = fe_values.finite_element_output.shape_gradients[snc][q_point];
+          = fe_values->finite_element_output.shape_gradients[snc][q_point];
         return return_value;
       }
     else
@@ -3370,7 +3370,7 @@ namespace FEValuesViews
         for (unsigned int d=0; d<dim; ++d)
           if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
             return_value[d]
-              = fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[d]][q_point];
+              = fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[d]][q_point];
 
         return return_value;
       }
@@ -3387,9 +3387,9 @@ namespace FEValuesViews
     // this function works like in
     // the case above
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (shape_function < fe_values.fe->dofs_per_cell,
-            ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
-    Assert (fe_values.update_flags & update_gradients,
+    Assert (shape_function < fe_values->fe->dofs_per_cell,
+            ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
+    Assert (fe_values->update_flags & update_gradients,
             typename FVB::ExcAccessToUninitializedField("update_gradients"));
 
     // same as for the scalar case except
@@ -3399,14 +3399,14 @@ namespace FEValuesViews
       return divergence_type();
     else if (snc != -1)
       return
-        fe_values.finite_element_output.shape_gradients[snc][q_point][shape_function_data[shape_function].single_nonzero_component_index];
+        fe_values->finite_element_output.shape_gradients[snc][q_point][shape_function_data[shape_function].single_nonzero_component_index];
     else
       {
         divergence_type return_value = 0;
         for (unsigned int d=0; d<dim; ++d)
           if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
             return_value
-            += fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[d]][q_point][d];
+            += fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[d]][q_point][d];
 
         return return_value;
       }
@@ -3422,9 +3422,9 @@ namespace FEValuesViews
     // this function works like in the case above
     typedef FEValuesBase<dim,spacedim> FVB;
 
-    Assert (shape_function < fe_values.fe->dofs_per_cell,
-            ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
-    Assert (fe_values.update_flags & update_gradients,
+    Assert (shape_function < fe_values->fe->dofs_per_cell,
+            ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
+    Assert (fe_values->update_flags & update_gradients,
             typename FVB::ExcAccessToUninitializedField("update_gradients"));
     // same as for the scalar case except that we have one more index
     const int snc = shape_function_data[shape_function].single_nonzero_component;
@@ -3452,9 +3452,9 @@ namespace FEValuesViews
               // can only be zero
               // or one in 2d
               if (shape_function_data[shape_function].single_nonzero_component_index == 0)
-                return_value[0] = -1.0 * fe_values.finite_element_output.shape_gradients[snc][q_point][1];
+                return_value[0] = -1.0 * fe_values->finite_element_output.shape_gradients[snc][q_point][1];
               else
-                return_value[0] = fe_values.finite_element_output.shape_gradients[snc][q_point][0];
+                return_value[0] = fe_values->finite_element_output.shape_gradients[snc][q_point][0];
 
               return return_value;
             }
@@ -3467,11 +3467,11 @@ namespace FEValuesViews
 
               if (shape_function_data[shape_function].is_nonzero_shape_function_component[0])
                 return_value[0]
-                -= fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][1];
+                -= fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][1];
 
               if (shape_function_data[shape_function].is_nonzero_shape_function_component[1])
                 return_value[0]
-                += fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][0];
+                += fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][0];
 
               return return_value;
             }
@@ -3488,23 +3488,23 @@ namespace FEValuesViews
                 case 0:
                 {
                   return_value[0] = 0;
-                  return_value[1] = fe_values.finite_element_output.shape_gradients[snc][q_point][2];
-                  return_value[2] = -1.0 * fe_values.finite_element_output.shape_gradients[snc][q_point][1];
+                  return_value[1] = fe_values->finite_element_output.shape_gradients[snc][q_point][2];
+                  return_value[2] = -1.0 * fe_values->finite_element_output.shape_gradients[snc][q_point][1];
                   return return_value;
                 }
 
                 case 1:
                 {
-                  return_value[0] = -1.0 * fe_values.finite_element_output.shape_gradients[snc][q_point][2];
+                  return_value[0] = -1.0 * fe_values->finite_element_output.shape_gradients[snc][q_point][2];
                   return_value[1] = 0;
-                  return_value[2] = fe_values.finite_element_output.shape_gradients[snc][q_point][0];
+                  return_value[2] = fe_values->finite_element_output.shape_gradients[snc][q_point][0];
                   return return_value;
                 }
 
                 default:
                 {
-                  return_value[0] = fe_values.finite_element_output.shape_gradients[snc][q_point][1];
-                  return_value[1] = -1.0 * fe_values.finite_element_output.shape_gradients[snc][q_point][0];
+                  return_value[0] = fe_values->finite_element_output.shape_gradients[snc][q_point][1];
+                  return_value[1] = -1.0 * fe_values->finite_element_output.shape_gradients[snc][q_point][0];
                   return_value[2] = 0;
                   return return_value;
                 }
@@ -3521,25 +3521,25 @@ namespace FEValuesViews
               if (shape_function_data[shape_function].is_nonzero_shape_function_component[0])
                 {
                   return_value[1]
-                  += fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][2];
+                  += fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][2];
                   return_value[2]
-                  -= fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][1];
+                  -= fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[0]][q_point][1];
                 }
 
               if (shape_function_data[shape_function].is_nonzero_shape_function_component[1])
                 {
                   return_value[0]
-                  -= fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][2];
+                  -= fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][2];
                   return_value[2]
-                  += fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][0];
+                  += fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[1]][q_point][0];
                 }
 
               if (shape_function_data[shape_function].is_nonzero_shape_function_component[2])
                 {
                   return_value[0]
-                  += fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[2]][q_point][1];
+                  += fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[2]][q_point][1];
                   return_value[1]
-                  -= fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[2]][q_point][0];
+                  -= fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[2]][q_point][0];
                 }
 
               return return_value;
@@ -3560,9 +3560,9 @@ namespace FEValuesViews
     // this function works like in
     // the case above
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (shape_function < fe_values.fe->dofs_per_cell,
-            ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
-    Assert (fe_values.update_flags & update_hessians,
+    Assert (shape_function < fe_values->fe->dofs_per_cell,
+            ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
+    Assert (fe_values->update_flags & update_hessians,
             typename FVB::ExcAccessToUninitializedField("update_hessians"));
 
     // same as for the scalar case except
@@ -3574,7 +3574,7 @@ namespace FEValuesViews
       {
         hessian_type return_value;
         return_value[shape_function_data[shape_function].single_nonzero_component_index]
-          = fe_values.finite_element_output.shape_hessians[snc][q_point];
+          = fe_values->finite_element_output.shape_hessians[snc][q_point];
         return return_value;
       }
     else
@@ -3583,7 +3583,7 @@ namespace FEValuesViews
         for (unsigned int d=0; d<dim; ++d)
           if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
             return_value[d]
-              = fe_values.finite_element_output.shape_hessians[shape_function_data[shape_function].row_index[d]][q_point];
+              = fe_values->finite_element_output.shape_hessians[shape_function_data[shape_function].row_index[d]][q_point];
 
         return return_value;
       }
@@ -3598,9 +3598,9 @@ namespace FEValuesViews
     // this function works like in
     // the case above
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (shape_function < fe_values.fe->dofs_per_cell,
-            ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
-    Assert (fe_values.update_flags & update_3rd_derivatives,
+    Assert (shape_function < fe_values->fe->dofs_per_cell,
+            ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
+    Assert (fe_values->update_flags & update_3rd_derivatives,
             typename FVB::ExcAccessToUninitializedField("update_3rd_derivatives"));
 
     // same as for the scalar case except
@@ -3612,7 +3612,7 @@ namespace FEValuesViews
       {
         third_derivative_type return_value;
         return_value[shape_function_data[shape_function].single_nonzero_component_index]
-          = fe_values.finite_element_output.shape_3rd_derivatives[snc][q_point];
+          = fe_values->finite_element_output.shape_3rd_derivatives[snc][q_point];
         return return_value;
       }
     else
@@ -3621,7 +3621,7 @@ namespace FEValuesViews
         for (unsigned int d=0; d<dim; ++d)
           if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
             return_value[d]
-              = fe_values.finite_element_output.shape_3rd_derivatives[shape_function_data[shape_function].row_index[d]][q_point];
+              = fe_values->finite_element_output.shape_3rd_derivatives[shape_function_data[shape_function].row_index[d]][q_point];
 
         return return_value;
       }
@@ -3712,9 +3712,9 @@ namespace FEValuesViews
                                             const unsigned int q_point) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (shape_function < fe_values.fe->dofs_per_cell,
-            ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
-    Assert (fe_values.update_flags & update_gradients,
+    Assert (shape_function < fe_values->fe->dofs_per_cell,
+            ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
+    Assert (fe_values->update_flags & update_gradients,
             typename FVB::ExcAccessToUninitializedField("update_gradients"));
 
     // same as for the scalar case except
@@ -3724,14 +3724,14 @@ namespace FEValuesViews
       return symmetric_gradient_type();
     else if (snc != -1)
       return symmetrize_single_row (shape_function_data[shape_function].single_nonzero_component_index,
-                                    fe_values.finite_element_output.shape_gradients[snc][q_point]);
+                                    fe_values->finite_element_output.shape_gradients[snc][q_point]);
     else
       {
         gradient_type return_value;
         for (unsigned int d=0; d<dim; ++d)
           if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
             return_value[d]
-              = fe_values.finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[d]][q_point];
+              = fe_values->finite_element_output.shape_gradients[shape_function_data[shape_function].row_index[d]][q_point];
 
         return symmetrize(return_value);
       }
@@ -3746,9 +3746,9 @@ namespace FEValuesViews
                                             const unsigned int q_point) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (shape_function < fe_values.fe->dofs_per_cell,
-            ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
-    Assert (fe_values.update_flags & update_values,
+    Assert (shape_function < fe_values->fe->dofs_per_cell,
+            ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
+    Assert (fe_values->update_flags & update_values,
             typename FVB::ExcAccessToUninitializedField("update_values"));
 
     // similar to the vector case where we
@@ -3771,7 +3771,7 @@ namespace FEValuesViews
         const unsigned int comp =
           shape_function_data[shape_function].single_nonzero_component_index;
         return_value[value_type::unrolled_to_component_indices(comp)]
-          = fe_values.finite_element_output.shape_values(snc,q_point);
+          = fe_values->finite_element_output.shape_values(snc,q_point);
         return return_value;
       }
     else
@@ -3780,7 +3780,7 @@ namespace FEValuesViews
         for (unsigned int d = 0; d < value_type::n_independent_components; ++d)
           if (shape_function_data[shape_function].is_nonzero_shape_function_component[d])
             return_value[value_type::unrolled_to_component_indices(d)]
-              = fe_values.finite_element_output.shape_values(shape_function_data[shape_function].row_index[d],q_point);
+              = fe_values->finite_element_output.shape_values(shape_function_data[shape_function].row_index[d],q_point);
         return return_value;
       }
   }
@@ -3793,9 +3793,9 @@ namespace FEValuesViews
                                                 const unsigned int q_point) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (shape_function < fe_values.fe->dofs_per_cell,
-            ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
-    Assert (fe_values.update_flags & update_gradients,
+    Assert (shape_function < fe_values->fe->dofs_per_cell,
+            ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
+    Assert (fe_values->update_flags & update_gradients,
             typename FVB::ExcAccessToUninitializedField("update_gradients"));
 
     const int snc = shape_function_data[shape_function].single_nonzero_component;
@@ -3857,7 +3857,7 @@ namespace FEValuesViews
         // b_jj := \dfrac{\partial phi_{ii,jj}}{\partial x_jj}.
         // again, all other entries of 'b' are
         // zero
-        const dealii::Tensor<1, spacedim> phi_grad = fe_values.finite_element_output.shape_gradients[snc][q_point];
+        const dealii::Tensor<1, spacedim> phi_grad = fe_values->finite_element_output.shape_gradients[snc][q_point];
 
         divergence_type return_value;
         return_value[ii] = phi_grad[jj];
@@ -3883,9 +3883,9 @@ namespace FEValuesViews
                                    const unsigned int q_point) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (shape_function < fe_values.fe->dofs_per_cell,
-            ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
-    Assert (fe_values.update_flags & update_values,
+    Assert (shape_function < fe_values->fe->dofs_per_cell,
+            ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
+    Assert (fe_values->update_flags & update_values,
             typename FVB::ExcAccessToUninitializedField("update_values"));
 
     // similar to the vector case where we
@@ -3908,7 +3908,7 @@ namespace FEValuesViews
         const unsigned int comp =
           shape_function_data[shape_function].single_nonzero_component_index;
         const TableIndices<2> indices = dealii::Tensor<2,spacedim>::unrolled_to_component_indices(comp);
-        return_value[indices] = fe_values.finite_element_output.shape_values(snc,q_point);
+        return_value[indices] = fe_values->finite_element_output.shape_values(snc,q_point);
         return return_value;
       }
     else
@@ -3919,7 +3919,7 @@ namespace FEValuesViews
             {
               const TableIndices<2> indices = dealii::Tensor<2,spacedim>::unrolled_to_component_indices(d);
               return_value[indices]
-                = fe_values.finite_element_output.shape_values(shape_function_data[shape_function].row_index[d],q_point);
+                = fe_values->finite_element_output.shape_values(shape_function_data[shape_function].row_index[d],q_point);
             }
         return return_value;
       }
@@ -3933,9 +3933,9 @@ namespace FEValuesViews
                                        const unsigned int q_point) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (shape_function < fe_values.fe->dofs_per_cell,
-            ExcIndexRange (shape_function, 0, fe_values.fe->dofs_per_cell));
-    Assert (fe_values.update_flags & update_gradients,
+    Assert (shape_function < fe_values->fe->dofs_per_cell,
+            ExcIndexRange (shape_function, 0, fe_values->fe->dofs_per_cell));
+    Assert (fe_values->update_flags & update_gradients,
             typename FVB::ExcAccessToUninitializedField("update_gradients"));
 
     const int snc = shape_function_data[shape_function].single_nonzero_component;
@@ -3970,7 +3970,7 @@ namespace FEValuesViews
         const unsigned int ii = indices[0];
         const unsigned int jj = indices[1];
 
-        const dealii::Tensor<1, spacedim> phi_grad = fe_values.finite_element_output.shape_gradients[snc][q_point];
+        const dealii::Tensor<1, spacedim> phi_grad = fe_values->finite_element_output.shape_gradients[snc][q_point];
 
         divergence_type return_value;
         return_value[jj] = phi_grad[ii];
index ace4f054164c9f1595baca4e400844f8c2a8fd0c..0a42759a6319bebdf22caa7e373e4c25f848d3b7 100644 (file)
@@ -98,36 +98,36 @@ namespace FEValuesViews
   Scalar<dim,spacedim>::Scalar (const FEValuesBase<dim,spacedim> &fe_values,
                                 const unsigned int                component)
     :
-    fe_values (fe_values),
+    fe_values (&fe_values),
     component (component),
-    shape_function_data (fe_values.fe->dofs_per_cell)
+    shape_function_data (this->fe_values->fe->dofs_per_cell)
   {
-    Assert (component < fe_values.fe->n_components(),
-            ExcIndexRange(component, 0, fe_values.fe->n_components()));
+    const FiniteElement<dim, spacedim> &fe = *this->fe_values->fe;
+    Assert (component < fe.n_components(),
+            ExcIndexRange(component, 0, fe.n_components()));
 
 //TODO: we'd like to use the fields with the same name as these
 // variables from FEValuesBase, but they aren't initialized yet
 // at the time we get here, so re-create it all
     const std::vector<unsigned int> shape_function_to_row_table
-      = make_shape_function_to_row_table (*fe_values.fe);
+      = make_shape_function_to_row_table (fe);
 
-    for (unsigned int i=0; i<fe_values.fe->dofs_per_cell; ++i)
+    for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
       {
-        const bool is_primitive = (fe_values.fe->is_primitive() ||
-                                   fe_values.fe->is_primitive(i));
+        const bool is_primitive = fe.is_primitive() || fe.is_primitive(i);
 
         if (is_primitive == true)
           shape_function_data[i].is_nonzero_shape_function_component
             = (component ==
-               fe_values.fe->system_to_component_index(i).first);
+               fe.system_to_component_index(i).first);
         else
           shape_function_data[i].is_nonzero_shape_function_component
-            = (fe_values.fe->get_nonzero_components(i)[component]
+            = (fe.get_nonzero_components(i)[component]
                == true);
 
         if (shape_function_data[i].is_nonzero_shape_function_component == true)
           shape_function_data[i].row_index
-            = shape_function_to_row_table[i*fe_values.fe->n_components()+component];
+            = shape_function_to_row_table[i*fe.n_components()+component];
         else
           shape_function_data[i].row_index = numbers::invalid_unsigned_int;
       }
@@ -138,7 +138,7 @@ namespace FEValuesViews
   template <int dim, int spacedim>
   Scalar<dim,spacedim>::Scalar ()
     :
-    fe_values (*static_cast<dealii::FEValuesBase<dim,spacedim>*>(0)),
+    fe_values (NULL),
     component (numbers::invalid_unsigned_int)
   {}
 
@@ -158,49 +158,49 @@ namespace FEValuesViews
   Vector<dim,spacedim>::Vector (const FEValuesBase<dim,spacedim> &fe_values,
                                 const unsigned int       first_vector_component)
     :
-    fe_values (fe_values),
+    fe_values (&fe_values),
     first_vector_component (first_vector_component),
-    shape_function_data (fe_values.fe->dofs_per_cell)
+    shape_function_data (this->fe_values->fe->dofs_per_cell)
   {
-    Assert (first_vector_component+spacedim-1 < fe_values.fe->n_components(),
+    const FiniteElement<dim, spacedim> &fe = *this->fe_values->fe;
+    Assert (first_vector_component+spacedim-1 < fe.n_components(),
             ExcIndexRange(first_vector_component+spacedim-1, 0,
-                          fe_values.fe->n_components()));
+                          fe.n_components()));
 
 //TODO: we'd like to use the fields with the same name as these
 // variables from FEValuesBase, but they aren't initialized yet
 // at the time we get here, so re-create it all
     const std::vector<unsigned int> shape_function_to_row_table
-      = make_shape_function_to_row_table (*fe_values.fe);
+      = make_shape_function_to_row_table (fe);
 
     for (unsigned int d=0; d<spacedim; ++d)
       {
         const unsigned int component = first_vector_component + d;
 
-        for (unsigned int i=0; i<fe_values.fe->dofs_per_cell; ++i)
+        for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
           {
-            const bool is_primitive = (fe_values.fe->is_primitive() ||
-                                       fe_values.fe->is_primitive(i));
+            const bool is_primitive = fe.is_primitive() || fe.is_primitive(i);
 
             if (is_primitive == true)
               shape_function_data[i].is_nonzero_shape_function_component[d]
                 = (component ==
-                   fe_values.fe->system_to_component_index(i).first);
+                   fe.system_to_component_index(i).first);
             else
               shape_function_data[i].is_nonzero_shape_function_component[d]
-                = (fe_values.fe->get_nonzero_components(i)[component]
+                = (fe.get_nonzero_components(i)[component]
                    == true);
 
             if (shape_function_data[i].is_nonzero_shape_function_component[d]
                 == true)
               shape_function_data[i].row_index[d]
-                = shape_function_to_row_table[i*fe_values.fe->n_components()+component];
+                = shape_function_to_row_table[i*fe.n_components()+component];
             else
               shape_function_data[i].row_index[d]
                 = numbers::invalid_unsigned_int;
           }
       }
 
-    for (unsigned int i=0; i<fe_values.fe->dofs_per_cell; ++i)
+    for (unsigned int i=0; i<fe.dofs_per_cell; ++i)
       {
         unsigned int n_nonzero_components = 0;
         for (unsigned int d=0; d<spacedim; ++d)
@@ -232,7 +232,7 @@ namespace FEValuesViews
   template <int dim, int spacedim>
   Vector<dim,spacedim>::Vector ()
     :
-    fe_values (*static_cast<dealii::FEValuesBase<dim,spacedim>*>(0)),
+    fe_values (NULL),
     first_vector_component (numbers::invalid_unsigned_int)
   {}
 
@@ -253,52 +253,52 @@ namespace FEValuesViews
   SymmetricTensor(const FEValuesBase<dim, spacedim> &fe_values,
                   const unsigned int first_tensor_component)
     :
-    fe_values(fe_values),
+    fe_values(&fe_values),
     first_tensor_component(first_tensor_component),
-    shape_function_data(fe_values.fe->dofs_per_cell)
+    shape_function_data(this->fe_values->fe->dofs_per_cell)
   {
+    const FiniteElement<dim, spacedim> &fe = *this->fe_values->fe;
     Assert(first_tensor_component + (dim*dim+dim)/2 - 1
            <
-           fe_values.fe->n_components(),
+           fe.n_components(),
            ExcIndexRange(first_tensor_component +
                          dealii::SymmetricTensor<2,dim>::n_independent_components - 1,
                          0,
-                         fe_values.fe->n_components()));
+                         fe.n_components()));
 //TODO: we'd like to use the fields with the same name as these
 // variables from FEValuesBase, but they aren't initialized yet
 // at the time we get here, so re-create it all
     const std::vector<unsigned int> shape_function_to_row_table
-      = make_shape_function_to_row_table (*fe_values.fe);
+      = make_shape_function_to_row_table (fe);
 
     for (unsigned int d = 0; d < dealii::SymmetricTensor<2,dim>::n_independent_components; ++d)
       {
         const unsigned int component = first_tensor_component + d;
 
-        for (unsigned int i = 0; i < fe_values.fe->dofs_per_cell; ++i)
+        for (unsigned int i = 0; i < fe.dofs_per_cell; ++i)
           {
-            const bool is_primitive = (fe_values.fe->is_primitive() ||
-                                       fe_values.fe->is_primitive(i));
+            const bool is_primitive = fe.is_primitive() || fe.is_primitive(i);
 
             if (is_primitive == true)
               shape_function_data[i].is_nonzero_shape_function_component[d]
                 = (component ==
-                   fe_values.fe->system_to_component_index(i).first);
+                   fe.system_to_component_index(i).first);
             else
               shape_function_data[i].is_nonzero_shape_function_component[d]
-                = (fe_values.fe->get_nonzero_components(i)[component]
+                = (fe.get_nonzero_components(i)[component]
                    == true);
 
             if (shape_function_data[i].is_nonzero_shape_function_component[d]
                 == true)
               shape_function_data[i].row_index[d]
-                = shape_function_to_row_table[i*fe_values.fe->n_components()+component];
+                = shape_function_to_row_table[i*fe.n_components()+component];
             else
               shape_function_data[i].row_index[d]
                 = numbers::invalid_unsigned_int;
           }
       }
 
-    for (unsigned int i = 0; i < fe_values.fe->dofs_per_cell; ++i)
+    for (unsigned int i = 0; i < fe.dofs_per_cell; ++i)
       {
         unsigned int n_nonzero_components = 0;
         for (unsigned int d = 0; d < dealii::SymmetricTensor<2,dim>::n_independent_components; ++d)
@@ -331,7 +331,7 @@ namespace FEValuesViews
   template <int dim, int spacedim>
   SymmetricTensor<2, dim, spacedim>::SymmetricTensor()
     :
-    fe_values(*static_cast<dealii::FEValuesBase<dim, spacedim>*> (0)),
+    fe_values(NULL),
     first_tensor_component(numbers::invalid_unsigned_int)
   {}
 
@@ -352,52 +352,52 @@ namespace FEValuesViews
   Tensor(const FEValuesBase<dim, spacedim> &fe_values,
          const unsigned int first_tensor_component)
     :
-    fe_values(fe_values),
+    fe_values(&fe_values),
     first_tensor_component(first_tensor_component),
-    shape_function_data(fe_values.fe->dofs_per_cell)
+    shape_function_data(this->fe_values->fe->dofs_per_cell)
   {
+    const FiniteElement<dim, spacedim> &fe = *this->fe_values->fe;
     Assert(first_tensor_component + dim*dim - 1
            <
-           fe_values.fe->n_components(),
+           fe.n_components(),
            ExcIndexRange(first_tensor_component +
                          dim*dim - 1,
                          0,
-                         fe_values.fe->n_components()));
+                         fe.n_components()));
 //TODO: we'd like to use the fields with the same name as these
 // variables from FEValuesBase, but they aren't initialized yet
 // at the time we get here, so re-create it all
     const std::vector<unsigned int> shape_function_to_row_table
-      = make_shape_function_to_row_table (*fe_values.fe);
+      = make_shape_function_to_row_table (fe);
 
     for (unsigned int d = 0; d < dim*dim; ++d)
       {
         const unsigned int component = first_tensor_component + d;
 
-        for (unsigned int i = 0; i < fe_values.fe->dofs_per_cell; ++i)
+        for (unsigned int i = 0; i < fe.dofs_per_cell; ++i)
           {
-            const bool is_primitive = (fe_values.fe->is_primitive() ||
-                                       fe_values.fe->is_primitive(i));
+            const bool is_primitive = fe.is_primitive() || fe.is_primitive(i);
 
             if (is_primitive == true)
               shape_function_data[i].is_nonzero_shape_function_component[d]
                 = (component ==
-                   fe_values.fe->system_to_component_index(i).first);
+                   fe.system_to_component_index(i).first);
             else
               shape_function_data[i].is_nonzero_shape_function_component[d]
-                = (fe_values.fe->get_nonzero_components(i)[component]
+                = (fe.get_nonzero_components(i)[component]
                    == true);
 
             if (shape_function_data[i].is_nonzero_shape_function_component[d]
                 == true)
               shape_function_data[i].row_index[d]
-                = shape_function_to_row_table[i*fe_values.fe->n_components()+component];
+                = shape_function_to_row_table[i*fe.n_components()+component];
             else
               shape_function_data[i].row_index[d]
                 = numbers::invalid_unsigned_int;
           }
       }
 
-    for (unsigned int i = 0; i < fe_values.fe->dofs_per_cell; ++i)
+    for (unsigned int i = 0; i < fe.dofs_per_cell; ++i)
       {
         unsigned int n_nonzero_components = 0;
         for (unsigned int d = 0; d < dim*dim; ++d)
@@ -430,7 +430,7 @@ namespace FEValuesViews
   template <int dim, int spacedim>
   Tensor<2, dim, spacedim>::Tensor()
     :
-    fe_values(*static_cast<dealii::FEValuesBase<dim, spacedim>*> (0)),
+    fe_values(NULL),
     first_tensor_component(numbers::invalid_unsigned_int)
   {}
 
@@ -1270,18 +1270,18 @@ namespace FEValuesViews
                        std::vector<typename ProductType<value_type,typename InputVector::value_type>::type> &values) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (fe_values.update_flags & update_values,
+    Assert (fe_values->update_flags & update_values,
             typename FVB::ExcAccessToUninitializedField("update_values"));
-    Assert (fe_values.present_cell.get() != 0,
+    Assert (fe_values->present_cell.get() != 0,
             ExcMessage ("FEValues object is not reinit'ed to any cell"));
     AssertDimension (fe_function.size(),
-                     fe_values.present_cell->n_dofs_for_dof_handler());
+                     fe_values->present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs on this cell and call internal worker function
-    dealii::Vector<typename InputVector::value_type> dof_values(fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values(fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_values<dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_values, shape_function_data, values);
+    (dof_values, fe_values->finite_element_output.shape_values, shape_function_data, values);
   }
 
 
@@ -1294,18 +1294,18 @@ namespace FEValuesViews
                           std::vector<typename ProductType<gradient_type,typename InputVector::value_type>::type> &gradients) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (fe_values.update_flags & update_gradients,
+    Assert (fe_values->update_flags & update_gradients,
             typename FVB::ExcAccessToUninitializedField("update_gradients"));
-    Assert (fe_values.present_cell.get() != 0,
+    Assert (fe_values->present_cell.get() != 0,
             ExcMessage ("FEValues object is not reinit'ed to any cell"));
     AssertDimension (fe_function.size(),
-                     fe_values.present_cell->n_dofs_for_dof_handler());
+                     fe_values->present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs on this cell
-    dealii::Vector<typename InputVector::value_type> dof_values (fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values (fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_derivatives<1,dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_gradients, shape_function_data, gradients);
+    (dof_values, fe_values->finite_element_output.shape_gradients, shape_function_data, gradients);
   }
 
 
@@ -1318,18 +1318,18 @@ namespace FEValuesViews
                          std::vector<typename ProductType<hessian_type,typename InputVector::value_type>::type> &hessians) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (fe_values.update_flags & update_hessians,
+    Assert (fe_values->update_flags & update_hessians,
             typename FVB::ExcAccessToUninitializedField("update_hessians"));
-    Assert (fe_values.present_cell.get() != 0,
+    Assert (fe_values->present_cell.get() != 0,
             ExcMessage ("FEValues object is not reinit'ed to any cell"));
     AssertDimension (fe_function.size(),
-                     fe_values.present_cell->n_dofs_for_dof_handler());
+                     fe_values->present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs on this cell
-    dealii::Vector<typename InputVector::value_type> dof_values (fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values (fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_derivatives<2,dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_hessians, shape_function_data, hessians);
+    (dof_values, fe_values->finite_element_output.shape_hessians, shape_function_data, hessians);
   }
 
 
@@ -1342,18 +1342,18 @@ namespace FEValuesViews
                            std::vector<typename ProductType<value_type,typename InputVector::value_type>::type> &laplacians) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (fe_values.update_flags & update_hessians,
+    Assert (fe_values->update_flags & update_hessians,
             typename FVB::ExcAccessToUninitializedField("update_hessians"));
-    Assert (fe_values.present_cell.get() != 0,
+    Assert (fe_values->present_cell.get() != 0,
             ExcMessage ("FEValues object is not reinit'ed to any cell"));
     AssertDimension (fe_function.size(),
-                     fe_values.present_cell->n_dofs_for_dof_handler());
+                     fe_values->present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs on this cell
-    dealii::Vector<typename InputVector::value_type> dof_values (fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values (fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_laplacians<dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_hessians, shape_function_data, laplacians);
+    (dof_values, fe_values->finite_element_output.shape_hessians, shape_function_data, laplacians);
   }
 
 
@@ -1366,18 +1366,18 @@ namespace FEValuesViews
                                   std::vector<typename ProductType<third_derivative_type,typename InputVector::value_type>::type> &third_derivatives) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (fe_values.update_flags & update_3rd_derivatives,
+    Assert (fe_values->update_flags & update_3rd_derivatives,
             typename FVB::ExcAccessToUninitializedField("update_3rd_derivatives"));
-    Assert (fe_values.present_cell.get() != 0,
+    Assert (fe_values->present_cell.get() != 0,
             ExcMessage ("FEValues object is not reinit'ed to any cell"));
     AssertDimension (fe_function.size(),
-                     fe_values.present_cell->n_dofs_for_dof_handler());
+                     fe_values->present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs on this cell
-    dealii::Vector<typename InputVector::value_type> dof_values (fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values (fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_derivatives<3,dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_3rd_derivatives, shape_function_data, third_derivatives);
+    (dof_values, fe_values->finite_element_output.shape_3rd_derivatives, shape_function_data, third_derivatives);
   }
 
 
@@ -1390,18 +1390,18 @@ namespace FEValuesViews
                        std::vector<typename ProductType<value_type,typename InputVector::value_type>::type> &values) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (fe_values.update_flags & update_values,
+    Assert (fe_values->update_flags & update_values,
             typename FVB::ExcAccessToUninitializedField("update_values"));
-    Assert (fe_values.present_cell.get() != 0,
+    Assert (fe_values->present_cell.get() != 0,
             ExcMessage ("FEValues object is not reinit'ed to any cell"));
     AssertDimension (fe_function.size(),
-                     fe_values.present_cell->n_dofs_for_dof_handler());
+                     fe_values->present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs on this cell
-    dealii::Vector<typename InputVector::value_type> dof_values (fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values (fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_values<dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_values, shape_function_data, values);
+    (dof_values, fe_values->finite_element_output.shape_values, shape_function_data, values);
   }
 
 
@@ -1415,18 +1415,18 @@ namespace FEValuesViews
                           std::vector<typename ProductType<gradient_type,typename InputVector::value_type>::type> &gradients) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (fe_values.update_flags & update_gradients,
+    Assert (fe_values->update_flags & update_gradients,
             typename FVB::ExcAccessToUninitializedField("update_gradients"));
-    Assert (fe_values.present_cell.get() != 0,
+    Assert (fe_values->present_cell.get() != 0,
             ExcMessage ("FEValues object is not reinit'ed to any cell"));
     AssertDimension (fe_function.size(),
-                     fe_values.present_cell->n_dofs_for_dof_handler());
+                     fe_values->present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs on this cell
-    dealii::Vector<typename InputVector::value_type> dof_values (fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values (fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_derivatives<1,dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_gradients, shape_function_data, gradients);
+    (dof_values, fe_values->finite_element_output.shape_gradients, shape_function_data, gradients);
   }
 
 
@@ -1439,18 +1439,18 @@ namespace FEValuesViews
                                     std::vector<typename ProductType<symmetric_gradient_type,typename InputVector::value_type>::type> &symmetric_gradients) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (fe_values.update_flags & update_gradients,
+    Assert (fe_values->update_flags & update_gradients,
             typename FVB::ExcAccessToUninitializedField("update_gradients"));
-    Assert (fe_values.present_cell.get() != 0,
+    Assert (fe_values->present_cell.get() != 0,
             ExcMessage ("FEValues object is not reinit'ed to any cell"));
     AssertDimension (fe_function.size(),
-                     fe_values.present_cell->n_dofs_for_dof_handler());
+                     fe_values->present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs on this cell
-    dealii::Vector<typename InputVector::value_type> dof_values (fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values (fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_symmetric_gradients<dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_gradients, shape_function_data,
+    (dof_values, fe_values->finite_element_output.shape_gradients, shape_function_data,
      symmetric_gradients);
   }
 
@@ -1464,19 +1464,19 @@ namespace FEValuesViews
                             std::vector<typename ProductType<divergence_type,typename InputVector::value_type>::type> &divergences) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (fe_values.update_flags & update_gradients,
+    Assert (fe_values->update_flags & update_gradients,
             typename FVB::ExcAccessToUninitializedField("update_gradients"));
-    Assert (fe_values.present_cell.get() != 0,
+    Assert (fe_values->present_cell.get() != 0,
             ExcMessage ("FEValues object is not reinit'ed to any cell"));
     AssertDimension (fe_function.size(),
-                     fe_values.present_cell->n_dofs_for_dof_handler());
+                     fe_values->present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs
     // on this cell
-    dealii::Vector<typename InputVector::value_type> dof_values (fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values (fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_divergences<dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_gradients, shape_function_data, divergences);
+    (dof_values, fe_values->finite_element_output.shape_gradients, shape_function_data, divergences);
   }
 
   template <int dim, int spacedim>
@@ -1488,18 +1488,18 @@ namespace FEValuesViews
   {
     typedef FEValuesBase<dim,spacedim> FVB;
 
-    Assert (fe_values.update_flags & update_gradients,
+    Assert (fe_values->update_flags & update_gradients,
             typename FVB::ExcAccessToUninitializedField("update_gradients"));
-    Assert (fe_values.present_cell.get () != 0,
+    Assert (fe_values->present_cell.get () != 0,
             ExcMessage ("FEValues object is not reinited to any cell"));
     AssertDimension (fe_function.size (),
-                     fe_values.present_cell->n_dofs_for_dof_handler ());
+                     fe_values->present_cell->n_dofs_for_dof_handler ());
 
     // get function values of dofs on this cell
-    dealii::Vector<typename InputVector::value_type> dof_values (fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values (fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values (fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values (fe_function, dof_values);
     internal::do_function_curls<dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_gradients, shape_function_data, curls);
+    (dof_values, fe_values->finite_element_output.shape_gradients, shape_function_data, curls);
   }
 
 
@@ -1511,18 +1511,18 @@ namespace FEValuesViews
                          std::vector<typename ProductType<hessian_type,typename InputVector::value_type>::type> &hessians) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (fe_values.update_flags & update_hessians,
+    Assert (fe_values->update_flags & update_hessians,
             typename FVB::ExcAccessToUninitializedField("update_hessians"));
-    Assert (fe_values.present_cell.get() != 0,
+    Assert (fe_values->present_cell.get() != 0,
             ExcMessage ("FEValues object is not reinit'ed to any cell"));
     AssertDimension (fe_function.size(),
-                     fe_values.present_cell->n_dofs_for_dof_handler());
+                     fe_values->present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs on this cell
-    dealii::Vector<typename InputVector::value_type> dof_values (fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values (fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_derivatives<2,dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_hessians, shape_function_data, hessians);
+    (dof_values, fe_values->finite_element_output.shape_hessians, shape_function_data, hessians);
   }
 
 
@@ -1535,21 +1535,21 @@ namespace FEValuesViews
                            std::vector<typename ProductType<value_type,typename InputVector::value_type>::type> &laplacians) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (fe_values.update_flags & update_hessians,
+    Assert (fe_values->update_flags & update_hessians,
             typename FVB::ExcAccessToUninitializedField("update_hessians"));
-    Assert (laplacians.size() == fe_values.n_quadrature_points,
-            ExcDimensionMismatch(laplacians.size(), fe_values.n_quadrature_points));
-    Assert (fe_values.present_cell.get() != 0,
+    Assert (laplacians.size() == fe_values->n_quadrature_points,
+            ExcDimensionMismatch(laplacians.size(), fe_values->n_quadrature_points));
+    Assert (fe_values->present_cell.get() != 0,
             ExcMessage ("FEValues object is not reinit'ed to any cell"));
-    Assert (fe_function.size() == fe_values.present_cell->n_dofs_for_dof_handler(),
+    Assert (fe_function.size() == fe_values->present_cell->n_dofs_for_dof_handler(),
             ExcDimensionMismatch(fe_function.size(),
-                                 fe_values.present_cell->n_dofs_for_dof_handler()));
+                                 fe_values->present_cell->n_dofs_for_dof_handler()));
 
     // get function values of dofs on this cell
-    dealii::Vector<typename InputVector::value_type> dof_values (fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values (fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_laplacians<dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_hessians, shape_function_data, laplacians);
+    (dof_values, fe_values->finite_element_output.shape_hessians, shape_function_data, laplacians);
   }
 
 
@@ -1561,18 +1561,18 @@ namespace FEValuesViews
                                   std::vector<typename ProductType<third_derivative_type,typename InputVector::value_type>::type> &third_derivatives) const
   {
     typedef FEValuesBase<dim,spacedim> FVB;
-    Assert (fe_values.update_flags & update_3rd_derivatives,
+    Assert (fe_values->update_flags & update_3rd_derivatives,
             typename FVB::ExcAccessToUninitializedField("update_3rd_derivatives"));
-    Assert (fe_values.present_cell.get() != 0,
+    Assert (fe_values->present_cell.get() != 0,
             ExcMessage ("FEValues object is not reinit'ed to any cell"));
     AssertDimension (fe_function.size(),
-                     fe_values.present_cell->n_dofs_for_dof_handler());
+                     fe_values->present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs on this cell
-    dealii::Vector<typename InputVector::value_type> dof_values (fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values (fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_derivatives<3,dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_3rd_derivatives, shape_function_data, third_derivatives);
+    (dof_values, fe_values->finite_element_output.shape_3rd_derivatives, shape_function_data, third_derivatives);
   }
 
 
@@ -1585,18 +1585,18 @@ namespace FEValuesViews
                       std::vector<typename ProductType<value_type,typename InputVector::value_type>::type> &values) const
   {
     typedef FEValuesBase<dim, spacedim> FVB;
-    Assert(fe_values.update_flags & update_values,
+    Assert(fe_values->update_flags & update_values,
            typename FVB::ExcAccessToUninitializedField("update_values"));
-    Assert(fe_values.present_cell.get() != 0,
+    Assert(fe_values->present_cell.get() != 0,
            ExcMessage("FEValues object is not reinit'ed to any cell"));
     AssertDimension(fe_function.size(),
-                    fe_values.present_cell->n_dofs_for_dof_handler());
+                    fe_values->present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs on this cell
-    dealii::Vector<typename InputVector::value_type> dof_values(fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values(fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_values<dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_values, shape_function_data, values);
+    (dof_values, fe_values->finite_element_output.shape_values, shape_function_data, values);
   }
 
 
@@ -1609,19 +1609,19 @@ namespace FEValuesViews
                            std::vector<typename ProductType<divergence_type,typename InputVector::value_type>::type> &divergences) const
   {
     typedef FEValuesBase<dim, spacedim> FVB;
-    Assert(fe_values.update_flags & update_gradients,
+    Assert(fe_values->update_flags & update_gradients,
            typename FVB::ExcAccessToUninitializedField("update_gradients"));
-    Assert(fe_values.present_cell.get() != 0,
+    Assert(fe_values->present_cell.get() != 0,
            ExcMessage("FEValues object is not reinit'ed to any cell"));
     AssertDimension(fe_function.size(),
-                    fe_values.present_cell->n_dofs_for_dof_handler());
+                    fe_values->present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs
     // on this cell
-    dealii::Vector<typename InputVector::value_type> dof_values(fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values(fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_divergences<dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_gradients, shape_function_data, divergences);
+    (dof_values, fe_values->finite_element_output.shape_gradients, shape_function_data, divergences);
   }
 
   template <int dim, int spacedim>
@@ -1632,18 +1632,18 @@ namespace FEValuesViews
                       std::vector<typename ProductType<value_type,typename InputVector::value_type>::type> &values) const
   {
     typedef FEValuesBase<dim, spacedim> FVB;
-    Assert(fe_values.update_flags & update_values,
+    Assert(fe_values->update_flags & update_values,
            typename FVB::ExcAccessToUninitializedField("update_values"));
-    Assert(fe_values.present_cell.get() != 0,
+    Assert(fe_values->present_cell.get() != 0,
            ExcMessage("FEValues object is not reinit'ed to any cell"));
     AssertDimension(fe_function.size(),
-                    fe_values.present_cell->n_dofs_for_dof_handler());
+                    fe_values->present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs on this cell
-    dealii::Vector<typename InputVector::value_type> dof_values(fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values(fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_values<dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_values, shape_function_data, values);
+    (dof_values, fe_values->finite_element_output.shape_values, shape_function_data, values);
   }
 
 
@@ -1656,19 +1656,19 @@ namespace FEValuesViews
                            std::vector<typename ProductType<divergence_type,typename InputVector::value_type>::type> &divergences) const
   {
     typedef FEValuesBase<dim, spacedim> FVB;
-    Assert(fe_values.update_flags & update_gradients,
+    Assert(fe_values->update_flags & update_gradients,
            typename FVB::ExcAccessToUninitializedField("update_gradients"));
-    Assert(fe_values.present_cell.get() != 0,
+    Assert(fe_values->present_cell.get() != 0,
            ExcMessage("FEValues object is not reinit'ed to any cell"));
     AssertDimension(fe_function.size(),
-                    fe_values.present_cell->n_dofs_for_dof_handler());
+                    fe_values->present_cell->n_dofs_for_dof_handler());
 
     // get function values of dofs
     // on this cell
-    dealii::Vector<typename InputVector::value_type> dof_values(fe_values.dofs_per_cell);
-    fe_values.present_cell->get_interpolated_dof_values(fe_function, dof_values);
+    dealii::Vector<typename InputVector::value_type> dof_values(fe_values->dofs_per_cell);
+    fe_values->present_cell->get_interpolated_dof_values(fe_function, dof_values);
     internal::do_function_divergences<dim,spacedim>
-    (dof_values, fe_values.finite_element_output.shape_gradients, shape_function_data, divergences);
+    (dof_values, fe_values->finite_element_output.shape_gradients, shape_function_data, divergences);
   }
 }
 
@@ -1682,14 +1682,9 @@ namespace internal
     {
       const FiniteElement<dim,spacedim> &fe = fe_values.get_fe();
 
-      // create the views objects. allocate a
-      // bunch of default-constructed ones
-      // then destroy them again and do
-      // in-place construction of those we
-      // actually want to use (copying stuff
-      // is wasteful and we can't do that
-      // anyway because the class has
-      // reference members)
+      // create the views objects: Allocate a bunch of default-constructed ones
+      // then destroy them again and do in-place construction of those we
+      // actually want to use.
       const unsigned int n_scalars = fe.n_components();
       scalars.resize (n_scalars);
       for (unsigned int component=0; component<n_scalars; ++component)

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.