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