Skip to content
Snippets Groups Projects
Commit ce5e2e7e authored by sergio's avatar sergio Committed by Sergio Ferraris
Browse files

ENH: volSurfaceMapping: new mapToField and mapToSurface funcs

parent db9460f0
1 merge request!475finiteArea: new models and solvers for films
...@@ -6,7 +6,7 @@ ...@@ -6,7 +6,7 @@
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
Copyright (C) 2016-2017 Wikki Ltd Copyright (C) 2016-2017 Wikki Ltd
Copyright (C) 2020 OpenCFD Ltd. Copyright (C) 2020-2021 OpenCFD Ltd.
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
This file is part of OpenFOAM. This file is part of OpenFOAM.
...@@ -65,6 +65,36 @@ Foam::tmp<Foam::Field<Type>> Foam::volSurfaceMapping::mapToSurface ...@@ -65,6 +65,36 @@ Foam::tmp<Foam::Field<Type>> Foam::volSurfaceMapping::mapToSurface
} }
template<class Type>
Foam::tmp<Foam::Field<Type>> Foam::volSurfaceMapping::mapToSurface
(
const Field<Type>& f
) const
{
const labelList& faceLabels = mesh_.faceLabels();
auto tresult = tmp<Field<Type>>::New(faceLabels.size(), Zero);
auto& result = tresult.ref();
const polyMesh& pMesh = mesh_();
const polyBoundaryMesh& bm = pMesh.boundaryMesh();
label patchID, faceID;
forAll(faceLabels, i)
{
if (faceLabels[i] < pMesh.nFaces())
{
patchID = bm.whichPatch(faceLabels[i]);
faceID = bm[patchID].whichFace(faceLabels[i]);
result[i] = f[faceID];
}
}
return tresult;
}
template<class Type> template<class Type>
Foam::tmp<Foam::Field<Type>> Foam::volSurfaceMapping::mapInternalToSurface Foam::tmp<Foam::Field<Type>> Foam::volSurfaceMapping::mapInternalToSurface
( (
...@@ -151,6 +181,19 @@ void Foam::volSurfaceMapping::mapToField ...@@ -151,6 +181,19 @@ void Foam::volSurfaceMapping::mapToField
const GeometricField<Type, faPatchField, areaMesh>& af, const GeometricField<Type, faPatchField, areaMesh>& af,
Field<Type>& f Field<Type>& f
) const ) const
{
const Field<Type>& afi = af.internalField();
mapToField(afi, f);
}
template<class Type>
void Foam::volSurfaceMapping::mapToField
(
const Field<Type>& af,
Field<Type>& f
) const
{ {
const labelList& faceLabels = mesh_.faceLabels(); const labelList& faceLabels = mesh_.faceLabels();
...@@ -158,15 +201,13 @@ void Foam::volSurfaceMapping::mapToField ...@@ -158,15 +201,13 @@ void Foam::volSurfaceMapping::mapToField
const polyBoundaryMesh& bm = pMesh.boundaryMesh(); const polyBoundaryMesh& bm = pMesh.boundaryMesh();
label patchID, faceID; label patchID, faceID;
const Field<Type>& afi = af.internalField();
forAll(faceLabels, i) forAll(faceLabels, i)
{ {
if (faceLabels[i] < pMesh.nFaces()) if (faceLabels[i] < pMesh.nFaces())
{ {
patchID = bm.whichPatch(faceLabels[i]); patchID = bm.whichPatch(faceLabels[i]);
faceID = bm[patchID].whichFace(faceLabels[i]); faceID = bm[patchID].whichFace(faceLabels[i]);
f[faceID] = afi[i]; f[faceID] = af[i];
} }
} }
} }
......
...@@ -6,7 +6,7 @@ ...@@ -6,7 +6,7 @@
\\/ M anipulation | \\/ M anipulation |
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
Copyright (C) 2016-2017 Wikki Ltd Copyright (C) 2016-2017 Wikki Ltd
Copyright (C) 2019-2020 OpenCFD Ltd. Copyright (C) 2019-2021 OpenCFD Ltd.
------------------------------------------------------------------------------- -------------------------------------------------------------------------------
License License
This file is part of OpenFOAM. This file is part of OpenFOAM.
...@@ -42,6 +42,7 @@ SourceFiles ...@@ -42,6 +42,7 @@ SourceFiles
#define volSurfaceMapping_H #define volSurfaceMapping_H
#include "faMesh.H" #include "faMesh.H"
#include "volMesh.H"
// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * // // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
...@@ -85,7 +86,7 @@ public: ...@@ -85,7 +86,7 @@ public:
// Member Functions // Member Functions
//- Map droplet cloud sources to surface //- Map Boundary field to surface
template<class Type> template<class Type>
tmp<Field<Type>> mapToSurface tmp<Field<Type>> mapToSurface
( (
...@@ -93,6 +94,9 @@ public: ...@@ -93,6 +94,9 @@ public:
GeometricField<Type, fvPatchField, volMesh>::Boundary& df GeometricField<Type, fvPatchField, volMesh>::Boundary& df
) const; ) const;
//- Map vol Field to surface Field
template<class Type>
tmp<Field<Type>> mapToSurface(const Field<Type>& f) const;
//- Map patch internal field to surface //- Map patch internal field to surface
template<class Type> template<class Type>
...@@ -118,14 +122,23 @@ public: ...@@ -118,14 +122,23 @@ public:
typename GeometricField<Type, fvPatchField, volMesh>::Boundary& bf typename GeometricField<Type, fvPatchField, volMesh>::Boundary& bf
) const; ) const;
//- Map surface field to field assumes Field //- Map surface field to field
//- faces in the same order as Boundary // Assumes Field faces in the same order as Boundary
template<class Type> template<class Type>
void mapToField void mapToField
( (
const GeometricField<Type, faPatchField, areaMesh>& af, const GeometricField<Type, faPatchField, areaMesh>& af,
Field<Type>& f Field<Type>& f
) const; ) const;
//- Map surface field to volume field
// Assumes Field faces in the same order as Boundary
template<class Type>
void mapToField
(
const Field<Type>& af,
Field<Type>& f
) const;
}; };
......
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment