sampledIsoSurface.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2021 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
29#include "sampledIsoSurface.H"
30#include "dictionary.H"
31#include "fvMesh.H"
32#include "volFields.H"
35#include "PtrList.H"
36
37// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38
39namespace Foam
40{
43 (
46 word,
47 isoSurface
48 );
49}
50
51// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
52
53void Foam::sampledIsoSurface::getIsoFields() const
54{
55 const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
56
57 // Get volField
58 // ~~~~~~~~~~~~
59
60 volFieldPtr_ = fvm.getObjectPtr<volScalarField>(isoField_);
61
62 if (volFieldPtr_)
63 {
65 << "Lookup volField " << isoField_ << endl;
66
67 storedVolFieldPtr_.clear();
68 }
69 else
70 {
71 // Bit of a hack. Read field and store.
72
74 << "Checking " << isoField_
75 << " for same time " << fvm.time().timeName() << endl;
76
77 if
78 (
79 !storedVolFieldPtr_
80 || (fvm.time().timeName() != storedVolFieldPtr_().instance())
81 )
82 {
84 << "Reading volField " << isoField_
85 << " from time " << fvm.time().timeName() << endl;
86
87 IOobject vfHeader
88 (
89 isoField_,
90 fvm.time().timeName(),
91 fvm,
94 false
95 );
96
97 if (vfHeader.typeHeaderOk<volScalarField>(true))
98 {
99 storedVolFieldPtr_.reset
100 (
102 (
103 vfHeader,
104 fvm
105 )
106 );
107 volFieldPtr_ = storedVolFieldPtr_.get(); // get(), not release()
108 }
109 else
110 {
112 << "Cannot find isosurface field " << isoField_
113 << " in database or directory " << vfHeader.path()
114 << exit(FatalError);
115 }
116 }
117 }
118
119
120 // Get pointField
121 // ~~~~~~~~~~~~~~
122
123 // In case of multiple iso values we don't want to calculate multiple e.g.
124 // "volPointInterpolate(p)" so register it and re-use it. This is the
125 // same as the 'cache' functionality from volPointInterpolate but
126 // unfortunately that one does not guarantee that the field pointer
127 // remain: e.g. some other functionObject might delete the cached version.
128 // (volPointInterpolation::interpolate with cache=false deletes any
129 // registered one or if mesh.changing())
130
131 if (!subMeshPtr_)
132 {
133 const word pointFldName =
134 "volPointInterpolate_"
135 + type()
136 + "("
137 + isoField_
138 + ')';
139
140
141 pointFieldPtr_ = fvm.getObjectPtr<pointScalarField>(pointFldName);
142
143 if (pointFieldPtr_)
144 {
146 << "lookup pointField " << pointFldName << endl;
147
148 if (!pointFieldPtr_->upToDate(*volFieldPtr_))
149 {
151 << "updating pointField " << pointFldName << endl;
152
153 // Update the interpolated value
154 volPointInterpolation::New(fvm).interpolate
155 (
156 *volFieldPtr_,
157 const_cast<pointScalarField&>(*pointFieldPtr_)
158 );
159 }
160 }
161 else
162 {
163 // Not in registry. Interpolate.
164
166 << "creating pointField " << pointFldName << endl;
167
168 // Interpolate without cache. Note that we're registering it
169 // below so next time round it goes into the condition
170 // above.
171 pointFieldPtr_ =
172 volPointInterpolation::New(fvm).interpolate
173 (
174 *volFieldPtr_,
175 pointFldName,
176 false
177 ).ptr();
178
179 const_cast<pointScalarField*>(pointFieldPtr_)->store();
180 }
181
182
183 // If averaging redo the volField.
184 // Can only be done now since needs the point field.
185 if (average_)
186 {
187 storedVolFieldPtr_.reset
188 (
189 pointAverage(*pointFieldPtr_).ptr()
190 );
191 volFieldPtr_ = storedVolFieldPtr_.get(); // get(), not release()
192 }
193
194
196 << "volField " << volFieldPtr_->name()
197 << " min:" << min(*volFieldPtr_).value()
198 << " max:" << max(*volFieldPtr_).value() << nl
199 << "pointField " << pointFieldPtr_->name()
200 << " min:" << gMin(pointFieldPtr_->primitiveField())
201 << " max:" << gMax(pointFieldPtr_->primitiveField()) << endl;
202 }
203 else
204 {
205 // Get subMesh variants
206 const fvMesh& subFvm = subMeshPtr_().subMesh();
207
208 // Either lookup on the submesh or subset the whole-mesh volField
209
210 volSubFieldPtr_ = subFvm.getObjectPtr<volScalarField>(isoField_);
211
212 if (volSubFieldPtr_)
213 {
215 << "Sub-mesh lookup volField " << isoField_ << endl;
216
217 storedVolSubFieldPtr_.clear();
218 }
219 else
220 {
222 << "Sub-setting volField " << isoField_ << endl;
223
224 storedVolSubFieldPtr_.reset
225 (
226 subMeshPtr_().interpolate
227 (
228 *volFieldPtr_
229 ).ptr()
230 );
231 storedVolSubFieldPtr_->checkOut();
232 volSubFieldPtr_ = storedVolSubFieldPtr_.get(); // not release()
233 }
234
235
236 // The point field on subMesh
237
238 const word pointFldName =
239 "volPointInterpolate_"
240 + type()
241 + "("
242 + volSubFieldPtr_->name()
243 + ')';
244
245
246 pointFieldPtr_ = subFvm.getObjectPtr<pointScalarField>(pointFldName);
247
248 if (pointFieldPtr_)
249 {
251 << "Sub-mesh lookup pointField " << pointFldName << endl;
252
253 if (!pointFieldPtr_->upToDate(*volSubFieldPtr_))
254 {
256 << "Updating submesh pointField " << pointFldName << endl;
257
258 // Update the interpolated value
259 volPointInterpolation::New(subFvm).interpolate
260 (
261 *volSubFieldPtr_,
262 const_cast<pointScalarField&>(*pointFieldPtr_)
263 );
264 }
265 }
266 else
267 {
269 << "Interpolating submesh volField "
270 << volSubFieldPtr_->name()
271 << " to get submesh pointField " << pointFldName << endl;
272
273 pointSubFieldPtr_ =
275 (
276 subFvm
277 ).interpolate(*volSubFieldPtr_).ptr();
278
279 const_cast<pointScalarField*>(pointSubFieldPtr_)->store();
280 }
281
282
283 // If averaging redo the volField. Can only be done now since needs the
284 // point field.
285 if (average_)
286 {
287 storedVolSubFieldPtr_.reset
288 (
289 pointAverage(*pointSubFieldPtr_).ptr()
290 );
291 volSubFieldPtr_ = storedVolSubFieldPtr_.get(); // not release()
292 }
293
294
296 << "volSubField "
297 << volSubFieldPtr_->name()
298 << " min:" << min(*volSubFieldPtr_).value()
299 << " max:" << max(*volSubFieldPtr_).value() << nl
300 << "pointSubField "
301 << pointSubFieldPtr_->name()
302 << " min:" << gMin(pointSubFieldPtr_->primitiveField())
303 << " max:" << gMax(pointSubFieldPtr_->primitiveField()) << endl;
304 }
305}
306
307
308void Foam::sampledIsoSurface::combineSurfaces
309(
310 PtrList<isoSurfaceBase>& isoSurfPtrs
311)
312{
313 isoSurfacePtr_.reset(nullptr);
314
315 // Already checked previously for ALGO_POINT, but do it again
316 // - ALGO_POINT still needs fields (for interpolate)
317 // The others can do straight transfer
318 if
319 (
320 isoParams_.algorithm() == isoSurfaceParams::ALGO_POINT
321 && isoSurfPtrs.size() == 1
322 )
323 {
324 // Shift from list to autoPtr
325 isoSurfacePtr_.reset(isoSurfPtrs.release(0));
326 }
327 else if (isoSurfPtrs.size() == 1)
328 {
329 autoPtr<isoSurfaceBase> ptr(isoSurfPtrs.release(0));
330 auto& surf = *ptr;
331
332 surface_.transfer(static_cast<meshedSurface&>(surf));
333 meshCells_.transfer(surf.meshCells());
334 }
335 else
336 {
337 // Combine faces with point offsets
338 //
339 // Note: use points().size() from surface, not nPoints()
340 // since there may be uncompacted dangling nodes
341
342 label nFaces = 0, nPoints = 0;
343
344 for (const auto& surf : isoSurfPtrs)
345 {
346 nFaces += surf.size();
347 nPoints += surf.points().size();
348 }
349
350 faceList newFaces(nFaces);
351 pointField newPoints(nPoints);
352 meshCells_.resize(nFaces);
353
354 surfZoneList newZones(isoSurfPtrs.size());
355
356 nFaces = 0;
357 nPoints = 0;
358 forAll(isoSurfPtrs, surfi)
359 {
360 autoPtr<isoSurfaceBase> ptr(isoSurfPtrs.release(surfi));
361 auto& surf = *ptr;
362
363 SubList<face> subFaces(newFaces, surf.size(), nFaces);
364 SubList<point> subPoints(newPoints, surf.points().size(), nPoints);
365 SubList<label> subCells(meshCells_, surf.size(), nFaces);
366
367 newZones[surfi] = surfZone
368 (
370 subFaces.size(), // size
371 nFaces, // start
372 surfi // index
373 );
374
375 subFaces = surf.surfFaces();
376 subPoints = surf.points();
377 subCells = surf.meshCells();
378
379 if (nPoints)
380 {
381 for (face& f : subFaces)
382 {
383 for (label& pointi : f)
384 {
385 pointi += nPoints;
386 }
387 }
388 }
389
390 nFaces += subFaces.size();
391 nPoints += subPoints.size();
392 }
393
394 meshedSurface combined
395 (
396 std::move(newPoints),
397 std::move(newFaces),
398 newZones
399 );
400
401 surface_.transfer(combined);
402 }
403
404 // Addressing into the full mesh
405 if (subMeshPtr_ && meshCells_.size())
406 {
407 meshCells_ =
408 UIndirectList<label>(subMeshPtr_->cellMap(), meshCells_);
409 }
410}
411
412
413bool Foam::sampledIsoSurface::updateGeometry() const
414{
415 const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
416
417 // No update needed
418 if (fvm.time().timeIndex() == prevTimeIndex_)
419 {
420 return false;
421 }
422
423 prevTimeIndex_ = fvm.time().timeIndex();
424
425 // Clear any previously stored topologies
426 surface_.clear();
427 meshCells_.clear();
428 isoSurfacePtr_.reset(nullptr);
429
430 // Clear derived data
432
433 const bool hasCellZones =
434 (-1 != mesh().cellZones().findIndex(zoneNames_));
435
436 // Geometry
437 if
438 (
439 simpleSubMesh_
440 && isoParams_.algorithm() != isoSurfaceParams::ALGO_POINT
441 )
442 {
443 subMeshPtr_.reset(nullptr);
444
445 // Handle cell zones as inverse (blocked) selection
446 if (!ignoreCellsPtr_)
447 {
448 ignoreCellsPtr_.reset(new bitSet);
449
450 if (hasCellZones)
451 {
452 bitSet select(mesh().cellZones().selection(zoneNames_));
453
454 if (select.any() && !select.all())
455 {
456 // From selection to blocking
457 select.flip();
458
459 *ignoreCellsPtr_ = std::move(select);
460 }
461 }
462 }
463 }
464 else
465 {
466 // A standard subMesh treatment
467
468 if (ignoreCellsPtr_)
469 {
470 ignoreCellsPtr_->clearStorage();
471 }
472 else
473 {
474 ignoreCellsPtr_.reset(new bitSet);
475 }
476
477 // Get sub-mesh if any
478 if (!subMeshPtr_ && hasCellZones)
479 {
480 const label exposedPatchi =
481 mesh().boundaryMesh().findPatchID(exposedPatchName_);
482
484 << "Allocating subset of size "
485 << mesh().cellZones().selection(zoneNames_).count()
486 << " with exposed faces into patch "
487 << exposedPatchi << endl;
488
489 subMeshPtr_.reset
490 (
491 new fvMeshSubset
492 (
493 fvm,
494 mesh().cellZones().selection(zoneNames_),
495 exposedPatchi
496 )
497 );
498 }
499 }
500
501
502 // The fields
503 getIsoFields();
504
505 refPtr<volScalarField> tvolFld(*volFieldPtr_);
506 refPtr<pointScalarField> tpointFld(*pointFieldPtr_);
507
508 if (subMeshPtr_)
509 {
510 tvolFld.cref(*volSubFieldPtr_);
511 tpointFld.cref(*pointSubFieldPtr_);
512 }
513
514
515 // Create surfaces for each iso level
516
517 PtrList<isoSurfaceBase> isoSurfPtrs(isoValues_.size());
518
519 forAll(isoValues_, surfi)
520 {
521 isoSurfPtrs.set
522 (
523 surfi,
525 (
526 isoParams_,
527 tvolFld(),
528 tpointFld().primitiveField(),
529 isoValues_[surfi],
530 *ignoreCellsPtr_
531 )
532 );
533 }
534
535 // And flatten
536 const_cast<sampledIsoSurface&>(*this)
537 .combineSurfaces(isoSurfPtrs);
538
539
540 // triangulate uses remapFaces()
541 // - this is somewhat less efficient since it recopies the faces
542 // that we just created, but we probably don't want to do this
543 // too often anyhow.
544 if
545 (
546 triangulate_
547 && surface_.size()
548 && (isoParams_.algorithm() == isoSurfaceParams::ALGO_TOPO)
549 )
550 {
552 surface_.triangulate(faceMap);
553 meshCells_ = UIndirectList<label>(meshCells_, faceMap)();
554 }
555
556 if (debug)
557 {
558 Pout<< "isoSurface::updateGeometry() : constructed iso:" << nl
559 << " field : " << isoField_ << nl
560 << " value : " << flatOutput(isoValues_) << nl
561 << " average : " << Switch(average_) << nl
562 << " filter : "
563 << Switch(bool(isoParams_.filter())) << nl
564 << " bounds : " << isoParams_.getClipBounds() << nl;
565 if (subMeshPtr_)
566 {
567 Pout<< " zone size : "
568 << subMeshPtr_->subMesh().nCells() << nl;
569 }
570 Pout<< " points : " << points().size() << nl
571 << " faces : " << surface().size() << nl
572 << " cut cells : " << meshCells().size()
573 << endl;
574 }
575
576 return true;
577}
578
579
580// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
581
583(
584 const isoSurfaceParams& params,
585 const word& name,
586 const polyMesh& mesh,
587 const dictionary& dict
588)
589:
591 isoField_(dict.get<word>("isoField")),
592 isoValues_(),
593 isoParams_(dict, params),
594 average_(dict.getOrDefault("average", false)),
595 triangulate_(dict.getOrDefault("triangulate", false)),
596 simpleSubMesh_(dict.getOrDefault("simpleSubMesh", false)),
597 zoneNames_(),
598 exposedPatchName_(),
599 prevTimeIndex_(-1),
600
601 surface_(),
602 meshCells_(),
603 isoSurfacePtr_(nullptr),
604
605 subMeshPtr_(nullptr),
606 ignoreCellsPtr_(nullptr),
607
608 storedVolFieldPtr_(nullptr),
609 volFieldPtr_(nullptr),
610 pointFieldPtr_(nullptr),
611
612 storedVolSubFieldPtr_(nullptr),
613 volSubFieldPtr_(nullptr),
614 pointSubFieldPtr_(nullptr)
615{
617 {
618 // Forced use of specified algorithm (ignore dictionary entry)
619 isoParams_.algorithm(params.algorithm());
620 }
621
622 // The isoValues or isoValue
623
624 if (!dict.readIfPresent("isoValues", isoValues_))
625 {
626 isoValues_.resize(1);
627 dict.readEntry("isoValue", isoValues_.first());
628 }
629
630 if (isoValues_.empty())
631 {
633 << "No isoValue or isoValues specified." << nl
634 << exit(FatalIOError);
635 }
636
637 if (isoValues_.size() > 1)
638 {
639 const label nOrig = isoValues_.size();
640
641 inplaceUniqueSort(isoValues_);
642
643 if (nOrig != isoValues_.size())
644 {
646 << "Removed non-unique isoValues" << nl;
647 }
648 }
649
650 if (isoParams_.algorithm() == isoSurfaceParams::ALGO_POINT)
651 {
652 // Not possible for ALGO_POINT
653 simpleSubMesh_ = false;
654
655 // Previously emitted an error about using ALGO_POINT with
656 // non-interpolated, but that was before we had "sampleScheme"
657 // at the top level
658
659 if (isoValues_.size() > 1)
660 {
662 << "Multiple values on iso-surface (point) not supported"
663 << " since needs original interpolators." << nl
664 << exit(FatalIOError);
665 }
666 }
667
668 if (isoParams_.algorithm() == isoSurfaceParams::ALGO_TOPO)
669 {
670 if
671 (
672 triangulate_
673 && (isoParams_.filter() == isoSurfaceParams::filterType::NONE)
674 )
675 {
677 << "Cannot triangulate without a regularise filter" << nl
678 << exit(FatalIOError);
679 }
680 }
681
682
683 // Zones
684
685 if (!dict.readIfPresent("zones", zoneNames_) && dict.found("zone"))
686 {
687 zoneNames_.resize(1);
688 dict.readEntry("zone", zoneNames_.first());
689 }
690
691 if (-1 != mesh.cellZones().findIndex(zoneNames_))
692 {
693 dict.readIfPresent("exposedPatchName", exposedPatchName_);
694
696 << "Restricting to cellZone(s) " << flatOutput(zoneNames_)
697 << " with exposed internal faces into patch "
698 << mesh.boundaryMesh().findPatchID(exposedPatchName_) << endl;
699 }
700}
701
702
704(
705 const word& name,
706 const polyMesh& mesh,
707 const dictionary& dict
708)
709:
711{}
712
713
714// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
715
717{}
718
719
720// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
721
723{
724 const fvMesh& fvm = static_cast<const fvMesh&>(mesh());
725
726 return fvm.time().timeIndex() != prevTimeIndex_;
727}
728
729
731{
732 surface_.clear();
733 meshCells_.clear();
734 isoSurfacePtr_.reset(nullptr);
735 subMeshPtr_.reset(nullptr);
736
737 // Clear derived data
739
740 // Already marked as expired
741 if (prevTimeIndex_ == -1)
742 {
743 return false;
744 }
745
746 // Force update
747 prevTimeIndex_ = -1;
748 return true;
749}
750
751
753{
754 return updateGeometry();
755}
756
757
760(
761 const interpolation<scalar>& sampler
762) const
763{
764 return sampleOnFaces(sampler);
765}
766
767
770(
771 const interpolation<vector>& sampler
772) const
773{
774 return sampleOnFaces(sampler);
775}
776
777
780(
781 const interpolation<sphericalTensor>& sampler
782) const
783{
784 return sampleOnFaces(sampler);
785}
786
787
790(
791 const interpolation<symmTensor>& sampler
792) const
793{
794 return sampleOnFaces(sampler);
795}
796
797
800(
801 const interpolation<tensor>& sampler
802) const
803{
804 return sampleOnFaces(sampler);
805}
806
807
810(
811 const interpolation<scalar>& interpolator
812) const
813{
814 return sampleOnPoints(interpolator);
815}
816
817
820(
821 const interpolation<vector>& interpolator
822) const
823{
824 return sampleOnPoints(interpolator);
825}
826
829(
830 const interpolation<sphericalTensor>& interpolator
831) const
832{
833 return sampleOnPoints(interpolator);
834}
835
836
839(
840 const interpolation<symmTensor>& interpolator
841) const
842{
843 return sampleOnPoints(interpolator);
844}
845
846
849(
850 const interpolation<tensor>& interpolator
851) const
852{
853 return sampleOnPoints(interpolator);
854}
855
856
858{
859 os << "isoSurface: " << name() << " :";
860 isoParams_.print(os);
861 os << " field:" << isoField_
862 << " value:" << flatOutput(isoValues_);
863}
864
865
866// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addNamedToRunTimeSelectionTable(baseType, thisType, argNames, lookupName)
Add to construction table with 'lookupName' as the key.
Minimal example by using system/controlDict.functions:
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
label timeIndex() const noexcept
Return current time index.
Definition: TimeStateI.H:37
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
T & first()
Return the first element of the list.
Definition: UListI.H:202
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
bitSet selection(const labelUList &zoneIds) const
Definition: ZoneMesh.C:605
label findIndex(const wordRe &key) const
Zone index for the first match, return -1 if not found.
Definition: ZoneMesh.C:497
unsigned int count(const bool on=true) const
Count number of bits set.
Definition: bitSetI.H:500
static word defaultName
The default cloud name: defaultCloud.
Definition: cloud.H:90
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, bool mandatory=true) const
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Abstract base class for volume field interpolation.
Definition: interpolation.H:60
Preferences for controlling iso-surface algorithms.
filterType filter() const noexcept
Get current filter type.
@ ALGO_DEFAULT
Use current 'standard' algorithm.
algorithmType algorithm() const noexcept
Get current algorithm.
scalar print()
Print to screen.
label findPatchID(const word &patchName, const bool allowNotFound=true) const
Find patch index given a name, return -1 if not found.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:504
bool upToDate(const regIOobject &) const
Return true if up-to-date with respect to given object.
Definition: regIOobject.C:339
A sampledSurface defined by a surface of iso value. It only recalculates the iso-surface if time chan...
virtual ~sampledIsoSurface()
Destructor.
virtual bool expire()
Mark the surface as needing an update.
virtual bool needsUpdate() const
Does the surface need an update?
virtual bool update()
Update the surface as required.
An abstract class for surfaces with sampling.
static tmp< VolumeField< Type > > pointAverage(const PointField< Type > &pfld)
Create cell values by averaging the point values.
virtual void clearGeom() const
Additional cleanup when clearing the geometry.
const polyMesh & mesh() const noexcept
Access to the underlying mesh.
bool interpolate() const noexcept
Same as isPointData()
A class for managing temporary objects.
Definition: tmp.H:65
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
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
const pointField & points
label nPoints
#define DebugInfo
Report an information message using Foam::Info.
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
#define DebugInFunction
Report an information message using Foam::Info.
List< bool > select(const label n, const labelUList &locations)
Definition: BitOps.C:142
const wordList surface
Standard surface field types (scalar, vector, tensor, etc)
Namespace for OpenFOAM.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
GeometricField< scalar, pointPatchField, pointMesh > pointScalarField
List< label > labelList
A List of labels.
Definition: List.H:66
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
MeshedSurface< face > meshedSurface
void inplaceUniqueSort(ListType &input)
Inplace sorting and removal of duplicates.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:215
IOerror FatalIOError
error FatalError
Type gMin(const FieldField< Field, Type > &f)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
List< surfZone > surfZoneList
Definition: surfZoneList.H:47
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Type gMax(const FieldField< Field, Type > &f)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333