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