static char help[] = "Tests LU, Cholesky factorization and MatMatSolve() for an Elemental dense matrix.\n\n"; #include int main(int argc,char **argv) { Mat A,F,B,X,C,Aher,G; Vec b,x,c,d,e; PetscErrorCode ierr; PetscInt m = 5,n,p,i,j,nrows,ncols; PetscScalar *v,*barray,rval; PetscReal norm,tol=1.e-11; PetscMPIInt size,rank; PetscRandom rand; const PetscInt *rows,*cols; IS isrows,iscols; PetscBool mats_view=PETSC_FALSE; MatFactorInfo finfo; ierr = PetscInitialize(&argc,&argv,(char*) 0,help);if (ierr) return ierr; CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); CHKERRQ(PetscRandomCreate(PETSC_COMM_WORLD,&rand)); CHKERRQ(PetscRandomSetFromOptions(rand)); /* Get local dimensions of matrices */ CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL)); n = m; CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL)); p = m/2; CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL)); CHKERRQ(PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view)); /* Create matrix A */ CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," Create Elemental matrix A\n")); CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); CHKERRQ(MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE)); CHKERRQ(MatSetType(A,MATELEMENTAL)); CHKERRQ(MatSetFromOptions(A)); CHKERRQ(MatSetUp(A)); /* Set local matrix entries */ CHKERRQ(MatGetOwnershipIS(A,&isrows,&iscols)); CHKERRQ(ISGetLocalSize(isrows,&nrows)); CHKERRQ(ISGetIndices(isrows,&rows)); CHKERRQ(ISGetLocalSize(iscols,&ncols)); CHKERRQ(ISGetIndices(iscols,&cols)); CHKERRQ(PetscMalloc1(nrows*ncols,&v)); for (i=0; i tol) { CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Warning: |Aher*x - b| for Cholesky %g\n",(double)norm)); } /* Check norm(Aher*X - B) */ CHKERRQ(MatMatMult(Aher,X,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&C)); CHKERRQ(MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN)); CHKERRQ(MatNorm(C,NORM_1,&norm)); if (norm > tol) { CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Warning: |Aher*X - B| for Cholesky %g\n",(double)norm)); } /* LU factorization */ /*------------------*/ CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," Test LU Solver \n")); /* In-place LU */ /* Create matrix factor F, then copy A to F */ CHKERRQ(MatCreate(PETSC_COMM_WORLD,&F)); CHKERRQ(MatSetSizes(F,m,n,PETSC_DECIDE,PETSC_DECIDE)); CHKERRQ(MatSetType(F,MATELEMENTAL)); CHKERRQ(MatSetFromOptions(F)); CHKERRQ(MatSetUp(F)); CHKERRQ(MatAssemblyBegin(F,MAT_FINAL_ASSEMBLY)); CHKERRQ(MatAssemblyEnd(F,MAT_FINAL_ASSEMBLY)); CHKERRQ(MatCopy(A,F,SAME_NONZERO_PATTERN)); /* Create vector d to test MatSolveAdd() */ CHKERRQ(VecDuplicate(x,&d)); CHKERRQ(VecCopy(x,d)); /* PF=LU or F=LU factorization - perms is ignored by Elemental; set finfo.dtcol !0 or 0 to enable/disable partial pivoting */ finfo.dtcol = 0.1; CHKERRQ(MatLUFactor(F,0,0,&finfo)); /* Solve LUX = PB or LUX = B */ CHKERRQ(MatSolveAdd(F,b,d,x)); CHKERRQ(MatMatSolve(F,B,X)); CHKERRQ(MatDestroy(&F)); /* Check norm(A*X - B) */ CHKERRQ(VecCreate(PETSC_COMM_WORLD,&e)); CHKERRQ(VecSetSizes(e,m,PETSC_DECIDE)); CHKERRQ(VecSetFromOptions(e)); CHKERRQ(MatMult(A,x,c)); CHKERRQ(MatMult(A,d,e)); CHKERRQ(VecAXPY(c,-1.0,e)); CHKERRQ(VecAXPY(c,-1.0,b)); CHKERRQ(VecNorm(c,NORM_1,&norm)); if (norm > tol) { CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Warning: |A*x - b| for LU %g\n",(double)norm)); } /* Reuse product C; replace Aher with A */ CHKERRQ(MatProductReplaceMats(A,NULL,NULL,C)); CHKERRQ(MatMatMult(A,X,MAT_REUSE_MATRIX,PETSC_DEFAULT,&C)); CHKERRQ(MatAXPY(C,-1.0,B,SAME_NONZERO_PATTERN)); CHKERRQ(MatNorm(C,NORM_1,&norm)); if (norm > tol) { CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Warning: |A*X - B| for LU %g\n",(double)norm)); } /* Out-place LU */ CHKERRQ(MatGetFactor(A,MATSOLVERELEMENTAL,MAT_FACTOR_LU,&F)); CHKERRQ(MatLUFactorSymbolic(F,A,0,0,&finfo)); CHKERRQ(MatLUFactorNumeric(F,A,&finfo)); if (mats_view) { CHKERRQ(MatView(F,PETSC_VIEWER_STDOUT_WORLD)); } CHKERRQ(MatSolve(F,b,x)); CHKERRQ(MatMatSolve(F,B,X)); CHKERRQ(MatDestroy(&F)); /* Free space */ CHKERRQ(MatDestroy(&A)); CHKERRQ(MatDestroy(&Aher)); CHKERRQ(MatDestroy(&B)); CHKERRQ(MatDestroy(&C)); CHKERRQ(MatDestroy(&X)); CHKERRQ(VecDestroy(&b)); CHKERRQ(VecDestroy(&c)); CHKERRQ(VecDestroy(&d)); CHKERRQ(VecDestroy(&e)); CHKERRQ(VecDestroy(&x)); CHKERRQ(PetscRandomDestroy(&rand)); ierr = PetscFinalize(); return ierr; } /*TEST build: requires: elemental test: nsize: 2 output_file: output/ex145.out test: suffix: 2 nsize: 6 output_file: output/ex145.out TEST*/