00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027 #if defined(BL_OLD_STL)
00028 #include <stdio.h>
00029 #else
00030 #include <cstdio>
00031 #endif
00032
00033 #include <Utility.H>
00034 #include <ParallelDescriptor.H>
00035 #include <ParmParse.H>
00036
00037 FabComTag::FabComTag ()
00038 {
00039 fromProc = 0;
00040 toProc = 0;
00041 fabIndex = 0;
00042 fineIndex = 0;
00043 srcComp = 0;
00044 destComp = 0;
00045 nComp = 0;
00046 face = 0;
00047 fabArrayId = 0;
00048 fillBoxId = 0;
00049 procThatNeedsData = 0;
00050 procThatHasData = 0;
00051 }
00052
00053
00054
00055 namespace ParallelDescriptor
00056 {
00057
00058
00059
00060 int m_MyId = -1;
00061
00062
00063
00064 int m_nProcs = -1;
00065
00066
00067
00068 int m_nProcsCFD = -1;
00069
00070
00071
00072 MPI_Comm m_comm;
00073
00074 const int ioProcessor = 0;
00075
00076 namespace util
00077 {
00078
00079
00080
00081 void DoAllReduceReal (Real& r, MPI_Op op);
00082 void DoAllReduceLong (long& r, MPI_Op op);
00083 void DoAllReduceInt (int& r, MPI_Op op);
00084 void DoReduceReal (Real& r, MPI_Op op, int cpu);
00085 void DoReduceLong (long& r, MPI_Op op, int cpu);
00086 void DoReduceInt (int& r, MPI_Op op, int cpu);
00087
00088
00089
00090 void SetNProcsCFD ();
00091 }
00092 }
00093
00094 int
00095 ParallelDescriptor::MyProc ()
00096 {
00097 BL_ASSERT(m_MyId != -1);
00098 return m_MyId;
00099 }
00100
00101 int
00102 ParallelDescriptor::NProcs ()
00103 {
00104 BL_ASSERT(m_nProcs != -1);
00105 return m_nProcs;
00106 }
00107
00108 int
00109 ParallelDescriptor::NProcsCFD ()
00110 {
00111 if (m_nProcsCFD == -1)
00112 util::SetNProcsCFD();
00113
00114 BL_ASSERT(m_nProcsCFD != -1);
00115
00116 return m_nProcsCFD;
00117 }
00118
00119 bool
00120 ParallelDescriptor::IOProcessor ()
00121 {
00122 BL_ASSERT(m_MyId != -1);
00123 return m_MyId == ioProcessor;
00124 }
00125
00126 int
00127 ParallelDescriptor::IOProcessorNumber ()
00128 {
00129 return ioProcessor;
00130 }
00131
00132 MPI_Comm
00133 ParallelDescriptor::Communicator ()
00134 {
00135 return m_comm;
00136 }
00137
00138 CommData::CommData ()
00139 {
00140 for (int i = 0; i < length(); i++)
00141 m_data[i] = 0;
00142 }
00143
00144 CommData::CommData (int face,
00145 int fabindex,
00146 int fromproc,
00147 int id,
00148 int ncomp,
00149 int srccomp,
00150 int fabarrayid,
00151 const Box& box)
00152 {
00153 m_data[0] = face;
00154 m_data[1] = fabindex;
00155 m_data[2] = fromproc;
00156 m_data[3] = id;
00157 m_data[4] = ncomp;
00158 m_data[5] = srccomp;
00159 m_data[6] = fabarrayid;
00160
00161 const IntVect& sm = box.smallEnd();
00162
00163 for (int i = 0; i < BL_SPACEDIM; i++)
00164 m_data[7+i] = sm[i];
00165
00166 const IntVect& bg = box.bigEnd();
00167
00168 for (int i = 0; i < BL_SPACEDIM; i++)
00169 m_data[7+BL_SPACEDIM+i] = bg[i];
00170
00171 IntVect typ = box.type();
00172
00173 for (int i = 0; i < BL_SPACEDIM; i++)
00174 m_data[7+2*BL_SPACEDIM+i] = typ[i];
00175 }
00176
00177 CommData::CommData (const CommData& rhs)
00178 {
00179 for (int i = 0; i < length(); i++)
00180 m_data[i] = rhs.m_data[i];
00181 }
00182
00183 CommData&
00184 CommData::operator= (const CommData& rhs)
00185 {
00186 if (!(this == &rhs))
00187 {
00188 for (int i = 0; i < length(); i++)
00189 m_data[i] = rhs.m_data[i];
00190 }
00191 return *this;
00192 }
00193
00194 bool
00195 CommData::operator== (const CommData& rhs) const
00196 {
00197 for (int i = 0; i < length(); i++)
00198 if (!(m_data[i] == rhs.m_data[i]))
00199 return false;
00200
00201 return true;
00202 }
00203
00204 std::ostream&
00205 operator<< (std::ostream& os,
00206 const CommData& cd)
00207 {
00208 os << cd.face() << ' '
00209 << cd.fabindex() << ' '
00210 << cd.fromproc() << ' '
00211 << cd.id() << ' '
00212 << cd.nComp() << ' '
00213 << cd.srcComp() << ' '
00214 << cd.fabarrayid() << ' '
00215 << cd.box();
00216
00217 return os;
00218 }
00219
00220 std::ostream&
00221 operator<< (std::ostream& os,
00222 const Array<CommData>& cd)
00223 {
00224 for (int i = 0; i < cd.size(); i++)
00225 os << cd[i] << '\n';
00226 return os;
00227 }
00228
00229 #ifdef BL_USE_MPI
00230
00231 #include <ccse-mpi.H>
00232
00233 namespace
00234 {
00235 const char*
00236 the_message_string (const char* file,
00237 int line,
00238 const char* call,
00239 int status)
00240 {
00241
00242
00243
00244 const int DIM = 1024;
00245 static char buf[DIM];
00246 if ( status )
00247 {
00248 std::sprintf(buf, "BoxLib MPI Error: File %s, line %d, %s: %s",
00249 file, line, call, ParallelDescriptor::ErrorString(status));
00250 }
00251 else
00252 {
00253 std::sprintf(buf, "BoxLib MPI Error: File %s, line %d, %s",
00254 file, line, call);
00255 }
00256 buf[DIM-1] = '\0';
00257 return buf;
00258 }
00259 }
00260
00261 namespace ParallelDescriptor
00262 {
00263 void
00264 MPI_Error (const char* file, int line, const char* str, int rc)
00265 {
00266 BoxLib::Error(the_message_string(file, line, str, rc));
00267 }
00268 }
00269
00270
00271 std::vector<size_t>
00272 ParallelDescriptor::Gather (const size_t& t, int root)
00273 {
00274 std::vector<size_t> resl;
00275 if ( root == MyProc() ) resl.resize(NProcs());
00276 BL_MPI_REQUIRE( MPI_Gather(const_cast<size_t*>(&t),
00277 1,
00278 Mpi_typemap<size_t>::type(),
00279 &resl[0],
00280 1,
00281 Mpi_typemap<size_t>::type(),
00282 root,
00283 Communicator()) );
00284 return resl;
00285 }
00286
00287 ParallelDescriptor::Message
00288 ParallelDescriptor::Send (const std::vector<size_t>& buf,
00289 int dst_pid,
00290 int tag)
00291 {
00292 BL_MPI_REQUIRE( MPI_Send(const_cast<size_t*>(&buf[0]),
00293 buf.size(),
00294 Mpi_typemap<size_t>::type(),
00295 dst_pid,
00296 tag,
00297 Communicator()) );
00298 return Message();
00299 }
00300
00301 ParallelDescriptor::Message
00302 ParallelDescriptor::Recv (std::vector<size_t>& buf,
00303 int src_pid,
00304 int tag)
00305 {
00306 MPI_Status stat;
00307 BL_MPI_REQUIRE( MPI_Recv(&buf[0],
00308 buf.size(),
00309 Mpi_typemap<size_t>::type(),
00310 src_pid,
00311 tag,
00312 Communicator(),
00313 &stat) );
00314 return Message(stat, Mpi_typemap<size_t>::type());
00315 }
00316
00317 ParallelDescriptor::Message
00318 ParallelDescriptor::Recv (std::vector<char>& buf,
00319 int src_pid,
00320 int tag)
00321 {
00322 MPI_Status stat;
00323 BL_MPI_REQUIRE( MPI_Recv(&buf[0],
00324 buf.size(),
00325 Mpi_typemap<char>::type(),
00326 src_pid,
00327 tag,
00328 Communicator(),
00329 &stat) );
00330 return Message(stat, Mpi_typemap<char>::type());
00331 }
00332
00333
00334 void
00335 ParallelDescriptor::Abort ()
00336 {
00337 #ifdef WIN32
00338 throw;
00339 #endif
00340 MPI_Abort(Communicator(), -1);
00341 }
00342
00343 void
00344 ParallelDescriptor::Abort (int errorcode)
00345 {
00346 BoxLib::Abort(ErrorString(errorcode));
00347 }
00348
00349 const char*
00350 ParallelDescriptor::ErrorString (int errorcode)
00351 {
00352 BL_ASSERT(errorcode > 0 && errorcode <= MPI_ERR_LASTCODE);
00353
00354 int len = 0;
00355
00356 static char msg[MPI_MAX_ERROR_STRING+1];
00357
00358 MPI_Error_string(errorcode, msg, &len);
00359
00360 BL_ASSERT(len <= MPI_MAX_ERROR_STRING);
00361
00362 return msg;
00363 }
00364
00365 ParallelDescriptor::Message::Message ()
00366 :
00367 m_finished(true),
00368 m_type(MPI_DATATYPE_NULL),
00369 m_req(MPI_REQUEST_NULL)
00370 {}
00371
00372 ParallelDescriptor::Message::Message (MPI_Request req_, MPI_Datatype type_)
00373 :
00374 m_finished(false),
00375 m_type(type_),
00376 m_req(req_)
00377 {}
00378
00379 ParallelDescriptor::Message::Message (MPI_Status stat_, MPI_Datatype type_)
00380 :
00381 m_finished(true),
00382 m_type(type_),
00383 m_req(MPI_REQUEST_NULL), m_stat(stat_)
00384 {}
00385
00386 void
00387 ParallelDescriptor::Message::wait ()
00388 {
00389 BL_MPI_REQUIRE( MPI_Wait(&m_req, &m_stat) );
00390 }
00391
00392 bool
00393 ParallelDescriptor::Message::test ()
00394 {
00395 int flag;
00396 BL_MPI_REQUIRE( MPI_Test(&m_req, &flag, &m_stat) );
00397 m_finished = flag != 0;
00398 return m_finished;
00399 }
00400
00401 int
00402 ParallelDescriptor::Message::tag () const
00403 {
00404 if ( !m_finished ) BoxLib::Error("Message::tag: Not Finished!");
00405 return m_stat.MPI_TAG;
00406 }
00407
00408 int
00409 ParallelDescriptor::Message::pid () const
00410 {
00411 if ( !m_finished ) BoxLib::Error("Message::pid: Not Finished!");
00412 return m_stat.MPI_SOURCE;
00413 }
00414
00415 size_t
00416 ParallelDescriptor::Message::count () const
00417 {
00418 if ( m_type == MPI_DATATYPE_NULL ) BoxLib::Error("Message::count: Bad Type!");
00419 if ( !m_finished ) BoxLib::Error("Message::count: Not Finished!");
00420 int cnt;
00421 BL_MPI_REQUIRE( MPI_Get_count(&m_stat, m_type, &cnt) );
00422 return cnt;
00423 }
00424
00425 MPI_Datatype
00426 ParallelDescriptor::Message::type () const
00427 {
00428 return m_type;
00429 }
00430
00431 MPI_Request
00432 ParallelDescriptor::Message::req () const
00433 {
00434 return m_req;
00435 }
00436
00437 void
00438 ParallelDescriptor::util::SetNProcsCFD ()
00439 {
00440 BL_ASSERT(m_nProcs != -1);
00441 BL_ASSERT(m_nProcsCFD == -1);
00442
00443 m_nProcsCFD = m_nProcs;
00444
00445 ParmParse pp("ParallelDescriptor");
00446
00447 if (pp.query("nProcsCFD",m_nProcsCFD))
00448 {
00449 if (!(m_nProcsCFD > 0 && m_nProcsCFD <= m_nProcs))
00450 m_nProcsCFD = m_nProcs;
00451
00452 if (ParallelDescriptor::IOProcessor())
00453 std::cout << "--> Running job with NProcsCFD = "
00454 << m_nProcsCFD
00455 << std::endl;
00456 }
00457 }
00458
00459 void
00460 ParallelDescriptor::StartParallel (int* argc,
00461 char*** argv)
00462 {
00463 BL_ASSERT(m_MyId == -1);
00464 BL_ASSERT(m_nProcs == -1);
00465 BL_ASSERT(m_nProcsCFD == -1);
00466
00467 m_comm = MPI_COMM_WORLD;
00468
00469 int sflag;
00470
00471 BL_MPI_REQUIRE( MPI_Initialized(&sflag) );
00472
00473 if (!sflag)
00474 BL_MPI_REQUIRE( MPI_Init(argc, argv) );
00475 BL_MPI_REQUIRE( MPI_Barrier(Communicator()) );
00476
00477 BL_MPI_REQUIRE( MPI_Comm_size(Communicator(), &m_nProcs) );
00478
00479 BL_MPI_REQUIRE( MPI_Comm_rank(Communicator(), &m_MyId) );
00480
00481
00482
00483
00484 BL_MPI_REQUIRE( MPI_Barrier(Communicator()) );
00485 }
00486
00487 void
00488 ParallelDescriptor::EndParallel ()
00489 {
00490 BL_ASSERT(m_MyId != -1);
00491 BL_ASSERT(m_nProcs != -1);
00492
00493 BL_MPI_REQUIRE( MPI_Finalize() );
00494 }
00495
00496 double
00497 ParallelDescriptor::second ()
00498 {
00499 return MPI_Wtime();
00500 }
00501
00502 void
00503 ParallelDescriptor::Barrier ()
00504 {
00505 BL_MPI_REQUIRE( MPI_Barrier(ParallelDescriptor::Communicator()) );
00506 }
00507
00508 void
00509 ParallelDescriptor::Barrier (MPI_Comm comm)
00510 {
00511 BL_MPI_REQUIRE( MPI_Barrier(comm) );
00512 }
00513
00514 void
00515 ParallelDescriptor::Test (MPI_Request& request, int& flag, MPI_Status& status)
00516 {
00517 BL_MPI_REQUIRE( MPI_Test(&request,&flag,&status) );
00518 }
00519
00520 void
00521 ParallelDescriptor::Comm_dup (MPI_Comm comm, MPI_Comm& newcomm)
00522 {
00523 BL_MPI_REQUIRE( MPI_Comm_dup(comm, &newcomm) );
00524 }
00525
00526 void
00527 ParallelDescriptor::util::DoAllReduceReal (Real& r,
00528 MPI_Op op)
00529 {
00530 Real recv;
00531
00532 BL_MPI_REQUIRE( MPI_Allreduce(&r,
00533 &recv,
00534 1,
00535 Mpi_typemap<Real>::type(),
00536 op,
00537 Communicator()) );
00538 r = recv;
00539 }
00540
00541 void
00542 ParallelDescriptor::util::DoReduceReal (Real& r,
00543 MPI_Op op,
00544 int cpu)
00545 {
00546 Real recv;
00547
00548 BL_MPI_REQUIRE( MPI_Reduce(&r,
00549 &recv,
00550 1,
00551 Mpi_typemap<Real>::type(),
00552 op,
00553 cpu,
00554 Communicator()) );
00555
00556 if (ParallelDescriptor::MyProc() == cpu)
00557 r = recv;
00558 }
00559
00560 void
00561 ParallelDescriptor::ReduceRealMax (Real& r)
00562 {
00563 util::DoAllReduceReal(r,MPI_MAX);
00564 }
00565
00566 void
00567 ParallelDescriptor::ReduceRealMin (Real& r)
00568 {
00569 util::DoAllReduceReal(r,MPI_MIN);
00570 }
00571
00572 void
00573 ParallelDescriptor::ReduceRealSum (Real& r)
00574 {
00575 util::DoAllReduceReal(r,MPI_SUM);
00576 }
00577
00578 void
00579 ParallelDescriptor::ReduceRealMax (Real& r, int cpu)
00580 {
00581 util::DoReduceReal(r,MPI_MAX,cpu);
00582 }
00583
00584 void
00585 ParallelDescriptor::ReduceRealMin (Real& r, int cpu)
00586 {
00587 util::DoReduceReal(r,MPI_MIN,cpu);
00588 }
00589
00590 void
00591 ParallelDescriptor::ReduceRealSum (Real& r, int cpu)
00592 {
00593 util::DoReduceReal(r,MPI_SUM,cpu);
00594 }
00595
00596 void
00597 ParallelDescriptor::util::DoAllReduceLong (long& r,
00598 MPI_Op op)
00599 {
00600 long recv;
00601
00602 BL_MPI_REQUIRE( MPI_Allreduce(&r,
00603 &recv,
00604 1,
00605 MPI_LONG,
00606 op,
00607 Communicator()) );
00608 r = recv;
00609 }
00610
00611 void
00612 ParallelDescriptor::util::DoReduceLong (long& r,
00613 MPI_Op op,
00614 int cpu)
00615 {
00616 long recv;
00617
00618 BL_MPI_REQUIRE( MPI_Reduce(&r,
00619 &recv,
00620 1,
00621 MPI_LONG,
00622 op,
00623 cpu,
00624 Communicator()));
00625
00626 if (ParallelDescriptor::MyProc() == cpu)
00627 r = recv;
00628 }
00629
00630 void
00631 ParallelDescriptor::ReduceLongAnd (long& r)
00632 {
00633 util::DoAllReduceLong(r,MPI_LAND);
00634 }
00635
00636 void
00637 ParallelDescriptor::ReduceLongSum (long& r)
00638 {
00639 util::DoAllReduceLong(r,MPI_SUM);
00640 }
00641
00642 void
00643 ParallelDescriptor::ReduceLongMax (long& r)
00644 {
00645 util::DoAllReduceLong(r,MPI_MAX);
00646 }
00647
00648 void
00649 ParallelDescriptor::ReduceLongMin (long& r)
00650 {
00651 util::DoAllReduceLong(r,MPI_MIN);
00652 }
00653
00654 void
00655 ParallelDescriptor::ReduceLongAnd (long& r, int cpu)
00656 {
00657 util::DoReduceLong(r,MPI_LAND,cpu);
00658 }
00659
00660 void
00661 ParallelDescriptor::ReduceLongSum (long& r, int cpu)
00662 {
00663 util::DoReduceLong(r,MPI_SUM,cpu);
00664 }
00665
00666 void
00667 ParallelDescriptor::ReduceLongMax (long& r, int cpu)
00668 {
00669 util::DoReduceLong(r,MPI_MAX,cpu);
00670 }
00671
00672 void
00673 ParallelDescriptor::ReduceLongMin (long& r, int cpu)
00674 {
00675 util::DoReduceLong(r,MPI_MIN,cpu);
00676 }
00677
00678 void
00679 ParallelDescriptor::util::DoAllReduceInt (int& r,
00680 MPI_Op op)
00681 {
00682 int recv;
00683
00684 BL_MPI_REQUIRE( MPI_Allreduce(&r,
00685 &recv,
00686 1,
00687 MPI_INT,
00688 op,
00689 Communicator()));
00690 r = recv;
00691 }
00692
00693 void
00694 ParallelDescriptor::util::DoReduceInt (int& r,
00695 MPI_Op op,
00696 int cpu)
00697 {
00698 int recv;
00699
00700 BL_MPI_REQUIRE( MPI_Reduce(&r,
00701 &recv,
00702 1,
00703 MPI_INT,
00704 op,
00705 cpu,
00706 Communicator()));
00707
00708 if (ParallelDescriptor::MyProc() == cpu)
00709 r = recv;
00710 }
00711
00712 void
00713 ParallelDescriptor::ReduceIntSum (int& r)
00714 {
00715 util::DoAllReduceInt(r,MPI_SUM);
00716 }
00717
00718 void
00719 ParallelDescriptor::ReduceIntMax (int& r)
00720 {
00721 util::DoAllReduceInt(r,MPI_MAX);
00722 }
00723
00724 void
00725 ParallelDescriptor::ReduceIntMin (int& r)
00726 {
00727 util::DoAllReduceInt(r,MPI_MIN);
00728 }
00729
00730 void
00731 ParallelDescriptor::ReduceIntSum (int& r, int cpu)
00732 {
00733 util::DoReduceInt(r,MPI_SUM,cpu);
00734 }
00735
00736 void
00737 ParallelDescriptor::ReduceIntMax (int& r, int cpu)
00738 {
00739 util::DoReduceInt(r,MPI_MAX,cpu);
00740 }
00741
00742 void
00743 ParallelDescriptor::ReduceIntMin (int& r, int cpu)
00744 {
00745 util::DoReduceInt(r,MPI_MIN,cpu);
00746 }
00747
00748 void
00749 ParallelDescriptor::ReduceBoolAnd (bool& r)
00750 {
00751 int src = r;
00752
00753 util::DoAllReduceInt(src,MPI_SUM);
00754
00755 r = (src == ParallelDescriptor::NProcs()) ? true : false;
00756 }
00757
00758 void
00759 ParallelDescriptor::ReduceBoolOr (bool& r)
00760 {
00761 int src = r;
00762
00763 util::DoAllReduceInt(src,MPI_SUM);
00764
00765 r = (src == 0) ? false : true;
00766 }
00767
00768 void
00769 ParallelDescriptor::ReduceBoolAnd (bool& r, int cpu)
00770 {
00771 int src = r;
00772
00773 util::DoReduceInt(src,MPI_SUM,cpu);
00774
00775 if (ParallelDescriptor::MyProc() == cpu)
00776 r = (src == ParallelDescriptor::NProcs()) ? true : false;
00777 }
00778
00779 void
00780 ParallelDescriptor::ReduceBoolOr (bool& r, int cpu)
00781 {
00782 int src = r;
00783
00784 util::DoReduceInt(src,MPI_SUM,cpu);
00785
00786 if (ParallelDescriptor::MyProc() == cpu)
00787 r = (src == 0) ? false : true;
00788 }
00789
00790 void
00791 ParallelDescriptor::Gather (Real* sendbuf,
00792 int nsend,
00793 Real* recvbuf,
00794 int root)
00795 {
00796 BL_ASSERT(root >= 0);
00797 BL_ASSERT(nsend > 0);
00798 BL_ASSERT(!(sendbuf == 0));
00799 BL_ASSERT(!(recvbuf == 0));
00800
00801 MPI_Datatype typ = Mpi_typemap<Real>::type();
00802
00803 BL_MPI_REQUIRE( MPI_Gather(sendbuf,
00804 nsend,
00805 typ,
00806 recvbuf,
00807 nsend,
00808 typ,
00809 root,
00810 Communicator()));
00811 }
00812
00813 MPI_Op
00814 ParallelDescriptor::Max::op ()
00815 {
00816 return MPI_MAX;
00817 }
00818
00819 MPI_Op
00820 ParallelDescriptor::Min::op ()
00821 {
00822 return MPI_MIN;
00823 }
00824
00825 MPI_Op
00826 ParallelDescriptor::Sum::op ()
00827 {
00828 return MPI_SUM;
00829 }
00830
00831 MPI_Op
00832 ParallelDescriptor::Prod::op ()
00833 {
00834 return MPI_PROD;
00835 }
00836
00837 MPI_Op
00838 ParallelDescriptor::Logical_And::op ()
00839 {
00840 return MPI_LAND;
00841 }
00842
00843 MPI_Op
00844 ParallelDescriptor::Boolean_And::op ()
00845 {
00846 return MPI_BAND;
00847 }
00848
00849 MPI_Op
00850 ParallelDescriptor::Logical_Or::op ()
00851 {
00852 return MPI_LOR;
00853 }
00854
00855 MPI_Op
00856 ParallelDescriptor::Boolean_Or::op ()
00857 {
00858 return MPI_BOR;
00859 }
00860
00861 MPI_Op
00862 ParallelDescriptor::Logical_XOr::op ()
00863 {
00864 return MPI_LXOR;
00865 }
00866
00867 MPI_Op
00868 ParallelDescriptor::Boolean_XOr::op ()
00869 {
00870 return MPI_BXOR;
00871 }
00872
00873 template <>
00874 MPI_Datatype
00875 ParallelDescriptor::Mpi_typemap<char>::type ()
00876 {
00877 return MPI_CHAR;
00878 }
00879
00880 template <>
00881 MPI_Datatype
00882 ParallelDescriptor::Mpi_typemap<short>::type ()
00883 {
00884 return MPI_SHORT;
00885 }
00886
00887 template <>
00888 MPI_Datatype
00889 ParallelDescriptor::Mpi_typemap<int>::type ()
00890 {
00891 return MPI_INT;
00892 }
00893
00894 template <>
00895 MPI_Datatype
00896 ParallelDescriptor::Mpi_typemap<long>::type ()
00897 {
00898 return MPI_LONG;
00899 }
00900
00901 template <>
00902 MPI_Datatype
00903 ParallelDescriptor::Mpi_typemap<unsigned char>::type ()
00904 {
00905 return MPI_UNSIGNED_CHAR;
00906 }
00907
00908 template <>
00909 MPI_Datatype
00910 ParallelDescriptor::Mpi_typemap<unsigned short>::type ()
00911 {
00912 return MPI_UNSIGNED_SHORT;
00913 }
00914
00915 template <>
00916 MPI_Datatype
00917 ParallelDescriptor::Mpi_typemap<unsigned int>::type ()
00918 {
00919 return MPI_UNSIGNED;
00920 }
00921
00922 template <>
00923 MPI_Datatype
00924 ParallelDescriptor::Mpi_typemap<unsigned long>::type ()
00925 {
00926 return MPI_UNSIGNED_LONG;
00927 }
00928
00929 template <>
00930 MPI_Datatype
00931 ParallelDescriptor::Mpi_typemap<float>::type ()
00932 {
00933 return MPI_FLOAT;
00934 }
00935
00936 template <>
00937 MPI_Datatype
00938 ParallelDescriptor::Mpi_typemap<double>::type ()
00939 {
00940 return MPI_DOUBLE;
00941 }
00942
00943 void
00944 ParallelDescriptor::Waitsome (Array<MPI_Request>& reqs,
00945 int& completed,
00946 Array<int>& indx,
00947 Array<MPI_Status>& status)
00948 {
00949 BL_MPI_REQUIRE( MPI_Waitsome(reqs.size(),
00950 reqs.dataPtr(),
00951 &completed,
00952 indx.dataPtr(),
00953 status.dataPtr()));
00954 }
00955
00956 #else
00958 void ParallelDescriptor::StartParallel (int*, char***)
00959 {
00960 m_nProcs = 1;
00961 m_nProcsCFD = 1;
00962 m_MyId = 0;
00963 m_comm = 0;
00964 }
00965
00966 void
00967 ParallelDescriptor::Gather (Real* sendbuf,
00968 int nsend,
00969 Real* recvbuf,
00970 int root)
00971 {
00972 BL_ASSERT(root == 0);
00973 BL_ASSERT(nsend > 0);
00974 BL_ASSERT(!(sendbuf == 0));
00975 BL_ASSERT(!(recvbuf == 0));
00976
00977 for (int i = 0; i < nsend; ++i)
00978 recvbuf[i] = sendbuf[i];
00979 }
00980
00981 std::vector<size_t>
00982 ParallelDescriptor::Gather(const size_t& t, int root)
00983 {
00984 std::vector<size_t> resl(1);
00985 resl[0] = t;
00986 return resl;
00987 }
00988
00989 ParallelDescriptor::Message
00990 ParallelDescriptor::Send(const std::vector<size_t>& buf, int dst_pid, int tag)
00991 {
00992 return Message();
00993 }
00994
00995 ParallelDescriptor::Message
00996 ParallelDescriptor::Recv(std::vector<size_t>& buf, int src_pid, int tag)
00997 {
00998 return Message();
00999 }
01000
01001 ParallelDescriptor::Message
01002 ParallelDescriptor::Recv(std::vector<char>& buf, int src_pid, int tag)
01003 {
01004 return Message();
01005 }
01006
01007
01008 ParallelDescriptor::Message::Message ()
01009 :
01010 m_finished(true)
01011 {}
01012
01013 void
01014 ParallelDescriptor::Message::wait ()
01015 {}
01016
01017 bool
01018 ParallelDescriptor::Message::test ()
01019 {
01020 return m_finished;
01021 }
01022
01023 MPI_Request
01024 ParallelDescriptor::Message::req () const
01025 {
01026 return m_req;
01027 }
01028
01029 void ParallelDescriptor::util::SetNProcsCFD () {}
01030
01031 void ParallelDescriptor::EndParallel () {}
01032
01033 void ParallelDescriptor::Abort ()
01034 {
01035 #ifdef WIN32
01036 throw;
01037 #else
01038 std::abort();
01039 #endif
01040 }
01041 void ParallelDescriptor::Abort (int)
01042 {
01043 #ifdef WIN32
01044 throw;
01045 #else
01046 std::abort();
01047 #endif
01048 }
01049
01050 const char* ParallelDescriptor::ErrorString (int) { return ""; }
01051
01052 void ParallelDescriptor::Barrier () {}
01053 void ParallelDescriptor::Barrier (MPI_Comm) {}
01054
01055 void ParallelDescriptor::Test (MPI_Request&, int&, MPI_Status&) {}
01056
01057 void ParallelDescriptor::Comm_dup (MPI_Comm, MPI_Comm&) {}
01058
01059 void ParallelDescriptor::ReduceRealMax (Real&) {}
01060 void ParallelDescriptor::ReduceRealMin (Real&) {}
01061 void ParallelDescriptor::ReduceRealSum (Real&) {}
01062
01063 void ParallelDescriptor::ReduceRealMax (Real&,int) {}
01064 void ParallelDescriptor::ReduceRealMin (Real&,int) {}
01065 void ParallelDescriptor::ReduceRealSum (Real&,int) {}
01066
01067 void ParallelDescriptor::ReduceLongAnd (long&) {}
01068 void ParallelDescriptor::ReduceLongSum (long&) {}
01069 void ParallelDescriptor::ReduceLongMax (long&) {}
01070 void ParallelDescriptor::ReduceLongMin (long&) {}
01071
01072 void ParallelDescriptor::ReduceLongAnd (long&,int) {}
01073 void ParallelDescriptor::ReduceLongSum (long&,int) {}
01074 void ParallelDescriptor::ReduceLongMax (long&,int) {}
01075 void ParallelDescriptor::ReduceLongMin (long&,int) {}
01076
01077 void ParallelDescriptor::ReduceIntSum (int&) {}
01078 void ParallelDescriptor::ReduceIntMax (int&) {}
01079 void ParallelDescriptor::ReduceIntMin (int&) {}
01080
01081 void ParallelDescriptor::ReduceIntSum (int&,int) {}
01082 void ParallelDescriptor::ReduceIntMax (int&,int) {}
01083 void ParallelDescriptor::ReduceIntMin (int&,int) {}
01084
01085 void ParallelDescriptor::ReduceBoolAnd (bool&) {}
01086 void ParallelDescriptor::ReduceBoolOr (bool&) {}
01087
01088 void ParallelDescriptor::ReduceBoolAnd (bool&,int) {}
01089 void ParallelDescriptor::ReduceBoolOr (bool&,int) {}
01090
01091
01092
01093 double
01094 ParallelDescriptor::second ()
01095 {
01096 return BoxLib::wsecond();
01097 }
01098
01099 void
01100 ParallelDescriptor::Waitsome (Array<MPI_Request>& reqs,
01101 int& completed,
01102 Array<int>& indx,
01103 Array<MPI_Status>& status)
01104 {}
01105
01106 #endif
01107
01108
01109
01110 int
01111 ParallelDescriptor::SeqNum ()
01112 {
01113 const int BEG = 1000;
01114 const int END = 9000;
01115
01116 static int seqno = BEG;
01117
01118 int result = seqno++;
01119
01120 if (seqno > END) seqno = BEG;
01121
01122 return result;
01123 }
01124
01125 #if defined(BL_FORT_USE_UPPERCASE)
01126 #define FORT_BL_PD_BARRIER BL_PD_BARRIER
01127 #define FORT_BL_PD_COMMUNICATOR BL_PD_COMMUNICATOR
01128 #define FORT_BL_PD_MYPROC BL_PD_MYPROC
01129 #define FORT_BL_PD_NPROCS BL_PD_NPROCS
01130 #define FORT_BL_PD_IOPROC BL_PD_IOPROC
01131 #define FORT_BL_PD_ABORT BL_PD_ABORT
01132 #elif defined(BL_FORT_USE_LOWERCASE)
01133 #define FORT_BL_PD_BARRIER bl_pd_barrier
01134 #define FORT_BL_PD_COMMUNICATOR bl_pd_communicator
01135 #define FORT_BL_PD_MYPROC bl_pd_myproc
01136 #define FORT_BL_PD_NPROCS bl_pd_nprocs
01137 #define FORT_BL_PD_IOPROC bl_pd_ioproc
01138 #define FORT_BL_PD_ABORT bl_pd_abort
01139 #elif defined(BL_FORT_USE_UNDERSCORE)
01140 #define FORT_BL_PD_BARRIER bl_pd_barrier_
01141 #define FORT_BL_PD_COMMUNICATOR bl_pd_communicator_
01142 #define FORT_BL_PD_MYPROC bl_pd_myproc_
01143 #define FORT_BL_PD_NPROCS bl_pd_nprocs_
01144 #define FORT_BL_PD_IOPROC bl_pd_ioproc_
01145 #define FORT_BL_PD_ABORT bl_pd_abort_
01146 #endif
01147
01148 extern "C"
01149 void
01150 FORT_BL_PD_BARRIER ()
01151 {
01152 ParallelDescriptor::Barrier();
01153 }
01154
01155 extern "C"
01156 void
01157 FORT_BL_PD_COMMUNICATOR (void* vcomm)
01158 {
01159 MPI_Comm* comm = reinterpret_cast<MPI_Comm*>(vcomm);
01160
01161 *comm = ParallelDescriptor::Communicator();
01162 }
01163
01164 extern "C"
01165 void
01166 FORT_BL_PD_MYPROC (int* myproc)
01167 {
01168 *myproc = ParallelDescriptor::MyProc();
01169 }
01170
01171 extern "C"
01172 void
01173 FORT_BL_PD_NPROCS (int* nprocs)
01174 {
01175 *nprocs = ParallelDescriptor::NProcs();
01176 }
01177
01178 extern "C"
01179 void
01180 FORT_BL_PD_IOPROC (int* ioproc)
01181 {
01182 *ioproc = ParallelDescriptor::IOProcessorNumber();
01183 }
01184
01185 extern "C"
01186 void
01187 FORT_BL_PD_ABORT ()
01188 {
01189 ParallelDescriptor::Abort();
01190 }
01191
01192 #ifdef BL_USE_MPI
01193 namespace ParallelDescriptor
01194 {
01195 template <> MPI_Datatype Mpi_typemap<IntVect>::type()
01196 {
01197 static MPI_Datatype mine(MPI_DATATYPE_NULL);
01198 if ( mine == MPI_DATATYPE_NULL )
01199 {
01200 IntVect iv[2];
01201 MPI_Datatype types[] = {
01202 MPI_LB,
01203 MPI_INT,
01204 MPI_UB};
01205 int blocklens[] = { 1, BL_SPACEDIM, 1};
01206 MPI_Aint disp[3];
01207 int n = 0;
01208 BL_MPI_REQUIRE( MPI_Address(&iv[0], &disp[n++]) );
01209 BL_MPI_REQUIRE( MPI_Address(&iv[0].vect, &disp[n++]) );
01210 BL_MPI_REQUIRE( MPI_Address(&iv[1], &disp[n++]) );
01211 for ( int i = n-1; i >= 0; i-- )
01212 {
01213 disp[i] -= disp[0];
01214 }
01215 BL_MPI_REQUIRE( MPI_Type_struct(n, blocklens, disp, types, &mine) );
01216 BL_MPI_REQUIRE( MPI_Type_commit( &mine ) );
01217 }
01218 return mine;
01219 }
01220
01221 template <> MPI_Datatype Mpi_typemap<IndexType>::type()
01222 {
01223 static MPI_Datatype mine(MPI_DATATYPE_NULL);
01224 if ( mine == MPI_DATATYPE_NULL )
01225 {
01226 IndexType iv[2];
01227 MPI_Datatype types[] = {
01228 MPI_LB,
01229 MPI_UNSIGNED,
01230 MPI_UB};
01231 int blocklens[] = { 1, 1, 1};
01232 MPI_Aint disp[3];
01233 int n = 0;
01234 BL_MPI_REQUIRE( MPI_Address(&iv[0], &disp[n++]) );
01235 BL_MPI_REQUIRE( MPI_Address(&iv[0].itype, &disp[n++]) );
01236 BL_MPI_REQUIRE( MPI_Address(&iv[1], &disp[n++]) );
01237 for ( int i = n-1; i >= 0; i-- )
01238 {
01239 disp[i] -= disp[0];
01240 }
01241 BL_MPI_REQUIRE( MPI_Type_struct(n, blocklens, disp, types, &mine) );
01242 BL_MPI_REQUIRE( MPI_Type_commit( &mine ) );
01243 }
01244 return mine;
01245 }
01246
01247 template <> MPI_Datatype Mpi_typemap<Box>::type()
01248 {
01249 static MPI_Datatype mine(MPI_DATATYPE_NULL);
01250 if ( mine == MPI_DATATYPE_NULL )
01251 {
01252 Box iv[2];
01253 MPI_Datatype types[] = {
01254 MPI_LB,
01255 Mpi_typemap<IntVect>::type(),
01256 Mpi_typemap<IntVect>::type(),
01257 Mpi_typemap<IndexType>::type(),
01258 MPI_UB};
01259 int blocklens[] = { 1, 1, 1, 1, 1};
01260 MPI_Aint disp[5];
01261 int n = 0;
01262 BL_MPI_REQUIRE( MPI_Address(&iv[0], &disp[n++]) );
01263 BL_MPI_REQUIRE( MPI_Address(&iv[0].smallend, &disp[n++]) );
01264 BL_MPI_REQUIRE( MPI_Address(&iv[0].bigend, &disp[n++]) );
01265 BL_MPI_REQUIRE( MPI_Address(&iv[0].btype, &disp[n++]) );
01266 BL_MPI_REQUIRE( MPI_Address(&iv[1], &disp[n++]) );
01267 for ( int i = n-1; i >= 0; i-- )
01268 {
01269 disp[i] -= disp[0];
01270 }
01271 BL_MPI_REQUIRE( MPI_Type_struct(n, blocklens, disp, types, &mine) );
01272 BL_MPI_REQUIRE( MPI_Type_commit( &mine ) );
01273 }
01274 return mine;
01275 }
01276 }
01277 #endif
01278