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 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 
40  (
41  functionObject,
42  surfaceDistance,
43  dictionary
44  );
45 }
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
50 
51 Foam::functionObjects::surfaceDistance::surfaceDistance
52 (
53  const word& name,
54  const Time& runTime,
55  const dictionary& dict
56 )
57 :
59 {
60  read(dict);
61 
62  volScalarField* procFieldPtr
63  (
64  new volScalarField
65  (
66  IOobject
67  (
68  "surfaceDistance",
69  mesh_.time().timeName(),
70  mesh_,
71  IOobject::NO_READ,
72  IOobject::NO_WRITE
73  ),
74  mesh_,
76  )
77  );
78 
79  mesh_.objectRegistry::store(procFieldPtr);
80 }
81 
82 
83 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
84 
86 {}
87 
88 
89 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
90 
92 (
93  const dictionary& dict
94 )
95 {
97 
98  doCells_ = dict.lookupOrDefault("calculateCells", true);
99 
100  geomPtr_.clear();
101  geomPtr_.set
102  (
104  (
105  IOobject
106  (
107  "abc", // dummy name
108  mesh_.time().constant(), // directory
109  "triSurface", // instance
110  mesh_.time(), // registry
113  ),
114  dict.subDict("geometry"),
115  true // allow single-region shortcut
116  )
117  );
118 
119  return true;
120 }
121 
122 
124 {
125  volScalarField& distance = mesh_.lookupObjectRef<volScalarField>
126  (
127  "surfaceDistance"
128  );
129 
130  volScalarField::Boundary& bfld = distance.boundaryFieldRef();
131  forAll(bfld, patchi)
132  {
133  if (!polyPatch::constraintType(bfld[patchi].patch().type()))
134  {
135  const pointField& fc = mesh_.C().boundaryField()[patchi];
136 
137  labelList surfaces;
138  List<pointIndexHit> nearestInfo;
139  geomPtr_().findNearest
140  (
141  fc,
142  scalarField(fc.size(), GREAT),
143  surfaces,
144  nearestInfo
145  );
146 
147  scalarField dist(fc.size());
148  forAll(nearestInfo, i)
149  {
150  dist[i] = mag(nearestInfo[i].hitPoint()-fc[i]);
151  }
152  bfld[patchi] == dist;
153  }
154  }
155 
156  if (doCells_)
157  {
158  const pointField& cc = mesh_.C();
159 
160  labelList surfaces;
161  List<pointIndexHit> nearestInfo;
162  geomPtr_().findNearest
163  (
164  cc,
165  scalarField(cc.size(), GREAT),
166  surfaces,
167  nearestInfo
168  );
169 
170  forAll(nearestInfo, celli)
171  {
172  distance[celli] = mag(nearestInfo[celli].hitPoint()-cc[celli]);
173  }
174  }
175  distance.correctBoundaryConditions();
176 
177  return true;
178 }
179 
180 
182 {
183  Log << " functionObjects::" << type() << " " << name()
184  << " writing distance-to-surface field" << endl;
185 
186  const volScalarField& distance =
187  mesh_.lookupObject<volScalarField>("surfaceDistance");
188 
189 // volScalarField::Boundary& bfld = distance.boundaryFieldRef();
190 // forAll(bfld, patchi)
191 // {
192 // if (!polyPatch::constraintType(bfld[patchi].patch().type()))
193 // {
194 // const pointField& fc = mesh_.C().boundaryField()[patchi];
195 //
196 // labelList surfaces;
197 // List<pointIndexHit> nearestInfo;
198 // geomPtr_().findNearest
199 // (
200 // fc,
201 // scalarField(fc.size(), GREAT),
202 // surfaces,
203 // nearestInfo
204 // );
205 //
206 // scalarField dist(fc.size());
207 // forAll(nearestInfo, i)
208 // {
209 // dist[i] = mag(nearestInfo[i].hitPoint()-fc[i]);
210 // }
211 // bfld[patchi] == dist;
212 // }
213 // }
214 //
215 // if (doCells_)
216 // {
217 // const pointField& cc = mesh_.C();
218 //
219 // labelList surfaces;
220 // List<pointIndexHit> nearestInfo;
221 // geomPtr_().findNearest
222 // (
223 // cc,
224 // scalarField(cc.size(), GREAT),
225 // surfaces,
226 // nearestInfo
227 // );
228 //
229 // forAll(nearestInfo, celli)
230 // {
231 // distance[celli] = mag(nearestInfo[celli].hitPoint()-cc[celli]);
232 // }
233 // }
234 // distance.correctBoundaryConditions();
235  distance.write();
236 
237  return true;
238 }
239 
240 
241 // ************************************************************************* //
Foam::IOobject::NO_WRITE
Definition: IOobject.H:130
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::scalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Definition: primitiveFieldsFwd.H:52
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::functionObjects::surfaceDistance::write
virtual bool write()
Write the interpolated fields.
Definition: surfaceDistance.C:181
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
surfaceDistance.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::functionObjects::fvMeshFunctionObject
Specialization of Foam::functionObject for an Foam::fvMesh, providing a reference to the Foam::fvMesh...
Definition: fvMeshFunctionObject.H:64
Foam::functionObjects::addToRunTimeSelectionTable
addToRunTimeSelectionTable(functionObject, add, dictionary)
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::Field< vector >
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::functionObjects::surfaceDistance::execute
virtual bool execute()
Calculate the interpolated fields.
Definition: surfaceDistance.C:123
Foam::functionObjects::surfaceDistance::~surfaceDistance
virtual ~surfaceDistance()
Destructor.
Definition: surfaceDistance.C:85
Foam::polyPatch::constraintType
static bool constraintType(const word &pt)
Return true if the given type is a constraint type.
Definition: polyPatch.C:273
Foam::blockMeshTools::read
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
Definition: blockMeshTools.C:33
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::functionObjects::regionFunctionObject::read
virtual bool read(const dictionary &dict)
Read optional controls.
Definition: regionFunctionObject.C:166
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::distance
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Foam::functionObjects::defineTypeNameAndDebug
defineTypeNameAndDebug(fvMeshFunctionObject, 0)
Foam::foamVersion::patch
const std::string patch
OpenFOAM patch number as a std::string.
Foam::List< label >
Foam::searchableSurfaces
Container for searchableSurfaces. The collection is specified as a dictionary. For example,...
Definition: searchableSurfaces.H:92
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::functionObjects::surfaceDistance::read
virtual bool read(const dictionary &)
Read the controls.
Definition: surfaceDistance.C:92
Foam::GeometricField< scalar, fvPatchField, volMesh >
Log
#define Log
Report write to Foam::Info if the local log switch is true.
Definition: messageStream.H:332
Foam::IOobject::MUST_READ
Definition: IOobject.H:120