23#include <deal.II/grid/tria_iterator.templates.h>
44template <
typename 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.");
61 template <
typename Number>
62 std::enable_if_t<!std::is_unsigned<Number>::value,
69 template <
typename Number>
70 std::enable_if_t<std::is_unsigned<Number>::value, Number>
79 template <
typename VectorType>
81 std::is_same<VectorType,
83 std::is_same<VectorType,
88 std::is_same<VectorType,
90 typename VectorType::value_type>>
::value ||
91 std::is_same<VectorType,
93 typename VectorType::value_type>>
::value;
101 decltype(std::declval<T const>().set_ghost_state(std::declval<bool>()));
103 template <
typename T>
105 is_supported_operation<set_ghost_state_t, T>;
107 template <
typename VectorType,
108 typename std::enable_if<has_set_ghost_state<VectorType>,
109 VectorType>::type * =
nullptr>
113 vector.set_ghost_state(ghosted);
116 template <
typename VectorType,
117 typename std::enable_if<!has_set_ghost_state<VectorType>,
118 VectorType>::type * =
nullptr>
138 OutputVector & values,
139 const bool perform_check)
144 if (perform_check && is_dealii_vector<OutputVector>)
146 const bool old_ghost_state = values.has_ghost_elements();
152 for (
unsigned int i = 0; i < cell.
get_fe().n_dofs_per_cell(); ++i)
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]) *
161 number>::real_type>::epsilon(),
162 ExcNonMatchingElementsSetDofValuesByInterpolation<number>(
163 local_values[i], local_values_old[i]));
183 OutputVector & values,
184 const unsigned int fe_index_,
187 OutputVector & values)> &processor)
189 const unsigned int fe_index =
195 if (cell.is_active() && !cell.is_artificial())
204 processor(cell, local_values, values);
209 ExcMessage(
"Incorrect size of local_values vector."));
212 cell.
get_fe().n_dofs_per_cell(),
215 cell.
get_fe().get_interpolation_matrix(
222 if ((tmp.
size() > 0) && (local_values.
size() > 0))
223 interpolation.
vmult(tmp, local_values);
226 processor(cell, tmp, values);
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."));
248 ExcVectorDoesNotMatch()));
251 ExcVectorDoesNotMatch()));
255 for (
unsigned int child = 0; child < cell.n_children(); ++child)
259 .
vmult(tmp, local_values);
261 *cell.
child(child), tmp, values, fe_index, processor);
270template <
int dim,
int spacedim,
bool lda>
271template <
class OutputVector,
typename number>
275 OutputVector & values,
276 const unsigned int fe_index_,
277 const bool perform_check)
const
279 internal::process_by_interpolation<dim, spacedim, lda, OutputVector, number>(
286 OutputVector & values) {
292template <
int dim,
int spacedim,
bool lda>
293template <
class OutputVector,
typename number>
298 OutputVector & values,
299 const unsigned int fe_index_)
const
301 internal::process_by_interpolation<dim, spacedim, lda, OutputVector, number>(
308 OutputVector & values) {
309 std::vector<types::global_dof_index> dof_indices(
310 cell.
get_fe().n_dofs_per_cell());
321#include "dof_accessor_set.inst"
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
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define Assert(cond, exc)
#define DeclException2(Exception2, type1, type2, outsequence)
static ::ExceptionBase & ExcNonMatchingElementsSetDofValuesByInterpolation(Number arg1, Number arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
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 > &)