16#ifndef dealii_mesh_worker_local_results_h
17#define dealii_mesh_worker_local_results_h
80 template <
typename number>
150 matrix(
const unsigned int i,
const bool external =
false);
159 matrix(
const unsigned int i,
const bool external =
false)
const;
214 template <
typename MatrixType>
226 template <
typename MatrixType>
246 template <
typename StreamType>
260 std::vector<number>
J;
266 std::vector<BlockVector<number>>
R;
272 std::vector<MatrixBlock<FullMatrix<number>>>
M1;
280 std::vector<MatrixBlock<FullMatrix<number>>>
M2;
290 template <
typename number>
298 template <
typename number>
306 template <
typename number>
307 template <
typename MatrixType>
313 M1.resize(matrices.
size());
315 M2.resize(matrices.
size());
316 for (
unsigned int i = 0; i < matrices.
size(); ++i)
318 const unsigned int row = matrices.
block(i).
row;
332 template <
typename number>
333 template <
typename MatrixType>
339 M1.resize(matrices.
size());
341 M2.resize(matrices.
size());
342 for (
unsigned int i = 0; i < matrices.
size(); ++i)
345 const unsigned int row = o[o.
min_level()].row;
346 const unsigned int col = o[o.
min_level()].column;
359 template <
typename number>
367 for (
unsigned int i = 0; i < n; ++i)
380 template <
typename number>
383 const unsigned int nv)
385 quadrature_data.reinit(np, nv);
389 template <
typename number>
397 template <
typename number>
405 template <
typename number>
413 template <
typename number>
417 return quadrature_data.n_rows();
421 template <
typename number>
425 return quadrature_data.n_cols();
429 template <
typename number>
438 template <
typename number>
447 template <
typename number>
461 template <
typename number>
464 const unsigned int i)
466 return quadrature_data(k, i);
470 template <
typename number>
474 return quadrature_data;
478 template <
typename number>
487 template <
typename number>
496 template <
typename number>
510 template <
typename number>
513 const unsigned int i)
const
515 return quadrature_data(k, i);
519 template <
typename number>
520 template <
typename StreamType>
524 os <<
"J: " << J.size() << std::endl;
525 os <<
"R: " << R.size() << std::endl;
526 for (
unsigned int i = 0; i < R.size(); ++i)
528 os <<
" " << R[i].n_blocks() <<
" -";
529 for (
unsigned int j = 0; j < R[i].n_blocks(); ++j)
530 os <<
' ' << R[i].block(j).size();
533 os <<
"M: " << M1.size() <<
" face " << M2.size() << std::endl;
534 for (
unsigned int i = 0; i < M1.size(); ++i)
536 os <<
" " << M1[i].row <<
',' << M1[i].column <<
' '
537 << M1[i].matrix.m() <<
'x' << M1[i].matrix.n();
539 os <<
" face " << M2[i].row <<
',' << M2[i].column <<
' '
540 << M2[i].matrix.m() <<
'x' << M2[i].matrix.n();
unsigned int min_level() const
const value_type & block(size_type i) const
unsigned int size() const
unsigned int size() const
Number of stored data objects.
const value_type & block(size_type i) const
MatrixBlock< FullMatrix< number > > & matrix(const unsigned int i, const bool external=false)
unsigned int n_quadrature_values() const
BlockVector< number > & vector(const unsigned int i)
unsigned int n_vectors() const
number & quadrature_value(const unsigned int k, const unsigned int i)
number value(const unsigned int i) const
std::vector< BlockVector< number > > R
std::size_t memory_consumption() const
const BlockVector< number > & vector(const unsigned int i) const
const MatrixBlock< FullMatrix< number > > & matrix(const unsigned int i, const bool external=false) const
number quadrature_value(const unsigned int k, const unsigned int i) const
Table< 2, number > quadrature_data
void initialize_matrices(const MatrixBlockVector< MatrixType > &matrices, bool both)
void initialize_matrices(const MGMatrixBlockVector< MatrixType > &matrices, bool both)
std::vector< MatrixBlock< FullMatrix< number > > > M2
void reinit(const BlockIndices &local_sizes)
void initialize_matrices(const unsigned int n, bool both)
std::vector< MatrixBlock< FullMatrix< number > > > M1
unsigned int n_values() const
void initialize_quadrature(const unsigned int np, const unsigned int nv)
Table< 2, number > & quadrature_values()
void initialize_numbers(const unsigned int n)
void initialize_vectors(const unsigned int n)
unsigned int n_matrices() const
void print_debug(StreamType &os) const
unsigned int n_quadrature_points() const
number & value(const unsigned int i)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define AssertIndexRange(index, range)