void Test2dVariableFibres() throw(Exception) { PlaneStimulusCellFactory<CellLuoRudy1991FromCellML, 2> cell_factory(-1000*1000); TetrahedralMesh<2,2> electrics_mesh; electrics_mesh.ConstructRegularSlabMesh(1.0/96.0/*stepsize*/, 1.0/*length*/, 1.0/*width*/, 1.0/*depth*/); QuadraticMesh<2> mechanics_mesh; mechanics_mesh.ConstructRegularSlabMesh(0.2, 1.0, 1.0, 1.0 /*as above with a different stepsize*/); std::vector<unsigned> fixed_nodes = NonlinearElasticityTools<2>::GetNodesByComponentValue(mechanics_mesh, 0, 0.0); // all the X=0.0 nodes ElectroMechanicsProblemDefinition<2> problem_defn(mechanics_mesh); problem_defn.SetContractionModel(NHS,1.0); problem_defn.SetUseDefaultCardiacMaterialLaw(INCOMPRESSIBLE); problem_defn.SetZeroDisplacementNodes(fixed_nodes); problem_defn.SetMechanicsSolveTimestep(1.0); FileFinder finder("heart/test/data/fibre_tests/5by5mesh_curving_fibres.ortho",RelativeTo::ChasteSourceRoot); problem_defn.SetVariableFibreSheetDirectionsFile(finder, false); HeartConfig::Instance()->SetSimulationDuration(125.0); CardiacElectroMechanicsProblem<2,1> problem(INCOMPRESSIBLE, MONODOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "TestCardiacEmVaryingFibres"); // // fibres going from (1,0) at X=0 to (1,1)-direction at X=1 // // the fibres file was created with the code (inside a class that owns a mesh) // for(unsigned elem_index=0; elem_index<mechanics_mesh.GetNumElements(); elem_index++) // { // double X = mechanics_mesh.GetElement(elem_index)->CalculateCentroid()[0]; // double theta = M_PI*X/4; // std::cout << cos(theta) << " " << sin(theta) << " " << -sin(theta) << " " << cos(theta) << "\n" << std::flush; // } // assert(0); // problem.SetNoElectricsOutput(); problem.Solve(); // test by checking the length of the tissue against hardcoded value std::vector<c_vector<double,2> >& r_deformed_position = problem.rGetDeformedPosition(); // visualised, looks good - contracts in X-direction near the fixed surface, // but on the other side the fibres are in the (1,1) direction, so contraction // pulls the tissue downward a bit TS_ASSERT_DELTA(r_deformed_position[5](0), 0.9055, 2e-3); //IntelProduction differs by about 1.6e-3... MechanicsEventHandler::Headings(); MechanicsEventHandler::Report(); }
// with stretch (and stretch-rate) independent contraction models the implicit and explicit schemes // are identical void TestCompareImplicitAndExplicitWithStretchIndependentContractionModel() throw(Exception) { QuadraticMesh<2> mesh(0.25, 1.0, 1.0); MooneyRivlinMaterialLaw<2> law(1); std::vector<unsigned> fixed_nodes = NonlinearElasticityTools<2>::GetNodesByComponentValue(mesh,0,0.0); ElectroMechanicsProblemDefinition<2> problem_defn(mesh); problem_defn.SetMaterialLaw(INCOMPRESSIBLE,&law); problem_defn.SetZeroDisplacementNodes(fixed_nodes); problem_defn.SetContractionModel(NONPHYSIOL1,0.01); problem_defn.SetMechanicsSolveTimestep(0.01); //This is only set to make ElectroMechanicsProblemDefinition::Validate pass //The following lines are not relevant to this test but need to be there TetrahedralMesh<2,2>* p_fine_mesh = new TetrahedralMesh<2,2>();//unused in this test p_fine_mesh->ConstructRegularSlabMesh(0.25, 1.0, 1.0); TetrahedralMesh<2,2>* p_coarse_mesh = new TetrahedralMesh<2,2>();//unused in this test p_coarse_mesh->ConstructRegularSlabMesh(0.25, 1.0, 1.0); FineCoarseMeshPair<2>* p_pair = new FineCoarseMeshPair<2>(*p_fine_mesh, *p_coarse_mesh);//also unused in this test p_pair->SetUpBoxesOnFineMesh(); ///////////////////////////////////////////////////////////////////// // NONPHYSIOL 1 - contraction model is of the form sin(t) IncompressibleExplicitSolver2d expl_solver(mesh,problem_defn,""/*"TestCompareExplAndImplCardiacSolvers_Exp"*/); p_pair->ComputeFineElementsAndWeightsForCoarseQuadPoints(*(expl_solver.GetQuadratureRule()), false); p_pair->DeleteFineBoxCollection(); expl_solver.SetFineCoarseMeshPair(p_pair); expl_solver.Initialise(); IncompressibleImplicitSolver2d impl_solver(mesh,problem_defn,""/*"TestCompareExplAndImplCardiacSolvers_Imp"*/); impl_solver.SetFineCoarseMeshPair(p_pair); impl_solver.Initialise(); double dt = 0.25; for(double t=0; t<3; t+=dt) { expl_solver.Solve(t,t+dt,dt); impl_solver.Solve(t,t+dt,dt); // computations should be identical TS_ASSERT_EQUALS(expl_solver.GetNumNewtonIterations(), impl_solver.GetNumNewtonIterations()); for(unsigned i=0; i<mesh.GetNumNodes(); i++) { TS_ASSERT_DELTA(expl_solver.rGetDeformedPosition()[i](0), impl_solver.rGetDeformedPosition()[i](0), 1e-9); TS_ASSERT_DELTA(expl_solver.rGetDeformedPosition()[i](1), impl_solver.rGetDeformedPosition()[i](1), 1e-9); } } //in need of deletion even if all these 3 have no influence at all on this test delete p_fine_mesh; delete p_coarse_mesh; delete p_pair; }
void TestWithSimpleContractionModel() throw(Exception) { QuadraticMesh<2> mesh(0.25, 1.0, 1.0); MooneyRivlinMaterialLaw<2> law(1); std::vector<unsigned> fixed_nodes = NonlinearElasticityTools<2>::GetNodesByComponentValue(mesh,0,0.0); ElectroMechanicsProblemDefinition<2> problem_defn(mesh); problem_defn.SetMaterialLaw(INCOMPRESSIBLE,&law); problem_defn.SetZeroDisplacementNodes(fixed_nodes); problem_defn.SetContractionModel(NONPHYSIOL1,0.01); problem_defn.SetMechanicsSolveTimestep(0.01); //This is only set to make ElectroMechanicsProblemDefinition::Validate pass // NONPHYSIOL1 => NonphysiologicalContractionModel 1 IncompressibleExplicitSolver2d solver(mesh,problem_defn,"TestExplicitCardiacMech"); //The following lines are not relevant to this test but need to be there TetrahedralMesh<2,2>* p_fine_mesh = new TetrahedralMesh<2,2>();//unused in this test p_fine_mesh->ConstructRegularSlabMesh(0.25, 1.0, 1.0); TetrahedralMesh<2,2>* p_coarse_mesh = new TetrahedralMesh<2,2>();//unused in this test p_coarse_mesh->ConstructRegularSlabMesh(0.25, 1.0, 1.0); FineCoarseMeshPair<2>* p_pair = new FineCoarseMeshPair<2>(*p_fine_mesh, *p_coarse_mesh);//also unused in this test p_pair->SetUpBoxesOnFineMesh(); p_pair->ComputeFineElementsAndWeightsForCoarseQuadPoints(*(solver.GetQuadratureRule()), false); p_pair->DeleteFineBoxCollection(); solver.SetFineCoarseMeshPair(p_pair); /////////////////////////////////////////////////////////////////////////// solver.Initialise(); // coverage QuadraturePointsGroup<2> quad_points(mesh, *(solver.GetQuadratureRule())); std::vector<double> calcium_conc(solver.GetTotalNumQuadPoints(), 0.0); std::vector<double> voltages(solver.GetTotalNumQuadPoints(), 0.0); solver.SetCalciumAndVoltage(calcium_conc, voltages); // solve UP TO t=0. So Ta(lam_n,t_{n+1})=5*sin(0)=0, ie no deformation solver.Solve(-0.01,0.0,0.01); TS_ASSERT_EQUALS(solver.GetNumNewtonIterations(),0u); solver.Solve(0.24,0.25,0.01); TS_ASSERT_DELTA(solver.rGetDeformedPosition()[4](0), 0.9732, 1e-2); TS_ASSERT_DELTA(solver.rGetDeformedPosition()[4](1), -0.0156, 1e-2); //in need of deletion even if all these 3 have no influence at all on this test delete p_fine_mesh; delete p_coarse_mesh; delete p_pair; }
void TestAssemblerMeshType() { TetrahedralMesh<2,2> mesh; ContinuumMechanicsProblemDefinition<2> problem_defn(mesh); TS_ASSERT_THROWS_CONTAINS(ContinuumMechanicsNeumannBcsAssembler<2>(&mesh, &problem_defn), "Continuum mechanics solvers require a quadratic mesh"); TetrahedralMesh<3,3> mesh3d; ContinuumMechanicsProblemDefinition<3> problem_defn3d(mesh3d); TS_ASSERT_THROWS_CONTAINS(ContinuumMechanicsNeumannBcsAssembler<3>(&mesh3d, &problem_defn3d), "Continuum mechanics solvers require a quadratic mesh"); }
void TestExceptions() throw(Exception) { EntirelyStimulatedTissueCellFactory cell_factory; TetrahedralMesh<2,2> electrics_mesh; electrics_mesh.ConstructRegularSlabMesh(0.01, 0.05, 0.05); QuadraticMesh<2> mechanics_mesh; mechanics_mesh.ConstructRegularSlabMesh(0.025, 0.05, 0.05); HeartConfig::Instance()->SetSimulationDuration(10.0); ElectroMechanicsProblemDefinition<2> problem_defn(mechanics_mesh); problem_defn.SetContractionModel(KERCHOFFS2003,1.0); problem_defn.SetUseDefaultCardiacMaterialLaw(COMPRESSIBLE); problem_defn.SetMechanicsSolveTimestep(1.0); // Note that TS_ASSERT_THROWS_THIS won't compile if you just construct // an object with two templates. This is because the macro thinks // that the comma between the template parameters separates arguments, so it thinks // there are 3 arguments and it is expecting 2 try { CardiacElectroMechanicsProblem<2,2> problem (COMPRESSIBLE, MONODOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "blahblah"); } catch (Exception& e) { TS_ASSERT_EQUALS(e.GetShortMessage(), "The second template parameter should be 1 when a monodomain problem is chosen"); } try { CardiacElectroMechanicsProblem<2,1> problem (COMPRESSIBLE, BIDOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "blahblah"); } catch (Exception& e) { TS_ASSERT_EQUALS(e.GetShortMessage(), "The second template parameter should be 2 when a bidomain problem is chosen"); } }
/* == 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); }
void Test3d() throw(Exception) { PlaneStimulusCellFactory<CellLuoRudy1991FromCellML, 3> cell_factory(-1000*1000); // set up two meshes of 1mm by 1mm by 1mm TetrahedralMesh<3,3> electrics_mesh; electrics_mesh.ConstructRegularSlabMesh(0.01, 0.1, 0.1, 0.1); QuadraticMesh<3> mechanics_mesh(0.1, 0.1, 0.1, 0.1); // fix the nodes on x=0 std::vector<unsigned> fixed_nodes = NonlinearElasticityTools<3>::GetNodesByComponentValue(mechanics_mesh,0,0); HeartConfig::Instance()->SetSimulationDuration(50.0); ElectroMechanicsProblemDefinition<3> problem_defn(mechanics_mesh); problem_defn.SetContractionModel(NHS,1.0); problem_defn.SetUseDefaultCardiacMaterialLaw(INCOMPRESSIBLE); problem_defn.SetZeroDisplacementNodes(fixed_nodes); problem_defn.SetMechanicsSolveTimestep(1.0); CardiacElectroMechanicsProblem<3,1> problem(INCOMPRESSIBLE, MONODOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "TestCardiacElectroMech3d"); problem.SetNoElectricsOutput(); problem.Solve(); // test by checking the length of the tissue against hardcoded value c_vector<double,3> undeformed_node_1 = mechanics_mesh.GetNode(1)->rGetLocation(); TS_ASSERT_DELTA(undeformed_node_1(0), 0.1, 1e-6); TS_ASSERT_DELTA(undeformed_node_1(1), 0.0, 1e-6); TS_ASSERT_DELTA(undeformed_node_1(2), 0.0, 1e-6); std::vector<c_vector<double,3> >& r_deformed_position = problem.rGetDeformedPosition(); TS_ASSERT_DELTA(r_deformed_position[1](0), 0.0879, 1e-3); MechanicsEventHandler::Headings(); MechanicsEventHandler::Report(); }
/* HOW_TO_TAG Cardiac/Electro-mechanics * Run electro-mechanical simulations using bidomain instead of monodomain * * This test is the same as above but with bidomain instead of monodomain. * Extracellular conductivities are set very high so the results should be the same. */ void TestWithHomogeneousEverythingCompressibleBidomain() throw(Exception) { EntirelyStimulatedTissueCellFactory cell_factory; TetrahedralMesh<2,2> electrics_mesh; electrics_mesh.ConstructRegularSlabMesh(0.01, 0.05, 0.05); QuadraticMesh<2> mechanics_mesh; mechanics_mesh.ConstructRegularSlabMesh(0.025, 0.05, 0.05); std::vector<unsigned> fixed_nodes; std::vector<c_vector<double,2> > fixed_node_locations; // fix the node at the origin so that the solution is well-defined (ie unique) fixed_nodes.push_back(0); fixed_node_locations.push_back(zero_vector<double>(2)); // for the rest of the nodes, if they lie on X=0, fix x=0 but leave y free. for(unsigned i=1 /*not 0*/; i<mechanics_mesh.GetNumNodes(); i++) { if(fabs(mechanics_mesh.GetNode(i)->rGetLocation()[0])<1e-6) { c_vector<double,2> new_position; new_position(0) = 0.0; new_position(1) = SolidMechanicsProblemDefinition<2>::FREE; fixed_nodes.push_back(i); fixed_node_locations.push_back(new_position); } } ElectroMechanicsProblemDefinition<2> problem_defn(mechanics_mesh); problem_defn.SetContractionModel(KERCHOFFS2003,1.0); problem_defn.SetUseDefaultCardiacMaterialLaw(COMPRESSIBLE); problem_defn.SetFixedNodes(fixed_nodes, fixed_node_locations); problem_defn.SetMechanicsSolveTimestep(1.0); // the following is just for coverage - applying a zero pressure so has no effect on deformation std::vector<BoundaryElement<1,2>*> boundary_elems; boundary_elems.push_back(* (mechanics_mesh.GetBoundaryElementIteratorBegin())); problem_defn.SetApplyNormalPressureOnDeformedSurface(boundary_elems, 0.0); HeartConfig::Instance()->SetSimulationDuration(10.0); HeartConfig::Instance()->SetExtracellularConductivities(Create_c_vector(1500,1500,1500)); //creates the EM problem with ELEC_PROB_DIM=2 CardiacElectroMechanicsProblem<2,2> problem(COMPRESSIBLE, BIDOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "TestCardiacEmHomogeneousEverythingCompressibleBidomain"); problem.Solve(); std::vector<c_vector<double,2> >& r_deformed_position = problem.rGetDeformedPosition(); // not sure how easy is would be determine what the deformation should be // exactly, but it certainly should be constant squash in X direction, constant // stretch in Y. // first, check node 8 starts is the far corner assert(fabs(mechanics_mesh.GetNode(8)->rGetLocation()[0] - 0.05)<1e-8); assert(fabs(mechanics_mesh.GetNode(8)->rGetLocation()[1] - 0.05)<1e-8); double X_scale_factor = r_deformed_position[8](0)/0.05; double Y_scale_factor = r_deformed_position[8](1)/0.05; std::cout << "Scale_factors = " << X_scale_factor << " " << Y_scale_factor << ", product = " << X_scale_factor*Y_scale_factor<<"\n"; for(unsigned i=0; i<mechanics_mesh.GetNumNodes(); i++) { double X = mechanics_mesh.GetNode(i)->rGetLocation()[0]; double Y = mechanics_mesh.GetNode(i)->rGetLocation()[1]; TS_ASSERT_DELTA( r_deformed_position[i](0), X * X_scale_factor, 1e-6); TS_ASSERT_DELTA( r_deformed_position[i](1), Y * Y_scale_factor, 1e-6); } //check interpolated voltages and calcium unsigned quad_points = problem.mpCardiacMechSolver->GetTotalNumQuadPoints(); TS_ASSERT_EQUALS(problem.mInterpolatedVoltages.size(), quad_points); TS_ASSERT_EQUALS(problem.mInterpolatedCalciumConcs.size(), quad_points); //two hardcoded values TS_ASSERT_DELTA(problem.mInterpolatedVoltages[0],9.267,1e-3); TS_ASSERT_DELTA(problem.mInterpolatedCalciumConcs[0],0.001464,1e-6); //for the rest, we check that, at the end of this simulation, all quad nodes have V and Ca above a certain threshold for(unsigned i = 0; i < quad_points; i++) { TS_ASSERT_LESS_THAN(9.2,problem.mInterpolatedVoltages[i]); TS_ASSERT_LESS_THAN(0.0014,problem.mInterpolatedCalciumConcs[i]); } //check default value of whether there is a bath or not TS_ASSERT_EQUALS(problem.mpElectricsProblem->GetHasBath(), false); //test the functionality of having phi_e on the mechanics mesh (values are tested somewhere else) Hdf5DataReader data_reader("TestCardiacEmHomogeneousEverythingCompressibleBidomain/electrics","voltage_mechanics_mesh"); TS_ASSERT_THROWS_NOTHING(data_reader.GetVariableOverTime("Phi_e",0u)); }
// In this test all nodes are stimulated at the same time, hence // exactly the same active tension is generated at all nodes, // so the internal force in entirely homogeneous. // We fix the x-coordinate of the X=0, but leave the y-coordinate // free - ie sliding boundary conditions. With a homogeneous // force this means the solution should be a perfect rectangle, // ie x=alpha*X, y=beta*Y for some alpha, beta. void TestWithHomogeneousEverythingCompressible() throw(Exception) { EntirelyStimulatedTissueCellFactory cell_factory; TetrahedralMesh<2,2> electrics_mesh; electrics_mesh.ConstructRegularSlabMesh(0.01, 0.05, 0.05); QuadraticMesh<2> mechanics_mesh; mechanics_mesh.ConstructRegularSlabMesh(0.025, 0.05, 0.05); std::vector<unsigned> fixed_nodes; std::vector<c_vector<double,2> > fixed_node_locations; // fix the node at the origin so that the solution is well-defined (ie unique) fixed_nodes.push_back(0); fixed_node_locations.push_back(zero_vector<double>(2)); // for the rest of the nodes, if they lie on X=0, fix x=0 but leave y free. for(unsigned i=1 /*not 0*/; i<mechanics_mesh.GetNumNodes(); i++) { if(fabs(mechanics_mesh.GetNode(i)->rGetLocation()[0])<1e-6) { c_vector<double,2> new_position; new_position(0) = 0.0; new_position(1) = SolidMechanicsProblemDefinition<2>::FREE; fixed_nodes.push_back(i); fixed_node_locations.push_back(new_position); } } ElectroMechanicsProblemDefinition<2> problem_defn(mechanics_mesh); problem_defn.SetContractionModel(KERCHOFFS2003,1.0); problem_defn.SetUseDefaultCardiacMaterialLaw(COMPRESSIBLE); problem_defn.SetFixedNodes(fixed_nodes, fixed_node_locations); problem_defn.SetMechanicsSolveTimestep(1.0); // the following is just for coverage - applying a zero pressure so has no effect on deformation std::vector<BoundaryElement<1,2>*> boundary_elems; boundary_elems.push_back(* (mechanics_mesh.GetBoundaryElementIteratorBegin())); problem_defn.SetApplyNormalPressureOnDeformedSurface(boundary_elems, 0.0); HeartConfig::Instance()->SetSimulationDuration(10.0); CardiacElectroMechanicsProblem<2,1> problem(COMPRESSIBLE, MONODOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "TestCardiacEmHomogeneousEverythingCompressible"); problem.Solve(); std::vector<c_vector<double,2> >& r_deformed_position = problem.rGetDeformedPosition(); // not sure how easy is would be determine what the deformation should be // exactly, but it certainly should be constant squash in X direction, constant // stretch in Y. // first, check node 8 starts is the far corner assert(fabs(mechanics_mesh.GetNode(8)->rGetLocation()[0] - 0.05)<1e-8); assert(fabs(mechanics_mesh.GetNode(8)->rGetLocation()[1] - 0.05)<1e-8); double X_scale_factor = r_deformed_position[8](0)/0.05; double Y_scale_factor = r_deformed_position[8](1)/0.05; std::cout << "Scale_factors = " << X_scale_factor << " " << Y_scale_factor << ", product = " << X_scale_factor*Y_scale_factor<<"\n"; for(unsigned i=0; i<mechanics_mesh.GetNumNodes(); i++) { double X = mechanics_mesh.GetNode(i)->rGetLocation()[0]; double Y = mechanics_mesh.GetNode(i)->rGetLocation()[1]; TS_ASSERT_DELTA( r_deformed_position[i](0), X * X_scale_factor, 1e-6); TS_ASSERT_DELTA( r_deformed_position[i](1), Y * Y_scale_factor, 1e-6); } // coverage TS_ASSERT(problem.GetMechanicsSolver()!=NULL); }
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); }
/** * Solve a problem in which the traction boundary condition is a normal pressure * applied to a surface of the deformed body. * * The initial body is the unit cube. The deformation is given by x=Rx', where * R is a rotation matrix and x'=[X/(lambda^2) lambda*Y lambda*Z], ie simple stretching * followed by a rotation. * * Ignoring R for the time being, the deformation gradient would be * F=diag(1/lambda^2, lambda, lambda). * Assuming a Neo-Hookean material law, the first Piola-Kirchoff tensor is * * S = 2c F^T - p F^{-1} * = diag ( 2c lam^{-2} - p lam^2, 2c lam - p lam^{-1}, 2c lam - p lam^{-1} * * We choose the internal pressure (p) to be 2clam^2, so that * * S = diag( 2c (lam^{-2} - lam^4), 0, 0) * * To obtain this deformation as the solution, we can specify fixed location boundary conditions * on X=0 and specify a traction of t=[2c(lam^{-2} - lam^4) 0 0] on the X=1 surface. * * Now, including the rotation matrix, we can specify appropriate fixed location boundary conditions * on the X=0 surface, and specify that a normal pressure should act on the *deformed* X=1 surface, * with value P = 2c(lam^{-2} - lam^4)/(lambda^2) [divided by lambda^2 as the deformed surface * has a smaller area than the undeformed surface. * */ void TestWithPressureOnDeformedSurfaceBoundaryCondition3d() throw (Exception) { double lambda = 0.85; unsigned num_elem_each_dir = 5; QuadraticMesh<3> mesh(1.0/num_elem_each_dir, 1.0, 1.0, 1.0); // Neo-Hookean material law double c1 = 1.0; MooneyRivlinMaterialLaw<3> law(c1, 0.0); // anything will do here c_matrix<double,3,3> rotation_matrix = identity_matrix<double>(3); rotation_matrix(0,0)=1.0/sqrt(2.0); rotation_matrix(0,1)=-1.0/sqrt(2.0); rotation_matrix(1,0)=1.0/sqrt(2.0); rotation_matrix(1,1)=1.0/sqrt(2.0); // Define displacement boundary conditions std::vector<unsigned> fixed_nodes; std::vector<c_vector<double,3> > locations; for (unsigned i=0; i<mesh.GetNumNodes(); i++) { double X = mesh.GetNode(i)->rGetLocation()[0]; double Y = mesh.GetNode(i)->rGetLocation()[1]; double Z = mesh.GetNode(i)->rGetLocation()[2]; // if X=0 if ( fabs(X)<1e-6) { fixed_nodes.push_back(i); c_vector<double,3> new_position_before_rotation; new_position_before_rotation(0) = 0.0; new_position_before_rotation(1) = lambda*Y; new_position_before_rotation(2) = lambda*Z; locations.push_back(prod(rotation_matrix, new_position_before_rotation)); } } assert(fixed_nodes.size()==(2*num_elem_each_dir+1)*(2*num_elem_each_dir+1)); // Define pressure boundary conditions X=1 surface std::vector<BoundaryElement<2,3>*> boundary_elems; double pressure_bc = 2*c1*((pow(lambda,-2.0)-pow(lambda,4.0)))/(lambda*lambda); for (TetrahedralMesh<3,3>::BoundaryElementIterator iter = mesh.GetBoundaryElementIteratorBegin(); iter != mesh.GetBoundaryElementIteratorEnd(); ++iter) { if (fabs((*iter)->CalculateCentroid()[0]-1)<1e-6) { BoundaryElement<2,3>* p_element = *iter; boundary_elems.push_back(p_element); } } assert(boundary_elems.size()==2*num_elem_each_dir*num_elem_each_dir); SolidMechanicsProblemDefinition<3> problem_defn(mesh); problem_defn.SetMaterialLaw(INCOMPRESSIBLE,&law); problem_defn.SetFixedNodes(fixed_nodes, locations); problem_defn.SetApplyNormalPressureOnDeformedSurface(boundary_elems, pressure_bc); IncompressibleNonlinearElasticitySolver<3> solver(mesh, problem_defn, "nonlin_elas_3d_pressure_on_deformed"); solver.Solve(); // Compare std::vector<c_vector<double,3> >& r_solution = solver.rGetDeformedPosition(); for (unsigned i=0; i<mesh.GetNumNodes(); i++) { double X = mesh.GetNode(i)->rGetLocation()[0]; double Y = mesh.GetNode(i)->rGetLocation()[1]; double Z = mesh.GetNode(i)->rGetLocation()[2]; c_vector<double,3> new_position_before_rotation; new_position_before_rotation(0) = X/(lambda*lambda); new_position_before_rotation(1) = lambda*Y; new_position_before_rotation(2) = lambda*Z; c_vector<double,3> new_position = prod(rotation_matrix, new_position_before_rotation); TS_ASSERT_DELTA(r_solution[i](0), new_position(0), 1e-2); TS_ASSERT_DELTA(r_solution[i](1), new_position(1), 1e-2); TS_ASSERT_DELTA(r_solution[i](2), new_position(2), 1e-2); } // test the pressures std::vector<double>& r_pressures = solver.rGetPressures(); for (unsigned i=0; i<r_pressures.size(); i++) { TS_ASSERT_DELTA(r_pressures[i], 2*c1*lambda*lambda, 0.05); } }
void TestAssembler3d() throw (Exception) { QuadraticMesh<3> mesh(1.0, 1.0, 1.0, 1.0); ContinuumMechanicsProblemDefinition<3> problem_defn(mesh); double t1 = 2.6854233; double t2 = 3.2574578; double t3 = 4.5342308; // for the boundary element on the z=0 surface, create a traction std::vector<BoundaryElement<2,3>*> boundary_elems; std::vector<c_vector<double,3> > tractions; c_vector<double,3> traction; traction(0) = t1; traction(1) = t2; traction(2) = t3; for (TetrahedralMesh<3,3>::BoundaryElementIterator iter = mesh.GetBoundaryElementIteratorBegin(); iter != mesh.GetBoundaryElementIteratorEnd(); ++iter) { if (fabs((*iter)->CalculateCentroid()[2])<1e-4) { BoundaryElement<2,3>* p_element = *iter; boundary_elems.push_back(p_element); tractions.push_back(traction); } } assert(boundary_elems.size()==2); problem_defn.SetTractionBoundaryConditions(boundary_elems, tractions); Vec vec = PetscTools::CreateVec(4*mesh.GetNumNodes()); ContinuumMechanicsNeumannBcsAssembler<3> assembler(&mesh, &problem_defn); assembler.SetVectorToAssemble(vec, true); assembler.AssembleVector(); ReplicatableVector vec_repl(vec); // For boundary elements in 3d - ie for 2d elements - and with QUADRATIC basis functions, we have // \intgl_{canonical element} \phi_i dV = 0.0 for i=0,1,2 (ie bases on vertices) // and // \intgl_{canonical element} \phi_i dV = 1/6 for i=3,4,5 (ie bases on mid-nodes) for(unsigned i=0; i<mesh.GetNumNodes(); i++) { double x = mesh.GetNode(i)->rGetLocation()[0]; double y = mesh.GetNode(i)->rGetLocation()[1]; double z = mesh.GetNode(i)->rGetLocation()[2]; if(fabs(z)<1e-8) // ie if z=0 { if( fabs(x-0.5)<1e-8 || fabs(y-0.5)<1e-8 ) // if x=0.5 or y=0.5 { unsigned num_surf_elems_contained_in = 1; if( fabs(x+y-1.0)<1e-8 ) // if x=0.5 AND y=0.5 { num_surf_elems_contained_in = 2; } // interior node on traction surface TS_ASSERT_DELTA(vec_repl[4*i], num_surf_elems_contained_in * t1/6.0, 1e-8); TS_ASSERT_DELTA(vec_repl[4*i+1], num_surf_elems_contained_in * t2/6.0, 1e-8); TS_ASSERT_DELTA(vec_repl[4*i+2], num_surf_elems_contained_in * t3/6.0, 1e-8); } else { TS_ASSERT_DELTA(vec_repl[4*i], 0.0, 1e-8); TS_ASSERT_DELTA(vec_repl[4*i+1], 0.0, 1e-8); TS_ASSERT_DELTA(vec_repl[4*i+2], 0.0, 1e-8); } } else { TS_ASSERT_DELTA(vec_repl[4*i], 0.0, 1e-8); TS_ASSERT_DELTA(vec_repl[4*i+1], 0.0, 1e-8); TS_ASSERT_DELTA(vec_repl[4*i+2], 0.0, 1e-8); } } for(unsigned i=0; i<mesh.GetNumNodes(); i++) { TS_ASSERT_DELTA(vec_repl[4*i+3], 0.0, 1e-8); } PetscTools::Destroy(vec); }
/* NOTE: This test has a twin in heart/test/tutorials/TestCardiacElectroMechanicsTutorial.hpp * TestCardiacElectroMechanicsTutorial::dontTestTwistingCube() * * If you need to re-generate the fibres for this test * * Remove "dont" from the tutorial * * Rerun it * * Copy output cp /tmp/$USER/testoutput/TutorialFibreFiles/5by5by5_fibres.orthoquad heart/test/data/fibre_tests/5by5by5_fibres_by_quadpt.orthoquad */ void TestTwistingCube() throw(Exception) { PlaneStimulusCellFactory<CellLuoRudy1991FromCellML, 3> cell_factory(-1000*1000); // set up two meshes of 1mm by 1mm by 1mm TetrahedralMesh<3,3> electrics_mesh; electrics_mesh.ConstructRegularSlabMesh(0.01, 0.1, 0.1, 0.1); QuadraticMesh<3> mechanics_mesh(0.02, 0.1, 0.1, 0.1); // fix the nodes on Z=0 std::vector<unsigned> fixed_nodes = NonlinearElasticityTools<3>::GetNodesByComponentValue(mechanics_mesh,2,0.0); HeartConfig::Instance()->SetSimulationDuration(50.0); ElectroMechanicsProblemDefinition<3> problem_defn(mechanics_mesh); problem_defn.SetContractionModel(KERCHOFFS2003,1.0); problem_defn.SetUseDefaultCardiacMaterialLaw(INCOMPRESSIBLE); problem_defn.SetZeroDisplacementNodes(fixed_nodes); problem_defn.SetMechanicsSolveTimestep(1.0); FileFinder finder("heart/test/data/fibre_tests/5by5by5_fibres_by_quadpt.orthoquad",RelativeTo::ChasteSourceRoot); problem_defn.SetVariableFibreSheetDirectionsFile(finder, true); CardiacElectroMechanicsProblem<3,1> problem(INCOMPRESSIBLE, MONODOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "TestCardiacElectroMech3dTwistingCube"); problem.Solve(); // verified that it twists by visualising, some hardcoded values here.. { //Check that we are indexing the relevant corners of the cube c_vector<double, 3> undeformed_node1 = mechanics_mesh.GetNode(6*6*5)->rGetLocation(); TS_ASSERT_DELTA(undeformed_node1(0), 0.0, 1e-6); TS_ASSERT_DELTA(undeformed_node1(1), 0.0, 1e-6); TS_ASSERT_DELTA(undeformed_node1(2), 0.1, 1e-6); c_vector<double, 3> undeformed_node2 = mechanics_mesh.GetNode(6*6*6-1)->rGetLocation(); TS_ASSERT_DELTA(undeformed_node2(0), 0.1, 1e-6); TS_ASSERT_DELTA(undeformed_node2(1), 0.1, 1e-6); TS_ASSERT_DELTA(undeformed_node2(2), 0.1, 1e-6); } std::vector<c_vector<double,3> >& r_deformed_position = problem.rGetDeformedPosition(); TS_ASSERT_DELTA(r_deformed_position[6*6*5](0), 0.0116, 1e-3); TS_ASSERT_DELTA(r_deformed_position[6*6*5](1), -0.0141, 1e-3); TS_ASSERT_DELTA(r_deformed_position[6*6*5](2), 0.1007, 1e-3); TS_ASSERT_DELTA(r_deformed_position[6*6*6-1](0), 0.0872, 1e-3); TS_ASSERT_DELTA(r_deformed_position[6*6*6-1](1), 0.1138, 1e-3); TS_ASSERT_DELTA(r_deformed_position[6*6*6-1](2), 0.1015, 1e-3); MechanicsEventHandler::Headings(); MechanicsEventHandler::Report(); }
/** * Solve 3D nonlinear elasticity problem with an exact solution. * * For full details see Pathmanathan, Gavaghan, Whiteley "A comparison of numerical * methods used in finite element modelling of soft tissue deformation", J. Strain * Analysis, to appear. * * We solve a 3d problem on a cube with a Neo-Hookean material, assuming the solution * will be * x = X+aX^2/2 * y = Y+bY^2/2 * z = Z/(1+aX)(1+bY) * with p=2c (c the material parameter), * which, note, has been constructed to be an incompressible. We assume displacement * boundary conditions on X=0 and traction boundary conditions on the remaining 5 surfaces. * It is then possible to determine the body force and surface tractions required for * this deformation, and they are defined in the above class. */ void TestSolve3d() throw(Exception) { unsigned num_runs=4; double l2_errors[4]; unsigned num_elem_each_dir[4] = {1,2,5,10}; for(unsigned run=0; run<num_runs; run++) { QuadraticMesh<3> mesh(1.0/num_elem_each_dir[run], 1.0, 1.0, 1.0); // Neo-Hookean material law MooneyRivlinMaterialLaw<3> law(ThreeDimensionalModelProblem::c1, 0.0); // Define displacement boundary conditions std::vector<unsigned> fixed_nodes; std::vector<c_vector<double,3> > locations; for (unsigned i=0; i<mesh.GetNumNodes(); i++) { double X = mesh.GetNode(i)->rGetLocation()[0]; double Y = mesh.GetNode(i)->rGetLocation()[1]; double Z = mesh.GetNode(i)->rGetLocation()[2]; // if X=0 if ( fabs(X)<1e-6) { fixed_nodes.push_back(i); c_vector<double,3> new_position; new_position(0) = 0.0; new_position(1) = Y + Y*Y*ThreeDimensionalModelProblem::b/2.0; new_position(2) = Z/((1+X*ThreeDimensionalModelProblem::a)*(1+Y*ThreeDimensionalModelProblem::b)); locations.push_back(new_position); } } assert(fixed_nodes.size()==(2*num_elem_each_dir[run]+1)*(2*num_elem_each_dir[run]+1)); // Define traction boundary conditions on all boundary elems that are not on X=0 std::vector<BoundaryElement<2,3>*> boundary_elems; for (TetrahedralMesh<3,3>::BoundaryElementIterator iter = mesh.GetBoundaryElementIteratorBegin(); iter != mesh.GetBoundaryElementIteratorEnd(); ++iter) { if (fabs((*iter)->CalculateCentroid()[0])>1e-6) { BoundaryElement<2,3>* p_element = *iter; boundary_elems.push_back(p_element); } } assert(boundary_elems.size()==10*num_elem_each_dir[run]*num_elem_each_dir[run]); SolidMechanicsProblemDefinition<3> problem_defn(mesh); problem_defn.SetMaterialLaw(INCOMPRESSIBLE,&law); problem_defn.SetFixedNodes(fixed_nodes, locations); problem_defn.SetBodyForce(ThreeDimensionalModelProblem::GetBodyForce); problem_defn.SetTractionBoundaryConditions(boundary_elems, ThreeDimensionalModelProblem::GetTraction); IncompressibleNonlinearElasticitySolver<3> solver(mesh, problem_defn, "nonlin_elas_3d_10"); solver.Solve(1e-12); // note newton tolerance of 1e-12 needed for convergence to occur solver.CreateCmguiOutput(); // Compare l2_errors[run] = 0; std::vector<c_vector<double,3> >& r_solution = solver.rGetDeformedPosition(); for (unsigned i=0; i<mesh.GetNumNodes(); i++) { double X = mesh.GetNode(i)->rGetLocation()[0]; double Y = mesh.GetNode(i)->rGetLocation()[1]; double Z = mesh.GetNode(i)->rGetLocation()[2]; double exact_x = X + X*X*ThreeDimensionalModelProblem::a/2.0; double exact_y = Y + Y*Y*ThreeDimensionalModelProblem::b/2.0; double exact_z = Z/((1+X*ThreeDimensionalModelProblem::a)*(1+Y*ThreeDimensionalModelProblem::b)); c_vector<double,3> error; error(0) = r_solution[i](0) - exact_x; error(1) = r_solution[i](1) - exact_y; error(2) = r_solution[i](2) - exact_z; l2_errors[run] += norm_2(error); if(num_elem_each_dir[run]==5u) { TS_ASSERT_DELTA(r_solution[i](0), exact_x, 1e-2); TS_ASSERT_DELTA(r_solution[i](1), exact_y, 1e-2); TS_ASSERT_DELTA(r_solution[i](2), exact_z, 1e-2); } } l2_errors[run] /= mesh.GetNumNodes(); std::vector<double>& r_pressures = solver.rGetPressures(); for (unsigned i=0; i<r_pressures.size(); i++) { TS_ASSERT_DELTA(r_pressures[i]/(2*ThreeDimensionalModelProblem::c1), 1.0, 2e-1); } } //for(unsigned run=0; run<num_runs; run++) //{ // std::cout << 1.0/num_elem_each_dir[run] << " " << l2_errors[run] << "\n"; //} TS_ASSERT_LESS_THAN(l2_errors[0], 0.0005) TS_ASSERT_LESS_THAN(l2_errors[1], 5e-5); TS_ASSERT_LESS_THAN(l2_errors[2], 2.1e-6); TS_ASSERT_LESS_THAN(l2_errors[3], 4e-7); }
/* == Mechano-electric feedback, and alternative boundary conditions == * * Let us now run a simulation with mechano-electric feedback (MEF), and with different boundary conditions. */ void TestWithMef() throw(Exception) { /* If we want to use MEF, where the stretch (in the fibre-direction) couples back to the cell * model and is used in stretch-activated channels (SACs), we can't just let Chaste convert * from cellml to C++ code as usual (see electro-physiology tutorials on how cell model files * are autogenerated from CellML during compilation), since these files don't use stretch and don't * have SACs. We have to use pycml to create a cell model class for us, rename and save it, and * manually add the SAC current. * * There is one example of this already in the code-base, which we will use it the following * simulation. It is the Noble 98 model, with a SAC added that depends linearly on stretches (>1). * It is found in the file !NobleVargheseKohlNoble1998WithSac.hpp, and is called * `CML_noble_varghese_kohl_noble_1998_basic_with_sac`. * * To add a SAC current to (or otherwise alter) your favourite cell model, you have to do the following. * Auto-generate the C++ code, by running the following on the cellml file: * `./python/ConvertCellModel.py heart/src/odes/cellml/LuoRudy1991.cellml` * (see [wiki:ChasteGuides/CodeGenerationFromCellML#Usingthehelperscript ChasteGuides/CodeGenerationFromCellML#Usingthehelperscript] * if you want further documentation on this script). * * Copy and rename the resultant .hpp and .cpp files (which can be found in the same folder as the * input cellml). For example, rename everything to `LuoRudy1991WithSac`. Then alter the class * to overload the method `AbstractCardiacCell::SetStretch(double stretch)` to store the stretch, * and then implement the SAC in the `GetIIonic()` method. `CML_noble_varghese_kohl_noble_1998_basic_with_sac` * provides an example of the changes that need to be made. * * Let us create a cell factory returning these Noble98 SAC cells, but with no stimulus - the * SAC switching on will lead be to activation. */ ZeroStimulusCellFactory<CML_noble_varghese_kohl_noble_1998_basic_with_sac, 2> cell_factory; /* Construct two meshes are before, in 2D */ TetrahedralMesh<2,2> electrics_mesh; electrics_mesh.ConstructRegularSlabMesh(0.01/*stepsize*/, 0.1/*length*/, 0.1/*width*/, 0.1/*depth*/); QuadraticMesh<2> mechanics_mesh; mechanics_mesh.ConstructRegularSlabMesh(0.02, 0.1, 0.1, 0.1 /*as above with a different stepsize*/); /* Collect the fixed nodes. This time we directly specify the new locations. We say the * nodes on X=0 are to be fixed, setting the deformed x=0, but leaving y to be free * (sliding boundary conditions). This functionality is described in more detail in the * solid mechanics tutorials. */ std::vector<unsigned> fixed_nodes; std::vector<c_vector<double,2> > fixed_node_locations; fixed_nodes.push_back(0); fixed_node_locations.push_back(zero_vector<double>(2)); for(unsigned i=1; i<mechanics_mesh.GetNumNodes(); i++) { double X = mechanics_mesh.GetNode(i)->rGetLocation()[0]; if(fabs(X) < 1e-6) // ie, if X==0 { c_vector<double,2> new_position; new_position(0) = 0.0; new_position(1) = ElectroMechanicsProblemDefinition<2>::FREE; fixed_nodes.push_back(i); fixed_node_locations.push_back(new_position); } } /* Now specify tractions on the top and bottom surfaces. For full descriptions of how * to apply tractions see the solid mechanics tutorials. Here, we collect the boundary * elements on the bottom and top surfaces, and apply inward tractions - this will have the * effect of stretching the tissue in the X-direction. */ std::vector<BoundaryElement<1,2>*> boundary_elems; std::vector<c_vector<double,2> > tractions; c_vector<double,2> traction; for (TetrahedralMesh<2,2>::BoundaryElementIterator iter = mechanics_mesh.GetBoundaryElementIteratorBegin(); iter != mechanics_mesh.GetBoundaryElementIteratorEnd(); ++iter) { if (fabs((*iter)->CalculateCentroid()[1] - 0.0) < 1e-6) // if Y=0 { BoundaryElement<1,2>* p_element = *iter; boundary_elems.push_back(p_element); traction(0) = 0.0; // kPa, since the contraction model and material law use kPa for stiffnesses traction(1) = 1.0; // kPa, since the contraction model and material law use kPa for stiffnesses tractions.push_back(traction); } if (fabs((*iter)->CalculateCentroid()[1] - 0.1) < 1e-6) // if Y=0.1 { BoundaryElement<1,2>* p_element = *iter; boundary_elems.push_back(p_element); traction(0) = 0.0; traction(1) = -1.0; tractions.push_back(traction); } } /* Now set up the problem. We will use a compressible approach. */ ElectroMechanicsProblemDefinition<2> problem_defn(mechanics_mesh); problem_defn.SetContractionModel(KERCHOFFS2003,0.01/*contraction model ODE timestep*/); problem_defn.SetUseDefaultCardiacMaterialLaw(INCOMPRESSIBLE); problem_defn.SetMechanicsSolveTimestep(1.0); /* Set the fixed node and traction info. */ problem_defn.SetFixedNodes(fixed_nodes, fixed_node_locations); problem_defn.SetTractionBoundaryConditions(boundary_elems, tractions); /* Now say that the deformation should affect the electro-physiology */ problem_defn.SetDeformationAffectsElectrophysiology(false /*deformation affects conductivity*/, true /*deformation affects cell models*/); /* Set the end time, create the problem, and solve */ HeartConfig::Instance()->SetSimulationDuration(50.0); CardiacElectroMechanicsProblem<2,1> problem(INCOMPRESSIBLE, MONODOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "TestCardiacElectroMechanicsWithMef"); problem.Solve(); /* Nothing exciting happens in the simulation as it is currently written. To get some interesting occurring, * alter the SAC conductance in the cell model from 0.035 to 0.35 (micro-Siemens). * (look for the line `const double g_sac = 0.035` in `NobleVargheseKohlNoble1998WithSac.hpp`). * * Rerun and visualise as usual, using Cmgui. By visualising the voltage on the deforming mesh, you can see that the * voltage gradually increases due to the SAC, since the tissue is stretched, until the threshold is reached * and activation occurs. * * For MEF simulations, we may want to visualise the electrical results on the electrics mesh using * Meshalyzer, for example to more easily visualise action potentials. This isn't (and currently * can't be) created by `CardiacElectroMechanicsProblem`. We can use a converter as follows * to post-process: */ FileFinder test_output_folder("TestCardiacElectroMechanicsWithMef/electrics", RelativeTo::ChasteTestOutput); Hdf5ToMeshalyzerConverter<2,2> converter(test_output_folder, "voltage", &electrics_mesh, false, HeartConfig::Instance()->GetVisualizerOutputPrecision()); /* Some other notes. If you want to apply time-dependent traction boundary conditions, this is possible by * specifying the traction in functional form - see solid mechanics tutorials. Similarly, more natural * 'pressure acting on the deformed body' boundary conditions are possible - see below tutorial. * * '''Robustness:''' Sometimes the nonlinear solver doesn't converge, and will give an error. This can be due to either * a non-physical (or not very physical) situation, or just because the current guess is quite far * from the solution and the solver can't find the solution (this can occur in nonlinear elasticity * problems if the loading is large, for example). Current work in progress is on making the solver * more robust, and also on parallelising the solver. One option when a solve fails is to decrease the * mechanics timestep. */ /* Ignore the following, it is just to check nothing has changed. */ Hdf5DataReader reader("TestCardiacElectroMechanicsWithMef/electrics", "voltage"); unsigned num_timesteps = reader.GetUnlimitedDimensionValues().size(); Vec voltage = PetscTools::CreateVec(electrics_mesh.GetNumNodes()); reader.GetVariableOverNodes(voltage, "V", num_timesteps-1); ReplicatableVector voltage_repl(voltage); for(unsigned i=0; i<voltage_repl.GetSize(); i++) { TS_ASSERT_DELTA(voltage_repl[i], -81.9080, 1e-3); } PetscTools::Destroy(voltage); }
/* == 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(); }
/* == 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 TestDynamicExpansionNoElectroMechanics() throw (Exception) { TetrahedralMesh<2,2> electrics_mesh; QuadraticMesh<2> mechanics_mesh; // could (should?) use finer electrics mesh (if there was electrical activity occurring) // but keeping electrics simulation time down TrianglesMeshReader<2,2> reader1("mesh/test/data/annuli/circular_annulus_960_elements"); electrics_mesh.ConstructFromMeshReader(reader1); TrianglesMeshReader<2,2> reader2("mesh/test/data/annuli/circular_annulus_960_elements_quad",2 /*quadratic elements*/); mechanics_mesh.ConstructFromMeshReader(reader2); ZeroStimulusCellFactory<CellLuoRudy1991FromCellML,2> cell_factory; std::vector<unsigned> fixed_nodes; std::vector<c_vector<double,2> > fixed_node_locations; for(unsigned i=0; i<mechanics_mesh.GetNumNodes(); i++) { double x = mechanics_mesh.GetNode(i)->rGetLocation()[0]; double y = mechanics_mesh.GetNode(i)->rGetLocation()[1]; if (fabs(x)<1e-6 && fabs(y+0.5)<1e-6) // fixed point (0.0,-0.5) at bottom of mesh { fixed_nodes.push_back(i); c_vector<double,2> new_position; new_position(0) = x; new_position(1) = y; fixed_node_locations.push_back(new_position); } if (fabs(x)<1e-6 && fabs(y-0.5)<1e-6) // constrained point (0.0,0.5) at top of mesh { fixed_nodes.push_back(i); c_vector<double,2> new_position; new_position(0) = x; new_position(1) = ElectroMechanicsProblemDefinition<2>::FREE; fixed_node_locations.push_back(new_position); } } HeartConfig::Instance()->SetSimulationDuration(110.0); ElectroMechanicsProblemDefinition<2> problem_defn(mechanics_mesh); problem_defn.SetContractionModel(KERCHOFFS2003,0.1); problem_defn.SetUseDefaultCardiacMaterialLaw(COMPRESSIBLE); //problem_defn.SetZeroDisplacementNodes(fixed_nodes); problem_defn.SetFixedNodes(fixed_nodes, fixed_node_locations); problem_defn.SetMechanicsSolveTimestep(1.0); FileFinder finder("heart/test/data/fibre_tests/circular_annulus_960_elements.ortho", RelativeTo::ChasteSourceRoot); problem_defn.SetVariableFibreSheetDirectionsFile(finder, false); // The snes solver seems more robust... problem_defn.SetSolveUsingSnes(); problem_defn.SetVerboseDuringSolve(); // #2091 // This is a 2d problem, so a direct solve (LU factorisation) is possible and will speed things up // markedly (might be able to remove this line after #2057 is done..) //PetscTools::SetOption("-pc_type", "lu"); // removed - see comments at the end of #1818. std::vector<BoundaryElement<1,2>*> boundary_elems; for (TetrahedralMesh<2,2>::BoundaryElementIterator iter = mechanics_mesh.GetBoundaryElementIteratorBegin(); iter != mechanics_mesh.GetBoundaryElementIteratorEnd(); ++iter) { ChastePoint<2> centroid = (*iter)->CalculateCentroid(); double r = sqrt( centroid[0]*centroid[0] + centroid[1]*centroid[1] ); if (r < 0.4) { BoundaryElement<1,2>* p_element = *iter; boundary_elems.push_back(p_element); } } problem_defn.SetApplyNormalPressureOnDeformedSurface(boundary_elems, LinearPressureFunction); CardiacElectroMechanicsProblem<2,1> problem(COMPRESSIBLE, MONODOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "TestEmOnAnnulusDiastolicFilling"); problem.Solve(); // Hardcoded test of deformed position of (partially constrained) top node of the annulus, to check nothing has changed and that // different systems give the same result. TS_ASSERT_DELTA(problem.rGetDeformedPosition()[2](0), 0.0, 1e-15); TS_ASSERT_DELTA(problem.rGetDeformedPosition()[2](1), 0.5753, 1e-4); //Node 0 is on the righthand side, initially at x=0.5 y=0.0 TS_ASSERT_DELTA(problem.rGetDeformedPosition()[0](0), 0.5376, 1e-4); TS_ASSERT_DELTA(problem.rGetDeformedPosition()[0](1), 0.0376, 1e-4); MechanicsEventHandler::Headings(); MechanicsEventHandler::Report(); }
void TestStaticExpansionAndElectroMechanics() throw (Exception) { TetrahedralMesh<2,2> electrics_mesh; QuadraticMesh<2> mechanics_mesh; // could (should?) use finer electrics mesh, but keeping electrics simulation time down TrianglesMeshReader<2,2> reader1("mesh/test/data/annuli/circular_annulus_960_elements"); electrics_mesh.ConstructFromMeshReader(reader1); TrianglesMeshReader<2,2> reader2("mesh/test/data/annuli/circular_annulus_960_elements_quad",2 /*quadratic elements*/); mechanics_mesh.ConstructFromMeshReader(reader2); PointStimulus2dCellFactory cell_factory; std::vector<unsigned> fixed_nodes; std::vector<c_vector<double,2> > fixed_node_locations; for(unsigned i=0; i<mechanics_mesh.GetNumNodes(); i++) { double x = mechanics_mesh.GetNode(i)->rGetLocation()[0]; double y = mechanics_mesh.GetNode(i)->rGetLocation()[1]; if (fabs(x)<1e-6 && fabs(y+0.5)<1e-6) // fixed point (0.0,-0.5) at bottom of mesh { fixed_nodes.push_back(i); c_vector<double,2> new_position; new_position(0) = x; new_position(1) = y; fixed_node_locations.push_back(new_position); } if (fabs(x)<1e-6 && fabs(y-0.5)<1e-6) // constrained point (0.0,0.5) at top of mesh { fixed_nodes.push_back(i); c_vector<double,2> new_position; new_position(0) = x; new_position(1) = ElectroMechanicsProblemDefinition<2>::FREE; fixed_node_locations.push_back(new_position); } } // NOTE: sometimes the solver fails to converge (nonlinear solver failing to find a solution) during repolarisation. // This occurs with this test (on some configurations?) when the end time is set to 400ms. This needs to be // investigated - often surprisingly large numbers of newton iterations are needed in part of the repolarisation stage, and // then nonlinear solve failure occurs. Sometimes the solver can get past this if the mechanics timestep is decreased. // // See ticket #2193 HeartConfig::Instance()->SetSimulationDuration(400.0); ElectroMechanicsProblemDefinition<2> problem_defn(mechanics_mesh); problem_defn.SetContractionModel(KERCHOFFS2003,0.1); problem_defn.SetUseDefaultCardiacMaterialLaw(COMPRESSIBLE); //problem_defn.SetZeroDisplacementNodes(fixed_nodes); problem_defn.SetFixedNodes(fixed_nodes, fixed_node_locations); problem_defn.SetMechanicsSolveTimestep(1.0); FileFinder finder("heart/test/data/fibre_tests/circular_annulus_960_elements.ortho", RelativeTo::ChasteSourceRoot); problem_defn.SetVariableFibreSheetDirectionsFile(finder, false); // The snes solver seems more robust... problem_defn.SetSolveUsingSnes(); problem_defn.SetVerboseDuringSolve(); // #2091 // This is a 2d problem, so a direct solve (LU factorisation) is possible and will speed things up // markedly (might be able to remove this line after #2057 is done..) //PetscTools::SetOption("-pc_type", "lu"); // removed - see comments at the end of #1818. std::vector<BoundaryElement<1,2>*> boundary_elems; for (TetrahedralMesh<2,2>::BoundaryElementIterator iter = mechanics_mesh.GetBoundaryElementIteratorBegin(); iter != mechanics_mesh.GetBoundaryElementIteratorEnd(); ++iter) { ChastePoint<2> centroid = (*iter)->CalculateCentroid(); double r = sqrt( centroid[0]*centroid[0] + centroid[1]*centroid[1] ); if (r < 0.4) { BoundaryElement<1,2>* p_element = *iter; boundary_elems.push_back(p_element); } } problem_defn.SetApplyNormalPressureOnDeformedSurface(boundary_elems, -1.0 /*1 KPa is about 8mmHg*/); problem_defn.SetNumIncrementsForInitialDeformation(10); CardiacElectroMechanicsProblem<2,1> problem(COMPRESSIBLE, MONODOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "TestEmOnAnnulus"); problem.Solve(); // We don't really test anything, we mainly just want to verify it solves OK past the initial jump and through // the cycle. Have visualised. // Hardcoded test of deformed position of (partially constrained) top node of the annulus, to check nothing has changed and that // different systems give the same result. TS_ASSERT_DELTA(problem.rGetDeformedPosition()[2](0), 0.0, 1e-16); TS_ASSERT_DELTA(problem.rGetDeformedPosition()[2](1), 0.5753, 1e-4); //Node 0 is on the righthand side, initially at x=0.5 y=0.0 TS_ASSERT_DELTA(problem.rGetDeformedPosition()[0](0), 0.5376, 1e-4); TS_ASSERT_DELTA(problem.rGetDeformedPosition()[0](1), 0.0376, 1e-4); MechanicsEventHandler::Headings(); MechanicsEventHandler::Report(); }
// 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(); }
// This test is virtually identical to one of the the above tests (TestWithHomogeneousEverythingCompressible), // except it uses incompressible solid mechanics. Since the solution should be // x=alpha*X, y=beta*Y for some alpha, beta (see comments for above test), // we also test that alpha*beta = 1.0 void TestWithHomogeneousEverythingIncompressible() throw(Exception) { EntirelyStimulatedTissueCellFactory cell_factory; TetrahedralMesh<2,2> electrics_mesh; electrics_mesh.ConstructRegularSlabMesh(0.01, 0.05, 0.05); QuadraticMesh<2> mechanics_mesh; mechanics_mesh.ConstructRegularSlabMesh(0.025, 0.05, 0.05); std::vector<unsigned> fixed_nodes; std::vector<c_vector<double,2> > fixed_node_locations; fixed_nodes.push_back(0); fixed_node_locations.push_back(zero_vector<double>(2)); for(unsigned i=1 /*not 0*/; i<mechanics_mesh.GetNumNodes(); i++) { if(fabs(mechanics_mesh.GetNode(i)->rGetLocation()[0])<1e-6) { c_vector<double,2> new_position; new_position(0) = 0.0; new_position(1) = SolidMechanicsProblemDefinition<2>::FREE; fixed_nodes.push_back(i); fixed_node_locations.push_back(new_position); } } ElectroMechanicsProblemDefinition<2> problem_defn(mechanics_mesh); problem_defn.SetContractionModel(KERCHOFFS2003,1.0); problem_defn.SetUseDefaultCardiacMaterialLaw(INCOMPRESSIBLE); problem_defn.SetFixedNodes(fixed_nodes, fixed_node_locations); problem_defn.SetMechanicsSolveTimestep(1.0); // coverage, this file is just X-direction fibres FileFinder fibre_file("heart/test/data/fibre_tests/2by2mesh_fibres.ortho", RelativeTo::ChasteSourceRoot); problem_defn.SetVariableFibreSheetDirectionsFile(fibre_file, false); HeartConfig::Instance()->SetSimulationDuration(10.0); CardiacElectroMechanicsProblem<2,1> problem(INCOMPRESSIBLE, MONODOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "TestCardiacEmHomogeneousEverythingIncompressible"); problem.Solve(); std::vector<c_vector<double,2> >& r_deformed_position = problem.rGetDeformedPosition(); // not sure how easy is would be determine what the deformation should be // exactly, but it certainly should be constant squash in X direction, constant // stretch in Y. // first, check node 8 starts is the far corner assert(fabs(mechanics_mesh.GetNode(8)->rGetLocation()[0] - 0.05)<1e-8); assert(fabs(mechanics_mesh.GetNode(8)->rGetLocation()[1] - 0.05)<1e-8); double X_scale_factor = r_deformed_position[8](0)/0.05; double Y_scale_factor = r_deformed_position[8](1)/0.05; std::cout << "Scale_factors = " << X_scale_factor << " " << Y_scale_factor << ", product = " << X_scale_factor*Y_scale_factor<<"\n"; for(unsigned i=0; i<mechanics_mesh.GetNumNodes(); i++) { double X = mechanics_mesh.GetNode(i)->rGetLocation()[0]; double Y = mechanics_mesh.GetNode(i)->rGetLocation()[1]; TS_ASSERT_DELTA( r_deformed_position[i](0), X * X_scale_factor, 1e-6); TS_ASSERT_DELTA( r_deformed_position[i](1), Y * Y_scale_factor, 1e-6); } // here, we also test the deformation is incompressible TS_ASSERT_DELTA(X_scale_factor * Y_scale_factor, 1.0, 1e-6); }
//This test is identical to the test above, just that we use a model that requires an explicit solver void TestMechanicsWithBidomainAndBathExplicit() throw(Exception) { EntirelyStimulatedTissueCellFactory cell_factory; TetrahedralMesh<2,2> electrics_mesh; electrics_mesh.ConstructRegularSlabMesh(0.01, 0.05, 0.05); //make everything a bath node except for the x=0 line for (TetrahedralMesh<2,2>::ElementIterator iter=electrics_mesh.GetElementIteratorBegin(); iter != electrics_mesh.GetElementIteratorEnd(); ++iter) { if ((*iter).CalculateCentroid()[0] > 0.001) { (*iter).SetAttribute(HeartRegionCode::GetValidBathId()); } else { (*iter).SetAttribute(HeartRegionCode::GetValidTissueId()); } } QuadraticMesh<2> mechanics_mesh; mechanics_mesh.ConstructRegularSlabMesh(0.025, 0.05, 0.05); //store the original node positions std::vector<c_vector<double,2> > original_node_position; c_vector<double,2> pos = zero_vector<double>(2); for(unsigned i=0; i<mechanics_mesh.GetNumNodes(); i++) { pos(0) = mechanics_mesh.GetNode(i)->rGetLocation()[0]; pos(1) = mechanics_mesh.GetNode(i)->rGetLocation()[1]; original_node_position.push_back(pos); } std::vector<unsigned> fixed_nodes; std::vector<c_vector<double,2> > fixed_node_locations; // fix the node at the origin so that the solution is well-defined (ie unique) fixed_nodes.push_back(0); fixed_node_locations.push_back(zero_vector<double>(2)); // for the rest of the nodes, if they lie on X=0, fix x=0 but leave y free. for(unsigned i=1 /*not 0*/; i<mechanics_mesh.GetNumNodes(); i++) { if(fabs(mechanics_mesh.GetNode(i)->rGetLocation()[0])<1e-6) { c_vector<double,2> new_position; new_position(0) = 0.0; new_position(1) = SolidMechanicsProblemDefinition<2>::FREE; fixed_nodes.push_back(i); fixed_node_locations.push_back(new_position); } } ElectroMechanicsProblemDefinition<2> problem_defn(mechanics_mesh); problem_defn.SetContractionModel(NASH2004,1.0); problem_defn.SetUseDefaultCardiacMaterialLaw(INCOMPRESSIBLE); problem_defn.SetFixedNodes(fixed_nodes, fixed_node_locations); HeartConfig::Instance()->SetOdePdeAndPrintingTimeSteps(0.01,0.1,1.0); problem_defn.SetMechanicsSolveTimestep(1.0); HeartConfig::Instance()->SetSimulationDuration(10.0); HeartConfig::Instance()->SetExtracellularConductivities(Create_c_vector(1500,1500,1500)); CardiacElectroMechanicsProblem<2,2> problem(INCOMPRESSIBLE, BIDOMAIN_WITH_BATH, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "TestCardiacEmWithBath"); problem.Solve(); std::vector<c_vector<double,2> >& r_deformed_position = problem.rGetDeformedPosition(); // first, check node 8 starts is the far corner assert(fabs(mechanics_mesh.GetNode(8)->rGetLocation()[0] - 0.05)<1e-8); assert(fabs(mechanics_mesh.GetNode(8)->rGetLocation()[1] - 0.05)<1e-8); for(unsigned i=0; i<mechanics_mesh.GetNumNodes(); i++) { TS_ASSERT_DELTA( r_deformed_position[i](0), original_node_position[i](0), 1e-6); TS_ASSERT_DELTA( r_deformed_position[i](1), original_node_position[i](1), 1e-6); } }
// 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); }
// runs 5 different 3d tests with fibres read to be in various directions, and // checks contraction occurs in the right directions (and bulging occurs in the // other directions) void TestFibreRead() throw(Exception) { PlaneStimulusCellFactory<CellLuoRudy1991FromCellML,3> cell_factory(-5000*1000); double tissue_initial_size = 0.05; TetrahedralMesh<3,3> electrics_mesh; electrics_mesh.ConstructRegularSlabMesh(0.01/*stepsize*/, tissue_initial_size/*length*/, tissue_initial_size/*width*/, tissue_initial_size); QuadraticMesh<3> mechanics_mesh; mechanics_mesh.ConstructRegularSlabMesh(0.05, tissue_initial_size, tissue_initial_size, tissue_initial_size /*as above with a different stepsize*/); std::vector<unsigned> fixed_nodes = NonlinearElasticityTools<3>::GetNodesByComponentValue(mechanics_mesh, 2, 0.0); HeartConfig::Instance()->SetSimulationDuration(10.0); ElectroMechanicsProblemDefinition<3> problem_defn(mechanics_mesh); problem_defn.SetContractionModel(KERCHOFFS2003,1.0); problem_defn.SetZeroDisplacementNodes(fixed_nodes); problem_defn.SetMechanicsSolveTimestep(1.0); problem_defn.SetUseDefaultCardiacMaterialLaw(COMPRESSIBLE); // This is how the fibres files that are used in the simulations are created. // alongX.ortho is: fibres in X axis, sheet in XY // -- same as default directions: should give shortening in X direction // and identical results to default fibre setting // alongY1.ortho is: fibres in Y axis, sheet in YX // alongY2.ortho is: fibres in Y axis, sheet in YZ // -- as material law is transversely isotropic, these two should // give identical results, shortening in Y-direction // alongZ.ortho is: fibres in Z axis, (sheet in XZ) // -- should shorten in Z-direction OutputFileHandler handler("FibreFiles"); out_stream p_X_file = handler.OpenOutputFile("alongX.ortho"); out_stream p_Y1_file = handler.OpenOutputFile("alongY1.ortho"); out_stream p_Y2_file = handler.OpenOutputFile("alongY2.ortho"); out_stream p_Z_file = handler.OpenOutputFile("alongZ.ortho"); *p_X_file << mechanics_mesh.GetNumElements() << "\n"; *p_Y1_file << mechanics_mesh.GetNumElements() << "\n"; *p_Y2_file << mechanics_mesh.GetNumElements() << "\n"; *p_Z_file << mechanics_mesh.GetNumElements() << "\n"; for(unsigned i=0; i<mechanics_mesh.GetNumElements(); i++) { //double X = mechanics_mesh.GetElement(i)->CalculateCentroid()(0); *p_X_file << "1 0 0 0 1 0 0 0 1\n"; *p_Y1_file << "0 1 0 1 0 0 0 0 1\n"; *p_Y2_file << "0 1 0 0 0 1 1 0 0\n"; *p_Z_file << "0 0 1 1 0 0 0 1 0\n"; } p_X_file->close(); p_Y1_file->close(); p_Y2_file->close(); p_Z_file->close(); ////////////////////////////////////////////////////////////////// // Solve with no fibres read. ////////////////////////////////////////////////////////////////// std::vector<c_vector<double,3> > r_deformed_position_no_fibres; { HeartConfig::Instance()->SetSimulationDuration(20.0); CardiacElectroMechanicsProblem<3,1> problem(COMPRESSIBLE, MONODOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "TestCardiacEmFibreRead"); problem.Solve(); r_deformed_position_no_fibres = problem.rGetDeformedPosition(); } ////////////////////////////////////////////////////////////////// // Solve with fibres read: fibres in X-direction ////////////////////////////////////////////////////////////////// std::vector<c_vector<double,3> > r_deformed_position_fibres_alongX; { HeartConfig::Instance()->SetSimulationDuration(20.0); FileFinder finder("heart/test/data/fibre_tests/alongX.ortho", RelativeTo::ChasteSourceRoot); problem_defn.SetVariableFibreSheetDirectionsFile(finder, false); CardiacElectroMechanicsProblem<3,1> problem(COMPRESSIBLE, MONODOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "TestCardiacEmFibreRead"); problem.Solve(); r_deformed_position_fibres_alongX = problem.rGetDeformedPosition(); } // test the two results are identical for(unsigned i=0; i<r_deformed_position_no_fibres.size(); i++) { for(unsigned j=0; j<3; j++) { TS_ASSERT_DELTA(r_deformed_position_no_fibres[i](j), r_deformed_position_fibres_alongX[i](j), 1e-8); } } // check node 7 starts is the far corner assert(fabs(mechanics_mesh.GetNode(7)->rGetLocation()[0] - tissue_initial_size)<1e-8); assert(fabs(mechanics_mesh.GetNode(7)->rGetLocation()[1] - tissue_initial_size)<1e-8); assert(fabs(mechanics_mesh.GetNode(7)->rGetLocation()[2] - tissue_initial_size)<1e-8); // Test that contraction occurred in the X-direction TS_ASSERT_LESS_THAN(r_deformed_position_fibres_alongX[7](0), tissue_initial_size); TS_ASSERT_LESS_THAN(tissue_initial_size, r_deformed_position_fibres_alongX[7](1)); TS_ASSERT_LESS_THAN(tissue_initial_size, r_deformed_position_fibres_alongX[7](2)); // hardcoded test to check nothing changes TS_ASSERT_DELTA(r_deformed_position_fibres_alongX[7](0), 0.0487, 1e-4); TS_ASSERT_DELTA(r_deformed_position_fibres_alongX[7](1), 0.0506, 1e-4); //////////////////////////////////////////////////////////////////// // Solve with fibres read: fibres in Y-direction, sheet in YX plane //////////////////////////////////////////////////////////////////// std::vector<c_vector<double,3> > r_deformed_position_fibres_alongY1; { HeartConfig::Instance()->SetSimulationDuration(20.0); FileFinder finder("heart/test/data/fibre_tests/alongY1.ortho", RelativeTo::ChasteSourceRoot); problem_defn.SetVariableFibreSheetDirectionsFile(finder, false); CardiacElectroMechanicsProblem<3,1> problem(COMPRESSIBLE, MONODOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "TestCardiacEmFibreRead"); problem.Solve(); r_deformed_position_fibres_alongY1 = problem.rGetDeformedPosition(); } //////////////////////////////////////////////////////////////////// // Solve with fibres read: fibres in Y-direction, sheet in YZ plane //////////////////////////////////////////////////////////////////// std::vector<c_vector<double,3> > r_deformed_position_fibres_alongY2; { HeartConfig::Instance()->SetSimulationDuration(20.0); FileFinder fibres_file("heart/test/data/fibre_tests/alongY2.ortho", RelativeTo::ChasteSourceRoot); problem_defn.SetVariableFibreSheetDirectionsFile(fibres_file, false); CardiacElectroMechanicsProblem<3,1> problem(COMPRESSIBLE, MONODOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "TestCardiacEmFibreRead"); problem.Solve(); r_deformed_position_fibres_alongY2 = problem.rGetDeformedPosition(); } // test the two results are identical for(unsigned i=0; i<r_deformed_position_no_fibres.size(); i++) { for(unsigned j=0; j<3; j++) { TS_ASSERT_DELTA(r_deformed_position_fibres_alongY1[i](j), r_deformed_position_fibres_alongY2[i](j), 1e-8); } } // Test that contraction occurred in the Y-direction TS_ASSERT_LESS_THAN(tissue_initial_size, r_deformed_position_fibres_alongY1[7](0)); TS_ASSERT_LESS_THAN(r_deformed_position_fibres_alongY1[7](1), tissue_initial_size); TS_ASSERT_LESS_THAN(tissue_initial_size, r_deformed_position_fibres_alongY1[7](2)); // hardcoded test to check nothing changes TS_ASSERT_DELTA(r_deformed_position_fibres_alongY1[7](1), 0.0487, 1e-4); TS_ASSERT_DELTA(r_deformed_position_fibres_alongY1[7](0), 0.0506, 1e-4); ////////////////////////////////////////////////////////////////// // Solve with fibres read: fibres in Z-direction ////////////////////////////////////////////////////////////////// std::vector<c_vector<double,3> > r_deformed_position_fibres_alongZ; { HeartConfig::Instance()->SetSimulationDuration(20.0); FileFinder finder("heart/test/data/fibre_tests/alongZ.ortho", RelativeTo::ChasteSourceRoot); problem_defn.SetVariableFibreSheetDirectionsFile(finder, false); CardiacElectroMechanicsProblem<3,1> problem(COMPRESSIBLE, MONODOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "TestCardiacEmFibreRead"); problem.Solve(); r_deformed_position_fibres_alongZ = problem.rGetDeformedPosition(); } // Test that contraction occurred in the X-direction TS_ASSERT_LESS_THAN(tissue_initial_size, r_deformed_position_fibres_alongZ[7](0)); TS_ASSERT_LESS_THAN(tissue_initial_size, r_deformed_position_fibres_alongZ[7](1)); TS_ASSERT_LESS_THAN(r_deformed_position_fibres_alongZ[7](2), tissue_initial_size); // hardcoded test to check nothing changes TS_ASSERT_DELTA(r_deformed_position_fibres_alongZ[7](2), 0.0466, 1e-4); TS_ASSERT_DELTA(r_deformed_position_fibres_alongZ[7](0), 0.0504, 1e-4); }
/* == 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); }
// Test the functionality specific to SolidMechanicsProblemDefinition void TestSolidMechanicsProblemDefinition() throw(Exception) { TS_ASSERT_EQUALS(SolidMechanicsProblemDefinition<2>::FREE, DBL_MAX); TS_ASSERT_LESS_THAN(0, SolidMechanicsProblemDefinition<2>::FREE); QuadraticMesh<2> mesh(0.5, 1.0, 1.0); SolidMechanicsProblemDefinition<2> problem_defn(mesh); ////////////////////////////////// // Fixed nodes ////////////////////////////////// std::vector<unsigned> fixed_nodes; fixed_nodes.push_back(0); fixed_nodes.push_back(4); problem_defn.SetZeroDisplacementNodes(fixed_nodes); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes().size(), 2u); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes()[0], 0u); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes()[1], 4u); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodeValues().size(), 2u); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[0](0), 0.0, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[0](1), 0.0, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[1](0), 0.0, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[1](1), 0.0, 1e-12); fixed_nodes.push_back(8); fixed_nodes.push_back(9); fixed_nodes.push_back(10); std::vector<c_vector<double,2> > locations; c_vector<double,2> location = zero_vector<double>(2); // Node 0 is to be placed at (0,0) locations.push_back(location); // Node 4 is to be placed at (0,0.1) location(1)=0.1; locations.push_back(location); // Node 8 is to be placed at (0.1,0.1) location(0)=0.1; locations.push_back(location); // Node 9 is to be placed at (0.5,FREE) location(0) = 0.5; location(1) = SolidMechanicsProblemDefinition<2>::FREE; locations.push_back(location); // Node 9 is to be placed at (FREE,1.5) location(0) = SolidMechanicsProblemDefinition<2>::FREE; location(1) = 1.5; locations.push_back(location); problem_defn.SetFixedNodes(fixed_nodes, locations); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes().size(), 5u); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodeValues().size(), 5u); // the fully fixed nodes TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes()[0], 0u); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes()[1], 4u); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes()[2], 8u); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[0](0), 0.0 - mesh.GetNode(0)->rGetLocation()[0], 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[0](1), 0.0 - mesh.GetNode(0)->rGetLocation()[1], 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[1](0), 0.0 - mesh.GetNode(4)->rGetLocation()[0], 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[1](1), 0.1 - mesh.GetNode(4)->rGetLocation()[1], 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[2](0), 0.1 - mesh.GetNode(8)->rGetLocation()[0], 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[2](1), 0.1 - mesh.GetNode(8)->rGetLocation()[1], 1e-12); // the partial fixed nodes TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes()[3], 9u); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes()[4], 10u); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[3](0), 0.5 - mesh.GetNode(9)->rGetLocation()[0], 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[3](1), SolidMechanicsProblemDefinition<2>::FREE, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[4](0), SolidMechanicsProblemDefinition<2>::FREE, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[4](1), 1.5 - mesh.GetNode(10)->rGetLocation()[1], 1e-12); /////////////////////////////////////// // Set an incompressible material law /////////////////////////////////////// TS_ASSERT_THROWS_THIS(problem_defn.Validate(), "No material law has been set"); // set a homogeneous law MooneyRivlinMaterialLaw<2> incomp_mooney_rivlin_law(1.0); problem_defn.SetMaterialLaw(INCOMPRESSIBLE,&incomp_mooney_rivlin_law); TS_ASSERT_EQUALS(problem_defn.IsHomogeneousMaterial(), true); TS_ASSERT_EQUALS(problem_defn.GetCompressibilityType(), INCOMPRESSIBLE); TS_ASSERT_EQUALS(problem_defn.GetIncompressibleMaterialLaw(0), &incomp_mooney_rivlin_law); // set a heterogeneous law MooneyRivlinMaterialLaw<2> incomp_mooney_rivlin_law_2(2.0); std::vector<AbstractMaterialLaw<2>*> laws; for(unsigned i=0; i<mesh.GetNumElements()/2; i++) { laws.push_back(&incomp_mooney_rivlin_law); } for(unsigned i=mesh.GetNumElements()/2; i<mesh.GetNumElements(); i++) { laws.push_back(&incomp_mooney_rivlin_law_2); } problem_defn.SetMaterialLaw(INCOMPRESSIBLE,laws); TS_ASSERT_EQUALS(problem_defn.IsHomogeneousMaterial(), false); for(unsigned i=0; i<mesh.GetNumElements()/2; i++) { TS_ASSERT_EQUALS(problem_defn.GetIncompressibleMaterialLaw(i), &incomp_mooney_rivlin_law); } for(unsigned i=mesh.GetNumElements()/2; i<mesh.GetNumElements(); i++) { TS_ASSERT_EQUALS(problem_defn.GetIncompressibleMaterialLaw(i), &incomp_mooney_rivlin_law_2); } ///////////////////////////////////////////////////////////////////////// // Set a compressible material law (clears the incompressible laws) ///////////////////////////////////////////////////////////////////////// CompressibleMooneyRivlinMaterialLaw<2> comp_mooney_rivlin_law(2.0, 1.0); problem_defn.SetMaterialLaw(COMPRESSIBLE,&comp_mooney_rivlin_law); TS_ASSERT_EQUALS(problem_defn.IsHomogeneousMaterial(), true); TS_ASSERT_EQUALS(problem_defn.GetCompressibilityType(), COMPRESSIBLE); TS_ASSERT_EQUALS(problem_defn.GetCompressibleMaterialLaw(0), &comp_mooney_rivlin_law); // set a heterogeneous law CompressibleMooneyRivlinMaterialLaw<2> comp_mooney_rivlin_law_2(4.0, 1.0); std::vector<AbstractMaterialLaw<2>*> comp_laws; for(unsigned i=0; i<mesh.GetNumElements()/2; i++) { comp_laws.push_back(&comp_mooney_rivlin_law); } for(unsigned i=mesh.GetNumElements()/2; i<mesh.GetNumElements(); i++) { comp_laws.push_back(&comp_mooney_rivlin_law_2); } problem_defn.SetMaterialLaw(COMPRESSIBLE,comp_laws); TS_ASSERT_EQUALS(problem_defn.IsHomogeneousMaterial(), false); for(unsigned i=0; i<mesh.GetNumElements()/2; i++) { TS_ASSERT_EQUALS(problem_defn.GetCompressibleMaterialLaw(i), &comp_mooney_rivlin_law); } for(unsigned i=mesh.GetNumElements()/2; i<mesh.GetNumElements(); i++) { TS_ASSERT_EQUALS(problem_defn.GetCompressibleMaterialLaw(i), &comp_mooney_rivlin_law_2); } // should not throw anything problem_defn.Validate(); TS_ASSERT_THROWS_THIS(problem_defn.SetMaterialLaw(INCOMPRESSIBLE,&comp_mooney_rivlin_law),"Compressibility type was declared as INCOMPRESSIBLE but a compressible material law was given"); TS_ASSERT_THROWS_THIS(problem_defn.SetMaterialLaw(COMPRESSIBLE,&incomp_mooney_rivlin_law),"Incompressibility type was declared as COMPRESSIBLE but an incompressible material law was given"); /////////////////////////////// // solver stuff /////////////////////////////// TS_ASSERT_EQUALS(problem_defn.GetSolveUsingSnes(), false); TS_ASSERT_EQUALS(problem_defn.GetVerboseDuringSolve(), false); problem_defn.SetSolveUsingSnes(); problem_defn.SetVerboseDuringSolve(); TS_ASSERT_EQUALS(problem_defn.GetSolveUsingSnes(), true); TS_ASSERT_EQUALS(problem_defn.GetVerboseDuringSolve(), true); problem_defn.SetSolveUsingSnes(false); problem_defn.SetVerboseDuringSolve(false); TS_ASSERT_EQUALS(problem_defn.GetSolveUsingSnes(), false); TS_ASSERT_EQUALS(problem_defn.GetVerboseDuringSolve(), false); }
/* == Internal pressures == * * Next, we run a simulation on a 2d annulus, with an internal pressure applied. */ void TestAnnulusWithInternalPressure() throw (Exception) { /* The following should require little explanation now */ TetrahedralMesh<2,2> electrics_mesh; QuadraticMesh<2> mechanics_mesh; // could (should?) use finer electrics mesh, but keeping electrics simulation time down TrianglesMeshReader<2,2> reader1("mesh/test/data/annuli/circular_annulus_960_elements"); electrics_mesh.ConstructFromMeshReader(reader1); TrianglesMeshReader<2,2> reader2("mesh/test/data/annuli/circular_annulus_960_elements_quad",2 /*quadratic elements*/); mechanics_mesh.ConstructFromMeshReader(reader2); PointStimulus2dCellFactory cell_factory; std::vector<unsigned> fixed_nodes; std::vector<c_vector<double,2> > fixed_node_locations; for(unsigned i=0; i<mechanics_mesh.GetNumNodes(); i++) { double x = mechanics_mesh.GetNode(i)->rGetLocation()[0]; double y = mechanics_mesh.GetNode(i)->rGetLocation()[1]; if (fabs(x)<1e-6 && fabs(y+0.5)<1e-6) // fixed point (0.0,-0.5) at bottom of mesh { fixed_nodes.push_back(i); c_vector<double,2> new_position; new_position(0) = x; new_position(1) = y; fixed_node_locations.push_back(new_position); } if (fabs(x)<1e-6 && fabs(y-0.5)<1e-6) // constrained point (0.0,0.5) at top of mesh { fixed_nodes.push_back(i); c_vector<double,2> new_position; new_position(0) = x; new_position(1) = ElectroMechanicsProblemDefinition<2>::FREE; fixed_node_locations.push_back(new_position); } } /* Increase this end time to see more contraction */ HeartConfig::Instance()->SetSimulationDuration(30.0); ElectroMechanicsProblemDefinition<2> problem_defn(mechanics_mesh); problem_defn.SetContractionModel(KERCHOFFS2003,0.1); problem_defn.SetUseDefaultCardiacMaterialLaw(COMPRESSIBLE); //problem_defn.SetZeroDisplacementNodes(fixed_nodes); problem_defn.SetFixedNodes(fixed_nodes, fixed_node_locations); problem_defn.SetMechanicsSolveTimestep(1.0); FileFinder finder("heart/test/data/fibre_tests/circular_annulus_960_elements.ortho",RelativeTo::ChasteSourceRoot); problem_defn.SetVariableFibreSheetDirectionsFile(finder, false); /* 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(); /* Now let us collect all the boundary elements on the inner (endocardial) surface. The following * uses knowledge about the geometry - the inner surface is r=0.3, the outer is r=0.5. */ std::vector<BoundaryElement<1,2>*> boundary_elems; for (TetrahedralMesh<2,2>::BoundaryElementIterator iter = mechanics_mesh.GetBoundaryElementIteratorBegin(); iter != mechanics_mesh.GetBoundaryElementIteratorEnd(); ++iter) { ChastePoint<2> centroid = (*iter)->CalculateCentroid(); double r = sqrt( centroid[0]*centroid[0] + centroid[1]*centroid[1] ); if (r < 0.4) { BoundaryElement<1,2>* p_element = *iter; boundary_elems.push_back(p_element); } } /* This is how to set the pressure to be applied to these boundary elements. The negative sign implies * inward pressure. */ problem_defn.SetApplyNormalPressureOnDeformedSurface(boundary_elems, -1.0 /*1 KPa is about 8mmHg*/); /* The solver computes the equilibrium solution (given the pressure loading) before the first timestep. * As there is a big deformation from the undeformed state to this loaded state, the nonlinear solver may * not converge. The following increments the loading (solves with p=-1/3, then p=-2/3, then p=-1), which * allows convergence to occur. */ problem_defn.SetNumIncrementsForInitialDeformation(3); CardiacElectroMechanicsProblem<2,1> problem(COMPRESSIBLE, MONODOMAIN, &electrics_mesh, &mechanics_mesh, &cell_factory, &problem_defn, "TestAnnulusWithInternalPressure"); /* If we want stresses and strains output, we can do the following. The deformation gradients and 2nd PK stresses * for each element will be written at the requested times. */ problem.SetOutputDeformationGradientsAndStress(10.0 /*how often (in ms) to write - should be a multiple of mechanics timestep*/); /* Since this test involves a large deformation at t=0, several Newton iterations are required. To see how the nonlinear * solve is progressing, you can run from the binary from the command line with the command line argument "-mech_verbose". */ problem.Solve(); /* Visualise using cmgui, and note the different shapes at t=-1 (undeformed) and t=0 (loaded) * * Note: if you want to have a time-dependent pressure, you can replace the second parameter (the pressure) * in `SetApplyNormalPressureOnDeformedSurface()` with a function pointer (the name of a function) which returns * the pressure as a function of time. */ }
void TestStokesFlowProblemDefinition() throw(Exception) { QuadraticMesh<2> mesh(0.5, 1.0, 1.0); StokesFlowProblemDefinition<2> problem_defn(mesh); TS_ASSERT_THROWS_THIS(problem_defn.GetViscosity(), "Viscosity hasn't been set yet (for the Stokes' flow problem)"); problem_defn.SetViscosity(1.3423423); TS_ASSERT_EQUALS(problem_defn.GetViscosity(),1.3423423); ////////////////////////////////// // Fixed nodes ////////////////////////////////// std::vector<unsigned> fixed_flow_nodes; fixed_flow_nodes.push_back(0); fixed_flow_nodes.push_back(4); problem_defn.SetZeroFlowNodes(fixed_flow_nodes); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes().size(), 2u); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes()[0], 0u); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes()[1], 4u); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodeValues().size(), 2u); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[0](0), 0.0, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[0](1), 0.0, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[1](0), 0.0, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[1](1), 0.0, 1e-12); fixed_flow_nodes.push_back(8); fixed_flow_nodes.push_back(9); fixed_flow_nodes.push_back(10); std::vector<c_vector<double,2> > fixed_flows; c_vector<double,2> flow = zero_vector<double>(2); fixed_flows.push_back(flow); flow(1)=0.1; fixed_flows.push_back(flow); flow(0)=0.1; fixed_flows.push_back(flow); flow(0) = 0.5; flow(1) = StokesFlowProblemDefinition<2>::FREE; fixed_flows.push_back(flow); flow(0) = StokesFlowProblemDefinition<2>::FREE; flow(1) = 1.5; fixed_flows.push_back(flow); problem_defn.SetPrescribedFlowNodes(fixed_flow_nodes, fixed_flows); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes().size(), 5u); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodeValues().size(), 5u); // the fully fixed nodes TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes()[0], 0u); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes()[1], 4u); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes()[2], 8u); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[0](0), 0.0, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[0](1), 0.0, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[1](0), 0.0, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[1](1), 0.1, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[2](0), 0.1, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[2](1), 0.1, 1e-12); // the partial fixed nodes TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes()[3], 9u); TS_ASSERT_EQUALS(problem_defn.rGetDirichletNodes()[4], 10u); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[3](0), 0.5, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[3](1), StokesFlowProblemDefinition<2>::FREE, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[4](0), StokesFlowProblemDefinition<2>::FREE, 1e-12); TS_ASSERT_DELTA(problem_defn.rGetDirichletNodeValues()[4](1), 1.5, 1e-12); // should not throw anything problem_defn.Validate(); }