From 9b3a70a7391a5c8a10f0bc4557b9e381c784159c Mon Sep 17 00:00:00 2001 From: "dakshina.ilangova" Date: Sat, 5 Apr 2014 19:33:09 +0000 Subject: [PATCH] Added tests for lac/ShiftedMatrix git-svn-id: https://svn.dealii.org/trunk@32720 0785d39b-7218-0410-832d-ea1e28bc413d --- tests/lac/shifted_matrix_01.cc | 73 ++++++++++++++++++++++ tests/lac/shifted_matrix_01.output | 5 ++ tests/lac/shifted_matrix_02.cc | 71 ++++++++++++++++++++++ tests/lac/shifted_matrix_02.output | 4 ++ tests/lac/shifted_matrix_03.cc | 81 +++++++++++++++++++++++++ tests/lac/shifted_matrix_03.output | 7 +++ tests/lac/shifted_matrix_04.cc | 88 +++++++++++++++++++++++++++ tests/lac/shifted_matrix_04.output | 10 +++ tests/lac/shifted_matrix_05.cc | 97 ++++++++++++++++++++++++++++++ tests/lac/shifted_matrix_05.output | 14 +++++ 10 files changed, 450 insertions(+) create mode 100644 tests/lac/shifted_matrix_01.cc create mode 100644 tests/lac/shifted_matrix_01.output create mode 100644 tests/lac/shifted_matrix_02.cc create mode 100644 tests/lac/shifted_matrix_02.output create mode 100644 tests/lac/shifted_matrix_03.cc create mode 100644 tests/lac/shifted_matrix_03.output create mode 100644 tests/lac/shifted_matrix_04.cc create mode 100644 tests/lac/shifted_matrix_04.output create mode 100644 tests/lac/shifted_matrix_05.cc create mode 100644 tests/lac/shifted_matrix_05.output diff --git a/tests/lac/shifted_matrix_01.cc b/tests/lac/shifted_matrix_01.cc new file mode 100644 index 0000000000..2d414c4148 --- /dev/null +++ b/tests/lac/shifted_matrix_01.cc @@ -0,0 +1,73 @@ +// --------------------------------------------------------------------- +// $Id: shifted_matrix.cc 32491 2014-03-13 dilangov $ +// +// Copyright (C) 2014 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. +// +// --------------------------------------------------------------------- + +// check ShiftedMatrix::checkConstructor + +#include "../tests.h" +#include +#include +#include +#include + +#include +#include +#include + +template + void + checkConstructor(FullMatrix &A, double sigma) + { + deallog << "constructor" << std::endl; + + ShiftedMatrix < FullMatrix > S(A, sigma); + + deallog << "Multiplying with all ones vector" << std::endl; + Vector V(A.n()); + for (unsigned int i = 0; i < V.size(); ++i) + V(i) = 1; + + Vector O(A.m()); + + S.vmult(O, V); + + // Check the dimensions of the result matrix + Assert(A.m() == O.size(), ExcInternalError()); + deallog << "Dimensions of result vector verified" << std::endl; + + for (unsigned int i = 0; i < O.size(); ++i) + deallog << O(i) << '\t'; + deallog << std::endl; + } + +int +main() +{ + std::ofstream logfile("output"); + deallog << std::fixed; + deallog << std::setprecision(4); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + const double Adata[] = + { 2, 3, 4, 5 }; + + FullMatrix A(2, 2); + + A.fill(Adata); + + checkConstructor(A, 2); +} diff --git a/tests/lac/shifted_matrix_01.output b/tests/lac/shifted_matrix_01.output new file mode 100644 index 0000000000..b90077284d --- /dev/null +++ b/tests/lac/shifted_matrix_01.output @@ -0,0 +1,5 @@ + +DEAL::constructor +DEAL::Multiplying with all ones vector +DEAL::Dimensions of result vector verified +DEAL::7.0000 11.0000 diff --git a/tests/lac/shifted_matrix_02.cc b/tests/lac/shifted_matrix_02.cc new file mode 100644 index 0000000000..1d7c2965ef --- /dev/null +++ b/tests/lac/shifted_matrix_02.cc @@ -0,0 +1,71 @@ +// --------------------------------------------------------------------- +// $Id: shifted_matrix.cc 32491 2014-03-13 dilangov $ +// +// Copyright (C) 2014 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. +// +// --------------------------------------------------------------------- + +// check ShiftedMatrix::checkVmult + +#include "../tests.h" +#include +#include +#include +#include + +#include +#include +#include + +template + void + checkVmult(FullMatrix &A, double sigma, Vector &V) + { + deallog << "vmult" << std::endl; + + ShiftedMatrix < FullMatrix > S(A, sigma); + Vector O(A.m()); + + S.vmult(O, V); + + // Check the dimensions of the result matrix + Assert(A.m() == O.size(), ExcInternalError()); + deallog << "Dimensions of result vector verified" << std::endl; + + for (unsigned int i = 0; i < O.size(); ++i) + deallog << O(i) << '\t'; + deallog << std::endl; + } + +int +main() +{ + std::ofstream logfile("output"); + deallog << std::fixed; + deallog << std::setprecision(4); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + const double Adata[] = + { 2, 3, 4, 5 }; + + FullMatrix A(2, 2); + + A.fill(Adata); + + Vector V(2); + V(0) = 1; + V(1) = 2; + + checkVmult(A, 2, V); +} diff --git a/tests/lac/shifted_matrix_02.output b/tests/lac/shifted_matrix_02.output new file mode 100644 index 0000000000..5c9d5fcf3b --- /dev/null +++ b/tests/lac/shifted_matrix_02.output @@ -0,0 +1,4 @@ + +DEAL::vmult +DEAL::Dimensions of result vector verified +DEAL::10.0000 18.0000 diff --git a/tests/lac/shifted_matrix_03.cc b/tests/lac/shifted_matrix_03.cc new file mode 100644 index 0000000000..fc4c1dab5d --- /dev/null +++ b/tests/lac/shifted_matrix_03.cc @@ -0,0 +1,81 @@ +// --------------------------------------------------------------------- +// $Id: shifted_matrix.cc 32491 2014-03-13 dilangov $ +// +// Copyright (C) 2014 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. +// +// --------------------------------------------------------------------- + +// check ShiftedMatrix::checkResidual + +#include "../tests.h" +#include +#include +#include +#include + +#include +#include +#include + +template + void + checkResidual(FullMatrix &A, double sigma, Vector &V, + Vector &R) + { + deallog << "residual" << std::endl; + + ShiftedMatrix < FullMatrix > S(A, sigma); + Vector O(A.m()); + double residual; + + residual = S.residual(O, V, R); + + // Check the dimensions of the result matrix + Assert(A.m() == O.size(), ExcInternalError()); + deallog << "Dimensions of result vector verified" << std::endl; + + deallog << "Residual vector" << std::endl; + for (unsigned int i = 0; i < O.size(); ++i) + deallog << O(i) << '\t'; + deallog << std::endl; + + deallog << "Residual value" << std::endl; + deallog << residual << std::endl; + } + +int +main() +{ + std::ofstream logfile("output"); + deallog << std::fixed; + deallog << std::setprecision(4); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + const double Adata[] = + { 2, 3, 4, 5 }; + + FullMatrix A(2, 2); + + A.fill(Adata); + + Vector V(2); + V(0) = 1; + V(1) = 2; + + Vector R(2); + R(0) = 1; + R(1) = 1; + + checkResidual(A, 2, V, R); +} diff --git a/tests/lac/shifted_matrix_03.output b/tests/lac/shifted_matrix_03.output new file mode 100644 index 0000000000..64cf220694 --- /dev/null +++ b/tests/lac/shifted_matrix_03.output @@ -0,0 +1,7 @@ + +DEAL::residual +DEAL::Dimensions of result vector verified +DEAL::Residual vector +DEAL::-9.0000 -17.0000 +DEAL::Residual value +DEAL::19.2354 diff --git a/tests/lac/shifted_matrix_04.cc b/tests/lac/shifted_matrix_04.cc new file mode 100644 index 0000000000..32355a044a --- /dev/null +++ b/tests/lac/shifted_matrix_04.cc @@ -0,0 +1,88 @@ +// --------------------------------------------------------------------- +// $Id: shifted_matrix.cc 32491 2014-03-13 dilangov $ +// +// Copyright (C) 2014 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. +// +// --------------------------------------------------------------------- + +// check ShiftedMatrix::checkSetSigma + +#include "../tests.h" +#include +#include +#include +#include + +#include +#include +#include + +template + void + checkSetSigma(FullMatrix &A, double sigma) + { + deallog << "shift(sigma)" << std::endl; + + deallog << "Init ShiftedMatrix with sigma=0" << std::endl; + ShiftedMatrix < FullMatrix > S(A, 0); + + deallog << "Multiplying with all ones vector" << std::endl; + Vector V(A.n()); + for (unsigned int i = 0; i < V.size(); ++i) + V(i) = 1; + + Vector O(A.m()); + + S.vmult(O, V); + + // Check the dimensions of the result matrix + Assert(A.m() == O.size(), ExcInternalError()); + deallog << "Dimensions of result vector verified" << std::endl; + + for (unsigned int i = 0; i < O.size(); ++i) + deallog << O(i) << '\t'; + deallog << std::endl; + + deallog << "Setting new sigma value" << std::endl; + S.shift(sigma); + + deallog << "Multiplying with all ones vector" << std::endl; + S.vmult(O, V); + + // Check the dimensions of the result matrix + Assert(A.m() == O.size(), ExcInternalError()); + deallog << "Dimensions of result vector verified" << std::endl; + + for (unsigned int i = 0; i < O.size(); ++i) + deallog << O(i) << '\t'; + deallog << std::endl; + } + +int +main() +{ + std::ofstream logfile("output"); + deallog << std::fixed; + deallog << std::setprecision(4); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + const double Adata[] = + { 2, 3, 4, 5 }; + + FullMatrix A(2, 2); + + A.fill(Adata); + + checkSetSigma(A, 2); +} diff --git a/tests/lac/shifted_matrix_04.output b/tests/lac/shifted_matrix_04.output new file mode 100644 index 0000000000..822d2ddc4e --- /dev/null +++ b/tests/lac/shifted_matrix_04.output @@ -0,0 +1,10 @@ + +DEAL::shift(sigma) +DEAL::Init ShiftedMatrix with sigma=0 +DEAL::Multiplying with all ones vector +DEAL::Dimensions of result vector verified +DEAL::5.0000 9.0000 +DEAL::Setting new sigma value +DEAL::Multiplying with all ones vector +DEAL::Dimensions of result vector verified +DEAL::7.0000 11.0000 diff --git a/tests/lac/shifted_matrix_05.cc b/tests/lac/shifted_matrix_05.cc new file mode 100644 index 0000000000..8af92c3de7 --- /dev/null +++ b/tests/lac/shifted_matrix_05.cc @@ -0,0 +1,97 @@ +// --------------------------------------------------------------------- +// $Id: shifted_matrix.cc 32491 2014-03-13 dilangov $ +// +// Copyright (C) 2014 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. +// +// --------------------------------------------------------------------- + +// check ShiftedMatrix::checkGetSigma + +#include "../tests.h" +#include +#include +#include +#include + +#include +#include +#include + +template + void + checkGetSigma(FullMatrix &A) + { + deallog << "shift()" << std::endl; + + deallog << "Init ShiftedMatrix with sigma=0" << std::endl; + ShiftedMatrix < FullMatrix > S(A, 0); + + deallog << "Multiplying with all ones vector" << std::endl; + Vector V(A.n()); + for (unsigned int i = 0; i < V.size(); ++i) + V(i) = 1; + + Vector O(A.m()); + + S.vmult(O, V); + + // Check the dimensions of the result matrix + Assert(A.m() == O.size(), ExcInternalError()); + deallog << "Dimensions of result vector verified" << std::endl; + + for (unsigned int i = 0; i < O.size(); ++i) + deallog << O(i) << '\t'; + deallog << std::endl; + + deallog << "Setting new sigma value by incrementing old value by 1" + << std::endl; + + deallog << "Old sigma value" << std::endl; + double sigma = S.shift(); + deallog << sigma << std::endl; + sigma = sigma + 1; + deallog << "New sigma value" << std::endl; + deallog << sigma << std::endl; + + S.shift(sigma); + + deallog << "Multiplying with all ones vector" << std::endl; + S.vmult(O, V); + + // Check the dimensions of the result matrix + Assert(A.m() == O.size(), ExcInternalError()); + deallog << "Dimensions of result vector verified" << std::endl; + + for (unsigned int i = 0; i < O.size(); ++i) + deallog << O(i) << '\t'; + deallog << std::endl; + } + +int +main() +{ + std::ofstream logfile("output"); + deallog << std::fixed; + deallog << std::setprecision(4); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + const double Adata[] = + { 2, 3, 4, 5 }; + + FullMatrix A(2, 2); + + A.fill(Adata); + + checkGetSigma(A); +} diff --git a/tests/lac/shifted_matrix_05.output b/tests/lac/shifted_matrix_05.output new file mode 100644 index 0000000000..27f92db5f1 --- /dev/null +++ b/tests/lac/shifted_matrix_05.output @@ -0,0 +1,14 @@ + +DEAL::shift() +DEAL::Init ShiftedMatrix with sigma=0 +DEAL::Multiplying with all ones vector +DEAL::Dimensions of result vector verified +DEAL::5.0000 9.0000 +DEAL::Setting new sigma value by incrementing old value by 1 +DEAL::Old sigma value +DEAL::0 +DEAL::New sigma value +DEAL::1.0000 +DEAL::Multiplying with all ones vector +DEAL::Dimensions of result vector verified +DEAL::6.0000 10.0000 -- 2.39.5