Skip to content

Commit e27c6bb

Browse files
Merge pull request #806 from kishore-nori/AD-for-SkeletonIntegration
Added a dummy tag for ForwardDiff configs being constructed in `Gridap` `Autodiff.jl` to fix issue #805
2 parents 02d96ad + c7f12dc commit e27c6bb

File tree

3 files changed

+42
-21
lines changed

3 files changed

+42
-21
lines changed

NEWS.md

+3
Original file line numberDiff line numberDiff line change
@@ -11,6 +11,9 @@ and this project adheres to [Semantic Versioning](https://semver.org/spec/v2.0.0
1111
- Extended the `MultiField` functional gradient (with respect to degrees-of-freedom of `MultiFieldFEFunction`) to functionals involving Skeleton integration. The interface remains the same `gradient(f,xh)`. Since PR [#799](https://github.com/gridap/Gridap.jl/pull/799)
1212
- Added functionality to take jacobian of functional involving integration (`DomainContribution`) over Skeleton faces (obtained from testing bilinear form with the whole set of test `fe_basis`), with respect to the degrees-of-freedom of `FEFunction`. The interface remains the same - `jacobian(f,uh)`. Since PR [#803](https://github.com/gridap/Gridap.jl/pull/803)
1313

14+
### Changed
15+
- Added a dummy tag for ForwardDiff configs being constructed in Gridap at `src/Arrays/Autodiff.jl` to fix issue [#805](https://github.com/gridap/Gridap.jl/issues/805). Since PR [#806](https://github.com/gridap/Gridap.jl/pull/806)
16+
1417
### Fixed
1518
- Fixed the behavior of `gradient` for functionals involving operations of `CellFields` inside `mean` and `jump` terms of Skeleton Integration terms. Since PR [#800](https://github.com/gridap/Gridap.jl/pull/800)
1619
- Fixed the behavior `SkeletonCellFieldPair` at the Boundary integration terms. Since PR [#800](https://github.com/gridap/Gridap.jl/pull/800)

src/Arrays/Autodiff.jl

+31-15
Original file line numberDiff line numberDiff line change
@@ -1,15 +1,17 @@
11

22
function autodiff_array_gradient(a,i_to_x)
3-
i_to_cfg = lazy_map(ConfigMap(ForwardDiff.gradient),i_to_x)
4-
i_to_xdual = lazy_map(DualizeMap(ForwardDiff.gradient),i_to_x)
3+
dummy_forwarddiff_tag = ()->()
4+
i_to_cfg = lazy_map(ConfigMap(ForwardDiff.gradient,dummy_forwarddiff_tag),i_to_x)
5+
i_to_xdual = lazy_map(DualizeMap(ForwardDiff.gradient,dummy_forwarddiff_tag),i_to_x)
56
i_to_ydual = a(i_to_xdual)
67
i_to_result = lazy_map(AutoDiffMap(ForwardDiff.gradient),i_to_ydual,i_to_x,i_to_cfg)
78
i_to_result
89
end
910

1011
function autodiff_array_jacobian(a,i_to_x)
11-
i_to_cfg = lazy_map(ConfigMap(ForwardDiff.jacobian),i_to_x)
12-
i_to_xdual = lazy_map(DualizeMap(ForwardDiff.jacobian),i_to_x)
12+
dummy_forwarddiff_tag = ()->()
13+
i_to_cfg = lazy_map(ConfigMap(ForwardDiff.jacobian,dummy_forwarddiff_tag),i_to_x)
14+
i_to_xdual = lazy_map(DualizeMap(ForwardDiff.jacobian,dummy_forwarddiff_tag),i_to_x)
1315
i_to_ydual = a(i_to_xdual)
1416
i_to_result = lazy_map(AutoDiffMap(ForwardDiff.jacobian),i_to_ydual,i_to_x,i_to_cfg)
1517
i_to_result
@@ -21,19 +23,21 @@ function autodiff_array_hessian(a,i_to_x)
2123
end
2224

2325
function autodiff_array_gradient(a,i_to_x,j_to_i)
24-
i_to_xdual = lazy_map(DualizeMap(ForwardDiff.gradient),i_to_x)
26+
dummy_forwarddiff_tag = ()->()
27+
i_to_xdual = lazy_map(DualizeMap(ForwardDiff.gradient,dummy_forwarddiff_tag),i_to_x)
2528
j_to_ydual = a(i_to_xdual)
2629
j_to_x = lazy_map(Reindex(i_to_x),j_to_i)
27-
j_to_cfg = lazy_map(ConfigMap(ForwardDiff.gradient),j_to_x)
30+
j_to_cfg = lazy_map(ConfigMap(ForwardDiff.gradient,dummy_forwarddiff_tag),j_to_x)
2831
j_to_result = lazy_map(AutoDiffMap(ForwardDiff.gradient),j_to_ydual,j_to_x,j_to_cfg)
2932
j_to_result
3033
end
3134

3235
function autodiff_array_jacobian(a,i_to_x,j_to_i)
33-
i_to_xdual = lazy_map(DualizeMap(ForwardDiff.jacobian),i_to_x)
36+
dummy_forwarddiff_tag = ()->()
37+
i_to_xdual = lazy_map(DualizeMap(ForwardDiff.jacobian,dummy_forwarddiff_tag),i_to_x)
3438
j_to_ydual = a(i_to_xdual)
3539
j_to_x = lazy_map(Reindex(i_to_x),j_to_i)
36-
j_to_cfg = lazy_map(ConfigMap(ForwardDiff.jacobian),j_to_x)
40+
j_to_cfg = lazy_map(ConfigMap(ForwardDiff.jacobian,dummy_forwarddiff_tag),j_to_x)
3741
j_to_result = lazy_map(AutoDiffMap(ForwardDiff.jacobian),j_to_ydual,j_to_x,j_to_cfg)
3842
j_to_result
3943
end
@@ -43,31 +47,43 @@ function autodiff_array_hessian(a,i_to_x,i_to_j)
4347
autodiff_array_jacobian(agrad,i_to_x,i_to_j)
4448
end
4549

46-
struct ConfigMap{F} <: Map
47-
f::F
50+
struct ConfigMap{
51+
F <: Union{typeof(ForwardDiff.gradient),typeof(ForwardDiff.jacobian)},
52+
T <: Union{<:Function,Nothing}} <: Map
53+
54+
f::F # ForwardDiff operation
55+
tag::T # function for config tag name
4856
end
57+
# constructor with nothing as the tag, for backwards compatibility
58+
ConfigMap(f) = ConfigMap(f,nothing)
4959

5060
# TODO Prescribing long chunk size can lead to slow compilation times!
5161
function return_cache(k::ConfigMap{typeof(ForwardDiff.gradient)},x)
52-
cfg = ForwardDiff.GradientConfig(nothing,x,ForwardDiff.Chunk{length(x)}())
62+
cfg = ForwardDiff.GradientConfig(k.tag,x,ForwardDiff.Chunk{length(x)}())
5363
cfg
5464
end
5565

5666
function return_cache(k::ConfigMap{typeof(ForwardDiff.jacobian)},x)
57-
cfg = ForwardDiff.JacobianConfig(nothing,x,ForwardDiff.Chunk{length(x)}())
67+
cfg = ForwardDiff.JacobianConfig(k.tag,x,ForwardDiff.Chunk{length(x)}())
5868
cfg
5969
end
6070

6171
function evaluate!(cfg,k::ConfigMap,x)
6272
cfg
6373
end
6474

65-
struct DualizeMap{F} <: Map
66-
f::F
75+
struct DualizeMap{
76+
F <: Union{typeof(ForwardDiff.gradient),typeof(ForwardDiff.jacobian)},
77+
T <: Union{<:Function,Nothing}} <: Map
78+
79+
f::F # ForwardDiff operation
80+
tag::T # function for config tag name
6781
end
82+
# constructor with nothing as the tag, for backwards compatibility
83+
DualizeMap(f) = DualizeMap(f,nothing)
6884

6985
function return_cache(k::DualizeMap,x)
70-
return_cache(ConfigMap(k.f),x)
86+
return_cache(ConfigMap(k.f,k.tag),x)
7187
end
7288

7389
function evaluate!(cfg,k::DualizeMap,x)

src/FESpaces/FEAutodiff.jl

+8-6
Original file line numberDiff line numberDiff line change
@@ -166,20 +166,21 @@ end
166166
modules without any circular dependency
167167
=#
168168
function autodiff_array_gradient(a, i_to_x, j_to_i::SkeletonPair)
169-
i_to_xdual = lazy_map(DualizeMap(ForwardDiff.gradient),i_to_x)
169+
dummy_forwarddiff_tag = ()->()
170+
i_to_xdual = lazy_map(DualizeMap(ForwardDiff.gradient,dummy_forwarddiff_tag),i_to_x)
170171

171172
# dual output of both sides at once
172173
j_to_ydual_plus, j_to_ydual_minus = a(i_to_xdual)
173174

174175
# Work for plus side
175176
j_to_x_plus = lazy_map(Reindex(i_to_x),j_to_i.plus)
176-
j_to_cfg_plus = lazy_map(ConfigMap(ForwardDiff.gradient),j_to_x_plus)
177+
j_to_cfg_plus = lazy_map(ConfigMap(ForwardDiff.gradient,dummy_forwarddiff_tag),j_to_x_plus)
177178
j_to_result_plus = lazy_map(AutoDiffMap(ForwardDiff.gradient),
178179
j_to_ydual_plus,j_to_x_plus,j_to_cfg_plus)
179180

180181
# Work for minus side
181182
j_to_x_minus = lazy_map(Reindex(i_to_x),j_to_i.minus)
182-
j_to_cfg_minus = lazy_map(ConfigMap(ForwardDiff.gradient),j_to_x_minus)
183+
j_to_cfg_minus = lazy_map(ConfigMap(ForwardDiff.gradient,dummy_forwarddiff_tag),j_to_x_minus)
183184
j_to_result_minus = lazy_map(AutoDiffMap(ForwardDiff.gradient),
184185
j_to_ydual_minus,j_to_x_minus,j_to_cfg_minus)
185186

@@ -190,7 +191,8 @@ function autodiff_array_gradient(a, i_to_x, j_to_i::SkeletonPair)
190191
end
191192

192193
function autodiff_array_jacobian(a,i_to_x,j_to_i::SkeletonPair)
193-
i_to_xdual = lazy_map(DualizeMap(ForwardDiff.jacobian),i_to_x)
194+
dummy_forwarddiff_tag = ()->()
195+
i_to_xdual = lazy_map(DualizeMap(ForwardDiff.jacobian,dummy_forwarddiff_tag),i_to_x)
194196
j_to_ydual_plus, j_to_ydual_minus = a(i_to_xdual)
195197

196198
# Bilinear form when tested with test basis functions dv results in
@@ -205,13 +207,13 @@ function autodiff_array_jacobian(a,i_to_x,j_to_i::SkeletonPair)
205207

206208
# Work for plus side
207209
j_to_x_plus = lazy_map(Reindex(i_to_x),j_to_i.plus)
208-
j_to_cfg_plus = lazy_map(ConfigMap(ForwardDiff.jacobian),j_to_x_plus)
210+
j_to_cfg_plus = lazy_map(ConfigMap(ForwardDiff.jacobian,dummy_forwarddiff_tag),j_to_x_plus)
209211
j_to_result_plus_dense = lazy_map(AutoDiffMap(ForwardDiff.jacobian),
210212
j_to_ydual_plus_dense,j_to_x_plus,j_to_cfg_plus)
211213

212214
# Work for minus side
213215
j_to_x_minus = lazy_map(Reindex(i_to_x),j_to_i.minus)
214-
j_to_cfg_minus = lazy_map(ConfigMap(ForwardDiff.jacobian),j_to_x_minus)
216+
j_to_cfg_minus = lazy_map(ConfigMap(ForwardDiff.jacobian,dummy_forwarddiff_tag),j_to_x_minus)
215217
j_to_result_minus_dense = lazy_map(AutoDiffMap(ForwardDiff.jacobian),
216218
j_to_ydual_minus_dense,j_to_x_minus,j_to_cfg_minus)
217219

0 commit comments

Comments
 (0)