5 #ifndef _compressible_neo_hookean_h_ 
    6 #define _compressible_neo_hookean_h_ 
    9 #include <deal2lkit/parsed_function.h> 
   26 template <
int dim, 
int spacedim>
 
   29   typedef FEValuesCache<dim,spacedim> Scratch;
 
   30   typedef Assembly::CopyData::piDoMUSPreconditioner<dim,dim> CopyPreconditioner;
 
   31   typedef Assembly::CopyData::piDoMUSSystem<dim,dim> CopySystem;
 
   50   template<
typename Number>
 
   54                              Number &energy) 
const;
 
   56   template<
typename Number>
 
   60                      Number &energy) 
const;
 
   65                                 const std::vector<shared_ptr<MAT> >,
 
   76   mutable shared_ptr<TrilinosWrappers::PreconditionAMG>    Amg_preconditioner;
 
   77   mutable shared_ptr<TrilinosWrappers::PreconditionJacobi> Mp_preconditioner;
 
   78   mutable shared_ptr<TrilinosWrappers::PreconditionJacobi> T_preconditioner;
 
   82 template <
int dim, 
int spacedim>
 
   85       "FESystem[FE_Q(1)^d]",
 
   86       "u,u,u", 
"1", 
"1",
"1")
 
   90 template <
int dim, 
int spacedim>
 
   91 template<
typename Number>
 
   94     CopyPreconditioner &data,
 
   97   Number alpha = this->alpha;
 
   98   this->reinit(alpha, cell, fe_cache);
 
  101   auto &us = fe_cache.get_values(
"solution", 
"u", displacement, alpha);
 
  103   const unsigned int n_q_points = us.size();
 
  105   auto &JxW = fe_cache.get_JxW_values();
 
  107   for (
unsigned int q=0; q<n_q_points; ++q)
 
  111       energy += 0.5*(u*u)*JxW[q];
 
  115 template <
int dim, 
int spacedim>
 
  116 template<
typename Number>
 
  120     Number &energy)
 const 
  122   Number alpha = this->alpha;
 
  123   this->reinit(alpha, cell, fe_cache);
 
  126   auto &us = fe_cache.get_values(
"solution", 
"u", displacement, alpha);
 
  127   auto &us_dot = fe_cache.get_values(
"solution_dot", 
"u_dot", displacement, alpha);
 
  128   auto &Fs = fe_cache.get_deformation_gradients(
"solution", 
"Fu", displacement, alpha);
 
  130   const unsigned int n_q_points = us.size();
 
  132   auto &JxW = fe_cache.get_JxW_values();
 
  134   for (
unsigned int q=0; q<n_q_points; ++q)
 
  142       Number Ic = trace(C);
 
  143       Number J = determinant(F);
 
  144       Number lnJ = std::log (J);
 
  146       Number psi = (
mu/2.)*(Ic-dim) - 
mu*lnJ + (lambda/2.)*(lnJ)*(lnJ);
 
  148       energy += (psi)*JxW[q];
 
  155 template <
int dim, 
int spacedim>
 
  159   this->add_parameter(prm, &E, 
"Young's modulus", 
"10.0", 
Patterns::Double(0.0));
 
  160   this->add_parameter(prm, &nu, 
"Poisson's ratio", 
"0.3", 
Patterns::Double(0.0));
 
  163 template <
int dim, 
int spacedim>
 
  167   mu = E/(2.0*(1.+nu));
 
  168   lambda = (E *nu)/((1.+nu)*(1.-2.*nu));
 
  171 template <
int dim, 
int spacedim>
 
  176     const std::vector<shared_ptr<MAT> >,
 
  181   std::vector<std::vector<bool> > constant_modes;
 
  197   Amg_preconditioner->initialize (preconditioner_matrix.
block(0,0),
 
  201   auto A  = linear_operator< TrilinosWrappers::MPI::Vector >( matrix.
block(0,0) );
 
  210   system_op  = block_operator<1, 1, VEC >({{
 
  218   prec_op = block_operator<1, 1, VEC >({{
 
Definition: conservative.h:20
double aggregation_threshold
std::vector< std::vector< bool > > constant_modes
void declare_parameters(ParameterHandler &prm)
Definition: compressible_neo_hookean.h:63
ActiveSelector::active_cell_iterator active_cell_iterator
unsigned int smoother_sweeps
Definition: compressible_neo_hookean.h:14
CompressibleNeoHookeanInterface()
Definition: compressible_neo_hookean.h:54
void preconditioner_energy(const typename DoFHandler< dim, spacedim >::active_cell_iterator &, Scratch &, CopyPreconditioner &, Number &energy) const 
Definition: compressible_neo_hookean.h:92
LinearOperator< typename Solver::vector_type, typename Solver::vector_type > inverse_operator(const LinearOperator< typename Solver::vector_type, typename Solver::vector_type > &op, Solver &solver, const Preconditioner &preconditioner)
void compute_system_operators(const DoFHandler< dim, spacedim > &, const std::vector< shared_ptr< LATrilinos::BlockMatrix > >, LinearOperator< LATrilinos::VectorType > &, LinearOperator< LATrilinos::VectorType > &) const 
Definition: compressible_neo_hookean.h:130
BlockType & block(const unsigned int row, const unsigned int column)
~CompressibleNeoHookeanInterface()
Definition: compressible_neo_hookean.h:38
void parse_parameters_call_back()
Definition: compressible_neo_hookean.h:71
void system_energy(const typename DoFHandler< dim, spacedim >::active_cell_iterator &, Scratch &, CopySystem &, Number &energy) const 
Definition: compressible_neo_hookean.h:117
const FiniteElement< dim, spacedim > & get_fe() const 
bool higher_order_elements