1d24d4204SJose E. Roman 2d24d4204SJose E. Roman static char help[] = "Tests ScaLAPACK interface.\n\n"; 3d24d4204SJose E. Roman 4d24d4204SJose E. Roman #include <petscmat.h> 5d24d4204SJose E. Roman 6d24d4204SJose E. Roman int main(int argc,char **args) 7d24d4204SJose E. Roman { 8f7ec113fSDamian Marek Mat Cdense,Caij,B,C,Ct,Asub; 9d24d4204SJose E. Roman Vec d; 10d24d4204SJose E. Roman PetscInt i,j,M = 5,N,mb = 2,nb,nrows,ncols,mloc,nloc; 11d24d4204SJose E. Roman const PetscInt *rows,*cols; 12d24d4204SJose E. Roman IS isrows,iscols; 13d24d4204SJose E. Roman PetscErrorCode ierr; 14d24d4204SJose E. Roman PetscScalar *v; 15f7ec113fSDamian Marek PetscMPIInt rank,color; 16d24d4204SJose E. Roman PetscReal Cnorm; 17d24d4204SJose E. Roman PetscBool flg,mats_view=PETSC_FALSE; 18f7ec113fSDamian Marek MPI_Comm subcomm; 19d24d4204SJose E. Roman 20d24d4204SJose E. Roman ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 21*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 23d24d4204SJose E. Roman N = M; 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL)); 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mb",&mb,NULL)); 26d24d4204SJose E. Roman nb = mb; 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-nb",&nb,NULL)); 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view)); 29d24d4204SJose E. Roman 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(C,MATSCALAPACK)); 32d24d4204SJose E. Roman mloc = PETSC_DECIDE; 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&mloc,&M)); 34d24d4204SJose E. Roman nloc = PETSC_DECIDE; 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSplitOwnershipEqual(PETSC_COMM_WORLD,&nloc,&N)); 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,mloc,nloc,M,N)); 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatScaLAPACKSetBlockSizes(C,mb,nb)); 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 40*5f80ce2aSJacob Faibussowitsch /*CHKERRQ(MatCreateScaLAPACK(PETSC_COMM_WORLD,mb,nb,M,N,0,0,&C)); */ 41d24d4204SJose E. Roman 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipIS(C,&isrows,&iscols)); 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(isrows,&nrows)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isrows,&rows)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(iscols,&ncols)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(iscols,&cols)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nrows*ncols,&v)); 48d24d4204SJose E. Roman for (i=0;i<nrows;i++) { 49d24d4204SJose E. Roman for (j=0;j<ncols;j++) v[i*ncols+j] = (PetscReal)(rows[i]+1+(cols[j]+1)*0.01); 50d24d4204SJose E. Roman } 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,nrows,rows,ncols,cols,v,INSERT_VALUES)); 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 53*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isrows,&rows)); 54*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(iscols,&cols)); 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrows)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscols)); 59d24d4204SJose E. Roman 60d24d4204SJose E. Roman /* Test MatView(), MatDuplicate() and out-of-place MatConvert() */ 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(C,MAT_COPY_VALUES,&B)); 62d24d4204SJose E. Roman if (mats_view) { 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Duplicated C:\n")); 64*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 65d24d4204SJose E. Roman } 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(C,MATDENSE,MAT_INITIAL_MATRIX,&Cdense)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(C,Cdense,5,&flg)); 692c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Check fails: Cdense != C"); 70d24d4204SJose E. Roman 71d24d4204SJose E. Roman /* Test MatNorm() */ 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C,NORM_1,&Cnorm)); 73d24d4204SJose E. Roman 74d24d4204SJose E. Roman /* Test MatTranspose(), MatZeroEntries() and MatGetDiagonal() */ 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatTranspose(C,MAT_INITIAL_MATRIX,&Ct)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConjugate(Ct)); 77d24d4204SJose E. Roman if (mats_view) { 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"C's Transpose Conjugate:\n")); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(Ct,PETSC_VIEWER_STDOUT_WORLD)); 80d24d4204SJose E. Roman } 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroEntries(Ct)); 82*5f80ce2aSJacob Faibussowitsch if (M>N) CHKERRQ(MatCreateVecs(C,&d,NULL)); 83*5f80ce2aSJacob Faibussowitsch else CHKERRQ(MatCreateVecs(C,NULL,&d)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(C,d)); 85d24d4204SJose E. Roman if (mats_view) { 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Diagonal of C:\n")); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(d,PETSC_VIEWER_STDOUT_WORLD)); 88d24d4204SJose E. Roman } 89d24d4204SJose E. Roman if (M>N) { 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(C,NULL,d)); 91d24d4204SJose E. Roman } else { 92*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(C,d,NULL)); 93d24d4204SJose E. Roman } 94d24d4204SJose E. Roman if (mats_view) { 95*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Diagonal Scaled C:\n")); 96*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 97d24d4204SJose E. Roman } 98d24d4204SJose E. Roman 99d24d4204SJose E. Roman /* Test MatAXPY(), MatAYPX() and in-place MatConvert() */ 100*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B,MATSCALAPACK)); 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B,mloc,nloc,PETSC_DECIDE,PETSC_DECIDE)); 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatScaLAPACKSetBlockSizes(B,mb,nb)); 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(B)); 106*5f80ce2aSJacob Faibussowitsch /* CHKERRQ(MatCreateScaLAPACK(PETSC_COMM_WORLD,mb,nb,M,N,0,0,&B)); */ 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipIS(B,&isrows,&iscols)); 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(isrows,&nrows)); 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isrows,&rows)); 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(iscols,&ncols)); 111*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(iscols,&cols)); 112*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nrows*ncols,&v)); 113d24d4204SJose E. Roman for (i=0;i<nrows;i++) { 114d24d4204SJose E. Roman for (j=0;j<ncols;j++) v[i*ncols+j] = (PetscReal)(1000*rows[i]+cols[j]); 115d24d4204SJose E. Roman } 116*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES)); 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isrows,&rows)); 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(iscols,&cols)); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 122d24d4204SJose E. Roman if (mats_view) { 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"B:\n")); 124*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 125d24d4204SJose E. Roman } 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(B,2.5,C,SAME_NONZERO_PATTERN)); 127*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAYPX(B,3.75,C,SAME_NONZERO_PATTERN)); 128*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatConvert(B,MATDENSE,MAT_INPLACE_MATRIX,&B)); 129d24d4204SJose E. Roman if (mats_view) { 130*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"B after MatAXPY and MatAYPX:\n")); 131*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 132d24d4204SJose E. Roman } 133*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrows)); 134*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscols)); 135*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 136d24d4204SJose E. Roman 137d24d4204SJose E. Roman /* Test MatMatTransposeMult(): B = C*C^T */ 138*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMult(C,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B)); 139*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatScale(C,2.0)); 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMult(C,C,MAT_REUSE_MATRIX,PETSC_DEFAULT,&B)); 141*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatTransposeMultEqual(C,C,B,10,&flg)); 1422c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Check fails: B != C*C^T"); 143d24d4204SJose E. Roman 144d24d4204SJose E. Roman if (mats_view) { 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"C MatMatTransposeMult C:\n")); 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 147d24d4204SJose E. Roman } 148d24d4204SJose E. Roman 149d24d4204SJose E. Roman /* Test MatMult() */ 150*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatComputeOperator(C,MATAIJ,&Caij)); 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(C,Caij,5,&flg)); 1522c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultEqual() fails"); 153*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeEqual(C,Caij,5,&flg)); 1542c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultTransposeEqual() fails"); 155d24d4204SJose E. Roman 156d24d4204SJose E. Roman /* Test MatMultAdd() and MatMultTransposeAddEqual() */ 157*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultAddEqual(C,Caij,5,&flg)); 1582c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultAddEqual() fails"); 159*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultTransposeAddEqual(C,Caij,5,&flg)); 1602c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"C != Caij. MatMultTransposeAddEqual() fails"); 161d24d4204SJose E. Roman 162d24d4204SJose E. Roman /* Test MatMatMult() */ 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-test_matmatmult",&flg)); 164d24d4204SJose E. Roman if (flg) { 165d24d4204SJose E. Roman Mat CC,CCaij; 166*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(C,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CC)); 167*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(Caij,Caij,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&CCaij)); 168*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatMultEqual(CC,CCaij,5,&flg)); 1692c71b3e2SJacob Faibussowitsch PetscCheckFalse(!flg,PETSC_COMM_WORLD,PETSC_ERR_ARG_NOTSAMETYPE,"CC != CCaij. MatMatMult() fails"); 170*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&CCaij)); 171*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&CC)); 172d24d4204SJose E. Roman } 173d24d4204SJose E. Roman 174f7ec113fSDamian Marek /* Test MatCreate() on subcomm */ 175f7ec113fSDamian Marek color = rank%2; 176*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_split(PETSC_COMM_WORLD,color,0,&subcomm)); 177f7ec113fSDamian Marek if (color==0) { 178*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(subcomm,&Asub)); 179*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(Asub,MATSCALAPACK)); 180f7ec113fSDamian Marek mloc = PETSC_DECIDE; 181*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSplitOwnershipEqual(subcomm,&mloc,&M)); 182f7ec113fSDamian Marek nloc = PETSC_DECIDE; 183*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSplitOwnershipEqual(subcomm,&nloc,&N)); 184*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(Asub,mloc,nloc,M,N)); 185*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatScaLAPACKSetBlockSizes(Asub,mb,nb)); 186*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(Asub)); 187*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(Asub)); 188*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Asub)); 189f7ec113fSDamian Marek } 190f7ec113fSDamian Marek 191*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Cdense)); 192*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Caij)); 193*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 194*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 195*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Ct)); 196*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&d)); 197*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_free(&subcomm)); 198d24d4204SJose E. Roman ierr = PetscFinalize(); 199d24d4204SJose E. Roman return ierr; 200d24d4204SJose E. Roman } 201d24d4204SJose E. Roman 202d24d4204SJose E. Roman /*TEST 203d24d4204SJose E. Roman 204d24d4204SJose E. Roman build: 205d24d4204SJose E. Roman requires: scalapack 206d24d4204SJose E. Roman 207d24d4204SJose E. Roman test: 208d24d4204SJose E. Roman nsize: 2 209d24d4204SJose E. Roman args: -mb 5 -nb 5 -M 12 -N 10 210d24d4204SJose E. Roman requires: scalapack 211d24d4204SJose E. Roman 212d24d4204SJose E. Roman test: 213d24d4204SJose E. Roman suffix: 2 214d24d4204SJose E. Roman nsize: 6 215d24d4204SJose E. Roman args: -mb 8 -nb 6 -M 20 -N 50 216d24d4204SJose E. Roman requires: scalapack 217d24d4204SJose E. Roman output_file: output/ex242_1.out 218d24d4204SJose E. Roman 219d24d4204SJose E. Roman test: 220d24d4204SJose E. Roman suffix: 3 221d24d4204SJose E. Roman nsize: 3 222d24d4204SJose E. Roman args: -mb 2 -nb 2 -M 20 -N 20 -test_matmatmult 223d24d4204SJose E. Roman requires: scalapack 224d24d4204SJose E. Roman output_file: output/ex242_1.out 225d24d4204SJose E. Roman 226d24d4204SJose E. Roman TEST*/ 227