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*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*) 0,help)); 205f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 215f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 22d24d4204SJose E. Roman 235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(rand)); 25d24d4204SJose E. Roman 26d24d4204SJose E. Roman /* Get local dimensions of matrices */ 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); 28d24d4204SJose E. Roman n = m; 295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); 30d24d4204SJose E. Roman p = m/2; 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL)); 325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view)); 33d24d4204SJose E. Roman 34d24d4204SJose E. Roman /* Create matrix A */ 355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," Create ScaLAPACK matrix A\n")); 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATSCALAPACK)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 41d24d4204SJose E. Roman /* Set local matrix entries */ 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipIS(A,&isrows,&iscols)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(isrows,&nrows)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isrows,&rows)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(iscols,&ncols)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(iscols,&cols)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nrows*ncols,&v)); 48d24d4204SJose E. Roman for (i=0;i<nrows;i++) { 49d24d4204SJose E. Roman for (j=0;j<ncols;j++) { 505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rand,&rval)); 51d24d4204SJose E. Roman v[i*ncols+j] = rval; 52d24d4204SJose E. Roman } 53d24d4204SJose E. Roman } 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,nrows,rows,ncols,cols,v,INSERT_VALUES)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isrows,&rows)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(iscols,&cols)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrows)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscols)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 62d24d4204SJose E. Roman if (mats_view) { 635f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "A: nrows %" PetscInt_FMT ", m %" PetscInt_FMT "; ncols %" PetscInt_FMT ", n %" PetscInt_FMT "\n",nrows,m,ncols,n)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 65d24d4204SJose E. Roman } 66d24d4204SJose E. Roman 67d24d4204SJose E. Roman /* Create rhs matrix B */ 685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," Create rhs matrix B\n")); 695f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(B,m,p,PETSC_DECIDE,PETSC_DECIDE)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B,MATSCALAPACK)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(B)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(B)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipIS(B,&isrows,&iscols)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(isrows,&nrows)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isrows,&rows)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(iscols,&ncols)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(iscols,&cols)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nrows*ncols,&v)); 80d24d4204SJose E. Roman for (i=0;i<nrows;i++) { 81d24d4204SJose E. Roman for (j=0;j<ncols;j++) { 825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rand,&rval)); 83d24d4204SJose E. Roman v[i*ncols+j] = rval; 84d24d4204SJose E. Roman } 85d24d4204SJose E. Roman } 865f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isrows,&rows)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(iscols,&cols)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrows)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscols)); 935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(v)); 94d24d4204SJose E. Roman if (mats_view) { 955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "B: nrows %" PetscInt_FMT ", m %" PetscInt_FMT "; ncols %" PetscInt_FMT ", p %" PetscInt_FMT "\n",nrows,m,ncols,p)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(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) */ 1005f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&b)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(b,m,PETSC_DECIDE)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(b)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(b,&barray)); 104d24d4204SJose E. Roman for (j=0;j<m;j++) { 1055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(rand,&rval)); 106d24d4204SJose E. Roman barray[j] = rval; 107d24d4204SJose E. Roman } 1085f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(b,&barray)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyBegin(b)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(VecAssemblyEnd(b)); 111d24d4204SJose E. Roman if (mats_view) { 1125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] b: m %" PetscInt_FMT "\n",rank,m)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSynchronizedFlush(PETSC_COMM_WORLD,PETSC_STDOUT)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(b,PETSC_VIEWER_STDOUT_WORLD)); 115d24d4204SJose E. Roman } 1165f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(b,&x)); 117d24d4204SJose E. Roman 118d24d4204SJose E. Roman /* Create matrix X - same size as B */ 1195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," Create solution matrix X\n")); 1205f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(B,MAT_DO_NOT_COPY_VALUES,&X)); 121d24d4204SJose E. Roman 122d24d4204SJose E. Roman /* Cholesky factorization */ 123d24d4204SJose E. Roman /*------------------------*/ 1245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," Create ScaLAPACK matrix Aher\n")); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(MatHermitianTranspose(A,MAT_INITIAL_MATRIX,&Aher)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(Aher,1.0,A,SAME_NONZERO_PATTERN)); /* Aher = A + A^T */ 1275f80ce2aSJacob Faibussowitsch CHKERRQ(MatShift(Aher,100.0)); /* add 100.0 to diagonals of Aher to make it spd */ 128d24d4204SJose E. Roman if (mats_view) { 1295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Aher:\n")); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(Aher,PETSC_VIEWER_STDOUT_WORLD)); 131d24d4204SJose E. Roman } 132d24d4204SJose E. Roman 133d24d4204SJose E. Roman /* Cholesky factorization */ 134d24d4204SJose E. Roman /*------------------------*/ 1355f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 1385f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(Aher,MAT_COPY_VALUES,&G)); 139d24d4204SJose E. Roman 140d24d4204SJose E. Roman /* G = L * L^T */ 1415f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactor(G,0,0)); 142d24d4204SJose E. Roman if (mats_view) { 1435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD, "Cholesky Factor G:\n")); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 1485f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(G,b,x)); 1495f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatSolve(G,B,X)); 1505f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&G)); 151d24d4204SJose E. Roman 152d24d4204SJose E. Roman /* Out-place Cholesky */ 1535f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(Aher,MATSOLVERSCALAPACK,MAT_FACTOR_CHOLESKY,&G)); 1545f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorSymbolic(G,Aher,0,NULL)); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorNumeric(G,Aher,NULL)); 156d24d4204SJose E. Roman if (mats_view) { 1575f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(G,PETSC_VIEWER_STDOUT_WORLD)); 158d24d4204SJose E. Roman } 1595f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(G,b,x)); 1605f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatSolve(G,B,X)); 1615f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&G)); 162d24d4204SJose E. Roman 163d24d4204SJose E. Roman /* Check norm(Aher*x - b) */ 1645f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&c)); 1655f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(c,m,PETSC_DECIDE)); 1665f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(c)); 1675f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(Aher,x,c)); 1685f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(c,-1.0,b)); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(c,NORM_1,&norm)); 170d24d4204SJose E. Roman if (norm > tol) { 1715f80ce2aSJacob Faibussowitsch CHKERRQ(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) */ 1755f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(Aher,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); 1765f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN)); 1775f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C,NORM_1,&norm)); 178d24d4204SJose E. Roman if (norm > tol) { 1795f80ce2aSJacob Faibussowitsch CHKERRQ(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 /*------------------*/ 1845f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 1875f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&F)); 188d24d4204SJose E. Roman /* Create vector d to test MatSolveAdd() */ 1895f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&d)); 1905f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(x,d)); 191d24d4204SJose E. Roman 192d24d4204SJose E. Roman /* PF=LU factorization */ 1935f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactor(F,0,0,NULL)); 194d24d4204SJose E. Roman 195d24d4204SJose E. Roman /* Solve LUX = PB */ 1965f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolveAdd(F,b,d,x)); 1975f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatSolve(F,B,X)); 1985f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&F)); 199d24d4204SJose E. Roman 200d24d4204SJose E. Roman /* Check norm(A*X - B) */ 2015f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&e)); 2025f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(e,m,PETSC_DECIDE)); 2035f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(e)); 2045f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,x,c)); 2055f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,d,e)); 2065f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(c,-1.0,e)); 2075f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(c,-1.0,b)); 2085f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(c,NORM_1,&norm)); 209d24d4204SJose E. Roman if (norm > tol) { 2105f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 2135f80ce2aSJacob Faibussowitsch CHKERRQ(MatProductReplaceMats(A,NULL,NULL,C)); 2145f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C)); 2155f80ce2aSJacob Faibussowitsch CHKERRQ(MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN)); 2165f80ce2aSJacob Faibussowitsch CHKERRQ(MatNorm(C,NORM_1,&norm)); 217d24d4204SJose E. Roman if (norm > tol) { 2185f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 2225f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(A,MATSOLVERSCALAPACK,MAT_FACTOR_LU,&F)); 2235f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorSymbolic(F,A,0,0,NULL)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorNumeric(F,A,NULL)); 225d24d4204SJose E. Roman if (mats_view) { 2265f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(F,PETSC_VIEWER_STDOUT_WORLD)); 227d24d4204SJose E. Roman } 2285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(F,b,x)); 2295f80ce2aSJacob Faibussowitsch CHKERRQ(MatMatSolve(F,B,X)); 2305f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&F)); 231d24d4204SJose E. Roman 232d24d4204SJose E. Roman /* Free space */ 2335f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 2345f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Aher)); 2355f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 2365f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 2375f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&X)); 2385f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&c)); 2405f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&d)); 2415f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&e)); 2425f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 2435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&rand)); 244*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 245*b122ec5aSJacob 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