1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 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 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463 464 465 466 467 468 469 470 471 472 473 474 475 476 477 478 479 480 481 482 483 484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816 817 818 819 820 821 822 823 824 825 826 827 828 829 830 831 832 833 834 835 836 837 838 839 840 841 842 843
|
// for automatic compilation with ff-c++
// FFCS - 23/5/12 - remove metis dependency because it interfers with identically-named libmetis.a from parmetis
//ff-c++-LIBRARY-dep: superlu_dist ptscotchparmetis ptscotch scotchmetis scotch mpi blas fc
//ff-c++-cpp-dep:
// ORIG-DATE: 02/2009
// -*- Mode : c++ -*-
//
// SUMMARY :
// USAGE : LGPL
// ORG : LJLL Universite Pierre et Marie Curie, Paris, FRANCE
// AUTHOR : Jacques Morice
// E-MAIL : jacques.morice@ann.jussieu.fr
//
/*
This file is part of Freefem++
Freefem++ is free software; you can redistribute it and/or modify
it under the terms of the GNU Lesser General Public License as published by
the Free Software Foundation; either version 2.1 of the License, or
(at your option) any later version.
Freefem++ is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU Lesser General Public License for more details.
You should have received a copy of the GNU Lesser General Public License
along with Freefem++; if not, write to the Free Software
Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
Thank to the ARN () FF2A3 grant
ref:ANR-07-CIS7-002-01
*/
/*
Interface freefem++ et SuperLU_DIST_2.3
/bin/sh ff-mpic++ zSuperLU_DIST.cpp -I/Users/morice/librairie/SuperLU_DIST_2.3/SRC/ -L/Users/morice/librairie/openmpi/lib/ -lmpi -lopen-pal -lopen-rte -L/Users/morice/librairie/PATCHVECLIB/ -lwrapperdotblas -framework veclib -L/Users/morice/librairie/ParMetis-3.1/ -lparmetis -lmetis -L/Users/morice/librairie/SuperLU_DIST_2.3/lib/ -lsuperlu_dist_2.3
*/
// FFCS - required to define __int64 for MSMPI
#include <stdint.h>
#include <mpi.h>
#include <iostream>
using namespace std;
#include "rgraph.hpp"
#include "error.hpp"
#include "AFunction.hpp"
//#include "lex.hpp"
#include "MatriceCreuse_tpl.hpp"
#include "superlu_zdefs.h"
#include "ffsuperludistoption.hpp"
template <class R> struct SuperLUmpiDISTDriver
{
};
template <> struct SuperLUmpiDISTDriver<Complex>
{
/* Driver routines */
/* Driver routines */
static Dtype_t R_SLU_T() { return SLU_Z;}
static doublecomplex *dc(Complex *p) { return (doublecomplex *) (void *) p;}
static doublecomplex **dc(Complex **p) { return (doublecomplex **) (void *) p;}
// Remplacement doublecomplex par Complex
static void
pgssvx(superlu_options_t *p1, SuperMatrix *p2,
ScalePermstruct_t *p3,
Complex *p4, int p5, int p6, gridinfo_t *p7,
LUstruct_t *p8, SOLVEstruct_t *p9, double *p10,
SuperLUStat_t *p11, int *p12)
{
pzgssvx( p1,p2,p3,
dc(p4),p5,p6,p7,p8,p9,p10,p11,p12 );
}
static void
pgssvx_ABglobal(superlu_options_t *p1, SuperMatrix *p2,
ScalePermstruct_t *p3,
Complex *p4, int p5, int p6, gridinfo_t *p7,
LUstruct_t *p8, double *p9,
SuperLUStat_t *p10, int *p11)
{
pzgssvx_ABglobal( p1,p2,p3,
dc(p4),p5,p6,p7,p8,p9,p10,p11 );
}
static void
Create_CompCol_Matrix_dist(SuperMatrix *p1, int_t p2, int_t p3, int_t p4,
Complex *p5, int_t *p6, int_t *p7,
Stype_t p8, Dtype_t p9, Mtype_t p10)
{
zCreate_CompCol_Matrix_dist( p1,p2,p3,
p4,dc(p5),p6,p7,p8,p9,p10 );
}
static void
Create_CompRowLoc_Matrix_dist(SuperMatrix *p1, int_t p2, int_t p3,
int_t p4, int_t p5, int_t p6,
Complex *p7, int_t *p8, int_t *p9,
Stype_t p10, Dtype_t p11, Mtype_t p12)
{
zCreate_CompRowLoc_Matrix_dist( p1,p2,p3,
p4,p5,p6,dc(p7),p8,p9,p10,p11,p12 );
}
static void
CompRow_to_CompCol_dist(int_t p1, int_t p2, int_t p3,
Complex *p4, int_t *p5, int_t *p6,
Complex **p7, int_t **p8, int_t **p9)
{
zCompRow_to_CompCol_dist( p1,p2,p3,
dc(p4),p5,p6,dc(p7),p8,p9 );
}
static void
Create_Dense_Matrix_dist(SuperMatrix *p1, int_t p2, int_t p3, Complex *p4,
int_t p5, Stype_t p6, Dtype_t p7,Mtype_t p8)
{
zCreate_Dense_Matrix_dist( p1,p2,p3,
dc(p4),p5,p6,p7,p8 );
}
static void
Create_SuperNode_Matrix_dist(SuperMatrix *p1, int_t p2, int_t p3, int_t p4,
Complex *p5, int_t *p6,
int_t *p7, int_t *p8,
int_t *p9, int_t *p10,
Stype_t p11, Dtype_t p12, Mtype_t p13)
{
zCreate_SuperNode_Matrix_dist( p1,p2,p3,
p4,dc(p5),p6,p7,p8,p9,p10,p11,p12,p13 );
}
static void
Print_CompRowLoc_Matrix_dist(SuperMatrix *p1)
{
zPrint_CompRowLoc_Matrix_dist(p1);
}
};
template<class R>
class ZSolveSuperLUmpi : public MatriceMorse<R>::VirtualSolver, public SuperLUmpiDISTDriver<R> {
double eps;
mutable double epsr;
double tgv;
double tol_pivot_sym,tol_pivot; //Add 31 oct 2005
//mutable char equed[1];
//yes_no_t equil;
mutable SuperMatrix A;
NCformat *Astore;
//NCformat *Ustore;
//SCformat *Lstore;
mutable superlu_options_t options;
mutable mem_usage_t mem_usage;
mutable ScalePermstruct_t ScalePermstruct;
mutable LUstruct_t LUstruct;
mutable SOLVEstruct_t SOLVEstruct;
mutable gridinfo_t grid;
string string_option;
string data_option;
R *a;
int *asub, *xa;
int_t m, n, nnz;
// rajout pour //
int_t nprow,npcol; /* process rows and process columns*/
int matrixdist; // type of distributed matrix
MPI_Comm commworld ;
static const int assembled =0;
static const int distributedglobal =1;
static const int distributed =2;
public:
ZSolveSuperLUmpi(const MatriceMorse<R> &AA,int strategy,double ttgv, double epsilon,
double pivot,double pivot_sym, string datafile,
string param_char, KN<long> &pperm_r, KN<long> &pperm_c, void * ccommworld ) :
eps(epsilon),epsr(0),
tgv(ttgv),string_option(param_char),data_option(datafile),
tol_pivot_sym(pivot_sym),tol_pivot(pivot)
{
commworld = ccommworld ? *static_cast<MPI_Comm*>( ccommworld) : MPI_COMM_WORLD;
R* B;
//R* X;
SuperLUStat_t stat;
int info, ldb, nrhs=0;
int i;
double* berr;
int iam;
// Add for distributed matrix
int_t m_loc, m_loc_fst, fst_row, nnz_loc, fst_nnz;
R *aloc;
int *asubloc, *xaloc;
// End Add for distributed matrix
// time variables
long int starttime,finishtime;
long int timeused;
if(verbosity) starttime = clock();
A.Store=0;
/* Defaults */
nrhs = 0;
/* lecture de nprow and npcol */
// Cas max deux procs
nprow = 1;
MPI_Comm_size(commworld,&npcol);
matrixdist=0;
/* set the default options */
set_default_options_dist(&options);
DiagScale_t optionDiagScale;
//if(verbosity > 10) print_options_dist(&options);
if(!string_option.empty()) read_nprow_npcol_freefem( &string_option, &nprow, &npcol, &matrixdist);
if(!string_option.empty()) read_options_freefem(&string_option,&options,&optionDiagScale);
if(!data_option.empty()) read_options_superlu_datafile(&data_option,&options,&nprow, &npcol, &matrixdist,&optionDiagScale);
//if(verbosity > 10) print_options_dist(&options);
/* ------------------------------------------------------------
INITIALIZE THE SUPERLU PROCESS GRID.
------------------------------------------------------------*/
cout << "Complex superlu_gridinit " << commworld << " "<< ccommworld << " : " << nprow << "X" << npcol <<endl;
superlu_gridinit(commworld , nprow, npcol, &grid);
cout << " --\n";
/* Bail out if I do not belong in the grid. */
iam = grid.iam;
if ( iam >= nprow * npcol ){
//superlu_gridexit(&grid);
printf("this process is not used in superlu %d \n",iam);
}
else
{
// matrix to procs and vectors
if( matrixdist == assembled ){
if(!iam){
cout << "iam=" << iam << endl;
printf("\tProcess grid\t%d X %d\n", grid.nprow, grid.npcol);
/* create the matrix for superlu_dist */
n=AA.n;
m=AA.m;
nnz=AA.nbcoef;
assert( AA.lg[n] == nnz );
printf("\tDimension\t%dx%d\t # nonzeros %d\n", m, n, nnz);
/* transform Row to Col */
// cela coute cher comme fonction //
// dallocateA_dist(n, nnz, &a, &asub, &xa);
// dCompRow_to_CompCol_dist(m,n,nnz,arow,asubrow,xarow,&a,&asub,&xa);
// FFCS - "this->" required by g++ 4.7
this->CompRow_to_CompCol_dist(m,n,nnz,AA.a,AA.cl,AA.lg,&a,&asub,&xa);
/* Broadcast matrix A to the other PEs. */
MPI_Bcast( &m, 1, mpi_int_t, 0, grid.comm );
MPI_Bcast( &n, 1, mpi_int_t, 0, grid.comm );
MPI_Bcast( &nnz, 1, mpi_int_t, 0, grid.comm );
MPI_Bcast( a, nnz, SuperLU_MPI_DOUBLE_COMPLEX, 0, grid.comm );
MPI_Bcast( asub, nnz, mpi_int_t, 0, grid.comm );
MPI_Bcast( xa, n+1, mpi_int_t, 0, grid.comm );
}
else{
/*
printf("\tProcess grid\t%d X %d\n", grid.nprow, grid.npcol);
Receive matrix A from PE 0. */
MPI_Bcast( &m, 1, mpi_int_t, 0, grid.comm );
MPI_Bcast( &n, 1, mpi_int_t, 0, grid.comm );
MPI_Bcast( &nnz, 1, mpi_int_t, 0, grid.comm );
/* Allocate storage for compressed column representation. */
// FFCS - "this->" required by g++ 4.7
zallocateA_dist(n, nnz, this->dc(&a), &asub, &xa);
MPI_Bcast( a, nnz, SuperLU_MPI_DOUBLE_COMPLEX, 0, grid.comm );
MPI_Bcast( asub, nnz, mpi_int_t, 0, grid.comm );
MPI_Bcast( xa, n+1, mpi_int_t, 0, grid.comm );
}
Dtype_t R_SLU = SuperLUmpiDISTDriver<R>::R_SLU_T();
cout << "Debut: Create_CompCol_Matrix_dist" <<endl;
// FFCS - "this->" required by g++ 4.7
this->Create_CompCol_Matrix_dist(&A, m, n, nnz, a, asub, xa, SLU_NC, R_SLU, SLU_GE);
cout << "Fin: Create_CompCol_Matrix_dist" <<endl;
/* creation of pseudo solution + second member */
if ( !(B = new Complex[m] ) ){
//if ( !(B = doublecomplexMalloc_dist(m )) ){
printf("probleme d allocation\n");
exit(1);
}
if(verbosity)
printf("Dimension %dx%d; # nonzeros %d\n", A.nrow, A.ncol, nnz);
// /* set the default options */
// set_default_options_dist(&options);
// DiagScale_t optionDiagScale;
// if(!string_option.empty()) read_options_freefem(&string_option,&options,&optionDiagScale);
/* Initialize ScalePermstruct and LUstruct. */
ScalePermstructInit(m, n, &ScalePermstruct);
if(!(pperm_r==NULL) || !(pperm_c==NULL) ) ScalePermstruct.DiagScale=optionDiagScale;
if( !(pperm_r==NULL) )
for(int ii=0; ii<m; ii++) ScalePermstruct.perm_r[ii] = pperm_r[ii];
if( !(pperm_c==NULL) )
for(int ii=0; ii<n; ii++) ScalePermstruct.perm_c[ii]= pperm_c[ii];
if( ScalePermstruct.DiagScale != NOEQUIL ){
printf("FreeFem++ doesn't support change of the original matrix");
exit(1);
}
LUstructInit(m, n, &LUstruct);
/* Initialize the statistics variables. */
PStatInit(&stat);
ldb = m;
nrhs=1;
if ( !(berr = doubleMalloc_dist(nrhs )) ){
printf("probleme d allocation\n");
exit(1);
}
berr[0]=0.;
if(verbosity)
printf("Dimension %dx%d; # nonzeros %d\n", A.nrow, A.ncol, nnz);
/* INIT LU struct*/
/* ONLY PERFORM THE LU DECOMPOSITION */
//B.ncol = 0; /* Indicate not to solve the system */
nrhs=0;
SuperLUmpiDISTDriver<R>::pgssvx_ABglobal(&options, &A, &ScalePermstruct, B, ldb, nrhs, &grid,
&LUstruct, berr, &stat, &info);
if(verbosity)
printf("LU factorization: pdgssvx()/p returns info %d\n", info);
if ( verbosity) PStatPrint(&options,&stat,&grid);
PStatFree(&stat);
}
//##########################################################
//
// matrix distributed with matrix global given
//
//##########################################################
else if( matrixdist == distributedglobal) {
if(!iam){
printf("\tProcess grid\t%d X %d\n", grid.nprow, grid.npcol);
/* create the matrix for superlu_dist */
n=AA.n;
m=AA.m;
nnz=AA.nbcoef;
a=AA.a;
asub=AA.cl;
xa=AA.lg;
xa[n] = nnz;
printf("\tDimension\t%dx%d\t # nonzeros %d\n", m, n, nnz);
/* Broadcast matrix A to the other PEs. */
MPI_Bcast( &m, 1, mpi_int_t, 0, grid.comm );
MPI_Bcast( &n, 1, mpi_int_t, 0, grid.comm );
MPI_Bcast( &nnz, 1, mpi_int_t, 0, grid.comm );
MPI_Bcast( AA.a, nnz, SuperLU_MPI_DOUBLE_COMPLEX, 0, grid.comm );
MPI_Bcast( AA.cl, nnz, mpi_int_t, 0, grid.comm );
MPI_Bcast( AA.lg, n+1, mpi_int_t, 0, grid.comm );
}
else
{
printf("\tProcess grid\t%d X %d\n", grid.nprow, grid.npcol);
/* Receive matrix A from PE 0. */
MPI_Bcast( &m, 1, mpi_int_t, 0, grid.comm );
MPI_Bcast( &n, 1, mpi_int_t, 0, grid.comm );
MPI_Bcast( &nnz, 1, mpi_int_t, 0, grid.comm );
/* Allocate storage for compressed column representation. */
// FFCS - "this->" required by g++ 4.7
zallocateA_dist(n, nnz, this->dc(&a), &asub, &xa);
MPI_Bcast( a, nnz, SuperLU_MPI_DOUBLE_COMPLEX, 0, grid.comm );
MPI_Bcast( asub, nnz, mpi_int_t, 0, grid.comm );
MPI_Bcast( xa, n+1, mpi_int_t, 0, grid.comm );
}
/* Compute the number of rows to be distributed to local process */
m_loc = m / (grid.nprow * grid.npcol);
m_loc_fst = m_loc;
/* When m / procs is not an integer */
if ((m_loc * grid.nprow * grid.npcol) != m) {
/*m_loc = m_loc+1;
m_loc_fst = m_loc;*/
if (iam == (grid.nprow * grid.npcol - 1)) /* last proc. gets all*/
m_loc = m - m_loc * (grid.nprow * grid.npcol - 1);
}
fst_row = iam * m_loc_fst;
nnz_loc = xa[fst_row+m_loc]-xa[fst_row];
// FFCS - "this->" required by g++ 4.7
zallocateA_dist(m_loc, nnz_loc, this->dc(&aloc), &asubloc, &xaloc);
//xaloc = (int_t*) intMalloc_dist(m_loc+1);
for(int ii=0; ii < m_loc; ii++){
xaloc[ii] = xa[fst_row+ii]-xa[fst_row];
}
xaloc[m_loc]=nnz_loc;
fst_nnz = xa[fst_row];
//aloc = new R[nnz_loc];
//aloc = (Complex*) doubleMalloc_dist(nnz_loc);
//asubloc = (int_t*) intMalloc_dist(nnz_loc);
for(int ii=0; ii < nnz_loc; ii++){
aloc[ii] = a[fst_nnz+ii];
asubloc[ii] = asub[fst_nnz+ii];
}
if( iam ){
SUPERLU_FREE( a );
SUPERLU_FREE( asub );
SUPERLU_FREE( xa );
}
Dtype_t R_SLU = SuperLUmpiDISTDriver<R>::R_SLU_T();
if(verbosity) cout << "Debut: Create_CompRowCol_Matrix_dist" <<endl;
// FFCS - "this->" required by g++ 4.7
if(verbosity) this->Create_CompRowLoc_Matrix_dist(&A, m, n, nnz_loc, m_loc, fst_row, aloc, asubloc, xaloc, SLU_NR_loc, R_SLU, SLU_GE);
cout << "Fin: Create_CompRowCol_Matrix_dist" <<endl;
/* creation of pseudo solution + second member */
if ( !(B = new Complex[m_loc]) ){
printf("probleme d allocation\n");
exit(1);
}
for(int ii=0; ii < m_loc; ii++){
B[ii] = 1.; //BB[fst_row+ii];
}
if(verbosity)
printf("Dimension %dx%d; # nonzeros %d\n", A.nrow, A.ncol, nnz);
/* set the default options */
set_default_options_dist(&options);
DiagScale_t optionDiagScale;
if(!string_option.empty()) read_options_freefem(&string_option,&options,&optionDiagScale);
m=A.nrow;
n=A.ncol;
printf("Dimension %dx%d; # nonzeros %d\n", A.nrow, A.ncol, nnz);
/* Initialize ScalePermstruct and LUstruct. */
ScalePermstructInit(m, n, &ScalePermstruct);
if(pperm_r || pperm_c ) ScalePermstruct.DiagScale=optionDiagScale;
if(pperm_r)
for(int ii=0; ii<m; ii++) ScalePermstruct.perm_r[ii] = pperm_r[fst_row+ii];
if(pperm_c)
for(int ii=0; ii<n; ii++) ScalePermstruct.perm_c[ii] = pperm_c[ii];
LUstructInit(m, n, &LUstruct);
/* Initialize the statistics variables. */
PStatInit(&stat);
ldb = m_loc;
//ldx = m_loc;
nrhs=1;
if ( !(berr = doubleMalloc_dist(nrhs )) ){
printf("probleme d allocation\n");
exit(1);
}
berr[0]=0.;
/* ONLY PERFORM THE LU DECOMPOSITION */
nrhs=0;
SuperLUmpiDISTDriver<R>::pgssvx(&options, &A, &ScalePermstruct, B, ldb, nrhs, &grid,
&LUstruct, &SOLVEstruct, berr, &stat, &info);
if(verbosity)
printf("LU factorization: pdgssvx()/p returns info %d\n", info);
if ( verbosity) PStatPrint(&options,&stat,&grid);
PStatFree(&stat);
}
else if( matrixdist == distributed) {
printf("in construction\n");
exit(1);
}
else{
printf("matrix choice for SuperLU_DIST is assembled, distributedglobal and distributed \n");
exit(1);
}
delete [] B;
options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
nrhs=1;
SUPERLU_FREE(berr);
if(iam==0){
finishtime = clock();
timeused= (finishtime-starttime)/(1000 );
printf("=====================================================\n");
cout << "SuperLU_DIST : time factorisation :: " << timeused << " ms" <<endl;
printf("=====================================================\n");
}
}
}
void Solver(const MatriceMorse<R> &AA,KN_<R> &x,const KN_<R> &b) const {
R* B;
SuperLUStat_t stat;
int iam;
int info=0, ldb=m, nrhs=1;
int i;
double* berr;
double ferr;
double rpg, rcond;
int_t m_loc,m_loc_fst,fst_row;
// time variables
long int starttime,finishtime;
long int timeused;
iam = grid.iam;
if( iam < nprow*npcol){
if(verbosity) starttime = clock();
if(n != m) exit(1);
ffassert ( &x[0] != &b[0]);
epsr = (eps < 0) ? (epsr >0 ? -epsr : -eps ) : eps ;
Dtype_t R_SLU = SuperLUmpiDISTDriver<R>::R_SLU_T();
nrhs= 1;
/* Initialize the statistics variables. */
PStatInit(&stat);
/* cas matrix assembled */
if( matrixdist == assembled ){
if( !(B = new R[m*nrhs] ) ){
printf("probleme d allocation\n");
exit(1);
}
for(int ii=0; ii<n; ii++){
B[ii]=b[ii];
}
if ( !(berr = doubleMalloc_dist(nrhs )) ){
printf("probleme d allocation\n");
exit(1);
}
berr[0]=0.;
options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
ldb = m;
//nrhs= 1;
SuperLUmpiDISTDriver<R>::pgssvx_ABglobal (&options, &A, &ScalePermstruct, B, ldb, nrhs, &grid,
&LUstruct, berr, &stat, &info );
if(verbosity)
printf("Triangular solve: dgssvx() returns info %d\n", info);
if(verbosity) PStatPrint(&options, &stat, &grid);
for(int ii=0; ii<n; ii++){
x[ii] = B[ii];
}
if(verbosity) cout << " x min max " << x.min() << " " <<x.max() << endl;
}
else if( matrixdist == distributedglobal) {
R* xtemp;
iam = grid.iam;
/* Compute the number of rows to be distributed to local process */
m_loc = m / (grid.nprow * grid.npcol);
m_loc_fst = m_loc;
/* When m / procs is not an integer */
if ((m_loc * grid.nprow * grid.npcol) != m) {
/*m_loc = m_loc+1;
m_loc_fst = m_loc;*/
if (iam == (grid.nprow * grid.npcol - 1)) /* last proc. gets all*/
m_loc = m - m_loc * (grid.nprow * grid.npcol - 1);
}
fst_row = iam * m_loc_fst;
if ( !(B = new R[m_loc] ) ){
printf("probleme d allocation\n");
exit(1);
}
//printf("initilisation B:");
for(int ii=0; ii<m_loc; ++ii){
B[ii] = b[ii+fst_row];
//printf(" B[%d]= %f ",ii,B[ii]);
}
//printf(" :: fin \n");
//fflush(stdout);
if ( !(berr = doubleMalloc_dist(nrhs )) ){
printf("probleme d allocation\n");
exit(1);
}
berr[0]=0.;
options.Fact = FACTORED; /* Indicate the factored form of A is supplied. */
//options.Equil = YES;
//options.Trans = TRANS;
ldb = m;
SuperLUmpiDISTDriver<R>::pgssvx(&options, &A, &ScalePermstruct, B, ldb, nrhs, &grid,
&LUstruct, &SOLVEstruct, berr, &stat, &info );
if(verbosity)
printf("Triangular solve: dgssvx() returns info %d\n", info);
if ( !(xtemp = new R[AA.n]) ){
printf("probleme d allocation de xtemp\n");
exit(1);
}
int disp[nprow*npcol];
MPI_Allgather(&fst_row, 1, MPI_INT, disp, 1, MPI_INT, grid.comm);
int recv[nprow*npcol];
MPI_Allgather(&m_loc, 1, MPI_INT, recv, 1, MPI_INT, grid.comm);
MPI_Allgatherv(B, m_loc, SuperLU_MPI_DOUBLE_COMPLEX, xtemp, recv, disp, SuperLU_MPI_DOUBLE_COMPLEX, grid.comm);
for(int ii= 0; ii< AA.n ; ii++)
x[ii] = xtemp[ii];
if(verbosity) cout << " x min max " << x.min() << " " <<x.max() << endl;
delete [] xtemp;
}
else if( matrixdist == distributed) {
printf("in construction\n");
exit(1);
}
else{
printf("matrix choice for SuperLU_DIST is assembled, distributedglobal and distributed \n");
exit(1);
}
delete [] B;
SUPERLU_FREE( berr );
PStatFree(&stat);
if(iam==0){
finishtime = clock();
timeused= (finishtime-starttime)/(1000 );
printf("=====================================================\n");
cout << " SuperLU_DIST : time solve :: " << timeused << " ms" <<endl;
printf("=====================================================\n");
}
}
}
~ZSolveSuperLUmpi() {
int iam;
iam = grid.iam;
if(iam < nprow*npcol){
if(verbosity)
cout << "~SolveSuperLUmpi Complex:" << endl;
if( matrixdist == assembled) {
//if( A.Store) Destroy_CompCol_Matrix_dist(&A);
//if( L.Store && U.Store ) {
Destroy_LU(n, &grid, &LUstruct);
ScalePermstructFree(&ScalePermstruct);
LUstructFree(&LUstruct);
//}
if ( options.SolveInitialized ) {
zSolveFinalize(&options, &SOLVEstruct);
}
}
else if( matrixdist == distributedglobal) {
if( A.Store) Destroy_CompRowLoc_Matrix_dist(&A);
Destroy_LU(n, &grid, &LUstruct);
ScalePermstructFree(&ScalePermstruct);
LUstructFree(&LUstruct);
if ( options.SolveInitialized ) {
zSolveFinalize(&options, &SOLVEstruct);
}
}
else if( matrixdist == distributed) {
printf("in construction\n");
exit(1);
}
else{
printf("matrix choice for SuperLU_DIST is assembled, distributedglobal and distributed \n");
exit(1);
}
}
printf("superlu_gridexit(&grid), %d\n",iam);
superlu_gridexit(&grid);
}
void addMatMul(const KN_<R> & x, KN_<R> & Ax) const
{
ffassert(x.N()==Ax.N());
Ax += (const MatriceMorse<R> &) (*this) * x;
}
};
MatriceMorse<Complex>::VirtualSolver *
BuildSolverSuperLUmpi(DCL_ARG_SPARSE_SOLVER(Complex,A))
{
if(verbosity>9)
cout << " BuildSolverSuperLUmpi<double>" << endl;
return new ZSolveSuperLUmpi<Complex>(*A,ds.strategy,ds.tgv,ds.epsilon,ds.tol_pivot,ds.tol_pivot_sym,
ds.data_filename, ds.sparams, ds.perm_r, ds.perm_c, ds.commworld);
}
/* --FH: class Init { public:
Init();
};*/
// the 2 default sparse solver double and complex
//DefSparseSolver<double>::SparseMatSolver SparseMatSolver_R ; ;
DefSparseSolver<Complex>::SparseMatSolver SparseMatSolver_C;
// the default probleme solver
TypeSolveMat::TSolveMat TypeSolveMatdefaultvalue=TypeSolveMat::defaultvalue;
bool SetDefault()
{
if(verbosity)
cout << " SetDefault sparse to default" << endl;
//DefSparseSolver<double>::solver =SparseMatSolver_R;
DefSparseSolver<Complex>::solver =SparseMatSolver_C;
TypeSolveMat::defaultvalue =TypeSolveMat::SparseSolver;
return false;
}
bool SetSuperLUmpi()
{
if(verbosity)
cout << " SetDefault sparse solver to SuperLUmpi double" << endl;
//DefSparseSolver<double>::solver =BuildSolverSuperLUmpi;
DefSparseSolver<Complex>::solver =BuildSolverSuperLUmpi;
TypeSolveMat::defaultvalue = TypeSolveMatdefaultvalue;
return false;
}
static void Load_Init()
{
//SparseMatSolver_R= DefSparseSolver<double>::solver;
SparseMatSolver_C= DefSparseSolver<Complex>::solver;
if(verbosity>1)
cout << "\n Add: Complex SuperLU_DIST, defaultsolver defaultsolverSuperLUdist" << endl;
TypeSolveMat::defaultvalue=TypeSolveMat::SparseSolver;
//DefSparseSolver<double>::solver =BuildSolverSuperLUmpi;
DefSparseSolver<Complex>::solver =BuildSolverSuperLUmpi;
if(! Global.Find("defaultsolver").NotNull() )
Global.Add("defaultsolver","(",new OneOperator0<bool>(SetDefault));
Global.Add("complexdefaulttoSuperLUdist","(",new OneOperator0<bool>(SetSuperLUmpi));
}
LOADFUNC(Load_Init)
|