QUGaR 0.0.9
Loading...
Searching...
No Matches
cut_quadrature.hpp
Go to the documentation of this file.
1// --------------------------------------------------------------------------
2//
3// Copyright (C) 2025-present by Pablo Antolin
4//
5// This file is part of the QUGaR library.
6//
7// SPDX-License-Identifier: MIT
8//
9// --------------------------------------------------------------------------
10
11#ifndef QUGAR_LIBRARY_CUT_QUADRATURE_HPP
12#define QUGAR_LIBRARY_CUT_QUADRATURE_HPP
13
20
21
22#include <qugar/point.hpp>
23#include <qugar/types.hpp>
25
26#include <cstddef>
27#include <cstdint>
28#include <memory>
29#include <vector>
30
31
32namespace qugar {
33
34enum ImmersedStatus : std::uint8_t { cut, full, empty };
35
36template<int dim> struct CutCellsQuad
37{
39 std::vector<std::int64_t> cells;
41 std::vector<int> n_pts_per_cell;
42
44 std::vector<Point<dim>> points;
46 std::vector<real> weights;
47
48 void reserve(const int n_cells, const int n_tot_pts);
49};
50
51template<int dim> struct CutIsoBoundsQuad
52{
54 std::vector<std::int64_t> cells;
56 std::vector<int> local_facet_ids;
58 std::vector<int> n_pts_per_facet;
59
61 std::vector<Point<dim>> points;
63 std::vector<real> weights;
64
65 void reserve(const int n_cells, const int n_tot_pts);
66};
67
68template<int dim> struct CutUnfBoundsQuad
69{
71 std::vector<std::int64_t> cells;
73 std::vector<int> n_pts_per_cell;
74
76 std::vector<Point<dim>> points;
78 std::vector<real> weights;
80 std::vector<Point<dim>> normals;
81
82 void reserve(const int n_cells, const int n_tot_pts);
83};
84
98template<int dim>
99std::shared_ptr<const CutCellsQuad<dim>>
100 create_quadrature(const UnfittedDomain<dim> &unf_domain, const std::vector<std::int64_t> &cells, int n_pts_dir);
101
114template<int dim>
115std::shared_ptr<const CutUnfBoundsQuad<dim>> create_unfitted_bound_quadrature(const UnfittedDomain<dim> &unf_domain,
116 const std::vector<std::int64_t> &cells,
117 int n_pts_dir,
118 bool include_facet_unf_bdry,
119 bool exclude_ext_bdry);
120
138template<int dim>
140 const UnfittedDomain<dim> &unf_domain,
141 const std::vector<std::int64_t> &cells,
142 const std::vector<int> &facets,
143 int n_pts_dir);
144
159template<int dim>
161 const UnfittedDomain<dim> &unf_domain,
162 const std::vector<std::int64_t> &cells,
163 const std::vector<int> &facets,
164 int n_pts_dir);
165
166
167}// namespace qugar
168
169#endif// QUGAR_LIBRARY_CUT_QUADRATURE_HPP
Definition unfitted_domain.hpp:51
QUGaR's main namespace.
Definition affine_transf.hpp:28
ImmersedStatus
Definition cut_quadrature.hpp:34
@ cut
Definition cut_quadrature.hpp:34
@ empty
Definition cut_quadrature.hpp:34
@ full
Definition cut_quadrature.hpp:34
std::shared_ptr< const CutIsoBoundsQuad< dim - 1 > > create_facets_quadrature_exterior_integral(const UnfittedDomain< dim > &unf_domain, const std::vector< std::int64_t > &cells, const std::vector< int > &facets, int n_pts_dir)
Creates quadrature for exterior integrals.
std::shared_ptr< const CutUnfBoundsQuad< dim > > create_unfitted_bound_quadrature(const UnfittedDomain< dim > &unf_domain, const std::vector< std::int64_t > &cells, int n_pts_dir, bool include_facet_unf_bdry, bool exclude_ext_bdry)
Creates a quadrature for the unfitted boundary.
std::shared_ptr< const CutCellsQuad< dim > > create_quadrature(const UnfittedDomain< dim > &unf_domain, const std::vector< std::int64_t > &cells, int n_pts_dir)
Creates quadrature for cells.
std::shared_ptr< const CutIsoBoundsQuad< dim - 1 > > create_facets_quadrature_interior_integral(const UnfittedDomain< dim > &unf_domain, const std::vector< std::int64_t > &cells, const std::vector< int > &facets, int n_pts_dir)
Creates quadrature for interior integrals.
Definition and implementation of Point class.
Definition cut_quadrature.hpp:37
std::vector< int > n_pts_per_cell
List of number of quadrature points per cut cell.
Definition cut_quadrature.hpp:41
std::vector< std::int64_t > cells
List of cut cells.
Definition cut_quadrature.hpp:39
std::vector< real > weights
Vector of weights.
Definition cut_quadrature.hpp:46
std::vector< Point< dim > > points
Vector of points.
Definition cut_quadrature.hpp:44
void reserve(const int n_cells, const int n_tot_pts)
Definition cut_quadrature.hpp:52
std::vector< int > local_facet_ids
List of (local) facet ids associated to the vector cells.
Definition cut_quadrature.hpp:56
void reserve(const int n_cells, const int n_tot_pts)
std::vector< Point< dim > > points
Vector of points.
Definition cut_quadrature.hpp:61
std::vector< std::int64_t > cells
List of cut cells.
Definition cut_quadrature.hpp:54
std::vector< int > n_pts_per_facet
List of number of quadrature points per cut facet.
Definition cut_quadrature.hpp:58
std::vector< real > weights
Vector of weights.
Definition cut_quadrature.hpp:63
Definition cut_quadrature.hpp:69
std::vector< Point< dim > > points
Vector of points.
Definition cut_quadrature.hpp:76
std::vector< std::int64_t > cells
List of cut cells.
Definition cut_quadrature.hpp:71
std::vector< real > weights
Vector of weights.
Definition cut_quadrature.hpp:78
std::vector< Point< dim > > normals
Vector of normals.
Definition cut_quadrature.hpp:80
std::vector< int > n_pts_per_cell
List of number of quadrature points per cut cell.
Definition cut_quadrature.hpp:73
void reserve(const int n_cells, const int n_tot_pts)
Declaration of types.
Declaration of of UnfittedDomain class.