Loading...
Searching...
No Matches
distributedDILUPreconditioner.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) 2023-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
30#include "cyclicLduInterface.H"
32#include "processorColour.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36namespace Foam
37{
39
40 lduMatrix::preconditioner::
41 addasymMatrixConstructorToTable<distributedDILUPreconditioner>
43
47}
48
49
50// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
51
53(
54 const bool add,
55 const FieldField<Field, scalar>& coupleCoeffs,
56 const labelList& selectedInterfaces,
57 const solveScalarField& psiif,
58 solveScalarField& result,
59 const direction cmpt
60) const
61{
62 const auto& matrix = solver_.matrix();
63 const auto& lduAddr = matrix.lduAddr();
64 const auto& interfaces = solver_.interfaces();
65
66 const label startOfRequests = Pstream::nRequests();
67
68 // From lduMatrix::initMatrixInterfaces
69
70 // Avoid any conflicts with inter-processor comms
71 const int oldTag = UPstream::incrMsgType(321);
72
73 for (const label interfacei : selectedInterfaces)
74 {
75 interfaces[interfacei].initInterfaceMatrixUpdate
76 (
77 result,
78 add,
79 lduAddr,
80 interfacei,
81 psiif,
82 coupleCoeffs[interfacei],
83 cmpt,
85 );
86 }
87
88 UPstream::waitRequests(startOfRequests);
89
90 // Consume
91 for (const label interfacei : selectedInterfaces)
92 {
93 interfaces[interfacei].updateInterfaceMatrix
94 (
95 result,
96 add,
97 lduAddr,
98 interfacei,
99 psiif,
100 coupleCoeffs[interfacei],
101 cmpt,
103 );
104 }
105
106 UPstream::msgType(oldTag); // Restore tag
107}
108
109
111(
112 const labelList& selectedInterfaces,
114 const label colouri
115) const
116{
117 const auto& interfaces = solver_.interfaces();
118
119 if (selectedInterfaces.size())
120 {
121 // Save old data
122 FieldField<Field, scalar> one(interfaces.size());
123 FieldField<Field, solveScalar> old(interfaces.size());
124 for (const label inti : selectedInterfaces)
125 {
126 const auto& intf = interfaces[inti].interface();
127 const auto& fc = intf.faceCells();
128 one.set(inti, new scalarField(fc.size(), 1.0));
129 old.set(inti, new solveScalarField(psi, intf.faceCells()));
130 }
131 updateMatrixInterfaces
132 (
133 false, // add to psi
134 one,
135 selectedInterfaces,
136 psi, // send data
137 psi, // result
138 0 // cmpt
139 );
140
141 if (!colourBufs_.set(colouri))
142 {
143 colourBufs_.set
144 (
145 colouri,
147 (
148 interfaces.size()
149 )
150 );
151 }
152 auto& colourBuf = colourBufs_[colouri];
153 colourBuf.setSize(interfaces.size());
154 for (const label inti : selectedInterfaces)
155 {
156 const auto& intf = interfaces[inti].interface();
157 const auto& fc = intf.faceCells();
158 if (!colourBuf.set(inti))
159 {
160 colourBuf.set(inti, new Field<solveScalar>(fc.size()));
161 }
162 auto& cb = colourBuf[inti];
163 //cb.resize_nocopy(fc.size());
164 auto& oldValues = old[inti];
165
166 forAll(cb, face)
167 {
168 const label cell = fc[face];
169 // Store change in boundary value
170 cb[face] = psi[cell]-oldValues[face];
171 // Restore old value
172 std::swap(psi[cell], oldValues[face]);
173 }
174 }
175 }
176}
177
178
180(
181 const labelList& selectedInterfaces,
183) const
184{
185 const auto& interfaces = solver_.interfaces();
186 const auto& interfaceBouCoeffs = solver_.interfaceBouCoeffs();
187
188 // Start reads
189 for (const label inti : selectedInterfaces)
190 {
191 const auto& intf = interfaces[inti].interface();
192 const auto* ppp = isA<processorLduInterface>(intf);
193
194 auto& recvBuf = recvBufs_[inti];
195 recvBuf.resize_nocopy(interfaceBouCoeffs[inti].size());
196
198 (
199 requests.emplace_back(),
200 ppp->neighbProcNo(),
201 recvBuf,
202 ppp->tag()+70, // random offset
203 ppp->comm()
204 );
205 }
206}
207
208
210(
211 const labelList& selectedInterfaces,
212 const solveScalarField& psiInternal,
213 DynamicList<UPstream::Request>& requests
214) const
215{
216 const auto& interfaces = solver_.interfaces();
217
218 // Start writes
219 for (const label inti : selectedInterfaces)
220 {
221 const auto& intf = interfaces[inti].interface();
222 const auto* ppp = isA<processorLduInterface>(intf);
223 const auto& faceCells = intf.faceCells();
224
225 auto& sendBuf = sendBufs_[inti];
226
227 sendBuf.resize_nocopy(faceCells.size());
228 forAll(faceCells, face)
229 {
230 sendBuf[face] = psiInternal[faceCells[face]];
231 }
232
234 (
235 requests.emplace_back(),
236 ppp->neighbProcNo(),
237 sendBuf,
238 ppp->tag()+70, // random offset
239 ppp->comm()
240 );
241 }
242}
243
244
246(
247 DynamicList<UPstream::Request>& requests,
248 const bool cancel
249) const
250{
251 if (cancel)
252 {
253 UPstream::cancelRequests(requests);
254 }
255 else
256 {
257 UPstream::waitRequests(requests);
258 }
259 requests.clear();
260}
261
262
263// Diagonal calculation
264// ~~~~~~~~~~~~~~~~~~~~
265
267(
269 const label inti,
270 const Field<solveScalar>& recvBuf
271) const
272{
273 const auto& interfaces = solver_.interfaces();
274 const auto& interfaceBouCoeffs = solver_.interfaceBouCoeffs();
275 const auto& interfaceIntCoeffs = solver_.interfaceIntCoeffs();
276
277 const auto& intf = interfaces[inti].interface();
278 // TBD: do not use patch faceCells but passed-in addressing?
279 const auto& faceCells = intf.faceCells();
280 const auto& bc = interfaceBouCoeffs[inti];
281 const auto& ic = interfaceIntCoeffs[inti];
282
283 forAll(recvBuf, face)
284 {
285 // Note:interfaceBouCoeffs, interfaceIntCoeffs on the receiving
286 // side are negated
287 rD[faceCells[face]] -= bc[face]*ic[face]/recvBuf[face];
288 }
289}
290
291
293(
295 const label colouri
296) const
297{
298 // Add forward constributions to diagonal
299
300 const auto& matrix = solver_.matrix();
301 const auto& lduAddr = matrix.lduAddr();
302
303 const label* const __restrict__ uPtr = lduAddr.upperAddr().begin();
304 const label* const __restrict__ lPtr = lduAddr.lowerAddr().begin();
305 const scalar* const __restrict__ upperPtr = matrix.upper().begin();
306 const scalar* const __restrict__ lowerPtr = matrix.lower().begin();
307
308 const label nFaces = matrix.upper().size();
309 if (cellColourPtr_)
310 {
311 const auto& cellColour = *cellColourPtr_;
312 for (label face=0; face<nFaces; face++)
313 {
314 const label cell = lPtr[face];
315 if (cellColour[cell] == colouri)
316 {
317 rD[uPtr[face]] -= upperPtr[face]*lowerPtr[face]/rD[cell];
318 }
319 }
320 }
321 else
322 {
323 for (label face=0; face<nFaces; face++)
325 rD[uPtr[face]] -= upperPtr[face]*lowerPtr[face]/rD[lPtr[face]];
326 }
327 }
328}
329
330
332(
334 const label inti,
335 const Field<solveScalar>& recvBuf
336) const
337{
338 const auto& interfaces = solver_.interfaces();
339 const auto& interfaceBouCoeffs = solver_.interfaceBouCoeffs();
340
341 const auto& intf = interfaces[inti].interface();
342 const auto& faceCells = intf.faceCells();
343 const auto& bc = interfaceBouCoeffs[inti];
344 const solveScalar* __restrict__ rDPtr = rD_.begin();
345 solveScalar* __restrict__ wAPtr = wA.begin();
346
347 forAll(recvBuf, face)
348 {
349 // Note: interfaceBouCoeffs is -upperPtr
350 const label cell = faceCells[face];
351 wAPtr[cell] += rDPtr[cell]*bc[face]*recvBuf[face];
352 }
353}
354
355
357(
359 const label colouri
360) const
361{
362 const auto& matrix = solver_.matrix();
363 const auto& lduAddr = matrix.lduAddr();
364
365 solveScalar* __restrict__ wAPtr = wA.begin();
366 const solveScalar* __restrict__ rDPtr = rD_.begin();
367
368 const label* const __restrict__ uPtr = lduAddr.upperAddr().begin();
369 const label* const __restrict__ lPtr = lduAddr.lowerAddr().begin();
370 const label* const __restrict__ losortPtr = lduAddr.losortAddr().begin();
371 const scalar* const __restrict__ lowerPtr = matrix.lower().begin();
372
373 const label nFaces = matrix.upper().size();
374 if (cellColourPtr_)
375 {
376 const auto& cellColour = *cellColourPtr_;
377 for (label face=0; face<nFaces; face++)
378 {
379 const label sface = losortPtr[face];
380 const label cell = lPtr[sface];
381 if (cellColour[cell] == colouri)
382 {
383 wAPtr[uPtr[sface]] -=
384 rDPtr[uPtr[sface]]*lowerPtr[sface]*wAPtr[cell];
385 }
386 }
387 }
388 else
389 {
390 for (label face=0; face<nFaces; face++)
391 {
392 const label sface = losortPtr[face];
393 wAPtr[uPtr[sface]] -=
394 rDPtr[uPtr[sface]]*lowerPtr[sface]*wAPtr[lPtr[sface]];
395 }
396 }
397}
398
399
401(
403 const label colouri
404) const
405{
406 const auto& matrix = solver_.matrix();
407 const auto& lduAddr = matrix.lduAddr();
408
409 solveScalar* __restrict__ wAPtr = wA.begin();
410 const solveScalar* __restrict__ rDPtr = rD_.begin();
411
412 const label* const __restrict__ uPtr = lduAddr.upperAddr().begin();
413 const label* const __restrict__ lPtr = lduAddr.lowerAddr().begin();
414 const scalar* const __restrict__ upperPtr = matrix.upper().begin();
415
416 const label nFacesM1 = matrix.upper().size() - 1;
417
418 if (cellColourPtr_)
419 {
420 const auto& cellColour = *cellColourPtr_;
421 for (label face=nFacesM1; face>=0; face--)
422 {
423 const label cell = uPtr[face];
424 if (cellColour[cell] == colouri)
425 {
426 // Note: lower cell guaranteed in same colour
427 wAPtr[lPtr[face]] -=
428 rDPtr[lPtr[face]]*upperPtr[face]*wAPtr[cell];
429 }
430 }
431 }
432 else
433 {
434 for (label face=nFacesM1; face>=0; face--)
435 {
436 wAPtr[lPtr[face]] -=
437 rDPtr[lPtr[face]]*upperPtr[face]*wAPtr[uPtr[face]];
438 }
439 }
440}
441
442
444(
446) const
447{
448 const auto& matrix = solver_.matrix();
449
450 // Make sure no outstanding receives
451 wait(lowerRecvRequests_);
452
453 // Start reads (into recvBufs_)
454 receive(lowerNbrs_, lowerRecvRequests_);
455
456 // Start with diagonal
457 const scalarField& diag = matrix.diag();
458 rD.resize_nocopy(diag.size());
459 std::copy(diag.begin(), diag.end(), rD.begin());
460
461
462 // Subtract coupled contributions
463 {
464 // Wait for finish. Received result in recvBufs
465 wait(lowerRecvRequests_);
466
467 for (const label inti : lowerNbrs_)
468 {
469 addInterfaceDiag(rD, inti, recvBufs_[inti]);
470 }
471 }
472
473
474 // - divide the internal mesh into domains/colour, similar to domain
475 // decomposition
476 // - we could use subsetted bits of the mesh or just loop over only
477 // the cells of the current domain
478 // - a domain only uses boundary values of lower numbered domains
479 // - this also limits the interfaces that need to be evaluated since
480 // we assume an interface only changes local faceCells and these
481 // are all of the same colour
482
483 for (label colouri = 0; colouri < nColours_; colouri++)
484 {
485 if (cellColourPtr_) // && colourBufs_.set(colouri))
486 {
487 for (const label inti : lowerGlobalRecv_[colouri])
488 {
489 addInterfaceDiag(rD, inti, colourBufs_[colouri][inti]);
490 }
491 }
492
493 forwardInternalDiag(rD, colouri);
494
495 // Store effect of exchanging rD to higher interfaces in colourBufs_
496 if (cellColourPtr_)
497 {
498 sendGlobal(higherGlobalSend_[colouri], rD, higherColour_[colouri]);
499 }
500 }
501
502 // Make sure no outstanding sends
503 wait(higherSendRequests_);
504
505 // Start writes of rD (using sendBufs)
506 send(higherNbrs_, rD, higherSendRequests_);
507
508 // Calculate the reciprocal of the preconditioned diagonal
509 const label nCells = rD.size();
510
511 for (label cell=0; cell<nCells; cell++)
512 {
513 rD[cell] = 1.0/rD[cell];
514 }
515}
516
517
518// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
519
521(
522 const lduMatrix::solver& sol,
523 const dictionary& dict
524)
525:
526 lduMatrix::preconditioner(sol),
527 coupled_(dict.getOrDefault<bool>("coupled", true, keyType::LITERAL)),
528 rD_(sol.matrix().diag().size())
529{
530 const lduMesh& mesh = sol.matrix().mesh();
531 const auto& interfaces = sol.interfaces();
532 const auto& interfaceBouCoeffs = sol.interfaceBouCoeffs();
533
534 // Allocate buffers
535 // ~~~~~~~~~~~~~~~~
536
537 sendBufs_.setSize(interfaces.size());
538 recvBufs_.setSize(interfaces.size());
539 forAll(interfaceBouCoeffs, inti)
540 {
541 if (interfaces.set(inti))
542 {
543 const auto& bc = interfaceBouCoeffs[inti];
544 sendBufs_.set(inti, new solveScalarField(bc.size(), Zero));
545 recvBufs_.set(inti, new solveScalarField(bc.size(), Zero));
546 }
547 }
548
549
550 // Determine processor colouring
551 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
552
553 //const processorColour& colours = processorColour::New(mesh);
554 //const label myColour = colours[Pstream::myProcNo()];
555 if (&mesh != meshPtr_)
556 {
557 procColoursPtr_.reset(new labelList(Pstream::nProcs(mesh.comm())));
558 const label nProcColours =
560
561 if (debug)
562 {
563 Info<< typeName << " : calculated processor colours:"
564 << nProcColours << endl;
565 }
566
567 meshPtr_ = &mesh;
568 }
569 const labelList& procColours = procColoursPtr_();
570
571 const label myColour = procColours[Pstream::myProcNo(mesh.comm())];
572
573 bool haveGlobalCoupled = false;
574 bool haveDistributedAMI = false;
575 forAll(interfaces, inti)
576 {
577 if (interfaces.set(inti))
578 {
579 const auto& intf = interfaces[inti].interface();
580 const auto* ppp = isA<processorLduInterface>(intf);
581 if (ppp)
582 {
583 const label nbrColour = procColours[ppp->neighbProcNo()];
584 //const label nbrColour = ppp->neighbProcNo();
585 if (nbrColour < myColour)
586 {
587 lowerNbrs_.append(inti);
588 }
589 else if (nbrColour > myColour)
590 {
591 higherNbrs_.append(inti);
592 }
593 else
594 {
596 << "weird processorLduInterface"
597 << " from " << ppp->myProcNo()
598 << " to " << ppp->neighbProcNo()
599 << endl;
600 }
601 }
602 else //if (isA<cyclicAMILduInterface>(intf))
603 {
604 haveGlobalCoupled = true;
605
606 const auto* AMIpp = isA<cyclicAMILduInterface>(intf);
607 if (AMIpp)
608 {
609 //const auto& AMI =
610 //(
611 // AMIpp->owner()
612 // ? AMIpp->AMI()
613 // : AMIpp->neighbPatch().AMI()
614 //);
615 //haveDistributedAMI = AMI.distributed();
616
617 if (debug)
618 {
619 Info<< typeName << " : interface " << inti
620 << " of type " << intf.type()
621 << " is distributed:" << haveDistributedAMI << endl;
622 }
623 }
624 }
625 }
626 }
627
628
629
630 if (debug)
631 {
632 Info<< typeName
633 << " : lower proc interfaces:" << flatOutput(lowerNbrs_)
634 << " higher proc interfaces:" << flatOutput(higherNbrs_)
635 << endl;
636 }
637
638
639 // Determine local colouring/zoneing
640 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
641
642 // Start off with all cells colour 0
643 nColours_ = 1;
644 cellColourPtr_.reset();
645 if (coupled_ && haveGlobalCoupled)
646 {
647 labelList patchToColour;
648 cellColourPtr_.reset(new labelList(mesh.lduAddr().size()));
649 //if (haveDistributedAMI)
650 //{
651 // nColours_ = 1;
652 // *cellColourPtr_ = 0;
653 // patchToColour = labelList(interfaces.size(), 0);
654 //}
655 //else
656 {
658 (
659 mesh,
661 patchToColour
662 );
663 }
664
665 lowerGlobalRecv_.resize_nocopy(nColours_);
666 lowerGlobalSend_.resize_nocopy(nColours_);
667 lowerColour_.setSize(nColours_, -1);
668 higherGlobalRecv_.resize_nocopy(nColours_);
669 higherGlobalSend_.resize_nocopy(nColours_);
670 higherColour_.setSize(nColours_, -1);
671 colourBufs_.setSize(nColours_);
672
673 forAll(interfaces, inti)
674 {
675 if (interfaces.set(inti))
676 {
677 const auto& intf = interfaces[inti].interface();
678
679 label nbrInti = -1;
680 const auto* AMIpp = isA<cyclicAMILduInterface>(intf);
681 if (AMIpp)
682 {
683 nbrInti = AMIpp->neighbPatchID();
684 }
685 const auto* cycpp = isA<cyclicLduInterface>(intf);
686 if (cycpp)
687 {
688 nbrInti = cycpp->neighbPatchID();
689 }
690
691 if (nbrInti != -1)
692 {
693 const label colouri = patchToColour[inti];
694 const label nbrColouri = patchToColour[nbrInti];
695
696 if
697 (
698 (haveDistributedAMI && (inti < nbrInti))
699 || (!haveDistributedAMI && (colouri < nbrColouri))
700 )
701 {
702 // Send to higher
703 higherGlobalSend_[colouri].append(nbrInti);
704 higherColour_[colouri] = nbrColouri;
705 // Receive from higher
706 higherGlobalRecv_[colouri].append(inti);
707 }
708 else
709 {
710 // Send to lower
711 lowerGlobalSend_[colouri].append(nbrInti);
712 lowerColour_[colouri] = nbrColouri;
713 // Receive from lower
714 lowerGlobalRecv_[colouri].append(inti);
715 }
716 }
717 }
718 }
719 }
720
721
722 if (debug)
723 {
724 Info<< typeName << " : local colours:" << nColours_ << endl;
726 forAll(lowerGlobalRecv_, colouri)
727 {
728 const auto& lowerRecv = lowerGlobalRecv_[colouri];
729 if (lowerRecv.size())
730 {
731 const auto& lowerSend = lowerGlobalSend_[colouri];
732 const auto& lowerCol = lowerColour_[colouri];
733 Info<< indent
734 << "Lower global interfaces for colour:" << colouri
735 << " receiving from:" << flatOutput(lowerRecv)
736 << " sending to:" << flatOutput(lowerSend)
737 << " with colour:" << lowerCol << endl;
738 }
739 }
740 forAll(higherGlobalRecv_, colouri)
741 {
742 const auto& higherRecv = higherGlobalRecv_[colouri];
743 if (higherRecv.size())
744 {
745 const auto& higherSend = higherGlobalSend_[colouri];
746 const auto& higherCol = higherColour_[colouri];
747 Info<< indent
748 << "Higher global interfaces for colour:" << colouri
749 << " receiving from:" << flatOutput(higherRecv)
750 << " sending to:" << flatOutput(higherSend)
751 << " with colour:" << higherCol << endl;
752 }
753 }
754 }
757}
758
759
760// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
761
763{
764 DebugPout<< "~distributedDILUPreconditioner()" << endl;
765
766 // Wait on all requests before storage is being taken down
767
768 wait(lowerSendRequests_);
769 wait(lowerRecvRequests_);
770 wait(higherSendRequests_);
771 wait(higherRecvRequests_);
772
773 // TBD: cancel/ignore outstanding messages
774 //wait(lowerSendRequests_, true);
775 //wait(lowerRecvRequests_, true);
776 //wait(higherSendRequests_, true);
777 //wait(higherRecvRequests_, true);
778}
779
780
781// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
782
784(
786 const solveScalarField& rA,
787 const direction
788) const
789{
790 solveScalar* __restrict__ wAPtr = wA.begin();
791 const solveScalar* __restrict__ rAPtr = rA.begin();
792 const solveScalar* __restrict__ rDPtr = rD_.begin();
793
794 const label nCells = wA.size();
795
796
797 // Forward sweep
798 // ~~~~~~~~~~~~~
799
800 // Make sure no receives are still in flight (should not happen)
801 wait(lowerRecvRequests_);
802
803 // Start reads (into recvBufs)
804 receive(lowerNbrs_, lowerRecvRequests_);
805
806 // Initialise 'internal' cells
807 for (label cell=0; cell<nCells; cell++)
808 {
809 wAPtr[cell] = rDPtr[cell]*rAPtr[cell];
810 }
811
812 // Do 'halo' contributions from lower numbered procs
813 {
814 // Wait for finish. Received result in recvBufs
815 wait(lowerRecvRequests_);
816
817 for (const label inti : lowerNbrs_)
818 {
819 addInterface(wA, inti, recvBufs_[inti]);
820 }
821
822 for (label colouri = 0; colouri < nColours_; colouri++)
823 {
824 // Do non-processor boundaries for this colour
825 if (cellColourPtr_) // && colourBufs_.set(colouri))
826 {
827 for (const label inti : lowerGlobalRecv_[colouri])
828 {
829 addInterface(wA, inti, colourBufs_[colouri][inti]);
830 }
831 }
832
833 forwardInternal(wA, colouri);
834
835 // Store effect of exchanging rD to higher interfaces
836 // in colourBufs_
837 if (cellColourPtr_)
838 {
839 sendGlobal
840 (
841 higherGlobalSend_[colouri],
842 wA,
843 higherColour_[colouri]
844 );
845 }
846 }
847 }
848
849 // Make sure no outstanding sends from previous iteration
850 wait(higherSendRequests_);
851
852 // Start writes of wA (using sendBufs)
853 send(higherNbrs_, wA, higherSendRequests_);
854
855
856 // Backward sweep
857 // ~~~~~~~~~~~~~~
858
859 // Make sure no outstanding receives
860 wait(higherRecvRequests_);
861
862 // Start receives
863 receive(higherNbrs_, higherRecvRequests_);
864
865 {
866 // Wait until receives finished
867 wait(higherRecvRequests_);
868
869 for (const label inti : higherNbrs_)
870 {
871 addInterface(wA, inti, recvBufs_[inti]);
872 }
873
874 for (label colouri = nColours_-1; colouri >= 0; colouri--)
875 {
876 // Do non-processor boundaries for this colour
877 if (cellColourPtr_) // && colourBufs_.set(colouri))
878 {
879 for (const label inti : higherGlobalRecv_[colouri])
880 {
881 addInterface(wA, inti, colourBufs_[colouri][inti]);
882 }
883 }
884
885 backwardInternal(wA, colouri);
886
887 // Store effect of exchanging rD to higher interfaces
888 // in colourBufs_
889 if (cellColourPtr_)
890 {
891 sendGlobal
892 (
893 lowerGlobalSend_[colouri],
894 wA,
895 lowerColour_[colouri]
896 );
897 }
898 }
899 }
900
901 // Make sure no outstanding sends
903
904 // Start writes of wA (using sendBufs)
906}
907
908
910(
911 const solverPerformance& s
912) const
913{
914 DebugPout<< "setFinished fieldName:" << s.fieldName() << endl;
915
916 // Wait on all requests before storage is being taken down
917 // (could rely on construction order?)
918 wait(lowerSendRequests_);
919 wait(lowerRecvRequests_);
920 wait(higherSendRequests_);
921 wait(higherRecvRequests_);
922
923 // TBD: cancel/ignore outstanding messages
924 //wait(lowerSendRequests_, true);
925 //wait(lowerRecvRequests_, true);
926 //wait(higherSendRequests_, true);
927 //wait(higherRecvRequests_, true);
928}
929
930
931// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void clear() noexcept
Clear the addressed list, i.e. set the size to zero.
T & emplace_back(Args &&... args)
Construct an element at the end of the list, return reference to the new list element.
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 resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition ListI.H:171
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
Definition PtrList.H:171
static std::streamsize read(const UPstream::commsTypes commsType, const int fromProcNo, Type *buffer, std::streamsize count, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm, UPstream::Request *req=nullptr)
Receive buffer contents (contiguous types) from given processor.
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition UListI.H:410
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static bool write(const UPstream::commsTypes commsType, const int toProcNo, const Type *buffer, std::streamsize count, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm, UPstream::Request *req=nullptr, const UPstream::sendModes sendMode=UPstream::sendModes::normal)
Write buffer contents (contiguous types only) to given processor.
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests).
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
Definition UPstream.H:84
static int & msgType() noexcept
Message tag of standard messages.
Definition UPstream.H:1926
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
static void cancelRequests(UList< UPstream::Request > &requests)
Non-blocking comms: cancel and free outstanding requests. Corresponds to MPI_Cancel() + MPI_Request_f...
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition UPstream.H:1948
static void waitRequests()
Wait for all requests to finish.
Definition UPstream.H:2497
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A cell is defined as a list of faces with extra functionality.
Definition cell.H:56
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Version of DILUpreconditioner that uses preconditioning across processor (and coupled) boundaries....
virtual void addInterface(solveScalarField &wA, const label inti, const Field< solveScalar > &recvBuf) const
Update preconditioned variable from interface.
virtual void backwardInternal(solveScalarField &wA, const label colouri) const
Update preconditioned variable walking backward on internal faces.
void send(const labelList &selectedInterfaces, const solveScalarField &psiInternal, DynamicList< UPstream::Request > &requests) const
Start sending sendBufs_.
List< DynamicList< label > > higherGlobalSend_
Interfaces to non-processor higher coupled interfaces.
PtrList< FieldField< Field, solveScalar > > colourBufs_
Global interfaces. Per colour the interfaces that (might) influence it.
static autoPtr< labelList > procColoursPtr_
Previous processor colours.
virtual void setFinished(const solverPerformance &perf) const
Signal end of solver.
List< DynamicList< label > > lowerGlobalSend_
Interfaces to non-processor lower coupled interfaces.
const bool coupled_
Precondition across global coupled bc.
void updateMatrixInterfaces(const bool add, const FieldField< Field, scalar > &coupleCoeffs, const labelList &selectedInterfaces, const solveScalarField &psiif, solveScalarField &result, const direction cmpt) const
Variant of lduMatrix::updateMatrixInterfaces on selected interfaces.
autoPtr< labelList > cellColourPtr_
Local (cell) colouring from global interfaces.
DynamicList< label > lowerNbrs_
Interfaces to lower coloured processors.
List< DynamicList< label > > lowerGlobalRecv_
Interfaces to non-processor lower coupled interfaces.
virtual void forwardInternalDiag(solveScalarField &rD, const label colouri) const
Update diagonal for all faces of a certain colour.
label nColours_
Number of colours (in case of multiple disconnected regions.
virtual void precondition(solveScalarField &wA, const solveScalarField &rA, const direction cmpt=0) const
Return wA the preconditioned form of residual rA.
DynamicList< UPstream::Request > higherSendRequests_
void sendGlobal(const labelList &selectedInterfaces, solveScalarField &psi, const label colouri) const
Send (and store in colourBufs_[colouri]) the effect of.
DynamicList< UPstream::Request > lowerSendRequests_
DynamicList< UPstream::Request > lowerRecvRequests_
List< label > lowerColour_
Corresponding destination colour (for lowerGlobal).
distributedDILUPreconditioner(const lduMatrix::solver &, const dictionary &solverControlsUnused)
Construct from matrix components and preconditioner solver controls.
void receive(const labelList &selectedInterfaces, DynamicList< UPstream::Request > &requests) const
Start receiving in recvBufs_.
virtual void addInterfaceDiag(solveScalarField &rD, const label inti, const Field< solveScalar > &recvBuf) const
Update diagonal for interface.
List< label > higherColour_
Corresponding destination colour (for higherGlobal).
DynamicList< label > higherNbrs_
Interfaces to higher coloured processors.
virtual void calcReciprocalD(solveScalarField &rD) const
Calculate reciprocal of diagonal.
static const lduMesh * meshPtr_
Processor interface buffers and colouring.
virtual void forwardInternal(solveScalarField &wA, const label colouri) const
Update preconditioned variable walking forward on internal faces.
FieldField< Field, solveScalar > sendBufs_
Buffers for sending and receiving data.
solveScalarField rD_
The reciprocal preconditioned diagonal.
DynamicList< UPstream::Request > higherRecvRequests_
List< DynamicList< label > > higherGlobalRecv_
Interfaces to non-processor higher coupled interfaces.
void wait(DynamicList< UPstream::Request > &requests, const bool cancel=false) const
Wait for requests or cancel/free requests.
Smooth ATC in cells next to a set of patches supplied by type.
Definition faceCells.H:55
A face is a list of labels corresponding to mesh vertices.
Definition face.H:71
A class for handling keywords in dictionaries.
Definition keyType.H:69
const solver & solver_
Reference to the base-solver this preconditioner is used with.
Definition lduMatrix.H:586
preconditioner(const solver &sol)
Construct for given solver.
Definition lduMatrix.H:634
Abstract base-class for lduMatrix solvers.
Definition lduMatrix.H:152
const lduInterfaceFieldPtrsList & interfaces() const noexcept
Definition lduMatrix.H:343
const lduMatrix & matrix() const noexcept
Definition lduMatrix.H:328
const FieldField< Field, scalar > & interfaceBouCoeffs() const noexcept
Definition lduMatrix.H:333
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
Definition lduMatrix.H:81
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition lduMatrix.H:769
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition lduMesh.H:54
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
Definition one.H:57
static label cellColour(const lduMesh &mesh, labelList &cellColour, labelList &patchToColour)
Calculate (locally) per cell colour according to walking from global patches. Returns number of colou...
static label colour(const lduMesh &mesh, labelList &procColour)
Calculate processor colouring from processor connectivity. Sets colour per processor such that no nei...
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
const volScalarField & psi
dynamicFvMesh & mesh
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))
#define WarningInFunction
Report a warning using Foam::Warning.
#define DebugPout
Report an information message using Foam::Pout.
Namespace for handling debugging switches.
Definition debug.C:45
Namespace for OpenFOAM.
List< label > labelList
A List of labels.
Definition List.H:62
messageStream Info
Information stream (stdout output on master, null elsewhere).
lduMatrix::preconditioner::addasymMatrixConstructorToTable< distributedDILUPreconditioner > adddistributedDILUPreconditionerAsymMatrixConstructorToTable_
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition Ostream.H:490
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
Ostream & indent(Ostream &os)
Indent stream.
Definition Ostream.H:481
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
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition FlatOutput.H:217
Field< solveScalar > solveScalarField
uint8_t direction
Definition direction.H:49
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
SolverPerformance< scalar > solverPerformance
SolverPerformance instantiated for a scalar.
dict add("bounds", meshBb)
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299