Loading...
Searching...
No Matches
fusedLeastSquaresGrad.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2016 OpenFOAM Foundation
9 Copyright (C) 2018-2021 OpenCFD Ltd.
10 Copyright (C) 2024 M. Janssens
11-------------------------------------------------------------------------------
12License
13 This file is part of OpenFOAM.
14
15 OpenFOAM is free software: you can redistribute it and/or modify it
16 under the terms of the GNU General Public License as published by
17 the Free Software Foundation, either version 3 of the License, or
18 (at your option) any later version.
19
20 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
21 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
22 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
23 for more details.
24
25 You should have received a copy of the GNU General Public License
26 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
27
28\*---------------------------------------------------------------------------*/
29
31#include "leastSquaresVectors.H"
32#include "gaussGrad.H"
33#include "fvMesh.H"
34#include "volMesh.H"
35#include "surfaceMesh.H"
36#include "GeometricField.H"
38#include "fvcSurfaceOps.H"
39
40// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41
42template<class Type>
44<
46 <
50 >
51>
53(
55 const word& name
56) const
57{
58 typedef typename outerProduct<vector, Type>::type GradType;
60
61 const fvMesh& mesh = vf.mesh();
62
63 DebugPout<< "fusedLeastSquaresGrad<Type>::calcGrad on " << vf.name()
64 << " with name " << name << endl;
65
66 tmp<GradFieldType> tlsGrad
67 (
68 new GradFieldType
69 (
71 (
72 name,
73 vf.instance(),
74 mesh,
77 ),
78 mesh,
81 )
82 );
83 GradFieldType& lsGrad = tlsGrad.ref();
84
85 // Get reference to least square vectors
87 const surfaceVectorField& ownLs = lsv.pVectors();
88 const surfaceVectorField& neiLs = lsv.nVectors();
89
90 const auto differenceOp = [&]
91 (
92 const vector& area,
93 const Type& ownVal,
94 const Type& neiVal
95 ) -> Type
96 {
97 return neiVal - ownVal;
98 };
99
101 (
102 vf,
103 ownLs,
104 neiLs,
105 differenceOp,
106 lsGrad
107 );
108
110
111 return tlsGrad;
112}
113
114
115// ************************************************************************* //
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
Generic GeometricField class.
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
const fileName & instance() const noexcept
Read access to instance path component.
Definition IOobjectI.H:289
static FOAM_NO_DANGLING_REFERENCE const leastSquaresVectors & New(const fvMesh &mesh, Args &&... args)
Generic dimensioned Type class.
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
static const word & extrapolatedCalculatedType() noexcept
The type name for extrapolatedCalculated patch fields combines zero-gradient and calculated.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
virtual tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > calcGrad(const GeometricField< Type, fvPatchField, volMesh > &vsf, const word &name) const
Return the gradient of the given field to the gradScheme::grad for optional caching.
static void correctBoundaryConditions(const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > &)
Correct the boundary values of the gradient using the patchField snGrad functions.
Definition gaussGrad.C:199
const fvMesh & mesh() const
Return const reference to mesh.
Definition gradScheme.H:138
Least-squares gradient scheme vectors.
const surfaceVectorField & pVectors() const
Return const reference to owner least square vectors.
const surfaceVectorField & nVectors() const
Return const reference to neighbour least square vectors.
typeOfRank< typenamepTraits< arg1 >::cmptType, direction(pTraits< arg1 >::rank)+direction(pTraits< arg2 >::rank)>::type type
Definition products.H:118
A class for managing temporary objects.
Definition tmp.H:75
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
Mesh data needed to do the Finite Volume discretisation.
Definition volMesh.H:47
A class for handling words, derived from Foam::string.
Definition word.H:66
Surface integrate surfaceField creating a volField. Surface sum a surfaceField creating a volField.
#define DebugPout
Report an information message using Foam::Pout.
void surfaceOp(const GeometricField< Type, fvPatchField, volMesh > &vf, const surfaceVectorField &ownLs, const surfaceVectorField &neiLs, const CombineOp &cop, GeometricField< ResultType, fvPatchField, volMesh > &result)
GeometricField< vector, fvsPatchField, surfaceMesh > surfaceVectorField
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
Vector< scalar > vector
Definition vector.H:57