
Public Member Functions | |
| unsigned int | n_values () const |
| unsigned int | n_vectors () const |
| unsigned int | n_matrices () const |
| unsigned int | n_quadrature_points () const |
| unsigned int | n_quadrature_values () const |
| number & | value (unsigned int i) |
| number | value (unsigned int i) const |
| BlockVector< number > & | vector (unsigned int i) |
| const BlockVector< number > & | vector (unsigned int i) const |
| MatrixBlock< FullMatrix < number > > & | matrix (unsigned int i, bool external=false) |
| const MatrixBlock< FullMatrix < number > > & | matrix (unsigned int i, bool external=false) const |
| Table< 2, number > & | quadrature_values () |
| number & | quadrature_value (unsigned int k, unsigned int i) |
| number | quadrature_value (unsigned int k, unsigned int i) const |
| void | initialize_numbers (const unsigned int n) |
| void | initialize_vectors (const unsigned int n) |
| void | initialize_matrices (unsigned int n, bool both) |
| template<class MATRIX > | |
| void | initialize_matrices (const MatrixBlockVector< MATRIX > &matrices, bool both) |
| template<class MATRIX > | |
| void | initialize_matrices (const MGMatrixBlockVector< MATRIX > &matrices, bool both) |
| void | initialize_quadrature (unsigned int np, unsigned int nv) |
| void | reinit (const BlockIndices &local_sizes) |
| std::size_t | memory_consumption () const |
Private Member Functions | |
| void | initialize_local (MatrixBlock< FullMatrix< number > > &M, const unsigned int row, const unsigned int col) |
Private Attributes | |
| std::vector< number > | J |
| std::vector< BlockVector < number > > | R |
| std::vector< MatrixBlock < FullMatrix< number > > > | M1 |
| std::vector< MatrixBlock < FullMatrix< number > > > | M2 |
| Table< 2, number > | quadrature_data |
The class providing the scrapbook to fill with results of local integration. Depending on the tesk the mesh worker loop is performing, local results can be of different types. They have in common that they are the result of local integration over a cell or face. Their actual type is determined by the Assember using them. It is also the assembler setting the arrays of local results to the sizes needed. Here is a list of the provided data types and the assembers using them:
n_values() numbers accessed with value(), and stored in the data member J.
false. These are stored in M1, and they are the matrices coupling degrees of freedom in the same cell. For fluxes across faces, there is an additional set M2 of matrices of the same size, but the dimension of the matrices being according to the degrees of freedom on both cells. These are accessed with matrix(), using the second argument true. The local matrices initialized by reinit() of the info object and then assembled into the global system by Assembler classes.
Definition at line 209 of file local_results.h.
| unsigned int MeshWorker::LocalResults< number >::n_values | ( | ) | const [inline] |
The number of scalar values.
This number is set to a nonzero value by Assember::CellsAndFaces
Definition at line 573 of file local_results.h.
| unsigned int MeshWorker::LocalResults< number >::n_vectors | ( | ) | const [inline] |
The number of vectors.
This number is set to a nonzero value by Assember::ResidualSimple and Assember::ResidualLocalBlocksToGlobalBlocks.
Definition at line 582 of file local_results.h.
| unsigned int MeshWorker::LocalResults< number >::n_matrices | ( | ) | const [inline] |
The number of matrices.
Definition at line 591 of file local_results.h.
| unsigned int MeshWorker::LocalResults< number >::n_quadrature_points | ( | ) | const [inline] |
The number of quadrature points in quadrature_values().
Definition at line 600 of file local_results.h.
Referenced by MeshWorker::Assembler::GnuplotPatch::assemble().
| unsigned int MeshWorker::LocalResults< number >::n_quadrature_values | ( | ) | const [inline] |
The number of values in each quadrature point in quadrature_values().
Definition at line 609 of file local_results.h.
Referenced by MeshWorker::Assembler::GnuplotPatch::assemble().
| number & MeshWorker::LocalResults< number >::value | ( | unsigned int | i ) | [inline] |
Access scalar value at index i.
Definition at line 618 of file local_results.h.
References AssertIndexRange.
| number MeshWorker::LocalResults< number >::value | ( | unsigned int | i ) | const [inline] |
Read scalar value at index i.
Definition at line 671 of file local_results.h.
References AssertIndexRange.
| BlockVector< number > & MeshWorker::LocalResults< number >::vector | ( | unsigned int | i ) | [inline] |
Access vector at index i.
Definition at line 628 of file local_results.h.
References AssertIndexRange.
| const BlockVector< number > & MeshWorker::LocalResults< number >::vector | ( | unsigned int | i ) | const [inline] |
Read vector at index i.
Definition at line 681 of file local_results.h.
References AssertIndexRange.
| MatrixBlock< FullMatrix< number > > & MeshWorker::LocalResults< number >::matrix | ( | unsigned int | i, |
| bool | external = false |
||
| ) | [inline] |
Access matrix at index i. For results on internal faces, a true value for external refers to the flux between cells, while false refers to entries coupling inside the cell.
Definition at line 638 of file local_results.h.
References AssertIndexRange.
| const MatrixBlock< FullMatrix< number > > & MeshWorker::LocalResults< number >::matrix | ( | unsigned int | i, |
| bool | external = false |
||
| ) | const [inline] |
Read matrix at index i. For results on internal faces, a true value for external refers to the flux between cells, while false refers to entries coupling inside the cell.
Definition at line 691 of file local_results.h.
References AssertIndexRange.
| Table< 2, number > & MeshWorker::LocalResults< number >::quadrature_values | ( | ) | [inline] |
Access to the vector quadrature_data of data in quadrature points, organized such that there is a vector for each point, containing one entry for each component.
Definition at line 662 of file local_results.h.
| number & MeshWorker::LocalResults< number >::quadrature_value | ( | unsigned int | k, |
| unsigned int | i | ||
| ) | [inline] |
Access the ith value at quadrature point k
Definition at line 653 of file local_results.h.
Referenced by MeshWorker::Assembler::GnuplotPatch::assemble().
| number MeshWorker::LocalResults< number >::quadrature_value | ( | unsigned int | k, |
| unsigned int | i | ||
| ) | const [inline] |
Read the ith value at quadrature point k
Definition at line 706 of file local_results.h.
| void MeshWorker::LocalResults< number >::initialize_numbers | ( | const unsigned int | n ) | [inline] |
Initialize the vector with scalar values.
Definition at line 474 of file local_results.h.
| void MeshWorker::LocalResults< number >::initialize_vectors | ( | const unsigned int | n ) | [inline] |
Initialize the vector with vector values.
Definition at line 482 of file local_results.h.
| void MeshWorker::LocalResults< number >::initialize_matrices | ( | unsigned int | n, |
| bool | both | ||
| ) | [inline] |
Allocate n local matrices. Additionally, set their block row and column coordinates to zero. The matrices themselves are resized by reinit().
The template parameter MatrixPtr should point to a MatrixBlock instantiation in order to provide row and column info.
Definition at line 543 of file local_results.h.
| void MeshWorker::LocalResults< number >::initialize_matrices | ( | const MatrixBlockVector< MATRIX > & | matrices, |
| bool | both | ||
| ) | [inline] |
Allocate a local matrix for each of the global ones in matrices. Additionally, set their block row and column coordinates. The matrices themselves are resized by reinit().
Definition at line 491 of file local_results.h.
References MatrixBlockVector< MATRIX >::block(), and NamedData< std_cxx1x::shared_ptr< MatrixBlock< MATRIX > > >::size().
| void MeshWorker::LocalResults< number >::initialize_matrices | ( | const MGMatrixBlockVector< MATRIX > & | matrices, |
| bool | both | ||
| ) | [inline] |
Allocate a local matrix for each of the global level objects in matrices. Additionally, set their block row and column coordinates. The matrices themselves are resized by reinit().
Definition at line 517 of file local_results.h.
References MGMatrixBlockVector< MATRIX >::block(), MGLevelObject< Object >::min_level(), and MGMatrixBlockVector< MATRIX >::size().
| void MeshWorker::LocalResults< number >::initialize_quadrature | ( | unsigned int | np, |
| unsigned int | nv | ||
| ) | [inline] |
Initialize quadrature values to nv values in np quadrature points.
Definition at line 564 of file local_results.h.
Referenced by MeshWorker::Assembler::GnuplotPatch::initialize_info().
| void MeshWorker::LocalResults< number >::reinit | ( | const BlockIndices & | local_sizes ) |
Reinitialize matrices for new cell. Does not resize any of the data vectors stored in this object, but resizes the vectors in R and the matrices in M1 and M2 for hp and sets them to zero.
Referenced by MeshWorker::DoFInfo< dim, spacedim, number >::reinit().
| std::size_t MeshWorker::LocalResults< number >::memory_consumption | ( | ) | const |
The memory used by this object.
| void MeshWorker::LocalResults< number >::initialize_local | ( | MatrixBlock< FullMatrix< number > > & | M, |
| const unsigned int | row, | ||
| const unsigned int | col | ||
| ) | [private] |
Initialize a single local matrix block. A helper function for initialize()
std::vector<number> MeshWorker::LocalResults< number >::J [private] |
The local numbers, computed on a cell or on a face.
Definition at line 430 of file local_results.h.
std::vector<BlockVector<number> > MeshWorker::LocalResults< number >::R [private] |
The local vectors. This field is public, so that local integrators can write to it.
Definition at line 438 of file local_results.h.
std::vector<MatrixBlock<FullMatrix<number> > > MeshWorker::LocalResults< number >::M1 [private] |
The local matrices coupling degrees of freedom in the cell itself or within the first cell on a face.
Definition at line 447 of file local_results.h.
std::vector<MatrixBlock<FullMatrix<number> > > MeshWorker::LocalResults< number >::M2 [private] |
The local matrices coupling test functions on the cell with trial functions on the other cell.
Only used on interior faces.
Definition at line 459 of file local_results.h.
Table<2, number> MeshWorker::LocalResults< number >::quadrature_data [private] |
Values in quadrature points.
Definition at line 467 of file local_results.h.
documentation generated on Fri Feb 3 2012 06:04:15 by
doxygen
1.7.2