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  {
280  const_cast<VolField<Type>&>
281  (
282  *fieldIter()
283  );
284 
285  if (fieldsToAdd.found(fld.name()))
286  {
287  const VolField<Type>& fldToAdd =
288  *fieldsToAdd[fld.name()];
289 
290  if (debug)
291  {
292  Pout<< "MapVolFields : mapping " << fld.name()
293  << " and " << fldToAdd.name() << endl;
294  }
295 
296  MapVolField<Type>(meshMap, fld, fldToAdd);
297  }
298  else
299  {
301  << "Not mapping field " << fld.name()
302  << " since not present on mesh to add"
303  << endl;
304  }
305  }
306 }
307 
308 
309 template<class Type>
310 void Foam::fvMeshAdder::MapSurfaceField
311 (
312  const mapAddedPolyMesh& meshMap,
313 
315  const SurfaceField<Type>& fldToAdd
316 )
317 {
318  const fvMesh& mesh = fld.mesh();
319  const labelList& oldPatchStarts = meshMap.oldPatchStarts();
320 
321  typename SurfaceField<Type>::
322  Boundary& bfld = fld.boundaryFieldRef();
323 
324  // Internal field
325  // ~~~~~~~~~~~~~~
326 
327  // Store old internal field
328  {
329  Field<Type> oldField(fld);
330 
331  // Modify internal field
332  Field<Type>& intFld = fld.primitiveFieldRef();
333 
334  intFld.setSize(mesh.nInternalFaces());
335 
336  intFld.rmap(oldField, meshMap.oldFaceMap());
337  intFld.rmap(fldToAdd, meshMap.addedFaceMap());
338 
339 
340  // Faces that were boundary faces but are not anymore.
341  // Use owner value (so lowest numbered cell, i.e. from 'old' not 'added'
342  // mesh)
343  forAll(bfld, patchi)
344  {
345  const fvsPatchField<Type>& pf = bfld[patchi];
346 
347  const label start = oldPatchStarts[patchi];
348 
349  forAll(pf, i)
350  {
351  const label newFacei = meshMap.oldFaceMap()[start + i];
352 
353  if (newFacei >= 0 && newFacei < mesh.nInternalFaces())
354  {
355  intFld[newFacei] = pf[i];
356  }
357  }
358  }
359  }
360 
361 
362  // Patch fields from old mesh
363  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
364 
365  {
366  const labelList& oldPatchMap = meshMap.oldPatchMap();
367  const labelList& oldPatchSizes = meshMap.oldPatchSizes();
368 
369  // Reorder old patches in order of new ones. Put removed patches at end.
370 
371  label unusedPatchi = 0;
372 
373  forAll(oldPatchMap, patchi)
374  {
375  const label newPatchi = oldPatchMap[patchi];
376 
377  if (newPatchi != -1)
378  {
379  unusedPatchi++;
380  }
381  }
382 
383  label nUsedPatches = unusedPatchi;
384 
385  // Reorder list for patchFields
386  labelList oldToNew(oldPatchMap.size());
387 
388  forAll(oldPatchMap, patchi)
389  {
390  const label newPatchi = oldPatchMap[patchi];
391 
392  if (newPatchi != -1)
393  {
394  oldToNew[patchi] = newPatchi;
395  }
396  else
397  {
398  oldToNew[patchi] = unusedPatchi++;
399  }
400  }
401 
402 
403  // Sort deleted ones last so is now in newPatch ordering
404  bfld.reorder(oldToNew);
405  // Extend to covers all patches
406  bfld.setSize(mesh.boundaryMesh().size());
407  // Delete unused patches
408  for
409  (
410  label newPatchi = nUsedPatches;
411  newPatchi < bfld.size();
412  newPatchi++
413  )
414  {
415  bfld.set(newPatchi, nullptr);
416  }
417 
418 
419  // Map old values
420  // ~~~~~~~~~~~~~~
421 
422  forAll(oldPatchMap, patchi)
423  {
424  const label newPatchi = oldPatchMap[patchi];
425 
426  if (newPatchi != -1)
427  {
428  const labelList newToOld
429  (
430  calcPatchMap
431  (
432  oldPatchStarts[patchi],
433  oldPatchSizes[patchi],
434  meshMap.oldFaceMap(),
435  mesh.boundaryMesh()[newPatchi],
436  -1 // unmapped value
437  )
438  );
439 
440  directFvPatchFieldMapper patchMapper(newToOld);
441 
442  // Create new patchField with same type as existing one.
443  // Note:
444  // - boundaryField already in new order so access with newPatchi
445  // - bfld[newPatchi] both used for type and old
446  // value
447  // - hope that field mapping allows aliasing since old and new
448  // are same memory!
449  bfld.set
450  (
451  newPatchi,
453  (
454  bfld[newPatchi], // old field
455  mesh.boundary()[newPatchi], // new fvPatch
456  fld(), // new internal field
457  patchMapper // mapper (new to old)
458  )
459  );
460  }
461  }
462  }
463 
464 
465 
466  // Patch fields from added mesh
467  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
468 
469  {
470  const labelList& addedPatchMap = meshMap.addedPatchMap();
471 
472  // Add addedMesh patches
473  forAll(addedPatchMap, patchi)
474  {
475  const label newPatchi = addedPatchMap[patchi];
476 
477  if (newPatchi != -1)
478  {
479  const polyPatch& newPatch = mesh.boundaryMesh()[newPatchi];
480  const polyPatch& oldPatch =
481  fldToAdd.mesh().boundaryMesh()[patchi];
482 
483  if (!bfld(newPatchi))
484  {
485  // First occurrence of newPatchi. Map from existing
486  // patchField
487 
488  // From new patch faces to patch faces on added mesh.
489  const labelList newToAdded
490  (
491  calcPatchMap
492  (
493  oldPatch.start(),
494  oldPatch.size(),
495  meshMap.addedFaceMap(),
496  newPatch,
497  -1 // unmapped values
498  )
499  );
500 
501  directFvPatchFieldMapper patchMapper(newToAdded);
502 
503  bfld.set
504  (
505  newPatchi,
507  (
508  fldToAdd.boundaryField()[patchi],// added field
509  mesh.boundary()[newPatchi], // new fvPatch
510  fld(), // new int. field
511  patchMapper // mapper
512  )
513  );
514  }
515  else
516  {
517  // PatchField will have correct size already. Just slot in
518  // my elements.
519 
520  labelList addedToNew(oldPatch.size(), -1);
521  forAll(addedToNew, i)
522  {
523  const label addedFacei = oldPatch.start() + i;
524  const label newFacei =
525  meshMap.addedFaceMap()[addedFacei];
526  const label patchFacei = newFacei-newPatch.start();
527 
528  if (patchFacei >= 0 && patchFacei < newPatch.size())
529  {
530  addedToNew[i] = patchFacei;
531  }
532  }
533 
534  bfld[newPatchi].map
535  (
536  fldToAdd.boundaryField()[patchi],
537  reverseFvPatchFieldMapper(addedToNew)
538  );
539  }
540  }
541  }
542  }
543 }
544 
545 
546 template<class Type>
548 (
549  const mapAddedPolyMesh& meshMap,
550  const fvMesh& mesh,
551  const fvMesh& meshToAdd
552 )
553 {
555  (
556  mesh.objectRegistry::lookupClass<SurfaceField<Type>>()
557  );
558 
560  (
561  meshToAdd.objectRegistry::lookupClass<SurfaceField<Type>>()
562  );
563 
564  for
565  (
566  typename HashTable<const SurfaceField<Type>*>::
567  iterator fieldIter = fields.begin();
568  fieldIter != fields.end();
569  ++fieldIter
570  )
571  {
572  SurfaceField<Type>& fld = const_cast<SurfaceField<Type>&>(*fieldIter());
573 
574  if (fieldsToAdd.found(fld.name()))
575  {
576  const SurfaceField<Type>& fldToAdd = *fieldsToAdd[fld.name()];
577 
578  if (debug)
579  {
580  Pout<< "MapSurfaceFields : mapping " << fld.name()
581  << " and " << fldToAdd.name() << endl;
582  }
583 
584  MapSurfaceField<Type>(meshMap, fld, fldToAdd);
585  }
586  else
587  {
589  << "Not mapping field " << fld.name()
590  << " since not present on mesh to add"
591  << endl;
592  }
593  }
594 }
595 
596 
597 template<class Type>
598 void Foam::fvMeshAdder::MapPointField
599 (
600  const pointMesh& mesh,
601  const mapAddedPolyMesh& meshMap,
602  const labelListList& oldMeshPoints,
603 
605  const PointField<Type>& fldToAdd
606 )
607 {
608  // This is a bit tricky:
609  // - mesh pointed to by fld is invalid
610  // - pointPatches pointed to be fld are invalid
611 
612  typename PointField<Type>::
613  Boundary& bfld = fld.boundaryFieldRef();
614 
615  // Internal field
616  // ~~~~~~~~~~~~~~
617 
618  // Store old internal field
619  {
620  Field<Type> oldField(fld);
621 
622  // Modify internal field
623  Field<Type>& intFld = fld.primitiveFieldRef();
624 
625  intFld.setSize(mesh.size());
626 
627  intFld.rmap(oldField, meshMap.oldPointMap());
628  intFld.rmap(fldToAdd, meshMap.addedPointMap());
629  }
630 
631 
632  // Patch fields from old mesh
633  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
634 
635  {
636  const labelList& oldPatchMap = meshMap.oldPatchMap();
637 
638  // Reorder old patches in order of new ones. Put removed patches at end.
639 
640  label unusedPatchi = 0;
641 
642  forAll(oldPatchMap, patchi)
643  {
644  const label newPatchi = oldPatchMap[patchi];
645 
646  if (newPatchi != -1)
647  {
648  unusedPatchi++;
649  }
650  }
651 
652  label nUsedPatches = unusedPatchi;
653 
654  // Reorder list for patchFields
655  labelList oldToNew(oldPatchMap.size());
656 
657  forAll(oldPatchMap, patchi)
658  {
659  const label newPatchi = oldPatchMap[patchi];
660 
661  if (newPatchi != -1)
662  {
663  oldToNew[patchi] = newPatchi;
664  }
665  else
666  {
667  oldToNew[patchi] = unusedPatchi++;
668  }
669  }
670 
671 
672  // Sort deleted ones last so is now in newPatch ordering
673  bfld.reorder(oldToNew);
674 
675  // Extend to covers all patches
676  bfld.setSize(mesh.boundary().size());
677 
678  // Delete unused patches
679  for
680  (
681  label newPatchi = nUsedPatches;
682  newPatchi < bfld.size();
683  newPatchi++
684  )
685  {
686  bfld.set(newPatchi, nullptr);
687  }
688 
689 
690  // Map old values
691  // ~~~~~~~~~~~~~~
692 
693  forAll(oldPatchMap, patchi)
694  {
695  const label newPatchi = oldPatchMap[patchi];
696 
697  if (newPatchi != -1)
698  {
699  const labelList& oldMp = oldMeshPoints[patchi];
700  const pointPatch& newPp = mesh.boundary()[newPatchi];
701  const labelList& newMeshPoints = newPp.meshPoints();
702  const labelList& oldPointMap = meshMap.oldPointMap();
703 
704  Map<label> newMeshPointMap(2*newMeshPoints.size());
705  forAll(newMeshPoints, ppi)
706  {
707  newMeshPointMap.insert(newMeshPoints[ppi], ppi);
708  }
709 
710  labelList newToOld(newPp.size(), -1);
711  forAll(oldMp, oldPointi)
712  {
713  const label newPointi = oldPointMap[oldMp[oldPointi]];
714 
716  newMeshPointMap.find(newPointi);
717 
718  if (fnd == newMeshPointMap.end())
719  {
720  // Possibly an internal point
721  }
722  else
723  {
724  // Part of new patch
725  newToOld[fnd()] = oldPointi;
726  }
727  }
728 
729  directPointPatchFieldMapper patchMapper(newToOld);
730 
731  // Create new patchField with same type as existing one.
732  // Note:
733  // - boundaryField already in new order so access with newPatchi
734  // - bfld[newPatchi] both used for type and old
735  // value
736  // - hope that field mapping allows aliasing since old and new
737  // are same memory!
738  bfld.set
739  (
740  newPatchi,
742  (
743  bfld[newPatchi], // old field
744  mesh.boundary()[newPatchi], // new pointPatch
745  fld(), // new internal field
746  patchMapper // mapper (new to old)
747  )
748  );
749  }
750  }
751  }
752 
753 
754 
755  // Patch fields from added mesh
756  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
757 
758  {
759  const labelList& addedPatchMap = meshMap.addedPatchMap();
760 
761  // Add addedMesh patches
762  forAll(addedPatchMap, patchi)
763  {
764  const label newPatchi = addedPatchMap[patchi];
765 
766  if (newPatchi != -1)
767  {
768  const pointPatch& oldPatch = fldToAdd.mesh().boundary()[patchi];
769  const labelList& oldMp = oldPatch.meshPoints();
770  const pointPatch& newPp = mesh.boundary()[newPatchi];
771  const labelList& newMeshPoints = newPp.meshPoints();
772  const labelList& addedPointMap = meshMap.addedPointMap();
773 
774  Map<label> newMpm(2*newMeshPoints.size());
775  forAll(newMeshPoints, ppi)
776  {
777  newMpm.insert(newMeshPoints[ppi], ppi);
778  }
779 
780  if (!bfld(newPatchi))
781  {
782  // First occurrence of newPatchi. Map from existing
783  // patchField
784 
785  labelList newToAdded(newPp.size(), -1);
786  forAll(oldMp, oldPointi)
787  {
788  const label newPointi = addedPointMap[oldMp[oldPointi]];
789 
790  Map<label>::const_iterator fnd = newMpm.find(newPointi);
791  if (fnd == newMpm.end())
792  {
793  // Possibly an internal point
794  }
795  else
796  {
797  // Part of new patch
798  newToAdded[fnd()] = oldPointi;
799  }
800  }
801 
802  directPointPatchFieldMapper patchMapper(newToAdded);
803 
804  bfld.set
805  (
806  newPatchi,
808  (
809  fldToAdd.boundaryField()[patchi],// added field
810  mesh.boundary()[newPatchi], // new pointPatch
811  fld(), // new int. field
812  patchMapper // mapper
813  )
814  );
815  }
816  else
817  {
818  // PatchField will have correct size already. Just slot in
819  // my elements.
820 
821  labelList oldToNew(oldPatch.size(), -1);
822  forAll(oldMp, oldPointi)
823  {
824  const label newPointi = addedPointMap[oldMp[oldPointi]];
825 
826  Map<label>::const_iterator fnd = newMpm.find(newPointi);
827  if (fnd != newMpm.end())
828  {
829  // Part of new patch
830  oldToNew[oldPointi] = fnd();
831  }
832  }
833 
834  bfld[newPatchi].map
835  (
836  fldToAdd.boundaryField()[patchi],
837  reversePointPatchFieldMapper(oldToNew)
838  );
839  }
840  }
841  }
842  }
843 }
844 
845 
846 template<class Type>
848 (
849  const mapAddedPolyMesh& meshMap,
850  const pointMesh& mesh,
851  const labelListList& oldMeshPoints,
852  const objectRegistry& meshToAdd
853 )
854 {
856  (
858  );
859 
861  (
862  meshToAdd.lookupClass<PointField<Type>>()
863  );
864 
865  for
866  (
867  typename HashTable<const PointField<Type>*>::
868  iterator fieldIter = fields.begin();
869  fieldIter != fields.end();
870  ++fieldIter
871  )
872  {
873  PointField<Type>& fld = const_cast<PointField<Type>&>(*fieldIter());
874 
875  if (fieldsToAdd.found(fld.name()))
876  {
877  const PointField<Type>& fldToAdd = *fieldsToAdd[fld.name()];
878 
879  if (debug)
880  {
881  Pout<< "MapPointFields : mapping " << fld.name()
882  << " and " << fldToAdd.name() << endl;
883  }
884 
885  MapPointField<Type>(mesh, meshMap, oldMeshPoints, fld, fldToAdd);
886  }
887  else
888  {
890  << "Not mapping field " << fld.name()
891  << " since not present on mesh to add"
892  << endl;
893  }
894  }
895 }
896 
897 
898 template<class Type>
899 void Foam::fvMeshAdder::MapDimField
900 (
901  const mapAddedPolyMesh& meshMap,
902 
904  const DimensionedField<Type, volMesh>& fldToAdd
905 )
906 {
907  const fvMesh& mesh = fld.mesh();
908 
909  // Store old field
910  Field<Type> oldField(fld);
911 
912  fld.setSize(mesh.nCells());
913 
914  fld.rmap(oldField, meshMap.oldCellMap());
915  fld.rmap(fldToAdd, meshMap.addedCellMap());
916 }
917 
918 
919 template<class Type>
921 (
922  const mapAddedPolyMesh& meshMap,
923  const fvMesh& mesh,
924  const fvMesh& meshToAdd
925 )
926 {
927  // Note: use strict flag on lookupClass to avoid picking up volFields
929  (
930  mesh.objectRegistry::lookupClass<VolInternalField<Type>>(true)
931  );
932 
934  (
935  meshToAdd.objectRegistry::lookupClass<VolInternalField<Type>>(true)
936  );
937 
938  for
939  (
940  typename HashTable<const VolInternalField<Type>*>::
941  iterator fieldIter = fields.begin();
942  fieldIter != fields.end();
943  ++fieldIter
944  )
945  {
947  const_cast<VolInternalField<Type>&>(*fieldIter());
948 
949  if (fieldsToAdd.found(fld.name()))
950  {
951  const VolInternalField<Type>& fldToAdd = *fieldsToAdd[fld.name()];
952 
953  if (debug)
954  {
955  Pout<< "MapDimFields : mapping " << fld.name()
956  << " and " << fldToAdd.name() << endl;
957  }
958 
959  MapDimField<Type>(meshMap, fld, fldToAdd);
960  }
961  else
962  {
963  WarningIn("fvMeshAdder::MapDimFields(..)")
964  << "Not mapping field " << fld.name()
965  << " since not present on mesh to add"
966  << endl;
967  }
968  }
969 }
970 
971 
972 // ************************************************************************* //
#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:893
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:116
label size() const
Return number of points.
Definition: pointMesh.H:92
const pointBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: pointMesh.H:104
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:403
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
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.