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; 162dce792eSToby Isaac PetscBool setupcalled; 1720cf1dd8SToby Isaac } PetscSpace_Subspace; 1820cf1dd8SToby Isaac 19d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceDestroy_Subspace(PetscSpace sp) 20d71ae5a4SJacob Faibussowitsch { 2120cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 2220cf1dd8SToby Isaac 2320cf1dd8SToby Isaac PetscFunctionBegin; 2420cf1dd8SToby Isaac subsp = (PetscSpace_Subspace *)sp->data; 2520cf1dd8SToby Isaac subsp->x = NULL; 269566063dSJacob Faibussowitsch PetscCall(PetscFree(subsp->x_alloc)); 2720cf1dd8SToby Isaac subsp->Jx = NULL; 289566063dSJacob Faibussowitsch PetscCall(PetscFree(subsp->Jx_alloc)); 2920cf1dd8SToby Isaac subsp->u = NULL; 309566063dSJacob Faibussowitsch PetscCall(PetscFree(subsp->u_alloc)); 3120cf1dd8SToby Isaac subsp->Ju = NULL; 329566063dSJacob Faibussowitsch PetscCall(PetscFree(subsp->Ju_alloc)); 339566063dSJacob Faibussowitsch PetscCall(PetscFree(subsp->Q)); 349566063dSJacob Faibussowitsch PetscCall(PetscSpaceDestroy(&subsp->origSpace)); 359566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceDestroy(&subsp->dualSubspace)); 369566063dSJacob Faibussowitsch PetscCall(PetscFree(subsp)); 3720cf1dd8SToby Isaac sp->data = NULL; 389566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePolynomialGetTensor_C", NULL)); 393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4020cf1dd8SToby Isaac } 4120cf1dd8SToby Isaac 42d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceView_Subspace(PetscSpace sp, PetscViewer viewer) 43d71ae5a4SJacob Faibussowitsch { 44*9f196a02SMartin Diehl PetscBool isascii; 4520cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 4620cf1dd8SToby Isaac 4720cf1dd8SToby Isaac PetscFunctionBegin; 4820cf1dd8SToby Isaac subsp = (PetscSpace_Subspace *)sp->data; 49*9f196a02SMartin Diehl PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii)); 50*9f196a02SMartin Diehl if (isascii) { 5120cf1dd8SToby Isaac PetscInt origDim, subDim, origNc, subNc, o, s; 5220cf1dd8SToby Isaac 539566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumVariables(subsp->origSpace, &origDim)); 549566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumComponents(subsp->origSpace, &origNc)); 559566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumVariables(sp, &subDim)); 569566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumComponents(sp, &subNc)); 5720cf1dd8SToby Isaac if (subsp->x) { 589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Subspace-to-space domain shift:\n\n")); 5948a46eb9SPierre Jolivet for (o = 0; o < origDim; o++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g\n", (double)subsp->x[o])); 6020cf1dd8SToby Isaac } 6120cf1dd8SToby Isaac if (subsp->Jx) { 629566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Subspace-to-space domain transform:\n\n")); 6320cf1dd8SToby Isaac for (o = 0; o < origDim; o++) { 649566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)subsp->Jx[o * subDim + 0])); 659566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 6648a46eb9SPierre Jolivet for (s = 1; s < subDim; s++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)subsp->Jx[o * subDim + s])); 679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 6920cf1dd8SToby Isaac } 7020cf1dd8SToby Isaac } 7120cf1dd8SToby Isaac if (subsp->u) { 729566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Space-to-subspace range shift:\n\n")); 7348a46eb9SPierre Jolivet for (o = 0; o < origNc; o++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g\n", (double)subsp->u[o])); 7420cf1dd8SToby Isaac } 7520cf1dd8SToby Isaac if (subsp->Ju) { 769566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Space-to-subsace domain transform:\n")); 7720cf1dd8SToby Isaac for (o = 0; o < origNc; o++) { 789566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_FALSE)); 7948a46eb9SPierre Jolivet for (s = 0; s < subNc; s++) PetscCall(PetscViewerASCIIPrintf(viewer, " %g", (double)subsp->Ju[o * subNc + s])); 809566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIUseTabs(viewer, PETSC_TRUE)); 8120cf1dd8SToby Isaac } 829566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 8320cf1dd8SToby Isaac } 849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Original space:\n")); 8520cf1dd8SToby Isaac } 869566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 879566063dSJacob Faibussowitsch PetscCall(PetscSpaceView(subsp->origSpace, viewer)); 889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9020cf1dd8SToby Isaac } 9120cf1dd8SToby Isaac 92d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceEvaluate_Subspace(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) 93d71ae5a4SJacob Faibussowitsch { 9420cf1dd8SToby Isaac PetscSpace_Subspace *subsp = (PetscSpace_Subspace *)sp->data; 9520cf1dd8SToby Isaac PetscSpace origsp; 9620cf1dd8SToby Isaac PetscInt origDim, subDim, origNc, subNc, subNb, origNb, i, j, k, l, m, n, o; 9720cf1dd8SToby Isaac PetscReal *inpoints, *inB = NULL, *inD = NULL, *inH = NULL; 9820cf1dd8SToby Isaac 9920cf1dd8SToby Isaac PetscFunctionBegin; 10020cf1dd8SToby Isaac origsp = subsp->origSpace; 1019566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumVariables(sp, &subDim)); 1029566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumVariables(origsp, &origDim)); 1039566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumComponents(sp, &subNc)); 1049566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumComponents(origsp, &origNc)); 1059566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDimension(sp, &subNb)); 1069566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDimension(origsp, &origNb)); 1079566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(sp->dm, npoints * origDim, MPIU_REAL, &inpoints)); 10820cf1dd8SToby Isaac for (i = 0; i < npoints; i++) { 10920cf1dd8SToby Isaac if (subsp->x) { 11020cf1dd8SToby Isaac for (j = 0; j < origDim; j++) inpoints[i * origDim + j] = subsp->x[j]; 11120cf1dd8SToby Isaac } else { 11220cf1dd8SToby Isaac for (j = 0; j < origDim; j++) inpoints[i * origDim + j] = 0.0; 11320cf1dd8SToby Isaac } 11420cf1dd8SToby Isaac if (subsp->Jx) { 11520cf1dd8SToby Isaac for (j = 0; j < origDim; j++) { 116ad540459SPierre Jolivet for (k = 0; k < subDim; k++) inpoints[i * origDim + j] += subsp->Jx[j * subDim + k] * points[i * subDim + k]; 11720cf1dd8SToby Isaac } 11820cf1dd8SToby Isaac } else { 119ad540459SPierre Jolivet for (j = 0; j < PetscMin(subDim, origDim); j++) inpoints[i * origDim + j] += points[i * subDim + j]; 12020cf1dd8SToby Isaac } 12120cf1dd8SToby Isaac } 12248a46eb9SPierre Jolivet if (B) PetscCall(DMGetWorkArray(sp->dm, npoints * origNb * origNc, MPIU_REAL, &inB)); 12348a46eb9SPierre Jolivet if (D) PetscCall(DMGetWorkArray(sp->dm, npoints * origNb * origNc * origDim, MPIU_REAL, &inD)); 12448a46eb9SPierre Jolivet if (H) PetscCall(DMGetWorkArray(sp->dm, npoints * origNb * origNc * origDim * origDim, MPIU_REAL, &inH)); 1259566063dSJacob Faibussowitsch PetscCall(PetscSpaceEvaluate(origsp, npoints, inpoints, inB, inD, inH)); 12620cf1dd8SToby Isaac if (H) { 12720cf1dd8SToby Isaac PetscReal *phi, *psi; 12820cf1dd8SToby Isaac 1299566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(sp->dm, origNc * origDim * origDim, MPIU_REAL, &phi)); 1309566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(sp->dm, origNc * subDim * subDim, MPIU_REAL, &psi)); 13120cf1dd8SToby Isaac for (i = 0; i < npoints * subNb * subNc * subDim; i++) D[i] = 0.0; 13220cf1dd8SToby Isaac for (i = 0; i < subNb; i++) { 13320cf1dd8SToby Isaac const PetscReal *subq = &subsp->Q[i * origNb]; 13420cf1dd8SToby Isaac 13520cf1dd8SToby Isaac for (j = 0; j < npoints; j++) { 13620cf1dd8SToby Isaac for (k = 0; k < origNc * origDim; k++) phi[k] = 0.; 13720cf1dd8SToby Isaac for (k = 0; k < origNc * subDim; k++) psi[k] = 0.; 13820cf1dd8SToby Isaac for (k = 0; k < origNb; k++) { 139ad540459SPierre Jolivet for (l = 0; l < origNc * origDim * origDim; l++) phi[l] += inH[(j * origNb + k) * origNc * origDim * origDim + l] * subq[k]; 14020cf1dd8SToby Isaac } 14120cf1dd8SToby Isaac if (subsp->Jx) { 14220cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 14320cf1dd8SToby Isaac for (l = 0; l < subDim; l++) { 14420cf1dd8SToby Isaac for (m = 0; m < origDim; m++) { 14520cf1dd8SToby Isaac for (n = 0; n < subDim; n++) { 146ad540459SPierre 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]; 14720cf1dd8SToby Isaac } 14820cf1dd8SToby Isaac } 14920cf1dd8SToby Isaac } 15020cf1dd8SToby Isaac } 15120cf1dd8SToby Isaac } else { 15220cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 15320cf1dd8SToby Isaac for (l = 0; l < PetscMin(subDim, origDim); l++) { 154ad540459SPierre Jolivet for (m = 0; m < PetscMin(subDim, origDim); m++) psi[(k * subDim + l) * subDim + m] += phi[(k * origDim + l) * origDim + m]; 15520cf1dd8SToby Isaac } 15620cf1dd8SToby Isaac } 15720cf1dd8SToby Isaac } 15820cf1dd8SToby Isaac if (subsp->Ju) { 15920cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 16020cf1dd8SToby Isaac for (l = 0; l < origNc; l++) { 161ad540459SPierre 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]; 16220cf1dd8SToby Isaac } 16320cf1dd8SToby Isaac } 1649371c9d4SSatish Balay } else { 16520cf1dd8SToby Isaac for (k = 0; k < PetscMin(subNc, origNc); k++) { 166ad540459SPierre Jolivet for (l = 0; l < subDim * subDim; l++) H[((j * subNb + i) * subNc + k) * subDim * subDim + l] += psi[k * subDim * subDim + l]; 16720cf1dd8SToby Isaac } 16820cf1dd8SToby Isaac } 16920cf1dd8SToby Isaac } 17020cf1dd8SToby Isaac } 1719566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(sp->dm, subNc * origDim, MPIU_REAL, &psi)); 1729566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(sp->dm, origNc * origDim, MPIU_REAL, &phi)); 1739566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(sp->dm, npoints * origNb * origNc * origDim, MPIU_REAL, &inH)); 17420cf1dd8SToby Isaac } 17520cf1dd8SToby Isaac if (D) { 17620cf1dd8SToby Isaac PetscReal *phi, *psi; 17720cf1dd8SToby Isaac 1789566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(sp->dm, origNc * origDim, MPIU_REAL, &phi)); 1799566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(sp->dm, origNc * subDim, MPIU_REAL, &psi)); 18020cf1dd8SToby Isaac for (i = 0; i < npoints * subNb * subNc * subDim; i++) D[i] = 0.0; 18120cf1dd8SToby Isaac for (i = 0; i < subNb; i++) { 18220cf1dd8SToby Isaac const PetscReal *subq = &subsp->Q[i * origNb]; 18320cf1dd8SToby Isaac 18420cf1dd8SToby Isaac for (j = 0; j < npoints; j++) { 18520cf1dd8SToby Isaac for (k = 0; k < origNc * origDim; k++) phi[k] = 0.; 18620cf1dd8SToby Isaac for (k = 0; k < origNc * subDim; k++) psi[k] = 0.; 18720cf1dd8SToby Isaac for (k = 0; k < origNb; k++) { 188ad540459SPierre Jolivet for (l = 0; l < origNc * origDim; l++) phi[l] += inD[(j * origNb + k) * origNc * origDim + l] * subq[k]; 18920cf1dd8SToby Isaac } 19020cf1dd8SToby Isaac if (subsp->Jx) { 19120cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 19220cf1dd8SToby Isaac for (l = 0; l < subDim; l++) { 193ad540459SPierre Jolivet for (m = 0; m < origDim; m++) psi[k * subDim + l] += subsp->Jx[m * subDim + l] * phi[k * origDim + m]; 19420cf1dd8SToby Isaac } 19520cf1dd8SToby Isaac } 19620cf1dd8SToby Isaac } else { 19720cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 198ad540459SPierre Jolivet for (l = 0; l < PetscMin(subDim, origDim); l++) psi[k * subDim + l] += phi[k * origDim + l]; 19920cf1dd8SToby Isaac } 20020cf1dd8SToby Isaac } 20120cf1dd8SToby Isaac if (subsp->Ju) { 20220cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 20320cf1dd8SToby Isaac for (l = 0; l < origNc; l++) { 204ad540459SPierre Jolivet for (m = 0; m < subDim; m++) D[((j * subNb + i) * subNc + k) * subDim + m] += subsp->Ju[k * origNc + l] * psi[l * subDim + m]; 20520cf1dd8SToby Isaac } 20620cf1dd8SToby Isaac } 2079371c9d4SSatish Balay } else { 20820cf1dd8SToby Isaac for (k = 0; k < PetscMin(subNc, origNc); k++) { 209ad540459SPierre Jolivet for (l = 0; l < subDim; l++) D[((j * subNb + i) * subNc + k) * subDim + l] += psi[k * subDim + l]; 21020cf1dd8SToby Isaac } 21120cf1dd8SToby Isaac } 21220cf1dd8SToby Isaac } 21320cf1dd8SToby Isaac } 2149566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(sp->dm, subNc * origDim, MPIU_REAL, &psi)); 2159566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(sp->dm, origNc * origDim, MPIU_REAL, &phi)); 2169566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(sp->dm, npoints * origNb * origNc * origDim, MPIU_REAL, &inD)); 21720cf1dd8SToby Isaac } 21820cf1dd8SToby Isaac if (B) { 21920cf1dd8SToby Isaac PetscReal *phi; 22020cf1dd8SToby Isaac 2219566063dSJacob Faibussowitsch PetscCall(DMGetWorkArray(sp->dm, origNc, MPIU_REAL, &phi)); 22220cf1dd8SToby Isaac if (subsp->u) { 22320cf1dd8SToby Isaac for (i = 0; i < npoints * subNb; i++) { 22420cf1dd8SToby Isaac for (j = 0; j < subNc; j++) B[i * subNc + j] = subsp->u[j]; 22520cf1dd8SToby Isaac } 22620cf1dd8SToby Isaac } else { 22720cf1dd8SToby Isaac for (i = 0; i < npoints * subNb * subNc; i++) B[i] = 0.0; 22820cf1dd8SToby Isaac } 22920cf1dd8SToby Isaac for (i = 0; i < subNb; i++) { 23020cf1dd8SToby Isaac const PetscReal *subq = &subsp->Q[i * origNb]; 23120cf1dd8SToby Isaac 23220cf1dd8SToby Isaac for (j = 0; j < npoints; j++) { 23320cf1dd8SToby Isaac for (k = 0; k < origNc; k++) phi[k] = 0.; 23420cf1dd8SToby Isaac for (k = 0; k < origNb; k++) { 235ad540459SPierre Jolivet for (l = 0; l < origNc; l++) phi[l] += inB[(j * origNb + k) * origNc + l] * subq[k]; 23620cf1dd8SToby Isaac } 23720cf1dd8SToby Isaac if (subsp->Ju) { 23820cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 239ad540459SPierre Jolivet for (l = 0; l < origNc; l++) B[(j * subNb + i) * subNc + k] += subsp->Ju[k * origNc + l] * phi[l]; 24020cf1dd8SToby Isaac } 2419371c9d4SSatish Balay } else { 242ad540459SPierre Jolivet for (k = 0; k < PetscMin(subNc, origNc); k++) B[(j * subNb + i) * subNc + k] += phi[k]; 24320cf1dd8SToby Isaac } 24420cf1dd8SToby Isaac } 24520cf1dd8SToby Isaac } 2469566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(sp->dm, origNc, MPIU_REAL, &phi)); 2479566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(sp->dm, npoints * origNb * origNc, MPIU_REAL, &inB)); 24820cf1dd8SToby Isaac } 2499566063dSJacob Faibussowitsch PetscCall(DMRestoreWorkArray(sp->dm, npoints * origDim, MPIU_REAL, &inpoints)); 2503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 25120cf1dd8SToby Isaac } 25220cf1dd8SToby Isaac 253d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Subspace(PetscSpace sp) 254d71ae5a4SJacob Faibussowitsch { 25520cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 256362febeeSStefano Zampini 257362febeeSStefano Zampini PetscFunctionBegin; 2584dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&subsp)); 25920cf1dd8SToby Isaac sp->data = (void *)subsp; 2603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 26120cf1dd8SToby Isaac } 26220cf1dd8SToby Isaac 263d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceGetDimension_Subspace(PetscSpace sp, PetscInt *dim) 264d71ae5a4SJacob Faibussowitsch { 26520cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 26620cf1dd8SToby Isaac 26720cf1dd8SToby Isaac PetscFunctionBegin; 26820cf1dd8SToby Isaac subsp = (PetscSpace_Subspace *)sp->data; 26920cf1dd8SToby Isaac *dim = subsp->Nb; 2703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 27120cf1dd8SToby Isaac } 27220cf1dd8SToby Isaac 273d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceSetUp_Subspace(PetscSpace sp) 274d71ae5a4SJacob Faibussowitsch { 27520cf1dd8SToby Isaac const PetscReal *x; 27620cf1dd8SToby Isaac const PetscReal *Jx; 27720cf1dd8SToby Isaac const PetscReal *u; 27820cf1dd8SToby Isaac const PetscReal *Ju; 27920cf1dd8SToby Isaac PetscDualSpace dualSubspace; 28020cf1dd8SToby Isaac PetscSpace origSpace; 28120cf1dd8SToby Isaac PetscInt origDim, subDim, origNc, subNc, origNb, subNb, f, i, j, numPoints, offset; 28220cf1dd8SToby Isaac PetscReal *allPoints, *allWeights, *B, *V; 28320cf1dd8SToby Isaac DM dm; 2842dce792eSToby Isaac PetscSpace_Subspace *subsp = (PetscSpace_Subspace *)sp->data; 28520cf1dd8SToby Isaac 28620cf1dd8SToby Isaac PetscFunctionBegin; 2872dce792eSToby Isaac if (subsp->setupcalled) PetscFunctionReturn(PETSC_SUCCESS); 2882dce792eSToby Isaac subsp->setupcalled = PETSC_TRUE; 28920cf1dd8SToby Isaac x = subsp->x; 29020cf1dd8SToby Isaac Jx = subsp->Jx; 29120cf1dd8SToby Isaac u = subsp->u; 29220cf1dd8SToby Isaac Ju = subsp->Ju; 29320cf1dd8SToby Isaac origSpace = subsp->origSpace; 29420cf1dd8SToby Isaac dualSubspace = subsp->dualSubspace; 2959566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumComponents(origSpace, &origNc)); 2969566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumVariables(origSpace, &origDim)); 2979566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(dualSubspace, &dm)); 2989566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &subDim)); 2999566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDimension(origSpace, &origNb)); 3009566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(dualSubspace, &subNb)); 3019566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetNumComponents(dualSubspace, &subNc)); 30220cf1dd8SToby Isaac 30320cf1dd8SToby Isaac for (f = 0, numPoints = 0; f < subNb; f++) { 30420cf1dd8SToby Isaac PetscQuadrature q; 30520cf1dd8SToby Isaac PetscInt qNp; 30620cf1dd8SToby Isaac 3079566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetFunctional(dualSubspace, f, &q)); 3089566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(q, NULL, NULL, &qNp, NULL, NULL)); 30920cf1dd8SToby Isaac numPoints += qNp; 31020cf1dd8SToby Isaac } 3119566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(subNb * origNb, &V)); 3129566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(numPoints * origDim, &allPoints, numPoints * origNc, &allWeights, numPoints * origNb * origNc, &B)); 31320cf1dd8SToby Isaac for (f = 0, offset = 0; f < subNb; f++) { 31420cf1dd8SToby Isaac PetscQuadrature q; 31520cf1dd8SToby Isaac PetscInt qNp, p; 31620cf1dd8SToby Isaac const PetscReal *qp; 31720cf1dd8SToby Isaac const PetscReal *qw; 31820cf1dd8SToby Isaac 3199566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetFunctional(dualSubspace, f, &q)); 3209566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(q, NULL, NULL, &qNp, &qp, &qw)); 32120cf1dd8SToby Isaac for (p = 0; p < qNp; p++, offset++) { 32220cf1dd8SToby Isaac if (x) { 32320cf1dd8SToby Isaac for (i = 0; i < origDim; i++) allPoints[origDim * offset + i] = x[i]; 32420cf1dd8SToby Isaac } else { 32520cf1dd8SToby Isaac for (i = 0; i < origDim; i++) allPoints[origDim * offset + i] = 0.0; 32620cf1dd8SToby Isaac } 32720cf1dd8SToby Isaac if (Jx) { 32820cf1dd8SToby Isaac for (i = 0; i < origDim; i++) { 329ad540459SPierre Jolivet for (j = 0; j < subDim; j++) allPoints[origDim * offset + i] += Jx[i * subDim + j] * qp[j]; 33020cf1dd8SToby Isaac } 33120cf1dd8SToby Isaac } else { 33220cf1dd8SToby Isaac for (i = 0; i < PetscMin(subDim, origDim); i++) allPoints[origDim * offset + i] += qp[i]; 33320cf1dd8SToby Isaac } 33420cf1dd8SToby Isaac for (i = 0; i < origNc; i++) allWeights[origNc * offset + i] = 0.0; 33520cf1dd8SToby Isaac if (Ju) { 33620cf1dd8SToby Isaac for (i = 0; i < origNc; i++) { 337ad540459SPierre Jolivet for (j = 0; j < subNc; j++) allWeights[offset * origNc + i] += qw[j] * Ju[j * origNc + i]; 33820cf1dd8SToby Isaac } 33920cf1dd8SToby Isaac } else { 34020cf1dd8SToby Isaac for (i = 0; i < PetscMin(subNc, origNc); i++) allWeights[offset * origNc + i] += qw[i]; 34120cf1dd8SToby Isaac } 34220cf1dd8SToby Isaac } 34320cf1dd8SToby Isaac } 3449566063dSJacob Faibussowitsch PetscCall(PetscSpaceEvaluate(origSpace, numPoints, allPoints, B, NULL, NULL)); 34520cf1dd8SToby Isaac for (f = 0, offset = 0; f < subNb; f++) { 34620cf1dd8SToby Isaac PetscInt b, p, s, qNp; 34720cf1dd8SToby Isaac PetscQuadrature q; 34820cf1dd8SToby Isaac const PetscReal *qw; 34920cf1dd8SToby Isaac 3509566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetFunctional(dualSubspace, f, &q)); 3519566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(q, NULL, NULL, &qNp, NULL, &qw)); 35220cf1dd8SToby Isaac if (u) { 35320cf1dd8SToby Isaac for (b = 0; b < origNb; b++) { 354ad540459SPierre Jolivet for (s = 0; s < subNc; s++) V[f * origNb + b] += qw[s] * u[s]; 35520cf1dd8SToby Isaac } 35620cf1dd8SToby Isaac } else { 35720cf1dd8SToby Isaac for (b = 0; b < origNb; b++) V[f * origNb + b] = 0.0; 35820cf1dd8SToby Isaac } 35920cf1dd8SToby Isaac for (p = 0; p < qNp; p++, offset++) { 36020cf1dd8SToby Isaac for (b = 0; b < origNb; b++) { 361ad540459SPierre Jolivet for (s = 0; s < origNc; s++) V[f * origNb + b] += B[(offset * origNb + b) * origNc + s] * allWeights[offset * origNc + s]; 36220cf1dd8SToby Isaac } 36320cf1dd8SToby Isaac } 36420cf1dd8SToby Isaac } 36520cf1dd8SToby Isaac /* orthnormalize rows of V */ 36620cf1dd8SToby Isaac for (f = 0; f < subNb; f++) { 36720cf1dd8SToby Isaac PetscReal rho = 0.0, scal; 36820cf1dd8SToby Isaac 36920cf1dd8SToby Isaac for (i = 0; i < origNb; i++) rho += PetscSqr(V[f * origNb + i]); 37020cf1dd8SToby Isaac 37120cf1dd8SToby Isaac scal = 1. / PetscSqrtReal(rho); 37220cf1dd8SToby Isaac 37320cf1dd8SToby Isaac for (i = 0; i < origNb; i++) V[f * origNb + i] *= scal; 37420cf1dd8SToby Isaac for (j = f + 1; j < subNb; j++) { 37520cf1dd8SToby Isaac for (i = 0, scal = 0.; i < origNb; i++) scal += V[f * origNb + i] * V[j * origNb + i]; 37620cf1dd8SToby Isaac for (i = 0; i < origNb; i++) V[j * origNb + i] -= V[f * origNb + i] * scal; 37720cf1dd8SToby Isaac } 37820cf1dd8SToby Isaac } 3799566063dSJacob Faibussowitsch PetscCall(PetscFree3(allPoints, allWeights, B)); 38020cf1dd8SToby Isaac subsp->Q = V; 3813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38220cf1dd8SToby Isaac } 38320cf1dd8SToby Isaac 384d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpacePolynomialGetTensor_Subspace(PetscSpace sp, PetscBool *poly) 385d71ae5a4SJacob Faibussowitsch { 38620cf1dd8SToby Isaac PetscSpace_Subspace *subsp = (PetscSpace_Subspace *)sp->data; 38720cf1dd8SToby Isaac 38820cf1dd8SToby Isaac PetscFunctionBegin; 38920cf1dd8SToby Isaac *poly = PETSC_FALSE; 3909566063dSJacob Faibussowitsch PetscCall(PetscSpacePolynomialGetTensor(subsp->origSpace, poly)); 39120cf1dd8SToby Isaac if (*poly) { 39220cf1dd8SToby Isaac if (subsp->Jx) { 39320cf1dd8SToby Isaac PetscInt subDim, origDim, i, j; 39420cf1dd8SToby Isaac PetscInt maxnnz; 39520cf1dd8SToby Isaac 3969566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumVariables(subsp->origSpace, &origDim)); 3979566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumVariables(sp, &subDim)); 39820cf1dd8SToby Isaac maxnnz = 0; 39920cf1dd8SToby Isaac for (i = 0; i < origDim; i++) { 40020cf1dd8SToby Isaac PetscInt nnz = 0; 40120cf1dd8SToby Isaac 40220cf1dd8SToby Isaac for (j = 0; j < subDim; j++) nnz += (subsp->Jx[i * subDim + j] != 0.); 40320cf1dd8SToby Isaac maxnnz = PetscMax(maxnnz, nnz); 40420cf1dd8SToby Isaac } 40520cf1dd8SToby Isaac for (j = 0; j < subDim; j++) { 40620cf1dd8SToby Isaac PetscInt nnz = 0; 40720cf1dd8SToby Isaac 40820cf1dd8SToby Isaac for (i = 0; i < origDim; i++) nnz += (subsp->Jx[i * subDim + j] != 0.); 40920cf1dd8SToby Isaac maxnnz = PetscMax(maxnnz, nnz); 41020cf1dd8SToby Isaac } 41120cf1dd8SToby Isaac if (maxnnz > 1) *poly = PETSC_FALSE; 41220cf1dd8SToby Isaac } 41320cf1dd8SToby Isaac } 4143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 41520cf1dd8SToby Isaac } 41620cf1dd8SToby Isaac 417d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscSpaceInitialize_Subspace(PetscSpace sp) 418d71ae5a4SJacob Faibussowitsch { 41920cf1dd8SToby Isaac PetscFunctionBegin; 42020cf1dd8SToby Isaac sp->ops->setup = PetscSpaceSetUp_Subspace; 42120cf1dd8SToby Isaac sp->ops->view = PetscSpaceView_Subspace; 42220cf1dd8SToby Isaac sp->ops->destroy = PetscSpaceDestroy_Subspace; 42320cf1dd8SToby Isaac sp->ops->getdimension = PetscSpaceGetDimension_Subspace; 42420cf1dd8SToby Isaac sp->ops->evaluate = PetscSpaceEvaluate_Subspace; 4259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Subspace)); 4263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 42720cf1dd8SToby Isaac } 428b24fb147SBarry Smith /*@ 429b24fb147SBarry Smith PetscSpaceCreateSubspace - creates a subspace from a an `origSpace` and its dual `dualSubspace` 43020cf1dd8SToby Isaac 431b24fb147SBarry Smith Input Parameters: 432b24fb147SBarry Smith + origSpace - the original `PetscSpace` 433b24fb147SBarry Smith . dualSubspace - no idea 434b24fb147SBarry Smith . x - no idea 435b24fb147SBarry Smith . Jx - no idea 436b24fb147SBarry Smith . u - no idea 437b24fb147SBarry Smith . Ju - no idea 438b24fb147SBarry Smith - copymode - whether to copy, borrow, or own some of the input arrays I guess 439b24fb147SBarry Smith 440b24fb147SBarry Smith Output Parameter: 441b24fb147SBarry Smith . subspace - the subspace 442b24fb147SBarry Smith 443b24fb147SBarry Smith Level: advanced 444b24fb147SBarry Smith 445b24fb147SBarry Smith .seealso: `PetscSpace`, `PetscDualSpace`, `PetscCopyMode`, `PetscSpaceType` 446b24fb147SBarry Smith @*/ 447d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscSpaceCreateSubspace(PetscSpace origSpace, PetscDualSpace dualSubspace, PetscReal *x, PetscReal *Jx, PetscReal *u, PetscReal *Ju, PetscCopyMode copymode, PetscSpace *subspace) 448d71ae5a4SJacob Faibussowitsch { 44920cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 45020cf1dd8SToby Isaac PetscInt origDim, subDim, origNc, subNc, subNb; 45120cf1dd8SToby Isaac PetscInt order; 45220cf1dd8SToby Isaac DM dm; 45320cf1dd8SToby Isaac 45420cf1dd8SToby Isaac PetscFunctionBegin; 45520cf1dd8SToby Isaac PetscValidHeaderSpecific(origSpace, PETSCSPACE_CLASSID, 1); 45620cf1dd8SToby Isaac PetscValidHeaderSpecific(dualSubspace, PETSCDUALSPACE_CLASSID, 2); 4574f572ea9SToby Isaac if (x) PetscAssertPointer(x, 3); 4584f572ea9SToby Isaac if (Jx) PetscAssertPointer(Jx, 4); 4594f572ea9SToby Isaac if (u) PetscAssertPointer(u, 5); 4604f572ea9SToby Isaac if (Ju) PetscAssertPointer(Ju, 6); 4614f572ea9SToby Isaac PetscAssertPointer(subspace, 8); 4629566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumComponents(origSpace, &origNc)); 4639566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetNumVariables(origSpace, &origDim)); 4649566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDM(dualSubspace, &dm)); 4659566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &subDim)); 4669566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetDimension(dualSubspace, &subNb)); 4679566063dSJacob Faibussowitsch PetscCall(PetscDualSpaceGetNumComponents(dualSubspace, &subNc)); 4689566063dSJacob Faibussowitsch PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)origSpace), subspace)); 4699566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetType(*subspace, PETSCSPACESUBSPACE)); 4709566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumVariables(*subspace, subDim)); 4719566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetNumComponents(*subspace, subNc)); 4729566063dSJacob Faibussowitsch PetscCall(PetscSpaceGetDegree(origSpace, &order, NULL)); 4739566063dSJacob Faibussowitsch PetscCall(PetscSpaceSetDegree(*subspace, order, PETSC_DETERMINE)); 47420cf1dd8SToby Isaac subsp = (PetscSpace_Subspace *)(*subspace)->data; 47520cf1dd8SToby Isaac subsp->Nb = subNb; 47620cf1dd8SToby Isaac switch (copymode) { 47720cf1dd8SToby Isaac case PETSC_OWN_POINTER: 47820cf1dd8SToby Isaac if (x) subsp->x_alloc = x; 47920cf1dd8SToby Isaac if (Jx) subsp->Jx_alloc = Jx; 48020cf1dd8SToby Isaac if (u) subsp->u_alloc = u; 48120cf1dd8SToby Isaac if (Ju) subsp->Ju_alloc = Ju; 482f4d061e9SPierre Jolivet /* fall through */ 48320cf1dd8SToby Isaac case PETSC_USE_POINTER: 48420cf1dd8SToby Isaac if (x) subsp->x = x; 48520cf1dd8SToby Isaac if (Jx) subsp->Jx = Jx; 48620cf1dd8SToby Isaac if (u) subsp->u = u; 48720cf1dd8SToby Isaac if (Ju) subsp->Ju = Ju; 48820cf1dd8SToby Isaac break; 48920cf1dd8SToby Isaac case PETSC_COPY_VALUES: 49020cf1dd8SToby Isaac if (x) { 4919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(origDim, &subsp->x_alloc)); 4929566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(subsp->x_alloc, x, origDim)); 49320cf1dd8SToby Isaac subsp->x = subsp->x_alloc; 49420cf1dd8SToby Isaac } 49520cf1dd8SToby Isaac if (Jx) { 4969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(origDim * subDim, &subsp->Jx_alloc)); 4979566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(subsp->Jx_alloc, Jx, origDim * subDim)); 49820cf1dd8SToby Isaac subsp->Jx = subsp->Jx_alloc; 49920cf1dd8SToby Isaac } 50020cf1dd8SToby Isaac if (u) { 5019566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(subNc, &subsp->u_alloc)); 5029566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(subsp->u_alloc, u, subNc)); 50320cf1dd8SToby Isaac subsp->u = subsp->u_alloc; 50420cf1dd8SToby Isaac } 50520cf1dd8SToby Isaac if (Ju) { 5069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(origNc * subNc, &subsp->Ju_alloc)); 5079566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(subsp->Ju_alloc, Ju, origNc * subNc)); 50820cf1dd8SToby Isaac subsp->Ju = subsp->Ju_alloc; 50920cf1dd8SToby Isaac } 51020cf1dd8SToby Isaac break; 511d71ae5a4SJacob Faibussowitsch default: 512d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)origSpace), PETSC_ERR_ARG_OUTOFRANGE, "Unknown copy mode"); 51320cf1dd8SToby Isaac } 5149566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)origSpace)); 51520cf1dd8SToby Isaac subsp->origSpace = origSpace; 5169566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)dualSubspace)); 51720cf1dd8SToby Isaac subsp->dualSubspace = dualSubspace; 5189566063dSJacob Faibussowitsch PetscCall(PetscSpaceInitialize_Subspace(*subspace)); 5193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 52020cf1dd8SToby Isaac } 521