7 #include <gtest/gtest.h>
8 #include <gmock/gmock.h>
9 using namespace ::testing;
13 #define PFASST_UNIT_TESTING
14 #include "../examples/advection_diffusion/mpi_pfasst.cpp"
15 #undef PFASST_UNIT_TESTING
22 :
public TestWithParam<tuple<size_t, pfasst::quadrature::QuadratureType>>
27 vector<size_t> nsteps;
29 vector<double> convrate;
35 this->nnodes = get<0>(GetParam());
36 this->nodetype = get<1>(GetParam());
38 switch (this->nodetype)
41 this->niters = 2 * this->nnodes - 2;
42 this->nsteps = { 4, 8, 16, 32 };
46 this->niters = 2 * this->nnodes;
47 this->nsteps = { 4, 8, 16, 32 };
51 this->niters = this->nnodes;
52 this->nsteps = { 4, 8, 16, 32 };
56 this->niters = this->nnodes;
57 this->nsteps = { 4, 8, 16, 32 };
65 for (
size_t i = 0; i < this->nsteps.size(); ++i) {
66 auto dt = 0.5 / double(this->nsteps[i]);
81 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
82 MPI_Comm_size(MPI_COMM_WORLD, &size);
85 this->err.push_back(errors.at(last_error));
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])));
106 MPI_Comm_rank(MPI_COMM_WORLD, &rank);
107 MPI_Comm_size(MPI_COMM_WORLD, &size);
108 if (rank != size-1) {
113 switch (this->nodetype)
116 order = 2 * this->nnodes - 2;
117 quad =
"Gauss-Lobatto";
121 order = 2 * this->nnodes;
122 quad =
"Gauss-Legendre";
126 order = 2 * this->nnodes - 1;
127 quad =
"Gauss-Radau";
131 order = this->nnodes;
132 quad =
"Clenshaw-Curtis";
136 order = this->nnodes;
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.";
154 Combine(Range<size_t>(5, 6),
162 int main(
int argc,
char** argv)
164 testing::InitGoogleTest(&argc, argv);
165 MPI_Init(&argc, &argv);
169 int result = 1, max_result;
170 result = RUN_ALL_TESTS();
171 MPI_Allreduce(&result, &max_result, 1, MPI_INT, MPI_MAX, MPI_COMM_WORLD);
QuadratureType
Quadrature type descriptors.
int main(int argc, char **argv)
static void init(int argc, char **argv, std::function< void()> opts=nullptr, std::function< void()> logs=nullptr)
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)))
Clenshaw-Curtis quadrature
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.
Gauss-Legendre quadrature
TEST_P(ConvergenceTest, AllNodes)
advection-diffusion sweeper with semi-implicit time-integration.
error_map::key_type ktype