Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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
29#include "fvMesh.H"
30#include "wallFvPatch.H"
31#include "surfaceFields.H"
32
33// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34
35void 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
79Foam::nearWallDistNoSearch::nearWallDistNoSearch(const Foam::fvMesh& mesh)
80:
81 volScalarField::Boundary
82 (
83 mesh.boundary(),
84 mesh.V(), // Dummy internal field
86 ),
87 mesh_(mesh)
89 doAll();
90}
91
92
93// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
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// ************************************************************************* //
const Type & operator[](const labelPair &index) const
Definition FieldField.C:381
GeometricBoundaryField< Type, PatchField, GeoMesh > Boundary
This boundary field type.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Template invariant parts for fvPatchField.
virtual void correct()
Correct for mesh geom/topo changes.
virtual ~nearWallDistNoSearch()
Destructor.
const polyBoundaryMesh & patches
faceListList boundary
dynamicFvMesh & mesh
const expr V(m.psi().mesh().V())
const word calculatedType
A calculated patch field type.
GeometricField< vector, fvPatchField, volMesh > volVectorField
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
Definition fvPatch.H:64
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
UList< label > labelUList
A UList of labels.
Definition UList.H:75
fvsPatchField< vector > fvsPatchVectorField
fvPatchField< vector > fvPatchVectorField
fvsPatchField< scalar > fvsPatchScalarField
fvPatchField< scalar > fvPatchScalarField
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.