xref: /petsc/src/mat/tests/ex19.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
109371c9d4SSatish Balay int FormElementStiffness(PetscReal H, PetscScalar *Ke) {
11c4762a1bSJed Brown   PetscFunctionBegin;
129371c9d4SSatish Balay   Ke[0]  = H / 6.0;
139371c9d4SSatish Balay   Ke[1]  = -.125 * H;
149371c9d4SSatish Balay   Ke[2]  = H / 12.0;
159371c9d4SSatish Balay   Ke[3]  = -.125 * H;
169371c9d4SSatish Balay   Ke[4]  = -.125 * H;
179371c9d4SSatish Balay   Ke[5]  = H / 6.0;
189371c9d4SSatish Balay   Ke[6]  = -.125 * H;
199371c9d4SSatish Balay   Ke[7]  = H / 12.0;
209371c9d4SSatish Balay   Ke[8]  = H / 12.0;
219371c9d4SSatish Balay   Ke[9]  = -.125 * H;
229371c9d4SSatish Balay   Ke[10] = H / 6.0;
239371c9d4SSatish Balay   Ke[11] = -.125 * H;
249371c9d4SSatish Balay   Ke[12] = -.125 * H;
259371c9d4SSatish Balay   Ke[13] = H / 12.0;
269371c9d4SSatish Balay   Ke[14] = -.125 * H;
279371c9d4SSatish Balay   Ke[15] = H / 6.0;
28c4762a1bSJed Brown   PetscFunctionReturn(0);
29c4762a1bSJed Brown }
30c4762a1bSJed Brown 
319371c9d4SSatish Balay int main(int argc, char **args) {
32c4762a1bSJed Brown   Mat         C;
33c4762a1bSJed Brown   Vec         u, b;
34c4762a1bSJed Brown   PetscMPIInt size, rank;
35c4762a1bSJed Brown   PetscInt    i, m = 5, N, start, end, M, idx[4];
36c4762a1bSJed Brown   PetscInt    j, nrsub, ncsub, *rsub, *csub, mystart, myend;
37c4762a1bSJed Brown   PetscBool   flg;
38c4762a1bSJed Brown   PetscScalar one = 1.0, Ke[16], *vals;
39c4762a1bSJed Brown   PetscReal   h, norm;
40c4762a1bSJed Brown 
41327415f7SBarry Smith   PetscFunctionBeginUser;
429566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   N = (m + 1) * (m + 1); /* dimension of matrix */
46c4762a1bSJed Brown   M = m * m;             /* number of elements */
47c4762a1bSJed Brown   h = 1.0 / m;           /* mesh width */
489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
499566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* Create stiffness matrix */
529566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
539566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, N, N));
549566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(C));
559566063dSJacob Faibussowitsch   PetscCall(MatSetUp(C));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   start = rank * (M / size) + ((M % size) < rank ? (M % size) : rank);
58c4762a1bSJed Brown   end   = start + M / size + ((M % size) > rank);
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   /* Form the element stiffness for the Laplacian */
619566063dSJacob Faibussowitsch   PetscCall(FormElementStiffness(h * h, Ke));
62c4762a1bSJed Brown   for (i = start; i < end; i++) {
63c4762a1bSJed Brown     /* location of lower left corner of element */
64c4762a1bSJed Brown     /* node numbers for the four corners of element */
65c4762a1bSJed Brown     idx[0] = (m + 1) * (i / m) + (i % m);
669371c9d4SSatish Balay     idx[1] = idx[0] + 1;
679371c9d4SSatish Balay     idx[2] = idx[1] + m + 1;
689371c9d4SSatish Balay     idx[3] = idx[2] - 1;
699566063dSJacob Faibussowitsch     PetscCall(MatSetValues(C, 4, idx, 4, idx, Ke, ADD_VALUES));
70c4762a1bSJed Brown   }
719566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
729566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   /* Assemble the matrix again */
759566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(C));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   for (i = start; i < end; i++) {
78c4762a1bSJed Brown     /* location of lower left corner of element */
79c4762a1bSJed Brown     /* node numbers for the four corners of element */
80c4762a1bSJed Brown     idx[0] = (m + 1) * (i / m) + (i % m);
819371c9d4SSatish Balay     idx[1] = idx[0] + 1;
829371c9d4SSatish Balay     idx[2] = idx[1] + m + 1;
839371c9d4SSatish Balay     idx[3] = idx[2] - 1;
849566063dSJacob Faibussowitsch     PetscCall(MatSetValues(C, 4, idx, 4, idx, Ke, ADD_VALUES));
85c4762a1bSJed Brown   }
869566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
879566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   /* Create test vectors */
909566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &u));
919566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(u, PETSC_DECIDE, N));
929566063dSJacob Faibussowitsch   PetscCall(VecSetFromOptions(u));
939566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(u, &b));
949566063dSJacob Faibussowitsch   PetscCall(VecSet(u, one));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* Check error */
979566063dSJacob Faibussowitsch   PetscCall(MatMult(C, u, b));
989566063dSJacob Faibussowitsch   PetscCall(VecNorm(b, NORM_2, &norm));
99*48a46eb9SPierre Jolivet   if (norm > PETSC_SQRT_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error b %g should be near 0\n", (double)norm));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /* Now test MatGetValues() */
1029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-get_values", &flg));
103c4762a1bSJed Brown   if (flg) {
1049566063dSJacob Faibussowitsch     PetscCall(MatGetOwnershipRange(C, &mystart, &myend));
1059371c9d4SSatish Balay     nrsub = myend - mystart;
1069371c9d4SSatish Balay     ncsub = 4;
1079566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrsub * ncsub, &vals));
1089566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nrsub, &rsub));
1099566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(ncsub, &csub));
110c4762a1bSJed Brown     for (i = myend - 1; i >= mystart; i--) rsub[myend - i - 1] = i;
111c4762a1bSJed Brown     for (i = 0; i < ncsub; i++) csub[i] = 2 * (ncsub - i) + mystart;
1129566063dSJacob Faibussowitsch     PetscCall(MatGetValues(C, nrsub, rsub, ncsub, csub, vals));
1139566063dSJacob Faibussowitsch     PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
1149566063dSJacob 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));
115c4762a1bSJed Brown     for (i = 0; i < nrsub; i++) {
116c4762a1bSJed Brown       for (j = 0; j < ncsub; j++) {
117c4762a1bSJed Brown         if (PetscImaginaryPart(vals[i * ncsub + j]) != 0.0) {
1189566063dSJacob 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])));
119c4762a1bSJed Brown         } else {
1209566063dSJacob Faibussowitsch           PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "  C[%" PetscInt_FMT ", %" PetscInt_FMT "] = %g\n", rsub[i], csub[j], (double)PetscRealPart(vals[i * ncsub + j])));
121c4762a1bSJed Brown         }
122c4762a1bSJed Brown       }
123c4762a1bSJed Brown     }
1249566063dSJacob Faibussowitsch     PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
1259566063dSJacob Faibussowitsch     PetscCall(PetscFree(rsub));
1269566063dSJacob Faibussowitsch     PetscCall(PetscFree(csub));
1279566063dSJacob Faibussowitsch     PetscCall(PetscFree(vals));
128c4762a1bSJed Brown   }
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* Free data structures */
1319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
1329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
1339566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&C));
1349566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
135b122ec5aSJacob Faibussowitsch   return 0;
136c4762a1bSJed Brown }
137c4762a1bSJed Brown 
138c4762a1bSJed Brown /*TEST
139c4762a1bSJed Brown 
140c4762a1bSJed Brown    test:
141c4762a1bSJed Brown       nsize: 4
142c4762a1bSJed Brown 
143c4762a1bSJed Brown TEST*/
144