KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
mpi_data_communicator.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main author: Jordi Cotela
11 //
12 
13 #pragma once
14 
15 // System includes
16 #include <string>
17 #include <iostream>
18 #include <mpi.h>
19 
20 // External includes
21 
22 // Project includes
23 #include "includes/define.h"
25 
26 #ifndef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_SYNC_SHAPE_INTERFACE_FOR_TYPE
27 #define KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_SYNC_SHAPE_INTERFACE_FOR_TYPE(...) \
28 bool SynchronizeShape(__VA_ARGS__&) const override; \
29 bool SynchronizeShape( \
30  const __VA_ARGS__& rSendValue, const int SendDestination, const int SendTag, \
31  __VA_ARGS__& rRecvValue, const int RecvSource, const int RecvTag) const override; \
32 
33 #endif
34 
35 #ifndef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_REDUCE_INTERFACE_FOR_TYPE
36 #define KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_REDUCE_INTERFACE_FOR_TYPE(...) \
37 __VA_ARGS__ Sum(const __VA_ARGS__& rLocalValue, const int Root) const override; \
38 std::vector<__VA_ARGS__> Sum(const std::vector<__VA_ARGS__>& rLocalValues, const int Root) const override; \
39 void Sum(const std::vector<__VA_ARGS__>& rLocalValues, std::vector<__VA_ARGS__>& rGlobalValues, const int Root) const override; \
40 __VA_ARGS__ Min(const __VA_ARGS__& rLocalValue, const int Root) const override; \
41 std::vector<__VA_ARGS__> Min(const std::vector<__VA_ARGS__>& rLocalValues, const int Root) const override; \
42 void Min(const std::vector<__VA_ARGS__>& rLocalValues, std::vector<__VA_ARGS__>& rGlobalValues, const int Root) const override; \
43 __VA_ARGS__ Max(const __VA_ARGS__& rLocalValue, const int Root) const override; \
44 std::vector<__VA_ARGS__> Max(const std::vector<__VA_ARGS__>& rLocalValues, const int Root) const override; \
45 void Max(const std::vector<__VA_ARGS__>& rLocalValues, std::vector<__VA_ARGS__>& rGlobalValues, const int Root) const override; \
46 
47 #endif
48 
49 #ifndef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_ALLREDUCE_INTERFACE_FOR_TYPE
50 #define KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_ALLREDUCE_INTERFACE_FOR_TYPE(...) \
51 __VA_ARGS__ SumAll(const __VA_ARGS__& rLocalValue) const override; \
52 std::vector<__VA_ARGS__> SumAll(const std::vector<__VA_ARGS__>& rLocalValues) const override; \
53 void SumAll(const std::vector<__VA_ARGS__>& rLocalValues, std::vector<__VA_ARGS__>& rGlobalValues) const override; \
54 __VA_ARGS__ MinAll(const __VA_ARGS__& rLocalValue) const override; \
55 std::vector<__VA_ARGS__> MinAll(const std::vector<__VA_ARGS__>& rLocalValues) const override; \
56 void MinAll(const std::vector<__VA_ARGS__>& rLocalValues, std::vector<__VA_ARGS__>& rGlobalValues) const override; \
57 __VA_ARGS__ MaxAll(const __VA_ARGS__& rLocalValue) const override; \
58 std::vector<__VA_ARGS__> MaxAll(const std::vector<__VA_ARGS__>& rLocalValues) const override; \
59 void MaxAll(const std::vector<__VA_ARGS__>& rLocalValues, std::vector<__VA_ARGS__>& rGlobalValues) const override; \
60 
61 #endif
62 
63 #ifndef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_ALLREDUCE_LOC_INTERFACE_FOR_TYPE
64 #define KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_ALLREDUCE_LOC_INTERFACE_FOR_TYPE(...) \
65 std::pair<__VA_ARGS__, int> MinLocAll(const __VA_ARGS__& rLocalValue) const override; \
66 std::pair<__VA_ARGS__, int> MaxLocAll(const __VA_ARGS__& rLocalValue) const override;
67 #endif
68 
69 #ifndef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_SCANSUM_INTERFACE_FOR_TYPE
70 #define KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_SCANSUM_INTERFACE_FOR_TYPE(...) \
71 __VA_ARGS__ ScanSum(const __VA_ARGS__& rLocalValue) const override; \
72 std::vector<__VA_ARGS__> ScanSum(const std::vector<__VA_ARGS__>& rLocalValues) const override; \
73 void ScanSum(const std::vector<__VA_ARGS__>& rLocalValues, std::vector<__VA_ARGS__>& rGlobalValues) const override; \
74 
75 #endif
76 
77 #ifndef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_SENDRECV_INTERFACE_FOR_TYPE
78 #define KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_SENDRECV_INTERFACE_FOR_TYPE(...) \
79 __VA_ARGS__ SendRecvImpl( \
80  const __VA_ARGS__& SendValue, const int SendDestination, const int SendTag, \
81  const int RecvSource, const int RecvTag) const override; \
82 std::vector<__VA_ARGS__> SendRecvImpl(const std::vector<__VA_ARGS__>& rSendValues, \
83  const int SendDestination, const int SendTag, \
84  const int RecvSource, const int RecvTag) const override; \
85 void SendRecvImpl( \
86  const __VA_ARGS__& SendValue, const int SendDestination, const int SendTag, \
87  __VA_ARGS__& RecvValue, const int RecvSource, const int RecvTag) const override; \
88 void SendRecvImpl( \
89  const std::vector<__VA_ARGS__>& rSendValues, const int SendDestination, const int SendTag, \
90  std::vector<__VA_ARGS__>& rRecvValues, const int RecvSource, const int RecvTag) const override; \
91 void SendImpl(const __VA_ARGS__& rSendValues, \
92  const int SendDestination, const int SendTag = 0) const override; \
93 void SendImpl(const std::vector<__VA_ARGS__>& rSendValues, \
94  const int SendDestination, const int SendTag = 0) const override; \
95 void RecvImpl(__VA_ARGS__& rRecvValues, \
96  const int RecvSource, const int RecvTag = 0) const override; \
97 void RecvImpl(std::vector<__VA_ARGS__>& rRecvValues, \
98  const int RecvSource, const int RecvTag = 0) const override; \
99 
100 #endif
101 
102 #ifndef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_BROADCAST_INTERFACE_FOR_TYPE
103 #define KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_BROADCAST_INTERFACE_FOR_TYPE(...) \
104 void BroadcastImpl(__VA_ARGS__& rBuffer, const int SourceRank) const override; \
105 void BroadcastImpl(std::vector<__VA_ARGS__>& rBuffer, const int SourceRank) const override; \
106 
107 #endif
108 
109 #ifndef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_SCATTER_INTERFACE_FOR_TYPE
110 #define KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_SCATTER_INTERFACE_FOR_TYPE(...) \
111 std::vector<__VA_ARGS__> Scatter( \
112  const std::vector<__VA_ARGS__>& rSendValues, const int SourceRank) const override; \
113 void Scatter( \
114  const std::vector<__VA_ARGS__>& rSendValues, std::vector<__VA_ARGS__>& rRecvValues, \
115  const int SourceRank) const override; \
116 std::vector<__VA_ARGS__> Scatterv( \
117  const std::vector<std::vector<__VA_ARGS__>>& rSendValues, const int SourceRank) const override; \
118 void Scatterv( \
119  const std::vector<__VA_ARGS__>& rSendValues, \
120  const std::vector<int>& rSendCounts, const std::vector<int>& rSendOffsets, \
121  std::vector<__VA_ARGS__>& rRecvValues, const int SourceRank) const override; \
122 
123 #endif
124 
125 #ifndef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_GATHER_INTERFACE_FOR_TYPE
126 #define KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_GATHER_INTERFACE_FOR_TYPE(...) \
127 std::vector<__VA_ARGS__> Gather(const std::vector<__VA_ARGS__>& rSendValues, const int DestinationRank) const override; \
128 void Gather( \
129  const std::vector<__VA_ARGS__>& rSendValues, std::vector<__VA_ARGS__>& rRecvValues, \
130  const int DestinationRank) const override; \
131 std::vector<std::vector<__VA_ARGS__>> Gatherv( \
132  const std::vector<__VA_ARGS__>& rSendValues, const int DestinationRank) const override; \
133 void Gatherv(const std::vector<__VA_ARGS__>& rSendValues, \
134  std::vector<__VA_ARGS__>& rRecvValues, \
135  const std::vector<int>& rRecvCounts, \
136  const std::vector<int>& rRecvOffsets, \
137  const int DestinationRank) const override; \
138 std::vector<__VA_ARGS__> AllGather(const std::vector<__VA_ARGS__>& rSendValues) const override; \
139 void AllGather(const std::vector<__VA_ARGS__>& rSendValues, std::vector<__VA_ARGS__>& rRecvValues) const override; \
140 std::vector<std::vector<__VA_ARGS__>> AllGatherv(const std::vector<__VA_ARGS__>& rSendValues) const override; \
141 void AllGatherv(const std::vector<__VA_ARGS__>& rSendValues, std::vector<__VA_ARGS__>& rRecvValues, \
142  const std::vector<int>& rRecvCounts, const std::vector<int>& rRecvOffsets) const override;
143 #endif
144 
145 #ifndef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_PUBLIC_INTERFACE_FOR_TYPE
146 #define KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_PUBLIC_INTERFACE_FOR_TYPE(...) \
147 KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_REDUCE_INTERFACE_FOR_TYPE(__VA_ARGS__) \
148 KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_ALLREDUCE_INTERFACE_FOR_TYPE(__VA_ARGS__) \
149 KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_SCANSUM_INTERFACE_FOR_TYPE(__VA_ARGS__) \
150 KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_SCATTER_INTERFACE_FOR_TYPE(__VA_ARGS__) \
151 KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_GATHER_INTERFACE_FOR_TYPE(__VA_ARGS__) \
152 KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_SYNC_SHAPE_INTERFACE_FOR_TYPE(__VA_ARGS__) \
153 
154 #endif
155 
156 #ifndef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_IMPLEMENTATION_FOR_TYPE
157 #define KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_IMPLEMENTATION_FOR_TYPE(...) \
158 KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_SENDRECV_INTERFACE_FOR_TYPE(__VA_ARGS__) \
159 KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_BROADCAST_INTERFACE_FOR_TYPE(__VA_ARGS__) \
160 
161 #endif
162 
163 namespace Kratos
164 {
167 
170 
172 
181 class KRATOS_API(KRATOS_MPI_CORE) MPIDataCommunicator: public DataCommunicator
182 {
183  public:
186 
189 
193 
195  explicit MPIDataCommunicator(MPI_Comm MPIComm);
196 
198  ~MPIDataCommunicator() override;
199 
203 
205 
210  static MPIDataCommunicator::UniquePointer Create(MPI_Comm MPIComm);
211 
212  void Barrier() const override;
213 
225 
226  // MinLoc and MaxLoc AllReduce operations
232 
233  // Reduce operations
234 
235  bool AndReduce(
236  const bool Value,
237  const int Root) const override;
238 
239  Kratos::Flags AndReduce(
240  const Kratos::Flags Values,
241  const Kratos::Flags Mask,
242  const int Root) const override;
243 
244  bool OrReduce(
245  const bool Value,
246  const int Root) const override;
247 
248  Kratos::Flags OrReduce(
249  const Kratos::Flags Values,
250  const Kratos::Flags Mask,
251  const int Root) const override;
252 
253  // Allreduce operations
254 
255  bool AndReduceAll(const bool Value) const override;
256 
257  Kratos::Flags AndReduceAll(const Kratos::Flags Values, const Kratos::Flags Mask) const override;
258 
259  bool OrReduceAll(const bool Value) const override;
260 
261  Kratos::Flags OrReduceAll(const Kratos::Flags Values, const Kratos::Flags Mask) const override;
262 
266 
268 
271  static MPI_Comm GetMPICommunicator(const DataCommunicator& rDataCommunicator);
272 
276 
282  int Rank() const override;
283 
289  int Size() const override;
290 
295  bool IsDistributed() const override;
296 
302  bool IsDefinedOnThisRank() const override;
303 
309  bool IsNullOnThisRank() const override;
310 
318  const DataCommunicator& GetSubDataCommunicator(
319  const std::vector<int>& rRanks,
320  const std::string& rNewCommunicatorName
321  ) const override;
322 
326 
327  bool BroadcastErrorIfTrue(bool Condition, const int SourceRank) const override;
328 
329  bool BroadcastErrorIfFalse(bool Condition, const int SourceRank) const override;
330 
331  bool ErrorIfTrueOnAnyRank(bool Condition) const override;
332 
333  bool ErrorIfFalseOnAnyRank(bool Condition) const override;
334 
338 
340  std::string Info() const override;
341 
343  void PrintInfo(std::ostream &rOStream) const override;
344 
346  void PrintData(std::ostream &rOStream) const override;
347 
349 
350  protected:
351 
363 
364  // Broadcast operations
365 
366  void BroadcastImpl(std::string& rBroadcastValues, const int SourceRank) const override;
367 
368  // Sendrecv operations
369 
370  std::string SendRecvImpl(
371  const std::string& rSendValues, const int SendDestination, const int SendTag,
372  const int RecvSource, const int RecvTag) const override;
373 
374  void SendRecvImpl(
375  const std::string& rSendValues, const int SendDestination, const int SendTag,
376  std::string& rRecvValues, const int RecvSource, const int RecvTag) const override;
377 
378  void SendImpl(const std::string& rSendValues, const int SendDestination, const int SendTag = 0) const override;
379 
380  void RecvImpl(std::string& rRecvValues, const int RecvSource, const int RecvTag = 0) const override;
381 
382  private:
385 
386  MPI_Comm mComm;
387 
391 
392  void CheckMPIErrorCode(const int ierr, const std::string& MPICallName) const;
393 
394  template<class TDataType> bool SynchronizeShapeDetail(TDataType& rValue) const;
395 
396  template<class TDataType> bool SynchronizeShapeDetail(
397  const TDataType& rSendValue,
398  const int SendDestination,
399  const int SendTag,
400  TDataType& rRecvValue,
401  const int RecvSource,
402  const int RecvTag) const;
403 
404  template<class TDataType> void ReduceDetail(
405  const TDataType& rLocalValues,
406  TDataType& rReducedValues,
407  MPI_Op Operation,
408  const int Root) const;
409 
410  template<class TDataType> TDataType ReduceDetail(
411  const TDataType& rLocalValues,
412  MPI_Op Operation,
413  const int Root) const;
414 
415  template<class TDataType> std::vector<TDataType> ReduceDetailVector(
416  const std::vector<TDataType>& rLocalValues,
417  MPI_Op Operation,
418  const int Root) const;
419 
420  template<class TDataType> void AllReduceDetail(
421  const TDataType& rLocalValues,
422  TDataType& rReducedValues,
423  MPI_Op Operation) const;
424 
425  template<class TDataType> TDataType AllReduceDetail(
426  const TDataType& rLocalValues, MPI_Op Operation) const;
427 
428  template<class TDataType> std::vector<TDataType> AllReduceDetailVector(
429  const std::vector<TDataType>& rLocalValues,
430  MPI_Op Operation) const;
431 
440  template<class TDataType>
441  std::pair<TDataType, int> AllReduceDetailWithLocation(
442  const std::pair<TDataType, int>& rLocalValues,
443  MPI_Op Operation
444  ) const;
445 
446  template<class TDataType> void ScanDetail(
447  const TDataType& rLocalValues,
448  TDataType& rReducedValues,
449  MPI_Op Operation) const;
450 
451  template<class TDataType> TDataType ScanDetail(
452  const TDataType rLocalValues,
453  MPI_Op Operation) const;
454 
455  template<class TDataType> std::vector<TDataType> ScanDetail(
456  const std::vector<TDataType>& rLocalValues,
457  MPI_Op Operation) const;
458 
459  template<class TDataType> void SendRecvDetail(
460  const TDataType& rSendMessage, const int SendDestination, const int SendTag,
461  TDataType& rRecvMessage, const int RecvSource, const int RecvTag) const;
462 
463  template<class TDataType> TDataType SendRecvDetail(
464  const TDataType& rSendMessage,
465  const int SendDestination, const int SendTag,
466  const int RecvSource, const int RecvTag) const;
467 
468  template<class TDataType> std::vector<TDataType> SendRecvDetail(
469  const std::vector<TDataType>& rSendMessage,
470  const int SendDestination, const int SendTag,
471  const int RecvSource, const int RecvTag) const;
472 
473  template<class TDataType> void SendDetail(
474  const TDataType& rSendValues, const int SendDestination, const int SendTag) const;
475 
476  template<class TDataType> void RecvDetail(
477  TDataType& rRecvValues, const int RecvSource, const int RecvTag) const;
478 
479  template<class TDataType> void BroadcastDetail(
480  TDataType& rBuffer, const int SourceRank) const;
481 
482  template<class TSendDataType, class TRecvDataType> void ScatterDetail(
483  const TSendDataType& rSendValues, TRecvDataType& rRecvValues, const int SourceRank) const;
484 
485  template<class TDataType> std::vector<TDataType> ScatterDetail(
486  const std::vector<TDataType>& rSendValues, const int SourceRank) const;
487 
488  template<class TDataType> void ScattervDetail(
489  const TDataType& rSendValues,
490  const std::vector<int>& rSendCounts, const std::vector<int>& rSendOffsets,
491  TDataType& rRecvValues, const int SourceRank) const;
492 
493  template<class TDataType> std::vector<TDataType> ScattervDetail(
494  const std::vector<std::vector<TDataType>>& rSendValues,const int SourceRank) const;
495 
496  template<class TSendDataType, class TRecvDataType> void GatherDetail(
497  const TSendDataType& rSendValues, TRecvDataType& rRecvValues, const int RecvRank) const;
498 
499  template<class TDataType> std::vector<TDataType> GatherDetail(
500  const std::vector<TDataType>& rSendValues, const int DestinationRank) const;
501 
502  template<class TDataType> void GathervDetail(
503  const TDataType& rSendValues, TDataType& rRecvValues,
504  const std::vector<int>& rRecvCounts, const std::vector<int>& rRecvOffsets,
505  const int RecvRank) const;
506 
507  template<class TDataType> std::vector<std::vector<TDataType>>
508  GathervDetail(const std::vector<TDataType>& rSendValues, const int DestinationRank) const;
509 
510  template<class TDataType> void AllGatherDetail(
511  const TDataType& rSendValues, TDataType& rRecvValues) const;
512 
513  template<class TDataType> std::vector<TDataType> AllGatherDetail(
514  const std::vector<TDataType>& rSendValues) const;
515 
516  template<class TDataType>
517  void AllGathervDetail(
518  const TDataType& rSendValues, TDataType& rRecvValues,
519  const std::vector<int>& rRecvCounts, const std::vector<int>& rRecvOffsets) const;
520 
521  template<class TDataType>
522  std::vector<std::vector<TDataType>> AllGathervDetail(
523  const std::vector<TDataType>& rSendValues) const;
524 
525  bool IsEqualOnAllRanks(const int LocalValue) const;
526 
527  bool IsValidRank(const int Rank) const;
528 
529  template<class TDataType> void ValidateScattervInput(
530  const TDataType& rSendValues,
531  const std::vector<int>& rSendCounts, const std::vector<int>& rSendOffsets,
532  TDataType& rRecvValues, const int SourceRank) const;
533 
534  template<class TDataType> void ValidateGathervInput(
535  const TDataType& rSendValues, TDataType& rRecvValues,
536  const std::vector<int>& rRecvCounts, const std::vector<int>& rRecvOffsets,
537  const int RecvRank) const;
538 
539  template<class TDataType>
540  void ValidateAllGathervInput(
541  const TDataType& rSendValues, TDataType& rRecvValues,
542  const std::vector<int>& rRecvCounts, const std::vector<int>& rRecvOffsets) const;
543 
544  template<class TDataType> void PrepareScattervBuffers(
545  const std::vector<std::vector<TDataType>>& rInputMessage,
546  std::vector<TDataType>& rScattervMessage,
547  std::vector<int>& rMessageLengths,
548  std::vector<int>& rMessageDistances,
549  std::vector<TDataType>& rResult,
550  const int SourceRank) const;
551 
552  template<class TDataType> void PrepareGathervBuffers(
553  const std::vector<TDataType>& rGathervInput,
554  std::vector<TDataType>& rGathervMessage,
555  std::vector<int>& rMessageLengths,
556  std::vector<int>& rMessageDistances,
557  const int DestinationRank) const;
558 
559  template<class TDataType>
560  void PrepareAllGathervBuffers(
561  const std::vector<TDataType>& rGathervInput,
562  std::vector<TDataType>& rGathervMessage,
563  std::vector<int>& rMessageLengths,
564  std::vector<int>& rMessageDistances) const;
565 
566  template<class TDataType> void PrepareGathervReturn(
567  const std::vector<TDataType>& rGathervMessage,
568  const std::vector<int>& rMessageLengths,
569  const std::vector<int>& rMessageDistances,
570  std::vector<std::vector<TDataType>>& rOutputMessage,
571  const int DestinationRank) const;
572 
573  template<class TDataType>
574  void PrepareAllGathervReturn(
575  const std::vector<TDataType>& rGathervMessage,
576  const std::vector<int>& rMessageLengths,
577  const std::vector<int>& rMessageDistances,
578  std::vector<std::vector<TDataType>>& rOutputMessage) const;
579 
580  template<class TValue> inline MPI_Datatype MPIDatatype(const TValue&) const;
581 
582  template<class TContainer> inline void* MPIBuffer(TContainer& rValues) const;
583 
584  template<class TContainer> inline const void* MPIBuffer(const TContainer& rValues) const;
585 
586  template<class TContainer> inline int MPIMessageSize(const TContainer& rValues) const;
587 
591 
593  MPIDataCommunicator(MPIDataCommunicator const &rOther) = delete;
594 
596  MPIDataCommunicator &operator=(MPIDataCommunicator const &rOther) = delete;
597 
599 
600 }; // Class MPIDataCommunicator
601 
603 
606 
608 inline std::istream &operator>>(std::istream &rIStream,
609  MPIDataCommunicator &rThis)
610 {
611  return rIStream;
612 }
613 
615 inline std::ostream &operator<<(std::ostream &rOStream,
616  const MPIDataCommunicator &rThis)
617 {
618  rThis.PrintInfo(rOStream);
619  rOStream << std::endl;
620  rThis.PrintData(rOStream);
621 
622  return rOStream;
623 }
625 
627 
628 } // namespace Kratos.
629 
630 #undef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_REDUCE_INTERFACE_FOR_TYPE
631 #undef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_ALLREDUCE_INTERFACE_FOR_TYPE
632 #undef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_SCANSUM_INTERFACE_FOR_TYPE
633 #undef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_SENDRECV_INTERFACE_FOR_TYPE
634 #undef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_BROADCAST_INTERFACE_FOR_TYPE
635 #undef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_SCATTER_INTERFACE_FOR_TYPE
636 #undef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_GATHER_INTERFACE_FOR_TYPE
637 #undef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_PUBLIC_INTERFACE_FOR_TYPE
638 #undef KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_IMPLEMENTATION_FOR_TYPE
PeriodicInterfaceProcess & operator=(const PeriodicInterfaceProcess &)=delete
std::string Info() const override
Turn back information as a string.
Definition: periodic_interface_process.hpp:93
Base class for all Conditions.
Definition: condition.h:59
Serial (do-nothing) version of a wrapper class for MPI communication.
Definition: data_communicator.h:318
Definition: flags.h:58
Wrapper for common MPI calls within Kratos.
Definition: mpi_data_communicator.h:182
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: mpi_data_communicator.cpp:512
KRATOS_CLASS_POINTER_DEFINITION(MPIDataCommunicator)
Pointer definition of MPIDataCommunicator.
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: mpi_data_communicator.cpp:517
The base class for all operations in Kratos.
Definition: operation.h:43
#define KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_ALLREDUCE_LOC_INTERFACE_FOR_TYPE(...)
Definition: mpi_data_communicator.h:64
#define KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_PUBLIC_INTERFACE_FOR_TYPE(...)
Definition: mpi_data_communicator.h:146
#define KRATOS_MPI_DATA_COMMUNICATOR_DECLARE_IMPLEMENTATION_FOR_TYPE(...)
Definition: mpi_data_communicator.h:157
Modeler::Pointer Create(const std::string &ModelerName, Model &rModel, const Parameters ModelParameters)
Checks if the modeler is registered.
Definition: modeler_factory.cpp:30
TSpaceType::IndexType Size(TSpaceType &dummy, typename TSpaceType::VectorType const &rV)
Definition: add_strategies_to_python.cpp:111
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432