Closed tamiko closed 1 year ago
Stacktrace: (matrix_free/categorize_by_boundary_ids_01
)
#0 dealii::TriaAccessor<1, 2, 2>::boundary_id() const (this=0x7fffffffc3c0) at /home/tamiko/workspace/dealii/include/deal.II/grid/tria_accessor.templates.h:2005
#1 dealii::TriaAccessor<1, 2, 2>::at_boundary() const (this=0x7fffffffc3c0) at /home/tamiko/workspace/dealii/include/deal.II/grid/tria_accessor.templates.h:2080
#2 dealii::MatrixFreeTools::categorize_by_boundary_ids<2, dealii::MatrixFree<2, double, dealii::VectorizedArray<double, 2ul> >::AdditionalData>(dealii::Triangulation<2, 2> const&, dealii::MatrixFree<2, double, dealii::VectorizedArray<double, 2ul> >::AdditionalData&)::{lambda(auto:1 const&)#1}::operator()<dealii::TriaActiveIterator<dealii::CellAccessor<2, 2> > >(dealii::TriaActiveIterator<dealii::CellAccessor<2, 2> > const&) const (cell=<optimized out>, __closure=<optimized out>)
at /home/tamiko/workspace/dealii/include/deal.II/matrix_free/tools.h:382
#3 dealii::MatrixFreeTools::categorize_by_boundary_ids<2, dealii::MatrixFree<2, double, dealii::VectorizedArray<double, 2ul> >::AdditionalData>(dealii::Triangulation<2, 2> const&, dealii::MatrixFree<2, double, dealii::VectorizedArray<double, 2ul> >::AdditionalData&) (tria=..., additional_data=...)
at /home/tamiko/workspace/dealii/include/deal.II/matrix_free/tools.h:399
#4 0x00005555555ba6d3 in test<2>(unsigned int) (n_refinements=n_refinements@entry=3)
at /home/tamiko/workspace/dealii/tests/matrix_free/categorize_by_boundary_ids_01.cc:105
#5 0x00005555555ba48c in main() () at /home/tamiko/workspace/dealii/tests/matrix_free/categorize_by_boundary_ids_01.cc:119
Stacktrace: (matrix_free/pbc_orientation_02
)
#0 dealii::TriaAccessor<2, 3, 3>::boundary_id() const (this=0x7fffffffbc00) at /home/tamiko/workspace/dealii/include/deal.II/grid/tria_accessor.templates.h:2005
#1 dealii::TriaAccessor<2, 3, 3>::at_boundary() const (this=0x7fffffffbc00) at /home/tamiko/workspace/dealii/include/deal.II/grid/tria_accessor.templates.h:2080
#2 dealii::MatrixFreeTools::categorize_by_boundary_ids<3, dealii::MatrixFree<3, double, dealii::VectorizedArray<double, 2ul> >::AdditionalData>(dealii::Triangulation<3, 3> const&, dealii::MatrixFree<3, double, dealii::VectorizedArray<double, 2ul> >::AdditionalData&)::{lambda(auto:1 const&)#1}::operator()<dealii::TriaActiveIterator<dealii::CellAccessor<3, 3> > >(dealii::TriaActiveIterator<dealii::CellAccessor<3, 3> > const&) const (cell=<optimized out>, __closure=<optimized out>)
at /home/tamiko/workspace/dealii/include/deal.II/matrix_free/tools.h:382
#3 dealii::MatrixFreeTools::categorize_by_boundary_ids<3, dealii::MatrixFree<3, double, dealii::VectorizedArray<double, 2ul> >::AdditionalData>(dealii::Triangulation<3, 3> const&, dealii::MatrixFree<3, double, dealii::VectorizedArray<double, 2ul> >::AdditionalData&) (tria=..., additional_data=...)
at /home/tamiko/workspace/dealii/include/deal.II/matrix_free/tools.h:399
#4 0x00005555555cbc80 in test<3, 3>() () at /home/tamiko/workspace/dealii/include/deal.II/base/smartpointer.h:369
#5 0x00005555555cb539 in main() () at /home/tamiko/workspace/dealii/tests/matrix_free/pbc_orientation_02.cc:217
This looks odd:
I don't remember why I derefernced the iterator (and should not assigned to a reference); but one does not need the accessor but can work directly with the iterator...
@peterrum Would you mind to open a PR that replaces the expression by the iterator using face->boundary_id()
? :smiley: Thanks for looking into this!
I agree that it looks funny, but I don't think it's wrong -- binding a temporary to a reference is allowed. But it's worth fixing anyway. Give it a try!
@tamiko Do any of these tests fail in debug mode?
@bangerth No, they do not fail in debug mode.
I think the mistake might be that a non-const lvalue cannot bind to a temporary.
Oh, the const
! You may be right!
But this should be a compiler error... I am confused.
But here we go:
=================================================================
==2034964==ERROR: AddressSanitizer: stack-use-after-scope on address 0x7f8f5d74a474 at pc 0x55b47321dcbc bp 0x7fffcf5e20f0 sp 0x7fffcf5e20e8
READ of size 4 at 0x7f8f5d74a474 thread T0
#0 0x55b47321dcbb in dealii::TriaAccessorBase<1, 2, 2>::state() const /srv/testsuite/dealii/include/deal.II/grid/tria_accessor.templates.h:261:32
#1 0x55b47321dcbb in dealii::TriaAccessor<1, 2, 2>::used() const /srv/testsuite/dealii/include/deal.II/grid/tria_accessor.templates.h:1241:3
#2 0x55b47321dcbb in dealii::TriaAccessor<1, 2, 2>::boundary_id() const /srv/testsuite/dealii/include/deal.II/grid/tria_accessor.templates.h:2003:3
#3 0x55b47321680b in dealii::TriaAccessor<1, 2, 2>::at_boundary() const /srv/testsuite/dealii/include/deal.II/grid/tria_accessor.templates.h:2080:11
#4 0x55b47321680b in auto void dealii::MatrixFreeTools::categorize_by_boundary_ids<2, dealii::MatrixFree<2, double, dealii::VectorizedArray<double, 2ul>>::AdditionalData
>(dealii::Triangulation<2, 2> const&, dealii::MatrixFree<2, double, dealii::VectorizedArray<double, 2ul>>::AdditionalData&)::'lambda'(auto const&)::operator()<dealii::TriaAc
tiveIterator<dealii::CellAccessor<2, 2>>>(auto const&) const /srv/testsuite/dealii/include/deal.II/matrix_free/tools.h:382:20
#5 0x55b473211d5a in void dealii::MatrixFreeTools::categorize_by_boundary_ids<2, dealii::MatrixFree<2, double, dealii::VectorizedArray<double, 2ul>>::AdditionalData>(dea
lii::Triangulation<2, 2> const&, dealii::MatrixFree<2, double, dealii::VectorizedArray<double, 2ul>>::AdditionalData&) /srv/testsuite/dealii/include/deal.II/matrix_free/tool
s.h:399:17
#6 0x55b4732099ae in void test<2>(unsigned int) /srv/testsuite/dealii/tests/matrix_free/categorize_by_boundary_ids_01.cc:105:3
#7 0x55b473207cb4 in main /srv/testsuite/dealii/tests/matrix_free/categorize_by_boundary_ids_01.cc:119:5
#8 0x7f8f6a450989 (/usr/lib64/libc.so.6+0x23989)
#9 0x7f8f6a450a44 in __libc_start_main (/usr/lib64/libc.so.6+0x23a44)
#10 0x55b4730f0d50 in _start (/srv/temp/build-clang/tests/matrix_free/categorize_by_boundary_ids_01.debug/categorize_by_boundary_ids_01.debug+0xa1d50)
Address 0x7f8f5d74a474 is located in stack of thread T0 at offset 1140 in frame
#0 0x55b4732164bf in auto void dealii::MatrixFreeTools::categorize_by_boundary_ids<2, dealii::MatrixFree<2, double, dealii::VectorizedArray<double, 2ul>>::AdditionalData
>(dealii::Triangulation<2, 2> const&, dealii::MatrixFree<2, double, dealii::VectorizedArray<double, 2ul>>::AdditionalData&)::'lambda'(auto const&)::operator()<dealii::TriaAc
tiveIterator<dealii::CellAccessor<2, 2>>>(auto const&) const /srv/testsuite/dealii/include/deal.II/matrix_free/tools.h:377
This frame has 4 object(s):
[32, 328) 'agg.tmp.i.i36'
[400, 696) 'agg.tmp.i'
[768, 1064) 'agg.tmp.i.i'
[1136, 1152) 'ref.tmp' (line 381) <== Memory access at offset 1140 is inside this variable
HINT: this may be a false positive if your program uses some custom stack unwind mechanism, swapcontext or vfork
(longjmp and C++ exceptions *are* supported)
SUMMARY: AddressSanitizer: stack-use-after-scope /srv/testsuite/dealii/include/deal.II/grid/tria_accessor.templates.h:261:32 in dealii::TriaAccessorBase<1, 2, 2>::state() co
nst
Shadow bytes around the buggy address:
0x7f8f5d74a180: f2 f2 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8
0x7f8f5d74a200: f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8
0x7f8f5d74a280: f8 f8 f8 f8 f8 f8 f8 f2 f2 f2 f2 f2 f2 f2 f2 f2
0x7f8f5d74a300: f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8
0x7f8f5d74a380: f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8 f8
=>0x7f8f5d74a400: f8 f8 f8 f8 f8 f2 f2 f2 f2 f2 f2 f2 f2 f2[f8]f8
0x7f8f5d74a480: f3 f3 f3 f3 00 00 00 00 00 00 00 00 00 00 00 00
0x7f8f5d74a500: 00 00 00 00 00 00 00 00 00 00 00 00 00 00 00 00
0x7f8f5d74a580: 00 00 00 00 00 00 00 00 00 00 00 00 00 00 00 00
0x7f8f5d74a600: 00 00 00 00 00 00 00 00 00 00 00 00 00 00 00 00
0x7f8f5d74a680: 00 00 00 00 00 00 00 00 00 00 00 00 00 00 00 00
Shadow byte legend (one shadow byte represents 8 application bytes):
Addressable: 00
Partially addressable: 01 02 03 04 05 06 07
Heap left redzone: fa
Freed heap region: fd
Stack left redzone: f1
Stack mid redzone: f2
Stack right redzone: f3
Stack after return: f5
Stack use after scope: f8
Global redzone: f9
Global init order: f6
Poisoned by user: f7
Container overflow: fc
Array cookie: ac
Intra object redzone: bb
ASan internal: fe
Left alloca redzone: ca
Right alloca redzone: cb
==2034964==ABORTING
According to the address sanitizer the issue is the reference capture in line 377: (confusion)const auto to_category = [&](const auto &cell) {
of matrix_free/tools.h
.
Ah, this is a really subtle point:
cell->face(i)
creates a temporary which is our TriaIterator
.
(Note: auto &iterator = cell->face(i)
throws an error as expected)TriaAccessor
and binding to this object with a reference is fine.TriaIterator
temporary and we are left with a use after scope because the Accessor was stored in the temporary Iterator.Minimal reproducer:
#include <deal.II/grid/tria.h>
#include <deal.II/grid/grid_generator.h>
using namespace dealii;
int main()
{
Triangulation<2> triangulation;
GridGenerator::hyper_cube(triangulation);
const auto cell = triangulation.begin();
auto &face = *cell->face(0);
face.boundary_id();
}
Keeping the iterator temporary alive "fixes" the problem:
int main()
{
Triangulation<2> triangulation;
GridGenerator::hyper_cube(triangulation);
const auto cell = triangulation.begin();
const auto iterator = cell->face(0);
auto &face = *iterator;
face.boundary_id();
}
Ah, you're right. *cell->face()
returns a reference to the accessor, which is a sub-object of the temporary. We bind to that sub-object, which does not extend the lifetime of the surrounding iterator object. So it really isn't about the const
or not, but about the sub-object.
I think the only proper fix is to return a value and not a reference.
Currently the syntax *cell->face()
is just ill formed. And semantically it is hard to explain why *cell->face()
should be a disallowed operation :-)
I am baffled that we never spotted this. (Actually thinking about this a bit more: the unavailability of auto
saved us. No user would be able to figure out the type of this thing in the first place :stuck_out_tongue_winking_eye: )
That's a problem for after the release :-)
I see a segmentation fault in release mode in the test
matrix_free/pbc_orientation_02.release
:I am working on producing a better stack trace. Related failing tests (failing in optimized release mode)
In reference to #15383