From: Peter Munch Date: Tue, 16 Aug 2022 07:32:42 +0000 (+0200) Subject: FEEvaluationData: add virtual destructor X-Git-Tag: v9.5.0-rc1~1029^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F14190%2Fhead;p=dealii.git FEEvaluationData: add virtual destructor --- diff --git a/include/deal.II/matrix_free/fe_evaluation_data.h b/include/deal.II/matrix_free/fe_evaluation_data.h index f6854fb0c6..03be313a7a 100644 --- a/include/deal.II/matrix_free/fe_evaluation_data.h +++ b/include/deal.II/matrix_free/fe_evaluation_data.h @@ -149,6 +149,11 @@ public: FEEvaluationData & operator=(const FEEvaluationData &other); + /** + * Destructor. + */ + virtual ~FEEvaluationData() = default; + /** * Sets the pointers for values, gradients, hessians to the central * scratch_data_array inside the given scratch array, for a given number of diff --git a/tests/matrix_free/fe_evaluation_dynamic_cast.cc b/tests/matrix_free/fe_evaluation_dynamic_cast.cc new file mode 100644 index 0000000000..ed60607a57 --- /dev/null +++ b/tests/matrix_free/fe_evaluation_dynamic_cast.cc @@ -0,0 +1,69 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2019 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + + + +// Test that FEEvaluationData is virtual. + +#include + +#include + +#include + +#include + +#include +#include + +#include "../tests.h" + +int +main() +{ + initlog(); + + const unsigned int dim = 2; + using Number = double; + using VectorizedArrayType = VectorizedArray; + + Triangulation tria; + GridGenerator::hyper_cube(tria); + DoFHandler dof(tria); + dof.distribute_dofs(FE_Q(1)); + + MatrixFree dummy; + + dummy.reinit(MappingQ1{}, + dof, + AffineConstraints(), + QGauss<1>(2)); + + std::unique_ptr> phi_base = + std::make_unique>( + dummy); + + if (dynamic_cast *>( + phi_base.get())) + deallog << "OK!" << std::endl; + else + deallog << "Fail!" << std::endl; + + if (dynamic_cast *>( + phi_base.get())) + deallog << "Fail!" << std::endl; + else + deallog << "OK!" << std::endl; +} diff --git a/tests/matrix_free/fe_evaluation_dynamic_cast.output b/tests/matrix_free/fe_evaluation_dynamic_cast.output new file mode 100644 index 0000000000..043986ad9a --- /dev/null +++ b/tests/matrix_free/fe_evaluation_dynamic_cast.output @@ -0,0 +1,3 @@ + +DEAL::OK! +DEAL::OK!