PFASST++
test_mpi_advection_diffusion_conv.cpp
Go to the documentation of this file.
1 /*
2  * Convergence tests for the advection-diffusion examples.
3  */
4 #include <cmath>
5 using namespace std;
6 
7 #include <gtest/gtest.h>
8 #include <gmock/gmock.h>
9 using namespace ::testing;
10 
11 #include <pfasst/quadrature.hpp>
12 
13 #define PFASST_UNIT_TESTING
14 #include "../examples/advection_diffusion/mpi_pfasst.cpp"
15 #undef PFASST_UNIT_TESTING
17 
18 /*
19  * parameterized test fixture with number of nodes as parameter
20  */
21 class ConvergenceTest
22  : public TestWithParam<tuple<size_t, pfasst::quadrature::QuadratureType>>
23 {
24  protected:
25  size_t nnodes;
26  size_t niters;
27  vector<size_t> nsteps;
28  vector<double> err;
29  vector<double> convrate;
31 
32  public:
33  virtual void SetUp()
34  {
35  this->nnodes = get<0>(GetParam());
36  this->nodetype = get<1>(GetParam());
37 
38  switch (this->nodetype)
39  {
41  this->niters = 2 * this->nnodes - 2;
42  this->nsteps = { 4, 8, 16, 32 };
43  break;
44 
46  this->niters = 2 * this->nnodes;
47  this->nsteps = { 4, 8, 16, 32 };
48  break;
49 
51  this->niters = this->nnodes;
52  this->nsteps = { 4, 8, 16, 32 };
53  break;
54 
56  this->niters = this->nnodes;
57  this->nsteps = { 4, 8, 16, 32 };
58  break;
59 
60  default:
61  break;
62  }
63 
64  // run to compute errors
65  for (size_t i = 0; i < this->nsteps.size(); ++i) {
66  auto dt = 0.5 / double(this->nsteps[i]);
67 
68  auto errors = run_mpi_pfasst(
69  0.0, // abs_res_tol
70  0.0, // rel_res_tol
71  this->niters, // niters
72  this->nsteps[i], // nsteps
73  dt, // dt
74  128, // ndofs_f
75  64, // ndofs_c
76  this->nnodes, // nnodes_f
77  (this->nnodes+1)/2-1 // nnodes_c
78  );
79 
80  int rank, size;
81  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
82  MPI_Comm_size(MPI_COMM_WORLD, &size);
83 
84  auto last_error = pfasst::examples::advection_diffusion::ktype(this->nsteps[i]-size+rank, this->niters-1);
85  this->err.push_back(errors.at(last_error));
86  }
87 
88  // compute convergence rates
89  for (size_t i = 0; i < this->nsteps.size() - 1; ++i) {
90  this->convrate.push_back(log10(this->err[i+1] / this->err[i]) /
91  log10(double(this->nsteps[i]) / double(this->nsteps[i + 1])));
92  }
93  }
94 };
95 
96 /*
97  * The test below verifies that the code reproduces the theoretically expected rate of convergence
98  * (or better) ON THE LAST PROCESSOR.
99  */
101 {
102  int order;
103  int rank, size;
104  string quad;
105 
106  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
107  MPI_Comm_size(MPI_COMM_WORLD, &size);
108  if (rank != size-1) {
109  EXPECT_TRUE(true);
110  return;
111  }
112 
113  switch (this->nodetype)
114  {
116  order = 2 * this->nnodes - 2;
117  quad = "Gauss-Lobatto";
118  break;
119 
121  order = 2 * this->nnodes;
122  quad = "Gauss-Legendre";
123  break;
124 
126  order = 2 * this->nnodes - 1;
127  quad = "Gauss-Radau";
128  break;
129 
131  order = this->nnodes;
132  quad = "Clenshaw-Curtis";
133  break;
134 
136  order = this->nnodes;
137  quad = "Uniform";
138  break;
139 
140  default:
141  EXPECT_TRUE(false);
142  break;
143  }
144 
145  for (size_t i = 0; i < this->nsteps.size() - 1; ++i) {
146  EXPECT_THAT(convrate[i], Ge<double>(order)) << "Convergence rate for "
147  << this->nnodes << " " << quad << " nodes"
148  << " for nsteps " << this->nsteps[i]
149  << " not within expected range.";
150  }
151 }
152 
153 INSTANTIATE_TEST_CASE_P(AdvectionDiffusionPFASST, ConvergenceTest,
154  Combine(Range<size_t>(5, 6),
156  // pfasst::quadrature::QuadratureType::GaussLegendre,
157  // pfasst::quadrature::QuadratureType::GaussRadau,
160 );
161 
162 int main(int argc, char** argv)
163 {
164  testing::InitGoogleTest(&argc, argv);
165  MPI_Init(&argc, &argv);
166  pfasst::init(argc, argv,
169  int result = 1, max_result; // GTest return value 1 (failure), 0 (success)
170  result = RUN_ALL_TESTS();
171  MPI_Allreduce(&result, &max_result, 1, MPI_INT, MPI_MAX, MPI_COMM_WORLD);
172  MPI_Finalize();
173  return max_result;
174 }
QuadratureType
Quadrature type descriptors.
Definition: interface.hpp:32
STL namespace.
int main(int argc, char **argv)
static void init(int argc, char **argv, std::function< void()> opts=nullptr, std::function< void()> logs=nullptr)
Definition: pfasst.hpp:13
INSTANTIATE_TEST_CASE_P(AdvectionDiffusionPFASST, ConvergenceTest, Combine(Range< size_t >(5, 6), Values(pfasst::quadrature::QuadratureType::GaussLobatto, pfasst::quadrature::QuadratureType::ClenshawCurtis, pfasst::quadrature::QuadratureType::Uniform)))
error_map run_mpi_pfasst(const double abs_res_tol, const double rel_res_tol, const size_t niters, const size_t nsteps, const double dt, const size_t ndofs_f, const size_t ndofs_c, const size_t nnodes_f, const size_t nnodes_c)
Advection/diffusion example using an encapsulated IMEX sweeper.
Definition: mpi_pfasst.cpp:43
TEST_P(ConvergenceTest, AllNodes)
advection-diffusion sweeper with semi-implicit time-integration.
float dt
Definition: plot.py:10