deal.II version GIT relicensing-2017-g7ae82fad8b 2024-10-22 19:20: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\}}\)
Loading...
Searching...
No Matches
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
23
24#include <map>
25
27
28template <typename number>
30
31template <int dim, int spacedim>
33class DoFHandler;
34
35template <typename number>
36class FullMatrix;
37template <int dim, typename Number>
38class Function;
39template <typename MeshType>
41class InterGridMap;
42template <int dim, int spacedim>
43class Mapping;
44
45namespace hp
46{
47 template <int dim, int spacedim>
48 class MappingCollection;
49}
50
51namespace VectorTools
52{
73 template <int dim, int spacedim, typename VectorType>
76 const Mapping<dim, spacedim> &mapping,
77 const DoFHandler<dim, spacedim> &dof,
78 const Function<spacedim, typename VectorType::value_type> &function,
79 VectorType &vec,
80 const ComponentMask &component_mask = {},
81 const unsigned int level = numbers::invalid_unsigned_int);
82
88 template <int dim, int spacedim, typename VectorType>
91 const hp::MappingCollection<dim, spacedim> &mapping,
92 const DoFHandler<dim, spacedim> &dof,
93 const Function<spacedim, typename VectorType::value_type> &function,
94 VectorType &vec,
95 const ComponentMask &component_mask = {},
96 const unsigned int level = numbers::invalid_unsigned_int);
97
98
105 template <int dim, int spacedim, typename VectorType>
108 const DoFHandler<dim, spacedim> &dof,
109 const Function<spacedim, typename VectorType::value_type> &function,
110 VectorType &vec,
111 const ComponentMask &component_mask = {},
112 const unsigned int level = numbers::invalid_unsigned_int);
113
135 template <int dim, class InVector, class OutVector, int spacedim>
138 void interpolate(const DoFHandler<dim, spacedim> &dof_1,
139 const DoFHandler<dim, spacedim> &dof_2,
140 const FullMatrix<double> &transfer,
141 const InVector &data_1,
142 OutVector &data_2);
143
192 template <int dim, int spacedim, typename VectorType>
193 DEAL_II_CXX20_REQUIRES(concepts::is_writable_dealii_vector_type<VectorType>)
195 const Mapping<dim, spacedim> &mapping,
196 const DoFHandler<dim, spacedim> &dof_handler,
197 const std::map<types::material_id,
198 const Function<spacedim, typename VectorType::value_type> *>
199 &function_map,
200 VectorType &dst,
201 const ComponentMask &component_mask = {});
202
241 template <int dim, int spacedim, typename VectorType>
244 const DoFHandler<dim, spacedim> &dof_handler_1,
245 const VectorType &u1,
246 const DoFHandler<dim, spacedim> &dof_handler_2,
247 VectorType &u2);
248
261 template <int dim, int spacedim, typename VectorType>
262 DEAL_II_CXX20_REQUIRES(concepts::is_writable_dealii_vector_type<VectorType>)
264 const DoFHandler<dim, spacedim> &dof_handler_1,
265 const VectorType &u1,
266 const DoFHandler<dim, spacedim> &dof_handler_2,
267 const AffineConstraints<typename VectorType::value_type> &constraints,
268 VectorType &u2);
269
280 template <int dim, int spacedim, typename VectorType>
281 DEAL_II_CXX20_REQUIRES(concepts::is_writable_dealii_vector_type<VectorType>)
283 const InterGridMap<DoFHandler<dim, spacedim>> &intergridmap,
284 const VectorType &u1,
285 const AffineConstraints<typename VectorType::value_type> &constraints,
286 VectorType &u2);
287
311 template <int dim, int spacedim, typename VectorType>
312 DEAL_II_CXX20_REQUIRES(concepts::is_writable_dealii_vector_type<VectorType>)
314 const DoFHandler<dim, spacedim> &dof_handler_fine,
315 const VectorType &u_fine,
316 const DoFHandler<dim, spacedim> &dof_handler_coarse,
317 const AffineConstraints<typename VectorType::value_type>
318 &constraints_coarse,
319 VectorType &u_coarse);
320
344 template <int dim, int spacedim, typename VectorType>
345 DEAL_II_CXX20_REQUIRES(concepts::is_writable_dealii_vector_type<VectorType>)
347 const DoFHandler<dim, spacedim> &dof_handler_coarse,
348 const VectorType &u_coarse,
349 const DoFHandler<dim, spacedim> &dof_handler_fine,
350 const AffineConstraints<typename VectorType::value_type> &constraints_fine,
351 VectorType &u_fine);
352
386 template <int dim, int spacedim, typename VectorType>
387 DEAL_II_CXX20_REQUIRES(concepts::is_writable_dealii_vector_type<VectorType>)
388 void get_position_vector(const DoFHandler<dim, spacedim> &dh,
389 VectorType &vector,
390 const ComponentMask &mask = {});
391
402 template <int dim, int spacedim, typename VectorType>
404 void get_position_vector(const Mapping<dim, spacedim> &mapping,
405 const DoFHandler<dim, spacedim> &dh,
406 VectorType &vector,
407 const ComponentMask &mask = {});
408
410} // namespace VectorTools
411
413
414#endif // dealii_vector_tools_interpolate_h
Abstract base class for mapping classes.
Definition mapping.h:318
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:175
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
unsigned int level
Definition grid_out.cc:4626
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
static const unsigned int invalid_unsigned_int
Definition types.h:220
STL namespace.
Definition types.h:32