40 lduMatrix::preconditioner::
41 addasymMatrixConstructorToTable<distributedDILUPreconditioner>
63 const auto& lduAddr = matrix.
lduAddr();
73 for (
const label interfacei : selectedInterfaces)
75 interfaces[interfacei].initInterfaceMatrixUpdate
82 coupleCoeffs[interfacei],
91 for (
const label interfacei : selectedInterfaces)
93 interfaces[interfacei].updateInterfaceMatrix
100 coupleCoeffs[interfacei],
117 const auto& interfaces = solver_.interfaces();
119 if (selectedInterfaces.size())
122 FieldField<Field, scalar> one(interfaces.size());
123 FieldField<Field, solveScalar> old(interfaces.size());
124 for (
const label inti : selectedInterfaces)
126 const auto& intf = interfaces[inti].interface();
127 const auto& fc = intf.faceCells();
131 updateMatrixInterfaces
141 if (!colourBufs_.set(colouri))
152 auto& colourBuf = colourBufs_[colouri];
153 colourBuf.setSize(interfaces.size());
154 for (
const label inti : selectedInterfaces)
156 const auto& intf = interfaces[inti].interface();
157 const auto& fc = intf.faceCells();
158 if (!colourBuf.set(inti))
162 auto& cb = colourBuf[inti];
164 auto& oldValues = old[inti];
185 const auto& interfaces = solver_.interfaces();
186 const auto& interfaceBouCoeffs = solver_.interfaceBouCoeffs();
189 for (
const label inti : selectedInterfaces)
191 const auto& intf = interfaces[inti].interface();
194 auto& recvBuf = recvBufs_[inti];
195 recvBuf.resize_nocopy(interfaceBouCoeffs[inti].size());
199 requests.emplace_back(),
213 DynamicList<UPstream::Request>& requests
216 const auto& interfaces = solver_.interfaces();
219 for (
const label inti : selectedInterfaces)
221 const auto& intf = interfaces[inti].interface();
223 const auto& faceCells = intf.faceCells();
225 auto& sendBuf = sendBufs_[inti];
227 sendBuf.resize_nocopy(faceCells.size());
230 sendBuf[face] = psiInternal[faceCells[face]];
235 requests.emplace_back(),
247 DynamicList<UPstream::Request>& requests,
273 const auto& interfaces = solver_.interfaces();
274 const auto& interfaceBouCoeffs = solver_.interfaceBouCoeffs();
275 const auto& interfaceIntCoeffs = solver_.interfaceIntCoeffs();
277 const auto& intf = interfaces[inti].interface();
279 const auto&
faceCells = intf.faceCells();
280 const auto& bc = interfaceBouCoeffs[inti];
281 const auto& ic = interfaceIntCoeffs[inti];
300 const auto& matrix = solver_.matrix();
301 const auto& lduAddr = matrix.lduAddr();
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();
308 const label nFaces = matrix.upper().
size();
311 const auto& cellColour = *cellColourPtr_;
312 for (label face=0; face<nFaces; face++)
314 const label cell = lPtr[face];
315 if (cellColour[cell] == colouri)
317 rD[uPtr[face]] -= upperPtr[face]*lowerPtr[face]/rD[cell];
338 const auto& interfaces = solver_.interfaces();
339 const auto& interfaceBouCoeffs = solver_.interfaceBouCoeffs();
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();
362 const auto& matrix = solver_.matrix();
363 const auto& lduAddr = matrix.lduAddr();
365 solveScalar* __restrict__ wAPtr = wA.begin();
366 const solveScalar* __restrict__ rDPtr = rD_.begin();
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();
373 const label nFaces = matrix.upper().size();
376 const auto& cellColour = *cellColourPtr_;
377 for (label face=0; face<nFaces; face++)
379 const label sface = losortPtr[face];
380 const label cell = lPtr[sface];
381 if (cellColour[cell] == colouri)
383 wAPtr[uPtr[sface]] -=
384 rDPtr[uPtr[sface]]*lowerPtr[sface]*wAPtr[cell];
392 const label sface = losortPtr[
face];
393 wAPtr[uPtr[sface]] -=
394 rDPtr[uPtr[sface]]*lowerPtr[sface]*wAPtr[lPtr[sface]];
406 const auto& matrix = solver_.matrix();
407 const auto& lduAddr = matrix.lduAddr();
409 solveScalar* __restrict__ wAPtr = wA.begin();
410 const solveScalar* __restrict__ rDPtr = rD_.begin();
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();
416 const label nFacesM1 = matrix.upper().size() - 1;
420 const auto& cellColour = *cellColourPtr_;
424 if (cellColour[
cell] == colouri)
448 const auto& matrix = solver_.matrix();
451 wait(lowerRecvRequests_);
454 receive(lowerNbrs_, lowerRecvRequests_);
458 rD.resize_nocopy(
diag.size());
459 std::copy(
diag.begin(),
diag.end(), rD.begin());
465 wait(lowerRecvRequests_);
467 for (
const label inti : lowerNbrs_)
469 addInterfaceDiag(rD, inti, recvBufs_[inti]);
483 for (label colouri = 0; colouri < nColours_; colouri++)
487 for (
const label inti : lowerGlobalRecv_[colouri])
489 addInterfaceDiag(rD, inti, colourBufs_[colouri][inti]);
493 forwardInternalDiag(rD, colouri);
498 sendGlobal(higherGlobalSend_[colouri], rD, higherColour_[colouri]);
503 wait(higherSendRequests_);
506 send(higherNbrs_, rD, higherSendRequests_);
509 const label nCells = rD.size();
522 const lduMatrix::solver& sol,
526 lduMatrix::preconditioner(sol),
527 coupled_(
dict.getOrDefault<bool>(
"coupled", true, keyType::LITERAL)),
528 rD_(sol.matrix().
diag().size())
530 const lduMesh&
mesh = sol.matrix().
mesh();
531 const auto& interfaces = sol.interfaces();
532 const auto& interfaceBouCoeffs = sol.interfaceBouCoeffs();
539 forAll(interfaceBouCoeffs, inti)
541 if (interfaces.set(inti))
543 const auto& bc = interfaceBouCoeffs[inti];
558 const label nProcColours =
564 << nProcColours <<
endl;
573 bool haveGlobalCoupled =
false;
574 bool haveDistributedAMI =
false;
577 if (interfaces.set(inti))
579 const auto& intf = interfaces[inti].interface();
583 const label nbrColour = procColours[ppp->neighbProcNo()];
585 if (nbrColour < myColour)
589 else if (nbrColour > myColour)
596 <<
"weird processorLduInterface"
597 <<
" from " << ppp->myProcNo()
598 <<
" to " << ppp->neighbProcNo()
604 haveGlobalCoupled =
true;
620 <<
" of type " << intf.type()
621 <<
" is distributed:" << haveDistributedAMI <<
endl;
675 if (interfaces.set(inti))
677 const auto& intf = interfaces[inti].interface();
683 nbrInti = AMIpp->neighbPatchID();
688 nbrInti = cycpp->neighbPatchID();
693 const label colouri = patchToColour[inti];
694 const label nbrColouri = patchToColour[nbrInti];
698 (haveDistributedAMI && (inti < nbrInti))
699 || (!haveDistributedAMI && (colouri < nbrColouri))
729 if (lowerRecv.size())
734 <<
"Lower global interfaces for colour:" << colouri
735 <<
" receiving from:" <<
flatOutput(lowerRecv)
737 <<
" with colour:" << lowerCol <<
endl;
743 if (higherRecv.size())
748 <<
"Higher global interfaces for colour:" << colouri
749 <<
" receiving from:" <<
flatOutput(higherRecv)
751 <<
" with colour:" << higherCol <<
endl;
768 wait(lowerSendRequests_);
769 wait(lowerRecvRequests_);
770 wait(higherSendRequests_);
771 wait(higherRecvRequests_);
790 solveScalar* __restrict__ wAPtr = wA.
begin();
791 const solveScalar* __restrict__ rAPtr = rA.
begin();
792 const solveScalar* __restrict__ rDPtr = rD_.
begin();
794 const label nCells = wA.
size();
801 wait(lowerRecvRequests_);
804 receive(lowerNbrs_, lowerRecvRequests_);
815 wait(lowerRecvRequests_);
817 for (
const label inti : lowerNbrs_)
819 addInterface(wA, inti, recvBufs_[inti]);
822 for (label colouri = 0; colouri < nColours_; colouri++)
827 for (
const label inti : lowerGlobalRecv_[colouri])
829 addInterface(wA, inti, colourBufs_[colouri][inti]);
833 forwardInternal(wA, colouri);
841 higherGlobalSend_[colouri],
843 higherColour_[colouri]
850 wait(higherSendRequests_);
853 send(higherNbrs_, wA, higherSendRequests_);
860 wait(higherRecvRequests_);
863 receive(higherNbrs_, higherRecvRequests_);
867 wait(higherRecvRequests_);
869 for (
const label inti : higherNbrs_)
871 addInterface(wA, inti, recvBufs_[inti]);
874 for (label colouri = nColours_-1; colouri >= 0; colouri--)
879 for (
const label inti : higherGlobalRecv_[colouri])
881 addInterface(wA, inti, colourBufs_[colouri][inti]);
885 backwardInternal(wA, colouri);
893 lowerGlobalSend_[colouri],
895 lowerColour_[colouri]
918 wait(lowerSendRequests_);
919 wait(lowerRecvRequests_);
920 wait(higherSendRequests_);
921 wait(higherRecvRequests_);
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
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.
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
const T * set(const label i) const
Return const pointer to element (can be nullptr), or nullptr for out-of-range access (ie,...
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.
void size(const label n)
Older name for setAddressableSize.
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...
static label nRequests() noexcept
Number of outstanding requests (on the internal list of requests).
@ nonBlocking
"nonBlocking" (immediate) : (MPI_Isend, MPI_Irecv)
static int & msgType() noexcept
Message tag of standard messages.
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
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.
static void waitRequests()
Wait for all requests to finish.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
A cell is defined as a list of faces with extra functionality.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
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.
virtual ~distributedDILUPreconditioner()
Destructor.
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.
FieldField< Field, solveScalar > recvBufs_
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.
A face is a list of labels corresponding to mesh vertices.
A class for handling keywords in dictionaries.
const solver & solver_
Reference to the base-solver this preconditioner is used with.
preconditioner(const solver &sol)
Construct for given solver.
Abstract base-class for lduMatrix solvers.
const lduInterfaceFieldPtrsList & interfaces() const noexcept
const lduMatrix & matrix() const noexcept
const FieldField< Field, scalar > & interfaceBouCoeffs() const noexcept
lduMatrix is a general matrix class in which the coefficients are stored as three arrays,...
const lduAddressing & lduAddr() const
Return the LDU addressing.
Abstract base class for meshes which provide LDU addressing for the construction of lduMatrix and LDU...
A class representing the concept of 1 (one) that can be used to avoid manipulating objects known to b...
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.
const volScalarField & psi
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.
List< label > labelList
A List of labels.
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.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Ostream & indent(Ostream &os)
Indent stream.
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.
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Field< solveScalar > solveScalarField
static constexpr const zero Zero
Global zero (0).
void diag(pointPatchField< vector > &, const pointPatchField< tensor > &)
SolverPerformance< scalar > solverPerformance
SolverPerformance instantiated for a scalar.
dict add("bounds", meshBb)
#define forAll(list, i)
Loop across all elements in list.