dolfin/adaptivity

Documentation for C++ code found in dolfin/adaptivity/*.h

Functions

adapt

C++ documentation for adapt from dolfin/adaptivity/adapt.h:

std::shared_ptr<DirichletBC> dolfin::adapt(const DirichletBC &bc, std::shared_ptr<const Mesh> adapted_mesh, const FunctionSpace &S)

Refine Dirichlet bc based on refined mesh.

Parameters:
  • bc
  • adapted_mesh
  • S

C++ documentation for adapt from dolfin/adaptivity/adapt.h:

std::shared_ptr<ErrorControl> dolfin::adapt(const ErrorControl &ec, std::shared_ptr<const Mesh> adapted_mesh, bool adapt_coefficients = true)

Adapt error control object based on adapted mesh

Parameters:
  • ec – (ErrorControl ) The error control object to be adapted
  • adapted_mesh – (Mesh ) The new mesh
  • adapt_coefficients – (bool) Optional argument, default is true. If false, any form coefficients are not explicitly adapted, but pre-adapted coefficients will be transferred.
Returns:

ErrorControl The adapted error control object

C++ documentation for adapt from dolfin/adaptivity/adapt.h:

std::shared_ptr<Form> dolfin::adapt(const Form &form, std::shared_ptr<const Mesh> adapted_mesh, bool adapt_coefficients = true)

Adapt form based on adapted mesh

Parameters:
  • form – (Form ) The form that should be adapted [direction=in]
  • adapted_mesh – (Mesh ) The new mesh [direction=in]
  • adapt_coefficients – (bool) Optional argument, default is true. If false, the form coefficients are not explicitly adapted, but pre-adapted coefficients will be transferred. [direction=in]
Returns:

Form The adapted form

C++ documentation for adapt from dolfin/adaptivity/adapt.h:

std::shared_ptr<Function> dolfin::adapt(const Function &function, std::shared_ptr<const Mesh> adapted_mesh, bool interpolate = true)

Adapt Function based on adapted mesh

Parameters:
  • function – (Function &) The function that should be adapted [direction=in]
  • adapted_mesh – (std::shared_ptr<const Mesh>) The new mesh [direction=in]
  • interpolate – (bool) Optional argument, default is true. If false, the function’s function space is adapted, but the values are not interpolated. [direction=in]
Returns:

Function The adapted function

C++ documentation for adapt from dolfin/adaptivity/adapt.h:

std::shared_ptr<FunctionSpace> dolfin::adapt(const FunctionSpace &space)

Refine function space uniformly

Parameters:space – (FunctionSpace ) [direction=in]
Returns:FunctionSpace

C++ documentation for adapt from dolfin/adaptivity/adapt.h:

std::shared_ptr<FunctionSpace> dolfin::adapt(const FunctionSpace &space, const MeshFunction<bool> &cell_markers)

Refine function space based on cell markers

Parameters:
  • space – (FunctionSpace &) [direction=in]
  • cell_markers – (MehsFunction<bool>&) [direction=in]
Returns:

FunctionSpace

C++ documentation for adapt from dolfin/adaptivity/adapt.h:

std::shared_ptr<FunctionSpace> dolfin::adapt(const FunctionSpace &space, std::shared_ptr<const Mesh> adapted_mesh)

Refine function space based on refined mesh

Parameters:
  • space – (FunctionSpace &) [direction=in]
  • adapted_mesh – (std::sahred_ptr<const Mesh>) [direction=in]
Returns:

FunctionSpace

C++ documentation for adapt from dolfin/adaptivity/adapt.h:

std::shared_ptr<LinearVariationalProblem> dolfin::adapt(const LinearVariationalProblem &problem, std::shared_ptr<const Mesh> adapted_mesh)

Refine linear variational problem based on mesh.

Parameters:
  • problem
  • adapted_mesh

C++ documentation for adapt from dolfin/adaptivity/adapt.h:

std::shared_ptr<Mesh> dolfin::adapt(const Mesh &mesh)

Refine mesh uniformly

Parameters:mesh – (Mesh ) Input mesh [direction=in]
Returns:std::shared_ptr<Mesh> adapted mesh

C++ documentation for adapt from dolfin/adaptivity/adapt.h:

std::shared_ptr<Mesh> dolfin::adapt(const Mesh &mesh, const MeshFunction<bool> &cell_markers)

Refine mesh based on cell markers

Parameters:
  • mesh – (Mesh ) Input mesh [direction=in]
  • cell_markers – (MeshFunction<bool>) Markers denoting cells to be refined [direction=in]

C++ documentation for adapt from dolfin/adaptivity/adapt.h:

std::shared_ptr<MeshFunction<std::size_t>> dolfin::adapt(const MeshFunction<std::size_t> &mesh_function, std::shared_ptr<const Mesh> adapted_mesh)

Refine mesh function<std::size_t> based on mesh.

Parameters:
  • mesh_function
  • adapted_mesh

C++ documentation for adapt from dolfin/adaptivity/adapt.h:

std::shared_ptr<NonlinearVariationalProblem> dolfin::adapt(const NonlinearVariationalProblem &problem, std::shared_ptr<const Mesh> adapted_mesh)

Refine nonlinear variational problem based on mesh.

Parameters:
  • problem
  • adapted_mesh

C++ documentation for adapt from dolfin/adaptivity/adapt.h:

std::shared_ptr<GenericFunction> dolfin::adapt(std::shared_ptr<const GenericFunction> function, std::shared_ptr<const Mesh> adapted_mesh)

Refine GenericFunction based on refined mesh

Parameters:
  • function – (GeericFunction) The function that should be adapted [direction=in]
  • adapted_mesh – (Mehs) The new mesh [direction=in]
Returns:

GenericFunction The adapted function

adapt_markers

C++ documentation for adapt_markers from dolfin/adaptivity/adapt.h:

void dolfin::adapt_markers(std::vector<std::size_t> &refined_markers, const Mesh &adapted_mesh, const std::vector<std::size_t> &markers, const Mesh &mesh)

Helper function for refinement of boundary conditions.

Parameters:
  • refined_markers
  • adapted_mesh
  • markers
  • mesh

dorfler_mark

C++ documentation for dorfler_mark from dolfin/adaptivity/marking.h:

void dolfin::dorfler_mark(MeshFunction<bool> &markers, const dolfin::MeshFunction<double> &indicators, const double fraction)

Mark cells using Dorfler marking

Parameters:
  • markers – (MeshFunction<bool>) the cell markers (to be computed)
  • indicators – (MeshFunction<double>) error indicators (one per cell)
  • fraction – (double) the marking fraction

mark

C++ documentation for mark from dolfin/adaptivity/marking.h:

void dolfin::mark(MeshFunction<bool> &markers, const dolfin::MeshFunction<double> &indicators, const std::string strategy, const double fraction)

Mark cells based on indicators and given marking strategy

Parameters:
  • markers – (MeshFunction<bool>) the cell markers (to be computed)
  • indicators – (MeshFunction<double>) error indicators (one per cell)
  • strategy – (std::string) the marking strategy
  • fraction – (double) the marking fraction

solve

C++ documentation for solve from dolfin/adaptivity/adaptivesolve.h:

void dolfin::solve(const Equation &equation, Function &u, const DirichletBC &bc, const Form &J, const double tol, GoalFunctional &M)

Solve linear variational problem F(u; v) = 0 with single boundary condition

Parameters:
  • equation
  • u
  • bc
  • J
  • tol
  • M

C++ documentation for solve from dolfin/adaptivity/adaptivesolve.h:

void dolfin::solve(const Equation &equation, Function &u, const DirichletBC &bc, const double tol, GoalFunctional &M)

Solve linear variational problem a(u, v) == L(v) with single boundary condition

Parameters:
  • equation
  • u
  • bc
  • tol
  • M

C++ documentation for solve from dolfin/adaptivity/adaptivesolve.h:

void dolfin::solve(const Equation &equation, Function &u, const Form &J, const double tol, GoalFunctional &M)

Solve nonlinear variational problem F(u; v) = 0 without essential boundary conditions

Parameters:
  • equation
  • u
  • J
  • tol
  • M

C++ documentation for solve from dolfin/adaptivity/adaptivesolve.h:

void dolfin::solve(const Equation &equation, Function &u, const double tol, GoalFunctional &M)

Solve linear variational problem a(u, v) == L(v) without essential boundary conditions

Parameters:
  • equation
  • u
  • tol
  • M

C++ documentation for solve from dolfin/adaptivity/adaptivesolve.h:

void dolfin::solve(const Equation &equation, Function &u, std::vector<const DirichletBC *> bcs, const Form &J, const double tol, GoalFunctional &M)

Solve linear variational problem F(u; v) = 0 with list of boundary conditions

Parameters:
  • equation
  • u
  • bcs
  • J
  • tol
  • M

C++ documentation for solve from dolfin/adaptivity/adaptivesolve.h:

void dolfin::solve(const Equation &equation, Function &u, std::vector<const DirichletBC *> bcs, const double tol, GoalFunctional &M)

Solve linear variational problem a(u, v) == L(v) with list of boundary conditions

Parameters:
  • equation
  • u
  • bcs
  • tol
  • M

Classes

AdaptiveLinearVariationalSolver

C++ documentation for AdaptiveLinearVariationalSolver from dolfin/adaptivity/AdaptiveLinearVariationalSolver.h:

class dolfin::AdaptiveLinearVariationalSolver : public dolfin::GenericAdaptiveVariationalSolver

A class for goal-oriented adaptive solution of linear variational problems. For a linear variational problem of the form: find u in V satisfying

a(u, v) = L(v) for all v in :math:`\hat V`

and a corresponding conforming discrete problem: find u_h in V_h satisfying

a(u_h, v) = L(v) for all v in :math:`\hat V_h`

and a given goal functional M and tolerance tol, the aim is to find a V_H and a u_H in V_H satisfying the discrete problem such that

\|M(u) - M(u_H)\| < tol

This strategy is based on dual-weighted residual error estimators designed and automatically generated for the primal problem and subsequent h-adaptivity.

dolfin::AdaptiveLinearVariationalSolver::AdaptiveLinearVariationalSolver(std::shared_ptr<LinearVariationalProblem> problem, std::shared_ptr<Form> goal, std::shared_ptr<ErrorControl> control)

Create AdaptiveLinearVariationalSolver() from variational problem, goal form and error control instance Arguments problem (LinearVariationalProblem ) The primal problem goal (Form ) The goal functional control (ErrorControl ) An error controller object

Parameters:
  • problem
  • goal
  • control
dolfin::AdaptiveLinearVariationalSolver::AdaptiveLinearVariationalSolver(std::shared_ptr<LinearVariationalProblem> problem, std::shared_ptr<GoalFunctional> goal)

Create AdaptiveLinearVariationalSolver() (shared ptr version) Arguments problem (LinearVariationalProblem ) The primal problem goal (GoalFunctional ) The goal functional

Parameters:
  • problem
  • goal
void dolfin::AdaptiveLinearVariationalSolver::adapt_problem(std::shared_ptr<const Mesh> mesh)

Adapt the problem to other mesh. Arguments mesh (Mesh ) The other mesh

Parameters:mesh
double dolfin::AdaptiveLinearVariationalSolver::evaluate_goal(Form &M, std::shared_ptr<const Function> u) const

Evaluate the goal functional. Arguments M (Form ) The functional to be evaluated u (Function ) The function at which to evaluate the functional Returns double The value of M evaluated at u

Parameters:
  • M
  • u
std::vector<std::shared_ptr<const DirichletBC>> dolfin::AdaptiveLinearVariationalSolver::extract_bcs() const

Extract the boundary conditions for the primal problem. Returns std::vector<DirichletBC > The primal boundary conditions

void dolfin::AdaptiveLinearVariationalSolver::init(std::shared_ptr<LinearVariationalProblem> problem, std::shared_ptr<GoalFunctional> goal)

Helper function for instance initialization Arguments problem (LinearVariationalProblem ) The primal problem u (GoalFunctional ) The goal functional

Parameters:
  • problem
  • goal
std::size_t dolfin::AdaptiveLinearVariationalSolver::num_dofs_primal()

Return the number of degrees of freedom for primal problem Returnsstd::size_t The number of degrees of freedom

std::shared_ptr<const Function> dolfin::AdaptiveLinearVariationalSolver::solve_primal()

Solve the primal problem. Returns:cpp:any:Function The solution to the primal problem

dolfin::AdaptiveLinearVariationalSolver::~AdaptiveLinearVariationalSolver()

Destructor.

AdaptiveNonlinearVariationalSolver

C++ documentation for AdaptiveNonlinearVariationalSolver from dolfin/adaptivity/AdaptiveNonlinearVariationalSolver.h:

class dolfin::AdaptiveNonlinearVariationalSolver : public dolfin::GenericAdaptiveVariationalSolver

A class for goal-oriented adaptive solution of nonlinear variational problems. For a nonlinear variational problem of the form: find u in V satisfying

F(u; v) = 0 for all v in :math:`\hat V`

and a corresponding conforming discrete problem: find u_h in V_h satisfying (at least approximately)

F(u_h; v) = 0 for all v in :math:`\hat V_h`

and a given goal functional M and tolerance tol, the aim is to find a V_H and a u_H in V_H satisfying the discrete problem such that

\|M(u) - M(u_H)\| < tol

This strategy is based on dual-weighted residual error estimators designed and automatically generated for the primal problem and subsequent h-adaptivity.

dolfin::AdaptiveNonlinearVariationalSolver::AdaptiveNonlinearVariationalSolver(std::shared_ptr<NonlinearVariationalProblem> problem, std::shared_ptr<Form> goal, std::shared_ptr<ErrorControl> control)

Create AdaptiveLinearVariationalSolver from variational problem, goal form and error control instance Arguments problem (NonlinearVariationalProblem ) The primal problem goal (Form ) The goal functional control (ErrorControl ) An error controller object

Parameters:
  • problem
  • goal
  • control
dolfin::AdaptiveNonlinearVariationalSolver::AdaptiveNonlinearVariationalSolver(std::shared_ptr<NonlinearVariationalProblem> problem, std::shared_ptr<GoalFunctional> goal)

Create AdaptiveNonlinearVariationalSolver() (shared ptr version) Arguments problem (NonlinearVariationalProblem ) The primal problem goal (GoalFunctional ) The goal functional

Parameters:
  • problem
  • goal
void dolfin::AdaptiveNonlinearVariationalSolver::adapt_problem(std::shared_ptr<const Mesh> mesh)

Adapt the problem to other mesh. Arguments mesh (Mesh ) The other mesh

Parameters:mesh
double dolfin::AdaptiveNonlinearVariationalSolver::evaluate_goal(Form &M, std::shared_ptr<const Function> u) const

Evaluate the goal functional. Arguments M (Form ) The functional to be evaluated u (Function ) The function at which to evaluate the functional Returns double The value of M evaluated at u

Parameters:
  • M
  • u
std::vector<std::shared_ptr<const DirichletBC>> dolfin::AdaptiveNonlinearVariationalSolver::extract_bcs() const

Extract the boundary conditions for the primal problem. Returns std::vector<DirichletBC > The primal boundary conditions

void dolfin::AdaptiveNonlinearVariationalSolver::init(std::shared_ptr<NonlinearVariationalProblem> problem, std::shared_ptr<GoalFunctional> goal)

Helper function for instance initialization Arguments problem (NonlinearVariationalProblem ) The primal problem u (GoalFunctional ) The goal functional

Parameters:
  • problem
  • goal
std::size_t dolfin::AdaptiveNonlinearVariationalSolver::num_dofs_primal()

Return the number of degrees of freedom for primal problem Returnsstd::size_t The number of degrees of freedom

std::shared_ptr<const Function> dolfin::AdaptiveNonlinearVariationalSolver::solve_primal()

Solve the primal problem. Returns:cpp:any:Function The solution to the primal problem

dolfin::AdaptiveNonlinearVariationalSolver::~AdaptiveNonlinearVariationalSolver()

Destructor.

ErrorControl

C++ documentation for ErrorControl from dolfin/adaptivity/ErrorControl.h:

class dolfin::ErrorControl

(Goal-oriented) Error Control class. The notation used here follows the notation in “Automated goal-oriented error control I: stationary variational problems”, ME Rognes and A Logg, 2010-2011.

Friends: adapt().

dolfin::ErrorControl::ErrorControl(std::shared_ptr<Form> a_star, std::shared_ptr<Form> L_star, std::shared_ptr<Form> residual, std::shared_ptr<Form> a_R_T, std::shared_ptr<Form> L_R_T, std::shared_ptr<Form> a_R_dT, std::shared_ptr<Form> L_R_dT, std::shared_ptr<Form> eta_T, bool is_linear)

Create error control object

Parameters:
  • a_star – (Form ) the bilinear form for the dual problem
  • L_star – (Form ) the linear form for the dual problem
  • residual – (Form ) a functional for the residual (error estimate)
  • a_R_T – (Form ) the bilinear form for the strong cell residual problem
  • L_R_T – (Form ) the linear form for the strong cell residual problem
  • a_R_dT – (Form ) the bilinear form for the strong facet residual problem
  • L_R_dT – (Form ) the linear form for the strong facet residual problem
  • eta_T – (Form ) a linear form over DG_0 for error indicators
  • is_linear – (bool) true iff primal problem is linear
void dolfin::ErrorControl::apply_bcs_to_extrapolation(const std::vector<std::shared_ptr<const DirichletBC>> bcs)
Parameters:bcs
void dolfin::ErrorControl::compute_cell_residual(Function &R_T, const Function &u)

Compute representation for the strong cell residual from the weak residual

Parameters:
  • R_T – (Function ) the strong cell residual (to be computed)
  • u – (Function ) the primal approximation
void dolfin::ErrorControl::compute_dual(Function &z, const std::vector<std::shared_ptr<const DirichletBC>> bcs)

Compute dual approximation defined by dual variational problem and dual boundary conditions given by homogenized primal boundary conditions.

Parameters:
  • z – (Function ) the dual approximation (to be computed)
  • bcs – (std::vector<DirichletBC>) the primal boundary conditions
void dolfin::ErrorControl::compute_extrapolation(const Function &z, const std::vector<std::shared_ptr<const DirichletBC>> bcs)

Compute extrapolation with boundary conditions

Parameters:
  • z – (Function ) the extrapolated function (to be computed)
  • bcs – (std::vector<DirichletBC >) the dual boundary conditions
void dolfin::ErrorControl::compute_facet_residual(SpecialFacetFunction &R_dT, const Function &u, const Function &R_T)

Compute representation for the strong facet residual from the weak residual and the strong cell residual

Parameters:
void dolfin::ErrorControl::compute_indicators(MeshFunction<double> &indicators, const Function &u)

Compute error indicators

Parameters:
  • indicators – (MeshFunction<double>) the error indicators (to be computed)
  • u – (Function ) the primal approximation
double dolfin::ErrorControl::estimate_error(const Function &u, const std::vector<std::shared_ptr<const DirichletBC>> bcs)

Estimate the error relative to the goal M of the discrete approximation ‘u’ relative to the variational formulation by evaluating the weak residual at an approximation to the dual solution.

Parameters:
  • u – (Function ) the primal approximation
  • bcs – (std::vector<DirichletBC >) the primal boundary conditions
Returns:

double error estimate

void dolfin::ErrorControl::residual_representation(Function &R_T, SpecialFacetFunction &R_dT, const Function &u)

Compute strong representation (strong cell and facet residuals) of the weak residual.

Parameters:
  • R_T – (Function ) the strong cell residual (to be computed)
  • R_dT – (SpecialFacetFunction ) the strong facet residual (to be computed)
  • u – (Function ) the primal approximation
dolfin::ErrorControl::~ErrorControl()

Destructor.

Extrapolation

C++ documentation for Extrapolation from dolfin/adaptivity/Extrapolation.h:

class dolfin::Extrapolation

This class implements an algorithm for extrapolating a function on a given function space from an approximation of that function on a possibly lower-order function space. This can be used to obtain a higher-order approximation of a computed dual solution, which is necessary when the computed dual approximation is in the test space of the primal problem, thereby being orthogonal to the residual. It is assumed that the extrapolation is computed on the same mesh as the original function.

void dolfin::Extrapolation::add_cell_equations(Eigen::MatrixXd &A, Eigen::VectorXd &b, const Cell &cell0, const Cell &cell1, const std::vector<double> &coordinate_dofs0, const std::vector<double> &coordinate_dofs1, const ufc::cell &c0, const ufc::cell &c1, const FunctionSpace &V, const FunctionSpace &W, const Function &v, std::map<std::size_t, std::size_t> &dof2row)
Parameters:
  • A
  • b
  • cell0
  • cell1
  • coordinate_dofs0
  • coordinate_dofs1
  • c0
  • c1
  • V
  • W
  • v
  • dof2row
void dolfin::Extrapolation::average_coefficients(Function &w, std::vector<std::vector<double>> &coefficients)
Parameters:
  • w
  • coefficients
void dolfin::Extrapolation::build_unique_dofs(std::set<std::size_t> &unique_dofs, std::map<std::size_t, std::map<std::size_t, std::size_t>> &cell2dof2row, const Cell &cell0, const FunctionSpace &V)
Parameters:
  • unique_dofs
  • cell2dof2row
  • cell0
  • V
void dolfin::Extrapolation::compute_coefficients(std::vector<std::vector<double>> &coefficients, const Function &v, const FunctionSpace &V, const FunctionSpace &W, const Cell &cell0, const std::vector<double> &coordinate_dofs0, const ufc::cell &c0, const Eigen::Ref<const Eigen::Matrix<dolfin::la_index, Eigen::Dynamic, 1>> dofs, std::size_t &offset)
Parameters:
  • coefficients
  • v
  • V
  • W
  • cell0
  • coordinate_dofs0
  • c0
  • dofs
  • offset
std::map<std::size_t, std::size_t> dolfin::Extrapolation::compute_unique_dofs(const Cell &cell, const FunctionSpace &V, std::size_t &row, std::set<std::size_t> &unique_dofs)
Parameters:
  • cell
  • V
  • row
  • unique_dofs
void dolfin::Extrapolation::extrapolate(Function &w, const Function &v)

Compute extrapolation w from v.

Parameters:
  • w
  • v

GenericAdaptiveVariationalSolver

C++ documentation for GenericAdaptiveVariationalSolver from dolfin/adaptivity/GenericAdaptiveVariationalSolver.h:

class dolfin::GenericAdaptiveVariationalSolver

An abstract class for goal-oriented adaptive solution of variational problems.

void dolfin::GenericAdaptiveVariationalSolver::adapt_problem(std::shared_ptr<const Mesh> mesh) = 0

Adapt the problem to other mesh. Must be overloaded in subclass. Arguments mesh (Mesh ) The other mesh

Parameters:mesh
std::vector<std::shared_ptr<Parameters>> dolfin::GenericAdaptiveVariationalSolver::adaptive_data() const

Return stored adaptive data Returns std::vector<Parameters > The data stored in the adaptive loop

std::shared_ptr<ErrorControl> dolfin::GenericAdaptiveVariationalSolver::control

Error control object.

double dolfin::GenericAdaptiveVariationalSolver::evaluate_goal(Form &M, std::shared_ptr<const Function> u) const = 0

Evaluate the goal functional. Must be overloaded in subclass. Arguments M (Form ) The functional to be evaluated u (Function ) The function of which to evaluate the functional Returns double The value of M evaluated at u

Parameters:
  • M
  • u
std::vector<std::shared_ptr<const DirichletBC>> dolfin::GenericAdaptiveVariationalSolver::extract_bcs() const = 0

Extract the boundary conditions for the primal problem. Must be overloaded in subclass. Returns std::vector<DirichletBC > The primal boundary conditions

std::shared_ptr<Form> dolfin::GenericAdaptiveVariationalSolver::goal

The goal functional.

std::size_t dolfin::GenericAdaptiveVariationalSolver::num_dofs_primal() = 0

Return the number of degrees of freedom for primal problem Returnsstd::size_t The number of degrees of freedom

void dolfin::GenericAdaptiveVariationalSolver::solve(const double tol)

Solve such that the functional error is less than the given tolerance. Note that each call to solve is based on the leaf-node of the variational problem Arguments tol (double) The error tolerance

Parameters:tol
std::shared_ptr<const Function> dolfin::GenericAdaptiveVariationalSolver::solve_primal() = 0

Solve the primal problem. Must be overloaded in subclass. Returns:cpp:any:Function The solution to the primal problem

void dolfin::GenericAdaptiveVariationalSolver::summary()

Present summary of all adaptive data and parameters.

dolfin::GenericAdaptiveVariationalSolver::~GenericAdaptiveVariationalSolver()

GoalFunctional

C++ documentation for GoalFunctional from dolfin/adaptivity/GoalFunctional.h:

class dolfin::GoalFunctional : public dolfin::Form

A GoalFunctional is a Form of rank 0 with an associated ErrorControl .

dolfin::GoalFunctional::GoalFunctional(std::size_t rank, std::size_t num_coefficients)

Create GoalFunctional() Arguments rank (int) the rank of the functional (should be 0) num_coefficients (int) the number of coefficients in functional

Parameters:
  • rank
  • num_coefficients
void dolfin::GoalFunctional::update_ec(const Form &a, const Form &L) = 0

Update error control instance with given forms Arguments a (Form ) a bilinear form L (Form ) a linear form

Parameters:
  • a
  • L

MeshFunction

C++ documentation for MeshFunction from dolfin/adaptivity/adapt.h:

class dolfin::MeshFunction : public dolfin::MeshFunction<T>

A MeshFunction is a function that can be evaluated at a set of mesh entities. A MeshFunction is discrete and is only defined at the set of mesh entities of a fixed topological dimension. A MeshFunction may for example be used to store a global numbering scheme for the entities of a (parallel) mesh, marking sub domains or boolean markers for mesh refinement.

dolfin::MeshFunction::MeshFunction()

Create empty mesh function.

dolfin::MeshFunction::MeshFunction(const MeshFunction<T> &f)

Copy constructor

Parameters:f – (MeshFunction() ) The object to be copied.
dolfin::MeshFunction::MeshFunction(std::shared_ptr<const Mesh> mesh)

Create empty mesh function on given mesh

Parameters:mesh – (Mesh ) The mesh to create mesh function on.
dolfin::MeshFunction::MeshFunction(std::shared_ptr<const Mesh> mesh, const MeshValueCollection<T> &value_collection)

Create function from a MeshValueCollecion (shared_ptr version)

Parameters:
  • mesh – (Mesh ) The mesh to create mesh function on.
  • value_collection – (MeshValueCollection ) The mesh value collection for the mesh function data.
dolfin::MeshFunction::MeshFunction(std::shared_ptr<const Mesh> mesh, const std::string filename)

Create function from data file (shared_ptr version)

Parameters:
  • mesh – (Mesh ) The mesh to create mesh function on.
  • filename – (std::string) The filename to create mesh function from.
dolfin::MeshFunction::MeshFunction(std::shared_ptr<const Mesh> mesh, std::size_t dim)

Create mesh function of given dimension on given mesh

Parameters:
  • mesh – (Mesh ) The mesh to create mesh function on.
  • dim – (std::size_t) The mesh entity dimension for the mesh function.
dolfin::MeshFunction::MeshFunction(std::shared_ptr<const Mesh> mesh, std::size_t dim, const MeshDomains &domains)

Create function from MeshDomains

Parameters:
  • mesh – (Mesh ) The mesh to create mesh function on.
  • dim – (std::size_t) The dimension of the MeshFunction()
  • domains – (_MeshDomains) The domains from which to extract the domain markers
dolfin::MeshFunction::MeshFunction(std::shared_ptr<const Mesh> mesh, std::size_t dim, const T &value)

Create mesh of given dimension on given mesh and initialize to a value

Parameters:
  • mesh – (Mesh ) The mesh to create mesh function on.
  • dim – (std::size_t) The mesh entity dimension.
  • value
    1. The value.
std::size_t dolfin::MeshFunction::dim() const

Return topological dimension

Returns:std::size_t The dimension.
bool dolfin::MeshFunction::empty() const

Return true if empty

Returns:bool True if empty.
void dolfin::MeshFunction::init(std::shared_ptr<const Mesh> mesh, std::size_t dim)

Initialize mesh function for given topological dimension

Parameters:
  • mesh – (Mesh ) The mesh.
  • dim – (std::size_t) The dimension.
void dolfin::MeshFunction::init(std::shared_ptr<const Mesh> mesh, std::size_t dim, std::size_t size)

Initialize mesh function for given topological dimension of given size (shared_ptr version)

Parameters:
  • mesh – (Mesh ) The mesh.
  • dim – (std::size_t) The dimension.
  • size – (std::size_t) The size.
void dolfin::MeshFunction::init(std::size_t dim)

Initialize mesh function for given topological dimension

Parameters:dim – (std::size_t) The dimension.
void dolfin::MeshFunction::init(std::size_t dim, std::size_t size)

Initialize mesh function for given topological dimension of given size

Parameters:
  • dim – (std::size_t) The dimension.
  • size – (std::size_t) The size.
std::shared_ptr<const Mesh> dolfin::MeshFunction::mesh() const

Return mesh associated with mesh function

Returns:Mesh The mesh.
MeshFunction<T> &dolfin::MeshFunction::operator=(const MeshFunction<T> &f)

Assign mesh function to other mesh function Assignment operator

Parameters:f – (MeshFunction ) A MeshFunction object to assign to another MeshFunction .
MeshFunction<T> &dolfin::MeshFunction::operator=(const MeshValueCollection<T> &mesh)

Assignment operator

Parameters:mesh – (MeshValueCollection ) A MeshValueCollection object used to construct a MeshFunction .
const MeshFunction<T> &dolfin::MeshFunction::operator=(const T &value)

Set all values to given value

Parameters:value
T &dolfin::MeshFunction::operator[](const MeshEntity &entity)

Return value at given mesh entity

return T The value at the given entity.

Parameters:entity – (MeshEntity ) The mesh entity.
const T &dolfin::MeshFunction::operator[](const MeshEntity &entity) const

Return value at given mesh entity (const version)

Parameters:entity – (MeshEntity ) The mesh entity.
Returns:T The value at the given entity.
T &dolfin::MeshFunction::operator[](std::size_t index)

Return value at given index

Parameters:index – (std::size_t) The index.
Returns:T The value at the given index.
const T &dolfin::MeshFunction::operator[](std::size_t index) const

Return value at given index (const version)

Parameters:index – (std::size_t) The index.
Returns:T The value at the given index.
void dolfin::MeshFunction::set_all(const T &value)

Set all values to given value

Parameters:value
  1. The value to set all values to.
void dolfin::MeshFunction::set_value(std::size_t index, const T &value)

Set value at given index

Parameters:
  • index – (std::size_t) The index.
  • value
    1. The value.
void dolfin::MeshFunction::set_value(std::size_t index, const T &value, const Mesh &mesh)

Compatibility function for use in SubDomains.

Parameters:
  • index
  • value
  • mesh
void dolfin::MeshFunction::set_values(const std::vector<T> &values)

Set values

Parameters:values – (std::vector<T>) The values.
std::size_t dolfin::MeshFunction::size() const

Return size (number of entities)

Returns:std::size_t The size.
std::string dolfin::MeshFunction::str(bool verbose) const

Return informal string representation (pretty-print)

Parameters:verbose
std::string dolfin::MeshFunction::str(bool verbose) const

Return informal string representation (pretty-print)

Parameters:verbose – (bool) Flag to turn on additional output.
Returns:std::string An informal representation.
T *dolfin::MeshFunction::values()

Return array of values return T The values.

const T *dolfin::MeshFunction::values() const

Return array of values (const. version) return T The values.

std::vector<std::size_t> dolfin::MeshFunction::where_equal(T value)

Get indices where meshfunction is equal to given value Arguments value (T) The value. Returns std::vector<T> The indices.

Parameters:value
dolfin::MeshFunction::~MeshFunction()

Destructor.