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();
61 //Linear system and spectrum
62 void setLinearSolver(linearSolver kspType, preconditioner pcType);
63 double getConditionNumber(bool isSingular=false, double tol=1e-6) const;
64 std::vector< double > getEigenvalues (int nev, EPSWhich which=EPS_SMALLEST_MAGNITUDE, double tol=1e-6) const;
65 std::vector< Vector > getEigenvectors(int nev, EPSWhich which=EPS_SMALLEST_MAGNITUDE, double tol=1e-6) const;
66 Field getEigenvectorsField(int nev, EPSWhich which=EPS_SMALLEST_MAGNITUDE, double tol=1e-6) const;
70 void terminate();//vide la mémoire et enregistre le résultat final
71 double computeDiffusionMatrix(bool & stop);
72 double computeTimeStep(bool & stop);//For coupling calculations
73 bool iterateNewtonStep(bool &ok);
76 /* Boundary conditions */
77 void setBoundaryFields(map<string, LimitFieldStationaryDiffusion> boundaryFields){
78 _limitField = boundaryFields;
80 /** \fn setDirichletBoundaryCondition
81 * \brief adds a new boundary condition of type Dirichlet
83 * \param [in] string : the name of the boundary
84 * \param [in] double : the value of the temperature at the boundary
87 void setDirichletBoundaryCondition(string groupName,double Temperature){
88 _limitField[groupName]=LimitFieldStationaryDiffusion(DirichletStationaryDiffusion,Temperature,-1);
91 /** \fn setNeumannBoundaryCondition
92 * \brief adds a new boundary condition of type Neumann
94 * \param [in] string : the name of the boundary
95 * \param [in] double : outward normal flux
98 void setNeumannBoundaryCondition(string groupName, double normalFlux=0){
99 _limitField[groupName]=LimitFieldStationaryDiffusion(NeumannStationaryDiffusion,-1, normalFlux);
102 void setDirichletValues(map< int, double> dirichletBoundaryValues);
103 void setNeumannValues (map< int, double> neumannBoundaryValues);
105 void setConductivity(double conductivite){
106 _conductivity=conductivite;
108 void setDiffusiontensor(Matrix DiffusionTensor){
109 _DiffusionTensor=DiffusionTensor;
113 //get input fields to prepare the simulation or coupling
114 vector<string> getInputFieldsNames();
115 void setInputField(const string& nameField, Field& inputField );//supply of a required input field
117 void setFluidTemperatureField(Field coupledTemperatureField);
118 void setFluidTemperature(double fluidTemperature){ _fluidTemperature=fluidTemperature; }
119 Field& getFluidTemperatureField(){ return _fluidTemperatureField; }
121 /** \fn setHeatPowerField
122 * \brief set the heat power field (variable in space)
127 void setHeatPowerField(Field heatPower);
129 /** \fn setHeatPowerField
130 * \brief set the heat power field (variable in space)
132 * \param [in] string fileName (including file path)
133 * \param [in] string fieldName
136 void setHeatPowerField(string fileName, string fieldName, int iteration = 0, int order = 0, int meshLevel=0);
138 /** \fn getHeatPowerField
139 * \brief returns the heat power field
144 Field getHeatPowerField(){
145 return _heatPowerField;
147 //get output fields names for postprocessing or coupling
148 vector<string> getOutputFieldsNames() ;//liste tous les champs que peut fournir le code pour le postraitement
149 Field& getOutputField(const string& nameField );//Renvoie un champs pour le postraitement
151 Field& getOutputTemperatureField();
152 Field& getRodTemperatureField();
155 * \brief Updates display options
161 void setVerbose(bool verbose, bool system=false)
171 int _Ndim;//space dimension
172 int _nVar;//Number of equations to solve=1
177 bool _initializedMemory;
178 int _Nmailles;//number of cells for FV calculation
179 int _neibMaxNbCells;//maximum number of cells around a cell
182 double _precision_Newton;
183 double _erreur_rel;//norme(Uk+1-Uk)
184 bool _computationCompletedSuccessfully;
186 //Linear solver and petsc
192 int _maxPetscIts;//nombre maximum d'iteration gmres autorisé au cours d'une resolution de système lineaire
193 int _maxNewtonIts;//nombre maximum d'iteration de Newton autorise au cours de la resolution d'un pas de temps
194 int _PetscIts;//the number of iterations of the linear solver
196 Mat _A;//Linear system matrix
197 Vec _b;//Linear system right hand side
198 double _MaxIterLinearSolver;//nombre maximum d'iteration gmres obtenu au cours par les resolution de systemes lineaires au cours d'un pas de tmeps
199 bool _conditionNumber;//computes an estimate of the condition number
201 map<string, LimitFieldStationaryDiffusion> _limitField;
202 bool _onlyNeumannBC;//if true then the linear system is singular and should be solved up to a constant vector
204 bool _diffusionMatrixSet;
206 Matrix _DiffusionTensor;
207 Vec _deltaT, _Tk, _Tkm1, _b0;
208 Vec _Tk_seq; // Local sequential copy of the parallel vector _Tk, used for saving result files
212 //Heat transfert variables
213 double _conductivity, _fluidTemperature;
214 Field _heatPowerField, _fluidTemperatureField;
215 bool _heatPowerFieldSet, _fluidTemperatureFieldSet;
216 double _heatTransfertCoeff, _heatSource;
219 bool _verbose, _system;
220 ofstream * _runLogFile;//for creation of a log file to save the history of the simulation
222 string _fileName;//name of the calculation
223 string _path;//path to execution directory used for saving results
224 saveFormat _saveFormat;//file saving format : MED, VTK or CSV
226 double computeRHS(bool & stop);
227 double computeDiffusionMatrixFV(bool & stop);
228 double computeDiffusionMatrixFE(bool & stop);
230 /************ Data for FE calculation *************/
232 int _Nnodes;/* number of nodes for FE calculation */
233 int _neibMaxNbNodes;/* maximum number of nodes around a node */
234 int _NunknownNodes;/* number of unknown nodes for FE calculation */
235 int _NboundaryNodes;/* total number of boundary nodes */
236 int _NdirichletNodes;/* number of boundary nodes with Dirichlet BC for FE calculation */
237 std::vector< int > _boundaryNodeIds;/* List of boundary nodes */
238 std::vector< int > _dirichletNodeIds;/* List of boundary nodes with Dirichlet BC */
240 /********* Possibility to set a boundary field as DirichletNeumann boundary condition *********/
241 bool _dirichletValuesSet;
242 bool _neumannValuesSet;
243 std::map< int, double> _dirichletBoundaryValues;
244 std::map< int, double> _neumannBoundaryValues;
246 /**** MPI related variables ***/
247 PetscMPIInt _mpi_size; /* size of communicator */
248 PetscMPIInt _mpi_rank; /* processor rank */
249 VecScatter _scat; /* For the distribution of a local vector */
250 int _globalNbUnknowns, _localNbUnknowns;
251 int _d_nnz, _o_nnz; /* local and "non local" numbers of non zeros corfficients */
254 #endif /* StationaryDiffusionEquation_HXX_ */