DOLFIN
DOLFIN C++ interface
Loading...
Searching...
No Matches
dolfin Namespace Reference

Namespaces

namespace  Encoder

Classes

class  MeshFunction
class  AdaptiveLinearVariationalSolver
class  AdaptiveNonlinearVariationalSolver
class  ErrorControl
 (Goal-oriented) Error Control class. More...
class  Extrapolation
class  GenericAdaptiveVariationalSolver
class  GoalFunctional
class  TimeSeries
class  ALE
class  HarmonicSmoothing
class  MeshDisplacement
class  Array
class  ArrayView
class  Hierarchical
class  IndexSet
class  MPIInfo
class  MPI
class  NoDeleter
 NoDeleter is a customised deleter intended for use with smart pointers. More...
class  RangedIndexSet
class  Set
class  SubSystemsManager
class  Timer
class  UniqueIdGenerator
class  Variable
 Common base class for DOLFIN variables. More...
class  Assembler
class  AssemblerBase
 Provide some common functions used in assembler classes. More...
class  BasisFunction
 Represention of a finite element basis function. More...
class  DirichletBC
 Interface for setting (strong) Dirichlet boundary conditions. More...
class  DiscreteOperators
 Discrete gradient operators providing derivatives of functions. More...
class  DofMap
 Degree-of-freedom map. More...
class  DofMapBuilder
 Builds a DofMap on a Mesh. More...
class  Equation
class  FiniteElement
 This is a wrapper for a UFC finite element (ufc::finite_element). More...
class  Form
 Base class for UFC code generated by FFC for DOLFIN with option -l. More...
class  GenericDofMap
 This class provides a generic interface for dof maps. More...
class  LinearTimeDependentProblem
class  LinearVariationalProblem
class  LinearVariationalSolver
 This class implements a solver for linear variational problems. More...
class  LocalAssembler
class  LocalSolver
 Solve problems cell-wise. More...
class  MixedAssembler
class  MixedLinearVariationalProblem
class  MixedLinearVariationalSolver
 This class implements a solver for mixed linear variational problems. More...
class  MixedNonlinearVariationalProblem
class  MixedNonlinearVariationalSolver
 This class implements a solver for mixed nonlinear variational problems. More...
class  MultiMeshAssembler
class  MultiMeshDirichletBC
class  MultiMeshDofMap
class  MultiMeshForm
class  NonlinearVariationalProblem
class  NonlinearVariationalSolver
class  PETScDMCollection
class  PointSource
class  SparsityPatternBuilder
class  SystemAssembler
class  UFC
class  CoefficientAssigner
class  Constant
 This class represents a constant-valued expression. More...
class  Expression
class  Function
class  FunctionAssigner
class  FunctionAXPY
class  FunctionSpace
class  GenericFunction
class  LagrangeInterpolator
class  MultiMeshCoefficientAssigner
class  MultiMeshFunction
class  MultiMeshFunctionSpace
class  MultiMeshSubSpace
class  SpecialFacetFunction
class  MeshCoordinates
 This Function represents the mesh coordinates on a given mesh. More...
class  FacetArea
class  BoxMesh
class  IntervalMesh
class  RectangleMesh
class  SphericalShellMesh
class  UnitCubeMesh
class  UnitDiscMesh
 A unit disc mesh in 2D or 3D geometry. More...
class  UnitIntervalMesh
class  UnitSquareMesh
class  UnitTetrahedronMesh
class  UnitTriangleMesh
class  BoundingBoxTree
class  BoundingBoxTree1D
 Specialization of bounding box implementation to 1D. More...
class  BoundingBoxTree2D
 Specialization of bounding box implementation to 2D. More...
class  BoundingBoxTree3D
 Specialization of bounding box implementation to 3D. More...
class  CollisionPredicates
class  ConvexTriangulation
class  GenericBoundingBoxTree
class  GeometryDebugging
class  GeometryPredicates
class  GeometryTools
 This class provides useful tools (functions) for computational geometry. More...
class  IntersectionConstruction
class  MeshPointIntersection
class  Point
class  PredicateInitialization
class  SimplexQuadrature
 This class defines quadrature rules for simplices. More...
class  BoostGraphColoring
 This class colors a graph using the Boost Graph Library. More...
class  BoostGraphOrdering
 This class computes graph re-orderings. It uses Boost Graph. More...
class  CSRGraph
 Compressed Sparse Row graph. More...
class  GraphBuilder
 This class builds a Graph corresponding to various objects. More...
class  GraphColoring
 This class provides a common interface to graph coloring libraries. More...
class  ParMETIS
 This class provides an interface to ParMETIS. More...
class  SCOTCH
 This class provides an interface to SCOTCH-PT (parallel version). More...
class  ZoltanInterface
class  File
class  MeshValueCollection
class  GenericFile
 Base class for file I/O objects. More...
class  HDF5Attribute
class  HDF5File
class  HDF5Interface
class  HDF5Utility
class  RAWFile
 Output of data in raw binary format. More...
class  SVGFile
class  VTKFile
 Output of meshes and functions in VTK format. More...
class  VTKWriter
 Write VTK Mesh representation. More...
class  X3DFile
class  X3DOMParameters
 Class data to store X3DOM view parameters. More...
class  X3DOM
class  XDMFFile
 Read and write Mesh, Function, MeshFunction and other objects in XDMF. More...
class  XMLArray
 I/O of array data in XML format. More...
class  XMLFile
 I/O of DOLFIN objects in XML format. More...
class  XMLFunctionData
 I/O for XML representation of Function. More...
class  XMLMesh
 I/O of XML representation of a Mesh. More...
class  XMLMeshFunction
 I/O of XML representation of MeshFunction. More...
class  XMLMeshValueCollection
 I/O of XML representation of a MeshValueCollection. More...
class  XMLParameters
 I/O of Parameters in XML format. More...
class  XMLTable
 Output of XML representation of DOLFIN Table. More...
class  XMLVector
 I/O of XML representation of GenericVector. More...
class  XYZFile
 Simple and light file format for use with Xd3d. More...
class  Amesos2LUSolver
class  BelosKrylovSolver
class  BlockMatrix
 Block Matrix. More...
class  BlockVector
 Block vector. More...
class  CoordinateMatrix
 Coordinate sparse matrix. More...
class  DefaultFactory
 Default linear algebra factory based on global parameter "linear_algebra_backend". More...
class  EigenFactory
 Eigen linear algebra factory. More...
class  EigenKrylovSolver
class  EigenLUSolver
class  EigenMatrix
class  EigenVector
class  GenericLinearAlgebraFactory
 Base class for LinearAlgebra factories. More...
class  GenericLinearOperator
class  GenericLinearSolver
 This class provides a general solver for linear systems Ax = b. More...
class  GenericMatrix
 This class defines a common interface for matrices. More...
class  GenericTensor
 A common interface for arbitrary rank tensors. More...
class  GenericVector
 This class defines a common interface for vectors. More...
class  Ifpack2Preconditioner
 Implements preconditioners using Ifpack2 from Trilinos. More...
class  IndexMap
class  KrylovSolver
class  LinearAlgebraObject
class  LinearOperator
class  LinearSolver
 This class provides a general solver for linear systems Ax = b. More...
class  LUSolver
 LU solver for the built-in LA backends. More...
class  Matrix
class  MueluPreconditioner
 Implements Muelu preconditioner from Trilinos. More...
class  PETScBaseMatrix
class  PETScFactory
 PETSc linear algebra factory. More...
class  PETScKrylovSolver
class  PETScLinearOperator
 PETSc version of the GenericLinearOperator. More...
class  PETScLUSolver
class  PETScMatrix
class  PETScNestMatrix
class  PETScObject
class  PETScOptions
class  PETScPreconditioner
class  PETScVector
class  Scalar
class  SLEPcEigenSolver
class  SparsityPattern
class  SUNDIALSNVector
class  TensorLayout
class  TpetraFactory
 Tpetra linear algebra factory. More...
class  TpetraMatrix
class  TpetraVector
class  TrilinosParameters
class  TrilinosPreconditioner
 This class provides a common base for Trilinos preconditioners. More...
class  Vector
class  VectorSpaceBasis
class  Event
class  Logger
 Handling of error messages, logging and informational display. More...
class  LogManager
 Logger initialisation. More...
class  LogStream
class  Progress
class  Table
class  TableEntry
 This class represents an entry in a Table. More...
class  Lagrange
class  Legendre
 Interface for computing Legendre polynomials via Boost. More...
class  BoundaryComputation
 Provide a set of basic algorithms for the computation of boundaries. More...
class  BoundaryMesh
class  Cell
 A Cell is a MeshEntity of topological codimension 0. More...
class  CellType
class  DistributedMeshTools
class  DomainBoundary
class  DynamicMeshEditor
class  Edge
 An Edge is a MeshEntity of topological dimension 1. More...
class  Face
 A Face is a MeshEntity of topological dimension 2. More...
class  Facet
 A Facet is a MeshEntity of topological codimension 1. More...
class  FacetCell
class  HexahedronCell
 This class implements functionality for hexahedral cell meshes. More...
class  IntervalCell
 This class implements functionality for interval cell meshes. More...
class  LocalMeshData
 This class stores mesh data on a local processor corresponding to a portion of a (larger) global mesh. More...
class  LocalMeshValueCollection
class  Mesh
class  MeshColoring
class  MeshConnectivity
class  MeshData
class  MeshDomains
class  MeshEditor
class  MeshEntity
class  MeshEntityIterator
class  MeshEntityIteratorBase
 Base class for MeshEntityIterators. More...
class  MeshGeometry
 MeshGeometry stores the geometry imposed on a mesh. More...
class  MeshHierarchy
 Experimental implementation of a list of Meshes as a hierarchy. More...
class  MeshOrdering
class  MeshPartitioning
class  MeshQuality
 The class provides functions to quantify mesh quality. More...
class  MeshRelation
class  MeshRenumbering
 This class implements renumbering algorithms for meshes. More...
class  MeshSmoothing
 This class implements various mesh smoothing algorithms. More...
class  MeshTopology
class  MeshTransformation
class  MeshView
class  MultiMesh
class  PeriodicBoundaryComputation
 This class computes map from slave entity to master entity. More...
class  PointCell
 This class implements functionality for point cell meshes. More...
class  QuadrilateralCell
 This class implements functionality for quadrilaterial cells. More...
class  SubDomain
class  SubMesh
class  SubsetIterator
class  TetrahedronCell
 This class implements functionality for tetrahedral cell meshes. More...
class  TopologyComputation
class  TriangleCell
 This class implements functionality for triangular meshes. More...
class  Vertex
 A Vertex is a MeshEntity of topological dimension 0. More...
class  MultiStageScheme
 Place-holder for forms and solutions for a multi-stage Butcher tableau based method. More...
class  PointIntegralSolver
 This class is a time integrator for general Runge Kutta forms. More...
class  RKSolver
 This class is a time integrator for general Runge Kutta problems. More...
class  NewtonSolver
class  NonlinearProblem
class  OptimisationProblem
class  PETScSNESSolver
class  PETScTAOSolver
class  TAOLinearBoundSolver
class  GlobalParameters
 This class defines the global DOLFIN parameter database. More...
class  Parameter
 Base class for parameters. More...
class  Parameters
class  BisectionRefinement1D
 This class implements mesh refinement in 1D. More...
class  LocalMeshCoarsening
 This class implements local mesh coarsening for different mesh types. More...
class  ParallelRefinement
 Data structure and methods for refining meshes in parallel. More...
class  PlazaRefinementND
class  RegularCutRefinement
class  CVode
 Wrapper class to SUNDIALS CVODE. More...

Typedefs

typedef PetscInt la_index
 Index type for compatibility with linear algebra backend(s).
typedef Form TensorProductForm
 FIXME: Temporary fix.
typedef dolfin::Set< int > graph_set_type
 Typedefs for simple graph data structures.
typedef std::vector< graph_set_typeGraph
 Vector of unordered Sets.
typedef MeshEntityIteratorBase< CellCellIterator
 A CellIterator is a MeshEntityIterator of topological codimension 0.
typedef MeshEntityIteratorBase< EdgeEdgeIterator
 An EdgeIterator is a MeshEntityIterator of topological dimension 1.
typedef MeshEntityIteratorBase< FaceFaceIterator
 A FaceIterator is a MeshEntityIterator of topological dimension 2.
typedef MeshEntityIteratorBase< FacetFacetIterator
typedef MeshEntityIteratorBase< VertexVertexIterator
 A VertexIterator is a MeshEntityIterator of topological dimension 0.

Enumerations

enum class  TimingClear : bool { keep = false , clear = true }
enum class  TimingType : int { wall = 0 , user = 1 , system = 2 }
enum  LogLevel {
  CRITICAL = 50 , ERROR = 40 , WARNING = 30 , INFO = 20 ,
  PROGRESS = 16 , TRACE = 13 , DBG = 10
}

Functions

std::shared_ptr< Meshadapt (const Mesh &mesh)
std::shared_ptr< Meshadapt (const Mesh &mesh, const MeshFunction< bool > &cell_markers)
std::shared_ptr< FunctionSpaceadapt (const FunctionSpace &space)
std::shared_ptr< FunctionSpaceadapt (const FunctionSpace &space, const MeshFunction< bool > &cell_markers)
std::shared_ptr< FunctionSpaceadapt (const FunctionSpace &space, std::shared_ptr< const Mesh > adapted_mesh)
std::shared_ptr< Functionadapt (const Function &function, std::shared_ptr< const Mesh > adapted_mesh, bool interpolate=true)
std::shared_ptr< GenericFunctionadapt (std::shared_ptr< const GenericFunction > function, std::shared_ptr< const Mesh > adapted_mesh)
std::shared_ptr< MeshFunction< std::size_t > > adapt (const MeshFunction< std::size_t > &mesh_function, std::shared_ptr< const Mesh > adapted_mesh)
 Refine mesh function<std::size_t> based on mesh.
std::shared_ptr< DirichletBCadapt (const DirichletBC &bc, std::shared_ptr< const Mesh > adapted_mesh, const FunctionSpace &S)
 Refine Dirichlet bc based on refined mesh.
void 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.
std::shared_ptr< Formadapt (const Form &form, std::shared_ptr< const Mesh > adapted_mesh, bool adapt_coefficients=true)
std::shared_ptr< LinearVariationalProblemadapt (const LinearVariationalProblem &problem, std::shared_ptr< const Mesh > adapted_mesh)
 Refine linear variational problem based on mesh.
std::shared_ptr< NonlinearVariationalProblemadapt (const NonlinearVariationalProblem &problem, std::shared_ptr< const Mesh > adapted_mesh)
 Refine nonlinear variational problem based on mesh.
std::shared_ptr< ErrorControladapt (const ErrorControl &ec, std::shared_ptr< const Mesh > adapted_mesh, bool adapt_coefficients=true)
void solve (const Equation &equation, Function &u, const double tol, GoalFunctional &M)
void solve (const Equation &equation, Function &u, const DirichletBC &bc, const double tol, GoalFunctional &M)
void solve (const Equation &equation, Function &u, std::vector< const DirichletBC * > bcs, const double tol, GoalFunctional &M)
void solve (const Equation &equation, Function &u, const Form &J, const double tol, GoalFunctional &M)
void solve (const Equation &equation, Function &u, const DirichletBC &bc, const Form &J, const double tol, GoalFunctional &M)
void solve (const Equation &equation, Function &u, std::vector< const DirichletBC * > bcs, const Form &J, const double tol, GoalFunctional &M)
void mark (MeshFunction< bool > &markers, const dolfin::MeshFunction< double > &indicators, const std::string strategy, const double fraction)
void dorfler_mark (MeshFunction< bool > &markers, const dolfin::MeshFunction< double > &indicators, const double fraction)
std::string dolfin_version ()
 Return DOLFIN version string.
std::string ufc_signature ()
 Return UFC signature string.
std::string git_commit_hash ()
std::size_t sizeof_la_index ()
 Return sizeof the dolfin::la_index type.
bool has_debug ()
bool has_openmp ()
 Return true if DOLFIN is compiled with OpenMP.
bool has_mpi ()
 Return true if DOLFIN is compiled with MPI.
bool has_petsc ()
 Return true if DOLFIN is compiled with PETSc.
bool has_slepc ()
 Return true if DOLFIN is compiled with SLEPc.
bool has_scotch ()
 Return true if DOLFIN is compiled with Scotch.
bool has_sundials ()
 Return true if DOLFIN is compiled with SUNDIALS.
bool has_umfpack ()
 Return true if DOLFIN is compiled with Umfpack.
bool has_cholmod ()
 Return true if DOLFIN is compiled with Cholmod.
bool has_parmetis ()
 Return true if DOLFIN is compiled with ParMETIS.
bool has_zlib ()
 Return true if DOLFIN is compiled with ZLIB.
bool has_hdf5 ()
 Return true if DOLFIN is compiled with HDF5.
bool has_hdf5_parallel ()
 Return true if DOLFIN is compiled with Parallel HDF5.
void init (int argc, char *argv[])
template<typename T>
std::shared_ptr< T > reference_to_no_delete_pointer (T &r)
 Helper function to construct shared pointer with NoDeleter with cleaner syntax.
void tic ()
 Start timing (should not be used internally in DOLFIN!).
double toc ()
 Return elapsed wall time (should not be used internally in DOLFIN!).
double time ()
 Return wall time elapsed since some implementation dependent epoch.
Table timings (TimingClear clear, std::set< TimingType > type)
void list_timings (TimingClear clear, std::set< TimingType > type)
void dump_timings_to_xml (std::string filename, TimingClear clear)
std::tuple< std::size_t, double, double, double > timing (std::string task, TimingClear clear)
std::string indent (std::string block)
 Indent string block.
template<typename T>
std::string container_to_string (const T &x, std::string delimiter, int precision, int linebreak=0)
std::string to_string (const double *x, std::size_t n)
 Return string representation of given array.
template<class T>
std::size_t hash_local (const T &x)
 Return a hash of a given object.
template<class T>
std::size_t hash_global (const MPI_Comm mpi_comm, const T &x)
void assemble (GenericTensor &A, const Form &a)
 Assemble tensor.
void assemble_mixed (GenericTensor &A, const Form &a, bool add=false)
 Assemble tensor.
void assemble_system (GenericMatrix &A, GenericVector &b, const Form &a, const Form &L, std::vector< std::shared_ptr< const DirichletBC > > bcs)
 Assemble system (A, b) and apply Dirichlet boundary conditions.
void assemble_system (GenericMatrix &A, GenericVector &b, const Form &a, const Form &L, std::vector< std::shared_ptr< const DirichletBC > > bcs, const GenericVector &x0)
void assemble_multimesh (GenericTensor &A, const MultiMeshForm &a)
 Assemble tensor from multimesh form.
double assemble (const Form &a)
 Assemble scalar.
double assemble_mixed (const Form &a, bool add=false)
 Assemble scalar.
double assemble_multimesh (const MultiMeshForm &a)
 Assemble scalar from multimesh form.
void assemble_local (Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &A_e, const Form &a, const Cell &cell)
 Assemble form to local tensor on a cell (Eigen version for pybind11).
void assemble_local (const Form &a, const Cell &cell, std::vector< double > &tensor)
std::vector< std::size_t > dof_to_vertex_map (const FunctionSpace &space)
std::vector< dolfin::la_indexvertex_to_dof_map (const FunctionSpace &space)
void set_coordinates (MeshGeometry &geometry, const Function &position)
void get_coordinates (Function &position, const MeshGeometry &geometry)
Mesh create_mesh (Function &coordinates)
void solve (const Equation &equation, Function &u, Parameters parameters=empty_parameters)
void solve (const Equation &equation, Function &u, const DirichletBC &bc, Parameters parameters=empty_parameters)
void solve (const Equation &equation, Function &u, std::vector< const DirichletBC * > bcs, Parameters parameters=empty_parameters)
void solve (const Equation &equation, Function &u, const Form &J, Parameters parameters=empty_parameters)
void solve (const Equation &equation, Function &u, const DirichletBC &bc, const Form &J, Parameters parameters=empty_parameters)
void solve (const Equation &equation, Function &u, std::vector< const DirichletBC * > bcs, const Form &J, Parameters parameters=empty_parameters)
void assign (std::shared_ptr< Function > receiving_func, std::shared_ptr< const Function > assigning_func)
void assign (std::shared_ptr< Function > receiving_func, std::vector< std::shared_ptr< const Function > > assigning_funcs)
void assign (std::vector< std::shared_ptr< Function > > receiving_funcs, std::shared_ptr< const Function > assigning_func)
bool check_cgal (bool result_dolfin, bool result_cgal, const std::string &function)
std::vector< Pointcgal_intersection_check (const std::vector< Point > &dolfin_result, const std::vector< Point > &cgal_result, const std::string &function)
bool cgal_collides_segment_point_2d (const Point &q0, const Point &q1, const Point &p, bool only_interior=false)
bool cgal_collides_segment_point_3d (const Point &q0, const Point &q1, const Point &p, bool only_interior=false)
bool cgal_collides_segment_segment_2d (const Point &p0, const Point &p1, const Point &q0, const Point &q1)
bool cgal_collides_segment_segment_3d (const Point &p0, const Point &p1, const Point &q0, const Point &q1)
bool cgal_collides_triangle_point_2d (const Point &p0, const Point &p1, const Point &p2, const Point &point)
bool cgal_collides_triangle_point_3d (const Point &p0, const Point &p1, const Point &p2, const Point &point)
bool cgal_collides_triangle_segment_2d (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1)
bool cgal_collides_triangle_segment_3d (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1)
bool cgal_collides_triangle_triangle_2d (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2)
bool cgal_collides_triangle_triangle_3d (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2)
bool cgal_collides_tetrahedron_point_3d (const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0)
bool cgal_collides_tetrahedron_segment_3d (const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1)
bool cgal_collides_tetrahedron_triangle_3d (const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1, const Point &q2)
bool cgal_collides_tetrahedron_tetrahedron_3d (const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1, const Point &q2, const Point &q3)
std::vector< Pointcgal_intersection_segment_segment_2d (const Point &p0, const Point &p1, const Point &q0, const Point &q1)
std::vector< Pointcgal_intersection_segment_segment_3d (const Point &p0, const Point &p1, const Point &q0, const Point &q1)
std::vector< Pointcgal_triangulate_segment_segment_2d (const Point &p0, const Point &p1, const Point &q0, const Point &q1)
std::vector< Pointcgal_triangulate_segment_segment_3d (const Point &p0, const Point &p1, const Point &q0, const Point &q1)
std::vector< Pointcgal_intersection_triangle_segment_2d (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1)
std::vector< Pointcgal_intersection_triangle_segment_3d (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1)
std::vector< Pointcgal_triangulate_triangle_segment_2d (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1)
std::vector< Pointcgal_triangulate_triangle_segment_3d (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1)
std::vector< Pointcgal_intersection_triangle_triangle_2d (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2)
std::vector< Pointcgal_intersection_triangle_triangle_3d (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2)
std::vector< std::vector< Point > > cgal_triangulate_triangle_triangle_2d (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2)
std::vector< std::vector< Point > > cgal_triangulate_triangle_triangle_3d (const Point &p0, const Point &p1, const Point &p2, const Point &q0, const Point &q1, const Point &q2)
std::vector< Pointcgal_intersection_tetrahedron_triangle (const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1, const Point &q2)
std::vector< std::vector< Point > > cgal_triangulate_tetrahedron_triangle (const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1, const Point &q2)
std::vector< Pointcgal_intersection_tetrahedron_tetrahedron_3d (const Point &p0, const Point &p1, const Point &p2, const Point &p3, const Point &q0, const Point &q1, const Point &q2, const Point &q3)
bool cgal_is_degenerate_2d (const std::vector< Point > &s)
bool cgal_is_degenerate_3d (const std::vector< Point > &s)
double cgal_polyhedron_volume (const std::vector< Point > &ch)
double cgal_tet_volume (const std::vector< Point > &ch)
bool cgal_tet_is_degenerate (const std::vector< Point > &t)
bool cgal_triangulation_has_degenerate (std::vector< std::vector< Point > > triangulation)
bool cgal_triangulation_overlap (std::vector< std::vector< Point > > triangulation)
std::shared_ptr< const MeshPointIntersectionintersect (const Mesh &mesh, const Point &point)
Point operator* (double a, const Point &p)
 Multiplication with scalar.
std::ostream & operator<< (std::ostream &stream, const Point &point)
 Output of Point to stream.
void exactinit ()
 Initialize tolerances for exact arithmetic.
double orient1d (double a, double b, double x)
 Compute relative orientation of point x wrt segment [a, b].
double _orient2d (const double *a, const double *b, const double *c)
double orient2d (const Point &a, const Point &b, const Point &c)
 Convenience function using dolfin::Point.
double _orient3d (const double *a, const double *b, const double *c, const double *d)
double orient3d (const Point &a, const Point &b, const Point &c, const Point &d)
 Convenience function using dolfin::Point.
template<>
hid_t HDF5Interface::hdf5_type< std::int64_t > ()
template<>
hid_t HDF5Interface::hdf5_type< std::size_t > ()
template<typename Y, typename X>
Y & as_type (X &x)
template<typename Y, typename X>
std::shared_ptr< Y > as_type (std::shared_ptr< X > x)
template<typename Y, typename X>
bool has_type (const X &x)
 Check whether the object matches a specific type.
int usermult (Mat A, Vec x, Vec y)
 Callback function for PETSc mult function.
std::size_t solve (const GenericLinearOperator &A, GenericVector &x, const GenericVector &b, std::string method="lu", std::string preconditioner="none")
 Solve linear system Ax = b.
void list_linear_algebra_backends ()
 List available linear algebra backends.
void list_linear_solver_methods ()
 List available solver methods for current linear algebra backend.
void list_lu_solver_methods ()
 List available LU methods for current linear algebra backend.
void list_krylov_solver_methods ()
 List available Krylov methods for current linear algebra backend.
void list_krylov_solver_preconditioners ()
bool has_linear_algebra_backend (std::string backend)
 Return true if a specific linear algebra backend is supported.
bool has_lu_solver_method (std::string method)
bool has_krylov_solver_method (std::string method)
bool has_krylov_solver_preconditioner (std::string preconditioner)
std::map< std::string, std::string > linear_algebra_backends ()
 Return available linear algebra backends.
std::map< std::string, std::string > linear_solver_methods ()
std::map< std::string, std::string > lu_solver_methods ()
std::map< std::string, std::string > krylov_solver_methods ()
std::map< std::string, std::string > krylov_solver_preconditioners ()
double residual (const GenericLinearOperator &A, const GenericVector &x, const GenericVector &b)
 Compute residual ||Ax - b||.
double norm (const GenericVector &x, std::string norm_type="l2")
double normalize (GenericVector &x, std::string normalization_type="average")
 Normalize vector according to given normalization type.
bool in_nullspace (const GenericLinearOperator &A, const VectorSpaceBasis &x, std::string type="right")
void info (std::string msg,...)
 Print message.
void info (const Parameters &parameters, bool verbose=false)
 Print parameter (using output of str() method).
void info (const Variable &variable, bool verbose=false)
 Print variable (using output of str() method).
void info_stream (std::ostream &out, std::string msg)
 Print message to stream.
void info_underline (std::string msg,...)
 Print underlined message.
void warning (std::string msg,...)
 Print warning.
void error (std::string msg,...)
void dolfin_error (std::string location, std::string task, std::string reason,...)
void deprecation (std::string feature, std::string version_deprecated, std::string message,...)
void log (int debug_level, std::string msg,...)
 Print message at given debug level.
void begin (std::string msg,...)
 Begin task (increase indentation level).
void begin (int debug_level, std::string msg,...)
 Begin task (increase indentation level).
void end ()
 End task (decrease indentation level).
void set_log_active (bool active=true)
 Turn logging on or off.
void set_log_level (int level)
 Set log level.
void set_indentation_level (std::size_t indentation_level)
 Set indentation level.
void set_output_stream (std::ostream &out)
 Set output stream.
int get_log_level ()
 Get log level.
void monitor_memory_usage ()
void not_working_in_parallel (std::string what)
void __debug (std::string file, unsigned long line, std::string function, std::string format,...)
void __dolfin_assert (std::string file, unsigned long line, std::string function, std::string check)
std::size_t ipow (std::size_t a, std::size_t n)
double rand ()
void seed (std::size_t s)
bool near (double x, double x0, double eps=DOLFIN_EPS)
bool between (double x, std::pair< double, double > range)
Mesh refine (const Mesh &mesh, bool redistribute=true)
std::shared_ptr< const MeshHierarchyrefine (const MeshHierarchy &hierarchy, const MeshFunction< bool > &markers)
 Refine a MeshHierarchy.
void refine (Mesh &refined_mesh, const Mesh &mesh, bool redistribute=true)
Mesh refine (const Mesh &mesh, const MeshFunction< bool > &cell_markers, bool redistribute=true)
void refine (Mesh &refined_mesh, const Mesh &mesh, const MeshFunction< bool > &cell_markers, bool redistribute=true)
void p_refine (Mesh &refined_mesh, const Mesh &mesh)
Mesh p_refine (const Mesh &mesh)

Variables

Timer __global_timer
Timer __tic_timer
PredicateInitialization predicate_initialization
 Initialize the predicate.
LogStream cout
 dolfin::cout
LogStream endl
 dolfin::endl;
bool rand_seeded = false
GlobalParameters parameters
 The global parameter database.
Parameters empty_parameters ("empty")
 Default empty parameters.

Detailed Description

This comment in in timing.h but I think it is providing a doxygen docstring for the whole dolfin namespace... FIXME.

Typedef Documentation

◆ FacetIterator

◆ graph_set_type

Typedefs for simple graph data structures.

DOLFIN container for graphs

Enumeration Type Documentation

◆ LogLevel

These log levels match the levels in the Python 'logging' module (and adds trace/progress).

◆ TimingClear

enum class dolfin::TimingClear : bool
strong

Parameter specifying whether to clear timing(s):

  • TimingClear::keep
  • TimingClear::clear

◆ TimingType

enum class dolfin::TimingType : int
strong

Timing types:

  • TimingType::wall wall-clock time
  • TimingType::user user (cpu) time
  • TimingType::system system (kernel) time

Precision of wall is around 1 microsecond, user and system are around 10 millisecond (on Linux).

Function Documentation

◆ _orient2d()

double dolfin::_orient2d ( const double * a,
const double * b,
const double * c )

Compute relative orientation of points a, b, c. The orientation is such that orient2d(a, b, c) > 0 if a, b, c are ordered counter-clockwise.

◆ _orient3d()

double dolfin::_orient3d ( const double * a,
const double * b,
const double * c,
const double * d )

Compute relative orientation of points a, b, c, d. The orientation is such that orient3d(a, b, c, d) > 0 if a, b, c, d are oriented according to the left hand rule.

◆ adapt() [1/9]

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

◆ adapt() [2/9]

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
[in]form(Form) The form that should be adapted
[in]adapted_mesh(Mesh) The new mesh
[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.
Returns
Form The adapted form

◆ adapt() [3/9]

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
[in]function(Function&) The function that should be adapted
[in]adapted_mesh(std::shared_ptr<const Mesh>) The new mesh
[in]interpolate(bool) Optional argument, default is true. If false, the function's function space is adapted, but the values are not interpolated.
Returns
Function The adapted function

◆ adapt() [4/9]

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

Refine function space uniformly

Parameters
[in]space(FunctionSpace)
Returns
FunctionSpace

◆ adapt() [5/9]

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

Refine function space based on cell markers

Parameters
[in]space(FunctionSpace&)
[in]cell_markers(MehsFunction<bool>&)
Returns
FunctionSpace

◆ adapt() [6/9]

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

Refine function space based on refined mesh

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

◆ adapt() [7/9]

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

Refine mesh uniformly

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

◆ adapt() [8/9]

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

Refine mesh based on cell markers

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

◆ adapt() [9/9]

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
[in]function(GeericFunction) The function that should be adapted
[in]adapted_mesh(Mehs) The new mesh
Returns
GenericFunction The adapted function

◆ as_type() [1/2]

template<typename Y, typename X>
std::shared_ptr< Y > dolfin::as_type ( std::shared_ptr< X > x)

Cast shared pointer object to its derived class, if possible. Caller must check for success (returns null if cast fails).

◆ as_type() [2/2]

template<typename Y, typename X>
Y & dolfin::as_type ( X & x)

Cast object to its derived class, if possible (non-const version). An error is thrown if the cast is unsuccessful.

◆ assemble_local()

void dolfin::assemble_local ( const Form & a,
const Cell & cell,
std::vector< double > & tensor )

Assemble form to local tensor on a cell (Legacy version for SWIG)

◆ assemble_system()

void dolfin::assemble_system ( GenericMatrix & A,
GenericVector & b,
const Form & a,
const Form & L,
std::vector< std::shared_ptr< const DirichletBC > > bcs,
const GenericVector & x0 )

Assemble system (A, b) on sub domains and apply Dirichlet boundary conditions

◆ assign() [1/3]

void dolfin::assign ( std::shared_ptr< Function > receiving_func,
std::shared_ptr< const Function > assigning_func )

Assign one function to another. The functions must reside in the same type of FunctionSpace. One or both functions can be sub functions.

Parameters
receiving_func(std::shared_ptr<Function>) The receiving function
assigning_func(std::shared_ptr<Function>) The assigning function

◆ assign() [2/3]

void dolfin::assign ( std::shared_ptr< Function > receiving_func,
std::vector< std::shared_ptr< const Function > > assigning_funcs )

Assign several functions to sub functions of a mixed receiving function. The number of receiving functions must sum up to the number of sub functions in the assigning mixed function. The sub spaces of the assigning mixed space must be of the same type ans size as the receiving spaces.

Parameters
receiving_func(std::shared_ptr<Function>) The receiving function
assigning_funcs(std::vector<std::shared_ptr<Function>>) The assigning functions

◆ assign() [3/3]

void dolfin::assign ( std::vector< std::shared_ptr< Function > > receiving_funcs,
std::shared_ptr< const Function > assigning_func )

Assign sub functions of a single mixed function to single receiving functions. The number of sub functions in the assigning mixed function must sum up to the number of receiving functions. The sub spaces of the receiving mixed space must be of the same type ans size as the assigning spaces.

Parameters
receiving_funcs(std::vector<std::shared_ptr<Function>>) The receiving functions
assigning_func(std::shared_ptr<Function>) The assigning function

◆ between()

bool dolfin::between ( double x,
std::pair< double, double > range )

Check whether x is between x0 and x1 (inclusive, to within DOLFIN_EPS)

Parameters
x(double) Value to check
range(std::pair<double, double>) Range to check
Returns
bool

◆ container_to_string()

template<typename T>
std::string dolfin::container_to_string ( const T & x,
std::string delimiter,
int precision,
int linebreak = 0 )

Return string representation of given container of ints, floats, etc.

◆ create_mesh()

Mesh dolfin::create_mesh ( Function & coordinates)

Creates mesh from coordinate function

Topology is given by underlying mesh of the function space and geometry is given by function values. Hence resulting mesh geometry has a degree of the function space degree. Geometry of function mesh is ignored.

Mesh connectivities d-0, d-1, ..., d-r are built on function mesh (where d is topological dimension of the mesh and r is maximal dimension of entity associated with any coordinate node). Consider clearing unneeded connectivities when finished.

Parameters
coordinates(Function) Vector Lagrange function of any degree
Returns
Mesh The mesh

◆ deprecation()

void dolfin::deprecation ( std::string feature,
std::string version_deprecated,
std::string message,
... )

Issue deprecation warning for removed feature

Arguments feature (std::string) Name of the feature that has been removed. version_deprecated (std::string) Version number of the release in which the feature is deprecated. message (std::string) A format string explaining the deprecation.

◆ dof_to_vertex_map()

std::vector< std::size_t > dolfin::dof_to_vertex_map ( const FunctionSpace & space)

Return a map between dof indices and vertex indices

Only works for FunctionSpace with dofs exclusively on vertices. For mixed FunctionSpaces vertex index is offset with the number of dofs per vertex.

In parallel the returned map maps both owned and unowned dofs (using local indices) thus covering all the vertices. Hence the returned map is an inversion of vertex_to_dof_map.

Parameters
space(FunctionSpace) The FunctionSpace for what the dof to vertex map should be computed for
Returns
std::vector<std::size_t> The dof to vertex map

◆ dolfin_error()

void dolfin::dolfin_error ( std::string location,
std::string task,
std::string reason,
... )

Print error message. Prefer this to the above generic error message.

Arguments location (std::string) Name of the file from which the error message was generated. task (std::string) Name of the task that failed. Note that this string should begin with lowercase. Note that this string should not be punctuated. reason (std::string) A format string explaining the reason for the failure. Note that this string should begin with uppercase. Note that this string should not be punctuated. Note that this string may contain printf style formatting. ... (primitive types like int, std::size_t, double, bool) Optional arguments for the format string.

Developers should read the file dolfin/log/README in the DOLFIN source tree for further notes about the use of this function.

◆ dorfler_mark()

void dolfin::dorfler_mark ( dolfin::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

◆ dump_timings_to_xml()

void dolfin::dump_timings_to_xml ( std::string filename,
TimingClear clear )

Dump a summary of timings and tasks to XML file, optionally clearing stored timings. MPI_MAX, MPI_MIN and MPI_AVG reductions are stored. Collective on MPI_COMM_WORLD.

Arguments filename (std::string) output filename; must have .xml suffix; existing file is silently overwritten clear (TimingClear)

  • TimingClear::clear resets stored timings
  • TimingClear::keep leaves stored timings intact

◆ error()

void dolfin::error ( std::string msg,
... )

Print error message and throw an exception. Note to developers: this function should not be used internally in DOLFIN. Use the more informative dolfin_error instead.

◆ get_coordinates()

void dolfin::get_coordinates ( Function & position,
const MeshGeometry & geometry )

Stores mesh coordinates into function

Mesh connectivities d-0, d-1, ..., d-r are built on function mesh (where d is topological dimension of the mesh and r is maximal dimension of entity associated with any coordinate node). Consider clearing unneeded connectivities when finished.

Parameters
position(Function) Vectorial Lagrange function with matching degree and mesh
geometry(MeshGeometry) Mesh geometry to be stored

◆ git_commit_hash()

std::string dolfin::git_commit_hash ( )

Return git changeset hash (returns "unknown" if changeset is not known)

◆ has_debug()

bool dolfin::has_debug ( )

Return true if DOLFIN is compiled in debugging mode, i.e., with assertions on

◆ has_krylov_solver_method()

bool dolfin::has_krylov_solver_method ( std::string method)

Return true if Krylov method for the current linear algebra backend is available

◆ has_krylov_solver_preconditioner()

bool dolfin::has_krylov_solver_preconditioner ( std::string preconditioner)

Return true if Preconditioner for the current linear algebra backend is available

◆ has_lu_solver_method()

bool dolfin::has_lu_solver_method ( std::string method)

Return true if LU method for the current linear algebra backend is available

◆ hash_global()

template<class T>
std::size_t dolfin::hash_global ( const MPI_Comm mpi_comm,
const T & x )

Return a hash for a distributed (MPI) object. A hash is computed on each process, and the hash of the std::vector of all local hash keys is returned. This function is collective.

◆ in_nullspace()

bool dolfin::in_nullspace ( const GenericLinearOperator & A,
const VectorSpaceBasis & x,
std::string type = "right" )

Check whether a vector space basis is in the nullspace of a given operator. The string option 'type' can be "right" for the right nullspace (Ax=0) or "left" for the left nullspace (A^Tx = 0). To test the left nullspace, A must also be of type GenericMatrix.

◆ info()

void dolfin::info ( std::string msg,
... )

Print message.

The DOLFIN log system provides the following set of functions for uniform handling of log messages, warnings and errors. In addition, macros are provided for debug messages and dolfin_assertions.

Only messages with a debug level higher than or equal to the current log level are printed (the default being zero). Logging may also be turned off by calling set_log_active(false).

◆ init()

void dolfin::init ( int argc,
char * argv[] )

Initialize DOLFIN (and PETSc) with command-line arguments. This should not be needed in most cases since the initialization is otherwise handled automatically.

◆ intersect()

std::shared_ptr< const MeshPointIntersection > dolfin::intersect ( const Mesh & mesh,
const Point & point )

Compute and return intersection between Mesh and Point.

Arguments mesh (Mesh) The mesh to be intersected. point (Point) The point to be intersected.

Returns MeshPointIntersection The intersection data.

◆ ipow()

std::size_t dolfin::ipow ( std::size_t a,
std::size_t n )

Return a to the power n. NOTE: Overflow is not checked!

Parameters
a(std::size_t) Value
n(std::size_t) Power
Returns
std::size_t

◆ krylov_solver_methods()

std::map< std::string, std::string > dolfin::krylov_solver_methods ( )

Return a list of available Krylov methods for current linear algebra backend

◆ krylov_solver_preconditioners()

std::map< std::string, std::string > dolfin::krylov_solver_preconditioners ( )

Return a list of available preconditioners for current linear algebra backend

◆ linear_solver_methods()

std::map< std::string, std::string > dolfin::linear_solver_methods ( )

Return a list of available solver methods for current linear algebra backend

◆ list_krylov_solver_preconditioners()

void dolfin::list_krylov_solver_preconditioners ( )

List available preconditioners for current linear algebra backend

◆ list_timings()

void dolfin::list_timings ( TimingClear clear,
std::set< TimingType > type )

List a summary of timings and tasks, optionally clearing stored timings. MPI_AVG reduction is printed. Collective on MPI_COMM_WORLD.

Arguments clear (TimingClear)

  • TimingClear::clear resets stored timings
  • TimingClear::keep leaves stored timings intact type (std::set<TimingType>) subset of { TimingType::wall, TimingType::user, TimingType::system }

◆ lu_solver_methods()

std::map< std::string, std::string > dolfin::lu_solver_methods ( )

Return a list of available LU methods for current linear algebra backend

◆ mark()

void dolfin::mark ( dolfin::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

◆ monitor_memory_usage()

void dolfin::monitor_memory_usage ( )

Monitor memory usage. Call this function at the start of a program to continuously monitor the memory usage of the process.

◆ near()

bool dolfin::near ( double x,
double x0,
double eps = DOLFIN_EPS )

Check whether x is close to x0 (to within DOLFIN_EPS)

Parameters
x(double) First value
x0(double) Second value
eps(double) Tolerance
Returns
bool

◆ norm()

double dolfin::norm ( const GenericVector & x,
std::string norm_type = "l2" )

Compute norm of vector. Valid norm types are "l2", "l1" and "linf".

◆ not_working_in_parallel()

void dolfin::not_working_in_parallel ( std::string what)

Report that functionality has not (yet) been implemented to work in parallel

◆ p_refine() [1/2]

dolfin::Mesh dolfin::p_refine ( const Mesh & mesh)

Return a p_refined mesh Increase the polynomial order of the mesh from 1 to 2, i.e. add points at the Edge midpoints, to make a quadratic mesh.

Parameters
mesh(Mesh) The original linear mesh.
Returns
Mesh

◆ p_refine() [2/2]

void dolfin::p_refine ( Mesh & refined_mesh,
const Mesh & mesh )

Increase the polynomial order of the mesh from 1 to 2, i.e. add points at the Edge midpoints, to make a quadratic mesh.

Parameters
refined_mesh(Mesh) The mesh that will be the quadratic mesh.
mesh(Mesh) The original linear mesh.

◆ rand()

double dolfin::rand ( )

Return a random number, uniformly distributed between [0.0, 1.0)

Returns
double

◆ refine() [1/4]

dolfin::Mesh dolfin::refine ( const Mesh & mesh,
bool redistribute = true )

Create uniformly refined mesh

Parameters
mesh(Mesh) The mesh to refine.
redistribute(bool) Optional argument to redistribute the refined mesh if mesh is a distributed mesh.
Returns
Mesh The refined mesh.
mesh = refine(mesh);
Mesh refine(const Mesh &mesh, bool redistribute=true)
Definition refine.cpp:39

◆ refine() [2/4]

dolfin::Mesh dolfin::refine ( const Mesh & mesh,
const MeshFunction< bool > & cell_markers,
bool redistribute = true )

Create locally refined mesh

Parameters
mesh(Mesh) The mesh to refine.
cell_markers(MeshFunction<bool>) A mesh function over booleans specifying which cells that should be refined (and which should not).
redistribute(bool) Optional argument to redistribute the refined mesh if mesh is a distributed mesh.
Returns
Mesh The locally refined mesh.
MeshFunction<bool> cell_markers(mesh, mesh->topology().dim());
cell_markers.set_all(false);
Point origin(0.0, 0.0, 0.0);
for (CellIterator cell(mesh); !cell.end(); ++cell)
{
Point p = cell->midpoint();
if (p.distance(origin) < 0.1)
cell_markers[*cell] = true;
}
mesh = refine(mesh, cell_markers);
Definition MeshFunction.h:58
Definition Point.h:41
double distance(const Point &p) const
Definition Point.h:229
MeshEntityIteratorBase< Cell > CellIterator
A CellIterator is a MeshEntityIterator of topological codimension 0.
Definition Cell.h:414

◆ refine() [3/4]

void dolfin::refine ( Mesh & refined_mesh,
const Mesh & mesh,
bool redistribute = true )

Create uniformly refined mesh

Parameters
refined_mesh(Mesh) The mesh that will be the refined mesh.
mesh(Mesh) The original mesh.
redistribute(bool) Optional argument to redistribute the refined mesh if mesh is a distributed mesh.

◆ refine() [4/4]

void dolfin::refine ( Mesh & refined_mesh,
const Mesh & mesh,
const MeshFunction< bool > & cell_markers,
bool redistribute = true )

Create locally refined mesh

Parameters
refined_mesh(Mesh) The mesh that will be the refined mesh.
mesh(Mesh) The original mesh.
cell_markers(MeshFunction<bool>) A mesh function over booleans specifying which cells that should be refined (and which should not).
redistribute(bool) Optional argument to redistribute the refined mesh if mesh is a distributed mesh.

◆ seed()

void dolfin::seed ( std::size_t s)

Seed random number generator

Parameters
s(std::size_t) Seed value

◆ set_coordinates()

void dolfin::set_coordinates ( MeshGeometry & geometry,
const Function & position )

Sets mesh coordinates from function

Mesh connectivities d-0, d-1, ..., d-r are built on function mesh (where d is topological dimension of the mesh and r is maximal dimension of entity associated with any coordinate node). Consider clearing unneeded connectivities when finished.

Parameters
geometry(MeshGeometry) Mesh geometry to be set
position(Function) Vectorial Lagrange function with matching degree and mesh

◆ solve() [1/12]

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

◆ solve() [2/12]

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

◆ solve() [3/12]

void dolfin::solve ( const Equation & equation,
Function & u,
const DirichletBC & bc,
const Form & J,
Parameters parameters = empty_parameters )

Solve nonlinear variational problem F(u; v) == 0 with a single boundary condition. The argument J should provide the Jacobian bilinear form J = dF/du.

Optional parameters can be passed to the LinearVariationalSolver or NonlinearVariationalSolver classes.

◆ solve() [4/12]

void dolfin::solve ( const Equation & equation,
Function & u,
const DirichletBC & bc,
Parameters parameters = empty_parameters )

Solve linear variational problem a(u, v) == L(v) or nonlinear variational problem F(u; v) = 0 with a single boundary condition.

Optional parameters can be passed to the LinearVariationalSolver or NonlinearVariationalSolver classes.

◆ solve() [5/12]

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

◆ solve() [6/12]

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

◆ solve() [7/12]

void dolfin::solve ( const Equation & equation,
Function & u,
const Form & J,
Parameters parameters = empty_parameters )

Solve nonlinear variational problem F(u; v) == 0 without boundary conditions. The argument J should provide the Jacobian bilinear form J = dF/du.

Optional parameters can be passed to the LinearVariationalSolver or NonlinearVariationalSolver classes.

◆ solve() [8/12]

void dolfin::solve ( const Equation & equation,
Function & u,
Parameters parameters = empty_parameters )

Solve linear variational problem a(u, v) == L(v) or nonlinear variational problem F(u; v) = 0 without boundary conditions.

Optional parameters can be passed to the LinearVariationalSolver or NonlinearVariationalSolver classes.

◆ solve() [9/12]

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

◆ solve() [10/12]

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

◆ solve() [11/12]

void dolfin::solve ( const Equation & equation,
Function & u,
std::vector< const DirichletBC * > bcs,
const Form & J,
Parameters parameters = empty_parameters )

Solve nonlinear variational problem F(u; v) == 0 with a list of boundary conditions. The argument J should provide the Jacobian bilinear form J = dF/du.

Optional parameters can be passed to the LinearVariationalSolver or NonlinearVariationalSolver classes.

◆ solve() [12/12]

void dolfin::solve ( const Equation & equation,
Function & u,
std::vector< const DirichletBC * > bcs,
Parameters parameters = empty_parameters )

Solve linear variational problem a(u, v) == L(v) or nonlinear variational problem F(u; v) = 0 with a list of boundary conditions.

Optional parameters can be passed to the LinearVariationalSolver or NonlinearVariationalSolver classes.

◆ timing()

std::tuple< std::size_t, double, double, double > dolfin::timing ( std::string task,
TimingClear clear )

Return timing (count, total wall time, total user time, total system time) for given task, optionally clearing all timings for the task

Arguments task (std::string) name of a task clear (TimingClear)

  • TimingClear::clear resets stored timings
  • TimingClear::keep leaves stored timings intact

Returns std::tuple<std::size_t, double, double, double> (count, total wall time, total user time, total system time)

◆ timings()

Table dolfin::timings ( TimingClear clear,
std::set< TimingType > type )

Return a summary of timings and tasks in a Table, optionally clearing stored timings

Arguments clear (TimingClear)

  • TimingClear::clear resets stored timings
  • TimingClear::keep leaves stored timings intact type (std::set<TimingType>) subset of { TimingType::wall, TimingType::user, TimingType::system }

Returns Table Table with timings

◆ vertex_to_dof_map()

std::vector< dolfin::la_index > dolfin::vertex_to_dof_map ( const FunctionSpace & space)

Return a map between vertex indices and dof indices

Only works for FunctionSpace with dofs exclusively on vertices. For mixed FunctionSpaces dof index is offset with the number of dofs per vertex.

Parameters
space(FunctionSpace) The FunctionSpace for what the vertex to dof map should be computed for
Returns
std::vector<dolfin::la_index> The vertex to dof map

Variable Documentation

◆ rand_seeded

bool dolfin::rand_seeded = false

Flag to determine whether to reseed dolfin::rand(). Normally on first call.