nearWallDistNoSearch.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 "nearWallDistNoSearch.H"
27 #include "fvMesh.H"
28 #include "wallFvPatch.H"
29 #include "surfaceFields.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
34 {
35  const volVectorField& cellCentres = mesh_.C();
36  const fvPatchList& patches = mesh_.boundary();
37 
38  forAll(patches, patchI)
39  {
40  fvPatchScalarField& ypatch = operator[](patchI);
41 
42  if (isA<wallFvPatch>(patches[patchI]))
43  {
44  const labelUList& faceCells = patches[patchI].faceCells();
45 
46  const fvPatchVectorField& patchCentres
47  = cellCentres.boundaryField()[patchI];
48 
49  const fvsPatchVectorField& Apatch
50  = mesh_.Sf().boundaryField()[patchI];
51 
52  const fvsPatchScalarField& magApatch
53  = mesh_.magSf().boundaryField()[patchI];
54 
55  forAll(patchCentres, facei)
56  {
57  ypatch[facei] =
58  (
59  Apatch[facei] &
60  (
61  patchCentres[facei]
62  - cellCentres[faceCells[facei]]
63  )
64  )/magApatch[facei];
65  }
66  }
67  else
68  {
69  ypatch = 0.0;
70  }
71  }
72 }
73 
74 
75 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
76 
78 :
79  volScalarField::GeometricBoundaryField
80  (
81  mesh.boundary(),
82  mesh.V(), // Dummy internal field
83  calculatedFvPatchScalarField::typeName
84  ),
85  mesh_(mesh)
86 {
87  doAll();
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
92 
94 {}
95 
96 
97 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
98 
100 {
101  if (mesh_.changing())
102  {
103  // Update size of GeometricBoundaryField
104  forAll(mesh_.boundary(), patchI)
105  {
106  operator[](patchI).setSize(mesh_.boundary()[patchI].size());
107  }
108  }
109 
110  doAll();
111 }
112 
113 
114 // ************************************************************************* //
Foam::fvPatchField
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:65
Foam::nearWallDistNoSearch::correct
virtual void correct()
Correct for mesh geom/topo changes.
Definition: nearWallDistNoSearch.C:99
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
boundary
faceListList boundary(nPatches)
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
wallFvPatch.H
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
surfaceFields.H
Foam::surfaceFields.
Foam::nearWallDistNoSearch::~nearWallDistNoSearch
virtual ~nearWallDistNoSearch()
Destructor.
Definition: nearWallDistNoSearch.C:93
Foam::fvMesh::magSf
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Definition: fvMeshGeometry.C:358
nearWallDistNoSearch.H
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::fvMesh::Sf
const surfaceVectorField & Sf() const
Return cell face area vectors.
Definition: fvMeshGeometry.C:347
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam::nearWallDistNoSearch::doAll
void doAll()
Do all calculations.
Definition: nearWallDistNoSearch.C:33
Foam::nearWallDistNoSearch::mesh_
const fvMesh & mesh_
Reference to mesh.
Definition: nearWallDistNoSearch.H:59
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:550
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:369
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
patches
patches[0]
Definition: createSingleCellMesh.H:36
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::nearWallDistNoSearch::nearWallDistNoSearch
nearWallDistNoSearch(const nearWallDistNoSearch &)
Disallow default bitwise copy construct.