Skip to content
GitLab
Explore
Sign in
Primary navigation
Search or go to…
Project
S
SemiSmoothNewton.jl
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Snippets
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Package registry
Container registry
Model registry
Operate
Environments
Terraform modules
Monitor
Incidents
Service Desk
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
GitLab community forum
Contribute to GitLab
Provide feedback
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
Stephan Hilb
SemiSmoothNewton.jl
Commits
20736d4b
Commit
20736d4b
authored
Aug 20, 2021
by
Stephan Hilb
Browse files
Options
Downloads
Patches
Plain Diff
implement l1-stable quasi-interpolation
parent
f31f6cb3
No related branches found
No related tags found
No related merge requests found
Changes
1
Show whitespace changes
Inline
Side-by-side
Showing
1 changed file
src/operator.jl
+74
-22
74 additions, 22 deletions
src/operator.jl
with
74 additions
and
22 deletions
src/operator.jl
+
74
−
22
View file @
20736d4b
...
...
@@ -3,7 +3,7 @@ using LinearAlgebra: det, dot
using
StaticArrays
:
SA
,
SArray
,
MArray
using
ForwardDiff
:
jacobian
export
Poisson
,
L2Projection
,
init_point!
,
assemble
,
project_img
export
Poisson
,
L2Projection
,
init_point!
,
assemble
,
project_img
,
project_img2
abstract type
Operator
end
...
...
@@ -119,22 +119,8 @@ end
project_img
(
space
::
FeSpace
,
img
)
=
(
u
=
FeFunction
(
space
);
project_img!
(
u
,
img
))
function
project_img!
(
u
::
FeFunction
,
img
)
d
=
2
# domain dimension
space
=
u
.
space
mesh
=
space
.
mesh
f
=
ImageFunction
(
mesh
,
img
)
opparams
=
(;
f
)
nrdims
=
prod
(
space
.
size
)
# number of element dofs (i.e. local dofs not counting range dimensions)
nldofs
=
ndofs
(
space
.
element
)
a
(
xloc
,
u
,
du
,
v
,
dv
;
f
)
=
dot
(
u
,
v
)
l
(
xloc
,
v
,
dv
;
f
)
=
dot
(
f
,
v
)
# composite midpoint quadrature on lagrange point lattice
function
quadrature
(
p
)
function
quadrature
_composite_lagrange_midpoint
(
p
)
d_
=
2
n
=
binomial
(
p
+
2
,
2
)
weights
=
Vector
{
Float64
}(
undef
,
n
)
...
...
@@ -152,11 +138,24 @@ function project_img!(u::FeFunction, img)
return
weights
,
points
end
function
project_img!
(
u
::
FeFunction
,
img
)
d
=
2
# domain dimension
space
=
u
.
space
mesh
=
space
.
mesh
f
=
ImageFunction
(
mesh
,
img
)
opparams
=
(;
f
)
nrdims
=
prod
(
space
.
size
)
# number of element dofs (i.e. local dofs not counting range dimensions)
nldofs
=
ndofs
(
space
.
element
)
a
(
xloc
,
u
,
du
,
v
,
dv
;
f
)
=
dot
(
u
,
v
)
l
(
xloc
,
v
,
dv
;
f
)
=
dot
(
f
,
v
)
I
=
Float64
[]
J
=
Float64
[]
V
=
Float64
[]
b
=
zeros
(
ndofs
(
space
))
gdof
=
LinearIndices
((
nrdims
,
ndofs
(
space
)))
# mesh cells
for
cell
in
cells
(
mesh
)
...
...
@@ -167,7 +166,7 @@ function project_img!(u::FeFunction, img)
intel
=
abs
(
det
(
delmap
))
p
=
ceil
(
Int
,
diam
(
mesh
,
cell
))
qw
,
qx
=
quadrature
(
p
)
qw
,
qx
=
quadrature
_composite_lagrange_midpoint
(
p
)
nqpts
=
length
(
qw
)
# number of quadrature points
qphi
=
zeros
(
nrdims
,
nrdims
,
nldofs
,
nqpts
)
...
...
@@ -219,3 +218,56 @@ function project_img!(u::FeFunction, img)
u
.
data
.=
A
\
b
return
u
end
project_img2
(
space
::
FeSpace
,
img
)
=
(
u
=
FeFunction
(
space
);
project_img2!
(
u
,
img
))
# L1-stable quasi-interpolation operator
# (https://arxiv.org/pdf/1505.06931.pdf)
#
# FIXME: currently only approximate quadrature by sampling on lagrange lattice
# based on element size. Exact evaluation is tricky to implement (lots of
# intersection handling).
function
project_img2!
(
u
::
FeFunction
,
img
)
u
.
space
.
element
==
P1
()
throw
(
ArgumentError
(
"element unsupported"
))
d
=
2
# domain dimension
space
=
u
.
space
mesh
=
space
.
mesh
f
=
ImageFunction
(
mesh
,
img
)
nrdims
=
prod
(
space
.
size
)
# number of element dofs (i.e. local dofs not counting range dimensions)
nldofs
=
ndofs
(
space
.
element
)
# count contributions to respective dof for subsequent averaging
gdofcount
=
zeros
(
Int
,
size
(
u
.
data
))
# mesh cells
for
cell
in
cells
(
mesh
)
bind!
(
f
,
cell
)
p
=
ceil
(
Int
,
diam
(
mesh
,
cell
))
# we actually only use the lagrange lattice points
qw
,
qx
=
quadrature_composite_lagrange_midpoint
(
p
)
nqpts
=
length
(
qw
)
# number of quadrature points
qphi
=
zeros
(
nldofs
,
nqpts
)
for
k
in
axes
(
qx
,
2
)
xhat
=
SVector
{
d
}(
view
(
qx
,
:
,
k
))
x
=
elmap
(
mesh
,
cell
)(
xhat
)
# size: nrdims
value
=
evaluate
(
f
,
xhat
)
# size: nldofs
phix
=
evaluate_basis
(
space
.
element
,
SVector
{
d
}(
view
(
qx
,
:
,
k
)))
# size: nrdims × nldofs
gdofs
=
u
.
space
.
dofmap
[
:
,
:
,
cell
]
u
.
data
[
gdofs
]
.+=
value
*
(
12
.*
phix
'
.-
3
)
gdofcount
[
gdofs
]
.+=
1
end
end
u
.
data
./=
gdofcount
return
u
end
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment