Skip to content
GitLab
Projects
Groups
Snippets
/
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in / Register
Toggle navigation
Menu
Open sidebar
Development
openfoam
Commits
c1cbfe7a
Commit
c1cbfe7a
authored
May 18, 2017
by
Andrew Heather
Browse files
Merge branch 'develop' of develop.openfoam.com:Development/OpenFOAM-plus into develop
parents
99f31a75
b6dec586
Changes
99
Hide whitespace changes
Inline
Side-by-side
applications/test/DLList/Test-DLList.C
View file @
c1cbfe7a
...
...
@@ -56,7 +56,6 @@ int main(int argc, char *argv[])
Info
<<
"element:"
<<
*
iter
<<
endl
;
}
Info
<<
nl
<<
"And again using the same STL iterator: "
<<
nl
<<
endl
;
forAllIters
(
myList
,
iter
)
...
...
applications/test/HashPtrTable/Test-hashPtrTable.C
View file @
c1cbfe7a
...
...
@@ -36,14 +36,9 @@ void printTable(const HashPtrTable<T>& table)
{
Info
<<
table
.
size
()
<<
nl
<<
"("
<<
nl
;
for
(
typename
HashPtrTable
<
T
>::
const_iterator
iter
=
table
.
cbegin
();
iter
!=
table
.
cend
();
++
iter
)
forAllConstIters
(
table
,
iter
)
{
const
T
*
ptr
=
*
iter
;
const
T
*
ptr
=
iter
.
object
()
;
Info
<<
iter
.
key
()
<<
" = "
;
if
(
ptr
)
{
...
...
@@ -57,6 +52,22 @@ void printTable(const HashPtrTable<T>& table)
}
Info
<<
")"
<<
endl
;
// Values only, with for-range
Info
<<
"values ("
;
for
(
auto
val
:
table
)
{
Info
<<
' '
;
if
(
val
)
{
Info
<<
*
val
;
}
else
{
Info
<<
"nullptr"
;
}
}
Info
<<
" )"
<<
nl
;
}
...
...
@@ -68,7 +79,9 @@ int main()
HashPtrTable
<
double
>
myTable
;
myTable
.
insert
(
"abc"
,
new
double
(
42
.
1
));
myTable
.
insert
(
"def"
,
nullptr
);
myTable
.
insert
(
"ghi"
,
new
double
(
3
.
14159
));
myTable
.
insert
(
"pi"
,
new
double
(
3
.
14159
));
myTable
.
insert
(
"natlog"
,
new
double
(
2
.
718282
));
myTable
.
insert
(
"sqrt2"
,
new
double
(
1
.
414214
));
// Info<< myTable << endl;
printTable
(
myTable
);
...
...
@@ -79,8 +92,20 @@ int main()
printTable
(
copy
);
Info
<<
copy
<<
endl
;
Info
<<
"
\n
erase some existing and non-existing entries"
<<
nl
;
auto
iter
=
myTable
.
find
(
"pi"
);
myTable
.
erase
(
iter
);
iter
=
myTable
.
find
(
"unknownKey"
);
myTable
.
erase
(
iter
);
myTable
.
erase
(
"abc"
);
myTable
.
erase
(
"unknownKey"
);
printTable
(
myTable
);
return
0
;
}
// ************************************************************************* //
applications/test/HashSet/Test-hashSet.C
View file @
c1cbfe7a
...
...
@@ -197,8 +197,12 @@ int main(int argc, char *argv[])
Info
<<
"setD has no 11"
<<
endl
;
}
Info
<<
"setB : "
<<
flatOutput
(
setB
)
<<
endl
;
Info
<<
"setD : "
<<
flatOutput
(
setD
)
<<
endl
;
setD
-=
setB
;
Info
<<
"setD -= setB : "
<<
flatOutput
(
setD
)
<<
endl
;
// This should not work (yet?)
// setD[12] = true;
...
...
applications/test/HashTable/Test-hashTable.C
View file @
c1cbfe7a
...
...
@@ -25,6 +25,9 @@ License
#include
"HashTable.H"
#include
"List.H"
#include
"SortableList.H"
#include
"DynamicList.H"
#include
"FlatOutput.H"
#include
"IOstreams.H"
#include
"IStringStream.H"
#include
"OStringStream.H"
...
...
@@ -163,15 +166,15 @@ int main()
<<
"
\n
table2"
<<
table2
<<
nl
;
Info
<<
"
\n
table3"
<<
table
3
<<
"
\n
clearStorage table
3
... "
;
table
3
.
clearStorage
();
Info
<<
table
3
<<
nl
;
Info
<<
"
\n
table3"
<<
table
2
<<
"
\n
clearStorage table
2
... "
;
table
2
.
clearStorage
();
Info
<<
table
2
<<
nl
;
table1
=
{
{
"ac
a
"
,
3
.
0
},
{
"
aaw
"
,
6
.
0
},
{
"a
b
c"
,
3
.
0
},
{
"
def
"
,
6
.
0
},
{
"acr"
,
8
.
0
},
{
"aec"
,
10
.
0
}
};
...
...
@@ -193,7 +196,43 @@ int main()
// These do not yet work. Issues resolving the distance.
//
// List<scalar> table1vals(table1.begin(), table1.end());
// wordList table1keys(table1.begin(), table1.end());
{
Info
<<
"distance/size: "
<<
std
::
distance
(
table1
.
begin
(),
table1
.
end
())
<<
"/"
<<
table1
.
size
()
<<
" and "
<<
std
::
distance
(
table1
.
keys
().
begin
(),
table1
.
keys
().
end
())
<<
"/"
<<
table1
.
keys
().
size
()
<<
nl
;
SortableList
<
word
>
sortKeys
// DynamicList<word> sortKeys
(
table1
.
keys
().
begin
(),
table1
.
keys
().
end
()
);
Info
<<
"sortKeys: "
<<
flatOutput
(
sortKeys
)
<<
nl
;
}
Info
<<
"
\n
From table1: "
<<
flatOutput
(
table1
.
sortedToc
())
<<
nl
<<
"retain keys: "
<<
flatOutput
(
table3
.
sortedToc
())
<<
nl
;
table1
.
retain
(
table3
);
Info
<<
"-> "
<<
flatOutput
(
table1
.
sortedToc
())
<<
nl
;
Info
<<
"Lookup non-existent"
<<
nl
;
Info
<<
table1
.
lookup
(
"missing-const"
,
1.2345e+6
)
<<
" // const-access"
<<
nl
;
Info
<<
table1
(
"missing-inadvertent"
,
3
.
14159
)
<<
" // (inadvertent?) non-const access"
<<
nl
;
Info
<<
table1
(
"missing-autovivify"
)
<<
" // Known auto-vivification (non-const access)"
<<
nl
;
Info
<<
"
\n
table1: "
<<
table1
<<
endl
;
Info
<<
"
\n
Done
\n
"
;
...
...
applications/test/List/Test-List.C
View file @
c1cbfe7a
...
...
@@ -42,6 +42,7 @@ See also
#include
"vector.H"
#include
"labelRange.H"
#include
"scalarList.H"
#include
"ListOps.H"
#include
"SubList.H"
...
...
@@ -144,6 +145,18 @@ int main(int argc, char *argv[])
labelList
longLabelList
=
identity
(
15
);
// This does not work:
// scalarList slist = identity(15);
//
// More writing, but does work:
scalarList
slist
(
labelRange
::
null
.
begin
(),
labelRange
::
identity
(
15
).
end
()
);
Info
<<
"scalar identity:"
<<
flatOutput
(
slist
)
<<
endl
;
Info
<<
"labels (contiguous="
<<
contiguous
<
label
>
()
<<
")"
<<
nl
;
Info
<<
"normal: "
<<
longLabelList
<<
nl
;
...
...
@@ -220,7 +233,32 @@ int main(int argc, char *argv[])
}
Info
<<
"sub-sorted: "
<<
flatOutput
(
longLabelList
)
<<
nl
;
// Info<<"Slice=" << longLabelList[labelRange(23,5)] << nl;
// construct from a label-range
labelRange
range
(
25
,
15
);
labelList
ident
(
range
.
begin
(),
range
.
end
());
Info
<<
"range-list (label)="
<<
ident
<<
nl
;
List
<
scalar
>
sident
(
range
.
begin
(),
range
.
end
());
Info
<<
"range-list (scalar)="
<<
sident
<<
nl
;
// Sub-ranges also work
List
<
scalar
>
sident2
(
range
(
3
),
range
(
10
));
Info
<<
"range-list (scalar)="
<<
sident2
<<
nl
;
// VERY BAD IDEA: List<scalar> sident3(range(10), range(3));
// This doesn't work, and don't know what it should do anyhow
// List<vector> vident(range.begin(), range.end());
// Info<<"range-list (vector)=" << vident << nl;
// Even weird things like this
List
<
scalar
>
sident4
(
labelRange
().
begin
(),
labelRange
::
identity
(
8
).
end
()
);
Info
<<
"range-list (scalar)="
<<
sident4
<<
nl
;
}
wordReList
reLst
;
...
...
applications/test/labelRanges/Test-labelRanges.C
View file @
c1cbfe7a
...
...
@@ -58,7 +58,7 @@ int main(int argc, char *argv[])
Info
<<
"test sorting"
<<
endl
;
DynamicList
<
labelRange
>
list1
(
10
);
list1
.
append
(
labelRange
(
25
,
8
));
list1
.
append
(
labelRange
(
0
,
10
));
list1
.
append
(
labelRange
::
identity
(
8
));
list1
.
append
(
labelRange
(
15
,
5
));
list1
.
append
(
labelRange
(
50
,
-
10
));
...
...
applications/utilities/mesh/conversion/ideasUnvToFoam/ideasUnvToFoam.C
View file @
c1cbfe7a
...
...
@@ -924,8 +924,8 @@ int main(int argc, char *argv[])
}
}
HashTable
<
labelList
,
word
>
cellZones
;
HashTable
<
labelList
,
word
>
faceZones
;
HashTable
<
labelList
>
cellZones
;
HashTable
<
labelList
>
faceZones
;
List
<
bool
>
isAPatch
(
patchNames
.
size
(),
true
);
if
(
dofVertIndices
.
size
())
...
...
applications/utilities/mesh/manipulation/checkMesh/Make/files
View file @
c1cbfe7a
writeFields.C
checkTools.C
checkTopology.C
checkGeometry.C
...
...
applications/utilities/mesh/manipulation/checkMesh/checkMesh.C
View file @
c1cbfe7a
...
...
@@ -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
-2017
OpenCFD Ltd.
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
...
...
@@ -57,12 +57,18 @@ Usage
Reconstruct all cellSets and faceSets geometry and write to postProcessing/
directory according to surfaceFormat (e.g. vtk or ensight)
\param -writeAllFields \n
Writes all mesh quality measures as fields.
\param -writeFields '(\<fieldName\>)' \n
Writes selected mesh quality measures as fields.
\*---------------------------------------------------------------------------*/
#include
"argList.H"
#include
"timeSelector.H"
#include
"Time.H"
#include
"
poly
Mesh.H"
#include
"
fv
Mesh.H"
#include
"globalMeshData.H"
#include
"surfaceWriter.H"
#include
"vtkSetWriter.H"
...
...
@@ -71,6 +77,7 @@ Usage
#include
"checkTopology.H"
#include
"checkGeometry.H"
#include
"checkMeshQuality.H"
#include
"writeFields.H"
using
namespace
Foam
;
...
...
@@ -96,6 +103,17 @@ int main(int argc, char *argv[])
"include extra topology checks"
);
argList
::
addBoolOption
(
"writeAllFields"
,
"write volFields with mesh quality parameters"
);
argList
::
addOption
(
"writeFields"
,
"wordList"
,
"write volFields with selected mesh quality parameters"
);
argList
::
addBoolOption
(
"meshQuality"
,
"read user-defined mesh quality criterions from system/meshQualityDict"
...
...
@@ -110,7 +128,7 @@ int main(int argc, char *argv[])
#include
"setRootCase.H"
#include
"createTime.H"
instantList
timeDirs
=
timeSelector
::
select0
(
runTime
,
args
);
#include
"createNamed
Poly
Mesh.H"
#include
"createNamedMesh.H"
const
bool
noTopology
=
args
.
optionFound
(
"noTopology"
);
const
bool
allGeometry
=
args
.
optionFound
(
"allGeometry"
);
...
...
@@ -119,6 +137,24 @@ int main(int argc, char *argv[])
word
surfaceFormat
;
const
bool
writeSets
=
args
.
optionReadIfPresent
(
"writeSets"
,
surfaceFormat
);
HashSet
<
word
>
selectedFields
;
bool
writeFields
=
args
.
optionReadIfPresent
(
"writeFields"
,
selectedFields
);
if
(
!
writeFields
&&
args
.
optionFound
(
"writeAllFields"
))
{
selectedFields
.
insert
(
"nonOrthoAngle"
);
selectedFields
.
insert
(
"faceWeight"
);
selectedFields
.
insert
(
"skewness"
);
selectedFields
.
insert
(
"cellDeterminant"
);
selectedFields
.
insert
(
"aspectRatio"
);
selectedFields
.
insert
(
"cellShapes"
);
selectedFields
.
insert
(
"cellVolume"
);
selectedFields
.
insert
(
"cellVolumeRatio"
);
}
if
(
noTopology
)
{
...
...
@@ -143,6 +179,11 @@ int main(int argc, char *argv[])
<<
" representation"
<<
" of all faceSets and cellSets."
<<
nl
<<
endl
;
}
if
(
selectedFields
.
size
())
{
Info
<<
"Writing mesh quality as fields "
<<
selectedFields
<<
nl
<<
endl
;
}
autoPtr
<
IOdictionary
>
qualDict
;
...
...
@@ -234,6 +275,10 @@ int main(int argc, char *argv[])
Info
<<
"
\n
Failed "
<<
nFailedChecks
<<
" mesh checks.
\n
"
<<
endl
;
}
// Write selected fields
Foam
::
writeFields
(
mesh
,
selectedFields
);
}
else
if
(
state
==
polyMesh
::
POINTS_MOVED
)
{
...
...
@@ -262,6 +307,10 @@ int main(int argc, char *argv[])
{
Info
<<
"
\n
Mesh OK.
\n
"
<<
endl
;
}
// Write selected fields
Foam
::
writeFields
(
mesh
,
selectedFields
);
}
}
...
...
applications/utilities/mesh/manipulation/checkMesh/writeFields.C
0 → 100644
View file @
c1cbfe7a
#include
"writeFields.H"
#include
"volFields.H"
#include
"polyMeshTools.H"
#include
"zeroGradientFvPatchFields.H"
#include
"syncTools.H"
using
namespace
Foam
;
void
maxFaceToCell
(
const
scalarField
&
faceData
,
volScalarField
&
cellData
)
{
const
cellList
&
cells
=
cellData
.
mesh
().
cells
();
scalarField
&
cellFld
=
cellData
.
ref
();
cellFld
=
-
GREAT
;
forAll
(
cells
,
cellI
)
{
const
cell
&
cFaces
=
cells
[
cellI
];
forAll
(
cFaces
,
i
)
{
cellFld
[
cellI
]
=
max
(
cellFld
[
cellI
],
faceData
[
cFaces
[
i
]]);
}
}
forAll
(
cellData
.
boundaryField
(),
patchI
)
{
fvPatchScalarField
&
fvp
=
cellData
.
boundaryFieldRef
()[
patchI
];
fvp
=
fvp
.
patch
().
patchSlice
(
faceData
);
}
cellData
.
correctBoundaryConditions
();
}
void
minFaceToCell
(
const
scalarField
&
faceData
,
volScalarField
&
cellData
)
{
const
cellList
&
cells
=
cellData
.
mesh
().
cells
();
scalarField
&
cellFld
=
cellData
.
ref
();
cellFld
=
GREAT
;
forAll
(
cells
,
cellI
)
{
const
cell
&
cFaces
=
cells
[
cellI
];
forAll
(
cFaces
,
i
)
{
cellFld
[
cellI
]
=
min
(
cellFld
[
cellI
],
faceData
[
cFaces
[
i
]]);
}
}
forAll
(
cellData
.
boundaryField
(),
patchI
)
{
fvPatchScalarField
&
fvp
=
cellData
.
boundaryFieldRef
()[
patchI
];
fvp
=
fvp
.
patch
().
patchSlice
(
faceData
);
}
cellData
.
correctBoundaryConditions
();
}
void
Foam
::
writeFields
(
const
fvMesh
&
mesh
,
const
HashSet
<
word
>&
selectedFields
)
{
if
(
selectedFields
.
empty
())
{
return
;
}
Info
<<
"Writing fields with mesh quality parameters"
<<
endl
;
if
(
selectedFields
.
found
(
"nonOrthoAngle"
))
{
//- Face based orthogonality
const
scalarField
faceOrthogonality
(
polyMeshTools
::
faceOrthogonality
(
mesh
,
mesh
.
faceAreas
(),
mesh
.
cellCentres
()
)
);
//- Face based angle
const
scalarField
nonOrthoAngle
(
radToDeg
(
Foam
::
acos
(
min
(
1
.
0
,
faceOrthogonality
))
)
);
//- Cell field - max of either face
volScalarField
cellNonOrthoAngle
(
IOobject
(
"nonOrthoAngle"
,
mesh
.
time
().
timeName
(),
mesh
,
IOobject
::
NO_READ
,
IOobject
::
AUTO_WRITE
),
mesh
,
dimensionedScalar
(
"angle"
,
dimless
,
0
),
calculatedFvPatchScalarField
::
typeName
);
//- Take max
maxFaceToCell
(
nonOrthoAngle
,
cellNonOrthoAngle
);
Info
<<
" Writing non-orthogonality (angle) to "
<<
cellNonOrthoAngle
.
name
()
<<
endl
;
cellNonOrthoAngle
.
write
();
}
if
(
selectedFields
.
found
(
"faceWeight"
))
{
const
scalarField
faceWeights
(
polyMeshTools
::
faceWeights
(
mesh
,
mesh
.
faceCentres
(),
mesh
.
faceAreas
(),
mesh
.
cellCentres
()
)
);
volScalarField
cellWeights
(
IOobject
(
"faceWeight"
,
mesh
.
time
().
timeName
(),
mesh
,
IOobject
::
NO_READ
,
IOobject
::
AUTO_WRITE
),
mesh
,
dimensionedScalar
(
"weight"
,
dimless
,
0
),
calculatedFvPatchScalarField
::
typeName
);
//- Take min
minFaceToCell
(
faceWeights
,
cellWeights
);
Info
<<
" Writing face interpolation weights (0..0.5) to "
<<
cellWeights
.
name
()
<<
endl
;
cellWeights
.
write
();
}
// Skewness
// ~~~~~~~~
if
(
selectedFields
.
found
(
"skewness"
))
{
//- Face based skewness
const
scalarField
faceSkewness
(
polyMeshTools
::
faceSkewness
(
mesh
,
mesh
.
points
(),
mesh
.
faceCentres
(),
mesh
.
faceAreas
(),
mesh
.
cellCentres
()
)
);
//- Cell field - max of either face
volScalarField
cellSkewness
(
IOobject
(
"skewness"
,
mesh
.
time
().
timeName
(),
mesh
,
IOobject
::
NO_READ
,
IOobject
::
AUTO_WRITE
),
mesh
,
dimensionedScalar
(
"skewness"
,
dimless
,
0
),
calculatedFvPatchScalarField
::
typeName
);
//- Take max
maxFaceToCell
(
faceSkewness
,
cellSkewness
);
Info
<<
" Writing face skewness to "
<<
cellSkewness
.
name
()
<<
endl
;
cellSkewness
.
write
();
}
// cellDeterminant
// ~~~~~~~~~~~~~~~
if
(
selectedFields
.
found
(
"cellDeterminant"
))
{
volScalarField
cellDeterminant
(
IOobject
(
"cellDeterminant"
,
mesh
.
time
().
timeName
(),
mesh
,
IOobject
::
NO_READ
,
IOobject
::
AUTO_WRITE
,
false
),
mesh
,
dimensionedScalar
(
"cellDeterminant"
,
dimless
,
0
),
zeroGradientFvPatchScalarField
::
typeName
);
cellDeterminant
.
primitiveFieldRef
()
=
primitiveMeshTools
::
cellDeterminant
(
mesh
,
mesh
.
geometricD
(),
mesh
.
faceAreas
(),