View Issue Details
ID | Project | Category | View Status | Date Submitted | Last Update |
---|---|---|---|---|---|
0001817 | OpenFOAM | Bug | public | 2015-08-07 14:35 | 2015-08-07 15:54 |
Reporter | wyldckat | Assigned To | henry | ||
Priority | normal | Severity | minor | Reproducibility | always |
Status | resolved | Resolution | fixed | ||
Platform | GNU/Linux | OS | Ubuntu | OS Version | 12.04 |
Product Version | dev | ||||
Summary | 0001817: Building with 64-bit labels triggers error messages in blockMeshMergeFast.C | ||||
Description | The error messages should be self-explanatory ;) blockMesh/blockMeshMergeFast.C: In function ‘void Foam::genFaceFaceRotMap()’: blockMesh/blockMeshMergeFast.C:78:31: error: call of overloaded ‘mag(int&)’ is ambiguous ... blockMesh/blockMeshMergeFast.C: In function ‘Foam::label Foam::mapij(int, Foam::label, Foam::label)’: blockMesh/blockMeshMergeFast.C:200:34: error: call of overloaded ‘mag(const int&)’ is ambiguous ... blockMesh/blockMeshMergeFast.C: In member function ‘void Foam::blockMesh::calcMergeInfoFast()’: blockMesh/blockMeshMergeFast.C:408:40: error: call of overloaded ‘mag(int&)’ is ambiguous ... blockMesh/blockMeshMergeFast.C:409:40: error: call of overloaded ‘mag(int&)’ is ambiguous Attached is the updated file "blockMeshMergeFast.C" for the path "src/mesh/blockMesh/blockMesh/", as well as a patch for visual inspection here on the bug tracker. | ||||
Tags | No tags attached. | ||||
|
blockMeshMergeFast.C (16,226 bytes)
/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2015 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License This file is part of OpenFOAM. OpenFOAM is free software: you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation, either version 3 of the License, or (at your option) any later version. OpenFOAM is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details. You should have received a copy of the GNU General Public License along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>. \*---------------------------------------------------------------------------*/ #include "blockMesh.H" // * * * * * * * * * * * * * * * Private Functions * * * * * * * * * * * * * // namespace Foam { // faces // 6 // ( // 4(0 4 7 3) // x-min // 4(1 2 6 5) // x-max // 4(0 1 5 4) // y-min // 4(3 7 6 2) // y-max // 4(0 3 2 1) // z-min // 4(4 5 6 7) // z-max // ); // Face-edge directions static const label faceEdgeDirs[6][4] = { {2, 1, -2, -1}, {1, 2, -1, -2}, {1, 2, -1, -2}, {2, 1, -2, -1}, {2, 1, -2, -1}, {1, 2, -1, -2} }; // The face-face-rotation direction correspondence map static Pair<label> faceFaceRotMap[6][6][4]; // Generate the face-face-rotation direction correspondence map void genFaceFaceRotMap() { for(label facePi=0; facePi<6; facePi++) { for(label faceNi=0; faceNi<6; faceNi++) { for(label rot=0; rot<4; rot++) { Pair<label>& map = faceFaceRotMap[facePi][faceNi][rot]; for(label Pp=0; Pp<2; Pp++) { label Pdir = faceEdgeDirs[facePi][Pp]; label Np = (3 - Pp + rot)%4; label Ndir = faceEdgeDirs[faceNi][Np]; map[Pdir-1] = -Ndir; } // Handle sign change due to match-face transpose if (mag(map[0]) == 2 && map[0]*map[1] < 0) { map[0] = -map[0]; map[1] = -map[1]; } } } } } // Return the direction map for the merge-faces Pair<label> faceMap ( const label facePi, const face& faceP, const label faceNi, const face& faceN ) { // Search for the point on faceN corresponding to the 0-point on faceP for(label rot=0; rot<4; rot++) { if (faceN[rot] == faceP[0]) { return faceFaceRotMap[facePi][faceNi][rot]; } } FatalErrorIn ( "faceMap(const label facePi, const face& faceP, " "const label faceNi, const face& faceN)" ) << "Cannot find point correspondance for faces " << faceP << " and " << faceN << exit(FatalError); return Pair<label>(0, 0); } // Set the block and face indices for all the merge faces void setBlockFaceCorrespondence ( const cellList& topoCells, const faceList::subList& topoInternalFaces, const labelList& topoFaceCell, List<Pair<label> >& mergeBlock ) { forAll(topoInternalFaces, topoFacei) { label topoPi = topoFaceCell[topoFacei]; const labelList& topoPfaces = topoCells[topoPi]; bool foundFace = false; label topoPfacei; for ( topoPfacei = 0; topoPfacei < topoPfaces.size(); topoPfacei++ ) { if (topoPfaces[topoPfacei] == topoFacei) { foundFace = true; break; } } if (!foundFace) { FatalErrorIn("setBlockFaceCorrespondence()") << "Cannot find merge face for block " << topoPi << exit(FatalError); } mergeBlock[topoFacei].first() = topoPi; mergeBlock[topoFacei].second() = topoPfacei; } } // Return the number of divisions in each direction for the face Pair<label> faceNij(const label facei, const block& block) { Pair<label> fnij; label i = facei/2; if (i == 0) { fnij.first() = block.meshDensity().y() + 1; fnij.second() = block.meshDensity().z() + 1; } else if (i == 1) { fnij.first() = block.meshDensity().x() + 1; fnij.second() = block.meshDensity().z() + 1; } else if (i == 2) { fnij.first() = block.meshDensity().x() + 1; fnij.second() = block.meshDensity().y() + 1; } return fnij; } // Sign the index corresponding to the map inline label signIndex(const label map, const label i) { return map < 0 ? -i-1 : i; } // Reverse a signed index with the number of divisions inline label unsignIndex(const label i, const label ni) { return i >= 0 ? i : ni + i + 1; } // Return the mapped index inline label mapij(const label map, const label i, const label j) { return signIndex(map, mag(map) == 1 ? i : j); } // Return the face point index inline label facePoint ( const label facei, const block& block, const label i, const label j ) { switch (facei) { case 0: return block.vtxLabel(0, i, j); case 1: return block.vtxLabel(block.meshDensity().x(), i, j); case 2: return block.vtxLabel(i, 0, j); case 3: return block.vtxLabel(i, block.meshDensity().y(), j); case 4: return block.vtxLabel(i, j, 0); case 5: return block.vtxLabel(i, j, block.meshDensity().z()); default: return -1; } } // Return the neighbour face point from the signed indices inline label facePointN ( const block& block, const label i, const label j, const label k ) { return block.vtxLabel ( unsignIndex(i, block.meshDensity().x()), unsignIndex(j, block.meshDensity().y()), unsignIndex(k, block.meshDensity().z()) ); } // Return the neighbour face point from the mapped indices inline label facePointN ( const label facei, const block& block, const label i, const label j ) { switch (facei) { case 0: return facePointN(block, 0, i, j); case 1: return facePointN(block, block.meshDensity().x(), i, j); case 2: return facePointN(block, i, 0, j); case 3: return facePointN(block, i, block.meshDensity().y(), j); case 4: return facePointN(block, i, j, 0); case 5: return facePointN(block, i, j, block.meshDensity().z()); default: return -1; } } // Return the neighbour face point using the map inline label facePointN ( const label facei, const Pair<label>& fmap, const block& block, const label i, const label j ) { return facePointN(facei, block, mapij(fmap[0], i, j), mapij(fmap[1], i, j)); } } // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // void Foam::blockMesh::calcMergeInfoFast() { // Generate the static face-face map genFaceFaceRotMap(); const blockList& blocks = *this; if (verboseOutput) { Info<< "Creating block offsets" << endl; } blockOffsets_.setSize(blocks.size()); nPoints_ = 0; nCells_ = 0; forAll(blocks, blockI) { blockOffsets_[blockI] = nPoints_; nPoints_ += blocks[blockI].nPoints(); nCells_ += blocks[blockI].nCells(); } if (verboseOutput) { Info<< "Creating merge list using the fast topological search" << flush; } // Size merge list and initialize to -1 mergeList_.setSize(nPoints_, -1); // Block mesh topology const pointField& topoPoints = topology().points(); const cellList& topoCells = topology().cells(); const faceList& topoFaces = topology().faces(); const labelList& topoFaceOwn = topology().faceOwner(); const labelList& topoFaceNei = topology().faceNeighbour(); // Topological merging is only necessary for internal block faces // Note edge and face collapse may apply to boundary faces // but is not yet supported in the "fast" algorithm const faceList::subList topoInternalFaces ( topoFaces, topology().nInternalFaces() ); List<Pair<label> > mergeBlockP(topoInternalFaces.size()); setBlockFaceCorrespondence ( topoCells, topoInternalFaces, topoFaceOwn, mergeBlockP ); List<Pair<label> > mergeBlockN(topoInternalFaces.size()); setBlockFaceCorrespondence ( topoCells, topoInternalFaces, topoFaceNei, mergeBlockN ); if (debug) { Info<< endl; } forAll(topoInternalFaces, topoFacei) { if (debug) { Info<< "Processing face " << topoFacei << endl; } label blockPi = mergeBlockP[topoFacei].first(); label blockPfacei = mergeBlockP[topoFacei].second(); label blockNi = mergeBlockN[topoFacei].first(); label blockNfacei = mergeBlockN[topoFacei].second(); Pair<label> fmap ( faceMap ( blockPfacei, blocks[blockPi].blockShape().faces()[blockPfacei], blockNfacei, blocks[blockNi].blockShape().faces()[blockNfacei] ) ); if (debug) { Info<< " Face map for faces " << blocks[blockPi].blockShape().faces()[blockPfacei] << " " << blocks[blockNi].blockShape().faces()[blockNfacei] << ": " << fmap << endl; } const pointField& blockPpoints = blocks[blockPi].points(); const pointField& blockNpoints = blocks[blockNi].points(); Pair<label> Pnij(faceNij(blockPfacei, blocks[blockPi])); // Check block subdivision correspondence { Pair<label> Nnij(faceNij(blockNfacei, blocks[blockNi])); Pair<label> NPnij; NPnij[0] = Nnij[mag(fmap[0]) - 1]; NPnij[1] = Nnij[mag(fmap[1]) - 1]; if (Pnij != NPnij) { FatalErrorIn("blockMesh::calcMergeInfoFast()") << "Sub-division mismatch between face " << blockPfacei << " of block " << blockPi << Pnij << " and face " << blockNfacei << " of block " << blockNi << Nnij << exit(FatalError); } } // Calculate a suitable test distance from the bounding box of the face. // Note this is used only as a sanity check and for diagnostics in // case there is a grading inconsistency. const boundBox bb(topoCells[blockPi].points(topoFaces, topoPoints)); const scalar testSqrDist = magSqr(1e-6*bb.span()); // Accumulate the maximum merge distance for diagnostics scalar maxSqrDist = 0; for (label j=0; j<Pnij.second(); j++) { for (label i=0; i<Pnij.first(); i++) { label blockPpointi = facePoint(blockPfacei, blocks[blockPi], i, j); label blockNpointi = facePointN(blockNfacei, fmap, blocks[blockNi], i, j); scalar sqrDist ( magSqr ( blockPpoints[blockPpointi] - blockNpoints[blockNpointi] ) ); if (sqrDist > testSqrDist) { FatalErrorIn("blockMesh::calcMergeInfoFast()") << "Point merge failure between face " << blockPfacei << " of block " << blockPi << " and face " << blockNfacei << " of block " << blockNi << endl << " This may be due to inconsistent grading." << exit(FatalError); } maxSqrDist = max(maxSqrDist, sqrDist); label Ppointi = blockPpointi + blockOffsets_[blockPi]; label Npointi = blockNpointi + blockOffsets_[blockNi]; label minPNi = min(Ppointi, Npointi); if (mergeList_[Ppointi] != -1) { minPNi = min(minPNi, mergeList_[Ppointi]); } if (mergeList_[Npointi] != -1) { minPNi = min(minPNi, mergeList_[Npointi]); } mergeList_[Ppointi] = mergeList_[Npointi] = minPNi; } } if (debug) { Info<< " Max distance between merge points: " << sqrt(maxSqrDist) << endl; } } bool changedPointMerge = false; label nPasses = 0; do { changedPointMerge = false; nPasses++; forAll(topoInternalFaces, topoFacei) { label blockPi = mergeBlockP[topoFacei].first(); label blockPfacei = mergeBlockP[topoFacei].second(); label blockNi = mergeBlockN[topoFacei].first(); label blockNfacei = mergeBlockN[topoFacei].second(); Pair<label> fmap ( faceMap ( blockPfacei, blocks[blockPi].blockShape().faces()[blockPfacei], blockNfacei, blocks[blockNi].blockShape().faces()[blockNfacei] ) ); Pair<label> Pnij(faceNij(blockPfacei, blocks[blockPi])); for (label j=0; j<Pnij.second(); j++) { for (label i=0; i<Pnij.first(); i++) { label blockPpointi = facePoint(blockPfacei, blocks[blockPi], i, j); label blockNpointi = facePointN(blockNfacei, fmap, blocks[blockNi], i, j); label Ppointi = blockPpointi + blockOffsets_[blockPi]; label Npointi = blockNpointi + blockOffsets_[blockNi]; if (mergeList_[Ppointi] != mergeList_[Npointi]) { changedPointMerge = true; mergeList_[Ppointi] = mergeList_[Npointi] = min(mergeList_[Ppointi], mergeList_[Npointi]); } } } } if (verboseOutput) { Info<< "." << flush; } if (nPasses > 100) { FatalErrorIn("blockMesh::calcMergeInfoFast()") << "Point merging failed after 100 passes." << exit(FatalError); } } while (changedPointMerge); // Sort merge list and count number of unique points label nUniqPoints = 0; forAll(mergeList_, pointi) { if (mergeList_[pointi] > pointi) { FatalErrorIn("blockMesh::calcMergeInfoFast()") << "Merge list contains point index out of range" << exit(FatalError); } if (mergeList_[pointi] == -1 || mergeList_[pointi] == pointi) { mergeList_[pointi] = nUniqPoints++; } else { mergeList_[pointi] = mergeList_[mergeList_[pointi]]; } } // Correct number of points in mesh nPoints_ = nUniqPoints; } // ************************************************************************* // |
|
blockMeshMergeFast.C.a15ee7de5e09.patch (4,496 bytes)
diff --git a/src/mesh/blockMesh/blockMesh/blockMeshMergeFast.C b/src/mesh/blockMesh/blockMesh/blockMeshMergeFast.C index d6d2a2d..ba93f1f 100644 --- a/src/mesh/blockMesh/blockMesh/blockMeshMergeFast.C +++ b/src/mesh/blockMesh/blockMesh/blockMeshMergeFast.C @@ -42,7 +42,7 @@ namespace Foam // ); // Face-edge directions -static const int faceEdgeDirs[6][4] = +static const label faceEdgeDirs[6][4] = { {2, 1, -2, -1}, {1, 2, -1, -2}, @@ -53,24 +53,24 @@ static const int faceEdgeDirs[6][4] = }; // The face-face-rotation direction correspondence map -static Pair<int> faceFaceRotMap[6][6][4]; +static Pair<label> faceFaceRotMap[6][6][4]; // Generate the face-face-rotation direction correspondence map void genFaceFaceRotMap() { - for(int facePi=0; facePi<6; facePi++) + for(label facePi=0; facePi<6; facePi++) { - for(int faceNi=0; faceNi<6; faceNi++) + for(label faceNi=0; faceNi<6; faceNi++) { - for(int rot=0; rot<4; rot++) + for(label rot=0; rot<4; rot++) { - Pair<int>& map = faceFaceRotMap[facePi][faceNi][rot]; + Pair<label>& map = faceFaceRotMap[facePi][faceNi][rot]; - for(int Pp=0; Pp<2; Pp++) + for(label Pp=0; Pp<2; Pp++) { - int Pdir = faceEdgeDirs[facePi][Pp]; - int Np = (3 - Pp + rot)%4; - int Ndir = faceEdgeDirs[faceNi][Np]; + label Pdir = faceEdgeDirs[facePi][Pp]; + label Np = (3 - Pp + rot)%4; + label Ndir = faceEdgeDirs[faceNi][Np]; map[Pdir-1] = -Ndir; } @@ -86,7 +86,7 @@ void genFaceFaceRotMap() } // Return the direction map for the merge-faces -Pair<int> faceMap +Pair<label> faceMap ( const label facePi, const face& faceP, @@ -95,7 +95,7 @@ Pair<int> faceMap ) { // Search for the point on faceN corresponding to the 0-point on faceP - for(int rot=0; rot<4; rot++) + for(label rot=0; rot<4; rot++) { if (faceN[rot] == faceP[0]) { @@ -111,7 +111,7 @@ Pair<int> faceMap << faceP << " and " << faceN << exit(FatalError); - return Pair<int>(0, 0); + return Pair<label>(0, 0); } // Set the block and face indices for all the merge faces @@ -161,7 +161,7 @@ Pair<label> faceNij(const label facei, const block& block) { Pair<label> fnij; - int i = facei/2; + label i = facei/2; if (i == 0) { @@ -183,7 +183,7 @@ Pair<label> faceNij(const label facei, const block& block) } // Sign the index corresponding to the map -inline label signIndex(const int map, const label i) +inline label signIndex(const label map, const label i) { return map < 0 ? -i-1 : i; } @@ -195,7 +195,7 @@ inline label unsignIndex(const label i, const label ni) } // Return the mapped index -inline label mapij(const int map, const label i, const label j) +inline label mapij(const label map, const label i, const label j) { return signIndex(map, mag(map) == 1 ? i : j); } @@ -203,7 +203,7 @@ inline label mapij(const int map, const label i, const label j) // Return the face point index inline label facePoint ( - const int facei, + const label facei, const block& block, const label i, const label j @@ -248,7 +248,7 @@ inline label facePointN // Return the neighbour face point from the mapped indices inline label facePointN ( - const int facei, + const label facei, const block& block, const label i, const label j @@ -276,8 +276,8 @@ inline label facePointN // Return the neighbour face point using the map inline label facePointN ( - const int facei, - const Pair<int>& fmap, + const label facei, + const Pair<label>& fmap, const block& block, const label i, const label j @@ -377,7 +377,7 @@ void Foam::blockMesh::calcMergeInfoFast() label blockNi = mergeBlockN[topoFacei].first(); label blockNfacei = mergeBlockN[topoFacei].second(); - Pair<int> fmap + Pair<label> fmap ( faceMap ( @@ -504,7 +504,7 @@ void Foam::blockMesh::calcMergeInfoFast() label blockNi = mergeBlockN[topoFacei].first(); label blockNfacei = mergeBlockN[topoFacei].second(); - Pair<int> fmap + Pair<label> fmap ( faceMap ( |
|
Thanks for the bug-report. I am testing 64bit label support now. I chose int rather than label for the addressing in blockMeshMergeFast.C specifically because 64bit addressing is not need for looping from 0 to 6. The ambiguity in mag is simply because mag is currently defined for label rather than int32_t and int64_t; I will rectify this. |
|
Resolved by commit 44e7fc454df95f529aeff77b11d0d6a2d1669bb0 |
Date Modified | Username | Field | Change |
---|---|---|---|
2015-08-07 14:35 | wyldckat | New Issue | |
2015-08-07 14:35 | wyldckat | File Added: blockMeshMergeFast.C | |
2015-08-07 14:36 | wyldckat | File Added: blockMeshMergeFast.C.a15ee7de5e09.patch | |
2015-08-07 14:56 | henry | Note Added: 0005205 | |
2015-08-07 15:54 | henry | Note Added: 0005207 | |
2015-08-07 15:54 | henry | Status | new => resolved |
2015-08-07 15:54 | henry | Resolution | open => fixed |
2015-08-07 15:54 | henry | Assigned To | => henry |