Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
10License
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"
31#include "oversetFvMeshBase.H"
32#include "syncTools.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36template<class Type>
38(
39 const fvPatch& p,
41)
42:
43 coupledFvPatchField<Type>(p, iF),
45 setHoleCellValue_(false),
46 fluxCorrection_(false),
52 zoneId_(-1)
53{}
54
55
56template<class Type>
58(
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_),
74 zoneId_(ptf.zoneId_)
75{}
76
77
78template<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))
117 }
118}
119
120
121template<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_),
136 zoneId_(ptf.zoneId_)
137{}
138
139
140template<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_),
157{}
158
159
160// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
161
162template<class Type>
164(
165 const fvMatrix<Type>& matrix
166)
167{
168 const fvMesh& mesh = this->internalField().mesh();
169
171 const labelUList& cellTypes = overlap.cellTypes();
172 const labelUList& zoneID = overlap.zoneID();
173 const labelUList& own = mesh.owner();
174 const labelUList& nei = mesh.neighbour();
175
176 const polyBoundaryMesh& boundaryMesh = mesh.boundaryMesh();
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 =
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 =
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
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] =
307 (
308 matrix.internalCoeffs()[patchi][facei],
309 0
310 );
311
312 fringeUpperCoeffs_[fringesFaces] =
314 (
315 matrix.boundaryCoeffs()[patchi][facei],
316 0
317 );
318
319 fringeFaces_[fringesFaces] = facei;
320
321 fringesFaces++;
322 }
324 }
325 }
326}
327
328
329template<class Type>
331(
332 const fvMatrix<Type>& matrix,
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();
349 const labelUList& 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 =
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>());
408 Info << " gSum(phi) on fringes " << phiIn << endl;
409 Info << " gSum(p.flux) on fringes " << massIn << endl;
410}
411
412
413template<class Type>
415(
417 const lduAddressing& lduAddr,
418 solveScalarField& result
419) const
420{
421 const fvMesh& mesh = this->internalField().mesh();
422
424 const labelUList& cellTypes = overlap.cellTypes();
425 const labelUList& 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 =
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
502
503 forAll(boundaryMesh, patchi)
504 {
505 const polyPatch& p = mesh.boundaryMesh()[patchi];
506
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
584template<class Type>
586(
587 const Pstream::commsTypes commsType
588)
589{
590 if (this->oversetPatch_.master())
591 {
592 // Trigger interpolation
593 const fvMesh& mesh = this->internalField().mesh();
595 const word& fldName = this->internalField().name();
596
597 if (&mesh.lduAddr() != &mesh.fvMesh::lduAddr())
598 {
599 // Running extended addressing. Flux correction already done
600 // in the linear solver (through the initUpdateInterfaceMatrix
601 // method below)
602 if (debug)
603 {
604 Info<< "Skipping overset interpolation for solved-for field "
605 << fldName << endl;
606 }
607 }
608 else if (!fvSchemes.found("oversetInterpolation"))
609 {
611 << "Missing required dictionary entry"
612 << " 'oversetInterpolation'"
613 << ". Skipping overset interpolation for field "
614 << fldName << endl;
615 }
616 else if (fvSchemes.found("oversetInterpolationRequired"))
617 {
618 // Backwards compatibility mode: only interpolate what is
619 // explicitly mentioned
620
621 if (fvSchemes.found("oversetInterpolationSuppressed"))
622 {
624 << "Cannot have both dictionary entry"
625 << " 'oversetInterpolationSuppresed' and "
626 << " 'oversetInterpolationRequired' for field "
627 << fldName << exit(FatalIOError);
628 }
629 const dictionary& intDict = fvSchemes.subDict
630 (
631 "oversetInterpolationRequired"
632 );
633 if (intDict.found(fldName))
634 {
635 if (debug)
636 {
637 Info<< "Interpolating field " << fldName << endl;
638 }
639
640 // Interpolate without bc update (since would trigger infinite
641 // recursion back to oversetFvPatchField<Type>::evaluate)
642 // The only problem is bcs that use the new cell values
643 // (e.g. zeroGradient, processor). These need to appear -after-
644 // the 'overset' bc.
646 (
647 const_cast<Field<Type>&>
648 (
649 this->primitiveField()
650 )
651 );
652 }
653 else if (debug)
654 {
655 Info<< "Skipping overset interpolation for field "
656 << fldName << endl;
657 }
658 }
659 else
660 {
661 const dictionary* dictPtr
662 (
663 fvSchemes.findDict("oversetInterpolationSuppressed")
664 );
665
666 const wordHashSet& suppress =
668
669 bool skipInterpolate = suppress.found(fldName);
670
671 if (dictPtr)
672 {
673 skipInterpolate =
674 skipInterpolate
675 || dictPtr->found(fldName);
676 }
677
678 if (skipInterpolate)
679 {
680 if (debug)
681 {
682 Info<< "Skipping suppressed overset interpolation"
683 << " for field " << fldName << endl;
684 }
685 }
686 else
687 {
688 if (debug)
689 {
690 Info<< "Interpolating non-suppressed field " << fldName
691 << endl;
692 }
693
694 // Interpolate without bc update (since would trigger infinite
695 // recursion back to oversetFvPatchField<Type>::evaluate)
696 // The only problem is bcs that use the new cell values
697 // (e.g. zeroGradient, processor). These need to appear -after-
698 // the 'overset' bc.
699
701
703 const_cast<Field<Type>&>(this->primitiveField());
704
705 // tbd: different weights for different variables ...
707 (
708 fld,
709 mesh,
710 overlap,
712 );
713
714 if (this->setHoleCellValue_)
715 {
716 const labelUList& types = overlap.cellTypes();
717 label nConstrained = 0;
718 forAll(types, celli)
719 {
720 const label cType = types[celli];
721 if
722 (
723 cType == cellCellStencil::HOLE
724 || cType == cellCellStencil::SPECIAL
725 )
726 {
727 fld[celli] = this->holeCellValue_;
728 nConstrained++;
729 }
730 }
731
732 if (debug)
733 {
734 Pout<< FUNCTION_NAME << " field:" << fldName
735 << " patch:" << this->oversetPatch_.name()
736 << " set:" << nConstrained << " cells to:"
737 << this->holeCellValue_ << endl;
738 }
739 }
740 }
741 }
743
745}
746
747
748template<class Type>
750(
751 fvMatrix<Type>& matrix
752)
753{
754 if (this->manipulatedMatrix())
755 {
756 return;
757 }
758
759 const oversetFvPatch& ovp = this->oversetPatch_;
760
761 if (ovp.master())
762 {
763 if (fluxCorrection_ || (debug & 2))
764 {
765 storeFringeCoefficients(matrix);
766 }
767
768 // Trigger interpolation
769 const fvMesh& mesh = this->internalField().mesh();
770 const word& fldName = this->internalField().name();
771
772 // Try to find out if the solve routine comes from the unadapted mesh
773 // TBD. This should be cleaner.
774 if (&mesh.lduAddr() == &mesh.fvMesh::lduAddr())
775 {
776 return;
777 }
778
779 if (debug)
780 {
781 Pout<< FUNCTION_NAME << " field:" << fldName
782 << " patch:" << ovp.name() << endl;
783 }
784
785
786 // Calculate stabilised diagonal as normalisation for interpolation
787 const scalarField norm
788 (
789 dynamic_cast<const oversetFvMeshBase&>(mesh).normalisation(matrix)
790 );
791
793 const labelUList& types = overlap.cellTypes();
794 const labelListList& stencil = overlap.cellStencil();
795
796 dynamic_cast<const oversetFvMeshBase&>(mesh).addInterpolation
797 (
798 matrix,
799 norm,
800 this->setHoleCellValue_,
801 this->holeCellValue_
802 );
803
804 if (debug)
805 {
806 pointField allCs(mesh.cellCentres());
808 map.distribute(allCs, false, UPstream::msgType()+1);
809
810 // Make sure we don't interpolate from a hole
811
812 scalarField marker(this->primitiveField().size(), 0);
813 forAll(types, celli)
814 {
815 if (types[celli] == cellCellStencil::HOLE)
816 {
817 marker[celli] = 1.0;
818 }
819 }
821 (
822 marker,
823 mesh,
824 overlap,
826 );
827
828 forAll(marker, celli)
829 {
830 if
831 (
832 types[celli] == cellCellStencil::INTERPOLATED
833 && marker[celli] > SMALL
834 )
835 {
837 << " field:" << fldName
838 << " patch:" << ovp.name()
839 << " found:" << celli
840 << " at:" << mesh.cellCentres()[celli]
841 << " donorSlots:" << stencil[celli]
842 << " at:"
843 << UIndirectList<point>(allCs, stencil[celli])
844 << " amount-of-hole:" << marker[celli]
845 << endl;
846 }
847 }
848
849 // Make sure we don't have matrix coefficients for interpolated
850 // or hole cells
851
852 const lduAddressing& addr = mesh.lduAddr();
853 const labelUList& upperAddr = addr.upperAddr();
854 const labelUList& lowerAddr = addr.lowerAddr();
855 const scalarField& lower = matrix.lower();
856 const scalarField& upper = matrix.upper();
857
858 forAll(lowerAddr, facei)
859 {
860 const label l = lowerAddr[facei];
861 const bool lHole = (types[l] == cellCellStencil::HOLE);
862 const label u = upperAddr[facei];
863 const bool uHole = (types[u] == cellCellStencil::HOLE);
864
865 if
866 (
867 (lHole && upper[facei] != 0.0)
868 || (uHole && lower[facei] != 0.0)
869 )
870 {
872 << "Hole-neighbouring face:" << facei
873 << " lower:" << l
874 << " type:" << types[l]
875 << " coeff:" << lower[facei]
876 << " upper:" << upperAddr[facei]
877 << " type:" << types[u]
878 << " coeff:" << upper[facei]
879 << exit(FatalError);
880 }
881
882
883 // Empty donor list: treat like hole but still allow to
884 // influence neighbouring domains
885 const bool lEmpty =
886 (
888 && stencil[l].empty()
889 );
890 const bool uEmpty =
891 (
893 && stencil[u].empty()
894 );
895
896 if
897 (
898 (lEmpty && upper[facei] != 0.0)
899 || (uEmpty && lower[facei] != 0.0)
900 )
901 {
903 << "Still connected face:" << facei << " lower:" << l
904 << " type:" << types[l]
905 << " coeff:" << lower[facei]
906 << " upper:" << u
907 << " type:" << types[u]
908 << " coeff:" << upper[facei]
909 << exit(FatalError);
910 }
911 }
912
913 forAll(matrix.internalCoeffs(), patchi)
914 {
915 const labelUList& fc = addr.patchAddr(patchi);
916 const Field<Type>& bouCoeffs = matrix.boundaryCoeffs()[patchi];
917
918 forAll(fc, i)
919 {
920 const label celli = fc[i];
921 const bool lHole = (types[celli] == cellCellStencil::HOLE);
922
923 if (lHole && bouCoeffs[i] != pTraits<Type>::zero)
924 {
926 << "Patch:" << patchi
927 << " patchFace:" << i
928 << " lower:" << celli
929 << " type:" << types[celli]
930 << " bouCoeff:" << bouCoeffs[i]
931 << exit(FatalError);
932 }
933
934 // Check whether this is influenced by neighbouring domains
935 const bool lEmpty =
936 (
937 types[celli] == cellCellStencil::INTERPOLATED
938 && stencil[celli].empty()
939 );
940
941 if (lEmpty && bouCoeffs[i] != pTraits<Type>::zero)
942 {
944 << "Patch:" << patchi
945 << " patchFace:" << i
946 << " lower:" << celli
947 << " type:" << types[celli]
948 << " bouCoeff:" << bouCoeffs[i]
949 << exit(FatalError);
950 }
951 }
952 }
953
954
955 // Make sure that diagonal is non-zero. Note: should add
956 // boundaryCoeff ...
957 const FieldField<Field, Type>& internalCoeffs =
958 matrix.internalCoeffs();
959
960 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; ++cmpt)
961 {
962 // Replacement for m.addBoundaryDiag(norm, cmpt);
963 scalarField diag(matrix.diag());
964 forAll(internalCoeffs, patchi)
965 {
966 const labelUList& fc = addr.patchAddr(patchi);
967 const Field<Type>& intCoeffs = internalCoeffs[patchi];
968 const scalarField cmptCoeffs(intCoeffs.component(cmpt));
969 forAll(fc, i)
970 {
971 diag[fc[i]] += cmptCoeffs[i];
972 }
973 }
974
975 forAll(diag, celli)
976 {
977 if (mag(diag[celli]) < SMALL)
978 {
980 << "Patch:" << ovp.name()
981 << " cell:" << celli
982 << " at:" << mesh.cellCentres()[celli]
983 << " diag:" << diag[celli]
984 << exit(FatalError);
985 }
986 }
987 }
988 }
990
992}
993
994
995template<class Type>
997(
998 solveScalarField& result,
999 const bool add,
1000 const lduAddressing& lduAddr,
1001 const label interfacei,
1002 const solveScalarField& psiInternal,
1003 const scalarField& coeffs,
1004 const direction cmpt,
1005 const Pstream::commsTypes commsType
1006) const
1007{
1008}
1009
1010
1011template<class Type>
1013(
1014 solveScalarField& result,
1015 const bool add,
1016 const lduAddressing& lduAddr,
1017 const label patchId,
1018 const solveScalarField& psiInternal,
1019 const scalarField& coeffs,
1020 const direction,
1021 const Pstream::commsTypes commsType
1022) const
1023{
1024 auto& psi = const_cast<solveScalarField&>(psiInternal);
1025
1026 if (fluxCorrection_ && this->oversetPatch_.master())
1028 adjustPsi(psi, lduAddr, result);
1029 }
1030}
1031
1032
1033template<class Type>
1034void Foam::oversetFvPatchField<Type>::write(Ostream& os) const
1035{
1037
1038 if (this->setHoleCellValue_)
1039 {
1040 os.writeEntry("setHoleCellValue", setHoleCellValue_);
1041 os.writeEntry("holeCellValue", holeCellValue_);
1042 os.writeEntryIfDifferent
1043 (
1044 "interpolateHoleCellValue",
1045 false,
1046 interpolateHoleCellValue_
1047 );
1048 }
1050 (
1051 "fluxCorrection",
1052 false,
1053 fluxCorrection_
1054 );
1056 (
1057 "zone",
1058 label(-1),
1059 zoneId_
1060 );
1061}
1062
1063
1064// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
void map(const UList< Type > &mapF, const labelUList &mapAddressing)
1 to 1 map from the given field
Definition Field.C:302
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition Field.C:607
bool found(const Key &key) const
Same as contains().
Definition HashTable.H:1370
A simple container of IOobject preferences. Can also be used for general handling of read/no-read/rea...
static FOAM_NO_DANGLING_REFERENCE const cellCellStencilObject & New(const fvMesh &mesh, Args &&... args)
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
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:346
A List with indirect addressing. Like IndirectList but does not store addressing.
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
commsTypes
Communications types.
Definition UPstream.H:81
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
Addressing for all faces on surface of mesh. Can either be read from polyMesh or from triSurface....
virtual const wordHashSet & nonInterpolatedFields() const
Return the names of any (stencil or mesh specific) fields that.
virtual const mapDistribute & cellInterpolationMap() const
Return a communication schedule.
virtual const labelUList & cellTypes() const
Return the cell type list.
virtual const List< scalarList > & cellInterpolationWeights() const
Weights for cellStencil.
virtual const labelListList & cellStencil() const
Per interpolated cell the neighbour cells (in terms of slots as.
Calculation of interpolation stencils.
static void interpolate(Field< T > &psi, const fvMesh &mesh, const cellCellStencil &overlap, const List< scalarList > &wghts)
Interpolation of acceptor cells from donor cells.
static const labelIOList & zoneID(const fvMesh &)
Helper: get reference to registered zoneID. Loads volScalarField.
Abstract base class for coupled patches.
virtual void write(Ostream &) const
Write includes "value" entry.
virtual void initEvaluate(const Pstream::commsTypes commsType)
Initialise the evaluation of the patch field.
coupledFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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.
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...
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition dictionary.C:441
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...
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition fvMatrix.H:118
const FieldField< Field, Type > & internalCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for internal cells
Definition fvMatrix.H:549
const FieldField< Field, Type > & boundaryCoeffs() const noexcept
fvBoundary scalar field containing pseudo-matrix coeffs for boundary cells
Definition fvMatrix.H:567
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
virtual const lduAddressing & lduAddr() const
Return ldu addressing.
Definition fvMesh.C:707
virtual void interpolate(volScalarField &) const
Interpolate interpolationCells only.
Definition fvMesh.H:459
const labelUList & owner() const
Internal face owner. Note bypassing virtual mechanism so.
Definition fvMesh.H:572
const labelUList & neighbour() const
Internal face neighbour.
Definition fvMesh.H:580
const fvBoundaryMesh & boundary() const noexcept
Return reference to boundary mesh.
Definition fvMesh.H:395
bool manipulatedMatrix() const noexcept
True if the matrix has already been manipulated.
A FieldMapper for finite-volume patch fields.
const Field< Type > & primitiveField() const noexcept
Return const-reference to the internal field values.
const DimensionedField< Type, volMesh > & internalField() const noexcept
Return const-reference to the dimensioned internal field.
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
void extrapolateInternal()
Assign the patch field from the internal field.
A finiteVolume patch using a polyPatch and a fvBoundaryMesh.
Definition fvPatch.H:71
virtual const word & name() const
Return name.
Definition fvPatch.H:210
label start() const noexcept
The patch start within the polyMesh face list.
Definition fvPatch.H:226
virtual const labelUList & faceCells() const
Return faceCells.
Definition fvPatch.C:107
Selector class for finite volume differencing schemes. fvMesh is derived from fvSchemes so that all f...
Definition fvSchemes.H:54
The class contains the addressing required by the lduMatrix: upper, lower and losort.
virtual const labelUList & upperAddr() const =0
Return upper addressing.
virtual const labelUList & lowerAddr() const =0
Return lower addressing.
virtual const labelUList & patchAddr(const label patchNo) const =0
Return patch to internal addressing given patch number.
const scalarField & diag() const
Definition lduMatrix.C:195
const scalarField & upper() const
Definition lduMatrix.C:235
const scalarField & lower() const
Definition lduMatrix.C:306
Class containing processor-to-processor mapping information.
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.
Support for overset functionality.
Boundary condition for use on overset patches. To be run in combination with special dynamicFvMesh ty...
const bool setHoleCellValue_
Flag to set hole cell values.
void storeFringeCoefficients(const fvMatrix< Type > &matrix)
Store fringe coefficients and faces.
const oversetFvPatch & oversetPatch_
Local reference cast into the overset patch.
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.
const bool fluxCorrection_
Flag to correct fluxes.
oversetFvPatchField(const fvPatch &, const DimensionedField< Type, volMesh > &)
Construct from patch and internal field.
virtual void initEvaluate(const Pstream::commsTypes commsType)
Initialise the evaluation of the patch field.
scalarField fringeUpperCoeffs_
Fringe upper coefficients.
virtual void write(Ostream &os) const
Write.
void adjustPsi(solveScalarField &psi, const lduAddressing &lduAddr, solveScalarField &result) const
Adjust psi for mass correction. Requires storeFringeCoefficients.
const Type holeCellValue_
Hole cell value.
void fringeFlux(const fvMatrix< Type > &m, const surfaceScalarField &phi) const
Calculate patch flux (helper function). Requires.
virtual void manipulateMatrix(fvMatrix< Type > &matrix)
Manipulate matrix.
labelField fringeFaces_
Fringe faces.
const bool interpolateHoleCellValue_
Flag to interpolate hole cell values (from nearby non-hole cell).
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 zoneId_
Zone to sum flux for mass conservation.
scalarField fringeLowerCoeffs_
Fringe lower coefficients.
Patch for indicating interpolated boundaries (in overset meshes).
virtual bool master() const
Am I the master interface.
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh,...
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition polyMesh.H:609
A patch is a list of labels that address the faces in the global face list.
Definition polyPatch.H:73
label nInternalFaces() const noexcept
Number of internal faces.
const vectorField & cellCentres() const
const dictionary & schemesDict() const
The entire dictionary or the optional "select" sub-dictionary.
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.
A class for handling words, derived from Foam::string.
Definition word.H:66
volScalarField & p
const volScalarField & psi
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition error.H:629
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
const cellCellStencilObject & overlap
Definition correctPhi.H:57
label patchId(-1)
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
#define WarningInFunction
Report a warning using Foam::Warning.
#define FUNCTION_NAME
Namespace for handling debugging switches.
Definition debug.C:45
string upper(const std::string &s)
Return string copy transformed with std::toupper on each character.
string lower(const std::string &s)
Return string copy transformed with std::tolower on each character.
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition HashSet.H:80
PtrList< fvPatch > fvPatchList
Store lists of fvPatch as a PtrList.
Definition fvPatch.H:64
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
void component(FieldField< Field, typename FieldField< Field, Type >::cmptType > &sf, const FieldField< Field, Type > &f, const direction d)
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void add(DimensionedField< scalar, GeoMesh > &result, const dimensioned< scalar > &dt1, const DimensionedField< scalar, GeoMesh > &f2)
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
void reduce(T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce).
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< solveScalar > solveScalarField
uint8_t direction
Definition direction.H:49
Field< label > labelField
Specialisation of Field<T> for label.
Definition labelField.H:48
IOerror FatalIOError
Error stream (stdout output on all processes), with additional 'FOAM FATAL IO ERROR' header text and ...
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
vectorField pointField
pointField is a vectorField.
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
UList< label > labelUList
A UList of labels.
Definition UList.H:75
dict add("bounds", meshBb)
dictionary dict
const labelUList & cellTypes
Definition setCellMask.H:27
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299