Loading [MathJax]/extensions/TeX/AMSsymbols.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
copy_data.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2019 - 2023 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_meshworker_copy_data_h
16#define dealii_meshworker_copy_data_h
17
18#include <deal.II/base/config.h>
19
21#include <deal.II/base/types.h>
22
24#include <deal.II/lac/vector.h>
25
26#include <algorithm>
27#include <vector>
28
30
31namespace MeshWorker
32{
50 template <int n_matrices = 1,
51 int n_vectors = n_matrices,
52 int n_dof_indices = n_matrices,
53 typename ScalarType = double>
54 struct CopyData
55 {
61 CopyData() = default;
62
67 explicit CopyData(const unsigned int size);
68
72 explicit CopyData(
73 const ndarray<unsigned int, n_matrices, 2> &matrix_sizes,
74 const std::array<unsigned int, n_vectors> &vector_sizes,
75 const std::array<unsigned int, n_dof_indices> &dof_indices_sizes);
76
80 CopyData(const CopyData &other) = default;
81
85 CopyData &
86 operator=(const CopyData &other) = default;
87
92 void
93 reinit(const unsigned int size);
94
106 void
107 reinit(const unsigned int index, const unsigned int size);
108
122 void
123 reinit(const unsigned int index,
124 const unsigned int size_rows,
125 const unsigned int size_columns);
126
130 std::array<FullMatrix<ScalarType>, n_matrices> matrices;
131
135 std::array<Vector<ScalarType>, n_vectors> vectors;
136
140 std::array<std::vector<types::global_dof_index>, n_dof_indices>
142 };
143
144
145#ifndef DOXYGEN
146 //
147 // Template definitions
148 //
149 template <int n_matrices,
150 int n_vectors,
151 int n_dof_indices,
152 typename ScalarType>
154 const unsigned int size)
155 {
156 reinit(size);
157 }
158
159
160
161 template <int n_matrices,
162 int n_vectors,
163 int n_dof_indices,
164 typename ScalarType>
166 const ndarray<unsigned int, n_matrices, 2> &matrix_sizes,
167 const std::array<unsigned int, n_vectors> &vector_sizes,
168 const std::array<unsigned int, n_dof_indices> &dof_indices_sizes)
169 {
170 for (unsigned int i = 0; i < n_matrices; ++i)
171 matrices[i].reinit(matrix_sizes[i++]);
172
173 for (unsigned int i = 0; i < n_vectors; ++i)
174 vectors[i].reinit(vector_sizes[i++]);
175
176 for (unsigned int i = 0; i < n_dof_indices; ++i)
177 local_dof_indices[i].resize(dof_indices_sizes[i++]);
178 }
179
180
181
182 template <int n_matrices,
183 int n_vectors,
184 int n_dof_indices,
185 typename ScalarType>
186 void
188 const unsigned int size)
189 {
190 for (auto &m : matrices)
191 m.reinit({size, size});
192 for (auto &v : vectors)
193 v.reinit(size);
194 for (auto &d : local_dof_indices)
195 d.resize(size);
196 }
197
198
199
200 template <int n_matrices,
201 int n_vectors,
202 int n_dof_indices,
203 typename ScalarType>
204 void
206 const unsigned int index,
207 const unsigned int size)
208 {
209 reinit(index, size, size);
210 }
211
212
213
214 template <int n_matrices,
215 int n_vectors,
216 int n_dof_indices,
217 typename ScalarType>
218 void
220 const unsigned int index,
221 const unsigned int size_rows,
222 const unsigned int size_columns)
223 {
224 // We permit different numbers of matrices, vectors and DoF index vectors.
225 // So we have to be a bit permissive here.
226 constexpr int max_index = std::max({n_matrices, n_vectors, n_dof_indices});
227 (void)max_index;
228 Assert(index < max_index, ExcIndexRange(index, 0, max_index));
229
230 if (index < n_matrices)
231 matrices[index].reinit({size_rows, size_columns});
232
233 if (index < n_vectors)
234 vectors[index].reinit(size_columns);
235
236 if (index < n_dof_indices)
237 local_dof_indices[index].resize(size_columns);
238 }
239
240#endif // DOXYGEN
241} // namespace MeshWorker
242
244
245#endif
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:35
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:36
#define Assert(cond, exc)
static ::ExceptionBase & ExcIndexRange(std::size_t arg1, std::size_t arg2, std::size_t arg3)
std::size_t size
Definition mpi.cc:745
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
void reinit(MatrixBlock< MatrixType > &v, const BlockSparsityPattern &p)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
typename internal::ndarray::HelperArray< T, Ns... >::type ndarray
Definition ndarray.h:107
CopyData(const CopyData &other)=default
void reinit(const unsigned int index, const unsigned int size_rows, const unsigned int size_columns)
CopyData & operator=(const CopyData &other)=default
std::array< FullMatrix< ScalarType >, n_matrices > matrices
Definition copy_data.h:130
CopyData(const unsigned int size)
std::array< Vector< ScalarType >, n_vectors > vectors
Definition copy_data.h:135
void reinit(const unsigned int size)
void reinit(const unsigned int index, const unsigned int size)
CopyData(const ndarray< unsigned int, n_matrices, 2 > &matrix_sizes, const std::array< unsigned int, n_vectors > &vector_sizes, const std::array< unsigned int, n_dof_indices > &dof_indices_sizes)
std::array< std::vector< types::global_dof_index >, n_dof_indices > local_dof_indices
Definition copy_data.h:141