15#ifndef dealii_dof_accessor_h
16#define dealii_dof_accessor_h
33#include <boost/container/small_vector.hpp>
45template <
typename number>
47template <
typename number>
49template <
typename number>
52template <
typename Accessor>
105 template <
int structdim,
int dim,
int spacedim>
121 template <
int dim,
int spacedim>
212template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
233 using BaseClass = typename ::internal::DoFAccessorImplementation::
234 Inheritance<structdim, dimension, space_dimension>::BaseClass;
303 template <
int structdim2,
int dim2,
int spacedim2>
310 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
317 template <
bool level_dof_access2>
354 template <
bool level_dof_access2>
390 typename ::internal::DoFHandlerImplementation::
391 Iterators<dim, spacedim, level_dof_access>::line_iterator
392 line(
const unsigned int i)
const;
399 typename ::internal::DoFHandlerImplementation::
400 Iterators<dim, spacedim, level_dof_access>::quad_iterator
401 quad(
const unsigned int i)
const;
450 std::vector<types::global_dof_index> &dof_indices,
462 std::vector<types::global_dof_index> &dof_indices,
471 const std::vector<types::global_dof_index> &dof_indices,
497 const unsigned int vertex,
498 const unsigned int i,
509 const unsigned int vertex,
510 const unsigned int i,
590 std::set<types::fe_index>
623 "This accessor object has not been "
624 "associated with any DoFHandler object.");
678 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
686 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
718 const unsigned int i,
724 const unsigned int i,
730 const unsigned int vertex,
731 const unsigned int i,
739 template <
int,
int,
int,
bool>
747 friend struct ::internal::DoFHandlerImplementation::Policy::
749 friend struct ::internal::DoFHandlerImplementation::Implementation;
750 friend struct ::internal::hp::DoFHandlerImplementation::Implementation;
751 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
752 friend struct ::internal::DoFAccessorImplementation::Implementation;
764template <
int spacedim,
bool level_dof_access>
849 template <
int structdim2,
int dim2,
int spacedim2>
856 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
908 template <
bool level_dof_access2>
932 child(const
unsigned int c) const;
940 typename ::
internal::DoFHandlerImplementation::
941 Iterators<1, spacedim, level_dof_access>::line_iterator
942 line(const
unsigned int i) const;
950 typename ::
internal::DoFHandlerImplementation::
951 Iterators<1, spacedim, level_dof_access>::quad_iterator
952 quad(const
unsigned int i) const;
999 std::vector<
types::global_dof_index> &dof_indices,
1000 const
types::fe_index fe_index =
numbers::invalid_fe_index) const;
1011 std::vector<
types::global_dof_index> &dof_indices,
1012 const
types::fe_index fe_index =
numbers::invalid_fe_index) const;
1032 types::global_dof_index
1034 const
unsigned int vertex,
1035 const
unsigned int i,
1036 const
types::fe_index fe_index =
numbers::invalid_fe_index) const;
1055 types::global_dof_index
1057 const
types::fe_index fe_index =
numbers::invalid_fe_index) const;
1119 "This accessor object has not been "
1120 "associated with any DoFHandler object.");
1163 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
1166 const
DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1171 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
1174 const
DoFAccessor<structdim2, dim2, spacedim2, level_dof_access2> &) const;
1202 const
unsigned int i,
1204 const
types::fe_index fe_index =
numbers::invalid_fe_index) const;
1216 friend struct ::
internal::DoFHandlerImplementation::Policy::
1248template <
int structdim,
int dim,
int spacedim = dim>
1266 const int level = -1,
1267 const int index = -1,
1283 template <
typename OtherAccessor>
1303 const unsigned int i,
1324template <
int dimension_,
int space_dimension_,
bool level_dof_access>
1334 static const unsigned int dim = dimension_;
1395 template <
int structdim2,
int dim2,
int spacedim2>
1402 template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
1505 boost::container::small_vector<
1535 const unsigned int subface_no)
const;
1545 const unsigned int subface_no)
const;
1575 template <
class InputVector,
typename number>
1596 template <
typename Number,
typename ForwardIterator>
1599 ForwardIterator local_values_begin,
1600 ForwardIterator local_values_end)
const;
1622 template <
class InputVector,
typename ForwardIterator>
1626 const InputVector &values,
1627 ForwardIterator local_values_begin,
1628 ForwardIterator local_values_end)
const;
1654 template <
class OutputVector,
typename number>
1657 OutputVector &values)
const;
1690 template <
typename Number>
1758 template <
class OutputVector,
typename number>
1762 OutputVector &values,
1764 const bool perform_check =
false)
const;
1775 template <
class OutputVector,
typename number>
1779 OutputVector &values,
1796 template <
typename number,
typename OutputVector>
1799 OutputVector &global_destination)
const;
1815 template <
typename ForwardIterator,
typename OutputVector>
1818 ForwardIterator local_source_end,
1819 OutputVector &global_destination)
const;
1834 template <
typename ForwardIterator,
typename OutputVector>
1838 ForwardIterator local_source_begin,
1839 ForwardIterator local_source_end,
1840 OutputVector &global_destination)
const;
1848 template <
typename number,
typename OutputMatrix>
1851 OutputMatrix &global_destination)
const;
1857 template <
typename number,
typename OutputMatrix,
typename OutputVector>
1861 OutputMatrix &global_matrix,
1862 OutputVector &global_vector)
const;
1890 std::vector<types::global_dof_index> &dof_indices)
const;
2117 friend struct ::internal::DoFCellAccessorImplementation::Implementation;
2121template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
2125 return level_dof_access;
2130template <
int structdim,
int dim,
int spacedim>
2131template <
typename OtherAccessor>
2133 const OtherAccessor &)
2136 ExcMessage(
"You are attempting an illegal conversion between "
2137 "iterator/accessor types. The constructor you call "
2138 "only exists to make certain template constructs "
2139 "easier to write as dimension independent code but "
2140 "the conversion is not valid in the current context."));
2148template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
2156template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
2162 : ::
internal::DoFAccessorImplementation::
2163 Inheritance<structdim, dim, spacedim>::
BaseClass(tria, level, index)
2167 tria ==
nullptr || &
dof_handler->get_triangulation() == tria,
2169 "You can't create a DoF accessor in which the DoFHandler object "
2170 "uses a different triangulation than the one you pass as argument."));
2175template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
2176template <
int structdim2,
int dim2,
int spacedim2>
2185template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
2186template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
2194 "You are trying to assign iterators that are incompatible. "
2195 "The reason for incompatibility is that they refer to objects of "
2196 "different dimensionality (e.g., assigning a line iterator "
2197 "to a quad iterator)."));
2202template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
2203template <
bool level_dof_access2>
2212template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
2223template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
2233template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
2239 BaseClass::copy_from(da);
2244template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
2245template <
bool level_dof_access2>
2250 BaseClass::copy_from(a);
2256template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
2257template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
2264 return (BaseClass::operator==(a));
2269template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
2270template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
2277 return (BaseClass::operator!=(a));
2282template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
2285 const unsigned int i)
const
2287 Assert(
static_cast<unsigned int>(this->present_level) <
2302 namespace DoFAccessorImplementation
2308 template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
2322 "It is not possible to specify a FE index if no hp support is used!"));
2342 "You need to specify a FE index if hp support is used!"));
2363 boost::container::small_vector<::types::global_dof_index, 27>;
2378 typename GlobalIndexType,
2379 typename DoFPProcessor>
2382 const unsigned int obj_level,
2383 const unsigned int obj_index,
2385 const unsigned int local_index,
2386 const std::integral_constant<int, structdim> &,
2387 GlobalIndexType &global_index,
2388 const DoFPProcessor &process)
2400 [obj_level][structdim]
2409 if (structdim == dim)
2413 [obj_level][structdim]
2441 "You are requesting an active FE index that is not assigned "
2442 "to any of the cells connected to this entity."));
2457 [obj_level][structdim]
2466 [obj_level][structdim]
2468 [obj_level][structdim]
2486 template <
int dim,
int spacedim,
int structdim>
2487 static std::pair<unsigned int, unsigned int>
2489 const unsigned int obj_level,
2490 const unsigned int obj_index,
2492 const std::integral_constant<int, structdim> &)
2498 if (structdim == dim)
2500 const unsigned int ptr_0 =
2502 const unsigned int length =
2503 dof_handler.
get_fe(fe_index).template n_dofs_per_object<dim>(0);
2505 return {ptr_0, length};
2514 const unsigned int ptr_0 =
2516 const unsigned int length =
2517 dof_handler.
object_dof_ptr[obj_level][structdim][obj_index + 1] -
2520 return {ptr_0, length};
2532 const auto fe_index_local_ptr =
2539 Assert(fe_index_local_ptr !=
2543 "You tried to call a function accessing DoF indices, but "
2544 "they appear not be available (yet) or inconsistent. "
2545 "Did you call distribute_dofs() first? Alternatively, if "
2546 "you are using different elements on different cells (i.e., "
2547 "you are using the hp capabilities of deal.II), did you "
2548 "change the active_fe_index of a cell since you last "
2549 "called distribute_dofs()?"));
2554 fe_index_local_ptr);
2560 const unsigned int ptr_0 =
2565 const unsigned int ptr_1 =
2569 fe_index_local + 1];
2571 return {ptr_0, ptr_1 - ptr_0};
2574 template <
int dim,
int spacedim,
int structdim,
bool level_dof_access>
2575 static std::pair<unsigned int, unsigned int>
2577 const ::DoFAccessor<structdim, dim, spacedim, level_dof_access>
2585 std::integral_constant<int, structdim>());
2588 template <
int dim,
int spacedim,
int structdim>
2589 static std::pair<unsigned int, unsigned int>
2611 typename DoFProcessor,
2612 typename DoFMapping>
2615 const unsigned int obj_level,
2616 const unsigned int obj_index,
2619 const std::integral_constant<int, structdim> &dd,
2621 const DoFProcessor &process)
2628 if (range.second == 0)
2631 std::vector<types::global_dof_index> &object_dof_indices =
2636 object_dof_indices.data() + range.first;
2639 for (
unsigned int i = 0; i < range.second; ++i)
2642 stored_indices[(structdim == 0 || structdim == dim) ? i :
2645 if (dof_indices_ptr !=
nullptr)
2662 template <
int dim,
int spacedim,
int structdim>
2665 const unsigned int obj_level,
2666 const unsigned int obj_index,
2668 const unsigned int local_index,
2669 const std::integral_constant<int, structdim> &dd,
2679 [](
auto &ptr,
const auto &
value) { ptr =
value; });
2693 template <
int dim,
int spacedim,
int structdim>
2696 const unsigned int obj_level,
2697 const unsigned int obj_index,
2699 const unsigned int local_index,
2700 const std::integral_constant<int, structdim> &dd)
2710 [](
const auto &ptr,
auto &
value) {
value = ptr; });
2711 return global_index;
2715 template <
int dim,
int spacedim>
2719 const unsigned int vertex_index,
2720 const unsigned int i)
2724 "DoFHandler in hp-mode does not implement multilevel DoFs."));
2741 template <
int dim,
int spacedim,
int structdim>
2744 const unsigned int obj_level,
2745 const unsigned int obj_index,
2746 const std::integral_constant<int, structdim> &)
2757 if (structdim == dim)
2780 template <
int dim,
int spacedim,
int structdim>
2783 const unsigned int obj_level,
2784 const unsigned int obj_index,
2785 const unsigned int local_index,
2786 const std::integral_constant<int, structdim> &)
2791 Assert(((structdim == dim) &&
2801 if (structdim == dim)
2832 template <
int dim,
int spacedim,
int structdim>
2833 static std::set<types::fe_index>
2835 const unsigned int obj_level,
2836 const unsigned int obj_index,
2837 const std::integral_constant<int, structdim> &t)
2846 if (structdim == dim)
2850 std::set<types::fe_index> active_fe_indices;
2851 for (
unsigned int i = 0;
2854 active_fe_indices.insert(
2856 return active_fe_indices;
2861 template <
int dim,
int spacedim,
int structdim>
2864 const unsigned int obj_level,
2865 const unsigned int obj_index,
2867 const std::integral_constant<int, structdim> &)
2876 if (structdim == dim)
2893 template <
typename InputVector,
typename ForwardIterator>
2898 ForwardIterator local_values_begin)
2900 values.extract_subvector_to(cache, cache_end, local_values_begin);
2905#ifdef DEAL_II_WITH_TRILINOS
2906 static std::vector<unsigned int>
2911 std::vector<unsigned int> idx(v_end - v_begin);
2912 std::iota(idx.begin(), idx.end(), 0u);
2915 std::sort(idx.begin(),
2917 [&v_begin](
unsigned int i1,
unsigned int i2) {
2918 return *(v_begin + i1) < *(v_begin + i2);
2926# ifdef DEAL_II_TRILINOS_WITH_TPETRA
2927 template <
typename ForwardIterator,
typename Number,
typename MemorySpace>
2934 ForwardIterator local_values_begin)
2936 std::vector<unsigned int> sorted_indices_pos =
2938 const unsigned int cache_size = cache_end - cache_begin;
2939 std::vector<types::global_dof_index> cache_indices(cache_size);
2940 for (
unsigned int i = 0; i < cache_size; ++i)
2941 cache_indices[i] = *(cache_begin + sorted_indices_pos[i]);
2943 IndexSet index_set(cache_indices.back() + 1);
2944 index_set.
add_indices(cache_indices.begin(), cache_indices.end());
2950 for (
unsigned int i = 0; i < cache_size; ++i, ++local_values_begin)
2951 *local_values_begin = read_write_vector[sorted_indices_pos[i]];
2957 template <
typename ForwardIterator>
2962 ForwardIterator local_values_begin)
2964 std::vector<unsigned int> sorted_indices_pos =
2966 const unsigned int cache_size = cache_end - cache_begin;
2967 std::vector<types::global_dof_index> cache_indices(cache_size);
2968 for (
unsigned int i = 0; i < cache_size; ++i)
2969 cache_indices[i] = *(cache_begin + sorted_indices_pos[i]);
2971 IndexSet index_set(cache_indices.back() + 1);
2972 index_set.
add_indices(cache_indices.begin(), cache_indices.end());
2978 for (
unsigned int i = 0; i < cache_size; ++i, ++local_values_begin)
2979 *local_values_begin = read_write_vector[sorted_indices_pos[i]];
2987 template <
int dim,
int spacedim,
bool level_dof_access,
int structdim>
2990 const ::DoFAccessor<structdim, dim, spacedim, level_dof_access>
2993 const bool count_level_dofs)
2998 if (count_level_dofs)
3000 const auto &fe = accessor.get_fe(fe_index_);
3003 dofs_per_vertex = fe.n_dofs_per_vertex(),
3004 dofs_per_line = fe.n_dofs_per_line(),
3005 dofs_per_quad = fe.n_dofs_per_quad(0 ),
3006 dofs_per_hex = fe.n_dofs_per_hex();
3008 unsigned int index = 0;
3011 index += dofs_per_vertex * accessor.n_vertices();
3014 if (structdim == 2 || structdim == 3)
3015 index += dofs_per_line * accessor.n_lines();
3019 index += dofs_per_quad * accessor.n_faces();
3022 const unsigned int interior_dofs =
3023 structdim == 1 ? dofs_per_line :
3024 (structdim == 2 ? dofs_per_quad : dofs_per_hex);
3026 index += interior_dofs;
3032 const auto fe_index =
3034 accessor, fe_index_);
3036 unsigned int index = 0;
3039 for (
const auto vertex : accessor.vertex_indices())
3042 accessor.vertex_index(vertex),
3044 std::integral_constant<int, 0>())
3048 if (structdim == 2 || structdim == 3)
3049 for (
const auto line : accessor.line_indices())
3055 for (
const auto face : accessor.face_indices())
3072 template <
typename ArrayType>
3076 return array.size();
3085 template <
typename ArrayType>
3107 bool level_dof_access,
3109 typename DoFIndicesType,
3110 typename DoFOperation,
3111 typename DoFProcessor>
3114 const ::DoFAccessor<structdim, dim, spacedim, level_dof_access>
3116 const DoFIndicesType &const_dof_indices,
3118 const DoFOperation &dof_operation,
3119 const DoFProcessor &dof_processor,
3120 const bool count_level_dofs)
3124 accessor, fe_index_);
3129 (void)count_level_dofs;
3131 const auto &fe = accessor.get_fe(fe_index);
3145 if (fe.n_dofs_per_vertex() > 0)
3146 for (
const auto vertex : accessor.vertex_indices())
3147 dof_operation.process_vertex_dofs(*accessor.dof_handler,
3148 accessor.vertex_index(vertex),
3155 if (structdim > 1 && fe.n_dofs_per_line() > 0)
3157 const auto line_indices = internal::TriaAccessorImplementation::
3158 Implementation::get_line_indices_of_cell(accessor);
3159 const auto line_orientations =
3160 internal::TriaAccessorImplementation::Implementation::
3161 get_line_orientations_of_cell(accessor);
3163 for (
const auto line : accessor.line_indices())
3165 const auto line_orientation = line_orientations[line];
3167 dof_operation.process_dofs(
3168 accessor.get_dof_handler(),
3172 [](
const auto d) { return d; },
3173 std::integral_constant<int, 1>(),
3178 Assert(line_orientation ==
3181 dof_operation.process_dofs(
3182 accessor.get_dof_handler(),
3186 [&fe, line_orientation](
const auto d) {
3187 return fe.adjust_line_dof_index_for_line_orientation(
3188 d, line_orientation);
3190 std::integral_constant<int, 1>(),
3200 if (structdim == 3 && fe.max_dofs_per_quad() > 0)
3201 for (
const auto face_no : accessor.face_indices())
3203 const auto combined_orientation =
3204 accessor.combined_face_orientation(face_no);
3205 const unsigned int quad_index = accessor.quad_index(face_no);
3206 if (combined_orientation ==
3208 dof_operation.process_dofs(
3209 accessor.get_dof_handler(),
3213 [](
const auto d) { return d; },
3214 std::integral_constant<int, 2>(),
3218 dof_operation.process_dofs(
3219 accessor.get_dof_handler(),
3224 return fe.adjust_quad_dof_index_for_face_orientation(
3225 d, face_no, combined_orientation);
3227 std::integral_constant<int, 2>(),
3236 if (((dim == 3 && structdim == 2) ?
3237 fe.max_dofs_per_quad() :
3238 fe.template n_dofs_per_object<structdim>()) > 0)
3239 dof_operation.process_dofs(
3240 accessor.get_dof_handler(),
3244 [&](
const auto d) { return d; },
3245 std::integral_constant<int, structdim>(),
3249 if (dof_indices_ptr !=
nullptr)
3262 for (; dof_indices_ptr < end_dof_indices; ++dof_indices_ptr)
3263 dof_processor(invalid_index, dof_indices_ptr);
3272 template <
int dim,
int spacedim>
3278 template <
typename DoFProcessor>
3281 const unsigned int vertex_index,
3284 const DoFProcessor &dof_processor)
const
3295 std::integral_constant<int, 0>(),
3303 template <
int structdim,
typename DoFMapping,
typename DoFProcessor>
3306 const unsigned int obj_level,
3307 const unsigned int obj_index,
3310 const std::integral_constant<int, structdim>,
3312 const DoFProcessor &dof_processor)
const
3320 std::integral_constant<
int,
std::min(structdim, dim)>(),
3332 template <
int dim,
int spacedim>
3345 template <
typename DoFProcessor>
3348 const unsigned int vertex_index,
3351 const DoFProcessor &dof_processor)
const
3353 const unsigned int n_indices =
3354 dof_handler.
get_fe(0).template n_dofs_per_object<0>();
3359 for (
unsigned int d = 0; d < n_indices; ++d, ++dof_indices_ptr)
3360 dof_processor(stored_indices[d], dof_indices_ptr);
3366 template <
int structdim,
typename DoFMapping,
typename DoFProcessor>
3370 const unsigned int obj_index,
3373 const std::integral_constant<int, structdim>,
3375 const DoFProcessor &dof_processor)
const
3377 const unsigned int n_indices =
3378 dof_handler.
get_fe(0).template n_dofs_per_object<structdim>();
3386 std::integral_constant<
int,
std::min(structdim, dim)>());
3387 for (
unsigned int d = 0; d < n_indices; ++d, ++dof_indices_ptr)
3388 dof_processor(stored_indices[structdim < dim ?
mapping(d) : d],
3398 template <
int dim,
int spacedim,
bool level_dof_access,
int structdim>
3401 const ::DoFAccessor<structdim, dim, spacedim, level_dof_access>
3403 std::vector<types::global_dof_index> &dof_indices,
3411 [](
auto stored_index,
auto dof_ptr) { *dof_ptr = stored_index; },
3417 template <
int dim,
int spacedim,
bool level_dof_access,
int structdim>
3420 const ::DoFAccessor<structdim, dim, spacedim, level_dof_access>
3422 const std::vector<types::global_dof_index> &dof_indices,
3433 "This function is intended to be used for DoFCellAccessor, i.e., "
3434 "dimension == structdim."));
3441 [](
auto &stored_index,
auto dof_ptr) { stored_index = *dof_ptr; },
3447 template <
int dim,
int spacedim,
bool level_dof_access,
int structdim>
3450 const ::DoFAccessor<structdim, dim, spacedim, level_dof_access>
3453 std::vector<types::global_dof_index> &dof_indices,
3457 ExcMessage(
"MG DoF indices cannot be queried in hp case"));
3463 [](
auto stored_index,
auto dof_ptr) { *dof_ptr = stored_index; },
3469 template <
int dim,
int spacedim,
bool level_dof_access,
int structdim>
3472 const ::DoFAccessor<structdim, dim, spacedim, level_dof_access>
3475 const std::vector<types::global_dof_index> &dof_indices,
3479 ExcMessage(
"MG DoF indices cannot be queried in hp case"));
3487 ExcMessage(
"This function is intended to be used for "
3488 "DoFCellAccessor, i.e., dimension == structdim."));
3495 [](
auto &stored_index,
auto dof_ptr) { stored_index = *dof_ptr; },
3501 template <
int dim,
int spacedim>
3509 const unsigned int obj_index,
3511 const unsigned int local_index,
3512 const std::integral_constant<int, dim>)
3517 return mg_level->dof_object.access_dof_index(
3526 template <
int dim,
int spacedim, std::enable_if_t<(dim > 1),
int> = 0>
3534 const unsigned int obj_index,
3536 const unsigned int local_index,
3537 const std::integral_constant<int, 1>)
3539 return mg_faces->lines.access_dof_index(
3548 template <
int spacedim>
3556 const unsigned int obj_index,
3558 const unsigned int local_index,
3559 const std::integral_constant<int, 2>)
3563 return mg_faces->quads.access_dof_index(
3573 template <
int dim,
int spacedim,
bool level_dof_access>
3576 const ::DoFCellAccessor<dim, spacedim, level_dof_access> &accessor,
3578 const unsigned int fe_index);
3584template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
3587 const unsigned int i,
3590 const auto fe_index =
3595 return ::internal::DoFAccessorImplementation::Implementation::
3601 std::integral_constant<int, structdim>());
3605template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
3609 const unsigned int i)
const
3614 this->dof_handler->mg_faces,
3618 std::integral_constant<int, structdim>());
3622template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
3625 const unsigned int i,
3629 const auto fe_index =
3640 std::integral_constant<int, structdim>(),
3646template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
3650 const unsigned int i,
3656 this->dof_handler->mg_faces,
3660 std::integral_constant<int, structdim>()) = index;
3665template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
3671 return ::internal::DoFAccessorImplementation::Implementation::
3675 std::integral_constant<int, structdim>());
3680template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
3683 const unsigned int n)
const
3686 return ::internal::DoFAccessorImplementation::Implementation::
3691 std::integral_constant<int, structdim>());
3696template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
3697inline std::set<types::fe_index>
3701 std::set<types::fe_index> active_fe_indices;
3704 return active_fe_indices;
3709template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
3715 return ::internal::DoFAccessorImplementation::Implementation::
3720 std::integral_constant<int, structdim>());
3725template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
3728 const unsigned int vertex,
3729 const unsigned int i,
3733 (((this->
dof_handler->hp_capability_enabled ==
false) &&
3747 return ::internal::DoFAccessorImplementation::Implementation::
3750 this->vertex_index(vertex),
3753 std::integral_constant<int, 0>());
3757template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
3761 const unsigned int vertex,
3762 const unsigned int i,
3765 const auto fe_index =
3771 ExcMessage(
"Multigrid DoF indices can only be accessed after "
3772 "DoFHandler::distribute_mg_dofs() has been called!"));
3778 "DoFHandler in hp-mode does not implement multilevel DoFs."));
3780 return this->
dof_handler->mg_vertex_dofs[this->vertex_index(vertex)]
3781 .access_index(level, i, this->
dof_handler->get_fe().n_dofs_per_vertex());
3786template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
3790 const unsigned int vertex,
3791 const unsigned int i,
3795 const auto fe_index =
3805 "DoFHandler in hp-mode does not implement multilevel DoFs."));
3807 this->
dof_handler->mg_vertex_dofs[this->vertex_index(vertex)].access_index(
3808 level, i, this->
dof_handler->get_fe().n_dofs_per_vertex()) = index;
3813template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
3819 ExcMessage(
"This function can only be called for active FE indices"));
3826template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
3827inline typename ::internal::DoFHandlerImplementation::
3828 Iterators<dim, spacedim, level_dof_access>::line_iterator
3830 const unsigned int i)
const
3838 ExcMessage(
"You can only ask for line zero if the "
3839 "current object is a line itself."));
3840 return typename ::internal::DoFHandlerImplementation::
3841 Iterators<dim, spacedim, level_dof_access>::cell_iterator(
3842 &this->get_triangulation(),
3853 return typename ::internal::DoFHandlerImplementation::
3854 Iterators<dim, spacedim, level_dof_access>::line_iterator(
3857 this->line_index(i),
3862template <
int structdim,
int dim,
int spacedim,
bool level_dof_access>
3863inline typename ::internal::DoFHandlerImplementation::
3864 Iterators<dim, spacedim, level_dof_access>::quad_iterator
3866 const unsigned int i)
const
3876 ExcMessage(
"You can only ask for quad zero if the "
3877 "current object is a quad itself."));
3878 return typename ::internal::DoFHandlerImplementation::
3879 Iterators<dim, spacedim>::cell_iterator(&this->get_triangulation(),
3890 return typename ::internal::DoFHandlerImplementation::
3891 Iterators<dim, spacedim, level_dof_access>::quad_iterator(
3894 this->quad_index(i),
3902template <
int spacedim,
bool level_dof_access>
3910template <
int spacedim,
bool level_dof_access>
3922template <
int spacedim,
bool level_dof_access>
3947 "This constructor can not be called for face iterators in 1d, "
3948 "except to default-construct iterator objects."));
3953template <
int spacedim,
bool level_dof_access>
3954template <
int structdim2,
int dim2,
int spacedim2>
3963template <
int spacedim,
bool level_dof_access>
3964template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
3973template <
int spacedim,
bool level_dof_access>
3984template <
int spacedim,
bool level_dof_access>
3987 const unsigned int ,
3996template <
int spacedim,
bool level_dof_access>
4005template <
int spacedim,
bool level_dof_access>
4008 std::vector<types::global_dof_index> &dof_indices,
4011 const auto fe_index =
4015 for (
unsigned int i = 0; i < dof_indices.size(); ++i)
4019 this->global_vertex_index,
4022 std::integral_constant<int, 0>());
4027template <
int spacedim,
bool level_dof_access>
4031 std::vector<types::global_dof_index> &dof_indices,
4034 const auto fe_index =
4040 for (
unsigned int i = 0; i < dof_indices.size(); ++i)
4048template <
int spacedim,
bool level_dof_access>
4051 const unsigned int vertex,
4052 const unsigned int i,
4055 const auto fe_index =
4061 return ::internal::DoFAccessorImplementation::Implementation::
4064 this->global_vertex_index,
4067 std::integral_constant<int, 0>());
4072template <
int spacedim,
bool level_dof_access>
4075 const unsigned int i,
4078 const auto fe_index =
4082 return ::internal::DoFAccessorImplementation::Implementation::
4088 std::integral_constant<int, 0>());
4093template <
int spacedim,
bool level_dof_access>
4102template <
int spacedim,
bool level_dof_access>
4105 const unsigned int )
const
4112template <
int spacedim,
bool level_dof_access>
4123template <
int spacedim,
bool level_dof_access>
4134template <
int spacedim,
bool level_dof_access>
4140 BaseClass::copy_from(da);
4145template <
int spacedim,
bool level_dof_access>
4146template <
bool level_dof_access2>
4151 BaseClass::copy_from(a);
4157template <
int spacedim,
bool level_dof_access>
4160 const unsigned int )
const
4167template <
int spacedim,
bool level_dof_access>
4168inline typename ::internal::DoFHandlerImplementation::
4169 Iterators<1, spacedim, level_dof_access>::line_iterator
4171 const unsigned int )
const
4174 return typename ::internal::DoFHandlerImplementation::
4175 Iterators<1, spacedim, level_dof_access>::line_iterator();
4180template <
int spacedim,
bool level_dof_access>
4181inline typename ::internal::DoFHandlerImplementation::
4182 Iterators<1, spacedim, level_dof_access>::quad_iterator
4184 const unsigned int )
const
4187 return typename ::internal::DoFHandlerImplementation::
4188 Iterators<1, spacedim, level_dof_access>::quad_iterator();
4193template <
int spacedim,
bool level_dof_access>
4194template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
4201 return (BaseClass::operator==(a));
4206template <
int spacedim,
bool level_dof_access>
4207template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
4214 return (BaseClass::operator!=(a));
4236 template <
int dim,
int spacedim,
bool level_dof_access>
4245 (
typename std::decay_t<
decltype(accessor)>::ExcInvalidObject()));
4246 Assert(
static_cast<unsigned int>(accessor.level()) <
4260 template <
int dim,
int spacedim,
bool level_dof_access>
4273 (
typename std::decay_t<
decltype(accessor)>::ExcInvalidObject()));
4274 Assert(
static_cast<unsigned int>(accessor.level()) <
4278 ExcMessage(
"Invalid finite element index."));
4282 [accessor.present_index] = i;
4291 template <
int dim,
int spacedim,
bool level_dof_access>
4300 (
typename std::decay_t<
decltype(accessor)>::ExcInvalidObject()));
4301 Assert(
static_cast<unsigned int>(accessor.level()) <
4308 [accessor.present_index];
4312 [accessor.present_index];
4320 template <
int dim,
int spacedim,
bool level_dof_access>
4333 (
typename std::decay_t<
decltype(accessor)>::ExcInvalidObject()));
4334 Assert(
static_cast<unsigned int>(accessor.level()) <
4338 ExcMessage(
"Invalid finite element index."));
4342 [accessor.present_index] = i;
4351 template <
int dim,
int spacedim,
bool level_dof_access>
4360 (
typename std::decay_t<
decltype(accessor)>::ExcInvalidObject()));
4361 Assert(
static_cast<unsigned int>(accessor.level()) <
4367 [accessor.present_index] !=
4377 template <
int dim,
int spacedim,
bool level_dof_access>
4386 (
typename std::decay_t<
decltype(accessor)>::ExcInvalidObject()));
4387 Assert(
static_cast<unsigned int>(accessor.level()) <
4393 [accessor.present_index] =
4402template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4408 :
DoFAccessor<dimension_, dimension_, space_dimension_, level_dof_access>(
4417template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4418template <
int structdim2,
int dim2,
int spacedim2>
4427template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4428template <
int structdim2,
int dim2,
int spacedim2,
bool level_dof_access2>
4437template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4441 const unsigned int i)
const
4445 this->neighbor_level(i),
4446 this->neighbor_index(i),
4459template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4463 const unsigned int i)
const
4466 q(this->tria, this->level() + 1, this->child_index(i), this->
dof_handler);
4478template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4479inline boost::container::small_vector<
4485 boost::container::small_vector<
4491 for (
unsigned int i = 0; i < this->n_children(); ++i)
4499template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4505 q(this->tria, this->level() - 1, this->parent_index(), this->
dof_handler);
4514 namespace DoFCellAccessorImplementation
4516 template <
int dim,
int spacedim,
bool level_dof_access>
4520 const ::DoFCellAccessor<dim, spacedim, level_dof_access> &cell,
4521 const unsigned int i,
4522 const std::integral_constant<int, 1>)
4525 &cell.get_triangulation(),
4526 ((i == 0) && cell.at_boundary(0) ?
4528 ((i == 1) && cell.at_boundary(1) ?
4531 cell.vertex_index(i),
4532 &cell.get_dof_handler());
4533 return ::TriaIterator<
4538 template <
int dim,
int spacedim,
bool level_dof_access>
4542 const ::DoFCellAccessor<dim, spacedim, level_dof_access> &cell,
4543 const unsigned int i,
4544 const std::integral_constant<int, 2>)
4546 return cell.
line(i);
4550 template <
int dim,
int spacedim,
bool level_dof_access>
4554 const ::DoFCellAccessor<dim, spacedim, level_dof_access> &cell,
4555 const unsigned int i,
4556 const std::integral_constant<int, 3>)
4558 return cell.
quad(i);
4565template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4568 level_dof_access>::face_iterator
4570 const unsigned int i)
const
4574 return ::internal::DoFCellAccessorImplementation::get_face(
4575 *
this, i, std::integral_constant<int, dimension_>());
4580template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4581inline boost::container::small_vector<
4588 boost::container::small_vector<
4594 for (
const unsigned int i : this->face_indices())
4597 *
this, i, std::integral_constant<int, dimension_>());
4604template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4608 std::vector<types::global_dof_index> &dof_indices)
const
4610 if (level_dof_access)
4618template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4619template <
class InputVector,
typename number>
4622 const InputVector &values,
4630template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4631template <
typename Number,
typename ForwardIterator>
4635 ForwardIterator local_values_begin,
4636 ForwardIterator local_values_end)
const
4638 (void)local_values_end;
4639 Assert(this->is_artificial() ==
false,
4640 ExcMessage(
"Can't ask for DoF indices on artificial cells."));
4644 Assert(
static_cast<unsigned int>(local_values_end - local_values_begin) ==
4645 this->
get_fe().n_dofs_per_cell(),
4647 Assert(values.size() == this->get_dof_handler().n_dofs(),
4651 dof_indices(this->
get_fe().n_dofs_per_cell());
4655 boost::container::small_vector<Number, 27> values_temp(local_values_end -
4656 local_values_begin);
4661 using view_type = std::remove_reference_t<
decltype(*local_values_begin)>;
4663 local_values_end - local_values_begin);
4664 std::copy(values_temp.begin(), values_temp.end(), values_view2.
begin());
4669template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4670template <
class InputVector,
typename ForwardIterator>
4674 const InputVector &values,
4675 ForwardIterator local_values_begin,
4676 ForwardIterator local_values_end)
const
4678 Assert(this->is_artificial() ==
false,
4679 ExcMessage(
"Can't ask for DoF indices on artificial cells."));
4682 Assert(
static_cast<unsigned int>(local_values_end - local_values_begin) ==
4683 this->
get_fe().n_dofs_per_cell(),
4685 Assert(values.size() == this->get_dof_handler().n_dofs(),
4690 dof_indices(this->
get_fe().n_dofs_per_cell());
4702template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4703template <
class OutputVector,
typename number>
4707 OutputVector &values)
const
4709 Assert(this->is_artificial() ==
false,
4710 ExcMessage(
"Can't ask for DoF indices on artificial cells."));
4713 Assert(
static_cast<unsigned int>(local_values.
size()) ==
4714 this->get_fe().n_dofs_per_cell(),
4716 Assert(values.size() == this->get_dof_handler().n_dofs(),
4722 dof_indices(this->
get_fe().n_dofs_per_cell());
4726 for (
unsigned int i = 0; i < this->
get_fe().n_dofs_per_cell(); ++i)
4734template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4742 "For DoFHandler objects in hp-mode, finite elements are only "
4743 "associated with active cells. Consequently, you can not ask "
4744 "for the active finite element on cells with children."));
4748 Assert(this->reference_cell() == fe.reference_cell(),
4750 fe.reference_cell()));
4757template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4765 "You can not ask for the active FE index on a cell that has "
4766 "children because no degrees of freedom are assigned "
4767 "to this cell and, consequently, no finite element "
4768 "is associated with it."));
4770 (this->is_locally_owned() || this->is_ghost()),
4771 ExcMessage(
"You can only query active FE index information on cells "
4772 "that are either locally owned or (after distributing "
4773 "degrees of freedom) are ghost cells."));
4775 return ::internal::DoFCellAccessorImplementation::Implementation::
4776 active_fe_index(*
this);
4781template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4788 ExcMessage(
"You can not set the active FE index on a cell that has "
4789 "children because no degrees of freedom will be assigned "
4793 this->is_locally_owned(),
4794 ExcMessage(
"You can only set active FE index information on cells "
4795 "that are locally owned. On ghost cells, this information "
4796 "will automatically be propagated from the owning process "
4797 "of that cell, and there is no information at all on "
4798 "artificial cells."));
4806template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4815 "For DoFHandler objects in hp-mode, finite elements are only "
4816 "associated with active cells. Consequently, you can not ask "
4817 "for the future finite element on cells with children."));
4824template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4830 (this->has_children() ==
false),
4832 "You can not ask for the future FE index on a cell that has "
4833 "children because no degrees of freedom are assigned "
4834 "to this cell and, consequently, no finite element "
4835 "is associated with it."));
4837 (this->is_locally_owned()),
4838 ExcMessage(
"You can only query future FE index information on cells "
4839 "that are locally owned."));
4841 return ::internal::DoFCellAccessorImplementation::Implementation::
4842 future_fe_index(*
this);
4847template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4853 (this->has_children() ==
false),
4854 ExcMessage(
"You can not set the future FE index on a cell that has "
4855 "children because no degrees of freedom will be assigned "
4859 this->is_locally_owned(),
4860 ExcMessage(
"You can only set future FE index information on cells "
4861 "that are locally owned."));
4869template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4875 (this->has_children() ==
false),
4877 "You can not ask for the future FE index on a cell that has "
4878 "children because no degrees of freedom are assigned "
4879 "to this cell and, consequently, no finite element "
4880 "is associated with it."));
4882 (this->is_locally_owned()),
4883 ExcMessage(
"You can only query future FE index information on cells "
4884 "that are locally owned."));
4886 return ::internal::DoFCellAccessorImplementation::Implementation::
4887 future_fe_index_set(*
this);
4892template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4898 (this->has_children() ==
false),
4900 "You can not ask for the future FE index on a cell that has "
4901 "children because no degrees of freedom are assigned "
4902 "to this cell and, consequently, no finite element "
4903 "is associated with it."));
4905 (this->is_locally_owned()),
4906 ExcMessage(
"You can only query future FE index information on cells "
4907 "that are locally owned."));
4915template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4916template <
typename number,
typename OutputVector>
4920 OutputVector &global_destination)
const
4924 global_destination);
4929template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4930template <
typename ForwardIterator,
typename OutputVector>
4934 ForwardIterator local_source_end,
4935 OutputVector &global_destination)
const
4939 Assert(
static_cast<unsigned int>(local_source_end - local_source_begin) ==
4940 this->
get_fe().n_dofs_per_cell(),
4948 internal::ArrayViewHelper::is_contiguous(local_source_begin,
4951 "This function can not be called with iterator types that do not point to contiguous memory."));
4953 const unsigned int n_dofs = local_source_end - local_source_begin;
4956 dof_indices(n_dofs);
4961 global_destination.add(n_dofs, dof_indices.data(), &(*local_source_begin));
4966template <
int dimension_,
int space_dimension_,
bool level_dof_access>
4967template <
typename ForwardIterator,
typename OutputVector>
4972 ForwardIterator local_source_begin,
4973 ForwardIterator local_source_end,
4974 OutputVector &global_destination)
const
4978 Assert(local_source_end - local_source_begin ==
4979 this->
get_fe().n_dofs_per_cell(),
4987 dof_indices(this->
get_fe().n_dofs_per_cell());
4995 global_destination);
5000template <
int dimension_,
int space_dimension_,
bool level_dof_access>
5001template <
typename number,
typename OutputMatrix>
5005 OutputMatrix &global_destination)
const
5009 Assert(local_source.
m() == this->get_fe().n_dofs_per_cell(),
5011 Assert(local_source.
n() == this->get_fe().n_dofs_per_cell(),
5020 const unsigned int n_dofs = local_source.
m();
5023 dof_indices(n_dofs);
5028 for (
unsigned int i = 0; i < n_dofs; ++i)
5029 global_destination.add(dof_indices[i],
5032 &local_source(i, 0));
5037template <
int dimension_,
int space_dimension_,
bool level_dof_access>
5038template <
typename number,
typename OutputMatrix,
typename OutputVector>
5043 OutputMatrix &global_matrix,
5044 OutputVector &global_vector)
const
5048 Assert(local_matrix.
m() == this->get_fe().n_dofs_per_cell(),
5050 Assert(local_matrix.
n() == this->get_fe().n_dofs_per_cell(),
5056 Assert(local_vector.
size() == this->get_fe().n_dofs_per_cell(),
5063 const unsigned int n_dofs = this->
get_fe().n_dofs_per_cell();
5065 dof_indices(n_dofs);
5070 for (
unsigned int i = 0; i < n_dofs; ++i)
5072 global_matrix.add(dof_indices[i],
5075 &local_matrix(i, 0));
5076 global_vector(dof_indices[i]) += local_vector(i);
ArrayView< std::remove_reference_t< typename std::iterator_traits< Iterator >::reference >, MemorySpaceType > make_array_view(const Iterator begin, const Iterator end)
void get_dof_values(const VectorType &global_vector, ForwardIteratorInd local_indices_begin, ForwardIteratorVec local_vector_begin, ForwardIteratorVec local_vector_end) const
void distribute_local_to_global(const InVector &local_vector, const std::vector< size_type > &local_dof_indices, OutVector &global_vector) const
friend struct ::internal::DoFHandlerImplementation::Policy::Implementation
void get_mg_dof_indices(const int level, std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index) const
DoFAccessor< 0, 1, spacedim, level_dof_access > & operator=(const DoFAccessor< 0, 1, spacedim, level_dof_access > &da)=delete
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index) const
typename::internal::DoFHandlerImplementation::Iterators< 1, spacedim, level_dof_access >::line_iterator line(const unsigned int i) const
DoFAccessor< 0, 1, spacedim, level_dof_access > & operator=(DoFAccessor< 0, 1, spacedim, level_dof_access > &&) noexcept=default
types::global_dof_index vertex_dof_index(const unsigned int vertex, const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
TriaAccessor< 0, 1, spacedim > BaseClass
DoFAccessor(DoFAccessor< 0, 1, spacedim, level_dof_access > &&)=default
const DoFHandler< 1, spacedim > & get_dof_handler() const
DoFHandler< 1, spacedim > * dof_handler
types::global_dof_index dof_index(const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
typename::internal::DoFHandlerImplementation::Iterators< 1, spacedim, level_dof_access >::quad_iterator quad(const unsigned int i) const
DoFAccessor(const DoFAccessor< 0, 1, spacedim, level_dof_access > &)=default
types::fe_index nth_active_fe_index(const unsigned int n) const
const FiniteElement< 1, spacedim > & get_fe(const types::fe_index fe_index) const
friend class TriaRawIterator
void set_dof_index(const unsigned int i, const types::global_dof_index index, const types::fe_index fe_index=numbers::invalid_fe_index) const
static constexpr unsigned int space_dimension
void set_dof_handler(DoFHandler< 1, spacedim > *dh)
unsigned int n_active_fe_indices() const
DoFHandler< 1, spacedim > AccessorData
TriaIterator< DoFAccessor< 0, 1, spacedim, level_dof_access > > child(const unsigned int c) const
bool fe_index_is_active(const types::fe_index fe_index) const
void copy_from(const DoFAccessor< 0, 1, spacedim, level_dof_access2 > &a)
static constexpr unsigned int dimension
typename::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::line_iterator line(const unsigned int i) const
const DoFHandler< dim, spacedim > & get_dof_handler() const
void set_mg_dof_index(const int level, const unsigned int i, const types::global_dof_index index) const
static constexpr unsigned int space_dimension
DoFAccessor< structdim, dim, spacedim, level_dof_access > & operator=(const DoFAccessor< structdim, dim, spacedim, level_dof_access > &da)=delete
types::fe_index nth_active_fe_index(const unsigned int n) const
DoFAccessor(DoFAccessor< structdim, dim, spacedim, level_dof_access > &&)
DoFAccessor(const DoFAccessor< structdim, dim, spacedim, level_dof_access2 > &)
static constexpr unsigned int dimension
TriaIterator< DoFAccessor< structdim, dim, spacedim, level_dof_access > > child(const unsigned int c) const
void set_mg_vertex_dof_index(const int level, const unsigned int vertex, const unsigned int i, const types::global_dof_index index, const types::fe_index fe_index=numbers::invalid_fe_index) const
void set_dof_handler(DoFHandler< dim, spacedim > *dh)
DoFAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
DoFAccessor< structdim, dim, spacedim, level_dof_access > & operator=(DoFAccessor< structdim, dim, spacedim, level_dof_access > &&)
std::set< types::fe_index > get_active_fe_indices() const
DoFAccessor(const DoFAccessor< structdim, dim, spacedim, level_dof_access > &)=default
void set_mg_dof_indices(const int level, const std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index)
void get_mg_dof_indices(const int level, std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index) const
types::global_dof_index mg_vertex_dof_index(const int level, const unsigned int vertex, const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
typename ::internal::DoFAccessorImplementation:: Inheritance< structdim, dimension, space_dimension >::BaseClass BaseClass
DoFAccessor(const Triangulation< dim, spacedim > *tria, const int level, const int index, const DoFHandler< dim, spacedim > *dof_handler)
bool operator!=(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &) const
void copy_from(const TriaAccessorBase< structdim, dim, spacedim > &da)
types::global_dof_index dof_index(const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
typename::internal::DoFHandlerImplementation::Iterators< dim, spacedim, level_dof_access >::quad_iterator quad(const unsigned int i) const
friend class TriaRawIterator
DoFHandler< dimension, space_dimension > AccessorData
void copy_from(const DoFAccessor< structdim, dim, spacedim, level_dof_access2 > &a)
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index fe_index) const
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index=numbers::invalid_fe_index) const
bool operator==(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &) const
types::global_dof_index mg_dof_index(const int level, const unsigned int i) const
void set_dof_index(const unsigned int i, const types::global_dof_index index, const types::fe_index fe_index=numbers::invalid_fe_index) const
unsigned int n_active_fe_indices() const
static bool is_level_cell()
types::global_dof_index vertex_dof_index(const unsigned int vertex, const unsigned int i, const types::fe_index fe_index=numbers::invalid_fe_index) const
DoFHandler< dim, spacedim > * dof_handler
DoFAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
bool fe_index_is_active(const types::fe_index fe_index) const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > periodic_neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
boost::container::small_vector< TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > >, GeometryInfo< dimension_ >::max_children_per_cell > child_iterators() const
DoFCellAccessor(DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &&)=default
DoFHandler< dimension_, space_dimension_ > AccessorData
void get_active_or_mg_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
void get_dof_values(const InputVector &values, Vector< number > &local_values) const
boost::container::small_vector< face_iterator, GeometryInfo< dimension_ >::faces_per_cell > face_iterators() const
types::fe_index future_fe_index() const
DoFCellAccessor< dimension_, space_dimension_, level_dof_access > & operator=(const DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &da)=delete
const FiniteElement< dimension_, space_dimension_ > & get_fe() const
DoFAccessor< dimension_, dimension_, space_dimension_, level_dof_access > BaseClass
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > child(const unsigned int i) const
TriaIterator< DoFAccessor< dimension_ - 1, dimension_, space_dimension_, level_dof_access > > face_iterator
DoFCellAccessor(const Triangulation< dimension_, space_dimension_ > *tria, const int level, const int index, const AccessorData *local_data)
void set_future_fe_index(const types::fe_index i) const
void distribute_local_to_global(const Vector< number > &local_source, OutputVector &global_destination) const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > neighbor_or_periodic_neighbor(const unsigned int i) const
const FiniteElement< dimension_, space_dimension_ > & get_future_fe() const
void get_dof_values(const AffineConstraints< typename InputVector::value_type > &constraints, const InputVector &values, ForwardIterator local_values_begin, ForwardIterator local_values_end) const
void distribute_local_to_global(const FullMatrix< number > &local_source, OutputMatrix &global_destination) const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > neighbor(const unsigned int i) const
static const unsigned int dim
void set_active_fe_index(const types::fe_index i) const
void distribute_local_to_global_by_interpolation(const Vector< number > &local_values, OutputVector &values, const types::fe_index fe_index=numbers::invalid_fe_index) const
DoFCellAccessor(const DoFAccessor< structdim2, dim2, spacedim2, level_dof_access2 > &)
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > periodic_neighbor(const unsigned int i) const
DoFCellAccessor(const InvalidAccessor< structdim2, dim2, spacedim2 > &)
static const unsigned int spacedim
void get_interpolated_dof_values(const ReadVector< Number > &values, Vector< Number > &interpolated_values, const types::fe_index fe_index=numbers::invalid_fe_index) const
face_iterator face(const unsigned int i) const
void clear_future_fe_index() const
void distribute_local_to_global(const FullMatrix< number > &local_matrix, const Vector< number > &local_vector, OutputMatrix &global_matrix, OutputVector &global_vector) const
void set_dof_indices(const std::vector< types::global_dof_index > &dof_indices)
void distribute_local_to_global(const AffineConstraints< typename OutputVector::value_type > &constraints, ForwardIterator local_source_begin, ForwardIterator local_source_end, OutputVector &global_destination) const
void set_dof_values(const Vector< number > &local_values, OutputVector &values) const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > neighbor_child_on_subface(const unsigned int face_no, const unsigned int subface_no) const
void set_dof_values_by_interpolation(const Vector< number > &local_values, OutputVector &values, const types::fe_index fe_index=numbers::invalid_fe_index, const bool perform_check=false) const
void set_mg_dof_indices(const std::vector< types::global_dof_index > &dof_indices)
DoFCellAccessor(const DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &)=default
void get_mg_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
void distribute_local_to_global(ForwardIterator local_source_begin, ForwardIterator local_source_end, OutputVector &global_destination) const
void get_dof_values(const ReadVector< Number > &values, ForwardIterator local_values_begin, ForwardIterator local_values_end) const
~DoFCellAccessor()=default
DoFHandler< dimension_, space_dimension_ > Container
bool future_fe_index_set() const
DoFCellAccessor< dimension_, space_dimension_, level_dof_access > & operator=(DoFCellAccessor< dimension_, space_dimension_, level_dof_access > &&)=default
void get_dof_indices(std::vector< types::global_dof_index > &dof_indices) const
types::fe_index active_fe_index() const
TriaIterator< DoFCellAccessor< dimension_, space_dimension_, level_dof_access > > parent() const
std::vector< std::unique_ptr<::internal::DoFHandlerImplementation::DoFLevel< dim > > > mg_levels
static const types::fe_index default_fe_index
std::vector< std::vector< types::fe_index > > hp_cell_future_fe_indices
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
std::vector< std::vector< types::fe_index > > hp_cell_active_fe_indices
std::vector< MGVertexDoFs > mg_vertex_dofs
std::vector< std::array< std::vector< types::global_dof_index >, dim+1 > > object_dof_indices
std::vector< std::array< std::vector< offset_type >, dim+1 > > object_dof_ptr
std::unique_ptr<::internal::DoFHandlerImplementation::DoFFaces< dim > > mg_faces
bool hp_capability_enabled
bool has_hp_capabilities() const
std::array< std::vector< offset_type >, dim+1 > hp_object_fe_ptr
std::array< std::vector< types::fe_index >, dim+1 > hp_object_fe_indices
typename InvalidAccessor< structdim, dim, spacedim >::AccessorData AccessorData
DoFInvalidAccessor(const void *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
unsigned int n_dofs_per_vertex() const
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
InvalidAccessor(const void *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
void import_elements(const ::Vector< Number > &vec, VectorOperation::values operation, const std::shared_ptr< const Utilities::MPI::CommunicationPatternBase > &communication_pattern={})
friend class TriaIterator
TriaAccessorBase(const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *=nullptr)
const Triangulation< dim, spacedim > * tria
TriaIterator< TriaAccessor< structdim, dim, spacedim > > child(const unsigned int i) const
unsigned int vertex_index(const unsigned int i) const
Point< spacedim > & vertex(const unsigned int i) const
TriaAccessor(const Triangulation< dim, spacedim > *parent=nullptr, const int level=-1, const int index=-1, const AccessorData *local_data=nullptr)
IteratorState::IteratorStates state() const
virtual size_type size() const override
#define DEAL_II_ALWAYS_INLINE
#define DEAL_II_NAMESPACE_OPEN
constexpr bool running_in_debug_mode()
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ASSERT_UNREACHABLE()
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcVectorNotEmpty()
#define DeclException0(Exception0)
static ::ExceptionBase & ExcVectorDoesNotMatch()
static ::ExceptionBase & ExcInvalidObject()
static ::ExceptionBase & ExcNonMatchingReferenceCellTypes(ReferenceCell arg1, ReferenceCell arg2)
static ::ExceptionBase & ExcVectorNotEmpty()
static ::ExceptionBase & ExcMatrixDoesNotMatch()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
static ::ExceptionBase & ExcNotImplementedWithHP()
static ::ExceptionBase & ExcCantCompareIterators()
static ::ExceptionBase & ExcNotActive()
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcMatrixDoesNotMatch()
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcInvalidObject()
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcNotActive()
static ::ExceptionBase & ExcCantCompareIterators()
static ::ExceptionBase & ExcVectorDoesNotMatch()
static ::ExceptionBase & ExcMessage(std::string arg1)
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
@ past_the_end
Iterator reached end of container.
types::fe_index get_fe_index_or_default(const DoFAccessor< structdim, dim, spacedim, level_dof_access > &cell, const types::fe_index fe_index)
void get_cell_dof_indices(const ::DoFCellAccessor< dim, spacedim, level_dof_access > &accessor, Implementation::dof_index_vector_type &dof_indices, const unsigned int fe_index)
TriaIterator< ::DoFAccessor< dim - 1, dim, spacedim, level_dof_access > > get_face(const ::DoFCellAccessor< dim, spacedim, level_dof_access > &cell, const unsigned int i, const std::integral_constant< int, 1 >)
constexpr types::global_dof_index invalid_dof_index
constexpr types::geometric_orientation reverse_line_orientation
constexpr types::geometric_orientation default_geometric_orientation
constexpr types::fe_index invalid_fe_index
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
unsigned short int fe_index
unsigned int global_dof_index
static constexpr unsigned int faces_per_cell
static constexpr unsigned int max_children_per_cell
void process_vertex_dofs(DoFHandler< dim, spacedim > &dof_handler, const unsigned int vertex_index, const types::fe_index fe_index, types::global_dof_index *&dof_indices_ptr, const DoFProcessor &dof_processor) const
void process_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const types::fe_index fe_index, const DoFMapping &mapping, const std::integral_constant< int, structdim >, types::global_dof_index *&dof_indices_ptr, const DoFProcessor &dof_processor) const
void process_dofs(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int, const unsigned int obj_index, const types::fe_index fe_index, const DoFMapping &mapping, const std::integral_constant< int, structdim >, types::global_dof_index *&dof_indices_ptr, const DoFProcessor &dof_processor) const
MGDoFIndexProcessor(const unsigned int level)
void process_vertex_dofs(DoFHandler< dim, spacedim > &dof_handler, const unsigned int vertex_index, const types::fe_index, types::global_dof_index *&dof_indices_ptr, const DoFProcessor &dof_processor) const
static types::global_dof_index * get_array_ptr(const std::tuple<> &)
static void process_dof_indices(const ::DoFAccessor< structdim, dim, spacedim, level_dof_access > &accessor, const DoFIndicesType &const_dof_indices, const types::fe_index fe_index_, const DoFOperation &dof_operation, const DoFProcessor &dof_processor, const bool count_level_dofs)
static void get_dof_indices(const ::DoFAccessor< structdim, dim, spacedim, level_dof_access > &accessor, std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index)
static void extract_subvector_to(const LinearAlgebra::TpetraWrappers::Vector< Number, MemorySpace > &values, const types::global_dof_index *cache_begin, const types::global_dof_index *cache_end, ForwardIterator local_values_begin)
static std::set< types::fe_index > get_active_fe_indices(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const std::integral_constant< int, structdim > &t)
static types::global_dof_index * get_array_ptr(const ArrayType &array)
static unsigned int get_array_length(const std::tuple<> &)
static unsigned int n_active_fe_indices(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const std::integral_constant< int, structdim > &)
static void set_dof_index(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const types::fe_index fe_index, const unsigned int local_index, const std::integral_constant< int, structdim > &dd, const types::global_dof_index global_index)
static void extract_subvector_to(const LinearAlgebra::EpetraWrappers::Vector &values, const types::global_dof_index *cache_begin, const types::global_dof_index *cache_end, ForwardIterator local_values_begin)
static void get_mg_dof_indices(const ::DoFAccessor< structdim, dim, spacedim, level_dof_access > &accessor, const int level, std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index)
static bool fe_index_is_active(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const types::fe_index fe_index, const std::integral_constant< int, structdim > &)
static types::global_dof_index & get_mg_dof_index(const DoFHandler< dim, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< dim > > &, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< dim > > &mg_faces, const unsigned int obj_index, const types::fe_index fe_index, const unsigned int local_index, const std::integral_constant< int, 1 >)
static types::global_dof_index & get_mg_dof_index(const DoFHandler< dim, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< dim > > &mg_level, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< dim > > &, const unsigned int obj_index, const types::fe_index fe_index, const unsigned int local_index, const std::integral_constant< int, dim >)
static std::pair< unsigned int, unsigned int > process_object_range(const ::DoFAccessor< structdim, dim, spacedim, level_dof_access > accessor, const types::fe_index fe_index)
static types::global_dof_index & mg_vertex_dof_index(DoFHandler< dim, spacedim > &dof_handler, const int level, const unsigned int vertex_index, const unsigned int i)
boost::container::small_vector<::types::global_dof_index, 27 > dof_index_vector_type
static void set_mg_dof_indices(const ::DoFAccessor< structdim, dim, spacedim, level_dof_access > &accessor, const int level, const std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index)
static void set_dof_indices(const ::DoFAccessor< structdim, dim, spacedim, level_dof_access > &accessor, const std::vector< types::global_dof_index > &dof_indices, const types::fe_index fe_index)
static void process_dof_index(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const types::fe_index fe_index, const unsigned int local_index, const std::integral_constant< int, structdim > &, GlobalIndexType &global_index, const DoFPProcessor &process)
static types::global_dof_index & get_mg_dof_index(const DoFHandler< 3, spacedim > &dof_handler, const std::unique_ptr< internal::DoFHandlerImplementation::DoFLevel< 3 > > &, const std::unique_ptr< internal::DoFHandlerImplementation::DoFFaces< 3 > > &mg_faces, const unsigned int obj_index, const types::fe_index fe_index, const unsigned int local_index, const std::integral_constant< int, 2 >)
static std::pair< unsigned int, unsigned int > process_object_range(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const types::fe_index fe_index, const std::integral_constant< int, structdim > &)
static void extract_subvector_to(const InputVector &values, const types::global_dof_index *cache, const types::global_dof_index *cache_end, ForwardIterator local_values_begin)
static std::vector< unsigned int > sort_indices(const types::global_dof_index *v_begin, const types::global_dof_index *v_end)
static types::fe_index nth_active_fe_index(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const unsigned int local_index, const std::integral_constant< int, structdim > &)
static void process_object(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const types::fe_index fe_index, const DoFMapping &mapping, const std::integral_constant< int, structdim > &dd, types::global_dof_index *&dof_indices_ptr, const DoFProcessor &process)
static types::global_dof_index get_dof_index(const DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const types::fe_index fe_index, const unsigned int local_index, const std::integral_constant< int, structdim > &dd)
static unsigned int n_dof_indices(const ::DoFAccessor< structdim, dim, spacedim, level_dof_access > &accessor, const types::fe_index fe_index_, const bool count_level_dofs)
static std::pair< unsigned int, unsigned int > process_object_range(::DoFInvalidAccessor< structdim, dim, spacedim >, const unsigned int)
static unsigned int get_array_length(const ArrayType &array)
::CellAccessor< dim, spacedim > BaseClass
::TriaAccessor< structdim, dim, spacedim > BaseClass
static types::fe_index future_fe_index(const DoFCellAccessor< dim, spacedim, level_dof_access > &accessor)
static void clear_future_fe_index(const DoFCellAccessor< dim, spacedim, level_dof_access > &accessor)
static bool future_fe_index_set(const DoFCellAccessor< dim, spacedim, level_dof_access > &accessor)
static void set_active_fe_index(const DoFCellAccessor< dim, spacedim, level_dof_access > &accessor, const types::fe_index i)
static void set_future_fe_index(const DoFCellAccessor< dim, spacedim, level_dof_access > &accessor, const types::fe_index i)
static types::fe_index active_fe_index(const DoFCellAccessor< dim, spacedim, level_dof_access > &accessor)
static void set(typename VectorType::value_type value, const types::global_dof_index i, VectorType &V)