Go to the documentation of this file.    1 #include <libmesh/dof_map.h> 
    2 #include <libmesh/fem_system.h> 
    3 #include <libmesh/newton_solver.h> 
    4 #include <libmesh/enum_solver_type.h> 
    5 #include <libmesh/enum_preconditioner_type.h> 
    6 #include <libmesh/parallel.h> 
    7 #include <libmesh/auto_ptr.h>  
   15 template<
typename TimeSolverType>
 
   26   template<
typename SystemType>
 
   32     SystemType & system = es.
add_system<SystemType>(
"ScalarSystem");
 
   34     system.time_solver = libmesh_make_unique<TimeSolverType>(system);
 
   38     DiffSolver & solver = *(system.time_solver->diff_solver().get());
 
   43     NewtonSolver & newton = cast_ref<NewtonSolver &>(solver);
 
   50     system.deltat = deltat;
 
   52     TimeSolverType * time_solver = cast_ptr<TimeSolverType *>(system.time_solver.get());
 
   53     this->aux_time_solver_init(*time_solver);
 
   59     std::vector<dof_id_type> solution_index;
 
   60     solution_index.push_back(0);
 
   61     const bool has_solution = system.get_dof_map().all_semilocal_indices(solution_index);
 
   63     for (
unsigned int t_step=0; t_step != n_timesteps; ++t_step)
 
   66         system.time_solver->advance_timestep();
 
   72             Number exact_soln = system.u(system.time);
 
   73             rel_error =  
std::abs((exact_soln - (*system.solution)(0))/exact_soln);
 
   75         system.comm().max(rel_error);
 
   78         LIBMESH_ASSERT_FP_EQUAL( 0.0,
 
   80                                  std::numeric_limits<Real>::epsilon()*10 );
 
   95                              const std::string & name_in,
 
   96                              const unsigned int number_in)
 
  113     this->time_evolving(_u_var,1);
 
  121     FEMContext & c = cast_ref<FEMContext &>(context);
 
  124     const unsigned int n_u_dofs =
 
  128     for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
  130         Number Fval = this->F(c,qp);
 
  132         for (
unsigned int i=0; i != n_u_dofs; i++)
 
  136     return request_jacobian;
 
  143     FEMContext & c = cast_ref<FEMContext &>(context);
 
  147     const unsigned int n_u_dofs =
 
  151     for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
  156         Number Mval = this->M(c,qp);
 
  158         for (
unsigned int i=0; i != n_u_dofs; i++)
 
  162             if (request_jacobian)
 
  163               for (
unsigned int j=0; j != n_u_dofs; j++)
 
  168     return request_jacobian;
 
  185                                                    const std::string & name_in,
 
  186                                                    const unsigned int number_in)
 
  193     this->time_evolving(_u_var,2);
 
  204     FEMContext & c = cast_ref<FEMContext &>(context);
 
  208     const unsigned int n_u_dofs =
 
  212     for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
  217         Number Cval = this->C(c,qp);
 
  219         for (
unsigned int i=0; i != n_u_dofs; i++)
 
  223             if (request_jacobian)
 
  224               for (
unsigned int j=0; j != n_u_dofs; j++)
 
  229     return request_jacobian;
 
  235     FEMContext & c = cast_ref<FEMContext &>(context);
 
  239     const unsigned int n_u_dofs =
 
  243     for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
  248         Number Mval = this->M(c,qp);
 
  250         for (
unsigned int i=0; i != n_u_dofs; i++)
 
  254             if (request_jacobian)
 
  255               for (
unsigned int j=0; j != n_u_dofs; j++)
 
  260     return request_jacobian;
 
  275                                                   const std::string & name_in,
 
  276                                                   const unsigned int number_in)
 
  284     FEMContext & c = cast_ref<FEMContext &>(context);
 
  286     unsigned int v_var = this->get_second_order_dot_var(_u_var);
 
  290     const unsigned int n_u_dofs =
 
  294     for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
  296         Number Fval = this->F(c,qp);
 
  298         for (
unsigned int i=0; i != n_u_dofs; i++)
 
  302     return request_jacobian;
 
  309     FEMContext & c = cast_ref<FEMContext &>(context);
 
  311     unsigned int v_var = this->get_second_order_dot_var(_u_var);
 
  317     const unsigned int n_u_dofs =
 
  321     for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
  326         Number Cval = this->C(c,qp);
 
  328         for (
unsigned int i=0; i != n_u_dofs; i++)
 
  332             if (request_jacobian)
 
  333               for (
unsigned int j=0; j != n_u_dofs; j++)
 
  338     return request_jacobian;
 
  344     FEMContext & c = cast_ref<FEMContext &>(context);
 
  346     unsigned int v_var = this->get_second_order_dot_var(_u_var);
 
  351     const unsigned int n_u_dofs =
 
  355     for (
unsigned int qp=0; qp != n_qpoints; qp++)
 
  360         Number Mval = this->M(c,qp);
 
  362         for (
unsigned int i=0; i != n_u_dofs; i++)
 
  366             if (request_jacobian)
 
  367               for (
unsigned int j=0; j != n_u_dofs; j++)
 
  372     return request_jacobian;
 
 
const DenseVector< Number > & get_elem_residual() const
Const accessor for element residual.
 
const QBase & get_element_qrule() const
Accessor for element interior quadrature rule for the dimension of the current _elem.
 
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
 
virtual bool mass_residual(bool request_jacobian, DiffContext &context) override
Note the nonlinear residual is F(u)-M(u)
 
virtual System & add_system(const std::string &system_type, const std::string &name)
Add the system of type system_type named name to the systems array.
 
virtual bool element_time_derivative(bool request_jacobian, DiffContext &context) override
Note the nonlinear residual is M(u)
 
const std::vector< dof_id_type > & get_dof_indices() const
Accessor for element dof indices.
 
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
 
Defines a dense submatrix for use in Finite Element-type computations.
 
unsigned int n_points() const
 
FirstOrderScalarSystemBase(EquationSystems &es, const std::string &name_in, const unsigned int number_in)
 
The libMesh namespace provides an interface to certain functionality in the library.
 
FEMSystem-based class for testing of TimeSolvers using first order SCALARs.
 
This is a generic class that defines a solver to handle ImplicitSystem classes, including NonlinearIm...
 
void interior_accel(unsigned int var, unsigned int qp, OutputType &u) const
 
This class provides a specific system class.
 
SecondOrderScalarSystemFirstOrderTimeSolverBase(EquationSystems &es, const std::string &name_in, const unsigned int number_in)
 
void interior_rate(unsigned int var, unsigned int qp, OutputType &u) const
 
SecondOrderScalarSystemSecondOrderTimeSolverBase(EquationSystems &es, const std::string &name_in, const unsigned int number_in)
 
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
 
void set_solver_type(const SolverType st)
Sets the type of solver to use.
 
Real relative_step_tolerance
 
libMesh::Parallel::Communicator * TestCommWorld
 
virtual void init()
Initialize all the systems.
 
virtual void aux_time_solver_init(TimeSolverType &)
 
FEMSystem-based class for testing of TimeSolvers using second order SCALARs.
 
void run_test_with_exact_soln(Real deltat, unsigned int n_timesteps)
 
This class defines a solver which uses the default libMesh linear solver in a quasiNewton method to h...
 
Real get_elem_solution_rate_derivative() const
The derivative of the current elem_solution_rate w.r.t.
 
This class provides all data required for a physics package (e.g.
 
This is the EquationSystems class.
 
Real relative_residual_tolerance
 
virtual bool damping_residual(bool request_jacobian, DiffContext &context) override
Note the nonlinear residual is M(u)
 
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
 
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
 
Real absolute_residual_tolerance
The DiffSolver should exit after the residual is reduced to either less than absolute_residual_tolera...
 
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
 
virtual bool mass_residual(bool request_jacobian, DiffContext &context) override
Note the nonlinear residual is F(u)-M(u)
 
LinearSolver< Number > & get_linear_solver()
 
virtual bool damping_residual(bool request_jacobian, DiffContext &context) override
Note the nonlinear residual is M(u)\ddot{u} + C(u)
 
void set_preconditioner_type(const PreconditionerType pct)
Sets the type of preconditioner to use.
 
FEMSystem-based class for testing of TimeSolvers using second order SCALARs.
 
virtual bool element_time_derivative(bool request_jacobian, DiffContext &context) override
Note the nonlinear residual is F(u)-M(u)
 
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
 
virtual bool mass_residual(bool request_jacobian, DiffContext &context) override
Note the nonlinear residual is F(u)-M(u)
 
Real get_elem_solution_accel_derivative() const
The derivative of the current elem_solution_accel w.r.t.
 
This class provides all data required for a physics package (e.g.