Skip to content

Commit

Permalink
correct deficiency in vector cache definition
Browse files Browse the repository at this point in the history
  • Loading branch information
PetrKryslUCSD committed Apr 22, 2023
1 parent eb2eba2 commit 40a5f76
Show file tree
Hide file tree
Showing 6 changed files with 47 additions and 34 deletions.
2 changes: 1 addition & 1 deletion Project.toml
Original file line number Diff line number Diff line change
@@ -1,7 +1,7 @@
name = "FinEtools"
uuid = "91bb5406-6c9a-523d-811d-0644c4229550"
authors = ["Petr Krysl <[email protected]>"]
version = "6.0.18"
version = "6.1.0"

[deps]
ChunkSplitters = "ae650224-84b6-46f8-82ea-d812ca08434e"
Expand Down
2 changes: 1 addition & 1 deletion src/FEMMBaseModule.jl
Original file line number Diff line number Diff line change
Expand Up @@ -539,7 +539,7 @@ function distribloads(
startassembly!(assembler, P.nfreedofs)
for i in 1:count(fes) # Loop over elements
gathervalues_asmat!(geom, ecoords, fes.conn[i])
fill!(Fe, 0.0)
fill!(Fe, T(0.0))
for j in 1:npts
locjac!(loc, J, ecoords, Ns[j], gradNparams[j])
Jac = Jacobianmdim(self.integdomain, J, loc, fes.conn[i], Ns[j], m)
Expand Down
22 changes: 11 additions & 11 deletions src/ForceIntensityModule.jl
Original file line number Diff line number Diff line change
Expand Up @@ -31,8 +31,8 @@ getforce!(forceout::Vector{T}, XYZ::Matrix{T}, tangents::Matrix{T}, fe_label) wh
The buffer `forceout` is filled with the value of the force. The vector
`forceout` is also returned for convenience.
"""
struct ForceIntensity{T<:Number, F<:Function}
_cache::VectorCache{T, F} # vector cache where the current value of the force can be retrieved
struct ForceIntensity{CT<:Number, T<:Number, F<:Function}
_cache::VectorCache{CT, T, F} # vector cache where the current value of the force can be retrieved
end

"""
Expand Down Expand Up @@ -63,12 +63,12 @@ location `XYZ`, using if appropriate the information supplied in the Jacobian
matrix `tangents`, and the label of the finite element, `fe_label`.
"""
function ForceIntensity(
::Type{T},
::Type{CT},
ndofn,
computeforce!::F,
) where {T<:Number, F<:Function}
) where {CT<:Number, F<:Function}
# Allocate the buffer to be ready for the first call
return ForceIntensity(VectorCache(T, ndofn, computeforce!))
return ForceIntensity(VectorCache(CT, ndofn, computeforce!))
end

"""
Expand Down Expand Up @@ -120,21 +120,21 @@ v = updateforce!(fi, XYZ, tangents, fe_label)
```
"""
function ForceIntensity(
::Type{T},
::Type{CT},
ndofn,
computeforce!::F,
time::T,
) where {T<:Number, F<:Function}
) where {CT<:Number, F<:Function, T}
# Allocate the buffer to be ready for the first call
return ForceIntensity(VectorCache(T, ndofn, computeforce!, time))
return ForceIntensity(VectorCache(CT, ndofn, computeforce!, time))
end

"""
ForceIntensity(force::Vector{T}) where {T<:Number}
Construct force intensity when the constant `force` vector is given.
"""
function ForceIntensity(force::Vector{T}) where {T<:Number}
function ForceIntensity(force::Vector{CT}) where {CT<:Number}
return ForceIntensity(VectorCache(deepcopy(force)))
end

Expand All @@ -145,8 +145,8 @@ Construct force intensity when the force is given as a scalar value.
The dimension of the force vector in this case is 1.
"""
function ForceIntensity(force::T) where {T<:Number}
return ForceIntensity(T[force])
function ForceIntensity(force::CT) where {CT<:Number}
return ForceIntensity(CT[force])
end

"""
Expand Down
4 changes: 2 additions & 2 deletions src/SurfaceNormalModule.jl
Original file line number Diff line number Diff line change
Expand Up @@ -29,9 +29,9 @@ computenormal!(normalout::Vector{T}, XYZ::Matrix{T}, tangents::Matrix{T}, fe_lab
The buffer `normalout` is filled with the value of the normal vector.
"""
struct SurfaceNormal{T<:Number, F<:Function}
struct SurfaceNormal{CT<:Number, T<:Number, F<:Function}
# Cache of the current value of the normal
_cache::VectorCache{T, F}
_cache::VectorCache{CT, T, F}
end

"""
Expand Down
40 changes: 21 additions & 19 deletions src/VectorCacheModule.jl
Original file line number Diff line number Diff line change
Expand Up @@ -27,24 +27,24 @@ fillcache!(cacheout::Vector{T},
```
The cache `cacheout` is filled with the value of the vector.
"""
struct VectorCache{T<:Number, F<:Function}
struct VectorCache{CT<:Number, T<:Number, F<:Function}
# Function to update and retrieve the vector
_fillcache!::F
# Cache where the current value of the vector can be retrieved
_cache::Vector{T}
_cache::Vector{CT}
# Current time (or current load factor). Do not set directly. Use `settime!`
_time::Ref{T}
end

"""
VectorCache(::Type{T}, nentries::IT, fillcache!::F) where {T<:Number,F<:Function, IT}
VectorCache(::Type{CT}, nentries::IT, fillcache!::F) where {CT<:Number, F<:Function, IT}
Construct vector cache. The function to fill the vector cache is given.
This constructor is intended for *time-independent* vector caches.
This function needs to have a signature of
```
fillcache!(cacheout::Vector{T},
fillcache!(cacheout::Vector{CT},
XYZ::Matrix{T}, tangents::Matrix{T},
fe_label; time::T = 0.0) where {T}
Calculate the vector and copy it into the cache....
Expand All @@ -56,36 +56,37 @@ the location `XYZ`, using the information supplied in the Jacobian
matrix `tangents`, and the label of the finite element, `fe_label`, if
appropriate.
"""
function VectorCache(::Type{T}, nentries::IT, fillcache!::F) where {T<:Number,F<:Function, IT}
function VectorCache(::Type{CT}, nentries::IT, fillcache!::F) where {CT<:Number, F<:Function, IT}
T = real(CT(0.0))
function fillcachenotime!(
cacheout::Vector{T},
cacheout::Vector{CT},
XYZ::Matrix{T},
tangents::Matrix{T},
fe_label;
time::T = 0.0,
) where {T}
) where {CT, T}
return fillcache!(cacheout, XYZ, tangents, fe_label)
end
# Allocate the cache to be ready for the first call
return VectorCache(fillcachenotime!, zeros(T, nentries), Ref(0.0))
return VectorCache(fillcachenotime!, zeros(CT, nentries), Ref(0.0))
end

"""
VectorCache(
::Type{T},
::Type{CT},
nentries::IT,
fillcache!::F,
time::T,
) where {T<:Number, F<:Function, IT}
) where {CT<:Number, F<:Function, IT}
Construct vector cache. The function to fill the vector cache is given.
This constructor is intended for *time-dependent* vector caches.
This function needs to have a signature of
```
fillcache!(cacheout::Vector{T},
fillcache!(cacheout::Vector{CT},
XYZ::Matrix{T}, tangents::Matrix{T},
fe_label; time::T = 0.0) where {T}
fe_label; time::T = 0.0) where {CT, T}
Calculate the vector and copy it into the cache....
return forceout
end
Expand All @@ -96,28 +97,29 @@ matrix `tangents`, and the label of the finite element, `fe_label`, if
appropriate. The time can also be supplied (keyword argument `time`).
"""
function VectorCache(
::Type{T},
::Type{CT},
nentries::IT,
fillcache!::F,
time::T,
) where {T<:Number, F<:Function, IT}
) where {CT<:Number, IT, F<:Function, T}
# Allocate the cache to be ready for the first call
return VectorCache(fillcache!, zeros(T, nentries), Ref(time))
return VectorCache(fillcache!, zeros(CT, nentries), Ref(time))
end

"""
VectorCache(vector::Vector{T}) where {T<:Number}
VectorCache(vector::Vector{CT}) where {CT<:Number}
Construct vector cache. The *constant* vector is given.
"""
function VectorCache(vector::Vector{T}) where {T<:Number}
function VectorCache(vector::Vector{CT}) where {CT<:Number}
T = real(CT(0.0))
function fillcache!(
cacheout::Vector{T},
cacheout::Vector{CT},
XYZ::Matrix{T},
tangents::Matrix{T},
fe_label;
time::T = 0.0,
)
) where {CT, T}
# do nothing: the vector is already in the cache
return cacheout
end
Expand Down
11 changes: 11 additions & 0 deletions test/test_basics.jl
Original file line number Diff line number Diff line change
Expand Up @@ -939,6 +939,17 @@ function test()
fi = ForceIntensity(11.0)
F = distribloads(femm, geom, psi, fi, 3)
@test abs(sum(F) - (*).(L, W, t, 11.0)) / 667 <= 1.0e-5

psi = NodalField(fill(1.0 + 1.0im, count(fens), 1))
numberdofs!(psi)

femm = FEMMBase(IntegDomain(fes, GaussRule(3, 2)))
fi = ForceIntensity([11.0im])
F = distribloads(femm, geom, psi, fi, 3)
@test abs(sum(F) - (*).(L, W, t, 11.0im)) / 667 <= 1.0e-5
fi = ForceIntensity(11.0im)
F = distribloads(femm, geom, psi, fi, 3)
@test abs(sum(F) - (*).(L, W, t, 11.0im)) / 667 <= 1.0e-5
true
end
end
Expand Down

2 comments on commit 40a5f76

@PetrKryslUCSD
Copy link
Owner Author

Choose a reason for hiding this comment

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

@JuliaRegistrator
Copy link

Choose a reason for hiding this comment

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

Registration pull request created: JuliaRegistries/General/82116

After the above pull request is merged, it is recommended that a tag is created on this repository for the registered package version.

This will be done automatically if the Julia TagBot GitHub Action is installed, or can be done manually through the github interface, or via:

git tag -a v6.1.0 -m "<description of version>" 40a5f76e5947feee0c670314b874c9e5593ec315
git push origin v6.1.0

Please sign in to comment.