Skip to content

Introduce isconstrained bitvector in ConstraintHandler#1310

Merged
KnutAM merged 4 commits intomasterfrom
kam/fasterconstraints
Mar 24, 2026
Merged

Introduce isconstrained bitvector in ConstraintHandler#1310
KnutAM merged 4 commits intomasterfrom
kam/fasterconstraints

Conversation

@KnutAM
Copy link
Copy Markdown
Member

@KnutAM KnutAM commented Mar 21, 2026

This makes it easier and faster to check if a dof is constrained.
Previous check was haskey(ch.dof_mapping, dofnr), now just ch.isconstrained[dofnr].

This speeds up zero_out_rows! and makes creating free_dofs easier (and a tiny bit faster).

Benchmarks

# close!
4.963 ms (59 allocations: 25.78 MiB) # PR
5.075 ms (57 allocations: 25.41 MiB) # master

# apply!
293.966 ms (0 allocations: 0 bytes) # PR
1.333 s (0 allocations: 0 bytes)    # master
benchmark code
using Ferrite, BenchmarkTools
grid = generate_grid(Hexahedron, (100,100,100));
dh = close!(add!(DofHandler(grid), :u, Lagrange{RefHexahedron, 1}()^3));
@btime close!(ch) setup = (ch = add!(ConstraintHandler(dh), Dirichlet(:u, getfacetset(grid,"left"), Returns(zero(Vec{3}))))) evals = 1;

ch = close!(add!(ConstraintHandler(dh), Dirichlet(:u, getfacetset(grid,"left"), Returns(zero(Vec{3})))));
K = allocate_matrix(dh); 
f = zeros(ndofs(dh));
@btime apply!($K, $f, $ch);

Extra runs for close!:

# PR 
  5.007 ms (59 allocations: 25.78 MiB)
  5.025 ms (59 allocations: 25.78 MiB)
  5.052 ms (59 allocations: 25.78 MiB)
  
# master
  5.100 ms (57 allocations: 25.41 MiB)
  5.121 ms (57 allocations: 25.41 MiB)
  5.111 ms (57 allocations: 25.41 MiB)

Also used to make close simpler and slightly faster, and make it easier to check if dof is constrained
@codecov
Copy link
Copy Markdown

codecov Bot commented Mar 21, 2026

Codecov Report

✅ All modified and coverable lines are covered by tests.
✅ Project coverage is 94.31%. Comparing base (8eb6a07) to head (f244f12).
⚠️ Report is 1 commits behind head on master.

Additional details and impacted files
@@           Coverage Diff           @@
##           master    #1310   +/-   ##
=======================================
  Coverage   94.30%   94.31%           
=======================================
  Files          40       40           
  Lines        6763     6768    +5     
=======================================
+ Hits         6378     6383    +5     
  Misses        385      385           

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

🚀 New features to boost your workflow:
  • ❄️ Test Analytics: Detect flaky tests, report on failures, and find test suite problems.

@KristofferC
Copy link
Copy Markdown
Collaborator

There is also BitSet but I haven't used it much so I don't really know the trade off.

@KnutAM
Copy link
Copy Markdown
Member Author

KnutAM commented Mar 21, 2026

There is also BitSet but I haven't used it much so I don't really know the trade off.

  BitSet([itr])

  Construct a sorted set of Ints generated by the given iterable object, 
  or an empty set. Implemented as a bit string, and therefore designed 
  for dense integer sets. If the set will be sparse (for example, holding 
  a few very large integers), use Set instead.

I thought of it but discarded when reading this, since for the constrained dofs we have very sparse sets for large grids. But is there a better way to use this than my naive thinking of prescribed = BitSet(ch.prescribed_dofs) and checking haskey(prescribed, dofnr)?

@KristofferC
Copy link
Copy Markdown
Collaborator

Yeah, probably not worthwhile. I don't think it matters in practice though, this should be a lot smaller than many of the other data structures used already.

Comment thread src/Dofs/ConstraintHandler.jl
Comment thread src/Dofs/ConstraintHandler.jl Outdated
nzval = K.nzval
return @inbounds for i in eachindex(colval, nzval)
if haskey(ch.dofmapping, colval[i])
return @inbounds for (i, col) in pairs(colval)
Copy link
Copy Markdown
Collaborator

Choose a reason for hiding this comment

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

pairs here feels odd. Should it be enumerate?

Copy link
Copy Markdown
Member Author

Choose a reason for hiding this comment

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

My interpretation here is that colval and nzval have the same indices, hence it should be pairs. Of course, in practice for SparseMatrixCSR both colval and nzval are Vectors, so it doesn't matter.

Copy link
Copy Markdown
Collaborator

Choose a reason for hiding this comment

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

Ok, I think that is fair.

@inbounds for i in eachindex(rowval, nzval)
if haskey(ch.dofmapping, rowval[i])
nzval[i] = 0
@inbounds for (i, row) in pairs(rowval)
Copy link
Copy Markdown
Collaborator

Choose a reason for hiding this comment

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

Same w.r.t pairs

Copy link
Copy Markdown
Member Author

Choose a reason for hiding this comment

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

As above, but here the storage types in AbstractSparseMatrixCSC could be offset arrays, and then I assume that it should be pairs since we assume rowval and nzval have the same indexing?

Copy link
Copy Markdown
Collaborator

Choose a reason for hiding this comment

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

👍

Copy link
Copy Markdown
Collaborator

@KristofferC KristofferC left a comment

Choose a reason for hiding this comment

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

Should be good to go I think?

@KnutAM KnutAM merged commit 9134b54 into master Mar 24, 2026
17 of 18 checks passed
@KnutAM KnutAM deleted the kam/fasterconstraints branch March 24, 2026 16:58
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Projects

None yet

Development

Successfully merging this pull request may close these issues.

3 participants