xref: /petsc/src/dm/interface/dmi.c (revision 9e47a1b0f33d9bc426ec8fd9a0e0035e6f7196c5)
1af0996ceSBarry Smith #include <petsc/private/dmimpl.h> /*I      "petscdm.h"     I*/
22764a2aaSMatthew G. Knepley #include <petscds.h>
311689aebSJed Brown 
45b8243e1SJed Brown // Greatest common divisor of two nonnegative integers
5d71ae5a4SJacob Faibussowitsch PetscInt PetscGCD(PetscInt a, PetscInt b)
6d71ae5a4SJacob Faibussowitsch {
75b8243e1SJed Brown   while (b != 0) {
85b8243e1SJed Brown     PetscInt tmp = b;
95b8243e1SJed Brown     b            = a % b;
105b8243e1SJed Brown     a            = tmp;
115b8243e1SJed Brown   }
125b8243e1SJed Brown   return a;
135b8243e1SJed Brown }
145b8243e1SJed Brown 
15d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateGlobalVector_Section_Private(DM dm, Vec *vec)
16d71ae5a4SJacob Faibussowitsch {
1711689aebSJed Brown   PetscSection gSection;
18cc30b04dSMatthew G Knepley   PetscInt     localSize, bs, blockSize = -1, pStart, pEnd, p;
195d1c0279SHong Zhang   PetscInt     in[2], out[2];
2011689aebSJed Brown 
2111689aebSJed Brown   PetscFunctionBegin;
229566063dSJacob Faibussowitsch   PetscCall(DMGetGlobalSection(dm, &gSection));
239566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(gSection, &pStart, &pEnd));
2411689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
2511689aebSJed Brown     PetscInt dof, cdof;
2611689aebSJed Brown 
279566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(gSection, p, &dof));
289566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetConstraintDof(gSection, p, &cdof));
290680ce0aSHong Zhang 
305b8243e1SJed Brown     if (dof - cdof > 0) {
315b8243e1SJed Brown       if (blockSize < 0) {
320680ce0aSHong Zhang         /* set blockSize */
330680ce0aSHong Zhang         blockSize = dof - cdof;
345b8243e1SJed Brown       } else {
355b8243e1SJed Brown         blockSize = PetscGCD(dof - cdof, blockSize);
3611689aebSJed Brown       }
3711689aebSJed Brown     }
380680ce0aSHong Zhang   }
395d1c0279SHong Zhang 
402c8be454SMatthew G. Knepley   in[0] = blockSize < 0 ? PETSC_MIN_INT : -blockSize;
415d1c0279SHong Zhang   in[1] = blockSize;
421c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(in, out, 2, MPIU_INT, MPI_MAX, PetscObjectComm((PetscObject)dm)));
435d1c0279SHong Zhang   /* -out[0] = min(blockSize), out[1] = max(blockSize) */
445d1c0279SHong Zhang   if (-out[0] == out[1]) {
455d1c0279SHong Zhang     bs = out[1];
465d1c0279SHong Zhang   } else bs = 1;
475d1c0279SHong Zhang 
485d1c0279SHong Zhang   if (bs < 0) { /* Everyone was empty */
495d1c0279SHong Zhang     blockSize = 1;
505d1c0279SHong Zhang     bs        = 1;
515d1c0279SHong Zhang   }
525d1c0279SHong Zhang 
539566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetConstrainedStorageSize(gSection, &localSize));
547a8be351SBarry Smith   PetscCheck(localSize % blockSize == 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Mismatch between blocksize %" PetscInt_FMT " and local storage size %" PetscInt_FMT, blockSize, localSize);
559566063dSJacob Faibussowitsch   PetscCall(VecCreate(PetscObjectComm((PetscObject)dm), vec));
569566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*vec, localSize, PETSC_DETERMINE));
579566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(*vec, bs));
589566063dSJacob Faibussowitsch   PetscCall(VecSetType(*vec, dm->vectype));
599566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
609566063dSJacob Faibussowitsch   /* PetscCall(VecSetLocalToGlobalMapping(*vec, dm->ltogmap)); */
613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6211689aebSJed Brown }
6311689aebSJed Brown 
64d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateLocalVector_Section_Private(DM dm, Vec *vec)
65d71ae5a4SJacob Faibussowitsch {
66fad22124SMatthew G Knepley   PetscSection section;
6711689aebSJed Brown   PetscInt     localSize, blockSize = -1, pStart, pEnd, p;
6811689aebSJed Brown 
6911689aebSJed Brown   PetscFunctionBegin;
709566063dSJacob Faibussowitsch   PetscCall(DMGetLocalSection(dm, &section));
719566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
7211689aebSJed Brown   for (p = pStart; p < pEnd; ++p) {
7311689aebSJed Brown     PetscInt dof;
7411689aebSJed Brown 
759566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(section, p, &dof));
7611689aebSJed Brown     if ((blockSize < 0) && (dof > 0)) blockSize = dof;
775b8243e1SJed Brown     if (dof > 0) blockSize = PetscGCD(dof, blockSize);
7811689aebSJed Brown   }
799566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(section, &localSize));
809566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, vec));
819566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(*vec, localSize, localSize));
829566063dSJacob Faibussowitsch   PetscCall(VecSetBlockSize(*vec, blockSize));
839566063dSJacob Faibussowitsch   PetscCall(VecSetType(*vec, dm->vectype));
849566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*vec, dm));
853ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8611689aebSJed Brown }
874d9407bcSMatthew G. Knepley 
88*9e47a1b0SMatthew G. Knepley static PetscErrorCode PetscSectionSelectFields_Private(PetscSection s, PetscSection gs, PetscBool useComp, PetscInt numFields, const PetscInt fields[], IS *is)
89d71ae5a4SJacob Faibussowitsch {
904d9407bcSMatthew G. Knepley   PetscInt *subIndices;
91*9e47a1b0SMatthew G. Knepley   PetscInt  bs, bsLocal[2], bsMinMax[2];
92*9e47a1b0SMatthew G. Knepley   PetscInt  pStart, pEnd, subSize = 0, subOff = 0;
934d9407bcSMatthew G. Knepley 
944d9407bcSMatthew G. Knepley   PetscFunctionBegin;
95*9e47a1b0SMatthew G. Knepley   PetscCall(PetscSectionGetChart(gs, &pStart, &pEnd));
96*9e47a1b0SMatthew G. Knepley   if (useComp) {
973a544194SStefano Zampini     PetscInt Nc;
983a544194SStefano Zampini 
99*9e47a1b0SMatthew G. Knepley     PetscCall(PetscSectionGetFieldComponents(s, 0, &Nc));
100*9e47a1b0SMatthew G. Knepley     for (PetscInt f = 0; f < numFields; ++f) PetscCheck(fields[f] < Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "FV component[%" PetscInt_FMT "]: %" PetscInt_FMT " not in [0, %" PetscInt_FMT ")", f, fields[f], Nc);
101*9e47a1b0SMatthew G. Knepley     bs = numFields;
102*9e47a1b0SMatthew G. Knepley     for (PetscInt p = pStart; p < pEnd; ++p) {
103*9e47a1b0SMatthew G. Knepley       PetscInt gdof, fcdof;
104*9e47a1b0SMatthew G. Knepley 
105*9e47a1b0SMatthew G. Knepley       PetscCall(PetscSectionGetDof(gs, p, &gdof));
106*9e47a1b0SMatthew G. Knepley       if (gdof > 0) {
107*9e47a1b0SMatthew G. Knepley         PetscCheck(gdof == Nc, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "This is collocated FV, so there should be %" PetscInt_FMT " dofs in each cell, not %" PetscInt_FMT, Nc, gdof);
108*9e47a1b0SMatthew G. Knepley         PetscCall(PetscSectionGetFieldConstraintDof(s, p, 0, &fcdof));
109*9e47a1b0SMatthew G. Knepley         PetscCheck(!fcdof, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "This is collocated FV, so there should be 0 constraints in each cell, not %" PetscInt_FMT, fcdof);
110*9e47a1b0SMatthew G. Knepley         subSize += numFields;
111*9e47a1b0SMatthew G. Knepley       }
112*9e47a1b0SMatthew G. Knepley     }
113*9e47a1b0SMatthew G. Knepley     PetscCall(PetscMalloc1(subSize, &subIndices));
114*9e47a1b0SMatthew G. Knepley     for (PetscInt p = pStart; p < pEnd; ++p) {
115*9e47a1b0SMatthew G. Knepley       PetscInt gdof, goff;
116*9e47a1b0SMatthew G. Knepley 
117*9e47a1b0SMatthew G. Knepley       PetscCall(PetscSectionGetDof(gs, p, &gdof));
118*9e47a1b0SMatthew G. Knepley       if (gdof > 0) {
119*9e47a1b0SMatthew G. Knepley         PetscCall(PetscSectionGetOffset(gs, p, &goff));
120*9e47a1b0SMatthew G. Knepley         for (PetscInt f = 0; f < numFields; ++f, ++subOff) subIndices[subOff] = goff + fields[f];
121*9e47a1b0SMatthew G. Knepley       }
122*9e47a1b0SMatthew G. Knepley     }
123*9e47a1b0SMatthew G. Knepley     PetscCheck(subSize == subOff, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "The offset array size %" PetscInt_FMT " != %" PetscInt_FMT " the number of indices", subSize, subOff);
124*9e47a1b0SMatthew G. Knepley     PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)gs), subSize, subIndices, PETSC_OWN_POINTER, is));
125*9e47a1b0SMatthew G. Knepley     PetscCall(ISSetBlockSize(*is, bs));
126*9e47a1b0SMatthew G. Knepley     PetscFunctionReturn(PETSC_SUCCESS);
127*9e47a1b0SMatthew G. Knepley   }
128*9e47a1b0SMatthew G. Knepley   bs = 0;
129*9e47a1b0SMatthew G. Knepley   for (PetscInt f = 0; f < numFields; ++f) {
130*9e47a1b0SMatthew G. Knepley     PetscInt Nc;
131*9e47a1b0SMatthew G. Knepley 
132*9e47a1b0SMatthew G. Knepley     PetscCall(PetscSectionGetFieldComponents(s, fields[f], &Nc));
1333a544194SStefano Zampini     bs += Nc;
1343a544194SStefano Zampini   }
135*9e47a1b0SMatthew G. Knepley   for (PetscInt p = pStart; p < pEnd; ++p) {
136b9eca265Ssarens     PetscInt gdof, pSubSize = 0;
1374d9407bcSMatthew G. Knepley 
138*9e47a1b0SMatthew G. Knepley     PetscCall(PetscSectionGetDof(gs, p, &gdof));
1394d9407bcSMatthew G. Knepley     if (gdof > 0) {
140*9e47a1b0SMatthew G. Knepley       for (PetscInt f = 0; f < numFields; ++f) {
1414d9407bcSMatthew G. Knepley         PetscInt fdof, fcdof;
1424d9407bcSMatthew G. Knepley 
143*9e47a1b0SMatthew G. Knepley         PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof));
144*9e47a1b0SMatthew G. Knepley         PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &fcdof));
145b9eca265Ssarens         pSubSize += fdof - fcdof;
146b9eca265Ssarens       }
147b9eca265Ssarens       subSize += pSubSize;
1483a544194SStefano Zampini       if (pSubSize && bs != pSubSize) {
149*9e47a1b0SMatthew G. Knepley         // Layout does not admit a pointwise block size
150b9eca265Ssarens         bs = 1;
1514d9407bcSMatthew G. Knepley       }
1524d9407bcSMatthew G. Knepley     }
1534d9407bcSMatthew G. Knepley   }
154*9e47a1b0SMatthew G. Knepley   // Must have same blocksize on all procs (some might have no points)
1559371c9d4SSatish Balay   bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs;
1569371c9d4SSatish Balay   bsLocal[1] = bs;
157*9e47a1b0SMatthew G. Knepley   PetscCall(PetscGlobalMinMaxInt(PetscObjectComm((PetscObject)gs), bsLocal, bsMinMax));
1589371c9d4SSatish Balay   if (bsMinMax[0] != bsMinMax[1]) {
1599371c9d4SSatish Balay     bs = 1;
1609371c9d4SSatish Balay   } else {
1619371c9d4SSatish Balay     bs = bsMinMax[0];
1629371c9d4SSatish Balay   }
1639566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(subSize, &subIndices));
164*9e47a1b0SMatthew G. Knepley   for (PetscInt p = pStart; p < pEnd; ++p) {
1654d9407bcSMatthew G. Knepley     PetscInt gdof, goff;
1664d9407bcSMatthew G. Knepley 
167*9e47a1b0SMatthew G. Knepley     PetscCall(PetscSectionGetDof(gs, p, &gdof));
1684d9407bcSMatthew G. Knepley     if (gdof > 0) {
169*9e47a1b0SMatthew G. Knepley       PetscCall(PetscSectionGetOffset(gs, p, &goff));
170*9e47a1b0SMatthew G. Knepley       for (PetscInt f = 0; f < numFields; ++f) {
171*9e47a1b0SMatthew G. Knepley         PetscInt fdof, fcdof, poff = 0;
1724d9407bcSMatthew G. Knepley 
1734d9407bcSMatthew G. Knepley         /* Can get rid of this loop by storing field information in the global section */
174*9e47a1b0SMatthew G. Knepley         for (PetscInt f2 = 0; f2 < fields[f]; ++f2) {
175*9e47a1b0SMatthew G. Knepley           PetscCall(PetscSectionGetFieldDof(s, p, f2, &fdof));
176*9e47a1b0SMatthew G. Knepley           PetscCall(PetscSectionGetFieldConstraintDof(s, p, f2, &fcdof));
1774d9407bcSMatthew G. Knepley           poff += fdof - fcdof;
1784d9407bcSMatthew G. Knepley         }
179*9e47a1b0SMatthew G. Knepley         PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof));
180*9e47a1b0SMatthew G. Knepley         PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &fcdof));
181*9e47a1b0SMatthew G. Knepley         for (PetscInt fc = 0; fc < fdof - fcdof; ++fc, ++subOff) subIndices[subOff] = goff + poff + fc;
1824d9407bcSMatthew G. Knepley       }
1834d9407bcSMatthew G. Knepley     }
1844d9407bcSMatthew G. Knepley   }
185*9e47a1b0SMatthew G. Knepley   PetscCheck(subSize == subOff, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "The offset array size %" PetscInt_FMT " != %" PetscInt_FMT " the number of indices", subSize, subOff);
186*9e47a1b0SMatthew G. Knepley   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)gs), subSize, subIndices, PETSC_OWN_POINTER, is));
187627349f5SMatthew G. Knepley   if (bs > 1) {
188*9e47a1b0SMatthew G. Knepley     // We need to check that the block size does not come from non-contiguous fields
189*9e47a1b0SMatthew G. Knepley     PetscInt set = 1, rset = 1;
190*9e47a1b0SMatthew G. Knepley     for (PetscInt i = 0; i < subSize; i += bs) {
191*9e47a1b0SMatthew G. Knepley       for (PetscInt j = 0; j < bs; ++j) {
1929371c9d4SSatish Balay         if (subIndices[i + j] != subIndices[i] + j) {
1939371c9d4SSatish Balay           set = 0;
1949371c9d4SSatish Balay           break;
1959371c9d4SSatish Balay         }
196627349f5SMatthew G. Knepley       }
197627349f5SMatthew G. Knepley     }
198*9e47a1b0SMatthew G. Knepley     PetscCall(MPIU_Allreduce(&set, &rset, 1, MPIU_INT, MPI_PROD, PetscObjectComm((PetscObject)gs)));
19988206212SAlexis Marboeuf     if (rset) PetscCall(ISSetBlockSize(*is, bs));
200627349f5SMatthew G. Knepley   }
201*9e47a1b0SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
2024d9407bcSMatthew G. Knepley }
203*9e47a1b0SMatthew G. Knepley 
204*9e47a1b0SMatthew G. Knepley static PetscErrorCode DMSelectFields_Private(DM dm, PetscSection section, PetscBool useComp, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
205*9e47a1b0SMatthew G. Knepley {
2064d9407bcSMatthew G. Knepley   PetscSection subsection;
2074d9407bcSMatthew G. Knepley   PetscBool    haveNull = PETSC_FALSE;
208*9e47a1b0SMatthew G. Knepley   PetscInt     nf = 0, of = 0;
2094d9407bcSMatthew G. Knepley 
210*9e47a1b0SMatthew G. Knepley   PetscFunctionBegin;
211*9e47a1b0SMatthew G. Knepley   if (useComp) {
212*9e47a1b0SMatthew G. Knepley     PetscCall(PetscSectionCreateComponentSubsection(section, numFields, fields, &subsection));
2139566063dSJacob Faibussowitsch     PetscCall(DMSetLocalSection(*subdm, subsection));
2149566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&subsection));
215*9e47a1b0SMatthew G. Knepley     (*subdm)->nullspaceConstructors[0] = dm->nullspaceConstructors[0];
216e5e52638SMatthew G. Knepley     if (dm->probs) {
217*9e47a1b0SMatthew G. Knepley       PetscFV  fv, fvNew;
218*9e47a1b0SMatthew G. Knepley       PetscInt fnum[1] = {0};
2190f21e855SMatthew G. Knepley 
220*9e47a1b0SMatthew G. Knepley       PetscCall(DMSetNumFields(*subdm, 1));
221*9e47a1b0SMatthew G. Knepley       PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fv));
222*9e47a1b0SMatthew G. Knepley       PetscCall(PetscFVClone(fv, &fvNew));
223*9e47a1b0SMatthew G. Knepley       PetscCall(PetscFVSetNumComponents(fvNew, numFields));
224*9e47a1b0SMatthew G. Knepley       PetscCall(DMSetField(*subdm, 0, NULL, (PetscObject)fvNew));
225*9e47a1b0SMatthew G. Knepley       PetscCall(PetscFVDestroy(&fvNew));
2269566063dSJacob Faibussowitsch       PetscCall(DMCreateDS(*subdm));
227f646a522SMatthew G. Knepley       if (numFields == 1 && is) {
2280f21e855SMatthew G. Knepley         PetscObject disc, space, pmat;
2294d9407bcSMatthew G. Knepley 
2309566063dSJacob Faibussowitsch         PetscCall(DMGetField(*subdm, 0, NULL, &disc));
2319566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nullspace", &space));
2329566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space));
2339566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "nearnullspace", &space));
2349566063dSJacob Faibussowitsch         if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space));
2359566063dSJacob Faibussowitsch         PetscCall(PetscObjectQuery(disc, "pmat", &pmat));
2369566063dSJacob Faibussowitsch         if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat));
2374d9407bcSMatthew G. Knepley       }
238*9e47a1b0SMatthew G. Knepley       PetscCall(PetscDSCopyConstants(dm->probs[0].ds, (*subdm)->probs[0].ds));
239*9e47a1b0SMatthew G. Knepley       PetscCall(PetscDSCopyBoundary(dm->probs[0].ds, 1, fnum, (*subdm)->probs[0].ds));
240*9e47a1b0SMatthew G. Knepley       PetscCall(PetscDSSelectEquations(dm->probs[0].ds, 1, fnum, (*subdm)->probs[0].ds));
241*9e47a1b0SMatthew G. Knepley     }
242*9e47a1b0SMatthew G. Knepley     if ((*subdm)->nullspaceConstructors[0] && is) {
243*9e47a1b0SMatthew G. Knepley       MatNullSpace nullSpace;
244*9e47a1b0SMatthew G. Knepley 
245*9e47a1b0SMatthew G. Knepley       PetscCall((*(*subdm)->nullspaceConstructors[0])(*subdm, 0, 0, &nullSpace));
246*9e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace));
247*9e47a1b0SMatthew G. Knepley       PetscCall(MatNullSpaceDestroy(&nullSpace));
248*9e47a1b0SMatthew G. Knepley     }
249*9e47a1b0SMatthew G. Knepley     if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh));
250*9e47a1b0SMatthew G. Knepley     PetscFunctionReturn(PETSC_SUCCESS);
251*9e47a1b0SMatthew G. Knepley   }
252*9e47a1b0SMatthew G. Knepley   PetscCall(PetscSectionCreateSubsection(section, numFields, fields, &subsection));
253*9e47a1b0SMatthew G. Knepley   PetscCall(DMSetLocalSection(*subdm, subsection));
254*9e47a1b0SMatthew G. Knepley   PetscCall(PetscSectionDestroy(&subsection));
255*9e47a1b0SMatthew G. Knepley   for (PetscInt f = 0; f < numFields; ++f) {
256*9e47a1b0SMatthew G. Knepley     (*subdm)->nullspaceConstructors[f] = dm->nullspaceConstructors[fields[f]];
257*9e47a1b0SMatthew G. Knepley     if ((*subdm)->nullspaceConstructors[f]) {
258*9e47a1b0SMatthew G. Knepley       haveNull = PETSC_TRUE;
259*9e47a1b0SMatthew G. Knepley       nf       = f;
260*9e47a1b0SMatthew G. Knepley       of       = fields[f];
261*9e47a1b0SMatthew G. Knepley     }
262*9e47a1b0SMatthew G. Knepley   }
263*9e47a1b0SMatthew G. Knepley   if (dm->probs) {
264*9e47a1b0SMatthew G. Knepley     PetscCall(DMSetNumFields(*subdm, numFields));
265*9e47a1b0SMatthew G. Knepley     for (PetscInt f = 0; f < numFields; ++f) {
266*9e47a1b0SMatthew G. Knepley       PetscObject disc;
267*9e47a1b0SMatthew G. Knepley 
268*9e47a1b0SMatthew G. Knepley       PetscCall(DMGetField(dm, fields[f], NULL, &disc));
269*9e47a1b0SMatthew G. Knepley       PetscCall(DMSetField(*subdm, f, NULL, disc));
270*9e47a1b0SMatthew G. Knepley     }
271*9e47a1b0SMatthew G. Knepley     // TODO: if only FV, then cut down the components
272*9e47a1b0SMatthew G. Knepley     PetscCall(DMCreateDS(*subdm));
273*9e47a1b0SMatthew G. Knepley     if (numFields == 1 && is) {
274*9e47a1b0SMatthew G. Knepley       PetscObject disc, space, pmat;
275*9e47a1b0SMatthew G. Knepley 
276*9e47a1b0SMatthew G. Knepley       PetscCall(DMGetField(*subdm, 0, NULL, &disc));
277*9e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectQuery(disc, "nullspace", &space));
278*9e47a1b0SMatthew G. Knepley       if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", space));
279*9e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectQuery(disc, "nearnullspace", &space));
280*9e47a1b0SMatthew G. Knepley       if (space) PetscCall(PetscObjectCompose((PetscObject)*is, "nearnullspace", space));
281*9e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectQuery(disc, "pmat", &pmat));
282*9e47a1b0SMatthew G. Knepley       if (pmat) PetscCall(PetscObjectCompose((PetscObject)*is, "pmat", pmat));
283*9e47a1b0SMatthew G. Knepley     }
284*9e47a1b0SMatthew G. Knepley     // Check if DSes record their DM fields
285e21d936aSMatthew G. Knepley     if (dm->probs[0].fields) {
2869310035eSMatthew G. Knepley       PetscInt d, e;
2879310035eSMatthew G. Knepley 
2889310035eSMatthew G. Knepley       for (d = 0, e = 0; d < dm->Nds && e < (*subdm)->Nds; ++d) {
2899310035eSMatthew G. Knepley         const PetscInt  Nf = dm->probs[d].ds->Nf;
2909310035eSMatthew G. Knepley         const PetscInt *fld;
2919310035eSMatthew G. Knepley         PetscInt        f, g;
2929310035eSMatthew G. Knepley 
2939566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(dm->probs[d].fields, &fld));
2949310035eSMatthew G. Knepley         for (f = 0; f < Nf; ++f) {
2959371c9d4SSatish Balay           for (g = 0; g < numFields; ++g)
2969371c9d4SSatish Balay             if (fld[f] == fields[g]) break;
2979310035eSMatthew G. Knepley           if (g < numFields) break;
2989310035eSMatthew G. Knepley         }
2999566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(dm->probs[d].fields, &fld));
3009310035eSMatthew G. Knepley         if (f == Nf) continue;
3019566063dSJacob Faibussowitsch         PetscCall(PetscDSCopyConstants(dm->probs[d].ds, (*subdm)->probs[e].ds));
3029566063dSJacob Faibussowitsch         PetscCall(PetscDSCopyBoundary(dm->probs[d].ds, numFields, fields, (*subdm)->probs[e].ds));
303*9e47a1b0SMatthew G. Knepley         // Translate DM fields to DS fields
3049310035eSMatthew G. Knepley         {
305e21d936aSMatthew G. Knepley           IS              infields, dsfields;
30674a61aaaSMatthew G. Knepley           const PetscInt *fld, *ofld;
30774a61aaaSMatthew G. Knepley           PetscInt       *fidx;
308*9e47a1b0SMatthew G. Knepley           PetscInt        onf, nf;
309e21d936aSMatthew G. Knepley 
3109566063dSJacob Faibussowitsch           PetscCall(ISCreateGeneral(PETSC_COMM_SELF, numFields, fields, PETSC_USE_POINTER, &infields));
3119566063dSJacob Faibussowitsch           PetscCall(ISIntersect(infields, dm->probs[d].fields, &dsfields));
3129566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&infields));
3139566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(dsfields, &nf));
3147a8be351SBarry Smith           PetscCheck(nf, PETSC_COMM_SELF, PETSC_ERR_PLIB, "DS cannot be supported on 0 fields");
3159566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(dsfields, &fld));
3169566063dSJacob Faibussowitsch           PetscCall(ISGetLocalSize(dm->probs[d].fields, &onf));
3179566063dSJacob Faibussowitsch           PetscCall(ISGetIndices(dm->probs[d].fields, &ofld));
3189566063dSJacob Faibussowitsch           PetscCall(PetscMalloc1(nf, &fidx));
319*9e47a1b0SMatthew G. Knepley           for (PetscInt f = 0, g = 0; f < onf && g < nf; ++f) {
32074a61aaaSMatthew G. Knepley             if (ofld[f] == fld[g]) fidx[g++] = f;
32174a61aaaSMatthew G. Knepley           }
3229566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(dm->probs[d].fields, &ofld));
3239566063dSJacob Faibussowitsch           PetscCall(ISRestoreIndices(dsfields, &fld));
3249566063dSJacob Faibussowitsch           PetscCall(ISDestroy(&dsfields));
3259566063dSJacob Faibussowitsch           PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds));
3269566063dSJacob Faibussowitsch           PetscCall(PetscDSSelectEquations(dm->probs[0].ds, nf, fidx, (*subdm)->probs[0].ds));
3279566063dSJacob Faibussowitsch           PetscCall(PetscFree(fidx));
3289310035eSMatthew G. Knepley         }
3299310035eSMatthew G. Knepley         ++e;
3309310035eSMatthew G. Knepley       }
331e21d936aSMatthew G. Knepley     } else {
3329566063dSJacob Faibussowitsch       PetscCall(PetscDSCopyConstants(dm->probs[0].ds, (*subdm)->probs[0].ds));
3339566063dSJacob Faibussowitsch       PetscCall(PetscDSCopyBoundary(dm->probs[0].ds, PETSC_DETERMINE, NULL, (*subdm)->probs[0].ds));
3349566063dSJacob Faibussowitsch       PetscCall(PetscDSSelectDiscretizations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds));
3359566063dSJacob Faibussowitsch       PetscCall(PetscDSSelectEquations(dm->probs[0].ds, numFields, fields, (*subdm)->probs[0].ds));
336d5af271aSMatthew G. Knepley     }
337e21d936aSMatthew G. Knepley   }
3388cda7954SMatthew G. Knepley   if (haveNull && is) {
3398cda7954SMatthew G. Knepley     MatNullSpace nullSpace;
3408cda7954SMatthew G. Knepley 
3419566063dSJacob Faibussowitsch     PetscCall((*(*subdm)->nullspaceConstructors[nf])(*subdm, of, nf, &nullSpace));
3429566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)*is, "nullspace", (PetscObject)nullSpace));
3439566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&nullSpace));
3448cda7954SMatthew G. Knepley   }
34548a46eb9SPierre Jolivet   if (dm->coarseMesh) PetscCall(DMCreateSubDM(dm->coarseMesh, numFields, fields, NULL, &(*subdm)->coarseMesh));
346*9e47a1b0SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
3474d9407bcSMatthew G. Knepley }
348*9e47a1b0SMatthew G. Knepley 
349*9e47a1b0SMatthew G. Knepley /*@C
350*9e47a1b0SMatthew G. Knepley   DMCreateSectionSubDM - Returns an `IS` and `subDM` containing a `PetscSection` that encapsulates a subproblem defined by a subset of the fields in a `PetscSection` in the `DM`.
351*9e47a1b0SMatthew G. Knepley 
352*9e47a1b0SMatthew G. Knepley   Not Collective
353*9e47a1b0SMatthew G. Knepley 
354*9e47a1b0SMatthew G. Knepley   Input Parameters:
355*9e47a1b0SMatthew G. Knepley + dm        - The `DM` object
356*9e47a1b0SMatthew G. Knepley . numFields - The number of fields to incorporate into `subdm`
357*9e47a1b0SMatthew G. Knepley - fields    - The field numbers of the selected fields
358*9e47a1b0SMatthew G. Knepley 
359*9e47a1b0SMatthew G. Knepley   Output Parameters:
360*9e47a1b0SMatthew G. Knepley + is    - The global indices for the subproblem or `NULL`
361*9e47a1b0SMatthew G. Knepley - subdm - The `DM` for the subproblem, which must already have be cloned from `dm` or `NULL`
362*9e47a1b0SMatthew G. Knepley 
363*9e47a1b0SMatthew G. Knepley   Level: intermediate
364*9e47a1b0SMatthew G. Knepley 
365*9e47a1b0SMatthew G. Knepley   Notes:
366*9e47a1b0SMatthew G. Knepley   If `is` and `subdm` are both `NULL` this does nothing
367*9e47a1b0SMatthew G. Knepley 
368*9e47a1b0SMatthew G. Knepley   If there is only one field, and it is a `PetscFV`, then the selected fields will be interpreted as selected components.
369*9e47a1b0SMatthew G. Knepley 
370*9e47a1b0SMatthew G. Knepley .seealso: `DMCreateSubDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()`
371*9e47a1b0SMatthew G. Knepley @*/
372*9e47a1b0SMatthew G. Knepley PetscErrorCode DMCreateSectionSubDM(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
373*9e47a1b0SMatthew G. Knepley {
374*9e47a1b0SMatthew G. Knepley   PetscSection section, sectionGlobal;
375*9e47a1b0SMatthew G. Knepley   PetscInt     Nf;
376*9e47a1b0SMatthew G. Knepley   PetscBool    useComp = PETSC_FALSE;
377*9e47a1b0SMatthew G. Knepley 
378*9e47a1b0SMatthew G. Knepley   PetscFunctionBegin;
379*9e47a1b0SMatthew G. Knepley   if (!numFields) PetscFunctionReturn(PETSC_SUCCESS);
380*9e47a1b0SMatthew G. Knepley   PetscCall(DMGetLocalSection(dm, &section));
381*9e47a1b0SMatthew G. Knepley   PetscCall(DMGetGlobalSection(dm, &sectionGlobal));
382*9e47a1b0SMatthew G. Knepley   PetscCheck(section, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
383*9e47a1b0SMatthew G. Knepley   PetscCheck(sectionGlobal, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
384*9e47a1b0SMatthew G. Knepley   PetscCall(PetscSectionGetNumFields(section, &Nf));
385*9e47a1b0SMatthew G. Knepley   PetscCheck(numFields <= Nf, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Number of requested fields %" PetscInt_FMT " greater than number of DM fields %" PetscInt_FMT, numFields, Nf);
386*9e47a1b0SMatthew G. Knepley 
387*9e47a1b0SMatthew G. Knepley   // If there is only one field and it is a PetscFV, then interpret selection fields as components
388*9e47a1b0SMatthew G. Knepley   {
389*9e47a1b0SMatthew G. Knepley     PetscInt Nf;
390*9e47a1b0SMatthew G. Knepley 
391*9e47a1b0SMatthew G. Knepley     PetscCall(DMGetNumFields(dm, &Nf));
392*9e47a1b0SMatthew G. Knepley     if (Nf == 1) {
393*9e47a1b0SMatthew G. Knepley       PetscObject  obj;
394*9e47a1b0SMatthew G. Knepley       PetscClassId id;
395*9e47a1b0SMatthew G. Knepley 
396*9e47a1b0SMatthew G. Knepley       PetscCall(DMGetField(dm, 0, NULL, &obj));
397*9e47a1b0SMatthew G. Knepley       PetscCall(PetscObjectGetClassId(obj, &id));
398*9e47a1b0SMatthew G. Knepley       if (id == PETSCFV_CLASSID) useComp = PETSC_TRUE;
399*9e47a1b0SMatthew G. Knepley     }
400*9e47a1b0SMatthew G. Knepley   }
401*9e47a1b0SMatthew G. Knepley   if (is) PetscCall(PetscSectionSelectFields_Private(section, sectionGlobal, useComp, numFields, fields, is));
402*9e47a1b0SMatthew G. Knepley   if (subdm) PetscCall(DMSelectFields_Private(dm, section, useComp, numFields, fields, is, subdm));
4033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4044d9407bcSMatthew G. Knepley }
4052adcc780SMatthew G. Knepley 
406792b654fSMatthew G. Knepley /*@C
4070b3275a6SBarry Smith   DMCreateSectionSuperDM - Returns an arrays of `IS` and a `DM` containing a `PetscSection` that encapsulates a superproblem defined by the array of `DM` and their `PetscSection`
408792b654fSMatthew G. Knepley 
40920f4b53cSBarry Smith   Not Collective
410792b654fSMatthew G. Knepley 
411d8d19677SJose E. Roman   Input Parameters:
4120b3275a6SBarry Smith + dms - The `DM` objects, the must all have the same topology; for example obtained with `DMClone()`
4130b3275a6SBarry Smith - len - The number of `DM` in `dms`
414792b654fSMatthew G. Knepley 
415792b654fSMatthew G. Knepley   Output Parameters:
41620f4b53cSBarry Smith + is      - The global indices for the subproblem, or `NULL`
4170b3275a6SBarry Smith - superdm - The `DM` for the superproblem, which must already have be cloned and contain the same topology as the `dms`
418792b654fSMatthew G. Knepley 
419792b654fSMatthew G. Knepley   Level: intermediate
420792b654fSMatthew G. Knepley 
421a4e35b19SJacob Faibussowitsch .seealso: `DMCreateSuperDM()`, `DMGetLocalSection()`, `DMPlexSetMigrationSF()`, `DMView()`
422792b654fSMatthew G. Knepley @*/
423d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateSectionSuperDM(DM dms[], PetscInt len, IS **is, DM *superdm)
424d71ae5a4SJacob Faibussowitsch {
4259796d432SMatthew G. Knepley   MPI_Comm     comm;
4269796d432SMatthew G. Knepley   PetscSection supersection, *sections, *sectionGlobals;
4278cda7954SMatthew G. Knepley   PetscInt    *Nfs, Nf = 0, f, supf, oldf = -1, nullf = -1, i;
4289796d432SMatthew G. Knepley   PetscBool    haveNull = PETSC_FALSE;
4292adcc780SMatthew G. Knepley 
4302adcc780SMatthew G. Knepley   PetscFunctionBegin;
4319566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dms[0], &comm));
4329796d432SMatthew G. Knepley   /* Pull out local and global sections */
4339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(len, &Nfs, len, &sections, len, &sectionGlobals));
4342adcc780SMatthew G. Knepley   for (i = 0; i < len; ++i) {
4359566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dms[i], &sections[i]));
4369566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(dms[i], &sectionGlobals[i]));
4377a8be351SBarry Smith     PetscCheck(sections[i], comm, PETSC_ERR_ARG_WRONG, "Must set default section for DM before splitting fields");
4387a8be351SBarry Smith     PetscCheck(sectionGlobals[i], comm, PETSC_ERR_ARG_WRONG, "Must set default global section for DM before splitting fields");
4399566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetNumFields(sections[i], &Nfs[i]));
4402adcc780SMatthew G. Knepley     Nf += Nfs[i];
4412adcc780SMatthew G. Knepley   }
4429796d432SMatthew G. Knepley   /* Create the supersection */
4439566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreateSupersection(sections, len, &supersection));
4449566063dSJacob Faibussowitsch   PetscCall(DMSetLocalSection(*superdm, supersection));
4459796d432SMatthew G. Knepley   /* Create ISes */
4462adcc780SMatthew G. Knepley   if (is) {
4479796d432SMatthew G. Knepley     PetscSection supersectionGlobal;
4489796d432SMatthew G. Knepley     PetscInt     bs = -1, startf = 0;
4492adcc780SMatthew G. Knepley 
4509566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(len, is));
4519566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalSection(*superdm, &supersectionGlobal));
4529796d432SMatthew G. Knepley     for (i = 0; i < len; startf += Nfs[i], ++i) {
4539796d432SMatthew G. Knepley       PetscInt *subIndices;
454ec4c761aSMatthew G. Knepley       PetscInt  subSize, subOff, pStart, pEnd, p, start, end, dummy;
4552adcc780SMatthew G. Knepley 
4569566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetChart(sectionGlobals[i], &pStart, &pEnd));
4579566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobals[i], &subSize));
4589566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(subSize, &subIndices));
4599796d432SMatthew G. Knepley       for (p = pStart, subOff = 0; p < pEnd; ++p) {
4609796d432SMatthew G. Knepley         PetscInt gdof, gcdof, gtdof, d;
4612adcc780SMatthew G. Knepley 
4629566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(sectionGlobals[i], p, &gdof));
4639566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetConstraintDof(sections[i], p, &gcdof));
4642adcc780SMatthew G. Knepley         gtdof = gdof - gcdof;
4652adcc780SMatthew G. Knepley         if (gdof > 0 && gtdof) {
4669371c9d4SSatish Balay           if (bs < 0) {
4679371c9d4SSatish Balay             bs = gtdof;
4689371c9d4SSatish Balay           } else if (bs != gtdof) {
4699371c9d4SSatish Balay             bs = 1;
4709371c9d4SSatish Balay           }
4719566063dSJacob Faibussowitsch           PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf, &start, &dummy));
4729566063dSJacob Faibussowitsch           PetscCall(DMGetGlobalFieldOffset_Private(*superdm, p, startf + Nfs[i] - 1, &dummy, &end));
47363a3b9bcSJacob Faibussowitsch           PetscCheck(end - start == gtdof, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Invalid number of global dofs %" PetscInt_FMT " != %" PetscInt_FMT " for dm %" PetscInt_FMT " on point %" PetscInt_FMT, end - start, gtdof, i, p);
4749796d432SMatthew G. Knepley           for (d = start; d < end; ++d, ++subOff) subIndices[subOff] = d;
4752adcc780SMatthew G. Knepley         }
4762adcc780SMatthew G. Knepley       }
4779566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(comm, subSize, subIndices, PETSC_OWN_POINTER, &(*is)[i]));
4782adcc780SMatthew G. Knepley       /* Must have same blocksize on all procs (some might have no points) */
4799796d432SMatthew G. Knepley       {
4809796d432SMatthew G. Knepley         PetscInt bs = -1, bsLocal[2], bsMinMax[2];
4819796d432SMatthew G. Knepley 
4829371c9d4SSatish Balay         bsLocal[0] = bs < 0 ? PETSC_MAX_INT : bs;
4839371c9d4SSatish Balay         bsLocal[1] = bs;
4849566063dSJacob Faibussowitsch         PetscCall(PetscGlobalMinMaxInt(comm, bsLocal, bsMinMax));
4859371c9d4SSatish Balay         if (bsMinMax[0] != bsMinMax[1]) {
4869371c9d4SSatish Balay           bs = 1;
4879371c9d4SSatish Balay         } else {
4889371c9d4SSatish Balay           bs = bsMinMax[0];
4899371c9d4SSatish Balay         }
4909566063dSJacob Faibussowitsch         PetscCall(ISSetBlockSize((*is)[i], bs));
4912adcc780SMatthew G. Knepley       }
4922adcc780SMatthew G. Knepley     }
4932adcc780SMatthew G. Knepley   }
4949796d432SMatthew G. Knepley   /* Preserve discretizations */
495e5e52638SMatthew G. Knepley   if (len && dms[0]->probs) {
4969566063dSJacob Faibussowitsch     PetscCall(DMSetNumFields(*superdm, Nf));
4979796d432SMatthew G. Knepley     for (i = 0, supf = 0; i < len; ++i) {
4989796d432SMatthew G. Knepley       for (f = 0; f < Nfs[i]; ++f, ++supf) {
4992adcc780SMatthew G. Knepley         PetscObject disc;
5002adcc780SMatthew G. Knepley 
5019566063dSJacob Faibussowitsch         PetscCall(DMGetField(dms[i], f, NULL, &disc));
5029566063dSJacob Faibussowitsch         PetscCall(DMSetField(*superdm, supf, NULL, disc));
5032adcc780SMatthew G. Knepley       }
5042adcc780SMatthew G. Knepley     }
5059566063dSJacob Faibussowitsch     PetscCall(DMCreateDS(*superdm));
5062adcc780SMatthew G. Knepley   }
5079796d432SMatthew G. Knepley   /* Preserve nullspaces */
5089796d432SMatthew G. Knepley   for (i = 0, supf = 0; i < len; ++i) {
5099796d432SMatthew G. Knepley     for (f = 0; f < Nfs[i]; ++f, ++supf) {
5109796d432SMatthew G. Knepley       (*superdm)->nullspaceConstructors[supf] = dms[i]->nullspaceConstructors[f];
5119796d432SMatthew G. Knepley       if ((*superdm)->nullspaceConstructors[supf]) {
5129796d432SMatthew G. Knepley         haveNull = PETSC_TRUE;
5139796d432SMatthew G. Knepley         nullf    = supf;
5148cda7954SMatthew G. Knepley         oldf     = f;
5152adcc780SMatthew G. Knepley       }
5169796d432SMatthew G. Knepley     }
5179796d432SMatthew G. Knepley   }
5189796d432SMatthew G. Knepley   /* Attach nullspace to IS */
5199796d432SMatthew G. Knepley   if (haveNull && is) {
5209796d432SMatthew G. Knepley     MatNullSpace nullSpace;
5219796d432SMatthew G. Knepley 
5229566063dSJacob Faibussowitsch     PetscCall((*(*superdm)->nullspaceConstructors[nullf])(*superdm, oldf, nullf, &nullSpace));
5239566063dSJacob Faibussowitsch     PetscCall(PetscObjectCompose((PetscObject)(*is)[nullf], "nullspace", (PetscObject)nullSpace));
5249566063dSJacob Faibussowitsch     PetscCall(MatNullSpaceDestroy(&nullSpace));
5259796d432SMatthew G. Knepley   }
5269566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&supersection));
5279566063dSJacob Faibussowitsch   PetscCall(PetscFree3(Nfs, sections, sectionGlobals));
5283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5292adcc780SMatthew G. Knepley }
530