$extrastylesheet
fem_context.C
Go to the documentation of this file.
00001 // The libMesh Finite Element Library.
00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
00003 
00004 // This library is free software; you can redistribute it and/or
00005 // modify it under the terms of the GNU Lesser General Public
00006 // License as published by the Free Software Foundation; either
00007 // version 2.1 of the License, or (at your option) any later version.
00008 
00009 // This library is distributed in the hope that it will be useful,
00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00012 // Lesser General Public License for more details.
00013 
00014 // You should have received a copy of the GNU Lesser General Public
00015 // License along with this library; if not, write to the Free Software
00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00017 
00018 
00019 
00020 #include "libmesh/boundary_info.h"
00021 #include "libmesh/dof_map.h"
00022 #include "libmesh/elem.h"
00023 #include "libmesh/fe_base.h"
00024 #include "libmesh/fe_interface.h"
00025 #include "libmesh/fem_context.h"
00026 #include "libmesh/libmesh_logging.h"
00027 #include "libmesh/mesh_base.h"
00028 #include "libmesh/quadrature.h"
00029 #include "libmesh/system.h"
00030 #include "libmesh/time_solver.h"
00031 #include "libmesh/unsteady_solver.h" // For euler_residual
00032 
00033 namespace libMesh
00034 {
00035 
00036 FEMContext::FEMContext (const System &sys)
00037   : DiffContext(sys),
00038     _mesh_sys(NULL),
00039     _mesh_x_var(0),
00040     _mesh_y_var(0),
00041     _mesh_z_var(0),
00042     side(0), edge(0),
00043     _boundary_info(sys.get_mesh().get_boundary_info()),
00044     _elem(NULL),
00045     _dim(sys.get_mesh().mesh_dimension()),
00046     _elem_dim(0), /* This will be reset in set_elem(). */
00047     _element_qrule(4,NULL),
00048     _side_qrule(4,NULL),
00049     _edge_qrule(NULL)
00050 {
00051   // Reserve space for the FEAbstract and QBase objects for each
00052   // element dimension possiblity (0,1,2,3)
00053   _element_fe.resize(4);
00054   _side_fe.resize(4);
00055   _element_fe_var.resize(4);
00056   _side_fe_var.resize(4);
00057 
00058   // We need to know which of our variables has the hardest
00059   // shape functions to numerically integrate.
00060 
00061   unsigned int nv = sys.n_vars();
00062 
00063   libmesh_assert (nv);
00064   FEType hardest_fe_type = sys.variable_type(0);
00065 
00066   bool have_scalar = false;
00067 
00068   for (unsigned int i=0; i != nv; ++i)
00069     {
00070       FEType fe_type = sys.variable_type(i);
00071 
00072       // FIXME - we don't yet handle mixed finite elements from
00073       // different families which require different quadrature rules
00074       // libmesh_assert_equal_to (fe_type.family, hardest_fe_type.family);
00075 
00076       if (fe_type.order > hardest_fe_type.order)
00077         hardest_fe_type = fe_type;
00078 
00079       // We need to detect SCALAR's so we can prepare FE objects for them.
00080       if( fe_type.family == SCALAR )
00081         have_scalar = true;
00082     }
00083 
00084   std::set<unsigned char> elem_dims( sys.get_mesh().elem_dimensions() );
00085   if(have_scalar)
00086     // SCALAR FEs have dimension 0 by assumption
00087     elem_dims.insert(0);
00088 
00089   for( std::set<unsigned char>::const_iterator dim_it = elem_dims.begin();
00090        dim_it != elem_dims.end(); ++dim_it )
00091     {
00092       const unsigned char dim = *dim_it;
00093 
00094       // Create an adequate quadrature rule
00095       _element_qrule[dim] = hardest_fe_type.default_quadrature_rule
00096         (dim, sys.extra_quadrature_order).release();
00097       _side_qrule[dim] = hardest_fe_type.default_quadrature_rule
00098         (dim-1, sys.extra_quadrature_order).release();
00099       if (dim == 3)
00100         _edge_qrule = hardest_fe_type.default_quadrature_rule
00101           (1, sys.extra_quadrature_order).release();
00102 
00103       // Next, create finite element objects
00104       _element_fe_var[dim].resize(nv);
00105       _side_fe_var[dim].resize(nv);
00106       if (dim == 3)
00107         _edge_fe_var.resize(nv);
00108 
00109 
00110       for (unsigned int i=0; i != nv; ++i)
00111         {
00112           FEType fe_type = sys.variable_type(i);
00113 
00114           if ( _element_fe[dim][fe_type] == NULL )
00115             {
00116               _element_fe[dim][fe_type] = FEAbstract::build(dim, fe_type).release();
00117               _element_fe[dim][fe_type]->attach_quadrature_rule(_element_qrule[dim]);
00118               _side_fe[dim][fe_type] = FEAbstract::build(dim, fe_type).release();
00119               _side_fe[dim][fe_type]->attach_quadrature_rule(_side_qrule[dim]);
00120 
00121               if (this->_dim == 3)
00122                 {
00123                   _edge_fe[fe_type] = FEAbstract::build(dim, fe_type).release();
00124                   _edge_fe[fe_type]->attach_quadrature_rule(_edge_qrule);
00125                 }
00126             }
00127 
00128           _element_fe_var[dim][i] = _element_fe[dim][fe_type];
00129           _side_fe_var[dim][i] = _side_fe[dim][fe_type];
00130           if ((dim) == 3)
00131             _edge_fe_var[i] = _edge_fe[fe_type];
00132 
00133         }
00134     }
00135 }
00136 
00137 FEMContext::~FEMContext()
00138 {
00139   // We don't want to store UniquePtrs in STL containers, but we don't
00140   // want to leak memory either
00141   for (std::vector<std::map<FEType, FEAbstract *> >::iterator d = _element_fe.begin();
00142        d != _element_fe.end(); ++d)
00143     for (std::map<FEType, FEAbstract *>::iterator i = d->begin();
00144          i != d->end(); ++i)
00145       delete i->second;
00146 
00147   for (std::vector<std::map<FEType, FEAbstract *> >::iterator d = _side_fe.begin();
00148        d != _side_fe.end(); ++d)
00149     for (std::map<FEType, FEAbstract *>::iterator i = d->begin();
00150          i != d->end(); ++i)
00151       delete i->second;
00152 
00153   for (std::map<FEType, FEAbstract *>::iterator i = _edge_fe.begin();
00154        i != _edge_fe.end(); ++i)
00155     delete i->second;
00156   _edge_fe.clear();
00157 
00158   for (std::vector<QBase*>::iterator i = _element_qrule.begin();
00159        i != _element_qrule.end(); ++i)
00160     delete *i;
00161   _element_qrule.clear();
00162 
00163   for (std::vector<QBase*>::iterator i = _side_qrule.begin();
00164        i != _side_qrule.end(); ++i)
00165     delete *i;
00166   _side_qrule.clear();
00167 
00168   delete _edge_qrule;
00169   _edge_qrule = NULL;
00170 }
00171 
00172 
00173 
00174 bool FEMContext::has_side_boundary_id(boundary_id_type id) const
00175 {
00176   return _boundary_info.has_boundary_id(&(this->get_elem()), side, id);
00177 }
00178 
00179 
00180 std::vector<boundary_id_type> FEMContext::side_boundary_ids() const
00181 {
00182   return _boundary_info.boundary_ids(&(this->get_elem()), side);
00183 }
00184 
00185 
00186 
00187 template<typename OutputType,
00188          FEMContext::diff_subsolution_getter subsolution_getter>
00189 void FEMContext::some_interior_value(unsigned int var, unsigned int qp, OutputType& u) const
00190 {
00191   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
00192 
00193   // Get local-to-global dof index lookup
00194   libmesh_assert_greater (this->get_dof_indices().size(), var);
00195   const unsigned int n_dofs = cast_int<unsigned int>
00196     (this->get_dof_indices(var).size());
00197 
00198   // Get current local coefficients
00199   const DenseSubVector<Number> &coef = (this->*subsolution_getter)(var);
00200 
00201   // Get finite element object
00202   FEGenericBase<OutputShape>* fe = NULL;
00203   this->get_element_fe<OutputShape>( var, fe );
00204 
00205   // Get shape function values at quadrature point
00206   const std::vector<std::vector<OutputShape> > &phi = fe->get_phi();
00207 
00208   // Accumulate solution value
00209   u = 0.;
00210 
00211   for (unsigned int l=0; l != n_dofs; l++)
00212     u += phi[l][qp] * coef(l);
00213 }
00214 
00215 
00216 
00217 template<typename OutputType,
00218          FEMContext::diff_subsolution_getter subsolution_getter>
00219 void FEMContext::some_interior_gradient(unsigned int var, unsigned int qp, OutputType& du) const
00220 {
00221   typedef typename TensorTools::MakeReal
00222     <typename TensorTools::DecrementRank<OutputType>::type>::type
00223     OutputShape;
00224 
00225   // Get local-to-global dof index lookup
00226   libmesh_assert_greater (this->get_dof_indices().size(), var);
00227   const unsigned int n_dofs = cast_int<unsigned int>
00228     (this->get_dof_indices(var).size());
00229 
00230   // Get current local coefficients
00231   const DenseSubVector<Number> &coef = (this->*subsolution_getter)(var);
00232 
00233   // Get finite element object
00234   FEGenericBase<OutputShape>* fe = NULL;
00235   this->get_element_fe<OutputShape>( var, fe );
00236 
00237   // Get shape function values at quadrature point
00238   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> > &dphi = fe->get_dphi();
00239 
00240   // Accumulate solution derivatives
00241   du = 0;
00242 
00243   for (unsigned int l=0; l != n_dofs; l++)
00244     du.add_scaled(dphi[l][qp], coef(l));
00245 
00246   return;
00247 }
00248 
00249 
00250 
00251 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
00252 template<typename OutputType,
00253          FEMContext::diff_subsolution_getter subsolution_getter>
00254 void FEMContext::some_interior_hessian(unsigned int var, unsigned int qp, OutputType& d2u) const
00255 {
00256   typedef typename TensorTools::MakeReal<
00257     typename TensorTools::DecrementRank<
00258     typename TensorTools::DecrementRank<
00259     OutputType>::type>::type>::type
00260     OutputShape;
00261 
00262   // Get local-to-global dof index lookup
00263   libmesh_assert_greater (this->get_dof_indices().size(), var);
00264   const unsigned int n_dofs = cast_int<unsigned int>
00265     (this->get_dof_indices(var).size());
00266 
00267   // Get current local coefficients
00268   const DenseSubVector<Number> &coef = (this->*subsolution_getter)(var);
00269 
00270   // Get finite element object
00271   FEGenericBase<OutputShape>* fe = NULL;
00272   this->get_element_fe<OutputShape>( var, fe );
00273 
00274   // Get shape function values at quadrature point
00275   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> > &d2phi = fe->get_d2phi();
00276 
00277   // Accumulate solution second derivatives
00278   d2u = 0.0;
00279 
00280   for (unsigned int l=0; l != n_dofs; l++)
00281     d2u.add_scaled(d2phi[l][qp], coef(l));
00282 
00283   return;
00284 }
00285 #endif
00286 
00287 
00288 
00289 template<typename OutputType,
00290          FEMContext::diff_subsolution_getter subsolution_getter>
00291 void FEMContext::some_side_value(unsigned int var, unsigned int qp, OutputType& u) const
00292 {
00293   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
00294 
00295   // Get local-to-global dof index lookup
00296   libmesh_assert_greater (this->get_dof_indices().size(), var);
00297   const unsigned int n_dofs = cast_int<unsigned int>
00298     (this->get_dof_indices(var).size());
00299 
00300   // Get current local coefficients
00301   const DenseSubVector<Number> &coef = (this->*subsolution_getter)(var);
00302 
00303   // Get finite element object
00304   FEGenericBase<OutputShape>* the_side_fe = NULL;
00305   this->get_side_fe<OutputShape>( var, the_side_fe );
00306 
00307   // Get shape function values at quadrature point
00308   const std::vector<std::vector<OutputShape> > &phi = the_side_fe->get_phi();
00309 
00310   // Accumulate solution value
00311   u = 0.;
00312 
00313   for (unsigned int l=0; l != n_dofs; l++)
00314     u += phi[l][qp] * coef(l);
00315 
00316   return;
00317 }
00318 
00319 Number FEMContext::interior_value(unsigned int var, unsigned int qp) const
00320 {
00321   Number u;
00322 
00323   this->interior_value( var, qp, u );
00324 
00325   return u;
00326 }
00327 
00328 template<typename OutputType>
00329 void FEMContext::interior_value(unsigned int var, unsigned int qp,
00330                                 OutputType& u) const
00331 {
00332   this->some_interior_value <OutputType,&DiffContext::get_elem_solution>
00333     (var, qp, u);
00334 }
00335 
00336 
00337 template<typename OutputType>
00338 void FEMContext::interior_values (unsigned int var,
00339                                   const NumericVector<Number> & _system_vector,
00340                                   std::vector<OutputType>& u_vals) const
00341 {
00342   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
00343 
00344   // Get local-to-global dof index lookup
00345   libmesh_assert_greater (this->get_dof_indices().size(), var);
00346   const unsigned int n_dofs = cast_int<unsigned int>
00347     (this->get_dof_indices(var).size());
00348 
00349   // Get current local coefficients
00350   const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);
00351 
00352   // Get the finite element object
00353   FEGenericBase<OutputShape>* fe = NULL;
00354   this->get_element_fe<OutputShape>( var, fe );
00355 
00356   // Get shape function values at quadrature point
00357   const std::vector<std::vector<OutputShape> > &phi = fe->get_phi();
00358 
00359   // Loop over all the q_points on this element
00360   for (unsigned int qp=0; qp != u_vals.size(); qp++)
00361     {
00362       OutputType &u = u_vals[qp];
00363 
00364       // Compute the value at this q_point
00365       u = 0.;
00366 
00367       for (unsigned int l=0; l != n_dofs; l++)
00368         u += phi[l][qp] * coef(l);
00369     }
00370 
00371   return;
00372 }
00373 
00374 Gradient FEMContext::interior_gradient(unsigned int var, unsigned int qp) const
00375 {
00376   Gradient du;
00377 
00378   this->interior_gradient( var, qp, du );
00379 
00380   return du;
00381 }
00382 
00383 
00384 
00385 template<typename OutputType>
00386 void FEMContext::interior_gradient(unsigned int var, unsigned int qp,
00387                                    OutputType& du) const
00388 {
00389   this->some_interior_gradient <OutputType,&DiffContext::get_elem_solution>
00390     (var, qp, du);
00391 }
00392 
00393 
00394 
00395 template<typename OutputType>
00396 void FEMContext::interior_gradients
00397 (unsigned int var,
00398  const NumericVector<Number> & _system_vector,
00399  std::vector<OutputType>& du_vals) const
00400 {
00401   typedef typename TensorTools::MakeReal
00402     <typename TensorTools::DecrementRank<OutputType>::type>::type
00403     OutputShape;
00404 
00405   // Get local-to-global dof index lookup
00406   libmesh_assert_greater (this->get_dof_indices().size(), var);
00407   const unsigned int n_dofs = cast_int<unsigned int>
00408     (this->get_dof_indices(var).size());
00409 
00410   // Get current local coefficients
00411   const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);
00412 
00413   // Get finite element object
00414   FEGenericBase<OutputShape>* fe = NULL;
00415   this->get_element_fe<OutputShape>( var, fe );
00416 
00417   // Get shape function values at quadrature point
00418   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> > &dphi = fe->get_dphi();
00419 
00420   // Loop over all the q_points in this finite element
00421   for (unsigned int qp=0; qp != du_vals.size(); qp++)
00422     {
00423       OutputType &du = du_vals[qp];
00424 
00425       // Compute the gradient at this q_point
00426       du = 0;
00427 
00428       for (unsigned int l=0; l != n_dofs; l++)
00429         du.add_scaled(dphi[l][qp], coef(l));
00430     }
00431 
00432   return;
00433 }
00434 
00435 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
00436 Tensor FEMContext::interior_hessian(unsigned int var, unsigned int qp) const
00437 {
00438   Tensor d2u;
00439 
00440   this->interior_hessian( var, qp, d2u );
00441 
00442   return d2u;
00443 }
00444 
00445 template<typename OutputType>
00446 void FEMContext::interior_hessian(unsigned int var, unsigned int qp,
00447                                   OutputType& d2u) const
00448 {
00449   this->some_interior_hessian <OutputType,&DiffContext::get_elem_solution>
00450     (var, qp, d2u);
00451 }
00452 
00453 
00454 template<typename OutputType>
00455 void FEMContext::interior_hessians
00456 (unsigned int var,
00457  const NumericVector<Number> & _system_vector,
00458  std::vector<OutputType>& d2u_vals) const
00459 {
00460   typedef typename TensorTools::MakeReal<
00461     typename TensorTools::DecrementRank<
00462     typename TensorTools::DecrementRank<
00463     OutputType>::type>::type>::type
00464     OutputShape;
00465 
00466   // Get local-to-global dof index lookup
00467   libmesh_assert_greater (this->get_dof_indices().size(), var);
00468   const unsigned int n_dofs = cast_int<unsigned int>
00469     (this->get_dof_indices(var).size());
00470 
00471   // Get current local coefficients
00472   const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);
00473 
00474   // Get finite element object
00475   FEGenericBase<OutputShape>* fe = NULL;
00476   this->get_element_fe<OutputShape>( var, fe );
00477 
00478   // Get shape function values at quadrature point
00479   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> > &d2phi = fe->get_d2phi();
00480 
00481   // Loop over all the q_points in this finite element
00482   for (unsigned int qp=0; qp != d2u_vals.size(); qp++)
00483     {
00484       OutputType &d2u = d2u_vals[qp];
00485 
00486       // Compute the gradient at this q_point
00487       d2u = 0;
00488 
00489       for (unsigned int l=0; l != n_dofs; l++)
00490         d2u.add_scaled(d2phi[l][qp], coef(l));
00491     }
00492 
00493   return;
00494 }
00495 
00496 
00497 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
00498 
00499 
00500 template<typename OutputType>
00501 void FEMContext::interior_curl(unsigned int var, unsigned int qp,
00502                                OutputType& curl_u) const
00503 {
00504   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
00505 
00506   // Get local-to-global dof index lookup
00507   libmesh_assert_greater (this->get_dof_indices().size(), var);
00508   const unsigned int n_dofs = cast_int<unsigned int>
00509     (this->get_dof_indices(var).size());
00510 
00511   // Get current local coefficients
00512   libmesh_assert_greater (this->_elem_subsolutions.size(), var);
00513   libmesh_assert(&(this->get_elem_solution(var)));
00514   const DenseSubVector<Number> &coef = this->get_elem_solution(var);
00515 
00516   // Get finite element object
00517   FEGenericBase<OutputShape>* fe = NULL;
00518   this->get_element_fe<OutputShape>( var, fe );
00519 
00520   // Get shape function values at quadrature point
00521   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputShape> > &curl_phi = fe->get_curl_phi();
00522 
00523   // Accumulate solution curl
00524   curl_u = 0.;
00525 
00526   for (unsigned int l=0; l != n_dofs; l++)
00527     curl_u.add_scaled(curl_phi[l][qp], coef(l));
00528 
00529   return;
00530 }
00531 
00532 
00533 template<typename OutputType>
00534 void FEMContext::interior_div(unsigned int var, unsigned int qp,
00535                               OutputType& div_u) const
00536 {
00537   typedef typename
00538     TensorTools::IncrementRank
00539     <typename TensorTools::MakeReal<OutputType>::type>::type OutputShape;
00540 
00541   // Get local-to-global dof index lookup
00542   libmesh_assert_greater (this->get_dof_indices().size(), var);
00543   const unsigned int n_dofs = cast_int<unsigned int>
00544     (this->get_dof_indices(var).size());
00545 
00546   // Get current local coefficients
00547   libmesh_assert_greater (this->_elem_subsolutions.size(), var);
00548   libmesh_assert(&(this->get_elem_solution(var)));
00549   const DenseSubVector<Number> &coef = this->get_elem_solution(var);
00550 
00551   // Get finite element object
00552   FEGenericBase<OutputShape>* fe = NULL;
00553   this->get_element_fe<OutputShape>( var, fe );
00554 
00555   // Get shape function values at quadrature point
00556   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputDivergence> > &div_phi = fe->get_div_phi();
00557 
00558   // Accumulate solution curl
00559   div_u = 0.;
00560 
00561   for (unsigned int l=0; l != n_dofs; l++)
00562     div_u += div_phi[l][qp] * coef(l);
00563 
00564   return;
00565 }
00566 
00567 
00568 Number FEMContext::side_value(unsigned int var, unsigned int qp) const
00569 {
00570   Number u = 0.;
00571 
00572   this->side_value( var, qp, u );
00573 
00574   return u;
00575 }
00576 
00577 
00578 template<typename OutputType>
00579 void FEMContext::side_value(unsigned int var, unsigned int qp,
00580                             OutputType& u) const
00581 {
00582   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
00583 
00584   // Get local-to-global dof index lookup
00585   libmesh_assert_greater (this->get_dof_indices().size(), var);
00586   const unsigned int n_dofs = cast_int<unsigned int>
00587     (this->get_dof_indices(var).size());
00588 
00589   // Get current local coefficients
00590   libmesh_assert_greater (this->_elem_subsolutions.size(), var);
00591   libmesh_assert(&(this->get_elem_solution(var)));
00592   const DenseSubVector<Number> &coef = this->get_elem_solution(var);
00593 
00594   // Get finite element object
00595   FEGenericBase<OutputShape>* the_side_fe = NULL;
00596   this->get_side_fe<OutputShape>( var, the_side_fe );
00597 
00598   // Get shape function values at quadrature point
00599   const std::vector<std::vector<OutputShape> > &phi = the_side_fe->get_phi();
00600 
00601   // Accumulate solution value
00602   u = 0.;
00603 
00604   for (unsigned int l=0; l != n_dofs; l++)
00605     u += phi[l][qp] * coef(l);
00606 
00607   return;
00608 }
00609 
00610 
00611 template<typename OutputType>
00612 void FEMContext::side_values
00613 (unsigned int var,
00614  const NumericVector<Number> & _system_vector,
00615  std::vector<OutputType>& u_vals) const
00616 {
00617   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
00618 
00619   // Get local-to-global dof index lookup
00620   libmesh_assert_greater (this->get_dof_indices().size(), var);
00621   const unsigned int n_dofs = cast_int<unsigned int>
00622     (this->get_dof_indices(var).size());
00623 
00624   // Get current local coefficients
00625   const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);
00626 
00627   // Get the finite element object
00628   FEGenericBase<OutputShape>* the_side_fe = NULL;
00629   this->get_side_fe<OutputShape>( var, the_side_fe );
00630 
00631   // Get shape function values at quadrature point
00632   const std::vector<std::vector<OutputShape> > &phi = the_side_fe->get_phi();
00633 
00634   // Loop over all the q_points on this element
00635   for (unsigned int qp=0; qp != u_vals.size(); qp++)
00636     {
00637       OutputType &u = u_vals[qp];
00638 
00639       // Compute the value at this q_point
00640       u = 0.;
00641 
00642       for (unsigned int l=0; l != n_dofs; l++)
00643         u += phi[l][qp] * coef(l);
00644     }
00645 
00646   return;
00647 }
00648 
00649 Gradient FEMContext::side_gradient(unsigned int var, unsigned int qp) const
00650 {
00651   Gradient du;
00652 
00653   this->side_gradient( var, qp, du );
00654 
00655   return du;
00656 }
00657 
00658 
00659 template<typename OutputType>
00660 void FEMContext::side_gradient(unsigned int var, unsigned int qp,
00661                                OutputType& du) const
00662 {
00663   typedef typename TensorTools::MakeReal
00664     <typename TensorTools::DecrementRank<OutputType>::type>::type
00665     OutputShape;
00666 
00667   // Get local-to-global dof index lookup
00668   libmesh_assert_greater (this->get_dof_indices().size(), var);
00669   const unsigned int n_dofs = cast_int<unsigned int>
00670     (this->get_dof_indices(var).size());
00671 
00672   // Get current local coefficients
00673   libmesh_assert_greater (this->_elem_subsolutions.size(), var);
00674   libmesh_assert(&(this->get_elem_solution(var)));
00675   const DenseSubVector<Number> &coef = this->get_elem_solution(var);
00676 
00677   // Get finite element object
00678   FEGenericBase<OutputShape>* the_side_fe = NULL;
00679   this->get_side_fe<OutputShape>( var, the_side_fe );
00680 
00681   // Get shape function values at quadrature point
00682   const std::vector<std::vector< typename FEGenericBase<OutputShape>::OutputGradient> > &dphi = the_side_fe->get_dphi();
00683 
00684   // Accumulate solution derivatives
00685   du = 0.;
00686 
00687   for (unsigned int l=0; l != n_dofs; l++)
00688     du.add_scaled(dphi[l][qp], coef(l));
00689 
00690   return;
00691 }
00692 
00693 
00694 
00695 template<typename OutputType>
00696 void FEMContext::side_gradients
00697 (unsigned int var,
00698  const NumericVector<Number> & _system_vector,
00699  std::vector<OutputType>& du_vals) const
00700 {
00701   typedef typename TensorTools::MakeReal
00702     <typename TensorTools::DecrementRank<OutputType>::type>::type
00703     OutputShape;
00704 
00705   // Get local-to-global dof index lookup
00706   libmesh_assert_greater (this->get_dof_indices().size(), var);
00707   const unsigned int n_dofs = cast_int<unsigned int>
00708     (this->get_dof_indices(var).size());
00709 
00710   // Get current local coefficients
00711   const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);
00712 
00713   // Get finite element object
00714   FEGenericBase<OutputShape>* the_side_fe = NULL;
00715   this->get_side_fe<OutputShape>( var, the_side_fe );
00716 
00717   // Get shape function values at quadrature point
00718   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> > &dphi = the_side_fe->get_dphi();
00719 
00720   // Loop over all the q_points in this finite element
00721   for (unsigned int qp=0; qp != du_vals.size(); qp++)
00722     {
00723       OutputType &du = du_vals[qp];
00724 
00725       du = 0;
00726 
00727       // Compute the gradient at this q_point
00728       for (unsigned int l=0; l != n_dofs; l++)
00729         du.add_scaled(dphi[l][qp], coef(l));
00730     }
00731 
00732   return;
00733 }
00734 
00735 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
00736 Tensor FEMContext::side_hessian(unsigned int var, unsigned int qp) const
00737 {
00738   Tensor d2u;
00739 
00740   this->side_hessian( var, qp, d2u );
00741 
00742   return d2u;
00743 }
00744 
00745 template<typename OutputType>
00746 void FEMContext::side_hessian(unsigned int var, unsigned int qp,
00747                               OutputType& d2u) const
00748 {
00749   typedef typename TensorTools::MakeReal<
00750     typename TensorTools::DecrementRank<
00751     typename TensorTools::DecrementRank<
00752     OutputType>::type>::type>::type
00753     OutputShape;
00754 
00755   // Get local-to-global dof index lookup
00756   libmesh_assert_greater (this->get_dof_indices().size(), var);
00757   const unsigned int n_dofs = cast_int<unsigned int>
00758     (this->get_dof_indices(var).size());
00759 
00760   // Get current local coefficients
00761   libmesh_assert_greater (this->_elem_subsolutions.size(), var);
00762   libmesh_assert(&(this->get_elem_solution(var)));
00763   const DenseSubVector<Number> &coef = this->get_elem_solution(var);
00764 
00765   // Get finite element object
00766   FEGenericBase<OutputShape>* the_side_fe = NULL;
00767   this->get_side_fe<OutputShape>( var, the_side_fe );
00768 
00769   // Get shape function values at quadrature point
00770   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> > &d2phi = the_side_fe->get_d2phi();
00771 
00772   // Accumulate solution second derivatives
00773   d2u = 0.0;
00774 
00775   for (unsigned int l=0; l != n_dofs; l++)
00776     d2u.add_scaled(d2phi[l][qp], coef(l));
00777 
00778   return;
00779 }
00780 
00781 
00782 template<typename OutputType>
00783 void FEMContext::side_hessians
00784 (unsigned int var,
00785  const NumericVector<Number> & _system_vector,
00786  std::vector<OutputType>& d2u_vals) const
00787 {
00788   typedef typename TensorTools::MakeReal<
00789     typename TensorTools::DecrementRank<
00790     typename TensorTools::DecrementRank<
00791     OutputType>::type>::type>::type
00792     OutputShape;
00793 
00794   // Get local-to-global dof index lookup
00795   libmesh_assert_greater (this->get_dof_indices().size(), var);
00796   const unsigned int n_dofs = cast_int<unsigned int>
00797     (this->get_dof_indices(var).size());
00798 
00799   // Get current local coefficients
00800   const DenseSubVector<Number> &coef = get_localized_subvector(_system_vector, var);
00801 
00802   // Get finite element object
00803   FEGenericBase<OutputShape>* the_side_fe = NULL;
00804   this->get_side_fe<OutputShape>( var, the_side_fe );
00805 
00806   // Get shape function values at quadrature point
00807   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> > &d2phi = the_side_fe->get_d2phi();
00808 
00809   // Loop over all the q_points in this finite element
00810   for (unsigned int qp=0; qp != d2u_vals.size(); qp++)
00811     {
00812       OutputType &d2u = d2u_vals[qp];
00813 
00814       // Compute the gradient at this q_point
00815       d2u = 0;
00816 
00817       for (unsigned int l=0; l != n_dofs; l++)
00818         d2u.add_scaled(d2phi[l][qp], coef(l));
00819     }
00820 
00821   return;
00822 }
00823 
00824 
00825 
00826 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
00827 
00828 
00829 
00830 Number FEMContext::point_value(unsigned int var, const Point &p) const
00831 {
00832   Number u = 0.;
00833 
00834   this->point_value( var, p, u );
00835 
00836   return u;
00837 }
00838 
00839 template<typename OutputType>
00840 void FEMContext::point_value(unsigned int var, const Point &p,
00841                              OutputType& u) const
00842 {
00843   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
00844 
00845   // Get local-to-global dof index lookup
00846   libmesh_assert_greater (this->get_dof_indices().size(), var);
00847   const unsigned int n_dofs = cast_int<unsigned int>
00848     (this->get_dof_indices(var).size());
00849 
00850   // Get current local coefficients
00851   libmesh_assert_greater (this->_elem_subsolutions.size(), var);
00852   libmesh_assert(&(this->get_elem_solution(var)));
00853   const DenseSubVector<Number> &coef = this->get_elem_solution(var);
00854 
00855   // Get finite element object
00856   FEGenericBase<OutputShape>* fe = NULL;
00857   this->get_element_fe<OutputShape>( var, fe );
00858 
00859   // Build a FE for calculating u(p)
00860   UniquePtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );
00861 
00862   // Get the values of the shape function derivatives
00863   const std::vector<std::vector<OutputShape> >&  phi = fe_new->get_phi();
00864 
00865   u = 0.;
00866 
00867   for (unsigned int l=0; l != n_dofs; l++)
00868     u += phi[l][0] * coef(l);
00869 
00870   return;
00871 }
00872 
00873 
00874 
00875 Gradient FEMContext::point_gradient(unsigned int var, const Point &p) const
00876 {
00877   Gradient grad_u;
00878 
00879   this->point_gradient( var, p, grad_u );
00880 
00881   return grad_u;
00882 }
00883 
00884 
00885 
00886 template<typename OutputType>
00887 void FEMContext::point_gradient(unsigned int var, const Point &p,
00888                                 OutputType& grad_u) const
00889 {
00890   typedef typename TensorTools::MakeReal
00891     <typename TensorTools::DecrementRank<OutputType>::type>::type
00892     OutputShape;
00893 
00894   // Get local-to-global dof index lookup
00895   libmesh_assert_greater (this->get_dof_indices().size(), var);
00896   const unsigned int n_dofs = cast_int<unsigned int>
00897     (this->get_dof_indices(var).size());
00898 
00899   // Get current local coefficients
00900   libmesh_assert_greater (this->_elem_subsolutions.size(), var);
00901   libmesh_assert(&(this->get_elem_solution(var)));
00902   const DenseSubVector<Number> &coef = this->get_elem_solution(var);
00903 
00904   // Get finite element object
00905   FEGenericBase<OutputShape>* fe = NULL;
00906   this->get_element_fe<OutputShape>( var, fe );
00907 
00908   // Build a FE for calculating u(p)
00909   UniquePtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );
00910 
00911   // Get the values of the shape function derivatives
00912   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> >&  dphi = fe_new->get_dphi();
00913 
00914   grad_u = 0.0;
00915 
00916   for (unsigned int l=0; l != n_dofs; l++)
00917     grad_u.add_scaled(dphi[l][0], coef(l));
00918 
00919   return;
00920 }
00921 
00922 
00923 
00924 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
00925 
00926 Tensor FEMContext::point_hessian(unsigned int var, const Point &p) const
00927 {
00928   Tensor hess_u;
00929 
00930   this->point_hessian( var, p, hess_u );
00931 
00932   return hess_u;
00933 }
00934 
00935 
00936 template<typename OutputType>
00937 void FEMContext::point_hessian(unsigned int var, const Point &p,
00938                                OutputType& hess_u) const
00939 {
00940   typedef typename TensorTools::MakeReal<
00941     typename TensorTools::DecrementRank<
00942     typename TensorTools::DecrementRank<
00943     OutputType>::type>::type>::type
00944     OutputShape;
00945 
00946   // Get local-to-global dof index lookup
00947   libmesh_assert_greater (this->get_dof_indices().size(), var);
00948   const unsigned int n_dofs = cast_int<unsigned int>
00949     (this->get_dof_indices(var).size());
00950 
00951   // Get current local coefficients
00952   libmesh_assert_greater (this->_elem_subsolutions.size(), var);
00953   libmesh_assert(&(this->get_elem_solution(var)));
00954   const DenseSubVector<Number> &coef = this->get_elem_solution(var);
00955 
00956   // Get finite element object
00957   FEGenericBase<OutputShape>* fe = NULL;
00958   this->get_element_fe<OutputShape>( var, fe );
00959 
00960   // Build a FE for calculating u(p)
00961   UniquePtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );
00962 
00963   // Get the values of the shape function derivatives
00964   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> >&  d2phi = fe_new->get_d2phi();
00965 
00966   hess_u = 0.0;
00967 
00968   for (unsigned int l=0; l != n_dofs; l++)
00969     hess_u.add_scaled(d2phi[l][0], coef(l));
00970 
00971   return;
00972 }
00973 
00974 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
00975 
00976 
00977 template<typename OutputType>
00978 void FEMContext::point_curl(unsigned int var, const Point &p,
00979                             OutputType& curl_u) const
00980 {
00981   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
00982 
00983   // Get local-to-global dof index lookup
00984   libmesh_assert_greater (this->get_dof_indices().size(), var);
00985   const unsigned int n_dofs = cast_int<unsigned int>
00986     (this->get_dof_indices(var).size());
00987 
00988   // Get current local coefficients
00989   libmesh_assert_greater (this->_elem_subsolutions.size(), var);
00990   libmesh_assert(&(this->get_elem_solution(var)));
00991   const DenseSubVector<Number> &coef = this->get_elem_solution(var);
00992 
00993   // Get finite element object
00994   FEGenericBase<OutputShape>* fe = NULL;
00995   this->get_element_fe<OutputShape>( var, fe );
00996 
00997   // Build a FE for calculating u(p)
00998   UniquePtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );
00999 
01000   // Get the values of the shape function derivatives
01001   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputShape> >&  curl_phi = fe_new->get_curl_phi();
01002 
01003   curl_u = 0.0;
01004 
01005   for (unsigned int l=0; l != n_dofs; l++)
01006     curl_u.add_scaled(curl_phi[l][0], coef(l));
01007 
01008   return;
01009 }
01010 
01011 
01012 
01013 Number FEMContext::fixed_interior_value(unsigned int var, unsigned int qp) const
01014 {
01015   Number u = 0.;
01016 
01017   this->fixed_interior_value( var, qp, u );
01018 
01019   return u;
01020 }
01021 
01022 
01023 
01024 template<typename OutputType>
01025 void FEMContext::fixed_interior_value(unsigned int var, unsigned int qp,
01026                                       OutputType& u) const
01027 {
01028   this->some_interior_value <OutputType,&DiffContext::get_elem_fixed_solution>
01029     (var, qp, u);
01030 }
01031 
01032 
01033 
01034 Gradient FEMContext::fixed_interior_gradient(unsigned int var, unsigned int qp) const
01035 {
01036   Gradient du;
01037 
01038   this->fixed_interior_gradient( var, qp, du );
01039 
01040   return du;
01041 }
01042 
01043 
01044 template<typename OutputType>
01045 void FEMContext::FEMContext::fixed_interior_gradient(unsigned int var, unsigned int qp,
01046                                                      OutputType& du) const
01047 {
01048   this->some_interior_gradient <OutputType,&DiffContext::get_elem_fixed_solution>
01049     (var, qp, du);
01050 }
01051 
01052 
01053 
01054 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
01055 Tensor FEMContext::fixed_interior_hessian(unsigned int var, unsigned int qp) const
01056 {
01057   Tensor d2u;
01058 
01059   this->fixed_interior_hessian( var, qp, d2u );
01060 
01061   return d2u;
01062 }
01063 
01064 
01065 template<typename OutputType>
01066 void FEMContext::fixed_interior_hessian(unsigned int var, unsigned int qp,
01067                                         OutputType& d2u) const
01068 {
01069   this->some_interior_hessian <OutputType,&DiffContext::get_elem_fixed_solution>
01070     (var, qp, d2u);
01071 }
01072 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
01073 
01074 
01075 
01076 Number FEMContext::fixed_side_value(unsigned int var, unsigned int qp) const
01077 {
01078   Number u = 0.;
01079 
01080   this->fixed_side_value( var, qp, u );
01081 
01082   return u;
01083 }
01084 
01085 
01086 template<typename OutputType>
01087 void FEMContext::fixed_side_value(unsigned int var, unsigned int qp,
01088                                   OutputType& u) const
01089 {
01090   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
01091 
01092   // Get local-to-global dof index lookup
01093   libmesh_assert_greater (this->get_dof_indices().size(), var);
01094   const unsigned int n_dofs = cast_int<unsigned int>
01095     (this->get_dof_indices(var).size());
01096 
01097   // Get current local coefficients
01098   libmesh_assert_greater (_elem_fixed_subsolutions.size(), var);
01099   libmesh_assert(&(this->get_elem_fixed_solution(var)));
01100   const DenseSubVector<Number> &coef = this->get_elem_fixed_solution(var);
01101 
01102   // Get finite element object
01103   FEGenericBase<OutputShape>* the_side_fe = NULL;
01104   this->get_side_fe<OutputShape>( var, the_side_fe );
01105 
01106   // Get shape function values at quadrature point
01107   const std::vector<std::vector<OutputShape> > &phi = the_side_fe->get_phi();
01108 
01109   // Accumulate solution value
01110   u = 0.0;
01111 
01112   for (unsigned int l=0; l != n_dofs; l++)
01113     u += phi[l][qp] * coef(l);
01114 
01115   return;
01116 }
01117 
01118 
01119 
01120 Gradient FEMContext::fixed_side_gradient(unsigned int var, unsigned int qp) const
01121 {
01122   Gradient du;
01123 
01124   this->fixed_side_gradient( var, qp, du );
01125 
01126   return du;
01127 }
01128 
01129 
01130 template<typename OutputType>
01131 void FEMContext::fixed_side_gradient(unsigned int var, unsigned int qp,
01132                                      OutputType& du) const
01133 {
01134   typedef typename TensorTools::MakeReal
01135     <typename TensorTools::DecrementRank<OutputType>::type>::type
01136     OutputShape;
01137 
01138   // Get local-to-global dof index lookup
01139   libmesh_assert_greater (this->get_dof_indices().size(), var);
01140   const unsigned int n_dofs = cast_int<unsigned int>
01141     (this->get_dof_indices(var).size());
01142 
01143   // Get current local coefficients
01144   libmesh_assert_greater (_elem_fixed_subsolutions.size(), var);
01145   libmesh_assert(&(this->get_elem_fixed_solution(var)));
01146   const DenseSubVector<Number> &coef = this->get_elem_fixed_solution(var);
01147 
01148   // Get finite element object
01149   FEGenericBase<OutputShape>* the_side_fe = NULL;
01150   this->get_side_fe<OutputShape>( var, the_side_fe );
01151 
01152   // Get shape function values at quadrature point
01153   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> > &dphi = the_side_fe->get_dphi();
01154 
01155   // Accumulate solution derivatives
01156   du = 0.0;
01157 
01158   for (unsigned int l=0; l != n_dofs; l++)
01159     du.add_scaled(dphi[l][qp], coef(l));
01160 
01161   return;
01162 }
01163 
01164 
01165 
01166 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
01167 Tensor FEMContext::fixed_side_hessian(unsigned int var, unsigned int qp) const
01168 {
01169   Tensor d2u;
01170 
01171   this->fixed_side_hessian( var, qp, d2u );
01172 
01173   return d2u;
01174 }
01175 
01176 template<typename OutputType>
01177 void FEMContext::fixed_side_hessian(unsigned int var, unsigned int qp,
01178                                     OutputType& d2u) const
01179 {
01180   typedef typename TensorTools::MakeReal<
01181     typename TensorTools::DecrementRank<
01182     typename TensorTools::DecrementRank<
01183     OutputType>::type>::type>::type
01184     OutputShape;
01185 
01186   // Get local-to-global dof index lookup
01187   libmesh_assert_greater (this->get_dof_indices().size(), var);
01188   const unsigned int n_dofs = cast_int<unsigned int>
01189     (this->get_dof_indices(var).size());
01190 
01191   // Get current local coefficients
01192   libmesh_assert_greater (_elem_fixed_subsolutions.size(), var);
01193   libmesh_assert(&(this->get_elem_fixed_solution(var)));
01194   const DenseSubVector<Number> &coef = this->get_elem_fixed_solution(var);
01195 
01196   // Get finite element object
01197   FEGenericBase<OutputShape>* the_side_fe = NULL;
01198   this->get_side_fe<OutputShape>( var, the_side_fe );
01199 
01200   // Get shape function values at quadrature point
01201   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> > &d2phi = the_side_fe->get_d2phi();
01202 
01203   // Accumulate solution second derivatives
01204   d2u = 0.0;
01205 
01206   for (unsigned int l=0; l != n_dofs; l++)
01207     d2u.add_scaled(d2phi[l][qp], coef(l));
01208 
01209   return;
01210 }
01211 #endif // ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
01212 
01213 
01214 
01215 Number FEMContext::fixed_point_value(unsigned int var, const Point &p) const
01216 {
01217   Number u = 0.;
01218 
01219   this->fixed_point_value( var, p, u );
01220 
01221   return u;
01222 }
01223 
01224 template<typename OutputType>
01225 void FEMContext::fixed_point_value(unsigned int var, const Point &p,
01226                                    OutputType& u) const
01227 {
01228   typedef typename TensorTools::MakeReal<OutputType>::type OutputShape;
01229 
01230   // Get local-to-global dof index lookup
01231   libmesh_assert_greater (this->get_dof_indices().size(), var);
01232   const unsigned int n_dofs = cast_int<unsigned int>
01233     (this->get_dof_indices(var).size());
01234 
01235   // Get current local coefficients
01236   libmesh_assert_greater (_elem_fixed_subsolutions.size(), var);
01237   libmesh_assert(&(this->get_elem_fixed_solution(var)));
01238   const DenseSubVector<Number> &coef = this->get_elem_fixed_solution(var);
01239 
01240   // Get finite element object
01241   FEGenericBase<OutputShape>* fe = NULL;
01242   this->get_element_fe<OutputShape>( var, fe );
01243 
01244   // Build a FE for calculating u(p)
01245   UniquePtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );
01246 
01247   // Get the values of the shape function derivatives
01248   const std::vector<std::vector<OutputShape> >&  phi = fe_new->get_phi();
01249 
01250   u = 0.;
01251 
01252   for (unsigned int l=0; l != n_dofs; l++)
01253     u += phi[l][0] * coef(l);
01254 
01255   return;
01256 }
01257 
01258 
01259 
01260 Gradient FEMContext::fixed_point_gradient(unsigned int var, const Point &p) const
01261 {
01262   Gradient grad_u;
01263 
01264   this->fixed_point_gradient( var, p, grad_u );
01265 
01266   return grad_u;
01267 }
01268 
01269 
01270 
01271 template<typename OutputType>
01272 void FEMContext::fixed_point_gradient(unsigned int var, const Point &p,
01273                                       OutputType& grad_u) const
01274 {
01275   typedef typename TensorTools::MakeReal
01276     <typename TensorTools::DecrementRank<OutputType>::type>::type
01277     OutputShape;
01278 
01279   // Get local-to-global dof index lookup
01280   libmesh_assert_greater (this->get_dof_indices().size(), var);
01281   const unsigned int n_dofs = cast_int<unsigned int>
01282     (this->get_dof_indices(var).size());
01283 
01284   // Get current local coefficients
01285   libmesh_assert_greater (_elem_fixed_subsolutions.size(), var);
01286   libmesh_assert(&(this->get_elem_fixed_solution(var)));
01287   const DenseSubVector<Number> &coef = this->get_elem_fixed_solution(var);
01288 
01289   // Get finite element object
01290   FEGenericBase<OutputShape>* fe = NULL;
01291   this->get_element_fe<OutputShape>( var, fe );
01292 
01293   // Build a FE for calculating u(p)
01294   UniquePtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );
01295 
01296   // Get the values of the shape function derivatives
01297   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputGradient> >&  dphi = fe_new->get_dphi();
01298 
01299   grad_u = 0.0;
01300 
01301   for (unsigned int l=0; l != n_dofs; l++)
01302     grad_u.add_scaled(dphi[l][0], coef(l));
01303 
01304   return;
01305 }
01306 
01307 
01308 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
01309 
01310 Tensor FEMContext::fixed_point_hessian(unsigned int var, const Point &p) const
01311 {
01312   Tensor hess_u;
01313 
01314   this->fixed_point_hessian( var, p, hess_u );
01315 
01316   return hess_u;
01317 }
01318 
01319 
01320 
01321 template<typename OutputType>
01322 void FEMContext::fixed_point_hessian(unsigned int var, const Point &p,
01323                                      OutputType& hess_u) const
01324 {
01325   typedef typename TensorTools::MakeReal<
01326     typename TensorTools::DecrementRank<
01327     typename TensorTools::DecrementRank<
01328     OutputType>::type>::type>::type
01329     OutputShape;
01330 
01331   // Get local-to-global dof index lookup
01332   libmesh_assert_greater (this->get_dof_indices().size(), var);
01333   const unsigned int n_dofs = cast_int<unsigned int>
01334     (this->get_dof_indices(var).size());
01335 
01336   // Get current local coefficients
01337   libmesh_assert_greater (_elem_fixed_subsolutions.size(), var);
01338   libmesh_assert(&(this->get_elem_fixed_solution(var)));
01339   const DenseSubVector<Number> &coef = this->get_elem_fixed_solution(var);
01340 
01341   // Get finite element object
01342   FEGenericBase<OutputShape>* fe = NULL;
01343   this->get_element_fe<OutputShape>( var, fe );
01344 
01345   // Build a FE for calculating u(p)
01346   UniquePtr<FEGenericBase<OutputShape> > fe_new = this->build_new_fe( fe, p );
01347 
01348   // Get the values of the shape function derivatives
01349   const std::vector<std::vector<typename FEGenericBase<OutputShape>::OutputTensor> >&  d2phi = fe_new->get_d2phi();
01350 
01351   hess_u = 0.0;
01352 
01353   for (unsigned int l=0; l != n_dofs; l++)
01354     hess_u.add_scaled(d2phi[l][0], coef(l));
01355 
01356   return;
01357 }
01358 
01359 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVES
01360 
01361 
01362 
01363 template<typename OutputType>
01364 void FEMContext::interior_rate(unsigned int var, unsigned int qp,
01365                                OutputType& u) const
01366 {
01367   this->some_interior_value <OutputType,&DiffContext::get_elem_solution_rate>
01368     (var, qp, u);
01369 }
01370 
01371 
01372 
01373 void FEMContext::elem_reinit(Real theta)
01374 {
01375   // Update the "time" variable of this context object
01376   this->_update_time_from_system(theta);
01377 
01378   // Handle a moving element if necessary.
01379   if (_mesh_sys)
01380     // We assume that the ``default'' state
01381     // of the mesh is its final, theta=1.0
01382     // position, so we don't bother with
01383     // mesh motion in that case.
01384     if (theta != 1.0)
01385       {
01386         // FIXME - ALE is not threadsafe yet!
01387         libmesh_assert_equal_to (libMesh::n_threads(), 1);
01388 
01389         elem_position_set(theta);
01390       }
01391   elem_fe_reinit();
01392 }
01393 
01394 
01395 void FEMContext::elem_side_reinit(Real theta)
01396 {
01397   // Update the "time" variable of this context object
01398   this->_update_time_from_system(theta);
01399 
01400   // Handle a moving element if necessary
01401   if (_mesh_sys)
01402     {
01403       // FIXME - not threadsafe yet!
01404       elem_position_set(theta);
01405       side_fe_reinit();
01406     }
01407 }
01408 
01409 
01410 void FEMContext::elem_edge_reinit(Real theta)
01411 {
01412   // Update the "time" variable of this context object
01413   this->_update_time_from_system(theta);
01414 
01415   // Handle a moving element if necessary
01416   if (_mesh_sys)
01417     {
01418       // FIXME - not threadsafe yet!
01419       elem_position_set(theta);
01420       edge_fe_reinit();
01421     }
01422 }
01423 
01424 
01425 void FEMContext::nonlocal_reinit(Real theta)
01426 {
01427   // Update the "time" variable of this context object
01428   this->_update_time_from_system(theta);
01429 
01430   // We can reuse the Elem FE safely here.
01431   elem_fe_reinit();
01432 }
01433 
01434 
01435 void FEMContext::elem_fe_reinit ()
01436 {
01437   // Initialize all the interior FE objects on elem.
01438   // Logging of FE::reinit is done in the FE functions
01439   // We only reinit the FE objects for the current element
01440   // dimension
01441   const unsigned char dim = this->get_elem_dim();
01442 
01443   libmesh_assert( !_element_fe[dim].empty() );
01444 
01445   std::map<FEType, FEAbstract *>::iterator local_fe_end = _element_fe[dim].end();
01446   for (std::map<FEType, FEAbstract *>::iterator i = _element_fe[dim].begin();
01447        i != local_fe_end; ++i)
01448     {
01449       if(this->has_elem())
01450         i->second->reinit(&(this->get_elem()));
01451       else
01452         // If !this->has_elem(), then we assume we are dealing with a SCALAR variable
01453         i->second->reinit(NULL);
01454     }
01455 }
01456 
01457 
01458 void FEMContext::side_fe_reinit ()
01459 {
01460   // Initialize all the side FE objects on elem/side.
01461   // Logging of FE::reinit is done in the FE functions
01462   // We only reinit the FE objects for the current element
01463   // dimension
01464   const unsigned char dim = this->get_elem_dim();
01465 
01466   libmesh_assert( !_side_fe[dim].empty() );
01467 
01468   std::map<FEType, FEAbstract *>::iterator local_fe_end = _side_fe[dim].end();
01469   for (std::map<FEType, FEAbstract *>::iterator i = _side_fe[dim].begin();
01470        i != local_fe_end; ++i)
01471     {
01472       i->second->reinit(&(this->get_elem()), this->get_side());
01473     }
01474 }
01475 
01476 
01477 
01478 void FEMContext::edge_fe_reinit ()
01479 {
01480   libmesh_assert_equal_to (this->get_elem_dim(), 3);
01481 
01482   // Initialize all the interior FE objects on elem/edge.
01483   // Logging of FE::reinit is done in the FE functions
01484   std::map<FEType, FEAbstract *>::iterator local_fe_end = _edge_fe.end();
01485   for (std::map<FEType, FEAbstract *>::iterator i = _edge_fe.begin();
01486        i != local_fe_end; ++i)
01487     {
01488       i->second->edge_reinit(&(this->get_elem()), this->get_edge());
01489     }
01490 }
01491 
01492 
01493 
01494 void FEMContext::elem_position_get()
01495 {
01496   // This is too expensive to call unless we've been asked to move the mesh
01497   libmesh_assert (_mesh_sys);
01498 
01499   // This will probably break with threading when two contexts are
01500   // operating on elements which share a node
01501   libmesh_assert_equal_to (libMesh::n_threads(), 1);
01502 
01503   // If the coordinate data is in our own system, it's already
01504   // been set up for us
01505   //  if (_mesh_sys == this->number())
01506   //    {
01507   unsigned int n_nodes = this->get_elem().n_nodes();
01508 
01509   const unsigned char dim = this->get_elem_dim();
01510 
01511   // For simplicity we demand that mesh coordinates be stored
01512   // in a format that allows a direct copy
01513   libmesh_assert(this->get_mesh_x_var() == libMesh::invalid_uint ||
01514                  (this->get_element_fe(this->get_mesh_x_var(), dim)->get_fe_type().family
01515                   == LAGRANGE &&
01516                   this->get_element_fe(this->get_mesh_x_var(), dim)->get_fe_type().order
01517                   == this->get_elem().default_order()));
01518   libmesh_assert(this->get_mesh_y_var() == libMesh::invalid_uint ||
01519                  (this->get_element_fe(this->get_mesh_y_var(), dim)->get_fe_type().family
01520                   == LAGRANGE &&
01521                   this->get_element_fe(this->get_mesh_y_var(), dim)->get_fe_type().order
01522                   == this->get_elem().default_order()));
01523   libmesh_assert(this->get_mesh_z_var() == libMesh::invalid_uint ||
01524                  (this->get_element_fe(this->get_mesh_z_var(), dim)->get_fe_type().family
01525                   == LAGRANGE &&
01526                   this->get_element_fe(this->get_mesh_z_var(), dim)->get_fe_type().order
01527                   == this->get_elem().default_order()));
01528 
01529   // Get degree of freedom coefficients from point coordinates
01530   if (this->get_mesh_x_var() != libMesh::invalid_uint)
01531     for (unsigned int i=0; i != n_nodes; ++i)
01532       (this->get_elem_solution(this->get_mesh_x_var()))(i) = this->get_elem().point(i)(0);
01533 
01534   if (this->get_mesh_y_var() != libMesh::invalid_uint)
01535     for (unsigned int i=0; i != n_nodes; ++i)
01536       (this->get_elem_solution(this->get_mesh_y_var()))(i) = this->get_elem().point(i)(1);
01537 
01538   if (this->get_mesh_z_var() != libMesh::invalid_uint)
01539     for (unsigned int i=0; i != n_nodes; ++i)
01540       (this->get_elem_solution(this->get_mesh_z_var()))(i) = this->get_elem().point(i)(2);
01541   //    }
01542   // FIXME - If the coordinate data is not in our own system, someone
01543   // had better get around to implementing that... - RHS
01544   //  else
01545   //    {
01546   //      libmesh_not_implemented();
01547   //    }
01548 }
01549 
01550 
01551 
01552 // We can ignore the theta argument in the current use of this
01553 // function, because elem_subsolutions will already have been set to
01554 // the theta value.
01555 //
01556 // To enable loose mesh movement coupling things will need to change.
01557 void FEMContext::_do_elem_position_set(Real)
01558 {
01559   // This is too expensive to call unless we've been asked to move the mesh
01560   libmesh_assert (_mesh_sys);
01561 
01562   // This will probably break with threading when two contexts are
01563   // operating on elements which share a node
01564   libmesh_assert_equal_to (libMesh::n_threads(), 1);
01565 
01566   const unsigned char dim = this->get_elem_dim();
01567 
01568   // If the coordinate data is in our own system, it's already
01569   // been set up for us, and we can ignore our input parameter theta
01570   //  if (_mesh_sys == this->number())
01571   //    {
01572   unsigned int n_nodes = this->get_elem().n_nodes();
01573   // For simplicity we demand that mesh coordinates be stored
01574   // in a format that allows a direct copy
01575   libmesh_assert(this->get_mesh_x_var() == libMesh::invalid_uint ||
01576                  (this->get_element_fe(this->get_mesh_x_var(), dim)->get_fe_type().family
01577                   == LAGRANGE &&
01578                   this->get_elem_solution(this->get_mesh_x_var()).size() == n_nodes));
01579   libmesh_assert(this->get_mesh_y_var() == libMesh::invalid_uint ||
01580                  (this->get_element_fe(this->get_mesh_y_var(), dim)->get_fe_type().family
01581                   == LAGRANGE &&
01582                   this->get_elem_solution(this->get_mesh_y_var()).size() == n_nodes));
01583   libmesh_assert(this->get_mesh_z_var() == libMesh::invalid_uint ||
01584                  (this->get_element_fe(this->get_mesh_z_var(), dim)->get_fe_type().family
01585                   == LAGRANGE &&
01586                   this->get_elem_solution(this->get_mesh_z_var()).size() == n_nodes));
01587 
01588   // Set the new point coordinates
01589   if (this->get_mesh_x_var() != libMesh::invalid_uint)
01590     for (unsigned int i=0; i != n_nodes; ++i)
01591       const_cast<Elem&>(this->get_elem()).point(i)(0) =
01592         libmesh_real(this->get_elem_solution(this->get_mesh_x_var())(i));
01593 
01594   if (this->get_mesh_y_var() != libMesh::invalid_uint)
01595     for (unsigned int i=0; i != n_nodes; ++i)
01596       const_cast<Elem&>(this->get_elem()).point(i)(1) =
01597         libmesh_real(this->get_elem_solution(this->get_mesh_y_var())(i));
01598 
01599   if (this->get_mesh_z_var() != libMesh::invalid_uint)
01600     for (unsigned int i=0; i != n_nodes; ++i)
01601       const_cast<Elem&>(this->get_elem()).point(i)(2) =
01602         libmesh_real(this->get_elem_solution(this->get_mesh_z_var())(i));
01603   //    }
01604   // FIXME - If the coordinate data is not in our own system, someone
01605   // had better get around to implementing that... - RHS
01606   //  else
01607   //    {
01608   //      libmesh_not_implemented();
01609   //    }
01610 }
01611 
01612 
01613 
01614 
01615 
01616 /*
01617   void FEMContext::reinit(const FEMSystem &sys, Elem *e)
01618   {
01619   // Initialize our elem pointer, algebraic objects
01620   this->pre_fe_reinit(e);
01621 
01622   // Moving the mesh may be necessary
01623   // Reinitializing the FE objects is definitely necessary
01624   this->elem_reinit(1.);
01625   }
01626 */
01627 
01628 
01629 
01630 void FEMContext::pre_fe_reinit(const System &sys, const Elem *e)
01631 {
01632   this->set_elem(e);
01633 
01634   // Initialize the per-element data for elem.
01635   if(this->has_elem())
01636     sys.get_dof_map().dof_indices (&(this->get_elem()), this->get_dof_indices());
01637   else
01638     // If !this->has_elem(), then we assume we are dealing with a SCALAR variable
01639     sys.get_dof_map().dof_indices (NULL, this->get_dof_indices());
01640 
01641   const unsigned int n_dofs = cast_int<unsigned int>
01642     (this->get_dof_indices().size());
01643   const std::size_t n_qoi = sys.qoi.size();
01644 
01645   // This also resizes elem_solution
01646   sys.current_local_solution->get(this->get_dof_indices(), this->get_elem_solution().get_values());
01647 
01648   if (sys.use_fixed_solution)
01649     this->get_elem_fixed_solution().resize(n_dofs);
01650   this->get_elem_solution_rate().resize(n_dofs);
01651 
01652   // These resize calls also zero out the residual and jacobian
01653   this->get_elem_residual().resize(n_dofs);
01654   this->get_elem_jacobian().resize(n_dofs, n_dofs);
01655 
01656   this->get_qoi_derivatives().resize(n_qoi);
01657   this->_elem_qoi_subderivatives.resize(n_qoi);
01658   for (std::size_t q=0; q != n_qoi; ++q)
01659     (this->get_qoi_derivatives())[q].resize(n_dofs);
01660 
01661   // Initialize the per-variable data for elem.
01662   {
01663     unsigned int sub_dofs = 0;
01664     for (unsigned int i=0; i != sys.n_vars(); ++i)
01665       {
01666         if(this->has_elem())
01667           sys.get_dof_map().dof_indices (&(this->get_elem()), this->get_dof_indices(i), i);
01668         else
01669           // If !this->has_elem(), then we assume we are dealing with a SCALAR variable
01670           sys.get_dof_map().dof_indices (NULL, this->get_dof_indices(i), i);
01671 
01672 
01673         const unsigned int n_dofs_var = cast_int<unsigned int>
01674           (this->get_dof_indices(i).size());
01675 
01676         this->get_elem_solution(i).reposition
01677           (sub_dofs, n_dofs_var);
01678 
01679         this->get_elem_solution_rate(i).reposition
01680           (sub_dofs, n_dofs_var);
01681 
01682         if (sys.use_fixed_solution)
01683           this->get_elem_fixed_solution(i).reposition
01684             (sub_dofs, n_dofs_var);
01685 
01686         this->get_elem_residual(i).reposition
01687           (sub_dofs, n_dofs_var);
01688 
01689         for (std::size_t q=0; q != n_qoi; ++q)
01690           this->get_qoi_derivatives(q,i).reposition
01691             (sub_dofs, n_dofs_var);
01692 
01693         for (unsigned int j=0; j != i; ++j)
01694           {
01695             const unsigned int n_dofs_var_j =
01696               cast_int<unsigned int>
01697               (this->get_dof_indices(j).size());
01698 
01699             this->get_elem_jacobian(i,j).reposition
01700               (sub_dofs, this->get_elem_residual(j).i_off(),
01701                n_dofs_var, n_dofs_var_j);
01702             this->get_elem_jacobian(j,i).reposition
01703               (this->get_elem_residual(j).i_off(), sub_dofs,
01704                n_dofs_var_j, n_dofs_var);
01705           }
01706         this->get_elem_jacobian(i,i).reposition
01707           (sub_dofs, sub_dofs,
01708            n_dofs_var,
01709            n_dofs_var);
01710         sub_dofs += n_dofs_var;
01711       }
01712     libmesh_assert_equal_to (sub_dofs, n_dofs);
01713   }
01714 
01715   // Now do the localization for the user requested vectors
01716   DiffContext::localized_vectors_iterator localized_vec_it = this->_localized_vectors.begin();
01717   const DiffContext::localized_vectors_iterator localized_vec_end = this->_localized_vectors.end();
01718 
01719   for(; localized_vec_it != localized_vec_end; ++localized_vec_it)
01720     {
01721       const NumericVector<Number>& current_localized_vector = *localized_vec_it->first;
01722       DenseVector<Number>& target_vector = localized_vec_it->second.first;
01723 
01724       current_localized_vector.get(this->get_dof_indices(), target_vector.get_values());
01725 
01726       // Initialize the per-variable data for elem.
01727       unsigned int sub_dofs = 0;
01728       for (unsigned int i=0; i != sys.n_vars(); ++i)
01729         {
01730           const unsigned int n_dofs_var = cast_int<unsigned int>
01731             (this->get_dof_indices(i).size());
01732           sys.get_dof_map().dof_indices (&(this->get_elem()), this->get_dof_indices(i), i);
01733 
01734           localized_vec_it->second.second[i]->reposition
01735             (sub_dofs, n_dofs_var);
01736 
01737           sub_dofs += n_dofs_var;
01738         }
01739       libmesh_assert_equal_to (sub_dofs, n_dofs);
01740     }
01741 }
01742 
01743 void FEMContext::set_elem( const Elem* e )
01744 {
01745   this->_elem = e;
01746 
01747   // If e is NULL, we assume it's SCALAR and set _elem_dim to 0.
01748   this->_elem_dim = this->_elem ? this->_elem->dim() : 0;
01749 }
01750 
01751 void FEMContext::_update_time_from_system(Real theta)
01752 {
01753   // Update the "time" variable based on the value of theta.  For this
01754   // to work, we need to know the value of deltat, a pointer to which is now
01755   // stored by our parent DiffContext class.  Note: get_deltat_value() will
01756   // assert in debug mode if the requested pointer is NULL.
01757   const Real deltat = this->get_deltat_value();
01758 
01759   this->set_time(theta*(this->get_system_time() + deltat) + (1.-theta)*this->get_system_time());
01760 }
01761 
01762 
01763 
01764 template<typename OutputShape>
01765 UniquePtr<FEGenericBase<OutputShape> > FEMContext::build_new_fe( const FEGenericBase<OutputShape>* fe,
01766                                                                  const Point &p ) const
01767 {
01768   FEType fe_type = fe->get_fe_type();
01769 
01770   // If we don't have an Elem to evaluate on, then the only functions
01771   // we can sensibly evaluate are the scalar dofs which are the same
01772   // everywhere.
01773   libmesh_assert(this->has_elem() || fe_type.family == SCALAR);
01774 
01775   unsigned int elem_dim = this->has_elem() ? this->get_elem().dim() : 0;
01776 
01777   FEGenericBase<OutputShape>* fe_new =
01778     FEGenericBase<OutputShape>::build(elem_dim, fe_type).release();
01779 
01780   // Map the physical co-ordinates to the master co-ordinates using the inverse_map from fe_interface.h
01781   // Build a vector of point co-ordinates to send to reinit
01782   Point master_point = this->has_elem() ?
01783     FEInterface::inverse_map(elem_dim, fe_type, &this->get_elem(), p) :
01784     Point(0);
01785 
01786   std::vector<Point> coor(1, master_point);
01787 
01788   // Reinitialize the element and compute the shape function values at coor
01789   if(this->has_elem())
01790     fe_new->reinit (&this->get_elem(), &coor);
01791   else
01792     // If !this->has_elem(), then we assume we are dealing with a SCALAR variable
01793     fe_new->reinit (NULL, &coor);
01794 
01795   return UniquePtr<FEGenericBase<OutputShape> >(fe_new);
01796 }
01797 
01798 
01799 
01800 
01801 
01802 // Instantiate member function templates
01803 template void FEMContext::interior_value<Number>(unsigned int, unsigned int, Number&) const;
01804 template void FEMContext::interior_values<Number>(unsigned int, const NumericVector<Number> &,
01805                                                   std::vector<Number>&) const;
01806 template void FEMContext::interior_value<Gradient>(unsigned int, unsigned int, Gradient&) const;
01807 template void FEMContext::interior_values<Gradient>(unsigned int, const NumericVector<Number> &,
01808                                                     std::vector<Gradient>&) const;
01809 
01810 template void FEMContext::interior_gradient<Gradient>(unsigned int, unsigned int, Gradient&) const;
01811 template void FEMContext::interior_gradients<Gradient>(unsigned int, const NumericVector<Number> &,
01812                                                        std::vector<Gradient>&) const;
01813 template void FEMContext::interior_gradient<Tensor>(unsigned int, unsigned int, Tensor&) const;
01814 template void FEMContext::interior_gradients<Tensor>(unsigned int, const NumericVector<Number> &,
01815                                                      std::vector<Tensor>&) const;
01816 
01817 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
01818 template void FEMContext::interior_hessian<Tensor>(unsigned int, unsigned int, Tensor&) const;
01819 template void FEMContext::interior_hessians<Tensor>(unsigned int, const NumericVector<Number> &,
01820                                                     std::vector<Tensor>&) const;
01821 //FIXME: Not everything is implemented yet for second derivatives of RealGradients
01822 //template void FEMContext::interior_hessian<??>(unsigned int, unsigned int, ??&) const;
01823 //template void FEMContext::interior_hessians<??>(unsigned int, const NumericVector<Number> &,
01824 //                                                std::vector<??>&) const;
01825 #endif
01826 
01827 template void FEMContext::interior_curl<Gradient>(unsigned int, unsigned int, Gradient&) const;
01828 
01829 template void FEMContext::interior_div<Number>(unsigned int, unsigned int, Number&) const;
01830 
01831 template void FEMContext::side_value<Number>(unsigned int, unsigned int, Number&) const;
01832 template void FEMContext::side_value<Gradient>(unsigned int, unsigned int, Gradient&) const;
01833 template void FEMContext::side_values<Number>(unsigned int, const NumericVector<Number> &,
01834                                               std::vector<Number>&) const;
01835 template void FEMContext::side_values<Gradient>(unsigned int, const NumericVector<Number> &,
01836                                                 std::vector<Gradient>&) const;
01837 
01838 template void FEMContext::side_gradient<Gradient>(unsigned int, unsigned int, Gradient&) const;
01839 template void FEMContext::side_gradients<Gradient>(unsigned int, const NumericVector<Number> &,
01840                                                    std::vector<Gradient>&) const;
01841 template void FEMContext::side_gradient<Tensor>(unsigned int, unsigned int, Tensor&) const;
01842 template void FEMContext::side_gradients<Tensor>(unsigned int, const NumericVector<Number> &,
01843                                                  std::vector<Tensor>&) const;
01844 
01845 
01846 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
01847 template void FEMContext::side_hessian<Tensor>(unsigned int, unsigned int, Tensor&) const;
01848 template void FEMContext::side_hessians<Tensor>(unsigned int, const NumericVector<Number> &,
01849                                                 std::vector<Tensor>&) const;
01850 //FIXME: Not everything is implemented yet for second derivatives of RealGradients
01851 //template void FEMContext::side_hessian<??>(unsigned int, unsigned int,
01852 //                                           ??&) const;
01853 //template void FEMContext::side_hessians<??>(unsigned int, const NumericVector<Number> &,
01854 //                                            std::vector<??>&) const;
01855 #endif
01856 
01857 template void FEMContext::point_value<Number>(unsigned int, const Point&, Number&) const;
01858 template void FEMContext::point_value<Gradient>(unsigned int, const Point&, Gradient&) const;
01859 
01860 template void FEMContext::point_gradient<Gradient>(unsigned int, const Point&, Gradient&) const;
01861 template void FEMContext::point_gradient<Tensor>(unsigned int, const Point&, Tensor&) const;
01862 
01863 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
01864 template void FEMContext::point_hessian<Tensor>(unsigned int, const Point&, Tensor&) const;
01865 //FIXME: Not everything is implemented yet for second derivatives of RealGradients
01866 //template void FEMContext::point_hessian<??>(unsigned int, const Point&, ??&) const;
01867 #endif
01868 
01869 template void FEMContext::point_curl<Gradient>(unsigned int, const Point&, Gradient&) const;
01870 
01871 template void FEMContext::fixed_interior_value<Number>(unsigned int, unsigned int, Number&) const;
01872 template void FEMContext::fixed_interior_value<Gradient>(unsigned int, unsigned int, Gradient&) const;
01873 
01874 template void FEMContext::fixed_interior_gradient<Gradient>(unsigned int, unsigned int, Gradient&) const;
01875 template void FEMContext::fixed_interior_gradient<Tensor>(unsigned int, unsigned int, Tensor&) const;
01876 
01877 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
01878 template void FEMContext::fixed_interior_hessian<Tensor>(unsigned int, unsigned int, Tensor&) const;
01879 //FIXME: Not everything is implemented yet for second derivatives of RealGradients
01880 //template void FEMContext::fixed_interior_hessian<??>(unsigned int, unsigned int, ??&) const;
01881 #endif
01882 
01883 template void FEMContext::fixed_side_value<Number>(unsigned int, unsigned int, Number&) const;
01884 template void FEMContext::fixed_side_value<Gradient>(unsigned int, unsigned int, Gradient&) const;
01885 
01886 template void FEMContext::fixed_side_gradient<Gradient>(unsigned int, unsigned int, Gradient&) const;
01887 template void FEMContext::fixed_side_gradient<Tensor>(unsigned int, unsigned int, Tensor&) const;
01888 
01889 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
01890 template void FEMContext::fixed_side_hessian<Tensor>(unsigned int, unsigned int, Tensor&) const;
01891 //FIXME: Not everything is implemented yet for second derivatives of RealGradients
01892 //template void FEMContext::fixed_side_hessian<??>(unsigned int, unsigned int, ??&) const;
01893 #endif
01894 
01895 template void FEMContext::fixed_point_value<Number>(unsigned int, const Point&, Number&) const;
01896 template void FEMContext::fixed_point_value<Gradient>(unsigned int, const Point&, Gradient&) const;
01897 
01898 template void FEMContext::fixed_point_gradient<Gradient>(unsigned int, const Point&, Gradient&) const;
01899 template void FEMContext::fixed_point_gradient<Tensor>(unsigned int, const Point&, Tensor&) const;
01900 
01901 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
01902 template void FEMContext::fixed_point_hessian<Tensor>(unsigned int, const Point&, Tensor&) const;
01903 //FIXME: Not everything is implemented yet for second derivatives of RealGradients
01904 //template void FEMContext::fixed_point_hessian<??>(unsigned int, const Point&, ??&) const;
01905 #endif
01906 
01907 template void FEMContext::interior_rate<Number>(unsigned int, unsigned int, Number&) const;
01908 template void FEMContext::interior_rate<Gradient>(unsigned int, unsigned int, Gradient&) const;
01909 
01910 template UniquePtr<FEGenericBase<Real> > FEMContext::build_new_fe( const FEGenericBase<Real>*, const Point & ) const;
01911 template UniquePtr<FEGenericBase<RealGradient> > FEMContext::build_new_fe( const FEGenericBase<RealGradient>*, const Point & ) const;
01912 
01913 
01914 } // namespace libMesh