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*327415f7SBarry Smith PetscFunctionBeginUser; 209566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*) 0,help)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 229566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 23d24d4204SJose E. Roman 249566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 259566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 26d24d4204SJose E. Roman 27d24d4204SJose E. Roman /* Get local dimensions of matrices */ 289566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 29d24d4204SJose E. Roman n = m; 309566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 31d24d4204SJose E. Roman p = m/2; 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL)); 339566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view)); 34d24d4204SJose E. Roman 35d24d4204SJose E. Roman /* Create matrix A */ 369566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Create ScaLAPACK matrix A\n")); 379566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 389566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE)); 399566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATSCALAPACK)); 409566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 419566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 42d24d4204SJose E. Roman /* Set local matrix entries */ 439566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(A,&isrows,&iscols)); 449566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows,&nrows)); 459566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows,&rows)); 469566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols,&ncols)); 479566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscols,&cols)); 489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows*ncols,&v)); 49d24d4204SJose E. Roman for (i=0;i<nrows;i++) { 50d24d4204SJose E. Roman for (j=0;j<ncols;j++) { 519566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 52d24d4204SJose E. Roman v[i*ncols+j] = rval; 53d24d4204SJose E. Roman } 54d24d4204SJose E. Roman } 559566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES)); 569566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 579566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 589566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows,&rows)); 599566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscols,&cols)); 609566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 619566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 629566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 63d24d4204SJose E. Roman if (mats_view) { 649566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "A: nrows %" PetscInt_FMT ", m %" PetscInt_FMT "; ncols %" PetscInt_FMT ", n %" PetscInt_FMT "\n",nrows,m,ncols,n)); 659566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 66d24d4204SJose E. Roman } 67d24d4204SJose E. Roman 68d24d4204SJose E. Roman /* Create rhs matrix B */ 699566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Create rhs matrix B\n")); 709566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); 719566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,m,p,PETSC_DECIDE,PETSC_DECIDE)); 729566063dSJacob Faibussowitsch PetscCall(MatSetType(B,MATSCALAPACK)); 739566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 749566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 759566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(B,&isrows,&iscols)); 769566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows,&nrows)); 779566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows,&rows)); 789566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols,&ncols)); 799566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscols,&cols)); 809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows*ncols,&v)); 81d24d4204SJose E. Roman for (i=0;i<nrows;i++) { 82d24d4204SJose E. Roman for (j=0;j<ncols;j++) { 839566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 84d24d4204SJose E. Roman v[i*ncols+j] = rval; 85d24d4204SJose E. Roman } 86d24d4204SJose E. Roman } 879566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES)); 889566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 899566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 909566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows,&rows)); 919566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscols,&cols)); 929566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 939566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 949566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 95d24d4204SJose E. Roman if (mats_view) { 969566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "B: nrows %" PetscInt_FMT ", m %" PetscInt_FMT "; ncols %" PetscInt_FMT ", p %" PetscInt_FMT "\n",nrows,m,ncols,p)); 979566063dSJacob Faibussowitsch PetscCall(MatView(B,PETSC_VIEWER_STDOUT_WORLD)); 98d24d4204SJose E. Roman } 99d24d4204SJose E. Roman 100d24d4204SJose E. Roman /* Create rhs vector b and solution x (same size as b) */ 1019566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&b)); 1029566063dSJacob Faibussowitsch PetscCall(VecSetSizes(b,m,PETSC_DECIDE)); 1039566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(b)); 1049566063dSJacob Faibussowitsch PetscCall(VecGetArray(b,&barray)); 105d24d4204SJose E. Roman for (j=0;j<m;j++) { 1069566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 107d24d4204SJose E. Roman barray[j] = rval; 108d24d4204SJose E. Roman } 1099566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(b,&barray)); 1109566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(b)); 1119566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(b)); 112d24d4204SJose E. Roman if (mats_view) { 1139566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] b: m %" PetscInt_FMT "\n",rank,m)); 1149566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 1159566063dSJacob Faibussowitsch PetscCall(VecView(b,PETSC_VIEWER_STDOUT_WORLD)); 116d24d4204SJose E. Roman } 1179566063dSJacob Faibussowitsch PetscCall(VecDuplicate(b,&x)); 118d24d4204SJose E. Roman 119d24d4204SJose E. Roman /* Create matrix X - same size as B */ 1209566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Create solution matrix X\n")); 1219566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&X)); 122d24d4204SJose E. Roman 123d24d4204SJose E. Roman /* Cholesky factorization */ 124d24d4204SJose E. Roman /*------------------------*/ 1259566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Create ScaLAPACK matrix Aher\n")); 1269566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&Aher)); 1279566063dSJacob Faibussowitsch PetscCall(MatAXPY(Aher,1.0,A,SAME_NONZERO_PATTERN)); /* Aher = A + A^T */ 1289566063dSJacob Faibussowitsch PetscCall(MatShift(Aher,100.0)); /* add 100.0 to diagonals of Aher to make it spd */ 129d24d4204SJose E. Roman if (mats_view) { 1309566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Aher:\n")); 1319566063dSJacob Faibussowitsch PetscCall(MatView(Aher,PETSC_VIEWER_STDOUT_WORLD)); 132d24d4204SJose E. Roman } 133d24d4204SJose E. Roman 134d24d4204SJose E. Roman /* Cholesky factorization */ 135d24d4204SJose E. Roman /*------------------------*/ 1369566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Test Cholesky Solver \n")); 137d24d4204SJose E. Roman /* In-place Cholesky */ 138d24d4204SJose E. Roman /* Create matrix factor G, with a copy of Aher */ 1399566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Aher,MAT_COPY_VALUES,&G)); 140d24d4204SJose E. Roman 141d24d4204SJose E. Roman /* G = L * L^T */ 1429566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactor(G,0,0)); 143d24d4204SJose E. Roman if (mats_view) { 1449566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Cholesky Factor G:\n")); 1459566063dSJacob Faibussowitsch PetscCall(MatView(G,PETSC_VIEWER_STDOUT_WORLD)); 146d24d4204SJose E. Roman } 147d24d4204SJose E. Roman 148d24d4204SJose E. Roman /* Solve L * L^T x = b and L * L^T * X = B */ 1499566063dSJacob Faibussowitsch PetscCall(MatSolve(G,b,x)); 1509566063dSJacob Faibussowitsch PetscCall(MatMatSolve(G,B,X)); 1519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&G)); 152d24d4204SJose E. Roman 153d24d4204SJose E. Roman /* Out-place Cholesky */ 1549566063dSJacob Faibussowitsch PetscCall(MatGetFactor(Aher,MATSOLVERSCALAPACK,MAT_FACTOR_CHOLESKY,&G)); 1559566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(G,Aher,0,NULL)); 1569566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(G,Aher,NULL)); 1571baa6e33SBarry Smith if (mats_view) PetscCall(MatView(G,PETSC_VIEWER_STDOUT_WORLD)); 1589566063dSJacob Faibussowitsch PetscCall(MatSolve(G,b,x)); 1599566063dSJacob Faibussowitsch PetscCall(MatMatSolve(G,B,X)); 1609566063dSJacob Faibussowitsch PetscCall(MatDestroy(&G)); 161d24d4204SJose E. Roman 162d24d4204SJose E. Roman /* Check norm(Aher*x - b) */ 1639566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&c)); 1649566063dSJacob Faibussowitsch PetscCall(VecSetSizes(c,m,PETSC_DECIDE)); 1659566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(c)); 1669566063dSJacob Faibussowitsch PetscCall(MatMult(Aher,x,c)); 1679566063dSJacob Faibussowitsch PetscCall(VecAXPY(c,-1.0,b)); 1689566063dSJacob Faibussowitsch PetscCall(VecNorm(c,NORM_1,&norm)); 169d24d4204SJose E. Roman if (norm > tol) { 1709566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: ||Aher*x - b||=%g for Cholesky\n",(double)norm)); 171d24d4204SJose E. Roman } 172d24d4204SJose E. Roman 173d24d4204SJose E. Roman /* Check norm(Aher*X - B) */ 1749566063dSJacob Faibussowitsch PetscCall(MatMatMult(Aher,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); 1759566063dSJacob Faibussowitsch PetscCall(MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN)); 1769566063dSJacob Faibussowitsch PetscCall(MatNorm(C,NORM_1,&norm)); 177d24d4204SJose E. Roman if (norm > tol) { 1789566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: ||Aher*X - B||=%g for Cholesky\n",(double)norm)); 179d24d4204SJose E. Roman } 180d24d4204SJose E. Roman 181d24d4204SJose E. Roman /* LU factorization */ 182d24d4204SJose E. Roman /*------------------*/ 1839566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Test LU Solver \n")); 184d24d4204SJose E. Roman /* In-place LU */ 185d24d4204SJose E. Roman /* Create matrix factor F, with a copy of A */ 1869566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&F)); 187d24d4204SJose E. Roman /* Create vector d to test MatSolveAdd() */ 1889566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&d)); 1899566063dSJacob Faibussowitsch PetscCall(VecCopy(x,d)); 190d24d4204SJose E. Roman 191d24d4204SJose E. Roman /* PF=LU factorization */ 1929566063dSJacob Faibussowitsch PetscCall(MatLUFactor(F,0,0,NULL)); 193d24d4204SJose E. Roman 194d24d4204SJose E. Roman /* Solve LUX = PB */ 1959566063dSJacob Faibussowitsch PetscCall(MatSolveAdd(F,b,d,x)); 1969566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F,B,X)); 1979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 198d24d4204SJose E. Roman 199d24d4204SJose E. Roman /* Check norm(A*X - B) */ 2009566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&e)); 2019566063dSJacob Faibussowitsch PetscCall(VecSetSizes(e,m,PETSC_DECIDE)); 2029566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(e)); 2039566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,c)); 2049566063dSJacob Faibussowitsch PetscCall(MatMult(A,d,e)); 2059566063dSJacob Faibussowitsch PetscCall(VecAXPY(c,-1.0,e)); 2069566063dSJacob Faibussowitsch PetscCall(VecAXPY(c,-1.0,b)); 2079566063dSJacob Faibussowitsch PetscCall(VecNorm(c,NORM_1,&norm)); 208d24d4204SJose E. Roman if (norm > tol) { 2099566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: ||A*x - b||=%g for LU\n",(double)norm)); 210d24d4204SJose E. Roman } 211d24d4204SJose E. Roman /* Reuse product C; replace Aher with A */ 2129566063dSJacob Faibussowitsch PetscCall(MatProductReplaceMats(A,NULL,NULL,C)); 2139566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C)); 2149566063dSJacob Faibussowitsch PetscCall(MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN)); 2159566063dSJacob Faibussowitsch PetscCall(MatNorm(C,NORM_1,&norm)); 216d24d4204SJose E. Roman if (norm > tol) { 2179566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: ||A*X - B||=%g for LU\n",(double)norm)); 218d24d4204SJose E. Roman } 219d24d4204SJose E. Roman 220d24d4204SJose E. Roman /* Out-place LU */ 2219566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A,MATSOLVERSCALAPACK,MAT_FACTOR_LU,&F)); 2229566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F,A,0,0,NULL)); 2239566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F,A,NULL)); 2241baa6e33SBarry Smith if (mats_view) PetscCall(MatView(F,PETSC_VIEWER_STDOUT_WORLD)); 2259566063dSJacob Faibussowitsch PetscCall(MatSolve(F,b,x)); 2269566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F,B,X)); 2279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 228d24d4204SJose E. Roman 229d24d4204SJose E. Roman /* Free space */ 2309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Aher)); 2329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 2349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&X)); 2359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 2369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&c)); 2379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d)); 2389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&e)); 2399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2409566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 2419566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 242b122ec5aSJacob Faibussowitsch return 0; 243d24d4204SJose E. Roman } 244d24d4204SJose E. Roman 245d24d4204SJose E. Roman /*TEST 246d24d4204SJose E. Roman 247d24d4204SJose E. Roman build: 248d24d4204SJose E. Roman requires: scalapack 249d24d4204SJose E. Roman 250d24d4204SJose E. Roman test: 251d24d4204SJose E. Roman nsize: 2 252d24d4204SJose E. Roman output_file: output/ex245.out 253d24d4204SJose E. Roman 254d24d4204SJose E. Roman test: 255d24d4204SJose E. Roman suffix: 2 256d24d4204SJose E. Roman nsize: 6 257d24d4204SJose E. Roman output_file: output/ex245.out 258d24d4204SJose E. Roman 259d24d4204SJose E. Roman TEST*/ 260