Loading...
Searching...
No Matches
createFields.H
Go to the documentation of this file.
1Info<< "Reading velocity field U\n" << endl;
2volVectorField U
3(
4 IOobject
5 (
6 "U",
7 runTime.timeName(),
8 mesh,
9 IOobject::MUST_READ,
10 IOobject::AUTO_WRITE
11 ),
12 mesh
13);
14
15// Initialise the velocity internal field to zero
16// Note: explicitly bypass evaluation of contraint patch overrides
17// (e.g. swirlFanVelocity might lookup phi,rho)
18//U = Zero;
19{
20 U.internalFieldRef() = Zero;
21 U.boundaryFieldRef() = Zero;
22}
23
25(
26 IOobject
27 (
28 "phi",
29 runTime.timeName(),
30 mesh,
31 IOobject::NO_READ,
32 IOobject::AUTO_WRITE
33 ),
34 fvc::flux(U)
35);
36
37if (args.found("initialiseUBCs"))
38{
39 U.correctBoundaryConditions();
40 phi = fvc::flux(U);
41}
42
43
44// Construct a pressure field
45// If it is available read it otherwise construct from the velocity BCs
46// converting fixed-value BCs to zero-gradient and vice versa.
47
48// Allow override from command-line -pName option
49const word pName = args.getOrDefault<word>("pName", "p");
50
51// Infer the pressure BCs from the velocity
52wordList pBCTypes
53(
54 U.boundaryField().size(),
55 fixedValueFvPatchScalarField::typeName
56);
57
58forAll(U.boundaryField(), patchi)
59{
60 if (U.boundaryField()[patchi].fixesValue())
61 {
62 pBCTypes[patchi] = zeroGradientFvPatchScalarField::typeName;
63 }
64}
65
66Info<< "Constructing pressure field " << pName << nl << endl;
68(
69 IOobject
70 (
71 pName,
72 runTime.timeName(),
73 mesh,
74 IOobject::READ_IF_PRESENT,
75 IOobject::NO_WRITE
76 ),
77 mesh,
78 dimensionedScalar(sqr(dimVelocity), Zero),
80);
81
82// Infer the velocity potential BCs from the pressure
83wordList PhiBCTypes
84(
85 p.boundaryField().size(),
86 zeroGradientFvPatchScalarField::typeName
87);
88
89forAll(p.boundaryField(), patchi)
90{
91 if (p.boundaryField()[patchi].fixesValue())
92 {
93 PhiBCTypes[patchi] = fixedValueFvPatchScalarField::typeName;
94 }
95}
96
97Info<< "Constructing velocity potential field Phi\n" << endl;
99(
100 IOobject
101 (
102 "Phi",
103 runTime.timeName(),
104 mesh,
105 IOobject::READ_IF_PRESENT,
106 IOobject::NO_WRITE
107 ),
108 mesh,
109 dimensionedScalar(dimLength*dimVelocity, Zero),
110 PhiBCTypes
111);
112
113label PhiRefCell = 0;
114scalar PhiRefValue = 0;
116(
117 Phi,
118 potentialFlow.dict(),
119 PhiRefCell,
120 PhiRefValue
121);
122mesh.setFluxRequired(Phi.name());
123
124#include "createMRF.H"
const dictionary & potentialFlow(mesh.solutionDict().subDict("potentialFlow"))
U
Definition pEqn.H:72
volScalarField & p
dynamicFvMesh & mesh
engineTime & runTime
List< word > wordList
List of word.
Definition fileName.H:60
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
Foam::argList args(argc, argv)
wordList pBCTypes(U.boundaryField().size(), fixedValueFvPatchScalarField::typeName)
setRefCell(p, pimple.dict(), pRefCell, pRefValue)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299