correctPhi.H
Go to the documentation of this file.
1 {
2  if (mesh.changing())
3  {
4  volVectorField::Boundary& bfld = U.boundaryFieldRef();
5  forAll(bfld, patchi)
6  {
7  if (bfld[patchi].fixesValue())
8  {
9  bfld[patchi].initEvaluate();
10  }
11  }
12 
13  surfaceScalarField::Boundary& phiBfld = phi.boundaryFieldRef();
14  forAll(bfld, patchi)
15  {
16  if (bfld[patchi].fixesValue())
17  {
18  bfld[patchi].evaluate();
19 
20  phiBfld[patchi] =
21  bfld[patchi]
22  & mesh.Sf().boundaryField()[patchi];
23  }
24  }
25  }
26 
28  (
29  p_rgh.boundaryField().size(),
30  zeroGradientFvPatchScalarField::typeName
31  );
32 
33  for (label i=0; i<p_rgh.boundaryField().size(); i++)
34  {
35  if (p_rgh.boundaryField()[i].fixesValue())
36  {
37  pcorrTypes[i] = fixedValueFvPatchScalarField::typeName;
38  }
39  }
40 
42  (
43  IOobject
44  (
45  "pcorr",
46  runTime.timeName(),
47  mesh,
48  IOobject::NO_READ,
49  IOobject::NO_WRITE
50  ),
51  mesh,
52  dimensionedScalar(p_rgh.dimensions(), Zero),
54  );
55 
56  if (pcorr.needReference())
57  {
59  adjustPhi(phi, U, pcorr);
61  }
62 
63  mesh.setFluxRequired(pcorr.name());
64 
65  dimensionedScalar rAUf("rAUf", dimTime/rho.dimensions(), 1.0);
66 
67  const cellCellStencilObject& overlap = Stencil::New(mesh);
68  const labelList& cellTypes = overlap.cellTypes();
69  const labelIOList& zoneIDs = overlap.zoneID();
70 
71  while (pimple.correctNonOrthogonal())
72  {
73  label nZones = gMax(zoneIDs)+1;
74  //label refCellI2 = -1;
75 
76  labelList refCells(nZones, -1);
77  labelList refZones(nZones, -1);
78 
79  forAll(zoneIDs, cellI)
80  {
81  label zoneId = zoneIDs[cellI];
82  if
83  (
84  refCells[zoneId] == -1
85  && cellTypes[cellI] == cellCellStencil::CALCULATED
86  && refZones[zoneId] == -1
87  )
88  {
89  refCells[zoneId] = cellI;
90  refZones[zoneId] = zoneId;
91  }
92  }
93 
94  fvScalarMatrix pcorrEqn
95  (
97  );
98 
99  // Only set reference for cells that are CALCULATED
100  {
101  DynamicList<label> validCells(refCells.size());
102  forAll(refCells, zoneId)
103  {
104  if (refCells[zoneId] != -1)
105  {
106  validCells.append(refCells[zoneId]);
107  }
108  }
109 
110  pcorrEqn.setReferences
111  (
112  validCells,
113  scalar(0),
114  true
115  );
116  }
117 
118 
119  const dictionary& d = mesh.solver
120  (
121  pcorr.select
122  (
123  pimple.finalInnerIter()
124  )
125  );
126 
127  //Bypass virtual layer
128  mesh.fvMesh::solve(pcorrEqn, d);
129 
130  if (pimple.finalNonOrthogonalIter())
131  {
132  phi -= pcorrEqn.flux();
133  }
134  }
135 }
nZones
label nZones
Definition: interpolatedFaces.H:24
rAUf
dimensionedScalar rAUf("rAUf", dimTime, 1.0)
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::adjustPhi
bool adjustPhi(surfaceScalarField &phi, const volVectorField &U, volScalarField &p)
Adjust the balance of fluxes to obey continuity.
Definition: adjustPhi.C:37
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::labelIOList
IOList< label > labelIOList
Label container classes.
Definition: labelIOList.H:44
Foam::fac::div
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:50
rho
rho
Definition: readInitialConditions.H:88
zoneIDs
const labelIOList & zoneIDs
Definition: correctPhi.H:59
pcorr
volScalarField pcorr(IOobject("pcorr", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedScalar(p.dimensions(), Zero), pcorrTypes)
Foam::fvc::makeRelative
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:77
Foam::wordList
List< word > wordList
A List of words.
Definition: fileName.H:62
pimple
pimpleControl & pimple
Definition: setRegionFluidFields.H:56
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:44
p_rgh
volScalarField & p_rgh
Definition: setRegionFluidFields.H:15
pcorrTypes
wordList pcorrTypes(p.boundaryField().size(), zeroGradientFvPatchScalarField::typeName)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
U
U
Definition: pEqn.H:72
overlap
const cellCellStencilObject & overlap
Definition: correctPhi.H:57
cellTypes
const labelList & cellTypes
Definition: correctPhi.H:58
Foam::fac::laplacian
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:47
Foam::fvc::makeAbsolute
void makeAbsolute(surfaceScalarField &phi, const volVectorField &U)
Make the given flux absolute.
Definition: fvcMeshPhi.C:116
forAll
forAll(p.boundaryField(), patchi)
Definition: correctPhi.H:36
Foam::gMax
Type gMax(const FieldField< Field, Type > &f)
Definition: FieldFieldFunctions.C:592