ENH: paraFoam: catch read errors. Fixes #798.

This commit is contained in:
mattijs
2018-04-12 15:12:04 +01:00
parent cb14a2020b
commit d2b1b1cdc0
4 changed files with 271 additions and 149 deletions

View File

@ -2,7 +2,7 @@
========= | ========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | \\ / O peration |
\\ / A nd | Copyright (C) 2017 OpenCFD Ltd. \\ / A nd | Copyright (C) 2017-2018 OpenCFD Ltd.
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
@ -276,11 +276,33 @@ void Foam::vtkPVFoam::convertVolFields
if (ioobj.headerClassName() == FieldType::typeName) if (ioobj.headerClassName() == FieldType::typeName)
{ {
// Load field // Throw FatalError, FatalIOError as exceptions
FieldType fld(ioobj, mesh); const bool throwingError = FatalError.throwExceptions();
const bool throwingIOerr = FatalIOError.throwExceptions();
// Convert try
convertVolField(patchInterpList, fld); {
// Load field
FieldType fld(ioobj, mesh);
// Convert
convertVolField(patchInterpList, fld);
}
catch (Foam::IOerror& ioErr)
{
ioErr.write(Warning, false);
Info << nl << endl;
}
catch (Foam::error& err)
{
// Bit of trickery to get the original message
err.write(Warning, false);
Info << nl << endl;
}
// Restore previous exception throwing state
FatalError.throwExceptions(throwingError);
FatalIOError.throwExceptions(throwingIOerr);
} }
} }
} }
@ -307,40 +329,62 @@ void Foam::vtkPVFoam::convertDimFields
continue; continue;
} }
// Load field // Throw FatalError, FatalIOError as exceptions
FieldType dimFld(ioobj, mesh); const bool throwingError = FatalError.throwExceptions();
const bool throwingIOerr = FatalIOError.throwExceptions();
// Construct volField with zero-gradient patch fields try
IOobject io(dimFld);
io.readOpt() = IOobject::NO_READ;
PtrList<fvPatchField<Type>> patchFields(mesh.boundary().size());
forAll(patchFields, patchI)
{ {
patchFields.set // Load field
( FieldType dimFld(ioobj, mesh);
patchI,
fvPatchField<Type>::New // Construct volField with zero-gradient patch fields
IOobject io(dimFld);
io.readOpt() = IOobject::NO_READ;
PtrList<fvPatchField<Type>> patchFields(mesh.boundary().size());
forAll(patchFields, patchI)
{
patchFields.set
( (
zeroGradientFvPatchField<Type>::typeName, patchI,
mesh.boundary()[patchI], fvPatchField<Type>::New
dimFld (
) zeroGradientFvPatchField<Type>::typeName,
mesh.boundary()[patchI],
dimFld
)
);
}
VolFieldType volFld
(
io,
dimFld.mesh(),
dimFld.dimensions(),
dimFld,
patchFields
); );
volFld.correctBoundaryConditions();
convertVolField(patchInterpList, volFld);
}
catch (Foam::IOerror& ioErr)
{
ioErr.write(Warning, false);
Info << nl << endl;
}
catch (Foam::error& err)
{
// Bit of trickery to get the original message
err.write(Warning, false);
Info << nl << endl;
} }
VolFieldType volFld // Restore previous exception throwing state
( FatalError.throwExceptions(throwingError);
io, FatalIOError.throwExceptions(throwingIOerr);
dimFld.mesh(),
dimFld.dimensions(),
dimFld,
patchFields
);
volFld.correctBoundaryConditions();
convertVolField(patchInterpList, volFld);
} }
} }
@ -421,31 +465,54 @@ void Foam::vtkPVFoam::convertAreaFields
if (ioobj.headerClassName() == FieldType::typeName) if (ioobj.headerClassName() == FieldType::typeName)
{ {
// Load field // Throw FatalError, FatalIOError as exceptions
FieldType fld(ioobj, mesh); const bool throwingError = FatalError.throwExceptions();
const bool throwingIOerr = FatalIOError.throwExceptions();
// Convert try
for (const auto partId : partIds)
{ {
const auto& longName = selectedPartIds_[partId]; // Load field
FieldType fld(ioobj, mesh);
auto iter = cachedVtp_.find(longName); // Convert
if (!iter.found() || !iter.object().dataset) for (const auto partId : partIds)
{ {
// Should not happen, but for safety require a vtk geometry const auto& longName = selectedPartIds_[partId];
continue;
auto iter = cachedVtp_.find(longName);
if (!iter.found() || !iter.object().dataset)
{
// Should not happen, but for safety require a vtk
// geometry
continue;
}
foamVtpData& vtpData = iter.object();
auto dataset = vtpData.dataset;
vtkSmartPointer<vtkFloatArray> cdata = convertFieldToVTK
(
fld.name(),
fld
);
dataset->GetCellData()->AddArray(cdata);
} }
foamVtpData& vtpData = iter.object();
auto dataset = vtpData.dataset;
vtkSmartPointer<vtkFloatArray> cdata = convertFieldToVTK
(
fld.name(),
fld
);
dataset->GetCellData()->AddArray(cdata);
} }
catch (Foam::IOerror& ioErr)
{
ioErr.write(Warning, false);
Info << nl << endl;
}
catch (Foam::error& err)
{
// Bit of trickery to get the original message
err.write(Warning, false);
Info << nl << endl;
}
// Restore previous exception throwing state
FatalError.throwExceptions(throwingError);
FatalIOError.throwExceptions(throwingIOerr);
} }
} }
} }
@ -483,91 +550,113 @@ void Foam::vtkPVFoam::convertPointFields
Info<< "convertPointFields : " << fieldName << nl; Info<< "convertPointFields : " << fieldName << nl;
} }
FieldType pfld(ioobj, pMesh); // Throw FatalError, FatalIOError as exceptions
const bool throwingError = FatalError.throwExceptions();
const bool throwingIOerr = FatalIOError.throwExceptions();
convertPointFieldBlock(pfld, rangeVolume_); // internalMesh try
convertPointFieldBlock(pfld, rangeCellZones_); // cellZones
convertPointFieldBlock(pfld, rangeCellSets_); // cellSets
// Patches - currently skip field conversion for groups
for
(
const auto partId
: rangePatches_.intersection(selectedPartIds_)
)
{ {
const auto& longName = selectedPartIds_[partId]; FieldType pfld(ioobj, pMesh);
auto iter = cachedVtp_.find(longName); convertPointFieldBlock(pfld, rangeVolume_); // internalMesh
if (!iter.found() || !iter.object().dataset) convertPointFieldBlock(pfld, rangeCellZones_); // cellZones
{ convertPointFieldBlock(pfld, rangeCellSets_); // cellSets
// Should not happen, but for safety require a vtk geometry
continue;
}
foamVtpData& vtpData = iter.object(); // Patches - currently skip field conversion for groups
auto dataset = vtpData.dataset; for
const labelList& patchIds = vtpData.additionalIds();
if (patchIds.size() != 1)
{
continue;
}
const label patchId = patchIds[0];
vtkSmartPointer<vtkFloatArray> pdata = convertFieldToVTK
( (
fieldName, const auto partId
pfld.boundaryField()[patchId].patchInternalField()() : rangePatches_.intersection(selectedPartIds_)
); )
dataset->GetPointData()->AddArray(pdata);
}
// Face Zones
for
(
const auto partId
: rangeFaceZones_.intersection(selectedPartIds_)
)
{
const auto& longName = selectedPartIds_[partId];
const word zoneName = getFoamName(longName);
auto iter = cachedVtp_.find(longName);
if (!iter.found() || !iter.object().dataset)
{ {
// Should not happen, but for safety require a vtk geometry const auto& longName = selectedPartIds_[partId];
continue;
}
foamVtpData& vtpData = iter.object(); auto iter = cachedVtp_.find(longName);
auto dataset = vtpData.dataset; if (!iter.found() || !iter.object().dataset)
{
// Should not happen, but for safety require a vtk geometry
continue;
}
const label zoneId = mesh.faceZones().findZoneID(zoneName); foamVtpData& vtpData = iter.object();
auto dataset = vtpData.dataset;
if (zoneId < 0) const labelList& patchIds = vtpData.additionalIds();
{ if (patchIds.size() != 1)
continue; {
} continue;
}
// Extract the field on the zone const label patchId = patchIds[0];
Field<Type> znfld
(
pfld.primitiveField(),
mesh.faceZones()[zoneId]().meshPoints()
);
vtkSmartPointer<vtkFloatArray> pdata = vtkSmartPointer<vtkFloatArray> pdata = convertFieldToVTK
convertFieldToVTK
( (
fieldName, fieldName,
znfld pfld.boundaryField()[patchId].patchInternalField()()
); );
dataset->GetPointData()->AddArray(pdata); dataset->GetPointData()->AddArray(pdata);
}
// Face Zones
for
(
const auto partId
: rangeFaceZones_.intersection(selectedPartIds_)
)
{
const auto& longName = selectedPartIds_[partId];
const word zoneName = getFoamName(longName);
auto iter = cachedVtp_.find(longName);
if (!iter.found() || !iter.object().dataset)
{
// Should not happen, but for safety require a vtk geometry
continue;
}
foamVtpData& vtpData = iter.object();
auto dataset = vtpData.dataset;
const label zoneId = mesh.faceZones().findZoneID(zoneName);
if (zoneId < 0)
{
continue;
}
// Extract the field on the zone
Field<Type> znfld
(
pfld.primitiveField(),
mesh.faceZones()[zoneId]().meshPoints()
);
vtkSmartPointer<vtkFloatArray> pdata =
convertFieldToVTK
(
fieldName,
znfld
);
dataset->GetPointData()->AddArray(pdata);
}
} }
catch (Foam::IOerror& ioErr)
{
ioErr.write(Warning, false);
Info << nl << endl;
}
catch (Foam::error& err)
{
// Bit of trickery to get the original message
err.write(Warning, false);
Info << nl << endl;
}
// Restore previous exception throwing state
FatalError.throwExceptions(throwingError);
FatalIOError.throwExceptions(throwingIOerr);
} }
} }
@ -735,18 +824,40 @@ void Foam::vtkPVFoam::convertLagrangianFields
if (ioobj.headerClassName() == IOField<Type>::typeName) if (ioobj.headerClassName() == IOField<Type>::typeName)
{ {
IOField<Type> fld(ioobj); // Throw FatalError, FatalIOError as exceptions
const bool throwingError = FatalError.throwExceptions();
const bool throwingIOerr = FatalIOError.throwExceptions();
vtkSmartPointer<vtkFloatArray> data = try
convertFieldToVTK {
( IOField<Type> fld(ioobj);
ioobj.name(),
fld
);
// Provide identical data as cell and as point data vtkSmartPointer<vtkFloatArray> data =
vtkmesh->GetCellData()->AddArray(data); convertFieldToVTK
vtkmesh->GetPointData()->AddArray(data); (
ioobj.name(),
fld
);
// Provide identical data as cell and as point data
vtkmesh->GetCellData()->AddArray(data);
vtkmesh->GetPointData()->AddArray(data);
}
catch (Foam::IOerror& ioErr)
{
ioErr.write(Warning, false);
Info << nl << endl;
}
catch (Foam::error& err)
{
// Bit of trickery to get the original message
err.write(Warning, false);
Info << nl << endl;
}
// Restore previous exception throwing state
FatalError.throwExceptions(throwingError);
FatalIOError.throwExceptions(throwingIOerr);
} }
} }
} }

View File

@ -46,7 +46,6 @@ wmake $targetType genericPatchFields
conversion/Allwmake $targetType $* conversion/Allwmake $targetType $*
wmake $targetType mesh/extrudeModel wmake $targetType mesh/extrudeModel
wmake $targetType dynamicMesh wmake $targetType dynamicMesh
wmake $targetType sampling
wmake $targetType dynamicFvMesh wmake $targetType dynamicFvMesh
wmake $targetType sampling wmake $targetType sampling
wmake $targetType topoChangerFvMesh wmake $targetType topoChangerFvMesh

View File

@ -260,34 +260,43 @@ void Foam::IOerror::abort()
} }
Foam::Ostream& Foam::operator<<(Ostream& os, const IOerror& err) void Foam::IOerror::write(Ostream& os, const bool includeTitle) const
{ {
if (!os.bad()) if (!os.bad())
{ {
os << nl os << nl;
<< err.title().c_str() << nl if (includeTitle)
<< err.message().c_str() << nl << endl;
os << "file: " << err.ioFileName().c_str();
if (err.ioStartLineNumber() >= 0 && err.ioEndLineNumber() >= 0)
{ {
os << " from line " << err.ioStartLineNumber() os << title().c_str() << nl;
<< " to line " << err.ioEndLineNumber() << '.';
} }
else if (err.ioStartLineNumber() >= 0) os << message().c_str() << nl << endl;
os << "file: " << ioFileName().c_str();
if (ioStartLineNumber() >= 0 && ioEndLineNumber() >= 0)
{ {
os << " at line " << err.ioStartLineNumber() << '.'; os << " from line " << ioStartLineNumber()
<< " to line " << ioEndLineNumber() << '.';
}
else if (ioStartLineNumber() >= 0)
{
os << " at line " << ioStartLineNumber() << '.';
} }
if (IOerror::level >= 2 && err.sourceFileLineNumber()) if (IOerror::level >= 2 && sourceFileLineNumber())
{ {
os << nl << nl os << nl << nl
<< " From function " << err.functionName().c_str() << endl << " From function " << functionName().c_str() << endl
<< " in file " << err.sourceFileName().c_str() << " in file " << sourceFileName().c_str()
<< " at line " << err.sourceFileLineNumber() << '.'; << " at line " << sourceFileLineNumber() << '.';
} }
} }
}
Foam::Ostream& Foam::operator<<(Ostream& os, const IOerror& err)
{
err.write(os);
return os; return os;
} }

View File

@ -3,7 +3,7 @@
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | \\ / O peration |
\\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation \\ / A nd | Copyright (C) 2011-2015 OpenFOAM Foundation
\\/ M anipulation | Copyright (C) 2015-2017 OpenCFD Ltd. \\/ M anipulation | Copyright (C) 2015-2018 OpenCFD Ltd.
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
This file is part of OpenFOAM. This file is part of OpenFOAM.
@ -320,6 +320,9 @@ public:
//- Abort : used to stop code for fatal errors //- Abort : used to stop code for fatal errors
void abort(); void abort();
//- Print error message
void write(Ostream& os, const bool includeTitle = true) const;
// Ostream operator // Ostream operator