createFields.H
Go to the documentation of this file.
1 Info<< "Reading velocity field U\n" << endl;
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 U = dimensionedVector(U.dimensions(), Zero);
17 
19 (
20  IOobject
21  (
22  "phi",
23  runTime.timeName(),
24  mesh,
25  IOobject::NO_READ,
26  IOobject::AUTO_WRITE
27  ),
28  fvc::flux(U)
29 );
30 
31 if (args.found("initialiseUBCs"))
32 {
33  U.correctBoundaryConditions();
34  phi = fvc::flux(U);
35 }
36 
37 
38 // Construct a pressure field
39 // If it is available read it otherwise construct from the velocity BCs
40 // converting fixed-value BCs to zero-gradient and vice versa.
41 
42 // Allow override from command-line -pName option
43 const word pName = args.getOrDefault<word>("pName", "p");
44 
45 // Infer the pressure BCs from the velocity
47 (
48  U.boundaryField().size(),
49  fixedValueFvPatchScalarField::typeName
50 );
51 
52 forAll(U.boundaryField(), patchi)
53 {
54  if (U.boundaryField()[patchi].fixesValue())
55  {
56  pBCTypes[patchi] = zeroGradientFvPatchScalarField::typeName;
57  }
58 }
59 
60 // Note that registerObject is false for the pressure field. The pressure
61 // field in this solver doesn't have a physical value during the solution.
62 // It shouldn't be looked up and used by sub models or boundary conditions.
63 Info<< "Constructing pressure field " << pName << nl << endl;
65 (
66  IOobject
67  (
68  pName,
69  runTime.timeName(),
70  mesh,
71  IOobject::READ_IF_PRESENT,
72  IOobject::NO_WRITE,
73  false
74  ),
75  mesh,
77  pBCTypes
78 );
79 
80 // Infer the velocity potential BCs from the pressure
81 wordList PhiBCTypes
82 (
83  p.boundaryField().size(),
84  zeroGradientFvPatchScalarField::typeName
85 );
86 
87 forAll(p.boundaryField(), patchi)
88 {
89  if (p.boundaryField()[patchi].fixesValue())
90  {
91  PhiBCTypes[patchi] = fixedValueFvPatchScalarField::typeName;
92  }
93 }
94 
95 Info<< "Constructing velocity potential field Phi\n" << endl;
97 (
98  IOobject
99  (
100  "Phi",
101  runTime.timeName(),
102  mesh,
103  IOobject::READ_IF_PRESENT,
104  IOobject::NO_WRITE
105  ),
106  mesh,
108  PhiBCTypes
109 );
110 
111 label PhiRefCell = 0;
112 scalar PhiRefValue = 0;
114 (
115  Phi,
116  potentialFlow.dict(),
117  PhiRefCell,
118  PhiRefValue
119 );
120 mesh.setFluxRequired(Phi.name());
121 
122 #include "createMRF.H"
123 
124 // Add solver-specific interpolations
125 {
126  wordHashSet& nonInt =
127  const_cast<wordHashSet&>(Stencil::New(mesh).nonInterpolatedFields());
128 
129  nonInt.insert("cellMask");
130  nonInt.insert("interpolatedCells");
131 
132 }
133 
134 // Mask field for zeroing out contributions on hole cells
135 #include "createCellMask.H"
136 
137 // Create bool field with interpolated cells
138 #include "createInterpolatedCells.H"
forAll(U.boundaryField(), patchi)
Definition: createFields.H:52
tmp< surfaceScalarField > flux(const volVectorField &vvf)
Return the face-flux field obtained from the given volVectorField.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
const dictionary & potentialFlow(mesh.solutionDict().subDict("potentialFlow"))
setRefCell(p, pimple.dict(), pRefCell, pRefValue)
T getOrDefault(const word &optName, const T &deflt) const
Get a value from the named option if present, or return default.
Definition: argListI.H:300
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:85
bool insert(const Key &key)
Insert a new entry, not overwriting existing entries.
Definition: HashSet.H:227
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:84
wordList pBCTypes(U.boundaryField().size(), fixedValueFvPatchScalarField::typeName)
dynamicFvMesh & mesh
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition: HashSet.H:73
Creates mask for interpolated cells.
volScalarField & p
Definition: createFields.H:23
volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedVector(dimVelocity, Zero))
List< word > wordList
A List of words.
Definition: fileName.H:58
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Creates mask for blocked out cells.
messageStream Info
Information stream (stdout output on master, null elsewhere)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Foam::argList args(argc, argv)
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:157
const dimensionSet dimVelocity