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

#include <Assembler.h>

Inheritance diagram for dolfin::Assembler:
Inheritance graph
[legend]
Collaboration diagram for dolfin::Assembler:
Collaboration graph
[legend]

Public Member Functions

 Assembler ()
 Constructor.
 
void assemble (GenericTensor &A, const Form &a)
 
void assemble_cells (GenericTensor &A, const Form &a, UFC &ufc, std::shared_ptr< const MeshFunction< std::size_t >> domains, std::vector< double > *values)
 
void assemble_exterior_facets (GenericTensor &A, const Form &a, UFC &ufc, std::shared_ptr< const MeshFunction< std::size_t >> domains, std::vector< double > *values)
 
void assemble_interior_facets (GenericTensor &A, const Form &a, UFC &ufc, std::shared_ptr< const MeshFunction< std::size_t >> domains, std::shared_ptr< const MeshFunction< std::size_t >> cell_domains, std::vector< double > *values)
 
void assemble_vertices (GenericTensor &A, const Form &a, UFC &ufc, std::shared_ptr< const MeshFunction< std::size_t >> domains)
 
- Public Member Functions inherited from dolfin::AssemblerBase
 AssemblerBase ()
 Constructor.
 
void init_global_tensor (GenericTensor &A, const Form &a)
 

Additional Inherited Members

- Public Attributes inherited from dolfin::AssemblerBase
bool add_values
 
bool finalize_tensor
 
bool keep_diagonal
 
- Static Protected Member Functions inherited from dolfin::AssemblerBase
static void check (const Form &a)
 Check form.
 
static std::string progress_message (std::size_t rank, std::string integral_type)
 Pretty-printing for progress bar.
 

Detailed Description

This class provides automated assembly of linear systems, or more generally, assembly of a sparse tensor from a given variational form.

Subdomains for cells and facets may be specified by assigning subdomain indicators specified by MeshFunction to the Form being assembled:

form.dx = cell_domains
form.ds = exterior_facet_domains
form.dS = interior_facet_domains

Member Function Documentation

void Assembler::assemble ( GenericTensor A,
const Form a 
)

Assemble tensor from given form

Parameters
[out]A(GenericTensor) The tensor to assemble.
[in]a(Form&) The form to assemble the tensor from.
void Assembler::assemble_cells ( GenericTensor A,
const Form a,
UFC ufc,
std::shared_ptr< const MeshFunction< std::size_t >>  domains,
std::vector< double > *  values 
)

Assemble tensor from given form over cells. This function is provided for users who wish to build a customized assembler.

Parameters
[out]A(GenericTensor&) The tensor to assemble.
[in]a(Form&) The form to assemble the tensor from.
[in]ufc(UFC&)
[in]domains(MeshFunction<std::size_t>)
[in]values(std::vector<double>*)
void Assembler::assemble_exterior_facets ( GenericTensor A,
const Form a,
UFC ufc,
std::shared_ptr< const MeshFunction< std::size_t >>  domains,
std::vector< double > *  values 
)

Assemble tensor from given form over exterior facets. This function is provided for users who wish to build a customized assembler.

Parameters
[out]A(GenericTensor) The tensor to assemble.
[in]a(Form) The form to assemble the tensor from.
[in]ufc(UFC)
[in]domains(MeshFunction)
[in]values(std::vector<double>*)
void Assembler::assemble_interior_facets ( GenericTensor A,
const Form a,
UFC ufc,
std::shared_ptr< const MeshFunction< std::size_t >>  domains,
std::shared_ptr< const MeshFunction< std::size_t >>  cell_domains,
std::vector< double > *  values 
)

Assemble tensor from given form over interior facets. This function is provided for users who wish to build a customized assembler.

Parameters
[out]A(GenericTensor) The tensor to assemble.
[in]a(Form) The form to assemble the tensor from.
[in]ufc(UFC)
[in]domains(MeshFunction<std::size_t>)
[in]cell_domains(MeshFunction<std::size_t>)
[in]values(std::vector<double>*)
void Assembler::assemble_vertices ( GenericTensor A,
const Form a,
UFC ufc,
std::shared_ptr< const MeshFunction< std::size_t >>  domains 
)

Assemble tensor from given form over vertices. This function is provided for users who wish to build a customized assembler.

Parameters
[out]A(GenericTersn) The tensor to assemble.
[in]a(Form) The form to assemble the tensor from.
[in]ufc(UFC)
[in]domains(MeshFunction)

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