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         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 = 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                 this->mpLinearSystem->AddLhsMultipleValues(p_indices, a_elem);
00130                 this->mpPreconditionMatrixLinearSystem->AddLhsMultipleValues(p_indices, a_elem_precond);
00131             }
00132 
00133             if (assembleResidual)
00134             {
00135                 this->mpLinearSystem->AddRhsMultipleValues(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                 this->mpLinearSystem->AddLhsMultipleValues(p_indices, a_boundary_elem);
00170                 this->mpPreconditionMatrixLinearSystem->AddLhsMultipleValues(p_indices, a_boundary_elem);
00171             }
00172 
00173             if (assembleResidual)
00174             {
00175                 this->mpLinearSystem->AddRhsMultipleValues(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     if (assembleResidual)
00189     {
00190         this->mpLinearSystem->AssembleRhsVector();
00191     }
00192     if (assembleJacobian)
00193     {
00194         this->mpLinearSystem->AssembleIntermediateLhsMatrix();
00195         this->mpPreconditionMatrixLinearSystem->AssembleIntermediateLhsMatrix();
00196     }
00197 
00198     // Apply Dirichlet boundary conditions
00199     this->ApplyBoundaryConditions(assembleJacobian);
00200 
00201     if (assembleResidual)
00202     {
00203         this->mpLinearSystem->AssembleRhsVector();
00204     }
00205     if (assembleJacobian)
00206     {
00207         this->mpLinearSystem->AssembleFinalLhsMatrix();
00208         this->mpPreconditionMatrixLinearSystem->AssembleFinalLhsMatrix();
00209     }
00210 }
00211 
00212 template<size_t DIM>
00213 void NonlinearElasticitySolver<DIM>::AssembleOnElement(
00214             Element<DIM, DIM>& rElement,
00215             c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElem,
00216             c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElemPrecond,
00217             c_vector<double, STENCIL_SIZE>& rBElem,
00218             bool assembleResidual,
00219             bool assembleJacobian)
00220 {
00221     static c_matrix<double,DIM,DIM> jacobian;
00222     static c_matrix<double,DIM,DIM> inverse_jacobian;
00223     double jacobian_determinant;
00224     
00225     this->mpQuadMesh->GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
00226 
00227     if (assembleJacobian)
00228     {
00229         rAElem.clear();
00230         rAElemPrecond.clear();
00231     }
00232 
00233     if (assembleResidual)
00234     {
00235         rBElem.clear();
00236     }
00237 
00239     // Get the current displacement at the nodes
00241     static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> element_current_displacements;
00242     static c_vector<double,NUM_VERTICES_PER_ELEMENT> element_current_pressures;
00243     for (unsigned II=0; II<NUM_NODES_PER_ELEMENT; II++)
00244     {
00245         for (unsigned JJ=0; JJ<DIM; JJ++)
00246         {
00247             element_current_displacements(JJ,II) = this->mCurrentSolution[DIM*rElement.GetNodeGlobalIndex(II) + JJ];
00248         }
00249     }
00250 
00252     // Get the current pressure at the vertices
00254     for (unsigned II=0; II<NUM_VERTICES_PER_ELEMENT; II++)
00255     {
00256 
00257         // At the moment we assume the vertices are the first num_vertices nodes in the list of nodes
00258         // in the mesh. Hence:
00259         unsigned vertex_index = rElement.GetNodeGlobalIndex(II);
00260         assert(vertex_index < this->mpQuadMesh->GetNumVertices());
00261 
00262         // In the future (or currently with AlexW's adaptive quadratic mesh class (project work)), we
00263         // will want to use this instead:
00264         //unsigned vertex_index = this->mpQuadMesh->GetVertexIndexOfNode( rElement.GetNodeGlobalIndex(II) );
00265 
00266         element_current_pressures(II) = this->mCurrentSolution[DIM*this->mpQuadMesh->GetNumNodes() + vertex_index];
00267     }
00268 
00269     // allocate memory for the basis functions values and derivative values
00270     static c_vector<double, NUM_VERTICES_PER_ELEMENT> linear_phi;
00271     static c_vector<double, NUM_NODES_PER_ELEMENT> quad_phi;
00272     static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> grad_quad_phi;
00273     static c_matrix<double, NUM_NODES_PER_ELEMENT, DIM> trans_grad_quad_phi;
00274 
00275 
00276     // get the material law
00277     AbstractIncompressibleMaterialLaw<DIM>* p_material_law;
00278     if (this->mMaterialLaws.size()==1)
00279     {
00280         // homogeneous
00281         p_material_law = this->mMaterialLaws[0];
00282     }
00283     else
00284     {
00285         // heterogeneous
00286         #define COVERAGE_IGNORE // not going to have tests in cts for everything
00287         p_material_law = this->mMaterialLaws[rElement.GetIndex()];
00288         #undef COVERAGE_IGNORE
00289     }
00290     
00291     
00292     static c_matrix<double,DIM,DIM> grad_u; // grad_u = (du_i/dX_M)
00293             
00294     static c_matrix<double,DIM,DIM> F;      // the deformation gradient, F = dx/dX, F_{iM} = dx_i/dX_M
00295     static c_matrix<double,DIM,DIM> C;      // Green deformation tensor, C = F^T F
00296     static c_matrix<double,DIM,DIM> inv_C;  // inverse(C)
00297     static c_matrix<double,DIM,DIM> inv_F;  // inverse(F)
00298     static c_matrix<double,DIM,DIM> T;      // Second Piola-Kirchoff stress tensor (= dW/dE = 2dW/dC) 
00299 
00300     static c_matrix<double,DIM,DIM> F_T;    // F*T
00301     static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> F_T_grad_quad_phi; // F*T*grad_quad_phi
00302 
00303     c_vector<double,DIM> body_force;
00304 
00305     static FourthOrderTensor<DIM,DIM,DIM,DIM> dTdE;    // dTdE(M,N,P,Q) = dT_{MN}/dE_{PQ} 
00306     static FourthOrderTensor<DIM,DIM,DIM,DIM> dSdF;    // dSdF(M,i,N,j) = dS_{Mi}/dF_{jN} 
00307 
00308     static FourthOrderTensor<NUM_NODES_PER_ELEMENT,DIM,DIM,DIM> temp_tensor;
00309     static FourthOrderTensor<NUM_NODES_PER_ELEMENT,DIM,NUM_NODES_PER_ELEMENT,DIM> dSdF_quad_quad;
00310 
00311     static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> temp_matrix;
00312     static c_matrix<double,NUM_NODES_PER_ELEMENT,DIM> grad_quad_phi_times_invF;
00313 
00314 
00315 
00321     for (unsigned quadrature_index=0; quadrature_index < this->mpQuadratureRule->GetNumQuadPoints(); quadrature_index++)
00322     {
00323         // this is needed by the cardiac mechanics solver
00324         unsigned current_quad_point_global_index =   rElement.GetIndex()*this->mpQuadratureRule->GetNumQuadPoints()
00325                                                    + quadrature_index;
00326 
00327         double wJ = jacobian_determinant * this->mpQuadratureRule->GetWeight(quadrature_index);
00328 
00329         const ChastePoint<DIM>& quadrature_point = this->mpQuadratureRule->rGetQuadPoint(quadrature_index);
00330 
00332         // set up basis function info
00334         LinearBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, linear_phi);
00335         QuadraticBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, quad_phi);
00336         QuadraticBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(quadrature_point, inverse_jacobian, grad_quad_phi);
00337         trans_grad_quad_phi = trans(grad_quad_phi);
00338         
00339         
00341         // get the body force, interpolating X if necessary
00343 
00344         if(assembleResidual)
00345         {
00346             if (this->mUsingBodyForceFunction)
00347             {
00348                 c_vector<double,DIM> X = zero_vector<double>(DIM);
00349                 // interpolate X (using the vertices and the /linear/ bases, as no curvilinear elements
00350                 for (unsigned node_index=0; node_index<NUM_VERTICES_PER_ELEMENT; node_index++)
00351                 {
00352                     X += linear_phi(node_index)*this->mpQuadMesh->GetNode( rElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
00353                 }
00354                 body_force = (*(this->mpBodyForceFunction))(X, this->mCurrentTime);
00355             }
00356             else
00357             {
00358                 body_force = this->mBodyForce;
00359             }
00360         }
00361 
00363         // interpolate grad_u and p
00365         grad_u = zero_matrix<double>(DIM,DIM); 
00366 
00367         for (unsigned node_index=0; node_index<NUM_NODES_PER_ELEMENT; node_index++)
00368         {
00369             for (unsigned i=0; i<DIM; i++)
00370             {
00371                 for (unsigned M=0; M<DIM; M++)
00372                 {
00373                     grad_u(i,M) += grad_quad_phi(M,node_index)*element_current_displacements(i,node_index);
00374                 }
00375             }
00376         }
00377 
00378         double pressure = 0;
00379         for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00380         {
00381             pressure += linear_phi(vertex_index)*element_current_pressures(vertex_index);
00382         }
00383 
00384 
00386         // calculate C, inv(C) and T
00388         for (unsigned i=0; i<DIM; i++)
00389         {
00390             for (unsigned M=0; M<DIM; M++)
00391             {
00392                 F(i,M) = (i==M?1:0) + grad_u(i,M);
00393             }
00394         }
00395 
00396         C = prod(trans(F),F);
00397         inv_C = Inverse(C);
00398         inv_F = Inverse(F);
00399 
00400         double detF = Determinant(F);
00401 
00402         ComputeStressAndStressDerivative(p_material_law, C, inv_C, pressure, rElement.GetIndex(), current_quad_point_global_index,
00403                                          T, dTdE, assembleJacobian);
00404 
00405 
00407         // residual vector
00409         if (assembleResidual)
00410         {
00411             F_T = prod(F,T);
00412             F_T_grad_quad_phi = prod(F_T, grad_quad_phi);
00413             
00414             for (unsigned index=0; index<NUM_NODES_PER_ELEMENT*DIM; index++)
00415             {
00416                 unsigned spatial_dim = index%DIM;
00417                 unsigned node_index = (index-spatial_dim)/DIM;
00418 
00419                 rBElem(index) +=  - this->mDensity
00420                                   * body_force(spatial_dim)
00421                                   * quad_phi(node_index)
00422                                   * wJ;
00423                                   
00424                 // the  T(M,N)*F(spatial_dim,M)*grad_quad_phi(N,node_index)  term                  
00425                 rBElem(index) +=   F_T_grad_quad_phi(spatial_dim,node_index)
00426                                  * wJ;
00427             }
00428 
00429             for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00430             {
00431                 rBElem( NUM_NODES_PER_ELEMENT*DIM + vertex_index ) +=   linear_phi(vertex_index)
00432                                                                       * (detF - 1)
00433                                                                       * wJ;
00434             }
00435         }
00436         
00437         
00439         // Jacobian matrix
00441         if (assembleJacobian) 
00442         {
00443             // save trans(grad_quad_phi) * invF
00444             grad_quad_phi_times_invF = prod(trans_grad_quad_phi, inv_F);
00445             
00447             // Set up the tensor dSdF
00448             //
00449             // dSdF as a function of T and dTdE (which is what the material law returns) is given by:
00450             //
00451             // dS_{Mi}/dF_{jN} = (dT_{MN}/dC_{PQ}+dT_{MN}/dC_{PQ}) F{iP} F_{jQ}  + T_{MN} delta_{ij}
00452             //
00453             // todo1: this should probably move into the material law (but need to make sure 
00454             // memory is handled efficiently
00455             // todo2: get material law to return this immediately, not dTdE
00457 
00458             // set up the tensor 0.5(dTdE(M,N,P,Q) + dTdE(M,N,Q,P))
00459             for (unsigned M=0; M<DIM; M++)
00460             {
00461                 for (unsigned N=0; N<DIM; N++)
00462                 {
00463                     for (unsigned P=0; P<DIM; P++)
00464                     {
00465                         for (unsigned Q=0; Q<DIM; Q++)
00466                         {
00467                             // this is NOT dSdF, just using this as storage space
00468                             dSdF(M,N,P,Q) = 0.5*(dTdE(M,N,P,Q) + dTdE(M,N,Q,P));
00469                         }
00470                     }
00471                 }
00472             }
00473             // This is NOT dTdE, just reusing memory. A^{MdPQ}  = F^d_N * dTdE_sym^{MNPQ}
00474             dTdE.template SetAsContractionOnSecondDimension<DIM>(F, dSdF);  
00475             // dSdF{MdPe} := F^d_N * F^e_Q * dTdE_sym^{MNPQ}
00476             dSdF.template SetAsContractionOnFourthDimension<DIM>(F, dTdE);  
00477             
00478             // now add the T_{MN} delta_{ij} term
00479             for(unsigned M=0; M<DIM; M++)
00480             {
00481                 for(unsigned N=0; N<DIM; N++)
00482                 {
00483                     for(unsigned i=0; i<DIM; i++)
00484                     {
00485                         dSdF(M,i,N,i) += T(M,N);
00486                     }
00487                 }
00488             }
00489             
00490 
00492             // Set up the tensor  
00493             //   dSdF_quad_quad(node_index1, spatial_dim1, node_index2, spatial_dim2)
00494             //            =    dS_{M,spatial_dim1}/d_F{spatial_dim2,N}      
00495             //               * grad_quad_phi(M,node_index1) 
00496             //               * grad_quad_phi(P,node_index2)
00497             //
00498             //            =    dSdF(M,spatial_index1,N,spatial_index2)
00499             //               * grad_quad_phi(M,node_index1) 
00500             //               * grad_quad_phi(P,node_index2)
00501             //            
00503             temp_tensor.template SetAsContractionOnFirstDimension<DIM>( trans_grad_quad_phi, dSdF );
00504             dSdF_quad_quad.template SetAsContractionOnThirdDimension<DIM>( trans_grad_quad_phi, temp_tensor );
00505 
00506 
00507 
00508 
00509             for (unsigned index1=0; index1<NUM_NODES_PER_ELEMENT*DIM; index1++)
00510             {
00511                 unsigned spatial_dim1 = index1%DIM;
00512                 unsigned node_index1 = (index1-spatial_dim1)/DIM;
00513 
00514 
00515                 for (unsigned index2=0; index2<NUM_NODES_PER_ELEMENT*DIM; index2++)
00516                 {
00517                     unsigned spatial_dim2 = index2%DIM;
00518                     unsigned node_index2 = (index2-spatial_dim2)/DIM;
00519 
00520                     // the dSdF*grad_quad_phi*grad_quad_phi term
00521                     rAElem(index1,index2)  +=   dSdF_quad_quad(node_index1,spatial_dim1,node_index2,spatial_dim2)
00522                                               * wJ;
00523                 }
00524 
00525                 for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00526                 {
00527                     unsigned index2 = NUM_NODES_PER_ELEMENT*DIM + vertex_index;
00528         
00529                     // the -invF(M,spatial_dim1)*grad_quad_phi(M,node_index1)*linear_phi(vertex_index) term
00530                     rAElem(index1,index2)  +=  - grad_quad_phi_times_invF(node_index1,spatial_dim1)
00531                                                * linear_phi(vertex_index)
00532                                                * wJ;
00533                 }
00534             }
00535 
00536             for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00537             {
00538                 unsigned index1 = NUM_NODES_PER_ELEMENT*DIM + vertex_index;
00539 
00540                 for (unsigned index2=0; index2<NUM_NODES_PER_ELEMENT*DIM; index2++)
00541                 {
00542                     unsigned spatial_dim2 = index2%DIM;
00543                     unsigned node_index2 = (index2-spatial_dim2)/DIM;
00544 
00545                     // same as (negative of) the opposite block (ie a few lines up), except for detF
00546                     rAElem(index1,index2) +=   detF
00547                                              * grad_quad_phi_times_invF(node_index2,spatial_dim2)
00548                                              * linear_phi(vertex_index)
00549                                              * wJ;
00550                 }
00551 
00553                 // Preconditioner matrix
00554                 // Fill the mass matrix (ie \intgl phi_i phi_j) in the
00555                 // pressure-pressure block. Note, the rest of the
00556                 // entries are filled in at the end
00558                 for (unsigned vertex_index2=0; vertex_index2<NUM_VERTICES_PER_ELEMENT; vertex_index2++)
00559                 {
00560                     unsigned index2 = NUM_NODES_PER_ELEMENT*DIM + vertex_index2;
00561                     rAElemPrecond(index1,index2) +=   linear_phi(vertex_index)
00562                                                     * linear_phi(vertex_index2)
00563                                                     * wJ;
00564                 }
00565             }
00566         }
00567     }
00568 
00569 
00570     if (assembleJacobian)
00571     {
00572         // Fill in the other blocks of the preconditioner matrix, by adding 
00573         // the jacobian matrix (this doesn't effect the pressure-pressure block 
00574         // of rAElemPrecond as the pressure-pressure block of rAElem is zero),
00575         // and the zero a block.
00576         //
00577         // The following altogether gives the preconditioner  [ A  B1^T ]
00578         //                                                    [ 0  M    ]
00579 
00580         rAElemPrecond = rAElemPrecond + rAElem;
00581         for(unsigned i=NUM_NODES_PER_ELEMENT*DIM; i<STENCIL_SIZE; i++)
00582         {
00583             for(unsigned j=0; j<NUM_NODES_PER_ELEMENT*DIM; j++)
00584             {
00585                 rAElemPrecond(i,j) = 0.0;
00586             }
00587         }
00588     }
00589 }
00590 
00591 
00592 template<size_t DIM>
00593 void NonlinearElasticitySolver<DIM>::ComputeStressAndStressDerivative(AbstractIncompressibleMaterialLaw<DIM>* pMaterialLaw,
00594                                                                       c_matrix<double,DIM,DIM>& rC, 
00595                                                                       c_matrix<double,DIM,DIM>& rInvC, 
00596                                                                       double pressure, 
00597                                                                       unsigned elementIndex,
00598                                                                       unsigned currentQuadPointGlobalIndex,
00599                                                                       c_matrix<double,DIM,DIM>& rT,
00600                                                                       FourthOrderTensor<DIM,DIM,DIM,DIM>& rDTdE,
00601                                                                       bool computeDTdE)
00602 {
00603     // just call the method on the material law
00604     pMaterialLaw->ComputeStressAndStressDerivative(rC,rInvC,pressure,rT,rDTdE,computeDTdE);
00605 }
00606 
00607 
00608 template<size_t DIM>
00609 void NonlinearElasticitySolver<DIM>::AssembleOnBoundaryElement(
00610             BoundaryElement<DIM-1,DIM>& rBoundaryElement,
00611             c_matrix<double,BOUNDARY_STENCIL_SIZE,BOUNDARY_STENCIL_SIZE>& rAelem,
00612             c_vector<double,BOUNDARY_STENCIL_SIZE>& rBelem,
00613             c_vector<double,DIM>& rTraction,
00614             bool assembleResidual,
00615             bool assembleJacobian)
00616 {
00617     rAelem.clear();
00618     rBelem.clear();
00619 
00620     if (assembleJacobian && !assembleResidual)
00621     {
00622         // nothing to do
00623         return;
00624     }
00625 
00626     c_vector<double, DIM> weighted_direction;
00627     double jacobian_determinant;
00628     this->mpQuadMesh->GetWeightedDirectionForBoundaryElement(rBoundaryElement.GetIndex(), weighted_direction, jacobian_determinant);
00629 
00630     c_vector<double,NUM_NODES_PER_BOUNDARY_ELEMENT> phi;
00631 
00632     for (unsigned quad_index=0; quad_index<this->mpBoundaryQuadratureRule->GetNumQuadPoints(); quad_index++)
00633     {
00634         double wJ = jacobian_determinant * this->mpBoundaryQuadratureRule->GetWeight(quad_index);
00635 
00636         const ChastePoint<DIM-1>& quad_point = this->mpBoundaryQuadratureRule->rGetQuadPoint(quad_index);
00637 
00638         QuadraticBasisFunction<DIM-1>::ComputeBasisFunctions(quad_point, phi);
00639 
00640         // get the required traction, interpolating X (slightly inefficiently, as interpolating
00641         // using quad bases) if necessary.
00642         c_vector<double,DIM> traction = zero_vector<double>(DIM);
00643         if (this->mUsingTractionBoundaryConditionFunction)
00644         {
00645             c_vector<double,DIM> X = zero_vector<double>(DIM);
00646             for (unsigned node_index=0; node_index<NUM_NODES_PER_BOUNDARY_ELEMENT; node_index++)
00647             {
00648                 X += phi(node_index)*this->mpQuadMesh->GetNode( rBoundaryElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
00649             }
00650             traction = (*(this->mpTractionBoundaryConditionFunction))(X, this->mCurrentTime);
00651         }
00652         else
00653         {
00654             traction = rTraction;
00655         }
00656 
00657 
00658         for (unsigned index=0; index<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index++)
00659         {
00660             unsigned spatial_dim = index%DIM;
00661             unsigned node_index = (index-spatial_dim)/DIM;
00662 
00663             assert(node_index < NUM_NODES_PER_BOUNDARY_ELEMENT);
00664 
00665             rBelem(index) -=    traction(spatial_dim)
00666                               * phi(node_index)
00667                               * wJ;
00668         }
00669     }
00670 }
00671 
00672 template<size_t DIM>
00673 void NonlinearElasticitySolver<DIM>::FormInitialGuess()
00674 {
00675     this->mCurrentSolution.resize(this->mNumDofs, 0.0);
00676 
00677     for (unsigned i=0; i<this->mpQuadMesh->GetNumElements(); i++)
00678     {
00679         double zero_strain_pressure;
00680         if (this->mMaterialLaws.size()==1)
00681         {
00682             // homogeneous
00683             zero_strain_pressure = this->mMaterialLaws[0]->GetZeroStrainPressure();
00684         }
00685         else
00686         {
00687             // heterogeneous
00688             zero_strain_pressure = this->mMaterialLaws[i]->GetZeroStrainPressure();
00689         }
00690 
00691         // loop over vertices and set pressure solution to be zero-strain-pressure
00692         for (unsigned j=0; j<NUM_VERTICES_PER_ELEMENT; j++)
00693         {
00694             // At the moment we assume the vertices are the first num_vertices nodes in the list of nodes
00695             // in the mesh. Hence:
00696             unsigned vertex_index = this->mpQuadMesh->GetElement(i)->GetNodeGlobalIndex(j);
00697             assert(vertex_index < this->mpQuadMesh->GetNumVertices());
00698 
00699             // In the future (or currently with AlexW's adaptive quadratic mesh class (project work)), we
00700             // will want to use this instead:
00701             //unsigned vertex_index = this->mpQuadMesh->GetVertexIndexOfNode(this->mpQuadMesh->GetElement(i)->GetNodeGlobalIndex(j));
00702             
00703             this->mCurrentSolution[ DIM*this->mpQuadMesh->GetNumNodes() + vertex_index ] = zero_strain_pressure;
00704         }
00705     }
00706 }
00707 
00708 
00709 template<size_t DIM>
00710 NonlinearElasticitySolver<DIM>::NonlinearElasticitySolver(
00711             QuadraticMesh<DIM>* pQuadMesh,
00712             AbstractIncompressibleMaterialLaw<DIM>* pMaterialLaw,
00713             c_vector<double,DIM> bodyForce,
00714             double density,
00715             std::string outputDirectory,
00716             std::vector<unsigned>& fixedNodes,
00717             std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00718     : AbstractNonlinearElasticitySolver<INCOMPRESSIBLE,DIM>(pQuadMesh,
00719                                                             bodyForce, density,
00720                                                             outputDirectory, fixedNodes)
00721 {
00722     assert(pMaterialLaw != NULL);
00723     mMaterialLaws.push_back(pMaterialLaw);
00724 
00725     Initialise(pFixedNodeLocations);
00726     FormInitialGuess();
00727 }
00728 
00729 
00730 template<size_t DIM>
00731 NonlinearElasticitySolver<DIM>::NonlinearElasticitySolver(
00732             QuadraticMesh<DIM>* pQuadMesh,
00733             std::vector<AbstractIncompressibleMaterialLaw<DIM>*>& rMaterialLaws,
00734             c_vector<double,DIM> bodyForce,
00735             double density,
00736             std::string outputDirectory,
00737             std::vector<unsigned>& fixedNodes,
00738             std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00739     : AbstractNonlinearElasticitySolver<INCOMPRESSIBLE,DIM>(pQuadMesh,
00740                                                             bodyForce, density,
00741                                                             outputDirectory, fixedNodes)
00742 {
00743     mMaterialLaws.resize(rMaterialLaws.size(), NULL);
00744     for (unsigned i=0; i<mMaterialLaws.size(); i++)
00745     {
00746         assert(rMaterialLaws[i] != NULL);
00747         mMaterialLaws[i] = rMaterialLaws[i];
00748     }
00749 
00750     assert(rMaterialLaws.size()==pQuadMesh->GetNumElements());
00751     Initialise(pFixedNodeLocations);
00752     FormInitialGuess();
00753 }
00754 
00755 
00756 template<size_t DIM>
00757 NonlinearElasticitySolver<DIM>::~NonlinearElasticitySolver()
00758 {
00759 //    //Post-hoc debugging
00760 //    Mat matrix = this->mpLinearSystem->rGetLhsMatrix();
00761 //    for(unsigned i=0; i<this->mpQuadMesh->GetNumNodes(); i++)
00762 //    {
00763 //        unsigned elements = this->mpQuadMesh->GetNode(i)->GetNumContainingElements();
00764 //        if(i<this->mpQuadMesh->GetNumVertices()) // then this is a vertex
00765 //        {
00766 //            TRACE("VERT");
00767 //        }
00768 //        else
00769 //        {
00770 //            TRACE("INT");
00771 //        }
00772 //        PRINT_3_VARIABLES(DIM, i, elements);
00773 //        for (unsigned dim=0; dim<DIM; dim++)
00774 //        {
00775 //            int row=DIM*i+dim;
00776 //            int ncols;
00777 //            MatGetRow(matrix, row, &ncols, NULL, NULL);
00778 //            PRINT_3_VARIABLES(i,row,ncols);
00779 //        }
00780 // 
00781 //        if(i<this->mpQuadMesh->GetNumVertices()) // then this is a vertex
00782 //        {
00783 //            int row=DIM*this->mpQuadMesh->GetNumNodes() + i;
00784 //            int ncols;
00785 //            MatGetRow(matrix, row, &ncols, NULL, NULL);
00786 //            PRINT_3_VARIABLES(i,row,ncols);
00787 //        }
00788 //    }
00789 
00790 
00791 }
00792 
00793 
00794 
00795 template<size_t DIM>
00796 std::vector<double>& NonlinearElasticitySolver<DIM>::rGetPressures()
00797 {
00798     mPressures.clear();
00799     mPressures.resize(this->mpQuadMesh->GetNumVertices());
00800 
00801     for (unsigned i=0; i<this->mpQuadMesh->GetNumVertices(); i++)
00802     {
00803         mPressures[i] = this->mCurrentSolution[DIM*this->mpQuadMesh->GetNumNodes() + i];
00804     }
00805     return mPressures;
00806 }
00807 
00808 
00809 
00811 // Explicit instantiation
00813 
00814 //template class NonlinearElasticitySolver<1>;
00815 template class NonlinearElasticitySolver<2>;
00816 template class NonlinearElasticitySolver<3>;

Generated on Mon Apr 18 11:35:36 2011 for Chaste by  doxygen 1.5.5