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 "NonlinearElasticityAssembler.hpp"
00042 #include "LinearBasisFunction.hpp"
00043 #include "QuadraticBasisFunction.hpp"
00044
00045
00046 template<size_t DIM>
00047 void NonlinearElasticityAssembler<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 this->mpLinearSystem->ZeroRhsVector();
00058 }
00059 if (assembleJacobian)
00060 {
00061 this->mpLinearSystem->ZeroLhsMatrix();
00062 this->mpPreconditionMatrixLinearSystem->ZeroLhsMatrix();
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
00072
00074 for (typename AbstractTetrahedralMesh<DIM, DIM>::ElementIterator iter = mpQuadMesh->GetElementIteratorBegin();
00075 iter != mpQuadMesh->GetElementIteratorEnd();
00076 ++iter)
00077 {
00078
00079
00080
00081
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
00089 unsigned p_indices[STENCIL_SIZE];
00090 for (unsigned i=0; i<NUM_NODES_PER_ELEMENT; i++)
00091 {
00092 for (unsigned j=0; j<DIM; j++)
00093 {
00094 p_indices[DIM*i+j] = DIM*element.GetNodeGlobalIndex(i) + j;
00095 }
00096 }
00097
00098 for (unsigned i=0; i<NUM_VERTICES_PER_ELEMENT; i++)
00099 {
00100 p_indices[DIM*NUM_NODES_PER_ELEMENT + i] = DIM*mpQuadMesh->GetNumNodes() + element.GetNodeGlobalIndex(i);
00101 }
00102
00103 if (assembleJacobian)
00104 {
00105 this->mpLinearSystem->AddLhsMultipleValues(p_indices, a_elem);
00106 this->mpPreconditionMatrixLinearSystem->AddLhsMultipleValues(p_indices, a_elem_precond);
00107 }
00108
00109 if (assembleResidual)
00110 {
00111 this->mpLinearSystem->AddRhsMultipleValues(p_indices, b_elem);
00112 }
00113 }
00114 }
00115
00117
00118
00120 c_vector<double, BOUNDARY_STENCIL_SIZE> b_boundary_elem;
00121 c_matrix<double, BOUNDARY_STENCIL_SIZE, BOUNDARY_STENCIL_SIZE> a_boundary_elem;
00122 if (mBoundaryElements.size()>0)
00123 {
00124 for (unsigned i=0; i<mBoundaryElements.size(); i++)
00125 {
00126 BoundaryElement<DIM-1,DIM>& r_boundary_element = *(mBoundaryElements[i]);
00127 AssembleOnBoundaryElement(r_boundary_element, a_boundary_elem, b_boundary_elem, this->mSurfaceTractions[i], assembleResidual, assembleJacobian);
00128
00129 unsigned p_indices[BOUNDARY_STENCIL_SIZE];
00130 for (unsigned i=0; i<NUM_NODES_PER_BOUNDARY_ELEMENT; i++)
00131 {
00132 for (unsigned j=0; j<DIM; j++)
00133 {
00134 p_indices[DIM*i+j] = DIM*r_boundary_element.GetNodeGlobalIndex(i) + j;
00135 }
00136 }
00137
00138 for (unsigned i=0; i<DIM ; i++)
00139 {
00140 p_indices[DIM*NUM_NODES_PER_BOUNDARY_ELEMENT + i] = DIM*mpQuadMesh->GetNumNodes() + r_boundary_element.GetNodeGlobalIndex(i);
00141 }
00142
00143 if (assembleJacobian)
00144 {
00145 this->mpLinearSystem->AddLhsMultipleValues(p_indices, a_boundary_elem);
00146 this->mpPreconditionMatrixLinearSystem->AddLhsMultipleValues(p_indices, a_boundary_elem);
00147 }
00148
00149 if (assembleResidual)
00150 {
00151 this->mpLinearSystem->AddRhsMultipleValues(p_indices, b_boundary_elem);
00152 }
00153
00154
00155 if (DIM==2)
00156 {
00157 assert(8==BOUNDARY_STENCIL_SIZE);
00158 assert(b_boundary_elem(6)==0);
00159 assert(b_boundary_elem(7)==0);
00160 }
00161 }
00162 }
00163
00164 if (assembleResidual)
00165 {
00166 this->mpLinearSystem->AssembleRhsVector();
00167 }
00168 if (assembleJacobian)
00169 {
00170 this->mpLinearSystem->AssembleIntermediateLhsMatrix();
00171 this->mpPreconditionMatrixLinearSystem->AssembleIntermediateLhsMatrix();
00172 }
00173
00174
00175
00176 this->ApplyBoundaryConditions(assembleJacobian);
00177
00178 if (assembleResidual)
00179 {
00180 this->mpLinearSystem->AssembleRhsVector();
00181 }
00182 if (assembleJacobian)
00183 {
00184 this->mpLinearSystem->AssembleFinalLhsMatrix();
00185 this->mpPreconditionMatrixLinearSystem->AssembleFinalLhsMatrix();
00186 }
00187 }
00188
00189 template<size_t DIM>
00190 void NonlinearElasticityAssembler<DIM>::AssembleOnElement(
00191 Element<DIM, DIM>& rElement,
00192 c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElem,
00193 c_matrix<double, STENCIL_SIZE, STENCIL_SIZE >& rAElemPrecond,
00194 c_vector<double, STENCIL_SIZE>& rBElem,
00195 bool assembleResidual,
00196 bool assembleJacobian)
00197 {
00198 c_matrix<double, DIM, DIM> jacobian, inverse_jacobian;
00199 double jacobian_determinant;
00200
00201 mpQuadMesh->GetInverseJacobianForElement(rElement.GetIndex(), jacobian, jacobian_determinant, inverse_jacobian);
00202
00203 if (assembleJacobian)
00204 {
00205 rAElem.clear();
00206 rAElemPrecond.clear();
00207 }
00208
00209 if (assembleResidual)
00210 {
00211 rBElem.clear();
00212 }
00213
00215
00217 static c_matrix<double,DIM,NUM_NODES_PER_ELEMENT> element_current_displacements;
00218 static c_vector<double,NUM_VERTICES_PER_ELEMENT> element_current_pressures;
00219 for (unsigned II=0; II<NUM_NODES_PER_ELEMENT; II++)
00220 {
00221 for (unsigned JJ=0; JJ<DIM; JJ++)
00222 {
00223 element_current_displacements(JJ,II) = this->mCurrentSolution[DIM*rElement.GetNodeGlobalIndex(II) + JJ];
00224 }
00225 }
00226
00228
00230 for (unsigned II=0; II<NUM_VERTICES_PER_ELEMENT; II++)
00231 {
00232 element_current_pressures(II) = this->mCurrentSolution[DIM*mpQuadMesh->GetNumNodes() + rElement.GetNodeGlobalIndex(II)];
00233 }
00234
00235
00236 c_vector<double, NUM_VERTICES_PER_ELEMENT> linear_phi;
00237 c_vector<double, NUM_NODES_PER_ELEMENT> quad_phi;
00238 c_matrix<double, DIM, NUM_NODES_PER_ELEMENT> grad_quad_phi;
00239
00240
00241 AbstractIncompressibleMaterialLaw<DIM> *p_material_law;
00242 if (this->mMaterialLaws.size()==1)
00243 {
00244
00245 p_material_law = this->mMaterialLaws[0];
00246 }
00247 else
00248 {
00249
00250 #define COVERAGE_IGNORE // not going to have tests in cts for everything
00251 p_material_law = this->mMaterialLaws[rElement.GetIndex()];
00252 #undef COVERAGE_IGNORE
00253 }
00254
00255
00261 for (unsigned quadrature_index=0; quadrature_index < this->mpQuadratureRule->GetNumQuadPoints(); quadrature_index++)
00262 {
00263 double wJ = jacobian_determinant * this->mpQuadratureRule->GetWeight(quadrature_index);
00264
00265 const ChastePoint<DIM>& quadrature_point = this->mpQuadratureRule->rGetQuadPoint(quadrature_index);
00266
00268
00270 LinearBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, linear_phi);
00271 QuadraticBasisFunction<DIM>::ComputeBasisFunctions(quadrature_point, quad_phi);
00272 QuadraticBasisFunction<DIM>::ComputeTransformedBasisFunctionDerivatives(quadrature_point, inverse_jacobian, grad_quad_phi);
00273
00274
00276
00278 c_vector<double,DIM> body_force;
00279
00280 if (this->mUsingBodyForceFunction)
00281 {
00282 c_vector<double,DIM> X = zero_vector<double>(DIM);
00283
00284 for (unsigned node_index=0; node_index<NUM_VERTICES_PER_ELEMENT; node_index++)
00285 {
00286 X += linear_phi(node_index)*mpQuadMesh->GetNode( rElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
00287 }
00288 body_force = (*(this->mpBodyForceFunction))(X);
00289 }
00290 else
00291 {
00292 body_force = this->mBodyForce;
00293 }
00294
00296
00298 static c_matrix<double,DIM,DIM> grad_u;
00299 grad_u = zero_matrix<double>(DIM,DIM);
00300
00301 for (unsigned node_index=0; node_index<NUM_NODES_PER_ELEMENT; node_index++)
00302 {
00303 for (unsigned i=0; i<DIM; i++)
00304 {
00305 for (unsigned M=0; M<DIM; M++)
00306 {
00307 grad_u(i,M) += grad_quad_phi(M,node_index)*element_current_displacements(i,node_index);
00308 }
00309 }
00310 }
00311
00312 double pressure = 0;
00313 for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00314 {
00315 pressure += linear_phi(vertex_index)*element_current_pressures(vertex_index);
00316 }
00317
00318
00320
00322 static c_matrix<double,DIM,DIM> F;
00323 static c_matrix<double,DIM,DIM> C;
00324 static c_matrix<double,DIM,DIM> inv_C;
00325 static c_matrix<double,DIM,DIM> inv_F;
00326 static c_matrix<double,DIM,DIM> T;
00327
00328 for (unsigned i=0; i<DIM; i++)
00329 {
00330 for (unsigned M=0; M<DIM; M++)
00331 {
00332 F(i,M) = (i==M?1:0) + grad_u(i,M);
00333 }
00334 }
00335
00336 C = prod(trans(F),F);
00337 inv_C = Inverse(C);
00338 inv_F = Inverse(F);
00339
00340 double detF = Determinant(F);
00341
00342 p_material_law->ComputeStressAndStressDerivative(C,inv_C,pressure,T,this->dTdE,assembleJacobian);
00343
00344 static FourthOrderTensor<DIM> dTdE_F;
00345 static FourthOrderTensor<DIM> dTdE_FF1;
00346 static FourthOrderTensor<DIM> dTdE_FF2;
00347
00348 dTdE_F.SetAsProduct(this->dTdE, F, 0);
00349 dTdE_FF1.SetAsProduct(dTdE_F, F, 3);
00350 dTdE_FF2.SetAsProduct(dTdE_F, F, 2);
00351
00352
00354
00356 if (assembleResidual)
00357 {
00358 for (unsigned index=0; index<NUM_NODES_PER_ELEMENT*DIM; index++)
00359 {
00360 unsigned spatial_dim = index%DIM;
00361 unsigned node_index = (index-spatial_dim)/DIM;
00362
00363 assert(node_index < NUM_NODES_PER_ELEMENT);
00364
00365 rBElem(index) += - this->mDensity
00366 * body_force(spatial_dim)
00367 * quad_phi(node_index)
00368 * wJ;
00369
00370 for (unsigned M=0; M<DIM; M++)
00371 {
00372 for (unsigned N=0; N<DIM; N++)
00373 {
00374 rBElem(index) += T(M,N)
00375 * F(spatial_dim,M)
00376 * grad_quad_phi(N,node_index)
00377 * wJ;
00378 }
00379 }
00380 }
00381
00382 for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00383 {
00384 rBElem( NUM_NODES_PER_ELEMENT*DIM + vertex_index ) += linear_phi(vertex_index)
00385 * (detF - 1)
00386 * wJ;
00387 }
00388 }
00389
00391
00393 if (assembleJacobian)
00394 {
00395 for (unsigned index1=0; index1<NUM_NODES_PER_ELEMENT*DIM; index1++)
00396 {
00397 unsigned spatial_dim1 = index1%DIM;
00398 unsigned node_index1 = (index1-spatial_dim1)/DIM;
00399
00400
00401 for (unsigned index2=0; index2<NUM_NODES_PER_ELEMENT*DIM; index2++)
00402 {
00403 unsigned spatial_dim2 = index2%DIM;
00404 unsigned node_index2 = (index2-spatial_dim2)/DIM;
00405
00406 for (unsigned M=0; M<DIM; M++)
00407 {
00408 for (unsigned N=0; N<DIM; N++)
00409 {
00410 rAElem(index1,index2) += T(M,N)
00411 * grad_quad_phi(N,node_index1)
00412 * grad_quad_phi(M,node_index2)
00413 * (spatial_dim1==spatial_dim2?1:0)
00414 * wJ;
00415
00417
00418
00419
00420
00421
00422
00423
00424
00425
00426
00427
00428
00429
00430
00431
00432
00433
00434
00435 }
00436 }
00437
00438 for (unsigned N=0; N<DIM; N++)
00439 {
00440 for (unsigned P=0; P<DIM; P++)
00441 {
00442 rAElem(index1,index2) += 0.5
00443 * dTdE_FF1(spatial_dim1,N,P,spatial_dim2)
00444 * grad_quad_phi(P,node_index2)
00445 * grad_quad_phi(N,node_index1)
00446 * wJ;
00447 }
00448
00449 for (unsigned Q=0; Q<DIM; Q++)
00450 {
00451 rAElem(index1,index2) += 0.5
00452 * dTdE_FF2(spatial_dim1,N,spatial_dim2,Q)
00453 * grad_quad_phi(Q,node_index2)
00454 * grad_quad_phi(N,node_index1)
00455 * wJ;
00456 }
00457 }
00458 }
00459
00460 for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00461 {
00462 unsigned index2 = NUM_NODES_PER_ELEMENT*DIM + vertex_index;
00463
00464 for (unsigned M=0; M<DIM; M++)
00465 {
00466 for (unsigned N=0; N<DIM; N++)
00467 {
00468 rAElem(index1,index2) += - F(spatial_dim1,M)
00469 * inv_C(M,N)
00470 * grad_quad_phi(N,node_index1)
00471 * linear_phi(vertex_index)
00472 * wJ;
00473 }
00474 }
00475 }
00476 }
00477
00478 for (unsigned vertex_index=0; vertex_index<NUM_VERTICES_PER_ELEMENT; vertex_index++)
00479 {
00480 unsigned index1 = NUM_NODES_PER_ELEMENT*DIM + vertex_index;
00481
00482 for (unsigned index2=0; index2<NUM_NODES_PER_ELEMENT*DIM; index2++)
00483 {
00484 unsigned spatial_dim2 = index2%DIM;
00485 unsigned node_index2 = (index2-spatial_dim2)/DIM;
00486
00487 for (unsigned M=0; M<DIM; M++)
00488 {
00489 rAElem(index1,index2) += linear_phi(vertex_index)
00490 * detF
00491 * inv_F(M,spatial_dim2)
00492 * grad_quad_phi(M,node_index2)
00493 * wJ;
00494 }
00495 }
00496
00498
00499
00500
00501
00503 for (unsigned vertex_index2=0; vertex_index2<NUM_VERTICES_PER_ELEMENT; vertex_index2++)
00504 {
00505 unsigned index2 = NUM_NODES_PER_ELEMENT*DIM + vertex_index2;
00506 rAElemPrecond(index1,index2) += linear_phi(vertex_index)
00507 * linear_phi(vertex_index2)
00508 * wJ;
00509 }
00510 }
00511 }
00512 }
00513
00514
00515 if (assembleJacobian)
00516 {
00517
00518
00519
00520 rAElemPrecond = rAElemPrecond + rAElem;
00521 }
00522 }
00523
00524 template<size_t DIM>
00525 void NonlinearElasticityAssembler<DIM>::AssembleOnBoundaryElement(
00526 BoundaryElement<DIM-1,DIM>& rBoundaryElement,
00527 c_matrix<double,BOUNDARY_STENCIL_SIZE,BOUNDARY_STENCIL_SIZE>& rAelem,
00528 c_vector<double,BOUNDARY_STENCIL_SIZE>& rBelem,
00529 c_vector<double,DIM>& rTraction,
00530 bool assembleResidual,
00531 bool assembleJacobian)
00532 {
00533 rAelem.clear();
00534 rBelem.clear();
00535
00536 if (assembleJacobian && !assembleResidual)
00537 {
00538
00539 return;
00540 }
00541
00542 c_vector<double, DIM> weighted_direction;
00543 double jacobian_determinant;
00544 mpQuadMesh->GetWeightedDirectionForBoundaryElement(rBoundaryElement.GetIndex(), weighted_direction, jacobian_determinant);
00545
00546 c_vector<double,NUM_NODES_PER_BOUNDARY_ELEMENT> phi;
00547
00548 for (unsigned quad_index=0; quad_index<this->mpBoundaryQuadratureRule->GetNumQuadPoints(); quad_index++)
00549 {
00550 double wJ = jacobian_determinant * this->mpBoundaryQuadratureRule->GetWeight(quad_index);
00551
00552 const ChastePoint<DIM-1>& quad_point = this->mpBoundaryQuadratureRule->rGetQuadPoint(quad_index);
00553
00554 QuadraticBasisFunction<DIM-1>::ComputeBasisFunctions(quad_point, phi);
00555
00556
00557
00558 c_vector<double,DIM> traction = zero_vector<double>(DIM);
00559 if (this->mUsingTractionBoundaryConditionFunction)
00560 {
00561 c_vector<double,DIM> X = zero_vector<double>(DIM);
00562 for (unsigned node_index=0; node_index<NUM_NODES_PER_BOUNDARY_ELEMENT; node_index++)
00563 {
00564 X += phi(node_index)*mpQuadMesh->GetNode( rBoundaryElement.GetNodeGlobalIndex(node_index) )->rGetLocation();
00565 }
00566 traction = (*(this->mpTractionBoundaryConditionFunction))(X);
00567 }
00568 else
00569 {
00570 traction = rTraction;
00571 }
00572
00573
00574 for (unsigned index=0; index<NUM_NODES_PER_BOUNDARY_ELEMENT*DIM; index++)
00575 {
00576 unsigned spatial_dim = index%DIM;
00577 unsigned node_index = (index-spatial_dim)/DIM;
00578
00579 assert(node_index < NUM_NODES_PER_BOUNDARY_ELEMENT);
00580
00581 rBelem(index) -= traction(spatial_dim)
00582 * phi(node_index)
00583 * wJ;
00584 }
00585 }
00586 }
00587
00588 template<size_t DIM>
00589 void NonlinearElasticityAssembler<DIM>::FormInitialGuess()
00590 {
00591 this->mCurrentSolution.resize(this->mNumDofs, 0.0);
00592
00593 for (unsigned i=0; i<mpQuadMesh->GetNumElements(); i++)
00594 {
00595 double zero_strain_pressure;
00596 if (this->mMaterialLaws.size()==1)
00597 {
00598
00599 zero_strain_pressure = this->mMaterialLaws[0]->GetZeroStrainPressure();
00600 }
00601 else
00602 {
00603
00604 zero_strain_pressure = this->mMaterialLaws[i]->GetZeroStrainPressure();
00605 }
00606
00607
00608 for (unsigned j=0; j<NUM_VERTICES_PER_ELEMENT; j++)
00609 {
00610 unsigned index = mpQuadMesh->GetElement(i)->GetNodeGlobalIndex(j);
00611 this->mCurrentSolution[ DIM*mpQuadMesh->GetNumNodes() + index ] = zero_strain_pressure;
00612 }
00613 }
00614 }
00615
00616 template<size_t DIM>
00617 void NonlinearElasticityAssembler<DIM>::Initialise(std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00618 {
00619 assert(mpQuadMesh);
00620
00621 for (unsigned i=0; i<this->mFixedNodes.size(); i++)
00622 {
00623 assert(this->mFixedNodes[i] < mpQuadMesh->GetNumNodes());
00624 }
00625
00626 this->mpQuadratureRule = new GaussianQuadratureRule<DIM>(3);
00627 this->mpBoundaryQuadratureRule = new GaussianQuadratureRule<DIM-1>(3);
00628
00629 FormInitialGuess();
00630
00631
00632
00633 if (pFixedNodeLocations == NULL)
00634 {
00635 this->mFixedNodeDisplacements.clear();
00636 for (unsigned i=0; i<this->mFixedNodes.size(); i++)
00637 {
00638 this->mFixedNodeDisplacements.push_back(zero_vector<double>(DIM));
00639 }
00640 }
00641 else
00642 {
00643 assert(pFixedNodeLocations->size()==this->mFixedNodes.size());
00644 for (unsigned i=0; i<this->mFixedNodes.size(); i++)
00645 {
00646 unsigned index = this->mFixedNodes[i];
00647 c_vector<double,DIM> displacement = (*pFixedNodeLocations)[i]-mpQuadMesh->GetNode(index)->rGetLocation();
00648 this->mFixedNodeDisplacements.push_back(displacement);
00649 }
00650 }
00651 assert(this->mFixedNodeDisplacements.size()==this->mFixedNodes.size());
00652 }
00653
00654 template<size_t DIM>
00655 NonlinearElasticityAssembler<DIM>::NonlinearElasticityAssembler(
00656 QuadraticMesh<DIM>* pQuadMesh,
00657 AbstractIncompressibleMaterialLaw<DIM>* pMaterialLaw,
00658 c_vector<double,DIM> bodyForce,
00659 double density,
00660 std::string outputDirectory,
00661 std::vector<unsigned>& fixedNodes,
00662 std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00663 : AbstractNonlinearElasticityAssembler<DIM>(DIM*pQuadMesh->GetNumNodes()+pQuadMesh->GetNumVertices(),
00664 pMaterialLaw, bodyForce, density,
00665 outputDirectory, fixedNodes),
00666 mpQuadMesh(pQuadMesh)
00667 {
00668 Initialise(pFixedNodeLocations);
00669 }
00670
00671
00672 template<size_t DIM>
00673 NonlinearElasticityAssembler<DIM>::NonlinearElasticityAssembler(
00674 QuadraticMesh<DIM>* pQuadMesh,
00675 std::vector<AbstractIncompressibleMaterialLaw<DIM>*>& rMaterialLaws,
00676 c_vector<double,DIM> bodyForce,
00677 double density,
00678 std::string outputDirectory,
00679 std::vector<unsigned>& fixedNodes,
00680 std::vector<c_vector<double,DIM> >* pFixedNodeLocations)
00681 : AbstractNonlinearElasticityAssembler<DIM>(DIM*pQuadMesh->GetNumNodes()+pQuadMesh->GetNumVertices(),
00682 rMaterialLaws, bodyForce, density,
00683 outputDirectory, fixedNodes),
00684 mpQuadMesh(pQuadMesh)
00685 {
00686 assert(rMaterialLaws.size()==pQuadMesh->GetNumElements());
00687 Initialise(pFixedNodeLocations);
00688 }
00689
00690
00691 template<size_t DIM>
00692 NonlinearElasticityAssembler<DIM>::~NonlinearElasticityAssembler()
00693 {
00694 delete mpQuadratureRule;
00695 delete mpBoundaryQuadratureRule;
00696 }
00697
00698 template<size_t DIM>
00699 void NonlinearElasticityAssembler<DIM>::SetSurfaceTractionBoundaryConditions(
00700 std::vector<BoundaryElement<DIM-1,DIM>*>& rBoundaryElements,
00701 std::vector<c_vector<double,DIM> >& rSurfaceTractions)
00702 {
00703 assert(rBoundaryElements.size()==rSurfaceTractions.size());
00704 mBoundaryElements = rBoundaryElements;
00705 this->mSurfaceTractions = rSurfaceTractions;
00706 }
00707
00708 template<size_t DIM>
00709 void NonlinearElasticityAssembler<DIM>::SetFunctionalTractionBoundaryCondition(
00710 std::vector<BoundaryElement<DIM-1,DIM>*> rBoundaryElements,
00711 c_vector<double,DIM> (*pFunction)(c_vector<double,DIM>&))
00712 {
00713 mBoundaryElements = rBoundaryElements;
00714 this->mUsingTractionBoundaryConditionFunction = true;
00715 this->mpTractionBoundaryConditionFunction = pFunction;
00716 }
00717
00718 template<size_t DIM>
00719 std::vector<double>& NonlinearElasticityAssembler<DIM>::rGetPressures()
00720 {
00721 this->mPressures.clear();
00722 this->mPressures.resize(mpQuadMesh->GetNumVertices());
00723
00724 for (unsigned i=0; i<mpQuadMesh->GetNumVertices(); i++)
00725 {
00726 this->mPressures[i] = this->mCurrentSolution[DIM*mpQuadMesh->GetNumNodes() + i];
00727 }
00728 return this->mPressures;
00729 }
00730
00731 template<size_t DIM>
00732 std::vector<c_vector<double,DIM> >& NonlinearElasticityAssembler<DIM>::rGetDeformedPosition()
00733 {
00734 this->mDeformedPosition.resize(mpQuadMesh->GetNumNodes(), zero_vector<double>(DIM));
00735 for (unsigned i=0; i<mpQuadMesh->GetNumNodes(); i++)
00736 {
00737 for (unsigned j=0; j<DIM; j++)
00738 {
00739 this->mDeformedPosition[i](j) = mpQuadMesh->GetNode(i)->rGetLocation()[j] + this->mCurrentSolution[DIM*i+j];
00740 }
00741 }
00742 return this->mDeformedPosition;
00743 }
00744
00745
00747
00749
00750
00751 template class NonlinearElasticityAssembler<2>;
00752 template class NonlinearElasticityAssembler<3>;