Arcane  4.1.12.0
Developer documentation
Loading...
Searching...
No Matches
MpiAdapter.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2026 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
4// See the top-level COPYRIGHT file for details.
5// SPDX-License-Identifier: Apache-2.0
6//-----------------------------------------------------------------------------
7/*---------------------------------------------------------------------------*/
8/* MpiAdapter.cc (C) 2000-2026 */
9/* */
10/* Parallelism manager using MPI. */
11/*---------------------------------------------------------------------------*/
12/*---------------------------------------------------------------------------*/
13
14#include "arccore/message_passing_mpi/internal/MpiAdapter.h"
15
16#include "arccore/trace/ITraceMng.h"
17
18#include "arccore/collections/Array.h"
19
20#include "arccore/message_passing/Request.h"
21#include "arccore/message_passing/IStat.h"
22#include "arccore/message_passing/internal/SubRequestCompletionInfo.h"
23
24#include "arccore/base/IStackTraceService.h"
25#include "arccore/base/TimeoutException.h"
26#include "arccore/base/String.h"
27#include "arccore/base/NotImplementedException.h"
28#include "arccore/base/PlatformUtils.h"
29#include "arccore/base/FatalErrorException.h"
30#include "arccore/base/TraceInfo.h"
31
32#include "arccore/message_passing_mpi/StandaloneMpiMessagePassingMng.h"
33#include "arccore/message_passing_mpi/internal/MpiLock.h"
34#include "arccore/message_passing_mpi/internal/NoMpiProfiling.h"
35#include "arccore/message_passing_mpi/internal/MpiRequest.h"
36#include "arccore/message_passing_mpi/internal/MpiMachineShMemWinBaseInternalCreator.h"
37
38#include <cstdint>
39
40/*---------------------------------------------------------------------------*/
41/*---------------------------------------------------------------------------*/
42
43namespace Arcane::MessagePassing::Mpi
44{
45
46/*---------------------------------------------------------------------------*/
47/*---------------------------------------------------------------------------*/
48
50: public TraceAccessor
51{
52 public:
53
55 {
56 TraceInfo m_trace;
57 String m_stack_trace;
58 };
59
60 public:
61
62 typedef std::map<MPI_Request, RequestInfo>::iterator Iterator;
63
64 public:
65
66 explicit RequestSet(ITraceMng* tm)
67 : TraceAccessor(tm)
68 {
69 m_trace_mng_ref = makeRef(tm);
70 if (arccoreIsCheck()) {
71 m_no_check_request = false;
72 m_request_error_is_fatal = true;
73 }
74 if (Platform::getEnvironmentVariable("ARCCORE_NOREPORT_ERROR_MPIREQUEST") == "TRUE")
75 m_is_report_error_in_request = false;
76 if (Platform::getEnvironmentVariable("ARCCORE_MPIREQUEST_STACKTRACE") == "TRUE")
77 m_use_trace_full_stack = true;
78 if (Platform::getEnvironmentVariable("ARCCORE_TRACE_MPIREQUEST") == "TRUE")
79 m_trace_mpirequest = true;
80 }
81
82 public:
83
84 void addRequest(MPI_Request request)
85 {
87 return;
88 if (m_trace_mpirequest)
89 info() << "MpiAdapter: AddRequest r=" << request;
90 _addRequest(request, TraceInfo());
91 }
92 void addRequest(MPI_Request request, const TraceInfo& ti)
93 {
95 return;
96 if (m_trace_mpirequest)
97 info() << "MpiAdapter: AddRequest r=" << request;
98 _addRequest(request, ti);
99 }
100 void removeRequest(MPI_Request request)
101 {
103 return;
104 if (m_trace_mpirequest)
105 info() << "MpiAdapter: RemoveRequest r=" << request;
106 _removeRequest(request);
107 }
108 void removeRequest(Iterator request_iter)
109 {
111 return;
112 if (request_iter == m_allocated_requests.end()) {
113 if (m_trace_mpirequest)
114 info() << "MpiAdapter: RemoveRequestIter null iterator";
115 return;
116 }
117 if (m_trace_mpirequest)
118 info() << "MpiAdapter: RemoveRequestIter r=" << request_iter->first;
119 m_allocated_requests.erase(request_iter);
120 }
122 Iterator findRequest(MPI_Request request)
123 {
125 return m_allocated_requests.end();
126
127 if (_isEmptyRequest(request))
128 return m_allocated_requests.end();
129 auto ireq = m_allocated_requests.find(request);
130 if (ireq == m_allocated_requests.end()) {
131 if (m_is_report_error_in_request || m_request_error_is_fatal) {
132 error() << "MpiAdapter::testRequest() request not referenced "
133 << " id=" << request;
134 _checkFatalInRequest();
135 }
136 }
137 return ireq;
138 }
139
140 private:
141
145 void _addRequest(MPI_Request request, const TraceInfo& trace_info)
146 {
147 if (request == MPI_REQUEST_NULL) {
148 if (m_is_report_error_in_request || m_request_error_is_fatal) {
149 error() << "MpiAdapter::_addRequest() trying to add null request";
150 _checkFatalInRequest();
151 }
152 return;
153 }
154 if (_isEmptyRequest(request))
155 return;
156 ++m_total_added_request;
157 //info() << "MPI_ADAPTER:ADD REQUEST " << request;
158 auto i = m_allocated_requests.find(request);
159 if (i != m_allocated_requests.end()) {
160 if (m_is_report_error_in_request || m_request_error_is_fatal) {
161 error() << "MpiAdapter::_addRequest() request already referenced "
162 << " id=" << request;
163 _checkFatalInRequest();
164 }
165 return;
166 }
167 RequestInfo rinfo;
168 rinfo.m_trace = trace_info;
169 if (m_use_trace_full_stack)
170 rinfo.m_stack_trace = Platform::getStackTrace();
171 m_allocated_requests.insert(std::make_pair(request, rinfo));
172 }
173
177 void _removeRequest(MPI_Request request)
178 {
179 //info() << "MPI_ADAPTER:REMOVE REQUEST " << request;
180 if (request == MPI_REQUEST_NULL) {
181 if (m_is_report_error_in_request || m_request_error_is_fatal) {
182 error() << "MpiAdapter::_removeRequest() null request (" << MPI_REQUEST_NULL << ")";
183 _checkFatalInRequest();
184 }
185 return;
186 }
187 if (_isEmptyRequest(request))
188 return;
189 auto i = m_allocated_requests.find(request);
190 if (i == m_allocated_requests.end()) {
191 if (m_is_report_error_in_request || m_request_error_is_fatal) {
192 error() << "MpiAdapter::_removeRequest() request not referenced "
193 << " id=" << request;
194 _checkFatalInRequest();
195 }
196 }
197 else
198 m_allocated_requests.erase(i);
199 }
200
201 public:
202
203 void _checkFatalInRequest()
204 {
205 if (m_request_error_is_fatal)
206 ARCCORE_FATAL("Error in requests management");
207 }
208 Int64 nbRequest() const { return m_allocated_requests.size(); }
209 Int64 totalAddedRequest() const { return m_total_added_request; }
210 void printRequests() const
211 {
212 info() << "PRINT REQUESTS\n";
213 for (auto& x : m_allocated_requests) {
214 info() << "Request id=" << x.first << " trace=" << x.second.m_trace
215 << " stack=" << x.second.m_stack_trace;
216 }
217 }
218 void setEmptyRequests(MPI_Request r1, MPI_Request r2)
219 {
220 m_empty_request1 = r1;
221 m_empty_request2 = r2;
222 }
223
224 public:
225
226 bool m_request_error_is_fatal = false;
227 bool m_is_report_error_in_request = true;
228 bool m_trace_mpirequest = false;
231
232 private:
233
234 std::map<MPI_Request, RequestInfo> m_allocated_requests;
235 bool m_use_trace_full_stack = false;
236 MPI_Request m_empty_request1 = MPI_REQUEST_NULL;
237 MPI_Request m_empty_request2 = MPI_REQUEST_NULL;
238 Int64 m_total_added_request = 0;
239 Ref<ITraceMng> m_trace_mng_ref;
240
241 private:
242
243 bool _isEmptyRequest(MPI_Request r) const
244 {
245 return (r == m_empty_request1 || r == m_empty_request2);
246 }
247};
248
249#define ARCCORE_ADD_REQUEST(request) \
250 m_request_set->addRequest(request, A_FUNCINFO);
251
252/*---------------------------------------------------------------------------*/
253/*---------------------------------------------------------------------------*/
254
255namespace
256{
257 int _checkSize(Int64 i64_size)
258 {
259 if (i64_size > INT32_MAX)
260 ARCCORE_FATAL("Can not convert '{0}' to type integer", i64_size);
261 return (int)i64_size;
262 }
263} // namespace
264
265/*---------------------------------------------------------------------------*/
266/*---------------------------------------------------------------------------*/
267
269MpiAdapter(ITraceMng* trace, IStat* stat, MPI_Comm comm,
270 MpiLock* mpi_lock, IMpiProfiling* mpi_op)
271: TraceAccessor(trace)
272, m_stat(stat)
273, m_mpi_lock(mpi_lock)
274, m_mpi_prof(mpi_op)
275, m_communicator(comm)
276, m_comm_rank(0)
277, m_comm_size(0)
278, m_empty_request1(MPI_REQUEST_NULL)
279, m_empty_request2(MPI_REQUEST_NULL)
280{
281 m_request_set = new RequestSet(trace);
282
283 if (Platform::getEnvironmentVariable("ARCCORE_TRACE_MPI") == "TRUE")
284 m_is_trace = true;
285 {
286 String s = Platform::getEnvironmentVariable("ARCCORE_ALLOW_NULL_RANK_FOR_MPI_ANY_SOURCE");
287 if (s == "1" || s == "TRUE")
288 m_is_allow_null_rank_for_any_source = true;
289 if (s == "0" || s == "FALSE")
290 m_is_allow_null_rank_for_any_source = false;
291 }
292
293 ::MPI_Comm_rank(m_communicator, &m_comm_rank);
294 ::MPI_Comm_size(m_communicator, &m_comm_size);
295
296 // By default, we do not do MPI profiling; we will use the appropriate set
297 // method to change it
298 if (!m_mpi_prof)
299 m_mpi_prof = new NoMpiProfiling();
300
311 MPI_Irecv(m_recv_buffer_for_empty_request, 1, MPI_CHAR, MPI_PROC_NULL,
313
314 /*
315 * Starting from version 4 of openmpi, it also seems that sends with small
316 * buffers always generate the same request. Therefore, it must also be
317 * removed from the requests to test. We also post the corresponding
318 * MPI_Recv to prevent MPI_ISend from being unintentionally used in a
319 * user MPI_Recv (e.g., via MPI_Recv(MPI_ANY_TAG)).
320 */
321 m_send_buffer_for_empty_request2[0] = 0;
322 MPI_Isend(m_send_buffer_for_empty_request2, 1, MPI_CHAR, m_comm_rank,
323 50505, m_communicator, &m_empty_request2);
324
325 MPI_Recv(m_recv_buffer_for_empty_request2, 1, MPI_CHAR, m_comm_rank,
326 50505, m_communicator, MPI_STATUS_IGNORE);
327
328 m_request_set->setEmptyRequests(m_empty_request1, m_empty_request2);
329}
330
331/*---------------------------------------------------------------------------*/
332/*---------------------------------------------------------------------------*/
333
334MpiAdapter::
335~MpiAdapter()
336{
337 if (m_empty_request1 != MPI_REQUEST_NULL)
338 MPI_Request_free(&m_empty_request1);
339 if (m_empty_request2 != MPI_REQUEST_NULL)
340 MPI_Request_free(&m_empty_request2);
341
342 delete m_request_set;
343 delete m_mpi_prof;
344}
345
346/*---------------------------------------------------------------------------*/
347/*---------------------------------------------------------------------------*/
348
350buildRequest(int ret, MPI_Request mpi_request)
351{
352 return MpiRequest(ret, this, mpi_request);
353}
354
355/*---------------------------------------------------------------------------*/
356/*---------------------------------------------------------------------------*/
357
358void MpiAdapter::
359_checkHasNoRequests()
360{
361 Int64 nb_request = m_request_set->nbRequest();
362 // We cannot perform this test in the destructor because it could
363 // potentially throw an exception, and this should not be done in a destructor.
364 if (nb_request != 0) {
365 warning() << " Pending mpi requests size=" << nb_request;
366 m_request_set->printRequests();
367 _checkFatalInRequest();
368 }
369}
370
371/*---------------------------------------------------------------------------*/
372/*---------------------------------------------------------------------------*/
373
375destroy()
376{
377 _checkHasNoRequests();
378 delete this;
379}
380
381/*---------------------------------------------------------------------------*/
382/*---------------------------------------------------------------------------*/
383
386{
387 m_request_set->m_request_error_is_fatal = v;
388}
389bool MpiAdapter::
390isRequestErrorAreFatal() const
391{
392 return m_request_set->m_request_error_is_fatal;
393}
394
397{
398 m_request_set->m_is_report_error_in_request = v;
399}
400bool MpiAdapter::
401isPrintRequestError() const
402{
403 return m_request_set->m_is_report_error_in_request;
404}
405
407setCheckRequest(bool v)
408{
409 m_request_set->m_no_check_request = !v;
410}
411
412bool MpiAdapter::
413isCheckRequest() const
414{
415 return !m_request_set->m_no_check_request;
416}
417
418/*---------------------------------------------------------------------------*/
419/*---------------------------------------------------------------------------*/
420
421int MpiAdapter::
422toMPISize(Int64 count)
423{
424 return _checkSize(count);
425}
426
427/*---------------------------------------------------------------------------*/
428/*---------------------------------------------------------------------------*/
429
430void MpiAdapter::
431_trace(const char* function)
432{
433 if (m_is_trace) {
435 if (stack_service)
436 info() << "MPI_TRACE: " << function << "\n"
437 << stack_service->stackTrace().toString();
438 else
439 info() << "MPI_TRACE: " << function;
440 }
441}
442
443/*---------------------------------------------------------------------------*/
444/*---------------------------------------------------------------------------*/
445
446void MpiAdapter::
447broadcast(void* buf, Int64 nb_elem, Int32 root, MPI_Datatype datatype)
448{
449 int _nb_elem = _checkSize(nb_elem);
450 _trace(MpiInfo(eMpiName::Bcast).name().localstr());
451 double begin_time = MPI_Wtime();
452 if (m_is_trace)
453 info() << "MPI_TRACE: MPI broadcast: before"
454 << " buf=" << buf
455 << " nb_elem=" << nb_elem
456 << " root=" << root
457 << " datatype=" << datatype;
458
459 m_mpi_prof->broadcast(buf, _nb_elem, datatype, root, m_communicator);
460 double end_time = MPI_Wtime();
461 double sr_time = (end_time - begin_time);
462 //TODO determine the message size
463 m_stat->add(MpiInfo(eMpiName::Bcast).name(), sr_time, 0);
464}
465
466/*---------------------------------------------------------------------------*/
467/*---------------------------------------------------------------------------*/
468
469Request MpiAdapter::
470nonBlockingBroadcast(void* buf, Int64 nb_elem, Int32 root, MPI_Datatype datatype)
471{
472 MPI_Request mpi_request = MPI_REQUEST_NULL;
473 int ret = -1;
474 int _nb_elem = _checkSize(nb_elem);
475 _trace(" MPI_Bcast");
476 double begin_time = MPI_Wtime();
477 ret = MPI_Ibcast(buf, _nb_elem, datatype, root, m_communicator, &mpi_request);
478 double end_time = MPI_Wtime();
479 double sr_time = (end_time - begin_time);
480 //TODO determine the message size
481 m_stat->add("IBroadcast", sr_time, 0);
482 ARCCORE_ADD_REQUEST(mpi_request);
483 return buildRequest(ret, mpi_request);
484}
485
486/*---------------------------------------------------------------------------*/
487/*---------------------------------------------------------------------------*/
488
489void MpiAdapter::
490gather(const void* send_buf, void* recv_buf, Int64 nb_elem, Int32 root, MPI_Datatype datatype)
491{
492 void* _sbuf = const_cast<void*>(send_buf);
493 int _nb_elem = _checkSize(nb_elem);
494 int _root = static_cast<int>(root);
495 _trace(MpiInfo(eMpiName::Gather).name().localstr());
496 double begin_time = MPI_Wtime();
497 m_mpi_prof->gather(_sbuf, _nb_elem, datatype, recv_buf, _nb_elem, datatype, _root, m_communicator);
498 double end_time = MPI_Wtime();
499 double sr_time = (end_time - begin_time);
500 //TODO determine the message size
501 m_stat->add(MpiInfo(eMpiName::Gather).name(), sr_time, 0);
502}
503
504/*---------------------------------------------------------------------------*/
505/*---------------------------------------------------------------------------*/
506
507Request MpiAdapter::
508nonBlockingGather(const void* send_buf, void* recv_buf,
509 Int64 nb_elem, Int32 root, MPI_Datatype datatype)
510{
511 MPI_Request mpi_request = MPI_REQUEST_NULL;
512 int ret = -1;
513 void* _sbuf = const_cast<void*>(send_buf);
514 int _nb_elem = _checkSize(nb_elem);
515 int _root = static_cast<int>(root);
516 _trace("MPI_Igather");
517 double begin_time = MPI_Wtime();
518 ret = MPI_Igather(_sbuf, _nb_elem, datatype, recv_buf, _nb_elem, datatype, _root,
519 m_communicator, &mpi_request);
520 double end_time = MPI_Wtime();
521 double sr_time = (end_time - begin_time);
522 //TODO determine the message size
523 m_stat->add("IGather", sr_time, 0);
524 ARCCORE_ADD_REQUEST(mpi_request);
525 return buildRequest(ret, mpi_request);
526}
527
528/*---------------------------------------------------------------------------*/
529/*---------------------------------------------------------------------------*/
530
531void MpiAdapter::
532allGather(const void* send_buf, void* recv_buf,
533 Int64 nb_elem, MPI_Datatype datatype)
534{
535 void* _sbuf = const_cast<void*>(send_buf);
536 int _nb_elem = _checkSize(nb_elem);
537 _trace(MpiInfo(eMpiName::Allgather).name().localstr());
538 double begin_time = MPI_Wtime();
539 m_mpi_prof->allGather(_sbuf, _nb_elem, datatype, recv_buf, _nb_elem, datatype, m_communicator);
540 double end_time = MPI_Wtime();
541 double sr_time = (end_time - begin_time);
542 //TODO determine the message size
543 m_stat->add(MpiInfo(eMpiName::Allgather).name(), sr_time, 0);
544}
545
546/*---------------------------------------------------------------------------*/
547/*---------------------------------------------------------------------------*/
548
549Request MpiAdapter::
550nonBlockingAllGather(const void* send_buf, void* recv_buf,
551 Int64 nb_elem, MPI_Datatype datatype)
552{
553 MPI_Request mpi_request = MPI_REQUEST_NULL;
554 int ret = -1;
555 void* _sbuf = const_cast<void*>(send_buf);
556 int _nb_elem = _checkSize(nb_elem);
557 _trace("MPI_Iallgather");
558 double begin_time = MPI_Wtime();
559 ret = MPI_Iallgather(_sbuf, _nb_elem, datatype, recv_buf, _nb_elem, datatype,
560 m_communicator, &mpi_request);
561 double end_time = MPI_Wtime();
562 double sr_time = (end_time - begin_time);
563 //TODO determine the message size
564 m_stat->add("IAllGather", sr_time, 0);
565 ARCCORE_ADD_REQUEST(mpi_request);
566 return buildRequest(ret, mpi_request);
567}
568
569/*---------------------------------------------------------------------------*/
570/*---------------------------------------------------------------------------*/
571
572void MpiAdapter::
573gatherVariable(const void* send_buf, void* recv_buf, const int* recv_counts,
574 const int* recv_indexes, Int64 nb_elem, Int32 root, MPI_Datatype datatype)
575{
576 void* _sbuf = const_cast<void*>(send_buf);
577 int _nb_elem = _checkSize(nb_elem);
578 int _root = static_cast<int>(root);
579 _trace(MpiInfo(eMpiName::Gatherv).name().localstr());
580 double begin_time = MPI_Wtime();
581 m_mpi_prof->gatherVariable(_sbuf, _nb_elem, datatype, recv_buf, recv_counts, recv_indexes, datatype, _root, m_communicator);
582 double end_time = MPI_Wtime();
583 double sr_time = (end_time - begin_time);
584 //TODO determine the message size
585 m_stat->add(MpiInfo(eMpiName::Gatherv).name().localstr(), sr_time, 0);
586}
587
588/*---------------------------------------------------------------------------*/
589/*---------------------------------------------------------------------------*/
590
591void MpiAdapter::
592allGatherVariable(const void* send_buf, void* recv_buf, const int* recv_counts,
593 const int* recv_indexes, Int64 nb_elem, MPI_Datatype datatype)
594{
595 void* _sbuf = const_cast<void*>(send_buf);
596 int _nb_elem = _checkSize(nb_elem);
597 _trace(MpiInfo(eMpiName::Allgatherv).name().localstr());
598 //info() << " ALLGATHERV N=" << _nb_elem;
599 //for( int i=0; i<m_comm_size; ++i )
600 //info() << " ALLGATHERV I=" << i << " recv_count=" << recv_counts[i]
601 // << " recv_indexes=" << recv_indexes[i];
602 double begin_time = MPI_Wtime();
603 m_mpi_prof->allGatherVariable(_sbuf, _nb_elem, datatype, recv_buf, recv_counts, recv_indexes, datatype, m_communicator);
604 double end_time = MPI_Wtime();
605 double sr_time = (end_time - begin_time);
606 //TODO determine the message size
607 m_stat->add(MpiInfo(eMpiName::Allgatherv).name().localstr(), sr_time, 0);
608}
609
610/*---------------------------------------------------------------------------*/
611/*---------------------------------------------------------------------------*/
612
613void MpiAdapter::
614scatterVariable(const void* send_buf, const int* send_count, const int* send_indexes,
615 void* recv_buf, Int64 nb_elem, Int32 root, MPI_Datatype datatype)
616{
617 void* _sbuf = const_cast<void*>(send_buf);
618 int* _send_count = const_cast<int*>(send_count);
619 int* _send_indexes = const_cast<int*>(send_indexes);
620 int _nb_elem = _checkSize(nb_elem);
621 _trace(MpiInfo(eMpiName::Scatterv).name().localstr());
622 double begin_time = MPI_Wtime();
623 m_mpi_prof->scatterVariable(_sbuf,
624 _send_count,
625 _send_indexes,
626 datatype,
627 recv_buf,
628 _nb_elem,
629 datatype,
630 root,
632 double end_time = MPI_Wtime();
633 double sr_time = (end_time - begin_time);
634 //TODO determine the message size
635 m_stat->add(MpiInfo(eMpiName::Scatterv).name(), sr_time, 0);
636}
637
638/*---------------------------------------------------------------------------*/
639/*---------------------------------------------------------------------------*/
640
641void MpiAdapter::
642allToAll(const void* send_buf, void* recv_buf, Integer count, MPI_Datatype datatype)
643{
644 void* _sbuf = const_cast<void*>(send_buf);
645 int icount = _checkSize(count);
646 _trace(MpiInfo(eMpiName::Alltoall).name().localstr());
647 double begin_time = MPI_Wtime();
648 m_mpi_prof->allToAll(_sbuf, icount, datatype, recv_buf, icount, datatype, m_communicator);
649 double end_time = MPI_Wtime();
650 double sr_time = (end_time - begin_time);
651 //TODO determine the message size
652 m_stat->add(MpiInfo(eMpiName::Alltoall).name().localstr(), sr_time, 0);
653}
654
655/*---------------------------------------------------------------------------*/
656/*---------------------------------------------------------------------------*/
657
658Request MpiAdapter::
659nonBlockingAllToAll(const void* send_buf, void* recv_buf, Integer count, MPI_Datatype datatype)
660{
661 MPI_Request mpi_request = MPI_REQUEST_NULL;
662 int ret = -1;
663 void* _sbuf = const_cast<void*>(send_buf);
664 int icount = _checkSize(count);
665 _trace("MPI_IAlltoall");
666 double begin_time = MPI_Wtime();
667 ret = MPI_Ialltoall(_sbuf, icount, datatype, recv_buf, icount, datatype, m_communicator, &mpi_request);
668 double end_time = MPI_Wtime();
669 double sr_time = (end_time - begin_time);
670 //TODO determine the message size
671 m_stat->add("IAllToAll", sr_time, 0);
672 ARCCORE_ADD_REQUEST(mpi_request);
673 return buildRequest(ret, mpi_request);
674}
675
676/*---------------------------------------------------------------------------*/
677/*---------------------------------------------------------------------------*/
678
679void MpiAdapter::
680allToAllVariable(const void* send_buf, const int* send_counts,
681 const int* send_indexes, void* recv_buf, const int* recv_counts,
682 const int* recv_indexes, MPI_Datatype datatype)
683{
684 void* _sbuf = const_cast<void*>(send_buf);
685 int* _send_counts = const_cast<int*>(send_counts);
686 int* _send_indexes = const_cast<int*>(send_indexes);
687 int* _recv_counts = const_cast<int*>(recv_counts);
688 int* _recv_indexes = const_cast<int*>(recv_indexes);
689
690 _trace(MpiInfo(eMpiName::Alltoallv).name().localstr());
691 double begin_time = MPI_Wtime();
692 m_mpi_prof->allToAllVariable(_sbuf, _send_counts, _send_indexes, datatype,
693 recv_buf, _recv_counts, _recv_indexes, datatype, m_communicator);
694 double end_time = MPI_Wtime();
695 double sr_time = (end_time - begin_time);
696 //TODO determine the message size
697 m_stat->add(MpiInfo(eMpiName::Alltoallv).name(), sr_time, 0);
698}
699
700/*---------------------------------------------------------------------------*/
701/*---------------------------------------------------------------------------*/
702
703Request MpiAdapter::
704nonBlockingAllToAllVariable(const void* send_buf, const int* send_counts,
705 const int* send_indexes, void* recv_buf, const int* recv_counts,
706 const int* recv_indexes, MPI_Datatype datatype)
707{
708 MPI_Request mpi_request = MPI_REQUEST_NULL;
709 int ret = -1;
710 void* _sbuf = const_cast<void*>(send_buf);
711 int* _send_counts = const_cast<int*>(send_counts);
712 int* _send_indexes = const_cast<int*>(send_indexes);
713 int* _recv_counts = const_cast<int*>(recv_counts);
714 int* _recv_indexes = const_cast<int*>(recv_indexes);
715
716 _trace("MPI_Ialltoallv");
717 double begin_time = MPI_Wtime();
718 ret = MPI_Ialltoallv(_sbuf, _send_counts, _send_indexes, datatype,
719 recv_buf, _recv_counts, _recv_indexes, datatype,
720 m_communicator, &mpi_request);
721 double end_time = MPI_Wtime();
722 double sr_time = (end_time - begin_time);
723 //TODO determine the message size
724 m_stat->add("IAllToAll", sr_time, 0);
725 ARCCORE_ADD_REQUEST(mpi_request);
726 return buildRequest(ret, mpi_request);
727}
728
729/*---------------------------------------------------------------------------*/
730/*---------------------------------------------------------------------------*/
731
732void MpiAdapter::
733barrier()
734{
735 // TODO: theoretically there should not be any pending requests
736 // between two barriers to avoid any issues.
737 // _checkHasNoRequests();
738 // TODO add corresponding trace for profiling.
739 MPI_Barrier(m_communicator);
740}
741
742/*---------------------------------------------------------------------------*/
743/*---------------------------------------------------------------------------*/
744
745Request MpiAdapter::
746nonBlockingBarrier()
747{
748 MPI_Request mpi_request = MPI_REQUEST_NULL;
749 int ret = -1;
750 ret = MPI_Ibarrier(m_communicator, &mpi_request);
751 ARCCORE_ADD_REQUEST(mpi_request);
752 return buildRequest(ret, mpi_request);
753}
754
755/*---------------------------------------------------------------------------*/
756/*---------------------------------------------------------------------------*/
757
758void MpiAdapter::
759allReduce(const void* send_buf, void* recv_buf, Int64 count, MPI_Datatype datatype, MPI_Op op)
760{
761 void* _sbuf = const_cast<void*>(send_buf);
762 int _n = _checkSize(count);
763 double begin_time = MPI_Wtime();
764 _trace(MpiInfo(eMpiName::Allreduce).name().localstr());
765 try {
766 ++m_nb_all_reduce;
767 m_mpi_prof->allReduce(_sbuf, recv_buf, _n, datatype, op, m_communicator);
768 }
769 catch (TimeoutException& ex) {
770 std::ostringstream ostr;
771 ostr << "MPI_Allreduce"
772 << " send_buf=" << send_buf
773 << " recv_buf=" << recv_buf
774 << " n=" << count
775 << " datatype=" << datatype
776 << " op=" << op
777 << " NB=" << m_nb_all_reduce;
778 ex.setAdditionalInfo(ostr.str());
779 throw;
780 }
781 double end_time = MPI_Wtime();
782 m_stat->add(MpiInfo(eMpiName::Allreduce).name(), end_time - begin_time, count);
783}
784
785/*---------------------------------------------------------------------------*/
786/*---------------------------------------------------------------------------*/
787
788Request MpiAdapter::
789nonBlockingAllReduce(const void* send_buf, void* recv_buf, Int64 count, MPI_Datatype datatype, MPI_Op op)
790{
791 MPI_Request mpi_request = MPI_REQUEST_NULL;
792 int ret = -1;
793 void* _sbuf = const_cast<void*>(send_buf);
794 int _n = _checkSize(count);
795 double begin_time = MPI_Wtime();
796 _trace("MPI_IAllreduce");
797 ret = MPI_Iallreduce(_sbuf, recv_buf, _n, datatype, op, m_communicator, &mpi_request);
798 double end_time = MPI_Wtime();
799 m_stat->add("IReduce", end_time - begin_time, _n);
800 ARCCORE_ADD_REQUEST(mpi_request);
801 return buildRequest(ret, mpi_request);
802}
803
804/*---------------------------------------------------------------------------*/
805/*---------------------------------------------------------------------------*/
806
807void MpiAdapter::
808reduce(const void* send_buf, void* recv_buf, Int64 count, MPI_Datatype datatype, MPI_Op op, Integer root)
809{
810 void* _sbuf = const_cast<void*>(send_buf);
811 int _n = _checkSize(count);
812 int _root = static_cast<int>(root);
813 double begin_time = MPI_Wtime();
814 _trace(MpiInfo(eMpiName::Reduce).name().localstr());
815 try {
816 ++m_nb_reduce;
817 m_mpi_prof->reduce(_sbuf, recv_buf, _n, datatype, op, _root, m_communicator);
818 }
819 catch (TimeoutException& ex) {
820 std::ostringstream ostr;
821 ostr << "MPI_reduce"
822 << " send_buf=" << send_buf
823 << " recv_buf=" << recv_buf
824 << " n=" << count
825 << " datatype=" << datatype
826 << " op=" << op
827 << " root=" << root
828 << " NB=" << m_nb_reduce;
829 ex.setAdditionalInfo(ostr.str());
830 throw;
831 }
832
833 double end_time = MPI_Wtime();
834 m_stat->add(MpiInfo(eMpiName::Reduce).name(), end_time - begin_time, 0);
835}
836
837/*---------------------------------------------------------------------------*/
838/*---------------------------------------------------------------------------*/
839
840void MpiAdapter::
841scan(const void* send_buf, void* recv_buf, Int64 count, MPI_Datatype datatype, MPI_Op op)
842{
843 void* _sbuf = const_cast<void*>(send_buf);
844 int _n = _checkSize(count);
845 double begin_time = MPI_Wtime();
846 _trace(MpiInfo(eMpiName::Scan).name().localstr());
847 m_mpi_prof->scan(_sbuf, recv_buf, _n, datatype, op, m_communicator);
848 double end_time = MPI_Wtime();
849 m_stat->add(MpiInfo(eMpiName::Scan).name(), end_time - begin_time, count);
850}
851
852/*---------------------------------------------------------------------------*/
853/*---------------------------------------------------------------------------*/
854
855void MpiAdapter::
856directSendRecv(const void* send_buffer, Int64 send_buffer_size,
857 void* recv_buffer, Int64 recv_buffer_size,
858 Int32 proc, Int64 elem_size, MPI_Datatype data_type)
859{
860 void* v_send_buffer = const_cast<void*>(send_buffer);
861 MPI_Status mpi_status;
862 double begin_time = MPI_Wtime();
863 _trace(MpiInfo(eMpiName::Sendrecv).name().localstr());
864 int sbuf_size = _checkSize(send_buffer_size);
865 int rbuf_size = _checkSize(recv_buffer_size);
866 m_mpi_prof->sendRecv(v_send_buffer, sbuf_size, data_type, proc, 99,
867 recv_buffer, rbuf_size, data_type, proc, 99,
868 m_communicator, &mpi_status);
869 double end_time = MPI_Wtime();
870 Int64 send_size = send_buffer_size * elem_size;
871 Int64 recv_size = recv_buffer_size * elem_size;
872 double sr_time = (end_time - begin_time);
873
874 //debug(Trace::High) << "MPI SendRecv: send " << send_size << " recv "
875 // << recv_size << " time " << sr_time ;
876 m_stat->add(MpiInfo(eMpiName::Sendrecv).name(), sr_time, send_size + recv_size);
877}
878
879/*---------------------------------------------------------------------------*/
880/*---------------------------------------------------------------------------*/
881
883sendNonBlockingNoStat(const void* send_buffer, Int64 send_buffer_size,
884 Int32 dest_rank, MPI_Datatype data_type, int mpi_tag)
885{
886 void* v_send_buffer = const_cast<void*>(send_buffer);
887 MPI_Request mpi_request = MPI_REQUEST_NULL;
888 int sbuf_size = _checkSize(send_buffer_size);
889 int ret = 0;
890 m_mpi_prof->iSend(v_send_buffer, sbuf_size, data_type, dest_rank, mpi_tag, m_communicator, &mpi_request);
891 if (m_is_trace)
892 info() << " ISend ret=" << ret << " proc=" << dest_rank << " tag=" << mpi_tag << " request=" << mpi_request;
893 ARCCORE_ADD_REQUEST(mpi_request);
894 return buildRequest(ret, mpi_request);
895}
896
897/*---------------------------------------------------------------------------*/
898/*---------------------------------------------------------------------------*/
899
900Request MpiAdapter::
901directSend(const void* send_buffer, Int64 send_buffer_size,
902 Int32 proc, Int64 elem_size, MPI_Datatype data_type,
903 int mpi_tag, bool is_blocked)
904{
905 void* v_send_buffer = const_cast<void*>(send_buffer);
906 MPI_Request mpi_request = MPI_REQUEST_NULL;
907
908 double begin_time = 0.0;
909 double end_time = 0.0;
910 Int64 send_size = send_buffer_size * elem_size;
911 int ret = 0;
912 if (m_is_trace)
913 info() << "MPI_TRACE: MPI Send: send before"
914 << " size=" << send_size
915 << " dest=" << proc
916 << " tag=" << mpi_tag
917 << " datatype=" << data_type
918 << " blocking " << is_blocked;
919 if (is_blocked) {
920 // if m_mpi_lock is not null, we must
921 // use an MPI_ISend followed by an
922 // active MPI_Test loop to avoid any
923 // dead lock issues.
924 if (m_mpi_lock) {
925 {
926 MpiLock::Section mls(m_mpi_lock);
927 begin_time = MPI_Wtime();
928 int sbuf_size = _checkSize(send_buffer_size);
929 m_mpi_prof->iSend(v_send_buffer, sbuf_size, data_type, proc, mpi_tag, m_communicator, &mpi_request);
930 }
931 int is_finished = 0;
932 MPI_Status mpi_status;
933 while (is_finished == 0) {
934 MpiLock::Section mls(m_mpi_lock);
935 MPI_Request_get_status(mpi_request, &is_finished, &mpi_status);
936 if (is_finished != 0) {
937 m_mpi_prof->wait(&mpi_request, (MPI_Status*)MPI_STATUS_IGNORE);
938 end_time = MPI_Wtime();
939 mpi_request = MPI_REQUEST_NULL;
940 }
941 }
942 }
943 else {
944 MpiLock::Section mls(m_mpi_lock);
945 begin_time = MPI_Wtime();
946 int sbuf_size = _checkSize(send_buffer_size);
947 m_mpi_prof->send(v_send_buffer, sbuf_size, data_type, proc, mpi_tag, m_communicator);
948 end_time = MPI_Wtime();
949 }
950 }
951 else {
952 {
953 MpiLock::Section mls(m_mpi_lock);
954 begin_time = MPI_Wtime();
955 int sbuf_size = _checkSize(send_buffer_size);
956 m_mpi_prof->iSend(v_send_buffer, sbuf_size, data_type, proc, mpi_tag, m_communicator, &mpi_request);
957 if (m_is_trace)
958 info() << " ISend ret=" << ret << " proc=" << proc << " tag=" << mpi_tag << " request=" << mpi_request;
959 end_time = MPI_Wtime();
960 ARCCORE_ADD_REQUEST(mpi_request);
961 }
962 if (m_is_trace) {
963 info() << "MPI Send: send after"
964 << " request=" << mpi_request;
965 }
966 }
967 double sr_time = (end_time - begin_time);
968
969 debug(Trace::High) << "MPI Send: send " << send_size
970 << " time " << sr_time << " blocking " << is_blocked;
971 // TODO(FL): look into how to profile Isend
972 m_stat->add(MpiInfo(eMpiName::Send).name(), end_time - begin_time, send_size);
973 return buildRequest(ret, mpi_request);
974}
975
976/*---------------------------------------------------------------------------*/
977/*---------------------------------------------------------------------------*/
978
979Request MpiAdapter::
980directSendPack(const void* send_buffer, Int64 send_buffer_size,
981 Int32 proc, int mpi_tag, bool is_blocked)
982{
983 return directSend(send_buffer, send_buffer_size, proc, 1, MPI_PACKED, mpi_tag, is_blocked);
984}
985
986/*---------------------------------------------------------------------------*/
987/*---------------------------------------------------------------------------*/
988
989MpiMessagePassingMng* MpiAdapter::
990commSplit(bool keep)
991{
992 MPI_Comm new_comm;
993
994 MPI_Comm_split(m_communicator, (keep) ? 1 : MPI_UNDEFINED, commRank(), &new_comm);
995 if (keep) {
996 // Failed if new_comm is MPI_COMM_NULL
997 return StandaloneMpiMessagePassingMng::create(new_comm, true);
998 }
999 return nullptr;
1000}
1001
1002/*---------------------------------------------------------------------------*/
1003/*---------------------------------------------------------------------------*/
1004
1006receiveNonBlockingNoStat(void* recv_buffer, Int64 recv_buffer_size,
1007 Int32 source_rank, MPI_Datatype data_type, int mpi_tag)
1008{
1009 int rbuf_size = _checkSize(recv_buffer_size);
1010 int ret = 0;
1011 MPI_Request mpi_request = MPI_REQUEST_NULL;
1012 m_mpi_prof->iRecv(recv_buffer, rbuf_size, data_type, source_rank, mpi_tag, m_communicator, &mpi_request);
1013 ARCCORE_ADD_REQUEST(mpi_request);
1014 return buildRequest(ret, mpi_request);
1015}
1016
1017/*---------------------------------------------------------------------------*/
1018/*---------------------------------------------------------------------------*/
1019
1020Request MpiAdapter::
1021directRecv(void* recv_buffer, Int64 recv_buffer_size,
1022 Int32 proc, Int64 elem_size, MPI_Datatype data_type,
1023 int mpi_tag, bool is_blocked)
1024{
1025 MPI_Status mpi_status;
1026 MPI_Request mpi_request = MPI_REQUEST_NULL;
1027 int ret = 0;
1028 double begin_time = 0.0;
1029 double end_time = 0.0;
1030
1031 int i_proc = 0;
1032 if (proc == A_PROC_NULL_RANK)
1033 ARCCORE_THROW(NotImplementedException, "Receive with MPI_PROC_NULL");
1034 if (proc == A_NULL_RANK && !m_is_allow_null_rank_for_any_source)
1035 ARCCORE_FATAL("Can not use A_NULL_RANK for any source. Use A_ANY_SOURCE_RANK instead");
1036 if (proc == A_NULL_RANK || proc == A_ANY_SOURCE_RANK)
1037 i_proc = MPI_ANY_SOURCE;
1038 else
1039 i_proc = static_cast<int>(proc);
1040
1041 Int64 recv_size = recv_buffer_size * elem_size;
1042 if (m_is_trace) {
1043 info() << "MPI_TRACE: MPI Recv: recv before "
1044 << " size=" << recv_size
1045 << " from=" << i_proc
1046 << " tag=" << mpi_tag
1047 << " datatype=" << data_type
1048 << " blocking=" << is_blocked;
1049 }
1050 if (is_blocked) {
1051 // if m_mpi_lock is not null, we must
1052 // use an MPI_IRecv followed by an
1053 // active MPI_Test loop to avoid any
1054 // dead lock issues.
1055 if (m_mpi_lock) {
1056 {
1057 MpiLock::Section mls(m_mpi_lock);
1058 begin_time = MPI_Wtime();
1059 int rbuf_size = _checkSize(recv_buffer_size);
1060 m_mpi_prof->iRecv(recv_buffer, rbuf_size, data_type, i_proc, mpi_tag, m_communicator, &mpi_request);
1061 }
1062 int is_finished = 0;
1063 MPI_Status mpi_status;
1064 while (is_finished == 0) {
1065 MpiLock::Section mls(m_mpi_lock);
1066 MPI_Request_get_status(mpi_request, &is_finished, &mpi_status);
1067 if (is_finished != 0) {
1068 end_time = MPI_Wtime();
1069 m_mpi_prof->wait(&mpi_request, (MPI_Status*)MPI_STATUS_IGNORE);
1070 mpi_request = MPI_REQUEST_NULL;
1071 }
1072 }
1073 }
1074 else {
1075 MpiLock::Section mls(m_mpi_lock);
1076 begin_time = MPI_Wtime();
1077 int rbuf_size = _checkSize(recv_buffer_size);
1078 m_mpi_prof->recv(recv_buffer, rbuf_size, data_type, i_proc, mpi_tag, m_communicator, &mpi_status);
1079 end_time = MPI_Wtime();
1080 }
1081 }
1082 else {
1083 {
1084 MpiLock::Section mls(m_mpi_lock);
1085 begin_time = MPI_Wtime();
1086 int rbuf_size = _checkSize(recv_buffer_size);
1087 m_mpi_prof->iRecv(recv_buffer, rbuf_size, data_type, i_proc, mpi_tag, m_communicator, &mpi_request);
1088 end_time = MPI_Wtime();
1089 ARCCORE_ADD_REQUEST(mpi_request);
1090 }
1091 if (m_is_trace) {
1092 info() << "MPI Recv: recv after "
1093 << " request=" << mpi_request;
1094 }
1095 }
1096 double sr_time = (end_time - begin_time);
1097
1098 debug(Trace::High) << "MPI Recv: recv after " << recv_size
1099 << " time " << sr_time << " blocking " << is_blocked;
1100 m_stat->add(MpiInfo(eMpiName::Recv).name(), end_time - begin_time, recv_size);
1101 return buildRequest(ret, mpi_request);
1102}
1103
1104/*---------------------------------------------------------------------------*/
1105/*---------------------------------------------------------------------------*/
1106
1107void MpiAdapter::
1108probeRecvPack(UniqueArray<Byte>& recv_buffer, Int32 proc)
1109{
1110 double begin_time = MPI_Wtime();
1111 MPI_Status status;
1112 int recv_buffer_size = 0;
1113 _trace("MPI_Probe");
1114 m_mpi_prof->probe(proc, 101, m_communicator, &status);
1115 m_mpi_prof->getCount(&status, MPI_PACKED, &recv_buffer_size);
1116
1117 recv_buffer.resize(recv_buffer_size);
1118 m_mpi_prof->recv(recv_buffer.data(), recv_buffer_size, MPI_PACKED, proc, 101, m_communicator, &status);
1119
1120 double end_time = MPI_Wtime();
1121 Int64 recv_size = recv_buffer_size;
1122 double sr_time = (end_time - begin_time);
1123 debug(Trace::High) << "MPI probeRecvPack " << recv_size
1124 << " time " << sr_time;
1125 m_stat->add(MpiInfo(eMpiName::Recv).name(), end_time - begin_time, recv_size);
1126}
1127
1128/*---------------------------------------------------------------------------*/
1129/*---------------------------------------------------------------------------*/
1130
1131MessageSourceInfo MpiAdapter::
1132_buildSourceInfoFromStatus(const MPI_Status& mpi_status)
1133{
1134 // Retrieves the message size in bytes.
1135 MPI_Count message_size = 0;
1136 MPI_Get_elements_x(&mpi_status, MPI_BYTE, &message_size);
1137 MessageTag tag(mpi_status.MPI_TAG);
1138 MessageRank rank(mpi_status.MPI_SOURCE);
1139 return MessageSourceInfo(rank, tag, message_size);
1140}
1141
1142/*---------------------------------------------------------------------------*/
1143/*---------------------------------------------------------------------------*/
1144
1145MessageId MpiAdapter::
1146_probeMessage(MessageRank source, MessageTag tag, bool is_blocking)
1147{
1148 MPI_Status mpi_status;
1149 int has_message = 0;
1150 MPI_Message message;
1151 int ret = 0;
1152 int mpi_source = source.value();
1153 if (source.isProcNull())
1154 ARCCORE_THROW(NotImplementedException, "Probe with MPI_PROC_NULL");
1155 if (source.isNull() && !m_is_allow_null_rank_for_any_source)
1156 ARCCORE_FATAL("Can not use MPI_Mprobe with null rank. Use MessageRank::anySourceRank() instead");
1157 if (source.isNull() || source.isAnySource())
1158 mpi_source = MPI_ANY_SOURCE;
1159 int mpi_tag = tag.value();
1160 if (tag.isNull())
1161 mpi_tag = MPI_ANY_TAG;
1162 if (is_blocking) {
1163 ret = MPI_Mprobe(mpi_source, mpi_tag, m_communicator, &message, &mpi_status);
1164 has_message = true;
1165 }
1166 else {
1167 ret = MPI_Improbe(mpi_source, mpi_tag, m_communicator, &has_message, &message, &mpi_status);
1168 }
1169 if (ret != 0)
1170 ARCCORE_FATAL("Error during call to MPI_Mprobe r={0}", ret);
1171 MessageId ret_message;
1172 if (has_message != 0) {
1173 MessageSourceInfo si(_buildSourceInfoFromStatus(mpi_status));
1174 ret_message = MessageId(si, message);
1175 }
1176 return ret_message;
1177}
1178
1179/*---------------------------------------------------------------------------*/
1180/*---------------------------------------------------------------------------*/
1181
1182MessageId MpiAdapter::
1183probeMessage(PointToPointMessageInfo message)
1184{
1185 if (!message.isValid())
1186 return MessageId();
1187
1188 // The message must be initialized with a (rank/tag) pair.
1189 if (!message.isRankTag())
1190 ARCCORE_FATAL("Invalid message_info: message.isRankTag() is false");
1191
1192 return _probeMessage(message.destinationRank(), message.tag(), message.isBlocking());
1193}
1194
1195/*---------------------------------------------------------------------------*/
1196/*---------------------------------------------------------------------------*/
1197
1198MessageSourceInfo MpiAdapter::
1199_legacyProbeMessage(MessageRank source, MessageTag tag, bool is_blocking)
1200{
1201 MPI_Status mpi_status;
1202 int has_message = 0;
1203 int ret = 0;
1204 int mpi_source = source.value();
1205 if (source.isProcNull())
1206 ARCCORE_THROW(NotImplementedException, "Probe with MPI_PROC_NULL");
1207 if (source.isNull() && !m_is_allow_null_rank_for_any_source)
1208 ARCCORE_FATAL("Can not use MPI_Probe with null rank. Use MessageRank::anySourceRank() instead");
1209 if (source.isNull() || source.isAnySource())
1210 mpi_source = MPI_ANY_SOURCE;
1211 int mpi_tag = tag.value();
1212 if (tag.isNull())
1213 mpi_tag = MPI_ANY_TAG;
1214 if (is_blocking) {
1215 ret = MPI_Probe(mpi_source, mpi_tag, m_communicator, &mpi_status);
1216 has_message = true;
1217 }
1218 else
1219 ret = MPI_Iprobe(mpi_source, mpi_tag, m_communicator, &has_message, &mpi_status);
1220 if (ret != 0)
1221 ARCCORE_FATAL("Error during call to MPI_Mprobe r={0}", ret);
1222 if (has_message != 0)
1223 return _buildSourceInfoFromStatus(mpi_status);
1224 return {};
1225}
1226
1227/*---------------------------------------------------------------------------*/
1228/*---------------------------------------------------------------------------*/
1229
1230MessageSourceInfo MpiAdapter::
1231legacyProbeMessage(PointToPointMessageInfo message)
1232{
1233 if (!message.isValid())
1234 return {};
1235
1236 // The message must be initialized with a (rank/tag) pair.
1237 if (!message.isRankTag())
1238 ARCCORE_FATAL("Invalid message_info: message.isRankTag() is false");
1239
1240 return _legacyProbeMessage(message.destinationRank(), message.tag(), message.isBlocking());
1241}
1242
1243/*---------------------------------------------------------------------------*/
1244/*---------------------------------------------------------------------------*/
1245
1247Request MpiAdapter::
1248directRecv(void* recv_buffer, Int64 recv_buffer_size,
1249 MessageId message, Int64 elem_size, MPI_Datatype data_type,
1250 bool is_blocked)
1251{
1252 MPI_Status mpi_status;
1253 MPI_Request mpi_request = MPI_REQUEST_NULL;
1254 MPI_Message mpi_message = (MPI_Message)message;
1255 int ret = 0;
1256 double begin_time = 0.0;
1257 double end_time = 0.0;
1258
1259 Int64 recv_size = recv_buffer_size * elem_size;
1260 if (m_is_trace) {
1261 info() << "MPI_TRACE: MPI Mrecv: recv before "
1262 << " size=" << recv_size
1263 << " from_msg=" << message
1264 << " datatype=" << data_type
1265 << " blocking=" << is_blocked;
1266 }
1267 if (is_blocked) {
1268 // if m_mpi_lock is not null, we must
1269 // use an MPI_IRecv followed by an
1270 // active MPI_Test loop to avoid any
1271 // dead lock issues.
1272 if (m_mpi_lock) {
1273 {
1274 MpiLock::Section mls(m_mpi_lock);
1275 begin_time = MPI_Wtime();
1276 int rbuf_size = _checkSize(recv_buffer_size);
1277 MPI_Imrecv(recv_buffer, rbuf_size, data_type, &mpi_message, &mpi_request);
1278 //m_mpi_prof->iRecv(recv_buffer, rbuf_size, data_type, i_proc, mpi_tag, m_communicator, &mpi_request);
1279 }
1280 int is_finished = 0;
1281 MPI_Status mpi_status;
1282 while (is_finished == 0) {
1283 MpiLock::Section mls(m_mpi_lock);
1284 MPI_Request_get_status(mpi_request, &is_finished, &mpi_status);
1285 if (is_finished != 0) {
1286 end_time = MPI_Wtime();
1287 m_mpi_prof->wait(&mpi_request, (MPI_Status*)MPI_STATUS_IGNORE);
1288 mpi_request = MPI_REQUEST_NULL;
1289 }
1290 }
1291 }
1292 else {
1293 MpiLock::Section mls(m_mpi_lock);
1294 begin_time = MPI_Wtime();
1295 int rbuf_size = _checkSize(recv_buffer_size);
1296 MPI_Mrecv(recv_buffer, rbuf_size, data_type, &mpi_message, &mpi_status);
1297 //m_mpi_prof->recv(recv_buffer, rbuf_size, data_type, i_proc, mpi_tag, m_communicator, &mpi_status);
1298 end_time = MPI_Wtime();
1299 }
1300 }
1301 else {
1302 {
1303 MpiLock::Section mls(m_mpi_lock);
1304 begin_time = MPI_Wtime();
1305 int rbuf_size = _checkSize(recv_buffer_size);
1306 //m_mpi_prof->iRecv(recv_buffer, rbuf_size, data_type, i_proc, mpi_tag, m_communicator, &mpi_request);
1307 ret = MPI_Imrecv(recv_buffer, rbuf_size, data_type, &mpi_message, &mpi_request);
1308 //m_mpi_prof->iRecv(recv_buffer, rbuf_size, data_type, i_proc, mpi_tag, m_communicator, &mpi_request);
1309 end_time = MPI_Wtime();
1310 ARCCORE_ADD_REQUEST(mpi_request);
1311 }
1312 if (m_is_trace) {
1313 info() << "MPI Recv: recv after "
1314 << " request=" << mpi_request;
1315 }
1316 }
1317 double sr_time = (end_time - begin_time);
1318
1319 debug(Trace::High) << "MPI Recv: recv after " << recv_size
1320 << " time " << sr_time << " blocking " << is_blocked;
1321 m_stat->add(MpiInfo(eMpiName::Recv).name(), end_time - begin_time, recv_size);
1322 return buildRequest(ret, mpi_request);
1323}
1324
1325/*---------------------------------------------------------------------------*/
1326/*---------------------------------------------------------------------------*/
1327
1328Request MpiAdapter::
1329directRecvPack(void* recv_buffer, Int64 recv_buffer_size,
1330 Int32 proc, int mpi_tag, bool is_blocking)
1331{
1332 return directRecv(recv_buffer, recv_buffer_size, proc, 1, MPI_PACKED, mpi_tag, is_blocking);
1333}
1334
1335/*---------------------------------------------------------------------------*/
1336/*---------------------------------------------------------------------------*/
1337
1338// FIXME: Implement direct method with MPI_STATUS_IGNORE
1339void MpiAdapter::
1340waitAllRequests(ArrayView<Request> requests)
1341{
1342 UniqueArray<bool> indexes(requests.size());
1343 UniqueArray<MPI_Status> mpi_status(requests.size());
1344 while (_waitAllRequestsMPI(requests, indexes, mpi_status)) {
1345 ; // Continue as long as there are requests.
1346 }
1347}
1348
1349/*---------------------------------------------------------------------------*/
1350/*---------------------------------------------------------------------------*/
1351
1352// FIXME: Implement direct method with MPI_STATUS_IGNORE
1353void MpiAdapter::
1354waitSomeRequests(ArrayView<Request> requests,
1355 ArrayView<bool> indexes,
1356 bool is_non_blocking)
1357{
1358 UniqueArray<MPI_Status> mpi_status(requests.size());
1359 waitSomeRequestsMPI(requests, indexes, mpi_status, is_non_blocking);
1360}
1361
1362/*---------------------------------------------------------------------------*/
1363/*---------------------------------------------------------------------------*/
1364
1366{
1367 SubRequestInfo(Ref<ISubRequest> sr, Integer i, int source_rank, int source_tag)
1368 : sub_request(sr)
1369 , index(i)
1370 , mpi_source_rank(source_rank)
1371 , mpi_source_tag(source_tag)
1372 {}
1373
1374 Ref<ISubRequest> sub_request;
1375 Integer index = -1;
1376 int mpi_source_rank = MPI_PROC_NULL;
1377 int mpi_source_tag = 0;
1378};
1379
1380/*---------------------------------------------------------------------------*/
1381/*---------------------------------------------------------------------------*/
1382
1383bool MpiAdapter::
1384_handleEndRequests(ArrayView<Request> requests, ArrayView<bool> done_indexes,
1385 ArrayView<MPI_Status> status)
1386{
1387 UniqueArray<SubRequestInfo> new_requests;
1388 Integer size = requests.size();
1389 {
1390 MpiLock::Section mls(m_mpi_lock);
1391 for (Integer i = 0; i < size; ++i) {
1392 if (done_indexes[i]) {
1393 // Be careful to use a reference, otherwise the reset won't
1394 // apply to the correct variable
1395 Request& r = requests[i];
1396 // Note: the request might not be valid (for example, if it is)
1397 // a blocking request but still have a sub-request.
1398 if (r.hasSubRequest()) {
1399 if (m_is_trace)
1400 info() << "Done request with sub-request r=" << r << " mpi_r=" << r << " i=" << i
1401 << " source_rank=" << status[i].MPI_SOURCE
1402 << " source_tag=" << status[i].MPI_TAG;
1403 new_requests.add(SubRequestInfo(r.subRequest(), i, status[i].MPI_SOURCE, status[i].MPI_TAG));
1404 }
1405 if (r.isValid()) {
1406 _removeRequest((MPI_Request)(r));
1407 r.reset();
1408 }
1409 }
1410 }
1411 }
1412
1413 // NOTE: calls to sub-requests can generate other requests.
1414 // Care must be taken not to use sub-requests while the lock is active.
1415 bool has_new_request = false;
1416 if (!new_requests.empty()) {
1417 // Contains the status of the i-th request
1418 UniqueArray<MPI_Status> old_status(size);
1419 {
1420 Integer index = 0;
1421 for (Integer i = 0; i < size; ++i) {
1422 if (done_indexes[i]) {
1423 old_status[i] = status[index];
1424 ++index;
1425 }
1426 }
1427 }
1428 // If there are new requests, the values in 'status' must be shifted
1429 for (SubRequestInfo& sri : new_requests) {
1430 Integer index = sri.index;
1431 if (m_is_trace)
1432 info() << "Before handle new request index=" << index
1433 << " sri.source_rank=" << sri.mpi_source_rank
1434 << " sri.source_tag=" << sri.mpi_source_tag;
1435 SubRequestCompletionInfo completion_info(MessageRank(old_status[index].MPI_SOURCE), MessageTag(old_status[index].MPI_TAG));
1436 Request r = sri.sub_request->executeOnCompletion(completion_info);
1437 if (m_is_trace)
1438 info() << "Handle new request index=" << index << " old_r=" << requests[index] << " new_r=" << r;
1439 // If there is a new request, it replaces
1440 // the old one, so we must act as if
1441 // the original request is not finished.
1442 if (r.isValid()) {
1443 has_new_request = true;
1444 requests[index] = r;
1445 done_indexes[index] = false;
1446 }
1447 }
1448 {
1449 Integer index = 0;
1450 for (Integer i = 0; i < size; ++i) {
1451 if (done_indexes[i]) {
1452 status[index] = old_status[i];
1453 ++index;
1454 }
1455 }
1456 }
1457 }
1458 return has_new_request;
1459}
1460
1461/*---------------------------------------------------------------------------*/
1462/*---------------------------------------------------------------------------*/
1463
1464bool MpiAdapter::
1465_waitAllRequestsMPI(ArrayView<Request> requests,
1466 ArrayView<bool> indexes,
1467 ArrayView<MPI_Status> mpi_status)
1468{
1469 Integer size = requests.size();
1470 if (size == 0)
1471 return false;
1472 //ATTENTION: Mpi modifies this array upon return from MPI_Waitall
1473 UniqueArray<MPI_Request> mpi_request(size);
1474 for (Integer i = 0; i < size; ++i) {
1475 mpi_request[i] = (MPI_Request)(requests[i]);
1476 }
1477 if (m_is_trace)
1478 info() << " MPI_waitall begin size=" << size;
1479 double diff_time = 0.0;
1480 if (m_mpi_lock) {
1481 double begin_time = MPI_Wtime();
1482 for (Integer i = 0; i < size; ++i) {
1483 MPI_Request request = (MPI_Request)(mpi_request[i]);
1484 int is_finished = 0;
1485 while (is_finished == 0) {
1486 MpiLock::Section mls(m_mpi_lock);
1487 m_mpi_prof->test(&request, &is_finished, (MPI_Status*)MPI_STATUS_IGNORE);
1488 }
1489 }
1490 double end_time = MPI_Wtime();
1491 diff_time = end_time - begin_time;
1492 }
1493 else {
1494 //TODO: transform into a while loop and MPI_Testall if m_mpi_lock is non-null
1495 MpiLock::Section mls(m_mpi_lock);
1496 double begin_time = MPI_Wtime();
1497 m_mpi_prof->waitAll(size, mpi_request.data(), mpi_status.data());
1498 double end_time = MPI_Wtime();
1499 diff_time = end_time - begin_time;
1500 }
1501
1502 // Indicates that each request has been processed because we performed a waitall.
1503 for (Integer i = 0; i < size; ++i) {
1504 indexes[i] = true;
1505 }
1506
1507 bool has_new_request = _handleEndRequests(requests, indexes, mpi_status);
1508 if (m_is_trace)
1509 info() << " MPI_waitall end size=" << size;
1510 m_stat->add(MpiInfo(eMpiName::Waitall).name(), diff_time, size);
1511 return has_new_request;
1512}
1513
1514/*---------------------------------------------------------------------------*/
1515/*---------------------------------------------------------------------------*/
1516
1517void MpiAdapter::
1518waitSomeRequestsMPI(ArrayView<Request> requests, ArrayView<bool> indexes,
1519 ArrayView<MPI_Status> mpi_status, bool is_non_blocking)
1520{
1521 Integer size = requests.size();
1522 if (size == 0)
1523 return;
1524 //TODO: use StackArray (when they become available...)
1525 UniqueArray<MPI_Request> mpi_request(size);
1526 UniqueArray<MPI_Request> saved_mpi_request(size);
1527 UniqueArray<int> completed_requests(size);
1528 int nb_completed_request = 0;
1529
1530 // Save the request to deallocate it in m_allocated_requests,
1531 // because its value will no longer be valid after calling MPI_Wait*
1532 for (Integer i = 0; i < size; ++i) {
1533 // In the case where this method is called multiple times
1534 // with the same requests array, there may be invalid
1535 // requests that will crash the MPI call.
1536 if (!requests[i].isValid()) {
1537 saved_mpi_request[i] = MPI_REQUEST_NULL;
1538 }
1539 else {
1540 saved_mpi_request[i] = static_cast<MPI_Request>(requests[i]);
1541 }
1542 }
1543
1544 // Only display debug in blocking mode or if explicitly requested to
1545 // avoid too many messages
1546 bool is_print_debug = m_is_trace || (!is_non_blocking);
1547 if (is_print_debug)
1548 debug() << "WaitRequestBegin is_non_blocking=" << is_non_blocking << " n=" << size;
1549
1550 double begin_time = MPI_Wtime();
1551
1552 try {
1553 if (is_non_blocking) {
1554 _trace(MpiInfo(eMpiName::Testsome).name().localstr());
1555 {
1556 MpiLock::Section mls(m_mpi_lock);
1557 m_mpi_prof->testSome(size, saved_mpi_request.data(), &nb_completed_request,
1558 completed_requests.data(), mpi_status.data());
1559 }
1560 //If there is no active handle in the list, it returns outcount = MPI_UNDEFINED.
1561 if (nb_completed_request == MPI_UNDEFINED) // If no requests were valid.
1562 nb_completed_request = 0;
1563 if (is_print_debug)
1564 debug() << "WaitSomeRequestMPI: TestSome nb_completed=" << nb_completed_request;
1565 }
1566 else {
1567 _trace(MpiInfo(eMpiName::Waitsome).name().localstr());
1568 {
1569 // TODO: if the lock exists, a testSome() loop must be performed
1570 // so as not to block.
1571 MpiLock::Section mls(m_mpi_lock);
1572 m_mpi_prof->waitSome(size, saved_mpi_request.data(), &nb_completed_request,
1573 completed_requests.data(), mpi_status.data());
1574 }
1575 // One must not use mpi_request[i] because it is modified by Mpi
1576 // mpi_request[i] == MPI_REQUEST_NULL
1577 if (nb_completed_request == MPI_UNDEFINED) // If no requests were valid.
1578 nb_completed_request = 0;
1579 if (is_print_debug)
1580 debug() << "WaitSomeRequest nb_completed=" << nb_completed_request;
1581 }
1582 }
1583 catch (TimeoutException& ex) {
1584 std::ostringstream ostr;
1585 if (is_non_blocking)
1586 ostr << MpiInfo(eMpiName::Testsome).name();
1587 else
1588 ostr << MpiInfo(eMpiName::Waitsome).name();
1589 ostr << " size=" << size
1590 << " is_non_blocking=" << is_non_blocking;
1591 ex.setAdditionalInfo(ostr.str());
1592 throw;
1593 }
1594
1595 for (int z = 0; z < nb_completed_request; ++z) {
1596 int index = completed_requests[z];
1597 if (is_print_debug)
1598 debug() << "Completed my_rank=" << m_comm_rank << " z=" << z
1599 << " index=" << index
1600 << " tag=" << mpi_status[z].MPI_TAG
1601 << " source=" << mpi_status[z].MPI_SOURCE;
1602
1603 indexes[index] = true;
1604 }
1605
1606 bool has_new_request = _handleEndRequests(requests, indexes, mpi_status);
1607 if (has_new_request) {
1608 // If there are new requests, it is possible that no
1609 // request has completed. In the case of testSome, this is not serious.
1610 // In the case of waitSome, this means that we must wait again.
1611 }
1612 double end_time = MPI_Wtime();
1613 m_stat->add(MpiInfo(eMpiName::Waitsome).name(), end_time - begin_time, size);
1614}
1615
1616/*---------------------------------------------------------------------------*/
1617/*---------------------------------------------------------------------------*/
1618
1619void MpiAdapter::
1620freeRequest(Request& request)
1621{
1622 if (!request.isValid()) {
1623 warning() << "MpiAdapter::freeRequest() null request r=" << (MPI_Request)request;
1624 _checkFatalInRequest();
1625 return;
1626 }
1627 {
1628 MpiLock::Section mls(m_mpi_lock);
1629
1630 auto mr = (MPI_Request)request;
1631 _removeRequest(mr);
1632 MPI_Request_free(&mr);
1633 }
1634 request.reset();
1635}
1636
1637/*---------------------------------------------------------------------------*/
1638/*---------------------------------------------------------------------------*/
1639
1640bool MpiAdapter::
1641testRequest(Request& request)
1642{
1643 // It is allowed by MPI to perform a test with a null request.
1644 if (!request.isValid())
1645 return true;
1646
1647 auto mr = (MPI_Request)request;
1648 int is_finished = 0;
1649
1650 {
1651 MpiLock::Section mls(m_mpi_lock);
1652
1653 // First, we must retrieve the request location because if it
1654 // is finished, it will be automatically freed by MPI during the test
1655 // and thus we will no longer be able to remove it
1656 RequestSet::Iterator request_iter = m_request_set->findRequest(mr);
1657
1658 m_mpi_prof->test(&mr, &is_finished, (MPI_Status*)MPI_STATUS_IGNORE);
1659 //info() << "** TEST REQUEST r=" << mr << " is_finished=" << is_finished;
1660 if (is_finished != 0) {
1661 m_request_set->removeRequest(request_iter);
1662 if (request.hasSubRequest())
1663 ARCCORE_THROW(NotImplementedException, "SubRequest support");
1664 request.reset();
1665 return true;
1666 }
1667 }
1668
1669 return false;
1670}
1671
1672/*---------------------------------------------------------------------------*/
1673/*---------------------------------------------------------------------------*/
1674
1679_addRequest(MPI_Request request)
1680{
1681 m_request_set->addRequest(request);
1682}
1683
1684/*---------------------------------------------------------------------------*/
1685/*---------------------------------------------------------------------------*/
1686
1691_removeRequest(MPI_Request request)
1692{
1693 m_request_set->removeRequest(request);
1694}
1695
1696/*---------------------------------------------------------------------------*/
1697/*---------------------------------------------------------------------------*/
1698
1699void MpiAdapter::
1700enableDebugRequest(bool enable_debug_request)
1701{
1702 m_stat->enable(enable_debug_request);
1703 //if (!m_enable_debug_request)
1704 //info() << "WARNING: Mpi adpater debug request is disabled (multi-threading)";
1705}
1706
1707/*---------------------------------------------------------------------------*/
1708/*---------------------------------------------------------------------------*/
1709
1710void MpiAdapter::
1711_checkFatalInRequest()
1712{
1713 if (isRequestErrorAreFatal())
1714 ARCCORE_FATAL("Error in requests management");
1715}
1716
1717/*---------------------------------------------------------------------------*/
1718/*---------------------------------------------------------------------------*/
1719
1720void MpiAdapter::
1721setMpiProfiling(IMpiProfiling* mpi_profiling)
1722{
1723 m_mpi_prof = mpi_profiling;
1724}
1725
1726/*---------------------------------------------------------------------------*/
1727/*---------------------------------------------------------------------------*/
1728
1729IMpiProfiling* MpiAdapter::
1730getMpiProfiling() const
1731{
1732 return m_mpi_prof;
1733}
1734
1735/*---------------------------------------------------------------------------*/
1736/*---------------------------------------------------------------------------*/
1737
1738void MpiAdapter::
1739setProfiler(IProfiler* profiler)
1740{
1741 if (!profiler) {
1742 m_mpi_prof = nullptr;
1743 return;
1744 }
1745
1746 IMpiProfiling* p = dynamic_cast<IMpiProfiling*>(profiler);
1747 if (!p)
1748 ARCCORE_FATAL("Invalid profiler. Profiler has to implemented interface 'IMpiProfiling'");
1749 m_mpi_prof = p;
1750}
1751
1752/*---------------------------------------------------------------------------*/
1753/*---------------------------------------------------------------------------*/
1754
1755IProfiler* MpiAdapter::
1756profiler() const
1757{
1758 return m_mpi_prof;
1759}
1760
1761/*---------------------------------------------------------------------------*/
1762/*---------------------------------------------------------------------------*/
1763
1764// Attention : Non thread-safe !
1765void MpiAdapter::
1766initializeWindowCreator(MPI_Comm comm_machine)
1767{
1768 if (m_window_creator.isNull()) {
1769 Integer machine_comm_rank = 0;
1770 Integer machine_comm_size = 0;
1771 ::MPI_Comm_rank(comm_machine, &machine_comm_rank);
1772 ::MPI_Comm_size(comm_machine, &machine_comm_size);
1773 m_window_creator = makeRef(new MpiMachineShMemWinBaseInternalCreator(comm_machine, machine_comm_rank, machine_comm_size, m_communicator, m_comm_size));
1774 }
1775}
1776
1777/*---------------------------------------------------------------------------*/
1778/*---------------------------------------------------------------------------*/
1779
1781windowCreator() const
1782{
1783 return m_window_creator.get();
1784}
1785
1786/*---------------------------------------------------------------------------*/
1787/*---------------------------------------------------------------------------*/
1788
1789} // End namespace Arcane::MessagePassing::Mpi
1790
1791/*---------------------------------------------------------------------------*/
1792/*---------------------------------------------------------------------------*/
#define ARCCORE_FATAL(...)
Macro throwing a FatalErrorException.
#define ARCCORE_THROW(exception_class,...)
Macro to throw an exception with formatting.
bool empty() const
Capacity (number of allocated elements) of the vector.
Modifiable view of an array of type T.
constexpr Integer size() const noexcept
Returns the size of the array.
void add(ConstReferenceType val)
Adds element val to the end of the array.
Interface of a function call tracing service.
virtual StackTrace stackTrace(int first_function=0)=0
Character string indicating the call stack.
Interface of a profiler for message exchanges.
Definition IProfiler.h:32
virtual void enable(bool is_enabled)=0
Enables or disables statistics.
Information about the source of a message.
Abstraction interface for MPI operations. Primarily used to employ a decorator for MPI functions in o...
Iterator findRequest(MPI_Request request)
Checks that the request is in the list.
bool m_no_check_request
True if requests are not checked.
void _addRequest(MPI_Request request, const TraceInfo &trace_info)
MpiAdapter(ITraceMng *msg, IStat *stat, MPI_Comm comm, MpiLock *mpi_lock, IMpiProfiling *mpi_prof=nullptr)
void setRequestErrorAreFatal(bool v)
Indicates if errors in the list of requests are fatal.
int commRank() const
Rank of this instance in the communicator.
Definition MpiAdapter.h:148
void destroy()
Destroys the instance. It should no longer be used afterward.
MPI_Comm m_communicator
MPI Communicator.
Definition MpiAdapter.h:227
Request sendNonBlockingNoStat(const void *send_buffer, Int64 send_buffer_size, Int32 proc, MPI_Datatype data_type, int mpi_tag)
Non-blocking version of send without temporal statistics.
void setCheckRequest(bool v)
Indicates if requests are checked.
void setPrintRequestError(bool v)
Indicates if messages are displayed for errors in the requests.
Request buildRequest(int ret, MPI_Request request)
Constructs an Arccore request from an MPI request.
MPI_Request m_empty_request1
Empty requests. See MpiAdapter.cc for more information.
Definition MpiAdapter.h:235
void _addRequest(MPI_Request request)
Request receiveNonBlockingNoStat(void *recv_buffer, Int64 recv_buffer_size, Int32 source_rank, MPI_Datatype data_type, int mpi_tag)
Non-blocking version of receive without temporal statistics.
void _removeRequest(MPI_Request request)
Informative structure linked to the enumerations for MPI operations. Provides the name associated wit...
MPI implementation of the message exchange manager.
MPI specialization of a 'Request'.
Definition MpiRequest.h:36
Implementation of the MPI operations interface. Corresponds to a simple call to MPI functions of the ...
static MpiMessagePassingMng * create(MPI_Comm comm, bool clean_comm=false)
Creates a manager associated with the communicator comm.
Information for sending/receiving a point-to-point message.
Reference to an instance.
const String & toString() const
String indicating the call stack.
TraceAccessor(ITraceMng *m)
Constructs an accessor via the trace manager m.
TraceMessageDbg debug(Trace::eDebugLevel=Trace::Medium) const
Flow for a debug message.
TraceMessage info() const
Flow for an information message.
TraceMessage error() const
Flow for an error message.
TraceMessage warning() const
Flow for a warning message.
1D data vector with value semantics (STL style).
String getStackTrace()
Returns a string containing the call stack.
IStackTraceService * getStackTraceService()
Service used to obtain the call stack.
String getEnvironmentVariable(const String &name)
Environment variable named name.
std::int64_t Int64
Signed integer type of 64 bits.
Int32 Integer
Type representing an integer.
auto makeRef(InstanceType *t) -> Ref< InstanceType >
Creates a reference on a pointer.
bool arccoreIsCheck()
True if in check mode.
std::int32_t Int32
Signed integer type of 32 bits.