dynamicOversetFvMeshTemplates.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-2019 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"
31#include "oversetFvPatchField.H"
34#include "processorFvPatch.H"
35#include "syncTools.H"
36
37// * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * * //
38
39template<class T>
41{
42 const cellCellStencil& overlap = Stencil::New(*this);
43 const labelListList& stencil = overlap.cellStencil();
44
45 if (stencil.size() != nCells())
46 {
47 return;
48 }
49
54
55 Field<T> work(psi);
56 map.mapDistributeBase::distribute(work, UPstream::msgType()+1);
57
58 forAll(cellIDs, i)
59 {
60 label celli = cellIDs[i];
61
62 const scalarList& w = wghts[celli];
63 const labelList& nbrs = stencil[celli];
64 const scalar f = factor[celli];
65
67 forAll(nbrs, nbrI)
68 {
69 s += w[nbrI]*work[nbrs[nbrI]];
70 }
71 //Pout<< "Interpolated value:" << s << endl;
72 //T oldPsi = psi[celli];
73 psi[celli] = (1.0-f)*psi[celli] + f*s;
74 //Pout<< "psi was:" << oldPsi << " now:" << psi[celli] << endl;
75 }
76}
77
78
79template<class GeoField>
81{
84}
85
86
87template<class GeoField>
89{
90 auto flds(this->lookupClass<GeoField>());
91 for (auto fldPtr : flds)
92 {
93 const word& name = fldPtr->name();
94 if (!suppressed.found(baseName(name)))
95 {
96 if (debug)
97 {
98 Pout<< "dynamicOversetFvMesh::interpolate: interpolating : "
99 << name << endl;
100 }
101 interpolate(fldPtr->primitiveFieldRef());
102 }
103 else
104 {
105 if (debug)
106 {
107 Pout<< "dynamicOversetFvMesh::interpolate: skipping : " << name
108 << endl;
109 }
110 }
111 }
112}
113
114
115template<class GeoField, class PatchType>
117(
118 typename GeoField::Boundary& bfld,
119 const bool typeOnly
120)
121{
122 const label startOfRequests = UPstream::nRequests();
123
124 forAll(bfld, patchi)
125 {
126 if (typeOnly == (isA<PatchType>(bfld[patchi]) != nullptr))
127 {
128 bfld[patchi].initEvaluate(UPstream::defaultCommsType);
129 }
130 }
131
132 // Wait for outstanding requests
133 if
134 (
137 )
138 {
139 UPstream::waitRequests(startOfRequests);
140 }
141
142 forAll(bfld, patchi)
143 {
144 if (typeOnly == (isA<PatchType>(bfld[patchi]) != nullptr))
145 {
146 bfld[patchi].evaluate(UPstream::defaultCommsType);
147 }
148 }
149}
150
151
152template<class Type>
154(
155 const fvMatrix<Type>& m
156) const
157{
158 // Determine normalisation. This is normally the original diagonal plus
159 // remote contributions. This needs to be stabilised for hole cells
160 // which can have a zero diagonal. Assume that if any component has
161 // a non-zero diagonal the cell does not need stabilisation.
163 scalarField& norm = tnorm.ref();
164
165 // Add remote coeffs to duplicate behaviour of fvMatrix::addBoundaryDiag
166 const FieldField<Field, Type>& internalCoeffs = m.internalCoeffs();
167 for (direction cmpt=0; cmpt<pTraits<Type>::nComponents; cmpt++)
168 {
169 forAll(internalCoeffs, patchi)
170 {
171 const labelUList& fc = lduAddr().patchAddr(patchi);
172 const Field<Type>& intCoeffs = internalCoeffs[patchi];
173 const scalarField cmptCoeffs(intCoeffs.component(cmpt));
174 forAll(fc, i)
175 {
176 norm[fc[i]] += cmptCoeffs[i];
177 }
178 }
179 }
180
181 // Count number of problematic cells
182 label nZeroDiag = 0;
183 forAll(norm, celli)
184 {
185 const scalar& n = norm[celli];
186 if (magSqr(n) < sqr(SMALL))
187 {
188 //Pout<< "For field " << m.psi().name()
189 // << " have diagonal " << n << " for cell " << celli
190 // << " at:" << cellCentres()[celli] << endl;
191 nZeroDiag++;
192 }
193 }
194
195 reduce(nZeroDiag, sumOp<label>());
196
197 if (debug)
198 {
199 Pout<< "For field " << m.psi().name() << " have zero diagonals for "
200 << nZeroDiag << " cells" << endl;
201 }
202
203 if (nZeroDiag > 0)
204 {
205 // Walk out the norm across hole cells
206
207 const labelList& own = faceOwner();
208 const labelList& nei = faceNeighbour();
210 const labelUList& types = overlap.cellTypes();
211
212 label nHoles = 0;
213 scalarField extrapolatedNorm(norm);
214 forAll(types, celli)
215 {
216 if (types[celli] == cellCellStencil::HOLE)
217 {
218 extrapolatedNorm[celli] = -GREAT;
219 nHoles++;
220 }
221 }
222
223 bitSet isFront(nFaces());
224 for (label facei = 0; facei < nInternalFaces(); facei++)
225 {
226 label ownType = types[own[facei]];
227 label neiType = types[nei[facei]];
228 if
229 (
230 (ownType == cellCellStencil::HOLE)
231 != (neiType == cellCellStencil::HOLE)
232 )
233 {
234 isFront.set(facei);
235 }
236 }
237 labelList nbrTypes;
238 syncTools::swapBoundaryCellList(*this, types, nbrTypes);
239 for (label facei = nInternalFaces(); facei < nFaces(); facei++)
240 {
241 label ownType = types[own[facei]];
242 label neiType = nbrTypes[facei-nInternalFaces()];
243 if
244 (
245 (ownType == cellCellStencil::HOLE)
246 != (neiType == cellCellStencil::HOLE)
247 )
248 {
249 isFront.set(facei);
250 }
251 }
252
253
254 while (true)
255 {
256 scalarField nbrNorm;
257 syncTools::swapBoundaryCellList(*this, extrapolatedNorm, nbrNorm);
258
259 bitSet newIsFront(nFaces());
260 scalarField newNorm(extrapolatedNorm);
261
262 label nChanged = 0;
263 for (const label facei : isFront)
264 {
265 if (extrapolatedNorm[own[facei]] == -GREAT)
266 {
267 // Average owner cell, add faces to newFront
268 newNorm[own[facei]] = cellAverage
269 (
270 types,
271 nbrTypes,
272 extrapolatedNorm,
273 nbrNorm,
274 own[facei],
275 newIsFront
276 );
277 nChanged++;
278 }
279 if
280 (
281 isInternalFace(facei)
282 && extrapolatedNorm[nei[facei]] == -GREAT
283 )
284 {
285 // Average nei cell, add faces to newFront
286 newNorm[nei[facei]] = cellAverage
287 (
288 types,
289 nbrTypes,
290 extrapolatedNorm,
291 nbrNorm,
292 nei[facei],
293 newIsFront
294 );
295 nChanged++;
296 }
297 }
298
299 reduce(nChanged, sumOp<label>());
300 if (nChanged == 0)
301 {
302 break;
303 }
304
305 // Transfer new front
306 extrapolatedNorm.transfer(newNorm);
307 isFront.transfer(newIsFront);
309 }
310
311
312 forAll(norm, celli)
313 {
314 scalar& n = norm[celli];
315 if (magSqr(n) < sqr(SMALL))
316 {
317 //Pout<< "For field " << m.psi().name()
318 // << " for cell " << celli
319 // << " at:" << cellCentres()[celli]
320 // << " have norm " << n
321 // << " have extrapolated norm " << extrapolatedNorm[celli]
322 // << endl;
323 // Override the norm
324 n = extrapolatedNorm[celli];
325 }
326 }
327 }
328 return tnorm;
329}
330
331
332template<class Type>
334(
336 const scalarField& normalisation
337) const
338{
341 const labelListList& stencil = overlap.cellStencil();
344 const labelUList& types = overlap.cellTypes();
345
346
347 // Force asymmetric matrix (if it wasn't already)
348 scalarField& lower = m.lower();
349 scalarField& upper = m.upper();
350 Field<Type>& source = m.source();
351 scalarField& diag = m.diag();
352
353
354 // Get the addressing. Note that the addressing is now extended with
355 // any interpolation faces.
356 const lduAddressing& addr = lduAddr();
357 const labelUList& upperAddr = addr.upperAddr();
358 const labelUList& lowerAddr = addr.lowerAddr();
359 const labelUList& ownerStartAddr = addr.ownerStartAddr();
360 const labelUList& losortAddr = addr.losortAddr();
361 const lduInterfacePtrsList& interfaces = allInterfaces_;
362
363 if (!isA<fvMeshPrimitiveLduAddressing>(addr))
364 {
366 << "Problem : addressing is not fvMeshPrimitiveLduAddressing"
367 << exit(FatalError);
368 }
369
370
371
372 // 1. Adapt lduMatrix for additional faces and new ordering
373 upper.setSize(upperAddr.size(), 0.0);
374 inplaceReorder(reverseFaceMap_, upper);
375 lower.setSize(lowerAddr.size(), 0.0);
376 inplaceReorder(reverseFaceMap_, lower);
377
378
379 //const label nOldInterfaces = dynamicMotionSolverFvMesh::interfaces().size();
380 const label nOldInterfaces = dynamicFvMesh::interfaces().size();
381
382
383 if (interfaces.size() > nOldInterfaces)
384 {
385 // Extend matrix coefficients
386 m.internalCoeffs().setSize(interfaces.size());
387 m.boundaryCoeffs().setSize(interfaces.size());
388
389 // 1b. Adapt for additional interfaces
390 for
391 (
392 label patchi = nOldInterfaces;
393 patchi < interfaces.size();
394 patchi++
395 )
396 {
397 const labelUList& fc = interfaces[patchi].faceCells();
398 m.internalCoeffs().set(patchi, new Field<Type>(fc.size(), Zero));
399 m.boundaryCoeffs().set(patchi, new Field<Type>(fc.size(), Zero));
400 }
401
402 // 1c. Adapt field for additional interfaceFields (note: solver uses
403 // GeometricField::scalarInterfaces() to get hold of interfaces)
405
406 typename GeoField::Boundary& bfld =
407 const_cast<GeoField&>(m.psi()).boundaryFieldRef();
408
409 bfld.setSize(interfaces.size());
410
411
412 // This gets quite interesting: we do not want to add additional
413 // fvPatches (since direct correspondence to polyMesh) so instead
414 // add a reference to an existing processor patch
415 label addPatchi = 0;
416 for (label patchi = 0; patchi < nOldInterfaces; patchi++)
417 {
418 if (isA<processorFvPatch>(bfld[patchi].patch()))
419 {
420 addPatchi = patchi;
421 break;
422 }
423 }
424
425 for
426 (
427 label patchi = nOldInterfaces;
428 patchi < interfaces.size();
429 patchi++
430 )
431 {
432 bfld.set
433 (
434 patchi,
436 (
437 interfaces[patchi],
438 bfld[addPatchi].patch(), // dummy processorFvPatch
439 m.psi()
440 )
441 );
442 }
443 }
444
445
446 // 2. Adapt fvMatrix level: faceFluxCorrectionPtr
447 // Question: do we need to do this?
448 // This seems to be set/used only by the gaussLaplacianScheme and
449 // fvMatrix:correction, both of which are outside the linear solver.
450
451
452
453 // Clear out existing connections on cells to be interpolated
454 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
455 // Note: could avoid doing the zeroing of the new faces since these
456 // are set to zero anyway.
457
458 forAll(upperAddr, facei)
459 {
460 if (types[upperAddr[facei]] == cellCellStencil::INTERPOLATED)
461 {
462 // Disconnect upper from lower
463 label celli = upperAddr[facei];
464 lower[facei] *= 1.0-factor[celli];
465 }
466 if (types[lowerAddr[facei]] == cellCellStencil::INTERPOLATED)
467 {
468 // Disconnect lower from upper
469 label celli = lowerAddr[facei];
470 upper[facei] *= 1.0-factor[celli];
471 }
472 }
473
474 for (label patchi = 0; patchi < nOldInterfaces; ++patchi)
475 {
476 const labelUList& fc = addr.patchAddr(patchi);
477 Field<Type>& intCoeffs = m.internalCoeffs()[patchi];
478 Field<Type>& bouCoeffs = m.boundaryCoeffs()[patchi];
479 forAll(fc, i)
480 {
481 label celli = fc[i];
482 {
483 if (types[celli] == cellCellStencil::INTERPOLATED)
484 {
485 scalar f = factor[celli];
486 intCoeffs[i] *= 1.0-f;
487 bouCoeffs[i] *= 1.0-f;
488 }
489 else if (types[celli] == cellCellStencil::HOLE)
490 {
491 intCoeffs[i] = pTraits<Type>::zero;
492 bouCoeffs[i] = pTraits<Type>::zero;
493 }
494 }
495 }
496 }
497
498
499
500 // Modify matrix
501 // ~~~~~~~~~~~~~
502
503 // Do hole cells. Note: maybe put into interpolationCells() loop above?
504 forAll(types, celli)
505 {
506 if (types[celli] == cellCellStencil::HOLE)
507 {
508 label startLabel = ownerStartAddr[celli];
509 label endLabel = ownerStartAddr[celli + 1];
510
511 for (label facei = startLabel; facei < endLabel; facei++)
512 {
513 upper[facei] = 0.0;
514 }
515
516 startLabel = addr.losortStartAddr()[celli];
517 endLabel = addr.losortStartAddr()[celli + 1];
518
519 for (label i = startLabel; i < endLabel; i++)
520 {
521 label facei = losortAddr[i];
522 lower[facei] = 0.0;
523 }
524
525 diag[celli] = normalisation[celli];
526 source[celli] = normalisation[celli]*m.psi()[celli];
527 }
528 }
529
530
531 //const globalIndex globalNumbering(V().size());
532 //labelList globalCellIDs(overlap.cellInterpolationMap().constructSize());
533 //forAll(V(), cellI)
534 //{
535 // globalCellIDs[cellI] = globalNumbering.toGlobal(cellI);
536 //}
537 //overlap.cellInterpolationMap().distribute(globalCellIDs);
538
539
540 forAll(cellIDs, i)
541 {
542 label celli = cellIDs[i];
543
544 const scalar f = factor[celli];
545 const scalarList& w = wghts[celli];
546 const labelList& nbrs = stencil[celli];
547 const labelList& nbrFaces = stencilFaces_[celli];
548 const labelList& nbrPatches = stencilPatches_[celli];
549
550 if (types[celli] == cellCellStencil::HOLE)
551 {
552 FatalErrorInFunction << "Found HOLE cell " << celli
553 << " at:" << C()[celli]
554 << " . Should this be in interpolationCells()????"
555 << abort(FatalError);
556 }
557 else
558 {
559 // Create interpolation stencil
560
561 diag[celli] *= (1.0-f);
562 source[celli] *= (1.0-f);
563
564 forAll(nbrs, nbri)
565 {
566 label patchi = nbrPatches[nbri];
567 label facei = nbrFaces[nbri];
568
569 if (patchi == -1)
570 {
571 label nbrCelli = nbrs[nbri];
572
573 // Add the coefficients
574 const scalar s = normalisation[celli]*f*w[nbri];
575
576 scalar& u = upper[facei];
577 scalar& l = lower[facei];
578 if (celli < nbrCelli)
579 {
580 diag[celli] += s;
581 u += -s;
582 }
583 else
584 {
585 diag[celli] += s;
586 l += -s;
587 }
588 }
589 else
590 {
591 // Patch face. Store in boundaryCoeffs. Note sign change.
592 //const label globalCelli = globalCellIDs[nbrs[nbri]];
593 //const label proci =
594 // globalNumbering.whichProcID(globalCelli);
595 //const label remoteCelli =
596 // globalNumbering.toLocal(proci, globalCelli);
597 //
598 //Pout<< "for cell:" << celli
599 // << " need weight from remote slot:" << nbrs[nbri]
600 // << " proc:" << proci << " remote cell:" << remoteCelli
601 // << " patch:" << patchi
602 // << " patchFace:" << facei
603 // << " weight:" << w[nbri]
604 // << endl;
605
606 const scalar s = normalisation[celli]*f*w[nbri];
607 m.boundaryCoeffs()[patchi][facei] += pTraits<Type>::one*s;
608 m.internalCoeffs()[patchi][facei] += pTraits<Type>::one*s;
609
610 // Note: do NOT add to diagonal - this is in the
611 // internalCoeffs and gets added to the diagonal
612 // inside fvMatrix::solve
613 }
614 }
615
616 //if (mag(diag[celli]) < SMALL)
617 //{
618 // Pout<< "for cell:" << celli
619 // << " at:" << this->C()[celli]
620 // << " diag:" << diag[celli] << endl;
621 //
622 // forAll(nbrs, nbri)
623 // {
624 // label patchi = nbrPatches[nbri];
625 // label facei = nbrFaces[nbri];
626 //
627 // const label globalCelli = globalCellIDs[nbrs[nbri]];
628 // const label proci =
629 // globalNumbering.whichProcID(globalCelli);
630 // const label remoteCelli =
631 // globalNumbering.toLocal(proci, globalCelli);
632 //
633 // Pout<< " need weight from slot:" << nbrs[nbri]
634 // << " proc:" << proci << " remote cell:"
635 // << remoteCelli
636 // << " patch:" << patchi
637 // << " patchFace:" << facei
638 // << " weight:" << w[nbri]
639 // << endl;
640 // }
641 // Pout<< endl;
642 //}
643 }
644 }
645}
646
647
648template<class Type>
650(
652 const dictionary& dict
653) const
654{
656 // Check we're running with bcs that can handle implicit matrix manipulation
657 typename GeoField::Boundary& bpsi =
658 const_cast<GeoField&>(m.psi()).boundaryFieldRef();
659
660 bool hasOverset = false;
661 forAll(bpsi, patchi)
662 {
663 if (isA<oversetFvPatchField<Type>>(bpsi[patchi]))
664 {
665 hasOverset = true;
666 break;
667 }
668 }
669
670 if (!hasOverset)
671 {
672 if (debug)
673 {
674 Pout<< "dynamicOversetFvMesh::solve() :"
675 << " bypassing matrix adjustment for field " << m.psi().name()
676 << endl;
677 }
678 //return dynamicMotionSolverFvMesh::solve(m, dict);
679 return dynamicFvMesh::solve(m, dict);
680 }
681
682 if (debug)
683 {
684 Pout<< "dynamicOversetFvMesh::solve() :"
685 << " adjusting matrix for interpolation for field "
686 << m.psi().name() << endl;
687 }
688
689 // Calculate stabilised diagonal as normalisation for interpolation
690 const scalarField norm(normalisation(m));
691
692 if (debug)
693 {
694 volScalarField scale
695 (
697 (
698 m.psi().name() + "_normalisation",
699 this->time().timeName(),
700 *this,
703 false
704 ),
705 *this,
707 );
708 scale.ref().field() = norm;
710 <
713 >(scale.boundaryFieldRef(), false);
714 scale.write();
715
716 if (debug)
717 {
718 Pout<< "dynamicOversetFvMesh::solve() :"
719 << " writing matrix normalisation for field " << m.psi().name()
720 << " to " << scale.name() << endl;
721 }
722 }
723
724
725 // Switch to extended addressing (requires mesh::update() having been
726 // called)
727 active(true);
728
729 // Adapt matrix
730 scalarField oldUpper(m.upper());
731 scalarField oldLower(m.lower());
734 const label oldSize = bpsi.size();
735
736 addInterpolation(m, norm);
737
738 // Swap psi values so added patches have patchNeighbourField
739 correctBoundaryConditions<GeoField, calculatedProcessorFvPatchField<Type>>
740 (
741 bpsi,
742 true
743 );
744
745
746 // Print a bit
747 //write(Pout, m, lduAddr(), interfaces());
748 //{
749 // const fvSolution& sol = static_cast<const fvSolution&>(*this);
750 // const dictionary& pDict = sol.subDict("solvers").subDict("p");
751 // writeAgglomeration(GAMGAgglomeration::New(m, pDict));
752 //}
753
754 // Use lower level solver
755 //SolverPerformance<Type> s(dynamicMotionSolverFvMesh::solve(m, dict));
757
758 // Restore boundary field
759 bpsi.setSize(oldSize);
760
761 // Restore matrix
762 m.upper().transfer(oldUpper);
763 m.lower().transfer(oldLower);
764 m.internalCoeffs().transfer(oldInt);
765 m.boundaryCoeffs().transfer(oldBou);
766
767 // Switch to original addressing
768 active(false);
769
770 return s;
771}
772
773
774template<class Type>
776(
777 Ostream& os,
778 const fvMatrix<Type>& m,
779 const lduAddressing& addr,
780 const lduInterfacePtrsList& interfaces
781) const
782{
783 os << "** Matrix **" << endl;
784 const labelUList& upperAddr = addr.upperAddr();
785 const labelUList& lowerAddr = addr.lowerAddr();
786 const labelUList& ownerStartAddr = addr.ownerStartAddr();
787 const labelUList& losortAddr = addr.losortAddr();
788
789 const scalarField& lower = m.lower();
790 const scalarField& upper = m.upper();
791 const Field<Type>& source = m.source();
792 const scalarField& diag = m.diag();
793
794
795 // Invert patch addressing
796 labelListList cellToPatch(addr.size());
797 labelListList cellToPatchFace(addr.size());
798 {
799 forAll(interfaces, patchi)
800 {
801 if (interfaces.set(patchi))
802 {
803 const labelUList& fc = interfaces[patchi].faceCells();
804
805 forAll(fc, i)
806 {
807 cellToPatch[fc[i]].append(patchi);
808 cellToPatchFace[fc[i]].append(i);
809 }
810 }
811 }
812 }
813
814 forAll(source, celli)
815 {
816 os << "cell:" << celli << " diag:" << diag[celli]
817 << " source:" << source[celli] << endl;
818
819 label startLabel = ownerStartAddr[celli];
820 label endLabel = ownerStartAddr[celli + 1];
821
822 for (label facei = startLabel; facei < endLabel; facei++)
823 {
824 os << " face:" << facei
825 << " nbr:" << upperAddr[facei] << " coeff:" << upper[facei]
826 << endl;
827 }
828
829 startLabel = addr.losortStartAddr()[celli];
830 endLabel = addr.losortStartAddr()[celli + 1];
831
832 for (label i = startLabel; i < endLabel; i++)
833 {
834 label facei = losortAddr[i];
835 os << " face:" << facei
836 << " nbr:" << lowerAddr[facei] << " coeff:" << lower[facei]
837 << endl;
838 }
839
840 forAll(cellToPatch[celli], i)
841 {
842 label patchi = cellToPatch[celli][i];
843 label patchFacei = cellToPatchFace[celli][i];
844
845 os << " patch:" << patchi
846 << " patchface:" << patchFacei
847 << " intcoeff:" << m.internalCoeffs()[patchi][patchFacei]
848 << " boucoeff:" << m.boundaryCoeffs()[patchi][patchFacei]
849 << endl;
850 }
851 }
852 forAll(m.internalCoeffs(), patchi)
853 {
854 if (m.internalCoeffs().set(patchi))
855 {
856 const labelUList& fc = addr.patchAddr(patchi);
857
858 os << "patch:" << patchi
859 //<< " type:" << interfaces[patchi].type()
860 << " size:" << fc.size() << endl;
861 if
862 (
863 interfaces.set(patchi)
864 && isA<processorLduInterface>(interfaces[patchi])
865 )
866 {
867 const processorLduInterface& ppp =
868 refCast<const processorLduInterface>(interfaces[patchi]);
869 os << "(processor with my:" << ppp.myProcNo()
870 << " nbr:" << ppp.neighbProcNo()
871 << ")" << endl;
872 }
873
874 forAll(fc, i)
875 {
876 os << " cell:" << fc[i]
877 << " int:" << m.internalCoeffs()[patchi][i]
878 << " boun:" << m.boundaryCoeffs()[patchi][i]
879 << endl;
880 }
881 }
882 }
883
884
885 lduInterfaceFieldPtrsList interfaceFields =
886 m.psi().boundaryField().scalarInterfaces();
887 forAll(interfaceFields, inti)
888 {
889 if (interfaceFields.set(inti))
890 {
891 os << "interface:" << inti
892 << " if type:" << interfaceFields[inti].interface().type()
893 << " fld type:" << interfaceFields[inti].type() << endl;
894 }
895 }
896
897 os << "** End of Matrix **" << endl;
898}
899
900
901template<class GeoField>
903{
904 typename GeoField::Boundary& bfld = fld.boundaryFieldRef();
905
906 const label startOfRequests = UPstream::nRequests();
907
908 forAll(bfld, patchi)
909 {
910 if (bfld[patchi].coupled())
911 {
912 //Pout<< "initEval of " << bfld[patchi].patch().name() << endl;
913 bfld[patchi].initEvaluate(Pstream::defaultCommsType);
914 }
915 }
916
917 // Wait for outstanding requests
918 if
919 (
922 )
923 {
924 UPstream::waitRequests(startOfRequests);
925 }
926
927 forAll(bfld, patchi)
928 {
929 if (bfld[patchi].coupled())
930 {
931 //Pout<< "eval of " << bfld[patchi].patch().name() << endl;
932 bfld[patchi].evaluate(Pstream::defaultCommsType);
933 }
934 }
935}
936
937
938template<class GeoField>
940{
941 Pout<< "** starting checking of " << fld.name() << endl;
942
943 GeoField fldCorr(fld.name()+"_correct", fld);
944 correctCoupledBoundaryConditions(fldCorr);
945
946 const typename GeoField::Boundary& bfld = fld.boundaryField();
947 const typename GeoField::Boundary& bfldCorr = fldCorr.boundaryField();
948
949 forAll(bfld, patchi)
950 {
951 const typename GeoField::Patch& pfld = bfld[patchi];
952 const typename GeoField::Patch& pfldCorr = bfldCorr[patchi];
953
954 Pout<< "Patch:" << pfld.patch().name() << endl;
955
956 forAll(pfld, i)
957 {
958 if (pfld[i] != pfldCorr[i])
959 {
960 Pout<< " " << i << " orig:" << pfld[i]
961 << " corrected:" << pfldCorr[i] << endl;
962 }
963 }
964 }
965 Pout<< "** end of " << fld.name() << endl;
966}
967
968
969// ************************************************************************* //
label n
Info<< nl<< "Wrote faMesh in vtk format: "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.mesh().time().globalPath()/"finiteArea-edges"));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
labelList cellIDs
Graphite solid properties.
Definition: C.H:53
const Field< Type > & field() const
Return field.
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:80
Generic templated field type.
Definition: Field.H:82
tmp< Field< cmptType > > component(const direction) const
Return a component field of the field.
Definition: Field.C:545
Generic GeometricField class.
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
void correctBoundaryConditions()
Correct boundary field.
bool found(const Key &key) const
Return true if hashed entry is found in table.
Definition: HashTableI.H:100
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:65
void transfer(List< T > &list)
Definition: List.C:447
void setSize(const label n)
Alias for resize()
Definition: List.H:218
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
virtual const fileName & name() const
Get the name of the stream.
Definition: OSstream.H:107
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics.
static autoPtr< Time > New()
Construct (dummy) Time - no functionObjects or libraries.
Definition: Time.C:717
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
@ nonBlocking
"nonBlocking"
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:556
static label nRequests()
Get number of outstanding requests.
Definition: UPstream.C:90
static void waitRequests(const label start=0)
Wait until all requests (from start onwards) have finished.
Definition: UPstream.C:100
static commsTypes defaultCommsType
Default commsType.
Definition: UPstream.H:281
static bool & parRun() noexcept
Test if this a parallel run.
Definition: UPstream.H:433
const T * set(const label i) const
Definition: UPtrList.H:248
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:590
A processorFvPatchField type bypassing fvPatch.
virtual const labelUList & interpolationCells() const
Indices of interpolated cells.
virtual const mapDistribute & cellInterpolationMap() const
Return a communication schedule.
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.
Calculation of interpolation stencils.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
static void checkCoupledBC(const GeoField &fld)
Debug: check halo swap is ok.
void addInterpolation(fvMatrix< Type > &, const scalarField &norm) const
Add interpolation to matrix (coefficients)
tmp< scalarField > normalisation(const fvMatrix< Type > &m) const
Freeze values at holes.
static void correctCoupledBoundaryConditions(GeoField &fld)
Debug: correct coupled bc.
lduInterfacePtrsList interfaces() const
virtual bool write()
Write the output fields.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:121
const FieldField< Field, Type > & internalCoeffs() const noexcept
Definition: fvMatrix.H:470
const FieldField< Field, Type > & boundaryCoeffs() const noexcept
Definition: fvMatrix.H:484
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:412
Field< Type > & source() noexcept
Definition: fvMatrix.H:458
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.
label size() const
Return number of equations.
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 labelUList & losortAddr() const
Return losort addressing.
scalarField & upper()
Definition: lduMatrix.C:203
scalarField & lower()
Definition: lduMatrix.C:174
scalarField & diag()
Definition: lduMatrix.C:192
Class containing processor-to-processor mapping information.
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition: one.H:62
Boundary condition for use on overset patches. To be run in combination with special dynamicFvMesh ty...
label nCells() const noexcept
Number of mesh cells.
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 valid=true) const
Write using setting from DB.
bool interpolate() const noexcept
Same as isPointData()
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:396
static void swapBoundaryCellList(const polyMesh &mesh, const UList< T > &cellData, List< T > &neighbourCellData)
Swap to obtain neighbour cell values for all boundary faces.
A class for managing temporary objects.
Definition: tmp.H:65
T & ref() const
Definition: tmpI.H:227
A class for handling words, derived from Foam::string.
Definition: word.H:68
const volScalarField & psi
const volScalarField & T
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
OBJstream os(runTime.globalPath()/outputName)
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
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
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
const TargetType * isA(const Type &t)
Check if dynamic_cast to TargetType is possible.
Definition: typeInfo.H:197
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
void inplaceReorder(const labelUList &oldToNew, ListType &input, const bool prune=false)
Inplace reorder the elements of a list.
void reduce(const List< UPstream::commsStruct > &comms, T &value, const BinaryOp &bop, const int tag, const label comm)
errorManip< error > abort(error &err)
Definition: errorManip.H:144
uint8_t direction
Definition: direction.H:56
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
bool interpolate(const vector &p1, const vector &p2, const vector &o, vector &n, scalar l)
Definition: curveTools.C:75
error FatalError
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
labelList f(nPoints)
dictionary dict
cellMask correctBoundaryConditions()
CEqn solve()
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
A non-counting (dummy) refCount.
Definition: refCount.H:59