16 #ifndef dealii_hp_dof_level_h 17 #define dealii_hp_dof_level_h 44 namespace DoFHandlerImplementation
46 struct Implementation;
49 namespace DoFCellAccessorImplementation
51 struct Implementation;
224 set_dof_index(
const unsigned int obj_index,
225 const unsigned int fe_index,
226 const unsigned int local_index,
241 get_dof_index(
const unsigned int obj_index,
242 const unsigned int fe_index,
243 const unsigned int local_index)
const;
249 active_fe_index(
const unsigned int obj_index)
const;
256 fe_index_is_active(
const unsigned int obj_index,
257 const unsigned int fe_index)
const;
263 set_active_fe_index(
const unsigned int obj_index,
264 const unsigned int fe_index);
271 future_fe_index(
const unsigned int obj_index)
const;
277 set_future_fe_index(
const unsigned int obj_index,
278 const unsigned int fe_index);
284 future_fe_index_set(
const unsigned int obj_index)
const;
290 clear_future_fe_index(
const unsigned int obj_index);
304 get_cell_cache_start(
const unsigned int obj_index,
305 const unsigned int dofs_per_cell)
const;
318 template <
class Archive>
320 serialize(Archive &ar,
const unsigned int version);
331 template <
int dim,
int spacedim>
334 const ::hp::FECollection<dim, spacedim> &fe_collection);
344 template <
int dim,
int spacedim>
347 const ::hp::FECollection<dim, spacedim> &fe_collection);
368 normalize_active_fe_indices();
374 friend class ::hp::DoFHandler;
375 friend struct ::internal::hp::DoFHandlerImplementation::
377 friend struct ::internal::DoFCellAccessorImplementation::
388 return (static_cast<signed_active_fe_index_type>(active_fe_index) < 0);
394 DoFLevel::get_toggled_compression_state(
406 DoFLevel::get_dof_index(
const unsigned int obj_index,
407 const unsigned int fe_index,
408 const unsigned int local_index)
const 415 Assert(dof_offsets[obj_index] != static_cast<offset_type>(-1),
416 ExcMessage(
"You are trying to access degree of freedom " 417 "information for an object on which no such " 418 "information is available"));
421 (is_compressed_entry(active_fe_indices[obj_index]) ==
false ?
422 active_fe_indices[obj_index] :
423 get_toggled_compression_state(active_fe_indices[obj_index])),
424 ExcMessage(
"FE index does not match that of the present cell"));
428 if (is_compressed_entry(active_fe_indices[obj_index]) ==
false)
429 return dof_indices[dof_offsets[obj_index] + local_index];
431 return dof_indices[dof_offsets[obj_index]] + local_index;
437 DoFLevel::set_dof_index(
const unsigned int obj_index,
438 const unsigned int fe_index,
439 const unsigned int local_index,
448 Assert(dof_offsets[obj_index] != static_cast<offset_type>(-1),
449 ExcMessage(
"You are trying to access degree of freedom " 450 "information for an object on which no such " 451 "information is available"));
453 is_compressed_entry(active_fe_indices[obj_index]) ==
false,
455 "This function can no longer be called after compressing the dof_indices array"));
456 Assert(fe_index == active_fe_indices[obj_index],
457 ExcMessage(
"FE index does not match that of the present cell"));
458 dof_indices[dof_offsets[obj_index] + local_index] =
global_index;
464 DoFLevel::active_fe_index(
const unsigned int obj_index)
const 468 if (is_compressed_entry(active_fe_indices[obj_index]) ==
false)
469 return active_fe_indices[obj_index];
471 return get_toggled_compression_state(active_fe_indices[obj_index]);
477 DoFLevel::fe_index_is_active(
const unsigned int obj_index,
478 const unsigned int fe_index)
const 480 return (fe_index == active_fe_index(obj_index));
486 DoFLevel::set_active_fe_index(
const unsigned int obj_index,
487 const unsigned int fe_index)
496 Assert(is_compressed_entry(fe_index) ==
false,
498 "You are using an active_fe_index that is larger than an " 499 "internal limitation for these objects. Try to work with " 500 "hp::FECollection objects that have a more modest size."));
501 Assert(fe_index != invalid_active_fe_index,
503 "You are using an active_fe_index that is reserved for " 504 "internal purposes for these objects. Try to work with " 505 "hp::FECollection objects that have a more modest size."));
507 active_fe_indices[obj_index] = fe_index;
513 DoFLevel::future_fe_index(
const unsigned int obj_index)
const 517 if (future_fe_index_set(obj_index))
518 return future_fe_indices[obj_index];
520 return active_fe_index(obj_index);
526 DoFLevel::set_future_fe_index(
const unsigned int obj_index,
527 const unsigned int fe_index)
536 Assert(is_compressed_entry(fe_index) ==
false,
538 "You are using a future_fe_index that is larger than an " 539 "internal limitation for these objects. Try to work with " 540 "hp::FECollection objects that have a more modest size."));
541 Assert(fe_index != invalid_active_fe_index,
543 "You are using a future_fe_index that is reserved for " 544 "internal purposes for these objects. Try to work with " 545 "hp::FECollection objects that have a more modest size."));
547 future_fe_indices[obj_index] = fe_index;
553 DoFLevel::future_fe_index_set(
const unsigned int obj_index)
const 557 return (future_fe_indices[obj_index] != invalid_active_fe_index);
563 DoFLevel::clear_future_fe_index(
const unsigned int obj_index)
567 future_fe_indices[obj_index] = invalid_active_fe_index;
573 DoFLevel::get_cell_cache_start(
const unsigned int obj_index,
574 const unsigned int dofs_per_cell)
const 578 (obj_index < cell_cache_offsets.size()) &&
579 (cell_cache_offsets[obj_index] + dofs_per_cell <=
580 cell_dof_indices_cache.size()),
582 "You are trying to access an element of the cache that stores " 583 "the indices of all degrees of freedom that live on one cell. " 584 "However, this element does not exist. Did you forget to call " 585 "DoFHandler::distribute_dofs(), or did you forget to call it " 586 "again after changing the active_fe_index of one of the cells?"));
588 return &cell_dof_indices_cache[cell_cache_offsets[obj_index]];
593 template <
class Archive>
595 DoFLevel::serialize(Archive &ar,
const unsigned int)
597 ar & this->active_fe_indices;
598 ar & this->cell_cache_offsets;
599 ar & this->cell_dof_indices_cache;
600 ar & this->dof_indices;
601 ar & this->dof_offsets;
602 ar & this->future_fe_indices;
#define AssertIndexRange(index, range)
std::vector< offset_type > cell_cache_offsets
unsigned short int active_fe_index_type
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
#define DEAL_II_NAMESPACE_CLOSE
std::vector< offset_type > dof_offsets
std::vector< types::global_dof_index > dof_indices
signed short int signed_active_fe_index_type
#define DEAL_II_NAMESPACE_OPEN
std::vector< types::global_dof_index > cell_dof_indices_cache
static const active_fe_index_type invalid_active_fe_index
std::vector< active_fe_index_type > future_fe_indices
TrilinosWrappers::types::int_type global_index(const Epetra_BlockMap &map, const ::types::global_dof_index i)
std::vector< active_fe_index_type > active_fe_indices
std::enable_if< std::is_fundamental< T >::value, std::size_t >::type memory_consumption(const T &t)