From 297495d7e6df2f20271b60326a39aad925b2f43b Mon Sep 17 00:00:00 2001 From: Constantine Khroulev Date: Wed, 23 Sep 2015 13:02:24 -0800 Subject: [PATCH] Remove some raw pointers. --- src/base/stressbalance/sia/siafd_test.cc | 19 +++++++------------ src/base/stressbalance/ssa/SSAFEM.cc | 15 +++++---------- src/earth/greens.cc | 5 ++--- 3 files changed, 14 insertions(+), 25 deletions(-) diff --git a/src/base/stressbalance/sia/siafd_test.cc b/src/base/stressbalance/sia/siafd_test.cc index 183d3a5c18..7b79a2a059 100644 --- a/src/base/stressbalance/sia/siafd_test.cc +++ b/src/base/stressbalance/sia/siafd_test.cc @@ -205,22 +205,18 @@ static void setInitStateF(IceGrid &grid, IceModelVec2S *surface, IceModelVec2S *thickness, IceModelVec3 *enthalpy) { - int Mz=grid.Mz(); - double *dummy1, *dummy2, *dummy3, *dummy4, *dummy5; + int Mz = grid.Mz(); double ST = 1.67e-5, Tmin = 223.15, // K LforFG = 750000; // m - dummy1=new double[Mz]; dummy2=new double[Mz]; - dummy3=new double[Mz]; dummy4=new double[Mz]; - dummy5=new double[Mz]; + std::vector dummy1(Mz), dummy2(Mz), dummy3(Mz), dummy4(Mz), dummy5(Mz); bed->set(0); mask->set(MASK_GROUNDED); - double *T; - T = new double[grid.Mz()]; + std::vector T(Mz); IceModelVec::AccessList list; list.add(*thickness); @@ -231,6 +227,7 @@ static void setInitStateF(IceGrid &grid, double r = radius(grid, i, j), Ts = Tmin + ST * r; + if (r > LforFG - 1.0) { // if (essentially) outside of sheet (*thickness)(i, j) = 0.0; for (int k = 0; k < Mz; k++) { @@ -239,9 +236,9 @@ static void setInitStateF(IceGrid &grid, } else { r = std::max(r, 1.0); // avoid singularity at origin bothexact(0.0, r, &(grid.z()[0]), Mz, 0.0, - &(*thickness)(i, j), dummy5, T, dummy1, dummy2, dummy3, dummy4); + &(*thickness)(i, j), &dummy5[0], &T[0], &dummy1[0], &dummy2[0], &dummy3[0], &dummy4[0]); } - enthalpy->set_column(i, j, T); + enthalpy->set_column(i, j, &T[0]); } @@ -249,9 +246,6 @@ static void setInitStateF(IceGrid &grid, surface->copy_from(*thickness); - delete [] dummy1; delete [] dummy2; delete [] dummy3; delete [] dummy4; - delete [] T; delete [] dummy5; - enthalpy_from_temperature_cold(EC, grid, *thickness, *enthalpy, *enthalpy); @@ -343,6 +337,7 @@ int main(int argc, char *argv[]) { // create grid and set defaults IceGrid::Ptr grid(new IceGrid(ctx, P)); + grid->report_parameters(); setVerbosityLevel(5); diff --git a/src/base/stressbalance/ssa/SSAFEM.cc b/src/base/stressbalance/ssa/SSAFEM.cc index 462833ff06..0f2ac815b3 100644 --- a/src/base/stressbalance/ssa/SSAFEM.cc +++ b/src/base/stressbalance/ssa/SSAFEM.cc @@ -242,15 +242,13 @@ void SSAFEM::cacheQuadPtValues() { using fem::Quadrature; using fem::FunctionGerm; - const double - *Enth_e[4]; - double - *Enth_q[4]; + std::vector Enth_q[Quadrature::Nq]; + const double *Enth_e[4]; double ice_density = m_config->get_double("ice_density"); for (unsigned int q=0; qMz()]; + Enth_q[q].resize(m_grid->Mz()); } GeometryCalculator gc(sea_level, *m_config); @@ -343,7 +341,8 @@ void SSAFEM::cacheQuadPtValues() { // Evaluate column integrals in flow law at every quadrature point's column coefficients[q].B = m_flow_law->averaged_hardness(coefficients[q].H, m_grid->kBelowHeight(coefficients[q].H), - &(m_grid->z()[0]), Enth_q[q]); + &(m_grid->z()[0]), + &(Enth_q[q])[0]); } } // j-loop @@ -352,10 +351,6 @@ void SSAFEM::cacheQuadPtValues() { loop.failed(); } loop.check(); - - for (unsigned int q = 0; q < Quadrature::Nq; q++) { - delete [] Enth_q[q]; - } } /** @brief Compute the "(effective viscosity) x (ice thickness)" diff --git a/src/earth/greens.cc b/src/earth/greens.cc index 440a351703..24f4648f48 100644 --- a/src/earth/greens.cc +++ b/src/earth/greens.cc @@ -16,6 +16,7 @@ // along with PISM; if not, write to the Free Software // Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA +#include #include #include #include @@ -97,11 +98,10 @@ double viscDisc(double t, double H0, double R0, double r, const int N_gsl_workspace = 1000; gsl_integration_workspace* w = gsl_integration_workspace_alloc(N_gsl_workspace); - double* pts; const int lengthpts = 142; // Matlab: pts=[10.^(-3:-0.05:-10) 1.0e-14]; - pts = new double[lengthpts]; + std::vector pts(lengthpts); for (int j=0; j < lengthpts-1; j++) { pts[j] = pow(10.0,-3.0 - 0.05 * (double) j); } @@ -127,7 +127,6 @@ double viscDisc(double t, double H0, double R0, double r, sum += result; } - delete [] pts; gsl_integration_workspace_free(w); // u(k)=rhoi*g*H0*R0*result; return rho * grav * H0 * R0 * sum;