oversetFvPatchField.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-2022 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 "volFields.H"
29 #include "cellCellStencil.H"
30 #include "cellCellStencilObject.H"
32 #include "syncTools.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 template<class Type>
38 (
39  const fvPatch& p,
41 )
42 :
43  coupledFvPatchField<Type>(p, iF),
44  oversetPatch_(refCast<const oversetFvPatch>(p)),
45  setHoleCellValue_(false),
46  fluxCorrection_(false),
47  interpolateHoleCellValue_(false),
48  holeCellValue_(pTraits<Type>::min),
49  fringeUpperCoeffs_(),
50  fringeLowerCoeffs_(),
51  fringeFaces_(),
52  zoneId_(-1)
53 {}
54 
55 
56 template<class Type>
58 (
59  const oversetFvPatchField<Type>& ptf,
60  const fvPatch& p,
62  const fvPatchFieldMapper& mapper
63 )
64 :
65  coupledFvPatchField<Type>(ptf, p, iF, mapper),
66  oversetPatch_(refCast<const oversetFvPatch>(p)),
67  setHoleCellValue_(ptf.setHoleCellValue_),
68  fluxCorrection_(ptf.fluxCorrection_),
69  interpolateHoleCellValue_(ptf.interpolateHoleCellValue_),
70  holeCellValue_(ptf.holeCellValue_),
71  fringeUpperCoeffs_(ptf.fringeUpperCoeffs_),
72  fringeLowerCoeffs_(ptf.fringeLowerCoeffs_),
73  fringeFaces_(ptf.fringeFaces_),
74  zoneId_(ptf.zoneId_)
75 {}
76 
77 
78 template<class Type>
80 (
81  const fvPatch& p,
83  const dictionary& dict
84 )
85 :
86  coupledFvPatchField<Type>(p, iF, dict, IOobjectOption::NO_READ),
87  oversetPatch_(refCast<const oversetFvPatch>(p, dict)),
88  setHoleCellValue_(dict.getOrDefault("setHoleCellValue", false)),
89  fluxCorrection_
90  (
91  dict.getOrDefaultCompat
92  (
93  "fluxCorrection",
94  {{"massCorrection", 2206}},
95  false
96  )
97  ),
98  interpolateHoleCellValue_
99  (
100  dict.getOrDefault("interpolateHoleCellValue", false)
101  ),
102  holeCellValue_
103  (
104  setHoleCellValue_
105  ? dict.get<Type>("holeCellValue")
107  ),
108  fringeUpperCoeffs_(),
109  fringeLowerCoeffs_(),
110  fringeFaces_(),
111  zoneId_(dict.getOrDefault<label>("zone", -1))
112 {
113  // Use 'value' supplied, or set to internal field
114  if (!this->readValueEntry(dict))
115  {
117  }
118 }
119 
120 
121 template<class Type>
123 (
124  const oversetFvPatchField<Type>& ptf
125 )
126 :
127  coupledFvPatchField<Type>(ptf),
128  oversetPatch_(ptf.oversetPatch_),
129  setHoleCellValue_(ptf.setHoleCellValue_),
130  fluxCorrection_(ptf.fluxCorrection_),
131  interpolateHoleCellValue_(ptf.interpolateHoleCellValue_),
132  holeCellValue_(ptf.holeCellValue_),
133  fringeUpperCoeffs_(ptf.fringeUpperCoeffs_),
134  fringeLowerCoeffs_(ptf.fringeLowerCoeffs_),
135  fringeFaces_(ptf.fringeFaces_),
136  zoneId_(ptf.zoneId_)
137 {}
138 
139 
140 template<class Type>
142 (
143  const oversetFvPatchField<Type>& ptf,
145 )
146 :
147  coupledFvPatchField<Type>(ptf, iF),
148  oversetPatch_(ptf.oversetPatch_),
149  setHoleCellValue_(ptf.setHoleCellValue_),
150  fluxCorrection_(ptf.fluxCorrection_),
151  interpolateHoleCellValue_(ptf.interpolateHoleCellValue_),
152  holeCellValue_(ptf.holeCellValue_),
153  fringeUpperCoeffs_(ptf.fringeUpperCoeffs_),
154  fringeLowerCoeffs_(ptf.fringeLowerCoeffs_),
155  fringeFaces_(ptf.fringeFaces_),
156  zoneId_(ptf.zoneId_)
157 {}
158 
159 
160 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
161 
162 template<class Type>
164 (
165  const fvMatrix<Type>& matrix
166 )
167 {
168  const fvMesh& mesh = this->internalField().mesh();
169 
172  const labelList& zoneID = overlap.zoneID();
173  const labelUList& own = mesh.owner();
174  const labelUList& nei = mesh.neighbour();
175 
177 
178  label fringesFaces = 0;
179 
180  forAll(own, facei)
181  {
182  const label zonei = zoneID[own[facei]];
183 
184  const label ownType = cellTypes[own[facei]];
185  const label neiType = cellTypes[nei[facei]];
186 
187  const bool ownCalc =
188  (ownType == cellCellStencil::CALCULATED)
189  && (neiType == cellCellStencil::INTERPOLATED);
190 
191  const bool neiCalc =
192  (ownType == cellCellStencil::INTERPOLATED)
193  && (neiType == cellCellStencil::CALCULATED);
194 
195  const bool ownNei = (ownCalc || neiCalc);
196 
197  if
198  (
199  (ownNei && (zonei == zoneId_)) || (ownNei && (zoneId_ == -1))
200  )
201  {
202  fringesFaces++;
203  }
204  }
205 
206  const fvPatchList& patches = mesh.boundary();
207 
208  labelList neiCellTypes;
210  {
211  forAll(patches, patchi)
212  {
213  const fvPatch& curPatch = patches[patchi];
214 
215  const labelUList& fc = curPatch.faceCells();
216 
217  const label start = curPatch.start();
218 
219  forAll(fc, i)
220  {
221  const label facei = start + i;
222  const label celli = fc[i];
223  const label ownType = cellTypes[celli];
224  const label neiType = neiCellTypes[facei-mesh.nInternalFaces()];
225 
226  const label zonei = zoneID[celli];
227 
228  const bool ownCalc =
229  (ownType == cellCellStencil::CALCULATED)
230  && (neiType == cellCellStencil::INTERPOLATED);
231 
232 
233  if (ownCalc && (zonei == zoneId_))
234  {
235  fringesFaces++;
236  }
237  }
238  }
239  }
240  fringeUpperCoeffs_.setSize(fringesFaces, Zero);
241  fringeLowerCoeffs_.setSize(fringesFaces, Zero);
242  fringeFaces_.setSize(fringesFaces, -1);
243 
244  const scalarField& upper = matrix.upper();
245  const scalarField& lower = matrix.lower();
246 
247  fringesFaces = 0;
248  forAll(own, facei)
249  {
250  const label zonei = zoneID[own[facei]];
251 
252  const label ownType = cellTypes[own[facei]];
253  const label neiType = cellTypes[nei[facei]];
254 
255  const bool ownCalc =
256  (ownType == cellCellStencil::CALCULATED)
257  && (neiType == cellCellStencil::INTERPOLATED);
258 
259  const bool neiCalc =
260  (ownType == cellCellStencil::INTERPOLATED)
261  && (neiType == cellCellStencil::CALCULATED);
262 
263  const bool ownNei = (ownCalc || neiCalc);
264 
265  if
266  (
267  (ownNei && (zonei == zoneId_)) || (ownNei && (zoneId_ == -1))
268  )
269  {
270  fringeUpperCoeffs_[fringesFaces] = upper[facei];
271  fringeLowerCoeffs_[fringesFaces] = lower[facei];
272  fringeFaces_[fringesFaces] = facei;
273  fringesFaces++;
274  }
275  }
276 
277  forAll(boundaryMesh, patchi)
278  {
279  const polyPatch& p = boundaryMesh[patchi];
280 
281  if (isA<coupledPolyPatch>(p))
282  {
283  const labelUList& fc = p.faceCells();
284  const label start = p.start();
285 
286  forAll(fc, i)
287  {
288  const label facei = start + i;
289  const label celli = fc[i];
290  const label ownType = cellTypes[celli];
291  const label neiType = neiCellTypes[facei-mesh.nInternalFaces()];
292 
293  const label zonei = zoneID[celli];
294 
295  const bool ownCalc =
296  (ownType == cellCellStencil::CALCULATED)
297  && (neiType == cellCellStencil::INTERPOLATED);
298 
299  const bool neiCalc =
300  (neiType == cellCellStencil::CALCULATED)
301  && (ownType == cellCellStencil::INTERPOLATED);
302 
303  if ((ownCalc||neiCalc) && (zonei == zoneId_))
304  {
305  fringeLowerCoeffs_[fringesFaces] =
306  component
307  (
308  matrix.internalCoeffs()[patchi][facei],
309  0
310  );
311 
312  fringeUpperCoeffs_[fringesFaces] =
313  component
314  (
315  matrix.boundaryCoeffs()[patchi][facei],
316  0
317  );
318 
319  fringeFaces_[fringesFaces] = facei;
320 
321  fringesFaces++;
322  }
323  }
324  }
325  }
326 }
327 
328 
329 template<class Type>
331 (
332  const fvMatrix<Type>& matrix,
333  const surfaceScalarField& phi
334 ) const
335 {
336  scalar massIn = 0;
337  scalar phiIn = 0;
338 
339  const Field<Type>& psi = matrix.psi();
340 
341  const scalarField& upper = matrix.upper();
342  const scalarField& lower = matrix.lower();
343 
344  if (this->oversetPatch_.master())
345  {
346  const fvMesh& mesh = this->internalField().mesh();
347  const cellCellStencilObject& overlap = Stencil::New(mesh);
349  const labelList& zoneID = overlap.zoneID();
350 
351  // Check all faces on the outside of interpolated cells
352  const labelUList& own = mesh.owner();
353  const labelUList& nei = mesh.neighbour();
354 
355  label fringesFaces = 0;
356  forAll(own, facei)
357  {
358  const label zonei = zoneID[own[facei]];
359 
360  const label ownType = cellTypes[own[facei]];
361  const label neiType = cellTypes[nei[facei]];
362 
363  const bool ownCalc =
364  (ownType == cellCellStencil::CALCULATED)
365  && (neiType == cellCellStencil::INTERPOLATED);
366 
367  const bool neiCalc =
368  (ownType == cellCellStencil::INTERPOLATED)
369  && (neiType == cellCellStencil::CALCULATED);
370 
371  const bool ownNei = (ownCalc || neiCalc);
372 
373  if
374  (
375  (ownNei && (zonei == zoneId_)) || (ownNei && (zoneId_ == -1))
376  )
377  {
378  const label fringei = fringeFaces_[fringesFaces];
379 
380  // Get fringe upper/lower coeffs
381  //const scalar& ufc = fringeUpperCoeffs_[fringesFaces];
382  //const scalar& lfc = fringeLowerCoeffs_[fringesFaces];
383 
384  const scalar& ufc = upper[fringei];
385  const scalar& lfc = lower[fringei];
386 
387  const Type curFlux =
388  ufc*psi[nei[fringei]] - lfc*psi[own[fringei]];
389 
390  if (neiCalc)
391  {
392  phiIn -= phi[fringei];
393  massIn -= curFlux;
394  }
395  else
396  {
397  phiIn += phi[fringei];
398  massIn += curFlux;
399  }
400  fringesFaces++;
401  }
402  }
403  }
404 
405  reduce(massIn, sumOp<scalar>());
406  reduce(phiIn, sumOp<scalar>());
407 
408  Info << " gSum(phi) on fringes " << phiIn << endl;
409  Info << " gSum(p.flux) on fringes " << massIn << endl;
410 }
411 
412 
413 template<class Type>
415 (
417  const lduAddressing& lduAddr,
418  solveScalarField& result
419 ) const
420 {
421  const fvMesh& mesh = this->internalField().mesh();
422 
425  const labelList& zoneID = overlap.zoneID();
426 
427  // Pass-1: accumulate all fluxes, calculate correction factor
428 
429  scalarField interpolatedCoeffs(fringeUpperCoeffs_.size(), Zero);
430 
431  // Options for scaling corrections
432  scalar massIn = 0;
433  scalar offDiagCoeffs = 0;
434  labelField facePerCell(cellTypes.size(), 0);
435 
436  // Check all faces on the outside of interpolated cells
437  const labelUList& own = mesh.owner();
438  const labelUList& nei = mesh.neighbour();
439 
440  label fringesFaces = 0;
441  {
442  forAll(own, facei)
443  {
444  const label zonei = zoneID[own[facei]];
445 
446  const label ownType = cellTypes[own[facei]];
447  const label neiType = cellTypes[nei[facei]];
448 
449  const bool ownCalc =
450  (ownType == cellCellStencil::CALCULATED)
451  && (neiType == cellCellStencil::INTERPOLATED);
452 
453  const bool neiCalc =
454  (ownType == cellCellStencil::INTERPOLATED)
455  && (neiType == cellCellStencil::CALCULATED);
456 
457  const bool ownNei = (ownCalc || neiCalc);
458 
459  if
460  (
461  (ownNei && (zonei == zoneId_)) || (ownNei && (zoneId_ == -1))
462  )
463  {
464  // Get fringe upper/lower coeffs
465  const scalar& ufc = fringeUpperCoeffs_[fringesFaces];
466  const scalar& lfc = fringeLowerCoeffs_[fringesFaces];
467 
468  const scalar curFlux =
469  ufc*psi[nei[facei]] - lfc*psi[own[facei]];
470 
471  if (neiCalc) // interpolated is owner
472  {
473  massIn -= curFlux;
474  offDiagCoeffs += lfc;
475  facePerCell[own[facei]]++;
476  }
477  else
478  {
479  massIn += curFlux;
480  offDiagCoeffs += ufc;
481  facePerCell[nei[facei]]++;
482  }
483  fringesFaces++;
484  }
485  }
486  }
487 
488  scalarField weights(facePerCell.size(), scalar(1));
489  forAll (weights, celli)
490  {
491  if (facePerCell[celli] > 1)
492  {
493  weights[celli] = scalar(1)/facePerCell[celli];
494  }
495  }
496 
497  // Check all coupled faces on the outside of interpolated cells
498  labelList neiCellTypes;
500 
501  const fvPatchList& boundaryMesh = mesh.boundary();
502 
503  forAll(boundaryMesh, patchi)
504  {
505  const polyPatch& p = mesh.boundaryMesh()[patchi];
506 
507  if (isA<coupledPolyPatch>(p))
508  {
509  const auto& coupledPatch = refCast<const coupledPolyPatch>(p);
510 
511  const labelUList& fc = p.faceCells();
512  const label start = p.start();
513 
514  forAll(fc, i)
515  {
516  const label facei = start + i;
517  const label celli = fc[i];
518  const label ownType = cellTypes[celli];
519  const label neiType = neiCellTypes[facei-mesh.nInternalFaces()];
520 
521  const label zonei = zoneID[celli];
522 
523  const bool ownCalc =
524  (ownType == cellCellStencil::CALCULATED)
525  && (neiType == cellCellStencil::INTERPOLATED);
526 
527  const bool neiCalc =
528  (neiType == cellCellStencil::CALCULATED)
529  && (ownType == cellCellStencil::INTERPOLATED);
530 
531  if ((ownCalc||neiCalc) && (zonei == zoneId_))
532  {
533  const scalar psiOwn = psi[celli];
534  const scalar& lfc = fringeLowerCoeffs_[fringesFaces];
535  const scalar curFlux = lfc*psiOwn;
536 
537  if (ownCalc)
538  {
539  massIn -= curFlux;
540 
541  if (coupledPatch.owner())
542  {
543  offDiagCoeffs -= lfc;
544  }
545 
546  fringesFaces++;
547  }
548  else
549  {
550  massIn += curFlux;
551 
552  if (coupledPatch.owner())
553  {
554  offDiagCoeffs -= lfc;
555  }
556 
557  fringesFaces++;
558  }
559  }
560  }
561  }
562  }
563 
564  reduce(massIn, sumOp<scalar>());
565  reduce(offDiagCoeffs, sumOp<scalar>());
566 
567  scalar psiCorr = -massIn/offDiagCoeffs;
568 
569  forAll (cellTypes, celli)
570  {
571  const bool bInter = (cellTypes[celli] == cellCellStencil::INTERPOLATED);
572  const label zonei = zoneID[celli];
573  if
574  (
575  (bInter && (zonei == zoneId_)) ||(bInter && (zoneId_ == -1))
576  )
577  {
578  psi[celli] += psiCorr;
579  }
580  }
581 }
582 
583 
584 template<class Type>
586 patchNeighbourField() const
587 {
588  return tmp<Field<Type>>::New(this->size(), Zero);
589 }
590 
591 
592 template<class Type>
594 (
595  const Pstream::commsTypes commsType
596 )
597 {
598  if (this->oversetPatch_.master())
599  {
600  // Trigger interpolation
601  const fvMesh& mesh = this->internalField().mesh();
603  const word& fldName = this->internalField().name();
604 
605  if (&mesh.lduAddr() != &mesh.fvMesh::lduAddr())
606  {
607  // Running extended addressing. Flux correction already done
608  // in the linear solver (through the initUpdateInterfaceMatrix
609  // method below)
610  if (debug)
611  {
612  Info<< "Skipping overset interpolation for solved-for field "
613  << fldName << endl;
614  }
615  }
616  else if (!fvSchemes.found("oversetInterpolation"))
617  {
618  IOWarningInFunction(fvSchemes)
619  << "Missing required dictionary entry"
620  << " 'oversetInterpolation'"
621  << ". Skipping overset interpolation for field "
622  << fldName << endl;
623  }
624  else if (fvSchemes.found("oversetInterpolationRequired"))
625  {
626  // Backwards compatibility mode: only interpolate what is
627  // explicitly mentioned
628 
629  if (fvSchemes.found("oversetInterpolationSuppressed"))
630  {
631  FatalIOErrorInFunction(fvSchemes)
632  << "Cannot have both dictionary entry"
633  << " 'oversetInterpolationSuppresed' and "
634  << " 'oversetInterpolationRequired' for field "
635  << fldName << exit(FatalIOError);
636  }
637  const dictionary& intDict = fvSchemes.subDict
638  (
639  "oversetInterpolationRequired"
640  );
641  if (intDict.found(fldName))
642  {
643  if (debug)
644  {
645  Info<< "Interpolating field " << fldName << endl;
646  }
647 
648  // Interpolate without bc update (since would trigger infinite
649  // recursion back to oversetFvPatchField<Type>::evaluate)
650  // The only problem is bcs that use the new cell values
651  // (e.g. zeroGradient, processor). These need to appear -after-
652  // the 'overset' bc.
654  (
655  const_cast<Field<Type>&>
656  (
657  this->primitiveField()
658  )
659  );
660  }
661  else if (debug)
662  {
663  Info<< "Skipping overset interpolation for field "
664  << fldName << endl;
665  }
666  }
667  else
668  {
669  const dictionary* dictPtr
670  (
671  fvSchemes.findDict("oversetInterpolationSuppressed")
672  );
673 
674  const wordHashSet& suppress =
676 
677  bool skipInterpolate = suppress.found(fldName);
678 
679  if (dictPtr)
680  {
681  skipInterpolate =
682  skipInterpolate
683  || dictPtr->found(fldName);
684  }
685 
686  if (skipInterpolate)
687  {
688  if (debug)
689  {
690  Info<< "Skipping suppressed overset interpolation"
691  << " for field " << fldName << endl;
692  }
693  }
694  else
695  {
696  if (debug)
697  {
698  Info<< "Interpolating non-suppressed field " << fldName
699  << endl;
700  }
701 
702  // Interpolate without bc update (since would trigger infinite
703  // recursion back to oversetFvPatchField<Type>::evaluate)
704  // The only problem is bcs that use the new cell values
705  // (e.g. zeroGradient, processor). These need to appear -after-
706  // the 'overset' bc.
707 
708  const cellCellStencil& overlap = Stencil::New(mesh);
709 
710  Field<Type>& fld =
711  const_cast<Field<Type>&>(this->primitiveField());
712 
713  // tbd: different weights for different variables ...
715  (
716  fld,
717  mesh,
718  overlap,
720  );
721 
722  if (this->setHoleCellValue_)
723  {
724  const labelUList& types = overlap.cellTypes();
725  label nConstrained = 0;
726  forAll(types, celli)
727  {
728  const label cType = types[celli];
729  if
730  (
731  cType == cellCellStencil::HOLE
732  || cType == cellCellStencil::SPECIAL
733  )
734  {
735  fld[celli] = this->holeCellValue_;
736  nConstrained++;
737  }
738  }
739 
740  if (debug)
741  {
742  Pout<< FUNCTION_NAME << " field:" << fldName
743  << " patch:" << this->oversetPatch_.name()
744  << " set:" << nConstrained << " cells to:"
745  << this->holeCellValue_ << endl;
746  }
747  }
748  }
749  }
750  }
753 }
754 
755 
756 template<class Type>
758 (
759  fvMatrix<Type>& matrix
760 )
761 {
762  if (this->manipulatedMatrix())
763  {
764  return;
765  }
766 
767  const oversetFvPatch& ovp = this->oversetPatch_;
768 
769  if (ovp.master())
770  {
771  if (fluxCorrection_ || (debug & 2))
772  {
773  storeFringeCoefficients(matrix);
774  }
775 
776  // Trigger interpolation
777  const fvMesh& mesh = this->internalField().mesh();
778  const word& fldName = this->internalField().name();
779 
780  // Try to find out if the solve routine comes from the unadapted mesh
781  // TBD. This should be cleaner.
782  if (&mesh.lduAddr() == &mesh.fvMesh::lduAddr())
783  {
784  return;
785  }
786 
787  if (debug)
788  {
789  Pout<< FUNCTION_NAME << " field:" << fldName
790  << " patch:" << ovp.name() << endl;
791  }
792 
793 
794  // Calculate stabilised diagonal as normalisation for interpolation
795  const scalarField norm
796  (
797  dynamic_cast<const oversetFvMeshBase&>(mesh).normalisation(matrix)
798  );
799 
800  const cellCellStencil& overlap = Stencil::New(mesh);
801  const labelUList& types = overlap.cellTypes();
802  const labelListList& stencil = overlap.cellStencil();
803 
804  dynamic_cast<const oversetFvMeshBase&>(mesh).addInterpolation
805  (
806  matrix,
807  norm,
808  this->setHoleCellValue_,
809  this->holeCellValue_
810  );
811 
812  if (debug)
813  {
814  pointField allCs(mesh.cellCentres());
815  const mapDistribute& map = overlap.cellInterpolationMap();
816  map.distribute(allCs, false, UPstream::msgType()+1);
817 
818  // Make sure we don't interpolate from a hole
819 
820  scalarField marker(this->primitiveField().size(), 0);
821  forAll(types, celli)
822  {
823  if (types[celli] == cellCellStencil::HOLE)
824  {
825  marker[celli] = 1.0;
826  }
827  }
829  (
830  marker,
831  mesh,
832  overlap,
834  );
835 
836  forAll(marker, celli)
837  {
838  if
839  (
840  types[celli] == cellCellStencil::INTERPOLATED
841  && marker[celli] > SMALL
842  )
843  {
845  << " field:" << fldName
846  << " patch:" << ovp.name()
847  << " found:" << celli
848  << " at:" << mesh.cellCentres()[celli]
849  << " donorSlots:" << stencil[celli]
850  << " at:"
851  << UIndirectList<point>(allCs, stencil[celli])
852  << " amount-of-hole:" << marker[celli]
853  << endl;
854  }
855  }
856 
857  // Make sure we don't have matrix coefficients for interpolated
858  // or hole cells
859 
860  const lduAddressing& addr = mesh.lduAddr();
861  const labelUList& upperAddr = addr.upperAddr();
862  const labelUList& lowerAddr = addr.lowerAddr();
863  const scalarField& lower = matrix.lower();
864  const scalarField& upper = matrix.upper();
865 
866  forAll(lowerAddr, facei)
867  {
868  const label l = lowerAddr[facei];
869  const bool lHole = (types[l] == cellCellStencil::HOLE);
870  const label u = upperAddr[facei];
871  const bool uHole = (types[u] == cellCellStencil::HOLE);
872 
873  if
874  (
875  (lHole && upper[facei] != 0.0)
876  || (uHole && lower[facei] != 0.0)
877  )
878  {
880  << "Hole-neighbouring face:" << facei
881  << " lower:" << l
882  << " type:" << types[l]
883  << " coeff:" << lower[facei]
884  << " upper:" << upperAddr[facei]
885  << " type:" << types[u]
886  << " coeff:" << upper[facei]
887  << exit(FatalError);
888  }
889 
890 
891  // Empty donor list: treat like hole but still allow to
892  // influence neighbouring domains
893  const bool lEmpty =
894  (
896  && stencil[l].empty()
897  );
898  const bool uEmpty =
899  (
901  && stencil[u].empty()
902  );
903 
904  if
905  (
906  (lEmpty && upper[facei] != 0.0)
907  || (uEmpty && lower[facei] != 0.0)
908  )
909  {
911  << "Still connected face:" << facei << " lower:" << l
912  << " type:" << types[l]
913  << " coeff:" << lower[facei]
914  << " upper:" << u
915  << " type:" << types[u]
916  << " coeff:" << upper[facei]
917  << exit(FatalError);
918  }
919  }
920 
921  forAll(matrix.internalCoeffs(), patchi)
922  {
923  const labelUList& fc = addr.patchAddr(patchi);
924  const Field<Type>& bouCoeffs = matrix.boundaryCoeffs()[patchi];
925 
926  forAll(fc, i)
927  {
928  const label celli = fc[i];
929  const bool lHole = (types[celli] == cellCellStencil::HOLE);
930 
931  if (lHole && bouCoeffs[i] != pTraits<Type>::zero)
932  {
934  << "Patch:" << patchi
935  << " patchFace:" << i
936  << " lower:" << celli
937  << " type:" << types[celli]
938  << " bouCoeff:" << bouCoeffs[i]
939  << exit(FatalError);
940  }
941 
942  // Check whether this is influenced by neighbouring domains
943  const bool lEmpty =
944  (
945  types[celli] == cellCellStencil::INTERPOLATED
946  && stencil[celli].empty()
947  );
948 
949  if (lEmpty && bouCoeffs[i] != pTraits<Type>::zero)
950  {
952  << "Patch:" << patchi
953  << " patchFace:" << i
954  << " lower:" << celli
955  << " type:" << types[celli]
956  << " bouCoeff:" << bouCoeffs[i]
957  << exit(FatalError);
958  }
959  }
960  }
961 
962 
963  // Make sure that diagonal is non-zero. Note: should add
964  // boundaryCoeff ...
965  const FieldField<Field, Type>& internalCoeffs =
966  matrix.internalCoeffs();
967 
968  for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; ++cmpt)
969  {
970  // Replacement for m.addBoundaryDiag(norm, cmpt);
971  scalarField diag(matrix.diag());
972  forAll(internalCoeffs, patchi)
973  {
974  const labelUList& fc = addr.patchAddr(patchi);
975  const Field<Type>& intCoeffs = internalCoeffs[patchi];
976  const scalarField cmptCoeffs(intCoeffs.component(cmpt));
977  forAll(fc, i)
978  {
979  diag[fc[i]] += cmptCoeffs[i];
980  }
981  }
982 
983  forAll(diag, celli)
984  {
985  if (mag(diag[celli]) < SMALL)
986  {
988  << "Patch:" << ovp.name()
989  << " cell:" << celli
990  << " at:" << mesh.cellCentres()[celli]
991  << " diag:" << diag[celli]
992  << exit(FatalError);
993  }
994  }
995  }
996  }
997  }
1000 }
1001 
1002 
1003 template<class Type>
1005 (
1006  solveScalarField& result,
1007  const bool add,
1008  const lduAddressing& lduAddr,
1009  const label interfacei,
1010  const solveScalarField& psiInternal,
1011  const scalarField& coeffs,
1012  const direction cmpt,
1013  const Pstream::commsTypes commsType
1014 ) const
1015 {
1016 }
1017 
1018 
1019 template<class Type>
1021 (
1022  solveScalarField& result,
1023  const bool add,
1024  const lduAddressing& lduAddr,
1025  const label patchId,
1026  const solveScalarField& psiInternal,
1027  const scalarField& coeffs,
1028  const direction,
1029  const Pstream::commsTypes commsType
1030 ) const
1031 {
1032  auto& psi = const_cast<solveScalarField&>(psiInternal);
1033 
1034  if (fluxCorrection_ && this->oversetPatch_.master())
1035  {
1036  adjustPsi(psi, lduAddr, result);
1037  }
1038 }
1039 
1040 
1041 template<class Type>
1042 void Foam::oversetFvPatchField<Type>::write(Ostream& os) const
1043 {
1045 
1046  if (this->setHoleCellValue_)
1047  {
1048  os.writeEntry("setHoleCellValue", setHoleCellValue_);
1049  os.writeEntry("holeCellValue", holeCellValue_);
1051  (
1052  "interpolateHoleCellValue",
1053  false,
1054  interpolateHoleCellValue_
1055  );
1056  }
1058  (
1059  "fluxCorrection",
1060  false,
1061  fluxCorrection_
1062  );
1064  (
1065  "zone",
1066  label(-1),
1067  zoneId_
1068  );
1069 }
1070 
1071 
1072 // ************************************************************************* //
label patchId(-1)
const scalarField & diag() const
Definition: lduMatrix.C:163
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
uint8_t direction
Definition: direction.H:46
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
commsTypes
Communications types.
Definition: UPstream.H:77
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
void fringeFlux(const fvMatrix< Type > &m, const surfaceScalarField &phi) const
Calculate patch flux (helper function). Requires.
bool found(const Key &key) const
Same as contains()
Definition: HashTable.H:1370
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute List data using default commsType, default flip/negate operator.
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData, const bool parRun=UPstream::parRun())
Extract and swap to obtain neighbour cell values for all boundary faces.
static const cellCellStencilObject & New(const fvMesh &mesh, Args &&... args)
Get existing or create MeshObject registered with typeName.
Definition: MeshObject.C:53
virtual void updateInterfaceMatrix(solveScalarField &result, const bool add, const lduAddressing &lduAddr, const label patchId, const solveScalarField &psiInternal, const scalarField &coeffs, const direction, const Pstream::commsTypes commsType) const
Update result field based on interface functionality.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition: typeInfo.H:172
virtual void initEvaluate(const Pstream::commsTypes commsType)
Initialise the evaluation of the patch field.
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface...
Definition: boundaryMesh.H:58
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition: fvPatch.H:70
virtual tmp< Field< Type > > patchNeighbourField() const
Return neighbour field. Dummy.
A traits class, which is primarily used for primitives and vector-space.
Definition: pTraits.H: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.
Patch for indicating interpolated boundaries (in overset meshes).
void extrapolateInternal()
Assign the patch field from the internal field.
Definition: fvPatchField.C:62
string upper(const std::string &s)
Return string copy transformed with std::toupper on each character.
Definition: stringOps.C:1187
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1252
::Foam::direction nComponents(const expressions::valueTypeCode) noexcept
The number of components associated with given valueTypeCode.
Definition: exprTraits.C:40
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
virtual void interpolate(volScalarField &) const
Interpolate interpolationCells only.
Definition: fvMesh.H:459
static void interpolate(Field< T > &psi, const fvMesh &mesh, const cellCellStencil &overlap, const List< scalarList > &wghts)
Interpolation of acceptor cells from donor cells.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
virtual void initEvaluate(const Pstream::commsTypes commsType)
Initialise the evaluation of the patch field.
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
UList< label > labelUList
A UList of labels.
Definition: UList.H:78
virtual const labelUList & cellTypes() const
Return the cell type list.
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:104
const labelUList & neighbour() const
Internal face neighbour.
Definition: fvMesh.H:580
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
virtual void write(Ostream &) const
Write includes "value" entry.
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
void storeFringeCoefficients(const fvMatrix< Type > &matrix)
Store fringe coefficients and faces.
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
void setSize(const label n)
Alias for resize()
Definition: List.H:320
dynamicFvMesh & mesh
virtual const labelListList & cellStencil() const
Per interpolated cell the neighbour cells (in terms of slots as.
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
A class for handling words, derived from Foam::string.
Definition: word.H:63
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition: fvMesh.C:752
A FieldMapper for finite-volume patch fields.
const cellCellStencilObject & overlap
Definition: correctPhi.H:57
virtual const List< scalarList > & cellInterpolationWeights() const
Weights for cellStencil.
const scalarField & lower() const
Definition: lduMatrix.C:268
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:64
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
Definition: fvPatchField.C:342
const dictionary & schemesDict() const
The entire dictionary or the optional "select" sub-dictionary.
virtual void initInterfaceMatrixUpdate(Field< Type > &, const bool add, const lduAddressing &, const label interfacei, const Field< Type > &, const scalarField &, const Pstream::commsTypes commsType) const
Initialise neighbour matrix update.
label nInternalFaces() const noexcept
Number of internal faces.
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition: HashSet.H:73
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
const vectorField & cellCentres() const
Ostream & writeEntryIfDifferent(const word &key, const T &value1, const T &value2)
Write a keyword/value entry only when the two values differ.
Definition: Ostream.H:336
Abstract base class for coupled patches.
const labelList & cellTypes
Definition: setCellMask.H:27
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO...
oversetFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
void adjustPsi(solveScalarField &psi, const lduAddressing &lduAddr, solveScalarField &result) const
Adjust psi for mass correction. Requires storeFringeCoefficients.
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
int debug
Static debugging option.
OBJstream os(runTime.globalPath()/outputName)
#define FUNCTION_NAME
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition: fvMesh.H:572
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))
const scalarField & upper() const
Definition: lduMatrix.C:203
string lower(const std::string &s)
Return string copy transformed with std::tolower on each character.
Definition: stringOps.C:1171
#define WarningInFunction
Report a warning using Foam::Warning.
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:637
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition: fvSchemes.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
virtual const mapDistribute & cellInterpolationMap() const
Return a communication schedule.
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
const polyBoundaryMesh & patches
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
messageStream Info
Information stream (stdout output on master, null elsewhere)
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition: fvMesh.H:395
const FieldField< Field, Type > & internalCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for internal cells
Definition: fvMatrix.H:547
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
const volScalarField & psi
The class contains the addressing required by the lduMatrix: upper, lower and losort.
const FieldField< Field, Type > & boundaryCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for boundary cells
Definition: fvMatrix.H:565
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
virtual void write(Ostream &os) const
Write.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
List< label > labelList
A List of labels.
Definition: List.H:62
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
static const labelIOList & zoneID(const fvMesh &)
Helper: get reference to registered zoneID. Loads volScalarField.
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
virtual const wordHashSet & nonInterpolatedFields() const
Return the names of any (stencil or mesh specific) fields that.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
void reduce(T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) using linear/tree communication schedule.
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
Definition: fvPatch.H:59
Boundary condition for use on overset patches. To be run in combination with special dynamicFvMesh ty...
const dictionary * findDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary pointer if present (and it is a dictionary) otherwise return nullptr...
Definition: dictionaryI.H:124
IOerror FatalIOError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL IO ERROR&#39; header text and ...
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127