diff --git a/Applications/FileIO/TetGenInterface.cpp b/Applications/FileIO/TetGenInterface.cpp
index 7726fa4518c1fb4632219260f414d9b863a96d1a..ef9ece69767adb70502546a416f9c1b9187ded6d 100644
--- a/Applications/FileIO/TetGenInterface.cpp
+++ b/Applications/FileIO/TetGenInterface.cpp
@@ -430,8 +430,8 @@ bool TetGenInterface::parseElements(std::ifstream& ins,
                                     bool region_attribute) const
 {
     std::string line;
-    auto* ids(static_cast<std::size_t*>(
-        alloca(sizeof(std::size_t) * n_nodes_per_tet)));
+    std::vector<std::size_t> ids(n_nodes_per_tet);
+
     elements.reserve(n_tets);
     materials.reserve(n_tets);
 
@@ -492,8 +492,9 @@ bool TetGenInterface::parseElements(std::ifstream& ins,
         }
         // insert new element into vector
         auto** tet_nodes = new MeshLib::Node*[4];
-        for (unsigned k(0); k<4; k++) {
-            tet_nodes[k] = nodes[ids[k]];
+        for (int n = 0; n < 4; n++)
+        {
+            tet_nodes[n] = nodes[ids[n]];
         }
         elements.push_back (new MeshLib::Tet(tet_nodes));
         materials.push_back(region);
diff --git a/GeoLib/Polyline.cpp b/GeoLib/Polyline.cpp
index 2c0525dc5066460d452c9b1db26f9f1e6b5617ad..1b47e65fbfc456f6d46b9b3e5bfd4ac790de5f02 100644
--- a/GeoLib/Polyline.cpp
+++ b/GeoLib/Polyline.cpp
@@ -347,6 +347,7 @@ Polyline* Polyline::constructPolylineFromSegments(const std::vector<Polyline*> &
         if (!ply_found)
         {
             ERR("Error in Polyline::contructPolylineFromSegments() - Not all segments are connected.");
+            delete new_ply;
             new_ply = nullptr;
             break;
         }
diff --git a/MeshLib/MeshEditing/MeshRevision.cpp b/MeshLib/MeshEditing/MeshRevision.cpp
index f53b67c7f431f64ecb111cc7a80eb15fdb6fbb14..8da790ff7c703cdbe5d0ed4a2bffaf1fc39aba85 100644
--- a/MeshLib/MeshEditing/MeshRevision.cpp
+++ b/MeshLib/MeshEditing/MeshRevision.cpp
@@ -487,6 +487,7 @@ unsigned MeshRevision::reduceHex(MeshLib::Element const*const org_elem,
                 prism_nodes[3] = nodes[org_elem->getNode(this->lutHexDiametralNode(org_elem->getNodeIDinElement(face->getNode(1))))->getID()];
                 prism_nodes[4] = nodes[org_elem->getNode(this->lutHexDiametralNode(org_elem->getNodeIDinElement(face->getNode(2))))->getID()];
                 prism_nodes[5] = nodes[org_elem->getNode(org_elem->getNodeIDinElement(face->getNode(0)))->getID()];
+                new_elements.push_back(new MeshLib::Prism(prism_nodes));
                 delete face;
                 return 1;
             }
diff --git a/MeshLib/MeshGenerators/MeshLayerMapper.cpp b/MeshLib/MeshGenerators/MeshLayerMapper.cpp
index d7c316e577d5dc79c8102ac35132681ca2698555..47500682008809ddad7c93b16c5407597dca314c 100644
--- a/MeshLib/MeshGenerators/MeshLayerMapper.cpp
+++ b/MeshLib/MeshGenerators/MeshLayerMapper.cpp
@@ -102,12 +102,20 @@ MeshLib::Mesh* MeshLayerMapper::createStaticLayers(MeshLib::Mesh const& mesh, st
                 e_nodes[j] = new_nodes[node_id+nNodes];
                 e_nodes[j+nElemNodes] = new_nodes[node_id];
             }
-            // extrude triangles to prism
             if (sfc_elem->getGeomType() == MeshLib::MeshElemType::TRIANGLE)
-                new_elems.push_back (new MeshLib::Prism(e_nodes));
-            // extrude quads to hexes
+            {
+                // extrude triangles to prism
+                new_elems.push_back(new MeshLib::Prism(e_nodes));
+            }
             else if (sfc_elem->getGeomType() == MeshLib::MeshElemType::QUAD)
-                new_elems.push_back (new MeshLib::Hex(e_nodes));
+            {
+                // extrude quads to hexes
+                new_elems.push_back(new MeshLib::Hex(e_nodes));
+            }
+            else
+            {
+                OGS_FATAL("MeshLayerMapper: Unknown element type to extrude.");
+            }
             materials->push_back(mat_id);
         }
     }
@@ -127,14 +135,13 @@ bool MeshLayerMapper::createRasterLayers(
         return false;
     }
 
-    auto* top(new MeshLib::Mesh(mesh));
+    auto top = std::make_unique<MeshLib::Mesh>(mesh);
     if (!layerMapping(*top, *rasters.back(), noDataReplacementValue))
         return false;
 
-    auto* bottom(new MeshLib::Mesh(mesh));
+    auto bottom = std::make_unique<MeshLib::Mesh>(mesh);
     if (!layerMapping(*bottom, *rasters[0], 0))
     {
-        delete top;
         return false;
     }
 
@@ -153,13 +160,11 @@ bool MeshLayerMapper::createRasterLayers(
     std::vector<MeshLib::Node*> const& nodes = bottom->getNodes();
     for (MeshLib::Node* node : nodes)
         _nodes.push_back(new MeshLib::Node(*node));
-    delete bottom;
 
     // add the other layers
     for (std::size_t i=0; i<nLayers-1; ++i)
         addLayerToMesh(*top, i, *rasters[i+1]);
 
-    delete top;
     return true;
 }
 
diff --git a/ProcessLib/BoundaryCondition/NonuniformNeumannBoundaryCondition.cpp b/ProcessLib/BoundaryCondition/NonuniformNeumannBoundaryCondition.cpp
index 6aad8fd1a3494c393baeb01998e69f1b08c75cae..96615b8de369a7949d4fb363d02ec82806a2b602 100644
--- a/ProcessLib/BoundaryCondition/NonuniformNeumannBoundaryCondition.cpp
+++ b/ProcessLib/BoundaryCondition/NonuniformNeumannBoundaryCondition.cpp
@@ -57,9 +57,10 @@ createNonuniformNeumannBoundaryCondition(
         boundary_mesh.getProperties().getPropertyVector<std::size_t>(
             mapping_to_bulk_nodes_property);
 
-    if (!(mapping_to_bulk_nodes && mapping_to_bulk_nodes->getMeshItemType() ==
-                                       MeshLib::MeshItemType::Node) &&
-        mapping_to_bulk_nodes->getNumberOfComponents() == 1)
+    if (mapping_to_bulk_nodes == nullptr ||
+        mapping_to_bulk_nodes->getMeshItemType() !=
+            MeshLib::MeshItemType::Node ||
+        mapping_to_bulk_nodes->getNumberOfComponents() != 1)
     {
         OGS_FATAL("Field `%s' is not set up properly.",
                   mapping_to_bulk_nodes_property.c_str());
diff --git a/Tests/GeoLib/TestGrid.cpp b/Tests/GeoLib/TestGrid.cpp
index 64b9c79c857195fc032edb97ec5fa223de7ab2e7..11b0fc144c130949bdb8ee8a49156b1725196dbe 100644
--- a/Tests/GeoLib/TestGrid.cpp
+++ b/Tests/GeoLib/TestGrid.cpp
@@ -82,6 +82,7 @@ TEST(GeoLib, SearchNearestPointInGrid)
     pnts.push_back(new GeoLib::Point(0.0,0.0,0.0));
     GeoLib::Grid<GeoLib::Point> *grid(nullptr);
     ASSERT_NO_THROW(grid = new GeoLib::Grid<GeoLib::Point>(pnts.begin(), pnts.end()));
+    ASSERT_NE(grid, nullptr);
 
     GeoLib::Point p0(0,10,10);
     GeoLib::Point* res(grid->getNearestPoint(p0));
diff --git a/Tests/GeoLib/TestPolygon.cpp b/Tests/GeoLib/TestPolygon.cpp
index 94327578f5b43fa9bbc981c1909d90d4cf34d874..886475657eb1e99afafb8cd353e4a24c516f6b48 100644
--- a/Tests/GeoLib/TestPolygon.cpp
+++ b/Tests/GeoLib/TestPolygon.cpp
@@ -85,14 +85,14 @@ TEST_F(PolygonTest, isPntInPolygonCheckCorners)
 TEST_F(PolygonTest, isPntInPolygonCheckPointsRestOnPolygonEdges)
 {
     for (std::size_t k(0); k<_pnts.size()-1; k++) {
-        for (double t(0.0); t<1.0; t+=0.001) {
-            EXPECT_TRUE(_polygon->isPntInPolygon(
-                GeoLib::Point(
-                    (*_pnts[k])[0] + t*((*_pnts[k+1])[0]-(*_pnts[k])[0]),
-                    (*_pnts[k])[1] + t*((*_pnts[k+1])[1]-(*_pnts[k])[1]),
-                    (*_pnts[k])[2] + t*((*_pnts[k+1])[2]-(*_pnts[k])[2])
-                ))
-            );
+        double t = 0;
+        while (t < 1.0)
+        {
+            EXPECT_TRUE(_polygon->isPntInPolygon(GeoLib::Point(
+                (*_pnts[k])[0] + t * ((*_pnts[k + 1])[0] - (*_pnts[k])[0]),
+                (*_pnts[k])[1] + t * ((*_pnts[k + 1])[1] - (*_pnts[k])[1]),
+                (*_pnts[k])[2] + t * ((*_pnts[k + 1])[2] - (*_pnts[k])[2]))));
+            t += 0.001;
         }
     }
 }
diff --git a/Tests/MeshLib/TestUniqueMeshId.cpp b/Tests/MeshLib/TestUniqueMeshId.cpp
index 3fa6d590f2a87eb71a48b3506e827632c5b83230..020d7f9c87eccf04178cc9b9230c663645d79d6d 100644
--- a/Tests/MeshLib/TestUniqueMeshId.cpp
+++ b/Tests/MeshLib/TestUniqueMeshId.cpp
@@ -23,13 +23,14 @@ TEST(MeshLib, UniqueMeshId)
     //
     // Test mesh counter increments.
     //
-    Mesh* m1 = new Mesh("second", std::vector<Node*>(), std::vector<Element*>());
+    auto m1 = std::make_unique<Mesh>("second", std::vector<Node*>(),
+                                     std::vector<Element*>());
     ASSERT_EQ(counter_value + std::size_t(1), m1->getID());
 
     Mesh m2("third", std::vector<Node*>(), std::vector<Element*>());
     ASSERT_EQ(counter_value + std::size_t(2), m2.getID());
 
-    delete m1;
+    m1.reset();
     ASSERT_EQ(counter_value + std::size_t(2), m2.getID());
 
     Mesh m3("fourth", std::vector<Node*>(), std::vector<Element*>());