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

ENH: Correct SuSp dimensions and adding function to mapping util

parent b63721f8
Branches
Tags
No related merge requests found
......@@ -165,9 +165,19 @@ SuSp
);
faMatrix<Type>& fam = tfam.ref();
fam.diag() += mesh.S()*max(sp.internalField(), scalar(0));
fam.diag() +=
mesh.S()*max
(
sp.internalField(),
dimensionedScalar("0", sp.dimensions(), Zero)
);
fam.source() -= mesh.S()*min(sp.internalField(), scalar(0))
fam.source() -=
mesh.S()*min
(
sp.internalField(),
dimensionedScalar("0", sp.dimensions(), Zero)
)
*vf.internalField();
return tfam;
......
......@@ -67,6 +67,44 @@ Foam::tmp<Foam::Field<Type>> Foam::volSurfaceMapping::mapToSurface
}
template<class Type>
Foam::tmp<Foam::Field<Type>> Foam::volSurfaceMapping::mapInternalToSurface
(
const typename GeometricField<Type, fvPatchField, volMesh>::Boundary& df
) const
{
// Grab labels for all faces in faMesh
const labelList& faceLabels = mesh_.faceLabels();
tmp<Field<Type>> tresult
(
new Field<Type>(faceLabels.size(), Zero)
);
Field<Type>& result = tresult.ref();
// Get reference to volume mesh
const polyMesh& pMesh = mesh_();
const polyBoundaryMesh& bm = pMesh.boundaryMesh();
label patchID, faceID;
// Grab droplet cloud source by identifying patch and face
forAll(faceLabels, i)
{
// Escape if face is beyond active faces, eg belongs to a face zone
if (faceLabels[i] < pMesh.nFaces())
{
patchID = bm.whichPatch(faceLabels[i]);
faceID = bm[patchID].whichFace(faceLabels[i]);
result[i] = df[patchID].patchInternalField()()[faceID];
}
}
return tresult;
}
template<class Type>
void Foam::volSurfaceMapping::mapToVolume
(
......
......@@ -96,6 +96,15 @@ public:
GeometricField<Type, fvPatchField, volMesh>::Boundary& df
) const;
//- Map patch internal field to surface
template<class Type>
tmp<Field<Type>> mapInternalToSurface
(
const typename
GeometricField<Type, fvPatchField, volMesh>::Boundary& df
) const;
//- Map surface field to volume boundary field
template<class Type>
void mapToVolume
......
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