Thanks to visit codestin.com
Credit goes to github.com

Skip to content
Draft
Next Next commit
start Scaled Contravariant Piola Map
  • Loading branch information
amboschman committed Jun 26, 2025
commit d491b71f13f07c7be628b41f6bf267259c50598e
14 changes: 11 additions & 3 deletions src/FESpaces/DivConformingFESpaces.jl
Original file line number Diff line number Diff line change
Expand Up @@ -29,7 +29,8 @@ struct TransformRTDofBasis{Dc,Dp} <: Map end ;
function get_cell_dof_basis(model::DiscreteModel,
cell_reffe::AbstractArray{<:GenericRefFE{<:DivConforming}},
::DivConformity,
sign_flip=get_sign_flip(model, cell_reffe))
map_type=nothing)
sign_flip = get_sign_flip(model, cell_reffe)
cell_map = get_cell_map(Triangulation(model))
phi = cell_map[1]
Jt = lazy_map(Broadcasting(∇),cell_map)
Expand All @@ -48,9 +49,16 @@ end
function get_cell_shapefuns(model::DiscreteModel,
cell_reffe::AbstractArray{<:GenericRefFE{<:DivConforming}},
::DivConformity,
sign_flip=get_sign_flip(model, cell_reffe))
map_type=nothing)
sign_flip=get_sign_flip(model, cell_reffe)
cell_reffe_shapefuns=lazy_map(get_shapefuns,cell_reffe)
k=ContraVariantPiolaMap()

if map_type == nothing
k=ContraVariantPiolaMap()
else
k=map_type
end

lazy_map(k,
cell_reffe_shapefuns,
get_cell_map(Triangulation(model)),
Expand Down
10 changes: 8 additions & 2 deletions src/FESpaces/FESpaceFactories.jl
Original file line number Diff line number Diff line change
Expand Up @@ -71,7 +71,8 @@ function FESpace(
dirichlet_tags=Int[],
dirichlet_masks=nothing,
constraint=nothing,
vector_type=nothing)
vector_type=nothing,
map_type=nothing)

conf = Conformity(testitem(cell_reffe),conformity)

Expand All @@ -89,8 +90,13 @@ function FESpace(
trian)
return V
end

if map_type == nothing
cell_fe = CellFE(model,cell_reffe,conf)
else
cell_fe = CellFE(model,cell_reffe,conf,map_type)
end

cell_fe = CellFE(model,cell_reffe,conf)
_vector_type = _get_vector_type(vector_type,cell_fe,trian)
if conformity in (L2Conformity(),:L2) && dirichlet_tags == Int[]
F = _DiscontinuousFESpace(_vector_type,trian,cell_fe)
Expand Down
34 changes: 32 additions & 2 deletions src/ReferenceFEs/RaviartThomasRefFEs.jl
Original file line number Diff line number Diff line change
Expand Up @@ -375,6 +375,7 @@ function _eval_moment_dof_basis!(dofs,vals::AbstractMatrix,b)
end

struct ContraVariantPiolaMap <: Map end
struct ScaledContraVariantPiolaMap <: Map end

function evaluate!(
cache,
Expand All @@ -387,6 +388,17 @@ function evaluate!(
Broadcasting(Operation(k))(∇v,Jt,detJ,sign_flip)
end

function evaluate!(
cache,
::Broadcasting{typeof(∇)},
a::Fields.BroadcastOpFieldArray{ScaledContraVariantPiolaMap})
v, Jt, detJ,sign_flip = a.args
# Assuming J comes from an affine map
∇v = Broadcasting(∇)(v)
k = ContraVariantPiolaMap()
Broadcasting(Operation(k))(∇v,Jt,detJ,sign_flip)
end

function lazy_map(
::Broadcasting{typeof(gradient)},
a::LazyArray{<:Fill{Broadcasting{Operation{ContraVariantPiolaMap}}}})
Expand All @@ -396,6 +408,15 @@ function lazy_map(
lazy_map(Broadcasting(Operation(k)),∇v,Jt,detJ,sign_flip)
end

function lazy_map(
::Broadcasting{typeof(gradient)},
a::LazyArray{<:Fill{Broadcasting{Operation{ScaledContraVariantPiolaMap}}}})
v, Jt, detJ,sign_flip = a.args
∇v = lazy_map(Broadcasting(∇),v)
k = ScaledContraVariantPiolaMap()
lazy_map(Broadcasting(Operation(k)),∇v,Jt,detJ,sign_flip)
end

function evaluate!(cache,::ContraVariantPiolaMap,
v::Number,
Jt::Number,
Expand All @@ -404,6 +425,15 @@ function evaluate!(cache,::ContraVariantPiolaMap,
((-1)^sign_flip*v)⋅((1/detJ)*Jt)
end

function evaluate!(cache,::ScaledContraVariantPiolaMap,
v::Number,
Jt::Number,
detJ::Number,
sign_flip::Bool)
D = size(Jt)[1]
((-1)^sign_flip*v)⋅((1/(detJ^(1/D)))*Jt)
end

function evaluate!(cache,
k::ContraVariantPiolaMap,
v::AbstractVector{<:Field},
Expand All @@ -415,7 +445,7 @@ function evaluate!(cache,
end

function lazy_map(
k::ContraVariantPiolaMap,
k::Union{ContraVariantPiolaMap,ScaledContraVariantPiolaMap},
cell_ref_shapefuns::AbstractArray{<:AbstractArray{<:Field}},
cell_map::AbstractArray{<:Field},
sign_flip::AbstractArray{<:AbstractArray{<:Field}})
Expand All @@ -424,4 +454,4 @@ function lazy_map(
cell_detJ = lazy_map(Operation(meas),cell_Jt)

lazy_map(Broadcasting(Operation(k)),cell_ref_shapefuns,cell_Jt,cell_detJ,sign_flip)
end
end
1 change: 1 addition & 0 deletions src/ReferenceFEs/ReferenceFEs.jl
Original file line number Diff line number Diff line change
Expand Up @@ -97,6 +97,7 @@ export TET_AXIS
export INVALID_PERM

export ContraVariantPiolaMap
export ScaledContraVariantPiolaMap

export Dof
export get_nodes
Expand Down