$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_NEWMARK_SYSTEM_H 00021 #define LIBMESH_NEWMARK_SYSTEM_H 00022 00023 // Local Includes 00024 #include "libmesh/linear_implicit_system.h" 00025 00026 // C++ includes 00027 00028 namespace libMesh 00029 { 00030 00031 // Forward Declarations 00032 00033 00051 // ------------------------------------------------------------ 00052 // NewmarkSystem class definition 00053 00054 class NewmarkSystem : public LinearImplicitSystem 00055 { 00056 public: 00057 00062 NewmarkSystem (EquationSystems& es, 00063 const std::string& name, 00064 const unsigned int number); 00065 00069 ~NewmarkSystem (); 00070 00071 00075 typedef NewmarkSystem sys_type; 00076 00081 virtual void clear (); 00082 00087 virtual void reinit (); 00088 00093 virtual void assemble (); 00094 00099 virtual std::string system_type () const { return "Newmark"; } 00100 00101 00102 //--------------------------------------------------------- 00103 // These members are specific to the Newmark system 00104 // 00105 00109 void initial_conditions (); 00110 00115 void compute_matrix (); 00116 00120 void update_rhs (); 00121 00125 void update_u_v_a (); 00126 00132 void set_newmark_parameters (const Real delta_T = _default_timestep, 00133 const Real alpha = _default_alpha, 00134 const Real delta = _default_delta); 00135 00136 protected: 00137 00138 00139 private: 00140 00144 Real _a_0; 00145 Real _a_1; 00146 Real _a_2; 00147 Real _a_3; 00148 Real _a_4; 00149 Real _a_5; 00150 Real _a_6; 00151 Real _a_7; 00152 00156 bool _finished_assemble; 00157 00161 static const Real _default_alpha; 00162 00166 static const Real _default_delta; 00167 00171 static const Real _default_timestep; 00172 00173 }; 00174 00175 00176 } // namespace libMesh 00177 00178 00179 // ------------------------------------------------------------ 00180 // NewmarkSystem inline methods 00181 00182 00183 00184 #endif // LIBMESH_NEWMARK_SYSTEM_H