Loading [MathJax]/extensions/TeX/newcommand.js
 deal.II version GIT relicensing-3075-gc235bd4825 2025-04-15 08:10:00+00:00
\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}} \newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=} \newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]} \newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages Concepts
vector_tools_interpolate.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2020 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_vector_tools_interpolate_h
16#define dealii_vector_tools_interpolate_h
17
18#include <deal.II/base/config.h>
19
21#include <deal.II/base/types.h>
22
24
25#include <map>
26
28
29template <typename number>
31
32template <int dim, int spacedim>
34class DoFHandler;
35
36template <typename number>
37class FullMatrix;
38template <int dim, typename Number>
39class Function;
40template <typename MeshType>
42class InterGridMap;
43template <int dim, int spacedim>
44class Mapping;
45
46namespace hp
47{
48 template <int dim, int spacedim>
49 class MappingCollection;
50}
51
52namespace VectorTools
53{
74 template <int dim, int spacedim, typename VectorType>
77 const Mapping<dim, spacedim> &mapping,
78 const DoFHandler<dim, spacedim> &dof,
79 const Function<spacedim, typename VectorType::value_type> &function,
80 VectorType &vec,
81 const ComponentMask &component_mask = {},
82 const unsigned int level = numbers::invalid_unsigned_int);
83
89 template <int dim, int spacedim, typename VectorType>
92 const hp::MappingCollection<dim, spacedim> &mapping,
93 const DoFHandler<dim, spacedim> &dof,
94 const Function<spacedim, typename VectorType::value_type> &function,
95 VectorType &vec,
96 const ComponentMask &component_mask = {},
97 const unsigned int level = numbers::invalid_unsigned_int);
98
99
106 template <int dim, int spacedim, typename VectorType>
109 const DoFHandler<dim, spacedim> &dof,
110 const Function<spacedim, typename VectorType::value_type> &function,
111 VectorType &vec,
112 const ComponentMask &component_mask = {},
113 const unsigned int level = numbers::invalid_unsigned_int);
114
136 template <int dim, class InVector, class OutVector, int spacedim>
139 void interpolate(const DoFHandler<dim, spacedim> &dof_1,
140 const DoFHandler<dim, spacedim> &dof_2,
141 const FullMatrix<double> &transfer,
142 const InVector &data_1,
143 OutVector &data_2);
144
193 template <int dim, int spacedim, typename VectorType>
194 DEAL_II_CXX20_REQUIRES(concepts::is_writable_dealii_vector_type<VectorType>)
196 const Mapping<dim, spacedim> &mapping,
197 const DoFHandler<dim, spacedim> &dof_handler,
198 const std::map<types::material_id,
199 const Function<spacedim, typename VectorType::value_type> *>
200 &function_map,
201 VectorType &dst,
202 const ComponentMask &component_mask = {});
203
242 template <int dim, int spacedim, typename VectorType>
245 const DoFHandler<dim, spacedim> &dof_handler_1,
246 const VectorType &u1,
247 const DoFHandler<dim, spacedim> &dof_handler_2,
248 VectorType &u2);
249
262 template <int dim, int spacedim, typename VectorType>
263 DEAL_II_CXX20_REQUIRES(concepts::is_writable_dealii_vector_type<VectorType>)
265 const DoFHandler<dim, spacedim> &dof_handler_1,
266 const VectorType &u1,
267 const DoFHandler<dim, spacedim> &dof_handler_2,
268 const AffineConstraints<typename VectorType::value_type> &constraints,
269 VectorType &u2);
270
281 template <int dim, int spacedim, typename VectorType>
282 DEAL_II_CXX20_REQUIRES(concepts::is_writable_dealii_vector_type<VectorType>)
284 const InterGridMap<DoFHandler<dim, spacedim>> &intergridmap,
285 const VectorType &u1,
286 const AffineConstraints<typename VectorType::value_type> &constraints,
287 VectorType &u2);
288
312 template <int dim, int spacedim, typename VectorType>
313 DEAL_II_CXX20_REQUIRES(concepts::is_writable_dealii_vector_type<VectorType>)
315 const DoFHandler<dim, spacedim> &dof_handler_fine,
316 const VectorType &u_fine,
317 const DoFHandler<dim, spacedim> &dof_handler_coarse,
318 const AffineConstraints<typename VectorType::value_type>
319 &constraints_coarse,
320 VectorType &u_coarse);
321
345 template <int dim, int spacedim, typename VectorType>
346 DEAL_II_CXX20_REQUIRES(concepts::is_writable_dealii_vector_type<VectorType>)
348 const DoFHandler<dim, spacedim> &dof_handler_coarse,
349 const VectorType &u_coarse,
350 const DoFHandler<dim, spacedim> &dof_handler_fine,
351 const AffineConstraints<typename VectorType::value_type> &constraints_fine,
352 VectorType &u_fine);
353
387 template <int dim, int spacedim, typename VectorType>
388 DEAL_II_CXX20_REQUIRES(concepts::is_writable_dealii_vector_type<VectorType>)
389 void get_position_vector(const DoFHandler<dim, spacedim> &dh,
390 VectorType &vector,
391 const ComponentMask &mask = {});
392
403 template <int dim, int spacedim, typename VectorType>
405 void get_position_vector(const Mapping<dim, spacedim> &mapping,
406 const DoFHandler<dim, spacedim> &dh,
407 VectorType &vector,
408 const ComponentMask &mask = {});
409
411} // namespace VectorTools
412
414
415#endif // dealii_vector_tools_interpolate_h
Abstract base class for mapping classes.
Definition mapping.h:320
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:35
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:242
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:36
unsigned int level
Definition grid_out.cc:4635
void interpolate_to_coarser_mesh(const DoFHandler< dim, spacedim > &dof_handler_fine, const VectorType &u_fine, const DoFHandler< dim, spacedim > &dof_handler_coarse, const AffineConstraints< typename VectorType::value_type > &constraints_coarse, VectorType &u_coarse)
void interpolate(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Function< spacedim, typename VectorType::value_type > &function, VectorType &vec, const ComponentMask &component_mask={}, const unsigned int level=numbers::invalid_unsigned_int)
void interpolate_based_on_material_id(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof_handler, const std::map< types::material_id, const Function< spacedim, typename VectorType::value_type > * > &function_map, VectorType &dst, const ComponentMask &component_mask={})
void get_position_vector(const DoFHandler< dim, spacedim > &dh, VectorType &vector, const ComponentMask &mask={})
void interpolate_to_finer_mesh(const DoFHandler< dim, spacedim > &dof_handler_coarse, const VectorType &u_coarse, const DoFHandler< dim, spacedim > &dof_handler_fine, const AffineConstraints< typename VectorType::value_type > &constraints_fine, VectorType &u_fine)
void interpolate_to_different_mesh(const DoFHandler< dim, spacedim > &dof_handler_1, const VectorType &u1, const DoFHandler< dim, spacedim > &dof_handler_2, VectorType &u2)
Definition hp.h:117
constexpr unsigned int invalid_unsigned_int
Definition types.h:238
STL namespace.
Definition types.h:32