xref: /petsc/src/mat/tests/ex19.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests reusing MPI parallel matrices and MatGetValues().\n\
3c4762a1bSJed Brown To test the parallel matrix assembly, this example intentionally lays out\n\
4c4762a1bSJed Brown the matrix across processors differently from the way it is assembled.\n\
5c4762a1bSJed Brown This example uses bilinear elements on the unit square.  Input arguments are:\n\
6c4762a1bSJed Brown   -m <size> : problem size\n\n";
7c4762a1bSJed Brown 
8c4762a1bSJed Brown #include <petscmat.h>
9c4762a1bSJed Brown 
10*3ba16761SJacob Faibussowitsch PetscErrorCode FormElementStiffness(PetscReal H, PetscScalar *Ke)
11d71ae5a4SJacob Faibussowitsch {
12c4762a1bSJed Brown   PetscFunctionBegin;
139371c9d4SSatish Balay   Ke[0]  = H / 6.0;
149371c9d4SSatish Balay   Ke[1]  = -.125 * H;
159371c9d4SSatish Balay   Ke[2]  = H / 12.0;
169371c9d4SSatish Balay   Ke[3]  = -.125 * H;
179371c9d4SSatish Balay   Ke[4]  = -.125 * H;
189371c9d4SSatish Balay   Ke[5]  = H / 6.0;
199371c9d4SSatish Balay   Ke[6]  = -.125 * H;
209371c9d4SSatish Balay   Ke[7]  = H / 12.0;
219371c9d4SSatish Balay   Ke[8]  = H / 12.0;
229371c9d4SSatish Balay   Ke[9]  = -.125 * H;
239371c9d4SSatish Balay   Ke[10] = H / 6.0;
249371c9d4SSatish Balay   Ke[11] = -.125 * H;
259371c9d4SSatish Balay   Ke[12] = -.125 * H;
269371c9d4SSatish Balay   Ke[13] = H / 12.0;
279371c9d4SSatish Balay   Ke[14] = -.125 * H;
289371c9d4SSatish Balay   Ke[15] = H / 6.0;
29*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
30c4762a1bSJed Brown }
31c4762a1bSJed Brown 
32d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
33d71ae5a4SJacob Faibussowitsch {
34c4762a1bSJed Brown   Mat         C;
35c4762a1bSJed Brown   Vec         u, b;
36c4762a1bSJed Brown   PetscMPIInt size, rank;
37c4762a1bSJed Brown   PetscInt    i, m = 5, N, start, end, M, idx[4];
38c4762a1bSJed Brown   PetscInt    j, nrsub, ncsub, *rsub, *csub, mystart, myend;
39c4762a1bSJed Brown   PetscBool   flg;
40c4762a1bSJed Brown   PetscScalar one = 1.0, Ke[16], *vals;
41c4762a1bSJed Brown   PetscReal   h, norm;
42c4762a1bSJed Brown 
43327415f7SBarry Smith   PetscFunctionBeginUser;
449566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
46c4762a1bSJed Brown 
47c4762a1bSJed Brown   N = (m + 1) * (m + 1); /* dimension of matrix */
48c4762a1bSJed Brown   M = m * m;             /* number of elements */
49c4762a1bSJed Brown   h = 1.0 / m;           /* mesh width */
509566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
519566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   /* Create stiffness matrix */
549566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
559566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, N, N));
569566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
579566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   start = rank * (M / size) + ((M % size) < rank ? (M % size) : rank);
60c4762a1bSJed Brown   end   = start + M / size + ((M % size) > rank);
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /* Form the element stiffness for the Laplacian */
639566063dSJacob Faibussowitsch   PetscCall(FormElementStiffness(h * h, Ke));
64c4762a1bSJed Brown   for (i = start; i < end; i++) {
65c4762a1bSJed Brown     /* location of lower left corner of element */
66c4762a1bSJed Brown     /* node numbers for the four corners of element */
67c4762a1bSJed Brown     idx[0] = (m + 1) * (i / m) + (i % m);
689371c9d4SSatish Balay     idx[1] = idx[0] + 1;
699371c9d4SSatish Balay     idx[2] = idx[1] + m + 1;
709371c9d4SSatish Balay     idx[3] = idx[2] - 1;
719566063dSJacob Faibussowitsch     PetscCall(MatSetValues(C, 4, idx, 4, idx, Ke, ADD_VALUES));
72c4762a1bSJed Brown   }
739566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* Assemble the matrix again */
779566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(C));
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   for (i = start; i < end; i++) {
80c4762a1bSJed Brown     /* location of lower left corner of element */
81c4762a1bSJed Brown     /* node numbers for the four corners of element */
82c4762a1bSJed Brown     idx[0] = (m + 1) * (i / m) + (i % m);
839371c9d4SSatish Balay     idx[1] = idx[0] + 1;
849371c9d4SSatish Balay     idx[2] = idx[1] + m + 1;
859371c9d4SSatish Balay     idx[3] = idx[2] - 1;
869566063dSJacob Faibussowitsch     PetscCall(MatSetValues(C, 4, idx, 4, idx, Ke, ADD_VALUES));
87c4762a1bSJed Brown   }
889566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
899566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /* Create test vectors */
929566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &u));
939566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(u, PETSC_DECIDE, N));
949566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(u));
959566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &b));
969566063dSJacob Faibussowitsch   PetscCall(VecSet(u, one));
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   /* Check error */
999566063dSJacob Faibussowitsch   PetscCall(MatMult(C, u, b));
1009566063dSJacob Faibussowitsch   PetscCall(VecNorm(b, NORM_2, &norm));
10148a46eb9SPierre Jolivet   if (norm > PETSC_SQRT_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error b %g should be near 0\n", (double)norm));
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   /* Now test MatGetValues() */
1049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-get_values", &flg));
105c4762a1bSJed Brown   if (flg) {
1069566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(C, &mystart, &myend));
1079371c9d4SSatish Balay     nrsub = myend - mystart;
1089371c9d4SSatish Balay     ncsub = 4;
1099566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrsub * ncsub, &vals));
1109566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrsub, &rsub));
1119566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(ncsub, &csub));
112c4762a1bSJed Brown     for (i = myend - 1; i >= mystart; i--) rsub[myend - i - 1] = i;
113c4762a1bSJed Brown     for (i = 0; i < ncsub; i++) csub[i] = 2 * (ncsub - i) + mystart;
1149566063dSJacob Faibussowitsch     PetscCall(MatGetValues(C, nrsub, rsub, ncsub, csub, vals));
1159566063dSJacob Faibussowitsch     PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
1169566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "processor number %d: start=%" PetscInt_FMT ", end=%" PetscInt_FMT ", mystart=%" PetscInt_FMT ", myend=%" PetscInt_FMT "\n", rank, start, end, mystart, myend));
117c4762a1bSJed Brown     for (i = 0; i < nrsub; i++) {
118c4762a1bSJed Brown       for (j = 0; j < ncsub; j++) {
119c4762a1bSJed Brown         if (PetscImaginaryPart(vals[i * ncsub + j]) != 0.0) {
1209566063dSJacob Faibussowitsch           PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "  C[%" PetscInt_FMT ", %" PetscInt_FMT "] = %g + %g i\n", rsub[i], csub[j], (double)PetscRealPart(vals[i * ncsub + j]), (double)PetscImaginaryPart(vals[i * ncsub + j])));
121c4762a1bSJed Brown         } else {
1229566063dSJacob Faibussowitsch           PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "  C[%" PetscInt_FMT ", %" PetscInt_FMT "] = %g\n", rsub[i], csub[j], (double)PetscRealPart(vals[i * ncsub + j])));
123c4762a1bSJed Brown         }
124c4762a1bSJed Brown       }
125c4762a1bSJed Brown     }
1269566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
1279566063dSJacob Faibussowitsch     PetscCall(PetscFree(rsub));
1289566063dSJacob Faibussowitsch     PetscCall(PetscFree(csub));
1299566063dSJacob Faibussowitsch     PetscCall(PetscFree(vals));
130c4762a1bSJed Brown   }
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   /* Free data structures */
1339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
1349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
1359566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1369566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
137b122ec5aSJacob Faibussowitsch   return 0;
138c4762a1bSJed Brown }
139c4762a1bSJed Brown 
140c4762a1bSJed Brown /*TEST
141c4762a1bSJed Brown 
142c4762a1bSJed Brown    test:
143c4762a1bSJed Brown       nsize: 4
144c4762a1bSJed Brown 
145c4762a1bSJed Brown TEST*/
146