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
@ -275,6 +275,12 @@ void Foam::vtkPVFoam::convertVolFields
const auto& ioobj = *(iter.object()); const auto& ioobj = *(iter.object());
if (ioobj.headerClassName() == FieldType::typeName) if (ioobj.headerClassName() == FieldType::typeName)
{
// Throw FatalError, FatalIOError as exceptions
const bool throwingError = FatalError.throwExceptions();
const bool throwingIOerr = FatalIOError.throwExceptions();
try
{ {
// Load field // Load field
FieldType fld(ioobj, mesh); FieldType fld(ioobj, mesh);
@ -282,6 +288,22 @@ void Foam::vtkPVFoam::convertVolFields
// Convert // Convert
convertVolField(patchInterpList, fld); 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,6 +329,12 @@ void Foam::vtkPVFoam::convertDimFields
continue; continue;
} }
// Throw FatalError, FatalIOError as exceptions
const bool throwingError = FatalError.throwExceptions();
const bool throwingIOerr = FatalIOError.throwExceptions();
try
{
// Load field // Load field
FieldType dimFld(ioobj, mesh); FieldType dimFld(ioobj, mesh);
@ -342,6 +370,22 @@ void Foam::vtkPVFoam::convertDimFields
convertVolField(patchInterpList, volFld); 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;
}
// Restore previous exception throwing state
FatalError.throwExceptions(throwingError);
FatalIOError.throwExceptions(throwingIOerr);
}
} }
@ -420,6 +464,12 @@ void Foam::vtkPVFoam::convertAreaFields
const auto& ioobj = *(iter.object()); const auto& ioobj = *(iter.object());
if (ioobj.headerClassName() == FieldType::typeName) if (ioobj.headerClassName() == FieldType::typeName)
{
// Throw FatalError, FatalIOError as exceptions
const bool throwingError = FatalError.throwExceptions();
const bool throwingIOerr = FatalIOError.throwExceptions();
try
{ {
// Load field // Load field
FieldType fld(ioobj, mesh); FieldType fld(ioobj, mesh);
@ -432,7 +482,8 @@ void Foam::vtkPVFoam::convertAreaFields
auto iter = cachedVtp_.find(longName); auto iter = cachedVtp_.find(longName);
if (!iter.found() || !iter.object().dataset) if (!iter.found() || !iter.object().dataset)
{ {
// Should not happen, but for safety require a vtk geometry // Should not happen, but for safety require a vtk
// geometry
continue; continue;
} }
@ -447,6 +498,22 @@ void Foam::vtkPVFoam::convertAreaFields
dataset->GetCellData()->AddArray(cdata); 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,6 +550,12 @@ void Foam::vtkPVFoam::convertPointFields
Info<< "convertPointFields : " << fieldName << nl; Info<< "convertPointFields : " << fieldName << nl;
} }
// Throw FatalError, FatalIOError as exceptions
const bool throwingError = FatalError.throwExceptions();
const bool throwingIOerr = FatalIOError.throwExceptions();
try
{
FieldType pfld(ioobj, pMesh); FieldType pfld(ioobj, pMesh);
convertPointFieldBlock(pfld, rangeVolume_); // internalMesh convertPointFieldBlock(pfld, rangeVolume_); // internalMesh
@ -569,6 +642,22 @@ void Foam::vtkPVFoam::convertPointFields
dataset->GetPointData()->AddArray(pdata); 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);
}
} }
@ -734,6 +823,12 @@ void Foam::vtkPVFoam::convertLagrangianFields
const auto& ioobj = *(iter.object()); const auto& ioobj = *(iter.object());
if (ioobj.headerClassName() == IOField<Type>::typeName) if (ioobj.headerClassName() == IOField<Type>::typeName)
{
// Throw FatalError, FatalIOError as exceptions
const bool throwingError = FatalError.throwExceptions();
const bool throwingIOerr = FatalIOError.throwExceptions();
try
{ {
IOField<Type> fld(ioobj); IOField<Type> fld(ioobj);
@ -748,6 +843,22 @@ void Foam::vtkPVFoam::convertLagrangianFields
vtkmesh->GetCellData()->AddArray(data); vtkmesh->GetCellData()->AddArray(data);
vtkmesh->GetPointData()->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