$extrastylesheet
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 #ifndef LIBMESH_FEM_SYSTEM_H 00021 #define LIBMESH_FEM_SYSTEM_H 00022 00023 // Local Includes 00024 #include "libmesh/diff_system.h" 00025 #include "libmesh/fem_physics.h" 00026 00027 // C++ includes 00028 #include <cstddef> 00029 00030 namespace libMesh 00031 { 00032 00033 // Forward Declarations 00034 class DiffContext; 00035 class FEMContext; 00036 00037 00052 // ------------------------------------------------------------ 00053 // FEMSystem class definition 00054 00055 class FEMSystem : public DifferentiableSystem, 00056 public FEMPhysics 00057 { 00058 public: 00059 00064 FEMSystem (EquationSystems& es, 00065 const std::string& name, 00066 const unsigned int number); 00067 00071 virtual ~FEMSystem (); 00072 00076 typedef FEMSystem sys_type; 00077 00081 typedef DifferentiableSystem Parent; 00082 00087 virtual void clear (); 00088 00094 virtual void assembly (bool get_residual, bool get_jacobian, 00095 bool apply_heterogeneous_constraints = false); 00096 00105 virtual void solve (); 00106 00111 void mesh_position_get(); 00112 00117 void mesh_position_set(); 00118 00127 virtual UniquePtr<DiffContext> build_context(); 00128 00129 /* 00130 * Prepares the result of a build_context() call for use. 00131 * 00132 * Most FEMSystem-based problems will need to reimplement this in order to 00133 * call FE::get_*() as their particular physics requires. 00134 */ 00135 virtual void init_context(DiffContext &); 00136 00141 virtual void postprocess (); 00142 00151 virtual void assemble_qoi 00152 (const QoISet& indices = QoISet()); 00153 00161 virtual void assemble_qoi_derivative 00162 (const QoISet &qoi_indices = QoISet(), 00163 bool include_liftfunc = true, 00164 bool apply_constraints = true); 00165 00172 bool fe_reinit_during_postprocess; 00173 00179 Real numerical_jacobian_h; 00180 00194 Real verify_analytic_jacobians; 00195 00199 typedef bool (TimeSolver::*TimeSolverResPtr)(bool, DiffContext&); 00200 00205 void numerical_jacobian (TimeSolverResPtr res, 00206 FEMContext &context) const; 00207 00213 void numerical_elem_jacobian (FEMContext &context) const; 00214 00220 void numerical_side_jacobian (FEMContext &context) const; 00221 00227 void numerical_nonlocal_jacobian (FEMContext &context) const; 00228 00229 protected: 00234 virtual void init_data (); 00235 }; 00236 00237 00238 } // namespace libMesh 00239 00240 00241 #endif // LIBMESH_FEM_SYSTEM_H