44 Foam::label Foam::mergePolyMesh::patchIndex(
const polyPatch& p)
48 const word& pType = p.type();
49 const word& pName = p.name();
51 bool nameFound =
false;
53 forAll(patchNames_, patchi)
55 if (patchNames_[patchi] == pName)
57 if (word(patchDicts_[patchi][
"type"]) == pType)
74 patchDicts_.append(dictionary(IStringStream(os.str())()));
81 const word& caseName = p.boundaryMesh().mesh().time().caseName();
83 patchNames_.append(pName +
"_" + caseName);
85 Info<<
"label patchIndex(const polyPatch& p) : " 86 <<
"Patch " << p.index() <<
" named " 87 << pName <<
" in mesh " << caseName
88 <<
" already exists, but patch types " 89 <<
" do not match.\nCreating a composite name as " 90 << patchNames_.last() <<
endl;
94 patchNames_.append(pName);
97 return patchNames_.size() - 1;
103 DynamicList<word>& names,
109 if (names[zoneI] == curName)
116 names.append(curName);
118 return names.size() - 1;
128 patchNames_(2*boundaryMesh().size()),
129 patchDicts_(2*boundaryMesh().size()),
135 wordList curPatchNames = boundaryMesh().names();
137 forAll(boundaryMesh(), patchi)
139 patchNames_.append(boundaryMesh()[patchi].
name());
142 boundaryMesh()[
patchi].write(os);
143 patchDicts_.append(dictionary(IStringStream(os.str())()));
149 wordList curPointZoneNames = pointZones().names();
150 if (curPointZoneNames.size())
152 pointZoneNames_.setCapacity(2*curPointZoneNames.size());
155 forAll(curPointZoneNames, zoneI)
157 pointZoneNames_.append(curPointZoneNames[zoneI]);
161 wordList curFaceZoneNames = faceZones().names();
163 if (curFaceZoneNames.size())
165 faceZoneNames_.setCapacity(2*curFaceZoneNames.size());
167 forAll(curFaceZoneNames, zoneI)
169 faceZoneNames_.append(curFaceZoneNames[zoneI]);
173 wordList curCellZoneNames = cellZones().names();
175 if (curCellZoneNames.size())
177 cellZoneNames_.setCapacity(2*curCellZoneNames.size());
179 forAll(curCellZoneNames, zoneI)
181 cellZoneNames_.append(curCellZoneNames[zoneI]);
207 pointZoneIndices[zoneI] = zoneIndex(pointZoneNames_, pz[zoneI].
name());
213 zoneID = pz.whichZone(pointi);
218 zoneID = pointZoneIndices[zoneID];
221 renumberPoints[pointi] =
244 cellZoneIndices[zoneI] = zoneIndex(cellZoneNames_, cz[zoneI].
name());
250 zoneID = cz.whichZone(celli);
255 zoneID = cellZoneIndices[zoneID];
258 renumberCells[celli] =
273 const polyBoundaryMesh& bm = m.boundaryMesh();
278 forAll(patchIndices, patchi)
280 patchIndices[
patchi] = patchIndex(bm[patchi]);
285 meshMod_.setNumPatches(patchNames_.size());
294 faceZoneIndices[zoneI] = zoneIndex(faceZoneNames_, fz[zoneI].
name());
301 const labelList& nei = m.faceNeighbour();
303 label newOwn, newNei, newPatch, newZone;
308 const face& curFace = f[facei];
310 face newFace(curFace.size());
314 newFace[pointi] = renumberPoints[curFace[pointi]];
320 if (
min(newFace) < 0)
323 <<
"Error in point mapping for face " << facei
324 <<
". Old face: " << curFace <<
" New face: " << newFace
329 if (facei < m.nInternalFaces() || facei >= m.nFaces())
335 newPatch = patchIndices[bm.whichPatch(facei)];
339 if (newOwn > -1) newOwn = renumberCells[newOwn];
348 newNei = renumberCells[newNei];
352 newZone = fz.whichZone(facei);
357 newZoneFlip = fz[newZone].flipMap()[fz[newZone].whichFace(facei)];
360 newZone = faceZoneIndices[newZone];
363 renumberFaces[facei] =
386 Info<<
"patch names: " << patchNames_ <<
nl 387 <<
"patch dicts: " << patchDicts_ <<
nl 388 <<
"point zone names: " << pointZoneNames_ <<
nl 389 <<
"face zone names: " << faceZoneNames_ <<
nl 390 <<
"cell zone names: " << cellZoneNames_ <<
endl;
393 if (patchNames_.size() != boundaryMesh().size())
395 Info<<
"Copying old patches" <<
endl;
397 List<polyPatch*> newPatches(patchNames_.size());
399 const polyBoundaryMesh& oldPatches = boundaryMesh();
404 for (patchi = 0; patchi < oldPatches.size(); patchi++)
406 newPatches[
patchi] = oldPatches[
patchi].clone(oldPatches).ptr();
409 Info<<
"Adding new patches. " <<
endl;
411 label endOfLastPatch =
414 : oldPatches[patchi - 1].start() + oldPatches[patchi - 1].size();
416 for (; patchi < patchNames_.size(); patchi++)
419 dictionary
dict(patchDicts_[patchi]);
421 dict.
set(
"startFace", endOfLastPatch);
436 addPatches(newPatches);
440 if (pointZoneNames_.size() > pointZones().size())
442 Info<<
"Adding new pointZones. " <<
endl;
443 label nZones = pointZones().size();
445 pointZones().setSize(pointZoneNames_.size());
447 for (
label zoneI = nZones; zoneI < pointZoneNames_.size(); zoneI++)
454 pointZoneNames_[zoneI],
462 if (cellZoneNames_.size() > cellZones().size())
464 Info<<
"Adding new cellZones. " <<
endl;
466 label nZones = cellZones().size();
468 cellZones().setSize(cellZoneNames_.size());
470 for (
label zoneI = nZones; zoneI < cellZoneNames_.size(); zoneI++)
477 cellZoneNames_[zoneI],
485 if (faceZoneNames_.size() > faceZones().size())
487 Info<<
"Adding new faceZones. " <<
endl;
489 label nZones = faceZones().size();
491 faceZones().setSize(faceZoneNames_.size());
493 for (
label zoneI = nZones; zoneI < faceZoneNames_.size(); zoneI++)
500 faceZoneNames_[zoneI],
511 meshMod_.changeMesh(*
this,
false);
#define forAll(list, i)
Loop across all elements in list.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const dictionary & dict() const
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const word & name() const
Ostream & endl(Ostream &os)
Add newline and flush stream.
mergePolyMesh(const IOobject &io)
Construct from IOobject.
void addMesh(const polyMesh &m)
Add a mesh.
List< bool > boolList
Bool container classes.
vectorField pointField
pointField is a vectorField.
List< label > labelList
A List of labels.
MeshZones< pointZone, polyMesh > meshPointZones
A MeshZones with the type pointZone.
errorManip< error > abort(error &err)
MeshZones< cellZone, polyMesh > meshCellZones
A MeshZones with the type cellZone.
void merge()
Merge meshes.
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
MeshZones< faceZone, polyMesh > meshFaceZones
A MeshZones with the type faceZone.
List< word > wordList
A List of words.
void set(entry *)
Assign a new entry, overwrite any existing entry.
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.
List< cell > cellList
list of cells