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 199566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*) 0,help)); 209566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 22d24d4204SJose E. Roman 239566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 249566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rand)); 25d24d4204SJose E. Roman 26d24d4204SJose E. Roman /* Get local dimensions of matrices */ 279566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 28d24d4204SJose E. Roman n = m; 299566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 30d24d4204SJose E. Roman p = m/2; 319566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL)); 329566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view)); 33d24d4204SJose E. Roman 34d24d4204SJose E. Roman /* Create matrix A */ 359566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Create ScaLAPACK matrix A\n")); 369566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 379566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE)); 389566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATSCALAPACK)); 399566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 409566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 41d24d4204SJose E. Roman /* Set local matrix entries */ 429566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(A,&isrows,&iscols)); 439566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows,&nrows)); 449566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows,&rows)); 459566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols,&ncols)); 469566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscols,&cols)); 479566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows*ncols,&v)); 48d24d4204SJose E. Roman for (i=0;i<nrows;i++) { 49d24d4204SJose E. Roman for (j=0;j<ncols;j++) { 509566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 51d24d4204SJose E. Roman v[i*ncols+j] = rval; 52d24d4204SJose E. Roman } 53d24d4204SJose E. Roman } 549566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES)); 559566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 569566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 579566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows,&rows)); 589566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscols,&cols)); 599566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 609566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 619566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 62d24d4204SJose E. Roman if (mats_view) { 639566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "A: nrows %" PetscInt_FMT ", m %" PetscInt_FMT "; ncols %" PetscInt_FMT ", n %" PetscInt_FMT "\n",nrows,m,ncols,n)); 649566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 65d24d4204SJose E. Roman } 66d24d4204SJose E. Roman 67d24d4204SJose E. Roman /* Create rhs matrix B */ 689566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Create rhs matrix B\n")); 699566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); 709566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B,m,p,PETSC_DECIDE,PETSC_DECIDE)); 719566063dSJacob Faibussowitsch PetscCall(MatSetType(B,MATSCALAPACK)); 729566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(B)); 739566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 749566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipIS(B,&isrows,&iscols)); 759566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isrows,&nrows)); 769566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isrows,&rows)); 779566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(iscols,&ncols)); 789566063dSJacob Faibussowitsch PetscCall(ISGetIndices(iscols,&cols)); 799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrows*ncols,&v)); 80d24d4204SJose E. Roman for (i=0;i<nrows;i++) { 81d24d4204SJose E. Roman for (j=0;j<ncols;j++) { 829566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 83d24d4204SJose E. Roman v[i*ncols+j] = rval; 84d24d4204SJose E. Roman } 85d24d4204SJose E. Roman } 869566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES)); 879566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 889566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 899566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isrows,&rows)); 909566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(iscols,&cols)); 919566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrows)); 929566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscols)); 939566063dSJacob Faibussowitsch PetscCall(PetscFree(v)); 94d24d4204SJose E. Roman if (mats_view) { 959566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "B: nrows %" PetscInt_FMT ", m %" PetscInt_FMT "; ncols %" PetscInt_FMT ", p %" PetscInt_FMT "\n",nrows,m,ncols,p)); 969566063dSJacob 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) */ 1009566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&b)); 1019566063dSJacob Faibussowitsch PetscCall(VecSetSizes(b,m,PETSC_DECIDE)); 1029566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(b)); 1039566063dSJacob Faibussowitsch PetscCall(VecGetArray(b,&barray)); 104d24d4204SJose E. Roman for (j=0;j<m;j++) { 1059566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(rand,&rval)); 106d24d4204SJose E. Roman barray[j] = rval; 107d24d4204SJose E. Roman } 1089566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(b,&barray)); 1099566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(b)); 1109566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(b)); 111d24d4204SJose E. Roman if (mats_view) { 1129566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] b: m %" PetscInt_FMT "\n",rank,m)); 1139566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 1149566063dSJacob Faibussowitsch PetscCall(VecView(b,PETSC_VIEWER_STDOUT_WORLD)); 115d24d4204SJose E. Roman } 1169566063dSJacob Faibussowitsch PetscCall(VecDuplicate(b,&x)); 117d24d4204SJose E. Roman 118d24d4204SJose E. Roman /* Create matrix X - same size as B */ 1199566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Create solution matrix X\n")); 1209566063dSJacob Faibussowitsch PetscCall(MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&X)); 121d24d4204SJose E. Roman 122d24d4204SJose E. Roman /* Cholesky factorization */ 123d24d4204SJose E. Roman /*------------------------*/ 1249566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Create ScaLAPACK matrix Aher\n")); 1259566063dSJacob Faibussowitsch PetscCall(MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&Aher)); 1269566063dSJacob Faibussowitsch PetscCall(MatAXPY(Aher,1.0,A,SAME_NONZERO_PATTERN)); /* Aher = A + A^T */ 1279566063dSJacob Faibussowitsch PetscCall(MatShift(Aher,100.0)); /* add 100.0 to diagonals of Aher to make it spd */ 128d24d4204SJose E. Roman if (mats_view) { 1299566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Aher:\n")); 1309566063dSJacob Faibussowitsch PetscCall(MatView(Aher,PETSC_VIEWER_STDOUT_WORLD)); 131d24d4204SJose E. Roman } 132d24d4204SJose E. Roman 133d24d4204SJose E. Roman /* Cholesky factorization */ 134d24d4204SJose E. Roman /*------------------------*/ 1359566063dSJacob 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 */ 1389566063dSJacob Faibussowitsch PetscCall(MatDuplicate(Aher,MAT_COPY_VALUES,&G)); 139d24d4204SJose E. Roman 140d24d4204SJose E. Roman /* G = L * L^T */ 1419566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactor(G,0,0)); 142d24d4204SJose E. Roman if (mats_view) { 1439566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Cholesky Factor G:\n")); 1449566063dSJacob 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 */ 1489566063dSJacob Faibussowitsch PetscCall(MatSolve(G,b,x)); 1499566063dSJacob Faibussowitsch PetscCall(MatMatSolve(G,B,X)); 1509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&G)); 151d24d4204SJose E. Roman 152d24d4204SJose E. Roman /* Out-place Cholesky */ 1539566063dSJacob Faibussowitsch PetscCall(MatGetFactor(Aher,MATSOLVERSCALAPACK,MAT_FACTOR_CHOLESKY,&G)); 1549566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorSymbolic(G,Aher,0,NULL)); 1559566063dSJacob Faibussowitsch PetscCall(MatCholeskyFactorNumeric(G,Aher,NULL)); 156*1baa6e33SBarry Smith if (mats_view) PetscCall(MatView(G,PETSC_VIEWER_STDOUT_WORLD)); 1579566063dSJacob Faibussowitsch PetscCall(MatSolve(G,b,x)); 1589566063dSJacob Faibussowitsch PetscCall(MatMatSolve(G,B,X)); 1599566063dSJacob Faibussowitsch PetscCall(MatDestroy(&G)); 160d24d4204SJose E. Roman 161d24d4204SJose E. Roman /* Check norm(Aher*x - b) */ 1629566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&c)); 1639566063dSJacob Faibussowitsch PetscCall(VecSetSizes(c,m,PETSC_DECIDE)); 1649566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(c)); 1659566063dSJacob Faibussowitsch PetscCall(MatMult(Aher,x,c)); 1669566063dSJacob Faibussowitsch PetscCall(VecAXPY(c,-1.0,b)); 1679566063dSJacob Faibussowitsch PetscCall(VecNorm(c,NORM_1,&norm)); 168d24d4204SJose E. Roman if (norm > tol) { 1699566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: ||Aher*x - b||=%g for Cholesky\n",(double)norm)); 170d24d4204SJose E. Roman } 171d24d4204SJose E. Roman 172d24d4204SJose E. Roman /* Check norm(Aher*X - B) */ 1739566063dSJacob Faibussowitsch PetscCall(MatMatMult(Aher,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); 1749566063dSJacob Faibussowitsch PetscCall(MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN)); 1759566063dSJacob Faibussowitsch PetscCall(MatNorm(C,NORM_1,&norm)); 176d24d4204SJose E. Roman if (norm > tol) { 1779566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: ||Aher*X - B||=%g for Cholesky\n",(double)norm)); 178d24d4204SJose E. Roman } 179d24d4204SJose E. Roman 180d24d4204SJose E. Roman /* LU factorization */ 181d24d4204SJose E. Roman /*------------------*/ 1829566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD," Test LU Solver \n")); 183d24d4204SJose E. Roman /* In-place LU */ 184d24d4204SJose E. Roman /* Create matrix factor F, with a copy of A */ 1859566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&F)); 186d24d4204SJose E. Roman /* Create vector d to test MatSolveAdd() */ 1879566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&d)); 1889566063dSJacob Faibussowitsch PetscCall(VecCopy(x,d)); 189d24d4204SJose E. Roman 190d24d4204SJose E. Roman /* PF=LU factorization */ 1919566063dSJacob Faibussowitsch PetscCall(MatLUFactor(F,0,0,NULL)); 192d24d4204SJose E. Roman 193d24d4204SJose E. Roman /* Solve LUX = PB */ 1949566063dSJacob Faibussowitsch PetscCall(MatSolveAdd(F,b,d,x)); 1959566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F,B,X)); 1969566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 197d24d4204SJose E. Roman 198d24d4204SJose E. Roman /* Check norm(A*X - B) */ 1999566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&e)); 2009566063dSJacob Faibussowitsch PetscCall(VecSetSizes(e,m,PETSC_DECIDE)); 2019566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(e)); 2029566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,c)); 2039566063dSJacob Faibussowitsch PetscCall(MatMult(A,d,e)); 2049566063dSJacob Faibussowitsch PetscCall(VecAXPY(c,-1.0,e)); 2059566063dSJacob Faibussowitsch PetscCall(VecAXPY(c,-1.0,b)); 2069566063dSJacob Faibussowitsch PetscCall(VecNorm(c,NORM_1,&norm)); 207d24d4204SJose E. Roman if (norm > tol) { 2089566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: ||A*x - b||=%g for LU\n",(double)norm)); 209d24d4204SJose E. Roman } 210d24d4204SJose E. Roman /* Reuse product C; replace Aher with A */ 2119566063dSJacob Faibussowitsch PetscCall(MatProductReplaceMats(A,NULL,NULL,C)); 2129566063dSJacob Faibussowitsch PetscCall(MatMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C)); 2139566063dSJacob Faibussowitsch PetscCall(MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN)); 2149566063dSJacob Faibussowitsch PetscCall(MatNorm(C,NORM_1,&norm)); 215d24d4204SJose E. Roman if (norm > tol) { 2169566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Warning: ||A*X - B||=%g for LU\n",(double)norm)); 217d24d4204SJose E. Roman } 218d24d4204SJose E. Roman 219d24d4204SJose E. Roman /* Out-place LU */ 2209566063dSJacob Faibussowitsch PetscCall(MatGetFactor(A,MATSOLVERSCALAPACK,MAT_FACTOR_LU,&F)); 2219566063dSJacob Faibussowitsch PetscCall(MatLUFactorSymbolic(F,A,0,0,NULL)); 2229566063dSJacob Faibussowitsch PetscCall(MatLUFactorNumeric(F,A,NULL)); 223*1baa6e33SBarry Smith if (mats_view) PetscCall(MatView(F,PETSC_VIEWER_STDOUT_WORLD)); 2249566063dSJacob Faibussowitsch PetscCall(MatSolve(F,b,x)); 2259566063dSJacob Faibussowitsch PetscCall(MatMatSolve(F,B,X)); 2269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&F)); 227d24d4204SJose E. Roman 228d24d4204SJose E. Roman /* Free space */ 2299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Aher)); 2319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 2329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 2339566063dSJacob Faibussowitsch PetscCall(MatDestroy(&X)); 2349566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b)); 2359566063dSJacob Faibussowitsch PetscCall(VecDestroy(&c)); 2369566063dSJacob Faibussowitsch PetscCall(VecDestroy(&d)); 2379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&e)); 2389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 2399566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rand)); 2409566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 241b122ec5aSJacob Faibussowitsch return 0; 242d24d4204SJose E. Roman } 243d24d4204SJose E. Roman 244d24d4204SJose E. Roman /*TEST 245d24d4204SJose E. Roman 246d24d4204SJose E. Roman build: 247d24d4204SJose E. Roman requires: scalapack 248d24d4204SJose E. Roman 249d24d4204SJose E. Roman test: 250d24d4204SJose E. Roman nsize: 2 251d24d4204SJose E. Roman output_file: output/ex245.out 252d24d4204SJose E. Roman 253d24d4204SJose E. Roman test: 254d24d4204SJose E. Roman suffix: 2 255d24d4204SJose E. Roman nsize: 6 256d24d4204SJose E. Roman output_file: output/ex245.out 257d24d4204SJose E. Roman 258d24d4204SJose E. Roman TEST*/ 259