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 { 7dc6ed827SStefano 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 PetscMPIInt rank; 11c4762a1bSJed Brown PetscScalar value[3]; 12c4762a1bSJed Brown PetscReal norm,tol=100*PETSC_MACHINE_EPSILON; 13c4762a1bSJed Brown IS perm,iperm; 14c4762a1bSJed Brown MatFactorInfo info; 15dc6ed827SStefano Zampini MatFactorType facttype = MAT_FACTOR_LU; 16dc6ed827SStefano Zampini char solvertype[64]; 17dc6ed827SStefano Zampini char factortype[64]; 18c4762a1bSJed Brown 19*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 205f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* Create and assemble matrices, all have same data structure */ 23c4762a1bSJed Brown for (k=0; k<num_numfac; k++) { 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A[k])); 255f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A[k],PETSC_DECIDE,PETSC_DECIDE,N,N)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A[k])); 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A[k])); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A[k],&rstart,&rend)); 29c4762a1bSJed Brown 30dc6ed827SStefano Zampini value[0] = -1.0*(k+1); 31dc6ed827SStefano Zampini value[1] = 2.0*(k+1); 32dc6ed827SStefano Zampini value[2] = -1.0*(k+1); 33c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 34c4762a1bSJed Brown col[0] = i-1; col[1] = i; col[2] = i+1; 35c4762a1bSJed Brown if (i == 0) { 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A[k],1,&i,2,col+1,value+1,INSERT_VALUES)); 37c4762a1bSJed Brown } else if (i == N-1) { 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A[k],1,&i,2,col,value,INSERT_VALUES)); 39c4762a1bSJed Brown } else { 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A[k],1,&i,3,col,value,INSERT_VALUES)); 41c4762a1bSJed Brown } 42c4762a1bSJed Brown } 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A[k],MAT_FINAL_ASSEMBLY)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A[k],MAT_FINAL_ASSEMBLY)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(A[k],MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE)); 46c4762a1bSJed Brown } 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* Create vectors */ 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A[0],&x,&b)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&u)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Set rhs vector b */ 535f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(b,1.0)); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Get a symbolic factor F from A[0] */ 565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncpy(solvertype,"petsc",sizeof(solvertype))); 575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL, NULL, "-mat_solver_type",solvertype,sizeof(solvertype),NULL)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetEnum(NULL,NULL,"-mat_factor_type",MatFactorTypes,(PetscEnum*)&facttype,NULL)); 59c4762a1bSJed Brown 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetFactor(A[0],solvertype,facttype,&F)); 61c4762a1bSJed Brown /* test mumps options */ 62dc6ed827SStefano Zampini #if defined(PETSC_HAVE_MUMPS) 635f80ce2aSJacob Faibussowitsch CHKERRQ(MatMumpsSetIcntl(F,7,5)); 64c4762a1bSJed Brown #endif 655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrncpy(factortype,MatFactorTypes[facttype],sizeof(factortype))); 665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrtoupper(solvertype)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscStrtoupper(factortype)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD," %s %s:\n",solvertype,factortype)); 69c4762a1bSJed Brown 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatFactorInfoInitialize(&info)); 71c4762a1bSJed Brown info.fill = 5.0; 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOrdering(A[0],MATORDERINGNATURAL,&perm,&iperm)); 73dc6ed827SStefano Zampini switch (facttype) { 74dc6ed827SStefano Zampini case MAT_FACTOR_LU: 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorSymbolic(F,A[0],perm,iperm,&info)); 76dc6ed827SStefano Zampini break; 77dc6ed827SStefano Zampini case MAT_FACTOR_ILU: 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatILUFactorSymbolic(F,A[0],perm,iperm,&info)); 79dc6ed827SStefano Zampini break; 80dc6ed827SStefano Zampini case MAT_FACTOR_ICC: 815f80ce2aSJacob Faibussowitsch CHKERRQ(MatICCFactorSymbolic(F,A[0],perm,&info)); 82dc6ed827SStefano Zampini break; 83dc6ed827SStefano Zampini case MAT_FACTOR_CHOLESKY: 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorSymbolic(F,A[0],perm,&info)); 85dc6ed827SStefano Zampini break; 86dc6ed827SStefano Zampini default: 8798921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Not for factor type %s",factortype); 88dc6ed827SStefano Zampini } 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* Compute numeric factors using same F, then solve */ 91c4762a1bSJed Brown for (k = 0; k < num_numfac; k++) { 92dc6ed827SStefano Zampini switch (facttype) { 93dc6ed827SStefano Zampini case MAT_FACTOR_LU: 94dc6ed827SStefano Zampini case MAT_FACTOR_ILU: 955f80ce2aSJacob Faibussowitsch CHKERRQ(MatLUFactorNumeric(F,A[k],&info)); 96dc6ed827SStefano Zampini break; 97dc6ed827SStefano Zampini case MAT_FACTOR_ICC: 98dc6ed827SStefano Zampini case MAT_FACTOR_CHOLESKY: 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatCholeskyFactorNumeric(F,A[k],&info)); 100dc6ed827SStefano Zampini break; 101dc6ed827SStefano Zampini default: 10298921bdaSJacob Faibussowitsch SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Not for factor type %s",factortype); 103dc6ed827SStefano Zampini } 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* Solve A[k] * x = b */ 1065f80ce2aSJacob Faibussowitsch CHKERRQ(MatSolve(F,b,x)); 107c4762a1bSJed Brown 108c4762a1bSJed Brown /* Check the residual */ 1095f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A[k],x,u)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(u,-1.0,b)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(u,NORM_INFINITY,&norm)); 112c4762a1bSJed Brown if (norm > tol) { 1135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%" PetscInt_FMT "-the %s numfact and solve: residual %g\n",k,factortype,(double)norm)); 114c4762a1bSJed Brown } 115c4762a1bSJed Brown } 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* Free data structures */ 118c4762a1bSJed Brown for (k=0; k<num_numfac; k++) { 1195f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A[k])); 120c4762a1bSJed Brown } 1215f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&F)); 1225f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&perm)); 1235f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iperm)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&b)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&u)); 127*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 128*b122ec5aSJacob Faibussowitsch return 0; 129c4762a1bSJed Brown } 130c4762a1bSJed Brown 131c4762a1bSJed Brown /*TEST 132c4762a1bSJed Brown 133c4762a1bSJed Brown test: 134c4762a1bSJed Brown 135c4762a1bSJed Brown test: 136c4762a1bSJed Brown suffix: 2 137c4762a1bSJed Brown args: -mat_solver_type superlu 138c4762a1bSJed Brown requires: superlu 139c4762a1bSJed Brown 140c4762a1bSJed Brown test: 141c4762a1bSJed Brown suffix: 3 142c4762a1bSJed Brown nsize: 2 143c4762a1bSJed Brown requires: mumps 144c4762a1bSJed Brown args: -mat_solver_type mumps 145c4762a1bSJed Brown 146dc6ed827SStefano Zampini test: 147dc6ed827SStefano Zampini suffix: 4 148dc6ed827SStefano Zampini args: -mat_solver_type cusparse -mat_type aijcusparse -mat_factor_type {{lu cholesky ilu icc}separate output} 149dc6ed827SStefano Zampini requires: cuda 150dc6ed827SStefano Zampini 151c4762a1bSJed Brown TEST*/ 152