File physicsmodel.hxx¶
Base class for Physics Models.
Changelog:
2013-08 Ben Dudson benjamin.dudson@york.ac.uk
- Initial version
Copyright 2013 B.D.Dudson
Contact: Ben Dudson, bd512@york.ac.uk
This file is part of BOUT++.
BOUT++ is free software: you can redistribute it and/or modify it under the terms of the GNU Lesser General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version.
BOUT++ is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License along with BOUT++. If not, see http://www.gnu.org/licenses/.
Defines
-
BOUTMAIN
(ModelClass)¶ Macro to define a simple main() which creates the given model and runs it. This should be sufficient for most use cases, but a user can define their own main() function if needed.
Example
class MyModel : public PhysicsModel { .. };
BOUTMAIN(MyModel);
-
SOLVE_FOR1
(var)¶ Macro to replace solver->add, passing variable name.
-
SOLVE_FOR2
(var1, var2)¶
-
SOLVE_FOR3
(var1, var2, var3)¶
-
SOLVE_FOR4
(var1, var2, var3, var4)¶
-
SOLVE_FOR5
(var1, var2, var3, var4, var5)¶
-
SOLVE_FOR6
(var1, var2, var3, var4, var5, var6)¶
-
class
PhysicsModel
¶ - #include <physicsmodel.hxx>
Base class for physics models
Subclassed by LegacyModel
Public Types
Public Functions
-
PhysicsModel
()¶
-
virtual
~PhysicsModel
()¶
-
void
initialise
(Solver *s)¶ Initialse the model, calling the init() and postInit() methods
Note: this is usually only called by the Solver
-
int
runRHS
(BoutReal time)¶ Run the RHS function, to calculate the time derivatives
Input
The system state should be in the evolving variables
- Parameters
time
: The simulation time
The time derivatives will be put in the ddt() variables
Returns a flag: 0 indicates success, non-zero an error flag
-
bool
splitOperator
()¶ True if this model uses split operators
-
int
runDiffusive
(BoutReal time, bool linear)¶ Run the diffusive (usually implicit) part of the model
-
bool
hasPrecon
()¶ True if a preconditioner has been defined
-
int
runPrecon
(BoutReal t, BoutReal gamma, BoutReal delta)¶ Run the preconditioner. The system state should be in the evolving variables, and the vector to be solved in the ddt() variables. The result will be put in the ddt() variables.
Note: this is usually only called by the Solver
-
bool
hasJacobian
()¶ True if a Jacobian function has been defined
Protected Functions
-
virtual int
init
(bool restarting) = 0¶ This function is called once by the solver at the start of a simulation.
A valid PhysicsModel must implement this function
Variables should be read from the inputs, and the variables to be evolved should be specified.
-
int
postInit
(bool restarting)¶ Post-initialise. This reads the restart file
- Parameters
restarting
: If true, will load state from restart file
-
virtual int
rhs
(BoutReal t)¶ This function is called by the time integration solver at least once per time step.
Variables being evolved will be set by the solver before the call, and this function must calculate and set the time-derivatives.
By default this function just returns an error, which will stop the simulation.
-
virtual int
outputMonitor
(BoutReal simtime, int iter, int NOUT)¶ Implemented by user code to monitor solution at output times
-
virtual int
timestepMonitor
(BoutReal simtime, BoutReal dt)¶ Timestep monitor. If enabled by setting solver:monitor_timestep=true then this function is called every internal timestep.
-
void
setSplitOperator
(bool split = true)¶ Specify that this model is split into a convective and diffusive part.
-
void
setPrecon
(preconfunc pset)¶ Specify a preconditioner function.
-
void
setJacobian
(jacobianfunc jset)¶ Specify a Jacobian-vector multiply function.
-
void
bout_solve
(Field2D &var, const char *name)¶ Specify a variable for the solver to evolve
Note that the variable must not be destroyed (e.g. go out of scope) after this call, since a pointer to
var
is stored in the solver.- Parameters
var
: The variable to evolvename
: The name to use for variable initialisation and output
To evolve the state, the solver will set
var
, and the user-supplied rhs() function should calculate ddt(var).
-
bool
bout_constrain
(Field3D &var, Field3D &F_var, const char *name)¶ Specify a constrained variable
var
, which will be adjusted to makeF_var
equal to zero. If the solver does not support constraints then this will throw an exception- Parameters
var
: The variable the solver should modifyF_var
: The control variable, which the user will setname
: The name to use for initialisation and output
Protected Attributes
-
Solver *
solver
= {nullptr}¶ This is set by a call to initialise, and can be used by models to specify evolving variables.
-
PhysicsModelMonitor
modelMonitor
¶ write restarts and pass outputMonitor method inside a Monitor subclass
Private Members
-
bool
splitop
= {false}¶ Split operator model?
-
preconfunc
userprecon
= {nullptr}¶ Pointer to user-supplied preconditioner function.
-
jacobianfunc
userjacobian
= {nullptr}¶ Pointer to user-supplied Jacobian-vector multiply function.
-
bool
initialised
= {false}¶ True if model already initialised.
-
class
PhysicsModelMonitor
: public Monitor¶ - #include <physicsmodel.hxx>
Monitor class for PhysicsModel
Public Functions
-
PhysicsModelMonitor
()¶
-
PhysicsModelMonitor
(PhysicsModel *model)¶
-
int
call
(Solver *solver, BoutReal time, int iter, int nout)¶ Callback function for the solver, called after timestep_ has passed
- Return
- non-zero if simulation should be stopped
- Parameters
solver
: The solver calling this monitortime
: The current simulation timeiter
: The current simulation iterationnout
: The total number of iterations for this simulation
Private Members
-
PhysicsModel *
model
¶
-
-