From 97b65b811d881d32efdd87589563b2c63d44364b Mon Sep 17 00:00:00 2001 From: Henry Date: Mon, 11 Feb 2013 17:02:41 +0000 Subject: [PATCH] linearUpwindV: calculate interpolated field on processor patches --- .../schemes/linearUpwind/linearUpwindV.C | 89 +++++++++++++++++-- 1 file changed, 84 insertions(+), 5 deletions(-) diff --git a/src/finiteVolume/interpolation/surfaceInterpolation/schemes/linearUpwind/linearUpwindV.C b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/linearUpwind/linearUpwindV.C index 9e5ad7f3a0..6368d2accd 100644 --- a/src/finiteVolume/interpolation/surfaceInterpolation/schemes/linearUpwind/linearUpwindV.C +++ b/src/finiteVolume/interpolation/surfaceInterpolation/schemes/linearUpwind/linearUpwindV.C @@ -2,7 +2,7 @@ ========= | \\ / F ield | OpenFOAM: The Open Source CFD Toolbox \\ / O peration | - \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation + \\ / A nd | Copyright (C) 2011-2013 OpenFOAM Foundation \\/ M anipulation | ------------------------------------------------------------------------------- License @@ -70,8 +70,8 @@ Foam::linearUpwindV::correction const labelList& own = mesh.owner(); const labelList& nei = mesh.neighbour(); - const vectorField& C = mesh.C(); - const vectorField& Cf = mesh.Cf(); + const volVectorField& C = mesh.C(); + const surfaceVectorField& Cf = mesh.Cf(); tmp < @@ -97,8 +97,7 @@ Foam::linearUpwindV::correction if (faceFlux[facei] > 0.0) { maxCorr = - (1.0 - w[facei]) - *(vf[nei[facei]] - vf[own[facei]]); + (1.0 - w[facei])*(vf[nei[facei]] - vf[own[facei]]); sfCorr[facei] = (Cf[facei] - C[own[facei]]) & gradVf[own[facei]]; @@ -139,6 +138,86 @@ Foam::linearUpwindV::correction } } + + typename GeometricField:: + GeometricBoundaryField& bSfCorr = sfCorr.boundaryField(); + + forAll(bSfCorr, patchi) + { + fvsPatchField& pSfCorr = bSfCorr[patchi]; + + if (pSfCorr.coupled()) + { + const labelUList& pOwner = + mesh.boundary()[patchi].faceCells(); + + const vectorField& pCf = Cf.boundaryField()[patchi]; + const scalarField& pW = w.boundaryField()[patchi]; + + const scalarField& pFaceFlux = faceFlux.boundaryField()[patchi]; + + const Field::type> pGradVfNei + ( + gradVf.boundaryField()[patchi].patchNeighbourField() + ); + + const Field pVfNei + ( + vf.boundaryField()[patchi].patchNeighbourField() + ); + + // Build the d-vectors + vectorField pd(Cf.boundaryField()[patchi].patch().delta()); + + forAll(pOwner, facei) + { + label own = pOwner[facei]; + + vector maxCorr; + + if (pFaceFlux[facei] > 0) + { + pSfCorr[facei] = (pCf[facei] - C[own]) & gradVf[own]; + + maxCorr = (1.0 - pW[facei])*(pVfNei[facei] - vf[own]); + } + else + { + pSfCorr[facei] = + (pCf[facei] - pd[facei] - C[own]) & pGradVfNei[facei]; + + maxCorr = pW[facei]*(vf[own] - pVfNei[facei]); + } + + scalar pSfCorrs = magSqr(pSfCorr[facei]); + scalar maxCorrs = pSfCorr[facei] & maxCorr; + + if (pSfCorrs > 0) + { + if (maxCorrs < 0) + { + pSfCorr[facei] = vector::zero; + } + else if (pSfCorrs > maxCorrs) + { + pSfCorr[facei] *= maxCorrs/(pSfCorrs + VSMALL); + } + } + else if (pSfCorrs < 0) + { + if (maxCorrs > 0) + { + pSfCorr[facei] = vector::zero; + } + else if (pSfCorrs < maxCorrs) + { + pSfCorr[facei] *= maxCorrs/(pSfCorrs - VSMALL); + } + } + } + } + } + return tsfCorr; }