1d24d4204SJose E. Roman 2d24d4204SJose E. Roman static char help[] = "Tests LU, Cholesky factorization and MatMatSolve() for a ScaLAPACK dense matrix.\n\n"; 3d24d4204SJose E. Roman 4d24d4204SJose E. Roman #include <petscmat.h> 5d24d4204SJose E. Roman 6d24d4204SJose E. Roman int main(int argc,char **argv) 7d24d4204SJose E. Roman { 8d24d4204SJose E. Roman Mat A,F,B,X,C,Aher,G; 9d24d4204SJose E. Roman Vec b,x,c,d,e; 10d24d4204SJose E. Roman PetscInt m=5,n,p,i,j,nrows,ncols; 11d24d4204SJose E. Roman PetscScalar *v,*barray,rval; 12d24d4204SJose E. Roman PetscReal norm,tol=1.e5*PETSC_MACHINE_EPSILON; 13d24d4204SJose E. Roman PetscMPIInt size,rank; 14d24d4204SJose E. Roman PetscRandom rand; 15d24d4204SJose E. Roman const PetscInt *rows,*cols; 16d24d4204SJose E. Roman IS isrows,iscols; 17d24d4204SJose E. Roman PetscBool mats_view=PETSC_FALSE; 18d24d4204SJose E. Roman 19*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*) 0,help)); 20*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 21*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 22d24d4204SJose E. Roman 23*9566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 24*9566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 25d24d4204SJose E. Roman 26d24d4204SJose E. Roman /* Get local dimensions of matrices */ 27*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 28d24d4204SJose E. Roman n = m; 29*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 30d24d4204SJose E. Roman p = m/2; 31*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL)); 32*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view)); 33d24d4204SJose E. Roman 34d24d4204SJose E. Roman /* Create matrix A */ 35*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Create ScaLAPACK matrix A\n")); 36*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 37*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE)); 38*9566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATSCALAPACK)); 39*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 40*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 41d24d4204SJose E. Roman /* Set local matrix entries */ 42*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(A,&isrows,&iscols)); 43*9566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows,&nrows)); 44*9566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows,&rows)); 45*9566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols,&ncols)); 46*9566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscols,&cols)); 47*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows*ncols,&v)); 48d24d4204SJose E. Roman for (i=0;i<nrows;i++) { 49d24d4204SJose E. Roman for (j=0;j<ncols;j++) { 50*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 51d24d4204SJose E. Roman v[i*ncols+j] = rval; 52d24d4204SJose E. Roman } 53d24d4204SJose E. Roman } 54*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES)); 55*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 56*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 57*9566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows,&rows)); 58*9566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscols,&cols)); 59*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 60*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 61*9566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 62d24d4204SJose E. Roman if (mats_view) { 63*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "A: nrows %" PetscInt_FMT ", m %" PetscInt_FMT "; ncols %" PetscInt_FMT ", n %" PetscInt_FMT "\n",nrows,m,ncols,n)); 64*9566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 65d24d4204SJose E. Roman } 66d24d4204SJose E. Roman 67d24d4204SJose E. Roman /* Create rhs matrix B */ 68*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Create rhs matrix B\n")); 69*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); 70*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,m,p,PETSC_DECIDE,PETSC_DECIDE)); 71*9566063dSJacob Faibussowitsch PetscCall(MatSetType(B,MATSCALAPACK)); 72*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 73*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 74*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(B,&isrows,&iscols)); 75*9566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows,&nrows)); 76*9566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows,&rows)); 77*9566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols,&ncols)); 78*9566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscols,&cols)); 79*9566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows*ncols,&v)); 80d24d4204SJose E. Roman for (i=0;i<nrows;i++) { 81d24d4204SJose E. Roman for (j=0;j<ncols;j++) { 82*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 83d24d4204SJose E. Roman v[i*ncols+j] = rval; 84d24d4204SJose E. Roman } 85d24d4204SJose E. Roman } 86*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES)); 87*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 88*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 89*9566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows,&rows)); 90*9566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscols,&cols)); 91*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 92*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 93*9566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 94d24d4204SJose E. Roman if (mats_view) { 95*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "B: nrows %" PetscInt_FMT ", m %" PetscInt_FMT "; ncols %" PetscInt_FMT ", p %" PetscInt_FMT "\n",nrows,m,ncols,p)); 96*9566063dSJacob Faibussowitsch PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 97d24d4204SJose E. Roman } 98d24d4204SJose E. Roman 99d24d4204SJose E. Roman /* Create rhs vector b and solution x (same size as b) */ 100*9566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&b)); 101*9566063dSJacob Faibussowitsch PetscCall(VecSetSizes(b,m,PETSC_DECIDE)); 102*9566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(b)); 103*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(b,&barray)); 104d24d4204SJose E. Roman for (j=0;j<m;j++) { 105*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 106d24d4204SJose E. Roman barray[j] = rval; 107d24d4204SJose E. Roman } 108*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(b,&barray)); 109*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(b)); 110*9566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(b)); 111d24d4204SJose E. Roman if (mats_view) { 112*9566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] b: m %" PetscInt_FMT "\n",rank,m)); 113*9566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 114*9566063dSJacob Faibussowitsch PetscCall(VecView(b,PETSC_VIEWER_STDOUT_WORLD)); 115d24d4204SJose E. Roman } 116*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(b,&x)); 117d24d4204SJose E. Roman 118d24d4204SJose E. Roman /* Create matrix X - same size as B */ 119*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Create solution matrix X\n")); 120*9566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&X)); 121d24d4204SJose E. Roman 122d24d4204SJose E. Roman /* Cholesky factorization */ 123d24d4204SJose E. Roman /*------------------------*/ 124*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Create ScaLAPACK matrix Aher\n")); 125*9566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&Aher)); 126*9566063dSJacob Faibussowitsch PetscCall(MatAXPY(Aher,1.0,A,SAME_NONZERO_PATTERN)); /* Aher = A + A^T */ 127*9566063dSJacob Faibussowitsch PetscCall(MatShift(Aher,100.0)); /* add 100.0 to diagonals of Aher to make it spd */ 128d24d4204SJose E. Roman if (mats_view) { 129*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Aher:\n")); 130*9566063dSJacob Faibussowitsch PetscCall(MatView(Aher,PETSC_VIEWER_STDOUT_WORLD)); 131d24d4204SJose E. Roman } 132d24d4204SJose E. Roman 133d24d4204SJose E. Roman /* Cholesky factorization */ 134d24d4204SJose E. Roman /*------------------------*/ 135*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Test Cholesky Solver \n")); 136d24d4204SJose E. Roman /* In-place Cholesky */ 137d24d4204SJose E. Roman /* Create matrix factor G, with a copy of Aher */ 138*9566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Aher,MAT_COPY_VALUES,&G)); 139d24d4204SJose E. Roman 140d24d4204SJose E. Roman /* G = L * L^T */ 141*9566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactor(G,0,0)); 142d24d4204SJose E. Roman if (mats_view) { 143*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Cholesky Factor G:\n")); 144*9566063dSJacob Faibussowitsch PetscCall(MatView(G,PETSC_VIEWER_STDOUT_WORLD)); 145d24d4204SJose E. Roman } 146d24d4204SJose E. Roman 147d24d4204SJose E. Roman /* Solve L * L^T x = b and L * L^T * X = B */ 148*9566063dSJacob Faibussowitsch PetscCall(MatSolve(G,b,x)); 149*9566063dSJacob Faibussowitsch PetscCall(MatMatSolve(G,B,X)); 150*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&G)); 151d24d4204SJose E. Roman 152d24d4204SJose E. Roman /* Out-place Cholesky */ 153*9566063dSJacob Faibussowitsch PetscCall(MatGetFactor(Aher,MATSOLVERSCALAPACK,MAT_FACTOR_CHOLESKY,&G)); 154*9566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(G,Aher,0,NULL)); 155*9566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(G,Aher,NULL)); 156d24d4204SJose E. Roman if (mats_view) { 157*9566063dSJacob Faibussowitsch PetscCall(MatView(G,PETSC_VIEWER_STDOUT_WORLD)); 158d24d4204SJose E. Roman } 159*9566063dSJacob Faibussowitsch PetscCall(MatSolve(G,b,x)); 160*9566063dSJacob Faibussowitsch PetscCall(MatMatSolve(G,B,X)); 161*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&G)); 162d24d4204SJose E. Roman 163d24d4204SJose E. Roman /* Check norm(Aher*x - b) */ 164*9566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&c)); 165*9566063dSJacob Faibussowitsch PetscCall(VecSetSizes(c,m,PETSC_DECIDE)); 166*9566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(c)); 167*9566063dSJacob Faibussowitsch PetscCall(MatMult(Aher,x,c)); 168*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(c,-1.0,b)); 169*9566063dSJacob Faibussowitsch PetscCall(VecNorm(c,NORM_1,&norm)); 170d24d4204SJose E. Roman if (norm > tol) { 171*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: ||Aher*x - b||=%g for Cholesky\n",(double)norm)); 172d24d4204SJose E. Roman } 173d24d4204SJose E. Roman 174d24d4204SJose E. Roman /* Check norm(Aher*X - B) */ 175*9566063dSJacob Faibussowitsch PetscCall(MatMatMult(Aher,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); 176*9566063dSJacob Faibussowitsch PetscCall(MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN)); 177*9566063dSJacob Faibussowitsch PetscCall(MatNorm(C,NORM_1,&norm)); 178d24d4204SJose E. Roman if (norm > tol) { 179*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: ||Aher*X - B||=%g for Cholesky\n",(double)norm)); 180d24d4204SJose E. Roman } 181d24d4204SJose E. Roman 182d24d4204SJose E. Roman /* LU factorization */ 183d24d4204SJose E. Roman /*------------------*/ 184*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Test LU Solver \n")); 185d24d4204SJose E. Roman /* In-place LU */ 186d24d4204SJose E. Roman /* Create matrix factor F, with a copy of A */ 187*9566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&F)); 188d24d4204SJose E. Roman /* Create vector d to test MatSolveAdd() */ 189*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&d)); 190*9566063dSJacob Faibussowitsch PetscCall(VecCopy(x,d)); 191d24d4204SJose E. Roman 192d24d4204SJose E. Roman /* PF=LU factorization */ 193*9566063dSJacob Faibussowitsch PetscCall(MatLUFactor(F,0,0,NULL)); 194d24d4204SJose E. Roman 195d24d4204SJose E. Roman /* Solve LUX = PB */ 196*9566063dSJacob Faibussowitsch PetscCall(MatSolveAdd(F,b,d,x)); 197*9566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F,B,X)); 198*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 199d24d4204SJose E. Roman 200d24d4204SJose E. Roman /* Check norm(A*X - B) */ 201*9566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&e)); 202*9566063dSJacob Faibussowitsch PetscCall(VecSetSizes(e,m,PETSC_DECIDE)); 203*9566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(e)); 204*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,c)); 205*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,d,e)); 206*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(c,-1.0,e)); 207*9566063dSJacob Faibussowitsch PetscCall(VecAXPY(c,-1.0,b)); 208*9566063dSJacob Faibussowitsch PetscCall(VecNorm(c,NORM_1,&norm)); 209d24d4204SJose E. Roman if (norm > tol) { 210*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: ||A*x - b||=%g for LU\n",(double)norm)); 211d24d4204SJose E. Roman } 212d24d4204SJose E. Roman /* Reuse product C; replace Aher with A */ 213*9566063dSJacob Faibussowitsch PetscCall(MatProductReplaceMats(A,NULL,NULL,C)); 214*9566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C)); 215*9566063dSJacob Faibussowitsch PetscCall(MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN)); 216*9566063dSJacob Faibussowitsch PetscCall(MatNorm(C,NORM_1,&norm)); 217d24d4204SJose E. Roman if (norm > tol) { 218*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: ||A*X - B||=%g for LU\n",(double)norm)); 219d24d4204SJose E. Roman } 220d24d4204SJose E. Roman 221d24d4204SJose E. Roman /* Out-place LU */ 222*9566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A,MATSOLVERSCALAPACK,MAT_FACTOR_LU,&F)); 223*9566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F,A,0,0,NULL)); 224*9566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F,A,NULL)); 225d24d4204SJose E. Roman if (mats_view) { 226*9566063dSJacob Faibussowitsch PetscCall(MatView(F,PETSC_VIEWER_STDOUT_WORLD)); 227d24d4204SJose E. Roman } 228*9566063dSJacob Faibussowitsch PetscCall(MatSolve(F,b,x)); 229*9566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F,B,X)); 230*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 231d24d4204SJose E. Roman 232d24d4204SJose E. Roman /* Free space */ 233*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 234*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Aher)); 235*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 236*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 237*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&X)); 238*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 239*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&c)); 240*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d)); 241*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&e)); 242*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 243*9566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 244*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 245b122ec5aSJacob Faibussowitsch return 0; 246d24d4204SJose E. Roman } 247d24d4204SJose E. Roman 248d24d4204SJose E. Roman /*TEST 249d24d4204SJose E. Roman 250d24d4204SJose E. Roman build: 251d24d4204SJose E. Roman requires: scalapack 252d24d4204SJose E. Roman 253d24d4204SJose E. Roman test: 254d24d4204SJose E. Roman nsize: 2 255d24d4204SJose E. Roman output_file: output/ex245.out 256d24d4204SJose E. Roman 257d24d4204SJose E. Roman test: 258d24d4204SJose E. Roman suffix: 2 259d24d4204SJose E. Roman nsize: 6 260d24d4204SJose E. Roman output_file: output/ex245.out 261d24d4204SJose E. Roman 262d24d4204SJose E. Roman TEST*/ 263