Reference documentation for deal.II version 8.4.2
simple.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2015 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_simple_h
18 #define dealii__mesh_worker_simple_h
19 
20 #include <deal.II/algorithms/any_data.h>
21 #include <deal.II/base/smartpointer.h>
22 #include <deal.II/base/mg_level_object.h>
23 #include <deal.II/lac/block_vector.h>
24 #include <deal.II/meshworker/dof_info.h>
25 #include <deal.II/meshworker/functional.h>
26 #include <deal.II/multigrid/mg_constrained_dofs.h>
27 
28 /*
29  * The header containing the classes MeshWorker::Assember::MatrixSimple,
30  * MeshWorker::Assember::MGMatrixSimple, MeshWorker::Assember::ResidualSimple,
31  * and MeshWorker::Assember::SystemSimple.
32  */
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 namespace MeshWorker
37 {
38  namespace Assembler
39  {
53  template <typename VectorType>
55  {
56  public:
63  void initialize(AnyData &results);
64 
69 
84 
92  template <class DOFINFO>
93  void initialize_info(DOFINFO &info, bool face) const;
94 
101  template <class DOFINFO>
102  void assemble(const DOFINFO &info);
103 
107  template<class DOFINFO>
108  void assemble(const DOFINFO &info1,
109  const DOFINFO &info2);
110  private:
119  };
120 
121 
153  template <typename MatrixType>
155  {
156  public:
161  MatrixSimple(double threshold = 1.e-12);
162 
166  void initialize(MatrixType &m);
167 
171  void initialize(std::vector<MatrixType> &m);
172 
181 
189  template <class DOFINFO>
190  void initialize_info(DOFINFO &info, bool face) const;
191 
196  template<class DOFINFO>
197  void assemble(const DOFINFO &info);
198 
203  template<class DOFINFO>
204  void assemble(const DOFINFO &info1,
205  const DOFINFO &info2);
206  private:
211  void assemble(const FullMatrix<double> &M,
212  const unsigned int index,
213  const std::vector<types::global_dof_index> &i1,
214  const std::vector<types::global_dof_index> &i2);
215 
219  std::vector<SmartPointer<MatrixType,MatrixSimple<MatrixType> > > matrix;
224 
230  const double threshold;
231 
232  };
233 
234 
245  template <typename MatrixType>
247  {
248  public:
253  MGMatrixSimple(double threshold = 1.e-12);
254 
259 
263  void initialize(const MGConstrainedDoFs &mg_constrained_dofs);
264 
269  void initialize_fluxes(MGLevelObject<MatrixType> &flux_up,
270  MGLevelObject<MatrixType> &flux_down);
271 
277  void initialize_interfaces(MGLevelObject<MatrixType> &interface_in,
278  MGLevelObject<MatrixType> &interface_out);
286  template <class DOFINFO>
287  void initialize_info(DOFINFO &info, bool face) const;
288 
292  template<class DOFINFO>
293  void assemble(const DOFINFO &info);
294 
299  template<class DOFINFO>
300  void assemble(const DOFINFO &info1,
301  const DOFINFO &info2);
302  private:
306  void assemble(MatrixType &G,
307  const FullMatrix<double> &M,
308  const std::vector<types::global_dof_index> &i1,
309  const std::vector<types::global_dof_index> &i2);
310 
314  void assemble(MatrixType &G,
315  const FullMatrix<double> &M,
316  const std::vector<types::global_dof_index> &i1,
317  const std::vector<types::global_dof_index> &i2,
318  const unsigned int level);
319 
324  void assemble_up(MatrixType &G,
325  const FullMatrix<double> &M,
326  const std::vector<types::global_dof_index> &i1,
327  const std::vector<types::global_dof_index> &i2,
328  const unsigned int level = numbers::invalid_unsigned_int);
333  void assemble_down(MatrixType &G,
334  const FullMatrix<double> &M,
335  const std::vector<types::global_dof_index> &i1,
336  const std::vector<types::global_dof_index> &i2,
337  const unsigned int level = numbers::invalid_unsigned_int);
338 
343  void assemble_in(MatrixType &G,
344  const FullMatrix<double> &M,
345  const std::vector<types::global_dof_index> &i1,
346  const std::vector<types::global_dof_index> &i2,
347  const unsigned int level = numbers::invalid_unsigned_int);
348 
353  void assemble_out(MatrixType &G,
354  const FullMatrix<double> &M,
355  const std::vector<types::global_dof_index> &i1,
356  const std::vector<types::global_dof_index> &i2,
357  const unsigned int level = numbers::invalid_unsigned_int);
358 
363 
369 
375 
381 
391 
397  const double threshold;
398 
399  };
400 
401 
412  template <typename MatrixType, typename VectorType>
413  class SystemSimple :
414  private MatrixSimple<MatrixType>,
415  private ResidualSimple<VectorType>
416  {
417  public:
421  SystemSimple(double threshold = 1.e-12);
422 
426  void initialize(MatrixType &m, VectorType &rhs);
427 
436 
444  template <class DOFINFO>
445  void initialize_info(DOFINFO &info, bool face) const;
446 
450  template<class DOFINFO>
451  void assemble(const DOFINFO &info);
452 
457  template<class DOFINFO>
458  void assemble(const DOFINFO &info1,
459  const DOFINFO &info2);
460  };
461 
462 
463 //----------------------------------------------------------------------//
464 
465  template <typename VectorType>
466  inline void
468  {
469  residuals = results;
470  }
471 
472  template <typename VectorType>
473  inline void
475  {
476  constraints = &c;
477  }
478 
479 
480  template <typename MatrixType>
481  inline void
483  {}
484 
485 
486  template <typename VectorType>
487  template <class DOFINFO>
488  inline void
490  {
491  info.initialize_vectors(residuals.size());
492  }
493 
494 
495  template <typename VectorType>
496  template <class DOFINFO>
497  inline void
499  {
500  for (unsigned int k=0; k<residuals.size(); ++k)
501  {
502  VectorType *v = residuals.entry<VectorType *>(k);
503  if (constraints == 0)
504  {
505  for (unsigned int i=0; i<info.vector(k).block(0).size(); ++i)
506  (*v)(info.indices[i]) += info.vector(k).block(0)(i);
507  }
508  else
509  {
510  if (info.indices_by_block.size() == 0)
511  constraints->distribute_local_to_global(info.vector(k).block(0), info.indices, *v);
512  else
513  for (unsigned int i=0; i != info.vector(k).n_blocks(); ++i)
514  constraints->distribute_local_to_global(info.vector(k).block(i), info.indices_by_block[i], *v);
515  }
516  }
517  }
518 
519  template <typename VectorType>
520  template <class DOFINFO>
521  inline void
523  const DOFINFO &info2)
524  {
525  for (unsigned int k=0; k<residuals.size(); ++k)
526  {
527  VectorType *v = residuals.entry<VectorType *>(k);
528  if (constraints == 0)
529  {
530  for (unsigned int i=0; i<info1.vector(k).block(0).size(); ++i)
531  (*v)(info1.indices[i]) += info1.vector(k).block(0)(i);
532  for (unsigned int i=0; i<info2.vector(k).block(0).size(); ++i)
533  (*v)(info2.indices[i]) += info2.vector(k).block(0)(i);
534  }
535  else
536  {
537  if (info1.indices_by_block.size() == 0 && info2.indices_by_block.size() == 0)
538  {
539  constraints->distribute_local_to_global
540  (info1.vector(k).block(0), info1.indices, *v);
541  constraints->distribute_local_to_global
542  (info2.vector(k).block(0), info2.indices, *v);
543  }
544  else if (info1.indices_by_block.size() != 0 && info2.indices_by_block.size() != 0)
545  {
546  for (unsigned int i=0; i<info1.vector(k).n_blocks(); ++i)
547  {
548  constraints->distribute_local_to_global
549  (info1.vector(k).block(i), info1.indices_by_block[i], *v);
550  constraints->distribute_local_to_global
551  (info2.vector(k).block(i), info2.indices_by_block[i], *v);
552  }
553  }
554  else
555  {
556  Assert(false, ExcNotImplemented());
557  }
558  }
559  }
560  }
561 
562 
563 //----------------------------------------------------------------------//
564 
565  template <typename MatrixType>
566  inline
568  :
569  threshold(threshold)
570  {}
571 
572 
573  template <typename MatrixType>
574  inline void
576  {
577  matrix.resize(1);
578  matrix[0] = &m;
579  }
580 
581 
582  template <typename MatrixType>
583  inline void
584  MatrixSimple<MatrixType>::initialize(std::vector<MatrixType> &m)
585  {
586  matrix.resize(m.size());
587  for (unsigned int i=0; i<m.size(); ++i)
588  matrix[i] = &m[i];
589  }
590 
591 
592  template <typename MatrixType>
593  inline void
595  {
596  constraints = &c;
597  }
598 
599 
600  template <typename MatrixType >
601  template <class DOFINFO>
602  inline void
603  MatrixSimple<MatrixType>::initialize_info(DOFINFO &info, bool face) const
604  {
605  Assert(matrix.size() != 0, ExcNotInitialized());
606 
607  const unsigned int n = info.indices_by_block.size();
608 
609  if (n == 0)
610  info.initialize_matrices(matrix.size(), face);
611  else
612  {
613  info.initialize_matrices(matrix.size()*n*n, face);
614  unsigned int k=0;
615  for (unsigned int m=0; m<matrix.size(); ++m)
616  for (unsigned int i=0; i<n; ++i)
617  for (unsigned int j=0; j<n; ++j,++k)
618  {
619  info.matrix(k,false).row = i;
620  info.matrix(k,false).column = j;
621  if (face)
622  {
623  info.matrix(k,true).row = i;
624  info.matrix(k,true).column = j;
625  }
626  }
627  }
628  }
629 
630 
631 
632  template <typename MatrixType>
633  inline void
635  const unsigned int index,
636  const std::vector<types::global_dof_index> &i1,
637  const std::vector<types::global_dof_index> &i2)
638  {
639  AssertDimension(M.m(), i1.size());
640  AssertDimension(M.n(), i2.size());
641 
642  if (constraints == 0)
643  {
644  for (unsigned int j=0; j<i1.size(); ++j)
645  for (unsigned int k=0; k<i2.size(); ++k)
646  if (std::fabs(M(j,k)) >= threshold)
647  matrix[index]->add(i1[j], i2[k], M(j,k));
648  }
649  else
650  constraints->distribute_local_to_global(M, i1, i2, *matrix[index]);
651  }
652 
653 
654  template <typename MatrixType>
655  template <class DOFINFO>
656  inline void
658  {
659  Assert(!info.level_cell, ExcMessage("Cell may not access level dofs"));
660  const unsigned int n = info.indices_by_block.size();
661 
662  if (n == 0)
663  for (unsigned int m=0; m<matrix.size(); ++m)
664  assemble(info.matrix(m,false).matrix, m, info.indices, info.indices);
665  else
666  {
667  for (unsigned int m=0; m<matrix.size(); ++m)
668  for (unsigned int k=0; k<n*n; ++k)
669  {
670  assemble(info.matrix(k+m*n*n,false).matrix, m,
671  info.indices_by_block[info.matrix(k+m*n*n,false).row],
672  info.indices_by_block[info.matrix(k+m*n*n,false).column]);
673  }
674  }
675  }
676 
677 
678  template <typename MatrixType>
679  template <class DOFINFO>
680  inline void
681  MatrixSimple<MatrixType>::assemble(const DOFINFO &info1, const DOFINFO &info2)
682  {
683  Assert(!info1.level_cell, ExcMessage("Cell may not access level dofs"));
684  Assert(!info2.level_cell, ExcMessage("Cell may not access level dofs"));
685  AssertDimension(info1.indices_by_block.size(),info2.indices_by_block.size());
686 
687  const unsigned int n = info1.indices_by_block.size();
688 
689  if (n == 0)
690  {
691  for (unsigned int m=0; m<matrix.size(); ++m)
692  {
693  assemble(info1.matrix(m,false).matrix, m, info1.indices, info1.indices);
694  assemble(info1.matrix(m,true).matrix, m, info1.indices, info2.indices);
695  assemble(info2.matrix(m,false).matrix, m, info2.indices, info2.indices);
696  assemble(info2.matrix(m,true).matrix, m, info2.indices, info1.indices);
697  }
698  }
699  else
700  {
701  for (unsigned int m=0; m<matrix.size(); ++m)
702  for (unsigned int k=0; k<n*n; ++k)
703  {
704  const unsigned int row = info1.matrix(k+m*n*n,false).row;
705  const unsigned int column = info1.matrix(k+m*n*n,false).column;
706 
707  assemble(info1.matrix(k+m*n*n,false).matrix, m,
708  info1.indices_by_block[row], info1.indices_by_block[column]);
709  assemble(info1.matrix(k+m*n*n,true).matrix, m,
710  info1.indices_by_block[row], info2.indices_by_block[column]);
711  assemble(info2.matrix(k+m*n*n,false).matrix, m,
712  info2.indices_by_block[row], info2.indices_by_block[column]);
713  assemble(info2.matrix(k+m*n*n,true).matrix, m,
714  info2.indices_by_block[row], info1.indices_by_block[column]);
715  }
716  }
717  }
718 
719 
720 //----------------------------------------------------------------------//
721 
722  template <typename MatrixType>
723  inline
725  :
726  threshold(threshold)
727  {}
728 
729 
730  template <typename MatrixType>
731  inline void
733  {
734  matrix = &m;
735  }
736 
737  template <typename MatrixType>
738  inline void
740  {
741  mg_constrained_dofs = &c;
742  }
743 
744 
745  template <typename MatrixType>
746  inline void
749  {
750  flux_up = &up;
751  flux_down = &down;
752  }
753 
754 
755  template <typename MatrixType>
756  inline void
759  {
760  interface_in = &in;
761  interface_out = &out;
762  }
763 
764 
765  template <typename MatrixType >
766  template <class DOFINFO>
767  inline void
768  MGMatrixSimple<MatrixType>::initialize_info(DOFINFO &info, bool face) const
769  {
770  const unsigned int n = info.indices_by_block.size();
771 
772  if (n == 0)
773  info.initialize_matrices(1, face);
774  else
775  {
776  info.initialize_matrices(n*n, face);
777  unsigned int k=0;
778  for (unsigned int i=0; i<n; ++i)
779  for (unsigned int j=0; j<n; ++j,++k)
780  {
781  info.matrix(k,false).row = i;
782  info.matrix(k,false).column = j;
783  if (face)
784  {
785  info.matrix(k,true).row = i;
786  info.matrix(k,true).column = j;
787  }
788  }
789  }
790  }
791 
792 
793  template <typename MatrixType>
794  inline void
796  (MatrixType &G,
797  const FullMatrix<double> &M,
798  const std::vector<types::global_dof_index> &i1,
799  const std::vector<types::global_dof_index> &i2)
800  {
801  AssertDimension(M.m(), i1.size());
802  AssertDimension(M.n(), i2.size());
803  Assert(mg_constrained_dofs == 0, ExcInternalError());
804 //TODO: Possibly remove this function all together
805 
806  for (unsigned int j=0; j<i1.size(); ++j)
807  for (unsigned int k=0; k<i2.size(); ++k)
808  if (std::fabs(M(j,k)) >= threshold)
809  G.add(i1[j], i2[k], M(j,k));
810  }
811 
812 
813  template <typename MatrixType>
814  inline void
816  (MatrixType &G,
817  const FullMatrix<double> &M,
818  const std::vector<types::global_dof_index> &i1,
819  const std::vector<types::global_dof_index> &i2,
820  const unsigned int level)
821  {
822  AssertDimension(M.m(), i1.size());
823  AssertDimension(M.n(), i2.size());
824 
825  if (mg_constrained_dofs == 0)
826  {
827  for (unsigned int j=0; j<i1.size(); ++j)
828  for (unsigned int k=0; k<i2.size(); ++k)
829  if (std::fabs(M(j,k)) >= threshold)
830  G.add(i1[j], i2[k], M(j,k));
831  }
832  else
833  {
834  for (unsigned int j=0; j<i1.size(); ++j)
835  for (unsigned int k=0; k<i2.size(); ++k)
836  {
837  // Only enter the local values into the global matrix,
838  // if the value is larger than the threshold
839  if (std::fabs(M(j,k)) < threshold)
840  continue;
841 
842  // Do not enter, if either the row or the column
843  // corresponds to an index on the refinement edge. The
844  // level problems are solved with homogeneous
845  // Dirichlet boundary conditions, therefore we
846  // eliminate these rows and columns. The corresponding
847  // matrix entries are entered by assemble_in() and
848  // assemble_out().
849  if (mg_constrained_dofs->at_refinement_edge(level, i1[j]) ||
850  mg_constrained_dofs->at_refinement_edge(level, i2[k]))
851  continue;
852 
853  // At the boundary, only enter the term on the
854  // diagonal, but not the coupling terms
855  if ((mg_constrained_dofs->is_boundary_index(level, i1[j]) ||
856  mg_constrained_dofs->is_boundary_index(level, i2[k])) &&
857  (i1[j] != i2[k]))
858  continue;
859 
860  G.add(i1[j], i2[k], M(j,k));
861  }
862  }
863  }
864 
865 
866  template <typename MatrixType>
867  inline void
869  (MatrixType &G,
870  const FullMatrix<double> &M,
871  const std::vector<types::global_dof_index> &i1,
872  const std::vector<types::global_dof_index> &i2,
873  const unsigned int level)
874  {
875  AssertDimension(M.n(), i1.size());
876  AssertDimension(M.m(), i2.size());
877 
878  if (mg_constrained_dofs == 0)
879  {
880  for (unsigned int j=0; j<i1.size(); ++j)
881  for (unsigned int k=0; k<i2.size(); ++k)
882  if (std::fabs(M(k,j)) >= threshold)
883  G.add(i1[j], i2[k], M(k,j));
884  }
885  else
886  {
887  for (unsigned int j=0; j<i1.size(); ++j)
888  for (unsigned int k=0; k<i2.size(); ++k)
889  if (std::fabs(M(k,j)) >= threshold)
890  if (!mg_constrained_dofs->at_refinement_edge(level, i2[k]))
891  G.add(i1[j], i2[k], M(k,j));
892  }
893  }
894 
895  template <typename MatrixType>
896  inline void
898  (MatrixType &G,
899  const FullMatrix<double> &M,
900  const std::vector<types::global_dof_index> &i1,
901  const std::vector<types::global_dof_index> &i2,
902  const unsigned int level)
903  {
904  AssertDimension(M.m(), i1.size());
905  AssertDimension(M.n(), i2.size());
906 
907  if (mg_constrained_dofs == 0)
908  {
909  for (unsigned int j=0; j<i1.size(); ++j)
910  for (unsigned int k=0; k<i2.size(); ++k)
911  if (std::fabs(M(j,k)) >= threshold)
912  G.add(i1[j], i2[k], M(j,k));
913  }
914  else
915  {
916  for (unsigned int j=0; j<i1.size(); ++j)
917  for (unsigned int k=0; k<i2.size(); ++k)
918  if (std::fabs(M(j,k)) >= threshold)
919  if (!mg_constrained_dofs->at_refinement_edge(level, i2[k]))
920  G.add(i1[j], i2[k], M(j,k));
921  }
922  }
923 
924  template <typename MatrixType>
925  inline void
927  (MatrixType &G,
928  const FullMatrix<double> &M,
929  const std::vector<types::global_dof_index> &i1,
930  const std::vector<types::global_dof_index> &i2,
931  const unsigned int level)
932  {
933  AssertDimension(M.m(), i1.size());
934  AssertDimension(M.n(), i2.size());
935  Assert(mg_constrained_dofs != 0, ExcInternalError());
936 
937  for (unsigned int j=0; j<i1.size(); ++j)
938  for (unsigned int k=0; k<i2.size(); ++k)
939  if (std::fabs(M(j,k)) >= threshold)
940  // Enter values into matrix only if j corresponds to a
941  // degree of freedom on the refinemenent edge, k does
942  // not, and both are not on the boundary. This is part
943  // the difference between the complete matrix with no
944  // boundary condition at the refinement edge and and
945  // the matrix assembled above by assemble().
946 
947  // Thus the logic is: enter the row if it is
948  // constrained by hanging node constraints (actually,
949  // the whole refinement edge), but not if it is
950  // constrained by a boundary constraint.
951  if (mg_constrained_dofs->at_refinement_edge(level, i1[j]) &&
952  !mg_constrained_dofs->at_refinement_edge(level, i2[k]))
953  {
954  if ((!mg_constrained_dofs->is_boundary_index(level, i1[j]) &&
955  !mg_constrained_dofs->is_boundary_index(level, i2[k]))
956  ||
957  (mg_constrained_dofs->is_boundary_index(level, i1[j]) &&
958  mg_constrained_dofs->is_boundary_index(level, i2[k]) &&
959  i1[j] == i2[k]))
960  G.add(i1[j], i2[k], M(j,k));
961  }
962  }
963 
964 
965  template <typename MatrixType>
966  inline void
968  (MatrixType &G,
969  const FullMatrix<double> &M,
970  const std::vector<types::global_dof_index> &i1,
971  const std::vector<types::global_dof_index> &i2,
972  const unsigned int level)
973  {
974  AssertDimension(M.n(), i1.size());
975  AssertDimension(M.m(), i2.size());
976  Assert(mg_constrained_dofs != 0, ExcInternalError());
977 
978  for (unsigned int j=0; j<i1.size(); ++j)
979  for (unsigned int k=0; k<i2.size(); ++k)
980  if (std::fabs(M(k,j)) >= threshold)
981  if (mg_constrained_dofs->at_refinement_edge(level, i1[j]) &&
982  !mg_constrained_dofs->at_refinement_edge(level, i2[k]))
983  {
984  if ((!mg_constrained_dofs->is_boundary_index(level, i1[j]) &&
985  !mg_constrained_dofs->is_boundary_index(level, i2[k]))
986  ||
987  (mg_constrained_dofs->is_boundary_index(level, i1[j]) &&
988  mg_constrained_dofs->is_boundary_index(level, i2[k]) &&
989  i1[j] == i2[k]))
990  G.add(i1[j], i2[k], M(k,j));
991  }
992  }
993 
994 
995  template <typename MatrixType>
996  template <class DOFINFO>
997  inline void
999  {
1000  Assert(info.level_cell, ExcMessage("Cell must access level dofs"));
1001  const unsigned int level = info.cell->level();
1002 
1003  if (info.indices_by_block.size() == 0)
1004  {
1005  assemble((*matrix)[level], info.matrix(0,false).matrix,
1006  info.indices, info.indices, level);
1007  if (mg_constrained_dofs != 0)
1008  {
1009  assemble_in((*interface_in)[level], info.matrix(0,false).matrix,
1010  info.indices, info.indices, level);
1011  assemble_out((*interface_out)[level],info.matrix(0,false).matrix,
1012  info.indices, info.indices, level);
1013  }
1014  }
1015  else
1016  for (unsigned int k=0; k<info.n_matrices(); ++k)
1017  {
1018  const unsigned int row = info.matrix(k,false).row;
1019  const unsigned int column = info.matrix(k,false).column;
1020 
1021  assemble((*matrix)[level], info.matrix(k,false).matrix,
1022  info.indices_by_block[row], info.indices_by_block[column], level);
1023 
1024  if (mg_constrained_dofs != 0)
1025  {
1026  assemble_in((*interface_in)[level], info.matrix(k,false).matrix,
1027  info.indices_by_block[row], info.indices_by_block[column], level);
1028  assemble_out((*interface_out)[level],info.matrix(k,false).matrix,
1029  info.indices_by_block[column], info.indices_by_block[row], level);
1030  }
1031  }
1032  }
1033 
1034 
1035  template <typename MatrixType>
1036  template <class DOFINFO>
1037  inline void
1039  const DOFINFO &info2)
1040  {
1041  Assert(info1.level_cell, ExcMessage("Cell must access level dofs"));
1042  Assert(info2.level_cell, ExcMessage("Cell must access level dofs"));
1043  const unsigned int level1 = info1.cell->level();
1044  const unsigned int level2 = info2.cell->level();
1045 
1046  if (info1.indices_by_block.size() == 0)
1047  {
1048  if (level1 == level2)
1049  {
1050  assemble((*matrix)[level1], info1.matrix(0,false).matrix, info1.indices, info1.indices, level1);
1051  assemble((*matrix)[level1], info1.matrix(0,true).matrix, info1.indices, info2.indices, level1);
1052  assemble((*matrix)[level1], info2.matrix(0,false).matrix, info2.indices, info2.indices, level1);
1053  assemble((*matrix)[level1], info2.matrix(0,true).matrix, info2.indices, info1.indices, level1);
1054  }
1055  else
1056  {
1057  Assert(level1 > level2, ExcInternalError());
1058  // Do not add info2.M1,
1059  // which is done by
1060  // the coarser cell
1061  assemble((*matrix)[level1], info1.matrix(0,false).matrix, info1.indices, info1.indices, level1);
1062  if (level1>0)
1063  {
1064  assemble_up((*flux_up)[level1],info1.matrix(0,true).matrix, info2.indices, info1.indices, level1);
1065  assemble_down((*flux_down)[level1], info2.matrix(0,true).matrix, info2.indices, info1.indices, level1);
1066  }
1067  }
1068  }
1069  else
1070  for (unsigned int k=0; k<info1.n_matrices(); ++k)
1071  {
1072  const unsigned int row = info1.matrix(k,false).row;
1073  const unsigned int column = info1.matrix(k,false).column;
1074 
1075  if (level1 == level2)
1076  {
1077  assemble((*matrix)[level1], info1.matrix(k,false).matrix, info1.indices_by_block[row], info1.indices_by_block[column], level1);
1078  assemble((*matrix)[level1], info1.matrix(k,true).matrix, info1.indices_by_block[row], info2.indices_by_block[column], level1);
1079  assemble((*matrix)[level1], info2.matrix(k,false).matrix, info2.indices_by_block[row], info2.indices_by_block[column], level1);
1080  assemble((*matrix)[level1], info2.matrix(k,true).matrix, info2.indices_by_block[row], info1.indices_by_block[column], level1);
1081  }
1082  else
1083  {
1084  Assert(level1 > level2, ExcInternalError());
1085  // Do not add info2.M1,
1086  // which is done by
1087  // the coarser cell
1088  assemble((*matrix)[level1], info1.matrix(k,false).matrix, info1.indices_by_block[row], info1.indices_by_block[column], level1);
1089  if (level1>0)
1090  {
1091  assemble_up((*flux_up)[level1],info1.matrix(k,true).matrix, info2.indices_by_block[row], info1.indices_by_block[column], level1);
1092  assemble_down((*flux_down)[level1], info2.matrix(k,true).matrix, info2.indices_by_block[row], info1.indices_by_block[column], level1);
1093  }
1094  }
1095  }
1096  }
1097 
1098 //----------------------------------------------------------------------//
1099 
1100  template <typename MatrixType, typename VectorType>
1102  :
1103  MatrixSimple<MatrixType>(t)
1104  {}
1105 
1106 
1107  template <typename MatrixType, typename VectorType>
1108  inline void
1109  SystemSimple<MatrixType,VectorType>::initialize(MatrixType &m, VectorType &rhs)
1110  {
1111  AnyData data;
1112  VectorType *p = &rhs;
1113  data.add(p, "right hand side");
1114 
1117  }
1118 
1119  template <typename MatrixType, typename VectorType>
1120  inline void
1122  {
1125  }
1126 
1127 
1128  template <typename MatrixType, typename VectorType>
1129  template <class DOFINFO>
1130  inline void
1132  bool face) const
1133  {
1136  }
1137 
1138 
1139  template <typename MatrixType, typename VectorType>
1140  template <class DOFINFO>
1141  inline void
1143  {
1146  }
1147 
1148 
1149  template <typename MatrixType, typename VectorType>
1150  template <class DOFINFO>
1151  inline void
1153  const DOFINFO &info2)
1154  {
1155  MatrixSimple<MatrixType>::assemble(info1, info2);
1157  }
1158  }
1159 }
1160 
1161 DEAL_II_NAMESPACE_CLOSE
1162 
1163 #endif
size_type m() const
static const unsigned int invalid_unsigned_int
Definition: types.h:164
type entry(const std::string &name)
Access to stored data object by name.
Definition: any_data.h:341
void initialize_interfaces(MGLevelObject< MatrixType > &interface_in, MGLevelObject< MatrixType > &interface_out)
Definition: simple.h:758
SystemSimple(double threshold=1.e-12)
Definition: simple.h:1101
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1052
void assemble_out(MatrixType &G, const FullMatrix< double > &M, const std::vector< types::global_dof_index > &i1, const std::vector< types::global_dof_index > &i2, const unsigned int level=numbers::invalid_unsigned_int)
Definition: simple.h:968
::ExceptionBase & ExcMessage(std::string arg1)
SmartPointer< const MGConstrainedDoFs, MGMatrixSimple< MatrixType > > mg_constrained_dofs
Definition: simple.h:390
Auxiliary class aiding in the handling of block structures like in BlockVector or FESystem...
Definition: block_indices.h:54
void initialize(MGLevelObject< MatrixType > &m)
Definition: simple.h:732
void assemble_up(MatrixType &G, const FullMatrix< double > &M, const std::vector< types::global_dof_index > &i1, const std::vector< types::global_dof_index > &i2, const unsigned int level=numbers::invalid_unsigned_int)
Definition: simple.h:869
void initialize_info(DOFINFO &info, bool face) const
Definition: simple.h:489
void initialize(MatrixType &m, VectorType &rhs)
Definition: simple.h:1109
void assemble_in(MatrixType &G, const FullMatrix< double > &M, const std::vector< types::global_dof_index > &i1, const std::vector< types::global_dof_index > &i2, const unsigned int level=numbers::invalid_unsigned_int)
Definition: simple.h:927
size_type n() const
#define Assert(cond, exc)
Definition: exceptions.h:294
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > interface_in
Definition: simple.h:380
void assemble(const DOFINFO &info)
Definition: simple.h:498
void initialize_fluxes(MGLevelObject< MatrixType > &flux_up, MGLevelObject< MatrixType > &flux_down)
Definition: simple.h:747
std::vector< SmartPointer< MatrixType, MatrixSimple< MatrixType > > > matrix
Definition: simple.h:219
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > flux_up
Definition: simple.h:368
void initialize(AnyData &results)
Definition: simple.h:467
MGMatrixSimple(double threshold=1.e-12)
Definition: simple.h:724
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > matrix
Definition: simple.h:362
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > interface_out
Definition: simple.h:386
void initialize_local_blocks(const BlockIndices &)
Definition: simple.h:482
void assemble_down(MatrixType &G, const FullMatrix< double > &M, const std::vector< types::global_dof_index > &i1, const std::vector< types::global_dof_index > &i2, const unsigned int level=numbers::invalid_unsigned_int)
Definition: simple.h:898
SmartPointer< const ConstraintMatrix, ResidualSimple< VectorType > > constraints
Definition: simple.h:118
void add(type entry, const std::string &name)
Add a new data object.
Definition: any_data.h:430
unsigned int size() const
Number of stored data objects.
Definition: any_data.h:204
void assemble(const DOFINFO &info)
Definition: simple.h:1142
void initialize_info(DOFINFO &info, bool face) const
Definition: simple.h:768
void initialize_info(DOFINFO &info, bool face) const
Definition: simple.h:1131
void assemble(const DOFINFO &info)
Definition: simple.h:998
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > flux_down
Definition: simple.h:374
void initialize(MatrixType &m)
Definition: simple.h:575
MatrixSimple(double threshold=1.e-12)
Definition: simple.h:567
SmartPointer< const ConstraintMatrix, MatrixSimple< MatrixType > > constraints
Definition: simple.h:223
void assemble(const DOFINFO &info)
Definition: simple.h:657
void initialize_info(DOFINFO &info, bool face) const
Definition: simple.h:603