surfaceDistance.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) 2019-2022 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 "surfaceDistance.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace functionObjects
36 {
37  defineTypeNameAndDebug(surfaceDistance, 0);
38  addToRunTimeSelectionTable(functionObject, surfaceDistance, dictionary);
39 }
40 }
41 
42 
43 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
44 
46 (
47  const word& name,
48  const Time& runTime,
49  const dictionary& dict
50 )
51 :
52  fvMeshFunctionObject(name, runTime, dict)
53 {
54  read(dict);
55 
56  volScalarField* distPtr
57  (
58  new volScalarField
59  (
60  IOobject
61  (
62  "surfaceDistance",
63  mesh_.time().timeName(),
64  mesh_,
68  ),
69  mesh_,
71  )
72  );
73 
74  mesh_.objectRegistry::store(distPtr);
75 }
76 
77 
78 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
79 
81 (
82  const dictionary& dict
83 )
84 {
86 
87  doCells_ = dict.getOrDefault("calculateCells", true);
88 
89  geomPtr_.reset(nullptr);
90  geomPtr_.reset
91  (
93  (
94  IOobject
95  (
96  "abc", // dummy name
97  mesh_.time().constant(), // directory
98  "triSurface", // instance
99  mesh_.time(), // registry
102  ),
103  dict.subDict("geometry"),
104  true // allow single-region shortcut
105  )
106  );
107 
108  return true;
109 }
110 
111 
113 {
114  volScalarField& distance = mesh_.lookupObjectRef<volScalarField>
115  (
116  "surfaceDistance"
117  );
118 
119  volScalarField::Boundary& bfld = distance.boundaryFieldRef();
120  forAll(bfld, patchi)
121  {
122  if (!polyPatch::constraintType(bfld[patchi].patch().type()))
123  {
124  const pointField& fc = mesh_.C().boundaryField()[patchi];
125 
126  labelList surfaces;
127  List<pointIndexHit> nearestInfo;
128  geomPtr_().findNearest
129  (
130  fc,
131  scalarField(fc.size(), GREAT),
132  surfaces,
133  nearestInfo
134  );
135 
136  scalarField dist(fc.size());
137  forAll(nearestInfo, i)
138  {
139  dist[i] = nearestInfo[i].hitPoint().dist(fc[i]);
140  }
141  bfld[patchi] == dist;
142  }
143  }
144 
145  if (doCells_)
146  {
147  const pointField& cc = mesh_.C().internalField();
148 
149  labelList surfaces;
150  List<pointIndexHit> nearestInfo;
151  geomPtr_().findNearest
152  (
153  cc,
154  scalarField(cc.size(), GREAT),
155  surfaces,
156  nearestInfo
157  );
158 
159  forAll(nearestInfo, celli)
160  {
161  distance[celli] = nearestInfo[celli].hitPoint().dist(cc[celli]);
162  }
163  }
164  distance.correctBoundaryConditions();
165 
166  return true;
167 }
168 
169 
171 {
172  Log << " functionObjects::" << type() << " " << name()
173  << " writing distance-to-surface field" << endl;
174 
175  const volScalarField& distance =
176  mesh_.lookupObject<volScalarField>("surfaceDistance");
177 
178 // volScalarField::Boundary& bfld = distance.boundaryFieldRef();
179 // forAll(bfld, patchi)
180 // {
181 // if (!polyPatch::constraintType(bfld[patchi].patch().type()))
182 // {
183 // const pointField& fc = mesh_.C().boundaryField()[patchi];
184 //
185 // labelList surfaces;
186 // List<pointIndexHit> nearestInfo;
187 // geomPtr_().findNearest
188 // (
189 // fc,
190 // scalarField(fc.size(), GREAT),
191 // surfaces,
192 // nearestInfo
193 // );
194 //
195 // scalarField dist(fc.size());
196 // forAll(nearestInfo, i)
197 // {
198 // dist[i] = nearestInfo[i].hitPoint().dist(fc[i]);
199 // }
200 // bfld[patchi] == dist;
201 // }
202 // }
203 //
204 // if (doCells_)
205 // {
206 // const pointField& cc = mesh_.C().internalField();
207 //
208 // labelList surfaces;
209 // List<pointIndexHit> nearestInfo;
210 // geomPtr_().findNearest
211 // (
212 // cc,
213 // scalarField(cc.size(), GREAT),
214 // surfaces,
215 // nearestInfo
216 // );
217 //
218 // forAll(nearestInfo, celli)
219 // {
220 // distance[celli] = nearestInfo[celli].hitPoint().dist(cc[celli]);
221 // }
222 // }
223 // distance.correctBoundaryConditions();
224  distance.write();
225 
226  return true;
227 }
228 
229 
230 // ************************************************************************* //
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
virtual bool read(const dictionary &)
Read the controls.
defineTypeNameAndDebug(ObukhovLength, 0)
virtual bool write()
Write the interpolated fields.
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Ignore writing from objectRegistry::writeObject()
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Macros for easy insertion into run-time selection tables.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
surfaceDistance(const word &name, const Time &runTime, const dictionary &dict)
Construct for given objectRegistry and dictionary.
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Container for searchableSurfaces. The collection is specified as a dictionary. For example...
virtual bool execute()
Calculate the interpolated fields.
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
addToRunTimeSelectionTable(functionObject, ObukhovLength, dictionary)
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition: polyPatch.C:269
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Nothing to be read.
#define Log
Definition: PDRblock.C:28
const std::string patch
OpenFOAM patch number as a std::string.
virtual bool read(const dictionary &dict)
Read optional controls.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Request registration (bool: true)
const fvMesh & mesh_
Reference to the fvMesh.
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127