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