]> SALOME platform Git repositories - tools/solverlab.git/blob - CoreFlows/Models/inc/StationaryDiffusionEquation.hxx
Salome HOME
Implemented ICoCo v2
[tools/solverlab.git] / CoreFlows / Models / inc / StationaryDiffusionEquation.hxx
1 //============================================================================
2 /**
3  * \file StationaryDiffusionEquation.hxx
4  * \author Michael NDJINGA
5  * \version 1.0
6  * \date June 2019
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.
10  * */
11 //============================================================================
12
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)
17  */
18  
19 #ifndef StationaryDiffusionEquation_HXX_
20 #define StationaryDiffusionEquation_HXX_
21
22 #include "ProblemCoreFlows.hxx"
23
24 using namespace std;
25
26 /*! Boundary condition type  */
27 enum BoundaryTypeStationaryDiffusion    { NeumannStationaryDiffusion, DirichletStationaryDiffusion, NoneBCStationaryDiffusion};
28
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;
35         }
36
37         BoundaryTypeStationaryDiffusion bcType;
38         double T; //for Dirichlet
39         double normalFlux; //for Neumann
40 };
41
42 class StationaryDiffusionEquation
43 {
44
45 public :
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
51                          *  */
52
53         StationaryDiffusionEquation( int dim,bool FECalculation=true,double lambda=1,MPI_Comm comm = MPI_COMM_WORLD);
54
55     void setMesh(const Mesh &M);
56     void setFileName(string fileName){
57         _fileName = fileName;
58     }
59     bool solveStationaryProblem();
60     
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;
67
68         //Gestion du calcul
69         void initialize();
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);
74         void save();
75
76     /* Boundary conditions */
77         void setBoundaryFields(map<string, LimitFieldStationaryDiffusion> boundaryFields){
78                 _limitField = boundaryFields;
79     };
80         /** \fn setDirichletBoundaryCondition
81                          * \brief adds a new boundary condition of type Dirichlet
82                          * \details
83                          * \param [in] string : the name of the boundary
84                          * \param [in] double : the value of the temperature at the boundary
85                          * \param [out] void
86                          *  */
87         void setDirichletBoundaryCondition(string groupName,double Temperature){
88                 _limitField[groupName]=LimitFieldStationaryDiffusion(DirichletStationaryDiffusion,Temperature,-1);
89         };
90
91         /** \fn setNeumannBoundaryCondition
92                          * \brief adds a new boundary condition of type Neumann
93                          * \details
94                          * \param [in] string : the name of the boundary
95                          * \param [in] double : outward normal flux
96                          * \param [out] void
97                          *  */
98         void setNeumannBoundaryCondition(string groupName, double normalFlux=0){
99                 _limitField[groupName]=LimitFieldStationaryDiffusion(NeumannStationaryDiffusion,-1, normalFlux);
100         };
101
102         void setDirichletValues(map< int, double> dirichletBoundaryValues);
103         void setNeumannValues  (map< int, double>   neumannBoundaryValues);
104         
105         void setConductivity(double conductivite){
106                 _conductivity=conductivite;
107         };
108         void setDiffusiontensor(Matrix DiffusionTensor){
109                 _DiffusionTensor=DiffusionTensor;
110         };
111
112
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
116         
117         void setFluidTemperatureField(Field coupledTemperatureField){
118                 _fluidTemperatureField=coupledTemperatureField;
119                 _fluidTemperatureFieldSet=true;
120         };
121         void setFluidTemperature(double fluidTemperature){
122         _fluidTemperature=fluidTemperature;
123         }
124         Field& getFluidTemperatureField(){
125                 return _fluidTemperatureField;
126         }
127         /** \fn setHeatPowerField
128          * \brief set the heat power field (variable in space)
129          * \details
130          * \param [in] Field
131          * \param [out] void
132          *  */
133         void setHeatPowerField(Field heatPower){
134                 _heatPowerField=heatPower;
135                 _heatPowerFieldSet=true;
136         }
137
138         /** \fn setHeatPowerField
139          * \brief set the heat power field (variable in space)
140          * \details
141          * \param [in] string fileName (including file path)
142          * \param [in] string fieldName
143          * \param [out] void
144          *  */
145         void setHeatPowerField(string fileName, string fieldName){
146                 _heatPowerField=Field(fileName, CELLS,fieldName);
147                 _heatPowerFieldSet=true;
148         }
149
150         /** \fn getHeatPowerField
151          * \brief returns the heat power field
152          * \details
153          * \param [in] void
154          * \param [out] Field
155          *  */
156         Field getHeatPowerField(){
157                 return _heatPowerField;
158         }
159         //get output fields names for postprocessing or coupling
160         vector<string> getOutputFieldsNames() ;//liste tous les champs que peut fournir le code pour le postraitement
161         Field&         getOutputField(const string& nameField );//Renvoie un champs pour le postraitement
162
163     Field& getOutputTemperatureField();
164         Field& getRodTemperatureField();
165
166         /** \fn setVerbose
167          * \brief Updates display options
168          * \details
169          * \param [in] bool
170          * \param [in] bool
171          * \param [out] void
172          *  */
173         void setVerbose(bool verbose,  bool system=false)
174         {
175                 _verbose = verbose;
176                 _system = system;
177         };
178
179 protected :
180         //Main unknown field
181         Field _VV;
182
183         int _Ndim;//space dimension
184         int _nVar;//Number of equations to solve=1
185
186     //Mesh data
187         Mesh _mesh;
188     bool _meshSet;
189         bool _initializedMemory;
190         int _Nmailles;//number of cells for FV calculation
191         int _neibMaxNbCells;//maximum number of cells around a cell
192     
193         double _precision;
194         double _precision_Newton;
195         double _erreur_rel;//norme(Uk+1-Uk)
196     bool _computationCompletedSuccessfully;
197     
198         //Linear solver and petsc
199         KSP _ksp;
200         KSPType _ksptype;
201         PC _pc;
202         PCType _pctype;
203         string _pc_hypre;
204         int _maxPetscIts;//nombre maximum d'iteration gmres autorisé au cours d'une resolution de système lineaire
205         int _maxNewtonIts;//nombre maximum d'iteration de Newton autorise au cours de la resolution d'un pas de temps
206         int _PetscIts;//the number of iterations of the linear solver
207         int _NEWTON_its;
208         Mat  _A;//Linear system matrix
209         Vec _b;//Linear system right hand side
210         double _MaxIterLinearSolver;//nombre maximum d'iteration gmres obtenu au cours par les resolution de systemes lineaires au cours d'un pas de tmeps
211         bool _conditionNumber;//computes an estimate of the condition number
212
213         map<string, LimitFieldStationaryDiffusion> _limitField;
214     bool _onlyNeumannBC;//if true then the linear system is singular and should be solved up to a constant vector
215     
216         bool _diffusionMatrixSet;
217         Vector _normale;
218         Matrix _DiffusionTensor;
219         Vec _deltaT, _Tk, _Tkm1, _b0;
220         Vec _Tk_seq; // Local sequential copy of the parallel vector _Tk, used for saving result files
221
222         double _dt_src;
223     
224         //Heat transfert variables
225         double _conductivity, _fluidTemperature;
226         Field _heatPowerField, _fluidTemperatureField;
227         bool _heatPowerFieldSet, _fluidTemperatureFieldSet;
228         double _heatTransfertCoeff, _heatSource;
229
230         //Display variables
231         bool _verbose, _system;
232         ofstream * _runLogFile;//for creation of a log file to save the history of the simulation
233     //saving parameters
234         string _fileName;//name of the calculation
235         string _path;//path to execution directory used for saving results
236         saveFormat _saveFormat;//file saving format : MED, VTK or CSV
237
238         double computeRHS(bool & stop);
239         double computeDiffusionMatrixFV(bool & stop);
240         double computeDiffusionMatrixFE(bool & stop);
241
242     /************ Data for FE calculation *************/
243     bool _FECalculation;
244         int _Nnodes;/* number of nodes for FE calculation */
245         int _neibMaxNbNodes;/* maximum number of nodes around a node */
246         int _NunknownNodes;/* number of unknown nodes for FE calculation */
247         int _NboundaryNodes;/* total number of boundary nodes */
248         int _NdirichletNodes;/* number of boundary nodes with Dirichlet BC for FE calculation */
249     std::vector< int > _boundaryNodeIds;/* List of boundary nodes */
250     std::vector< int > _dirichletNodeIds;/* List of boundary nodes with Dirichlet BC */
251
252     /********* Possibility to set a boundary field as DirichletNeumann boundary condition *********/
253     bool _dirichletValuesSet;
254     bool _neumannValuesSet;
255     std::map< int, double> _dirichletBoundaryValues;
256     std::map< int, double> _neumannBoundaryValues;
257
258         /**** MPI related variables ***/
259         PetscMPIInt    _mpi_size;        /* size of communicator */
260         PetscMPIInt    _mpi_rank;        /* processor rank */
261         VecScatter _scat;                       /* For the distribution of a local vector */
262         int _globalNbUnknowns, _localNbUnknowns;
263         int _d_nnz, _o_nnz;                     /* local and "non local" numbers of non zeros corfficients */
264 };
265
266 #endif /* StationaryDiffusionEquation_HXX_ */