xref: /petsc/src/mat/tests/ex106.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Test repeated LU factorizations. Used for checking memory leak\n\
3c4762a1bSJed Brown   -m <size> : problem size\n\
4c4762a1bSJed Brown   -mat_nonsym : use nonsymmetric matrix (default is symmetric)\n\n";
5c4762a1bSJed Brown 
6c4762a1bSJed Brown #include <petscmat.h>
7c4762a1bSJed Brown int main(int argc,char **args)
8c4762a1bSJed Brown {
9c4762a1bSJed Brown   Mat            C,F;                /* matrix */
10c4762a1bSJed Brown   Vec            x,u,b;          /* approx solution, RHS, exact solution */
11c4762a1bSJed Brown   PetscReal      norm;             /* norm of solution error */
12c4762a1bSJed Brown   PetscScalar    v,none = -1.0;
13c4762a1bSJed Brown   PetscInt       I,J,ldim,low,high,iglobal,Istart,Iend;
14c4762a1bSJed Brown   PetscInt       i,j,m = 3,n = 2,its;
15c4762a1bSJed Brown   PetscMPIInt    size,rank;
16c4762a1bSJed Brown   PetscBool      mat_nonsymmetric;
17c4762a1bSJed Brown   PetscInt       its_max;
18c4762a1bSJed Brown   MatFactorInfo  factinfo;
19c4762a1bSJed Brown   IS             perm,iperm;
20c4762a1bSJed Brown 
21*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
235f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
245f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
25c4762a1bSJed Brown   n    = 2*size;
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   /*
28c4762a1bSJed Brown      Set flag if we are doing a nonsymmetric problem; the default is symmetric.
29c4762a1bSJed Brown   */
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-mat_nonsym",&mat_nonsymmetric));
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   /*
33c4762a1bSJed Brown      Create parallel matrix, specifying only its global dimensions.
34c4762a1bSJed Brown      When using MatCreate(), the matrix format can be specified at
35c4762a1bSJed Brown      runtime. Also, the parallel partitioning of the matrix is
36c4762a1bSJed Brown      determined by PETSc at runtime.
37c4762a1bSJed Brown   */
385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(C,&Istart,&Iend));
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /*
44c4762a1bSJed Brown      Set matrix entries matrix in parallel.
45c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
46c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
47c4762a1bSJed Brown         appropriate processor during matrix assembly).
48c4762a1bSJed Brown       - Always specify global row and columns of matrix entries.
49c4762a1bSJed Brown   */
50c4762a1bSJed Brown   for (I=Istart; I<Iend; I++) {
51c4762a1bSJed Brown     v = -1.0; i = I/n; j = I - i*n;
525f80ce2aSJacob Faibussowitsch     if (i>0)   {J = I - n; CHKERRQ(MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES));}
535f80ce2aSJacob Faibussowitsch     if (i<m-1) {J = I + n; CHKERRQ(MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES));}
545f80ce2aSJacob Faibussowitsch     if (j>0)   {J = I - 1; CHKERRQ(MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES));}
555f80ce2aSJacob Faibussowitsch     if (j<n-1) {J = I + 1; CHKERRQ(MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES));}
565f80ce2aSJacob Faibussowitsch     v = 4.0; CHKERRQ(MatSetValues(C,1,&I,1,&I,&v,ADD_VALUES));
57c4762a1bSJed Brown   }
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /*
60c4762a1bSJed Brown      Make the matrix nonsymmetric if desired
61c4762a1bSJed Brown   */
62c4762a1bSJed Brown   if (mat_nonsymmetric) {
63c4762a1bSJed Brown     for (I=Istart; I<Iend; I++) {
64c4762a1bSJed Brown       v = -1.5; i = I/n;
655f80ce2aSJacob Faibussowitsch       if (i>1)   {J = I-n-1; CHKERRQ(MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES));}
66c4762a1bSJed Brown     }
67c4762a1bSJed Brown   } else {
685f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOption(C,MAT_SYMMETRIC,PETSC_TRUE));
695f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetOption(C,MAT_SYMMETRY_ETERNAL,PETSC_TRUE));
70c4762a1bSJed Brown   }
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /*
73c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
74c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
75c4762a1bSJed Brown      Computations can be done while messages are in transition
76c4762a1bSJed Brown      by placing code between these two statements.
77c4762a1bSJed Brown   */
785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   its_max=1000;
82c4762a1bSJed Brown   /*
83c4762a1bSJed Brown      Create parallel vectors.
84c4762a1bSJed Brown       - When using VecSetSizes(), we specify only the vector's global
85c4762a1bSJed Brown         dimension; the parallel partitioning is determined at runtime.
86c4762a1bSJed Brown       - Note: We form 1 vector from scratch and then duplicate as needed.
87c4762a1bSJed Brown   */
885f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&u));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetSizes(u,PETSC_DECIDE,m*n));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetFromOptions(u));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(u,&b));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(b,&x));
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   /*
95c4762a1bSJed Brown      Currently, all parallel PETSc vectors are partitioned by
96c4762a1bSJed Brown      contiguous chunks across the processors.  Determine which
97c4762a1bSJed Brown      range of entries are locally owned.
98c4762a1bSJed Brown   */
995f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetOwnershipRange(x,&low,&high));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /*
102c4762a1bSJed Brown     Set elements within the exact solution vector in parallel.
103c4762a1bSJed Brown      - Each processor needs to insert only elements that it owns
104c4762a1bSJed Brown        locally (but any non-local entries will be sent to the
105c4762a1bSJed Brown        appropriate processor during vector assembly).
106c4762a1bSJed Brown      - Always specify global locations of vector entries.
107c4762a1bSJed Brown   */
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetLocalSize(x,&ldim));
109c4762a1bSJed Brown   for (i=0; i<ldim; i++) {
110c4762a1bSJed Brown     iglobal = i + low;
111c4762a1bSJed Brown     v       = (PetscScalar)(i + 100*rank);
1125f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetValues(u,1,&iglobal,&v,INSERT_VALUES));
113c4762a1bSJed Brown   }
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   /*
116c4762a1bSJed Brown      Assemble vector, using the 2-step process:
117c4762a1bSJed Brown        VecAssemblyBegin(), VecAssemblyEnd()
118c4762a1bSJed Brown      Computations can be done while messages are in transition,
119c4762a1bSJed Brown      by placing code between these two statements.
120c4762a1bSJed Brown   */
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyBegin(u));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAssemblyEnd(u));
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   /* Compute right-hand-side vector */
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(C,u,b));
126c4762a1bSJed Brown 
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOrdering(C,MATORDERINGNATURAL,&perm,&iperm));
128c4762a1bSJed Brown   its_max = 2000;
129c4762a1bSJed Brown   for (i=0; i<its_max; i++) {
1305f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&F));
1315f80ce2aSJacob Faibussowitsch     CHKERRQ(MatLUFactorSymbolic(F,C,perm,iperm,&factinfo));
132c4762a1bSJed Brown     for (j=0; j<1; j++) {
1335f80ce2aSJacob Faibussowitsch       CHKERRQ(MatLUFactorNumeric(F,C,&factinfo));
134c4762a1bSJed Brown     }
1355f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSolve(F,b,x));
1365f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&F));
137c4762a1bSJed Brown   }
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&perm));
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&iperm));
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* Check the error */
1425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(x,none,u));
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(x,NORM_2,&norm));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %t\n",(double)norm));
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   /* Free work space. */
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&u));
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
151*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
152*b122ec5aSJacob Faibussowitsch   return 0;
153c4762a1bSJed Brown }
154