MRG: Incremental foundation merge

This commit is contained in:
Andrew Heather
2016-10-03 08:27:06 +01:00
23 changed files with 94 additions and 64 deletions

View File

@ -97,7 +97,7 @@ Foam::scalarField Foam::cellShapeControl::cellSize
Foam::scalar Foam::cellShapeControl::cellSize(const point& pt) const
{
scalarList bary;
FixedList<scalar, 4> bary;
cellShapeControlMesh::Cell_handle ch;
shapeControlMesh_.barycentricCoords(pt, bary, ch);
@ -172,7 +172,7 @@ Foam::scalar Foam::cellShapeControl::cellSize(const point& pt) const
Foam::tensor Foam::cellShapeControl::cellAlignment(const point& pt) const
{
scalarList bary;
FixedList<scalar, 4> bary;
cellShapeControlMesh::Cell_handle ch;
shapeControlMesh_.barycentricCoords(pt, bary, ch);
@ -244,7 +244,7 @@ void Foam::cellShapeControl::cellSizeAndAlignment
tensor& alignment
) const
{
scalarList bary;
FixedList<scalar, 4> bary;
cellShapeControlMesh::Cell_handle ch;
shapeControlMesh_.barycentricCoords(pt, bary, ch);

View File

@ -450,7 +450,7 @@ Foam::cellShapeControlMesh::~cellShapeControlMesh()
void Foam::cellShapeControlMesh::barycentricCoords
(
const Foam::point& pt,
scalarList& bary,
FixedList<scalar, 4>& bary,
Cell_handle& ch
) const
{

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012-2015 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2012-2016 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -122,7 +122,7 @@ public:
void barycentricCoords
(
const Foam::point& pt,
scalarList& bary,
FixedList<scalar, 4>& bary,
Cell_handle& ch
) const;

View File

@ -94,7 +94,7 @@ Foam::fileControl::~fileControl()
//
//Foam::scalar Foam::fileControl::cellSize(const point& pt) const
//{
// scalarList bary;
// FixedList<scalar, 4> bary;
// Cell_handle ch;
//
// triangulatedMesh_.barycentricCoords(pt, bary, ch);
@ -112,7 +112,7 @@ Foam::fileControl::~fileControl()
////- Return the cell alignment at the given location
//Foam::tensor Foam::fileControl::cellAlignment(const point& pt) const
//{
// scalarList bary;
// FixedList<scalar, 4> bary;
// Cell_handle ch;
//
// triangulatedMesh_.barycentricCoords(pt, bary, ch);
@ -144,7 +144,7 @@ Foam::fileControl::~fileControl()
// tensor& alignment
//) const
//{
// scalarList bary;
// FixedList<scalar, 4> bary;
// Cell_handle ch;
//
// triangulatedMesh_.barycentricCoords(pt, bary, ch);

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2012-2015 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2012-2016 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -118,8 +118,7 @@ Foam::scalar Foam::nonUniformField::interpolate
pts[faceHitByPt[2]]
);
scalarList bary(3, 0.0);
FixedList<scalar, 3> bary;
tri.barycentric(pt, bary);
// return pointCellSize_[pMap[faceHitByPt[0]]]*bary[0]

View File

@ -91,6 +91,11 @@ public:
HashTable<T, label, Hash<label>>(map)
{}
//- Construct from an initializer list
Map(std::initializer_list<Tuple2<label, T>> map)
:
HashTable<T, label, Hash<label>>(map)
{}
};

View File

@ -72,6 +72,15 @@ Foam::SortableList<T>::SortableList(const SortableList<T>& lst)
{}
template<class T>
Foam::SortableList<T>::SortableList(std::initializer_list<T> values)
:
List<T>(values)
{
sort();
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
@ -138,20 +147,27 @@ inline void Foam::SortableList<T>::operator=(const T& t)
template<class T>
inline void Foam::SortableList<T>::operator=(const UList<T>& rhs)
inline void Foam::SortableList<T>::operator=(const UList<T>& lst)
{
List<T>::operator=(rhs);
List<T>::operator=(lst);
indices_.clear();
}
template<class T>
inline void Foam::SortableList<T>::operator=(const SortableList<T>& rhs)
inline void Foam::SortableList<T>::operator=(const SortableList<T>& lst)
{
List<T>::operator=(rhs);
indices_ = rhs.indices();
List<T>::operator=(lst);
indices_ = lst.indices();
}
template<class T>
inline void Foam::SortableList<T>::operator=(std::initializer_list<T> lst)
{
List<T>::operator=(lst);
sort();
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
// ************************************************************************* //

View File

@ -67,33 +67,36 @@ public:
//- Null constructor, sort later (eg, after assignment or transfer)
SortableList();
//- Construct from UList, sorting immediately.
//- Construct from UList, sorting immediately
explicit SortableList(const UList<T>&);
//- Construct from transferred List, sorting immediately.
//- Construct from transferred List, sorting immediately
explicit SortableList(const Xfer<List<T>>&);
//- Construct given size. Sort later on.
//- Construct given size. Sort later on
// The indices remain empty until the list is sorted
explicit SortableList(const label size);
//- Construct given size and initial value. Sort later on.
//- Construct given size and initial value. Sort later on
// The indices remain empty until the list is sorted
SortableList(const label size, const T&);
//- Construct as copy.
//- Construct as copy
SortableList(const SortableList<T>&);
//- Construct from an initializer list, sorting immediately
SortableList(std::initializer_list<T>);
// Member Functions
//- Return the list of sorted indices. Updated every sort.
//- Return the list of sorted indices. Updated every sort
const labelList& indices() const
{
return indices_;
}
//- Return non-const access to the sorted indices. Updated every sort.
//- Return non-const access to the sorted indices. Updated every sort
labelList& indices()
{
return indices_;
@ -121,12 +124,14 @@ public:
//- Assignment of all entries to the given value
inline void operator=(const T&);
//- Assignment to UList operator. Takes linear time.
//- Assignment to UList operator. Takes linear time
inline void operator=(const UList<T>&);
//- Assignment operator. Takes linear time.
//- Assignment operator. Takes linear time
inline void operator=(const SortableList<T>&);
//- Assignment to an initializer list
void operator=(std::initializer_list<T>);
};

View File

@ -244,6 +244,13 @@ Foam::Field<Type>::Field(const UList<Type>& list)
{}
template<class Type>
Foam::Field<Type>::Field(const UIndirectList<Type>& list)
:
List<Type>(list)
{}
#ifndef NoConstructFromTmp
template<class Type>
Foam::Field<Type>::Field(const tmp<Field<Type>>& tf)

View File

@ -123,6 +123,9 @@ public:
//- Construct as copy of a UList\<Type\>
explicit Field(const UList<Type>&);
//- Construct as copy of a UIndirectList\<Type\>
explicit Field(const UIndirectList<Type>&);
//- Construct by transferring the List contents
explicit Field(const Xfer<List<Type>>&);

View File

@ -249,7 +249,7 @@ public:
inline scalar barycentric
(
const point& pt,
List<scalar>& bary
FixedList<scalar, 4>& bary
) const;
//- Return nearest point to p on tetrahedron. Is p itself

View File

@ -316,7 +316,7 @@ template<class Point, class PointRef>
Foam::scalar Foam::tetrahedron<Point, PointRef>::barycentric
(
const point& pt,
List<scalar>& bary
FixedList<scalar, 4>& bary
) const
{
// Reference:
@ -346,8 +346,6 @@ Foam::scalar Foam::tetrahedron<Point, PointRef>::barycentric
vector res = inv(t, detT) & (pt - d_);
bary.setSize(4);
bary[0] = res.x();
bary[1] = res.y();
bary[2] = res.z();

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -248,7 +248,7 @@ public:
inline scalar barycentric
(
const point& pt,
List<scalar>& bary
FixedList<scalar, 3>& bary
) const;
//- Return point intersection with a ray.

View File

@ -280,7 +280,7 @@ template<class Point, class PointRef>
Foam::scalar Foam::triangle<Point, PointRef>::barycentric
(
const point& pt,
List<scalar>& bary
FixedList<scalar, 3>& bary
) const
{
// Reference:
@ -302,13 +302,11 @@ Foam::scalar Foam::triangle<Point, PointRef>::barycentric
{
// Degenerate triangle, returning 1/3 barycentric coordinates.
bary = List<scalar>(3, 1.0/3.0);
bary = FixedList<scalar, 3>(1.0/3.0);
return denom;
}
bary.setSize(3);
bary[1] = (d11*d20 - d01*d21)/denom;
bary[2] = (d00*d21 - d01*d20)/denom;
bary[0] = 1.0 - bary[1] - bary[2];

View File

@ -166,7 +166,7 @@ void Foam::cellPointWeight::findTriangle
{
const tetIndices& tetIs = faceTets[tetI];
List<scalar> triWeights(3);
FixedList<scalar, 3> triWeights;
// Barycentric coordinates of the position
scalar det = tetIs.faceTri(mesh).barycentric(position, triWeights);
@ -234,7 +234,7 @@ void Foam::cellPointWeight::findTriangle
// determinant is suitable. If not, the return from barycentric
// to triWeights is safe.
List<scalar> triWeights(3);
FixedList<scalar, 3> triWeights;
tetIs.faceTri(mesh).barycentric(position, triWeights);
@ -260,9 +260,7 @@ Foam::cellPointWeight::cellPointWeight
const label facei
)
:
celli_(celli),
weights_(4),
faceVertices_(3)
celli_(celli)
{
if (facei < 0)
{

View File

@ -58,10 +58,10 @@ protected:
const label celli_;
//- Weights applied to tet vertices
List<scalar> weights_;
FixedList<scalar, 4> weights_;
//- Face vertex indices
List<label> faceVertices_;
FixedList<label, 3> faceVertices_;
// Protected Member Functions
@ -112,13 +112,13 @@ public:
}
//- Interpolation weights
inline const List<scalar>& weights() const
inline const FixedList<scalar, 4>& weights() const
{
return weights_;
}
//- Interpolation addressing for points on face
inline const List<label>& faceVertices() const
inline const FixedList<label, 3>& faceVertices() const
{
return faceVertices_;
}

View File

@ -31,8 +31,8 @@ inline Type Foam::interpolationCellPoint<Type>::interpolate
const cellPointWeight& cpw
) const
{
const List<scalar>& weights = cpw.weights();
const List<label>& faceVertices = cpw.faceVertices();
const FixedList<scalar, 4>& weights = cpw.weights();
const FixedList<label, 3>& faceVertices = cpw.faceVertices();
Type t = this->psi_[cpw.cell()]*weights[0];
t += psip_[faceVertices[0]]*weights[1];
@ -79,7 +79,7 @@ inline Type Foam::interpolationCellPoint<Type>::interpolate
}
}
List<scalar> weights;
FixedList<scalar, 4> weights;
tetIs.tet(this->pMesh_).barycentric(position, weights);

View File

@ -31,8 +31,8 @@ inline Type Foam::interpolationCellPointWallModified<Type>::interpolate
const cellPointWeightWallModified& cpw
) const
{
const List<scalar>& weights = cpw.weights();
const List<label>& faceVertices = cpw.faceVertices();
const FixedList<scalar, 4>& weights = cpw.weights();
const FixedList<label, 3>& faceVertices = cpw.faceVertices();
Type t = this->psi_[cpw.cell()]*weights[0];
t += this->psip_[faceVertices[0]]*weights[1];

View File

@ -55,9 +55,7 @@ Foam::AveragingMethods::Dual<Type>::Dual
volumeCell_(mesh.V()),
volumeDual_(mesh.nPoints(), 0.0),
dataCell_(FieldField<Field, Type>::operator[](0)),
dataDual_(FieldField<Field, Type>::operator[](1)),
tetVertices_(3),
tetCoordinates_(4)
dataDual_(FieldField<Field, Type>::operator[](1))
{
forAll(this->mesh_.C(), celli)
{
@ -123,7 +121,10 @@ void Foam::AveragingMethods::Dual<Type>::tetGeometry
tetIs.tet(this->mesh_).barycentric(position, tetCoordinates_);
tetCoordinates_ = max(tetCoordinates_, scalar(0));
forAll(tetCoordinates_, i)
{
tetCoordinates_[i] = max(tetCoordinates_[i], scalar(0));
}
}

View File

@ -94,10 +94,10 @@ private:
Field<Type>& dataDual_;
//- Tet vertex labels
mutable List<label> tetVertices_;
mutable FixedList<label, 3> tetVertices_;
//- Tet barycentric coordinates
mutable List<scalar> tetCoordinates_;
mutable FixedList<scalar, 4> tetCoordinates_;
//- Private static member functions

View File

@ -340,7 +340,7 @@ Foam::vector Foam::PackingModels::Implicit<CloudType>::velocityCorrection
const label celli = p.cell();
const label facei = p.tetFace();
const tetIndices tetIs(celli, facei, p.tetPt(), mesh);
List<scalar> tetCoordinates(4);
FixedList<scalar, 4> tetCoordinates;
tetIs.tet(mesh).barycentric(p.position(), tetCoordinates);
// cell velocity

View File

@ -3,7 +3,7 @@
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2014 OpenFOAM Foundation
\\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
\\/ M anipulation | Copyright (C) 2015-2016 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -116,7 +116,7 @@ point offsetSurface::operator()
const triSurface& base = baseSurfPtr_();
const triPointRef baseTri(base[triI].tri(base.points()));
List<scalar> bary;
FixedList<scalar, 3> bary;
baseTri.barycentric(surfacePoint, bary);
const triSurface& offset = offsetSurfPtr_();
@ -135,8 +135,8 @@ point offsetSurface::operator()
);
//- Either return interpolatedPoint or re-project onto surface (since
// snapping might not have do so exactly)
// Either return interpolatedPoint or re-project onto surface (since
// snapping might not have do so exactly)
if (project_)
{

View File

@ -3,7 +3,7 @@
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | Copyright (C) 2011-2016 OpenFOAM Foundation
\\/ M anipulation | Copyright (C) 2015 OpenCFD Ltd.
\\/ M anipulation | Copyright (C) 2015-2016 OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
@ -1240,7 +1240,7 @@ void Foam::isoSurface::trimToBox
dynInterpolatedOldPoints.append(oldPoints);
triPointRef tri(oldTriPoints, oldPoints);
scalarList bary;
FixedList<scalar, 3> bary;
tri.barycentric(pt, bary);
FixedList<scalar, 3> weights;
weights[0] = bary[0];