16 #ifndef dealii__hp_dof_faces_h 17 #define dealii__hp_dof_faces_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/hp/fe_collection.h> 25 DEAL_II_NAMESPACE_OPEN
29 template <
int dim,
int spacedim>
100 template <
int structdim>
117 std::vector<types::global_dof_index>
dofs;
131 template <
int dim,
int spacedim>
133 set_dof_index (const ::hp::DoFHandler<dim,spacedim> &dof_handler,
134 const unsigned int obj_index,
135 const unsigned int fe_index,
136 const unsigned int local_index,
138 const unsigned int obj_level);
151 template <
int dim,
int spacedim>
153 get_dof_index (const ::hp::DoFHandler<dim,spacedim> &dof_handler,
154 const unsigned int obj_index,
155 const unsigned int fe_index,
156 const unsigned int local_index,
157 const unsigned int obj_level)
const;
170 template <
int dim,
int spacedim>
173 const unsigned int obj_index)
const;
178 template <
int dim,
int spacedim>
181 const unsigned int obj_level,
182 const unsigned int obj_index,
183 const unsigned int n)
const;
189 template <
int dim,
int spacedim>
192 const unsigned int obj_index,
193 const unsigned int fe_index,
194 const unsigned int obj_level)
const;
206 template <
class Archive>
208 const unsigned int version);
259 std::size_t memory_consumption ()
const;
265 template <
class Archive>
266 void serialize(Archive &ar,
267 const unsigned int version);
290 std::size_t memory_consumption ()
const;
296 template <
class Archive>
297 void serialize(Archive &ar,
298 const unsigned int version);
326 std::size_t memory_consumption ()
const;
332 template <
class Archive>
333 void serialize(Archive &ar,
334 const unsigned int version);
339 template <
class Archive>
345 template <
class Archive>
353 template <
class Archive>
360 template <
int structdim>
361 template <
int dim,
int spacedim>
366 const unsigned int obj_index,
367 const unsigned int fe_index,
368 const unsigned int local_index,
369 const unsigned int )
const 372 ExcMessage (
"You need to specify a FE index when working " 373 "with hp DoFHandlers"));
374 Assert (fe_index < dof_handler.get_fe().size(),
377 dof_handler.get_fe()[fe_index].template n_dofs_per_object<structdim>(),
379 dof_handler.get_fe()[fe_index]
380 .template n_dofs_per_object<structdim>()));
381 Assert (obj_index < dof_offsets.size(),
388 ExcMessage (
"You are trying to access degree of freedom " 389 "information for an object on which no such " 390 "information is available"));
392 Assert (structdim<dim,
ExcMessage (
"This object can not be used for cells."));
405 if (*pointer == fe_index)
406 return *(pointer + 1 + local_index);
409 dof_handler.get_fe()[*pointer]
410 .template n_dofs_per_object<structdim>() + 1);
416 template <
int structdim>
417 template <
int dim,
int spacedim>
422 const unsigned int obj_index,
423 const unsigned int fe_index,
424 const unsigned int local_index,
429 ExcMessage (
"You need to specify a FE index when working " 430 "with hp DoFHandlers"));
431 Assert (fe_index < dof_handler.get_fe().size(),
434 dof_handler.get_fe()[fe_index].template n_dofs_per_object<structdim>(),
436 dof_handler.get_fe()[fe_index]
437 .template n_dofs_per_object<structdim>()));
438 Assert (obj_index < dof_offsets.size(),
445 ExcMessage (
"You are trying to access degree of freedom " 446 "information for an object on which no such " 447 "information is available"));
449 Assert (structdim<dim,
ExcMessage (
"This object can not be used for cells."));
462 if (*pointer == fe_index)
464 *(pointer + 1 + local_index) = global_index;
468 pointer += dof_handler.get_fe()[*pointer]
469 .template n_dofs_per_object<structdim>() + 1;
475 template <
int structdim>
476 template <
int dim,
int spacedim>
481 const unsigned int obj_index)
const 483 Assert (obj_index < dof_offsets.size(),
492 Assert (structdim<dim,
ExcMessage (
"This object can not be used for cells."));
499 const unsigned int starting_offset = dof_offsets[obj_index];
501 unsigned int counter = 0;
510 pointer += dof_handler.get_fe()[*pointer]
511 .template n_dofs_per_object<structdim>() + 1;
518 template <
int structdim>
519 template <
int dim,
int spacedim>
525 const unsigned int obj_index,
526 const unsigned int n)
const 528 Assert (obj_index < dof_offsets.size(),
535 ExcMessage (
"You are trying to access degree of freedom " 536 "information for an object on which no such " 537 "information is available"));
539 Assert (structdim<dim,
ExcMessage (
"This object can not be used for cells."));
541 Assert (n < n_active_fe_indices(dof_handler, obj_index),
543 n_active_fe_indices(dof_handler, obj_index)));
550 const unsigned int starting_offset = dof_offsets[obj_index];
552 unsigned int counter = 0;
558 const unsigned int fe_index = *pointer;
560 Assert (fe_index < dof_handler.get_fe().size(),
567 pointer += dof_handler.get_fe()[fe_index]
568 .template n_dofs_per_object<structdim>() + 1;
574 template <
int structdim>
575 template <
int dim,
int spacedim>
580 const unsigned int obj_index,
581 const unsigned int fe_index,
582 const unsigned int )
const 584 Assert (obj_index < dof_offsets.size(),
585 ExcIndexRange (obj_index, 0, static_cast<unsigned int>(dof_offsets.size())));
587 ExcMessage (
"You need to specify a FE index when working " 588 "with hp DoFHandlers"));
589 Assert (fe_index < dof_handler.get_fe().size(),
596 ExcMessage (
"You are trying to access degree of freedom " 597 "information for an object on which no such " 598 "information is available"));
600 Assert (structdim<dim,
ExcMessage (
"This object can not be used for cells."));
614 else if (*pointer == fe_index)
618 dof_handler.get_fe()[*pointer]
619 .template n_dofs_per_object<structdim>()+1);
623 template <
int structdim>
624 template <
class Archive>
637 DEAL_II_NAMESPACE_CLOSE
static const unsigned int invalid_unsigned_int
void set_dof_index(const ::hp::DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const types::global_dof_index global_index, const unsigned int obj_level)
types::global_dof_index get_dof_index(const ::hp::DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_index, const unsigned int fe_index, const unsigned int local_index, const unsigned int obj_level) const
unsigned int n_active_fe_indices(const ::hp::DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_index) const
internal::hp::DoFIndicesOnFacesOrEdges< 1 > lines
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
static ::ExceptionBase & ExcMessage(std::string arg1)
unsigned int global_dof_index
#define Assert(cond, exc)
void serialize(Archive &ar, const unsigned int version)
std::size_t memory_consumption() const
std::vector< types::global_dof_index > dofs
std::vector< unsigned int > dof_offsets
internal::hp::DoFIndicesOnFacesOrEdges< 2 > quads
bool fe_index_is_active(const ::hp::DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_index, const unsigned int fe_index, const unsigned int obj_level) const
types::global_dof_index nth_active_fe_index(const ::hp::DoFHandler< dim, spacedim > &dof_handler, const unsigned int obj_level, const unsigned int obj_index, const unsigned int n) const
const types::global_dof_index invalid_dof_index
internal::hp::DoFIndicesOnFacesOrEdges< 1 > lines
static ::ExceptionBase & ExcInternalError()