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
00042
00043 #include "CompressibleNonlinearElasticitySolver.hpp"
00044 #include "LinearBasisFunction.hpp"
00045 #include "QuadraticBasisFunction.hpp"
00046 #include <algorithm>
00047
00048 template<size_t DIM>
00049 void CompressibleNonlinearElasticitySolver<DIM>::AssembleSystem(bool assembleResidual,
00050 bool assembleJacobian)
00051 {
00052
00053 assert(assembleResidual || assembleJacobian);
00054 assert(this->mCurrentSolution.size()==this->mNumDofs);
00055
00056
00057 if (assembleResidual)
00058 {
00059 this->mpLinearSystem->ZeroRhsVector();
00060 }
00061 if (assembleJacobian)
00062 {
00063 this->mpLinearSystem->ZeroLhsMatrix();
00064 this->mpPreconditionMatrixLinearSystem->ZeroLhsMatrix();
00065 }
00066
00067 c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem;
00068
00069
00070 c_matrix<double, STENCIL_SIZE, STENCIL_SIZE> a_elem_precond;
00071 c_vector<double, STENCIL_SIZE> b_elem;
00072
00074
00076 for (typename AbstractTetrahedralMesh<DIM, DIM>::ElementIterator iter = this->mpQuadMesh->GetElementIteratorBegin();
00077 iter != this->mpQuadMesh->GetElementIteratorEnd();
00078 ++iter)
00079 {
00080 #ifdef MECH_VERY_VERBOSE
00081 if (assembleJacobian)
00082 {
00083 std::cout << "\r[" << PetscTools::GetMyRank() << "]: Element " << (*iter).GetIndex() << " of " << this->mpQuadMesh->GetNumElements() << std::flush;
00084 }
00085 #endif
00086
00087 Element<DIM, DIM>& element = *iter;
00088
00089 if (element.GetOwnership() == true)
00090 {
00091 AssembleOnElement(element, a_elem, a_elem_precond, b_elem, assembleResidual, assembleJacobian);
00092
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105 unsigned p_indices[STENCIL_SIZE];
00106 for (unsigned i=0; i<NUM_NODES_PER_ELEMENT; i++)
00107 {
00108 for (unsigned j=0; j<DIM; j++)
00109 {
00110 p_indices[DIM*i+j] = DIM*element.GetNodeGlobalIndex(i) + j;
00111 }
00112 }
00113
00114 if (assembleJacobian)
00115 {
00116 this->mpLinearSystem->AddLhsMultipleValues(p_indices, a_elem);
00117 this->mpPreconditionMatrixLinearSystem->AddLhsMultipleValues(p_indices, a_elem_precond);
00118 }
00119
00120 if (assembleResidual)
00121 {
00122 this->mpLinearSystem->AddRhsMultipleValues(p_indices, b_elem);
00123 }
00124 }
00125 }
00126
00128
00129
00131 c_vector<double, BOUNDARY_STENCIL_SIZE> b_boundary_elem;
00132 c_matrix<double, BOUNDARY_STENCIL_SIZE, BOUNDARY_STENCIL_SIZE> a_boundary_elem;
00133 if (this->mBoundaryElements.size()>0)
00134 {
00135 for (unsigned i=0; i<this->mBoundaryElements.size(); i++)
00136 {
00137 BoundaryElement<DIM-1,DIM>& r_boundary_element = *(this->mBoundaryElements[i]);
00138 AssembleOnBoundaryElement(r_boundary_element, a_boundary_elem, b_boundary_elem, this->mSurfaceTractions[i], assembleResidual, assembleJacobian);
00139
00140 unsigned p_indices[BOUNDARY_STENCIL_SIZE];
00141 for (unsigned i=0; i<NUM_NODES_PER_BOUNDARY_ELEMENT; i++)
00142 {
00143 for (unsigned j=0; j<DIM; j++)
00144 {
00145 p_indices[DIM*i+j] = DIM*r_boundary_element.GetNodeGlobalIndex(i) + j;
00146 }
00147 }
00148
00149 if (assembleJacobian)
00150 {
00151 this->mpLinearSystem->AddLhsMultipleValues(p_indices, a_boundary_elem);
00152 this->mpPreconditionMatrixLinearSystem->AddLhsMultipleValues(p_indices, a_boundary_elem);
00153 }
00154
00155 if (assembleResidual)
00156 {
00157 this->mpLinearSystem->AddRhsMultipleValues(p_indices, b_boundary_elem);
00158 }
00159 }
00160 }
00161
00162 if (assembleResidual)
00163 {
00164 this->mpLinearSystem->AssembleRhsVector();
00165 }
00166 if (assembleJacobian)
00167 {
00168 this->mpLinearSystem->AssembleIntermediateLhsMatrix();
00169 this->mpPreconditionMatrixLinearSystem->AssembleIntermediateLhsMatrix();
00170 }
00171
00172
00173 this->ApplyBoundaryConditions(assembleJacobian);
00174
00175 if (assembleResidual)
00176 {
00177 this->mpLinearSystem->AssembleRhsVector();
00178 }
00179 if (assembleJacobian)
00180 {
00181 this->mpLinearSystem->AssembleFinalLhsMatrix();
00182 this->mpPreconditionMatrixLinearSystem->AssembleFinalLhsMatrix();
00183 }
00184 }
00185
00186 template<size_t DIM>
00187 void CompressibleNonlinearElasticitySolver<DIM>::AssembleOnElement(
00188 Element<DIM, DIM>& rElement,
00189 c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElem,
00190 c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElemPrecond,
00191 c_vector<double, STENCIL_SIZE>& rBElem,
00192 bool assembleResidual,
00193 bool assembleJacobian)
00194 {
00195 static c_matrix<double,DIM,DIM> jacobian;
00196 static c_matrix<double,DIM,DIM> inverse_jacobian;
00197 double jacobian_determinant;
00198
00199 this->mpQuadMesh->GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
00200
00201 if (assembleJacobian)
00202 {
00203 rAElem.clear();
00204 rAElemPrecond.clear();
00205 }
00206
00207 if (assembleResidual)
00208 {
00209 rBElem.clear();
00210 }
00211
00213
00215 static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> element_current_displacements;
00216 for (unsigned II=0; II<NUM_NODES_PER_ELEMENT; II++)
00217 {
00218 for (unsigned JJ=0; JJ<DIM; JJ++)
00219 {
00220 element_current_displacements(JJ,II) = this->mCurrentSolution[DIM*rElement.GetNodeGlobalIndex(II) + JJ];
00221 }
00222 }
00223
00224
00225 static c_vector<double, NUM_VERTICES_PER_ELEMENT> linear_phi;
00226 static c_vector<double, NUM_NODES_PER_ELEMENT> quad_phi;
00227 static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> grad_quad_phi;
00228 static c_matrix<double, NUM_NODES_PER_ELEMENT, DIM> trans_grad_quad_phi;
00229
00230
00231
00232 AbstractMaterialLaw<DIM>* p_material_law;
00233 if (this->mMaterialLaws.size()==1)
00234 {
00235
00236 p_material_law = this->mMaterialLaws[0];
00237 }
00238 else
00239 {
00240
00241 #define COVERAGE_IGNORE // not going to have tests in cts for everything
00242 p_material_law = this->mMaterialLaws[rElement.GetIndex()];
00243 #undef COVERAGE_IGNORE
00244 }
00245
00246
00247 static c_matrix<double,DIM,DIM> grad_u;
00248
00249 static c_matrix<double,DIM,DIM> F;
00250 static c_matrix<double,DIM,DIM> C;
00251 static c_matrix<double,DIM,DIM> inv_C;
00252 static c_matrix<double,DIM,DIM> inv_F;
00253 static c_matrix<double,DIM,DIM> T;
00254
00255 static c_matrix<double,DIM,DIM> F_T;
00256 static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> F_T_grad_quad_phi;
00257
00258 c_vector<double,DIM> body_force;
00259
00260 static FourthOrderTensor<DIM,DIM,DIM,DIM> dTdE;
00261 static FourthOrderTensor<DIM,DIM,DIM,DIM> dSdF;
00262
00263 static FourthOrderTensor<NUM_NODES_PER_ELEMENT,DIM,DIM,DIM> temp_tensor;
00264 static FourthOrderTensor<NUM_NODES_PER_ELEMENT,DIM,NUM_NODES_PER_ELEMENT,DIM> dSdF_quad_quad;
00265
00266 static c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> temp_matrix;
00267 static c_matrix<double,NUM_NODES_PER_ELEMENT,DIM> grad_quad_phi_times_invF;
00268
00269
00270
00276 for (unsigned quadrature_index=0; quadrature_index < this->mpQuadratureRule->GetNumQuadPoints(); quadrature_index++)
00277 {
00278
00279 unsigned current_quad_point_global_index = rElement.GetIndex()*this->mpQuadratureRule->GetNumQuadPoints()
00280 + quadrature_index;
00281
00282 double wJ = jacobian_determinant * this->mpQuadratureRule->GetWeight(quadrature_index);
00283
00284 const ChastePoint<DIM>& quadrature_point = this->mpQuadratureRule->rGetQuadPoint(quadrature_index);
00285
00287
00289 LinearBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, linear_phi);
00290 QuadraticBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, quad_phi);
00291 QuadraticBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(quadrature_point, inverse_jacobian, grad_quad_phi);
00292 trans_grad_quad_phi = trans(grad_quad_phi);
00293
00294
00296
00298
00299 if(assembleResidual)
00300 {
00301 if (this->mUsingBodyForceFunction)
00302 {
00303 c_vector<double,DIM> X = zero_vector<double>(DIM);
00304
00305 for (unsigned node_index=0; node_index<NUM_VERTICES_PER_ELEMENT; node_index++)
00306 {
00307 X += linear_phi(node_index)*this->mpQuadMesh->GetNode( rElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
00308 }
00309 body_force = (*(this->mpBodyForceFunction))(X, this->mCurrentTime);
00310 }
00311 else
00312 {
00313 body_force = this->mBodyForce;
00314 }
00315 }
00316
00318
00320 grad_u = zero_matrix<double>(DIM,DIM);
00321
00322 for (unsigned node_index=0; node_index<NUM_NODES_PER_ELEMENT; node_index++)
00323 {
00324 for (unsigned i=0; i<DIM; i++)
00325 {
00326 for (unsigned M=0; M<DIM; M++)
00327 {
00328 grad_u(i,M) += grad_quad_phi(M,node_index)*element_current_displacements(i,node_index);
00329 }
00330 }
00331 }
00332
00333
00335
00337 for (unsigned i=0; i<DIM; i++)
00338 {
00339 for (unsigned M=0; M<DIM; M++)
00340 {
00341 F(i,M) = (i==M?1:0) + grad_u(i,M);
00342 }
00343 }
00344
00345 C = prod(trans(F),F);
00346 inv_C = Inverse(C);
00347 inv_F = Inverse(F);
00348
00349 ComputeStressAndStressDerivative(p_material_law, C, inv_C, 0.0, rElement.GetIndex(), current_quad_point_global_index,
00350 T, dTdE, assembleJacobian);
00351
00352
00354
00356 if (assembleResidual)
00357 {
00358 F_T = prod(F,T);
00359 F_T_grad_quad_phi = prod(F_T, grad_quad_phi);
00360
00361 for (unsigned index=0; index<NUM_NODES_PER_ELEMENT*DIM; index++)
00362 {
00363 unsigned spatial_dim = index%DIM;
00364 unsigned node_index = (index-spatial_dim)/DIM;
00365
00366 rBElem(index) += - this->mDensity
00367 * body_force(spatial_dim)
00368 * quad_phi(node_index)
00369 * wJ;
00370
00371
00372 rBElem(index) += F_T_grad_quad_phi(spatial_dim,node_index)
00373 * wJ;
00374 }
00375 }
00376
00377
00379
00381 if (assembleJacobian)
00382 {
00383
00384 grad_quad_phi_times_invF = prod(trans_grad_quad_phi, inv_F);
00385
00387
00388
00389
00390
00391
00392
00393
00394
00395
00397
00398
00399 for (unsigned M=0; M<DIM; M++)
00400 {
00401 for (unsigned N=0; N<DIM; N++)
00402 {
00403 for (unsigned P=0; P<DIM; P++)
00404 {
00405 for (unsigned Q=0; Q<DIM; Q++)
00406 {
00407
00408 dSdF(M,N,P,Q) = 0.5*(dTdE(M,N,P,Q) + dTdE(M,N,Q,P));
00409 }
00410 }
00411 }
00412 }
00413
00414 dTdE.template SetAsContractionOnSecondDimension<DIM>(F, dSdF);
00415
00416 dSdF.template SetAsContractionOnFourthDimension<DIM>(F, dTdE);
00417
00418
00419 for(unsigned M=0; M<DIM; M++)
00420 {
00421 for(unsigned N=0; N<DIM; N++)
00422 {
00423 for(unsigned i=0; i<DIM; i++)
00424 {
00425 dSdF(M,i,N,i) += T(M,N);
00426 }
00427 }
00428 }
00429
00430
00432
00433
00434
00435
00436
00437
00438
00439
00440
00441
00443 temp_tensor.template SetAsContractionOnFirstDimension<DIM>( trans_grad_quad_phi, dSdF );
00444 dSdF_quad_quad.template SetAsContractionOnThirdDimension<DIM>( trans_grad_quad_phi, temp_tensor );
00445
00446
00447
00448
00449 for (unsigned index1=0; index1<NUM_NODES_PER_ELEMENT*DIM; index1++)
00450 {
00451 unsigned spatial_dim1 = index1%DIM;
00452 unsigned node_index1 = (index1-spatial_dim1)/DIM;
00453
00454
00455 for (unsigned index2=0; index2<NUM_NODES_PER_ELEMENT*DIM; index2++)
00456 {
00457 unsigned spatial_dim2 = index2%DIM;
00458 unsigned node_index2 = (index2-spatial_dim2)/DIM;
00459
00460
00461 rAElem(index1,index2) += dSdF_quad_quad(node_index1,spatial_dim1,node_index2,spatial_dim2)
00462 * wJ;
00463 }
00464 }
00465 }
00466 }
00467
00468
00469 if (assembleJacobian)
00470 {
00471 rAElemPrecond = rAElem;
00472 }
00473 }
00474
00475
00476 template<size_t DIM>
00477 void CompressibleNonlinearElasticitySolver<DIM>::ComputeStressAndStressDerivative(AbstractMaterialLaw<DIM>* pMaterialLaw,
00478 c_matrix<double,DIM,DIM>& rC,
00479 c_matrix<double,DIM,DIM>& rInvC,
00480 double pressure,
00481 unsigned elementIndex,
00482 unsigned currentQuadPointGlobalIndex,
00483 c_matrix<double,DIM,DIM>& rT,
00484 FourthOrderTensor<DIM,DIM,DIM,DIM>& rDTdE,
00485 bool computeDTdE)
00486 {
00487
00488 pMaterialLaw->ComputeStressAndStressDerivative(rC,rInvC,pressure,rT,rDTdE,computeDTdE);
00489 }
00490
00491
00492 template<size_t DIM>
00493 void CompressibleNonlinearElasticitySolver<DIM>::AssembleOnBoundaryElement(
00494 BoundaryElement<DIM-1,DIM>& rBoundaryElement,
00495 c_matrix<double,BOUNDARY_STENCIL_SIZE,BOUNDARY_STENCIL_SIZE>& rAelem,
00496 c_vector<double,BOUNDARY_STENCIL_SIZE>& rBelem,
00497 c_vector<double,DIM>& rTraction,
00498 bool assembleResidual,
00499 bool assembleJacobian)
00500 {
00501 rAelem.clear();
00502 rBelem.clear();
00503
00504 if (assembleJacobian && !assembleResidual)
00505 {
00506
00507 return;
00508 }
00509
00510 c_vector<double, DIM> weighted_direction;
00511 double jacobian_determinant;
00512 this->mpQuadMesh->GetWeightedDirectionForBoundaryElement(rBoundaryElement.GetIndex(), weighted_direction, jacobian_determinant);
00513
00514 c_vector<double,NUM_NODES_PER_BOUNDARY_ELEMENT> phi;
00515
00516 for (unsigned quad_index=0; quad_index<this->mpBoundaryQuadratureRule->GetNumQuadPoints(); quad_index++)
00517 {
00518 double wJ = jacobian_determinant * this->mpBoundaryQuadratureRule->GetWeight(quad_index);
00519
00520 const ChastePoint<DIM-1>& quad_point = this->mpBoundaryQuadratureRule->rGetQuadPoint(quad_index);
00521
00522 QuadraticBasisFunction<DIM-1>::ComputeBasisFunctions(quad_point, phi);
00523
00524
00525
00526 c_vector<double,DIM> traction = zero_vector<double>(DIM);
00527 if (this->mUsingTractionBoundaryConditionFunction)
00528 {
00529 c_vector<double,DIM> X = zero_vector<double>(DIM);
00530 for (unsigned node_index=0; node_index<NUM_NODES_PER_BOUNDARY_ELEMENT; node_index++)
00531 {
00532 X += phi(node_index)*this->mpQuadMesh->GetNode( rBoundaryElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
00533 }
00534 traction = (*(this->mpTractionBoundaryConditionFunction))(X, this->mCurrentTime);
00535 }
00536 else
00537 {
00538 traction = rTraction;
00539 }
00540
00541
00542 for (unsigned index=0; index<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index++)
00543 {
00544 unsigned spatial_dim = index%DIM;
00545 unsigned node_index = (index-spatial_dim)/DIM;
00546
00547 assert(node_index < NUM_NODES_PER_BOUNDARY_ELEMENT);
00548
00549 rBelem(index) -= traction(spatial_dim)
00550 * phi(node_index)
00551 * wJ;
00552 }
00553 }
00554 }
00555
00556
00557
00558 template<size_t DIM>
00559 CompressibleNonlinearElasticitySolver<DIM>::CompressibleNonlinearElasticitySolver(
00560 QuadraticMesh<DIM>* pQuadMesh,
00561 AbstractMaterialLaw<DIM>* pMaterialLaw,
00562 c_vector<double,DIM> bodyForce,
00563 double density,
00564 std::string outputDirectory,
00565 std::vector<unsigned>& fixedNodes,
00566 std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00567 : AbstractNonlinearElasticitySolver<COMPRESSIBLE,DIM>(pQuadMesh,
00568 bodyForce, density,
00569 outputDirectory, fixedNodes)
00570 {
00571 assert(pMaterialLaw != NULL);
00572 mMaterialLaws.push_back(pMaterialLaw);
00573
00574 Initialise(pFixedNodeLocations);
00575 }
00576
00577
00578 template<size_t DIM>
00579 CompressibleNonlinearElasticitySolver<DIM>::CompressibleNonlinearElasticitySolver(
00580 QuadraticMesh<DIM>* pQuadMesh,
00581 std::vector<AbstractMaterialLaw<DIM>*>& rMaterialLaws,
00582 c_vector<double,DIM> bodyForce,
00583 double density,
00584 std::string outputDirectory,
00585 std::vector<unsigned>& fixedNodes,
00586 std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00587 : AbstractNonlinearElasticitySolver<COMPRESSIBLE,DIM>(pQuadMesh,
00588 bodyForce, density,
00589 outputDirectory, fixedNodes)
00590 {
00591 mMaterialLaws.resize(rMaterialLaws.size(), NULL);
00592 for (unsigned i=0; i<mMaterialLaws.size(); i++)
00593 {
00594 assert(rMaterialLaws[i] != NULL);
00595 mMaterialLaws[i] = rMaterialLaws[i];
00596 }
00597
00598 assert(rMaterialLaws.size()==pQuadMesh->GetNumElements());
00599 Initialise(pFixedNodeLocations);
00600 }
00601
00602
00603 template<size_t DIM>
00604 CompressibleNonlinearElasticitySolver<DIM>::~CompressibleNonlinearElasticitySolver()
00605 {
00606 }
00608
00610
00611
00612 template class CompressibleNonlinearElasticitySolver<2>;
00613 template class CompressibleNonlinearElasticitySolver<3>;