/usr/lib/python2.7/dist-packages/ufl/geometry.py is in python-ufl 2017.2.0.0-2.
This file is owned by root:root, with mode 0o644.
The actual contents of the file can be viewed below.
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 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843 844 845 846 847 848 | # -*- coding: utf-8 -*-
"Types for representing symbolic expressions for geometric quantities."
# Copyright (C) 2008-2016 Martin Sandve Alnæs
#
# This file is part of UFL.
#
# UFL is free software: you can redistribute it and/or modify
# it under the terms of the GNU Lesser General Public License as published by
# the Free Software Foundation, either version 3 of the License, or
# (at your option) any later version.
#
# UFL is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
# GNU Lesser General Public License for more details.
#
# You should have received a copy of the GNU Lesser General Public License
# along with UFL. If not, see <http://www.gnu.org/licenses/>.
#
# Modified by Anders Logg, 2009.
# Modified by Kristian B. Oelgaard, 2009
# Modified by Marie E. Rognes 2012
# Modified by Massimiliano Leoni, 2016
from ufl.utils.py23 import as_native_str
from ufl.utils.py23 import as_native_strings
from ufl.log import error
from ufl.core.ufl_type import ufl_type
from ufl.core.terminal import Terminal
from ufl.domain import as_domain
"""
Possible coordinate bootstrapping:
Xf = Xf[q]
FacetCoordinate = quadrature point on facet (ds,dS)
X = X[q]
CellCoordinate = quadrature point on cell (dx)
x = x[q]
SpatialCoordinate = quadrature point from input array (dc)
Jacobians of mappings between coordinates:
Jcf = dX/dXf = grad_Xf X(Xf)
CellFacetJacobian
Jxc = dx/dX = grad_X x(X)
Jacobian
Jxf = dx/dXf = grad_Xf x(Xf) = Jxc Jcf = dx/dX dX/dXf = grad_X x(X) grad_Xf X(Xf)
FacetJacobian = Jacobian * CellFacetJacobian
Possible computation of X from Xf:
X = Jcf Xf + X0f
CellCoordinate = CellFacetJacobian * FacetCoordinate + CellFacetOrigin
Possible computation of x from X:
x = f(X)
SpatialCoordinate = sum_k xdofs_k * xphi_k(X)
x = Jxc X + x0
SpatialCoordinate = Jacobian * CellCoordinate + CellOrigin
Possible computation of x from Xf:
x = x(X(Xf))
x = Jxf Xf + x0f
SpatialCoordinate = FacetJacobian * FacetCoordinate + FacetOrigin
Inverse relations:
X = K * (x - x0)
CellCoordinate = JacobianInverse * (SpatialCoordinate - CellOrigin)
Xf = FK * (x - x0f)
FacetCoordinate = FacetJacobianInverse * (SpatialCoordinate - FacetOrigin)
Xf = CFK * (X - X0f)
FacetCoordinate = CellFacetJacobianInverse * (CellCoordinate - CellFacetOrigin)
"""
# --- Expression node types
@ufl_type(is_abstract=True)
class GeometricQuantity(Terminal):
__slots__ = as_native_strings(("_domain",))
def __init__(self, domain):
Terminal.__init__(self)
self._domain = as_domain(domain)
def ufl_domains(self):
return (self._domain,)
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell (or over each facet for facet quantities)."
# NB! Geometric quantities are piecewise constant by
# default. Override if needed.
return True
# NB! Geometric quantities are scalar by default. Override if
# needed.
ufl_shape = ()
def _ufl_signature_data_(self, renumbering):
"Signature data of geometric quantities depend on the domain numbering."
return (self._ufl_class_.__name__,) + self._domain._ufl_signature_data_(renumbering)
def __str__(self):
return self._ufl_class_.name
def __repr__(self):
r = "%s(%s)" % (self._ufl_class_.__name__, repr(self._domain))
return as_native_str(r)
def _ufl_compute_hash_(self):
return hash((type(self).__name__,) + self._domain._ufl_hash_data_())
def __eq__(self, other):
return isinstance(other, self._ufl_class_) and other._domain == self._domain
@ufl_type(is_abstract=True)
class GeometricCellQuantity(GeometricQuantity):
__slots__ = ()
@ufl_type(is_abstract=True)
class GeometricFacetQuantity(GeometricQuantity):
__slots__ = ()
# --- Coordinate represented in different coordinate systems
@ufl_type()
class SpatialCoordinate(GeometricCellQuantity):
"""UFL geometry representation: The coordinate in a domain.
In the context of expression integration,
represents the domain coordinate of each quadrature point.
In the context of expression evaluation in a point,
represents the value of that point.
"""
__slots__ = ()
name = "x"
@property
def ufl_shape(self):
"""Return the number of coordinates defined (i.e. the geometric dimension
of the domain)."""
g = self._domain.geometric_dimension()
return (g,)
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# Only case this is true is if the domain is a vertex cell.
t = self._domain.topological_dimension()
return t == 0
def evaluate(self, x, mapping, component, index_values):
"Return the value of the coordinate."
if component == ():
if isinstance(x, (tuple, list)):
return float(x[0])
else:
return float(x)
else:
return float(x[component[0]])
@ufl_type()
class CellCoordinate(GeometricCellQuantity):
"""UFL geometry representation: The coordinate in a reference cell.
In the context of expression integration,
represents the reference cell coordinate of each quadrature point.
In the context of expression evaluation in a point in a cell,
represents that point in the reference coordinate system of the cell.
"""
__slots__ = ()
name = "X"
@property
def ufl_shape(self):
t = self._domain.topological_dimension()
return (t,)
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# Only case this is true is if the domain is a vertex cell.
t = self._domain.topological_dimension()
return t == 0
@ufl_type()
class FacetCoordinate(GeometricFacetQuantity):
"""UFL geometry representation: The coordinate in a reference cell of a facet.
In the context of expression integration over a facet,
represents the reference facet coordinate of each quadrature point.
In the context of expression evaluation in a point on a facet,
represents that point in the reference coordinate system of the facet.
"""
__slots__ = ()
name = "Xf"
def __init__(self, domain):
GeometricFacetQuantity.__init__(self, domain)
t = self._domain.topological_dimension()
if t < 2:
error("FacetCoordinate is only defined for topological dimensions >= 2.")
@property
def ufl_shape(self):
t = self._domain.topological_dimension()
return (t-1,)
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# Only case this is true is if the domain is an interval cell
# (with a vertex facet).
t = self._domain.topological_dimension()
return t <= 1
# --- Origin of coordinate systems in larger coordinate systems
@ufl_type()
class CellOrigin(GeometricCellQuantity):
"""UFL geometry representation: The spatial coordinate corresponding to origin of a reference cell."""
__slots__ = ()
name = "x0"
@property
def ufl_shape(self):
g = self._domain.geometric_dimension()
return (g,)
def is_cellwise_constant(self):
return True
@ufl_type()
class FacetOrigin(GeometricFacetQuantity):
"""UFL geometry representation: The spatial coordinate corresponding to origin of a reference facet."""
__slots__ = ()
name = "x0f"
@property
def ufl_shape(self):
g = self._domain.geometric_dimension()
return (g,)
@ufl_type()
class CellFacetOrigin(GeometricFacetQuantity):
"""UFL geometry representation: The reference cell coordinate corresponding to origin of a reference facet."""
__slots__ = ()
name = "X0f"
@property
def ufl_shape(self):
t = self._domain.topological_dimension()
return (t,)
# --- Jacobians of mappings between coordinate systems
@ufl_type()
class Jacobian(GeometricCellQuantity):
"""UFL geometry representation: The Jacobian of the mapping from reference cell to spatial coordinates.
.. math:: J_{ij} = \\frac{dx_i}{dX_j}
"""
__slots__ = ()
name = "J"
@property
def ufl_shape(self):
"""Return the number of coordinates defined (i.e. the geometric dimension
of the domain)."""
g = self._domain.geometric_dimension()
t = self._domain.topological_dimension()
return (g, t)
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# Only true for a piecewise linear coordinate field in simplex cells
return self._domain.is_piecewise_linear_simplex_domain()
@ufl_type()
class FacetJacobian(GeometricFacetQuantity):
"""UFL geometry representation: The Jacobian of the mapping from reference facet to spatial coordinates.
FJ_ij = dx_i/dXf_j
The FacetJacobian is the product of the Jacobian and CellFacetJacobian:
FJ = dx/dXf = dx/dX dX/dXf = J * CFJ
"""
__slots__ = ()
name = "FJ"
def __init__(self, domain):
GeometricFacetQuantity.__init__(self, domain)
t = self._domain.topological_dimension()
if t < 2:
error("FacetJacobian is only defined for topological dimensions >= 2.")
@property
def ufl_shape(self):
g = self._domain.geometric_dimension()
t = self._domain.topological_dimension()
return (g, t-1)
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# Only true for a piecewise linear coordinate field in simplex
# cells
return self._domain.is_piecewise_linear_simplex_domain()
@ufl_type()
class CellFacetJacobian(GeometricFacetQuantity): # dX/dXf
"""UFL geometry representation: The Jacobian of the mapping from reference facet to reference cell coordinates.
CFJ_ij = dX_i/dXf_j
"""
__slots__ = ()
name = "CFJ"
def __init__(self, domain):
GeometricFacetQuantity.__init__(self, domain)
t = self._domain.topological_dimension()
if t < 2:
error("CellFacetJacobian is only defined for topological dimensions >= 2.")
@property
def ufl_shape(self):
t = self._domain.topological_dimension()
return (t, t-1)
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# This is always a constant mapping between two reference
# coordinate systems.
return True
@ufl_type()
class ReferenceCellEdgeVectors(GeometricCellQuantity):
"""UFL geometry representation: The vectors between reference cell vertices for each edge in cell."""
__slots__ = ()
name = "RCEV"
def __init__(self, domain):
GeometricCellQuantity.__init__(self, domain)
t = self._domain.topological_dimension()
if t < 2:
error("CellEdgeVectors is only defined for topological dimensions >= 2.")
@property
def ufl_shape(self):
cell = self.ufl_domain().ufl_cell()
ne = cell.num_edges()
t = cell.topological_dimension()
return (ne, t)
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# This is always constant for a given cell type
return True
@ufl_type()
class ReferenceFacetEdgeVectors(GeometricFacetQuantity):
"""UFL geometry representation: The vectors between reference cell vertices for each edge in current facet."""
__slots__ = ()
name = "RFEV"
def __init__(self, domain):
GeometricFacetQuantity.__init__(self, domain)
t = self._domain.topological_dimension()
if t < 3:
error("FacetEdgeVectors is only defined for topological dimensions >= 3.")
@property
def ufl_shape(self):
cell = self.ufl_domain().ufl_cell()
nfe = cell.num_facet_edges()
t = cell.topological_dimension()
return (nfe, t)
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# This is always constant for a given cell type
return True
@ufl_type()
class CellVertices(GeometricCellQuantity):
"""UFL geometry representation: Physical cell vertices."""
__slots__ = ()
name = "CV"
def __init__(self, domain):
GeometricCellQuantity.__init__(self, domain)
@property
def ufl_shape(self):
cell = self.ufl_domain().ufl_cell()
nv = cell.num_vertices()
g = cell.geometric_dimension()
return (nv, g)
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# This is always constant for a given cell type
return True
@ufl_type()
class CellEdgeVectors(GeometricCellQuantity):
"""UFL geometry representation: The vectors between physical cell vertices for each edge in cell."""
__slots__ = ()
name = "CEV"
def __init__(self, domain):
GeometricCellQuantity.__init__(self, domain)
t = self._domain.topological_dimension()
if t < 2:
error("CellEdgeVectors is only defined for topological dimensions >= 2.")
@property
def ufl_shape(self):
cell = self.ufl_domain().ufl_cell()
ne = cell.num_edges()
g = cell.geometric_dimension()
return (ne, g)
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# This is always constant for a given cell type
return True
@ufl_type()
class FacetEdgeVectors(GeometricFacetQuantity):
"""UFL geometry representation: The vectors between physical cell vertices for each edge in current facet."""
__slots__ = ()
name = "FEV"
def __init__(self, domain):
GeometricFacetQuantity.__init__(self, domain)
t = self._domain.topological_dimension()
if t < 3:
error("FacetEdgeVectors is only defined for topological dimensions >= 3.")
@property
def ufl_shape(self):
cell = self.ufl_domain().ufl_cell()
nfe = cell.num_facet_edges()
g = cell.geometric_dimension()
return (nfe, g)
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# This is always constant for a given cell type
return True
# --- Determinants (signed or pseudo) of geometry mapping Jacobians
@ufl_type()
class JacobianDeterminant(GeometricCellQuantity):
"""UFL geometry representation: The determinant of the Jacobian.
Represents the signed determinant of a square Jacobian or the pseudo-determinant of a non-square Jacobian.
"""
__slots__ = ()
name = "detJ"
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# Only true for a piecewise linear coordinate field in simplex
# cells
return self._domain.is_piecewise_linear_simplex_domain()
@ufl_type()
class FacetJacobianDeterminant(GeometricFacetQuantity):
"""UFL geometry representation: The pseudo-determinant of the FacetJacobian."""
__slots__ = ()
name = "detFJ"
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# Only true for a piecewise linear coordinate field in simplex cells
return self._domain.is_piecewise_linear_simplex_domain()
@ufl_type()
class CellFacetJacobianDeterminant(GeometricFacetQuantity):
"""UFL geometry representation: The pseudo-determinant of the CellFacetJacobian."""
__slots__ = ()
name = "detCFJ"
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# Only true for a piecewise linear coordinate field in simplex
# cells
return self._domain.is_piecewise_linear_simplex_domain()
# --- Inverses (signed or pseudo) of geometry mapping Jacobians
@ufl_type()
class JacobianInverse(GeometricCellQuantity):
"""UFL geometry representation: The inverse of the Jacobian.
Represents the inverse of a square Jacobian or the pseudo-inverse of a non-square Jacobian.
"""
__slots__ = ()
name = "K"
@property
def ufl_shape(self):
"""Return the number of coordinates defined (i.e. the geometric dimension
of the domain)."""
g = self._domain.geometric_dimension()
t = self._domain.topological_dimension()
return (t, g)
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# Only true for a piecewise linear coordinate field in simplex
# cells
return self._domain.is_piecewise_linear_simplex_domain()
@ufl_type()
class FacetJacobianInverse(GeometricFacetQuantity):
"""UFL geometry representation: The pseudo-inverse of the FacetJacobian."""
__slots__ = ()
name = "FK"
def __init__(self, domain):
GeometricFacetQuantity.__init__(self, domain)
t = self._domain.topological_dimension()
if t < 2:
error("FacetJacobianInverse is only defined for topological dimensions >= 2.")
@property
def ufl_shape(self):
g = self._domain.geometric_dimension()
t = self._domain.topological_dimension()
return (t-1, g)
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# Only true for a piecewise linear coordinate field in simplex
# cells
return self._domain.is_piecewise_linear_simplex_domain()
@ufl_type()
class CellFacetJacobianInverse(GeometricFacetQuantity):
"""UFL geometry representation: The pseudo-inverse of the CellFacetJacobian."""
__slots__ = ()
name = "CFK"
def __init__(self, domain):
GeometricFacetQuantity.__init__(self, domain)
t = self._domain.topological_dimension()
if t < 2:
error("CellFacetJacobianInverse is only defined for topological dimensions >= 2.")
@property
def ufl_shape(self):
t = self._domain.topological_dimension()
return (t-1, t)
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# Only true for a piecewise linear coordinate field in simplex cells
return self._domain.is_piecewise_linear_simplex_domain()
# --- Types representing normal or tangent vectors
@ufl_type()
class FacetNormal(GeometricFacetQuantity):
"""UFL geometry representation: The outwards pointing normal vector of the current facet."""
__slots__ = ()
name = "n"
@property
def ufl_shape(self):
"""Return the number of coordinates defined (i.e. the geometric dimension
of the domain)."""
g = self._domain.geometric_dimension()
return (g,)
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# For product cells, this is only true for some but not all
# facets. Seems like too much work to fix right now. Only
# true for a piecewise linear coordinate field with simplex
# _facets_.
is_piecewise_linear = self._domain.ufl_coordinate_element().degree() == 1
return is_piecewise_linear and self._domain.ufl_cell().has_simplex_facets()
@ufl_type()
class CellNormal(GeometricCellQuantity):
"""UFL geometry representation: The upwards pointing normal vector of the current manifold cell."""
__slots__ = ()
name = "cell_normal"
@property
def ufl_shape(self):
"""Return the number of coordinates defined (i.e. the geometric dimension
of the domain)."""
g = self._domain.geometric_dimension()
# t = self._domain.topological_dimension()
# return (g-t,g) # TODO: Should it be CellNormals? For interval in 3D we have two!
return (g,)
@ufl_type()
class ReferenceNormal(GeometricFacetQuantity):
"""UFL geometry representation: The outwards pointing normal vector of the current facet on the reference cell"""
__slots__ = ()
name = "reference_normal"
@property
def ufl_shape(self):
t = self._domain.topological_dimension()
return (t,)
# TODO: Implement in apply_geometry_lowering and enable
# @ufl_type()
# class FacetTangents(GeometricFacetQuantity):
# """UFL geometry representation: The tangent vectors of the current facet."""
# __slots__ = ()
# name = "t"
#
# def __init__(self, domain):
# GeometricFacetQuantity.__init__(self, domain)
# t = self._domain.topological_dimension()
# if t < 2:
# error("FacetTangents is only defined for topological dimensions >= 2.")
#
# @property
# def ufl_shape(self):
# g = self._domain.geometric_dimension()
# t = self._domain.topological_dimension()
# return (t-1,g)
#
# def is_cellwise_constant(self): # NB! Copied from FacetNormal
# "Return whether this expression is spatially constant over each cell."
# # For product cells, this is only true for some but not all facets. Seems like too much work to fix right now.
# # Only true for a piecewise linear coordinate field with simplex _facets_.
# is_piecewise_linear = self._domain.ufl_coordinate_element().degree() == 1
# return is_piecewise_linear and self._domain.ufl_cell().has_simplex_facets()
# TODO: Implement in apply_geometry_lowering and enable
# @ufl_type()
# class CellTangents(GeometricCellQuantity):
# """UFL geometry representation: The tangent vectors of the current manifold cell."""
# __slots__ = ()
# name = "cell_tangents"
#
# @property
# def ufl_shape(self):
# g = self._domain.geometric_dimension()
# t = self._domain.topological_dimension()
# return (t,g)
# --- Types representing midpoint coordinates
# TODO: Implement in the rest of fenics
# @ufl_type()
# class CellMidpoint(GeometricCellQuantity):
# """UFL geometry representation: The midpoint coordinate of the current cell."""
# __slots__ = ()
# name = "cell_midpoint"
#
# @property
# def ufl_shape(self):
# g = self._domain.geometric_dimension()
# return (g,)
# TODO: Implement in the rest of fenics
# @ufl_type()
# class FacetMidpoint(GeometricFacetQuantity):
# """UFL geometry representation: The midpoint coordinate of the current facet."""
# __slots__ = ()
# name = "facet_midpoint"
#
# @property
# def ufl_shape(self):
# g = self._domain.geometric_dimension()
# return (g,)
# --- Types representing measures of the cell and entities of the cell, typically used for stabilisation terms
# TODO: Clean up this set of types? Document!
@ufl_type()
class ReferenceCellVolume(GeometricCellQuantity):
"""UFL geometry representation: The volume of the reference cell."""
__slots__ = ()
name = "reference_cell_volume"
@ufl_type()
class ReferenceFacetVolume(GeometricFacetQuantity):
"""UFL geometry representation: The volume of the reference cell of the current facet."""
__slots__ = ()
name = "reference_facet_volume"
@ufl_type()
class CellVolume(GeometricCellQuantity):
"""UFL geometry representation: The volume of the cell."""
__slots__ = ()
name = "volume"
@ufl_type()
class Circumradius(GeometricCellQuantity):
"""UFL geometry representation: The circumradius of the cell."""
__slots__ = ()
name = "circumradius"
@ufl_type()
class CellDiameter(GeometricCellQuantity):
"""UFL geometry representation: The diameter of the cell, i.e.,
maximal distance of two points in the cell."""
__slots__ = ()
name = "diameter"
# @ufl_type()
# class CellSurfaceArea(GeometricCellQuantity):
# """UFL geometry representation: The total surface area of the cell."""
# __slots__ = ()
# name = "surfacearea"
@ufl_type()
class FacetArea(GeometricFacetQuantity): # FIXME: Should this be allowed for interval domain?
"""UFL geometry representation: The area of the facet."""
__slots__ = ()
name = "facetarea"
@ufl_type()
class MinCellEdgeLength(GeometricCellQuantity):
"""UFL geometry representation: The minimum edge length of the cell."""
__slots__ = ()
name = "mincelledgelength"
@ufl_type()
class MaxCellEdgeLength(GeometricCellQuantity):
"""UFL geometry representation: The maximum edge length of the cell."""
__slots__ = ()
name = "maxcelledgelength"
@ufl_type()
class MinFacetEdgeLength(GeometricFacetQuantity):
"""UFL geometry representation: The minimum edge length of the facet."""
__slots__ = ()
name = "minfacetedgelength"
@ufl_type()
class MaxFacetEdgeLength(GeometricFacetQuantity):
"""UFL geometry representation: The maximum edge length of the facet."""
__slots__ = ()
name = "maxfacetedgelength"
# --- Types representing other stuff
@ufl_type()
class CellOrientation(GeometricCellQuantity):
"""UFL geometry representation: The orientation (+1/-1) of the current cell.
For non-manifold cells (tdim == gdim), this equals the sign
of the Jacobian determinant, i.e. +1 if the physical cell is
oriented the same way as the reference cell and -1 otherwise.
For manifold cells of tdim==gdim-1 this is input data belonging
to the mesh, used to distinguish between the sides of the manifold.
"""
__slots__ = ()
name = "cell_orientation"
@ufl_type()
class FacetOrientation(GeometricFacetQuantity):
"""UFL geometry representation: The orientation (+1/-1) of the current facet relative to the reference cell."""
__slots__ = ()
name = "facet_orientation"
# This doesn't quite fit anywhere. Make a special set of symbolic
# terminal types instead?
@ufl_type()
class QuadratureWeight(GeometricQuantity):
"""UFL geometry representation: The current quadrature weight.
Only used inside a quadrature context.
"""
__slots__ = ()
name = "weight"
def is_cellwise_constant(self):
"Return whether this expression is spatially constant over each cell."
# The weight usually varies with the quadrature points
return False
|