wallDistData.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2018 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "wallDistData.H"
27 #include "patchDataWave.H"
28 #include "wallPolyPatch.H"
29 #include "emptyFvPatchFields.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class TransferType>
35 (
36  const Foam::fvMesh& mesh,
38  const bool correctWalls
39 )
40 :
42  (
43  IOobject
44  (
45  "y",
46  mesh.time().timeName(),
47  mesh
48  ),
49  mesh,
50  dimensionedScalar("y", dimLength, great)
51  ),
52  cellDistFuncs(mesh),
53  field_(field),
54  correctWalls_(correctWalls),
55  nUnset_(0)
56 {
57  correct();
58 }
59 
60 
61 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
62 
63 template<class TransferType>
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
70 template<class TransferType>
72 {
73  const polyMesh& mesh = cellDistFuncs::mesh();
74 
75  //
76  // Fill data on wall patches with initial values
77  //
78 
79  // Get patchids of walls
80  labelHashSet wallPatchIDs(getPatchIDs<wallPolyPatch>());
81 
82  // Collect pointers to data on patches
83  UPtrList<Field<Type>> patchData(mesh.boundaryMesh().size());
84 
86  Boundary& fieldBf = field_.boundaryFieldRef();
87 
88  forAll(fieldBf, patchi)
89  {
90  patchData.set(patchi, &fieldBf[patchi]);
91  }
92 
93  // Do mesh wave
95  (
96  mesh,
97  wallPatchIDs,
98  patchData,
99  correctWalls_
100  );
101 
102  // Transfer cell values from wave into *this and field_
103  transfer(wave.distance());
104 
105  field_.transfer(wave.cellData());
106 
108  Boundary& bf = boundaryFieldRef();
109 
110  // Transfer values on patches into boundaryField of *this and field_
111  forAll(bf, patchi)
112  {
113  scalarField& waveFld = wave.patchDistance()[patchi];
114 
115  if (!isA<emptyFvPatchScalarField>(boundaryField()[patchi]))
116  {
117  bf[patchi].transfer(waveFld);
118  Field<Type>& wavePatchData = wave.patchData()[patchi];
119  fieldBf[patchi].transfer(wavePatchData);
120  }
121  }
122 
123  // Transfer number of unset values
124  nUnset_ = wave.nUnset();
125 }
126 
127 
128 // ************************************************************************* //
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:424
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
virtual ~wallDistData()
Destructor.
Definition: wallDistData.C:64
label nUnset() const
Generic GeometricField class.
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
Takes a set of patches to start MeshWave from.
Definition: patchDataWave.H:63
Collection of functions used in wall distance calculation.
Definition: cellDistFuncs.H:61
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dynamicFvMesh & mesh
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:54
const FieldField< Field, Type > & patchData() const
virtual void correct()
Correct for mesh geom/topo changes.
Definition: wallDistData.C:71
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/ubuntu/OpenFOAM-6/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:72
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
const scalarField & distance() const
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
void transfer(UPtrList< T > &)
Transfer the contents of the argument UPtrList into this.
Definition: UPtrList.C:94
const FieldField< Field, scalar > & patchDistance() const
Wall distance calculation. Like wallDist but also transports extra data (template argument)...
Definition: wallDistData.H:57
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
void transfer(List< Type > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
const Field< Type > & cellData() const