conformedFvPatchField.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 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 
26 #include "conformedFvPatchField.H"
27 #include "fvMeshStitcherTools.H"
28 #include "nonConformalBoundary.H"
29 #include "nonConformalFvPatch.H"
32 #include "forwardFieldMapper.H"
36 #include "volFields.H"
37 
38 // * * * * * * * * * * Private Static Member Functions * * * * * * * * * * * //
39 
40 template<class Type>
42 (
43  const fvBoundaryMesh& fvbm
44 )
45 {
46  labelList result(fvbm.size(), -1);
47  labelList origNcPatchCount(fvbm.size(), 0);
48 
49  // Global non-conformal patches
50  forAll(fvbm, ncPatchi)
51  {
52  const fvPatch& fvp = fvbm[ncPatchi];
53 
54  if (!isA<nonConformalFvPatch>(fvp)) continue;
55  if (isA<processorCyclicFvPatch>(fvp)) continue;
56 
57  const nonConformalFvPatch& ncFvp =
58  refCast<const nonConformalFvPatch>(fvp);
59 
60  result[ncPatchi] = origNcPatchCount[ncFvp.origPatchIndex()] ++;
61  }
62 
63  // Non-conformal processor cyclics
64  forAll(fvbm, ncPatchi)
65  {
66  const fvPatch& fvp = fvbm[ncPatchi];
67 
68  if (!isA<nonConformalFvPatch>(fvp)) continue;
69  if (!isA<processorCyclicFvPatch>(fvp)) continue;
70 
71  const nonConformalProcessorCyclicFvPatch& ncpcFvp =
72  refCast<const nonConformalProcessorCyclicFvPatch>(fvp);
73 
74  result[ncPatchi] = result[ncpcFvp.referPatchIndex()];
75  }
76 
77  return result;
78 }
79 
80 
81 // * * * * * * * * * * * * * Private Constructors * * * * * * * * * * * * * //
82 
83 template<class Type>
85 (
86  const fvPatch& p,
87  const DimensionedField<Type, volMesh>& iF,
88  autoPtr<fvPatchField<Type>>&& origFieldPtr,
89  PtrList<fvPatchField<Type>>&& ncFieldPtrs,
90  PtrList<scalarField>&& ncCoverages
91 )
92 :
93  fvPatchField<Type>(p, iF),
94  origFieldPtr_(origFieldPtr),
95  ncFieldPtrs_(ncFieldPtrs),
96  ncCoverages_(ncCoverages)
97 {}
98 
99 
100 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
101 
102 template<class Type>
104 {
105  const fvPatch& origFvp = this->patch();
106  const fvBoundaryMesh& fvbm = origFvp.boundaryMesh();
107 
108  // Determine the indices of the non-conformal patches associated with this
109  // original patch. This assumes that all meshes that we might be mapping
110  // between have their (non-processor) patches in the same order.
111  DynamicList<label> result(ncFieldPtrs_.size());
112  forAll(fvbm, ncPatchi)
113  {
114  const fvPatch& fvp = fvbm[ncPatchi];
115 
116  if (!isA<nonConformalFvPatch>(fvp)) continue;
117  if (isA<processorCyclicFvPatch>(fvp)) continue;
118 
119  const nonConformalFvPatch& ncFvp =
120  refCast<const nonConformalFvPatch>(fvp);
121 
122  if (ncFvp.origPatchName() != origFvp.name()) continue;
123 
124  result.append(ncPatchi);
125  }
126 
127  // Check the number of identified indices corresponds to the number of
128  // stored patch fields
129  if (result.size() != ncFieldPtrs_.size())
130  {
131  const DimensionedField<Type, volMesh>& iF = this->internalField();
132 
134  << "Conformed field " << (isNull(iF) ? "" : iF.name() + ' ')
135  << "on patch " << origFvp.name() << " has " << ncFieldPtrs_.size()
136  << " non-conformal fields stored, but the patch has "
137  << result.size() << " associated non-conformal patches"
138  << exit(FatalError);
139  }
140 
141  return result;
142 }
143 
144 
145 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
146 
147 template<class Type>
149 (
150  typename VolField<Type>::Boundary& bF
151 )
152 {
153  const DimensionedField<Type, volMesh>& iF = bF[0].internalField();
154 
155  const fvBoundaryMesh& fvbm = iF.mesh().boundary();
156 
157  const labelList origPatchIndices =
159 
160  const labelList ncOrigNcField =
162 
163  // Replace every patch field on an orig boundary with a conformed patch
164  // field, and move the previous orig boundary field into the member data of
165  // the conformed patch field
166  forAll(origPatchIndices, i)
167  {
168  const label origPatchi = origPatchIndices[i];
169  const fvPatch& origFvp = fvbm[origPatchi];
170 
172  (
174  (
175  origFvp,
176  iF,
177  bF.set(origPatchi, nullptr),
180  )
181  );
182 
183  pF() == pF->origFieldPtr_();
184 
185  bF.set(origPatchi, pF.ptr());
186  }
187 
188  // Copy every (non-processor) non-conformal patch field into the conformed
189  // patch field on the associated orig patch.
190  forAll(fvbm, ncPatchi)
191  {
192  const fvPatch& fvp = fvbm[ncPatchi];
193 
194  if (!isA<nonConformalFvPatch>(fvp)) continue;
195  if (isA<processorCyclicFvPatch>(fvp)) continue;
196 
197  const nonConformalFvPatch& ncFvp =
198  refCast<const nonConformalFvPatch>(fvp);
199 
200  const label origPatchi = ncFvp.origPatchIndex();
201  const fvPatch& origFvp = fvbm[origPatchi];
202 
204  refCast<conformedFvPatchField<Type>>(bF[origPatchi]);
205 
206  const label n =
207  max(cpF.ncFieldPtrs_.size(), ncOrigNcField[ncPatchi] + 1);
208 
209  cpF.ncFieldPtrs_.resize(n);
210  cpF.ncFieldPtrs_.set
211  (
212  ncOrigNcField[ncPatchi],
214  (
215  bF[ncPatchi],
216  fvp,
217  iF,
219  )
220  );
221 
222  cpF.ncCoverages_.resize(n);
223  cpF.ncCoverages_.set
224  (
225  ncOrigNcField[ncPatchi],
226  scalarField(origFvp.size(), 0)
227  );
228  }
229 
230  // Compute the coverage for every non-conformal patch
231  forAll(fvbm, ncPatchi)
232  {
233  const fvPatch& fvp = fvbm[ncPatchi];
234 
235  if (!isA<nonConformalFvPatch>(fvp)) continue;
236 
237  const nonConformalFvPatch& ncFvp =
238  refCast<const nonConformalFvPatch>(fvp);
239 
240  const label origPatchi = ncFvp.origPatchIndex();
241  const fvPatch& origFvp = fvbm[origPatchi];
242 
243  const labelList ncOrigPatchFace =
244  ncFvp.polyFaces() - origFvp.start();
245 
246  const scalarField origNcMagSf
247  (
249  (
250  ncFvp.patch().magSf(),
251  origFvp.size(),
252  ncOrigPatchFace
253  )
254  );
255 
257  refCast<conformedFvPatchField<Type>>(bF[origPatchi]);
258 
259  cpF.ncCoverages_[ncOrigNcField[ncPatchi]] +=
260  origNcMagSf/origFvp.magSf();
261  }
262 
263  // Map the data from the non-conformal patch fields to the conformed fields
264  // by doing an area-weighted average
265  forAll(fvbm, ncPatchi)
266  {
267  const fvPatch& fvp = fvbm[ncPatchi];
268 
269  if (!isA<nonConformalFvPatch>(fvp)) continue;
270 
271  const nonConformalFvPatch& ncFvp =
272  refCast<const nonConformalFvPatch>(fvp);
273 
274  const label origPatchi = ncFvp.origPatchIndex();
275  const fvPatch& origFvp = fvbm[origPatchi];
276 
277  const labelList ncOrigPatchFace =
278  ncFvp.polyFaces() - origFvp.start();
279 
281  refCast<conformedFvPatchField<Type>>(bF[origPatchi]);
282 
283  cpF.ncFieldPtrs_[ncOrigNcField[ncPatchi]].map
284  (
285  bF[ncPatchi],
287  (
288  ncOrigPatchFace,
289  ncFvp.patch().magSf()
290  /max
291  (
293  (
294  cpF.ncCoverages_[ncOrigNcField[ncPatchi]]
295  *origFvp.magSf(),
296  ncOrigPatchFace
297  ),
298  rootVSmall
299  )
300  )
301  );
302  }
303 }
304 
305 
306 template<class Type>
308 (
309  typename VolField<Type>::Boundary& bF
310 )
311 {
312  const DimensionedField<Type, volMesh>& iF = bF[0].internalField();
313 
314  const fvBoundaryMesh& fvbm = iF.mesh().boundary();
315 
316  const labelList origPatchIndices =
318 
319  const labelList ncOrigNcField =
321 
322  // Map out of the conformed fields into the non-conformal patches
323  forAllReverse(fvbm, ncPatchi)
324  {
325  const fvPatch& fvp = fvbm[ncPatchi];
326 
327  if (!isA<nonConformalFvPatch>(fvp)) continue;
328 
329  const nonConformalFvPatch& ncFvp =
330  refCast<const nonConformalFvPatch>(fvp);
331 
332  const label origPatchi = ncFvp.origPatchIndex();
333  const fvPatch& origFvp = fvbm[origPatchi];
334 
336  refCast<conformedFvPatchField<Type>>(bF[origPatchi]);
337 
338  labelList ncOrigPatchFace =
339  ncFvp.polyFaces() - origFvp.start();
340 
341  forAll(ncOrigPatchFace, ncPatchFacei)
342  {
343  const label origPatchFacei = ncOrigPatchFace[ncPatchFacei];
344 
345  if (cpF.ncCoverages_[ncOrigNcField[ncPatchi]][origPatchFacei] == 0)
346  {
347  ncOrigPatchFace[ncPatchFacei] = -1;
348  }
349  }
350 
351  // If this is a calculated condition, then set unmapped face values to
352  // the value in the adjacent cell. This is a tolerable approximation
353  // and prevents errors in calculations that occur before the values are
354  // re-computed as part of the solution.
355  if (isType<calculatedFvPatchField<Type>>(bF[ncPatchi]))
356  {
357  bF[ncPatchi] = bF[ncPatchi].patchInternalField();
358 
359  bF[ncPatchi].map
360  (
361  cpF.ncFieldPtrs_[ncOrigNcField[ncPatchi]],
362  forwardFieldMapper(ncOrigPatchFace)
363  );
364  }
365  // Otherwise let the boundary condition determine what an appropriate
366  // value is for an unmapped face
367  else
368  {
369  bF[ncPatchi].map
370  (
371  cpF.ncFieldPtrs_[ncOrigNcField[ncPatchi]],
372  forwardOrAssignPatchFieldMapper(bF[ncPatchi], ncOrigPatchFace)
373  );
374  }
375  }
376 
377  // Move the orig patch fields into the boundary field
378  forAll(origPatchIndices, i)
379  {
380  const label origPatchi = origPatchIndices[i];
381 
383  refCast<conformedFvPatchField<Type>>(bF[origPatchi]);
384 
385  bF.set(origPatchi, cpF.origFieldPtr_.ptr());
386  }
387 }
388 
389 
390 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
391 
392 template<class Type>
394 (
395  const fvPatch& p,
397 )
398 :
399  fvPatchField<Type>(p, iF)
400 {
402 }
403 
404 
405 template<class Type>
407 (
408  const fvPatch& p,
410  const dictionary& dict
411 )
412 :
413  fvPatchField<Type>(p, iF, dict),
414  origFieldPtr_
415  (
416  fvPatchField<Type>::New(p, iF, dict.subDict("origField")).ptr()
417  ),
418  ncFieldPtrs_(dict.subDict("ncFields").size()),
419  ncCoverages_(dict.subDict("ncFields").size())
420 {
421  label i = 0;
422 
423  forAllConstIter(dictionary, dict.subDict("ncFields"), iter)
424  {
425  const word& ncPName = iter().keyword();
426  const fvPatch& ncP = p.boundaryMesh()[ncPName];
427 
428  ncFieldPtrs_.set
429  (
430  i,
431  fvPatchField<Type>::New(ncP, iF, iter().dict())
432  );
433  ncCoverages_.set
434  (
435  i,
436  new scalarField("coverage", iter().dict(), p.size())
437  );
438 
439  ++ i;
440  }
441 }
442 
443 
444 template<class Type>
446 (
447  const conformedFvPatchField<Type>& ptf,
448  const fvPatch& p,
450  const fieldMapper& mapper
451 )
452 :
453  fvPatchField<Type>(ptf, p, iF, mapper),
454  origFieldPtr_
455  (
456  fvPatchField<Type>::New(ptf.origFieldPtr_(), p, iF, mapper).ptr()
457  ),
458  ncFieldPtrs_(ptf.ncFieldPtrs_.size()),
459  ncCoverages_(ptf.ncCoverages_.size())
460 {
461  const labelList ncPatchIndices(this->ncPatchIndices());
462 
463  forAll(ptf.ncFieldPtrs_, i)
464  {
465  ncFieldPtrs_.set
466  (
467  i,
469  (
470  ptf.ncFieldPtrs_[i],
471  p.boundaryMesh()[ncPatchIndices[i]],
472  iF,
473  mapper
474  )
475  );
476  ncCoverages_.set
477  (
478  i,
479  mapper(ptf.ncCoverages_[i]).ptr()
480  );
481  }
482 }
483 
484 
485 template<class Type>
487 (
488  const conformedFvPatchField<Type>& ptf,
490 )
491 :
492  fvPatchField<Type>(ptf, iF),
493  origFieldPtr_(ptf.origFieldPtr_->clone(iF).ptr()),
494  ncFieldPtrs_(ptf.ncFieldPtrs_, iF),
495  ncCoverages_(ptf.ncCoverages_)
496 {}
497 
498 
499 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
500 
501 template<class Type>
503 (
504  const fvPatchField<Type>& ptf,
505  const fieldMapper& mapper
506 )
507 {
508  fvPatchField<Type>::map(ptf, mapper);
509 
511  {
512  const conformedFvPatchField<Type>& cptf =
513  refCast<const conformedFvPatchField<Type>>(ptf);
514 
515  origFieldPtr_->map(cptf.origFieldPtr_(), mapper);
516  forAll(ncFieldPtrs_, i)
517  {
518  ncFieldPtrs_[i].map(cptf.ncFieldPtrs_[i], mapper);
519  mapper(ncCoverages_[i], cptf.ncCoverages_[i]);
520  }
521  }
522  else
523  {
525  }
526 }
527 
528 
529 template<class Type>
531 {
533 
535  {
536  const conformedFvPatchField<Type>& cptf =
537  refCast<const conformedFvPatchField<Type>>(ptf);
538 
539  origFieldPtr_->reset(cptf.origFieldPtr_());
540  forAll(ncFieldPtrs_, i)
541  {
542  ncFieldPtrs_[i].reset(cptf.ncFieldPtrs_[i]);
543  ncCoverages_[i] = cptf.ncCoverages_[i];
544  }
545  }
546  else
547  {
549  }
550 }
551 
552 
553 template<class Type>
555 {
557  writeEntry(os, "value", *this);
558 
559  writeKeyword(os, "origField") << nl;
560  os << indent << token::BEGIN_BLOCK << incrIndent << nl;
561  origFieldPtr_->write(os);
562  os << decrIndent << indent << token::END_BLOCK << nl;
563 
564  writeKeyword(os, "ncFields") << nl;
565  os << indent << token::BEGIN_BLOCK << incrIndent << nl;
566  forAll(ncFieldPtrs_, i)
567  {
568  writeKeyword(os, ncFieldPtrs_[i].patch().name()) << nl;
569  os << indent << token::BEGIN_BLOCK << incrIndent << nl;
570  ncFieldPtrs_[i].write(os);
571  writeEntry(os, "coverage", ncCoverages_[i]);
572  os << decrIndent << indent << token::END_BLOCK << nl;
573  }
574  os << decrIndent << indent << token::END_BLOCK << nl;
575 }
576 
577 
578 // ************************************************************************* //
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllReverse(list, i)
Reverse loop across all elements in list.
Definition: UList.H:446
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
static nonConformalBoundary & New(const word &name, const polyMesh &mesh)
Construct and return the named DemandDrivenMeshObject.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const Mesh & mesh() const
Return mesh.
Generic GeometricBoundaryField class.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
T * ptr()
Return object pointer for reuse.
Definition: autoPtrI.H:90
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
This vol field boundary condition holds data from both the original faces and any associated non-conf...
virtual void write(Ostream &) const
Write.
virtual void reset(const fvPatchField< Type > &)
Reset the fvPatchField to the given fvPatchField.
static void conform(typename VolField< Type >::Boundary &bF)
Conform the given boundary field.
static void unconform(typename VolField< Type >::Boundary &bF)
Un-conform the given boundary field.
virtual void map(const fvPatchField< Type > &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Abstract base class for field mapping.
Definition: fieldMapper.H:48
Forward field mapper.
Foam::fvBoundaryMesh.
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Definition: fvPatchField.H:88
virtual void write(Ostream &) const
Write.
Definition: fvPatchField.C:229
virtual void reset(const fvPatchField< Type > &)
Reset the fvPatchField to the given fvPatchField.
Definition: fvPatchField.C:195
virtual void map(const fvPatchField< Type > &, const fieldMapper &)
Map the given fvPatchField onto this fvPatchField.
Definition: fvPatchField.C:185
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual label size() const
Return size.
Definition: fvPatch.H:138
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:132
const scalarField & magSf() const
Return face area magnitudes.
Definition: fvPatch.C:142
labelList allOrigPatchIndices() const
Return a list of the orig patch indices.
Non-conformal FV patch. Provides the necessary interface for a FV patch which does not conform to the...
const fvPatch & patch() const
Reference to the fvPatch.
const labelList & polyFaces() const
Return face face-poly-faces.
label origPatchIndex() const
Original patch ID.
Mapper which sets the field size and initialises all values to zero. It does not actually map values.
@ BEGIN_BLOCK
Definition: token.H:113
@ END_BLOCK
Definition: token.H:114
A class for handling words, derived from string.
Definition: word.H:62
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
tmp< Field< Type > > fieldMap(const Field< Type > &f, const labelUList &addr, const label addr0=0)
Map a field with an (optional) addressing offset.
tmp< Field< Type > > fieldRMapSum(const Field< Type > &f, const label size, const labelUList &addr, const label addr0=0)
Reverse map a field with an (optional) addressing offset, initialising the.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:241
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
bool isType(const Type &t)
Check the typeid.
Definition: typeInfo.H:158
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:234
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
bool isA(const Type &t)
Check if a dynamic_cast to typeid is possible.
Definition: typeInfo.H:166
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
T clone(const T &t)
Definition: List.H:55
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
bool isNull(const T &t)
Return true if t is a reference to the nullObject of type T.
Definition: nullObjectI.H:58
error FatalError
Ostream & writeKeyword(Foam::Ostream &os, const keyType &kw)
Write the keyword to the Ostream with the current level of indentation.
Definition: keyType.C:155
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:227
static const char nl
Definition: Ostream.H:266
dictionary dict
volScalarField & p