From: Wolfgang Bangerth Date: Tue, 14 Apr 2020 03:27:38 +0000 (-0600) Subject: Add a test case. X-Git-Tag: v9.2.0-rc1~237^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=0935c219d0bfb0e4c15930d214a7fb1ae4d9e095;p=dealii.git Add a test case. --- diff --git a/tests/data_out/data_out_complex_postprocessor_03.cc b/tests/data_out/data_out_complex_postprocessor_03.cc new file mode 100644 index 0000000000..ba0e39d34e --- /dev/null +++ b/tests/data_out/data_out_complex_postprocessor_03.cc @@ -0,0 +1,205 @@ +// --------------------------------------------------------------------- +// +// 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. +// +// --------------------------------------------------------------------- + + +// Check DataOut for complex vectors, using two postprocessors. The +// testcase is extracted from the step-58 tutorial. There was an +// indexing error when using more than one postprocessor. + +#include + +#include +#include +#include + +#include +#include + +#include +#include +#include + +#include + +#include + +#include +#include +#include + +#include "../tests.h" + + +namespace DataPostprocessors +{ + template + class ComplexMagnitude : public DataPostprocessorScalar + { + public: + ComplexMagnitude(); + + virtual void + evaluate_vector_field( + const DataPostprocessorInputs::Vector &inputs, + std::vector> & computed_quantities) const; + }; + + template + ComplexMagnitude::ComplexMagnitude() + : DataPostprocessorScalar("Magnitude", update_values) + {} + + + template + void + ComplexMagnitude::evaluate_vector_field( + const DataPostprocessorInputs::Vector &inputs, + std::vector> & computed_quantities) const + { + Assert(computed_quantities.size() == inputs.solution_values.size(), + ExcDimensionMismatch(computed_quantities.size(), + inputs.solution_values.size())); + + for (unsigned int q = 0; q < computed_quantities.size(); ++q) + { + Assert(computed_quantities[q].size() == 1, + ExcDimensionMismatch(computed_quantities[q].size(), 1)); + Assert(inputs.solution_values[q].size() == 2, + ExcDimensionMismatch(inputs.solution_values[q].size(), 2)); + + computed_quantities[q](0) = + std::abs(std::complex(inputs.solution_values[q](0), + inputs.solution_values[q](1))); + } + } + + + + template + class ComplexPhase : public DataPostprocessorScalar + { + public: + ComplexPhase(); + + virtual void + evaluate_vector_field( + const DataPostprocessorInputs::Vector &inputs, + std::vector> & computed_quantities) const; + }; + + template + ComplexPhase::ComplexPhase() + : DataPostprocessorScalar("Phase", update_values) + {} + + + template + void + ComplexPhase::evaluate_vector_field( + const DataPostprocessorInputs::Vector &inputs, + std::vector> & computed_quantities) const + { + Assert(computed_quantities.size() == inputs.solution_values.size(), + ExcDimensionMismatch(computed_quantities.size(), + inputs.solution_values.size())); + + for (unsigned int q = 0; q < computed_quantities.size(); ++q) + { + Assert(computed_quantities[q].size() == 1, + ExcDimensionMismatch(computed_quantities[q].size(), 1)); + Assert(inputs.solution_values[q].size() == 2, + ExcDimensionMismatch(inputs.solution_values[q].size(), 2)); + + computed_quantities[q](0) = + std::arg(std::complex(inputs.solution_values[q](0), + inputs.solution_values[q](1))); + } + } + +} // namespace DataPostprocessors + + + +template +void +check() +{ + Triangulation tria; + GridGenerator::hyper_cube(tria, 0., 1.); + + FE_Q fe(1); + DoFHandler dof_handler(tria); + dof_handler.distribute_dofs(fe); + + // Initialize the vector with a constant 0+1i. This will lead to + // easily recognizable magnitudes and phases in the output file. + Vector> v(dof_handler.n_dofs()); + for (unsigned int i = 0; i < v.size(); ++i) + v(i) = std::complex(0, 1); + + const DataPostprocessors::ComplexMagnitude complex_magnitude; + const DataPostprocessors::ComplexPhase complex_phase; + + DataOut data_out; + data_out.attach_dof_handler(dof_handler); + + data_out.add_data_vector(v, complex_magnitude); + data_out.add_data_vector(v, complex_phase); + data_out.add_data_vector(v, "Psi"); + + data_out.build_patches(); + + data_out.write_gnuplot(deallog.get_file_stream()); +} + + + +int +main() +{ + initlog(); + + try + { + check<1>(); + check<2>(); + check<3>(); + } + catch (std::exception &exc) + { + deallog << std::endl + << std::endl + << "----------------------------------------------------" + << std::endl; + deallog << "Exception on processing: " << std::endl + << exc.what() << std::endl + << "Aborting!" << std::endl + << "----------------------------------------------------" + << std::endl; + return 1; + } + catch (...) + { + deallog << std::endl + << std::endl + << "----------------------------------------------------" + << std::endl; + deallog << "Unknown exception!" << std::endl + << "Aborting!" << std::endl + << "----------------------------------------------------" + << std::endl; + return 1; + } +} diff --git a/tests/data_out/data_out_complex_postprocessor_03.with_complex_values=on.output b/tests/data_out/data_out_complex_postprocessor_03.with_complex_values=on.output new file mode 100644 index 0000000000..37ca714dcd --- /dev/null +++ b/tests/data_out/data_out_complex_postprocessor_03.with_complex_values=on.output @@ -0,0 +1,81 @@ + +# This file was generated by the deal.II library. + + +# +# For a description of the GNUPLOT format see the GNUPLOT manual. +# +# +0.00000 1.00000 1.57080 0.00000 1.00000 +1.00000 1.00000 1.57080 0.00000 1.00000 + + +# This file was generated by the deal.II library. + + +# +# For a description of the GNUPLOT format see the GNUPLOT manual. +# +# +0.00000 0.00000 1.00000 1.57080 0.00000 1.00000 +1.00000 0.00000 1.00000 1.57080 0.00000 1.00000 + +0.00000 1.00000 1.00000 1.57080 0.00000 1.00000 +1.00000 1.00000 1.00000 1.57080 0.00000 1.00000 + + +# This file was generated by the deal.II library. + + +# +# For a description of the GNUPLOT format see the GNUPLOT manual. +# +# +0.00000 0.00000 0.00000 1.00000 1.57080 0.00000 1.00000 +1.00000 0.00000 0.00000 1.00000 1.57080 0.00000 1.00000 + + +0.00000 0.00000 0.00000 1.00000 1.57080 0.00000 1.00000 +0.00000 1.00000 0.00000 1.00000 1.57080 0.00000 1.00000 + + +0.00000 0.00000 0.00000 1.00000 1.57080 0.00000 1.00000 +0.00000 0.00000 1.00000 1.00000 1.57080 0.00000 1.00000 + + +1.00000 0.00000 0.00000 1.00000 1.57080 0.00000 1.00000 +1.00000 1.00000 0.00000 1.00000 1.57080 0.00000 1.00000 + + +1.00000 0.00000 0.00000 1.00000 1.57080 0.00000 1.00000 +1.00000 0.00000 1.00000 1.00000 1.57080 0.00000 1.00000 + + +0.00000 1.00000 0.00000 1.00000 1.57080 0.00000 1.00000 +1.00000 1.00000 0.00000 1.00000 1.57080 0.00000 1.00000 + + +0.00000 1.00000 0.00000 1.00000 1.57080 0.00000 1.00000 +0.00000 1.00000 1.00000 1.00000 1.57080 0.00000 1.00000 + + +1.00000 1.00000 0.00000 1.00000 1.57080 0.00000 1.00000 +1.00000 1.00000 1.00000 1.00000 1.57080 0.00000 1.00000 + + +0.00000 0.00000 1.00000 1.00000 1.57080 0.00000 1.00000 +1.00000 0.00000 1.00000 1.00000 1.57080 0.00000 1.00000 + + +0.00000 0.00000 1.00000 1.00000 1.57080 0.00000 1.00000 +0.00000 1.00000 1.00000 1.00000 1.57080 0.00000 1.00000 + + +1.00000 0.00000 1.00000 1.00000 1.57080 0.00000 1.00000 +1.00000 1.00000 1.00000 1.00000 1.57080 0.00000 1.00000 + + +0.00000 1.00000 1.00000 1.00000 1.57080 0.00000 1.00000 +1.00000 1.00000 1.00000 1.00000 1.57080 0.00000 1.00000 + +