Loading [MathJax]/extensions/tex2jax.js
PISM, A Parallel Ice Sheet Model 2.2.2-d6b3a29ca committed by Constantine Khrulev on 2025-03-28
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Quadrature.cc
Go to the documentation of this file.
1/* Copyright (C) 2020, 2023, 2025 PISM Authors
2 *
3 * This file is part of PISM.
4 *
5 * PISM is free software; you can redistribute it and/or modify it under the
6 * terms of the GNU General Public License as published by the Free Software
7 * Foundation; either version 3 of the License, or (at your option) any later
8 * version.
9 *
10 * PISM is distributed in the hope that it will be useful, but WITHOUT ANY
11 * WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
12 * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more
13 * details.
14 *
15 * You should have received a copy of the GNU General Public License
16 * along with PISM; if not, write to the Free Software
17 * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
18 */
19
20#include "pism/util/fem/Quadrature.hh"
21#include <cmath>
22
23namespace pism {
24namespace fem {
25
26const std::vector<QuadPoint>& Quadrature::points() const {
27 return m_points;
28}
29
30const std::vector<double>& Quadrature::weights() const {
31 return m_weights;
32}
33
34//! Build quadrature points and weights for a tensor product quadrature based on a 1D quadrature
35//! rule. Uses the same 1D quadrature in both directions.
36/**
37 @param[in] n 1D quadrature size (the resulting quadrature has size n*n)
38 @param[in] points1 1D quadrature points
39 @param[in] weights1 1D quadrature weights
40 @param[out] points resulting 2D quadrature points
41 @param[out] weights resulting 2D quadrature weights
42 */
43static void tensor_product_quadrature(unsigned int n,
44 const double *points1,
45 const double *weights1,
46 std::vector<QuadPoint>& points,
47 std::vector<double>& weights) {
48 points.resize(n * n);
49 weights.resize(n * n);
50
51 unsigned int q = 0;
52 for (unsigned int j = 0; j < n; ++j) {
53 for (unsigned int i = 0; i < n; ++i) {
54 points[q] = {points1[i], points1[j], 0.0};
55 weights[q] = weights1[i] * weights1[j];
56
57 ++q;
58 }
59 }
60}
62
63 // coordinates and weights of the 2-point 1D Gaussian quadrature
64 double A = 1.0 / std::sqrt(3.0);
65
66 m_points = {{-A, 0.0, 0.0}, {A, 0.0, 0.0}};
67 m_weights = {0.5 * D, 0.5 * D};
68}
69
70//! One-point quadrature on a rectangle.
72 m_points = {{0.0, 0.0, 0.0}};
73 m_weights = {4.0};
74}
75
76//! Two-by-two Gaussian quadrature on a rectangle.
78
79 // coordinates and weights of the 2-point 1D Gaussian quadrature
80 const double
81 A = 1.0 / sqrt(3.0),
82 points2[2] = {-A, A},
83 weights2[2] = {1.0, 1.0};
84
85 tensor_product_quadrature(2, points2, weights2, m_points, m_weights);
86}
87
89 const double
90 A = 0.0,
91 B = sqrt(0.6),
92 points3[3] = {-B, A, B};
93
94 const double
95 w1 = 5.0 / 9.0,
96 w2 = 8.0 / 9.0,
97 weights3[3] = {w1, w2, w1};
98
99 tensor_product_quadrature(3, points3, weights3, m_points, m_weights);
100}
101
103 const double
104 A = sqrt(3.0 / 7.0 - (2.0 / 7.0) * sqrt(6.0 / 5.0)), // smaller magnitude
105 B = sqrt(3.0 / 7.0 + (2.0 / 7.0) * sqrt(6.0 / 5.0)), // larger magnitude
106 points4[4] = {-B, -A, A, B};
107
108 // The weights w_i for Gaussian quadrature on the reference element with these
109 // quadrature points
110 const double
111 w1 = (18.0 + sqrt(30.0)) / 36.0, // larger
112 w2 = (18.0 - sqrt(30.0)) / 36.0, // smaller
113 weights4[4] = {w2, w1, w1, w2};
114
115 tensor_product_quadrature(4, points4, weights4, m_points, m_weights);
116}
117
118
119//! @brief N*N-point uniform (*not* Gaussian) quadrature for integrating discontinuous
120//! functions.
122
123 std::vector<double> xi(N), w(N);
124 const double dxi = 2.0 / N;
125 for (unsigned int k = 0; k < N; ++k) {
126 xi[k] = -1.0 + dxi*(k + 0.5);
127 w[k] = 2.0 / N;
128 }
129
130 tensor_product_quadrature(N, xi.data(), w.data(), m_points, m_weights);
131}
132
133/*!
134 * 3-point Gaussian quadrature on the triangle (0,0)-(1,0)-(0,1)
135 */
137
138 const double
139 one_over_six = 1.0 / 6.0,
140 two_over_three = 2.0 / 3.0;
141
142 m_points = {{two_over_three, one_over_six, 0.0},
143 {one_over_six, two_over_three, 0.0},
144 {one_over_six, one_over_six, 0.0}};
145
146 m_weights = {one_over_six, one_over_six, one_over_six};
147}
148
150 double xis[8] = {-1.0, 1.0, 1.0, -1.0, -1.0, 1.0, 1.0, -1.0};
151 double etas[8] = {-1.0, -1.0, 1.0, 1.0, -1.0, -1.0, 1.0, 1.0};
152 double zetas[8] = {-1.0, -1.0, -1.0, -1.0, 1.0, 1.0, 1.0, 1.0};
153
154 m_points.resize(8);
155 m_weights.resize(8);
156
157 double C = 1.0 / sqrt(3.0);
158
159 for (int k = 0; k < 8; ++k) {
160 m_points[k] = {C * xis[k], C * etas[k], C * zetas[k]};
161 m_weights[k] = 1.0;
162 }
163}
164
166 m_points = {{0.0, 0.0, 0.0}};
167 m_weights = {8.0};
168}
169
171 const double
172 A = sqrt(3.0 / 7.0 - (2.0 / 7.0) * sqrt(6.0 / 5.0)), // smaller magnitude
173 B = sqrt(3.0 / 7.0 + (2.0 / 7.0) * sqrt(6.0 / 5.0)), // larger magnitude
174 pt[4] = {-B, -A, A, B};
175
176 // The weights w_i for Gaussian quadrature on the reference element with these
177 // quadrature points
178 const double
179 w1 = (18.0 + sqrt(30.0)) / 36.0, // larger
180 w2 = (18.0 - sqrt(30.0)) / 36.0, // smaller
181 w[4] = {w2, w1, w1, w2};
182
183 m_points.resize(64);
184 m_weights.resize(64);
185 int q = 0;
186 for (int i = 0; i < 4; ++i) {
187 for (int j = 0; j < 4; ++j) {
188 for (int k = 0; k < 4; ++k) {
189 m_points[q] = {pt[i], pt[j], pt[k]};
190 m_weights[q] = w[i] * w[j] * w[k];
191 ++q;
192 }
193 }
194 }
195}
196
197} // end of namespace fem
198} // end of namespace pism
Q1Quadrature1()
One-point quadrature on a rectangle.
Definition Quadrature.cc:71
Q1Quadrature4()
Two-by-two Gaussian quadrature on a rectangle.
Definition Quadrature.cc:77
Q1QuadratureN(unsigned int n)
N*N-point uniform (not Gaussian) quadrature for integrating discontinuous functions.
std::vector< QuadPoint > m_points
Definition Quadrature.hh:76
const std::vector< double > & weights() const
Definition Quadrature.cc:30
std::vector< double > m_weights
Definition Quadrature.hh:77
const std::vector< QuadPoint > & points() const
Definition Quadrature.cc:26
#define n
Definition exactTestM.c:37
static void tensor_product_quadrature(unsigned int n, const double *points1, const double *weights1, std::vector< QuadPoint > &points, std::vector< double > &weights)
Definition Quadrature.cc:43
static const double k
Definition exactTestP.cc:42