domainDecompositionDistribute.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "domainDecomposition.H"
27 #include "decompositionMethod.H"
28 #include "cpuTime.H"
29 #include "cellSet.H"
30 #include "regionSplit.H"
31 #include "Tuple2.H"
32 #include "faceSet.H"
33 #include "decompositionModel.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 void Foam::domainDecomposition::distributeCells()
38 {
39  Info<< "\nCalculating distribution of cells" << endl;
40 
41  cpuTime decompositionTime;
42 
43  const decompositionModel& method = decompositionModel::New(*this);
44 
45  scalarField cellWeights;
46  if (method.found("weightField"))
47  {
48  word weightName = method.lookup("weightField");
49 
51  (
52  IOobject
53  (
54  weightName,
55  time().timeName(),
56  *this,
59  ),
60  *this
61  );
62  cellWeights = weights.primitiveField();
63  }
64 
65  cellToProc_ = method.decomposer().decompose(*this, cellWeights);
66 
67  Info<< "\nFinished decomposition in "
68  << decompositionTime.elapsedCpuTime()
69  << " s" << endl;
70 }
71 
72 
73 // ************************************************************************* //
static const decompositionModel & New(const polyMesh &mesh, const fileName &decompDictFile="")
Read (optionallly from absolute path) & register on mesh.
IOobject(const word &name, const fileName &instance, const objectRegistry &registry, readOption r=NO_READ, writeOption w=NO_WRITE, bool registerObject=true)
Construct from name, instance, registry, io options.
Definition: IOobject.C:114
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
const Internal::FieldType & primitiveField() const
Return a const-reference to the internal field.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
word timeName
Definition: getTimeIndex.H:3
const surfaceScalarField & weights() const
Return reference to linear difference weighting factors.
messageStream Info
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:243