xref: /petsc/src/dm/dt/dualspace/impls/sum/dualspacesum.c (revision 834855d6effb0d027771461c8e947ee1ce5a1e17)
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 @*/
PetscDualSpaceSumGetNumSubspaces(PetscDualSpace sp,PetscInt * numSumSpaces)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 @*/
PetscDualSpaceSumSetNumSubspaces(PetscDualSpace sp,PetscInt numSumSpaces)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 @*/
PetscDualSpaceSumGetConcatenate(PetscDualSpace sp,PetscBool * concatenate)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 @*/
PetscDualSpaceSumSetConcatenate(PetscDualSpace sp,PetscBool concatenate)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 @*/
PetscDualSpaceSumGetSubspace(PetscDualSpace sp,PetscInt s,PetscDualSpace * subsp)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 @*/
PetscDualSpaceSumSetSubspace(PetscDualSpace sp,PetscInt s,PetscDualSpace subsp)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 
PetscDualSpaceSumGetNumSubspaces_Sum(PetscDualSpace space,PetscInt * numSumSpaces)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 
PetscDualSpaceSumSetNumSubspaces_Sum(PetscDualSpace space,PetscInt numSumSpaces)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;
165371d2eb7SMartin Diehl   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 
PetscDualSpaceSumGetConcatenate_Sum(PetscDualSpace sp,PetscBool * concatenate)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 
PetscDualSpaceSumSetConcatenate_Sum(PetscDualSpace sp,PetscBool concatenate)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;
189371d2eb7SMartin Diehl   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 
PetscDualSpaceSumGetSubspace_Sum(PetscDualSpace space,PetscInt s,PetscDualSpace * subspace)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 
PetscDualSpaceSumSetSubspace_Sum(PetscDualSpace space,PetscInt s,PetscDualSpace subspace)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;
214371d2eb7SMartin Diehl   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 
PetscDualSpaceDuplicate_Sum(PetscDualSpace sp,PetscDualSpace spNew)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 
PetscDualSpaceSumCreateQuadrature(PetscInt Ns,PetscInt cdim,PetscBool uniform_points,PetscQuadrature subquads[],PetscQuadrature * fullquad)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 
PetscDualSpaceSumCreateMatrix(PetscDualSpace sp,Mat submats[],ISLocalToGlobalMapping map_rows[],ISLocalToGlobalMapping map_cols[],PetscQuadrature fullquad,Mat * fullmat)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;
304dd460d27SBarry Smith   PetscInt     nrows = 0, ncols, nnz = 0;
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   PetscCall(PetscQuadratureGetData(fullquad, NULL, NULL, &npoints, NULL, NULL));
3112dce792eSToby Isaac   ncols = Nc * npoints;
3122dce792eSToby Isaac   for (PetscInt s = 0; s < num_subspaces; s++) {
3132dce792eSToby Isaac     // Get the COO data for each matrix, map the is and js, and append to growing COO data
3142dce792eSToby Isaac     PetscInt        sNrows, sNcols;
3152dce792eSToby Isaac     Mat             smat;
3162dce792eSToby Isaac     const PetscInt *si;
3172dce792eSToby Isaac     const PetscInt *sj;
3182dce792eSToby Isaac     PetscScalar    *sv;
3192dce792eSToby Isaac     PetscMemType    memtype;
3202dce792eSToby Isaac     PetscInt        snz;
3212dce792eSToby Isaac     PetscInt        snz_actual;
3222dce792eSToby Isaac     PetscInt       *cooi;
3232dce792eSToby Isaac     PetscInt       *cooj;
3242dce792eSToby Isaac     PetscScalar    *coov;
3252dce792eSToby Isaac     PetscScalar    *v_new;
3262dce792eSToby Isaac     PetscInt       *i_new;
3272dce792eSToby Isaac     PetscInt       *j_new;
3282dce792eSToby Isaac 
3292dce792eSToby Isaac     if (!submats[s]) continue;
3302dce792eSToby Isaac     PetscCall(MatGetSize(submats[s], &sNrows, &sNcols));
3312dce792eSToby Isaac     nrows += sNrows;
3322dce792eSToby Isaac     PetscCall(MatConvert(submats[s], MATSEQAIJ, MAT_INITIAL_MATRIX, &smat));
3332dce792eSToby Isaac     PetscCall(MatBindToCPU(smat, PETSC_TRUE));
3342dce792eSToby Isaac     PetscCall(MatSeqAIJGetCSRAndMemType(smat, &si, &sj, &sv, &memtype));
3352dce792eSToby Isaac     PetscCheck(memtype == PETSC_MEMTYPE_HOST, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Could not convert matrix to host memory");
3362dce792eSToby Isaac     snz = si[sNrows];
3372dce792eSToby Isaac 
3382dce792eSToby Isaac     PetscCall(PetscMalloc1(nnz + snz, &v_new));
3392dce792eSToby Isaac     PetscCall(PetscArraycpy(v_new, v, nnz));
3402dce792eSToby Isaac     PetscCall(PetscFree(v));
3412dce792eSToby Isaac     v = v_new;
3422dce792eSToby Isaac 
3432dce792eSToby Isaac     PetscCall(PetscMalloc1(nnz + snz, &i_new));
3442dce792eSToby Isaac     PetscCall(PetscArraycpy(i_new, i, nnz));
3452dce792eSToby Isaac     PetscCall(PetscFree(i));
3462dce792eSToby Isaac     i = i_new;
3472dce792eSToby Isaac 
3482dce792eSToby Isaac     PetscCall(PetscMalloc1(nnz + snz, &j_new));
3492dce792eSToby Isaac     PetscCall(PetscArraycpy(j_new, j, nnz));
3502dce792eSToby Isaac     PetscCall(PetscFree(j));
3512dce792eSToby Isaac     j = j_new;
3522dce792eSToby Isaac 
3532dce792eSToby Isaac     PetscCall(PetscMalloc2(snz, &cooi, snz, &cooj));
3542dce792eSToby Isaac 
3552dce792eSToby Isaac     coov = &v[nnz];
3562dce792eSToby Isaac 
3572dce792eSToby Isaac     snz_actual = 0;
3582dce792eSToby Isaac     for (PetscInt row = 0; row < sNrows; row++) {
3592dce792eSToby Isaac       for (PetscInt k = si[row]; k < si[row + 1]; k++, snz_actual++) {
3602dce792eSToby Isaac         cooi[snz_actual] = row;
3612dce792eSToby Isaac         cooj[snz_actual] = sj[k];
3622dce792eSToby Isaac         coov[snz_actual] = sv[k];
3632dce792eSToby Isaac       }
3642dce792eSToby Isaac     }
3652dce792eSToby Isaac     PetscCall(MatDestroy(&smat));
3662dce792eSToby Isaac 
3672dce792eSToby Isaac     if (snz_actual > 0) {
3682dce792eSToby Isaac       PetscCall(ISLocalToGlobalMappingApply(map_rows[s], snz_actual, cooi, &i[nnz]));
3692dce792eSToby Isaac       PetscCall(ISLocalToGlobalMappingApply(map_cols[s], snz_actual, cooj, &j[nnz]));
3702dce792eSToby Isaac       nnz += snz_actual;
3712dce792eSToby Isaac     }
3722dce792eSToby Isaac     PetscCall(PetscFree2(cooi, cooj));
3732dce792eSToby Isaac   }
3742dce792eSToby Isaac   PetscCall(MatCreate(PETSC_COMM_SELF, fullmat));
3752dce792eSToby Isaac   mat = *fullmat;
3762dce792eSToby Isaac   PetscCall(MatSetSizes(mat, nrows, ncols, nrows, ncols));
3772dce792eSToby Isaac   PetscCall(MatSetType(mat, MATSEQAIJ));
3782dce792eSToby Isaac   PetscCall(MatSetPreallocationCOO(mat, nnz, i, j));
3792dce792eSToby Isaac   PetscCall(MatSetValuesCOO(mat, v, INSERT_VALUES));
3802dce792eSToby Isaac   PetscCall(PetscFree(i));
3812dce792eSToby Isaac   PetscCall(PetscFree(j));
3822dce792eSToby Isaac   PetscCall(PetscFree(v));
3832dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
3842dce792eSToby Isaac }
3852dce792eSToby Isaac 
PetscDualSpaceSumCreateMappings(PetscDualSpace sp,PetscBool interior,PetscBool uniform_points,ISLocalToGlobalMapping map_row[],ISLocalToGlobalMapping map_col[])3862dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumCreateMappings(PetscDualSpace sp, PetscBool interior, PetscBool uniform_points, ISLocalToGlobalMapping map_row[], ISLocalToGlobalMapping map_col[])
3872dce792eSToby Isaac {
3882dce792eSToby Isaac   PetscSection section;
3892dce792eSToby Isaac   PetscInt     pStart, pEnd, Ns, Nc;
3902dce792eSToby Isaac   PetscInt     num_points = 0;
3912dce792eSToby Isaac   PetscBool    interleave_basis, interleave_components, concatenate;
3922dce792eSToby Isaac   PetscInt    *roffset;
3932dce792eSToby Isaac 
3942dce792eSToby Isaac   PetscFunctionBegin;
3952dce792eSToby Isaac   if (interior) {
3962dce792eSToby Isaac     PetscCall(PetscDualSpaceGetInteriorSection(sp, &section));
3972dce792eSToby Isaac   } else {
3982dce792eSToby Isaac     PetscCall(PetscDualSpaceGetSection(sp, &section));
3992dce792eSToby Isaac   }
4002dce792eSToby Isaac   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
4012dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetNumSubspaces(sp, &Ns));
4022dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetInterleave(sp, &interleave_basis, &interleave_components));
4032dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetConcatenate(sp, &concatenate));
4042dce792eSToby Isaac   PetscCall(PetscSectionGetChart(section, &pStart, &pEnd));
4052dce792eSToby Isaac   for (PetscInt p = pStart; p < pEnd; p++) {
4062dce792eSToby Isaac     PetscInt dof;
4072dce792eSToby Isaac 
4082dce792eSToby Isaac     PetscCall(PetscSectionGetDof(section, p, &dof));
4092dce792eSToby Isaac     num_points += (dof > 0);
4102dce792eSToby Isaac   }
4112dce792eSToby Isaac   PetscCall(PetscCalloc1(pEnd - pStart, &roffset));
4122dce792eSToby Isaac   for (PetscInt s = 0, coffset = 0; s < Ns; s++) {
4132dce792eSToby Isaac     PetscDualSpace  subsp;
4142dce792eSToby Isaac     PetscSection    subsection;
4152dce792eSToby Isaac     IS              is_row, is_col;
4162dce792eSToby Isaac     PetscInt        subNb;
4172dce792eSToby Isaac     PetscInt        sNc, sNpoints, sNcols;
4182dce792eSToby Isaac     PetscQuadrature quad;
4192dce792eSToby Isaac 
4202dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
4212dce792eSToby Isaac     PetscCall(PetscDualSpaceGetNumComponents(subsp, &sNc));
4222dce792eSToby Isaac     if (interior) {
4232dce792eSToby Isaac       PetscCall(PetscDualSpaceGetInteriorSection(subsp, &subsection));
4242dce792eSToby Isaac       PetscCall(PetscDualSpaceGetInteriorData(subsp, &quad, NULL));
4252dce792eSToby Isaac     } else {
4262dce792eSToby Isaac       PetscCall(PetscDualSpaceGetSection(subsp, &subsection));
4272dce792eSToby Isaac       PetscCall(PetscDualSpaceGetAllData(subsp, &quad, NULL));
4282dce792eSToby Isaac     }
4292dce792eSToby Isaac     PetscCall(PetscSectionGetStorageSize(subsection, &subNb));
4302dce792eSToby Isaac     if (num_points == 1) {
4312dce792eSToby Isaac       PetscInt rstride;
4322dce792eSToby Isaac 
4332dce792eSToby Isaac       rstride = interleave_basis ? Ns : 1;
4342dce792eSToby Isaac       PetscCall(ISCreateStride(PETSC_COMM_SELF, subNb, roffset[0], rstride, &is_row));
4352dce792eSToby Isaac       roffset[0] += interleave_basis ? 1 : subNb;
4362dce792eSToby Isaac     } else {
4372dce792eSToby Isaac       PetscInt *rows;
4382dce792eSToby Isaac 
4392dce792eSToby Isaac       PetscCall(PetscMalloc1(subNb, &rows));
4402dce792eSToby Isaac       for (PetscInt p = pStart; p < pEnd; p++) {
4412dce792eSToby Isaac         PetscInt subdof, dof, off, suboff, stride;
4422dce792eSToby Isaac 
4432dce792eSToby Isaac         PetscCall(PetscSectionGetOffset(subsection, p, &suboff));
4442dce792eSToby Isaac         PetscCall(PetscSectionGetDof(subsection, p, &subdof));
4452dce792eSToby Isaac         PetscCall(PetscSectionGetOffset(section, p, &off));
4462dce792eSToby Isaac         PetscCall(PetscSectionGetDof(section, p, &dof));
44700045ab3SPierre Jolivet         PetscCheck(subdof * Ns == dof || !interleave_basis, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Basis cannot be interleaved");
4482dce792eSToby Isaac         stride = interleave_basis ? Ns : 1;
449*ac530a7eSPierre Jolivet         for (PetscInt k = 0; k < subdof; k++) rows[suboff + k] = off + roffset[p - pStart] + k * stride;
4502dce792eSToby Isaac         roffset[p - pStart] += interleave_basis ? 1 : subdof;
4512dce792eSToby Isaac       }
4522dce792eSToby Isaac       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, subNb, rows, PETSC_OWN_POINTER, &is_row));
4532dce792eSToby Isaac     }
4542dce792eSToby Isaac 
4552dce792eSToby Isaac     sNpoints = 0;
4562dce792eSToby Isaac     if (quad) PetscCall(PetscQuadratureGetData(quad, NULL, NULL, &sNpoints, NULL, NULL));
4572dce792eSToby Isaac     sNcols = sNpoints * sNc;
4582dce792eSToby Isaac 
4592dce792eSToby Isaac     if (!concatenate) {
4602dce792eSToby Isaac       PetscCall(ISCreateStride(PETSC_COMM_SELF, sNcols, coffset, 1, &is_col));
4612dce792eSToby Isaac       coffset += uniform_points ? 0 : sNcols;
4622dce792eSToby Isaac     } else {
4632dce792eSToby Isaac       if (uniform_points && interleave_components) {
4642dce792eSToby Isaac         PetscCall(ISCreateStride(PETSC_COMM_SELF, sNcols, coffset, Ns, &is_col));
4652dce792eSToby Isaac         coffset += 1;
4662dce792eSToby Isaac       } else {
4672dce792eSToby Isaac         PetscInt *cols;
4682dce792eSToby Isaac 
4692dce792eSToby Isaac         PetscCall(PetscMalloc1(sNcols, &cols));
4702dce792eSToby Isaac         for (PetscInt p = 0, r = 0; p < sNpoints; p++) {
471*ac530a7eSPierre Jolivet           for (PetscInt c = 0; c < sNc; c++, r++) cols[r] = coffset + p * Nc + c;
4722dce792eSToby Isaac         }
4732dce792eSToby Isaac         coffset += uniform_points ? sNc : Nc * sNpoints + sNc;
4742dce792eSToby Isaac         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, sNcols, cols, PETSC_OWN_POINTER, &is_col));
4752dce792eSToby Isaac       }
4762dce792eSToby Isaac     }
4772dce792eSToby Isaac     PetscCall(ISLocalToGlobalMappingCreateIS(is_row, &map_row[s]));
4782dce792eSToby Isaac     PetscCall(ISLocalToGlobalMappingCreateIS(is_col, &map_col[s]));
4792dce792eSToby Isaac     PetscCall(ISDestroy(&is_row));
4802dce792eSToby Isaac     PetscCall(ISDestroy(&is_col));
4812dce792eSToby Isaac   }
4822dce792eSToby Isaac   PetscCall(PetscFree(roffset));
4832dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
4842dce792eSToby Isaac }
4852dce792eSToby Isaac 
PetscDualSpaceCreateFacetSubspace_Sum(PetscDualSpace sp,PetscInt f,PetscDualSpace * bdsp)4862dce792eSToby Isaac static PetscErrorCode PetscDualSpaceCreateFacetSubspace_Sum(PetscDualSpace sp, PetscInt f, PetscDualSpace *bdsp)
4872dce792eSToby Isaac {
4882dce792eSToby Isaac   PetscInt       k, Nc, Nk, Nknew, Ncnew, Ns;
4892dce792eSToby Isaac   PetscInt       dim, pointDim = -1;
4902dce792eSToby Isaac   PetscInt       depth, Ncopies;
4912dce792eSToby Isaac   PetscBool      any;
4922dce792eSToby Isaac   DM             dm, K;
4932dce792eSToby Isaac   DMPolytopeType ct;
4942dce792eSToby Isaac 
4952dce792eSToby Isaac   PetscFunctionBegin;
4962dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetNumSubspaces(sp, &Ns));
4972dce792eSToby Isaac   any = PETSC_FALSE;
4982dce792eSToby Isaac   for (PetscInt s = 0; s < Ns; s++) {
4992dce792eSToby Isaac     PetscDualSpace subsp, subpointsp;
5002dce792eSToby Isaac 
5012dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
5022dce792eSToby Isaac     PetscCall(PetscDualSpaceGetPointSubspace(subsp, f, &subpointsp));
5032dce792eSToby Isaac     if (subpointsp) any = PETSC_TRUE;
5042dce792eSToby Isaac   }
5052dce792eSToby Isaac   if (!any) {
5062dce792eSToby Isaac     *bdsp = NULL;
5072dce792eSToby Isaac     PetscFunctionReturn(PETSC_SUCCESS);
5082dce792eSToby Isaac   }
5092dce792eSToby Isaac   PetscCall(PetscDualSpaceGetDM(sp, &dm));
5102dce792eSToby Isaac   PetscCall(DMGetDimension(dm, &dim));
5112dce792eSToby Isaac   PetscCall(DMPlexGetDepth(dm, &depth));
5122dce792eSToby Isaac   PetscCall(PetscDualSpaceDuplicate(sp, bdsp));
5132dce792eSToby Isaac   PetscCheck((depth == dim) || (depth == 1), PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported interpolation state of reference element");
5142dce792eSToby Isaac   pointDim = (depth == dim) ? (dim - 1) : 0;
5152dce792eSToby Isaac   PetscCall(DMPlexGetCellType(dm, f, &ct));
5162dce792eSToby Isaac   PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, ct, &K));
5172dce792eSToby Isaac   PetscCall(PetscDualSpaceSetDM(*bdsp, K));
5182dce792eSToby Isaac   PetscCall(PetscDualSpaceGetFormDegree(sp, &k));
5192dce792eSToby Isaac   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
5202dce792eSToby Isaac   PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(k), &Nk));
5212dce792eSToby Isaac   Ncopies = Nc / Nk;
5222dce792eSToby Isaac   PetscCall(PetscDTBinomialInt(pointDim, PetscAbsInt(k), &Nknew));
5232dce792eSToby Isaac   Ncnew = Nknew * Ncopies;
5242dce792eSToby Isaac   PetscCall(PetscDualSpaceSetNumComponents(*bdsp, Ncnew));
5252dce792eSToby Isaac   for (PetscInt s = 0; s < Ns; s++) {
5262dce792eSToby Isaac     PetscDualSpace subsp, subpointsp;
5272dce792eSToby Isaac 
5282dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
5292dce792eSToby Isaac     PetscCall(PetscDualSpaceGetPointSubspace(subsp, f, &subpointsp));
5302dce792eSToby Isaac     if (subpointsp) {
5312dce792eSToby Isaac       PetscCall(PetscObjectReference((PetscObject)subpointsp));
5322dce792eSToby Isaac     } else {
5332dce792eSToby Isaac       // make an empty dual space
5342dce792eSToby Isaac       PetscInt subNc, subNcopies, subpointNc;
5352dce792eSToby Isaac 
5362dce792eSToby Isaac       PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)subsp), &subpointsp));
5372dce792eSToby Isaac       PetscCall(PetscDualSpaceGetNumComponents(subsp, &subNc));
5382dce792eSToby Isaac       subNcopies = subNc / Nk;
5392dce792eSToby Isaac       subpointNc = subNcopies * Nknew;
5402dce792eSToby Isaac       PetscCall(PetscDualSpaceSetType(subpointsp, PETSCDUALSPACESIMPLE));
5412dce792eSToby Isaac       PetscCall(PetscDualSpaceSimpleSetDimension(subpointsp, 0));
5422dce792eSToby Isaac       PetscCall(PetscDualSpaceSetFormDegree(subpointsp, k));
5432dce792eSToby Isaac       PetscCall(PetscDualSpaceSetNumComponents(subpointsp, subpointNc));
5442dce792eSToby Isaac     }
5452dce792eSToby Isaac     // K should be equal to subpointsp->dm, but we want it to be exactly the same
5462dce792eSToby Isaac     PetscCall(PetscObjectReference((PetscObject)K));
5472dce792eSToby Isaac     PetscCall(DMDestroy(&subpointsp->dm));
5482dce792eSToby Isaac     subpointsp->dm = K;
5492dce792eSToby Isaac     PetscCall(PetscDualSpaceSetUp(subpointsp));
5502dce792eSToby Isaac     PetscCall(PetscDualSpaceSumSetSubspace(*bdsp, s, subpointsp));
5512dce792eSToby Isaac     PetscCall(PetscDualSpaceDestroy(&subpointsp));
5522dce792eSToby Isaac   }
5532dce792eSToby Isaac   PetscCall(DMDestroy(&K));
5542dce792eSToby Isaac   PetscCall(PetscDualSpaceSetUp(*bdsp));
5552dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
5562dce792eSToby Isaac }
5572dce792eSToby Isaac 
PetscDualSpaceSumIsUniform(PetscDualSpace sp,PetscBool * is_uniform)5582dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumIsUniform(PetscDualSpace sp, PetscBool *is_uniform)
5592dce792eSToby Isaac {
5602dce792eSToby Isaac   PetscDualSpace_Sum *sum     = (PetscDualSpace_Sum *)sp->data;
5612dce792eSToby Isaac   PetscBool           uniform = PETSC_TRUE;
5622dce792eSToby Isaac 
5632dce792eSToby Isaac   PetscFunctionBegin;
5642dce792eSToby Isaac   for (PetscInt s = 1; s < sum->numSumSpaces; s++) {
5652dce792eSToby Isaac     if (sum->sumspaces[s] != sum->sumspaces[0]) {
5662dce792eSToby Isaac       uniform = PETSC_FALSE;
5672dce792eSToby Isaac       break;
5682dce792eSToby Isaac     }
5692dce792eSToby Isaac   }
5702dce792eSToby Isaac   *is_uniform = uniform;
5712dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
5722dce792eSToby Isaac }
5732dce792eSToby Isaac 
PetscDualSpaceGetSymmetries_Sum(PetscDualSpace sp,const PetscInt **** perms,const PetscScalar **** flips)5742dce792eSToby Isaac static PetscErrorCode PetscDualSpaceGetSymmetries_Sum(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips)
5752dce792eSToby Isaac {
5762dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data;
5772dce792eSToby Isaac 
5782dce792eSToby Isaac   PetscFunctionBegin;
5792dce792eSToby Isaac   if (!sum->symComputed) {
5802dce792eSToby Isaac     PetscInt       Ns;
5812dce792eSToby Isaac     PetscBool      any_perms = PETSC_FALSE;
5822dce792eSToby Isaac     PetscBool      any_flips = PETSC_FALSE;
5832dce792eSToby Isaac     PetscInt    ***symperms  = NULL;
5842dce792eSToby Isaac     PetscScalar ***symflips  = NULL;
5852dce792eSToby Isaac 
5862dce792eSToby Isaac     sum->symComputed = PETSC_TRUE;
5872dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetNumSubspaces(sp, &Ns));
5882dce792eSToby Isaac     for (PetscInt s = 0; s < Ns; s++) {
5892dce792eSToby Isaac       PetscDualSpace       subsp;
5902dce792eSToby Isaac       const PetscInt    ***sub_perms;
5912dce792eSToby Isaac       const PetscScalar ***sub_flips;
5922dce792eSToby Isaac 
5932dce792eSToby Isaac       PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
5942dce792eSToby Isaac       PetscCall(PetscDualSpaceGetSymmetries(subsp, &sub_perms, &sub_flips));
5952dce792eSToby Isaac       if (sub_perms) any_perms = PETSC_TRUE;
5962dce792eSToby Isaac       if (sub_flips) any_flips = PETSC_TRUE;
5972dce792eSToby Isaac     }
5982dce792eSToby Isaac     if (any_perms || any_flips) {
5992dce792eSToby Isaac       DM       K;
6002dce792eSToby Isaac       PetscInt pStart, pEnd, numPoints;
6012dce792eSToby Isaac       PetscInt spintdim;
6022dce792eSToby Isaac 
6032dce792eSToby Isaac       PetscCall(PetscDualSpaceGetDM(sp, &K));
6042dce792eSToby Isaac       PetscCall(DMPlexGetChart(K, &pStart, &pEnd));
6052dce792eSToby Isaac       numPoints = pEnd - pStart;
6062dce792eSToby Isaac       PetscCall(PetscCalloc1(numPoints, &symperms));
6072dce792eSToby Isaac       PetscCall(PetscCalloc1(numPoints, &symflips));
6082dce792eSToby Isaac       PetscCall(PetscDualSpaceGetBoundarySymmetries_Internal(sp, symperms, symflips));
6092dce792eSToby Isaac       // get interior symmetries
6102dce792eSToby Isaac       PetscCall(PetscDualSpaceGetInteriorDimension(sp, &spintdim));
6112dce792eSToby Isaac       if (spintdim) {
6122dce792eSToby Isaac         PetscInt       groupSize;
6132dce792eSToby Isaac         PetscInt     **cellPerms;
6142dce792eSToby Isaac         PetscScalar  **cellFlips;
6152dce792eSToby Isaac         DMPolytopeType ct;
6162dce792eSToby Isaac 
6172dce792eSToby Isaac         PetscCall(DMPlexGetCellType(K, 0, &ct));
61885036b15SMatthew G. Knepley         groupSize       = DMPolytopeTypeGetNumArrangements(ct);
6192dce792eSToby Isaac         sum->numSelfSym = groupSize;
6202dce792eSToby Isaac         sum->selfSymOff = groupSize / 2;
6212dce792eSToby Isaac         PetscCall(PetscCalloc1(groupSize, &cellPerms));
6222dce792eSToby Isaac         PetscCall(PetscCalloc1(groupSize, &cellFlips));
6232dce792eSToby Isaac         symperms[0] = &cellPerms[groupSize / 2];
6242dce792eSToby Isaac         symflips[0] = &cellFlips[groupSize / 2];
6252dce792eSToby Isaac         for (PetscInt o = -groupSize / 2; o < groupSize / 2; o++) {
6262dce792eSToby Isaac           PetscBool any_o_perms = PETSC_FALSE;
6272dce792eSToby Isaac           PetscBool any_o_flips = PETSC_FALSE;
6282dce792eSToby Isaac 
6292dce792eSToby Isaac           for (PetscInt s = 0; s < Ns; s++) {
6302dce792eSToby Isaac             PetscDualSpace       subsp;
6312dce792eSToby Isaac             const PetscInt    ***sub_perms;
6322dce792eSToby Isaac             const PetscScalar ***sub_flips;
6332dce792eSToby Isaac 
6342dce792eSToby Isaac             PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
6352dce792eSToby Isaac             PetscCall(PetscDualSpaceGetSymmetries(subsp, &sub_perms, &sub_flips));
6362dce792eSToby Isaac             if (sub_perms && sub_perms[0] && sub_perms[0][o]) any_o_perms = PETSC_TRUE;
6372dce792eSToby Isaac             if (sub_flips && sub_flips[0] && sub_flips[0][o]) any_o_flips = PETSC_TRUE;
6382dce792eSToby Isaac           }
6392dce792eSToby Isaac           if (any_o_perms) {
6402dce792eSToby Isaac             PetscInt *o_perm;
6412dce792eSToby Isaac 
6422dce792eSToby Isaac             PetscCall(PetscMalloc1(spintdim, &o_perm));
6432dce792eSToby Isaac             for (PetscInt i = 0; i < spintdim; i++) o_perm[i] = i;
6442dce792eSToby Isaac             for (PetscInt s = 0; s < Ns; s++) {
6452dce792eSToby Isaac               PetscDualSpace       subsp;
6462dce792eSToby Isaac               const PetscInt    ***sub_perms;
6472dce792eSToby Isaac               const PetscScalar ***sub_flips;
6482dce792eSToby Isaac 
6492dce792eSToby Isaac               PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
6502dce792eSToby Isaac               PetscCall(PetscDualSpaceGetSymmetries(subsp, &sub_perms, &sub_flips));
6512dce792eSToby Isaac               if (sub_perms && sub_perms[0] && sub_perms[0][o]) {
6522dce792eSToby Isaac                 PetscInt  subspdim;
6532dce792eSToby Isaac                 PetscInt *range, *domain;
6542dce792eSToby Isaac                 PetscInt *range_mapped, *domain_mapped;
6552dce792eSToby Isaac 
6562dce792eSToby Isaac                 PetscCall(PetscDualSpaceGetInteriorDimension(subsp, &subspdim));
6572dce792eSToby Isaac                 PetscCall(PetscMalloc4(subspdim, &range, subspdim, &range_mapped, subspdim, &domain, subspdim, &domain_mapped));
6582dce792eSToby Isaac                 for (PetscInt i = 0; i < subspdim; i++) domain[i] = i;
6592dce792eSToby Isaac                 PetscCall(PetscArraycpy(range, sub_perms[0][o], subspdim));
6602dce792eSToby Isaac                 PetscCall(ISLocalToGlobalMappingApply(sum->int_rows[s], subspdim, domain, domain_mapped));
6612dce792eSToby Isaac                 PetscCall(ISLocalToGlobalMappingApply(sum->int_rows[s], subspdim, range, range_mapped));
6622dce792eSToby Isaac                 for (PetscInt i = 0; i < subspdim; i++) o_perm[domain_mapped[i]] = range_mapped[i];
6632dce792eSToby Isaac                 PetscCall(PetscFree4(range, range_mapped, domain, domain_mapped));
6642dce792eSToby Isaac               }
6652dce792eSToby Isaac             }
6662dce792eSToby Isaac             symperms[0][o] = o_perm;
6672dce792eSToby Isaac           }
6682dce792eSToby Isaac           if (any_o_flips) {
6692dce792eSToby Isaac             PetscScalar *o_flip;
6702dce792eSToby Isaac 
6712dce792eSToby Isaac             PetscCall(PetscMalloc1(spintdim, &o_flip));
6722dce792eSToby Isaac             for (PetscInt i = 0; i < spintdim; i++) o_flip[i] = 1.0;
6732dce792eSToby Isaac             for (PetscInt s = 0; s < Ns; s++) {
6742dce792eSToby Isaac               PetscDualSpace       subsp;
6752dce792eSToby Isaac               const PetscInt    ***sub_perms;
6762dce792eSToby Isaac               const PetscScalar ***sub_flips;
6772dce792eSToby Isaac 
6782dce792eSToby Isaac               PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
6792dce792eSToby Isaac               PetscCall(PetscDualSpaceGetSymmetries(subsp, &sub_perms, &sub_flips));
6802dce792eSToby Isaac               if (sub_perms && sub_perms[0] && sub_perms[0][o]) {
6812dce792eSToby Isaac                 PetscInt  subspdim;
6822dce792eSToby Isaac                 PetscInt *domain;
6832dce792eSToby Isaac                 PetscInt *domain_mapped;
6842dce792eSToby Isaac 
6852dce792eSToby Isaac                 PetscCall(PetscDualSpaceGetInteriorDimension(subsp, &subspdim));
6862dce792eSToby Isaac                 PetscCall(PetscMalloc2(subspdim, &domain, subspdim, &domain_mapped));
6872dce792eSToby Isaac                 for (PetscInt i = 0; i < subspdim; i++) domain[i] = i;
6882dce792eSToby Isaac                 PetscCall(ISLocalToGlobalMappingApply(sum->int_rows[s], subspdim, domain, domain_mapped));
6892dce792eSToby Isaac                 for (PetscInt i = 0; i < subspdim; i++) o_flip[domain_mapped[i]] = sub_perms[0][o][i];
6902dce792eSToby Isaac                 PetscCall(PetscFree2(domain, domain_mapped));
6912dce792eSToby Isaac               }
6922dce792eSToby Isaac             }
6932dce792eSToby Isaac             symflips[0][o] = o_flip;
6942dce792eSToby Isaac           }
6952dce792eSToby Isaac         }
6962dce792eSToby Isaac         {
6972dce792eSToby Isaac           PetscBool any_perms = PETSC_FALSE;
6982dce792eSToby Isaac           PetscBool any_flips = PETSC_FALSE;
6992dce792eSToby Isaac           for (PetscInt o = -groupSize / 2; o < groupSize / 2; o++) {
7002dce792eSToby Isaac             if (symperms[0][o]) any_perms = PETSC_TRUE;
7012dce792eSToby Isaac             if (symflips[0][o]) any_flips = PETSC_TRUE;
7022dce792eSToby Isaac           }
7032dce792eSToby Isaac           if (!any_perms) {
7042dce792eSToby Isaac             PetscCall(PetscFree(cellPerms));
7052dce792eSToby Isaac             symperms[0] = NULL;
7062dce792eSToby Isaac           }
7072dce792eSToby Isaac           if (!any_flips) {
7082dce792eSToby Isaac             PetscCall(PetscFree(cellFlips));
7092dce792eSToby Isaac             symflips[0] = NULL;
7102dce792eSToby Isaac           }
7112dce792eSToby Isaac         }
7122dce792eSToby Isaac       }
7132dce792eSToby Isaac       if (!any_perms) {
7142dce792eSToby Isaac         PetscCall(PetscFree(symperms));
7152dce792eSToby Isaac         symperms = NULL;
7162dce792eSToby Isaac       }
7172dce792eSToby Isaac       if (!any_flips) {
7182dce792eSToby Isaac         PetscCall(PetscFree(symflips));
7192dce792eSToby Isaac         symflips = NULL;
7202dce792eSToby Isaac       }
7212dce792eSToby Isaac     }
7222dce792eSToby Isaac     sum->symperms = symperms;
7232dce792eSToby Isaac     sum->symflips = symflips;
7242dce792eSToby Isaac   }
7252dce792eSToby Isaac   if (perms) *perms = (const PetscInt ***)sum->symperms;
7262dce792eSToby Isaac   if (flips) *flips = (const PetscScalar ***)sum->symflips;
7272dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
7282dce792eSToby Isaac }
7292dce792eSToby Isaac 
PetscDualSpaceSetUp_Sum(PetscDualSpace sp)7302dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSetUp_Sum(PetscDualSpace sp)
7312dce792eSToby Isaac {
7322dce792eSToby Isaac   PetscDualSpace_Sum *sum         = (PetscDualSpace_Sum *)sp->data;
7332dce792eSToby Isaac   PetscBool           concatenate = PETSC_TRUE;
7342dce792eSToby Isaac   PetscBool           uniform;
7352dce792eSToby Isaac   PetscInt            Ns, Nc, i, sum_Nc = 0;
7362dce792eSToby Isaac   PetscInt            minNc, maxNc;
7372dce792eSToby Isaac   PetscInt            minForm, maxForm, cdim, depth;
7382dce792eSToby Isaac   DM                  K;
7392dce792eSToby Isaac   PetscQuadrature    *all_quads = NULL;
7402dce792eSToby Isaac   PetscQuadrature    *int_quads = NULL;
7412dce792eSToby Isaac   Mat                *all_mats  = NULL;
7422dce792eSToby Isaac   Mat                *int_mats  = NULL;
7432dce792eSToby Isaac 
7442dce792eSToby Isaac   PetscFunctionBegin;
745371d2eb7SMartin Diehl   if (sum->setupcalled) PetscFunctionReturn(PETSC_SUCCESS);
746371d2eb7SMartin Diehl   sum->setupcalled = PETSC_TRUE;
7472dce792eSToby Isaac 
7482dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetNumSubspaces(sp, &Ns));
7492dce792eSToby Isaac   PetscCheck(Ns >= 0, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Cannot have %" PetscInt_FMT " subspaces", Ns);
7502dce792eSToby Isaac 
7512dce792eSToby Isaac   // step 1: make sure they share a DM
7522dce792eSToby Isaac   PetscCall(PetscDualSpaceGetDM(sp, &K));
7532dce792eSToby Isaac   PetscCall(DMGetCoordinateDim(K, &cdim));
7542dce792eSToby Isaac   PetscCall(DMPlexGetDepth(K, &depth));
7552dce792eSToby Isaac   PetscCall(PetscDualSpaceSumIsUniform(sp, &sp->uniform));
7562dce792eSToby Isaac   uniform = sp->uniform;
7572dce792eSToby Isaac   {
7582dce792eSToby Isaac     for (PetscInt s = 0; s < Ns; s++) {
7592dce792eSToby Isaac       PetscDualSpace subsp;
7602dce792eSToby Isaac       DM             sub_K;
7612dce792eSToby Isaac 
7622dce792eSToby Isaac       PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
7632dce792eSToby Isaac       PetscCall(PetscDualSpaceSetUp(subsp));
7642dce792eSToby Isaac       PetscCall(PetscDualSpaceGetDM(subsp, &sub_K));
7652dce792eSToby Isaac       if (s == 0 && K == NULL) {
7662dce792eSToby Isaac         PetscCall(PetscDualSpaceSetDM(sp, sub_K));
7672dce792eSToby Isaac         K = sub_K;
7682dce792eSToby Isaac       }
769835f2295SStefano Zampini       PetscCheck(sub_K == K, PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Subspace %" PetscInt_FMT " does not have the same DM as the sum space", s);
7702dce792eSToby Isaac     }
7712dce792eSToby Isaac   }
7722dce792eSToby Isaac 
7732dce792eSToby Isaac   // step 2: count components
7742dce792eSToby Isaac   PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc));
7752dce792eSToby Isaac   PetscCall(PetscDualSpaceSumGetConcatenate(sp, &concatenate));
7762dce792eSToby Isaac   minNc   = Nc;
7772dce792eSToby Isaac   maxNc   = Nc;
7782dce792eSToby Isaac   minForm = sp->k == PETSC_FORM_DEGREE_UNDEFINED ? PETSC_INT_MAX : sp->k;
7792dce792eSToby Isaac   maxForm = sp->k == PETSC_FORM_DEGREE_UNDEFINED ? PETSC_INT_MIN : sp->k;
7802dce792eSToby Isaac   for (i = 0; i < Ns; ++i) {
7812dce792eSToby Isaac     PetscInt       sNc, formDegree;
7822dce792eSToby Isaac     PetscDualSpace si;
7832dce792eSToby Isaac 
7842dce792eSToby Isaac     PetscCall(PetscDualSpaceSumGetSubspace(sp, i, &si));
7852dce792eSToby Isaac     PetscCall(PetscDualSpaceSetUp(si));
7862dce792eSToby Isaac     PetscCall(PetscDualSpaceGetNumComponents(si, &sNc));
7872dce792eSToby 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");
7882dce792eSToby Isaac     minNc = PetscMin(minNc, sNc);
7892dce792eSToby Isaac     maxNc = PetscMax(maxNc, sNc);
7902dce792eSToby Isaac     sum_Nc += sNc;
7912dce792eSToby Isaac     PetscCall(PetscDualSpaceGetFormDegree(si, &formDegree));
7922dce792eSToby Isaac     minForm = PetscMin(minForm, formDegree);
7932dce792eSToby Isaac     maxForm = PetscMax(maxForm, formDegree);
7942dce792eSToby Isaac   }
7952dce792eSToby Isaac   sp->k = (minForm != maxForm) ? PETSC_FORM_DEGREE_UNDEFINED : minForm;
7962dce792eSToby Isaac 
7972dce792eSToby 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);
7982dce792eSToby 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.");
7992dce792eSToby Isaac 
8002dce792eSToby Isaac   /* A PetscDualSpace should have a fixed number of components, but
8012dce792eSToby Isaac      if the spaces we are combining have different form degrees, they will not
8022dce792eSToby Isaac      have the same number of components on subcomponents of the boundary,
8032dce792eSToby Isaac      so we do not try to create boundary dual spaces in this case */
8042dce792eSToby Isaac   if (sp->k != PETSC_FORM_DEGREE_UNDEFINED && depth > 0) {
8052dce792eSToby Isaac     PetscInt  pStart, pEnd;
8062dce792eSToby Isaac     PetscInt *pStratStart, *pStratEnd;
8072dce792eSToby Isaac 
8082dce792eSToby Isaac     PetscCall(DMPlexGetDepth(K, &depth));
8092dce792eSToby Isaac     PetscCall(DMPlexGetChart(K, &pStart, &pEnd));
810f4f49eeaSPierre Jolivet     PetscCall(PetscCalloc1(pEnd, &sp->pointSpaces));
8112dce792eSToby Isaac     PetscCall(PetscMalloc2(depth + 1, &pStratStart, depth + 1, &pStratEnd));
8122dce792eSToby Isaac     for (PetscInt d = 0; d <= depth; ++d) PetscCall(DMPlexGetDepthStratum(K, d, &pStratStart[d], &pStratEnd[d]));
8132dce792eSToby Isaac 
8142dce792eSToby Isaac     for (PetscInt p = pStratStart[depth - 1]; p < pStratEnd[depth - 1]; p++) { /* calculate the facet dual spaces */
8152dce792eSToby Isaac       PetscReal      v0[3];
8162dce792eSToby Isaac       DMPolytopeType ptype;
8172dce792eSToby Isaac       PetscReal      J[9], detJ;
8182dce792eSToby Isaac       PetscInt       q;
8192dce792eSToby Isaac 
8202dce792eSToby Isaac       PetscCall(DMPlexComputeCellGeometryAffineFEM(K, p, v0, J, NULL, &detJ));
8212dce792eSToby Isaac       PetscCall(DMPlexGetCellType(K, p, &ptype));
8222dce792eSToby Isaac 
8232dce792eSToby Isaac       /* compare to previous facets: if computed, reference that dualspace */
8242dce792eSToby Isaac       for (q = pStratStart[depth - 1]; q < p; q++) {
8252dce792eSToby Isaac         DMPolytopeType qtype;
8262dce792eSToby Isaac 
8272dce792eSToby Isaac         PetscCall(DMPlexGetCellType(K, q, &qtype));
8282dce792eSToby Isaac         if (qtype == ptype) break;
8292dce792eSToby Isaac       }
8302dce792eSToby Isaac       if (q < p) { /* this facet has the same dual space as that one */
8312dce792eSToby Isaac         PetscCall(PetscObjectReference((PetscObject)sp->pointSpaces[q]));
8322dce792eSToby Isaac         sp->pointSpaces[p] = sp->pointSpaces[q];
8332dce792eSToby Isaac         continue;
8342dce792eSToby Isaac       }
8352dce792eSToby Isaac       /* if not, recursively compute this dual space */
8362dce792eSToby Isaac       PetscCall(PetscDualSpaceCreateFacetSubspace_Sum(sp, p, &sp->pointSpaces[p]));
8372dce792eSToby Isaac     }
8382dce792eSToby Isaac     for (PetscInt h = 2; h <= depth; h++) { /* get the higher subspaces from the facet subspaces */
8392dce792eSToby Isaac       PetscInt hd = depth - h;
8402dce792eSToby Isaac 
8412dce792eSToby Isaac       for (PetscInt p = pStratStart[hd]; p < pStratEnd[hd]; p++) {
8422dce792eSToby Isaac         PetscInt        suppSize;
8432dce792eSToby Isaac         const PetscInt *supp;
8442dce792eSToby Isaac         DM              qdm;
8452dce792eSToby Isaac         PetscDualSpace  qsp, psp;
8462dce792eSToby Isaac         PetscInt        c, coneSize, q;
8472dce792eSToby Isaac         const PetscInt *cone;
8482dce792eSToby Isaac         const PetscInt *refCone;
8492dce792eSToby Isaac 
8502dce792eSToby Isaac         PetscCall(DMPlexGetSupportSize(K, p, &suppSize));
8512dce792eSToby Isaac         PetscCall(DMPlexGetSupport(K, p, &supp));
8522dce792eSToby Isaac         q   = supp[0];
8532dce792eSToby Isaac         qsp = sp->pointSpaces[q];
8542dce792eSToby Isaac         PetscCall(DMPlexGetConeSize(K, q, &coneSize));
8552dce792eSToby Isaac         PetscCall(DMPlexGetCone(K, q, &cone));
8562dce792eSToby Isaac         for (c = 0; c < coneSize; c++)
8572dce792eSToby Isaac           if (cone[c] == p) break;
8582dce792eSToby Isaac         PetscCheck(c != coneSize, PetscObjectComm((PetscObject)K), PETSC_ERR_PLIB, "cone/support mismatch");
8592dce792eSToby Isaac         if (!qsp) {
8602dce792eSToby Isaac           sp->pointSpaces[p] = NULL;
8612dce792eSToby Isaac           continue;
8622dce792eSToby Isaac         }
8632dce792eSToby Isaac         PetscCall(PetscDualSpaceGetDM(qsp, &qdm));
8642dce792eSToby Isaac         PetscCall(DMPlexGetCone(qdm, 0, &refCone));
8652dce792eSToby Isaac         /* get the equivalent dual space from the support dual space */
8662dce792eSToby Isaac         PetscCall(PetscDualSpaceGetPointSubspace(qsp, refCone[c], &psp));
8672dce792eSToby Isaac         PetscCall(PetscObjectReference((PetscObject)psp));
8682dce792eSToby Isaac         sp->pointSpaces[p] = psp;
8692dce792eSToby Isaac       }
8702dce792eSToby Isaac     }
8712dce792eSToby Isaac     PetscCall(PetscFree2(pStratStart, pStratEnd));
8722dce792eSToby Isaac   }
8732dce792eSToby Isaac 
8742dce792eSToby Isaac   sum->uniform = uniform;
8752dce792eSToby Isaac   PetscCall(PetscCalloc1(Ns, &sum->all_rows));
8762dce792eSToby Isaac   PetscCall(PetscCalloc1(Ns, &sum->all_cols));
8772dce792eSToby Isaac   PetscCall(PetscCalloc1(Ns, &sum->int_rows));
8782dce792eSToby Isaac   PetscCall(PetscCalloc1(Ns, &sum->int_cols));
8792dce792eSToby Isaac   PetscCall(PetscMalloc4(Ns, &all_quads, Ns, &all_mats, Ns, &int_quads, Ns, &int_mats));
8802dce792eSToby Isaac   {
8812dce792eSToby Isaac     // test for uniform all points and uniform interior points
8822dce792eSToby Isaac     PetscBool       uniform_all         = PETSC_TRUE;
8832dce792eSToby Isaac     PetscBool       uniform_interior    = PETSC_TRUE;
8842dce792eSToby Isaac     PetscQuadrature quad_all_first      = NULL;
8852dce792eSToby Isaac     PetscQuadrature quad_interior_first = NULL;
8862dce792eSToby Isaac     PetscInt        pStart, pEnd;
8872dce792eSToby Isaac 
8882dce792eSToby Isaac     PetscCall(DMPlexGetChart(K, &pStart, &pEnd));
8892dce792eSToby Isaac     PetscCall(PetscDualSpaceSectionCreate_Internal(sp, &sp->pointSection));
8902dce792eSToby Isaac 
8912dce792eSToby Isaac     for (PetscInt p = pStart; p < pEnd; p++) {
8922dce792eSToby Isaac       PetscInt full_dof = 0;
8932dce792eSToby Isaac 
8942dce792eSToby Isaac       for (PetscInt s = 0; s < Ns; s++) {
8952dce792eSToby Isaac         PetscDualSpace subsp;
8962dce792eSToby Isaac         PetscSection   subsection;
8972dce792eSToby Isaac         PetscInt       dof;
8982dce792eSToby Isaac 
8992dce792eSToby Isaac         PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
9002dce792eSToby Isaac         PetscCall(PetscDualSpaceGetSection(subsp, &subsection));
9012dce792eSToby Isaac         PetscCall(PetscSectionGetDof(subsection, p, &dof));
9022dce792eSToby Isaac         full_dof += dof;
9032dce792eSToby Isaac       }
9042dce792eSToby Isaac       PetscCall(PetscSectionSetDof(sp->pointSection, p, full_dof));
9052dce792eSToby Isaac     }
9062dce792eSToby Isaac     PetscCall(PetscDualSpaceSectionSetUp_Internal(sp, sp->pointSection));
9072dce792eSToby Isaac 
9082dce792eSToby Isaac     for (PetscInt s = 0; s < Ns; s++) {
9092dce792eSToby Isaac       PetscDualSpace  subsp;
9102dce792eSToby Isaac       PetscQuadrature subquad_all;
9112dce792eSToby Isaac       PetscQuadrature subquad_interior;
9122dce792eSToby Isaac       Mat             submat_all;
9132dce792eSToby Isaac       Mat             submat_interior;
9142dce792eSToby Isaac 
9152dce792eSToby Isaac       PetscCall(PetscDualSpaceSumGetSubspace(sp, s, &subsp));
9162dce792eSToby Isaac       PetscCall(PetscDualSpaceGetAllData(subsp, &subquad_all, &submat_all));
9172dce792eSToby Isaac       PetscCall(PetscDualSpaceGetInteriorData(subsp, &subquad_interior, &submat_interior));
9182dce792eSToby Isaac       if (!s) {
9192dce792eSToby Isaac         quad_all_first      = subquad_all;
9202dce792eSToby Isaac         quad_interior_first = subquad_interior;
9212dce792eSToby Isaac       } else {
9222dce792eSToby Isaac         if (subquad_all != quad_all_first) uniform_all = PETSC_FALSE;
9232dce792eSToby Isaac         if (subquad_interior != quad_interior_first) uniform_interior = PETSC_FALSE;
9242dce792eSToby Isaac       }
9252dce792eSToby Isaac       all_quads[s] = subquad_all;
9262dce792eSToby Isaac       int_quads[s] = subquad_interior;
9272dce792eSToby Isaac       all_mats[s]  = submat_all;
9282dce792eSToby Isaac       int_mats[s]  = submat_interior;
9292dce792eSToby Isaac     }
9302dce792eSToby Isaac     sum->uniform_all_points      = uniform_all;
9312dce792eSToby Isaac     sum->uniform_interior_points = uniform_interior;
9322dce792eSToby Isaac 
9332dce792eSToby Isaac     PetscCall(PetscDualSpaceSumCreateMappings(sp, PETSC_TRUE, uniform_interior, sum->int_rows, sum->int_cols));
9342dce792eSToby Isaac     PetscCall(PetscDualSpaceSumCreateQuadrature(Ns, cdim, uniform_interior, int_quads, &sp->intNodes));
9352dce792eSToby Isaac     if (sp->intNodes) PetscCall(PetscDualSpaceSumCreateMatrix(sp, int_mats, sum->int_rows, sum->int_cols, sp->intNodes, &sp->intMat));
9362dce792eSToby Isaac 
9372dce792eSToby Isaac     PetscCall(PetscDualSpaceSumCreateMappings(sp, PETSC_FALSE, uniform_all, sum->all_rows, sum->all_cols));
9382dce792eSToby Isaac     PetscCall(PetscDualSpaceSumCreateQuadrature(Ns, cdim, uniform_all, all_quads, &sp->allNodes));
9392dce792eSToby Isaac     if (sp->allNodes) PetscCall(PetscDualSpaceSumCreateMatrix(sp, all_mats, sum->all_rows, sum->all_cols, sp->allNodes, &sp->allMat));
9402dce792eSToby Isaac   }
9412dce792eSToby Isaac   PetscCall(PetscFree4(all_quads, all_mats, int_quads, int_mats));
9422dce792eSToby Isaac   PetscCall(PetscDualSpaceComputeFunctionalsFromAllData(sp));
9432dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
9442dce792eSToby Isaac }
9452dce792eSToby Isaac 
PetscDualSpaceSumView_Ascii(PetscDualSpace sp,PetscViewer v)9462dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumView_Ascii(PetscDualSpace sp, PetscViewer v)
9472dce792eSToby Isaac {
9482dce792eSToby Isaac   PetscDualSpace_Sum *sum         = (PetscDualSpace_Sum *)sp->data;
9492dce792eSToby Isaac   PetscBool           concatenate = sum->concatenate;
9502dce792eSToby Isaac   PetscInt            i, Ns = sum->numSumSpaces;
9512dce792eSToby Isaac 
9522dce792eSToby Isaac   PetscFunctionBegin;
9532dce792eSToby Isaac   if (concatenate) PetscCall(PetscViewerASCIIPrintf(v, "Sum dual space of %" PetscInt_FMT " concatenated subspaces%s\n", Ns, sum->uniform ? " (all identical)" : ""));
9542dce792eSToby Isaac   else PetscCall(PetscViewerASCIIPrintf(v, "Sum dual space of %" PetscInt_FMT " subspaces%s\n", Ns, sum->uniform ? " (all identical)" : ""));
9552dce792eSToby Isaac   for (i = 0; i < (sum->uniform ? (Ns > 0 ? 1 : 0) : Ns); ++i) {
9562dce792eSToby Isaac     PetscCall(PetscViewerASCIIPushTab(v));
9572dce792eSToby Isaac     PetscCall(PetscDualSpaceView(sum->sumspaces[i], v));
9582dce792eSToby Isaac     PetscCall(PetscViewerASCIIPopTab(v));
9592dce792eSToby Isaac   }
9602dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
9612dce792eSToby Isaac }
9622dce792eSToby Isaac 
PetscDualSpaceView_Sum(PetscDualSpace sp,PetscViewer viewer)9632dce792eSToby Isaac static PetscErrorCode PetscDualSpaceView_Sum(PetscDualSpace sp, PetscViewer viewer)
9642dce792eSToby Isaac {
9659f196a02SMartin Diehl   PetscBool isascii;
9662dce792eSToby Isaac 
9672dce792eSToby Isaac   PetscFunctionBegin;
9689f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
9699f196a02SMartin Diehl   if (isascii) PetscCall(PetscDualSpaceSumView_Ascii(sp, viewer));
9702dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
9712dce792eSToby Isaac }
9722dce792eSToby Isaac 
PetscDualSpaceDestroy_Sum(PetscDualSpace sp)9732dce792eSToby Isaac static PetscErrorCode PetscDualSpaceDestroy_Sum(PetscDualSpace sp)
9742dce792eSToby Isaac {
9752dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data;
9762dce792eSToby Isaac   PetscInt            i, Ns = sum->numSumSpaces;
9772dce792eSToby Isaac 
9782dce792eSToby Isaac   PetscFunctionBegin;
9792dce792eSToby Isaac   if (sum->symperms) {
9802dce792eSToby Isaac     PetscInt **selfSyms = sum->symperms[0];
9812dce792eSToby Isaac 
9822dce792eSToby Isaac     if (selfSyms) {
9832dce792eSToby Isaac       PetscInt i, **allocated = &selfSyms[-sum->selfSymOff];
9842dce792eSToby Isaac 
9852dce792eSToby Isaac       for (i = 0; i < sum->numSelfSym; i++) PetscCall(PetscFree(allocated[i]));
9862dce792eSToby Isaac       PetscCall(PetscFree(allocated));
9872dce792eSToby Isaac     }
9882dce792eSToby Isaac     PetscCall(PetscFree(sum->symperms));
9892dce792eSToby Isaac   }
9902dce792eSToby Isaac   if (sum->symflips) {
9912dce792eSToby Isaac     PetscScalar **selfSyms = sum->symflips[0];
9922dce792eSToby Isaac 
9932dce792eSToby Isaac     if (selfSyms) {
9942dce792eSToby Isaac       PetscInt      i;
9952dce792eSToby Isaac       PetscScalar **allocated = &selfSyms[-sum->selfSymOff];
9962dce792eSToby Isaac 
9972dce792eSToby Isaac       for (i = 0; i < sum->numSelfSym; i++) PetscCall(PetscFree(allocated[i]));
9982dce792eSToby Isaac       PetscCall(PetscFree(allocated));
9992dce792eSToby Isaac     }
10002dce792eSToby Isaac     PetscCall(PetscFree(sum->symflips));
10012dce792eSToby Isaac   }
10022dce792eSToby Isaac   for (i = 0; i < Ns; ++i) {
10032dce792eSToby Isaac     PetscCall(PetscDualSpaceDestroy(&sum->sumspaces[i]));
10042dce792eSToby Isaac     if (sum->all_rows) PetscCall(ISLocalToGlobalMappingDestroy(&sum->all_rows[i]));
10052dce792eSToby Isaac     if (sum->all_cols) PetscCall(ISLocalToGlobalMappingDestroy(&sum->all_cols[i]));
10062dce792eSToby Isaac     if (sum->int_rows) PetscCall(ISLocalToGlobalMappingDestroy(&sum->int_rows[i]));
10072dce792eSToby Isaac     if (sum->int_cols) PetscCall(ISLocalToGlobalMappingDestroy(&sum->int_cols[i]));
10082dce792eSToby Isaac   }
10092dce792eSToby Isaac   PetscCall(PetscFree(sum->sumspaces));
10102dce792eSToby Isaac   PetscCall(PetscFree(sum->all_rows));
10112dce792eSToby Isaac   PetscCall(PetscFree(sum->all_cols));
10122dce792eSToby Isaac   PetscCall(PetscFree(sum->int_rows));
10132dce792eSToby Isaac   PetscCall(PetscFree(sum->int_cols));
10142dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetSubspace_C", NULL));
10152dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetSubspace_C", NULL));
10162dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetNumSubspaces_C", NULL));
10172dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetNumSubspaces_C", NULL));
10182dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetConcatenate_C", NULL));
10192dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetConcatenate_C", NULL));
10202dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetInterleave_C", NULL));
10212dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetInterleave_C", NULL));
10222dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetContinuity_C", NULL));
10232dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetContinuity_C", NULL));
10242dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetMomentOrder_C", NULL));
10252dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetMomentOrder_C", NULL));
10262dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetNodeType_C", NULL));
10272dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetNodeType_C", NULL));
10282dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTensor_C", NULL));
10292dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTensor_C", NULL));
10302dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTrimmed_C", NULL));
10312dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTrimmed_C", NULL));
10322dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetUseMoments_C", NULL));
10332dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetUseMoments_C", NULL));
10342dce792eSToby Isaac   PetscCall(PetscFree(sum));
10352dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
10362dce792eSToby Isaac }
10372dce792eSToby Isaac 
10382dce792eSToby Isaac /*@
10392dce792eSToby Isaac   PetscDualSpaceSumSetInterleave - Set whether the basis functions and components of a uniform sum are interleaved
10402dce792eSToby Isaac 
10412dce792eSToby Isaac   Logically collective
10422dce792eSToby Isaac 
10432dce792eSToby Isaac   Input Parameters:
10442dce792eSToby Isaac + sp                    - a `PetscDualSpace` of type `PETSCDUALSPACESUM`
10452dce792eSToby Isaac . interleave_basis      - if `PETSC_TRUE`, the basis vectors of the subspaces are interleaved
10462dce792eSToby Isaac - interleave_components - if `PETSC_TRUE` and the space concatenates components (`PetscDualSpaceSumGetConcatenate()`),
10472dce792eSToby Isaac                           interleave the concatenated components
10482dce792eSToby Isaac 
10492dce792eSToby Isaac   Level: developer
10502dce792eSToby Isaac 
10512dce792eSToby Isaac .seealso: `PetscDualSpace`, `PETSCDUALSPACESUM`, `PETSCFEVECTOR`, `PetscDualSpaceSumGetInterleave()`
10522dce792eSToby Isaac @*/
PetscDualSpaceSumSetInterleave(PetscDualSpace sp,PetscBool interleave_basis,PetscBool interleave_components)10532dce792eSToby Isaac PetscErrorCode PetscDualSpaceSumSetInterleave(PetscDualSpace sp, PetscBool interleave_basis, PetscBool interleave_components)
10542dce792eSToby Isaac {
10552dce792eSToby Isaac   PetscFunctionBegin;
10562dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
10572dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumSetInterleave_C", (PetscDualSpace, PetscBool, PetscBool), (sp, interleave_basis, interleave_components));
10582dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
10592dce792eSToby Isaac }
10602dce792eSToby Isaac 
PetscDualSpaceSumSetInterleave_Sum(PetscDualSpace sp,PetscBool interleave_basis,PetscBool interleave_components)10612dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumSetInterleave_Sum(PetscDualSpace sp, PetscBool interleave_basis, PetscBool interleave_components)
10622dce792eSToby Isaac {
10632dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data;
10644d86920dSPierre Jolivet 
10652dce792eSToby Isaac   PetscFunctionBegin;
10662dce792eSToby Isaac   sum->interleave_basis      = interleave_basis;
10672dce792eSToby Isaac   sum->interleave_components = interleave_components;
10682dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
10692dce792eSToby Isaac }
10702dce792eSToby Isaac 
10712dce792eSToby Isaac /*@
10722dce792eSToby Isaac   PetscDualSpaceSumGetInterleave - Get whether the basis functions and components of a uniform sum are interleaved
10732dce792eSToby Isaac 
10742dce792eSToby Isaac   Logically collective
10752dce792eSToby Isaac 
10762dce792eSToby Isaac   Input Parameter:
10772dce792eSToby Isaac . sp - a `PetscDualSpace` of type `PETSCDUALSPACESUM`
10782dce792eSToby Isaac 
10792dce792eSToby Isaac   Output Parameters:
10802dce792eSToby Isaac + interleave_basis      - if `PETSC_TRUE`, the basis vectors of the subspaces are interleaved
10812dce792eSToby Isaac - interleave_components - if `PETSC_TRUE` and the space concatenates components (`PetscDualSpaceSumGetConcatenate()`),
10822dce792eSToby Isaac                           interleave the concatenated components
10832dce792eSToby Isaac 
10842dce792eSToby Isaac   Level: developer
10852dce792eSToby Isaac 
10862dce792eSToby Isaac .seealso: `PetscDualSpace`, `PETSCDUALSPACESUM`, `PETSCFEVECTOR`, `PetscDualSpaceSumSetInterleave()`
10872dce792eSToby Isaac @*/
PetscDualSpaceSumGetInterleave(PetscDualSpace sp,PeOp PetscBool * interleave_basis,PeOp PetscBool * interleave_components)1088ce78bad3SBarry Smith PetscErrorCode PetscDualSpaceSumGetInterleave(PetscDualSpace sp, PeOp PetscBool *interleave_basis, PeOp PetscBool *interleave_components)
10892dce792eSToby Isaac {
10902dce792eSToby Isaac   PetscFunctionBegin;
10912dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
10922dce792eSToby Isaac   if (interleave_basis) PetscAssertPointer(interleave_basis, 2);
10932dce792eSToby Isaac   if (interleave_components) PetscAssertPointer(interleave_components, 3);
10942dce792eSToby Isaac   PetscTryMethod(sp, "PetscDualSpaceSumGetInterleave_C", (PetscDualSpace, PetscBool *, PetscBool *), (sp, interleave_basis, interleave_components));
10952dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
10962dce792eSToby Isaac }
10972dce792eSToby Isaac 
PetscDualSpaceSumGetInterleave_Sum(PetscDualSpace sp,PetscBool * interleave_basis,PetscBool * interleave_components)10982dce792eSToby Isaac static PetscErrorCode PetscDualSpaceSumGetInterleave_Sum(PetscDualSpace sp, PetscBool *interleave_basis, PetscBool *interleave_components)
10992dce792eSToby Isaac {
11002dce792eSToby Isaac   PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data;
11014d86920dSPierre Jolivet 
11022dce792eSToby Isaac   PetscFunctionBegin;
11032dce792eSToby Isaac   if (interleave_basis) *interleave_basis = sum->interleave_basis;
11042dce792eSToby Isaac   if (interleave_components) *interleave_components = sum->interleave_components;
11052dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11062dce792eSToby Isaac }
11072dce792eSToby Isaac 
11082dce792eSToby Isaac #define PetscDualSpaceSumPassthrough(sp, func, ...) \
11092dce792eSToby Isaac   do { \
11102dce792eSToby Isaac     PetscDualSpace_Sum *sum = (PetscDualSpace_Sum *)sp->data; \
11112dce792eSToby Isaac     PetscBool           is_uniform; \
11122dce792eSToby Isaac     PetscCall(PetscDualSpaceSumIsUniform(sp, &is_uniform)); \
11132dce792eSToby Isaac     if (is_uniform && sum->numSumSpaces > 0) { \
11142dce792eSToby Isaac       PetscDualSpace subsp; \
11152dce792eSToby Isaac       PetscCall(PetscDualSpaceSumGetSubspace(sp, 0, &subsp)); \
11162dce792eSToby Isaac       PetscCall(func(subsp, __VA_ARGS__)); \
11172dce792eSToby Isaac     } \
11182dce792eSToby Isaac   } while (0)
11192dce792eSToby Isaac 
PetscDualSpaceLagrangeGetContinuity_Sum(PetscDualSpace sp,PetscBool * value)11202dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetContinuity_Sum(PetscDualSpace sp, PetscBool *value)
11212dce792eSToby Isaac {
11222dce792eSToby Isaac   PetscFunctionBegin;
11232dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetContinuity, value);
11242dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11252dce792eSToby Isaac }
11262dce792eSToby Isaac 
PetscDualSpaceLagrangeSetContinuity_Sum(PetscDualSpace sp,PetscBool value)11272dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetContinuity_Sum(PetscDualSpace sp, PetscBool value)
11282dce792eSToby Isaac {
11292dce792eSToby Isaac   PetscFunctionBegin;
11302dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetContinuity, value);
11312dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11322dce792eSToby Isaac }
11332dce792eSToby Isaac 
PetscDualSpaceLagrangeGetTensor_Sum(PetscDualSpace sp,PetscBool * value)11342dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetTensor_Sum(PetscDualSpace sp, PetscBool *value)
11352dce792eSToby Isaac {
11362dce792eSToby Isaac   PetscFunctionBegin;
11372dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetTensor, value);
11382dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11392dce792eSToby Isaac }
11402dce792eSToby Isaac 
PetscDualSpaceLagrangeSetTensor_Sum(PetscDualSpace sp,PetscBool value)11412dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetTensor_Sum(PetscDualSpace sp, PetscBool value)
11422dce792eSToby Isaac {
11432dce792eSToby Isaac   PetscFunctionBegin;
11442dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetTensor, value);
11452dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11462dce792eSToby Isaac }
11472dce792eSToby Isaac 
PetscDualSpaceLagrangeGetTrimmed_Sum(PetscDualSpace sp,PetscBool * value)11482dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetTrimmed_Sum(PetscDualSpace sp, PetscBool *value)
11492dce792eSToby Isaac {
11502dce792eSToby Isaac   PetscFunctionBegin;
11512dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetTrimmed, value);
11522dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11532dce792eSToby Isaac }
11542dce792eSToby Isaac 
PetscDualSpaceLagrangeSetTrimmed_Sum(PetscDualSpace sp,PetscBool value)11552dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetTrimmed_Sum(PetscDualSpace sp, PetscBool value)
11562dce792eSToby Isaac {
11572dce792eSToby Isaac   PetscFunctionBegin;
11582dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetTrimmed, value);
11592dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11602dce792eSToby Isaac }
11612dce792eSToby Isaac 
PetscDualSpaceLagrangeGetUseMoments_Sum(PetscDualSpace sp,PetscBool * value)11622dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetUseMoments_Sum(PetscDualSpace sp, PetscBool *value)
11632dce792eSToby Isaac {
11642dce792eSToby Isaac   PetscFunctionBegin;
11652dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetUseMoments, value);
11662dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11672dce792eSToby Isaac }
11682dce792eSToby Isaac 
PetscDualSpaceLagrangeSetUseMoments_Sum(PetscDualSpace sp,PetscBool value)11692dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetUseMoments_Sum(PetscDualSpace sp, PetscBool value)
11702dce792eSToby Isaac {
11712dce792eSToby Isaac   PetscFunctionBegin;
11722dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetUseMoments, value);
11732dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11742dce792eSToby Isaac }
11752dce792eSToby Isaac 
PetscDualSpaceLagrangeGetMomentOrder_Sum(PetscDualSpace sp,PetscInt * value)11762dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetMomentOrder_Sum(PetscDualSpace sp, PetscInt *value)
11772dce792eSToby Isaac {
11782dce792eSToby Isaac   PetscFunctionBegin;
11792dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetMomentOrder, value);
11802dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11812dce792eSToby Isaac }
11822dce792eSToby Isaac 
PetscDualSpaceLagrangeSetMomentOrder_Sum(PetscDualSpace sp,PetscInt value)11832dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetMomentOrder_Sum(PetscDualSpace sp, PetscInt value)
11842dce792eSToby Isaac {
11852dce792eSToby Isaac   PetscFunctionBegin;
11862dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetMomentOrder, value);
11872dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11882dce792eSToby Isaac }
11892dce792eSToby Isaac 
PetscDualSpaceLagrangeGetNodeType_Sum(PetscDualSpace sp,PetscDTNodeType * node_type,PetscBool * include_endpoints,PetscReal * exponent)11902dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeGetNodeType_Sum(PetscDualSpace sp, PetscDTNodeType *node_type, PetscBool *include_endpoints, PetscReal *exponent)
11912dce792eSToby Isaac {
11922dce792eSToby Isaac   PetscFunctionBegin;
11932dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeGetNodeType, node_type, include_endpoints, exponent);
11942dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
11952dce792eSToby Isaac }
11962dce792eSToby Isaac 
PetscDualSpaceLagrangeSetNodeType_Sum(PetscDualSpace sp,PetscDTNodeType node_type,PetscBool include_endpoints,PetscReal exponent)11972dce792eSToby Isaac static PetscErrorCode PetscDualSpaceLagrangeSetNodeType_Sum(PetscDualSpace sp, PetscDTNodeType node_type, PetscBool include_endpoints, PetscReal exponent)
11982dce792eSToby Isaac {
11992dce792eSToby Isaac   PetscFunctionBegin;
12002dce792eSToby Isaac   PetscDualSpaceSumPassthrough(sp, PetscDualSpaceLagrangeSetNodeType, node_type, include_endpoints, exponent);
12012dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
12022dce792eSToby Isaac }
12032dce792eSToby Isaac 
PetscDualSpaceInitialize_Sum(PetscDualSpace sp)12042dce792eSToby Isaac static PetscErrorCode PetscDualSpaceInitialize_Sum(PetscDualSpace sp)
12052dce792eSToby Isaac {
12062dce792eSToby Isaac   PetscFunctionBegin;
12072dce792eSToby Isaac   sp->ops->destroy              = PetscDualSpaceDestroy_Sum;
12082dce792eSToby Isaac   sp->ops->view                 = PetscDualSpaceView_Sum;
12092dce792eSToby Isaac   sp->ops->setfromoptions       = NULL;
12102dce792eSToby Isaac   sp->ops->duplicate            = PetscDualSpaceDuplicate_Sum;
12112dce792eSToby Isaac   sp->ops->setup                = PetscDualSpaceSetUp_Sum;
12122dce792eSToby Isaac   sp->ops->createheightsubspace = NULL;
12132dce792eSToby Isaac   sp->ops->createpointsubspace  = NULL;
12142dce792eSToby Isaac   sp->ops->getsymmetries        = PetscDualSpaceGetSymmetries_Sum;
12152dce792eSToby Isaac   sp->ops->apply                = PetscDualSpaceApplyDefault;
12162dce792eSToby Isaac   sp->ops->applyall             = PetscDualSpaceApplyAllDefault;
12172dce792eSToby Isaac   sp->ops->applyint             = PetscDualSpaceApplyInteriorDefault;
12182dce792eSToby Isaac   sp->ops->createalldata        = PetscDualSpaceCreateAllDataDefault;
12192dce792eSToby Isaac   sp->ops->createintdata        = PetscDualSpaceCreateInteriorDataDefault;
12202dce792eSToby Isaac 
12212dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetNumSubspaces_C", PetscDualSpaceSumGetNumSubspaces_Sum));
12222dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetNumSubspaces_C", PetscDualSpaceSumSetNumSubspaces_Sum));
12232dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetSubspace_C", PetscDualSpaceSumGetSubspace_Sum));
12242dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetSubspace_C", PetscDualSpaceSumSetSubspace_Sum));
12252dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetConcatenate_C", PetscDualSpaceSumGetConcatenate_Sum));
12262dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetConcatenate_C", PetscDualSpaceSumSetConcatenate_Sum));
12272dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumGetInterleave_C", PetscDualSpaceSumGetInterleave_Sum));
12282dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceSumSetInterleave_C", PetscDualSpaceSumSetInterleave_Sum));
12292dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetContinuity_C", PetscDualSpaceLagrangeGetContinuity_Sum));
12302dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetContinuity_C", PetscDualSpaceLagrangeSetContinuity_Sum));
12312dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetMomentOrder_C", PetscDualSpaceLagrangeGetMomentOrder_Sum));
12322dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetMomentOrder_C", PetscDualSpaceLagrangeSetMomentOrder_Sum));
12332dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetNodeType_C", PetscDualSpaceLagrangeGetNodeType_Sum));
12342dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetNodeType_C", PetscDualSpaceLagrangeSetNodeType_Sum));
12352dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTensor_C", PetscDualSpaceLagrangeGetTensor_Sum));
12362dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTensor_C", PetscDualSpaceLagrangeSetTensor_Sum));
12372dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetTrimmed_C", PetscDualSpaceLagrangeGetTrimmed_Sum));
12382dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetTrimmed_C", PetscDualSpaceLagrangeSetTrimmed_Sum));
12392dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeGetUseMoments_C", PetscDualSpaceLagrangeGetUseMoments_Sum));
12402dce792eSToby Isaac   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscDualSpaceLagrangeSetUseMoments_C", PetscDualSpaceLagrangeSetUseMoments_Sum));
12412dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
12422dce792eSToby Isaac }
12432dce792eSToby Isaac 
12442dce792eSToby Isaac /*MC
12452dce792eSToby Isaac   PETSCDUALSPACESUM = "sum" - A `PetscDualSpace` object that encapsulates a sum of subspaces.
12462dce792eSToby Isaac 
12472dce792eSToby Isaac   Level: intermediate
12482dce792eSToby Isaac 
12492dce792eSToby Isaac   Note:
12502dce792eSToby Isaac   That sum can either be direct or a concatenation. For example if A and B are spaces each with 2 components,
12512dce792eSToby 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
12522dce792eSToby Isaac   same reference element.
12532dce792eSToby Isaac 
12542dce792eSToby Isaac .seealso: `PetscDualSpace`, `PetscDualSpaceType`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`, `PetscDualSpaceSumGetNumSubspaces()`, `PetscDualSpaceSumSetNumSubspaces()`,
12552dce792eSToby Isaac           `PetscDualSpaceSumGetConcatenate()`, `PetscDualSpaceSumSetConcatenate()`, `PetscDualSpaceSumSetInterleave()`, `PetscDualSpaceSumGetInterleave()`
12562dce792eSToby Isaac M*/
PetscDualSpaceCreate_Sum(PetscDualSpace sp)12572dce792eSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate_Sum(PetscDualSpace sp)
12582dce792eSToby Isaac {
12592dce792eSToby Isaac   PetscDualSpace_Sum *sum;
12602dce792eSToby Isaac 
12612dce792eSToby Isaac   PetscFunctionBegin;
12622dce792eSToby Isaac   PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1);
12632dce792eSToby Isaac   PetscCall(PetscNew(&sum));
12642dce792eSToby Isaac   sum->numSumSpaces = PETSC_DEFAULT;
12652dce792eSToby Isaac   sp->data          = sum;
12662dce792eSToby Isaac   sp->k             = PETSC_FORM_DEGREE_UNDEFINED;
12672dce792eSToby Isaac   PetscCall(PetscDualSpaceInitialize_Sum(sp));
12682dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
12692dce792eSToby Isaac }
12702dce792eSToby Isaac 
12712dce792eSToby Isaac /*@
12722dce792eSToby Isaac   PetscDualSpaceCreateSum - Create a finite element dual basis that is the sum of other dual bases
12732dce792eSToby Isaac 
12742dce792eSToby Isaac   Collective
12752dce792eSToby Isaac 
12762dce792eSToby Isaac   Input Parameters:
12772dce792eSToby Isaac + numSubspaces - the number of spaces that will be added together
12782dce792eSToby Isaac . subspaces    - an array of length `numSubspaces` of spaces
12792dce792eSToby 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
12802dce792eSToby Isaac 
12812dce792eSToby Isaac   Output Parameter:
12822dce792eSToby Isaac . sumSpace - a `PetscDualSpace` of type `PETSCDUALSPACESUM`
12832dce792eSToby Isaac 
12842dce792eSToby Isaac   Level: advanced
12852dce792eSToby Isaac 
12862dce792eSToby Isaac .seealso: `PetscDualSpace`, `PETSCDUALSPACESUM`, `PETSCSPACESUM`
12872dce792eSToby Isaac @*/
PetscDualSpaceCreateSum(PetscInt numSubspaces,const PetscDualSpace subspaces[],PetscBool concatenate,PetscDualSpace * sumSpace)12882dce792eSToby Isaac PetscErrorCode PetscDualSpaceCreateSum(PetscInt numSubspaces, const PetscDualSpace subspaces[], PetscBool concatenate, PetscDualSpace *sumSpace)
12892dce792eSToby Isaac {
12902dce792eSToby Isaac   PetscInt i, Nc = 0;
12912dce792eSToby Isaac 
12922dce792eSToby Isaac   PetscFunctionBegin;
12932dce792eSToby Isaac   PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)subspaces[0]), sumSpace));
12942dce792eSToby Isaac   PetscCall(PetscDualSpaceSetType(*sumSpace, PETSCDUALSPACESUM));
12952dce792eSToby Isaac   PetscCall(PetscDualSpaceSumSetNumSubspaces(*sumSpace, numSubspaces));
12962dce792eSToby Isaac   PetscCall(PetscDualSpaceSumSetConcatenate(*sumSpace, concatenate));
12972dce792eSToby Isaac   for (i = 0; i < numSubspaces; ++i) {
12982dce792eSToby Isaac     PetscInt sNc;
12992dce792eSToby Isaac 
13002dce792eSToby Isaac     PetscCall(PetscDualSpaceSumSetSubspace(*sumSpace, i, subspaces[i]));
13012dce792eSToby Isaac     PetscCall(PetscDualSpaceGetNumComponents(subspaces[i], &sNc));
13022dce792eSToby Isaac     if (concatenate) Nc += sNc;
13032dce792eSToby Isaac     else Nc = sNc;
13042dce792eSToby Isaac 
13052dce792eSToby Isaac     if (i == 0) {
13062dce792eSToby Isaac       DM dm;
13072dce792eSToby Isaac 
13082dce792eSToby Isaac       PetscCall(PetscDualSpaceGetDM(subspaces[i], &dm));
13092dce792eSToby Isaac       PetscCall(PetscDualSpaceSetDM(*sumSpace, dm));
13102dce792eSToby Isaac     }
13112dce792eSToby Isaac   }
13122dce792eSToby Isaac   PetscCall(PetscDualSpaceSetNumComponents(*sumSpace, Nc));
13132dce792eSToby Isaac   PetscFunctionReturn(PETSC_SUCCESS);
13142dce792eSToby Isaac }
1315