Fluid structure interaction suite
components.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2022 by Luca Heltai
4 //
5 // This file is part of the FSI-suite platform, based on the deal.II library.
6 //
7 // The FSI-suite platform is free software; you can use it, redistribute it,
8 // and/or modify it under the terms of the GNU Lesser General Public License as
9 // published by the Free Software Foundation; either version 3.0 of the License,
10 // or (at your option) any later version. The full text of the license can be
11 // found in the file LICENSE at the top level of the FSI-suite platform
12 // distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef parsed_tools_components_h
17 #define parsed_tools_components_h
18 
19 #include <deal.II/base/config.h>
20 
23 #include <deal.II/base/utilities.h>
24 
26 #include <deal.II/fe/fe_dgq.h>
27 #include <deal.II/fe/fe_q.h>
29 #include <deal.II/fe/fe_system.h>
30 
31 #include <deal.II/grid/tria.h>
32 
33 #include "parsed_tools/enum.h"
34 
35 namespace ParsedTools
36 {
37  namespace Components
38  {
43  enum class Type
44  {
45  all = 0,
46  scalar = 1,
47  vector = 2,
48  normal = 3,
49  tangential = 4,
50  };
51 
66  template <int dim, int spacedim>
67  std::unique_ptr<dealii::FiniteElement<dim, spacedim>>
69  const dealii::Triangulation<dim, spacedim> &tria,
70  const unsigned int degree = 1,
71  const int continuity = 0);
72 
82  template <int dim, int spacedim>
83  dealii::Quadrature<dim>
84  get_cell_quadrature(const dealii::Triangulation<dim, spacedim> &tria,
85  const unsigned int degree);
86 
96  template <int dim, int spacedim>
97  dealii::Quadrature<dim - 1>
98  get_face_quadrature(const dealii::Triangulation<dim, spacedim> &tria,
99  const unsigned int degree);
100 
105  template <typename Container>
106  std::string
107  join(const Container &strings, const std::string &separator);
108 
113  unsigned int
114  n_components(const std::string &component_names);
115 
120  unsigned int
121  n_blocks(const std::string &component_names);
122 
123 
127  std::string
128  blocks_to_names(const std::vector<std::string> &components,
129  const std::vector<unsigned int> &multiplicities);
130 
134  std::pair<std::vector<std::string>, std::vector<unsigned int>>
135  names_to_blocks(const std::string &component_names);
136 
163  std::vector<unsigned int>
164  block_indices(const std::string &component_names,
165  const std::string &selected_components);
166 
197  std::vector<unsigned int>
198  component_indices(const std::string &component_names,
199  const std::string &selected_components);
200 
233  std::pair<unsigned int, unsigned int>
234  component_to_indices(const std::string &component_names,
235  const std::string &selected_component);
236 
261  std::string
262  component_name(const std::string &component_names,
263  const std::string &selected_component);
264 
292  Type
293  type(const std::string &component_name,
294  const std::string &selected_component);
295 
300  dealii::ComponentMask
301  mask(const std::string &component_names,
302  const std::string &selected_component);
303 
304 
305 
306 #ifndef DOXYGEN
307  // Template implementation
308  template <typename Container>
309  std::string
310  join(const Container &strings, const std::string &separator)
311  {
312  std::string result;
313  std::string sep = "";
314  for (const auto &s : strings)
315  {
316  result += sep + s;
317  sep = separator;
318  }
319  return result;
320  }
321 
322 
323 
324  template <int dim, int spacedim>
325  std::unique_ptr<dealii::FiniteElement<dim, spacedim>>
327  const dealii::Triangulation<dim, spacedim> &tria,
328  const unsigned int degree,
329  const int continuity)
330  {
331  const auto ref_cells = tria.get_reference_cells();
332  AssertThrow(
333  ref_cells.size() == 1,
335  "This function does nots support mixed simplx/hex grid types."));
336  AssertThrow(continuity == -1 || continuity == 0,
338  "only -1 and 0 are supported for continuity"));
339  std::unique_ptr<dealii::FiniteElement<dim, spacedim>> result;
340  if (ref_cells[0].is_simplex())
341  {
342  if (continuity == 0)
343  result.reset(new dealii::FE_SimplexP<dim, spacedim>(degree));
344  else
345  result.reset(new dealii::FE_SimplexDGP<dim, spacedim>(degree));
346  }
347  else
348  {
349  Assert(ref_cells[0].is_hyper_cube(),
351  "Only simplex and hex cells are supported"));
352  if (continuity == 0)
353  result.reset(new dealii::FE_Q<dim, spacedim>(degree));
354  else
355  result.reset(new dealii::FE_DGQ<dim, spacedim>(degree));
356  }
357  return result;
358  }
359 
360 
361 
362  template <int dim, int spacedim>
363  dealii::Quadrature<dim>
364  get_cell_quadrature(const dealii::Triangulation<dim, spacedim> &tria,
365  const unsigned int degree)
366  {
367  const auto ref_cells = tria.get_reference_cells();
368 
369  AssertThrow(
370  ref_cells.size() == 1,
372  "This function does nots support mixed simplx/hex grid types."));
373  return ref_cells[0].template get_gauss_type_quadrature<dim>(degree);
374  }
375 
376 
377 
378  template <int dim, int spacedim>
379  dealii::Quadrature<dim - 1>
380  get_face_quadrature(const dealii::Triangulation<dim, spacedim> &tria,
381  const unsigned int degree)
382  {
383  if constexpr (dim == 1)
384  {
385  return dealii::QGauss<dim - 1>(degree);
386  }
387  else
388  {
389  const auto ref_cells = tria.get_reference_cells();
390 
391  AssertThrow(
392  ref_cells.size() == 1,
394  "This function does nots support mixed simplx/hex grid types."));
395 
396  const dealii::ReferenceCell face_type =
397  ref_cells[0].face_reference_cell(0);
398  return face_type.template get_gauss_type_quadrature<dim - 1>(degree);
399  }
400  }
401 #endif
402  } // namespace Components
403 } // namespace ParsedTools
404 #endif
#define Assert(cond, exc)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
std::pair< unsigned int, unsigned int > component_to_indices(const std::string &component_names, const std::string &selected_component)
Return the indices within the components and within the blocks corresponding to first component of th...
Definition: components.cc:170
std::unique_ptr< FiniteElement< dim, spacedim > > get_lagrangian_finite_element(const Triangulation< dim, spacedim > &tria, const unsigned int degree=1, const int continuity=0)
Get a Lagrangian FiniteElement object compatible with the given Triangulation object.
std::pair< std::vector< std::string >, std::vector< unsigned int > > names_to_blocks(const std::string &component_names)
Build block names and multiplicities from component names.
Definition: components.cc:56
Quadrature< dim - 1 > get_face_quadrature(const Triangulation< dim, spacedim > &tria, const unsigned int degree)
Return a Quadrature object that can be used on the given Triangulation faces.
std::string component_name(const std::string &component_names, const std::string &selected_component)
Return the canonical component name for the given selected component.
Definition: components.cc:72
Type
Enumerator used to identify patterns of components and their size in the a block system.
Definition: components.h:44
@ tangential
Tangential component of a vector.
@ normal
Normal component of a vector.
unsigned int n_blocks(const std::string &component_names)
Count the number of components in the given list of comma separated components.
Definition: components.cc:217
std::string blocks_to_names(const std::vector< std::string > &components, const std::vector< unsigned int > &multiplicities)
Build component names from block names and multiplicities.
Definition: components.cc:40
Type type(const std::string &component_name, const std::string &selected_component)
Return the component type for the given selected component.
Definition: components.cc:186
std::vector< unsigned int > block_indices(const std::string &component_names, const std::string &selected_components)
Return the indices within the block corresponding to the given selected components,...
Definition: components.cc:134
std::vector< unsigned int > component_indices(const std::string &component_names, const std::string &selected_components)
Return the indices within the components corresponding to first component of the given selected compo...
Definition: components.cc:156
Quadrature< dim > get_cell_quadrature(const Triangulation< dim, spacedim > &tria, const unsigned int degree)
Return a Quadrature object that can be used on the given Triangulation cells.
std::string join(const Container &strings, const std::string &separator)
Join strings in a container together using a given separator.
ComponentMask mask(const std::string &component_names, const std::string &selected_component)
Return a component mask corresponding to a given selected component, from a list of comma separated c...
Definition: components.cc:226
unsigned int n_components(const std::string &component_names)
Count the number of components in the given list of comma separated components.
Definition: components.cc:32
We collect in this namespace some wrappers around commonly used deal.II classes, derived from the Par...