Reference documentation for deal.II version GIT relicensing-245-g36f19064f7 2024-03-29 07:20:02+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
component_mask.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2012 - 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#ifndef dealii_fe_component_mask_h
16#define dealii_fe_component_mask_h
17
18#include <deal.II/base/config.h>
19
22
23#include <algorithm>
24#include <iosfwd>
25#include <vector>
26
28
29
30
81{
82public:
89 ComponentMask() = default;
90
94 template <typename = void>
95 DEAL_II_DEPRECATED_EARLY_WITH_COMMENT(
96 "Implicit conversions from std::vector<bool> to ComponentMask are deprecated!")
97 ComponentMask(const std::vector<bool> &block_mask);
98
108 explicit ComponentMask(const std::vector<bool> &component_mask);
109
118 ComponentMask(const unsigned int n_components, const bool initializer);
119
123 void
124 set(const unsigned int index, const bool value);
125
134 unsigned int
135 size() const;
136
150 bool
151 operator[](const unsigned int component_index) const;
152
162 bool
163 represents_n_components(const unsigned int n) const;
164
178 unsigned int
179 n_selected_components(const unsigned int overall_number_of_components =
180 numbers::invalid_unsigned_int) const;
181
188 unsigned int
189 first_selected_component(const unsigned int overall_number_of_components =
190 numbers::invalid_unsigned_int) const;
191
197 bool
199
205 operator|(const ComponentMask &mask) const;
206
212 operator&(const ComponentMask &mask) const;
213
217 bool
218 operator==(const ComponentMask &mask) const;
219
223 bool
224 operator!=(const ComponentMask &mask) const;
225
230 std::size_t
231 memory_consumption() const;
232
237 "The number of selected components in a mask "
238 "must be greater than zero.");
239
240private:
244 std::vector<bool> component_mask;
245
246 // make the output operator a friend so it can access
247 // the component_mask array
248 friend std::ostream &
249 operator<<(std::ostream &out, const ComponentMask &mask);
250};
251
252
263std::ostream &
264operator<<(std::ostream &out, const ComponentMask &mask);
265
266#ifndef DOXYGEN
267// -------------------- inline functions ---------------------
268
269template <typename>
270inline ComponentMask::ComponentMask(const std::vector<bool> &component_mask)
272{}
273
274
275inline ComponentMask::ComponentMask(const std::vector<bool> &component_mask)
276 : component_mask(component_mask)
277{}
278
279
280inline ComponentMask::ComponentMask(const unsigned int n_components,
281 const bool initializer)
282 : component_mask(n_components, initializer)
283{}
284
285
286inline unsigned int
288{
289 return component_mask.size();
290}
291
292
293inline void
294ComponentMask::set(const unsigned int index, const bool value)
295{
296 AssertIndexRange(index, component_mask.size());
298}
299
300
301inline bool
302ComponentMask::operator[](const unsigned int component_index) const
303{
304 // if the mask represents the all-component mask
305 // then always return true
306 if (component_mask.empty())
307 return true;
308 else
309 {
310 // otherwise check the validity of the index and
311 // return whatever is appropriate
312 AssertIndexRange(component_index, component_mask.size());
313 return component_mask[component_index];
314 }
315}
316
317
318inline bool
319ComponentMask::represents_n_components(const unsigned int n) const
320{
321 return ((component_mask.empty()) || (component_mask.size() == n));
322}
323
324
325inline unsigned int
326ComponentMask::n_selected_components(const unsigned int n) const
327{
328 if ((n != numbers::invalid_unsigned_int) && (size() > 0))
329 AssertDimension(n, size());
330
331 const unsigned int real_n = (n != numbers::invalid_unsigned_int ? n : size());
332 if (component_mask.empty())
333 return real_n;
334 else
335 {
336 AssertDimension(real_n, component_mask.size());
337 return std::count_if(component_mask.begin(),
338 component_mask.end(),
339 [](const bool selected) { return selected; });
340 }
341}
342
343
344inline unsigned int
345ComponentMask::first_selected_component(const unsigned int n) const
346{
347 if ((n != numbers::invalid_unsigned_int) && (size() > 0))
348 AssertDimension(n, size());
349
350 if (component_mask.empty())
351 return 0;
352 else
353 {
354 for (unsigned int c = 0; c < component_mask.size(); ++c)
355 if (component_mask[c] == true)
356 return c;
357
358 Assert(false, ExcMessage("No component is selected at all!"));
360 }
361}
362
363
364
365inline bool
367{
368 return (component_mask.empty());
369}
370
371
372
373inline ComponentMask
375{
376 // if one of the two masks denotes the all-component mask,
377 // then return the other one
378 if (component_mask.empty())
379 return mask;
380 else if (mask.component_mask.empty())
381 return *this;
382 else
383 {
384 // if both masks have individual entries set, form
385 // the combination of the two
386 AssertDimension(component_mask.size(), mask.component_mask.size());
387 std::vector<bool> new_mask(component_mask.size());
388 for (unsigned int i = 0; i < component_mask.size(); ++i)
389 new_mask[i] = (component_mask[i] || mask.component_mask[i]);
390
391 return ComponentMask(new_mask);
392 }
393}
394
395
396inline ComponentMask
398{
399 // if one of the two masks denotes the all-component mask,
400 // then return the other one
401 if (component_mask.empty())
402 return mask;
403 else if (mask.component_mask.empty())
404 return *this;
405 else
406 {
407 // if both masks have individual entries set, form
408 // the combination of the two
409 AssertDimension(component_mask.size(), mask.component_mask.size());
410 std::vector<bool> new_mask(component_mask.size());
411 for (unsigned int i = 0; i < component_mask.size(); ++i)
412 new_mask[i] = (component_mask[i] && mask.component_mask[i]);
413
414 return ComponentMask(new_mask);
415 }
416}
417
418
419inline bool
421{
422 return component_mask == mask.component_mask;
423}
424
425
426inline bool
428{
429 return component_mask != mask.component_mask;
430}
431#endif // DOXYGEN
432
433
435
436#endif
bool operator[](const unsigned int component_index) const
bool represents_n_components(const unsigned int n) const
ComponentMask()=default
ComponentMask operator&(const ComponentMask &mask) const
void set(const unsigned int index, const bool value)
std::size_t memory_consumption() const
bool represents_the_all_selected_mask() const
bool operator!=(const ComponentMask &mask) const
std::vector< bool > component_mask
unsigned int size() const
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
ComponentMask operator|(const ComponentMask &mask) const
unsigned int first_selected_component(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
bool operator==(const ComponentMask &mask) const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:494
static ::ExceptionBase & ExcNoComponentSelected()
static ::ExceptionBase & ExcMessage(std::string arg1)
static const unsigned int invalid_unsigned_int
Definition types.h:220
STL namespace.