wallDistData.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2015 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 // Construct from components
34 template<class TransferType>
36 (
37  const Foam::fvMesh& mesh,
39  const bool correctWalls
40 )
41 :
43  (
44  IOobject
45  (
46  "y",
47  mesh.time().timeName(),
48  mesh
49  ),
50  mesh,
51  dimensionedScalar("y", dimLength, GREAT)
52  ),
53  cellDistFuncs(mesh),
54  field_(field),
55  correctWalls_(correctWalls),
56  nUnset_(0)
57 {
58  correct();
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
64 template<class TransferType>
66 {}
67 
68 
69 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70 
71 // Correct for mesh geom/topo changes
72 template<class TransferType>
74 {
75  const polyMesh& mesh = cellDistFuncs::mesh();
76 
77  //
78  // Fill data on wall patches with initial values
79  //
80 
81  // Get patchids of walls
82  labelHashSet wallPatchIDs(getPatchIDs<wallPolyPatch>());
83 
84  // Collect pointers to data on patches
85  UPtrList<Field<Type> > patchData(mesh.boundaryMesh().size());
86 
87  forAll(field_.boundaryField(), patchI)
88  {
89  patchData.set(patchI, &field_.boundaryField()[patchI]);
90  }
91 
92  // Do mesh wave
94  (
95  mesh,
96  wallPatchIDs,
97  patchData,
98  correctWalls_
99  );
100 
101  // Transfer cell values from wave into *this and field_
102  transfer(wave.distance());
103 
104  field_.transfer(wave.cellData());
105 
106  // Transfer values on patches into boundaryField of *this and field_
107  forAll(boundaryField(), patchI)
108  {
109  scalarField& waveFld = wave.patchDistance()[patchI];
110 
111  if (!isA<emptyFvPatchScalarField>(boundaryField()[patchI]))
112  {
113  boundaryField()[patchI].transfer(waveFld);
114 
115  Field<Type>& wavePatchData = wave.patchData()[patchI];
116 
117  field_.boundaryField()[patchI].transfer(wavePatchData);
118  }
119  }
120 
121  // Transfer number of unset values
122  nUnset_ = wave.nUnset();
123 }
124 
125 
126 // ************************************************************************* //
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:52
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label size() const
Return the number of elements in the PtrList.
Definition: PtrListI.H:32
label nUnset() const
const FieldField< Field, Type > & patchData() const
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/dm2/henry/OpenFOAM/OpenFOAM-3.0.x/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:74
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:741
const scalarField & distance() const
const Field< Type > & cellData() const
dynamicFvMesh & mesh
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
#define forAll(list, i)
Definition: UList.H:421
const FieldField< Field, scalar > & patchDistance() const
virtual ~wallDistData()
Destructor.
Definition: wallDistData.C:65
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Generic GeometricField class.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Wall distance calculation. Like wallDist but also transports extra data (template argument)...
Definition: wallDistData.H:57
Collection of functions used in wall distance calculation.
Definition: cellDistFuncs.H:61
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Takes a set of patches to start MeshWave from.
Definition: patchDataWave.H:63
virtual void correct()
Correct for mesh geom/topo changes.
Definition: wallDistData.C:73