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
FEM.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#include <cassert>
20
21#include "pism/util/fem/FEM.hh"
22#include "pism/util/node_types.hh"
23
24namespace pism {
25namespace fem {
26
27namespace linear {
28
29//! Linear basis functions on the interval [-1, -1]
30Germ chi(unsigned int k, const QuadPoint &pt) {
31 // coordinates of reference element nodes
32 static const double xis[n_chi] = {-1.0, 1.0};
33
34 assert(k < linear::n_chi);
35
36 return {0.5 * (1.0 + xis[k] * pt.xi),
37 0.5 * xis[k],
38 0.0, 0.0}; // unused
39}
40
41} // end of namespace linear
42
43namespace q0 {
44/*!
45 * Piecewise-constant shape functions.
46 */
47Germ chi(unsigned int k, const QuadPoint &pt) {
48 assert(k < q0::n_chi);
49
50 Germ result;
51
52 if ((k == 0 and pt.xi <= 0.0 and pt.eta <= 0.0) or
53 (k == 1 and pt.xi > 0.0 and pt.eta <= 0.0) or
54 (k == 2 and pt.xi > 0.0 and pt.eta > 0.0) or
55 (k == 3 and pt.xi <= 0.0 and pt.eta > 0.0)) {
56 result.val = 1.0;
57 } else {
58 result.val = 0.0;
59 }
60
61 result.dx = 0.0;
62 result.dy = 0.0;
63 result.dz = 0.0;
64
65 return result;
66}
67} // end of namespace q0
68
69namespace q1 {
70
71//! Q1 basis functions on the reference element with nodes (-1,-1), (1,-1), (1,1), (-1,1).
72Germ chi(unsigned int k, const QuadPoint &pt) {
73 // coordinates of reference element nodes
74 static const double xi[n_chi] = {-1.0, 1.0, 1.0, -1.0};
75 static const double eta[n_chi] = {-1.0, -1.0, 1.0, 1.0};
76
77 assert(k < q1::n_chi);
78
79 return {0.25 * (1.0 + xi[k] * pt.xi) * (1.0 + eta[k] * pt.eta),
80 0.25 * xi[k] * (1.0 + eta[k] * pt.eta),
81 0.25 * eta[k] * (1.0 + xi[k] * pt.xi),
82 0.0}; // unused
83}
84
85} // end of namespace q1
86
87namespace p1 {
88
89//! P1 basis functions on the reference element with nodes (0,0), (1,0), (0,1).
90Germ chi(unsigned int k, const QuadPoint &pt) {
91 assert(k < q1::n_chi);
92
93 switch (k) {
94 default:
95 case 0:
96 return {1.0 - pt.xi - pt.eta, -1.0, -1.0, 0.0};
97 case 1:
98 return {pt.xi, 1.0, 0.0, 0.0};
99 case 2:
100 return {pt.eta, 0.0, 1.0, 0.0};
101 }
102}
103
104} // end of namespace p1
105
106ElementType element_type(const int node_type[q1::n_chi]) {
107
108 // number of exterior nodes in this element
109 const int n_exterior_nodes = (static_cast<int>(node_type[0] == NODE_EXTERIOR) +
110 static_cast<int>(node_type[1] == NODE_EXTERIOR) +
111 static_cast<int>(node_type[2] == NODE_EXTERIOR) +
112 static_cast<int>(node_type[3] == NODE_EXTERIOR));
113
114 // an element is a "Q1 interior" if all its nodes are interior or boundary
115 if (n_exterior_nodes == 0) {
116 return ELEMENT_Q;
117 }
118
119 if (n_exterior_nodes == 1) {
120 // an element is a "P1 interior" if it has exactly 1 exterior node
121
122 for (unsigned int k = 0; k < q1::n_chi; ++k) {
123 // Consider a Q1 element with one exterior node and three interior (or boundary)
124 // nodes. This is not an interior element by itself, but it contains an embedded P1
125 // element that *is* interior and should contribute. We need to find which of the four
126 // types of embedded P1 elements to use, but P1 elements are numbered using the node
127 // at the right angle of the reference element, not the "missing" node. Here we map
128 // "missing" nodes to "opposite" nodes, i.e. nodes used to choose P1 element types.
129 if (node_type[k] == NODE_EXTERIOR) {
130 return ElementType((k + 2) % 4);
131 }
132 }
133 }
134
135 return ELEMENT_EXTERIOR;
136}
137
138namespace q13d {
139
140Germ chi(unsigned int k, const QuadPoint &p) {
141 /* Coordinated of the nodes of the reference element: */
142 double xis[8] = {-1.0, 1.0, 1.0, -1.0, -1.0, 1.0, 1.0, -1.0};
143 double etas[8] = {-1.0, -1.0, 1.0, 1.0, -1.0, -1.0, 1.0, 1.0};
144 double zetas[8] = {-1.0, -1.0, -1.0, -1.0, 1.0, 1.0, 1.0, 1.0};
145
146 return {0.125 * (1.0 + xis[k]*p.xi) * (1.0 + etas[k]*p.eta) * (1.0 + zetas[k]*p.zeta),
147 0.125 * xis[k] * (1.0 + etas[k] * p.eta) * (1.0 + zetas[k] * p.zeta),
148 0.125 * etas[k] * (1.0 + xis[k] * p.xi) * (1.0 + zetas[k] * p.zeta),
149 0.125 * zetas[k] * (1.0 + xis[k] * p.xi) * (1.0 + etas[k] * p.eta)};
150}
151
152} // end of namespace q13d
153
154} // end of namespace fem
155} // end of namespace pism
Germ chi(unsigned int k, const QuadPoint &pt)
Linear basis functions on the interval [-1, -1].
Definition FEM.cc:30
const int n_chi
Definition FEM.hh:177
Germ chi(unsigned int k, const QuadPoint &pt)
P1 basis functions on the reference element with nodes (0,0), (1,0), (0,1).
Definition FEM.cc:90
Germ chi(unsigned int k, const QuadPoint &pt)
Definition FEM.cc:47
const int n_chi
Definition FEM.hh:184
Germ chi(unsigned int k, const QuadPoint &p)
Evaluate a Q1 shape function and its derivatives with respect to xi and eta.
Definition FEM.cc:140
const int n_chi
Definition FEM.hh:191
Germ chi(unsigned int k, const QuadPoint &pt)
Q1 basis functions on the reference element with nodes (-1,-1), (1,-1), (1,1), (-1,...
Definition FEM.cc:72
ElementType element_type(const int node_type[q1::n_chi])
Definition FEM.cc:106
ElementType
Definition FEM.hh:203
@ ELEMENT_Q
Definition FEM.hh:203
@ ELEMENT_EXTERIOR
Definition FEM.hh:205
@ NODE_EXTERIOR
Definition node_types.hh:36
static const double k
Definition exactTestP.cc:42
double dy
Function derivative with respect to y.
Definition FEM.hh:157
double val
Function value.
Definition FEM.hh:153
double dz
Function derivative with respect to z.
Definition FEM.hh:159
double dx
Function deriviative with respect to x.
Definition FEM.hh:155
Struct for gathering the value and derivative of a function at a point.
Definition FEM.hh:151