sampledPatchInternalFieldTemplates.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 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
32
33// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
34
35template<class Type>
37Foam::sampledPatchInternalField::sampleOnFaces
38(
39 const interpolation<Type>& sampler
40) const
41{
42 const auto& vField = sampler.psi();
43
44 // One value per face
45 auto tvalues = tmp<Field<Type>>::New(patchFaceLabels().size());
46 auto& values = tvalues.ref();
47
48 forAll(patchStart(), i)
49 {
50 // Get patchface wise data by sampling internal field
51 Field<Type> interpVals = vField.primitiveField();
52 mappers_[i].map().distribute(interpVals);
53
54 // Store at correct position in values
55 label end =
56 (
57 i < patchStart().size()-1
58 ? patchStart()[i+1]
60 );
61
62 for (label triI = patchStart()[i]; triI < end; triI++)
63 {
64 values[triI] = interpVals[patchFaceLabels()[triI]];
65 }
66 }
67
68 return tvalues;
69}
70
71
72template<class Type>
74Foam::sampledPatchInternalField::sampleOnPoints
75(
76 const interpolation<Type>& interpolator
77) const
78{
79 label sz = 0;
80 forAll(patchIDs(), i)
81 {
82 sz += mesh().boundaryMesh()[patchIDs()[i]].size();
83 }
84
85 Field<Type> allPatchVals(sz);
86 sz = 0;
87
88 forAll(patchIDs(), i)
89 {
90 // See mappedFixedValueFvPatchField
91 const mapDistribute& distMap = mappers_[i].map();
92
93 // Send back sample points to processor that holds the cell.
94 // Mark cells with point::max so we know which ones we need
95 // to interpolate (since expensive).
96 vectorField samples(mappers_[i].samplePoints());
97 distMap.reverseDistribute(mesh().nCells(), point::max, samples);
98
99 Field<Type> patchVals(mesh().nCells());
100
101 forAll(samples, celli)
102 {
103 if (samples[celli] != point::max)
104 {
105 patchVals[celli] = interpolator.interpolate
106 (
107 samples[celli],
108 celli
109 );
110 }
111 }
112
113 distMap.distribute(patchVals);
114
115 // Now patchVals holds the interpolated data in patch face order.
116 // Collect.
117 SubList<Type>(allPatchVals, patchVals.size(), sz) = patchVals;
118 sz += patchVals.size();
119 }
120
121 // Interpolate to points. Reconstruct the patch of all faces to aid
122 // interpolation.
123
124 labelList meshFaceLabels(allPatchVals.size());
125 sz = 0;
126 for (const label patchId : patchIDs())
127 {
128 const polyPatch& pp = mesh().boundaryMesh()[patchId];
129 forAll(pp, i)
130 {
131 meshFaceLabels[sz++] = pp.start()+i;
132 }
133 }
134
135 indirectPrimitivePatch allPatches
136 (
137 IndirectList<face>(mesh().faces(), meshFaceLabels),
138 mesh().points()
139 );
140
142 (
143 allPatches
144 ).faceToPointInterpolate(allPatchVals);
145}
146
147
148// ************************************************************************* //
Generic templated field type.
Definition: Field.H:82
A List with indirect addressing.
Definition: IndirectList.H:119
label size() const
The surface size is the number of faces.
Interpolation class within a primitive patch. Allows interpolation from points to faces and vice vers...
tmp< Field< Type > > faceToPointInterpolate(const Field< Type > &ff) const
Interpolate from faces to points.
A list of faces which address into the list of points.
A List obtained as a section of another List.
Definition: SubList.H:70
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
Abstract base class for volume field interpolation.
Definition: interpolation.H:60
const GeometricField< Type, fvPatchField, volMesh > & psi() const noexcept
Return the field to be interpolated.
virtual Type interpolate(const vector &position, const label celli, const label facei=-1) const =0
Interpolate field to the given point in the given cell.
Class containing processor-to-processor mapping information.
void reverseDistribute(const label constructSize, List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Reverse distribute data using default commsType.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute data using default commsType.
static const complex max
complex (VGREAT,VGREAT)
Definition: complex.H:293
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:456
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:75
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:364
const labelList & patchStart() const
The offset into patchIndex, patchFaceLabels.
Definition: sampledPatch.H:168
const labelList & patchFaceLabels() const
For each face, the patch local face ID.
Definition: sampledPatch.H:180
A class for managing temporary objects.
Definition: tmp.H:65
dynamicFvMesh & mesh
const pointField & points
label patchId(-1)
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
scalarField samples(nIntervals, Zero)