Reference documentation for deal.II version 8.5.1
assembler.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 2017 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_assembler_h
18 #define dealii__mesh_worker_assembler_h
19 
20 #include <deal.II/base/smartpointer.h>
21 #include <deal.II/base/mg_level_object.h>
22 #include <deal.II/lac/block_vector.h>
23 #include <deal.II/meshworker/dof_info.h>
24 #include <deal.II/meshworker/functional.h>
25 #include <deal.II/meshworker/simple.h>
26 #include <deal.II/multigrid/mg_constrained_dofs.h>
27 #include <deal.II/algorithms/any_data.h>
28 
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 namespace MeshWorker
33 {
83  namespace Assembler
84  {
106  template <typename VectorType>
108  {
109  public:
113  void initialize(const BlockInfo *block_info,
114  AnyData &residuals);
126  template <class DOFINFO>
127  void initialize_info(DOFINFO &info, bool face) const;
128 
129 
133  template<class DOFINFO>
134  void assemble(const DOFINFO &info);
135 
139  template<class DOFINFO>
140  void assemble(const DOFINFO &info1,
141  const DOFINFO &info2);
142  private:
146  void assemble(VectorType &global,
147  const BlockVector<double> &local,
148  const std::vector<types::global_dof_index> &dof);
149 
154 
158  SmartPointer<const BlockInfo,
165  };
166 
167 
195  template <typename MatrixType, typename number = double>
197  {
198  public:
204 
209  void initialize(const BlockInfo *block_info,
211 
223  template <class DOFINFO>
224  void initialize_info(DOFINFO &info, bool face) const;
225 
226 
230  template<class DOFINFO>
231  void assemble(const DOFINFO &info);
232 
236  template<class DOFINFO>
237  void assemble(const DOFINFO &info1,
238  const DOFINFO &info2);
239 
240  private:
244  void assemble(MatrixBlock<MatrixType> &global,
245  const FullMatrix<number> &local,
246  const unsigned int block_row,
247  const unsigned int block_col,
248  const std::vector<types::global_dof_index> &dof1,
249  const std::vector<types::global_dof_index> &dof2);
250 
256 
260  SmartPointer<const BlockInfo,
267 
273  const double threshold;
274 
275  };
276 
301  template <typename MatrixType, typename number = double>
303  {
304  public:
308 
314 
319  void initialize(const BlockInfo *block_info,
321 
326 
333 
347  template <class DOFINFO>
348  void initialize_info(DOFINFO &info, bool face) const;
349 
350 
354  template<class DOFINFO>
355  void assemble(const DOFINFO &info);
356 
360  template<class DOFINFO>
361  void assemble(const DOFINFO &info1,
362  const DOFINFO &info2);
363 
364  private:
368  void assemble(MatrixType &global,
369  const FullMatrix<number> &local,
370  const unsigned int block_row,
371  const unsigned int block_col,
372  const std::vector<types::global_dof_index> &dof1,
373  const std::vector<types::global_dof_index> &dof2,
374  const unsigned int level1,
375  const unsigned int level2,
376  bool transpose = false);
377 
381  void assemble_fluxes(MatrixType &global,
382  const FullMatrix<number> &local,
383  const unsigned int block_row,
384  const unsigned int block_col,
385  const std::vector<types::global_dof_index> &dof1,
386  const std::vector<types::global_dof_index> &dof2,
387  const unsigned int level1,
388  const unsigned int level2);
389 
393  void assemble_up(MatrixType &global,
394  const FullMatrix<number> &local,
395  const unsigned int block_row,
396  const unsigned int block_col,
397  const std::vector<types::global_dof_index> &dof1,
398  const std::vector<types::global_dof_index> &dof2,
399  const unsigned int level1,
400  const unsigned int level2);
401 
405  void assemble_down(MatrixType &global,
406  const FullMatrix<number> &local,
407  const unsigned int block_row,
408  const unsigned int block_col,
409  const std::vector<types::global_dof_index> &dof1,
410  const std::vector<types::global_dof_index> &dof2,
411  const unsigned int level1,
412  const unsigned int level2);
413 
417  void assemble_in(MatrixType &global,
418  const FullMatrix<number> &local,
419  const unsigned int block_row,
420  const unsigned int block_col,
421  const std::vector<types::global_dof_index> &dof1,
422  const std::vector<types::global_dof_index> &dof2,
423  const unsigned int level1,
424  const unsigned int level2);
425 
429  void assemble_out(MatrixType &global,
430  const FullMatrix<number> &local,
431  const unsigned int block_row,
432  const unsigned int block_col,
433  const std::vector<types::global_dof_index> &dof1,
434  const std::vector<types::global_dof_index> &dof2,
435  const unsigned int level1,
436  const unsigned int level2);
437 
442 
448 
454 
460 
466 
471 
476 
477 
483  const double threshold;
484 
485  };
486 
487 //----------------------------------------------------------------------//
488 
489  template <typename VectorType>
490  inline void
492  AnyData &m)
493  {
494  block_info = b;
495  residuals = m;
496  }
497 
498  template <typename VectorType>
499  inline void
501  (const ConstraintMatrix &c)
502  {
503  constraints = &c;
504  }
505 
506 
507  template <typename VectorType>
508  template <class DOFINFO>
509  inline void
511  (DOFINFO &info, bool) const
512  {
513  info.initialize_vectors(residuals.size());
514  }
515 
516  template <typename VectorType>
517  inline void
519  (VectorType &global,
520  const BlockVector<double> &local,
521  const std::vector<types::global_dof_index> &dof)
522  {
523  if (constraints == 0)
524  {
525  for (unsigned int b=0; b<local.n_blocks(); ++b)
526  for (unsigned int j=0; j<local.block(b).size(); ++j)
527  {
528  // The coordinates of
529  // the current entry in
530  // DoFHandler
531  // numbering, which
532  // differs from the
533  // block-wise local
534  // numbering we use in
535  // our local vectors
536  const unsigned int jcell = this->block_info->local().local_to_global(b, j);
537  global(dof[jcell]) += local.block(b)(j);
538  }
539  }
540  else
541  constraints->distribute_local_to_global(local, dof, global);
542  }
543 
544 
545  template <typename VectorType>
546  template <class DOFINFO>
547  inline void
549  (const DOFINFO &info)
550  {
551  for (unsigned int i=0; i<residuals.size(); ++i)
552  assemble(*(residuals.entry<VectorType>(i)), info.vector(i), info.indices);
553  }
554 
555 
556  template <typename VectorType>
557  template <class DOFINFO>
558  inline void
560  (const DOFINFO &info1,
561  const DOFINFO &info2)
562  {
563  for (unsigned int i=0; i<residuals.size(); ++i)
564  {
565  assemble(*(residuals.entry<VectorType>(i)), info1.vector(i), info1.indices);
566  assemble(*(residuals.entry<VectorType>(i)), info2.vector(i), info2.indices);
567  }
568  }
569 
570 
571 //----------------------------------------------------------------------//
572 
573  template <typename MatrixType, typename number>
574  inline
576  (double threshold)
577  :
578  threshold(threshold)
579  {}
580 
581 
582  template <typename MatrixType, typename number>
583  inline void
585  (const BlockInfo *b,
587  {
588  block_info = b;
589  matrices = &m;
590  }
591 
592 
593 
594  template <typename MatrixType, typename number>
595  inline void
597  (const ConstraintMatrix &c)
598  {
599  constraints = &c;
600  }
601 
602 
603 
604  template <typename MatrixType ,typename number>
605  template <class DOFINFO>
606  inline void
608  (DOFINFO &info,
609  bool face) const
610  {
611  info.initialize_matrices(*matrices, face);
612  }
613 
614 
615 
616  template <typename MatrixType, typename number>
617  inline void
620  const FullMatrix<number> &local,
621  const unsigned int block_row,
622  const unsigned int block_col,
623  const std::vector<types::global_dof_index> &dof1,
624  const std::vector<types::global_dof_index> &dof2)
625  {
626  if (constraints == 0)
627  {
628  for (unsigned int j=0; j<local.n_rows(); ++j)
629  for (unsigned int k=0; k<local.n_cols(); ++k)
630  if (std::fabs(local(j,k)) >= threshold)
631  {
632  // The coordinates of
633  // the current entry in
634  // DoFHandler
635  // numbering, which
636  // differs from the
637  // block-wise local
638  // numbering we use in
639  // our local matrices
640  const unsigned int jcell = this->block_info->local().local_to_global(block_row, j);
641  const unsigned int kcell = this->block_info->local().local_to_global(block_col, k);
642 
643  global.add(dof1[jcell], dof2[kcell], local(j,k));
644  }
645  }
646  else
647  {
648  const BlockIndices &bi = this->block_info->local();
649  std::vector<types::global_dof_index> sliced_row_indices (bi.block_size(block_row));
650  for (unsigned int i=0; i<sliced_row_indices.size(); ++i)
651  sliced_row_indices[i] = dof1[bi.block_start(block_row)+i];
652 
653  std::vector<types::global_dof_index> sliced_col_indices (bi.block_size(block_col));
654  for (unsigned int i=0; i<sliced_col_indices.size(); ++i)
655  sliced_col_indices[i] = dof2[bi.block_start(block_col)+i];
656 
657  constraints->distribute_local_to_global(local,
658  sliced_row_indices, sliced_col_indices, global);
659  }
660  }
661 
662 
663  template <typename MatrixType, typename number>
664  template <class DOFINFO>
665  inline void
667  {
668  for (unsigned int i=0; i<matrices->size(); ++i)
669  {
670  // Row and column index of
671  // the block we are dealing with
672  const types::global_dof_index row = matrices->block(i).row;
673  const types::global_dof_index col = matrices->block(i).column;
674 
675  assemble(matrices->block(i), info.matrix(i,false).matrix, row, col, info.indices, info.indices);
676  }
677  }
678 
679 
680  template <typename MatrixType, typename number>
681  template <class DOFINFO>
682  inline void
684  const DOFINFO &info2)
685  {
686  for (unsigned int i=0; i<matrices->size(); ++i)
687  {
688  // Row and column index of
689  // the block we are dealing with
690  const types::global_dof_index row = matrices->block(i).row;
691  const types::global_dof_index col = matrices->block(i).column;
692 
693  assemble(matrices->block(i), info1.matrix(i,false).matrix, row, col, info1.indices, info1.indices);
694  assemble(matrices->block(i), info1.matrix(i,true).matrix, row, col, info1.indices, info2.indices);
695  assemble(matrices->block(i), info2.matrix(i,false).matrix, row, col, info2.indices, info2.indices);
696  assemble(matrices->block(i), info2.matrix(i,true).matrix, row, col, info2.indices, info1.indices);
697  }
698  }
699 
700 
701 // ----------------------------------------------------------------------//
702 
703  template <typename MatrixType, typename number>
704  inline
706  (double threshold)
707  :
708  threshold(threshold)
709  {}
710 
711 
712  template <typename MatrixType, typename number>
713  inline void
715  (const BlockInfo *b,
716  MatrixPtrVector &m)
717  {
718  block_info = b;
719  AssertDimension(block_info->local().size(), block_info->global().size());
720  matrices = &m;
721  }
722 
723 
724  template <typename MatrixType, typename number>
725  inline void
727  (const MGConstrainedDoFs &mg_c)
728  {
729  mg_constrained_dofs = &mg_c;
730  }
731 
732 
733  template <typename MatrixType ,typename number>
734  template <class DOFINFO>
735  inline void
737  (DOFINFO &info,
738  bool face) const
739  {
740  info.initialize_matrices(*matrices, face);
741  }
742 
743 
744 
745  template <typename MatrixType, typename number>
746  inline void
749  MatrixPtrVector &down)
750  {
751  flux_up = up;
752  flux_down = down;
753  }
754 
755 
756  template <typename MatrixType, typename number>
757  inline void
760  MatrixPtrVector &out)
761  {
762  interface_in = in;
763  interface_out = out;
764  }
765 
766 
767  template <typename MatrixType, typename number>
768  inline void
770  (MatrixType &global,
771  const FullMatrix<number> &local,
772  const unsigned int block_row,
773  const unsigned int block_col,
774  const std::vector<types::global_dof_index> &dof1,
775  const std::vector<types::global_dof_index> &dof2,
776  const unsigned int level1,
777  const unsigned int level2,
778  bool transpose)
779  {
780  for (unsigned int j=0; j<local.n_rows(); ++j)
781  for (unsigned int k=0; k<local.n_cols(); ++k)
782  if (std::fabs(local(j,k)) >= threshold)
783  {
784  // The coordinates of
785  // the current entry in
786  // DoFHandler
787  // numbering, which
788  // differs from the
789  // block-wise local
790  // numbering we use in
791  // our local matrices
792  const unsigned int jcell = this->block_info->local().local_to_global(block_row, j);
793  const unsigned int kcell = this->block_info->local().local_to_global(block_col, k);
794 
795  // The global dof
796  // indices to assemble
797  // in. Since we may
798  // have face matrices
799  // coupling two
800  // different cells, we
801  // provide two sets of
802  // dof indices.
803  const unsigned int jglobal = this->block_info->level(level1).global_to_local(dof1[jcell]).second;
804  const unsigned int kglobal = this->block_info->level(level2).global_to_local(dof2[kcell]).second;
805 
806  if (mg_constrained_dofs == 0)
807  {
808  if (transpose)
809  global.add(kglobal, jglobal, local(j,k));
810  else
811  global.add(jglobal, kglobal, local(j,k));
812  }
813  else
814  {
815  if (!mg_constrained_dofs->at_refinement_edge(level1, jglobal) &&
816  !mg_constrained_dofs->at_refinement_edge(level2, kglobal))
817  {
818  if (mg_constrained_dofs->set_boundary_values())
819  {
820  if ((!mg_constrained_dofs->is_boundary_index(level1, jglobal) &&
821  !mg_constrained_dofs->is_boundary_index(level2, kglobal))
822  ||
823  (mg_constrained_dofs->is_boundary_index(level1, jglobal) &&
824  mg_constrained_dofs->is_boundary_index(level2, kglobal) &&
825  jglobal == kglobal))
826  {
827  if (transpose)
828  global.add(kglobal, jglobal, local(j,k));
829  else
830  global.add(jglobal, kglobal, local(j,k));
831  }
832  }
833  else
834  {
835  if (transpose)
836  global.add(kglobal, jglobal, local(j,k));
837  else
838  global.add(jglobal, kglobal, local(j,k));
839  }
840  }
841  }
842  }
843  }
844 
845 
846  template <typename MatrixType, typename number>
847  inline void
849  (MatrixType &global,
850  const FullMatrix<number> &local,
851  const unsigned int block_row,
852  const unsigned int block_col,
853  const std::vector<types::global_dof_index> &dof1,
854  const std::vector<types::global_dof_index> &dof2,
855  const unsigned int level1,
856  const unsigned int level2)
857  {
858  for (unsigned int j=0; j<local.n_rows(); ++j)
859  for (unsigned int k=0; k<local.n_cols(); ++k)
860  if (std::fabs(local(j,k)) >= threshold)
861  {
862  // The coordinates of
863  // the current entry in
864  // DoFHandler
865  // numbering, which
866  // differs from the
867  // block-wise local
868  // numbering we use in
869  // our local matrices
870  const unsigned int jcell = this->block_info->local().local_to_global(block_row, j);
871  const unsigned int kcell = this->block_info->local().local_to_global(block_col, k);
872 
873  // The global dof
874  // indices to assemble
875  // in. Since we may
876  // have face matrices
877  // coupling two
878  // different cells, we
879  // provide two sets of
880  // dof indices.
881  const unsigned int jglobal = this->block_info->level(level1).global_to_local(dof1[jcell]).second;
882  const unsigned int kglobal = this->block_info->level(level2).global_to_local(dof2[kcell]).second;
883 
884  if (mg_constrained_dofs == 0)
885  global.add(jglobal, kglobal, local(j,k));
886  else
887  {
888  if (!mg_constrained_dofs->non_refinement_edge_index(level1, jglobal) &&
889  !mg_constrained_dofs->non_refinement_edge_index(level2, kglobal))
890  {
891  if (!mg_constrained_dofs->at_refinement_edge(level1, jglobal) &&
892  !mg_constrained_dofs->at_refinement_edge(level2, kglobal))
893  global.add(jglobal, kglobal, local(j,k));
894  }
895  }
896  }
897  }
898 
899  template <typename MatrixType, typename number>
900  inline void
902  (MatrixType &global,
903  const FullMatrix<number> &local,
904  const unsigned int block_row,
905  const unsigned int block_col,
906  const std::vector<types::global_dof_index> &dof1,
907  const std::vector<types::global_dof_index> &dof2,
908  const unsigned int level1,
909  const unsigned int level2)
910  {
911  for (unsigned int j=0; j<local.n_rows(); ++j)
912  for (unsigned int k=0; k<local.n_cols(); ++k)
913  if (std::fabs(local(j,k)) >= threshold)
914  {
915  // The coordinates of
916  // the current entry in
917  // DoFHandler
918  // numbering, which
919  // differs from the
920  // block-wise local
921  // numbering we use in
922  // our local matrices
923  const unsigned int jcell = this->block_info->local().local_to_global(block_row, j);
924  const unsigned int kcell = this->block_info->local().local_to_global(block_col, k);
925 
926  // The global dof
927  // indices to assemble
928  // in. Since we may
929  // have face matrices
930  // coupling two
931  // different cells, we
932  // provide two sets of
933  // dof indices.
934  const unsigned int jglobal = this->block_info->level(level1).global_to_local(dof1[jcell]).second;
935  const unsigned int kglobal = this->block_info->level(level2).global_to_local(dof2[kcell]).second;
936 
937  if (mg_constrained_dofs == 0)
938  global.add(jglobal, kglobal, local(j,k));
939  else
940  {
941  if (!mg_constrained_dofs->non_refinement_edge_index(level1, jglobal) &&
942  !mg_constrained_dofs->non_refinement_edge_index(level2, kglobal))
943  {
944  if (!mg_constrained_dofs->at_refinement_edge(level1, jglobal) &&
945  !mg_constrained_dofs->at_refinement_edge(level2, kglobal))
946  global.add(jglobal, kglobal, local(j,k));
947  }
948  }
949  }
950  }
951 
952  template <typename MatrixType, typename number>
953  inline void
955  (MatrixType &global,
956  const FullMatrix<number> &local,
957  const unsigned int block_row,
958  const unsigned int block_col,
959  const std::vector<types::global_dof_index> &dof1,
960  const std::vector<types::global_dof_index> &dof2,
961  const unsigned int level1,
962  const unsigned int level2)
963  {
964  for (unsigned int j=0; j<local.n_rows(); ++j)
965  for (unsigned int k=0; k<local.n_cols(); ++k)
966  if (std::fabs(local(k,j)) >= threshold)
967  {
968  // The coordinates of
969  // the current entry in
970  // DoFHandler
971  // numbering, which
972  // differs from the
973  // block-wise local
974  // numbering we use in
975  // our local matrices
976  const unsigned int jcell = this->block_info->local().local_to_global(block_row, j);
977  const unsigned int kcell = this->block_info->local().local_to_global(block_col, k);
978 
979  // The global dof
980  // indices to assemble
981  // in. Since we may
982  // have face matrices
983  // coupling two
984  // different cells, we
985  // provide two sets of
986  // dof indices.
987  const unsigned int jglobal = this->block_info->level(level1).global_to_local(dof1[jcell]).second;
988  const unsigned int kglobal = this->block_info->level(level2).global_to_local(dof2[kcell]).second;
989 
990  if (mg_constrained_dofs == 0)
991  global.add(jglobal, kglobal, local(k,j));
992  else
993  {
994  if (!mg_constrained_dofs->non_refinement_edge_index(level1, jglobal) &&
995  !mg_constrained_dofs->non_refinement_edge_index(level2, kglobal))
996  {
997  if (!mg_constrained_dofs->at_refinement_edge(level1, jglobal) &&
998  !mg_constrained_dofs->at_refinement_edge(level2, kglobal))
999  global.add(jglobal, kglobal, local(k,j));
1000  }
1001  }
1002  }
1003  }
1004 
1005  template <typename MatrixType, typename number>
1006  inline void
1008  (MatrixType &global,
1009  const FullMatrix<number> &local,
1010  const unsigned int block_row,
1011  const unsigned int block_col,
1012  const std::vector<types::global_dof_index> &dof1,
1013  const std::vector<types::global_dof_index> &dof2,
1014  const unsigned int level1,
1015  const unsigned int level2)
1016  {
1017 // AssertDimension(local.n(), dof1.size());
1018 // AssertDimension(local.m(), dof2.size());
1019 
1020  for (unsigned int j=0; j<local.n_rows(); ++j)
1021  for (unsigned int k=0; k<local.n_cols(); ++k)
1022  if (std::fabs(local(j,k)) >= threshold)
1023  {
1024  // The coordinates of
1025  // the current entry in
1026  // DoFHandler
1027  // numbering, which
1028  // differs from the
1029  // block-wise local
1030  // numbering we use in
1031  // our local matrices
1032  const unsigned int jcell = this->block_info->local().local_to_global(block_row, j);
1033  const unsigned int kcell = this->block_info->local().local_to_global(block_col, k);
1034 
1035  // The global dof
1036  // indices to assemble
1037  // in. Since we may
1038  // have face matrices
1039  // coupling two
1040  // different cells, we
1041  // provide two sets of
1042  // dof indices.
1043  const unsigned int jglobal = this->block_info->level(level1).global_to_local(dof1[jcell]).second;
1044  const unsigned int kglobal = this->block_info->level(level2).global_to_local(dof2[kcell]).second;
1045 
1046  if (mg_constrained_dofs == 0)
1047  global.add(jglobal, kglobal, local(j,k));
1048  else
1049  {
1050  if (mg_constrained_dofs->at_refinement_edge(level1, jglobal) &&
1051  !mg_constrained_dofs->at_refinement_edge(level2, kglobal))
1052  {
1053  if (mg_constrained_dofs->set_boundary_values())
1054  {
1055  if ((!mg_constrained_dofs->is_boundary_index(level1, jglobal) &&
1056  !mg_constrained_dofs->is_boundary_index(level2, kglobal))
1057  ||
1058  (mg_constrained_dofs->is_boundary_index(level1, jglobal) &&
1059  mg_constrained_dofs->is_boundary_index(level2, kglobal) &&
1060  jglobal == kglobal))
1061  global.add(jglobal, kglobal, local(j,k));
1062  }
1063  else
1064  global.add(jglobal, kglobal, local(j,k));
1065  }
1066  }
1067  }
1068  }
1069 
1070  template <typename MatrixType, typename number>
1071  inline void
1073  (MatrixType &global,
1074  const FullMatrix<number> &local,
1075  const unsigned int block_row,
1076  const unsigned int block_col,
1077  const std::vector<types::global_dof_index> &dof1,
1078  const std::vector<types::global_dof_index> &dof2,
1079  const unsigned int level1,
1080  const unsigned int level2)
1081  {
1082 // AssertDimension(local.n(), dof1.size());
1083 // AssertDimension(local.m(), dof2.size());
1084 
1085  for (unsigned int j=0; j<local.n_rows(); ++j)
1086  for (unsigned int k=0; k<local.n_cols(); ++k)
1087  if (std::fabs(local(k,j)) >= threshold)
1088  {
1089  // The coordinates of
1090  // the current entry in
1091  // DoFHandler
1092  // numbering, which
1093  // differs from the
1094  // block-wise local
1095  // numbering we use in
1096  // our local matrices
1097  const unsigned int jcell = this->block_info->local().local_to_global(block_row, j);
1098  const unsigned int kcell = this->block_info->local().local_to_global(block_col, k);
1099 
1100  // The global dof
1101  // indices to assemble
1102  // in. Since we may
1103  // have face matrices
1104  // coupling two
1105  // different cells, we
1106  // provide two sets of
1107  // dof indices.
1108  const unsigned int jglobal = this->block_info->level(level1).global_to_local(dof1[jcell]).second;
1109  const unsigned int kglobal = this->block_info->level(level2).global_to_local(dof2[kcell]).second;
1110 
1111  if (mg_constrained_dofs == 0)
1112  global.add(jglobal, kglobal, local(k,j));
1113  else
1114  {
1115  if (mg_constrained_dofs->at_refinement_edge(level1, jglobal) &&
1116  !mg_constrained_dofs->at_refinement_edge(level2, kglobal))
1117  {
1118  if (mg_constrained_dofs->set_boundary_values())
1119  {
1120  if ((!mg_constrained_dofs->is_boundary_index(level1, jglobal) &&
1121  !mg_constrained_dofs->is_boundary_index(level2, kglobal))
1122  ||
1123  (mg_constrained_dofs->is_boundary_index(level1, jglobal) &&
1124  mg_constrained_dofs->is_boundary_index(level2, kglobal) &&
1125  jglobal == kglobal))
1126  global.add(jglobal, kglobal, local(k,j));
1127  }
1128  else
1129  global.add(jglobal, kglobal, local(k,j));
1130  }
1131  }
1132  }
1133  }
1134 
1135 
1136  template <typename MatrixType, typename number>
1137  template <class DOFINFO>
1138  inline void
1140  {
1141  const unsigned int level = info.cell->level();
1142 
1143  for (unsigned int i=0; i<matrices->size(); ++i)
1144  {
1145  // Row and column index of
1146  // the block we are dealing with
1147  const unsigned int row = matrices->block(i)[level].row;
1148  const unsigned int col = matrices->block(i)[level].column;
1149 
1150  assemble(matrices->block(i)[level].matrix, info.matrix(i,false).matrix, row, col,
1151  info.indices, info.indices, level, level);
1152  if (mg_constrained_dofs != 0)
1153  {
1154  if (interface_in != 0)
1155  assemble_in(interface_in->block(i)[level], info.matrix(i,false).matrix, row, col,
1156  info.indices, info.indices, level, level);
1157  if (interface_out != 0)
1158  assemble_in(interface_out->block(i)[level], info.matrix(i,false).matrix, row, col,
1159  info.indices, info.indices, level, level);
1160 
1161  assemble_in(matrices->block_in(i)[level], info.matrix(i,false).matrix, row, col,
1162  info.indices, info.indices, level, level);
1163  assemble_out(matrices->block_out(i)[level], info.matrix(i,false).matrix, row, col,
1164  info.indices, info.indices, level, level);
1165  }
1166  }
1167  }
1168 
1169 
1170  template <typename MatrixType, typename number>
1171  template <class DOFINFO>
1172  inline void
1174  (const DOFINFO &info1,
1175  const DOFINFO &info2)
1176  {
1177  const unsigned int level1 = info1.cell->level();
1178  const unsigned int level2 = info2.cell->level();
1179 
1180  for (unsigned int i=0; i<matrices->size(); ++i)
1181  {
1182  MGLevelObject<MatrixBlock<MatrixType> > &o = matrices->block(i);
1183 
1184  // Row and column index of
1185  // the block we are dealing with
1186  const unsigned int row = o[level1].row;
1187  const unsigned int col = o[level1].column;
1188 
1189  if (level1 == level2)
1190  {
1191  if (mg_constrained_dofs == 0)
1192  {
1193  assemble(o[level1].matrix, info1.matrix(i,false).matrix, row, col,
1194  info1.indices, info1.indices, level1, level1);
1195  assemble(o[level1].matrix, info1.matrix(i,true).matrix, row, col,
1196  info1.indices, info2.indices, level1, level2);
1197  assemble(o[level1].matrix, info2.matrix(i,false).matrix, row, col,
1198  info2.indices, info2.indices, level2, level2);
1199  assemble(o[level1].matrix, info2.matrix(i,true).matrix, row, col,
1200  info2.indices, info1.indices, level2, level1);
1201  }
1202  else
1203  {
1204  assemble_fluxes(o[level1].matrix, info1.matrix(i,false).matrix, row, col,
1205  info1.indices, info1.indices, level1, level1);
1206  assemble_fluxes(o[level1].matrix, info1.matrix(i,true).matrix, row, col,
1207  info1.indices, info2.indices, level1, level2);
1208  assemble_fluxes(o[level1].matrix, info2.matrix(i,false).matrix, row, col,
1209  info2.indices, info2.indices, level2, level2);
1210  assemble_fluxes(o[level1].matrix, info2.matrix(i,true).matrix, row, col,
1211  info2.indices, info1.indices, level2, level1);
1212  }
1213  }
1214  else
1215  {
1216  Assert(level1 > level2, ExcNotImplemented());
1217  if (flux_up->size() != 0)
1218  {
1219  // Do not add M22,
1220  // which is done by
1221  // the coarser cell
1222  assemble_fluxes(o[level1].matrix, info1.matrix(i,false).matrix, row, col,
1223  info1.indices, info1.indices, level1, level1);
1224  assemble_up(flux_up->block(i)[level1].matrix, info1.matrix(i,true).matrix, row, col,
1225  info1.indices, info2.indices, level1, level2);
1226  assemble_down(flux_down->block(i)[level1].matrix, info2.matrix(i,true).matrix, row, col,
1227  info2.indices, info1.indices, level2, level1);
1228  }
1229  }
1230  }
1231  }
1232  }
1233 }
1234 
1235 DEAL_II_NAMESPACE_CLOSE
1236 
1237 #endif
SmartPointer< const BlockInfo, MatrixLocalBlocksToGlobalBlocks< MatrixType, number > > block_info
Definition: assembler.h:261
void assemble_down(MatrixType &global, const FullMatrix< number > &local, const unsigned int block_row, const unsigned int block_col, const std::vector< types::global_dof_index > &dof1, const std::vector< types::global_dof_index > &dof2, const unsigned int level1, const unsigned int level2)
Definition: assembler.h:955
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
void initialize_edge_flux(MatrixPtrVector &up, MatrixPtrVector &down)
Definition: assembler.h:748
void initialize(const BlockInfo *block_info, AnyData &residuals)
Definition: assembler.h:491
void initialize_info(DOFINFO &info, bool face) const
Definition: assembler.h:608
void add(const size_type i, const size_type j, const typename MatrixType::value_type value)
Definition: matrix_block.h:649
size_type block_size(const unsigned int i) const
void initialize_info(DOFINFO &info, bool face) const
Definition: assembler.h:511
SmartPointer< const ConstraintMatrix, MatrixLocalBlocksToGlobalBlocks< MatrixType, number > > constraints
Definition: assembler.h:266
SmartPointer< const MGConstrainedDoFs, MGMatrixLocalBlocksToGlobalBlocks< MatrixType, number > > mg_constrained_dofs
Definition: assembler.h:475
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:313
void assemble_fluxes(MatrixType &global, const FullMatrix< number > &local, const unsigned int block_row, const unsigned int block_col, const std::vector< types::global_dof_index > &dof1, const std::vector< types::global_dof_index > &dof2, const unsigned int level1, const unsigned int level2)
Definition: assembler.h:849
SmartPointer< const BlockInfo, MGMatrixLocalBlocksToGlobalBlocks< MatrixType, number > > block_info
Definition: assembler.h:470
void initialize(const BlockInfo *block_info, MatrixBlockVector< MatrixType > &matrices)
Definition: assembler.h:585
SmartPointer< const BlockInfo, ResidualLocalBlocksToGlobalBlocks< VectorType > > block_info
Definition: assembler.h:159
void assemble_out(MatrixType &global, const FullMatrix< number > &local, const unsigned int block_row, const unsigned int block_col, const std::vector< types::global_dof_index > &dof1, const std::vector< types::global_dof_index > &dof2, const unsigned int level1, const unsigned int level2)
Definition: assembler.h:1073
void initialize_interfaces(MatrixPtrVector &interface_in, MatrixPtrVector &interface_out)
Definition: assembler.h:759
SmartPointer< const ConstraintMatrix, ResidualLocalBlocksToGlobalBlocks< VectorType > > constraints
Definition: assembler.h:164
void assemble_up(MatrixType &global, const FullMatrix< number > &local, const unsigned int block_row, const unsigned int block_col, const std::vector< types::global_dof_index > &dof1, const std::vector< types::global_dof_index > &dof2, const unsigned int level1, const unsigned int level2)
Definition: assembler.h:902
size_type block_start(const unsigned int i) const
A small class collecting the different BlockIndices involved in global, multilevel and local computat...
Definition: block_info.h:93
static ::ExceptionBase & ExcNotImplemented()
SmartPointer< MatrixBlockVector< MatrixType >, MatrixLocalBlocksToGlobalBlocks< MatrixType, number > > matrices
Definition: assembler.h:255
BlockType & block(const unsigned int i)
void initialize_info(DOFINFO &info, bool face) const
Definition: assembler.h:737
MatrixLocalBlocksToGlobalBlocks(double threshold=1.e-12)
Definition: assembler.h:576
void initialize(const BlockInfo *block_info, MatrixPtrVector &matrices)
Definition: assembler.h:715
void assemble_in(MatrixType &global, const FullMatrix< number > &local, const unsigned int block_row, const unsigned int block_col, const std::vector< types::global_dof_index > &dof1, const std::vector< types::global_dof_index > &dof2, const unsigned int level1, const unsigned int level2)
Definition: assembler.h:1008