Skip to content

Commit

Permalink
Merge pull request #14190 from peterrum/feevaldata_vd
Browse files Browse the repository at this point in the history
FEEvaluationData: add virtual destructor
  • Loading branch information
kronbichler committed Aug 23, 2022
2 parents 9cb3474 + d67a019 commit 55ce20b
Show file tree
Hide file tree
Showing 3 changed files with 77 additions and 0 deletions.
5 changes: 5 additions & 0 deletions include/deal.II/matrix_free/fe_evaluation_data.h
Original file line number Diff line number Diff line change
Expand Up @@ -149,6 +149,11 @@ class FEEvaluationData
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
Expand Down
69 changes: 69 additions & 0 deletions tests/matrix_free/fe_evaluation_dynamic_cast.cc
Original file line number Diff line number Diff line change
@@ -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 <deal.II/dofs/dof_handler.h>

#include <deal.II/fe/fe_q.h>

#include <deal.II/grid/grid_generator.h>

#include <deal.II/lac/affine_constraints.h>

#include <deal.II/matrix_free/fe_evaluation.h>
#include <deal.II/matrix_free/matrix_free.h>

#include "../tests.h"

int
main()
{
initlog();

const unsigned int dim = 2;
using Number = double;
using VectorizedArrayType = VectorizedArray<Number>;

Triangulation<dim> tria;
GridGenerator::hyper_cube(tria);
DoFHandler<dim> dof(tria);
dof.distribute_dofs(FE_Q<dim>(1));

MatrixFree<dim, Number, VectorizedArrayType> dummy;

dummy.reinit(MappingQ1<dim>{},
dof,
AffineConstraints<double>(),
QGauss<1>(2));

std::unique_ptr<FEEvaluationData<dim, VectorizedArrayType, false>> phi_base =
std::make_unique<FEEvaluation<dim, -1, 0, 1, Number, VectorizedArrayType>>(
dummy);

if (dynamic_cast<FEEvaluation<dim, -1, 0, 1, Number, VectorizedArrayType> *>(
phi_base.get()))
deallog << "OK!" << std::endl;
else
deallog << "Fail!" << std::endl;

if (dynamic_cast<FEEvaluation<dim, -1, 0, 2, Number, VectorizedArrayType> *>(
phi_base.get()))
deallog << "Fail!" << std::endl;
else
deallog << "OK!" << std::endl;
}
3 changes: 3 additions & 0 deletions tests/matrix_free/fe_evaluation_dynamic_cast.output
Original file line number Diff line number Diff line change
@@ -0,0 +1,3 @@

DEAL::OK!
DEAL::OK!

0 comments on commit 55ce20b

Please sign in to comment.