1c4762a1bSJed Brown static char help[] = "Illustrate how to do one symbolic factorization and multiple numeric factorizations using same matrix structure. \n\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include <petscmat.h> 4c4762a1bSJed Brown 5c4762a1bSJed Brown int main(int argc,char **args) 6c4762a1bSJed Brown { 7*dc6ed827SStefano Zampini PetscInt i,rstart,rend,N=10,num_numfac=5,col[3],k; 8c4762a1bSJed Brown Mat A[5],F; 9c4762a1bSJed Brown Vec u,x,b; 10c4762a1bSJed Brown PetscErrorCode ierr; 11c4762a1bSJed Brown PetscMPIInt rank; 12c4762a1bSJed Brown PetscScalar value[3]; 13c4762a1bSJed Brown PetscReal norm,tol=100*PETSC_MACHINE_EPSILON; 14c4762a1bSJed Brown IS perm,iperm; 15c4762a1bSJed Brown MatFactorInfo info; 16*dc6ed827SStefano Zampini MatFactorType facttype = MAT_FACTOR_LU; 17*dc6ed827SStefano Zampini char solvertype[64]; 18*dc6ed827SStefano Zampini char factortype[64]; 19c4762a1bSJed Brown 20c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 21c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRQ(ierr); 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* Create and assemble matrices, all have same data structure */ 24c4762a1bSJed Brown for (k=0; k<num_numfac; k++) { 25c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A[k]);CHKERRQ(ierr); 26c4762a1bSJed Brown ierr = MatSetSizes(A[k],PETSC_DECIDE,PETSC_DECIDE,N,N);CHKERRQ(ierr); 27c4762a1bSJed Brown ierr = MatSetFromOptions(A[k]);CHKERRQ(ierr); 28c4762a1bSJed Brown ierr = MatSetUp(A[k]);CHKERRQ(ierr); 29c4762a1bSJed Brown ierr = MatGetOwnershipRange(A[k],&rstart,&rend);CHKERRQ(ierr); 30c4762a1bSJed Brown 31*dc6ed827SStefano Zampini value[0] = -1.0*(k+1); 32*dc6ed827SStefano Zampini value[1] = 2.0*(k+1); 33*dc6ed827SStefano Zampini value[2] = -1.0*(k+1); 34c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 35c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 36c4762a1bSJed Brown if (i == 0) { 37c4762a1bSJed Brown ierr = MatSetValues(A[k],1,&i,2,col+1,value+1,INSERT_VALUES);CHKERRQ(ierr); 38c4762a1bSJed Brown } else if (i == N-1) { 39c4762a1bSJed Brown ierr = MatSetValues(A[k],1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr); 40c4762a1bSJed Brown } else { 41c4762a1bSJed Brown ierr = MatSetValues(A[k],1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr); 42c4762a1bSJed Brown } 43c4762a1bSJed Brown } 44c4762a1bSJed Brown ierr = MatAssemblyBegin(A[k],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 45c4762a1bSJed Brown ierr = MatAssemblyEnd(A[k],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 46c4762a1bSJed Brown ierr = MatSetOption(A[k],MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 47c4762a1bSJed Brown } 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* Create vectors */ 50c4762a1bSJed Brown ierr = MatCreateVecs(A[0],&x,&b);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = VecDuplicate(x,&u);CHKERRQ(ierr); 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* Set rhs vector b */ 54c4762a1bSJed Brown ierr = VecSet(b,1.0);CHKERRQ(ierr); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Get a symbolic factor F from A[0] */ 57*dc6ed827SStefano Zampini ierr = PetscStrncpy(solvertype,"petsc",sizeof(solvertype));CHKERRQ(ierr); 58*dc6ed827SStefano Zampini ierr = PetscOptionsGetString(NULL, NULL, "-mat_solver_type",solvertype,sizeof(solvertype),NULL);CHKERRQ(ierr); 59*dc6ed827SStefano Zampini ierr = PetscOptionsGetEnum(NULL,NULL,"-mat_factor_type",MatFactorTypes,(PetscEnum*)&facttype,NULL);CHKERRQ(ierr); 60c4762a1bSJed Brown 61*dc6ed827SStefano Zampini ierr = MatGetFactor(A[0],solvertype,facttype,&F);CHKERRQ(ierr); 62c4762a1bSJed Brown /* test mumps options */ 63*dc6ed827SStefano Zampini #if defined(PETSC_HAVE_MUMPS) 64*dc6ed827SStefano Zampini ierr = MatMumpsSetIcntl(F,7,5);CHKERRQ(ierr); 65c4762a1bSJed Brown #endif 66*dc6ed827SStefano Zampini ierr = PetscStrncpy(factortype,MatFactorTypes[facttype],sizeof(factortype));CHKERRQ(ierr); 67*dc6ed827SStefano Zampini ierr = PetscStrtoupper(solvertype);CHKERRQ(ierr); 68*dc6ed827SStefano Zampini ierr = PetscStrtoupper(factortype);CHKERRQ(ierr); 69*dc6ed827SStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD," %s %s:\n",solvertype,factortype);CHKERRQ(ierr); 70c4762a1bSJed Brown 71c4762a1bSJed Brown ierr = MatFactorInfoInitialize(&info);CHKERRQ(ierr); 72c4762a1bSJed Brown info.fill = 5.0; 73c4762a1bSJed Brown ierr = MatGetOrdering(A[0],MATORDERINGNATURAL,&perm,&iperm);CHKERRQ(ierr); 74*dc6ed827SStefano Zampini switch (facttype) { 75*dc6ed827SStefano Zampini case MAT_FACTOR_LU: 76c4762a1bSJed Brown ierr = MatLUFactorSymbolic(F,A[0],perm,iperm,&info);CHKERRQ(ierr); 77*dc6ed827SStefano Zampini break; 78*dc6ed827SStefano Zampini case MAT_FACTOR_ILU: 79*dc6ed827SStefano Zampini ierr = MatILUFactorSymbolic(F,A[0],perm,iperm,&info);CHKERRQ(ierr); 80*dc6ed827SStefano Zampini break; 81*dc6ed827SStefano Zampini case MAT_FACTOR_ICC: 82*dc6ed827SStefano Zampini ierr = MatICCFactorSymbolic(F,A[0],perm,&info);CHKERRQ(ierr); 83*dc6ed827SStefano Zampini break; 84*dc6ed827SStefano Zampini case MAT_FACTOR_CHOLESKY: 85*dc6ed827SStefano Zampini ierr = MatCholeskyFactorSymbolic(F,A[0],perm,&info);CHKERRQ(ierr); 86*dc6ed827SStefano Zampini break; 87*dc6ed827SStefano Zampini default: 88*dc6ed827SStefano Zampini SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Not for factor type %s\n",factortype); 89*dc6ed827SStefano Zampini break; 90*dc6ed827SStefano Zampini } 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* Compute numeric factors using same F, then solve */ 93c4762a1bSJed Brown for (k = 0; k < num_numfac; k++) { 94*dc6ed827SStefano Zampini switch (facttype) { 95*dc6ed827SStefano Zampini case MAT_FACTOR_LU: 96*dc6ed827SStefano Zampini case MAT_FACTOR_ILU: 97c4762a1bSJed Brown ierr = MatLUFactorNumeric(F,A[k],&info);CHKERRQ(ierr); 98*dc6ed827SStefano Zampini break; 99*dc6ed827SStefano Zampini case MAT_FACTOR_ICC: 100*dc6ed827SStefano Zampini case MAT_FACTOR_CHOLESKY: 101*dc6ed827SStefano Zampini ierr = MatCholeskyFactorNumeric(F,A[k],&info);CHKERRQ(ierr); 102*dc6ed827SStefano Zampini break; 103*dc6ed827SStefano Zampini default: 104*dc6ed827SStefano Zampini SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Not for factor type %s\n",factortype); 105*dc6ed827SStefano Zampini break; 106*dc6ed827SStefano Zampini } 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* Solve A[k] * x = b */ 109c4762a1bSJed Brown ierr = MatSolve(F,b,x);CHKERRQ(ierr); 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* Check the residual */ 112c4762a1bSJed Brown ierr = MatMult(A[k],x,u);CHKERRQ(ierr); 113c4762a1bSJed Brown ierr = VecAXPY(u,-1.0,b);CHKERRQ(ierr); 114c4762a1bSJed Brown ierr = VecNorm(u,NORM_INFINITY,&norm);CHKERRQ(ierr); 115c4762a1bSJed Brown if (norm > tol) { 116*dc6ed827SStefano Zampini ierr = PetscPrintf(PETSC_COMM_WORLD,"%D-the %s numfact and solve: residual %g\n",k,factortype,(double)norm);CHKERRQ(ierr); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown } 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* Free data structures */ 121c4762a1bSJed Brown for (k=0; k<num_numfac; k++) { 122c4762a1bSJed Brown ierr = MatDestroy(&A[k]);CHKERRQ(ierr); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown ierr = MatDestroy(&F);CHKERRQ(ierr); 125c4762a1bSJed Brown ierr = ISDestroy(&perm);CHKERRQ(ierr); 126c4762a1bSJed Brown ierr = ISDestroy(&iperm);CHKERRQ(ierr); 127c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 128c4762a1bSJed Brown ierr = VecDestroy(&b);CHKERRQ(ierr); 129c4762a1bSJed Brown ierr = VecDestroy(&u);CHKERRQ(ierr); 130c4762a1bSJed Brown ierr = PetscFinalize(); 131c4762a1bSJed Brown return ierr; 132c4762a1bSJed Brown } 133c4762a1bSJed Brown 134c4762a1bSJed Brown /*TEST 135c4762a1bSJed Brown 136c4762a1bSJed Brown test: 137c4762a1bSJed Brown 138c4762a1bSJed Brown test: 139c4762a1bSJed Brown suffix: 2 140c4762a1bSJed Brown args: -mat_solver_type superlu 141c4762a1bSJed Brown requires: superlu 142c4762a1bSJed Brown 143c4762a1bSJed Brown test: 144c4762a1bSJed Brown suffix: 3 145c4762a1bSJed Brown nsize: 2 146c4762a1bSJed Brown requires: mumps 147c4762a1bSJed Brown args: -mat_solver_type mumps 148c4762a1bSJed Brown 149*dc6ed827SStefano Zampini test: 150*dc6ed827SStefano Zampini suffix: 4 151*dc6ed827SStefano Zampini args: -mat_solver_type cusparse -mat_type aijcusparse -mat_factor_type {{lu cholesky ilu icc}separate output} 152*dc6ed827SStefano Zampini requires: cuda 153*dc6ed827SStefano Zampini 154c4762a1bSJed Brown TEST*/ 155