Loading...
Searching...
No Matches
oversetFvMeshBaseTemplates.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) 2014-2025 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 "fvMatrix.H"
30#include "oversetFvPatchField.H"
33#include "processorFvPatch.H"
34#include "syncTools.H"
35
36// * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
37
38template<class Type>
40(
41 Field<Type>& coeffs,
42 const labelUList& types,
43 const scalarList& factor,
44 const bool setHoleCellValue,
45 const label celli,
46 const label facei
47) const
48{
49 const label cType = types[celli];
50 const scalar f = factor[celli];
51
53 {
54 coeffs[facei] *= 1.0-f;
55 }
56 else if (cType == cellCellStencil::HOLE)
57 {
58 // Disconnect hole cell from influence of neighbour
59 coeffs[facei] = pTraits<Type>::zero;
60 }
61 else if (cType == cellCellStencil::SPECIAL)
62 {
63 if (setHoleCellValue)
64 {
65 // Behave like hole
66 coeffs[facei] = pTraits<Type>::zero;
67 }
68 else
69 {
70 // Behave like interpolated
71 coeffs[facei] *= 1.0-f;
72 }
73 }
74}
75
76template<class GeoField, class PatchType, bool TypeOnly>
78(
79 typename GeoField::Boundary& bfld
80)
81{
82 if constexpr (TypeOnly)
83 {
84 bfld.evaluate_if
85 (
86 [](const auto& pfld)
87 {
88 return (bool(isA<PatchType>(pfld)));
89 }
90 );
91 }
92 else
93 {
94 bfld.evaluate_if
95 (
96 [](const auto& pfld)
97 {
98 return (!bool(isA<PatchType>(pfld)));
99 }
100 );
101 }
102}
103
104
105template<class Type>
106Foam::tmp<Foam::scalarField> Foam::oversetFvMeshBase::normalisation
107(
108 const fvMatrix<Type>& m
109) const
110{
111 // Determine normalisation. This is normally the original diagonal plus
112 // remote contributions. This needs to be stabilised for hole cells
113 // which can have a zero diagonal. Assume that if any component has
114 // a non-zero diagonal the cell does not need stabilisation.
116 scalarField& norm = tnorm.ref();
117
118 // Add remote coeffs to duplicate behaviour of fvMatrix::addBoundaryDiag
119 const FieldField<Field, Type>& internalCoeffs = m.internalCoeffs();
120 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
121 {
122 forAll(internalCoeffs, patchi)
123 {
124 const labelUList& fc = mesh_.lduAddr().patchAddr(patchi);
125 const Field<Type>& intCoeffs = internalCoeffs[patchi];
126 const scalarField cmptCoeffs(intCoeffs.component(cmpt));
127 forAll(fc, i)
128 {
129 norm[fc[i]] += cmptCoeffs[i];
130 }
131 }
132 }
133
134 // Count number of problematic cells
135 label nZeroDiag = 0;
136 forAll(norm, celli)
137 {
138 const scalar& n = norm[celli];
139 if (magSqr(n) < sqr(SMALL))
140 {
141 //Pout<< "For field " << m.psi().name()
142 // << " have diagonal " << n << " for cell " << celli
143 // << " at:" << cellCentres()[celli] << endl;
144 nZeroDiag++;
145 }
146 }
147
148 if (debug)
149 {
150 Pout<< "For field " << m.psi().name() << " have zero diagonals for "
151 << returnReduce(nZeroDiag, sumOp<label>()) << " cells" << endl;
152 }
153
154 if (returnReduceOr(nZeroDiag))
155 {
156 // Walk out the norm across hole cells
157
158 const labelList& own = mesh_.faceOwner();
159 const labelList& nei = mesh_.faceNeighbour();
161 const labelUList& types = overlap.cellTypes();
162
163 // label nHoles = 0;
164 scalarField extrapolatedNorm(norm);
165 forAll(types, celli)
166 {
167 if (types[celli] == cellCellStencil::HOLE)
168 {
169 extrapolatedNorm[celli] = -GREAT;
170 // ++nHoles;
171 }
172 }
173
174 bitSet isFront(mesh_.nFaces());
175 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
176 {
177 label ownType = types[own[facei]];
178 label neiType = types[nei[facei]];
179 if
180 (
181 (ownType == cellCellStencil::HOLE)
182 != (neiType == cellCellStencil::HOLE)
183 )
184 {
185 isFront.set(facei);
186 }
187 }
188 labelList nbrTypes;
189 syncTools::swapBoundaryCellList(mesh_, types, nbrTypes);
190 for
191 (
192 label facei = mesh_.nInternalFaces();
193 facei < mesh_.nFaces();
194 ++facei
195 )
196 {
197 const label ownType = types[own[facei]];
198 const label neiType = nbrTypes[facei-mesh_.nInternalFaces()];
199 if
200 (
201 (ownType == cellCellStencil::HOLE)
202 != (neiType == cellCellStencil::HOLE)
203 )
204 {
205 isFront.set(facei);
206 }
207 }
208
209
210 while (true)
211 {
212 scalarField nbrNorm;
213 syncTools::swapBoundaryCellList(mesh_, extrapolatedNorm, nbrNorm);
214
215 bitSet newIsFront(mesh_.nFaces());
216 scalarField newNorm(extrapolatedNorm);
217
218 label nChanged = 0;
219 for (const label facei : isFront)
220 {
221 if (extrapolatedNorm[own[facei]] == -GREAT)
222 {
223 // Average owner cell, add faces to newFront
224 newNorm[own[facei]] = cellAverage
225 (
226 types,
227 nbrTypes,
228 extrapolatedNorm,
229 nbrNorm,
230 own[facei],
231 newIsFront
232 );
233 nChanged++;
234 }
235 if
236 (
237 mesh_.isInternalFace(facei)
238 && extrapolatedNorm[nei[facei]] == -GREAT
239 )
240 {
241 // Average nei cell, add faces to newFront
242 newNorm[nei[facei]] = cellAverage
243 (
244 types,
245 nbrTypes,
246 extrapolatedNorm,
247 nbrNorm,
248 nei[facei],
249 newIsFront
250 );
251 nChanged++;
252 }
253 }
254
255 if (!returnReduceOr(nChanged))
256 {
257 break;
258 }
259
260 // Transfer new front
261 extrapolatedNorm.transfer(newNorm);
262 isFront.transfer(newIsFront);
264 }
265
266
267 forAll(norm, celli)
268 {
269 scalar& n = norm[celli];
270 if (magSqr(n) < sqr(SMALL))
271 {
272 //Pout<< "For field " << m.psi().name()
273 // << " for cell " << celli
274 // << " at:" << cellCentres()[celli]
275 // << " have norm " << n
276 // << " have extrapolated norm " << extrapolatedNorm[celli]
277 // << endl;
278 // Override the norm
279 n = extrapolatedNorm[celli];
280 }
282 }
283 return tnorm;
284}
285
286
287template<class Type>
289(
291 const scalarField& normalisation,
292 const bool setHoleCellValue,
293 const Type& holeCellValue
294) const
295{
298 const labelListList& stencil = overlap.cellStencil();
300 const labelUList& types = overlap.cellTypes();
301
302
303 // Force asymmetric matrix (if it wasn't already)
304 scalarField& lower = m.lower();
305 scalarField& upper = m.upper();
306 Field<Type>& source = m.source();
307 scalarField& diag = m.diag();
308
309
310 // Get the addressing. Note that the addressing is now extended with
311 // any interpolation faces.
312 const lduAddressing& addr = lduAddr();
313 const labelUList& upperAddr = addr.upperAddr();
314 const labelUList& lowerAddr = addr.lowerAddr();
315 const lduInterfacePtrsList& interfaces = allInterfaces_;
316
318 {
320 << "Problem : addressing is not fvMeshPrimitiveLduAddressing"
321 << exit(FatalError);
322 }
323
324
325 // 1. Adapt lduMatrix for additional faces and new ordering
326 upper.setSize(upperAddr.size(), 0.0);
327 inplaceReorder(reverseFaceMap_, upper);
328 lower.setSize(lowerAddr.size(), 0.0);
329 inplaceReorder(reverseFaceMap_, lower);
330
331
332 //const label nOldInterfaces = dynamicMotionSolverFvMesh::interfaces().size();
333 const label nOldInterfaces = mesh_.fvMesh::interfaces().size();
334
335
336 if (interfaces.size() > nOldInterfaces)
337 {
338 // Extend matrix coefficients
339 m.internalCoeffs().setSize(interfaces.size());
340 m.boundaryCoeffs().setSize(interfaces.size());
341
342 // 1b. Adapt for additional interfaces
343 for
344 (
345 label patchi = nOldInterfaces;
346 patchi < interfaces.size();
347 patchi++
348 )
349 {
350 const labelUList& fc = interfaces[patchi].faceCells();
351 m.internalCoeffs().set(patchi, new Field<Type>(fc.size(), Zero));
352 m.boundaryCoeffs().set(patchi, new Field<Type>(fc.size(), Zero));
353 }
354
355 // 1c. Adapt field for additional interfaceFields (note: solver uses
356 // GeometricField::scalarInterfaces() to get hold of interfaces)
357
358 auto& bfld = m.psi().constCast().boundaryFieldRef();
359
360 bfld.setSize(interfaces.size());
361
362
363 // This gets quite interesting: we do not want to add additional
364 // fvPatches (since direct correspondence to polyMesh) so instead
365 // add a reference to an existing processor patch
366 label addPatchi = 0;
367 for (label patchi = 0; patchi < nOldInterfaces; patchi++)
368 {
369 if (isA<processorFvPatch>(bfld[patchi].patch()))
370 {
371 addPatchi = patchi;
372 break;
373 }
374 }
375
376 for
377 (
378 label patchi = nOldInterfaces;
379 patchi < interfaces.size();
380 patchi++
381 )
382 {
383 bfld.set
384 (
385 patchi,
387 (
388 interfaces[patchi],
389 bfld[addPatchi].patch(), // dummy processorFvPatch
390 m.psi()
391 )
392 );
393 }
394 }
395
396
397 // 2. Adapt fvMatrix level: faceFluxCorrectionPtr
398 // Question: do we need to do this?
399 // This seems to be set/used only by the gaussLaplacianScheme and
400 // fvMatrix:correction, both of which are outside the linear solver.
401
402
403 // Clear out existing connections on cells to be interpolated
404 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
405 // Note: could avoid doing the zeroing of the new faces since these
406 // are set to zero anyway.
407
408 forAll(upperAddr, facei)
409 {
410 const label l = lowerAddr[facei];
411 scaleConnection(upper, types, factor, setHoleCellValue, l, facei);
412 const label u = upperAddr[facei];
413 scaleConnection(lower, types, factor, setHoleCellValue, u, facei);
414 /*
415 if (types[upperAddr[facei]] == cellCellStencil::INTERPOLATED)
416 {
417 // Disconnect upper from lower
418 label celli = upperAddr[facei];
419 lower[facei] *= 1.0-factor[celli];
420 }
421 if (types[lowerAddr[facei]] == cellCellStencil::INTERPOLATED)
422 {
423 // Disconnect lower from upper
424 label celli = lowerAddr[facei];
425 upper[facei] *= 1.0-factor[celli];
426 }
427 */
428 }
429
430 for (label patchi = 0; patchi < nOldInterfaces; ++patchi)
431 {
432 const labelUList& fc = addr.patchAddr(patchi);
433 Field<Type>& bCoeffs = m.boundaryCoeffs()[patchi];
434 Field<Type>& iCoeffs = m.internalCoeffs()[patchi];
435 forAll(fc, i)
436 {
437 scaleConnection(bCoeffs, types, factor, setHoleCellValue, fc[i], i);
438
439 scaleConnection(iCoeffs, types, factor, setHoleCellValue, fc[i], i);
440 }
441 /*
442 const labelUList& fc = addr.patchAddr(patchi);
443 Field<Type>& intCoeffs = m.internalCoeffs()[patchi];
444 Field<Type>& bouCoeffs = m.boundaryCoeffs()[patchi];
445 forAll(fc, i)
446 {
447 label celli = fc[i];
448 {
449 if (types[celli] == cellCellStencil::INTERPOLATED)
450 {
451 scalar f = factor[celli];
452 intCoeffs[i] *= 1.0-f;
453 bouCoeffs[i] *= 1.0-f;
454 }
455 else if (types[celli] == cellCellStencil::HOLE)
456 {
457 intCoeffs[i] = pTraits<Type>::zero;
458 bouCoeffs[i] = pTraits<Type>::zero;
459 }
460 }
461 }
462 */
463 }
464
465
466
467 // Modify matrix
468 // ~~~~~~~~~~~~~
469
470 // Do hole cells. Note: maybe put into interpolationCells() loop above?
471 forAll(types, celli)
472 {
473 if
474 (
475 types[celli] == cellCellStencil::HOLE
476 || (setHoleCellValue && types[celli] == cellCellStencil::SPECIAL)
477 )
478 {
479 const Type wantedValue
480 (
481 setHoleCellValue
482 ? holeCellValue
483 : m.psi()[celli]
484 );
485 diag[celli] = normalisation[celli];
486 source[celli] = normalisation[celli]*wantedValue;
487 }
488 else if
489 (
490 types[celli] == cellCellStencil::INTERPOLATED
491 || (!setHoleCellValue && types[celli] == cellCellStencil::SPECIAL)
492 )
493 {
494 const scalar f = factor[celli];
495 const scalarList& w = wghts[celli];
496 const labelList& nbrs = stencil[celli];
497 const labelList& nbrFaces = stencilFaces_[celli];
498 const labelList& nbrPatches = stencilPatches_[celli];
499
500 diag[celli] *= (1.0-f);
501 source[celli] *= (1.0-f);
502
503 forAll(nbrs, nbri)
504 {
505 const label patchi = nbrPatches[nbri];
506 const label facei = nbrFaces[nbri];
507
508 if (patchi == -1)
509 {
510 const label nbrCelli = nbrs[nbri];
511 // Add the coefficients
512 const scalar s = normalisation[celli]*f*w[nbri];
513
514 scalar& u = upper[facei];
515 scalar& l = lower[facei];
516 if (celli < nbrCelli)
517 {
518 diag[celli] += s;
519 u += -s;
520 }
521 else
522 {
523 diag[celli] += s;
524 l += -s;
525 }
526 }
527 else
528 {
529 // Patch face. Store in boundaryCoeffs. Note sign change.
530 //const label globalCelli = globalCellIDs[nbrs[nbri]];
531 //const label proci =
532 // globalNumbering.whichProcID(globalCelli);
533 //const label remoteCelli =
534 // globalNumbering.toLocal(proci, globalCelli);
535 //
536 //Pout<< "for cell:" << celli
537 // << " need weight from remote slot:" << nbrs[nbri]
538 // << " proc:" << proci << " remote cell:" << remoteCelli
539 // << " patch:" << patchi
540 // << " patchFace:" << facei
541 // << " weight:" << w[nbri]
542 // << endl;
543
544 const scalar s = normalisation[celli]*f*w[nbri];
545 m.boundaryCoeffs()[patchi][facei] += pTraits<Type>::one*s;
546 m.internalCoeffs()[patchi][facei] += pTraits<Type>::one*s;
547
548 // Note: do NOT add to diagonal - this is in the
549 // internalCoeffs and gets added to the diagonal
550 // inside fvMatrix::solve
551 }
552 }
553 }
554 /*
555 label startLabel = ownerStartAddr[celli];
556 label endLabel = ownerStartAddr[celli + 1];
557
558 for (label facei = startLabel; facei < endLabel; facei++)
559 {
560 upper[facei] = 0.0;
561 }
562
563 startLabel = addr.losortStartAddr()[celli];
564 endLabel = addr.losortStartAddr()[celli + 1];
565
566 for (label i = startLabel; i < endLabel; i++)
567 {
568 label facei = losortAddr[i];
569 lower[facei] = 0.0;
570 }
571
572 diag[celli] = normalisation[celli];
573 source[celli] = normalisation[celli]*m.psi()[celli];
574 */
575 }
576
577
578 //const globalIndex globalNumbering(V().size());
579 //labelList globalCellIDs(overlap.cellInterpolationMap().constructSize());
580 //forAll(V(), cellI)
581 //{
582 // globalCellIDs[cellI] = globalNumbering.toGlobal(cellI);
583 //}
584 //overlap.cellInterpolationMap().distribute(globalCellIDs);
585
586/*
587 forAll(cellIDs, i)
588 {
589 label celli = cellIDs[i];
590
591 const scalar f = factor[celli];
592 const scalarList& w = wghts[celli];
593 const labelList& nbrs = stencil[celli];
594 const labelList& nbrFaces = stencilFaces_[celli];
595 const labelList& nbrPatches = stencilPatches_[celli];
596
597 if (types[celli] == cellCellStencil::HOLE)
598 {
599 FatalErrorInFunction << "Found HOLE cell " << celli
600 << " at:" << C()[celli]
601 << " . Should this be in interpolationCells()????"
602 << abort(FatalError);
603 }
604 else
605 {
606 // Create interpolation stencil
607
608 diag[celli] *= (1.0-f);
609 source[celli] *= (1.0-f);
610
611 forAll(nbrs, nbri)
612 {
613 label patchi = nbrPatches[nbri];
614 label facei = nbrFaces[nbri];
615
616 if (patchi == -1)
617 {
618 label nbrCelli = nbrs[nbri];
619
620 // Add the coefficients
621 const scalar s = normalisation[celli]*f*w[nbri];
622
623 scalar& u = upper[facei];
624 scalar& l = lower[facei];
625 if (celli < nbrCelli)
626 {
627 diag[celli] += s;
628 u += -s;
629 }
630 else
631 {
632 diag[celli] += s;
633 l += -s;
634 }
635 }
636 else
637 {
638 // Patch face. Store in boundaryCoeffs. Note sign change.
639 //const label globalCelli = globalCellIDs[nbrs[nbri]];
640 //const label proci =
641 // globalNumbering.whichProcID(globalCelli);
642 //const label remoteCelli =
643 // globalNumbering.toLocal(proci, globalCelli);
644 //
645 //Pout<< "for cell:" << celli
646 // << " need weight from remote slot:" << nbrs[nbri]
647 // << " proc:" << proci << " remote cell:" << remoteCelli
648 // << " patch:" << patchi
649 // << " patchFace:" << facei
650 // << " weight:" << w[nbri]
651 // << endl;
652
653 const scalar s = normalisation[celli]*f*w[nbri];
654 m.boundaryCoeffs()[patchi][facei] += pTraits<Type>::one*s;
655 m.internalCoeffs()[patchi][facei] += pTraits<Type>::one*s;
656
657 // Note: do NOT add to diagonal - this is in the
658 // internalCoeffs and gets added to the diagonal
659 // inside fvMatrix::solve
660 }
661 }
662
663 //if (mag(diag[celli]) < SMALL)
664 //{
665 // Pout<< "for cell:" << celli
666 // << " at:" << this->C()[celli]
667 // << " diag:" << diag[celli] << endl;
668 //
669 // forAll(nbrs, nbri)
670 // {
671 // label patchi = nbrPatches[nbri];
672 // label facei = nbrFaces[nbri];
673 //
674 // const label globalCelli = globalCellIDs[nbrs[nbri]];
675 // const label proci =
676 // globalNumbering.whichProcID(globalCelli);
677 // const label remoteCelli =
678 // globalNumbering.toLocal(proci, globalCelli);
679 //
680 // Pout<< " need weight from slot:" << nbrs[nbri]
681 // << " proc:" << proci << " remote cell:"
682 // << remoteCelli
683 // << " patch:" << patchi
684 // << " patchFace:" << facei
685 // << " weight:" << w[nbri]
686 // << endl;
687 // }
688 // Pout<< endl;
689 //}
691 }
692*/
693}
694
695
696template<class Type>
698(
700 const dictionary& dict
701) const
702{
704 // Check we're running with bcs that can handle implicit matrix manipulation
705 auto& bpsi = m.psi().constCast().boundaryFieldRef();
706
707 bool hasOverset = false;
708 forAll(bpsi, patchi)
709 {
710 if (isA<oversetFvPatchField<Type>>(bpsi[patchi]))
711 {
712 hasOverset = true;
713 break;
714 }
715 }
716
717 if (!hasOverset)
718 {
719 if (debug)
720 {
721 Pout<< "oversetFvMeshBase::solveOverset() :"
722 << " bypassing matrix adjustment for field " << m.psi().name()
723 << endl;
724 }
725 //return dynamicMotionSolverFvMesh::solve(m, dict);
726 return mesh_.fvMesh::solve(m, dict);
727 }
728
729 if (debug)
730 {
731 Pout<< "oversetFvMeshBase::solveOverset() :"
732 << " adjusting matrix for interpolation for field "
733 << m.psi().name() << endl;
734 }
735
736 // Calculate stabilised diagonal as normalisation for interpolation
737 const scalarField norm(normalisation(m));
738
739 if (debug && mesh_.time().writeTime())
740 {
741 volScalarField scale
742 (
744 (
745 m.psi().name() + "_normalisation",
746 mesh_.time().timeName(),
747 mesh_,
751 ),
752 mesh_,
754 );
755 scale.ref().field() = norm;
757 <
760 false
761 >(scale.boundaryFieldRef());
762 scale.write();
763
764 if (debug)
765 {
766 Pout<< "oversetFvMeshBase::solveOverset() :"
767 << " writing matrix normalisation for field " << m.psi().name()
768 << " to " << scale.name() << endl;
769 }
770 }
771
772
773 // Switch to extended addressing (requires mesh::update() having been
774 // called)
775 active(true);
776
777 // Adapt matrix
778 scalarField oldUpper(m.upper());
779 scalarField oldLower(m.lower());
782
783 Field<Type> oldSource(m.source());
784 scalarField oldDiag(m.diag());
785
786 const label oldSize = bpsi.size();
787
788 // Insert the interpolation into the matrix (done inside
789 // oversetFvPatchField<Type>::manipulateMatrix)
790 m.boundaryManipulate(bpsi);
791
792 // Swap psi values so added patches have patchNeighbourField
794 <
795 GeoField,
797 true
798 >(bpsi);
799
800 // Use lower level solver
801 //SolverPerformance<Type> s(dynamicMotionSolverFvMesh::solve(m, dict));
802 SolverPerformance<Type> s(mesh_.fvMesh::solve(m, dict));
803
804 // Restore boundary field
805 bpsi.setSize(oldSize);
806
807 // Restore matrix
808 m.upper().transfer(oldUpper);
809 m.lower().transfer(oldLower);
810
811 m.source().transfer(oldSource);
812 m.diag().transfer(oldDiag);
813
814 m.internalCoeffs().transfer(oldInt);
815 m.boundaryCoeffs().transfer(oldBou);
816
817
819 const labelUList& types = overlap.cellTypes();
820 const labelList& own = mesh_.faceOwner();
821 const labelList& nei = mesh_.faceNeighbour();
822
823 auto& psi = m.psi().constCast();
824
825 // Mirror hole cell values next to calculated cells
826 for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
827 {
828 const label ownType = types[own[facei]];
829 const label neiType = types[nei[facei]];
830
831 if
832 (
833 ownType == cellCellStencil::HOLE
834 && neiType == cellCellStencil::CALCULATED)
835 {
836 psi[own[facei]] = psi[nei[facei]];
837 }
838 else if
839 (
841 && neiType == cellCellStencil::HOLE
842 )
843 {
844 psi[nei[facei]] = psi[own[facei]];
845 }
846 }
847
848 // Switch to original addressing
849 active(false);
850
851 return s;
852}
853
854
855template<class Type>
857(
858 Ostream& os,
859 const fvMatrix<Type>& m,
860 const lduAddressing& addr,
861 const lduInterfacePtrsList& interfaces
862) const
863{
864 os << "** Matrix **" << endl;
865 const labelUList& upperAddr = addr.upperAddr();
866 const labelUList& lowerAddr = addr.lowerAddr();
867 const labelUList& ownerStartAddr = addr.ownerStartAddr();
868 const labelUList& losortAddr = addr.losortAddr();
869
870 const scalarField& lower = m.lower();
871 const scalarField& upper = m.upper();
872 const Field<Type>& source = m.source();
873 const scalarField& diag = m.diag();
874
875
876 // Invert patch addressing
877 labelListList cellToPatch(addr.size());
878 labelListList cellToPatchFace(addr.size());
879 {
880 forAll(interfaces, patchi)
881 {
882 if (interfaces.set(patchi))
883 {
884 const labelUList& fc = interfaces[patchi].faceCells();
885
886 forAll(fc, i)
887 {
888 cellToPatch[fc[i]].append(patchi);
889 cellToPatchFace[fc[i]].append(i);
890 }
891 }
892 }
893 }
894
895 forAll(source, celli)
896 {
897 os << "cell:" << celli << " diag:" << diag[celli]
898 << " source:" << source[celli] << endl;
899
900 label startLabel = ownerStartAddr[celli];
901 label endLabel = ownerStartAddr[celli + 1];
902
903 for (label facei = startLabel; facei < endLabel; facei++)
904 {
905 os << " face:" << facei
906 << " nbr:" << upperAddr[facei] << " coeff:" << upper[facei]
907 << endl;
908 }
909
910 startLabel = addr.losortStartAddr()[celli];
911 endLabel = addr.losortStartAddr()[celli + 1];
912
913 for (label i = startLabel; i < endLabel; i++)
914 {
915 label facei = losortAddr[i];
916 os << " face:" << facei
917 << " nbr:" << lowerAddr[facei] << " coeff:" << lower[facei]
918 << endl;
919 }
920
921 forAll(cellToPatch[celli], i)
922 {
923 label patchi = cellToPatch[celli][i];
924 label patchFacei = cellToPatchFace[celli][i];
925
926 os << " patch:" << patchi
927 << " patchface:" << patchFacei
928 << " intcoeff:" << m.internalCoeffs()[patchi][patchFacei]
929 << " boucoeff:" << m.boundaryCoeffs()[patchi][patchFacei]
930 << endl;
931 }
932 }
933 forAll(m.internalCoeffs(), patchi)
934 {
935 if (m.internalCoeffs().set(patchi))
936 {
937 const labelUList& fc = addr.patchAddr(patchi);
938
939 os << "patch:" << patchi
940 //<< " type:" << interfaces[patchi].type()
941 << " size:" << fc.size() << endl;
942 if
943 (
944 interfaces.set(patchi)
945 && isA<processorLduInterface>(interfaces[patchi])
946 )
947 {
948 const processorLduInterface& ppp =
949 refCast<const processorLduInterface>(interfaces[patchi]);
950 os << "(processor with my:" << ppp.myProcNo()
951 << " nbr:" << ppp.neighbProcNo()
952 << ")" << endl;
953 }
954
955 forAll(fc, i)
956 {
957 os << " cell:" << fc[i]
958 << " int:" << m.internalCoeffs()[patchi][i]
959 << " boun:" << m.boundaryCoeffs()[patchi][i]
960 << endl;
961 }
962 }
963 }
964
965
966 lduInterfaceFieldPtrsList interfaceFields =
968 forAll(interfaceFields, inti)
969 {
970 if (interfaceFields.set(inti))
971 {
972 os << "interface:" << inti
973 << " if type:" << interfaceFields[inti].interface().type()
974 << " fld type:" << interfaceFields[inti].type() << endl;
975 }
977
978 os << "** End of Matrix **" << endl;
979}
980
981
982template<class GeoField>
985 // Evaluate all coupled fields
986 fld.boundaryFieldRef().template evaluateCoupled<void>();
987}
988
989
990template<class GeoField>
992{
993 Pout<< "** starting checking of " << fld.name() << endl;
994
995 GeoField fldCorr(fld.name()+"_correct", fld);
996 correctCoupledBoundaryConditions(fldCorr);
997
998 const auto& bfld = fld.boundaryField();
999 const auto& bfldCorr = fldCorr.boundaryField();
1000
1001 forAll(bfld, patchi)
1002 {
1003 const auto& pfld = bfld[patchi];
1004 const auto& pfldCorr = bfldCorr[patchi];
1005
1006 Pout<< "Patch:" << pfld.patch().name() << endl;
1007
1008 forAll(pfld, i)
1009 {
1010 if (pfld[i] != pfldCorr[i])
1011 {
1012 Pout<< " " << i << " orig:" << pfld[i]
1013 << " corrected:" << pfldCorr[i] << endl;
1014 }
1015 }
1016 }
1017 Pout<< "** end of " << fld.name() << endl;
1018}
1019
1020
1021// ************************************************************************* //
label n
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))
const DynamicField< Type > & field() const noexcept
Return const-reference to the primitive field values.
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
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition Field.C:607
lduInterfaceFieldPtrsList scalarInterfaces() const
Return a list of pointers for each patch field with only those pointing to interfaces being set.
Generic GeometricField class.
this_type & constCast() const noexcept
Return non-const reference to this field.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
@ NO_REGISTER
Do not request registration (bool: false).
@ NO_READ
Nothing to be read.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
void append(const T &val)
Append an element at the end of the list.
Definition List.H:497
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
void setSize(const label n)
Same as resize().
Definition PtrList.H:357
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition UPtrList.H:366
void append(T *ptr)
Append an element to the end of the list.
Definition UPtrList.H:877
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition bitSet.H:61
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition bitSetI.H:502
A processorFvPatchField type bypassing fvPatch.
virtual const scalarList & cellInterpolationWeight() const
Per interpolated cell the interpolation factor. (0 = use.
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.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
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
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition fvMatrix.C:1260
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition fvMatrix.H:487
Field< Type > & source() noexcept
Definition fvMatrix.H:535
The class contains the addressing required by the lduMatrix: upper, lower and losort.
const labelUList & ownerStartAddr() const
Return owner start addressing.
const labelUList & losortStartAddr() const
Return losort start addressing.
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.
label size() const noexcept
Return number of equations.
const labelUList & losortAddr() const
Return losort addressing.
const scalarField & diag() const
Definition lduMatrix.C:195
const scalarField & upper() const
Definition lduMatrix.C:235
const scalarField & lower() const
Definition lduMatrix.C:306
void scaleConnection(Field< Type > &coeffs, const labelUList &types, const scalarList &factor, const bool setHoleCellValue, const label celli, const label facei) const
Freeze values at holes.
static void correctBoundaryConditions(typename GeoField::Boundary &bfld)
Correct boundary conditions of certain type (TypeOnly = true) or explicitly not of the type (TypeOnly...
virtual lduInterfacePtrsList interfaces() const
Return a list of pointers for each patch.
const fvMesh & mesh_
Reference to mesh.
SolverPerformance< Type > solveOverset(fvMatrix< Type > &m, const dictionary &) const
Solve given dictionary with settings.
virtual const lduAddressing & lduAddr() const
Return ldu addressing. If active: is (extended).
static void checkCoupledBC(const GeoField &fld)
Debug: check halo swap is ok.
labelListList stencilFaces_
Corresponding faces (in above lduPtr) to the stencil.
labelListList stencilPatches_
Corresponding patches (in above lduPtr) to the stencil.
void addInterpolation(fvMatrix< Type > &m, const scalarField &normalisation, const bool setHoleCellValue, const Type &holeCellValue) const
Manipulate the matrix to add the interpolation/set hole.
void write(Ostream &, const fvMatrix< Type > &, const lduAddressing &, const lduInterfacePtrsList &) const
Debug: print matrix.
lduInterfacePtrsList allInterfaces_
Interfaces for above mesh. Contains both original and above added processorLduInterfaces.
tmp< scalarField > normalisation(const fvMatrix< Type > &m) const
Determine normalisation for interpolation. This equals the original diagonal except stabilised for ze...
static void correctCoupledBoundaryConditions(GeoField &fld)
Debug: correct coupled bc.
scalar cellAverage(const labelUList &types, const labelUList &nbrTypes, const scalarField &norm, const scalarField &nbrNorm, const label celli, bitSet &isFront) const
Average norm of valid neighbours.
bool active() const
Return true if using extended addressing.
labelList reverseFaceMap_
From old to new face labels.
Boundary condition for use on overset patches. To be run in combination with special dynamicFvMesh ty...
virtual void write(Ostream &os) const
Write.
A traits class, which is primarily used for primitives and vector-space.
Definition pTraits.H:64
An abstract base class for processor coupled interfaces.
virtual int neighbProcNo() const =0
Return neighbour processor number (rank in communicator).
virtual int myProcNo() const =0
Return processor number (rank in communicator).
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
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 void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop, const bool parRun=UPstream::parRun())
Synchronize values on all mesh faces.
Definition syncTools.H:465
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
const volScalarField & psi
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
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;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
const cellCellStencilObject & overlap
Definition correctPhi.H:57
Namespace for handling debugging switches.
Definition debug.C:45
const std::string patch
OpenFOAM patch number as a std::string.
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
bool returnReduceOr(const bool value, const int communicator=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
const dimensionSet dimless
Dimensionless.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
UPtrList< const lduInterface > lduInterfacePtrsList
Store lists of lduInterface as a UPtrList.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
uint8_t direction
Definition direction.H:49
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.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
UPtrList< const lduInterfaceField > lduInterfaceFieldPtrsList
List of coupled interface fields to be used in coupling.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
labelList f(nPoints)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299