xref: /petsc/src/mat/tests/ex106.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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*327415f7SBarry Smith   PetscFunctionBeginUser;
229566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
249566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
26c4762a1bSJed Brown   n    = 2*size;
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /*
29c4762a1bSJed Brown      Set flag if we are doing a nonsymmetric problem; the default is symmetric.
30c4762a1bSJed Brown   */
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL,NULL,"-mat_nonsym",&mat_nonsymmetric));
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   /*
34c4762a1bSJed Brown      Create parallel matrix, specifying only its global dimensions.
35c4762a1bSJed Brown      When using MatCreate(), the matrix format can be specified at
36c4762a1bSJed Brown      runtime. Also, the parallel partitioning of the matrix is
37c4762a1bSJed Brown      determined by PETSc at runtime.
38c4762a1bSJed Brown   */
399566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
409566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n));
419566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
429566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /*
45c4762a1bSJed Brown      Set matrix entries matrix in parallel.
46c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
47c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
48c4762a1bSJed Brown         appropriate processor during matrix assembly).
49c4762a1bSJed Brown       - Always specify global row and columns of matrix entries.
50c4762a1bSJed Brown   */
51c4762a1bSJed Brown   for (I=Istart; I<Iend; I++) {
52c4762a1bSJed Brown     v = -1.0; i = I/n; j = I - i*n;
539566063dSJacob Faibussowitsch     if (i>0)   {J = I - n; PetscCall(MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES));}
549566063dSJacob Faibussowitsch     if (i<m-1) {J = I + n; PetscCall(MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES));}
559566063dSJacob Faibussowitsch     if (j>0)   {J = I - 1; PetscCall(MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES));}
569566063dSJacob Faibussowitsch     if (j<n-1) {J = I + 1; PetscCall(MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES));}
579566063dSJacob Faibussowitsch     v = 4.0; PetscCall(MatSetValues(C,1,&I,1,&I,&v,ADD_VALUES));
58c4762a1bSJed Brown   }
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /*
61c4762a1bSJed Brown      Make the matrix nonsymmetric if desired
62c4762a1bSJed Brown   */
63c4762a1bSJed Brown   if (mat_nonsymmetric) {
64c4762a1bSJed Brown     for (I=Istart; I<Iend; I++) {
65c4762a1bSJed Brown       v = -1.5; i = I/n;
669566063dSJacob Faibussowitsch       if (i>1)   {J = I-n-1; PetscCall(MatSetValues(C,1,&I,1,&J,&v,ADD_VALUES));}
67c4762a1bSJed Brown     }
68c4762a1bSJed Brown   } else {
699566063dSJacob Faibussowitsch     PetscCall(MatSetOption(C,MAT_SYMMETRIC,PETSC_TRUE));
709566063dSJacob Faibussowitsch     PetscCall(MatSetOption(C,MAT_SYMMETRY_ETERNAL,PETSC_TRUE));
71c4762a1bSJed Brown   }
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   /*
74c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
75c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd()
76c4762a1bSJed Brown      Computations can be done while messages are in transition
77c4762a1bSJed Brown      by placing code between these two statements.
78c4762a1bSJed Brown   */
799566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
809566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   its_max=1000;
83c4762a1bSJed Brown   /*
84c4762a1bSJed Brown      Create parallel vectors.
85c4762a1bSJed Brown       - When using VecSetSizes(), we specify only the vector's global
86c4762a1bSJed Brown         dimension; the parallel partitioning is determined at runtime.
87c4762a1bSJed Brown       - Note: We form 1 vector from scratch and then duplicate as needed.
88c4762a1bSJed Brown   */
899566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&u));
909566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(u,PETSC_DECIDE,m*n));
919566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(u));
929566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u,&b));
939566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(b,&x));
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /*
96c4762a1bSJed Brown      Currently, all parallel PETSc vectors are partitioned by
97c4762a1bSJed Brown      contiguous chunks across the processors.  Determine which
98c4762a1bSJed Brown      range of entries are locally owned.
99c4762a1bSJed Brown   */
1009566063dSJacob Faibussowitsch   PetscCall(VecGetOwnershipRange(x,&low,&high));
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   /*
103c4762a1bSJed Brown     Set elements within the exact solution vector in parallel.
104c4762a1bSJed Brown      - Each processor needs to insert only elements that it owns
105c4762a1bSJed Brown        locally (but any non-local entries will be sent to the
106c4762a1bSJed Brown        appropriate processor during vector assembly).
107c4762a1bSJed Brown      - Always specify global locations of vector entries.
108c4762a1bSJed Brown   */
1099566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(x,&ldim));
110c4762a1bSJed Brown   for (i=0; i<ldim; i++) {
111c4762a1bSJed Brown     iglobal = i + low;
112c4762a1bSJed Brown     v       = (PetscScalar)(i + 100*rank);
1139566063dSJacob Faibussowitsch     PetscCall(VecSetValues(u,1,&iglobal,&v,INSERT_VALUES));
114c4762a1bSJed Brown   }
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /*
117c4762a1bSJed Brown      Assemble vector, using the 2-step process:
118c4762a1bSJed Brown        VecAssemblyBegin(), VecAssemblyEnd()
119c4762a1bSJed Brown      Computations can be done while messages are in transition,
120c4762a1bSJed Brown      by placing code between these two statements.
121c4762a1bSJed Brown   */
1229566063dSJacob Faibussowitsch   PetscCall(VecAssemblyBegin(u));
1239566063dSJacob Faibussowitsch   PetscCall(VecAssemblyEnd(u));
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* Compute right-hand-side vector */
1269566063dSJacob Faibussowitsch   PetscCall(MatMult(C,u,b));
127c4762a1bSJed Brown 
1289566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(C,MATORDERINGNATURAL,&perm,&iperm));
129c4762a1bSJed Brown   its_max = 2000;
130c4762a1bSJed Brown   for (i=0; i<its_max; i++) {
1319566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(C,MATSOLVERPETSC,MAT_FACTOR_LU,&F));
1329566063dSJacob Faibussowitsch     PetscCall(MatLUFactorSymbolic(F,C,perm,iperm,&factinfo));
133c4762a1bSJed Brown     for (j=0; j<1; j++) {
1349566063dSJacob Faibussowitsch       PetscCall(MatLUFactorNumeric(F,C,&factinfo));
135c4762a1bSJed Brown     }
1369566063dSJacob Faibussowitsch     PetscCall(MatSolve(F,b,x));
1379566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&F));
138c4762a1bSJed Brown   }
1399566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
1409566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iperm));
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   /* Check the error */
1439566063dSJacob Faibussowitsch   PetscCall(VecAXPY(x,none,u));
1449566063dSJacob Faibussowitsch   PetscCall(VecNorm(x,NORM_2,&norm));
1459566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Norm of error %t\n",(double)norm));
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /* Free work space. */
1489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
1499566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1509566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
1519566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1529566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
153b122ec5aSJacob Faibussowitsch   return 0;
154c4762a1bSJed Brown }
155