deal.II version 9.7.0
\(\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
fe_values.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2021 - 2025 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
17
28#include <deal.II/lac/vector.h>
29
31
33
34namespace NonMatching
35{
41
42
43
44 template <int dim>
45 template <typename Number>
47 const Quadrature<1> &quadrature,
50 const DoFHandler<dim> &dof_handler,
51 const ReadVector<Number> &level_set,
52 const AdditionalData &additional_data)
55 , q_collection_1D(quadrature)
59 dof_handler,
60 level_set,
61 additional_data)
62 {
63 // Tensor products of each quadrature in q_collection_1d. Used on the
64 // non-intersected cells.
65 hp::QCollection<dim> q_collection;
66 for (const auto &quadrature : q_collection_1D)
67 q_collection.push_back(Quadrature<dim>(quadrature));
68
69 initialize(q_collection);
70 }
71
72
73
74 template <int dim>
75 template <typename Number>
97
98
99
100 template <int dim>
101 void
103 {
106
107 Assert(fe_collection->size() > 0,
108 ExcMessage("Incoming hp::FECollection can not be empty."));
109 Assert(mapping_collection->size() == fe_collection->size() ||
110 mapping_collection->size() == 1,
111 ExcMessage("Size of hp::MappingCollection must be "
112 "the same as hp::FECollection or 1."));
113 Assert(q_collection.size() == fe_collection->size() ||
114 q_collection.size() == 1,
115 ExcMessage("Size of hp::QCollection<dim> must be the "
116 "same as hp::FECollection or 1."));
117 Assert(q_collection_1D.size() == fe_collection->size() ||
118 q_collection_1D.size() == 1,
119 ExcMessage("Size of hp::QCollection<1> must be the "
120 "same as hp::FECollection or 1."));
121
122 // For each element in fe_collection, create ::FEValues objects to use
123 // on the non-intersected cells.
126 for (unsigned int fe_index = 0; fe_index < fe_collection->size();
127 ++fe_index)
128 {
129 const unsigned int mapping_index =
130 mapping_collection->size() > 1 ? fe_index : 0;
131 const unsigned int q_index = q_collection.size() > 1 ? fe_index : 0;
132
133 fe_values_inside_full_quadrature[fe_index].emplace(
134 (*mapping_collection)[mapping_index],
135 (*fe_collection)[fe_index],
136 q_collection[q_index],
137 region_update_flags.inside);
138 fe_values_outside_full_quadrature[fe_index].emplace(
139 (*mapping_collection)[mapping_index],
140 (*fe_collection)[fe_index],
141 q_collection[q_index],
142 region_update_flags.outside);
143 }
144 }
145
146
147
148 template <int dim>
149 template <bool level_dof_access>
150 void
153 const unsigned int q_index,
154 const unsigned int mapping_index)
155 {
156 this->reinit_internal(cell,
157 q_index,
158 mapping_index,
159 cell->active_fe_index());
160 }
161
162
163
164 template <int dim>
165 void
167 const unsigned int q_index,
168 const unsigned int mapping_index,
169 const unsigned int fe_index)
170 {
171 this->reinit_internal(cell, q_index, mapping_index, fe_index);
172 }
173
174
175
176 template <int dim>
177 template <typename CellIteratorType>
178 void
179 FEValues<dim>::reinit_internal(const CellIteratorType &cell,
180 const unsigned int q_index_in,
181 const unsigned int mapping_index_in,
182 const unsigned int fe_index_in)
183 {
184 current_cell_location = mesh_classifier->location_to_level_set(cell);
185
186 if (fe_index_in == numbers::invalid_unsigned_int)
187 this->active_fe_index = 0;
188 else
189 this->active_fe_index = fe_index_in;
190
191 unsigned int mapping_index = mapping_index_in;
192 unsigned int q_index = q_index_in;
193 unsigned int q_index_1D = q_index_in;
194
195 if (mapping_index == numbers::invalid_unsigned_int)
196 {
197 if (mapping_collection->size() > 1)
198 mapping_index = active_fe_index;
199 else
200 mapping_index = 0;
201 }
202
203 if (q_index == numbers::invalid_unsigned_int)
204 {
206 q_index = active_fe_index;
207 else
208 q_index = 0;
209 }
210
211 if (q_index_1D == numbers::invalid_unsigned_int)
212 {
213 if (q_collection_1D.size() > 1)
214 q_index_1D = active_fe_index;
215 else
216 q_index_1D = 0;
217 }
218
219 // These objects were created with a quadrature based on the previous cell
220 // and are thus no longer valid.
221 fe_values_inside.reset();
222 fe_values_surface.reset();
223 fe_values_outside.reset();
224
225 switch (current_cell_location)
226 {
228 {
229 Assert((active_fe_index == mapping_index) ||
230 ((mapping_collection->size() == 1) &&
231 (mapping_index == 0)),
234
235 fe_values_inside_full_quadrature.at(q_index)->reinit(cell);
236 break;
237 }
239 {
240 Assert((active_fe_index == mapping_index) ||
241 ((mapping_collection->size() == 1) &&
242 (mapping_index == 0)),
245
246 fe_values_outside_full_quadrature.at(q_index)->reinit(cell);
247 break;
248 }
250 {
251 quadrature_generator.set_1D_quadrature(q_index_1D);
252 quadrature_generator.generate(cell);
253
254 const Quadrature<dim> &inside_quadrature =
255 quadrature_generator.get_inside_quadrature();
256 const Quadrature<dim> &outside_quadrature =
257 quadrature_generator.get_outside_quadrature();
258 const ImmersedSurfaceQuadrature<dim> &surface_quadrature =
259 quadrature_generator.get_surface_quadrature();
260
261 // Even if a cell is formally intersected the number of created
262 // quadrature points can be 0. Avoid creating an FEValues object
263 // if that is the case.
264 if (inside_quadrature.size() > 0)
265 {
266 fe_values_inside.emplace((*mapping_collection)[mapping_index],
268 inside_quadrature,
269 region_update_flags.inside);
270
271 fe_values_inside->reinit(cell);
272 }
273
274 if (outside_quadrature.size() > 0)
275 {
276 fe_values_outside.emplace((*mapping_collection)[mapping_index],
278 outside_quadrature,
279 region_update_flags.outside);
280
281 fe_values_outside->reinit(cell);
282 }
283
284 if (surface_quadrature.size() > 0)
285 {
286 fe_values_surface.emplace((*mapping_collection)[mapping_index],
288 surface_quadrature,
289 region_update_flags.surface);
290 fe_values_surface->reinit(cell);
291 }
292
293 break;
294 }
295 default:
296 {
298 break;
299 }
300 }
301 }
302
303
304
305 template <int dim>
306 const std::optional<::FEValues<dim>> &
314
315
316
317 template <int dim>
318 const std::optional<::FEValues<dim>> &
326
327
328
329 template <int dim>
330 const std::optional<FEImmersedSurfaceValues<dim>> &
335
336
337
338 template <int dim>
339 template <typename Number>
342 const Quadrature<1> &quadrature,
345 const DoFHandler<dim> &dof_handler,
346 const ReadVector<Number> &level_set,
347 const AdditionalData &additional_data)
350 , q_collection_1D(quadrature)
354 dof_handler,
355 level_set,
356 additional_data)
357 {
358 // Tensor products of each quadrature in q_collection_1d. Used on the
359 // non-intersected cells.
360 hp::QCollection<dim - 1> q_collection;
361 for (const auto &quadrature : q_collection_1D)
362 q_collection.push_back(quadrature);
363
364 initialize(q_collection);
365 }
366
367
368
369 template <int dim>
370 template <typename Number>
393
394
395
396 template <int dim>
397 void
399 const hp::QCollection<dim - 1> &q_collection)
400 {
402
403 Assert(fe_collection->size() > 0,
404 ExcMessage("Incoming hp::FECollection can not be empty."));
405 Assert(
406 mapping_collection->size() == fe_collection->size() ||
407 mapping_collection->size() == 1,
409 "Size of hp::MappingCollection must be the same as hp::FECollection or 1."));
410 Assert(
411 q_collection.size() == fe_collection->size() || q_collection.size() == 1,
413 "Size of hp::QCollection<dim> must be the same as hp::FECollection or 1."));
414 Assert(
415 q_collection_1D.size() == fe_collection->size() ||
416 q_collection_1D.size() == 1,
418 "Size of hp::QCollection<1> must be the same as hp::FECollection or 1."));
419
422 q_collection,
423 region_update_flags.inside);
426 q_collection,
427 region_update_flags.outside);
428 }
429
430
431
432 template <int dim>
433 template <typename CellAccessorType>
434 void
437 const unsigned int face_no,
438 const unsigned int q_index_in,
439 const unsigned int active_fe_index_in,
440 const std::function<void(::FEInterfaceValues<dim> &,
441 const unsigned int)> &call_reinit)
442 {
444 mesh_classifier->location_to_level_set(cell, face_no);
445
446 // These objects were created with a quadrature based on the previous cell
447 // and are thus no longer valid.
448 fe_values_inside.reset();
449 fe_values_outside.reset();
450
451 switch (current_face_location)
452 {
454 {
455 call_reinit(*fe_values_inside_full_quadrature, q_index_in);
456 break;
457 }
459 {
460 call_reinit(*fe_values_outside_full_quadrature, q_index_in);
461 break;
462 }
464 {
465 unsigned int q_index = q_index_in;
466
467 if (q_index == numbers::invalid_unsigned_int)
468 {
469 unsigned int active_fe_index = active_fe_index_in;
470
471 if (active_fe_index == numbers::invalid_unsigned_int)
472 {
473 if constexpr (std::is_same_v<
475 CellAccessorType> ||
476 std::is_same_v<
478 CellAccessorType>)
479 active_fe_index = cell->active_fe_index();
480 else
481 active_fe_index = 0;
482 }
483
484 if (q_collection_1D.size() > 1)
485 q_index = active_fe_index;
486 else
487 q_index = 0;
488 }
489
490 AssertIndexRange(q_index, q_collection_1D.size());
491
492 face_quadrature_generator.set_1D_quadrature(q_index);
493 face_quadrature_generator.generate(cell, face_no);
494
495 const Quadrature<dim - 1> &inside_quadrature =
496 face_quadrature_generator.get_inside_quadrature();
497 const Quadrature<dim - 1> &outside_quadrature =
498 face_quadrature_generator.get_outside_quadrature();
499
500 // Even if a cell is formally intersected the number of created
501 // quadrature points can be 0. Avoid creating an FEInterfaceValues
502 // object if that is the case.
503 if (inside_quadrature.size() > 0)
504 {
508 inside_quadrature),
509 region_update_flags.inside);
510
511 call_reinit(*fe_values_inside, /*q_index=*/0);
512 }
513
514 if (outside_quadrature.size() > 0)
515 {
519 outside_quadrature),
520 region_update_flags.outside);
521
522 call_reinit(*fe_values_outside, /*q_index=*/0);
523 }
524 break;
525 }
526 default:
527 {
529 break;
530 }
531 }
532 }
533
534
535
536 template <int dim>
537 const std::optional<::FEInterfaceValues<dim>> &
545
546
547
548 template <int dim>
549 const std::optional<::FEInterfaceValues<dim>> &
557
558
559#include "non_matching/fe_values.inst"
560
561} // namespace NonMatching
const std::optional<::FEInterfaceValues< dim > > & get_outside_fe_values() const
Definition fe_values.cc:550
void initialize(const hp::QCollection< dim - 1 > &q_collection)
Definition fe_values.cc:398
const ObserverPointer< const hp::MappingCollection< dim > > mapping_collection
Definition fe_values.h:653
const hp::QCollection< 1 > q_collection_1D
Definition fe_values.h:664
const RegionUpdateFlags region_update_flags
Definition fe_values.h:674
void do_reinit(const TriaIterator< CellAccessorType > &cell, const unsigned int face_no, const unsigned int q_index, const unsigned int active_fe_index, const std::function< void(::FEInterfaceValues< dim > &, const unsigned int)> &call_reinit)
Definition fe_values.cc:435
LocationToLevelSet current_face_location
Definition fe_values.h:669
const ObserverPointer< const MeshClassifier< dim > > mesh_classifier
Definition fe_values.h:679
std::optional<::FEInterfaceValues< dim > > fe_values_inside_full_quadrature
Definition fe_values.h:686
const std::optional<::FEInterfaceValues< dim > > & get_inside_fe_values() const
Definition fe_values.cc:538
const ObserverPointer< const hp::FECollection< dim > > fe_collection
Definition fe_values.h:658
std::optional<::FEInterfaceValues< dim > > fe_values_outside_full_quadrature
Definition fe_values.h:693
std::optional<::FEInterfaceValues< dim > > fe_values_inside
Definition fe_values.h:702
DiscreteFaceQuadratureGenerator< dim > face_quadrature_generator
Definition fe_values.h:716
FEInterfaceValues(const hp::FECollection< dim > &fe_collection, const Quadrature< 1 > &quadrature, const RegionUpdateFlags region_update_flags, const MeshClassifier< dim > &mesh_classifier, const DoFHandler< dim > &dof_handler, const ReadVector< Number > &level_set, const AdditionalData &additional_data=AdditionalData())
Definition fe_values.cc:340
typename FaceQuadratureGenerator< dim >::AdditionalData AdditionalData
Definition fe_values.h:491
std::optional<::FEInterfaceValues< dim > > fe_values_outside
Definition fe_values.h:711
std::optional< NonMatching::FEImmersedSurfaceValues< dim > > fe_values_surface
Definition fe_values.h:415
std::deque< std::optional<::FEValues< dim > > > fe_values_inside_full_quadrature
Definition fe_values.h:372
std::optional<::FEValues< dim > > fe_values_outside
Definition fe_values.h:406
LocationToLevelSet current_cell_location
Definition fe_values.h:341
const ObserverPointer< const hp::FECollection< dim > > fe_collection
Definition fe_values.h:330
const std::optional<::FEValues< dim > > & get_outside_fe_values() const
Definition fe_values.cc:319
const ObserverPointer< const MeshClassifier< dim > > mesh_classifier
Definition fe_values.h:356
void initialize(const hp::QCollection< dim > &q_collection)
Definition fe_values.cc:102
const std::optional< FEImmersedSurfaceValues< dim > > & get_surface_fe_values() const
Definition fe_values.cc:331
void reinit_internal(const CellIteratorType &cell, const unsigned int q_index, const unsigned int mapping_index, const unsigned int fe_index)
Definition fe_values.cc:179
const hp::QCollection< 1 > q_collection_1D
Definition fe_values.h:336
const ObserverPointer< const hp::MappingCollection< dim > > mapping_collection
Definition fe_values.h:325
const RegionUpdateFlags region_update_flags
Definition fe_values.h:351
void reinit(const TriaIterator< DoFCellAccessor< dim, dim, level_dof_access > > &cell, const unsigned int q_index=numbers::invalid_unsigned_int, const unsigned int mapping_index=numbers::invalid_unsigned_int)
Definition fe_values.cc:151
std::optional<::FEValues< dim > > fe_values_inside
Definition fe_values.h:397
FEValues(const hp::FECollection< dim > &fe_collection, const Quadrature< 1 > &quadrature, const RegionUpdateFlags region_update_flags, const MeshClassifier< dim > &mesh_classifier, const DoFHandler< dim > &dof_handler, const ReadVector< Number > &level_set, const AdditionalData &additional_data=AdditionalData())
Definition fe_values.cc:46
typename QuadratureGenerator< dim >::AdditionalData AdditionalData
Definition fe_values.h:146
DiscreteQuadratureGenerator< dim > quadrature_generator
Definition fe_values.h:420
std::deque< std::optional<::FEValues< dim > > > fe_values_outside_full_quadrature
Definition fe_values.h:388
unsigned int active_fe_index
Definition fe_values.h:346
const std::optional<::FEValues< dim > > & get_inside_fe_values() const
Definition fe_values.cc:307
unsigned int size() const
unsigned int size() const
Definition collection.h:316
void push_back(const Quadrature< dim_in > &new_quadrature)
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_ASSERT_UNREACHABLE()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcMessage(std::string arg1)
@ update_default
No update.
Definition hp.h:117
constexpr unsigned int invalid_unsigned_int
Definition types.h:238