NonlinearElasticitySolver.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2010
00004 
00005 University of Oxford means the Chancellor, Masters and Scholars of the
00006 University of Oxford, having an administrative office at Wellington
00007 Square, Oxford OX1 2JD, UK.
00008 
00009 This file is part of Chaste.
00010 
00011 Chaste is free software: you can redistribute it and/or modify it
00012 under the terms of the GNU Lesser General Public License as published
00013 by the Free Software Foundation, either version 2.1 of the License, or
00014 (at your option) any later version.
00015 
00016 Chaste is distributed in the hope that it will be useful, but WITHOUT
00017 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018 FITNESS FOR A PARTICULAR PURPOSE.  See the GNU Lesser General Public
00019 License for more details. The offer of Chaste under the terms of the
00020 License is subject to the License being interpreted in accordance with
00021 English Law and subject to any action against the University of Oxford
00022 being under the jurisdiction of the English Courts.
00023 
00024 You should have received a copy of the GNU Lesser General Public License
00025 along with Chaste. If not, see <http://www.gnu.org/licenses/>.
00026 
00027 */
00028 
00029 
00030 
00031 /*
00032  * NOTE ON COMPILATION ERRORS:
00033  *
00034  * This file won't compile with Intel icpc version 9.1.039, with error message:
00035  * "Terminate with:
00036   (0): internal error: backend signals"
00037  *
00038  * Try recompiling with icpc version 10.0.025.
00039  */
00040 
00041 #include "NonlinearElasticitySolver.hpp"
00042 #include "LinearBasisFunction.hpp"
00043 #include "QuadraticBasisFunction.hpp"
00044 #include <algorithm>
00045 
00046 //#include "Debug.hpp"
00047 template<size_t DIM>
00048 void NonlinearElasticitySolver<DIM>::AssembleSystem(bool assembleResidual,
00049                                                     bool assembleJacobian)
00050 {
00051     // Check we've actually been asked to do something!
00052     assert(assembleResidual || assembleJacobian);
00053     assert(this->mCurrentSolution.size()==this->mNumDofs);
00054 
00055     // Zero the matrix/vector if it is to be assembled
00056     if (assembleResidual)
00057     {
00058         this->mpLinearSystem->ZeroRhsVector();
00059     }
00060     if (assembleJacobian)
00061     {
00062         this->mpLinearSystem->ZeroLhsMatrix();
00063         this->mpPreconditionMatrixLinearSystem->ZeroLhsMatrix();
00064     }
00065 
00066     c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem;
00067     // the (element) preconditioner matrix: this is the same as the jacobian, but
00068     // with the mass matrix (ie \intgl phi_i phi_j) in the pressure-pressure block.
00069     c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem_precond;
00070     c_vector<double, STENCIL_SIZE> b_elem;
00071 
00073     // loop over elements
00075     for (typename AbstractTetrahedralMesh<DIM, DIM>::ElementIterator iter = mpQuadMesh->GetElementIteratorBegin();
00076          iter != mpQuadMesh->GetElementIteratorEnd();
00077          ++iter)
00078     {
00079         #ifdef MECH_VERY_VERBOSE
00080         if (assembleJacobian) // && ((*iter).GetIndex()%500==0))
00081         {
00082             std::cout << "\r[" << PetscTools::GetMyRank() << "]: Element " << (*iter).GetIndex() << " of " << mpQuadMesh->GetNumElements() << std::flush;
00083         }
00084         #endif
00085 
00086         Element<DIM, DIM>& element = *iter;
00087 
00088         if (element.GetOwnership() == true)
00089         {
00090             AssembleOnElement(element, a_elem, a_elem_precond, b_elem, assembleResidual, assembleJacobian);
00091             
00095             //for(unsigned i=0; i<STENCIL_SIZE; i++)
00096             //{
00097             //    for(unsigned j=0; j<STENCIL_SIZE; j++)
00098             //    {
00099             //        a_elem(i,j)=1.0;
00100             //    }
00101             //}
00102             
00103 
00104             unsigned p_indices[STENCIL_SIZE];
00105             for (unsigned i=0; i<NUM_NODES_PER_ELEMENT; i++)
00106             {
00107                 for (unsigned j=0; j<DIM; j++)
00108                 {
00109                     p_indices[DIM*i+j] = DIM*element.GetNodeGlobalIndex(i) + j;
00110                 }
00111             }
00112 
00113             for (unsigned i=0; i<NUM_VERTICES_PER_ELEMENT; i++)
00114             {
00115                 p_indices[DIM*NUM_NODES_PER_ELEMENT + i] = DIM*mpQuadMesh->GetNumNodes() + element.GetNodeGlobalIndex(i);
00116             }
00117 
00118             if (assembleJacobian)
00119             {
00120                 this->mpLinearSystem->AddLhsMultipleValues(p_indices, a_elem);
00121                 this->mpPreconditionMatrixLinearSystem->AddLhsMultipleValues(p_indices, a_elem_precond);
00122             }
00123 
00124             if (assembleResidual)
00125             {
00126                 this->mpLinearSystem->AddRhsMultipleValues(p_indices, b_elem);
00127             }
00128         }
00129     }
00130 
00132     // loop over specified boundary elements and compute
00133     // surface traction terms
00135     c_vector<double, BOUNDARY_STENCIL_SIZE> b_boundary_elem;
00136     c_matrix<double, BOUNDARY_STENCIL_SIZE, BOUNDARY_STENCIL_SIZE> a_boundary_elem;
00137     if (mBoundaryElements.size()>0)
00138     {
00139         for (unsigned i=0; i<mBoundaryElements.size(); i++)
00140         {
00141             BoundaryElement<DIM-1,DIM>& r_boundary_element = *(mBoundaryElements[i]);
00142             AssembleOnBoundaryElement(r_boundary_element, a_boundary_elem, b_boundary_elem, this->mSurfaceTractions[i], assembleResidual, assembleJacobian);
00143 
00144             unsigned p_indices[BOUNDARY_STENCIL_SIZE];
00145             for (unsigned i=0; i<NUM_NODES_PER_BOUNDARY_ELEMENT; i++)
00146             {
00147                 for (unsigned j=0; j<DIM; j++)
00148                 {
00149                     p_indices[DIM*i+j] = DIM*r_boundary_element.GetNodeGlobalIndex(i) + j;
00150                 }
00151             }
00152 
00153             for (unsigned i=0; i<DIM /*vertices per boundary elem */; i++)
00154             {
00155                 p_indices[DIM*NUM_NODES_PER_BOUNDARY_ELEMENT + i] = DIM*mpQuadMesh->GetNumNodes() + r_boundary_element.GetNodeGlobalIndex(i);
00156             }
00157 
00158             if (assembleJacobian)
00159             {
00160                 this->mpLinearSystem->AddLhsMultipleValues(p_indices, a_boundary_elem);
00161                 this->mpPreconditionMatrixLinearSystem->AddLhsMultipleValues(p_indices, a_boundary_elem);
00162             }
00163 
00164             if (assembleResidual)
00165             {
00166                 this->mpLinearSystem->AddRhsMultipleValues(p_indices, b_boundary_elem);
00167             }
00168 
00169             // some extra checking
00170             if (DIM==2)
00171             {
00172                 assert(8==BOUNDARY_STENCIL_SIZE);
00173                 assert(b_boundary_elem(6)==0);
00174                 assert(b_boundary_elem(7)==0);
00175             }
00176         }
00177     }
00178 
00179     if (assembleResidual)
00180     {
00181         this->mpLinearSystem->AssembleRhsVector();
00182     }
00183     if (assembleJacobian)
00184     {
00185         this->mpLinearSystem->AssembleIntermediateLhsMatrix();
00186         this->mpPreconditionMatrixLinearSystem->AssembleIntermediateLhsMatrix();
00187     }
00188 
00189     // Apply Dirichlet boundary conditions
00190     this->ApplyBoundaryConditions(assembleJacobian);
00191 
00192     if (assembleResidual)
00193     {
00194         this->mpLinearSystem->AssembleRhsVector();
00195     }
00196     if (assembleJacobian)
00197     {
00198         this->mpLinearSystem->AssembleFinalLhsMatrix();
00199         this->mpPreconditionMatrixLinearSystem->AssembleFinalLhsMatrix();
00200     }
00201 }
00202 
00203 template<size_t DIM>
00204 void NonlinearElasticitySolver<DIM>::AssembleOnElement(
00205             Element<DIM, DIM>& rElement,
00206             c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElem,
00207             c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElemPrecond,
00208             c_vector<double, STENCIL_SIZE>& rBElem,
00209             bool assembleResidual,
00210             bool assembleJacobian)
00211 {
00212     static c_matrix<double,DIM,DIM> jacobian;
00213     static c_matrix<double,DIM,DIM> inverse_jacobian;
00214     double jacobian_determinant;
00215     
00216     mpQuadMesh->GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
00217 
00218     if (assembleJacobian)
00219     {
00220         rAElem.clear();
00221         rAElemPrecond.clear();
00222     }
00223 
00224     if (assembleResidual)
00225     {
00226         rBElem.clear();
00227     }
00228 
00230     // Get the current displacement at the nodes
00232     static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> element_current_displacements;
00233     static c_vector<double,NUM_VERTICES_PER_ELEMENT> element_current_pressures;
00234     for (unsigned II=0; II<NUM_NODES_PER_ELEMENT; II++)
00235     {
00236         for (unsigned JJ=0; JJ<DIM; JJ++)
00237         {
00238             element_current_displacements(JJ,II) = this->mCurrentSolution[DIM*rElement.GetNodeGlobalIndex(II) + JJ];
00239         }
00240     }
00241 
00243     // Get the current pressure at the vertices
00245     for (unsigned II=0; II<NUM_VERTICES_PER_ELEMENT; II++)
00246     {
00247         element_current_pressures(II) = this->mCurrentSolution[DIM*mpQuadMesh->GetNumNodes() + rElement.GetNodeGlobalIndex(II)];
00248     }
00249 
00250     // allocate memory for the basis functions values and derivative values
00251     static c_vector<double, NUM_VERTICES_PER_ELEMENT> linear_phi;
00252     static c_vector<double, NUM_NODES_PER_ELEMENT> quad_phi;
00253     static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> grad_quad_phi;
00254     static c_matrix<double, NUM_NODES_PER_ELEMENT, DIM> trans_grad_quad_phi;
00255 
00256 
00257     // get the material law
00258     AbstractIncompressibleMaterialLaw<DIM>* p_material_law;
00259     if (this->mMaterialLaws.size()==1)
00260     {
00261         // homogeneous
00262         p_material_law = this->mMaterialLaws[0];
00263     }
00264     else
00265     {
00266         // heterogeneous
00267         #define COVERAGE_IGNORE // not going to have tests in cts for everything
00268         p_material_law = this->mMaterialLaws[rElement.GetIndex()];
00269         #undef COVERAGE_IGNORE
00270     }
00271     
00272     
00273     static c_matrix<double,DIM,DIM> grad_u; // grad_u = (du_i/dX_M)
00274             
00275     static c_matrix<double,DIM,DIM> F;      // the deformation gradient, F = dx/dX, F_{iM} = dx_i/dX_M
00276     static c_matrix<double,DIM,DIM> C;      // Green deformation tensor, C = F^T F
00277     static c_matrix<double,DIM,DIM> inv_C;  // inverse(C)
00278     static c_matrix<double,DIM,DIM> inv_F;  // inverse(F)
00279     static c_matrix<double,DIM,DIM> T;      // Second Piola-Kirchoff stress tensor (= dW/dE = 2dW/dC) 
00280 
00281     static c_matrix<double,DIM,DIM> F_T;    // F*T
00282     static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> F_T_grad_quad_phi; // F*T*grad_quad_phi
00283 
00284     c_vector<double,DIM> body_force;
00285 
00286     static FourthOrderTensor<DIM,DIM,DIM,DIM> dTdE;    // dTdE(M,N,P,Q) = dT_{MN}/dE_{PQ} 
00287     static FourthOrderTensor<DIM,DIM,DIM,DIM> dSdF;    // dSdF(M,i,N,j) = dS_{Mi}/dF_{jN} 
00288 
00289     static FourthOrderTensor<NUM_NODES_PER_ELEMENT,DIM,DIM,DIM> temp_tensor;
00290     static FourthOrderTensor<NUM_NODES_PER_ELEMENT,DIM,NUM_NODES_PER_ELEMENT,DIM> dSdF_quad_quad;
00291 
00292     static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> temp_matrix;
00293     static c_matrix<double,NUM_NODES_PER_ELEMENT,DIM> grad_quad_phi_times_invF;
00294 
00295 
00296 
00302     for (unsigned quadrature_index=0; quadrature_index < this->mpQuadratureRule->GetNumQuadPoints(); quadrature_index++)
00303     {
00304         // this is needed by the cardiac mechanics solver
00305         unsigned current_quad_point_global_index =   rElement.GetIndex()*this->mpQuadratureRule->GetNumQuadPoints()
00306                                                    + quadrature_index;
00307 
00308         double wJ = jacobian_determinant * this->mpQuadratureRule->GetWeight(quadrature_index);
00309 
00310         const ChastePoint<DIM>& quadrature_point = this->mpQuadratureRule->rGetQuadPoint(quadrature_index);
00311 
00313         // set up basis function info
00315         LinearBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, linear_phi);
00316         QuadraticBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, quad_phi);
00317         QuadraticBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(quadrature_point, inverse_jacobian, grad_quad_phi);
00318         trans_grad_quad_phi = trans(grad_quad_phi);
00319         
00320         
00322         // get the body force, interpolating X if necessary
00324 
00325         if(assembleResidual)
00326         {
00327             if (this->mUsingBodyForceFunction)
00328             {
00329                 c_vector<double,DIM> X = zero_vector<double>(DIM);
00330                 // interpolate X (using the vertices and the /linear/ bases, as no curvilinear elements
00331                 for (unsigned node_index=0; node_index<NUM_VERTICES_PER_ELEMENT; node_index++)
00332                 {
00333                     X += linear_phi(node_index)*mpQuadMesh->GetNode( rElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
00334                 }
00335                 body_force = (*(this->mpBodyForceFunction))(X);
00336             }
00337             else
00338             {
00339                 body_force = this->mBodyForce;
00340             }
00341         }
00342 
00344         // interpolate grad_u and p
00346         grad_u = zero_matrix<double>(DIM,DIM); 
00347 
00348         for (unsigned node_index=0; node_index<NUM_NODES_PER_ELEMENT; node_index++)
00349         {
00350             for (unsigned i=0; i<DIM; i++)
00351             {
00352                 for (unsigned M=0; M<DIM; M++)
00353                 {
00354                     grad_u(i,M) += grad_quad_phi(M,node_index)*element_current_displacements(i,node_index);
00355                 }
00356             }
00357         }
00358 
00359         double pressure = 0;
00360         for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00361         {
00362             pressure += linear_phi(vertex_index)*element_current_pressures(vertex_index);
00363         }
00364 
00365 
00367         // calculate C, inv(C) and T
00369         for (unsigned i=0; i<DIM; i++)
00370         {
00371             for (unsigned M=0; M<DIM; M++)
00372             {
00373                 F(i,M) = (i==M?1:0) + grad_u(i,M);
00374             }
00375         }
00376 
00377         C = prod(trans(F),F);
00378         inv_C = Inverse(C);
00379         inv_F = Inverse(F);
00380 
00381         double detF = Determinant(F);
00382 
00383         ComputeStressAndStressDerivative(p_material_law, C, inv_C, pressure, rElement.GetIndex(), current_quad_point_global_index,
00384                                          T, dTdE, assembleJacobian);
00385 
00386 
00388         // residual vector
00390         if (assembleResidual)
00391         {
00392             F_T = prod(F,T);
00393             F_T_grad_quad_phi = prod(F_T, grad_quad_phi);
00394             
00395             for (unsigned index=0; index<NUM_NODES_PER_ELEMENT*DIM; index++)
00396             {
00397                 unsigned spatial_dim = index%DIM;
00398                 unsigned node_index = (index-spatial_dim)/DIM;
00399 
00400                 rBElem(index) +=  - this->mDensity
00401                                   * body_force(spatial_dim)
00402                                   * quad_phi(node_index)
00403                                   * wJ;
00404                                   
00405                 // the  T(M,N)*F(spatial_dim,M)*grad_quad_phi(N,node_index)  term                  
00406                 rBElem(index) +=   F_T_grad_quad_phi(spatial_dim,node_index)
00407                                  * wJ;
00408             }
00409 
00410             for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00411             {
00412                 rBElem( NUM_NODES_PER_ELEMENT*DIM + vertex_index ) +=   linear_phi(vertex_index)
00413                                                                       * (detF - 1)
00414                                                                       * wJ;
00415             }
00416         }
00417         
00418         
00420         // Jacobian matrix
00422         if (assembleJacobian) 
00423         {
00424             // save trans(grad_quad_phi) * invF
00425             grad_quad_phi_times_invF = prod(trans_grad_quad_phi, inv_F);
00426             
00428             // Set up the tensor dSdF
00429             //
00430             // dSdF as a function of T and dTdE (which is what the material law returns) is given by:
00431             //
00432             // dS_{Mi}/dF_{jN} = (dT_{MN}/dC_{PQ}+dT_{MN}/dC_{PQ}) F{iP} F_{jQ}  + T_{MN} delta_{ij}
00433             //
00434             // todo1: this should probably move into the material law (but need to make sure 
00435             // memory is handled efficiently
00436             // todo2: get material law to return this immediately, not dTdE
00438 
00439             // set up the tensor 0.5(dTdE(M,N,P,Q) + dTdE(M,N,Q,P))
00440             for (unsigned M=0; M<DIM; M++)
00441             {
00442                 for (unsigned N=0; N<DIM; N++)
00443                 {
00444                     for (unsigned P=0; P<DIM; P++)
00445                     {
00446                         for (unsigned Q=0; Q<DIM; Q++)
00447                         {
00448                             // this is NOT dSdF, just using this as storage space
00449                             dSdF(M,N,P,Q) = 0.5*(dTdE(M,N,P,Q) + dTdE(M,N,Q,P));
00450                         }
00451                     }
00452                 }
00453             }
00454             // This is NOT dTdE, just reusing memory. A^{MdPQ}  = F^d_N * dTdE_sym^{MNPQ}
00455             dTdE.template SetAsContractionOnSecondDimension<DIM>(F, dSdF);  
00456             // dSdF{MdPe} := F^d_N * F^e_Q * dTdE_sym^{MNPQ}
00457             dSdF.template SetAsContractionOnFourthDimension<DIM>(F, dTdE);  
00458             
00459             // now add the T_{MN} delta_{ij} term
00460             for(unsigned M=0; M<DIM; M++)
00461             {
00462                 for(unsigned N=0; N<DIM; N++)
00463                 {
00464                     for(unsigned i=0; i<DIM; i++)
00465                     {
00466                         dSdF(M,i,N,i) += T(M,N);
00467                     }
00468                 }
00469             }
00470             
00471 
00473             // Set up the tensor  
00474             //   dSdF_quad_quad(node_index1, spatial_dim1, node_index2, spatial_dim2)
00475             //            =    dS_{M,spatial_dim1}/d_F{spatial_dim2,N}      
00476             //               * grad_quad_phi(M,node_index1) 
00477             //               * grad_quad_phi(P,node_index2)
00478             //
00479             //            =    dSdF(M,spatial_index1,N,spatial_index2)
00480             //               * grad_quad_phi(M,node_index1) 
00481             //               * grad_quad_phi(P,node_index2)
00482             //            
00484             temp_tensor.template SetAsContractionOnFirstDimension<DIM>( trans_grad_quad_phi, dSdF );
00485             dSdF_quad_quad.template SetAsContractionOnThirdDimension<DIM>( trans_grad_quad_phi, temp_tensor );
00486 
00487 
00488 
00489 
00490             for (unsigned index1=0; index1<NUM_NODES_PER_ELEMENT*DIM; index1++)
00491             {
00492                 unsigned spatial_dim1 = index1%DIM;
00493                 unsigned node_index1 = (index1-spatial_dim1)/DIM;
00494 
00495 
00496                 for (unsigned index2=0; index2<NUM_NODES_PER_ELEMENT*DIM; index2++)
00497                 {
00498                     unsigned spatial_dim2 = index2%DIM;
00499                     unsigned node_index2 = (index2-spatial_dim2)/DIM;
00500 
00501                     // the dSdF*grad_quad_phi*grad_quad_phi term
00502                     rAElem(index1,index2)  +=   dSdF_quad_quad(node_index1,spatial_dim1,node_index2,spatial_dim2)
00503                                               * wJ;
00504                 }
00505 
00506                 for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00507                 {
00508                     unsigned index2 = NUM_NODES_PER_ELEMENT*DIM + vertex_index;
00509         
00510                     // the -invF(M,spatial_dim1)*grad_quad_phi(M,node_index1)*linear_phi(vertex_index) term
00511                     rAElem(index1,index2)  +=  - grad_quad_phi_times_invF(node_index1,spatial_dim1)
00512                                                * linear_phi(vertex_index)
00513                                                * wJ;
00514                 }
00515             }
00516 
00517             for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00518             {
00519                 unsigned index1 = NUM_NODES_PER_ELEMENT*DIM + vertex_index;
00520 
00521                 for (unsigned index2=0; index2<NUM_NODES_PER_ELEMENT*DIM; index2++)
00522                 {
00523                     unsigned spatial_dim2 = index2%DIM;
00524                     unsigned node_index2 = (index2-spatial_dim2)/DIM;
00525 
00526                     // same as (negative of) the opposite block (ie a few lines up), except for detF
00527                     rAElem(index1,index2) +=   detF
00528                                              * grad_quad_phi_times_invF(node_index2,spatial_dim2)
00529                                              * linear_phi(vertex_index)
00530                                              * wJ;
00531                 }
00532 
00534                 // Preconditioner matrix
00535                 // Fill the mass matrix (ie \intgl phi_i phi_j) in the
00536                 // pressure-pressure block. Note, the rest of the
00537                 // entries are filled in at the end
00539                 for (unsigned vertex_index2=0; vertex_index2<NUM_VERTICES_PER_ELEMENT; vertex_index2++)
00540                 {
00541                     unsigned index2 = NUM_NODES_PER_ELEMENT*DIM + vertex_index2;
00542                     rAElemPrecond(index1,index2) +=   linear_phi(vertex_index)
00543                                                     * linear_phi(vertex_index2)
00544                                                     * wJ;
00545                 }
00546             }
00547         }
00548     }
00549 
00550 
00551     if (assembleJacobian)
00552     {
00553         // Fill in the other blocks of the preconditioner matrix, by adding 
00554         // the jacobian matrix (this doesn't effect the pressure-pressure block 
00555         // of rAElemPrecond as the pressure-pressure block of rAElem is zero),
00556         // and the zero a block.
00557         //
00558         // The following altogether gives the preconditioner  [ A  B1^T ]
00559         //                                                    [ 0  M    ]
00560 
00561         rAElemPrecond = rAElemPrecond + rAElem;
00562         for(unsigned i=NUM_NODES_PER_ELEMENT*DIM; i<STENCIL_SIZE; i++)
00563         {
00564             for(unsigned j=0; j<NUM_NODES_PER_ELEMENT*DIM; j++)
00565             {
00566                 rAElemPrecond(i,j) = 0.0;
00567             }
00568         }
00569     }
00570 }
00571 
00572 
00573 template<size_t DIM>
00574 void NonlinearElasticitySolver<DIM>::ComputeStressAndStressDerivative(AbstractIncompressibleMaterialLaw<DIM>* pMaterialLaw,
00575                                                                       c_matrix<double,DIM,DIM>& rC, 
00576                                                                       c_matrix<double,DIM,DIM>& rInvC, 
00577                                                                       double pressure, 
00578                                                                       unsigned elementIndex,
00579                                                                       unsigned currentQuadPointGlobalIndex,
00580                                                                       c_matrix<double,DIM,DIM>& rT,
00581                                                                       FourthOrderTensor<DIM,DIM,DIM,DIM>& rDTdE,
00582                                                                       bool computeDTdE)
00583 {
00584     // just call the method on the material law
00585     pMaterialLaw->ComputeStressAndStressDerivative(rC,rInvC,pressure,rT,rDTdE,computeDTdE);
00586 }
00587 
00588 
00589 template<size_t DIM>
00590 void NonlinearElasticitySolver<DIM>::AssembleOnBoundaryElement(
00591             BoundaryElement<DIM-1,DIM>& rBoundaryElement,
00592             c_matrix<double,BOUNDARY_STENCIL_SIZE,BOUNDARY_STENCIL_SIZE>& rAelem,
00593             c_vector<double,BOUNDARY_STENCIL_SIZE>& rBelem,
00594             c_vector<double,DIM>& rTraction,
00595             bool assembleResidual,
00596             bool assembleJacobian)
00597 {
00598     rAelem.clear();
00599     rBelem.clear();
00600 
00601     if (assembleJacobian && !assembleResidual)
00602     {
00603         // nothing to do
00604         return;
00605     }
00606 
00607     c_vector<double, DIM> weighted_direction;
00608     double jacobian_determinant;
00609     mpQuadMesh->GetWeightedDirectionForBoundaryElement(rBoundaryElement.GetIndex(), weighted_direction, jacobian_determinant);
00610 
00611     c_vector<double,NUM_NODES_PER_BOUNDARY_ELEMENT> phi;
00612 
00613     for (unsigned quad_index=0; quad_index<this->mpBoundaryQuadratureRule->GetNumQuadPoints(); quad_index++)
00614     {
00615         double wJ = jacobian_determinant * this->mpBoundaryQuadratureRule->GetWeight(quad_index);
00616 
00617         const ChastePoint<DIM-1>& quad_point = this->mpBoundaryQuadratureRule->rGetQuadPoint(quad_index);
00618 
00619         QuadraticBasisFunction<DIM-1>::ComputeBasisFunctions(quad_point, phi);
00620 
00621         // get the required traction, interpolating X (slightly inefficiently, as interpolating
00622         // using quad bases) if necessary.
00623         c_vector<double,DIM> traction = zero_vector<double>(DIM);
00624         if (this->mUsingTractionBoundaryConditionFunction)
00625         {
00626             c_vector<double,DIM> X = zero_vector<double>(DIM);
00627             for (unsigned node_index=0; node_index<NUM_NODES_PER_BOUNDARY_ELEMENT; node_index++)
00628             {
00629                 X += phi(node_index)*mpQuadMesh->GetNode( rBoundaryElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
00630             }
00631             traction = (*(this->mpTractionBoundaryConditionFunction))(X);
00632         }
00633         else
00634         {
00635             traction = rTraction;
00636         }
00637 
00638 
00639         for (unsigned index=0; index<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index++)
00640         {
00641             unsigned spatial_dim = index%DIM;
00642             unsigned node_index = (index-spatial_dim)/DIM;
00643 
00644             assert(node_index < NUM_NODES_PER_BOUNDARY_ELEMENT);
00645 
00646             rBelem(index) -=    traction(spatial_dim)
00647                               * phi(node_index)
00648                               * wJ;
00649         }
00650     }
00651 }
00652 
00653 template<size_t DIM>
00654 void NonlinearElasticitySolver<DIM>::FormInitialGuess()
00655 {
00656     this->mCurrentSolution.resize(this->mNumDofs, 0.0);
00657 
00658     for (unsigned i=0; i<mpQuadMesh->GetNumElements(); i++)
00659     {
00660         double zero_strain_pressure;
00661         if (this->mMaterialLaws.size()==1)
00662         {
00663             // homogeneous
00664             zero_strain_pressure = this->mMaterialLaws[0]->GetZeroStrainPressure();
00665         }
00666         else
00667         {
00668             // heterogeneous
00669             zero_strain_pressure = this->mMaterialLaws[i]->GetZeroStrainPressure();
00670         }
00671 
00672         // loop over vertices and set pressure solution to be zero-strain-pressure
00673         for (unsigned j=0; j<NUM_VERTICES_PER_ELEMENT; j++)
00674         {
00675             unsigned index = mpQuadMesh->GetElement(i)->GetNodeGlobalIndex(j);
00676             this->mCurrentSolution[ DIM*mpQuadMesh->GetNumNodes() + index ] = zero_strain_pressure;
00677         }
00678     }
00679 }
00680 
00681 template<size_t DIM>
00682 void NonlinearElasticitySolver<DIM>::Initialise(std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00683 {
00684     assert(mpQuadMesh);
00685 
00686     AllocateMatrixMemory();
00687 
00688     for (unsigned i=0; i<this->mFixedNodes.size(); i++)
00689     {
00690         assert(this->mFixedNodes[i] < mpQuadMesh->GetNumNodes());
00691     }
00692 
00693     this->mpQuadratureRule = new GaussianQuadratureRule<DIM>(3);
00694     this->mpBoundaryQuadratureRule = new GaussianQuadratureRule<DIM-1>(3);
00695 
00696     FormInitialGuess();
00697 
00698     // compute the displacements at each of the fixed nodes, given the
00699     // fixed nodes locations.
00700     if (pFixedNodeLocations == NULL)
00701     {
00702         this->mFixedNodeDisplacements.clear();
00703         for (unsigned i=0; i<this->mFixedNodes.size(); i++)
00704         {
00705             this->mFixedNodeDisplacements.push_back(zero_vector<double>(DIM));
00706         }
00707     }
00708     else
00709     {
00710         assert(pFixedNodeLocations->size()==this->mFixedNodes.size());
00711         for (unsigned i=0; i<this->mFixedNodes.size(); i++)
00712         {
00713             unsigned index = this->mFixedNodes[i];
00714             c_vector<double,DIM> displacement = (*pFixedNodeLocations)[i]-mpQuadMesh->GetNode(index)->rGetLocation();
00715             this->mFixedNodeDisplacements.push_back(displacement);
00716         }
00717     }
00718     assert(this->mFixedNodeDisplacements.size()==this->mFixedNodes.size());
00719 }
00720 
00721 template<size_t DIM>
00722 void NonlinearElasticitySolver<DIM>::AllocateMatrixMemory()
00723 {
00724     if(DIM==2)
00725     {
00726         // 2D: N elements around a point => 7N+3 non-zeros in that row? Assume N<=10 (structured mesh would have N_max=6) => 73.
00727         unsigned num_non_zeros = std::min(75u, this->mNumDofs);
00728 
00729         this->mpLinearSystem = new LinearSystem(this->mNumDofs, num_non_zeros);
00730         this->mpPreconditionMatrixLinearSystem = new LinearSystem(this->mNumDofs, num_non_zeros);
00731     }
00732     else
00733     {
00734         assert(DIM==3);
00735 
00736         // pass in 0 as the preallocation number, and the false says don't preallocate
00737         this->mpLinearSystem = new LinearSystem(this->mNumDofs, 0);
00738         this->mpPreconditionMatrixLinearSystem = new LinearSystem(this->mNumDofs, 0);
00739 
00740         // 3D: N elements around a point. nz < (3*10+6)N (lazy estimate). Better estimate is 23N+4?. Assume N<20 => 500ish
00741 
00742         // in 3d we get the number of containing elements for each node and use that to obtain an upper bound
00743         // for the number of non-zeros for each DOF associated with that node.
00744 
00745         int* num_non_zeros_each_row = new int[this->mNumDofs];
00746         for(unsigned i=0; i<this->mNumDofs; i++)
00747         {
00748             num_non_zeros_each_row[i] = 0;
00749         }
00750 
00751         for(unsigned i=0; i<mpQuadMesh->GetNumNodes(); i++)
00752         {
00753             // this upper bound neglects the fact that two containing elements will share the same nodes..
00754             // 4 = max num dofs associated with this node
00755             // 30 = 3*9+3 = 3 dimensions x 9 other nodes on this element   +  3 vertices with a pressure unknown
00756             unsigned num_non_zeros_upper_bound = 4 + 30*mpQuadMesh->GetNode(i)->GetNumContainingElements();
00757             
00758             num_non_zeros_upper_bound = std::min(num_non_zeros_upper_bound, this->mNumDofs);
00759 
00760             num_non_zeros_each_row[DIM*i + 0] = num_non_zeros_upper_bound;
00761             num_non_zeros_each_row[DIM*i + 1] = num_non_zeros_upper_bound;
00762             num_non_zeros_each_row[DIM*i + 2] = num_non_zeros_upper_bound;
00763 
00764             //Could do !mpQuadMesh->GetNode(i)->IsInternal()
00765             if(i<mpQuadMesh->GetNumVertices()) // then this is a vertex
00766             {
00767                 num_non_zeros_each_row[DIM*mpQuadMesh->GetNumNodes() + i] = num_non_zeros_upper_bound;
00768             }
00769         }
00770         
00771         // NOTE: PetscTools::SetupMat() creates a MATAIJ matrix, which means the matrix will
00772         // be of type MATSEQAIJ if num_procs=1 and MATMPIAIJ otherwise. In the former case
00773         // MatSeqAIJSetPreallocation MUST be called [MatMPIAIJSetPreallocation will have 
00774         // no effect (silently)], and vice versa in the latter case
00775         if(PetscTools::GetNumProcs()==1)
00776         {
00777             MatSeqAIJSetPreallocation(this->mpLinearSystem->rGetLhsMatrix(),                   PETSC_NULL, num_non_zeros_each_row);
00778             MatSeqAIJSetPreallocation(this->mpPreconditionMatrixLinearSystem->rGetLhsMatrix(), PETSC_NULL, num_non_zeros_each_row);
00779         }
00780         else
00781         {
00782             PetscInt lo, hi;
00783             this->mpLinearSystem->GetOwnershipRange(lo, hi);
00784             int* num_non_zeros_each_row_this_proc = new int[hi-lo];
00785             int* zero = new int[hi-lo];
00786             for(unsigned i=0; i<unsigned(hi-lo); i++)
00787             {
00788                 num_non_zeros_each_row_this_proc[i] = num_non_zeros_each_row[lo+i];
00789                 zero[i] = 0;
00790             }
00791 
00792             MatMPIAIJSetPreallocation(this->mpLinearSystem->rGetLhsMatrix(), PETSC_NULL, num_non_zeros_each_row_this_proc, PETSC_NULL, num_non_zeros_each_row_this_proc);
00793             MatMPIAIJSetPreallocation(this->mpPreconditionMatrixLinearSystem->rGetLhsMatrix(), PETSC_NULL, num_non_zeros_each_row_this_proc, PETSC_NULL, num_non_zeros_each_row_this_proc);
00794         }
00795 
00796         //unsigned total_non_zeros = 0;
00797         //for(unsigned i=0; i<this->mNumDofs; i++)
00798         //{
00799         //   total_non_zeros += num_non_zeros_each_row[i];
00800         //}
00801         //std::cout << total_non_zeros << " versus " << 500*this->mNumDofs << "\n" << std::flush;
00802 
00803         delete [] num_non_zeros_each_row;
00804     }
00805 }
00806 
00807 
00808 
00809 template<size_t DIM>
00810 NonlinearElasticitySolver<DIM>::NonlinearElasticitySolver(
00811             QuadraticMesh<DIM>* pQuadMesh,
00812             AbstractIncompressibleMaterialLaw<DIM>* pMaterialLaw,
00813             c_vector<double,DIM> bodyForce,
00814             double density,
00815             std::string outputDirectory,
00816             std::vector<unsigned>& fixedNodes,
00817             std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00818     : AbstractNonlinearElasticitySolver<DIM>(DIM*pQuadMesh->GetNumNodes()+pQuadMesh->GetNumVertices(),
00819                                              pMaterialLaw, bodyForce, density,
00820                                              outputDirectory, fixedNodes),
00821       mpQuadMesh(pQuadMesh)
00822 {
00823     Initialise(pFixedNodeLocations);
00824 }
00825 
00826 
00827 template<size_t DIM>
00828 NonlinearElasticitySolver<DIM>::NonlinearElasticitySolver(
00829             QuadraticMesh<DIM>* pQuadMesh,
00830             std::vector<AbstractIncompressibleMaterialLaw<DIM>*>& rMaterialLaws,
00831             c_vector<double,DIM> bodyForce,
00832             double density,
00833             std::string outputDirectory,
00834             std::vector<unsigned>& fixedNodes,
00835             std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00836     : AbstractNonlinearElasticitySolver<DIM>(DIM*pQuadMesh->GetNumNodes()+pQuadMesh->GetNumVertices(),
00837                                              rMaterialLaws, bodyForce, density,
00838                                              outputDirectory, fixedNodes),
00839       mpQuadMesh(pQuadMesh)
00840 {
00841     assert(rMaterialLaws.size()==pQuadMesh->GetNumElements());
00842     Initialise(pFixedNodeLocations);
00843 }
00844 
00845 
00846 template<size_t DIM>
00847 NonlinearElasticitySolver<DIM>::~NonlinearElasticitySolver()
00848 {
00849 //    //Post-hoc debugging
00850 //    Mat matrix = this->mpLinearSystem->rGetLhsMatrix();
00851 //    for(unsigned i=0; i<mpQuadMesh->GetNumNodes(); i++)
00852 //    {
00853 //        unsigned elements = mpQuadMesh->GetNode(i)->GetNumContainingElements();
00854 //        if(i<mpQuadMesh->GetNumVertices()) // then this is a vertex
00855 //        {
00856 //            TRACE("VERT");
00857 //        }
00858 //        else
00859 //        {
00860 //            TRACE("INT");
00861 //        }
00862 //        PRINT_3_VARIABLES(DIM, i, elements);
00863 //        for (unsigned dim=0; dim<DIM; dim++)
00864 //        {
00865 //            int row=DIM*i+dim;
00866 //            int ncols;
00867 //            MatGetRow(matrix, row, &ncols, NULL, NULL);
00868 //            PRINT_3_VARIABLES(i,row,ncols);
00869 //        }
00870 // 
00871 //        if(i<mpQuadMesh->GetNumVertices()) // then this is a vertex
00872 //        {
00873 //            int row=DIM*mpQuadMesh->GetNumNodes() + i;
00874 //            int ncols;
00875 //            MatGetRow(matrix, row, &ncols, NULL, NULL);
00876 //            PRINT_3_VARIABLES(i,row,ncols);
00877 //        }
00878 //    }
00879 
00880     delete mpQuadratureRule;
00881     delete mpBoundaryQuadratureRule;
00882 }
00883 
00884 template<size_t DIM>
00885 void NonlinearElasticitySolver<DIM>::SetSurfaceTractionBoundaryConditions(
00886             std::vector<BoundaryElement<DIM-1,DIM>*>& rBoundaryElements,
00887             std::vector<c_vector<double,DIM> >& rSurfaceTractions)
00888 {
00889     assert(rBoundaryElements.size()==rSurfaceTractions.size());
00890     mBoundaryElements = rBoundaryElements;
00891     this->mSurfaceTractions = rSurfaceTractions;
00892 }
00893 
00894 template<size_t DIM>
00895 void NonlinearElasticitySolver<DIM>::SetFunctionalTractionBoundaryCondition(
00896             std::vector<BoundaryElement<DIM-1,DIM>*> rBoundaryElements,
00897             c_vector<double,DIM> (*pFunction)(c_vector<double,DIM>&))
00898 {
00899     mBoundaryElements = rBoundaryElements;
00900     this->mUsingTractionBoundaryConditionFunction = true;
00901     this->mpTractionBoundaryConditionFunction = pFunction;
00902 }
00903 
00904 template<size_t DIM>
00905 std::vector<double>& NonlinearElasticitySolver<DIM>::rGetPressures()
00906 {
00907     this->mPressures.clear();
00908     this->mPressures.resize(mpQuadMesh->GetNumVertices());
00909 
00910     for (unsigned i=0; i<mpQuadMesh->GetNumVertices(); i++)
00911     {
00912         this->mPressures[i] = this->mCurrentSolution[DIM*mpQuadMesh->GetNumNodes() + i];
00913     }
00914     return this->mPressures;
00915 }
00916 
00917 template<size_t DIM>
00918 std::vector<c_vector<double,DIM> >& NonlinearElasticitySolver<DIM>::rGetDeformedPosition()
00919 {
00920     this->mDeformedPosition.resize(mpQuadMesh->GetNumNodes(), zero_vector<double>(DIM));
00921     for (unsigned i=0; i<mpQuadMesh->GetNumNodes(); i++)
00922     {
00923         for (unsigned j=0; j<DIM; j++)
00924         {
00925             this->mDeformedPosition[i](j) = mpQuadMesh->GetNode(i)->rGetLocation()[j] + this->mCurrentSolution[DIM*i+j];
00926         }
00927     }
00928     return this->mDeformedPosition;
00929 }
00930 
00931 
00933 // Explicit instantiation
00935 
00936 //template class NonlinearElasticitySolver<1>;
00937 template class NonlinearElasticitySolver<2>;
00938 template class NonlinearElasticitySolver<3>;

Generated on Mon Nov 1 12:35:24 2010 for Chaste by  doxygen 1.5.5