Chaste Release::3.1
CompressibleNonlinearElasticitySolver.cpp
00001 /*
00002 
00003 Copyright (c) 2005-2012, University of Oxford.
00004 All rights reserved.
00005 
00006 University of Oxford means the Chancellor, Masters and Scholars of the
00007 University of Oxford, having an administrative office at Wellington
00008 Square, Oxford OX1 2JD, UK.
00009 
00010 This file is part of Chaste.
00011 
00012 Redistribution and use in source and binary forms, with or without
00013 modification, are permitted provided that the following conditions are met:
00014  * Redistributions of source code must retain the above copyright notice,
00015    this list of conditions and the following disclaimer.
00016  * Redistributions in binary form must reproduce the above copyright notice,
00017    this list of conditions and the following disclaimer in the documentation
00018    and/or other materials provided with the distribution.
00019  * Neither the name of the University of Oxford nor the names of its
00020    contributors may be used to endorse or promote products derived from this
00021    software without specific prior written permission.
00022 
00023 THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
00024 AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
00025 IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
00026 ARE DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE
00027 LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR
00028 CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE
00029 GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
00030 HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
00031 LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT
00032 OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
00033 
00034 */
00035 
00036 /*
00037  * NOTE ON COMPILATION ERRORS:
00038  *
00039  * (The following applies to IncompressibleNonlinearElasticityAssembler; possibly/probably holds for this class too).
00040  *
00041  * This file won't compile with Intel icpc version 9.1.039, with error message:
00042  * "Terminate with:
00043   (0): internal error: backend signals"
00044  *
00045  * Try recompiling with icpc version 10.0.025.
00046  */
00047 
00048 #include "CompressibleNonlinearElasticitySolver.hpp"
00049 #include "LinearBasisFunction.hpp"
00050 #include "QuadraticBasisFunction.hpp"
00051 #include <algorithm>
00052 
00053 template<size_t DIM>
00054 void CompressibleNonlinearElasticitySolver<DIM>::AssembleSystem(bool assembleResidual,
00055                                                                 bool assembleJacobian)
00056 {
00057     // Check we've actually been asked to do something!
00058     assert(assembleResidual || assembleJacobian);
00059     assert(this->mCurrentSolution.size()==this->mNumDofs);
00060 
00061     // Zero the matrix/vector if it is to be assembled
00062     if (assembleResidual)
00063     {
00064         PetscVecTools::Finalise(this->mResidualVector);
00065         PetscVecTools::Zero(this->mResidualVector);
00066     }
00067     if (assembleJacobian)
00068     {
00069         PetscMatTools::Zero(this->mrJacobianMatrix);
00070         PetscMatTools::Zero(this->mPreconditionMatrix);
00071     }
00072 
00073     c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem;
00074     // The (element) preconditioner matrix: this is the same as the jacobian, but
00075     // with the mass matrix (ie \intgl phi_i phi_j) in the pressure-pressure block.
00076     c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem_precond;
00077     c_vector<double, STENCIL_SIZE> b_elem;
00078 
00079     // Loop over elements
00080     for (typename AbstractTetrahedralMesh<DIM, DIM>::ElementIterator iter = this->mrQuadMesh.GetElementIteratorBegin();
00081          iter != this->mrQuadMesh.GetElementIteratorEnd();
00082          ++iter)
00083     {
00084         Element<DIM, DIM>& element = *iter;
00085 
00086         if (element.GetOwnership() == true)
00087         {
00088             #define COVERAGE_IGNORE
00089             // note: if assembleJacobian only
00090             if(CommandLineArguments::Instance()->OptionExists("-mech_very_verbose") && assembleJacobian)
00091             {
00092                 std::cout << "\r[" << PetscTools::GetMyRank() << "]: Element " << (*iter).GetIndex() << " of " << this->mrQuadMesh.GetNumElements() << std::flush;
00093             }
00094             #undef COVERAGE_IGNORE
00095 
00096             AssembleOnElement(element, a_elem, a_elem_precond, b_elem, assembleResidual, assembleJacobian);
00097 
00101             //for (unsigned i=0; i<STENCIL_SIZE; i++)
00102             //{
00103             //    for (unsigned j=0; j<STENCIL_SIZE; j++)
00104             //    {
00105             //        a_elem(i,j)=1.0;
00106             //    }
00107             //}
00108 
00109             unsigned p_indices[STENCIL_SIZE];
00110             for (unsigned i=0; i<NUM_NODES_PER_ELEMENT; i++)
00111             {
00112                 for (unsigned j=0; j<DIM; j++)
00113                 {
00114                     p_indices[DIM*i+j] = DIM*element.GetNodeGlobalIndex(i) + j;
00115                 }
00116             }
00117 
00118             if (assembleJacobian)
00119             {
00120                 PetscMatTools::AddMultipleValues<STENCIL_SIZE>(this->mrJacobianMatrix, p_indices, a_elem);
00121                 PetscMatTools::AddMultipleValues<STENCIL_SIZE>(this->mPreconditionMatrix, p_indices, a_elem_precond);
00122             }
00123 
00124             if (assembleResidual)
00125             {
00126                 PetscVecTools::AddMultipleValues<STENCIL_SIZE>(this->mResidualVector, p_indices, b_elem);
00127             }
00128         }
00129     }
00130 
00131     // Loop over specified boundary elements and compute surface traction terms
00132     c_vector<double, BOUNDARY_STENCIL_SIZE> b_boundary_elem;
00133     c_matrix<double, BOUNDARY_STENCIL_SIZE, BOUNDARY_STENCIL_SIZE> a_boundary_elem;
00134     if (this->mrProblemDefinition.GetTractionBoundaryConditionType() != NO_TRACTIONS)
00135     {
00136         for (unsigned bc_index=0; bc_index<this->mrProblemDefinition.rGetTractionBoundaryElements().size(); bc_index++)
00137         {
00138             BoundaryElement<DIM-1,DIM>& r_boundary_element = *(this->mrProblemDefinition.rGetTractionBoundaryElements()[bc_index]);
00139 
00140             // If the BCs are tractions applied on a given surface, the boundary integral is independent of u,
00141             // so a_boundary_elem will be zero (no contribution to jacobian).
00142             // If the BCs are normal pressure applied to the deformed body, the boundary depends on the deformation,
00143             // so there is a contribution to the jacobian, and a_boundary_elem is non-zero. Note however that
00144             // the AssembleOnBoundaryElement() method might decide not to include this, as it can actually
00145             // cause divergence if the current guess is not close to the true solution
00146 
00147             this->AssembleOnBoundaryElement(r_boundary_element, a_boundary_elem, b_boundary_elem, assembleResidual, assembleJacobian, bc_index);
00148 
00149             unsigned p_indices[BOUNDARY_STENCIL_SIZE];
00150             for (unsigned i=0; i<NUM_NODES_PER_BOUNDARY_ELEMENT; i++)
00151             {
00152                 for (unsigned j=0; j<DIM; j++)
00153                 {
00154                     p_indices[DIM*i+j] = DIM*r_boundary_element.GetNodeGlobalIndex(i) + j;
00155                 }
00156             }
00157 
00158             if (assembleJacobian)
00159             {
00160                 PetscMatTools::AddMultipleValues<BOUNDARY_STENCIL_SIZE>(this->mrJacobianMatrix, p_indices, a_boundary_elem);
00161                 PetscMatTools::AddMultipleValues<BOUNDARY_STENCIL_SIZE>(this->mPreconditionMatrix, p_indices, a_boundary_elem);
00162             }
00163 
00164             if (assembleResidual)
00165             {
00166                 PetscVecTools::AddMultipleValues<BOUNDARY_STENCIL_SIZE>(this->mResidualVector, p_indices, b_boundary_elem);
00167             }
00168         }
00169     }
00170 
00171     this->FinishAssembleSystem(assembleResidual, assembleJacobian);
00172 }
00173 
00174 template<size_t DIM>
00175 void CompressibleNonlinearElasticitySolver<DIM>::AssembleOnElement(
00176             Element<DIM, DIM>& rElement,
00177             c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElem,
00178             c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElemPrecond,
00179             c_vector<double, STENCIL_SIZE>& rBElem,
00180             bool assembleResidual,
00181             bool assembleJacobian)
00182 {
00183     static c_matrix<double,DIM,DIM> jacobian;
00184     static c_matrix<double,DIM,DIM> inverse_jacobian;
00185     double jacobian_determinant;
00186 
00187     this->mrQuadMesh.GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
00188 
00189     if (assembleJacobian)
00190     {
00191         rAElem.clear();
00192         rAElemPrecond.clear();
00193     }
00194 
00195     if (assembleResidual)
00196     {
00197         rBElem.clear();
00198     }
00199 
00200     // Get the current displacement at the nodes
00201     static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> element_current_displacements;
00202     for (unsigned II=0; II<NUM_NODES_PER_ELEMENT; II++)
00203     {
00204         for (unsigned JJ=0; JJ<DIM; JJ++)
00205         {
00206             element_current_displacements(JJ,II) = this->mCurrentSolution[DIM*rElement.GetNodeGlobalIndex(II) + JJ];
00207         }
00208     }
00209 
00210     // Allocate memory for the basis functions values and derivative values
00211     static c_vector<double, NUM_VERTICES_PER_ELEMENT> linear_phi;
00212     static c_vector<double, NUM_NODES_PER_ELEMENT> quad_phi;
00213     static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> grad_quad_phi;
00214     static c_matrix<double, NUM_NODES_PER_ELEMENT, DIM> trans_grad_quad_phi;
00215 
00216     // Get the material law
00217     AbstractCompressibleMaterialLaw<DIM>* p_material_law
00218        = this->mrProblemDefinition.GetCompressibleMaterialLaw(rElement.GetIndex());
00219 
00220 
00221     static c_matrix<double,DIM,DIM> grad_u; // grad_u = (du_i/dX_M)
00222 
00223     static c_matrix<double,DIM,DIM> F;      // the deformation gradient, F = dx/dX, F_{iM} = dx_i/dX_M
00224     static c_matrix<double,DIM,DIM> C;      // Green deformation tensor, C = F^T F
00225     static c_matrix<double,DIM,DIM> inv_C;  // inverse(C)
00226     static c_matrix<double,DIM,DIM> inv_F;  // inverse(F)
00227     static c_matrix<double,DIM,DIM> T;      // Second Piola-Kirchoff stress tensor (= dW/dE = 2dW/dC)
00228 
00229     static c_matrix<double,DIM,DIM> F_T;    // F*T
00230     static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> F_T_grad_quad_phi; // F*T*grad_quad_phi
00231 
00232     c_vector<double,DIM> body_force;
00233 
00234     static FourthOrderTensor<DIM,DIM,DIM,DIM> dTdE;    // dTdE(M,N,P,Q) = dT_{MN}/dE_{PQ}
00235     static FourthOrderTensor<DIM,DIM,DIM,DIM> dSdF;    // dSdF(M,i,N,j) = dS_{Mi}/dF_{jN}
00236 
00237     static FourthOrderTensor<NUM_NODES_PER_ELEMENT,DIM,DIM,DIM> temp_tensor;
00238     static FourthOrderTensor<NUM_NODES_PER_ELEMENT,DIM,NUM_NODES_PER_ELEMENT,DIM> dSdF_quad_quad;
00239 
00240     static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> temp_matrix;
00241     static c_matrix<double,NUM_NODES_PER_ELEMENT,DIM> grad_quad_phi_times_invF;
00242 
00243     if(this->mSetComputeAverageStressPerElement)
00244     {
00245         this->mAverageStressesPerElement[rElement.GetIndex()] = zero_vector<double>(DIM*(DIM+1)/2);
00246     }
00247 
00248     // Loop over Gauss points
00249     for (unsigned quadrature_index=0; quadrature_index < this->mpQuadratureRule->GetNumQuadPoints(); quadrature_index++)
00250     {
00251         // This is needed by the cardiac mechanics solver
00252         unsigned current_quad_point_global_index =   rElement.GetIndex()*this->mpQuadratureRule->GetNumQuadPoints()
00253                                                    + quadrature_index;
00254 
00255         double wJ = jacobian_determinant * this->mpQuadratureRule->GetWeight(quadrature_index);
00256 
00257         const ChastePoint<DIM>& quadrature_point = this->mpQuadratureRule->rGetQuadPoint(quadrature_index);
00258 
00259         // Set up basis function information
00260         LinearBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, linear_phi);
00261         QuadraticBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, quad_phi);
00262         QuadraticBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(quadrature_point, inverse_jacobian, grad_quad_phi);
00263         trans_grad_quad_phi = trans(grad_quad_phi);
00264 
00265         // Get the body force, interpolating X if necessary
00266         if (assembleResidual)
00267         {
00268             switch (this->mrProblemDefinition.GetBodyForceType())
00269             {
00270                 case FUNCTIONAL_BODY_FORCE:
00271                 {
00272                     c_vector<double,DIM> X = zero_vector<double>(DIM);
00273                     // interpolate X (using the vertices and the /linear/ bases, as no curvilinear elements
00274                     for (unsigned node_index=0; node_index<NUM_VERTICES_PER_ELEMENT; node_index++)
00275                     {
00276                         X += linear_phi(node_index)*this->mrQuadMesh.GetNode( rElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
00277                     }
00278                     body_force = this->mrProblemDefinition.EvaluateBodyForceFunction(X, this->mCurrentTime);
00279                     break;
00280                 }
00281                 case CONSTANT_BODY_FORCE:
00282                 {
00283                     body_force = this->mrProblemDefinition.GetConstantBodyForce();
00284                     break;
00285                 }
00286                 default:
00287                     NEVER_REACHED;
00288             }
00289         }
00290 
00291         // Interpolate grad_u
00292         grad_u = zero_matrix<double>(DIM,DIM);
00293         for (unsigned node_index=0; node_index<NUM_NODES_PER_ELEMENT; node_index++)
00294         {
00295             for (unsigned i=0; i<DIM; i++)
00296             {
00297                 for (unsigned M=0; M<DIM; M++)
00298                 {
00299                     grad_u(i,M) += grad_quad_phi(M,node_index)*element_current_displacements(i,node_index);
00300                 }
00301             }
00302         }
00303 
00304         // Calculate C, inv(C) and T
00305         for (unsigned i=0; i<DIM; i++)
00306         {
00307             for (unsigned M=0; M<DIM; M++)
00308             {
00309                 F(i,M) = (i==M?1:0) + grad_u(i,M);
00310             }
00311         }
00312 
00313         C = prod(trans(F),F);
00314         inv_C = Inverse(C);
00315         inv_F = Inverse(F);
00316 
00317         // Compute the passive stress, and dTdE corresponding to passive stress
00318         this->SetupChangeOfBasisMatrix(rElement.GetIndex(), current_quad_point_global_index);
00319         p_material_law->SetChangeOfBasisMatrix(this->mChangeOfBasisMatrix);
00320         p_material_law->ComputeStressAndStressDerivative(C, inv_C, 0.0, T, dTdE, assembleJacobian);
00321 
00322         if(this->mIncludeActiveTension)
00323         {
00324             // Add any active stresses, if there are any. Requires subclasses to overload this method,
00325             // see for example the cardiac mechanics assemblers.
00326             this->AddActiveStressAndStressDerivative(C, rElement.GetIndex(), current_quad_point_global_index,
00327                                                      T, dTdE, assembleJacobian);
00328         }
00329 
00330         if(this->mSetComputeAverageStressPerElement)
00331         {
00332             this->AddStressToAverageStressPerElement(T,rElement.GetIndex());
00333         }
00334 
00335         // Residual vector
00336         if (assembleResidual)
00337         {
00338             F_T = prod(F,T);
00339             F_T_grad_quad_phi = prod(F_T, grad_quad_phi);
00340 
00341             for (unsigned index=0; index<NUM_NODES_PER_ELEMENT*DIM; index++)
00342             {
00343                 unsigned spatial_dim = index%DIM;
00344                 unsigned node_index = (index-spatial_dim)/DIM;
00345 
00346                 rBElem(index) += - this->mrProblemDefinition.GetDensity()
00347                                  * body_force(spatial_dim)
00348                                  * quad_phi(node_index)
00349                                  * wJ;
00350 
00351                 // The T(M,N)*F(spatial_dim,M)*grad_quad_phi(N,node_index) term
00352                 rBElem(index) +=   F_T_grad_quad_phi(spatial_dim,node_index)
00353                                  * wJ;
00354             }
00355         }
00356 
00357         // Jacobian matrix
00358         if (assembleJacobian)
00359         {
00360             // Save trans(grad_quad_phi) * invF
00361             grad_quad_phi_times_invF = prod(trans_grad_quad_phi, inv_F);
00362 
00364             // Set up the tensor dSdF
00365             //
00366             // dSdF as a function of T and dTdE (which is what the material law returns) is given by:
00367             //
00368             // dS_{Mi}/dF_{jN} = (dT_{MN}/dC_{PQ}+dT_{MN}/dC_{PQ}) F{iP} F_{jQ}  + T_{MN} delta_{ij}
00369             //
00370             // todo1: this should probably move into the material law (but need to make sure
00371             // memory is handled efficiently
00372             // todo2: get material law to return this immediately, not dTdE
00374 
00375             // Set up the tensor 0.5(dTdE(M,N,P,Q) + dTdE(M,N,Q,P))
00376             for (unsigned M=0; M<DIM; M++)
00377             {
00378                 for (unsigned N=0; N<DIM; N++)
00379                 {
00380                     for (unsigned P=0; P<DIM; P++)
00381                     {
00382                         for (unsigned Q=0; Q<DIM; Q++)
00383                         {
00384                             // this is NOT dSdF, just using this as storage space
00385                             dSdF(M,N,P,Q) = 0.5*(dTdE(M,N,P,Q) + dTdE(M,N,Q,P));
00386                         }
00387                     }
00388                 }
00389             }
00390 
00391             // This is NOT dTdE, just reusing memory. A^{MdPQ}  = F^d_N * dTdE_sym^{MNPQ}
00392             dTdE.template SetAsContractionOnSecondDimension<DIM>(F, dSdF);
00393 
00394             // dSdF{MdPe} := F^d_N * F^e_Q * dTdE_sym^{MNPQ}
00395             dSdF.template SetAsContractionOnFourthDimension<DIM>(F, dTdE);
00396 
00397             // Now add the T_{MN} delta_{ij} term
00398             for (unsigned M=0; M<DIM; M++)
00399             {
00400                 for (unsigned N=0; N<DIM; N++)
00401                 {
00402                     for (unsigned i=0; i<DIM; i++)
00403                     {
00404                         dSdF(M,i,N,i) += T(M,N);
00405                     }
00406                 }
00407             }
00408 
00410             // Set up the tensor
00411             //   dSdF_quad_quad(node_index1, spatial_dim1, node_index2, spatial_dim2)
00412             //            =    dS_{M,spatial_dim1}/d_F{spatial_dim2,N}
00413             //               * grad_quad_phi(M,node_index1)
00414             //               * grad_quad_phi(P,node_index2)
00415             //
00416             //            =    dSdF(M,spatial_index1,N,spatial_index2)
00417             //               * grad_quad_phi(M,node_index1)
00418             //               * grad_quad_phi(P,node_index2)
00419             //
00421             temp_tensor.template SetAsContractionOnFirstDimension<DIM>(trans_grad_quad_phi, dSdF);
00422             dSdF_quad_quad.template SetAsContractionOnThirdDimension<DIM>(trans_grad_quad_phi, temp_tensor);
00423 
00424             for (unsigned index1=0; index1<NUM_NODES_PER_ELEMENT*DIM; index1++)
00425             {
00426                 unsigned spatial_dim1 = index1%DIM;
00427                 unsigned node_index1 = (index1-spatial_dim1)/DIM;
00428 
00429                 for (unsigned index2=0; index2<NUM_NODES_PER_ELEMENT*DIM; index2++)
00430                 {
00431                     unsigned spatial_dim2 = index2%DIM;
00432                     unsigned node_index2 = (index2-spatial_dim2)/DIM;
00433 
00434                     // The dSdF*grad_quad_phi*grad_quad_phi term
00435                     rAElem(index1,index2) +=   dSdF_quad_quad(node_index1,spatial_dim1,node_index2,spatial_dim2)
00436                                              * wJ;
00437                 }
00438             }
00439         }
00440     }
00441 
00442     rAElemPrecond.clear();
00443     if (assembleJacobian)
00444     {
00445         rAElemPrecond = rAElem;
00446     }
00447 
00448     if(this->mSetComputeAverageStressPerElement)
00449     {
00450         for(unsigned i=0; i<DIM*(DIM+1)/2; i++)
00451         {
00452             this->mAverageStressesPerElement[rElement.GetIndex()](i) /= this->mpQuadratureRule->GetNumQuadPoints();
00453         }
00454     }
00455 
00456 }
00457 
00458 
00459 template<size_t DIM>
00460 CompressibleNonlinearElasticitySolver<DIM>::CompressibleNonlinearElasticitySolver(QuadraticMesh<DIM>& rQuadMesh,
00461                                                                                   SolidMechanicsProblemDefinition<DIM>& rProblemDefinition,
00462                                                                                   std::string outputDirectory)
00463     : AbstractNonlinearElasticitySolver<DIM>(rQuadMesh,
00464                                              rProblemDefinition,
00465                                              outputDirectory,
00466                                              COMPRESSIBLE)
00467 {
00468     if(rProblemDefinition.GetCompressibilityType() != COMPRESSIBLE)
00469     {
00470         EXCEPTION("SolidMechanicsProblemDefinition object contains incompressible material laws");
00471     }
00472 }
00473 
00474 
00475 template<size_t DIM>
00476 CompressibleNonlinearElasticitySolver<DIM>::~CompressibleNonlinearElasticitySolver()
00477 {
00478 }
00479 
00481 // Explicit instantiation
00483 
00484 template class CompressibleNonlinearElasticitySolver<2>;
00485 template class CompressibleNonlinearElasticitySolver<3>;