Skip to content

Commit

Permalink
Trying to optimize
Browse files Browse the repository at this point in the history
function get_cell_dof_basis(model::DiscreteModel,
                            cell_reffe::AbstractArray{<:GenericRefFE{RaviartThomas}},
                            ::DivConformity)

To this end:
  (1) Replaced a function by a new Map instance.
  (2) I early evaluated the Jacobian at all integration points
      (instead of on a cell by cell basis).

Did not observe signficant performance improvements, though.
  • Loading branch information
amartinhuertas committed Aug 14, 2021
1 parent 2a3874f commit a250431
Showing 1 changed file with 30 additions and 36 deletions.
66 changes: 30 additions & 36 deletions src/FESpaces/DivConformingFESpaces.jl
Original file line number Diff line number Diff line change
Expand Up @@ -24,14 +24,24 @@
# shape function corresponding to the global DoF.
# * We do NOT have to use the signed determinant, but its absolute value, in the Piola Map.

struct TransformRTDofBasis{Dc,Dp} <: Map end ;

function get_cell_dof_basis(model::DiscreteModel,
cell_reffe::AbstractArray{<:GenericRefFE{RaviartThomas}},
::DivConformity)
cell_reffe::AbstractArray{<:GenericRefFE{RaviartThomas}},
::DivConformity)
sign_flip = get_sign_flip(model, cell_reffe)
lazy_map(_transform_rt_dof_basis,
cell_reffe,
get_cell_map(Triangulation(model)),
sign_flip)
cell_map = get_cell_map(Triangulation(model))
phi = cell_map[1]
Jt = lazy_map(Broadcasting(∇),cell_map)
x = lazy_map(get_nodes,lazy_map(get_dof_basis,cell_reffe))
Jtx = lazy_map(evaluate,Jt,x)
reffe = cell_reffe[1]
Dc = num_dims(reffe)
et = return_type(get_prebasis(reffe))
pt = Point{Dc,et}
Dp = first(size(return_type(phi,zero(pt))))
k = TransformRTDofBasis{Dc,Dp}()
lazy_map(k,cell_reffe,Jtx,sign_flip)
end

function get_cell_shapefuns(model::DiscreteModel,
Expand Down Expand Up @@ -112,17 +122,10 @@ function get_sign_flip(model::DiscreteModel,
get_cell_to_bgcell(model))
end

function _transform_rt_dof_basis(reffe::GenericRefFE{RaviartThomas},
phi::Field,
sign_flip::AbstractVector{Bool})
cache = return_cache(_transform_rt_dof_basis,reffe,phi,sign_flip)
evaluate!(cache,_transform_rt_dof_basis,reffe,phi,sign_flip)
end

function return_cache(::typeof(_transform_rt_dof_basis),
function return_cache(::TransformRTDofBasis{Dc,Dp},
reffe::GenericRefFE{RaviartThomas},
phi::Field,
::AbstractVector{Bool})
Jtx,
::AbstractVector{Bool}) where {Dc,Dp}
p = get_polytope(reffe)
prebasis = get_prebasis(reffe)
order = get_order(prebasis)
Expand All @@ -132,47 +135,38 @@ function return_cache(::typeof(_transform_rt_dof_basis),
get_face_nodes_dofs(dofs),
get_face_moments(dofs)
db = MomentBasedDofBasis(nodes,nf_moments,nf_nodes)
face_moments = [ similar(i,VectorValue{Dp,et}) for i in nf_moments ]

# Determine the type of the elements
# in the range of Field. This is the
# element type of the face_moments array.
D = num_dims(reffe)
pt = Point{D,et}
T=return_type(phi,zero(pt))
face_moments = [ similar(i,T) for i in nf_moments ]

Jt_q_cache = return_cache((phi),db.nodes)
cache = (db.nodes, db.face_nodes, nf_moments, face_moments, Jt_q_cache)
cache = (db.nodes, db.face_nodes, nf_moments, face_moments)
cache
end


function evaluate!(cache,
::typeof(_transform_rt_dof_basis),
::TransformRTDofBasis,
reffe::GenericRefFE{RaviartThomas},
phi::Field,
Jt_q,
sign_flip::AbstractVector{Bool})
nodes, nf_nodes, nf_moments, face_moments, Jt_q_cache = cache

Jt_q = evaluate!(Jt_q_cache,(phi),nodes)
nodes, nf_nodes, nf_moments, face_moments = cache
face_own_dofs=get_face_own_dofs(reffe)
for face in 1:length(face_moments)
moments = nf_moments[face]
if length(moments) > 0
nf_moments_face = nf_moments[face]
face_moments_face = face_moments[face]
if length(nf_moments_face) > 0
sign = (-1)^sign_flip[face_own_dofs[face][1]]
num_qpoints, num_moments = size(moments)
num_qpoints, num_moments = size(nf_moments_face)
for i in 1:num_qpoints
Jt_q_i = Jt_q[nf_nodes[face][i]]
change = sign * meas(Jt_q_i) * pinvJt(Jt_q_i)
for j in 1:num_moments
face_moments[face][i,j] = change nf_moments[face][i,j]
face_moments_face[i,j] = change nf_moments_face[i,j]
end
end
end
end
MomentBasedDofBasis(nodes,face_moments,nf_nodes)
end


# Support for DIV operator
function _DIV(f::LazyArray{<:Fill{T}}) where T
df=_DIV(f.args[1])
Expand Down

0 comments on commit a250431

Please sign in to comment.