From 27e7e8be660ca87c210f779bfc2f59a8fccb5de0 Mon Sep 17 00:00:00 2001 From: Martin Kronbichler Date: Mon, 30 May 2022 22:34:59 +0200 Subject: [PATCH] New test case by Pascal Kraft --- tests/numerics/project_complex.cc | 137 ++++++++++++++++++++++++++ tests/numerics/project_complex.output | 2 + 2 files changed, 139 insertions(+) create mode 100644 tests/numerics/project_complex.cc create mode 100644 tests/numerics/project_complex.output diff --git a/tests/numerics/project_complex.cc b/tests/numerics/project_complex.cc new file mode 100644 index 0000000000..9254ff0ef9 --- /dev/null +++ b/tests/numerics/project_complex.cc @@ -0,0 +1,137 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2006 - 2018 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. +// +// --------------------------------------------------------------------- + + + +// Check that VectorTools::project can be instantiated for std::complex. + +#include +#include +#include +#include + +#include + +#include + +#include +#include +#include + +#include +#include +#include +#include + +#include +#include +#include +#include + +#include "../tests.h" + +using ComplexNumber = std::complex; +using Constraints = dealii::AffineConstraints; +using Position = dealii::Point<3, double>; + + +template void +dealii::VectorTools::project<3, dealii::Vector, 3>( + const dealii::Mapping<3, 3> &, + const dealii::DoFHandler<3, 3> &, + const Constraints &, + const dealii::Quadrature<3> &, + const dealii::Function<3, ComplexNumber> &, + dealii::Vector &, + const bool, + const dealii::Quadrature<2> &, + const bool); + +class Field : public dealii::Function<3, ComplexNumber> +{ +public: + Field() + {} + + ComplexNumber + value(const Position &, const unsigned int) const + { + return ComplexNumber(0, 0); + } + + void + vector_value(const Position &, dealii::Vector &value) const + { + value.reinit(3); + for (unsigned int i = 0; i < 3; i++) + { + value[i] = {0, 0}; + } + } + + dealii::Tensor<1, 3, ComplexNumber> + curl(const Position &) const + { + dealii::Tensor<1, 3, ComplexNumber> ret; + for (unsigned int i = 0; i < 3; i++) + { + ret[i] = {0, 0}; + } + return ret; + } + + dealii::Tensor<1, 3, ComplexNumber> + val(const Position &) const + { + dealii::Tensor<1, 3, ComplexNumber> ret; + for (unsigned int i = 0; i < 3; i++) + { + ret[i] = {0, 0}; + } + return ret; + } +}; + + +int +main() +{ + initlog(); + + Triangulation<3> triangulation; + Field f; + std::vector repetitions; + repetitions.push_back(5); + repetitions.push_back(5); + repetitions.push_back(5); + dealii::Point<3, double> l(0, 0, 0); + dealii::Point<3, double> r(1, 1, 1); + dealii::GridGenerator::subdivided_hyper_rectangle( + triangulation, repetitions, l, r, false); + + FE_Q<3> fe(1); + dealii::DoFHandler<3> dof_handler(triangulation); + dof_handler.distribute_dofs(fe); + Constraints local_constraints; + local_constraints.close(); + + dealii::QGauss<3> quadrature(2); + dealii::Vector interpolated_exact_solution( + dof_handler.n_dofs()); + VectorTools::project( + dof_handler, local_constraints, quadrature, f, interpolated_exact_solution); + + deallog << "OK" << std::endl; +} diff --git a/tests/numerics/project_complex.output b/tests/numerics/project_complex.output new file mode 100644 index 0000000000..0fd8fc12f0 --- /dev/null +++ b/tests/numerics/project_complex.output @@ -0,0 +1,2 @@ + +DEAL::OK -- 2.39.5