createFaFields.H
Go to the documentation of this file.
1  Info<< "Reading field h" << endl;
3  (
4  IOobject
5  (
6  "h",
7  runTime.timeName(),
8  mesh,
9  IOobject::MUST_READ,
10  IOobject::AUTO_WRITE
11  ),
12  aMesh
13  );
14 
15 
16  Info<< "Reading field Us" << endl;
18  (
19  IOobject
20  (
21  "Us",
22  runTime.timeName(),
23  mesh,
24  IOobject::MUST_READ,
25  IOobject::AUTO_WRITE
26  ),
27  aMesh
28  );
29 
30 
32  (
33  IOobject
34  (
35  "phis",
36  runTime.timeName(),
37  mesh,
38  IOobject::READ_IF_PRESENT,
39  IOobject::AUTO_WRITE
40  ),
41  fac::interpolate(Us) & aMesh.Le()
42  );
43 
44 
45  edgeScalarField phi2s
46  (
47  IOobject
48  (
49  "phi2s",
50  runTime.timeName(),
51  mesh,
52  IOobject::READ_IF_PRESENT,
53  IOobject::AUTO_WRITE
54  ),
55  fac::interpolate(h*Us) & aMesh.Le()
56  );
57 
58 
59  const areaVectorField& Ns = aMesh.faceAreaNormals();
60  areaVectorField Gs(g - Ns*(Ns & g));
61  areaScalarField Gn(mag(g - Gs));
62 
63  // Mass source
65  (
66  IOobject
67  (
68  "Sm",
69  runTime.timeName(),
70  mesh,
71  IOobject::NO_READ,
72  IOobject::NO_WRITE
73  ),
74  aMesh,
76  );
77 
78  // Mass sink
80  (
81  IOobject
82  (
83  "Sd",
84  runTime.timeName(),
85  mesh,
86  IOobject::NO_READ,
87  IOobject::NO_WRITE
88  ),
89  aMesh,
91  );
92 
94  (
95  IOobject
96  (
97  "Sg",
98  runTime.timeName(),
99  mesh,
100  IOobject::NO_READ,
101  IOobject::NO_WRITE
102  ),
103  aMesh,
105  );
106 
107 
108  // Surface pressure
109  areaScalarField ps
110  (
111  IOobject
112  (
113  "ps",
114  runTime.timeName(),
115  mesh,
116  IOobject::NO_READ,
117  IOobject::AUTO_WRITE
118  ),
119  rhol*Gn*h - sigma*fac::laplacian(h)
120  );
121 
122  // Friction factor
123  areaScalarField dotProduct
124  (
125  aMesh.faceAreaNormals() & (g/mag(g))
126  );
127 
128  Info<< "View factor: min = " << min(dotProduct.internalField())
129  << " max = " << max(dotProduct.internalField()) << endl;
130 
131  areaScalarField manningField
132  (
133  IOobject
134  (
135  "manningField",
136  runTime.timeName(),
137  mesh,
138  IOobject::MUST_READ,
139  IOobject::AUTO_WRITE
140  ),
141  aMesh
142  );
143 
144  areaScalarField frictionFactor
145  (
146  IOobject
147  (
148  "frictionFactor",
149  runTime.timeName(),
150  mesh,
151  IOobject::NO_READ,
152  IOobject::NO_WRITE
153  ),
154  aMesh,
155  dimensionedScalar("one", dimless, 0.01)
156  );
157 
158  aMesh.setFluxRequired("h");
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimVelocity
const dimensionSet dimVelocity
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
rhol
dimensionedScalar rhol("rhol", dimDensity, transportProperties)
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
aMesh
faMesh aMesh(mesh)
Foam::areaVectorField
GeometricField< vector, faPatchField, areaMesh > areaVectorField
Definition: areaFieldsFwd.H:58
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::constant::universal::h
const dimensionedScalar h
Planck constant.
Definition: setRegionSolidFields.H:33
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::edgeScalarField
GeometricField< scalar, faePatchField, edgeMesh > edgeScalarField
Definition: edgeFieldsFwd.H:52
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
g
const uniformDimensionedVectorField & g
Definition: createFluidFields.H:26
Foam::areaScalarField
GeometricField< scalar, faPatchField, areaMesh > areaScalarField
Definition: areaFieldsFwd.H:53
Us
Us
Definition: createFaFields.H:51
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
sigma
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47
Foam::fac::interpolate
static tmp< GeometricField< Type, faePatchField, edgeMesh > > interpolate(const GeometricField< Type, faPatchField, areaMesh > &tvf, const edgeScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
phis
edgeScalarField phis(IOobject("phis", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), linearEdgeInterpolate(Us) &aMesh.Le())
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189