]> SALOME platform Git repositories - tools/solverlab.git/blob - CoreFlows/examples/C/DiffusionEquation_1DHeatedRod.cxx
Salome HOME
Updated GUI documentation
[tools/solverlab.git] / CoreFlows / examples / C / DiffusionEquation_1DHeatedRod.cxx
1 #include "DiffusionEquation.hxx"\r
2 \r
3 using namespace std;\r
4 \r
5 #define PI 3.14159265\r
6 \r
7 void power_field_diffusionTest(Field & Phi){\r
8         double L=4.2;\r
9         double lambda=0.2;\r
10         double phi=1e5;\r
11         double x;\r
12         Mesh M = Phi.getMesh();\r
13         int nbCells = M.getNumberOfCells();\r
14         for (int j = 0; j < nbCells; j++) {\r
15                 x=M.getCell(j).x();\r
16                 Phi(j) = phi*cos(PI*(x-L/2)/(L+lambda));\r
17         }\r
18 }\r
19 \r
20 int main(int argc, char** argv)\r
21 {\r
22         //Preprocessing: mesh and group creation\r
23         double xinf=0.0;\r
24         double xsup=4.2;\r
25         int nx=10;\r
26         cout << "Building of a 1D mesh with "<<nx<<" cells" << endl;\r
27         Mesh M(xinf,xsup,nx);\r
28         double eps=1.E-6;\r
29         M.setGroupAtPlan(xsup,0,eps,"Neumann");\r
30         M.setGroupAtPlan(xinf,0,eps,"Neumann");\r
31         int spaceDim = M.getSpaceDimension();\r
32 \r
33 \r
34         //Solid parameters\r
35         double cp_ur=300;//Uranium specific heat\r
36         double rho_ur=10000;//Uranium density\r
37         double lambda_ur=5;\r
38  \r
39     bool FEcalculation=false;\r
40         DiffusionEquation  myProblem(spaceDim,FEcalculation,rho_ur,cp_ur,lambda_ur);\r
41 \r
42         //Setting initial field\r
43         Vector VV_Constant(1);\r
44         VV_Constant(0) = 623;//Rod clad temperature\r
45 \r
46         cout << "Building initial data" << endl;\r
47         myProblem.setInitialFieldConstant(M,VV_Constant);\r
48 \r
49         //Set fluid temperature (temperature du fluide)\r
50         double fluidTemp=573;//fluid mean temperature\r
51         double heatTransfertCoeff=1000;//fluid/solid exchange coefficient\r
52         myProblem.setFluidTemperature(fluidTemp);\r
53         myProblem.setHeatTransfertCoeff(heatTransfertCoeff);\r
54         //Set heat source\r
55         Field Phi("Heat power field", CELLS, M, 1);\r
56         power_field_diffusionTest(Phi);\r
57         myProblem.setHeatPowerField(Phi);\r
58         Phi.writeVTK("1DheatPowerField");\r
59 \r
60         //set the boundary conditions\r
61         myProblem.setNeumannBoundaryCondition("Neumann");\r
62 \r
63         // set the numerical method\r
64         myProblem.setTimeScheme( Explicit);\r
65 \r
66         // name result file\r
67         string fileName = "1DRodTemperature_FV";\r
68 \r
69         // parameters calculation\r
70         unsigned MaxNbOfTimeStep =3;\r
71         int freqSave = 1;\r
72         double cfl = 0.5;\r
73         double maxTime = 1000000;\r
74         double precision = 1e-6;\r
75 \r
76         myProblem.setCFL(cfl);\r
77         myProblem.setPrecision(precision);\r
78         myProblem.setMaxNbOfTimeStep(MaxNbOfTimeStep);\r
79         myProblem.setTimeMax(maxTime);\r
80         myProblem.setFreqSave(freqSave);\r
81         myProblem.setFileName(fileName);\r
82 \r
83         // set display option to monitor the calculation\r
84         myProblem.setVerbose( true);\r
85         //set file saving format\r
86         myProblem.setSaveFileFormat(CSV);\r
87 \r
88         // evolution\r
89         myProblem.initialize();\r
90         bool ok = myProblem.run();\r
91         if (ok)\r
92                 cout << "Simulation "<<fileName<<" is successful !" << endl;\r
93         else\r
94                 cout << "Simulation "<<fileName<<"  failed ! " << endl;\r
95 \r
96         cout << "------------ End of simulation -----------" << endl;\r
97         myProblem.terminate();\r
98 \r
99         return EXIT_SUCCESS;\r
100 }\r