Reference documentation for deal.II version 9.4.1
\(\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
dof_accessor_set.cc
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1998 - 2022 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.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
19
20#include <deal.II/fe/fe.h>
21
23#include <deal.II/grid/tria_iterator.templates.h>
24
26
38#include <deal.II/lac/vector.h>
39
40#include <vector>
41
43
44template <typename Number>
46 Number,
47 Number,
48 << "Called set_dof_values_by_interpolation(), but"
49 << " the element to be set, value " << std::setprecision(16)
50 << arg1 << ", does not match with the non-zero value "
51 << std::setprecision(16) << arg2 << " already set before.");
52
53namespace internal
54{
55#ifdef DEBUG
61 template <typename Number>
62 std::enable_if_t<!std::is_unsigned<Number>::value,
64 get_abs(const Number a)
65 {
66 return std::abs(a);
67 }
68
69 template <typename Number>
70 std::enable_if_t<std::is_unsigned<Number>::value, Number>
71 get_abs(const Number a)
72 {
73 return a;
74 }
75
79 template <typename VectorType>
80 constexpr bool is_dealii_vector =
81 std::is_same<VectorType,
83 std::is_same<VectorType,
85 std::is_same<
86 VectorType,
88 std::is_same<VectorType,
90 typename VectorType::value_type>>::value ||
91 std::is_same<VectorType,
93 typename VectorType::value_type>>::value;
94
99 template <typename T>
101 decltype(std::declval<T const>().set_ghost_state(std::declval<bool>()));
102
103 template <typename T>
104 constexpr bool has_set_ghost_state =
105 is_supported_operation<set_ghost_state_t, T>;
106
107 template <typename VectorType,
108 typename std::enable_if<has_set_ghost_state<VectorType>,
109 VectorType>::type * = nullptr>
110 void
111 set_ghost_state(VectorType &vector, const bool ghosted)
112 {
113 vector.set_ghost_state(ghosted);
114 }
115
116 template <typename VectorType,
117 typename std::enable_if<!has_set_ghost_state<VectorType>,
118 VectorType>::type * = nullptr>
119 void
120 set_ghost_state(VectorType &, const bool)
121 {
122 // serial vector: nothing to do
123 }
124#endif
125
130 template <int dim,
131 int spacedim,
132 bool lda,
133 class OutputVector,
134 typename number>
135 void
137 const Vector<number> & local_values,
138 OutputVector & values,
139 const bool perform_check)
140 {
141 (void)perform_check;
142
143#ifdef DEBUG
144 if (perform_check && is_dealii_vector<OutputVector>)
145 {
146 const bool old_ghost_state = values.has_ghost_elements();
147 set_ghost_state(values, true);
148
149 Vector<number> local_values_old(cell.get_fe().n_dofs_per_cell());
150 cell.get_dof_values(values, local_values_old);
151
152 for (unsigned int i = 0; i < cell.get_fe().n_dofs_per_cell(); ++i)
153 {
154 // a check consistent with the one in
155 // Utilities::MPI::Partitioner::import_from_ghosted_array_finish()
156 Assert(local_values_old[i] == number() ||
157 get_abs(local_values_old[i] - local_values[i]) <=
158 get_abs(local_values_old[i] + local_values[i]) *
159 100000. *
160 std::numeric_limits<typename numbers::NumberTraits<
161 number>::real_type>::epsilon(),
162 ExcNonMatchingElementsSetDofValuesByInterpolation<number>(
163 local_values[i], local_values_old[i]));
164 }
165
166 set_ghost_state(values, old_ghost_state);
167 }
168#endif
169
170 cell.set_dof_values(local_values, values);
171 }
172
173
174 template <int dim,
175 int spacedim,
176 bool lda,
177 class OutputVector,
178 typename number>
179 void
182 const Vector<number> & local_values,
183 OutputVector & values,
184 const unsigned int fe_index_,
185 const std::function<void(const DoFCellAccessor<dim, spacedim, lda> &cell,
186 const Vector<number> &local_values,
187 OutputVector & values)> &processor)
188 {
189 const unsigned int fe_index =
190 (cell.get_dof_handler().has_hp_capabilities() == false &&
193 fe_index_;
194
195 if (cell.is_active() && !cell.is_artificial())
196 {
197 if ((cell.get_dof_handler().has_hp_capabilities() == false) ||
198 // for hp-DoFHandlers, we need to require that on
199 // active cells, you either don't specify an fe_index,
200 // or that you specify the correct one
201 (fe_index == cell.active_fe_index()) ||
203 // simply set the values on this cell
204 processor(cell, local_values, values);
205 else
206 {
207 Assert(local_values.size() ==
208 cell.get_dof_handler().get_fe(fe_index).n_dofs_per_cell(),
209 ExcMessage("Incorrect size of local_values vector."));
210
211 FullMatrix<double> interpolation(
212 cell.get_fe().n_dofs_per_cell(),
213 cell.get_dof_handler().get_fe(fe_index).n_dofs_per_cell());
214
215 cell.get_fe().get_interpolation_matrix(
216 cell.get_dof_handler().get_fe(fe_index), interpolation);
217
218 // do the interpolation to the target space. for historical
219 // reasons, matrices are set to size 0x0 internally even if
220 // we reinit as 4x0, so we have to treat this case specially
221 Vector<number> tmp(cell.get_fe().n_dofs_per_cell());
222 if ((tmp.size() > 0) && (local_values.size() > 0))
223 interpolation.vmult(tmp, local_values);
224
225 // now set the dof values in the global vector
226 processor(cell, tmp, values);
227 }
228 }
229 else
230 // otherwise distribute them to the children
231 {
232 Assert((cell.get_dof_handler().has_hp_capabilities() == false) ||
235 "You cannot call this function on non-active cells "
236 "of DoFHandler objects unless you provide an explicit "
237 "finite element index because they do not have naturally "
238 "associated finite element spaces associated: degrees "
239 "of freedom are only distributed on active cells for which "
240 "the active FE index has been set."));
241
243 cell.get_dof_handler().get_fe(fe_index);
244 const unsigned int dofs_per_cell = fe.n_dofs_per_cell();
245
246 Assert(local_values.size() == dofs_per_cell,
248 ExcVectorDoesNotMatch()));
249 Assert(values.size() == cell.get_dof_handler().n_dofs(),
251 ExcVectorDoesNotMatch()));
252
253 Vector<number> tmp(dofs_per_cell);
254
255 for (unsigned int child = 0; child < cell.n_children(); ++child)
256 {
257 if (tmp.size() > 0)
258 fe.get_prolongation_matrix(child, cell.refinement_case())
259 .vmult(tmp, local_values);
261 *cell.child(child), tmp, values, fe_index, processor);
262 }
263 }
264 }
265
266} // namespace internal
267
268
269
270template <int dim, int spacedim, bool lda>
271template <class OutputVector, typename number>
272void
274 const Vector<number> &local_values,
275 OutputVector & values,
276 const unsigned int fe_index_,
277 const bool perform_check) const
278{
279 internal::process_by_interpolation<dim, spacedim, lda, OutputVector, number>(
280 *this,
281 local_values,
282 values,
283 fe_index_,
284 [perform_check](const DoFCellAccessor<dim, spacedim, lda> &cell,
285 const Vector<number> & local_values,
286 OutputVector & values) {
287 internal::set_dof_values(cell, local_values, values, perform_check);
288 });
289}
290
291
292template <int dim, int spacedim, bool lda>
293template <class OutputVector, typename number>
294void
297 const Vector<number> &local_values,
298 OutputVector & values,
299 const unsigned int fe_index_) const
300{
301 internal::process_by_interpolation<dim, spacedim, lda, OutputVector, number>(
302 *this,
303 local_values,
304 values,
305 fe_index_,
307 const Vector<number> & local_values,
308 OutputVector & values) {
309 std::vector<types::global_dof_index> dof_indices(
310 cell.get_fe().n_dofs_per_cell());
311 cell.get_dof_indices(dof_indices);
313 dof_indices,
314 values);
315 });
316}
317
318
319// --------------------------------------------------------------------------
320// explicit instantiations
321#include "dof_accessor_set.inst"
322
void distribute_local_to_global(const InVector &local_vector, const std::vector< size_type > &local_dof_indices, OutVector &global_vector) const
const DoFHandler< dim, spacedim > & get_dof_handler() const
void get_dof_values(const InputVector &values, Vector< number > &local_values) const
const FiniteElement< dimension_, space_dimension_ > & get_fe() const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > child(const unsigned int i) const
void distribute_local_to_global_by_interpolation(const Vector< number > &local_values, OutputVector &values, const unsigned int fe_index=DoFHandler< dimension_, space_dimension_ >::invalid_fe_index) const
unsigned int active_fe_index() const
void set_dof_values(const Vector< number > &local_values, OutputVector &values) const
void set_dof_values_by_interpolation(const Vector< number > &local_values, OutputVector &values, const unsigned int fe_index=DoFHandler< dimension_, space_dimension_ >::invalid_fe_index, const bool perform_check=false) const
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
bool has_hp_capabilities() const
types::global_dof_index n_dofs() const
unsigned int n_dofs_per_cell() const
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
void vmult(Vector< number2 > &w, const Vector< number2 > &v, const bool adding=false) const
Definition: vector.h:109
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
#define Assert(cond, exc)
Definition: exceptions.h:1473
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:532
static ::ExceptionBase & ExcNonMatchingElementsSetDofValuesByInterpolation(Number arg1, Number arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
size_type size() const
void process_by_interpolation(const DoFCellAccessor< dim, spacedim, lda > &cell, const Vector< number > &local_values, OutputVector &values, const unsigned int fe_index_, const std::function< void(const DoFCellAccessor< dim, spacedim, lda > &cell, const Vector< number > &local_values, OutputVector &values)> &processor)
constexpr bool has_set_ghost_state
decltype(std::declval< T const >().set_ghost_state(std::declval< bool >())) set_ghost_state_t
constexpr bool is_dealii_vector
void set_ghost_state(VectorType &vector, const bool ghosted)
void set_dof_values(const DoFCellAccessor< dim, spacedim, lda > &cell, const Vector< number > &local_values, OutputVector &values, const bool perform_check)
std::enable_if_t<!std::is_unsigned< Number >::value, typename numbers::NumberTraits< Number >::real_type > get_abs(const Number a)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)