/* == 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);
    }
    /* == 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);
    }
    /* 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);
    }
Пример #4
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;
            }
        }
    }
}