DOLFIN
DOLFIN C++ interface
Static Public Member Functions | List of all members
dolfin::LocalAssembler Class Reference

#include <LocalAssembler.h>

Static Public Member Functions

static void assemble (Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &A, UFC &ufc, const std::vector< double > &coordinate_dofs, ufc::cell &ufc_cell, const Cell &cell, const MeshFunction< std::size_t > *cell_domains, const MeshFunction< std::size_t > *exterior_facet_domains, const MeshFunction< std::size_t > *interior_facet_domains)
 
static void assemble_cell (Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &A, UFC &ufc, const std::vector< double > &coordinate_dofs, const ufc::cell &ufc_cell, const Cell &cell, const MeshFunction< std::size_t > *cell_domains)
 
static void assemble_exterior_facet (Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &A, UFC &ufc, const std::vector< double > &coordinate_dofs, const ufc::cell &ufc_cell, const Cell &cell, const Facet &facet, const std::size_t local_facet, const MeshFunction< std::size_t > *exterior_facet_domains)
 
static void assemble_interior_facet (Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &A, UFC &ufc, const std::vector< double > &coordinate_dofs, const ufc::cell &ufc_cell, const Cell &cell, const Facet &facet, const std::size_t local_facet, const MeshFunction< std::size_t > *interior_facet_domains, const MeshFunction< std::size_t > *cell_domains)
 

Detailed Description

Assembly of local cell tensors. Used by the adaptivity and LocalSolver functionality in dolfin. The local assembly functionality provided here is also wrapped as a free function assemble_local(form_a, cell) in Python for easier usage. Use from the C++ interface defined below will be faster than the free function as fewer objects need to be created and destroyed.

Member Function Documentation

void LocalAssembler::assemble ( Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &  A,
UFC ufc,
const std::vector< double > &  coordinate_dofs,
ufc::cell &  ufc_cell,
const Cell cell,
const MeshFunction< std::size_t > *  cell_domains,
const MeshFunction< std::size_t > *  exterior_facet_domains,
const MeshFunction< std::size_t > *  interior_facet_domains 
)
static

Assemble a local tensor on a cell. Internally calls assemble_cell(), assemble_exterior_facet(), assemble_interior_facet().

Parameters
[out]AThe tensor to assemble.
[in]ufc
[in]coordinate_dofs
[in]ufc_cell
[in]cell
[in]cell_domains
[in]exterior_facet_domains
[in]interior_facet_domains
void LocalAssembler::assemble_cell ( Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &  A,
UFC ufc,
const std::vector< double > &  coordinate_dofs,
const ufc::cell &  ufc_cell,
const Cell cell,
const MeshFunction< std::size_t > *  cell_domains 
)
static

Worker method called by assemble() to perform assembly of volume integrals (UFL measure dx).

Parameters
[out]AThe tensor to assemble.
[in]ufc
[in]coordinate_dofs
[in]ufc_cell
[in]cell
[in]cell_domains
void LocalAssembler::assemble_exterior_facet ( Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &  A,
UFC ufc,
const std::vector< double > &  coordinate_dofs,
const ufc::cell &  ufc_cell,
const Cell cell,
const Facet facet,
const std::size_t  local_facet,
const MeshFunction< std::size_t > *  exterior_facet_domains 
)
static

Worker method called by assemble() for each of the cell's external facets to perform assembly of external facet integrals (UFL measure ds).

Parameters
[out]AThe tensor to assemble.
[in]ufc
[in]coordinate_dofs
[in]ufc_cell
[in]cell
[in]facet
[in]local_facet
[in]exterior_facet_domains
void LocalAssembler::assemble_interior_facet ( Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic, Eigen::RowMajor > &  A,
UFC ufc,
const std::vector< double > &  coordinate_dofs,
const ufc::cell &  ufc_cell,
const Cell cell,
const Facet facet,
const std::size_t  local_facet,
const MeshFunction< std::size_t > *  interior_facet_domains,
const MeshFunction< std::size_t > *  cell_domains 
)
static

Worker method called by assemble() for each of the cell's internal facets to perform assembly of internal facet integrals (UFL measure dS)

Parameters
[out]AThe tensor to assemble.
[in]ufc
[in]coordinate_dofs
[in]ufc_cell
[in]cell
[in]facet
[in]local_facet
[in]interior_facet_domains
[in]cell_domains

The documentation for this class was generated from the following files: