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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "nearWallDistNoSearch.H"
29 #include "fvMesh.H"
30 #include "wallFvPatch.H"
31 #include "surfaceFields.H"
32 
33 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34 
35 void Foam::nearWallDistNoSearch::doAll()
36 {
37  const volVectorField& cellCentres = mesh_.C();
38  const fvPatchList& patches = mesh_.boundary();
39 
40  forAll(patches, patchi)
41  {
42  fvPatchScalarField& ypatch = operator[](patchi);
43 
44  if (isA<wallFvPatch>(patches[patchi]))
45  {
46  const labelUList& faceCells = patches[patchi].faceCells();
47 
48  const fvPatchVectorField& patchCentres
49  = cellCentres.boundaryField()[patchi];
50 
51  const fvsPatchVectorField& Apatch
52  = mesh_.Sf().boundaryField()[patchi];
53 
54  const fvsPatchScalarField& magApatch
55  = mesh_.magSf().boundaryField()[patchi];
56 
57  forAll(patchCentres, facei)
58  {
59  ypatch[facei] =
60  (
61  Apatch[facei] &
62  (
63  patchCentres[facei]
64  - cellCentres[faceCells[facei]]
65  )
66  )/magApatch[facei];
67  }
68  }
69  else
70  {
71  ypatch = 0.0;
72  }
73  }
74 }
75 
76 
77 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
78 
79 Foam::nearWallDistNoSearch::nearWallDistNoSearch(const Foam::fvMesh& mesh)
80 :
81  volScalarField::Boundary
82  (
83  mesh.boundary(),
84  mesh.V(), // Dummy internal field
85  calculatedFvPatchScalarField::typeName
86  ),
87  mesh_(mesh)
88 {
89  doAll();
90 }
91 
92 
93 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
94 
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
102 {
103  if (mesh_.changing())
104  {
105  // Update size of Boundary
106  forAll(mesh_.boundary(), patchi)
107  {
108  operator[](patchi).setSize(mesh_.boundary()[patchi].size());
109  }
110  }
111 
112  doAll();
113 }
114 
115 
116 // ************************************************************************* //
Foam::fvPatchScalarField
fvPatchField< scalar > fvPatchScalarField
Definition: fvPatchFieldsFwd.H:40
Foam::nearWallDistNoSearch::correct
virtual void correct()
Correct for mesh geom/topo changes.
Definition: nearWallDistNoSearch.C:101
Foam::fvsPatchScalarField
fvsPatchField< scalar > fvsPatchScalarField
Definition: fvsPatchFieldsFwd.H:45
wallFvPatch.H
surfaceFields.H
Foam::surfaceFields.
Foam::nearWallDistNoSearch::~nearWallDistNoSearch
virtual ~nearWallDistNoSearch()
Destructor.
Definition: nearWallDistNoSearch.C:95
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::fvPatchList
PtrList< fvPatch > fvPatchList
container classes for fvPatch
Definition: fvPatchList.H:47
nearWallDistNoSearch.H
Foam::fvsPatchVectorField
fvsPatchField< vector > fvsPatchVectorField
Definition: fvsPatchFieldsFwd.H:48
Foam::fvMesh::magSf
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
Definition: fvMeshGeometry.C:330
Foam::fvMesh::C
const volVectorField & C() const
Return cell centres as volVectorField.
Definition: fvMeshGeometry.C:341
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
fvMesh.H
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:62
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:685
Foam::fvPatchVectorField
fvPatchField< vector > fvPatchVectorField
Definition: fvPatchFieldsFwd.H:43
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::labelUList
UList< label > labelUList
A UList of labels.
Definition: UList.H:85
Foam::GeometricField< scalar, fvPatchField, volMesh >
boundary
faceListList boundary
Definition: createBlockMesh.H:4
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::fvMesh::Sf
const surfaceVectorField & Sf() const
Return cell face area vectors.
Definition: fvMeshGeometry.C:319