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 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "processorLduInterface.H"
30 #include "cyclicLduInterface.H"
31 #include "cyclicAMILduInterface.H"
32 #include "processorColour.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
39 
40  lduMatrix::preconditioner::
41  addasymMatrixConstructorToTable<distributedDILUPreconditioner>
43 
45 
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,
146  new FieldField<Field, solveScalar>
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,
182  DynamicList<UPstream::Request>& requests
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<const processorLduInterface>(intf);
193 
194  auto& recvBuf = recvBufs_[inti];
195  recvBuf.resize_nocopy(interfaceBouCoeffs[inti].size());
196 
197  requests.push_back(UPstream::Request());
199  (
200  requests.back(),
201  ppp->neighbProcNo(),
202  recvBuf.data_bytes(),
203  recvBuf.size_bytes(),
204  ppp->tag()+70, // random offset
205  ppp->comm()
206  );
207  }
208 }
209 
210 
212 (
213  const labelList& selectedInterfaces,
214  const solveScalarField& psiInternal,
215  DynamicList<UPstream::Request>& requests
216 ) const
217 {
218  const auto& interfaces = solver_.interfaces();
219 
220  // Start writes
221  for (const label inti : selectedInterfaces)
222  {
223  const auto& intf = interfaces[inti].interface();
224  const auto* ppp = isA<const processorLduInterface>(intf);
225  const auto& faceCells = intf.faceCells();
226 
227  auto& sendBuf = sendBufs_[inti];
228 
229  sendBuf.resize_nocopy(faceCells.size());
230  forAll(faceCells, face)
231  {
232  sendBuf[face] = psiInternal[faceCells[face]];
233  }
234 
235  requests.push_back(UPstream::Request());
237  (
238  requests.back(),
239  ppp->neighbProcNo(),
240  sendBuf.cdata_bytes(),
241  sendBuf.size_bytes(),
242  ppp->tag()+70, // random offset
243  ppp->comm()
244  );
245  }
246 }
247 
248 
250 (
251  DynamicList<UPstream::Request>& requests,
252  const bool cancel
253 ) const
254 {
255  if (cancel)
256  {
257  UPstream::cancelRequests(requests);
258  }
259  else
260  {
261  UPstream::waitRequests(requests);
262  }
263  requests.clear();
264 }
265 
266 
267 // Diagonal calculation
268 // ~~~~~~~~~~~~~~~~~~~~
269 
271 (
272  solveScalarField& rD,
273  const label inti,
274  const Field<solveScalar>& recvBuf
275 ) const
276 {
277  const auto& interfaces = solver_.interfaces();
278  const auto& interfaceBouCoeffs = solver_.interfaceBouCoeffs();
279  const auto& interfaceIntCoeffs = solver_.interfaceIntCoeffs();
280 
281  const auto& intf = interfaces[inti].interface();
282  // TBD: do not use patch faceCells but passed-in addressing?
283  const auto& faceCells = intf.faceCells();
284  const auto& bc = interfaceBouCoeffs[inti];
285  const auto& ic = interfaceIntCoeffs[inti];
286 
287  forAll(recvBuf, face)
288  {
289  // Note:interfaceBouCoeffs, interfaceIntCoeffs on the receiving
290  // side are negated
291  rD[faceCells[face]] -= bc[face]*ic[face]/recvBuf[face];
292  }
293 }
294 
295 
297 (
298  solveScalarField& rD,
299  const label colouri
300 ) const
301 {
302  // Add forward constributions to diagonal
303 
304  const auto& matrix = solver_.matrix();
305  const auto& lduAddr = matrix.lduAddr();
306 
307  const label* const __restrict__ uPtr = lduAddr.upperAddr().begin();
308  const label* const __restrict__ lPtr = lduAddr.lowerAddr().begin();
309  const scalar* const __restrict__ upperPtr = matrix.upper().begin();
310  const scalar* const __restrict__ lowerPtr = matrix.lower().begin();
311 
312  const label nFaces = matrix.upper().size();
313  if (cellColourPtr_)
314  {
315  const auto& cellColour = *cellColourPtr_;
316  for (label face=0; face<nFaces; face++)
317  {
318  const label cell = lPtr[face];
319  if (cellColour[cell] == colouri)
320  {
321  rD[uPtr[face]] -= upperPtr[face]*lowerPtr[face]/rD[cell];
322  }
323  }
324  }
325  else
326  {
327  for (label face=0; face<nFaces; face++)
328  {
329  rD[uPtr[face]] -= upperPtr[face]*lowerPtr[face]/rD[lPtr[face]];
330  }
331  }
332 }
333 
334 
336 (
337  solveScalarField& wA,
338  const label inti,
339  const Field<solveScalar>& recvBuf
340 ) const
341 {
342  const auto& interfaces = solver_.interfaces();
343  const auto& interfaceBouCoeffs = solver_.interfaceBouCoeffs();
344 
345  const auto& intf = interfaces[inti].interface();
346  const auto& faceCells = intf.faceCells();
347  const auto& bc = interfaceBouCoeffs[inti];
348  const solveScalar* __restrict__ rDPtr = rD_.begin();
349  solveScalar* __restrict__ wAPtr = wA.begin();
350 
351  forAll(recvBuf, face)
352  {
353  // Note: interfaceBouCoeffs is -upperPtr
354  const label cell = faceCells[face];
355  wAPtr[cell] += rDPtr[cell]*bc[face]*recvBuf[face];
356  }
357 }
358 
359 
361 (
362  solveScalarField& wA,
363  const label colouri
364 ) const
365 {
366  const auto& matrix = solver_.matrix();
367  const auto& lduAddr = matrix.lduAddr();
368 
369  solveScalar* __restrict__ wAPtr = wA.begin();
370  const solveScalar* __restrict__ rDPtr = rD_.begin();
371 
372  const label* const __restrict__ uPtr = lduAddr.upperAddr().begin();
373  const label* const __restrict__ lPtr = lduAddr.lowerAddr().begin();
374  const label* const __restrict__ losortPtr = lduAddr.losortAddr().begin();
375  const scalar* const __restrict__ lowerPtr = matrix.lower().begin();
376 
377  const label nFaces = matrix.upper().size();
378  if (cellColourPtr_)
379  {
380  const auto& cellColour = *cellColourPtr_;
381  for (label face=0; face<nFaces; face++)
382  {
383  const label sface = losortPtr[face];
384  const label cell = lPtr[sface];
385  if (cellColour[cell] == colouri)
386  {
387  wAPtr[uPtr[sface]] -=
388  rDPtr[uPtr[sface]]*lowerPtr[sface]*wAPtr[cell];
389  }
390  }
391  }
392  else
393  {
394  for (label face=0; face<nFaces; face++)
395  {
396  const label sface = losortPtr[face];
397  wAPtr[uPtr[sface]] -=
398  rDPtr[uPtr[sface]]*lowerPtr[sface]*wAPtr[lPtr[sface]];
399  }
400  }
401 }
402 
403 
405 (
406  solveScalarField& wA,
407  const label colouri
408 ) const
409 {
410  const auto& matrix = solver_.matrix();
411  const auto& lduAddr = matrix.lduAddr();
412 
413  solveScalar* __restrict__ wAPtr = wA.begin();
414  const solveScalar* __restrict__ rDPtr = rD_.begin();
415 
416  const label* const __restrict__ uPtr = lduAddr.upperAddr().begin();
417  const label* const __restrict__ lPtr = lduAddr.lowerAddr().begin();
418  const scalar* const __restrict__ upperPtr = matrix.upper().begin();
419 
420  const label nFacesM1 = matrix.upper().size() - 1;
421 
422  if (cellColourPtr_)
423  {
424  const auto& cellColour = *cellColourPtr_;
425  for (label face=nFacesM1; face>=0; face--)
426  {
427  const label cell = uPtr[face];
428  if (cellColour[cell] == colouri)
429  {
430  // Note: lower cell guaranteed in same colour
431  wAPtr[lPtr[face]] -=
432  rDPtr[lPtr[face]]*upperPtr[face]*wAPtr[cell];
433  }
434  }
435  }
436  else
437  {
438  for (label face=nFacesM1; face>=0; face--)
439  {
440  wAPtr[lPtr[face]] -=
441  rDPtr[lPtr[face]]*upperPtr[face]*wAPtr[uPtr[face]];
442  }
443  }
444 }
445 
446 
448 (
449  solveScalarField& rD
450 ) const
451 {
452  const auto& matrix = solver_.matrix();
453 
454  // Make sure no outstanding receives
455  wait(lowerRecvRequests_);
456 
457  // Start reads (into recvBufs_)
458  receive(lowerNbrs_, lowerRecvRequests_);
459 
460  // Start with diagonal
461  const scalarField& diag = matrix.diag();
462  rD.resize_nocopy(diag.size());
463  std::copy(diag.begin(), diag.end(), rD.begin());
464 
465 
466  // Subtract coupled contributions
467  {
468  // Wait for finish. Received result in recvBufs
469  wait(lowerRecvRequests_);
470 
471  for (const label inti : lowerNbrs_)
472  {
473  addInterfaceDiag(rD, inti, recvBufs_[inti]);
474  }
475  }
476 
477 
478  // - divide the internal mesh into domains/colour, similar to domain
479  // decomposition
480  // - we could use subsetted bits of the mesh or just loop over only
481  // the cells of the current domain
482  // - a domain only uses boundary values of lower numbered domains
483  // - this also limits the interfaces that need to be evaluated since
484  // we assume an interface only changes local faceCells and these
485  // are all of the same colour
486 
487  for (label colouri = 0; colouri < nColours_; colouri++)
488  {
489  if (cellColourPtr_) // && colourBufs_.set(colouri))
490  {
491  for (const label inti : lowerGlobalRecv_[colouri])
492  {
493  addInterfaceDiag(rD, inti, colourBufs_[colouri][inti]);
494  }
495  }
496 
497  forwardInternalDiag(rD, colouri);
498 
499  // Store effect of exchanging rD to higher interfaces in colourBufs_
500  if (cellColourPtr_)
501  {
502  sendGlobal(higherGlobalSend_[colouri], rD, higherColour_[colouri]);
503  }
504  }
505 
506  // Make sure no outstanding sends
507  wait(higherSendRequests_);
508 
509  // Start writes of rD (using sendBufs)
510  send(higherNbrs_, rD, higherSendRequests_);
511 
512  // Calculate the reciprocal of the preconditioned diagonal
513  const label nCells = rD.size();
514 
515  for (label cell=0; cell<nCells; cell++)
516  {
517  rD[cell] = 1.0/rD[cell];
518  }
519 }
520 
521 
522 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
523 
525 (
526  const lduMatrix::solver& sol,
527  const dictionary& dict
528 )
529 :
530  lduMatrix::preconditioner(sol),
531  coupled_(dict.getOrDefault<bool>("coupled", true, keyType::LITERAL)),
532  rD_(sol.matrix().diag().size())
533 {
534  const lduMesh& mesh = sol.matrix().mesh();
535  const auto& interfaces = sol.interfaces();
536  const auto& interfaceBouCoeffs = sol.interfaceBouCoeffs();
537 
538  // Allocate buffers
539  // ~~~~~~~~~~~~~~~~
540 
541  sendBufs_.setSize(interfaces.size());
542  recvBufs_.setSize(interfaces.size());
543  forAll(interfaceBouCoeffs, inti)
544  {
545  if (interfaces.set(inti))
546  {
547  const auto& bc = interfaceBouCoeffs[inti];
548  sendBufs_.set(inti, new solveScalarField(bc.size(), Zero));
549  recvBufs_.set(inti, new solveScalarField(bc.size(), Zero));
550  }
551  }
552 
553 
554  // Determine processor colouring
555  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
556 
557  //const processorColour& colours = processorColour::New(mesh);
558  //const label myColour = colours[Pstream::myProcNo()];
559  if (&mesh != meshPtr_)
560  {
561  procColoursPtr_.reset(new labelList(Pstream::nProcs(mesh.comm())));
562  const label nProcColours =
564 
565  if (debug)
566  {
567  Info<< typeName << " : calculated processor colours:"
568  << nProcColours << endl;
569  }
570 
571  meshPtr_ = &mesh;
572  }
573  const labelList& procColours = procColoursPtr_();
574 
575  const label myColour = procColours[Pstream::myProcNo(mesh.comm())];
576 
577  bool haveGlobalCoupled = false;
578  bool haveDistributedAMI = false;
579  forAll(interfaces, inti)
580  {
581  if (interfaces.set(inti))
582  {
583  const auto& intf = interfaces[inti].interface();
584  const auto* ppp = isA<const processorLduInterface>(intf);
585  if (ppp)
586  {
587  const label nbrColour = procColours[ppp->neighbProcNo()];
588  //const label nbrColour = ppp->neighbProcNo();
589  if (nbrColour < myColour)
590  {
591  lowerNbrs_.append(inti);
592  }
593  else if (nbrColour > myColour)
594  {
595  higherNbrs_.append(inti);
596  }
597  else
598  {
600  << "weird processorLduInterface"
601  << " from " << ppp->myProcNo()
602  << " to " << ppp->neighbProcNo()
603  << endl;
604  }
605  }
606  else //if (isA<const cyclicAMILduInterface>(intf))
607  {
608  haveGlobalCoupled = true;
609 
610  const auto* AMIpp = isA<const cyclicAMILduInterface>(intf);
611  if (AMIpp)
612  {
613  //const auto& AMI =
614  //(
615  // AMIpp->owner()
616  // ? AMIpp->AMI()
617  // : AMIpp->neighbPatch().AMI()
618  //);
619  //haveDistributedAMI = AMI.distributed();
620 
621  if (debug)
622  {
623  Info<< typeName << " : interface " << inti
624  << " of type " << intf.type()
625  << " is distributed:" << haveDistributedAMI << endl;
626  }
627  }
628  }
629  }
630  }
631 
632 
633 
634  if (debug)
635  {
636  Info<< typeName
637  << " : lower proc interfaces:" << flatOutput(lowerNbrs_)
638  << " higher proc interfaces:" << flatOutput(higherNbrs_)
639  << endl;
640  }
641 
642 
643  // Determine local colouring/zoneing
644  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
645 
646  // Start off with all cells colour 0
647  nColours_ = 1;
648  cellColourPtr_.reset();
649  if (coupled_ && haveGlobalCoupled)
650  {
651  labelList patchToColour;
652  cellColourPtr_.reset(new labelList(mesh.lduAddr().size()));
653  //if (haveDistributedAMI)
654  //{
655  // nColours_ = 1;
656  // *cellColourPtr_ = 0;
657  // patchToColour = labelList(interfaces.size(), 0);
658  //}
659  //else
660  {
662  (
663  mesh,
664  cellColourPtr_(),
665  patchToColour
666  );
667  }
668 
675  colourBufs_.setSize(nColours_);
676 
677  forAll(interfaces, inti)
678  {
679  if (interfaces.set(inti))
680  {
681  const auto& intf = interfaces[inti].interface();
682 
683  label nbrInti = -1;
684  const auto* AMIpp = isA<const cyclicAMILduInterface>(intf);
685  if (AMIpp)
686  {
687  nbrInti = AMIpp->neighbPatchID();
688  }
689  const auto* cycpp = isA<const cyclicLduInterface>(intf);
690  if (cycpp)
691  {
692  nbrInti = cycpp->neighbPatchID();
693  }
694 
695  if (nbrInti != -1)
696  {
697  const label colouri = patchToColour[inti];
698  const label nbrColouri = patchToColour[nbrInti];
699 
700  if
701  (
702  (haveDistributedAMI && (inti < nbrInti))
703  || (!haveDistributedAMI && (colouri < nbrColouri))
704  )
705  {
706  // Send to higher
707  higherGlobalSend_[colouri].append(nbrInti);
708  higherColour_[colouri] = nbrColouri;
709  // Receive from higher
710  higherGlobalRecv_[colouri].append(inti);
711  }
712  else
713  {
714  // Send to lower
715  lowerGlobalSend_[colouri].append(nbrInti);
716  lowerColour_[colouri] = nbrColouri;
717  // Receive from lower
718  lowerGlobalRecv_[colouri].append(inti);
719  }
720  }
721  }
722  }
723  }
724 
725 
726  if (debug)
727  {
728  Info<< typeName << " : local colours:" << nColours_ << endl;
729  Info<< incrIndent;
730  forAll(lowerGlobalRecv_, colouri)
731  {
732  const auto& lowerRecv = lowerGlobalRecv_[colouri];
733  if (lowerRecv.size())
734  {
735  const auto& lowerSend = lowerGlobalSend_[colouri];
736  const auto& lowerCol = lowerColour_[colouri];
737  Info<< indent
738  << "Lower global interfaces for colour:" << colouri
739  << " receiving from:" << flatOutput(lowerRecv)
740  << " sending to:" << flatOutput(lowerSend)
741  << " with colour:" << lowerCol << endl;
742  }
743  }
744  forAll(higherGlobalRecv_, colouri)
745  {
746  const auto& higherRecv = higherGlobalRecv_[colouri];
747  if (higherRecv.size())
748  {
749  const auto& higherSend = higherGlobalSend_[colouri];
750  const auto& higherCol = higherColour_[colouri];
751  Info<< indent
752  << "Higher global interfaces for colour:" << colouri
753  << " receiving from:" << flatOutput(higherRecv)
754  << " sending to:" << flatOutput(higherSend)
755  << " with colour:" << higherCol << endl;
756  }
757  }
758  }
761 }
762 
763 
764 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
765 
767 {
768  DebugPout<< "~distributedDILUPreconditioner()" << endl;
769 
770  // Wait on all requests before storage is being taken down
771 
772  wait(lowerSendRequests_);
773  wait(lowerRecvRequests_);
774  wait(higherSendRequests_);
775  wait(higherRecvRequests_);
776 
777  // TBD: cancel/ignore outstanding messages
778  //wait(lowerSendRequests_, true);
779  //wait(lowerRecvRequests_, true);
780  //wait(higherSendRequests_, true);
781  //wait(higherRecvRequests_, true);
782 }
783 
784 
785 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
786 
788 (
789  solveScalarField& wA,
790  const solveScalarField& rA,
791  const direction
792 ) const
793 {
794  solveScalar* __restrict__ wAPtr = wA.begin();
795  const solveScalar* __restrict__ rAPtr = rA.begin();
796  const solveScalar* __restrict__ rDPtr = rD_.begin();
797 
798  const label nCells = wA.size();
799 
800 
801  // Forward sweep
802  // ~~~~~~~~~~~~~
803 
804  // Make sure no receives are still in flight (should not happen)
805  wait(lowerRecvRequests_);
806 
807  // Start reads (into recvBufs)
808  receive(lowerNbrs_, lowerRecvRequests_);
809 
810  // Initialise 'internal' cells
811  for (label cell=0; cell<nCells; cell++)
812  {
813  wAPtr[cell] = rDPtr[cell]*rAPtr[cell];
814  }
815 
816  // Do 'halo' contributions from lower numbered procs
817  {
818  // Wait for finish. Received result in recvBufs
819  wait(lowerRecvRequests_);
820 
821  for (const label inti : lowerNbrs_)
822  {
823  addInterface(wA, inti, recvBufs_[inti]);
824  }
825 
826  for (label colouri = 0; colouri < nColours_; colouri++)
827  {
828  // Do non-processor boundaries for this colour
829  if (cellColourPtr_) // && colourBufs_.set(colouri))
830  {
831  for (const label inti : lowerGlobalRecv_[colouri])
832  {
833  addInterface(wA, inti, colourBufs_[colouri][inti]);
834  }
835  }
836 
837  forwardInternal(wA, colouri);
838 
839  // Store effect of exchanging rD to higher interfaces
840  // in colourBufs_
841  if (cellColourPtr_)
842  {
843  sendGlobal
844  (
845  higherGlobalSend_[colouri],
846  wA,
847  higherColour_[colouri]
848  );
849  }
850  }
851  }
852 
853  // Make sure no outstanding sends from previous iteration
854  wait(higherSendRequests_);
855 
856  // Start writes of wA (using sendBufs)
857  send(higherNbrs_, wA, higherSendRequests_);
858 
859 
860  // Backward sweep
861  // ~~~~~~~~~~~~~~
862 
863  // Make sure no outstanding receives
864  wait(higherRecvRequests_);
865 
866  // Start receives
867  receive(higherNbrs_, higherRecvRequests_);
868 
869  {
870  // Wait until receives finished
871  wait(higherRecvRequests_);
872 
873  for (const label inti : higherNbrs_)
874  {
875  addInterface(wA, inti, recvBufs_[inti]);
876  }
877 
878  for (label colouri = nColours_-1; colouri >= 0; colouri--)
879  {
880  // Do non-processor boundaries for this colour
881  if (cellColourPtr_) // && colourBufs_.set(colouri))
882  {
883  for (const label inti : higherGlobalRecv_[colouri])
884  {
885  addInterface(wA, inti, colourBufs_[colouri][inti]);
886  }
887  }
888 
889  backwardInternal(wA, colouri);
890 
891  // Store effect of exchanging rD to higher interfaces
892  // in colourBufs_
893  if (cellColourPtr_)
894  {
895  sendGlobal
896  (
897  lowerGlobalSend_[colouri],
898  wA,
899  lowerColour_[colouri]
900  );
901  }
902  }
903  }
904 
905  // Make sure no outstanding sends
906  wait(lowerSendRequests_);
908  // Start writes of wA (using sendBufs)
909  send(lowerNbrs_, wA, lowerSendRequests_);
910 }
911 
912 
914 (
915  const solverPerformance& s
916 ) const
917 {
918  DebugPout<< "setFinished fieldName:" << s.fieldName() << endl;
919 
920  // Wait on all requests before storage is being taken down
921  // (could rely on construction order?)
922  wait(lowerSendRequests_);
923  wait(lowerRecvRequests_);
924  wait(higherSendRequests_);
925  wait(higherRecvRequests_);
926 
927  // TBD: cancel/ignore outstanding messages
928  //wait(lowerSendRequests_, true);
929  //wait(lowerRecvRequests_, true);
930  //wait(higherSendRequests_, true);
931  //wait(higherRecvRequests_, true);
932 }
933 
934 
935 // ************************************************************************* //
virtual void addInterfaceDiag(solveScalarField &rD, const label inti, const Field< solveScalar > &recvBuf) const
Update diagonal for interface.
FieldField< Field, solveScalar > sendBufs_
Buffers for sending and receiving data.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
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...
virtual void precondition(solveScalarField &wA, const solveScalarField &rA, const direction cmpt=0) const
Return wA the preconditioned form of residual rA.
uint8_t direction
Definition: direction.H:46
Ostream & indent(Ostream &os)
Indent stream.
Definition: Ostream.H:493
void send(const labelList &selectedInterfaces, const solveScalarField &psiInternal, DynamicList< UPstream::Request > &requests) const
Start sending sendBufs_.
static label read(const UPstream::commsTypes commsType, const int fromProcNo, char *buf, const std::streamsize bufSize, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm, UPstream::Request *req=nullptr)
Read buffer contents from given processor.
Definition: UIPstreamRead.C:35
static int incrMsgType(int val=1) noexcept
Increment the message tag for standard messages.
Definition: UPstream.H:1251
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
Field< solveScalar > solveScalarField
distributedDILUPreconditioner(const lduMatrix::solver &, const dictionary &solverControlsUnused)
Construct from matrix components and preconditioner solver controls.
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:517
const solver & solver_
Reference to the base-solver this preconditioner is used with.
Definition: lduMatrix.H:558
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests)
autoPtr< labelList > cellColourPtr_
Local (cell) colouring from global interfaces.
const lduMatrix & matrix() const noexcept
Definition: lduMatrix.H:306
virtual void forwardInternalDiag(solveScalarField &rD, const label colouri) const
Update diagonal for all faces of a certain colour.
solveScalarField rD_
The reciprocal preconditioned diagonal.
List< DynamicList< label > > lowerGlobalSend_
Interfaces to non-processor lower coupled interfaces.
virtual void setFinished(const solverPerformance &perf) const
Signal end of solver.
static const lduMesh * meshPtr_
Processor interface buffers and colouring.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
Definition: lduMesh.H:53
List< DynamicList< label > > higherGlobalSend_
Interfaces to non-processor higher coupled interfaces.
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition: ListI.H:175
static int & msgType() noexcept
Message tag of standard messages.
Definition: UPstream.H:1229
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1074
Smooth ATC in cells next to a set of patches supplied by type.
Definition: faceCells.H:52
FieldField< Field, solveScalar > recvBufs_
static void waitRequests()
Wait for all requests to finish.
Definition: UPstream.H:1538
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
List< label > higherColour_
Corresponding destination colour (for higherGlobal)
static autoPtr< labelList > procColoursPtr_
Previous processor colours.
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
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:1065
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
virtual void forwardInternal(solveScalarField &wA, const label colouri) const
Update preconditioned variable walking forward on internal faces.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
SolverPerformance is the class returned by the LduMatrix solver containing performance statistics...
virtual void addInterface(solveScalarField &wA, const label inti, const Field< solveScalar > &recvBuf) const
Update preconditioned variable from interface.
const lduInterfaceFieldPtrsList & interfaces() const noexcept
Definition: lduMatrix.H:321
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicList.H:584
Version of DILUpreconditioner that uses preconditioning across processor (and coupled) boundaries...
virtual void backwardInternal(solveScalarField &wA, const label colouri) const
Update preconditioned variable walking backward on internal faces.
iterator begin() noexcept
Return an iterator to begin traversing the UList.
Definition: UListI.H:391
virtual void calcReciprocalD(solveScalarField &rD) const
Calculate reciprocal of diagonal.
void receive(const labelList &selectedInterfaces, DynamicList< UPstream::Request > &requests) const
Start receiving in recvBufs_.
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.
void add(FieldField< Field1, typename typeOfSum< Type1, Type2 >::type > &f, const FieldField< Field1, Type1 > &f1, const FieldField< Field2, Type2 > &f2)
DynamicList< label > lowerNbrs_
Interfaces to lower coloured processors.
int debug
Static debugging option.
PtrList< FieldField< Field, solveScalar > > colourBufs_
Global interfaces. Per colour the interfaces that (might) influence it.
defineTypeNameAndDebug(combustionModel, 0)
static label colour(const lduMesh &mesh, labelList &procColour)
Calculate processor colouring from processor connectivity. Sets colour per processor such that no nei...
List< DynamicList< label > > lowerGlobalRecv_
Interfaces to non-processor lower coupled interfaces.
static void cancelRequests(UList< UPstream::Request > &requests)
Non-blocking comms: cancel and free outstanding requests. Corresponds to MPI_Cancel() + MPI_Request_f...
void sendGlobal(const labelList &selectedInterfaces, solveScalarField &psi, const label colouri) const
Send (and store in colourBufs_[colouri]) the effect of.
DynamicList< label > higherNbrs_
Interfaces to higher coloured processors.
label nColours_
Number of colours (in case of multiple disconnected regions.
#define WarningInFunction
Report a warning using Foam::Warning.
const lduAddressing & lduAddr() const
Return the LDU addressing.
Definition: lduMatrix.H:734
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:53
static bool write(const UPstream::commsTypes commsType, const int toProcNo, const char *buf, const std::streamsize bufSize, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm, UPstream::Request *req=nullptr, const UPstream::sendModes sendMode=UPstream::sendModes::normal)
Write buffer contents to given processor.
#define DebugPout
Report an information message using Foam::Pout.
"nonBlocking" : (MPI_Isend, MPI_Irecv)
messageStream Info
Information stream (stdout output on master, null elsewhere)
void wait(DynamicList< UPstream::Request > &requests, const bool cancel=false) const
Wait for requests or cancel/free requests.
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
const volScalarField & psi
List< label > labelList
A List of labels.
Definition: List.H:62
List< label > lowerColour_
Corresponding destination colour (for lowerGlobal)
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))
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:502
List< DynamicList< label > > higherGlobalRecv_
Interfaces to non-processor higher coupled interfaces.
lduMatrix::preconditioner::addasymMatrixConstructorToTable< distributedDILUPreconditioner > adddistributedDILUPreconditionerAsymMatrixConstructorToTable_
Namespace for OpenFOAM.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:225
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127