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

hierarchGeomDecomp.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 "hierarchGeomDecomp.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <OpenFOAM/PstreamReduceOps.H>
00029 #include <OpenFOAM/SortableList.H>
00030 
00031 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00032 
00033 namespace Foam
00034 {
00035     defineTypeNameAndDebug(hierarchGeomDecomp, 0);
00036 
00037     addToRunTimeSelectionTable
00038     (
00039         decompositionMethod,
00040         hierarchGeomDecomp,
00041         dictionary
00042     );
00043 
00044     addToRunTimeSelectionTable
00045     (
00046         decompositionMethod,
00047         hierarchGeomDecomp,
00048         dictionaryMesh
00049     );
00050 }
00051 
00052 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00053 
00054 void Foam::hierarchGeomDecomp::setDecompOrder()
00055 {
00056     word order(geomDecomDict_.lookup("order"));
00057 
00058     if (order.size() != 3)
00059     {
00060         FatalIOErrorIn
00061         (
00062             "hierarchGeomDecomp::hierarchGeomDecomp"
00063             "(const dictionary& decompositionDict)",
00064             decompositionDict_
00065         )   << "number of characters in order (" << order << ") != 3"
00066             << exit(FatalIOError);
00067     }
00068 
00069     for (label i = 0; i < 3; i++)
00070     {
00071         if (order[i] == 'x')
00072         {
00073             decompOrder_[i] = 0;
00074         }
00075         else if (order[i] == 'y')
00076         {
00077             decompOrder_[i] = 1;
00078         }
00079         else if (order[i] == 'z')
00080         {
00081             decompOrder_[i] = 2;
00082         }
00083         else
00084         {
00085             FatalIOErrorIn
00086             (
00087                 "hierarchGeomDecomp::hierarchGeomDecomp"
00088                 "(const dictionary& decompositionDict)",
00089                 decompositionDict_
00090             )   << "Illegal decomposition order " << order << endl
00091                 << "It should only contain x, y or z" << exit(FatalError);
00092         }
00093     }
00094 }
00095 
00096 
00097 Foam::label Foam::hierarchGeomDecomp::findLower
00098 (
00099     const List<scalar>& l,
00100     const scalar t,
00101     const label initLow,
00102     const label initHigh
00103 )
00104 {
00105     if (initHigh <= initLow)
00106     {
00107         return initLow;
00108     }
00109 
00110     label low = initLow;
00111     label high = initHigh;
00112 
00113     while ((high - low) > 1)
00114     {
00115         label mid = (low + high)/2;
00116 
00117         if (l[mid] < t)
00118         {
00119             low = mid;
00120         }
00121         else
00122         {
00123             high = mid;
00124         }
00125     }
00126 
00127     // high and low can still differ by one. Choose best.
00128 
00129     label tIndex = -1;
00130 
00131     if (l[high-1] < t)
00132     {
00133         tIndex = high;
00134     }
00135     else
00136     {
00137         tIndex = low;
00138     }
00139 
00140     return tIndex;
00141 }
00142 
00143 
00144 // Create a mapping between the index and the weighted size.
00145 // For convenience, sortedWeightedSize is one size bigger than current. This
00146 // avoids extra tests.
00147 void Foam::hierarchGeomDecomp::calculateSortedWeightedSizes
00148 (
00149     const labelList& current,
00150     const labelList& indices,
00151     const scalarField& weights,
00152     const label globalCurrentSize,
00153 
00154     scalarField& sortedWeightedSizes
00155 )
00156 {
00157     // Evaluate cumulative weights.
00158     sortedWeightedSizes[0] = 0;
00159     forAll(current, i)
00160     {
00161         label pointI = current[indices[i]];
00162         sortedWeightedSizes[i + 1] = sortedWeightedSizes[i] + weights[pointI];
00163     }
00164     // Non-dimensionalise and multiply by size.
00165     scalar globalCurrentLength = returnReduce
00166     (
00167         sortedWeightedSizes[current.size()],
00168         sumOp<scalar>()
00169     );
00170     // Normalise weights by global sum of weights and multiply through
00171     // by global size.
00172     sortedWeightedSizes *= (globalCurrentSize/globalCurrentLength);
00173 }
00174 
00175 
00176 // Find position in values so between minIndex and this position there
00177 // are wantedSize elements.
00178 void Foam::hierarchGeomDecomp::findBinary
00179 (
00180     const label sizeTol,
00181     const List<scalar>& values,
00182     const label minIndex,       // index of previous value
00183     const scalar minValue,      // value at minIndex
00184     const scalar maxValue,      // global max of values
00185     const scalar wantedSize,    // wanted size
00186 
00187     label& mid,                 // index where size of bin is
00188                                 // wantedSize (to within sizeTol)
00189     scalar& midValue            // value at mid
00190 )
00191 {
00192     label low = minIndex;
00193     scalar lowValue = minValue;
00194 
00195     scalar highValue = maxValue;
00196     // (one beyond) index of highValue
00197     label high = values.size();
00198 
00199     // Safeguards to avoid infinite loop.
00200     scalar midValuePrev = VGREAT;
00201 
00202     while (true)
00203     {
00204         label size = returnReduce(mid-minIndex, sumOp<label>());
00205 
00206         if (debug)
00207         {
00208             Pout<< "    low:" << low << " lowValue:" << lowValue
00209                 << " high:" << high << " highValue:" << highValue
00210                 << " mid:" << mid << " midValue:" << midValue << endl
00211                 << "    globalSize:" << size << " wantedSize:" << wantedSize
00212                 << " sizeTol:" << sizeTol << endl;
00213         }
00214 
00215         if (wantedSize < size - sizeTol)
00216         {
00217             high = mid;
00218             highValue = midValue;
00219         }
00220         else if (wantedSize > size + sizeTol)
00221         {
00222             low = mid;
00223             lowValue = midValue;
00224         }
00225         else
00226         {
00227             break;
00228         }
00229 
00230         // Update mid, midValue
00231         midValue = 0.5*(lowValue+highValue);
00232         mid = findLower(values, midValue, low, high);
00233 
00234         // Safeguard if same as previous.
00235         bool hasNotChanged = (mag(midValue-midValuePrev) < SMALL);
00236 
00237         if (returnReduce(hasNotChanged, andOp<bool>()))
00238         {
00239             WarningIn("hierarchGeomDecomp::findBinary(..)")
00240                 << "unable to find desired deomposition split, making do!"
00241                 << endl;
00242             break;
00243         }
00244 
00245         midValuePrev = midValue;
00246     }
00247 }
00248 
00249 
00250 // Find position in values so between minIndex and this position there
00251 // are wantedSize elements.
00252 void Foam::hierarchGeomDecomp::findBinary
00253 (
00254     const label sizeTol,
00255     const List<scalar>& sortedWeightedSizes,
00256     const List<scalar>& values,
00257     const label minIndex,       // index of previous value
00258     const scalar minValue,      // value at minIndex
00259     const scalar maxValue,      // global max of values
00260     const scalar wantedSize,    // wanted size
00261 
00262     label& mid,                 // index where size of bin is
00263                                 // wantedSize (to within sizeTol)
00264     scalar& midValue            // value at mid
00265 )
00266 {
00267     label low = minIndex;
00268     scalar lowValue = minValue;
00269 
00270     scalar highValue = maxValue;
00271     // (one beyond) index of highValue
00272     label high = values.size();
00273 
00274     // Safeguards to avoid infinite loop.
00275     scalar midValuePrev = VGREAT;
00276 
00277     while (true)
00278     {
00279         scalar weightedSize = returnReduce
00280         (
00281             sortedWeightedSizes[mid] - sortedWeightedSizes[minIndex],
00282             sumOp<scalar>()
00283         );
00284 
00285         if (debug)
00286         {
00287             Pout<< "    low:" << low << " lowValue:" << lowValue
00288                 << " high:" << high << " highValue:" << highValue
00289                 << " mid:" << mid << " midValue:" << midValue << endl
00290                 << "    globalSize:" << weightedSize
00291                 << " wantedSize:" << wantedSize
00292                 << " sizeTol:" << sizeTol << endl;
00293         }
00294 
00295         if (wantedSize < weightedSize - sizeTol)
00296         {
00297             high = mid;
00298             highValue = midValue;
00299         }
00300         else if (wantedSize > weightedSize + sizeTol)
00301         {
00302             low = mid;
00303             lowValue = midValue;
00304         }
00305         else
00306         {
00307             break;
00308         }
00309 
00310         // Update mid, midValue
00311         midValue = 0.5*(lowValue+highValue);
00312         mid = findLower(values, midValue, low, high);
00313 
00314         // Safeguard if same as previous.
00315         bool hasNotChanged = (mag(midValue-midValuePrev) < SMALL);
00316 
00317         if (returnReduce(hasNotChanged, andOp<bool>()))
00318         {
00319             WarningIn("hierarchGeomDecomp::findBinary(..)")
00320                 << "unable to find desired deomposition split, making do!"
00321                 << endl;
00322             break;
00323         }
00324 
00325         midValuePrev = midValue;
00326     }
00327 }
00328 
00329 
00330 // Sort points into bins according to one component. Recurses to next component.
00331 void Foam::hierarchGeomDecomp::sortComponent
00332 (
00333     const label sizeTol,
00334     const pointField& points,
00335     const labelList& current,       // slice of points to decompose
00336     const direction componentIndex, // index in decompOrder_
00337     const label mult,               // multiplication factor for finalDecomp
00338     labelList& finalDecomp
00339 )
00340 {
00341     // Current component
00342     label compI = decompOrder_[componentIndex];
00343 
00344     if (debug)
00345     {
00346         Pout<< "sortComponent : Sorting slice of size " << current.size()
00347             << " in component " << compI << endl;
00348     }
00349 
00350     // Storage for sorted component compI
00351     SortableList<scalar> sortedCoord(current.size());
00352 
00353     forAll(current, i)
00354     {
00355         label pointI = current[i];
00356 
00357         sortedCoord[i] = points[pointI][compI];
00358     }
00359     sortedCoord.sort();
00360 
00361     label globalCurrentSize = returnReduce(current.size(), sumOp<label>());
00362 
00363     scalar minCoord = returnReduce
00364     (
00365         (
00366             sortedCoord.size()
00367           ? sortedCoord[0]
00368           : GREAT
00369         ),
00370         minOp<scalar>()
00371     );
00372 
00373     scalar maxCoord = returnReduce
00374     (
00375         (
00376             sortedCoord.size()
00377           ? sortedCoord[sortedCoord.size()-1]
00378           : -GREAT
00379         ),
00380         maxOp<scalar>()
00381     );
00382 
00383     if (debug)
00384     {
00385         Pout<< "sortComponent : minCoord:" << minCoord
00386             << " maxCoord:" << maxCoord << endl;
00387     }
00388 
00389     // starting index (in sortedCoord) of bin (= local)
00390     label leftIndex = 0;
00391     // starting value of bin (= global since coordinate)
00392     scalar leftCoord = minCoord;
00393 
00394     // Sort bins of size n
00395     for (label bin = 0; bin < n_[compI]; bin++)
00396     {
00397         // Now we need to determine the size of the bin (dx). This is
00398         // determined by the 'pivot' values - everything to the left of this
00399         // value goes in the current bin, everything to the right into the next
00400         // bins.
00401 
00402         // Local number of elements
00403         label localSize = -1;     // offset from leftOffset
00404 
00405         // Value at right of bin (leftIndex+localSize-1)
00406         scalar rightCoord = -GREAT;
00407 
00408         if (bin == n_[compI]-1)
00409         {
00410             // Last bin. Copy all.
00411             localSize = current.size()-leftIndex;
00412             rightCoord = maxCoord;                  // note: not used anymore
00413         }
00414         else if (Pstream::nProcs() == 1)
00415         {
00416             // No need for binary searching of bin size
00417             localSize = label(current.size()/n_[compI]);
00418             rightCoord = sortedCoord[leftIndex+localSize];
00419         }
00420         else
00421         {
00422             // For the current bin (starting at leftCoord) we want a rightCoord
00423             // such that the sum of all sizes are globalCurrentSize/n_[compI].
00424             // We have to iterate to obtain this.
00425 
00426             label rightIndex = current.size();
00427             rightCoord = maxCoord;
00428 
00429             // Calculate rightIndex/rightCoord to have wanted size
00430             findBinary
00431             (
00432                 sizeTol,
00433                 sortedCoord,
00434                 leftIndex,
00435                 leftCoord,
00436                 maxCoord,
00437                 globalCurrentSize/n_[compI],  // wanted size
00438 
00439                 rightIndex,
00440                 rightCoord
00441             );
00442             localSize = rightIndex - leftIndex;
00443         }
00444 
00445         if (debug)
00446         {
00447             Pout<< "For component " << compI << ", bin " << bin
00448                 << " copying" << endl
00449                 << "from " << leftCoord << " at local index "
00450                 << leftIndex << endl
00451                 << "to " << rightCoord << " localSize:"
00452                 << localSize << endl
00453                 << endl;
00454         }
00455 
00456 
00457         // Copy localSize elements starting from leftIndex.
00458         labelList slice(localSize);
00459 
00460         forAll(slice, i)
00461         {
00462             label pointI = current[sortedCoord.indices()[leftIndex+i]];
00463 
00464             // Mark point into correct bin
00465             finalDecomp[pointI] += bin*mult;
00466 
00467             // And collect for next sorting action
00468             slice[i] = pointI;
00469         }
00470 
00471         // Sort slice in next component
00472         if (componentIndex < 2)
00473         {
00474             string oldPrefix;
00475             if (debug)
00476             {
00477                 oldPrefix = Pout.prefix();
00478                 Pout.prefix() = "  " + oldPrefix;
00479             }
00480 
00481             sortComponent
00482             (
00483                 sizeTol,
00484                 points,
00485                 slice,
00486                 componentIndex+1,
00487                 mult*n_[compI],     // Multiplier to apply to decomposition.
00488                 finalDecomp
00489             );
00490 
00491             if (debug)
00492             {
00493                 Pout.prefix() = oldPrefix;
00494             }
00495         }
00496 
00497         // Step to next bin.
00498         leftIndex += localSize;
00499         leftCoord = rightCoord;
00500     }
00501 }
00502 
00503 
00504 // Sort points into bins according to one component. Recurses to next component.
00505 void Foam::hierarchGeomDecomp::sortComponent
00506 (
00507     const label sizeTol,
00508     const scalarField& weights,
00509     const pointField& points,
00510     const labelList& current,       // slice of points to decompose
00511     const direction componentIndex, // index in decompOrder_
00512     const label mult,               // multiplication factor for finalDecomp
00513     labelList& finalDecomp
00514 )
00515 {
00516     // Current component
00517     label compI = decompOrder_[componentIndex];
00518 
00519     if (debug)
00520     {
00521         Pout<< "sortComponent : Sorting slice of size " << current.size()
00522             << " in component " << compI << endl;
00523     }
00524 
00525     // Storage for sorted component compI
00526     SortableList<scalar> sortedCoord(current.size());
00527 
00528     forAll(current, i)
00529     {
00530         label pointI = current[i];
00531 
00532         sortedCoord[i] = points[pointI][compI];
00533     }
00534     sortedCoord.sort();
00535 
00536     label globalCurrentSize = returnReduce(current.size(), sumOp<label>());
00537 
00538     // Now evaluate local cumulative weights, based on the sorting.
00539     // Make one bigger than the nodes.
00540     scalarField sortedWeightedSizes(current.size()+1, 0);
00541     calculateSortedWeightedSizes
00542     (
00543         current,
00544         sortedCoord.indices(),
00545         weights,
00546         globalCurrentSize,
00547         sortedWeightedSizes
00548     );
00549 
00550     scalar minCoord = returnReduce
00551     (
00552         (
00553             sortedCoord.size()
00554           ? sortedCoord[0]
00555           : GREAT
00556         ),
00557         minOp<scalar>()
00558     );
00559 
00560     scalar maxCoord = returnReduce
00561     (
00562         (
00563             sortedCoord.size()
00564           ? sortedCoord[sortedCoord.size()-1]
00565           : -GREAT
00566         ),
00567         maxOp<scalar>()
00568     );
00569 
00570     if (debug)
00571     {
00572         Pout<< "sortComponent : minCoord:" << minCoord
00573             << " maxCoord:" << maxCoord << endl;
00574     }
00575 
00576     // starting index (in sortedCoord) of bin (= local)
00577     label leftIndex = 0;
00578     // starting value of bin (= global since coordinate)
00579     scalar leftCoord = minCoord;
00580 
00581     // Sort bins of size n
00582     for (label bin = 0; bin < n_[compI]; bin++)
00583     {
00584         // Now we need to determine the size of the bin (dx). This is
00585         // determined by the 'pivot' values - everything to the left of this
00586         // value goes in the current bin, everything to the right into the next
00587         // bins.
00588 
00589         // Local number of elements
00590         label localSize = -1;     // offset from leftOffset
00591 
00592         // Value at right of bin (leftIndex+localSize-1)
00593         scalar rightCoord = -GREAT;
00594 
00595         if (bin == n_[compI]-1)
00596         {
00597             // Last bin. Copy all.
00598             localSize = current.size()-leftIndex;
00599             rightCoord = maxCoord;                  // note: not used anymore
00600         }
00601         else
00602         {
00603             // For the current bin (starting at leftCoord) we want a rightCoord
00604             // such that the sum of all weighted sizes are
00605             // globalCurrentLength/n_[compI].
00606             // We have to iterate to obtain this.
00607 
00608             label rightIndex = current.size();
00609             rightCoord = maxCoord;
00610 
00611             // Calculate rightIndex/rightCoord to have wanted size
00612             findBinary
00613             (
00614                 sizeTol,
00615                 sortedWeightedSizes,
00616                 sortedCoord,
00617                 leftIndex,
00618                 leftCoord,
00619                 maxCoord,
00620                 globalCurrentSize/n_[compI],  // wanted size
00621 
00622                 rightIndex,
00623                 rightCoord
00624             );
00625             localSize = rightIndex - leftIndex;
00626         }
00627 
00628         if (debug)
00629         {
00630             Pout<< "For component " << compI << ", bin " << bin
00631                 << " copying" << endl
00632                 << "from " << leftCoord << " at local index "
00633                 << leftIndex << endl
00634                 << "to " << rightCoord << " localSize:"
00635                 << localSize << endl
00636                 << endl;
00637         }
00638 
00639 
00640         // Copy localSize elements starting from leftIndex.
00641         labelList slice(localSize);
00642 
00643         forAll(slice, i)
00644         {
00645             label pointI = current[sortedCoord.indices()[leftIndex+i]];
00646 
00647             // Mark point into correct bin
00648             finalDecomp[pointI] += bin*mult;
00649 
00650             // And collect for next sorting action
00651             slice[i] = pointI;
00652         }
00653 
00654         // Sort slice in next component
00655         if (componentIndex < 2)
00656         {
00657             string oldPrefix;
00658             if (debug)
00659             {
00660                 oldPrefix = Pout.prefix();
00661                 Pout.prefix() = "  " + oldPrefix;
00662             }
00663 
00664             sortComponent
00665             (
00666                 sizeTol,
00667                 weights,
00668                 points,
00669                 slice,
00670                 componentIndex+1,
00671                 mult*n_[compI],     // Multiplier to apply to decomposition.
00672                 finalDecomp
00673             );
00674 
00675             if (debug)
00676             {
00677                 Pout.prefix() = oldPrefix;
00678             }
00679         }
00680 
00681         // Step to next bin.
00682         leftIndex += localSize;
00683         leftCoord = rightCoord;
00684     }
00685 }
00686 
00687 
00688 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00689 
00690 Foam::hierarchGeomDecomp::hierarchGeomDecomp
00691 (
00692     const dictionary& decompositionDict
00693 )
00694 :
00695     geomDecomp(decompositionDict, typeName),
00696     decompOrder_()
00697 {
00698     setDecompOrder();
00699 }
00700 
00701 
00702 Foam::hierarchGeomDecomp::hierarchGeomDecomp
00703 (
00704     const dictionary& decompositionDict,
00705     const polyMesh&
00706 )
00707 :
00708     geomDecomp(decompositionDict, hierarchGeomDecomp::typeName),
00709     decompOrder_()
00710 {
00711     setDecompOrder();
00712 }
00713 
00714 
00715 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00716 
00717 Foam::labelList Foam::hierarchGeomDecomp::decompose
00718 (
00719     const pointField& points
00720 )
00721 {
00722     // construct a list for the final result
00723     labelList finalDecomp(points.size(), 0);
00724 
00725     // Start off with every point sorted onto itself.
00726     labelList slice(points.size());
00727     forAll(slice, i)
00728     {
00729         slice[i] = i;
00730     }
00731 
00732     pointField rotatedPoints = rotDelta_ & points;
00733 
00734 
00735     // Calculate tolerance of cell distribution. For large cases finding
00736     // distibution to the cell exact would cause too many iterations so allow
00737     // some slack.
00738     label allSize = points.size();
00739     reduce(allSize, sumOp<label>());
00740 
00741     const label sizeTol = max(1, label(1E-3*allSize/nProcessors_));
00742 
00743     // Sort recursive
00744     sortComponent
00745     (
00746         sizeTol,
00747         rotatedPoints,
00748         slice,
00749         0,              // Sort first component in decompOrder.
00750         1,              // Offset for different x bins.
00751         finalDecomp
00752     );
00753 
00754     return finalDecomp;
00755 }
00756 
00757 
00758 Foam::labelList Foam::hierarchGeomDecomp::decompose
00759 (
00760     const pointField& points,
00761     const scalarField& weights
00762 )
00763 {
00764     // construct a list for the final result
00765     labelList finalDecomp(points.size(), 0);
00766 
00767     // Start off with every point sorted onto itself.
00768     labelList slice(points.size());
00769     forAll(slice, i)
00770     {
00771         slice[i] = i;
00772     }
00773 
00774     pointField rotatedPoints = rotDelta_ & points;
00775 
00776     // Calculate tolerance of cell distribution. For large cases finding
00777     // distibution to the cell exact would cause too many iterations so allow
00778     // some slack.
00779     label allSize = points.size();
00780     reduce(allSize, sumOp<label>());
00781 
00782     const label sizeTol = max(1, label(1E-3*allSize/nProcessors_));
00783 
00784     // Sort recursive
00785     sortComponent
00786     (
00787         sizeTol,
00788         weights,
00789         rotatedPoints,
00790         slice,
00791         0,              // Sort first component in decompOrder.
00792         1,              // Offset for different x bins.
00793         finalDecomp
00794     );
00795 
00796     return finalDecomp;
00797 }
00798 
00799 
00800 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines