-
Notifications
You must be signed in to change notification settings - Fork 86
/
Copy pathboundary_facet_basis.py
225 lines (190 loc) · 8.21 KB
/
boundary_facet_basis.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
import logging
from typing import Callable, Optional, Tuple, Any
import numpy as np
from numpy import ndarray
from skfem.element import (BOUNDARY_ELEMENT_MAP, DiscreteField, Element,
ElementHex0, ElementQuad0, ElementTetP0,
ElementTriP0)
from skfem.mapping import Mapping
from skfem.mesh import Mesh, MeshHex, MeshLine, MeshQuad, MeshTet, MeshTri
from .abstract_basis import AbstractBasis
from .cell_basis import CellBasis
from ..dofs import Dofs
logger = logging.getLogger(__name__)
class BoundaryFacetBasis(AbstractBasis):
"""For fields defined on the boundary of the domain."""
def __init__(self,
mesh: Mesh,
elem: Element,
mapping: Optional[Mapping] = None,
intorder: Optional[int] = None,
quadrature: Optional[Tuple[ndarray, ndarray]] = None,
facets: Optional[Any] = None,
dofs: Optional[Dofs] = None):
"""Precomputed global basis on boundary facets.
Parameters
----------
mesh
An object of type :class:`~skfem.mesh.Mesh`.
elem
An object of type :class:`~skfem.element.Element`.
mapping
An object of type :class:`skfem.mapping.Mapping`. If `None`, uses
`mesh.mapping`.
intorder
Optional integration order, i.e. the degree of polynomials that are
integrated exactly by the used quadrature. Not used if `quadrature`
is specified.
quadrature
Optional tuple of quadrature points and weights.
facets
Optional subset of facet indices.
dofs
Optional :class:`~skfem.assembly.Dofs` object.
"""
typestr = ("{}({}, {})".format(type(self).__name__,
type(mesh).__name__,
type(elem).__name__))
logger.info("Initializing {}".format(typestr))
super(BoundaryFacetBasis, self).__init__(
mesh,
elem,
mapping,
intorder,
quadrature,
mesh.brefdom,
dofs,
)
# initialize indices
if isinstance(facets, tuple):
self.find, self.tind, self.tind_normals = facets
else:
if facets is None:
self.find = np.nonzero(self.mesh.f2t[1] == -1)[0]
else:
self.find = mesh.normalize_facets(facets)
self.tind = self.mesh.f2t[0, self.find]
self.tind_normals = self.mesh.f2t[0, self.find]
if len(self.find) == 0:
logger.warning("Initializing {} with no facets.".format(typestr))
# boundary refdom to global facet
x = self.mapping.G(self.X, find=self.find)
# global facet to refdom facet
Y = self.mapping.invF(x, tind=self.tind)
# calculate normals
Y0 = self.mapping.invF(x, tind=self.tind_normals)
self.normals = DiscreteField(
value=self.mapping.normals(Y0,
self.tind_normals,
self.find,
self.mesh.t2f)
)
self.nelems = len(self.find)
self.basis = [self.elem.gbasis(self.mapping, Y, j, tind=self.tind)
for j in range(self.Nbfun)]
self.dx = (np.abs(self.mapping.detDG(self.X, find=self.find))
* np.tile(self.W, (self.nelems, 1)))
logger.info("Initializing finished.")
@classmethod
def init_subdomain(cls, mesh, elem, elements=None, side=0, **kwargs):
assert elements is not None
elements = mesh.normalize_elements(elements)
subdomain_facets, sub_counts = np.unique(mesh.t2f[:, elements],
return_counts=True)
facets = subdomain_facets[sub_counts == 1]
tind = mesh.f2t[:, facets].flatten('F')
ix = np.isin(tind, elements)
tind = tind[ix] if side == 0 else tind[~ix]
return cls(mesh, elem, facets=(facets, tind, tind), **kwargs)
def default_parameters(self):
"""Return default parameters for `~skfem.assembly.asm`."""
return {
'x': self.global_coordinates(),
'h': self.mesh_parameters(),
'n': self.normals,
}
def global_coordinates(self) -> DiscreteField:
return DiscreteField(self.mapping.G(self.X, find=self.find))
def mesh_parameters(self) -> DiscreteField:
return DiscreteField((np.abs(self.mapping.detDG(self.X, self.find))
** (1. / (self.mesh.dim() - 1.)))
if self.mesh.dim() != 1 else np.array([0.]))
def with_element(self, elem: Element) -> 'BoundaryFacetBasis':
"""Return a similar basis using a different element."""
return type(self)(
self.mesh,
elem,
mapping=self.mapping,
quadrature=self.quadrature,
facets=(self.find, self.tind, self.tind_normals),
)
def project(self, interp, facets=None):
"""Perform :math:`L^2` projection onto the basis on the boundary.
The values of the interior DOFs remain zero.
Parameters
----------
interp
An object of type :class:`~skfem.element.DiscreteField` which is a
function (to be projected) evaluated at global quadrature points at
the boundary of the domain. If a function is given, then
:class:`~skfem.element.DiscreteField` is created by passing
an array of global quadrature point locations to the function.
facets
Optionally perform the projection on a subset of facets. The
values of the remaining DOFs are zero.
"""
from skfem.utils import solve, condense
M, f = self._projection(interp)
if facets is not None:
return solve(*condense(M, f, I=self.get_dofs(facets=facets)))
return solve(*condense(M, f, I=self.get_dofs(facets=self.find)))
def _trace_project(self,
x: ndarray,
elem: Element) -> ndarray:
from skfem.utils import projection
fbasis = BoundaryFacetBasis(self.mesh,
elem,
facets=self.find,
quadrature=(self.X, self.W))
I = fbasis.get_dofs(self.find).all()
if len(I) == 0: # special case: no facet DOFs
if fbasis.dofs.interior_dofs is not None:
if fbasis.dofs.interior_dofs.shape[0] > 1:
# no one-to-one restriction: requires interpolation
raise NotImplementedError
# special case: piecewise constant elem
I = fbasis.dofs.interior_dofs[:, self.tind].flatten()
else:
raise ValueError
return projection(x, fbasis, self, I=I)
def trace(self,
x: ndarray,
projection: Callable[[ndarray], ndarray],
target_elem: Optional[Element] = None) -> Tuple[CellBasis,
ndarray]:
DEFAULT_TARGET = {
MeshTri: ElementTriP0,
MeshQuad: ElementQuad0,
MeshTet: ElementTetP0,
MeshHex: ElementHex0,
}
meshcls = type(self.mesh)
if meshcls not in DEFAULT_TARGET:
raise NotImplementedError("Mesh type not supported.")
if target_elem is None:
target_elem = DEFAULT_TARGET[meshcls]()
if type(target_elem) not in BOUNDARY_ELEMENT_MAP:
raise Exception("The specified element not supported.")
elemcls = BOUNDARY_ELEMENT_MAP[type(target_elem)]
target_meshcls = {
MeshTri: MeshLine,
MeshQuad: MeshLine,
MeshTet: MeshTri,
MeshHex: MeshQuad,
}[meshcls]
assert callable(target_meshcls) # to satisfy mypy
p, t = self.mesh._reix(self.mesh.facets[:, self.find])
return (
CellBasis(target_meshcls(projection(p), t), elemcls()),
self._trace_project(x, target_elem)
)