xref: /petsc/src/dm/dt/dualspace/impls/sum/dualspacesum.c (revision f4f49eeac7efa77fffa46b7ff95a3ed169f659ed)
12dce792eSToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
22dce792eSToby Isaac 
32dce792eSToby Isaac /*@
42dce792eSToby Isaac   PetscDualSpaceSumGetNumSubspaces - Get the number of spaces in the sum space
52dce792eSToby Isaac 
62dce792eSToby Isaac   Input Parameter:
72dce792eSToby Isaac . sp - the dual space object
82dce792eSToby Isaac 
92dce792eSToby Isaac   Output Parameter:
102dce792eSToby Isaac . numSumSpaces - the number of spaces
112dce792eSToby Isaac 
122dce792eSToby Isaac   Level: intermediate
132dce792eSToby Isaac 
142dce792eSToby Isaac   Note:
152dce792eSToby Isaac   The name NumSubspaces is slightly misleading because it is actually getting the number of defining spaces of the sum, not a number of Subspaces of it
162dce792eSToby Isaac 
172dce792eSToby Isaac .seealso: `PETSCDUALSPACESUM`, `PetscDualSpace`, `PetscDualSpaceSumSetNumSubspaces()`
182dce792eSToby Isaac @*/
192dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumGetNumSubspaces(PetscDualSpace sp, PetscInt *numSumSpaces)
202dce792eSToby Isaac {
212dce792eSToby Isaac   PetscFunctionBegin;
222dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
232dce792eSToby Isaac   PetscAssertPointer(numSumSpaces, 2);
242dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumGetNumSubspaces_C", (PetscDualSpace, PetscInt *), (sp, numSumSpaces));
252dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
262dce792eSToby Isaac }
272dce792eSToby Isaac 
282dce792eSToby Isaac /*@
292dce792eSToby Isaac   PetscDualSpaceSumSetNumSubspaces - Set the number of spaces in the sum space
302dce792eSToby Isaac 
312dce792eSToby Isaac   Input Parameters:
322dce792eSToby Isaac + sp           - the dual space object
332dce792eSToby Isaac - numSumSpaces - the number of spaces
342dce792eSToby Isaac 
352dce792eSToby Isaac   Level: intermediate
362dce792eSToby Isaac 
372dce792eSToby Isaac   Note:
382dce792eSToby Isaac   The name NumSubspaces is slightly misleading because it is actually setting the number of defining spaces of the sum, not a number of Subspaces of it
392dce792eSToby Isaac 
402dce792eSToby Isaac .seealso: `PETSCDUALSPACESUM`, `PetscDualSpace`, `PetscDualSpaceSumGetNumSubspaces()`
412dce792eSToby Isaac @*/
422dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumSetNumSubspaces(PetscDualSpace sp, PetscInt numSumSpaces)
432dce792eSToby Isaac {
442dce792eSToby Isaac   PetscFunctionBegin;
452dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
462dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumSetNumSubspaces_C", (PetscDualSpace, PetscInt), (sp, numSumSpaces));
472dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
482dce792eSToby Isaac }
492dce792eSToby Isaac 
502dce792eSToby Isaac /*@
512dce792eSToby Isaac   PetscDualSpaceSumGetConcatenate - Get the concatenate flag for this space.
522dce792eSToby Isaac 
532dce792eSToby Isaac   Input Parameter:
542dce792eSToby Isaac . sp - the dual space object
552dce792eSToby Isaac 
562dce792eSToby Isaac   Output Parameter:
572dce792eSToby Isaac . concatenate - flag indicating whether subspaces are concatenated.
582dce792eSToby Isaac 
592dce792eSToby Isaac   Level: intermediate
602dce792eSToby Isaac 
612dce792eSToby Isaac   Note:
622dce792eSToby Isaac   A concatenated sum space will have the number of components equal to the sum of the number of
632dce792eSToby Isaac   components of all subspaces. A non-concatenated, or direct sum space will have the same
642dce792eSToby Isaac   number of components as its subspaces.
652dce792eSToby Isaac 
662dce792eSToby Isaac .seealso: `PETSCDUALSPACESUM`, `PetscDualSpace`, `PetscDualSpaceSumSetConcatenate()`
672dce792eSToby Isaac @*/
682dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumGetConcatenate(PetscDualSpace sp, PetscBool *concatenate)
692dce792eSToby Isaac {
702dce792eSToby Isaac   PetscFunctionBegin;
712dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
722dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumGetConcatenate_C", (PetscDualSpace, PetscBool *), (sp, concatenate));
732dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
742dce792eSToby Isaac }
752dce792eSToby Isaac 
762dce792eSToby Isaac /*@
772dce792eSToby Isaac   PetscDualSpaceSumSetConcatenate - Sets the concatenate flag for this space.
782dce792eSToby Isaac 
792dce792eSToby Isaac   Input Parameters:
802dce792eSToby Isaac + sp          - the dual space object
812dce792eSToby Isaac - concatenate - are subspaces concatenated components (true) or direct summands (false)
822dce792eSToby Isaac 
832dce792eSToby Isaac   Level: intermediate
842dce792eSToby Isaac 
852dce792eSToby Isaac   Notes:
862dce792eSToby Isaac   A concatenated sum space will have the number of components equal to the sum of the number of
872dce792eSToby Isaac   components of all subspaces. A non-concatenated, or direct sum space will have the same
882dce792eSToby Isaac   number of components as its subspaces .
892dce792eSToby Isaac 
902dce792eSToby Isaac .seealso: `PETSCDUALSPACESUM`, `PetscDualSpace`, `PetscDualSpaceSumGetConcatenate()`
912dce792eSToby Isaac @*/
922dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumSetConcatenate(PetscDualSpace sp, PetscBool concatenate)
932dce792eSToby Isaac {
942dce792eSToby Isaac   PetscFunctionBegin;
952dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
962dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumSetConcatenate_C", (PetscDualSpace, PetscBool), (sp, concatenate));
972dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
982dce792eSToby Isaac }
992dce792eSToby Isaac 
1002dce792eSToby Isaac /*@
1012dce792eSToby Isaac   PetscDualSpaceSumGetSubspace - Get a space in the sum space
1022dce792eSToby Isaac 
1032dce792eSToby Isaac   Input Parameters:
1042dce792eSToby Isaac + sp - the dual space object
1052dce792eSToby Isaac - s  - The space number
1062dce792eSToby Isaac 
1072dce792eSToby Isaac   Output Parameter:
1082dce792eSToby Isaac . subsp - the `PetscDualSpace`
1092dce792eSToby Isaac 
1102dce792eSToby Isaac   Level: intermediate
1112dce792eSToby Isaac 
1122dce792eSToby Isaac   Note:
1132dce792eSToby Isaac   The name GetSubspace is slightly misleading because it is actually getting one of the defining spaces of the sum, not a Subspace of it
1142dce792eSToby Isaac 
1152dce792eSToby Isaac .seealso: `PETSCDUALSPACESUM`, `PetscDualSpace`, `PetscDualSpaceSumSetSubspace()`
1162dce792eSToby Isaac @*/
1172dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumGetSubspace(PetscDualSpace sp, PetscInt s, PetscDualSpace *subsp)
1182dce792eSToby Isaac {
1192dce792eSToby Isaac   PetscFunctionBegin;
1202dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1212dce792eSToby Isaac   PetscAssertPointer(subsp, 3);
1222dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumGetSubspace_C", (PetscDualSpace, PetscInt, PetscDualSpace *), (sp, s, subsp));
1232dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1242dce792eSToby Isaac }
1252dce792eSToby Isaac 
1262dce792eSToby Isaac /*@
1272dce792eSToby Isaac   PetscDualSpaceSumSetSubspace - Set a space in the sum space
1282dce792eSToby Isaac 
1292dce792eSToby Isaac   Input Parameters:
1302dce792eSToby Isaac + sp    - the dual space object
1312dce792eSToby Isaac . s     - The space number
1322dce792eSToby Isaac - subsp - the number of spaces
1332dce792eSToby Isaac 
1342dce792eSToby Isaac   Level: intermediate
1352dce792eSToby Isaac 
1362dce792eSToby Isaac   Note:
1372dce792eSToby Isaac   The name SetSubspace is slightly misleading because it is actually setting one of the defining spaces of the sum, not a Subspace of it
1382dce792eSToby Isaac 
1392dce792eSToby Isaac .seealso: `PETSCDUALSPACESUM`, `PetscDualSpace`, `PetscDualSpaceSumGetSubspace()`
1402dce792eSToby Isaac @*/
1412dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumSetSubspace(PetscDualSpace sp, PetscInt s, PetscDualSpace subsp)
1422dce792eSToby Isaac {
1432dce792eSToby Isaac   PetscFunctionBegin;
1442dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
1452dce792eSToby Isaac   if (subsp) PetscValidHeaderSpecific(subsp, PETSCDUALSPACE_CLASSID, 3);
1462dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumSetSubspace_C", (PetscDualSpace, PetscInt, PetscDualSpace), (sp, s, subsp));
1472dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1482dce792eSToby Isaac }
1492dce792eSToby Isaac 
1502dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumGetNumSubspaces_Sum(PetscDualSpace space, PetscInt *numSumSpaces)
1512dce792eSToby Isaac {
1522dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)space->data;
1532dce792eSToby Isaac 
1542dce792eSToby Isaac   PetscFunctionBegin;
1552dce792eSToby Isaac   *numSumSpaces = sum->numSumSpaces;
1562dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1572dce792eSToby Isaac }
1582dce792eSToby Isaac 
1592dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumSetNumSubspaces_Sum(PetscDualSpace space, PetscInt numSumSpaces)
1602dce792eSToby Isaac {
1612dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)space->data;
1622dce792eSToby Isaac   PetscInt            Ns  = sum->numSumSpaces;
1632dce792eSToby Isaac 
1642dce792eSToby Isaac   PetscFunctionBegin;
1652dce792eSToby Isaac   PetscCheck(!sum->setupCalled, PetscObjectComm((PetscObject)space), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of subspaces after setup called");
1662dce792eSToby Isaac   if (numSumSpaces == Ns) PetscFunctionReturn(PETSC_SUCCESS);
1672dce792eSToby Isaac   for (PetscInt s = 0; s < Ns; ++s) PetscCall(PetscDualSpaceDestroy(&sum->sumspaces[s]));
1682dce792eSToby Isaac   PetscCall(PetscFree(sum->sumspaces));
1692dce792eSToby Isaac 
1702dce792eSToby Isaac   Ns = sum->numSumSpaces = numSumSpaces;
1712dce792eSToby Isaac   PetscCall(PetscCalloc1(Ns, &sum->sumspaces));
1722dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1732dce792eSToby Isaac }
1742dce792eSToby Isaac 
1752dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumGetConcatenate_Sum(PetscDualSpace sp, PetscBool *concatenate)
1762dce792eSToby Isaac {
1772dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data;
1782dce792eSToby Isaac 
1792dce792eSToby Isaac   PetscFunctionBegin;
1802dce792eSToby Isaac   *concatenate = sum->concatenate;
1812dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1822dce792eSToby Isaac }
1832dce792eSToby Isaac 
1842dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumSetConcatenate_Sum(PetscDualSpace sp, PetscBool concatenate)
1852dce792eSToby Isaac {
1862dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data;
1872dce792eSToby Isaac 
1882dce792eSToby Isaac   PetscFunctionBegin;
1892dce792eSToby Isaac   PetscCheck(!sum->setupCalled, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change space concatenation after setup called.");
1902dce792eSToby Isaac 
1912dce792eSToby Isaac   sum->concatenate = concatenate;
1922dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
1932dce792eSToby Isaac }
1942dce792eSToby Isaac 
1952dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumGetSubspace_Sum(PetscDualSpace space, PetscInt s, PetscDualSpace *subspace)
1962dce792eSToby Isaac {
1972dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)space->data;
1982dce792eSToby Isaac   PetscInt            Ns  = sum->numSumSpaces;
1992dce792eSToby Isaac 
2002dce792eSToby Isaac   PetscFunctionBegin;
2012dce792eSToby Isaac   PetscCheck(Ns >= 0, PetscObjectComm((PetscObject)space), PETSC_ERR_ARG_WRONGSTATE, "Must call PetscDualSpaceSumSetNumSubspaces() first");
2022dce792eSToby Isaac   PetscCheck(s >= 0 && s < Ns, PetscObjectComm((PetscObject)space), PETSC_ERR_ARG_OUTOFRANGE, "Invalid subspace number %" PetscInt_FMT, s);
2032dce792eSToby Isaac 
2042dce792eSToby Isaac   *subspace = sum->sumspaces[s];
2052dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
2062dce792eSToby Isaac }
2072dce792eSToby Isaac 
2082dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumSetSubspace_Sum(PetscDualSpace space, PetscInt s, PetscDualSpace subspace)
2092dce792eSToby Isaac {
2102dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)space->data;
2112dce792eSToby Isaac   PetscInt            Ns  = sum->numSumSpaces;
2122dce792eSToby Isaac 
2132dce792eSToby Isaac   PetscFunctionBegin;
2142dce792eSToby Isaac   PetscCheck(!sum->setupCalled, PetscObjectComm((PetscObject)space), PETSC_ERR_ARG_WRONGSTATE, "Cannot change subspace after setup called");
2152dce792eSToby Isaac   PetscCheck(Ns >= 0, PetscObjectComm((PetscObject)space), PETSC_ERR_ARG_WRONGSTATE, "Must call PetscDualSpaceSumSetNumSubspaces() first");
2162dce792eSToby Isaac   PetscCheck(s >= 0 && s < Ns, PetscObjectComm((PetscObject)space), PETSC_ERR_ARG_OUTOFRANGE, "Invalid subspace number %" PetscInt_FMT, s);
2172dce792eSToby Isaac 
2182dce792eSToby Isaac   PetscCall(PetscObjectReference((PetscObject)subspace));
2192dce792eSToby Isaac   PetscCall(PetscDualSpaceDestroy(&sum->sumspaces[s]));
2202dce792eSToby Isaac   sum->sumspaces[s] = subspace;
2212dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
2222dce792eSToby Isaac }
2232dce792eSToby Isaac 
2242dce792eSToby Isaac static PetscErrorCode PetscDualSpaceDuplicate_Sum(PetscDualSpace sp, PetscDualSpace spNew)
2252dce792eSToby Isaac {
2262dce792eSToby Isaac   PetscInt       num_subspaces, Nc;
2272dce792eSToby Isaac   PetscBool      concatenate, interleave_basis, interleave_components;
2282dce792eSToby Isaac   PetscDualSpace subsp_first     = NULL;
2292dce792eSToby Isaac   PetscDualSpace subsp_dup_first = NULL;
2302dce792eSToby Isaac   DM             K;
2312dce792eSToby Isaac 
2322dce792eSToby Isaac   PetscFunctionBegin;
2332dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetNumSubspaces(sp, &num_subspaces));
2342dce792eSToby Isaac   PetscCall(PetscDualSpaceSumSetNumSubspaces(spNew, num_subspaces));
2352dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetConcatenate(sp, &concatenate));
2362dce792eSToby Isaac   PetscCall(PetscDualSpaceSumSetConcatenate(spNew, concatenate));
2372dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetInterleave(sp, &interleave_basis, &interleave_components));
2382dce792eSToby Isaac   PetscCall(PetscDualSpaceSumSetInterleave(spNew, interleave_basis, interleave_components));
2392dce792eSToby Isaac   PetscCall(PetscDualSpaceGetDM(sp, &K));
2402dce792eSToby Isaac   PetscCall(PetscDualSpaceSetDM(spNew, K));
2412dce792eSToby Isaac   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
2422dce792eSToby Isaac   PetscCall(PetscDualSpaceSetNumComponents(spNew, Nc));
2432dce792eSToby Isaac   for (PetscInt s = 0; s < num_subspaces; s++) {
2442dce792eSToby Isaac     PetscDualSpace subsp, subspNew;
2452dce792eSToby Isaac 
2462dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
2472dce792eSToby Isaac     if (s == 0) {
2482dce792eSToby Isaac       subsp_first = subsp;
2492dce792eSToby Isaac       PetscCall(PetscDualSpaceDuplicate(subsp, &subsp_dup_first));
2502dce792eSToby Isaac       subspNew = subsp_dup_first;
2512dce792eSToby Isaac     } else if (subsp == subsp_first) {
2522dce792eSToby Isaac       PetscCall(PetscObjectReference((PetscObject)subsp_dup_first));
2532dce792eSToby Isaac       subspNew = subsp_dup_first;
2542dce792eSToby Isaac     } else {
2552dce792eSToby Isaac       PetscCall(PetscDualSpaceDuplicate(subsp, &subspNew));
2562dce792eSToby Isaac     }
2572dce792eSToby Isaac     PetscCall(PetscDualSpaceSumSetSubspace(spNew, s, subspNew));
2582dce792eSToby Isaac     PetscCall(PetscDualSpaceDestroy(&subspNew));
2592dce792eSToby Isaac   }
2602dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
2612dce792eSToby Isaac }
2622dce792eSToby Isaac 
2632dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumCreateQuadrature(PetscInt Ns, PetscInt cdim, PetscBool uniform_points, PetscQuadrature subquads[], PetscQuadrature *fullquad)
2642dce792eSToby Isaac {
2652dce792eSToby Isaac   PetscReal *points;
2662dce792eSToby Isaac   PetscInt   Npoints;
2672dce792eSToby Isaac 
2682dce792eSToby Isaac   PetscFunctionBegin;
2692dce792eSToby Isaac   if (uniform_points) {
2702dce792eSToby Isaac     PetscCall(PetscObjectReference((PetscObject)subquads[0]));
2712dce792eSToby Isaac     *fullquad = subquads[0];
2722dce792eSToby Isaac     PetscFunctionReturn(PETSC_SUCCESS);
2732dce792eSToby Isaac   }
2742dce792eSToby Isaac   Npoints = 0;
2752dce792eSToby Isaac   for (PetscInt s = 0; s < Ns; s++) {
2762dce792eSToby Isaac     PetscInt subNpoints;
2772dce792eSToby Isaac 
2782dce792eSToby Isaac     if (!subquads[s]) continue;
2792dce792eSToby Isaac     PetscCall(PetscQuadratureGetData(subquads[s], NULL, NULL, &subNpoints, NULL, NULL));
2802dce792eSToby Isaac     Npoints += subNpoints;
2812dce792eSToby Isaac   }
2822dce792eSToby Isaac   *fullquad = NULL;
2832dce792eSToby Isaac   if (!Npoints) PetscFunctionReturn(PETSC_SUCCESS);
2842dce792eSToby Isaac   PetscCall(PetscMalloc1(Npoints * cdim, &points));
2852dce792eSToby Isaac   for (PetscInt s = 0, offset = 0; s < Ns; s++) {
2862dce792eSToby Isaac     PetscInt         subNpoints;
2872dce792eSToby Isaac     const PetscReal *subpoints;
2882dce792eSToby Isaac 
2892dce792eSToby Isaac     PetscCall(PetscQuadratureGetData(subquads[s], NULL, NULL, &subNpoints, &subpoints, NULL));
2902dce792eSToby Isaac     PetscCall(PetscArraycpy(&points[offset], subpoints, cdim * subNpoints));
2912dce792eSToby Isaac     offset += cdim * subNpoints;
2922dce792eSToby Isaac   }
2932dce792eSToby Isaac   PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, fullquad));
2942dce792eSToby Isaac   PetscCall(PetscQuadratureSetData(*fullquad, cdim, 0, Npoints, points, NULL));
2952dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
2962dce792eSToby Isaac }
2972dce792eSToby Isaac 
2982dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumCreateMatrix(PetscDualSpace sp, Mat submats[], ISLocalToGlobalMapping map_rows[], ISLocalToGlobalMapping map_cols[], PetscQuadrature fullquad, Mat *fullmat)
2992dce792eSToby Isaac {
3002dce792eSToby Isaac   Mat          mat;
3012dce792eSToby Isaac   PetscInt    *i = NULL, *j = NULL;
3022dce792eSToby Isaac   PetscScalar *v = NULL;
3032dce792eSToby Isaac   PetscInt     npoints;
3042dce792eSToby Isaac   PetscInt     nrows, ncols, nnz;
3052dce792eSToby Isaac   PetscInt     Nc, num_subspaces;
3062dce792eSToby Isaac 
3072dce792eSToby Isaac   PetscFunctionBegin;
3082dce792eSToby Isaac   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
3092dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetNumSubspaces(sp, &num_subspaces));
3102dce792eSToby Isaac   nrows = 0;
3112dce792eSToby Isaac   ncols = 0;
3122dce792eSToby Isaac   nnz   = 0;
3132dce792eSToby Isaac   PetscCall(PetscQuadratureGetData(fullquad, NULL, NULL, &npoints, NULL, NULL));
3142dce792eSToby Isaac   ncols = Nc * npoints;
3152dce792eSToby Isaac   for (PetscInt s = 0; s < num_subspaces; s++) {
3162dce792eSToby Isaac     // Get the COO data for each matrix, map the is and js, and append to growing COO data
3172dce792eSToby Isaac     PetscInt        sNrows, sNcols;
3182dce792eSToby Isaac     Mat             smat;
3192dce792eSToby Isaac     const PetscInt *si;
3202dce792eSToby Isaac     const PetscInt *sj;
3212dce792eSToby Isaac     PetscScalar    *sv;
3222dce792eSToby Isaac     PetscMemType    memtype;
3232dce792eSToby Isaac     PetscInt        snz;
3242dce792eSToby Isaac     PetscInt        snz_actual;
3252dce792eSToby Isaac     PetscInt       *cooi;
3262dce792eSToby Isaac     PetscInt       *cooj;
3272dce792eSToby Isaac     PetscScalar    *coov;
3282dce792eSToby Isaac     PetscScalar    *v_new;
3292dce792eSToby Isaac     PetscInt       *i_new;
3302dce792eSToby Isaac     PetscInt       *j_new;
3312dce792eSToby Isaac 
3322dce792eSToby Isaac     if (!submats[s]) continue;
3332dce792eSToby Isaac     PetscCall(MatGetSize(submats[s], &sNrows, &sNcols));
3342dce792eSToby Isaac     nrows += sNrows;
3352dce792eSToby Isaac     PetscCall(MatConvert(submats[s], MATSEQAIJ, MAT_INITIAL_MATRIX, &smat));
3362dce792eSToby Isaac     PetscCall(MatBindToCPU(smat, PETSC_TRUE));
3372dce792eSToby Isaac     PetscCall(MatSeqAIJGetCSRAndMemType(smat, &si, &sj, &sv, &memtype));
3382dce792eSToby Isaac     PetscCheck(memtype == PETSC_MEMTYPE_HOST, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not convert matrix to host memory");
3392dce792eSToby Isaac     snz = si[sNrows];
3402dce792eSToby Isaac 
3412dce792eSToby Isaac     PetscCall(PetscMalloc1(nnz + snz, &v_new));
3422dce792eSToby Isaac     PetscCall(PetscArraycpy(v_new, v, nnz));
3432dce792eSToby Isaac     PetscCall(PetscFree(v));
3442dce792eSToby Isaac     v = v_new;
3452dce792eSToby Isaac 
3462dce792eSToby Isaac     PetscCall(PetscMalloc1(nnz + snz, &i_new));
3472dce792eSToby Isaac     PetscCall(PetscArraycpy(i_new, i, nnz));
3482dce792eSToby Isaac     PetscCall(PetscFree(i));
3492dce792eSToby Isaac     i = i_new;
3502dce792eSToby Isaac 
3512dce792eSToby Isaac     PetscCall(PetscMalloc1(nnz + snz, &j_new));
3522dce792eSToby Isaac     PetscCall(PetscArraycpy(j_new, j, nnz));
3532dce792eSToby Isaac     PetscCall(PetscFree(j));
3542dce792eSToby Isaac     j = j_new;
3552dce792eSToby Isaac 
3562dce792eSToby Isaac     PetscCall(PetscMalloc2(snz, &cooi, snz, &cooj));
3572dce792eSToby Isaac 
3582dce792eSToby Isaac     coov = &v[nnz];
3592dce792eSToby Isaac 
3602dce792eSToby Isaac     snz_actual = 0;
3612dce792eSToby Isaac     for (PetscInt row = 0; row < sNrows; row++) {
3622dce792eSToby Isaac       for (PetscInt k = si[row]; k < si[row + 1]; k++, snz_actual++) {
3632dce792eSToby Isaac         cooi[snz_actual] = row;
3642dce792eSToby Isaac         cooj[snz_actual] = sj[k];
3652dce792eSToby Isaac         coov[snz_actual] = sv[k];
3662dce792eSToby Isaac       }
3672dce792eSToby Isaac     }
3682dce792eSToby Isaac     PetscCall(MatDestroy(&smat));
3692dce792eSToby Isaac 
3702dce792eSToby Isaac     if (snz_actual > 0) {
3712dce792eSToby Isaac       PetscCall(ISLocalToGlobalMappingApply(map_rows[s], snz_actual, cooi, &i[nnz]));
3722dce792eSToby Isaac       PetscCall(ISLocalToGlobalMappingApply(map_cols[s], snz_actual, cooj, &j[nnz]));
3732dce792eSToby Isaac       nnz += snz_actual;
3742dce792eSToby Isaac     }
3752dce792eSToby Isaac     PetscCall(PetscFree2(cooi, cooj));
3762dce792eSToby Isaac   }
3772dce792eSToby Isaac   PetscCall(MatCreate(PETSC_COMM_SELF, fullmat));
3782dce792eSToby Isaac   mat = *fullmat;
3792dce792eSToby Isaac   PetscCall(MatSetSizes(mat, nrows, ncols, nrows, ncols));
3802dce792eSToby Isaac   PetscCall(MatSetType(mat, MATSEQAIJ));
3812dce792eSToby Isaac   PetscCall(MatSetPreallocationCOO(mat, nnz, i, j));
3822dce792eSToby Isaac   PetscCall(MatSetValuesCOO(mat, v, INSERT_VALUES));
3832dce792eSToby Isaac   PetscCall(PetscFree(i));
3842dce792eSToby Isaac   PetscCall(PetscFree(j));
3852dce792eSToby Isaac   PetscCall(PetscFree(v));
3862dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
3872dce792eSToby Isaac }
3882dce792eSToby Isaac 
3892dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumCreateMappings(PetscDualSpace sp, PetscBool interior, PetscBool uniform_points, ISLocalToGlobalMapping map_row[], ISLocalToGlobalMapping map_col[])
3902dce792eSToby Isaac {
3912dce792eSToby Isaac   PetscSection section;
3922dce792eSToby Isaac   PetscInt     pStart, pEnd, Ns, Nc;
3932dce792eSToby Isaac   PetscInt     num_points = 0;
3942dce792eSToby Isaac   PetscBool    interleave_basis, interleave_components, concatenate;
3952dce792eSToby Isaac   PetscInt    *roffset;
3962dce792eSToby Isaac 
3972dce792eSToby Isaac   PetscFunctionBegin;
3982dce792eSToby Isaac   if (interior) {
3992dce792eSToby Isaac     PetscCall(PetscDualSpaceGetInteriorSection(sp, &section));
4002dce792eSToby Isaac   } else {
4012dce792eSToby Isaac     PetscCall(PetscDualSpaceGetSection(sp, &section));
4022dce792eSToby Isaac   }
4032dce792eSToby Isaac   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
4042dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetNumSubspaces(sp, &Ns));
4052dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetInterleave(sp, &interleave_basis, &interleave_components));
4062dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetConcatenate(sp, &concatenate));
4072dce792eSToby Isaac   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
4082dce792eSToby Isaac   for (PetscInt p = pStart; p < pEnd; p++) {
4092dce792eSToby Isaac     PetscInt dof;
4102dce792eSToby Isaac 
4112dce792eSToby Isaac     PetscCall(PetscSectionGetDof(section, p, &dof));
4122dce792eSToby Isaac     num_points += (dof > 0);
4132dce792eSToby Isaac   }
4142dce792eSToby Isaac   PetscCall(PetscCalloc1(pEnd - pStart, &roffset));
4152dce792eSToby Isaac   for (PetscInt s = 0, coffset = 0; s < Ns; s++) {
4162dce792eSToby Isaac     PetscDualSpace  subsp;
4172dce792eSToby Isaac     PetscSection    subsection;
4182dce792eSToby Isaac     IS              is_row, is_col;
4192dce792eSToby Isaac     PetscInt        subNb;
4202dce792eSToby Isaac     PetscInt        sNc, sNpoints, sNcols;
4212dce792eSToby Isaac     PetscQuadrature quad;
4222dce792eSToby Isaac 
4232dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
4242dce792eSToby Isaac     PetscCall(PetscDualSpaceGetNumComponents(subsp, &sNc));
4252dce792eSToby Isaac     if (interior) {
4262dce792eSToby Isaac       PetscCall(PetscDualSpaceGetInteriorSection(subsp, &subsection));
4272dce792eSToby Isaac       PetscCall(PetscDualSpaceGetInteriorData(subsp, &quad, NULL));
4282dce792eSToby Isaac     } else {
4292dce792eSToby Isaac       PetscCall(PetscDualSpaceGetSection(subsp, &subsection));
4302dce792eSToby Isaac       PetscCall(PetscDualSpaceGetAllData(subsp, &quad, NULL));
4312dce792eSToby Isaac     }
4322dce792eSToby Isaac     PetscCall(PetscSectionGetStorageSize(subsection, &subNb));
4332dce792eSToby Isaac     if (num_points == 1) {
4342dce792eSToby Isaac       PetscInt rstride;
4352dce792eSToby Isaac 
4362dce792eSToby Isaac       rstride = interleave_basis ? Ns : 1;
4372dce792eSToby Isaac       PetscCall(ISCreateStride(PETSC_COMM_SELF, subNb, roffset[0], rstride, &is_row));
4382dce792eSToby Isaac       roffset[0] += interleave_basis ? 1 : subNb;
4392dce792eSToby Isaac     } else {
4402dce792eSToby Isaac       PetscInt *rows;
4412dce792eSToby Isaac 
4422dce792eSToby Isaac       PetscCall(PetscMalloc1(subNb, &rows));
4432dce792eSToby Isaac       for (PetscInt p = pStart; p < pEnd; p++) {
4442dce792eSToby Isaac         PetscInt subdof, dof, off, suboff, stride;
4452dce792eSToby Isaac 
4462dce792eSToby Isaac         PetscCall(PetscSectionGetOffset(subsection, p, &suboff));
4472dce792eSToby Isaac         PetscCall(PetscSectionGetDof(subsection, p, &subdof));
4482dce792eSToby Isaac         PetscCall(PetscSectionGetOffset(section, p, &off));
4492dce792eSToby Isaac         PetscCall(PetscSectionGetDof(section, p, &dof));
4502dce792eSToby Isaac         PetscCheck(subdof * Ns == dof || !interleave_basis, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Basis cannot be interleaved\n");
4512dce792eSToby Isaac         stride = interleave_basis ? Ns : 1;
4522dce792eSToby Isaac         for (PetscInt k = 0; k < subdof; k++) { rows[suboff + k] = off + roffset[p - pStart] + k * stride; }
4532dce792eSToby Isaac         roffset[p - pStart] += interleave_basis ? 1 : subdof;
4542dce792eSToby Isaac       }
4552dce792eSToby Isaac       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, subNb, rows, PETSC_OWN_POINTER, &is_row));
4562dce792eSToby Isaac     }
4572dce792eSToby Isaac 
4582dce792eSToby Isaac     sNpoints = 0;
4592dce792eSToby Isaac     if (quad) PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &sNpoints, NULL, NULL));
4602dce792eSToby Isaac     sNcols = sNpoints * sNc;
4612dce792eSToby Isaac 
4622dce792eSToby Isaac     if (!concatenate) {
4632dce792eSToby Isaac       PetscCall(ISCreateStride(PETSC_COMM_SELF, sNcols, coffset, 1, &is_col));
4642dce792eSToby Isaac       coffset += uniform_points ? 0 : sNcols;
4652dce792eSToby Isaac     } else {
4662dce792eSToby Isaac       if (uniform_points && interleave_components) {
4672dce792eSToby Isaac         PetscCall(ISCreateStride(PETSC_COMM_SELF, sNcols, coffset, Ns, &is_col));
4682dce792eSToby Isaac         coffset += 1;
4692dce792eSToby Isaac       } else {
4702dce792eSToby Isaac         PetscInt *cols;
4712dce792eSToby Isaac 
4722dce792eSToby Isaac         PetscCall(PetscMalloc1(sNcols, &cols));
4732dce792eSToby Isaac         for (PetscInt p = 0, r = 0; p < sNpoints; p++) {
4742dce792eSToby Isaac           for (PetscInt c = 0; c < sNc; c++, r++) { cols[r] = coffset + p * Nc + c; }
4752dce792eSToby Isaac         }
4762dce792eSToby Isaac         coffset += uniform_points ? sNc : Nc * sNpoints + sNc;
4772dce792eSToby Isaac         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, sNcols, cols, PETSC_OWN_POINTER, &is_col));
4782dce792eSToby Isaac       }
4792dce792eSToby Isaac     }
4802dce792eSToby Isaac     PetscCall(ISLocalToGlobalMappingCreateIS(is_row, &map_row[s]));
4812dce792eSToby Isaac     PetscCall(ISLocalToGlobalMappingCreateIS(is_col, &map_col[s]));
4822dce792eSToby Isaac     PetscCall(ISDestroy(&is_row));
4832dce792eSToby Isaac     PetscCall(ISDestroy(&is_col));
4842dce792eSToby Isaac   }
4852dce792eSToby Isaac   PetscCall(PetscFree(roffset));
4862dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
4872dce792eSToby Isaac }
4882dce792eSToby Isaac 
4892dce792eSToby Isaac static PetscErrorCode PetscDualSpaceCreateFacetSubspace_Sum(PetscDualSpace sp, PetscInt f, PetscDualSpace *bdsp)
4902dce792eSToby Isaac {
4912dce792eSToby Isaac   PetscInt       k, Nc, Nk, Nknew, Ncnew, Ns;
4922dce792eSToby Isaac   PetscInt       dim, pointDim = -1;
4932dce792eSToby Isaac   PetscInt       depth, Ncopies;
4942dce792eSToby Isaac   PetscBool      any;
4952dce792eSToby Isaac   DM             dm, K;
4962dce792eSToby Isaac   DMPolytopeType ct;
4972dce792eSToby Isaac 
4982dce792eSToby Isaac   PetscFunctionBegin;
4992dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetNumSubspaces(sp, &Ns));
5002dce792eSToby Isaac   any = PETSC_FALSE;
5012dce792eSToby Isaac   for (PetscInt s = 0; s < Ns; s++) {
5022dce792eSToby Isaac     PetscDualSpace subsp, subpointsp;
5032dce792eSToby Isaac 
5042dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
5052dce792eSToby Isaac     PetscCall(PetscDualSpaceGetPointSubspace(subsp, f, &subpointsp));
5062dce792eSToby Isaac     if (subpointsp) any = PETSC_TRUE;
5072dce792eSToby Isaac   }
5082dce792eSToby Isaac   if (!any) {
5092dce792eSToby Isaac     *bdsp = NULL;
5102dce792eSToby Isaac     PetscFunctionReturn(PETSC_SUCCESS);
5112dce792eSToby Isaac   }
5122dce792eSToby Isaac   PetscCall(PetscDualSpaceGetDM(sp, &dm));
5132dce792eSToby Isaac   PetscCall(DMGetDimension(dm, &dim));
5142dce792eSToby Isaac   PetscCall(DMPlexGetDepth(dm, &depth));
5152dce792eSToby Isaac   PetscCall(PetscDualSpaceDuplicate(sp, bdsp));
5162dce792eSToby Isaac   PetscCheck((depth == dim) || (depth == 1), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported interpolation state of reference element");
5172dce792eSToby Isaac   pointDim = (depth == dim) ? (dim - 1) : 0;
5182dce792eSToby Isaac   PetscCall(DMPlexGetCellType(dm, f, &ct));
5192dce792eSToby Isaac   PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &K));
5202dce792eSToby Isaac   PetscCall(PetscDualSpaceSetDM(*bdsp, K));
5212dce792eSToby Isaac   PetscCall(PetscDualSpaceGetFormDegree(sp, &k));
5222dce792eSToby Isaac   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
5232dce792eSToby Isaac   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk));
5242dce792eSToby Isaac   Ncopies = Nc / Nk;
5252dce792eSToby Isaac   PetscCall(PetscDTBinomialInt(pointDim, PetscAbsInt(k), &Nknew));
5262dce792eSToby Isaac   Ncnew = Nknew * Ncopies;
5272dce792eSToby Isaac   PetscCall(PetscDualSpaceSetNumComponents(*bdsp, Ncnew));
5282dce792eSToby Isaac   for (PetscInt s = 0; s < Ns; s++) {
5292dce792eSToby Isaac     PetscDualSpace subsp, subpointsp;
5302dce792eSToby Isaac 
5312dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
5322dce792eSToby Isaac     PetscCall(PetscDualSpaceGetPointSubspace(subsp, f, &subpointsp));
5332dce792eSToby Isaac     if (subpointsp) {
5342dce792eSToby Isaac       PetscCall(PetscObjectReference((PetscObject)subpointsp));
5352dce792eSToby Isaac     } else {
5362dce792eSToby Isaac       // make an empty dual space
5372dce792eSToby Isaac       PetscInt subNc, subNcopies, subpointNc;
5382dce792eSToby Isaac 
5392dce792eSToby Isaac       PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)subsp), &subpointsp));
5402dce792eSToby Isaac       PetscCall(PetscDualSpaceGetNumComponents(subsp, &subNc));
5412dce792eSToby Isaac       subNcopies = subNc / Nk;
5422dce792eSToby Isaac       subpointNc = subNcopies * Nknew;
5432dce792eSToby Isaac       PetscCall(PetscDualSpaceSetType(subpointsp, PETSCDUALSPACESIMPLE));
5442dce792eSToby Isaac       PetscCall(PetscDualSpaceSimpleSetDimension(subpointsp, 0));
5452dce792eSToby Isaac       PetscCall(PetscDualSpaceSetFormDegree(subpointsp, k));
5462dce792eSToby Isaac       PetscCall(PetscDualSpaceSetNumComponents(subpointsp, subpointNc));
5472dce792eSToby Isaac     }
5482dce792eSToby Isaac     // K should be equal to subpointsp->dm, but we want it to be exactly the same
5492dce792eSToby Isaac     PetscCall(PetscObjectReference((PetscObject)K));
5502dce792eSToby Isaac     PetscCall(DMDestroy(&subpointsp->dm));
5512dce792eSToby Isaac     subpointsp->dm = K;
5522dce792eSToby Isaac     PetscCall(PetscDualSpaceSetUp(subpointsp));
5532dce792eSToby Isaac     PetscCall(PetscDualSpaceSumSetSubspace(*bdsp, s, subpointsp));
5542dce792eSToby Isaac     PetscCall(PetscDualSpaceDestroy(&subpointsp));
5552dce792eSToby Isaac   }
5562dce792eSToby Isaac   PetscCall(DMDestroy(&K));
5572dce792eSToby Isaac   PetscCall(PetscDualSpaceSetUp(*bdsp));
5582dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
5592dce792eSToby Isaac }
5602dce792eSToby Isaac 
5612dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumIsUniform(PetscDualSpace sp, PetscBool *is_uniform)
5622dce792eSToby Isaac {
5632dce792eSToby Isaac   PetscDualSpace_Sum *sum     = (PetscDualSpace_Sum *)sp->data;
5642dce792eSToby Isaac   PetscBool           uniform = PETSC_TRUE;
5652dce792eSToby Isaac 
5662dce792eSToby Isaac   PetscFunctionBegin;
5672dce792eSToby Isaac   for (PetscInt s = 1; s < sum->numSumSpaces; s++) {
5682dce792eSToby Isaac     if (sum->sumspaces[s] != sum->sumspaces[0]) {
5692dce792eSToby Isaac       uniform = PETSC_FALSE;
5702dce792eSToby Isaac       break;
5712dce792eSToby Isaac     }
5722dce792eSToby Isaac   }
5732dce792eSToby Isaac   *is_uniform = uniform;
5742dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
5752dce792eSToby Isaac }
5762dce792eSToby Isaac 
5772dce792eSToby Isaac static PetscErrorCode PetscDualSpaceGetSymmetries_Sum(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
5782dce792eSToby Isaac {
5792dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data;
5802dce792eSToby Isaac 
5812dce792eSToby Isaac   PetscFunctionBegin;
5822dce792eSToby Isaac   if (!sum->symComputed) {
5832dce792eSToby Isaac     PetscInt       Ns;
5842dce792eSToby Isaac     PetscBool      any_perms = PETSC_FALSE;
5852dce792eSToby Isaac     PetscBool      any_flips = PETSC_FALSE;
5862dce792eSToby Isaac     PetscInt    ***symperms  = NULL;
5872dce792eSToby Isaac     PetscScalar ***symflips  = NULL;
5882dce792eSToby Isaac 
5892dce792eSToby Isaac     sum->symComputed = PETSC_TRUE;
5902dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetNumSubspaces(sp, &Ns));
5912dce792eSToby Isaac     for (PetscInt s = 0; s < Ns; s++) {
5922dce792eSToby Isaac       PetscDualSpace       subsp;
5932dce792eSToby Isaac       const PetscInt    ***sub_perms;
5942dce792eSToby Isaac       const PetscScalar ***sub_flips;
5952dce792eSToby Isaac 
5962dce792eSToby Isaac       PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
5972dce792eSToby Isaac       PetscCall(PetscDualSpaceGetSymmetries(subsp, &sub_perms, &sub_flips));
5982dce792eSToby Isaac       if (sub_perms) any_perms = PETSC_TRUE;
5992dce792eSToby Isaac       if (sub_flips) any_flips = PETSC_TRUE;
6002dce792eSToby Isaac     }
6012dce792eSToby Isaac     if (any_perms || any_flips) {
6022dce792eSToby Isaac       DM       K;
6032dce792eSToby Isaac       PetscInt pStart, pEnd, numPoints;
6042dce792eSToby Isaac       PetscInt spintdim;
6052dce792eSToby Isaac 
6062dce792eSToby Isaac       PetscCall(PetscDualSpaceGetDM(sp, &K));
6072dce792eSToby Isaac       PetscCall(DMPlexGetChart(K, &pStart, &pEnd));
6082dce792eSToby Isaac       numPoints = pEnd - pStart;
6092dce792eSToby Isaac       PetscCall(PetscCalloc1(numPoints, &symperms));
6102dce792eSToby Isaac       PetscCall(PetscCalloc1(numPoints, &symflips));
6112dce792eSToby Isaac       PetscCall(PetscDualSpaceGetBoundarySymmetries_Internal(sp, symperms, symflips));
6122dce792eSToby Isaac       // get interior symmetries
6132dce792eSToby Isaac       PetscCall(PetscDualSpaceGetInteriorDimension(sp, &spintdim));
6142dce792eSToby Isaac       if (spintdim) {
6152dce792eSToby Isaac         PetscInt       groupSize;
6162dce792eSToby Isaac         PetscInt     **cellPerms;
6172dce792eSToby Isaac         PetscScalar  **cellFlips;
6182dce792eSToby Isaac         DMPolytopeType ct;
6192dce792eSToby Isaac 
6202dce792eSToby Isaac         PetscCall(DMPlexGetCellType(K, 0, &ct));
62185036b15SMatthew G. Knepley         groupSize       = DMPolytopeTypeGetNumArrangements(ct);
6222dce792eSToby Isaac         sum->numSelfSym = groupSize;
6232dce792eSToby Isaac         sum->selfSymOff = groupSize / 2;
6242dce792eSToby Isaac         PetscCall(PetscCalloc1(groupSize, &cellPerms));
6252dce792eSToby Isaac         PetscCall(PetscCalloc1(groupSize, &cellFlips));
6262dce792eSToby Isaac         symperms[0] = &cellPerms[groupSize / 2];
6272dce792eSToby Isaac         symflips[0] = &cellFlips[groupSize / 2];
6282dce792eSToby Isaac         for (PetscInt o = -groupSize / 2; o < groupSize / 2; o++) {
6292dce792eSToby Isaac           PetscBool any_o_perms = PETSC_FALSE;
6302dce792eSToby Isaac           PetscBool any_o_flips = PETSC_FALSE;
6312dce792eSToby Isaac 
6322dce792eSToby Isaac           for (PetscInt s = 0; s < Ns; s++) {
6332dce792eSToby Isaac             PetscDualSpace       subsp;
6342dce792eSToby Isaac             const PetscInt    ***sub_perms;
6352dce792eSToby Isaac             const PetscScalar ***sub_flips;
6362dce792eSToby Isaac 
6372dce792eSToby Isaac             PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
6382dce792eSToby Isaac             PetscCall(PetscDualSpaceGetSymmetries(subsp, &sub_perms, &sub_flips));
6392dce792eSToby Isaac             if (sub_perms && sub_perms[0] && sub_perms[0][o]) any_o_perms = PETSC_TRUE;
6402dce792eSToby Isaac             if (sub_flips && sub_flips[0] && sub_flips[0][o]) any_o_flips = PETSC_TRUE;
6412dce792eSToby Isaac           }
6422dce792eSToby Isaac           if (any_o_perms) {
6432dce792eSToby Isaac             PetscInt *o_perm;
6442dce792eSToby Isaac 
6452dce792eSToby Isaac             PetscCall(PetscMalloc1(spintdim, &o_perm));
6462dce792eSToby Isaac             for (PetscInt i = 0; i < spintdim; i++) o_perm[i] = i;
6472dce792eSToby Isaac             for (PetscInt s = 0; s < Ns; s++) {
6482dce792eSToby Isaac               PetscDualSpace       subsp;
6492dce792eSToby Isaac               const PetscInt    ***sub_perms;
6502dce792eSToby Isaac               const PetscScalar ***sub_flips;
6512dce792eSToby Isaac 
6522dce792eSToby Isaac               PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
6532dce792eSToby Isaac               PetscCall(PetscDualSpaceGetSymmetries(subsp, &sub_perms, &sub_flips));
6542dce792eSToby Isaac               if (sub_perms && sub_perms[0] && sub_perms[0][o]) {
6552dce792eSToby Isaac                 PetscInt  subspdim;
6562dce792eSToby Isaac                 PetscInt *range, *domain;
6572dce792eSToby Isaac                 PetscInt *range_mapped, *domain_mapped;
6582dce792eSToby Isaac 
6592dce792eSToby Isaac                 PetscCall(PetscDualSpaceGetInteriorDimension(subsp, &subspdim));
6602dce792eSToby Isaac                 PetscCall(PetscMalloc4(subspdim, &range, subspdim, &range_mapped, subspdim, &domain, subspdim, &domain_mapped));
6612dce792eSToby Isaac                 for (PetscInt i = 0; i < subspdim; i++) domain[i] = i;
6622dce792eSToby Isaac                 PetscCall(PetscArraycpy(range, sub_perms[0][o], subspdim));
6632dce792eSToby Isaac                 PetscCall(ISLocalToGlobalMappingApply(sum->int_rows[s], subspdim, domain, domain_mapped));
6642dce792eSToby Isaac                 PetscCall(ISLocalToGlobalMappingApply(sum->int_rows[s], subspdim, range, range_mapped));
6652dce792eSToby Isaac                 for (PetscInt i = 0; i < subspdim; i++) o_perm[domain_mapped[i]] = range_mapped[i];
6662dce792eSToby Isaac                 PetscCall(PetscFree4(range, range_mapped, domain, domain_mapped));
6672dce792eSToby Isaac               }
6682dce792eSToby Isaac             }
6692dce792eSToby Isaac             symperms[0][o] = o_perm;
6702dce792eSToby Isaac           }
6712dce792eSToby Isaac           if (any_o_flips) {
6722dce792eSToby Isaac             PetscScalar *o_flip;
6732dce792eSToby Isaac 
6742dce792eSToby Isaac             PetscCall(PetscMalloc1(spintdim, &o_flip));
6752dce792eSToby Isaac             for (PetscInt i = 0; i < spintdim; i++) o_flip[i] = 1.0;
6762dce792eSToby Isaac             for (PetscInt s = 0; s < Ns; s++) {
6772dce792eSToby Isaac               PetscDualSpace       subsp;
6782dce792eSToby Isaac               const PetscInt    ***sub_perms;
6792dce792eSToby Isaac               const PetscScalar ***sub_flips;
6802dce792eSToby Isaac 
6812dce792eSToby Isaac               PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
6822dce792eSToby Isaac               PetscCall(PetscDualSpaceGetSymmetries(subsp, &sub_perms, &sub_flips));
6832dce792eSToby Isaac               if (sub_perms && sub_perms[0] && sub_perms[0][o]) {
6842dce792eSToby Isaac                 PetscInt  subspdim;
6852dce792eSToby Isaac                 PetscInt *domain;
6862dce792eSToby Isaac                 PetscInt *domain_mapped;
6872dce792eSToby Isaac 
6882dce792eSToby Isaac                 PetscCall(PetscDualSpaceGetInteriorDimension(subsp, &subspdim));
6892dce792eSToby Isaac                 PetscCall(PetscMalloc2(subspdim, &domain, subspdim, &domain_mapped));
6902dce792eSToby Isaac                 for (PetscInt i = 0; i < subspdim; i++) domain[i] = i;
6912dce792eSToby Isaac                 PetscCall(ISLocalToGlobalMappingApply(sum->int_rows[s], subspdim, domain, domain_mapped));
6922dce792eSToby Isaac                 for (PetscInt i = 0; i < subspdim; i++) o_flip[domain_mapped[i]] = sub_perms[0][o][i];
6932dce792eSToby Isaac                 PetscCall(PetscFree2(domain, domain_mapped));
6942dce792eSToby Isaac               }
6952dce792eSToby Isaac             }
6962dce792eSToby Isaac             symflips[0][o] = o_flip;
6972dce792eSToby Isaac           }
6982dce792eSToby Isaac         }
6992dce792eSToby Isaac         {
7002dce792eSToby Isaac           PetscBool any_perms = PETSC_FALSE;
7012dce792eSToby Isaac           PetscBool any_flips = PETSC_FALSE;
7022dce792eSToby Isaac           for (PetscInt o = -groupSize / 2; o < groupSize / 2; o++) {
7032dce792eSToby Isaac             if (symperms[0][o]) any_perms = PETSC_TRUE;
7042dce792eSToby Isaac             if (symflips[0][o]) any_flips = PETSC_TRUE;
7052dce792eSToby Isaac           }
7062dce792eSToby Isaac           if (!any_perms) {
7072dce792eSToby Isaac             PetscCall(PetscFree(cellPerms));
7082dce792eSToby Isaac             symperms[0] = NULL;
7092dce792eSToby Isaac           }
7102dce792eSToby Isaac           if (!any_flips) {
7112dce792eSToby Isaac             PetscCall(PetscFree(cellFlips));
7122dce792eSToby Isaac             symflips[0] = NULL;
7132dce792eSToby Isaac           }
7142dce792eSToby Isaac         }
7152dce792eSToby Isaac       }
7162dce792eSToby Isaac       if (!any_perms) {
7172dce792eSToby Isaac         PetscCall(PetscFree(symperms));
7182dce792eSToby Isaac         symperms = NULL;
7192dce792eSToby Isaac       }
7202dce792eSToby Isaac       if (!any_flips) {
7212dce792eSToby Isaac         PetscCall(PetscFree(symflips));
7222dce792eSToby Isaac         symflips = NULL;
7232dce792eSToby Isaac       }
7242dce792eSToby Isaac     }
7252dce792eSToby Isaac     sum->symperms = symperms;
7262dce792eSToby Isaac     sum->symflips = symflips;
7272dce792eSToby Isaac   }
7282dce792eSToby Isaac   if (perms) *perms = (const PetscInt ***)sum->symperms;
7292dce792eSToby Isaac   if (flips) *flips = (const PetscScalar ***)sum->symflips;
7302dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
7312dce792eSToby Isaac }
7322dce792eSToby Isaac 
7332dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSetUp_Sum(PetscDualSpace sp)
7342dce792eSToby Isaac {
7352dce792eSToby Isaac   PetscDualSpace_Sum *sum         = (PetscDualSpace_Sum *)sp->data;
7362dce792eSToby Isaac   PetscBool           concatenate = PETSC_TRUE;
7372dce792eSToby Isaac   PetscBool           uniform;
7382dce792eSToby Isaac   PetscInt            Ns, Nc, i, sum_Nc = 0;
7392dce792eSToby Isaac   PetscInt            minNc, maxNc;
7402dce792eSToby Isaac   PetscInt            minForm, maxForm, cdim, depth;
7412dce792eSToby Isaac   DM                  K;
7422dce792eSToby Isaac   PetscQuadrature    *all_quads = NULL;
7432dce792eSToby Isaac   PetscQuadrature    *int_quads = NULL;
7442dce792eSToby Isaac   Mat                *all_mats  = NULL;
7452dce792eSToby Isaac   Mat                *int_mats  = NULL;
7462dce792eSToby Isaac 
7472dce792eSToby Isaac   PetscFunctionBegin;
7482dce792eSToby Isaac   if (sum->setupCalled) PetscFunctionReturn(PETSC_SUCCESS);
7492dce792eSToby Isaac   sum->setupCalled = PETSC_TRUE;
7502dce792eSToby Isaac 
7512dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetNumSubspaces(sp, &Ns));
7522dce792eSToby Isaac   PetscCheck(Ns >= 0, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Cannot have %" PetscInt_FMT " subspaces", Ns);
7532dce792eSToby Isaac 
7542dce792eSToby Isaac   // step 1: make sure they share a DM
7552dce792eSToby Isaac   PetscCall(PetscDualSpaceGetDM(sp, &K));
7562dce792eSToby Isaac   PetscCall(DMGetCoordinateDim(K, &cdim));
7572dce792eSToby Isaac   PetscCall(DMPlexGetDepth(K, &depth));
7582dce792eSToby Isaac   PetscCall(PetscDualSpaceSumIsUniform(sp, &sp->uniform));
7592dce792eSToby Isaac   uniform = sp->uniform;
7602dce792eSToby Isaac   {
7612dce792eSToby Isaac     for (PetscInt s = 0; s < Ns; s++) {
7622dce792eSToby Isaac       PetscDualSpace subsp;
7632dce792eSToby Isaac       DM             sub_K;
7642dce792eSToby Isaac 
7652dce792eSToby Isaac       PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
7662dce792eSToby Isaac       PetscCall(PetscDualSpaceSetUp(subsp));
7672dce792eSToby Isaac       PetscCall(PetscDualSpaceGetDM(subsp, &sub_K));
7682dce792eSToby Isaac       if (s == 0 && K == NULL) {
7692dce792eSToby Isaac         PetscCall(PetscDualSpaceSetDM(sp, sub_K));
7702dce792eSToby Isaac         K = sub_K;
7712dce792eSToby Isaac       }
7722dce792eSToby Isaac       PetscCheck(sub_K == K, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Subspace %d does not have the same DM as the sum space", (int)s);
7732dce792eSToby Isaac     }
7742dce792eSToby Isaac   }
7752dce792eSToby Isaac 
7762dce792eSToby Isaac   // step 2: count components
7772dce792eSToby Isaac   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
7782dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetConcatenate(sp, &concatenate));
7792dce792eSToby Isaac   minNc   = Nc;
7802dce792eSToby Isaac   maxNc   = Nc;
7812dce792eSToby Isaac   minForm = sp->k == PETSC_FORM_DEGREE_UNDEFINED ? PETSC_INT_MAX : sp->k;
7822dce792eSToby Isaac   maxForm = sp->k == PETSC_FORM_DEGREE_UNDEFINED ? PETSC_INT_MIN : sp->k;
7832dce792eSToby Isaac   for (i = 0; i < Ns; ++i) {
7842dce792eSToby Isaac     PetscInt       sNc, formDegree;
7852dce792eSToby Isaac     PetscDualSpace si;
7862dce792eSToby Isaac 
7872dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetSubspace(sp, i, &si));
7882dce792eSToby Isaac     PetscCall(PetscDualSpaceSetUp(si));
7892dce792eSToby Isaac     PetscCall(PetscDualSpaceGetNumComponents(si, &sNc));
7902dce792eSToby Isaac     if (sNc != Nc) PetscCheck(concatenate, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Subspace as a different number of components but space does not concatenate components");
7912dce792eSToby Isaac     minNc = PetscMin(minNc, sNc);
7922dce792eSToby Isaac     maxNc = PetscMax(maxNc, sNc);
7932dce792eSToby Isaac     sum_Nc += sNc;
7942dce792eSToby Isaac     PetscCall(PetscDualSpaceGetFormDegree(si, &formDegree));
7952dce792eSToby Isaac     minForm = PetscMin(minForm, formDegree);
7962dce792eSToby Isaac     maxForm = PetscMax(maxForm, formDegree);
7972dce792eSToby Isaac   }
7982dce792eSToby Isaac   sp->k = (minForm != maxForm) ? PETSC_FORM_DEGREE_UNDEFINED : minForm;
7992dce792eSToby Isaac 
8002dce792eSToby Isaac   if (concatenate) PetscCheck(sum_Nc == Nc, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Total number of subspace components (%" PetscInt_FMT ") does not match number of target space components (%" PetscInt_FMT ").", sum_Nc, Nc);
8012dce792eSToby Isaac   else PetscCheck(minNc == Nc && maxNc == Nc, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Subspaces must have same number of components as the target space.");
8022dce792eSToby Isaac 
8032dce792eSToby Isaac   /* A PetscDualSpace should have a fixed number of components, but
8042dce792eSToby Isaac      if the spaces we are combining have different form degrees, they will not
8052dce792eSToby Isaac      have the same number of components on subcomponents of the boundary,
8062dce792eSToby Isaac      so we do not try to create boundary dual spaces in this case */
8072dce792eSToby Isaac   if (sp->k != PETSC_FORM_DEGREE_UNDEFINED && depth > 0) {
8082dce792eSToby Isaac     PetscInt  pStart, pEnd;
8092dce792eSToby Isaac     PetscInt *pStratStart, *pStratEnd;
8102dce792eSToby Isaac 
8112dce792eSToby Isaac     PetscCall(DMPlexGetDepth(K, &depth));
8122dce792eSToby Isaac     PetscCall(DMPlexGetChart(K, &pStart, &pEnd));
813*f4f49eeaSPierre Jolivet     PetscCall(PetscCalloc1(pEnd, &sp->pointSpaces));
8142dce792eSToby Isaac     PetscCall(PetscMalloc2(depth + 1, &pStratStart, depth + 1, &pStratEnd));
8152dce792eSToby Isaac     for (PetscInt d = 0; d <= depth; ++d) PetscCall(DMPlexGetDepthStratum(K, d, &pStratStart[d], &pStratEnd[d]));
8162dce792eSToby Isaac 
8172dce792eSToby Isaac     for (PetscInt p = pStratStart[depth - 1]; p < pStratEnd[depth - 1]; p++) { /* calculate the facet dual spaces */
8182dce792eSToby Isaac       PetscReal      v0[3];
8192dce792eSToby Isaac       DMPolytopeType ptype;
8202dce792eSToby Isaac       PetscReal      J[9], detJ;
8212dce792eSToby Isaac       PetscInt       q;
8222dce792eSToby Isaac 
8232dce792eSToby Isaac       PetscCall(DMPlexComputeCellGeometryAffineFEM(K, p, v0, J, NULL, &detJ));
8242dce792eSToby Isaac       PetscCall(DMPlexGetCellType(K, p, &ptype));
8252dce792eSToby Isaac 
8262dce792eSToby Isaac       /* compare to previous facets: if computed, reference that dualspace */
8272dce792eSToby Isaac       for (q = pStratStart[depth - 1]; q < p; q++) {
8282dce792eSToby Isaac         DMPolytopeType qtype;
8292dce792eSToby Isaac 
8302dce792eSToby Isaac         PetscCall(DMPlexGetCellType(K, q, &qtype));
8312dce792eSToby Isaac         if (qtype == ptype) break;
8322dce792eSToby Isaac       }
8332dce792eSToby Isaac       if (q < p) { /* this facet has the same dual space as that one */
8342dce792eSToby Isaac         PetscCall(PetscObjectReference((PetscObject)sp->pointSpaces[q]));
8352dce792eSToby Isaac         sp->pointSpaces[p] = sp->pointSpaces[q];
8362dce792eSToby Isaac         continue;
8372dce792eSToby Isaac       }
8382dce792eSToby Isaac       /* if not, recursively compute this dual space */
8392dce792eSToby Isaac       PetscCall(PetscDualSpaceCreateFacetSubspace_Sum(sp, p, &sp->pointSpaces[p]));
8402dce792eSToby Isaac     }
8412dce792eSToby Isaac     for (PetscInt h = 2; h <= depth; h++) { /* get the higher subspaces from the facet subspaces */
8422dce792eSToby Isaac       PetscInt hd = depth - h;
8432dce792eSToby Isaac 
8442dce792eSToby Isaac       for (PetscInt p = pStratStart[hd]; p < pStratEnd[hd]; p++) {
8452dce792eSToby Isaac         PetscInt        suppSize;
8462dce792eSToby Isaac         const PetscInt *supp;
8472dce792eSToby Isaac         DM              qdm;
8482dce792eSToby Isaac         PetscDualSpace  qsp, psp;
8492dce792eSToby Isaac         PetscInt        c, coneSize, q;
8502dce792eSToby Isaac         const PetscInt *cone;
8512dce792eSToby Isaac         const PetscInt *refCone;
8522dce792eSToby Isaac 
8532dce792eSToby Isaac         PetscCall(DMPlexGetSupportSize(K, p, &suppSize));
8542dce792eSToby Isaac         PetscCall(DMPlexGetSupport(K, p, &supp));
8552dce792eSToby Isaac         q   = supp[0];
8562dce792eSToby Isaac         qsp = sp->pointSpaces[q];
8572dce792eSToby Isaac         PetscCall(DMPlexGetConeSize(K, q, &coneSize));
8582dce792eSToby Isaac         PetscCall(DMPlexGetCone(K, q, &cone));
8592dce792eSToby Isaac         for (c = 0; c < coneSize; c++)
8602dce792eSToby Isaac           if (cone[c] == p) break;
8612dce792eSToby Isaac         PetscCheck(c != coneSize, PetscObjectComm((PetscObject)K), PETSC_ERR_PLIB, "cone/support mismatch");
8622dce792eSToby Isaac         if (!qsp) {
8632dce792eSToby Isaac           sp->pointSpaces[p] = NULL;
8642dce792eSToby Isaac           continue;
8652dce792eSToby Isaac         }
8662dce792eSToby Isaac         PetscCall(PetscDualSpaceGetDM(qsp, &qdm));
8672dce792eSToby Isaac         PetscCall(DMPlexGetCone(qdm, 0, &refCone));
8682dce792eSToby Isaac         /* get the equivalent dual space from the support dual space */
8692dce792eSToby Isaac         PetscCall(PetscDualSpaceGetPointSubspace(qsp, refCone[c], &psp));
8702dce792eSToby Isaac         PetscCall(PetscObjectReference((PetscObject)psp));
8712dce792eSToby Isaac         sp->pointSpaces[p] = psp;
8722dce792eSToby Isaac       }
8732dce792eSToby Isaac     }
8742dce792eSToby Isaac     PetscCall(PetscFree2(pStratStart, pStratEnd));
8752dce792eSToby Isaac   }
8762dce792eSToby Isaac 
8772dce792eSToby Isaac   sum->uniform = uniform;
8782dce792eSToby Isaac   PetscCall(PetscCalloc1(Ns, &sum->all_rows));
8792dce792eSToby Isaac   PetscCall(PetscCalloc1(Ns, &sum->all_cols));
8802dce792eSToby Isaac   PetscCall(PetscCalloc1(Ns, &sum->int_rows));
8812dce792eSToby Isaac   PetscCall(PetscCalloc1(Ns, &sum->int_cols));
8822dce792eSToby Isaac   PetscCall(PetscMalloc4(Ns, &all_quads, Ns, &all_mats, Ns, &int_quads, Ns, &int_mats));
8832dce792eSToby Isaac   {
8842dce792eSToby Isaac     // test for uniform all points and uniform interior points
8852dce792eSToby Isaac     PetscBool       uniform_all         = PETSC_TRUE;
8862dce792eSToby Isaac     PetscBool       uniform_interior    = PETSC_TRUE;
8872dce792eSToby Isaac     PetscQuadrature quad_all_first      = NULL;
8882dce792eSToby Isaac     PetscQuadrature quad_interior_first = NULL;
8892dce792eSToby Isaac     PetscInt        pStart, pEnd;
8902dce792eSToby Isaac 
8912dce792eSToby Isaac     PetscCall(DMPlexGetChart(K, &pStart, &pEnd));
8922dce792eSToby Isaac     PetscCall(PetscDualSpaceSectionCreate_Internal(sp, &sp->pointSection));
8932dce792eSToby Isaac 
8942dce792eSToby Isaac     for (PetscInt p = pStart; p < pEnd; p++) {
8952dce792eSToby Isaac       PetscInt full_dof = 0;
8962dce792eSToby Isaac 
8972dce792eSToby Isaac       for (PetscInt s = 0; s < Ns; s++) {
8982dce792eSToby Isaac         PetscDualSpace subsp;
8992dce792eSToby Isaac         PetscSection   subsection;
9002dce792eSToby Isaac         PetscInt       dof;
9012dce792eSToby Isaac 
9022dce792eSToby Isaac         PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
9032dce792eSToby Isaac         PetscCall(PetscDualSpaceGetSection(subsp, &subsection));
9042dce792eSToby Isaac         PetscCall(PetscSectionGetDof(subsection, p, &dof));
9052dce792eSToby Isaac         full_dof += dof;
9062dce792eSToby Isaac       }
9072dce792eSToby Isaac       PetscCall(PetscSectionSetDof(sp->pointSection, p, full_dof));
9082dce792eSToby Isaac     }
9092dce792eSToby Isaac     PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, sp->pointSection));
9102dce792eSToby Isaac 
9112dce792eSToby Isaac     for (PetscInt s = 0; s < Ns; s++) {
9122dce792eSToby Isaac       PetscDualSpace  subsp;
9132dce792eSToby Isaac       PetscQuadrature subquad_all;
9142dce792eSToby Isaac       PetscQuadrature subquad_interior;
9152dce792eSToby Isaac       Mat             submat_all;
9162dce792eSToby Isaac       Mat             submat_interior;
9172dce792eSToby Isaac 
9182dce792eSToby Isaac       PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
9192dce792eSToby Isaac       PetscCall(PetscDualSpaceGetAllData(subsp, &subquad_all, &submat_all));
9202dce792eSToby Isaac       PetscCall(PetscDualSpaceGetInteriorData(subsp, &subquad_interior, &submat_interior));
9212dce792eSToby Isaac       if (!s) {
9222dce792eSToby Isaac         quad_all_first      = subquad_all;
9232dce792eSToby Isaac         quad_interior_first = subquad_interior;
9242dce792eSToby Isaac       } else {
9252dce792eSToby Isaac         if (subquad_all != quad_all_first) uniform_all = PETSC_FALSE;
9262dce792eSToby Isaac         if (subquad_interior != quad_interior_first) uniform_interior = PETSC_FALSE;
9272dce792eSToby Isaac       }
9282dce792eSToby Isaac       all_quads[s] = subquad_all;
9292dce792eSToby Isaac       int_quads[s] = subquad_interior;
9302dce792eSToby Isaac       all_mats[s]  = submat_all;
9312dce792eSToby Isaac       int_mats[s]  = submat_interior;
9322dce792eSToby Isaac     }
9332dce792eSToby Isaac     sum->uniform_all_points      = uniform_all;
9342dce792eSToby Isaac     sum->uniform_interior_points = uniform_interior;
9352dce792eSToby Isaac 
9362dce792eSToby Isaac     PetscCall(PetscDualSpaceSumCreateMappings(sp, PETSC_TRUE, uniform_interior, sum->int_rows, sum->int_cols));
9372dce792eSToby Isaac     PetscCall(PetscDualSpaceSumCreateQuadrature(Ns, cdim, uniform_interior, int_quads, &sp->intNodes));
9382dce792eSToby Isaac     if (sp->intNodes) PetscCall(PetscDualSpaceSumCreateMatrix(sp, int_mats, sum->int_rows, sum->int_cols, sp->intNodes, &sp->intMat));
9392dce792eSToby Isaac 
9402dce792eSToby Isaac     PetscCall(PetscDualSpaceSumCreateMappings(sp, PETSC_FALSE, uniform_all, sum->all_rows, sum->all_cols));
9412dce792eSToby Isaac     PetscCall(PetscDualSpaceSumCreateQuadrature(Ns, cdim, uniform_all, all_quads, &sp->allNodes));
9422dce792eSToby Isaac     if (sp->allNodes) PetscCall(PetscDualSpaceSumCreateMatrix(sp, all_mats, sum->all_rows, sum->all_cols, sp->allNodes, &sp->allMat));
9432dce792eSToby Isaac   }
9442dce792eSToby Isaac   PetscCall(PetscFree4(all_quads, all_mats, int_quads, int_mats));
9452dce792eSToby Isaac   PetscCall(PetscDualSpaceComputeFunctionalsFromAllData(sp));
9462dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
9472dce792eSToby Isaac }
9482dce792eSToby Isaac 
9492dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumView_Ascii(PetscDualSpace sp, PetscViewer v)
9502dce792eSToby Isaac {
9512dce792eSToby Isaac   PetscDualSpace_Sum *sum         = (PetscDualSpace_Sum *)sp->data;
9522dce792eSToby Isaac   PetscBool           concatenate = sum->concatenate;
9532dce792eSToby Isaac   PetscInt            i, Ns = sum->numSumSpaces;
9542dce792eSToby Isaac 
9552dce792eSToby Isaac   PetscFunctionBegin;
9562dce792eSToby Isaac   if (concatenate) PetscCall(PetscViewerASCIIPrintf(v, "Sum dual space of %" PetscInt_FMT " concatenated subspaces%s\n", Ns, sum->uniform ? " (all identical)" : ""));
9572dce792eSToby Isaac   else PetscCall(PetscViewerASCIIPrintf(v, "Sum dual space of %" PetscInt_FMT " subspaces%s\n", Ns, sum->uniform ? " (all identical)" : ""));
9582dce792eSToby Isaac   for (i = 0; i < (sum->uniform ? (Ns > 0 ? 1 : 0) : Ns); ++i) {
9592dce792eSToby Isaac     PetscCall(PetscViewerASCIIPushTab(v));
9602dce792eSToby Isaac     PetscCall(PetscDualSpaceView(sum->sumspaces[i], v));
9612dce792eSToby Isaac     PetscCall(PetscViewerASCIIPopTab(v));
9622dce792eSToby Isaac   }
9632dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
9642dce792eSToby Isaac }
9652dce792eSToby Isaac 
9662dce792eSToby Isaac static PetscErrorCode PetscDualSpaceView_Sum(PetscDualSpace sp, PetscViewer viewer)
9672dce792eSToby Isaac {
9682dce792eSToby Isaac   PetscBool iascii;
9692dce792eSToby Isaac 
9702dce792eSToby Isaac   PetscFunctionBegin;
9712dce792eSToby Isaac   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
9722dce792eSToby Isaac   if (iascii) PetscCall(PetscDualSpaceSumView_Ascii(sp, viewer));
9732dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
9742dce792eSToby Isaac }
9752dce792eSToby Isaac 
9762dce792eSToby Isaac static PetscErrorCode PetscDualSpaceDestroy_Sum(PetscDualSpace sp)
9772dce792eSToby Isaac {
9782dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data;
9792dce792eSToby Isaac   PetscInt            i, Ns = sum->numSumSpaces;
9802dce792eSToby Isaac 
9812dce792eSToby Isaac   PetscFunctionBegin;
9822dce792eSToby Isaac   if (sum->symperms) {
9832dce792eSToby Isaac     PetscInt **selfSyms = sum->symperms[0];
9842dce792eSToby Isaac 
9852dce792eSToby Isaac     if (selfSyms) {
9862dce792eSToby Isaac       PetscInt i, **allocated = &selfSyms[-sum->selfSymOff];
9872dce792eSToby Isaac 
9882dce792eSToby Isaac       for (i = 0; i < sum->numSelfSym; i++) PetscCall(PetscFree(allocated[i]));
9892dce792eSToby Isaac       PetscCall(PetscFree(allocated));
9902dce792eSToby Isaac     }
9912dce792eSToby Isaac     PetscCall(PetscFree(sum->symperms));
9922dce792eSToby Isaac   }
9932dce792eSToby Isaac   if (sum->symflips) {
9942dce792eSToby Isaac     PetscScalar **selfSyms = sum->symflips[0];
9952dce792eSToby Isaac 
9962dce792eSToby Isaac     if (selfSyms) {
9972dce792eSToby Isaac       PetscInt      i;
9982dce792eSToby Isaac       PetscScalar **allocated = &selfSyms[-sum->selfSymOff];
9992dce792eSToby Isaac 
10002dce792eSToby Isaac       for (i = 0; i < sum->numSelfSym; i++) PetscCall(PetscFree(allocated[i]));
10012dce792eSToby Isaac       PetscCall(PetscFree(allocated));
10022dce792eSToby Isaac     }
10032dce792eSToby Isaac     PetscCall(PetscFree(sum->symflips));
10042dce792eSToby Isaac   }
10052dce792eSToby Isaac   for (i = 0; i < Ns; ++i) {
10062dce792eSToby Isaac     PetscCall(PetscDualSpaceDestroy(&sum->sumspaces[i]));
10072dce792eSToby Isaac     if (sum->all_rows) PetscCall(ISLocalToGlobalMappingDestroy(&sum->all_rows[i]));
10082dce792eSToby Isaac     if (sum->all_cols) PetscCall(ISLocalToGlobalMappingDestroy(&sum->all_cols[i]));
10092dce792eSToby Isaac     if (sum->int_rows) PetscCall(ISLocalToGlobalMappingDestroy(&sum->int_rows[i]));
10102dce792eSToby Isaac     if (sum->int_cols) PetscCall(ISLocalToGlobalMappingDestroy(&sum->int_cols[i]));
10112dce792eSToby Isaac   }
10122dce792eSToby Isaac   PetscCall(PetscFree(sum->sumspaces));
10132dce792eSToby Isaac   PetscCall(PetscFree(sum->all_rows));
10142dce792eSToby Isaac   PetscCall(PetscFree(sum->all_cols));
10152dce792eSToby Isaac   PetscCall(PetscFree(sum->int_rows));
10162dce792eSToby Isaac   PetscCall(PetscFree(sum->int_cols));
10172dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetSubspace_C", NULL));
10182dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetSubspace_C", NULL));
10192dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetNumSubspaces_C", NULL));
10202dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetNumSubspaces_C", NULL));
10212dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetConcatenate_C", NULL));
10222dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetConcatenate_C", NULL));
10232dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetInterleave_C", NULL));
10242dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetInterleave_C", NULL));
10252dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetContinuity_C", NULL));
10262dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetContinuity_C", NULL));
10272dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetMomentOrder_C", NULL));
10282dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetMomentOrder_C", NULL));
10292dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetNodeType_C", NULL));
10302dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetNodeType_C", NULL));
10312dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTensor_C", NULL));
10322dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTensor_C", NULL));
10332dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTrimmed_C", NULL));
10342dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTrimmed_C", NULL));
10352dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetUseMoments_C", NULL));
10362dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetUseMoments_C", NULL));
10372dce792eSToby Isaac   PetscCall(PetscFree(sum));
10382dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
10392dce792eSToby Isaac }
10402dce792eSToby Isaac 
10412dce792eSToby Isaac /*@
10422dce792eSToby Isaac   PetscDualSpaceSumSetInterleave - Set whether the basis functions and components of a uniform sum are interleaved
10432dce792eSToby Isaac 
10442dce792eSToby Isaac   Logically collective
10452dce792eSToby Isaac 
10462dce792eSToby Isaac   Input Parameters:
10472dce792eSToby Isaac + sp                    - a `PetscDualSpace` of type `PETSCDUALSPACESUM`
10482dce792eSToby Isaac . interleave_basis      - if `PETSC_TRUE`, the basis vectors of the subspaces are interleaved
10492dce792eSToby Isaac - interleave_components - if `PETSC_TRUE` and the space concatenates components (`PetscDualSpaceSumGetConcatenate()`),
10502dce792eSToby Isaac                           interleave the concatenated components
10512dce792eSToby Isaac 
10522dce792eSToby Isaac   Level: developer
10532dce792eSToby Isaac 
10542dce792eSToby Isaac .seealso: `PetscDualSpace`, `PETSCDUALSPACESUM`, `PETSCFEVECTOR`, `PetscDualSpaceSumGetInterleave()`
10552dce792eSToby Isaac @*/
10562dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumSetInterleave(PetscDualSpace sp, PetscBool interleave_basis, PetscBool interleave_components)
10572dce792eSToby Isaac {
10582dce792eSToby Isaac   PetscFunctionBegin;
10592dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
10602dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumSetInterleave_C", (PetscDualSpace, PetscBool, PetscBool), (sp, interleave_basis, interleave_components));
10612dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
10622dce792eSToby Isaac }
10632dce792eSToby Isaac 
10642dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumSetInterleave_Sum(PetscDualSpace sp, PetscBool interleave_basis, PetscBool interleave_components)
10652dce792eSToby Isaac {
10662dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data;
10674d86920dSPierre Jolivet 
10682dce792eSToby Isaac   PetscFunctionBegin;
10692dce792eSToby Isaac   sum->interleave_basis      = interleave_basis;
10702dce792eSToby Isaac   sum->interleave_components = interleave_components;
10712dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
10722dce792eSToby Isaac }
10732dce792eSToby Isaac 
10742dce792eSToby Isaac /*@
10752dce792eSToby Isaac   PetscDualSpaceSumGetInterleave - Get whether the basis functions and components of a uniform sum are interleaved
10762dce792eSToby Isaac 
10772dce792eSToby Isaac   Logically collective
10782dce792eSToby Isaac 
10792dce792eSToby Isaac   Input Parameter:
10802dce792eSToby Isaac . sp - a `PetscDualSpace` of type `PETSCDUALSPACESUM`
10812dce792eSToby Isaac 
10822dce792eSToby Isaac   Output Parameters:
10832dce792eSToby Isaac + interleave_basis      - if `PETSC_TRUE`, the basis vectors of the subspaces are interleaved
10842dce792eSToby Isaac - interleave_components - if `PETSC_TRUE` and the space concatenates components (`PetscDualSpaceSumGetConcatenate()`),
10852dce792eSToby Isaac                           interleave the concatenated components
10862dce792eSToby Isaac 
10872dce792eSToby Isaac   Level: developer
10882dce792eSToby Isaac 
10892dce792eSToby Isaac .seealso: `PetscDualSpace`, `PETSCDUALSPACESUM`, `PETSCFEVECTOR`, `PetscDualSpaceSumSetInterleave()`
10902dce792eSToby Isaac @*/
10912dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumGetInterleave(PetscDualSpace sp, PetscBool *interleave_basis, PetscBool *interleave_components)
10922dce792eSToby Isaac {
10932dce792eSToby Isaac   PetscFunctionBegin;
10942dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
10952dce792eSToby Isaac   if (interleave_basis) PetscAssertPointer(interleave_basis, 2);
10962dce792eSToby Isaac   if (interleave_components) PetscAssertPointer(interleave_components, 3);
10972dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumGetInterleave_C", (PetscDualSpace, PetscBool *, PetscBool *), (sp, interleave_basis, interleave_components));
10982dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
10992dce792eSToby Isaac }
11002dce792eSToby Isaac 
11012dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumGetInterleave_Sum(PetscDualSpace sp, PetscBool *interleave_basis, PetscBool *interleave_components)
11022dce792eSToby Isaac {
11032dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data;
11044d86920dSPierre Jolivet 
11052dce792eSToby Isaac   PetscFunctionBegin;
11062dce792eSToby Isaac   if (interleave_basis) *interleave_basis = sum->interleave_basis;
11072dce792eSToby Isaac   if (interleave_components) *interleave_components = sum->interleave_components;
11082dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11092dce792eSToby Isaac }
11102dce792eSToby Isaac 
11112dce792eSToby Isaac #define PetscDualSpaceSumPassthrough(sp, func, ...) \
11122dce792eSToby Isaac   do { \
11132dce792eSToby Isaac     PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data; \
11142dce792eSToby Isaac     PetscBool           is_uniform; \
11152dce792eSToby Isaac     PetscCall(PetscDualSpaceSumIsUniform(sp, &is_uniform)); \
11162dce792eSToby Isaac     if (is_uniform && sum->numSumSpaces > 0) { \
11172dce792eSToby Isaac       PetscDualSpace subsp; \
11182dce792eSToby Isaac       PetscCall(PetscDualSpaceSumGetSubspace(sp, 0, &subsp)); \
11192dce792eSToby Isaac       PetscCall(func(subsp, __VA_ARGS__)); \
11202dce792eSToby Isaac     } \
11212dce792eSToby Isaac   } while (0)
11222dce792eSToby Isaac 
11232dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetContinuity_Sum(PetscDualSpace sp, PetscBool *value)
11242dce792eSToby Isaac {
11252dce792eSToby Isaac   PetscFunctionBegin;
11262dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetContinuity, value);
11272dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11282dce792eSToby Isaac }
11292dce792eSToby Isaac 
11302dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetContinuity_Sum(PetscDualSpace sp, PetscBool value)
11312dce792eSToby Isaac {
11322dce792eSToby Isaac   PetscFunctionBegin;
11332dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetContinuity, value);
11342dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11352dce792eSToby Isaac }
11362dce792eSToby Isaac 
11372dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetTensor_Sum(PetscDualSpace sp, PetscBool *value)
11382dce792eSToby Isaac {
11392dce792eSToby Isaac   PetscFunctionBegin;
11402dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetTensor, value);
11412dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11422dce792eSToby Isaac }
11432dce792eSToby Isaac 
11442dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetTensor_Sum(PetscDualSpace sp, PetscBool value)
11452dce792eSToby Isaac {
11462dce792eSToby Isaac   PetscFunctionBegin;
11472dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetTensor, value);
11482dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11492dce792eSToby Isaac }
11502dce792eSToby Isaac 
11512dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetTrimmed_Sum(PetscDualSpace sp, PetscBool *value)
11522dce792eSToby Isaac {
11532dce792eSToby Isaac   PetscFunctionBegin;
11542dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetTrimmed, value);
11552dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11562dce792eSToby Isaac }
11572dce792eSToby Isaac 
11582dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetTrimmed_Sum(PetscDualSpace sp, PetscBool value)
11592dce792eSToby Isaac {
11602dce792eSToby Isaac   PetscFunctionBegin;
11612dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetTrimmed, value);
11622dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11632dce792eSToby Isaac }
11642dce792eSToby Isaac 
11652dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetUseMoments_Sum(PetscDualSpace sp, PetscBool *value)
11662dce792eSToby Isaac {
11672dce792eSToby Isaac   PetscFunctionBegin;
11682dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetUseMoments, value);
11692dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11702dce792eSToby Isaac }
11712dce792eSToby Isaac 
11722dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetUseMoments_Sum(PetscDualSpace sp, PetscBool value)
11732dce792eSToby Isaac {
11742dce792eSToby Isaac   PetscFunctionBegin;
11752dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetUseMoments, value);
11762dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11772dce792eSToby Isaac }
11782dce792eSToby Isaac 
11792dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetMomentOrder_Sum(PetscDualSpace sp, PetscInt *value)
11802dce792eSToby Isaac {
11812dce792eSToby Isaac   PetscFunctionBegin;
11822dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetMomentOrder, value);
11832dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11842dce792eSToby Isaac }
11852dce792eSToby Isaac 
11862dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetMomentOrder_Sum(PetscDualSpace sp, PetscInt value)
11872dce792eSToby Isaac {
11882dce792eSToby Isaac   PetscFunctionBegin;
11892dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetMomentOrder, value);
11902dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11912dce792eSToby Isaac }
11922dce792eSToby Isaac 
11932dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetNodeType_Sum(PetscDualSpace sp, PetscDTNodeType *node_type, PetscBool *include_endpoints, PetscReal *exponent)
11942dce792eSToby Isaac {
11952dce792eSToby Isaac   PetscFunctionBegin;
11962dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetNodeType, node_type, include_endpoints, exponent);
11972dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11982dce792eSToby Isaac }
11992dce792eSToby Isaac 
12002dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetNodeType_Sum(PetscDualSpace sp, PetscDTNodeType node_type, PetscBool include_endpoints, PetscReal exponent)
12012dce792eSToby Isaac {
12022dce792eSToby Isaac   PetscFunctionBegin;
12032dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetNodeType, node_type, include_endpoints, exponent);
12042dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
12052dce792eSToby Isaac }
12062dce792eSToby Isaac 
12072dce792eSToby Isaac static PetscErrorCode PetscDualSpaceInitialize_Sum(PetscDualSpace sp)
12082dce792eSToby Isaac {
12092dce792eSToby Isaac   PetscFunctionBegin;
12102dce792eSToby Isaac   sp->ops->destroy              = PetscDualSpaceDestroy_Sum;
12112dce792eSToby Isaac   sp->ops->view                 = PetscDualSpaceView_Sum;
12122dce792eSToby Isaac   sp->ops->setfromoptions       = NULL;
12132dce792eSToby Isaac   sp->ops->duplicate            = PetscDualSpaceDuplicate_Sum;
12142dce792eSToby Isaac   sp->ops->setup                = PetscDualSpaceSetUp_Sum;
12152dce792eSToby Isaac   sp->ops->createheightsubspace = NULL;
12162dce792eSToby Isaac   sp->ops->createpointsubspace  = NULL;
12172dce792eSToby Isaac   sp->ops->getsymmetries        = PetscDualSpaceGetSymmetries_Sum;
12182dce792eSToby Isaac   sp->ops->apply                = PetscDualSpaceApplyDefault;
12192dce792eSToby Isaac   sp->ops->applyall             = PetscDualSpaceApplyAllDefault;
12202dce792eSToby Isaac   sp->ops->applyint             = PetscDualSpaceApplyInteriorDefault;
12212dce792eSToby Isaac   sp->ops->createalldata        = PetscDualSpaceCreateAllDataDefault;
12222dce792eSToby Isaac   sp->ops->createintdata        = PetscDualSpaceCreateInteriorDataDefault;
12232dce792eSToby Isaac 
12242dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetNumSubspaces_C", PetscDualSpaceSumGetNumSubspaces_Sum));
12252dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetNumSubspaces_C", PetscDualSpaceSumSetNumSubspaces_Sum));
12262dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetSubspace_C", PetscDualSpaceSumGetSubspace_Sum));
12272dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetSubspace_C", PetscDualSpaceSumSetSubspace_Sum));
12282dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetConcatenate_C", PetscDualSpaceSumGetConcatenate_Sum));
12292dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetConcatenate_C", PetscDualSpaceSumSetConcatenate_Sum));
12302dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetInterleave_C", PetscDualSpaceSumGetInterleave_Sum));
12312dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetInterleave_C", PetscDualSpaceSumSetInterleave_Sum));
12322dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetContinuity_C", PetscDualSpaceLagrangeGetContinuity_Sum));
12332dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetContinuity_C", PetscDualSpaceLagrangeSetContinuity_Sum));
12342dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetMomentOrder_C", PetscDualSpaceLagrangeGetMomentOrder_Sum));
12352dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetMomentOrder_C", PetscDualSpaceLagrangeSetMomentOrder_Sum));
12362dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetNodeType_C", PetscDualSpaceLagrangeGetNodeType_Sum));
12372dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetNodeType_C", PetscDualSpaceLagrangeSetNodeType_Sum));
12382dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTensor_C", PetscDualSpaceLagrangeGetTensor_Sum));
12392dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTensor_C", PetscDualSpaceLagrangeSetTensor_Sum));
12402dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTrimmed_C", PetscDualSpaceLagrangeGetTrimmed_Sum));
12412dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTrimmed_C", PetscDualSpaceLagrangeSetTrimmed_Sum));
12422dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetUseMoments_C", PetscDualSpaceLagrangeGetUseMoments_Sum));
12432dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetUseMoments_C", PetscDualSpaceLagrangeSetUseMoments_Sum));
12442dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
12452dce792eSToby Isaac }
12462dce792eSToby Isaac 
12472dce792eSToby Isaac /*MC
12482dce792eSToby Isaac   PETSCDUALSPACESUM = "sum" - A `PetscDualSpace` object that encapsulates a sum of subspaces.
12492dce792eSToby Isaac 
12502dce792eSToby Isaac   Level: intermediate
12512dce792eSToby Isaac 
12522dce792eSToby Isaac   Note:
12532dce792eSToby Isaac   That sum can either be direct or a concatenation. For example if A and B are spaces each with 2 components,
12542dce792eSToby Isaac   the direct sum of A and B will also have 2 components while the concatenated sum will have 4 components. In both cases A and B must be defined over the
12552dce792eSToby Isaac   same reference element.
12562dce792eSToby Isaac 
12572dce792eSToby Isaac .seealso: `PetscDualSpace`, `PetscDualSpaceType`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`, `PetscDualSpaceSumGetNumSubspaces()`, `PetscDualSpaceSumSetNumSubspaces()`,
12582dce792eSToby Isaac           `PetscDualSpaceSumGetConcatenate()`, `PetscDualSpaceSumSetConcatenate()`, `PetscDualSpaceSumSetInterleave()`, `PetscDualSpaceSumGetInterleave()`
12592dce792eSToby Isaac M*/
12602dce792eSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Sum(PetscDualSpace sp)
12612dce792eSToby Isaac {
12622dce792eSToby Isaac   PetscDualSpace_Sum *sum;
12632dce792eSToby Isaac 
12642dce792eSToby Isaac   PetscFunctionBegin;
12652dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
12662dce792eSToby Isaac   PetscCall(PetscNew(&sum));
12672dce792eSToby Isaac   sum->numSumSpaces = PETSC_DEFAULT;
12682dce792eSToby Isaac   sp->data          = sum;
12692dce792eSToby Isaac   sp->k             = PETSC_FORM_DEGREE_UNDEFINED;
12702dce792eSToby Isaac   PetscCall(PetscDualSpaceInitialize_Sum(sp));
12712dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
12722dce792eSToby Isaac }
12732dce792eSToby Isaac 
12742dce792eSToby Isaac /*@
12752dce792eSToby Isaac   PetscDualSpaceCreateSum - Create a finite element dual basis that is the sum of other dual bases
12762dce792eSToby Isaac 
12772dce792eSToby Isaac   Collective
12782dce792eSToby Isaac 
12792dce792eSToby Isaac   Input Parameters:
12802dce792eSToby Isaac + numSubspaces - the number of spaces that will be added together
12812dce792eSToby Isaac . subspaces    - an array of length `numSubspaces` of spaces
12822dce792eSToby Isaac - concatenate  - if `PETSC_FALSE`, the sum-space has the same components as the individual dual spaces (`PetscDualSpaceGetNumComponents()`); if `PETSC_TRUE`, the individual components are concatenated to create a dual space with more components
12832dce792eSToby Isaac 
12842dce792eSToby Isaac   Output Parameter:
12852dce792eSToby Isaac . sumSpace - a `PetscDualSpace` of type `PETSCDUALSPACESUM`
12862dce792eSToby Isaac 
12872dce792eSToby Isaac   Level: advanced
12882dce792eSToby Isaac 
12892dce792eSToby Isaac .seealso: `PetscDualSpace`, `PETSCDUALSPACESUM`, `PETSCSPACESUM`
12902dce792eSToby Isaac @*/
12912dce792eSToby Isaac PetscErrorCode PetscDualSpaceCreateSum(PetscInt numSubspaces, const PetscDualSpace subspaces[], PetscBool concatenate, PetscDualSpace *sumSpace)
12922dce792eSToby Isaac {
12932dce792eSToby Isaac   PetscInt i, Nc = 0;
12942dce792eSToby Isaac 
12952dce792eSToby Isaac   PetscFunctionBegin;
12962dce792eSToby Isaac   PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)subspaces[0]), sumSpace));
12972dce792eSToby Isaac   PetscCall(PetscDualSpaceSetType(*sumSpace, PETSCDUALSPACESUM));
12982dce792eSToby Isaac   PetscCall(PetscDualSpaceSumSetNumSubspaces(*sumSpace, numSubspaces));
12992dce792eSToby Isaac   PetscCall(PetscDualSpaceSumSetConcatenate(*sumSpace, concatenate));
13002dce792eSToby Isaac   for (i = 0; i < numSubspaces; ++i) {
13012dce792eSToby Isaac     PetscInt sNc;
13022dce792eSToby Isaac 
13032dce792eSToby Isaac     PetscCall(PetscDualSpaceSumSetSubspace(*sumSpace, i, subspaces[i]));
13042dce792eSToby Isaac     PetscCall(PetscDualSpaceGetNumComponents(subspaces[i], &sNc));
13052dce792eSToby Isaac     if (concatenate) Nc += sNc;
13062dce792eSToby Isaac     else Nc = sNc;
13072dce792eSToby Isaac 
13082dce792eSToby Isaac     if (i == 0) {
13092dce792eSToby Isaac       DM dm;
13102dce792eSToby Isaac 
13112dce792eSToby Isaac       PetscCall(PetscDualSpaceGetDM(subspaces[i], &dm));
13122dce792eSToby Isaac       PetscCall(PetscDualSpaceSetDM(*sumSpace, dm));
13132dce792eSToby Isaac     }
13142dce792eSToby Isaac   }
13152dce792eSToby Isaac   PetscCall(PetscDualSpaceSetNumComponents(*sumSpace, Nc));
13162dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
13172dce792eSToby Isaac }
1318