]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Implement the second FEValues::get_function_2nd_derivatives function, for vector...
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 13 Nov 2000 13:33:56 +0000 (13:33 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 13 Nov 2000 13:33:56 +0000 (13:33 +0000)
git-svn-id: https://svn.dealii.org/trunk@3483 0785d39b-7218-0410-832d-ea1e28bc413d

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

index d85ab1d6744d8fd20fa3d030088a995d0889890e..e32fdc8ba660e1bb25b9ac0aaed03de3312284bc 100644 (file)
@@ -468,6 +468,33 @@ class FEValuesBase
     void get_function_2nd_derivatives (const InputVector      &fe_function,
                                       vector<Tensor<2,dim> > &second_derivatives) const;
 
+    
+                                    /**
+                                     * Return the tensor of second
+                                     * derivatives of the finite
+                                     * element function characterized
+                                     * by @p{fe_function} restricted to
+                                     * @p{cell} at the quadrature points.
+                                     *
+                                     * The function assumes that the
+                                     * @p{second_derivatives} object already has
+                                     * the right size.
+                                     *
+                                     * This function does the same as
+                                     * the other one with the same
+                                     * name, but applies to
+                                     * vector-valued finite elements.
+                                     *
+                                     * The actual data type of the
+                                     * input vector may be either a
+                                     * @p{Vector<double>},
+                                     * @p{Vector<float>}, or
+                                     * @p{BlockVector<double,...>}.
+                                     */
+    template <class InputVector>
+    void get_function_2nd_derivatives (const InputVector      &fe_function,
+                                      vector<vector<Tensor<2,dim> > > &second_derivatives) const;
+    
                                     /**
                                      * Return the position of the @p{i}th
                                      * quadrature point in real space.
index 29ed57a3cbcc6344c97e89e041f84d5e53c48f84..3e02dbb2f1b3b344ec27d22aada34c756afa274c 100644 (file)
@@ -286,6 +286,48 @@ void FEValuesBase<dim>::get_function_2nd_derivatives (const InputVector      &fe
 
 
 
+template <int dim>
+template <class InputVector>
+void
+FEValuesBase<dim>::
+get_function_2nd_derivatives (const InputVector               &fe_function,
+                             vector<vector<Tensor<2,dim> > > &second_derivs) const
+{
+  Assert (n_quadrature_points == second_derivs.size(),
+         ExcWrongNoOfComponents());
+  Assert (selected_dataset<shape_values.size(),
+         ExcIndexRange (selected_dataset, 0, shape_values.size()));
+  for (unsigned i=0;i<second_derivs.size();++i)
+    Assert (second_derivs[i].size() == fe->n_components(),
+           ExcWrongVectorSize(second_derivs[i].size(), fe->n_components()));
+  Assert (update_flags & update_second_derivatives, ExcAccessToUninitializedField());
+
+                                  // get function values of dofs
+                                  // on this cell
+  Vector<typename InputVector::value_type> dof_values (dofs_per_cell);
+  if (present_cell->active())
+    present_cell->get_dof_values (fe_function, dof_values);
+  else
+    present_cell->get_interpolated_dof_values(fe_function, dof_values);
+
+                                  // initialize with zero
+  for (unsigned i=0;i<second_derivs.size();++i)
+    fill_n (second_derivs[i].begin(), second_derivs[i].size(), Tensor<2,dim>());
+
+                                  // add up contributions of trial
+                                  // functions
+  for (unsigned int point=0; point<n_quadrature_points; ++point)
+    for (unsigned int shape_func=0; shape_func<dofs_per_cell; ++shape_func)
+      {
+       Tensor<2,dim> tmp(shape_2nd_derivatives[shape_func][point]);
+       tmp *= dof_values(shape_func);
+       second_derivs[point][fe->system_to_component_index(shape_func).first]
+         += tmp;
+      };
+};
+
+
+
 template <int dim>
 const Point<dim> &
 FEValuesBase<dim>::quadrature_point (const unsigned int i) const
@@ -929,3 +971,14 @@ void FEValuesBase<deal_II_dimension>::get_function_2nd_derivatives (const Vector
 template
 void FEValuesBase<deal_II_dimension>::get_function_2nd_derivatives (const BlockVector<double> &,
                                             vector<Tensor<2,deal_II_dimension> > &) const;
+//-----------------------------------------------------------------------------
+
+template
+void FEValuesBase<deal_II_dimension>::get_function_2nd_derivatives (const Vector<double> &,
+                                           vector<vector<Tensor<2,deal_II_dimension> > > &) const;
+template
+void FEValuesBase<deal_II_dimension>::get_function_2nd_derivatives (const Vector<float> &,
+                                            vector<vector<Tensor<2,deal_II_dimension> > > &) const;
+template
+void FEValuesBase<deal_II_dimension>::get_function_2nd_derivatives (const BlockVector<double> &,
+                                            vector<vector<Tensor<2,deal_II_dimension> > > &) const;

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.