1 //============================================================================
3 * \file StationaryDiffusionEquation.hxx
4 * \author Michael NDJINGA
7 * \brief Stationary heat diffusion equation solved with either finite elements or finite volume method.
8 * -\lambda\Delta T=\Phi + \lambda_{sf} (T_{fluid}-T)
9 * Dirichlet (imposed temperature) or Neumann (imposed normal flux) boundary conditions.
11 //============================================================================
13 /*! \class StationaryDiffusionEquation StationaryDiffusionEquation.hxx "StationaryDiffusionEquation.hxx"
14 * \brief Scalar stationary heat equation solved with either finite elements or finite volume method.
15 * \details see \ref StationaryDiffusionEqPage for more details
16 * -\lambda\Delta T=\Phi(T) + \lambda_{sf} (T_{fluid}-T)
19 #ifndef StationaryDiffusionEquation_HXX_
20 #define StationaryDiffusionEquation_HXX_
22 #include "ProblemCoreFlows.hxx"
26 /*! Boundary condition type */
27 enum BoundaryTypeStationaryDiffusion { NeumannStationaryDiffusion, DirichletStationaryDiffusion, NoneBCStationaryDiffusion};
29 /** \struct LimitField
30 * \brief value of some fields on the boundary */
31 struct LimitFieldStationaryDiffusion{
32 LimitFieldStationaryDiffusion(){bcType=NoneBCStationaryDiffusion; T=0; normalFlux=0;}
33 LimitFieldStationaryDiffusion(BoundaryTypeStationaryDiffusion _bcType, double _T, double _normalFlux){
34 bcType=_bcType; T=_T; normalFlux=_normalFlux;
37 BoundaryTypeStationaryDiffusion bcType;
38 double T; //for Dirichlet
39 double normalFlux; //for Neumann
42 class StationaryDiffusionEquation
46 /** \fn StationaryDiffusionEquation
47 * \brief Constructor for the temperature diffusion in a solid
48 * \param [in] int : space dimension
49 * \param [in] bool : numerical method
50 * \param [in] double : solid conductivity
53 StationaryDiffusionEquation( int dim,bool FECalculation=true,double lambda=1,MPI_Comm comm = MPI_COMM_WORLD);
55 void setMesh(const Mesh &M);
56 void setFileName(string fileName){
59 bool solveStationaryProblem();
60 Field getOutputTemperatureField();
62 //Linear system and spectrum
63 void setLinearSolver(linearSolver kspType, preconditioner pcType);
64 double getConditionNumber(bool isSingular=false, double tol=1e-6) const;
65 std::vector< double > getEigenvalues (int nev, EPSWhich which=EPS_SMALLEST_MAGNITUDE, double tol=1e-6) const;
66 std::vector< Vector > getEigenvectors(int nev, EPSWhich which=EPS_SMALLEST_MAGNITUDE, double tol=1e-6) const;
67 Field getEigenvectorsField(int nev, EPSWhich which=EPS_SMALLEST_MAGNITUDE, double tol=1e-6) const;
71 void terminate();//vide la mémoire et enregistre le résultat final
72 double computeDiffusionMatrix(bool & stop);
73 double computeTimeStep(bool & stop);//For coupling calculations
74 bool iterateNewtonStep(bool &ok);
77 /* Boundary conditions */
78 void setBoundaryFields(map<string, LimitFieldStationaryDiffusion> boundaryFields){
79 _limitField = boundaryFields;
81 /** \fn setDirichletBoundaryCondition
82 * \brief adds a new boundary condition of type Dirichlet
84 * \param [in] string : the name of the boundary
85 * \param [in] double : the value of the temperature at the boundary
88 void setDirichletBoundaryCondition(string groupName,double Temperature){
89 _limitField[groupName]=LimitFieldStationaryDiffusion(DirichletStationaryDiffusion,Temperature,-1);
92 /** \fn setNeumannBoundaryCondition
93 * \brief adds a new boundary condition of type Neumann
95 * \param [in] string : the name of the boundary
96 * \param [in] double : outward normal flux
99 void setNeumannBoundaryCondition(string groupName, double normalFlux=0){
100 _limitField[groupName]=LimitFieldStationaryDiffusion(NeumannStationaryDiffusion,-1, normalFlux);
103 void setDirichletValues(map< int, double> dirichletBoundaryValues);
104 void setNeumannValues (map< int, double> neumannBoundaryValues);
106 void setConductivity(double conductivite){
107 _conductivity=conductivite;
109 void setFluidTemperatureField(Field coupledTemperatureField){
110 _fluidTemperatureField=coupledTemperatureField;
111 _fluidTemperatureFieldSet=true;
113 void setFluidTemperature(double fluidTemperature){
114 _fluidTemperature=fluidTemperature;
116 Field& getRodTemperatureField(){
119 Field& getFluidTemperatureField(){
120 return _fluidTemperatureField;
122 void setDiffusiontensor(Matrix DiffusionTensor){
123 _DiffusionTensor=DiffusionTensor;
126 /** \fn setHeatPowerField
127 * \brief set the heat power field (variable in space)
132 void setHeatPowerField(Field heatPower){
133 _heatPowerField=heatPower;
134 _heatPowerFieldSet=true;
137 /** \fn setHeatPowerField
138 * \brief set the heat power field (variable in space)
140 * \param [in] string fileName (including file path)
141 * \param [in] string fieldName
144 void setHeatPowerField(string fileName, string fieldName){
145 _heatPowerField=Field(fileName, CELLS,fieldName);
146 _heatPowerFieldSet=true;
150 * \brief Updates display options
156 void setVerbose(bool verbose, bool system=false)
166 int _Ndim;//space dimension
167 int _nVar;//Number of equations to solve=1
172 bool _initializedMemory;
173 int _Nmailles;//number of cells for FV calculation
174 int _neibMaxNbCells;//maximum number of cells around a cell
177 double _precision_Newton;
178 double _erreur_rel;//norme(Uk+1-Uk)
179 bool _computationCompletedSuccessfully;
181 //Linear solver and petsc
187 int _maxPetscIts;//nombre maximum d'iteration gmres autorisé au cours d'une resolution de système lineaire
188 int _maxNewtonIts;//nombre maximum d'iteration de Newton autorise au cours de la resolution d'un pas de temps
189 int _PetscIts;//the number of iterations of the linear solver
191 Mat _A;//Linear system matrix
192 Vec _b;//Linear system right hand side
193 double _MaxIterLinearSolver;//nombre maximum d'iteration gmres obtenu au cours par les resolution de systemes lineaires au cours d'un pas de tmeps
194 bool _conditionNumber;//computes an estimate of the condition number
196 map<string, LimitFieldStationaryDiffusion> _limitField;
197 bool _onlyNeumannBC;//if true then the linear system is singular and should be solved up to a constant vector
199 bool _diffusionMatrixSet;
201 Matrix _DiffusionTensor;
202 Vec _deltaT, _Tk, _Tkm1, _b0;
205 //Heat transfert variables
206 double _conductivity, _fluidTemperature;
207 Field _heatPowerField, _fluidTemperatureField;
208 bool _heatPowerFieldSet, _fluidTemperatureFieldSet;
209 double _heatTransfertCoeff, _heatSource;
212 bool _verbose, _system;
213 ofstream * _runLogFile;//for creation of a log file to save the history of the simulation
215 string _fileName;//name of the calculation
216 string _path;//path to execution directory used for saving results
217 saveFormat _saveFormat;//file saving format : MED, VTK or CSV
219 double computeRHS(bool & stop);
220 double computeDiffusionMatrixFV(bool & stop);
221 double computeDiffusionMatrixFE(bool & stop);
223 /************ Data for FE calculation *************/
225 int _Nnodes;/* number of nodes for FE calculation */
226 int _neibMaxNbNodes;/* maximum number of nodes around a node */
227 int _NunknownNodes;/* number of unknown nodes for FE calculation */
228 int _NboundaryNodes;/* total number of boundary nodes */
229 int _NdirichletNodes;/* number of boundary nodes with Dirichlet BC for FE calculation */
230 std::vector< int > _boundaryNodeIds;/* List of boundary nodes */
231 std::vector< int > _dirichletNodeIds;/* List of boundary nodes with Dirichlet BC */
233 /********* Possibility to set a boundary field as DirichletNeumann boundary condition *********/
234 bool _dirichletValuesSet;
235 bool _neumannValuesSet;
236 std::map< int, double> _dirichletBoundaryValues;
237 std::map< int, double> _neumannBoundaryValues;
240 PetscMPIInt _size; /* size of communicator */
241 PetscMPIInt _rank; /* processor rank */
244 #endif /* StationaryDiffusionEquation_HXX_ */