$extrastylesheet
fem_system.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/dof_map.h"
00021 #include "libmesh/elem.h"
00022 #include "libmesh/equation_systems.h"
00023 #include "libmesh/fe_base.h"
00024 #include "libmesh/fem_context.h"
00025 #include "libmesh/fem_system.h"
00026 #include "libmesh/libmesh_logging.h"
00027 #include "libmesh/mesh_base.h"
00028 #include "libmesh/numeric_vector.h"
00029 #include "libmesh/parallel.h"
00030 #include "libmesh/parallel_algebra.h"
00031 #include "libmesh/parallel_ghost_sync.h"
00032 #include "libmesh/quadrature.h"
00033 #include "libmesh/sparse_matrix.h"
00034 #include "libmesh/time_solver.h"
00035 #include "libmesh/unsteady_solver.h" // For eulerian_residual
00036 #include "libmesh/fe_interface.h"
00037 
00038 namespace {
00039 using namespace libMesh;
00040 
00041 // give this guy some scope since there
00042 // is underlying vector allocation upon
00043 // creation/deletion
00044 ConstElemRange elem_range;
00045 
00046 typedef Threads::spin_mutex femsystem_mutex;
00047 femsystem_mutex assembly_mutex;
00048 
00049 void assemble_unconstrained_element_system
00050 (const FEMSystem& _sys,
00051  const bool _get_jacobian,
00052  const bool _constrain_heterogeneously,
00053  FEMContext &_femcontext)
00054 {
00055   if (_sys.print_element_solutions)
00056     {
00057       std::streamsize old_precision = libMesh::out.precision();
00058       libMesh::out.precision(16);
00059       if (_femcontext.has_elem())
00060         libMesh::out << "U_elem " << _femcontext.get_elem().id();
00061       else
00062         libMesh::out << "U_scalar ";
00063       libMesh::out << " = " << _femcontext.get_elem_solution() << std::endl;
00064 
00065       if (_sys.use_fixed_solution)
00066         {
00067           if (_femcontext.has_elem())
00068             libMesh::out << "Ufixed_elem " << _femcontext.get_elem().id();
00069           else
00070             libMesh::out << "Ufixed_scalar ";
00071           libMesh::out << " = " << _femcontext.get_elem_fixed_solution() << std::endl;
00072           libMesh::out.precision(old_precision);
00073         }
00074     }
00075 
00076   // We need jacobians to do heterogeneous residual constraints
00077   const bool need_jacobian =
00078     _get_jacobian || _constrain_heterogeneously;
00079 
00080   bool jacobian_computed =
00081     _sys.time_solver->element_residual(need_jacobian, _femcontext);
00082 
00083   // Compute a numeric jacobian if we have to
00084   if (need_jacobian && !jacobian_computed)
00085     {
00086       // Make sure we didn't compute a jacobian and lie about it
00087       libmesh_assert_equal_to (_femcontext.get_elem_jacobian().l1_norm(), 0.0);
00088       // Logging of numerical jacobians is done separately
00089       _sys.numerical_elem_jacobian(_femcontext);
00090     }
00091 
00092   // Compute a numeric jacobian if we're asked to verify the
00093   // analytic jacobian we got
00094   if (need_jacobian && jacobian_computed &&
00095       _sys.verify_analytic_jacobians != 0.0)
00096     {
00097       DenseMatrix<Number> analytic_jacobian(_femcontext.get_elem_jacobian());
00098 
00099       _femcontext.get_elem_jacobian().zero();
00100       // Logging of numerical jacobians is done separately
00101       _sys.numerical_elem_jacobian(_femcontext);
00102 
00103       Real analytic_norm = analytic_jacobian.l1_norm();
00104       Real numerical_norm = _femcontext.get_elem_jacobian().l1_norm();
00105 
00106       // If we can continue, we'll probably prefer the analytic jacobian
00107       analytic_jacobian.swap(_femcontext.get_elem_jacobian());
00108 
00109       // The matrix "analytic_jacobian" will now hold the error matrix
00110       analytic_jacobian.add(-1.0, _femcontext.get_elem_jacobian());
00111       Real error_norm = analytic_jacobian.l1_norm();
00112 
00113       Real relative_error = error_norm /
00114         std::max(analytic_norm, numerical_norm);
00115 
00116       if (relative_error > _sys.verify_analytic_jacobians)
00117         {
00118           libMesh::err << "Relative error " << relative_error
00119                        << " detected in analytic jacobian on element "
00120                        << _femcontext.get_elem().id() << '!' << std::endl;
00121 
00122           std::streamsize old_precision = libMesh::out.precision();
00123           libMesh::out.precision(16);
00124           libMesh::out << "J_analytic " << _femcontext.get_elem().id() << " = "
00125                        << _femcontext.get_elem_jacobian() << std::endl;
00126           analytic_jacobian.add(1.0, _femcontext.get_elem_jacobian());
00127           libMesh::out << "J_numeric " << _femcontext.get_elem().id() << " = "
00128                        << analytic_jacobian << std::endl;
00129 
00130           libMesh::out.precision(old_precision);
00131 
00132           libmesh_error_msg("Relative error too large, exiting!");
00133         }
00134     }
00135 
00136   for (_femcontext.side = 0;
00137        _femcontext.side != _femcontext.get_elem().n_sides();
00138        ++_femcontext.side)
00139     {
00140       // Don't compute on non-boundary sides unless requested
00141       if (!_sys.get_physics()->compute_internal_sides &&
00142           _femcontext.get_elem().neighbor(_femcontext.side) != NULL)
00143         continue;
00144 
00145       // Any mesh movement has already been done (and restored,
00146       // if the TimeSolver isn't broken), but
00147       // reinitializing the side FE objects is still necessary
00148       _femcontext.side_fe_reinit();
00149 
00150       DenseMatrix<Number> old_jacobian;
00151       // If we're in DEBUG mode, we should always verify that the
00152       // user's side_residual function doesn't alter our existing
00153       // jacobian and then lie about it
00154 #ifndef DEBUG
00155       // Even if we're not in DEBUG mode, when we're verifying
00156       // analytic jacobians we'll want to verify each side's
00157       // jacobian contribution separately.
00158       /* PB: We also need to account for the case when the user wants to
00159          use numerical Jacobians and not analytic Jacobians */
00160       if ( (_sys.verify_analytic_jacobians != 0.0 && need_jacobian) ||
00161            (!jacobian_computed && need_jacobian) )
00162 #endif // ifndef DEBUG
00163         {
00164           old_jacobian = _femcontext.get_elem_jacobian();
00165           _femcontext.get_elem_jacobian().zero();
00166         }
00167       jacobian_computed =
00168         _sys.time_solver->side_residual(need_jacobian, _femcontext);
00169 
00170       // Compute a numeric jacobian if we have to
00171       if (need_jacobian && !jacobian_computed)
00172         {
00173           // In DEBUG mode, we've already set elem_jacobian == 0,
00174           // so we can make sure side_residual didn't compute a
00175           // jacobian and lie about it
00176           libmesh_assert_equal_to (_femcontext.get_elem_jacobian().l1_norm(), 0.0);
00177 
00178           // Logging of numerical jacobians is done separately
00179           _sys.numerical_side_jacobian(_femcontext);
00180 
00181           // Add back in element interior numerical Jacobian
00182           _femcontext.get_elem_jacobian() += old_jacobian;
00183         }
00184 
00185       // Compute a numeric jacobian if we're asked to verify the
00186       // analytic jacobian we got
00187       if (need_jacobian && jacobian_computed &&
00188           _sys.verify_analytic_jacobians != 0.0)
00189         {
00190           DenseMatrix<Number> analytic_jacobian(_femcontext.get_elem_jacobian());
00191 
00192           _femcontext.get_elem_jacobian().zero();
00193           // Logging of numerical jacobians is done separately
00194           _sys.numerical_side_jacobian(_femcontext);
00195 
00196           Real analytic_norm = analytic_jacobian.l1_norm();
00197           Real numerical_norm = _femcontext.get_elem_jacobian().l1_norm();
00198 
00199           // If we can continue, we'll probably prefer the analytic jacobian
00200           analytic_jacobian.swap(_femcontext.get_elem_jacobian());
00201 
00202           // The matrix "analytic_jacobian" will now hold the error matrix
00203           analytic_jacobian.add(-1.0, _femcontext.get_elem_jacobian());
00204           Real error_norm = analytic_jacobian.l1_norm();
00205 
00206           Real relative_error = error_norm /
00207             std::max(analytic_norm, numerical_norm);
00208 
00209           if (relative_error > _sys.verify_analytic_jacobians)
00210             {
00211               libMesh::err << "Relative error " << relative_error
00212                            << " detected in analytic jacobian on element "
00213                            << _femcontext.get_elem().id()
00214                            << ", side "
00215                            << static_cast<unsigned int>(_femcontext.side) << '!' << std::endl;
00216 
00217               std::streamsize old_precision = libMesh::out.precision();
00218               libMesh::out.precision(16);
00219               libMesh::out << "J_analytic " << _femcontext.get_elem().id() << " = "
00220                            << _femcontext.get_elem_jacobian() << std::endl;
00221               analytic_jacobian.add(1.0, _femcontext.get_elem_jacobian());
00222               libMesh::out << "J_numeric " << _femcontext.get_elem().id() << " = "
00223                            << analytic_jacobian << std::endl;
00224               libMesh::out.precision(old_precision);
00225 
00226               libmesh_error_msg("Relative error too large, exiting!");
00227             }
00228           // Once we've verified a side, we'll want to add back the
00229           // rest of the accumulated jacobian
00230           _femcontext.get_elem_jacobian() += old_jacobian;
00231         }
00232       // In DEBUG mode, we've set elem_jacobian == 0, and we
00233       // may still need to add the old jacobian back
00234 #ifdef DEBUG
00235       if (need_jacobian && jacobian_computed &&
00236           _sys.verify_analytic_jacobians == 0.0)
00237         {
00238           _femcontext.get_elem_jacobian() += old_jacobian;
00239         }
00240 #endif // ifdef DEBUG
00241     }
00242 }
00243 
00244 void add_element_system
00245 (const FEMSystem& _sys,
00246  const bool _get_residual,
00247  const bool _get_jacobian,
00248  const bool _constrain_heterogeneously,
00249  FEMContext &_femcontext)
00250 {
00251 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00252   if (_get_residual && _sys.print_element_residuals)
00253     {
00254       std::streamsize old_precision = libMesh::out.precision();
00255       libMesh::out.precision(16);
00256       if (_femcontext.has_elem())
00257         libMesh::out << "Rraw_elem " << _femcontext.get_elem().id();
00258       else
00259         libMesh::out << "Rraw_scalar ";
00260       libMesh::out << " = " << _femcontext.get_elem_residual() << std::endl;
00261       libMesh::out.precision(old_precision);
00262     }
00263 
00264   // We turn off the asymmetric constraint application;
00265   // enforce_constraints_exactly() should be called in the solver
00266   if (_get_residual && _get_jacobian)
00267     {
00268       if (_constrain_heterogeneously)
00269         _sys.get_dof_map().heterogenously_constrain_element_matrix_and_vector
00270           (_femcontext.get_elem_jacobian(),
00271            _femcontext.get_elem_residual(),
00272            _femcontext.get_dof_indices(), false);
00273       else
00274         _sys.get_dof_map().constrain_element_matrix_and_vector
00275           (_femcontext.get_elem_jacobian(),
00276            _femcontext.get_elem_residual(),
00277            _femcontext.get_dof_indices(), false);
00278     }
00279   else if (_get_residual)
00280     {
00281       if (_constrain_heterogeneously)
00282         _sys.get_dof_map().heterogenously_constrain_element_vector
00283           (_femcontext.get_elem_jacobian(),
00284            _femcontext.get_elem_residual(),
00285            _femcontext.get_dof_indices(), false);
00286       else
00287         _sys.get_dof_map().constrain_element_vector
00288           (_femcontext.get_elem_residual(), _femcontext.get_dof_indices(), false);
00289     }
00290   else if (_get_jacobian)
00291     {
00292       // Heterogeneous and homogeneous constraints are the same on the
00293       // matrix
00294         _sys.get_dof_map().constrain_element_matrix
00295           (_femcontext.get_elem_jacobian(), _femcontext.get_dof_indices(), false);
00296     }
00297 #endif // #ifdef LIBMESH_ENABLE_CONSTRAINTS
00298 
00299   if (_get_residual && _sys.print_element_residuals)
00300     {
00301       std::streamsize old_precision = libMesh::out.precision();
00302       libMesh::out.precision(16);
00303       if (_femcontext.has_elem())
00304         libMesh::out << "R_elem " << _femcontext.get_elem().id();
00305       else
00306         libMesh::out << "R_scalar ";
00307       libMesh::out << " = " << _femcontext.get_elem_residual() << std::endl;
00308       libMesh::out.precision(old_precision);
00309     }
00310 
00311   if (_get_jacobian && _sys.print_element_jacobians)
00312     {
00313       std::streamsize old_precision = libMesh::out.precision();
00314       libMesh::out.precision(16);
00315       if (_femcontext.has_elem())
00316         libMesh::out << "J_elem " << _femcontext.get_elem().id();
00317       else
00318         libMesh::out << "J_scalar ";
00319       libMesh::out << " = " << _femcontext.get_elem_jacobian() << std::endl;
00320       libMesh::out.precision(old_precision);
00321     }
00322 
00323   { // A lock is necessary around access to the global system
00324     femsystem_mutex::scoped_lock lock(assembly_mutex);
00325 
00326     if (_get_jacobian)
00327       _sys.matrix->add_matrix (_femcontext.get_elem_jacobian(),
00328                                _femcontext.get_dof_indices());
00329     if (_get_residual)
00330       _sys.rhs->add_vector (_femcontext.get_elem_residual(),
00331                             _femcontext.get_dof_indices());
00332   } // Scope for assembly mutex
00333 }
00334 
00335 
00336 
00337 class AssemblyContributions
00338 {
00339 public:
00343   AssemblyContributions(FEMSystem &sys,
00344                         bool get_residual,
00345                         bool get_jacobian,
00346                         bool constrain_heterogeneously) :
00347     _sys(sys),
00348     _get_residual(get_residual),
00349     _get_jacobian(get_jacobian),
00350     _constrain_heterogeneously(constrain_heterogeneously) {}
00351 
00355   void operator()(const ConstElemRange &range) const
00356   {
00357     UniquePtr<DiffContext> con = _sys.build_context();
00358     FEMContext &_femcontext = cast_ref<FEMContext&>(*con);
00359     _sys.init_context(_femcontext);
00360 
00361     for (ConstElemRange::const_iterator elem_it = range.begin();
00362          elem_it != range.end(); ++elem_it)
00363       {
00364         Elem *el = const_cast<Elem *>(*elem_it);
00365 
00366         _femcontext.pre_fe_reinit(_sys, el);
00367         _femcontext.elem_fe_reinit();
00368 
00369         assemble_unconstrained_element_system
00370           (_sys, _get_jacobian, _constrain_heterogeneously,
00371            _femcontext);
00372 
00373         add_element_system
00374           (_sys, _get_residual, _get_jacobian,
00375            _constrain_heterogeneously, _femcontext);
00376       }
00377   }
00378 
00379 private:
00380 
00381   FEMSystem& _sys;
00382 
00383   const bool _get_residual, _get_jacobian, _constrain_heterogeneously;
00384 };
00385 
00386 class PostprocessContributions
00387 {
00388 public:
00392   explicit
00393   PostprocessContributions(FEMSystem &sys) : _sys(sys) {}
00394 
00398   void operator()(const ConstElemRange &range) const
00399   {
00400     UniquePtr<DiffContext> con = _sys.build_context();
00401     FEMContext &_femcontext = cast_ref<FEMContext&>(*con);
00402     _sys.init_context(_femcontext);
00403 
00404     for (ConstElemRange::const_iterator elem_it = range.begin();
00405          elem_it != range.end(); ++elem_it)
00406       {
00407         Elem *el = const_cast<Elem *>(*elem_it);
00408         _femcontext.pre_fe_reinit(_sys, el);
00409 
00410         // Optionally initialize all the interior FE objects on elem.
00411         if (_sys.fe_reinit_during_postprocess)
00412           _femcontext.elem_fe_reinit();
00413 
00414         _sys.element_postprocess(_femcontext);
00415 
00416         for (_femcontext.side = 0;
00417              _femcontext.side != _femcontext.get_elem().n_sides();
00418              ++_femcontext.side)
00419           {
00420             // Don't compute on non-boundary sides unless requested
00421             if (!_sys.postprocess_sides ||
00422                 (!_sys.get_physics()->compute_internal_sides &&
00423                  _femcontext.get_elem().neighbor(_femcontext.side) != NULL))
00424               continue;
00425 
00426             // Optionally initialize all the FE objects on this side.
00427             if (_sys.fe_reinit_during_postprocess)
00428               _femcontext.side_fe_reinit();
00429 
00430             _sys.side_postprocess(_femcontext);
00431           }
00432       }
00433   }
00434 
00435 private:
00436 
00437   FEMSystem& _sys;
00438 };
00439 
00440 class QoIContributions
00441 {
00442 public:
00446   explicit
00447   QoIContributions(FEMSystem &sys, DifferentiableQoI &diff_qoi,
00448                    const QoISet &qoi_indices) :
00449     qoi(sys.qoi.size(), 0.), _sys(sys), _diff_qoi(diff_qoi),_qoi_indices(qoi_indices) {}
00450 
00454   QoIContributions(const QoIContributions &other,
00455                    Threads::split) :
00456     qoi(other._sys.qoi.size(), 0.), _sys(other._sys), _diff_qoi(other._diff_qoi) {}
00457 
00461   void operator()(const ConstElemRange &range)
00462   {
00463     UniquePtr<DiffContext> con = _sys.build_context();
00464     FEMContext &_femcontext = cast_ref<FEMContext&>(*con);
00465     _diff_qoi.init_context(_femcontext);
00466 
00467     std::vector<bool> have_heterogenous_qoi_bc(_sys.qoi.size(), false);
00468     bool have_some_heterogenous_qoi_bc = false;
00469 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00470     for (unsigned int q=0; q != _sys.qoi.size(); ++q)
00471       if (_qoi_indices.has_index(q) &&
00472           _sys.get_dof_map().has_heterogenous_adjoint_constraints(q))
00473         {
00474           have_heterogenous_qoi_bc[q] = true;
00475           have_some_heterogenous_qoi_bc = true;
00476         }
00477 #endif
00478 
00479     for (ConstElemRange::const_iterator elem_it = range.begin();
00480          elem_it != range.end(); ++elem_it)
00481       {
00482         Elem *el = const_cast<Elem *>(*elem_it);
00483 
00484         _femcontext.pre_fe_reinit(_sys, el);
00485 
00486         // We might have some heterogenous dofs here; let's see for
00487         // certain
00488 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00489         if (have_some_heterogenous_qoi_bc)
00490           {
00491             have_some_heterogenous_qoi_bc = false;
00492             for (unsigned int q=0; q != _sys.qoi.size(); ++q)
00493               {
00494                 if (have_heterogenous_qoi_bc[q])
00495                   {
00496                     have_heterogenous_qoi_bc[q] = false;
00497                     for (unsigned int d=0;
00498                          d != _femcontext.get_dof_indices().size(); ++d)
00499                       if (_sys.get_dof_map().has_heterogenous_adjoint_constraint
00500                           (q, _femcontext.get_dof_indices()[d]) != Number(0))
00501                         {
00502                           have_some_heterogenous_qoi_bc = true;
00503                           have_heterogenous_qoi_bc[q] = true;
00504                           break;
00505                         }
00506                   }
00507               }
00508           }
00509 #endif
00510 
00511         if (_diff_qoi.assemble_qoi_elements ||
00512             have_some_heterogenous_qoi_bc)
00513           _femcontext.elem_fe_reinit();
00514 
00515         if (_diff_qoi.assemble_qoi_elements)
00516           _diff_qoi.element_qoi(_femcontext, _qoi_indices);
00517 
00518         // If we have some heterogenous dofs here, those are
00519         // themselves part of a regularized flux QoI which the library
00520         // handles integrating
00521 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00522         if (have_some_heterogenous_qoi_bc)
00523           {
00524             _sys.time_solver->element_residual(false, _femcontext);
00525 
00526             for (unsigned int q=0; q != _sys.qoi.size(); ++q)
00527               {
00528                 if (have_heterogenous_qoi_bc[q])
00529                   {
00530                     for (unsigned int d=0;
00531                          d != _femcontext.get_dof_indices().size(); ++d)
00532                       this->qoi[q] -= _femcontext.get_elem_residual()(d) *
00533                         _sys.get_dof_map().has_heterogenous_adjoint_constraint(q, _femcontext.get_dof_indices()[d]);
00534                   }
00535               }
00536 
00537           }
00538 #endif
00539 
00540         for (_femcontext.side = 0;
00541              _femcontext.side != _femcontext.get_elem().n_sides();
00542              ++_femcontext.side)
00543           {
00544             // Don't compute on non-boundary sides unless requested
00545             if (!_diff_qoi.assemble_qoi_sides ||
00546                 (!_diff_qoi.assemble_qoi_internal_sides &&
00547                  _femcontext.get_elem().neighbor(_femcontext.side) != NULL))
00548               continue;
00549 
00550             _femcontext.side_fe_reinit();
00551 
00552             _diff_qoi.side_qoi(_femcontext, _qoi_indices);
00553           }
00554       }
00555 
00556     this->_diff_qoi.thread_join( this->qoi, _femcontext.get_qois(), _qoi_indices );
00557   }
00558 
00559   void join (const QoIContributions& other)
00560   {
00561     libmesh_assert_equal_to (this->qoi.size(), other.qoi.size());
00562     this->_diff_qoi.thread_join( this->qoi, other.qoi, _qoi_indices );
00563   }
00564 
00565   std::vector<Number> qoi;
00566 
00567 private:
00568 
00569   FEMSystem& _sys;
00570   DifferentiableQoI& _diff_qoi;
00571 
00572   const QoISet _qoi_indices;
00573 };
00574 
00575 class QoIDerivativeContributions
00576 {
00577 public:
00581   QoIDerivativeContributions(FEMSystem &sys, const QoISet& qoi_indices,
00582                              DifferentiableQoI &qoi,
00583                              bool include_liftfunc,
00584                              bool apply_constraints ) :
00585     _sys(sys),
00586     _qoi_indices(qoi_indices),
00587     _qoi(qoi),
00588     _include_liftfunc(include_liftfunc),
00589     _apply_constraints(apply_constraints) {}
00590 
00594   void operator()(const ConstElemRange &range) const
00595   {
00596     UniquePtr<DiffContext> con = _sys.build_context();
00597     FEMContext &_femcontext = cast_ref<FEMContext&>(*con);
00598     _qoi.init_context(_femcontext);
00599 
00600     bool have_some_heterogenous_qoi_bc = false;
00601 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00602     std::vector<bool> have_heterogenous_qoi_bc(_sys.qoi.size(), false);
00603     if (_include_liftfunc || _apply_constraints)
00604       for (unsigned int q=0; q != _sys.qoi.size(); ++q)
00605         if (_qoi_indices.has_index(q) &&
00606             _sys.get_dof_map().has_heterogenous_adjoint_constraints(q))
00607           {
00608             have_heterogenous_qoi_bc[q] = true;
00609             have_some_heterogenous_qoi_bc = true;
00610           }
00611 #endif
00612 
00613     for (ConstElemRange::const_iterator elem_it = range.begin();
00614          elem_it != range.end(); ++elem_it)
00615       {
00616         Elem *el = const_cast<Elem *>(*elem_it);
00617 
00618         _femcontext.pre_fe_reinit(_sys, el);
00619 
00620         // We might have some heterogenous dofs here; let's see for
00621         // certain
00622 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00623         bool elem_has_some_heterogenous_qoi_bc = false;
00624         std::vector<bool> elem_has_heterogenous_qoi_bc(_sys.qoi.size(), false);
00625         if (have_some_heterogenous_qoi_bc)
00626           {
00627             for (unsigned int q=0; q != _sys.qoi.size(); ++q)
00628               {
00629                 if (have_heterogenous_qoi_bc[q])
00630                   {
00631                     for (unsigned int d=0;
00632                          d != _femcontext.get_dof_indices().size(); ++d)
00633                       if (_sys.get_dof_map().has_heterogenous_adjoint_constraint
00634                           (q, _femcontext.get_dof_indices()[d]) != Number(0))
00635                         {
00636                           elem_has_some_heterogenous_qoi_bc = true;
00637                           elem_has_heterogenous_qoi_bc[q] = true;
00638                           break;
00639                         }
00640                   }
00641               }
00642           }
00643 #endif
00644 
00645         // If we're going to call a user integral, then we need FE
00646         // information to call element_qoi.
00647         // If we're going to evaluate lift-function-based components
00648         // of a QoI, then we need FE information to assemble the
00649         // element residual.
00650         if (_qoi.assemble_qoi_elements ||
00651             (_include_liftfunc &&
00652              elem_has_some_heterogenous_qoi_bc))
00653           _femcontext.elem_fe_reinit();
00654 
00655         if (_qoi.assemble_qoi_elements)
00656           _qoi.element_qoi_derivative(_femcontext, _qoi_indices);
00657 
00658 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00659         // If we need to use heterogenous dofs here, we need the
00660         // Jacobian either for the regularized flux QoI integration
00661         // and/or for constraint application.
00662         if ((_include_liftfunc || _apply_constraints) &&
00663             elem_has_some_heterogenous_qoi_bc)
00664           _sys.time_solver->element_residual(true, _femcontext);
00665 
00666         // If we have some heterogenous dofs here, those are
00667         // themselves part of a regularized flux QoI which the library
00668         // may handle integrating
00669         if (_include_liftfunc && elem_has_some_heterogenous_qoi_bc)
00670           {
00671             for (unsigned int q=0; q != _sys.qoi.size(); ++q)
00672               {
00673                 if (elem_has_heterogenous_qoi_bc[q])
00674                   {
00675                     for (unsigned int i=0;
00676                          i != _femcontext.get_dof_indices().size(); ++i)
00677                       {
00678                         Number liftfunc_val =
00679                           _sys.get_dof_map().has_heterogenous_adjoint_constraint(q, _femcontext.get_dof_indices()[i]);
00680 
00681                         if (liftfunc_val != Number(0))
00682                           {
00683                             for (unsigned int j=0;
00684                                  j != _femcontext.get_dof_indices().size(); ++j)
00685                               _femcontext.get_qoi_derivatives()[q](j) -=
00686                                 _femcontext.get_elem_jacobian()(i,j) *
00687                                 liftfunc_val;
00688                           }
00689                       }
00690                   }
00691               }
00692           }
00693 #endif
00694 
00695 
00696         for (_femcontext.side = 0;
00697              _femcontext.side != _femcontext.get_elem().n_sides();
00698              ++_femcontext.side)
00699           {
00700             // Don't compute on non-boundary sides unless requested
00701             if (!_qoi.assemble_qoi_sides ||
00702                 (!_qoi.assemble_qoi_internal_sides &&
00703                  _femcontext.get_elem().neighbor(_femcontext.side) != NULL))
00704               continue;
00705 
00706             _femcontext.side_fe_reinit();
00707 
00708             _qoi.side_qoi_derivative(_femcontext, _qoi_indices);
00709           }
00710 
00711         // We need some unmodified indices to use for constraining
00712         // multiple vector
00713         // FIXME - there should be a DofMap::constrain_element_vectors
00714         // to do this more efficiently
00715 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00716         std::vector<dof_id_type> original_dofs = _femcontext.get_dof_indices();
00717 #endif
00718 
00719         { // A lock is necessary around access to the global system
00720           femsystem_mutex::scoped_lock lock(assembly_mutex);
00721 
00722 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00723           // We'll need to see if any heterogenous constraints apply
00724           // to the QoI dofs on this element *or* to any of the dofs
00725           // they depend on, so let's get those dependencies
00726           if (_apply_constraints)
00727             _sys.get_dof_map().constrain_nothing(_femcontext.get_dof_indices());
00728 #endif
00729 
00730           for (unsigned int i=0; i != _sys.qoi.size(); ++i)
00731             if (_qoi_indices.has_index(i))
00732               {
00733 #ifdef LIBMESH_ENABLE_CONSTRAINTS
00734                 if (_apply_constraints)
00735                   {
00736 #ifndef NDEBUG
00737                     bool has_heterogenous_constraint = false;
00738                     for (unsigned int d=0;
00739                          d != _femcontext.get_dof_indices().size(); ++d)
00740                       if (_sys.get_dof_map().has_heterogenous_adjoint_constraint
00741                           (i, _femcontext.get_dof_indices()[d]) != Number(0))
00742                         {
00743                           has_heterogenous_constraint = true;
00744                           libmesh_assert(elem_has_heterogenous_qoi_bc[i]);
00745                           libmesh_assert(elem_has_some_heterogenous_qoi_bc);
00746                           break;
00747                         }
00748 #else
00749                     bool has_heterogenous_constraint =
00750                       elem_has_heterogenous_qoi_bc[i];
00751 #endif
00752 
00753                     _femcontext.get_dof_indices() = original_dofs;
00754 
00755                     // If we're going to need K to impose a heterogenous
00756                     // constraint, we may not have already computed it above.
00757                     if (has_heterogenous_constraint)
00758                       {
00759                         _sys.get_dof_map().heterogenously_constrain_element_vector
00760                           (_femcontext.get_elem_jacobian(),
00761                            _femcontext.get_qoi_derivatives()[i],
00762                            _femcontext.get_dof_indices(), false, i);
00763                       }
00764                     else
00765                       {
00766                         _sys.get_dof_map().constrain_element_vector
00767                           (_femcontext.get_qoi_derivatives()[i],
00768                            _femcontext.get_dof_indices(), false);
00769                       }
00770                   }
00771 #endif
00772 
00773                 _sys.get_adjoint_rhs(i).add_vector
00774                   (_femcontext.get_qoi_derivatives()[i], _femcontext.get_dof_indices());
00775               }
00776         }
00777       }
00778   }
00779 
00780 private:
00781 
00782   FEMSystem& _sys;
00783   const QoISet& _qoi_indices;
00784   DifferentiableQoI& _qoi;
00785   bool _include_liftfunc, _apply_constraints;
00786 };
00787 
00788 
00789 }
00790 
00791 
00792 namespace libMesh
00793 {
00794 
00795 
00796 
00797 
00798 
00799 FEMSystem::FEMSystem (EquationSystems& es,
00800                       const std::string& name_in,
00801                       const unsigned int number_in)
00802   : Parent(es, name_in, number_in),
00803     fe_reinit_during_postprocess(true),
00804     numerical_jacobian_h(TOLERANCE),
00805     verify_analytic_jacobians(0.0)
00806 {
00807 }
00808 
00809 
00810 FEMSystem::~FEMSystem ()
00811 {
00812   this->clear();
00813 }
00814 
00815 
00816 
00817 void FEMSystem::clear()
00818 {
00819   Parent::clear();
00820 }
00821 
00822 
00823 
00824 void FEMSystem::init_data ()
00825 {
00826   // First initialize LinearImplicitSystem data
00827   Parent::init_data();
00828 }
00829 
00830 
00831 void FEMSystem::assembly (bool get_residual, bool get_jacobian,
00832                           bool apply_heterogeneous_constraints)
00833 {
00834   libmesh_assert(get_residual || get_jacobian);
00835   std::string log_name;
00836   if (get_residual && get_jacobian)
00837     log_name = "assembly()";
00838   else if (get_residual)
00839     log_name = "assembly(get_residual)";
00840   else
00841     log_name = "assembly(get_jacobian)";
00842 
00843   START_LOG(log_name, "FEMSystem");
00844 
00845   const MeshBase& mesh = this->get_mesh();
00846 
00847   //  this->get_vector("_nonlinear_solution").localize
00848   //    (*current_local_nonlinear_solution,
00849   //     dof_map.get_send_list());
00850   this->update();
00851 
00852   if (print_solution_norms)
00853     {
00854       //      this->get_vector("_nonlinear_solution").close();
00855       this->solution->close();
00856 
00857       std::streamsize old_precision = libMesh::out.precision();
00858       libMesh::out.precision(16);
00859       libMesh::out << "|U| = "
00860         //                    << this->get_vector("_nonlinear_solution").l1_norm()
00861                    << this->solution->l1_norm()
00862                    << std::endl;
00863       libMesh::out.precision(old_precision);
00864     }
00865   if (print_solutions)
00866     {
00867       std::streamsize old_precision = libMesh::out.precision();
00868       libMesh::out.precision(16);
00869       //      libMesh::out << "U = [" << this->get_vector("_nonlinear_solution")
00870       libMesh::out << "U = [" << *(this->solution)
00871                    << "];" << std::endl;
00872       libMesh::out.precision(old_precision);
00873     }
00874 
00875   // Is this definitely necessary? [RHS]
00876   // Yes. [RHS 2012]
00877   if (get_jacobian)
00878     matrix->zero();
00879   if (get_residual)
00880     rhs->zero();
00881 
00882   // Stupid C++ lets you set *Real* verify_analytic_jacobians = true!
00883   if (verify_analytic_jacobians > 0.5)
00884     {
00885       libMesh::err << "WARNING!  verify_analytic_jacobians was set "
00886                    << "to absurdly large value of "
00887                    << verify_analytic_jacobians << std::endl;
00888       libMesh::err << "Resetting to 1e-6!" << std::endl;
00889       verify_analytic_jacobians = 1e-6;
00890     }
00891 
00892   // In time-dependent problems, the nonlinear function we're trying
00893   // to solve at each timestep may depend on the particular solver
00894   // we're using
00895   libmesh_assert(time_solver.get());
00896 
00897   // Build the residual and jacobian contributions on every active
00898   // mesh element on this processor
00899   Threads::parallel_for
00900     (elem_range.reset(mesh.active_local_elements_begin(),
00901                       mesh.active_local_elements_end()),
00902      AssemblyContributions(*this, get_residual, get_jacobian,
00903                            apply_heterogeneous_constraints));
00904 
00905   // Check and see if we have SCALAR variables
00906   bool have_scalar = false;
00907   for(unsigned int i=0; i != this->n_variable_groups(); ++i)
00908     {
00909       if( this->variable_group(i).type().family == SCALAR )
00910         {
00911           have_scalar = true;
00912           break;
00913         }
00914     }
00915 
00916   // SCALAR dofs are stored on the last processor, so we'll evaluate
00917   // their equation terms there and only if we have a SCALAR variable
00918   if ( this->processor_id() == (this->n_processors()-1) && have_scalar )
00919     {
00920       UniquePtr<DiffContext> con = this->build_context();
00921       FEMContext &_femcontext = cast_ref<FEMContext&>(*con);
00922       this->init_context(_femcontext);
00923       _femcontext.pre_fe_reinit(*this, NULL);
00924 
00925       bool jacobian_computed =
00926         this->time_solver->nonlocal_residual(get_jacobian, _femcontext);
00927 
00928       // Nonlocal residuals are likely to be length 0, in which case we
00929       // don't need to do any more.  And we shouldn't try to do any
00930       // more; lots of DenseVector/DenseMatrix code assumes rank>0.
00931       if (_femcontext.get_elem_residual().size())
00932         {
00933           // Compute a numeric jacobian if we have to
00934           if (get_jacobian && !jacobian_computed)
00935             {
00936               // Make sure we didn't compute a jacobian and lie about it
00937               libmesh_assert_equal_to (_femcontext.get_elem_jacobian().l1_norm(), 0.0);
00938               // Logging of numerical jacobians is done separately
00939               this->numerical_nonlocal_jacobian(_femcontext);
00940             }
00941 
00942           // Compute a numeric jacobian if we're asked to verify the
00943           // analytic jacobian we got
00944           if (get_jacobian && jacobian_computed &&
00945               this->verify_analytic_jacobians != 0.0)
00946             {
00947               DenseMatrix<Number> analytic_jacobian(_femcontext.get_elem_jacobian());
00948 
00949               _femcontext.get_elem_jacobian().zero();
00950               // Logging of numerical jacobians is done separately
00951               this->numerical_nonlocal_jacobian(_femcontext);
00952 
00953               Real analytic_norm = analytic_jacobian.l1_norm();
00954               Real numerical_norm = _femcontext.get_elem_jacobian().l1_norm();
00955 
00956               // If we can continue, we'll probably prefer the analytic jacobian
00957               analytic_jacobian.swap(_femcontext.get_elem_jacobian());
00958 
00959               // The matrix "analytic_jacobian" will now hold the error matrix
00960               analytic_jacobian.add(-1.0, _femcontext.get_elem_jacobian());
00961               Real error_norm = analytic_jacobian.l1_norm();
00962 
00963               Real relative_error = error_norm /
00964                 std::max(analytic_norm, numerical_norm);
00965 
00966               if (relative_error > this->verify_analytic_jacobians)
00967                 {
00968                   libMesh::err << "Relative error " << relative_error
00969                                << " detected in analytic jacobian on nonlocal dofs!"
00970                                << std::endl;
00971 
00972                   std::streamsize old_precision = libMesh::out.precision();
00973                   libMesh::out.precision(16);
00974                   libMesh::out << "J_analytic nonlocal = "
00975                                << _femcontext.get_elem_jacobian() << std::endl;
00976                   analytic_jacobian.add(1.0, _femcontext.get_elem_jacobian());
00977                   libMesh::out << "J_numeric nonlocal = "
00978                                << analytic_jacobian << std::endl;
00979 
00980                   libMesh::out.precision(old_precision);
00981 
00982                   libmesh_error_msg("Relative error too large, exiting!");
00983                 }
00984             }
00985 
00986           add_element_system
00987             (*this, get_residual, get_jacobian,
00988              apply_heterogeneous_constraints, _femcontext);
00989         }
00990     }
00991 
00992   if (get_residual && (print_residual_norms || print_residuals))
00993     this->rhs->close();
00994   if (get_residual && print_residual_norms)
00995     {
00996       std::streamsize old_precision = libMesh::out.precision();
00997       libMesh::out.precision(16);
00998       libMesh::out << "|F| = " << this->rhs->l1_norm() << std::endl;
00999       libMesh::out.precision(old_precision);
01000     }
01001   if (get_residual && print_residuals)
01002     {
01003       std::streamsize old_precision = libMesh::out.precision();
01004       libMesh::out.precision(16);
01005       libMesh::out << "F = [" << *(this->rhs) << "];" << std::endl;
01006       libMesh::out.precision(old_precision);
01007     }
01008 
01009   if (get_jacobian && (print_jacobian_norms || print_jacobians))
01010     this->matrix->close();
01011   if (get_jacobian && print_jacobian_norms)
01012     {
01013       std::streamsize old_precision = libMesh::out.precision();
01014       libMesh::out.precision(16);
01015       libMesh::out << "|J| = " << this->matrix->l1_norm() << std::endl;
01016       libMesh::out.precision(old_precision);
01017     }
01018   if (get_jacobian && print_jacobians)
01019     {
01020       std::streamsize old_precision = libMesh::out.precision();
01021       libMesh::out.precision(16);
01022       libMesh::out << "J = [" << *(this->matrix) << "];" << std::endl;
01023       libMesh::out.precision(old_precision);
01024     }
01025   STOP_LOG(log_name, "FEMSystem");
01026 }
01027 
01028 
01029 
01030 void FEMSystem::solve()
01031 {
01032   // We are solving the primal problem
01033   Parent::solve();
01034 
01035   // On a moving mesh we want the mesh to reflect the new solution
01036   this->mesh_position_set();
01037 }
01038 
01039 
01040 
01041 void FEMSystem::mesh_position_set()
01042 {
01043   // If we don't need to move the mesh, we're done
01044   if (_mesh_sys != this)
01045     return;
01046 
01047   MeshBase& mesh = this->get_mesh();
01048 
01049   UniquePtr<DiffContext> con = this->build_context();
01050   FEMContext &_femcontext = cast_ref<FEMContext&>(*con);
01051   this->init_context(_femcontext);
01052 
01053   // Move every mesh element we can
01054   MeshBase::const_element_iterator el =
01055     mesh.active_local_elements_begin();
01056   const MeshBase::const_element_iterator end_el =
01057     mesh.active_local_elements_end();
01058 
01059   for ( ; el != end_el; ++el)
01060     {
01061       // We need the algebraic data
01062       _femcontext.pre_fe_reinit(*this, *el);
01063       // And when asserts are on, we also need the FE so
01064       // we can assert that the mesh data is of the right type.
01065 #ifndef NDEBUG
01066       _femcontext.elem_fe_reinit();
01067 #endif
01068 
01069       // This code won't handle moving subactive elements
01070       libmesh_assert(!_femcontext.get_elem().has_children());
01071 
01072       _femcontext.elem_position_set(1.);
01073     }
01074 
01075   // We've now got positions set on all local nodes (and some
01076   // semilocal nodes); let's request positions for non-local nodes
01077   // from their processors.
01078 
01079   SyncNodalPositions sync_object(mesh);
01080   Parallel::sync_dofobject_data_by_id
01081     (this->comm(), mesh.nodes_begin(), mesh.nodes_end(), sync_object);
01082 }
01083 
01084 
01085 
01086 void FEMSystem::postprocess ()
01087 {
01088   START_LOG("postprocess()", "FEMSystem");
01089 
01090   const MeshBase& mesh = this->get_mesh();
01091 
01092   this->update();
01093 
01094   // Get the time solver object associated with the system, and tell it that
01095   // we are not solving the adjoint problem
01096   this->get_time_solver().set_is_adjoint(false);
01097 
01098   // Loop over every active mesh element on this processor
01099   Threads::parallel_for(elem_range.reset(mesh.active_local_elements_begin(),
01100                                          mesh.active_local_elements_end()),
01101                         PostprocessContributions(*this));
01102 
01103   STOP_LOG("postprocess()", "FEMSystem");
01104 }
01105 
01106 
01107 
01108 void FEMSystem::assemble_qoi (const QoISet &qoi_indices)
01109 {
01110   START_LOG("assemble_qoi()", "FEMSystem");
01111 
01112   const MeshBase& mesh = this->get_mesh();
01113 
01114   this->update();
01115 
01116   const unsigned int Nq = cast_int<unsigned int>(qoi.size());
01117 
01118   // the quantity of interest is assumed to be a sum of element and
01119   // side terms
01120   for (unsigned int i=0; i != Nq; ++i)
01121     if (qoi_indices.has_index(i))
01122       qoi[i] = 0;
01123 
01124   // Create a non-temporary qoi_contributions object, so we can query
01125   // its results after the reduction
01126   QoIContributions qoi_contributions(*this, *(this->diff_qoi), qoi_indices);
01127 
01128   // Loop over every active mesh element on this processor
01129   Threads::parallel_reduce(elem_range.reset(mesh.active_local_elements_begin(),
01130                                             mesh.active_local_elements_end()),
01131                            qoi_contributions);
01132 
01133   this->diff_qoi->parallel_op( this->comm(), this->qoi, qoi_contributions.qoi, qoi_indices );
01134 
01135   STOP_LOG("assemble_qoi()", "FEMSystem");
01136 }
01137 
01138 
01139 
01140 void FEMSystem::assemble_qoi_derivative (const QoISet& qoi_indices,
01141                                          bool include_liftfunc,
01142                                          bool apply_constraints)
01143 {
01144   START_LOG("assemble_qoi_derivative()", "FEMSystem");
01145 
01146   const MeshBase& mesh = this->get_mesh();
01147 
01148   this->update();
01149 
01150   // The quantity of interest derivative assembly accumulates on
01151   // initially zero vectors
01152   for (unsigned int i=0; i != qoi.size(); ++i)
01153     if (qoi_indices.has_index(i))
01154       this->add_adjoint_rhs(i).zero();
01155 
01156   // Loop over every active mesh element on this processor
01157   Threads::parallel_for(elem_range.reset(mesh.active_local_elements_begin(),
01158                                          mesh.active_local_elements_end()),
01159                         QoIDerivativeContributions(*this, qoi_indices,
01160                                                    *(this->diff_qoi),
01161                                                    include_liftfunc,
01162                                                    apply_constraints));
01163 
01164   STOP_LOG("assemble_qoi_derivative()", "FEMSystem");
01165 }
01166 
01167 
01168 
01169 void FEMSystem::numerical_jacobian (TimeSolverResPtr res,
01170                                     FEMContext &context) const
01171 {
01172   // Logging is done by numerical_elem_jacobian
01173   // or numerical_side_jacobian
01174 
01175   DenseVector<Number> original_residual(context.get_elem_residual());
01176   DenseVector<Number> backwards_residual(context.get_elem_residual());
01177   DenseMatrix<Number> numeric_jacobian(context.get_elem_jacobian());
01178 #ifdef DEBUG
01179   DenseMatrix<Number> old_jacobian(context.get_elem_jacobian());
01180 #endif
01181 
01182   Real numerical_point_h = 0.;
01183   if (_mesh_sys == this)
01184     numerical_point_h = numerical_jacobian_h * context.get_elem().hmin();
01185 
01186   for (unsigned int j = 0; j != context.get_dof_indices().size(); ++j)
01187     {
01188       // Take the "minus" side of a central differenced first derivative
01189       Number original_solution = context.get_elem_solution()(j);
01190       context.get_elem_solution()(j) -= numerical_jacobian_h;
01191 
01192       // Make sure to catch any moving mesh terms
01193       // FIXME - this could be less ugly
01194       Real *coord = NULL;
01195       if (_mesh_sys == this)
01196         {
01197           if (_mesh_x_var != libMesh::invalid_uint)
01198             for (unsigned int k = 0;
01199                  k != context.get_dof_indices( _mesh_x_var ).size(); ++k)
01200               if (context.get_dof_indices( _mesh_x_var )[k] ==
01201                   context.get_dof_indices()[j])
01202                 coord = &(context.get_elem().point(k)(0));
01203           if (_mesh_y_var != libMesh::invalid_uint)
01204             for (unsigned int k = 0;
01205                  k != context.get_dof_indices( _mesh_y_var ).size(); ++k)
01206               if (context.get_dof_indices( _mesh_y_var )[k] ==
01207                   context.get_dof_indices()[j])
01208                 coord = &(context.get_elem().point(k)(1));
01209           if (_mesh_z_var != libMesh::invalid_uint)
01210             for (unsigned int k = 0;
01211                  k != context.get_dof_indices( _mesh_z_var ).size(); ++k)
01212               if (context.get_dof_indices( _mesh_z_var )[k] ==
01213                   context.get_dof_indices()[j])
01214                 coord = &(context.get_elem().point(k)(2));
01215         }
01216       if (coord)
01217         {
01218           // We have enough information to scale the perturbations
01219           // here appropriately
01220           context.get_elem_solution()(j) = original_solution - numerical_point_h;
01221           *coord = libmesh_real(context.get_elem_solution()(j));
01222         }
01223 
01224       context.get_elem_residual().zero();
01225       ((*time_solver).*(res))(false, context);
01226 #ifdef DEBUG
01227       libmesh_assert_equal_to (old_jacobian, context.get_elem_jacobian());
01228 #endif
01229       backwards_residual = context.get_elem_residual();
01230 
01231       // Take the "plus" side of a central differenced first derivative
01232       context.get_elem_solution()(j) = original_solution + numerical_jacobian_h;
01233       if (coord)
01234         {
01235           context.get_elem_solution()(j) = original_solution + numerical_point_h;
01236           *coord = libmesh_real(context.get_elem_solution()(j));
01237         }
01238       context.get_elem_residual().zero();
01239       ((*time_solver).*(res))(false, context);
01240 #ifdef DEBUG
01241       libmesh_assert_equal_to (old_jacobian, context.get_elem_jacobian());
01242 #endif
01243 
01244       context.get_elem_solution()(j) = original_solution;
01245       if (coord)
01246         {
01247           *coord = libmesh_real(context.get_elem_solution()(j));
01248           for (unsigned int i = 0; i != context.get_dof_indices().size(); ++i)
01249             {
01250               numeric_jacobian(i,j) =
01251                 (context.get_elem_residual()(i) - backwards_residual(i)) /
01252                 2. / numerical_point_h;
01253             }
01254         }
01255       else
01256         {
01257           for (unsigned int i = 0; i != context.get_dof_indices().size(); ++i)
01258             {
01259               numeric_jacobian(i,j) =
01260                 (context.get_elem_residual()(i) - backwards_residual(i)) /
01261                 2. / numerical_jacobian_h;
01262             }
01263         }
01264     }
01265 
01266   context.get_elem_residual() = original_residual;
01267   context.get_elem_jacobian() = numeric_jacobian;
01268 }
01269 
01270 
01271 
01272 void FEMSystem::numerical_elem_jacobian (FEMContext &context) const
01273 {
01274   START_LOG("numerical_elem_jacobian()", "FEMSystem");
01275   this->numerical_jacobian(&TimeSolver::element_residual, context);
01276   STOP_LOG("numerical_elem_jacobian()", "FEMSystem");
01277 }
01278 
01279 
01280 
01281 void FEMSystem::numerical_side_jacobian (FEMContext &context) const
01282 {
01283   START_LOG("numerical_side_jacobian()", "FEMSystem");
01284   this->numerical_jacobian(&TimeSolver::side_residual, context);
01285   STOP_LOG("numerical_side_jacobian()", "FEMSystem");
01286 }
01287 
01288 
01289 
01290 void FEMSystem::numerical_nonlocal_jacobian (FEMContext &context) const
01291 {
01292   START_LOG("numerical_nonlocal_jacobian()", "FEMSystem");
01293   this->numerical_jacobian(&TimeSolver::nonlocal_residual, context);
01294   STOP_LOG("numerical_nonlocal_jacobian()", "FEMSystem");
01295 }
01296 
01297 
01298 
01299 UniquePtr<DiffContext> FEMSystem::build_context ()
01300 {
01301   FEMContext* fc = new FEMContext(*this);
01302 
01303   DifferentiablePhysics* phys = this->get_physics();
01304 
01305   libmesh_assert (phys);
01306 
01307   // If we are solving a moving mesh problem, tell that to the Context
01308   fc->set_mesh_system(phys->get_mesh_system());
01309   fc->set_mesh_x_var(phys->get_mesh_x_var());
01310   fc->set_mesh_y_var(phys->get_mesh_y_var());
01311   fc->set_mesh_z_var(phys->get_mesh_z_var());
01312 
01313   fc->set_deltat_pointer( &deltat );
01314 
01315   // If we are solving the adjoint problem, tell that to the Context
01316   fc->is_adjoint() = this->get_time_solver().is_adjoint();
01317 
01318   return UniquePtr<DiffContext>(fc);
01319 }
01320 
01321 
01322 
01323 void FEMSystem::init_context(DiffContext &c)
01324 {
01325   // Parent::init_context(c);  // may be a good idea in derived classes
01326 
01327   // Although we do this in DiffSystem::build_context() and
01328   // FEMSystem::build_context() as well, we do it here just to be
01329   // extra sure that the deltat pointer gets set.  Since the
01330   // intended behavior is for classes derived from FEMSystem to
01331   // call Parent::init_context() in their own init_context()
01332   // overloads, we can ensure that those classes get the correct
01333   // deltat pointers even if they have different build_context()
01334   // overloads.
01335   c.set_deltat_pointer ( &deltat );
01336 
01337   FEMContext &context = cast_ref<FEMContext&>(c);
01338 
01339   // Make sure we're prepared to do mass integration
01340   for (unsigned int var = 0; var != this->n_vars(); ++var)
01341     if (this->get_physics()->is_time_evolving(var))
01342       {
01343         // Request shape functions based on FEType
01344         switch( FEInterface::field_type( this->variable_type(var) ) )
01345           {
01346           case( TYPE_SCALAR ):
01347             {
01348               FEBase* elem_fe = NULL;
01349               context.get_element_fe(var, elem_fe);
01350               elem_fe->get_JxW();
01351               elem_fe->get_phi();
01352             }
01353             break;
01354           case( TYPE_VECTOR ):
01355             {
01356               FEGenericBase<RealGradient>* elem_fe = NULL;
01357               context.get_element_fe(var, elem_fe);
01358               elem_fe->get_JxW();
01359               elem_fe->get_phi();
01360             }
01361             break;
01362           default:
01363             libmesh_error_msg("Unrecognized field type!");
01364           }
01365       }
01366 }
01367 
01368 
01369 
01370 void FEMSystem::mesh_position_get()
01371 {
01372   // This function makes no sense unless we've already picked out some
01373   // variable(s) to reflect mesh position coordinates
01374   if (!_mesh_sys)
01375     libmesh_error_msg("_mesh_sys was NULL!");
01376 
01377   // We currently assume mesh variables are in our own system
01378   if (_mesh_sys != this)
01379     libmesh_not_implemented();
01380 
01381   // Loop over every active mesh element on this processor
01382   const MeshBase& mesh = this->get_mesh();
01383 
01384   MeshBase::const_element_iterator el =
01385     mesh.active_local_elements_begin();
01386   const MeshBase::const_element_iterator end_el =
01387     mesh.active_local_elements_end();
01388 
01389   UniquePtr<DiffContext> con = this->build_context();
01390   FEMContext &_femcontext = cast_ref<FEMContext&>(*con);
01391   this->init_context(_femcontext);
01392 
01393   // Get the solution's mesh variables from every element
01394   for ( ; el != end_el; ++el)
01395     {
01396       _femcontext.pre_fe_reinit(*this, *el);
01397 
01398       _femcontext.elem_position_get();
01399 
01400       if (_mesh_x_var != libMesh::invalid_uint)
01401         this->solution->insert(_femcontext.get_elem_solution(_mesh_x_var),
01402                                _femcontext.get_dof_indices(_mesh_x_var) );
01403       if (_mesh_y_var != libMesh::invalid_uint)
01404         this->solution->insert(_femcontext.get_elem_solution(_mesh_y_var),
01405                                _femcontext.get_dof_indices(_mesh_y_var));
01406       if (_mesh_z_var != libMesh::invalid_uint)
01407         this->solution->insert(_femcontext.get_elem_solution(_mesh_z_var),
01408                                _femcontext.get_dof_indices(_mesh_z_var));
01409     }
01410 
01411   this->solution->close();
01412 
01413   // And make sure the current_local_solution is up to date too
01414   this->System::update();
01415 }
01416 
01417 } // namespace libMesh