Skip to content
Snippets Groups Projects
Commit 980ce9de authored by Max Kahnt's avatar Max Kahnt
Browse files

Hold transfer operators with shared pointers instead of raw pointers.

parent c2444907
Branches
No related tags found
No related merge requests found
Pipeline #
......@@ -67,7 +67,7 @@ class ObstacleTNNMGStep
typedef typename Dune::BitSetVector<VectorType::block_type::dimension> BitVector;
typedef MultigridTransfer<Vector, BitVector, Matrix> Transfer;
typedef typename std::vector<Transfer*> TransferPointerVector;
typedef typename std::vector<std::shared_ptr<Transfer>> TransferPointerVector;
protected:
......
......@@ -55,10 +55,13 @@ struct MultigridTestSuite
typedef MultigridTransfer<Vector, BitVector, Matrix> Transfer;
typedef CompressedMultigridTransfer<Vector, BitVector, Matrix> TransferImplementation;
// create transfer operators from level i to i+1 (note that this will only work for either uniformly refined grids or adaptive grids with RefinementType=COPY)
std::vector<Transfer> transfer(grid.maxLevel());
for (size_t i = 0; i < transfer.size(); ++i) {
transfer[i].setup(grid, i, i+1);
std::vector<std::shared_ptr<Transfer>> transfer(grid.maxLevel());
for (size_t i = 0; i < transfer.size(); ++i)
{
// create transfer operators from level i to i+1 (note that this will only work for either uniformly refined grids or adaptive grids with RefinementType=COPY)
auto t = std::make_shared<TransferImplementation>();
t->setup(grid, i, i+1);
transfer[i] = t;
}
// set up smoothers and basesolver
......
......@@ -53,12 +53,11 @@ void solveObstacleProblemByTNNMG(const GridType& grid, const MatrixType& mat, Ve
const int blockSize = VectorType::block_type::dimension;
// we need a vector of pointers to the transfer operator base class
std::vector<Transfer*> transfer(grid.maxLevel());
std::vector<std::shared_ptr<Transfer>> transfer(grid.maxLevel());
for (size_t i = 0; i < transfer.size(); ++i)
{
// create transfer operator from level i to i+1
TransferImplementation* t = new TransferImplementation;
// create transfer operators from level i to i+1
auto t = std::make_shared<TransferImplementation>();
t->setup(grid, i, i+1);
transfer[i] = t;
}
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment