1828beda2SMark Adams static char help[] = "Test DMPlexGetLETKFLocalizationMatrix.\n\n";
2828beda2SMark Adams
3828beda2SMark Adams #include <petscdmplex.h>
4*ee102026SMark Adams #include <petscdmda.h>
5828beda2SMark Adams
main(int argc,char ** argv)6828beda2SMark Adams int main(int argc, char **argv)
7828beda2SMark Adams {
8828beda2SMark Adams DM dm;
9*ee102026SMark Adams Mat H, Q = NULL;
10*ee102026SMark Adams PetscInt nvertexobs, ndof = 1, n_state_global;
11*ee102026SMark Adams PetscInt dim = 1, n, vStart, vEnd;
12*ee102026SMark Adams PetscInt faces[3] = {1, 1, 1};
13828beda2SMark Adams PetscReal lower[3] = {0.0, 0.0, 0.0};
14828beda2SMark Adams PetscReal upper[3] = {1.0, 1.0, 1.0};
15*ee102026SMark Adams Vec Vecxyz[3] = {NULL, NULL, NULL};
16*ee102026SMark Adams PetscBool isda, isplex, print = PETSC_FALSE;
17*ee102026SMark Adams char type[256] = DMPLEX;
18828beda2SMark Adams
19828beda2SMark Adams PetscFunctionBeginUser;
20828beda2SMark Adams PetscCall(PetscInitialize(&argc, &argv, NULL, help));
21828beda2SMark Adams
22828beda2SMark Adams /* Get dimension and from options. We need the data here and Plex does not have access functions */
23*ee102026SMark Adams PetscOptionsBegin(PETSC_COMM_WORLD, "", "DMField Tutorial Options", "DM");
24*ee102026SMark Adams PetscCall(PetscOptionsFList("-dm_type", "DM implementation on which to define field", "ex20.c", DMList, type, type, 256, NULL));
25*ee102026SMark Adams PetscCall(PetscStrncmp(type, DMPLEX, 256, &isplex));
26*ee102026SMark Adams PetscCall(PetscStrncmp(type, DMDA, 256, &isda));
27*ee102026SMark Adams PetscCall(PetscOptionsGetBool(NULL, NULL, "-ex20_print", &print, NULL));
28828beda2SMark Adams PetscCall(PetscOptionsGetInt(NULL, NULL, "-dm_plex_dim", &dim, NULL));
29*ee102026SMark Adams PetscCheck(dim <= 3 && dim >= 1, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "dm_plex_dim = %" PetscInt_FMT, dim);
30828beda2SMark Adams n = 3;
31828beda2SMark Adams PetscCall(PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &n, NULL));
32*ee102026SMark Adams PetscCheck(n == 0 || n == dim, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "dm_plex_box_faces dimension %" PetscInt_FMT " does not match requested dimension %" PetscInt_FMT, n, dim);
33828beda2SMark Adams n = 3;
34828beda2SMark Adams PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-dm_plex_box_lower", lower, &n, NULL));
35*ee102026SMark Adams PetscCheck(n == 0 || n == dim, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "dm_plex_box_lower dimension %" PetscInt_FMT " does not match requested dimension %" PetscInt_FMT, n, dim);
36828beda2SMark Adams n = 3;
37828beda2SMark Adams PetscCall(PetscOptionsGetRealArray(NULL, NULL, "-dm_plex_box_upper", upper, &n, NULL));
38*ee102026SMark Adams PetscCheck(n == 0 || n == dim, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "dm_plex_box_upper dimension %" PetscInt_FMT " does not match requested dimension %" PetscInt_FMT, n, dim);
39*ee102026SMark Adams PetscOptionsEnd();
40828beda2SMark Adams
41*ee102026SMark Adams if (isplex) {
42828beda2SMark Adams PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
43828beda2SMark Adams PetscCall(DMSetType(dm, DMPLEX));
44828beda2SMark Adams PetscCall(DMSetFromOptions(dm));
45*ee102026SMark Adams {
46828beda2SMark Adams PetscSection section;
47828beda2SMark Adams PetscInt pStart, pEnd;
48*ee102026SMark Adams
49828beda2SMark Adams PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
50*ee102026SMark Adams PetscCall(DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd));
51828beda2SMark Adams PetscCall(PetscSectionCreate(PETSC_COMM_WORLD, §ion));
52828beda2SMark Adams PetscCall(PetscSectionSetNumFields(section, 1));
53828beda2SMark Adams PetscCall(PetscSectionSetChart(section, pStart, pEnd));
54828beda2SMark Adams for (PetscInt v = vStart; v < vEnd; ++v) PetscCall(PetscSectionSetDof(section, v, 1));
55828beda2SMark Adams PetscCall(PetscSectionSetUp(section));
56828beda2SMark Adams PetscCall(DMSetLocalSection(dm, section));
57828beda2SMark Adams PetscCall(PetscSectionDestroy(§ion));
58828beda2SMark Adams
59*ee102026SMark Adams for (PetscInt d = 0; d < dim; d++) {
60*ee102026SMark Adams Vec loc_vec;
61*ee102026SMark Adams Vec coordinates;
62*ee102026SMark Adams PetscSection coordSection, s;
63*ee102026SMark Adams const PetscScalar *coordArray;
64*ee102026SMark Adams PetscScalar *xArray;
65*ee102026SMark Adams
66*ee102026SMark Adams PetscCall(DMCreateGlobalVector(dm, &Vecxyz[d]));
67*ee102026SMark Adams PetscCall(PetscObjectSetName((PetscObject)Vecxyz[d], d == 0 ? "x_coordinate" : (d == 1 ? "y_coordinate" : "z_coordinate")));
68*ee102026SMark Adams PetscCall(DMGetLocalVector(dm, &loc_vec));
69*ee102026SMark Adams
70*ee102026SMark Adams PetscCall(DMGetLocalSection(dm, &s));
71*ee102026SMark Adams PetscCall(DMGetCoordinatesLocal(dm, &coordinates));
72*ee102026SMark Adams PetscCall(DMGetCoordinateSection(dm, &coordSection));
73*ee102026SMark Adams PetscCall(VecGetArrayRead(coordinates, &coordArray));
74*ee102026SMark Adams PetscCall(VecGetArray(loc_vec, &xArray));
75*ee102026SMark Adams
76*ee102026SMark Adams for (PetscInt v = vStart; v < vEnd; v++) {
77*ee102026SMark Adams PetscInt cOff, sOff;
78*ee102026SMark Adams
79*ee102026SMark Adams PetscCall(PetscSectionGetOffset(coordSection, v, &cOff));
80*ee102026SMark Adams PetscCall(PetscSectionGetOffset(s, v, &sOff));
81*ee102026SMark Adams xArray[sOff] = coordArray[cOff + d];
82*ee102026SMark Adams }
83*ee102026SMark Adams PetscCall(VecRestoreArrayRead(coordinates, &coordArray));
84*ee102026SMark Adams PetscCall(VecRestoreArray(loc_vec, &xArray));
85*ee102026SMark Adams
86*ee102026SMark Adams PetscCall(DMLocalToGlobal(dm, loc_vec, INSERT_VALUES, Vecxyz[d]));
87*ee102026SMark Adams PetscCall(DMRestoreLocalVector(dm, &loc_vec));
88*ee102026SMark Adams PetscCall(VecGetSize(Vecxyz[d], &n_state_global));
89*ee102026SMark Adams }
90*ee102026SMark Adams }
91*ee102026SMark Adams
92*ee102026SMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Created DMPlex in %" PetscInt_FMT "D with faces (%" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "), global vector size %" PetscInt_FMT "\n", dim, faces[0], faces[1], faces[2], n_state_global));
93*ee102026SMark Adams } else if (isda) {
94*ee102026SMark Adams switch (dim) {
95*ee102026SMark Adams case 1:
96*ee102026SMark Adams PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, faces[0], ndof, 1, NULL, &dm));
97*ee102026SMark Adams break;
98*ee102026SMark Adams case 2:
99*ee102026SMark Adams PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, faces[0], faces[1] + 1, PETSC_DETERMINE, PETSC_DETERMINE, ndof, 1, NULL, NULL, &dm));
100*ee102026SMark Adams break;
101*ee102026SMark Adams default:
102*ee102026SMark Adams PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_PERIODIC, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, faces[0], faces[1] + 1, faces[2] + 1, PETSC_DETERMINE, PETSC_DETERMINE, PETSC_DETERMINE, ndof, 1, NULL, NULL, NULL, &dm));
103*ee102026SMark Adams break;
104*ee102026SMark Adams }
105*ee102026SMark Adams PetscCall(DMSetUp(dm));
106*ee102026SMark Adams PetscCall(DMDASetUniformCoordinates(dm, lower[0], upper[0], lower[1], upper[1], lower[2], upper[2]));
107*ee102026SMark Adams {
108*ee102026SMark Adams Vec coord;
109*ee102026SMark Adams PetscCall(DMGetCoordinates(dm, &coord));
110*ee102026SMark Adams for (PetscInt d = 0; d < dim; d++) {
111*ee102026SMark Adams PetscCall(DMCreateGlobalVector(dm, &Vecxyz[d]));
112*ee102026SMark Adams PetscCall(VecSetFromOptions(Vecxyz[d]));
113*ee102026SMark Adams PetscCall(PetscObjectSetName((PetscObject)Vecxyz[d], d == 0 ? "x_coordinate" : (d == 1 ? "y_coordinate" : "z_coordinate")));
114*ee102026SMark Adams PetscCall(VecStrideGather(coord, d, Vecxyz[d], INSERT_VALUES));
115*ee102026SMark Adams PetscCall(VecGetSize(Vecxyz[d], &n_state_global));
116*ee102026SMark Adams }
117*ee102026SMark Adams }
118*ee102026SMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Created DMDA of type %s in %" PetscInt_FMT "D with faces (%" PetscInt_FMT ", %" PetscInt_FMT ", %" PetscInt_FMT "), global vector size %" PetscInt_FMT "\n", type, dim, faces[0], faces[1], faces[2], n_state_global));
119*ee102026SMark Adams } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This test does not run for DM type %s", type);
120*ee102026SMark Adams PetscCall(DMViewFromOptions(dm, NULL, "-ex20_dm_view")); // PetscSleep(10);
121*ee102026SMark Adams
122*ee102026SMark Adams /* Set number of local observations to use: 3^dim */
123*ee102026SMark Adams nvertexobs = 1;
124*ee102026SMark Adams for (PetscInt d = 0; d < dim && d < 2; d++) nvertexobs *= 3;
125*ee102026SMark Adams PetscCheck(nvertexobs > 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "nvertexobs %" PetscInt_FMT " must be > 0 locally for now", nvertexobs);
126828beda2SMark Adams
127828beda2SMark Adams /* Count observations (every other vertex in each dimension) */
128*ee102026SMark Adams PetscInt nobs_local = 0;
129*ee102026SMark Adams PetscBool *isObs;
130*ee102026SMark Adams PetscInt nloc;
131*ee102026SMark Adams
132*ee102026SMark Adams PetscCall(VecGetLocalSize(Vecxyz[0], &nloc));
133*ee102026SMark Adams PetscCall(PetscMalloc1(nloc, &isObs));
134828beda2SMark Adams {
135*ee102026SMark Adams const PetscScalar *coords[3];
136*ee102026SMark Adams PetscReal gridSpacing[3];
137*ee102026SMark Adams for (PetscInt d = 0; d < dim; d++) PetscCall(VecGetArrayRead(Vecxyz[d], &coords[d]));
138*ee102026SMark Adams for (PetscInt d = 0; d < dim; d++) gridSpacing[d] = (upper[d] - lower[d]) / faces[d];
139828beda2SMark Adams
140*ee102026SMark Adams for (PetscInt v = 0; v < nloc; v++) {
141*ee102026SMark Adams PetscReal c[3] = {0.0, 0.0, 0.0};
142828beda2SMark Adams
143*ee102026SMark Adams isObs[v] = PETSC_TRUE;
144*ee102026SMark Adams for (PetscInt d = 0; d < dim; d++) c[d] = PetscRealPart(coords[d][v]);
145828beda2SMark Adams /* Check if this vertex is at an observation location (every other grid point) */
146828beda2SMark Adams for (PetscInt d = 0; d < dim; d++) {
147*ee102026SMark Adams PetscReal relCoord = c[d] - lower[d];
148*ee102026SMark Adams PetscInt gridIdx = (PetscInt)PetscFloorReal(relCoord / gridSpacing[d] + 0.5);
149*ee102026SMark Adams PetscCheck(PetscAbsReal(relCoord - gridIdx * gridSpacing[d]) < PETSC_SMALL, PETSC_COMM_WORLD, PETSC_ERR_ARG_INCOMP, "Error vertex v %" PetscInt_FMT " (dim %" PetscInt_FMT "): %g not on grid (h= %g, distance to grid %g)", v, d, (double)c[d], (double)gridSpacing[d], (double)PetscAbsReal(relCoord - gridIdx * gridSpacing[d]));
150828beda2SMark Adams if (gridIdx % 2 != 0) {
151*ee102026SMark Adams isObs[v] = PETSC_FALSE;
152828beda2SMark Adams break;
153828beda2SMark Adams }
154828beda2SMark Adams }
155*ee102026SMark Adams if (isObs[v]) nobs_local++;
156828beda2SMark Adams }
157*ee102026SMark Adams for (PetscInt d = 0; d < dim; d++) PetscCall(VecRestoreArrayRead(Vecxyz[d], &coords[d]));
158828beda2SMark Adams }
159828beda2SMark Adams
160*ee102026SMark Adams /* Create H matrix n_obs X n_state */
161828beda2SMark Adams PetscCall(MatCreate(PETSC_COMM_WORLD, &H));
162*ee102026SMark Adams PetscCall(MatSetSizes(H, nobs_local, PETSC_DECIDE, PETSC_DECIDE, n_state_global)); //
163*ee102026SMark Adams PetscCall(MatSetBlockSizes(H, 1, ndof));
164828beda2SMark Adams PetscCall(MatSetType(H, MATAIJ));
165828beda2SMark Adams PetscCall(MatSeqAIJSetPreallocation(H, 1, NULL));
166*ee102026SMark Adams PetscCall(MatMPIAIJSetPreallocation(H, 1, NULL, 1, NULL)); // assumes boolean observation operator, could use interpolation
167828beda2SMark Adams PetscCall(PetscObjectSetName((PetscObject)H, "H_observation_operator"));
168*ee102026SMark Adams PetscCall(MatSetFromOptions(H));
169828beda2SMark Adams
170828beda2SMark Adams /* Fill H matrix */
171*ee102026SMark Adams PetscInt globalRowIdx, globalColIdx, obsIdx = 0;
172*ee102026SMark Adams PetscCall(VecGetOwnershipRange(Vecxyz[0], &globalColIdx, NULL));
173*ee102026SMark Adams PetscCall(MatGetOwnershipRange(H, &globalRowIdx, NULL));
174*ee102026SMark Adams for (PetscInt v = 0; v < nloc; v++) {
175*ee102026SMark Adams if (isObs[v]) {
176*ee102026SMark Adams PetscInt grow = globalRowIdx + obsIdx++, gcol = globalColIdx + v;
177*ee102026SMark Adams PetscCall(MatSetValue(H, grow, gcol, 1.0, INSERT_VALUES));
178828beda2SMark Adams }
179828beda2SMark Adams }
180*ee102026SMark Adams PetscCall(PetscFree(isObs));
181828beda2SMark Adams PetscCall(MatAssemblyBegin(H, MAT_FINAL_ASSEMBLY));
182828beda2SMark Adams PetscCall(MatAssemblyEnd(H, MAT_FINAL_ASSEMBLY));
183828beda2SMark Adams
184828beda2SMark Adams /* View H */
185828beda2SMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Observation Operator H:\n"));
186*ee102026SMark Adams if (print) PetscCall(MatView(H, PETSC_VIEWER_STDOUT_WORLD));
187828beda2SMark Adams
188828beda2SMark Adams /* Perturb interior vertex coordinates */
189828beda2SMark Adams {
190*ee102026SMark Adams PetscScalar *coords[3];
191*ee102026SMark Adams PetscInt nloc;
192828beda2SMark Adams unsigned long seed = 123456789;
193828beda2SMark Adams
194*ee102026SMark Adams PetscCall(VecGetLocalSize(Vecxyz[0], &nloc));
195*ee102026SMark Adams for (PetscInt d = 0; d < dim; d++) PetscCall(VecGetArray(Vecxyz[d], &coords[d]));
196828beda2SMark Adams
197*ee102026SMark Adams for (PetscInt v = 0; v < nloc; v++) {
198828beda2SMark Adams for (PetscInt d = 0; d < dim; d++) {
199828beda2SMark Adams PetscReal noise, gridSpacing = (upper[d] - lower[d]) / faces[d];
200828beda2SMark Adams
201828beda2SMark Adams seed = (1103515245 * seed + 12345) % 2147483648;
202828beda2SMark Adams noise = (PetscReal)seed / 2147483648.0;
203*ee102026SMark Adams coords[d][v] += (noise - 0.5) * 0.001 * gridSpacing;
204828beda2SMark Adams }
205828beda2SMark Adams }
206*ee102026SMark Adams for (PetscInt d = 0; d < dim; d++) PetscCall(VecRestoreArray(Vecxyz[d], &coords[d]));
207828beda2SMark Adams }
208828beda2SMark Adams
209828beda2SMark Adams /* Call the function */
210*ee102026SMark Adams PetscCall(DMPlexGetLETKFLocalizationMatrix(nvertexobs, nobs_local, ndof, Vecxyz, H, &Q));
211828beda2SMark Adams PetscCall(PetscObjectSetName((PetscObject)Q, "Q_localization"));
212828beda2SMark Adams
213*ee102026SMark Adams // View Q
214828beda2SMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Localization Matrix Q:\n"));
215*ee102026SMark Adams if (print) PetscCall(MatView(Q, PETSC_VIEWER_STDOUT_WORLD));
216828beda2SMark Adams
217828beda2SMark Adams /* Cleanup */
218*ee102026SMark Adams for (PetscInt d = 0; d < dim; d++) PetscCall(VecDestroy(&Vecxyz[d]));
219828beda2SMark Adams PetscCall(MatDestroy(&H));
220828beda2SMark Adams PetscCall(MatDestroy(&Q));
221828beda2SMark Adams PetscCall(DMDestroy(&dm));
222828beda2SMark Adams PetscCall(PetscFinalize());
223828beda2SMark Adams return 0;
224828beda2SMark Adams }
225828beda2SMark Adams
226828beda2SMark Adams /*TEST
227828beda2SMark Adams
228828beda2SMark Adams test:
229*ee102026SMark Adams requires: kokkos_kernels
230828beda2SMark Adams suffix: 1
231828beda2SMark Adams diff_args: -j
232*ee102026SMark Adams args: -dm_plex_dim 1 -dm_plex_box_faces 16 -dm_plex_simplex 0 -dm_plex_box_bd periodic -dm_plex_box_upper 5 -ex20_print -ex20_dm_view -ex20_dm_view -mat_type aijkokkos -dm_vec_type kokkos
233828beda2SMark Adams
234828beda2SMark Adams test:
235*ee102026SMark Adams requires: kokkos_kernels
236828beda2SMark Adams suffix: 2
237828beda2SMark Adams diff_args: -j
238*ee102026SMark Adams args: -dm_plex_dim 2 -dm_plex_box_faces 7,7 -dm_plex_simplex 0 -dm_plex_box_bd periodic,none -dm_plex_box_upper 5,5 -ex20_print -ex20_dm_view -ex20_dm_view -mat_type aijkokkos -dm_vec_type kokkos
239828beda2SMark Adams
240828beda2SMark Adams test:
241*ee102026SMark Adams requires: kokkos_kernels
242*ee102026SMark Adams suffix: da2
243*ee102026SMark Adams diff_args: -j
244*ee102026SMark Adams args: -dm_type da -dm_plex_dim 2 -dm_plex_box_faces 7,7 -dm_plex_box_upper 5,5 -ex20_print -ex20_dm_view -ex20_dm_view -mat_type aijkokkos -vec_type kokkos
245*ee102026SMark Adams
246*ee102026SMark Adams test:
247*ee102026SMark Adams requires: kokkos_kernels
248828beda2SMark Adams suffix: 3
249828beda2SMark Adams diff_args: -j
250*ee102026SMark Adams args: -dm_plex_dim 3 -dm_plex_box_faces 5,5,5 -dm_plex_simplex 0 -dm_plex_box_bd periodic,none,none -dm_plex_box_upper 5,5,5 -ex20_print -ex20_dm_view -mat_type aijkokkos -dm_vec_type kokkos
251828beda2SMark Adams
252828beda2SMark Adams TEST*/
253