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