00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041 #include "CompressibleNonlinearElasticitySolver.hpp"
00042 #include "LinearBasisFunction.hpp"
00043 #include "QuadraticBasisFunction.hpp"
00044 #include <algorithm>
00045
00046 template<size_t DIM>
00047 void CompressibleNonlinearElasticitySolver<DIM>::AssembleSystem(bool assembleResidual,
00048 bool assembleJacobian)
00049 {
00050
00051 assert(assembleResidual || assembleJacobian);
00052 assert(this->mCurrentSolution.size()==this->mNumDofs);
00053
00054
00055 if (assembleResidual)
00056 {
00057 PetscVecTools::Zero(this->mResidualVector);
00058 }
00059 if (assembleJacobian)
00060 {
00061 PetscMatTools::Zero(this->mJacobianMatrix);
00062 PetscMatTools::Zero(this->mPreconditionMatrix);
00063 }
00064
00065 c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem;
00066
00067
00068 c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem_precond;
00069 c_vector<double, STENCIL_SIZE> b_elem;
00070
00071
00072 for (typename AbstractTetrahedralMesh<DIM, DIM>::ElementIterator iter = this->mrQuadMesh.GetElementIteratorBegin();
00073 iter != this->mrQuadMesh.GetElementIteratorEnd();
00074 ++iter)
00075 {
00076 #ifdef MECH_VERY_VERBOSE
00077 if (assembleJacobian)
00078 {
00079 std::cout << "\r[" << PetscTools::GetMyRank() << "]: Element " << (*iter).GetIndex() << " of " << this->mrQuadMesh.GetNumElements() << std::flush;
00080 }
00081 #endif
00082
00083 Element<DIM, DIM>& element = *iter;
00084
00085 if (element.GetOwnership() == true)
00086 {
00087 AssembleOnElement(element, a_elem, a_elem_precond, b_elem, assembleResidual, assembleJacobian);
00088
00092
00093
00094
00095
00096
00097
00098
00099
00100 unsigned p_indices[STENCIL_SIZE];
00101 for (unsigned i=0; i<NUM_NODES_PER_ELEMENT; i++)
00102 {
00103 for (unsigned j=0; j<DIM; j++)
00104 {
00105 p_indices[DIM*i+j] = DIM*element.GetNodeGlobalIndex(i) + j;
00106 }
00107 }
00108
00109 if (assembleJacobian)
00110 {
00111 PetscMatTools::AddMultipleValues<STENCIL_SIZE>(this->mJacobianMatrix, p_indices, a_elem);
00112 PetscMatTools::AddMultipleValues<STENCIL_SIZE>(this->mPreconditionMatrix, p_indices, a_elem_precond);
00113 }
00114
00115 if (assembleResidual)
00116 {
00117 PetscVecTools::AddMultipleValues<STENCIL_SIZE>(this->mResidualVector, p_indices, b_elem);
00118 }
00119 }
00120 }
00121
00122
00123 c_vector<double, BOUNDARY_STENCIL_SIZE> b_boundary_elem;
00124 c_matrix<double, BOUNDARY_STENCIL_SIZE, BOUNDARY_STENCIL_SIZE> a_boundary_elem;
00125 if (this->mrProblemDefinition.GetTractionBoundaryConditionType() != NO_TRACTIONS)
00126 {
00127 for (unsigned bc_index=0; bc_index<this->mrProblemDefinition.rGetTractionBoundaryElements().size(); bc_index++)
00128 {
00129 BoundaryElement<DIM-1,DIM>& r_boundary_element = *(this->mrProblemDefinition.rGetTractionBoundaryElements()[bc_index]);
00130 AssembleOnBoundaryElement(r_boundary_element, a_boundary_elem, b_boundary_elem, assembleResidual, assembleJacobian, bc_index);
00131
00132 unsigned p_indices[BOUNDARY_STENCIL_SIZE];
00133 for (unsigned i=0; i<NUM_NODES_PER_BOUNDARY_ELEMENT; i++)
00134 {
00135 for (unsigned j=0; j<DIM; j++)
00136 {
00137 p_indices[DIM*i+j] = DIM*r_boundary_element.GetNodeGlobalIndex(i) + j;
00138 }
00139 }
00140
00141 if (assembleJacobian)
00142 {
00143 PetscMatTools::AddMultipleValues<BOUNDARY_STENCIL_SIZE>(this->mJacobianMatrix, p_indices, a_boundary_elem);
00144 PetscMatTools::AddMultipleValues<BOUNDARY_STENCIL_SIZE>(this->mPreconditionMatrix, p_indices, a_boundary_elem);
00145 }
00146
00147 if (assembleResidual)
00148 {
00149 PetscVecTools::AddMultipleValues<BOUNDARY_STENCIL_SIZE>(this->mResidualVector, p_indices, b_boundary_elem);
00150 }
00151 }
00152 }
00153
00154 this->FinishAssembleSystem(assembleResidual, assembleJacobian);
00155 }
00156
00157 template<size_t DIM>
00158 void CompressibleNonlinearElasticitySolver<DIM>::AssembleOnElement(
00159 Element<DIM, DIM>& rElement,
00160 c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElem,
00161 c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElemPrecond,
00162 c_vector<double, STENCIL_SIZE>& rBElem,
00163 bool assembleResidual,
00164 bool assembleJacobian)
00165 {
00166 static c_matrix<double,DIM,DIM> jacobian;
00167 static c_matrix<double,DIM,DIM> inverse_jacobian;
00168 double jacobian_determinant;
00169
00170 this->mrQuadMesh.GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
00171
00172 if (assembleJacobian)
00173 {
00174 rAElem.clear();
00175 rAElemPrecond.clear();
00176 }
00177
00178 if (assembleResidual)
00179 {
00180 rBElem.clear();
00181 }
00182
00183
00184 static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> element_current_displacements;
00185 for (unsigned II=0; II<NUM_NODES_PER_ELEMENT; II++)
00186 {
00187 for (unsigned JJ=0; JJ<DIM; JJ++)
00188 {
00189 element_current_displacements(JJ,II) = this->mCurrentSolution[DIM*rElement.GetNodeGlobalIndex(II) + JJ];
00190 }
00191 }
00192
00193
00194 static c_vector<double, NUM_VERTICES_PER_ELEMENT> linear_phi;
00195 static c_vector<double, NUM_NODES_PER_ELEMENT> quad_phi;
00196 static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> grad_quad_phi;
00197 static c_matrix<double, NUM_NODES_PER_ELEMENT, DIM> trans_grad_quad_phi;
00198
00199
00200 AbstractCompressibleMaterialLaw<DIM>* p_material_law
00201 = this->mrProblemDefinition.GetCompressibleMaterialLaw(rElement.GetIndex());
00202
00203
00204 static c_matrix<double,DIM,DIM> grad_u;
00205
00206 static c_matrix<double,DIM,DIM> F;
00207 static c_matrix<double,DIM,DIM> C;
00208 static c_matrix<double,DIM,DIM> inv_C;
00209 static c_matrix<double,DIM,DIM> inv_F;
00210 static c_matrix<double,DIM,DIM> T;
00211
00212 static c_matrix<double,DIM,DIM> F_T;
00213 static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> F_T_grad_quad_phi;
00214
00215 c_vector<double,DIM> body_force;
00216
00217 static FourthOrderTensor<DIM,DIM,DIM,DIM> dTdE;
00218 static FourthOrderTensor<DIM,DIM,DIM,DIM> dSdF;
00219
00220 static FourthOrderTensor<NUM_NODES_PER_ELEMENT,DIM,DIM,DIM> temp_tensor;
00221 static FourthOrderTensor<NUM_NODES_PER_ELEMENT,DIM,NUM_NODES_PER_ELEMENT,DIM> dSdF_quad_quad;
00222
00223 static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> temp_matrix;
00224 static c_matrix<double,NUM_NODES_PER_ELEMENT,DIM> grad_quad_phi_times_invF;
00225
00226
00227 for (unsigned quadrature_index=0; quadrature_index < this->mpQuadratureRule->GetNumQuadPoints(); quadrature_index++)
00228 {
00229
00230 unsigned current_quad_point_global_index = rElement.GetIndex()*this->mpQuadratureRule->GetNumQuadPoints()
00231 + quadrature_index;
00232
00233 double wJ = jacobian_determinant * this->mpQuadratureRule->GetWeight(quadrature_index);
00234
00235 const ChastePoint<DIM>& quadrature_point = this->mpQuadratureRule->rGetQuadPoint(quadrature_index);
00236
00237
00238 LinearBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, linear_phi);
00239 QuadraticBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, quad_phi);
00240 QuadraticBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(quadrature_point, inverse_jacobian, grad_quad_phi);
00241 trans_grad_quad_phi = trans(grad_quad_phi);
00242
00243
00244 if (assembleResidual)
00245 {
00246 switch (this->mrProblemDefinition.GetBodyForceType())
00247 {
00248 case FUNCTIONAL_BODY_FORCE:
00249 {
00250 c_vector<double,DIM> X = zero_vector<double>(DIM);
00251
00252 for (unsigned node_index=0; node_index<NUM_VERTICES_PER_ELEMENT; node_index++)
00253 {
00254 X += linear_phi(node_index)*this->mrQuadMesh.GetNode( rElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
00255 }
00256 body_force = this->mrProblemDefinition.EvaluateBodyForceFunction(X, this->mCurrentTime);
00257 break;
00258 }
00259 case CONSTANT_BODY_FORCE:
00260 {
00261 body_force = this->mrProblemDefinition.GetConstantBodyForce();
00262 break;
00263 }
00264 default:
00265 NEVER_REACHED;
00266 }
00267 }
00268
00269
00270 grad_u = zero_matrix<double>(DIM,DIM);
00271 for (unsigned node_index=0; node_index<NUM_NODES_PER_ELEMENT; node_index++)
00272 {
00273 for (unsigned i=0; i<DIM; i++)
00274 {
00275 for (unsigned M=0; M<DIM; M++)
00276 {
00277 grad_u(i,M) += grad_quad_phi(M,node_index)*element_current_displacements(i,node_index);
00278 }
00279 }
00280 }
00281
00282
00283 for (unsigned i=0; i<DIM; i++)
00284 {
00285 for (unsigned M=0; M<DIM; M++)
00286 {
00287 F(i,M) = (i==M?1:0) + grad_u(i,M);
00288 }
00289 }
00290
00291 C = prod(trans(F),F);
00292 inv_C = Inverse(C);
00293 inv_F = Inverse(F);
00294
00295 this->ComputeStressAndStressDerivative(p_material_law, C, inv_C, 0.0, rElement.GetIndex(), current_quad_point_global_index,
00296 T, dTdE, assembleJacobian);
00297
00298
00299 if (assembleResidual)
00300 {
00301 F_T = prod(F,T);
00302 F_T_grad_quad_phi = prod(F_T, grad_quad_phi);
00303
00304 for (unsigned index=0; index<NUM_NODES_PER_ELEMENT*DIM; index++)
00305 {
00306 unsigned spatial_dim = index%DIM;
00307 unsigned node_index = (index-spatial_dim)/DIM;
00308
00309 rBElem(index) += - this->mrProblemDefinition.GetDensity()
00310 * body_force(spatial_dim)
00311 * quad_phi(node_index)
00312 * wJ;
00313
00314
00315 rBElem(index) += F_T_grad_quad_phi(spatial_dim,node_index)
00316 * wJ;
00317 }
00318 }
00319
00320
00321 if (assembleJacobian)
00322 {
00323
00324 grad_quad_phi_times_invF = prod(trans_grad_quad_phi, inv_F);
00325
00327
00328
00329
00330
00331
00332
00333
00334
00335
00337
00338
00339 for (unsigned M=0; M<DIM; M++)
00340 {
00341 for (unsigned N=0; N<DIM; N++)
00342 {
00343 for (unsigned P=0; P<DIM; P++)
00344 {
00345 for (unsigned Q=0; Q<DIM; Q++)
00346 {
00347
00348 dSdF(M,N,P,Q) = 0.5*(dTdE(M,N,P,Q) + dTdE(M,N,Q,P));
00349 }
00350 }
00351 }
00352 }
00353
00354
00355 dTdE.template SetAsContractionOnSecondDimension<DIM>(F, dSdF);
00356
00357
00358 dSdF.template SetAsContractionOnFourthDimension<DIM>(F, dTdE);
00359
00360
00361 for (unsigned M=0; M<DIM; M++)
00362 {
00363 for (unsigned N=0; N<DIM; N++)
00364 {
00365 for (unsigned i=0; i<DIM; i++)
00366 {
00367 dSdF(M,i,N,i) += T(M,N);
00368 }
00369 }
00370 }
00371
00373
00374
00375
00376
00377
00378
00379
00380
00381
00382
00384 temp_tensor.template SetAsContractionOnFirstDimension<DIM>(trans_grad_quad_phi, dSdF);
00385 dSdF_quad_quad.template SetAsContractionOnThirdDimension<DIM>(trans_grad_quad_phi, temp_tensor);
00386
00387 for (unsigned index1=0; index1<NUM_NODES_PER_ELEMENT*DIM; index1++)
00388 {
00389 unsigned spatial_dim1 = index1%DIM;
00390 unsigned node_index1 = (index1-spatial_dim1)/DIM;
00391
00392 for (unsigned index2=0; index2<NUM_NODES_PER_ELEMENT*DIM; index2++)
00393 {
00394 unsigned spatial_dim2 = index2%DIM;
00395 unsigned node_index2 = (index2-spatial_dim2)/DIM;
00396
00397
00398 rAElem(index1,index2) += dSdF_quad_quad(node_index1,spatial_dim1,node_index2,spatial_dim2)
00399 * wJ;
00400 }
00401 }
00402 }
00403 }
00404
00405 rAElemPrecond.clear();
00406 if (assembleJacobian)
00407 {
00408
00409
00410
00411
00412
00413
00414
00415
00416
00417
00418
00419 rAElemPrecond = rAElem;
00420 }
00421 }
00422
00423 template<size_t DIM>
00424 void CompressibleNonlinearElasticitySolver<DIM>::AssembleOnBoundaryElement(
00425 BoundaryElement<DIM-1,DIM>& rBoundaryElement,
00426 c_matrix<double,BOUNDARY_STENCIL_SIZE,BOUNDARY_STENCIL_SIZE>& rAelem,
00427 c_vector<double,BOUNDARY_STENCIL_SIZE>& rBelem,
00428 bool assembleResidual,
00429 bool assembleJacobian,
00430 unsigned boundaryConditionIndex)
00431 {
00432 rAelem.clear();
00433 rBelem.clear();
00434
00435 if (assembleJacobian && !assembleResidual)
00436 {
00437
00438 return;
00439 }
00440
00441 c_vector<double, DIM> weighted_direction;
00442 double jacobian_determinant;
00443
00444 this->mrQuadMesh.GetWeightedDirectionForBoundaryElement(rBoundaryElement.GetIndex(), weighted_direction, jacobian_determinant);
00445
00446
00447
00448
00449
00450 c_vector<double,DIM> deformed_normal;
00451 if (this->mrProblemDefinition.GetTractionBoundaryConditionType()==PRESSURE_ON_DEFORMED)
00452 {
00453
00454 static std::vector<c_vector<double,DIM> > element_current_displacements(DIM);
00455 for (unsigned II=0; II<DIM; II++)
00456 {
00457 for (unsigned JJ=0; JJ<DIM; JJ++)
00458 {
00459 element_current_displacements[II](JJ) = this->mCurrentSolution[DIM*rBoundaryElement.GetNodeGlobalIndex(II) + JJ];
00460 }
00461 }
00462
00463 this->mDeformedBoundaryElement.ApplyUndeformedElementAndDisplacement(&rBoundaryElement, element_current_displacements);
00464
00465 this->mDeformedBoundaryElement.CalculateWeightedDirection(weighted_direction, jacobian_determinant);
00466
00467 deformed_normal = this->mDeformedBoundaryElement.CalculateNormal();
00468 }
00469
00470 c_vector<double,NUM_NODES_PER_BOUNDARY_ELEMENT> phi;
00471
00472 for (unsigned quad_index=0; quad_index<this->mpBoundaryQuadratureRule->GetNumQuadPoints(); quad_index++)
00473 {
00474 double wJ = jacobian_determinant * this->mpBoundaryQuadratureRule->GetWeight(quad_index);
00475
00476 const ChastePoint<DIM-1>& quad_point = this->mpBoundaryQuadratureRule->rGetQuadPoint(quad_index);
00477
00478 QuadraticBasisFunction<DIM-1>::ComputeBasisFunctions(quad_point, phi);
00479
00480
00481
00482 c_vector<double,DIM> traction = zero_vector<double>(DIM);
00483 switch (this->mrProblemDefinition.GetTractionBoundaryConditionType())
00484 {
00485 case FUNCTIONAL_TRACTION:
00486 {
00487 c_vector<double,DIM> X = zero_vector<double>(DIM);
00488 for (unsigned node_index=0; node_index<NUM_NODES_PER_BOUNDARY_ELEMENT; node_index++)
00489 {
00490 X += phi(node_index)*this->mrQuadMesh.GetNode( rBoundaryElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
00491 }
00492 traction = this->mrProblemDefinition.EvaluateTractionFunction(X, this->mCurrentTime);
00493 break;
00494 }
00495 case ELEMENTWISE_TRACTION:
00496 {
00497 traction = this->mrProblemDefinition.rGetElementwiseTractions()[boundaryConditionIndex];
00498 break;
00499 }
00500 case PRESSURE_ON_DEFORMED:
00501 {
00502
00503 traction = this->mrProblemDefinition.GetNormalPressure()*deformed_normal;
00504 break;
00505 }
00506 default:
00507 NEVER_REACHED;
00508 }
00509
00510
00511 for (unsigned index=0; index<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index++)
00512 {
00513 unsigned spatial_dim = index%DIM;
00514 unsigned node_index = (index-spatial_dim)/DIM;
00515
00516 assert(node_index < NUM_NODES_PER_BOUNDARY_ELEMENT);
00517
00518 rBelem(index) -= traction(spatial_dim)
00519 * phi(node_index)
00520 * wJ;
00521 }
00522 }
00523 }
00524
00525 template<size_t DIM>
00526 CompressibleNonlinearElasticitySolver<DIM>::CompressibleNonlinearElasticitySolver(QuadraticMesh<DIM>& rQuadMesh,
00527 SolidMechanicsProblemDefinition<DIM>& rProblemDefinition,
00528 std::string outputDirectory)
00529 : AbstractNonlinearElasticitySolver<DIM>(rQuadMesh,
00530 rProblemDefinition,
00531 outputDirectory,
00532 COMPRESSIBLE)
00533 {
00534 if(rProblemDefinition.GetCompressibilityType() != COMPRESSIBLE)
00535 {
00536 EXCEPTION("SolidMechanicsProblemDefinition object contains incompressible material laws");
00537 }
00538
00539 this->Initialise();
00540 }
00541
00542
00543 template<size_t DIM>
00544 CompressibleNonlinearElasticitySolver<DIM>::~CompressibleNonlinearElasticitySolver()
00545 {
00546 }
00547
00549
00551
00552 template class CompressibleNonlinearElasticitySolver<2>;
00553 template class CompressibleNonlinearElasticitySolver<3>;