libMesh
biharmonic_jr.h
Go to the documentation of this file.
1 #ifndef BIHARMONIC_JR_H
2 #define BIHARMONIC_JR_H
3 
4 // LibMesh includes
5 #include "libmesh/transient_system.h"
6 #include "libmesh/nonlinear_solver.h"
7 
8 
9 // Example includes
10 #include "biharmonic.h"
11 
12 // Bring in bits from the libMesh namespace.
13 // Just the bits we're using, since this is a header.
15 using libMesh::Gradient;
17 using libMesh::Number;
20 using libMesh::Point;
22 using libMesh::System;
24 
25 
30  public NonlinearImplicitSystem::ComputeResidual,
31  public NonlinearImplicitSystem::ComputeResidualandJacobian,
32  public NonlinearImplicitSystem::ComputeBounds,
33  public System::Initialization
34 {
35 public:
39  JR(EquationSystems & eqSys,
40  const std::string & name,
41  const unsigned int number);
42 
43  void initialize() override;
44 
48  static Number InitialDensityBall(const Point & p,
49  const Parameters & params,
50  const std::string &,
51  const std::string &);
52 
53  static Number InitialDensityRod(const Point & p,
54  const Parameters & params,
55  const std::string &,
56  const std::string &);
57 
58  static Number InitialDensityStrip(const Point & p,
59  const Parameters & params,
60  const std::string &,
61  const std::string &);
62 
63  static Gradient InitialGradientZero(const Point &,
64  const Parameters &,
65  const std::string &,
66  const std::string &);
67 
71  void residual(const NumericVector<Number> & u,
73  NonlinearImplicitSystem & sys) override;
74 
81  NonlinearImplicitSystem &) override;
82 
86  void bounds(NumericVector<Number> & XL,
88  NonlinearImplicitSystem &) override;
89 
90 private:
92 };
93 
94 
95 #endif // BIHARMONIC_JR_H
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
Biharmonic & _biharmonic
Definition: biharmonic_jr.h:91
This is the EquationSystems class.
void initialize() override
Definition: biharmonic_jr.C:97
static Number InitialDensityStrip(const Point &p, const Parameters &params, const std::string &, const std::string &)
void residual_and_jacobian(const NumericVector< Number > &u, NumericVector< Number > *R, SparseMatrix< Number > *J, NonlinearImplicitSystem &) override
The residual and Jacobian assembly function for the Biharmonic system.
JR(EquationSystems &eqSys, const std::string &name, const unsigned int number)
Constructor.
Definition: biharmonic_jr.C:19
This class provides the ability to map between arbitrary, user-defined strings and several data types...
Definition: parameters.h:67
TransientSystem< NonlinearImplicitSystem > TransientNonlinearImplicitSystem
Biharmonic&#39;s friend class definition.
Definition: biharmonic_jr.h:29
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
Definition: vector_fe_ex5.C:44
The Biharmonic class encapsulates most of the data structures necessary to calculate the biharmonic r...
Definition: biharmonic.h:45
Generic sparse matrix.
Definition: vector_fe_ex5.C:46
static Gradient InitialGradientZero(const Point &, const Parameters &, const std::string &, const std::string &)
void bounds(NumericVector< Number > &XL, NumericVector< Number > &XU, NonlinearImplicitSystem &) override
Function defining the bounds of the Biharmonic system.
Manages consistently variables, degrees of freedom, and coefficient vectors.
Definition: system.h:96
NumberVectorValue Gradient
Manages consistently variables, degrees of freedom, coefficient vectors, matrices and non-linear solv...
static Number InitialDensityRod(const Point &p, const Parameters &params, const std::string &, const std::string &)
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
void residual(const NumericVector< Number > &u, NumericVector< Number > &R, NonlinearImplicitSystem &sys) override
The residual assembly function for the Biharmonic system.
static Number InitialDensityBall(const Point &p, const Parameters &params, const std::string &, const std::string &)
Static functions to be used for initialization.