From 256dafeea511159e472ff5a755dec15cb79dff9d Mon Sep 17 00:00:00 2001 From: mattijs <mattijs> Date: Wed, 5 Mar 2025 09:11:47 +0000 Subject: [PATCH] ENH: fvMatrix: adding setValues for patch Speedup in offloading --- .../epsilonWallFunctionFvPatchScalarField.C | 2 +- .../omegaWallFunctionFvPatchScalarField.C | 2 +- .../fvMatrices/fvMatrix/fvMatrix.C | 107 +++++++++++++++++- .../fvMatrices/fvMatrix/fvMatrix.H | 11 +- 4 files changed, 118 insertions(+), 4 deletions(-) diff --git a/src/TurbulenceModels/turbulenceModels/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.C b/src/TurbulenceModels/turbulenceModels/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.C index 8c78d596140..9f9df980b13 100644 --- a/src/TurbulenceModels/turbulenceModels/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.C +++ b/src/TurbulenceModels/turbulenceModels/derivedFvPatchFields/wallFunctions/epsilonWallFunctions/epsilonWallFunction/epsilonWallFunctionFvPatchScalarField.C @@ -613,7 +613,7 @@ void Foam::epsilonWallFunctionFvPatchScalarField::manipulateMatrix return; } - matrix.setValues(patch().faceCells(), patchInternalField()); + matrix.setValues(patch().index(), patchInternalField()()); fvPatchField<scalar>::manipulateMatrix(matrix); } diff --git a/src/TurbulenceModels/turbulenceModels/derivedFvPatchFields/wallFunctions/omegaWallFunctions/omegaWallFunction/omegaWallFunctionFvPatchScalarField.C b/src/TurbulenceModels/turbulenceModels/derivedFvPatchFields/wallFunctions/omegaWallFunctions/omegaWallFunction/omegaWallFunctionFvPatchScalarField.C index 2a1770659a0..8b5a77c0890 100644 --- a/src/TurbulenceModels/turbulenceModels/derivedFvPatchFields/wallFunctions/omegaWallFunctions/omegaWallFunction/omegaWallFunctionFvPatchScalarField.C +++ b/src/TurbulenceModels/turbulenceModels/derivedFvPatchFields/wallFunctions/omegaWallFunctions/omegaWallFunction/omegaWallFunctionFvPatchScalarField.C @@ -603,7 +603,7 @@ void Foam::omegaWallFunctionFvPatchScalarField::manipulateMatrix return; } - matrix.setValues(patch().faceCells(), patchInternalField()); + matrix.setValues(patch().index(), patchInternalField()); fvPatchField<scalar>::manipulateMatrix(matrix); } diff --git a/src/finiteVolume/fvMatrices/fvMatrix/fvMatrix.C b/src/finiteVolume/fvMatrices/fvMatrix/fvMatrix.C index b2b64f0d13f..0aaeded6572 100644 --- a/src/finiteVolume/fvMatrices/fvMatrix/fvMatrix.C +++ b/src/finiteVolume/fvMatrices/fvMatrix/fvMatrix.C @@ -6,7 +6,7 @@ \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2011-2017 OpenFOAM Foundation - Copyright (C) 2016-2024 OpenCFD Ltd. + Copyright (C) 2016-2025 OpenCFD Ltd. ------------------------------------------------------------------------------- License This file is part of OpenFOAM. @@ -1006,6 +1006,111 @@ void Foam::fvMatrix<Type>::setValues } +template<class Type> +void Foam::fvMatrix<Type>::setValues +( + const label setPatchi, + const UList<Type>& values +) +{ +// this->setValuesFromList(cellLabels, values); + const fvMesh& mesh = psi_.mesh(); + const auto& cellLabels = mesh.boundaryMesh()[setPatchi].faceCells(); + + const cellList& cells = mesh.cells(); + const labelUList& own = mesh.owner(); + const labelUList& nei = mesh.neighbour(); + + scalarField& Diag = diag(); + Field<Type>& psi = + const_cast + < + GeometricField<Type, fvPatchField, volMesh>& + >(psi_).primitiveFieldRef(); + + + // Following actions: + // - adjust local field psi + // - set local matrix to be diagonal (so adjust source) + // - cut connections to neighbours + // - make (on non-adjusted cells) contribution explicit + + if (symmetric() || asymmetric()) + { + forAll(cellLabels, i) + { + const label celli = cellLabels[i]; + const Type& value = values[i]; + + for (const label facei : cells[celli]) + { + const label patchi = mesh.boundaryMesh().patchID(facei); + + if (patchi == -1) + { + if (symmetric()) + { + if (celli == own[facei]) + { + source_[nei[facei]] -= upper()[facei]*value; + } + else + { + source_[own[facei]] -= upper()[facei]*value; + } + + upper()[facei] = 0.0; + } + else + { + if (celli == own[facei]) + { + source_[nei[facei]] -= lower()[facei]*value; + } + else + { + source_[own[facei]] -= upper()[facei]*value; + } + + upper()[facei] = 0.0; + lower()[facei] = 0.0; + } + } + else if (patchi != setPatchi) + { + if (internalCoeffs_[patchi].size()) + { + const label patchFacei = + mesh.boundaryMesh()[patchi].whichFace(facei); + + internalCoeffs_[patchi][patchFacei] = Zero; + boundaryCoeffs_[patchi][patchFacei] = Zero; + } + } + } + } + + // Can explictly zero out contribution from setPatchi + if (internalCoeffs_[setPatchi].size()) + { + internalCoeffs_[setPatchi] = Zero; + boundaryCoeffs_[setPatchi] = Zero; + } + } + + // Note: above loop might have affected source terms on adjusted cells + // so make sure to adjust them afterwards + forAll(cellLabels, i) + { + const label celli = cellLabels[i]; + const Type& value = values[i]; + + psi[celli] = value; + source_[celli] = value*Diag[celli]; + } +} + + template<class Type> void Foam::fvMatrix<Type>::setReference ( diff --git a/src/finiteVolume/fvMatrices/fvMatrix/fvMatrix.H b/src/finiteVolume/fvMatrices/fvMatrix/fvMatrix.H index 78f18d22c3e..68405a70809 100644 --- a/src/finiteVolume/fvMatrices/fvMatrix/fvMatrix.H +++ b/src/finiteVolume/fvMatrices/fvMatrix/fvMatrix.H @@ -6,7 +6,7 @@ \\/ M anipulation | ------------------------------------------------------------------------------- Copyright (C) 2011-2017 OpenFOAM Foundation - Copyright (C) 2016-2023 OpenCFD Ltd. + Copyright (C) 2016-2025 OpenCFD Ltd. ------------------------------------------------------------------------------- License This file is part of OpenFOAM. @@ -541,6 +541,15 @@ public: const UIndirectList<Type>& values ); + //- Set solution in cells next to the given patch to the specified + //- values and eliminate the corresponding equations from the + //- matrix. + void setValues + ( + const label patchi, + const UList<Type>& values + ); + //- Set reference level for solution void setReference ( -- GitLab