Serac  0.1
Serac is an implicit thermal strucural mechanics simulation code.
Classes | Public Member Functions | Public Attributes | List of all members
serac::GradientAssemblyLookupTables Struct Reference

this object figures out the sparsity pattern associated with a finite element discretization of the given test and trial function spaces, and records which nonzero each element "stiffness" matrix maps to, to facilitate assembling the element matrices into the global sparse matrix. e.g. More...

#include <dof_numbering.hpp>

Classes

struct  Entry
 a type for representing a nonzero entry in a sparse matrix More...
 

Public Member Functions

 GradientAssemblyLookupTables (const serac::BlockElementRestriction &block_test_dofs, const serac::BlockElementRestriction &block_trial_dofs)
 create lookup tables describing which degrees of freedom correspond to each domain/boundary element More...
 
uint32_t operator() (int i, int j) const
 return the index (into the nonzero entries) corresponding to entry (i,j) More...
 

Public Attributes

uint32_t nnz
 how many nonzero entries appear in the sparse matrix
 
std::vector< int > row_ptr
 array holding the offsets for a given row of the sparse matrix i.e. row r corresponds to the indices [row_ptr[r], row_ptr[r+1])
 
std::vector< int > col_ind
 array holding the column associated with each nonzero entry
 
std::unordered_map< Entry, uint32_t, Entry::Hashernz_LUT
 nz_LUT returns the index of the col_ind / value CSR arrays corresponding to the (i,j) entry
 

Detailed Description

this object figures out the sparsity pattern associated with a finite element discretization of the given test and trial function spaces, and records which nonzero each element "stiffness" matrix maps to, to facilitate assembling the element matrices into the global sparse matrix. e.g.

element_nonzero_LUT(e, i, j) says where (in the global sparse matrix) to put the (i,j) component of the matrix associated with element element matrix e

Note: due to an internal inconsistency between mfem::FiniteElementSpace and mfem::FaceRestriction, we choose to use the Restriction operator as the "source of truth", since we are also using its convention for quadrature point numbering.

Definition at line 253 of file dof_numbering.hpp.

Constructor & Destructor Documentation

◆ GradientAssemblyLookupTables()

serac::GradientAssemblyLookupTables::GradientAssemblyLookupTables ( const serac::BlockElementRestriction block_test_dofs,
const serac::BlockElementRestriction block_trial_dofs 
)
inline

create lookup tables describing which degrees of freedom correspond to each domain/boundary element

Parameters
block_test_dofsobject containing information about dofs for the test space
block_trial_dofsobject containing information about dofs for the trial space

Definition at line 287 of file dof_numbering.hpp.

Member Function Documentation

◆ operator()()

uint32_t serac::GradientAssemblyLookupTables::operator() ( int  i,
int  j 
) const
inline

return the index (into the nonzero entries) corresponding to entry (i,j)

Parameters
ithe row
jthe column

Definition at line 352 of file dof_numbering.hpp.


The documentation for this struct was generated from the following file: