]> SALOME platform Git repositories - tools/solverlab.git/blob - CoreFlows/examples/C/SinglePhase_2DHeatDrivenCavity.cxx
Salome HOME
Updated GUI documentation
[tools/solverlab.git] / CoreFlows / examples / C / SinglePhase_2DHeatDrivenCavity.cxx
1 #include "SinglePhase.hxx"
2
3 using namespace std;
4
5 int main(int argc, char** argv)
6 {
7         /*Preprocessing: mesh and group creation*/
8         double xinf=0;
9         double xsup=1;
10         double yinf=0;
11         double ysup=1;
12         int nx=10;
13         int ny=10;
14         cout << "Building a regular mesh with "<<nx<<" times "<< ny<< " cells " << endl;
15         Mesh M(xinf,xsup,nx,yinf,ysup,ny);
16         double eps=1.E-6;
17         M.setGroupAtPlan(xinf,0,eps,"coldWall");
18         M.setGroupAtPlan(xsup,0,eps,"hotWall");
19         M.setGroupAtPlan(yinf,1,eps,"coldWall");
20         M.setGroupAtPlan(ysup,1,eps,"hotWall");
21         int spaceDim = M.getSpaceDimension();
22
23         // physical constants
24         vector<double> viscosite(1), conductivite(1);
25         viscosite[0]= 8.85e-5;
26         conductivite[0]=1000;//transfert de chaleur du à l'ébullition en paroi.
27         vector<double> gravite(spaceDim,0.) ;
28         gravite[1]=-10;
29         gravite[0]=0;
30
31         // set the limit field for each boundary
32         LimitField limitColdWall,  limitHotWall;
33         map<string, LimitField> boundaryFields;
34         limitColdWall.bcType=Wall;
35         limitColdWall.T = 590;//Temperature de la parois froide
36         limitColdWall.v_x = vector<double>(1,0);
37         limitColdWall.v_y = vector<double>(1,0);
38         boundaryFields["coldWall"]= limitColdWall;
39
40         limitHotWall.bcType=Wall;
41         limitHotWall.T = 560;//Temperature des parois chauffantes
42         limitHotWall.v_x = vector<double>(1,0);
43         limitHotWall.v_y = vector<double>(1,0);
44         boundaryFields["hotWall"]= limitHotWall;
45
46         SinglePhase  myProblem(Liquid,around155bars600K,spaceDim);
47         int nVar = myProblem.getNumberOfVariables();
48
49         //Initial field creation
50         cout << "Building initial data" << endl;
51         /* First case constant initial data */
52         Vector VV_Constant(nVar);
53         // constant vector
54         VV_Constant(0) = 155e5;
55         VV_Constant(1) = 0;
56         VV_Constant(2) = 0;
57         VV_Constant(3) = 573;
58
59         myProblem.setInitialFieldConstant(M,VV_Constant);
60
61         //set the boundary conditions
62         myProblem.setBoundaryFields(boundaryFields);
63
64         // physical parameters
65         myProblem.setViscosity(viscosite);
66         myProblem.setConductivity(conductivite);
67         myProblem.setGravity(gravite);
68
69         // set the numerical method
70         myProblem.setNumericalScheme(upwind, Implicit);
71
72         // set the Petsc resolution
73         myProblem.setLinearSolver(GMRES,LU,false);
74
75         // name result file
76         string fileName = "2DHeatDrivenCavity";
77
78         // parameters calculation
79         unsigned MaxNbOfTimeStep = 3;
80         int freqSave = 1;
81         double cfl = 10;
82         double maxTime = 50;
83         double precision = 1e-6;
84
85         myProblem.setCFL(cfl);
86         myProblem.setPrecision(precision);
87         myProblem.setMaxNbOfTimeStep(MaxNbOfTimeStep);
88         myProblem.setTimeMax(maxTime);
89         myProblem.setFreqSave(freqSave);
90         myProblem.setFileName(fileName);
91         myProblem.setNewtonSolver(precision,50);
92         myProblem.saveVelocity();
93
94         // evolution
95         myProblem.initialize();
96
97         bool ok = myProblem.run();
98         if (ok)
99                 cout << "Simulation "<<fileName<<" is successful !" << endl;
100         else
101                 cout << "Simulation "<<fileName<<"  failed ! " << endl;
102
103         cout << "------------ End of calculation !!! -----------" << endl;
104
105         myProblem.terminate();
106
107         return EXIT_SUCCESS;
108 }