$extrastylesheet
newmark_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 // C++ includes
00021 
00022 // Local includes
00023 #include "libmesh/newmark_system.h"
00024 #include "libmesh/equation_systems.h"
00025 #include "libmesh/sparse_matrix.h"
00026 #include "libmesh/libmesh_logging.h"
00027 #include "libmesh/numeric_vector.h"
00028 
00029 namespace libMesh
00030 {
00031 
00032 
00033 
00034 
00035 // ------------------------------------------------------------
00036 // NewmarkSystem static members
00037 const Real NewmarkSystem::_default_alpha    = .25;
00038 const Real NewmarkSystem::_default_delta    = .5;
00039 const Real NewmarkSystem::_default_timestep = 1.;
00040 
00041 
00042 
00043 // ------------------------------------------------------------
00044 // NewmarkSystem implementation
00045 NewmarkSystem::NewmarkSystem (EquationSystems& es,
00046                               const std::string& name_in,
00047                               const unsigned int number_in) :
00048   LinearImplicitSystem (es, name_in, number_in),
00049   _a_0                 (1./(_default_alpha*_default_timestep*_default_timestep)),
00050   _a_1                 (_default_delta/(_default_alpha*_default_timestep)),
00051   _a_2                 (1./(_default_alpha*_default_timestep)),
00052   _a_3                 (1./(2.*_default_alpha)-1.),
00053   _a_4                 (_default_delta/_default_alpha -1.),
00054   _a_5                 (_default_timestep/2.*(_default_delta/_default_alpha-2.)),
00055   _a_6                 (_default_timestep*(1.-_default_delta)),
00056   _a_7                 (_default_delta*_default_timestep),
00057   _finished_assemble   (false)
00058 
00059 {
00060   // default values of the newmark parameters
00061   es.parameters.set<Real>("Newmark alpha") = _default_alpha;
00062   es.parameters.set<Real>("Newmark delta") = _default_delta;
00063 
00064   // time step size.
00065   // should be handled at a later stage through EquationSystems?
00066   es.parameters.set<Real>("Newmark time step") = _default_timestep;
00067 
00068   // add additional matrices and vectors that will be used in the
00069   // newmark algorithm to the data structure
00070   // functions LinearImplicitSystem::add_matrix and LinearImplicitSystem::add_vector
00071   // are used so we do not have to bother about initialization and
00072   // dof mapping
00073 
00074   // system matrices
00075   this->add_matrix ("stiffness");
00076   this->add_matrix ("damping");
00077   this->add_matrix ("mass");
00078 
00079   // load vector
00080   this->add_vector ("force");
00081 
00082   // the displacement and the time derivatives
00083   this->add_vector ("displacement");
00084   this->add_vector ("velocity");
00085   this->add_vector ("acceleration");
00086 
00087   // contributions to the rhs
00088   this->add_vector ("rhs_m");
00089   this->add_vector ("rhs_c");
00090 
00091   // results from the previous time step
00092   this->add_vector ("old_solution");
00093   this->add_vector ("old_acceleration");
00094 }
00095 
00096 
00097 
00098 NewmarkSystem::~NewmarkSystem ()
00099 {
00100   this->clear();
00101 }
00102 
00103 
00104 
00105 void NewmarkSystem::clear ()
00106 {
00107   // use parent clear this will also clear the
00108   // matrices and vectors added in the constructor
00109   LinearImplicitSystem::clear();
00110 
00111   // Get a reference to the EquationSystems
00112   EquationSystems& es =
00113     this->get_equation_systems();
00114 
00115   // default values of the newmark parameters
00116   es.parameters.set<Real>("Newmark alpha") = _default_alpha;
00117   es.parameters.set<Real>("Newmark delta") = _default_delta;
00118 
00119   // time step size.  should be handled at a later stage through EquationSystems?
00120   es.parameters.set<Real>("Newmark time step") = _default_timestep;
00121 
00122   // set bool to false
00123   _finished_assemble = false;
00124 }
00125 
00126 
00127 
00128 void NewmarkSystem::reinit ()
00129 {
00130   libmesh_not_implemented();
00131 }
00132 
00133 
00134 
00135 void NewmarkSystem::assemble ()
00136 {
00137   if (!_finished_assemble)
00138     {
00139       // prepare matrix with the help of the _dof_map,
00140       // fill with sparsity pattern
00141       LinearImplicitSystem::assemble();
00142 
00143       // compute the effective system matrix
00144       this->compute_matrix();
00145 
00146       // apply initial conditions
00147       this->initial_conditions();
00148 
00149       _finished_assemble = true;
00150     }
00151 }
00152 
00153 
00154 void NewmarkSystem::initial_conditions ()
00155 {
00156   // libmesh_assert(init_cond_fptr);
00157 
00158   // Log how long the user's matrix assembly code takes
00159   START_LOG("initial_conditions ()", "NewmarkSystem");
00160 
00161   // Set all values to 0, then
00162   // call the user-specified function for initial conditions.
00163   this->get_vector("displacement").zero();
00164   this->get_vector("velocity").zero();
00165   this->get_vector("acceleration").zero();
00166   this->user_initialization();
00167 
00168   // Stop logging the user code
00169   STOP_LOG("initial_conditions ()", "NewmarkSystem");
00170 }
00171 
00172 
00173 
00174 void NewmarkSystem::compute_matrix ()
00175 {
00176   // close the component matrices
00177   this->get_matrix ("stiffness").close();
00178   this->get_matrix ("mass"     ).close();
00179   this->get_matrix ("damping"  ).close();
00180 
00181   // close & zero the system matrix
00182   this->matrix->close (); this->matrix->zero();
00183 
00184   // add up the matrices
00185   this->matrix->add (1.,   this->get_matrix ("stiffness"));
00186   this->matrix->add (_a_0, this->get_matrix ("mass"));
00187   this->matrix->add (_a_1, this->get_matrix ("damping"));
00188 
00189 }
00190 
00191 
00192 
00193 void NewmarkSystem::update_rhs ()
00194 {
00195   START_LOG("update_rhs ()", "NewmarkSystem");
00196 
00197   // zero the rhs-vector
00198   NumericVector<Number>& the_rhs = *this->rhs;
00199   the_rhs.zero();
00200 
00201   // get writable references to some vectors
00202   NumericVector<Number>& rhs_m = this->get_vector("rhs_m");
00203   NumericVector<Number>& rhs_c = this->get_vector("rhs_c");
00204 
00205 
00206   // zero the vectors for matrix-vector product
00207   rhs_m.zero();
00208   rhs_c.zero();
00209 
00210   // compute auxiliary vectors rhs_m and rhs_c
00211   rhs_m.add(_a_0, this->get_vector("displacement"));
00212   rhs_m.add(_a_2, this->get_vector("velocity"));
00213   rhs_m.add(_a_3, this->get_vector("acceleration"));
00214 
00215   rhs_c.add(_a_1, this->get_vector("displacement"));
00216   rhs_c.add(_a_4, this->get_vector("velocity"));
00217   rhs_c.add(_a_5, this->get_vector("acceleration"));
00218 
00219   // compute rhs
00220   the_rhs.add(this->get_vector("force"));
00221   the_rhs.add_vector(rhs_m, this->get_matrix("mass"));
00222   the_rhs.add_vector(rhs_c, this->get_matrix("damping"));
00223 
00224   STOP_LOG("update_rhs ()", "NewmarkSystem");
00225 }
00226 
00227 
00228 
00229 void NewmarkSystem::update_u_v_a ()
00230 {
00231   START_LOG("update_u_v_a ()", "NewmarkSystem");
00232 
00233   // get some references for convenience
00234   const NumericVector<Number>&  solu = *this->solution;
00235 
00236   NumericVector<Number>&  disp_vec   = this->get_vector("displacement");
00237   NumericVector<Number>&  vel_vec    = this->get_vector("velocity");
00238   NumericVector<Number>&  acc_vec    = this->get_vector("acceleration");
00239   NumericVector<Number>&  old_acc    = this->get_vector("old_acceleration");
00240   NumericVector<Number>&  old_solu   = this->get_vector("old_solution");
00241 
00242   // copy data
00243   old_solu = disp_vec;
00244   disp_vec = solu;
00245   old_acc  = acc_vec;
00246 
00247   // compute the new acceleration vector
00248   acc_vec.scale(-_a_3);
00249   acc_vec.add(_a_0, disp_vec);
00250   acc_vec.add(-_a_0, old_solu);
00251   acc_vec.add(-_a_2,vel_vec);
00252 
00253   // compute the new velocity vector
00254   vel_vec.add(_a_6,old_acc);
00255   vel_vec.add(_a_7,acc_vec);
00256 
00257   STOP_LOG("update_u_v_a ()", "NewmarkSystem");
00258 }
00259 
00260 
00261 
00262 void NewmarkSystem::set_newmark_parameters (const Real delta_T,
00263                                             const Real alpha,
00264                                             const Real delta)
00265 {
00266   libmesh_assert_not_equal_to (delta_T, 0.);
00267 
00268   // Get a reference to the EquationSystems
00269   EquationSystems& es =
00270     this->get_equation_systems();
00271 
00272   // the newmark parameters
00273   es.parameters.set<Real>("Newmark alpha") = alpha;
00274   es.parameters.set<Real>("Newmark delta") = delta;
00275 
00276   // time step size.
00277   // should be handled at a later stage through EquationSystems?
00278   es.parameters.set<Real>("Newmark time step") = delta_T;
00279 
00280   // the constants for time integration
00281   _a_0 = 1./(alpha*delta_T*delta_T);
00282   _a_1 = delta/(alpha*delta_T);
00283   _a_2 = 1./(alpha*delta_T);
00284   _a_3 = 1./(2.*alpha)-1.;
00285   _a_4 = delta/alpha -1.;
00286   _a_5 = delta_T/2.*(delta/alpha-2.);
00287   _a_6 = delta_T*(1.-delta);
00288   _a_7 = delta*delta_T;
00289 }
00290 
00291 } // namespace libMesh