xref: /petsc/src/dm/dt/space/impls/subspace/spacesubspace.c (revision b24fb147d2f783efb2f58813f80260c02fe8ea96)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac 
320cf1dd8SToby Isaac typedef struct {
420cf1dd8SToby Isaac   PetscDualSpace dualSubspace;
520cf1dd8SToby Isaac   PetscSpace     origSpace;
620cf1dd8SToby Isaac   PetscReal     *x;
720cf1dd8SToby Isaac   PetscReal     *x_alloc;
820cf1dd8SToby Isaac   PetscReal     *Jx;
920cf1dd8SToby Isaac   PetscReal     *Jx_alloc;
1020cf1dd8SToby Isaac   PetscReal     *u;
1120cf1dd8SToby Isaac   PetscReal     *u_alloc;
1220cf1dd8SToby Isaac   PetscReal     *Ju;
1320cf1dd8SToby Isaac   PetscReal     *Ju_alloc;
1420cf1dd8SToby Isaac   PetscReal     *Q;
1520cf1dd8SToby Isaac   PetscInt       Nb;
1620cf1dd8SToby Isaac } PetscSpace_Subspace;
1720cf1dd8SToby Isaac 
18d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceDestroy_Subspace(PetscSpace sp)
19d71ae5a4SJacob Faibussowitsch {
2020cf1dd8SToby Isaac   PetscSpace_Subspace *subsp;
2120cf1dd8SToby Isaac 
2220cf1dd8SToby Isaac   PetscFunctionBegin;
2320cf1dd8SToby Isaac   subsp    = (PetscSpace_Subspace *)sp->data;
2420cf1dd8SToby Isaac   subsp->x = NULL;
259566063dSJacob Faibussowitsch   PetscCall(PetscFree(subsp->x_alloc));
2620cf1dd8SToby Isaac   subsp->Jx = NULL;
279566063dSJacob Faibussowitsch   PetscCall(PetscFree(subsp->Jx_alloc));
2820cf1dd8SToby Isaac   subsp->u = NULL;
299566063dSJacob Faibussowitsch   PetscCall(PetscFree(subsp->u_alloc));
3020cf1dd8SToby Isaac   subsp->Ju = NULL;
319566063dSJacob Faibussowitsch   PetscCall(PetscFree(subsp->Ju_alloc));
329566063dSJacob Faibussowitsch   PetscCall(PetscFree(subsp->Q));
339566063dSJacob Faibussowitsch   PetscCall(PetscSpaceDestroy(&subsp->origSpace));
349566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceDestroy(&subsp->dualSubspace));
359566063dSJacob Faibussowitsch   PetscCall(PetscFree(subsp));
3620cf1dd8SToby Isaac   sp->data = NULL;
379566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePolynomialGetTensor_C", NULL));
383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3920cf1dd8SToby Isaac }
4020cf1dd8SToby Isaac 
41d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceView_Subspace(PetscSpace sp, PetscViewer viewer)
42d71ae5a4SJacob Faibussowitsch {
4320cf1dd8SToby Isaac   PetscBool            iascii;
4420cf1dd8SToby Isaac   PetscSpace_Subspace *subsp;
4520cf1dd8SToby Isaac 
4620cf1dd8SToby Isaac   PetscFunctionBegin;
4720cf1dd8SToby Isaac   subsp = (PetscSpace_Subspace *)sp->data;
489566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
4920cf1dd8SToby Isaac   if (iascii) {
5020cf1dd8SToby Isaac     PetscInt origDim, subDim, origNc, subNc, o, s;
5120cf1dd8SToby Isaac 
529566063dSJacob Faibussowitsch     PetscCall(PetscSpaceGetNumVariables(subsp->origSpace, &origDim));
539566063dSJacob Faibussowitsch     PetscCall(PetscSpaceGetNumComponents(subsp->origSpace, &origNc));
549566063dSJacob Faibussowitsch     PetscCall(PetscSpaceGetNumVariables(sp, &subDim));
559566063dSJacob Faibussowitsch     PetscCall(PetscSpaceGetNumComponents(sp, &subNc));
5620cf1dd8SToby Isaac     if (subsp->x) {
579566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Subspace-to-space domain shift:\n\n"));
5848a46eb9SPierre Jolivet       for (o = 0; o < origDim; o++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g\n", (double)subsp->x[o]));
5920cf1dd8SToby Isaac     }
6020cf1dd8SToby Isaac     if (subsp->Jx) {
619566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Subspace-to-space domain transform:\n\n"));
6220cf1dd8SToby Isaac       for (o = 0; o < origDim; o++) {
639566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)subsp->Jx[o * subDim + 0]));
649566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
6548a46eb9SPierre Jolivet         for (s = 1; s < subDim; s++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)subsp->Jx[o * subDim + s]));
669566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
679566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
6820cf1dd8SToby Isaac       }
6920cf1dd8SToby Isaac     }
7020cf1dd8SToby Isaac     if (subsp->u) {
719566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Space-to-subspace range shift:\n\n"));
7248a46eb9SPierre Jolivet       for (o = 0; o < origNc; o++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g\n", (double)subsp->u[o]));
7320cf1dd8SToby Isaac     }
7420cf1dd8SToby Isaac     if (subsp->Ju) {
759566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Space-to-subsace domain transform:\n"));
7620cf1dd8SToby Isaac       for (o = 0; o < origNc; o++) {
779566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
7848a46eb9SPierre Jolivet         for (s = 0; s < subNc; s++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)subsp->Ju[o * subNc + s]));
799566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
8020cf1dd8SToby Isaac       }
819566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
8220cf1dd8SToby Isaac     }
839566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Original space:\n"));
8420cf1dd8SToby Isaac   }
859566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
869566063dSJacob Faibussowitsch   PetscCall(PetscSpaceView(subsp->origSpace, viewer));
879566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8920cf1dd8SToby Isaac }
9020cf1dd8SToby Isaac 
91d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceEvaluate_Subspace(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
92d71ae5a4SJacob Faibussowitsch {
9320cf1dd8SToby Isaac   PetscSpace_Subspace *subsp = (PetscSpace_Subspace *)sp->data;
9420cf1dd8SToby Isaac   PetscSpace           origsp;
9520cf1dd8SToby Isaac   PetscInt             origDim, subDim, origNc, subNc, subNb, origNb, i, j, k, l, m, n, o;
9620cf1dd8SToby Isaac   PetscReal           *inpoints, *inB = NULL, *inD = NULL, *inH = NULL;
9720cf1dd8SToby Isaac 
9820cf1dd8SToby Isaac   PetscFunctionBegin;
9920cf1dd8SToby Isaac   origsp = subsp->origSpace;
1009566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumVariables(sp, &subDim));
1019566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumVariables(origsp, &origDim));
1029566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumComponents(sp, &subNc));
1039566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumComponents(origsp, &origNc));
1049566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDimension(sp, &subNb));
1059566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDimension(origsp, &origNb));
1069566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(sp->dm, npoints * origDim, MPIU_REAL, &inpoints));
10720cf1dd8SToby Isaac   for (i = 0; i < npoints; i++) {
10820cf1dd8SToby Isaac     if (subsp->x) {
10920cf1dd8SToby Isaac       for (j = 0; j < origDim; j++) inpoints[i * origDim + j] = subsp->x[j];
11020cf1dd8SToby Isaac     } else {
11120cf1dd8SToby Isaac       for (j = 0; j < origDim; j++) inpoints[i * origDim + j] = 0.0;
11220cf1dd8SToby Isaac     }
11320cf1dd8SToby Isaac     if (subsp->Jx) {
11420cf1dd8SToby Isaac       for (j = 0; j < origDim; j++) {
115ad540459SPierre Jolivet         for (k = 0; k < subDim; k++) inpoints[i * origDim + j] += subsp->Jx[j * subDim + k] * points[i * subDim + k];
11620cf1dd8SToby Isaac       }
11720cf1dd8SToby Isaac     } else {
118ad540459SPierre Jolivet       for (j = 0; j < PetscMin(subDim, origDim); j++) inpoints[i * origDim + j] += points[i * subDim + j];
11920cf1dd8SToby Isaac     }
12020cf1dd8SToby Isaac   }
12148a46eb9SPierre Jolivet   if (B) PetscCall(DMGetWorkArray(sp->dm, npoints * origNb * origNc, MPIU_REAL, &inB));
12248a46eb9SPierre Jolivet   if (D) PetscCall(DMGetWorkArray(sp->dm, npoints * origNb * origNc * origDim, MPIU_REAL, &inD));
12348a46eb9SPierre Jolivet   if (H) PetscCall(DMGetWorkArray(sp->dm, npoints * origNb * origNc * origDim * origDim, MPIU_REAL, &inH));
1249566063dSJacob Faibussowitsch   PetscCall(PetscSpaceEvaluate(origsp, npoints, inpoints, inB, inD, inH));
12520cf1dd8SToby Isaac   if (H) {
12620cf1dd8SToby Isaac     PetscReal *phi, *psi;
12720cf1dd8SToby Isaac 
1289566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(sp->dm, origNc * origDim * origDim, MPIU_REAL, &phi));
1299566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(sp->dm, origNc * subDim * subDim, MPIU_REAL, &psi));
13020cf1dd8SToby Isaac     for (i = 0; i < npoints * subNb * subNc * subDim; i++) D[i] = 0.0;
13120cf1dd8SToby Isaac     for (i = 0; i < subNb; i++) {
13220cf1dd8SToby Isaac       const PetscReal *subq = &subsp->Q[i * origNb];
13320cf1dd8SToby Isaac 
13420cf1dd8SToby Isaac       for (j = 0; j < npoints; j++) {
13520cf1dd8SToby Isaac         for (k = 0; k < origNc * origDim; k++) phi[k] = 0.;
13620cf1dd8SToby Isaac         for (k = 0; k < origNc * subDim; k++) psi[k] = 0.;
13720cf1dd8SToby Isaac         for (k = 0; k < origNb; k++) {
138ad540459SPierre Jolivet           for (l = 0; l < origNc * origDim * origDim; l++) phi[l] += inH[(j * origNb + k) * origNc * origDim * origDim + l] * subq[k];
13920cf1dd8SToby Isaac         }
14020cf1dd8SToby Isaac         if (subsp->Jx) {
14120cf1dd8SToby Isaac           for (k = 0; k < subNc; k++) {
14220cf1dd8SToby Isaac             for (l = 0; l < subDim; l++) {
14320cf1dd8SToby Isaac               for (m = 0; m < origDim; m++) {
14420cf1dd8SToby Isaac                 for (n = 0; n < subDim; n++) {
145ad540459SPierre Jolivet                   for (o = 0; o < origDim; o++) psi[(k * subDim + l) * subDim + n] += subsp->Jx[m * subDim + l] * subsp->Jx[o * subDim + n] * phi[(k * origDim + m) * origDim + o];
14620cf1dd8SToby Isaac                 }
14720cf1dd8SToby Isaac               }
14820cf1dd8SToby Isaac             }
14920cf1dd8SToby Isaac           }
15020cf1dd8SToby Isaac         } else {
15120cf1dd8SToby Isaac           for (k = 0; k < subNc; k++) {
15220cf1dd8SToby Isaac             for (l = 0; l < PetscMin(subDim, origDim); l++) {
153ad540459SPierre Jolivet               for (m = 0; m < PetscMin(subDim, origDim); m++) psi[(k * subDim + l) * subDim + m] += phi[(k * origDim + l) * origDim + m];
15420cf1dd8SToby Isaac             }
15520cf1dd8SToby Isaac           }
15620cf1dd8SToby Isaac         }
15720cf1dd8SToby Isaac         if (subsp->Ju) {
15820cf1dd8SToby Isaac           for (k = 0; k < subNc; k++) {
15920cf1dd8SToby Isaac             for (l = 0; l < origNc; l++) {
160ad540459SPierre Jolivet               for (m = 0; m < subDim * subDim; m++) H[((j * subNb + i) * subNc + k) * subDim * subDim + m] += subsp->Ju[k * origNc + l] * psi[l * subDim * subDim + m];
16120cf1dd8SToby Isaac             }
16220cf1dd8SToby Isaac           }
1639371c9d4SSatish Balay         } else {
16420cf1dd8SToby Isaac           for (k = 0; k < PetscMin(subNc, origNc); k++) {
165ad540459SPierre Jolivet             for (l = 0; l < subDim * subDim; l++) H[((j * subNb + i) * subNc + k) * subDim * subDim + l] += psi[k * subDim * subDim + l];
16620cf1dd8SToby Isaac           }
16720cf1dd8SToby Isaac         }
16820cf1dd8SToby Isaac       }
16920cf1dd8SToby Isaac     }
1709566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(sp->dm, subNc * origDim, MPIU_REAL, &psi));
1719566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(sp->dm, origNc * origDim, MPIU_REAL, &phi));
1729566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(sp->dm, npoints * origNb * origNc * origDim, MPIU_REAL, &inH));
17320cf1dd8SToby Isaac   }
17420cf1dd8SToby Isaac   if (D) {
17520cf1dd8SToby Isaac     PetscReal *phi, *psi;
17620cf1dd8SToby Isaac 
1779566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(sp->dm, origNc * origDim, MPIU_REAL, &phi));
1789566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(sp->dm, origNc * subDim, MPIU_REAL, &psi));
17920cf1dd8SToby Isaac     for (i = 0; i < npoints * subNb * subNc * subDim; i++) D[i] = 0.0;
18020cf1dd8SToby Isaac     for (i = 0; i < subNb; i++) {
18120cf1dd8SToby Isaac       const PetscReal *subq = &subsp->Q[i * origNb];
18220cf1dd8SToby Isaac 
18320cf1dd8SToby Isaac       for (j = 0; j < npoints; j++) {
18420cf1dd8SToby Isaac         for (k = 0; k < origNc * origDim; k++) phi[k] = 0.;
18520cf1dd8SToby Isaac         for (k = 0; k < origNc * subDim; k++) psi[k] = 0.;
18620cf1dd8SToby Isaac         for (k = 0; k < origNb; k++) {
187ad540459SPierre Jolivet           for (l = 0; l < origNc * origDim; l++) phi[l] += inD[(j * origNb + k) * origNc * origDim + l] * subq[k];
18820cf1dd8SToby Isaac         }
18920cf1dd8SToby Isaac         if (subsp->Jx) {
19020cf1dd8SToby Isaac           for (k = 0; k < subNc; k++) {
19120cf1dd8SToby Isaac             for (l = 0; l < subDim; l++) {
192ad540459SPierre Jolivet               for (m = 0; m < origDim; m++) psi[k * subDim + l] += subsp->Jx[m * subDim + l] * phi[k * origDim + m];
19320cf1dd8SToby Isaac             }
19420cf1dd8SToby Isaac           }
19520cf1dd8SToby Isaac         } else {
19620cf1dd8SToby Isaac           for (k = 0; k < subNc; k++) {
197ad540459SPierre Jolivet             for (l = 0; l < PetscMin(subDim, origDim); l++) psi[k * subDim + l] += phi[k * origDim + l];
19820cf1dd8SToby Isaac           }
19920cf1dd8SToby Isaac         }
20020cf1dd8SToby Isaac         if (subsp->Ju) {
20120cf1dd8SToby Isaac           for (k = 0; k < subNc; k++) {
20220cf1dd8SToby Isaac             for (l = 0; l < origNc; l++) {
203ad540459SPierre Jolivet               for (m = 0; m < subDim; m++) D[((j * subNb + i) * subNc + k) * subDim + m] += subsp->Ju[k * origNc + l] * psi[l * subDim + m];
20420cf1dd8SToby Isaac             }
20520cf1dd8SToby Isaac           }
2069371c9d4SSatish Balay         } else {
20720cf1dd8SToby Isaac           for (k = 0; k < PetscMin(subNc, origNc); k++) {
208ad540459SPierre Jolivet             for (l = 0; l < subDim; l++) D[((j * subNb + i) * subNc + k) * subDim + l] += psi[k * subDim + l];
20920cf1dd8SToby Isaac           }
21020cf1dd8SToby Isaac         }
21120cf1dd8SToby Isaac       }
21220cf1dd8SToby Isaac     }
2139566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(sp->dm, subNc * origDim, MPIU_REAL, &psi));
2149566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(sp->dm, origNc * origDim, MPIU_REAL, &phi));
2159566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(sp->dm, npoints * origNb * origNc * origDim, MPIU_REAL, &inD));
21620cf1dd8SToby Isaac   }
21720cf1dd8SToby Isaac   if (B) {
21820cf1dd8SToby Isaac     PetscReal *phi;
21920cf1dd8SToby Isaac 
2209566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(sp->dm, origNc, MPIU_REAL, &phi));
22120cf1dd8SToby Isaac     if (subsp->u) {
22220cf1dd8SToby Isaac       for (i = 0; i < npoints * subNb; i++) {
22320cf1dd8SToby Isaac         for (j = 0; j < subNc; j++) B[i * subNc + j] = subsp->u[j];
22420cf1dd8SToby Isaac       }
22520cf1dd8SToby Isaac     } else {
22620cf1dd8SToby Isaac       for (i = 0; i < npoints * subNb * subNc; i++) B[i] = 0.0;
22720cf1dd8SToby Isaac     }
22820cf1dd8SToby Isaac     for (i = 0; i < subNb; i++) {
22920cf1dd8SToby Isaac       const PetscReal *subq = &subsp->Q[i * origNb];
23020cf1dd8SToby Isaac 
23120cf1dd8SToby Isaac       for (j = 0; j < npoints; j++) {
23220cf1dd8SToby Isaac         for (k = 0; k < origNc; k++) phi[k] = 0.;
23320cf1dd8SToby Isaac         for (k = 0; k < origNb; k++) {
234ad540459SPierre Jolivet           for (l = 0; l < origNc; l++) phi[l] += inB[(j * origNb + k) * origNc + l] * subq[k];
23520cf1dd8SToby Isaac         }
23620cf1dd8SToby Isaac         if (subsp->Ju) {
23720cf1dd8SToby Isaac           for (k = 0; k < subNc; k++) {
238ad540459SPierre Jolivet             for (l = 0; l < origNc; l++) B[(j * subNb + i) * subNc + k] += subsp->Ju[k * origNc + l] * phi[l];
23920cf1dd8SToby Isaac           }
2409371c9d4SSatish Balay         } else {
241ad540459SPierre Jolivet           for (k = 0; k < PetscMin(subNc, origNc); k++) B[(j * subNb + i) * subNc + k] += phi[k];
24220cf1dd8SToby Isaac         }
24320cf1dd8SToby Isaac       }
24420cf1dd8SToby Isaac     }
2459566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(sp->dm, origNc, MPIU_REAL, &phi));
2469566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(sp->dm, npoints * origNb * origNc, MPIU_REAL, &inB));
24720cf1dd8SToby Isaac   }
2489566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(sp->dm, npoints * origDim, MPIU_REAL, &inpoints));
2493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
25020cf1dd8SToby Isaac }
25120cf1dd8SToby Isaac 
252d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Subspace(PetscSpace sp)
253d71ae5a4SJacob Faibussowitsch {
25420cf1dd8SToby Isaac   PetscSpace_Subspace *subsp;
255362febeeSStefano Zampini 
256362febeeSStefano Zampini   PetscFunctionBegin;
2574dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&subsp));
25820cf1dd8SToby Isaac   sp->data = (void *)subsp;
2593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
26020cf1dd8SToby Isaac }
26120cf1dd8SToby Isaac 
262d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceGetDimension_Subspace(PetscSpace sp, PetscInt *dim)
263d71ae5a4SJacob Faibussowitsch {
26420cf1dd8SToby Isaac   PetscSpace_Subspace *subsp;
26520cf1dd8SToby Isaac 
26620cf1dd8SToby Isaac   PetscFunctionBegin;
26720cf1dd8SToby Isaac   subsp = (PetscSpace_Subspace *)sp->data;
26820cf1dd8SToby Isaac   *dim  = subsp->Nb;
2693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27020cf1dd8SToby Isaac }
27120cf1dd8SToby Isaac 
272d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceSetUp_Subspace(PetscSpace sp)
273d71ae5a4SJacob Faibussowitsch {
27420cf1dd8SToby Isaac   const PetscReal     *x;
27520cf1dd8SToby Isaac   const PetscReal     *Jx;
27620cf1dd8SToby Isaac   const PetscReal     *u;
27720cf1dd8SToby Isaac   const PetscReal     *Ju;
27820cf1dd8SToby Isaac   PetscDualSpace       dualSubspace;
27920cf1dd8SToby Isaac   PetscSpace           origSpace;
28020cf1dd8SToby Isaac   PetscInt             origDim, subDim, origNc, subNc, origNb, subNb, f, i, j, numPoints, offset;
28120cf1dd8SToby Isaac   PetscReal           *allPoints, *allWeights, *B, *V;
28220cf1dd8SToby Isaac   DM                   dm;
28320cf1dd8SToby Isaac   PetscSpace_Subspace *subsp;
28420cf1dd8SToby Isaac 
28520cf1dd8SToby Isaac   PetscFunctionBegin;
28620cf1dd8SToby Isaac   subsp        = (PetscSpace_Subspace *)sp->data;
28720cf1dd8SToby Isaac   x            = subsp->x;
28820cf1dd8SToby Isaac   Jx           = subsp->Jx;
28920cf1dd8SToby Isaac   u            = subsp->u;
29020cf1dd8SToby Isaac   Ju           = subsp->Ju;
29120cf1dd8SToby Isaac   origSpace    = subsp->origSpace;
29220cf1dd8SToby Isaac   dualSubspace = subsp->dualSubspace;
2939566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumComponents(origSpace, &origNc));
2949566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumVariables(origSpace, &origDim));
2959566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(dualSubspace, &dm));
2969566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &subDim));
2979566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDimension(origSpace, &origNb));
2989566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(dualSubspace, &subNb));
2999566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetNumComponents(dualSubspace, &subNc));
30020cf1dd8SToby Isaac 
30120cf1dd8SToby Isaac   for (f = 0, numPoints = 0; f < subNb; f++) {
30220cf1dd8SToby Isaac     PetscQuadrature q;
30320cf1dd8SToby Isaac     PetscInt        qNp;
30420cf1dd8SToby Isaac 
3059566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetFunctional(dualSubspace, f, &q));
3069566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(q, NULL, NULL, &qNp, NULL, NULL));
30720cf1dd8SToby Isaac     numPoints += qNp;
30820cf1dd8SToby Isaac   }
3099566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(subNb * origNb, &V));
3109566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(numPoints * origDim, &allPoints, numPoints * origNc, &allWeights, numPoints * origNb * origNc, &B));
31120cf1dd8SToby Isaac   for (f = 0, offset = 0; f < subNb; f++) {
31220cf1dd8SToby Isaac     PetscQuadrature  q;
31320cf1dd8SToby Isaac     PetscInt         qNp, p;
31420cf1dd8SToby Isaac     const PetscReal *qp;
31520cf1dd8SToby Isaac     const PetscReal *qw;
31620cf1dd8SToby Isaac 
3179566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetFunctional(dualSubspace, f, &q));
3189566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(q, NULL, NULL, &qNp, &qp, &qw));
31920cf1dd8SToby Isaac     for (p = 0; p < qNp; p++, offset++) {
32020cf1dd8SToby Isaac       if (x) {
32120cf1dd8SToby Isaac         for (i = 0; i < origDim; i++) allPoints[origDim * offset + i] = x[i];
32220cf1dd8SToby Isaac       } else {
32320cf1dd8SToby Isaac         for (i = 0; i < origDim; i++) allPoints[origDim * offset + i] = 0.0;
32420cf1dd8SToby Isaac       }
32520cf1dd8SToby Isaac       if (Jx) {
32620cf1dd8SToby Isaac         for (i = 0; i < origDim; i++) {
327ad540459SPierre Jolivet           for (j = 0; j < subDim; j++) allPoints[origDim * offset + i] += Jx[i * subDim + j] * qp[j];
32820cf1dd8SToby Isaac         }
32920cf1dd8SToby Isaac       } else {
33020cf1dd8SToby Isaac         for (i = 0; i < PetscMin(subDim, origDim); i++) allPoints[origDim * offset + i] += qp[i];
33120cf1dd8SToby Isaac       }
33220cf1dd8SToby Isaac       for (i = 0; i < origNc; i++) allWeights[origNc * offset + i] = 0.0;
33320cf1dd8SToby Isaac       if (Ju) {
33420cf1dd8SToby Isaac         for (i = 0; i < origNc; i++) {
335ad540459SPierre Jolivet           for (j = 0; j < subNc; j++) allWeights[offset * origNc + i] += qw[j] * Ju[j * origNc + i];
33620cf1dd8SToby Isaac         }
33720cf1dd8SToby Isaac       } else {
33820cf1dd8SToby Isaac         for (i = 0; i < PetscMin(subNc, origNc); i++) allWeights[offset * origNc + i] += qw[i];
33920cf1dd8SToby Isaac       }
34020cf1dd8SToby Isaac     }
34120cf1dd8SToby Isaac   }
3429566063dSJacob Faibussowitsch   PetscCall(PetscSpaceEvaluate(origSpace, numPoints, allPoints, B, NULL, NULL));
34320cf1dd8SToby Isaac   for (f = 0, offset = 0; f < subNb; f++) {
34420cf1dd8SToby Isaac     PetscInt         b, p, s, qNp;
34520cf1dd8SToby Isaac     PetscQuadrature  q;
34620cf1dd8SToby Isaac     const PetscReal *qw;
34720cf1dd8SToby Isaac 
3489566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetFunctional(dualSubspace, f, &q));
3499566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(q, NULL, NULL, &qNp, NULL, &qw));
35020cf1dd8SToby Isaac     if (u) {
35120cf1dd8SToby Isaac       for (b = 0; b < origNb; b++) {
352ad540459SPierre Jolivet         for (s = 0; s < subNc; s++) V[f * origNb + b] += qw[s] * u[s];
35320cf1dd8SToby Isaac       }
35420cf1dd8SToby Isaac     } else {
35520cf1dd8SToby Isaac       for (b = 0; b < origNb; b++) V[f * origNb + b] = 0.0;
35620cf1dd8SToby Isaac     }
35720cf1dd8SToby Isaac     for (p = 0; p < qNp; p++, offset++) {
35820cf1dd8SToby Isaac       for (b = 0; b < origNb; b++) {
359ad540459SPierre Jolivet         for (s = 0; s < origNc; s++) V[f * origNb + b] += B[(offset * origNb + b) * origNc + s] * allWeights[offset * origNc + s];
36020cf1dd8SToby Isaac       }
36120cf1dd8SToby Isaac     }
36220cf1dd8SToby Isaac   }
36320cf1dd8SToby Isaac   /* orthnormalize rows of V */
36420cf1dd8SToby Isaac   for (f = 0; f < subNb; f++) {
36520cf1dd8SToby Isaac     PetscReal rho = 0.0, scal;
36620cf1dd8SToby Isaac 
36720cf1dd8SToby Isaac     for (i = 0; i < origNb; i++) rho += PetscSqr(V[f * origNb + i]);
36820cf1dd8SToby Isaac 
36920cf1dd8SToby Isaac     scal = 1. / PetscSqrtReal(rho);
37020cf1dd8SToby Isaac 
37120cf1dd8SToby Isaac     for (i = 0; i < origNb; i++) V[f * origNb + i] *= scal;
37220cf1dd8SToby Isaac     for (j = f + 1; j < subNb; j++) {
37320cf1dd8SToby Isaac       for (i = 0, scal = 0.; i < origNb; i++) scal += V[f * origNb + i] * V[j * origNb + i];
37420cf1dd8SToby Isaac       for (i = 0; i < origNb; i++) V[j * origNb + i] -= V[f * origNb + i] * scal;
37520cf1dd8SToby Isaac     }
37620cf1dd8SToby Isaac   }
3779566063dSJacob Faibussowitsch   PetscCall(PetscFree3(allPoints, allWeights, B));
37820cf1dd8SToby Isaac   subsp->Q = V;
3793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38020cf1dd8SToby Isaac }
38120cf1dd8SToby Isaac 
382d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpacePolynomialGetTensor_Subspace(PetscSpace sp, PetscBool *poly)
383d71ae5a4SJacob Faibussowitsch {
38420cf1dd8SToby Isaac   PetscSpace_Subspace *subsp = (PetscSpace_Subspace *)sp->data;
38520cf1dd8SToby Isaac 
38620cf1dd8SToby Isaac   PetscFunctionBegin;
38720cf1dd8SToby Isaac   *poly = PETSC_FALSE;
3889566063dSJacob Faibussowitsch   PetscCall(PetscSpacePolynomialGetTensor(subsp->origSpace, poly));
38920cf1dd8SToby Isaac   if (*poly) {
39020cf1dd8SToby Isaac     if (subsp->Jx) {
39120cf1dd8SToby Isaac       PetscInt subDim, origDim, i, j;
39220cf1dd8SToby Isaac       PetscInt maxnnz;
39320cf1dd8SToby Isaac 
3949566063dSJacob Faibussowitsch       PetscCall(PetscSpaceGetNumVariables(subsp->origSpace, &origDim));
3959566063dSJacob Faibussowitsch       PetscCall(PetscSpaceGetNumVariables(sp, &subDim));
39620cf1dd8SToby Isaac       maxnnz = 0;
39720cf1dd8SToby Isaac       for (i = 0; i < origDim; i++) {
39820cf1dd8SToby Isaac         PetscInt nnz = 0;
39920cf1dd8SToby Isaac 
40020cf1dd8SToby Isaac         for (j = 0; j < subDim; j++) nnz += (subsp->Jx[i * subDim + j] != 0.);
40120cf1dd8SToby Isaac         maxnnz = PetscMax(maxnnz, nnz);
40220cf1dd8SToby Isaac       }
40320cf1dd8SToby Isaac       for (j = 0; j < subDim; j++) {
40420cf1dd8SToby Isaac         PetscInt nnz = 0;
40520cf1dd8SToby Isaac 
40620cf1dd8SToby Isaac         for (i = 0; i < origDim; i++) nnz += (subsp->Jx[i * subDim + j] != 0.);
40720cf1dd8SToby Isaac         maxnnz = PetscMax(maxnnz, nnz);
40820cf1dd8SToby Isaac       }
40920cf1dd8SToby Isaac       if (maxnnz > 1) *poly = PETSC_FALSE;
41020cf1dd8SToby Isaac     }
41120cf1dd8SToby Isaac   }
4123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
41320cf1dd8SToby Isaac }
41420cf1dd8SToby Isaac 
415d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceInitialize_Subspace(PetscSpace sp)
416d71ae5a4SJacob Faibussowitsch {
41720cf1dd8SToby Isaac   PetscFunctionBegin;
41820cf1dd8SToby Isaac   sp->ops->setup        = PetscSpaceSetUp_Subspace;
41920cf1dd8SToby Isaac   sp->ops->view         = PetscSpaceView_Subspace;
42020cf1dd8SToby Isaac   sp->ops->destroy      = PetscSpaceDestroy_Subspace;
42120cf1dd8SToby Isaac   sp->ops->getdimension = PetscSpaceGetDimension_Subspace;
42220cf1dd8SToby Isaac   sp->ops->evaluate     = PetscSpaceEvaluate_Subspace;
4239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Subspace));
4243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
42520cf1dd8SToby Isaac }
426*b24fb147SBarry Smith /*@
427*b24fb147SBarry Smith     PetscSpaceCreateSubspace - creates a subspace from a an `origSpace` and its dual `dualSubspace`
42820cf1dd8SToby Isaac 
429*b24fb147SBarry Smith     Input Parameters:
430*b24fb147SBarry Smith +   origSpace - the original `PetscSpace`
431*b24fb147SBarry Smith .   dualSubspace - no idea
432*b24fb147SBarry Smith .   x - no idea
433*b24fb147SBarry Smith .   Jx - no idea
434*b24fb147SBarry Smith .   u - no idea
435*b24fb147SBarry Smith .   Ju - no idea
436*b24fb147SBarry Smith -   copymode - whether to copy, borrow, or own some of the input arrays I guess
437*b24fb147SBarry Smith 
438*b24fb147SBarry Smith     Output Parameter:
439*b24fb147SBarry Smith .   subspace - the subspace
440*b24fb147SBarry Smith 
441*b24fb147SBarry Smith     Level: advanced
442*b24fb147SBarry Smith 
443*b24fb147SBarry Smith .seealso: `PetscSpace`, `PetscDualSpace`, `PetscCopyMode`, `PetscSpaceType`
444*b24fb147SBarry Smith @*/
445d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSpaceCreateSubspace(PetscSpace origSpace, PetscDualSpace dualSubspace, PetscReal *x, PetscReal *Jx, PetscReal *u, PetscReal *Ju, PetscCopyMode copymode, PetscSpace *subspace)
446d71ae5a4SJacob Faibussowitsch {
44720cf1dd8SToby Isaac   PetscSpace_Subspace *subsp;
44820cf1dd8SToby Isaac   PetscInt             origDim, subDim, origNc, subNc, subNb;
44920cf1dd8SToby Isaac   PetscInt             order;
45020cf1dd8SToby Isaac   DM                   dm;
45120cf1dd8SToby Isaac 
45220cf1dd8SToby Isaac   PetscFunctionBegin;
45320cf1dd8SToby Isaac   PetscValidHeaderSpecific(origSpace, PETSCSPACE_CLASSID, 1);
45420cf1dd8SToby Isaac   PetscValidHeaderSpecific(dualSubspace, PETSCDUALSPACE_CLASSID, 2);
45520cf1dd8SToby Isaac   if (x) PetscValidRealPointer(x, 3);
45620cf1dd8SToby Isaac   if (Jx) PetscValidRealPointer(Jx, 4);
45720cf1dd8SToby Isaac   if (u) PetscValidRealPointer(u, 5);
45820cf1dd8SToby Isaac   if (Ju) PetscValidRealPointer(Ju, 6);
459064a246eSJacob Faibussowitsch   PetscValidPointer(subspace, 8);
4609566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumComponents(origSpace, &origNc));
4619566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumVariables(origSpace, &origDim));
4629566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(dualSubspace, &dm));
4639566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &subDim));
4649566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(dualSubspace, &subNb));
4659566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetNumComponents(dualSubspace, &subNc));
4669566063dSJacob Faibussowitsch   PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)origSpace), subspace));
4679566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetType(*subspace, PETSCSPACESUBSPACE));
4689566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetNumVariables(*subspace, subDim));
4699566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetNumComponents(*subspace, subNc));
4709566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDegree(origSpace, &order, NULL));
4719566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetDegree(*subspace, order, PETSC_DETERMINE));
47220cf1dd8SToby Isaac   subsp     = (PetscSpace_Subspace *)(*subspace)->data;
47320cf1dd8SToby Isaac   subsp->Nb = subNb;
47420cf1dd8SToby Isaac   switch (copymode) {
47520cf1dd8SToby Isaac   case PETSC_OWN_POINTER:
47620cf1dd8SToby Isaac     if (x) subsp->x_alloc = x;
47720cf1dd8SToby Isaac     if (Jx) subsp->Jx_alloc = Jx;
47820cf1dd8SToby Isaac     if (u) subsp->u_alloc = u;
47920cf1dd8SToby Isaac     if (Ju) subsp->Ju_alloc = Ju;
480f4d061e9SPierre Jolivet     /* fall through */
48120cf1dd8SToby Isaac   case PETSC_USE_POINTER:
48220cf1dd8SToby Isaac     if (x) subsp->x = x;
48320cf1dd8SToby Isaac     if (Jx) subsp->Jx = Jx;
48420cf1dd8SToby Isaac     if (u) subsp->u = u;
48520cf1dd8SToby Isaac     if (Ju) subsp->Ju = Ju;
48620cf1dd8SToby Isaac     break;
48720cf1dd8SToby Isaac   case PETSC_COPY_VALUES:
48820cf1dd8SToby Isaac     if (x) {
4899566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(origDim, &subsp->x_alloc));
4909566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(subsp->x_alloc, x, origDim));
49120cf1dd8SToby Isaac       subsp->x = subsp->x_alloc;
49220cf1dd8SToby Isaac     }
49320cf1dd8SToby Isaac     if (Jx) {
4949566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(origDim * subDim, &subsp->Jx_alloc));
4959566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(subsp->Jx_alloc, Jx, origDim * subDim));
49620cf1dd8SToby Isaac       subsp->Jx = subsp->Jx_alloc;
49720cf1dd8SToby Isaac     }
49820cf1dd8SToby Isaac     if (u) {
4999566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(subNc, &subsp->u_alloc));
5009566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(subsp->u_alloc, u, subNc));
50120cf1dd8SToby Isaac       subsp->u = subsp->u_alloc;
50220cf1dd8SToby Isaac     }
50320cf1dd8SToby Isaac     if (Ju) {
5049566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(origNc * subNc, &subsp->Ju_alloc));
5059566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(subsp->Ju_alloc, Ju, origNc * subNc));
50620cf1dd8SToby Isaac       subsp->Ju = subsp->Ju_alloc;
50720cf1dd8SToby Isaac     }
50820cf1dd8SToby Isaac     break;
509d71ae5a4SJacob Faibussowitsch   default:
510d71ae5a4SJacob Faibussowitsch     SETERRQ(PetscObjectComm((PetscObject)origSpace), PETSC_ERR_ARG_OUTOFRANGE, "Unknown copy mode");
51120cf1dd8SToby Isaac   }
5129566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)origSpace));
51320cf1dd8SToby Isaac   subsp->origSpace = origSpace;
5149566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dualSubspace));
51520cf1dd8SToby Isaac   subsp->dualSubspace = dualSubspace;
5169566063dSJacob Faibussowitsch   PetscCall(PetscSpaceInitialize_Subspace(*subspace));
5173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
51820cf1dd8SToby Isaac }
519