Salome HOME
Simplified PETSc includes
[tools/solverlab.git] / CoreFlows / Models / inc / ProblemCoreFlows.hxx
1 //============================================================================
2 // Name        : ProblemCoreFlows
3 // Author      : M. Ndjinga
4 // Version     :
5 // Copyright   : CEA Saclay 2014
6 // Description : Generic class for PDEs problems
7 //============================================================================
8 /* A ProblemCoreFlows class */
9
10 /*! \class ProblemCoreFlows ProblemCoreFlows.hxx "ProblemCoreFlows.hxx"
11  *  \brief Generic class for thermal hydraulics problems
12  *  \details  Common functions to CoreFlows models
13  */
14
15 #ifndef PROBLEMCOREFLOWS_HXX_
16 #define PROBLEMCOREFLOWS_HXX_
17
18 #include <iostream>
19 #include <fstream>
20 #include <unistd.h>
21 #include <vector>
22 #include <string>
23 #include <map>
24
25 #include <petscksp.h>
26 #include <slepceps.h>
27 #include <slepcsvd.h>
28
29 #include "Field.hxx"
30 #include "Mesh.hxx"
31 #include "Cell.hxx"
32 #include "Face.hxx"
33 #include "CdmathException.hxx"
34
35 using namespace std;
36
37 //! enumeration linearSolver
38 /*! the linearSolver can be GMRES or BiCGStab (see Petsc documentation) */
39 enum linearSolver
40 {
41         GMRES,/**< linearSolver is GMRES */
42         BCGS,/**< linearSolver is BiCGSstab */
43         CG/**< linearSolver is CG */
44 };
45
46 //! enumeration preconditioner
47 /*! the preconditioner can be ILU or LU  (see Petsc documentation) */
48 enum preconditioner
49 {
50         ILU,/**< preconditioner is ILU(0) */
51         LU,/**< preconditioner is actually a direct solver (LU factorisation)*/
52         NONE,/**< no preconditioner used */
53         ICC,/**< preconditioner is ICC(0) */
54         CHOLESKY/**< preconditioner is actually a direct solver for symmetric matrices (CHOLESKY factorisation)*/
55 };
56
57 //! enumeration saveFormat
58 /*! the numerical results are saved using MED, VTK or CSV format */
59 enum saveFormat
60 {
61         MED,/**< MED format is used  */
62         VTK,/**< VTK format is used */
63         CSV/**< CSV format is used */
64 };
65
66 //! enumeration TimeScheme
67 /*! The numerical method can be Explicit or Implicit  */
68 enum TimeScheme
69 {
70         Explicit,/**<  Explicit numerical scheme */
71         Implicit/**< Implicit numerical scheme */
72 };
73
74 class ProblemCoreFlows
75 {
76 public :
77         //! Constructeur par défaut
78         ProblemCoreFlows();
79         virtual ~ProblemCoreFlows();
80         
81         // -*-*-*- Gestion du calcul (interface ICoCo -*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
82
83         /** \fn initialize
84          * \brief Alloue la mémoire et vérifie que  le maillage et les conditions limites/initiales sont bien définis
85          * \Details c'est une fonction virtuelle pure, on la surcharge dans le problem Fluide .
86          * @param  void
87          *  */
88         virtual void initialize()=0;
89
90         /** \fn terminate
91          * \brief vide la mémoire et enregistre le résultat final
92          * \Details on la surcharge dans le problem fluid
93          *  @param void
94          *  */
95         virtual void terminate()=0;
96
97         /** \fn computeTimeStep
98          * \brief Propose un pas de temps pour le calcul
99          * \details Pour proposer un pas de temps la fonction nécessite de discrétiser les opérateurs
100          *  convection, diffusion et sources. Pour chacun on emploie la condition cfl.
101          *   En cas de problème durant ce calcul (exemple t=tmax), la fonction renvoie stop=true sinon stop = false
102          *   Cest une fonction virtuelle pure, on la surcharge dans Problème Fulide
103          *  @param Stop booléen correspond a True si ya un problème dans le code (exemple t=tmax) False sinon
104          *  \return  dt le pas de temps
105          *  */
106         virtual double computeTimeStep(bool & stop)=0;
107
108         /** \fn initTimeStep
109          * \brief Enregistre le nouveau pas de temps dt et l'intègre aux opérateurs
110          *  discrets (convection, diffusion, sources)
111          *  \details c'est une fonction virtuelle pure, on la surcharge dans le problemfluid
112          *  @param  dt est le nouvel pas de temps (double)
113          *  \return false si dt <0 et True sinon
114          *                         */
115         virtual bool initTimeStep(double dt)=0;
116
117         /** \fn solveTimeStep
118          * \brief calcule les valeurs inconnues au pas de temps +1 .
119          *  \details c'est une fonction virtuelle ,
120          *  @param  void
121          *  \return Renvoie false en cas de problème durant le calcul (valeurs non physiques..)
122          *  */
123         virtual bool solveTimeStep();//
124
125         /** \fn validateTimeStep
126          * \brief Valide le calcule au temps courant
127          * \details met à jour le temps présent t=t+dt, sauvegarde les champs inconnus
128          * et reinitialise dt à 0, teste la stationnarité .
129          * c'est une fonction virtuel , on la surchage dans chacun des modèles
130          * @param  void
131          * \return  Renvoie false en cas de problème durant le calcul
132          *  */
133         virtual void validateTimeStep()=0;//
134
135         /** \fn abortTimeStep
136          * \brief efface les inconnues calculées par solveTimeStep() et reinitialise dt à 0
137          * \details c'est une fonction virtuelle pure, elle est surchargée dans ProblemFluid, TransportEquation et DiffusionEquation
138          *  */
139         virtual void abortTimeStep()=0;
140
141         /** \fn run
142          * \brief Vérifie tous les paramètres et lance le code en supposant que la cfl ne changera pas
143          * \details  Renvoie "false" en cas de problème durant le calcul
144          * C'est une fonction virtuelle pure, on la surcharge dans problème Fluide
145          * @param void
146          * \return Renvoie "false" en cas de problème durant le calcul
147          *  */
148         virtual bool run();//
149
150         /** \fn iterateTimeStep
151          * \brief Calcul d'une sous-itération du pas de temps en cours, typiquement une itération de Newton pour un schéma implicite
152          * \details Deux paramètres booléen (converged et ok) sont retournés
153          * converged, Vaut true si on peut passer au pas de temps suivant (shéma explicite ou convergence du schéma de Newton dans le schéma implicite)
154          * ok  vaut true si le calcul n'a pas rencontré d'erreur ou de problème particulier dans la mise à jour des variables physiques
155          * \param [in] bool, passage par reférence.
156          * \param [out] bool
157          *  */
158         virtual bool iterateTimeStep(bool &converged) = 0; //??
159
160         /** \fn isStationary
161          * \brief vérifie la stationnairité du problème .
162          * \details Renvoie "true" si le problème atteint l'état stationnaire et "false" sinon
163          * \param [in] void
164          * \param [out] bool
165          *  */
166         virtual bool isStationary() const;
167
168         /** \fn presentTime
169          * \brief Calcule la valeur du temps courant .
170          * \details
171          * \param [in] void
172          * \param [out] double
173          *  */
174         virtual double presentTime() const;
175
176         /*
177         //Coupling interface
178         virtual vector<string> getInputFieldsNames()=0 ;//Renvoie les noms des champs dont le problème a besoin (données initiales)
179         virtual  Field& getInputFieldTemplate(const string& name)=0;//Renvoie le format de champs attendu (maillage, composantes etc)
180         virtual void setInputField(const string& name, const Field& afield)=0;//enregistre les valeurs d'une donnée initiale
181         virtual vector<string> getOutputFieldsNames()=0 ;//liste tous les champs que peut fournir le code pour le postraitement
182         virtual Field& getOutputField(const string& nameField )=0;//Renvoie un champs pour le postraitement
183          */
184
185         Field getUnknownField() const;
186         
187         //paramètres du calcul -*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
188
189         /** \fn setPresentTime
190          * \brief met à jour _time (le temps courant du calcul)
191          * \details
192          * \param [in] double
193          * \param [out] void
194          *  */
195         void setPresentTime (double time);
196
197         /** \fn setMaxNbOfTimeStep
198          * \brief met à jour _maxNbOfTimeStep ( le nombre maximum d'itération du calcul )
199          * \details
200          * \param [in] int
201          * \param [out] void
202          *  */
203         void setMaxNbOfTimeStep(int maxNbOfTimeStep);
204
205         /** \fn setTimeMax
206          * \brief met à jour _timeMax (Le temps maximum du calcul)
207          * \details
208          * \param [in] double
209          * \param [out] void
210          *  */
211         void setTimeMax(double timeMax);
212
213         /** \fn setCFL
214          * \brief met à jour la _CFL
215          * \details
216          * \param [in] double
217          * \param [out] void
218          *  */
219         void setCFL(double cfl);
220
221         /** \fn setPrecision
222          * \brief met à jour _precision (la précision du calcule)
223          * \details
224          * \param [in] double
225          * \param [out] void
226          *  */
227         void setPrecision(double precision);
228
229         /** \fn setInitialField
230          * \brief sets the initial field
231          * \details
232          * \param [in] Field
233          * \param [out] void
234          *  */
235         void setInitialField(const Field &VV);
236
237         /** \fn setInitialField
238          * \brief sets the initial field from a field in a med file
239          * \details
240          * \param [in] string : the file name
241          * \param [in] string : the field name
242          * \param [in] int : the time step number
243          * \param [out] void
244          *  */
245         void setInitialField(string fileName, string fieldName, int timeStepNumber);
246
247         /** \fn setInitialFieldConstant
248          * \brief sets a constant initial field on a mesh stored in a med file
249          * \details
250          * \param [in] string : the file name
251          * \param [in] vector<double> : the value in each cell
252          * \param [out] void
253          *  */
254         void setInitialFieldConstant(string fileName, const vector<double> Vconstant);
255
256         /** \fn setInitialFieldConstant
257          * \brief sets a constant initial field 
258          * \details
259          * \param [in] Mesh 
260          * \param [in] Vector
261          * \param [out] void
262          *  */
263         void setInitialFieldConstant(const Mesh& M, const Vector Vconstant);
264
265         /** \fn setInitialFieldConstant
266          * \brief sets a constant initial field
267          * \details
268          * \param [in] Mesh
269          * \param [in] vector<double>
270          * \param [out] void
271          *  */
272         void setInitialFieldConstant(const Mesh& M, const vector<double> Vconstant);
273
274         /** \fn setInitialFieldConstant
275          * \brief sets a constant initial field
276          * \details
277          * \param [in] int the space dimension
278          * \param [in] vector<double> the value in each cell
279          * \param [in] double the lowest value in the x direction
280          * \param [in] double the highest value in the x direction
281          * \param [in] string name of the left boundary
282          * \param [in] string name of the right boundary
283          * \param [in] double the lowest value in the y direction
284          * \param [in] double the highest value in the y direction
285          * \param [in] string name of the back boundary
286          * \param [in] string name of the front boundary
287          * \param [in] double the lowest value in the z direction
288          * \param [in] double the highest value in the z direction
289          * \param [in] string name of the bottom boundary
290          * \param [in] string name of the top boundary
291          * \param [out] void
292          *  */
293         void setInitialFieldConstant( int nDim, const vector<double> Vconstant, double xmin, double xmax,int nx, string leftSide, string rightSide,
294                         double ymin=0, double ymax=0, int ny=0, string backSide="", string frontSide="",
295                         double zmin=0, double zmax=0, int nz=0, string bottomSide="", string topSide="");
296
297         /** \fn setInitialFieldStepFunction
298          * \brief sets a step function initial field (Riemann problem)
299          * \details
300          * \param [in] Mesh
301          * \param [in] Vector
302          * \param [in] Vector
303          * \param [in] double position of the discontinuity on one of the three axis
304          * \param [in] int direction (axis carrying the discontinuity) : 0 for x, 1 for y, 2 for z
305          * \param [out] void
306          *  */
307         void setInitialFieldStepFunction(const Mesh M, const Vector Vleft, const Vector Vright, double disc_pos, int direction=0);
308
309         /** \fn setInitialFieldStepFunction
310          * \brief sets a constant initial field
311          * \details
312          * \param [in] int the space dimension
313          * \param [in] vector<double> the value left of the discontinuity
314          * \param [in] vector<double> the value right of the discontinuity
315          * \param [in] double the position of the discontinuity in the x direction
316          * \param [in] double the lowest value in the x direction
317          * \param [in] double the highest value in the x direction
318          * \param [in] string name of the left boundary
319          * \param [in] string name of the right boundary
320          * \param [in] double the lowest value in the y direction
321          * \param [in] double the highest value in the y direction
322          * \param [in] string name of the back boundary
323          * \param [in] string name of the front boundary
324          * \param [in] double the lowest value in the z direction
325          * \param [in] double the highest value in the z direction
326          * \param [in] string name of the bottom boundary
327          * \param [in] string name of the top boundary
328          * \param [out] void
329          *  */
330         void setInitialFieldStepFunction( int nDim, const vector<double> VV_Left, vector<double> VV_Right, double xstep,
331                         double xmin, double xmax,int nx, string leftSide, string rightSide,
332                         double ymin=0, double ymax=0, int ny=0, string backSide="", string frontSide="",
333                         double zmin=0, double zmax=0, int nz=0, string bottomSide="", string topSide="");
334
335         /** \fn setInitialFieldSphericalStepFunction
336          * \brief sets a step function initial field with value Vin inside the ball with radius Radius and Vout outside
337          * \details
338          * \param [in] Mesh
339          * \param [in] Vector Vin, value inside the ball
340          * \param [in] Vector Vout, value outside the ball
341          * \param [in] double radius of the ball
342          * \param [in] Vector Center, coordinates of the ball center
343          * \param [out] void
344          *  */
345         void setInitialFieldSphericalStepFunction(const Mesh M, const Vector Vin, const Vector Vout, double Radius, Vector Center);
346
347         /** \fn getTime
348          * \brief renvoie _time (le temps courant du calcul)
349          * \details
350          * \param [in] void
351          * \param [out] double
352          *  */
353         double getTime();
354
355         /** \fn getNbTimeStep
356          * \brief renvoie _nbTimeStep le Numéro d'itération courant
357          * \details
358          * \param [in] void
359          * \param [out] unsigned
360          *  */
361         unsigned getNbTimeStep();
362
363         /** \fn getCFL
364          * \brief renvoie la _CFL
365          * \details
366          * \param [in] void
367          * \param [out] double
368          *  */
369         double getCFL();
370
371         /** \fn getPrecision
372          * \brief renvoie _precision (la précision du calcul)
373          * \details
374          * \param [in] void
375          * \param [out] double
376          *  */
377         double getPrecision();
378
379         /** \fn getMesh
380          * \brief renvoie _Mesh (le maillage du problème)
381          * \details
382          * \param [in] void
383          * \param [out] Mesh
384          *  */
385         Mesh getMesh();
386
387         /** \fn setFileName
388          * \brief met à jour _fileName le nom du fichier
389          * \details
390          * \param [in]  string
391          * \param [out] void
392          *  */
393         void setFileName(string fileName);
394
395         /** \fn setFreqSave
396          * \brief met à jour _FreqSave (la fréquence du sauvgarde de la solution)
397          * \details
398          * \param [in] double
399          * \param [out] void
400          *  */
401         void setFreqSave(int freqSave);
402
403         /** \fn save
404          * \brief sauvgarde les données dans des fichiers MED ou VTK
405          * \details c'est une fonction virtuelle pure , on la surcharge
406          * dans chacun des modèles
407          * @param  void
408          */
409         virtual void save() = 0;
410
411         /** \fn getLinearSolver
412          * \brief renvoie _ksptype (le type du solveur linéaire utilisé)
413          * \details
414          * \param [in] void
415          * \param [out] string
416          *  */
417         string getLinearSolver() {
418                 return _ksptype;
419         };
420
421         /** \fn getNumberOfVariables
422          * \brief le nombre d'inconnues du problème
423          * \details
424          * @param void
425          * \return renvoie _nVar (le nombre d'inconnues du problème)
426          *  */
427         int getNumberOfVariables(){
428                 return _nVar;
429         };
430
431         /** \fn setWellBalancedCorrection
432          * \brief include a well balanced correction to treat stiff source terms
433          * @param boolean that is true if a well balanced correction should be applied
434          * */
435         void setWellBalancedCorrection(bool wellBalancedCorr){
436                 _wellBalancedCorrection=wellBalancedCorr;
437         }
438
439         /** \fn setLinearSolver
440          * \brief sets the linear solver and preconditioner
441          * \details virtual function overloaded by intanciable classes
442          * @param kspType linear solver type (GMRES or BICGSTAB)
443          * @param pcType preconditioner (ILU,LU or NONE)
444          */
445         void setLinearSolver(linearSolver solverName, preconditioner pcType);
446
447         /** \fn setNewtonSolver
448          * \brief set the Newton algorithm parameters
449          * \param [in] int maximum number of newton iterations
450          * \param [in] double precision required for the convergence of the newton scheme
451          * \param [out] void
452          *  */
453         void setNewtonSolver(double precision,int iterations=20)
454         {
455                 _maxNewtonIts=iterations;
456                 _precision_Newton=precision;
457         };
458
459         /** \fn displayConditionNumber
460          * \brief display the condition number of the preconditioned linear systems
461          */
462         void displayConditionNumber(bool display=true){
463                 _conditionNumber=display;
464         }
465
466         /** \fn setSaveFileFormat
467          * \brief sets the numerical results file format (MED, VTK or CSV)
468          * \details
469          * \param [in] saveFormat
470          * \param [out] void
471          *  */
472         void setSaveFileFormat(saveFormat saveFileFormat){
473                 _saveFormat=saveFileFormat;
474         }
475
476         //Couplages Thermohydraulique-thermique-neutronique *-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*-*
477
478         /** \fn setHeatPowerField
479          * \brief set the heat power field (variable in space)
480          * \details
481          * \param [in] Field
482          * \param [out] void
483          *  */
484         void setHeatPowerField(Field heatPower){
485                 _heatPowerField=heatPower;
486                 _heatPowerFieldSet=true;
487         }
488
489         /** \fn setHeatPowerField
490          * \brief set the heat power field (variable in space)
491          * \details
492          * \param [in] string fileName (including file path)
493          * \param [in] string fieldName
494          * \param [out] void
495          *  */
496         void setHeatPowerField(string fileName, string fieldName){
497                 _heatPowerField=Field(fileName, CELLS,fieldName);
498                 _heatPowerFieldSet=true;
499         }
500
501         /** \fn setHeatSource
502          * \brief met à jour la puissance thermique ( _phi )
503          * \details
504          * \param [in] double
505          * \param [out] void
506          *  */
507         void setHeatSource(double phi){
508                 _heatSource=phi;
509         }
510
511         /** \fn getHeatPowerField
512          * \brief renvoie le champs ?? ( _heatPowerField )
513          * \details
514          * \param [in] void
515          * \param [out] Field
516          *  */
517         Field getHeatPowerField(){
518                 return _heatPowerField;
519         }
520
521         /** \fn setRodTemperatureField ??
522          * \brief
523          * \details
524          * \param [in] Field
525          * \param [out] void
526          *  */
527         void setRodTemperatureField(Field rodTemperature){
528                 _rodTemperatureField=rodTemperature;
529                 _rodTemperatureFieldSet=true;
530         }
531
532         /** \fn setRodTemperature ??
533          * \brief
534          * \details
535          * \param [in] double
536          * \param [out] void
537          *  */
538         void setRodTemperature(double rodTemp){
539                 _rodTemperature=rodTemp;
540         }
541
542         /** \fn getRodTemperatureField
543          * \brief
544          * \details
545          * \param [in] void
546          * \param [out] Field
547          *  */
548         virtual Field& getRodTemperatureField(){ // ?? je ne retrouve pas cet attribut dans le file.cxx
549                 return _rodTemperatureField;
550         }
551
552         /** \fn setHeatTransfertCoeff
553          * \brief set the heat transfert coefficient for heat exchange between fluid and solid
554          * \details
555          * \param [in] double
556          * \param [out] void
557          *  */
558         void setHeatTransfertCoeff(double heatTransfertCoeff){
559                 _heatTransfertCoeff=heatTransfertCoeff;
560         }
561
562         /** \fn setDISPLAY
563          * \brief met à jour les paramètres de l'affichage
564          * \details
565          * \param [in] bool
566          * \param [in] bool
567          * \param [out] void
568          *  */
569         void setVerbose(bool verbose,  bool system=false)
570         {
571                 _verbose = verbose;
572                 _system = system;
573         };
574
575     //Spectrum analysis
576     double getConditionNumber(bool isSingular=false, double tol=1e-6) const;
577     std::vector< double > getEigenvalues (int nev, EPSWhich which=EPS_SMALLEST_MAGNITUDE, double tol=1e-6) const;
578     std::vector< Vector > getEigenvectors(int nev, EPSWhich which=EPS_SMALLEST_MAGNITUDE, double tol=1e-6) const;
579     Field getEigenvectorsField(int nev, EPSWhich which=EPS_SMALLEST_MAGNITUDE, double tol=1e-6) const;
580
581         //  some supplementary functions
582
583         /** \fn displayMatrix
584          * \brief displays a matrix of size "size x size" for profiling
585          * @param  matrix is a pointer of size "size"
586          * @param size, size of the matrix
587          * @param name, string, name or description of the matrix
588          * @return displays the matrix on the terminal
589          *  */
590         void displayMatrix(double *matrix, int size, string name);
591
592         /** \fn displayMatrix
593          * \brief displays a vector of size "size" for profiling
594          * @param  vector is a pointer of size "size"
595          * @param size, size of the vector
596          * @param name, string, name or description of the vector
597          * @return displays the vector on the terminal
598          *  */
599         void displayVector(double *vector, int size, string name);
600
601         /** \fn getTimeScheme
602          * \brief returns the  time scheme name
603          * \param [in] void
604          * \param [out] enum TimeScheme (explicit or implicit)
605          *  */
606         TimeScheme getTimeScheme();
607
608         /** \fn setNumericalScheme
609          * \brief sets the numerical method ( explicit vs implicit )
610          * \details
611          * \param [in] TimeScheme
612          * \param [out] void
613          *  */
614         void setTimeScheme( TimeScheme method);
615
616
617 protected :
618
619         int _Ndim;//space dimension
620         int _nVar;//Number of equations to sole
621         int _Nmailles;//number of cells
622         int _Nnodes;//number of nodes
623         int _Nfaces;//number of faces
624         int _neibMaxNb;//maximum number of neighbours around a cell
625         int _neibMaxNbNodes;/* maximum number of nodes around a node */
626         Mesh _mesh;
627         Field _perimeters;
628
629         //Main unknown field
630         Field _VV;
631
632         //Numerical method
633         double _dt;
634         double _cfl;
635         double _maxvp;//valeur propre max pour calcul cfl
636         double _minl;//minimum cell diameter
637     bool _FECalculation;
638         /** boolean used to specify that a well balanced correction should be used */
639         bool _wellBalancedCorrection;
640         TimeScheme _timeScheme;
641
642         //Linear solver and petsc
643         KSP _ksp;
644         KSPType _ksptype;
645         PC _pc;
646         PCType _pctype;
647         string _pc_hypre;
648         int _maxPetscIts;//nombre maximum d'iteration gmres autorise au cours d'une resolution de systeme lineaire
649         int _PetscIts;//the number of iterations of the linear solver
650         int _maxNewtonIts;//nombre maximum d'iteration de Newton autorise au cours de la resolution d'un pas de temps
651         int _NEWTON_its;
652         Mat  _A;//Linear system matrix
653         Vec _b;//Linear system right hand side
654         double _MaxIterLinearSolver;//nombre maximum d'iteration gmres obtenu au cours par les resolution de systemes lineaires au cours d'un pas de tmeps
655         bool _conditionNumber;//computes an estimate of the condition number
656
657         //simulation monitoring variables
658         bool _isStationary;
659         bool _initialDataSet;
660         bool _initializedMemory;
661         bool _restartWithNewTimeScheme;
662         bool _restartWithNewFileName;
663         double _timeMax,_time;
664         int _maxNbOfTimeStep,_nbTimeStep;
665         double _precision;
666         double _precision_Newton;
667         double _erreur_rel;//norme(Un+1-Un)
668         string _fileName;//name of the calculation
669         int _freqSave;
670         ofstream * _runLogFile;//for creation of a log file to save the history of the simulation
671
672         //Heat transfert variables
673         Field _heatPowerField, _rodTemperatureField;
674         bool _heatPowerFieldSet, _rodTemperatureFieldSet;
675         double _heatTransfertCoeff;
676         double _heatSource, _rodTemperature;
677         double _hsatv, _hsatl;//all models appart from DiffusionEquation will need this
678
679         //Display variables
680         bool _verbose, _system;
681
682         string _path;//path to execution directory used for saving results
683         saveFormat _saveFormat;//file saving format : MED, VTK or CSV
684 };
685
686 #endif /* PROBLEMCOREFLOWS_HXX_ */