fvMeshAdderTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2015-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 #include "emptyFvPatchField.H"
33 
34 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
35 
36 template<class Type>
37 void Foam::fvMeshAdder::MapVolField
38 (
39  const mapAddedPolyMesh& meshMap,
40 
41  GeometricField<Type, fvPatchField, volMesh>& fld,
42  const GeometricField<Type, fvPatchField, volMesh>& fldToAdd,
43  const bool fullyMapped
44 )
45 {
46  const fvMesh& mesh = fld.mesh();
47 
48  // Internal field
49  // ~~~~~~~~~~~~~~
50 
51  {
52  // Store old internal field
53  const Field<Type> oldInternalField(fld.primitiveField());
54 
55  // Modify internal field
56  Field<Type>& intFld = fld.primitiveFieldRef();
57 
58  intFld.setSize(mesh.nCells());
59 
60  intFld.rmap(oldInternalField, meshMap.oldCellMap());
61  intFld.rmap(fldToAdd.primitiveField(), meshMap.addedCellMap());
62  }
63 
64 
65  // Patch fields from old mesh
66  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
67 
68  auto& bfld = fld.boundaryFieldRef();
69 
70  {
71  const labelList& oldPatchMap = meshMap.oldPatchMap();
72  const labelList& oldPatchStarts = meshMap.oldPatchStarts();
73  const labelList& oldPatchSizes = meshMap.oldPatchSizes();
74 
75  // Reorder old patches in order of new ones. Put removed patches at end.
76 
77  label unusedPatchi = 0;
78 
79  forAll(oldPatchMap, patchi)
80  {
81  label newPatchi = oldPatchMap[patchi];
82 
83  if (newPatchi != -1)
84  {
85  unusedPatchi++;
86  }
87  }
88 
89  label nUsedPatches = unusedPatchi;
90 
91  // Reorder list for patchFields
92  labelList oldToNew(oldPatchMap.size());
93 
94  forAll(oldPatchMap, patchi)
95  {
96  const label newPatchi = oldPatchMap[patchi];
97 
98  if (newPatchi != -1)
99  {
100  oldToNew[patchi] = newPatchi;
101  }
102  else
103  {
104  oldToNew[patchi] = unusedPatchi++;
105  }
106  }
107 
108 
109  // Sort deleted ones last so is now in newPatch ordering
110  bfld.reorder(oldToNew);
111  // Extend to covers all patches
112  bfld.resize(mesh.boundaryMesh().size());
113  // Delete unused patches
114  for
115  (
116  label newPatchi = nUsedPatches;
117  newPatchi < bfld.size();
118  newPatchi++
119  )
120  {
121  bfld.release(newPatchi);
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  // Override mapping (for use in e.g. fvMeshDistribute where
149  // it sorts mapping out itself)
150  if (fullyMapped)
151  {
152  patchMapper.hasUnmapped() = false;
153  }
154 
155  // Create new patchField with same type as existing one.
156  // Note:
157  // - boundaryField already in new order so access with newPatchi
158  // - fld.boundaryField()[newPatchi] both used for type and old
159  // value
160  // - hope that field mapping allows aliasing since old and new
161  // are same memory!
162  bfld.set
163  (
164  newPatchi,
166  (
167  bfld[newPatchi], // old field
168  mesh.boundary()[newPatchi], // new fvPatch
169  fld(), // new internal field
170  patchMapper // mapper (new to old)
171  )
172  );
173  }
174  }
175  }
176 
177 
178 
179  // Patch fields from added mesh
180  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
181 
182  {
183  const labelList& addedPatchMap = meshMap.addedPatchMap();
184 
185  // Add addedMesh patches
186  forAll(addedPatchMap, patchi)
187  {
188  const label newPatchi = addedPatchMap[patchi];
189 
190  if (newPatchi != -1)
191  {
192  const polyPatch& newPatch = mesh.boundaryMesh()[newPatchi];
193  const polyPatch& oldPatch =
194  fldToAdd.mesh().boundaryMesh()[patchi];
195 
196  if (!bfld.set(newPatchi))
197  {
198  // First occurrence of newPatchi. Map from existing
199  // patchField
200 
201  // From new patch faces to patch faces on added mesh.
202  labelList newToAdded
203  (
204  calcPatchMap
205  (
206  oldPatch.start(),
207  oldPatch.size(),
208  meshMap.addedFaceMap(),
209  newPatch,
210  -1 // unmapped values
211  )
212  );
213 
214  directFvPatchFieldMapper patchMapper(newToAdded);
215 
216  // Override mapping (for use in e.g. fvMeshDistribute where
217  // it sorts mapping out itself)
218  if (fullyMapped)
219  {
220  patchMapper.hasUnmapped() = false;
221  }
222 
223  bfld.set
224  (
225  newPatchi,
227  (
228  fldToAdd.boundaryField()[patchi], // added field
229  mesh.boundary()[newPatchi], // new fvPatch
230  fld(), // new int. field
231  patchMapper // mapper
232  )
233  );
234  }
235  else
236  {
237  // PatchField will have correct size already. Just slot in
238  // my elements.
239 
240  labelList addedToNew(oldPatch.size(), -1);
241  forAll(addedToNew, i)
242  {
243  label addedFacei = oldPatch.start()+i;
244  label newFacei = meshMap.addedFaceMap()[addedFacei];
245  label patchFacei = newFacei-newPatch.start();
246  if (patchFacei >= 0 && patchFacei < newPatch.size())
247  {
248  addedToNew[i] = patchFacei;
249  }
250  }
251 
252  bfld[newPatchi].rmap
253  (
254  fldToAdd.boundaryField()[patchi],
255  addedToNew
256  );
257  }
258  }
259  }
260  }
261 }
262 
263 
264 template<class Type>
266 (
267  const mapAddedPolyMesh& meshMap,
268  const fvMesh& mesh,
269  const fvMesh& meshToAdd,
270  const bool fullyMapped
271 )
272 {
273  typedef GeometricField<Type, fvPatchField, volMesh> fldType;
274 
275  const UPtrList<const fldType> fields
276  (
277  mesh.objectRegistry::csorted<fldType>()
278  );
279 
280  // It is necessary to enforce that all old-time fields are stored
281  // before the mapping is performed. Otherwise, if the
282  // old-time-level field is mapped before the field itself, sizes
283  // will not match.
284 
285  for (const auto& field : fields)
286  {
287  DebugPout
288  << "MapVolFields : Storing old time for "
289  << field.name() << endl;
290 
291  const_cast<fldType&>(field).storeOldTimes();
292  }
293 
294 
295  for (const auto& field : fields)
296  {
297  fldType& fld = const_cast<fldType&>(field);
298 
299  const auto* toAdd =
300  meshToAdd.objectRegistry::cfindObject<fldType>(field.name());
301 
302  if (toAdd)
303  {
304  DebugPout
305  << "MapVolFields : mapping " << field.name() << endl;
306 
307  MapVolField<Type>(meshMap, fld, *toAdd, fullyMapped);
308  }
309  else
310  {
312  << "Not mapping field " << field.name()
313  << " - not present on mesh to add" << endl;
314  }
315  }
316 }
317 
318 
319 template<class Type>
320 void Foam::fvMeshAdder::MapSurfaceField
321 (
322  const mapAddedPolyMesh& meshMap,
323 
324  GeometricField<Type, fvsPatchField, surfaceMesh>& fld,
325  const GeometricField<Type, fvsPatchField, surfaceMesh>& fldToAdd,
326  const bool fullyMapped
327 )
328 {
329  const fvMesh& mesh = fld.mesh();
330  const labelList& oldPatchStarts = meshMap.oldPatchStarts();
331 
332  auto& bfld = fld.boundaryFieldRef();
333 
334  // Internal field
335  // ~~~~~~~~~~~~~~
336 
337  // Store old internal field
338  {
339  const Field<Type> oldField(fld);
340 
341  // Modify internal field
342  Field<Type>& intFld = fld.primitiveFieldRef();
343 
344  intFld.setSize(mesh.nInternalFaces());
345 
346  intFld.rmap(oldField, meshMap.oldFaceMap());
347  intFld.rmap(fldToAdd, meshMap.addedFaceMap());
348 
349 
350  // Faces that were boundary faces but are not anymore.
351  // Use owner value (so lowest numbered cell, i.e. from 'old' not 'added'
352  // mesh)
353  forAll(bfld, patchi)
354  {
355  const fvsPatchField<Type>& pf = bfld[patchi];
356 
357  label start = oldPatchStarts[patchi];
358 
359  forAll(pf, i)
360  {
361  label newFacei = meshMap.oldFaceMap()[start + i];
362 
363  if (newFacei >= 0 && newFacei < mesh.nInternalFaces())
364  {
365  intFld[newFacei] = pf[i];
366  }
367  }
368  }
369  }
370 
371 
372  // Patch fields from old mesh
373  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
374 
375  {
376  const labelList& oldPatchMap = meshMap.oldPatchMap();
377  const labelList& oldPatchSizes = meshMap.oldPatchSizes();
378 
379  // Reorder old patches in order of new ones. Put removed patches at end.
380 
381  label unusedPatchi = 0;
382 
383  forAll(oldPatchMap, patchi)
384  {
385  const label newPatchi = oldPatchMap[patchi];
386 
387  if (newPatchi != -1)
388  {
389  unusedPatchi++;
390  }
391  }
392 
393  label nUsedPatches = unusedPatchi;
394 
395  // Reorder list for patchFields
396  labelList oldToNew(oldPatchMap.size());
397 
398  forAll(oldPatchMap, patchi)
399  {
400  const label newPatchi = oldPatchMap[patchi];
401 
402  if (newPatchi != -1)
403  {
404  oldToNew[patchi] = newPatchi;
405  }
406  else
407  {
408  oldToNew[patchi] = unusedPatchi++;
409  }
410  }
411 
412 
413  // Sort deleted ones last so is now in newPatch ordering
414  bfld.reorder(oldToNew);
415  // Extend to covers all patches
416  bfld.resize(mesh.boundaryMesh().size());
417  // Delete unused patches
418  for
419  (
420  label newPatchi = nUsedPatches;
421  newPatchi < bfld.size();
422  newPatchi++
423  )
424  {
425  bfld.release(newPatchi);
426  }
427 
428 
429  // Map old values
430  // ~~~~~~~~~~~~~~
431 
432  forAll(oldPatchMap, patchi)
433  {
434  const label newPatchi = oldPatchMap[patchi];
435 
436  if (newPatchi != -1)
437  {
438  labelList newToOld
439  (
440  calcPatchMap
441  (
442  oldPatchStarts[patchi],
443  oldPatchSizes[patchi],
444  meshMap.oldFaceMap(),
445  mesh.boundaryMesh()[newPatchi],
446  -1 // unmapped value
447  )
448  );
449 
450  directFvPatchFieldMapper patchMapper(newToOld);
451 
452  // Override mapping (for use in e.g. fvMeshDistribute where
453  // it sorts mapping out itself)
454  if (fullyMapped)
455  {
456  patchMapper.hasUnmapped() = false;
457  }
458 
459  // Create new patchField with same type as existing one.
460  // Note:
461  // - boundaryField already in new order so access with newPatchi
462  // - bfld[newPatchi] both used for type and old
463  // value
464  // - hope that field mapping allows aliasing since old and new
465  // are same memory!
466  bfld.set
467  (
468  newPatchi,
470  (
471  bfld[newPatchi], // old field
472  mesh.boundary()[newPatchi], // new fvPatch
473  fld(), // new internal field
474  patchMapper // mapper (new to old)
475  )
476  );
477  }
478  }
479  }
480 
481 
482 
483  // Patch fields from added mesh
484  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
485 
486  {
487  const labelList& addedPatchMap = meshMap.addedPatchMap();
488 
489  // Add addedMesh patches
490  forAll(addedPatchMap, patchi)
491  {
492  const label newPatchi = addedPatchMap[patchi];
493 
494  if (newPatchi != -1)
495  {
496  const polyPatch& newPatch = mesh.boundaryMesh()[newPatchi];
497  const polyPatch& oldPatch =
498  fldToAdd.mesh().boundaryMesh()[patchi];
499 
500  if (!bfld.set(newPatchi))
501  {
502  // First occurrence of newPatchi. Map from existing
503  // patchField
504 
505  // From new patch faces to patch faces on added mesh.
506  labelList newToAdded
507  (
508  calcPatchMap
509  (
510  oldPatch.start(),
511  oldPatch.size(),
512  meshMap.addedFaceMap(),
513  newPatch,
514  -1 // unmapped values
515  )
516  );
517 
518  directFvPatchFieldMapper patchMapper(newToAdded);
519 
520  // Override mapping (for use in e.g. fvMeshDistribute where
521  // it sorts mapping out itself)
522  if (fullyMapped)
523  {
524  patchMapper.hasUnmapped() = false;
525  }
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  const bool fullyMapped
575 )
576 {
577  typedef GeometricField<Type, fvsPatchField, surfaceMesh> fldType;
578 
579  const UPtrList<const fldType> fields
580  (
581  mesh.objectRegistry::csorted<fldType>()
582  );
583 
584  // It is necessary to enforce that all old-time fields are stored
585  // before the mapping is performed. Otherwise, if the
586  // old-time-level field is mapped before the field itself, sizes
587  // will not match.
588 
589  for (const auto& field : fields)
590  {
591  DebugPout
592  << "MapSurfaceFields : Storing old time for "
593  << field.name() << endl;
594 
595  const_cast<fldType&>(field).storeOldTimes();
596  }
597 
598  for (const auto& field : fields)
599  {
600  fldType& fld = const_cast<fldType&>(field);
601 
602  const auto* toAdd =
603  meshToAdd.objectRegistry::cfindObject<fldType>(field.name());
604 
605  if (toAdd)
606  {
607  DebugPout
608  << "MapSurfaceFields : mapping " << field.name() << endl;
609 
610  MapSurfaceField<Type>(meshMap, fld, *toAdd, fullyMapped);
611  }
612  else
613  {
615  << "Not mapping field " << field.name()
616  << " - not present on mesh to add" << endl;
617  }
618  }
619 }
620 
621 
622 template<class Type>
623 void Foam::fvMeshAdder::MapDimField
624 (
625  const mapAddedPolyMesh& meshMap,
626 
627  DimensionedField<Type, volMesh>& fld,
628  const DimensionedField<Type, volMesh>& fldToAdd
629 )
630 {
631  const fvMesh& mesh = fld.mesh();
632 
633  // Store old field
634  const Field<Type> oldField(fld);
635 
636  fld.setSize(mesh.nCells());
637 
638  fld.rmap(oldField, meshMap.oldCellMap());
639  fld.rmap(fldToAdd, meshMap.addedCellMap());
640 }
641 
642 
643 template<class Type>
645 (
646  const mapAddedPolyMesh& meshMap,
647  const fvMesh& mesh,
648  const fvMesh& meshToAdd
649 )
650 {
651  typedef DimensionedField<Type, volMesh> fldType;
652 
653  // NB: strict=true to avoid picking up volFields
655  (
656  mesh.objectRegistry::csorted<fldType, true>()
657  );
658 
659  for (const auto& field : fields)
660  {
661  fldType& fld = const_cast<fldType&>(field);
662 
663  const auto* toAdd =
664  meshToAdd.objectRegistry::cfindObject<fldType>(field.name());
665 
666  // Apply strict check (to avoid picking volFields)
667  if (toAdd && Foam::isType<fldType>(*toAdd))
668  {
669  DebugPout
670  << "MapDimFields : mapping " << field.name() << endl;
671 
672  MapDimField<Type>(meshMap, fld, *toAdd);
673  }
674  else
675  {
677  << "Not mapping field " << field.name()
678  << " - not present on mesh to add" << endl;
679  }
680  }
681 }
682 
684 //- Multi-mesh mapping
685 
686 template<class Type>
687 void Foam::fvMeshAdder::MapDimField
688 (
690  const labelListList& cellProcAddressing,
691  const bool fullyMapped
692 )
693 {
694  // Add fields to fields[0] after adding the meshes to meshes[0].
695  // Mesh[0] is the sum of all meshes. Fields are not yet mapped.
696 
697  if (!flds.test(0) || cellProcAddressing.size() != flds.size())
698  {
700  << "Not valid field at element 0 in list of size "
701  << flds.size() << nl
702  << exit(FatalError);
703  }
704 
705 
706  // Internal field
707  // ~~~~~~~~~~~~~~
708 
709  {
710  // Store old internal field
711  const Field<Type> oldInternalField(flds[0]);
712 
713  // Modify internal field
714  Field<Type>& intFld = flds[0];
715 
716  // Set to new mesh size
717  intFld.setSize(flds[0].mesh().nCells());
718  // Add fld0
719  intFld.rmap(oldInternalField, cellProcAddressing[0]);
720 
721  for (label meshi = 1; meshi < flds.size(); meshi++)
722  {
723  if (flds.set(meshi))
724  {
725  const Field<Type>& addFld = flds[meshi];
726  intFld.rmap(addFld, cellProcAddressing[meshi]);
727  }
728  }
729  }
730 }
731 
732 
733 template<class Type>
734 void Foam::fvMeshAdder::MapVolField
735 (
737  const labelList& oldPatchStarts0,
738  const labelList& oldPatchSizes0,
739  const labelListList& patchProcAddressing,
740  const labelListList& cellProcAddressing,
741  const labelListList& faceProcAddressing,
742  const bool fullyMapped
743 )
744 {
745  // Add fields to fields[0] after adding the meshes to meshes[0].
746  // Mesh[0] is the sum of all meshes. Fields are not yet mapped.
747 
748  if (!flds.test(0))
749  {
751  << "Not valid field at element 0 in list of size "
752  << flds.size() << nl
753  << exit(FatalError);
754  }
755 
756 
757  // Internal field
758  // ~~~~~~~~~~~~~~
759 
760  {
761  // Store old internal field
762  const Field<Type> oldInternalField(flds[0].primitiveField());
763 
764  // Modify internal field
765  Field<Type>& intFld = flds[0].primitiveFieldRef();
766 
767  // Set to new mesh size
768  intFld.setSize(flds[0].mesh().nCells());
769  // Add fld0
770  intFld.rmap(oldInternalField, cellProcAddressing[0]);
771 
772  for (label meshi = 1; meshi < flds.size(); meshi++)
773  {
774  if (flds.set(meshi))
775  {
776  const Field<Type>& addFld = flds[meshi].primitiveFieldRef();
777  intFld.rmap(addFld, cellProcAddressing[meshi]);
778  }
779  }
780  }
781 
782 
783  // Patch fields from old mesh
784  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
785 
786  auto& bfld0 = flds[0].boundaryFieldRef();
787  //Pout<< " Adding patchfields for field " << flds[0].name() << endl;
788  forAll(bfld0, patchi)
789  {
790  //Pout<< " patch:" << patchi
791  // << " patch:" << flds[0].mesh().boundaryMesh()[patchi].name()
792  // << " size:" << flds[0].mesh().boundaryMesh()[patchi].size()
793  // << endl;
794  //Pout<< " patchField:" << bfld0[patchi].size()
795  // << endl;
796 
797  labelList newToOld
798  (
799  calcPatchMap
800  (
801  oldPatchStarts0[patchi],
802  oldPatchSizes0[patchi],
803  faceProcAddressing[0],
804  bfld0[patchi].patch().patch(),
805  -1 // unmapped value
806  )
807  );
808 
809  directFvPatchFieldMapper patchMapper(newToOld);
810 
811  // Override mapping (for use in e.g. fvMeshDistribute where
812  // it sorts mapping out itself)
813  if (fullyMapped)
814  {
815  patchMapper.hasUnmapped() = false;
816  }
817 
818  bfld0[patchi].autoMap(patchMapper);
819  }
820 
821  for (label meshi = 1; meshi < flds.size(); meshi++)
822  {
823  if (flds.set(meshi))
824  {
825  const auto& bfld = flds[meshi].boundaryFieldRef();
826 
827  const labelList& patchMap = patchProcAddressing[meshi];
828 
829  forAll(patchMap, oldPatchi)
830  {
831  const auto& fvp = bfld[oldPatchi].patch();
832  const label newPatchi = patchMap[oldPatchi];
833 
834  //Pout<< " oldPatch:" << oldPatchi
835  // << " newPatch:" << newPatchi << endl;
836 
837  if (newPatchi >= 0 && newPatchi < bfld0.size())
838  {
839  const auto& fvp0 = bfld0[newPatchi].patch();
840  labelList addedToNew(bfld[oldPatchi].size(), -1);
841  forAll(addedToNew, i)
842  {
843  const label newFacei =
844  faceProcAddressing[meshi][fvp.start()+i];
845  const label patchFacei = newFacei-fvp0.start();
846  if
847  (
848  patchFacei >= 0
849  && patchFacei < fvp0.size()
850  )
851  {
852  addedToNew[i] = patchFacei;
853  }
854  }
855 
856  bfld0[newPatchi].rmap(bfld[oldPatchi], addedToNew);
857  }
858  else
859  {
860  WarningInFunction << "Ignoring old patch "
861  << bfld[oldPatchi].patch().name() << " on field "
862  << flds[meshi].name() << endl; //exit(FatalError);
863  }
864  }
865  }
866  }
867 }
868 
869 
870 template<class Type>
871 void Foam::fvMeshAdder::MapSurfaceField
872 (
874  const labelList& oldFaceOwner0,
875  const labelList& oldPatchStarts0,
876  const labelList& oldPatchSizes0,
877  const labelListList& patchProcAddressing,
878  const labelListList& cellProcAddressing,
879  const labelListList& faceProcAddressing,
880  const bool fullyMapped
881 )
882 {
883  // Add fields to fields[0] after adding the meshes to meshes[0].
884  // Mesh[0] is the sum of all meshes. Fields are not yet mapped.
885 
886  if (!flds.test(0))
887  {
889  << "Not valid field at element 0 in list of size "
890  << flds.size() << nl
891  << exit(FatalError);
892  }
893 
894  const fvMesh& mesh0 = flds[0].mesh();
895 
896 
897  // Internal field
898  // ~~~~~~~~~~~~~~
899 
900  {
901  // Store old internal field
902  const Field<Type> oldInternalField(flds[0].primitiveField());
903 
904  // Modify internal field
905  Field<Type>& intFld = flds[0].primitiveFieldRef();
906 
907  // Set to new mesh size
908  intFld.setSize(mesh0.nInternalFaces());
909 
910  // Map
911  forAll(flds, meshi)
912  {
913  if (flds.set(meshi))
914  {
915  const labelList& faceMap = faceProcAddressing[meshi];
916  const auto& fld = flds[meshi];
917 
918  // Map internal field
919  if (meshi == 0)
920  {
921  intFld.rmap(oldInternalField, faceMap);
922  }
923  else
924  {
925  intFld.rmap(fld.primitiveField(), faceMap);
926  }
927 
928  // Map faces that were boundary faces but are not anymore.
929  // There will be two meshes that provide the same face. Use
930  // owner one.
931  const auto& bfld = flds[meshi].boundaryField();
932 
933  forAll(bfld, oldPatchi)
934  {
935  const fvsPatchField<Type>& pf = bfld[oldPatchi];
936  //Pout<< "patch:" << pf.patch().name() << endl;
937  forAll(pf, patchFacei)
938  {
939  // Get new face, mapped face owner
940  label newFacei;
941  label newOwn;
942  if (meshi == 0)
943  {
944  // Do not access mesh information since in-place
945  // modified
946  const label oldFacei =
947  oldPatchStarts0[oldPatchi]+patchFacei;
948  newFacei = faceProcAddressing[meshi][oldFacei];
949  const label oldOwn = oldFaceOwner0[oldFacei];
950  newOwn = cellProcAddressing[meshi][oldOwn];
951 
952  //Pout<< "MESH0: pfi:" << patchFacei
953  // << " old face:" << oldFacei
954  // << " new face:" << newFacei
955  // << " at:" << mesh0.faceCentres()[newFacei]
956  // << " oldOwn:" << oldOwn
957  // << " newOwn:" << newOwn << endl;
958  }
959  else
960  {
961  const label oldFacei =
962  pf.patch().start()+patchFacei;
963  newFacei = faceProcAddressing[meshi][oldFacei];
964  const label oldOwn =
965  fld.mesh().faceOwner()[oldFacei];
966  newOwn = cellProcAddressing[meshi][oldOwn];
967 
968  //Pout<< "MESH:" << meshi << " pfi:" << patchFacei
969  // << " old face:" << oldFacei
970  // << " new face:" << newFacei
971  // << " at:" << mesh0.faceCentres()[newFacei]
972  // << " oldOwn:" << oldOwn
973  // << " newOwn:" << newOwn << endl;
974  }
975 
976  if
977  (
978  newFacei >= 0
979  && newFacei < mesh0.nInternalFaces()
980  && (newOwn == mesh0.faceOwner()[newFacei])
981  )
982  {
983  intFld[newFacei] = pf[patchFacei];
984  }
985  //else
986  //{
987  // Pout<< " ignoring pfi:" << patchFacei
988  // << " value:" << pf[patchFacei]
989  // << " since newFacei:" << newFacei
990  // << " since newOwn:" << newOwn
991  // << " own:" << mesh0.faceOwner()[newFacei]
992  // << endl;
993  //}
994  }
995  }
996  }
997  }
998  }
999 
1000 
1001  // Patch fields from old mesh
1002  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
1003 
1004  auto& bfld0 = flds[0].boundaryFieldRef();
1005  forAll(bfld0, patchi)
1006  {
1007  labelList newToOld
1008  (
1009  calcPatchMap
1010  (
1011  oldPatchStarts0[patchi],
1012  oldPatchSizes0[patchi],
1013  faceProcAddressing[0],
1014  bfld0[patchi].patch().patch(),
1015  -1 // unmapped value
1016  )
1017  );
1018 
1019  directFvPatchFieldMapper patchMapper(newToOld);
1020 
1021  // Override mapping (for use in e.g. fvMeshDistribute where
1022  // it sorts mapping out itself)
1023  if (fullyMapped)
1024  {
1025  patchMapper.hasUnmapped() = false;
1026  }
1027 
1028  bfld0[patchi].autoMap(patchMapper);
1029  }
1030 
1031  for (label meshi = 1; meshi < flds.size(); meshi++)
1032  {
1033  if (flds.set(meshi))
1034  {
1035  const auto& bfld = flds[meshi].boundaryFieldRef();
1036 
1037  const labelList& patchMap = patchProcAddressing[meshi];
1038 
1039  forAll(patchMap, oldPatchi)
1040  {
1041  const auto& fvp = bfld[oldPatchi].patch();
1042  const label newPatchi = patchMap[oldPatchi];
1043  if (newPatchi >= 0 && newPatchi < bfld0.size())
1044  {
1045  const auto& fvp0 = bfld0[newPatchi].patch();
1046  labelList addedToNew(bfld[oldPatchi].size(), -1);
1047  forAll(addedToNew, i)
1048  {
1049  const label newFacei =
1050  faceProcAddressing[meshi][fvp.start()+i];
1051  const label patchFacei = newFacei-fvp0.start();
1052  if
1053  (
1054  patchFacei >= 0
1055  && patchFacei < fvp0.size()
1056  )
1057  {
1058  addedToNew[i] = patchFacei;
1059  }
1060  }
1061 
1062  bfld0[newPatchi].rmap(bfld[oldPatchi], addedToNew);
1063  }
1064  else
1065  {
1066  WarningInFunction << "Ignoring old patch "
1067  << bfld[oldPatchi].patch().name() << " on field "
1068  << flds[meshi].name() << endl; //exit(FatalError);
1069  }
1070  }
1071  }
1072  }
1074 
1075 
1076 template<class Type>
1078 (
1079  const UPtrList<fvMesh>& meshes,
1080  const labelList& oldPatchStarts0,
1081  const labelList& oldPatchSizes0,
1082  const labelListList& patchProcAddressing,
1083  const labelListList& cellProcAddressing,
1084  const labelListList& faceProcAddressing,
1085  const labelListList& pointProcAddressing,
1086  const bool fullyMapped
1087 )
1088 {
1090 
1091  if (!meshes.test(0))
1092  {
1094  << "Not valid mesh at element 0 in list of size "
1095  << meshes.size() << nl
1096  << exit(FatalError);
1097  }
1098  const auto& mesh0 = meshes[0];
1099 
1100  const UPtrList<const fldType> fields
1101  (
1102  mesh0.objectRegistry::csorted<fldType>()
1103  );
1104 
1105 
1106  // It is necessary to enforce that all old-time fields are stored
1107  // before the mapping is performed. Otherwise, if the
1108  // old-time-level field is mapped before the field itself, sizes
1109  // will not match.
1110 
1111  for (const auto& field : fields)
1112  {
1113  DebugPout
1114  << "MapVolFields : Storing old time for "
1115  << field.name() << endl;
1116 
1117  const_cast<fldType&>(field).storeOldTimes();
1118  }
1119 
1120 
1121  for (const auto& field : fields)
1122  {
1123  const word& name0 = field.name();
1124 
1125  DebugPout
1126  << "MapVolFields : mapping " << name0 << endl;
1127 
1128  UPtrList<fldType> meshToField(meshes.size());
1129  forAll(meshes, meshi)
1130  {
1131  if (meshes.set(meshi))
1132  {
1133  auto& meshFld = meshes[meshi].
1134  objectRegistry::lookupObjectRef<fldType>(name0);
1135  meshToField.set(meshi, &meshFld);
1136  }
1137  }
1138 
1139  MapVolField
1140  (
1141  meshToField,
1142  oldPatchStarts0,
1143  oldPatchSizes0,
1144  patchProcAddressing,
1145  cellProcAddressing,
1146  faceProcAddressing,
1147  fullyMapped
1148  );
1149  }
1151 
1152 
1153 template<class Type>
1155 (
1156  const UPtrList<fvMesh>& meshes,
1157  const labelListList& cellProcAddressing,
1158  const bool fullyMapped
1159 )
1160 {
1161  typedef DimensionedField<Type, volMesh> fldType;
1162 
1163  if (!meshes.test(0))
1164  {
1166  << "Not valid mesh at element 0 in list of size "
1167  << meshes.size() << nl
1168  << exit(FatalError);
1169  }
1170  const auto& mesh0 = meshes[0];
1171 
1172  // NB: strict=true to avoid picking up volFields
1173  const UPtrList<const fldType> fields
1174  (
1175  mesh0.objectRegistry::csorted<fldType, true>()
1176  );
1177 
1178  for (const auto& field : fields)
1179  {
1180  const word& name0 = field.name();
1181 
1182  DebugPout
1183  << "MapDimFields : mapping " << name0 << endl;
1184 
1185  UPtrList<fldType> meshToField(meshes.size());
1186  forAll(meshes, meshi)
1187  {
1188  if (meshes.set(meshi))
1189  {
1190  auto& meshFld = meshes[meshi].
1191  objectRegistry::lookupObjectRef<fldType>(name0);
1192  meshToField.set(meshi, &meshFld);
1193  }
1194  }
1195 
1196  MapDimField(meshToField, cellProcAddressing, fullyMapped);
1197  }
1199 
1200 
1201 template<class Type>
1203 (
1204  const UPtrList<fvMesh>& meshes,
1205  const labelList& oldFaceOwner0,
1206  const labelList& oldPatchStarts0,
1207  const labelList& oldPatchSizes0,
1208  const labelListList& patchProcAddressing,
1209  const labelListList& cellProcAddressing,
1210  const labelListList& faceProcAddressing,
1211  const labelListList& pointProcAddressing,
1212  const bool fullyMapped
1213 )
1214 {
1216 
1217  if (!meshes.test(0))
1218  {
1220  << "Not valid mesh at element 0 in list of size "
1221  << meshes.size() << nl
1222  << exit(FatalError);
1223  }
1224  const auto& mesh0 = meshes[0];
1225 
1226  const UPtrList<const fldType> fields
1227  (
1228  mesh0.objectRegistry::csorted<fldType>()
1229  );
1230 
1231 
1232  // It is necessary to enforce that all old-time fields are stored
1233  // before the mapping is performed. Otherwise, if the
1234  // old-time-level field is mapped before the field itself, sizes
1235  // will not match.
1236 
1237  for (const auto& field : fields)
1238  {
1239  DebugPout
1240  << "MapSurfaceFields : Storing old time for "
1241  << field.name() << endl;
1242 
1243  const_cast<fldType&>(field).storeOldTimes();
1244  }
1245 
1246 
1247  for (const auto& field : fields)
1248  {
1249  const word& name0 = field.name();
1250 
1251  DebugPout
1252  << "MapSurfaceFields : Mapping " << field.name() << endl;
1253 
1254  UPtrList<fldType> meshToField(meshes.size());
1255  forAll(meshes, meshi)
1256  {
1257  if (meshes.set(meshi))
1258  {
1259  auto& meshFld = meshes[meshi].
1260  objectRegistry::lookupObjectRef<fldType>(name0);
1261  meshToField.set(meshi, &meshFld);
1262  }
1263  }
1264 
1265  MapSurfaceField
1266  (
1267  meshToField,
1268  oldFaceOwner0,
1269  oldPatchStarts0,
1270  oldPatchSizes0,
1271  patchProcAddressing,
1272  cellProcAddressing,
1273  faceProcAddressing,
1274  fullyMapped
1275  );
1276  }
1277 }
1278 
1279 
1280 // ************************************************************************* //
Foam::surfaceFields.
const T * test(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: UPtrListI.H:127
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
rDeltaTY field()
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Class containing mesh-to-mesh mapping information after a mesh addition where we add a mesh (&#39;added m...
Generic GeometricField class.
Definition: areaFieldsFwd.H:50
DirectFieldMapper< fvPatchFieldMapper > directFvPatchFieldMapper
A fvPatchFieldMapper with direct mapping.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:97
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
const polyMesh & mesh() const noexcept
Return the mesh reference.
Generic templated field type.
Definition: Field.H:62
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
label size() const noexcept
The number of entries in the list.
Definition: UPtrListI.H:106
Foam::PtrList< Foam::fvMesh > meshes(regionNames.size())
label nInternalFaces() const noexcept
Number of internal faces.
A list of pointers to objects of type <T>, without allocation/deallocation management of the pointers...
Definition: HashTable.H:106
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie...
Definition: PtrList.H:159
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
static void MapVolFields(const mapAddedPolyMesh &, const fvMesh &mesh, const fvMesh &meshToAdd, const bool fullyMapped=false)
Map all volFields of Type.
void rmap(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 reverse-map from the given field
Definition: Field.C:529
static void MapSurfaceFields(const mapAddedPolyMesh &, const fvMesh &mesh, const fvMesh &meshToAdd, const bool fullyMapped=false)
Map all surfaceFields of Type.
#define WarningInFunction
Report a warning using Foam::Warning.
label nCells() const noexcept
Number of mesh cells.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
const std::string patch
OpenFOAM patch number as a std::string.
#define DebugPout
Report an information message using Foam::Pout.
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
static void MapDimFields(const mapAddedPolyMesh &, const fvMesh &mesh, const fvMesh &meshToAdd)
Map all DimensionedFields of Type.
List< label > labelList
A List of labels.
Definition: List.H:62
static tmp< fvPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
static tmp< fvsPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, surfaceMesh > &)
Return a pointer to a new patchField created on freestore given.