directionalMeshWavePatchDistMethod.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) 2020,2023 OpenCFD Ltd.
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 "patchDataWave.H"
30 #include "fvMesh.H"
31 #include "volFields.H"
33 #include "emptyFvPatchFields.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace patchDistMethods
41 {
42  defineTypeNameAndDebug(directionalMeshWave, 0);
44  (
45  patchDistMethod,
46  directionalMeshWave,
47  dictionary
48  );
49 }
50 }
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
54 Foam::patchDistMethods::directionalMeshWave::directionalMeshWave
55 (
56  const dictionary& dict,
57  const fvMesh& mesh,
58  const labelHashSet& patchIDs
59 )
60 :
61  Foam::patchDistMethods::meshWave(dict, mesh, patchIDs),
62  n_(dict.get<vector>("normal"))
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
69 {
70  y = dimensionedScalar(dimLength, GREAT);
71 
73  (
74  IOobject
75  (
76  "nWall",
77  y.time().timeName(),
78  mesh_,
82  ),
83  mesh_,
85  patchDistMethod::patchTypes<scalar>(mesh_, patchIDs_)
86  );
87 
88  const fvPatchList& patches = mesh_.boundary();
89 
90  volVectorField::Boundary& nbf = n.boundaryFieldRef();
91 
92  for (const label patchi : patchIDs_)
93  {
94  nbf[patchi] == patches[patchi].nf();
95  }
96 
97  return correct(y, n);
98 }
99 
100 
102 (
103  volScalarField& y,
105 )
106 {
107  y = dimensionedScalar(dimLength, GREAT);
108 
109  // Collect pointers to data on patches
110  UPtrList<vectorField> patchData(mesh_.boundaryMesh().size());
111 
112  volVectorField::Boundary& nbf = n.boundaryFieldRef();
113 
114  forAll(nbf, patchi)
115  {
116  patchData.set(patchi, &nbf[patchi]);
117  }
118 
119  // Do mesh wave
120  vector testDirection(n_);
121 
122  patchDataWave<directionalWallPointData<vector>, vector> wave
123  (
124  mesh_,
125  patchIDs_,
126  patchData,
127  correctWalls_,
128  testDirection
129  );
130 
131  // Transfer cell values from wave into y and n
132  y.transfer(wave.distance());
133 
134  n.transfer(wave.cellData());
135 
136  // Transfer values on patches into boundaryField of y and n
137  volScalarField::Boundary& ybf = y.boundaryFieldRef();
138 
139  forAll(ybf, patchi)
140  {
141  scalarField& waveFld = wave.patchDistance()[patchi];
142 
143  if (!isA<emptyFvPatchScalarField>(ybf[patchi]))
144  {
145  ybf[patchi].transfer(waveFld);
146 
147  vectorField& wavePatchData = wave.patchData()[patchi];
148 
149  nbf[patchi].transfer(wavePatchData);
150  }
151  }
152 
153  // Make sure boundary values are up-to-date
154  y.correctBoundaryConditions();
155  n.correctBoundaryConditions();
156 
157  // Transfer number of unset values
158  this->nUnset_ = wave.nUnset();
159 
160  return this->nUnset_ > 0;
161 }
162 
163 
164 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
const labelList patchIDs(pbm.indices(polyPatchNames, true))
virtual bool correct(volScalarField &y)
Correct the given distance-to-patch field.
Fast topological mesh-wave method for calculating the distance to nearest patch for all cells and bou...
dictionary dict
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition: List.C:326
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
addToRunTimeSelectionTable(patchDistMethod, advectionDiffusion, dictionary)
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
Macros for easy insertion into run-time selection tables.
const fvMesh & mesh_
Reference to the mesh.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
scalar y
dynamicFvMesh & mesh
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void transfer(UPtrList< T > &list)
Transfer contents into this list and annul the argument.
Definition: UPtrListI.H:216
Vector< scalar > vector
Definition: vector.H:57
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: HashTable.H:106
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
defineTypeNameAndDebug(advectionDiffusion, 0)
const polyBoundaryMesh & patches
Nothing to be read.
const labelHashSet patchIDs_
Set of patch IDs.
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
label n
Field< vector > vectorField
Specialisation of Field<T> for vector.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Do not request registration (bool: false)
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127