Chaste  Release::3.4
AbstractNonlinearElasticitySolver.hpp
1 /*
2 
3 Copyright (c) 2005-2016, University of Oxford.
4 All rights reserved.
5 
6 University of Oxford means the Chancellor, Masters and Scholars of the
7 University of Oxford, having an administrative office at Wellington
8 Square, Oxford OX1 2JD, UK.
9 
10 This file is part of Chaste.
11 
12 Redistribution and use in source and binary forms, with or without
13 modification, are permitted provided that the following conditions are met:
14  * Redistributions of source code must retain the above copyright notice,
15  this list of conditions and the following disclaimer.
16  * Redistributions in binary form must reproduce the above copyright notice,
17  this list of conditions and the following disclaimer in the documentation
18  and/or other materials provided with the distribution.
19  * Neither the name of the University of Oxford nor the names of its
20  contributors may be used to endorse or promote products derived from this
21  software without specific prior written permission.
22 
23 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
24 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
25 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
26 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
27 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
28 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
29 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
30 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
31 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
32 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
33 
34 */
35 
36 #ifndef ABSTRACTNONLINEARELASTICITYSOLVER_HPP_
37 #define ABSTRACTNONLINEARELASTICITYSOLVER_HPP_
38 
39 #include <vector>
40 #include <cmath>
41 #include "AbstractContinuumMechanicsSolver.hpp"
42 #include "LinearSystem.hpp"
43 #include "LogFile.hpp"
44 #include "MechanicsEventHandler.hpp"
45 #include "ReplicatableVector.hpp"
46 #include "FourthOrderTensor.hpp"
47 #include "CmguiDeformedSolutionsWriter.hpp"
48 #include "AbstractMaterialLaw.hpp"
49 #include "QuadraticBasisFunction.hpp"
50 #include "SolidMechanicsProblemDefinition.hpp"
51 #include "Timer.hpp"
52 #include "AbstractPerElementWriter.hpp"
53 #include "petscsnes.h"
54 
55 
56 
57 //#define MECH_USE_HYPRE // uses HYPRE to solve linear systems, requires PETSc to be installed with HYPRE
58 
59 
64 typedef enum StrainType_
65 {
66  DEFORMATION_GRADIENT_F = 0,
67  DEFORMATION_TENSOR_C,
68  LAGRANGE_STRAIN_E
69 } StrainType;
70 
71 
72 
73 // Bizarrely PETSc 2.2 has this, but doesn't put it in the petscksp.h header...
74 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
75 extern PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec);
76 #endif
77 
78 
80 // Globals functions used by the SNES solver
82 
92 template<unsigned DIM>
93 PetscErrorCode AbstractNonlinearElasticitySolver_ComputeResidual(SNES snes,
94  Vec currentGuess,
95  Vec residualVector,
96  void* pContext);
97 
98 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
99 
111  template<unsigned DIM>
112  PetscErrorCode AbstractNonlinearElasticitySolver_ComputeJacobian(SNES snes,
113  Vec currentGuess,
114  Mat globalJacobian,
115  Mat preconditioner,
116  void* pContext);
117 #else
118 
129 template<unsigned DIM>
130 PetscErrorCode AbstractNonlinearElasticitySolver_ComputeJacobian(SNES snes,
131  Vec currentGuess,
132  Mat* pGlobalJacobian,
133  Mat* pPreconditioner,
134  MatStructure* pMatStructure,
135  void* pContext);
136 #endif
137 
138 template <unsigned DIM>
139 class AbstractNonlinearElasticitySolver; //Forward declaration
140 
147 template<unsigned DIM>
148 class StressPerElementWriter : public AbstractPerElementWriter<DIM, DIM, DIM*DIM>
149 {
150 private:
152 public:
153 
160  : AbstractPerElementWriter<DIM, DIM, DIM*DIM>(pMesh),
161  mpSolver(pSolver)
162  {
163  }
164 
172  void Visit(Element<DIM, DIM>* pElement, unsigned localIndex, c_vector<double, DIM*DIM>& rData)
173  {
177  assert(localIndex == pElement->GetIndex()); //Will fail when we move to DistributedQuadraticMesh
178  //Flatten the matrix
179  c_matrix<double, DIM, DIM> data = mpSolver->GetAverageStressPerElement(localIndex);
180  for (unsigned i=0; i<DIM; i++)
181  {
182  for (unsigned j=0; j<DIM; j++)
183  {
184  rData[i*DIM+j] = data(i,j);
185  }
186  }
187  }
188 };
189 
201 template<unsigned DIM>
203 {
204  friend class StressRecoveror<DIM>;
205  friend class VtkNonlinearElasticitySolutionWriter<DIM>;
206 
207 
208 protected:
209 
211  static const size_t NUM_VERTICES_PER_ELEMENT = DIM+1;
212 
214  static const size_t NUM_NODES_PER_ELEMENT = (DIM+1)*(DIM+2)/2; // assuming quadratic
215 
217  static const size_t NUM_NODES_PER_BOUNDARY_ELEMENT = DIM*(DIM+1)/2; // assuming quadratic
218 
224 
232  static double MAX_NEWTON_ABS_TOL;
233 
235  static double MIN_NEWTON_ABS_TOL;
236 
238  static double NEWTON_REL_TOL;
239 
245 
251 
258  c_matrix<double,DIM,DIM> mChangeOfBasisMatrix;
259 
266 
273 
278 
281 
287  double mCurrentTime;
288 
289 
294 
300 
307 
312 
318 
323 
332 
338 
346  std::vector<c_vector<double,DIM*(DIM+1)/2> > mAverageStressesPerElement;
347 
355  void AddStressToAverageStressPerElement(c_matrix<double,DIM,DIM>& rT, unsigned elementIndex);
356 
365  virtual void SetKspSolverAndPcType(KSP solver);
366 
367 
380  virtual void AssembleSystem(bool assembleResidual, bool assembleLinearSystem)=0;
381 
382 
383 
391  virtual void FinishAssembleSystem(bool assembleResidual, bool assembleLinearSystem);
392 
393 
395  //
396  // Element level methods
397  //
399 
407  void GetElementCentroidStrain(StrainType strainType,
408  Element<DIM,DIM>& rElement,
409  c_matrix<double,DIM,DIM>& rDeformationGradient);
410 
429  virtual void AddActiveStressAndStressDerivative(c_matrix<double,DIM,DIM>& rC,
430  unsigned elementIndex,
431  unsigned currentQuadPointGlobalIndex,
432  c_matrix<double,DIM,DIM>& rT,
434  bool addToDTdE)
435  {
436  // does nothing - subclass needs to overload this if there are active stresses
437  }
438 
446  virtual void SetupChangeOfBasisMatrix(unsigned elementIndex, unsigned currentQuadPointGlobalIndex)
447  {
448  }
449 
450 
451 
471  c_matrix<double, BOUNDARY_STENCIL_SIZE, BOUNDARY_STENCIL_SIZE>& rAelem,
472  c_vector<double, BOUNDARY_STENCIL_SIZE>& rBelem,
473  bool assembleResidual,
474  bool assembleJacobian,
475  unsigned boundaryConditionIndex);
476 
477 
493 
514  c_matrix<double,BOUNDARY_STENCIL_SIZE,BOUNDARY_STENCIL_SIZE>& rAelem,
515  c_vector<double,BOUNDARY_STENCIL_SIZE>& rBelem,
516  bool assembleResidual,
517  bool assembleJacobian,
518  unsigned boundaryConditionIndex);
519 
521  //
522  // These methods form the non-SNES nonlinear solver
523  //
525 
537  double ComputeResidualAndGetNorm(bool allowException);
538 
540  double CalculateResidualNorm();
541 
551  void VectorSum(std::vector<double>& rX, ReplicatableVector& rY, double a, std::vector<double>& rZ);
552 
560  void PrintLineSearchResult(double s, double residNorm);
561 
569  double TakeNewtonStep();
570 
580  double UpdateSolutionUsingLineSearch(Vec solution);
581 
589  virtual void PostNewtonStep(unsigned counter, double normResidual);
590 
596  void SolveNonSnes(double tol=-1.0);
597 
598 
600  //
601  // These methods form the SNES nonlinear solver
602  //
604 public: // need to be public as are called by global functions
611  void ComputeResidual(Vec currentGuess, Vec residualVector);
612 
620  void ComputeJacobian(Vec currentGuess, Mat* pJacobian, Mat* pPreconditioner);
621 
622 private:
627  void SolveSnes();
628 
629 public:
630 
642  SolidMechanicsProblemDefinition<DIM>& rProblemDefinition,
643  std::string outputDirectory,
644  CompressibilityType compressibilityType);
645 
650 
656  void Solve(double tol=-1.0);
657 
667  void SetIncludeActiveTension(bool includeActiveTension = true)
668  {
669  mIncludeActiveTension = includeActiveTension;
670  }
671 
675  unsigned GetNumNewtonIterations();
676 
677 
685  void SetWriteOutputEachNewtonIteration(bool writeOutputEachNewtonIteration=true)
686  {
687  mWriteOutputEachNewtonIteration = writeOutputEachNewtonIteration;
688  }
689 
696  void SetKspAbsoluteTolerance(double kspAbsoluteTolerance)
697  {
698  assert(kspAbsoluteTolerance > 0);
699  mKspAbsoluteTol = kspAbsoluteTolerance;
700  }
701 
715  void SetTakeFullFirstNewtonStep(bool takeFullFirstStep = true)
716  {
717  mTakeFullFirstNewtonStep = takeFullFirstStep;
718  }
719 
733  void SetUsePetscDirectSolve(bool usePetscDirectSolve = true)
734  {
735  mPetscDirectSolve = usePetscDirectSolve;
736  }
737 
738 
745  void SetCurrentTime(double time)
746  {
747  mCurrentTime = time;
748  }
749 
754  void CreateCmguiOutput();
755 
756 
770  void WriteCurrentStrains(StrainType strainType, std::string fileName, int counterToAppend = -1);
771 
778  void SetComputeAverageStressPerElementDuringSolve(bool setComputeAverageStressPerElement = true);
779 
793  void WriteCurrentAverageElementStresses(std::string fileName, int counterToAppend = -1);
794 
799  std::vector<c_vector<double,DIM> >& rGetSpatialSolution();
800 
805  std::vector<c_vector<double,DIM> >& rGetDeformedPosition();
806 
815  c_matrix<double,DIM,DIM> GetAverageStressPerElement(unsigned elementIndex);
816 };
817 
819 // Implementation: first, the non-nonlinear-solve methods
821 
822 template<unsigned DIM>
824  SolidMechanicsProblemDefinition<DIM>& rProblemDefinition,
825  std::string outputDirectory,
826  CompressibilityType compressibilityType)
827  : AbstractContinuumMechanicsSolver<DIM>(rQuadMesh, rProblemDefinition, outputDirectory, compressibilityType),
828  mrProblemDefinition(rProblemDefinition),
829  mrJacobianMatrix(this->mSystemLhsMatrix),
830  mKspAbsoluteTol(-1),
831  mWriteOutputEachNewtonIteration(false),
832  mNumNewtonIterations(0),
833  mCurrentTime(0.0),
834  mCheckedOutwardNormals(false),
835  mLastDampingValue(0.0),
836  mIncludeActiveTension(true),
837  mSetComputeAverageStressPerElement(false)
838 {
839  mUseSnesSolver = (mrProblemDefinition.GetSolveUsingSnes() ||
840  CommandLineArguments::Instance()->OptionExists("-mech_use_snes") );
841 
842  mChangeOfBasisMatrix = identity_matrix<double>(DIM,DIM);
843 
844  mTakeFullFirstNewtonStep = CommandLineArguments::Instance()->OptionExists("-mech_full_first_newton_step");
845  mPetscDirectSolve = CommandLineArguments::Instance()->OptionExists("-mech_petsc_direct_solve");
846 }
847 
848 template<unsigned DIM>
850 {
851 }
852 
853 
854 template<unsigned DIM>
855 void AbstractNonlinearElasticitySolver<DIM>::FinishAssembleSystem(bool assembleResidual, bool assembleJacobian)
856 {
857  PetscVecTools::Finalise(this->mResidualVector);
858 
859  if (assembleJacobian)
860  {
861  PetscMatTools::SwitchWriteMode(mrJacobianMatrix);
862  PetscMatTools::SwitchWriteMode(this->mPreconditionMatrix);
863 
864  VecCopy(this->mResidualVector, this->mLinearSystemRhsVector);
865  }
866 
867  // Apply Dirichlet boundary conditions
868  if(assembleJacobian)
869  {
870  this->ApplyDirichletBoundaryConditions(NONLINEAR_PROBLEM_APPLY_TO_EVERYTHING, this->mCompressibilityType==COMPRESSIBLE);
871  }
872  else if (assembleResidual)
873  {
874  this->ApplyDirichletBoundaryConditions(NONLINEAR_PROBLEM_APPLY_TO_RESIDUAL_ONLY, this->mCompressibilityType==COMPRESSIBLE);
875  }
876 
877  if (assembleResidual)
878  {
879  PetscVecTools::Finalise(this->mResidualVector);
880  }
881  if (assembleJacobian)
882  {
883  PetscMatTools::Finalise(mrJacobianMatrix);
884  PetscMatTools::Finalise(this->mPreconditionMatrix);
885  PetscVecTools::Finalise(this->mLinearSystemRhsVector);
886  }
887 }
888 
889 
890 template<unsigned DIM>
892 {
893  this->mSpatialSolution.clear();
894  this->mSpatialSolution.resize(this->mrQuadMesh.GetNumNodes(), zero_vector<double>(DIM));
895  for (unsigned i=0; i<this->mrQuadMesh.GetNumNodes(); i++)
896  {
897  for (unsigned j=0; j<DIM; j++)
898  {
899  this->mSpatialSolution[i](j) = this->mrQuadMesh.GetNode(i)->rGetLocation()[j] + this->mCurrentSolution[this->mProblemDimension*i+j];
900  }
901  }
902  return this->mSpatialSolution;
903 }
904 
905 template<unsigned DIM>
907 {
908  return rGetSpatialSolution();
909 }
910 
911 
912 template<unsigned DIM>
913 void AbstractNonlinearElasticitySolver<DIM>::WriteCurrentStrains(StrainType strainType, std::string fileName, int counterToAppend)
914 {
915  if (!this->mWriteOutput)
916  {
917  return;
918  }
919 
920  std::stringstream file_name;
921  file_name << fileName;
922  if (counterToAppend >= 0)
923  {
924  file_name << "_" << counterToAppend;
925  }
926  file_name << ".strain";
927 
928  out_stream p_file = this->mpOutputFileHandler->OpenOutputFile(file_name.str());
929 
930  c_matrix<double,DIM,DIM> strain;
931 
932  for (typename AbstractTetrahedralMesh<DIM,DIM>::ElementIterator iter = this->mrQuadMesh.GetElementIteratorBegin();
933  iter != this->mrQuadMesh.GetElementIteratorEnd();
934  ++iter)
935  {
936  GetElementCentroidStrain(strainType, *iter, strain);
937  for(unsigned i=0; i<DIM; i++)
938  {
939  for(unsigned j=0; j<DIM; j++)
940  {
941  *p_file << strain(i,j) << " ";
942  }
943  }
944  *p_file << "\n";
945  }
946  p_file->close();
947 }
948 
949 
950 
951 template<unsigned DIM>
953 {
954  if (!this->mWriteOutput)
955  {
956  return;
957  }
958 
959  if(!mSetComputeAverageStressPerElement)
960  {
961  EXCEPTION("Call SetComputeAverageStressPerElementDuringSolve() before solve if calling WriteCurrentAverageElementStresses()");
962  }
963 
964  std::stringstream file_name;
965  file_name << fileName;
966  if (counterToAppend >= 0)
967  {
968  file_name << "_" << counterToAppend;
969  }
970  file_name << ".stress";
971  assert(mAverageStressesPerElement.size()==this->mrQuadMesh.GetNumElements());
972 
973  StressPerElementWriter<DIM> stress_writer(&(this->mrQuadMesh),this);
974  stress_writer.WriteData(*(this->mpOutputFileHandler), file_name.str());
975 }
976 
977 
978 template<unsigned DIM>
980 {
981  if (this->mOutputDirectory=="")
982  {
983  EXCEPTION("No output directory was given so no output was written, cannot convert to cmgui format");
984  }
985 
986  CmguiDeformedSolutionsWriter<DIM> writer(this->mOutputDirectory + "/cmgui",
987  "solution",
988  this->mrQuadMesh,
989  WRITE_QUADRATIC_MESH);
990 
991  std::vector<c_vector<double,DIM> >& r_deformed_positions = this->rGetDeformedPosition();
992  writer.WriteInitialMesh(); // this writes solution_0.exnode and .exelem
993  writer.WriteDeformationPositions(r_deformed_positions, 1); // this writes the final solution as solution_1.exnode
994  writer.WriteCmguiScript(); // writes LoadSolutions.com
995 }
996 
997 template<unsigned DIM>
999 {
1000  mSetComputeAverageStressPerElement = setComputeAverageStressPerElement;
1001  if(setComputeAverageStressPerElement && mAverageStressesPerElement.size()==0)
1002  {
1003  mAverageStressesPerElement.resize(this->mrQuadMesh.GetNumElements(), zero_vector<double>(DIM*(DIM+1)/2));
1004  }
1005 }
1006 
1007 template<unsigned DIM>
1008 void AbstractNonlinearElasticitySolver<DIM>::AddStressToAverageStressPerElement(c_matrix<double,DIM,DIM>& rT, unsigned elemIndex)
1009 {
1010  assert(mSetComputeAverageStressPerElement);
1011  assert(elemIndex<this->mrQuadMesh.GetNumElements());
1012 
1013  // In 2d the matrix is
1014  // [T00 T01]
1015  // [T10 T11]
1016  // where T01 = T10. We store this as a vector
1017  // [T00 T01 T11]
1018  //
1019  // Similarly, for 3d we store
1020  // [T00 T01 T02 T11 T12 T22]
1021  for(unsigned i=0; i<DIM*(DIM+1)/2; i++)
1022  {
1023  unsigned row;
1024  unsigned col;
1025  if(DIM==2)
1026  {
1027  row = i<=1 ? 0 : 1;
1028  col = i==0 ? 0 : 1;
1029  }
1030  else // DIM==3
1031  {
1032  row = i<=2 ? 0 : (i<=4? 1 : 2);
1033  col = i==0 ? 0 : (i==1 || i==3? 1 : 2);
1034  }
1035 
1036  this->mAverageStressesPerElement[elemIndex](i) += rT(row,col);
1037  }
1038 }
1039 
1040 
1041 template<unsigned DIM>
1042 c_matrix<double,DIM,DIM> AbstractNonlinearElasticitySolver<DIM>::GetAverageStressPerElement(unsigned elementIndex)
1043 {
1044  if(!mSetComputeAverageStressPerElement)
1045  {
1046  EXCEPTION("Call SetComputeAverageStressPerElementDuringSolve() before solve if calling GetAverageStressesPerElement()");
1047  }
1048  assert(elementIndex<this->mrQuadMesh.GetNumElements());
1049 
1050  c_matrix<double,DIM,DIM> stress;
1051 
1052  // In 2d the matrix is
1053  // [T00 T01]
1054  // [T10 T11]
1055  // where T01 = T10, and was stored as
1056  // [T00 T01 T11]
1057  //
1058  // Similarly, for 3d the matrix was stored as
1059  // [T00 T01 T02 T11 T12 T22]
1060  if(DIM==2)
1061  {
1062  stress(0,0) = mAverageStressesPerElement[elementIndex](0);
1063  stress(1,0) = stress(0,1) = mAverageStressesPerElement[elementIndex](1);
1064  stress(1,1) = mAverageStressesPerElement[elementIndex](2);
1065  }
1066  else
1067  {
1068  stress(0,0) = mAverageStressesPerElement[elementIndex](0);
1069  stress(1,0) = stress(0,1) = mAverageStressesPerElement[elementIndex](1);
1070  stress(2,0) = stress(0,2) = mAverageStressesPerElement[elementIndex](2);
1071  stress(1,1) = mAverageStressesPerElement[elementIndex](3);
1072  stress(2,1) = stress(1,2) = mAverageStressesPerElement[elementIndex](4);
1073  stress(2,2) = mAverageStressesPerElement[elementIndex](5);
1074  }
1075 
1076  return stress;
1077 }
1078 
1079 
1081 // Methods at the 'element level'.
1083 
1084 template<unsigned DIM>
1086  Element<DIM,DIM>& rElement,
1087  c_matrix<double,DIM,DIM>& rStrain)
1088 {
1089  static c_matrix<double,DIM,DIM> jacobian;
1090  static c_matrix<double,DIM,DIM> inverse_jacobian;
1091  double jacobian_determinant;
1092 
1093  this->mrQuadMesh.GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
1094 
1095  // Get the current displacement at the nodes
1096  static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> element_current_displacements;
1097  for (unsigned II=0; II<NUM_NODES_PER_ELEMENT; II++)
1098  {
1099  for (unsigned JJ=0; JJ<DIM; JJ++)
1100  {
1101  element_current_displacements(JJ,II) = this->mCurrentSolution[this->mProblemDimension*rElement.GetNodeGlobalIndex(II) + JJ];
1102  }
1103  }
1104 
1105  // Allocate memory for the basis functions values and derivative values
1106  static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> grad_quad_phi;
1107  static c_matrix<double,DIM,DIM> grad_u; // grad_u = (du_i/dX_M)
1108 
1109  // we need the point in the canonical element which corresponds to the centroid of the
1110  // version of the element in physical space. This point can be shown to be (1/3,1/3).
1111  ChastePoint<DIM> quadrature_point;
1112  if(DIM==2)
1113  {
1114  quadrature_point.rGetLocation()(0) = 1.0/3.0;
1115  quadrature_point.rGetLocation()(1) = 1.0/3.0;
1116  }
1117  else
1118  {
1119  assert(DIM==3);
1120  quadrature_point.rGetLocation()(0) = 1.0/4.0;
1121  quadrature_point.rGetLocation()(1) = 1.0/4.0;
1122  quadrature_point.rGetLocation()(2) = 1.0/4.0;
1123  }
1124 
1125  QuadraticBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(quadrature_point, inverse_jacobian, grad_quad_phi);
1126 
1127  // Interpolate grad_u
1128  grad_u = zero_matrix<double>(DIM,DIM);
1129  for (unsigned node_index=0; node_index<NUM_NODES_PER_ELEMENT; node_index++)
1130  {
1131  for (unsigned i=0; i<DIM; i++)
1132  {
1133  for (unsigned M=0; M<DIM; M++)
1134  {
1135  grad_u(i,M) += grad_quad_phi(M,node_index)*element_current_displacements(i,node_index);
1136  }
1137  }
1138  }
1139 
1140  c_matrix<double,DIM,DIM> deformation_gradient;
1141 
1142  for (unsigned i=0; i<DIM; i++)
1143  {
1144  for (unsigned M=0; M<DIM; M++)
1145  {
1146  deformation_gradient(i,M) = (i==M?1:0) + grad_u(i,M);
1147  }
1148  }
1149 
1150  switch(strainType)
1151  {
1152  case DEFORMATION_GRADIENT_F:
1153  {
1154  rStrain = deformation_gradient;
1155  break;
1156  }
1157  case DEFORMATION_TENSOR_C:
1158  {
1159  rStrain = prod(trans(deformation_gradient),deformation_gradient);
1160  break;
1161  }
1162  case LAGRANGE_STRAIN_E:
1163  {
1164  c_matrix<double,DIM,DIM> C = prod(trans(deformation_gradient),deformation_gradient);
1165  for (unsigned M=0; M<DIM; M++)
1166  {
1167  for (unsigned N=0; N<DIM; N++)
1168  {
1169  rStrain(M,N) = 0.5* ( C(M,N)-(M==N?1:0) );
1170  }
1171  }
1172  break;
1173  }
1174  default:
1175  {
1176  NEVER_REACHED;
1177  break;
1178  }
1179  }
1180 }
1181 
1182 
1183 
1184 template<unsigned DIM>
1186  BoundaryElement<DIM-1,DIM>& rBoundaryElement,
1187  c_matrix<double,BOUNDARY_STENCIL_SIZE,BOUNDARY_STENCIL_SIZE>& rAelem,
1188  c_vector<double,BOUNDARY_STENCIL_SIZE>& rBelem,
1189  bool assembleResidual,
1190  bool assembleJacobian,
1191  unsigned boundaryConditionIndex)
1192 {
1193  if( this->mrProblemDefinition.GetTractionBoundaryConditionType() == PRESSURE_ON_DEFORMED
1194  || this->mrProblemDefinition.GetTractionBoundaryConditionType() == FUNCTIONAL_PRESSURE_ON_DEFORMED)
1195  {
1196  AssembleOnBoundaryElementForPressureOnDeformedBc(rBoundaryElement, rAelem, rBelem,
1197  assembleResidual, assembleJacobian, boundaryConditionIndex);
1198  return;
1199  }
1200 
1201  rAelem.clear();
1202  rBelem.clear();
1203 
1204  if (assembleJacobian && !assembleResidual)
1205  {
1206  // Nothing to do
1207  return;
1208  }
1209 
1210  c_vector<double, DIM> weighted_direction;
1211  double jacobian_determinant;
1212  this->mrQuadMesh.GetWeightedDirectionForBoundaryElement(rBoundaryElement.GetIndex(), weighted_direction, jacobian_determinant);
1213 
1214  c_vector<double,NUM_NODES_PER_BOUNDARY_ELEMENT> phi;
1215 
1216  for (unsigned quad_index=0; quad_index<this->mpBoundaryQuadratureRule->GetNumQuadPoints(); quad_index++)
1217  {
1218  double wJ = jacobian_determinant * this->mpBoundaryQuadratureRule->GetWeight(quad_index);
1219 
1220  const ChastePoint<DIM-1>& quad_point = this->mpBoundaryQuadratureRule->rGetQuadPoint(quad_index);
1221 
1223 
1224  // Get the required traction, interpolating X (slightly inefficiently,
1225  // as interpolating using quad bases) if necessary
1226  c_vector<double,DIM> traction = zero_vector<double>(DIM);
1227  switch (this->mrProblemDefinition.GetTractionBoundaryConditionType())
1228  {
1229  case FUNCTIONAL_TRACTION:
1230  {
1231  c_vector<double,DIM> X = zero_vector<double>(DIM);
1232  for (unsigned node_index=0; node_index<NUM_NODES_PER_BOUNDARY_ELEMENT; node_index++)
1233  {
1234  X += phi(node_index)*this->mrQuadMesh.GetNode( rBoundaryElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
1235  }
1236  traction = this->mrProblemDefinition.EvaluateTractionFunction(X, this->mCurrentTime);
1237  break;
1238  }
1239  case ELEMENTWISE_TRACTION:
1240  {
1241  traction = this->mrProblemDefinition.rGetElementwiseTractions()[boundaryConditionIndex];
1242  break;
1243  }
1244  default:
1245  NEVER_REACHED;
1246  }
1247 
1248 
1249  for (unsigned index=0; index<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index++)
1250  {
1251  unsigned spatial_dim = index%DIM;
1252  unsigned node_index = (index-spatial_dim)/DIM;
1253 
1254  assert(node_index < NUM_NODES_PER_BOUNDARY_ELEMENT);
1255 
1256  rBelem(index) -= traction(spatial_dim)
1257  * phi(node_index)
1258  * wJ;
1259  }
1260  }
1261 }
1262 
1263 template<unsigned DIM>
1265 {
1266  if(mUseSnesSolver)
1267  {
1268  // although not using this in the first few steps might be useful when the deformation
1269  // is large, the snes solver is more robust, so we have this on all the time. (Also because
1270  // for cardiac problems and in timesteps after the initial large deformation, we want this on
1271  // in the first step
1272  return true;
1273 
1274  // could do something like this, if make the snes a member variable
1275  //PetscInt iteration_number;
1276  //SNESGetIterationNumber(mSnes,&iteration_number);
1277  //return (iteration_number >= 3);
1278  }
1279  else
1280  {
1281  return (mLastDampingValue >= 0.5);
1282  }
1283 
1284 }
1285 
1286 template<unsigned DIM>
1288  BoundaryElement<DIM-1,DIM>& rBoundaryElement,
1289  c_matrix<double,BOUNDARY_STENCIL_SIZE,BOUNDARY_STENCIL_SIZE>& rAelem,
1290  c_vector<double,BOUNDARY_STENCIL_SIZE>& rBelem,
1291  bool assembleResidual,
1292  bool assembleJacobian,
1293  unsigned boundaryConditionIndex)
1294 {
1295  assert( this->mrProblemDefinition.GetTractionBoundaryConditionType()==PRESSURE_ON_DEFORMED
1296  || this->mrProblemDefinition.GetTractionBoundaryConditionType()==FUNCTIONAL_PRESSURE_ON_DEFORMED);
1297 
1298  rAelem.clear();
1299  rBelem.clear();
1300 
1301  c_vector<double, DIM> weighted_direction;
1302  double jacobian_determinant;
1303  // note: jacobian determinant may be over-written below
1304  this->mrQuadMesh.GetWeightedDirectionForBoundaryElement(rBoundaryElement.GetIndex(), weighted_direction, jacobian_determinant);
1305 
1306 
1308  // Find the volume element of the mesh which
1309  // contains this boundary element
1311 
1312  Element<DIM,DIM>* p_containing_vol_element = NULL;
1313 
1314  std::set<unsigned> potential_elements = rBoundaryElement.GetNode(0)->rGetContainingElementIndices();
1315  for(std::set<unsigned>::iterator iter = potential_elements.begin();
1316  iter != potential_elements.end();
1317  iter++)
1318  {
1319  p_containing_vol_element = this->mrQuadMesh.GetElement(*iter);
1320 
1321  bool this_vol_ele_contains_surf_ele = true;
1322  // loop over the nodes of boundary element and see if they are in the volume element
1323  for(unsigned i=1; i<NUM_NODES_PER_BOUNDARY_ELEMENT; i++) // don't need to start at 0, given looping over contain elems of node 0
1324  {
1325  unsigned surf_element_node_index = rBoundaryElement.GetNodeGlobalIndex(i);
1326  bool found_this_node = false;
1327  for(unsigned j=0; j<p_containing_vol_element->GetNumNodes(); j++)
1328  {
1329  unsigned vol_element_node_index = p_containing_vol_element->GetNodeGlobalIndex(j);
1330  if(surf_element_node_index == vol_element_node_index)
1331  {
1332  found_this_node = true;
1333  break;
1334  }
1335  }
1336  if(!found_this_node)
1337  {
1338  this_vol_ele_contains_surf_ele = false;
1339  break;
1340  }
1341  }
1342  if(this_vol_ele_contains_surf_ele)
1343  {
1344  break;
1345  }
1346  }
1347 
1348  // Find the local node index in the volume element for each node in the boundary element
1349  std::vector<unsigned> surf_to_vol_map(NUM_NODES_PER_BOUNDARY_ELEMENT);
1350  for(unsigned i=0; i<NUM_NODES_PER_BOUNDARY_ELEMENT; i++)
1351  {
1352  unsigned index = rBoundaryElement.GetNodeGlobalIndex(i);
1353  for(unsigned j=0; j<NUM_NODES_PER_ELEMENT; j++)
1354  {
1355  if(p_containing_vol_element->GetNodeGlobalIndex(j)==index)
1356  {
1357  surf_to_vol_map[i] = j;
1358  break;
1359  }
1360  }
1361  }
1362 
1363 
1364  // We require the volume element to compute F, which requires grad_phi on the volume element. For this we will
1365  // need the inverse jacobian for the volume element
1366  static c_matrix<double,DIM,DIM> jacobian_vol_element;
1367  static c_matrix<double,DIM,DIM> inverse_jacobian_vol_element;
1368  double jacobian_determinant_vol_element;
1369  this->mrQuadMesh.GetInverseJacobianForElement(p_containing_vol_element->GetIndex(), jacobian_vol_element, jacobian_determinant_vol_element, inverse_jacobian_vol_element);
1370 
1371  // Get the current displacements at each node of the volume element, to be used in computing F
1372  static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> element_current_displacements;
1373  for (unsigned II=0; II<NUM_NODES_PER_ELEMENT; II++)
1374  {
1375  for (unsigned JJ=0; JJ<DIM; JJ++)
1376  {
1377  element_current_displacements(JJ,II) = this->mCurrentSolution[this->mProblemDimension*p_containing_vol_element->GetNodeGlobalIndex(II) + JJ];
1378  }
1379  }
1380 
1381 
1382  // We will need both {grad phi_i} for the quadratic bases of the volume element, for computing F..
1383  static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> grad_quad_phi_vol_element;
1384  // ..the phi_i for each of the quadratic bases of the surface element, for the standard FE assembly part.
1385  c_vector<double,NUM_NODES_PER_BOUNDARY_ELEMENT> quad_phi_surf_element;
1386  // We need this too, which is obtained by taking a subset of grad_quad_phi_vol_element
1387  static c_matrix<double, DIM, NUM_NODES_PER_BOUNDARY_ELEMENT> grad_quad_phi_surf_element;
1388 
1389  c_matrix<double,DIM,DIM> F;
1390  c_matrix<double,DIM,DIM> invF;
1391 
1392  c_vector<double,DIM> normal = rBoundaryElement.CalculateNormal();
1393  c_matrix<double,1,DIM> normal_as_mat;
1394  for(unsigned i=0; i<DIM; i++)
1395  {
1396  normal_as_mat(0,i) = normal(i);
1397  }
1398 
1399  double normal_pressure;
1400  switch (this->mrProblemDefinition.GetTractionBoundaryConditionType())
1401  {
1402  case PRESSURE_ON_DEFORMED:
1403  normal_pressure = this->mrProblemDefinition.GetNormalPressure();
1404  break;
1405  case FUNCTIONAL_PRESSURE_ON_DEFORMED:
1406  normal_pressure = this->mrProblemDefinition.EvaluateNormalPressureFunction(this->mCurrentTime);
1407  break;
1408  default:
1409  NEVER_REACHED;
1410  }
1411 
1412 
1413 
1414  for (unsigned quad_index=0; quad_index<this->mpBoundaryQuadratureRule->GetNumQuadPoints(); quad_index++)
1415  {
1416  double wJ = jacobian_determinant * this->mpBoundaryQuadratureRule->GetWeight(quad_index);
1417 
1418  // Get the quadrature point on this surface element (in canonical space) - so eg, for a 2D problem,
1419  // the quad point is in 1D space
1420  const ChastePoint<DIM-1>& quadrature_point = this->mpBoundaryQuadratureRule->rGetQuadPoint(quad_index);
1421  QuadraticBasisFunction<DIM-1>::ComputeBasisFunctions(quadrature_point, quad_phi_surf_element);
1422 
1423  // We will need the xi coordinates of this quad point in the volume element. We could do this by figuring
1424  // out how the nodes of the surface element are ordered in the list of nodes in the volume element,
1425  // however it is less fiddly to compute directly. Firstly, compute the corresponding physical location
1426  // of the quad point, by interpolating
1427  c_vector<double,DIM> X = zero_vector<double>(DIM);
1428  for (unsigned node_index=0; node_index<NUM_NODES_PER_BOUNDARY_ELEMENT; node_index++)
1429  {
1430  X += quad_phi_surf_element(node_index)*rBoundaryElement.GetNode(node_index)->rGetLocation();
1431  }
1432 
1433 
1434  // Now compute the xi coordinates of the quad point in the volume element
1435  c_vector<double,DIM+1> weight = p_containing_vol_element->CalculateInterpolationWeights(X);
1436  c_vector<double,DIM> xi;
1437  for(unsigned i=0; i<DIM; i++)
1438  {
1439  xi(i) = weight(i+1); // Note, in 2d say, weights = [1-xi(0)-xi(1), xi(0), xi(1)]
1440  }
1441 
1442  // check one of the weights was zero, as the quad point is on the boundary of the volume element
1443  if(DIM==2)
1444  {
1445  assert( DIM!=2 || (fabs(weight(0))<1e-6) || (fabs(weight(1))<1e-6) || (fabs(weight(2))<1e-6) );
1446  }
1447  else
1448  {
1449  #define COVERAGE_IGNORE
1450  assert( DIM!=3 || (fabs(weight(0))<1e-6) || (fabs(weight(1))<1e-6) || (fabs(weight(2))<1e-6) || (fabs(weight(3))<1e-6) );
1451  #undef COVERAGE_IGNORE
1452  }
1453 
1454  // Now we can compute the grad_phi and then interpolate F
1455  QuadraticBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(xi, inverse_jacobian_vol_element, grad_quad_phi_vol_element);
1456 
1457  F = identity_matrix<double>(DIM,DIM);
1458  for (unsigned node_index=0; node_index<NUM_NODES_PER_ELEMENT; node_index++)
1459  {
1460  for (unsigned i=0; i<DIM; i++)
1461  {
1462  for (unsigned M=0; M<DIM; M++)
1463  {
1464  F(i,M) += grad_quad_phi_vol_element(M,node_index)*element_current_displacements(i,node_index);
1465  }
1466  }
1467  }
1468 
1469  double detF = Determinant(F);
1470  invF = Inverse(F);
1471 
1472  if(assembleResidual)
1473  {
1474  c_vector<double,DIM> traction = detF*normal_pressure*prod(trans(invF),normal);
1475 
1476  // assemble
1477  for (unsigned index=0; index<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index++)
1478  {
1479  unsigned spatial_dim = index%DIM;
1480  unsigned node_index = (index-spatial_dim)/DIM;
1481 
1482  assert(node_index < NUM_NODES_PER_BOUNDARY_ELEMENT);
1483 
1484  rBelem(index) -= traction(spatial_dim)
1485  * quad_phi_surf_element(node_index)
1486  * wJ;
1487  }
1488  }
1489 
1490 
1491  // Sometimes we don't include the analytic jacobian for this integral
1492  // in the jacobian that is assembled - the ShouldAssembleMatrixTermForPressureOnDeformedBc()
1493  // bit below - see the documentation for this method to see why.
1494  if(assembleJacobian && ShouldAssembleMatrixTermForPressureOnDeformedBc())
1495  {
1496  for(unsigned II=0; II<NUM_NODES_PER_BOUNDARY_ELEMENT; II++)
1497  {
1498  for(unsigned N=0; N<DIM; N++)
1499  {
1500  grad_quad_phi_surf_element(N,II) = grad_quad_phi_vol_element(N,surf_to_vol_map[II]);
1501  }
1502  }
1503 
1504  static FourthOrderTensor<DIM,DIM,DIM,DIM> tensor1;
1505  for(unsigned N=0; N<DIM; N++)
1506  {
1507  for(unsigned e=0; e<DIM; e++)
1508  {
1509  for(unsigned M=0; M<DIM; M++)
1510  {
1511  for(unsigned d=0; d<DIM; d++)
1512  {
1513  tensor1(N,e,M,d) = invF(N,e)*invF(M,d) - invF(M,e)*invF(N,d);
1514  }
1515  }
1516  }
1517  }
1518 
1519  // tensor2(II,e,M,d) = tensor1(N,e,M,d)*grad_quad_phi_surf_element(N,II)
1521  tensor2.template SetAsContractionOnFirstDimension<DIM>( trans(grad_quad_phi_surf_element), tensor1);
1522 
1523  // tensor3 is really a third-order tensor
1524  // tensor3(II,e,0,d) = tensor2(II,e,M,d)*normal(M)
1526  tensor3.template SetAsContractionOnThirdDimension<DIM>( normal_as_mat, tensor2);
1527 
1528  for (unsigned index1=0; index1<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index1++)
1529  {
1530  unsigned spatial_dim1 = index1%DIM;
1531  unsigned node_index1 = (index1-spatial_dim1)/DIM;
1532 
1533  for (unsigned index2=0; index2<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index2++)
1534  {
1535  unsigned spatial_dim2 = index2%DIM;
1536  unsigned node_index2 = (index2-spatial_dim2)/DIM;
1537 
1538  rAelem(index1,index2) -= normal_pressure
1539  * detF
1540  * tensor3(node_index2,spatial_dim2,0,spatial_dim1)
1541  * quad_phi_surf_element(node_index1)
1542  * wJ;
1543  }
1544  }
1545  }
1546  }
1547 }
1548 
1549 
1550 
1551 
1552 template<unsigned DIM>
1554 {
1555 
1556  // check the problem definition is set up correctly (and fully).
1557  mrProblemDefinition.Validate();
1558 
1559  // If the problem includes specified pressures on deformed surfaces (as opposed
1560  // to specified tractions), the code needs to compute normals, and they need
1561  // to be consistently all facing outward (or all facing inward). Check the undeformed
1562  // mesh boundary elements has nodes that are ordered so that all normals are
1563  // outward-facing
1564  if(mrProblemDefinition.GetTractionBoundaryConditionType()==PRESSURE_ON_DEFORMED && mCheckedOutwardNormals==false)
1565  {
1566  this->mrQuadMesh.CheckOutwardNormals();
1567  mCheckedOutwardNormals = true;
1568  }
1569 
1570  // Write the initial solution
1571  this->WriteCurrentSpatialSolution("initial", "nodes");
1572 
1573  if(mUseSnesSolver)
1574  {
1575  SolveSnes();
1576  }
1577  else
1578  {
1579  SolveNonSnes(tol);
1580  }
1581 
1582  // Remove pressure dummy values (P=0 at internal nodes, which should have been
1583  // been the result of the solve above), by linear interpolating from vertices of
1584  // edges to the internal node of the edge
1585  if(this->mCompressibilityType==INCOMPRESSIBLE)
1586  {
1587  this->RemovePressureDummyValuesThroughLinearInterpolation();
1588  }
1589 
1590  // Write the final solution
1591  this->WriteCurrentSpatialSolution("solution", "nodes");
1592 
1593 }
1594 
1595 
1599 template<unsigned DIM>
1601 {
1602  // Four alternatives
1603  // (a) Petsc direct solve
1604  // Otherwise iterative solve with:
1605  // (b) Incompressible: GMRES with ILU preconditioner (or bjacobi=ILU on each process) [default]. Very poor on large problems.
1606  // (c) Incompressible: GMRES with AMG preconditioner. Uncomment #define MECH_USE_HYPRE above. Requires Petsc3 with HYPRE installed.
1607  // (d) Compressible: CG with ICC
1608 
1609  PC pc;
1610  KSPGetPC(solver, &pc);
1611 
1612 
1613  if (mPetscDirectSolve)
1614  {
1615  if(this->mCompressibilityType==COMPRESSIBLE)
1616  {
1617  KSPSetType(solver,KSPPREONLY);
1618 
1619  }
1620  PCSetType(pc, PCLU);
1621 
1622  // See #2057
1623  // PCFactorSetMatSolverPackage(pc,"mumps");
1624  }
1625  else
1626  {
1627  if (this->mCompressibilityType==COMPRESSIBLE)
1628  {
1629  KSPSetType(solver,KSPCG);
1631  {
1632  PCSetType(pc, PCICC);
1633  //Note that PetscOptionsSetValue is dangerous because we can't easily do
1634  //regression testing. If a name changes, then the behaviour of the code changes
1635  //because it won't recognise the old name. However, it won't fail to compile/run.
1636  #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 1) //PETSc 3.1 or later
1637  PetscOptionsSetValue("-pc_factor_shift_type", "positive_definite");
1638  #else
1639  PetscOptionsSetValue("-pc_factor_shift_positive_definite", "");
1640  #endif
1641  }
1642  else
1643  {
1644  PCSetType(pc, PCBJACOBI);
1645  }
1646  }
1647  else
1648  {
1649  unsigned num_restarts = 100;
1650  KSPSetType(solver,KSPGMRES);
1651  KSPGMRESSetRestart(solver,num_restarts);
1652 
1653  #ifndef MECH_USE_HYPRE
1654  PCSetType(pc, PCBJACOBI); // BJACOBI = ILU on each block (block = part of matrix on each process)
1655  #else
1656  // Speed up linear solve time massively for larger simulations (in fact GMRES may stagnate without
1658  // this for larger problems), by using a AMG preconditioner -- needs HYPRE installed
1660  PetscOptionsSetValue("-pc_hypre_type", "boomeramg");
1661  // PetscOptionsSetValue("-pc_hypre_boomeramg_max_iter", "1");
1662  // PetscOptionsSetValue("-pc_hypre_boomeramg_strong_threshold", "0.0");
1663 
1664  PCSetType(pc, PCHYPRE);
1665  #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >=2) //PETSc 3.2 or later
1666  KSPSetPCSide(solver, PC_RIGHT);
1667  #else
1668  KSPSetPreconditionerSide(solver, PC_RIGHT);
1669  #endif
1670 
1671  // other possible preconditioners..
1672  //PCBlockDiagonalMechanics* p_custom_pc = new PCBlockDiagonalMechanics(solver, this->mPreconditionMatrix, mBlock1Size, mBlock2Size);
1673  //PCLDUFactorisationMechanics* p_custom_pc = new PCLDUFactorisationMechanics(solver, this->mPreconditionMatrix, mBlock1Size, mBlock2Size);
1674  //remember to delete memory..
1675  //KSPSetPreconditionerSide(solver, PC_RIGHT);
1676  #endif
1677  }
1678  }
1679 }
1680 
1681 
1682 
1684 // The code for the non-SNES solver - maybe remove all this
1685 // as SNES solver appears better
1687 
1688 template<unsigned DIM>
1690 {
1691  if (!allowException)
1692  {
1693  // Assemble the residual
1694  AssembleSystem(true, false);
1695  }
1696  else
1697  {
1698  try
1699  {
1700  // Try to assemble the residual using this current solution
1701  AssembleSystem(true, false);
1702  }
1703  catch(Exception&)
1704  {
1705  // If fail (because e.g. ODEs fail to solve, or strains are too large for material law), return infinity
1706  return DBL_MAX;
1707  }
1708  }
1709 
1710  // Return the scaled norm of the residual
1711  return CalculateResidualNorm();
1712 }
1713 
1714 template<unsigned DIM>
1716 {
1717  double norm;
1718 
1719  //\todo Change to NORM_1 and remove the division by mNumDofs...
1720  VecNorm(this->mResidualVector, NORM_2, &norm);
1721  return norm/this->mNumDofs;
1722 }
1723 
1724 template<unsigned DIM>
1726  ReplicatableVector& rY,
1727  double a,
1728  std::vector<double>& rZ)
1729 {
1730  assert(rX.size()==rY.GetSize());
1731  assert(rY.GetSize()==rZ.size());
1732  for (unsigned i=0; i<rX.size(); i++)
1733  {
1734  rZ[i] = rX[i] + a*rY[i];
1735  }
1736 }
1737 
1738 
1739 template<unsigned DIM>
1741 {
1742  if(this->mVerbose)
1743  {
1744  Timer::Reset();
1745  }
1746 
1748  // Assemble Jacobian (and preconditioner)
1750  MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ASSEMBLE);
1751  AssembleSystem(true, true);
1752  MechanicsEventHandler::EndEvent(MechanicsEventHandler::ASSEMBLE);
1753  if(this->mVerbose)
1754  {
1755  Timer::PrintAndReset("AssembleSystem");
1756  }
1757 
1759  // Solve the linear system.
1761  MechanicsEventHandler::BeginEvent(MechanicsEventHandler::SOLVE);
1762 
1763  Vec solution;
1764  VecDuplicate(this->mResidualVector,&solution);
1765 
1766  KSP solver;
1767  KSPCreate(PETSC_COMM_WORLD,&solver);
1768 
1769 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
1770  KSPSetOperators(solver, mrJacobianMatrix, this->mPreconditionMatrix);
1771 #else
1772  KSPSetOperators(solver, mrJacobianMatrix, this->mPreconditionMatrix, DIFFERENT_NONZERO_PATTERN /*in precond between successive solves*/);
1773 #endif
1774 
1775  // Set the type of KSP solver (CG, GMRES etc) and preconditioner (ILU, HYPRE, etc)
1776  SetKspSolverAndPcType(solver);
1777 
1778  //PetscOptionsSetValue("-ksp_monitor","");
1779  //PetscOptionsSetValue("-ksp_norm_type","natural");
1780 
1781  KSPSetFromOptions(solver);
1782  KSPSetUp(solver);
1783 
1784 
1785  // Set the linear system absolute tolerance.
1786  // This is either the user provided value, or set to
1787  // max {1e-6 * initial_residual, 1e-12}
1788  if (mKspAbsoluteTol < 0)
1789  {
1790  Vec temp;
1791  VecDuplicate(this->mResidualVector, &temp);
1792  Vec temp2;
1793  VecDuplicate(this->mResidualVector, &temp2);
1794  Vec linsys_residual;
1795  VecDuplicate(this->mResidualVector, &linsys_residual);
1796 
1797  KSPInitialResidual(solver, solution, temp, temp2, linsys_residual, this->mLinearSystemRhsVector);
1798  double initial_resid_norm;
1799  VecNorm(linsys_residual, NORM_2, &initial_resid_norm);
1800 
1801  PetscTools::Destroy(temp);
1802  PetscTools::Destroy(temp2);
1803  PetscTools::Destroy(linsys_residual);
1804 
1805  double ksp_rel_tol = 1e-6;
1806  double absolute_tol = ksp_rel_tol * initial_resid_norm;
1807  if(absolute_tol < 1e-12)
1808  {
1809  absolute_tol = 1e-12;
1810  }
1811  KSPSetTolerances(solver, 1e-16, absolute_tol, PETSC_DEFAULT, 1000 /* max iters */); // Note: some machines - max iters seems to be 1000 whatever we give here
1812  }
1813  else
1814  {
1815  KSPSetTolerances(solver, 1e-16, mKspAbsoluteTol, PETSC_DEFAULT, 1000 /* max iters */); // Note: some machines - max iters seems to be 1000 whatever we give here
1816  }
1817 
1818  if(this->mVerbose)
1819  {
1820  Timer::PrintAndReset("KSP Setup");
1821  }
1822 
1823  KSPSolve(solver,this->mLinearSystemRhsVector,solution);
1824 
1825 // ///// For printing matrix when debugging
1826 // OutputFileHandler handler("TEMP",false);
1827 // std::stringstream ss;
1828 // static unsigned counter = 0;
1829 // ss << "all_" << counter++ << ".txt";
1830 // out_stream p_file = handler.OpenOutputFile(ss.str());
1831 // *p_file << std::setprecision(10);
1832 // for(unsigned i=0; i<this->mNumDofs; i++)
1833 // {
1834 // for(unsigned j=0; j<this->mNumDofs; j++)
1835 // {
1836 // *p_file << PetscMatTools::GetElement(mrJacobianMatrix, i, j) << " ";
1837 // }
1838 // *p_file << PetscVecTools::GetElement(this->mLinearSystemRhsVector, i) << " ";
1839 // *p_file << PetscVecTools::GetElement(solution, i) << "\n";
1840 // }
1841 // p_file->close();
1842 
1843 
1845  // Error checking for linear solve
1847 
1848  // warn if ksp reports failure
1849  KSPConvergedReason reason;
1850  KSPGetConvergedReason(solver,&reason);
1851 
1852  if(reason != KSP_DIVERGED_ITS)
1853  {
1854  // Throw an exception if the solver failed for any reason other than DIVERGED_ITS.
1855  // This is not covered as would be difficult to cover - requires a bad matrix to
1856  // assembled, for example.
1857  #define COVERAGE_IGNORE
1858  KSPEXCEPT(reason);
1859  #undef COVERAGE_IGNORE
1860  }
1861  else
1862  {
1863  // DIVERGED_ITS just means it didn't converge in the given maximum number of iterations,
1864  // which is potentially not a problem, as the nonlinear solver may (and often will) still converge.
1865  // Just warn once.
1866  // (Very difficult to cover in normal tests, requires relative and absolute ksp tols to be very small, there
1867  // is no interface for setting both of these. Could be covered by setting up a problem the solver
1868  // finds difficult to solve, but this would be overkill.)
1869  #define COVERAGE_IGNORE
1870  WARN_ONCE_ONLY("Linear solve (within a mechanics solve) didn't converge, but this may not stop nonlinear solve converging")
1871  #undef COVERAGE_IGNORE
1872  }
1873 
1874  // quit if no ksp iterations were done
1875  int num_iters;
1876  KSPGetIterationNumber(solver, &num_iters);
1877  if (num_iters==0)
1878  {
1879  PetscTools::Destroy(solution);
1880  KSPDestroy(PETSC_DESTROY_PARAM(solver));
1881  EXCEPTION("KSP Absolute tolerance was too high, linear system wasn't solved - there will be no decrease in Newton residual. Decrease KspAbsoluteTolerance");
1882  }
1883 
1884 
1885  if(this->mVerbose)
1886  {
1887  Timer::PrintAndReset("KSP Solve");
1888  std::cout << "[" << PetscTools::GetMyRank() << "]: Num iterations = " << num_iters << "\n" << std::flush;
1889  }
1890 
1891  MechanicsEventHandler::EndEvent(MechanicsEventHandler::SOLVE);
1892 
1894  // Update the solution
1895  // Newton method: sol = sol - update, where update=Jac^{-1}*residual
1896  // Newton with damping: sol = sol - s*update, where s is chosen
1897  // such that |residual(sol)| is minimised. Damping is important to
1898  // avoid initial divergence.
1899  //
1900  // Normally, finding the best s from say 0.05,0.1,0.2,..,1.0 is cheap,
1901  // but this is not the case in cardiac electromechanics calculations.
1902  // Therefore, we initially check s=1 (expected to be the best most of the
1903  // time, then s=0.9. If the norm of the residual increases, we assume
1904  // s=1 is the best. Otherwise, check s=0.8 to see if s=0.9 is a local min.
1906  MechanicsEventHandler::BeginEvent(MechanicsEventHandler::UPDATE);
1907  double new_norm_resid = UpdateSolutionUsingLineSearch(solution);
1908  MechanicsEventHandler::EndEvent(MechanicsEventHandler::UPDATE);
1909 
1910  PetscTools::Destroy(solution);
1911  KSPDestroy(PETSC_DESTROY_PARAM(solver));
1912 
1913  return new_norm_resid;
1914 }
1915 
1916 
1917 template<unsigned DIM>
1919 {
1920  if(this->mVerbose)
1921  {
1922  std::cout << "\tTesting s = " << s << ", |f| = " << residNorm << "\n" << std::flush;
1923  }
1924 }
1925 
1926 template<unsigned DIM>
1928 {
1929  double initial_norm_resid = CalculateResidualNorm();
1930  if(this->mVerbose)
1931  {
1932  std::cout << "\tInitial |f| [corresponding to s=0] is " << initial_norm_resid << "\n" << std::flush;
1933  }
1934 
1935  ReplicatableVector update(solution);
1936 
1937  std::vector<double> old_solution = this->mCurrentSolution;
1938 
1939  std::vector<double> damping_values; // = {1.0, 0.9, .., 0.2, 0.1, 0.05} ie size 11
1940  for (unsigned i=10; i>=1; i--)
1941  {
1942  damping_values.push_back((double)i/10.0);
1943  }
1944  damping_values.push_back(0.05);
1945  assert(damping_values.size()==11);
1946 
1948  // let mCurrentSolution = old_solution - damping_val[0]*update; and compute residual
1949  unsigned index = 0;
1950  VectorSum(old_solution, update, -damping_values[index], this->mCurrentSolution);
1951  double current_resid_norm = ComputeResidualAndGetNorm(true);
1952  PrintLineSearchResult(damping_values[index], current_resid_norm);
1953 
1955  // let mCurrentSolution = old_solution - damping_val[1]*update; and compute residual
1956  index = 1;
1957  VectorSum(old_solution, update, -damping_values[index], this->mCurrentSolution);
1958  double next_resid_norm = ComputeResidualAndGetNorm(true);
1959  PrintLineSearchResult(damping_values[index], next_resid_norm);
1960 
1961  index = 2;
1962  // While f(s_next) < f(s_current), [f = residnorm], keep trying new damping values,
1963  // ie exit thus loop when next norm of the residual first increases
1964  while ( (next_resid_norm==DBL_MAX) // the residual is returned as infinity if the deformation is so large to cause exceptions in the material law/EM contraction model
1965  || ( (next_resid_norm < current_resid_norm) && (index<damping_values.size()) ) )
1966  {
1967  current_resid_norm = next_resid_norm;
1968 
1969  // let mCurrentSolution = old_solution - damping_val*update; and compute residual
1970  VectorSum(old_solution, update, -damping_values[index], this->mCurrentSolution);
1971  next_resid_norm = ComputeResidualAndGetNorm(true);
1972  PrintLineSearchResult(damping_values[index], next_resid_norm);
1973 
1974  index++;
1975  }
1976 
1977  unsigned best_index;
1978 
1979  if (index==damping_values.size() && (next_resid_norm < current_resid_norm))
1980  {
1981  // Difficult to come up with large forces/tractions such that it had to
1982  // test right down to s=0.05, but overall doesn't fail.
1983  // The possible damping values have been manually temporarily altered to
1984  // get this code to be called, it appears to work correctly. Even if it
1985  // didn't tests wouldn't fail, they would just be v. slightly less efficient.
1986  #define COVERAGE_IGNORE
1987  // if we exited because we got to the end of the possible damping values, the
1988  // best one was the last one (excl the final index++ at the end)
1989  current_resid_norm = next_resid_norm;
1990  best_index = index-1;
1991  #undef COVERAGE_IGNORE
1992  }
1993  else
1994  {
1995  // else the best one must have been the second last one (excl the final index++ at the end)
1996  // (as we would have exited when the resid norm first increased)
1997  best_index = index-2;
1998  }
1999 
2000  // See documentation for SetTakeFullFirstNewtonStep()
2001  bool full_first_step = mTakeFullFirstNewtonStep && mFirstStep;
2002 
2003 
2004  // Check out best was better than the original residual-norm
2005  if (initial_norm_resid < current_resid_norm && !full_first_step)
2006  {
2007  #define COVERAGE_IGNORE
2008  EXCEPTION("Residual does not appear to decrease in newton direction, quitting");
2009  #undef COVERAGE_IGNORE
2010  }
2011 
2012  // See documentation for SetTakeFullFirstNewtonStep()
2013  if(full_first_step)
2014  {
2015  if(this->mVerbose)
2016  {
2017  std::cout << "\tTaking full first Newton step...\n";
2018  }
2019  best_index = 0;
2020  }
2021 
2022  if(this->mVerbose)
2023  {
2024  std::cout << "\tChoosing s = " << damping_values[best_index] << "\n" << std::flush;
2025  }
2026 
2027 
2029 //
2030 // double l_inf_disp = 0.0;
2031 // double l_inf_pressure = 0.0;
2032 //
2033 // if( this->mCompressibilityType==INCOMPRESSIBLE )
2034 // {
2035 // for(unsigned i=0; i<this->mrQuadMesh.GetNumNodes(); i++)
2036 // {
2037 // for(unsigned j=0; j<DIM; j++)
2038 // {
2039 // double value = update[(DIM+1)*i + j]*damping_values[best_index];
2040 // l_inf_disp = std::max(l_inf_disp, fabs(value));
2041 // }
2042 // l_inf_pressure = std::max(l_inf_pressure, fabs(update[(DIM+1)*i + DIM]*damping_values[best_index]));
2043 // }
2044 // std::cout << "l_inf_disp, l_inf_pressure = " << l_inf_disp << " " << l_inf_pressure << "\n";
2045 // }
2046 // else
2047 // {
2048 // for(unsigned i=0; i<this->mrQuadMesh.GetNumNodes(); i++)
2049 // {
2050 // for(unsigned j=0; j<DIM; j++)
2051 // {
2052 // double value = update[DIM*i + j]*damping_values[best_index];
2053 // l_inf_disp = std::max(l_inf_disp, fabs(value));
2054 // }
2055 // }
2056 // std::cout << "l_inf_disp = " << l_inf_disp << "\n";
2057 // }
2058 
2059 
2060 
2061  VectorSum(old_solution, update, -damping_values[best_index], this->mCurrentSolution);
2062 
2063  mLastDampingValue = damping_values[best_index];
2064  return current_resid_norm;
2065 }
2066 
2067 template<unsigned DIM>
2068 void AbstractNonlinearElasticitySolver<DIM>::PostNewtonStep(unsigned counter, double normResidual)
2069 {
2070 }
2071 
2072 
2073 template<unsigned DIM>
2075 {
2076  mLastDampingValue = 0;
2077 
2078  if (mWriteOutputEachNewtonIteration)
2079  {
2080  this->WriteCurrentSpatialSolution("newton_iteration", "nodes", 0);
2081  }
2082 
2083  // Compute residual
2084  double norm_resid = ComputeResidualAndGetNorm(false);
2085  if(this->mVerbose)
2086  {
2087  std::cout << "\nNorm of residual is " << norm_resid << "\n";
2088  }
2089 
2090  mNumNewtonIterations = 0;
2091  unsigned iteration_number = 1;
2092 
2093  if (tol < 0) // i.e. if wasn't passed in as a parameter
2094  {
2095  tol = NEWTON_REL_TOL*norm_resid;
2096 
2097  #define COVERAGE_IGNORE // not going to have tests in cts for everything
2098  if (tol > MAX_NEWTON_ABS_TOL)
2099  {
2100  tol = MAX_NEWTON_ABS_TOL;
2101  }
2102  if (tol < MIN_NEWTON_ABS_TOL)
2103  {
2104  tol = MIN_NEWTON_ABS_TOL;
2105  }
2106  #undef COVERAGE_IGNORE
2107  }
2108 
2109  if(this->mVerbose)
2110  {
2111  std::cout << "Solving with tolerance " << tol << "\n";
2112  }
2113 
2114  while (norm_resid > tol)
2115  {
2116  if(this->mVerbose)
2117  {
2118  std::cout << "\n-------------------\n"
2119  << "Newton iteration " << iteration_number
2120  << ":\n-------------------\n";
2121  }
2122 
2123  // take newton step (and get returned residual)
2124  mFirstStep = (iteration_number==1);
2125  norm_resid = TakeNewtonStep();
2126 
2127  if(this->mVerbose)
2128  {
2129  std::cout << "Norm of residual is " << norm_resid << "\n";
2130  }
2131 
2132  if (mWriteOutputEachNewtonIteration)
2133  {
2134  this->WriteCurrentSpatialSolution("newton_iteration", "nodes", iteration_number);
2135  }
2136 
2137  mNumNewtonIterations = iteration_number;
2138 
2139  PostNewtonStep(iteration_number,norm_resid);
2140 
2141  iteration_number++;
2142  if (iteration_number==20)
2143  {
2144  #define COVERAGE_IGNORE
2145  EXCEPTION("Not converged after 20 newton iterations, quitting");
2146  #undef COVERAGE_IGNORE
2147  }
2148  }
2149 
2150  if (norm_resid > tol)
2151  {
2152  #define COVERAGE_IGNORE
2153  EXCEPTION("Failed to converge");
2154  #undef COVERAGE_IGNORE
2155  }
2156 }
2157 
2158 
2159 
2160 template<unsigned DIM>
2162 {
2163  return mNumNewtonIterations;
2164 }
2165 
2166 
2167 
2169 // SNES version of the nonlinear solver
2171 
2172 
2173 template<unsigned DIM>
2175 {
2176  // Set up solution guess for residuals
2177  Vec initial_guess;
2178  VecDuplicate(this->mResidualVector, &initial_guess);
2179  double* p_initial_guess;
2180  VecGetArray(initial_guess, &p_initial_guess);
2181  int lo, hi;
2182  VecGetOwnershipRange(initial_guess, &lo, &hi);
2183  for (int global_index=lo; global_index<hi; global_index++)
2184  {
2185  int local_index = global_index - lo;
2186  p_initial_guess[local_index] = this->mCurrentSolution[global_index];
2187  }
2188  VecRestoreArray(initial_guess, &p_initial_guess);
2189  PetscVecTools::Finalise(initial_guess);
2190 
2191  Vec snes_residual_vec;
2192  VecDuplicate(this->mResidualVector, &snes_residual_vec);
2193 
2194  SNES snes;
2195 
2196  SNESCreate(PETSC_COMM_WORLD, &snes);
2197  SNESSetFunction(snes, snes_residual_vec, &AbstractNonlinearElasticitySolver_ComputeResidual<DIM>, this);
2198  SNESSetJacobian(snes, mrJacobianMatrix, this->mPreconditionMatrix, &AbstractNonlinearElasticitySolver_ComputeJacobian<DIM>, this);
2199 #if (PETSC_VERSION_MAJOR == 3 && PETSC_VERSION_MINOR >= 4) //PETSc 3.4 or later
2200  SNESSetType(snes, SNESNEWTONLS);
2201 #else
2202  SNESSetType(snes, SNESLS);
2203 #endif
2204  SNESSetTolerances(snes,1e-5,1e-5,1e-5,PETSC_DEFAULT,PETSC_DEFAULT);
2205  SNESSetMaxLinearSolveFailures(snes,100);
2206 
2207  KSP ksp;
2208  SNESGetKSP(snes, &ksp);
2209 
2210  KSPSetTolerances(ksp, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 1000 /* max iters */); // Note: some machines - max iters seems to be 1000 whatever we give here
2211 
2212  // Set the type of KSP solver (CG, GMRES etc) and preconditioner (ILU, HYPRE, etc)
2213  SetKspSolverAndPcType(ksp);
2214 
2215  if(this->mVerbose)
2216  {
2217  PetscOptionsSetValue("-snes_monitor","");
2218  }
2219  SNESSetFromOptions(snes);
2220 
2221 #if (PETSC_VERSION_MAJOR == 2 && PETSC_VERSION_MINOR == 2) //PETSc 2.2
2222  PetscErrorCode err = SNESSolve(snes, initial_guess);
2223 #else
2224  PetscErrorCode err = SNESSolve(snes, PETSC_NULL, initial_guess);
2225 #endif
2226 
2227  SNESConvergedReason reason;
2228  SNESGetConvergedReason(snes,&reason);
2229 
2230 #define COVERAGE_IGNORE
2231  if (err != 0)
2232  {
2233  std::stringstream err_stream;
2234  err_stream << err;
2235  PetscTools::Destroy(initial_guess);
2236  PetscTools::Destroy(snes_residual_vec);
2237  SNESDestroy(PETSC_DESTROY_PARAM(snes));
2238  EXCEPTION("Nonlinear Solver failed. PETSc error code: "+err_stream.str()+" .");
2239  }
2240 
2241  if (reason < 0)
2242  {
2243  std::stringstream reason_stream;
2244  reason_stream << reason;
2245  PetscTools::Destroy(initial_guess);
2246  PetscTools::Destroy(snes_residual_vec);
2247  SNESDestroy(PETSC_DESTROY_PARAM(snes));
2248  EXCEPTION("Nonlinear Solver did not converge. PETSc reason code: "+reason_stream.str()+" .");
2249  }
2250 #undef COVERAGE_IGNORE
2251 
2252  PetscInt num_iters;
2253  SNESGetIterationNumber(snes,&num_iters);
2254  mNumNewtonIterations = num_iters;
2255 
2256  PetscTools::Destroy(initial_guess);
2257  PetscTools::Destroy(snes_residual_vec);
2258  SNESDestroy(PETSC_DESTROY_PARAM(snes));
2259 }
2260 
2261 
2262 template<unsigned DIM>
2264 {
2265  // Note: AssembleSystem() assumes the current solution is in this->mCurrentSolution and assembles
2266  // this->mResiduaVector and/or this->mrJacobianMatrix. Since PETSc wants us to use the input
2267  // currentGuess, and write the output to residualVector, we have to copy do some copies below.
2268 
2269  ReplicatableVector guess_repl(currentGuess);
2270  for(unsigned i=0; i<guess_repl.GetSize(); i++)
2271  {
2272  this->mCurrentSolution[i] = guess_repl[i];
2273  }
2274  AssembleSystem(true,false);
2275  VecCopy(this->mResidualVector, residualVector);
2276 }
2277 
2278 template<unsigned DIM>
2279 void AbstractNonlinearElasticitySolver<DIM>::ComputeJacobian(Vec currentGuess, Mat* pJacobian, Mat* pPreconditioner)
2280 {
2281  // Note: AssembleSystem() assumes the current solution is in this->mCurrentSolution and assembles
2282  // this->mResiduaVector and/or this->mrJacobianMatrix.
2283  // We need to copy the input currentGuess into the local mCurrentGuess.
2284  // We don't have to copy mrJacobianMatrix to pJacobian, which would be expensive, as they will
2285  // point to the same memory.
2286 
2287  // check Petsc data corresponds to internal Mats
2288  assert(mrJacobianMatrix==*pJacobian);
2289  assert(this->mPreconditionMatrix==*pPreconditioner);
2290 
2291  MechanicsEventHandler::BeginEvent(MechanicsEventHandler::ASSEMBLE);
2292  ReplicatableVector guess_repl(currentGuess);
2293  for(unsigned i=0; i<guess_repl.GetSize(); i++)
2294  {
2295  this->mCurrentSolution[i] = guess_repl[i];
2296  }
2297 
2298  AssembleSystem(false,true);
2299  MechanicsEventHandler::EndEvent(MechanicsEventHandler::ASSEMBLE);
2300 }
2301 
2302 
2303 
2304 template<unsigned DIM>
2305 PetscErrorCode AbstractNonlinearElasticitySolver_ComputeResidual(SNES snes,
2306  Vec currentGuess,
2307  Vec residualVector,
2308  void* pContext)
2309 {
2310  // Extract the solver from the void*
2312  p_solver->ComputeResidual(currentGuess, residualVector);
2313  return 0;
2314 }
2315 
2316 template<unsigned DIM>
2317 #if ( PETSC_VERSION_MAJOR==3 && PETSC_VERSION_MINOR>=5 )
2318  PetscErrorCode AbstractNonlinearElasticitySolver_ComputeJacobian(SNES snes,
2319  Vec currentGuess,
2320  Mat globalJacobian,
2321  Mat preconditioner,
2322  void* pContext)
2323  {
2324  // Extract the solver from the void*
2326  p_solver->ComputeJacobian(currentGuess, &globalJacobian, &preconditioner);
2327  return 0;
2328  }
2329 #else
2330  PetscErrorCode AbstractNonlinearElasticitySolver_ComputeJacobian(SNES snes,
2331  Vec currentGuess,
2332  Mat* pGlobalJacobian,
2333  Mat* pPreconditioner,
2334  MatStructure* pMatStructure,
2335  void* pContext)
2336  {
2337  // Extract the solver from the void*
2339  p_solver->ComputeJacobian(currentGuess, pGlobalJacobian, pPreconditioner);
2340  return 0;
2341  }
2342 #endif
2343 
2344 
2345 // Constant setting definitions
2346 
2347 template<unsigned DIM>
2349 
2350 template<unsigned DIM>
2352 
2353 template<unsigned DIM>
2355 
2356 #endif /*ABSTRACTNONLINEARELASTICITYSOLVER_HPP_*/
void WriteInitialMesh(std::string fileName="")
void SetComputeAverageStressPerElementDuringSolve(bool setComputeAverageStressPerElement=true)
ElementIterator GetElementIteratorBegin(bool skipDeletedElements=true)
void GetElementCentroidStrain(StrainType strainType, Element< DIM, DIM > &rElement, c_matrix< double, DIM, DIM > &rDeformationGradient)
c_vector< double, DIM > & rGetLocation()
Definition: ChastePoint.cpp:76
boost::numeric::ublas::c_matrix< T, 1, 1 > Inverse(const boost::numeric::ublas::c_matrix< T, 1, 1 > &rM)
AbstractNonlinearElasticitySolver(AbstractTetrahedralMesh< DIM, DIM > &rQuadMesh, SolidMechanicsProblemDefinition< DIM > &rProblemDefinition, std::string outputDirectory, CompressibilityType compressibilityType)
unsigned GetNodeGlobalIndex(unsigned localIndex) const
#define PETSC_DESTROY_PARAM(x)
Definition: PetscTools.hpp:69
#define EXCEPTION(message)
Definition: Exception.hpp:143
virtual void AddActiveStressAndStressDerivative(c_matrix< double, DIM, DIM > &rC, unsigned elementIndex, unsigned currentQuadPointGlobalIndex, c_matrix< double, DIM, DIM > &rT, FourthOrderTensor< DIM, DIM, DIM, DIM > &rDTdE, bool addToDTdE)
std::vector< c_vector< double, DIM *(DIM+1)/2 > > mAverageStressesPerElement
static void ComputeTransformedBasisFunctionDerivatives(const ChastePoint< ELEMENT_DIM > &rPoint, const c_matrix< double, ELEMENT_DIM, ELEMENT_DIM > &rInverseJacobian, c_matrix< double, ELEMENT_DIM,(ELEMENT_DIM+1)*(ELEMENT_DIM+2)/2 > &rReturnValue)
T Determinant(const boost::numeric::ublas::c_matrix< T, 1, 1 > &rM)
c_matrix< double, DIM, DIM > GetAverageStressPerElement(unsigned elementIndex)
AbstractNonlinearElasticitySolver< DIM > * mpSolver
virtual void SetupChangeOfBasisMatrix(unsigned elementIndex, unsigned currentQuadPointGlobalIndex)
static void PrintAndReset(std::string message)
Definition: Timer.cpp:70
#define NEVER_REACHED
Definition: Exception.hpp:206
void PrintLineSearchResult(double s, double residNorm)
bool OptionExists(const std::string &rOption)
Node< SPACE_DIM > * GetNode(unsigned localIndex) const
void Visit(Element< DIM, DIM > *pElement, unsigned localIndex, c_vector< double, DIM *DIM > &rData)
static bool IsSequential()
Definition: PetscTools.cpp:91
void SetTakeFullFirstNewtonStep(bool takeFullFirstStep=true)
SolidMechanicsProblemDefinition< DIM > & mrProblemDefinition
virtual void AssembleSystem(bool assembleResidual, bool assembleLinearSystem)=0
StressPerElementWriter(AbstractTetrahedralMesh< DIM, DIM > *pMesh, AbstractNonlinearElasticitySolver< DIM > *pSolver)
virtual void PostNewtonStep(unsigned counter, double normResidual)
void WriteDeformationPositions(std::vector< c_vector< double, DIM > > &rDeformedPositions, unsigned counter)
void SetKspAbsoluteTolerance(double kspAbsoluteTolerance)
static void ComputeBasisFunctions(const ChastePoint< ELEMENT_DIM > &rPoint, c_vector< double,(ELEMENT_DIM+1)*(ELEMENT_DIM+2)/2 > &rReturnValue)
unsigned GetNumNodes() const
static void SwitchWriteMode(Mat matrix)
std::vector< c_vector< double, DIM > > & rGetSpatialSolution()
void AssembleOnBoundaryElementForPressureOnDeformedBc(BoundaryElement< DIM-1, DIM > &rBoundaryElement, c_matrix< double, BOUNDARY_STENCIL_SIZE, BOUNDARY_STENCIL_SIZE > &rAelem, c_vector< double, BOUNDARY_STENCIL_SIZE > &rBelem, bool assembleResidual, bool assembleJacobian, unsigned boundaryConditionIndex)
void ComputeResidual(Vec currentGuess, Vec residualVector)
static void Destroy(Vec &rVec)
Definition: PetscTools.hpp:351
void AssembleOnBoundaryElement(BoundaryElement< DIM-1, DIM > &rBoundaryElement, c_matrix< double, BOUNDARY_STENCIL_SIZE, BOUNDARY_STENCIL_SIZE > &rAelem, c_vector< double, BOUNDARY_STENCIL_SIZE > &rBelem, bool assembleResidual, bool assembleJacobian, unsigned boundaryConditionIndex)
void AddStressToAverageStressPerElement(c_matrix< double, DIM, DIM > &rT, unsigned elementIndex)
c_vector< double, SPACE_DIM > CalculateNormal()
void VectorSum(std::vector< double > &rX, ReplicatableVector &rY, double a, std::vector< double > &rZ)
void WriteCurrentAverageElementStresses(std::string fileName, int counterToAppend=-1)
void ComputeJacobian(Vec currentGuess, Mat *pJacobian, Mat *pPreconditioner)
unsigned GetIndex() const
void SetWriteOutputEachNewtonIteration(bool writeOutputEachNewtonIteration=true)
static void Finalise(Mat matrix)
static CommandLineArguments * Instance()
void WriteCmguiScript(std::string fieldBaseName="", std::string undeformedBaseName="")
virtual void FinishAssembleSystem(bool assembleResidual, bool assembleLinearSystem)
static void Reset()
Definition: Timer.cpp:44
void SetIncludeActiveTension(bool includeActiveTension=true)
static unsigned GetMyRank()
Definition: PetscTools.cpp:114
void WriteCurrentStrains(StrainType strainType, std::string fileName, int counterToAppend=-1)
c_vector< double, SPACE_DIM+1 > CalculateInterpolationWeights(const ChastePoint< SPACE_DIM > &rTestPoint)
Definition: Element.cpp:228
void SetUsePetscDirectSolve(bool usePetscDirectSolve=true)
FourthOrderTensor< DIM, DIM, DIM, DIM > dTdE
std::vector< c_vector< double, DIM > > & rGetDeformedPosition()
static void Finalise(Vec vector)