Reference documentation for deal.II version 8.5.1
operators.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2011 - 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__matrix_free_operators_h
18 #define dealii__matrix_free_operators_h
19 
20 
21 #include <deal.II/base/exceptions.h>
22 #include <deal.II/base/vectorization.h>
23 #include <deal.II/base/subscriptor.h>
24 
25 #include <deal.II/lac/diagonal_matrix.h>
26 #include <deal.II/lac/la_parallel_vector.h>
27 #include <deal.II/lac/vector_view.h>
28 #include <deal.II/multigrid/mg_constrained_dofs.h>
29 #include <deal.II/matrix_free/matrix_free.h>
30 #include <deal.II/matrix_free/fe_evaluation.h>
31 
32 
33 DEAL_II_NAMESPACE_OPEN
34 
35 
36 namespace MatrixFreeOperators
37 {
38 
39  namespace
40  {
41  // workaroud for unifying non-block vector and block vector implementations
42  // a non-block vector has one block and the only subblock is the vector itself
43  template <typename VectorType>
44  typename std_cxx11::enable_if<IsBlockVector<VectorType>::value,
45  unsigned int>::type
46  n_blocks(const VectorType &vector)
47  {
48  return vector.n_blocks();
49  }
50 
51  template <typename VectorType>
52  typename std_cxx11::enable_if<!IsBlockVector<VectorType>::value,
53  unsigned int>::type
54  n_blocks(const VectorType &)
55  {
56  return 1;
57  }
58 
59  template <typename VectorType>
60  typename std_cxx11::enable_if<IsBlockVector<VectorType>::value,
61  typename VectorType::BlockType &>::type
62  subblock(VectorType &vector, unsigned int block_no)
63  {
64  return vector.block(block_no);
65  }
66 
67  template <typename VectorType>
68  typename std_cxx11::enable_if<IsBlockVector<VectorType>::value,
69  const typename VectorType::BlockType &>::type
70  subblock(const VectorType &vector, unsigned int block_no)
71  {
72  AssertIndexRange (block_no, vector.n_blocks());
73  return vector.block(block_no);
74  }
75 
76  template <typename VectorType>
77  typename std_cxx11::enable_if<!IsBlockVector<VectorType>::value,
78  VectorType &>::type
79  subblock(VectorType &vector, unsigned int)
80  {
81  return vector;
82  }
83 
84  template <typename VectorType>
85  typename std_cxx11::enable_if<!IsBlockVector<VectorType>::value,
86  const VectorType &>::type
87  subblock(const VectorType &vector, unsigned int)
88  {
89  return vector;
90  }
91 
92  template <typename VectorType>
93  typename std_cxx11::enable_if<IsBlockVector<VectorType>::value,
94  void>::type
95  collect_sizes(VectorType &vector)
96  {
97  vector.collect_sizes();
98  }
99 
100  template <typename VectorType>
101  typename std_cxx11::enable_if<!IsBlockVector<VectorType>::value,
102  void>::type
103  collect_sizes(const VectorType &)
104  {}
105  }
106 
163  template <int dim, typename VectorType = LinearAlgebra::distributed::Vector<double> >
164  class Base : public Subscriptor
165  {
166  public:
170  typedef typename VectorType::value_type value_type;
171 
175  typedef typename VectorType::size_type size_type;
176 
180  Base ();
181 
185  virtual ~Base();
186 
191  virtual void clear();
192 
209  void initialize (std_cxx11::shared_ptr<const MatrixFree<dim,value_type> > data,
210  const std::vector<unsigned int> &selected_row_blocks = std::vector<unsigned int>(),
211  const std::vector<unsigned int> &selected_column_blocks = std::vector<unsigned int>());
212 
226  void initialize (std_cxx11::shared_ptr<const MatrixFree<dim,value_type> > data,
227  const MGConstrainedDoFs &mg_constrained_dofs,
228  const unsigned int level,
229  const std::vector<unsigned int> &selected_row_blocks = std::vector<unsigned int>());
230 
245  void
246  initialize(std_cxx11::shared_ptr<const MatrixFree<dim, value_type> > data_,
247  const std::vector<MGConstrainedDoFs> &mg_constrained_dofs,
248  const unsigned int level,
249  const std::vector<unsigned int> &selected_row_blocks = std::vector<unsigned int>());
250 
254  size_type m () const;
255 
259  size_type n () const;
260 
264  void vmult_interface_down(VectorType &dst,
265  const VectorType &src) const;
266 
270  void vmult_interface_up(VectorType &dst,
271  const VectorType &src) const;
272 
276  void vmult (VectorType &dst,
277  const VectorType &src) const;
278 
282  void Tvmult (VectorType &dst,
283  const VectorType &src) const;
284 
288  void vmult_add (VectorType &dst,
289  const VectorType &src) const;
290 
294  void Tvmult_add (VectorType &dst,
295  const VectorType &src) const;
296 
301  value_type el (const unsigned int row,
302  const unsigned int col) const;
303 
307  virtual std::size_t memory_consumption () const;
308 
312  void initialize_dof_vector (VectorType &vec) const;
313 
320  virtual void compute_diagonal () = 0;
321 
325  std_cxx11::shared_ptr<const MatrixFree<dim,value_type> >
326  get_matrix_free () const;
327 
331  const std_cxx11::shared_ptr<DiagonalMatrix<VectorType> > &
333 
339  void precondition_Jacobi(VectorType &dst,
340  const VectorType &src,
341  const value_type omega) const;
342 
343  protected:
344 
349  void set_constrained_entries_to_one (VectorType &dst) const;
350 
354  virtual void apply_add(VectorType &dst,
355  const VectorType &src) const = 0;
356 
362  virtual void Tapply_add(VectorType &dst,
363  const VectorType &src) const;
364 
368  std_cxx11::shared_ptr<const MatrixFree<dim,value_type> > data;
369 
374  std_cxx11::shared_ptr<DiagonalMatrix<VectorType > > inverse_diagonal_entries;
375 
380  std::vector<unsigned int> selected_rows;
381 
386  std::vector<unsigned int> selected_columns;
387 
388  private:
389 
393  std::vector<std::vector<unsigned int> > edge_constrained_indices;
394 
398  mutable std::vector<std::vector<std::pair<value_type,value_type> > >
400 
406 
411  void mult_add (VectorType &dst,
412  const VectorType &src,
413  const bool transpose) const;
414 
422  void adjust_ghost_range_if_necessary(const VectorType &vec,
423  const bool is_row) const;
424  };
425 
426 
427 
464  template <typename OperatorType>
466  {
467  public:
471  typedef typename OperatorType::value_type value_type;
472 
477 
481  void clear();
482 
486  void initialize (const OperatorType &operator_in);
487 
491  template <typename VectorType>
492  void vmult (VectorType &dst,
493  const VectorType &src) const;
494 
498  template <typename VectorType>
499  void Tvmult (VectorType &dst,
500  const VectorType &src) const;
501 
502  private:
507  };
508 
509 
510 
530  template <int dim, int fe_degree, int n_components = 1, typename Number = double>
532  {
533  public:
539 
548  void apply(const AlignedVector<VectorizedArray<Number> > &inverse_coefficient,
549  const unsigned int n_actual_components,
550  const VectorizedArray<Number> *in_array,
551  VectorizedArray<Number> *out_array) const;
552 
559 
560  private:
565 
570  };
571 
572 
573 
583  template <int dim, int fe_degree, int n_q_points_1d = fe_degree+1, int n_components = 1, typename VectorType = LinearAlgebra::distributed::Vector<double> >
584  class MassOperator : public Base<dim, VectorType>
585  {
586  public:
591 
596 
600  MassOperator ();
601 
605  virtual void compute_diagonal ();
606 
607  private:
613  virtual void apply_add (VectorType &dst,
614  const VectorType &src) const;
615 
620  VectorType &dst,
621  const VectorType &src,
622  const std::pair<unsigned int,unsigned int> &cell_range) const;
623  };
624 
625 
626 
638  template <int dim, int fe_degree, int n_q_points_1d = fe_degree+1, int n_components = 1, typename VectorType = LinearAlgebra::distributed::Vector<double> >
639  class LaplaceOperator : public Base<dim, VectorType>
640  {
641  public:
646 
651 
655  LaplaceOperator ();
656 
663  virtual void compute_diagonal ();
664 
708  void set_coefficient(const std_cxx11::shared_ptr<Table<2, VectorizedArray<value_type> > > &scalar_coefficient );
709 
710  virtual void clear();
711 
718  std_cxx11::shared_ptr< Table<2, VectorizedArray<value_type> > > get_coefficient();
719 
720  private:
726  virtual void apply_add (VectorType &dst,
727  const VectorType &src) const;
728 
733  VectorType &dst,
734  const VectorType &src,
735  const std::pair<unsigned int,unsigned int> &cell_range) const;
736 
741  VectorType &dst,
742  const unsigned int &,
743  const std::pair<unsigned int,unsigned int> &cell_range) const;
744 
749  const unsigned int cell) const;
750 
754  std_cxx11::shared_ptr< Table<2, VectorizedArray<value_type> > > scalar_coefficient;
755  };
756 
757 
758 
759  // ------------------------------------ inline functions ---------------------
760 
761  template <int dim, int fe_degree, int n_components, typename Number>
762  inline
765  :
766  fe_eval (fe_eval)
767  {
768  FullMatrix<double> shapes_1d(fe_degree+1, fe_degree+1);
769  for (unsigned int i=0, c=0; i<shapes_1d.m(); ++i)
770  for (unsigned int j=0; j<shapes_1d.n(); ++j, ++c)
771  shapes_1d(i,j) = fe_eval.get_shape_info().shape_values_number[c];
772  shapes_1d.gauss_jordan();
773  const unsigned int stride = (fe_degree+2)/2;
774  inverse_shape.resize(stride*(fe_degree+1));
775  for (unsigned int i=0; i<stride; ++i)
776  for (unsigned int q=0; q<(fe_degree+2)/2; ++q)
777  {
778  inverse_shape[i*stride+q] =
779  0.5 * (shapes_1d(i,q) + shapes_1d(i,fe_degree-q));
780  inverse_shape[(fe_degree-i)*stride+q] =
781  0.5 * (shapes_1d(i,q) - shapes_1d(i,fe_degree-q));
782  }
783  if (fe_degree % 2 == 0)
784  for (unsigned int q=0; q<(fe_degree+2)/2; ++q)
785  inverse_shape[fe_degree/2*stride+q] = shapes_1d(fe_degree/2,q);
786  }
787 
788 
789 
790  template <int dim, int fe_degree, int n_components, typename Number>
791  inline
792  void
795  {
796  const unsigned int dofs_per_cell = Utilities::fixed_int_power<fe_degree+1,dim>::value;
797  Assert(inverse_jxw.size() > 0 &&
798  inverse_jxw.size() % dofs_per_cell == 0,
799  ExcMessage("Expected diagonal to be a multiple of scalar dof per cells"));
800 
801  // temporarily reduce size of inverse_jxw to dofs_per_cell to get JxW values
802  // from fe_eval (will not reallocate any memory)
803  const unsigned int previous_size = inverse_jxw.size();
804  inverse_jxw.resize(dofs_per_cell);
805  fe_eval.fill_JxW_values(inverse_jxw);
806 
807  // invert
808  inverse_jxw.resize_fast(previous_size);
809  for (unsigned int q=0; q<dofs_per_cell; ++q)
810  inverse_jxw[q] = 1./inverse_jxw[q];
811  // copy values to rest of vector
812  for (unsigned int q=dofs_per_cell; q<previous_size; )
813  for (unsigned int i=0; i<dofs_per_cell; ++i, ++q)
814  inverse_jxw[q] = inverse_jxw[i];
815  }
816 
817 
818 
819  template <int dim, int fe_degree, int n_components, typename Number>
820  inline
821  void
823  ::apply(const AlignedVector<VectorizedArray<Number> > &inverse_coefficients,
824  const unsigned int n_actual_components,
825  const VectorizedArray<Number> *in_array,
826  VectorizedArray<Number> *out_array) const
827  {
828  const unsigned int dofs_per_cell = Utilities::fixed_int_power<fe_degree+1,dim>::value;
829  Assert(inverse_coefficients.size() > 0 &&
830  inverse_coefficients.size() % dofs_per_cell == 0,
831  ExcMessage("Expected diagonal to be a multiple of scalar dof per cells"));
832  if (inverse_coefficients.size() != dofs_per_cell)
833  AssertDimension(n_actual_components * dofs_per_cell, inverse_coefficients.size());
834 
835  Assert(dim == 2 || dim == 3, ExcNotImplemented());
836 
838  fe_degree+1, VectorizedArray<Number> >
839  evaluator(inverse_shape, inverse_shape, inverse_shape);
840 
841  const unsigned int shift_coefficient =
842  inverse_coefficients.size() > dofs_per_cell ? dofs_per_cell : 0;
843  const VectorizedArray<Number> *inv_coefficient = &inverse_coefficients[0];
844  VectorizedArray<Number> temp_data_field[dofs_per_cell];
845  for (unsigned int d=0; d<n_actual_components; ++d)
846  {
847  const VectorizedArray<Number> *in = in_array+d*dofs_per_cell;
848  VectorizedArray<Number> *out = out_array+d*dofs_per_cell;
849  // Need to select 'apply' method with hessian slot because values
850  // assume symmetries that do not exist in the inverse shapes
851  evaluator.template hessians<0,false,false> (in, temp_data_field);
852  evaluator.template hessians<1,false,false> (temp_data_field, out);
853 
854  if (dim == 3)
855  {
856  evaluator.template hessians<2,false,false> (out, temp_data_field);
857  for (unsigned int q=0; q<dofs_per_cell; ++q)
858  temp_data_field[q] *= inv_coefficient[q];
859  evaluator.template hessians<2,true,false> (temp_data_field, out);
860  }
861  else if (dim == 2)
862  for (unsigned int q=0; q<dofs_per_cell; ++q)
863  out[q] *= inv_coefficient[q];
864 
865  evaluator.template hessians<1,true,false>(out, temp_data_field);
866  evaluator.template hessians<0,true,false>(temp_data_field, out);
867 
868  inv_coefficient += shift_coefficient;
869  }
870  }
871 
872  //----------------- Base operator -----------------------------
873  template <int dim, typename VectorType>
875  {
876  }
877 
878 
879 
880  template <int dim, typename VectorType>
882  :
883  Subscriptor(),
884  have_interface_matrices(false)
885  {
886  // boost-1.62.0 doesn't allow initializing a shared_ptr
887  // with NULL. Make sure the default constructor does that.
888  Assert(data.get() == NULL, ExcInternalError());
889  }
890 
891 
892 
893  template <int dim, typename VectorType>
896  {
897  Assert(data.get() != NULL,
899  typename Base<dim, VectorType>::size_type total_size = 0;
900  for (unsigned int i=0; i<selected_rows.size(); ++i)
901  total_size += data->get_vector_partitioner(selected_rows[i])->size();
902  return total_size;
903  }
904 
905 
906 
907  template <int dim, typename VectorType>
910  {
911  Assert(data.get() != NULL,
913  typename Base<dim, VectorType>::size_type total_size = 0;
914  for (unsigned int i=0; i<selected_columns.size(); ++i)
915  total_size += data->get_vector_partitioner(selected_columns[i])->size();
916  return total_size;
917  }
918 
919 
920 
921  template <int dim, typename VectorType>
922  void
924  {
925  data.reset();
926  inverse_diagonal_entries.reset();
927  }
928 
929 
930 
931  template <int dim, typename VectorType>
933  Base<dim,VectorType>::el (const unsigned int row,
934  const unsigned int col) const
935  {
936  (void) col;
937  Assert (row == col, ExcNotImplemented());
938  Assert (inverse_diagonal_entries.get() != NULL &&
939  inverse_diagonal_entries->m() > 0, ExcNotInitialized());
940  return 1.0/(*inverse_diagonal_entries)(row,row);
941  }
942 
943 
944 
945  template <int dim, typename VectorType>
946  void
948  {
949  Assert(data.get() != NULL,
951  AssertDimension(n_blocks(vec), selected_rows.size());
952  for (unsigned int i = 0; i < n_blocks(vec); ++i)
953  {
954  const unsigned int index = selected_rows[i];
955  if (!subblock(vec,index).partitioners_are_compatible
956  (*data->get_dof_info(index).vector_partitioner))
957  data->initialize_dof_vector(subblock(vec,index));
958  Assert(subblock(vec,index).partitioners_are_globally_compatible
959  (*data->get_dof_info(index).vector_partitioner),
960  ExcInternalError());
961  }
962  collect_sizes(vec);
963  }
964 
965 
966 
967  template <int dim, typename VectorType>
968  void
970  initialize (std_cxx11::shared_ptr<const MatrixFree<dim,typename Base<dim,VectorType>::value_type> > data_,
971  const std::vector<unsigned int> &given_row_selection,
972  const std::vector<unsigned int> &given_column_selection)
973  {
974  data = data_;
975 
976  selected_rows.clear();
977  selected_columns.clear();
978  if (given_row_selection.empty())
979  for (unsigned int i=0; i<data_->n_components(); ++i)
980  selected_rows.push_back(i);
981  else
982  {
983  for (unsigned int i=0; i<given_row_selection.size(); ++i)
984  {
985  AssertIndexRange(given_row_selection[i], data_->n_components());
986  for (unsigned int j=0; j<given_row_selection.size(); ++j)
987  if (j!=i)
988  Assert(given_row_selection[j] != given_row_selection[i],
989  ExcMessage("Given row indices must be unique"));
990 
991  selected_rows.push_back(given_row_selection[i]);
992  }
993  }
994  if (given_column_selection.size() == 0)
995  selected_columns = selected_rows;
996  else
997  {
998  for (unsigned int i=0; i<given_column_selection.size(); ++i)
999  {
1000  AssertIndexRange(given_column_selection[i], data_->n_components());
1001  for (unsigned int j=0; j<given_column_selection.size(); ++j)
1002  if (j!=i)
1003  Assert(given_column_selection[j] != given_column_selection[i],
1004  ExcMessage("Given column indices must be unique"));
1005 
1006  selected_columns.push_back(given_column_selection[i]);
1007  }
1008  }
1009 
1010  edge_constrained_indices.clear();
1011  edge_constrained_indices.resize(selected_rows.size());
1012  edge_constrained_values.clear();
1013  edge_constrained_values.resize(selected_rows.size());
1014  have_interface_matrices = false;
1015  }
1016 
1017 
1018 
1019  template <int dim, typename VectorType>
1020  void
1022  initialize (std_cxx11::shared_ptr<const MatrixFree<dim,typename Base<dim,VectorType>::value_type> > data_,
1023  const MGConstrainedDoFs &mg_constrained_dofs,
1024  const unsigned int level,
1025  const std::vector<unsigned int> &given_row_selection)
1026  {
1027  std::vector<MGConstrainedDoFs> mg_constrained_dofs_vector(1, mg_constrained_dofs);
1028  initialize(data_, mg_constrained_dofs_vector, level, given_row_selection);
1029  }
1030 
1031 
1032 
1033  template <int dim, typename VectorType>
1034  void
1036  initialize (std_cxx11::shared_ptr<const MatrixFree<dim, value_type> > data_,
1037  const std::vector<MGConstrainedDoFs> &mg_constrained_dofs,
1038  const unsigned int level,
1039  const std::vector<unsigned int> &given_row_selection)
1040  {
1042  ExcMessage("level is not set"));
1043 
1044  selected_rows.clear();
1045  selected_columns.clear();
1046  if (given_row_selection.empty())
1047  for (unsigned int i=0; i<data_->n_components(); ++i)
1048  selected_rows.push_back(i);
1049  else
1050  {
1051  for (unsigned int i=0; i<given_row_selection.size(); ++i)
1052  {
1053  AssertIndexRange(given_row_selection[i], data_->n_components());
1054  for (unsigned int j=0; j<given_row_selection.size(); ++j)
1055  if (j!=i)
1056  Assert(given_row_selection[j] != given_row_selection[i],
1057  ExcMessage("Given row indices must be unique"));
1058 
1059  selected_rows.push_back(given_row_selection[i]);
1060  }
1061  }
1062  selected_columns = selected_rows;
1063 
1064  AssertDimension(mg_constrained_dofs.size(), selected_rows.size());
1065  edge_constrained_indices.clear();
1066  edge_constrained_indices.resize(selected_rows.size());
1067  edge_constrained_values.clear();
1068  edge_constrained_values.resize(selected_rows.size());
1069 
1070  data = data_;
1071 
1072  for (unsigned int j = 0; j < selected_rows.size(); ++j)
1073  {
1074  if (data_->n_macro_cells() > 0)
1075  {
1076  AssertDimension(static_cast<int>(level),
1077  data_->get_cell_iterator(0,0,j)->level());
1078  }
1079 
1080  // setup edge_constrained indices
1081  std::vector<types::global_dof_index> interface_indices;
1082  mg_constrained_dofs[j].get_refinement_edge_indices(level)
1083  .fill_index_vector(interface_indices);
1084  edge_constrained_indices[j].clear();
1085  edge_constrained_indices[j].reserve(interface_indices.size());
1086  edge_constrained_values[j].resize(interface_indices.size());
1087  const IndexSet &locally_owned = data->get_dof_handler(selected_rows[j]).locally_owned_mg_dofs(level);
1088  for (unsigned int i=0; i<interface_indices.size(); ++i)
1089  if (locally_owned.is_element(interface_indices[i]))
1090  edge_constrained_indices[j].push_back
1091  (locally_owned.index_within_set(interface_indices[i]));
1092  have_interface_matrices |=
1093  Utilities::MPI::max(static_cast<unsigned int>(edge_constrained_indices[j].size()),
1094  data->get_vector_partitioner()->get_mpi_communicator()) > 0;
1095  }
1096  }
1097 
1098 
1099 
1100  template <int dim, typename VectorType>
1101  void
1103  {
1104  for (unsigned int j=0; j<n_blocks(dst); ++j)
1105  {
1106  const std::vector<unsigned int> &
1107  constrained_dofs = data->get_constrained_dofs(selected_rows[j]);
1108  for (unsigned int i=0; i<constrained_dofs.size(); ++i)
1109  subblock(dst,j).local_element(constrained_dofs[i]) = 1.;
1110  for (unsigned int i=0; i<edge_constrained_indices[j].size(); ++i)
1111  subblock(dst,j).local_element(edge_constrained_indices[j][i]) = 1.;
1112  }
1113  }
1114 
1115 
1116 
1117  template <int dim, typename VectorType>
1118  void
1120  const VectorType &src) const
1121  {
1122  typedef typename Base<dim,VectorType>::value_type Number;
1123  dst = Number(0.);
1124  vmult_add (dst, src);
1125  }
1126 
1127 
1128 
1129  template <int dim, typename VectorType>
1130  void
1132  const VectorType &src) const
1133  {
1134  mult_add (dst, src, false);
1135  }
1136 
1137 
1138 
1139  template <int dim, typename VectorType>
1140  void
1142  const VectorType &src) const
1143  {
1144  mult_add (dst, src, true);
1145  }
1146 
1147 
1148 
1149  template <int dim, typename VectorType>
1150  void
1152  const bool is_row) const
1153  {
1154  typedef typename Base<dim,VectorType>::value_type Number;
1155  for (unsigned int i=0; i<n_blocks(src); ++i)
1156  {
1157  const unsigned int mf_component = is_row ? selected_rows[i] : selected_columns[i];
1158  // If both vectors use the same partitioner -> done
1159  if (subblock(src,i).get_partitioner().get() ==
1160  data->get_dof_info(mf_component).vector_partitioner.get())
1161  continue;
1162 
1163  // If not, assert that the local ranges are the same and reset to the
1164  // current partitioner
1165  Assert(subblock(src,i).get_partitioner()->local_size() ==
1166  data->get_dof_info(mf_component).vector_partitioner->local_size(),
1167  ExcMessage("The vector passed to the vmult() function does not have "
1168  "the correct size for compatibility with MatrixFree."));
1169 
1170  // copy the vector content to a temporary vector so that it does not get
1171  // lost
1172  VectorView<Number> view_src_in(subblock(src,i).local_size(),
1173  subblock(src,i).begin());
1174  const Vector<Number> copy_vec = view_src_in;
1175  subblock(const_cast<VectorType &>(src),i).
1176  reinit(data->get_dof_info(mf_component).vector_partitioner);
1177  VectorView<Number> view_src_out(subblock(src,i).local_size(),
1178  subblock(src,i).begin());
1179  static_cast<Vector<Number>&>(view_src_out) = copy_vec;
1180  }
1181  }
1182 
1183 
1184 
1185  template <int dim, typename VectorType>
1186  void
1188  const VectorType &src,
1189  const bool transpose) const
1190  {
1191  AssertDimension(dst.size(), src.size());
1192  AssertDimension(n_blocks(dst), n_blocks(src));
1193  AssertDimension(n_blocks(dst), selected_rows.size());
1194  typedef typename Base<dim,VectorType>::value_type Number;
1195  adjust_ghost_range_if_necessary(src, false);
1196  adjust_ghost_range_if_necessary(dst, true);
1197 
1198  // set zero Dirichlet values on the input vector (and remember the src and
1199  // dst values because we need to reset them at the end)
1200  for (unsigned int j=0; j<n_blocks(dst); ++j)
1201  {
1202  for (unsigned int i=0; i<edge_constrained_indices[j].size(); ++i)
1203  {
1204  edge_constrained_values[j][i] =
1205  std::pair<Number,Number>
1206  (subblock(src,j).local_element(edge_constrained_indices[j][i]),
1207  subblock(dst,j).local_element(edge_constrained_indices[j][i]));
1208  subblock(const_cast<VectorType &>(src),j).local_element(edge_constrained_indices[j][i]) = 0.;
1209  }
1210  }
1211 
1212  if (transpose)
1213  Tapply_add(dst,src);
1214  else
1215  apply_add(dst,src);
1216 
1217  for (unsigned int j=0; j<n_blocks(dst); ++j)
1218  {
1219  const std::vector<unsigned int> &
1220  constrained_dofs = data->get_constrained_dofs(selected_rows[j]);
1221  for (unsigned int i=0; i<constrained_dofs.size(); ++i)
1222  subblock(dst,j).local_element(constrained_dofs[i]) +=
1223  subblock(src,j).local_element(constrained_dofs[i]);
1224  }
1225 
1226  // reset edge constrained values, multiply by unit matrix and add into
1227  // destination
1228  for (unsigned int j=0; j<n_blocks(dst); ++j)
1229  {
1230  for (unsigned int i=0; i<edge_constrained_indices[j].size(); ++i)
1231  {
1232  subblock(const_cast<VectorType &>(src),j).local_element
1233  (edge_constrained_indices[j][i])
1234  = edge_constrained_values[j][i].first;
1235  subblock(dst,j).local_element(edge_constrained_indices[j][i])
1236  = edge_constrained_values[j][i].second +
1237  edge_constrained_values[j][i].first;
1238  }
1239  }
1240  }
1241 
1242 
1243 
1244  template <int dim, typename VectorType>
1245  void
1247  vmult_interface_down(VectorType &dst,
1248  const VectorType &src) const
1249  {
1250  typedef typename Base<dim,VectorType>::value_type Number;
1251  AssertDimension(dst.size(), src.size());
1252  adjust_ghost_range_if_necessary(src, false);
1253  adjust_ghost_range_if_necessary(dst, true);
1254 
1255  dst = Number(0.);
1256 
1257  if (!have_interface_matrices)
1258  return;
1259 
1260  // set zero Dirichlet values on the input vector (and remember the src and
1261  // dst values because we need to reset them at the end)
1262  for (unsigned int j=0; j<n_blocks(dst); ++j)
1263  for (unsigned int i=0; i<edge_constrained_indices[j].size(); ++i)
1264  {
1265  edge_constrained_values[j][i] =
1266  std::pair<Number,Number>
1267  (subblock(src,j).local_element(edge_constrained_indices[j][i]),
1268  subblock(dst,j).local_element(edge_constrained_indices[j][i]));
1269  subblock(const_cast<VectorType &>(src),j)
1270  .local_element(edge_constrained_indices[j][i]) = 0.;
1271  }
1272 
1273  apply_add(dst,src);
1274 
1275  for (unsigned int j=0; j<n_blocks(dst); ++j)
1276  {
1277  unsigned int c=0;
1278  for (unsigned int i=0; i<edge_constrained_indices[j].size(); ++i)
1279  {
1280  for ( ; c<edge_constrained_indices[j][i]; ++c)
1281  subblock(dst,j).local_element(c) = 0.;
1282  ++c;
1283 
1284  // reset the src values
1285  subblock(const_cast<VectorType &>(src),j)
1286  .local_element(edge_constrained_indices[j][i])
1287  = edge_constrained_values[j][i].first;
1288  }
1289  for ( ; c<dst.local_size(); ++c)
1290  subblock(dst,j).local_element(c) = 0.;
1291  }
1292  }
1293 
1294 
1295 
1296  template <int dim, typename VectorType>
1297  void
1299  vmult_interface_up(VectorType &dst,
1300  const VectorType &src) const
1301  {
1302  typedef typename Base<dim,VectorType>::value_type Number;
1303  AssertDimension(dst.size(), src.size());
1304  adjust_ghost_range_if_necessary(src, false);
1305  adjust_ghost_range_if_necessary(dst, true);
1306 
1307  dst = Number(0.);
1308 
1309  if (!have_interface_matrices)
1310  return;
1311 
1312  VectorType src_cpy (src);
1313  for (unsigned int j=0; j<n_blocks(dst); ++j)
1314  {
1315  unsigned int c=0;
1316  for (unsigned int i=0; i<edge_constrained_indices[j].size(); ++i)
1317  {
1318  for ( ; c<edge_constrained_indices[j][i]; ++c)
1319  subblock(src_cpy,j).local_element(c) = 0.;
1320  ++c;
1321  }
1322  for ( ; c<subblock(src_cpy,j).local_size(); ++c)
1323  subblock(src_cpy,j).local_element(c) = 0.;
1324  }
1325 
1326  apply_add(dst,src_cpy);
1327 
1328  for (unsigned int j=0; j<n_blocks(dst); ++j)
1329  for (unsigned int i=0; i<edge_constrained_indices[j].size(); ++i)
1330  subblock(dst,j).local_element(edge_constrained_indices[j][i]) = 0.;
1331  }
1332 
1333 
1334 
1335  template <int dim, typename VectorType>
1336  void
1338  const VectorType &src) const
1339  {
1340  typedef typename Base<dim,VectorType>::value_type Number;
1341  dst = Number(0.);
1342  Tvmult_add (dst,src);
1343  }
1344 
1345 
1346 
1347  template <int dim, typename VectorType>
1348  std::size_t
1350  {
1351  return inverse_diagonal_entries.get() != NULL ? inverse_diagonal_entries->memory_consumption() : sizeof(*this);
1352  }
1353 
1354 
1355 
1356  template <int dim, typename VectorType>
1357  std_cxx11::shared_ptr<const MatrixFree<dim,typename Base<dim,VectorType>::value_type> >
1359  {
1360  return data;
1361  }
1362 
1363 
1364 
1365  template <int dim, typename VectorType>
1366  const std_cxx11::shared_ptr<DiagonalMatrix<VectorType> > &
1368  {
1369  Assert(inverse_diagonal_entries.get() != NULL &&
1370  inverse_diagonal_entries->m() > 0, ExcNotInitialized());
1371  return inverse_diagonal_entries;
1372  }
1373 
1374 
1375 
1376  template <int dim, typename VectorType>
1377  void
1379  const VectorType &src) const
1380  {
1381  apply_add(dst,src);
1382  }
1383 
1384 
1385 
1386  template <int dim, typename VectorType>
1387  void
1389  const VectorType &src,
1390  const typename Base<dim,VectorType>::value_type omega) const
1391  {
1392  Assert(inverse_diagonal_entries.get() &&
1393  inverse_diagonal_entries->m() > 0, ExcNotInitialized());
1394  inverse_diagonal_entries->vmult(dst,src);
1395  dst *= omega;
1396  }
1397 
1398 
1399 
1400  //------------------------- MGInterfaceOperator ------------------------------
1401 
1402  template <typename OperatorType>
1404  :
1405  Subscriptor(),
1406  mf_base_operator(NULL)
1407  {
1408  }
1409 
1410 
1411 
1412  template <typename OperatorType>
1413  void
1415  {
1416  mf_base_operator = NULL;
1417  }
1418 
1419 
1420 
1421  template <typename OperatorType>
1422  void
1423  MGInterfaceOperator<OperatorType>::initialize (const OperatorType &operator_in)
1424  {
1425  mf_base_operator = &operator_in;
1426  }
1427 
1428 
1429 
1430  template <typename OperatorType>
1431  template <typename VectorType>
1432  void
1434  const VectorType &src) const
1435  {
1436 #ifdef DEAL_II_WITH_CXX11
1437 #ifndef DEAL_II_MSVC
1438  static_assert (std::is_same<typename VectorType::value_type,value_type>::value,
1439  "The vector type must be based on the same value type as this"
1440  "operator");
1441 #endif
1442 #endif
1443 
1444  Assert(mf_base_operator != NULL,
1445  ExcNotInitialized());
1446 
1447  mf_base_operator->vmult_interface_down(dst, src);
1448  }
1449 
1450 
1451 
1452  template <typename OperatorType>
1453  template <typename VectorType>
1454  void
1456  const VectorType &src) const
1457  {
1458 #ifdef DEAL_II_WITH_CXX11
1459 #ifndef DEAL_II_MSVC
1460  static_assert (std::is_same<typename VectorType::value_type,value_type>::value,
1461  "The vector type must be based on the same value type as this"
1462  "operator");
1463 #endif
1464 #endif
1465 
1466  Assert(mf_base_operator != NULL,
1467  ExcNotInitialized());
1468 
1469  mf_base_operator->vmult_interface_up(dst, src);
1470  }
1471 
1472 
1473 
1474  //-----------------------------MassOperator----------------------------------
1475 
1476  template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType>
1479  :
1480  Base<dim, VectorType>()
1481  {}
1482 
1483 
1484 
1485  template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType>
1486  void
1489  {
1490  typedef typename Base<dim,VectorType>::value_type Number;
1492 
1493  this->inverse_diagonal_entries.
1494  reset(new DiagonalMatrix<VectorType>());
1495  VectorType &inverse_diagonal_vector = this->inverse_diagonal_entries->get_vector();
1496  VectorType ones;
1497  this->initialize_dof_vector(ones);
1498  this->initialize_dof_vector(inverse_diagonal_vector);
1499  ones = Number(1.);
1500  apply_add(inverse_diagonal_vector, ones);
1501 
1502  this->set_constrained_entries_to_one(inverse_diagonal_vector);
1503 
1504  const unsigned int local_size = inverse_diagonal_vector.local_size();
1505  for (unsigned int i=0; i<local_size; ++i)
1506  inverse_diagonal_vector.local_element(i)
1507  =1./inverse_diagonal_vector.local_element(i);
1508 
1509  inverse_diagonal_vector.update_ghost_values();
1510  }
1511 
1512 
1513 
1514  template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType>
1515  void
1517  apply_add (VectorType &dst,
1518  const VectorType &src) const
1519  {
1521  this, dst, src);
1522  }
1523 
1524 
1525 
1526  template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType>
1527  void
1530  VectorType &dst,
1531  const VectorType &src,
1532  const std::pair<unsigned int,unsigned int> &cell_range) const
1533  {
1534  typedef typename Base<dim,VectorType>::value_type Number;
1536  for (unsigned int cell=cell_range.first; cell<cell_range.second; ++cell)
1537  {
1538  phi.reinit (cell);
1539  phi.read_dof_values(src);
1540  phi.evaluate (true,false,false);
1541  for (unsigned int q=0; q<phi.n_q_points; ++q)
1542  phi.submit_value (phi.get_value(q), q);
1543  phi.integrate (true,false);
1544  phi.distribute_local_to_global (dst);
1545  }
1546  }
1547 
1548 
1549  //-----------------------------LaplaceOperator----------------------------------
1550 
1551  template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType>
1554  :
1555  Base<dim, VectorType>()
1556  {
1557  }
1558 
1559 
1560 
1561  template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType>
1562  void
1565  {
1567  scalar_coefficient.reset();
1568  }
1569 
1570 
1571 
1572  template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType>
1573  void
1575  set_coefficient(const std_cxx11::shared_ptr<Table<2, VectorizedArray<typename Base<dim,VectorType>::value_type> > > &scalar_coefficient_ )
1576  {
1577  scalar_coefficient = scalar_coefficient_;
1578  }
1579 
1580 
1581 
1582  template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType>
1583  std_cxx11::shared_ptr< Table<2, VectorizedArray< typename LaplaceOperator<dim, fe_degree, n_q_points_1d, n_components, VectorType>::value_type> > >
1586  {
1587  Assert (scalar_coefficient.get(),
1588  ExcNotInitialized());
1589  return scalar_coefficient;
1590  }
1591 
1592 
1593 
1594  template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType>
1595  void
1598  {
1599  typedef typename Base<dim,VectorType>::value_type Number;
1601 
1602  unsigned int dummy = 0;
1603  this->inverse_diagonal_entries.
1604  reset(new DiagonalMatrix<VectorType>());
1605  VectorType &inverse_diagonal_vector = this->inverse_diagonal_entries->get_vector();
1606  this->initialize_dof_vector(inverse_diagonal_vector);
1607 
1608  this->data->cell_loop (&LaplaceOperator::local_diagonal_cell,
1609  this, inverse_diagonal_vector, dummy);
1610  this->set_constrained_entries_to_one(inverse_diagonal_vector);
1611 
1612  for (unsigned int i=0; i<inverse_diagonal_vector.local_size(); ++i)
1613  if (std::abs(inverse_diagonal_vector.local_element(i)) > std::sqrt(std::numeric_limits<Number>::epsilon()))
1614  inverse_diagonal_vector.local_element(i) = 1./inverse_diagonal_vector.local_element(i);
1615  else
1616  inverse_diagonal_vector.local_element(i) = 1.;
1617 
1618  inverse_diagonal_vector.update_ghost_values();
1619  }
1620 
1621 
1622 
1623  template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType>
1624  void
1626  apply_add (VectorType &dst,
1627  const VectorType &src) const
1628  {
1630  this, dst, src);
1631  }
1632 
1633  namespace
1634  {
1635  template<typename Number>
1636  bool
1637  non_negative(const VectorizedArray<Number> &n)
1638  {
1639  for (unsigned int v=0; v<VectorizedArray<Number>::n_array_elements; ++v)
1640  if (n[v] < 0.)
1641  return false;
1642 
1643  return true;
1644  }
1645  }
1646 
1647 
1648 
1649  template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType>
1650  void
1652  do_operation_on_cell(FEEvaluation<dim,fe_degree,n_q_points_1d,n_components,typename Base<dim,VectorType>::value_type> &phi,
1653  const unsigned int cell) const
1654  {
1655  phi.evaluate (false,true,false);
1656  if (scalar_coefficient.get())
1657  {
1658  for (unsigned int q=0; q<phi.n_q_points; ++q)
1659  {
1660  Assert (non_negative((*scalar_coefficient)(cell,q)),
1661  ExcMessage("Coefficient must be non-negative"));
1662  phi.submit_gradient ((*scalar_coefficient)(cell,q)*phi.get_gradient(q), q);
1663  }
1664  }
1665  else
1666  {
1667  for (unsigned int q=0; q<phi.n_q_points; ++q)
1668  {
1669  phi.submit_gradient (phi.get_gradient(q), q);
1670  }
1671  }
1672  phi.integrate (false,true);
1673  }
1674 
1675 
1676 
1677 
1678  template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType>
1679  void
1682  VectorType &dst,
1683  const VectorType &src,
1684  const std::pair<unsigned int,unsigned int> &cell_range) const
1685  {
1686  typedef typename Base<dim,VectorType>::value_type Number;
1687  FEEvaluation<dim,fe_degree,n_q_points_1d,n_components,Number> phi (data, this->selected_rows[0]);
1688  for (unsigned int cell=cell_range.first; cell<cell_range.second; ++cell)
1689  {
1690  phi.reinit (cell);
1691  phi.read_dof_values(src);
1692  do_operation_on_cell(phi,cell);
1693  phi.distribute_local_to_global (dst);
1694  }
1695  }
1696 
1697 
1698  template <int dim, int fe_degree, int n_q_points_1d, int n_components, typename VectorType>
1699  void
1702  VectorType &dst,
1703  const unsigned int &,
1704  const std::pair<unsigned int,unsigned int> &cell_range) const
1705  {
1706  typedef typename Base<dim,VectorType>::value_type Number;
1707  FEEvaluation<dim,fe_degree,n_q_points_1d,n_components,Number> phi (data, this->selected_rows[0]);
1708  for (unsigned int cell=cell_range.first; cell<cell_range.second; ++cell)
1709  {
1710  phi.reinit (cell);
1711  VectorizedArray<Number> local_diagonal_vector[phi.tensor_dofs_per_cell];
1712  for (unsigned int i=0; i<phi.dofs_per_cell; ++i)
1713  {
1714  for (unsigned int j=0; j<phi.dofs_per_cell; ++j)
1716  phi.begin_dof_values()[i] = 1.;
1717  do_operation_on_cell(phi,cell);
1718  local_diagonal_vector[i] = phi.begin_dof_values()[i];
1719  }
1720  for (unsigned int i=0; i<phi.tensor_dofs_per_cell; ++i)
1721  phi.begin_dof_values()[i] = local_diagonal_vector[i];
1722  phi.distribute_local_to_global (dst);
1723  }
1724  }
1725 
1726 
1727 } // end of namespace MatrixFreeOperators
1728 
1729 
1730 DEAL_II_NAMESPACE_CLOSE
1731 
1732 #endif
std_cxx11::shared_ptr< Table< 2, VectorizedArray< value_type > > > get_coefficient()
Definition: operators.h:1585
size_type m() const
void initialize(std_cxx11::shared_ptr< const MatrixFree< dim, value_type > > data, const std::vector< unsigned int > &selected_row_blocks=std::vector< unsigned int >(), const std::vector< unsigned int > &selected_column_blocks=std::vector< unsigned int >())
void do_operation_on_cell(FEEvaluation< dim, fe_degree, n_q_points_1d, n_components, value_type > &phi, const unsigned int cell) const
Definition: operators.h:1652
const std_cxx11::shared_ptr< DiagonalMatrix< VectorType > > & get_matrix_diagonal_inverse() const
Definition: operators.h:1367
static const unsigned int invalid_unsigned_int
Definition: types.h:170
void vmult_interface_down(VectorType &dst, const VectorType &src) const
Definition: operators.h:1247
OperatorType::value_type value_type
Definition: operators.h:471
std_cxx11::shared_ptr< DiagonalMatrix< VectorType > > inverse_diagonal_entries
Definition: operators.h:374
std::vector< unsigned int > selected_rows
Definition: operators.h:380
#define AssertDimension(dim1, dim2)
Definition: exceptions.h:1146
void reinit(const unsigned int cell)
virtual void apply_add(VectorType &dst, const VectorType &src) const
Definition: operators.h:1626
void local_apply_cell(const MatrixFree< dim, value_type > &data, VectorType &dst, const VectorType &src, const std::pair< unsigned int, unsigned int > &cell_range) const
Definition: operators.h:1681
void apply(const AlignedVector< VectorizedArray< Number > > &inverse_coefficient, const unsigned int n_actual_components, const VectorizedArray< Number > *in_array, VectorizedArray< Number > *out_array) const
Definition: operators.h:823
void gauss_jordan()
const unsigned int dofs_per_cell
void evaluate(const bool evaluate_values, const bool evaluate_gradients, const bool evaluate_hessians=false)
void read_dof_values(const VectorType &src, const unsigned int first_index=0)
void integrate(const bool integrate_values, const bool integrate_gradients)
void precondition_Jacobi(VectorType &dst, const VectorType &src, const value_type omega) const
Definition: operators.h:1388
#define AssertIndexRange(index, range)
Definition: exceptions.h:1170
const FEEvaluationBase< dim, n_components, Number > & fe_eval
Definition: operators.h:564
void mult_add(VectorType &dst, const VectorType &src, const bool transpose) const
Definition: operators.h:1187
static ::ExceptionBase & ExcNotInitialized()
void vmult_interface_up(VectorType &dst, const VectorType &src) const
Definition: operators.h:1299
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
void set_constrained_entries_to_one(VectorType &dst) const
Definition: operators.h:1102
std_cxx11::shared_ptr< const MatrixFree< dim, value_type > > get_matrix_free() const
Definition: operators.h:1358
const VectorizedArray< Number > * begin_dof_values() const
Base< dim, VectorType >::value_type value_type
Definition: operators.h:590
std::vector< unsigned int > selected_columns
Definition: operators.h:386
value_type el(const unsigned int row, const unsigned int col) const
Definition: operators.h:933
VectorType::value_type value_type
Definition: operators.h:170
static ::ExceptionBase & ExcMessage(std::string arg1)
const unsigned int n_q_points
size_type n() const
Definition: operators.h:909
size_type n() const
virtual void Tapply_add(VectorType &dst, const VectorType &src) const
Definition: operators.h:1378
void set_coefficient(const std_cxx11::shared_ptr< Table< 2, VectorizedArray< value_type > > > &scalar_coefficient)
Definition: operators.h:1575
void Tvmult(VectorType &dst, const VectorType &src) const
Definition: operators.h:1337
#define Assert(cond, exc)
Definition: exceptions.h:313
void local_diagonal_cell(const MatrixFree< dim, value_type > &data, VectorType &dst, const unsigned int &, const std::pair< unsigned int, unsigned int > &cell_range) const
Definition: operators.h:1701
size_type index_within_set(const size_type global_index) const
Definition: index_set.h:1653
virtual void compute_diagonal()=0
std_cxx11::shared_ptr< const MatrixFree< dim, value_type > > data
Definition: operators.h:368
VectorType::size_type size_type
Definition: operators.h:175
std::vector< std::vector< unsigned int > > edge_constrained_indices
Definition: operators.h:393
AlignedVector< VectorizedArray< Number > > inverse_shape
Definition: operators.h:569
void Tvmult_add(VectorType &dst, const VectorType &src) const
Definition: operators.h:1141
virtual void apply_add(VectorType &dst, const VectorType &src) const
Definition: operators.h:1517
Base< dim, VectorType >::size_type size_type
Definition: operators.h:650
void vmult(VectorType &dst, const VectorType &src) const
Definition: operators.h:1119
size_type m() const
Definition: operators.h:895
void distribute_local_to_global(VectorType &dst, const unsigned int first_index=0) const
virtual void apply_add(VectorType &dst, const VectorType &src) const =0
void vmult(VectorType &dst, const VectorType &src) const
Definition: operators.h:1433
void submit_value(const value_type val_in, const unsigned int q_point)
CellwiseInverseMassMatrix(const FEEvaluationBase< dim, n_components, Number > &fe_eval)
Definition: operators.h:764
std::vector< std::vector< std::pair< value_type, value_type > > > edge_constrained_values
Definition: operators.h:399
void local_apply_cell(const MatrixFree< dim, value_type > &data, VectorType &dst, const VectorType &src, const std::pair< unsigned int, unsigned int > &cell_range) const
Definition: operators.h:1529
virtual std::size_t memory_consumption() const
Definition: operators.h:1349
void fill_inverse_JxW_values(AlignedVector< VectorizedArray< Number > > &inverse_jxw) const
Definition: operators.h:794
const internal::MatrixFreeFunctions::ShapeInfo< Number > & get_shape_info() const
static ::ExceptionBase & ExcNotImplemented()
std_cxx11::shared_ptr< Table< 2, VectorizedArray< value_type > > > scalar_coefficient
Definition: operators.h:754
Base< dim, VectorType >::value_type value_type
Definition: operators.h:645
Definition: table.h:33
bool is_element(const size_type index) const
Definition: index_set.h:1489
void adjust_ghost_range_if_necessary(const VectorType &vec, const bool is_row) const
Definition: operators.h:1151
void vmult_add(VectorType &dst, const VectorType &src) const
Definition: operators.h:1131
void initialize(const OperatorType &operator_in)
Definition: operators.h:1423
virtual void clear()
Definition: operators.h:923
SmartPointer< const OperatorType > mf_base_operator
Definition: operators.h:506
value_type get_value(const unsigned int q_point) const
Base< dim, VectorType >::size_type size_type
Definition: operators.h:595
void Tvmult(VectorType &dst, const VectorType &src) const
Definition: operators.h:1455
void initialize_dof_vector(VectorType &vec) const
Definition: operators.h:947
T max(const T &t, const MPI_Comm &mpi_communicator)
static ::ExceptionBase & ExcInternalError()