16 #ifndef dealii__hp_dof_level_h 17 #define dealii__hp_dof_level_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/exceptions.h> 26 DEAL_II_NAMESPACE_OPEN
31 template <
int,
int>
class FECollection;
41 struct Implementation;
46 struct Implementation;
182 set_dof_index (
const unsigned int obj_index,
183 const unsigned int fe_index,
184 const unsigned int local_index,
199 get_dof_index (
const unsigned int obj_index,
200 const unsigned int fe_index,
201 const unsigned int local_index)
const;
207 active_fe_index (
const unsigned int obj_index)
const;
214 fe_index_is_active (
const unsigned int obj_index,
215 const unsigned int fe_index)
const;
221 set_active_fe_index (
const unsigned int obj_index,
222 const unsigned int fe_index);
236 get_cell_cache_start (
const unsigned int obj_index,
237 const unsigned int dofs_per_cell)
const;
243 std::size_t memory_consumption ()
const;
249 template <
class Archive>
250 void serialize(Archive &ar,
251 const unsigned int version);
262 template <
int dim,
int spacedim>
263 void compress_data (const ::hp::FECollection<dim,spacedim> &fe_collection);
273 template <
int dim,
int spacedim>
274 void uncompress_data (const ::hp::FECollection<dim,spacedim> &fe_collection);
280 template <
int,
int>
friend class ::hp::DoFHandler;
281 friend struct ::internal::hp::DoFHandler::Implementation;
282 friend struct ::internal::DoFCellAccessor::Implementation;
291 get_dof_index (
const unsigned int obj_index,
292 const unsigned int fe_index,
293 const unsigned int local_index)
const 296 Assert (obj_index < dof_offsets.size(),
297 ExcIndexRange (obj_index, 0, dof_offsets.size()));
303 ExcMessage (
"You are trying to access degree of freedom " 304 "information for an object on which no such " 305 "information is available"));
307 Assert (fe_index == active_fe_indices[obj_index],
308 ExcMessage (
"FE index does not match that of the present cell"));
313 return dof_indices[dof_offsets[obj_index]+local_index];
315 return dof_indices[dof_offsets[obj_index]]+local_index;
323 set_dof_index (
const unsigned int obj_index,
324 const unsigned int fe_index,
325 const unsigned int local_index,
329 Assert (obj_index < dof_offsets.size(),
330 ExcIndexRange (obj_index, 0, dof_offsets.size()));
336 ExcMessage (
"You are trying to access degree of freedom " 337 "information for an object on which no such " 338 "information is available"));
340 ExcMessage (
"This function can no longer be called after compressing the dof_indices array"));
341 Assert (fe_index == active_fe_indices[obj_index],
342 ExcMessage (
"FE index does not match that of the present cell"));
343 dof_indices[dof_offsets[obj_index]+local_index] = global_index;
351 active_fe_index (
const unsigned int obj_index)
const 353 Assert (obj_index < active_fe_indices.size(),
354 ExcIndexRange (obj_index, 0, active_fe_indices.size()));
357 return active_fe_indices[obj_index];
367 fe_index_is_active (
const unsigned int obj_index,
368 const unsigned int fe_index)
const 370 return (fe_index == active_fe_index(obj_index));
377 set_active_fe_index (
const unsigned int obj_index,
378 const unsigned int fe_index)
380 Assert (obj_index < active_fe_indices.size(),
381 ExcIndexRange (obj_index, 0, active_fe_indices.size()));
383 active_fe_indices[obj_index] = fe_index;
390 DoFLevel::get_cell_cache_start (
const unsigned int obj_index,
391 const unsigned int dofs_per_cell)
const 394 Assert (obj_index < cell_cache_offsets.size(),
396 Assert (cell_cache_offsets[obj_index]+dofs_per_cell
398 cell_dof_indices_cache.size(),
401 return &cell_dof_indices_cache[cell_cache_offsets[obj_index]];
404 template <
class Archive>
407 DoFLevel::serialize(Archive &ar,
410 ar &this->active_fe_indices;
411 ar &this->cell_cache_offsets;
412 ar &this->cell_dof_indices_cache;
413 ar &this->dof_indices;
414 ar &this->dof_offsets;
420 DEAL_II_NAMESPACE_CLOSE
::ExceptionBase & ExcMessage(std::string arg1)
std::vector< offset_type > cell_cache_offsets
unsigned short int active_fe_index_type
unsigned int global_dof_index
#define Assert(cond, exc)
std::vector< offset_type > dof_offsets
std::vector< types::global_dof_index > dof_indices
std::vector< types::global_dof_index > cell_dof_indices_cache
signed short int signed_active_fe_index_type
std::vector< active_fe_index_type > active_fe_indices