ENH: improvements for makeFaMesh, checkFaMesh

- added -dry-run, -write-vtk options.
  Additional mesh information after creation.

- add parallel reductions and more information for checkFaMesh

ENH: minor cleanup of faPatch internals

- align pointLabels and pointEdges creation more closely with coding
  patterns used in PrimitivePatch

- use fileHandler when loading "S0" field.
This commit is contained in:
Mark Olesen
2021-10-12 08:47:18 +02:00
committed by Andrew Heather
parent a95427c28a
commit 7fc943c178
17 changed files with 579 additions and 190 deletions

View File

@ -1,9 +1,11 @@
EXE_INC = \ EXE_INC = \
-I$(LIB_SRC)/finiteArea/lnInclude \ -I$(LIB_SRC)/finiteArea/lnInclude \
-I$(LIB_SRC)/finiteVolume/lnInclude \ -I$(LIB_SRC)/fileFormats/lnInclude \
-I$(LIB_SRC)/meshTools/lnInclude \ -I$(LIB_SRC)/surfMesh/lnInclude \
-I$(LIB_SRC)/meshTools/lnInclude
EXE_LIBS = \ EXE_LIBS = \
-lfiniteArea \ -lfiniteArea \
-lfiniteVolume \ -lfileFormats \
-lsurfMesh \
-lmeshTools -lmeshTools

View File

@ -6,6 +6,7 @@
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
Copyright (C) 2016-2017 Wikki Ltd Copyright (C) 2016-2017 Wikki Ltd
Copyright (C) 2021 OpenCFD Ltd.
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
This file is part of OpenFOAM. This file is part of OpenFOAM.
@ -29,14 +30,22 @@ Application
Description Description
Check a finiteArea mesh Check a finiteArea mesh
Author Original Authors
Zeljko Tukovic, FAMENA Zeljko Tukovic, FAMENA
Hrvoje Jasak, Wikki Ltd. Hrvoje Jasak, Wikki Ltd.
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
#include "fvCFD.H" #include "Time.H"
#include "faCFD.H" #include "argList.H"
#include "faMesh.H"
#include "polyMesh.H"
#include "areaFaMesh.H"
#include "edgeFaMesh.H"
#include "areaFields.H"
#include "edgeFields.H"
#include "processorFaPatch.H"
#include "foamVtkUIndPatchWriter.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
@ -51,37 +60,32 @@ int main(int argc, char *argv[])
"Check a finiteArea mesh" "Check a finiteArea mesh"
); );
#include "addRegionOption.H" argList::addBoolOption
(
"write-vtk",
"Write mesh as a vtp (vtk) file for display or debugging"
);
#include "addRegionOption.H"
#include "setRootCase.H" #include "setRootCase.H"
#include "createTime.H" #include "createTime.H"
#include "createNamedMesh.H" #include "createNamedPolyMesh.H"
#include "createFaMesh.H"
// Create
faMesh aMesh(mesh);
Info<< "Time = " << runTime.timeName() << nl << endl; Info<< "Time = " << runTime.timeName() << nl << endl;
// General mesh statistics #include "printMeshSummary.H"
Info<< "Number of points: " << aMesh.nPoints() << nl
<< "Number of internal edges: " << aMesh.nInternalEdges() << nl
<< "Number of edges: " << aMesh.nEdges() << nl
<< "Number of faces: " << aMesh.nFaces() << nl
<< endl;
// Check geometry if (args.found("write-vtk"))
Info<< "Face area: min = " << min(aMesh.S().field()) {
<< " max = " << max(aMesh.S().field()) << nl #include "faMeshWriteVTK.H"
<< "Internal edge length: min = " }
<< min(aMesh.magLe().internalField()) << nl
<< " max = " << max(aMesh.magLe().internalField()) << nl
<< "Edge length: min = "
<< min(aMesh.magLe()).value() << nl
<< " max = " << max(aMesh.magLe()).value() << nl
<< "Face area normals: min = " << min(aMesh.faceAreaNormals().field())
<< " max = " << max(aMesh.faceAreaNormals().field()) << nl
<< endl;
Info<< "\nEnd\n" << endl;
Info << "\nEnd" << endl;
return 0; return 0;
} }

View File

@ -0,0 +1,47 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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, distributed under GPL-3.0-or-later.
Description
VTK output of faMesh with some geometric or debug fields
\*---------------------------------------------------------------------------*/
{
// finiteArea
vtk::uindirectPatchWriter writer
(
aMesh.patch(),
// vtk::formatType::INLINE_ASCII,
fileName
(
aMesh.mesh().time().globalPath() / "finiteArea"
)
);
writer.writeGeometry();
// CellData
writer.beginCellData(3);
writer.writeProcIDs();
writer.write("area", aMesh.S().field());
writer.write("normal", aMesh.faceAreaNormals());
// PointData
writer.beginPointData(1);
writer.write("normal", aMesh.pointAreaNormals());
Info<< nl
<< "Wrote faMesh in vtk format: " << writer.output().name() << nl;
}
// ************************************************************************* //

View File

@ -0,0 +1,134 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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, distributed under GPL-3.0-or-later.
Description
Summary of faMesh information
\*---------------------------------------------------------------------------*/
{
const faBoundaryMesh& patches = aMesh.boundary();
const label nNonProcessor = patches.nNonProcessor();
const label nPatches = patches.size();
Info<< "----------------" << nl
<< "Mesh Information" << nl
<< "----------------" << nl
<< " " << "boundingBox: " << boundBox(aMesh.points()) << nl;
Info<< " Number of points: "
<< returnReduce(aMesh.nPoints(), sumOp<label>()) << nl
<< " Number of faces: "
<< returnReduce(aMesh.nFaces(), sumOp<label>()) << nl;
Info<< " Number of edges: "
<< returnReduce(aMesh.nEdges(), sumOp<label>()) << nl
<< " Number of internal edges: "
<< returnReduce(aMesh.nInternalEdges(), sumOp<label>()) << nl;
label nProcEdges = 0;
if (Pstream::parRun())
{
for (const faPatch& fap : patches)
{
const auto* cpp = isA<processorFaPatch>(fap);
if (cpp)
{
nProcEdges += fap.nEdges();
}
}
}
const label nBoundEdges = aMesh.nBoundaryEdges() - nProcEdges;
Info<< " Number of boundary edges: "
<< returnReduce(nBoundEdges - nProcEdges, sumOp<label>()) << nl;
if (Pstream::parRun())
{
Info<< " Number of processor edges: "
<< returnReduce(nProcEdges, sumOp<label>()) << nl;
}
Info<< "----------------" << nl
<< "Patches" << nl
<< "----------------" << nl;
for (label patchi = 0; patchi < nNonProcessor; ++patchi)
{
const faPatch& p = patches[patchi];
Info<< " " << "patch " << p.index()
<< " (size: " << returnReduce(p.size(), sumOp<label>())
<< ") name: " << p.name()
<< nl;
}
// Geometry information
Info<< nl;
{
scalarMinMax limit(gMinMax(aMesh.S().field()));
Info<< "Face area:" << nl
<< " min = " << limit.min() << " max = " << limit.max() << nl;
}
{
scalarMinMax limit(minMax(aMesh.magLe().primitiveField()));
// Include processor boundaries into 'internal' edges
if (Pstream::parRun())
{
for (label patchi = nNonProcessor; patchi < nPatches; ++patchi)
{
limit.add(minMax(aMesh.magLe().boundaryField()[patchi]));
}
reduce(limit, minMaxOp<scalar>());
}
Info<< "Edge length (internal):" << nl
<< " min = " << limit.min() << " max = " << limit.max() << nl;
// Include (non-processor) boundaries
for (label patchi = 0; patchi < nNonProcessor; ++patchi)
{
limit.add(minMax(aMesh.magLe().boundaryField()[patchi]));
}
if (Pstream::parRun())
{
reduce(limit, minMaxOp<scalar>());
}
Info<< "Edge length:" << nl
<< " min = " << limit.min() << " max = " << limit.max() << nl;
}
// Not particularly meaningful
#if 0
{
MinMax<vector> limit(gMinMax(aMesh.faceAreaNormals().field()));
Info<< "Face area normals:" << nl
<< " min = " << limit.min() << " max = " << limit.max() << nl;
}
#endif
}
// ************************************************************************* //

View File

@ -1,6 +1,7 @@
EXE_INC = \ EXE_INC = \
-I$(LIB_SRC)/finiteArea/lnInclude \ -I$(LIB_SRC)/finiteArea/lnInclude \
-I$(LIB_SRC)/fileFormats/lnInclude \ -I$(LIB_SRC)/fileFormats/lnInclude \
-I$(LIB_SRC)/surfMesh/lnInclude \
-I$(LIB_SRC)/meshTools/lnInclude \ -I$(LIB_SRC)/meshTools/lnInclude \
-I$(LIB_SRC)/parallel/decompose/faDecompose/lnInclude \ -I$(LIB_SRC)/parallel/decompose/faDecompose/lnInclude \
-I$(LIB_SRC)/parallel/reconstruct/faReconstruct/lnInclude -I$(LIB_SRC)/parallel/reconstruct/faReconstruct/lnInclude
@ -8,6 +9,7 @@ EXE_INC = \
EXE_LIBS = \ EXE_LIBS = \
-lfiniteArea \ -lfiniteArea \
-lfileFormats \ -lfileFormats \
-lsurfMesh \
-lmeshTools \ -lmeshTools \
-lfaDecompose \ -lfaDecompose \
-lfaReconstruct -lfaReconstruct

View File

@ -17,9 +17,11 @@ Description
if (Pstream::parRun()) if (Pstream::parRun())
{ {
faMeshReconstructor reconstructor(areaMesh); faMeshReconstructor reconstructor(aMesh);
reconstructor.writeAddressing(); reconstructor.writeAddressing();
Info<< "Wrote proc-addressing" << nl << endl;
// Handle area fields // Handle area fields
// ------------------ // ------------------
@ -72,22 +74,53 @@ if (Pstream::parRun())
if (nAreaFields) if (nAreaFields)
{ {
Info<< "Decomposing " << nAreaFields << " area fields" << nl << endl; Info<< "Decomposing " << nAreaFields << " area fields" << nl;
faFieldDecomposer fieldDecomposer faFieldDecomposer fieldDecomposer
( (
fullMesh, fullMesh,
areaMesh, aMesh,
reconstructor.edgeProcAddressing(), reconstructor.edgeProcAddressing(),
reconstructor.faceProcAddressing(), reconstructor.faceProcAddressing(),
reconstructor.boundaryProcAddressing() reconstructor.boundaryProcAddressing()
); );
fieldDecomposer.decomposeFields(areaScalarFields); if (areaScalarFields.size())
fieldDecomposer.decomposeFields(areaVectorFields); {
fieldDecomposer.decomposeFields(areaSphTensorFields); Info<< " scalars: "
fieldDecomposer.decomposeFields(areaSymmTensorFields); << flatOutput(PtrListOps::names(areaScalarFields)) << nl;
fieldDecomposer.decomposeFields(areaTensorFields); fieldDecomposer.decomposeFields(areaScalarFields);
}
if (areaVectorFields.size())
{
Info<< " vectors: "
<< flatOutput(PtrListOps::names(areaVectorFields)) << nl;
fieldDecomposer.decomposeFields(areaVectorFields);
}
if (areaSphTensorFields.size())
{
Info<< " sphTensors: "
<< flatOutput(PtrListOps::names(areaSphTensorFields)) << nl;
fieldDecomposer.decomposeFields(areaSphTensorFields);
}
if (areaSymmTensorFields.size())
{
Info<< " symmTensors: "
<< flatOutput(PtrListOps::names(areaSymmTensorFields)) << nl;
fieldDecomposer.decomposeFields(areaSymmTensorFields);
}
if (areaTensorFields.size())
{
Info<< " tensors: "
<< flatOutput(PtrListOps::names(areaTensorFields)) << nl;
fieldDecomposer.decomposeFields(areaTensorFields);
}
Info<< endl;
} }
} }

View File

@ -16,15 +16,18 @@ Description
\*---------------------------------------------------------------------------*/ \*---------------------------------------------------------------------------*/
{ {
Info<< "Writing edges in obj format" << endl; Info<< nl
<< "Writing edges in obj format" << endl;
word outputName("faMesh-edges.obj"); word outputName("finiteArea-edges.obj");
if (Pstream::parRun()) if (Pstream::parRun())
{ {
outputName = word outputName = word
( (
"faMesh-edges-" + Foam::name(Pstream::myProcNo()) + ".obj" "finiteArea-edges-proc"
+ Foam::name(Pstream::myProcNo())
+ ".obj"
); );
} }
@ -36,7 +39,7 @@ Description
false false
); );
os.write(areaMesh.patch().edges(), areaMesh.patch().localPoints()); os.write(aMesh.patch().edges(), aMesh.patch().localPoints());
} }

View File

@ -0,0 +1,47 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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, distributed under GPL-3.0-or-later.
Description
VTK output of faMesh with some geometric or debug fields
\*---------------------------------------------------------------------------*/
{
// finiteArea
vtk::uindirectPatchWriter writer
(
aMesh.patch(),
// vtk::formatType::INLINE_ASCII,
fileName
(
aMesh.mesh().time().globalPath() / "finiteArea"
)
);
writer.writeGeometry();
// CellData
writer.beginCellData(3);
writer.writeProcIDs();
writer.write("area", aMesh.S().field());
writer.write("normal", aMesh.faceAreaNormals());
// PointData
writer.beginPointData(1);
writer.write("normal", aMesh.pointAreaNormals());
Info<< nl
<< "Wrote faMesh in vtk format: " << writer.output().name() << nl;
}
// ************************************************************************* //

View File

@ -32,7 +32,7 @@ Description
When called in parallel, it will also try to act like decomposePar, When called in parallel, it will also try to act like decomposePar,
create procAddressing and decompose serial finite-area fields. create procAddressing and decompose serial finite-area fields.
Author Original Authors
Zeljko Tukovic, FAMENA Zeljko Tukovic, FAMENA
Hrvoje Jasak, Wikki Ltd. Hrvoje Jasak, Wikki Ltd.
@ -48,6 +48,8 @@ Author
#include "areaFields.H" #include "areaFields.H"
#include "faFieldDecomposer.H" #include "faFieldDecomposer.H"
#include "faMeshReconstructor.H" #include "faMeshReconstructor.H"
#include "PtrListOps.H"
#include "foamVtkUIndPatchWriter.H"
#include "OBJstream.H" #include "OBJstream.H"
using namespace Foam; using namespace Foam;
@ -69,11 +71,21 @@ int main(int argc, char *argv[])
); );
argList::addOption("dict", "file", "Alternative faMeshDefinition"); argList::addOption("dict", "file", "Alternative faMeshDefinition");
argList::addBoolOption
(
"dry-run",
"Create but do not write"
);
argList::addBoolOption
(
"write-vtk",
"Write mesh as a vtp (vtk) file for display or debugging"
);
argList::addBoolOption argList::addBoolOption
( (
"write-edges-obj", "write-edges-obj",
"Write mesh edges as obj files and exit", "Write mesh edges as obj files (one per processor)",
false // could make an advanced option true // advanced option (debugging only)
); );
#include "addRegionOption.H" #include "addRegionOption.H"
@ -81,6 +93,8 @@ int main(int argc, char *argv[])
#include "createTime.H" #include "createTime.H"
#include "createNamedPolyMesh.H" #include "createNamedPolyMesh.H"
const bool dryrun = args.found("dry-run");
// Reading faMeshDefinition dictionary // Reading faMeshDefinition dictionary
#include "findMeshDefinitionDict.H" #include "findMeshDefinitionDict.H"
@ -91,33 +105,40 @@ int main(int argc, char *argv[])
meshDefDict.add("emptyPatch", patchName, true); meshDefDict.add("emptyPatch", patchName, true);
} }
// Create
faMesh areaMesh(mesh, meshDefDict);
bool quickExit = false; // Create
faMesh aMesh(mesh, meshDefDict);
// Mesh information
#include "printMeshSummary.H"
if (args.found("write-edges-obj")) if (args.found("write-edges-obj"))
{ {
quickExit = true;
#include "faMeshWriteEdgesOBJ.H" #include "faMeshWriteEdgesOBJ.H"
} }
if (quickExit) if (args.found("write-vtk"))
{ {
Info<< "\nEnd\n" << endl; #include "faMeshWriteVTK.H"
return 0;
} }
// Set the precision of the points data to 10 if (dryrun)
IOstream::defaultPrecision(10); {
Info<< "\ndry-run: not writing mesh or decomposing fields\n" << nl;
}
else
{
// Set the precision of the points data to 10
IOstream::defaultPrecision(10);
Info<< nl << "Write finite area mesh." << nl; Info<< nl << "Write finite area mesh." << nl;
areaMesh.write(); aMesh.write();
Info<< endl;
#include "decomposeFaFields.H" Info<< endl;
#include "decomposeFaFields.H"
}
Info << "\nEnd\n" << endl; Info<< "\nEnd\n" << endl;
return 0; return 0;
} }

View File

@ -0,0 +1,99 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / 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, distributed under GPL-3.0-or-later.
Description
Summary of faMesh information
\*---------------------------------------------------------------------------*/
{
const faBoundaryMesh& patches = aMesh.boundary();
const label nNonProcessor = patches.nNonProcessor();
const label nPatches = patches.size();
Info<< "----------------" << nl
<< "Mesh Information" << nl
<< "----------------" << nl
<< " " << "boundingBox: " << boundBox(aMesh.points()) << nl
<< " " << "nFaces: " << returnReduce(aMesh.nFaces(), sumOp<label>())
<< nl;
Info<< "----------------" << nl
<< "Patches" << nl
<< "----------------" << nl;
for (label patchi = 0; patchi < nNonProcessor; ++patchi)
{
const faPatch& p = patches[patchi];
Info<< " " << "patch " << p.index()
<< " (size: " << returnReduce(p.size(), sumOp<label>())
<< ") name: " << p.name()
<< nl;
}
// Geometry information
Info<< nl;
{
scalarMinMax limit(gMinMax(aMesh.S().field()));
Info<< "Face area:" << nl
<< " min = " << limit.min() << " max = " << limit.max() << nl;
}
{
scalarMinMax limit(minMax(aMesh.magLe().primitiveField()));
// Include processor boundaries into 'internal' edges
if (Pstream::parRun())
{
for (label patchi = nNonProcessor; patchi < nPatches; ++patchi)
{
limit.add(minMax(aMesh.magLe().boundaryField()[patchi]));
}
reduce(limit, minMaxOp<scalar>());
}
Info<< "Edge length (internal):" << nl
<< " min = " << limit.min() << " max = " << limit.max() << nl;
// Include (non-processor) boundaries
for (label patchi = 0; patchi < nNonProcessor; ++patchi)
{
limit.add(minMax(aMesh.magLe().boundaryField()[patchi]));
}
if (Pstream::parRun())
{
reduce(limit, minMaxOp<scalar>());
}
Info<< "Edge length:" << nl
<< " min = " << limit.min()
<< " max = " << limit.max() << nl;
}
// Not particularly meaningful
#if 0
{
MinMax<vector> limit(gMinMax(aMesh.faceAreaNormals().field()));
Info<< "Face area normals:" << nl
<< " min = " << limit.min() << " max = " << limit.max() << nl;
}
#endif
}
// ************************************************************************* //

View File

@ -75,8 +75,15 @@ Foam::faBoundaryMesh::faBoundaryMesh
regIOobject(io), regIOobject(io),
mesh_(mesh) mesh_(mesh)
{ {
if (readOpt() == IOobject::MUST_READ) if
(
readOpt() == IOobject::MUST_READ
|| readOpt() == IOobject::MUST_READ_IF_MODIFIED
)
{ {
// Warn for MUST_READ_IF_MODIFIED
warnNoRereading<faBoundaryMesh>();
faPatchList& patches = *this; faPatchList& patches = *this;
// Read faPatch list // Read faPatch list

View File

@ -343,7 +343,7 @@ Foam::faMesh::faMesh(const polyMesh& pMesh)
// Calculate the geometry for the patches (transformation tensors etc.) // Calculate the geometry for the patches (transformation tensors etc.)
boundary_.calcGeometry(); boundary_.calcGeometry();
if (isFile(pMesh.time().timePath()/mesh().dbDir()/"S0")) if (fileHandler().isFile(pMesh.time().timePath()/"S0"))
{ {
S0Ptr_ = new DimensionedField<scalar, areaMesh> S0Ptr_ = new DimensionedField<scalar, areaMesh>
( (
@ -508,7 +508,7 @@ Foam::faMesh::faMesh
// Calculate the geometry for the patches (transformation tensors etc.) // Calculate the geometry for the patches (transformation tensors etc.)
boundary_.calcGeometry(); boundary_.calcGeometry();
if (isFile(mesh().time().timePath()/"S0")) if (fileHandler().isFile(pMesh.time().timePath()/"S0"))
{ {
S0Ptr_ = new DimensionedField<scalar, areaMesh> S0Ptr_ = new DimensionedField<scalar, areaMesh>
( (

View File

@ -6,7 +6,7 @@
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
Copyright (C) 2016-2017 Wikki Ltd Copyright (C) 2016-2017 Wikki Ltd
Copyright (C) 2019 OpenCFD Ltd. Copyright (C) 2019-2021 OpenCFD Ltd.
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
This file is part of OpenFOAM. This file is part of OpenFOAM.
@ -47,34 +47,33 @@ const Foam::scalar Foam::cyclicFaPatch::matchTol_ = 1e-3;
void Foam::cyclicFaPatch::calcTransforms() void Foam::cyclicFaPatch::calcTransforms()
{ {
const label sizeby2 = this->size()/2;
if (size() > 0) if (size() > 0)
{ {
// const label sizeby2 = this->size()/2; pointField half0Ctrs(sizeby2);
pointField half0Ctrs(size()/2); pointField half1Ctrs(sizeby2);
pointField half1Ctrs(size()/2); for (label edgei=0; edgei < sizeby2; ++edgei)
for (label i=0; i<size()/2; ++i)
{ {
half0Ctrs[i] = this->edgeCentres()[i]; half0Ctrs[edgei] = this->edgeCentres()[edgei];
half1Ctrs[i] = this->edgeCentres()[i+size()/2]; half1Ctrs[edgei] = this->edgeCentres()[edgei + sizeby2];
} }
vectorField half0Normals(size()/2); vectorField half0Normals(sizeby2);
vectorField half1Normals(size()/2); vectorField half1Normals(sizeby2);
const vectorField eN(edgeNormals()*magEdgeLengths()); const vectorField eN(edgeNormals()*magEdgeLengths());
scalar maxMatchError = 0; scalar maxMatchError = 0;
label errorEdge = -1; label errorEdge = -1;
for (label edgei = 0; edgei < size()/2; ++edgei) for (label edgei = 0; edgei < sizeby2; ++edgei)
{ {
half0Normals[edgei] = eN[edgei]; half0Normals[edgei] = eN[edgei];
label nbrEdgei = edgei + size()/2; half1Normals[edgei] = eN[edgei + sizeby2];
half1Normals[edgei] = eN[nbrEdgei];
scalar magLe = mag(half0Normals[edgei]); scalar magLe = mag(half0Normals[edgei]);
scalar nbrMagLe = mag(half1Normals[edgei]); scalar nbrMagLe = mag(half1Normals[edgei]);
scalar avLe = (magLe + nbrMagLe)/2.0; scalar avLe = 0.5*(magLe + nbrMagLe);
if (magLe < ROOTVSMALL && nbrMagLe < ROOTVSMALL) if (magLe < ROOTVSMALL && nbrMagLe < ROOTVSMALL)
{ {
@ -101,10 +100,10 @@ void Foam::cyclicFaPatch::calcTransforms()
// Check for error in edge matching // Check for error in edge matching
if (maxMatchError > matchTol_) if (maxMatchError > matchTol_)
{ {
label nbrEdgei = errorEdge + size()/2; label nbrEdgei = errorEdge + sizeby2;
scalar magLe = mag(half0Normals[errorEdge]); scalar magLe = mag(half0Normals[errorEdge]);
scalar nbrMagLe = mag(half1Normals[errorEdge]); scalar nbrMagLe = mag(half1Normals[errorEdge]);
scalar avLe = (magLe + nbrMagLe)/2.0; scalar avLe = 0.5*(magLe + nbrMagLe);
FatalErrorInFunction FatalErrorInFunction
<< "edge " << errorEdge << "edge " << errorEdge
@ -162,7 +161,7 @@ void Foam::cyclicFaPatch::makeWeights(scalarField& w) const
for (label edgei = 0; edgei < sizeby2; ++edgei) for (label edgei = 0; edgei < sizeby2; ++edgei)
{ {
scalar avL = (magL[edgei] + magL[edgei + sizeby2])/2.0; scalar avL = 0.5*(magL[edgei] + magL[edgei + sizeby2]);
if if
( (
@ -191,7 +190,7 @@ void Foam::cyclicFaPatch::makeWeights(scalarField& w) const
// Check for error in matching // Check for error in matching
if (maxMatchError > matchTol_) if (maxMatchError > matchTol_)
{ {
scalar avL = (magL[errorEdge] + magL[errorEdge + sizeby2])/2.0; scalar avL = 0.5*(magL[errorEdge] + magL[errorEdge + sizeby2]);
FatalErrorInFunction FatalErrorInFunction
<< "edge " << errorEdge << " and " << errorEdge + sizeby2 << "edge " << errorEdge << " and " << errorEdge + sizeby2
@ -313,7 +312,7 @@ Foam::tmp<Foam::labelField> Foam::cyclicFaPatch::transfer
const label sizeby2 = this->size()/2; const label sizeby2 = this->size()/2;
for (label edgei=0; edgei<sizeby2; ++edgei) for (label edgei=0; edgei < sizeby2; ++edgei)
{ {
pnf[edgei] = interfaceData[edgei + sizeby2]; pnf[edgei] = interfaceData[edgei + sizeby2];
pnf[edgei + sizeby2] = interfaceData[edgei]; pnf[edgei + sizeby2] = interfaceData[edgei];
@ -345,7 +344,7 @@ Foam::tmp<Foam::labelField> Foam::cyclicFaPatch::internalFieldTransfer
const label sizeby2 = this->size()/2; const label sizeby2 = this->size()/2;
for (label edgei=0; edgei<sizeby2; ++edgei) for (label edgei=0; edgei < sizeby2; ++edgei)
{ {
pnf[edgei] = iF[edgeCells[edgei + sizeby2]]; pnf[edgei] = iF[edgeCells[edgei + sizeby2]];
pnf[edgei + sizeby2] = iF[edgeCells[edgei]]; pnf[edgei + sizeby2] = iF[edgeCells[edgei]];

View File

@ -80,12 +80,7 @@ void Foam::processorFaPatch::makeNonGlobalPatchPoints() const
|| !boundaryMesh().mesh()().globalData().nGlobalPoints() || !boundaryMesh().mesh()().globalData().nGlobalPoints()
) )
{ {
nonGlobalPatchPointsPtr_ = new labelList(nPoints()); nonGlobalPatchPointsPtr_ = new labelList(identity(nPoints()));
labelList& ngpp = *nonGlobalPatchPointsPtr_;
forAll(ngpp, i)
{
ngpp[i] = i;
}
} }
else else
{ {

View File

@ -6,7 +6,7 @@
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
Copyright (C) 2016-2017 Wikki Ltd Copyright (C) 2016-2017 Wikki Ltd
Copyright (C) 2019 OpenCFD Ltd. Copyright (C) 2019-2021 OpenCFD Ltd.
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
This file is part of OpenFOAM. This file is part of OpenFOAM.
@ -177,19 +177,19 @@ public:
// Member Functions // Member Functions
//- Return interface size //- Return interface size
virtual label interfaceSize() const virtual label interfaceSize() const noexcept
{ {
return size(); return size();
} }
//- Return processor number //- Return processor number
int myProcNo() const int myProcNo() const noexcept
{ {
return myProcNo_; return myProcNo_;
} }
//- Return neighbour processor number //- Return neighbour processor number
int neighbProcNo() const int neighbProcNo() const noexcept
{ {
return neighbProcNo_; return neighbProcNo_;
} }
@ -197,18 +197,11 @@ public:
//- Return true if running parallel //- Return true if running parallel
virtual bool coupled() const virtual bool coupled() const
{ {
if (Pstream::parRun()) return Pstream::parRun();
{
return true;
}
else
{
return false;
}
} }
//- Is this the master side? //- Is this the master side?
virtual bool master() const virtual bool master() const noexcept
{ {
return (myProcNo_ < neighbProcNo_); return (myProcNo_ < neighbProcNo_);
} }

View File

@ -210,88 +210,6 @@ const Foam::labelList& Foam::faPatch::pointLabels() const
} }
void Foam::faPatch::calcPointLabels() const
{
SLList<label> labels;
const edgeList::subList edges = patchSlice(boundaryMesh().mesh().edges());
forAll(edges, edgeI)
{
bool existStart = false;
bool existEnd = false;
forAllIters(labels, iter)
{
if (*iter == edges[edgeI].start())
{
existStart = true;
}
if (*iter == edges[edgeI].end())
{
existEnd = true;
}
}
if (!existStart)
{
labels.append(edges[edgeI].start());
}
if (!existEnd)
{
labels.append(edges[edgeI].end());
}
}
pointLabelsPtr_ = new labelList(labels);
}
void Foam::faPatch::calcPointEdges() const
{
const labelList& points = pointLabels();
const edgeList::subList e = patchSlice(boundaryMesh().mesh().edges());
// set up storage for pointEdges
List<SLList<label>> pointEdgs(points.size());
forAll(e, edgeI)
{
const edge& curPoints = e[edgeI];
forAll(curPoints, pointI)
{
const label localPointIndex = points.find(curPoints[pointI]);
pointEdgs[localPointIndex].append(edgeI);
}
}
// sort out the list
pointEdgesPtr_ = new labelListList(pointEdgs.size());
labelListList& pEdges = *pointEdgesPtr_;
forAll(pointEdgs, pointI)
{
pEdges[pointI].setSize(pointEdgs[pointI].size());
label i = 0;
for
(
SLList<label>::iterator curEdgesIter = pointEdgs[pointI].begin();
curEdgesIter != pointEdgs[pointI].end();
++curEdgesIter, ++i
)
{
pEdges[pointI][i] = curEdgesIter();
}
}
}
const Foam::labelListList& Foam::faPatch::pointEdges() const const Foam::labelListList& Foam::faPatch::pointEdges() const
{ {
if (!pointEdgesPtr_) if (!pointEdgesPtr_)
@ -303,6 +221,94 @@ const Foam::labelListList& Foam::faPatch::pointEdges() const
} }
void Foam::faPatch::calcPointLabels() const
{
const edgeList::subList edges = patchSlice(boundaryMesh().mesh().edges());
// Walk boundary edges.
// The edge orientation corresponds to the face orientation
// (outwards normal).
// Note: could combine this with calcPointEdges for more efficiency
// Map<label> markedPoints(4*edges.size());
labelHashSet markedPoints(4*edges.size());
DynamicList<label> dynEdgePoints(2*edges.size());
for (const edge& e : edges)
{
// if (markedPoints.insert(e.first(), markedPoints.size()))
if (markedPoints.insert(e.first()))
{
dynEdgePoints.append(e.first());
}
// if (markedPoints.insert(e.second(), markedPoints.size()))
if (markedPoints.insert(e.second()))
{
dynEdgePoints.append(e.second());
}
}
// Transfer to plain list (reuse storage)
pointLabelsPtr_ = new labelList(std::move(dynEdgePoints));
/// const labelList& edgePoints = *pointLabelsPtr_;
///
/// // Cannot use invertManyToMany - we have non-local edge numbering
///
/// // Intermediate storage for pointEdges.
/// // Points on the boundary will normally connect 1 or 2 edges only.
/// List<DynamicList<label,2>> dynPointEdges(edgePoints.size());
///
/// forAll(edges, edgei)
/// {
/// const edge& e = edges[edgei];
///
/// dynPointEdges[markedPoints[e.first()]].append(edgei);
/// dynPointEdges[markedPoints[e.second()]].append(edgei);
/// }
///
/// // Flatten to regular list
/// pointEdgesPtr_ = new labelListList(edgePoints.size());
/// auto& pEdges = *pointEdgesPtr_;
///
/// forAll(pEdges, pointi)
/// {
/// pEdges[pointi] = std::move(dynPointEdges[pointi]);
/// }
}
void Foam::faPatch::calcPointEdges() const
{
const edgeList::subList edges = patchSlice(boundaryMesh().mesh().edges());
const labelList& edgePoints = pointLabels();
// Cannot use invertManyToMany - we have non-local edge numbering
// Intermediate storage for pointEdges.
// Points on the boundary will normally connect 1 or 2 edges only.
List<DynamicList<label,2>> dynPointEdges(edgePoints.size());
forAll(edges, edgei)
{
const edge& e = edges[edgei];
dynPointEdges[edgePoints.find(e.first())].append(edgei);
dynPointEdges[edgePoints.find(e.second())].append(edgei);
}
// Flatten to regular list
pointEdgesPtr_ = new labelListList(edgePoints.size());
auto& pEdges = *pointEdgesPtr_;
forAll(pEdges, pointi)
{
pEdges[pointi] = std::move(dynPointEdges[pointi]);
}
}
Foam::tmp<Foam::vectorField> Foam::faPatch::ngbPolyPatchFaceNormals() const Foam::tmp<Foam::vectorField> Foam::faPatch::ngbPolyPatchFaceNormals() const
{ {
if (nbrPolyPatchId_ < 0) if (nbrPolyPatchId_ < 0)
@ -383,7 +389,7 @@ const Foam::scalarField& Foam::faPatch::magEdgeLengths() const
Foam::tmp<Foam::vectorField> Foam::faPatch::edgeNormals() const Foam::tmp<Foam::vectorField> Foam::faPatch::edgeNormals() const
{ {
tmp<vectorField> tedgeNorm(edgeLengths()); auto tedgeNorm = tmp<vectorField>::New(edgeLengths());
for (vector& n : tedgeNorm.ref()) for (vector& n : tedgeNorm.ref())
{ {

View File

@ -80,17 +80,14 @@ class faPatch
//- Reference to boundary mesh //- Reference to boundary mesh
const faBoundaryMesh& boundaryMesh_; const faBoundaryMesh& boundaryMesh_;
//- Demand-driven: edge-face addressing
mutable labelList::subList* edgeFacesPtr_;
// Demand-driven private data //- Demand-driven: local points labels
mutable labelList* pointLabelsPtr_;
//- Edge-face addressing //- Demand-driven: point-edge addressing
mutable labelList::subList* edgeFacesPtr_; mutable labelListList* pointEdgesPtr_;
//- Local points labels
mutable labelList* pointLabelsPtr_;
//- Point-edge addressing
mutable labelListList* pointEdgesPtr_;
// Private Member Functions // Private Member Functions