/* == Sliding boundary conditions ==
     *
     * It is common to require a Dirichlet boundary condition where the displacement/position in one dimension
     * is fixed, but the displacement/position in the others are free. This can be easily done when
     * collecting the new locations for the fixed nodes, as shown in the following example. Here, we
     * take a unit square, apply gravity downward, and suppose the Y=0 surface is like a frictionless boundary,
     * so that, for the nodes on Y=0, we specify y=0 but leave x free (Here (X,Y)=old position, (x,y)=new position).
     * (Note though that this wouldn't be enough to uniquely specify the final solution - an arbitrary
     * translation in the Y direction could be added a solution to obtain another valid solution, so we
     * fully fix the node at the origin.)
     */
    void TestWithSlidingDirichletBoundaryConditions() throw(Exception)
    {
        QuadraticMesh<2> mesh;
        mesh.ConstructRegularSlabMesh(0.1 /*stepsize*/, 1.0 /*width*/, 1.0 /*height*/);

        ExponentialMaterialLaw<2> law(1.0, 0.5); // First parameter is 'a', second 'b', in W=a*exp(b(I1-3))

        /* Create fixed nodes and locations... */
        std::vector<unsigned> fixed_nodes;
        std::vector<c_vector<double,2> > locations;

        /* Fix node 0 (the node at the origin) */
        fixed_nodes.push_back(0);
        locations.push_back(zero_vector<double>(2));

        /* For the rest, if the node is on the Y=0 surface.. */
        for (unsigned i=1; i<mesh.GetNumNodes(); i++)
        {
            if ( fabs(mesh.GetNode(i)->rGetLocation()[1])<1e-6)
            {
                /* ..add it to the list of fixed nodes.. */
                fixed_nodes.push_back(i);
                /* ..and define y to be 0 but x is fixed */
                c_vector<double,2> new_location;
                new_location(0) = SolidMechanicsProblemDefinition<2>::FREE;
                new_location(1) = 0.0;
                locations.push_back(new_location);
            }
        }

        /* Set the material law and fixed nodes, add some gravity, and solve */
        SolidMechanicsProblemDefinition<2> problem_defn(mesh);
        problem_defn.SetMaterialLaw(INCOMPRESSIBLE,&law);
        problem_defn.SetFixedNodes(fixed_nodes, locations);
        c_vector<double,2> gravity = zero_vector<double>(2);
        gravity(1) = -0.5;
        problem_defn.SetBodyForce(gravity);

        IncompressibleNonlinearElasticitySolver<2> solver(mesh,
                                                          problem_defn,
                                                          "ElasticitySlidingBcsExample");
        solver.Solve();
        solver.CreateCmguiOutput();

        /* Check the node at (1,0) has moved but has stayed on Y=0 */
        TS_ASSERT_LESS_THAN(1.0, solver.rGetDeformedPosition()[10](0));
        TS_ASSERT_DELTA(solver.rGetDeformedPosition()[10](1), 0.0, 1e-3);
    }
Exemplo n.º 2
0
    void TestBarPressureOnUndersideCompressible() throw(Exception)
    {
        EXIT_IF_PARALLEL; // petsc's direct solve only runs one 1 proc - MUMPS may be the answer for direct solves in parallel

        std::string base = "BarPressureOnUndersideCompressible";
        std::stringstream output_dir;

        /////////////////////////////////////
        // Solve on coarse mesh..
        /////////////////////////////////////
        QuadraticMesh<3> mesh;
        double h = 1.0;
        mesh.ConstructRegularSlabMesh(h, 10.0, 1.0, 1.0);
        output_dir << base << "h_" << h;
        std::vector<double> solution;
        unsigned num_iters = SolvePressureOnUndersideCompressible(mesh, output_dir.str(), solution, false, true);
        TS_ASSERT_EQUALS(num_iters, 5u);
    }
Exemplo n.º 3
0
    void TestReadMeshes(void) throw(Exception)
    {
        {
            READER_2D reader("mesh/test/data/square_4_elements_gmsh.msh");
            TetrahedralMesh<2,2> mesh;
            mesh.ConstructFromMeshReader(reader);
            TS_ASSERT_EQUALS(mesh.GetNumNodes(), 5u);
            TS_ASSERT_EQUALS(mesh.GetNumElements(), 4u);
            TS_ASSERT_EQUALS(mesh.GetNumBoundaryElements(), 4u);
        }

        {
            READER_3D reader("mesh/test/data/simple_cube_gmsh.msh");
            TetrahedralMesh<3,3> mesh;
            mesh.ConstructFromMeshReader(reader);
            TS_ASSERT_EQUALS(mesh.GetNumNodes(), 14u);
            TS_ASSERT_EQUALS(mesh.GetNumElements(), 24u);
            TS_ASSERT_EQUALS(mesh.GetNumBoundaryElements(), 24u);
        }

        {
           READER_2D reader("mesh/test/data/quad_square_4_elements_gmsh.msh",2,2);
           QuadraticMesh<2> mesh;
           mesh.ConstructFromMeshReader(reader);
           TS_ASSERT_EQUALS(mesh.GetNumNodes(), 13u);
           TS_ASSERT_EQUALS(mesh.GetNumElements(), 4u);
           TS_ASSERT_EQUALS(mesh.GetNumBoundaryElements(), 4u);
        }

        {
           READER_3D reader("mesh/test/data/quad_cube_gmsh.msh",2,2);
           QuadraticMesh<3> mesh;
           mesh.ConstructFromMeshReader(reader);
           TS_ASSERT_EQUALS(mesh.GetNumNodes(), 63u);
           TS_ASSERT_EQUALS(mesh.GetNumElements(), 24u);
           TS_ASSERT_EQUALS(mesh.GetNumBoundaryElements(), 24u);
        }
    }
Exemplo n.º 4
0
    void TestBarPressureOnUnderside() throw(Exception)
    {
        EXIT_IF_PARALLEL; // petsc's direct solve only runs one 1 proc - MUMPS may be the answer for direct solves in parallel

        std::string base = "BarPressureOnUnderside";
        std::stringstream output_dir;


        /////////////////////////////////////
        // Solve on coarse mesh..
        /////////////////////////////////////
        QuadraticMesh<3> mesh;
        double h = 1.0;
        mesh.ConstructRegularSlabMesh(h, 10.0, 1.0, 1.0);
        output_dir << base << "h_" << h;
        std::vector<double> solution;
        unsigned num_iters = SolvePressureOnUnderside(mesh, output_dir.str(), solution, false, true);
        TS_ASSERT_EQUALS(num_iters, 6u);

        /////////////////////////////////////
        // Solve on slightly finer mesh..
        /////////////////////////////////////
        h = 0.5;
        QuadraticMesh<3> mesh2;
        mesh2.ConstructRegularSlabMesh(h, 10.0, 1.0, 1.0);
        std::vector<double> solution2;

        InterpolateMechanicsSolutionToNewMesh<3>(mesh, solution, mesh2, solution2, INCOMPRESSIBLE);

        output_dir.str("");
        output_dir << base << "h_" << h;
        num_iters = SolvePressureOnUnderside(mesh2, output_dir.str(), solution2, true, true);
        TS_ASSERT_EQUALS(num_iters, 2u);

        ///////////////////////////////////////
        // Could continue onto finer meshes..
        ///////////////////////////////////////
    }
    /* == Incompressible deformation: non-zero displacement boundary conditions, functional tractions ==
     *
     * We now consider a more complicated example. We prescribe particular new locations for the nodes
     * on the Dirichlet boundary, and also show how to prescribe a traction that is given in functional form
     * rather than prescribed for each boundary element.
     */
    void TestIncompressibleProblemMoreComplicatedExample() throw(Exception)
    {
        /* Create a mesh */
        QuadraticMesh<2> mesh;
        mesh.ConstructRegularSlabMesh(0.1 /*stepsize*/, 1.0 /*width*/, 1.0 /*height*/);

        /* Use a different material law this time, an exponential material law.
         * The material law needs to inherit from `AbstractIncompressibleMaterialLaw`,
         * and there are a few implemented, see `continuum_mechanics/src/problem/material_laws` */
        ExponentialMaterialLaw<2> law(1.0, 0.5); // First parameter is 'a', second 'b', in W=a*exp(b(I1-3))
        /* Now specify the fixed nodes, and their new locations. Create `std::vector`s for each. */
        std::vector<unsigned> fixed_nodes;
        std::vector<c_vector<double,2> > locations;
        /* Loop over the mesh nodes */
        for (unsigned i=0; i<mesh.GetNumNodes(); i++)
        {
            /* If the node is on the Y=0 surface (the LHS) */
            if ( fabs(mesh.GetNode(i)->rGetLocation()[1])<1e-6)
            {
                /* Add it to the list of fixed nodes */
                fixed_nodes.push_back(i);
                /* and define a new position x=(X,0.1*X^2^) */
                c_vector<double,2> new_location;
                double X = mesh.GetNode(i)->rGetLocation()[0];
                new_location(0) = X;
                new_location(1) = 0.1*X*X;
                locations.push_back(new_location);
            }
        }

        /* Now collect all the boundary elements on the top surface, as before, except
         * here we don't create the tractions for each element
         */
        std::vector<BoundaryElement<1,2>*> boundary_elems;
        for (TetrahedralMesh<2,2>::BoundaryElementIterator iter = mesh.GetBoundaryElementIteratorBegin();
             iter != mesh.GetBoundaryElementIteratorEnd();
             ++iter)
        {
            /* If Y=1, have found a boundary element */
            if (fabs((*iter)->CalculateCentroid()[1] - 1.0)<1e-6)
            {
                BoundaryElement<1,2>* p_element = *iter;
                boundary_elems.push_back(p_element);
            }
        }

        /* Create a problem definition object, and this time calling `SetFixedNodes`
         * which takes in the new locations of the fixed nodes.
         */
        SolidMechanicsProblemDefinition<2> problem_defn(mesh);
        problem_defn.SetMaterialLaw(INCOMPRESSIBLE,&law);
        problem_defn.SetFixedNodes(fixed_nodes, locations);
        /* Now call `SetTractionBoundaryConditions`, which takes in a vector of
         * boundary elements as in the previous test. However this time the second argument
         * is a ''function pointer'' (just the name of the function) to a
         * function returning traction in terms of position (and time [see below]).
         * This function is defined above, before the tests. It has to take in a `c_vector` (position)
         *  and a double (time), and returns a `c_vector` (traction), and will only be called
         * using points in the boundary elements being passed in. The function `MyTraction`
         * (defined above, before the tests) above defines a horizontal traction (ie a shear stress, since it is
         * applied to the top surface) which increases in magnitude across the object.
          */
        problem_defn.SetTractionBoundaryConditions(boundary_elems, MyTraction);
        /* Note: You can also call `problem_defn.SetBodyForce(MyBodyForce)`, passing in a function
         * instead of a vector, although isn't really physically useful, it is only really useful
         * for constructing problems with exact solutions.
         *
         * Create the solver as before */
        IncompressibleNonlinearElasticitySolver<2> solver(mesh,
                                                          problem_defn,
                                                          "IncompressibleElasticityMoreComplicatedExample");


        /* Call `Solve()` */
        solver.Solve();

        /* Another quick check */
        TS_ASSERT_EQUALS(solver.GetNumNewtonIterations(), 6u);
        /* Visualise as before.
         *
         * '''Advanced:''' Note that the function `MyTraction` takes in time, which it didn't use. In the above it would have been called
         * with t=0. The current time can be set using `SetCurrentTime()`. The idea is that the user may want to solve a
         * sequence of static problems with time-dependent tractions (say), for which they should allow `MyTraction` to
         * depend on time, and put the solve inside a time-loop, for example:
         */
        //for (double t=0; t<T; t+=dt)
        //{
        //    solver.SetCurrentTime(t);
        //    solver.Solve();
        //}
        /* In this the current time would be passed through to `MyTraction`
         *
         * Create Cmgui output
         */
        solver.CreateCmguiOutput();

        /* This is just to check that nothing has been accidentally changed in this test */
        TS_ASSERT_DELTA(solver.rGetDeformedPosition()[98](0), 1.4543, 1e-3);
        TS_ASSERT_DELTA(solver.rGetDeformedPosition()[98](1), 0.5638, 1e-3);
    }
    /* == Compressible deformation, and other bits and pieces ==
     *
     * In this test, we will show the (very minor) changes required to solve a compressible nonlinear
     * elasticity problem, we will describe and show how to specify 'pressure on deformed body'
     * boundary conditions, we illustrate how a quadratic mesh can be generated using a linear mesh
     * input files, and we also illustrate how `Solve()` can be called repeatedly, with loading
     * changing between the solves.
     *
     * Note: for other examples of compressible solves, including problems with an exact solution, see the
     * file `continuum_mechanics/test/TestCompressibleNonlinearElasticitySolver.hpp`
     */
    void TestSolvingCompressibleProblem() throw (Exception)
    {
        /* All mechanics problems must take in quadratic meshes, but the mesh files for
         * (standard) linear meshes in Triangles/Tetgen can be automatically converted
         * to quadratic meshes, by simply doing the following. (The mesh loaded here is a disk
         * centred at the origin with radius 1).
         */
        QuadraticMesh<2> mesh;
        TrianglesMeshReader<2,2> reader("mesh/test/data/disk_522_elements");
        mesh.ConstructFromLinearMeshReader(reader);

        /* Compressible problems require a compressible material law, ie one that
         * inherits from `AbstractCompressibleMaterialLaw`. The `CompressibleMooneyRivlinMaterialLaw`
         * is one such example; instantiate one of these
         */
        CompressibleMooneyRivlinMaterialLaw<2> law(1.0, 0.5);

        /* For this problem, we fix the nodes on the surface for which Y < -0.9 */
        std::vector<unsigned> fixed_nodes;
        for ( TetrahedralMesh<2,2>::BoundaryNodeIterator iter = mesh.GetBoundaryNodeIteratorBegin();
              iter != mesh.GetBoundaryNodeIteratorEnd();
              ++iter)
        {
            double Y = (*iter)->rGetLocation()[1];
            if(Y < -0.9)
            {
                fixed_nodes.push_back((*iter)->GetIndex());
            }
        }

        /* We will (later) apply Neumann boundary conditions to surface elements which lie below Y=0,
         * and these are collected here. (Minor, subtle, comment: we don't bother here checking Y>-0.9,
         * so the surface elements collected here include the ones on the Dirichlet boundary. This doesn't
         * matter as the Dirichlet boundary conditions to the nonlinear system essentially overwrite
         * an Neumann-related effects).
         */
        std::vector<BoundaryElement<1,2>*> boundary_elems;
        for (TetrahedralMesh<2,2>::BoundaryElementIterator iter
              = mesh.GetBoundaryElementIteratorBegin();
            iter != mesh.GetBoundaryElementIteratorEnd();
            ++iter)
        {
           BoundaryElement<1,2>* p_element = *iter;
           if(p_element->CalculateCentroid()[1]<0.0)
           {
               boundary_elems.push_back(p_element);
           }
        }
        assert(boundary_elems.size()>0);

        /* Create the problem definition class, and set the law again, this time
         * stating that the law is compressible
         */
        SolidMechanicsProblemDefinition<2> problem_defn(mesh);
        problem_defn.SetMaterialLaw(COMPRESSIBLE,&law);

        /* Set the fixed nodes and gravity */
        problem_defn.SetZeroDisplacementNodes(fixed_nodes);

        /* The elasticity solvers have two nonlinear solvers implemented, one hand-coded and one which uses PETSc's SNES
         * solver. The latter is not the default but can be more robust (and will probably be the default in later
         * versions). This is how it can be used. (This option can also be called if the compiled binary is run from
         * the command line (see ChasteGuides/RunningBinariesFromCommandLine) using the option "-mech_use_snes").
         */
        problem_defn.SetSolveUsingSnes();
        /* This line tells the solver to output info about the nonlinear solve as it progresses, and can be used with
         * or without the SNES option above. The corresponding command line option is "-mech_verbose"
         */
        problem_defn.SetVerboseDuringSolve();

        c_vector<double,2> gravity;
        gravity(0) = 0;
        gravity(1) = 0.1;
        problem_defn.SetBodyForce(gravity);

        /* Declare the compressible solver, which has the same interface as the incompressible
         * one, and call `Solve()`
         */
        CompressibleNonlinearElasticitySolver<2> solver(mesh,
                                                        problem_defn,
                                                        "CompressibleSolidMechanicsExample");
        solver.Solve();

        /* Now we call add additional boundary conditions, and call `Solve() again. Firstly: these
         * Neumann conditions here are not specified traction boundary conditions (such BCs are specified
         * on the undeformed body), but instead, the (more natural) specification of a pressure
         * exactly in the ''normal direction on the deformed body''. We have to provide a set of boundary
         * elements of the mesh, and a pressure to act on those elements. The solver will automatically
         * compute the deformed normal directions on which the pressure acts. Note: with this type of
         * BC, the ordering of the nodes on the boundary elements needs to be consistent, otherwise some
         * normals will be computed to be inward and others outward. The solver will check this on the
         * original mesh and throw an exception if this is not the case. (The required ordering is such that:
         * in 2D, surface element nodes are ordered anticlockwise (looking at the whole mesh); in 3D, looking
         * at any surface element from outside the mesh, the three nodes are ordered anticlockwise. (Triangle/tetgen
         * automatically create boundary elements like this)).
         */
        double external_pressure = -0.4; // negative sign => inward pressure
        problem_defn.SetApplyNormalPressureOnDeformedSurface(boundary_elems, external_pressure);
        /* Call `Solve()` again, so now solving with fixed nodes, gravity, and pressure. The solution from
         * the previous solve will be used as the initial guess. Although at the moment the solution from the
         * previous call to `Solve()` will be over-written, calling `Solve()` repeatedly may be useful for
         * some problems: sometimes, Newton's method will fail to converge for given force/pressures etc, and it can
         * be (very) helpful to ''increment'' the loading. For example, set the gravity to be (0,-9.81/3), solve,
         * then set it to be (0,-2*9.81/3), solve again, and finally set it to be (0,-9.81) and solve for a
         * final time
         */
        solver.Solve();
        solver.CreateCmguiOutput();
    }
    void TestAssembler2d() throw (Exception)
    {
        QuadraticMesh<2> mesh;
        TrianglesMeshReader<2,2> mesh_reader("mesh/test/data/canonical_triangle_quadratic", 2, 2, false);
        mesh.ConstructFromMeshReader(mesh_reader);

        ContinuumMechanicsProblemDefinition<2> problem_defn(mesh);

        double t1 = 2.6854233;
        double t2 = 3.2574578;

        // for the boundary element on the y=0 surface, create a traction
        std::vector<BoundaryElement<1,2>*> boundary_elems;
        std::vector<c_vector<double,2> > tractions;
        c_vector<double,2> traction;
        traction(0) = t1;
        traction(1) = t2;
        for (TetrahedralMesh<2,2>::BoundaryElementIterator iter
              = mesh.GetBoundaryElementIteratorBegin();
             iter != mesh.GetBoundaryElementIteratorEnd();
             ++iter)
        {
            if (fabs((*iter)->CalculateCentroid()[1])<1e-4)
            {
                BoundaryElement<1,2>* p_element = *iter;
                boundary_elems.push_back(p_element);
                tractions.push_back(traction);
            }
        }
        assert(boundary_elems.size()==1);
        problem_defn.SetTractionBoundaryConditions(boundary_elems, tractions);


        ContinuumMechanicsNeumannBcsAssembler<2> assembler(&mesh, &problem_defn);

        TS_ASSERT_THROWS_THIS(assembler.AssembleVector(), "Vector to be assembled has not been set");
        Vec bad_sized_vec = PetscTools::CreateVec(2);
        assembler.SetVectorToAssemble(bad_sized_vec, true);
        TS_ASSERT_THROWS_THIS(assembler.AssembleVector(), "Vector provided to be assembled has size 2, not expected size of 18 ((dim+1)*num_nodes)");


        Vec vec = PetscTools::CreateVec(3*mesh.GetNumNodes());

        assembler.SetVectorToAssemble(vec, true);
        assembler.AssembleVector();

        ReplicatableVector vec_repl(vec);

        // Note: on a 1d boundary element,  intgl phi_i dx = 1/6 for the bases on the vertices
        // and intgl phi_i dx = 4/6 for the basis at the interior node. (Here the integrals
        // are over the canonical 1d element, [0,1], which is also the physical element for this
        // mesh.

        // node 0 is on the surface, and is a vertex
        TS_ASSERT_DELTA(vec_repl[0], t1/6.0, 1e-8);
        TS_ASSERT_DELTA(vec_repl[1], t2/6.0, 1e-8);

        // node 1 is on the surface, and is a vertex
        TS_ASSERT_DELTA(vec_repl[3], t1/6.0, 1e-8);
        TS_ASSERT_DELTA(vec_repl[4], t2/6.0, 1e-8);

        // nodes 2, 3, 4 are not on the surface
        for(unsigned i=2; i<5; i++)
        {
            TS_ASSERT_DELTA(vec_repl[3*i], 0.0, 1e-8);
            TS_ASSERT_DELTA(vec_repl[3*i], 0.0, 1e-8);
        }

        // node 5 is on the surface, and is an interior node
        TS_ASSERT_DELTA(vec_repl[15], 4*t1/6.0, 1e-8);
        TS_ASSERT_DELTA(vec_repl[16], 4*t2/6.0, 1e-8);


        // the rest of the vector is the pressure block and should be zero.
        TS_ASSERT_DELTA(vec_repl[2], 0.0, 1e-8);
        TS_ASSERT_DELTA(vec_repl[5], 0.0, 1e-8);
        TS_ASSERT_DELTA(vec_repl[8], 0.0, 1e-8);
        TS_ASSERT_DELTA(vec_repl[11], 0.0, 1e-8);
        TS_ASSERT_DELTA(vec_repl[14], 0.0, 1e-8);
        TS_ASSERT_DELTA(vec_repl[17], 0.0, 1e-8);

        PetscTools::Destroy(vec);
        PetscTools::Destroy(bad_sized_vec);
    }
Exemplo n.º 8
0
    /*
     * Test that the matrix is calculated correctly on the cannonical triangle.
     * Tests against the analytical solution calculated by hand.
     */
    void TestAssembler()  throw(Exception)
    {
        QuadraticMesh<2> mesh;
        TrianglesMeshReader<2,2> mesh_reader("mesh/test/data/canonical_triangle_quadratic", 2, 2, false);
        mesh.ConstructFromMeshReader(mesh_reader);

        double mu = 2.0;
        c_vector<double,2> body_force;
        double g1 = 1.34254;
        double g2 = 75.3422;
        body_force(0) = g1;
        body_force(1) = g2;

        StokesFlowProblemDefinition<2> problem_defn(mesh);
        problem_defn.SetViscosity(mu);
        problem_defn.SetBodyForce(body_force);

        StokesFlowAssembler<2> assembler(&mesh, &problem_defn);

        // The tests below test the assembler against hand-calculated variables for
        // an OLD weak form (corresponding to different boundary conditions), not the
        // current Stokes weak form. This factor converts the assembler to use the old
        // weak form. See documentation for this variable for more details.
        assembler.mScaleFactor = 0.0;

        Vec vec = PetscTools::CreateVec(18);
        Mat mat;
        PetscTools::SetupMat(mat, 18, 18, 18);

        assembler.SetVectorToAssemble(vec, true);
        assembler.SetMatrixToAssemble(mat, true);
        assembler.Assemble();
        PetscMatTools::Finalise(mat);

        double A[6][6] = {
            {      1.0,  1.0/6.0,  1.0/6.0,      0.0, -2.0/3.0, -2.0/3.0},
            {  1.0/6.0,  1.0/2.0,      0.0,      0.0,      0.0, -2.0/3.0},
            {  1.0/6.0,      0.0,  1.0/2.0,      0.0, -2.0/3.0,      0.0},
            {      0.0,      0.0,      0.0,  8.0/3.0, -4.0/3.0, -4.0/3.0},
            { -2.0/3.0,      0.0, -2.0/3.0, -4.0/3.0,  8.0/3.0,      0.0},
            { -2.0/3.0, -2.0/3.0,      0.0, -4.0/3.0,      0.0,  8.0/3.0}
        };

        double Bx[6][3] = {
            { -1.0/6.0,      0.0,      0.0},
            {      0.0,  1.0/6.0,      0.0},
            {      0.0,      0.0,      0.0},
            {  1.0/6.0,  1.0/6.0,  1.0/3.0},
            { -1.0/6.0, -1.0/6.0, -1.0/3.0},
            {  1.0/6.0, -1.0/6.0,      0.0},
        };

        double By[6][3] = {
            { -1.0/6.0,      0.0,      0.0},
            {      0.0,      0.0,      0.0},
            {      0.0,      0.0,  1.0/6.0},
            {  1.0/6.0,  1.0/3.0,  1.0/6.0},
            {  1.0/6.0,      0.0, -1.0/6.0},
            { -1.0/6.0, -1.0/3.0, -1.0/6.0},
        };

        c_matrix<double,18,18> exact_A = zero_matrix<double>(18);

        // The diagonal 6x6 blocks
        for (unsigned i=0; i<6; i++)
        {
            for (unsigned j=0; j<6; j++)
            {
                exact_A(3*i,  3*j)   = mu*A[i][j];
                exact_A(3*i+1,3*j+1) = mu*A[i][j];
            }
        }


        // The 6x3 Blocks
        for (unsigned i=0; i<6; i++)
        {
            for (unsigned j=0; j<3; j++)
            {
                exact_A(3*i,3*j+2)   = -Bx[i][j];
                exact_A(3*i+1,3*j+2) = -By[i][j];
                //- as -Div(U)=0
                exact_A(3*j+2,3*i)   = -Bx[i][j];
                exact_A(3*j+2,3*i+1) = -By[i][j];
            }
        }

        int lo, hi;
        MatGetOwnershipRange(mat, &lo, &hi);
        for (unsigned i=lo; i<(unsigned)hi; i++)
        {
            for (unsigned j=0; j<18; j++)
            {
                TS_ASSERT_DELTA(PetscMatTools::GetElement(mat,i,j), exact_A(i,j), 1e-9);
            }
        }
        ReplicatableVector vec_repl(vec);

        // The first 6 entries in the vector correspond to nodes 0, 1, 2, i.e. the vertices.
        // For these nodes, it can be shown that the integral of the corresponding
        // basis function is zero, i.e. \intgl_{canonical element} \phi_i dV = 0.0  for i=0,1,2, phi_i the
        // i-th QUADRATIC basis.
        for(unsigned i=0; i<3; i++)
        {
            TS_ASSERT_DELTA(vec_repl[3*i],   g1*0.0, 1e-8);
            TS_ASSERT_DELTA(vec_repl[3*i+1], g2*0.0, 1e-8);
        }

        // The next 6 entries in the vector correspond to nodes 3, 4, 5, i.e. the internal edges.
        // For these nodes, it can be shown that the integral of the corresponding
        // basis function is 1/6, i.e. \intgl_{canonical element} \phi_i dV = 1/6  for i=3,4,5, phi_i the
        // i-th QUADRATIC basis.
        for(unsigned i=3; i<6; i++)
        {
            TS_ASSERT_DELTA(vec_repl[3*i],   g1/6.0, 1e-8);
            TS_ASSERT_DELTA(vec_repl[3*i+1], g2/6.0, 1e-8);
        }

        // The pressure-block of the RHS vector should be zero.
        TS_ASSERT_DELTA(vec_repl[2], 0.0, 1e-9);
        TS_ASSERT_DELTA(vec_repl[5], 0.0, 1e-9);
        TS_ASSERT_DELTA(vec_repl[8], 0.0, 1e-9);
        TS_ASSERT_DELTA(vec_repl[11], 0.0, 1e-9);
        TS_ASSERT_DELTA(vec_repl[14], 0.0, 1e-9);
        TS_ASSERT_DELTA(vec_repl[17], 0.0, 1e-9);

        // Replace the body force with a functional body force (see MyBodyForce) above, and
        // assemble the vector again. This bit isn't so much to test the vector, but
        // to test the physical location being integrated at is interpolated correctly
        // and passed into the force function - see asserts in MyBodyForce.
        mesh.Translate(0.5, 0.8);

        Vec vec2 = PetscTools::CreateVec(18);
        assembler.SetVectorToAssemble(vec2, true);
        problem_defn.SetBodyForce(MyBodyForce);
        assembler.Assemble();

        ReplicatableVector vec_repl2(vec2);
        for(unsigned i=3; i<6; i++)
        {
            TS_ASSERT_DELTA(vec_repl2[3*i],   10.0/6.0, 1e-8);
            TS_ASSERT_DELTA(vec_repl2[3*i+1], 20.0/6.0, 1e-8);
        }

        PetscTools::Destroy(vec);
        PetscTools::Destroy(vec2);
        PetscTools::Destroy(mat);
    }
    /* == Incompressible deformation: 2D shape hanging under gravity with a balancing traction ==
     *
     * We now repeat the above test but include a traction on the bottom surface (Y=0). We apply this
     * in the inward direction so that is counters (somewhat) the effect of gravity. We also show how stresses
     * and strains can be written to file.
     */
    void TestIncompressibleProblemWithTractions() throw(Exception)
    {
        /* All of this is exactly as above */
        QuadraticMesh<2> mesh;
        mesh.ConstructRegularSlabMesh(0.1 /*stepsize*/, 0.8 /*width*/, 1.0 /*height*/);

        MooneyRivlinMaterialLaw<2> law(1.0);

        c_vector<double,2> body_force;
        body_force(0) =  0.0;
        body_force(1) = -2.0;

        std::vector<unsigned> fixed_nodes = NonlinearElasticityTools<2>::GetNodesByComponentValue(mesh, 1, 1.0);

        /* Now the traction boundary conditions. We need to collect all the boundary elements on the surface which we want to
         * apply non-zero tractions, put them in a `std::vector`, and create a corresponding `std::vector` of the tractions
         * for each of the boundary elements. Note that the each traction is a 2D vector with dimensions of pressure.
         *
         * First, declare the data structures:
         */
        std::vector<BoundaryElement<1,2>*> boundary_elems;
        std::vector<c_vector<double,2> > tractions;
        /* Create a constant traction */
        c_vector<double,2> traction;
        traction(0) = 0;
        traction(1) = 1.0; // this choice of sign corresponds to an inward force (if applied to the bottom surface)
        /* Loop over boundary elements */
        for (TetrahedralMesh<2,2>::BoundaryElementIterator iter = mesh.GetBoundaryElementIteratorBegin();
             iter != mesh.GetBoundaryElementIteratorEnd();
             ++iter)
        {
            /* If the centre of the element has Y value of 0.0, it is on the surface we need */
            if (fabs((*iter)->CalculateCentroid()[1] - 0.0) < 1e-6)
            {
                /* Put the boundary element and the constant traction into the stores. */
                BoundaryElement<1,2>* p_element = *iter;
                boundary_elems.push_back(p_element);
                tractions.push_back(traction);
            }
        }
        /* A quick check */
        assert(boundary_elems.size() == 8u);

        /* Now create the problem definition object, setting the material law, fixed nodes and body force as
         * before (this time not calling `SetDensity()`, so using the default density of 1.0,
         * and also calling a method for setting tractions, which takes in the boundary elements
         * and tractions for each of those elements.
         */
        SolidMechanicsProblemDefinition<2> problem_defn(mesh);
        problem_defn.SetMaterialLaw(INCOMPRESSIBLE,&law);
        problem_defn.SetZeroDisplacementNodes(fixed_nodes);
        problem_defn.SetBodyForce(body_force);
        problem_defn.SetTractionBoundaryConditions(boundary_elems, tractions);

        /* Create solver as before */
        IncompressibleNonlinearElasticitySolver<2> solver(mesh,
                                                          problem_defn,
                                                          "IncompressibleElasticityWithTractionsTutorial");

        /* In this test we also output the stress and strain. For the former, we have to tell the solver to store
         * the stresses that are computed during the solve.
         */
        solver.SetComputeAverageStressPerElementDuringSolve();


        /* Call `Solve()` */
        solver.Solve();

        /* If VTK output is written (discussed above) strains can be visualised. Alternatively, we can create text files
         * for strains and stresses by doing the following.
         *
         * Write the final deformation gradients to file. The i-th line of this file provides the deformation gradient F,
         * written as 'F(0,0) F(0,1) F(1,0) F(1,1)', evaluated at the centroid of the i-th element. The first variable
         * can also be DEFORMATION_TENSOR_C or LAGRANGE_STRAIN_E to write C or E. The second parameter is the file name.
         */
        solver.WriteCurrentStrains(DEFORMATION_GRADIENT_F,"deformation_grad");

        /* Since we called `SetComputeAverageStressPerElementDuringSolve`, we can write the stresses to file too. However,
         * note that for each element this is not the stress evaluated at the centroid, but the mean average of the stresses
         * evaluated at the quadrature points - for technical cardiac electromechanics reasons, it is difficult to
         * define the stress at non-quadrature points.
         */
        solver.WriteCurrentAverageElementStresses("2nd_PK_stress");


        /* Another quick check */
        TS_ASSERT_EQUALS(solver.GetNumNewtonIterations(), 4u);

        /* Visualise as before by going to the output directory and doing
         * `x=load('solution.nodes'); plot(x(:,1),x(:,2),'m*')` in Matlab/octave, or by using Cmgui.
         * The effect of the traction should be clear (especially when compared to
         * the results of the first test).
         *
         * Create Cmgui output
         */
        solver.CreateCmguiOutput();

        /* This is just to check that nothing has been accidentally changed in this test */
        TS_ASSERT_DELTA(solver.rGetDeformedPosition()[8](0), 0.8561, 1e-3);
        TS_ASSERT_DELTA(solver.rGetDeformedPosition()[8](1), 0.0310, 1e-3);
    }
    /* In the first test we use INCOMPRESSIBLE nonlinear elasticity. For incompressible elasticity, there
     * is a constraint on the deformation, which results in a pressure field (a Lagrange multiplier)
     * which is solved for together with the deformation.
     *
     * All the mechanics solvers solve for the deformation using the finite element method with QUADRATIC
     * basis functions for the deformation. This necessitates the use of a `QuadraticMesh` - such meshes have
     * extra nodes that aren't vertices of elements, in this case midway along each edge. (The displacement
     * is solved for at ''each node'' in the mesh (including internal [non-vertex] nodes), whereas the pressure
     * is only solved for at each vertex - in FEM terms, quadratic interpolation for displacement, linear
     * interpolation for pressure, which is required for stability. The pressure at internal nodes is computed
     * by linear interpolation).
     *
     */
    void TestSimpleIncompressibleProblem() throw(Exception)
    {
        /* First, define the geometry. This should be specified using the `QuadraticMesh` class, which inherits from `TetrahedralMesh`
         * and has mostly the same interface. Here we define a 0.8 by 1 rectangle, with elements 0.1 wide.
         * (`QuadraticMesh`s can also be read in using `TrianglesMeshReader`; see next tutorial/rest of code base for examples of this).
         */
        QuadraticMesh<2> mesh;
        mesh.ConstructRegularSlabMesh(0.1 /*stepsize*/, 0.8 /*width*/, 1.0 /*height*/);

        /* We use a Mooney-Rivlin material law, which applies to isotropic materials and has two parameters.
         * Restricted to 2D however, it only has one parameter, which can be thought of as the total
         * stiffness. We declare a Mooney-Rivlin law, setting the parameter to 1.
         */
        MooneyRivlinMaterialLaw<2> law(1.0);

        /* Next, the body force density. In realistic problems this will either be
         * acceleration due to gravity (ie b=(0,-9.81)) or zero if the effect of gravity can be neglected.
         * In this problem we apply a gravity-like downward force.
         */
        c_vector<double,2> body_force;
        body_force(0) =  0.0;
        body_force(1) = -2.0;


        /* Two types of boundary condition are required: displacement and traction. As with the other PDE solvers,
         * the displacement (Dirichlet) boundary conditions are specified at nodes, whereas traction (Neumann) boundary
         * conditions are specified on boundary elements.
         *
         * In this test we apply displacement boundary conditions on one surface of the mesh, the upper (Y=1.0) surface.
         * We are going to specify zero-displacement at these nodes.
         * We do not specify any traction boundary conditions, which means that (effectively) zero-traction boundary
         * conditions (ie zero pressures) are applied on the three other surfaces.
         *
         * We need to get a `std::vector` of all the node indices that we want to fix. The `NonlinearElasticityTools`
         * has a static method for helping do this: the following gets all the nodes for which Y=1.0. The second
         * argument (the '1') indicates Y . (So, for example, `GetNodesByComponentValue(mesh, 0, 10)` would get the nodes on X=10).
         */
        std::vector<unsigned> fixed_nodes = NonlinearElasticityTools<2>::GetNodesByComponentValue(mesh, 1, 1.0);

        /*
         * Before creating the solver we create a `SolidMechanicsProblemDefinition` object,  which contains
         * everything that defines the problem: mesh, material law, body force,
         * the fixed nodes and their locations, any traction boundary conditions, and the density
         * (which multiplies the body force, otherwise isn't used).
         */
        SolidMechanicsProblemDefinition<2> problem_defn(mesh);

        /* Set the material problem on the problem definition object, saying that the problem, and
         * the material law, is incompressible. All material law files can be found in
         * `continuum_mechanics/src/problem/material_laws`. */
        problem_defn.SetMaterialLaw(INCOMPRESSIBLE,&law);

        /* Set the fixed nodes, choosing zero displacement for these nodes (see later for how
         * to provide locations for the fixed nodes). */
        problem_defn.SetZeroDisplacementNodes(fixed_nodes);
        /* Set the body force and the density. (Note that the second line isn't technically
         * needed, as internally the density is initialised to 1)
         */
        problem_defn.SetBodyForce(body_force);
        problem_defn.SetDensity(1.0);

        /* Now we create the (incompressible) solver, passing in the mesh, problem definition
         * and output directory
         */
        IncompressibleNonlinearElasticitySolver<2> solver(mesh,
                                                          problem_defn,
                                                          "SimpleIncompressibleElasticityTutorial");

        /* .. and to compute the solution, just call `Solve()` */
        solver.Solve();


        /* '''Visualisation'''. Go to the folder `SimpleIncompressibleElasticityTutorial` in your test-output directory.
         * There should be 2 files, initial.nodes and solution.nodes. These are the original nodal positions and the deformed
         * positions. Each file has two columns, the x and y locations of each node. To visualise the solution in say
         * Matlab or Octave, you could do: `x=load('solution.nodes'); plot(x(:,1),x(:,2),'k*')`. For Cmgui output, see below.
         *
         * To get the actual solution from the solver, use these two methods. Note that the first
         * gets the deformed position (ie the new location, not the displacement). They are both of size
         * num_total_nodes.
         */
        std::vector<c_vector<double,2> >& r_deformed_positions = solver.rGetDeformedPosition();
        std::vector<double>& r_pressures = solver.rGetPressures();
        /* Let us obtain the values of the new position, and the pressure, at the bottom right corner node. */
        unsigned node_index = 8;
        assert( fabs(mesh.GetNode(node_index)->rGetLocation()[0] - 0.8) < 1e-6); // check that X=0.8, ie that we have the correct node,
        assert( fabs(mesh.GetNode(node_index)->rGetLocation()[1] - 0.0) < 1e-6); // check that Y=0.0, ie that we have the correct node,
        std::cout << "New position: " << r_deformed_positions[node_index](0) << " " << r_deformed_positions[node_index](1) << "\n";
        std::cout << "Pressure: " << r_pressures[node_index] << "\n";

        /* HOW_TO_TAG Continuum mechanics
         * Visualise nonlinear elasticity problems solutions, including visualisng strains
         */

        /* One visualiser is Cmgui. This method can be used to convert all the output files to Cmgui format.
         * They are placed in `[OUTPUT_DIRECTORY]/cmgui`. A script is created to easily load the data: in a
         * terminal cd to this directory and call `cmgui LoadSolutions.com`. (In this directory, the initial position is given by
         * solution_0.exnode, the deformed by solution_1.exnode).
         */
        solver.CreateCmguiOutput();

        /* The recommended visualiser is Paraview, for which Chaste must be installed with VTK. With paraview, strains (and in the future
         * stresses) can be visualised on the undeformed/deformed geometry). We can create VTK output using
         * the `VtkNonlinearElasticitySolutionWriter` class. The undeformed geometry, solution displacement, and pressure (if incompressible
         * problem) are written to file, and below we also choose to write the deformation tensor C for each element.
         */
        VtkNonlinearElasticitySolutionWriter<2> vtk_writer(solver);
        vtk_writer.SetWriteElementWiseStrains(DEFORMATION_TENSOR_C); // other options are DEFORMATION_GRADIENT_F and LAGRANGE_STRAIN_E
        vtk_writer.Write();


        /* These are just to check that nothing has been accidentally changed in this test.
         * Newton's method (with damping) was used to solve the nonlinear problem, and we check that
         * 4 iterations were needed to converge.
         */
        TS_ASSERT_DELTA(r_deformed_positions[node_index](0),  0.7980, 1e-3);
        TS_ASSERT_DELTA(r_deformed_positions[node_index](1), -0.1129, 1e-3);
        TS_ASSERT_EQUALS(solver.GetNumNewtonIterations(), 4u);
    }
VoltageInterpolaterOntoMechanicsMesh<DIM>::VoltageInterpolaterOntoMechanicsMesh(
                                     TetrahedralMesh<DIM,DIM>& rElectricsMesh,
                                     QuadraticMesh<DIM>& rMechanicsMesh,
                                     std::vector<std::string>& rVariableNames,
                                     std::string directory,
                                     std::string inputFileNamePrefix)
{
    // Read the data from the HDF5 file
    Hdf5DataReader reader(directory,inputFileNamePrefix);

    unsigned num_timesteps = reader.GetUnlimitedDimensionValues().size();

    // set up the elements and weights for the coarse nodes in the fine mesh
    FineCoarseMeshPair<DIM> mesh_pair(rElectricsMesh, rMechanicsMesh);
    mesh_pair.SetUpBoxesOnFineMesh();
    mesh_pair.ComputeFineElementsAndWeightsForCoarseNodes(true);
    assert(mesh_pair.rGetElementsAndWeights().size()==rMechanicsMesh.GetNumNodes());

    // create and setup a writer
    Hdf5DataWriter* p_writer = new Hdf5DataWriter(*rMechanicsMesh.GetDistributedVectorFactory(),
                                                  directory,
                                                  "voltage_mechanics_mesh",
                                                  false, //don't clean
                                                  false);

    std::vector<int> columns_id;
    for (unsigned var_index = 0; var_index < rVariableNames.size(); var_index++)
    {
        std::string var_name = rVariableNames[var_index];
        columns_id.push_back( p_writer->DefineVariable(var_name,"mV") );
    }

    p_writer->DefineUnlimitedDimension("Time","msecs", num_timesteps);
    p_writer->DefineFixedDimension( rMechanicsMesh.GetNumNodes() );
    p_writer->EndDefineMode();

    assert(columns_id.size() == rVariableNames.size());

    // set up a vector to read into
    DistributedVectorFactory factory(rElectricsMesh.GetNumNodes());
    Vec voltage = factory.CreateVec();
    std::vector<double> interpolated_voltages(rMechanicsMesh.GetNumNodes());
    Vec voltage_coarse = NULL;

    for(unsigned time_step=0; time_step<num_timesteps; time_step++)
    {
        for (unsigned var_index = 0; var_index < rVariableNames.size(); var_index++)
        {
            std::string var_name = rVariableNames[var_index];
            // read
            reader.GetVariableOverNodes(voltage, var_name, time_step);
            ReplicatableVector voltage_repl(voltage);

            // interpolate
            for(unsigned i=0; i<mesh_pair.rGetElementsAndWeights().size(); i++)
            {
                double interpolated_voltage = 0;

                Element<DIM,DIM>& element = *(rElectricsMesh.GetElement(mesh_pair.rGetElementsAndWeights()[i].ElementNum));
                for(unsigned node_index = 0; node_index<element.GetNumNodes(); node_index++)
                {
                    unsigned global_node_index = element.GetNodeGlobalIndex(node_index);
                    interpolated_voltage += voltage_repl[global_node_index]*mesh_pair.rGetElementsAndWeights()[i].Weights(node_index);
                }

                interpolated_voltages[i] = interpolated_voltage;
            }

            if(voltage_coarse!=NULL)
            {
                PetscTools::Destroy(voltage_coarse);
            }
            voltage_coarse = PetscTools::CreateVec(interpolated_voltages);
            // write
            p_writer->PutVector(columns_id[var_index], voltage_coarse);
        }
        p_writer->PutUnlimitedVariable(time_step);
        p_writer->AdvanceAlongUnlimitedDimension();
    }

    if(voltage_coarse!=NULL)
    {
        PetscTools::Destroy(voltage);
        PetscTools::Destroy(voltage_coarse);
    }

    // delete to flush
    delete p_writer;

    // Convert the new data to CMGUI format.
    // alter the directory in HeartConfig as that is where Hdf5ToCmguiConverter decides
    // where to output
    std::string config_directory = HeartConfig::Instance()->GetOutputDirectory();
    HeartConfig::Instance()->SetOutputDirectory(directory);
    Hdf5ToCmguiConverter<DIM,DIM> converter(FileFinder(directory, RelativeTo::ChasteTestOutput),
                                            "voltage_mechanics_mesh",
                                            &rMechanicsMesh,
                                            false);
    HeartConfig::Instance()->SetOutputDirectory(config_directory);
}
Exemplo n.º 12
0
void InterpolateMechanicsSolutionToNewMesh(QuadraticMesh<3>& rCoarseMesh, std::vector<double>& rCoarseSolution,
                                           QuadraticMesh<3>& rFineMesh, std::vector<double>& rFineSolution,
                                           CompressibilityType compressibilityType)
{
    unsigned NUM_UNKNOWNS = (compressibilityType==INCOMPRESSIBLE ? DIM+1 : DIM);

    if(rCoarseSolution.size() != rCoarseMesh.GetNumNodes()*NUM_UNKNOWNS)
    {
        EXCEPTION("rCoarseSolution not correct size");
    }
    if(rFineSolution.size() != rFineMesh.GetNumNodes()*NUM_UNKNOWNS)
    {
        rFineSolution.resize(rFineMesh.GetNumNodes()*NUM_UNKNOWNS);
    }

    c_vector<double, (DIM+1)*(DIM+2)/2> quad_basis;

    for(unsigned i=0; i<rFineMesh.GetNumNodes(); i++)
    {
        // find containing elements and weights in coarse mesh
        ChastePoint<DIM> point = rFineMesh.GetNode(i)->GetPoint();
        unsigned coarse_element_index = rCoarseMesh.GetContainingElementIndex(point,
                                                                              false);
        Element<DIM,DIM>* p_coarse_element = rCoarseMesh.GetElement(coarse_element_index);
        c_vector<double,DIM+1> weight = p_coarse_element->CalculateInterpolationWeights(point);

        c_vector<double,DIM> xi;
        xi(0) = weight(1);
        xi(1) = weight(2);
        if(DIM==3)
        {
            xi(2) = weight(3);
        }

        QuadraticBasisFunction<DIM>::ComputeBasisFunctions(xi, quad_basis);

        // interpolate (u,p) (don't do anything for p if compressible)
        c_vector<double,DIM+1> fine_solution = zero_vector<double>(DIM+1);

        for(unsigned elem_node_index=0; elem_node_index<(DIM+1)*(DIM+2)/2; elem_node_index++)
        {
            unsigned coarse_node = p_coarse_element->GetNodeGlobalIndex(elem_node_index);
            c_vector<double,DIM+1> coarse_solution_at_node;

            for(unsigned j=0; j<DIM; j++)
            {
                coarse_solution_at_node(j) = rCoarseSolution[NUM_UNKNOWNS*coarse_node + j];
            }
            if(compressibilityType==INCOMPRESSIBLE)
            {
                coarse_solution_at_node(DIM) = rCoarseSolution[NUM_UNKNOWNS*coarse_node + DIM];
            }

            fine_solution += coarse_solution_at_node*quad_basis(elem_node_index);
        }

        for(unsigned j=0; j<DIM; j++)
        {
            rFineSolution[NUM_UNKNOWNS*i + j] = fine_solution(j);
        }

        if(compressibilityType==INCOMPRESSIBLE)
        {
            rFineSolution[NUM_UNKNOWNS*i + DIM] = fine_solution(DIM);

            // Whilst the returned p from a solve is defined properly at all nodes, during the solve linear basis functions are
            // used for p and therefore p not computed explicitly at internal nodes, and the solver solves for p=0 at these internal
            // nodes. (After the solve, p is interpolated from vertices to internal nodes)
            if(rFineMesh.GetNode(i)->IsInternal())
            {
                rFineSolution[NUM_UNKNOWNS*i + DIM] = 0.0;
            }
        }
    }
}
Exemplo n.º 13
0
// Similar to above but uses a compressible material
unsigned SolvePressureOnUndersideCompressible(QuadraticMesh<3>& rMesh, std::string outputDirectory,
                                              std::vector<double>& rSolution,
                                              bool useSolutionAsGuess,
                                              double scaleFactor = 1.0)
{
    CompressibleExponentialLaw<3> law;
    std::vector<unsigned> fixed_nodes = NonlinearElasticityTools<3>::GetNodesByComponentValue(rMesh, 0, 0.0);

    std::vector<BoundaryElement<2,3>*> boundary_elems;

    double pressure = 0.001;

    for (TetrahedralMesh<3,3>::BoundaryElementIterator iter
           = rMesh.GetBoundaryElementIteratorBegin();
         iter != rMesh.GetBoundaryElementIteratorEnd();
         ++iter)
    {
        BoundaryElement<2,3>* p_element = *iter;
        double Z = p_element->CalculateCentroid()[2];
        if(fabs(Z)<1e-6)
        {
            boundary_elems.push_back(p_element);
        }
    }

    SolidMechanicsProblemDefinition<3> problem_defn(rMesh);
    problem_defn.SetMaterialLaw(COMPRESSIBLE,&law);
    problem_defn.SetZeroDisplacementNodes(fixed_nodes);

    problem_defn.SetApplyNormalPressureOnDeformedSurface(boundary_elems, pressure*scaleFactor);


    problem_defn.SetVerboseDuringSolve();

    CompressibleNonlinearElasticitySolver<3> solver(rMesh,problem_defn,outputDirectory);
    solver.SetWriteOutputEachNewtonIteration();

    // cover the SetTakeFullFirstNewtonStep() method, and the SetUsePetscDirectSolve() method
    solver.SetTakeFullFirstNewtonStep();
    solver.SetUsePetscDirectSolve();

    if(useSolutionAsGuess)
    {
        if(solver.rGetCurrentSolution().size()!=rSolution.size())
        {
            EXCEPTION("Badly-sized input");
        }
        for(unsigned i=0; i<rSolution.size(); i++)
        {
            solver.rGetCurrentSolution()[i] = rSolution[i];
        }
    }

    if(scaleFactor < 1.0)
    {
        try
        {
            solver.Solve();
        }
        catch (Exception& e)
        {
            // not final Solve, so don't quit
            WARNING(e.GetMessage());
        }
    }
    else
    {
        solver.Solve();
        solver.CreateCmguiOutput();
        
        VtkNonlinearElasticitySolutionWriter<3> vtk_writer(solver);
        vtk_writer.SetWriteElementWiseStrains(DEFORMATION_TENSOR_C);
        vtk_writer.Write();        
    }

    rSolution.clear();
    rSolution.resize(solver.rGetCurrentSolution().size());
    for(unsigned i=0; i<rSolution.size(); i++)
    {
        rSolution[i] = solver.rGetCurrentSolution()[i];
    }

    return solver.GetNumNewtonIterations();
}