40 Foam::label Foam::mergePolyMesh::patchIndex(
const polyPatch&
p)
44 const word& pType =
p.type();
45 const word& pName =
p.name();
47 bool nameFound =
false;
51 if (patchNames_[
patchi] == pName)
53 if (word(patchDicts_[
patchi][
"type"]) == pType)
70 patchDicts_.append(dictionary(IStringStream(os.str())()));
77 const word& caseName =
p.boundaryMesh().mesh().time().caseName();
79 patchNames_.
append(pName +
"_" + caseName);
81 Info<<
"label patchIndex(const polyPatch& p) : "
82 <<
"Patch " <<
p.index() <<
" named "
83 << pName <<
" in mesh " << caseName
84 <<
" already exists, but patch types "
85 <<
" do not match.\nCreating a composite name as "
93 return patchNames_.
size() - 1;
99 DynamicList<word>& names,
105 if (names[zonei] == curName)
112 names.append(curName);
114 return names.size() - 1;
124 patchNames_(2*mesh_.boundaryMesh().size()),
125 patchDicts_(2*mesh_.boundaryMesh().size()),
131 wordList curPatchNames = mesh_.boundaryMesh().names();
135 patchNames_.append(mesh_.boundaryMesh()[
patchi].name());
138 mesh_.boundaryMesh()[
patchi].write(os);
139 patchDicts_.append(dictionary(IStringStream(os.str())()));
145 wordList curPointZoneNames = mesh_.pointZones().toc();
146 if (curPointZoneNames.size())
148 pointZoneNames_.setCapacity(2*curPointZoneNames.size());
151 forAll(curPointZoneNames, zonei)
153 pointZoneNames_.append(curPointZoneNames[zonei]);
156 pointZonesAddedPoints_.setSize(pointZoneNames_.size());
159 wordList curFaceZoneNames = mesh_.faceZones().toc();
161 if (curFaceZoneNames.size())
163 faceZoneNames_.setCapacity(2*curFaceZoneNames.size());
166 forAll(curFaceZoneNames, zonei)
168 faceZoneNames_.append(curFaceZoneNames[zonei]);
171 faceZonesAddedFaces_.setSize(faceZoneNames_.size());
174 wordList curCellZoneNames = mesh_.cellZones().toc();
176 if (curCellZoneNames.size())
178 cellZoneNames_.setCapacity(2*curCellZoneNames.size());
181 forAll(curCellZoneNames, zonei)
183 cellZoneNames_.append(curCellZoneNames[zonei]);
186 cellZonesAddedCells_.setSize(cellZoneNames_.size());
201 const pointZoneList& pz = m.pointZones();
206 pointZoneIndices[zonei] = zoneIndex(pointZoneNames_, pz[zonei].
name());
207 pointZonesAddedPoints_.setSize(pointZoneNames_.
size());
212 renumberPoints[pointi] = meshMod_.
addPoint
219 const labelList zones(pz.whichZones(pointi));
222 pointZonesAddedPoints_[pointZoneIndices[zonei]]
223 .insert(renumberPoints[pointi]);
232 const cellZoneList& cz = m.cellZones();
237 cellZoneIndices[zonei] = zoneIndex(cellZoneNames_, cz[zonei].
name());
238 cellZonesAddedCells_.setSize(cellZoneNames_.
size());
243 renumberCells[celli] = meshMod_.
addCell(-1);
245 const labelList zones(cz.whichZones(celli));
248 cellZonesAddedCells_[cellZoneIndices[zonei]]
249 .insert(renumberCells[celli]);
254 const polyBoundaryMesh& bm = m.boundaryMesh();
270 const faceZoneList& fz = m.faceZones();
275 faceZoneIndices[zonei] = zoneIndex(faceZoneNames_, fz[zonei].
name());
276 faceZonesAddedFaces_.setSize(faceZoneNames_.
size());
283 const labelList& nei = m.faceNeighbour();
285 label newOwn, newNei, newPatch;
289 const face& curFace =
f[facei];
291 face newFace(curFace.size());
295 newFace[pointi] = renumberPoints[curFace[pointi]];
301 if (
min(newFace) < 0)
304 <<
"Error in point mapping for face " << facei
305 <<
". Old face: " << curFace <<
" New face: " << newFace
310 if (facei < m.nInternalFaces() || facei >= m.nFaces())
316 newPatch = patchIndices[bm.whichPatch(facei)];
320 if (newOwn > -1) newOwn = renumberCells[newOwn];
329 newNei = renumberCells[newNei];
332 renumberFaces[facei] = meshMod_.
addFace
342 const labelList zones(fz.whichZones(facei));
345 const faceZone& fzi = fz[zones[zonei]];
346 const bool flip = fzi.flipMap()[fzi.localIndex(facei)];
348 faceZonesAddedFaces_[faceZoneIndices[zonei]]
349 .insert(renumberFaces[facei], flip);
359 Info<<
"patch names: " << patchNames_ <<
nl
360 <<
"patch dicts: " << patchDicts_ <<
nl
361 <<
"point zone names: " << pointZoneNames_ <<
nl
362 <<
"face zone names: " << faceZoneNames_ <<
nl
363 <<
"cell zone names: " << cellZoneNames_ <<
endl;
367 if (patchNames_.size() != mesh_.boundaryMesh().size())
371 Info<<
"Copying old patches" <<
endl;
374 List<polyPatch*> newPatches(patchNames_.size());
376 const polyBoundaryMesh& oldPatches = mesh_.boundaryMesh();
383 newPatches[
patchi] = oldPatches[
patchi].clone(oldPatches).ptr();
388 Info<<
"Adding new patches. " <<
endl;
391 label endOfLastPatch =
394 : oldPatches[
patchi - 1].start() + oldPatches[
patchi - 1].size();
401 dict.
set(
"startFace", endOfLastPatch);
415 mesh_.removeBoundary();
416 mesh_.addPatches(newPatches);
420 if (pointZoneNames_.size() > mesh_.pointZones().size())
424 Info<<
"Adding new pointZones. " <<
endl;
427 label nZones = mesh_.pointZones().size();
429 mesh_.pointZones().setSize(pointZoneNames_.size());
431 for (
label zonei = nZones; zonei < pointZoneNames_.size(); zonei++)
433 mesh_.pointZones().set
438 pointZoneNames_[zonei],
446 if (cellZoneNames_.size() > mesh_.cellZones().size())
450 Info<<
"Adding new cellZones. " <<
endl;
453 label nZones = mesh_.cellZones().size();
455 mesh_.cellZones().setSize(cellZoneNames_.size());
457 for (
label zonei = nZones; zonei < cellZoneNames_.size(); zonei++)
459 mesh_.cellZones().set
464 cellZoneNames_[zonei],
472 if (faceZoneNames_.size() > mesh_.faceZones().size())
476 Info<<
"Adding new faceZones. " <<
endl;
479 label nZones = mesh_.faceZones().size();
481 mesh_.faceZones().setSize(faceZoneNames_.size());
483 for (
label zonei = nZones; zonei < faceZoneNames_.size(); zonei++)
485 mesh_.faceZones().set
490 faceZoneNames_[zonei],
500 autoPtr<polyTopoChangeMap> map(meshMod_.changeMesh(mesh_));
503 mesh_.pointZones().insert(pointZonesAddedPoints_);
506 mesh_.faceZones().insert(faceZonesAddedFaces_);
509 mesh_.cellZones().insert(cellZonesAddedCells_);
511 mesh_.topoChange(map);
#define forAll(list, i)
Loop across all elements in list.
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
void size(const label)
Override size to be inconsistent with allocated storage.
T & last()
Return the last element of the list.
void set(entry *)
Assign a new entry, overwrite any existing entry.
void addMesh(const polyMesh &m)
Add a mesh.
mergePolyMesh(polyMesh &mesh)
Construct from polyMesh.
void merge()
Merge meshes.
static autoPtr< polyPatch > New(const word &patchType, const word &name, const label size, const label start, const label index, const polyBoundaryMesh &bm)
Return a pointer to a new patch created on freestore from.
void setNumPatches(const label nPatches)
Explicitly set the number of patches if construct-without-mesh.
label addCell(const label masterCellID)
Add cell and return new cell index.
label addPoint(const point &, const label masterPointID, const bool inCell)
Add point and return new point index.
label addFace(const face &f, const label own, const label nei, const label masterFaceID, const bool flipFaceFlux, const label patchID)
Add face to cells and return new face index.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const dimensionedScalar c
Speed of light in a vacuum.
List< word > wordList
A List of words.
List< label > labelList
A List of labels.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
List< cell > cellList
list of cells
Ostream & endl(Ostream &os)
Add newline and flush stream.
word name(const bool)
Return a word representation of a bool.
errorManip< error > abort(error &err)
vectorField pointField
pointField is a vectorField.
List< bool > boolList
Bool container classes.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
defineTypeNameAndDebug(combustionModel, 0)