sampledSurfaces.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-2024 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "sampledSurfaces.H"
30 #include "polySurface.H"
31 
32 #include "mapPolyMesh.H"
33 #include "volFields.H"
34 #include "surfaceFields.H"
35 #include "HashOps.H"
36 #include "ListOps.H"
37 #include "Time.H"
38 #include "IndirectList.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45  defineTypeNameAndDebug(sampledSurfaces, 0);
46 
48  (
49  functionObject,
50  sampledSurfaces,
51  dictionary
52  );
53 }
54 
55 Foam::scalar Foam::sampledSurfaces::mergeTol_ = 1e-10;
56 
57 
58 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
59 
60 void Foam::sampledSurfaces::storeRegistrySurface
61 (
62  const sampledSurface& s
63 )
64 {
65  s.storeRegistrySurface
66  (
67  storedObjects(),
68  IOobject::groupName(name(), s.name()) // surfaceName
69  );
70 }
71 
72 
73 Foam::IOobjectList Foam::sampledSurfaces::preCheckFields()
74 {
75  wordList allFields; // Just needed for warnings
76  HashTable<wordHashSet> selected;
77 
78  IOobjectList objects;
79 
80  if (loadFromFiles_)
81  {
82  // Check files for a particular time
83  objects = IOobjectList(obr_, obr_.time().timeName());
84 
85  allFields = objects.names();
86  selected = objects.classes(fieldSelection_);
87  }
88  else
89  {
90  // Check currently available fields
91  allFields = obr_.names();
92  selected = obr_.classes(fieldSelection_);
93  }
94 
95  // Parallel consistency (no-op in serial)
96  Pstream::mapCombineReduce(selected, HashSetOps::plusEqOp<word>());
97 
98 
99  DynamicList<label> missed(fieldSelection_.size());
100 
101  // Detect missing fields
102  forAll(fieldSelection_, i)
103  {
104  if (!ListOps::found(allFields, fieldSelection_[i]))
105  {
106  missed.push_back(i);
107  }
108  }
109 
110  if (missed.size())
111  {
113  << nl
114  << "Cannot find "
115  << (loadFromFiles_ ? "field file" : "registered field")
116  << " matching "
117  << UIndirectList<wordRe>(fieldSelection_, missed) << endl;
118  }
119 
120 
121  // Currently only support volume and surface field types
122  label nVolumeFields = 0;
123  label nSurfaceFields = 0;
124 
125  forAllConstIters(selected, iter)
126  {
127  const word& clsName = iter.key();
128  const label n = iter.val().size();
129 
130  if (fieldTypes::volume.contains(clsName))
131  {
132  nVolumeFields += n;
133  }
134  else if (fieldTypes::surface.contains(clsName))
135  {
136  nSurfaceFields += n;
137  }
138  }
139 
140  // Now propagate field counts (per surface)
141  // - can update writer even when not writing without problem
142 
143  forAll(*this, surfi)
144  {
145  const sampledSurface& s = (*this)[surfi];
146  surfaceWriter& outWriter = writers_[surfi];
147 
148  outWriter.nFields
149  (
150  nVolumeFields
151  + (s.withSurfaceFields() ? nSurfaceFields : 0)
152  +
153  (
154  // Face-id field, but not if the writer manages that itself
155  !s.isPointData() && s.hasFaceIds() && !outWriter.usesFaceIds()
156  ? 1 : 0
157  )
158  );
159  }
160 
161  return objects;
162 }
163 
164 
165 Foam::autoPtr<Foam::surfaceWriter> Foam::sampledSurfaces::newWriter
166 (
167  word writerType,
168  const dictionary& topDict,
169  const dictionary& surfDict
170 )
171 {
172  // Per-surface adjustment
173  surfDict.readIfPresent<word>("surfaceFormat", writerType);
174 
175  return surfaceWriter::New
176  (
177  writerType,
178  // Top-level/surface-specific "formatOptions"
179  surfaceWriter::formatOptions(topDict, surfDict, writerType)
180  );
181 }
182 
183 
184 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
185 
186 Foam::sampledSurfaces::sampledSurfaces
187 (
188  const word& name,
189  const Time& runTime,
190  const dictionary& dict
191 )
192 :
193  functionObjects::fvMeshFunctionObject(name, runTime, dict),
195  loadFromFiles_(false),
196  verbose_(false),
197  onExecute_(false),
198  outputPath_
199  (
200  time_.globalPath()/functionObject::outputPrefix/name
201  )
202 {
203  outputPath_.clean(); // Remove unneeded ".."
204 
205  read(dict);
206 }
207 
208 
209 Foam::sampledSurfaces::sampledSurfaces
210 (
211  const word& name,
212  const objectRegistry& obr,
213  const dictionary& dict,
214  const bool loadFromFiles
215 )
216 :
217  functionObjects::fvMeshFunctionObject(name, obr, dict),
219  loadFromFiles_(loadFromFiles),
220  verbose_(false),
221  onExecute_(false),
222  outputPath_
223  (
224  time_.globalPath()/functionObject::outputPrefix/name
225  )
226 {
227  outputPath_.clean(); // Remove unneeded ".."
229  read(dict);
230 }
231 
232 
233 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
234 
236 {
237  bool old(verbose_);
238  verbose_ = on;
239  return old;
240 }
241 
242 
244 {
246 
248  writers_.clear();
249  actions_.clear();
250  hasFaces_.clear();
251  fieldSelection_.clear();
252 
253  verbose_ = dict.getOrDefault("verbose", false);
254  onExecute_ = dict.getOrDefault("sampleOnExecute", false);
255 
256  sampleFaceScheme_ =
257  dict.getOrDefault<word>("sampleScheme", "cell");
258 
259  sampleNodeScheme_ =
260  dict.getOrDefault<word>("interpolationScheme", "cellPoint");
261 
262  const entry* eptr = dict.findEntry("surfaces", keyType::LITERAL);
263 
264  // Surface writer type and format options
265  const word writerType =
266  (eptr ? dict.get<word>("surfaceFormat") : word::null);
267 
268  // Store on registry?
269  const bool dfltStore = dict.getOrDefault("store", false);
270 
271  if (eptr && eptr->isDict())
272  {
273  PtrList<sampledSurface> surfs(eptr->dict().size());
274 
275  actions_.resize(surfs.size(), ACTION_WRITE); // Default action
276  writers_.resize(surfs.size());
277 
278  label surfi = 0;
279 
280  for (const entry& dEntry : eptr->dict())
281  {
282  if (!dEntry.isDict())
283  {
284  continue;
285  }
286 
287  const dictionary& surfDict = dEntry.dict();
288 
289  autoPtr<sampledSurface> surf =
291  (
292  dEntry.keyword(),
293  mesh_,
294  surfDict
295  );
296 
297  if (!surf || !surf->enabled())
298  {
299  continue;
300  }
301 
302  // Define the surface
303  surfs.set(surfi, surf);
304 
305  // Define additional action(s)
306  if (surfDict.getOrDefault("store", dfltStore))
307  {
308  actions_[surfi] |= ACTION_STORE;
309  }
310 
311  // Define surface writer, but do NOT yet attach a surface
312  writers_.set
313  (
314  surfi,
315  newWriter(writerType, dict, surfDict)
316  );
317 
318  writers_[surfi].isPointData(surfs[surfi].isPointData());
319 
320  // Use outputDir/TIME/surface-name
321  writers_[surfi].useTimeDir(true);
322  writers_[surfi].verbose(verbose_);
323 
324  ++surfi;
325  }
326 
327  surfs.resize(surfi);
328  actions_.resize(surfi);
329  writers_.resize(surfi);
330  surfaces().transfer(surfs);
331  }
332  else if (eptr)
333  {
334  // This is slightly trickier.
335  // We want access to the individual dictionaries used for construction
336 
337  DynamicList<dictionary> capture;
338 
339  PtrList<sampledSurface> input
340  (
341  eptr->stream(),
342  sampledSurface::iNewCapture(mesh_, capture)
343  );
344 
345  PtrList<sampledSurface> surfs(input.size());
346 
347  actions_.resize(surfs.size(), ACTION_WRITE); // Default action
348  writers_.resize(surfs.size());
349 
350  label surfi = 0;
351 
352  forAll(input, inputi)
353  {
354  const dictionary& surfDict = capture[inputi];
355 
356  autoPtr<sampledSurface> surf = input.release(inputi);
357 
358  if (!surf || !surf->enabled())
359  {
360  continue;
361  }
362 
363  // Define the surface
364  surfs.set(surfi, surf);
365 
366  // Define additional action(s)
367  if (surfDict.getOrDefault("store", dfltStore))
368  {
369  actions_[surfi] |= ACTION_STORE;
370  }
371 
372  // Define surface writer, but do NOT yet attach a surface
373  writers_.set
374  (
375  surfi,
376  newWriter(writerType, dict, surfDict)
377  );
378 
379  writers_[surfi].isPointData(surfs[surfi].isPointData());
380 
381  // Use outputDir/TIME/surface-name
382  writers_[surfi].useTimeDir(true);
383  writers_[surfi].verbose(verbose_);
384 
385  ++surfi;
386  }
387 
388  surfs.resize(surfi);
389  actions_.resize(surfi);
390  writers_.resize(surfi);
391  surfaces().transfer(surfs);
392  }
393 
394 
395  const auto& surfs = surfaces();
396 
397  // Have some surfaces, so sort out which fields are needed and report
398 
399  hasFaces_.resize_nocopy(surfs.size());
400  hasFaces_ = false;
401 
402  if (surfs.size())
403  {
404  dict.readEntry("fields", fieldSelection_);
405  fieldSelection_.uniq();
406 
407  forAll(*this, surfi)
408  {
409  const sampledSurface& s = (*this)[surfi];
410 
411  if (!surfi)
412  {
413  Info<< "Sampled surface:" << nl;
414  }
415 
416  Info<< " " << s.name() << " -> " << writers_[surfi].type();
417  if (actions_[surfi] & ACTION_STORE)
418  {
419  Info<< ", store on registry ("
420  << IOobject::groupName(name(), s.name()) << ')';
421  }
422  Info<< nl;
423  Info<< " ";
424  s.print(Info, 0);
425  Info<< nl;
426  }
427  Info<< nl;
428  }
429 
430  if (debug && UPstream::master())
431  {
432  Pout<< "sample fields:" << fieldSelection_ << nl
433  << "sample surfaces:" << nl << '(' << nl;
434 
435  for (const sampledSurface& s : surfaces())
436  {
437  Pout<< " " << s << nl;
438  }
439  Pout<< ')' << endl;
440  }
441 
442  // Ensure all surfaces and merge information are expired
443  expire(true);
444 
445  return true;
446 }
447 
448 
449 bool Foam::sampledSurfaces::performAction(unsigned request)
450 {
451  // Update surfaces and store
452  bool ok = false;
453 
454  forAll(*this, surfi)
455  {
456  sampledSurface& s = (*this)[surfi];
457 
458  if (request & actions_[surfi])
459  {
460  if (s.update())
461  {
462  writers_[surfi].expire();
463  hasFaces_[surfi] = false;
464  }
465 
466  if (!hasFaces_[surfi])
467  {
468  hasFaces_[surfi] = returnReduceOr(s.faces().size());
469  }
470 
471  ok = ok || hasFaces_[surfi];
472 
473  // Store surfaces (even empty ones) otherwise we miss geometry
474  // updates.
475  // Any associated fields will be removed if the size changes
476 
477  if ((request & actions_[surfi]) & ACTION_STORE)
478  {
479  storeRegistrySurface(s);
480  }
481  }
482  }
483 
484  if (!ok)
485  {
486  // No surface with an applicable action or with faces to sample
487  return true;
488  }
489 
490 
491  // Determine availability of fields.
492  // Count per-surface number of fields, including Ids etc
493  // which only seems to be needed for VTK legacy
494 
495  IOobjectList objects = preCheckFields();
496 
497  // Update writers
498 
499  forAll(*this, surfi)
500  {
501  const sampledSurface& s = (*this)[surfi];
502 
503  if (((request & actions_[surfi]) & ACTION_WRITE) && hasFaces_[surfi])
504  {
505  surfaceWriter& outWriter = writers_[surfi];
506 
507  if (outWriter.needsUpdate())
508  {
509  outWriter.setSurface(s);
510  }
511 
512  outWriter.open(outputPath_/s.name());
513 
514  outWriter.beginTime(obr_.time());
515 
516  // Face-id field, but not if the writer manages that itself
517  if (!s.isPointData() && s.hasFaceIds() && !outWriter.usesFaceIds())
518  {
519  // This is still quite horrible.
520 
521  Field<label> ids(s.faceIds());
522 
523  if
524  (
526  (
527  !ListOps::found(ids, lessOp1<label>(0))
528  )
529  )
530  {
531  // From 0-based to 1-based, provided there are no
532  // negative ids (eg, encoded solid/side)
533 
534  ids += label(1);
535  }
536 
537  writeSurface(outWriter, ids, "Ids");
538  }
539  }
540  }
541 
542  // Sample fields
543 
544  performAction<volScalarField>(objects, request);
545  performAction<volVectorField>(objects, request);
546  performAction<volSphericalTensorField>(objects, request);
547  performAction<volSymmTensorField>(objects, request);
548  performAction<volTensorField>(objects, request);
549 
550  // Only bother with surface fields if a sampler supports them
551  if
552  (
553  testAny
554  (
555  surfaces(),
556  [](const sampledSurface& s) { return s.withSurfaceFields(); }
557  )
558  )
559  {
560  performAction<surfaceScalarField>(objects, request);
561  performAction<surfaceVectorField>(objects, request);
562  performAction<surfaceSphericalTensorField>(objects, request);
563  performAction<surfaceSymmTensorField>(objects, request);
564  performAction<surfaceTensorField>(objects, request);
565  }
566 
567 
568  // Finish this time step
569  forAll(writers_, surfi)
570  {
571  if (((request & actions_[surfi]) & ACTION_WRITE) && hasFaces_[surfi])
572  {
573  // Write geometry if no fields were written so that we still
574  // can have something to look at
575 
576  if (!writers_[surfi].wroteData())
577  {
578  writers_[surfi].write();
579  }
580 
581  writers_[surfi].endTime();
582  }
583  }
584 
585  return true;
586 }
587 
588 
590 {
591  if (onExecute_)
592  {
593  return performAction(ACTION_ALL & ~ACTION_WRITE);
594  }
595 
596  return true;
597 }
598 
601 {
602  return performAction(ACTION_ALL);
603 }
604 
605 
607 {
608  if (&mpm.mesh() == &mesh_)
609  {
610  expire();
611  }
612 
613  // pointMesh and interpolation will have been reset in mesh.update
614 }
615 
616 
617 void Foam::sampledSurfaces::movePoints(const polyMesh& mesh)
618 {
619  if (&mesh == &mesh_)
620  {
621  expire();
622  }
623 }
624 
625 
627 {
628  if (state != polyMesh::UNCHANGED)
629  {
630  // May want to use force expiration here
631  expire();
632  }
633 }
634 
635 
636 bool Foam::sampledSurfaces::needsUpdate() const
637 {
638  for (const sampledSurface& s : surfaces())
639  {
640  if (s.needsUpdate())
641  {
642  return true;
643  }
644  }
645 
646  return false;
647 }
648 
649 
650 bool Foam::sampledSurfaces::expire(const bool force)
651 {
652  // Dimension as fraction of mesh bounding box
653  const scalar mergeDim = mergeTol_ * mesh_.bounds().mag();
654 
655  bool changed = false;
656 
657  forAll(*this, surfi)
658  {
659  sampledSurface& s = (*this)[surfi];
660 
661  if (s.invariant() && !force)
662  {
663  // 'Invariant' - does not change when geometry does
664  continue;
665  }
666  if (s.expire())
667  {
668  changed = true;
669  }
670 
671  writers_[surfi].expire();
672  writers_[surfi].mergeDim(mergeDim);
673  hasFaces_[surfi] = false;
674  }
675 
676  // True if any surfaces just expired
677  return changed;
678 }
679 
680 
681 bool Foam::sampledSurfaces::update()
682 {
683  if (!needsUpdate())
684  {
685  return false;
686  }
687 
688  bool changed = false;
689 
690  forAll(*this, surfi)
691  {
692  sampledSurface& s = (*this)[surfi];
693 
694  if (s.update())
695  {
696  changed = true;
697  writers_[surfi].expire();
698  hasFaces_[surfi] = returnReduceOr(s.faces().size());
699  }
700  }
701 
702  return changed;
703 }
704 
707 {
708  return mergeTol_;
709 }
710 
711 
712 Foam::scalar Foam::sampledSurfaces::mergeTol(scalar tol) noexcept
713 {
714  scalar old(mergeTol_);
715  mergeTol_ = tol;
716  return old;
717 }
718 
719 
720 // ************************************************************************* //
Foam::surfaceFields.
virtual bool write()
Sample and write.
dictionary dict
const polyMesh & mesh() const noexcept
Return polyMesh.
Definition: mapPolyMesh.H:440
objectRegistry & storedObjects()
Write access to the output objects ("functionObjectObjects") registered on Time.
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
An abstract class for surfaces with sampling.
virtual void movePoints(const polyMesh &mesh)
Update for mesh point-motion - expires the surfaces.
bool found(const ListType &input, const UnaryPredicate &pred, const label start=0)
Same as found_if.
Definition: ListOps.H:858
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
virtual const dictionary & dict() const =0
Return dictionary, if entry is a dictionary, otherwise Fatal.
virtual void readUpdate(const polyMesh::readUpdateState state)
Update for changes of mesh due to readUpdate - expires the surfaces.
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
Abstract base-class for Time/database function objects.
static dictionary formatOptions(const dictionary &dict, const word &formatName, const word &entryName="formatOptions")
Same as fileFormats::getFormatOptions.
Definition: surfaceWriter.C:57
static autoPtr< sampledSurface > New(const word &name, const polyMesh &mesh, const dictionary &dict)
Return a reference to the selected surface.
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Macros for easy insertion into run-time selection tables.
const word & name() const noexcept
Return the name of this functionObject.
static scalar mergeTol() noexcept
Get merge tolerance.
Various functions to operate on Lists.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
bool returnReduceAnd(const bool value, const label comm=UPstream::worldComm)
Perform logical (and) MPI Allreduce on a copy. Uses UPstream::reduceAnd.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
virtual bool read(const dictionary &dict)
Read the sampledSurfaces dictionary.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
static Istream & input(Istream &is, IntRange< T > &range)
Definition: IntRanges.C:33
static const word null
An empty word.
Definition: word.H:84
String literal.
Definition: keyType.H:82
const wordList surface
Standard surface field types (scalar, vector, tensor, etc)
const direction noexcept
Definition: Scalar.H:258
const wordList volume
Standard volume field types (scalar, vector, tensor, etc)
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
int debug
Static debugging option.
defineTypeNameAndDebug(combustionModel, 0)
static bool clean(std::string &str)
Cleanup filename string, possibly applies other transformations such as changing the path separator e...
Definition: fileName.C:192
virtual bool execute()
Sample and store if the sampleOnExecute is enabled.
void resize(const label newLen)
Adjust size of PtrList.
Definition: PtrList.C:95
bool verbose(bool on) noexcept
Enable/disable verbose output.
List< word > wordList
List of word.
Definition: fileName.H:59
#define WarningInFunction
Report a warning using Foam::Warning.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
static bool master(const label communicator=worldComm)
True if process corresponds to the master rank in the communicator.
Definition: UPstream.H:1094
const entry * findEntry(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:84
void clear()
Clear the PtrList. Delete allocated entries and set size to zero.
Definition: PtrListI.H:81
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
Registry of regIOobjects.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
static void mapCombineReduce(Container &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) applying cop to inplace combine map values from different processo...
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
static autoPtr< surfaceWriter > New(const word &writeType)
Select construct a surfaceWriter.
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
A keyword and a list of tokens is an &#39;entry&#39;.
Definition: entry.H:63
addToRunTimeSelectionTable(functionObject, pointHistory, dictionary)
virtual void updateMesh(const mapPolyMesh &mpm)
Update for changes of mesh - expires the surfaces.