PFASST++
test_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/serial_mlsdc.cpp"
15 #undef PFASST_UNIT_TESTING
17 
18 /*
19  * parameterized test fixture with number of nodes as parameter
20  */
22  : public TestWithParam<tuple<size_t, pfasst::quadrature::QuadratureType>>
23 {
24  protected:
25  size_t nnodes;
26  size_t niters;
27  size_t nlevs = 3;
28  vector<size_t> nsteps;
29  vector<double> err;
30  vector<double> convrate;
32 
33  public:
34  virtual void SetUp()
35  {
36  this->nnodes = get<0>(GetParam());
37  this->nodetype = get<1>(GetParam());
38 
39  switch (this->nodetype)
40  {
42  this->niters = 2 * this->nnodes - 2;
43  this->nsteps = { 4, 8, 16, 32 };
44  nlevs = 2;
45  break;
46 
48  this->niters = 2 * this->nnodes;
49  this->nsteps = { 4, 8, 16, 32 };
50  break;
51 
53  this->niters = this->nnodes;
54  this->nsteps = { 4, 8, 16, 32 };
55  break;
56 
58  this->niters = this->nnodes;
59  this->nsteps = { 4, 8, 16, 32 };
60  break;
61 
62  default:
63  break;
64  }
65 
66  // run to compute errors
67  for (size_t i = 0; i < this->nsteps.size(); ++i) {
68  auto dt = 0.5 / double(this->nsteps[i]);
69 
70  auto errors_and_residuals = run_serial_mlsdc(
71  nlevs, // nlevs
72  this->nsteps[i], // nsteps
73  dt, // dt
74  this->niters, // niters
75  this->nnodes, // nnodes
76  128 // ndofs
77  );
78 
79  auto errors = get<0>(errors_and_residuals);
80  auto last_error = pfasst::examples::advection_diffusion::ktype(this->nsteps[i]-1, this->niters-1);
81  this->err.push_back(errors.at(last_error));
82  }
83 
84  // compute convergence rates
85  for (size_t i = 0; i < this->nsteps.size() - 1; ++i) {
86  this->convrate.push_back(log10(this->err.at(i+1) / this->err.at(i)) /
87  log10(double(this->nsteps.at(i)) / double(this->nsteps.at(i + 1))));
88  }
89  }
90 };
91 
92 /*
93  * The test below verifies that the code reproduces the theoretically expected rate of convergence
94  * (or better).
95  */
97 {
98  string quad;
99  int order = 1;
100  double fudge = 1.0;
101 
102  switch (this->nodetype)
103  {
105  order = 2 * this->nnodes - 2;
106  quad = "Gauss-Lobatto";
107  fudge = 0.9;
108  break;
109 
111  order = 2 * this->nnodes;
112  quad = "Gauss-Legendre";
113  break;
114 
116  order = 2 * this->nnodes - 1;
117  quad = "Gauss-Radau";
118  break;
119 
121  order = this->nnodes;
122  quad = "Clenshaw-Curtis";
123  break;
124 
126  order = this->nnodes;
127  quad = "Uniform";
128  break;
129 
130  default:
131  EXPECT_TRUE(false);
132  break;
133  }
134 
135  for (size_t i = 0; i < this->nsteps.size() - 1; ++i) {
136  EXPECT_THAT(convrate[i], Ge<double>(fudge * order)) << "Convergence rate for "
137  << this->nnodes << " " << quad << " nodes"
138  << " for nsteps " << this->nsteps[i]
139  << " not within expected range.";
140  }
141 }
142 
143 INSTANTIATE_TEST_CASE_P(AdvectionDiffusionPFASST, ConvergenceTest,
144  Combine(Range<size_t>(5, 6),
148 );
149 
150 int main(int argc, char** argv)
151 {
152  testing::InitGoogleTest(&argc, argv);
153  pfasst::init(argc, argv,
156  return RUN_ALL_TESTS();
157 }
pfasst::quadrature::QuadratureType nodetype
QuadratureType
Quadrature type descriptors.
Definition: interface.hpp:32
STL namespace.
static void init(int argc, char **argv, std::function< void()> opts=nullptr, std::function< void()> logs=nullptr)
Definition: pfasst.hpp:13
int main(int argc, char **argv)
TEST_P(ConvergenceTest, AllNodes)
advection-diffusion sweeper with semi-implicit time-integration.
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)))
float dt
Definition: plot.py:10
tuple< error_map, residual_map > run_serial_mlsdc(size_t nlevs, size_t nsteps_in=4, double step_size_in=0.01, size_t num_iter_in=8, size_t nnodes_in=5, size_t ndofs_in=128)
Advection/diffusion example using an encapsulated IMEX sweeper.