Reference documentation for deal.II version 9.0.0
simple.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2018 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  protected:
115 
120  };
121 
122 
154  template <typename MatrixType>
156  {
157  public:
162  MatrixSimple(double threshold = 1.e-12);
163 
167  void initialize(MatrixType &m);
168 
172  void initialize(std::vector<MatrixType> &m);
173 
182 
190  template <class DOFINFO>
191  void initialize_info(DOFINFO &info, bool face) const;
192 
197  template <class DOFINFO>
198  void assemble(const DOFINFO &info);
199 
204  template <class DOFINFO>
205  void assemble(const DOFINFO &info1,
206  const DOFINFO &info2);
207  protected:
211  std::vector<SmartPointer<MatrixType,MatrixSimple<MatrixType> > > matrix;
212 
218  const double threshold;
219 
220  private:
225  void assemble(const FullMatrix<double> &M,
226  const unsigned int index,
227  const std::vector<types::global_dof_index> &i1,
228  const std::vector<types::global_dof_index> &i2);
229 
234 
235  };
236 
237 
248  template <typename MatrixType>
250  {
251  public:
256  MGMatrixSimple(double threshold = 1.e-12);
257 
262 
267 
274 
288  template <class DOFINFO>
289  void initialize_info(DOFINFO &info, bool face) const;
290 
294  template <class DOFINFO>
295  void assemble(const DOFINFO &info);
296 
301  template <class DOFINFO>
302  void assemble(const DOFINFO &info1,
303  const DOFINFO &info2);
304  private:
308  void assemble(MatrixType &G,
309  const FullMatrix<double> &M,
310  const std::vector<types::global_dof_index> &i1,
311  const std::vector<types::global_dof_index> &i2);
312 
316  void assemble(MatrixType &G,
317  const FullMatrix<double> &M,
318  const std::vector<types::global_dof_index> &i1,
319  const std::vector<types::global_dof_index> &i2,
320  const unsigned int level);
321 
326  void assemble_up(MatrixType &G,
327  const FullMatrix<double> &M,
328  const std::vector<types::global_dof_index> &i1,
329  const std::vector<types::global_dof_index> &i2,
330  const unsigned int level = numbers::invalid_unsigned_int);
335  void assemble_down(MatrixType &G,
336  const FullMatrix<double> &M,
337  const std::vector<types::global_dof_index> &i1,
338  const std::vector<types::global_dof_index> &i2,
339  const unsigned int level = numbers::invalid_unsigned_int);
340 
345  void assemble_in(MatrixType &G,
346  const FullMatrix<double> &M,
347  const std::vector<types::global_dof_index> &i1,
348  const std::vector<types::global_dof_index> &i2,
349  const unsigned int level = numbers::invalid_unsigned_int);
350 
355  void assemble_out(MatrixType &G,
356  const FullMatrix<double> &M,
357  const std::vector<types::global_dof_index> &i1,
358  const std::vector<types::global_dof_index> &i2,
359  const unsigned int level = numbers::invalid_unsigned_int);
360 
365 
371 
377 
383 
393 
399  const double threshold;
400 
401  };
402 
403 
414  template <typename MatrixType, typename VectorType>
415  class SystemSimple :
416  private MatrixSimple<MatrixType>,
417  private ResidualSimple<VectorType>
418  {
419  public:
423  SystemSimple(double threshold = 1.e-12);
424 
428  void initialize(MatrixType &m, VectorType &rhs);
429 
438 
446  template <class DOFINFO>
447  void initialize_info(DOFINFO &info, bool face) const;
448 
452  template <class DOFINFO>
453  void assemble(const DOFINFO &info);
454 
459  template <class DOFINFO>
460  void assemble(const DOFINFO &info1,
461  const DOFINFO &info2);
462 
463  private:
468  void assemble(const FullMatrix<double> &M,
469  const Vector<double> &vector,
470  const unsigned int index,
471  const std::vector<types::global_dof_index> &indices);
472 
473  void assemble(const FullMatrix<double> &M,
474  const Vector<double> &vector,
475  const unsigned int index,
476  const std::vector<types::global_dof_index> &i1,
477  const std::vector<types::global_dof_index> &i2);
478  };
479 
480 
481 //----------------------------------------------------------------------//
482 
483  template <typename VectorType>
484  inline void
486  {
487  residuals = results;
488  }
489 
490  template <typename VectorType>
491  inline void
493  {
494  constraints = &c;
495  }
496 
497 
498  template <typename MatrixType>
499  inline void
501  {}
502 
503 
504  template <typename VectorType>
505  template <class DOFINFO>
506  inline void
508  {
509  info.initialize_vectors(residuals.size());
510  }
511 
512 
513  template <typename VectorType>
514  template <class DOFINFO>
515  inline void
517  {
518  for (unsigned int k=0; k<residuals.size(); ++k)
519  {
520  VectorType *v = residuals.entry<VectorType *>(k);
521  for (unsigned int i=0; i != info.vector(k).n_blocks(); ++i)
522  {
523  const std::vector<types::global_dof_index> &ldi = info.vector(k).n_blocks()==1?
524  info.indices:
525  info.indices_by_block[i];
526 
527  if (constraints != nullptr)
528  constraints->distribute_local_to_global(info.vector(k).block(i), ldi, *v);
529  else
530  v->add(ldi, info.vector(k).block(i));
531  }
532  }
533  }
534 
535  template <typename VectorType>
536  template <class DOFINFO>
537  inline void
539  const DOFINFO &info2)
540  {
541  assemble(info1);
542  assemble(info2);
543  }
544 
545 
546 //----------------------------------------------------------------------//
547 
548  template <typename MatrixType>
549  inline
551  :
552  threshold(threshold)
553  {}
554 
555 
556  template <typename MatrixType>
557  inline void
559  {
560  matrix.resize(1);
561  matrix[0] = &m;
562  }
563 
564 
565  template <typename MatrixType>
566  inline void
567  MatrixSimple<MatrixType>::initialize(std::vector<MatrixType> &m)
568  {
569  matrix.resize(m.size());
570  for (unsigned int i=0; i<m.size(); ++i)
571  matrix[i] = &m[i];
572  }
573 
574 
575  template <typename MatrixType>
576  inline void
578  {
579  constraints = &c;
580  }
581 
582 
583  template <typename MatrixType >
584  template <class DOFINFO>
585  inline void
586  MatrixSimple<MatrixType>::initialize_info(DOFINFO &info, bool face) const
587  {
588  Assert(matrix.size() != 0, ExcNotInitialized());
589 
590  const unsigned int n = info.indices_by_block.size();
591 
592  if (n == 0)
593  info.initialize_matrices(matrix.size(), face);
594  else
595  {
596  info.initialize_matrices(matrix.size()*n*n, face);
597  unsigned int k=0;
598  for (unsigned int m=0; m<matrix.size(); ++m)
599  for (unsigned int i=0; i<n; ++i)
600  for (unsigned int j=0; j<n; ++j,++k)
601  {
602  info.matrix(k,false).row = i;
603  info.matrix(k,false).column = j;
604  if (face)
605  {
606  info.matrix(k,true).row = i;
607  info.matrix(k,true).column = j;
608  }
609  }
610  }
611  }
612 
613 
614 
615  template <typename MatrixType>
616  inline void
618  const unsigned int index,
619  const std::vector<types::global_dof_index> &i1,
620  const std::vector<types::global_dof_index> &i2)
621  {
622  AssertDimension(M.m(), i1.size());
623  AssertDimension(M.n(), i2.size());
624 
625  if (constraints == nullptr)
626  {
627  for (unsigned int j=0; j<i1.size(); ++j)
628  for (unsigned int k=0; k<i2.size(); ++k)
629  if (std::fabs(M(j,k)) >= threshold)
630  matrix[index]->add(i1[j], i2[k], M(j,k));
631  }
632  else
633  constraints->distribute_local_to_global(M, i1, i2, *matrix[index]);
634  }
635 
636 
637  template <typename MatrixType>
638  template <class DOFINFO>
639  inline void
641  {
642  Assert(!info.level_cell, ExcMessage("Cell may not access level dofs"));
643  const unsigned int n = info.indices_by_block.size();
644 
645  if (n == 0)
646  for (unsigned int m=0; m<matrix.size(); ++m)
647  assemble(info.matrix(m,false).matrix, m, info.indices, info.indices);
648  else
649  {
650  for (unsigned int m=0; m<matrix.size(); ++m)
651  for (unsigned int k=0; k<n*n; ++k)
652  {
653  assemble(info.matrix(k+m*n*n,false).matrix, m,
654  info.indices_by_block[info.matrix(k+m*n*n,false).row],
655  info.indices_by_block[info.matrix(k+m*n*n,false).column]);
656  }
657  }
658  }
659 
660 
661  template <typename MatrixType>
662  template <class DOFINFO>
663  inline void
664  MatrixSimple<MatrixType>::assemble(const DOFINFO &info1, const DOFINFO &info2)
665  {
666  Assert(!info1.level_cell, ExcMessage("Cell may not access level dofs"));
667  Assert(!info2.level_cell, ExcMessage("Cell may not access level dofs"));
668  AssertDimension(info1.indices_by_block.size(),info2.indices_by_block.size());
669 
670  const unsigned int n = info1.indices_by_block.size();
671 
672  if (n == 0)
673  {
674  for (unsigned int m=0; m<matrix.size(); ++m)
675  {
676  assemble(info1.matrix(m,false).matrix, m, info1.indices, info1.indices);
677  assemble(info1.matrix(m,true).matrix, m, info1.indices, info2.indices);
678  assemble(info2.matrix(m,false).matrix, m, info2.indices, info2.indices);
679  assemble(info2.matrix(m,true).matrix, m, info2.indices, info1.indices);
680  }
681  }
682  else
683  {
684  for (unsigned int m=0; m<matrix.size(); ++m)
685  for (unsigned int k=0; k<n*n; ++k)
686  {
687  const unsigned int row = info1.matrix(k+m*n*n,false).row;
688  const unsigned int column = info1.matrix(k+m*n*n,false).column;
689 
690  assemble(info1.matrix(k+m*n*n,false).matrix, m,
691  info1.indices_by_block[row], info1.indices_by_block[column]);
692  assemble(info1.matrix(k+m*n*n,true).matrix, m,
693  info1.indices_by_block[row], info2.indices_by_block[column]);
694  assemble(info2.matrix(k+m*n*n,false).matrix, m,
695  info2.indices_by_block[row], info2.indices_by_block[column]);
696  assemble(info2.matrix(k+m*n*n,true).matrix, m,
697  info2.indices_by_block[row], info1.indices_by_block[column]);
698  }
699  }
700  }
701 
702 
703 //----------------------------------------------------------------------//
704 
705  template <typename MatrixType>
706  inline
708  :
709  threshold(threshold)
710  {}
711 
712 
713  template <typename MatrixType>
714  inline void
716  {
717  matrix = &m;
718  }
719 
720  template <typename MatrixType>
721  inline void
723  {
724  mg_constrained_dofs = &c;
725  }
726 
727 
728  template <typename MatrixType>
729  inline void
732  {
733  flux_up = &up;
734  flux_down = &down;
735  }
736 
737 
738  template <typename MatrixType>
739  inline void
742  {
743  interface_in = &in;
744  interface_out = &out;
745  }
746 
747 
748  template <typename MatrixType >
749  template <class DOFINFO>
750  inline void
751  MGMatrixSimple<MatrixType>::initialize_info(DOFINFO &info, bool face) const
752  {
753  const unsigned int n = info.indices_by_block.size();
754 
755  if (n == 0)
756  info.initialize_matrices(1, face);
757  else
758  {
759  info.initialize_matrices(n*n, face);
760  unsigned int k=0;
761  for (unsigned int i=0; i<n; ++i)
762  for (unsigned int j=0; j<n; ++j,++k)
763  {
764  info.matrix(k,false).row = i;
765  info.matrix(k,false).column = j;
766  if (face)
767  {
768  info.matrix(k,true).row = i;
769  info.matrix(k,true).column = j;
770  }
771  }
772  }
773  }
774 
775 
776  template <typename MatrixType>
777  inline void
779  (MatrixType &G,
780  const FullMatrix<double> &M,
781  const std::vector<types::global_dof_index> &i1,
782  const std::vector<types::global_dof_index> &i2)
783  {
784  AssertDimension(M.m(), i1.size());
785  AssertDimension(M.n(), i2.size());
786  Assert(mg_constrained_dofs == 0, ExcInternalError());
787 //TODO: Possibly remove this function all together
788 
789  for (unsigned int j=0; j<i1.size(); ++j)
790  for (unsigned int k=0; k<i2.size(); ++k)
791  if (std::fabs(M(j,k)) >= threshold)
792  G.add(i1[j], i2[k], M(j,k));
793  }
794 
795 
796  template <typename MatrixType>
797  inline void
799  (MatrixType &G,
800  const FullMatrix<double> &M,
801  const std::vector<types::global_dof_index> &i1,
802  const std::vector<types::global_dof_index> &i2,
803  const unsigned int level)
804  {
805  AssertDimension(M.m(), i1.size());
806  AssertDimension(M.n(), i2.size());
807 
808  if (mg_constrained_dofs == nullptr)
809  {
810  for (unsigned int j=0; j<i1.size(); ++j)
811  for (unsigned int k=0; k<i2.size(); ++k)
812  if (std::fabs(M(j,k)) >= threshold)
813  G.add(i1[j], i2[k], M(j,k));
814  }
815  else
816  {
817  for (unsigned int j=0; j<i1.size(); ++j)
818  for (unsigned int k=0; k<i2.size(); ++k)
819  {
820  // Only enter the local values into the global matrix,
821  // if the value is larger than the threshold
822  if (std::fabs(M(j,k)) < threshold)
823  continue;
824 
825  // Do not enter, if either the row or the column
826  // corresponds to an index on the refinement edge. The
827  // level problems are solved with homogeneous
828  // Dirichlet boundary conditions, therefore we
829  // eliminate these rows and columns. The corresponding
830  // matrix entries are entered by assemble_in() and
831  // assemble_out().
832  if (mg_constrained_dofs->at_refinement_edge(level, i1[j]) ||
833  mg_constrained_dofs->at_refinement_edge(level, i2[k]))
834  continue;
835 
836  // At the boundary, only enter the term on the
837  // diagonal, but not the coupling terms
838  if ((mg_constrained_dofs->is_boundary_index(level, i1[j]) ||
839  mg_constrained_dofs->is_boundary_index(level, i2[k])) &&
840  (i1[j] != i2[k]))
841  continue;
842 
843  G.add(i1[j], i2[k], M(j,k));
844  }
845  }
846  }
847 
848 
849  template <typename MatrixType>
850  inline void
852  (MatrixType &G,
853  const FullMatrix<double> &M,
854  const std::vector<types::global_dof_index> &i1,
855  const std::vector<types::global_dof_index> &i2,
856  const unsigned int level)
857  {
858  AssertDimension(M.n(), i1.size());
859  AssertDimension(M.m(), i2.size());
860 
861  if (mg_constrained_dofs == nullptr)
862  {
863  for (unsigned int j=0; j<i1.size(); ++j)
864  for (unsigned int k=0; k<i2.size(); ++k)
865  if (std::fabs(M(k,j)) >= threshold)
866  G.add(i1[j], i2[k], M(k,j));
867  }
868  else
869  {
870  for (unsigned int j=0; j<i1.size(); ++j)
871  for (unsigned int k=0; k<i2.size(); ++k)
872  if (std::fabs(M(k,j)) >= threshold)
873  if (!mg_constrained_dofs->at_refinement_edge(level, i2[k]))
874  G.add(i1[j], i2[k], M(k,j));
875  }
876  }
877 
878  template <typename MatrixType>
879  inline void
881  (MatrixType &G,
882  const FullMatrix<double> &M,
883  const std::vector<types::global_dof_index> &i1,
884  const std::vector<types::global_dof_index> &i2,
885  const unsigned int level)
886  {
887  AssertDimension(M.m(), i1.size());
888  AssertDimension(M.n(), i2.size());
889 
890  if (mg_constrained_dofs == nullptr)
891  {
892  for (unsigned int j=0; j<i1.size(); ++j)
893  for (unsigned int k=0; k<i2.size(); ++k)
894  if (std::fabs(M(j,k)) >= threshold)
895  G.add(i1[j], i2[k], M(j,k));
896  }
897  else
898  {
899  for (unsigned int j=0; j<i1.size(); ++j)
900  for (unsigned int k=0; k<i2.size(); ++k)
901  if (std::fabs(M(j,k)) >= threshold)
902  if (!mg_constrained_dofs->at_refinement_edge(level, i2[k]))
903  G.add(i1[j], i2[k], M(j,k));
904  }
905  }
906 
907  template <typename MatrixType>
908  inline void
910  (MatrixType &G,
911  const FullMatrix<double> &M,
912  const std::vector<types::global_dof_index> &i1,
913  const std::vector<types::global_dof_index> &i2,
914  const unsigned int level)
915  {
916  AssertDimension(M.m(), i1.size());
917  AssertDimension(M.n(), i2.size());
918  Assert(mg_constrained_dofs != nullptr, ExcInternalError());
919 
920  for (unsigned int j=0; j<i1.size(); ++j)
921  for (unsigned int k=0; k<i2.size(); ++k)
922  if (std::fabs(M(j,k)) >= threshold)
923  // Enter values into matrix only if j corresponds to a
924  // degree of freedom on the refinement edge, k does
925  // not, and both are not on the boundary. This is part
926  // the difference between the complete matrix with no
927  // boundary condition at the refinement edge and and
928  // the matrix assembled above by assemble().
929 
930  // Thus the logic is: enter the row if it is
931  // constrained by hanging node constraints (actually,
932  // the whole refinement edge), but not if it is
933  // constrained by a boundary constraint.
934  if (mg_constrained_dofs->at_refinement_edge(level, i1[j]) &&
935  !mg_constrained_dofs->at_refinement_edge(level, i2[k]))
936  {
937  if ((!mg_constrained_dofs->is_boundary_index(level, i1[j]) &&
938  !mg_constrained_dofs->is_boundary_index(level, i2[k]))
939  ||
940  (mg_constrained_dofs->is_boundary_index(level, i1[j]) &&
941  mg_constrained_dofs->is_boundary_index(level, i2[k]) &&
942  i1[j] == i2[k]))
943  G.add(i1[j], i2[k], M(j,k));
944  }
945  }
946 
947 
948  template <typename MatrixType>
949  inline void
951  (MatrixType &G,
952  const FullMatrix<double> &M,
953  const std::vector<types::global_dof_index> &i1,
954  const std::vector<types::global_dof_index> &i2,
955  const unsigned int level)
956  {
957  AssertDimension(M.n(), i1.size());
958  AssertDimension(M.m(), i2.size());
959  Assert(mg_constrained_dofs != nullptr, ExcInternalError());
960 
961  for (unsigned int j=0; j<i1.size(); ++j)
962  for (unsigned int k=0; k<i2.size(); ++k)
963  if (std::fabs(M(k,j)) >= threshold)
964  if (mg_constrained_dofs->at_refinement_edge(level, i1[j]) &&
965  !mg_constrained_dofs->at_refinement_edge(level, i2[k]))
966  {
967  if ((!mg_constrained_dofs->is_boundary_index(level, i1[j]) &&
968  !mg_constrained_dofs->is_boundary_index(level, i2[k]))
969  ||
970  (mg_constrained_dofs->is_boundary_index(level, i1[j]) &&
971  mg_constrained_dofs->is_boundary_index(level, i2[k]) &&
972  i1[j] == i2[k]))
973  G.add(i1[j], i2[k], M(k,j));
974  }
975  }
976 
977 
978  template <typename MatrixType>
979  template <class DOFINFO>
980  inline void
982  {
983  Assert(info.level_cell, ExcMessage("Cell must access level dofs"));
984  const unsigned int level = info.cell->level();
985 
986  if (info.indices_by_block.size() == 0)
987  {
988  assemble((*matrix)[level], info.matrix(0,false).matrix,
989  info.indices, info.indices, level);
990  if (mg_constrained_dofs != nullptr)
991  {
992  assemble_in((*interface_in)[level], info.matrix(0,false).matrix,
993  info.indices, info.indices, level);
994  assemble_out((*interface_out)[level],info.matrix(0,false).matrix,
995  info.indices, info.indices, level);
996  }
997  }
998  else
999  for (unsigned int k=0; k<info.n_matrices(); ++k)
1000  {
1001  const unsigned int row = info.matrix(k,false).row;
1002  const unsigned int column = info.matrix(k,false).column;
1003 
1004  assemble((*matrix)[level], info.matrix(k,false).matrix,
1005  info.indices_by_block[row], info.indices_by_block[column], level);
1006 
1007  if (mg_constrained_dofs != nullptr)
1008  {
1009  assemble_in((*interface_in)[level], info.matrix(k,false).matrix,
1010  info.indices_by_block[row], info.indices_by_block[column], level);
1011  assemble_out((*interface_out)[level],info.matrix(k,false).matrix,
1012  info.indices_by_block[column], info.indices_by_block[row], level);
1013  }
1014  }
1015  }
1016 
1017 
1018  template <typename MatrixType>
1019  template <class DOFINFO>
1020  inline void
1022  const DOFINFO &info2)
1023  {
1024  Assert(info1.level_cell, ExcMessage("Cell must access level dofs"));
1025  Assert(info2.level_cell, ExcMessage("Cell must access level dofs"));
1026  const unsigned int level1 = info1.cell->level();
1027  const unsigned int level2 = info2.cell->level();
1028 
1029  if (info1.indices_by_block.size() == 0)
1030  {
1031  if (level1 == level2)
1032  {
1033  assemble((*matrix)[level1], info1.matrix(0,false).matrix, info1.indices, info1.indices, level1);
1034  assemble((*matrix)[level1], info1.matrix(0,true).matrix, info1.indices, info2.indices, level1);
1035  assemble((*matrix)[level1], info2.matrix(0,false).matrix, info2.indices, info2.indices, level1);
1036  assemble((*matrix)[level1], info2.matrix(0,true).matrix, info2.indices, info1.indices, level1);
1037  }
1038  else
1039  {
1040  Assert(level1 > level2, ExcInternalError());
1041  // Do not add info2.M1,
1042  // which is done by
1043  // the coarser cell
1044  assemble((*matrix)[level1], info1.matrix(0,false).matrix, info1.indices, info1.indices, level1);
1045  if (level1>0)
1046  {
1047  assemble_up((*flux_up)[level1],info1.matrix(0,true).matrix, info2.indices, info1.indices, level1);
1048  assemble_down((*flux_down)[level1], info2.matrix(0,true).matrix, info2.indices, info1.indices, level1);
1049  }
1050  }
1051  }
1052  else
1053  for (unsigned int k=0; k<info1.n_matrices(); ++k)
1054  {
1055  const unsigned int row = info1.matrix(k,false).row;
1056  const unsigned int column = info1.matrix(k,false).column;
1057 
1058  if (level1 == level2)
1059  {
1060  assemble((*matrix)[level1], info1.matrix(k,false).matrix, info1.indices_by_block[row], info1.indices_by_block[column], level1);
1061  assemble((*matrix)[level1], info1.matrix(k,true).matrix, info1.indices_by_block[row], info2.indices_by_block[column], level1);
1062  assemble((*matrix)[level1], info2.matrix(k,false).matrix, info2.indices_by_block[row], info2.indices_by_block[column], level1);
1063  assemble((*matrix)[level1], info2.matrix(k,true).matrix, info2.indices_by_block[row], info1.indices_by_block[column], level1);
1064  }
1065  else
1066  {
1067  Assert(level1 > level2, ExcInternalError());
1068  // Do not add info2.M1,
1069  // which is done by
1070  // the coarser cell
1071  assemble((*matrix)[level1], info1.matrix(k,false).matrix, info1.indices_by_block[row], info1.indices_by_block[column], level1);
1072  if (level1>0)
1073  {
1074  assemble_up((*flux_up)[level1],info1.matrix(k,true).matrix, info2.indices_by_block[column], info1.indices_by_block[row], level1);
1075  assemble_down((*flux_down)[level1], info2.matrix(k,true).matrix, info2.indices_by_block[row], info1.indices_by_block[column], level1);
1076  }
1077  }
1078  }
1079  }
1080 
1081 //----------------------------------------------------------------------//
1082 
1083  template <typename MatrixType, typename VectorType>
1085  :
1086  MatrixSimple<MatrixType>(t)
1087  {}
1088 
1089 
1090  template <typename MatrixType, typename VectorType>
1091  inline void
1092  SystemSimple<MatrixType,VectorType>::initialize(MatrixType &m, VectorType &rhs)
1093  {
1094  AnyData data;
1095  VectorType *p = &rhs;
1096  data.add(p, "right hand side");
1097 
1100  }
1101 
1102  template <typename MatrixType, typename VectorType>
1103  inline void
1105  {
1107  }
1108 
1109 
1110  template <typename MatrixType, typename VectorType>
1111  template <class DOFINFO>
1112  inline void
1114  bool face) const
1115  {
1118  }
1119 
1120  template <typename MatrixType,typename VectorType>
1121  inline void
1123  const Vector<double> &vector,
1124  const unsigned int index,
1125  const std::vector<types::global_dof_index> &indices)
1126  {
1127  AssertDimension(M.m(), indices.size());
1128  AssertDimension(M.n(), indices.size());
1129 
1131  VectorType *v = residuals.entry<VectorType *>(index);
1132 
1134  {
1135  for (unsigned int i=0; i<indices.size(); ++i)
1136  (*v)(indices[i]) += vector(i);
1137 
1138  for (unsigned int j=0; j<indices.size(); ++j)
1139  for (unsigned int k=0; k<indices.size(); ++k)
1140  if (std::fabs(M(j,k)) >= MatrixSimple<MatrixType>::threshold)
1141  MatrixSimple<MatrixType>::matrix[index]->add(indices[j], indices[k], M(j,k));
1142  }
1143  else
1144  {
1145  ResidualSimple<VectorType>::constraints->distribute_local_to_global(M,vector,indices,*MatrixSimple<MatrixType>::matrix[index],*v, true);
1146  }
1147  }
1148 
1149  template <typename MatrixType,typename VectorType>
1150  inline void
1152  const Vector<double> &vector,
1153  const unsigned int index,
1154  const std::vector<types::global_dof_index> &i1,
1155  const std::vector<types::global_dof_index> &i2)
1156  {
1157  AssertDimension(M.m(), i1.size());
1158  AssertDimension(M.n(), i2.size());
1159 
1161  VectorType *v = residuals.entry<VectorType *>(index);
1162 
1164  {
1165  for (unsigned int j=0; j<i1.size(); ++j)
1166  for (unsigned int k=0; k<i2.size(); ++k)
1167  if (std::fabs(M(j,k)) >= MatrixSimple<MatrixType>::threshold)
1168  MatrixSimple<MatrixType>::matrix[index]->add(i1[j], i2[k], M(j,k));
1169  }
1170  else
1171  {
1172  ResidualSimple<VectorType>::constraints->distribute_local_to_global(vector, i1, i2, *v, M, false);
1173  ResidualSimple<VectorType>::constraints->distribute_local_to_global(M, i1, i2, *MatrixSimple<MatrixType>::matrix[index]);
1174  }
1175  }
1176 
1177 
1178  template <typename MatrixType, typename VectorType>
1179  template <class DOFINFO>
1180  inline void
1182  {
1184  Assert(!info.level_cell, ExcMessage("Cell may not access level dofs"));
1185  const unsigned int n = info.indices_by_block.size();
1186 
1187  if (n == 0)
1188  {
1189  for (unsigned int m=0; m<MatrixSimple<MatrixType>::matrix.size(); ++m)
1190  assemble(info.matrix(m,false).matrix,info.vector(m).block(0), m, info.indices);
1191  }
1192  else
1193  {
1194  for (unsigned int m=0; m<MatrixSimple<MatrixType>::matrix.size(); ++m)
1195  for (unsigned int k=0; k<n*n; ++k)
1196  {
1197  const unsigned int row = info.matrix(k+m*n*n,false).row;
1198  const unsigned int column = info.matrix(k+m*n*n,false).column;
1199 
1200  if (row == column)
1201  assemble(info.matrix(k+m*n*n,false).matrix,
1202  info.vector(m).block(row), m,
1203  info.indices_by_block[row]);
1204  else
1205  assemble(info.matrix(k+m*n*n,false).matrix,
1206  info.vector(m).block(row), m,
1207  info.indices_by_block[row],
1208  info.indices_by_block[column]);
1209  }
1210  }
1211 
1212  }
1213 
1214 
1215  template <typename MatrixType, typename VectorType>
1216  template <class DOFINFO>
1217  inline void
1219  const DOFINFO &info2)
1220  {
1221  Assert(!info1.level_cell, ExcMessage("Cell may not access level dofs"));
1222  Assert(!info2.level_cell, ExcMessage("Cell may not access level dofs"));
1223  AssertDimension(info1.indices_by_block.size(),info2.indices_by_block.size());
1224 
1225  const unsigned int n = info1.indices_by_block.size();
1226 
1227  if (n == 0)
1228  {
1229  for (unsigned int m=0; m<MatrixSimple<MatrixType>::matrix.size(); ++m)
1230  {
1231  assemble(info1.matrix(m,false).matrix, info1.vector(m).block(0), m, info1.indices);
1232  assemble(info1.matrix(m,true).matrix, info1.vector(m).block(0), m, info1.indices, info2.indices);
1233  assemble(info2.matrix(m,false).matrix, info2.vector(m).block(0), m, info2.indices);
1234  assemble(info2.matrix(m,true).matrix, info2.vector(m).block(0), m, info2.indices, info1.indices);
1235  }
1236  }
1237  else
1238  {
1239  for (unsigned int m=0; m<MatrixSimple<MatrixType>::matrix.size(); ++m)
1240  for (unsigned int k=0; k<n*n; ++k)
1241  {
1242  const unsigned int row = info1.matrix(k+m*n*n,false).row;
1243  const unsigned int column = info1.matrix(k+m*n*n,false).column;
1244 
1245  if (row == column)
1246  {
1247  assemble(info1.matrix(k+m*n*n,false).matrix, info1.vector(m).block(row), m,
1248  info1.indices_by_block[row]);
1249  assemble(info2.matrix(k+m*n*n,false).matrix, info2.vector(m).block(row), m,
1250  info2.indices_by_block[row]);
1251  }
1252  else
1253  {
1254  assemble(info1.matrix(k+m*n*n,false).matrix, info1.vector(m).block(row), m,
1255  info1.indices_by_block[row], info1.indices_by_block[column]);
1256  assemble(info2.matrix(k+m*n*n,false).matrix, info2.vector(m).block(row), m,
1257  info2.indices_by_block[row], info2.indices_by_block[column]);
1258  }
1259  assemble(info1.matrix(k+m*n*n,true).matrix, info1.vector(m).block(row), m,
1260  info1.indices_by_block[row], info2.indices_by_block[column]);
1261  assemble(info2.matrix(k+m*n*n,true).matrix, info2.vector(m).block(row), m,
1262  info2.indices_by_block[row], info1.indices_by_block[column]);
1263  }
1264  }
1265  }
1266  }
1267 }
1268 
1269 DEAL_II_NAMESPACE_CLOSE
1270 
1271 #endif
size_type m() const
static const unsigned int invalid_unsigned_int
Definition: types.h:173
void initialize_interfaces(MGLevelObject< MatrixType > &interface_in, MGLevelObject< MatrixType > &interface_out)
Definition: simple.h:741
SystemSimple(double threshold=1.e-12)
Definition: simple.h:1084
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1248
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:951
SmartPointer< const MGConstrainedDoFs, MGMatrixSimple< MatrixType > > mg_constrained_dofs
Definition: simple.h:392
void initialize(MGLevelObject< MatrixType > &m)
Definition: simple.h:715
static ::ExceptionBase & ExcNotInitialized()
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:852
void initialize_info(DOFINFO &info, bool face) const
Definition: simple.h:507
void initialize(MatrixType &m, VectorType &rhs)
Definition: simple.h:1092
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:910
static ::ExceptionBase & ExcMessage(std::string arg1)
size_type n() const
#define Assert(cond, exc)
Definition: exceptions.h:1142
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > interface_in
Definition: simple.h:382
void assemble(const DOFINFO &info)
Definition: simple.h:516
void initialize_fluxes(MGLevelObject< MatrixType > &flux_up, MGLevelObject< MatrixType > &flux_down)
Definition: simple.h:730
std::vector< SmartPointer< MatrixType, MatrixSimple< MatrixType > > > matrix
Definition: simple.h:211
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > flux_up
Definition: simple.h:370
void initialize(AnyData &results)
Definition: simple.h:485
MGMatrixSimple(double threshold=1.e-12)
Definition: simple.h:707
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > matrix
Definition: simple.h:364
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > interface_out
Definition: simple.h:388
void initialize_local_blocks(const BlockIndices &)
Definition: simple.h:500
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:881
SmartPointer< const ConstraintMatrix, ResidualSimple< VectorType > > constraints
Definition: simple.h:119
void add(type entry, const std::string &name)
Add a new data object.
Definition: any_data.h:427
void assemble(const DOFINFO &info)
Definition: simple.h:1181
void initialize_info(DOFINFO &info, bool face) const
Definition: simple.h:751
void initialize_info(DOFINFO &info, bool face) const
Definition: simple.h:1113
void assemble(const DOFINFO &info)
Definition: simple.h:981
SmartPointer< MGLevelObject< MatrixType >, MGMatrixSimple< MatrixType > > flux_down
Definition: simple.h:376
void initialize(MatrixType &m)
Definition: simple.h:558
MatrixSimple(double threshold=1.e-12)
Definition: simple.h:550
SmartPointer< const ConstraintMatrix, MatrixSimple< MatrixType > > constraints
Definition: simple.h:233
void assemble(const DOFINFO &info)
Definition: simple.h:640
void initialize_info(DOFINFO &info, bool face) const
Definition: simple.h:586
static ::ExceptionBase & ExcInternalError()