fvMeshAdderTemplates.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-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2019 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 "volFields.H"
30 #include "surfaceFields.H"
31 #include "emptyFvPatchField.H"
33 
34 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35 
36 template<class Type>
37 void Foam::fvMeshAdder::MapVolField
38 (
39  const mapAddedPolyMesh& meshMap,
40 
41  GeometricField<Type, fvPatchField, volMesh>& fld,
42  const GeometricField<Type, fvPatchField, volMesh>& fldToAdd,
43  const bool fullyMapped
44 )
45 {
46  const fvMesh& mesh = fld.mesh();
47 
48  // Internal field
49  // ~~~~~~~~~~~~~~
50 
51  {
52  // Store old internal field
53  Field<Type> oldInternalField(fld.primitiveField());
54 
55  // Modify internal field
56  Field<Type>& intFld = fld.primitiveFieldRef();
57 
58  intFld.setSize(mesh.nCells());
59 
60  intFld.rmap(oldInternalField, meshMap.oldCellMap());
61  intFld.rmap(fldToAdd.primitiveField(), meshMap.addedCellMap());
62  }
63 
64 
65  // Patch fields from old mesh
66  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
67 
68  auto& bfld = fld.boundaryFieldRef();
69 
70  {
71  const labelList& oldPatchMap = meshMap.oldPatchMap();
72  const labelList& oldPatchStarts = meshMap.oldPatchStarts();
73  const labelList& oldPatchSizes = meshMap.oldPatchSizes();
74 
75  // Reorder old patches in order of new ones. Put removed patches at end.
76 
77  label unusedPatchi = 0;
78 
79  forAll(oldPatchMap, patchi)
80  {
81  label newPatchi = oldPatchMap[patchi];
82 
83  if (newPatchi != -1)
84  {
85  unusedPatchi++;
86  }
87  }
88 
89  label nUsedPatches = unusedPatchi;
90 
91  // Reorder list for patchFields
92  labelList oldToNew(oldPatchMap.size());
93 
94  forAll(oldPatchMap, patchi)
95  {
96  const label newPatchi = oldPatchMap[patchi];
97 
98  if (newPatchi != -1)
99  {
100  oldToNew[patchi] = newPatchi;
101  }
102  else
103  {
104  oldToNew[patchi] = unusedPatchi++;
105  }
106  }
107 
108 
109  // Sort deleted ones last so is now in newPatch ordering
110  bfld.reorder(oldToNew);
111  // Extend to covers all patches
112  bfld.setSize(mesh.boundaryMesh().size());
113  // Delete unused patches
114  for
115  (
116  label newPatchi = nUsedPatches;
117  newPatchi < bfld.size();
118  newPatchi++
119  )
120  {
121  bfld.set(newPatchi, nullptr);
122  }
123 
124 
125  // Map old values
126  // ~~~~~~~~~~~~~~
127 
128  forAll(oldPatchMap, patchi)
129  {
130  label newPatchi = oldPatchMap[patchi];
131 
132  if (newPatchi != -1)
133  {
134  labelList newToOld
135  (
136  calcPatchMap
137  (
138  oldPatchStarts[patchi],
139  oldPatchSizes[patchi],
140  meshMap.oldFaceMap(),
141  mesh.boundaryMesh()[newPatchi],
142  -1 // unmapped value
143  )
144  );
145 
146  directFvPatchFieldMapper patchMapper(newToOld);
147 
148  // Override mapping (for use in e.g. fvMeshDistribute where
149  // it sorts mapping out itself)
150  if (fullyMapped)
151  {
152  patchMapper.hasUnmapped() = false;
153  }
154 
155  // Create new patchField with same type as existing one.
156  // Note:
157  // - boundaryField already in new order so access with newPatchi
158  // - fld.boundaryField()[newPatchi] both used for type and old
159  // value
160  // - hope that field mapping allows aliasing since old and new
161  // are same memory!
162  bfld.set
163  (
164  newPatchi,
166  (
167  bfld[newPatchi], // old field
168  mesh.boundary()[newPatchi], // new fvPatch
169  fld(), // new internal field
170  patchMapper // mapper (new to old)
171  )
172  );
173  }
174  }
175  }
176 
177 
178 
179  // Patch fields from added mesh
180  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
181 
182  {
183  const labelList& addedPatchMap = meshMap.addedPatchMap();
184 
185  // Add addedMesh patches
186  forAll(addedPatchMap, patchi)
187  {
188  label newPatchi = addedPatchMap[patchi];
189 
190  if (newPatchi != -1)
191  {
192  const polyPatch& newPatch = mesh.boundaryMesh()[newPatchi];
193  const polyPatch& oldPatch =
194  fldToAdd.mesh().boundaryMesh()[patchi];
195 
196  if (!bfld(newPatchi))
197  {
198  // First occurrence of newPatchi. Map from existing
199  // patchField
200 
201  // From new patch faces to patch faces on added mesh.
202  labelList newToAdded
203  (
204  calcPatchMap
205  (
206  oldPatch.start(),
207  oldPatch.size(),
208  meshMap.addedFaceMap(),
209  newPatch,
210  -1 // unmapped values
211  )
212  );
213 
214  directFvPatchFieldMapper patchMapper(newToAdded);
215 
216  // Override mapping (for use in e.g. fvMeshDistribute where
217  // it sorts mapping out itself)
218  if (fullyMapped)
219  {
220  patchMapper.hasUnmapped() = false;
221  }
222 
223  bfld.set
224  (
225  newPatchi,
227  (
228  fldToAdd.boundaryField()[patchi], // added field
229  mesh.boundary()[newPatchi], // new fvPatch
230  fld(), // new int. field
231  patchMapper // mapper
232  )
233  );
234  }
235  else
236  {
237  // PatchField will have correct size already. Just slot in
238  // my elements.
239 
240  labelList addedToNew(oldPatch.size(), -1);
241  forAll(addedToNew, i)
242  {
243  label addedFacei = oldPatch.start()+i;
244  label newFacei = meshMap.addedFaceMap()[addedFacei];
245  label patchFacei = newFacei-newPatch.start();
246  if (patchFacei >= 0 && patchFacei < newPatch.size())
247  {
248  addedToNew[i] = patchFacei;
249  }
250  }
251 
252  bfld[newPatchi].rmap
253  (
254  fldToAdd.boundaryField()[patchi],
255  addedToNew
256  );
257  }
258  }
259  }
260  }
261 }
262 
263 
264 template<class Type>
266 (
267  const mapAddedPolyMesh& meshMap,
268  const fvMesh& mesh,
269  const fvMesh& meshToAdd,
270  const bool fullyMapped
271 )
272 {
274 
276  (
277  mesh.objectRegistry::lookupClass<fldType>()
278  );
279 
280  HashTable<const fldType*> fieldsToAdd
281  (
282  meshToAdd.objectRegistry::lookupClass<fldType>()
283  );
284 
285  // It is necessary to enforce that all old-time fields are stored
286  // before the mapping is performed. Otherwise, if the
287  // old-time-level field is mapped before the field itself, sizes
288  // will not match.
289 
290  forAllIters(fields, fieldIter)
291  {
292  fldType& fld = const_cast<fldType&>(*fieldIter());
293 
294  DebugPout
295  << "MapVolFields : Storing old time for " << fld.name()
296  << endl;
297 
298  fld.storeOldTimes();
299  }
300 
301 
302  forAllIters(fields, fieldIter)
303  {
304  fldType& fld = const_cast<fldType&>(*fieldIter());
305 
306  if (fieldsToAdd.found(fld.name()))
307  {
308  const fldType& fldToAdd = *fieldsToAdd[fld.name()];
309 
310  DebugPout
311  << "MapVolFields : mapping " << fld.name()
312  << " and " << fldToAdd.name() << endl;
313 
314  MapVolField<Type>(meshMap, fld, fldToAdd, fullyMapped);
315  }
316  else
317  {
319  << "Not mapping field " << fld.name()
320  << " since not present on mesh to add"
321  << endl;
322  }
323  }
324 }
325 
326 
327 template<class Type>
328 void Foam::fvMeshAdder::MapSurfaceField
329 (
330  const mapAddedPolyMesh& meshMap,
331 
334  const bool fullyMapped
335 )
336 {
337  const fvMesh& mesh = fld.mesh();
338  const labelList& oldPatchStarts = meshMap.oldPatchStarts();
339 
340  auto& bfld = fld.boundaryFieldRef();
341 
342  // Internal field
343  // ~~~~~~~~~~~~~~
344 
345  // Store old internal field
346  {
347  Field<Type> oldField(fld);
348 
349  // Modify internal field
350  Field<Type>& intFld = fld.primitiveFieldRef();
351 
352  intFld.setSize(mesh.nInternalFaces());
353 
354  intFld.rmap(oldField, meshMap.oldFaceMap());
355  intFld.rmap(fldToAdd, meshMap.addedFaceMap());
356 
357 
358  // Faces that were boundary faces but are not anymore.
359  // Use owner value (so lowest numbered cell, i.e. from 'old' not 'added'
360  // mesh)
361  forAll(bfld, patchi)
362  {
363  const fvsPatchField<Type>& pf = bfld[patchi];
364 
365  label start = oldPatchStarts[patchi];
366 
367  forAll(pf, i)
368  {
369  label newFacei = meshMap.oldFaceMap()[start + i];
370 
371  if (newFacei >= 0 && newFacei < mesh.nInternalFaces())
372  {
373  intFld[newFacei] = pf[i];
374  }
375  }
376  }
377  }
378 
379 
380  // Patch fields from old mesh
381  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
382 
383  {
384  const labelList& oldPatchMap = meshMap.oldPatchMap();
385  const labelList& oldPatchSizes = meshMap.oldPatchSizes();
386 
387  // Reorder old patches in order of new ones. Put removed patches at end.
388 
389  label unusedPatchi = 0;
390 
391  forAll(oldPatchMap, patchi)
392  {
393  label newPatchi = oldPatchMap[patchi];
394 
395  if (newPatchi != -1)
396  {
397  unusedPatchi++;
398  }
399  }
400 
401  label nUsedPatches = unusedPatchi;
402 
403  // Reorder list for patchFields
404  labelList oldToNew(oldPatchMap.size());
405 
406  forAll(oldPatchMap, patchi)
407  {
408  const label newPatchi = oldPatchMap[patchi];
409 
410  if (newPatchi != -1)
411  {
412  oldToNew[patchi] = newPatchi;
413  }
414  else
415  {
416  oldToNew[patchi] = unusedPatchi++;
417  }
418  }
419 
420 
421  // Sort deleted ones last so is now in newPatch ordering
422  bfld.reorder(oldToNew);
423  // Extend to covers all patches
424  bfld.setSize(mesh.boundaryMesh().size());
425  // Delete unused patches
426  for
427  (
428  label newPatchi = nUsedPatches;
429  newPatchi < bfld.size();
430  newPatchi++
431  )
432  {
433  bfld.set(newPatchi, nullptr);
434  }
435 
436 
437  // Map old values
438  // ~~~~~~~~~~~~~~
439 
440  forAll(oldPatchMap, patchi)
441  {
442  label newPatchi = oldPatchMap[patchi];
443 
444  if (newPatchi != -1)
445  {
446  labelList newToOld
447  (
448  calcPatchMap
449  (
450  oldPatchStarts[patchi],
451  oldPatchSizes[patchi],
452  meshMap.oldFaceMap(),
453  mesh.boundaryMesh()[newPatchi],
454  -1 // unmapped value
455  )
456  );
457 
458  directFvPatchFieldMapper patchMapper(newToOld);
459 
460  // Override mapping (for use in e.g. fvMeshDistribute where
461  // it sorts mapping out itself)
462  if (fullyMapped)
463  {
464  patchMapper.hasUnmapped() = false;
465  }
466 
467  // Create new patchField with same type as existing one.
468  // Note:
469  // - boundaryField already in new order so access with newPatchi
470  // - bfld[newPatchi] both used for type and old
471  // value
472  // - hope that field mapping allows aliasing since old and new
473  // are same memory!
474  bfld.set
475  (
476  newPatchi,
478  (
479  bfld[newPatchi], // old field
480  mesh.boundary()[newPatchi], // new fvPatch
481  fld(), // new internal field
482  patchMapper // mapper (new to old)
483  )
484  );
485  }
486  }
487  }
488 
489 
490 
491  // Patch fields from added mesh
492  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
493 
494  {
495  const labelList& addedPatchMap = meshMap.addedPatchMap();
496 
497  // Add addedMesh patches
498  forAll(addedPatchMap, patchi)
499  {
500  label newPatchi = addedPatchMap[patchi];
501 
502  if (newPatchi != -1)
503  {
504  const polyPatch& newPatch = mesh.boundaryMesh()[newPatchi];
505  const polyPatch& oldPatch =
506  fldToAdd.mesh().boundaryMesh()[patchi];
507 
508  if (!bfld(newPatchi))
509  {
510  // First occurrence of newPatchi. Map from existing
511  // patchField
512 
513  // From new patch faces to patch faces on added mesh.
514  labelList newToAdded
515  (
516  calcPatchMap
517  (
518  oldPatch.start(),
519  oldPatch.size(),
520  meshMap.addedFaceMap(),
521  newPatch,
522  -1 // unmapped values
523  )
524  );
525 
526  directFvPatchFieldMapper patchMapper(newToAdded);
527 
528  // Override mapping (for use in e.g. fvMeshDistribute where
529  // it sorts mapping out itself)
530  if (fullyMapped)
531  {
532  patchMapper.hasUnmapped() = false;
533  }
534 
535  bfld.set
536  (
537  newPatchi,
539  (
540  fldToAdd.boundaryField()[patchi],// added field
541  mesh.boundary()[newPatchi], // new fvPatch
542  fld(), // new int. field
543  patchMapper // mapper
544  )
545  );
546  }
547  else
548  {
549  // PatchField will have correct size already. Just slot in
550  // my elements.
551 
552  labelList addedToNew(oldPatch.size(), -1);
553  forAll(addedToNew, i)
554  {
555  label addedFacei = oldPatch.start()+i;
556  label newFacei = meshMap.addedFaceMap()[addedFacei];
557  label patchFacei = newFacei-newPatch.start();
558  if (patchFacei >= 0 && patchFacei < newPatch.size())
559  {
560  addedToNew[i] = patchFacei;
561  }
562  }
563 
564  bfld[newPatchi].rmap
565  (
566  fldToAdd.boundaryField()[patchi],
567  addedToNew
568  );
569  }
570  }
571  }
572  }
573 }
574 
575 
576 template<class Type>
578 (
579  const mapAddedPolyMesh& meshMap,
580  const fvMesh& mesh,
581  const fvMesh& meshToAdd,
582  const bool fullyMapped
583 )
584 {
586 
588  (
589  mesh.objectRegistry::lookupClass<fldType>()
590  );
591 
592  HashTable<const fldType*> fieldsToAdd
593  (
594  meshToAdd.objectRegistry::lookupClass<fldType>()
595  );
596 
597  // It is necessary to enforce that all old-time fields are stored
598  // before the mapping is performed. Otherwise, if the
599  // old-time-level field is mapped before the field itself, sizes
600  // will not match.
601 
602  forAllIters(fields, fieldIter)
603  {
604  fldType& fld = const_cast<fldType&>(*fieldIter());
605 
606  DebugPout
607  << "MapSurfaceFields : Storing old time for "
608  << fld.name() << endl;
609 
610  fld.storeOldTimes();
611  }
612 
613 
614  forAllIters(fields, fieldIter)
615  {
616  fldType& fld = const_cast<fldType&>(*fieldIter());
617 
618  if (fieldsToAdd.found(fld.name()))
619  {
620  const fldType& fldToAdd = *fieldsToAdd[fld.name()];
621 
622  DebugPout
623  << "MapSurfaceFields : mapping " << fld.name()
624  << " and " << fldToAdd.name() << endl;
625 
626  MapSurfaceField<Type>(meshMap, fld, fldToAdd, fullyMapped);
627  }
628  else
629  {
631  << "Not mapping field " << fld.name()
632  << " since not present on mesh to add"
633  << endl;
634  }
635  }
636 }
637 
638 
639 template<class Type>
640 void Foam::fvMeshAdder::MapDimField
641 (
642  const mapAddedPolyMesh& meshMap,
643 
645  const DimensionedField<Type, volMesh>& fldToAdd
646 )
647 {
648  const fvMesh& mesh = fld.mesh();
649 
650  // Store old field
651  Field<Type> oldField(fld);
652 
653  fld.setSize(mesh.nCells());
654 
655  fld.rmap(oldField, meshMap.oldCellMap());
656  fld.rmap(fldToAdd, meshMap.addedCellMap());
657 }
658 
659 
660 template<class Type>
662 (
663  const mapAddedPolyMesh& meshMap,
664  const fvMesh& mesh,
665  const fvMesh& meshToAdd
666 )
667 {
668  typedef DimensionedField<Type, volMesh> fldType;
669 
670  // Note: use strict flag on lookupClass to avoid picking up
671  // volFields
673  (
674  mesh.objectRegistry::lookupClass<fldType>(true)
675  );
676 
677  HashTable<const fldType*> fieldsToAdd
678  (
679  meshToAdd.objectRegistry::lookupClass<fldType>(true)
680  );
681 
682  forAllIters(fields, fieldIter)
683  {
684  fldType& fld = const_cast<fldType&>(*fieldIter());
685 
686  if (fieldsToAdd.found(fld.name()))
687  {
688  const fldType& fldToAdd = *fieldsToAdd[fld.name()];
689 
690  DebugPout
691  << "MapDimFields : mapping " << fld.name()
692  << " and " << fldToAdd.name() << endl;
693 
694  MapDimField<Type>(meshMap, fld, fldToAdd);
695  }
696  else
697  {
699  << "Not mapping field " << fld.name()
700  << " since not present on mesh to add"
701  << endl;
702  }
703  }
704 }
705 
706 
707 // ************************************************************************* //
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
volFields.H
Foam::mapAddedPolyMesh::oldPatchStarts
const labelList & oldPatchStarts() const
Return list of the old patch start labels.
Definition: mapAddedPolyMesh.H:178
Foam::mapAddedPolyMesh::addedCellMap
const labelList & addedCellMap() const
Definition: mapAddedPolyMesh.H:216
Foam::fvMeshAdder::MapSurfaceFields
static void MapSurfaceFields(const mapAddedPolyMesh &, const fvMesh &mesh, const fvMesh &meshToAdd, const bool fullyMapped=false)
Map all surfaceFields of Type.
Definition: fvMeshAdderTemplates.C:578
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Number of internal faces.
Definition: primitiveMeshI.H:78
Foam::fvMeshAdder::MapVolFields
static void MapVolFields(const mapAddedPolyMesh &, const fvMesh &mesh, const fvMesh &meshToAdd, const bool fullyMapped=false)
Map all volFields of Type.
Definition: fvMeshAdderTemplates.C:266
Foam::fvsPatchField
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:68
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
surfaceFields.H
Foam::surfaceFields.
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::mapAddedPolyMesh::addedFaceMap
const labelList & addedFaceMap() const
Definition: mapAddedPolyMesh.H:212
Foam::primitiveMesh::nCells
label nCells() const
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::Field
Generic templated field type.
Definition: Field.H:63
directFvPatchFieldMapper.H
fld
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;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputLagrangian.H:23
forAllIters
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:223
Foam::Field::rmap
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:478
Foam::fvsPatchField::New
static tmp< fvsPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, surfaceMesh > &)
Return a pointer to a new patchField created on freestore given.
Definition: fvsPatchFieldNew.C:74
Foam::mapAddedPolyMesh::oldCellMap
const labelList & oldCellMap() const
Definition: mapAddedPolyMesh.H:159
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
DebugPout
#define DebugPout
Report an information message using Foam::Pout.
Definition: messageStream.H:370
Foam::mapAddedPolyMesh::oldFaceMap
const labelList & oldFaceMap() const
Definition: mapAddedPolyMesh.H:155
Foam::mapAddedPolyMesh::addedPatchMap
const labelList & addedPatchMap() const
From added mesh patch index to new patch index or -1 if.
Definition: mapAddedPolyMesh.H:223
Foam::HashTable
A HashTable similar to std::unordered_map.
Definition: HashTable.H:105
Foam::fvMesh::boundary
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:679
Foam::fvPatchField::New
static tmp< fvPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
Definition: fvPatchFieldNew.C:88
emptyFvPatchField.H
Foam::List< label >
Foam::mapAddedPolyMesh::oldPatchSizes
const labelList & oldPatchSizes() const
Return list of the old patch sizes.
Definition: mapAddedPolyMesh.H:172
Foam::fvMeshAdder::MapDimFields
static void MapDimFields(const mapAddedPolyMesh &, const fvMesh &mesh, const fvMesh &meshToAdd)
Map all DimensionedFields of Type.
Definition: fvMeshAdderTemplates.C:662
Foam::HashTable::found
bool found(const Key &key) const
Return true if hashed entry is found in table.
Definition: HashTableI.H:100
Foam::mapAddedPolyMesh
Class containing mesh-to-mesh mapping information after a mesh addition where we add a mesh ('added m...
Definition: mapAddedPolyMesh.H:58
Foam::GeometricField< Type, fvPatchField, volMesh >
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:303
Foam::mapAddedPolyMesh::oldPatchMap
const labelList & oldPatchMap() const
From old patch index to new patch index or -1 if patch.
Definition: mapAddedPolyMesh.H:166
fields
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
Foam::GeometricField::boundaryField
const Boundary & boundaryField() const
Return const-reference to the boundary field.
Definition: GeometricFieldI.H:62
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54