Loading...
Searching...
No Matches
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-2024 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"
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:
52 fvMeshFunctionObject(name, runTime, dict)
53{
54 read(dict);
55
56 volScalarField* distPtr
57 (
59 (
61 (
62 "surfaceDistance",
64 mesh_.thisDb(),
68 ),
69 mesh_,
71 )
72 );
74 regIOobject::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 (
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
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 }
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// ************************************************************************* //
#define Log
Definition PDRblock.C:28
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition Time.H:75
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
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Abstract base-class for Time/database function objects.
virtual bool read(const dictionary &dict)
Read and set the function object if its data have changed.
const fvMesh & mesh_
Reference to the fvMesh.
fvMeshFunctionObject(const fvMeshFunctionObject &)=delete
No copy construct.
Computes the distance to the nearest surface from a given geometry.
autoPtr< searchableSurfaces > geomPtr_
Geometry.
surfaceDistance(const word &name, const Time &runTime, const dictionary &dict)
Construct for given objectRegistry and dictionary.
virtual bool read(const dictionary &dict)
Read the function-object dictionary.
Switch doCells_
Switch to calculate distance-to-cells.
virtual bool execute()
Execute the function-object operations.
virtual bool write()
Write the function-object results.
virtual const objectRegistry & thisDb() const
Return the object registry - resolve conflict polyMesh/lduMesh.
Definition fvMesh.H:376
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
static bool constraintType(const word &patchType)
Return true if the given type is a constraint type.
Definition polyPatch.C:255
bool store()
Register object with its registry and transfer ownership to the registry.
Container for searchableSurfaces. The collection is specified as a dictionary. For example,...
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
engineTime & runTime
auto & name
Function objects are OpenFOAM utilities to ease workflow configurations and enhance workflows.
Namespace for OpenFOAM.
scalar distance(const vector &p1, const vector &p2)
Definition curveTools.C:12
List< label > labelList
A List of labels.
Definition List.H:62
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299