OpenLB 1.7
Loading...
Searching...
No Matches
Classes | Public Member Functions | List of all members
olb::ConcreteBlockCommunicator< ConcreteBlockLattice< T, DESCRIPTOR, Platform::GPU_CUDA > >::RecvTask Class Reference

Wrapper for a non-blocking block propagation receive request. More...

#include <communicator.hh>

+ Collaboration diagram for olb::ConcreteBlockCommunicator< ConcreteBlockLattice< T, DESCRIPTOR, Platform::GPU_CUDA > >::RecvTask:

Classes

class  ref
 Manual replacement for std::reference_wrapper<RecvTask> More...
 

Public Member Functions

 RecvTask (MPI_Comm comm, int tag, int rank, const std::vector< std::type_index > &fields, const std::vector< CellID > &cells, ConcreteBlockLattice< T, DESCRIPTOR, Platform::GPU_CUDA > &block)
 
 ~RecvTask ()
 
bool operator< (const RecvTask &rhs) const
 
void receive ()
 
bool isDone ()
 
void unpack ()
 
void wait ()
 

Detailed Description

template<typename T, typename DESCRIPTOR>
class olb::ConcreteBlockCommunicator< ConcreteBlockLattice< T, DESCRIPTOR, Platform::GPU_CUDA > >::RecvTask

Wrapper for a non-blocking block propagation receive request.

Definition at line 415 of file communicator.hh.

Constructor & Destructor Documentation

◆ RecvTask()

template<typename T , typename DESCRIPTOR >
olb::ConcreteBlockCommunicator< ConcreteBlockLattice< T, DESCRIPTOR, Platform::GPU_CUDA > >::RecvTask::RecvTask ( MPI_Comm comm,
int tag,
int rank,
const std::vector< std::type_index > & fields,
const std::vector< CellID > & cells,
ConcreteBlockLattice< T, DESCRIPTOR, Platform::GPU_CUDA > & block )
inline

Definition at line 458 of file communicator.hh.

461 :
462 _tag(tag),
463 _rank(rank),
464 _fields(block.getDataRegistry().deviceFieldArrays(fields)),
465 _onlyPopulationField(fields.size() == 1 && fields[0] == typeid(descriptors::POPULATION)),
466 _cells(cells),
467 _target(block),
468 _stream(std::make_unique<gpu::cuda::device::Stream>(cudaStreamNonBlocking))
469 {
470 std::size_t size = 0;
471 for (auto& field : fields) {
472 size += _target.getCommunicatable(field).size(cells);
473 }
474 _buffer = gpu::cuda::device::malloc<std::uint8_t>(size);
475 _request = std::make_unique<MpiRecvRequest>(
476 _buffer.get(), size, _rank, _tag, comm);
477 }
auto & getDataRegistry()
Return reference to Data's FieldTypeRegistry.
Communicatable & getCommunicatable(std::type_index field) override
virtual std::size_t size(ConstSpan< CellID > indices) const =0

References olb::gpu::cuda::device::unique_ptr< T >::get(), olb::ConcreteBlockLattice< T, DESCRIPTOR, PLATFORM >::getCommunicatable(), and olb::Communicatable::size().

+ Here is the call graph for this function:

◆ ~RecvTask()

template<typename T , typename DESCRIPTOR >
olb::ConcreteBlockCommunicator< ConcreteBlockLattice< T, DESCRIPTOR, Platform::GPU_CUDA > >::RecvTask::~RecvTask ( )
inline

Member Function Documentation

◆ isDone()

template<typename T , typename DESCRIPTOR >
bool olb::ConcreteBlockCommunicator< ConcreteBlockLattice< T, DESCRIPTOR, Platform::GPU_CUDA > >::RecvTask::isDone ( )
inline

Definition at line 495 of file communicator.hh.

496 {
497 return _request->isDone();
498 }

◆ operator<()

template<typename T , typename DESCRIPTOR >
bool olb::ConcreteBlockCommunicator< ConcreteBlockLattice< T, DESCRIPTOR, Platform::GPU_CUDA > >::RecvTask::operator< ( const RecvTask & rhs) const
inline

Definition at line 484 of file communicator.hh.

485 {
486 return _rank < rhs._rank
487 || (_rank == rhs._rank && _tag < rhs._tag);
488 }

◆ receive()

template<typename T , typename DESCRIPTOR >
void olb::ConcreteBlockCommunicator< ConcreteBlockLattice< T, DESCRIPTOR, Platform::GPU_CUDA > >::RecvTask::receive ( )
inline

Definition at line 490 of file communicator.hh.

491 {
492 _request->start();
493 };

◆ unpack()

template<typename T , typename DESCRIPTOR >
void olb::ConcreteBlockCommunicator< ConcreteBlockLattice< T, DESCRIPTOR, Platform::GPU_CUDA > >::RecvTask::unpack ( )
inline

Definition at line 500 of file communicator.hh.

501 {
502 if (_onlyPopulationField) {
503 gpu::cuda::DeviceContext<T,DESCRIPTOR> lattice(_target);
504 gpu::cuda::async_scatter_field<descriptors::POPULATION>(_stream->get(), lattice, _cells, _buffer.get());
505 } else {
506 gpu::cuda::async_scatter_any_fields(_stream->get(), _fields, _cells, _buffer.get());
507 }
508 }
void async_scatter_any_fields(cudaStream_t stream, thrust::device_vector< AnyDeviceFieldArrayD > &fields, const thrust::device_vector< CellID > &indices, std::uint8_t *buffer)
Non-blocking scatter of fields data in buffer to given indices.

References olb::gpu::cuda::async_scatter_any_fields(), and olb::gpu::cuda::device::unique_ptr< T >::get().

+ Here is the call graph for this function:

◆ wait()

template<typename T , typename DESCRIPTOR >
void olb::ConcreteBlockCommunicator< ConcreteBlockLattice< T, DESCRIPTOR, Platform::GPU_CUDA > >::RecvTask::wait ( )
inline

Definition at line 510 of file communicator.hh.

511 {
512 _stream->synchronize();
513 }

The documentation for this class was generated from the following file: