fvMeshAdderTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "volFields.H"
27 #include "surfaceFields.H"
28 #include "pointFields.H"
33 
34 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35 
36 template<class Type>
37 void Foam::fvMeshAdder::MapVolField
38 (
39  const mapAddedPolyMesh& meshMap,
40 
41  VolField<Type>& fld,
42  const VolField<Type>& fldToAdd
43 )
44 {
45  const fvMesh& mesh = fld.mesh();
46 
47  // Internal field
48  // ~~~~~~~~~~~~~~
49 
50  {
51  // Store old internal field
52  Field<Type> oldInternalField(fld.primitiveField());
53 
54  // Modify internal field
55  Field<Type>& intFld = fld.primitiveFieldRef();
56 
57  intFld.setSize(mesh.nCells());
58 
59  intFld.rmap(oldInternalField, meshMap.oldCellMap());
60  intFld.rmap(fldToAdd.primitiveField(), meshMap.addedCellMap());
61  }
62 
63 
64  // Patch fields from old mesh
65  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
66 
67  typename VolField<Type>::
68  Boundary& 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  const 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  const 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 
149  // Create new patchField with same type as existing one.
150  // Note:
151  // - boundaryField already in new order so access with newPatchi
152  // - fld.boundaryField()[newPatchi] both used for type and old
153  // value
154  // - hope that field mapping allows aliasing since old and new
155  // are same memory!
156  bfld.set
157  (
158  newPatchi,
160  (
161  bfld[newPatchi], // old field
162  mesh.boundary()[newPatchi], // new fvPatch
163  fld(), // new internal field
164  patchMapper // mapper (new to old)
165  )
166  );
167  }
168  }
169  }
170 
171 
172 
173  // Patch fields from added mesh
174  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
175 
176  {
177  const labelList& addedPatchMap = meshMap.addedPatchMap();
178 
179  // Add addedMesh patches
180  forAll(addedPatchMap, patchi)
181  {
182  const label newPatchi = addedPatchMap[patchi];
183 
184  if (newPatchi != -1)
185  {
186  const polyPatch& newPatch = mesh.boundaryMesh()[newPatchi];
187  const polyPatch& oldPatch =
188  fldToAdd.mesh().boundaryMesh()[patchi];
189 
190  if (!bfld(newPatchi))
191  {
192  // First occurrence of newPatchi. Map from existing
193  // patchField
194 
195  // From new patch faces to patch faces on added mesh.
196  const labelList newToAdded
197  (
198  calcPatchMap
199  (
200  oldPatch.start(),
201  oldPatch.size(),
202  meshMap.addedFaceMap(),
203  newPatch,
204  -1 // unmapped values
205  )
206  );
207 
208  directFvPatchFieldMapper patchMapper(newToAdded);
209 
210  bfld.set
211  (
212  newPatchi,
214  (
215  fldToAdd.boundaryField()[patchi], // added field
216  mesh.boundary()[newPatchi], // new fvPatch
217  fld(), // new int. field
218  patchMapper // mapper
219  )
220  );
221  }
222  else
223  {
224  // PatchField will have correct size already. Just slot in
225  // my elements.
226 
227  labelList addedToNew(oldPatch.size(), -1);
228  forAll(addedToNew, i)
229  {
230  const label addedFacei = oldPatch.start() + i;
231  const label newFacei =
232  meshMap.addedFaceMap()[addedFacei];
233  const label patchFacei = newFacei-newPatch.start();
234 
235  if (patchFacei >= 0 && patchFacei < newPatch.size())
236  {
237  addedToNew[i] = patchFacei;
238  }
239  }
240 
241  bfld[newPatchi].map
242  (
243  fldToAdd.boundaryField()[patchi],
244  reverseFvPatchFieldMapper(addedToNew)
245  );
246  }
247  }
248  }
249  }
250 }
251 
252 
253 template<class Type>
255 (
256  const mapAddedPolyMesh& meshMap,
257  const fvMesh& mesh,
258  const fvMesh& meshToAdd
259 )
260 {
262  (
263  mesh.objectRegistry::lookupClass<VolField<Type>>()
264  );
265 
266  HashTable<const VolField<Type>*> fieldsToAdd
267  (
268  meshToAdd.objectRegistry::lookupClass<VolField<Type>>()
269  );
270 
271  for
272  (
273  typename HashTable<const VolField<Type>*>::
274  iterator fieldIter = fields.begin();
275  fieldIter != fields.end();
276  ++fieldIter
277  )
278  {
279  if (fvMesh::geometryFields.found(fieldIter()->name())) continue;
280 
282  const_cast<VolField<Type>&>
283  (
284  *fieldIter()
285  );
286 
287  if (fieldsToAdd.found(fld.name()))
288  {
289  const VolField<Type>& fldToAdd =
290  *fieldsToAdd[fld.name()];
291 
292  if (debug)
293  {
294  Pout<< "MapVolFields : mapping " << fld.name()
295  << " and " << fldToAdd.name() << endl;
296  }
297 
298  MapVolField<Type>(meshMap, fld, fldToAdd);
299  }
300  else
301  {
303  << "Not mapping field " << fld.name()
304  << " since not present on mesh to add"
305  << endl;
306  }
307  }
308 }
309 
310 
311 template<class Type>
312 void Foam::fvMeshAdder::MapSurfaceField
313 (
314  const mapAddedPolyMesh& meshMap,
315 
317  const SurfaceField<Type>& fldToAdd
318 )
319 {
320  const fvMesh& mesh = fld.mesh();
321  const labelList& oldPatchStarts = meshMap.oldPatchStarts();
322 
323  typename SurfaceField<Type>::
324  Boundary& bfld = fld.boundaryFieldRef();
325 
326  // Internal field
327  // ~~~~~~~~~~~~~~
328 
329  // Store old internal field
330  {
331  Field<Type> oldField(fld);
332 
333  // Modify internal field
334  Field<Type>& intFld = fld.primitiveFieldRef();
335 
336  intFld.setSize(mesh.nInternalFaces());
337 
338  intFld.rmap(oldField, meshMap.oldFaceMap());
339  intFld.rmap(fldToAdd, meshMap.addedFaceMap());
340 
341 
342  // Faces that were boundary faces but are not anymore.
343  // Use owner value (so lowest numbered cell, i.e. from 'old' not 'added'
344  // mesh)
345  forAll(bfld, patchi)
346  {
347  const fvsPatchField<Type>& pf = bfld[patchi];
348 
349  const label start = oldPatchStarts[patchi];
350 
351  forAll(pf, i)
352  {
353  const label newFacei = meshMap.oldFaceMap()[start + i];
354 
355  if (newFacei >= 0 && newFacei < mesh.nInternalFaces())
356  {
357  intFld[newFacei] = pf[i];
358  }
359  }
360  }
361  }
362 
363 
364  // Patch fields from old mesh
365  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
366 
367  {
368  const labelList& oldPatchMap = meshMap.oldPatchMap();
369  const labelList& oldPatchSizes = meshMap.oldPatchSizes();
370 
371  // Reorder old patches in order of new ones. Put removed patches at end.
372 
373  label unusedPatchi = 0;
374 
375  forAll(oldPatchMap, patchi)
376  {
377  const label newPatchi = oldPatchMap[patchi];
378 
379  if (newPatchi != -1)
380  {
381  unusedPatchi++;
382  }
383  }
384 
385  label nUsedPatches = unusedPatchi;
386 
387  // Reorder list for patchFields
388  labelList oldToNew(oldPatchMap.size());
389 
390  forAll(oldPatchMap, patchi)
391  {
392  const label newPatchi = oldPatchMap[patchi];
393 
394  if (newPatchi != -1)
395  {
396  oldToNew[patchi] = newPatchi;
397  }
398  else
399  {
400  oldToNew[patchi] = unusedPatchi++;
401  }
402  }
403 
404 
405  // Sort deleted ones last so is now in newPatch ordering
406  bfld.reorder(oldToNew);
407  // Extend to covers all patches
408  bfld.setSize(mesh.boundaryMesh().size());
409  // Delete unused patches
410  for
411  (
412  label newPatchi = nUsedPatches;
413  newPatchi < bfld.size();
414  newPatchi++
415  )
416  {
417  bfld.set(newPatchi, nullptr);
418  }
419 
420 
421  // Map old values
422  // ~~~~~~~~~~~~~~
423 
424  forAll(oldPatchMap, patchi)
425  {
426  const label newPatchi = oldPatchMap[patchi];
427 
428  if (newPatchi != -1)
429  {
430  const labelList newToOld
431  (
432  calcPatchMap
433  (
434  oldPatchStarts[patchi],
435  oldPatchSizes[patchi],
436  meshMap.oldFaceMap(),
437  mesh.boundaryMesh()[newPatchi],
438  -1 // unmapped value
439  )
440  );
441 
442  directFvPatchFieldMapper patchMapper(newToOld);
443 
444  // Create new patchField with same type as existing one.
445  // Note:
446  // - boundaryField already in new order so access with newPatchi
447  // - bfld[newPatchi] both used for type and old
448  // value
449  // - hope that field mapping allows aliasing since old and new
450  // are same memory!
451  bfld.set
452  (
453  newPatchi,
455  (
456  bfld[newPatchi], // old field
457  mesh.boundary()[newPatchi], // new fvPatch
458  fld(), // new internal field
459  patchMapper // mapper (new to old)
460  )
461  );
462  }
463  }
464  }
465 
466 
467 
468  // Patch fields from added mesh
469  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
470 
471  {
472  const labelList& addedPatchMap = meshMap.addedPatchMap();
473 
474  // Add addedMesh patches
475  forAll(addedPatchMap, patchi)
476  {
477  const label newPatchi = addedPatchMap[patchi];
478 
479  if (newPatchi != -1)
480  {
481  const polyPatch& newPatch = mesh.boundaryMesh()[newPatchi];
482  const polyPatch& oldPatch =
483  fldToAdd.mesh().boundaryMesh()[patchi];
484 
485  if (!bfld(newPatchi))
486  {
487  // First occurrence of newPatchi. Map from existing
488  // patchField
489 
490  // From new patch faces to patch faces on added mesh.
491  const labelList newToAdded
492  (
493  calcPatchMap
494  (
495  oldPatch.start(),
496  oldPatch.size(),
497  meshMap.addedFaceMap(),
498  newPatch,
499  -1 // unmapped values
500  )
501  );
502 
503  directFvPatchFieldMapper patchMapper(newToAdded);
504 
505  bfld.set
506  (
507  newPatchi,
509  (
510  fldToAdd.boundaryField()[patchi],// added field
511  mesh.boundary()[newPatchi], // new fvPatch
512  fld(), // new int. field
513  patchMapper // mapper
514  )
515  );
516  }
517  else
518  {
519  // PatchField will have correct size already. Just slot in
520  // my elements.
521 
522  labelList addedToNew(oldPatch.size(), -1);
523  forAll(addedToNew, i)
524  {
525  const label addedFacei = oldPatch.start() + i;
526  const label newFacei =
527  meshMap.addedFaceMap()[addedFacei];
528  const label patchFacei = newFacei-newPatch.start();
529 
530  if (patchFacei >= 0 && patchFacei < newPatch.size())
531  {
532  addedToNew[i] = patchFacei;
533  }
534  }
535 
536  bfld[newPatchi].map
537  (
538  fldToAdd.boundaryField()[patchi],
539  reverseFvPatchFieldMapper(addedToNew)
540  );
541  }
542  }
543  }
544  }
545 }
546 
547 
548 template<class Type>
550 (
551  const mapAddedPolyMesh& meshMap,
552  const fvMesh& mesh,
553  const fvMesh& meshToAdd
554 )
555 {
557  (
558  mesh.objectRegistry::lookupClass<SurfaceField<Type>>()
559  );
560 
562  (
563  meshToAdd.objectRegistry::lookupClass<SurfaceField<Type>>()
564  );
565 
566  for
567  (
568  typename HashTable<const SurfaceField<Type>*>::
569  iterator fieldIter = fields.begin();
570  fieldIter != fields.end();
571  ++fieldIter
572  )
573  {
574  if (fvMesh::geometryFields.found(fieldIter()->name())) continue;
575 
576  SurfaceField<Type>& fld = const_cast<SurfaceField<Type>&>(*fieldIter());
577 
578  if (fieldsToAdd.found(fld.name()))
579  {
580  const SurfaceField<Type>& fldToAdd = *fieldsToAdd[fld.name()];
581 
582  if (debug)
583  {
584  Pout<< "MapSurfaceFields : mapping " << fld.name()
585  << " and " << fldToAdd.name() << endl;
586  }
587 
588  MapSurfaceField<Type>(meshMap, fld, fldToAdd);
589  }
590  else
591  {
593  << "Not mapping field " << fld.name()
594  << " since not present on mesh to add"
595  << endl;
596  }
597  }
598 }
599 
600 
601 template<class Type>
602 void Foam::fvMeshAdder::MapPointField
603 (
604  const pointMesh& mesh,
605  const mapAddedPolyMesh& meshMap,
606  const labelListList& oldMeshPoints,
607 
609  const PointField<Type>& fldToAdd
610 )
611 {
612  // This is a bit tricky:
613  // - mesh pointed to by fld is invalid
614  // - pointPatches pointed to be fld are invalid
615 
616  typename PointField<Type>::
617  Boundary& bfld = fld.boundaryFieldRef();
618 
619  // Internal field
620  // ~~~~~~~~~~~~~~
621 
622  // Store old internal field
623  {
624  Field<Type> oldField(fld);
625 
626  // Modify internal field
627  Field<Type>& intFld = fld.primitiveFieldRef();
628 
629  intFld.setSize(mesh.size());
630 
631  intFld.rmap(oldField, meshMap.oldPointMap());
632  intFld.rmap(fldToAdd, meshMap.addedPointMap());
633  }
634 
635 
636  // Patch fields from old mesh
637  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
638 
639  {
640  const labelList& oldPatchMap = meshMap.oldPatchMap();
641 
642  // Reorder old patches in order of new ones. Put removed patches at end.
643 
644  label unusedPatchi = 0;
645 
646  forAll(oldPatchMap, patchi)
647  {
648  const label newPatchi = oldPatchMap[patchi];
649 
650  if (newPatchi != -1)
651  {
652  unusedPatchi++;
653  }
654  }
655 
656  label nUsedPatches = unusedPatchi;
657 
658  // Reorder list for patchFields
659  labelList oldToNew(oldPatchMap.size());
660 
661  forAll(oldPatchMap, patchi)
662  {
663  const label newPatchi = oldPatchMap[patchi];
664 
665  if (newPatchi != -1)
666  {
667  oldToNew[patchi] = newPatchi;
668  }
669  else
670  {
671  oldToNew[patchi] = unusedPatchi++;
672  }
673  }
674 
675 
676  // Sort deleted ones last so is now in newPatch ordering
677  bfld.reorder(oldToNew);
678 
679  // Extend to covers all patches
680  bfld.setSize(mesh.boundary().size());
681 
682  // Delete unused patches
683  for
684  (
685  label newPatchi = nUsedPatches;
686  newPatchi < bfld.size();
687  newPatchi++
688  )
689  {
690  bfld.set(newPatchi, nullptr);
691  }
692 
693 
694  // Map old values
695  // ~~~~~~~~~~~~~~
696 
697  forAll(oldPatchMap, patchi)
698  {
699  const label newPatchi = oldPatchMap[patchi];
700 
701  if (newPatchi != -1)
702  {
703  const labelList& oldMp = oldMeshPoints[patchi];
704  const pointPatch& newPp = mesh.boundary()[newPatchi];
705  const labelList& newMeshPoints = newPp.meshPoints();
706  const labelList& oldPointMap = meshMap.oldPointMap();
707 
708  Map<label> newMeshPointMap(2*newMeshPoints.size());
709  forAll(newMeshPoints, ppi)
710  {
711  newMeshPointMap.insert(newMeshPoints[ppi], ppi);
712  }
713 
714  labelList newToOld(newPp.size(), -1);
715  forAll(oldMp, oldPointi)
716  {
717  const label newPointi = oldPointMap[oldMp[oldPointi]];
718 
720  newMeshPointMap.find(newPointi);
721 
722  if (fnd == newMeshPointMap.end())
723  {
724  // Possibly an internal point
725  }
726  else
727  {
728  // Part of new patch
729  newToOld[fnd()] = oldPointi;
730  }
731  }
732 
733  directPointPatchFieldMapper patchMapper(newToOld);
734 
735  // Create new patchField with same type as existing one.
736  // Note:
737  // - boundaryField already in new order so access with newPatchi
738  // - bfld[newPatchi] both used for type and old
739  // value
740  // - hope that field mapping allows aliasing since old and new
741  // are same memory!
742  bfld.set
743  (
744  newPatchi,
746  (
747  bfld[newPatchi], // old field
748  mesh.boundary()[newPatchi], // new pointPatch
749  fld(), // new internal field
750  patchMapper // mapper (new to old)
751  )
752  );
753  }
754  }
755  }
756 
757 
758 
759  // Patch fields from added mesh
760  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
761 
762  {
763  const labelList& addedPatchMap = meshMap.addedPatchMap();
764 
765  // Add addedMesh patches
766  forAll(addedPatchMap, patchi)
767  {
768  const label newPatchi = addedPatchMap[patchi];
769 
770  if (newPatchi != -1)
771  {
772  const pointPatch& oldPatch = fldToAdd.mesh().boundary()[patchi];
773  const labelList& oldMp = oldPatch.meshPoints();
774  const pointPatch& newPp = mesh.boundary()[newPatchi];
775  const labelList& newMeshPoints = newPp.meshPoints();
776  const labelList& addedPointMap = meshMap.addedPointMap();
777 
778  Map<label> newMpm(2*newMeshPoints.size());
779  forAll(newMeshPoints, ppi)
780  {
781  newMpm.insert(newMeshPoints[ppi], ppi);
782  }
783 
784  if (!bfld(newPatchi))
785  {
786  // First occurrence of newPatchi. Map from existing
787  // patchField
788 
789  labelList newToAdded(newPp.size(), -1);
790  forAll(oldMp, oldPointi)
791  {
792  const label newPointi = addedPointMap[oldMp[oldPointi]];
793 
794  Map<label>::const_iterator fnd = newMpm.find(newPointi);
795  if (fnd == newMpm.end())
796  {
797  // Possibly an internal point
798  }
799  else
800  {
801  // Part of new patch
802  newToAdded[fnd()] = oldPointi;
803  }
804  }
805 
806  directPointPatchFieldMapper patchMapper(newToAdded);
807 
808  bfld.set
809  (
810  newPatchi,
812  (
813  fldToAdd.boundaryField()[patchi],// added field
814  mesh.boundary()[newPatchi], // new pointPatch
815  fld(), // new int. field
816  patchMapper // mapper
817  )
818  );
819  }
820  else
821  {
822  // PatchField will have correct size already. Just slot in
823  // my elements.
824 
825  labelList oldToNew(oldPatch.size(), -1);
826  forAll(oldMp, oldPointi)
827  {
828  const label newPointi = addedPointMap[oldMp[oldPointi]];
829 
830  Map<label>::const_iterator fnd = newMpm.find(newPointi);
831  if (fnd != newMpm.end())
832  {
833  // Part of new patch
834  oldToNew[oldPointi] = fnd();
835  }
836  }
837 
838  bfld[newPatchi].map
839  (
840  fldToAdd.boundaryField()[patchi],
841  reversePointPatchFieldMapper(oldToNew)
842  );
843  }
844  }
845  }
846  }
847 }
848 
849 
850 template<class Type>
852 (
853  const mapAddedPolyMesh& meshMap,
854  const pointMesh& mesh,
855  const labelListList& oldMeshPoints,
856  const objectRegistry& meshToAdd
857 )
858 {
860  (
862  );
863 
865  (
866  meshToAdd.lookupClass<PointField<Type>>()
867  );
868 
869  for
870  (
871  typename HashTable<const PointField<Type>*>::
872  iterator fieldIter = fields.begin();
873  fieldIter != fields.end();
874  ++fieldIter
875  )
876  {
877  if (fvMesh::geometryFields.found(fieldIter()->name())) continue;
878 
879  PointField<Type>& fld = const_cast<PointField<Type>&>(*fieldIter());
880 
881  if (fieldsToAdd.found(fld.name()))
882  {
883  const PointField<Type>& fldToAdd = *fieldsToAdd[fld.name()];
884 
885  if (debug)
886  {
887  Pout<< "MapPointFields : mapping " << fld.name()
888  << " and " << fldToAdd.name() << endl;
889  }
890 
891  MapPointField<Type>(mesh, meshMap, oldMeshPoints, fld, fldToAdd);
892  }
893  else
894  {
896  << "Not mapping field " << fld.name()
897  << " since not present on mesh to add"
898  << endl;
899  }
900  }
901 }
902 
903 
904 template<class Type>
905 void Foam::fvMeshAdder::MapDimField
906 (
907  const mapAddedPolyMesh& meshMap,
908 
910  const DimensionedField<Type, volMesh>& fldToAdd
911 )
912 {
913  const fvMesh& mesh = fld.mesh();
914 
915  // Store old field
916  Field<Type> oldField(fld);
917 
918  fld.setSize(mesh.nCells());
919 
920  fld.rmap(oldField, meshMap.oldCellMap());
921  fld.rmap(fldToAdd, meshMap.addedCellMap());
922 }
923 
924 
925 template<class Type>
927 (
928  const mapAddedPolyMesh& meshMap,
929  const fvMesh& mesh,
930  const fvMesh& meshToAdd
931 )
932 {
933  // Note: use strict flag on lookupClass to avoid picking up volFields
935  (
936  mesh.objectRegistry::lookupClass<VolInternalField<Type>>(true)
937  );
938 
940  (
941  meshToAdd.objectRegistry::lookupClass<VolInternalField<Type>>(true)
942  );
943 
944  for
945  (
946  typename HashTable<const VolInternalField<Type>*>::
947  iterator fieldIter = fields.begin();
948  fieldIter != fields.end();
949  ++fieldIter
950  )
951  {
953  const_cast<VolInternalField<Type>&>(*fieldIter());
954 
955  if (fieldsToAdd.found(fld.name()))
956  {
957  const VolInternalField<Type>& fldToAdd = *fieldsToAdd[fld.name()];
958 
959  if (debug)
960  {
961  Pout<< "MapDimFields : mapping " << fld.name()
962  << " and " << fldToAdd.name() << endl;
963  }
964 
965  MapDimField<Type>(meshMap, fld, fldToAdd);
966  }
967  else
968  {
969  WarningIn("fvMeshAdder::MapDimFields(..)")
970  << "Not mapping field " << fld.name()
971  << " since not present on mesh to add"
972  << endl;
973  }
974  }
975 }
976 
977 
978 // ************************************************************************* //
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const
Return mesh.
Pre-declare SubField and related Field type.
Definition: Field.H:82
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition: Field.C:265
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:362
Generic GeometricField class.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
An STL-conforming hash table.
Definition: HashTable.H:127
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
const word & name() const
Return name.
Definition: IOobject.H:310
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
Map point field on topology change. This is a partial template specialisation for GeoMesh=pointMesh.
HashTable< label, label, Hash< label > >::const_iterator const_iterator
Definition: Map.H:59
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
static void MapVolFields(const mapAddedPolyMesh &, const fvMesh &mesh, const fvMesh &meshToAdd)
Map all volFields of Type.
static void MapDimFields(const mapAddedPolyMesh &, const fvMesh &mesh, const fvMesh &meshToAdd)
Map all DimensionedFields of Type.
static void MapSurfaceFields(const mapAddedPolyMesh &, const fvMesh &mesh, const fvMesh &meshToAdd)
Map all surfaceFields of Type.
static void MapPointFields(const mapAddedPolyMesh &, const pointMesh &mesh, const labelListList &oldMeshPoints, const objectRegistry &meshToAdd)
Map all surfaceFields of Type.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:895
static const HashSet< word > geometryFields
Set of names of registered geometric fields.
Definition: fvMesh.H:317
static tmp< fvPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:82
static tmp< fvsPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, surfaceMesh > &)
Return a pointer to a new patchField created on freestore given.
Class containing mesh-to-mesh mapping information after a mesh addition where we add a mesh ('added m...
const labelList & addedFaceMap() const
const labelList & oldPatchMap() const
From old patch index to new patch index or -1 if patch.
const labelList & addedCellMap() const
const labelList & addedPatchMap() const
From added mesh patch index to new patch index or -1 if.
const labelList & oldPointMap() const
From old mesh point/face/cell to new mesh point/face/cell.
const labelList & oldCellMap() const
const labelList & oldPatchStarts() const
Return list of the old patch start labels.
const labelList & addedPointMap() const
From added mesh point/face/cell to new mesh point/face/cell.
const labelList & oldFaceMap() const
const labelList & oldPatchSizes() const
Return list of the old patch sizes.
Registry of regIOobjects.
HashTable< const Type * > lookupClass(const bool strict=false) const
Lookup and return all objects of the given Type.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:52
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:124
label size() const
Return number of points.
Definition: pointMesh.H:100
const pointBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: pointMesh.H:112
static autoPtr< pointPatchField< Type > > New(const word &, const pointPatch &, const DimensionedField< Type, pointMesh > &)
Return a pointer to a new patchField created on freestore given.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:405
label nCells() const
label nInternalFaces() const
label patchi
gmvFile<< "tracers "<< particles.size()<< nl;{ pointField positions(particles.size());label particlei=0;forAllConstIter(Cloud< passiveParticle >, particles, iter) { positions[particlei++]=iter().position(mesh);} for(i=0;i< pTraits< point >::nComponents;i++) { forAll(positions, particlei) { gmvFile<< component(positions[particlei], i)<< ' ';} gmvFile<< nl;}}forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
#define WarningInFunction
Report a warning using Foam::Warning.
#define WarningIn(functionName)
Report a warning using Foam::Warning.
List< label > labelList
A List of labels.
Definition: labelList.H:56
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
typename VolField< Type >::Internal VolInternalField
Definition: volFieldsFwd.H:58
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
label newPointi
Definition: readKivaGrid.H:495
Foam::surfaceFields.