17 #include <deal.II/base/quadrature_lib.h> 18 #include <deal.II/base/qprojector.h> 19 #include <deal.II/base/thread_management.h> 20 #include <deal.II/base/utilities.h> 21 #include <deal.II/lac/vector.h> 22 #include <deal.II/lac/la_vector.h> 23 #include <deal.II/lac/block_vector.h> 24 #include <deal.II/lac/la_parallel_vector.h> 25 #include <deal.II/lac/la_parallel_block_vector.h> 26 #include <deal.II/lac/petsc_parallel_vector.h> 27 #include <deal.II/lac/petsc_block_vector.h> 28 #include <deal.II/lac/petsc_parallel_block_vector.h> 29 #include <deal.II/lac/trilinos_vector.h> 30 #include <deal.II/lac/trilinos_block_vector.h> 31 #include <deal.II/lac/constraint_matrix.h> 32 #include <deal.II/grid/tria.h> 33 #include <deal.II/grid/tria_iterator.h> 34 #include <deal.II/grid/grid_generator.h> 35 #include <deal.II/fe/fe_tools.h> 36 #include <deal.II/fe/fe.h> 37 #include <deal.II/fe/fe_values.h> 38 #include <deal.II/dofs/dof_handler.h> 39 #include <deal.II/dofs/dof_accessor.h> 40 #include <deal.II/dofs/dof_tools.h> 41 #include <deal.II/hp/dof_handler.h> 43 #include <deal.II/base/std_cxx11/shared_ptr.h> 45 #include <deal.II/base/index_set.h> 50 DEAL_II_NAMESPACE_OPEN
54 template <
int dim,
int spacedim,
55 template <
int,
int>
class DoFHandlerType1,
56 template <
int,
int>
class DoFHandlerType2,
57 class InVector,
class OutVector>
61 const DoFHandlerType2<dim, spacedim> &dof2,
71 template <
int dim,
int spacedim,
72 template <
int,
int>
class DoFHandlerType1,
73 template <
int,
int>
class DoFHandlerType2,
74 class InVector,
class OutVector>
78 const DoFHandlerType2<dim, spacedim> &dof2,
84 Assert(u1.size()==dof1.n_dofs(),
86 Assert(u2.size()==dof2.n_dofs(),
90 const IndexSet u2_elements = u2.locally_owned_elements();
92 const IndexSet &dof1_local_dofs = dof1.locally_owned_dofs();
93 const IndexSet &dof2_local_dofs = dof2.locally_owned_dofs();
94 const IndexSet u1_elements = u1.locally_owned_elements();
95 Assert(u1_elements == dof1_local_dofs,
96 ExcMessage(
"The provided vector and DoF handler should have the same" 98 Assert(u2_elements == dof2_local_dofs,
99 ExcMessage(
"The provided vector and DoF handler should have the same" 114 std::map<const FiniteElement<dim,spacedim> *,
115 std::map<const FiniteElement<dim,spacedim> *,
116 std_cxx11::shared_ptr<FullMatrix<double> > > >
117 interpolation_matrices;
119 typename DoFHandlerType1<dim,spacedim>::active_cell_iterator cell1 = dof1.begin_active(),
121 typename DoFHandlerType2<dim,spacedim>::active_cell_iterator cell2 = dof2.begin_active(),
125 std::vector<types::global_dof_index> dofs;
126 dofs.reserve (DoFTools::max_dofs_per_cell (dof2));
129 OutVector touch_count(u2);
137 dof1.get_triangulation().locally_owned_subdomain();
139 for (; cell1!=endc1; ++cell1, ++cell2)
140 if ((cell1->subdomain_id() == subdomain_id)
144 Assert(cell1->get_fe().n_components() == cell2->get_fe().n_components(),
146 cell2->get_fe().n_components()));
155 const bool hanging_nodes_not_allowed
156 = ((cell2->get_fe().dofs_per_vertex != 0) &&
159 if (hanging_nodes_not_allowed)
160 for (
unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
161 Assert (cell1->at_boundary(face) ||
162 cell1->neighbor(face)->level() == cell1->level(),
166 const unsigned int dofs_per_cell1 = cell1->get_fe().dofs_per_cell;
167 const unsigned int dofs_per_cell2 = cell2->get_fe().dofs_per_cell;
168 u1_local.
reinit (dofs_per_cell1);
169 u2_local.
reinit (dofs_per_cell2);
175 if (interpolation_matrices[&cell1->get_fe()][&cell2->get_fe()].get() == 0)
177 std_cxx11::shared_ptr<FullMatrix<double> >
180 interpolation_matrices[&cell1->get_fe()][&cell2->get_fe()]
181 = interpolation_matrix;
185 *interpolation_matrix);
188 cell1->get_dof_values(u1, u1_local);
189 interpolation_matrices[&cell1->get_fe()][&cell2->get_fe()]
190 ->vmult(u2_local, u1_local);
192 dofs.resize (dofs_per_cell2);
193 cell2->get_dof_indices(dofs);
195 for (
unsigned int i=0; i<dofs_per_cell2; ++i)
201 u2(dofs[i])+=u2_local(i);
202 touch_count(dofs[i]) += 1;
216 IndexSet locally_owned_dofs = dof2.locally_owned_dofs();
227 Assert(static_cast<typename OutVector::value_type>(touch_count(i)) !=
typename OutVector::value_type(0),
229 u2(i) /= touch_count(i);
241 template <
int,
int>
class DoFHandlerType,
242 class InVector,
class OutVector,
int spacedim>
247 OutVector &u1_interpolated)
251 Assert(u1.size() == dof1.n_dofs(),
253 Assert(u1_interpolated.size() == dof1.n_dofs(),
257 const IndexSet &dof1_local_dofs = dof1.locally_owned_dofs();
258 const IndexSet u1_elements = u1.locally_owned_elements();
259 const IndexSet u1_interpolated_elements = u1_interpolated.locally_owned_elements();
260 Assert(u1_elements == dof1_local_dofs,
261 ExcMessage(
"The provided vector and DoF handler should have the same" 263 Assert(u1_interpolated_elements == dof1_local_dofs,
264 ExcMessage(
"The provided vector and DoF handler should have the same" 272 dof1.get_triangulation().locally_owned_subdomain();
274 typename DoFHandlerType<dim, spacedim>::active_cell_iterator cell = dof1.begin_active(),
280 std::map<const FiniteElement<dim> *,
281 std_cxx11::shared_ptr<FullMatrix<double> > > interpolation_matrices;
283 for (; cell!=endc; ++cell)
284 if ((cell->subdomain_id() == subdomain_id)
295 const bool hanging_nodes_not_allowed=
298 if (hanging_nodes_not_allowed)
299 for (
unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
300 Assert (cell->at_boundary(face) ||
301 cell->neighbor(face)->level() == cell->level(),
304 const unsigned int dofs_per_cell1 = cell->get_fe().dofs_per_cell;
308 if (interpolation_matrices[&cell->get_fe()] == 0)
310 interpolation_matrices[&cell->get_fe()] =
311 std_cxx11::shared_ptr<FullMatrix<double> >
314 *interpolation_matrices[&cell->get_fe()]);
317 u1_local.
reinit (dofs_per_cell1);
318 u1_int_local.
reinit (dofs_per_cell1);
320 cell->get_dof_values(u1, u1_local);
321 interpolation_matrices[&cell->get_fe()]->vmult(u1_int_local, u1_local);
322 cell->set_dof_values(u1_int_local, u1_interpolated);
335 template <
int dim,
int spacedim,
class InVector>
341 InVector &u1_interpolated)
345 interpolate(dof2, u2, dof1, constraints1, u1_interpolated);
349 #ifdef DEAL_II_WITH_PETSC 350 template <
int dim,
int spacedim>
362 IndexSet dof2_locally_relevant_dofs;
364 dof2_locally_relevant_dofs);
370 dof2_locally_relevant_dofs,
373 interpolate(dof2, u2, dof1, constraints1, u1_interpolated);
378 #ifdef DEAL_II_WITH_TRILINOS 379 template <
int dim,
int spacedim>
391 IndexSet dof2_locally_relevant_dofs;
393 dof2_locally_relevant_dofs);
399 dof2_locally_relevant_dofs,
402 interpolate(dof2, u2, dof1, constraints1, u1_interpolated);
407 template <
int dim,
int spacedim,
typename Number>
416 IndexSet dof2_locally_relevant_dofs;
418 dof2_locally_relevant_dofs);
421 u2 (dof2_locally_owned_dofs,
422 dof2_locally_relevant_dofs,
426 u2.update_ghost_values ();
427 interpolate(dof2, u2, dof1, constraints1, u1_interpolated);
434 template <
int dim,
class InVector,
class OutVector,
int spacedim>
440 OutVector &u1_interpolated)
444 if (dof1.
get_fe().dofs_per_vertex==0 && dof2.
get_fe().dofs_per_vertex==0
458 internal::back_interpolate(dof1, constraints1, u1, dof2, constraints2,
465 template <
int dim,
class InVector,
class OutVector,
int spacedim>
469 OutVector &u1_difference)
479 const IndexSet u1_elements = u1.locally_owned_elements();
480 const IndexSet u1_difference_elements = u1_difference.locally_owned_elements();
481 Assert(u1_elements == dof1_local_dofs,
482 ExcMessage(
"The provided vector and DoF handler should have the same" 484 Assert(u1_difference_elements == dof1_local_dofs,
485 ExcMessage(
"The provided vector and DoF handler should have the same" 496 const bool hanging_nodes_not_allowed=
499 const unsigned int dofs_per_cell=dof1.
get_fe().dofs_per_cell;
514 for (; cell!=endc; ++cell)
515 if ((cell->subdomain_id() == subdomain_id)
519 if (hanging_nodes_not_allowed)
520 for (
unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
521 Assert (cell->at_boundary(face) ||
522 cell->neighbor(face)->level() == cell->level(),
525 cell->get_dof_values(u1, u1_local);
526 difference_matrix.
vmult(u1_diff_local, u1_local);
527 cell->set_dof_values(u1_diff_local, u1_difference);
542 template <
int dim,
class InVector,
class OutVector,
int spacedim>
548 OutVector &u1_difference)
550 back_interpolate(dof1, constraints1, u1, dof2, constraints2, u1_difference);
551 u1_difference.sadd(-1., 1., u1);
555 #ifdef DEAL_II_WITH_TRILINOS 556 template <
int dim,
int spacedim>
564 back_interpolate(dof1, constraints1, u1, dof2, constraints2, u1_difference);
572 u1_completely_distributed = u1;
574 u1_difference.
sadd(-1, u1_completely_distributed);
582 template <
int dim,
class InVector,
class OutVector,
int spacedim>
588 OutVector &u1_difference)
594 if (dof1.
get_fe().dofs_per_vertex==0 && dof2.
get_fe().dofs_per_vertex==0
599 internal::interpolation_difference(dof1, constraints1, u1, dof2, constraints2, u1_difference);
605 template <
int dim,
class InVector,
class OutVector,
int spacedim>
621 const unsigned int n1 = dof1.
get_fe().dofs_per_cell;
622 const unsigned int n2 = dof2.
get_fe().dofs_per_cell;
626 std::vector<types::global_dof_index> dofs(n2);
634 cell1->get_dof_values(u1, u1_local);
635 matrix.vmult(u2_local, u1_local);
636 cell2->get_dof_indices(dofs);
637 for (
unsigned int i=0; i<n2; ++i)
639 u2(dofs[i])+=u2_local(i);
651 #include "fe_tools_interpolate.inst" 656 DEAL_II_NAMESPACE_CLOSE
const types::subdomain_id invalid_subdomain_id
void vmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
static ::ExceptionBase & ExcHangingNodesNotAllowed(int arg1)
void sadd(const TrilinosScalar s, const VectorBase &V)
cell_iterator end() const
ActiveSelector::active_cell_iterator active_cell_iterator
const FiniteElement< dim, spacedim > & get_fe() const
active_cell_iterator begin_active(const unsigned int level=0) const
static ::ExceptionBase & ExcMessage(std::string arg1)
unsigned int global_dof_index
#define Assert(cond, exc)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
const Epetra_Map & vector_partitioner() const
size_type n_constraints() const
types::global_dof_index n_dofs() const
const MPI_Comm & get_mpi_communicator() const
const MPI_Comm & get_mpi_communicator() const
unsigned int subdomain_id
void distribute(VectorType &vec) const
unsigned int n_components() const
const MPI_Comm & get_mpi_communicator() const
virtual void reinit(const size_type N, const bool omit_zeroing_entries=false)
const Triangulation< dim, spacedim > & get_triangulation() const
const unsigned int dofs_per_vertex
bool is_element(const size_type index) const
const IndexSet & locally_owned_dofs() const
static ::ExceptionBase & ExcTriangulationMismatch()
static ::ExceptionBase & ExcInternalError()