From: Wolfgang Bangerth Date: Mon, 9 Nov 2015 01:18:22 +0000 (-0600) Subject: Add test. X-Git-Tag: v8.4.0-rc2~240^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=refs%2Fpull%2F1850%2Fhead;p=dealii.git Add test. This is a simplified test of FESystems whose base elements are again FESystems. The test helped me debug a couple of issues, and might be of help again at a later time. It is a simpler version of fe/shapes_system. --- diff --git a/tests/fe/shapes_system_02.cc b/tests/fe/shapes_system_02.cc new file mode 100644 index 0000000000..5161356138 --- /dev/null +++ b/tests/fe/shapes_system_02.cc @@ -0,0 +1,64 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2013 - 2015 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 at +// the top level of the deal.II distribution. +// +// --------------------------------------------------------------------- + + +// a simplified version of shapes_system.cc + +#include "../tests.h" +#include "shapes.h" +#include +#include +#include +#include +#include +#include +#include +#include + +#define PRECISION 8 + +template +void plot_FE_System_shape_functions() +{ + MappingQGeneric m(1); + + FESystem p3(FE_Q(1), 1, + FESystem (FE_Q(1),2), 2); + test_compute_functions(m, p3, "System_1"); + + FESystem p4(p3, 2, + FESystem (p3,3), 1, + p3, 1); + test_compute_functions(m, p4, "System_2"); +} + + +int +main() +{ + std::ofstream logfile ("output"); + deallog << std::setprecision(PRECISION) << std::fixed; + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + deallog << "FE_System<1>" << std::endl; + plot_FE_System_shape_functions<1>(); + deallog << "FE_System<2>" << std::endl; + plot_FE_System_shape_functions<2>(); + deallog << "FE_System<3>" << std::endl; + plot_FE_System_shape_functions<3>(); + + return 0; +} diff --git a/tests/fe/shapes_system_02.output b/tests/fe/shapes_system_02.output new file mode 100644 index 0000000000..34dd77fdf5 --- /dev/null +++ b/tests/fe/shapes_system_02.output @@ -0,0 +1,4 @@ + +DEAL::FE_System<1> +DEAL::FE_System<2> +DEAL::FE_System<3>