/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
\\ / A nd | www.openfoam.com
\\/ M anipulation |
-------------------------------------------------------------------------------
Copyright (C) 2015 OpenFOAM Foundation
Copyright (C) 2016-2018 OpenCFD Ltd.
-------------------------------------------------------------------------------
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 .
\*---------------------------------------------------------------------------*/
#include "parFvFieldReconstructor.H"
#include "Time.H"
#include "PtrList.H"
#include "fvPatchFields.H"
#include "emptyFvPatch.H"
#include "emptyFvPatchField.H"
#include "emptyFvsPatchField.H"
#include "IOobjectList.H"
#include "mapDistributePolyMesh.H"
#include "processorFvPatch.H"
#include "directFvPatchFieldMapper.H"
#include "distributedUnallocatedDirectFieldMapper.H"
#include "distributedUnallocatedDirectFvPatchFieldMapper.H"
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template
Foam::tmp>
Foam::parFvFieldReconstructor::reconstructFvVolumeInternalField
(
const DimensionedField& fld
) const
{
distributedUnallocatedDirectFieldMapper mapper
(
labelUList::null(),
distMap_.cellMap()
);
Field internalField(fld, mapper);
// Construct a volField
IOobject baseIO
(
fld.name(),
baseMesh_.time().timeName(),
fld.local(),
baseMesh_,
IOobject::NO_READ,
IOobject::NO_WRITE
);
auto tfield = tmp>::New
(
baseIO,
baseMesh_,
fld.dimensions(),
internalField
);
tfield.ref().oriented() = fld.oriented();
return tfield;
}
template
Foam::tmp>
Foam::parFvFieldReconstructor::reconstructFvVolumeInternalField
(
const IOobject& fieldIoObject
) const
{
// Read the field
DimensionedField fld
(
fieldIoObject,
procMesh_
);
// Distribute onto baseMesh
return reconstructFvVolumeInternalField(fld);
}
// Reconstruct a field onto the baseMesh
template
Foam::tmp>
Foam::parFvFieldReconstructor::reconstructFvVolumeField
(
const GeometricField& fld
) const
{
// Create the internalField by remote mapping
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
distributedUnallocatedDirectFieldMapper mapper
(
labelUList::null(),
distMap_.cellMap()
);
Field internalField(fld.internalField(), mapper);
// Create the patchFields by remote mapping
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
// Note: patchFields still on mesh, not baseMesh
PtrList> patchFields(fld.mesh().boundary().size());
const typename GeometricField::Boundary&
bfld = fld.boundaryField();
forAll(bfld, patchI)
{
if (patchFaceMaps_.set(patchI))
{
// Clone local patch field
patchFields.set(patchI, bfld[patchI].clone());
distributedUnallocatedDirectFvPatchFieldMapper mapper
(
labelUList::null(),
patchFaceMaps_[patchI]
);
// Map into local copy
patchFields[patchI].autoMap(mapper);
}
}
PtrList> basePatchFields
(
baseMesh_.boundary().size()
);
// Clone the patchFields onto the base patches. This is just to reset
// the reference to the patch, size and content stay the same.
forAll(patchFields, patchI)
{
if (patchFields.set(patchI))
{
const fvPatch& basePatch = baseMesh_.boundary()[patchI];
const fvPatchField& pfld = patchFields[patchI];
labelList dummyMap(identity(pfld.size()));
directFvPatchFieldMapper dummyMapper(dummyMap);
basePatchFields.set
(
patchI,
fvPatchField::New
(
pfld,
basePatch,
DimensionedField::null(),
dummyMapper
)
);
}
}
// Add some empty patches on remaining patches (tbd.probably processor
// patches)
forAll(basePatchFields, patchI)
{
if (patchI >= patchFields.size() || !patchFields.set(patchI))
{
basePatchFields.set
(
patchI,
fvPatchField::New
(
emptyFvPatchField::typeName,
baseMesh_.boundary()[patchI],
DimensionedField::null()
)
);
}
}
// Construct a volField
IOobject baseIO
(
fld.name(),
baseMesh_.time().timeName(),
fld.local(),
baseMesh_,
IOobject::NO_READ,
IOobject::NO_WRITE
);
auto tfield = tmp>::New
(
baseIO,
baseMesh_,
fld.dimensions(),
internalField,
basePatchFields
);
tfield.ref().oriented()= fld.oriented();
return tfield;
}
template
Foam::tmp>
Foam::parFvFieldReconstructor::reconstructFvVolumeField
(
const IOobject& fieldIoObject
) const
{
// Read the field
GeometricField fld
(
fieldIoObject,
procMesh_
);
// Distribute onto baseMesh
return reconstructFvVolumeField(fld);
}
template
Foam::tmp>
Foam::parFvFieldReconstructor::reconstructFvSurfaceField
(
const GeometricField& fld
) const
{
// Create the internalField by remote mapping
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
distributedUnallocatedDirectFieldMapper mapper
(
labelUList::null(),
distMap_.faceMap()
);
// Create flat field of internalField + all patch fields
Field flatFld(fld.mesh().nFaces(), Type(Zero));
SubList(flatFld, fld.internalField().size()) = fld.internalField();
forAll(fld.boundaryField(), patchI)
{
const fvsPatchField& fvp = fld.boundaryField()[patchI];
SubList(flatFld, fvp.size(), fvp.patch().start()) = fvp;
}
// Map all faces
Field internalField(flatFld, mapper, fld.oriented()());
// Trim to internal faces (note: could also have special mapper)
internalField.setSize
(
min
(
internalField.size(),
baseMesh_.nInternalFaces()
)
);
// Create the patchFields by remote mapping
// ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
// Note: patchFields still on mesh, not baseMesh
PtrList> patchFields(fld.mesh().boundary().size());
const typename GeometricField::Boundary&
bfld = fld.boundaryField();
forAll(bfld, patchI)
{
if (patchFaceMaps_.set(patchI))
{
// Clone local patch field
patchFields.set(patchI, bfld[patchI].clone());
distributedUnallocatedDirectFvPatchFieldMapper mapper
(
labelUList::null(),
patchFaceMaps_[patchI]
);
// Map into local copy
patchFields[patchI].autoMap(mapper);
}
}
PtrList> basePatchFields
(
baseMesh_.boundary().size()
);
// Clone the patchFields onto the base patches. This is just to reset
// the reference to the patch, size and content stay the same.
forAll(patchFields, patchI)
{
if (patchFields.set(patchI))
{
const fvPatch& basePatch = baseMesh_.boundary()[patchI];
const fvsPatchField& pfld = patchFields[patchI];
labelList dummyMap(identity(pfld.size()));
directFvPatchFieldMapper dummyMapper(dummyMap);
basePatchFields.set
(
patchI,
fvsPatchField::New
(
pfld,
basePatch,
DimensionedField::null(),
dummyMapper
)
);
}
}
// Add some empty patches on remaining patches (tbd.probably processor
// patches)
forAll(basePatchFields, patchI)
{
if (patchI >= patchFields.size() || !patchFields.set(patchI))
{
basePatchFields.set
(
patchI,
fvsPatchField::New
(
emptyFvsPatchField::typeName,
baseMesh_.boundary()[patchI],
DimensionedField::null()
)
);
}
}
// Construct a volField
IOobject baseIO
(
fld.name(),
baseMesh_.time().timeName(),
fld.local(),
baseMesh_,
IOobject::NO_READ,
IOobject::NO_WRITE
);
auto tfield = tmp>::New
(
baseIO,
baseMesh_,
fld.dimensions(),
internalField,
basePatchFields
);
tfield.ref().oriented() = fld.oriented();
return tfield;
}
template
Foam::tmp>
Foam::parFvFieldReconstructor::reconstructFvSurfaceField
(
const IOobject& fieldIoObject
) const
{
// Read the field
GeometricField fld
(
fieldIoObject,
procMesh_
);
return reconstructFvSurfaceField(fld);
}
template
Foam::label Foam::parFvFieldReconstructor::reconstructFvVolumeInternalFields
(
const IOobjectList& objects,
const wordRes& selectedFields
) const
{
typedef DimensionedField fieldType;
// Available fields, sorted order
const wordList fieldNames =
(
selectedFields.empty()
? objects.sortedNames()
: objects.sortedNames(selectedFields)
);
label nFields = 0;
for (const word& fieldName : fieldNames)
{
if (!nFields++)
{
Info<< " Reconstructing "
<< fieldType::typeName << "s\n" << nl;
}
Info<< " " << fieldName << nl;
tmp tfld
(
reconstructFvVolumeInternalField(*(objects[fieldName]))
);
if (isWriteProc_)
{
tfld().write();
}
}
if (nFields) Info<< endl;
return nFields;
}
template
Foam::label Foam::parFvFieldReconstructor::reconstructFvVolumeFields
(
const IOobjectList& objects,
const wordRes& selectedFields
) const
{
typedef GeometricField fieldType;
// Available fields, sorted order
const wordList fieldNames =
(
selectedFields.empty()
? objects.sortedNames()
: objects.sortedNames(selectedFields)
);
label nFields = 0;
for (const word& fieldName : fieldNames)
{
if ("cellDist" == fieldName)
{
continue;
}
if (!nFields++)
{
Info<< " Reconstructing "
<< fieldType::typeName << "s\n" << nl;
}
Info<< " " << fieldName << nl;
tmp tfld
(
reconstructFvVolumeField(*(objects[fieldName]))
);
if (isWriteProc_)
{
tfld().write();
}
}
if (nFields) Info<< endl;
return nFields;
}
template
Foam::label Foam::parFvFieldReconstructor::reconstructFvSurfaceFields
(
const IOobjectList& objects,
const wordRes& selectedFields
) const
{
typedef GeometricField fieldType;
// Available fields, sorted order
const wordList fieldNames =
(
selectedFields.empty()
? objects.sortedNames()
: objects.sortedNames(selectedFields)
);
label nFields = 0;
for (const word& fieldName : fieldNames)
{
if (!nFields++)
{
Info<< " Reconstructing "
<< fieldType::typeName << "s\n" << nl;
}
Info<< " " << fieldName << nl;
tmp tfld
(
reconstructFvSurfaceField(*(objects[fieldName]))
);
if (isWriteProc_)
{
tfld().write();
}
}
if (nFields) Info<< endl;
return nFields;
}
// ************************************************************************* //