libMesh
time_solver_test_common.h
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> // libmesh_make_unique
8 
9 #include "test_comm.h"
10 #include "libmesh_cppunit.h"
11 
12 
13 using namespace libMesh;
14 
15 template<typename TimeSolverType>
17 {
18 protected:
19 
20  // Any specialized initialization that's needed for the test
21  virtual void aux_time_solver_init( TimeSolverType & /*time_solver*/ ){}
22 
23  // Implementation for solving ODE of SystemType
24  // Note this test assumes that the time integrator gets the *exact* solution
25  // to within floating point tolerance.
26  template<typename SystemType>
27  void run_test_with_exact_soln(Real deltat, unsigned int n_timesteps)
28  {
32  SystemType & system = es.add_system<SystemType>("ScalarSystem");
33 
34  system.time_solver = libmesh_make_unique<TimeSolverType>(system);
35 
36  es.init();
37 
38  DiffSolver & solver = *(system.time_solver->diff_solver().get());
39  solver.relative_step_tolerance = std::numeric_limits<Real>::epsilon()*10;
40  solver.relative_residual_tolerance = std::numeric_limits<Real>::epsilon()*10;
41  solver.absolute_residual_tolerance = std::numeric_limits<Real>::epsilon()*10;
42 
43  NewtonSolver & newton = cast_ref<NewtonSolver &>(solver);
44 
45  // LASPACK GMRES + ILU defaults don't like these problems, so
46  // we'll use a sophisticated "just divide the scalars" solver instead.
49 
50  system.deltat = deltat;
51 
52  TimeSolverType * time_solver = cast_ptr<TimeSolverType *>(system.time_solver.get());
53  this->aux_time_solver_init(*time_solver);
54 
55  // We're going to want to check our solution, and when we run
56  // "make check" with LIBMESH_RUN='mpirun -np N" for N>1 then we'll
57  // need to keep that check in sync with the processors that are just
58  // twiddling their thumbs, not owning our mesh point.
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);
62 
63  for (unsigned int t_step=0; t_step != n_timesteps; ++t_step)
64  {
65  system.solve();
66  system.time_solver->advance_timestep();
67 
68  Real rel_error = 0;
69 
70  if (has_solution)
71  {
72  Number exact_soln = system.u(system.time);
73  rel_error = std::abs((exact_soln - (*system.solution)(0))/exact_soln);
74  }
75  system.comm().max(rel_error);
76 
77  // Using relative error for comparison, so "exact" is 0
78  LIBMESH_ASSERT_FP_EQUAL( 0.0,
79  rel_error,
80  std::numeric_limits<Real>::epsilon()*10 );
81  }
82  }
83 
84 };
85 
87 
92 {
93 public:
95  const std::string & name_in,
96  const unsigned int number_in)
97  : FEMSystem(es, name_in, number_in)
98  {}
99 
100 
102  virtual Number F( FEMContext & context, unsigned int qp ) =0;
103 
105  virtual Number M( FEMContext & context, unsigned int qp ) =0;
106 
108  virtual Number u( Real t ) =0;
109 
110  virtual void init_data () override
111  {
112  _u_var = this->add_variable ("u", FIRST, LAGRANGE);
113  this->time_evolving(_u_var,1);
115  }
116 
118  virtual bool element_time_derivative (bool request_jacobian,
119  DiffContext & context) override
120  {
121  FEMContext & c = cast_ref<FEMContext &>(context);
122  DenseSubVector<Number> & Fu = c.get_elem_residual(_u_var);
123 
124  const unsigned int n_u_dofs =
125  cast_int<unsigned int>(c.get_dof_indices(_u_var).size());
126  unsigned int n_qpoints = c.get_element_qrule().n_points();
127 
128  for (unsigned int qp=0; qp != n_qpoints; qp++)
129  {
130  Number Fval = this->F(c,qp);
131 
132  for (unsigned int i=0; i != n_u_dofs; i++)
133  Fu(i) += Fval;
134  }
135 
136  return request_jacobian;
137  }
138 
140  virtual bool mass_residual (bool request_jacobian,
141  DiffContext & context) override
142  {
143  FEMContext & c = cast_ref<FEMContext &>(context);
144  DenseSubVector<Number> & Fu = c.get_elem_residual(_u_var);
145  DenseSubMatrix<Number> & Kuu = c.get_elem_jacobian(_u_var, _u_var);
146 
147  const unsigned int n_u_dofs =
148  cast_int<unsigned int>(c.get_dof_indices(_u_var).size());
149  unsigned int n_qpoints = c.get_element_qrule().n_points();
150 
151  for (unsigned int qp=0; qp != n_qpoints; qp++)
152  {
153  libMesh::Number udot;
154  c.interior_rate(_u_var,qp,udot);
155 
156  Number Mval = this->M(c,qp);
157 
158  for (unsigned int i=0; i != n_u_dofs; i++)
159  {
160  Fu(i) -= Mval*udot;
161 
162  if (request_jacobian)
163  for (unsigned int j=0; j != n_u_dofs; j++)
164  Kuu(i,j) -= Mval*context.get_elem_solution_rate_derivative();
165  }
166  }
167 
168  return request_jacobian;
169  }
170 
171 protected:
172  unsigned int _u_var;
173 };
174 
176 
182 {
183 public:
185  const std::string & name_in,
186  const unsigned int number_in)
187  : FirstOrderScalarSystemBase(es, name_in, number_in)
188  {}
189 
190  virtual void init_data () override
191  {
192  _u_var = this->add_variable ("u", FIRST, LAGRANGE);
193  this->time_evolving(_u_var,2);
195  }
196 
198  virtual Number C( FEMContext & context, unsigned int qp ) =0;
199 
201  virtual bool damping_residual (bool request_jacobian,
202  DiffContext & context) override
203  {
204  FEMContext & c = cast_ref<FEMContext &>(context);
205  DenseSubVector<Number> & Fu = c.get_elem_residual(_u_var);
206  DenseSubMatrix<Number> & Kuu = c.get_elem_jacobian(_u_var, _u_var);
207 
208  const unsigned int n_u_dofs =
209  cast_int<unsigned int>(c.get_dof_indices(_u_var).size());
210  unsigned int n_qpoints = c.get_element_qrule().n_points();
211 
212  for (unsigned int qp=0; qp != n_qpoints; qp++)
213  {
214  libMesh::Number udot;
215  c.interior_rate(_u_var,qp,udot);
216 
217  Number Cval = this->C(c,qp);
218 
219  for (unsigned int i=0; i != n_u_dofs; i++)
220  {
221  Fu(i) += Cval*udot;
222 
223  if (request_jacobian)
224  for (unsigned int j=0; j != n_u_dofs; j++)
225  Kuu(i,j) += Cval*context.get_elem_solution_rate_derivative();
226  }
227  }
228 
229  return request_jacobian;
230  }
231 
232  virtual bool mass_residual (bool request_jacobian,
233  DiffContext & context) override
234  {
235  FEMContext & c = cast_ref<FEMContext &>(context);
236  DenseSubVector<Number> & Fu = c.get_elem_residual(_u_var);
237  DenseSubMatrix<Number> & Kuu = c.get_elem_jacobian(_u_var, _u_var);
238 
239  const unsigned int n_u_dofs =
240  cast_int<unsigned int>(c.get_dof_indices(_u_var).size());
241  unsigned int n_qpoints = c.get_element_qrule().n_points();
242 
243  for (unsigned int qp=0; qp != n_qpoints; qp++)
244  {
245  libMesh::Number uddot;
246  c.interior_accel(_u_var,qp,uddot);
247 
248  Number Mval = this->M(c,qp);
249 
250  for (unsigned int i=0; i != n_u_dofs; i++)
251  {
252  Fu(i) += Mval*uddot;
253 
254  if (request_jacobian)
255  for (unsigned int j=0; j != n_u_dofs; j++)
256  Kuu(i,j) += Mval*context.get_elem_solution_accel_derivative();
257  }
258  }
259 
260  return request_jacobian;
261  }
262 };
263 
264 
266 
272 {
273 public:
275  const std::string & name_in,
276  const unsigned int number_in)
277  : SecondOrderScalarSystemSecondOrderTimeSolverBase(es, name_in, number_in)
278  {}
279 
281  virtual bool element_time_derivative (bool request_jacobian,
282  DiffContext & context) override
283  {
284  FEMContext & c = cast_ref<FEMContext &>(context);
285 
286  unsigned int v_var = this->get_second_order_dot_var(_u_var);
287 
289 
290  const unsigned int n_u_dofs =
291  cast_int<unsigned int>(c.get_dof_indices(_u_var).size());
292  unsigned int n_qpoints = c.get_element_qrule().n_points();
293 
294  for (unsigned int qp=0; qp != n_qpoints; qp++)
295  {
296  Number Fval = this->F(c,qp);
297 
298  for (unsigned int i=0; i != n_u_dofs; i++)
299  Fv(i) += Fval;
300  }
301 
302  return request_jacobian;
303  }
304 
306  virtual bool damping_residual (bool request_jacobian,
307  DiffContext & context) override
308  {
309  FEMContext & c = cast_ref<FEMContext &>(context);
310 
311  unsigned int v_var = this->get_second_order_dot_var(_u_var);
312 
314 
315  DenseSubMatrix<Number> & Kvv = c.get_elem_jacobian(v_var, v_var);
316 
317  const unsigned int n_u_dofs =
318  cast_int<unsigned int>(c.get_dof_indices(_u_var).size());
319  unsigned int n_qpoints = c.get_element_qrule().n_points();
320 
321  for (unsigned int qp=0; qp != n_qpoints; qp++)
322  {
323  libMesh::Number udot;
324  c.interior_rate(v_var,qp,udot);
325 
326  Number Cval = this->C(c,qp);
327 
328  for (unsigned int i=0; i != n_u_dofs; i++)
329  {
330  Fv(i) += Cval*udot;
331 
332  if (request_jacobian)
333  for (unsigned int j=0; j != n_u_dofs; j++)
334  Kvv(i,j) += Cval*context.get_elem_solution_rate_derivative();
335  }
336  }
337 
338  return request_jacobian;
339  }
340 
341  virtual bool mass_residual (bool request_jacobian,
342  DiffContext & context) override
343  {
344  FEMContext & c = cast_ref<FEMContext &>(context);
345 
346  unsigned int v_var = this->get_second_order_dot_var(_u_var);
347 
349  DenseSubMatrix<Number> & Kvv = c.get_elem_jacobian(v_var, v_var);
350 
351  const unsigned int n_u_dofs =
352  cast_int<unsigned int>(c.get_dof_indices(_u_var).size());
353  unsigned int n_qpoints = c.get_element_qrule().n_points();
354 
355  for (unsigned int qp=0; qp != n_qpoints; qp++)
356  {
357  libMesh::Number uddot;
358  c.interior_accel(v_var,qp,uddot);
359 
360  Number Mval = this->M(c,qp);
361 
362  for (unsigned int i=0; i != n_u_dofs; i++)
363  {
364  Fv(i) += Mval*uddot;
365 
366  if (request_jacobian)
367  for (unsigned int j=0; j != n_u_dofs; j++)
368  Kvv(i,j) += Mval*context.get_elem_solution_accel_derivative();
369  }
370  }
371 
372  return request_jacobian;
373  }
374 };
libMesh::Number
Real Number
Definition: libmesh_common.h:195
libMesh::DiffContext::get_elem_residual
const DenseVector< Number > & get_elem_residual() const
Const accessor for element residual.
Definition: diff_context.h:249
libMesh::FEMContext::get_element_qrule
const QBase & get_element_qrule() const
Accessor for element interior quadrature rule for the dimension of the current _elem.
Definition: fem_context.h:790
libMesh::Mesh
The Mesh class is a thin wrapper, around the ReplicatedMesh class by default.
Definition: mesh.h:50
FirstOrderScalarSystemBase::mass_residual
virtual bool mass_residual(bool request_jacobian, DiffContext &context) override
Note the nonlinear residual is F(u)-M(u)
Definition: time_solver_test_common.h:140
libMesh::EquationSystems::add_system
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.
Definition: equation_systems.C:345
SecondOrderScalarSystemFirstOrderTimeSolverBase::element_time_derivative
virtual bool element_time_derivative(bool request_jacobian, DiffContext &context) override
Note the nonlinear residual is M(u)
Definition: time_solver_test_common.h:281
libMesh::DiffContext::get_dof_indices
const std::vector< dof_id_type > & get_dof_indices() const
Accessor for element dof indices.
Definition: diff_context.h:367
libMesh::FEMSystem::init_data
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
Definition: fem_system.C:843
libMesh::DenseSubMatrix
Defines a dense submatrix for use in Finite Element-type computations.
Definition: dense_submatrix.h:45
libMesh::JACOBI
Definition: enum_solver_type.h:46
libMesh::QBase::n_points
unsigned int n_points() const
Definition: quadrature.h:126
FirstOrderScalarSystemBase::FirstOrderScalarSystemBase
FirstOrderScalarSystemBase(EquationSystems &es, const std::string &name_in, const unsigned int number_in)
Definition: time_solver_test_common.h:94
libMesh::IDENTITY_PRECOND
Definition: enum_preconditioner_type.h:34
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
FirstOrderScalarSystemBase
FEMSystem-based class for testing of TimeSolvers using first order SCALARs.
Definition: time_solver_test_common.h:91
libMesh::DiffSolver
This is a generic class that defines a solver to handle ImplicitSystem classes, including NonlinearIm...
Definition: diff_solver.h:70
mesh
MeshBase & mesh
Definition: mesh_communication.C:1257
libMesh::FEMContext::interior_accel
void interior_accel(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:1349
libMesh::FEMSystem
This class provides a specific system class.
Definition: fem_system.h:53
SecondOrderScalarSystemFirstOrderTimeSolverBase::SecondOrderScalarSystemFirstOrderTimeSolverBase
SecondOrderScalarSystemFirstOrderTimeSolverBase(EquationSystems &es, const std::string &name_in, const unsigned int number_in)
Definition: time_solver_test_common.h:274
TimeSolverTestImplementation
Definition: time_solver_test_common.h:16
libMesh::FEMContext::interior_rate
void interior_rate(unsigned int var, unsigned int qp, OutputType &u) const
Definition: fem_context.C:1318
SecondOrderScalarSystemSecondOrderTimeSolverBase::SecondOrderScalarSystemSecondOrderTimeSolverBase
SecondOrderScalarSystemSecondOrderTimeSolverBase(EquationSystems &es, const std::string &name_in, const unsigned int number_in)
Definition: time_solver_test_common.h:184
std::abs
MetaPhysicL::DualNumber< T, D > abs(const MetaPhysicL::DualNumber< T, D > &in)
libMesh::LinearSolver::set_solver_type
void set_solver_type(const SolverType st)
Sets the type of solver to use.
Definition: linear_solver.h:126
libMesh::DiffSolver::relative_step_tolerance
Real relative_step_tolerance
Definition: diff_solver.h:204
TestCommWorld
libMesh::Parallel::Communicator * TestCommWorld
Definition: driver.C:111
libMesh::EquationSystems::init
virtual void init()
Initialize all the systems.
Definition: equation_systems.C:96
TimeSolverTestImplementation::aux_time_solver_init
virtual void aux_time_solver_init(TimeSolverType &)
Definition: time_solver_test_common.h:21
SecondOrderScalarSystemFirstOrderTimeSolverBase
FEMSystem-based class for testing of TimeSolvers using second order SCALARs.
Definition: time_solver_test_common.h:271
TimeSolverTestImplementation::run_test_with_exact_soln
void run_test_with_exact_soln(Real deltat, unsigned int n_timesteps)
Definition: time_solver_test_common.h:27
libMesh::NewtonSolver
This class defines a solver which uses the default libMesh linear solver in a quasiNewton method to h...
Definition: newton_solver.h:46
libMesh::DiffContext::get_elem_solution_rate_derivative
Real get_elem_solution_rate_derivative() const
The derivative of the current elem_solution_rate w.r.t.
Definition: diff_context.h:445
libMesh::DiffContext
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
libMesh::EquationSystems
This is the EquationSystems class.
Definition: equation_systems.h:74
libMesh::DenseSubVector< Number >
FirstOrderScalarSystemBase::_u_var
unsigned int _u_var
Definition: time_solver_test_common.h:172
libMesh::DiffSolver::relative_residual_tolerance
Real relative_residual_tolerance
Definition: diff_solver.h:192
SecondOrderScalarSystemFirstOrderTimeSolverBase::damping_residual
virtual bool damping_residual(bool request_jacobian, DiffContext &context) override
Note the nonlinear residual is M(u)
Definition: time_solver_test_common.h:306
FirstOrderScalarSystemBase::init_data
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
Definition: time_solver_test_common.h:110
SecondOrderScalarSystemSecondOrderTimeSolverBase::init_data
virtual void init_data() override
Initializes the member data fields associated with the system, so that, e.g., assemble() may be used.
Definition: time_solver_test_common.h:190
libMesh::DiffSolver::absolute_residual_tolerance
Real absolute_residual_tolerance
The DiffSolver should exit after the residual is reduced to either less than absolute_residual_tolera...
Definition: diff_solver.h:191
libmesh_cppunit.h
libMesh::DiffContext::get_elem_jacobian
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
Definition: diff_context.h:283
SecondOrderScalarSystemSecondOrderTimeSolverBase::mass_residual
virtual bool mass_residual(bool request_jacobian, DiffContext &context) override
Note the nonlinear residual is F(u)-M(u)
Definition: time_solver_test_common.h:232
libMesh::NewtonSolver::get_linear_solver
LinearSolver< Number > & get_linear_solver()
Definition: newton_solver.h:81
SecondOrderScalarSystemSecondOrderTimeSolverBase::damping_residual
virtual bool damping_residual(bool request_jacobian, DiffContext &context) override
Note the nonlinear residual is M(u)\ddot{u} + C(u)
Definition: time_solver_test_common.h:201
libMesh::LinearSolver::set_preconditioner_type
void set_preconditioner_type(const PreconditionerType pct)
Sets the type of preconditioner to use.
Definition: linear_solver.C:108
SecondOrderScalarSystemSecondOrderTimeSolverBase
FEMSystem-based class for testing of TimeSolvers using second order SCALARs.
Definition: time_solver_test_common.h:181
FirstOrderScalarSystemBase::element_time_derivative
virtual bool element_time_derivative(bool request_jacobian, DiffContext &context) override
Note the nonlinear residual is F(u)-M(u)
Definition: time_solver_test_common.h:118
libMesh::LAGRANGE
Definition: enum_fe_family.h:36
libMesh::MeshTools::Generation::build_point
void build_point(UnstructuredMesh &mesh, const ElemType type=INVALID_ELEM, const bool gauss_lobatto_grid=false)
A specialized build_cube() for 0D meshes.
Definition: mesh_generation.C:1462
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
test_comm.h
SecondOrderScalarSystemFirstOrderTimeSolverBase::mass_residual
virtual bool mass_residual(bool request_jacobian, DiffContext &context) override
Note the nonlinear residual is F(u)-M(u)
Definition: time_solver_test_common.h:341
libMesh::DiffContext::get_elem_solution_accel_derivative
Real get_elem_solution_accel_derivative() const
The derivative of the current elem_solution_accel w.r.t.
Definition: diff_context.h:454
libMesh::FIRST
Definition: enum_order.h:42
libMesh::FEMContext
This class provides all data required for a physics package (e.g.
Definition: fem_context.h:62