Newer
Older
/*---------------------------------------------------------------------------*\
========= |
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
\\ / O peration |
-------------------------------------------------------------------------------
-------------------------------------------------------------------------------
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 "surfaceInterpolationScheme.H"
#include "volFields.H"
#include "surfaceFields.H"
#include "geometricOneField.H"
#include "coupledFvPatchField.H"
// * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
template<class Type>
Henry Weller
committed
Foam::tmp<Foam::surfaceInterpolationScheme<Type>>
Henry Weller
committed
Foam::surfaceInterpolationScheme<Type>::New
(
const fvMesh& mesh,
Istream& schemeData
)
{
if (schemeData.eof())
{
FatalIOErrorInFunction(schemeData)
<< "Discretisation scheme not specified\n\n"
<< "Valid schemes:\n"
<< MeshConstructorTablePtr_->sortedToc()
<< exit(FatalIOError);
}
Mark Olesen
committed
const word schemeName(schemeData);
if (surfaceInterpolation::debug || surfaceInterpolationScheme<Type>::debug)
{
InfoInFunction << "Discretisation scheme = " << schemeName << endl;
auto* ctorPtr = MeshConstructorTable(schemeName);
FatalIOErrorInLookup
(
schemeData,
"discretisation",
schemeName,
*MeshConstructorTablePtr_
) << exit(FatalIOError);
return ctorPtr(mesh, schemeData);
}
template<class Type>
Henry Weller
committed
Foam::tmp<Foam::surfaceInterpolationScheme<Type>>
Henry Weller
committed
Foam::surfaceInterpolationScheme<Type>::New
(
const fvMesh& mesh,
const surfaceScalarField& faceFlux,
Istream& schemeData
)
{
if (schemeData.eof())
{
FatalIOErrorInFunction(schemeData)
<< "Discretisation scheme not specified"
<< endl << endl
<< "Valid schemes are :" << endl
<< MeshConstructorTablePtr_->sortedToc()
<< exit(FatalIOError);
}
Mark Olesen
committed
const word schemeName(schemeData);
if (surfaceInterpolation::debug || surfaceInterpolationScheme<Type>::debug)
{
InfoInFunction << "Discretisation scheme = " << schemeName << endl;
auto* ctorPtr = MeshFluxConstructorTable(schemeName);
FatalIOErrorInLookup
(
schemeData,
"discretisation",
schemeName,
*MeshFluxConstructorTablePtr_
) << exit(FatalIOError);
return ctorPtr(mesh, faceFlux, schemeData);
}
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
template<class Type>
Henry Weller
committed
Foam::tmp<Foam::GeometricField<Type, Foam::fvsPatchField, Foam::surfaceMesh>>
Henry Weller
committed
Foam::surfaceInterpolationScheme<Type>::interpolate
(
const GeometricField<Type, fvPatchField, volMesh>& vf,
const tmp<surfaceScalarField>& tlambdas,
const tmp<surfaceScalarField>& tys
)
{
if (surfaceInterpolation::debug)
{
InfoInFunction
<< "Interpolating "
<< vf.type() << " "
<< vf.name()
<< " from cells to faces without explicit correction"
<< endl;
}
const surfaceScalarField& lambdas = tlambdas();
const surfaceScalarField& ys = tys();
const Field<Type>& vfi = vf;
const scalarField& lambda = lambdas;
const scalarField& y = ys;
const fvMesh& mesh = vf.mesh();
const labelUList& P = mesh.owner();
const labelUList& N = mesh.neighbour();
Henry Weller
committed
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> tsf
(
new GeometricField<Type, fvsPatchField, surfaceMesh>
(
IOobject
(
"interpolate("+vf.name()+')',
vf.instance(),
vf.db()
),
mesh,
vf.dimensions()
)
);
GeometricField<Type, fvsPatchField, surfaceMesh>& sf = tsf.ref();
Field<Type>& sfi = sf.primitiveFieldRef();
for (label fi=0; fi<P.size(); fi++)
{
sfi[fi] = lambda[fi]*vfi[P[fi]] + y[fi]*vfi[N[fi]];
}
// Interpolate across coupled patches using given lambdas and ys
typename GeometricField<Type, fvsPatchField, surfaceMesh>::
Boundary& sfbf = sf.boundaryFieldRef();
forAll(lambdas.boundaryField(), pi)
{
const fvsPatchScalarField& pLambda = lambdas.boundaryField()[pi];
const fvsPatchScalarField& pY = ys.boundaryField()[pi];
if (vf.boundaryField()[pi].coupled())
{
sfbf[pi] =
pLambda*vf.boundaryField()[pi].patchInternalField()
+ pY*vf.boundaryField()[pi].patchNeighbourField();
}
else
{
sfbf[pi] = vf.boundaryField()[pi];
}
}
tlambdas.clear();
tys.clear();
return tsf;
}
template<class Type>
template<class SFType>
Foam::tmp
<
Foam::GeometricField
<
typename Foam::innerProduct<typename SFType::value_type, Type>::type,
Foam::fvsPatchField,
Foam::surfaceMesh
>
>
Foam::surfaceInterpolationScheme<Type>::dotInterpolate
const SFType& Sf,
const GeometricField<Type, fvPatchField, volMesh>& vf,
const tmp<surfaceScalarField>& tlambdas
)
{
if (surfaceInterpolation::debug)
{
InfoInFunction
<< "Interpolating "
<< vf.type() << " "
<< vf.name()
<< " from cells to faces without explicit correction"
typedef typename Foam::innerProduct<typename SFType::value_type, Type>::type
RetType;
const surfaceScalarField& lambdas = tlambdas();
const Field<Type>& vfi = vf;
const scalarField& lambda = lambdas;
const fvMesh& mesh = vf.mesh();
const labelUList& P = mesh.owner();
const labelUList& N = mesh.neighbour();
tmp<GeometricField<RetType, fvsPatchField, surfaceMesh>> tsf
new GeometricField<RetType, fvsPatchField, surfaceMesh>
(
IOobject
(
"interpolate("+vf.name()+')',
vf.instance(),
vf.db()
),
mesh,
Sf.dimensions()*vf.dimensions()
GeometricField<RetType, fvsPatchField, surfaceMesh>& sf = tsf.ref();
Field<RetType>& sfi = sf.primitiveFieldRef();
const typename SFType::Internal& Sfi = Sf.internalField();
for (label fi=0; fi<P.size(); fi++)
// Same as:
// sfi[fi] = Sfi[fi] & lerp(vfi[N[fi]], vfi[P[fi]], lambda[fi]);
// but maybe the compiler notices the fused multiply add form
sfi[fi] = Sfi[fi] & (lambda[fi]*(vfi[P[fi]] - vfi[N[fi]]) + vfi[N[fi]]);
}
// Interpolate across coupled patches using given lambdas
typename GeometricField<RetType, fvsPatchField, surfaceMesh>::
Boundary& sfbf = sf.boundaryFieldRef();
forAll(lambdas.boundaryField(), pi)
{
const fvsPatchScalarField& pLambda = lambdas.boundaryField()[pi];
Henry Weller
committed
const typename SFType::Patch& pSf = Sf.boundaryField()[pi];
fvsPatchField<RetType>& psf = sfbf[pi];
if (vf.boundaryField()[pi].coupled())
{
psf =
pSf
& lerp
(
vf.boundaryField()[pi].patchNeighbourField(),
vf.boundaryField()[pi].patchInternalField(),
pLambda
psf = pSf & vf.boundaryField()[pi];
}
}
tlambdas.clear();
// tsf.ref().oriented() = Sf.oriented();
return tsf;
}
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
template<class Type>
Foam::tmp<Foam::GeometricField<Type, Foam::fvsPatchField, Foam::surfaceMesh>>
Foam::surfaceInterpolationScheme<Type>::interpolate
(
const GeometricField<Type, fvPatchField, volMesh>& vf,
const tmp<surfaceScalarField>& tlambdas
)
{
return dotInterpolate(geometricOneField(), vf, tlambdas);
}
template<class Type>
Foam::tmp
<
Foam::GeometricField
<
typename Foam::innerProduct<Foam::vector, Type>::type,
Foam::fvsPatchField,
Foam::surfaceMesh
>
>
Foam::surfaceInterpolationScheme<Type>::dotInterpolate
(
const surfaceVectorField& Sf,
const GeometricField<Type, fvPatchField, volMesh>& vf
) const
{
if (surfaceInterpolation::debug)
{
InfoInFunction
<< "Interpolating "
<< vf.type() << " "
<< vf.name()
<< " from cells to faces"
<< endl;
}
tmp
<
GeometricField
<
typename Foam::innerProduct<Foam::vector, Type>::type,
fvsPatchField,
surfaceMesh
>
> tsf = dotInterpolate(Sf, vf, weights(vf));
tsf.ref().oriented() = Sf.oriented();
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
if (corrected())
{
tsf.ref() += Sf & correction(vf);
}
return tsf;
}
template<class Type>
Foam::tmp
<
Foam::GeometricField
<
typename Foam::innerProduct<Foam::vector, Type>::type,
Foam::fvsPatchField,
Foam::surfaceMesh
>
>
Foam::surfaceInterpolationScheme<Type>::dotInterpolate
(
const surfaceVectorField& Sf,
const tmp<GeometricField<Type, fvPatchField, volMesh>>& tvf
) const
{
tmp
<
GeometricField
<
typename Foam::innerProduct<Foam::vector, Type>::type,
fvsPatchField,
surfaceMesh
>
> tSfDotinterpVf = dotInterpolate(Sf, tvf());
tvf.clear();
return tSfDotinterpVf;
}
Henry Weller
committed
Foam::tmp<Foam::GeometricField<Type, Foam::fvsPatchField, Foam::surfaceMesh>>
Henry Weller
committed
Foam::surfaceInterpolationScheme<Type>::interpolate
(
const GeometricField<Type, fvPatchField, volMesh>& vf
) const
{
if (surfaceInterpolation::debug)
{
InfoInFunction
<< "Interpolating "
<< vf.type() << " "
<< vf.name()
<< " from cells to faces"
Henry Weller
committed
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> tsf
= interpolate(vf, weights(vf));
if (corrected())
{
tsf.ref() += correction(vf);
}
return tsf;
}
template<class Type>
Henry Weller
committed
Foam::tmp<Foam::GeometricField<Type, Foam::fvsPatchField, Foam::surfaceMesh>>
Henry Weller
committed
Foam::surfaceInterpolationScheme<Type>::interpolate
Henry Weller
committed
const tmp<GeometricField<Type, fvPatchField, volMesh>>& tvf
Henry Weller
committed
tmp<GeometricField<Type, fvsPatchField, surfaceMesh>> tinterpVf
= interpolate(tvf());
tvf.clear();
return tinterpVf;
}
// ************************************************************************* //