Click here to read what the information on this page means.
This element is implemented in
FIAT ,
Symfem , and
(legacy) UFL.
↓ Show implementation detail ↓↑ Hide implementation detail ↑FIAT | FIAT.KongMulderVeldhuizen ↓ Show FIAT examples ↓↑ Hide FIAT examples ↑Before running this example, you must install FIAT: pip3 install git+https://github.com/firedrakeproject/fiat.git This element can then be created with the following lines of Python: import FIAT
# Create Kong-Mulder-Veldhuizen degree 1 element = FIAT.KongMulderVeldhuizen(FIAT.ufc_cell("triangle"), 1)
# Create Kong-Mulder-Veldhuizen degree 2 element = FIAT.KongMulderVeldhuizen(FIAT.ufc_cell("triangle"), 2)
# Create Kong-Mulder-Veldhuizen degree 1 element = FIAT.KongMulderVeldhuizen(FIAT.ufc_cell("tetrahedron"), 1) This implementation is correct for all the examples below. |
Symfem | "KMV" ↓ Show Symfem examples ↓↑ Hide Symfem examples ↑Before running this example, you must install Symfem: pip3 install symfem This element can then be created with the following lines of Python: import symfem
# Create Kong-Mulder-Veldhuizen degree 1 on a triangle element = symfem.create_element("triangle", "KMV", 1)
# Create Kong-Mulder-Veldhuizen degree 2 on a triangle element = symfem.create_element("triangle", "KMV", 2)
# Create Kong-Mulder-Veldhuizen degree 1 on a tetrahedron element = symfem.create_element("tetrahedron", "KMV", 1) This implementation is used to compute the examples below and verify other implementations. |
(legacy) UFL | "KMV" ↓ Show (legacy) UFL examples ↓↑ Hide (legacy) UFL examples ↑Before running this example, you must install (legacy) UFL: pip3 install setuptools pip3 install fenics-ufl-legacy This element can then be created with the following lines of Python: import ufl_legacy
# Create Kong-Mulder-Veldhuizen degree 1 on a triangle element = ufl_legacy.FiniteElement("KMV", "triangle", 1)
# Create Kong-Mulder-Veldhuizen degree 2 on a triangle element = ufl_legacy.FiniteElement("KMV", "triangle", 2)
# Create Kong-Mulder-Veldhuizen degree 1 on a tetrahedron element = ufl_legacy.FiniteElement("KMV", "tetrahedron", 1) |
- Chin-Joe-Kong, M. J. S., Mulder, Wim A., and Van Veldhuizen, M. Higher-order triangular and tetrahedral finite elements with mass lumping for solving the wave equation, Journal of Engineering Mathematics 35, 405–426, 1999. [DOI: 10.1023/A:1004420829610] [BibTeX]
Element added | 09 May 2021 |
Element last updated | 27 September 2024 |