xref: /petsc/src/dm/impls/plex/tutorials/ex20.c (revision 828beda2407cdd56b6d541804943d5bc46cf6ab9)
1*828beda2SMark Adams static char help[] = "Test DMPlexGetLETKFLocalizationMatrix.\n\n";
2*828beda2SMark Adams 
3*828beda2SMark Adams #include <petscdmplex.h>
4*828beda2SMark Adams 
5*828beda2SMark Adams int main(int argc, char **argv)
6*828beda2SMark Adams {
7*828beda2SMark Adams   DM        dm;
8*828beda2SMark Adams   Mat       H, Q;
9*828beda2SMark Adams   PetscInt  numobservations;
10*828beda2SMark Adams   PetscInt  dim = 1, n;
11*828beda2SMark Adams   PetscInt  faces[3];
12*828beda2SMark Adams   PetscReal lower[3] = {0.0, 0.0, 0.0};
13*828beda2SMark Adams   PetscReal upper[3] = {1.0, 1.0, 1.0};
14*828beda2SMark Adams 
15*828beda2SMark Adams   PetscFunctionBeginUser;
16*828beda2SMark Adams   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
17*828beda2SMark Adams 
18*828beda2SMark Adams   /* Get dimension and from options. We need the data here and Plex does not have access functions */
19*828beda2SMark Adams   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dm_plex_dim", &dim, NULL));
20*828beda2SMark Adams   n = 3;
21*828beda2SMark Adams   PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &n, NULL));
22*828beda2SMark Adams   n = 3;
23*828beda2SMark Adams   PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-dm_plex_box_lower", lower, &n, NULL));
24*828beda2SMark Adams   n = 3;
25*828beda2SMark Adams   PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-dm_plex_box_upper", upper, &n, NULL));
26*828beda2SMark Adams 
27*828beda2SMark Adams   /* Create the mesh using DMPlexCreateBoxMesh (could pass parameters) */
28*828beda2SMark Adams   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
29*828beda2SMark Adams   PetscCall(DMSetType(dm, DMPLEX));
30*828beda2SMark Adams   PetscCall(DMSetFromOptions(dm));
31*828beda2SMark Adams   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
32*828beda2SMark Adams 
33*828beda2SMark Adams   /* Verify dimension matches */
34*828beda2SMark Adams   PetscInt dmDim;
35*828beda2SMark Adams   PetscCall(DMGetDimension(dm, &dmDim));
36*828beda2SMark Adams   PetscCheck(dmDim == dim, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "DM dimension %" PetscInt_FMT " does not match requested dimension %" PetscInt_FMT, dmDim, dim);
37*828beda2SMark Adams 
38*828beda2SMark Adams   /* Set number of local observations to use: 3^dim */
39*828beda2SMark Adams   numobservations = 1;
40*828beda2SMark Adams   for (PetscInt d = 0; d < dim && d < 2; d++) numobservations *= 3;
41*828beda2SMark Adams 
42*828beda2SMark Adams   /* Get number of vertices */
43*828beda2SMark Adams   PetscInt vStart, vEnd, numVertices;
44*828beda2SMark Adams   PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
45*828beda2SMark Adams   numVertices = vEnd - vStart;
46*828beda2SMark Adams 
47*828beda2SMark Adams   /* Create a section for vertices (required for Global Point mapping) */
48*828beda2SMark Adams   PetscSection section;
49*828beda2SMark Adams   PetscInt     pStart, pEnd;
50*828beda2SMark Adams   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
51*828beda2SMark Adams   PetscCall(PetscSectionCreate(PETSC_COMM_WORLD, &section));
52*828beda2SMark Adams   PetscCall(PetscSectionSetNumFields(section, 1));
53*828beda2SMark Adams   PetscCall(PetscSectionSetChart(section, pStart, pEnd));
54*828beda2SMark Adams   for (PetscInt v = vStart; v < vEnd; ++v) PetscCall(PetscSectionSetDof(section, v, 1));
55*828beda2SMark Adams   PetscCall(PetscSectionSetUp(section));
56*828beda2SMark Adams   PetscCall(DMSetLocalSection(dm, section));
57*828beda2SMark Adams   PetscCall(PetscSectionDestroy(&section));
58*828beda2SMark Adams 
59*828beda2SMark Adams   /* Create global section */
60*828beda2SMark Adams   PetscSection globalSection;
61*828beda2SMark Adams   PetscCall(DMGetGlobalSection(dm, &globalSection));
62*828beda2SMark Adams 
63*828beda2SMark Adams   /* Count observations (every other vertex in each dimension) */
64*828beda2SMark Adams   PetscInt numlocalobs = 0;
65*828beda2SMark Adams   {
66*828beda2SMark Adams     Vec                coordinates;
67*828beda2SMark Adams     PetscSection       coordSection;
68*828beda2SMark Adams     const PetscScalar *coordArray;
69*828beda2SMark Adams     PetscInt           offset;
70*828beda2SMark Adams 
71*828beda2SMark Adams     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
72*828beda2SMark Adams     PetscCall(DMGetCoordinateSection(dm, &coordSection));
73*828beda2SMark Adams     PetscCall(VecGetArrayRead(coordinates, &coordArray));
74*828beda2SMark Adams 
75*828beda2SMark Adams     for (PetscInt v = vStart; v < vEnd; v++) {
76*828beda2SMark Adams       PetscReal coords[3] = {0.0, 0.0, 0.0};
77*828beda2SMark Adams       PetscBool isObs     = PETSC_TRUE;
78*828beda2SMark Adams 
79*828beda2SMark Adams       PetscCall(PetscSectionGetOffset(coordSection, v, &offset));
80*828beda2SMark Adams       for (PetscInt d = 0; d < dim; d++) coords[d] = PetscRealPart(coordArray[offset + d]);
81*828beda2SMark Adams 
82*828beda2SMark Adams       /* Check if this vertex is at an observation location (every other grid point) */
83*828beda2SMark Adams       for (PetscInt d = 0; d < dim; d++) {
84*828beda2SMark Adams         PetscReal gridSpacing = (upper[d] - lower[d]) / faces[d];
85*828beda2SMark Adams         PetscInt  gridIdx     = (PetscInt)((coords[d] - lower[d]) / gridSpacing + 0.5);
86*828beda2SMark Adams         if (gridIdx % 2 != 0) {
87*828beda2SMark Adams           isObs = PETSC_FALSE;
88*828beda2SMark Adams           break;
89*828beda2SMark Adams         }
90*828beda2SMark Adams       }
91*828beda2SMark Adams       if (isObs) numlocalobs++;
92*828beda2SMark Adams     }
93*828beda2SMark Adams     PetscCall(VecRestoreArrayRead(coordinates, &coordArray));
94*828beda2SMark Adams   }
95*828beda2SMark Adams 
96*828beda2SMark Adams   /* Create H matrix */
97*828beda2SMark Adams   PetscCall(MatCreate(PETSC_COMM_WORLD, &H));
98*828beda2SMark Adams   PetscCall(MatSetSizes(H, numlocalobs, PETSC_DECIDE, PETSC_DECIDE, numVertices));
99*828beda2SMark Adams   PetscCall(MatSetType(H, MATAIJ));
100*828beda2SMark Adams   PetscCall(MatSeqAIJSetPreallocation(H, 1, NULL));
101*828beda2SMark Adams   PetscCall(MatMPIAIJSetPreallocation(H, 1, NULL, 0, NULL));
102*828beda2SMark Adams   PetscCall(PetscObjectSetName((PetscObject)H, "H_observation_operator"));
103*828beda2SMark Adams 
104*828beda2SMark Adams   /* Fill H matrix */
105*828beda2SMark Adams   {
106*828beda2SMark Adams     Vec                coordinates;
107*828beda2SMark Adams     PetscSection       coordSection;
108*828beda2SMark Adams     const PetscScalar *coordArray;
109*828beda2SMark Adams     PetscInt           obsIdx = 0;
110*828beda2SMark Adams     PetscInt           offset;
111*828beda2SMark Adams 
112*828beda2SMark Adams     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
113*828beda2SMark Adams     PetscCall(DMGetCoordinateSection(dm, &coordSection));
114*828beda2SMark Adams     PetscCall(VecGetArrayRead(coordinates, &coordArray));
115*828beda2SMark Adams 
116*828beda2SMark Adams     for (PetscInt v = vStart; v < vEnd; v++) {
117*828beda2SMark Adams       PetscReal coords[3] = {0.0, 0.0, 0.0};
118*828beda2SMark Adams       PetscBool isObs     = PETSC_TRUE;
119*828beda2SMark Adams 
120*828beda2SMark Adams       PetscCall(PetscSectionGetOffset(coordSection, v, &offset));
121*828beda2SMark Adams       for (PetscInt d = 0; d < dim; d++) coords[d] = PetscRealPart(coordArray[offset + d]);
122*828beda2SMark Adams 
123*828beda2SMark Adams       /* Check if this vertex is at an observation location (every other grid point) */
124*828beda2SMark Adams       for (PetscInt d = 0; d < dim; d++) {
125*828beda2SMark Adams         PetscReal gridSpacing = (upper[d] - lower[d]) / faces[d];
126*828beda2SMark Adams         PetscInt  gridIdx     = (PetscInt)((coords[d] - lower[d]) / gridSpacing + 0.5);
127*828beda2SMark Adams         if (gridIdx % 2 != 0) {
128*828beda2SMark Adams           isObs = PETSC_FALSE;
129*828beda2SMark Adams           break;
130*828beda2SMark Adams         }
131*828beda2SMark Adams       }
132*828beda2SMark Adams 
133*828beda2SMark Adams       if (isObs) {
134*828beda2SMark Adams         PetscCall(MatSetValue(H, obsIdx, v - vStart, 1.0, INSERT_VALUES));
135*828beda2SMark Adams         obsIdx++;
136*828beda2SMark Adams       }
137*828beda2SMark Adams     }
138*828beda2SMark Adams     PetscCall(VecRestoreArrayRead(coordinates, &coordArray));
139*828beda2SMark Adams   }
140*828beda2SMark Adams   PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
141*828beda2SMark Adams   PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
142*828beda2SMark Adams 
143*828beda2SMark Adams   /* View H */
144*828beda2SMark Adams   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Observation Operator H:\n"));
145*828beda2SMark Adams   PetscCall(MatView(H, PETSC_VIEWER_STDOUT_WORLD));
146*828beda2SMark Adams 
147*828beda2SMark Adams   /* Perturb interior vertex coordinates */
148*828beda2SMark Adams   {
149*828beda2SMark Adams     Vec           coordinates;
150*828beda2SMark Adams     PetscSection  coordSection;
151*828beda2SMark Adams     PetscScalar  *coordArray;
152*828beda2SMark Adams     unsigned long seed = 123456789;
153*828beda2SMark Adams 
154*828beda2SMark Adams     PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
155*828beda2SMark Adams     PetscCall(DMGetCoordinateSection(dm, &coordSection));
156*828beda2SMark Adams     PetscCall(VecGetArray(coordinates, &coordArray));
157*828beda2SMark Adams 
158*828beda2SMark Adams     for (PetscInt v = vStart; v < vEnd; v++) {
159*828beda2SMark Adams       PetscInt offset;
160*828beda2SMark Adams 
161*828beda2SMark Adams       PetscCall(PetscSectionGetOffset(coordSection, v, &offset));
162*828beda2SMark Adams 
163*828beda2SMark Adams       for (PetscInt d = 0; d < dim; d++) {
164*828beda2SMark Adams         PetscReal noise, gridSpacing = (upper[d] - lower[d]) / faces[d];
165*828beda2SMark Adams 
166*828beda2SMark Adams         seed  = (1103515245 * seed + 12345) % 2147483648;
167*828beda2SMark Adams         noise = (PetscReal)seed / 2147483648.0;
168*828beda2SMark Adams         coordArray[offset + d] += (noise - 0.5) * 0.05 * gridSpacing;
169*828beda2SMark Adams       }
170*828beda2SMark Adams     }
171*828beda2SMark Adams     PetscCall(VecRestoreArray(coordinates, &coordArray));
172*828beda2SMark Adams   }
173*828beda2SMark Adams 
174*828beda2SMark Adams   /* Call the function */
175*828beda2SMark Adams   PetscCall(DMPlexGetLETKFLocalizationMatrix(dm, numobservations, numlocalobs, H, &Q));
176*828beda2SMark Adams   PetscCall(PetscObjectSetName((PetscObject)Q, "Q_localization"));
177*828beda2SMark Adams 
178*828beda2SMark Adams   /* View Q */
179*828beda2SMark Adams   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Localization Matrix Q:\n"));
180*828beda2SMark Adams   PetscCall(MatView(Q, PETSC_VIEWER_STDOUT_WORLD));
181*828beda2SMark Adams 
182*828beda2SMark Adams   /* Cleanup */
183*828beda2SMark Adams   PetscCall(MatDestroy(&H));
184*828beda2SMark Adams   PetscCall(MatDestroy(&Q));
185*828beda2SMark Adams   PetscCall(DMDestroy(&dm));
186*828beda2SMark Adams   PetscCall(PetscFinalize());
187*828beda2SMark Adams   return 0;
188*828beda2SMark Adams }
189*828beda2SMark Adams 
190*828beda2SMark Adams /*TEST
191*828beda2SMark Adams 
192*828beda2SMark Adams   test:
193*828beda2SMark Adams     requires: kokkos
194*828beda2SMark Adams     suffix: 1
195*828beda2SMark Adams     diff_args: -j
196*828beda2SMark Adams     args: -dm_plex_dim 1 -dm_plex_box_faces 16 -dm_plex_simplex 0
197*828beda2SMark Adams 
198*828beda2SMark Adams   test:
199*828beda2SMark Adams     requires: kokkos
200*828beda2SMark Adams     suffix: 2
201*828beda2SMark Adams     diff_args: -j
202*828beda2SMark Adams     args: -dm_plex_dim 2 -dm_plex_box_faces 8,8 -dm_plex_simplex 0
203*828beda2SMark Adams 
204*828beda2SMark Adams   test:
205*828beda2SMark Adams     requires: kokkos
206*828beda2SMark Adams     suffix: 3
207*828beda2SMark Adams     diff_args: -j
208*828beda2SMark Adams     args: -dm_plex_dim 3 -dm_plex_box_faces 5,5,5 -dm_plex_simplex 0
209*828beda2SMark Adams 
210*828beda2SMark Adams TEST*/
211