/*---------------------------------------------------------------------------*\ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | \\ / A nd | Copyright (C) 2004-2010 OpenCFD Ltd. \\/ M anipulation | ------------------------------------------------------------------------------- | Copyright (C) 2011-2017 OpenFOAM Foundation ------------------------------------------------------------------------------- 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 . \*---------------------------------------------------------------------------*/ #include "PointEdgeWave.H" #include "polyMesh.H" #include "processorPolyPatch.H" #include "cyclicPolyPatch.H" #include "OPstream.H" #include "IPstream.H" #include "PstreamCombineReduceOps.H" #include "debug.H" #include "typeInfo.H" #include "globalMeshData.H" #include "pointFields.H" // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * // template Foam::scalar Foam::PointEdgeWave::propagationTol_ = 0.01; template int Foam::PointEdgeWave::dummyTrackData_ = 12345; namespace Foam { //- Reduction class. If x and y are not equal assign value. template class combineEqOp { TrackingData& td_; public: combineEqOp(TrackingData& td) : td_(td) {} void operator()(Type& x, const Type& y) const { if (!x.valid(td_) && y.valid(td_)) { x = y; } } }; } // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * // // Handle leaving domain. Implementation referred to Type template void Foam::PointEdgeWave::leaveDomain ( const polyPatch& patch, const labelList& patchPointLabels, List& pointInfo ) const { const labelList& meshPoints = patch.meshPoints(); forAll(patchPointLabels, i) { label patchPointi = patchPointLabels[i]; const point& pt = patch.points()[meshPoints[patchPointi]]; pointInfo[i].leaveDomain(patch, patchPointi, pt, td_); } } // Handle entering domain. Implementation referred to Type template void Foam::PointEdgeWave::enterDomain ( const polyPatch& patch, const labelList& patchPointLabels, List& pointInfo ) const { const labelList& meshPoints = patch.meshPoints(); forAll(patchPointLabels, i) { label patchPointi = patchPointLabels[i]; const point& pt = patch.points()[meshPoints[patchPointi]]; pointInfo[i].enterDomain(patch, patchPointi, pt, td_); } } // Transform. Implementation referred to Type template void Foam::PointEdgeWave::transform ( const polyPatch& patch, const tensorField& rotTensor, List& pointInfo ) const { if (rotTensor.size() == 1) { const tensor& T = rotTensor[0]; forAll(pointInfo, i) { pointInfo[i].transform(T, td_); } } else { FatalErrorInFunction << "Non-uniform transformation on patch " << patch.name() << " of type " << patch.type() << " not supported for point fields" << abort(FatalError); forAll(pointInfo, i) { pointInfo[i].transform(rotTensor[i], td_); } } } // Update info for pointi, at position pt, with information from // neighbouring edge. // Updates: // - changedPoint_, changedPoints_, nChangedPoints_, // - statistics: nEvals_, nUnvisitedPoints_ template bool Foam::PointEdgeWave::updatePoint ( const label pointi, const label neighbourEdgeI, const Type& neighbourInfo, Type& pointInfo ) { nEvals_++; bool wasValid = pointInfo.valid(td_); bool propagate = pointInfo.updatePoint ( mesh_, pointi, neighbourEdgeI, neighbourInfo, propagationTol_, td_ ); if (propagate) { if (!changedPoint_[pointi]) { changedPoint_[pointi] = true; changedPoints_[nChangedPoints_++] = pointi; } } if (!wasValid && pointInfo.valid(td_)) { --nUnvisitedPoints_; } return propagate; } // Update info for pointi, at position pt, with information from // same point. // Updates: // - changedPoint_, changedPoints_, nChangedPoints_, // - statistics: nEvals_, nUnvisitedPoints_ template bool Foam::PointEdgeWave::updatePoint ( const label pointi, const Type& neighbourInfo, Type& pointInfo ) { nEvals_++; bool wasValid = pointInfo.valid(td_); bool propagate = pointInfo.updatePoint ( mesh_, pointi, neighbourInfo, propagationTol_, td_ ); if (propagate) { if (!changedPoint_[pointi]) { changedPoint_[pointi] = true; changedPoints_[nChangedPoints_++] = pointi; } } if (!wasValid && pointInfo.valid(td_)) { --nUnvisitedPoints_; } return propagate; } // Update info for edgeI, at position pt, with information from // neighbouring point. // Updates: // - changedEdge_, changedEdges_, nChangedEdges_, // - statistics: nEvals_, nUnvisitedEdge_ template bool Foam::PointEdgeWave::updateEdge ( const label edgeI, const label neighbourPointi, const Type& neighbourInfo, Type& edgeInfo ) { nEvals_++; bool wasValid = edgeInfo.valid(td_); bool propagate = edgeInfo.updateEdge ( mesh_, edgeI, neighbourPointi, neighbourInfo, propagationTol_, td_ ); if (propagate) { if (!changedEdge_[edgeI]) { changedEdge_[edgeI] = true; changedEdges_[nChangedEdges_++] = edgeI; } } if (!wasValid && edgeInfo.valid(td_)) { --nUnvisitedEdges_; } return propagate; } // Check if patches of given type name are present template template Foam::label Foam::PointEdgeWave::countPatchType() const { label nPatches = 0; forAll(mesh_.boundaryMesh(), patchi) { if (isA(mesh_.boundaryMesh()[patchi])) { nPatches++; } } return nPatches; } // Transfer all the information to/from neighbouring processors template void Foam::PointEdgeWave::handleProcPatches() { // 1. Send all point info on processor patches. PstreamBuffers pBufs(Pstream::commsTypes::nonBlocking); DynamicList patchInfo; DynamicList