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-2019 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"
31 
32 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
33 
34 template<class Type>
35 void Foam::fvMeshAdder::MapVolField
36 (
37  const mapAddedPolyMesh& meshMap,
38 
39  GeometricField<Type, fvPatchField, volMesh>& fld,
40  const GeometricField<Type, fvPatchField, volMesh>& fldToAdd
41 )
42 {
43  const fvMesh& mesh = fld.mesh();
44 
45  // Internal field
46  // ~~~~~~~~~~~~~~
47 
48  {
49  // Store old internal field
50  Field<Type> oldInternalField(fld.primitiveField());
51 
52  // Modify internal field
53  Field<Type>& intFld = fld.primitiveFieldRef();
54 
55  intFld.setSize(mesh.nCells());
56 
57  intFld.rmap(oldInternalField, meshMap.oldCellMap());
58  intFld.rmap(fldToAdd.primitiveField(), meshMap.addedCellMap());
59  }
60 
61 
62  // Patch fields from old mesh
63  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
64 
65  typename GeometricField<Type, fvPatchField, volMesh>::
66  Boundary& bfld = fld.boundaryFieldRef();
67 
68  {
69  const labelList& oldPatchMap = meshMap.oldPatchMap();
70  const labelList& oldPatchStarts = meshMap.oldPatchStarts();
71  const labelList& oldPatchSizes = meshMap.oldPatchSizes();
72 
73  // Reorder old patches in order of new ones. Put removed patches at end.
74 
75  label unusedPatchi = 0;
76 
77  forAll(oldPatchMap, patchi)
78  {
79  label newPatchi = oldPatchMap[patchi];
80 
81  if (newPatchi != -1)
82  {
83  unusedPatchi++;
84  }
85  }
86 
87  label nUsedPatches = unusedPatchi;
88 
89  // Reorder list for patchFields
90  labelList oldToNew(oldPatchMap.size());
91 
92  forAll(oldPatchMap, patchi)
93  {
94  label newPatchi = oldPatchMap[patchi];
95 
96  if (newPatchi != -1)
97  {
98  oldToNew[patchi] = newPatchi;
99  }
100  else
101  {
102  oldToNew[patchi] = unusedPatchi++;
103  }
104  }
105 
106 
107  // Sort deleted ones last so is now in newPatch ordering
108  bfld.reorder(oldToNew);
109  // Extend to covers all patches
110  bfld.setSize(mesh.boundaryMesh().size());
111  // Delete unused patches
112  for
113  (
114  label newPatchi = nUsedPatches;
115  newPatchi < bfld.size();
116  newPatchi++
117  )
118  {
119  bfld.set(newPatchi, nullptr);
120  }
121 
122 
123  // Map old values
124  // ~~~~~~~~~~~~~~
125 
126  forAll(oldPatchMap, patchi)
127  {
128  label newPatchi = oldPatchMap[patchi];
129 
130  if (newPatchi != -1)
131  {
132  labelList newToOld
133  (
134  calcPatchMap
135  (
136  oldPatchStarts[patchi],
137  oldPatchSizes[patchi],
138  meshMap.oldFaceMap(),
139  mesh.boundaryMesh()[newPatchi],
140  -1 // unmapped value
141  )
142  );
143 
144  directFvPatchFieldMapper patchMapper(newToOld);
145 
146 
147  // Create new patchField with same type as existing one.
148  // Note:
149  // - boundaryField already in new order so access with newPatchi
150  // - fld.boundaryField()[newPatchi] both used for type and old
151  // value
152  // - hope that field mapping allows aliasing since old and new
153  // are same memory!
154  bfld.set
155  (
156  newPatchi,
158  (
159  bfld[newPatchi], // old field
160  mesh.boundary()[newPatchi], // new fvPatch
161  fld(), // new internal field
162  patchMapper // mapper (new to old)
163  )
164  );
165  }
166  }
167  }
168 
169 
170 
171  // Patch fields from added mesh
172  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
173 
174  {
175  const labelList& addedPatchMap = meshMap.addedPatchMap();
176 
177  // Add addedMesh patches
178  forAll(addedPatchMap, patchi)
179  {
180  label newPatchi = addedPatchMap[patchi];
181 
182  if (newPatchi != -1)
183  {
184  const polyPatch& newPatch = mesh.boundaryMesh()[newPatchi];
185  const polyPatch& oldPatch =
186  fldToAdd.mesh().boundaryMesh()[patchi];
187 
188  if (!bfld(newPatchi))
189  {
190  // First occurrence of newPatchi. Map from existing
191  // patchField
192 
193  // From new patch faces to patch faces on added mesh.
194  labelList newToAdded
195  (
196  calcPatchMap
197  (
198  oldPatch.start(),
199  oldPatch.size(),
200  meshMap.addedFaceMap(),
201  newPatch,
202  -1 // unmapped values
203  )
204  );
205 
206  directFvPatchFieldMapper patchMapper(newToAdded);
207 
208  bfld.set
209  (
210  newPatchi,
212  (
213  fldToAdd.boundaryField()[patchi], // added field
214  mesh.boundary()[newPatchi], // new fvPatch
215  fld(), // new int. field
216  patchMapper // mapper
217  )
218  );
219  }
220  else
221  {
222  // PatchField will have correct size already. Just slot in
223  // my elements.
224 
225  labelList addedToNew(oldPatch.size(), -1);
226  forAll(addedToNew, i)
227  {
228  label addedFacei = oldPatch.start()+i;
229  label newFacei = meshMap.addedFaceMap()[addedFacei];
230  label patchFacei = newFacei-newPatch.start();
231  if (patchFacei >= 0 && patchFacei < newPatch.size())
232  {
233  addedToNew[i] = patchFacei;
234  }
235  }
236 
237  bfld[newPatchi].rmap
238  (
239  fldToAdd.boundaryField()[patchi],
240  addedToNew
241  );
242  }
243  }
244  }
245  }
246 }
247 
248 
249 template<class Type>
251 (
252  const mapAddedPolyMesh& meshMap,
253  const fvMesh& mesh,
254  const fvMesh& meshToAdd
255 )
256 {
258  (
259  mesh.objectRegistry::lookupClass
261  ()
262  );
263 
265  (
266  meshToAdd.objectRegistry::lookupClass
268  ()
269  );
270 
271  // It is necessary to enforce that all old-time fields are stored
272  // before the mapping is performed. Otherwise, if the
273  // old-time-level field is mapped before the field itself, sizes
274  // will not match.
275 
276  for
277  (
279  iterator fieldIter = fields.begin();
280  fieldIter != fields.end();
281  ++fieldIter
282  )
283  {
284  if (debug)
285  {
286  Pout<< "MapVolFields : Storing old time for " << fieldIter()->name()
287  << endl;
288  }
289 
290  const_cast<GeometricField<Type, fvPatchField, volMesh>*>(fieldIter())
291  ->storeOldTimes();
292  }
293 
294 
295  for
296  (
298  iterator fieldIter = fields.begin();
299  fieldIter != fields.end();
300  ++fieldIter
301  )
302  {
305  (
306  *fieldIter()
307  );
308 
309  if (fieldsToAdd.found(fld.name()))
310  {
312  *fieldsToAdd[fld.name()];
313 
314  if (debug)
315  {
316  Pout<< "MapVolFields : mapping " << fld.name()
317  << " and " << fldToAdd.name() << endl;
318  }
319 
320  MapVolField<Type>(meshMap, fld, fldToAdd);
321  }
322  else
323  {
325  << "Not mapping field " << fld.name()
326  << " since not present on mesh to add"
327  << endl;
328  }
329  }
330 }
331 
332 
333 template<class Type>
334 void Foam::fvMeshAdder::MapSurfaceField
335 (
336  const mapAddedPolyMesh& meshMap,
337 
340 )
341 {
342  const fvMesh& mesh = fld.mesh();
343  const labelList& oldPatchStarts = meshMap.oldPatchStarts();
344 
346  Boundary& bfld = fld.boundaryFieldRef();
347 
348  // Internal field
349  // ~~~~~~~~~~~~~~
350 
351  // Store old internal field
352  {
353  Field<Type> oldField(fld);
354 
355  // Modify internal field
356  Field<Type>& intFld = fld.primitiveFieldRef();
357 
358  intFld.setSize(mesh.nInternalFaces());
359 
360  intFld.rmap(oldField, meshMap.oldFaceMap());
361  intFld.rmap(fldToAdd, meshMap.addedFaceMap());
362 
363 
364  // Faces that were boundary faces but are not anymore.
365  // Use owner value (so lowest numbered cell, i.e. from 'old' not 'added'
366  // mesh)
367  forAll(bfld, patchi)
368  {
369  const fvsPatchField<Type>& pf = bfld[patchi];
370 
371  label start = oldPatchStarts[patchi];
372 
373  forAll(pf, i)
374  {
375  label newFacei = meshMap.oldFaceMap()[start + i];
376 
377  if (newFacei >= 0 && newFacei < mesh.nInternalFaces())
378  {
379  intFld[newFacei] = pf[i];
380  }
381  }
382  }
383  }
384 
385 
386  // Patch fields from old mesh
387  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
388 
389  {
390  const labelList& oldPatchMap = meshMap.oldPatchMap();
391  const labelList& oldPatchSizes = meshMap.oldPatchSizes();
392 
393  // Reorder old patches in order of new ones. Put removed patches at end.
394 
395  label unusedPatchi = 0;
396 
397  forAll(oldPatchMap, patchi)
398  {
399  label newPatchi = oldPatchMap[patchi];
400 
401  if (newPatchi != -1)
402  {
403  unusedPatchi++;
404  }
405  }
406 
407  label nUsedPatches = unusedPatchi;
408 
409  // Reorder list for patchFields
410  labelList oldToNew(oldPatchMap.size());
411 
412  forAll(oldPatchMap, patchi)
413  {
414  label newPatchi = oldPatchMap[patchi];
415 
416  if (newPatchi != -1)
417  {
418  oldToNew[patchi] = newPatchi;
419  }
420  else
421  {
422  oldToNew[patchi] = unusedPatchi++;
423  }
424  }
425 
426 
427  // Sort deleted ones last so is now in newPatch ordering
428  bfld.reorder(oldToNew);
429  // Extend to covers all patches
430  bfld.setSize(mesh.boundaryMesh().size());
431  // Delete unused patches
432  for
433  (
434  label newPatchi = nUsedPatches;
435  newPatchi < bfld.size();
436  newPatchi++
437  )
438  {
439  bfld.set(newPatchi, nullptr);
440  }
441 
442 
443  // Map old values
444  // ~~~~~~~~~~~~~~
445 
446  forAll(oldPatchMap, patchi)
447  {
448  label newPatchi = oldPatchMap[patchi];
449 
450  if (newPatchi != -1)
451  {
452  labelList newToOld
453  (
454  calcPatchMap
455  (
456  oldPatchStarts[patchi],
457  oldPatchSizes[patchi],
458  meshMap.oldFaceMap(),
459  mesh.boundaryMesh()[newPatchi],
460  -1 // unmapped value
461  )
462  );
463 
464  directFvPatchFieldMapper patchMapper(newToOld);
465 
466  // Create new patchField with same type as existing one.
467  // Note:
468  // - boundaryField already in new order so access with newPatchi
469  // - bfld[newPatchi] both used for type and old
470  // value
471  // - hope that field mapping allows aliasing since old and new
472  // are same memory!
473  bfld.set
474  (
475  newPatchi,
477  (
478  bfld[newPatchi], // old field
479  mesh.boundary()[newPatchi], // new fvPatch
480  fld(), // new internal field
481  patchMapper // mapper (new to old)
482  )
483  );
484  }
485  }
486  }
487 
488 
489 
490  // Patch fields from added mesh
491  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
492 
493  {
494  const labelList& addedPatchMap = meshMap.addedPatchMap();
495 
496  // Add addedMesh patches
497  forAll(addedPatchMap, patchi)
498  {
499  label newPatchi = addedPatchMap[patchi];
500 
501  if (newPatchi != -1)
502  {
503  const polyPatch& newPatch = mesh.boundaryMesh()[newPatchi];
504  const polyPatch& oldPatch =
505  fldToAdd.mesh().boundaryMesh()[patchi];
506 
507  if (!bfld(newPatchi))
508  {
509  // First occurrence of newPatchi. Map from existing
510  // patchField
511 
512  // From new patch faces to patch faces on added mesh.
513  labelList newToAdded
514  (
515  calcPatchMap
516  (
517  oldPatch.start(),
518  oldPatch.size(),
519  meshMap.addedFaceMap(),
520  newPatch,
521  -1 // unmapped values
522  )
523  );
524 
525  directFvPatchFieldMapper patchMapper(newToAdded);
526 
527  bfld.set
528  (
529  newPatchi,
531  (
532  fldToAdd.boundaryField()[patchi],// added field
533  mesh.boundary()[newPatchi], // new fvPatch
534  fld(), // new int. field
535  patchMapper // mapper
536  )
537  );
538  }
539  else
540  {
541  // PatchField will have correct size already. Just slot in
542  // my elements.
543 
544  labelList addedToNew(oldPatch.size(), -1);
545  forAll(addedToNew, i)
546  {
547  label addedFacei = oldPatch.start()+i;
548  label newFacei = meshMap.addedFaceMap()[addedFacei];
549  label patchFacei = newFacei-newPatch.start();
550  if (patchFacei >= 0 && patchFacei < newPatch.size())
551  {
552  addedToNew[i] = patchFacei;
553  }
554  }
555 
556  bfld[newPatchi].rmap
557  (
558  fldToAdd.boundaryField()[patchi],
559  addedToNew
560  );
561  }
562  }
563  }
564  }
565 }
566 
567 
568 template<class Type>
570 (
571  const mapAddedPolyMesh& meshMap,
572  const fvMesh& mesh,
573  const fvMesh& meshToAdd
574 )
575 {
577 
579  (
580  mesh.objectRegistry::lookupClass<fldType>()
581  );
582 
583  HashTable<const fldType*> fieldsToAdd
584  (
585  meshToAdd.objectRegistry::lookupClass<fldType>()
586  );
587 
588  // It is necessary to enforce that all old-time fields are stored
589  // before the mapping is performed. Otherwise, if the
590  // old-time-level field is mapped before the field itself, sizes
591  // will not match.
592 
593  for
594  (
595  typename HashTable<const fldType*>::
596  iterator fieldIter = fields.begin();
597  fieldIter != fields.end();
598  ++fieldIter
599  )
600  {
601  if (debug)
602  {
603  Pout<< "MapSurfaceFields : Storing old time for "
604  << fieldIter()->name() << endl;
605  }
606 
607  const_cast<fldType*>(fieldIter())->storeOldTimes();
608  }
609 
610 
611  for
612  (
613  typename HashTable<const fldType*>::
614  iterator fieldIter = fields.begin();
615  fieldIter != fields.end();
616  ++fieldIter
617  )
618  {
619  fldType& fld = const_cast<fldType&>(*fieldIter());
620 
621  if (fieldsToAdd.found(fld.name()))
622  {
623  const fldType& fldToAdd = *fieldsToAdd[fld.name()];
624 
625  if (debug)
626  {
627  Pout<< "MapSurfaceFields : mapping " << fld.name()
628  << " and " << fldToAdd.name() << endl;
629  }
630 
631  MapSurfaceField<Type>(meshMap, fld, fldToAdd);
632  }
633  else
634  {
636  << "Not mapping field " << fld.name()
637  << " since not present on mesh to add"
638  << endl;
639  }
640  }
641 }
642 
643 
644 template<class Type>
645 void Foam::fvMeshAdder::MapPointField
646 (
647  const pointMesh& mesh,
648  const mapAddedPolyMesh& meshMap,
649  const labelListList& oldMeshPoints,
650 
653 )
654 {
655  // This is a bit tricky:
656  // - mesh pointed to by fld is invalid
657  // - pointPatches pointed to be fld are invalid
658 
660  Boundary& bfld = fld.boundaryFieldRef();
661 
662  // Internal field
663  // ~~~~~~~~~~~~~~
664 
665  // Store old internal field
666  {
667  Field<Type> oldField(fld);
668 
669  // Modify internal field
670  Field<Type>& intFld = fld.primitiveFieldRef();
671 
672  intFld.setSize(mesh.size());
673 
674  intFld.rmap(oldField, meshMap.oldPointMap());
675  intFld.rmap(fldToAdd, meshMap.addedPointMap());
676  }
677 
678 
679  // Patch fields from old mesh
680  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
681 
682  {
683  const labelList& oldPatchMap = meshMap.oldPatchMap();
684 
685  // Reorder old patches in order of new ones. Put removed patches at end.
686 
687  label unusedPatchi = 0;
688 
689  forAll(oldPatchMap, patchi)
690  {
691  label newPatchi = oldPatchMap[patchi];
692 
693  if (newPatchi != -1)
694  {
695  unusedPatchi++;
696  }
697  }
698 
699  label nUsedPatches = unusedPatchi;
700 
701  // Reorder list for patchFields
702  labelList oldToNew(oldPatchMap.size());
703 
704  forAll(oldPatchMap, patchi)
705  {
706  label newPatchi = oldPatchMap[patchi];
707 
708  if (newPatchi != -1)
709  {
710  oldToNew[patchi] = newPatchi;
711  }
712  else
713  {
714  oldToNew[patchi] = unusedPatchi++;
715  }
716  }
717 
718 
719  // Sort deleted ones last so is now in newPatch ordering
720  bfld.reorder(oldToNew);
721  // Extend to covers all patches
722  bfld.setSize(mesh.boundary().size());
723  // Delete unused patches
724  for
725  (
726  label newPatchi = nUsedPatches;
727  newPatchi < bfld.size();
728  newPatchi++
729  )
730  {
731  bfld.set(newPatchi, nullptr);
732  }
733 
734 
735  // Map old values
736  // ~~~~~~~~~~~~~~
737 
738  forAll(oldPatchMap, patchi)
739  {
740  label newPatchi = oldPatchMap[patchi];
741 
742  if (newPatchi != -1)
743  {
744  const labelList& oldMp = oldMeshPoints[patchi];
745  const pointPatch& newPp = mesh.boundary()[newPatchi];
746  const labelList& newMeshPoints = newPp.meshPoints();
747  Map<label> newMeshPointMap(2*newMeshPoints.size());
748  forAll(newMeshPoints, ppi)
749  {
750  newMeshPointMap.insert(newMeshPoints[ppi], ppi);
751  }
752 
753  labelList newToOld(newPp.size(), -1);
754  forAll(oldMp, oldPointi)
755  {
756  label newPointi = oldMp[oldPointi];
757 
759  newMeshPointMap.find(newPointi);
760  if (fnd == newMeshPointMap.end())
761  {
762  // Possibly an internal point
763  }
764  else
765  {
766  // Part of new patch
767  newToOld[fnd()] = oldPointi;
768  }
769  }
770 
771  directPointPatchFieldMapper patchMapper(newToOld);
772 
773  // Create new patchField with same type as existing one.
774  // Note:
775  // - boundaryField already in new order so access with newPatchi
776  // - bfld[newPatchi] both used for type and old
777  // value
778  // - hope that field mapping allows aliasing since old and new
779  // are same memory!
780  bfld.set
781  (
782  newPatchi,
784  (
785  bfld[newPatchi], // old field
786  mesh.boundary()[newPatchi], // new pointPatch
787  fld(), // new internal field
788  patchMapper // mapper (new to old)
789  )
790  );
791  }
792  }
793  }
794 
795 
796 
797  // Patch fields from added mesh
798  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
799 
800  {
801  const labelList& addedPatchMap = meshMap.addedPatchMap();
802 
803  // Add addedMesh patches
804  forAll(addedPatchMap, patchi)
805  {
806  label newPatchi = addedPatchMap[patchi];
807 
808  if (newPatchi != -1)
809  {
810  const pointPatch& oldPatch = fldToAdd.mesh().boundary()[patchi];
811  const labelList& oldMp = oldPatch.meshPoints();
812  const pointPatch& newPp = mesh.boundary()[newPatchi];
813  const labelList& newMeshPoints = newPp.meshPoints();
814  Map<label> newMpm(2*newMeshPoints.size());
815  forAll(newMeshPoints, ppi)
816  {
817  newMpm.insert(newMeshPoints[ppi], ppi);
818  }
819 
820  if (!bfld(newPatchi))
821  {
822  // First occurrence of newPatchi. Map from existing
823  // patchField
824 
825  labelList newToAdded(newPp.size(), -1);
826  forAll(oldMp, oldPointi)
827  {
828  label newPointi = oldMp[oldPointi];
829  Map<label>::const_iterator fnd = newMpm.find(newPointi);
830  if (fnd == newMpm.end())
831  {
832  // Possibly an internal point
833  }
834  else
835  {
836  // Part of new patch
837  newToAdded[fnd()] = oldPointi;
838  }
839  }
840 
841  directPointPatchFieldMapper patchMapper(newToAdded);
842 
843  bfld.set
844  (
845  newPatchi,
847  (
848  fldToAdd.boundaryField()[patchi],// added field
849  mesh.boundary()[newPatchi], // new pointPatch
850  fld(), // new int. field
851  patchMapper // mapper
852  )
853  );
854  }
855  else
856  {
857  // PatchField will have correct size already. Just slot in
858  // my elements.
859 
860  labelList oldToNew(oldPatch.size(), -1);
861  forAll(oldMp, oldPointi)
862  {
863  label newPointi = oldMp[oldPointi];
864  Map<label>::const_iterator fnd = newMpm.find(newPointi);
865  if (fnd != newMpm.end())
866  {
867  // Part of new patch
868  oldToNew[oldPointi] = fnd();
869  }
870  }
871 
872  bfld[newPatchi].rmap
873  (
874  fldToAdd.boundaryField()[patchi],
875  oldToNew
876  );
877  }
878  }
879  }
880  }
881 }
882 
883 
884 template<class Type>
886 (
887  const mapAddedPolyMesh& meshMap,
888  const pointMesh& mesh,
889  const labelListList& oldMeshPoints,
890  const objectRegistry& meshToAdd
891 )
892 {
894 
896  HashTable<const fldType*> fieldsToAdd(meshToAdd.lookupClass<fldType>());
897 
898  // It is necessary to enforce that all old-time fields are stored
899  // before the mapping is performed. Otherwise, if the
900  // old-time-level field is mapped before the field itself, sizes
901  // will not match.
902 
903  for
904  (
905  typename HashTable<const fldType*>::
906  iterator fieldIter = fields.begin();
907  fieldIter != fields.end();
908  ++fieldIter
909  )
910  {
911  if (debug)
912  {
913  Pout<< "MapPointFields : Storing old time for "
914  << fieldIter()->name() << endl;
915  }
916 
917  const_cast<fldType*>(fieldIter())->storeOldTimes();
918  }
919 
920 
921  for
922  (
923  typename HashTable<const fldType*>::
924  iterator fieldIter = fields.begin();
925  fieldIter != fields.end();
926  ++fieldIter
927  )
928  {
929  fldType& fld = const_cast<fldType&>(*fieldIter());
930 
931  if (fieldsToAdd.found(fld.name()))
932  {
933  const fldType& fldToAdd = *fieldsToAdd[fld.name()];
934 
935  if (debug)
936  {
937  Pout<< "MapPointFields : mapping " << fld.name()
938  << " and " << fldToAdd.name() << endl;
939  }
940 
941  MapPointField<Type>(mesh, meshMap, oldMeshPoints, fld, fldToAdd);
942  }
943  else
944  {
946  << "Not mapping field " << fld.name()
947  << " since not present on mesh to add"
948  << endl;
949  }
950  }
951 }
952 
953 
954 template<class Type>
955 void Foam::fvMeshAdder::MapDimField
956 (
957  const mapAddedPolyMesh& meshMap,
958 
960  const DimensionedField<Type, volMesh>& fldToAdd
961 )
962 {
963  const fvMesh& mesh = fld.mesh();
964 
965  // Store old field
966  Field<Type> oldField(fld);
967 
968  fld.setSize(mesh.nCells());
969 
970  fld.rmap(oldField, meshMap.oldCellMap());
971  fld.rmap(fldToAdd, meshMap.addedCellMap());
972 }
973 
974 
975 template<class Type>
977 (
978  const mapAddedPolyMesh& meshMap,
979  const fvMesh& mesh,
980  const fvMesh& meshToAdd
981 )
982 {
983  typedef DimensionedField<Type, volMesh> fldType;
984 
985  // Note: use strict flag on lookupClass to avoid picking up
986  // volFields
988  (
989  mesh.objectRegistry::lookupClass<fldType>(true)
990  );
991 
992  HashTable<const fldType*> fieldsToAdd
993  (
994  meshToAdd.objectRegistry::lookupClass<fldType>(true)
995  );
996 
997  for
998  (
999  typename HashTable<const fldType*>::
1000  iterator fieldIter = fields.begin();
1001  fieldIter != fields.end();
1002  ++fieldIter
1003  )
1004  {
1005  fldType& fld = const_cast<fldType&>(*fieldIter());
1006 
1007  if (fieldsToAdd.found(fld.name()))
1008  {
1009  const fldType& fldToAdd = *fieldsToAdd[fld.name()];
1010 
1011  if (debug)
1012  {
1013  Pout<< "MapDimFields : mapping " << fld.name()
1014  << " and " << fldToAdd.name() << endl;
1015  }
1016 
1017  MapDimField<Type>(meshMap, fld, fldToAdd);
1018  }
1019  else
1020  {
1021  WarningIn("fvMeshAdder::MapDimFields(..)")
1022  << "Not mapping field " << fld.name()
1023  << " since not present on mesh to add"
1024  << endl;
1025  }
1026  }
1027 }
1028 
1029 
1030 // ************************************************************************* //
const labelList & oldPointMap() const
From old mesh point/face/cell to new mesh point/face/cell.
Foam::surfaceFields.
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:434
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:82
#define WarningIn(functionName)
Report a warning using Foam::Warning.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const labelList & oldPatchSizes() const
Return list of the old patch sizes.
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
const word & name() const
Return name.
Definition: IOobject.H:303
static iteratorEnd end()
iteratorEnd set to beyond the end of any HashTable
Definition: HashTable.H:112
Map point field on topology change. This is a partial template specialisation for GeoMesh=pointMesh...
label nInternalFaces() const
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
label nCells() const
const labelList & oldCellMap() const
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const labelList & addedPatchMap() const
From added mesh patch index to new patch index or -1 if.
const labelList & addedCellMap() const
Class containing mesh-to-mesh mapping information after a mesh addition where we add a mesh (&#39;added m...
const labelList & oldPatchStarts() const
Return list of the old patch start labels.
Generic GeometricField class.
Info<< "Calculating turbulent flame speed field St\"<< endl;volScalarField St(IOobject("St", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:228
const labelList & addedPointMap() const
From added mesh point/face/cell to new mesh point/face/cell.
Mesh representing a set of points created from polyMesh.
Definition: pointMesh.H:48
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
Abstract base class for point-mesh patch fields.
dynamicFvMesh & mesh
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< ' ';}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< ' ';}gmvFile<< nl;forAll(lagrangianScalarNames, i){ const word &name=lagrangianScalarNames[i];IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:113
Pre-declare SubField and related Field type.
Definition: Field.H:56
word name() const
Return file name (part beyond last /)
Definition: fileName.C:195
const labelList & oldPatchMap() const
From old patch index to new patch index or -1 if patch.
const objectRegistry & thisDb() const
Return database. For now is its polyMesh.
Definition: pointMesh.H:115
const labelList & oldFaceMap() const
static void MapPointFields(const mapAddedPolyMesh &, const pointMesh &mesh, const labelListList &oldMeshPoints, const objectRegistry &meshToAdd)
Map all surfaceFields of Type.
static tmp< fvPatchField< Type > > New(const word &, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
List< label > labelList
A List of labels.
Definition: labelList.H:56
iterator begin()
Iterator set to the beginning of the HashTable.
Definition: HashTableI.H:411
An STL-conforming hash table.
Definition: HashTable.H:61
label size() const
Return number of points.
Definition: pointMesh.H:91
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
const Mesh & mesh() const
Return mesh.
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
void setSize(const label)
Reset size of List.
Definition: List.C:281
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:362
label patchi
#define WarningInFunction
Report a warning using Foam::Warning.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
label newPointi
Definition: readKivaGrid.H:501
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:309
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static void MapSurfaceFields(const mapAddedPolyMesh &, const fvMesh &mesh, const fvMesh &meshToAdd)
Map all surfaceFields of Type.
static void MapDimFields(const mapAddedPolyMesh &, const fvMesh &mesh, const fvMesh &meshToAdd)
Map all DimensionedFields of Type.
static void MapVolFields(const mapAddedPolyMesh &, const fvMesh &mesh, const fvMesh &meshToAdd)
Map all volFields of Type.
Registry of regIOobjects.
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
const pointBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: pointMesh.H:103
HashTable< const Type * > lookupClass(const bool strict=false) const
Lookup and return all objects of the given Type.
const labelList & addedFaceMap() const
An abstract base class with a fat-interface to all derived classes covering all possible ways in whic...
Definition: fvsPatchField.H:65
virtual const labelList & meshPoints() const =0
Return mesh points.
virtual label size() const =0
Return size.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:540