NonlinearElasticitySolver.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2011
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 
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         PetscVecTools::Zero(this->mResidualVector);
00059     }
00060     if (assembleJacobian)
00061     {
00062         PetscMatTools::Zero(this->mJacobianMatrix);
00063         PetscMatTools::Zero(this->mPreconditionMatrix);
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 = this->mpQuadMesh->GetElementIteratorBegin();
00076          iter != this->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 " << this->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                 // At the moment we assume the vertices are the first num_vertices nodes in the list of nodes
00116                 // in the mesh. Hence:
00117                 unsigned vertex_index = element.GetNodeGlobalIndex(i);
00118                 assert(vertex_index < this->mpQuadMesh->GetNumVertices());
00119                 
00120                 // In the future (or currently with AlexW's adaptive quadratic mesh class (project work)), we
00121                 // will want to use this instead:
00122                 //unsigned vertex_index = this->mpQuadMesh->GetVertexIndexOfNode(element.GetNodeGlobalIndex(i));
00123 
00124                 p_indices[DIM*NUM_NODES_PER_ELEMENT + i] = DIM*this->mpQuadMesh->GetNumNodes() + vertex_index;
00125             }
00126 
00127             if (assembleJacobian)
00128             {
00129                 PetscMatTools::AddMultipleValues<STENCIL_SIZE>(this->mJacobianMatrix, p_indices, a_elem);
00130                 PetscMatTools::AddMultipleValues<STENCIL_SIZE>(this->mPreconditionMatrix, p_indices, a_elem_precond);
00131             }
00132 
00133             if (assembleResidual)
00134             {
00135                 PetscVecTools::AddMultipleValues<STENCIL_SIZE>(this->mResidualVector, p_indices, b_elem);
00136             }
00137         }
00138     }
00139 
00141     // loop over specified boundary elements and compute
00142     // surface traction terms
00144     c_vector<double, BOUNDARY_STENCIL_SIZE> b_boundary_elem;
00145     c_matrix<double, BOUNDARY_STENCIL_SIZE, BOUNDARY_STENCIL_SIZE> a_boundary_elem;
00146     if (this->mBoundaryElements.size()>0)
00147     {
00148         for (unsigned i=0; i<this->mBoundaryElements.size(); i++)
00149         {
00150             BoundaryElement<DIM-1,DIM>& r_boundary_element = *(this->mBoundaryElements[i]);
00151             AssembleOnBoundaryElement(r_boundary_element, a_boundary_elem, b_boundary_elem, this->mSurfaceTractions[i], assembleResidual, assembleJacobian);
00152 
00153             unsigned p_indices[BOUNDARY_STENCIL_SIZE];
00154             for (unsigned i=0; i<NUM_NODES_PER_BOUNDARY_ELEMENT; i++)
00155             {
00156                 for (unsigned j=0; j<DIM; j++)
00157                 {
00158                     p_indices[DIM*i+j] = DIM*r_boundary_element.GetNodeGlobalIndex(i) + j;
00159                 }
00160             }
00161 
00162             for (unsigned i=0; i<DIM /*vertices per boundary elem */; i++)
00163             {
00164                 p_indices[DIM*NUM_NODES_PER_BOUNDARY_ELEMENT + i] = DIM*this->mpQuadMesh->GetNumNodes() + r_boundary_element.GetNodeGlobalIndex(i);
00165             }
00166 
00167             if (assembleJacobian)
00168             {
00169                 PetscMatTools::AddMultipleValues<BOUNDARY_STENCIL_SIZE>(this->mJacobianMatrix, p_indices, a_boundary_elem);
00170                 PetscMatTools::AddMultipleValues<BOUNDARY_STENCIL_SIZE>(this->mPreconditionMatrix, p_indices, a_boundary_elem);
00171             }
00172 
00173             if (assembleResidual)
00174             {
00175                 PetscVecTools::AddMultipleValues<BOUNDARY_STENCIL_SIZE>(this->mResidualVector, p_indices, b_boundary_elem);
00176             }
00177 
00178             // some extra checking
00179             if (DIM==2)
00180             {
00181                 assert(8==BOUNDARY_STENCIL_SIZE);
00182                 assert(b_boundary_elem(6)==0);
00183                 assert(b_boundary_elem(7)==0);
00184             }
00185         }
00186     }
00187 
00188     this->FinishAssembleSystem(assembleResidual, assembleJacobian);
00189 }
00190 
00191 template<size_t DIM>
00192 void NonlinearElasticitySolver<DIM>::AssembleOnElement(
00193             Element<DIM, DIM>& rElement,
00194             c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElem,
00195             c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElemPrecond,
00196             c_vector<double, STENCIL_SIZE>& rBElem,
00197             bool assembleResidual,
00198             bool assembleJacobian)
00199 {
00200     static c_matrix<double,DIM,DIM> jacobian;
00201     static c_matrix<double,DIM,DIM> inverse_jacobian;
00202     double jacobian_determinant;
00203     
00204     this->mpQuadMesh->GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
00205 
00206     if (assembleJacobian)
00207     {
00208         rAElem.clear();
00209         rAElemPrecond.clear();
00210     }
00211 
00212     if (assembleResidual)
00213     {
00214         rBElem.clear();
00215     }
00216 
00218     // Get the current displacement at the nodes
00220     static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> element_current_displacements;
00221     static c_vector<double,NUM_VERTICES_PER_ELEMENT> element_current_pressures;
00222     for (unsigned II=0; II<NUM_NODES_PER_ELEMENT; II++)
00223     {
00224         for (unsigned JJ=0; JJ<DIM; JJ++)
00225         {
00226             element_current_displacements(JJ,II) = this->mCurrentSolution[DIM*rElement.GetNodeGlobalIndex(II) + JJ];
00227         }
00228     }
00229 
00231     // Get the current pressure at the vertices
00233     for (unsigned II=0; II<NUM_VERTICES_PER_ELEMENT; II++)
00234     {
00235 
00236         // At the moment we assume the vertices are the first num_vertices nodes in the list of nodes
00237         // in the mesh. Hence:
00238         unsigned vertex_index = rElement.GetNodeGlobalIndex(II);
00239         assert(vertex_index < this->mpQuadMesh->GetNumVertices());
00240 
00241         // In the future (or currently with AlexW's adaptive quadratic mesh class (project work)), we
00242         // will want to use this instead:
00243         //unsigned vertex_index = this->mpQuadMesh->GetVertexIndexOfNode( rElement.GetNodeGlobalIndex(II) );
00244 
00245         element_current_pressures(II) = this->mCurrentSolution[DIM*this->mpQuadMesh->GetNumNodes() + vertex_index];
00246     }
00247 
00248     // allocate memory for the basis functions values and derivative values
00249     static c_vector<double, NUM_VERTICES_PER_ELEMENT> linear_phi;
00250     static c_vector<double, NUM_NODES_PER_ELEMENT> quad_phi;
00251     static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> grad_quad_phi;
00252     static c_matrix<double, NUM_NODES_PER_ELEMENT, DIM> trans_grad_quad_phi;
00253 
00254 
00255     // get the material law
00256     AbstractIncompressibleMaterialLaw<DIM>* p_material_law;
00257     if (this->mMaterialLaws.size()==1)
00258     {
00259         // homogeneous
00260         p_material_law = this->mMaterialLaws[0];
00261     }
00262     else
00263     {
00264         // heterogeneous
00265         #define COVERAGE_IGNORE // not going to have tests in cts for everything
00266         p_material_law = this->mMaterialLaws[rElement.GetIndex()];
00267         #undef COVERAGE_IGNORE
00268     }
00269     
00270     
00271     static c_matrix<double,DIM,DIM> grad_u; // grad_u = (du_i/dX_M)
00272             
00273     static c_matrix<double,DIM,DIM> F;      // the deformation gradient, F = dx/dX, F_{iM} = dx_i/dX_M
00274     static c_matrix<double,DIM,DIM> C;      // Green deformation tensor, C = F^T F
00275     static c_matrix<double,DIM,DIM> inv_C;  // inverse(C)
00276     static c_matrix<double,DIM,DIM> inv_F;  // inverse(F)
00277     static c_matrix<double,DIM,DIM> T;      // Second Piola-Kirchoff stress tensor (= dW/dE = 2dW/dC) 
00278 
00279     static c_matrix<double,DIM,DIM> F_T;    // F*T
00280     static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> F_T_grad_quad_phi; // F*T*grad_quad_phi
00281 
00282     c_vector<double,DIM> body_force;
00283 
00284     static FourthOrderTensor<DIM,DIM,DIM,DIM> dTdE;    // dTdE(M,N,P,Q) = dT_{MN}/dE_{PQ} 
00285     static FourthOrderTensor<DIM,DIM,DIM,DIM> dSdF;    // dSdF(M,i,N,j) = dS_{Mi}/dF_{jN} 
00286 
00287     static FourthOrderTensor<NUM_NODES_PER_ELEMENT,DIM,DIM,DIM> temp_tensor;
00288     static FourthOrderTensor<NUM_NODES_PER_ELEMENT,DIM,NUM_NODES_PER_ELEMENT,DIM> dSdF_quad_quad;
00289 
00290     static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> temp_matrix;
00291     static c_matrix<double,NUM_NODES_PER_ELEMENT,DIM> grad_quad_phi_times_invF;
00292 
00293 
00294 
00300     for (unsigned quadrature_index=0; quadrature_index < this->mpQuadratureRule->GetNumQuadPoints(); quadrature_index++)
00301     {
00302         // this is needed by the cardiac mechanics solver
00303         unsigned current_quad_point_global_index =   rElement.GetIndex()*this->mpQuadratureRule->GetNumQuadPoints()
00304                                                    + quadrature_index;
00305 
00306         double wJ = jacobian_determinant * this->mpQuadratureRule->GetWeight(quadrature_index);
00307 
00308         const ChastePoint<DIM>& quadrature_point = this->mpQuadratureRule->rGetQuadPoint(quadrature_index);
00309 
00311         // set up basis function info
00313         LinearBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, linear_phi);
00314         QuadraticBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, quad_phi);
00315         QuadraticBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(quadrature_point, inverse_jacobian, grad_quad_phi);
00316         trans_grad_quad_phi = trans(grad_quad_phi);
00317         
00318         
00320         // get the body force, interpolating X if necessary
00322 
00323         if(assembleResidual)
00324         {
00325             if (this->mUsingBodyForceFunction)
00326             {
00327                 c_vector<double,DIM> X = zero_vector<double>(DIM);
00328                 // interpolate X (using the vertices and the /linear/ bases, as no curvilinear elements
00329                 for (unsigned node_index=0; node_index<NUM_VERTICES_PER_ELEMENT; node_index++)
00330                 {
00331                     X += linear_phi(node_index)*this->mpQuadMesh->GetNode( rElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
00332                 }
00333                 body_force = (*(this->mpBodyForceFunction))(X, this->mCurrentTime);
00334             }
00335             else
00336             {
00337                 body_force = this->mBodyForce;
00338             }
00339         }
00340 
00342         // interpolate grad_u and p
00344         grad_u = zero_matrix<double>(DIM,DIM); 
00345 
00346         for (unsigned node_index=0; node_index<NUM_NODES_PER_ELEMENT; node_index++)
00347         {
00348             for (unsigned i=0; i<DIM; i++)
00349             {
00350                 for (unsigned M=0; M<DIM; M++)
00351                 {
00352                     grad_u(i,M) += grad_quad_phi(M,node_index)*element_current_displacements(i,node_index);
00353                 }
00354             }
00355         }
00356 
00357         double pressure = 0;
00358         for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00359         {
00360             pressure += linear_phi(vertex_index)*element_current_pressures(vertex_index);
00361         }
00362 
00363 
00365         // calculate C, inv(C) and T
00367         for (unsigned i=0; i<DIM; i++)
00368         {
00369             for (unsigned M=0; M<DIM; M++)
00370             {
00371                 F(i,M) = (i==M?1:0) + grad_u(i,M);
00372             }
00373         }
00374 
00375         C = prod(trans(F),F);
00376         inv_C = Inverse(C);
00377         inv_F = Inverse(F);
00378 
00379         double detF = Determinant(F);
00380 
00381         this->ComputeStressAndStressDerivative(p_material_law, C, inv_C, pressure, rElement.GetIndex(), current_quad_point_global_index,
00382                                                T, dTdE, assembleJacobian);
00383 
00384 
00386         // residual vector
00388         if (assembleResidual)
00389         {
00390             F_T = prod(F,T);
00391             F_T_grad_quad_phi = prod(F_T, grad_quad_phi);
00392             
00393             for (unsigned index=0; index<NUM_NODES_PER_ELEMENT*DIM; index++)
00394             {
00395                 unsigned spatial_dim = index%DIM;
00396                 unsigned node_index = (index-spatial_dim)/DIM;
00397 
00398                 rBElem(index) +=  - this->mDensity
00399                                   * body_force(spatial_dim)
00400                                   * quad_phi(node_index)
00401                                   * wJ;
00402                                   
00403                 // the  T(M,N)*F(spatial_dim,M)*grad_quad_phi(N,node_index)  term                  
00404                 rBElem(index) +=   F_T_grad_quad_phi(spatial_dim,node_index)
00405                                  * wJ;
00406             }
00407 
00408             for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00409             {
00410                 rBElem( NUM_NODES_PER_ELEMENT*DIM + vertex_index ) +=   linear_phi(vertex_index)
00411                                                                       * (detF - 1)
00412                                                                       * wJ;
00413             }
00414         }
00415         
00416         
00418         // Jacobian matrix
00420         if (assembleJacobian) 
00421         {
00422             // save trans(grad_quad_phi) * invF
00423             grad_quad_phi_times_invF = prod(trans_grad_quad_phi, inv_F);
00424             
00426             // Set up the tensor dSdF
00427             //
00428             // dSdF as a function of T and dTdE (which is what the material law returns) is given by:
00429             //
00430             // dS_{Mi}/dF_{jN} = (dT_{MN}/dC_{PQ}+dT_{MN}/dC_{PQ}) F{iP} F_{jQ}  + T_{MN} delta_{ij}
00431             //
00432             // todo1: this should probably move into the material law (but need to make sure 
00433             // memory is handled efficiently
00434             // todo2: get material law to return this immediately, not dTdE
00436 
00437             // set up the tensor 0.5(dTdE(M,N,P,Q) + dTdE(M,N,Q,P))
00438             for (unsigned M=0; M<DIM; M++)
00439             {
00440                 for (unsigned N=0; N<DIM; N++)
00441                 {
00442                     for (unsigned P=0; P<DIM; P++)
00443                     {
00444                         for (unsigned Q=0; Q<DIM; Q++)
00445                         {
00446                             // this is NOT dSdF, just using this as storage space
00447                             dSdF(M,N,P,Q) = 0.5*(dTdE(M,N,P,Q) + dTdE(M,N,Q,P));
00448                         }
00449                     }
00450                 }
00451             }
00452             // This is NOT dTdE, just reusing memory. A^{MdPQ}  = F^d_N * dTdE_sym^{MNPQ}
00453             dTdE.template SetAsContractionOnSecondDimension<DIM>(F, dSdF);  
00454             // dSdF{MdPe} := F^d_N * F^e_Q * dTdE_sym^{MNPQ}
00455             dSdF.template SetAsContractionOnFourthDimension<DIM>(F, dTdE);  
00456             
00457             // now add the T_{MN} delta_{ij} term
00458             for(unsigned M=0; M<DIM; M++)
00459             {
00460                 for(unsigned N=0; N<DIM; N++)
00461                 {
00462                     for(unsigned i=0; i<DIM; i++)
00463                     {
00464                         dSdF(M,i,N,i) += T(M,N);
00465                     }
00466                 }
00467             }
00468             
00469 
00471             // Set up the tensor  
00472             //   dSdF_quad_quad(node_index1, spatial_dim1, node_index2, spatial_dim2)
00473             //            =    dS_{M,spatial_dim1}/d_F{spatial_dim2,N}      
00474             //               * grad_quad_phi(M,node_index1) 
00475             //               * grad_quad_phi(P,node_index2)
00476             //
00477             //            =    dSdF(M,spatial_index1,N,spatial_index2)
00478             //               * grad_quad_phi(M,node_index1) 
00479             //               * grad_quad_phi(P,node_index2)
00480             //            
00482             temp_tensor.template SetAsContractionOnFirstDimension<DIM>( trans_grad_quad_phi, dSdF );
00483             dSdF_quad_quad.template SetAsContractionOnThirdDimension<DIM>( trans_grad_quad_phi, temp_tensor );
00484 
00485 
00486 
00487 
00488             for (unsigned index1=0; index1<NUM_NODES_PER_ELEMENT*DIM; index1++)
00489             {
00490                 unsigned spatial_dim1 = index1%DIM;
00491                 unsigned node_index1 = (index1-spatial_dim1)/DIM;
00492 
00493 
00494                 for (unsigned index2=0; index2<NUM_NODES_PER_ELEMENT*DIM; index2++)
00495                 {
00496                     unsigned spatial_dim2 = index2%DIM;
00497                     unsigned node_index2 = (index2-spatial_dim2)/DIM;
00498 
00499                     // the dSdF*grad_quad_phi*grad_quad_phi term
00500                     rAElem(index1,index2)  +=   dSdF_quad_quad(node_index1,spatial_dim1,node_index2,spatial_dim2)
00501                                               * wJ;
00502                 }
00503 
00504                 for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00505                 {
00506                     unsigned index2 = NUM_NODES_PER_ELEMENT*DIM + vertex_index;
00507         
00508                     // the -invF(M,spatial_dim1)*grad_quad_phi(M,node_index1)*linear_phi(vertex_index) term
00509                     rAElem(index1,index2)  +=  - grad_quad_phi_times_invF(node_index1,spatial_dim1)
00510                                                * linear_phi(vertex_index)
00511                                                * wJ;
00512                 }
00513             }
00514 
00515             for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00516             {
00517                 unsigned index1 = NUM_NODES_PER_ELEMENT*DIM + vertex_index;
00518 
00519                 for (unsigned index2=0; index2<NUM_NODES_PER_ELEMENT*DIM; index2++)
00520                 {
00521                     unsigned spatial_dim2 = index2%DIM;
00522                     unsigned node_index2 = (index2-spatial_dim2)/DIM;
00523 
00524                     // same as (negative of) the opposite block (ie a few lines up), except for detF
00525                     rAElem(index1,index2) +=   detF
00526                                              * grad_quad_phi_times_invF(node_index2,spatial_dim2)
00527                                              * linear_phi(vertex_index)
00528                                              * wJ;
00529                 }
00530 
00532                 // Preconditioner matrix
00533                 // Fill the mass matrix (ie \intgl phi_i phi_j) in the
00534                 // pressure-pressure block. Note, the rest of the
00535                 // entries are filled in at the end
00537                 for (unsigned vertex_index2=0; vertex_index2<NUM_VERTICES_PER_ELEMENT; vertex_index2++)
00538                 {
00539                     unsigned index2 = NUM_NODES_PER_ELEMENT*DIM + vertex_index2;
00540                     rAElemPrecond(index1,index2) +=   linear_phi(vertex_index)
00541                                                     * linear_phi(vertex_index2)
00542                                                     * wJ;
00543                 }
00544             }
00545         }
00546     }
00547 
00548 
00549     if (assembleJacobian)
00550     {
00551         // Fill in the other blocks of the preconditioner matrix, by adding 
00552         // the jacobian matrix (this doesn't effect the pressure-pressure block 
00553         // of rAElemPrecond as the pressure-pressure block of rAElem is zero),
00554         // and the zero a block.
00555         //
00556         // The following altogether gives the preconditioner  [ A  B1^T ]
00557         //                                                    [ 0  M    ]
00558 
00559         rAElemPrecond = rAElemPrecond + rAElem;
00560         for(unsigned i=NUM_NODES_PER_ELEMENT*DIM; i<STENCIL_SIZE; i++)
00561         {
00562             for(unsigned j=0; j<NUM_NODES_PER_ELEMENT*DIM; j++)
00563             {
00564                 rAElemPrecond(i,j) = 0.0;
00565             }
00566         }
00567     }
00568 }
00569 
00570 
00571 
00572 
00573 
00574 template<size_t DIM>
00575 void NonlinearElasticitySolver<DIM>::AssembleOnBoundaryElement(
00576             BoundaryElement<DIM-1,DIM>& rBoundaryElement,
00577             c_matrix<double,BOUNDARY_STENCIL_SIZE,BOUNDARY_STENCIL_SIZE>& rAelem,
00578             c_vector<double,BOUNDARY_STENCIL_SIZE>& rBelem,
00579             c_vector<double,DIM>& rTraction,
00580             bool assembleResidual,
00581             bool assembleJacobian)
00582 {
00583     rAelem.clear();
00584     rBelem.clear();
00585 
00586     if (assembleJacobian && !assembleResidual)
00587     {
00588         // nothing to do
00589         return;
00590     }
00591 
00592     c_vector<double, DIM> weighted_direction;
00593     double jacobian_determinant;
00594     this->mpQuadMesh->GetWeightedDirectionForBoundaryElement(rBoundaryElement.GetIndex(), weighted_direction, jacobian_determinant);
00595 
00596     c_vector<double,NUM_NODES_PER_BOUNDARY_ELEMENT> phi;
00597 
00598     for (unsigned quad_index=0; quad_index<this->mpBoundaryQuadratureRule->GetNumQuadPoints(); quad_index++)
00599     {
00600         double wJ = jacobian_determinant * this->mpBoundaryQuadratureRule->GetWeight(quad_index);
00601 
00602         const ChastePoint<DIM-1>& quad_point = this->mpBoundaryQuadratureRule->rGetQuadPoint(quad_index);
00603 
00604         QuadraticBasisFunction<DIM-1>::ComputeBasisFunctions(quad_point, phi);
00605 
00606         // get the required traction, interpolating X (slightly inefficiently, as interpolating
00607         // using quad bases) if necessary.
00608         c_vector<double,DIM> traction = zero_vector<double>(DIM);
00609         if (this->mUsingTractionBoundaryConditionFunction)
00610         {
00611             c_vector<double,DIM> X = zero_vector<double>(DIM);
00612             for (unsigned node_index=0; node_index<NUM_NODES_PER_BOUNDARY_ELEMENT; node_index++)
00613             {
00614                 X += phi(node_index)*this->mpQuadMesh->GetNode( rBoundaryElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
00615             }
00616             traction = (*(this->mpTractionBoundaryConditionFunction))(X, this->mCurrentTime);
00617         }
00618         else
00619         {
00620             traction = rTraction;
00621         }
00622 
00623 
00624         for (unsigned index=0; index<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index++)
00625         {
00626             unsigned spatial_dim = index%DIM;
00627             unsigned node_index = (index-spatial_dim)/DIM;
00628 
00629             assert(node_index < NUM_NODES_PER_BOUNDARY_ELEMENT);
00630 
00631             rBelem(index) -=    traction(spatial_dim)
00632                               * phi(node_index)
00633                               * wJ;
00634         }
00635     }
00636 }
00637 
00638 template<size_t DIM>
00639 void NonlinearElasticitySolver<DIM>::FormInitialGuess()
00640 {
00641     this->mCurrentSolution.resize(this->mNumDofs, 0.0);
00642 
00643     for (unsigned i=0; i<this->mpQuadMesh->GetNumElements(); i++)
00644     {
00645         double zero_strain_pressure;
00646         if (this->mMaterialLaws.size()==1)
00647         {
00648             // homogeneous
00649             zero_strain_pressure = this->mMaterialLaws[0]->GetZeroStrainPressure();
00650         }
00651         else
00652         {
00653             // heterogeneous
00654             zero_strain_pressure = this->mMaterialLaws[i]->GetZeroStrainPressure();
00655         }
00656 
00657         // loop over vertices and set pressure solution to be zero-strain-pressure
00658         for (unsigned j=0; j<NUM_VERTICES_PER_ELEMENT; j++)
00659         {
00660             // At the moment we assume the vertices are the first num_vertices nodes in the list of nodes
00661             // in the mesh. Hence:
00662             unsigned vertex_index = this->mpQuadMesh->GetElement(i)->GetNodeGlobalIndex(j);
00663             assert(vertex_index < this->mpQuadMesh->GetNumVertices());
00664 
00665             // In the future (or currently with AlexW's adaptive quadratic mesh class (project work)), we
00666             // will want to use this instead:
00667             //unsigned vertex_index = this->mpQuadMesh->GetVertexIndexOfNode(this->mpQuadMesh->GetElement(i)->GetNodeGlobalIndex(j));
00668             
00669             this->mCurrentSolution[ DIM*this->mpQuadMesh->GetNumNodes() + vertex_index ] = zero_strain_pressure;
00670         }
00671     }
00672 }
00673 
00674 
00675 template<size_t DIM>
00676 NonlinearElasticitySolver<DIM>::NonlinearElasticitySolver(
00677             QuadraticMesh<DIM>* pQuadMesh,
00678             AbstractMaterialLaw<DIM>* pMaterialLaw,
00679             c_vector<double,DIM> bodyForce,
00680             double density,
00681             std::string outputDirectory,
00682             std::vector<unsigned>& fixedNodes,
00683             std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00684     : AbstractNonlinearElasticitySolver<DIM>(pQuadMesh,
00685                                              bodyForce, density,
00686                                              outputDirectory, fixedNodes,
00687                                              INCOMPRESSIBLE)
00688 {
00689     assert(pMaterialLaw != NULL);
00690 
00691     AbstractIncompressibleMaterialLaw<DIM>* p_law = dynamic_cast<AbstractIncompressibleMaterialLaw<DIM>*>(pMaterialLaw);
00692     if(!p_law)
00693     {
00694         EXCEPTION("NonlinearElasticitySolver must take in an incompressible material law (ie of type AbstractIncompressibleMaterialLaw)");
00695     }
00696     mMaterialLaws.push_back(p_law);
00697 
00698     Initialise(pFixedNodeLocations);
00699     FormInitialGuess();
00700 }
00701 
00702 
00703 template<size_t DIM>
00704 NonlinearElasticitySolver<DIM>::NonlinearElasticitySolver(
00705             QuadraticMesh<DIM>* pQuadMesh,
00706             std::vector<AbstractMaterialLaw<DIM>*>& rMaterialLaws,
00707             c_vector<double,DIM> bodyForce,
00708             double density,
00709             std::string outputDirectory,
00710             std::vector<unsigned>& fixedNodes,
00711             std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00712     : AbstractNonlinearElasticitySolver<DIM>(pQuadMesh,
00713                                              bodyForce, density,
00714                                              outputDirectory, fixedNodes,
00715                                              INCOMPRESSIBLE)
00716 {
00717     mMaterialLaws.resize(rMaterialLaws.size(), NULL);
00718     for (unsigned i=0; i<mMaterialLaws.size(); i++)
00719     {
00720         assert(rMaterialLaws[i] != NULL);
00721         AbstractIncompressibleMaterialLaw<DIM>* p_law = dynamic_cast<AbstractIncompressibleMaterialLaw<DIM>*>(rMaterialLaws[i]);
00722         if(!p_law)
00723         {
00724             EXCEPTION("NonlinearElasticitySolver must take in an incompressible material law (ie of type AbstractIncompressibleMaterialLaw)");
00725         }
00726         mMaterialLaws[i] = p_law;
00727     }
00728 
00729     assert(rMaterialLaws.size()==pQuadMesh->GetNumElements());
00730     Initialise(pFixedNodeLocations);
00731     FormInitialGuess();
00732 }
00733 
00734 
00735 template<size_t DIM>
00736 NonlinearElasticitySolver<DIM>::~NonlinearElasticitySolver()
00737 {
00738 }
00739 
00740 
00741 
00742 template<size_t DIM>
00743 std::vector<double>& NonlinearElasticitySolver<DIM>::rGetPressures()
00744 {
00745     mPressures.clear();
00746     mPressures.resize(this->mpQuadMesh->GetNumVertices());
00747 
00748     for (unsigned i=0; i<this->mpQuadMesh->GetNumVertices(); i++)
00749     {
00750         mPressures[i] = this->mCurrentSolution[DIM*this->mpQuadMesh->GetNumNodes() + i];
00751     }
00752     return mPressures;
00753 }
00754 
00755 
00756 
00758 // Explicit instantiation
00760 
00761 //template class NonlinearElasticitySolver<1>;
00762 template class NonlinearElasticitySolver<2>;
00763 template class NonlinearElasticitySolver<3>;

Generated on Tue May 31 14:31:49 2011 for Chaste by  doxygen 1.5.5