deal.II version GIT relicensing-2165-gc91f007519 2024-11-20 01:40:00+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
dof_info.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2009 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
16#ifndef dealii_mesh_worker_dof_info_h
17#define dealii_mesh_worker_dof_info_h
18
19#include <deal.II/base/config.h>
20
22
24
26
29
30#include <memory>
31
33
34namespace MeshWorker
35{
36 // Forward declaration
37#ifndef DOXYGEN
38 template <int dim, class DOFINFO>
39 class DoFInfoBox;
40#endif
41
73 template <int dim, int spacedim = dim, typename number = double>
74 class DoFInfo : public LocalResults<number>
75 {
76 public:
79
82
89 unsigned int face_number;
90
97 unsigned int sub_number;
98
103 std::vector<types::global_dof_index> indices;
104
109 std::vector<std::vector<types::global_dof_index>> indices_by_block;
110
115
120 DoFInfo(const DoFHandler<dim, spacedim> &dof_handler);
121
125 template <class DHCellIterator>
126 void
127 reinit(const DHCellIterator &c);
128
132 template <class DHCellIterator, class DHFaceIterator>
133 void
134 reinit(const DHCellIterator &c,
135 const DHFaceIterator &f,
136 const unsigned int face_no);
137
141 template <class DHCellIterator, class DHFaceIterator>
142 void
143 reinit(const DHCellIterator &c,
144 const DHFaceIterator &f,
145 const unsigned int face_no,
146 const unsigned int subface_no);
147
152 template <class DHFaceIterator>
153 void
154 set_face(const DHFaceIterator &f, const unsigned int face_no);
155
160 template <class DHFaceIterator>
161 void
162 set_subface(const DHFaceIterator &f,
163 const unsigned int face_no,
164 const unsigned int subface_no);
165
166 const BlockIndices &
167 local_indices() const;
168
169
172
177
178 private:
184 DoFInfo();
185
187 void
189
191 template <class DHCellIterator>
192 void
193 get_indices(const DHCellIterator &c);
194
196 std::vector<types::global_dof_index> indices_org;
197
204
205 friend class DoFInfoBox<dim, DoFInfo<dim, spacedim, number>>;
206 };
207
208
219 template <int dim, class DOFINFO>
221 {
222 public:
226 DoFInfoBox(const DOFINFO &seed);
227
233
237 DoFInfoBox &
239
243 void
244 reset();
245
251 template <class ASSEMBLER>
252 void
253 assemble(ASSEMBLER &ass) const;
254
255
259 DOFINFO cell;
268
274
280
285 };
286
287 //----------------------------------------------------------------------//
288
289 template <int dim, int spacedim, typename number>
291 : face_number(numbers::invalid_unsigned_int)
292 , sub_number(numbers::invalid_unsigned_int)
293 , level_cell(false)
294 {}
295
296
297
298 template <int dim, int spacedim, typename number>
300 const DoFHandler<dim, spacedim> &dof_handler)
301 : face_number(numbers::invalid_unsigned_int)
302 , sub_number(numbers::invalid_unsigned_int)
303 , level_cell(false)
304 {
305 std::vector<types::global_dof_index> aux(1);
306 aux[0] = dof_handler.get_fe().n_dofs_per_cell();
308 }
309
310
311 template <int dim, int spacedim, typename number>
312 template <class DHCellIterator>
313 inline void
315 {
316 indices.resize(c->get_fe().n_dofs_per_cell());
317 if (block_info == nullptr || block_info->local().size() == 0)
318 c->get_active_or_mg_dof_indices(indices);
319 else
320 {
321 indices_org.resize(c->get_fe().n_dofs_per_cell());
322 c->get_active_or_mg_dof_indices(indices_org);
323 set_block_indices();
324 }
325 }
326
327
328 template <int dim, int spacedim, typename number>
329 template <class DHCellIterator>
330 inline void
332 {
333 get_indices(c);
334 level_cell = c->is_level_cell();
335
337 face_number = numbers::invalid_unsigned_int;
339 if (block_info)
340 LocalResults<number>::reinit(block_info->local());
341 else
342 LocalResults<number>::reinit(aux_local_indices);
343 }
344
345
346 template <int dim, int spacedim, typename number>
347 template <class DHFaceIterator>
348 inline void
350 const unsigned int face_no)
351 {
352 face = static_cast<typename Triangulation<dim, spacedim>::face_iterator>(f);
353 face_number = face_no;
355 }
356
357
358 template <int dim, int spacedim, typename number>
359 template <class DHCellIterator, class DHFaceIterator>
360 inline void
362 const DHFaceIterator &f,
363 const unsigned int face_no)
364 {
365 if ((cell.state() != IteratorState::valid) ||
367 get_indices(c);
368 level_cell = c->is_level_cell();
369
371 set_face(f, face_no);
372
373 if (block_info)
374 LocalResults<number>::reinit(block_info->local());
375 else
376 LocalResults<number>::reinit(aux_local_indices);
377 }
378
379
380 template <int dim, int spacedim, typename number>
381 template <class DHFaceIterator>
382 inline void
384 const unsigned int face_no,
385 const unsigned int subface_no)
386 {
387 face = static_cast<typename Triangulation<dim, spacedim>::face_iterator>(f);
388 face_number = face_no;
389 sub_number = subface_no;
390 }
391
392
393 template <int dim, int spacedim, typename number>
394 template <class DHCellIterator, class DHFaceIterator>
395 inline void
397 const DHFaceIterator &f,
398 const unsigned int face_no,
399 const unsigned int subface_no)
400 {
401 if (cell.state() != IteratorState::valid ||
402 cell !=
403 static_cast<typename Triangulation<dim, spacedim>::cell_iterator>(c))
404 get_indices(c);
405 level_cell = c->is_level_cell();
406
407 cell = static_cast<typename Triangulation<dim, spacedim>::cell_iterator>(c);
408 set_subface(f, face_no, subface_no);
409
410 if (block_info)
411 LocalResults<number>::reinit(block_info->local());
412 else
413 LocalResults<number>::reinit(aux_local_indices);
414 }
415
416
417 template <int dim, int spacedim, typename number>
418 inline const BlockIndices &
420 {
421 if (block_info)
422 return block_info->local();
423 return aux_local_indices;
424 }
425
426 //----------------------------------------------------------------------//
427
428 template <int dim, class DOFINFO>
429 inline DoFInfoBox<dim, DOFINFO>::DoFInfoBox(const DOFINFO &seed)
430 : cell(seed)
431 , cell_valid(true)
432 {
433 for (const unsigned int i : GeometryInfo<dim>::face_indices())
434 {
435 exterior[i] = seed;
436 interior[i] = seed;
437 interior_face_available[i] = false;
438 exterior_face_available[i] = false;
439 }
440 }
441
442
443 template <int dim, class DOFINFO>
445 const DoFInfoBox<dim, DOFINFO> &other)
446 : cell(other.cell)
447 , cell_valid(other.cell_valid)
448 {
449 for (const unsigned int i : GeometryInfo<dim>::face_indices())
450 {
451 exterior[i] = other.exterior[i];
452 interior[i] = other.interior[i];
453 interior_face_available[i] = false;
454 exterior_face_available[i] = false;
455 }
456 }
457
458
459 template <int dim, class DOFINFO>
462 {
463 cell = other.cell;
464 cell_valid = other.cell_valid;
465 for (const unsigned int i : GeometryInfo<dim>::face_indices())
466 {
467 exterior[i] = other.exterior[i];
468 interior[i] = other.interior[i];
469 interior_face_available[i] = false;
470 exterior_face_available[i] = false;
471 }
472 return *this;
473 }
474
475
476 template <int dim, class DOFINFO>
477 inline void
479 {
480 cell_valid = false;
481 for (const unsigned int i : GeometryInfo<dim>::face_indices())
482 {
483 interior_face_available[i] = false;
484 exterior_face_available[i] = false;
485 }
486 }
487
488
489 template <int dim, class DOFINFO>
490 template <class ASSEMBLER>
491 inline void
492 DoFInfoBox<dim, DOFINFO>::assemble(ASSEMBLER &assembler) const
493 {
494 if (!cell_valid)
495 return;
496
497 assembler.assemble(cell);
498 for (const unsigned int i : GeometryInfo<dim>::face_indices())
499 {
500 // Only do something if data available
501 if (interior_face_available[i])
502 {
503 // If both data
504 // available, it is an
505 // interior face
506 if (exterior_face_available[i])
507 assembler.assemble(interior[i], exterior[i]);
508 else
509 assembler.assemble(interior[i]);
510 }
511 }
512 }
513} // namespace MeshWorker
514
516
517#endif
void reinit(const unsigned int n_blocks, const size_type n_elements_per_block)
A small class collecting the different BlockIndices involved in global, multilevel and local computat...
Definition block_info.h:95
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
unsigned int n_dofs_per_cell() const
DOFINFO interior[GeometryInfo< dim >::faces_per_cell]
Definition dof_info.h:263
void assemble(ASSEMBLER &ass) const
Definition dof_info.h:492
bool exterior_face_available[GeometryInfo< dim >::faces_per_cell]
Definition dof_info.h:279
bool interior_face_available[GeometryInfo< dim >::faces_per_cell]
Definition dof_info.h:273
DOFINFO exterior[GeometryInfo< dim >::faces_per_cell]
Definition dof_info.h:267
DoFInfoBox & operator=(const DoFInfoBox< dim, DOFINFO > &)
Definition dof_info.h:461
DoFInfoBox(const DOFINFO &seed)
Definition dof_info.h:429
unsigned int sub_number
Definition dof_info.h:97
void set_face(const DHFaceIterator &f, const unsigned int face_no)
Definition dof_info.h:349
std::vector< types::global_dof_index > indices_org
Auxiliary vector.
Definition dof_info.h:196
BlockIndices aux_local_indices
Definition dof_info.h:203
ObserverPointer< const BlockInfo, DoFInfo< dim, spacedim > > block_info
The block structure of the system.
Definition dof_info.h:171
DoFInfo(const BlockInfo &block_info)
std::vector< types::global_dof_index > indices
Definition dof_info.h:103
void set_subface(const DHFaceIterator &f, const unsigned int face_no, const unsigned int subface_no)
Definition dof_info.h:383
void set_block_indices()
Set up local block indices.
const BlockIndices & local_indices() const
Definition dof_info.h:419
void get_indices(const DHCellIterator &c)
Fill index vector with active indices.
Definition dof_info.h:314
Triangulation< dim, spacedim >::face_iterator face
The current face.
Definition dof_info.h:81
std::vector< std::vector< types::global_dof_index > > indices_by_block
Definition dof_info.h:109
unsigned int face_number
Definition dof_info.h:89
void reinit(const DHCellIterator &c)
Definition dof_info.h:331
Triangulation< dim, spacedim >::cell_iterator cell
The current cell.
Definition dof_info.h:78
void reinit(const BlockIndices &local_sizes)
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition tria.h:1556
@ valid
Iterator points to a valid object.
static const unsigned int invalid_unsigned_int
Definition types.h:220
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices()