1 #include "StationaryDiffusionEquation.hxx"
2 #include "DiffusionEquation.hxx"
4 #include "SparseMatrixPetsc.hxx"
12 StationaryDiffusionEquation::StationaryDiffusionEquation(int dim, bool FECalculation, double lambda,MPI_Comm comm){
13 /* Initialisation of PETSC */
14 //check if PETSC is already initialised
15 PetscBool petscInitialized;
16 PetscInitialized(&petscInitialized);
18 {//check if MPI is already initialised
20 MPI_Initialized(&mpiInitialized);
22 PETSC_COMM_WORLD = comm;
23 PetscInitialize(NULL,NULL,0,0);//Note this is ok if MPI has been been initialised independently from PETSC
25 MPI_Comm_rank(PETSC_COMM_WORLD,&_mpi_rank);
26 MPI_Comm_size(PETSC_COMM_WORLD,&_mpi_size);
27 PetscPrintf(PETSC_COMM_WORLD,"Simulation on %d processors\n",_mpi_size);//Prints to standard out, only from the first processor in the communicator. Calls from other processes are ignored.
28 PetscSynchronizedPrintf(PETSC_COMM_WORLD,"Processor [%d] ready for action\n",_mpi_rank);//Prints synchronized output from several processors. Output of the first processor is followed by that of the second, etc.
29 PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT);
33 std::cout<<"Conductivity="<<lambda<<endl;
34 throw CdmathException("Error : conductivity parameter lambda cannot be negative");
38 std::cout<<"Space dimension="<<dim<<endl;
39 throw CdmathException("Error : parameter dim cannot be negative");
42 _FECalculation=FECalculation;
46 _nVar=1;//scalar prolem
48 _diffusionMatrixSet=false;
49 _initializedMemory=false;
57 _boundaryNodeIds=std::vector< int >(0);
58 _dirichletNodeIds=std::vector< int >(0);
62 _dirichletValuesSet=false;
63 _neumannValuesSet=false;
67 _precision_Newton=_precision;
68 _MaxIterLinearSolver=0;//During several newton iterations, stores the max petssc interations
72 int _PetscIts=0;//the number of iterations of the linear solver
73 _ksptype = (char*)&KSPGMRES;
75 _pctype = (char*)&PCNONE;
77 _pctype = (char*)&PCLU;
79 _conditionNumber=false;
82 //parameters for monitoring simulation
85 _runLogFile=new ofstream;
87 //result save parameters
88 _fileName = "StationaryDiffusionProblem";
89 char result[ PATH_MAX ];//extracting current directory
90 getcwd(result, PATH_MAX );
91 _path=string( result );
93 _computationCompletedSuccessfully=false;
95 //heat transfer parameters
97 _fluidTemperatureFieldSet=false;
99 _heatPowerFieldSet=false;
100 _heatTransfertCoeff=0;
103 /* Default diffusion tensor is diagonal */
104 _DiffusionTensor=Matrix(_Ndim);
105 for(int idim=0;idim<_Ndim;idim++)
106 _DiffusionTensor(idim,idim)=_conductivity;
109 void StationaryDiffusionEquation::initialize()
113 _runLogFile->open((_fileName+".log").c_str(), ios::out | ios::trunc);;//for creation of a log file to save the history of the simulation
116 throw CdmathException("StationaryDiffusionEquation::initialize() set mesh first");
119 cout<<"!!!! Initialisation of the computation of the temperature diffusion in a solid using ";
120 *_runLogFile<<"!!!!! Initialisation of the computation of the temperature diffusion in a solid using ";
123 cout<< "Finite volumes method"<<endl<<endl;
124 *_runLogFile<< "Finite volumes method"<<endl<<endl;
128 cout<< "Finite elements method"<<endl<<endl;
129 *_runLogFile<< "Finite elements method"<<endl<<endl;
133 /**************** Field creation *********************/
135 if(!_heatPowerFieldSet){
136 _heatPowerField=Field("Heat power",_VV.getTypeOfField(),_mesh,1);
137 for(int i =0; i<_VV.getNumberOfElements(); i++)
138 _heatPowerField(i) = _heatSource;
139 _heatPowerFieldSet=true;
141 if(!_fluidTemperatureFieldSet){
142 _fluidTemperatureField=Field("Fluid temperature",_VV.getTypeOfField(),_mesh,1);
143 for(int i =0; i<_VV.getNumberOfElements(); i++)
144 _fluidTemperatureField(i) = _fluidTemperature;
145 _fluidTemperatureFieldSet=true;
148 /* Détection des noeuds frontière avec une condition limite de Dirichlet */
151 if(_NboundaryNodes==_Nnodes)
152 cout<<"!!!!! Warning : all nodes are boundary nodes !!!!!"<<endl<<endl;
154 for(int i=0; i<_NboundaryNodes; i++)
156 std::map<int,double>::iterator it=_dirichletBoundaryValues.find(_boundaryNodeIds[i]);
157 if( it != _dirichletBoundaryValues.end() )
158 _dirichletNodeIds.push_back(_boundaryNodeIds[i]);
159 else if( _mesh.getNode(_boundaryNodeIds[i]).getGroupNames().size()==0 )
161 cout<<"!!! No boundary group set for boundary node" << _boundaryNodeIds[i]<< endl;
162 *_runLogFile<< "!!! No boundary group set for boundary node" << _boundaryNodeIds[i]<<endl;
163 _runLogFile->close();
164 throw CdmathException("Missing boundary group");
166 else if(_limitField[_mesh.getNode(_boundaryNodeIds[i]).getGroupName()].bcType==NoneBCStationaryDiffusion)
168 cout<<"!!! No boundary condition set for boundary node " << _boundaryNodeIds[i]<< endl;
169 cout<<"!!! Accepted boundary conditions are DirichletStationaryDiffusion "<< DirichletStationaryDiffusion <<" and NeumannStationaryDiffusion "<< NeumannStationaryDiffusion << endl;
170 *_runLogFile<< "!!! No boundary condition set for boundary node " << _boundaryNodeIds[i]<<endl;
171 *_runLogFile<< "!!! Accepted boundary conditions are DirichletStationaryDiffusion "<< DirichletStationaryDiffusion <<" and NeumannStationaryDiffusion "<< NeumannStationaryDiffusion <<endl;
172 _runLogFile->close();
173 throw CdmathException("Missing boundary condition");
175 else if(_limitField[_mesh.getNode(_boundaryNodeIds[i]).getGroupName()].bcType==DirichletStationaryDiffusion)
176 _dirichletNodeIds.push_back(_boundaryNodeIds[i]);
177 else if(_limitField[_mesh.getNode(_boundaryNodeIds[i]).getGroupName()].bcType!=NeumannStationaryDiffusion)
179 cout<<"!!! Wrong boundary condition "<< _limitField[_mesh.getNode(_boundaryNodeIds[i]).getGroupName()].bcType<< " set for boundary node " << _boundaryNodeIds[i]<< endl;
180 cout<<"!!! Accepted boundary conditions are DirichletStationaryDiffusion "<< DirichletStationaryDiffusion <<" and NeumannStationaryDiffusion "<< NeumannStationaryDiffusion << endl;
181 *_runLogFile<< "!!! Wrong boundary condition "<< _limitField[_mesh.getNode(_boundaryNodeIds[i]).getGroupName()].bcType<< " set for boundary node " << _boundaryNodeIds[i]<<endl;
182 *_runLogFile<< "!!! Accepted boundary conditions are DirichletStationaryDiffusion "<< DirichletStationaryDiffusion <<" and NeumannStationaryDiffusion "<< NeumannStationaryDiffusion <<endl;
183 _runLogFile->close();
184 throw CdmathException("Wrong boundary condition");
187 _NdirichletNodes=_dirichletNodeIds.size();
188 _NunknownNodes=_Nnodes - _NdirichletNodes;
189 cout<<"Number of unknown nodes " << _NunknownNodes <<", Number of boundary nodes " << _NboundaryNodes<< ", Number of Dirichlet boundary nodes " << _NdirichletNodes <<endl<<endl;
190 *_runLogFile<<"Number of unknown nodes " << _NunknownNodes <<", Number of boundary nodes " << _NboundaryNodes<< ", Number of Dirichlet boundary nodes " << _NdirichletNodes <<endl<<endl;
195 _globalNbUnknowns = _Nmailles*_nVar;
197 MPI_Bcast(&_NunknownNodes, 1, MPI_INT, 0, PETSC_COMM_WORLD);
198 _globalNbUnknowns = _NunknownNodes*_nVar;
201 /* Vectors creations */
202 VecCreate(PETSC_COMM_WORLD, &_Tk);//main unknown
203 VecSetSizes(_Tk,PETSC_DECIDE,_globalNbUnknowns);
204 VecSetFromOptions(_Tk);
205 VecGetLocalSize(_Tk, &_localNbUnknowns);
207 VecDuplicate(_Tk, &_Tkm1);
208 VecDuplicate(_Tk, &_deltaT);
209 VecDuplicate(_Tk, &_b);//RHS of the linear system: _b=Tn/dt + _b0 + puisance volumique + couplage thermique avec le fluide
211 /* Matrix creation */
212 MatCreateAIJ(PETSC_COMM_WORLD, _localNbUnknowns, _localNbUnknowns, _globalNbUnknowns, _globalNbUnknowns, _d_nnz, PETSC_NULL, _o_nnz, PETSC_NULL, &_A);
214 /* Local sequential vector creation */
215 if(_mpi_size>1 && _mpi_rank == 0)
216 VecCreateSeq(PETSC_COMM_SELF,_globalNbUnknowns,&_Tk_seq);//For saving results on proc 0
219 VecScatterCreateToZero(_Tk,&_scat,&_Tk_seq);
222 KSPCreate(PETSC_COMM_WORLD, &_ksp);
223 KSPSetType(_ksp, _ksptype);
224 // if(_ksptype == KSPGMRES) KSPGMRESSetRestart(_ksp,10000);
225 KSPSetTolerances(_ksp,_precision,_precision,PETSC_DEFAULT,_maxPetscIts);
226 KSPGetPC(_ksp, &_pc);
227 PCSetType(_pc, _pctype);
229 //Checking whether all boundary conditions are Neumann boundary condition
230 //if(_FECalculation) _onlyNeumannBC = _NdirichletNodes==0;
231 if(!_neumannValuesSet)//Boundary conditions set via LimitField structure
233 map<string, LimitFieldStationaryDiffusion>::iterator it = _limitField.begin();
234 while(it != _limitField.end() and (it->second).bcType == NeumannStationaryDiffusion)
236 _onlyNeumannBC = (it == _limitField.end() && _limitField.size()>0);//what if _limitField.size()==0 ???
240 _onlyNeumannBC = _neumannBoundaryValues.size()==_NboundaryNodes;
242 _onlyNeumannBC = _neumannBoundaryValues.size()==_mesh.getBoundaryFaceIds().size();
244 //If only Neumann BC, then matrix is singular and solution should be sought in space of mean zero vectors
247 std::cout<<"### Warning : all boundary conditions are Neumann. System matrix is not invertible since constant vectors are in the kernel."<<std::endl;
248 std::cout<<"### Check the compatibility condition between the right hand side and the boundary data. For homogeneous Neumann BCs, the right hand side must have integral equal to zero."<<std::endl;
249 std::cout<<"### The system matrix being singular, we seek a zero sum solution, and exact (LU and CHOLESKY) and incomplete factorisations (ILU and ICC) may fail."<<std::endl<<endl;
250 *_runLogFile<<"### Warning : all boundary condition are Neumann. System matrix is not invertible since constant vectors are in the kernel."<<std::endl;
251 *_runLogFile<<"### The system matrix being singular, we seek a zero sum solution, and exact (LU and CHOLESKY) and incomplete factorisations (ILU and ICC) may fail."<<std::endl<<endl;
252 *_runLogFile<<"### Check the compatibility condition between the right hand side and the boundary data. For homogeneous Neumann BCs, the right hand side must have integral equal to zero."<<std::endl;
255 MatNullSpaceCreate(PETSC_COMM_WORLD, PETSC_TRUE, 0, PETSC_NULL, &nullsp);
256 MatSetNullSpace(_A, nullsp);
257 MatSetTransposeNullSpace(_A, nullsp);
258 MatNullSpaceDestroy(&nullsp);
259 //PCFactorSetShiftType(_pc,MAT_SHIFT_NONZERO);
260 //PCFactorSetShiftAmount(_pc,1e-10);
262 _initializedMemory=true;
265 double StationaryDiffusionEquation::computeTimeStep(bool & stop){
266 if(!_diffusionMatrixSet)//The diffusion matrix is computed once and for all time steps
267 computeDiffusionMatrix(stop);
269 _dt_src=computeRHS(stop);
274 double StationaryDiffusionEquation::computeDiffusionMatrix(bool & stop)
282 result=computeDiffusionMatrixFE(stop);
284 result=computeDiffusionMatrixFV(stop);
286 MatAssemblyBegin(_A, MAT_FINAL_ASSEMBLY);
287 MatAssemblyEnd( _A, MAT_FINAL_ASSEMBLY);
289 //Contribution from the solid/fluid heat exchange with assumption of constant heat transfer coefficient
290 //update value here if variable heat transfer coefficient
291 if(_heatTransfertCoeff>_precision)
292 MatShift(_A,_heatTransfertCoeff);//Contribution from the liquit/solid heat transfer
294 if(_verbose or _system)
295 MatView(_A,PETSC_VIEWER_STDOUT_WORLD);
300 double StationaryDiffusionEquation::computeDiffusionMatrixFE(bool & stop){
307 Matrix M(_Ndim+1,_Ndim+1);//cell geometry matrix
308 std::vector< Vector > GradShapeFuncs(_Ndim+1);//shape functions of cell nodes
309 std::vector< int > nodeIds(_Ndim+1);//cell node Ids
310 std::vector< Node > nodes(_Ndim+1);//cell nodes
311 int i_int, j_int; //index of nodes j and k considered as unknown nodes
312 bool dirichletCell_treated;
314 std::vector< vector< double > > values(_Ndim+1,vector< double >(_Ndim+1,0));//values of shape functions on cell node
315 for (int idim=0; idim<_Ndim+1;idim++)
316 values[idim][idim]=1;
318 /* parameters for boundary treatment */
319 vector< double > valuesBorder(_Ndim+1);
320 Vector GradShapeFuncBorder(_Ndim+1);
322 for (int j=0; j<_Nmailles;j++)
324 Cj = _mesh.getCell(j);
326 for (int idim=0; idim<_Ndim+1;idim++){
327 nodeIds[idim]=Cj.getNodeId(idim);
328 nodes[idim]=_mesh.getNode(nodeIds[idim]);
329 for (int jdim=0; jdim<_Ndim;jdim++)
330 M(idim,jdim)=nodes[idim].getPoint()[jdim];
333 for (int idim=0; idim<_Ndim+1;idim++)
334 GradShapeFuncs[idim]=DiffusionEquation::gradientNodal(M,values[idim])/DiffusionEquation::fact(_Ndim);
336 /* Loop on the edges of the cell */
337 for (int idim=0; idim<_Ndim+1;idim++)
339 if(find(_dirichletNodeIds.begin(),_dirichletNodeIds.end(),nodeIds[idim])==_dirichletNodeIds.end())//!_mesh.isBorderNode(nodeIds[idim])
340 {//First node of the edge is not Dirichlet node
341 i_int=DiffusionEquation::unknownNodeIndex(nodeIds[idim], _dirichletNodeIds);//assumes Dirichlet boundary node numbering is strictly increasing
342 dirichletCell_treated=false;
343 for (int jdim=0; jdim<_Ndim+1;jdim++)
345 if(find(_dirichletNodeIds.begin(),_dirichletNodeIds.end(),nodeIds[jdim])==_dirichletNodeIds.end())//!_mesh.isBorderNode(nodeIds[jdim])
346 {//Second node of the edge is not Dirichlet node
347 j_int= DiffusionEquation::unknownNodeIndex(nodeIds[jdim], _dirichletNodeIds);//assumes Dirichlet boundary node numbering is strictly increasing
348 MatSetValue(_A,i_int,j_int,(_DiffusionTensor*GradShapeFuncs[idim])*GradShapeFuncs[jdim]/Cj.getMeasure(), ADD_VALUES);
350 else if (!dirichletCell_treated)
351 {//Second node of the edge is a Dirichlet node
352 dirichletCell_treated=true;
353 for (int kdim=0; kdim<_Ndim+1;kdim++)
355 std::map<int,double>::iterator it=_dirichletBoundaryValues.find(nodeIds[kdim]);
356 if( it != _dirichletBoundaryValues.end() )
358 if( _dirichletValuesSet )
359 valuesBorder[kdim]=_dirichletBoundaryValues[it->second];
361 valuesBorder[kdim]=_limitField[_mesh.getNode(nodeIds[kdim]).getGroupName()].T;
364 valuesBorder[kdim]=0;
366 GradShapeFuncBorder=DiffusionEquation::gradientNodal(M,valuesBorder)/DiffusionEquation::fact(_Ndim);
367 coeff =-1.*(_DiffusionTensor*GradShapeFuncBorder)*GradShapeFuncs[idim]/Cj.getMeasure();
368 VecSetValue(_b,i_int,coeff, ADD_VALUES);
375 //Calcul de la contribution de la condition limite de Neumann au second membre
376 if( _NdirichletNodes !=_NboundaryNodes)
378 vector< int > boundaryFaces = _mesh.getBoundaryFaceIds();
379 int NboundaryFaces=boundaryFaces.size();
380 for(int i = 0; i< NboundaryFaces ; i++)//On parcourt les faces du bord
382 Face Fi = _mesh.getFace(boundaryFaces[i]);
383 for(int j = 0 ; j<_Ndim ; j++)//On parcourt les noeuds de la face
385 if(find(_dirichletNodeIds.begin(),_dirichletNodeIds.end(),Fi.getNodeId(j))==_dirichletNodeIds.end())//node j is a Neumann BC node (not a Dirichlet BC node)
387 j_int=DiffusionEquation::unknownNodeIndex(Fi.getNodeId(j), _dirichletNodeIds);//indice du noeud j en tant que noeud inconnu
388 if( _neumannValuesSet )
389 coeff =Fi.getMeasure()/_Ndim*_neumannBoundaryValues[Fi.getNodeId(j)];
391 coeff =Fi.getMeasure()/_Ndim*_limitField[_mesh.getNode(Fi.getNodeId(j)).getGroupName()].normalFlux;
392 VecSetValue(_b, j_int, coeff, ADD_VALUES);
399 if(_onlyNeumannBC) //Check that the matrix is symmetric
401 PetscBool isSymetric;
402 MatIsSymmetric(_A,_precision,&isSymetric);
405 cout<<"Singular matrix is not symmetric, tolerance= "<< _precision<<endl;
406 throw CdmathException("Singular matrix should be symmetric with kernel composed of constant vectors");
410 _diffusionMatrixSet=true;
416 double StationaryDiffusionEquation::computeDiffusionMatrixFV(bool & stop){
419 long nbFaces = _mesh.getNumberOfFaces();
423 double inv_dxi, inv_dxj;
424 double barycenterDistance;
425 Vector normale(_Ndim);
428 std::vector< int > idCells;
430 for (int j=0; j<nbFaces;j++){
431 Fj = _mesh.getFace(j);
433 // compute the normal vector corresponding to face j : from idCells[0] to idCells[1]
434 idCells = Fj.getCellsId();
435 Cell1 = _mesh.getCell(idCells[0]);
437 for(int l=0; l<Cell1.getNumberOfFaces(); l++){
438 if (j == Cell1.getFacesId()[l]){
439 for (int idim = 0; idim < _Ndim; ++idim)
440 normale[idim] = Cell1.getNormalVector(l,idim);
445 //Compute velocity at the face Fj
446 dn=(_DiffusionTensor*normale)*normale;
448 // compute 1/dxi = volume of Ci/area of Fj
449 inv_dxi = Fj.getMeasure()/Cell1.getMeasure();
451 // If Fj is on the boundary
452 if (Fj.getNumberOfCells()==1) {
455 cout << "face numero " << j << " cellule frontiere " << idCells[0] << " ; vecteur normal=(";
456 for(int p=0; p<_Ndim; p++)
457 cout << normale[p] << ",";
461 std::map<int,double>::iterator it=_dirichletBoundaryValues.find(j);
462 if( it != _dirichletBoundaryValues.end() )
464 barycenterDistance=Cell1.getBarryCenter().distance(Fj.getBarryCenter());
465 MatSetValue(_A,idm,idm,dn*inv_dxi/barycenterDistance , ADD_VALUES);
466 VecSetValue(_b,idm, dn*inv_dxi/barycenterDistance*it->second, ADD_VALUES);
470 nameOfGroup = Fj.getGroupName();
472 if (_limitField[nameOfGroup].bcType==NeumannStationaryDiffusion){
473 VecSetValue(_b,idm, -dn*inv_dxi*_limitField[nameOfGroup].normalFlux, ADD_VALUES);
475 else if(_limitField[nameOfGroup].bcType==DirichletStationaryDiffusion){
476 barycenterDistance=Cell1.getBarryCenter().distance(Fj.getBarryCenter());
477 MatSetValue(_A,idm,idm,dn*inv_dxi/barycenterDistance , ADD_VALUES);
478 VecSetValue(_b,idm, dn*inv_dxi/barycenterDistance*_limitField[nameOfGroup].T, ADD_VALUES);
482 cout<<"!!!!!!!!!!!!!!! Error StationaryDiffusionEquation::computeDiffusionMatrixFV !!!!!!!!!!"<<endl;
483 cout<<"!!!!!! No boundary condition set for boundary named "<<nameOfGroup<< "!!!!!!!!!! _limitField[nameOfGroup].bcType= "<<_limitField[nameOfGroup].bcType<<endl;
484 cout<<"Accepted boundary conditions are NeumannStationaryDiffusion "<<NeumannStationaryDiffusion<< " and DirichletStationaryDiffusion "<<DirichletStationaryDiffusion<<endl;
485 *_runLogFile<<"!!!!!! Boundary condition not set for boundary named "<<nameOfGroup<< "!!!!!!!!!! _limitField[nameOfGroup].bcType= "<<_limitField[nameOfGroup].bcType<<endl;
486 _runLogFile->close();
487 throw CdmathException("Boundary condition not set");
490 // if Fj is inside the domain
491 } else if (Fj.getNumberOfCells()==2 ){
494 cout << "face numero " << j << " cellule gauche " << idCells[0] << " cellule droite " << idCells[1];
495 cout << " ; vecteur normal=(";
496 for(int p=0; p<_Ndim; p++)
497 cout << normale[p] << ",";
500 Cell2 = _mesh.getCell(idCells[1]);
503 inv_dxj = Fj.getMeasure()/Cell2.getMeasure();
505 inv_dxj = 1/Cell2.getMeasure();
507 barycenterDistance=Cell1.getBarryCenter().distance(Cell2.getBarryCenter());
509 MatSetValue(_A,idm,idm, dn*inv_dxi/barycenterDistance, ADD_VALUES);
510 MatSetValue(_A,idm,idn,-dn*inv_dxi/barycenterDistance, ADD_VALUES);
511 MatSetValue(_A,idn,idn, dn*inv_dxj/barycenterDistance, ADD_VALUES);
512 MatSetValue(_A,idn,idm,-dn*inv_dxj/barycenterDistance, ADD_VALUES);
516 *_runLogFile<<"StationaryDiffusionEquation::computeDiffusionMatrixFV(): incompatible number of cells around a face"<<endl;
517 throw CdmathException("StationaryDiffusionEquation::computeDiffusionMatrixFV(): incompatible number of cells around a face");
522 if(_onlyNeumannBC) //Check that the matrix is symmetric
524 PetscBool isSymetric;
525 MatIsSymmetric(_A,_precision,&isSymetric);
528 cout<<"Singular matrix is not symmetric, tolerance= "<< _precision<<endl;
529 throw CdmathException("Singular matrix should be symmetric with kernel composed of constant vectors");
533 _diffusionMatrixSet=true;
539 double StationaryDiffusionEquation::computeRHS(bool & stop)//Contribution of the PDE RHS to the linear systemm RHS (boundary conditions do contribute to the system RHS via the function computeDiffusionMatrix
545 for (int i=0; i<_Nmailles;i++)
546 VecSetValue(_b,i, _heatTransfertCoeff*_fluidTemperatureField(i) + _heatPowerField(i) ,ADD_VALUES);
550 std::vector< int > nodesId;
551 for (int i=0; i<_Nmailles;i++)
554 nodesId=Ci.getNodesId();
555 for (int j=0; j<nodesId.size();j++)
556 if(!_mesh.isBorderNode(nodesId[j]))
558 double coeff = _heatTransfertCoeff*_fluidTemperatureField(nodesId[j]) + _heatPowerField(nodesId[j]);
559 VecSetValue(_b,DiffusionEquation::unknownNodeIndex(nodesId[j], _dirichletNodeIds), coeff*Ci.getMeasure()/(_Ndim+1),ADD_VALUES);
564 VecAssemblyBegin(_b);
567 if(_verbose or _system)
568 VecView(_b,PETSC_VIEWER_STDOUT_WORLD);
574 bool StationaryDiffusionEquation::iterateNewtonStep(bool &converged)
578 //Only implicit scheme considered
580 #if PETSC_VERSION_GREATER_3_5
581 KSPSetOperators(_ksp, _A, _A);
583 KSPSetOperators(_ksp, _A, _A,SAME_NONZERO_PATTERN);
587 KSPSetComputeEigenvalues(_ksp,PETSC_TRUE);
588 KSPSolve(_ksp, _b, _Tk);
590 KSPConvergedReason reason;
591 KSPGetConvergedReason(_ksp,&reason);
592 KSPGetIterationNumber(_ksp, &_PetscIts);
594 KSPGetResidualNorm(_ksp,&residu);
595 if (reason!=2 and reason!=3)
597 PetscPrintf(PETSC_COMM_WORLD,"!!!!!!!!!!!!! Erreur système linéaire : pas de convergence de Petsc.");
598 PetscPrintf(PETSC_COMM_WORLD,"!!!!!!!!!!!!! Itérations maximales %d atteintes, résidu = %1.2f, précision demandée= %1.2f",_maxPetscIts,residu,_precision);
599 PetscPrintf(PETSC_COMM_WORLD,"Solver used %s, preconditioner %s, Final number of iteration = %d",_ksptype,_pctype,_PetscIts);
600 *_runLogFile<<"!!!!!!!!!!!!! Erreur système linéaire : pas de convergence de Petsc."<<endl;
601 *_runLogFile<<"!!!!!!!!!!!!! Itérations maximales "<<_maxPetscIts<<" atteintes, résidu="<<residu<<", précision demandée= "<<_precision<<endl;
602 *_runLogFile<<"Solver used "<< _ksptype<<", preconditioner "<<_pctype<<", Final number of iteration= "<<_PetscIts<<endl;
603 _runLogFile->close();
608 if( _MaxIterLinearSolver < _PetscIts)
609 _MaxIterLinearSolver = _PetscIts;
610 PetscPrintf(PETSC_COMM_WORLD,"## Système linéaire résolu en %d itérations par le solveur %s et le preconditioneur %s, précision demandée = %1.2f",_PetscIts,_ksptype,_pctype,_precision);
611 *_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;
612 VecCopy(_Tk, _deltaT);//ici on a deltaT=Tk
613 VecAXPY(_deltaT, -1, _Tkm1);//On obtient deltaT=Tk-Tkm1
615 VecAssemblyBegin(_Tk);
616 VecAssemblyEnd( _Tk);
617 VecAssemblyBegin(_deltaT);
618 VecAssemblyEnd( _deltaT);
621 PetscPrintf(PETSC_COMM_WORLD,"Début calcul de l'erreur maximale");
623 VecNorm(_deltaT,NORM_INFINITY,&_erreur_rel);
626 PetscPrintf(PETSC_COMM_WORLD,"Fin calcul de la variation relative, erreur maximale : %1.2f", _erreur_rel );
628 converged = (_erreur_rel <= _precision) ;//converged=convergence des iterations de Newton
636 void StationaryDiffusionEquation::setMesh(const Mesh &M)
640 if(_Ndim != M.getSpaceDimension() or _Ndim!=M.getMeshDimension())//for the moment we must have space dim=mesh dim
642 cout<< "Problem : dimension defined is "<<_Ndim<< " but mesh dimension= "<<M.getMeshDimension()<<", and space dimension is "<<M.getSpaceDimension()<<endl;
643 *_runLogFile<< "Problem : dim = "<<_Ndim<< " but mesh dim= "<<M.getMeshDimension()<<", mesh space dim= "<<M.getSpaceDimension()<<endl;
644 *_runLogFile<<"StationaryDiffusionEquation::setMesh: mesh has incorrect dimension"<<endl;
645 _runLogFile->close();
646 throw CdmathException("StationaryDiffusionEquation::setMesh: mesh has incorrect dimension");
650 _Nmailles = _mesh.getNumberOfCells();
651 _Nnodes = _mesh.getNumberOfNodes();
653 cout<<"Mesh has "<< _Nmailles << " cells and " << _Nnodes << " nodes"<<endl<<endl;;
654 *_runLogFile<<"Mesh has "<< _Nmailles << " cells and " << _Nnodes << " nodes"<<endl<<endl;
656 // find maximum nb of neibourghs
659 _VV=Field ("Temperature", CELLS, _mesh, 1);
660 _neibMaxNbCells=_mesh.getMaxNbNeighbours(CELLS);
664 if(_Ndim==1 )//The 1D cdmath mesh is necessarily made of segments
665 cout<<"1D Finite element method on segments"<<endl;
668 if( _mesh.isTriangular() )//Mesh dim=2
669 cout<<"2D Finite element method on triangles"<<endl;
670 else if (_mesh.getMeshDimension()==1)//Mesh dim=1
671 cout<<"1D Finite element method on a 2D network : space dimension is "<<_Ndim<< ", mesh dimension is "<<_mesh.getMeshDimension()<<endl;
674 cout<<"Error Finite element with Space dimension "<<_Ndim<< ", and mesh dimension "<<_mesh.getMeshDimension()<< ", mesh should be either triangular either 1D network"<<endl;
675 *_runLogFile<<"StationaryDiffusionEquation::setMesh: mesh has incorrect dimension"<<endl;
676 _runLogFile->close();
677 throw CdmathException("StationaryDiffusionEquation::setMesh: mesh has incorrect cell types");
682 if( _mesh.isTetrahedral() )//Mesh dim=3
683 cout<<"3D Finite element method on tetrahedra"<<endl;
684 else if (_mesh.getMeshDimension()==2 and _mesh.isTriangular())//Mesh dim=2
685 cout<<"2D Finite element method on a 3D surface : space dimension is "<<_Ndim<< ", mesh dimension is "<<_mesh.getMeshDimension()<<endl;
686 else if (_mesh.getMeshDimension()==1)//Mesh dim=1
687 cout<<"1D Finite element method on a 3D network : space dimension is "<<_Ndim<< ", mesh dimension is "<<_mesh.getMeshDimension()<<endl;
690 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;
691 *_runLogFile<<"StationaryDiffusionEquation::setMesh: mesh has incorrect dimension"<<endl;
692 _runLogFile->close();
693 throw CdmathException("StationaryDiffusionEquation::setMesh: mesh has incorrect cell types");
697 _VV=Field ("Temperature", NODES, _mesh, 1);
699 _neibMaxNbNodes=_mesh.getMaxNbNeighbours(NODES);
700 _boundaryNodeIds = _mesh.getBoundaryNodeIds();
701 _NboundaryNodes=_boundaryNodeIds.size();
705 /* MPI distribution parameters */
706 int nbVoisinsMax;//Mettre en attribut ?
708 MPI_Bcast(&_Nmailles , 1, MPI_INT, 0, PETSC_COMM_WORLD);
709 MPI_Bcast(&_neibMaxNbCells, 1, MPI_INT, 0, PETSC_COMM_WORLD);
710 nbVoisinsMax = _neibMaxNbCells;
713 MPI_Bcast(&_Nnodes , 1, MPI_INT, 0, PETSC_COMM_WORLD);
714 MPI_Bcast(&_neibMaxNbNodes, 1, MPI_INT, 0, PETSC_COMM_WORLD);
715 nbVoisinsMax = _neibMaxNbNodes;
717 _d_nnz = (nbVoisinsMax+1)*_nVar;
718 _o_nnz = nbVoisinsMax *_nVar;
723 void StationaryDiffusionEquation::setLinearSolver(linearSolver kspType, preconditioner pcType)
725 //_maxPetscIts=maxIterationsPetsc;
726 // set linear solver algorithm
728 _ksptype = (char*)&KSPGMRES;
729 else if (kspType==CG)
730 _ksptype = (char*)&KSPCG;
731 else if (kspType==BCGS)
732 _ksptype = (char*)&KSPBCGS;
734 cout << "!!! Error : only 'GMRES', 'CG' or 'BCGS' is acceptable as a linear solver !!!" << endl;
735 *_runLogFile << "!!! Error : only 'GMRES', 'CG' or 'BCGS' is acceptable as a linear solver !!!" << endl;
736 _runLogFile->close();
737 throw CdmathException("!!! Error : only 'GMRES', 'CG' or 'BCGS' algorithm is acceptable !!!");
739 // set preconditioner
741 _pctype = (char*)&PCNONE;
742 else if (pcType ==LU)
743 _pctype = (char*)&PCLU;
744 else if (pcType == ILU)
745 _pctype = (char*)&PCILU;
746 else if (pcType ==CHOLESKY)
747 _pctype = (char*)&PCCHOLESKY;
748 else if (pcType == ICC)
749 _pctype = (char*)&PCICC;
751 cout << "!!! Error : only 'NOPC', 'LU', 'ILU', 'CHOLESKY' or 'ICC' preconditioners are acceptable !!!" << endl;
752 *_runLogFile << "!!! Error : only 'NOPC' or 'LU' or 'ILU' preconditioners are acceptable !!!" << endl;
753 _runLogFile->close();
754 throw CdmathException("!!! Error : only 'NOPC' or 'LU' or 'ILU' preconditioners are acceptable !!!" );
758 bool StationaryDiffusionEquation::solveStationaryProblem()
760 if(!_initializedMemory)
762 *_runLogFile<< "ProblemCoreFlows::run() call initialize() first"<< _fileName<<endl;
763 _runLogFile->close();
764 throw CdmathException("ProblemCoreFlows::run() call initialize() first");
766 bool stop=false; // Does the Problem want to stop (error) ?
767 bool converged=false; // has the newton scheme converged (end) ?
769 PetscPrintf(PETSC_COMM_WORLD,"!!! Running test case %s using ",_fileName);
770 *_runLogFile<< "!!! Running test case "<< _fileName<< " using ";
774 PetscPrintf(PETSC_COMM_WORLD,"Finite volumes method\n\n");
775 *_runLogFile<< "Finite volumes method"<<endl<<endl;
779 PetscPrintf(PETSC_COMM_WORLD,"Finite elements method\n\n");
780 *_runLogFile<< "Finite elements method"<< endl<<endl;
783 computeDiffusionMatrix( stop);//For the moment the conductivity does not depend on the temperature (linear LHS)
785 PetscPrintf(PETSC_COMM_WORLD,"Error : failed computing diffusion matrix, stopping calculation\n");
786 *_runLogFile << "Error : failed computing diffusion matrix, stopping calculation"<< endl;
787 _runLogFile->close();
788 throw CdmathException("Failed computing diffusion matrix");
790 computeRHS(stop);//For the moment the heat power does not depend on the unknown temperature (linear RHS)
792 PetscPrintf(PETSC_COMM_WORLD,"Error : failed computing right hand side, stopping calculation\n");
793 *_runLogFile << "Error : failed computing right hand side, stopping calculation"<< endl;
794 throw CdmathException("Failed computing right hand side");
796 stop = iterateNewtonStep(converged);
798 PetscPrintf(PETSC_COMM_WORLD,"Error : failed solving linear system, stopping calculation\n");
799 *_runLogFile << "Error : failed linear system, stopping calculation"<< endl;
800 _runLogFile->close();
801 throw CdmathException("Failed solving linear system");
804 _computationCompletedSuccessfully=true;
807 // Newton iteration loop for non linear problems
809 while(!stop and !converged and _NEWTON_its<_maxNewtonIts)
811 computeDiffusionMatrix( stop);//case when the conductivity depends on the temperature (nonlinear LHS)
812 computeRHS(stop);//case the heat power depends on the unknown temperature (nonlinear RHS)
813 stop = iterateNewtonStep(converged);
817 cout << "Error : failed solving Newton iteration "<<_NEWTON_its<<", stopping calculation"<< endl;
818 *_runLogFile << "Error : failed solving Newton iteration "<<_NEWTON_its<<", stopping calculation"<< endl;
819 throw CdmathException("Failed solving a Newton iteration");
821 else if(_NEWTON_its==_maxNewtonIts){
822 cout << "Error : no convergence of Newton scheme. Maximum Newton iterations "<<_maxNewtonIts<<" reached, stopping calculation"<< endl;
823 *_runLogFile << "Error : no convergence of Newton scheme. Maximum Newton iterations "<<_maxNewtonIts<<" reached, stopping calculation"<< endl;
824 throw CdmathException("No convergence of Newton scheme");
827 cout << "Convergence of Newton scheme at iteration "<<_NEWTON_its<<", end of calculation"<< endl;
828 *_runLogFile << "Convergence of Newton scheme at iteration "<<_NEWTON_its<<", end of calculation"<< endl;
833 *_runLogFile<< "!!!!!! Computation successful !!!!!!"<< endl;
834 _runLogFile->close();
839 void StationaryDiffusionEquation::save(){
840 PetscPrintf(PETSC_COMM_WORLD,"Saving numerical results\n\n");
841 *_runLogFile<< "Saving numerical results"<< endl<<endl;
843 string resultFile(_path+"/StationaryDiffusionEquation");//Results
846 resultFile+=_fileName;
848 // create mesh and component info
849 string suppress ="rm -rf "+resultFile+"_*";
850 system(suppress.c_str());//Nettoyage des précédents calculs identiques
853 VecScatterBegin(_scat,_Tk,_Tk_seq,INSERT_VALUES,SCATTER_FORWARD);
854 VecScatterEnd( _scat,_Tk,_Tk_seq,INSERT_VALUES,SCATTER_FORWARD);
857 if(_verbose or _system)
858 VecView(_Tk,PETSC_VIEWER_STDOUT_WORLD);
863 for(int i=0; i<_Nmailles; i++)
866 VecGetValues(_Tk_seq, 1, &i, &Ti);
868 VecGetValues(_Tk , 1, &i, &Ti);
874 for(int i=0; i<_NunknownNodes; i++)
877 VecGetValues(_Tk_seq, 1, &i, &Ti);
879 VecGetValues(_Tk , 1, &i, &Ti);
880 globalIndex = DiffusionEquation::globalNodeIndex(i, _dirichletNodeIds);
886 for(int i=0; i<_NdirichletNodes; i++)
888 Ni=_mesh.getNode(_dirichletNodeIds[i]);
889 nameOfGroup = Ni.getGroupName();
890 _VV(_dirichletNodeIds[i])=_limitField[nameOfGroup].T;
894 _VV.setInfoOnComponent(0,"Temperature_(K)");
898 _VV.writeVTK(resultFile);
901 _VV.writeMED(resultFile);
904 _VV.writeCSV(resultFile);
910 void StationaryDiffusionEquation::terminate()
914 VecDestroy(&_deltaT);
917 if(_mpi_size>1 && _mpi_rank == 0)
918 VecDestroy(&_Tk_seq);
920 PetscBool petscInitialized;
921 PetscInitialized(&petscInitialized);
928 StationaryDiffusionEquation::setDirichletValues(map< int, double> dirichletBoundaryValues)
930 _dirichletValuesSet=true;
931 _dirichletBoundaryValues=dirichletBoundaryValues;
935 StationaryDiffusionEquation::setNeumannValues(map< int, double> neumannBoundaryValues)
937 _neumannValuesSet=true;
938 _neumannBoundaryValues=neumannBoundaryValues;
942 StationaryDiffusionEquation::getConditionNumber(bool isSingular, double tol) const
944 SparseMatrixPetsc A = SparseMatrixPetsc(_A);
945 return A.getConditionNumber( isSingular, tol);
947 std::vector< double >
948 StationaryDiffusionEquation::getEigenvalues(int nev, EPSWhich which, double tol) const
950 SparseMatrixPetsc A = SparseMatrixPetsc(_A);
952 if(_FECalculation)//We need to scale the FE matrix, otherwise the eigenvalues go to zero as the mesh is refined
954 Vector nodal_volumes(_NunknownNodes);
956 for(int i = 0; i< _Nmailles ; i++)//On parcourt les cellules du maillage
958 Cell Ci = _mesh.getCell(i);
959 for(int j = 0 ; j<_Ndim+1 ; j++)//On parcourt les noeuds de la cellule
961 if(find(_dirichletNodeIds.begin(),_dirichletNodeIds.end(),Ci.getNodeId(j))==_dirichletNodeIds.end())//node j is an unknown node (not a Dirichlet node)
963 j_int=DiffusionEquation::unknownNodeIndex(Ci.getNodeId(j), _dirichletNodeIds);//indice du noeud j en tant que noeud inconnu
964 nodal_volumes[j_int]+=Ci.getMeasure()/(_Ndim+1);
968 for( j_int = 0; j_int< _NunknownNodes ; j_int++)
969 nodal_volumes[j_int]=1/nodal_volumes[j_int];
970 A.leftDiagonalScale(nodal_volumes);
973 return A.getEigenvalues( nev, which, tol);
975 std::vector< Vector >
976 StationaryDiffusionEquation::getEigenvectors(int nev, EPSWhich which, double tol) const
978 SparseMatrixPetsc A = SparseMatrixPetsc(_A);
979 return A.getEigenvectors( nev, which, tol);
982 StationaryDiffusionEquation::getEigenvectorsField(int nev, EPSWhich which, double tol) const
984 SparseMatrixPetsc A = SparseMatrixPetsc(_A);
985 MEDCoupling::DataArrayDouble * d = A.getEigenvectorsDataArrayDouble( nev, which, tol);
988 my_eigenfield = Field("Eigenvectors field", _VV.getTypeOfField(), _mesh, nev);
990 my_eigenfield.setFieldByDataArrayDouble(d);
992 return my_eigenfield;
996 StationaryDiffusionEquation::getOutputTemperatureField()
998 if(!_computationCompletedSuccessfully)
999 throw("Computation not performed yet or failed. No temperature field available");
1005 StationaryDiffusionEquation::getRodTemperatureField()
1007 return getOutputTemperatureField();
1011 StationaryDiffusionEquation::getInputFieldsNames()
1013 vector<string> result(2);
1015 result[0]="FluidTemperature";
1016 result[1]="HeatPower";
1021 StationaryDiffusionEquation::getOutputFieldsNames()
1023 vector<string> result(1);
1025 result[0]="RodTemperature";
1031 StationaryDiffusionEquation::getOutputField(const string& nameField )
1033 if(nameField=="RodTemperature" || nameField=="RODTEMPERATURE" || nameField=="TEMPERATURECOMBUSTIBLE" || nameField=="TemperatureCombustible" )
1034 return getRodTemperatureField();
1037 cout<<"Error : Field name "<< nameField << " does not exist, call getOutputFieldsNames first" << endl;
1038 throw CdmathException("DiffusionEquation::getOutputField error : Unknown Field name");
1043 StationaryDiffusionEquation::setInputField(const string& nameField, Field& inputField )
1045 if(nameField=="FluidTemperature" || nameField=="FLUIDTEMPERATURE" || nameField=="TemperatureFluide" || nameField=="TEMPERATUREFLUIDE")
1046 return setFluidTemperatureField( inputField) ;
1047 else if(nameField=="HeatPower" || nameField=="HEATPOWER" || nameField=="PuissanceThermique" || nameField=="PUISSANCETHERMIQUE" )
1048 return setHeatPowerField( inputField );
1051 cout<<"Error : Field name "<< nameField << " is not an input field name, call getInputFieldsNames first" << endl;
1052 throw CdmathException("StationaryDiffusionEquation::setInputField error : Unknown Field name");