libMesh
biharmonic.h
Go to the documentation of this file.
1 #ifndef BIHARMONIC_H
2 #define BIHARMONIC_H
3 
4 #include "libmesh/equation_systems.h"
5 #include "libmesh/replicated_mesh.h"
6 #include "libmesh/exodusII_io.h"
7 #include "libmesh/mesh_refinement.h"
8 
9 // Bring in bits from the libMesh namespace.
10 // Just the bits we're using, since this is a header.
13 using libMesh::Point;
14 using libMesh::Real;
16 using libMesh::UniquePtr;
17 
18 #ifdef LIBMESH_ENABLE_AMR
20 #endif
21 
47 {
48 public:
49  // Initial state enumeration
51  ROD = 1,
52  BALL = 2};
53 
54  // Free energy enumeration
59 
68 
69 
74  {
75  // delete _meshRefinement;
76  }
77 
78 
79  // Misc. getters
80  bool verbose() { return _verbose; }
81  Real dt0() { return _dt0; }
82  Real dt() { return _dt; }
83 
84 
85  // Public interface functions
86  void viewParameters();
87  void init();
88  void step(const Real & dt = -1.0);
89  void output(int timestep, const Real & t, Real & o_t, bool force = false);
90  void run();
91 
92 private:
93  unsigned int _dim, _N;
99  bool _verbose;
105  //
106  std::string _ofile_base, _ofile;
107  UniquePtr<ExodusII_IO> _exio;
109  int _o_count;
110  //
111  friend class JR;
112  class JR; // forward
114  JR * _jr;
115 };
116 
117 #endif // BIHARMONIC_H
void run()
Definition: biharmonic.C:286
This is the EquationSystems class.
The ReplicatedMesh class is derived from the MeshBase class, and is used to store identical copies of...
std::string _ofile
Definition: biharmonic.h:106
The ExodusII_IO class implements reading meshes in the ExodusII file format from Sandia National Labs...
Definition: exodusII_io.h:52
Biharmonic&#39;s friend class definition.
Definition: biharmonic_jr.h:29
bool _verbose
Definition: biharmonic.h:99
Biharmonic(ReplicatedMesh &mesh)
Constructor retrieves command-line options, setting defaults, if necessary.
Definition: biharmonic.C:13
Real _cnWeight
Definition: biharmonic.h:104
MeshBase & mesh
void init()
Initialize all the systems.
Definition: biharmonic.C:212
Real _kappa
Definition: biharmonic.h:94
Real _theta_c
Definition: biharmonic.h:94
void output(int timestep, const Real &t, Real &o_t, bool force=false)
Definition: biharmonic.C:250
The Biharmonic class encapsulates most of the data structures necessary to calculate the biharmonic r...
Definition: biharmonic.h:46
Point _initialCenter
Definition: biharmonic.h:101
Real _tol
Definition: biharmonic.h:95
Real dt0()
Definition: biharmonic.h:81
std::unique_ptr< T > UniquePtr
Definition: auto_ptr.h:46
bool _netforce
Definition: biharmonic.h:96
int _o_count
Definition: biharmonic.h:109
void viewParameters()
Definition: biharmonic.C:165
bool _cahn_hillard
Definition: biharmonic.h:96
~Biharmonic()
Destructor.
Definition: biharmonic.h:73
This is the MeshRefinement class.
bool verbose()
Definition: biharmonic.h:80
FreeEnergyEnum _energy
Definition: biharmonic.h:97
int _log_truncation
Definition: biharmonic.h:98
unsigned int _dim
Definition: biharmonic.h:93
Real dt()
Definition: biharmonic.h:82
void step(const Real &dt=-1.0)
Definition: biharmonic.C:229
Real _initialWidth
Definition: biharmonic.h:102
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool _growth
Definition: biharmonic.h:96
Real _o_dt
Definition: biharmonic.h:108
unsigned int _N
Definition: biharmonic.h:93
UniquePtr< ExodusII_IO > _exio
Definition: biharmonic.h:107
ReplicatedMesh & _mesh
Definition: biharmonic.h:112
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
std::string _ofile_base
Definition: biharmonic.h:106
InitialStateEnum _initialState
Definition: biharmonic.h:100
Real _theta
Definition: biharmonic.h:94
bool _degenerate
Definition: biharmonic.h:96