discreteSurfaceTemplates.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) 2016-2018 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "discreteSurface.H"
29 #include "surfFields.H"
30 
31 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
32 
33 template<class Type>
35 (
36  const objectRegistry& obr,
37  const word& fieldName,
38  const word& sampleScheme
39 ) const
40 {
42  typedef DimensionedField<Type, surfGeoMesh> SurfFieldType;
43  typedef IOField<Type> TmpFieldType;
44 
45  const auto* volFldPtr = mesh().findObject<VolFieldType>(fieldName);
46 
47  if (!volFldPtr)
48  {
49  return false;
50  }
51 
52  const auto& volFld = *volFldPtr;
53  auto samplerPtr = interpolation<Type>::New(sampleScheme, *volFldPtr);
54 
55  tmp<Field<Type>> tfield = sampleOnFaces(*samplerPtr);
56 
57  // The rest could be moved into a separate helper
58  if (isA<surfMesh>(obr))
59  {
60  const surfMesh& surf = dynamicCast<const surfMesh>(obr);
61 
62  SurfFieldType* ptr = surf.getObjectPtr<SurfFieldType>(fieldName);
63  if (!ptr)
64  {
65  // Doesn't exist or the wrong type
66  ptr = new SurfFieldType
67  (
68  IOobject
69  (
70  fieldName,
71  surf.time().timeName(),
72  surf,
73  IOobject::NO_READ,
74  IOobject::NO_WRITE
75  ),
76  surf,
77  dimensioned<Type>(volFld.dimensions(), Zero)
78  );
79  ptr->writeOpt() = IOobject::NO_WRITE;
80 
81  surf.store(ptr);
82  }
83 
84  ptr->field() = tfield;
85  }
86  else
87  {
88  TmpFieldType* ptr = obr.getObjectPtr<TmpFieldType>(fieldName);
89  if (!ptr)
90  {
91  // Doesn't exist or the wrong type
92  ptr = new TmpFieldType
93  (
94  IOobject
95  (
96  fieldName,
97  obr.time().timeName(),
98  obr,
99  IOobject::NO_READ,
100  IOobject::NO_WRITE
101  ),
102  tfield().size()
103  );
104  ptr->writeOpt() = IOobject::NO_WRITE;
105 
106  obr.store(ptr);
107  }
108 
109  *ptr = tfield;
110  }
111 
112  return true;
113 }
114 
115 
116 template<class Type>
119 (
120  const word& fieldName,
121  const word& sampleScheme
122 ) const
123 {
124  typedef GeometricField<Type, fvPatchField, volMesh> VolFieldType;
125 
126  const auto* volFldPtr = mesh().findObject<VolFieldType>(fieldName);
127 
128  if (volFldPtr)
129  {
130  auto samplerPtr = interpolation<Type>::New(sampleScheme, *volFldPtr);
131  return sampleOnFaces(*samplerPtr);
132  }
133 
134  return nullptr;
135 }
136 
137 
138 template<class Type>
141 (
142  const interpolation<Type>& sampler
143 ) const
144 {
145  const labelList& elements = sampleElements_;
146 
147  const auto& vField = sampler.psi();
148  const label len = elements.size();
149 
150  // One value per face
151  auto tvalues = tmp<Field<Type>>::New(len);
152  auto& values = tvalues.ref();
153 
154  if (!onBoundary())
155  {
156  // Sample cells
157 
158  const faceList& fcs = this->surfFaces();
159  const pointField& pts = points();
160 
161  for (label i=0; i < len; ++i)
162  {
163  const label celli = elements[i];
164  const point pt = fcs[i].centre(pts);
165 
166  values[i] = sampler.interpolate(pt, celli);
167  }
168  }
169  else
170  {
171  // Sample boundary faces
172 
173  const polyBoundaryMesh& pbm = mesh().boundaryMesh();
174  const label nBnd = mesh().nBoundaryFaces();
175 
176  // Create flat boundary field
177 
178  Field<Type> bVals(nBnd, Zero);
179 
180  const auto& bField = vField.boundaryField();
181 
182  forAll(bField, patchi)
183  {
184  const label bFacei = pbm[patchi].start() - mesh().nInternalFaces();
185 
187  (
188  bVals,
189  bField[patchi].size(),
190  bFacei
191  ) = bField[patchi];
192  }
193 
194  // Sample in flat boundary field
195 
196  for (label i=0; i < len; ++i)
197  {
198  const label bFacei = (elements[i] - mesh().nInternalFaces());
199  values[i] = bVals[bFacei];
200  }
201  }
202 
203  return tvalues;
204 }
205 
206 
207 template<class Type>
210 (
211  const interpolation<Type>& interpolator
212 ) const
213 {
214  // One value per vertex
215  auto tvalues = tmp<Field<Type>>::New(sampleElements_.size());
216  auto& values = tvalues.ref();
217 
218  if (!onBoundary())
219  {
220  // Sample cells.
221 
222  forAll(sampleElements_, pointi)
223  {
224  values[pointi] = interpolator.interpolate
225  (
226  samplePoints_[pointi],
227  sampleElements_[pointi]
228  );
229  }
230  }
231  else
232  {
233  // Sample boundary faces.
234 
235  forAll(samplePoints_, pointi)
236  {
237  const label facei = sampleElements_[pointi];
238 
239  values[pointi] = interpolator.interpolate
240  (
241  samplePoints_[pointi],
242  mesh().faceOwner()[facei],
243  facei
244  );
245  }
246  }
247 
248  return tvalues;
249 }
250 
251 
252 // ************************************************************************* //
Foam::objectRegistry::getObjectPtr
Type * getObjectPtr(const word &name, const bool recursive=false) const
Definition: objectRegistryTemplates.C:423
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::discreteSurface::sampleType
bool sampleType(const objectRegistry &store, const word &fieldName, const word &sampleScheme) const
Sample the volume field onto surface,.
Definition: discreteSurfaceTemplates.C:35
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:61
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:62
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:59
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::HashTableOps::values
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:149
Foam::SubList
A List obtained as a section of another List.
Definition: SubList.H:53
Foam::objectRegistry::time
const Time & time() const
Return time.
Definition: objectRegistry.H:186
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:764
Foam::surfMesh
A surface mesh consisting of general polygon faces.
Definition: surfMesh.H:64
Foam::polyBoundaryMesh::start
label start() const
The start label of the boundary faces in the polyMesh face list.
Definition: polyBoundaryMesh.C:638
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::interpolation::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Return the field to be interpolated.
Definition: interpolation.H:127
Foam::objectRegistry
Registry of regIOobjects.
Definition: objectRegistry.H:60
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:62
Foam::Field< vector >
Foam::interpolation
Abstract base class for interpolation.
Definition: mappedPatchFieldBase.H:95
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:43
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
discreteSurface.H
Foam::regIOobject::store
void store()
Transfer ownership of this object to its registry.
Definition: regIOobjectI.H:37
Foam::Vector< scalar >
Foam::List< label >
points
const pointField & points
Definition: gmvOutputHeader.H:1
surfFields.H
Foam::discreteSurface::sampleOnPoints
tmp< Field< Type > > sampleOnPoints(const interpolation< Type > &interpolator) const
Sample field on surface points.
Foam::GeometricField< Type, fvPatchField, volMesh >
Foam::interpolation::interpolate
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.
Foam::discreteSurface::sampleOnFaces
tmp< Field< Type > > sampleOnFaces(const interpolation< Type > &sampler) const
Sample field on surface faces.
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54