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