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

simpleGeomDecomp.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 "simpleGeomDecomp.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 #include <OpenFOAM/SortableList.H>
00029 
00030 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00031 
00032 namespace Foam
00033 {
00034     defineTypeNameAndDebug(simpleGeomDecomp, 0);
00035 
00036     addToRunTimeSelectionTable
00037     (
00038         decompositionMethod,
00039         simpleGeomDecomp,
00040         dictionary
00041     );
00042 
00043     addToRunTimeSelectionTable
00044     (
00045         decompositionMethod,
00046         simpleGeomDecomp,
00047         dictionaryMesh
00048     );
00049 }
00050 
00051 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00052 
00053 // assignToProcessorGroup : given nCells cells and nProcGroup processor
00054 // groups to share them, how do we share them out? Answer : each group
00055 // gets nCells/nProcGroup cells, and the first few get one
00056 // extra to make up the numbers. This should produce almost
00057 // perfect load balancing
00058 
00059 void Foam::simpleGeomDecomp::assignToProcessorGroup
00060 (
00061     labelList& processorGroup,
00062     const label nProcGroup
00063 )
00064 {
00065     label jump = processorGroup.size()/nProcGroup;
00066     label jumpb = jump + 1;
00067     label fstProcessorGroup = processorGroup.size() - jump*nProcGroup;
00068 
00069     label ind = 0;
00070     label j = 0;
00071 
00072     // assign cells to the first few processor groups (those with
00073     // one extra cell each
00074     for (j=0; j<fstProcessorGroup; j++)
00075     {
00076         for (register label k=0; k<jumpb; k++)
00077         {
00078             processorGroup[ind++] = j;
00079         }
00080     }
00081 
00082     // and now to the `normal' processor groups
00083     for (; j<nProcGroup; j++)
00084     {
00085         for (register label k=0; k<jump; k++)
00086         {
00087             processorGroup[ind++] = j;
00088         }
00089     }
00090 }
00091 
00092 
00093 void Foam::simpleGeomDecomp::assignToProcessorGroup
00094 (
00095     labelList& processorGroup,
00096     const label nProcGroup,
00097     const labelList& indices,
00098     const scalarField& weights,
00099     const scalar summedWeights
00100 )
00101 {
00102     // This routine gets the sorted points.
00103     // Easiest to explain with an example.
00104     // E.g. 400 points, sum of weights : 513.
00105     // Now with number of divisions in this direction (nProcGroup) : 4
00106     // gives the split at 513/4 = 128
00107     // So summed weight from 0..128 goes into bin 0,
00108     //     ,,              128..256 goes into bin 1
00109     //   etc.
00110     // Finally any remaining ones go into the last bin (3).
00111 
00112     const scalar jump = summedWeights/nProcGroup;
00113     const label nProcGroupM1 = nProcGroup - 1;
00114     scalar sumWeights = 0;
00115     label ind = 0;
00116     label j = 0;
00117 
00118     // assign cells to all except last group.
00119     for (j=0; j<nProcGroupM1; j++)
00120     {
00121         const scalar limit = jump*scalar(j + 1);
00122         while (sumWeights < limit)
00123         {
00124             sumWeights += weights[indices[ind]];
00125             processorGroup[ind++] = j;
00126         }
00127     }
00128     // Ensure last included.
00129     while (ind < processorGroup.size())
00130     {
00131        processorGroup[ind++] = nProcGroupM1;
00132     }
00133 }
00134 
00135 
00136 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00137 
00138 Foam::simpleGeomDecomp::simpleGeomDecomp(const dictionary& decompositionDict)
00139 :
00140     geomDecomp(decompositionDict, typeName)
00141 {}
00142 
00143 
00144 Foam::simpleGeomDecomp::simpleGeomDecomp
00145 (
00146     const dictionary& decompositionDict,
00147     const polyMesh&
00148 )
00149 :
00150     geomDecomp(decompositionDict, typeName)
00151 {}
00152 
00153 
00154 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00155 
00156 Foam::labelList Foam::simpleGeomDecomp::decompose(const pointField& points)
00157 {
00158     // construct a list for the final result
00159     labelList finalDecomp(points.size());
00160 
00161     labelList processorGroups(points.size());
00162 
00163     labelList pointIndices(points.size());
00164     forAll(pointIndices, i)
00165     {
00166         pointIndices[i] = i;
00167     }
00168 
00169     pointField rotatedPoints = rotDelta_ & points;
00170 
00171     // and one to take the processor group id's. For each direction.
00172     // we assign the processors to groups of processors labelled
00173     // 0..nX to give a banded structure on the mesh. Then we
00174     // construct the actual processor number by treating this as
00175     // the units part of the processor number.
00176     sort
00177     (
00178         pointIndices,
00179         UList<scalar>::less(rotatedPoints.component(vector::X))
00180     );
00181 
00182     assignToProcessorGroup(processorGroups, n_.x());
00183 
00184     forAll(points, i)
00185     {
00186         finalDecomp[pointIndices[i]] = processorGroups[i];
00187     }
00188 
00189 
00190     // now do the same thing in the Y direction. These processor group
00191     // numbers add multiples of nX to the proc. number (columns)
00192     sort
00193     (
00194         pointIndices,
00195         UList<scalar>::less(rotatedPoints.component(vector::Y))
00196     );
00197 
00198     assignToProcessorGroup(processorGroups, n_.y());
00199 
00200     forAll(points, i)
00201     {
00202         finalDecomp[pointIndices[i]] += n_.x()*processorGroups[i];
00203     }
00204 
00205 
00206     // finally in the Z direction. Now we add multiples of nX*nY to give
00207     // layers
00208     sort
00209     (
00210         pointIndices,
00211         UList<scalar>::less(rotatedPoints.component(vector::Z))
00212     );
00213 
00214     assignToProcessorGroup(processorGroups, n_.z());
00215 
00216     forAll(points, i)
00217     {
00218         finalDecomp[pointIndices[i]] += n_.x()*n_.y()*processorGroups[i];
00219     }
00220 
00221     return finalDecomp;
00222 }
00223 
00224 
00225 Foam::labelList Foam::simpleGeomDecomp::decompose
00226 (
00227     const pointField& points,
00228     const scalarField& weights
00229 )
00230 {
00231     // construct a list for the final result
00232     labelList finalDecomp(points.size());
00233 
00234     labelList processorGroups(points.size());
00235 
00236     labelList pointIndices(points.size());
00237     forAll(pointIndices, i)
00238     {
00239         pointIndices[i] = i;
00240     }
00241 
00242     pointField rotatedPoints = rotDelta_ & points;
00243 
00244     // and one to take the processor group id's. For each direction.
00245     // we assign the processors to groups of processors labelled
00246     // 0..nX to give a banded structure on the mesh. Then we
00247     // construct the actual processor number by treating this as
00248     // the units part of the processor number.
00249     sort
00250     (
00251         pointIndices,
00252         UList<scalar>::less(rotatedPoints.component(vector::X))
00253     );
00254 
00255     const scalar summedWeights = sum(weights);
00256     assignToProcessorGroup
00257     (
00258         processorGroups,
00259         n_.x(),
00260         pointIndices,
00261         weights,
00262         summedWeights
00263     );
00264 
00265     forAll(points, i)
00266     {
00267         finalDecomp[pointIndices[i]] = processorGroups[i];
00268     }
00269 
00270 
00271     // now do the same thing in the Y direction. These processor group
00272     // numbers add multiples of nX to the proc. number (columns)
00273     sort
00274     (
00275         pointIndices,
00276         UList<scalar>::less(rotatedPoints.component(vector::Y))
00277     );
00278 
00279     assignToProcessorGroup
00280     (
00281         processorGroups,
00282         n_.y(),
00283         pointIndices,
00284         weights,
00285         summedWeights
00286     );
00287 
00288     forAll(points, i)
00289     {
00290         finalDecomp[pointIndices[i]] += n_.x()*processorGroups[i];
00291     }
00292 
00293 
00294     // finally in the Z direction. Now we add multiples of nX*nY to give
00295     // layers
00296     sort
00297     (
00298         pointIndices,
00299         UList<scalar>::less(rotatedPoints.component(vector::Z))
00300     );
00301 
00302     assignToProcessorGroup
00303     (
00304         processorGroups,
00305         n_.z(),
00306         pointIndices,
00307         weights,
00308         summedWeights
00309     );
00310 
00311     forAll(points, i)
00312     {
00313         finalDecomp[pointIndices[i]] += n_.x()*n_.y()*processorGroups[i];
00314     }
00315 
00316     return finalDecomp;
00317 }
00318 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines