ENH: support 'transform' specification for geometric decomposition

- can be used for block-like meshes that are not aligned with the global
  coordinate directions. Alternatively, for general testing purposes.

  Example,

    method  simple;
    coeffs
    {
        n       ( 2 2 2 );
        transform
        {
            origin  (-0.15 0.15 0);
            e1      (1 1 0);
            e3      (0 0 1);
        }
    }
This commit is contained in:
Mark Olesen
2021-04-29 13:54:37 +02:00
parent 36266a7e7d
commit 492d5cb645
15 changed files with 477 additions and 140 deletions

View File

@ -71,6 +71,7 @@ $(crot)/axisAngleRotation.C
$(crot)/coordinateRotation.C $(crot)/coordinateRotation.C
$(crot)/cylindricalRotation.C $(crot)/cylindricalRotation.C
$(crot)/identityRotation.C $(crot)/identityRotation.C
$(crot)/specifiedRotation.C
$(crot)/EulerCoordinateRotation.C $(crot)/EulerCoordinateRotation.C
$(crot)/STARCDCoordinateRotation.C $(crot)/STARCDCoordinateRotation.C

View File

@ -0,0 +1,104 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2021 OpenCFD Ltd.
-------------------------------------------------------------------------------
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 <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "specifiedRotation.H"
#include "dictionary.H"
#include "addToRunTimeSelectionTable.H"
// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
namespace Foam
{
namespace coordinateRotations
{
defineTypeName(specified);
//FUTURE addToRunTimeSelectionTable(coordinateRotation, specified, dictionary);
}
}
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::coordinateRotations::specified::specified()
:
coordinateRotation(),
Rmatrix_(sphericalTensor::I)
{}
Foam::coordinateRotations::specified::specified(const tensor& rot)
:
coordinateRotation(),
Rmatrix_(rot)
{}
Foam::coordinateRotations::specified::specified(const dictionary& dict)
:
specified()
{
// Not yet implemented, pending definition
}
// * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
void Foam::coordinateRotations::specified::clear()
{
Rmatrix_ = sphericalTensor::I;
}
Foam::tensor Foam::coordinateRotations::specified::R() const
{
return Rmatrix_;
}
void Foam::coordinateRotations::specified::write(Ostream& os) const
{
os << "specified rotation";
}
void Foam::coordinateRotations::specified::writeEntry
(
const word& keyword,
Ostream& os
) const
{
os.beginBlock(keyword);
os.writeEntry("type", type());
os.endBlock();
}
// ************************************************************************* //

View File

@ -0,0 +1,127 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2021 OpenCFD Ltd.
-------------------------------------------------------------------------------
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 <http://www.gnu.org/licenses/>.
Class
Foam::coordinateRotations::specified
Description
An user-specified coordinateRotation, primarily to be used internally
within coding when the rotation matrix is already known but
needs to be wrapped as a coordinateRotation for use in a coordinate
system.
Note
Currently no runtime selection mechanism, since specifying a
rotation matrix via a dictionary can be fragile due to rounding
factors and uncertainty if a forward or reverse rotation matrix
is intended.
SourceFiles
specifiedRotation.C
\*---------------------------------------------------------------------------*/
#ifndef coordinateRotations_specified_H
#define coordinateRotations_specified_H
#include "coordinateRotation.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
namespace coordinateRotations
{
/*---------------------------------------------------------------------------*\
Class coordinateRotations::specified Declaration
\*---------------------------------------------------------------------------*/
class specified
:
public coordinateRotation
{
// Private Data
//- The rotation matrix
tensor Rmatrix_;
public:
//- Runtime type information
TypeNameNoDebug("specified");
// Constructors
//- Default construct - an identity matrix
specified();
//- Construct from transformation matrix
explicit specified(const tensor& rot);
//- Construct from dictionary
explicit specified(const dictionary& unused);
//- Return clone
autoPtr<coordinateRotation> clone() const
{
return
autoPtr<coordinateRotation>::NewFrom
<coordinateRotations::specified>(*this);
}
//- Destructor
virtual ~specified() = default;
// Member Functions
//- Reset specification
virtual void clear();
//- Return the rotation tensor
virtual tensor R() const;
//- Write information
virtual void write(Ostream& os) const;
//- Write dictionary entry
virtual void writeEntry(const word& keyword, Ostream& os) const;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace coordinateRotations
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -6,7 +6,7 @@
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
Copyright (C) 2011-2014 OpenFOAM Foundation Copyright (C) 2011-2014 OpenFOAM Foundation
Copyright (C) 2018-2019 OpenCFD Ltd. Copyright (C) 2018-2021 OpenCFD Ltd.
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
This file is part of OpenFOAM. This file is part of OpenFOAM.
@ -70,6 +70,12 @@ Foam::coordSystem::cartesian::cartesian(autoPtr<coordinateSystem>&& csys)
{} {}
Foam::coordSystem::cartesian::cartesian(const coordinateRotation& crot)
:
coordinateSystem(crot)
{}
Foam::coordSystem::cartesian::cartesian Foam::coordSystem::cartesian::cartesian
( (
const point& origin, const point& origin,

View File

@ -6,7 +6,7 @@
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
Copyright (C) 2011-2013 OpenFOAM Foundation Copyright (C) 2011-2013 OpenFOAM Foundation
Copyright (C) 2018-2019 OpenCFD Ltd. Copyright (C) 2018-2021 OpenCFD Ltd.
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
This file is part of OpenFOAM. This file is part of OpenFOAM.
@ -75,7 +75,7 @@ public:
// Constructors // Constructors
//- Construct null. This is an identity coordinateSystem. //- Default construct. This is an identity coordinate system
cartesian(); cartesian();
//- Copy construct //- Copy construct
@ -93,6 +93,12 @@ public:
//- Move construct from autoPtr of another coordinateSystem type //- Move construct from autoPtr of another coordinateSystem type
explicit cartesian(autoPtr<coordinateSystem>&& csys); explicit cartesian(autoPtr<coordinateSystem>&& csys);
//- Copy construct from rotation with origin=0
explicit cartesian(const coordinateRotation& crot);
//- Move construct from rotation with origin=0
explicit cartesian(coordinateRotation&& crot);
//- Construct from origin and rotation //- Construct from origin and rotation
cartesian(const point& origin, const coordinateRotation& crot); cartesian(const point& origin, const coordinateRotation& crot);

View File

@ -260,7 +260,7 @@ public:
// Constructors // Constructors
//- Construct null. This is an identity coordinateSystem. //- Default construct. This is an identity coordinate system
coordinateSystem(); coordinateSystem();
//- Copy construct from rotation with origin=0 //- Copy construct from rotation with origin=0

View File

@ -6,7 +6,7 @@
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
Copyright (C) 2011-2016 OpenFOAM Foundation Copyright (C) 2011-2016 OpenFOAM Foundation
Copyright (C) 2018-2020 OpenCFD Ltd. Copyright (C) 2018-2021 OpenCFD Ltd.
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
This file is part of OpenFOAM. This file is part of OpenFOAM.
@ -108,6 +108,12 @@ Foam::coordSystem::cylindrical::cylindrical(autoPtr<coordinateSystem>&& csys)
{} {}
Foam::coordSystem::cylindrical::cylindrical(const coordinateRotation& crot)
:
coordinateSystem(crot)
{}
Foam::coordSystem::cylindrical::cylindrical Foam::coordSystem::cylindrical::cylindrical
( (
const point& origin, const point& origin,

View File

@ -6,7 +6,7 @@
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
Copyright (C) 2011-2014 OpenFOAM Foundation Copyright (C) 2011-2014 OpenFOAM Foundation
Copyright (C) 2018 OpenCFD Ltd. Copyright (C) 2018-2021 OpenCFD Ltd.
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
This file is part of OpenFOAM. This file is part of OpenFOAM.
@ -114,7 +114,7 @@ public:
// Constructors // Constructors
//- Construct null (identity coordinateSystem) //- Default construct. This is an identity coordinate system
cylindrical(); cylindrical();
//- Copy construct //- Copy construct
@ -132,6 +132,9 @@ public:
//- Move construct from autoPtr of another coordinateSystem type //- Move construct from autoPtr of another coordinateSystem type
explicit cylindrical(autoPtr<coordinateSystem>&& csys); explicit cylindrical(autoPtr<coordinateSystem>&& csys);
//- Copy construct from rotation with origin=0
explicit cylindrical(const coordinateRotation& crot);
//- Construct from origin and rotation //- Construct from origin and rotation
cylindrical(const point& origin, const coordinateRotation& crot); cylindrical(const point& origin, const coordinateRotation& crot);

View File

@ -6,7 +6,7 @@
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
Copyright (C) 2011-2017 OpenFOAM Foundation Copyright (C) 2011-2017 OpenFOAM Foundation
Copyright (C) 2018 OpenCFD Ltd. Copyright (C) 2018-2021 OpenCFD Ltd.
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
This file is part of OpenFOAM. This file is part of OpenFOAM.
@ -27,36 +27,103 @@ License
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
#include "geomDecomp.H" #include "geomDecomp.H"
#include "specifiedRotation.H"
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * // // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
void Foam::geomDecomp::setOrder()
{
const word order(coeffsDict_.getOrDefault<word>("order", ""));
if (order.empty())
{
return;
}
else if (order.size() != 3)
{
FatalIOErrorInFunction(decompDict_)
<< "Number of characters in order (" << order << ") != 3"
<< exit(FatalIOError);
}
for (int i = 0; i < 3; ++i)
{
// Change [x-z] -> [0-2]
switch (order[i])
{
case 'x': order_[i] = 0; break;
case 'y': order_[i] = 1; break;
case 'z': order_[i] = 2; break;
default:
FatalIOErrorInFunction(decompDict_)
<< "Illegal decomposition order " << order << nl
<< "It should only contain x, y or z"
<< exit(FatalError);
break;
}
}
}
void Foam::geomDecomp::readCoeffs() void Foam::geomDecomp::readCoeffs()
{ {
coeffsDict_.readIfPresent("delta", delta_); coeffsDict_.readIfPresent("delta", delta_);
coeffsDict_.readEntry("n", n_); coeffsDict_.readEntry("n", n_);
// Verify that the input makes sense
if (nDomains_ != n_.x()*n_.y()*n_.z()) if (nDomains_ != n_.x()*n_.y()*n_.z())
{ {
// Verify that the input makes sense
FatalErrorInFunction FatalErrorInFunction
<< "Wrong number of domain divisions in geomDecomp:" << nl << "Wrong number of domain divisions in geomDecomp:" << nl
<< "Number of domains : " << nDomains_ << nl << "Number of domains : " << nDomains_ << nl
<< "Wanted decomposition : " << n_ << "Wanted decomposition : " << n_
<< exit(FatalError); << exit(FatalError);
} }
setOrder();
const scalar d = 1 - 0.5*delta_*delta_; const dictionary* transformDict =
const scalar d2 = sqr(d); coeffsDict_.findDict("transform", keyType::LITERAL);
const scalar a = delta_; if (transformDict)
const scalar a2 = sqr(a); {
csys_ = coordinateSystem(*transformDict);
}
else if (equal(delta_, 0))
{
csys_.clear(); // Reset to identity
}
else
{
const scalar d = 1 - 0.5*delta_*delta_;
const scalar d2 = sqr(d);
rotDelta_ = tensor const scalar a = delta_;
( const scalar a2 = sqr(a);
d2, -a*d, a,
a*d - a2*d, a*a2 + d2, -2*a*d, // Direction (forward/reverse) doesn't matter much
a*d2 + a2, a*d - a2*d, d2 - a2 tensor rot
); (
d2, -a*d, a,
a*d - a2*d, a*a2 + d2, -2*a*d,
a*d2 + a2, a*d - a2*d, d2 - a2
);
// origin=0
csys_ = coordinateSystem(coordinateRotations::specified(rot));
}
}
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
Foam::tmp<Foam::pointField> Foam::geomDecomp::adjustPoints
(
const pointField& points
) const
{
return csys_.localPosition(points);
} }
@ -90,10 +157,11 @@ Foam::geomDecomp::geomDecomp
) )
: :
decompositionMethod(decompDict), decompositionMethod(decompDict),
coeffsDict_(findCoeffsDict(derivedType + "Coeffs", select)),
n_(1,1,1),
delta_(0.001), delta_(0.001),
rotDelta_(I) csys_(),
n_(1,1,1),
order_(0,1,2),
coeffsDict_(findCoeffsDict(derivedType + "Coeffs", select))
{ {
readCoeffs(); readCoeffs();
} }
@ -108,10 +176,11 @@ Foam::geomDecomp::geomDecomp
) )
: :
decompositionMethod(decompDict, regionName), decompositionMethod(decompDict, regionName),
coeffsDict_(findCoeffsDict(derivedType + "Coeffs", select)),
n_(1,1,1),
delta_(0.001), delta_(0.001),
rotDelta_(I) csys_(),
n_(1,1,1),
order_(0,1,2),
coeffsDict_(findCoeffsDict(derivedType + "Coeffs", select))
{ {
readCoeffs(); readCoeffs();
} }

View File

@ -6,7 +6,7 @@
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
Copyright (C) 2011 OpenFOAM Foundation Copyright (C) 2011 OpenFOAM Foundation
Copyright (C) 2018 OpenCFD Ltd. Copyright (C) 2018-2021 OpenCFD Ltd.
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
This file is part of OpenFOAM. This file is part of OpenFOAM.
@ -32,9 +32,11 @@ Description
Base coefficients: Base coefficients:
\table \table
Property | Description | Required | Default Property | Description | Required | Default
n | (nx ny nz) | yes n | (nx ny nz) | yes |
delta | delta for rotation matrix | no | 0.001 order | order of operation | no | xyz
delta | delta (jitter) for rotation matrix | no | 0.001
rotation | coordinate rotation | no |
\endtable \endtable
SourceFiles SourceFiles
@ -46,7 +48,9 @@ SourceFiles
#define geomDecomp_H #define geomDecomp_H
#include "decompositionMethod.H" #include "decompositionMethod.H"
#include "cartesianCS.H"
#include "Vector.H" #include "Vector.H"
#include "tmp.H"
namespace Foam namespace Foam
{ {
@ -59,29 +63,45 @@ class geomDecomp
: :
public decompositionMethod public decompositionMethod
{ {
// Private Data
//- Small delta (default: 0.001) to avoid staggering when
//- mesh itself is aligned with x/y/z
scalar delta_;
//- Local coordinates, or delta rotation tensor
coordSystem::cartesian csys_;
// Private Member Functions // Private Member Functions
//- Read input values and initialize the rotDelta_ //- Convert ordering string ("xyz") into list of components.
// Checks for bad entries, but no check for duplicate entries.
void setOrder();
//- Read input values and initialize rotation matrix
void readCoeffs(); void readCoeffs();
protected: protected:
// Protected data // Protected Data
//- The divisions
Vector<label> n_;
//- Decomposition order in terms of components (optional)
Vector<direction> order_;
//- Coefficients for all derived methods //- Coefficients for all derived methods
const dictionary& coeffsDict_; const dictionary& coeffsDict_;
Vector<label> n_;
//- Default = 0.001
scalar delta_;
tensor rotDelta_;
// Protected Member Functions // Protected Member Functions
//- Apply delta (jitter) or rotation to coordinates
tmp<pointField> adjustPoints(const pointField&) const;
//- Check that mesh directions are compatible with decomposition //- Check that mesh directions are compatible with decomposition
void checkDecompositionDirections(const Vector<label>&) const; void checkDecompositionDirections(const Vector<label>&) const;
@ -120,7 +140,6 @@ public:
//- Decompose with uniform weights on the points //- Decompose with uniform weights on the points
virtual labelList decompose(const pointField& points) const = 0; virtual labelList decompose(const pointField& points) const = 0;
}; };

View File

@ -47,42 +47,6 @@ namespace Foam
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * // // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
void Foam::hierarchGeomDecomp::setOrder()
{
const word order(coeffsDict_.getOrDefault<word>("order", ""));
if (order.empty())
{
return;
}
else if (order.size() != 3)
{
FatalIOErrorInFunction(decompDict_)
<< "Number of characters in order (" << order << ") != 3"
<< exit(FatalIOError);
}
for (int i = 0; i < 3; ++i)
{
// Change [x-z] -> [0-2]
switch (order[i])
{
case 'x': order_[i] = 0; break;
case 'y': order_[i] = 1; break;
case 'z': order_[i] = 2; break;
default:
FatalIOErrorInFunction(decompDict_)
<< "Illegal decomposition order " << order << nl
<< "It should only contain x, y or z"
<< exit(FatalError);
break;
}
}
}
Foam::label Foam::hierarchGeomDecomp::findLower Foam::label Foam::hierarchGeomDecomp::findLower
( (
const UList<scalar>& list, const UList<scalar>& list,
@ -363,7 +327,7 @@ Foam::label Foam::hierarchGeomDecomp::sortComponent
( (
( (
sortedCoord.size() sortedCoord.size()
? sortedCoord[0] ? sortedCoord.first()
: GREAT : GREAT
), ),
minOp<scalar>() minOp<scalar>()
@ -567,7 +531,7 @@ Foam::label Foam::hierarchGeomDecomp::sortComponent
( (
( (
sortedCoord.size() sortedCoord.size()
? sortedCoord[0] ? sortedCoord.first()
: GREAT : GREAT
), ),
minOp<scalar>() minOp<scalar>()
@ -716,11 +680,8 @@ Foam::hierarchGeomDecomp::hierarchGeomDecomp
const word& regionName const word& regionName
) )
: :
geomDecomp(typeName, decompDict, regionName), geomDecomp(typeName, decompDict, regionName)
order_({0,1,2}) {}
{
setOrder();
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
@ -736,7 +697,7 @@ Foam::labelList Foam::hierarchGeomDecomp::decompose
// Start off with every point sorted onto itself. // Start off with every point sorted onto itself.
labelList slice(identity(points.size())); labelList slice(identity(points.size()));
pointField rotatedPoints(rotDelta_ & points); const pointField rotatedPoints(adjustPoints(points));
// Calculate tolerance of cell distribution. For large cases finding // Calculate tolerance of cell distribution. For large cases finding
// distribution to the cell exact would cause too many iterations so allow // distribution to the cell exact would cause too many iterations so allow
@ -778,7 +739,7 @@ Foam::labelList Foam::hierarchGeomDecomp::decompose
// Start off with every point sorted onto itself. // Start off with every point sorted onto itself.
labelList slice(identity(points.size())); labelList slice(identity(points.size()));
pointField rotatedPoints(rotDelta_ & points); const pointField rotatedPoints(adjustPoints(points));
// Calculate tolerance of cell distribution. For large cases finding // Calculate tolerance of cell distribution. For large cases finding
// distribution to the cell exact would cause too many iterations so allow // distribution to the cell exact would cause too many iterations so allow

View File

@ -52,8 +52,9 @@ Description
\table \table
Property | Description | Required | Default Property | Description | Required | Default
n | (nx ny nz) | yes | n | (nx ny nz) | yes |
delta | delta for rotation matrix | no | 0.001
order | order of operation | no | xyz order | order of operation | no | xyz
delta | delta (jitter) for rotation matrix | no | 0.001
transform | cartesian coordinate transformation | no |
\endtable \endtable
SourceFiles SourceFiles
@ -65,8 +66,6 @@ SourceFiles
#define hierarchGeomDecomp_H #define hierarchGeomDecomp_H
#include "geomDecomp.H" #include "geomDecomp.H"
#include "FixedList.H"
#include "direction.H"
namespace Foam namespace Foam
{ {
@ -79,18 +78,8 @@ class hierarchGeomDecomp
: :
public geomDecomp public geomDecomp
{ {
// Private Data
//- Decomposition order in terms of components.
FixedList<direction, 3> order_;
// Private Member Functions // Private Member Functions
//- Convert ordering string ("xyz") into list of components.
// Checks for bad entries, but no check for duplicate entries.
void setOrder();
//- Find index of value in list between //- Find index of value in list between
//- first (inclusive) and last (exclusive) //- first (inclusive) and last (exclusive)
static label findLower static label findLower
@ -172,15 +161,14 @@ class hierarchGeomDecomp
) const; ) const;
//- No copy construct
hierarchGeomDecomp(const hierarchGeomDecomp&) = delete;
//- No copy assignment
void operator=(const hierarchGeomDecomp&) = delete;
public: public:
//- No copy construct
hierarchGeomDecomp(const hierarchGeomDecomp&) = delete;
//- No copy assignment
void operator=(const hierarchGeomDecomp&) = delete;
//- Runtime type information //- Runtime type information
TypeName("hierarchical"); TypeName("hierarchical");
@ -255,7 +243,7 @@ public:
// - the connections are across coupled patches // - the connections are across coupled patches
virtual labelList decompose virtual labelList decompose
( (
const labelListList& globalCellCells, const labelListList& globalCellCells, // unused
const pointField& cc, const pointField& cc,
const scalarField& cWeights const scalarField& cWeights
) const ) const
@ -265,7 +253,7 @@ public:
virtual labelList decompose virtual labelList decompose
( (
const labelListList& globalCellCells, const labelListList& globalCellCells, // unused
const pointField& cc const pointField& cc
) const ) const
{ {

View File

@ -28,7 +28,6 @@ License
#include "simpleGeomDecomp.H" #include "simpleGeomDecomp.H"
#include "addToRunTimeSelectionTable.H" #include "addToRunTimeSelectionTable.H"
#include "SortableList.H"
#include "globalIndex.H" #include "globalIndex.H"
#include "SubField.H" #include "SubField.H"
@ -174,7 +173,7 @@ Foam::labelList Foam::simpleGeomDecomp::decomposeOneProc
labelList pointIndices(identity(points.size())); labelList pointIndices(identity(points.size()));
const pointField rotatedPoints(rotDelta_ & points); const pointField rotatedPoints(adjustPoints(points));
vectorLessOp sorter(rotatedPoints); vectorLessOp sorter(rotatedPoints);
@ -239,7 +238,7 @@ Foam::labelList Foam::simpleGeomDecomp::decomposeOneProc
labelList pointIndices(identity(points.size())); labelList pointIndices(identity(points.size()));
const pointField rotatedPoints(rotDelta_ & points); const pointField rotatedPoints(adjustPoints(points));
vectorLessOp sorter(rotatedPoints); vectorLessOp sorter(rotatedPoints);
@ -350,11 +349,11 @@ Foam::labelList Foam::simpleGeomDecomp::decompose
SubField<point>(allPoints, points.size()) = points; SubField<point>(allPoints, points.size()) = points;
nTotalPoints += points.size(); nTotalPoints += points.size();
// Add slaves // Add received
for (const int slave : Pstream::subProcs()) for (const int subproci : Pstream::subProcs())
{ {
IPstream fromSlave(Pstream::commsTypes::scheduled, slave); IPstream fromProc(Pstream::commsTypes::scheduled, subproci);
pointField nbrPoints(fromSlave); pointField nbrPoints(fromProc);
SubField<point> SubField<point>
( (
allPoints, allPoints,
@ -368,14 +367,14 @@ Foam::labelList Foam::simpleGeomDecomp::decompose
labelList finalDecomp(decomposeOneProc(allPoints)); labelList finalDecomp(decomposeOneProc(allPoints));
// Send back // Send back
for (const int slave : Pstream::subProcs()) for (const int subproci : Pstream::subProcs())
{ {
OPstream toSlave(Pstream::commsTypes::scheduled, slave); OPstream toProc(Pstream::commsTypes::scheduled, subproci);
toSlave << SubField<label> toProc << SubField<label>
( (
finalDecomp, finalDecomp,
globalNumbers.localSize(slave), globalNumbers.localSize(subproci),
globalNumbers.offset(slave) globalNumbers.offset(subproci)
); );
} }
// Get my own part // Get my own part
@ -435,12 +434,12 @@ Foam::labelList Foam::simpleGeomDecomp::decompose
SubField<scalar>(allWeights, points.size()) = weights; SubField<scalar>(allWeights, points.size()) = weights;
nTotalPoints += points.size(); nTotalPoints += points.size();
// Add slaves // Add received
for (const int slave : Pstream::subProcs()) for (const int subproci : Pstream::subProcs())
{ {
IPstream fromSlave(Pstream::commsTypes::scheduled, slave); IPstream fromProc(Pstream::commsTypes::scheduled, subproci);
pointField nbrPoints(fromSlave); pointField nbrPoints(fromProc);
scalarField nbrWeights(fromSlave); scalarField nbrWeights(fromProc);
SubField<point> SubField<point>
( (
allPoints, allPoints,
@ -460,14 +459,14 @@ Foam::labelList Foam::simpleGeomDecomp::decompose
labelList finalDecomp(decomposeOneProc(allPoints, allWeights)); labelList finalDecomp(decomposeOneProc(allPoints, allWeights));
// Send back // Send back
for (const int slave : Pstream::subProcs()) for (const int subproci : Pstream::subProcs())
{ {
OPstream toSlave(Pstream::commsTypes::scheduled, slave); OPstream toProc(Pstream::commsTypes::scheduled, subproci);
toSlave << SubField<label> toProc << SubField<label>
( (
finalDecomp, finalDecomp,
globalNumbers.localSize(slave), globalNumbers.localSize(subproci),
globalNumbers.offset(slave) globalNumbers.offset(subproci)
); );
} }
// Get my own part // Get my own part

View File

@ -32,9 +32,11 @@ Description
Method coefficients: Method coefficients:
\table \table
Property | Description | Required | Default Property | Description | Required | Default
n | (nx ny nz) | yes n | (nx ny nz) | yes |
delta | delta for rotation matrix | no | 0.001 order | order of operation (unused) | no | xyz
delta | delta (jitter) for rotation matrix | no | 0.001
transform | cartesian coordinate transformation | no |
\endtable \endtable
SourceFiles SourceFiles
@ -60,7 +62,11 @@ class simpleGeomDecomp
{ {
// Private Member Functions // Private Member Functions
static void assignToProcessorGroup(labelList&, const label); static void assignToProcessorGroup
(
labelList& processorGroup,
const label nProcGroup
);
static void assignToProcessorGroup static void assignToProcessorGroup
( (
@ -80,15 +86,14 @@ class simpleGeomDecomp
) const; ) const;
//- No copy construct
void operator=(const simpleGeomDecomp&) = delete;
//- No copy assignment
simpleGeomDecomp(const simpleGeomDecomp&) = delete;
public: public:
//- No copy construct
void operator=(const simpleGeomDecomp&) = delete;
//- No copy assignment
simpleGeomDecomp(const simpleGeomDecomp&) = delete;
//- Runtime type information //- Runtime type information
TypeName("simple"); TypeName("simple");
@ -151,7 +156,7 @@ public:
//- Explicitly provided connectivity //- Explicitly provided connectivity
virtual labelList decompose virtual labelList decompose
( (
const labelListList& globalCellCells, const labelListList& globalCellCells, // unused
const pointField& cc, const pointField& cc,
const scalarField& cWeights const scalarField& cWeights
) const ) const

View File

@ -0,0 +1,43 @@
/*--------------------------------*- C++ -*----------------------------------*\
| ========= | |
| \\ / F ield | OpenFOAM: The Open Source CFD Toolbox |
| \\ / O peration | Version: v2106 |
| \\ / A nd | Website: www.openfoam.com |
| \\/ M anipulation | |
\*---------------------------------------------------------------------------*/
FoamFile
{
version 2.0;
format ascii;
class dictionary;
location "system";
object decomposeParDict;
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
numberOfSubdomains 8;
method simple;
coeffs
{
n ( 2 2 2 );
// Optional coordinate transformation for sorting
transform
{
origin (-0.15 0.15 0);
rotation
{
type axisAngle;
axis (0 0 1);
angle 44.5;
// Or disabled
//type none;
}
}
}
// ************************************************************************* //