xref: /petsc/src/mat/tests/ex106.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Test repeated LU factorizations. Used for checking memory leak\n\
3*c4762a1bSJed Brown   -m <size> : problem size\n\
4*c4762a1bSJed Brown   -mat_nonsym : use nonsymmetric matrix (default is symmetric)\n\n";
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown #include <petscmat.h>
7*c4762a1bSJed Brown int main(int argc,char **args)
8*c4762a1bSJed Brown {
9*c4762a1bSJed Brown   Mat            C,F;                /* matrix */
10*c4762a1bSJed Brown   Vec            x,u,b;          /* approx solution, RHS, exact solution */
11*c4762a1bSJed Brown   PetscReal      norm;             /* norm of solution error */
12*c4762a1bSJed Brown   PetscScalar    v,none = -1.0;
13*c4762a1bSJed Brown   PetscInt       I,J,ldim,low,high,iglobal,Istart,Iend;
14*c4762a1bSJed Brown   PetscErrorCode ierr;
15*c4762a1bSJed Brown   PetscInt       i,j,m = 3,n = 2,its;
16*c4762a1bSJed Brown   PetscMPIInt    size,rank;
17*c4762a1bSJed Brown   PetscBool      mat_nonsymmetric;
18*c4762a1bSJed Brown   PetscInt       its_max;
19*c4762a1bSJed Brown   MatFactorInfo  factinfo;
20*c4762a1bSJed Brown   IS             perm,iperm;
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
23*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
24*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
25*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
26*c4762a1bSJed Brown   n    = 2*size;
27*c4762a1bSJed Brown 
28*c4762a1bSJed Brown   /*
29*c4762a1bSJed Brown      Set flag if we are doing a nonsymmetric problem; the default is symmetric.
30*c4762a1bSJed Brown   */
31*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-mat_nonsym",&mat_nonsymmetric);CHKERRQ(ierr);
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown   /*
34*c4762a1bSJed Brown      Create parallel matrix, specifying only its global dimensions.
35*c4762a1bSJed Brown      When using MatCreate(), the matrix format can be specified at
36*c4762a1bSJed Brown      runtime. Also, the parallel partitioning of the matrix is
37*c4762a1bSJed Brown      determined by PETSc at runtime.
38*c4762a1bSJed Brown   */
39*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(C,&Istart,&Iend);CHKERRQ(ierr);
43*c4762a1bSJed Brown 
44*c4762a1bSJed Brown   /*
45*c4762a1bSJed Brown      Set matrix entries matrix in parallel.
46*c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
47*c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
48*c4762a1bSJed Brown         appropriate processor during matrix assembly).
49*c4762a1bSJed Brown       - Always specify global row and columns of matrix entries.
50*c4762a1bSJed Brown   */
51*c4762a1bSJed Brown   for (I=Istart; I<Iend; I++) {
52*c4762a1bSJed Brown     v = -1.0; i = I/n; j = I - i*n;
53*c4762a1bSJed Brown     if (i>0)   {J = I - n; ierr = MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);}
54*c4762a1bSJed Brown     if (i<m-1) {J = I + n; ierr = MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);}
55*c4762a1bSJed Brown     if (j>0)   {J = I - 1; ierr = MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);}
56*c4762a1bSJed Brown     if (j<n-1) {J = I + 1; ierr = MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);}
57*c4762a1bSJed Brown     v = 4.0; ierr = MatSetValues(C,1,&I,1,&I,&v,ADD_VALUES);
58*c4762a1bSJed Brown   }
59*c4762a1bSJed Brown 
60*c4762a1bSJed Brown   /*
61*c4762a1bSJed Brown      Make the matrix nonsymmetric if desired
62*c4762a1bSJed Brown   */
63*c4762a1bSJed Brown   if (mat_nonsymmetric) {
64*c4762a1bSJed Brown     for (I=Istart; I<Iend; I++) {
65*c4762a1bSJed Brown       v = -1.5; i = I/n;
66*c4762a1bSJed Brown       if (i>1)   {J = I-n-1; ierr = MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES);CHKERRQ(ierr);}
67*c4762a1bSJed Brown     }
68*c4762a1bSJed Brown   } else {
69*c4762a1bSJed Brown     ierr = MatSetOption(C,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr);
70*c4762a1bSJed Brown     ierr = MatSetOption(C,MAT_SYMMETRY_ETERNAL,PETSC_TRUE);CHKERRQ(ierr);
71*c4762a1bSJed Brown   }
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown   /*
74*c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
75*c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
76*c4762a1bSJed Brown      Computations can be done while messages are in transition
77*c4762a1bSJed Brown      by placing code between these two statements.
78*c4762a1bSJed Brown   */
79*c4762a1bSJed Brown   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
80*c4762a1bSJed Brown   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
81*c4762a1bSJed Brown 
82*c4762a1bSJed Brown   its_max=1000;
83*c4762a1bSJed Brown   /*
84*c4762a1bSJed Brown      Create parallel vectors.
85*c4762a1bSJed Brown       - When using VecSetSizes(), we specify only the vector's global
86*c4762a1bSJed Brown         dimension; the parallel partitioning is determined at runtime.
87*c4762a1bSJed Brown       - Note: We form 1 vector from scratch and then duplicate as needed.
88*c4762a1bSJed Brown   */
89*c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&u);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = VecSetSizes(u,PETSC_DECIDE,m*n);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = VecSetFromOptions(u);CHKERRQ(ierr);
92*c4762a1bSJed Brown   ierr = VecDuplicate(u,&b);CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr = VecDuplicate(b,&x);CHKERRQ(ierr);
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown   /*
96*c4762a1bSJed Brown      Currently, all parallel PETSc vectors are partitioned by
97*c4762a1bSJed Brown      contiguous chunks across the processors.  Determine which
98*c4762a1bSJed Brown      range of entries are locally owned.
99*c4762a1bSJed Brown   */
100*c4762a1bSJed Brown   ierr = VecGetOwnershipRange(x,&low,&high);CHKERRQ(ierr);
101*c4762a1bSJed Brown 
102*c4762a1bSJed Brown   /*
103*c4762a1bSJed Brown     Set elements within the exact solution vector in parallel.
104*c4762a1bSJed Brown      - Each processor needs to insert only elements that it owns
105*c4762a1bSJed Brown        locally (but any non-local entries will be sent to the
106*c4762a1bSJed Brown        appropriate processor during vector assembly).
107*c4762a1bSJed Brown      - Always specify global locations of vector entries.
108*c4762a1bSJed Brown   */
109*c4762a1bSJed Brown   ierr = VecGetLocalSize(x,&ldim);CHKERRQ(ierr);
110*c4762a1bSJed Brown   for (i=0; i<ldim; i++) {
111*c4762a1bSJed Brown     iglobal = i + low;
112*c4762a1bSJed Brown     v       = (PetscScalar)(i + 100*rank);
113*c4762a1bSJed Brown     ierr    = VecSetValues(u,1,&iglobal,&v,INSERT_VALUES);CHKERRQ(ierr);
114*c4762a1bSJed Brown   }
115*c4762a1bSJed Brown 
116*c4762a1bSJed Brown   /*
117*c4762a1bSJed Brown      Assemble vector, using the 2-step process:
118*c4762a1bSJed Brown        VecAssemblyBegin(), VecAssemblyEnd()
119*c4762a1bSJed Brown      Computations can be done while messages are in transition,
120*c4762a1bSJed Brown      by placing code between these two statements.
121*c4762a1bSJed Brown   */
122*c4762a1bSJed Brown   ierr = VecAssemblyBegin(u);CHKERRQ(ierr);
123*c4762a1bSJed Brown   ierr = VecAssemblyEnd(u);CHKERRQ(ierr);
124*c4762a1bSJed Brown 
125*c4762a1bSJed Brown   /* Compute right-hand-side vector */
126*c4762a1bSJed Brown   ierr = MatMult(C,u,b);CHKERRQ(ierr);
127*c4762a1bSJed Brown 
128*c4762a1bSJed Brown   ierr    = MatGetOrdering(C,MATORDERINGNATURAL,&perm,&iperm);CHKERRQ(ierr);
129*c4762a1bSJed Brown   its_max = 2000;
130*c4762a1bSJed Brown   for (i=0; i<its_max; i++) {
131*c4762a1bSJed Brown     ierr = MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&F);CHKERRQ(ierr);
132*c4762a1bSJed Brown     ierr = MatLUFactorSymbolic(F,C,perm,iperm,&factinfo);CHKERRQ(ierr);
133*c4762a1bSJed Brown     for (j=0; j<1; j++) {
134*c4762a1bSJed Brown       ierr = MatLUFactorNumeric(F,C,&factinfo);CHKERRQ(ierr);
135*c4762a1bSJed Brown     }
136*c4762a1bSJed Brown     ierr = MatSolve(F,b,x);CHKERRQ(ierr);
137*c4762a1bSJed Brown     ierr = MatDestroy(&F);CHKERRQ(ierr);
138*c4762a1bSJed Brown   }
139*c4762a1bSJed Brown   ierr = ISDestroy(&perm);CHKERRQ(ierr);
140*c4762a1bSJed Brown   ierr = ISDestroy(&iperm);CHKERRQ(ierr);
141*c4762a1bSJed Brown 
142*c4762a1bSJed Brown   /* Check the error */
143*c4762a1bSJed Brown   ierr = VecAXPY(x,none,u);CHKERRQ(ierr);
144*c4762a1bSJed Brown   ierr = VecNorm(x,NORM_2,&norm);CHKERRQ(ierr);
145*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Norm of error %t\n",(double)norm);CHKERRQ(ierr);
146*c4762a1bSJed Brown 
147*c4762a1bSJed Brown   /* Free work space. */
148*c4762a1bSJed Brown   ierr = VecDestroy(&u);CHKERRQ(ierr);
149*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
150*c4762a1bSJed Brown   ierr = VecDestroy(&b);CHKERRQ(ierr);
151*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
152*c4762a1bSJed Brown   ierr = PetscFinalize();
153*c4762a1bSJed Brown   return ierr;
154*c4762a1bSJed Brown }
155*c4762a1bSJed Brown 
156*c4762a1bSJed Brown 
157