xref: /petsc/src/dm/impls/plex/tests/ex64.c (revision 061e922f3926be00487707c73b78dd3d40309129)
1a8db3e61SBlaise Bourdin static char help[] = "Test FEM layout and GlobalToNaturalSF\n\n";
2a8db3e61SBlaise Bourdin 
3a8db3e61SBlaise Bourdin /*
4a8db3e61SBlaise Bourdin   In order to see the vectors which are being tested, use
5a8db3e61SBlaise Bourdin 
6a8db3e61SBlaise Bourdin      -ua_vec_view -s_vec_view
7a8db3e61SBlaise Bourdin */
8a8db3e61SBlaise Bourdin 
9a8db3e61SBlaise Bourdin #include <petsc.h>
10a8db3e61SBlaise Bourdin #include <exodusII.h>
11a8db3e61SBlaise Bourdin 
12a8db3e61SBlaise Bourdin #include <petsc/private/dmpleximpl.h>
13a8db3e61SBlaise Bourdin 
main(int argc,char ** argv)14*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
15*d71ae5a4SJacob Faibussowitsch {
16a8db3e61SBlaise Bourdin   DM              dm, pdm, dmU, dmA, dmS, dmUA, dmUA2, *dmList;
17a8db3e61SBlaise Bourdin   DM              seqdmU, seqdmA, seqdmS, seqdmUA, seqdmUA2, *seqdmList;
18a8db3e61SBlaise Bourdin   Vec             X, U, A, S, UA, UA2;
19a8db3e61SBlaise Bourdin   Vec             seqX, seqU, seqA, seqS, seqUA, seqUA2;
20a8db3e61SBlaise Bourdin   IS              isU, isA, isS, isUA;
21a8db3e61SBlaise Bourdin   IS              seqisU, seqisA, seqisS, seqisUA;
22a8db3e61SBlaise Bourdin   PetscSection    section;
23a8db3e61SBlaise Bourdin   const PetscInt  fieldU     = 0;
24a8db3e61SBlaise Bourdin   const PetscInt  fieldA     = 2;
25a8db3e61SBlaise Bourdin   const PetscInt  fieldS     = 1;
26a8db3e61SBlaise Bourdin   const PetscInt  fieldUA[2] = {0, 2};
27a8db3e61SBlaise Bourdin   char            ifilename[PETSC_MAX_PATH_LEN];
28a8db3e61SBlaise Bourdin   IS              csIS;
29a8db3e61SBlaise Bourdin   const PetscInt *csID;
30a8db3e61SBlaise Bourdin   PetscInt       *pStartDepth, *pEndDepth;
31a8db3e61SBlaise Bourdin   PetscInt        order = 1;
32a8db3e61SBlaise Bourdin   PetscInt        sdim, d, pStart, pEnd, p, numCS, set;
33a8db3e61SBlaise Bourdin   PetscMPIInt     rank, size;
34a8db3e61SBlaise Bourdin   PetscSF         migrationSF;
35a8db3e61SBlaise Bourdin 
36a8db3e61SBlaise Bourdin   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
37a8db3e61SBlaise Bourdin   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
38a8db3e61SBlaise Bourdin   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
39a8db3e61SBlaise Bourdin   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "FEM Layout Options", "ex64");
40a8db3e61SBlaise Bourdin   PetscCall(PetscOptionsString("-i", "Filename to read", "ex64", ifilename, ifilename, sizeof(ifilename), NULL));
41a8db3e61SBlaise Bourdin   PetscOptionsEnd();
42a8db3e61SBlaise Bourdin 
43a8db3e61SBlaise Bourdin   /* Read the mesh from a file in any supported format */
44a8db3e61SBlaise Bourdin   PetscCall(DMPlexCreateFromFile(PETSC_COMM_WORLD, ifilename, "ex64_plex", PETSC_TRUE, &dm));
45a8db3e61SBlaise Bourdin   PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE));
46a8db3e61SBlaise Bourdin   PetscCall(DMSetFromOptions(dm));
47a8db3e61SBlaise Bourdin   PetscCall(DMViewFromOptions(dm, NULL, "-dm_view"));
48a8db3e61SBlaise Bourdin   PetscCall(DMGetDimension(dm, &sdim));
49a8db3e61SBlaise Bourdin 
50a8db3e61SBlaise Bourdin   /* Create the main section containing all fields */
51a8db3e61SBlaise Bourdin   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dm), &section));
52a8db3e61SBlaise Bourdin   PetscCall(PetscSectionSetNumFields(section, 3));
53a8db3e61SBlaise Bourdin   PetscCall(PetscSectionSetFieldName(section, fieldU, "U"));
54a8db3e61SBlaise Bourdin   PetscCall(PetscSectionSetFieldName(section, fieldA, "Alpha"));
55a8db3e61SBlaise Bourdin   PetscCall(PetscSectionSetFieldName(section, fieldS, "Sigma"));
56a8db3e61SBlaise Bourdin   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
57a8db3e61SBlaise Bourdin   PetscCall(PetscSectionSetChart(section, pStart, pEnd));
58a8db3e61SBlaise Bourdin   PetscCall(PetscMalloc2(sdim + 1, &pStartDepth, sdim + 1, &pEndDepth));
5948a46eb9SPierre Jolivet   for (d = 0; d <= sdim; ++d) PetscCall(DMPlexGetDepthStratum(dm, d, &pStartDepth[d], &pEndDepth[d]));
60a8db3e61SBlaise Bourdin   /* Vector field U, Scalar field Alpha, Tensor field Sigma */
61a8db3e61SBlaise Bourdin   PetscCall(PetscSectionSetFieldComponents(section, fieldU, sdim));
62a8db3e61SBlaise Bourdin   PetscCall(PetscSectionSetFieldComponents(section, fieldA, 1));
63a8db3e61SBlaise Bourdin   PetscCall(PetscSectionSetFieldComponents(section, fieldS, sdim * (sdim + 1) / 2));
64a8db3e61SBlaise Bourdin 
65a8db3e61SBlaise Bourdin   /* Going through cell sets then cells, and setting up storage for the sections */
66a8db3e61SBlaise Bourdin   PetscCall(DMGetLabelSize(dm, "Cell Sets", &numCS));
67a8db3e61SBlaise Bourdin   PetscCall(DMGetLabelIdIS(dm, "Cell Sets", &csIS));
6848a46eb9SPierre Jolivet   if (csIS) PetscCall(ISGetIndices(csIS, &csID));
69a8db3e61SBlaise Bourdin   for (set = 0; set < numCS; set++) {
70a8db3e61SBlaise Bourdin     IS              cellIS;
71a8db3e61SBlaise Bourdin     const PetscInt *cellID;
72a8db3e61SBlaise Bourdin     PetscInt        numCells, cell, closureSize, *closureA = NULL;
73a8db3e61SBlaise Bourdin 
74a8db3e61SBlaise Bourdin     PetscCall(DMGetStratumSize(dm, "Cell Sets", csID[set], &numCells));
75a8db3e61SBlaise Bourdin     PetscCall(DMGetStratumIS(dm, "Cell Sets", csID[set], &cellIS));
76a8db3e61SBlaise Bourdin     if (numCells > 0) {
77a8db3e61SBlaise Bourdin       /* dof layout ordered by increasing height in the DAG: cell, face, edge, vertex */
78a8db3e61SBlaise Bourdin       PetscInt  dofUP1Tri[]  = {2, 0, 0};
79a8db3e61SBlaise Bourdin       PetscInt  dofAP1Tri[]  = {1, 0, 0};
80a8db3e61SBlaise Bourdin       PetscInt  dofUP2Tri[]  = {2, 2, 0};
81a8db3e61SBlaise Bourdin       PetscInt  dofAP2Tri[]  = {1, 1, 0};
82a8db3e61SBlaise Bourdin       PetscInt  dofUP1Quad[] = {2, 0, 0};
83a8db3e61SBlaise Bourdin       PetscInt  dofAP1Quad[] = {1, 0, 0};
84a8db3e61SBlaise Bourdin       PetscInt  dofUP2Quad[] = {2, 2, 2};
85a8db3e61SBlaise Bourdin       PetscInt  dofAP2Quad[] = {1, 1, 1};
86a8db3e61SBlaise Bourdin       PetscInt  dofS2D[]     = {0, 0, 3};
87a8db3e61SBlaise Bourdin       PetscInt  dofUP1Tet[]  = {3, 0, 0, 0};
88a8db3e61SBlaise Bourdin       PetscInt  dofAP1Tet[]  = {1, 0, 0, 0};
89a8db3e61SBlaise Bourdin       PetscInt  dofUP2Tet[]  = {3, 3, 0, 0};
90a8db3e61SBlaise Bourdin       PetscInt  dofAP2Tet[]  = {1, 1, 0, 0};
91a8db3e61SBlaise Bourdin       PetscInt  dofUP1Hex[]  = {3, 0, 0, 0};
92a8db3e61SBlaise Bourdin       PetscInt  dofAP1Hex[]  = {1, 0, 0, 0};
93a8db3e61SBlaise Bourdin       PetscInt  dofUP2Hex[]  = {3, 3, 3, 3};
94a8db3e61SBlaise Bourdin       PetscInt  dofAP2Hex[]  = {1, 1, 1, 1};
95a8db3e61SBlaise Bourdin       PetscInt  dofS3D[]     = {0, 0, 0, 6};
96a8db3e61SBlaise Bourdin       PetscInt *dofU, *dofA, *dofS;
97a8db3e61SBlaise Bourdin 
98a8db3e61SBlaise Bourdin       switch (sdim) {
99*d71ae5a4SJacob Faibussowitsch       case 2:
100*d71ae5a4SJacob Faibussowitsch         dofS = dofS2D;
101*d71ae5a4SJacob Faibussowitsch         break;
102*d71ae5a4SJacob Faibussowitsch       case 3:
103*d71ae5a4SJacob Faibussowitsch         dofS = dofS3D;
104*d71ae5a4SJacob Faibussowitsch         break;
105*d71ae5a4SJacob Faibussowitsch       default:
106*d71ae5a4SJacob Faibussowitsch         SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "No layout for dimension %" PetscInt_FMT, sdim);
107a8db3e61SBlaise Bourdin       }
108a8db3e61SBlaise Bourdin 
109a8db3e61SBlaise Bourdin       /* Identify cell type based on closure size only. This works for Tri/Tet/Quad/Hex meshes
110a8db3e61SBlaise Bourdin          It will not be enough to identify more exotic elements like pyramid or prisms...  */
111a8db3e61SBlaise Bourdin       PetscCall(ISGetIndices(cellIS, &cellID));
112a8db3e61SBlaise Bourdin       PetscCall(DMPlexGetTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA));
113a8db3e61SBlaise Bourdin       switch (closureSize) {
114a8db3e61SBlaise Bourdin       case 7: /* Tri */
115a8db3e61SBlaise Bourdin         if (order == 1) {
116a8db3e61SBlaise Bourdin           dofU = dofUP1Tri;
117a8db3e61SBlaise Bourdin           dofA = dofAP1Tri;
118a8db3e61SBlaise Bourdin         } else {
119a8db3e61SBlaise Bourdin           dofU = dofUP2Tri;
120a8db3e61SBlaise Bourdin           dofA = dofAP2Tri;
121a8db3e61SBlaise Bourdin         }
122a8db3e61SBlaise Bourdin         break;
123a8db3e61SBlaise Bourdin       case 9: /* Quad */
124a8db3e61SBlaise Bourdin         if (order == 1) {
125a8db3e61SBlaise Bourdin           dofU = dofUP1Quad;
126a8db3e61SBlaise Bourdin           dofA = dofAP1Quad;
127a8db3e61SBlaise Bourdin         } else {
128a8db3e61SBlaise Bourdin           dofU = dofUP2Quad;
129a8db3e61SBlaise Bourdin           dofA = dofAP2Quad;
130a8db3e61SBlaise Bourdin         }
131a8db3e61SBlaise Bourdin         break;
132a8db3e61SBlaise Bourdin       case 15: /* Tet */
133a8db3e61SBlaise Bourdin         if (order == 1) {
134a8db3e61SBlaise Bourdin           dofU = dofUP1Tet;
135a8db3e61SBlaise Bourdin           dofA = dofAP1Tet;
136a8db3e61SBlaise Bourdin         } else {
137a8db3e61SBlaise Bourdin           dofU = dofUP2Tet;
138a8db3e61SBlaise Bourdin           dofA = dofAP2Tet;
139a8db3e61SBlaise Bourdin         }
140a8db3e61SBlaise Bourdin         break;
141a8db3e61SBlaise Bourdin       case 27: /* Hex */
142a8db3e61SBlaise Bourdin         if (order == 1) {
143a8db3e61SBlaise Bourdin           dofU = dofUP1Hex;
144a8db3e61SBlaise Bourdin           dofA = dofAP1Hex;
145a8db3e61SBlaise Bourdin         } else {
146a8db3e61SBlaise Bourdin           dofU = dofUP2Hex;
147a8db3e61SBlaise Bourdin           dofA = dofAP2Hex;
148a8db3e61SBlaise Bourdin         }
149a8db3e61SBlaise Bourdin         break;
150*d71ae5a4SJacob Faibussowitsch       default:
151*d71ae5a4SJacob Faibussowitsch         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Unknown element with closure size %" PetscInt_FMT, closureSize);
152a8db3e61SBlaise Bourdin       }
153a8db3e61SBlaise Bourdin       PetscCall(DMPlexRestoreTransitiveClosure(dm, cellID[0], PETSC_TRUE, &closureSize, &closureA));
154a8db3e61SBlaise Bourdin 
155a8db3e61SBlaise Bourdin       for (cell = 0; cell < numCells; cell++) {
156a8db3e61SBlaise Bourdin         PetscInt *closure = NULL;
157a8db3e61SBlaise Bourdin 
158a8db3e61SBlaise Bourdin         PetscCall(DMPlexGetTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure));
159a8db3e61SBlaise Bourdin         for (p = 0; p < closureSize; ++p) {
160a8db3e61SBlaise Bourdin           /* Find depth of p */
161a8db3e61SBlaise Bourdin           for (d = 0; d <= sdim; ++d) {
162a8db3e61SBlaise Bourdin             if ((closure[2 * p] >= pStartDepth[d]) && (closure[2 * p] < pEndDepth[d])) {
163a8db3e61SBlaise Bourdin               PetscCall(PetscSectionSetDof(section, closure[2 * p], dofU[d] + dofA[d] + dofS[d]));
164a8db3e61SBlaise Bourdin               PetscCall(PetscSectionSetFieldDof(section, closure[2 * p], fieldU, dofU[d]));
165a8db3e61SBlaise Bourdin               PetscCall(PetscSectionSetFieldDof(section, closure[2 * p], fieldA, dofA[d]));
166a8db3e61SBlaise Bourdin               PetscCall(PetscSectionSetFieldDof(section, closure[2 * p], fieldS, dofS[d]));
167a8db3e61SBlaise Bourdin             }
168a8db3e61SBlaise Bourdin           }
169a8db3e61SBlaise Bourdin         }
170a8db3e61SBlaise Bourdin         PetscCall(DMPlexRestoreTransitiveClosure(dm, cellID[cell], PETSC_TRUE, &closureSize, &closure));
171a8db3e61SBlaise Bourdin       }
172a8db3e61SBlaise Bourdin       PetscCall(ISRestoreIndices(cellIS, &cellID));
173a8db3e61SBlaise Bourdin       PetscCall(ISDestroy(&cellIS));
174a8db3e61SBlaise Bourdin     }
175a8db3e61SBlaise Bourdin   }
17648a46eb9SPierre Jolivet   if (csIS) PetscCall(ISRestoreIndices(csIS, &csID));
177a8db3e61SBlaise Bourdin   PetscCall(ISDestroy(&csIS));
178a8db3e61SBlaise Bourdin   PetscCall(PetscSectionSetUp(section));
179a8db3e61SBlaise Bourdin   PetscCall(DMSetLocalSection(dm, section));
180a8db3e61SBlaise Bourdin   PetscCall(PetscObjectViewFromOptions((PetscObject)section, NULL, "-dm_section_view"));
181a8db3e61SBlaise Bourdin   PetscCall(PetscSectionDestroy(&section));
182a8db3e61SBlaise Bourdin 
183a8db3e61SBlaise Bourdin   /* Get DM and IS for each field of dm */
184a8db3e61SBlaise Bourdin   PetscCall(DMCreateSubDM(dm, 1, &fieldU, &seqisU, &seqdmU));
185a8db3e61SBlaise Bourdin   PetscCall(DMCreateSubDM(dm, 1, &fieldA, &seqisA, &seqdmA));
186a8db3e61SBlaise Bourdin   PetscCall(DMCreateSubDM(dm, 1, &fieldS, &seqisS, &seqdmS));
187a8db3e61SBlaise Bourdin   PetscCall(DMCreateSubDM(dm, 2, fieldUA, &seqisUA, &seqdmUA));
188a8db3e61SBlaise Bourdin 
189a8db3e61SBlaise Bourdin   PetscCall(PetscMalloc1(2, &seqdmList));
190a8db3e61SBlaise Bourdin   seqdmList[0] = seqdmU;
191a8db3e61SBlaise Bourdin   seqdmList[1] = seqdmA;
192a8db3e61SBlaise Bourdin 
193a8db3e61SBlaise Bourdin   PetscCall(DMCreateSuperDM(seqdmList, 2, NULL, &seqdmUA2));
194a8db3e61SBlaise Bourdin   PetscCall(PetscFree(seqdmList));
195a8db3e61SBlaise Bourdin 
196a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(dm, &seqX));
197a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(seqdmU, &seqU));
198a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(seqdmA, &seqA));
199a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(seqdmS, &seqS));
200a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(seqdmUA, &seqUA));
201a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(seqdmUA2, &seqUA2));
202a8db3e61SBlaise Bourdin 
203a8db3e61SBlaise Bourdin   PetscCall(PetscObjectSetName((PetscObject)seqX, "seqX"));
204a8db3e61SBlaise Bourdin   PetscCall(PetscObjectSetName((PetscObject)seqU, "seqU"));
205a8db3e61SBlaise Bourdin   PetscCall(PetscObjectSetName((PetscObject)seqA, "seqAlpha"));
206a8db3e61SBlaise Bourdin   PetscCall(PetscObjectSetName((PetscObject)seqS, "seqSigma"));
207a8db3e61SBlaise Bourdin   PetscCall(PetscObjectSetName((PetscObject)seqUA, "seqUAlpha"));
208a8db3e61SBlaise Bourdin   PetscCall(PetscObjectSetName((PetscObject)seqUA2, "seqUAlpha2"));
209a8db3e61SBlaise Bourdin   PetscCall(VecSet(seqX, -111.));
210a8db3e61SBlaise Bourdin 
211a8db3e61SBlaise Bourdin   /* Setting u to [x,y,z]  and alpha to x^2+y^2+z^2 by writing in UAlpha then restricting to U and Alpha */
212a8db3e61SBlaise Bourdin   {
213a8db3e61SBlaise Bourdin     PetscSection sectionUA;
214a8db3e61SBlaise Bourdin     Vec          UALoc;
215a8db3e61SBlaise Bourdin     PetscSection coordSection;
216a8db3e61SBlaise Bourdin     Vec          coord;
217a8db3e61SBlaise Bourdin     PetscScalar *cval, *xyz;
218a8db3e61SBlaise Bourdin     PetscInt     clSize, i, j;
219a8db3e61SBlaise Bourdin 
220a8db3e61SBlaise Bourdin     PetscCall(DMGetLocalSection(seqdmUA, &sectionUA));
221a8db3e61SBlaise Bourdin     PetscCall(DMGetLocalVector(seqdmUA, &UALoc));
222a8db3e61SBlaise Bourdin     PetscCall(VecGetArray(UALoc, &cval));
223a8db3e61SBlaise Bourdin     PetscCall(DMGetCoordinateSection(seqdmUA, &coordSection));
224a8db3e61SBlaise Bourdin     PetscCall(DMGetCoordinatesLocal(seqdmUA, &coord));
225a8db3e61SBlaise Bourdin     PetscCall(DMPlexGetChart(seqdmUA, &pStart, &pEnd));
226a8db3e61SBlaise Bourdin 
227a8db3e61SBlaise Bourdin     for (p = pStart; p < pEnd; ++p) {
228a8db3e61SBlaise Bourdin       PetscInt dofUA, offUA;
229a8db3e61SBlaise Bourdin 
230a8db3e61SBlaise Bourdin       PetscCall(PetscSectionGetDof(sectionUA, p, &dofUA));
231a8db3e61SBlaise Bourdin       if (dofUA > 0) {
232a8db3e61SBlaise Bourdin         xyz = NULL;
233a8db3e61SBlaise Bourdin         PetscCall(PetscSectionGetOffset(sectionUA, p, &offUA));
234a8db3e61SBlaise Bourdin         PetscCall(DMPlexVecGetClosure(seqdmUA, coordSection, coord, p, &clSize, &xyz));
235a8db3e61SBlaise Bourdin         cval[offUA + sdim] = 0.;
236a8db3e61SBlaise Bourdin         for (i = 0; i < sdim; ++i) {
237a8db3e61SBlaise Bourdin           cval[offUA + i] = 0;
238ad540459SPierre Jolivet           for (j = 0; j < clSize / sdim; ++j) cval[offUA + i] += xyz[j * sdim + i];
239a8db3e61SBlaise Bourdin           cval[offUA + i] = cval[offUA + i] * sdim / clSize;
240a8db3e61SBlaise Bourdin           cval[offUA + sdim] += PetscSqr(cval[offUA + i]);
241a8db3e61SBlaise Bourdin         }
242a8db3e61SBlaise Bourdin         PetscCall(DMPlexVecRestoreClosure(seqdmUA, coordSection, coord, p, &clSize, &xyz));
243a8db3e61SBlaise Bourdin       }
244a8db3e61SBlaise Bourdin     }
245a8db3e61SBlaise Bourdin     PetscCall(VecRestoreArray(UALoc, &cval));
246a8db3e61SBlaise Bourdin     PetscCall(DMLocalToGlobalBegin(seqdmUA, UALoc, INSERT_VALUES, seqUA));
247a8db3e61SBlaise Bourdin     PetscCall(DMLocalToGlobalEnd(seqdmUA, UALoc, INSERT_VALUES, seqUA));
248a8db3e61SBlaise Bourdin     PetscCall(DMRestoreLocalVector(seqdmUA, &UALoc));
249a8db3e61SBlaise Bourdin 
250a8db3e61SBlaise Bourdin     /* Update X */
251a8db3e61SBlaise Bourdin     PetscCall(VecISCopy(seqX, seqisUA, SCATTER_FORWARD, seqUA));
252a8db3e61SBlaise Bourdin     PetscCall(VecViewFromOptions(seqUA, NULL, "-ua_vec_view"));
253a8db3e61SBlaise Bourdin 
254a8db3e61SBlaise Bourdin     /* Restrict to U and Alpha */
255a8db3e61SBlaise Bourdin     PetscCall(VecISCopy(seqX, seqisU, SCATTER_REVERSE, seqU));
256a8db3e61SBlaise Bourdin     PetscCall(VecISCopy(seqX, seqisA, SCATTER_REVERSE, seqA));
257a8db3e61SBlaise Bourdin 
258a8db3e61SBlaise Bourdin     /* restrict to UA2 */
259a8db3e61SBlaise Bourdin     PetscCall(VecISCopy(seqX, seqisUA, SCATTER_REVERSE, seqUA2));
260a8db3e61SBlaise Bourdin     PetscCall(VecViewFromOptions(seqUA2, NULL, "-ua2_vec_view"));
261a8db3e61SBlaise Bourdin   }
262a8db3e61SBlaise Bourdin 
263a8db3e61SBlaise Bourdin   {
264a8db3e61SBlaise Bourdin     PetscInt         ovlp = 0;
265a8db3e61SBlaise Bourdin     PetscPartitioner part;
266a8db3e61SBlaise Bourdin 
267a8db3e61SBlaise Bourdin     PetscCall(DMSetUseNatural(dm, PETSC_TRUE));
268a8db3e61SBlaise Bourdin     PetscCall(DMPlexGetPartitioner(dm, &part));
269a8db3e61SBlaise Bourdin     PetscCall(PetscPartitionerSetFromOptions(part));
270a8db3e61SBlaise Bourdin     PetscCall(DMPlexDistribute(dm, ovlp, &migrationSF, &pdm));
271a8db3e61SBlaise Bourdin     if (!pdm) pdm = dm;
272a8db3e61SBlaise Bourdin     if (migrationSF) {
273a8db3e61SBlaise Bourdin       PetscCall(DMPlexSetMigrationSF(pdm, migrationSF));
274a8db3e61SBlaise Bourdin       PetscCall(PetscSFDestroy(&migrationSF));
275a8db3e61SBlaise Bourdin     }
276a8db3e61SBlaise Bourdin     PetscCall(DMViewFromOptions(pdm, NULL, "-dm_view"));
277a8db3e61SBlaise Bourdin   }
278a8db3e61SBlaise Bourdin 
279a8db3e61SBlaise Bourdin   /* Get DM and IS for each field of dm */
280a8db3e61SBlaise Bourdin   PetscCall(DMCreateSubDM(pdm, 1, &fieldU, &isU, &dmU));
281a8db3e61SBlaise Bourdin   PetscCall(DMCreateSubDM(pdm, 1, &fieldA, &isA, &dmA));
282a8db3e61SBlaise Bourdin   PetscCall(DMCreateSubDM(pdm, 1, &fieldS, &isS, &dmS));
283a8db3e61SBlaise Bourdin   PetscCall(DMCreateSubDM(pdm, 2, fieldUA, &isUA, &dmUA));
284a8db3e61SBlaise Bourdin 
285a8db3e61SBlaise Bourdin   PetscCall(PetscMalloc1(2, &dmList));
286a8db3e61SBlaise Bourdin   dmList[0] = dmU;
287a8db3e61SBlaise Bourdin   dmList[1] = dmA;
288a8db3e61SBlaise Bourdin 
289a8db3e61SBlaise Bourdin   PetscCall(DMCreateSuperDM(dmList, 2, NULL, &dmUA2));
290a8db3e61SBlaise Bourdin   PetscCall(PetscFree(dmList));
291a8db3e61SBlaise Bourdin 
292a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(pdm, &X));
293a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(dmU, &U));
294a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(dmA, &A));
295a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(dmS, &S));
296a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(dmUA, &UA));
297a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(dmUA2, &UA2));
298a8db3e61SBlaise Bourdin 
299a8db3e61SBlaise Bourdin   PetscCall(PetscObjectSetName((PetscObject)X, "X"));
300a8db3e61SBlaise Bourdin   PetscCall(PetscObjectSetName((PetscObject)U, "U"));
301a8db3e61SBlaise Bourdin   PetscCall(PetscObjectSetName((PetscObject)A, "Alpha"));
302a8db3e61SBlaise Bourdin   PetscCall(PetscObjectSetName((PetscObject)S, "Sigma"));
303a8db3e61SBlaise Bourdin   PetscCall(PetscObjectSetName((PetscObject)UA, "UAlpha"));
304a8db3e61SBlaise Bourdin   PetscCall(PetscObjectSetName((PetscObject)UA2, "UAlpha2"));
305a8db3e61SBlaise Bourdin   PetscCall(VecSet(X, -111.));
306a8db3e61SBlaise Bourdin 
307a8db3e61SBlaise Bourdin   /* Setting u to [x,y,z]  and alpha to x^2+y^2+z^2 by writing in UAlpha then restricting to U and Alpha */
308a8db3e61SBlaise Bourdin   {
309a8db3e61SBlaise Bourdin     PetscSection sectionUA;
310a8db3e61SBlaise Bourdin     Vec          UALoc;
311a8db3e61SBlaise Bourdin     PetscSection coordSection;
312a8db3e61SBlaise Bourdin     Vec          coord;
313a8db3e61SBlaise Bourdin     PetscScalar *cval, *xyz;
314a8db3e61SBlaise Bourdin     PetscInt     clSize, i, j;
315a8db3e61SBlaise Bourdin 
316a8db3e61SBlaise Bourdin     PetscCall(DMGetLocalSection(dmUA, &sectionUA));
317a8db3e61SBlaise Bourdin     PetscCall(DMGetLocalVector(dmUA, &UALoc));
318a8db3e61SBlaise Bourdin     PetscCall(VecGetArray(UALoc, &cval));
319a8db3e61SBlaise Bourdin     PetscCall(DMGetCoordinateSection(dmUA, &coordSection));
320a8db3e61SBlaise Bourdin     PetscCall(DMGetCoordinatesLocal(dmUA, &coord));
321a8db3e61SBlaise Bourdin     PetscCall(DMPlexGetChart(dmUA, &pStart, &pEnd));
322a8db3e61SBlaise Bourdin 
323a8db3e61SBlaise Bourdin     for (p = pStart; p < pEnd; ++p) {
324a8db3e61SBlaise Bourdin       PetscInt dofUA, offUA;
325a8db3e61SBlaise Bourdin 
326a8db3e61SBlaise Bourdin       PetscCall(PetscSectionGetDof(sectionUA, p, &dofUA));
327a8db3e61SBlaise Bourdin       if (dofUA > 0) {
328a8db3e61SBlaise Bourdin         xyz = NULL;
329a8db3e61SBlaise Bourdin         PetscCall(PetscSectionGetOffset(sectionUA, p, &offUA));
330a8db3e61SBlaise Bourdin         PetscCall(DMPlexVecGetClosure(dmUA, coordSection, coord, p, &clSize, &xyz));
331a8db3e61SBlaise Bourdin         cval[offUA + sdim] = 0.;
332a8db3e61SBlaise Bourdin         for (i = 0; i < sdim; ++i) {
333a8db3e61SBlaise Bourdin           cval[offUA + i] = 0;
334ad540459SPierre Jolivet           for (j = 0; j < clSize / sdim; ++j) cval[offUA + i] += xyz[j * sdim + i];
335a8db3e61SBlaise Bourdin           cval[offUA + i] = cval[offUA + i] * sdim / clSize;
336a8db3e61SBlaise Bourdin           cval[offUA + sdim] += PetscSqr(cval[offUA + i]);
337a8db3e61SBlaise Bourdin         }
338a8db3e61SBlaise Bourdin         PetscCall(DMPlexVecRestoreClosure(dmUA, coordSection, coord, p, &clSize, &xyz));
339a8db3e61SBlaise Bourdin       }
340a8db3e61SBlaise Bourdin     }
341a8db3e61SBlaise Bourdin     PetscCall(VecRestoreArray(UALoc, &cval));
342a8db3e61SBlaise Bourdin     PetscCall(DMLocalToGlobalBegin(dmUA, UALoc, INSERT_VALUES, UA));
343a8db3e61SBlaise Bourdin     PetscCall(DMLocalToGlobalEnd(dmUA, UALoc, INSERT_VALUES, UA));
344a8db3e61SBlaise Bourdin     PetscCall(DMRestoreLocalVector(dmUA, &UALoc));
345a8db3e61SBlaise Bourdin 
346a8db3e61SBlaise Bourdin     /* Update X */
347a8db3e61SBlaise Bourdin     PetscCall(VecISCopy(X, isUA, SCATTER_FORWARD, UA));
348a8db3e61SBlaise Bourdin     PetscCall(VecViewFromOptions(UA, NULL, "-ua_vec_view"));
349a8db3e61SBlaise Bourdin 
350a8db3e61SBlaise Bourdin     /* Restrict to U and Alpha */
351a8db3e61SBlaise Bourdin     PetscCall(VecISCopy(X, isU, SCATTER_REVERSE, U));
352a8db3e61SBlaise Bourdin     PetscCall(VecISCopy(X, isA, SCATTER_REVERSE, A));
353a8db3e61SBlaise Bourdin 
354a8db3e61SBlaise Bourdin     /* restrict to UA2 */
355a8db3e61SBlaise Bourdin     PetscCall(VecISCopy(X, isUA, SCATTER_REVERSE, UA2));
356a8db3e61SBlaise Bourdin     PetscCall(VecViewFromOptions(UA2, NULL, "-ua2_vec_view"));
357a8db3e61SBlaise Bourdin   }
358a8db3e61SBlaise Bourdin 
359a8db3e61SBlaise Bourdin   /* Verification */
360a8db3e61SBlaise Bourdin 
361a8db3e61SBlaise Bourdin   Vec       Xnat, Unat, Anat, UAnat, Snat, UA2nat;
362a8db3e61SBlaise Bourdin   PetscReal norm;
363a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(dm, &Xnat));
364a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(seqdmU, &Unat));
365a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(seqdmA, &Anat));
366a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(seqdmUA, &UAnat));
367a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(seqdmS, &Snat));
368a8db3e61SBlaise Bourdin   PetscCall(DMGetGlobalVector(seqdmUA2, &UA2nat));
369a8db3e61SBlaise Bourdin 
370a8db3e61SBlaise Bourdin   PetscCall(DMPlexGlobalToNaturalBegin(pdm, X, Xnat));
371a8db3e61SBlaise Bourdin   PetscCall(DMPlexGlobalToNaturalEnd(pdm, X, Xnat));
372a8db3e61SBlaise Bourdin   PetscCall(VecAXPY(Xnat, -1.0, seqX));
373a8db3e61SBlaise Bourdin   PetscCall(VecNorm(Xnat, NORM_INFINITY, &norm));
374a8db3e61SBlaise Bourdin   PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "X ||Vin - Vout|| = %g", (double)norm);
375a8db3e61SBlaise Bourdin 
376a8db3e61SBlaise Bourdin   PetscCall(DMPlexGlobalToNaturalBegin(dmU, U, Unat));
377a8db3e61SBlaise Bourdin   PetscCall(DMPlexGlobalToNaturalEnd(dmU, U, Unat));
378a8db3e61SBlaise Bourdin   PetscCall(VecAXPY(Unat, -1.0, seqU));
379a8db3e61SBlaise Bourdin   PetscCall(VecNorm(Unat, NORM_INFINITY, &norm));
380a8db3e61SBlaise Bourdin   PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "U ||Vin - Vout|| = %g", (double)norm);
381a8db3e61SBlaise Bourdin 
382a8db3e61SBlaise Bourdin   PetscCall(DMPlexGlobalToNaturalBegin(dmA, A, Anat));
383a8db3e61SBlaise Bourdin   PetscCall(DMPlexGlobalToNaturalEnd(dmA, A, Anat));
384a8db3e61SBlaise Bourdin   PetscCall(VecAXPY(Anat, -1.0, seqA));
385a8db3e61SBlaise Bourdin   PetscCall(VecNorm(Anat, NORM_INFINITY, &norm));
386a8db3e61SBlaise Bourdin   PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "A ||Vin - Vout|| = %g", (double)norm);
387a8db3e61SBlaise Bourdin 
388a8db3e61SBlaise Bourdin   PetscCall(DMPlexGlobalToNaturalBegin(dmUA, UA, UAnat));
389a8db3e61SBlaise Bourdin   PetscCall(DMPlexGlobalToNaturalEnd(dmUA, UA, UAnat));
390a8db3e61SBlaise Bourdin   PetscCall(VecAXPY(UAnat, -1.0, seqUA));
391a8db3e61SBlaise Bourdin   PetscCall(VecNorm(UAnat, NORM_INFINITY, &norm));
392a8db3e61SBlaise Bourdin   PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "UA ||Vin - Vout|| = %g", (double)norm);
393a8db3e61SBlaise Bourdin 
394a8db3e61SBlaise Bourdin   PetscCall(DMPlexGlobalToNaturalBegin(dmS, S, Snat));
395a8db3e61SBlaise Bourdin   PetscCall(DMPlexGlobalToNaturalEnd(dmS, S, Snat));
396a8db3e61SBlaise Bourdin   PetscCall(VecAXPY(Snat, -1.0, seqS));
397a8db3e61SBlaise Bourdin   PetscCall(VecNorm(Snat, NORM_INFINITY, &norm));
398a8db3e61SBlaise Bourdin   PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "S ||Vin - Vout|| = %g", (double)norm);
399a8db3e61SBlaise Bourdin 
400a8db3e61SBlaise Bourdin   PetscCall(DMPlexGlobalToNaturalBegin(dmUA2, UA2, UA2nat));
401a8db3e61SBlaise Bourdin   PetscCall(DMPlexGlobalToNaturalEnd(dmUA2, UA2, UA2nat));
402a8db3e61SBlaise Bourdin   PetscCall(VecAXPY(UA2nat, -1.0, seqUA2));
403a8db3e61SBlaise Bourdin   PetscCall(VecNorm(UA2nat, NORM_INFINITY, &norm));
404a8db3e61SBlaise Bourdin   PetscCheck(norm <= PETSC_SQRT_MACHINE_EPSILON, PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "UAlpha2 ||Vin - Vout|| = %g", (double)norm);
405a8db3e61SBlaise Bourdin 
406a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(dm, &Xnat));
407a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(seqdmU, &Unat));
408a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(seqdmA, &Anat));
409a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(seqdmUA, &UAnat));
410a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(seqdmS, &Snat));
411a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(seqdmUA2, &UA2nat));
412a8db3e61SBlaise Bourdin 
413a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(seqdmUA2, &seqUA2));
414a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(seqdmUA, &seqUA));
415a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(seqdmS, &seqS));
416a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(seqdmA, &seqA));
417a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(seqdmU, &seqU));
418a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(dm, &seqX));
4199371c9d4SSatish Balay   PetscCall(DMDestroy(&seqdmU));
4209371c9d4SSatish Balay   PetscCall(ISDestroy(&seqisU));
4219371c9d4SSatish Balay   PetscCall(DMDestroy(&seqdmA));
4229371c9d4SSatish Balay   PetscCall(ISDestroy(&seqisA));
4239371c9d4SSatish Balay   PetscCall(DMDestroy(&seqdmS));
4249371c9d4SSatish Balay   PetscCall(ISDestroy(&seqisS));
4259371c9d4SSatish Balay   PetscCall(DMDestroy(&seqdmUA));
4269371c9d4SSatish Balay   PetscCall(ISDestroy(&seqisUA));
427a8db3e61SBlaise Bourdin   PetscCall(DMDestroy(&seqdmUA2));
428a8db3e61SBlaise Bourdin 
429a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(dmUA2, &UA2));
430a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(dmUA, &UA));
431a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(dmS, &S));
432a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(dmA, &A));
433a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(dmU, &U));
434a8db3e61SBlaise Bourdin   PetscCall(DMRestoreGlobalVector(pdm, &X));
4359371c9d4SSatish Balay   PetscCall(DMDestroy(&dmU));
4369371c9d4SSatish Balay   PetscCall(ISDestroy(&isU));
4379371c9d4SSatish Balay   PetscCall(DMDestroy(&dmA));
4389371c9d4SSatish Balay   PetscCall(ISDestroy(&isA));
4399371c9d4SSatish Balay   PetscCall(DMDestroy(&dmS));
4409371c9d4SSatish Balay   PetscCall(ISDestroy(&isS));
4419371c9d4SSatish Balay   PetscCall(DMDestroy(&dmUA));
4429371c9d4SSatish Balay   PetscCall(ISDestroy(&isUA));
443a8db3e61SBlaise Bourdin   PetscCall(DMDestroy(&dmUA2));
444a8db3e61SBlaise Bourdin   PetscCall(DMDestroy(&pdm));
445a8db3e61SBlaise Bourdin   if (size > 1) PetscCall(DMDestroy(&dm));
446a8db3e61SBlaise Bourdin   PetscCall(PetscFree2(pStartDepth, pEndDepth));
447a8db3e61SBlaise Bourdin   PetscCall(PetscFinalize());
448a8db3e61SBlaise Bourdin   return 0;
449a8db3e61SBlaise Bourdin }
450a8db3e61SBlaise Bourdin 
451a8db3e61SBlaise Bourdin /*TEST
452a8db3e61SBlaise Bourdin 
453a8db3e61SBlaise Bourdin   build:
454a8db3e61SBlaise Bourdin     requires: !complex parmetis exodusii pnetcdf
455a8db3e61SBlaise Bourdin   # 2D seq
456a8db3e61SBlaise Bourdin   test:
457a8db3e61SBlaise Bourdin     suffix: 0
458a8db3e61SBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -dm_view -petscpartitioner_type parmetis
459a8db3e61SBlaise Bourdin   test:
460a8db3e61SBlaise Bourdin     suffix: 1
461a8db3e61SBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -dm_view -petscpartitioner_type parmetis
462a8db3e61SBlaise Bourdin   test:
463a8db3e61SBlaise Bourdin     suffix: 2
464a8db3e61SBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -dm_view -petscpartitioner_type parmetis
465a8db3e61SBlaise Bourdin 
466a8db3e61SBlaise Bourdin   # 2D par
467a8db3e61SBlaise Bourdin   test:
468a8db3e61SBlaise Bourdin     suffix: 3
469a8db3e61SBlaise Bourdin     nsize: 2
470a8db3e61SBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareT-large.exo -dm_view -petscpartitioner_type parmetis
471a8db3e61SBlaise Bourdin   test:
472a8db3e61SBlaise Bourdin     suffix: 4
473a8db3e61SBlaise Bourdin     nsize: 2
474a8db3e61SBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareQ-large.exo -dm_view -petscpartitioner_type parmetis
475a8db3e61SBlaise Bourdin   test:
476a8db3e61SBlaise Bourdin     suffix: 5
477a8db3e61SBlaise Bourdin     nsize: 2
478a8db3e61SBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourSquareH-large.exo -dm_view -petscpartitioner_type parmetis
479a8db3e61SBlaise Bourdin 
480a8db3e61SBlaise Bourdin   #3d seq
481a8db3e61SBlaise Bourdin   test:
482a8db3e61SBlaise Bourdin     suffix: 6
483a8db3e61SBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -dm_view -petscpartitioner_type parmetis
484a8db3e61SBlaise Bourdin   test:
485a8db3e61SBlaise Bourdin     suffix: 7
486a8db3e61SBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -dm_view -petscpartitioner_type parmetis
487a8db3e61SBlaise Bourdin 
488a8db3e61SBlaise Bourdin   #3d par
489a8db3e61SBlaise Bourdin   test:
490a8db3e61SBlaise Bourdin     suffix: 8
491a8db3e61SBlaise Bourdin     nsize: 2
492a8db3e61SBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickHex-large.exo -dm_view -petscpartitioner_type parmetis
493a8db3e61SBlaise Bourdin   test:
494a8db3e61SBlaise Bourdin     suffix: 9
495a8db3e61SBlaise Bourdin     nsize: 2
496a8db3e61SBlaise Bourdin     args: -i ${wPETSC_DIR}/share/petsc/datafiles/meshes/FourBrickTet-large.exo -dm_view -petscpartitioner_type parmetis
497a8db3e61SBlaise Bourdin 
498a8db3e61SBlaise Bourdin TEST*/
499