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

Final touches phys space #218

Merged
merged 3 commits into from
Mar 23, 2020
Merged
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
2 changes: 1 addition & 1 deletion Project.toml
Original file line number Diff line number Diff line change
@@ -1,7 +1,7 @@
name = "Gridap"
uuid = "56d4f2e9-7ea1-5844-9cf6-b9c51ca7ce8e"
authors = ["Santiago Badia <[email protected]>", "Francesc Verdugo <[email protected]>"]
version = "0.8.1"
version = "0.9.0"

[deps]
BSON = "fbb218c0-5317-5bc6-957e-2ee96dd4b1f0"
Expand Down
14 changes: 9 additions & 5 deletions src/FESpaces/ConformingFESpaces.jl
Original file line number Diff line number Diff line change
Expand Up @@ -71,11 +71,7 @@ function GradConformingFESpace(
cell_to_ctype = get_cell_type(grid_topology)
cell_map = get_cell_map(grid)

if is_ref
cell_shapefuns, cell_dof_basis = compute_cell_space(reffes, cell_to_ctype, cell_map)
else
cell_shapefuns, cell_dof_basis = compute_cell_space_physical(reffes, cell_to_ctype, cell_map)
end
cell_shapefuns, cell_dof_basis = compute_cell_space(reffes, cell_to_ctype, cell_map, Val(is_ref))

UnconstrainedFESpace(
nfree,
Expand All @@ -90,6 +86,14 @@ function GradConformingFESpace(

end

function compute_cell_space(reffes, cell_to_ctype, cell_map,ref_style::Val{true})
compute_cell_space(reffes, cell_to_ctype, cell_map)
end

function compute_cell_space(reffes, cell_to_ctype, cell_map,ref_style::Val{false})
compute_cell_space_physical(reffes, cell_to_ctype, cell_map)
end

"""
"""
function compute_cell_space(reffes, cell_to_ctype, cell_map)
Expand Down
7 changes: 1 addition & 6 deletions src/FESpaces/CurlConformingFESpaces.jl
Original file line number Diff line number Diff line change
Expand Up @@ -25,12 +25,7 @@ function CurlConformingFESpace(
cell_to_ctype = get_cell_type(grid_topology)
cell_map = get_cell_map(grid)

# cell_shapefuns, cell_dof_basis = _compute_hcurl_cell_space(reffes, cell_to_ctype, cell_map)
if is_ref
cell_shapefuns, cell_dof_basis = compute_cell_space(reffes,cell_to_ctype,cell_map)
else
cell_shapefuns, cell_dof_basis = compute_cell_space_physical(reffes,cell_to_ctype,cell_map)
end
cell_shapefuns, cell_dof_basis = compute_cell_space(reffes,cell_to_ctype,cell_map,Val(is_ref))

UnconstrainedFESpace(
nfree,
Expand Down
6 changes: 1 addition & 5 deletions src/FESpaces/DiscontinuousFESpaces.jl
Original file line number Diff line number Diff line change
Expand Up @@ -14,11 +14,7 @@ function DiscontinuousFESpace(reffes::Vector{<:ReferenceFE}, trian::Triangulatio
dirichlet_cells = Int[]
ntags = 0

if is_ref
cell_shapefuns, cell_dof_basis = compute_cell_space(reffes, cell_to_ctype, cell_map)
else
cell_shapefuns, cell_dof_basis = compute_cell_space_physical(reffes, cell_to_ctype, cell_map)
end
cell_shapefuns, cell_dof_basis = compute_cell_space(reffes, cell_to_ctype, cell_map,Val(is_ref))

UnconstrainedFESpace(
nfree,
Expand Down
8 changes: 1 addition & 7 deletions src/FESpaces/DivConformingFESpaces.jl
Original file line number Diff line number Diff line change
Expand Up @@ -24,13 +24,7 @@ function DivConformingFESpace(
cell_to_ctype = get_cell_type(grid_topology)
cell_map = get_cell_map(grid)

# cell_shapefuns, cell_dof_basis = _compute_hdiv_cell_space(reffes, cell_to_ctype, cell_map)
if is_ref
cell_shapefuns, cell_dof_basis = compute_cell_space(reffes,cell_to_ctype,cell_map)
else
cell_shapefuns, cell_dof_basis = compute_cell_space_physical(reffes,cell_to_ctype,cell_map)
end

cell_shapefuns, cell_dof_basis = compute_cell_space(reffes,cell_to_ctype,cell_map,Val(is_ref))

UnconstrainedFESpace(
nfree,
Expand Down
4 changes: 2 additions & 2 deletions src/Geometry/CellFields.jl
Original file line number Diff line number Diff line change
Expand Up @@ -26,8 +26,8 @@ physical space
RefStyle(::Type{<:CellFieldLike}) = @notimplemented

# We use duck typing here for all types marked with the RefStyle
# @santiagobadia : That is dangerous, it creates overflows if called with
RefStyle(::T) where T = RefStyle(T)
RefStyle(::Type) = @notimplemented
RefStyle(a) = RefStyle(typeof(a))
is_in_ref_space(::Type{T}) where T = get_val_parameter(RefStyle(T))
is_in_ref_space(::T) where T = is_in_ref_space(T)
is_in_physical_space(::Type{T}) where T = !is_in_ref_space(T)
Expand Down
1 change: 1 addition & 0 deletions src/MultiField/MultiField.jl
Original file line number Diff line number Diff line change
Expand Up @@ -37,6 +37,7 @@ using Gridap.FESpaces: CellMatKernel
using Gridap.FESpaces: CellVecKernel
import Gridap.Geometry: get_cell_map
import Gridap.Geometry: similar_object
import Gridap.Geometry: change_ref_style
import Gridap.Geometry: restrict
import Gridap.Fields: integrate
import Gridap.Fields: evaluate
Expand Down
12 changes: 12 additions & 0 deletions src/MultiField/MultiFieldCellBases.jl
Original file line number Diff line number Diff line change
Expand Up @@ -21,6 +21,12 @@ TrialStyle(::Type{CellBasisWithFieldID{S,R}}) where {S,R} = Val{S}()

RefStyle(::Type{CellBasisWithFieldID{S,R}}) where {S,R} = Val{R}()

function change_ref_style(cf::CellBasisWithFieldID)
cell_basis = change_ref_style(cf.cell_basis)
field_id = cf.field_id
CellBasisWithFieldID(cell_basis,field_id)
end

function similar_object(cb::CellBasisWithFieldID,array::AbstractArray)
cell_basis = similar_object(cb.cell_basis,array)
field_id = cb.field_id
Expand Down Expand Up @@ -99,6 +105,12 @@ function _have_same_field_ids(a::CellMatrixFieldWithFieldIds,b::CellMatrixFieldW
(a.field_id_rows == b.field_id_rows) && (a.field_id_cols == b.field_id_cols)
end

function change_ref_style(cf::CellMatrixFieldWithFieldIds)
cell_matrix_field = change_ref_style(cf.cell_matrix_field)
CellMatrixFieldWithFieldIds(
cell_matrix_field,cf.field_id_rows,cf.field_id_cols,RefStyle(cell_matrix_field))
end

function similar_object(cf::CellMatrixFieldWithFieldIds,a::AbstractArray)
cell_matrix_field = similar_object(cf.cell_matrix_field,a)
field_id_rows = cf.field_id_rows
Expand Down
1 change: 0 additions & 1 deletion test/FESpacesTests/CellDofBasesTests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -28,7 +28,6 @@ cell_basis = get_cell_basis(V)

test_cell_dof_basis(cell_dof_basis,cell_basis)

# TODO not working: it should be the cell-wise identity...
#display(evaluate(cell_dof_basis,cell_basis))

V = TestFESpace(
Expand Down
4 changes: 4 additions & 0 deletions test/MultiFieldTests/MultiFieldCellBasesTests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -29,6 +29,7 @@ cell_basis_v = get_cell_basis(V)
field_id_v = 3
v = CellBasisWithFieldID(cell_basis_v,field_id_v)
vq = collect(evaluate(v,q))
test_cell_basis(v,q,vq)
@test is_in_ref_space(v)

field_id_a = 2
Expand All @@ -48,6 +49,7 @@ cell_basis_vp = get_cell_basis(Vp)
field_id_vp = 6
vp = CellBasisWithFieldID(cell_basis_vp,field_id_vp)
vpq = collect(evaluate(vp,q))
test_cell_basis(vp,q,vpq)
@test !is_in_ref_space(vp)

r = 2*v
Expand All @@ -74,11 +76,13 @@ r = ∇(u)
@test is_trial(r)

r = v * u
test_cell_matrix_field(r,q,collect(evaluate(r,q)))
@test isa(r,CellMatrixFieldWithFieldIds)
@test r.field_id_rows == v.field_id
@test r.field_id_cols == u.field_id

r = u * v
test_cell_matrix_field(r,q,collect(evaluate(r,q)))
@test isa(r,CellMatrixFieldWithFieldIds)
@test r.field_id_rows == v.field_id
@test r.field_id_cols == u.field_id
Expand Down