NonlinearElasticityAssembler.cpp

00001 /*
00002 
00003 Copyright (C) University of Oxford, 2005-2009
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 "NonlinearElasticityAssembler.hpp"
00042 
00043 #include "LinearBasisFunction.hpp"
00044 #include "QuadraticBasisFunction.hpp"
00045 
00046 
00047 template<size_t DIM>
00048 void NonlinearElasticityAssembler<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     // Get an iterator over the elements of the mesh
00067     typename AbstractMesh<DIM, DIM>::ElementIterator
00068         iter = mpQuadMesh->GetElementIteratorBegin();
00069 
00070     c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem;
00071     // the (element) preconditioner matrix: this is the same as the jacobian, but
00072     // with the mass matrix (ie \intgl phi_i phi_j) in the pressure-pressure block.
00073     c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem_precond; 
00074     c_vector<double, STENCIL_SIZE> b_elem;
00075 
00077     // loop over elements
00079     while (iter != mpQuadMesh->GetElementIteratorEnd())
00080     {
00081         //if(assembleJacobian)
00082         //{
00083         //    LOG_AND_COUT(1, "Element " << (*iter)->GetIndex() << " of " << this->mpQuadMesh->GetNumElements());
00084         //}
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 
00092             unsigned p_indices[STENCIL_SIZE];
00093             for(unsigned i=0; i<NUM_NODES_PER_ELEMENT; i++)
00094             {
00095                 for(unsigned j=0; j<DIM; j++)
00096                 {
00097                     p_indices[DIM*i+j] = DIM*element.GetNodeGlobalIndex(i) + j;
00098                 }
00099             }
00100 
00101             for(unsigned i=0; i<NUM_VERTICES_PER_ELEMENT; i++)
00102             {
00103                 p_indices[DIM*NUM_NODES_PER_ELEMENT + i] = DIM*mpQuadMesh->GetNumNodes() + element.GetNodeGlobalIndex(i);
00104             }
00105 
00106             if (assembleJacobian)
00107             {
00108                 this->mpLinearSystem->AddLhsMultipleValues(p_indices, a_elem);
00109                 this->mpPreconditionMatrixLinearSystem->AddLhsMultipleValues(p_indices, a_elem_precond); 
00110             }
00111 
00112             if (assembleResidual)
00113             {
00114                 this->mpLinearSystem->AddRhsMultipleValues(p_indices, b_elem);
00115             }
00116         }
00117 
00118         iter++;
00119     }
00120     
00122     // loop over specified boundary elements and compute
00123     // surface traction terms 
00125     c_vector<double, BOUNDARY_STENCIL_SIZE> b_boundary_elem;
00126     c_matrix<double, BOUNDARY_STENCIL_SIZE, BOUNDARY_STENCIL_SIZE> a_boundary_elem;
00127     if (mBoundaryElements.size()>0)
00128     {
00129         for (unsigned i=0; i<mBoundaryElements.size(); i++)
00130         {
00131             BoundaryElement<DIM-1,DIM>& r_boundary_element = *(mBoundaryElements[i]);
00132             AssembleOnBoundaryElement(r_boundary_element, a_boundary_elem, b_boundary_elem, this->mSurfaceTractions[i], assembleResidual, assembleJacobian);
00133 
00134             unsigned p_indices[BOUNDARY_STENCIL_SIZE];
00135             for (unsigned i=0; i<NUM_NODES_PER_BOUNDARY_ELEMENT; i++)
00136             {
00137                 for (unsigned j=0; j<DIM; j++)
00138                 {
00139                     p_indices[DIM*i+j] = DIM*r_boundary_element.GetNodeGlobalIndex(i) + j;
00140                 }
00141             }
00142 
00143             for (unsigned i=0; i<DIM /*vertices per boundary elem */; i++)
00144             {
00145                 p_indices[DIM*NUM_NODES_PER_BOUNDARY_ELEMENT + i] = DIM*mpQuadMesh->GetNumNodes() + r_boundary_element.GetNodeGlobalIndex(i);
00146             }
00147 
00148             if (assembleJacobian)
00149             {
00150                 this->mpLinearSystem->AddLhsMultipleValues(p_indices, a_boundary_elem);
00151                 this->mpPreconditionMatrixLinearSystem->AddLhsMultipleValues(p_indices, a_boundary_elem); 
00152             }
00153 
00154             if (assembleResidual)
00155             {
00156                 this->mpLinearSystem->AddRhsMultipleValues(p_indices, b_boundary_elem);
00157             }
00158 
00159             // some extra checking
00160             if (DIM==2)
00161             {   
00162                 assert(8==BOUNDARY_STENCIL_SIZE);
00163                 assert(b_boundary_elem(6)==0);
00164                 assert(b_boundary_elem(7)==0);
00165             }
00166         }
00167     }
00168 
00169 #ifndef ___USE_DEALII_LINEAR_SYSTEM___
00170     if (assembleResidual)
00171     {
00172         this->mpLinearSystem->AssembleRhsVector();
00173     }
00174     if (assembleJacobian)
00175     {
00176         this->mpLinearSystem->AssembleIntermediateLhsMatrix();
00177     }
00178 #endif
00179     if (assembleJacobian)
00180     {
00181         this->mpPreconditionMatrixLinearSystem->AssembleIntermediateLhsMatrix(); 
00182     }
00183 
00184 
00185     // Apply dirichlet boundary conditions
00186     this->ApplyBoundaryConditions(assembleJacobian);
00187 
00188 #ifndef ___USE_DEALII_LINEAR_SYSTEM___
00189     if (assembleResidual)
00190     {
00191         this->mpLinearSystem->AssembleRhsVector();
00192     }
00193     if (assembleJacobian)
00194     {
00195         this->mpLinearSystem->AssembleFinalLhsMatrix();
00196     }
00197 #endif
00198     if (assembleJacobian)
00199     {
00200         this->mpPreconditionMatrixLinearSystem->AssembleFinalLhsMatrix(); 
00201     }
00202 }
00203 
00204 template<size_t DIM>
00205 void NonlinearElasticityAssembler<DIM>::AssembleOnElement(
00206             Element<DIM, DIM>& rElement,
00207             c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElem,
00208             c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElemPrecond, 
00209             c_vector<double, STENCIL_SIZE>& rBElem,
00210             bool assembleResidual,
00211             bool assembleJacobian)
00212 {
00213     c_matrix<double, DIM, DIM> jacobian, 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     c_vector<double, NUM_VERTICES_PER_ELEMENT> linear_phi;
00252     c_vector<double, NUM_NODES_PER_ELEMENT> quad_phi;
00253     c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> grad_quad_phi;
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     
00276     for (unsigned quadrature_index=0; quadrature_index < this->mpQuadratureRule->GetNumQuadPoints(); quadrature_index++)
00277     {
00278         double wJ = jacobian_determinant * this->mpQuadratureRule->GetWeight(quadrature_index);
00279 
00280         const ChastePoint<DIM>& quadrature_point = this->mpQuadratureRule->rGetQuadPoint(quadrature_index);
00281 
00283         // set up basis function info
00285         LinearBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, linear_phi);
00286         QuadraticBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, quad_phi);
00287         QuadraticBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(quadrature_point, inverse_jacobian, grad_quad_phi);
00288         
00289 
00291         // get the body force, interpolating X if necessary
00293         c_vector<double,DIM> body_force;
00294 
00295         if(this->mUsingBodyForceFunction)
00296         {
00297             c_vector<double,DIM> X = zero_vector<double>(DIM);
00298             // interpolate X (using the vertices and the /linear/ bases, as no curvilinear elements
00299             for(unsigned node_index=0; node_index<NUM_VERTICES_PER_ELEMENT; node_index++)
00300             {
00301                 X += linear_phi(node_index)*mpQuadMesh->GetNode( rElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
00302             }
00303             body_force = (*(this->mpBodyForceFunction))(X);
00304         }
00305         else
00306         {
00307             body_force = this->mBodyForce;
00308         }
00309 
00311         // interpolate grad_u and p
00313         static c_matrix<double,DIM,DIM> grad_u; // grad_u = (du_i/dX_M)
00314         grad_u = zero_matrix<double>(DIM,DIM);  // must be on new line!!
00315 
00316         for(unsigned node_index=0; node_index<NUM_NODES_PER_ELEMENT; node_index++)
00317         {
00318             for (unsigned i=0; i<DIM; i++)
00319             {
00320                 for(unsigned M=0; M<DIM; M++)
00321                 {
00322                     grad_u(i,M) += grad_quad_phi(M,node_index)*element_current_displacements(i,node_index);
00323                 }
00324             }
00325         }
00326 
00327         double pressure = 0;
00328         for(unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00329         {
00330             pressure += linear_phi(vertex_index)*element_current_pressures(vertex_index);
00331         }
00332 
00333     
00335         // calculate C, inv(C) and T
00337         static c_matrix<double,DIM,DIM> F;
00338         static c_matrix<double,DIM,DIM> C;
00339         static c_matrix<double,DIM,DIM> inv_C;
00340         static c_matrix<double,DIM,DIM> inv_F;
00341         static c_matrix<double,DIM,DIM> T;
00342 
00343         for (unsigned i=0; i<DIM; i++)
00344         {
00345             for (unsigned M=0; M<DIM; M++)
00346             {
00347                 F(i,M) = (i==M?1:0) + grad_u(i,M);
00348             }
00349         }
00350 
00351         C = prod(trans(F),F);
00352         inv_C = Inverse(C);
00353         inv_F = Inverse(F);
00354 
00355         double detF = Determinant(F);
00356 
00357         p_material_law->ComputeStressAndStressDerivative(C,inv_C,pressure,T,this->dTdE,assembleJacobian);
00358 
00360         // residual vector
00362         if (assembleResidual)
00363         {
00364             for(unsigned index=0; index<NUM_NODES_PER_ELEMENT*DIM; index++)
00365             {
00366                 unsigned spatial_dim = index%DIM;
00367                 unsigned node_index = (index-spatial_dim)/DIM;
00368 
00369                 assert(node_index < NUM_NODES_PER_ELEMENT);
00370 
00371                 rBElem(index) +=  - this->mDensity 
00372                                   * body_force(spatial_dim)
00373                                   * quad_phi(node_index)
00374                                   * wJ;
00375 
00376                 for (unsigned M=0; M<DIM; M++)
00377                 {
00378                     for (unsigned N=0; N<DIM; N++)
00379                     {
00380                         rBElem(index) +=   T(M,N)
00381                                          * F(spatial_dim,M)
00382                                          * grad_quad_phi(N,node_index)
00383                                          * wJ;
00384                     }
00385                 }
00386             }
00387             
00388             for(unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00389             {
00390                 rBElem( NUM_NODES_PER_ELEMENT*DIM + vertex_index ) +=   linear_phi(vertex_index)
00391                                                                       * (detF - 1)
00392                                                                       * wJ;
00393             }
00394         }
00395 
00397         // Jacobian matrix
00399         if(assembleJacobian)
00400         {
00401             for(unsigned index1=0; index1<NUM_NODES_PER_ELEMENT*DIM; index1++)
00402             {
00403                 unsigned spatial_dim1 = index1%DIM;
00404                 unsigned node_index1 = (index1-spatial_dim1)/DIM;
00405                 
00406                 
00407                 for(unsigned index2=0; index2<NUM_NODES_PER_ELEMENT*DIM; index2++)
00408                 {
00409                     unsigned spatial_dim2 = index2%DIM;
00410                     unsigned node_index2 = (index2-spatial_dim2)/DIM;
00411 
00412                     for (unsigned M=0; M<DIM; M++)
00413                     {
00414                         for (unsigned N=0; N<DIM; N++)
00415                         {
00416                             rAElem(index1,index2) +=   T(M,N)
00417                                                      * grad_quad_phi(N,node_index1)
00418                                                      * grad_quad_phi(M,node_index2)
00419                                                      * (spatial_dim1==spatial_dim2?1:0)
00420                                                      * wJ;
00421 
00422                             for (unsigned P=0; P<DIM; P++)
00423                             {
00424                                 for (unsigned Q=0; Q<DIM; Q++)
00425                                 {
00426                                     rAElem(index1,index2)  +=   0.5
00427                                                               * this->dTdE(M,N,P,Q)
00428                                                               * (
00429                                                                   grad_quad_phi(P,node_index2)
00430                                                                 * F(spatial_dim2,Q)
00431                                                                    +
00432                                                                   grad_quad_phi(Q,node_index2)
00433                                                                 * F(spatial_dim2,P)
00434                                                                  )
00435                                                               * F(spatial_dim1,M)
00436                                                               * grad_quad_phi(N,node_index1)
00437                                                               * wJ;
00438                                 }
00439                             }
00440                         }
00441                     }
00442                 }
00443                 for(unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00444                 {
00445                     unsigned index2 = NUM_NODES_PER_ELEMENT*DIM + vertex_index;
00446                     
00447                     for (unsigned M=0; M<DIM; M++)
00448                     {
00449                         for (unsigned N=0; N<DIM; N++)
00450                         {
00451                             rAElem(index1,index2)  +=  - F(spatial_dim1,M)
00452                                                        * inv_C(M,N)
00453                                                        * grad_quad_phi(N,node_index1)
00454                                                        * linear_phi(vertex_index)
00455                                                        * wJ;
00456                         }
00457                     }                 
00458                 }
00459             }
00460 
00461             for(unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00462             {
00463                 unsigned index1 = NUM_NODES_PER_ELEMENT*DIM + vertex_index;
00464 
00465                 for(unsigned index2=0; index2<NUM_NODES_PER_ELEMENT*DIM; index2++)
00466                 {
00467                     unsigned spatial_dim2 = index2%DIM;
00468                     unsigned node_index2 = (index2-spatial_dim2)/DIM;
00469 
00470                     for (unsigned M=0; M<DIM; M++)
00471                     {
00472                         rAElem(index1,index2) +=   linear_phi(vertex_index)
00473                                                  * detF
00474                                                  * inv_F(M,spatial_dim2)
00475                                                  * grad_quad_phi(M,node_index2)
00476                                                  * wJ;
00477                     }
00478                 }
00479                 
00481                 // Preconditioner matrix
00482                 // Fill the mass matrix (ie \intgl phi_i phi_j) in the 
00483                 // pressure-pressure block. Note, the rest of the 
00484                 // entries are filled in at the end
00486                 for(unsigned vertex_index2=0; vertex_index2<NUM_VERTICES_PER_ELEMENT; vertex_index2++) 
00487                 { 
00488                     unsigned index2 = NUM_NODES_PER_ELEMENT*DIM + vertex_index2;
00489                     rAElemPrecond(index1,index2) +=   linear_phi(vertex_index)
00490                                                     * linear_phi(vertex_index2) 
00491                                                     * wJ; 
00492                 } 
00493             }
00494         }
00495     }
00496 
00497 
00498     if (assembleJacobian) 
00499     { 
00500         // Fill in the other blocks of the preconditioner matrix. (This doesn't
00501         // effect the pressure-pressure block of the rAElemPrecond but the 
00502         // pressure-pressure block of rAElem is zero
00503         rAElemPrecond = rAElemPrecond + rAElem; 
00504     } 
00505 }
00506 
00507 template<size_t DIM>
00508 void NonlinearElasticityAssembler<DIM>::AssembleOnBoundaryElement(
00509             BoundaryElement<DIM-1,DIM>& rBoundaryElement,
00510             c_matrix<double,BOUNDARY_STENCIL_SIZE,BOUNDARY_STENCIL_SIZE>& rAelem,
00511             c_vector<double,BOUNDARY_STENCIL_SIZE>& rBelem,
00512             c_vector<double,DIM>& rTraction,
00513             bool assembleResidual,
00514             bool assembleJacobian)
00515 {
00516     rAelem.clear();
00517     rBelem.clear();
00518 
00519     if(assembleJacobian && !assembleResidual)
00520     {
00521         // nothing to do
00522         return;
00523     }
00524 
00525     c_vector<double, DIM> weighted_direction;
00526     double jacobian_determinant;
00527     mpQuadMesh->GetWeightedDirectionForBoundaryElement(rBoundaryElement.GetIndex(), weighted_direction, jacobian_determinant);
00528 
00529     c_vector<double,NUM_NODES_PER_BOUNDARY_ELEMENT> phi;
00530 
00531     for (unsigned quad_index=0; quad_index<this->mpBoundaryQuadratureRule->GetNumQuadPoints(); quad_index++)
00532     {
00533         double wJ = jacobian_determinant * this->mpBoundaryQuadratureRule->GetWeight(quad_index);
00534 
00535         const ChastePoint<DIM-1>& quad_point = this->mpBoundaryQuadratureRule->rGetQuadPoint(quad_index);
00536 
00537         QuadraticBasisFunction<DIM-1>::ComputeBasisFunctions(quad_point, phi);
00538 
00539         // get the required traction, interpolating X (slightly inefficiently, as interpolating
00540         // using quad bases) if necessary.
00541         c_vector<double,DIM> traction = zero_vector<double>(DIM);
00542         if(this->mUsingTractionBoundaryConditionFunction)
00543         {
00544             c_vector<double,DIM> X = zero_vector<double>(DIM);
00545             for(unsigned node_index=0; node_index<NUM_NODES_PER_BOUNDARY_ELEMENT; node_index++)
00546             {
00547                 X += phi(node_index)*mpQuadMesh->GetNode( rBoundaryElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
00548             }
00549             traction = (*(this->mpTractionBoundaryConditionFunction))(X);
00550         }
00551         else
00552         {
00553             traction = rTraction;
00554         }
00555 
00556 
00557         for(unsigned index=0; index<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index++)
00558         {
00559             unsigned spatial_dim = index%DIM;
00560             unsigned node_index = (index-spatial_dim)/DIM;
00561 
00562             assert(node_index < NUM_NODES_PER_BOUNDARY_ELEMENT);
00563 
00564             rBelem(index) -=    traction(spatial_dim)
00565                               * phi(node_index)
00566                               * wJ;
00567         }
00568     }
00569 }
00570 
00571 template<size_t DIM>
00572 void NonlinearElasticityAssembler<DIM>::FormInitialGuess()
00573 {
00574     this->mCurrentSolution.resize(this->mNumDofs, 0.0);
00575     
00576     for(unsigned i=0; i<mpQuadMesh->GetNumElements(); i++)
00577     {
00578         double zero_strain_pressure;
00579         if(this->mMaterialLaws.size()==1)
00580         {
00581             // homogeneous
00582             zero_strain_pressure = this->mMaterialLaws[0]->GetZeroStrainPressure();
00583         }
00584         else
00585         {
00586             // heterogeneous
00587             zero_strain_pressure = this->mMaterialLaws[i]->GetZeroStrainPressure();
00588         }
00589         
00590         // loop over vertices and set pressure solution to be zero-strain-pressure
00591         for(unsigned j=0; j<NUM_VERTICES_PER_ELEMENT; j++)
00592         {
00593             unsigned index = mpQuadMesh->GetElement(i)->GetNodeGlobalIndex(j);
00594             this->mCurrentSolution[ DIM*mpQuadMesh->GetNumNodes() + index ] = zero_strain_pressure;
00595         }
00596     }
00597 }
00598 
00599 template<size_t DIM>
00600 void NonlinearElasticityAssembler<DIM>::Initialise(std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00601 {
00602     assert(mpQuadMesh);
00603 
00604     for(unsigned i=0; i<this->mFixedNodes.size(); i++)
00605     {
00606         assert(this->mFixedNodes[i] < mpQuadMesh->GetNumNodes());
00607     }
00608 
00609     this->mpQuadratureRule = new GaussianQuadratureRule<DIM>(3);
00610     this->mpBoundaryQuadratureRule = new GaussianQuadratureRule<DIM-1>(3);
00611     
00612     FormInitialGuess();
00613 
00614     // compute the displacements at each of the fixed nodes, given the
00615     // fixed nodes locations.
00616     if(pFixedNodeLocations == NULL)
00617     {
00618         this->mFixedNodeDisplacements.clear();
00619         for(unsigned i=0; i<this->mFixedNodes.size(); i++)
00620         {
00621             this->mFixedNodeDisplacements.push_back(zero_vector<double>(DIM));
00622         }
00623     }
00624     else
00625     {
00626         assert(pFixedNodeLocations->size()==this->mFixedNodes.size());
00627         for(unsigned i=0; i<this->mFixedNodes.size(); i++)
00628         {
00629             unsigned index = this->mFixedNodes[i];
00630             c_vector<double,DIM> displacement = (*pFixedNodeLocations)[i]-mpQuadMesh->GetNode(index)->rGetLocation();
00631             this->mFixedNodeDisplacements.push_back(displacement);
00632         }
00633     }
00634     assert(this->mFixedNodeDisplacements.size()==this->mFixedNodes.size());
00635 }
00636 
00637 template<size_t DIM>
00638 NonlinearElasticityAssembler<DIM>::NonlinearElasticityAssembler(
00639             QuadraticMesh<DIM>* pQuadMesh,
00640             AbstractIncompressibleMaterialLaw<DIM>* pMaterialLaw,
00641             c_vector<double,DIM> bodyForce,
00642             double density,
00643             std::string outputDirectory,
00644             std::vector<unsigned>& fixedNodes,
00645             std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00646     : AbstractNonlinearElasticityAssembler<DIM>(DIM*pQuadMesh->GetNumNodes()+pQuadMesh->GetNumVertices(),
00647                                                 pMaterialLaw, bodyForce, density,
00648                                                 outputDirectory, fixedNodes),
00649       mpQuadMesh(pQuadMesh)
00650 {
00651 #ifdef ___USE_DEALII_LINEAR_SYSTEM___
00652     // has to be done in parent as needs mesh
00653     this->mpLinearSystem = new DealiiLinearSystem(*pQuadMesh);
00654 #endif
00655 
00656     Initialise(pFixedNodeLocations);
00657 }
00658 
00659 
00660 template<size_t DIM>
00661 NonlinearElasticityAssembler<DIM>::NonlinearElasticityAssembler(
00662             QuadraticMesh<DIM>* pQuadMesh,
00663             std::vector<AbstractIncompressibleMaterialLaw<DIM>*>& rMaterialLaws,
00664             c_vector<double,DIM> bodyForce,
00665             double density,
00666             std::string outputDirectory,
00667             std::vector<unsigned>& fixedNodes,
00668             std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00669     : AbstractNonlinearElasticityAssembler<DIM>(DIM*pQuadMesh->GetNumNodes()+pQuadMesh->GetNumVertices(),
00670                                                 rMaterialLaws, bodyForce, density,
00671                                                 outputDirectory, fixedNodes),
00672       mpQuadMesh(pQuadMesh)
00673 {
00674 #ifdef ___USE_DEALII_LINEAR_SYSTEM___
00676     this->mpLinearSystem = new DealiiLinearSystem(*pQuadMesh);
00677 #endif
00678 
00679     assert(rMaterialLaws.size()==pQuadMesh->GetNumElements());
00680     Initialise(pFixedNodeLocations);
00681 }
00682 
00683 
00684 template<size_t DIM>
00685 NonlinearElasticityAssembler<DIM>::~NonlinearElasticityAssembler()
00686 {
00687     delete mpQuadratureRule;
00688     delete mpBoundaryQuadratureRule;
00689 }
00690 
00691 template<size_t DIM>
00692 void NonlinearElasticityAssembler<DIM>::SetSurfaceTractionBoundaryConditions(
00693             std::vector<BoundaryElement<DIM-1,DIM>*>& rBoundaryElements,
00694             std::vector<c_vector<double,DIM> >& rSurfaceTractions)
00695 {
00696     assert(rBoundaryElements.size()==rSurfaceTractions.size());
00697     mBoundaryElements = rBoundaryElements;
00698     this->mSurfaceTractions = rSurfaceTractions;
00699 }
00700 
00701 template<size_t DIM>
00702 void NonlinearElasticityAssembler<DIM>::SetFunctionalTractionBoundaryCondition(
00703             std::vector<BoundaryElement<DIM-1,DIM>*> rBoundaryElements,
00704             c_vector<double,DIM> (*pFunction)(c_vector<double,DIM>&))
00705 {
00706     mBoundaryElements = rBoundaryElements;
00707     this->mUsingTractionBoundaryConditionFunction = true;
00708     this->mpTractionBoundaryConditionFunction = pFunction;
00709 }
00710 
00711 template<size_t DIM>
00712 std::vector<double>& NonlinearElasticityAssembler<DIM>::rGetPressures()
00713 {
00714     this->mPressures.clear();
00715     this->mPressures.resize(mpQuadMesh->GetNumVertices());
00716     
00717     for (unsigned i=0; i<mpQuadMesh->GetNumVertices(); i++)
00718     {
00719         this->mPressures[i] = this->mCurrentSolution[DIM*mpQuadMesh->GetNumNodes() + i];
00720     }
00721     return this->mPressures;
00722 }
00723 
00724 template<size_t DIM>
00725 std::vector<c_vector<double,DIM> >& NonlinearElasticityAssembler<DIM>::rGetDeformedPosition()
00726 {
00727     this->mDeformedPosition.resize(mpQuadMesh->GetNumNodes(), zero_vector<double>(DIM));
00728     for (unsigned i=0; i<mpQuadMesh->GetNumNodes(); i++)
00729     {
00730         for (unsigned j=0; j<DIM; j++)
00731         {
00732             this->mDeformedPosition[i](j) = mpQuadMesh->GetNode(i)->rGetLocation()[j] + this->mCurrentSolution[DIM*i+j];
00733         }
00734     }
00735     return this->mDeformedPosition;
00736 }
00737 
00738 
00740 // Explicit instantiation
00742 
00743 //template class NonlinearElasticityAssembler<1>;
00744 template class NonlinearElasticityAssembler<2>;
00745 template class NonlinearElasticityAssembler<3>;

Generated on Wed Mar 18 12:51:56 2009 for Chaste by  doxygen 1.5.5