Loading...
Searching...
No Matches
pointVolInterpolate.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) Wikki Ltd
9 Copyright (C) 2019 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
30#include "volFields.H"
31#include "pointFields.H"
33#include "emptyFvPatch.H"
34#include "globalMeshData.H"
35
36// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
37
38template<class Type>
40(
43) const
44{
46 << "interpolating field from points to cells"
47 << endl;
48
49 const FieldField<Field, scalar>& weights = volWeights();
50 const labelListList& cellPoints = vf.mesh().cellPoints();
51
52 // Multiply pointField by weighting factor matrix to create volField
53 forAll(cellPoints, cellI)
54 {
55 vf[cellI] = pTraits<Type>::zero;
56
57 const labelList& curCellPoints = cellPoints[cellI];
58
59 forAll(curCellPoints, cellPointI)
60 {
61 vf[cellI] +=
62 weights[cellI][cellPointI]*pf[curCellPoints[cellPointI]];
63 }
64 }
65
66
67 // Interpolate patch values: over-ride the internal values for the points
68 // on the patch with the interpolated point values from the faces
69 const fvBoundaryMesh& bm = vMesh().boundary();
70
71 const PtrList<primitivePatchInterpolation>& pi = patchInterpolators();
72 forAll(bm, patchI)
73 {
74 // If the patch is empty, skip it
75 if (bm[patchI].type() != emptyFvPatch::typeName)
76 {
77 vf.boundaryFieldRef()[patchI] =
78 pi[patchI].pointToFaceInterpolate
79 (
80 pf.boundaryField()[patchI].patchInternalField()
81 );
82 }
83 }
84
86
88 << "finished interpolating field from points to cells"
89 << endl;
90}
91
92
93template<class Type>
96(
98) const
99{
100 // Construct tmp<pointField>
102 (
104 (
106 (
107 "pointVolInterpolate(" + pf.name() + ')',
108 pf.instance(),
109 pf.db()
110 ),
111 vMesh(),
112 pf.dimensions()
113 )
114 );
115
116 // Perform interpolation
117 interpolate(pf, tvf.ref());
129{
130 // Construct tmp<volField>
132 interpolate(tpf());
133 tpf.clear();
134 return tvf;
135}
136
137
138// ************************************************************************* //
constexpr scalar pi(M_PI)
const Mesh & mesh() const noexcept
Return const reference to mesh.
const dimensionSet & dimensions() const noexcept
Return dimensions.
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
Generic GeometricField class.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
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 objectRegistry & db() const noexcept
Return the local objectRegistry.
Definition IOobject.C:450
const fileName & instance() const noexcept
Read access to instance path component.
Definition IOobjectI.H:289
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers....
Definition PtrList.H:67
A fvBoundaryMesh is a fvPatch list with a reference to the associated fvMesh, with additional search ...
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
const FieldField< Field, scalar > & volWeights() const
Return reference to weights arrays.
const fvMesh & vMesh() const noexcept
void interpolate(const GeometricField< Type, pointPatchField, pointMesh > &, GeometricField< Type, fvPatchField, volMesh > &) const
Interpolate from pointField to volField.
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
#define DebugInFunction
Report an information message using Foam::Info.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition curveTools.C:75
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299