1 #include "LinearElasticityModel.hxx"
2 #include "SparseMatrixPetsc.hxx"
11 int LinearElasticityModel::fact(int n)
13 return (n == 1 || n == 0) ? 1 : fact(n - 1) * n;
15 int LinearElasticityModel::unknownNodeIndex(int globalIndex, std::vector< int > dirichletNodes)
16 {//assumes Dirichlet node numbering is strictly increasing
17 int j=0;//indice de parcours des noeuds frontière avec CL Dirichlet
18 int boundarySize=dirichletNodes.size();
19 while(j<boundarySize and dirichletNodes[j]<globalIndex)
22 return globalIndex-boundarySize;
23 else if (dirichletNodes[j]>globalIndex)
26 throw CdmathException("LinearElasticityModel::unknownNodeIndex : Error : node is a Dirichlet boundary node");
29 int LinearElasticityModel::globalNodeIndex(int unknownNodeIndex, std::vector< int > dirichletNodes)
30 {//assumes Dirichlet boundary node numbering is strictly increasing
31 int boundarySize=dirichletNodes.size();
32 /* trivial case where all boundary nodes are Neumann BC */
34 return unknownNodeIndex;
36 double unknownNodeMax=-1;//max unknown node number in the interval between jth and (j+1)th Dirichlet boundary nodes
37 int j=0;//indice de parcours des noeuds frontière
38 //On cherche l'intervale [j,j+1] qui contient le noeud de numéro interieur unknownNodeIndex
39 while(j+1<boundarySize and unknownNodeMax<unknownNodeIndex)
41 unknownNodeMax += dirichletNodes[j+1]-dirichletNodes[j]-1;
46 return unknownNodeIndex+boundarySize;
47 else //unknownNodeMax>=unknownNodeIndex) hence our node global number is between dirichletNodes[j-1] and dirichletNodes[j]
48 return unknownNodeIndex - unknownNodeMax + dirichletNodes[j]-1;
51 LinearElasticityModel::LinearElasticityModel(int dim, bool FECalculation, double rho, double lambda, double mu){
52 PetscBool petscInitialized;
53 PetscInitialized(&petscInitialized);
55 PetscInitialize(NULL,NULL,0,0);
59 std::cout<<"First Lamé coefficient="<<lambda<<endl;
60 throw CdmathException("Error : First Lamé coefficient lambda cannot be negative");
62 if(2*mu+dim*lambda < 0.)
64 std::cout<<"First Lamé coefficient="<<lambda<<", second Lamé coefficient="<<mu<<", 2*mu+dim*lambda= "<<2*mu+dim*lambda<<endl;
65 throw CdmathException("Error : 2*mu+dim*lambda cannot be negative");
69 std::cout<<"space dimension="<<dim<<endl;
70 throw CdmathException("Error : parameter dim cannot be negative");
73 _FECalculation=FECalculation;
78 _initializedMemory=false;
86 _boundaryNodeIds=std::vector< int >(0);
87 _dirichletNodeIds=std::vector< int >(0);
91 _dirichletValuesSet=false;
92 _neumannValuesSet=false;
96 _precision_Newton=_precision;
97 _MaxIterLinearSolver=0;//During several newton iterations, stores the max petssc interations
99 int _PetscIts=0;//the number of iterations of the linear solver
100 _ksptype = (char*)&KSPGMRES;
101 _pctype = (char*)&PCLU;
102 _conditionNumber=false;
105 //parameters for monitoring simulation
108 _runLogFile=new ofstream;
110 //result save parameters
111 _fileName = "LinearElasticityProblem";
112 char result[ PATH_MAX ];//extracting current directory
113 getcwd(result, PATH_MAX );
114 _path=string( result );
116 _computationCompletedSuccessfully=false;
118 //heat transfer parameters
122 _densityFieldSet=false;
125 void LinearElasticityModel::initialize()
127 _runLogFile->open((_fileName+".log").c_str(), ios::out | ios::trunc);;//for creation of a log file to save the history of the simulation
130 throw CdmathException("LinearElasticityModel::initialize() set mesh first");
133 cout<<"!!!! Initialisation of the computation of the elastic deformation of a solid using ";
134 *_runLogFile<<"!!!!! Initialisation of the computation of the elastic deformation of a solid using ";
137 cout<< "Finite volumes method"<<endl<<endl;
138 *_runLogFile<< "Finite volumes method"<<endl<<endl;
142 cout<< "Finite elements method"<<endl<<endl;
143 *_runLogFile<< "Finite elements method"<<endl<<endl;
146 /**************** Field creation *********************/
148 if(!_densityFieldSet){
150 _densityField=Field("Density",NODES,_mesh,1);
151 for(int i =0; i<_Nnodes; i++)
152 _densityField(i) = _rho;
155 _densityField=Field("Density",CELLS,_mesh,1);
156 for(int i =0; i<_Nmailles; i++)
157 _densityField(i) = _rho;
159 _densityFieldSet=true;
162 /* Détection des noeuds frontière avec une condition limite de Dirichlet */
165 if(_NboundaryNodes==_Nnodes)
166 cout<<"!!!!! Warning : all nodes are boundary nodes !!!!!"<<endl<<endl;
168 for(int i=0; i<_NboundaryNodes; i++)
170 std::map<int,double>::iterator it=_dirichletBoundaryValues.find(_boundaryNodeIds[i]);
171 if( it != _dirichletBoundaryValues.end() )
172 _dirichletNodeIds.push_back(_boundaryNodeIds[i]);
173 else if( _mesh.getNode(_boundaryNodeIds[i]).getGroupNames().size()==0 )
175 cout<<"!!! No boundary value set for boundary node" << _boundaryNodeIds[i]<< endl;
176 *_runLogFile<< "!!! No boundary value set for boundary node" << _boundaryNodeIds[i]<<endl;
177 _runLogFile->close();
178 throw CdmathException("Missing boundary value");
180 else if(_limitField[_mesh.getNode(_boundaryNodeIds[i]).getGroupName()].bcType==NoneBCLinearElasticity)
182 cout<<"!!! No boundary condition set for boundary node " << _boundaryNodeIds[i]<< endl;
183 *_runLogFile<< "!!!No boundary condition set for boundary node " << _boundaryNodeIds[i]<<endl;
184 _runLogFile->close();
185 throw CdmathException("Missing boundary condition");
187 else if(_limitField[_mesh.getNode(_boundaryNodeIds[i]).getGroupName()].bcType==DirichletLinearElasticity)
188 _dirichletNodeIds.push_back(_boundaryNodeIds[i]);
189 else if(_limitField[_mesh.getNode(_boundaryNodeIds[i]).getGroupName()].bcType!=NeumannLinearElasticity)
191 cout<<"!!! Wrong boundary condition "<< _limitField[_mesh.getNode(_boundaryNodeIds[i]).getGroupName()].bcType<< " set for boundary node " << _boundaryNodeIds[i]<< endl;
192 cout<<"!!! Accepted boundary conditions are Dirichlet "<< DirichletLinearElasticity <<" and Neumann "<< NeumannLinearElasticity << endl;
193 *_runLogFile<< "Wrong boundary condition "<< _limitField[_mesh.getNode(_boundaryNodeIds[i]).getGroupName()].bcType<< " set for boundary node " << _boundaryNodeIds[i]<<endl;
194 *_runLogFile<< "Accepted boundary conditions are Dirichlet "<< DirichletLinearElasticity <<" and Neumann "<< NeumannLinearElasticity <<endl;
195 _runLogFile->close();
196 throw CdmathException("Wrong boundary condition");
199 _NdirichletNodes=_dirichletNodeIds.size();
200 _NunknownNodes=_Nnodes - _NdirichletNodes;
201 cout<<"Number of unknown nodes " << _NunknownNodes <<", Number of boundary nodes " << _NboundaryNodes<< ", Number of Dirichlet boundary nodes " << _NdirichletNodes <<endl<<endl;
202 *_runLogFile<<"Number of unknown nodes " << _NunknownNodes <<", Number of boundary nodes " << _NboundaryNodes<< ", Number of Dirichlet boundary nodes " << _NdirichletNodes <<endl<<endl;
205 //creation de la matrice
207 MatCreateSeqAIJ(PETSC_COMM_SELF, _Nmailles*_nVar, _Nmailles*_nVar, (1+_neibMaxNbCells), PETSC_NULL, &_A);
209 MatCreateSeqAIJ(PETSC_COMM_SELF, _NunknownNodes*_nVar, _NunknownNodes*_nVar, (1+_neibMaxNbNodes), PETSC_NULL, &_A);
211 VecCreate(PETSC_COMM_SELF, &_displacements);
213 VecDuplicate(_displacements, &_b);//RHS of the linear system
216 KSPCreate(PETSC_COMM_SELF, &_ksp);
217 KSPSetType(_ksp, _ksptype);
218 // if(_ksptype == KSPGMRES) KSPGMRESSetRestart(_ksp,10000);
219 KSPSetTolerances(_ksp,_precision,_precision,PETSC_DEFAULT,_maxPetscIts);
220 KSPGetPC(_ksp, &_pc);
221 PCSetType(_pc, _pctype);
223 //Checking whether all boundary conditions are Neumann boundary condition
224 if(!_neumannValuesSet)//Boundary conditions set via LimitField structure
226 map<string, LimitFieldLinearElasticity>::iterator it = _limitField.begin();
227 while(it != _limitField.end() and (it->second).bcType == NeumannLinearElasticity)
229 _onlyNeumannBC = (it == _limitField.end() && _limitField.size()>0);//what if _limitField.size()==0 ???
233 _onlyNeumannBC = _neumannBoundaryValues.size()==_NboundaryNodes;
235 _onlyNeumannBC = _neumannBoundaryValues.size()==_mesh.getBoundaryFaceIds().size();
237 //If only Neumann BC, then matrix is singular and solution should be sought in space of mean zero vectors
240 std::cout<<"## Warning all boundary conditions are Neumann. System matrix is not invertible since constant vectors are in the kernel."<<std::endl;
241 std::cout<<"## As a consequence we seek a zero sum solution, and exact (LU and CHOLESKY) and incomplete factorisations (ILU and ICC) may fail."<<std::endl<<endl;
242 *_runLogFile<<"## Warning all boundary condition are Neumann. System matrix is not invertible since constant vectors are in the kernel."<<std::endl;
243 *_runLogFile<<"## As a consequence we seek a zero sum solution, and exact (LU and CHOLESKY) and incomplete factorisations (ILU and ICC) may fail."<<std::endl<<endl;
245 //Check that the matrix is symmetric
246 PetscBool isSymetric;
247 MatIsSymmetric(_A,_precision,&isSymetric);
250 cout<<"Singular matrix is not symmetric, tolerance= "<< _precision<<endl;
251 throw CdmathException("Singular matrix should be symmetric with kernel composed of constant vectors");
254 MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, PETSC_NULL, &nullsp);
255 MatSetNullSpace(_A, nullsp);
256 MatSetTransposeNullSpace(_A, nullsp);
257 MatNullSpaceDestroy(&nullsp);
258 //PCFactorSetShiftType(_pc,MAT_SHIFT_NONZERO);
259 //PCFactorSetShiftAmount(_pc,1e-10);
262 _initializedMemory=true;
266 Vector LinearElasticityModel::gradientNodal(Matrix M, vector< double > values){
267 vector< Matrix > matrices(_Ndim);
269 for (int idim=0; idim<_Ndim;idim++){
270 matrices[idim]=M.deepCopy();
271 for (int jdim=0; jdim<_Ndim+1;jdim++)
272 matrices[idim](jdim,idim) = values[jdim] ;
275 Vector result(_Ndim);
276 for (int idim=0; idim<_Ndim;idim++)
277 result[idim] = matrices[idim].determinant();
282 double LinearElasticityModel::computeStiffnessMatrix(bool & stop)
287 result=computeStiffnessMatrixFE(stop);
289 result=computeStiffnessMatrixFV(stop);
291 if(_verbose or _system)
292 MatView(_A,PETSC_VIEWER_STDOUT_SELF);
297 double LinearElasticityModel::computeStiffnessMatrixFE(bool & stop){
304 Matrix M(_Ndim+1,_Ndim+1);//cell geometry matrix
305 std::vector< Vector > GradShapeFuncs(_Ndim+1);//shape functions of cell nodes
306 std::vector< int > nodeIds(_Ndim+1);//cell node Ids
307 std::vector< Node > nodes(_Ndim+1);//cell nodes
308 int i_int, j_int; //index of nodes j and k considered as unknown nodes
309 bool dirichletCell_treated;
311 std::vector< vector< double > > values(_Ndim+1,vector< double >(_Ndim+1,0));//values of shape functions on cell node
312 for (int idim=0; idim<_Ndim+1;idim++)
313 values[idim][idim]=1;
315 /* parameters for boundary treatment */
316 vector< Vector > valuesBorder(_Ndim+1);
317 Vector GradShapeFuncBorder(_Ndim+1);
319 for (int j=0; j<_Nmailles;j++)
321 Cj = _mesh.getCell(j);
323 for (int idim=0; idim<_Ndim+1;idim++){
324 nodeIds[idim]=Cj.getNodeId(idim);
325 nodes[idim]=_mesh.getNode(nodeIds[idim]);
326 for (int jdim=0; jdim<_Ndim;jdim++)
327 M(idim,jdim)=nodes[idim].getPoint()[jdim];
330 for (int idim=0; idim<_Ndim+1;idim++)
331 GradShapeFuncs[idim]=gradientNodal(M,values[idim])/fact(_Ndim);
333 /* Loop on the edges of the cell */
334 for (int idim=0; idim<_Ndim+1;idim++)
336 if(find(_dirichletNodeIds.begin(),_dirichletNodeIds.end(),nodeIds[idim])==_dirichletNodeIds.end())//!_mesh.isBorderNode(nodeIds[idim])
337 {//First node of the edge is not Dirichlet node
338 i_int=unknownNodeIndex(nodeIds[idim], _dirichletNodeIds);//assumes Dirichlet boundary node numbering is strictly increasing
339 dirichletCell_treated=false;
340 for (int jdim=0; jdim<_Ndim+1;jdim++)
342 if(find(_dirichletNodeIds.begin(),_dirichletNodeIds.end(),nodeIds[jdim])==_dirichletNodeIds.end())//!_mesh.isBorderNode(nodeIds[jdim])
343 {//Second node of the edge is not Dirichlet node
344 j_int= unknownNodeIndex(nodeIds[jdim], _dirichletNodeIds);//assumes Dirichlet boundary node numbering is strictly increasing
345 MatSetValue(_A,i_int,j_int,_conductivity*(_DiffusionTensor*GradShapeFuncs[idim])*GradShapeFuncs[jdim]/Cj.getMeasure(), ADD_VALUES);
347 else if (!dirichletCell_treated)
348 {//Second node of the edge is a Dirichlet node
349 dirichletCell_treated=true;
350 for (int kdim=0; kdim<_Ndim+1;kdim++)
352 std::map<int,double>::iterator it=_dirichletBoundaryValues.find(nodeIds[kdim]);
353 if( it != _dirichletBoundaryValues.end() )
355 if( _dirichletValuesSet )//New way of storing BC
356 valuesBorder[kdim]=_dirichletBoundaryValues[it->second];
357 else //old way of storing BC
358 valuesBorder[kdim]=_limitField[_mesh.getNode(nodeIds[kdim]).getGroupName()].displacement;
361 valuesBorder[kdim]=Vector(_Ndim);
363 GradShapeFuncBorder=gradientNodal(M,valuesBorder)/fact(_Ndim);
364 double coeff =-_conductivity*(_DiffusionTensor*GradShapeFuncBorder)*GradShapeFuncs[idim]/Cj.getMeasure();
365 VecSetValue(_b,i_int,coeff, ADD_VALUES);
372 //Calcul de la contribution de la condition limite de Neumann au second membre
373 if( _NdirichletNodes !=_NboundaryNodes)
375 vector< int > boundaryFaces = _mesh.getBoundaryFaceIds();
376 int NboundaryFaces=boundaryFaces.size();
377 for(int i = 0; i< NboundaryFaces ; i++)//On parcourt les faces du bord
379 Face Fi = _mesh.getFace(i);
380 for(int j = 0 ; j<_Ndim ; j++)//On parcourt les noeuds de la face
382 if(find(_dirichletNodeIds.begin(),_dirichletNodeIds.end(),Fi.getNodeId(j))==_dirichletNodeIds.end())//node j is an Neumann BC node (not a Dirichlet BC node)
384 j_int=unknownNodeIndex(Fi.getNodeId(j), _dirichletNodeIds);//indice du noeud j en tant que noeud inconnu
385 if( _neumannValuesSet )
386 coeff =Fi.getMeasure()/(_Ndim+1)*_neumannBoundaryValues[Fi.getNodeId(j)];
388 coeff =Fi.getMeasure()/(_Ndim+1)*_limitField[_mesh.getNode(Fi.getNodeId(j)).getGroupName()].normalForce;
389 VecSetValue(_b, j_int, coeff, ADD_VALUES);
395 MatAssemblyBegin(_A, MAT_FINAL_ASSEMBLY);
396 MatAssemblyEnd(_A, MAT_FINAL_ASSEMBLY);
397 VecAssemblyBegin(_b);
405 double LinearElasticityModel::computeStiffnessMatrixFV(bool & stop){
406 long nbFaces = _mesh.getNumberOfFaces();
410 double inv_dxi, inv_dxj;
411 double barycenterDistance,coeff;
412 Vector normale(_Ndim);
415 std::vector< int > idCells;
419 for (int j=0; j<nbFaces;j++){
420 Fj = _mesh.getFace(j);
422 // compute the normal vector corresponding to face j : from idCells[0] to idCells[1]
423 idCells = Fj.getCellsId();
424 Cell1 = _mesh.getCell(idCells[0]);
426 for(int l=0; l<Cell1.getNumberOfFaces(); l++){
427 if (j == Cell1.getFacesId()[l]){
428 for (int idim = 0; idim < _Ndim; ++idim)
429 normale[idim] = Cell1.getNormalVector(l,idim);
434 //Compute velocity at the face Fj
435 dn=_conductivity*(_DiffusionTensor*normale)*normale;
437 // compute 1/dxi = volume of Ci/area of Fj
438 inv_dxi = Fj.getMeasure()/Cell1.getMeasure();
440 // If Fj is on the boundary
441 if (Fj.getNumberOfCells()==1) {
444 cout << "face numero " << j << " cellule frontiere " << idCells[0] << " ; vecteur normal=(";
445 for(int p=0; p<_Ndim; p++)
446 cout << normale[p] << ",";
450 std::map<int,double>::iterator it=_dirichletBoundaryValues.find(j);
451 if( it != _dirichletBoundaryValues.end() )
453 barycenterDistance=Cell1.getBarryCenter().distance(Fj.getBarryCenter());
454 MatSetValue(_A,idm,idm,dn*inv_dxi/barycenterDistance , ADD_VALUES);
455 VecSetValue(_b,idm, dn*inv_dxi/barycenterDistance*it->second, ADD_VALUES);
459 nameOfGroup = Fj.getGroupName();
461 if (_limitField[nameOfGroup].bcType==NeumannLinearElasticity){
462 if( _neumannValuesSet )
463 coeff =Fi.getMeasure()*_neumannBoundaryValues[j];
465 coeff =Fi.getMeasure()*_limitField[nameOfGroup].normalForce;
466 VecSetValue(_b,idm, coeff, ADD_VALUES);
468 else if(_limitField[nameOfGroup].bcType==DirichletLinearElasticity){
469 barycenterDistance=Cell1.getBarryCenter().distance(Fj.getBarryCenter());
470 MatSetValue(_A,idm,idm,dn*inv_dxi/barycenterDistance , ADD_VALUES);
471 VecSetValue(_b,idm, dn*inv_dxi/barycenterDistance*_limitField[nameOfGroup].displacement, ADD_VALUES);
475 cout<<"!!!!!!!!!!!!!!! Error LinearElasticityModel::computeStiffnessMatrixFV !!!!!!!!!!"<<endl;
476 cout<<"!!!!!! Boundary condition not accepted for boundary named !!!!!!!!!!"<<nameOfGroup<< ", _limitField[nameOfGroup].bcType= "<<_limitField[nameOfGroup].bcType<<endl;
477 cout<<"Accepted boundary conditions are Neumann "<<NeumannLinearElasticity<< " and Dirichlet "<<DirichletLinearElasticity<<endl;
478 *_runLogFile<<"!!!!!! Boundary condition not accepted for boundary named !!!!!!!!!!"<<nameOfGroup<< ", _limitField[nameOfGroup].bcType= "<<_limitField[nameOfGroup].bcType<<endl;
479 _runLogFile->close();
480 throw CdmathException("Boundary condition not accepted");
483 // if Fj is inside the domain
484 } else if (Fj.getNumberOfCells()==2 ){
487 cout << "face numero " << j << " cellule gauche " << idCells[0] << " cellule droite " << idCells[1];
488 cout << " ; vecteur normal=(";
489 for(int p=0; p<_Ndim; p++)
490 cout << normale[p] << ",";
493 Cell2 = _mesh.getCell(idCells[1]);
496 inv_dxj = Fj.getMeasure()/Cell2.getMeasure();
498 inv_dxj = 1/Cell2.getMeasure();
500 barycenterDistance=Cell1.getBarryCenter().distance(Cell2.getBarryCenter());
502 MatSetValue(_A,idm,idm, dn*inv_dxi/barycenterDistance, ADD_VALUES);
503 MatSetValue(_A,idm,idn,-dn*inv_dxi/barycenterDistance, ADD_VALUES);
504 MatSetValue(_A,idn,idn, dn*inv_dxj/barycenterDistance, ADD_VALUES);
505 MatSetValue(_A,idn,idm,-dn*inv_dxj/barycenterDistance, ADD_VALUES);
509 *_runLogFile<<"LinearElasticityModel::computeStiffnessMatrixFV(): incompatible number of cells around a face"<<endl;
510 throw CdmathException("LinearElasticityModel::computeStiffnessMatrixFV(): incompatible number of cells around a face");
514 MatAssemblyBegin(_A, MAT_FINAL_ASSEMBLY);
515 MatAssemblyEnd(_A, MAT_FINAL_ASSEMBLY);
516 VecAssemblyBegin(_b);
524 double LinearElasticityModel::computeRHS(bool & stop)//Contribution of the PDE RHS to the linear systemm RHS (boundary conditions do contribute to the system RHS via the function computeStiffnessMatrix)
526 VecAssemblyBegin(_b);
529 for (int i=0; i<_Nmailles;i++)
530 for (int j=0; j<_Ndim;j++)
531 VecSetValue(_b,i*nVar+j,_gravity(j)*_densityField(i),ADD_VALUES);
535 std::vector< int > nodesId;
536 for (int i=0; i<_Nmailles;i++)
539 nodesId=Ci.getNodesId();
540 for (int j=0; j<nodesId.size();j++)
541 if(!_mesh.isBorderNode(nodesId[j]))
542 for (int k=0; k<_Ndim; k++)
543 VecSetValue(_b,unknownNodeIndex(nodesId[j], _dirichletNodeIds)*nVar+k, _gravity(k)*_densityField(j)*Ci.getMeasure()/(_Ndim+1),ADD_VALUES);
549 if(_verbose or _system)
550 VecView(_b,PETSC_VIEWER_STDOUT_SELF);
556 bool LinearElasticityModel::solveLinearSystem()
560 //Only implicit scheme considered
561 MatAssemblyBegin(_A, MAT_FINAL_ASSEMBLY);
562 MatAssemblyEnd( _A, MAT_FINAL_ASSEMBLY);
564 #if PETSC_VERSION_GREATER_3_5
565 KSPSetOperators(_ksp, _A, _A);
567 KSPSetOperators(_ksp, _A, _A,SAME_NONZERO_PATTERN);
571 KSPSetComputeEigenvalues(_ksp,PETSC_TRUE);
572 KSPSolve(_ksp, _b, _displacements);
574 KSPConvergedReason reason;
575 KSPGetConvergedReason(_ksp,&reason);
576 KSPGetIterationNumber(_ksp, &_PetscIts);
578 KSPGetResidualNorm(_ksp,&residu);
579 if (reason!=2 and reason!=3)
581 cout<<"!!!!!!!!!!!!! Erreur système linéaire : pas de convergence de Petsc."<<endl;
582 cout<<"!!!!!!!!!!!!! Itérations maximales "<<_maxPetscIts<<" atteintes, résidu="<<residu<<", précision demandée= "<<_precision<<endl;
583 cout<<"Solver used "<< _ksptype<<", preconditioner "<<_pctype<<", Final number of iteration= "<<_PetscIts<<endl;
584 *_runLogFile<<"!!!!!!!!!!!!! Erreur système linéaire : pas de convergence de Petsc."<<endl;
585 *_runLogFile<<"!!!!!!!!!!!!! Itérations maximales "<<_maxPetscIts<<" atteintes, résidu="<<residu<<", précision demandée= "<<_precision<<endl;
586 *_runLogFile<<"Solver used "<< _ksptype<<", preconditioner "<<_pctype<<", Final number of iteration= "<<_PetscIts<<endl;
587 _runLogFile->close();
589 resolutionOK = false;
592 if( _MaxIterLinearSolver < _PetscIts)
593 _MaxIterLinearSolver = _PetscIts;
594 cout<<"## Système linéaire résolu en "<<_PetscIts<<" itérations par le solveur "<< _ksptype<<" et le preconditioneur "<<_pctype<<", précision demandée= "<<_precision<<endl<<endl;
595 *_runLogFile<<"## Système linéaire résolu en "<<_PetscIts<<" itérations par le solveur "<< _ksptype<<" et le preconditioneur "<<_pctype<<", précision demandée= "<<_precision<<endl<<endl;
603 void LinearElasticityModel::setMesh(const Mesh &M)
605 if(_Ndim != M.getSpaceDimension() or _Ndim!=M.getMeshDimension())//for the moment we must have space dim=mesh dim
607 cout<< "Problem : dimension defined is "<<_Ndim<< " but mesh dimension= "<<M.getMeshDimension()<<", and space dimension is "<<M.getSpaceDimension()<<endl;
608 *_runLogFile<< "Problem : dim = "<<_Ndim<< " but mesh dim= "<<M.getMeshDimension()<<", mesh space dim= "<<M.getSpaceDimension()<<endl;
609 *_runLogFile<<"LinearElasticityModel::setMesh: mesh has incorrect dimension"<<endl;
610 _runLogFile->close();
611 throw CdmathException("LinearElasticityModel::setMesh: mesh has incorrect dimension");
615 _Nmailles = _mesh.getNumberOfCells();
616 _Nnodes = _mesh.getNumberOfNodes();
618 cout<<"Mesh has "<< _Nmailles << " cells and " << _Nnodes << " nodes"<<endl<<endl;;
619 *_runLogFile<<"Mesh has "<< _Nmailles << " cells and " << _Nnodes << " nodes"<<endl<<endl;
621 // find maximum nb of neibourghs
624 _VV=Field ("Displacements", CELLS, _mesh, _Ndim);
625 _neibMaxNbCells=_mesh.getMaxNbNeighbours(CELLS);
629 if(_Ndim==1 )//The 1D cdmath mesh is necessarily made of segments
630 cout<<"1D Finite element method on segments"<<endl;
633 if( _mesh.isTriangular() )//Mesh dim=2
634 cout<<"2D Finite element method on triangles"<<endl;
635 else if (_mesh.getMeshDimension()==1)//Mesh dim=1
636 cout<<"1D Finite element method on a 2D network : space dimension is "<<_Ndim<< ", mesh dimension is "<<_mesh.getMeshDimension()<<endl;
639 cout<<"Error Finite element with Space dimension "<<_Ndim<< ", and mesh dimension "<<_mesh.getMeshDimension()<< ", mesh should be either triangular either 1D network"<<endl;
640 *_runLogFile<<"LinearElasticityModel::setMesh: mesh has incorrect dimension"<<endl;
641 _runLogFile->close();
642 throw CdmathException("LinearElasticityModel::setMesh: mesh has incorrect cell types");
647 if( _mesh.isTetrahedral() )//Mesh dim=3
648 cout<<"3D Finite element method on tetrahedra"<<endl;
649 else if (_mesh.getMeshDimension()==2 and _mesh.isTriangular())//Mesh dim=2
650 cout<<"2D Finite element method on a 3D surface : space dimension is "<<_Ndim<< ", mesh dimension is "<<_mesh.getMeshDimension()<<endl;
651 else if (_mesh.getMeshDimension()==1)//Mesh dim=1
652 cout<<"1D Finite element method on a 3D network : space dimension is "<<_Ndim<< ", mesh dimension is "<<_mesh.getMeshDimension()<<endl;
655 cout<<"Error Finite element with Space dimension "<<_Ndim<< ", and mesh dimension "<<_mesh.getMeshDimension()<< ", mesh should be either tetrahedral, either a triangularised surface or 1D network"<<endl;
656 *_runLogFile<<"LinearElasticityModel::setMesh: mesh has incorrect dimension"<<endl;
657 _runLogFile->close();
658 throw CdmathException("LinearElasticityModel::setMesh: mesh has incorrect cell types");
662 _VV=Field ("Temperature", NODES, _mesh, _Ndim);
664 _neibMaxNbNodes=_mesh.getMaxNbNeighbours(NODES);
665 _boundaryNodeIds = _mesh.getBoundaryNodeIds();
666 _NboundaryNodes=_boundaryNodeIds.size();
672 void LinearElasticityModel::setLinearSolver(linearSolver kspType, preconditioner pcType)
674 //_maxPetscIts=maxIterationsPetsc;
675 // set linear solver algorithm
677 _ksptype = (char*)&KSPGMRES;
678 else if (kspType==CG)
679 _ksptype = (char*)&KSPCG;
680 else if (kspType==BCGS)
681 _ksptype = (char*)&KSPBCGS;
683 cout << "!!! Error : only 'GMRES', 'CG' or 'BCGS' is acceptable as a linear solver !!!" << endl;
684 *_runLogFile << "!!! Error : only 'GMRES', 'CG' or 'BCGS' is acceptable as a linear solver !!!" << endl;
685 _runLogFile->close();
686 throw CdmathException("!!! Error : only 'GMRES', 'CG' or 'BCGS' algorithm is acceptable !!!");
688 // set preconditioner
690 _pctype = (char*)&PCNONE;
691 else if (pcType ==LU)
692 _pctype = (char*)&PCLU;
693 else if (pcType == ILU)
694 _pctype = (char*)&PCILU;
695 else if (pcType ==CHOLESKY)
696 _pctype = (char*)&PCCHOLESKY;
697 else if (pcType == ICC)
698 _pctype = (char*)&PCICC;
700 cout << "!!! Error : only 'NONE', 'LU', 'ILU', 'CHOLESKY' or 'ICC' preconditioners are acceptable !!!" << endl;
701 *_runLogFile << "!!! Error : only 'NONE' or 'LU' or 'ILU' preconditioners are acceptable !!!" << endl;
702 _runLogFile->close();
703 throw CdmathException("!!! Error : only 'NONE' or 'LU' or 'ILU' preconditioners are acceptable !!!" );
707 bool LinearElasticityModel::solveStationaryProblem()
709 if(!_initializedMemory)
711 *_runLogFile<< "ProblemCoreFlows::run() call initialize() first"<< _fileName<<endl;
712 _runLogFile->close();
713 throw CdmathException("ProblemCoreFlows::run() call initialize() first");
715 bool stop=false; // Does the Problem want to stop (error) ?
717 cout<< "!!! Running test case "<< _fileName << " using ";
718 *_runLogFile<< "!!! Running test case "<< _fileName<< " using ";
722 cout<< "Finite volumes method"<<endl<<endl;
723 *_runLogFile<< "Finite volumes method"<<endl<<endl;
727 cout<< "Finite elements method"<<endl<<endl;
728 *_runLogFile<< "Finite elements method"<< endl<<endl;
731 computeStiffnessMatrix( stop);
733 cout << "Error : failed computing diffusion matrix, stopping calculation"<< endl;
734 *_runLogFile << "Error : failed computing diffusion matrix, stopping calculation"<< endl;
735 _runLogFile->close();
736 throw CdmathException("Failed computing diffusion matrix");
740 cout << "Error : failed computing right hand side, stopping calculation"<< endl;
741 *_runLogFile << "Error : failed computing right hand side, stopping calculation"<< endl;
742 throw CdmathException("Failed computing right hand side");
744 stop = !solveLinearSystem();
746 cout << "Error : failed solving linear system, stopping calculation"<< endl;
747 *_runLogFile << "Error : failed linear system, stopping calculation"<< endl;
748 _runLogFile->close();
749 throw CdmathException("Failed solving linear system");
752 _computationCompletedSuccessfully=true;
755 *_runLogFile<< "!!!!!! Computation successful"<< endl;
756 _runLogFile->close();
761 void LinearElasticityModel::save(){
762 cout<< "Saving numerical results"<<endl<<endl;
763 *_runLogFile<< "Saving numerical results"<< endl<<endl;
765 string resultFile(_path+"/LinearElasticityModel");//Results
768 resultFile+=_fileName;
770 // create mesh and component info
771 string suppress ="rm -rf "+resultFile+"_*";
772 system(suppress.c_str());//Nettoyage des précédents calculs identiques
774 if(_verbose or _system)
775 VecView(_displacements,PETSC_VIEWER_STDOUT_SELF);
779 for(int i=0; i<_Nmailles; i++)
781 for(int j=0; j<_nVar; j++)
784 VecGetValues(_displacements, 1, &k, &uk);
791 for(int i=0; i<_NunknownNodes; i++)
793 globalIndex = globalNodeIndex(i, _dirichletNodeIds);
794 for(int j=0; j<_nVar; j++)
797 VecGetValues(_displacements, 1, &k, &uk);
798 _VV(globalIndex,j)=uk;
804 for(int i=0; i<_NdirichletNodes; i++)
806 Ni=_mesh.getNode(_dirichletNodeIds[i]);
807 nameOfGroup = Ni.getGroupName();
808 for(int j=0; j<_nVar; j++)
809 _VV(_dirichletNodeIds[i])=_limitField[nameOfGroup].displacement[i];
816 _VV.writeVTK(resultFile);
819 _VV.writeMED(resultFile);
822 _VV.writeCSV(resultFile);
827 LinearElasticityModel::getOutputDisplacementField()
829 if(!_computationCompletedSuccessfully)
830 throw("Computation not performed yet or failed. No displacement field available");
835 void LinearElasticityModel::terminate()
837 VecDestroy(&_displacements);
842 LinearElasticityModel::setDirichletValues(map< int, double> dirichletBoundaryValues)
844 _dirichletValuesSet=true;
845 _dirichletBoundaryValues=dirichletBoundaryValues;
848 LinearElasticityModel::setNeumannValues(map< int, double> neumannBoundaryValues)
850 _neumannValuesSet=true;
851 _neumannBoundaryValues=neumannBoundaryValues;
855 LinearElasticityModel::getConditionNumber(bool isSingular, double tol) const
857 SparseMatrixPetsc A = SparseMatrixPetsc(_A);
858 return A.getConditionNumber( isSingular, tol);