decomposePar, reconstructPar: Do all regions simultaneously

DecomposePar and reconstructPar now interleave the processing of
multiple regions. This means that with the -allRegions option, the
earlier times are completed in their entirety before later times are
considered. It also lets regions to access each other during
decomposition and reconstruction, which will be important for
non-conformal region interfaces.

To aid interpretation of the log, region prefixing is now used by both
utilities in the same way as is done by foamMultiRun.

DecomposePar has been overhauled so that it matches reconstructPar much
more closely, both in terms of output and of iteration sequence. All
meshes and addressing are loaded simultaneously and each field is
considered in turn. Previously, all the fields were loaded, and each
process and addressing set was considered in turn. This new strategy
optimises memory usage for cases with lots of fields.
This commit is contained in:
Will Bainbridge
2023-07-20 11:42:18 +01:00
parent 58e15f296c
commit 71ccf51ba5
43 changed files with 4246 additions and 4028 deletions

View File

@ -1,6 +1,6 @@
reconstructPar.C
fvFieldReconstructor.C
pointFieldReconstructor.C
reconstructLagrangianPositions.C
lagrangianFieldReconstructor.C
EXE = $(FOAM_APPBIN)/reconstructPar

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -66,8 +66,7 @@ Foam::fvFieldReconstructor::fvFieldReconstructor
procMeshes_(procMeshes),
faceProcAddressing_(faceProcAddressing),
cellProcAddressing_(cellProcAddressing),
faceProcAddressingBf_(faceProcAddressingBf),
nReconstructed_(0)
faceProcAddressingBf_(faceProcAddressingBf)
{
forAll(procMeshes_, proci)
{
@ -93,4 +92,26 @@ Foam::fvFieldReconstructor::fvFieldReconstructor
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
bool Foam::fvFieldReconstructor::reconstructs
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
)
{
bool result = false;
#define DO_FV_FIELDS_TYPE(Type, nullArg) \
result = result \
|| reconstructs<VolField<Type>::Internal>(objects, selectedFields) \
|| reconstructs<VolField<Type>>(objects, selectedFields) \
|| reconstructs<SurfaceField<Type>>(objects, selectedFields);
FOR_ALL_FIELD_TYPES(DO_FV_FIELDS_TYPE)
#undef DO_FV_FIELDS_TYPE
return result;
}
// ************************************************************************* //

View File

@ -29,7 +29,7 @@ Description
SourceFiles
fvFieldReconstructor.C
fvFieldReconstructorReconstructFields.C
fvFieldReconstructorTemplates.C
\*---------------------------------------------------------------------------*/
@ -71,12 +71,17 @@ class fvFieldReconstructor
//- Boundary field of face addressing
const PtrList<surfaceLabelField::Boundary>& faceProcAddressingBf_;
//- Number of fields reconstructed
label nReconstructed_;
// Private Member Functions
//- Return whether anything in the object list gets reconstructed
template<class FieldType>
static bool reconstructs
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
);
//- Convert a processor patch to the corresponding complete patch index
label completePatchID(const label proci, const label procPatchi) const;
@ -90,6 +95,21 @@ class fvFieldReconstructor
const bool isFlux
);
//- Read and reconstruct a volume internal field
template<class Type>
tmp<DimensionedField<Type, volMesh>>
reconstructVolInternalField(const IOobject& fieldIoObject) const;
//- Read and reconstruct a volume field
template<class Type>
tmp<VolField<Type>>
reconstructVolField(const IOobject& fieldIoObject) const;
//- Read and reconstruct a surface field
template<class Type>
tmp<SurfaceField<Type>>
reconstructFvSurfaceField(const IOobject& fieldIoObject) const;
public:
@ -111,59 +131,16 @@ public:
// Member Functions
//- Return number of fields reconstructed
label nReconstructed() const
{
return nReconstructed_;
}
//- Reconstruct volume internal field
template<class Type>
tmp<DimensionedField<Type, volMesh>>
reconstructFvVolumeInternalField
//- Return whether anything in the object list gets reconstructed
static bool reconstructs
(
const IOobject& fieldIoObject,
const PtrList<DimensionedField<Type, volMesh>>& procFields
) const;
//- Read and reconstruct volume internal field
template<class Type>
tmp<DimensionedField<Type, volMesh>>
reconstructFvVolumeInternalField(const IOobject& fieldIoObject) const;
//- Reconstruct volume field
template<class Type>
tmp<VolField<Type>>
reconstructFvVolumeField
(
const IOobject& fieldIoObject,
const PtrList<VolField<Type>>&
) const;
//- Read and reconstruct volume field
template<class Type>
tmp<VolField<Type>>
reconstructFvVolumeField(const IOobject& fieldIoObject) const;
//- Reconstruct surface field
template<class Type>
tmp<SurfaceField<Type>>
reconstructFvSurfaceField
(
const IOobject& fieldIoObject,
const PtrList<SurfaceField<Type>>&
) const;
//- Read and reconstruct surface field
template<class Type>
tmp<SurfaceField<Type>>
reconstructFvSurfaceField(const IOobject& fieldIoObject) const;
const IOobjectList& objects,
const HashSet<word>& selectedFields
);
//- Read, reconstruct and write all/selected volume internal fields
template<class Type>
void reconstructFvVolumeInternalFields
void reconstructVolInternalFields
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
@ -171,7 +148,7 @@ public:
//- Read, reconstruct and write all/selected volume fields
template<class Type>
void reconstructFvVolumeFields
void reconstructVolFields
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
@ -200,7 +177,7 @@ public:
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
#include "fvFieldReconstructorReconstructFields.C"
#include "fvFieldReconstructorTemplates.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

View File

@ -34,7 +34,33 @@ License
#include "reverseFvPatchFieldMapper.H"
#include "stringOps.H"
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
template<class FieldType>
bool Foam::fvFieldReconstructor::reconstructs
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
)
{
IOobjectList fields = objects.lookupClass(FieldType::typeName);
if (fields.size() && selectedFields.empty())
{
return true;
}
forAllConstIter(IOobjectList, fields, fieldIter)
{
if (selectedFields.found(fieldIter()->name()))
{
return true;
}
}
return false;
}
template<class Type>
void Foam::fvFieldReconstructor::rmapFaceToFace
@ -53,54 +79,15 @@ void Foam::fvFieldReconstructor::rmapFaceToFace
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class Type>
Foam::tmp<Foam::DimensionedField<Type, Foam::volMesh>>
Foam::fvFieldReconstructor::reconstructFvVolumeInternalField
(
const IOobject& fieldIoObject,
const PtrList<DimensionedField<Type, volMesh>>& procFields
) const
{
// Create the internalField
Field<Type> internalField(completeMesh_.nCells());
forAll(procMeshes_, proci)
{
const DimensionedField<Type, volMesh>& procField = procFields[proci];
// Set the cell values in the reconstructed field
internalField.rmap
(
procField.field(),
cellProcAddressing_[proci]
);
}
return tmp<DimensionedField<Type, volMesh>>
(
new DimensionedField<Type, volMesh>
(
fieldIoObject,
completeMesh_,
procFields[0].dimensions(),
internalField
)
);
}
template<class Type>
Foam::tmp<Foam::DimensionedField<Type, Foam::volMesh>>
Foam::fvFieldReconstructor::reconstructFvVolumeInternalField
Foam::fvFieldReconstructor::reconstructVolInternalField
(
const IOobject& fieldIoObject
) const
{
PtrList<DimensionedField<Type, volMesh>>
procFields(procMeshes_.size());
// Read the field for all the processors
PtrList<DimensionedField<Type, volMesh>> procFields(procMeshes_.size());
forAll(procMeshes_, proci)
{
procFields.set
@ -122,30 +109,72 @@ Foam::fvFieldReconstructor::reconstructFvVolumeInternalField
);
}
return reconstructFvVolumeInternalField
(
IOobject
// Create the internalField
Field<Type> internalField(completeMesh_.nCells());
forAll(procMeshes_, proci)
{
const DimensionedField<Type, volMesh>& procField = procFields[proci];
// Set the cell values in the reconstructed field
internalField.rmap
(
fieldIoObject.name(),
completeMesh_.time().name(),
procField.field(),
cellProcAddressing_[proci]
);
}
return tmp<DimensionedField<Type, volMesh>>
(
new DimensionedField<Type, volMesh>
(
IOobject
(
fieldIoObject.name(),
completeMesh_.time().name(),
completeMesh_,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
completeMesh_,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
procFields
procFields[0].dimensions(),
internalField
)
);
}
template<class Type>
Foam::tmp<Foam::VolField<Type>>
Foam::fvFieldReconstructor::reconstructFvVolumeField
Foam::fvFieldReconstructor::reconstructVolField
(
const IOobject& fieldIoObject,
const PtrList<VolField<Type>>& procFields
const IOobject& fieldIoObject
) const
{
// Read the field for all the processors
PtrList<VolField<Type>> procFields(procMeshes_.size());
forAll(procMeshes_, proci)
{
procFields.set
(
proci,
new VolField<Type>
(
IOobject
(
fieldIoObject.name(),
procMeshes_[proci].time().name(),
procMeshes_[proci],
IOobject::MUST_READ,
IOobject::NO_WRITE,
false
),
procMeshes_[proci]
)
);
}
// Create the internalField
Field<Type> internalField(completeMesh_.nCells());
@ -251,7 +280,15 @@ Foam::fvFieldReconstructor::reconstructFvVolumeField
(
new VolField<Type>
(
fieldIoObject,
IOobject
(
fieldIoObject.name(),
completeMesh_.time().name(),
completeMesh_,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
completeMesh_,
procFields[0].dimensions(),
internalField,
@ -262,21 +299,20 @@ Foam::fvFieldReconstructor::reconstructFvVolumeField
template<class Type>
Foam::tmp<Foam::VolField<Type>>
Foam::fvFieldReconstructor::reconstructFvVolumeField
Foam::tmp<Foam::SurfaceField<Type>>
Foam::fvFieldReconstructor::reconstructFvSurfaceField
(
const IOobject& fieldIoObject
) const
{
PtrList<VolField<Type>>
procFields(procMeshes_.size());
// Read the field for all the processors
PtrList<SurfaceField<Type>> procFields(procMeshes_.size());
forAll(procMeshes_, proci)
{
procFields.set
(
proci,
new VolField<Type>
new SurfaceField<Type>
(
IOobject
(
@ -292,30 +328,6 @@ Foam::fvFieldReconstructor::reconstructFvVolumeField
);
}
return reconstructFvVolumeField
(
IOobject
(
fieldIoObject.name(),
completeMesh_.time().name(),
completeMesh_,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
procFields
);
}
template<class Type>
Foam::tmp<Foam::SurfaceField<Type>>
Foam::fvFieldReconstructor::reconstructFvSurfaceField
(
const IOobject& fieldIoObject,
const PtrList<SurfaceField<Type>>& procFields
) const
{
// Create the internalField
Field<Type> internalField(completeMesh_.nInternalFaces());
@ -420,7 +432,15 @@ Foam::fvFieldReconstructor::reconstructFvSurfaceField
(
new SurfaceField<Type>
(
fieldIoObject,
IOobject
(
fieldIoObject.name(),
completeMesh_.time().name(),
completeMesh_,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
completeMesh_,
procFields[0].dimensions(),
internalField,
@ -430,55 +450,10 @@ Foam::fvFieldReconstructor::reconstructFvSurfaceField
}
template<class Type>
Foam::tmp<Foam::SurfaceField<Type>>
Foam::fvFieldReconstructor::reconstructFvSurfaceField
(
const IOobject& fieldIoObject
) const
{
PtrList<SurfaceField<Type>>
procFields(procMeshes_.size());
forAll(procMeshes_, proci)
{
procFields.set
(
proci,
new SurfaceField<Type>
(
IOobject
(
fieldIoObject.name(),
procMeshes_[proci].time().name(),
procMeshes_[proci],
IOobject::MUST_READ,
IOobject::NO_WRITE,
false
),
procMeshes_[proci]
)
);
}
return reconstructFvSurfaceField
(
IOobject
(
fieldIoObject.name(),
completeMesh_.time().name(),
completeMesh_,
IOobject::NO_READ,
IOobject::NO_WRITE,
false
),
procFields
);
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class Type>
void Foam::fvFieldReconstructor::reconstructFvVolumeInternalFields
void Foam::fvFieldReconstructor::reconstructVolInternalFields
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
@ -490,7 +465,8 @@ void Foam::fvFieldReconstructor::reconstructFvVolumeInternalFields
if (fields.size())
{
Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
Info<< nl << " Reconstructing " << fieldClassName << "s"
<< nl << endl;
forAllConstIter(IOobjectList, fields, fieldIter)
{
@ -502,18 +478,15 @@ void Foam::fvFieldReconstructor::reconstructFvVolumeInternalFields
{
Info<< " " << fieldIter()->name() << endl;
reconstructFvVolumeInternalField<Type>(*fieldIter())().write();
nReconstructed_++;
reconstructVolInternalField<Type>(*fieldIter())().write();
}
}
Info<< endl;
}
}
template<class Type>
void Foam::fvFieldReconstructor::reconstructFvVolumeFields
void Foam::fvFieldReconstructor::reconstructVolFields
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
@ -526,7 +499,8 @@ void Foam::fvFieldReconstructor::reconstructFvVolumeFields
if (fields.size())
{
Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
Info<< nl << " Reconstructing " << fieldClassName << "s"
<< nl << endl;
forAllConstIter(IOobjectList, fields, fieldIter)
{
@ -538,12 +512,9 @@ void Foam::fvFieldReconstructor::reconstructFvVolumeFields
{
Info<< " " << fieldIter()->name() << endl;
reconstructFvVolumeField<Type>(*fieldIter())().write();
nReconstructed_++;
reconstructVolField<Type>(*fieldIter())().write();
}
}
Info<< endl;
}
}
@ -562,7 +533,8 @@ void Foam::fvFieldReconstructor::reconstructFvSurfaceFields
if (fields.size())
{
Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
Info<< nl << " Reconstructing " << fieldClassName << "s"
<< nl << endl;
forAllConstIter(IOobjectList, fields, fieldIter)
{
@ -575,11 +547,8 @@ void Foam::fvFieldReconstructor::reconstructFvSurfaceFields
Info<< " " << fieldIter()->name() << endl;
reconstructFvSurfaceField<Type>(*fieldIter())().write();
nReconstructed_++;
}
}
Info<< endl;
}
}

View File

@ -0,0 +1,124 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "lagrangianFieldReconstructor.H"
#include "passiveParticleCloud.H"
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
Foam::lagrangianFieldReconstructor::lagrangianFieldReconstructor
(
const fvMesh& completeMesh,
const PtrList<fvMesh>& procMeshes,
const labelListList& faceProcAddressing,
const labelListList& cellProcAddressing,
const word& cloudName
)
:
completeMesh_(completeMesh),
procMeshes_(procMeshes),
cloudName_(cloudName)
{
// Construct and empty cloud for the complete positions
passiveParticleCloud completePositions
(
completeMesh_,
cloudName_,
IDLList<passiveParticle>()
);
forAll(procMeshes_, proci)
{
// Read the processor positions
Cloud<passiveParticle> procPositions
(
procMeshes_[proci],
cloudName_,
false
);
// Combine the processor's positions into the complete cloud
forAllConstIter(Cloud<passiveParticle>, procPositions, iter)
{
const passiveParticle& p = iter();
const label completeCelli = cellProcAddressing[proci][p.cell()];
const label completeFacei =
mag(faceProcAddressing[proci][p.tetFace()]) - 1;
completePositions.append
(
new passiveParticle
(
completeMesh_,
p.coordinates(),
completeCelli,
completeFacei,
p.procTetPt
(
procMeshes_[proci],
completeMesh_,
completeCelli,
completeFacei
)
)
);
}
}
// Write
IOPosition<Cloud<passiveParticle>>(completePositions).write();
}
// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
Foam::lagrangianFieldReconstructor::~lagrangianFieldReconstructor()
{}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
bool Foam::lagrangianFieldReconstructor::reconstructs
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
)
{
bool result = false;
#define DO_LAGRANGIAN_FIELDS_TYPE(Type, nullArg) \
result = result \
|| reconstructs<IOField<Type>>(objects, selectedFields) \
|| reconstructs<IOField<Field<Type>>>(objects, selectedFields) \
|| reconstructs<CompactIOField<Field<Type>>>(objects, selectedFields);
DO_LAGRANGIAN_FIELDS_TYPE(label, )
FOR_ALL_FIELD_TYPES(DO_LAGRANGIAN_FIELDS_TYPE)
#undef DO_LAGRANGIAN_FIELDS_TYPE
return result;
}
// ************************************************************************* //

View File

@ -0,0 +1,165 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
Class
Foam::lagrangianFieldReconstructor
Description
Lagrangian field reconstructor.
SourceFiles
lagrangianFieldReconstructor.C
lagrangianFieldReconstructorTemplates.C
\*---------------------------------------------------------------------------*/
#ifndef lagrangianFieldReconstructor_H
#define lagrangianFieldReconstructor_H
#include "cloud.H"
#include "fvMesh.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
class IOobjectList;
/*---------------------------------------------------------------------------*\
Class lagrangianFieldReconstructor Declaration
\*---------------------------------------------------------------------------*/
class lagrangianFieldReconstructor
{
// Private Data
//- Reference to complete mesh
const fvMesh& completeMesh_;
//- List of processor meshes
const PtrList<fvMesh>& procMeshes_;
//- The name of the cloud
const word cloudName_;
// Private Member Functions
//- Return whether anything in the object list gets reconstructed
template<class FieldType>
static bool reconstructs
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
);
//- Reconstruct a field
template
<
class Type,
template<class> class IOContainer,
template<class> class IOContainerType
>
tmp<IOContainer<Type>>
reconstructField(const IOobject& fieldIoObject) const;
//- Read, reconstruct and write all fields
template
<
class Type,
template<class> class IOContainer,
template<class> class IOContainerType
>
void reconstructFields
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
);
public:
// Constructors
//- Construct from components
lagrangianFieldReconstructor
(
const fvMesh& completeMesh,
const PtrList<fvMesh>& procMeshes,
const labelListList& faceProcAddressing,
const labelListList& cellProcAddressing,
const word& cloudName
);
//- Disallow default bitwise copy construction
lagrangianFieldReconstructor
(
const lagrangianFieldReconstructor&
) = delete;
//- Destructor
~lagrangianFieldReconstructor();
// Member Functions
//- Return whether anything in the object list gets reconstructed
static bool reconstructs
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
);
//- Read, reconstruct and write all fields
template<class Type>
void reconstructFields
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
);
// Member Operators
//- Disallow default bitwise assignment
void operator=(const lagrangianFieldReconstructor&) = delete;
};
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
#include "lagrangianFieldReconstructorTemplates.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -0,0 +1,179 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "lagrangianFieldReconstructor.H"
#include "IOobjectList.H"
#include "CompactIOField.H"
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
template<class FieldType>
bool Foam::lagrangianFieldReconstructor::reconstructs
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
)
{
IOobjectList fields = objects.lookupClass(FieldType::typeName);
if (fields.size() && selectedFields.empty())
{
return true;
}
forAllConstIter(IOobjectList, fields, fieldIter)
{
if (selectedFields.found(fieldIter()->name()))
{
return true;
}
}
return false;
}
template
<
class Type,
template<class> class IOContainer,
template<class> class IOContainerType
>
Foam::tmp<IOContainer<Type>>
Foam::lagrangianFieldReconstructor::reconstructField
(
const IOobject& fieldIoObject
) const
{
// Construct the complete field
tmp<IOContainer<Type>> tfield
(
new IOContainer<Type>
(
IOobject
(
fieldIoObject.name(),
completeMesh_.time().name(),
cloud::prefix/cloudName_,
completeMesh_,
IOobject::NO_READ,
IOobject::NO_WRITE
),
Field<Type>(0)
)
);
Field<Type>& field = tfield.ref();
// Combine the processor fields into the complete field
forAll(procMeshes_, proci)
{
typeIOobject<IOContainerType<Type>> localIOobject
(
fieldIoObject.name(),
procMeshes_[proci].time().name(),
cloud::prefix/cloudName_,
procMeshes_[proci],
IOobject::MUST_READ,
IOobject::NO_WRITE
);
if (localIOobject.headerOk())
{
IOContainer<Type> fieldi(localIOobject);
label offset = field.size();
field.setSize(offset + fieldi.size());
forAll(fieldi, j)
{
field[offset + j] = fieldi[j];
}
}
}
return tfield;
}
template
<
class Type,
template<class> class IOContainer,
template<class> class IOContainerType
>
void Foam::lagrangianFieldReconstructor::reconstructFields
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
)
{
const word& fieldClassName = IOContainerType<Type>::typeName;
IOobjectList fields = objects.lookupClass(fieldClassName);
if (fields.size())
{
Info<< nl << " Reconstructing lagrangian "
<< fieldClassName << "s" << nl << endl;
forAllConstIter(IOobjectList, fields, fieldIter)
{
if
(
selectedFields.empty()
|| selectedFields.found(fieldIter()->name())
)
{
Info<< " " << fieldIter()->name() << endl;
reconstructField<Type, IOContainer, IOContainerType>
(
*fieldIter()
)().write();
}
}
}
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class Type>
void Foam::lagrangianFieldReconstructor::reconstructFields
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
)
{
reconstructFields<Type, IOField, IOField>
(objects, selectedFields);
reconstructFields<Field<Type>, CompactIOField, IOField>
(objects, selectedFields);
reconstructFields<Field<Type>, CompactIOField, CompactIOField>
(objects, selectedFields);
}
// ************************************************************************* //

View File

@ -2,7 +2,7 @@
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
\\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
@ -38,8 +38,7 @@ Foam::pointFieldReconstructor::pointFieldReconstructor
completeMesh_(completeMesh),
procMeshes_(procMeshes),
pointProcAddressing_(pointProcAddressing),
patchPointAddressing_(procMeshes.size()),
nReconstructed_(0)
patchPointAddressing_(procMeshes.size())
{
// Inverse-addressing of the patch point labels.
labelList pointMap(completeMesh_.size(), -1);
@ -91,4 +90,24 @@ Foam::pointFieldReconstructor::pointFieldReconstructor
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
bool Foam::pointFieldReconstructor::reconstructs
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
)
{
bool result = false;
#define DO_POINT_FIELDS_TYPE(Type, nullArg) \
result = result \
|| reconstructs<PointField<Type>>(objects, selectedFields);
FOR_ALL_FIELD_TYPES(DO_POINT_FIELDS_TYPE)
#undef DO_POINT_FIELDS_TYPE
return result;
}
// ************************************************************************* //

View File

@ -29,6 +29,7 @@ Description
SourceFiles
pointFieldReconstructor.C
pointFieldReconstructorTemplates.C
\*---------------------------------------------------------------------------*/
@ -68,8 +69,21 @@ class pointFieldReconstructor
//- Point patch addressing
labelListListList patchPointAddressing_;
//- Number of fields reconstructed
label nReconstructed_;
// Private Member Functions
//- Return whether anything in the object list gets reconstructed
template<class FieldType>
static bool reconstructs
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
);
//- Reconstruct field
template<class Type>
tmp<PointField<Type>>
reconstructField(const IOobject& fieldIoObject);
public:
@ -90,16 +104,12 @@ public:
// Member Functions
//- Return number of fields reconstructed
label nReconstructed() const
{
return nReconstructed_;
}
//- Reconstruct field
template<class Type>
tmp<PointField<Type>>
reconstructField(const IOobject& fieldIoObject);
//- Return whether anything in the object list gets reconstructed
static bool reconstructs
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
);
//- Reconstruct and write all fields
template<class Type>
@ -124,7 +134,7 @@ public:
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
#include "pointFieldReconstructorReconstructFields.C"
#include "pointFieldReconstructorTemplates.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

View File

@ -27,17 +27,40 @@ License
#include "fvMesh.H"
#include "reversePointPatchFieldMapper.H"
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
template<class FieldType>
bool Foam::pointFieldReconstructor::reconstructs
(
const IOobjectList& objects,
const HashSet<word>& selectedFields
)
{
IOobjectList fields = objects.lookupClass(FieldType::typeName);
if (fields.size() && selectedFields.empty())
{
return true;
}
forAllConstIter(IOobjectList, fields, fieldIter)
{
if (selectedFields.found(fieldIter()->name()))
{
return true;
}
}
return false;
}
template<class Type>
Foam::tmp<Foam::PointField<Type>>
Foam::pointFieldReconstructor::reconstructField(const IOobject& fieldIoObject)
{
// Read the field for all the processors
PtrList<PointField<Type>> procFields
(
procMeshes_.size()
);
PtrList<PointField<Type>> procFields(procMeshes_.size());
forAll(procMeshes_, proci)
{
@ -59,14 +82,12 @@ Foam::pointFieldReconstructor::reconstructField(const IOobject& fieldIoObject)
);
}
// Create the internalField
Field<Type> internalField(completeMesh_.size());
// Create the patch fields
PtrList<pointPatchField<Type>> patchFields(completeMesh_.boundary().size());
forAll(procMeshes_, proci)
{
const PointField<Type>&
@ -122,8 +143,7 @@ Foam::pointFieldReconstructor::reconstructField(const IOobject& fieldIoObject)
}
}
// Construct and write the field
// setting the internalField and patchFields
// Construct and return the field
return tmp<PointField<Type>>
(
new PointField<Type>
@ -145,7 +165,8 @@ Foam::pointFieldReconstructor::reconstructField(const IOobject& fieldIoObject)
}
// Reconstruct and write all point fields
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class Type>
void Foam::pointFieldReconstructor::reconstructFields
(
@ -162,7 +183,8 @@ void Foam::pointFieldReconstructor::reconstructFields
if (fields.size())
{
Info<< " Reconstructing " << fieldClassName << "s\n" << endl;
Info<< nl << " Reconstructing " << fieldClassName << "s"
<< nl << endl;
forAllConstIter(IOobjectList, fields, fieldIter)
{
@ -175,12 +197,8 @@ void Foam::pointFieldReconstructor::reconstructFields
Info<< " " << fieldIter()->name() << endl;
reconstructField<Type>(*fieldIter())().write();
nReconstructed_++;
}
}
Info<< endl;
}
}

View File

@ -1,117 +0,0 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
InClass
Foam::reconstructLagrangian
Description
SourceFiles
reconstructLagrangianPositions.C
reconstructLagrangianFields.C
\*---------------------------------------------------------------------------*/
#ifndef reconstructLagrangian_H
#define reconstructLagrangian_H
#include "cloud.H"
#include "polyMesh.H"
#include "IOobjectList.H"
#include "CompactIOField.H"
#include "fvMesh.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
void reconstructLagrangianPositions
(
const polyMesh& mesh,
const word& cloudName,
const PtrList<fvMesh>& meshes,
const labelListList& faceProcAddressing,
const labelListList& cellProcAddressing
);
template<class Type>
tmp<IOField<Type>> reconstructLagrangianField
(
const word& cloudName,
const polyMesh& mesh,
const PtrList<fvMesh>& meshes,
const word& fieldName
);
template<class Type>
tmp<CompactIOField<Field<Type>>> reconstructLagrangianFieldField
(
const word& cloudName,
const polyMesh& mesh,
const PtrList<fvMesh>& meshes,
const word& fieldName
);
template<class Type>
void reconstructLagrangianFields
(
const word& cloudName,
const polyMesh& mesh,
const PtrList<fvMesh>& meshes,
const IOobjectList& objects,
const HashSet<word>& selectedFields
);
template<class Type>
void reconstructLagrangianFieldFields
(
const word& cloudName,
const polyMesh& mesh,
const PtrList<fvMesh>& meshes,
const IOobjectList& objects,
const HashSet<word>& selectedFields
);
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
} // End namespace Foam
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#ifdef NoRepository
#include "reconstructLagrangianFields.C"
#endif
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
#endif
// ************************************************************************* //

View File

@ -1,273 +0,0 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "IOField.H"
#include "CompactIOField.H"
#include "Time.H"
// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
template<class Type>
Foam::tmp<Foam::IOField<Type>> Foam::reconstructLagrangianField
(
const word& cloudName,
const polyMesh& mesh,
const PtrList<fvMesh>& meshes,
const word& fieldName
)
{
// Construct empty field on mesh
tmp<IOField<Type>> tfield
(
new IOField<Type>
(
IOobject
(
fieldName,
mesh.time().name(),
cloud::prefix/cloudName,
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
Field<Type>(0)
)
);
Field<Type>& field = tfield.ref();
forAll(meshes, i)
{
// Check object on local mesh
typeIOobject<IOField<Type>> localIOobject
(
fieldName,
meshes[i].time().name(),
cloud::prefix/cloudName,
meshes[i],
IOobject::MUST_READ,
IOobject::NO_WRITE
);
if (localIOobject.headerOk())
{
IOField<Type> fieldi(localIOobject);
label offset = field.size();
field.setSize(offset + fieldi.size());
forAll(fieldi, j)
{
field[offset + j] = fieldi[j];
}
}
}
return tfield;
}
template<class Type>
Foam::tmp<Foam::CompactIOField<Foam::Field<Type>>>
Foam::reconstructLagrangianFieldField
(
const word& cloudName,
const polyMesh& mesh,
const PtrList<fvMesh>& meshes,
const word& fieldName
)
{
// Construct empty field on mesh
tmp<CompactIOField<Field<Type>>> tfield
(
new CompactIOField<Field<Type>>
(
IOobject
(
fieldName,
mesh.time().name(),
cloud::prefix/cloudName,
mesh,
IOobject::NO_READ,
IOobject::NO_WRITE
),
Field<Field<Type>>(0)
)
);
Field<Field<Type>>& field = tfield.ref();
forAll(meshes, i)
{
// No type checking is done to handle CompactIOField and IOField
IOobject localIOobject
(
fieldName,
meshes[i].time().name(),
cloud::prefix/cloudName,
meshes[i],
IOobject::MUST_READ,
IOobject::NO_WRITE
);
if (localIOobject.headerOk())
{
CompactIOField<Field<Type>> fieldi(localIOobject);
label offset = field.size();
field.setSize(offset + fieldi.size());
forAll(fieldi, j)
{
field[offset + j] = fieldi[j];
}
}
}
return tfield;
}
template<class Type>
void Foam::reconstructLagrangianFields
(
const word& cloudName,
const polyMesh& mesh,
const PtrList<fvMesh>& meshes,
const IOobjectList& objects,
const HashSet<word>& selectedFields
)
{
const word fieldClassName(IOField<Type>::typeName);
IOobjectList fields = objects.lookupClass(fieldClassName);
if (fields.size())
{
Info<< " Reconstructing lagrangian "
<< fieldClassName << "s\n" << endl;
forAllConstIter(IOobjectList, fields, fieldIter)
{
if
(
selectedFields.empty()
|| selectedFields.found(fieldIter()->name())
)
{
Info<< " " << fieldIter()->name() << endl;
reconstructLagrangianField<Type>
(
cloudName,
mesh,
meshes,
fieldIter()->name()
)().write();
}
}
Info<< endl;
}
}
template<class Type>
void Foam::reconstructLagrangianFieldFields
(
const word& cloudName,
const polyMesh& mesh,
const PtrList<fvMesh>& meshes,
const IOobjectList& objects,
const HashSet<word>& selectedFields
)
{
{
const word fieldClassName(CompactIOField<Field<Type>>::typeName);
IOobjectList fields = objects.lookupClass(fieldClassName);
if (fields.size())
{
Info<< " Reconstructing lagrangian "
<< fieldClassName << "s\n" << endl;
forAllConstIter(IOobjectList, fields, fieldIter)
{
if
(
selectedFields.empty()
|| selectedFields.found(fieldIter()->name())
)
{
Info<< " " << fieldIter()->name() << endl;
reconstructLagrangianFieldField<Type>
(
cloudName,
mesh,
meshes,
fieldIter()->name()
)().write();
}
}
Info<< endl;
}
}
{
const word fieldClassName(IOField<Field<Type>>::typeName);
IOobjectList fields = objects.lookupClass(fieldClassName);
if (fields.size())
{
Info<< " Reconstructing lagrangian "
<< fieldClassName << "s\n" << endl;
forAllConstIter(IOobjectList, fields, fieldIter)
{
if
(
selectedFields.empty()
|| selectedFields.found(fieldIter()->name())
)
{
Info<< " " << fieldIter()->name() << endl;
reconstructLagrangianFieldField<Type>
(
cloudName,
mesh,
meshes,
fieldIter()->name()
)().write();
}
}
Info<< endl;
}
}
}
// ************************************************************************* //

View File

@ -1,83 +0,0 @@
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration | Website: https://openfoam.org
\\ / A nd | Copyright (C) 2011-2022 OpenFOAM Foundation
\\/ M anipulation |
-------------------------------------------------------------------------------
License
This file is part of OpenFOAM.
OpenFOAM is free software: you can redistribute it and/or modify it
under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
for more details.
You should have received a copy of the GNU General Public License
along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
\*---------------------------------------------------------------------------*/
#include "reconstructLagrangian.H"
#include "labelIOList.H"
#include "passiveParticleCloud.H"
// * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
void Foam::reconstructLagrangianPositions
(
const polyMesh& mesh,
const word& cloudName,
const PtrList<fvMesh>& meshes,
const labelListList& faceProcAddressing,
const labelListList& cellProcAddressing
)
{
passiveParticleCloud lagrangianPositions
(
mesh,
cloudName,
IDLList<passiveParticle>()
);
forAll(meshes, i)
{
const labelList& cellMap = cellProcAddressing[i];
const labelList& faceMap = faceProcAddressing[i];
Cloud<passiveParticle> lpi(meshes[i], cloudName, false);
forAllConstIter(Cloud<passiveParticle>, lpi, iter)
{
const passiveParticle& ppi = iter();
const label mappedCell = cellMap[ppi.cell()];
// Inverting sign if necessary and subtracting 1 from
// faceProcAddressing
label mappedTetFace = mag(faceMap[ppi.tetFace()]) - 1;
lagrangianPositions.append
(
new passiveParticle
(
mesh,
ppi.coordinates(),
mappedCell,
mappedTetFace,
ppi.procTetPt(meshes[i], mesh, mappedCell, mappedTetFace)
)
);
}
}
IOPosition<Cloud<passiveParticle>>(lagrangianPositions).write();
}
// ************************************************************************* //

View File

@ -34,10 +34,10 @@ Description
#include "timeSelector.H"
#include "IOobjectList.H"
#include "processorRunTimes.H"
#include "domainDecomposition.H"
#include "multiDomainDecomposition.H"
#include "fvFieldReconstructor.H"
#include "pointFieldReconstructor.H"
#include "reconstructLagrangian.H"
#include "lagrangianFieldReconstructor.H"
using namespace Foam;
@ -46,21 +46,37 @@ using namespace Foam;
namespace Foam
{
bool haveAllTimes
bool haveUniform
(
const HashSet<word>& masterTimeDirSet,
const instantList& timeDirs
const processorRunTimes& runTimes,
const word& regionDir = word::null
)
{
// Loop over all times
forAll(timeDirs, timei)
{
if (!masterTimeDirSet.found(timeDirs[timei].name()))
{
return false;
}
}
return true;
return
fileHandler().isDir
(
fileHandler().filePath
(
runTimes.procTimes()[0].timePath()/regionDir/"uniform"
)
);
}
void reconstructUniform
(
const processorRunTimes& runTimes,
const word& regionDir = word::null
)
{
fileHandler().cp
(
fileHandler().filePath
(
runTimes.procTimes()[0].timePath()/regionDir/"uniform"
),
runTimes.completeTime().timePath()/regionDir
);
}
@ -85,14 +101,41 @@ void writeDecomposition(const domainDecomposition& meshes)
cellProc.write();
Info<< "Wrote decomposition as volScalarField::Internal to "
<< cellProc.name() << " for use in postprocessing."
<< nl << endl;
<< cellProc.name() << " for use in postprocessing"
<< endl;
}
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
namespace Foam
{
class delayedNewLine
{
mutable bool first_;
public:
delayedNewLine()
:
first_(true)
{}
friend Ostream& operator<<(Ostream& os, const delayedNewLine& dnl)
{
if (!dnl.first_) os << nl;
dnl.first_ = false;
return os;
}
};
}
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
int main(int argc, char *argv[])
@ -102,9 +145,6 @@ int main(int argc, char *argv[])
"Reconstruct fields of a parallel case"
);
// Enable -constant ... if someone really wants it
// Enable -withZero to prevent accidentally trashing the initial fields
timeSelector::addOptions(true, true);
argList::noParallel();
#include "addRegionOption.H"
#include "addAllRegionsOption.H"
@ -112,7 +152,7 @@ int main(int argc, char *argv[])
(
"cellProc",
"write cell processor indices as a volScalarField::Internal for "
"post-processing."
"post-processing"
);
argList::addOption
(
@ -132,7 +172,7 @@ int main(int argc, char *argv[])
"list",
"specify a list of lagrangian fields to be reconstructed. Eg, '(U d)' -"
"regular expressions not currently supported, "
"positions always included."
"positions always included"
);
argList::addBoolOption
(
@ -150,6 +190,10 @@ int main(int argc, char *argv[])
"only reconstruct new times (i.e. that do not exist already)"
);
// Include explicit constant options, and explicit zero option (to prevent
// the user accidentally trashing the initial fields)
timeSelector::addOptions(true, true);
#include "setRootCase.H"
const bool writeCellProc = args.optionFound("cellProc");
@ -164,8 +208,7 @@ int main(int argc, char *argv[])
if (noFields)
{
Info<< "Skipping reconstructing fields"
<< nl << endl;
Info<< "Skipping reconstructing fields" << nl << endl;
}
const bool noLagrangian = args.optionFound("noLagrangian");
@ -191,7 +234,7 @@ int main(int argc, char *argv[])
{
FatalErrorInFunction
<< "Cannot specify noLagrangian and lagrangianFields "
<< "options together."
<< "options together"
<< exit(FatalError);
}
@ -199,11 +242,11 @@ int main(int argc, char *argv[])
}
// Set time from database
Info<< "Create time\n" << endl;
Info<< "Create time" << nl << endl;
processorRunTimes runTimes(Foam::Time::controlDictName, args);
// Allow override of time
const instantList times = runTimes.selectProc(args);
// Get the times to reconstruct
instantList times = runTimes.selectProc(args);
const Time& runTime = runTimes.procTimes()[0];
@ -228,484 +271,335 @@ int main(int argc, char *argv[])
// Warn fileHandler of number of processors
const_cast<fileOperation&>(fileHandler()).setNProcs(nProcs);
// Note that we do not set the runTime time so it is still the
// one set through the controlDict. The -time option
// only affects the selected set of times from processor0.
// - can be illogical
// + any point motion handled through mesh.readUpdate
// Quit if no times
if (times.empty())
{
WarningInFunction << "No times selected" << endl;
WarningInFunction << "No times selected" << nl << endl;
exit(1);
}
// Get current times if -newTimes
const bool newTimes = args.optionFound("newTimes");
instantList masterTimeDirs;
if (newTimes)
// If only reconstructing new times then filter out existing times
if (args.optionFound("newTimes"))
{
masterTimeDirs = runTimes.completeTime().times();
// Get all existing times
const instantList existingTimes = runTimes.completeTime().times();
// Put into a set
HashSet<word> existingTimesSet;
existingTimesSet.resize(2*existingTimes.size());
forAll(existingTimes, i)
{
existingTimesSet.insert(existingTimes[i].name());
}
// Remove times from the existing time set by shuffling up
label timei = 0;
forAll(times, timej)
{
if (!existingTimesSet.found(times[timej].name()))
{
times[timei ++] = times[timej];
}
}
times.resize(timei);
}
HashSet<word> masterTimeDirSet(2*masterTimeDirs.size());
forAll(masterTimeDirs, i)
// Quit if no times
if (times.empty())
{
masterTimeDirSet.insert(masterTimeDirs[i].name());
}
if
(
newTimes
&& regionNames.size() == 1
&& regionNames[0] == fvMesh::defaultRegion
&& haveAllTimes(masterTimeDirSet, times)
)
{
Info<< "All times already reconstructed.\n\nEnd\n" << endl;
Info<< "All times already reconstructed" << nl << nl
<< "End" << nl << endl;
return 0;
}
// Reconstruct all regions
forAll(regionNames, regioni)
// Create meshes
multiDomainDecomposition regionMeshes(runTimes, regionNames);
if (regionMeshes.readReconstruct(!noReconstructSets))
{
const word& regionName = regionNames[regioni];
Info<< endl;
const word& regionDir =
regionName == polyMesh::defaultRegion
? word::null
: regionName;
// Create meshes
Info<< "\n\nReconstructing mesh " << regionName << nl << endl;
domainDecomposition meshes(runTimes, regionName);
if (meshes.readReconstruct(!noReconstructSets) && writeCellProc)
if (writeCellProc)
{
writeDecomposition(meshes);
fileHandler().flush();
}
// Loop over all times
forAll(times, timei)
{
if (newTimes && masterTimeDirSet.found(times[timei].name()))
forAll(regionNames, regioni)
{
Info<< "Skipping time " << times[timei].name()
<< endl << endl;
continue;
}
// Set the time
runTimes.setTime(times[timei], timei);
Info<< "Time = " << runTimes.completeTime().userTimeName()
<< nl << endl;
// Update the meshes
const fvMesh::readUpdateState state =
meshes.readUpdateReconstruct();
// Write the mesh out, if necessary
if (state != fvMesh::UNCHANGED)
{
meshes.writeComplete(!noReconstructSets);
}
// Write the decomposition, if necessary
if
(
writeCellProc
&& meshes.completeMesh().facesInstance()
== runTimes.completeTime().name()
)
{
writeDecomposition(meshes);
writeDecomposition(regionMeshes.meshes(regioni)());
Info<< endl;
fileHandler().flush();
}
}
}
// Get list of objects from processor0 database
IOobjectList objects
(
meshes.procMeshes()[0],
runTimes.procTimes()[0].name()
);
// Loop over all times
forAll(times, timei)
{
// Set the time
runTimes.setTime(times[timei], timei);
if (!noFields)
Info<< "Time = " << runTimes.completeTime().userTimeName()
<< nl << endl;
// Update the meshes
const fvMesh::readUpdateState stat =
regionMeshes.readUpdateReconstruct();
if (stat >= fvMesh::TOPO_CHANGE) Info<< endl;
// Write the mesh out (if anything has changed)
regionMeshes.writeComplete(!noReconstructSets);
// Write the decomposition, if necessary
forAll(regionNames, regioni)
{
if (writeCellProc && stat >= fvMesh::TOPO_CHANGE)
{
// If there are any FV fields, reconstruct them
Info<< "Reconstructing FV fields" << nl << endl;
fvFieldReconstructor fvReconstructor
(
meshes.completeMesh(),
meshes.procMeshes(),
meshes.procFaceAddressing(),
meshes.procCellAddressing(),
meshes.procFaceAddressingBf()
);
fvReconstructor.reconstructFvVolumeInternalFields<scalar>
(
objects,
selectedFields
);
fvReconstructor.reconstructFvVolumeInternalFields<vector>
(
objects,
selectedFields
);
fvReconstructor.reconstructFvVolumeInternalFields
<sphericalTensor>
(
objects,
selectedFields
);
fvReconstructor.reconstructFvVolumeInternalFields<symmTensor>
(
objects,
selectedFields
);
fvReconstructor.reconstructFvVolumeInternalFields<tensor>
(
objects,
selectedFields
);
fvReconstructor.reconstructFvVolumeFields<scalar>
(
objects,
selectedFields
);
fvReconstructor.reconstructFvVolumeFields<vector>
(
objects,
selectedFields
);
fvReconstructor.reconstructFvVolumeFields<sphericalTensor>
(
objects,
selectedFields
);
fvReconstructor.reconstructFvVolumeFields<symmTensor>
(
objects,
selectedFields
);
fvReconstructor.reconstructFvVolumeFields<tensor>
(
objects,
selectedFields
);
fvReconstructor.reconstructFvSurfaceFields<scalar>
(
objects,
selectedFields
);
fvReconstructor.reconstructFvSurfaceFields<vector>
(
objects,
selectedFields
);
fvReconstructor.reconstructFvSurfaceFields<sphericalTensor>
(
objects,
selectedFields
);
fvReconstructor.reconstructFvSurfaceFields<symmTensor>
(
objects,
selectedFields
);
fvReconstructor.reconstructFvSurfaceFields<tensor>
(
objects,
selectedFields
);
if (fvReconstructor.nReconstructed() == 0)
{
Info<< "No FV fields" << nl << endl;
}
writeDecomposition(regionMeshes.meshes(regioni)());
Info<< endl;
fileHandler().flush();
}
}
if (!noFields)
// Do a region-by-region reconstruction of all the available fields
forAll(regionNames, regioni)
{
const word& regionName = regionNames[regioni];
const word regionDir =
regionName == polyMesh::defaultRegion ? word::null : regionName;
const delayedNewLine dnl;
// Prefixed scope
{
Info<< "Reconstructing point fields" << nl << endl;
const RegionConstRef<domainDecomposition> meshes =
regionMeshes.meshes(regioni);
const pointMesh& completePMesh =
pointMesh::New(meshes.completeMesh());
pointFieldReconstructor pointReconstructor
// Search for objects at this time
IOobjectList objects
(
completePMesh,
meshes.procMeshes(),
meshes.procPointAddressing()
meshes().procMeshes()[0],
runTimes.procTimes()[0].name()
);
pointReconstructor.reconstructFields<scalar>
(
objects,
selectedFields
);
pointReconstructor.reconstructFields<vector>
(
objects,
selectedFields
);
pointReconstructor.reconstructFields<sphericalTensor>
(
objects,
selectedFields
);
pointReconstructor.reconstructFields<symmTensor>
(
objects,
selectedFields
);
pointReconstructor.reconstructFields<tensor>
(
objects,
selectedFields
);
if (pointReconstructor.nReconstructed() == 0)
if (!noFields)
{
Info<< "No point fields" << nl << endl;
}
}
Info<< dnl << "Reconstructing FV fields" << endl;
// If there are any clouds, reconstruct them.
// The problem is that a cloud of size zero will not get written so
// in pass 1 we determine the cloud names and per cloud name the
// fields. Note that the fields are stored as IOobjectList from
// the first processor that has them. They are in pass2 only used
// for name and type (scalar, vector etc).
if (!noLagrangian)
{
HashTable<IOobjectList> cloudObjects;
forAll(runTimes.procTimes(), proci)
{
fileName lagrangianDir
if
(
fileHandler().filePath
fvFieldReconstructor::reconstructs
(
runTimes.procTimes()[proci].timePath()
/regionDir
/cloud::prefix
objects,
selectedFields
)
);
fileNameList cloudDirs;
if (!lagrangianDir.empty())
)
{
cloudDirs = fileHandler().readDir
fvFieldReconstructor fvReconstructor
(
lagrangianDir,
fileType::directory
meshes().completeMesh(),
meshes().procMeshes(),
meshes().procFaceAddressing(),
meshes().procCellAddressing(),
meshes().procFaceAddressingBf()
);
#define DO_FV_VOL_INTERNAL_FIELDS_TYPE(Type, nullArg) \
fvReconstructor.reconstructVolInternalFields<Type> \
(objects, selectedFields);
FOR_ALL_FIELD_TYPES(DO_FV_VOL_INTERNAL_FIELDS_TYPE)
#undef DO_FV_VOL_INTERNAL_FIELDS_TYPE
#define DO_FV_VOL_FIELDS_TYPE(Type, nullArg) \
fvReconstructor.reconstructVolFields<Type> \
(objects, selectedFields);
FOR_ALL_FIELD_TYPES(DO_FV_VOL_FIELDS_TYPE)
#undef DO_FV_VOL_FIELDS_TYPE
#define DO_FV_SURFACE_FIELDS_TYPE(Type, nullArg) \
fvReconstructor.reconstructFvSurfaceFields<Type> \
(objects, selectedFields);
FOR_ALL_FIELD_TYPES(DO_FV_SURFACE_FIELDS_TYPE)
#undef DO_FV_SURFACE_FIELDS_TYPE
}
forAll(cloudDirs, i)
else
{
// Check if we already have cloud objects for this
// cloudname
HashTable<IOobjectList>::const_iterator iter =
cloudObjects.find(cloudDirs[i]);
Info<< dnl << " (no FV fields)" << endl;
}
}
if (iter == cloudObjects.end())
{
// Do local scan for valid cloud objects
IOobjectList sprayObjs
if (!noFields)
{
Info<< dnl << "Reconstructing point fields" << endl;
if
(
pointFieldReconstructor::reconstructs
(
objects,
selectedFields
)
)
{
pointFieldReconstructor pointReconstructor
(
pointMesh::New(meshes().completeMesh()),
meshes().procMeshes(),
meshes().procPointAddressing()
);
#define DO_POINT_FIELDS_TYPE(Type, nullArg) \
pointReconstructor.reconstructFields<Type> \
(objects, selectedFields);
FOR_ALL_FIELD_TYPES(DO_POINT_FIELDS_TYPE)
#undef DO_POINT_FIELDS_TYPE
}
else
{
Info<< dnl << " (no point fields)" << endl;
}
}
if (!noLagrangian)
{
// Search for clouds that exist on any processor and add
// them into this table of cloud objects
HashTable<IOobjectList> cloudsObjects;
forAll(runTimes.procTimes(), proci)
{
// Find cloud directories
fileNameList cloudDirs
(
fileHandler().readDir
(
meshes.procMeshes()[proci],
fileHandler().filePath
(
runTimes.procTimes()[proci].timePath()
/regionDir
/cloud::prefix
),
fileType::directory
)
);
// Add objects in any found cloud directories
forAll(cloudDirs, i)
{
// Pass if we already have an objects for this name
HashTable<IOobjectList>::const_iterator iter =
cloudsObjects.find(cloudDirs[i]);
if (iter != cloudsObjects.end()) continue;
// Do local scan for valid cloud objects
IOobjectList cloudObjs
(
meshes().procMeshes()[proci],
runTimes.procTimes()[proci].name(),
cloud::prefix/cloudDirs[i]
cloud::prefix/cloudDirs[i],
IOobject::MUST_READ,
IOobject::NO_WRITE,
false
);
IOobject* positionsPtr =
sprayObjs.lookup(word("positions"));
if (positionsPtr)
// If "positions" is present, then add to the table
if (cloudObjs.lookup(word("positions")))
{
cloudObjects.insert(cloudDirs[i], sprayObjs);
cloudsObjects.insert(cloudDirs[i], cloudObjs);
}
}
}
// Reconstruct the objects found above
if (cloudsObjects.size())
{
forAllConstIter
(
HashTable<IOobjectList>,
cloudsObjects,
iter
)
{
const word cloudName =
string::validate<word>(iter.key());
const IOobjectList& cloudObjects = iter();
Info<< dnl << "Reconstructing lagrangian fields "
<< "for cloud " << cloudName << endl;
if
(
lagrangianFieldReconstructor::reconstructs
(
cloudObjects,
selectedLagrangianFields
)
)
{
lagrangianFieldReconstructor
lagrangianReconstructor
(
meshes().completeMesh(),
meshes().procMeshes(),
meshes().procFaceAddressing(),
meshes().procCellAddressing(),
cloudName
);
#define DO_CLOUD_FIELDS_TYPE(Type, nullArg) \
lagrangianReconstructor \
.reconstructFields<Type> \
(cloudObjects, selectedLagrangianFields);
DO_CLOUD_FIELDS_TYPE(label, );
FOR_ALL_FIELD_TYPES(DO_CLOUD_FIELDS_TYPE)
#undef DO_CLOUD_FIELDS_TYPE
}
else
{
Info<< dnl << " (no lagrangian fields)"
<< endl;
}
}
}
}
if (cloudObjects.size())
{
// Pass2: reconstruct the cloud
forAllConstIter(HashTable<IOobjectList>, cloudObjects, iter)
{
const word cloudName =
string::validate<word>(iter.key());
// Objects (on arbitrary processor)
const IOobjectList& sprayObjs = iter();
Info<< "Reconstructing lagrangian fields for cloud "
<< cloudName << nl << endl;
reconstructLagrangianPositions
(
meshes.completeMesh(),
cloudName,
meshes.procMeshes(),
meshes.procFaceAddressing(),
meshes.procCellAddressing()
);
reconstructLagrangianFields<label>
(
cloudName,
meshes.completeMesh(),
meshes.procMeshes(),
sprayObjs,
selectedLagrangianFields
);
reconstructLagrangianFieldFields<label>
(
cloudName,
meshes.completeMesh(),
meshes.procMeshes(),
sprayObjs,
selectedLagrangianFields
);
reconstructLagrangianFields<scalar>
(
cloudName,
meshes.completeMesh(),
meshes.procMeshes(),
sprayObjs,
selectedLagrangianFields
);
reconstructLagrangianFieldFields<scalar>
(
cloudName,
meshes.completeMesh(),
meshes.procMeshes(),
sprayObjs,
selectedLagrangianFields
);
reconstructLagrangianFields<vector>
(
cloudName,
meshes.completeMesh(),
meshes.procMeshes(),
sprayObjs,
selectedLagrangianFields
);
reconstructLagrangianFieldFields<vector>
(
cloudName,
meshes.completeMesh(),
meshes.procMeshes(),
sprayObjs,
selectedLagrangianFields
);
reconstructLagrangianFields<sphericalTensor>
(
cloudName,
meshes.completeMesh(),
meshes.procMeshes(),
sprayObjs,
selectedLagrangianFields
);
reconstructLagrangianFieldFields<sphericalTensor>
(
cloudName,
meshes.completeMesh(),
meshes.procMeshes(),
sprayObjs,
selectedLagrangianFields
);
reconstructLagrangianFields<symmTensor>
(
cloudName,
meshes.completeMesh(),
meshes.procMeshes(),
sprayObjs,
selectedLagrangianFields
);
reconstructLagrangianFieldFields<symmTensor>
(
cloudName,
meshes.completeMesh(),
meshes.procMeshes(),
sprayObjs,
selectedLagrangianFields
);
reconstructLagrangianFields<tensor>
(
cloudName,
meshes.completeMesh(),
meshes.procMeshes(),
sprayObjs,
selectedLagrangianFields
);
reconstructLagrangianFieldFields<tensor>
(
cloudName,
meshes.completeMesh(),
meshes.procMeshes(),
sprayObjs,
selectedLagrangianFields
);
}
}
else
{
Info<< "No lagrangian fields" << nl << endl;
}
}
// If there is a "uniform" directory in the time region
// directory copy from the master processor
Info<< dnl;
}
// Collect the uniform directory
if (haveUniform(runTimes))
{
Info<< "Collecting uniform files" << endl;
reconstructUniform(runTimes);
Info<< endl;
}
if (regionNames == wordList(1, polyMesh::defaultRegion)) continue;
// Collect the region uniform directories
forAll(regionNames, regioni)
{
const word& regionName = regionNames[regioni];
const word regionDir =
regionName == polyMesh::defaultRegion ? word::null : regionName;
if (haveUniform(runTimes, regionDir))
{
fileName uniformDir0
(
fileHandler().filePath
(
runTimes.procTimes()[0].timePath()/regionDir/"uniform"
)
);
if (!uniformDir0.empty() && fileHandler().isDir(uniformDir0))
// Prefixed scope
{
fileHandler().cp
(
uniformDir0,
runTimes.completeTime().timePath()/regionDir
);
}
}
const RegionConstRef<domainDecomposition> meshes =
regionMeshes.meshes(regioni);
// For the first region of a multi-region case additionally
// copy the "uniform" directory in the time directory
if (regioni == 0 && regionDir != word::null)
{
fileName uniformDir0
(
fileHandler().filePath
(
runTimes.procTimes()[0].timePath()/"uniform"
)
);
Info<< "Collecting uniform files" << endl;
if (!uniformDir0.empty() && fileHandler().isDir(uniformDir0))
{
fileHandler().cp
(
uniformDir0,
runTimes.completeTime().timePath()
);
reconstructUniform(runTimes, regionDir);
}
Info<< endl;
}
}
}
Info<< "\nEnd\n" << endl;
Info<< "End" << nl << endl;
return 0;
}