Alamo
CahnHilliard.H
Go to the documentation of this file.
1//
2// This implements a basic two-phase field model with Cahn-Hilliard kinetics.
3//
4// The free energy is
5//
6// .. math::
7//
8// F[\eta] = \int_\Omega \Big[\frac{1}{4}(\eta^2 - 1)^2 +
9// \frac{1}{2}\gamma |\nabla\eta|^2\Big] d\mathbf{x}
10//
11// The corresponding governing equation under conservative kinetics is
12//
13// .. math::
14//
15// \frac{\partial\eta}{\partial t} = L\nabla^2\Big(\eta^3 - \eta - \gamma\nabla^2\eta\Big)
16//
17// which is integrated using a forward Euler scheme.
18//
19// This is tested in :ref:`CahnHilliard`
20
21#ifndef INTEGRATOR_CAHNHILLIARD_H
22#define INTEGRATOR_CAHNHILLIARD_H
23
24#include <AMReX.H>
25#include <AMReX_MLMG.H>
26
27#include "IC/IC.H"
28#include "BC/BC.H"
29#include "IO/ParmParse.H"
31
32namespace Integrator
33{
34
36{
37public:
38 static constexpr const char* name = "cahnhilliard";
39
40 /// Basic constructor (don't use)
42
43 /// Destroy pointers defined in Parse
45
46 /// Use this constructor
48 { Parse(*this, pp); }
49
50 /// Scan input values and initialize fields
51 static void Parse(CahnHilliard& value, IO::ParmParse& pp);
52
53protected:
54
55 /// Set values in fields
56 void Initialize (int lev) override;
57 /// Integrate eta over one timestep on lev
58 void Advance (int lev, Set::Scalar time, Set::Scalar dt) override;
59 /// spectral version of advance
60 void AdvanceSpectral (int lev, Set::Scalar time, Set::Scalar dt);
61 /// realspace version of advance
62 void AdvanceReal (int lev, Set::Scalar time, Set::Scalar dt);
63 /// Mark any cells that need to be refined
64 void TagCellsForRefinement (int lev, amrex::TagBoxArray& tags, amrex::Real time, int ngrow) override;
65
66private:
67
68 Set::Field<Set::Scalar> etanew_mf; /// The new value for eta this timestep
69 Set::Field<Set::Scalar> etaold_mf; /// Last timestep's value for eta
70 Set::Field<Set::Scalar> intermediate; /// Intermediate field used for CH kinetics
71
72 BC::BC<Set::Scalar> *bc; /// eta's bc object
73 IC::IC<Set::Scalar> *ic; /// eta's ic object
74
76 Set::Scalar L = NAN;
78
79 std::string method;
80};
81}
82#endif
Definition BC.H:42
Pure abstract IC object from which all other IC objects inherit.
Definition IC.H:23
CahnHilliard()
Basic constructor (don't use)
IC::IC< Set::Scalar > * ic
eta's bc object
void AdvanceReal(int lev, Set::Scalar time, Set::Scalar dt)
realspace version of advance
Set::Scalar gamma
eta's ic object
Set::Field< Set::Scalar > intermediate
Last timestep's value for eta.
void Initialize(int lev) override
Set values in fields.
Set::Field< Set::Scalar > etaold_mf
The new value for eta this timestep.
Set::Field< Set::Scalar > etanew_mf
CahnHilliard(IO::ParmParse &pp)
Use this constructor.
void TagCellsForRefinement(int lev, amrex::TagBoxArray &tags, amrex::Real time, int ngrow) override
Mark any cells that need to be refined.
Set::Scalar refinement_threshold
static constexpr const char * name
BC::BC< Set::Scalar > * bc
Intermediate field used for CH kinetics.
void AdvanceSpectral(int lev, Set::Scalar time, Set::Scalar dt)
spectral version of advance
~CahnHilliard()
Destroy pointers defined in Parse.
void Advance(int lev, Set::Scalar time, Set::Scalar dt) override
Integrate eta over one timestep on lev.
static void Parse(CahnHilliard &value, IO::ParmParse &pp)
Scan input values and initialize fields.
amrex::Vector< amrex::Real > dt
Timesteps for each level of refinement.
Definition Integrator.H:381
Collection of numerical integrator objects.
Definition AllenCahn.H:41
amrex::Real Scalar
Definition Base.H:19
Eigen::Matrix< amrex::Real, AMREX_SPACEDIM, 1 > Vector
Definition Base.H:20