Reference documentation for deal.II version 8.5.1
local_results.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2006 - 2016 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #ifndef dealii__mesh_worker_local_results_h
18 #define dealii__mesh_worker_local_results_h
19 
20 #include <deal.II/base/config.h>
21 #include <deal.II/base/std_cxx11/function.h>
22 #include <deal.II/base/geometry_info.h>
23 #include <deal.II/lac/matrix_block.h>
24 #include <deal.II/lac/block_vector.h>
25 #include <deal.II/meshworker/vector_selector.h>
26 
27 DEAL_II_NAMESPACE_OPEN
28 
29 class BlockIndices;
30 
174 namespace MeshWorker
175 {
206  template <typename number>
208  {
209  public:
216  unsigned int n_values () const;
217 
224  unsigned int n_vectors () const;
225 
229  unsigned int n_matrices () const;
230 
234  unsigned int n_quadrature_points() const;
235 
239  unsigned int n_quadrature_values() const;
240 
244  number &value(unsigned int i);
245 
249  number value(unsigned int i) const;
250 
254  BlockVector<number> &vector(unsigned int i);
255 
259  const BlockVector<number> &vector(unsigned int i) const;
260 
266  MatrixBlock<FullMatrix<number> > &matrix(unsigned int i, bool external = false);
267 
273  const MatrixBlock<FullMatrix<number> > &matrix(unsigned int i, bool external = false) const;
274 
281 
285  number &quadrature_value(unsigned int k, unsigned int i);
286 
290  number quadrature_value(unsigned int k, unsigned int i) const;
291 
297  void initialize_numbers(const unsigned int n);
298 
304  void initialize_vectors(const unsigned int n);
305 
313  void initialize_matrices(unsigned int n, bool both);
314 
322  template <typename MatrixType>
324  bool both);
325 
333  template <typename MatrixType>
335  bool both);
336 
341  void initialize_quadrature(unsigned int np, unsigned int nv);
342 
348  void reinit(const BlockIndices &local_sizes);
349 
350  template <class StreamType>
351  void print_debug(StreamType &os) const;
352 
356  std::size_t memory_consumption () const;
357 
358  private:
364  const unsigned int row,
365  const unsigned int col);
366 
370  std::vector<number> J;
371 
376  std::vector<BlockVector<number> > R;
377 
382  std::vector<MatrixBlock<FullMatrix<number> > > M1;
383 
390  std::vector<MatrixBlock<FullMatrix<number> > > M2;
391 
396  };
397 
398 //----------------------------------------------------------------------//
399 
400  template <typename number>
401  inline void
403  {
404  J.resize(n);
405  }
406 
407 
408  template <typename number>
409  inline void
411  {
412  R.resize(n);
413  }
414 
415 
416  template <typename number>
417  template <typename MatrixType>
418  inline void
421  bool both)
422  {
423  M1.resize(matrices.size());
424  if (both)
425  M2.resize(matrices.size());
426  for (unsigned int i=0; i<matrices.size(); ++i)
427  {
428  const unsigned int row = matrices.block(i).row;
429  const unsigned int col = matrices.block(i).column;
430 
431  M1[i].row = row;
432  M1[i].column = col;
433  if (both)
434  {
435  M2[i].row = row;
436  M2[i].column = col;
437  }
438  }
439  }
440 
441 
442  template <typename number>
443  template <typename MatrixType>
444  inline void
447  bool both)
448  {
449  M1.resize(matrices.size());
450  if (both)
451  M2.resize(matrices.size());
452  for (unsigned int i=0; i<matrices.size(); ++i)
453  {
454  const MGLevelObject<MatrixBlock<MatrixType> > &o = matrices.block(i);
455  const unsigned int row = o[o.min_level()].row;
456  const unsigned int col = o[o.min_level()].column;
457 
458  M1[i].row = row;
459  M1[i].column = col;
460  if (both)
461  {
462  M2[i].row = row;
463  M2[i].column = col;
464  }
465  }
466  }
467 
468 
469  template <typename number>
470  inline void
472  const bool both)
473  {
474  M1.resize(n);
475  if (both)
476  M2.resize(n);
477  for (unsigned int i=0; i<n; ++i)
478  {
479  M1[i].row = 0;
480  M1[i].column = 0;
481  if (both)
482  {
483  M2[i].row = 0;
484  M2[i].column = 0;
485  }
486  }
487  }
488 
489 
490  template <typename number>
491  inline void
492  LocalResults<number>::initialize_quadrature(unsigned int np, unsigned int nv)
493  {
494  quadrature_data.reinit(np, nv);
495  }
496 
497 
498  template <typename number>
499  inline
500  unsigned int
502  {
503  return J.size();
504  }
505 
506 
507  template <typename number>
508  inline
509  unsigned int
511  {
512  return R.size();
513  }
514 
515 
516  template <typename number>
517  inline
518  unsigned int
520  {
521  return M1.size();
522  }
523 
524 
525  template <typename number>
526  inline
527  unsigned int
529  {
530  return quadrature_data.n_rows();
531  }
532 
533 
534  template <typename number>
535  inline
536  unsigned int
538  {
539  return quadrature_data.n_cols();
540  }
541 
542 
543  template <typename number>
544  inline
545  number &
547  {
548  AssertIndexRange(i,J.size());
549  return J[i];
550  }
551 
552 
553  template <typename number>
554  inline
557  {
558  AssertIndexRange(i,R.size());
559  return R[i];
560  }
561 
562 
563  template <typename number>
564  inline
566  LocalResults<number>::matrix(unsigned int i, bool external)
567  {
568  if (external)
569  {
570  AssertIndexRange(i,M2.size());
571  return M2[i];
572  }
573  AssertIndexRange(i,M1.size());
574  return M1[i];
575  }
576 
577 
578  template <typename number>
579  inline
580  number &
581  LocalResults<number>::quadrature_value(unsigned int k, unsigned int i)
582  {
583  return quadrature_data(k,i);
584  }
585 
586 
587  template <typename number>
588  inline
591  {
592  return quadrature_data;
593  }
594 
595 
596  template <typename number>
597  inline
598  number
599  LocalResults<number>::value(unsigned int i) const
600  {
601  AssertIndexRange(i,J.size());
602  return J[i];
603  }
604 
605 
606  template <typename number>
607  inline
608  const BlockVector<number> &
609  LocalResults<number>::vector(unsigned int i) const
610  {
611  AssertIndexRange(i,R.size());
612  return R[i];
613  }
614 
615 
616  template <typename number>
617  inline
619  LocalResults<number>::matrix(unsigned int i, bool external) const
620  {
621  if (external)
622  {
623  AssertIndexRange(i,M2.size());
624  return M2[i];
625  }
626  AssertIndexRange(i,M1.size());
627  return M1[i];
628  }
629 
630 
631  template <typename number>
632  inline
633  number
634  LocalResults<number>::quadrature_value(unsigned int k, unsigned int i) const
635  {
636  return quadrature_data(k,i);
637  }
638 
639 
640  template <typename number>
641  template <class StreamType>
642  void
643  LocalResults<number>::print_debug(StreamType &os) const
644  {
645  os << "J: " << J.size() << std::endl;
646  os << "R: " << R.size() << std::endl;
647  for (unsigned int i=0; i<R.size(); ++i)
648  {
649  os << " " << R[i].n_blocks() << " -";
650  for (unsigned int j=0; j<R[i].n_blocks(); ++j)
651  os << ' ' << R[i].block(j).size();
652  os << std::endl;
653  }
654  os << "M: " << M1.size() << " face " << M2.size() << std::endl;
655  for (unsigned int i=0; i<M1.size(); ++i)
656  {
657  os << " " << M1[i].row << "," << M1[i].column
658  << " " << M1[i].matrix.m() << 'x' << M1[i].matrix.n();
659  if (i < M2.size())
660  os << " face " << M2[i].row << "," << M2[i].column
661  << " " << M2[i].matrix.m() << 'x' << M2[i].matrix.n();
662  os << std::endl;
663  }
664  }
665 
666 }
667 
668 
669 DEAL_II_NAMESPACE_CLOSE
670 
671 #endif
unsigned int n_matrices() const
const value_type & block(size_type i) const
Definition: matrix_block.h:861
unsigned int n_quadrature_values() const
std::vector< BlockVector< number > > R
void initialize_quadrature(unsigned int np, unsigned int nv)
#define AssertIndexRange(index, range)
Definition: exceptions.h:1170
const value_type & block(size_type i) const
Definition: matrix_block.h:930
MatrixBlock< FullMatrix< number > > & matrix(unsigned int i, bool external=false)
void initialize_vectors(const unsigned int n)
size_type column
Definition: matrix_block.h:286
number & value(unsigned int i)
number & quadrature_value(unsigned int k, unsigned int i)
Table< 2, number > quadrature_data
unsigned int min_level() const
void initialize_matrices(unsigned int n, bool both)
unsigned int n_vectors() const
std::vector< number > J
std::vector< MatrixBlock< FullMatrix< number > > > M2
Table< 2, number > & quadrature_values()
unsigned int n_values() const
void initialize_local(MatrixBlock< FullMatrix< number > > &M, const unsigned int row, const unsigned int col)
MatrixType matrix
Definition: matrix_block.h:291
unsigned int size() const
Definition: matrix_block.h:898
void reinit(const BlockIndices &local_sizes)
Definition: mesh_worker.cc:27
std::size_t memory_consumption() const
Definition: mesh_worker.cc:45
std::vector< MatrixBlock< FullMatrix< number > > > M1
BlockVector< number > & vector(unsigned int i)
void initialize_numbers(const unsigned int n)
size_type row
Definition: matrix_block.h:281
unsigned int n_quadrature_points() const
unsigned int size() const
Number of stored data objects.
Definition: any_data.h:204