createRainFields.H
Go to the documentation of this file.
1 PtrList<volVectorField> Urain;
2 Urain.setSize(phases.size());
3 PtrList<surfaceScalarField> phirain;
4 phirain.setSize(phases.size());
5 PtrList<volScalarField> alpharain;
6 alpharain.setSize(phases.size());
7 PtrList<volScalarField> scr;
8 scr.setSize(phases.size());
9 
10 std::stringstream converter;
11 for (int phase_no = 0; phase_no < phases.size(); phase_no++)
12 {
13  string phase_Uname = "U";
14  string phase_phiname = "phi";
15  string phase_alphaname = "alpha";
16 
17  string phase_no_str;
18  std::stringstream convert;
19  convert << phase_no+1;
20  phase_no_str = convert.str();
21 
22  phase_Uname.append(phase_no_str);
23  phase_phiname.append(phase_no_str);
24  phase_alphaname.append(phase_no_str);
25 
26  tmp<volVectorField> Utemp(new volVectorField
27  (
28  IOobject
29  (
30  phase_Uname,
31  runTime.timeName(),
32  mesh,
33  IOobject::MUST_READ,
34  IOobject::AUTO_WRITE
35  ),
36  mesh
37  ));
38  Urain.set(phase_no,Utemp);
39 
40  tmp<surfaceScalarField> phitemp(new surfaceScalarField
41  (
42  IOobject
43  (
44  phase_phiname,
45  runTime.timeName(),
46  mesh,
47  IOobject::NO_READ,
48  IOobject::NO_WRITE
49  ),
50  linearInterpolate(Urain[phase_no]) & mesh.Sf()
51  ));
52  phirain.set(phase_no,phitemp);
53 
54  tmp<volScalarField> alphatemp(new volScalarField
55  (
56  IOobject
57  (
58  phase_alphaname,
59  runTime.timeName(),
60  mesh,
61  IOobject::MUST_READ,
62  IOobject::AUTO_WRITE
63  ),
64  mesh
65  ));
66  alpharain.set(phase_no,alphatemp);
67 
68 }
69 
70 
phirain
PtrList< surfaceScalarField > phirain
Definition: createRainFields.H:3
scr
PtrList< volScalarField > scr
Definition: createRainFields.H:7
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
alpharain
PtrList< volScalarField > alpharain
Definition: createRainFields.H:5
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:55
converter
std::stringstream converter
Definition: createRainFields.H:10
Foam::linearInterpolate
tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > linearInterpolate(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: linear.H:107
Urain
PtrList< volVectorField > Urain
Definition: createRainFields.H:1
Foam::surfaceScalarField
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Definition: surfaceFieldsFwd.H:52
phases
multiphaseSystem::phaseModelList & phases
Definition: createFields.H:11