FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

PointEdgeWave.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 \*---------------------------------------------------------------------------*/
00025 
00026 #include "PointEdgeWave.H"
00027 #include <OpenFOAM/polyMesh.H>
00028 #include <OpenFOAM/processorPolyPatch.H>
00029 #include <OpenFOAM/cyclicPolyPatch.H>
00030 #include <OpenFOAM/OPstream.H>
00031 #include <OpenFOAM/IPstream.H>
00032 #include <OpenFOAM/PstreamCombineReduceOps.H>
00033 #include <OpenFOAM/debug.H>
00034 #include <OpenFOAM/typeInfo.H>
00035 
00036 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00037 
00038 template <class Type>
00039 Foam::scalar Foam::PointEdgeWave<Type>::propagationTol_ = 0.01;
00040 
00041 
00042 // Offset labelList. Used for transferring from one cyclic half to the other.
00043 template <class Type>
00044 void Foam::PointEdgeWave<Type>::offset(const label val, labelList& elems)
00045 {
00046     forAll(elems, i)
00047     {
00048         elems[i] += val;
00049     }
00050 }
00051 
00052 
00053 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00054 
00055 // Gets point-point correspondence. Is 
00056 // - list of halfA points (in cyclic patch points)
00057 // - list of halfB points (can overlap with A!)
00058 // - for every patchPoint its corresponding point
00059 template <class Type>
00060 void Foam::PointEdgeWave<Type>::calcCyclicAddressing()
00061 {
00062     label cycHalf = 0;
00063 
00064     forAll(mesh_.boundaryMesh(), patchI)
00065     {
00066         const polyPatch& patch = mesh_.boundaryMesh()[patchI];
00067 
00068         if (isA<cyclicPolyPatch>(patch))
00069         {
00070             label halfSize = patch.size()/2;
00071 
00072             SubList<face> halfAFaces
00073             (
00074                 mesh_.faces(),
00075                 halfSize,
00076                 patch.start()
00077             );
00078 
00079             cycHalves_.set
00080             (
00081                 cycHalf++,
00082                 new primitivePatch(halfAFaces, mesh_.points())
00083             );
00084 
00085             SubList<face> halfBFaces
00086             (
00087                 mesh_.faces(),
00088                 halfSize,
00089                 patch.start() + halfSize
00090             );
00091 
00092             cycHalves_.set
00093             (
00094                 cycHalf++,
00095                 new primitivePatch(halfBFaces, mesh_.points())
00096             );
00097         }
00098     }
00099 }
00100 
00101 
00102 // Handle leaving domain. Implementation referred to Type
00103 template <class Type>
00104 void Foam::PointEdgeWave<Type>::leaveDomain
00105 (
00106     const polyPatch& meshPatch,
00107     const primitivePatch& patch,
00108     const labelList& patchPointLabels,
00109     List<Type>& pointInfo
00110 ) const
00111 {
00112     const labelList& meshPoints = patch.meshPoints();
00113     
00114     forAll(patchPointLabels, i)
00115     {
00116         label patchPointI = patchPointLabels[i];
00117 
00118         const point& pt = patch.points()[meshPoints[patchPointI]];
00119 
00120         pointInfo[i].leaveDomain(meshPatch, patchPointI, pt);
00121     }
00122 }
00123 
00124 
00125 // Handle entering domain. Implementation referred to Type
00126 template <class Type>
00127 void Foam::PointEdgeWave<Type>::enterDomain
00128 (
00129     const polyPatch& meshPatch,
00130     const primitivePatch& patch,
00131     const labelList& patchPointLabels,
00132     List<Type>& pointInfo
00133 ) const
00134 {
00135     const labelList& meshPoints = patch.meshPoints();
00136     
00137     forAll(patchPointLabels, i)
00138     {
00139         label patchPointI = patchPointLabels[i];
00140 
00141         const point& pt = patch.points()[meshPoints[patchPointI]];
00142 
00143         pointInfo[i].enterDomain(meshPatch, patchPointI, pt);
00144     }
00145 }
00146 
00147 
00148 // Transform. Implementation referred to Type
00149 template <class Type>
00150 void Foam::PointEdgeWave<Type>::transform
00151 (
00152     const tensorField& rotTensor,
00153     List<Type>& pointInfo
00154 ) const
00155 {
00156     if (rotTensor.size() == 1)
00157     {
00158         const tensor& T = rotTensor[0];
00159 
00160         forAll(pointInfo, i)
00161         {
00162             pointInfo[i].transform(T);
00163         }
00164     }
00165     else
00166     {
00167         FatalErrorIn
00168         (
00169             "PointEdgeWave<Type>::transform(const tensorField&, List<Type>&)"
00170         )   << "Parallel cyclics not supported" << abort(FatalError);
00171     
00172         forAll(pointInfo, i)
00173         {
00174             pointInfo[i].transform(rotTensor[i]);
00175         }
00176     }
00177 }
00178 
00179 
00180 // Update info for pointI, at position pt, with information from
00181 // neighbouring edge.
00182 // Updates:
00183 //      - changedPoint_, changedPoints_, nChangedPoints_,
00184 //      - statistics: nEvals_, nUnvisitedPoints_
00185 template <class Type>
00186 bool Foam::PointEdgeWave<Type>::updatePoint
00187 (
00188     const label pointI,
00189     const label neighbourEdgeI,
00190     const Type& neighbourInfo,
00191     const scalar tol,
00192     Type& pointInfo
00193 )
00194 {
00195     nEvals_++;
00196 
00197     bool wasValid = pointInfo.valid();
00198 
00199     bool propagate =
00200         pointInfo.updatePoint
00201         (
00202             mesh_,
00203             pointI,
00204             neighbourEdgeI,
00205             neighbourInfo,
00206             tol
00207         );
00208 
00209     if (propagate)
00210     {
00211         if (!changedPoint_[pointI])
00212         {
00213             changedPoint_[pointI] = true;
00214             changedPoints_[nChangedPoints_++] = pointI;
00215         }
00216     }
00217 
00218     if (!wasValid && pointInfo.valid())
00219     {
00220         --nUnvisitedPoints_;
00221     }
00222 
00223     return propagate;
00224 }
00225 
00226 
00227 // Update info for pointI, at position pt, with information from
00228 // same point.
00229 // Updates:
00230 //      - changedPoint_, changedPoints_, nChangedPoints_,
00231 //      - statistics: nEvals_, nUnvisitedPoints_
00232 template <class Type>
00233 bool Foam::PointEdgeWave<Type>::updatePoint
00234 (
00235     const label pointI,
00236     const Type& neighbourInfo,
00237     const scalar tol,
00238     Type& pointInfo
00239 )
00240 {
00241     nEvals_++;
00242 
00243     bool wasValid = pointInfo.valid();
00244 
00245     bool propagate =
00246         pointInfo.updatePoint
00247         (
00248             mesh_,
00249             pointI,
00250             neighbourInfo,
00251             tol
00252         );
00253 
00254     if (propagate)
00255     {
00256         if (!changedPoint_[pointI])
00257         {
00258             changedPoint_[pointI] = true;
00259             changedPoints_[nChangedPoints_++] = pointI;
00260         }
00261     }
00262 
00263     if (!wasValid && pointInfo.valid())
00264     {
00265         --nUnvisitedPoints_;
00266     }
00267 
00268     return propagate;
00269 }
00270 
00271 
00272 // Update info for edgeI, at position pt, with information from
00273 // neighbouring point.
00274 // Updates:
00275 //      - changedEdge_, changedEdges_, nChangedEdges_,
00276 //      - statistics: nEvals_, nUnvisitedEdge_
00277 template <class Type>
00278 bool Foam::PointEdgeWave<Type>::updateEdge
00279 (
00280     const label edgeI,
00281     const label neighbourPointI,
00282     const Type& neighbourInfo,
00283     const scalar tol,
00284     Type& edgeInfo
00285 )
00286 {
00287     nEvals_++;
00288 
00289     bool wasValid = edgeInfo.valid();
00290 
00291     bool propagate =
00292         edgeInfo.updateEdge
00293         (
00294             mesh_,
00295             edgeI,
00296             neighbourPointI,
00297             neighbourInfo,
00298             tol
00299         );
00300 
00301     if (propagate)
00302     {
00303         if (!changedEdge_[edgeI])
00304         {
00305             changedEdge_[edgeI] = true;
00306             changedEdges_[nChangedEdges_++] = edgeI;
00307         }
00308     }
00309 
00310     if (!wasValid && edgeInfo.valid())
00311     {
00312         --nUnvisitedEdges_;
00313     }
00314 
00315     return propagate;
00316 }
00317 
00318 
00319 // Check if patches of given type name are present
00320 template <class Type>
00321 template <class PatchType>
00322 Foam::label Foam::PointEdgeWave<Type>::countPatchType() const
00323 {
00324     label nPatches = 0;
00325 
00326     forAll(mesh_.boundaryMesh(), patchI)
00327     {
00328         if (isA<PatchType>(mesh_.boundaryMesh()[patchI]))
00329         {
00330             nPatches++;
00331         }
00332     }
00333     return nPatches;
00334 }
00335 
00336 
00337 // Collect changed patch points
00338 template <class Type>
00339 void Foam::PointEdgeWave<Type>::getChangedPatchPoints
00340 (
00341     const primitivePatch& patch,
00342 
00343     DynamicList<Type>& patchInfo,
00344     DynamicList<label>& patchPoints,
00345     DynamicList<label>& owner,
00346     DynamicList<label>& ownerIndex
00347 ) const
00348 {
00349     const labelList& meshPoints = patch.meshPoints();
00350     const faceList& localFaces = patch.localFaces();
00351     const labelListList& pointFaces = patch.pointFaces();
00352 
00353     forAll(meshPoints, patchPointI)
00354     {
00355         label meshPointI = meshPoints[patchPointI];
00356 
00357         if (changedPoint_[meshPointI])
00358         {
00359             patchInfo.append(allPointInfo_[meshPointI]);
00360             patchPoints.append(patchPointI);
00361 
00362             label patchFaceI = pointFaces[patchPointI][0];
00363 
00364             const face& f = localFaces[patchFaceI];
00365 
00366             label index = findIndex(f, patchPointI);
00367 
00368             owner.append(patchFaceI);
00369             ownerIndex.append(index);
00370         }
00371     }
00372 
00373     patchInfo.shrink();
00374     patchPoints.shrink();
00375     owner.shrink();
00376     ownerIndex.shrink();
00377 }
00378 
00379 
00380 // Update overall for changed patch points
00381 template <class Type>
00382 void Foam::PointEdgeWave<Type>::updateFromPatchInfo
00383 (
00384     const polyPatch& meshPatch,
00385     const primitivePatch& patch,
00386     const labelList& owner,
00387     const labelList& ownerIndex,
00388     List<Type>& patchInfo
00389 )
00390 {
00391     const faceList& localFaces = patch.localFaces();
00392     const labelList& meshPoints = patch.meshPoints();
00393 
00394     // Get patch and mesh points.
00395     labelList changedPatchPoints(patchInfo.size());
00396     labelList changedMeshPoints(patchInfo.size());
00397 
00398     forAll(owner, i)
00399     {
00400         label faceI = owner[i];
00401 
00402         const face& f = localFaces[faceI];
00403 
00404         label index = (f.size() - ownerIndex[i]) % f.size();
00405 
00406         changedPatchPoints[i] = f[index];
00407         changedMeshPoints[i] = meshPoints[f[index]];
00408     }
00409 
00410     // Adapt for entering domain
00411     enterDomain(meshPatch, patch, changedPatchPoints, patchInfo);
00412 
00413     // Merge received info
00414     forAll(patchInfo, i)
00415     {
00416         updatePoint
00417         (
00418             changedMeshPoints[i],
00419             patchInfo[i],
00420             propagationTol_,
00421             allPointInfo_[changedMeshPoints[i]]
00422         );
00423     }
00424 }
00425 
00426 
00427 
00428 // Transfer all the information to/from neighbouring processors
00429 template <class Type>
00430 void Foam::PointEdgeWave<Type>::handleProcPatches()
00431 {
00432     // 1. Send all point info on processor patches. Send as
00433     // face label + offset in face.
00434 
00435     forAll(mesh_.boundaryMesh(), patchI)
00436     {
00437         const polyPatch& patch = mesh_.boundaryMesh()[patchI];
00438 
00439         if (isA<processorPolyPatch>(patch))
00440         {
00441             // Get all changed points in relative addressing
00442 
00443             DynamicList<Type> patchInfo(patch.nPoints());
00444             DynamicList<label> patchPoints(patch.nPoints());
00445             DynamicList<label> owner(patch.nPoints());
00446             DynamicList<label> ownerIndex(patch.nPoints());
00447 
00448             getChangedPatchPoints
00449             (
00450                 patch,
00451                 patchInfo,
00452                 patchPoints,
00453                 owner,
00454                 ownerIndex
00455             );
00456 
00457             // Adapt for leaving domain
00458             leaveDomain(patch, patch, patchPoints, patchInfo);
00459 
00460             const processorPolyPatch& procPatch =
00461                 refCast<const processorPolyPatch>(patch);
00462 
00463             if (debug)
00464             {
00465                 Pout<< "Processor patch " << patchI << ' ' << patch.name()
00466                     << " communicating with " << procPatch.neighbProcNo()
00467                     << "  Sending:" << patchInfo.size() << endl;
00468             }
00469 
00470             {
00471                 OPstream toNeighbour
00472                 (   
00473                     Pstream::blocking,
00474                     procPatch.neighbProcNo()
00475                 );
00476 
00477                 toNeighbour << owner << ownerIndex << patchInfo;
00478             }
00479         }
00480     }
00481 
00482 
00483     //
00484     // 2. Receive all point info on processor patches.
00485     //
00486 
00487     forAll(mesh_.boundaryMesh(), patchI)
00488     {
00489         const polyPatch& patch = mesh_.boundaryMesh()[patchI];
00490 
00491         if (isA<processorPolyPatch>(patch))
00492         {
00493             const processorPolyPatch& procPatch =
00494                 refCast<const processorPolyPatch>(patch);
00495 
00496             List<Type> patchInfo;
00497             labelList owner;
00498             labelList ownerIndex;
00499             {
00500                 IPstream fromNeighbour
00501                 (
00502                     Pstream::blocking,
00503                     procPatch.neighbProcNo()
00504                 );
00505 
00506                 fromNeighbour >> owner >> ownerIndex >> patchInfo;
00507             }    
00508 
00509             if (debug)
00510             {
00511                 Pout<< "Processor patch " << patchI << ' ' << patch.name()
00512                     << " communicating with " << procPatch.neighbProcNo()
00513                     << "  Received:" << patchInfo.size() << endl;
00514             }
00515 
00516             // Apply transform to received data for non-parallel planes
00517             if (!procPatch.parallel())
00518             {
00519                 transform(procPatch.forwardT(), patchInfo);
00520             }
00521 
00522             updateFromPatchInfo
00523             (
00524                 patch,
00525                 patch,
00526                 owner,
00527                 ownerIndex,
00528                 patchInfo
00529             );
00530         }
00531     }
00532 
00533 
00534 
00535     //
00536     // 3. Handle all shared points
00537     //    (Note:irrespective if changed or not for now)
00538     //
00539 
00540     const globalMeshData& pd = mesh_.globalData();
00541 
00542     List<Type> sharedData(pd.nGlobalPoints());
00543 
00544     forAll(pd.sharedPointLabels(), i)
00545     {
00546         label meshPointI = pd.sharedPointLabels()[i];
00547 
00548         // Fill my entries in the shared points
00549         sharedData[pd.sharedPointAddr()[i]] = allPointInfo_[meshPointI];
00550     }
00551 
00552     // Combine on master. Reduce operator has to handle a list and call
00553     // Type.updatePoint for all elements
00554     combineReduce(sharedData, listUpdateOp<Type>());
00555 
00556     forAll(pd.sharedPointLabels(), i)
00557     {
00558         label meshPointI = pd.sharedPointLabels()[i];
00559 
00560         // Retrieve my entries from the shared points
00561         updatePoint
00562         (
00563             meshPointI,
00564             sharedData[pd.sharedPointAddr()[i]],
00565             propagationTol_,
00566             allPointInfo_[meshPointI]
00567         );
00568     }
00569 }
00570 
00571 
00572 template <class Type>
00573 void Foam::PointEdgeWave<Type>::handleCyclicPatches()
00574 {
00575     // 1. Send all point info on cyclic patches. Send as
00576     // face label + offset in face.
00577 
00578     label cycHalf = 0;
00579 
00580     forAll(mesh_.boundaryMesh(), patchI)
00581     {
00582         const polyPatch& patch = mesh_.boundaryMesh()[patchI];
00583 
00584         if (isA<cyclicPolyPatch>(patch))
00585         {
00586             const primitivePatch& halfA = cycHalves_[cycHalf++];
00587             const primitivePatch& halfB = cycHalves_[cycHalf++];
00588 
00589             // HalfA : get all changed points in relative addressing
00590 
00591             DynamicList<Type> halfAInfo(halfA.nPoints());
00592             DynamicList<label> halfAPoints(halfA.nPoints());
00593             DynamicList<label> halfAOwner(halfA.nPoints());
00594             DynamicList<label> halfAIndex(halfA.nPoints());
00595 
00596             getChangedPatchPoints
00597             (
00598                 halfA,
00599                 halfAInfo,
00600                 halfAPoints,
00601                 halfAOwner,
00602                 halfAIndex
00603             );
00604 
00605             // HalfB : get all changed points in relative addressing
00606 
00607             DynamicList<Type> halfBInfo(halfB.nPoints());
00608             DynamicList<label> halfBPoints(halfB.nPoints());
00609             DynamicList<label> halfBOwner(halfB.nPoints());
00610             DynamicList<label> halfBIndex(halfB.nPoints());
00611 
00612             getChangedPatchPoints
00613             (
00614                 halfB,
00615                 halfBInfo,
00616                 halfBPoints,
00617                 halfBOwner,
00618                 halfBIndex
00619             );
00620 
00621 
00622             // HalfA : adapt for leaving domain
00623             leaveDomain(patch, halfA, halfAPoints, halfAInfo);
00624 
00625             // HalfB : adapt for leaving domain
00626             leaveDomain(patch, halfB, halfBPoints, halfBInfo);
00627 
00628 
00629             // Apply rotation for non-parallel planes
00630             const cyclicPolyPatch& cycPatch =
00631                 refCast<const cyclicPolyPatch>(patch);
00632 
00633             if (!cycPatch.parallel())
00634             {
00635                 // received data from half1
00636                 transform(cycPatch.reverseT(), halfAInfo);
00637 
00638                 // received data from half2
00639                 transform(cycPatch.forwardT(), halfBInfo);
00640             }
00641 
00642             if (debug)
00643             {
00644                 Pout<< "Cyclic patch " << patchI << ' ' << patch.name() 
00645                     << "  Changed on first half : " << halfAInfo.size()
00646                     << "  Changed on second half : " << halfBInfo.size()
00647                     << endl;
00648             }
00649 
00650             // Half1: update with data from halfB
00651             updateFromPatchInfo
00652             (
00653                 patch,
00654                 halfA,
00655                 halfBOwner,
00656                 halfBIndex,
00657                 halfBInfo
00658             );
00659 
00660             // Half2: update with data from halfA
00661             updateFromPatchInfo
00662             (
00663                 patch,
00664                 halfB,
00665                 halfAOwner,
00666                 halfAIndex,
00667                 halfAInfo
00668             );
00669 
00670             if (debug)
00671             {
00672                 //checkCyclic(patch);
00673             }
00674         }
00675     }
00676 }
00677 
00678 
00679 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00680 
00681 // Iterate, propagating changedPointsInfo across mesh, until no change (or 
00682 // maxIter reached). Initial point values specified.
00683 template <class Type>
00684 Foam::PointEdgeWave<Type>::PointEdgeWave
00685 (
00686     const polyMesh& mesh,
00687     const labelList& changedPoints,
00688     const List<Type>& changedPointsInfo,
00689 
00690     List<Type>& allPointInfo,
00691     List<Type>& allEdgeInfo,
00692 
00693     const label maxIter
00694 )
00695 :
00696     mesh_(mesh),
00697     allPointInfo_(allPointInfo),
00698     allEdgeInfo_(allEdgeInfo),
00699     changedPoint_(mesh_.nPoints(), false),
00700     changedPoints_(mesh_.nPoints()),
00701     nChangedPoints_(0),
00702     changedEdge_(mesh_.nEdges(), false),
00703     changedEdges_(mesh_.nEdges()),
00704     nChangedEdges_(0),
00705     nCyclicPatches_(countPatchType<cyclicPolyPatch>()),
00706     cycHalves_(2*nCyclicPatches_),
00707     nEvals_(0),
00708     nUnvisitedPoints_(mesh_.nPoints()),
00709     nUnvisitedEdges_(mesh_.nEdges())
00710 {
00711     if (allPointInfo_.size() != mesh_.nPoints())
00712     {
00713         FatalErrorIn
00714         (
00715             "PointEdgeWave<Type>::PointEdgeWave"
00716             "(const polyMesh&, const labelList&, const List<Type>,"
00717             " List<Type>&, List<Type>&, const label maxIter)"
00718         )   << "size of pointInfo work array is not equal to the number"
00719             << " of points in the mesh" << endl
00720             << "    pointInfo   :" << allPointInfo_.size() << endl
00721             << "    mesh.nPoints:" << mesh_.nPoints()
00722             << exit(FatalError);
00723     }
00724     if (allEdgeInfo_.size() != mesh_.nEdges())
00725     {
00726         FatalErrorIn
00727         (
00728             "PointEdgeWave<Type>::PointEdgeWave"
00729             "(const polyMesh&, const labelList&, const List<Type>,"
00730             " List<Type>&, List<Type>&, const label maxIter)"
00731         )   << "size of edgeInfo work array is not equal to the number"
00732             << " of edges in the mesh" << endl
00733             << "    edgeInfo   :" << allEdgeInfo_.size() << endl
00734             << "    mesh.nEdges:" << mesh_.nEdges()
00735             << exit(FatalError);
00736     }
00737 
00738 
00739     // Calculate cyclic halves addressing.
00740     if (nCyclicPatches_ > 0)
00741     {
00742         calcCyclicAddressing();
00743     }
00744 
00745 
00746     // Set from initial changed points data
00747     setPointInfo(changedPoints, changedPointsInfo);
00748 
00749     if (debug)
00750     {
00751         Pout<< "Seed points               : " << nChangedPoints_ << endl;
00752     }
00753 
00754     // Iterate until nothing changes
00755     label iter = iterate(maxIter);
00756 
00757     if ((maxIter > 0) && (iter >= maxIter))
00758     {
00759         FatalErrorIn
00760         (
00761             "PointEdgeWave<Type>::PointEdgeWave"
00762             "(const polyMesh&, const labelList&, const List<Type>,"
00763             " List<Type>&, List<Type>&, const label maxIter)"
00764         )   << "Maximum number of iterations reached. Increase maxIter." << endl
00765             << "    maxIter:" << maxIter << endl
00766             << "    nChangedPoints:" << nChangedPoints_ << endl
00767             << "    nChangedEdges:" << nChangedEdges_ << endl
00768             << exit(FatalError);
00769     }
00770 }
00771 
00772 
00773 // * * * * * * * * * * * * * * * * Destructor  * * * * * * * * * * * * * * * //
00774 
00775 template <class Type>
00776 Foam::PointEdgeWave<Type>::~PointEdgeWave()
00777 {}
00778 
00779 
00780 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00781 
00782 
00783 template <class Type>
00784 Foam::label Foam::PointEdgeWave<Type>::getUnsetPoints() const
00785 {
00786     return nUnvisitedPoints_;
00787 }
00788 
00789 
00790 template <class Type>
00791 Foam::label Foam::PointEdgeWave<Type>::getUnsetEdges() const
00792 {
00793     return nUnvisitedEdges_;
00794 }
00795 
00796 
00797 // Copy point information into member data
00798 template <class Type>
00799 void Foam::PointEdgeWave<Type>::setPointInfo
00800 (
00801     const labelList& changedPoints,
00802     const List<Type>& changedPointsInfo
00803 )
00804 {
00805     forAll(changedPoints, changedPointI)
00806     {
00807         label pointI = changedPoints[changedPointI];
00808 
00809         bool wasValid = allPointInfo_[pointI].valid();
00810 
00811         // Copy info for pointI
00812         allPointInfo_[pointI] = changedPointsInfo[changedPointI];
00813 
00814         // Maintain count of unset points
00815         if (!wasValid && allPointInfo_[pointI].valid())
00816         {
00817             --nUnvisitedPoints_;
00818         }
00819 
00820         // Mark pointI as changed, both on list and on point itself.
00821 
00822         if (!changedPoint_[pointI])
00823         {
00824             changedPoint_[pointI] = true;
00825             changedPoints_[nChangedPoints_++] = pointI;
00826         }
00827     }
00828 }
00829 
00830 
00831 // Propagate information from edge to point. Return number of points changed.
00832 template <class Type>
00833 Foam::label Foam::PointEdgeWave<Type>::edgeToPoint()
00834 {
00835     for
00836     (
00837         label changedEdgeI = 0;
00838         changedEdgeI < nChangedEdges_;
00839         changedEdgeI++
00840     )
00841     {
00842         label edgeI = changedEdges_[changedEdgeI];
00843 
00844         if (!changedEdge_[edgeI])
00845         {
00846             FatalErrorIn("PointEdgeWave<Type>::edgeToPoint()")
00847                 << "edge " << edgeI
00848                 << " not marked as having been changed" << nl
00849                 << "This might be caused by multiple occurences of the same"
00850                 << " seed point." << abort(FatalError);
00851         }
00852 
00853 
00854         const Type& neighbourWallInfo = allEdgeInfo_[edgeI];
00855 
00856         // Evaluate all connected points (= edge endpoints)
00857         const edge& e = mesh_.edges()[edgeI];
00858 
00859         forAll(e, eI)
00860         {
00861             Type& currentWallInfo = allPointInfo_[e[eI]];
00862 
00863             if (currentWallInfo != neighbourWallInfo)
00864             {
00865                 updatePoint
00866                 (
00867                     e[eI],
00868                     edgeI,
00869                     neighbourWallInfo,
00870                     propagationTol_,
00871                     currentWallInfo
00872                 );
00873             }
00874         }
00875 
00876         // Reset status of edge
00877         changedEdge_[edgeI] = false;
00878     }
00879 
00880     // Handled all changed edges by now
00881     nChangedEdges_ = 0;
00882 
00883     if (nCyclicPatches_ > 0)
00884     {
00885         // Transfer changed points across cyclic halves
00886         handleCyclicPatches();
00887     }
00888     if (Pstream::parRun())
00889     {
00890         // Transfer changed points from neighbouring processors.
00891         handleProcPatches();
00892     }
00893 
00894     if (debug)
00895     {
00896         Pout<< "Changed points            : " << nChangedPoints_ << endl;
00897     }
00898 
00899     // Sum nChangedPoints over all procs
00900     label totNChanged = nChangedPoints_;
00901 
00902     reduce(totNChanged, sumOp<label>());
00903 
00904     return totNChanged;
00905 }
00906 
00907 
00908 // Propagate information from point to edge. Return number of edges changed.
00909 template <class Type>
00910 Foam::label Foam::PointEdgeWave<Type>::pointToEdge()
00911 {
00912     const labelListList& pointEdges = mesh_.pointEdges();
00913 
00914     for
00915     (
00916         label changedPointI = 0;
00917         changedPointI < nChangedPoints_;
00918         changedPointI++
00919     )
00920     {
00921         label pointI = changedPoints_[changedPointI];
00922 
00923         if (!changedPoint_[pointI])
00924         {
00925             FatalErrorIn("PointEdgeWave<Type>::pointToEdge()")
00926                 << "Point " << pointI
00927                 << " not marked as having been changed" << nl
00928                 << "This might be caused by multiple occurences of the same"
00929                 << " seed point." << abort(FatalError);
00930         }
00931 
00932         const Type& neighbourWallInfo = allPointInfo_[pointI];
00933 
00934         // Evaluate all connected edges
00935 
00936         const labelList& edgeLabels = pointEdges[pointI];
00937         forAll(edgeLabels, edgeLabelI)
00938         {
00939             label edgeI = edgeLabels[edgeLabelI];
00940 
00941             Type& currentWallInfo = allEdgeInfo_[edgeI];
00942 
00943             if (currentWallInfo != neighbourWallInfo)
00944             {
00945                 updateEdge
00946                 (
00947                     edgeI,
00948                     pointI,
00949                     neighbourWallInfo,
00950                     propagationTol_,
00951                     currentWallInfo
00952                 );
00953             }
00954         }
00955 
00956         // Reset status of point
00957         changedPoint_[pointI] = false;
00958     }
00959 
00960     // Handled all changed points by now
00961     nChangedPoints_ = 0;
00962 
00963     if (debug)
00964     {
00965         Pout<< "Changed edges             : " << nChangedEdges_ << endl;
00966     }
00967 
00968     // Sum nChangedPoints over all procs
00969     label totNChanged = nChangedEdges_;
00970 
00971     reduce(totNChanged, sumOp<label>());
00972 
00973     return totNChanged;
00974 }
00975 
00976 
00977 // Iterate
00978 template <class Type>
00979 Foam::label Foam::PointEdgeWave<Type>::iterate(const label maxIter)
00980 {
00981     if (nCyclicPatches_ > 0)
00982     {
00983         // Transfer changed points across cyclic halves
00984         handleCyclicPatches();
00985     }
00986     if (Pstream::parRun())
00987     {
00988         // Transfer changed points from neighbouring processors.
00989         handleProcPatches();
00990     }
00991 
00992     nEvals_ = 0;
00993 
00994     label iter = 0;
00995 
00996     while(iter < maxIter)
00997     {
00998         if (debug)
00999         {
01000             Pout<< "Iteration " << iter << endl;
01001         }
01002 
01003         label nEdges = pointToEdge();
01004 
01005         if (debug)
01006         {
01007             Pout<< "Total changed edges       : " << nEdges << endl;
01008         }
01009 
01010         if (nEdges == 0)
01011         {
01012             break;
01013         }
01014 
01015         label nPoints = edgeToPoint();
01016 
01017         if (debug)
01018         {
01019             Pout<< "Total changed points      : " << nPoints << endl;
01020 
01021             Pout<< "Total evaluations         : " << nEvals_ << endl;
01022 
01023             Pout<< "Remaining unvisited points: " << nUnvisitedPoints_ << endl;
01024 
01025             Pout<< "Remaining unvisited edges : " << nUnvisitedEdges_ << endl;
01026 
01027             Pout<< endl;
01028         }
01029 
01030         if (nPoints == 0)
01031         {
01032             break;
01033         }
01034 
01035         iter++;
01036     }
01037 
01038     return iter;
01039 }
01040 
01041 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines