DUNE-FEM (unstable)
A local matrix with a small array as storage. More...
#include <dune/fem/operator/common/temporarylocalmatrix.hh>
Public Types | |
| typedef std::vector< TemporaryLocalMatrix< DomainSpaceImp, RangeSpaceImp > > | RowType | 
| remember the value type  | |
| typedef ThisType | LocalMatrixInterfaceType | 
| type of this interface  | |
| typedef Traits::LocalMatrixType | LocalMatrixType | 
| type of local matrix implementation  | |
| typedef Traits::LittleBlockType | LittleBlockType | 
Public Member Functions | |
| template<class DomainEntityType , class RangeEntityType > | |
| void | init (const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity) | 
| initialize the local matrix to entities  More... | |
| template<class DomainEntityType , class RangeEntityType > | |
| void | bind (const DomainEntityType &domainEntity, const RangeEntityType &rangeEntity) | 
| initialize the local matrix to entities  More... | |
| void | unbind () | 
| clear local matrix from entities  More... | |
| void | add (const int localRow, const int localCol, const RangeFieldType &value) | 
| add value to matrix entry (row,col) where row and col are local row and local column  More... | |
| void | set (const int localRow, const int localCol, const RangeFieldType &value) | 
| set value of matrix entry (row,col) where row and col are local row and local column  More... | |
| const RangeFieldType | get (const int localRow, const int localCol) const | 
| get value of matrix entry (row,col) where row and col are local row and local column  More... | |
| void | scale (const RangeFieldType &value) | 
| scale matrix with scalar value  More... | |
| void | clear () | 
| set all entries of local matrix to zero  More... | |
| void | clearRow (const int localRow) | 
| set row to zero values  More... | |
| const DomainBasisFunctionSetType & | domainBasisFunctionSet () const | 
| access to the base function set within the domain space  More... | |
| const RangeBasisFunctionSetType & | rangeBasisFunctionSet () const | 
| access to the base function set within the range space  More... | |
| void | multiply_AT_A (const DenseMatrix &A) | 
| this = A^T * A  | |
| void | scale (const TemporaryLocalMatrix< DomainSpaceImp, RangeSpaceImp > &val) | 
| scale matrix with scalar  | |
| DenseMatrix< TemporaryLocalMatrix< DomainSpaceImp, RangeSpaceImp > > & | operator+= (const DenseMatrix &org) | 
| add matrix  | |
| DenseMatrix< TemporaryLocalMatrix< DomainSpaceImp, RangeSpaceImp > > & | operator-= (const DenseMatrix &org) | 
| substract matrix  | |
| void | resort () | 
| resort ordering in global matrix (if possible)  More... | |
| void | finalize () | 
| finalize local matrix setup and possibly add values to real matrix  More... | |
| const DomainSpaceType & | domainSpace () const | 
| access to the domain space  More... | |
| const RangeSpaceType & | rangeSpace () const | 
| access to the range space  More... | |
| void | multiplyAdd (const DomainLocalFunctionType &lhs, RangeLocalFunctionType &rhs) const | 
| multiply left hand side with local matrix and add to right hand side rhs += Matrix * lhs  More... | |
| void | clearCol (const int localCol) | 
| ser column entries to zero  More... | |
| MatrixColumnType | column (const unsigned int col) | 
| return column object for local matrix which contains axpy methods for convenience  More... | |
Detailed Description
class Dune::Fem::TemporaryLocalMatrix< DomainSpaceImp, RangeSpaceImp >
A local matrix with a small array as storage.
A TemporaryLocalMatrix is an implementation of the LocalMatrixInterface storing the matrix values in an array. It is useful when generating multiple local matrices that shall then be added together.
- Note
 - Due to the backing array, accesses to the matrix should be very fast.
 
- Parameters
 - 
  
DomainSpaceImp DiscreteFunctionSpace modelling the domain RangeSpaceImp DiscreteFunctionSpace modelling the range  
Member Typedef Documentation
◆ LittleBlockType
      
  | 
  inherited | 
type of block (i.e. FieldMatrix for BlockMatrices
Member Function Documentation
◆ add()
      
  | 
  inline | 
add value to matrix entry (row,col) where row and col are local row and local column
- Parameters
 - 
  
[in] localRow local row [in] localCol local column [in] value value to add  
◆ bind()
      
  | 
  inline | 
initialize the local matrix to entities
- Parameters
 - 
  
[in] domainEntity entity within grid of domain space, [in] rangeEntity entity within grid of range space  
References Dune::Fem::LocalMatrixDefault< LocalMatrixTraits >::bind().
Referenced by Dune::Fem::TemporaryLocalMatrix< DomainSpaceImp, RangeSpaceImp >::init().
◆ clear()
      
  | 
  inline | 
set all entries of local matrix to zero
◆ clearCol()
      
  | 
  inlineinherited | 
ser column entries to zero
- Parameters
 - 
  
[in] localCol local column that is set to zero  
◆ clearRow()
      
  | 
  inline | 
set row to zero values
- Parameters
 - 
  
[in] localRow local row that is set to zero  
◆ column()
      
  | 
  inlineinherited | 
return column object for local matrix which contains axpy methods for convenience
- Parameters
 - 
  
col local column number  
- Returns
 - object of type MatrixColumnObject
 
◆ domainBasisFunctionSet()
      
  | 
  inline | 
access to the base function set within the domain space
◆ domainSpace()
      
  | 
  inlineinherited | 
access to the domain space
◆ finalize()
      
  | 
  inlineinherited | 
finalize local matrix setup and possibly add values to real matrix
◆ get()
      
  | 
  inline | 
get value of matrix entry (row,col) where row and col are local row and local column
- Parameters
 - 
  
[in] localRow local row [in] localCol local column  
- Returns
 - value of matrix entry
 
◆ init()
      
  | 
  inline | 
initialize the local matrix to entities
- Parameters
 - 
  
[in] domainEntity entity within grid of domain space, [in] rangeEntity entity within grid of range space  
References Dune::Fem::TemporaryLocalMatrix< DomainSpaceImp, RangeSpaceImp >::bind().
◆ multiplyAdd()
      
  | 
  inlineinherited | 
multiply left hand side with local matrix and add to right hand side rhs += Matrix * lhs
- Parameters
 - 
  
[in] lhs left hand side [out] rhs right hand side  
◆ rangeBasisFunctionSet()
      
  | 
  inline | 
access to the base function set within the range space
◆ rangeSpace()
      
  | 
  inlineinherited | 
access to the range space
◆ resort()
      
  | 
  inlineinherited | 
resort ordering in global matrix (if possible)
◆ scale()
      
  | 
  inline | 
scale matrix with scalar value
- Parameters
 - 
  
[in] scalar scalar value that scales the matrix  
References Dune::size().
◆ set()
      
  | 
  inline | 
set value of matrix entry (row,col) where row and col are local row and local column
- Parameters
 - 
  
[in] localRow local row [in] localCol local column [in] value value to set  
◆ unbind()
      
  | 
  inline | 
clear local matrix from entities
References Dune::Fem::LocalMatrixDefault< LocalMatrixTraits >::unbind().
The documentation for this class was generated from the following file:
- dune/fem/operator/common/temporarylocalmatrix.hh
 
   | 
                                Legal Statements / Impressum  | 
                                Hosted by  TU Dresden & Uni Heidelberg  | 
				  generated with Hugo v0.111.3
								(Nov 3, 23:36, 2025)