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-2021 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  const 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  const 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  const 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  const 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  const 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  const label addedFacei = oldPatch.start() + i;
229  const label newFacei =
230  meshMap.addedFaceMap()[addedFacei];
231  const label patchFacei = newFacei-newPatch.start();
232 
233  if (patchFacei >= 0 && patchFacei < newPatch.size())
234  {
235  addedToNew[i] = patchFacei;
236  }
237  }
238 
239  bfld[newPatchi].rmap
240  (
241  fldToAdd.boundaryField()[patchi],
242  addedToNew
243  );
244  }
245  }
246  }
247  }
248 }
249 
250 
251 template<class Type>
253 (
254  const mapAddedPolyMesh& meshMap,
255  const fvMesh& mesh,
256  const fvMesh& meshToAdd
257 )
258 {
260  (
261  mesh.objectRegistry::lookupClass
263  ()
264  );
265 
267  (
268  meshToAdd.objectRegistry::lookupClass
270  ()
271  );
272 
273  // It is necessary to enforce that all old-time fields are stored
274  // before the mapping is performed. Otherwise, if the
275  // old-time-level field is mapped before the field itself, sizes
276  // will not match.
277 
278  for
279  (
281  iterator fieldIter = fields.begin();
282  fieldIter != fields.end();
283  ++fieldIter
284  )
285  {
286  if (debug)
287  {
288  Pout<< "MapVolFields : Storing old time for " << fieldIter()->name()
289  << endl;
290  }
291 
292  const_cast<GeometricField<Type, fvPatchField, volMesh>*>(fieldIter())
293  ->storeOldTimes();
294  }
295 
296 
297  for
298  (
300  iterator fieldIter = fields.begin();
301  fieldIter != fields.end();
302  ++fieldIter
303  )
304  {
307  (
308  *fieldIter()
309  );
310 
311  if (fieldsToAdd.found(fld.name()))
312  {
314  *fieldsToAdd[fld.name()];
315 
316  if (debug)
317  {
318  Pout<< "MapVolFields : mapping " << fld.name()
319  << " and " << fldToAdd.name() << endl;
320  }
321 
322  MapVolField<Type>(meshMap, fld, fldToAdd);
323  }
324  else
325  {
327  << "Not mapping field " << fld.name()
328  << " since not present on mesh to add"
329  << endl;
330  }
331  }
332 }
333 
334 
335 template<class Type>
336 void Foam::fvMeshAdder::MapSurfaceField
337 (
338  const mapAddedPolyMesh& meshMap,
339 
342 )
343 {
344  const fvMesh& mesh = fld.mesh();
345  const labelList& oldPatchStarts = meshMap.oldPatchStarts();
346 
348  Boundary& bfld = fld.boundaryFieldRef();
349 
350  // Internal field
351  // ~~~~~~~~~~~~~~
352 
353  // Store old internal field
354  {
355  Field<Type> oldField(fld);
356 
357  // Modify internal field
358  Field<Type>& intFld = fld.primitiveFieldRef();
359 
360  intFld.setSize(mesh.nInternalFaces());
361 
362  intFld.rmap(oldField, meshMap.oldFaceMap());
363  intFld.rmap(fldToAdd, meshMap.addedFaceMap());
364 
365 
366  // Faces that were boundary faces but are not anymore.
367  // Use owner value (so lowest numbered cell, i.e. from 'old' not 'added'
368  // mesh)
369  forAll(bfld, patchi)
370  {
371  const fvsPatchField<Type>& pf = bfld[patchi];
372 
373  const label start = oldPatchStarts[patchi];
374 
375  forAll(pf, i)
376  {
377  const label newFacei = meshMap.oldFaceMap()[start + i];
378 
379  if (newFacei >= 0 && newFacei < mesh.nInternalFaces())
380  {
381  intFld[newFacei] = pf[i];
382  }
383  }
384  }
385  }
386 
387 
388  // Patch fields from old mesh
389  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
390 
391  {
392  const labelList& oldPatchMap = meshMap.oldPatchMap();
393  const labelList& oldPatchSizes = meshMap.oldPatchSizes();
394 
395  // Reorder old patches in order of new ones. Put removed patches at end.
396 
397  label unusedPatchi = 0;
398 
399  forAll(oldPatchMap, patchi)
400  {
401  const label newPatchi = oldPatchMap[patchi];
402 
403  if (newPatchi != -1)
404  {
405  unusedPatchi++;
406  }
407  }
408 
409  label nUsedPatches = unusedPatchi;
410 
411  // Reorder list for patchFields
412  labelList oldToNew(oldPatchMap.size());
413 
414  forAll(oldPatchMap, patchi)
415  {
416  const label newPatchi = oldPatchMap[patchi];
417 
418  if (newPatchi != -1)
419  {
420  oldToNew[patchi] = newPatchi;
421  }
422  else
423  {
424  oldToNew[patchi] = unusedPatchi++;
425  }
426  }
427 
428 
429  // Sort deleted ones last so is now in newPatch ordering
430  bfld.reorder(oldToNew);
431  // Extend to covers all patches
432  bfld.setSize(mesh.boundaryMesh().size());
433  // Delete unused patches
434  for
435  (
436  label newPatchi = nUsedPatches;
437  newPatchi < bfld.size();
438  newPatchi++
439  )
440  {
441  bfld.set(newPatchi, nullptr);
442  }
443 
444 
445  // Map old values
446  // ~~~~~~~~~~~~~~
447 
448  forAll(oldPatchMap, patchi)
449  {
450  const label newPatchi = oldPatchMap[patchi];
451 
452  if (newPatchi != -1)
453  {
454  const labelList newToOld
455  (
456  calcPatchMap
457  (
458  oldPatchStarts[patchi],
459  oldPatchSizes[patchi],
460  meshMap.oldFaceMap(),
461  mesh.boundaryMesh()[newPatchi],
462  -1 // unmapped value
463  )
464  );
465 
466  directFvPatchFieldMapper patchMapper(newToOld);
467 
468  // Create new patchField with same type as existing one.
469  // Note:
470  // - boundaryField already in new order so access with newPatchi
471  // - bfld[newPatchi] both used for type and old
472  // value
473  // - hope that field mapping allows aliasing since old and new
474  // are same memory!
475  bfld.set
476  (
477  newPatchi,
479  (
480  bfld[newPatchi], // old field
481  mesh.boundary()[newPatchi], // new fvPatch
482  fld(), // new internal field
483  patchMapper // mapper (new to old)
484  )
485  );
486  }
487  }
488  }
489 
490 
491 
492  // Patch fields from added mesh
493  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
494 
495  {
496  const labelList& addedPatchMap = meshMap.addedPatchMap();
497 
498  // Add addedMesh patches
499  forAll(addedPatchMap, patchi)
500  {
501  const label newPatchi = addedPatchMap[patchi];
502 
503  if (newPatchi != -1)
504  {
505  const polyPatch& newPatch = mesh.boundaryMesh()[newPatchi];
506  const polyPatch& oldPatch =
507  fldToAdd.mesh().boundaryMesh()[patchi];
508 
509  if (!bfld(newPatchi))
510  {
511  // First occurrence of newPatchi. Map from existing
512  // patchField
513 
514  // From new patch faces to patch faces on added mesh.
515  const labelList newToAdded
516  (
517  calcPatchMap
518  (
519  oldPatch.start(),
520  oldPatch.size(),
521  meshMap.addedFaceMap(),
522  newPatch,
523  -1 // unmapped values
524  )
525  );
526 
527  directFvPatchFieldMapper patchMapper(newToAdded);
528 
529  bfld.set
530  (
531  newPatchi,
533  (
534  fldToAdd.boundaryField()[patchi],// added field
535  mesh.boundary()[newPatchi], // new fvPatch
536  fld(), // new int. field
537  patchMapper // mapper
538  )
539  );
540  }
541  else
542  {
543  // PatchField will have correct size already. Just slot in
544  // my elements.
545 
546  labelList addedToNew(oldPatch.size(), -1);
547  forAll(addedToNew, i)
548  {
549  const label addedFacei = oldPatch.start() + i;
550  const label newFacei =
551  meshMap.addedFaceMap()[addedFacei];
552  const label patchFacei = newFacei-newPatch.start();
553 
554  if (patchFacei >= 0 && patchFacei < newPatch.size())
555  {
556  addedToNew[i] = patchFacei;
557  }
558  }
559 
560  bfld[newPatchi].rmap
561  (
562  fldToAdd.boundaryField()[patchi],
563  addedToNew
564  );
565  }
566  }
567  }
568  }
569 }
570 
571 
572 template<class Type>
574 (
575  const mapAddedPolyMesh& meshMap,
576  const fvMesh& mesh,
577  const fvMesh& meshToAdd
578 )
579 {
581 
583  (
584  mesh.objectRegistry::lookupClass<fldType>()
585  );
586 
587  HashTable<const fldType*> fieldsToAdd
588  (
589  meshToAdd.objectRegistry::lookupClass<fldType>()
590  );
591 
592  // It is necessary to enforce that all old-time fields are stored
593  // before the mapping is performed. Otherwise, if the
594  // old-time-level field is mapped before the field itself, sizes
595  // will not match.
596 
597  for
598  (
599  typename HashTable<const fldType*>::
600  iterator fieldIter = fields.begin();
601  fieldIter != fields.end();
602  ++fieldIter
603  )
604  {
605  if (debug)
606  {
607  Pout<< "MapSurfaceFields : Storing old time for "
608  << fieldIter()->name() << endl;
609  }
610 
611  const_cast<fldType*>(fieldIter())->storeOldTimes();
612  }
613 
614 
615  for
616  (
617  typename HashTable<const fldType*>::
618  iterator fieldIter = fields.begin();
619  fieldIter != fields.end();
620  ++fieldIter
621  )
622  {
623  fldType& fld = const_cast<fldType&>(*fieldIter());
624 
625  if (fieldsToAdd.found(fld.name()))
626  {
627  const fldType& fldToAdd = *fieldsToAdd[fld.name()];
628 
629  if (debug)
630  {
631  Pout<< "MapSurfaceFields : mapping " << fld.name()
632  << " and " << fldToAdd.name() << endl;
633  }
634 
635  MapSurfaceField<Type>(meshMap, fld, fldToAdd);
636  }
637  else
638  {
640  << "Not mapping field " << fld.name()
641  << " since not present on mesh to add"
642  << endl;
643  }
644  }
645 }
646 
647 
648 template<class Type>
649 void Foam::fvMeshAdder::MapPointField
650 (
651  const pointMesh& mesh,
652  const mapAddedPolyMesh& meshMap,
653  const labelListList& oldMeshPoints,
654 
657 )
658 {
659  // This is a bit tricky:
660  // - mesh pointed to by fld is invalid
661  // - pointPatches pointed to be fld are invalid
662 
664  Boundary& bfld = fld.boundaryFieldRef();
665 
666  // Internal field
667  // ~~~~~~~~~~~~~~
668 
669  // Store old internal field
670  {
671  Field<Type> oldField(fld);
672 
673  // Modify internal field
674  Field<Type>& intFld = fld.primitiveFieldRef();
675 
676  intFld.setSize(mesh.size());
677 
678  intFld.rmap(oldField, meshMap.oldPointMap());
679  intFld.rmap(fldToAdd, meshMap.addedPointMap());
680  }
681 
682 
683  // Patch fields from old mesh
684  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
685 
686  {
687  const labelList& oldPatchMap = meshMap.oldPatchMap();
688 
689  // Reorder old patches in order of new ones. Put removed patches at end.
690 
691  label unusedPatchi = 0;
692 
693  forAll(oldPatchMap, patchi)
694  {
695  const label newPatchi = oldPatchMap[patchi];
696 
697  if (newPatchi != -1)
698  {
699  unusedPatchi++;
700  }
701  }
702 
703  label nUsedPatches = unusedPatchi;
704 
705  // Reorder list for patchFields
706  labelList oldToNew(oldPatchMap.size());
707 
708  forAll(oldPatchMap, patchi)
709  {
710  const label newPatchi = oldPatchMap[patchi];
711 
712  if (newPatchi != -1)
713  {
714  oldToNew[patchi] = newPatchi;
715  }
716  else
717  {
718  oldToNew[patchi] = unusedPatchi++;
719  }
720  }
721 
722 
723  // Sort deleted ones last so is now in newPatch ordering
724  bfld.reorder(oldToNew);
725 
726  // Extend to covers all patches
727  bfld.setSize(mesh.boundary().size());
728 
729  // Delete unused patches
730  for
731  (
732  label newPatchi = nUsedPatches;
733  newPatchi < bfld.size();
734  newPatchi++
735  )
736  {
737  bfld.set(newPatchi, nullptr);
738  }
739 
740 
741  // Map old values
742  // ~~~~~~~~~~~~~~
743 
744  forAll(oldPatchMap, patchi)
745  {
746  const label newPatchi = oldPatchMap[patchi];
747 
748  if (newPatchi != -1)
749  {
750  const labelList& oldMp = oldMeshPoints[patchi];
751  const pointPatch& newPp = mesh.boundary()[newPatchi];
752  const labelList& newMeshPoints = newPp.meshPoints();
753  const labelList& oldPointMap = meshMap.oldPointMap();
754 
755  Map<label> newMeshPointMap(2*newMeshPoints.size());
756  forAll(newMeshPoints, ppi)
757  {
758  newMeshPointMap.insert(newMeshPoints[ppi], ppi);
759  }
760 
761  labelList newToOld(newPp.size(), -1);
762  forAll(oldMp, oldPointi)
763  {
764  const label newPointi = oldPointMap[oldMp[oldPointi]];
765 
767  newMeshPointMap.find(newPointi);
768 
769  if (fnd == newMeshPointMap.end())
770  {
771  // Possibly an internal point
772  }
773  else
774  {
775  // Part of new patch
776  newToOld[fnd()] = oldPointi;
777  }
778  }
779 
780  directPointPatchFieldMapper patchMapper(newToOld);
781 
782  // Create new patchField with same type as existing one.
783  // Note:
784  // - boundaryField already in new order so access with newPatchi
785  // - bfld[newPatchi] both used for type and old
786  // value
787  // - hope that field mapping allows aliasing since old and new
788  // are same memory!
789  bfld.set
790  (
791  newPatchi,
793  (
794  bfld[newPatchi], // old field
795  mesh.boundary()[newPatchi], // new pointPatch
796  fld(), // new internal field
797  patchMapper // mapper (new to old)
798  )
799  );
800  }
801  }
802  }
803 
804 
805 
806  // Patch fields from added mesh
807  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
808 
809  {
810  const labelList& addedPatchMap = meshMap.addedPatchMap();
811 
812  // Add addedMesh patches
813  forAll(addedPatchMap, patchi)
814  {
815  const label newPatchi = addedPatchMap[patchi];
816 
817  if (newPatchi != -1)
818  {
819  const pointPatch& oldPatch = fldToAdd.mesh().boundary()[patchi];
820  const labelList& oldMp = oldPatch.meshPoints();
821  const pointPatch& newPp = mesh.boundary()[newPatchi];
822  const labelList& newMeshPoints = newPp.meshPoints();
823  const labelList& addedPointMap = meshMap.addedPointMap();
824 
825  Map<label> newMpm(2*newMeshPoints.size());
826  forAll(newMeshPoints, ppi)
827  {
828  newMpm.insert(newMeshPoints[ppi], ppi);
829  }
830 
831  if (!bfld(newPatchi))
832  {
833  // First occurrence of newPatchi. Map from existing
834  // patchField
835 
836  labelList newToAdded(newPp.size(), -1);
837  forAll(oldMp, oldPointi)
838  {
839  const label newPointi = addedPointMap[oldMp[oldPointi]];
840 
841  Map<label>::const_iterator fnd = newMpm.find(newPointi);
842  if (fnd == newMpm.end())
843  {
844  // Possibly an internal point
845  }
846  else
847  {
848  // Part of new patch
849  newToAdded[fnd()] = oldPointi;
850  }
851  }
852 
853  directPointPatchFieldMapper patchMapper(newToAdded);
854 
855  bfld.set
856  (
857  newPatchi,
859  (
860  fldToAdd.boundaryField()[patchi],// added field
861  mesh.boundary()[newPatchi], // new pointPatch
862  fld(), // new int. field
863  patchMapper // mapper
864  )
865  );
866  }
867  else
868  {
869  // PatchField will have correct size already. Just slot in
870  // my elements.
871 
872  labelList oldToNew(oldPatch.size(), -1);
873  forAll(oldMp, oldPointi)
874  {
875  const label newPointi = addedPointMap[oldMp[oldPointi]];
876 
877  Map<label>::const_iterator fnd = newMpm.find(newPointi);
878  if (fnd != newMpm.end())
879  {
880  // Part of new patch
881  oldToNew[oldPointi] = fnd();
882  }
883  }
884 
885  bfld[newPatchi].rmap
886  (
887  fldToAdd.boundaryField()[patchi],
888  oldToNew
889  );
890  }
891  }
892  }
893  }
894 }
895 
896 
897 template<class Type>
899 (
900  const mapAddedPolyMesh& meshMap,
901  const pointMesh& mesh,
902  const labelListList& oldMeshPoints,
903  const objectRegistry& meshToAdd
904 )
905 {
907 
909  HashTable<const fldType*> fieldsToAdd(meshToAdd.lookupClass<fldType>());
910 
911  // It is necessary to enforce that all old-time fields are stored
912  // before the mapping is performed. Otherwise, if the
913  // old-time-level field is mapped before the field itself, sizes
914  // will not match.
915 
916  for
917  (
918  typename HashTable<const fldType*>::
919  iterator fieldIter = fields.begin();
920  fieldIter != fields.end();
921  ++fieldIter
922  )
923  {
924  if (debug)
925  {
926  Pout<< "MapPointFields : Storing old time for "
927  << fieldIter()->name() << endl;
928  }
929 
930  const_cast<fldType*>(fieldIter())->storeOldTimes();
931  }
932 
933 
934  for
935  (
936  typename HashTable<const fldType*>::
937  iterator fieldIter = fields.begin();
938  fieldIter != fields.end();
939  ++fieldIter
940  )
941  {
942  fldType& fld = const_cast<fldType&>(*fieldIter());
943 
944  if (fieldsToAdd.found(fld.name()))
945  {
946  const fldType& fldToAdd = *fieldsToAdd[fld.name()];
947 
948  if (debug)
949  {
950  Pout<< "MapPointFields : mapping " << fld.name()
951  << " and " << fldToAdd.name() << endl;
952  }
953 
954  MapPointField<Type>(mesh, meshMap, oldMeshPoints, fld, fldToAdd);
955  }
956  else
957  {
959  << "Not mapping field " << fld.name()
960  << " since not present on mesh to add"
961  << endl;
962  }
963  }
964 }
965 
966 
967 template<class Type>
968 void Foam::fvMeshAdder::MapDimField
969 (
970  const mapAddedPolyMesh& meshMap,
971 
973  const DimensionedField<Type, volMesh>& fldToAdd
974 )
975 {
976  const fvMesh& mesh = fld.mesh();
977 
978  // Store old field
979  Field<Type> oldField(fld);
980 
981  fld.setSize(mesh.nCells());
982 
983  fld.rmap(oldField, meshMap.oldCellMap());
984  fld.rmap(fldToAdd, meshMap.addedCellMap());
985 }
986 
987 
988 template<class Type>
990 (
991  const mapAddedPolyMesh& meshMap,
992  const fvMesh& mesh,
993  const fvMesh& meshToAdd
994 )
995 {
996  typedef DimensionedField<Type, volMesh> fldType;
997 
998  // Note: use strict flag on lookupClass to avoid picking up
999  // volFields
1001  (
1002  mesh.objectRegistry::lookupClass<fldType>(true)
1003  );
1004 
1005  HashTable<const fldType*> fieldsToAdd
1006  (
1007  meshToAdd.objectRegistry::lookupClass<fldType>(true)
1008  );
1009 
1010  for
1011  (
1012  typename HashTable<const fldType*>::
1013  iterator fieldIter = fields.begin();
1014  fieldIter != fields.end();
1015  ++fieldIter
1016  )
1017  {
1018  fldType& fld = const_cast<fldType&>(*fieldIter());
1019 
1020  if (fieldsToAdd.found(fld.name()))
1021  {
1022  const fldType& fldToAdd = *fieldsToAdd[fld.name()];
1023 
1024  if (debug)
1025  {
1026  Pout<< "MapDimFields : mapping " << fld.name()
1027  << " and " << fldToAdd.name() << endl;
1028  }
1029 
1030  MapDimField<Type>(meshMap, fld, fldToAdd);
1031  }
1032  else
1033  {
1034  WarningIn("fvMeshAdder::MapDimFields(..)")
1035  << "Not mapping field " << fld.name()
1036  << " since not present on mesh to add"
1037  << endl;
1038  }
1039  }
1040 }
1041 
1042 
1043 // ************************************************************************* //
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:453
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.
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const word & name() const
Return name.
Definition: IOobject.H:315
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.
const labelList & addedPointMap() const
From added mesh point/face/cell to new mesh point/face/cell.
fvMesh & mesh
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.
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:109
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:85
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
Info<< "Reading field p_rgh\"<< endl;volScalarField p_rgh(IOobject("p_rgh", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, p_rgh, pimple.dict(), thermo.incompressible());mesh.schemes().setFluxRequired(p_rgh.name());hydrostaticInitialisation(p_rgh, p, rho, U, gh, ghf, pRef, thermo, pimple.dict());Info<< "Creating field dpdt\"<< endl;volScalarField dpdt(IOobject("dpdt", runTime.timeName(), mesh), mesh, dimensionedScalar(p.dimensions()/dimTime, 0));Info<< "Creating field kinetic energy K\"<< endl;volScalarField K("K", 0.5 *magSqr(U));dimensionedScalar initialMass=fvc::domainIntegrate(rho);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:131
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:95
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:306
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:97
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:800