fvMeshStitcherTemplates.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) 2022-2023 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 Description
25  Perform mapping of finite volume fields required by stitching.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "volFields.H"
30 #include "surfaceFields.H"
31 #include "fvPatchFieldMapper.H"
32 #include "fvMeshStitcher.H"
34 #include "nonConformalBoundary.H"
38 #include "processorFvPatch.H"
39 #include "surfaceToVolVelocity.H"
40 
41 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
42 
43 template<class Type, template<class> class GeoField>
44 void Foam::fvMeshStitcher::resizePatchFields()
45 {
46  HashTable<GeoField<Type>*> fields(mesh_.lookupClass<GeoField<Type>>());
47  forAllIter(typename HashTable<GeoField<Type>*>, fields, iter)
48  {
49  forAll(mesh_.boundary(), patchi)
50  {
51  typename GeoField<Type>::Patch& pf =
52  iter()->boundaryFieldRef()[patchi];
53 
54  if (isA<nonConformalFvPatch>(pf.patch()))
55  {
56  pf.map(pf, setSizeFvPatchFieldMapper(pf.patch().size()));
57  }
58  }
59  }
60 }
61 
62 
63 template<template<class> class GeoField>
64 void Foam::fvMeshStitcher::resizePatchFields()
65 {
66  #define ResizePatchFields(Type, nullArg) \
67  resizePatchFields<Type, GeoField>();
69  #undef ResizePatchFields
70 }
71 
72 
73 template<class Type>
74 void Foam::fvMeshStitcher::preConformSurfaceFields()
75 {
76  HashTable<SurfaceField<Type>*> fields
77  (
78  mesh_.lookupClass<SurfaceField<Type>>()
79  );
80 
81  const surfaceScalarField* phiPtr =
82  mesh_.foundObject<surfaceScalarField>("meshPhi")
83  ? &mesh_.lookupObject<surfaceScalarField>("meshPhi")
84  : nullptr;
85 
86  forAllIter(typename HashTable<SurfaceField<Type>*>, fields, iter)
87  {
88  SurfaceField<Type>& field = *iter();
89 
90  if
91  (
92  &field == static_cast<const regIOobject*>(phiPtr)
93  || field.isOldTime()
94  ) continue;
95 
96  autoPtr<SurfaceField<Type>> nccFieldPtr
97  (
98  new SurfaceField<Type>
99  (
100  IOobject
101  (
102  nccFieldPrefix_ + field.name(),
103  mesh_.time().name(),
104  mesh_
105  ),
106  field
107  )
108  );
109 
110  for (label i = 0; i <= field.nOldTimes(); ++ i)
111  {
112  SurfaceField<Type>& field0 = field.oldTime(i);
113 
114  nccFieldPtr->oldTime(i).boundaryFieldRef() =
115  conformalNccBoundaryField<Type>(field0.boundaryField());
116 
117  field0.boundaryFieldRef() =
118  conformalOrigBoundaryField<Type>(field0.boundaryField());
119  }
120 
121  nccFieldPtr.ptr()->store();
122  }
123 }
124 
125 
126 inline void Foam::fvMeshStitcher::preConformSurfaceFields()
127 {
128  #define PreConformSurfaceFields(Type, nullArg) \
129  preConformSurfaceFields<Type>();
131  #undef PreConformSurfaceFields
132 }
133 
134 
135 template<class Type>
136 void Foam::fvMeshStitcher::postNonConformSurfaceFields()
137 {
138  HashTable<SurfaceField<Type>*> fields
139  (
140  mesh_.lookupClass<SurfaceField<Type>>()
141  );
142 
143  const surfaceScalarField* phiPtr =
144  mesh_.foundObject<surfaceScalarField>("meshPhi")
145  ? &mesh_.lookupObject<surfaceScalarField>("meshPhi")
146  : nullptr;
147 
148  forAllIter(typename HashTable<SurfaceField<Type>*>, fields, iter)
149  {
150  if (iter.key()(nccFieldPrefix_.size()) == nccFieldPrefix_) continue;
151 
152  SurfaceField<Type>& field = *iter();
153 
154  if
155  (
156  &field == static_cast<const regIOobject*>(phiPtr)
157  || field.isOldTime()
158  || mesh_.topoChanged()
159  ) continue;
160 
161  const word nccFieldName = nccFieldPrefix_ + field.name();
162 
163  const SurfaceField<Type>& nccField =
164  mesh_.lookupObject<SurfaceField<Type>>(nccFieldName);
165 
166  for (label i = 0; i <= field.nOldTimes(); ++ i)
167  {
168  SurfaceField<Type>& field0 = field.oldTime(i);
169 
170  field0.boundaryFieldRef() =
171  nonConformalBoundaryField<Type>
172  (
173  nccField.oldTime(i).boundaryField(),
174  field0.boundaryField()
175  );
176 
177  field0.boundaryFieldRef() =
178  synchronisedBoundaryField<Type>(field0.boundaryField());
179  }
180  }
181 
182  // Remove NCC fields after all fields have been mapped. This is so that
183  // old-time fields aren't removed by current-time fields in advance of the
184  // old-time field being mapped.
185  forAllIter(typename HashTable<SurfaceField<Type>*>, fields, iter)
186  {
187  if (iter.key()(nccFieldPrefix_.size()) == nccFieldPrefix_) continue;
188 
189  SurfaceField<Type>& field = *iter();
190 
191  if
192  (
193  &field == static_cast<const regIOobject*>(phiPtr)
194  || field.isOldTime()
195  ) continue;
196 
197  const word nccFieldName = nccFieldPrefix_ + field.name();
198 
199  SurfaceField<Type>& nccField =
200  mesh_.lookupObjectRef<SurfaceField<Type>>(nccFieldName);
201 
202  nccField.checkOut();
203  }
204 
205  // Check there are no NCC fields left over
206  fields = mesh_.lookupClass<SurfaceField<Type>>();
207  forAllIter(typename HashTable<SurfaceField<Type>*>, fields, iter)
208  {
209  if (iter.key()(nccFieldPrefix_.size()) != nccFieldPrefix_) continue;
210 
212  << "Stitching mapping field \"" << iter()->name()
213  << "\" found, but the field it corresponds to no longer exists"
214  << exit(FatalError);
215  }
216 }
217 
218 
219 inline void Foam::fvMeshStitcher::postNonConformSurfaceFields()
220 {
221  #define PostNonConformSurfaceFields(Type, nullArg) \
222  postNonConformSurfaceFields<Type>();
224  #undef PostNonConformSurfaceFields
225 }
226 
227 
228 template<class Type>
229 void Foam::fvMeshStitcher::evaluateVolFields()
230 {
231  HashTable<VolField<Type>*> fields(mesh_.lookupClass<VolField<Type>>());
232  forAllIter(typename HashTable<VolField<Type>*>, fields, iter)
233  {
234  const label nReq = Pstream::nRequests();
235 
236  forAll(mesh_.boundary(), patchi)
237  {
238  typename VolField<Type>::Patch& pf =
239  iter()->boundaryFieldRef()[patchi];
240 
241  if (isA<nonConformalFvPatch>(pf.patch()))
242  {
243  pf.initEvaluate(Pstream::defaultCommsType);
244  }
245  }
246 
247  if
248  (
251  )
252  {
253  Pstream::waitRequests(nReq);
254  }
255 
256  forAll(mesh_.boundary(), patchi)
257  {
258  typename VolField<Type>::Patch& pf =
259  iter()->boundaryFieldRef()[patchi];
260 
261  if (isA<nonConformalFvPatch>(pf.patch()))
262  {
263  pf.evaluate(Pstream::defaultCommsType);
264  }
265  }
266  }
267 }
268 
269 
270 inline void Foam::fvMeshStitcher::evaluateVolFields()
271 {
272  #define EvaluateVolFields(Type, nullArg) \
273  evaluateVolFields<Type>();
275  #undef EvaluateVolFields
276 }
277 
278 
279 inline void Foam::fvMeshStitcher::postNonConformSurfaceVelocities()
280 {
281  if (mesh_.topoChanged())
282  {
283  HashTable<surfaceVectorField*> Ufs
284  (
285  mesh_.lookupClass<surfaceVectorField>()
286  );
287 
288  forAllIter(HashTable<surfaceVectorField*>, Ufs, iter)
289  {
290  surfaceVectorField& Uf = *iter();
291 
293 
294  if (!isNull(U))
295  {
296  forAll(Uf.boundaryField(), patchi)
297  {
298  if (isA<nonConformalFvPatch>(mesh_.boundary()[patchi]))
299  {
300  Uf.boundaryFieldRef()[patchi] ==
301  U.boundaryField()[patchi];
302  }
303  }
304  }
305  }
306  }
307 }
308 
309 
310 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
311 
312 template<class GeoBoundaryField>
314 (
315  const SurfaceFieldBoundary<label>& polyFacesBf,
316  GeoBoundaryField& fieldBf
317 )
318 {
319  forAll(polyFacesBf, nccPatchi)
320  {
321  if (isA<nonConformalFvPatch>(polyFacesBf[nccPatchi].patch()))
322  {
323  fieldBf[nccPatchi].map
324  (
325  fieldBf[nccPatchi],
326  setSizeFvPatchFieldMapper(polyFacesBf[nccPatchi].size())
327  );
328  }
329  }
330 }
331 
332 template<class GeoField>
334 (
335  const SurfaceFieldBoundary<label>& polyFacesBf,
336  GeoField& field
337 )
338 {
339  for (label i = 0; i <= field.nOldTimes(); ++ i)
340  {
341  resizeBoundaryFieldPatchFields
342  (
343  polyFacesBf,
344  field.oldTime(i).boundaryFieldRef()
345  );
346  }
347 }
348 
349 
350 template<class Type>
352 (
353  const Field<Type>& f,
354  const label size,
355  const labelUList& addr
356 )
357 {
358  tmp<Field<Type>> tresult(new Field<Type>(size, Zero));
359  forAll(addr, i)
360  {
361  tresult.ref()[addr[i]] += f[i];
362  }
363  return tresult;
364 }
365 
366 
367 template<class Type>
369 (
370  const tmp<Field<Type>>& f,
371  const label size,
372  const labelUList& addr
373 )
374 {
375  tmp<Field<Type>> tresult = fieldRMapSum(f(), size, addr);
376  f.clear();
377  return tresult;
378 }
379 
380 
383 {
384  const fvBoundaryMesh& fvbm = mesh_.boundary();
385 
386  const surfaceScalarField::Boundary& magSfb =
387  fvbm.mesh().magSf().boundaryField();
388 
389  tmp<surfaceScalarField::Boundary> tOrigNccMagSfb
390  (
392  (
393  fvbm,
396  )
397  );
398 
399  surfaceScalarField::Boundary& origNccMagSfb = tOrigNccMagSfb.ref();
400 
401  origNccMagSfb == 0;
402 
403  forAll(fvbm, nccPatchi)
404  {
405  const fvPatch& fvp = fvbm[nccPatchi];
406 
407  if (isA<nonConformalCoupledFvPatch>(fvp))
408  {
409  const nonConformalCoupledFvPatch& nccFvp =
410  refCast<const nonConformalCoupledFvPatch>(fvp);
411 
412  const label origPatchi = nccFvp.origPatchID();
413  const fvPatch& origFvp = nccFvp.origPatch();
414 
415  const labelList nccOrigPatchFace =
416  nccFvp.polyFaces() - origFvp.start();
417 
418  origNccMagSfb[origPatchi] +=
419  fieldRMapSum
420  (
421  magSfb[nccPatchi],
422  origFvp.size(),
423  nccOrigPatchFace
424  );
425  }
426  }
427 
428  return tOrigNccMagSfb;
429 }
430 
431 
432 template<class Type>
435 (
436  const SurfaceFieldBoundary<Type>& fieldb
437 ) const
438 {
439  const bool isFluxField = isFlux(fieldb[0].internalField());
440 
441  const fvBoundaryMesh& fvbm = fieldb[0].patch().boundaryMesh();
442 
443  const surfaceScalarField::Boundary& magSfb =
444  fvbm.mesh().magSf().boundaryField();
445 
447  (
449  (
451  fieldb
452  )
453  );
454 
455  SurfaceFieldBoundary<Type>& nccFieldb = tnccFieldb.ref();
456 
457  nccFieldb == pTraits<Type>::zero;
458 
459  const surfaceScalarField::Boundary origNccMagSfb
460  (
462  getOrigNccMagSfb()
463  );
464 
465  // Accumulate the non-conformal parts of the field into the original faces
466  forAll(fvbm, nccPatchi)
467  {
468  const fvPatch& fvp = fvbm[nccPatchi];
469 
470  if (isA<nonConformalCoupledFvPatch>(fvp))
471  {
472  const nonConformalCoupledFvPatch& nccFvp =
473  refCast<const nonConformalCoupledFvPatch>(fvbm[nccPatchi]);
474 
475  const label origPatchi = nccFvp.origPatchID();
476  const fvPatch& origFvp = nccFvp.origPatch();
477 
478  const labelList nccOrigPatchFace =
479  nccFvp.polyFaces() - origFvp.start();
480 
481  // If this is a flux then sum
482  if (isFluxField)
483  {
484  nccFieldb[origPatchi] +=
485  fieldRMapSum
486  (
487  fieldb[nccPatchi],
488  origFvp.size(),
489  nccOrigPatchFace
490  );
491  }
492 
493  // If not a flux then do an area-weighted sum
494  else
495  {
496  nccFieldb[origPatchi] +=
497  fieldRMapSum
498  (
499  fieldb[nccPatchi]*magSfb[nccPatchi],
500  origFvp.size(),
501  nccOrigPatchFace
502  );
503  }
504  }
505  }
506 
507  const labelList origPatchIDs =
509 
510  // Scale or average as appropriate
511  forAll(origPatchIDs, i)
512  {
513  const label origPatchi = origPatchIDs[i];
514 
515  // If this is a flux then scale to the total size of the face
516  if (isFluxField)
517  {
518  const scalarField origSumMagSf
519  (
520  magSfb[origPatchi] + origNccMagSfb[origPatchi]
521  );
522 
523  nccFieldb[origPatchi] *=
524  origSumMagSf
525  /max(origNccMagSfb[origPatchi], small*origSumMagSf);
526  }
527 
528  // If this is not a flux then convert to an area-weighted average
529  else
530  {
531  nccFieldb[origPatchi] /=
532  max(origNccMagSfb[origPatchi], vSmall);
533  }
534  }
535 
536  return tnccFieldb;
537 }
538 
539 
540 template<class Type>
543 (
544  const SurfaceFieldBoundary<Type>& fieldb
545 ) const
546 {
547  const bool isFluxField = isFlux(fieldb[0].internalField());
548 
549  const fvBoundaryMesh& fvbm = fieldb[0].patch().boundaryMesh();
550 
551  const surfaceScalarField::Boundary& magSfb =
552  fvbm.mesh().magSf().boundaryField();
553 
554  tmp<SurfaceFieldBoundary<Type>> torigFieldb
555  (
557  (
559  fieldb
560  )
561  );
562 
563  // If this is a flux then scale up to the total face areas
564  if (isFluxField)
565  {
566  const surfaceScalarField::Boundary origNccMagSfb
567  (
569  getOrigNccMagSfb()
570  );
571 
572  SurfaceFieldBoundary<Type>& origFieldb = torigFieldb.ref();
573 
574  const labelList origPatchIDs =
576 
577  forAll(origPatchIDs, i)
578  {
579  const label origPatchi = origPatchIDs[i];
580 
581  const scalarField origSumMagSf
582  (
583  magSfb[origPatchi] + origNccMagSfb[origPatchi]
584  );
585 
586  origFieldb[origPatchi] *=
587  origSumMagSf
588  /max(magSfb[origPatchi], small*origSumMagSf);
589  }
590  }
591 
592  return torigFieldb;
593 }
594 
595 
596 template<class Type>
599 (
600  const SurfaceFieldBoundary<Type>& nccFieldb,
601  const SurfaceFieldBoundary<Type>& origFieldb
602 ) const
603 {
604  const bool isFluxField = isFlux(origFieldb[0].internalField());
605 
606  const fvBoundaryMesh& fvbm = origFieldb[0].patch().boundaryMesh();
607 
608  const surfaceScalarField::Boundary& magSfb =
609  fvbm.mesh().magSf().boundaryField();
610 
612  (
614  (
616  origFieldb
617  )
618  );
619 
620  SurfaceFieldBoundary<Type>& fieldb = tfieldb.ref();
621 
622  // Set the coupled values
623  forAll(fvbm, nccPatchi)
624  {
625  const fvPatch& fvp = fvbm[nccPatchi];
626 
627  if (isA<nonConformalCoupledFvPatch>(fvp))
628  {
629  const nonConformalCoupledFvPatch& nccFvp =
630  refCast<const nonConformalCoupledFvPatch>(fvp);
631 
632  const label origPatchi = nccFvp.origPatchID();
633  const fvPatch& origFvp = nccFvp.origPatch();
634 
635  const labelList nccOrigPatchFace =
636  nccFvp.polyFaces() - origFvp.start();
637 
638  // Set the cyclic values
639  fieldb[nccPatchi] =
640  Field<Type>(nccFieldb[origPatchi], nccOrigPatchFace);
641  }
642  }
643 
644  // If a flux then scale down to the part face area
645  if (isFluxField)
646  {
647  const surfaceScalarField::Boundary origNccMagSfb
648  (
650  getOrigNccMagSfb()
651  );
652 
653  forAll(fvbm, nccPatchi)
654  {
655  const fvPatch& fvp = fvbm[nccPatchi];
656 
657  if (isA<nonConformalCoupledFvPatch>(fvp))
658  {
659  const nonConformalCoupledFvPatch& nccFvp =
660  refCast<const nonConformalCoupledFvPatch>(fvp);
661 
662  const label origPatchi = nccFvp.origPatchID();
663  const fvPatch& origFvp = nccFvp.origPatch();
664 
665  const labelList nccOrigPatchFace =
666  nccFvp.polyFaces() - origFvp.start();
667 
668  const scalarField origSumMagSf
669  (
670  magSfb[origPatchi] + origNccMagSfb[origPatchi]
671  );
672  const scalarField nccSumMagSf(origSumMagSf, nccOrigPatchFace);
673 
674  fieldb[nccPatchi] *= magSfb[nccPatchi]/nccSumMagSf;
675 
676  if (!isA<processorFvPatch>(fvp))
677  {
678  fieldb[origPatchi] *= magSfb[origPatchi]/origSumMagSf;
679  }
680  }
681  }
682  }
683 
684  // Set error values
685  forAll(fvbm, patchi)
686  {
687  const fvPatch& fvp = fvbm[patchi];
688 
689  if (isA<nonConformalErrorFvPatch>(fvp))
690  {
691  const label errorPatchi = patchi;
692 
693  const nonConformalErrorFvPatch& errorFvp =
694  refCast<const nonConformalErrorFvPatch>(fvp);
695 
696  const label origPatchi = errorFvp.origPatchID();
697  const fvPatch& origFvp = errorFvp.origPatch();
698 
699  const labelList errorOrigPatchFace =
700  errorFvp.polyFaces() - origFvp.start();
701 
702  if (isFluxField)
703  {
704  fieldb[errorPatchi] = Zero;
705  }
706  else
707  {
708  fieldb[errorPatchi] =
709  Field<Type>(origFieldb[origPatchi], errorOrigPatchFace);
710  }
711  }
712  }
713 
714  return tfieldb;
715 }
716 
717 
718 template<class Type>
721 (
722  const SurfaceFieldBoundary<Type>& fieldb,
723  const bool flip,
724  const scalar ownerWeight,
725  const scalar neighbourWeight
726 ) const
727 {
728  const fvBoundaryMesh& fvbm = fieldb[0].patch().boundaryMesh();
729 
730  tmp<SurfaceFieldBoundary<Type>> tsyncFieldb
731  (
733  (
735  fieldb
736  )
737  );
738 
739  SurfaceFieldBoundary<Type>& syncFieldb = tsyncFieldb.ref();
740 
742  (
744  fieldb.boundaryNeighbourField()
745  );
746 
747  forAll(fvbm, patchi)
748  {
749  const fvPatch& fvp = fvbm[patchi];
750 
751  if (fvp.coupled())
752  {
753  const coupledFvPatch& cfvp = refCast<const coupledFvPatch>(fvp);
754 
755  const scalar w = cfvp.owner() ? ownerWeight : neighbourWeight;
756  const scalar v = cfvp.owner() ? neighbourWeight : ownerWeight;
757 
758  syncFieldb[patchi] =
759  w*syncFieldb[patchi] + (flip ? -v : +v)*fieldbNbr[patchi];
760  }
761  }
762 
763  return tsyncFieldb;
764 }
765 
766 
767 template<class Type>
770 (
771  const SurfaceFieldBoundary<Type>& fieldb
772 ) const
773 {
774  const bool isFluxField = isFlux(fieldb[0].internalField());
775 
776  return synchronisedBoundaryField<Type>
777  (
778  fieldb,
779  isFluxField,
780  0.5,
781  0.5
782  );
783 }
784 
785 
786 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
static const DimensionedField< Type, GeoMesh > & null()
Return a null DimensionedField.
Pre-declare SubField and related Field type.
Definition: Field.H:82
Generic GeometricBoundaryField class.
tmp< GeometricBoundaryField > boundaryNeighbourField() const
Return BoundaryField of the values on the other side of couples.
Generic GeometricField class.
PatchField< Type > Patch
Type of the patch field of which the Boundary is composed.
const Boundary & boundaryField() const
Return const-reference to the boundary field.
void clear()
Clear the list, i.e. set size to zero.
Definition: ListI.H:125
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:137
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:147
static bool & parRun()
Is this a parallel run?
Definition: UPstream.H:399
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:272
This boundary condition is not designed to be evaluated; it is assumed that the value is assigned via...
An abstract base class for patches that couple regions of the computational domain e....
virtual bool owner() const
Does this side own the patch ?
Foam::fvBoundaryMesh.
const fvMesh & mesh() const
Return the mesh reference.
tmp< SurfaceFieldBoundary< Type > > conformalNccBoundaryField(const SurfaceFieldBoundary< Type > &fieldb) const
Extract the non-conformal-coupled parts of the boundary field.
tmp< SurfaceFieldBoundary< Type > > nonConformalBoundaryField(const SurfaceFieldBoundary< Type > &nccFieldb, const SurfaceFieldBoundary< Type > &origFieldb) const
Combine non-conformal-coupled and original parts of the.
static tmp< Field< Type > > fieldRMapSum(const Field< Type > &f, const label size, const labelUList &addr)
Reverse-map sum the values of a field.
tmp< SurfaceFieldBoundary< Type > > conformalOrigBoundaryField(const SurfaceFieldBoundary< Type > &fieldb) const
Extract the original parts of the boundary field and store it.
void resizeBoundaryFieldPatchFields(const SurfaceFieldBoundary< label > &polyFacesBf, GeoBoundaryField &fieldBf)
Resize the patch field of a given field to match the mesh.
tmp< SurfaceFieldBoundary< Type > > synchronisedBoundaryField(const SurfaceFieldBoundary< Type > &fieldb, const bool flip, const scalar ownerWeight, const scalar neighbourWeight) const
Synchronise the boundary field by combining corresponding.
tmp< SurfaceFieldBoundary< scalar > > getOrigNccMagSfb() const
Return the total non-conformal area associated with each.
void resizeFieldPatchFields(const SurfaceFieldBoundary< label > &polyFacesBf, GeoField &field)
Resize the patch field of a given field to match the mesh.
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:893
const surfaceScalarField & magSf() const
Return cell face area magnitudes.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:64
virtual label size() const
Return size.
Definition: fvPatch.H:157
virtual bool coupled() const
Return true if this patch is coupled.
Definition: fvPatch.H:163
virtual label start() const
Return start label of this patch in the polyMesh face list.
Definition: fvPatch.H:151
labelList allOrigPatchIDs() const
Return a list of the orig patch indices.
Non-conformal coupled FV patch. As nonConformalFvPatch, but is also coupled to another non-conformal ...
Non-conformal error FV patch. As nonConformalFvPatch. This patch is a non-coupled non-conformal patch...
const labelList & polyFaces() const
Return face face-poly-faces.
const labelList & polyFaces() const
Return face face-poly-faces.
const fvPatch & origPatch() const
Original patch.
label origPatchID() const
Original patch ID.
HashTable< const Type * > lookupClass(const bool strict=false) const
Lookup and return all objects of the given Type.
Traits class for primitives.
Definition: pTraits.H:53
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
#define PostNonConformSurfaceFields(Type, nullArg)
#define PreConformSurfaceFields(Type, nullArg)
#define ResizePatchFields(Type, nullArg)
#define EvaluateVolFields(Type, nullArg)
label patchi
Info<< "Calculating turbulent flame speed field St\n"<< endl;volScalarField St(IOobject("St", runTime.name(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), flameWrinkling->Xi() *Su);multivariateSurfaceInterpolationScheme< scalar >::fieldTable fields
Definition: createFields.H:230
U
Definition: pEqn.H:72
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
VolField< vector > volVectorField
Definition: volFieldsFwd.H:62
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
SurfaceField< scalar > surfaceScalarField
bool isFlux(const DimensionedField< Type, surfaceMesh > &df)
Check if surfaceField is a flux.
Definition: surfaceFields.H:55
FOR_ALL_FIELD_TYPES(DefineContiguousFvWallLocationDataType)
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:52
error FatalError
SurfaceField< vector > surfaceVectorField
const volVectorField & surfaceToVolVelocity(const surfaceVectorField &Uf)
Get the cell velocity field corresponding to a given face velocity, or a.
labelList f(nPoints)
Foam::surfaceFields.
Return the vol-field velocity corresponding to a given surface-field velocity.