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 
17 #ifdef LIBMESH_ENABLE_AMR
19 #endif
20 
46 {
47 public:
48  // Initial state enumeration
50  ROD = 1,
51  BALL = 2};
52 
53  // Free energy enumeration
58 
67 
68 
73  {
74  // delete _meshRefinement;
75  }
76 
77 
78  // Misc. getters
79  bool verbose() { return _verbose; }
80  Real dt0() { return _dt0; }
81  Real dt() { return _dt; }
82 
83 
84  // Public interface functions
85  void viewParameters();
86  void init();
87  void step(const Real & dt = -1.0);
88  void output(int timestep, const Real & t, Real & o_t, bool force = false);
89  void run();
90 
91 private:
92  unsigned int _dim, _N;
98  bool _verbose;
104  //
105  std::string _ofile_base, _ofile;
106  std::unique_ptr<ExodusII_IO> _exio;
108  int _o_count;
109  //
110  friend class JR;
111  class JR; // forward
113  JR * _jr;
114 };
115 
116 #endif // BIHARMONIC_H
void run()
Definition: biharmonic.C:289
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:105
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
std::unique_ptr< ExodusII_IO > _exio
Definition: biharmonic.h:106
bool _verbose
Definition: biharmonic.h:98
Biharmonic(ReplicatedMesh &mesh)
Constructor retrieves command-line options, setting defaults, if necessary.
Definition: biharmonic.C:16
Real _cnWeight
Definition: biharmonic.h:103
MeshBase & mesh
void init()
Initialize all the systems.
Definition: biharmonic.C:212
Real _kappa
Definition: biharmonic.h:93
Real _theta_c
Definition: biharmonic.h:93
void output(int timestep, const Real &t, Real &o_t, bool force=false)
Definition: biharmonic.C:253
The Biharmonic class encapsulates most of the data structures necessary to calculate the biharmonic r...
Definition: biharmonic.h:45
Point _initialCenter
Definition: biharmonic.h:100
Real _tol
Definition: biharmonic.h:94
Real dt0()
Definition: biharmonic.h:80
bool _netforce
Definition: biharmonic.h:95
int _o_count
Definition: biharmonic.h:108
void viewParameters()
Definition: biharmonic.C:165
bool _cahn_hillard
Definition: biharmonic.h:95
~Biharmonic()
Destructor.
Definition: biharmonic.h:72
Implements (adaptive) mesh refinement algorithms for a MeshBase.
bool verbose()
Definition: biharmonic.h:79
FreeEnergyEnum _energy
Definition: biharmonic.h:96
int _log_truncation
Definition: biharmonic.h:97
unsigned int _dim
Definition: biharmonic.h:92
Real dt()
Definition: biharmonic.h:81
void step(const Real &dt=-1.0)
Definition: biharmonic.C:232
Real _initialWidth
Definition: biharmonic.h:101
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool _growth
Definition: biharmonic.h:95
Real _o_dt
Definition: biharmonic.h:107
unsigned int _N
Definition: biharmonic.h:92
ReplicatedMesh & _mesh
Definition: biharmonic.h:111
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:39
std::string _ofile_base
Definition: biharmonic.h:105
InitialStateEnum _initialState
Definition: biharmonic.h:99
Real _theta
Definition: biharmonic.h:93
bool _degenerate
Definition: biharmonic.h:95