Context.C 21.1 KB
Newer Older
Francois Gygi committed
1 2
////////////////////////////////////////////////////////////////////////////////
//
Francois Gygi committed
3 4 5 6
// Copyright (c) 2008 The Regents of the University of California
//
// This file is part of Qbox
//
Francois Gygi committed
7 8
// Qbox is distributed under the terms of the GNU General Public License
// as published by the Free Software Foundation, either version 2 of
Francois Gygi committed
9 10 11 12 13 14
// the License, or (at your option) any later version.
// See the file COPYING in the root directory of this distribution
// or <http://www.gnu.org/licenses/>.
//
////////////////////////////////////////////////////////////////////////////////
//
Francois Gygi committed
15 16 17 18
// Context.C
//
////////////////////////////////////////////////////////////////////////////////

19
#include "Context.h"
Francois Gygi committed
20 21 22 23
#include <iostream>
#include <iomanip>
#include <cmath>
#include <cassert>
24
#include <string>
25
#include <vector>
Francois Gygi committed
26 27
#include "blacs.h"
#include <mpi.h>
28
using namespace std;
Francois Gygi committed
29 30

////////////////////////////////////////////////////////////////////////////////
31
ContextRep::ContextRep(MPI_Comm comm) : ictxt_(-1), myrow_(-1), mycol_(-1)
Francois Gygi committed
32
{
33
  // construct a single-row Context
Francois Gygi committed
34 35
  int nprocs;
  char order='R';
36 37 38 39
  MPI_Comm_dup(comm,&comm_);
  MPI_Comm_size(comm_,&nprocs);
  MPI_Comm_rank(comm_,&mype_);
  ictxt_ = Csys2blacs_handle(comm_);
Francois Gygi committed
40 41 42 43 44 45 46 47
  nprow_ = 1;
  npcol_ = nprocs;

  Cblacs_gridinit( &ictxt_, &order, nprow_, npcol_ );

  // get values of nprow_, npcol_, myrow_ and mycol_ in the new context
  if ( ictxt_ >= 0 )
    Cblacs_gridinfo(ictxt_, &nprow_, &npcol_, &myrow_, &mycol_);
48

Francois Gygi committed
49 50 51 52 53 54 55 56 57 58
  size_ = nprow_ * npcol_;
  myproc_ = myrow_ < 0 ? -1 : mycol_ + npcol_ * myrow_;
  onpe0_ = ( mype_ == 0 );
  active_ = ( ictxt_ >= 0 );

  pmap_.resize(size_);
  for ( int i = 0; i < size_; i++ )
    pmap_[i] = i;

  MPI_Group group_world, subgroup;
59
  MPI_Comm_group(comm,&group_world);
Francois Gygi committed
60
  MPI_Group_incl(group_world,size_,&pmap_[0],&subgroup);
61
  MPI_Comm_create(comm,subgroup,&comm_);
62 63
  MPI_Group_free(&group_world);
  MPI_Group_free(&subgroup);
Francois Gygi committed
64 65 66
}

////////////////////////////////////////////////////////////////////////////////
67
ContextRep::ContextRep(MPI_Comm comm, int nprow, int npcol) :
Francois Gygi committed
68 69
  ictxt_(-1), myrow_(-1), mycol_(-1), nprow_(nprow), npcol_(npcol)
{
70 71
  int nprocs;
  char order = 'C';
72 73 74 75
  MPI_Comm_dup(comm,&comm_);
  MPI_Comm_size(comm_,&nprocs);
  MPI_Comm_rank(comm_,&mype_);
  ictxt_ = Csys2blacs_handle(comm_);
Francois Gygi committed
76 77
  if ( nprocs < nprow * npcol )
  {
78 79 80
    cout << " nprocs=" << nprocs << endl;
    cout << " Context nprow*npcol > nprocs" << endl;
    Cblacs_abort(ictxt_, 1);
Francois Gygi committed
81 82 83 84 85 86
  }
  Cblacs_gridinit( &ictxt_, &order, nprow, npcol );

  // get values of nprow_, npcol_, myrow_ and mycol_ in the new context
  if ( ictxt_ >= 0 )
    Cblacs_gridinfo(ictxt_, &nprow_, &npcol_, &myrow_, &mycol_);
87

Francois Gygi committed
88
  size_ = nprow_ * npcol_;
89
  myproc_ = Cblacs_pnum(ictxt_,myrow_,mycol_);
Francois Gygi committed
90 91
  onpe0_ = ( mype_ == 0 );
  active_ = ( ictxt_ >= 0 );
92

Francois Gygi committed
93
  pmap_.resize(size_);
94 95 96 97 98 99 100 101
  // column-major order
  int i = 0;
  for ( int ic = 0; ic < npcol; ic++ )
    for ( int ir = 0; ir < nprow; ir++ )
    {
      pmap_[ir+nprow*ic] = i;
      i++;
    }
Francois Gygi committed
102 103

  MPI_Group group_world, subgroup;
104
  MPI_Comm_group(comm,&group_world);
Francois Gygi committed
105
  MPI_Group_incl(group_world,size_,&pmap_[0],&subgroup);
106
  MPI_Comm_create(comm,subgroup,&comm_);
107 108
  MPI_Group_free(&group_world);
  MPI_Group_free(&subgroup);
Francois Gygi committed
109
}
110

Francois Gygi committed
111
////////////////////////////////////////////////////////////////////////////////
112
ContextRep::~ContextRep()
Francois Gygi committed
113
{
114
  if ( myrow_ != -1 )
Francois Gygi committed
115
  {
116 117 118
    // cout << " ContextRep destructor called on ictxt = " << ictxt_ << endl;
    Cblacs_gridexit( ictxt_ );
    MPI_Comm_free(&comm_);
Francois Gygi committed
119
  }
120 121 122 123
}
////////////////////////////////////////////////////////////////////////////////
void ContextRep::dsend(int m, int n, double* a, int lda, int rdest, int cdest)
  const { Cdgesd2d(ictxt_,m,n,a,lda,rdest,cdest); }
124

125 126 127
////////////////////////////////////////////////////////////////////////////////
void ContextRep::drecv(int m, int n, double* a, int lda, int rsrc, int csrc)
  const { Cdgerv2d(ictxt_,m,n,a,lda,rsrc,csrc); }
128

129 130 131 132
////////////////////////////////////////////////////////////////////////////////
void ContextRep::dsum(char scope, char topology, int m, int n, double* a,
  int lda, int rdest, int cdest) const
{ Cdgsum2d(ictxt_,&scope,&topology,m,n,a,lda,rdest,cdest); }
133

134 135 136 137
////////////////////////////////////////////////////////////////////////////////
void ContextRep::dmax(char scope, char topology, int m, int n, double* a,
  int lda, int rdest, int cdest) const
{ Cdgamx2d(ictxt_,&scope,&topology,m,n,a,lda,(int*)0,(int*)0,-1,rdest,cdest); }
138

139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217
////////////////////////////////////////////////////////////////////////////////
void ContextRep::dmax(char scope, char topology, int m, int n, double* a,
  int lda, int* ra, int* ca, int rcflag, int rdest, int cdest) const
{ Cdgamx2d(ictxt_,&scope,&topology,m,n,a,lda,ra,ca,rcflag,rdest,cdest); }

////////////////////////////////////////////////////////////////////////////////
void ContextRep::dmin(char scope, char topology, int m, int n, double* a,
   int lda, int* ra, int* ca, int rcflag, int rdest, int cdest) const
{ Cdgamn2d(ictxt_,&scope,&topology,m,n,a,lda,ra,ca,rcflag,rdest,cdest); }

////////////////////////////////////////////////////////////////////////////////
void ContextRep::dmin(char scope, char topology, int m, int n, double* a,
  int lda, int rdest, int cdest) const
{ Cdgamn2d(ictxt_,&scope,&topology,m,n,a,lda,(int*)0,(int*)0,-1,rdest,cdest); }

////////////////////////////////////////////////////////////////////////////////
void ContextRep::dbcast_send(char scope, char topology, int m, int n,
  double* a,int lda) const
{ Cdgebs2d(ictxt_,&scope,&topology,m,n,a,lda); }

////////////////////////////////////////////////////////////////////////////////
void ContextRep::dbcast_recv(char scope, char topology, int m, int n,
  double* a, int lda, int rsrc, int csrc) const
{ Cdgebr2d(ictxt_,&scope,&topology,m,n,a,lda,rsrc,csrc); }

////////////////////////////////////////////////////////////////////////////////
void ContextRep::isend(int m, int n, int* a, int lda, int rdest, int cdest)
  const { Cigesd2d(ictxt_,m,n,a,lda,rdest,cdest); }

////////////////////////////////////////////////////////////////////////////////
void ContextRep::irecv(int m, int n, int* a, int lda, int rsrc, int csrc) const
{ Cigerv2d(ictxt_,m,n,a,lda,rsrc,csrc); }

////////////////////////////////////////////////////////////////////////////////
void ContextRep::isum(char scope, char topology, int m, int n, int* a, int lda,
  int rdest, int cdest) const
{ Cigsum2d(ictxt_,&scope,&topology,m,n,a,lda,rdest,cdest); }

////////////////////////////////////////////////////////////////////////////////
void ContextRep::imax(char scope, char topology, int m, int n, int* a, int lda,
  int* ra, int* ca, int rcflag, int rdest, int cdest) const
{ Cigamx2d(ictxt_,&scope,&topology,m,n,a,lda,ra,ca,rcflag,rdest,cdest); }

////////////////////////////////////////////////////////////////////////////////
void ContextRep::imax(char scope, char topology, int m, int n, int* a, int lda,
  int rdest, int cdest) const
{ Cigamx2d(ictxt_,&scope,&topology,m,n,a,lda,(int*)0,(int*)0,-1,rdest,cdest); }

////////////////////////////////////////////////////////////////////////////////
void ContextRep::imin(char scope, char topology, int m, int n, int* a, int lda,
  int* ra, int* ca, int rcflag, int rdest, int cdest) const
{ Cigamn2d(ictxt_,&scope,&topology,m,n,a,lda,ra,ca,rcflag,rdest,cdest); }

////////////////////////////////////////////////////////////////////////////////
void ContextRep::imin(char scope, char topology, int m, int n, int* a, int lda,
  int rdest, int cdest) const
{ Cigamn2d(ictxt_,&scope,&topology,m,n,a,lda,(int*)0,(int*)0,-1,rdest,cdest); }

////////////////////////////////////////////////////////////////////////////////
void ContextRep::ibcast_send(char scope, char topology,
  int m, int n, int* a,int lda) const
{ Cigebs2d(ictxt_,&scope,&topology,m,n,a,lda); }

////////////////////////////////////////////////////////////////////////////////
void ContextRep::ibcast_recv(char scope, char topology, int m, int n, int* a,
  int lda, int rsrc, int csrc) const
{ Cigebr2d(ictxt_,&scope,&topology,m,n,a,lda,rsrc,csrc); }

////////////////////////////////////////////////////////////////////////////////
void ContextRep::string_send(std::string& s, int rdest, int cdest) const
{
  int len = s.size();
  isend(1,1,&len,1,rdest,cdest);
  int ilen = len/(sizeof(int)/sizeof(char));
  if ( len%(sizeof(int)/sizeof(char)) != 0 ) ilen++;
  int* ibuf = new int[ilen];
  s.copy((char*)ibuf,std::string::npos);
  isend(ilen,1,ibuf,ilen,rdest,cdest);
  delete [] ibuf;
Francois Gygi committed
218 219 220
}

////////////////////////////////////////////////////////////////////////////////
221
void ContextRep::string_recv(std::string& s, int rsrc, int csrc) const
Francois Gygi committed
222
{
223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238
  int len = -1;
  irecv(1,1,&len,1,rsrc,csrc);
  int ilen = len/(sizeof(int)/sizeof(char));
  if ( len%(sizeof(int)/sizeof(char)) != 0 ) ilen++;
  int* ibuf = new int[ilen];
  irecv(ilen,1,ibuf,ilen,rsrc,csrc);
  s.resize(len);
  s.assign((char*)ibuf,len);
  delete [] ibuf;
}

////////////////////////////////////////////////////////////////////////////////
void ContextRep::string_bcast(std::string& s, int isrc) const
{
  int len;
  if ( mype() == isrc )
Francois Gygi committed
239
  {
240 241 242 243 244 245 246 247 248 249
    len = s.length();
  }
  MPI_Bcast(&len,1,MPI_INT,isrc,comm());
  char* buf = new char[len+1];
  // s is initialized only on task isrc
  if ( mype() == isrc )
  {
    s.copy(buf,std::string::npos);
    buf[len]=0;
    assert(buf[len]=='\0');
Francois Gygi committed
250
  }
251 252 253
  MPI_Bcast(buf,len+1,MPI_CHAR,isrc,comm());
  s = buf;
  delete [] buf;
Francois Gygi committed
254 255 256 257 258 259 260 261 262 263 264 265 266 267
}

////////////////////////////////////////////////////////////////////////////////
void ContextRep::print(ostream& os) const
{
  if ( active_ )
  {
    os << " " << nprow_ << "x" << npcol_ << " ";
    os << "{ ";
    for ( int ir = 0; ir < nprow_; ir++ )
    {
      os << "{ ";
      for ( int ic = 0; ic < npcol_; ic++ )
      {
268
        os << pmap_[ir+nprow_*ic] << " ";
Francois Gygi committed
269 270 271
      }
      os << "} ";
    }
272 273 274 275 276 277
    //os << "} (ictxt=" << ictxt_ << ")" << endl;
    os << "} (ictxt=" << ictxt_ << ")";
    int handle;
    Cblacs_get(ictxt_,10,&handle);
    os << "(handle=" << handle << ")";
    os << "(comm=" << comm_<< ")" << endl;
Francois Gygi committed
278 279 280 281 282 283 284 285
  }
  else
  {
    os << " (inactive)" << endl;
  }
}

////////////////////////////////////////////////////////////////////////////////
286
int Context::ictxt() const { return rep->ictxt(); }
Francois Gygi committed
287 288

////////////////////////////////////////////////////////////////////////////////
289
int Context::myrow() const { return rep->myrow(); }
Francois Gygi committed
290 291

////////////////////////////////////////////////////////////////////////////////
292
int Context::mycol() const { return rep->mycol(); }
Francois Gygi committed
293 294

////////////////////////////////////////////////////////////////////////////////
295
int Context::nprow() const { return rep->nprow(); }
Francois Gygi committed
296 297

////////////////////////////////////////////////////////////////////////////////
298
int Context::npcol() const { return rep->npcol(); }
Francois Gygi committed
299 300

////////////////////////////////////////////////////////////////////////////////
301
int Context::size() const { return rep->size(); }
Francois Gygi committed
302 303

////////////////////////////////////////////////////////////////////////////////
304
int Context::myproc() const { return rep->myproc(); }
Francois Gygi committed
305 306

////////////////////////////////////////////////////////////////////////////////
307
int Context::mype() const { return rep->mype(); }
Francois Gygi committed
308 309 310

////////////////////////////////////////////////////////////////////////////////
int Context::pmap(int irow, int icol) const
311
{ return rep->pmap(irow,icol); }
Francois Gygi committed
312 313

////////////////////////////////////////////////////////////////////////////////
314
bool Context::onpe0(void) const { return rep->onpe0(); }
Francois Gygi committed
315 316

////////////////////////////////////////////////////////////////////////////////
317
bool Context::active(void) const { return rep->active(); }
Francois Gygi committed
318 319

////////////////////////////////////////////////////////////////////////////////
320
void Context::abort(int ierr) const { rep->abort(ierr); }
Francois Gygi committed
321 322

////////////////////////////////////////////////////////////////////////////////
323
void Context::barrier(void) const { rep->barrier(); }
Francois Gygi committed
324 325

////////////////////////////////////////////////////////////////////////////////
326
void Context::barrier(char scope) const { rep->barrier(scope); }
Francois Gygi committed
327 328

////////////////////////////////////////////////////////////////////////////////
329
void Context::dsend(int m, int n, double* a,
Francois Gygi committed
330
  int lda, int rdest, int cdest) const
331
{ rep->dsend(m,n,a,lda,rdest,cdest); }
Francois Gygi committed
332

333
////////////////////////////////////////////////////////////////////////////////
334
void Context::drecv(int m, int n, double* a,
Francois Gygi committed
335
  int lda, int rsrc, int csrc) const
336
{ rep->drecv(m,n,a,lda,rsrc,csrc); }
Francois Gygi committed
337

338
////////////////////////////////////////////////////////////////////////////////
339
void Context::dsum(char scope, char topology,
Francois Gygi committed
340
  int m, int n, double* a, int lda, int rdest, int cdest) const
341
{ rep->dsum(scope,topology,m,n,a,lda,rdest,cdest); }
Francois Gygi committed
342

343
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
344
void Context::dsum(char scope, int m, int n, double* a, int lda) const
345
{ rep->dsum(scope,' ',m,n,a,lda,-1,-1); }
Francois Gygi committed
346

347
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
348
void Context::dsum(int m, int n, double* a, int lda) const
349
{ rep->dsum('A',' ',m,n,a,lda,-1,-1); }
Francois Gygi committed
350

351
////////////////////////////////////////////////////////////////////////////////
352 353
void Context::dmax(char scope, char topology,
  int m, int n, double* a, int lda, int* ra, int* ca, int rcflag,
Francois Gygi committed
354
  int rdest, int cdest) const
355
{ rep->dmax(scope,topology,m,n,a,lda,ra,ca,rcflag,rdest,cdest); }
Francois Gygi committed
356

357
////////////////////////////////////////////////////////////////////////////////
358
void Context::dmax(char scope, char topology,
Francois Gygi committed
359
  int m, int n, double* a, int lda, int rdest, int cdest) const
360
{ rep->dmax(scope,topology,m,n,a,lda,rdest,cdest); }
Francois Gygi committed
361

362
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
363
void Context::dmax(char scope, int m, int n, double* a, int lda) const
364
{ rep->dmax(scope,' ',m,n,a,lda,-1,-1); }
Francois Gygi committed
365

366
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
367
void Context::dmax(int m, int n, double* a, int lda) const
368
{ rep->dmax('A',' ',m,n,a,lda,-1,-1); }
Francois Gygi committed
369

370
////////////////////////////////////////////////////////////////////////////////
371 372
void Context::dmin(char scope, char topology,
  int m, int n, double* a, int lda, int* ra, int* ca, int rcflag,
Francois Gygi committed
373
  int rdest, int cdest) const
374
{ rep->dmin(scope,topology,m,n,a,lda,ra,ca,rcflag,rdest,cdest); }
Francois Gygi committed
375

376
////////////////////////////////////////////////////////////////////////////////
377
void Context::dmin(char scope, char topology,
Francois Gygi committed
378
  int m, int n, double* a, int lda, int rdest, int cdest) const
379
{ rep->dmin(scope,topology,m,n,a,lda,rdest,cdest); }
Francois Gygi committed
380

381
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
382
void Context::dmin(char scope, int m, int n, double* a, int lda) const
383
{ rep->dmin(scope,' ',m,n,a,lda,-1,-1); }
Francois Gygi committed
384

385
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
386
void Context::dmin(int m, int n, double* a, int lda) const
387
{ rep->dmin('A',' ',m,n,a,lda,-1,-1); }
Francois Gygi committed
388

389
////////////////////////////////////////////////////////////////////////////////
390
void Context::dbcast_send(char scope, char topology,
Francois Gygi committed
391
  int m, int n, double* a, int lda) const
392
{ rep->dbcast_send(scope,topology,m,n,a,lda); }
Francois Gygi committed
393

394
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
395
void Context::dbcast_send(char scope, int m, int n, double* a, int lda) const
396
{ rep->dbcast_send(scope,' ',m,n,a,lda); }
Francois Gygi committed
397

398
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
399
void Context::dbcast_send(int m, int n, double* a, int lda) const
400
{ rep->dbcast_send('A',' ',m,n,a,lda); }
Francois Gygi committed
401

402
////////////////////////////////////////////////////////////////////////////////
403
void Context::dbcast_recv(char scope, char topology,
Francois Gygi committed
404
  int m, int n, double* a, int lda, int rsrc,int csrc) const
405
{ rep->dbcast_recv(scope,topology,m,n,a,lda,rsrc,csrc); }
Francois Gygi committed
406

407
////////////////////////////////////////////////////////////////////////////////
408
void Context::dbcast_recv(char scope,
Francois Gygi committed
409
  int m, int n, double* a, int lda, int rsrc, int csrc) const
410
{ rep->dbcast_recv(scope,' ',m,n,a,lda,rsrc,csrc); }
Francois Gygi committed
411

412
////////////////////////////////////////////////////////////////////////////////
413
void Context::dbcast_recv(int m, int n, double* a, int lda,
Francois Gygi committed
414
  int rsrc,int csrc) const
415
{ rep->dbcast_recv('A',' ',m,n,a,lda,rsrc,csrc); }
Francois Gygi committed
416 417

////////////////////////////////////////////////////////////////////////////////
418
void Context::isend(int m, int n, int* a,
Francois Gygi committed
419
  int lda, int rdest, int cdest) const
420
{ rep->isend(m,n,a,lda,rdest,cdest); }
Francois Gygi committed
421

422
////////////////////////////////////////////////////////////////////////////////
423
void Context::irecv(int m, int n, int* a,
Francois Gygi committed
424
  int lda, int rsrc, int csrc) const
425
{ rep->irecv(m,n,a,lda,rsrc,csrc); }
Francois Gygi committed
426

427
////////////////////////////////////////////////////////////////////////////////
428
void Context::isum(char scope, char topology,
Francois Gygi committed
429
  int m, int n, int* a, int lda, int rdest, int cdest) const
430
{ rep->isum(scope,topology,m,n,a,lda,rdest,cdest); }
Francois Gygi committed
431

432
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
433
void Context::isum(char scope, int m, int n, int* a, int lda) const
434
{ rep->isum(scope,' ',m,n,a,lda,-1,-1); }
Francois Gygi committed
435

436
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
437
void Context::isum(int m, int n, int* a, int lda) const
438
{ rep->isum('A',' ',m,n,a,lda,-1,-1); }
Francois Gygi committed
439

440
////////////////////////////////////////////////////////////////////////////////
441 442
void Context::imax(char scope, char topology,
  int m, int n, int* a, int lda, int* ra, int* ca, int rcflag,
Francois Gygi committed
443
  int rdest, int cdest) const
444
{ rep->imax(scope,topology,m,n,a,lda,ra,ca,rcflag,rdest,cdest); }
Francois Gygi committed
445

446
////////////////////////////////////////////////////////////////////////////////
447
void Context::imax(char scope, char topology,
Francois Gygi committed
448
  int m, int n, int* a, int lda, int rdest, int cdest) const
449
{ rep->imax(scope,topology,m,n,a,lda,rdest,cdest); }
Francois Gygi committed
450

451
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
452
void Context::imax(char scope, int m, int n, int* a, int lda) const
453
{ rep->imax(scope,' ',m,n,a,lda,-1,-1); }
Francois Gygi committed
454

455
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
456
void Context::imax(int m, int n, int* a, int lda) const
457
{ rep->imax('A',' ',m,n,a,lda,-1,-1); }
Francois Gygi committed
458

459
////////////////////////////////////////////////////////////////////////////////
460 461
void Context::imin(char scope, char topology,
  int m, int n, int* a, int lda, int* ra, int* ca, int rcflag,
Francois Gygi committed
462
  int rdest, int cdest) const
463
{ rep->imin(scope,topology,m,n,a,lda,ra,ca,rcflag,rdest,cdest); }
Francois Gygi committed
464

465
////////////////////////////////////////////////////////////////////////////////
466
void Context::imin(char scope, char topology,
Francois Gygi committed
467
  int m, int n, int* a, int lda, int rdest, int cdest) const
468
{ rep->imin(scope,topology,m,n,a,lda,rdest,cdest); }
Francois Gygi committed
469

470
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
471
void Context::imin(char scope, int m, int n, int* a, int lda) const
472
{ rep->imin(scope,' ',m,n,a,lda,-1,-1); }
Francois Gygi committed
473

474
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
475
void Context::imin(int m, int n, int* a, int lda) const
476
{ rep->imin('A',' ',m,n,a,lda,-1,-1); }
Francois Gygi committed
477

478
////////////////////////////////////////////////////////////////////////////////
479
void Context::ibcast_send(char scope, char topology,
Francois Gygi committed
480
  int m, int n, int* a, int lda) const
481
{ rep->ibcast_send(scope,topology,m,n,a,lda); }
Francois Gygi committed
482

483
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
484
void Context::ibcast_send(char scope, int m, int n, int* a, int lda) const
485
{ rep->ibcast_send(scope,' ',m,n,a,lda); }
Francois Gygi committed
486

487
////////////////////////////////////////////////////////////////////////////////
Francois Gygi committed
488
void Context::ibcast_send(int m, int n, int* a, int lda) const
489
{ rep->ibcast_send('A',' ',m,n,a,lda); }
Francois Gygi committed
490

491
////////////////////////////////////////////////////////////////////////////////
492
void Context::ibcast_recv(char scope, char topology,
Francois Gygi committed
493
  int m, int n, int* a, int lda, int rsrc,int csrc) const
494
{ rep->ibcast_recv(scope,topology,m,n,a,lda,rsrc,csrc); }
Francois Gygi committed
495

496
////////////////////////////////////////////////////////////////////////////////
497
void Context::ibcast_recv(char scope,
Francois Gygi committed
498
  int m, int n, int* a, int lda, int rsrc, int csrc) const
499
{ rep->ibcast_recv(scope,' ',m,n,a,lda,rsrc,csrc); }
Francois Gygi committed
500

501
////////////////////////////////////////////////////////////////////////////////
502
void Context::ibcast_recv(int m, int n, int* a, int lda,
Francois Gygi committed
503
  int rsrc,int csrc) const
504
{ rep->ibcast_recv('A',' ',m,n,a,lda,rsrc,csrc); }
Francois Gygi committed
505

506 507 508
////////////////////////////////////////////////////////////////////////////////
void Context::string_send(std::string& s, int rdest, int cdest) const
{ rep->string_send(s,rdest,cdest); }
Francois Gygi committed
509

510 511 512
////////////////////////////////////////////////////////////////////////////////
void Context::string_recv(std::string& s, int rsrc, int csrc) const
{ rep->string_recv(s,rsrc,csrc); }
Francois Gygi committed
513

514 515 516
////////////////////////////////////////////////////////////////////////////////
void Context::string_bcast(std::string& s, int isrc) const
{ rep->string_bcast(s,isrc); }
Francois Gygi committed
517 518 519

////////////////////////////////////////////////////////////////////////////////
bool Context::operator==(const Context& ctxt) const
520
{ return ( rep->ictxt() == ctxt.ictxt() ); }
Francois Gygi committed
521 522

////////////////////////////////////////////////////////////////////////////////
523
MPI_Comm Context::comm(void) const { return rep->comm(); }
524

Francois Gygi committed
525
////////////////////////////////////////////////////////////////////////////////
526
void Context::print(ostream& os) const { rep->print(os);}
Francois Gygi committed
527 528 529 530 531 532 533

////////////////////////////////////////////////////////////////////////////////
ostream& operator<<(ostream& os, const Context& c)
{
  c.print(os);
  return os;
}