Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Add H(curl div) sobolev space and the covariant contravariant Piola mapping #312

Open
wants to merge 2 commits into
base: master
Choose a base branch
from
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
1 change: 1 addition & 0 deletions tsfc/finatinterface.py
Original file line number Diff line number Diff line change
Expand Up @@ -61,6 +61,7 @@
"Nedelec 2nd kind H(curl)": finat.NedelecSecondKind,
"Raviart-Thomas": finat.RaviartThomas,
"Regge": finat.Regge,
"Gopalakrishnan-Lederer-Schoberl": finat.GopalakrishnanLedererSchoberl,
"BDMCE": finat.BrezziDouglasMariniCubeEdge,
"BDMCF": finat.BrezziDouglasMariniCubeFace,
# These require special treatment below
Expand Down
11 changes: 11 additions & 0 deletions tsfc/ufl_utils.py
Original file line number Diff line number Diff line change
Expand Up @@ -403,6 +403,10 @@ def apply_mapping(expression, element, domain):

G(X) = det(J)^2 K g(x) K^T i.e. G_il(X)=(detJ)^2 K_ij g_jk K_lk

'covariant contravariant piola' mapping for g:

G(X) = det(J) J^T g(x) K^T i.e. G_il(X) = det(J) J_ji g_jk(x) K_lk

If 'contravariant piola' or 'covariant piola' (or their double
variants) are applied to a matrix-valued function, the appropriate
mappings are applied row-by-row.
Expand Down Expand Up @@ -442,6 +446,13 @@ def apply_mapping(expression, element, domain):
*k, i, j, m, n = indices(len(expression.ufl_shape) + 2)
kmn = (*k, m, n)
rexpression = as_tensor(detJ**2 * K[i, m] * expression[kmn] * K[j, n], (*k, i, j))
elif mapping == "covariant contravariant piola":
J = Jacobian(mesh)
K = JacobianInverse(mesh)
detJ = JacobianDeterminant(mesh)
*k, i, j, m, n = indices(len(expression.ufl_shape) + 2)
kmn = (*k, m, n)
rexpression = as_tensor(detJ * J[m, i] * expression[kmn] * K[j, n], (*k, i, j))
elif mapping == "symmetries":
# This tells us how to get from the pieces of the reference
# space expression to the physical space one.
Expand Down
Loading