1c4762a1bSJed Brown static char help[] = "Tests reusing MPI parallel matrices and MatGetValues().\n\
2c4762a1bSJed Brown To test the parallel matrix assembly, this example intentionally lays out\n\
3c4762a1bSJed Brown the matrix across processors differently from the way it is assembled.\n\
4c4762a1bSJed Brown This example uses bilinear elements on the unit square. Input arguments are:\n\
5c4762a1bSJed Brown -m <size> : problem size\n\n";
6c4762a1bSJed Brown
7c4762a1bSJed Brown #include <petscmat.h>
8c4762a1bSJed Brown
FormElementStiffness(PetscReal H,PetscScalar * Ke)93ba16761SJacob Faibussowitsch PetscErrorCode FormElementStiffness(PetscReal H, PetscScalar *Ke)
10d71ae5a4SJacob Faibussowitsch {
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;
283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
29c4762a1bSJed Brown }
30c4762a1bSJed Brown
main(int argc,char ** args)31d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
32d71ae5a4SJacob Faibussowitsch {
33c4762a1bSJed Brown Mat C;
34c4762a1bSJed Brown Vec u, b;
35c4762a1bSJed Brown PetscMPIInt size, rank;
36c4762a1bSJed Brown PetscInt i, m = 5, N, start, end, M, idx[4];
37c4762a1bSJed Brown PetscInt j, nrsub, ncsub, *rsub, *csub, mystart, myend;
38c4762a1bSJed Brown PetscBool flg;
39c4762a1bSJed Brown PetscScalar one = 1.0, Ke[16], *vals;
40c4762a1bSJed Brown PetscReal h, norm;
41c4762a1bSJed Brown
42327415f7SBarry Smith PetscFunctionBeginUser;
43c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
449566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-m", &m, NULL));
45c4762a1bSJed Brown
46c4762a1bSJed Brown N = (m + 1) * (m + 1); /* dimension of matrix */
47c4762a1bSJed Brown M = m * m; /* number of elements */
48c4762a1bSJed Brown h = 1.0 / m; /* mesh width */
499566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
51c4762a1bSJed Brown
52c4762a1bSJed Brown /* Create stiffness matrix */
539566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
549566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, N, N));
559566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C));
569566063dSJacob Faibussowitsch PetscCall(MatSetUp(C));
57c4762a1bSJed Brown
58c4762a1bSJed Brown start = rank * (M / size) + ((M % size) < rank ? (M % size) : rank);
59c4762a1bSJed Brown end = start + M / size + ((M % size) > rank);
60c4762a1bSJed Brown
61c4762a1bSJed Brown /* Form the element stiffness for the Laplacian */
629566063dSJacob Faibussowitsch PetscCall(FormElementStiffness(h * h, Ke));
63c4762a1bSJed Brown for (i = start; i < end; i++) {
64c4762a1bSJed Brown /* location of lower left corner of element */
65c4762a1bSJed Brown /* node numbers for the four corners of element */
66c4762a1bSJed Brown idx[0] = (m + 1) * (i / m) + (i % m);
679371c9d4SSatish Balay idx[1] = idx[0] + 1;
689371c9d4SSatish Balay idx[2] = idx[1] + m + 1;
699371c9d4SSatish Balay idx[3] = idx[2] - 1;
709566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 4, idx, 4, idx, Ke, ADD_VALUES));
71c4762a1bSJed Brown }
729566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
739566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
74c4762a1bSJed Brown
75c4762a1bSJed Brown /* Assemble the matrix again */
769566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(C));
77c4762a1bSJed Brown
78c4762a1bSJed Brown for (i = start; i < end; i++) {
79c4762a1bSJed Brown /* location of lower left corner of element */
80c4762a1bSJed Brown /* node numbers for the four corners of element */
81c4762a1bSJed Brown idx[0] = (m + 1) * (i / m) + (i % m);
829371c9d4SSatish Balay idx[1] = idx[0] + 1;
839371c9d4SSatish Balay idx[2] = idx[1] + m + 1;
849371c9d4SSatish Balay idx[3] = idx[2] - 1;
859566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 4, idx, 4, idx, Ke, ADD_VALUES));
86c4762a1bSJed Brown }
879566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
889566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
89c4762a1bSJed Brown
90c4762a1bSJed Brown /* Create test vectors */
919566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD, &u));
929566063dSJacob Faibussowitsch PetscCall(VecSetSizes(u, PETSC_DECIDE, N));
939566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(u));
949566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &b));
959566063dSJacob Faibussowitsch PetscCall(VecSet(u, one));
96c4762a1bSJed Brown
97c4762a1bSJed Brown /* Check error */
989566063dSJacob Faibussowitsch PetscCall(MatMult(C, u, b));
999566063dSJacob Faibussowitsch PetscCall(VecNorm(b, NORM_2, &norm));
10048a46eb9SPierre Jolivet if (norm > PETSC_SQRT_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Norm of error b %g should be near 0\n", (double)norm));
101c4762a1bSJed Brown
102c4762a1bSJed Brown /* Now test MatGetValues() */
1039566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-get_values", &flg));
104c4762a1bSJed Brown if (flg) {
1059566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C, &mystart, &myend));
1069371c9d4SSatish Balay nrsub = myend - mystart;
1079371c9d4SSatish Balay ncsub = 4;
1089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrsub * ncsub, &vals));
1099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nrsub, &rsub));
1109566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(ncsub, &csub));
111c4762a1bSJed Brown for (i = myend - 1; i >= mystart; i--) rsub[myend - i - 1] = i;
112c4762a1bSJed Brown for (i = 0; i < ncsub; i++) csub[i] = 2 * (ncsub - i) + mystart;
1139566063dSJacob Faibussowitsch PetscCall(MatGetValues(C, nrsub, rsub, ncsub, csub, vals));
1149566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
1159566063dSJacob 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));
116c4762a1bSJed Brown for (i = 0; i < nrsub; i++) {
117c4762a1bSJed Brown for (j = 0; j < ncsub; j++) {
118c4762a1bSJed Brown if (PetscImaginaryPart(vals[i * ncsub + j]) != 0.0) {
1199566063dSJacob 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])));
120c4762a1bSJed Brown } else {
1219566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, " C[%" PetscInt_FMT ", %" PetscInt_FMT "] = %g\n", rsub[i], csub[j], (double)PetscRealPart(vals[i * ncsub + j])));
122c4762a1bSJed Brown }
123c4762a1bSJed Brown }
124c4762a1bSJed Brown }
1259566063dSJacob Faibussowitsch PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
1269566063dSJacob Faibussowitsch PetscCall(PetscFree(rsub));
1279566063dSJacob Faibussowitsch PetscCall(PetscFree(csub));
1289566063dSJacob Faibussowitsch PetscCall(PetscFree(vals));
129c4762a1bSJed Brown }
130c4762a1bSJed Brown
131c4762a1bSJed Brown /* Free data structures */
1329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u));
1339566063dSJacob Faibussowitsch PetscCall(VecDestroy(&b));
1349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
1359566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
136b122ec5aSJacob Faibussowitsch return 0;
137c4762a1bSJed Brown }
138c4762a1bSJed Brown
139c4762a1bSJed Brown /*TEST
140c4762a1bSJed Brown
141c4762a1bSJed Brown test:
142c4762a1bSJed Brown nsize: 4
143*3886731fSPierre Jolivet output_file: output/empty.out
144c4762a1bSJed Brown
145c4762a1bSJed Brown TEST*/
146