fvMeshSubsetTemplates.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) 2019-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
29 #include "fvMeshSubset.H"
30 #include "emptyFvsPatchField.H"
31 #include "emptyPointPatchField.H"
32 #include "emptyFvPatchFields.H"
35 #include "flipOp.H"
36 
37 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
38 
39 template<class Type>
41 <
43 >
45 (
46  const GeometricField<Type, fvPatchField, volMesh>& vf,
47  const fvMesh& sMesh,
48  const labelUList& patchMap,
49  const labelUList& cellMap,
50  const labelUList& faceMap,
51  const bool allowUnmapped
52 )
53 {
54  // 1. Create the complete field with dummy patch fields
55  PtrList<fvPatchField<Type>> patchFields(patchMap.size());
56 
57  forAll(patchFields, patchi)
58  {
59  // Set the first one by hand as it corresponds to the
60  // exposed internal faces. Additional interpolation can be put here
61  // as necessary.
62  if (patchMap[patchi] == -1)
63  {
64  patchFields.set
65  (
66  patchi,
67  new emptyFvPatchField<Type>
68  (
69  sMesh.boundary()[patchi],
71  )
72  );
73  }
74  else
75  {
76  patchFields.set
77  (
78  patchi,
80  (
82  sMesh.boundary()[patchi],
84  )
85  );
86  }
87  }
88 
89  auto tresult = tmp<GeometricField<Type, fvPatchField, volMesh>>::New
90  (
91  IOobject
92  (
93  "subset"+vf.name(),
94  sMesh.time().timeName(),
95  sMesh,
98  ),
99  sMesh,
100  vf.dimensions(),
101  Field<Type>(vf.primitiveField(), cellMap),
102  patchFields
103  );
104  auto& result = tresult.ref();
105  result.oriented() = vf.oriented();
106 
107 
108  // 2. Change the fvPatchFields to the correct type using a mapper
109  // constructor (with reference to the now correct internal field)
110 
111  auto& bf = result.boundaryFieldRef();
112 
113  forAll(bf, patchi)
114  {
115  const label basePatchId = patchMap[patchi];
116 
117  if (basePatchId != -1)
118  {
119  // Construct addressing
120  const fvPatch& subPatch = sMesh.boundary()[patchi];
121  const fvPatch& basePatch = vf.mesh().boundary()[basePatchId];
122  const label baseStart = basePatch.start();
123  const label baseSize = basePatch.size();
124 
125  labelList directAddressing(subPatch.size());
126 
127  forAll(directAddressing, i)
128  {
129  const label baseFacei = faceMap[subPatch.start()+i];
130 
131  if (baseFacei >= baseStart && baseFacei < baseStart+baseSize)
132  {
133  directAddressing[i] = baseFacei-baseStart;
134  }
135  else
136  {
137  // Mapped from internal face. Do what? Leave up to
138  // fvPatchField
139  directAddressing[i] = -1;
140  }
141  }
142 
143 
144  directFvPatchFieldMapper mapper(directAddressing);
145 
146  // allowUnmapped : special mode for if we do not want to be
147  // warned for unmapped faces (e.g. from fvMeshDistribute).
148 
149  const bool hasUnmapped = mapper.hasUnmapped();
150  if (allowUnmapped)
151  {
152  mapper.hasUnmapped() = false;
153  }
154 
155  bf.set
156  (
157  patchi,
159  (
160  vf.boundaryField()[basePatchId],
161  subPatch,
162  result(),
163  mapper
164  )
165  );
166 
167  if (allowUnmapped && hasUnmapped)
168  {
169  // Set unmapped values to zeroGradient. This is the default
170  // action for unmapped fvPatchFields. Note that this bypasses
171  // any special logic for handling unmapped fvPatchFields but
172  // since this is only used inside fvMeshDistribute ...
173 
174  tmp<Field<Type>> tfld(bf[patchi].patchInternalField());
175  const Field<Type>& fld = tfld();
176 
177  Field<Type> value(bf[patchi]);
178  forAll(directAddressing, i)
179  {
180  if (directAddressing[i] == -1)
181  {
182  value[i] = fld[i];
183  }
184  }
185  bf[patchi].fvPatchField<Type>::operator=(value);
186  }
187  }
188  }
189 
190  return tresult;
191 }
192 
193 
194 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
195 
196 template<class Type>
197 Foam::tmp
198 <
200 >
202 (
203  const GeometricField<Type, fvPatchField, volMesh>& vf,
204  const bool allowUnmapped
205 ) const
206 {
207  return interpolate
208  (
209  vf,
210  subMesh(),
211  patchMap(),
212  cellMap(),
213  faceMap(),
214  allowUnmapped
215  );
216 }
217 
218 
219 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
220 
221 template<class Type>
222 Foam::tmp
223 <
225 >
227 (
228  const GeometricField<Type, fvsPatchField, surfaceMesh>& vf,
229  const fvMesh& sMesh,
230  const labelUList& patchMap,
231  const labelUList& cellMap,
232  const labelUList& faceMap
233 )
234 {
235  // 1. Create the complete field with dummy patch fields
236  PtrList<fvsPatchField<Type>> patchFields(patchMap.size());
237 
238  forAll(patchFields, patchi)
239  {
240  // Set the first one by hand as it corresponds to the
241  // exposed internal faces. Additional interpolation can be put here
242  // as necessary.
243  if (patchMap[patchi] == -1)
244  {
245  patchFields.set
246  (
247  patchi,
248  new emptyFvsPatchField<Type>
249  (
250  sMesh.boundary()[patchi],
252  )
253  );
254  }
255  else
256  {
257  patchFields.set
258  (
259  patchi,
261  (
263  sMesh.boundary()[patchi],
265  )
266  );
267  }
268  }
269 
270  // Create the complete field from the pieces
271  auto tresult = tmp<GeometricField<Type, fvsPatchField, surfaceMesh>>::New
272  (
273  IOobject
274  (
275  "subset"+vf.name(),
276  sMesh.time().timeName(),
277  sMesh,
280  ),
281  sMesh,
282  vf.dimensions(),
283  Field<Type>
284  (
285  vf.primitiveField(),
286  SubList<label>(faceMap, sMesh.nInternalFaces())
287  ),
288  patchFields
289  );
290  auto& result = tresult.ref();
291  result.oriented() = vf.oriented();
292 
293 
294  // 2. Change the fvsPatchFields to the correct type using a mapper
295  // constructor (with reference to the now correct internal field)
296 
297  auto& bf = result.boundaryFieldRef();
298 
299  forAll(bf, patchi)
300  {
301  if (patchMap[patchi] != -1)
302  {
303  // Construct addressing
304  const fvPatch& subPatch = sMesh.boundary()[patchi];
305  const fvPatch& basePatch = vf.mesh().boundary()[patchMap[patchi]];
306  const label baseStart = basePatch.start();
307  const label baseSize = basePatch.size();
308 
309  labelList directAddressing(subPatch.size());
310 
311  forAll(directAddressing, i)
312  {
313  label baseFacei = faceMap[subPatch.start()+i];
314 
315  if (baseFacei >= baseStart && baseFacei < baseStart+baseSize)
316  {
317  directAddressing[i] = baseFacei-baseStart;
318  }
319  else
320  {
321  // Mapped from internal face. Do what? Leave up to
322  // patchField. This would require also to pass in
323  // original internal field so for now do as post-processing
324  directAddressing[i] = -1;
325  }
326  }
327 
328  bf.set
329  (
330  patchi,
332  (
333  vf.boundaryField()[patchMap[patchi]],
334  subPatch,
335  result(),
336  directFvPatchFieldMapper(directAddressing)
337  )
338  );
339 
340 
341  // Post-process patch field for exposed faces
342 
343  fvsPatchField<Type>& pfld = bf[patchi];
344  const labelUList& fc = bf[patchi].patch().faceCells();
345  const labelList& own = vf.mesh().faceOwner();
346 
347  forAll(pfld, i)
348  {
349  label baseFacei = faceMap[subPatch.start()+i];
350  if (baseFacei < vf.primitiveField().size())
351  {
352  Type val = vf.internalField()[baseFacei];
353 
354  if (cellMap[fc[i]] == own[baseFacei] || !vf.is_oriented())
355  {
356  pfld[i] = val;
357  }
358  else
359  {
360  pfld[i] = flipOp()(val);
361  }
362  }
363  else
364  {
365  // Exposed face from other patch.
366  // Only possible in case of a coupled boundary
367  label patchi = vf.mesh().boundaryMesh().whichPatch
368  (
369  baseFacei
370  );
371  const fvPatch& otherPatch = vf.mesh().boundary()[patchi];
372  label patchFacei = otherPatch.patch().whichFace(baseFacei);
373  pfld[i] = vf.boundaryField()[patchi][patchFacei];
374  }
375  }
376  }
377  }
378 
379  return tresult;
380 }
381 
382 
383 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
384 
385 template<class Type>
386 Foam::tmp
387 <
389 >
391 (
392  const GeometricField<Type, fvsPatchField, surfaceMesh>& sf,
393  const bool allowUnmapped
394 ) const
395 {
396  return interpolate
397  (
398  sf,
399  subMesh(),
400  patchMap(),
401  cellMap(),
402  faceMap()
403  );
404 }
405 
406 
407 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
408 
409 template<class Type>
410 Foam::tmp
411 <
413 >
415 (
416  const GeometricField<Type, pointPatchField, pointMesh>& vf,
417  const pointMesh& sMesh,
418  const labelUList& patchMap,
419  const labelUList& pointMap
420 )
421 {
422  // 1. Create the complete field with dummy patch fields
423  PtrList<pointPatchField<Type>> patchFields(patchMap.size());
424 
425  forAll(patchFields, patchi)
426  {
427  // Set the first one by hand as it corresponds to the
428  // exposed internal faces. Additional interpolation can be put here
429  // as necessary.
430  if (patchMap[patchi] == -1)
431  {
432  patchFields.set
433  (
434  patchi,
435  new emptyPointPatchField<Type>
436  (
437  sMesh.boundary()[patchi],
439  )
440  );
441  }
442  else
443  {
444  patchFields.set
445  (
446  patchi,
448  (
450  sMesh.boundary()[patchi],
452  )
453  );
454  }
455  }
456 
457  // Create the complete field from the pieces
458  auto tresult = tmp<GeometricField<Type, pointPatchField, pointMesh>>::New
459  (
460  IOobject
461  (
462  "subset"+vf.name(),
463  sMesh.time().timeName(),
464  sMesh.thisDb(),
467  ),
468  sMesh,
469  vf.dimensions(),
470  Field<Type>(vf.primitiveField(), pointMap),
471  patchFields
472  );
473  auto& result = tresult.ref();
474  result.oriented() = vf.oriented();
475 
476 
477  // 2. Change the pointPatchFields to the correct type using a mapper
478  // constructor (with reference to the now correct internal field)
479 
480  auto& bf = result.boundaryFieldRef();
481 
482  forAll(bf, patchi)
483  {
484  // Set the first one by hand as it corresponds to the
485  // exposed internal faces. Additional interpolation can be put here
486  // as necessary.
487  if (patchMap[patchi] != -1)
488  {
489  // Construct addressing
490  const pointPatch& basePatch =
491  vf.mesh().boundary()[patchMap[patchi]];
492 
493  const labelList& meshPoints = basePatch.meshPoints();
494 
495  // Make addressing from mesh to patch point
496  Map<label> meshPointMap(2*meshPoints.size());
497  forAll(meshPoints, localI)
498  {
499  meshPointMap.insert(meshPoints[localI], localI);
500  }
501 
502  // Find which subpatch points originate from which patch point
503  const pointPatch& subPatch = sMesh.boundary()[patchi];
504  const labelList& subMeshPoints = subPatch.meshPoints();
505 
506  // If mapped from outside patch leave handling up to patchField
507  labelList directAddressing(subPatch.size(), -1);
508 
509  forAll(subMeshPoints, localI)
510  {
511  // Get mesh point on original mesh.
512  label meshPointi = pointMap[subMeshPoints[localI]];
513 
514  const auto iter = meshPointMap.cfind(meshPointi);
515 
516  if (iter.good())
517  {
518  directAddressing[localI] = *iter;
519  }
520  }
521 
522  bf.set
523  (
524  patchi,
526  (
527  vf.boundaryField()[patchMap[patchi]],
528  subPatch,
529  result(),
530  directPointPatchFieldMapper(directAddressing)
531  )
532  );
533  }
534  }
535 
536  return tresult;
537 }
538 
539 
540 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
541 
542 template<class Type>
543 Foam::tmp
544 <
546 >
548 (
549  const GeometricField<Type, pointPatchField, pointMesh>& sf,
550  const bool allowUnmapped
551 ) const
552 {
553  return interpolate
554  (
555  sf,
556  pointMesh::New(subMesh()), // subsetted point mesh
557  patchMap(),
558  pointMap()
559  );
560 }
561 
562 
563 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
564 
565 template<class Type>
566 Foam::tmp
567 <
569 >
571 (
572  const DimensionedField<Type, volMesh>& df,
573  const fvMesh& sMesh,
574  const labelUList& cellMap
575 )
576 {
577  auto tresult = tmp<DimensionedField<Type, volMesh>>::New
578  (
579  IOobject
580  (
581  "subset"+df.name(),
582  sMesh.time().timeName(),
583  sMesh,
586  ),
587  sMesh,
588  df.dimensions(),
589  Field<Type>(df, cellMap)
590  );
591  tresult.ref().oriented() = df.oriented();
592 
593  return tresult;
594 }
595 
596 
597 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
598 
599 template<class Type>
600 Foam::tmp
601 <
603 >
605 (
606  const DimensionedField<Type, volMesh>& df,
607  const bool allowUnmapped
608 ) const
609 {
610  return interpolate(df, subMesh(), cellMap());
611 }
612 
613 
614 // ************************************************************************* //
static tmp< DimensionedField< Type, volMesh > > interpolate(const DimensionedField< Type, volMesh > &, const fvMesh &sMesh, const labelUList &cellMap)
Map volume internal (dimensioned) field.
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
static const pointMesh & New(const polyMesh &mesh, Args &&... args)
Get existing or create a new MeshObject. Registered with typeName.
Definition: MeshObject.C:53
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Ignore writing from objectRegistry::writeObject()
DirectFieldMapper< fvPatchFieldMapper > directFvPatchFieldMapper
A fvPatchFieldMapper with direct mapping.
const labelList & faceMap() const
Return face map.
Definition: fvMeshSubsetI.H:65
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
const labelList & patchMap() const
Return patch map.
Definition: fvMeshSubsetI.H:92
static const DimensionedField< Type, GeoMesh > & null()
Return a NullObjectRef DimensionedField.
static autoPtr< pointPatchField< Type > > New(const word &patchFieldType, const pointPatch &p, const DimensionedField< Type, pointMesh > &iF)
Return a pointer to a new patchField created on freestore given.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
const labelList & cellMap() const
Return cell map.
Definition: fvMeshSubsetI.H:84
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< ' ';}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< ' ';}gmvFile<< nl;for(const word &name :lagrangianScalarNames){ IOField< scalar > fld(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Nothing to be read.
static const word & calculatedType() noexcept
The type name for calculated patch fields.
Definition: fvPatchField.H:204
List< label > labelList
A List of labels.
Definition: List.H:62
DirectFieldMapper< pointPatchFieldMapper > directPointPatchFieldMapper
A direct pointPatchFieldMapper.
A class for managing temporary objects.
Definition: HashPtrTable.H:50
static tmp< fvPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, volMesh > &)
Return a pointer to a new patchField created on freestore given.
static tmp< fvsPatchField< Type > > New(const word &patchFieldType, const fvPatch &, const DimensionedField< Type, surfaceMesh > &)
Return a pointer to a new patchField created on freestore given.