Loading...
Searching...
No Matches
searchableSurfaceControl.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) 2012-2017 OpenFOAM Foundation
9 Copyright (C) 2020-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
31#include "cellSizeFunction.H"
32#include "triSurfaceMesh.H"
33#include "searchableBox.H"
34#include "vectorTools.H"
35#include "quaternion.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
41
44(
48);
49
50}
51
52
53// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
54
55//Foam::tensor Foam::surfaceControl::requiredAlignment
56//(
57// const Foam::point& pt,
58// const vectorField& ptNormals
59//) const
60//{
94//
100//
101// vector np = ptNormals[0];
102//
103// const tensor Rp = rotationTensor(vector(0,0,1), np);
104//
105// vector na = Zero;
106//
107// scalar smallestAngle = GREAT;
108//
109// for (label pnI = 1; pnI < ptNormals.size(); ++pnI)
110// {
111// const vector& nextNormal = ptNormals[pnI];
112//
113// const scalar cosPhi = vectorTools::cosPhi(np, nextNormal);
114//
115// if (mag(cosPhi) < smallestAngle)
116// {
117// na = nextNormal;
118// smallestAngle = cosPhi;
119// }
120// }
121//
122// // Secondary alignment
123// vector ns = np ^ na;
124//
125// if (mag(ns) < SMALL)
126// {
127// WarningInFunction
128// << "Parallel normals detected in spoke search." << nl
129// << "point: " << pt << nl
130// << "np : " << np << nl
131// << "na : " << na << nl
132// << "ns : " << ns << nl
133// << endl;
134//
135// ns = np;
136// }
137//
138// ns /= mag(ns);
139//
140// tensor Rs = rotationTensor((Rp & vector(0,1,0)), ns);
141//
148//
149// return (Rs & Rp);
150//}
151
152
153// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
154
155Foam::searchableSurfaceControl::searchableSurfaceControl
156(
157 const Time& runTime,
158 const word& name,
159 const dictionary& controlFunctionDict,
160 const conformationSurfaces& geometryToConformTo,
161 const scalar& defaultCellSize
162)
163:
164 cellSizeAndAlignmentControl
165 (
166 runTime,
167 name,
168 controlFunctionDict,
169 geometryToConformTo,
170 defaultCellSize
171 ),
172 surfaceName_(controlFunctionDict.getOrDefault<word>("surface", name)),
173 searchableSurface_(geometryToConformTo.geometry()[surfaceName_]),
174 geometryToConformTo_(geometryToConformTo),
175 cellSizeFunctions_(1),
176 regionToCellSizeFunctions_(searchableSurface_.regions().size(), -1),
177 maxPriority_(-1)
178{
179 Info<< indent << "Master settings:" << endl;
181
182 cellSizeFunctions_.set
183 (
184 0,
185 cellSizeFunction::New
186 (
187 controlFunctionDict,
188 searchableSurface_,
189 defaultCellSize_,
190 labelList()
191 )
192 );
193
195
196 PtrList<cellSizeFunction> regionCellSizeFunctions;
197
198 DynamicList<label> defaultCellSizeRegions;
199
200 label nRegionCellSizeFunctions = 0;
201
202 // Loop over regions - if any entry is not specified they should
203 // inherit values from the parent surface.
204 if (controlFunctionDict.found("regions"))
205 {
206 const dictionary& regionsDict = controlFunctionDict.subDict("regions");
207 const wordList& regionNames = searchableSurface_.regions();
208
209 label nRegions = regionsDict.size();
210
211 regionCellSizeFunctions.setSize(nRegions);
212 defaultCellSizeRegions.setCapacity(nRegions);
213
214 forAll(regionNames, regionI)
215 {
216 const word& regionName = regionNames[regionI];
217
218 label regionID = geometryToConformTo_.geometry().findSurfaceRegionID
219 (
220 this->name(),
221 regionName
222 );
223
224 if (regionsDict.found(regionName))
225 {
226 // Get the dictionary for region
227 const dictionary& regionDict = regionsDict.subDict(regionName);
228
229 Info<< indent << "Region " << regionName
230 << " (ID = " << regionID << ")" << " settings:"
231 << endl;
232 Info<< incrIndent;
233
234 regionCellSizeFunctions.set
235 (
236 nRegionCellSizeFunctions,
237 cellSizeFunction::New
238 (
239 regionDict,
240 searchableSurface_,
241 defaultCellSize_,
242 labelList(1, regionID)
243 )
244 );
245 Info<< decrIndent;
246
247 regionToCellSizeFunctions_[regionID] = nRegionCellSizeFunctions;
248
249 nRegionCellSizeFunctions++;
250 }
251 else
252 {
253 // Add to default list
254 defaultCellSizeRegions.append(regionID);
255 }
256 }
257 }
258
259 if (defaultCellSizeRegions.empty() && !regionCellSizeFunctions.empty())
260 {
261 cellSizeFunctions_.transfer(regionCellSizeFunctions);
262 }
263 else if (nRegionCellSizeFunctions > 0)
264 {
265 regionCellSizeFunctions.setSize(nRegionCellSizeFunctions + 1);
266
267 regionCellSizeFunctions.set
268 (
269 nRegionCellSizeFunctions,
270 cellSizeFunction::New
271 (
272 controlFunctionDict,
273 searchableSurface_,
274 defaultCellSize_,
275 labelList()
276 )
277 );
278
279 const wordList& regionNames = searchableSurface_.regions();
280
281 forAll(regionNames, regionI)
282 {
283 if (regionToCellSizeFunctions_[regionI] == -1)
284 {
285 regionToCellSizeFunctions_[regionI] = nRegionCellSizeFunctions;
286 }
287 }
288
289 cellSizeFunctions_.transfer(regionCellSizeFunctions);
290 }
291 else
292 {
293 const wordList& regionNames = searchableSurface_.regions();
294
295 forAll(regionNames, regionI)
296 {
297 if (regionToCellSizeFunctions_[regionI] == -1)
298 {
299 regionToCellSizeFunctions_[regionI] = 0;
300 }
301 }
302 }
303
304
305 forAll(cellSizeFunctions_, funcI)
306 {
307 const label funcPriority = cellSizeFunctions_[funcI].priority();
308
309 if (funcPriority > maxPriority_)
310 {
311 maxPriority_ = funcPriority;
312 }
313 }
314
315 // Sort controlFunctions_ by maxPriority
316 SortableList<label> functionPriorities(cellSizeFunctions_.size());
317
318 forAll(cellSizeFunctions_, funcI)
319 {
320 functionPriorities[funcI] = cellSizeFunctions_[funcI].priority();
321 }
322
323 functionPriorities.reverseSort();
324
325 labelList invertedFunctionPriorities =
326 invert(functionPriorities.size(), functionPriorities.indices());
327
328 cellSizeFunctions_.reorder(invertedFunctionPriorities);
329
330 Info<< nl << indent << "There are " << cellSizeFunctions_.size()
331 << " region control functions" << endl;
332}
333
334
335// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
336
338{}
339
340
341// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
342
344(
346 scalarField& sizes,
347 triadField& alignments
348) const
349{
350 pts = searchableSurface_.points();
351 sizes.setSize(pts.size());
352 alignments.setSize(pts.size());
353
354 const scalar nearFeatDistSqrCoeff = 1e-8;
355
356 forAll(pts, pI)
357 {
358 // Is the point in the extendedFeatureEdgeMesh? If so get the
359 // point normal, otherwise get the surface normal from
360 // searchableSurface
361
362 pointIndexHit info;
363 label infoFeature;
364 geometryToConformTo_.findFeaturePointNearest
365 (
366 pts[pI],
367 nearFeatDistSqrCoeff,
368 info,
369 infoFeature
370 );
371
372 scalar limitedCellSize = GREAT;
373
374 autoPtr<triad> pointAlignment;
375
376 if (info.hit())
377 {
378 const extendedFeatureEdgeMesh& features =
379 geometryToConformTo_.features()[infoFeature];
380
381 vectorField norms = features.featurePointNormals(info.index());
382
383 // Create a triad from these norms.
384 pointAlignment.reset(new triad());
385 forAll(norms, nI)
386 {
387 pointAlignment() += norms[nI];
388 }
389
390 pointAlignment().normalize();
391 pointAlignment().orthogonalize();
392 }
393 else
394 {
395 geometryToConformTo_.findEdgeNearest
396 (
397 pts[pI],
398 nearFeatDistSqrCoeff,
399 info,
400 infoFeature
401 );
402
403 if (info.hit())
404 {
405 const extendedFeatureEdgeMesh& features =
406 geometryToConformTo_.features()[infoFeature];
407
408 vectorField norms = features.edgeNormals(info.index());
409
410 // Create a triad from these norms.
411 pointAlignment.reset(new triad());
412 forAll(norms, nI)
413 {
414 pointAlignment() += norms[nI];
415 }
416
417 pointAlignment().normalize();
418 pointAlignment().orthogonalize();
419 }
420 else
421 {
422 pointField ptField(1, pts[pI]);
423 scalarField distField(1, nearFeatDistSqrCoeff);
424 List<pointIndexHit> infoList(1, pointIndexHit());
425
426 searchableSurface_.findNearest(ptField, distField, infoList);
427
428 vectorField normals(1);
429 searchableSurface_.getNormal(infoList, normals);
430
431 if (mag(normals[0]) < SMALL)
432 {
433 normals[0] = vector::one;
434 }
435
436 pointAlignment.reset(new triad(normals[0]));
437
438 if (infoList[0].hit())
439 {
440 // Limit cell size
441 const vector vN =
442 infoList[0].point()
443 - 2.0*normals[0]*defaultCellSize_;
444
445 List<pointIndexHit> intersectionList;
446 searchableSurface_.findLineAny
447 (
448 ptField,
449 pointField(1, vN),
450 intersectionList
451 );
452 }
453
454// if (intersectionList[0].hit())
455// {
456// scalar dist = intersectionList[0].point().dist(pts[pI]);
457//
458// limitedCellSize = dist/2.0;
459// }
460 }
461 }
462
463 label priority = -1;
464 if (!cellSize(pts[pI], sizes[pI], priority))
465 {
466 sizes[pI] = defaultCellSize_;
467// FatalErrorInFunction
468// << "Could not calculate cell size"
469// << abort(FatalError);
470 }
471
472 sizes[pI] = min(limitedCellSize, sizes[pI]);
473
474 alignments[pI] = pointAlignment();
475 }
476}
477
478
480(
483) const
484{
485 const tmp<pointField> tpoints(searchableSurface_.points());
486 const pointField& points = tpoints();
487
488 const scalar nearFeatDistSqrCoeff = 1e-8;
489
490 pointField ptField(1, vector::min);
491 scalarField distField(1, nearFeatDistSqrCoeff);
492 List<pointIndexHit> infoList(1, pointIndexHit());
493
494 vectorField normals(1);
495 labelList region(1, label(-1));
496
497 forAll(points, pI)
498 {
499 // Is the point in the extendedFeatureEdgeMesh? If so get the
500 // point normal, otherwise get the surface normal from
501 // searchableSurface
502 ptField[0] = points[pI];
503
504 searchableSurface_.findNearest(ptField, distField, infoList);
505
506 if (infoList[0].hit())
507 {
508 searchableSurface_.getNormal(infoList, normals);
509 searchableSurface_.getRegion(infoList, region);
510
511 const cellSizeFunction& sizeFunc =
512 sizeFunctions()[regionToCellSizeFunctions_[region[0]]];
513
514 pointField extraPts;
515 scalarField extraSizes;
516 sizeFunc.sizeLocations
517 (
518 infoList[0],
519 normals[0],
520 extraPts,
521 extraSizes
522 );
523
524 pts.append(extraPts);
525 sizes.append(extraSizes);
526 }
527 }
528}
529
530
532(
533 const Foam::point& pt,
534 scalar& cellSize,
535 label& priority
536) const
537{
538 bool anyFunctionFound = false;
539
540 forAll(sizeFunctions(), funcI)
541 {
542 const cellSizeFunction& sizeFunc = sizeFunctions()[funcI];
543
544 if (sizeFunc.priority() < priority)
545 {
546 continue;
547 }
548
549 scalar sizeI = -1;
550
551 if (sizeFunc.cellSize(pt, sizeI))
552 {
553 anyFunctionFound = true;
554
555 if (sizeFunc.priority() == priority)
556 {
557 if (sizeI < cellSize)
558 {
559 cellSize = sizeI;
560 }
561 }
562 else
563 {
564 cellSize = sizeI;
565
566 priority = sizeFunc.priority();
567 }
568
569 if (debug)
570 {
571 Info<< " sizeI " << sizeI
572 <<" minSize " << cellSize << endl;
573 }
574 }
575 }
576
577 return anyFunctionFound;
578}
579
580
581// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
if(patchID !=-1)
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
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
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
void setSize(label n)
Alias for resize().
Definition List.H:536
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static const Form one
static const Form min
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
Abstract base class for specifying target cell sizes.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
virtual void cellSizeFunctionVertices(DynamicList< Foam::point > &pts, DynamicList< scalar > &sizes) const
~searchableSurfaceControl()
Destructor.
bool cellSize(const Foam::point &pt, scalar &cellSize, label &priority) const
virtual void initialVertices(pointField &pts, scalarField &sizes, triadField &alignments) const
A class for managing temporary objects.
Definition tmp.H:75
Representation of a 3D Cartesian coordinate system as a Vector of row vectors.
Definition triad.H:63
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
engineTime & runTime
auto & name
wordList regionNames
const pointField & points
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< word > wordList
List of word.
Definition fileName.H:60
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition Ostream.H:490
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
Field< triad > triadField
Specialisation of Field<T> for triad.
labelList invert(const label len, const labelUList &map)
Create an inverse one-to-one mapping.
Definition ListOps.C:28
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition Ostream.H:499
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
const pointField & pts
volScalarField & e
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299