c_vector<double,2> AnotherFunction(c_vector<double,2>& rX, double t) { c_vector<double,2> body_force; body_force(0) = rX(0)*t; body_force(1) = 10*rX(1)*t; return body_force; }
c_vector<double,2> SomeFunction(c_vector<double,2>& rX, double t) { c_vector<double,2> body_force; body_force(0) = rX(0)+t; body_force(1) = 2*(rX(1)+t); return body_force; }
c_vector<double,2> MyBodyForce(c_vector<double,2>& rX, double t) { // check the point rX has been interpolated correctly. The mesh // is the canonical triangle translated by (0.5,0.8). assert(rX(0)>0.0 + 0.5); assert(rX(1)>0.0 + 0.8); assert(rX(0)+rX(1)<1.0 + 0.5 + 0.8); c_vector<double,2> body_force; body_force(0) = 10.0; body_force(1) = 20.0; return body_force; }
static c_vector<double,3> GetBodyForce(c_vector<double,3>& rX, double t) { assert(rX(0)>=0 && rX(0)<=1 && rX(1)>=0 && rX(1)<=1 && rX(2)>=0 && rX(2)<=1); double lam1 = 1+a*rX(0); double lam2 = 1+b*rX(1); double invlam1 = 1.0/lam1; double invlam2 = 1.0/lam2; c_vector<double,3> body_force; body_force(0) = a; body_force(1) = b; body_force(2) = 2*rX(2)*invlam1*invlam2*( a*a*invlam1*invlam1 + b*b*invlam2*invlam2 ); return -2*c1*body_force; }
// Test all the functionality inside ContinuumMechanicsProblemDefinition, // which will be common to other problem definition classes void TestContinuumMechanicsProblemDefinition() throw(Exception) { QuadraticMesh<2> mesh(0.5, 1.0, 1.0); ContinuumMechanicsProblemDefinition<2> problem_defn(mesh); TS_ASSERT_THROWS_THIS(problem_defn.Validate(), "No Dirichlet boundary conditions (eg fixed displacement or fixed flow) have been set"); TS_ASSERT_DELTA(problem_defn.GetDensity(), 1.0, 1e-12); TS_ASSERT_EQUALS(problem_defn.GetBodyForceType(), CONSTANT_BODY_FORCE); TS_ASSERT_DELTA(problem_defn.GetConstantBodyForce()(0), 0.0, 1e-12); TS_ASSERT_DELTA(problem_defn.GetConstantBodyForce()(1), 0.0, 1e-12); TS_ASSERT_EQUALS(problem_defn.GetTractionBoundaryConditionType(), NO_TRACTIONS); problem_defn.SetDensity(2.0); TS_ASSERT_DELTA(problem_defn.GetDensity(), 2.0, 1e-12); ////////////////////////////////// // Body force ////////////////////////////////// problem_defn.SetBodyForce(SomeFunction); TS_ASSERT_EQUALS(problem_defn.GetBodyForceType(), FUNCTIONAL_BODY_FORCE); c_vector<double,2> X; X(0) = 10.0; X(1) = 11.0; double t = 0.5; TS_ASSERT_DELTA(problem_defn.EvaluateBodyForceFunction(X,t)(0), 10.5, 1e-12); TS_ASSERT_DELTA(problem_defn.EvaluateBodyForceFunction(X,t)(1), 23.0, 1e-12); TS_ASSERT_DELTA(problem_defn.GetBodyForce(X,t)(0), 10.5, 1e-12); TS_ASSERT_DELTA(problem_defn.GetBodyForce(X,t)(1), 23.0, 1e-12); c_vector<double,2> body_force; body_force(0) = -9.81; body_force(1) = 0.01; problem_defn.SetBodyForce(body_force); TS_ASSERT_EQUALS(problem_defn.GetBodyForceType(), CONSTANT_BODY_FORCE); TS_ASSERT_DELTA(problem_defn.GetConstantBodyForce()(0), -9.81, 1e-12); TS_ASSERT_DELTA(problem_defn.GetConstantBodyForce()(1), 0.01, 1e-12); TS_ASSERT_DELTA(problem_defn.GetBodyForce(X,t)(0), -9.81, 1e-12); TS_ASSERT_DELTA(problem_defn.GetBodyForce(X,t)(1), 0.01, 1e-12); ////////////////////////////////// // Traction ////////////////////////////////// std::vector<BoundaryElement<1,2>*> boundary_elements; std::vector<c_vector<double,2> > tractions; TetrahedralMesh<2,2>::BoundaryElementIterator iter = mesh.GetBoundaryElementIteratorBegin(); c_vector<double,2> vec = zero_vector<double>(2); vec(0)=1.0; boundary_elements.push_back(*iter); tractions.push_back(vec); ++iter; vec(1)=2.0; boundary_elements.push_back(*iter); tractions.push_back(vec); problem_defn.SetTractionBoundaryConditions(boundary_elements, tractions); TS_ASSERT_EQUALS(problem_defn.GetTractionBoundaryConditionType(), ELEMENTWISE_TRACTION); TS_ASSERT_EQUALS(problem_defn.rGetTractionBoundaryElements().size(), 2u); TS_ASSERT_EQUALS(problem_defn.rGetElementwiseTractions().size(), 2u); // comparing addresses TS_ASSERT_EQUALS(problem_defn.rGetTractionBoundaryElements()[0], boundary_elements[0]); TS_ASSERT_EQUALS(problem_defn.rGetTractionBoundaryElements()[1], boundary_elements[1]); TS_ASSERT_DELTA(problem_defn.rGetElementwiseTractions()[0](0), 1.0, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetElementwiseTractions()[0](1), 0.0, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetElementwiseTractions()[1](0), 1.0, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetElementwiseTractions()[1](1), 2.0, 1e-12); ++iter; boundary_elements.push_back(*iter); double pressure = 3423.342; problem_defn.SetApplyNormalPressureOnDeformedSurface(boundary_elements, pressure); TS_ASSERT_EQUALS(problem_defn.rGetTractionBoundaryElements().size(), 3u); // comparing addresses TS_ASSERT_EQUALS(problem_defn.rGetTractionBoundaryElements()[0], boundary_elements[0]); TS_ASSERT_EQUALS(problem_defn.rGetTractionBoundaryElements()[1], boundary_elements[1]); TS_ASSERT_EQUALS(problem_defn.rGetTractionBoundaryElements()[2], boundary_elements[2]); TS_ASSERT_DELTA(problem_defn.GetNormalPressure(), 3423.342, 1e-12); problem_defn.SetPressureScaling(10.0/43.0); TS_ASSERT_DELTA(problem_defn.GetNormalPressure(), 3423.342*10/43, 1e-12); problem_defn.SetPressureScaling(1); ++iter; boundary_elements.push_back(*iter); problem_defn.SetTractionBoundaryConditions(boundary_elements, AnotherFunction); TS_ASSERT_EQUALS(problem_defn.rGetTractionBoundaryElements().size(), 4u); TS_ASSERT_DELTA(problem_defn.EvaluateTractionFunction(X,t)(0), 5.0, 1e-12); TS_ASSERT_DELTA(problem_defn.EvaluateTractionFunction(X,t)(1), 55.0, 1e-12); problem_defn.SetApplyNormalPressureOnDeformedSurface(boundary_elements, PressureFunction); TS_ASSERT_EQUALS(problem_defn.GetTractionBoundaryConditionType(), FUNCTIONAL_PRESSURE_ON_DEFORMED); TS_ASSERT_DELTA(problem_defn.EvaluateNormalPressureFunction(3.64455), 3*3.64455, 1e-12); std::vector<unsigned> fixed_nodes; fixed_nodes.push_back(0); problem_defn.SetZeroDirichletNodes(fixed_nodes); // note, this functionality is all tested properly below // should not throw anything problem_defn.Validate(); }
/* * 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); }
void CompressibleNonlinearElasticitySolver<DIM>::AssembleOnElement( Element<DIM, DIM>& rElement, c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElem, c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElemPrecond, c_vector<double, STENCIL_SIZE>& rBElem, bool assembleResidual, bool assembleJacobian) { static c_matrix<double,DIM,DIM> jacobian; static c_matrix<double,DIM,DIM> inverse_jacobian; double jacobian_determinant; this->mrQuadMesh.GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian); if (assembleJacobian) { rAElem.clear(); rAElemPrecond.clear(); } if (assembleResidual) { rBElem.clear(); } // Get the current displacement at the nodes static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> element_current_displacements; for (unsigned II=0; II<NUM_NODES_PER_ELEMENT; II++) { for (unsigned JJ=0; JJ<DIM; JJ++) { element_current_displacements(JJ,II) = this->mCurrentSolution[DIM*rElement.GetNodeGlobalIndex(II) + JJ]; } } // Allocate memory for the basis functions values and derivative values static c_vector<double, NUM_VERTICES_PER_ELEMENT> linear_phi; static c_vector<double, NUM_NODES_PER_ELEMENT> quad_phi; static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> grad_quad_phi; static c_matrix<double, NUM_NODES_PER_ELEMENT, DIM> trans_grad_quad_phi; // Get the material law AbstractCompressibleMaterialLaw<DIM>* p_material_law = this->mrProblemDefinition.GetCompressibleMaterialLaw(rElement.GetIndex()); static c_matrix<double,DIM,DIM> grad_u; // grad_u = (du_i/dX_M) static c_matrix<double,DIM,DIM> F; // the deformation gradient, F = dx/dX, F_{iM} = dx_i/dX_M static c_matrix<double,DIM,DIM> C; // Green deformation tensor, C = F^T F static c_matrix<double,DIM,DIM> inv_C; // inverse(C) static c_matrix<double,DIM,DIM> inv_F; // inverse(F) static c_matrix<double,DIM,DIM> T; // Second Piola-Kirchoff stress tensor (= dW/dE = 2dW/dC) static c_matrix<double,DIM,DIM> F_T; // F*T static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> F_T_grad_quad_phi; // F*T*grad_quad_phi c_vector<double,DIM> body_force; static FourthOrderTensor<DIM,DIM,DIM,DIM> dTdE; // dTdE(M,N,P,Q) = dT_{MN}/dE_{PQ} static FourthOrderTensor<DIM,DIM,DIM,DIM> dSdF; // dSdF(M,i,N,j) = dS_{Mi}/dF_{jN} static FourthOrderTensor<NUM_NODES_PER_ELEMENT,DIM,DIM,DIM> temp_tensor; static FourthOrderTensor<NUM_NODES_PER_ELEMENT,DIM,NUM_NODES_PER_ELEMENT,DIM> dSdF_quad_quad; static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> temp_matrix; static c_matrix<double,NUM_NODES_PER_ELEMENT,DIM> grad_quad_phi_times_invF; if(this->mSetComputeAverageStressPerElement) { this->mAverageStressesPerElement[rElement.GetIndex()] = zero_vector<double>(DIM*(DIM+1)/2); } // Loop over Gauss points for (unsigned quadrature_index=0; quadrature_index < this->mpQuadratureRule->GetNumQuadPoints(); quadrature_index++) { // This is needed by the cardiac mechanics solver unsigned current_quad_point_global_index = rElement.GetIndex()*this->mpQuadratureRule->GetNumQuadPoints() + quadrature_index; double wJ = jacobian_determinant * this->mpQuadratureRule->GetWeight(quadrature_index); const ChastePoint<DIM>& quadrature_point = this->mpQuadratureRule->rGetQuadPoint(quadrature_index); // Set up basis function information LinearBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, linear_phi); QuadraticBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, quad_phi); QuadraticBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(quadrature_point, inverse_jacobian, grad_quad_phi); trans_grad_quad_phi = trans(grad_quad_phi); // Get the body force, interpolating X if necessary if (assembleResidual) { switch (this->mrProblemDefinition.GetBodyForceType()) { case FUNCTIONAL_BODY_FORCE: { c_vector<double,DIM> X = zero_vector<double>(DIM); // interpolate X (using the vertices and the /linear/ bases, as no curvilinear elements for (unsigned node_index=0; node_index<NUM_VERTICES_PER_ELEMENT; node_index++) { X += linear_phi(node_index)*this->mrQuadMesh.GetNode( rElement.GetNodeGlobalIndex(node_index) )->rGetLocation(); } body_force = this->mrProblemDefinition.EvaluateBodyForceFunction(X, this->mCurrentTime); break; } case CONSTANT_BODY_FORCE: { body_force = this->mrProblemDefinition.GetConstantBodyForce(); break; } default: NEVER_REACHED; } } // Interpolate grad_u grad_u = zero_matrix<double>(DIM,DIM); for (unsigned node_index=0; node_index<NUM_NODES_PER_ELEMENT; node_index++) { for (unsigned i=0; i<DIM; i++) { for (unsigned M=0; M<DIM; M++) { grad_u(i,M) += grad_quad_phi(M,node_index)*element_current_displacements(i,node_index); } } } // Calculate C, inv(C) and T for (unsigned i=0; i<DIM; i++) { for (unsigned M=0; M<DIM; M++) { F(i,M) = (i==M?1:0) + grad_u(i,M); } } C = prod(trans(F),F); inv_C = Inverse(C); inv_F = Inverse(F); // Compute the passive stress, and dTdE corresponding to passive stress this->SetupChangeOfBasisMatrix(rElement.GetIndex(), current_quad_point_global_index); p_material_law->SetChangeOfBasisMatrix(this->mChangeOfBasisMatrix); p_material_law->ComputeStressAndStressDerivative(C, inv_C, 0.0, T, dTdE, assembleJacobian); if(this->mIncludeActiveTension) { // Add any active stresses, if there are any. Requires subclasses to overload this method, // see for example the cardiac mechanics assemblers. this->AddActiveStressAndStressDerivative(C, rElement.GetIndex(), current_quad_point_global_index, T, dTdE, assembleJacobian); } if(this->mSetComputeAverageStressPerElement) { this->AddStressToAverageStressPerElement(T,rElement.GetIndex()); } // Residual vector if (assembleResidual) { F_T = prod(F,T); F_T_grad_quad_phi = prod(F_T, grad_quad_phi); for (unsigned index=0; index<NUM_NODES_PER_ELEMENT*DIM; index++) { unsigned spatial_dim = index%DIM; unsigned node_index = (index-spatial_dim)/DIM; rBElem(index) += - this->mrProblemDefinition.GetDensity() * body_force(spatial_dim) * quad_phi(node_index) * wJ; // The T(M,N)*F(spatial_dim,M)*grad_quad_phi(N,node_index) term rBElem(index) += F_T_grad_quad_phi(spatial_dim,node_index) * wJ; } } // Jacobian matrix if (assembleJacobian) { // Save trans(grad_quad_phi) * invF grad_quad_phi_times_invF = prod(trans_grad_quad_phi, inv_F); ///////////////////////////////////////////////////////////////////////////////////////////// // Set up the tensor dSdF // // dSdF as a function of T and dTdE (which is what the material law returns) is given by: // // dS_{Mi}/dF_{jN} = (dT_{MN}/dC_{PQ}+dT_{MN}/dC_{PQ}) F{iP} F_{jQ} + T_{MN} delta_{ij} // // todo1: this should probably move into the material law (but need to make sure // memory is handled efficiently // todo2: get material law to return this immediately, not dTdE ///////////////////////////////////////////////////////////////////////////////////////////// // Set up the tensor 0.5(dTdE(M,N,P,Q) + dTdE(M,N,Q,P)) for (unsigned M=0; M<DIM; M++) { for (unsigned N=0; N<DIM; N++) { for (unsigned P=0; P<DIM; P++) { for (unsigned Q=0; Q<DIM; Q++) { // this is NOT dSdF, just using this as storage space dSdF(M,N,P,Q) = 0.5*(dTdE(M,N,P,Q) + dTdE(M,N,Q,P)); } } } } // This is NOT dTdE, just reusing memory. A^{MdPQ} = F^d_N * dTdE_sym^{MNPQ} dTdE.template SetAsContractionOnSecondDimension<DIM>(F, dSdF); // dSdF{MdPe} := F^d_N * F^e_Q * dTdE_sym^{MNPQ} dSdF.template SetAsContractionOnFourthDimension<DIM>(F, dTdE); // Now add the T_{MN} delta_{ij} term for (unsigned M=0; M<DIM; M++) { for (unsigned N=0; N<DIM; N++) { for (unsigned i=0; i<DIM; i++) { dSdF(M,i,N,i) += T(M,N); } } } /////////////////////////////////////////////////////// // Set up the tensor // dSdF_quad_quad(node_index1, spatial_dim1, node_index2, spatial_dim2) // = dS_{M,spatial_dim1}/d_F{spatial_dim2,N} // * grad_quad_phi(M,node_index1) // * grad_quad_phi(P,node_index2) // // = dSdF(M,spatial_index1,N,spatial_index2) // * grad_quad_phi(M,node_index1) // * grad_quad_phi(P,node_index2) // /////////////////////////////////////////////////////// temp_tensor.template SetAsContractionOnFirstDimension<DIM>(trans_grad_quad_phi, dSdF); dSdF_quad_quad.template SetAsContractionOnThirdDimension<DIM>(trans_grad_quad_phi, temp_tensor); for (unsigned index1=0; index1<NUM_NODES_PER_ELEMENT*DIM; index1++) { unsigned spatial_dim1 = index1%DIM; unsigned node_index1 = (index1-spatial_dim1)/DIM; for (unsigned index2=0; index2<NUM_NODES_PER_ELEMENT*DIM; index2++) { unsigned spatial_dim2 = index2%DIM; unsigned node_index2 = (index2-spatial_dim2)/DIM; // The dSdF*grad_quad_phi*grad_quad_phi term rAElem(index1,index2) += dSdF_quad_quad(node_index1,spatial_dim1,node_index2,spatial_dim2) * wJ; } } } } rAElemPrecond.clear(); if (assembleJacobian) { rAElemPrecond = rAElem; } if(this->mSetComputeAverageStressPerElement) { for(unsigned i=0; i<DIM*(DIM+1)/2; i++) { this->mAverageStressesPerElement[rElement.GetIndex()](i) /= this->mpQuadratureRule->GetNumQuadPoints(); } } }