patchCorrectedInterpolationTemplates.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 | Copyright (C) 2015 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
27 #include "PointData.H"
28 #include "PointEdgeWave.H"
29 #include "volPointInterpolation.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 template <class Type>
36 (
37  const GeometricField<Type, fvPatchField, volMesh>& cellDisplacement,
39 ) const
40 {
41  // Create an uncorrected field
43  pointUncorrectedDisplacement
44  (
45  IOobject
46  (
47  "pointUncorrectedDisplacement",
48  mesh().time().timeName(),
49  mesh()
50  ),
51  pointDisplacement.mesh(),
52  pointDisplacement.dimensions(),
54  );
55 
56  // Interpolate to the uncorrected field, overwriting the fixed value
57  // boundary conditions
58  pointUncorrectedDisplacement ==
59  volPointInterpolation::New(mesh()).interpolate
60  (
61  cellDisplacement,
62  wordList
63  (
64  pointUncorrectedDisplacement.boundaryField().size(),
66  )
67  );
68 
69  // Set the point displacement to the uncorrected result everywhere except
70  // for on the boundary
71  pointDisplacement.internalField() =
72  pointUncorrectedDisplacement.internalField();
73  pointDisplacement.correctBoundaryConditions();
74 
75  // Set the residual displacement as the difference between the boundary
76  // specification and the uncorrected solution
77  // (this reuses the uncorrected displacement field as the residual)
78  pointUncorrectedDisplacement ==
79  pointDisplacement - pointUncorrectedDisplacement;
80 
81  // Interpolate the residual from the boundary into the field
82  interpolateDataFromPatchGroups(pointUncorrectedDisplacement);
83 
84  // Add the residual to the point displacement and correct the boundary
85  pointDisplacement += pointUncorrectedDisplacement;
86  pointDisplacement.correctBoundaryConditions();
87 }
88 
89 
90 template <class Type>
92 (
94 ) const
95 {
96  // Initialise
97  pointScalarField weight
98  (
99  IOobject
100  (
101  "weight",
102  mesh().time().timeName(),
103  mesh()
104  ),
105  data.mesh(),
106  dimensionedScalar("zero", dimless, 0),
108  );
109  data = dimensioned<Type>("zero", data.dimensions(), pTraits<Type>::zero);
110 
111  forAll(patchGroups_, patchGroupI)
112  {
113  // Distance and data for the current group
114  pointScalarField patchDistance
115  (
116  IOobject
117  (
118  "patchDistance",
119  mesh().time().timeName(),
120  mesh()
121  ),
122  data.mesh(),
123  dimensionedScalar("zero", data.dimensions(), 0),
125  );
127 
128  // Wave the data through the mesh
129  propagateDataFromPatchGroup
130  (
131  patchGroupI,
132  patchDistance,
133  patchData
134  );
135 
136  // Calculate the weight and add to weighted sum
137  const scalarField patchWeight
138  (
139  1/max(sqr(patchDistance.internalField()), SMALL)
140  );
141  data.internalField() += patchWeight*patchData.internalField();
142  weight.internalField() += patchWeight;
143  }
144 
145  // Complete the average
146  data /= weight;
147 }
148 
149 
150 template <class Type>
152 (
153  const label patchGroupI,
156 ) const
157 {
158  const labelList& patchGroup(patchGroups_[patchGroupI]);
159 
160  // Get the size of the seed info
161  label nSeedInfo(0);
162  forAll(patchGroup, patchGroupI)
163  {
164  const label patchI(patchGroup[patchGroupI]);
165 
166  nSeedInfo += data.boundaryField()[patchI].size();
167  }
168 
169  // Generate the seed labels and info
170  labelList seedLabels(nSeedInfo);
171  List<PointData<Type> > seedInfo(nSeedInfo);
172  nSeedInfo = 0;
173  forAll(patchGroup, patchGroupI)
174  {
175  const label patchI(patchGroup[patchGroupI]);
176 
177  pointPatchField<Type>& patchDataField(data.boundaryField()[patchI]);
178 
179  patchDataField.updateCoeffs();
180 
181  const pointPatch& patch(patchDataField.patch());
182  const Field<Type> patchData(patchDataField.patchInternalField());
183 
184  forAll(patch.meshPoints(), patchPointI)
185  {
186  const label pointI(patch.meshPoints()[patchPointI]);
187 
188  seedLabels[nSeedInfo] = pointI;
189 
190  seedInfo[nSeedInfo] =
192  (
193  mesh().points()[pointI],
194  0,
195  patchData[patchPointI]
196  );
197 
198  nSeedInfo ++;
199  }
200  }
201 
202  // Wave the data through the mesh
203  List<PointData<Type> > allPointInfo(mesh().nPoints());
204  List<PointData<Type> > allEdgeInfo(mesh().nEdges());
206  (
207  mesh(),
208  seedLabels,
209  seedInfo,
210  allPointInfo,
211  allEdgeInfo,
212  mesh().globalData().nTotalPoints()
213  );
214 
215  // Copy result into the fields
216  forAll(allPointInfo, pointI)
217  {
218  distance[pointI] = sqrt(allPointInfo[pointI].distSqr());
219  data[pointI] = allPointInfo[pointI].data();
220  }
221 }
222 
223 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::compressible::New
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
Definition: turbulentFluidThermoModel.C:36
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::patchCorrectedInterpolation::interpolateDataFromPatchGroups
void interpolateDataFromPatchGroups(GeometricField< Type, pointPatchField, pointMesh > &) const
Interpolate patch data by inverse distance weighting.
Definition: patchCorrectedInterpolationTemplates.C:92
Foam::PointEdgeWave
Wave propagation of information through grid. Every iteration information goes through one layer of e...
Definition: PointEdgeWave.H:85
Foam::patchCorrectedInterpolation::interpolateType
void interpolateType(const GeometricField< Type, fvPatchField, volMesh > &, GeometricField< Type, pointPatchField, pointMesh > &) const
Interpolate the given cell displacement.
Definition: patchCorrectedInterpolationTemplates.C:36
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::pointPatchField::patchInternalField
tmp< Field< Type > > patchInternalField() const
Return field created from appropriate internal field values.
Foam::pointPatch
Basic pointPatch represents a set of points from the mesh.
Definition: pointPatch.H:56
Foam::PointData
Variant of pointEdgePoint with some transported additional data. Templated on the transported data ty...
Definition: medialAxisMeshMover.H:57
Foam::pointPatchField< Type >
Foam::pointPatchField::patch
const pointPatch & patch() const
Return patch.
Definition: pointPatchField.H:271
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::fixedValuePointPatchField
A FixedValue boundary condition for pointField.
Definition: fixedValuePointPatchField.H:49
Foam::label
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
PointEdgeWave.H
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
fixedValuePointPatchField.H
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::PointEdgeWave::data
const TrackingData & data() const
Additional data to be passed into container.
Definition: PointEdgeWave.H:283
Foam::distance
scalar distance(const vector &p1, const vector &p2)
Definition: curveTools.C:12
Foam::GeometricField::correctBoundaryConditions
void correctBoundaryConditions()
Correct boundary field.
Definition: GeometricField.C:885
volPointInterpolation.H
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
PointData.H
Foam::pointPatch::meshPoints
virtual const labelList & meshPoints() const =0
Return mesh points.
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::patchCorrectedInterpolation::propagateDataFromPatchGroup
void propagateDataFromPatchGroup(const label, pointScalarField &, GeometricField< Type, pointPatchField, pointMesh > &) const
Propagate data from a number of patches into the field.
Definition: patchCorrectedInterpolationTemplates.C:152
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::pTraits
Traits class for primitives.
Definition: pTraits.H:50
points
const pointField & points
Definition: gmvOutputHeader.H:1
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::pointPatchField::updateCoeffs
virtual void updateCoeffs()
Update the coefficients associated with the patch field.
Definition: pointPatchField.H:409
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
zeroGradientPointPatchField.H
Foam::zeroGradientPointPatchField
Foam::zeroGradientPointPatchField.
Definition: zeroGradientPointPatchField.H:49
Foam::data
Database for solution data, solver performance and other reduced data.
Definition: data.H:52