xref: /petsc/src/dm/dt/space/impls/subspace/spacesubspace.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
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 
189371c9d4SSatish Balay static PetscErrorCode PetscSpaceDestroy_Subspace(PetscSpace sp) {
1920cf1dd8SToby Isaac   PetscSpace_Subspace *subsp;
2020cf1dd8SToby Isaac 
2120cf1dd8SToby Isaac   PetscFunctionBegin;
2220cf1dd8SToby Isaac   subsp    = (PetscSpace_Subspace *)sp->data;
2320cf1dd8SToby Isaac   subsp->x = NULL;
249566063dSJacob Faibussowitsch   PetscCall(PetscFree(subsp->x_alloc));
2520cf1dd8SToby Isaac   subsp->Jx = NULL;
269566063dSJacob Faibussowitsch   PetscCall(PetscFree(subsp->Jx_alloc));
2720cf1dd8SToby Isaac   subsp->u = NULL;
289566063dSJacob Faibussowitsch   PetscCall(PetscFree(subsp->u_alloc));
2920cf1dd8SToby Isaac   subsp->Ju = NULL;
309566063dSJacob Faibussowitsch   PetscCall(PetscFree(subsp->Ju_alloc));
319566063dSJacob Faibussowitsch   PetscCall(PetscFree(subsp->Q));
329566063dSJacob Faibussowitsch   PetscCall(PetscSpaceDestroy(&subsp->origSpace));
339566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceDestroy(&subsp->dualSubspace));
349566063dSJacob Faibussowitsch   PetscCall(PetscFree(subsp));
3520cf1dd8SToby Isaac   sp->data = NULL;
369566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePolynomialGetTensor_C", NULL));
3720cf1dd8SToby Isaac   PetscFunctionReturn(0);
3820cf1dd8SToby Isaac }
3920cf1dd8SToby Isaac 
409371c9d4SSatish Balay static PetscErrorCode PetscSpaceView_Subspace(PetscSpace sp, PetscViewer viewer) {
4120cf1dd8SToby Isaac   PetscBool            iascii;
4220cf1dd8SToby Isaac   PetscSpace_Subspace *subsp;
4320cf1dd8SToby Isaac 
4420cf1dd8SToby Isaac   PetscFunctionBegin;
4520cf1dd8SToby Isaac   subsp = (PetscSpace_Subspace *)sp->data;
469566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
4720cf1dd8SToby Isaac   if (iascii) {
4820cf1dd8SToby Isaac     PetscInt origDim, subDim, origNc, subNc, o, s;
4920cf1dd8SToby Isaac 
509566063dSJacob Faibussowitsch     PetscCall(PetscSpaceGetNumVariables(subsp->origSpace, &origDim));
519566063dSJacob Faibussowitsch     PetscCall(PetscSpaceGetNumComponents(subsp->origSpace, &origNc));
529566063dSJacob Faibussowitsch     PetscCall(PetscSpaceGetNumVariables(sp, &subDim));
539566063dSJacob Faibussowitsch     PetscCall(PetscSpaceGetNumComponents(sp, &subNc));
5420cf1dd8SToby Isaac     if (subsp->x) {
559566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Subspace-to-space domain shift:\n\n"));
5648a46eb9SPierre Jolivet       for (o = 0; o < origDim; o++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g\n", (double)subsp->x[o]));
5720cf1dd8SToby Isaac     }
5820cf1dd8SToby Isaac     if (subsp->Jx) {
599566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Subspace-to-space domain transform:\n\n"));
6020cf1dd8SToby Isaac       for (o = 0; o < origDim; o++) {
619566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)subsp->Jx[o * subDim + 0]));
629566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
6348a46eb9SPierre Jolivet         for (s = 1; s < subDim; s++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)subsp->Jx[o * subDim + s]));
649566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
659566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
6620cf1dd8SToby Isaac       }
6720cf1dd8SToby Isaac     }
6820cf1dd8SToby Isaac     if (subsp->u) {
699566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Space-to-subspace range shift:\n\n"));
7048a46eb9SPierre Jolivet       for (o = 0; o < origNc; o++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g\n", (double)subsp->u[o]));
7120cf1dd8SToby Isaac     }
7220cf1dd8SToby Isaac     if (subsp->Ju) {
739566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "Space-to-subsace domain transform:\n"));
7420cf1dd8SToby Isaac       for (o = 0; o < origNc; o++) {
759566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE));
7648a46eb9SPierre Jolivet         for (s = 0; s < subNc; s++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)subsp->Ju[o * subNc + s]));
779566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE));
7820cf1dd8SToby Isaac       }
799566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "\n"));
8020cf1dd8SToby Isaac     }
819566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Original space:\n"));
8220cf1dd8SToby Isaac   }
839566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
849566063dSJacob Faibussowitsch   PetscCall(PetscSpaceView(subsp->origSpace, viewer));
859566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
8620cf1dd8SToby Isaac   PetscFunctionReturn(0);
8720cf1dd8SToby Isaac }
8820cf1dd8SToby Isaac 
899371c9d4SSatish Balay static PetscErrorCode PetscSpaceEvaluate_Subspace(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) {
9020cf1dd8SToby Isaac   PetscSpace_Subspace *subsp = (PetscSpace_Subspace *)sp->data;
9120cf1dd8SToby Isaac   PetscSpace           origsp;
9220cf1dd8SToby Isaac   PetscInt             origDim, subDim, origNc, subNc, subNb, origNb, i, j, k, l, m, n, o;
9320cf1dd8SToby Isaac   PetscReal           *inpoints, *inB = NULL, *inD = NULL, *inH = NULL;
9420cf1dd8SToby Isaac 
9520cf1dd8SToby Isaac   PetscFunctionBegin;
9620cf1dd8SToby Isaac   origsp = subsp->origSpace;
979566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumVariables(sp, &subDim));
989566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumVariables(origsp, &origDim));
999566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumComponents(sp, &subNc));
1009566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumComponents(origsp, &origNc));
1019566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDimension(sp, &subNb));
1029566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDimension(origsp, &origNb));
1039566063dSJacob Faibussowitsch   PetscCall(DMGetWorkArray(sp->dm, npoints * origDim, MPIU_REAL, &inpoints));
10420cf1dd8SToby Isaac   for (i = 0; i < npoints; i++) {
10520cf1dd8SToby Isaac     if (subsp->x) {
10620cf1dd8SToby Isaac       for (j = 0; j < origDim; j++) inpoints[i * origDim + j] = subsp->x[j];
10720cf1dd8SToby Isaac     } else {
10820cf1dd8SToby Isaac       for (j = 0; j < origDim; j++) inpoints[i * origDim + j] = 0.0;
10920cf1dd8SToby Isaac     }
11020cf1dd8SToby Isaac     if (subsp->Jx) {
11120cf1dd8SToby Isaac       for (j = 0; j < origDim; j++) {
112ad540459SPierre Jolivet         for (k = 0; k < subDim; k++) inpoints[i * origDim + j] += subsp->Jx[j * subDim + k] * points[i * subDim + k];
11320cf1dd8SToby Isaac       }
11420cf1dd8SToby Isaac     } else {
115ad540459SPierre Jolivet       for (j = 0; j < PetscMin(subDim, origDim); j++) inpoints[i * origDim + j] += points[i * subDim + j];
11620cf1dd8SToby Isaac     }
11720cf1dd8SToby Isaac   }
11848a46eb9SPierre Jolivet   if (B) PetscCall(DMGetWorkArray(sp->dm, npoints * origNb * origNc, MPIU_REAL, &inB));
11948a46eb9SPierre Jolivet   if (D) PetscCall(DMGetWorkArray(sp->dm, npoints * origNb * origNc * origDim, MPIU_REAL, &inD));
12048a46eb9SPierre Jolivet   if (H) PetscCall(DMGetWorkArray(sp->dm, npoints * origNb * origNc * origDim * origDim, MPIU_REAL, &inH));
1219566063dSJacob Faibussowitsch   PetscCall(PetscSpaceEvaluate(origsp, npoints, inpoints, inB, inD, inH));
12220cf1dd8SToby Isaac   if (H) {
12320cf1dd8SToby Isaac     PetscReal *phi, *psi;
12420cf1dd8SToby Isaac 
1259566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(sp->dm, origNc * origDim * origDim, MPIU_REAL, &phi));
1269566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(sp->dm, origNc * subDim * subDim, MPIU_REAL, &psi));
12720cf1dd8SToby Isaac     for (i = 0; i < npoints * subNb * subNc * subDim; i++) D[i] = 0.0;
12820cf1dd8SToby Isaac     for (i = 0; i < subNb; i++) {
12920cf1dd8SToby Isaac       const PetscReal *subq = &subsp->Q[i * origNb];
13020cf1dd8SToby Isaac 
13120cf1dd8SToby Isaac       for (j = 0; j < npoints; j++) {
13220cf1dd8SToby Isaac         for (k = 0; k < origNc * origDim; k++) phi[k] = 0.;
13320cf1dd8SToby Isaac         for (k = 0; k < origNc * subDim; k++) psi[k] = 0.;
13420cf1dd8SToby Isaac         for (k = 0; k < origNb; k++) {
135ad540459SPierre Jolivet           for (l = 0; l < origNc * origDim * origDim; l++) phi[l] += inH[(j * origNb + k) * origNc * origDim * origDim + l] * subq[k];
13620cf1dd8SToby Isaac         }
13720cf1dd8SToby Isaac         if (subsp->Jx) {
13820cf1dd8SToby Isaac           for (k = 0; k < subNc; k++) {
13920cf1dd8SToby Isaac             for (l = 0; l < subDim; l++) {
14020cf1dd8SToby Isaac               for (m = 0; m < origDim; m++) {
14120cf1dd8SToby Isaac                 for (n = 0; n < subDim; n++) {
142ad540459SPierre 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];
14320cf1dd8SToby Isaac                 }
14420cf1dd8SToby Isaac               }
14520cf1dd8SToby Isaac             }
14620cf1dd8SToby Isaac           }
14720cf1dd8SToby Isaac         } else {
14820cf1dd8SToby Isaac           for (k = 0; k < subNc; k++) {
14920cf1dd8SToby Isaac             for (l = 0; l < PetscMin(subDim, origDim); l++) {
150ad540459SPierre Jolivet               for (m = 0; m < PetscMin(subDim, origDim); m++) psi[(k * subDim + l) * subDim + m] += phi[(k * origDim + l) * origDim + m];
15120cf1dd8SToby Isaac             }
15220cf1dd8SToby Isaac           }
15320cf1dd8SToby Isaac         }
15420cf1dd8SToby Isaac         if (subsp->Ju) {
15520cf1dd8SToby Isaac           for (k = 0; k < subNc; k++) {
15620cf1dd8SToby Isaac             for (l = 0; l < origNc; l++) {
157ad540459SPierre 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];
15820cf1dd8SToby Isaac             }
15920cf1dd8SToby Isaac           }
1609371c9d4SSatish Balay         } else {
16120cf1dd8SToby Isaac           for (k = 0; k < PetscMin(subNc, origNc); k++) {
162ad540459SPierre Jolivet             for (l = 0; l < subDim * subDim; l++) H[((j * subNb + i) * subNc + k) * subDim * subDim + l] += psi[k * subDim * subDim + l];
16320cf1dd8SToby Isaac           }
16420cf1dd8SToby Isaac         }
16520cf1dd8SToby Isaac       }
16620cf1dd8SToby Isaac     }
1679566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(sp->dm, subNc * origDim, MPIU_REAL, &psi));
1689566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(sp->dm, origNc * origDim, MPIU_REAL, &phi));
1699566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(sp->dm, npoints * origNb * origNc * origDim, MPIU_REAL, &inH));
17020cf1dd8SToby Isaac   }
17120cf1dd8SToby Isaac   if (D) {
17220cf1dd8SToby Isaac     PetscReal *phi, *psi;
17320cf1dd8SToby Isaac 
1749566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(sp->dm, origNc * origDim, MPIU_REAL, &phi));
1759566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(sp->dm, origNc * subDim, MPIU_REAL, &psi));
17620cf1dd8SToby Isaac     for (i = 0; i < npoints * subNb * subNc * subDim; i++) D[i] = 0.0;
17720cf1dd8SToby Isaac     for (i = 0; i < subNb; i++) {
17820cf1dd8SToby Isaac       const PetscReal *subq = &subsp->Q[i * origNb];
17920cf1dd8SToby Isaac 
18020cf1dd8SToby Isaac       for (j = 0; j < npoints; j++) {
18120cf1dd8SToby Isaac         for (k = 0; k < origNc * origDim; k++) phi[k] = 0.;
18220cf1dd8SToby Isaac         for (k = 0; k < origNc * subDim; k++) psi[k] = 0.;
18320cf1dd8SToby Isaac         for (k = 0; k < origNb; k++) {
184ad540459SPierre Jolivet           for (l = 0; l < origNc * origDim; l++) phi[l] += inD[(j * origNb + k) * origNc * origDim + l] * subq[k];
18520cf1dd8SToby Isaac         }
18620cf1dd8SToby Isaac         if (subsp->Jx) {
18720cf1dd8SToby Isaac           for (k = 0; k < subNc; k++) {
18820cf1dd8SToby Isaac             for (l = 0; l < subDim; l++) {
189ad540459SPierre Jolivet               for (m = 0; m < origDim; m++) psi[k * subDim + l] += subsp->Jx[m * subDim + l] * phi[k * origDim + m];
19020cf1dd8SToby Isaac             }
19120cf1dd8SToby Isaac           }
19220cf1dd8SToby Isaac         } else {
19320cf1dd8SToby Isaac           for (k = 0; k < subNc; k++) {
194ad540459SPierre Jolivet             for (l = 0; l < PetscMin(subDim, origDim); l++) psi[k * subDim + l] += phi[k * origDim + l];
19520cf1dd8SToby Isaac           }
19620cf1dd8SToby Isaac         }
19720cf1dd8SToby Isaac         if (subsp->Ju) {
19820cf1dd8SToby Isaac           for (k = 0; k < subNc; k++) {
19920cf1dd8SToby Isaac             for (l = 0; l < origNc; l++) {
200ad540459SPierre Jolivet               for (m = 0; m < subDim; m++) D[((j * subNb + i) * subNc + k) * subDim + m] += subsp->Ju[k * origNc + l] * psi[l * subDim + m];
20120cf1dd8SToby Isaac             }
20220cf1dd8SToby Isaac           }
2039371c9d4SSatish Balay         } else {
20420cf1dd8SToby Isaac           for (k = 0; k < PetscMin(subNc, origNc); k++) {
205ad540459SPierre Jolivet             for (l = 0; l < subDim; l++) D[((j * subNb + i) * subNc + k) * subDim + l] += psi[k * subDim + l];
20620cf1dd8SToby Isaac           }
20720cf1dd8SToby Isaac         }
20820cf1dd8SToby Isaac       }
20920cf1dd8SToby Isaac     }
2109566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(sp->dm, subNc * origDim, MPIU_REAL, &psi));
2119566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(sp->dm, origNc * origDim, MPIU_REAL, &phi));
2129566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(sp->dm, npoints * origNb * origNc * origDim, MPIU_REAL, &inD));
21320cf1dd8SToby Isaac   }
21420cf1dd8SToby Isaac   if (B) {
21520cf1dd8SToby Isaac     PetscReal *phi;
21620cf1dd8SToby Isaac 
2179566063dSJacob Faibussowitsch     PetscCall(DMGetWorkArray(sp->dm, origNc, MPIU_REAL, &phi));
21820cf1dd8SToby Isaac     if (subsp->u) {
21920cf1dd8SToby Isaac       for (i = 0; i < npoints * subNb; i++) {
22020cf1dd8SToby Isaac         for (j = 0; j < subNc; j++) B[i * subNc + j] = subsp->u[j];
22120cf1dd8SToby Isaac       }
22220cf1dd8SToby Isaac     } else {
22320cf1dd8SToby Isaac       for (i = 0; i < npoints * subNb * subNc; i++) B[i] = 0.0;
22420cf1dd8SToby Isaac     }
22520cf1dd8SToby Isaac     for (i = 0; i < subNb; i++) {
22620cf1dd8SToby Isaac       const PetscReal *subq = &subsp->Q[i * origNb];
22720cf1dd8SToby Isaac 
22820cf1dd8SToby Isaac       for (j = 0; j < npoints; j++) {
22920cf1dd8SToby Isaac         for (k = 0; k < origNc; k++) phi[k] = 0.;
23020cf1dd8SToby Isaac         for (k = 0; k < origNb; k++) {
231ad540459SPierre Jolivet           for (l = 0; l < origNc; l++) phi[l] += inB[(j * origNb + k) * origNc + l] * subq[k];
23220cf1dd8SToby Isaac         }
23320cf1dd8SToby Isaac         if (subsp->Ju) {
23420cf1dd8SToby Isaac           for (k = 0; k < subNc; k++) {
235ad540459SPierre Jolivet             for (l = 0; l < origNc; l++) B[(j * subNb + i) * subNc + k] += subsp->Ju[k * origNc + l] * phi[l];
23620cf1dd8SToby Isaac           }
2379371c9d4SSatish Balay         } else {
238ad540459SPierre Jolivet           for (k = 0; k < PetscMin(subNc, origNc); k++) B[(j * subNb + i) * subNc + k] += phi[k];
23920cf1dd8SToby Isaac         }
24020cf1dd8SToby Isaac       }
24120cf1dd8SToby Isaac     }
2429566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(sp->dm, origNc, MPIU_REAL, &phi));
2439566063dSJacob Faibussowitsch     PetscCall(DMRestoreWorkArray(sp->dm, npoints * origNb * origNc, MPIU_REAL, &inB));
24420cf1dd8SToby Isaac   }
2459566063dSJacob Faibussowitsch   PetscCall(DMRestoreWorkArray(sp->dm, npoints * origDim, MPIU_REAL, &inpoints));
24620cf1dd8SToby Isaac   PetscFunctionReturn(0);
24720cf1dd8SToby Isaac }
24820cf1dd8SToby Isaac 
2499371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Subspace(PetscSpace sp) {
25020cf1dd8SToby Isaac   PetscSpace_Subspace *subsp;
251362febeeSStefano Zampini 
252362febeeSStefano Zampini   PetscFunctionBegin;
253*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&subsp));
25420cf1dd8SToby Isaac   sp->data = (void *)subsp;
25520cf1dd8SToby Isaac   PetscFunctionReturn(0);
25620cf1dd8SToby Isaac }
25720cf1dd8SToby Isaac 
2589371c9d4SSatish Balay static PetscErrorCode PetscSpaceGetDimension_Subspace(PetscSpace sp, PetscInt *dim) {
25920cf1dd8SToby Isaac   PetscSpace_Subspace *subsp;
26020cf1dd8SToby Isaac 
26120cf1dd8SToby Isaac   PetscFunctionBegin;
26220cf1dd8SToby Isaac   subsp = (PetscSpace_Subspace *)sp->data;
26320cf1dd8SToby Isaac   *dim  = subsp->Nb;
26420cf1dd8SToby Isaac   PetscFunctionReturn(0);
26520cf1dd8SToby Isaac }
26620cf1dd8SToby Isaac 
2679371c9d4SSatish Balay static PetscErrorCode PetscSpaceSetUp_Subspace(PetscSpace sp) {
26820cf1dd8SToby Isaac   const PetscReal     *x;
26920cf1dd8SToby Isaac   const PetscReal     *Jx;
27020cf1dd8SToby Isaac   const PetscReal     *u;
27120cf1dd8SToby Isaac   const PetscReal     *Ju;
27220cf1dd8SToby Isaac   PetscDualSpace       dualSubspace;
27320cf1dd8SToby Isaac   PetscSpace           origSpace;
27420cf1dd8SToby Isaac   PetscInt             origDim, subDim, origNc, subNc, origNb, subNb, f, i, j, numPoints, offset;
27520cf1dd8SToby Isaac   PetscReal           *allPoints, *allWeights, *B, *V;
27620cf1dd8SToby Isaac   DM                   dm;
27720cf1dd8SToby Isaac   PetscSpace_Subspace *subsp;
27820cf1dd8SToby Isaac 
27920cf1dd8SToby Isaac   PetscFunctionBegin;
28020cf1dd8SToby Isaac   subsp        = (PetscSpace_Subspace *)sp->data;
28120cf1dd8SToby Isaac   x            = subsp->x;
28220cf1dd8SToby Isaac   Jx           = subsp->Jx;
28320cf1dd8SToby Isaac   u            = subsp->u;
28420cf1dd8SToby Isaac   Ju           = subsp->Ju;
28520cf1dd8SToby Isaac   origSpace    = subsp->origSpace;
28620cf1dd8SToby Isaac   dualSubspace = subsp->dualSubspace;
2879566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumComponents(origSpace, &origNc));
2889566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumVariables(origSpace, &origDim));
2899566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(dualSubspace, &dm));
2909566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &subDim));
2919566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDimension(origSpace, &origNb));
2929566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(dualSubspace, &subNb));
2939566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetNumComponents(dualSubspace, &subNc));
29420cf1dd8SToby Isaac 
29520cf1dd8SToby Isaac   for (f = 0, numPoints = 0; f < subNb; f++) {
29620cf1dd8SToby Isaac     PetscQuadrature q;
29720cf1dd8SToby Isaac     PetscInt        qNp;
29820cf1dd8SToby Isaac 
2999566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetFunctional(dualSubspace, f, &q));
3009566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(q, NULL, NULL, &qNp, NULL, NULL));
30120cf1dd8SToby Isaac     numPoints += qNp;
30220cf1dd8SToby Isaac   }
3039566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(subNb * origNb, &V));
3049566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(numPoints * origDim, &allPoints, numPoints * origNc, &allWeights, numPoints * origNb * origNc, &B));
30520cf1dd8SToby Isaac   for (f = 0, offset = 0; f < subNb; f++) {
30620cf1dd8SToby Isaac     PetscQuadrature  q;
30720cf1dd8SToby Isaac     PetscInt         qNp, p;
30820cf1dd8SToby Isaac     const PetscReal *qp;
30920cf1dd8SToby Isaac     const PetscReal *qw;
31020cf1dd8SToby Isaac 
3119566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetFunctional(dualSubspace, f, &q));
3129566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(q, NULL, NULL, &qNp, &qp, &qw));
31320cf1dd8SToby Isaac     for (p = 0; p < qNp; p++, offset++) {
31420cf1dd8SToby Isaac       if (x) {
31520cf1dd8SToby Isaac         for (i = 0; i < origDim; i++) allPoints[origDim * offset + i] = x[i];
31620cf1dd8SToby Isaac       } else {
31720cf1dd8SToby Isaac         for (i = 0; i < origDim; i++) allPoints[origDim * offset + i] = 0.0;
31820cf1dd8SToby Isaac       }
31920cf1dd8SToby Isaac       if (Jx) {
32020cf1dd8SToby Isaac         for (i = 0; i < origDim; i++) {
321ad540459SPierre Jolivet           for (j = 0; j < subDim; j++) allPoints[origDim * offset + i] += Jx[i * subDim + j] * qp[j];
32220cf1dd8SToby Isaac         }
32320cf1dd8SToby Isaac       } else {
32420cf1dd8SToby Isaac         for (i = 0; i < PetscMin(subDim, origDim); i++) allPoints[origDim * offset + i] += qp[i];
32520cf1dd8SToby Isaac       }
32620cf1dd8SToby Isaac       for (i = 0; i < origNc; i++) allWeights[origNc * offset + i] = 0.0;
32720cf1dd8SToby Isaac       if (Ju) {
32820cf1dd8SToby Isaac         for (i = 0; i < origNc; i++) {
329ad540459SPierre Jolivet           for (j = 0; j < subNc; j++) allWeights[offset * origNc + i] += qw[j] * Ju[j * origNc + i];
33020cf1dd8SToby Isaac         }
33120cf1dd8SToby Isaac       } else {
33220cf1dd8SToby Isaac         for (i = 0; i < PetscMin(subNc, origNc); i++) allWeights[offset * origNc + i] += qw[i];
33320cf1dd8SToby Isaac       }
33420cf1dd8SToby Isaac     }
33520cf1dd8SToby Isaac   }
3369566063dSJacob Faibussowitsch   PetscCall(PetscSpaceEvaluate(origSpace, numPoints, allPoints, B, NULL, NULL));
33720cf1dd8SToby Isaac   for (f = 0, offset = 0; f < subNb; f++) {
33820cf1dd8SToby Isaac     PetscInt         b, p, s, qNp;
33920cf1dd8SToby Isaac     PetscQuadrature  q;
34020cf1dd8SToby Isaac     const PetscReal *qw;
34120cf1dd8SToby Isaac 
3429566063dSJacob Faibussowitsch     PetscCall(PetscDualSpaceGetFunctional(dualSubspace, f, &q));
3439566063dSJacob Faibussowitsch     PetscCall(PetscQuadratureGetData(q, NULL, NULL, &qNp, NULL, &qw));
34420cf1dd8SToby Isaac     if (u) {
34520cf1dd8SToby Isaac       for (b = 0; b < origNb; b++) {
346ad540459SPierre Jolivet         for (s = 0; s < subNc; s++) V[f * origNb + b] += qw[s] * u[s];
34720cf1dd8SToby Isaac       }
34820cf1dd8SToby Isaac     } else {
34920cf1dd8SToby Isaac       for (b = 0; b < origNb; b++) V[f * origNb + b] = 0.0;
35020cf1dd8SToby Isaac     }
35120cf1dd8SToby Isaac     for (p = 0; p < qNp; p++, offset++) {
35220cf1dd8SToby Isaac       for (b = 0; b < origNb; b++) {
353ad540459SPierre Jolivet         for (s = 0; s < origNc; s++) V[f * origNb + b] += B[(offset * origNb + b) * origNc + s] * allWeights[offset * origNc + s];
35420cf1dd8SToby Isaac       }
35520cf1dd8SToby Isaac     }
35620cf1dd8SToby Isaac   }
35720cf1dd8SToby Isaac   /* orthnormalize rows of V */
35820cf1dd8SToby Isaac   for (f = 0; f < subNb; f++) {
35920cf1dd8SToby Isaac     PetscReal rho = 0.0, scal;
36020cf1dd8SToby Isaac 
36120cf1dd8SToby Isaac     for (i = 0; i < origNb; i++) rho += PetscSqr(V[f * origNb + i]);
36220cf1dd8SToby Isaac 
36320cf1dd8SToby Isaac     scal = 1. / PetscSqrtReal(rho);
36420cf1dd8SToby Isaac 
36520cf1dd8SToby Isaac     for (i = 0; i < origNb; i++) V[f * origNb + i] *= scal;
36620cf1dd8SToby Isaac     for (j = f + 1; j < subNb; j++) {
36720cf1dd8SToby Isaac       for (i = 0, scal = 0.; i < origNb; i++) scal += V[f * origNb + i] * V[j * origNb + i];
36820cf1dd8SToby Isaac       for (i = 0; i < origNb; i++) V[j * origNb + i] -= V[f * origNb + i] * scal;
36920cf1dd8SToby Isaac     }
37020cf1dd8SToby Isaac   }
3719566063dSJacob Faibussowitsch   PetscCall(PetscFree3(allPoints, allWeights, B));
37220cf1dd8SToby Isaac   subsp->Q = V;
37320cf1dd8SToby Isaac   PetscFunctionReturn(0);
37420cf1dd8SToby Isaac }
37520cf1dd8SToby Isaac 
3769371c9d4SSatish Balay static PetscErrorCode PetscSpacePolynomialGetTensor_Subspace(PetscSpace sp, PetscBool *poly) {
37720cf1dd8SToby Isaac   PetscSpace_Subspace *subsp = (PetscSpace_Subspace *)sp->data;
37820cf1dd8SToby Isaac 
37920cf1dd8SToby Isaac   PetscFunctionBegin;
38020cf1dd8SToby Isaac   *poly = PETSC_FALSE;
3819566063dSJacob Faibussowitsch   PetscCall(PetscSpacePolynomialGetTensor(subsp->origSpace, poly));
38220cf1dd8SToby Isaac   if (*poly) {
38320cf1dd8SToby Isaac     if (subsp->Jx) {
38420cf1dd8SToby Isaac       PetscInt subDim, origDim, i, j;
38520cf1dd8SToby Isaac       PetscInt maxnnz;
38620cf1dd8SToby Isaac 
3879566063dSJacob Faibussowitsch       PetscCall(PetscSpaceGetNumVariables(subsp->origSpace, &origDim));
3889566063dSJacob Faibussowitsch       PetscCall(PetscSpaceGetNumVariables(sp, &subDim));
38920cf1dd8SToby Isaac       maxnnz = 0;
39020cf1dd8SToby Isaac       for (i = 0; i < origDim; i++) {
39120cf1dd8SToby Isaac         PetscInt nnz = 0;
39220cf1dd8SToby Isaac 
39320cf1dd8SToby Isaac         for (j = 0; j < subDim; j++) nnz += (subsp->Jx[i * subDim + j] != 0.);
39420cf1dd8SToby Isaac         maxnnz = PetscMax(maxnnz, nnz);
39520cf1dd8SToby Isaac       }
39620cf1dd8SToby Isaac       for (j = 0; j < subDim; j++) {
39720cf1dd8SToby Isaac         PetscInt nnz = 0;
39820cf1dd8SToby Isaac 
39920cf1dd8SToby Isaac         for (i = 0; i < origDim; i++) nnz += (subsp->Jx[i * subDim + j] != 0.);
40020cf1dd8SToby Isaac         maxnnz = PetscMax(maxnnz, nnz);
40120cf1dd8SToby Isaac       }
40220cf1dd8SToby Isaac       if (maxnnz > 1) *poly = PETSC_FALSE;
40320cf1dd8SToby Isaac     }
40420cf1dd8SToby Isaac   }
40520cf1dd8SToby Isaac   PetscFunctionReturn(0);
40620cf1dd8SToby Isaac }
40720cf1dd8SToby Isaac 
4089371c9d4SSatish Balay static PetscErrorCode PetscSpaceInitialize_Subspace(PetscSpace sp) {
40920cf1dd8SToby Isaac   PetscFunctionBegin;
41020cf1dd8SToby Isaac   sp->ops->setup        = PetscSpaceSetUp_Subspace;
41120cf1dd8SToby Isaac   sp->ops->view         = PetscSpaceView_Subspace;
41220cf1dd8SToby Isaac   sp->ops->destroy      = PetscSpaceDestroy_Subspace;
41320cf1dd8SToby Isaac   sp->ops->getdimension = PetscSpaceGetDimension_Subspace;
41420cf1dd8SToby Isaac   sp->ops->evaluate     = PetscSpaceEvaluate_Subspace;
4159566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Subspace));
41620cf1dd8SToby Isaac   PetscFunctionReturn(0);
41720cf1dd8SToby Isaac }
41820cf1dd8SToby Isaac 
4199371c9d4SSatish Balay PetscErrorCode PetscSpaceCreateSubspace(PetscSpace origSpace, PetscDualSpace dualSubspace, PetscReal *x, PetscReal *Jx, PetscReal *u, PetscReal *Ju, PetscCopyMode copymode, PetscSpace *subspace) {
42020cf1dd8SToby Isaac   PetscSpace_Subspace *subsp;
42120cf1dd8SToby Isaac   PetscInt             origDim, subDim, origNc, subNc, subNb;
42220cf1dd8SToby Isaac   PetscInt             order;
42320cf1dd8SToby Isaac   DM                   dm;
42420cf1dd8SToby Isaac 
42520cf1dd8SToby Isaac   PetscFunctionBegin;
42620cf1dd8SToby Isaac   PetscValidHeaderSpecific(origSpace, PETSCSPACE_CLASSID, 1);
42720cf1dd8SToby Isaac   PetscValidHeaderSpecific(dualSubspace, PETSCDUALSPACE_CLASSID, 2);
42820cf1dd8SToby Isaac   if (x) PetscValidRealPointer(x, 3);
42920cf1dd8SToby Isaac   if (Jx) PetscValidRealPointer(Jx, 4);
43020cf1dd8SToby Isaac   if (u) PetscValidRealPointer(u, 5);
43120cf1dd8SToby Isaac   if (Ju) PetscValidRealPointer(Ju, 6);
432064a246eSJacob Faibussowitsch   PetscValidPointer(subspace, 8);
4339566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumComponents(origSpace, &origNc));
4349566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetNumVariables(origSpace, &origDim));
4359566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDM(dualSubspace, &dm));
4369566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &subDim));
4379566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetDimension(dualSubspace, &subNb));
4389566063dSJacob Faibussowitsch   PetscCall(PetscDualSpaceGetNumComponents(dualSubspace, &subNc));
4399566063dSJacob Faibussowitsch   PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)origSpace), subspace));
4409566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetType(*subspace, PETSCSPACESUBSPACE));
4419566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetNumVariables(*subspace, subDim));
4429566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetNumComponents(*subspace, subNc));
4439566063dSJacob Faibussowitsch   PetscCall(PetscSpaceGetDegree(origSpace, &order, NULL));
4449566063dSJacob Faibussowitsch   PetscCall(PetscSpaceSetDegree(*subspace, order, PETSC_DETERMINE));
44520cf1dd8SToby Isaac   subsp     = (PetscSpace_Subspace *)(*subspace)->data;
44620cf1dd8SToby Isaac   subsp->Nb = subNb;
44720cf1dd8SToby Isaac   switch (copymode) {
44820cf1dd8SToby Isaac   case PETSC_OWN_POINTER:
44920cf1dd8SToby Isaac     if (x) subsp->x_alloc = x;
45020cf1dd8SToby Isaac     if (Jx) subsp->Jx_alloc = Jx;
45120cf1dd8SToby Isaac     if (u) subsp->u_alloc = u;
45220cf1dd8SToby Isaac     if (Ju) subsp->Ju_alloc = Ju;
45320cf1dd8SToby Isaac   case PETSC_USE_POINTER:
45420cf1dd8SToby Isaac     if (x) subsp->x = x;
45520cf1dd8SToby Isaac     if (Jx) subsp->Jx = Jx;
45620cf1dd8SToby Isaac     if (u) subsp->u = u;
45720cf1dd8SToby Isaac     if (Ju) subsp->Ju = Ju;
45820cf1dd8SToby Isaac     break;
45920cf1dd8SToby Isaac   case PETSC_COPY_VALUES:
46020cf1dd8SToby Isaac     if (x) {
4619566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(origDim, &subsp->x_alloc));
4629566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(subsp->x_alloc, x, origDim));
46320cf1dd8SToby Isaac       subsp->x = subsp->x_alloc;
46420cf1dd8SToby Isaac     }
46520cf1dd8SToby Isaac     if (Jx) {
4669566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(origDim * subDim, &subsp->Jx_alloc));
4679566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(subsp->Jx_alloc, Jx, origDim * subDim));
46820cf1dd8SToby Isaac       subsp->Jx = subsp->Jx_alloc;
46920cf1dd8SToby Isaac     }
47020cf1dd8SToby Isaac     if (u) {
4719566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(subNc, &subsp->u_alloc));
4729566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(subsp->u_alloc, u, subNc));
47320cf1dd8SToby Isaac       subsp->u = subsp->u_alloc;
47420cf1dd8SToby Isaac     }
47520cf1dd8SToby Isaac     if (Ju) {
4769566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(origNc * subNc, &subsp->Ju_alloc));
4779566063dSJacob Faibussowitsch       PetscCall(PetscArraycpy(subsp->Ju_alloc, Ju, origNc * subNc));
47820cf1dd8SToby Isaac       subsp->Ju = subsp->Ju_alloc;
47920cf1dd8SToby Isaac     }
48020cf1dd8SToby Isaac     break;
4819371c9d4SSatish Balay   default: SETERRQ(PetscObjectComm((PetscObject)origSpace), PETSC_ERR_ARG_OUTOFRANGE, "Unknown copy mode");
48220cf1dd8SToby Isaac   }
4839566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)origSpace));
48420cf1dd8SToby Isaac   subsp->origSpace = origSpace;
4859566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)dualSubspace));
48620cf1dd8SToby Isaac   subsp->dualSubspace = dualSubspace;
4879566063dSJacob Faibussowitsch   PetscCall(PetscSpaceInitialize_Subspace(*subspace));
48820cf1dd8SToby Isaac   PetscFunctionReturn(0);
48920cf1dd8SToby Isaac }
490