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

Split values into geometric mapping and function values #764

Merged
merged 175 commits into from
Dec 6, 2023

Conversation

KnutAM
Copy link
Member

@KnutAM KnutAM commented Jul 12, 2023

This PR splits Values into geometric values and function values, the main advantages are

  • Reuse code for reinit! for FaceValues and CellValues
  • The new parameterization is free to allow e.g. GPU-arrays

For a more complete picture, please see #798.

For normal use, this is not breaking wrt. master (i.e. no direct field access).
However, CellValues{<:VectorInterpolation} cannot be dispatched on anymore. This can be added if required, but in that case perhaps (re-)introducing CellVectorValues and CellScalarValues as aliases is cleaner.

Description of internal changes

Description of internal changes

struct CellValues{FV, GV, QR, detT<:AbstractVector} <: AbstractCellValues
    fun_values::FV # FunctionValues
    geo_values::GV # GeometryValues
    qr::QR         # QuadratureRule
    detJdV::detT   # AbstractVector{<:Number}
end

struct FaceValues{FV, GV, QR, detT, nT, V_FV<:AbstractVector{FV}, V_GV<:AbstractVector{GV}} <: AbstractFaceValues
    fun_values::V_FV # AbstractVector{FunctionValues}
    geo_values::V_GV # AbstractVector{GeometryValues}
    qr::QR           # FaceQuadratureRule
    detJdV::detT     # AbstractVector{<:Number}
    normals::nT      # AbstractVector{<:Vec}
    current_face::ScalarWrapper{Int}
end

where

QuadratureRule/FaceQuadratureRule: No changes, contains the location and
weights for the quadrature on the reference element.

GeometryValues contain information about the geometric interpolation on
the reference element. These are used to calculate the required quantities
for (1) mapping mutable parts of FunctionValues from the reference element to the
actual element and (2) calculating properties for the current cell, such as the detJdV
for each quadrature point and the face normals.

struct GeometryValues{IP, M_t, dMdξ_t, d2Mdξ2_t}
    ip::IP             # ::Interpolation                Geometric interpolation 
    M::M_t             # ::AbstractVector{<:Number}     Values of geometric shape functions
    dMdξ::dMdξ_t       # ::AbstractVector{<:Vec}        Gradients of geometric shape functions in ref-domain
    d2Mdξ2::d2Mdξ2_t   # ::AbstractVector{<:Tensor{2}}  Hessians of geometric shape functions in ref-domain
                       # ::Nothing                      When hessians are not required
end

FunctionValues contain information about the function interpolation on
the reference element, as well as caches for the value and gradient on the
current element.

struct FunctionValues{IP, N_t, dNdx_t, dNdξ_t}
    ip::IP          # ::Interpolation
    N_x::N_t        # ::AbstractMatrix{Union{<:Tensor,<:Number}}
    N_ξ::N_t        # ::AbstractMatrix{Union{<:Tensor,<:Number}}
    dNdx::dNdx_t    # ::AbstractMatrix{Union{<:Tensor,<:StaticArray}}
    dNdξ::dNdξ_t    # ::AbstractMatrix{Union{<:Tensor,<:StaticArray}}
end

Construction

The function interpolation (type parameter to FunctionValues) decides how
it should be mapped, by defining get_mapping_type(::Interpolation), possible
options are currently IdentityMapping() (same as what currently is in Ferrite),
CovariantPiolaMapping() (H(curl)), and ContravariantPiolaMapping() (H(div)).

The type of mapping determines what information is required to be precalculated in
GeometryValues, specifically if the hessian, d²M/dξ² must be calculated or not.

reinit! (Calculation and application of mapping)

First, mapping::MappingValues is calculated for the current quadrature point,
q_point, by using the precalculated values for the geometric interpolation in
geometry_values and the cell's coordinates, x:

mapping = calculate_mapping(geometry_values, q_point, x)

This mapping object can be querried for information required to perform the mapping,
specifically the jacobian and the hessian. If the latter should be calculated is determined
by whether d²M/dξ² is calculated in GeometryValues (which in turn was decided on construction based on the needs of the mapping type requested by the function interpolation).

We use this mapping to first calculate the detJdV (+ normals in FaceValues)
and later in apply_mapping! to update the values in FunctionValues appropriately fulling the requested mapping type.

@termi-official termi-official self-requested a review July 12, 2023 10:35
@codecov-commenter
Copy link

codecov-commenter commented Jul 12, 2023

Codecov Report

Attention: 6 lines in your changes are missing coverage. Please review.

Comparison is base (2dfcaf9) 93.08% compared to head (a255614) 93.15%.

Files Patch % Lines
src/FEValues/face_values.jl 95.77% 3 Missing ⚠️
src/FEValues/cell_values.jl 96.66% 2 Missing ⚠️
src/FEValues/FunctionValues.jl 98.79% 1 Missing ⚠️

❗ Your organization needs to install the Codecov GitHub app to enable full functionality.

Additional details and impacted files
@@            Coverage Diff             @@
##           master     #764      +/-   ##
==========================================
+ Coverage   93.08%   93.15%   +0.06%     
==========================================
  Files          34       36       +2     
  Lines        5163     5229      +66     
==========================================
+ Hits         4806     4871      +65     
- Misses        357      358       +1     

☔ View full report in Codecov by Sentry.
📢 Have feedback on the report? Share it here.

@KnutAM
Copy link
Member Author

KnutAM commented Jul 25, 2023

Seems to fix #616 (At least the example provided by @termi-official is not erroring, but tests are needed)

@KnutAM KnutAM linked an issue Jul 25, 2023 that may be closed by this pull request
Copy link
Member

@termi-official termi-official left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Sorry for the late review! I think this looks really good already. I left some questions on the design.

src/FEValues/GeometryValues.jl Outdated Show resolved Hide resolved
docs/src/literate-tutorials/linear_shell.jl Show resolved Hide resolved
src/FEValues/CellValues.jl Outdated Show resolved Hide resolved
src/FEValues/CellValues.jl Outdated Show resolved Hide resolved
Copy link
Member

@fredrikekre fredrikekre left a comment

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

LGTM, just some stylistic nits that you don't have to fix if you don't agree. Can also merge and then do smaller clean ups later to make it easier to review.

Please squash-merge (and use a better commit message than the default message containing all commits in the branch).

src/FEValues/FunctionValues.jl Show resolved Hide resolved
src/FEValues/FunctionValues.jl Outdated Show resolved Hide resolved
src/FEValues/FunctionValues.jl Outdated Show resolved Hide resolved
src/FEValues/FunctionValues.jl Outdated Show resolved Hide resolved
src/FEValues/GeometryMapping.jl Outdated Show resolved Hide resolved
src/FEValues/cell_values.jl Outdated Show resolved Hide resolved
src/FEValues/cell_values.jl Outdated Show resolved Hide resolved
src/PointEval/point_values.jl Outdated Show resolved Hide resolved
@termi-official
Copy link
Member

termi-official commented Dec 6, 2023

I think @AbdAlazezAhmed just had a good point on slack. Should we also update the new InterfaceValues in this PR or in a separate one?

@KnutAM
Copy link
Member Author

KnutAM commented Dec 6, 2023

I think @AbdAlazezAhmed just had a good point on slack. Should we also update the new InterfaceValues in this PR or in a separate one?

Thanks - did it now with 2cbd2b9

@KnutAM KnutAM linked an issue Dec 6, 2023 that may be closed by this pull request
@KnutAM KnutAM merged commit 2d21665 into master Dec 6, 2023
7 checks passed
@KnutAM KnutAM deleted the kam/FunctionValues branch December 6, 2023 13:16
@AbdAlazezAhmed AbdAlazezAhmed mentioned this pull request Dec 6, 2023
@fredrikekre fredrikekre removed the awaiting review PR is finished from the authors POV, waiting for feedback label Dec 11, 2023
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
None yet
Projects
None yet
Development

Successfully merging this pull request may close these issues.

SimpleCellValues tutorial Bug in PointValues Field dim is not necessarily spatial dim
6 participants