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::ComputeResidualandJacobian,
31  public NonlinearImplicitSystem::ComputeBounds,
32  public System::Initialization
33 {
34 public:
38  JR(EquationSystems & eqSys,
39  const std::string & name,
40  const unsigned int number);
41 
42  void initialize();
43 
47  static Number InitialDensityBall(const Point & p,
48  const Parameters & parameters,
49  const std::string &,
50  const std::string &);
51 
52  static Number InitialDensityRod(const Point & p,
53  const Parameters & parameters,
54  const std::string &,
55  const std::string &);
56 
57  static Number InitialDensityStrip(const Point & p,
58  const Parameters & parameters,
59  const std::string &,
60  const std::string &);
61 
62  static Gradient InitialGradientZero(const Point &,
63  const Parameters &,
64  const std::string &,
65  const std::string &);
66 
74 
78  void bounds(NumericVector<Number> & XL,
81 
82 private:
84 };
85 
86 
87 #endif // BIHARMONIC_JR_H
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:39
Biharmonic & _biharmonic
Definition: biharmonic_jr.h:83
This is the EquationSystems class.
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:63
TransientSystem< NonlinearImplicitSystem > TransientNonlinearImplicitSystem
Biharmonic&#39;s friend class definition.
Definition: biharmonic_jr.h:29
Numeric vector.
Definition: dof_map.h:66
The Biharmonic class encapsulates most of the data structures necessary to calculate the biharmonic r...
Definition: biharmonic.h:46
void residual_and_jacobian(const NumericVector< Number > &u, NumericVector< Number > *R, SparseMatrix< Number > *J, NonlinearImplicitSystem &)
The residual and Jacobian assembly function for the Biharmonic system.
Generic sparse matrix.
Definition: dof_map.h:65
static Gradient InitialGradientZero(const Point &, const Parameters &, const std::string &, const std::string &)
This is the base class for classes which contain information related to any physical process that mig...
Definition: system.h:76
NumberVectorValue Gradient
This class provides a specific system class.
static Number InitialDensityRod(const Point &p, const Parameters &parameters, const std::string &, const std::string &)
static Number InitialDensityStrip(const Point &p, const Parameters &parameters, const std::string &, const std::string &)
void bounds(NumericVector< Number > &XL, NumericVector< Number > &XU, NonlinearImplicitSystem &)
Function defining the bounds of the Biharmonic system.
static Number InitialDensityBall(const Point &p, const Parameters &parameters, const std::string &, const std::string &)
Static functions to be used for initialization.
Parameters parameters
Data structure holding arbitrary parameters.
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38