Skip to content

Commit

Permalink
Merge pull request #728 from gridap/passing_cellfe_optargs
Browse files Browse the repository at this point in the history
CellFE constructor now gets optional arguments and pass them down
  • Loading branch information
amartinhuertas authored Dec 3, 2021
2 parents 7f5f15e + e8c1cd1 commit 8984d73
Show file tree
Hide file tree
Showing 3 changed files with 12 additions and 9 deletions.
7 changes: 4 additions & 3 deletions src/FESpaces/ConformingFESpaces.jl
Original file line number Diff line number Diff line change
Expand Up @@ -94,14 +94,15 @@ Generate a CellFE from a vector of reference fes
function CellFE(
model::DiscreteModel,
cell_reffe::AbstractArray{<:ReferenceFE},
conformity::Conformity
conformity::Conformity,
args...
)
cell_conformity = CellConformity(cell_reffe,conformity)
ctype_reffe, cell_ctype = compress_cell_data(cell_reffe)
ctype_num_dofs = map(num_dofs,ctype_reffe)
ctype_ldof_comp = map(reffe->get_dof_to_comp(reffe),ctype_reffe)
cell_shapefuns = get_cell_shapefuns(model,cell_reffe,conformity)
cell_dof_basis = get_cell_dof_basis(model,cell_reffe,conformity)
cell_shapefuns = get_cell_shapefuns(model,cell_reffe,conformity,args...)
cell_dof_basis = get_cell_dof_basis(model,cell_reffe,conformity,args...)
cell_shapefuns_domain = ReferenceDomain()
cell_dof_basis_domain = cell_shapefuns_domain
max_order = maximum(map(get_order,ctype_reffe))
Expand Down
8 changes: 4 additions & 4 deletions src/FESpaces/DivConformingFESpaces.jl
Original file line number Diff line number Diff line change
Expand Up @@ -28,8 +28,8 @@ struct TransformRTDofBasis{Dc,Dp} <: Map end ;

function get_cell_dof_basis(model::DiscreteModel,
cell_reffe::AbstractArray{<:GenericRefFE{RaviartThomas}},
::DivConformity)
sign_flip = get_sign_flip(model, cell_reffe)
::DivConformity,
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 @@ -46,8 +46,8 @@ end

function get_cell_shapefuns(model::DiscreteModel,
cell_reffe::AbstractArray{<:GenericRefFE{RaviartThomas}},
::DivConformity)
sign_flip = get_sign_flip(model, cell_reffe)
::DivConformity,
sign_flip=get_sign_flip(model, cell_reffe))
cell_reffe_shapefuns=lazy_map(get_shapefuns,cell_reffe)
k=ContraVariantPiolaMap()
lazy_map(k,
Expand Down
6 changes: 4 additions & 2 deletions src/FESpaces/FESpaceFactories.jl
Original file line number Diff line number Diff line change
Expand Up @@ -119,13 +119,15 @@ function FESpace(
return V
end

function FESpace(model::DiscreteModel, reffe::Tuple{<:ReferenceFEName,Any,Any}; kwargs...)
function FESpace(model::DiscreteModel,
reffe::Tuple{<:ReferenceFEName,Any,Any}; kwargs...)
basis, reffe_args,reffe_kwargs = reffe
cell_reffe = ReferenceFE(model,basis,reffe_args...;reffe_kwargs...)
FESpace(model,cell_reffe;kwargs...)
end

function FESpace(model::DiscreteModel, reffe::ReferenceFE; kwargs...)
function FESpace(model::DiscreteModel,
reffe::ReferenceFE; kwargs...)
cell_reffe = Fill(reffe,num_cells(model))
FESpace(model,cell_reffe;kwargs...)
end
Expand Down

0 comments on commit 8984d73

Please sign in to comment.