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 1820cf1dd8SToby Isaac static PetscErrorCode PetscSpaceDestroy_Subspace(PetscSpace sp) 1920cf1dd8SToby Isaac { 2020cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 2120cf1dd8SToby Isaac PetscErrorCode ierr; 2220cf1dd8SToby Isaac 2320cf1dd8SToby Isaac PetscFunctionBegin; 2420cf1dd8SToby Isaac subsp = (PetscSpace_Subspace *) sp->data; 2520cf1dd8SToby Isaac subsp->x = NULL; 2620cf1dd8SToby Isaac ierr = PetscFree(subsp->x_alloc);CHKERRQ(ierr); 2720cf1dd8SToby Isaac subsp->Jx = NULL; 2820cf1dd8SToby Isaac ierr = PetscFree(subsp->Jx_alloc);CHKERRQ(ierr); 2920cf1dd8SToby Isaac subsp->u = NULL; 3020cf1dd8SToby Isaac ierr = PetscFree(subsp->u_alloc);CHKERRQ(ierr); 3120cf1dd8SToby Isaac subsp->Ju = NULL; 3220cf1dd8SToby Isaac ierr = PetscFree(subsp->Ju_alloc);CHKERRQ(ierr); 3320cf1dd8SToby Isaac ierr = PetscFree(subsp->Q);CHKERRQ(ierr); 3420cf1dd8SToby Isaac ierr = PetscSpaceDestroy(&subsp->origSpace);CHKERRQ(ierr); 3520cf1dd8SToby Isaac ierr = PetscDualSpaceDestroy(&subsp->dualSubspace);CHKERRQ(ierr); 3620cf1dd8SToby Isaac ierr = PetscFree(subsp);CHKERRQ(ierr); 3720cf1dd8SToby Isaac sp->data = NULL; 3820cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", NULL);CHKERRQ(ierr); 3920cf1dd8SToby Isaac PetscFunctionReturn(0); 4020cf1dd8SToby Isaac } 4120cf1dd8SToby Isaac 4220cf1dd8SToby Isaac static PetscErrorCode PetscSpaceView_Subspace(PetscSpace sp, PetscViewer viewer) 4320cf1dd8SToby Isaac { 4420cf1dd8SToby Isaac PetscBool iascii; 4520cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 4620cf1dd8SToby Isaac PetscErrorCode ierr; 4720cf1dd8SToby Isaac 4820cf1dd8SToby Isaac PetscFunctionBegin; 4920cf1dd8SToby Isaac subsp = (PetscSpace_Subspace *) sp->data; 5020cf1dd8SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 5120cf1dd8SToby Isaac if (iascii) { 5220cf1dd8SToby Isaac PetscInt origDim, subDim, origNc, subNc, o, s; 5320cf1dd8SToby Isaac 5420cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(subsp->origSpace,&origDim);CHKERRQ(ierr); 5520cf1dd8SToby Isaac ierr = PetscSpaceGetNumComponents(subsp->origSpace,&origNc);CHKERRQ(ierr); 5620cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(sp,&subDim);CHKERRQ(ierr); 5720cf1dd8SToby Isaac ierr = PetscSpaceGetNumComponents(sp,&subNc);CHKERRQ(ierr); 5820cf1dd8SToby Isaac if (subsp->x) { 5920cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"Subspace-to-space domain shift:\n\n");CHKERRQ(ierr); 6020cf1dd8SToby Isaac for (o = 0; o < origDim; o++) { 6120cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer," %g\n", (double)subsp->x[o]);CHKERRQ(ierr); 6220cf1dd8SToby Isaac } 6320cf1dd8SToby Isaac } 6420cf1dd8SToby Isaac if (subsp->Jx) { 6520cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"Subspace-to-space domain transform:\n\n");CHKERRQ(ierr); 6620cf1dd8SToby Isaac for (o = 0; o < origDim; o++) { 6720cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer," %g", (double)subsp->Jx[o * subDim + 0]);CHKERRQ(ierr); 6820cf1dd8SToby Isaac ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 6920cf1dd8SToby Isaac for (s = 1; s < subDim; s++) { 7020cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer," %g", (double)subsp->Jx[o * subDim + s]);CHKERRQ(ierr); 7120cf1dd8SToby Isaac } 7220cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 7320cf1dd8SToby Isaac ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 7420cf1dd8SToby Isaac } 7520cf1dd8SToby Isaac } 7620cf1dd8SToby Isaac if (subsp->u) { 7720cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"Space-to-subspace range shift:\n\n");CHKERRQ(ierr); 7820cf1dd8SToby Isaac for (o = 0; o < origNc; o++) { 7920cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer," %d\n", subsp->u[o]);CHKERRQ(ierr); 8020cf1dd8SToby Isaac } 8120cf1dd8SToby Isaac } 8220cf1dd8SToby Isaac if (subsp->Ju) { 8320cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"Space-to-subsace domain transform:\n");CHKERRQ(ierr); 8420cf1dd8SToby Isaac for (o = 0; o < origNc; o++) { 8520cf1dd8SToby Isaac ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 8620cf1dd8SToby Isaac for (s = 0; s < subNc; s++) { 8720cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer," %d", subsp->Ju[o * subNc + s]);CHKERRQ(ierr); 8820cf1dd8SToby Isaac } 8920cf1dd8SToby Isaac ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 9020cf1dd8SToby Isaac } 9120cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 9220cf1dd8SToby Isaac } 9320cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"Original space:\n");CHKERRQ(ierr); 9420cf1dd8SToby Isaac } 9520cf1dd8SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 9620cf1dd8SToby Isaac ierr = PetscSpaceView(subsp->origSpace,viewer);CHKERRQ(ierr); 9720cf1dd8SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 9820cf1dd8SToby Isaac PetscFunctionReturn(0); 9920cf1dd8SToby Isaac } 10020cf1dd8SToby Isaac 10120cf1dd8SToby Isaac static PetscErrorCode PetscSpaceEvaluate_Subspace(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) 10220cf1dd8SToby Isaac { 10320cf1dd8SToby Isaac PetscSpace_Subspace *subsp = (PetscSpace_Subspace *) sp->data; 10420cf1dd8SToby Isaac PetscSpace origsp; 10520cf1dd8SToby Isaac PetscInt origDim, subDim, origNc, subNc, subNb, origNb, i, j, k, l, m, n, o; 10620cf1dd8SToby Isaac PetscReal *inpoints, *inB = NULL, *inD = NULL, *inH = NULL; 10720cf1dd8SToby Isaac PetscErrorCode ierr; 10820cf1dd8SToby Isaac 10920cf1dd8SToby Isaac PetscFunctionBegin; 11020cf1dd8SToby Isaac origsp = subsp->origSpace; 11120cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(sp,&subDim);CHKERRQ(ierr); 11220cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(origsp,&origDim);CHKERRQ(ierr); 11320cf1dd8SToby Isaac ierr = PetscSpaceGetNumComponents(sp,&subNc);CHKERRQ(ierr); 11420cf1dd8SToby Isaac ierr = PetscSpaceGetNumComponents(origsp,&origNc);CHKERRQ(ierr); 11520cf1dd8SToby Isaac ierr = PetscSpaceGetDimension(sp,&subNb);CHKERRQ(ierr); 11620cf1dd8SToby Isaac ierr = PetscSpaceGetDimension(origsp,&origNb);CHKERRQ(ierr); 11720cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,npoints*origDim,MPIU_REAL,&inpoints);CHKERRQ(ierr); 11820cf1dd8SToby Isaac for (i = 0; i < npoints; i++) { 11920cf1dd8SToby Isaac if (subsp->x) { 12020cf1dd8SToby Isaac for (j = 0; j < origDim; j++) inpoints[i * origDim + j] = subsp->x[j]; 12120cf1dd8SToby Isaac } else { 12220cf1dd8SToby Isaac for (j = 0; j < origDim; j++) inpoints[i * origDim + j] = 0.0; 12320cf1dd8SToby Isaac } 12420cf1dd8SToby Isaac if (subsp->Jx) { 12520cf1dd8SToby Isaac for (j = 0; j < origDim; j++) { 12620cf1dd8SToby Isaac for (k = 0; k < subDim; k++) { 12720cf1dd8SToby Isaac inpoints[i * origDim + j] += subsp->Jx[j * subDim + k] * points[i * subDim + k]; 12820cf1dd8SToby Isaac } 12920cf1dd8SToby Isaac } 13020cf1dd8SToby Isaac } else { 13120cf1dd8SToby Isaac for (j = 0; j < PetscMin(subDim, origDim); j++) { 13220cf1dd8SToby Isaac inpoints[i * origDim + j] += points[i * subDim + j]; 13320cf1dd8SToby Isaac } 13420cf1dd8SToby Isaac } 13520cf1dd8SToby Isaac } 13620cf1dd8SToby Isaac if (B) { 13720cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,npoints*origNb*origNc,MPIU_REAL,&inB);CHKERRQ(ierr); 13820cf1dd8SToby Isaac } 13920cf1dd8SToby Isaac if (D) { 14020cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,npoints*origNb*origNc*origDim,MPIU_REAL,&inD);CHKERRQ(ierr); 14120cf1dd8SToby Isaac } 14220cf1dd8SToby Isaac if (H) { 14320cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,npoints*origNb*origNc*origDim*origDim,MPIU_REAL,&inH);CHKERRQ(ierr); 14420cf1dd8SToby Isaac } 14520cf1dd8SToby Isaac ierr = PetscSpaceEvaluate(origsp,npoints,inpoints,inB,inD,inH);CHKERRQ(ierr); 14620cf1dd8SToby Isaac if (H) { 14720cf1dd8SToby Isaac PetscReal *phi, *psi; 14820cf1dd8SToby Isaac 14920cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,origNc*origDim*origDim,MPIU_REAL,&phi);CHKERRQ(ierr); 15020cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,origNc*subDim*subDim,MPIU_REAL,&psi);CHKERRQ(ierr); 15120cf1dd8SToby Isaac for (i = 0; i < npoints * subNb * subNc * subDim; i++) D[i] = 0.0; 15220cf1dd8SToby Isaac for (i = 0; i < subNb; i++) { 15320cf1dd8SToby Isaac const PetscReal *subq = &subsp->Q[i * origNb]; 15420cf1dd8SToby Isaac 15520cf1dd8SToby Isaac for (j = 0; j < npoints; j++) { 15620cf1dd8SToby Isaac for (k = 0; k < origNc * origDim; k++) phi[k] = 0.; 15720cf1dd8SToby Isaac for (k = 0; k < origNc * subDim; k++) psi[k] = 0.; 15820cf1dd8SToby Isaac for (k = 0; k < origNb; k++) { 15920cf1dd8SToby Isaac for (l = 0; l < origNc * origDim * origDim; l++) { 16020cf1dd8SToby Isaac phi[l] += inH[(j * origNb + k) * origNc * origDim * origDim + l] * subq[k]; 16120cf1dd8SToby Isaac } 16220cf1dd8SToby Isaac } 16320cf1dd8SToby Isaac if (subsp->Jx) { 16420cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 16520cf1dd8SToby Isaac for (l = 0; l < subDim; l++) { 16620cf1dd8SToby Isaac for (m = 0; m < origDim; m++) { 16720cf1dd8SToby Isaac for (n = 0; n < subDim; n++) { 16820cf1dd8SToby Isaac for (o = 0; o < origDim; o++) { 16920cf1dd8SToby Isaac psi[(k * subDim + l) * subDim + n] += subsp->Jx[m * subDim + l] * subsp->Jx[o * subDim + n] * phi[(k * origDim + m) * origDim + o]; 17020cf1dd8SToby Isaac } 17120cf1dd8SToby Isaac } 17220cf1dd8SToby Isaac } 17320cf1dd8SToby Isaac } 17420cf1dd8SToby Isaac } 17520cf1dd8SToby Isaac } else { 17620cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 17720cf1dd8SToby Isaac for (l = 0; l < PetscMin(subDim, origDim); l++) { 17820cf1dd8SToby Isaac for (m = 0; m < PetscMin(subDim, origDim); m++) { 17920cf1dd8SToby Isaac psi[(k * subDim + l) * subDim + m] += phi[(k * origDim + l) * origDim + m]; 18020cf1dd8SToby Isaac } 18120cf1dd8SToby Isaac } 18220cf1dd8SToby Isaac } 18320cf1dd8SToby Isaac } 18420cf1dd8SToby Isaac if (subsp->Ju) { 18520cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 18620cf1dd8SToby Isaac for (l = 0; l < origNc; l++) { 18720cf1dd8SToby Isaac for (m = 0; m < subDim * subDim; m++) { 18820cf1dd8SToby Isaac H[((j * subNb + i) * subNc + k) * subDim * subDim + m] += subsp->Ju[k * origNc + l] * psi[l * subDim * subDim + m]; 18920cf1dd8SToby Isaac } 19020cf1dd8SToby Isaac } 19120cf1dd8SToby Isaac } 19220cf1dd8SToby Isaac } 19320cf1dd8SToby Isaac else { 19420cf1dd8SToby Isaac for (k = 0; k < PetscMin(subNc, origNc); k++) { 19520cf1dd8SToby Isaac for (l = 0; l < subDim * subDim; l++) { 19620cf1dd8SToby Isaac H[((j * subNb + i) * subNc + k) * subDim * subDim + l] += psi[k * subDim * subDim + l]; 19720cf1dd8SToby Isaac } 19820cf1dd8SToby Isaac } 19920cf1dd8SToby Isaac } 20020cf1dd8SToby Isaac } 20120cf1dd8SToby Isaac } 20220cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,subNc*origDim,MPIU_REAL,&psi);CHKERRQ(ierr); 20320cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,origNc*origDim,MPIU_REAL,&phi);CHKERRQ(ierr); 20420cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,npoints*origNb*origNc*origDim,MPIU_REAL,&inH);CHKERRQ(ierr); 20520cf1dd8SToby Isaac } 20620cf1dd8SToby Isaac if (D) { 20720cf1dd8SToby Isaac PetscReal *phi, *psi; 20820cf1dd8SToby Isaac 20920cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,origNc*origDim,MPIU_REAL,&phi);CHKERRQ(ierr); 21020cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,origNc*subDim,MPIU_REAL,&psi);CHKERRQ(ierr); 21120cf1dd8SToby Isaac for (i = 0; i < npoints * subNb * subNc * subDim; i++) D[i] = 0.0; 21220cf1dd8SToby Isaac for (i = 0; i < subNb; i++) { 21320cf1dd8SToby Isaac const PetscReal *subq = &subsp->Q[i * origNb]; 21420cf1dd8SToby Isaac 21520cf1dd8SToby Isaac for (j = 0; j < npoints; j++) { 21620cf1dd8SToby Isaac for (k = 0; k < origNc * origDim; k++) phi[k] = 0.; 21720cf1dd8SToby Isaac for (k = 0; k < origNc * subDim; k++) psi[k] = 0.; 21820cf1dd8SToby Isaac for (k = 0; k < origNb; k++) { 21920cf1dd8SToby Isaac for (l = 0; l < origNc * origDim; l++) { 22020cf1dd8SToby Isaac phi[l] += inD[(j * origNb + k) * origNc * origDim + l] * subq[k]; 22120cf1dd8SToby Isaac } 22220cf1dd8SToby Isaac } 22320cf1dd8SToby Isaac if (subsp->Jx) { 22420cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 22520cf1dd8SToby Isaac for (l = 0; l < subDim; l++) { 22620cf1dd8SToby Isaac for (m = 0; m < origDim; m++) { 22720cf1dd8SToby Isaac psi[k * subDim + l] += subsp->Jx[m * subDim + l] * phi[k * origDim + m]; 22820cf1dd8SToby Isaac } 22920cf1dd8SToby Isaac } 23020cf1dd8SToby Isaac } 23120cf1dd8SToby Isaac } else { 23220cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 23320cf1dd8SToby Isaac for (l = 0; l < PetscMin(subDim, origDim); l++) { 23420cf1dd8SToby Isaac psi[k * subDim + l] += phi[k * origDim + l]; 23520cf1dd8SToby Isaac } 23620cf1dd8SToby Isaac } 23720cf1dd8SToby Isaac } 23820cf1dd8SToby Isaac if (subsp->Ju) { 23920cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 24020cf1dd8SToby Isaac for (l = 0; l < origNc; l++) { 24120cf1dd8SToby Isaac for (m = 0; m < subDim; m++) { 24220cf1dd8SToby Isaac D[((j * subNb + i) * subNc + k) * subDim + m] += subsp->Ju[k * origNc + l] * psi[l * subDim + m]; 24320cf1dd8SToby Isaac } 24420cf1dd8SToby Isaac } 24520cf1dd8SToby Isaac } 24620cf1dd8SToby Isaac } 24720cf1dd8SToby Isaac else { 24820cf1dd8SToby Isaac for (k = 0; k < PetscMin(subNc, origNc); k++) { 24920cf1dd8SToby Isaac for (l = 0; l < subDim; l++) { 25020cf1dd8SToby Isaac D[((j * subNb + i) * subNc + k) * subDim + l] += psi[k * subDim + l]; 25120cf1dd8SToby Isaac } 25220cf1dd8SToby Isaac } 25320cf1dd8SToby Isaac } 25420cf1dd8SToby Isaac } 25520cf1dd8SToby Isaac } 25620cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,subNc*origDim,MPIU_REAL,&psi);CHKERRQ(ierr); 25720cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,origNc*origDim,MPIU_REAL,&phi);CHKERRQ(ierr); 25820cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,npoints*origNb*origNc*origDim,MPIU_REAL,&inD);CHKERRQ(ierr); 25920cf1dd8SToby Isaac } 26020cf1dd8SToby Isaac if (B) { 26120cf1dd8SToby Isaac PetscReal *phi; 26220cf1dd8SToby Isaac 26320cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,origNc,MPIU_REAL,&phi);CHKERRQ(ierr); 26420cf1dd8SToby Isaac if (subsp->u) { 26520cf1dd8SToby Isaac for (i = 0; i < npoints * subNb; i++) { 26620cf1dd8SToby Isaac for (j = 0; j < subNc; j++) B[i * subNc + j] = subsp->u[j]; 26720cf1dd8SToby Isaac } 26820cf1dd8SToby Isaac } else { 26920cf1dd8SToby Isaac for (i = 0; i < npoints * subNb * subNc; i++) B[i] = 0.0; 27020cf1dd8SToby Isaac } 27120cf1dd8SToby Isaac for (i = 0; i < subNb; i++) { 27220cf1dd8SToby Isaac const PetscReal *subq = &subsp->Q[i * origNb]; 27320cf1dd8SToby Isaac 27420cf1dd8SToby Isaac for (j = 0; j < npoints; j++) { 27520cf1dd8SToby Isaac for (k = 0; k < origNc; k++) phi[k] = 0.; 27620cf1dd8SToby Isaac for (k = 0; k < origNb; k++) { 27720cf1dd8SToby Isaac for (l = 0; l < origNc; l++) { 27820cf1dd8SToby Isaac phi[l] += inB[(j * origNb + k) * origNc + l] * subq[k]; 27920cf1dd8SToby Isaac } 28020cf1dd8SToby Isaac } 28120cf1dd8SToby Isaac if (subsp->Ju) { 28220cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 28320cf1dd8SToby Isaac for (l = 0; l < origNc; l++) { 28420cf1dd8SToby Isaac B[(j * subNb + i) * subNc + k] += subsp->Ju[k * origNc + l] * phi[l]; 28520cf1dd8SToby Isaac } 28620cf1dd8SToby Isaac } 28720cf1dd8SToby Isaac } 28820cf1dd8SToby Isaac else { 28920cf1dd8SToby Isaac for (k = 0; k < PetscMin(subNc, origNc); k++) { 29020cf1dd8SToby Isaac B[(j * subNb + i) * subNc + k] += phi[k]; 29120cf1dd8SToby Isaac } 29220cf1dd8SToby Isaac } 29320cf1dd8SToby Isaac } 29420cf1dd8SToby Isaac } 29520cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,origNc,MPIU_REAL,&phi);CHKERRQ(ierr); 29620cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,npoints*origNb*origNc,MPIU_REAL,&inB);CHKERRQ(ierr); 29720cf1dd8SToby Isaac } 29820cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,npoints*origDim,MPIU_REAL,&inpoints);CHKERRQ(ierr); 29920cf1dd8SToby Isaac PetscFunctionReturn(0); 30020cf1dd8SToby Isaac } 30120cf1dd8SToby Isaac 30220cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Subspace(PetscSpace sp) 30320cf1dd8SToby Isaac { 30420cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 30520cf1dd8SToby Isaac PetscErrorCode ierr; 306*362febeeSStefano Zampini 307*362febeeSStefano Zampini PetscFunctionBegin; 30820cf1dd8SToby Isaac ierr = PetscNewLog(sp,&subsp);CHKERRQ(ierr); 30920cf1dd8SToby Isaac sp->data = (void *) subsp; 31020cf1dd8SToby Isaac PetscFunctionReturn(0); 31120cf1dd8SToby Isaac } 31220cf1dd8SToby Isaac 31320cf1dd8SToby Isaac static PetscErrorCode PetscSpaceGetDimension_Subspace(PetscSpace sp, PetscInt *dim) 31420cf1dd8SToby Isaac { 31520cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 31620cf1dd8SToby Isaac 31720cf1dd8SToby Isaac PetscFunctionBegin; 31820cf1dd8SToby Isaac subsp = (PetscSpace_Subspace *) sp->data; 31920cf1dd8SToby Isaac *dim = subsp->Nb; 32020cf1dd8SToby Isaac PetscFunctionReturn(0); 32120cf1dd8SToby Isaac } 32220cf1dd8SToby Isaac 32320cf1dd8SToby Isaac static PetscErrorCode PetscSpaceSetUp_Subspace(PetscSpace sp) 32420cf1dd8SToby Isaac { 32520cf1dd8SToby Isaac const PetscReal *x; 32620cf1dd8SToby Isaac const PetscReal *Jx; 32720cf1dd8SToby Isaac const PetscReal *u; 32820cf1dd8SToby Isaac const PetscReal *Ju; 32920cf1dd8SToby Isaac PetscDualSpace dualSubspace; 33020cf1dd8SToby Isaac PetscSpace origSpace; 33120cf1dd8SToby Isaac PetscInt origDim, subDim, origNc, subNc, origNb, subNb, f, i, j, numPoints, offset; 33220cf1dd8SToby Isaac PetscReal *allPoints, *allWeights, *B, *V; 33320cf1dd8SToby Isaac DM dm; 33420cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 33520cf1dd8SToby Isaac PetscErrorCode ierr; 33620cf1dd8SToby Isaac 33720cf1dd8SToby Isaac PetscFunctionBegin; 33820cf1dd8SToby Isaac subsp = (PetscSpace_Subspace *) sp->data; 33920cf1dd8SToby Isaac x = subsp->x; 34020cf1dd8SToby Isaac Jx = subsp->Jx; 34120cf1dd8SToby Isaac u = subsp->u; 34220cf1dd8SToby Isaac Ju = subsp->Ju; 34320cf1dd8SToby Isaac origSpace = subsp->origSpace; 34420cf1dd8SToby Isaac dualSubspace = subsp->dualSubspace; 34520cf1dd8SToby Isaac ierr = PetscSpaceGetNumComponents(origSpace,&origNc);CHKERRQ(ierr); 34620cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(origSpace,&origDim);CHKERRQ(ierr); 34720cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(dualSubspace,&dm);CHKERRQ(ierr); 34820cf1dd8SToby Isaac ierr = DMGetDimension(dm,&subDim);CHKERRQ(ierr); 34920cf1dd8SToby Isaac ierr = PetscSpaceGetDimension(origSpace,&origNb);CHKERRQ(ierr); 35020cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(dualSubspace,&subNb);CHKERRQ(ierr); 35120cf1dd8SToby Isaac ierr = PetscDualSpaceGetNumComponents(dualSubspace,&subNc);CHKERRQ(ierr); 35220cf1dd8SToby Isaac 35320cf1dd8SToby Isaac for (f = 0, numPoints = 0; f < subNb; f++) { 35420cf1dd8SToby Isaac PetscQuadrature q; 35520cf1dd8SToby Isaac PetscInt qNp; 35620cf1dd8SToby Isaac 35720cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(dualSubspace,f,&q);CHKERRQ(ierr); 35820cf1dd8SToby Isaac ierr = PetscQuadratureGetData(q,NULL,NULL,&qNp,NULL,NULL);CHKERRQ(ierr); 35920cf1dd8SToby Isaac numPoints += qNp; 36020cf1dd8SToby Isaac } 36120cf1dd8SToby Isaac ierr = PetscMalloc1(subNb*origNb,&V);CHKERRQ(ierr); 36220cf1dd8SToby Isaac ierr = PetscMalloc3(numPoints*origDim,&allPoints,numPoints*origNc,&allWeights,numPoints*origNb*origNc,&B);CHKERRQ(ierr); 36320cf1dd8SToby Isaac for (f = 0, offset = 0; f < subNb; f++) { 36420cf1dd8SToby Isaac PetscQuadrature q; 36520cf1dd8SToby Isaac PetscInt qNp, p; 36620cf1dd8SToby Isaac const PetscReal *qp; 36720cf1dd8SToby Isaac const PetscReal *qw; 36820cf1dd8SToby Isaac 36920cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(dualSubspace,f,&q);CHKERRQ(ierr); 37020cf1dd8SToby Isaac ierr = PetscQuadratureGetData(q,NULL,NULL,&qNp,&qp,&qw);CHKERRQ(ierr); 37120cf1dd8SToby Isaac for (p = 0; p < qNp; p++, offset++) { 37220cf1dd8SToby Isaac if (x) { 37320cf1dd8SToby Isaac for (i = 0; i < origDim; i++) allPoints[origDim * offset + i] = x[i]; 37420cf1dd8SToby Isaac } else { 37520cf1dd8SToby Isaac for (i = 0; i < origDim; i++) allPoints[origDim * offset + i] = 0.0; 37620cf1dd8SToby Isaac } 37720cf1dd8SToby Isaac if (Jx) { 37820cf1dd8SToby Isaac for (i = 0; i < origDim; i++) { 37920cf1dd8SToby Isaac for (j = 0; j < subDim; j++) { 38020cf1dd8SToby Isaac allPoints[origDim * offset + i] += Jx[i * subDim + j] * qp[j]; 38120cf1dd8SToby Isaac } 38220cf1dd8SToby Isaac } 38320cf1dd8SToby Isaac } else { 38420cf1dd8SToby Isaac for (i = 0; i < PetscMin(subDim, origDim); i++) allPoints[origDim * offset + i] += qp[i]; 38520cf1dd8SToby Isaac } 38620cf1dd8SToby Isaac for (i = 0; i < origNc; i++) allWeights[origNc * offset + i] = 0.0; 38720cf1dd8SToby Isaac if (Ju) { 38820cf1dd8SToby Isaac for (i = 0; i < origNc; i++) { 38920cf1dd8SToby Isaac for (j = 0; j < subNc; j++) { 39020cf1dd8SToby Isaac allWeights[offset * origNc + i] += qw[j] * Ju[j * origNc + i]; 39120cf1dd8SToby Isaac } 39220cf1dd8SToby Isaac } 39320cf1dd8SToby Isaac } else { 39420cf1dd8SToby Isaac for (i = 0; i < PetscMin(subNc, origNc); i++) allWeights[offset * origNc + i] += qw[i]; 39520cf1dd8SToby Isaac } 39620cf1dd8SToby Isaac } 39720cf1dd8SToby Isaac } 39820cf1dd8SToby Isaac ierr = PetscSpaceEvaluate(origSpace,numPoints,allPoints,B,NULL,NULL);CHKERRQ(ierr); 39920cf1dd8SToby Isaac for (f = 0, offset = 0; f < subNb; f++) { 40020cf1dd8SToby Isaac PetscInt b, p, s, qNp; 40120cf1dd8SToby Isaac PetscQuadrature q; 40220cf1dd8SToby Isaac const PetscReal *qw; 40320cf1dd8SToby Isaac 40420cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(dualSubspace,f,&q);CHKERRQ(ierr); 40520cf1dd8SToby Isaac ierr = PetscQuadratureGetData(q,NULL,NULL,&qNp,NULL,&qw);CHKERRQ(ierr); 40620cf1dd8SToby Isaac if (u) { 40720cf1dd8SToby Isaac for (b = 0; b < origNb; b++) { 40820cf1dd8SToby Isaac for (s = 0; s < subNc; s++) { 40920cf1dd8SToby Isaac V[f * origNb + b] += qw[s] * u[s]; 41020cf1dd8SToby Isaac } 41120cf1dd8SToby Isaac } 41220cf1dd8SToby Isaac } else { 41320cf1dd8SToby Isaac for (b = 0; b < origNb; b++) V[f * origNb + b] = 0.0; 41420cf1dd8SToby Isaac } 41520cf1dd8SToby Isaac for (p = 0; p < qNp; p++, offset++) { 41620cf1dd8SToby Isaac for (b = 0; b < origNb; b++) { 41720cf1dd8SToby Isaac for (s = 0; s < origNc; s++) { 41820cf1dd8SToby Isaac V[f * origNb + b] += B[(offset * origNb + b) * origNc + s] * allWeights[offset * origNc + s]; 41920cf1dd8SToby Isaac } 42020cf1dd8SToby Isaac } 42120cf1dd8SToby Isaac } 42220cf1dd8SToby Isaac } 42320cf1dd8SToby Isaac /* orthnormalize rows of V */ 42420cf1dd8SToby Isaac for (f = 0; f < subNb; f++) { 42520cf1dd8SToby Isaac PetscReal rho = 0.0, scal; 42620cf1dd8SToby Isaac 42720cf1dd8SToby Isaac for (i = 0; i < origNb; i++) rho += PetscSqr(V[f * origNb + i]); 42820cf1dd8SToby Isaac 42920cf1dd8SToby Isaac scal = 1. / PetscSqrtReal(rho); 43020cf1dd8SToby Isaac 43120cf1dd8SToby Isaac for (i = 0; i < origNb; i++) V[f * origNb + i] *= scal; 43220cf1dd8SToby Isaac for (j = f + 1; j < subNb; j++) { 43320cf1dd8SToby Isaac for (i = 0, scal = 0.; i < origNb; i++) scal += V[f * origNb + i] * V[j * origNb + i]; 43420cf1dd8SToby Isaac for (i = 0; i < origNb; i++) V[j * origNb + i] -= V[f * origNb + i] * scal; 43520cf1dd8SToby Isaac } 43620cf1dd8SToby Isaac } 43720cf1dd8SToby Isaac ierr = PetscFree3(allPoints,allWeights,B);CHKERRQ(ierr); 43820cf1dd8SToby Isaac subsp->Q = V; 43920cf1dd8SToby Isaac PetscFunctionReturn(0); 44020cf1dd8SToby Isaac } 44120cf1dd8SToby Isaac 44220cf1dd8SToby Isaac static PetscErrorCode PetscSpacePolynomialGetTensor_Subspace(PetscSpace sp, PetscBool *poly) 44320cf1dd8SToby Isaac { 44420cf1dd8SToby Isaac PetscSpace_Subspace *subsp = (PetscSpace_Subspace *) sp->data; 44520cf1dd8SToby Isaac PetscErrorCode ierr; 44620cf1dd8SToby Isaac 44720cf1dd8SToby Isaac PetscFunctionBegin; 44820cf1dd8SToby Isaac *poly = PETSC_FALSE; 44920cf1dd8SToby Isaac ierr = PetscSpacePolynomialGetTensor(subsp->origSpace,poly);CHKERRQ(ierr); 45020cf1dd8SToby Isaac if (*poly) { 45120cf1dd8SToby Isaac if (subsp->Jx) { 45220cf1dd8SToby Isaac PetscInt subDim, origDim, i, j; 45320cf1dd8SToby Isaac PetscInt maxnnz; 45420cf1dd8SToby Isaac 45520cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(subsp->origSpace,&origDim);CHKERRQ(ierr); 45620cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(sp,&subDim);CHKERRQ(ierr); 45720cf1dd8SToby Isaac maxnnz = 0; 45820cf1dd8SToby Isaac for (i = 0; i < origDim; i++) { 45920cf1dd8SToby Isaac PetscInt nnz = 0; 46020cf1dd8SToby Isaac 46120cf1dd8SToby Isaac for (j = 0; j < subDim; j++) nnz += (subsp->Jx[i * subDim + j] != 0.); 46220cf1dd8SToby Isaac maxnnz = PetscMax(maxnnz,nnz); 46320cf1dd8SToby Isaac } 46420cf1dd8SToby Isaac for (j = 0; j < subDim; j++) { 46520cf1dd8SToby Isaac PetscInt nnz = 0; 46620cf1dd8SToby Isaac 46720cf1dd8SToby Isaac for (i = 0; i < origDim; i++) nnz += (subsp->Jx[i * subDim + j] != 0.); 46820cf1dd8SToby Isaac maxnnz = PetscMax(maxnnz,nnz); 46920cf1dd8SToby Isaac } 47020cf1dd8SToby Isaac if (maxnnz > 1) *poly = PETSC_FALSE; 47120cf1dd8SToby Isaac } 47220cf1dd8SToby Isaac } 47320cf1dd8SToby Isaac PetscFunctionReturn(0); 47420cf1dd8SToby Isaac } 47520cf1dd8SToby Isaac 47620cf1dd8SToby Isaac static PetscErrorCode PetscSpaceInitialize_Subspace(PetscSpace sp) 47720cf1dd8SToby Isaac { 47820cf1dd8SToby Isaac PetscErrorCode ierr; 47920cf1dd8SToby Isaac 48020cf1dd8SToby Isaac PetscFunctionBegin; 48120cf1dd8SToby Isaac sp->ops->setup = PetscSpaceSetUp_Subspace; 48220cf1dd8SToby Isaac sp->ops->view = PetscSpaceView_Subspace; 48320cf1dd8SToby Isaac sp->ops->destroy = PetscSpaceDestroy_Subspace; 48420cf1dd8SToby Isaac sp->ops->getdimension = PetscSpaceGetDimension_Subspace; 48520cf1dd8SToby Isaac sp->ops->evaluate = PetscSpaceEvaluate_Subspace; 48620cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Subspace);CHKERRQ(ierr); 48720cf1dd8SToby Isaac PetscFunctionReturn(0); 48820cf1dd8SToby Isaac } 48920cf1dd8SToby Isaac 49020cf1dd8SToby Isaac PetscErrorCode PetscSpaceCreateSubspace(PetscSpace origSpace, PetscDualSpace dualSubspace, PetscReal *x, PetscReal *Jx, PetscReal *u, PetscReal *Ju, PetscCopyMode copymode, PetscSpace *subspace) 49120cf1dd8SToby Isaac { 49220cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 49320cf1dd8SToby Isaac PetscInt origDim, subDim, origNc, subNc, subNb; 49420cf1dd8SToby Isaac PetscInt order; 49520cf1dd8SToby Isaac DM dm; 49620cf1dd8SToby Isaac PetscErrorCode ierr; 49720cf1dd8SToby Isaac 49820cf1dd8SToby Isaac PetscFunctionBegin; 49920cf1dd8SToby Isaac PetscValidHeaderSpecific(origSpace,PETSCSPACE_CLASSID,1); 50020cf1dd8SToby Isaac PetscValidHeaderSpecific(dualSubspace,PETSCDUALSPACE_CLASSID,2); 50120cf1dd8SToby Isaac if (x) PetscValidRealPointer(x,3); 50220cf1dd8SToby Isaac if (Jx) PetscValidRealPointer(Jx,4); 50320cf1dd8SToby Isaac if (u) PetscValidRealPointer(u,5); 50420cf1dd8SToby Isaac if (Ju) PetscValidRealPointer(Ju,6); 505064a246eSJacob Faibussowitsch PetscValidPointer(subspace,8); 50620cf1dd8SToby Isaac ierr = PetscSpaceGetNumComponents(origSpace,&origNc);CHKERRQ(ierr); 50720cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(origSpace,&origDim);CHKERRQ(ierr); 50820cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(dualSubspace,&dm);CHKERRQ(ierr); 50920cf1dd8SToby Isaac ierr = DMGetDimension(dm,&subDim);CHKERRQ(ierr); 51020cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(dualSubspace,&subNb);CHKERRQ(ierr); 51120cf1dd8SToby Isaac ierr = PetscDualSpaceGetNumComponents(dualSubspace,&subNc);CHKERRQ(ierr); 51220cf1dd8SToby Isaac ierr = PetscSpaceCreate(PetscObjectComm((PetscObject)origSpace),subspace);CHKERRQ(ierr); 51320cf1dd8SToby Isaac ierr = PetscSpaceSetType(*subspace,PETSCSPACESUBSPACE);CHKERRQ(ierr); 51420cf1dd8SToby Isaac ierr = PetscSpaceSetNumVariables(*subspace,subDim);CHKERRQ(ierr); 51520cf1dd8SToby Isaac ierr = PetscSpaceSetNumComponents(*subspace,subNc);CHKERRQ(ierr); 51620cf1dd8SToby Isaac ierr = PetscSpaceGetDegree(origSpace,&order,NULL);CHKERRQ(ierr); 517d39dd5f5SToby Isaac ierr = PetscSpaceSetDegree(*subspace,order,PETSC_DETERMINE);CHKERRQ(ierr); 51820cf1dd8SToby Isaac subsp = (PetscSpace_Subspace *) (*subspace)->data; 51920cf1dd8SToby Isaac subsp->Nb = subNb; 52020cf1dd8SToby Isaac switch (copymode) { 52120cf1dd8SToby Isaac case PETSC_OWN_POINTER: 52220cf1dd8SToby Isaac if (x) subsp->x_alloc = x; 52320cf1dd8SToby Isaac if (Jx) subsp->Jx_alloc = Jx; 52420cf1dd8SToby Isaac if (u) subsp->u_alloc = u; 52520cf1dd8SToby Isaac if (Ju) subsp->Ju_alloc = Ju; 52620cf1dd8SToby Isaac case PETSC_USE_POINTER: 52720cf1dd8SToby Isaac if (x) subsp->x = x; 52820cf1dd8SToby Isaac if (Jx) subsp->Jx = Jx; 52920cf1dd8SToby Isaac if (u) subsp->u = u; 53020cf1dd8SToby Isaac if (Ju) subsp->Ju = Ju; 53120cf1dd8SToby Isaac break; 53220cf1dd8SToby Isaac case PETSC_COPY_VALUES: 53320cf1dd8SToby Isaac if (x) { 53420cf1dd8SToby Isaac ierr = PetscMalloc1(origDim,&subsp->x_alloc);CHKERRQ(ierr); 535580bdb30SBarry Smith ierr = PetscArraycpy(subsp->x_alloc,x,origDim);CHKERRQ(ierr); 53620cf1dd8SToby Isaac subsp->x = subsp->x_alloc; 53720cf1dd8SToby Isaac } 53820cf1dd8SToby Isaac if (Jx) { 53920cf1dd8SToby Isaac ierr = PetscMalloc1(origDim * subDim,&subsp->Jx_alloc);CHKERRQ(ierr); 540580bdb30SBarry Smith ierr = PetscArraycpy(subsp->Jx_alloc,Jx,origDim * subDim);CHKERRQ(ierr); 54120cf1dd8SToby Isaac subsp->Jx = subsp->Jx_alloc; 54220cf1dd8SToby Isaac } 54320cf1dd8SToby Isaac if (u) { 54420cf1dd8SToby Isaac ierr = PetscMalloc1(subNc,&subsp->u_alloc);CHKERRQ(ierr); 545580bdb30SBarry Smith ierr = PetscArraycpy(subsp->u_alloc,u,subNc);CHKERRQ(ierr); 54620cf1dd8SToby Isaac subsp->u = subsp->u_alloc; 54720cf1dd8SToby Isaac } 54820cf1dd8SToby Isaac if (Ju) { 54920cf1dd8SToby Isaac ierr = PetscMalloc1(origNc * subNc,&subsp->Ju_alloc);CHKERRQ(ierr); 550580bdb30SBarry Smith ierr = PetscArraycpy(subsp->Ju_alloc,Ju,origNc * subNc);CHKERRQ(ierr); 55120cf1dd8SToby Isaac subsp->Ju = subsp->Ju_alloc; 55220cf1dd8SToby Isaac } 55320cf1dd8SToby Isaac break; 55420cf1dd8SToby Isaac default: 55520cf1dd8SToby Isaac SETERRQ(PetscObjectComm((PetscObject)origSpace),PETSC_ERR_ARG_OUTOFRANGE,"Unknown copy mode"); 55620cf1dd8SToby Isaac } 55720cf1dd8SToby Isaac ierr = PetscObjectReference((PetscObject)origSpace);CHKERRQ(ierr); 55820cf1dd8SToby Isaac subsp->origSpace = origSpace; 55920cf1dd8SToby Isaac ierr = PetscObjectReference((PetscObject)dualSubspace);CHKERRQ(ierr); 56020cf1dd8SToby Isaac subsp->dualSubspace = dualSubspace; 56120cf1dd8SToby Isaac ierr = PetscSpaceInitialize_Subspace(*subspace);CHKERRQ(ierr); 56220cf1dd8SToby Isaac PetscFunctionReturn(0); 56320cf1dd8SToby Isaac } 56420cf1dd8SToby Isaac 565