Reference documentation for deal.II version 8.5.1
fe_tools_interpolate.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 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 #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>
42 
43 #include <deal.II/base/std_cxx11/shared_ptr.h>
44 
45 #include <deal.II/base/index_set.h>
46 
47 #include <iostream>
48 
49 
50 DEAL_II_NAMESPACE_OPEN
51 
52 namespace FETools
53 {
54  template <int dim, int spacedim,
55  template <int, int> class DoFHandlerType1,
56  template <int, int> class DoFHandlerType2,
57  class InVector, class OutVector>
58  void
59  interpolate(const DoFHandlerType1<dim, spacedim> &dof1,
60  const InVector &u1,
61  const DoFHandlerType2<dim, spacedim> &dof2,
62  OutVector &u2)
63  {
64  ConstraintMatrix dummy;
65  dummy.close();
66  interpolate(dof1, u1, dof2, dummy, u2);
67  }
68 
69 
70 
71  template <int dim, int spacedim,
72  template <int, int> class DoFHandlerType1,
73  template <int, int> class DoFHandlerType2,
74  class InVector, class OutVector>
75  void
76  interpolate (const DoFHandlerType1<dim, spacedim> &dof1,
77  const InVector &u1,
78  const DoFHandlerType2<dim, spacedim> &dof2,
79  const ConstraintMatrix &constraints,
80  OutVector &u2)
81  {
82  Assert(&dof1.get_triangulation()==&dof2.get_triangulation(), ExcTriangulationMismatch());
83 
84  Assert(u1.size()==dof1.n_dofs(),
85  ExcDimensionMismatch(u1.size(), dof1.n_dofs()));
86  Assert(u2.size()==dof2.n_dofs(),
87  ExcDimensionMismatch(u2.size(), dof2.n_dofs()));
88 
89 
90  const IndexSet u2_elements = u2.locally_owned_elements();
91 #ifdef DEBUG
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"
97  " index sets."));
98  Assert(u2_elements == dof2_local_dofs,
99  ExcMessage("The provided vector and DoF handler should have the same"
100  " index sets."));
101 #endif
102 
103  // allocate vectors at maximal
104  // size. will be reinited in inner
105  // cell, but Vector makes sure that
106  // this does not lead to
107  // reallocation of memory
108  Vector<typename OutVector::value_type> u1_local(DoFTools::max_dofs_per_cell(dof1));
109  Vector<typename OutVector::value_type> u2_local(DoFTools::max_dofs_per_cell(dof2));
110 
111  // have a map for interpolation
112  // matrices. shared_ptr make sure
113  // that memory is released again
114  std::map<const FiniteElement<dim,spacedim> *,
115  std::map<const FiniteElement<dim,spacedim> *,
116  std_cxx11::shared_ptr<FullMatrix<double> > > >
117  interpolation_matrices;
118 
119  typename DoFHandlerType1<dim,spacedim>::active_cell_iterator cell1 = dof1.begin_active(),
120  endc1 = dof1.end();
121  typename DoFHandlerType2<dim,spacedim>::active_cell_iterator cell2 = dof2.begin_active(),
122  endc2 = dof2.end();
123  (void)endc2;
124 
125  std::vector<types::global_dof_index> dofs;
126  dofs.reserve (DoFTools::max_dofs_per_cell (dof2));
127 
128  u2 = 0;
129  OutVector touch_count(u2);
130  touch_count = 0;
131 
132  // for distributed triangulations,
133  // we can only interpolate u1 on
134  // a cell, which this processor owns,
135  // so we have to know the subdomain_id
136  const types::subdomain_id subdomain_id =
137  dof1.get_triangulation().locally_owned_subdomain();
138 
139  for (; cell1!=endc1; ++cell1, ++cell2)
140  if ((cell1->subdomain_id() == subdomain_id)
141  ||
142  (subdomain_id == numbers::invalid_subdomain_id))
143  {
144  Assert(cell1->get_fe().n_components() == cell2->get_fe().n_components(),
145  ExcDimensionMismatch (cell1->get_fe().n_components(),
146  cell2->get_fe().n_components()));
147 
148  // for continuous elements on
149  // grids with hanging nodes we
150  // need hanging node
151  // constraints. Consequently,
152  // if there are no constraints
153  // then hanging nodes are not
154  // allowed.
155  const bool hanging_nodes_not_allowed
156  = ((cell2->get_fe().dofs_per_vertex != 0) &&
157  (constraints.n_constraints() == 0));
158 
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(),
164 
165 
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);
170 
171  // check if interpolation
172  // matrix for this particular
173  // pair of elements is already
174  // there
175  if (interpolation_matrices[&cell1->get_fe()][&cell2->get_fe()].get() == 0)
176  {
177  std_cxx11::shared_ptr<FullMatrix<double> >
178  interpolation_matrix (new FullMatrix<double> (dofs_per_cell2,
179  dofs_per_cell1));
180  interpolation_matrices[&cell1->get_fe()][&cell2->get_fe()]
181  = interpolation_matrix;
182 
183  get_interpolation_matrix(cell1->get_fe(),
184  cell2->get_fe(),
185  *interpolation_matrix);
186  }
187 
188  cell1->get_dof_values(u1, u1_local);
189  interpolation_matrices[&cell1->get_fe()][&cell2->get_fe()]
190  ->vmult(u2_local, u1_local);
191 
192  dofs.resize (dofs_per_cell2);
193  cell2->get_dof_indices(dofs);
194 
195  for (unsigned int i=0; i<dofs_per_cell2; ++i)
196  {
197  // if dof is locally_owned
198  const types::global_dof_index gdi = dofs[i];
199  if (u2_elements.is_element(gdi))
200  {
201  u2(dofs[i])+=u2_local(i);
202  touch_count(dofs[i]) += 1;
203  }
204  }
205  }
206  // cell1 is at the end, so should
207  // be cell2
208  Assert (cell2 == endc2, ExcInternalError());
209 
210  u2.compress(VectorOperation::add);
211  touch_count.compress(VectorOperation::add);
212 
213  // if we work on parallel distributed
214  // vectors, we have to ensure, that we only
215  // work on dofs this processor owns.
216  IndexSet locally_owned_dofs = dof2.locally_owned_dofs();
217 
218  // when a discontinuous element is
219  // interpolated to a continuous
220  // one, we take the mean values.
221  // for parallel vectors check,
222  // if this component is owned by
223  // this processor.
224  for (types::global_dof_index i=0; i<dof2.n_dofs(); ++i)
225  if (locally_owned_dofs.is_element(i))
226  {
227  Assert(static_cast<typename OutVector::value_type>(touch_count(i)) != typename OutVector::value_type(0),
228  ExcInternalError());
229  u2(i) /= touch_count(i);
230  }
231 
232  // finish the work on parallel vectors
233  u2.compress(VectorOperation::insert);
234  // Apply hanging node constraints.
235  constraints.distribute(u2);
236  }
237 
238 
239 
240  template <int dim,
241  template <int, int> class DoFHandlerType,
242  class InVector, class OutVector, int spacedim>
243  void
244  back_interpolate(const DoFHandlerType<dim, spacedim> &dof1,
245  const InVector &u1,
246  const FiniteElement<dim,spacedim> &fe2,
247  OutVector &u1_interpolated)
248  {
249  Assert(dof1.get_fe().n_components() == fe2.n_components(),
250  ExcDimensionMismatch(dof1.get_fe().n_components(), fe2.n_components()));
251  Assert(u1.size() == dof1.n_dofs(),
252  ExcDimensionMismatch(u1.size(), dof1.n_dofs()));
253  Assert(u1_interpolated.size() == dof1.n_dofs(),
254  ExcDimensionMismatch(u1_interpolated.size(), dof1.n_dofs()));
255 
256 #ifdef DEBUG
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"
262  " index sets."));
263  Assert(u1_interpolated_elements == dof1_local_dofs,
264  ExcMessage("The provided vector and DoF handler should have the same"
265  " index sets."));
266 #endif
267 
268  Vector<typename OutVector::value_type> u1_local(DoFTools::max_dofs_per_cell(dof1));
269  Vector<typename OutVector::value_type> u1_int_local(DoFTools::max_dofs_per_cell(dof1));
270 
271  const types::subdomain_id subdomain_id =
272  dof1.get_triangulation().locally_owned_subdomain();
273 
274  typename DoFHandlerType<dim, spacedim>::active_cell_iterator cell = dof1.begin_active(),
275  endc = dof1.end();
276 
277  // map from possible fe objects in
278  // dof1 to the back_interpolation
279  // matrices
280  std::map<const FiniteElement<dim> *,
281  std_cxx11::shared_ptr<FullMatrix<double> > > interpolation_matrices;
282 
283  for (; cell!=endc; ++cell)
284  if ((cell->subdomain_id() == subdomain_id)
285  ||
286  (subdomain_id == numbers::invalid_subdomain_id))
287  {
288  // For continuous elements on
289  // grids with hanging nodes we
290  // need hanging node
291  // constraints. Consequently,
292  // when the elements are
293  // continuous no hanging node
294  // constraints are allowed.
295  const bool hanging_nodes_not_allowed=
296  (cell->get_fe().dofs_per_vertex != 0) || (fe2.dofs_per_vertex != 0);
297 
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(),
303 
304  const unsigned int dofs_per_cell1 = cell->get_fe().dofs_per_cell;
305 
306  // make sure back_interpolation
307  // matrix is available
308  if (interpolation_matrices[&cell->get_fe()] == 0)
309  {
310  interpolation_matrices[&cell->get_fe()] =
311  std_cxx11::shared_ptr<FullMatrix<double> >
312  (new FullMatrix<double>(dofs_per_cell1, dofs_per_cell1));
313  get_back_interpolation_matrix(cell->get_fe(), fe2,
314  *interpolation_matrices[&cell->get_fe()]);
315  }
316 
317  u1_local.reinit (dofs_per_cell1);
318  u1_int_local.reinit (dofs_per_cell1);
319 
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);
323  };
324 
325  // if we work on a parallel vector, we have to finish the work
326  u1_interpolated.compress(VectorOperation::insert);
327  }
328 
329 
330 
331  namespace internal
332  {
333  namespace
334  {
335  template <int dim, int spacedim, class InVector>
336  void back_interpolate (const DoFHandler<dim,spacedim> &dof1,
337  const ConstraintMatrix &constraints1,
338  const InVector &u1,
339  const DoFHandler<dim,spacedim> &dof2,
340  const ConstraintMatrix &constraints2,
341  InVector &u1_interpolated)
342  {
344  interpolate(dof1, u1, dof2, constraints2, u2);
345  interpolate(dof2, u2, dof1, constraints1, u1_interpolated);
346  }
347 
348  // special version for PETSc
349 #ifdef DEAL_II_WITH_PETSC
350  template <int dim, int spacedim>
351  void back_interpolate (const DoFHandler<dim,spacedim> &dof1,
352  const ConstraintMatrix &constraints1,
353  const PETScWrappers::MPI::Vector &u1,
354  const DoFHandler<dim,spacedim> &dof2,
355  const ConstraintMatrix &constraints2,
356  PETScWrappers::MPI::Vector &u1_interpolated)
357  {
358  // if u1 is a parallel distributed PETSc vector, we create a
359  // vector u2 with based on the sets of locally owned and relevant
360  // dofs of dof2
361  IndexSet dof2_locally_owned_dofs = dof2.locally_owned_dofs();
362  IndexSet dof2_locally_relevant_dofs;
364  dof2_locally_relevant_dofs);
365 
366  PETScWrappers::MPI::Vector u2_out (dof2_locally_owned_dofs,
367  u1.get_mpi_communicator());
368  interpolate(dof1, u1, dof2, constraints2, u2_out);
369  PETScWrappers::MPI::Vector u2 (dof2_locally_owned_dofs,
370  dof2_locally_relevant_dofs,
371  u1.get_mpi_communicator());
372  u2 = u2_out;
373  interpolate(dof2, u2, dof1, constraints1, u1_interpolated);
374  }
375 #endif
376 
377  // special version for Trilinos
378 #ifdef DEAL_II_WITH_TRILINOS
379  template <int dim, int spacedim>
380  void back_interpolate (const DoFHandler<dim,spacedim> &dof1,
381  const ConstraintMatrix &constraints1,
383  const DoFHandler<dim,spacedim> &dof2,
384  const ConstraintMatrix &constraints2,
385  TrilinosWrappers::MPI::Vector &u1_interpolated)
386  {
387  // if u1 is a parallel distributed Trilinos vector, we create a
388  // vector u2 with based on the sets of locally owned and relevant
389  // dofs of dof2
390  IndexSet dof2_locally_owned_dofs = dof2.locally_owned_dofs();
391  IndexSet dof2_locally_relevant_dofs;
393  dof2_locally_relevant_dofs);
394 
395  TrilinosWrappers::MPI::Vector u2_out (dof2_locally_owned_dofs,
396  u1.get_mpi_communicator());
397  interpolate(dof1, u1, dof2, constraints2, u2_out);
398  TrilinosWrappers::MPI::Vector u2 (dof2_locally_owned_dofs,
399  dof2_locally_relevant_dofs,
400  u1.get_mpi_communicator());
401  u2 = u2_out;
402  interpolate(dof2, u2, dof1, constraints1, u1_interpolated);
403  }
404 #endif
405 
406  // special version for LinearAlgebra::distributed::Vector
407  template <int dim, int spacedim, typename Number>
408  void back_interpolate (const DoFHandler<dim,spacedim> &dof1,
409  const ConstraintMatrix &constraints1,
411  const DoFHandler<dim,spacedim> &dof2,
412  const ConstraintMatrix &constraints2,
414  {
415  IndexSet dof2_locally_owned_dofs = dof2.locally_owned_dofs();
416  IndexSet dof2_locally_relevant_dofs;
418  dof2_locally_relevant_dofs);
419 
421  u2 (dof2_locally_owned_dofs,
422  dof2_locally_relevant_dofs,
423  u1.get_mpi_communicator());
424 
425  interpolate(dof1, u1, dof2, constraints2, u2);
426  u2.update_ghost_values ();
427  interpolate(dof2, u2, dof1, constraints1, u1_interpolated);
428  }
429  }
430  }
431 
432 
433 
434  template <int dim, class InVector, class OutVector, int spacedim>
436  const ConstraintMatrix &constraints1,
437  const InVector &u1,
438  const DoFHandler<dim,spacedim> &dof2,
439  const ConstraintMatrix &constraints2,
440  OutVector &u1_interpolated)
441  {
442  // For discontinuous elements without constraints take the simpler version
443  // of the back_interpolate function.
444  if (dof1.get_fe().dofs_per_vertex==0 && dof2.get_fe().dofs_per_vertex==0
445  && constraints1.n_constraints()==0 && constraints2.n_constraints()==0)
446  back_interpolate(dof1, u1, dof2.get_fe(), u1_interpolated);
447  else
448  {
449  Assert(dof1.get_fe().n_components() == dof2.get_fe().n_components(),
450  ExcDimensionMismatch(dof1.get_fe().n_components(), dof2.get_fe().n_components()));
451  Assert(u1.size()==dof1.n_dofs(), ExcDimensionMismatch(u1.size(), dof1.n_dofs()));
452  Assert(u1_interpolated.size()==dof1.n_dofs(),
453  ExcDimensionMismatch(u1_interpolated.size(), dof1.n_dofs()));
454 
455  // For continuous elements first interpolate to dof2, taking into
456  // account constraints2, and then interpolate back to dof1 taking into
457  // account constraints1
458  internal::back_interpolate(dof1, constraints1, u1, dof2, constraints2,
459  u1_interpolated);
460  }
461  }
462 
463 
464 
465  template <int dim, class InVector, class OutVector, int spacedim>
467  const InVector &u1,
468  const FiniteElement<dim,spacedim> &fe2,
469  OutVector &u1_difference)
470  {
471  Assert(dof1.get_fe().n_components() == fe2.n_components(),
472  ExcDimensionMismatch(dof1.get_fe().n_components(), fe2.n_components()));
473  Assert(u1.size()==dof1.n_dofs(), ExcDimensionMismatch(u1.size(), dof1.n_dofs()));
474  Assert(u1_difference.size()==dof1.n_dofs(),
475  ExcDimensionMismatch(u1_difference.size(), dof1.n_dofs()));
476 
477 #ifdef DEBUG
478  const IndexSet &dof1_local_dofs = dof1.locally_owned_dofs();
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"
483  " index sets."));
484  Assert(u1_difference_elements == dof1_local_dofs,
485  ExcMessage("The provided vector and DoF handler should have the same"
486  " index sets."));
487 #endif
488 
489  // For continuous elements on grids
490  // with hanging nodes we need
491  // hanging node
492  // constraints. Consequently, when
493  // the elements are continuous no
494  // hanging node constraints are
495  // allowed.
496  const bool hanging_nodes_not_allowed=
497  (dof1.get_fe().dofs_per_vertex != 0) || (fe2.dofs_per_vertex != 0);
498 
499  const unsigned int dofs_per_cell=dof1.get_fe().dofs_per_cell;
500 
501  Vector<typename OutVector::value_type> u1_local(dofs_per_cell);
502  Vector<typename OutVector::value_type> u1_diff_local(dofs_per_cell);
503 
504  const types::subdomain_id subdomain_id =
505  dof1.get_triangulation().locally_owned_subdomain();
506 
507  FullMatrix<double> difference_matrix(dofs_per_cell, dofs_per_cell);
509  difference_matrix);
510 
512  endc = dof1.end();
513 
514  for (; cell!=endc; ++cell)
515  if ((cell->subdomain_id() == subdomain_id)
516  ||
517  (subdomain_id == numbers::invalid_subdomain_id))
518  {
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(),
524 
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);
528  }
529 
530  // if we work on a parallel PETSc vector
531  // we have to finish the work and
532  // update ghost values
533  u1_difference.compress(VectorOperation::insert);
534  }
535 
536 
537 
538  namespace internal
539  {
540  namespace
541  {
542  template <int dim, class InVector, class OutVector, int spacedim>
544  const ConstraintMatrix &constraints1,
545  const InVector &u1,
546  const DoFHandler<dim,spacedim> &dof2,
547  const ConstraintMatrix &constraints2,
548  OutVector &u1_difference)
549  {
550  back_interpolate(dof1, constraints1, u1, dof2, constraints2, u1_difference);
551  u1_difference.sadd(-1., 1., u1);
552  }
553 
554  // special version for Trilinos
555 #ifdef DEAL_II_WITH_TRILINOS
556  template <int dim, int spacedim>
558  const ConstraintMatrix &constraints1,
560  const DoFHandler<dim,spacedim> &dof2,
561  const ConstraintMatrix &constraints2,
562  TrilinosWrappers::MPI::Vector &u1_difference)
563  {
564  back_interpolate(dof1, constraints1, u1, dof2, constraints2, u1_difference);
565 
566  // Trilinos vectors with and without ghost entries are very different
567  // and we cannot use the sadd function directly, so we have to create
568  // a completely distributed vector first and copy the local entries
569  // from the vector with ghost entries
570  TrilinosWrappers::MPI::Vector u1_completely_distributed (u1_difference.vector_partitioner ());
571 
572  u1_completely_distributed = u1;
573 
574  u1_difference.sadd(-1, u1_completely_distributed);
575  }
576 #endif
577  }
578  }
579 
580 
581 
582  template <int dim, class InVector, class OutVector, int spacedim>
584  const ConstraintMatrix &constraints1,
585  const InVector &u1,
586  const DoFHandler<dim,spacedim> &dof2,
587  const ConstraintMatrix &constraints2,
588  OutVector &u1_difference)
589  {
590  // For discontinuous elements
591  // without constraints take the
592  // cheaper version of the
593  // interpolation_difference function.
594  if (dof1.get_fe().dofs_per_vertex==0 && dof2.get_fe().dofs_per_vertex==0
595  && constraints1.n_constraints()==0 && constraints2.n_constraints()==0)
596  interpolation_difference(dof1, u1, dof2.get_fe(), u1_difference);
597  else
598  {
599  internal::interpolation_difference(dof1, constraints1, u1, dof2, constraints2, u1_difference);
600  }
601  }
602 
603 
604 
605  template <int dim, class InVector, class OutVector, int spacedim>
607  const InVector &u1,
608  const DoFHandler<dim,spacedim> &dof2,
609  OutVector &u2)
610  {
612  Assert(dof1.get_fe().n_components() == dof2.get_fe().n_components(),
613  ExcDimensionMismatch(dof1.get_fe().n_components(), dof2.get_fe().n_components()));
614  Assert(u1.size()==dof1.n_dofs(), ExcDimensionMismatch(u1.size(), dof1.n_dofs()));
615  Assert(u2.size()==dof2.n_dofs(), ExcDimensionMismatch(u2.size(), dof2.n_dofs()));
616 
620 
621  const unsigned int n1 = dof1.get_fe().dofs_per_cell;
622  const unsigned int n2 = dof2.get_fe().dofs_per_cell;
623 
626  std::vector<types::global_dof_index> dofs(n2);
627 
628  FullMatrix<double> matrix(n2,n1);
629  get_projection_matrix(dof1.get_fe(), dof2.get_fe(), matrix);
630 
631  u2 = 0;
632  while (cell2 != end)
633  {
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)
638  {
639  u2(dofs[i])+=u2_local(i);
640  }
641 
642  ++cell1;
643  ++cell2;
644  }
645  }
646 } // end of namespace FETools
647 
648 
649 
650 /*-------------- Explicit Instantiations -------------------------------*/
651 #include "fe_tools_interpolate.inst"
652 
653 
654 /*---------------------------- fe_tools.cc ---------------------------*/
655 
656 DEAL_II_NAMESPACE_CLOSE
Definition: tria.h:67
void project_dg(const DoFHandler< dim, spacedim > &dof1, const InVector &u1, const DoFHandler< dim, spacedim > &dof2, OutVector &u2)
const types::subdomain_id invalid_subdomain_id
Definition: types.h:245
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
Definition: dof_handler.cc:756
void interpolation_difference(const DoFHandler< dim, spacedim > &dof1, const InVector &z1, const FiniteElement< dim, spacedim > &fe2, OutVector &z1_difference)
ActiveSelector::active_cell_iterator active_cell_iterator
Definition: dof_handler.h:228
const FiniteElement< dim, spacedim > & get_fe() const
void get_projection_matrix(const FiniteElement< dim, spacedim > &fe1, const FiniteElement< dim, spacedim > &fe2, FullMatrix< number > &matrix)
Definition: fe_tools.cc:1531
active_cell_iterator begin_active(const unsigned int level=0) const
Definition: dof_handler.cc:740
static ::ExceptionBase & ExcMessage(std::string arg1)
void get_interpolation_matrix(const FiniteElement< dim, spacedim > &fe1, const FiniteElement< dim, spacedim > &fe2, FullMatrix< number > &interpolation_matrix)
Definition: fe_tools.cc:1415
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:313
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
const Epetra_Map & vector_partitioner() const
size_type n_constraints() const
void get_back_interpolation_matrix(const FiniteElement< dim, spacedim > &fe1, const FiniteElement< dim, spacedim > &fe2, FullMatrix< number > &interpolation_matrix)
Definition: fe_tools.cc:1479
types::global_dof_index n_dofs() const
void extract_locally_relevant_dofs(const DoFHandlerType &dof_handler, IndexSet &dof_set)
Definition: dof_tools.cc:980
const MPI_Comm & get_mpi_communicator() const
const MPI_Comm & get_mpi_communicator() const
unsigned int subdomain_id
Definition: types.h:42
void distribute(VectorType &vec) const
unsigned int n_components() const
void get_interpolation_difference_matrix(const FiniteElement< dim, spacedim > &fe1, const FiniteElement< dim, spacedim > &fe2, FullMatrix< number > &difference_matrix)
Definition: fe_tools.cc:1505
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
Definition: fe_base.h:238
bool is_element(const size_type index) const
Definition: index_set.h:1489
const IndexSet & locally_owned_dofs() const
static ::ExceptionBase & ExcTriangulationMismatch()
void back_interpolate(const DoFHandlerType< dim, spacedim > &dof1, const InVector &u1, const FiniteElement< dim, spacedim > &fe2, OutVector &u1_interpolated)
void interpolate(const DoFHandlerType1< dim, spacedim > &dof1, const InVector &u1, const DoFHandlerType2< dim, spacedim > &dof2, OutVector &u2)
static ::ExceptionBase & ExcInternalError()