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