1*20cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 2*20cf1dd8SToby Isaac 3*20cf1dd8SToby Isaac typedef struct { 4*20cf1dd8SToby Isaac PetscDualSpace dualSubspace; 5*20cf1dd8SToby Isaac PetscSpace origSpace; 6*20cf1dd8SToby Isaac PetscReal *x; 7*20cf1dd8SToby Isaac PetscReal *x_alloc; 8*20cf1dd8SToby Isaac PetscReal *Jx; 9*20cf1dd8SToby Isaac PetscReal *Jx_alloc; 10*20cf1dd8SToby Isaac PetscReal *u; 11*20cf1dd8SToby Isaac PetscReal *u_alloc; 12*20cf1dd8SToby Isaac PetscReal *Ju; 13*20cf1dd8SToby Isaac PetscReal *Ju_alloc; 14*20cf1dd8SToby Isaac PetscReal *Q; 15*20cf1dd8SToby Isaac PetscInt Nb; 16*20cf1dd8SToby Isaac } PetscSpace_Subspace; 17*20cf1dd8SToby Isaac 18*20cf1dd8SToby Isaac static PetscErrorCode PetscSpaceDestroy_Subspace(PetscSpace sp) 19*20cf1dd8SToby Isaac { 20*20cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 21*20cf1dd8SToby Isaac PetscErrorCode ierr; 22*20cf1dd8SToby Isaac 23*20cf1dd8SToby Isaac PetscFunctionBegin; 24*20cf1dd8SToby Isaac subsp = (PetscSpace_Subspace *) sp->data; 25*20cf1dd8SToby Isaac subsp->x = NULL; 26*20cf1dd8SToby Isaac ierr = PetscFree(subsp->x_alloc);CHKERRQ(ierr); 27*20cf1dd8SToby Isaac subsp->Jx = NULL; 28*20cf1dd8SToby Isaac ierr = PetscFree(subsp->Jx_alloc);CHKERRQ(ierr); 29*20cf1dd8SToby Isaac subsp->u = NULL; 30*20cf1dd8SToby Isaac ierr = PetscFree(subsp->u_alloc);CHKERRQ(ierr); 31*20cf1dd8SToby Isaac subsp->Ju = NULL; 32*20cf1dd8SToby Isaac ierr = PetscFree(subsp->Ju_alloc);CHKERRQ(ierr); 33*20cf1dd8SToby Isaac ierr = PetscFree(subsp->Q);CHKERRQ(ierr); 34*20cf1dd8SToby Isaac ierr = PetscSpaceDestroy(&subsp->origSpace);CHKERRQ(ierr); 35*20cf1dd8SToby Isaac ierr = PetscDualSpaceDestroy(&subsp->dualSubspace);CHKERRQ(ierr); 36*20cf1dd8SToby Isaac ierr = PetscFree(subsp);CHKERRQ(ierr); 37*20cf1dd8SToby Isaac sp->data = NULL; 38*20cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", NULL);CHKERRQ(ierr); 39*20cf1dd8SToby Isaac PetscFunctionReturn(0); 40*20cf1dd8SToby Isaac } 41*20cf1dd8SToby Isaac 42*20cf1dd8SToby Isaac static PetscErrorCode PetscSpaceView_Subspace(PetscSpace sp, PetscViewer viewer) 43*20cf1dd8SToby Isaac { 44*20cf1dd8SToby Isaac PetscBool iascii; 45*20cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 46*20cf1dd8SToby Isaac PetscErrorCode ierr; 47*20cf1dd8SToby Isaac 48*20cf1dd8SToby Isaac PetscFunctionBegin; 49*20cf1dd8SToby Isaac subsp = (PetscSpace_Subspace *) sp->data; 50*20cf1dd8SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 51*20cf1dd8SToby Isaac if (iascii) { 52*20cf1dd8SToby Isaac PetscInt origDim, subDim, origNc, subNc, o, s; 53*20cf1dd8SToby Isaac 54*20cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(subsp->origSpace,&origDim);CHKERRQ(ierr); 55*20cf1dd8SToby Isaac ierr = PetscSpaceGetNumComponents(subsp->origSpace,&origNc);CHKERRQ(ierr); 56*20cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(sp,&subDim);CHKERRQ(ierr); 57*20cf1dd8SToby Isaac ierr = PetscSpaceGetNumComponents(sp,&subNc);CHKERRQ(ierr); 58*20cf1dd8SToby Isaac if (subsp->x) { 59*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"Subspace-to-space domain shift:\n\n");CHKERRQ(ierr); 60*20cf1dd8SToby Isaac for (o = 0; o < origDim; o++) { 61*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer," %g\n", (double)subsp->x[o]);CHKERRQ(ierr); 62*20cf1dd8SToby Isaac } 63*20cf1dd8SToby Isaac } 64*20cf1dd8SToby Isaac if (subsp->Jx) { 65*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"Subspace-to-space domain transform:\n\n");CHKERRQ(ierr); 66*20cf1dd8SToby Isaac for (o = 0; o < origDim; o++) { 67*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer," %g", (double)subsp->Jx[o * subDim + 0]);CHKERRQ(ierr); 68*20cf1dd8SToby Isaac ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 69*20cf1dd8SToby Isaac for (s = 1; s < subDim; s++) { 70*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer," %g", (double)subsp->Jx[o * subDim + s]);CHKERRQ(ierr); 71*20cf1dd8SToby Isaac } 72*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 73*20cf1dd8SToby Isaac ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 74*20cf1dd8SToby Isaac } 75*20cf1dd8SToby Isaac } 76*20cf1dd8SToby Isaac if (subsp->u) { 77*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"Space-to-subspace range shift:\n\n");CHKERRQ(ierr); 78*20cf1dd8SToby Isaac for (o = 0; o < origNc; o++) { 79*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer," %d\n", subsp->u[o]);CHKERRQ(ierr); 80*20cf1dd8SToby Isaac } 81*20cf1dd8SToby Isaac } 82*20cf1dd8SToby Isaac if (subsp->Ju) { 83*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"Space-to-subsace domain transform:\n");CHKERRQ(ierr); 84*20cf1dd8SToby Isaac for (o = 0; o < origNc; o++) { 85*20cf1dd8SToby Isaac ierr = PetscViewerASCIIUseTabs(viewer,PETSC_FALSE);CHKERRQ(ierr); 86*20cf1dd8SToby Isaac for (s = 0; s < subNc; s++) { 87*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer," %d", subsp->Ju[o * subNc + s]);CHKERRQ(ierr); 88*20cf1dd8SToby Isaac } 89*20cf1dd8SToby Isaac ierr = PetscViewerASCIIUseTabs(viewer,PETSC_TRUE);CHKERRQ(ierr); 90*20cf1dd8SToby Isaac } 91*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"\n");CHKERRQ(ierr); 92*20cf1dd8SToby Isaac } 93*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPrintf(viewer,"Original space:\n");CHKERRQ(ierr); 94*20cf1dd8SToby Isaac } 95*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 96*20cf1dd8SToby Isaac ierr = PetscSpaceView(subsp->origSpace,viewer);CHKERRQ(ierr); 97*20cf1dd8SToby Isaac ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 98*20cf1dd8SToby Isaac PetscFunctionReturn(0); 99*20cf1dd8SToby Isaac } 100*20cf1dd8SToby Isaac 101*20cf1dd8SToby Isaac static PetscErrorCode PetscSpaceEvaluate_Subspace(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) 102*20cf1dd8SToby Isaac { 103*20cf1dd8SToby Isaac PetscSpace_Subspace *subsp = (PetscSpace_Subspace *) sp->data; 104*20cf1dd8SToby Isaac PetscSpace origsp; 105*20cf1dd8SToby Isaac PetscInt origDim, subDim, origNc, subNc, subNb, origNb, i, j, k, l, m, n, o; 106*20cf1dd8SToby Isaac PetscReal *inpoints, *inB = NULL, *inD = NULL, *inH = NULL; 107*20cf1dd8SToby Isaac PetscErrorCode ierr; 108*20cf1dd8SToby Isaac 109*20cf1dd8SToby Isaac PetscFunctionBegin; 110*20cf1dd8SToby Isaac origsp = subsp->origSpace; 111*20cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(sp,&subDim);CHKERRQ(ierr); 112*20cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(origsp,&origDim);CHKERRQ(ierr); 113*20cf1dd8SToby Isaac ierr = PetscSpaceGetNumComponents(sp,&subNc);CHKERRQ(ierr); 114*20cf1dd8SToby Isaac ierr = PetscSpaceGetNumComponents(origsp,&origNc);CHKERRQ(ierr); 115*20cf1dd8SToby Isaac ierr = PetscSpaceGetDimension(sp,&subNb);CHKERRQ(ierr); 116*20cf1dd8SToby Isaac ierr = PetscSpaceGetDimension(origsp,&origNb);CHKERRQ(ierr); 117*20cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,npoints*origDim,MPIU_REAL,&inpoints);CHKERRQ(ierr); 118*20cf1dd8SToby Isaac for (i = 0; i < npoints; i++) { 119*20cf1dd8SToby Isaac if (subsp->x) { 120*20cf1dd8SToby Isaac for (j = 0; j < origDim; j++) inpoints[i * origDim + j] = subsp->x[j]; 121*20cf1dd8SToby Isaac } else { 122*20cf1dd8SToby Isaac for (j = 0; j < origDim; j++) inpoints[i * origDim + j] = 0.0; 123*20cf1dd8SToby Isaac } 124*20cf1dd8SToby Isaac if (subsp->Jx) { 125*20cf1dd8SToby Isaac for (j = 0; j < origDim; j++) { 126*20cf1dd8SToby Isaac for (k = 0; k < subDim; k++) { 127*20cf1dd8SToby Isaac inpoints[i * origDim + j] += subsp->Jx[j * subDim + k] * points[i * subDim + k]; 128*20cf1dd8SToby Isaac } 129*20cf1dd8SToby Isaac } 130*20cf1dd8SToby Isaac } else { 131*20cf1dd8SToby Isaac for (j = 0; j < PetscMin(subDim, origDim); j++) { 132*20cf1dd8SToby Isaac inpoints[i * origDim + j] += points[i * subDim + j]; 133*20cf1dd8SToby Isaac } 134*20cf1dd8SToby Isaac } 135*20cf1dd8SToby Isaac } 136*20cf1dd8SToby Isaac if (B) { 137*20cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,npoints*origNb*origNc,MPIU_REAL,&inB);CHKERRQ(ierr); 138*20cf1dd8SToby Isaac } 139*20cf1dd8SToby Isaac if (D) { 140*20cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,npoints*origNb*origNc*origDim,MPIU_REAL,&inD);CHKERRQ(ierr); 141*20cf1dd8SToby Isaac } 142*20cf1dd8SToby Isaac if (H) { 143*20cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,npoints*origNb*origNc*origDim*origDim,MPIU_REAL,&inH);CHKERRQ(ierr); 144*20cf1dd8SToby Isaac } 145*20cf1dd8SToby Isaac ierr = PetscSpaceEvaluate(origsp,npoints,inpoints,inB,inD,inH);CHKERRQ(ierr); 146*20cf1dd8SToby Isaac if (H) { 147*20cf1dd8SToby Isaac PetscReal *phi, *psi; 148*20cf1dd8SToby Isaac 149*20cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,origNc*origDim*origDim,MPIU_REAL,&phi);CHKERRQ(ierr); 150*20cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,origNc*subDim*subDim,MPIU_REAL,&psi);CHKERRQ(ierr); 151*20cf1dd8SToby Isaac for (i = 0; i < npoints * subNb * subNc * subDim; i++) D[i] = 0.0; 152*20cf1dd8SToby Isaac for (i = 0; i < subNb; i++) { 153*20cf1dd8SToby Isaac const PetscReal *subq = &subsp->Q[i * origNb]; 154*20cf1dd8SToby Isaac 155*20cf1dd8SToby Isaac for (j = 0; j < npoints; j++) { 156*20cf1dd8SToby Isaac for (k = 0; k < origNc * origDim; k++) phi[k] = 0.; 157*20cf1dd8SToby Isaac for (k = 0; k < origNc * subDim; k++) psi[k] = 0.; 158*20cf1dd8SToby Isaac for (k = 0; k < origNb; k++) { 159*20cf1dd8SToby Isaac for (l = 0; l < origNc * origDim * origDim; l++) { 160*20cf1dd8SToby Isaac phi[l] += inH[(j * origNb + k) * origNc * origDim * origDim + l] * subq[k]; 161*20cf1dd8SToby Isaac } 162*20cf1dd8SToby Isaac } 163*20cf1dd8SToby Isaac if (subsp->Jx) { 164*20cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 165*20cf1dd8SToby Isaac for (l = 0; l < subDim; l++) { 166*20cf1dd8SToby Isaac for (m = 0; m < origDim; m++) { 167*20cf1dd8SToby Isaac for (n = 0; n < subDim; n++) { 168*20cf1dd8SToby Isaac for (o = 0; o < origDim; o++) { 169*20cf1dd8SToby Isaac psi[(k * subDim + l) * subDim + n] += subsp->Jx[m * subDim + l] * subsp->Jx[o * subDim + n] * phi[(k * origDim + m) * origDim + o]; 170*20cf1dd8SToby Isaac } 171*20cf1dd8SToby Isaac } 172*20cf1dd8SToby Isaac } 173*20cf1dd8SToby Isaac } 174*20cf1dd8SToby Isaac } 175*20cf1dd8SToby Isaac } else { 176*20cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 177*20cf1dd8SToby Isaac for (l = 0; l < PetscMin(subDim, origDim); l++) { 178*20cf1dd8SToby Isaac for (m = 0; m < PetscMin(subDim, origDim); m++) { 179*20cf1dd8SToby Isaac psi[(k * subDim + l) * subDim + m] += phi[(k * origDim + l) * origDim + m]; 180*20cf1dd8SToby Isaac } 181*20cf1dd8SToby Isaac } 182*20cf1dd8SToby Isaac } 183*20cf1dd8SToby Isaac } 184*20cf1dd8SToby Isaac if (subsp->Ju) { 185*20cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 186*20cf1dd8SToby Isaac for (l = 0; l < origNc; l++) { 187*20cf1dd8SToby Isaac for (m = 0; m < subDim * subDim; m++) { 188*20cf1dd8SToby Isaac H[((j * subNb + i) * subNc + k) * subDim * subDim + m] += subsp->Ju[k * origNc + l] * psi[l * subDim * subDim + m]; 189*20cf1dd8SToby Isaac } 190*20cf1dd8SToby Isaac } 191*20cf1dd8SToby Isaac } 192*20cf1dd8SToby Isaac } 193*20cf1dd8SToby Isaac else { 194*20cf1dd8SToby Isaac for (k = 0; k < PetscMin(subNc, origNc); k++) { 195*20cf1dd8SToby Isaac for (l = 0; l < subDim * subDim; l++) { 196*20cf1dd8SToby Isaac H[((j * subNb + i) * subNc + k) * subDim * subDim + l] += psi[k * subDim * subDim + l]; 197*20cf1dd8SToby Isaac } 198*20cf1dd8SToby Isaac } 199*20cf1dd8SToby Isaac } 200*20cf1dd8SToby Isaac } 201*20cf1dd8SToby Isaac } 202*20cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,subNc*origDim,MPIU_REAL,&psi);CHKERRQ(ierr); 203*20cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,origNc*origDim,MPIU_REAL,&phi);CHKERRQ(ierr); 204*20cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,npoints*origNb*origNc*origDim,MPIU_REAL,&inH);CHKERRQ(ierr); 205*20cf1dd8SToby Isaac } 206*20cf1dd8SToby Isaac if (D) { 207*20cf1dd8SToby Isaac PetscReal *phi, *psi; 208*20cf1dd8SToby Isaac 209*20cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,origNc*origDim,MPIU_REAL,&phi);CHKERRQ(ierr); 210*20cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,origNc*subDim,MPIU_REAL,&psi);CHKERRQ(ierr); 211*20cf1dd8SToby Isaac for (i = 0; i < npoints * subNb * subNc * subDim; i++) D[i] = 0.0; 212*20cf1dd8SToby Isaac for (i = 0; i < subNb; i++) { 213*20cf1dd8SToby Isaac const PetscReal *subq = &subsp->Q[i * origNb]; 214*20cf1dd8SToby Isaac 215*20cf1dd8SToby Isaac for (j = 0; j < npoints; j++) { 216*20cf1dd8SToby Isaac for (k = 0; k < origNc * origDim; k++) phi[k] = 0.; 217*20cf1dd8SToby Isaac for (k = 0; k < origNc * subDim; k++) psi[k] = 0.; 218*20cf1dd8SToby Isaac for (k = 0; k < origNb; k++) { 219*20cf1dd8SToby Isaac for (l = 0; l < origNc * origDim; l++) { 220*20cf1dd8SToby Isaac phi[l] += inD[(j * origNb + k) * origNc * origDim + l] * subq[k]; 221*20cf1dd8SToby Isaac } 222*20cf1dd8SToby Isaac } 223*20cf1dd8SToby Isaac if (subsp->Jx) { 224*20cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 225*20cf1dd8SToby Isaac for (l = 0; l < subDim; l++) { 226*20cf1dd8SToby Isaac for (m = 0; m < origDim; m++) { 227*20cf1dd8SToby Isaac psi[k * subDim + l] += subsp->Jx[m * subDim + l] * phi[k * origDim + m]; 228*20cf1dd8SToby Isaac } 229*20cf1dd8SToby Isaac } 230*20cf1dd8SToby Isaac } 231*20cf1dd8SToby Isaac } else { 232*20cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 233*20cf1dd8SToby Isaac for (l = 0; l < PetscMin(subDim, origDim); l++) { 234*20cf1dd8SToby Isaac psi[k * subDim + l] += phi[k * origDim + l]; 235*20cf1dd8SToby Isaac } 236*20cf1dd8SToby Isaac } 237*20cf1dd8SToby Isaac } 238*20cf1dd8SToby Isaac if (subsp->Ju) { 239*20cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 240*20cf1dd8SToby Isaac for (l = 0; l < origNc; l++) { 241*20cf1dd8SToby Isaac for (m = 0; m < subDim; m++) { 242*20cf1dd8SToby Isaac D[((j * subNb + i) * subNc + k) * subDim + m] += subsp->Ju[k * origNc + l] * psi[l * subDim + m]; 243*20cf1dd8SToby Isaac } 244*20cf1dd8SToby Isaac } 245*20cf1dd8SToby Isaac } 246*20cf1dd8SToby Isaac } 247*20cf1dd8SToby Isaac else { 248*20cf1dd8SToby Isaac for (k = 0; k < PetscMin(subNc, origNc); k++) { 249*20cf1dd8SToby Isaac for (l = 0; l < subDim; l++) { 250*20cf1dd8SToby Isaac D[((j * subNb + i) * subNc + k) * subDim + l] += psi[k * subDim + l]; 251*20cf1dd8SToby Isaac } 252*20cf1dd8SToby Isaac } 253*20cf1dd8SToby Isaac } 254*20cf1dd8SToby Isaac } 255*20cf1dd8SToby Isaac } 256*20cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,subNc*origDim,MPIU_REAL,&psi);CHKERRQ(ierr); 257*20cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,origNc*origDim,MPIU_REAL,&phi);CHKERRQ(ierr); 258*20cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,npoints*origNb*origNc*origDim,MPIU_REAL,&inD);CHKERRQ(ierr); 259*20cf1dd8SToby Isaac } 260*20cf1dd8SToby Isaac if (B) { 261*20cf1dd8SToby Isaac PetscReal *phi; 262*20cf1dd8SToby Isaac 263*20cf1dd8SToby Isaac ierr = DMGetWorkArray(sp->dm,origNc,MPIU_REAL,&phi);CHKERRQ(ierr); 264*20cf1dd8SToby Isaac if (subsp->u) { 265*20cf1dd8SToby Isaac for (i = 0; i < npoints * subNb; i++) { 266*20cf1dd8SToby Isaac for (j = 0; j < subNc; j++) B[i * subNc + j] = subsp->u[j]; 267*20cf1dd8SToby Isaac } 268*20cf1dd8SToby Isaac } else { 269*20cf1dd8SToby Isaac for (i = 0; i < npoints * subNb * subNc; i++) B[i] = 0.0; 270*20cf1dd8SToby Isaac } 271*20cf1dd8SToby Isaac for (i = 0; i < subNb; i++) { 272*20cf1dd8SToby Isaac const PetscReal *subq = &subsp->Q[i * origNb]; 273*20cf1dd8SToby Isaac 274*20cf1dd8SToby Isaac for (j = 0; j < npoints; j++) { 275*20cf1dd8SToby Isaac for (k = 0; k < origNc; k++) phi[k] = 0.; 276*20cf1dd8SToby Isaac for (k = 0; k < origNb; k++) { 277*20cf1dd8SToby Isaac for (l = 0; l < origNc; l++) { 278*20cf1dd8SToby Isaac phi[l] += inB[(j * origNb + k) * origNc + l] * subq[k]; 279*20cf1dd8SToby Isaac } 280*20cf1dd8SToby Isaac } 281*20cf1dd8SToby Isaac if (subsp->Ju) { 282*20cf1dd8SToby Isaac for (k = 0; k < subNc; k++) { 283*20cf1dd8SToby Isaac for (l = 0; l < origNc; l++) { 284*20cf1dd8SToby Isaac B[(j * subNb + i) * subNc + k] += subsp->Ju[k * origNc + l] * phi[l]; 285*20cf1dd8SToby Isaac } 286*20cf1dd8SToby Isaac } 287*20cf1dd8SToby Isaac } 288*20cf1dd8SToby Isaac else { 289*20cf1dd8SToby Isaac for (k = 0; k < PetscMin(subNc, origNc); k++) { 290*20cf1dd8SToby Isaac B[(j * subNb + i) * subNc + k] += phi[k]; 291*20cf1dd8SToby Isaac } 292*20cf1dd8SToby Isaac } 293*20cf1dd8SToby Isaac } 294*20cf1dd8SToby Isaac } 295*20cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,origNc,MPIU_REAL,&phi);CHKERRQ(ierr); 296*20cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,npoints*origNb*origNc,MPIU_REAL,&inB);CHKERRQ(ierr); 297*20cf1dd8SToby Isaac } 298*20cf1dd8SToby Isaac ierr = DMRestoreWorkArray(sp->dm,npoints*origDim,MPIU_REAL,&inpoints);CHKERRQ(ierr); 299*20cf1dd8SToby Isaac PetscFunctionReturn(0); 300*20cf1dd8SToby Isaac } 301*20cf1dd8SToby Isaac 302*20cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Subspace(PetscSpace sp) 303*20cf1dd8SToby Isaac { 304*20cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 305*20cf1dd8SToby Isaac 306*20cf1dd8SToby Isaac PetscErrorCode ierr; 307*20cf1dd8SToby Isaac ierr = PetscNewLog(sp,&subsp);CHKERRQ(ierr); 308*20cf1dd8SToby Isaac sp->data = (void *) subsp; 309*20cf1dd8SToby Isaac PetscFunctionReturn(0); 310*20cf1dd8SToby Isaac } 311*20cf1dd8SToby Isaac 312*20cf1dd8SToby Isaac static PetscErrorCode PetscSpaceGetDimension_Subspace(PetscSpace sp, PetscInt *dim) 313*20cf1dd8SToby Isaac { 314*20cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 315*20cf1dd8SToby Isaac 316*20cf1dd8SToby Isaac PetscFunctionBegin; 317*20cf1dd8SToby Isaac subsp = (PetscSpace_Subspace *) sp->data; 318*20cf1dd8SToby Isaac *dim = subsp->Nb; 319*20cf1dd8SToby Isaac PetscFunctionReturn(0); 320*20cf1dd8SToby Isaac } 321*20cf1dd8SToby Isaac 322*20cf1dd8SToby Isaac static PetscErrorCode PetscSpaceSetUp_Subspace(PetscSpace sp) 323*20cf1dd8SToby Isaac { 324*20cf1dd8SToby Isaac const PetscReal *x; 325*20cf1dd8SToby Isaac const PetscReal *Jx; 326*20cf1dd8SToby Isaac const PetscReal *u; 327*20cf1dd8SToby Isaac const PetscReal *Ju; 328*20cf1dd8SToby Isaac PetscDualSpace dualSubspace; 329*20cf1dd8SToby Isaac PetscSpace origSpace; 330*20cf1dd8SToby Isaac PetscInt origDim, subDim, origNc, subNc, origNb, subNb, f, i, j, numPoints, offset; 331*20cf1dd8SToby Isaac PetscReal *allPoints, *allWeights, *B, *V; 332*20cf1dd8SToby Isaac DM dm; 333*20cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 334*20cf1dd8SToby Isaac PetscErrorCode ierr; 335*20cf1dd8SToby Isaac 336*20cf1dd8SToby Isaac PetscFunctionBegin; 337*20cf1dd8SToby Isaac subsp = (PetscSpace_Subspace *) sp->data; 338*20cf1dd8SToby Isaac x = subsp->x; 339*20cf1dd8SToby Isaac Jx = subsp->Jx; 340*20cf1dd8SToby Isaac u = subsp->u; 341*20cf1dd8SToby Isaac Ju = subsp->Ju; 342*20cf1dd8SToby Isaac origSpace = subsp->origSpace; 343*20cf1dd8SToby Isaac dualSubspace = subsp->dualSubspace; 344*20cf1dd8SToby Isaac ierr = PetscSpaceGetNumComponents(origSpace,&origNc);CHKERRQ(ierr); 345*20cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(origSpace,&origDim);CHKERRQ(ierr); 346*20cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(dualSubspace,&dm);CHKERRQ(ierr); 347*20cf1dd8SToby Isaac ierr = DMGetDimension(dm,&subDim);CHKERRQ(ierr); 348*20cf1dd8SToby Isaac ierr = PetscSpaceGetDimension(origSpace,&origNb);CHKERRQ(ierr); 349*20cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(dualSubspace,&subNb);CHKERRQ(ierr); 350*20cf1dd8SToby Isaac ierr = PetscDualSpaceGetNumComponents(dualSubspace,&subNc);CHKERRQ(ierr); 351*20cf1dd8SToby Isaac 352*20cf1dd8SToby Isaac for (f = 0, numPoints = 0; f < subNb; f++) { 353*20cf1dd8SToby Isaac PetscQuadrature q; 354*20cf1dd8SToby Isaac PetscInt qNp; 355*20cf1dd8SToby Isaac 356*20cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(dualSubspace,f,&q);CHKERRQ(ierr); 357*20cf1dd8SToby Isaac ierr = PetscQuadratureGetData(q,NULL,NULL,&qNp,NULL,NULL);CHKERRQ(ierr); 358*20cf1dd8SToby Isaac numPoints += qNp; 359*20cf1dd8SToby Isaac } 360*20cf1dd8SToby Isaac ierr = PetscMalloc1(subNb*origNb,&V);CHKERRQ(ierr); 361*20cf1dd8SToby Isaac ierr = PetscMalloc3(numPoints*origDim,&allPoints,numPoints*origNc,&allWeights,numPoints*origNb*origNc,&B);CHKERRQ(ierr); 362*20cf1dd8SToby Isaac for (f = 0, offset = 0; f < subNb; f++) { 363*20cf1dd8SToby Isaac PetscQuadrature q; 364*20cf1dd8SToby Isaac PetscInt qNp, p; 365*20cf1dd8SToby Isaac const PetscReal *qp; 366*20cf1dd8SToby Isaac const PetscReal *qw; 367*20cf1dd8SToby Isaac 368*20cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(dualSubspace,f,&q);CHKERRQ(ierr); 369*20cf1dd8SToby Isaac ierr = PetscQuadratureGetData(q,NULL,NULL,&qNp,&qp,&qw);CHKERRQ(ierr); 370*20cf1dd8SToby Isaac for (p = 0; p < qNp; p++, offset++) { 371*20cf1dd8SToby Isaac if (x) { 372*20cf1dd8SToby Isaac for (i = 0; i < origDim; i++) allPoints[origDim * offset + i] = x[i]; 373*20cf1dd8SToby Isaac } else { 374*20cf1dd8SToby Isaac for (i = 0; i < origDim; i++) allPoints[origDim * offset + i] = 0.0; 375*20cf1dd8SToby Isaac } 376*20cf1dd8SToby Isaac if (Jx) { 377*20cf1dd8SToby Isaac for (i = 0; i < origDim; i++) { 378*20cf1dd8SToby Isaac for (j = 0; j < subDim; j++) { 379*20cf1dd8SToby Isaac allPoints[origDim * offset + i] += Jx[i * subDim + j] * qp[j]; 380*20cf1dd8SToby Isaac } 381*20cf1dd8SToby Isaac } 382*20cf1dd8SToby Isaac } else { 383*20cf1dd8SToby Isaac for (i = 0; i < PetscMin(subDim, origDim); i++) allPoints[origDim * offset + i] += qp[i]; 384*20cf1dd8SToby Isaac } 385*20cf1dd8SToby Isaac for (i = 0; i < origNc; i++) allWeights[origNc * offset + i] = 0.0; 386*20cf1dd8SToby Isaac if (Ju) { 387*20cf1dd8SToby Isaac for (i = 0; i < origNc; i++) { 388*20cf1dd8SToby Isaac for (j = 0; j < subNc; j++) { 389*20cf1dd8SToby Isaac allWeights[offset * origNc + i] += qw[j] * Ju[j * origNc + i]; 390*20cf1dd8SToby Isaac } 391*20cf1dd8SToby Isaac } 392*20cf1dd8SToby Isaac } else { 393*20cf1dd8SToby Isaac for (i = 0; i < PetscMin(subNc, origNc); i++) allWeights[offset * origNc + i] += qw[i]; 394*20cf1dd8SToby Isaac } 395*20cf1dd8SToby Isaac } 396*20cf1dd8SToby Isaac } 397*20cf1dd8SToby Isaac ierr = PetscSpaceEvaluate(origSpace,numPoints,allPoints,B,NULL,NULL);CHKERRQ(ierr); 398*20cf1dd8SToby Isaac for (f = 0, offset = 0; f < subNb; f++) { 399*20cf1dd8SToby Isaac PetscInt b, p, s, qNp; 400*20cf1dd8SToby Isaac PetscQuadrature q; 401*20cf1dd8SToby Isaac const PetscReal *qw; 402*20cf1dd8SToby Isaac 403*20cf1dd8SToby Isaac ierr = PetscDualSpaceGetFunctional(dualSubspace,f,&q);CHKERRQ(ierr); 404*20cf1dd8SToby Isaac ierr = PetscQuadratureGetData(q,NULL,NULL,&qNp,NULL,&qw);CHKERRQ(ierr); 405*20cf1dd8SToby Isaac if (u) { 406*20cf1dd8SToby Isaac for (b = 0; b < origNb; b++) { 407*20cf1dd8SToby Isaac for (s = 0; s < subNc; s++) { 408*20cf1dd8SToby Isaac V[f * origNb + b] += qw[s] * u[s]; 409*20cf1dd8SToby Isaac } 410*20cf1dd8SToby Isaac } 411*20cf1dd8SToby Isaac } else { 412*20cf1dd8SToby Isaac for (b = 0; b < origNb; b++) V[f * origNb + b] = 0.0; 413*20cf1dd8SToby Isaac } 414*20cf1dd8SToby Isaac for (p = 0; p < qNp; p++, offset++) { 415*20cf1dd8SToby Isaac for (b = 0; b < origNb; b++) { 416*20cf1dd8SToby Isaac for (s = 0; s < origNc; s++) { 417*20cf1dd8SToby Isaac V[f * origNb + b] += B[(offset * origNb + b) * origNc + s] * allWeights[offset * origNc + s]; 418*20cf1dd8SToby Isaac } 419*20cf1dd8SToby Isaac } 420*20cf1dd8SToby Isaac } 421*20cf1dd8SToby Isaac } 422*20cf1dd8SToby Isaac /* orthnormalize rows of V */ 423*20cf1dd8SToby Isaac for (f = 0; f < subNb; f++) { 424*20cf1dd8SToby Isaac PetscReal rho = 0.0, scal; 425*20cf1dd8SToby Isaac 426*20cf1dd8SToby Isaac for (i = 0; i < origNb; i++) rho += PetscSqr(V[f * origNb + i]); 427*20cf1dd8SToby Isaac 428*20cf1dd8SToby Isaac scal = 1. / PetscSqrtReal(rho); 429*20cf1dd8SToby Isaac 430*20cf1dd8SToby Isaac for (i = 0; i < origNb; i++) V[f * origNb + i] *= scal; 431*20cf1dd8SToby Isaac for (j = f + 1; j < subNb; j++) { 432*20cf1dd8SToby Isaac for (i = 0, scal = 0.; i < origNb; i++) scal += V[f * origNb + i] * V[j * origNb + i]; 433*20cf1dd8SToby Isaac for (i = 0; i < origNb; i++) V[j * origNb + i] -= V[f * origNb + i] * scal; 434*20cf1dd8SToby Isaac } 435*20cf1dd8SToby Isaac } 436*20cf1dd8SToby Isaac ierr = PetscFree3(allPoints,allWeights,B);CHKERRQ(ierr); 437*20cf1dd8SToby Isaac subsp->Q = V; 438*20cf1dd8SToby Isaac PetscFunctionReturn(0); 439*20cf1dd8SToby Isaac } 440*20cf1dd8SToby Isaac 441*20cf1dd8SToby Isaac static PetscErrorCode PetscSpacePolynomialGetTensor_Subspace(PetscSpace sp, PetscBool *poly) 442*20cf1dd8SToby Isaac { 443*20cf1dd8SToby Isaac PetscSpace_Subspace *subsp = (PetscSpace_Subspace *) sp->data; 444*20cf1dd8SToby Isaac PetscErrorCode ierr; 445*20cf1dd8SToby Isaac 446*20cf1dd8SToby Isaac PetscFunctionBegin; 447*20cf1dd8SToby Isaac *poly = PETSC_FALSE; 448*20cf1dd8SToby Isaac ierr = PetscSpacePolynomialGetTensor(subsp->origSpace,poly);CHKERRQ(ierr); 449*20cf1dd8SToby Isaac if (*poly) { 450*20cf1dd8SToby Isaac if (subsp->Jx) { 451*20cf1dd8SToby Isaac PetscInt subDim, origDim, i, j; 452*20cf1dd8SToby Isaac PetscInt maxnnz; 453*20cf1dd8SToby Isaac 454*20cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(subsp->origSpace,&origDim);CHKERRQ(ierr); 455*20cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(sp,&subDim);CHKERRQ(ierr); 456*20cf1dd8SToby Isaac maxnnz = 0; 457*20cf1dd8SToby Isaac for (i = 0; i < origDim; i++) { 458*20cf1dd8SToby Isaac PetscInt nnz = 0; 459*20cf1dd8SToby Isaac 460*20cf1dd8SToby Isaac for (j = 0; j < subDim; j++) nnz += (subsp->Jx[i * subDim + j] != 0.); 461*20cf1dd8SToby Isaac maxnnz = PetscMax(maxnnz,nnz); 462*20cf1dd8SToby Isaac } 463*20cf1dd8SToby Isaac for (j = 0; j < subDim; j++) { 464*20cf1dd8SToby Isaac PetscInt nnz = 0; 465*20cf1dd8SToby Isaac 466*20cf1dd8SToby Isaac for (i = 0; i < origDim; i++) nnz += (subsp->Jx[i * subDim + j] != 0.); 467*20cf1dd8SToby Isaac maxnnz = PetscMax(maxnnz,nnz); 468*20cf1dd8SToby Isaac } 469*20cf1dd8SToby Isaac if (maxnnz > 1) *poly = PETSC_FALSE; 470*20cf1dd8SToby Isaac } 471*20cf1dd8SToby Isaac } 472*20cf1dd8SToby Isaac PetscFunctionReturn(0); 473*20cf1dd8SToby Isaac } 474*20cf1dd8SToby Isaac 475*20cf1dd8SToby Isaac static PetscErrorCode PetscSpaceInitialize_Subspace(PetscSpace sp) 476*20cf1dd8SToby Isaac { 477*20cf1dd8SToby Isaac PetscErrorCode ierr; 478*20cf1dd8SToby Isaac 479*20cf1dd8SToby Isaac PetscFunctionBegin; 480*20cf1dd8SToby Isaac sp->ops->setup = PetscSpaceSetUp_Subspace; 481*20cf1dd8SToby Isaac sp->ops->view = PetscSpaceView_Subspace; 482*20cf1dd8SToby Isaac sp->ops->destroy = PetscSpaceDestroy_Subspace; 483*20cf1dd8SToby Isaac sp->ops->getdimension = PetscSpaceGetDimension_Subspace; 484*20cf1dd8SToby Isaac sp->ops->evaluate = PetscSpaceEvaluate_Subspace; 485*20cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Subspace);CHKERRQ(ierr); 486*20cf1dd8SToby Isaac PetscFunctionReturn(0); 487*20cf1dd8SToby Isaac } 488*20cf1dd8SToby Isaac 489*20cf1dd8SToby Isaac PetscErrorCode PetscSpaceCreateSubspace(PetscSpace origSpace, PetscDualSpace dualSubspace, PetscReal *x, PetscReal *Jx, PetscReal *u, PetscReal *Ju, PetscCopyMode copymode, PetscSpace *subspace) 490*20cf1dd8SToby Isaac { 491*20cf1dd8SToby Isaac PetscSpace_Subspace *subsp; 492*20cf1dd8SToby Isaac PetscInt origDim, subDim, origNc, subNc, subNb; 493*20cf1dd8SToby Isaac PetscInt order; 494*20cf1dd8SToby Isaac DM dm; 495*20cf1dd8SToby Isaac PetscErrorCode ierr; 496*20cf1dd8SToby Isaac 497*20cf1dd8SToby Isaac PetscFunctionBegin; 498*20cf1dd8SToby Isaac PetscValidHeaderSpecific(origSpace,PETSCSPACE_CLASSID,1); 499*20cf1dd8SToby Isaac PetscValidHeaderSpecific(dualSubspace,PETSCDUALSPACE_CLASSID,2); 500*20cf1dd8SToby Isaac if (x) PetscValidRealPointer(x,3); 501*20cf1dd8SToby Isaac if (Jx) PetscValidRealPointer(Jx,4); 502*20cf1dd8SToby Isaac if (u) PetscValidRealPointer(u,5); 503*20cf1dd8SToby Isaac if (Ju) PetscValidRealPointer(Ju,6); 504*20cf1dd8SToby Isaac PetscValidPointer(subspace,7); 505*20cf1dd8SToby Isaac ierr = PetscSpaceGetNumComponents(origSpace,&origNc);CHKERRQ(ierr); 506*20cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(origSpace,&origDim);CHKERRQ(ierr); 507*20cf1dd8SToby Isaac ierr = PetscDualSpaceGetDM(dualSubspace,&dm);CHKERRQ(ierr); 508*20cf1dd8SToby Isaac ierr = DMGetDimension(dm,&subDim);CHKERRQ(ierr); 509*20cf1dd8SToby Isaac ierr = PetscDualSpaceGetDimension(dualSubspace,&subNb);CHKERRQ(ierr); 510*20cf1dd8SToby Isaac ierr = PetscDualSpaceGetNumComponents(dualSubspace,&subNc);CHKERRQ(ierr); 511*20cf1dd8SToby Isaac ierr = PetscSpaceCreate(PetscObjectComm((PetscObject)origSpace),subspace);CHKERRQ(ierr); 512*20cf1dd8SToby Isaac ierr = PetscSpaceSetType(*subspace,PETSCSPACESUBSPACE);CHKERRQ(ierr); 513*20cf1dd8SToby Isaac ierr = PetscSpaceSetNumVariables(*subspace,subDim);CHKERRQ(ierr); 514*20cf1dd8SToby Isaac ierr = PetscSpaceSetNumComponents(*subspace,subNc);CHKERRQ(ierr); 515*20cf1dd8SToby Isaac ierr = PetscSpaceGetDegree(origSpace,&order,NULL);CHKERRQ(ierr); 516*20cf1dd8SToby Isaac ierr = PetscSpaceSetDegree(*subspace,order);CHKERRQ(ierr); 517*20cf1dd8SToby Isaac subsp = (PetscSpace_Subspace *) (*subspace)->data; 518*20cf1dd8SToby Isaac subsp->Nb = subNb; 519*20cf1dd8SToby Isaac switch (copymode) { 520*20cf1dd8SToby Isaac case PETSC_OWN_POINTER: 521*20cf1dd8SToby Isaac if (x) subsp->x_alloc = x; 522*20cf1dd8SToby Isaac if (Jx) subsp->Jx_alloc = Jx; 523*20cf1dd8SToby Isaac if (u) subsp->u_alloc = u; 524*20cf1dd8SToby Isaac if (Ju) subsp->Ju_alloc = Ju; 525*20cf1dd8SToby Isaac case PETSC_USE_POINTER: 526*20cf1dd8SToby Isaac if (x) subsp->x = x; 527*20cf1dd8SToby Isaac if (Jx) subsp->Jx = Jx; 528*20cf1dd8SToby Isaac if (u) subsp->u = u; 529*20cf1dd8SToby Isaac if (Ju) subsp->Ju = Ju; 530*20cf1dd8SToby Isaac break; 531*20cf1dd8SToby Isaac case PETSC_COPY_VALUES: 532*20cf1dd8SToby Isaac if (x) { 533*20cf1dd8SToby Isaac ierr = PetscMalloc1(origDim,&subsp->x_alloc);CHKERRQ(ierr); 534*20cf1dd8SToby Isaac ierr = PetscMemcpy(subsp->x_alloc,x,origDim*sizeof(*subsp->x_alloc));CHKERRQ(ierr); 535*20cf1dd8SToby Isaac subsp->x = subsp->x_alloc; 536*20cf1dd8SToby Isaac } 537*20cf1dd8SToby Isaac if (Jx) { 538*20cf1dd8SToby Isaac ierr = PetscMalloc1(origDim * subDim,&subsp->Jx_alloc);CHKERRQ(ierr); 539*20cf1dd8SToby Isaac ierr = PetscMemcpy(subsp->Jx_alloc,Jx,origDim * subDim*sizeof(*subsp->Jx_alloc));CHKERRQ(ierr); 540*20cf1dd8SToby Isaac subsp->Jx = subsp->Jx_alloc; 541*20cf1dd8SToby Isaac } 542*20cf1dd8SToby Isaac if (u) { 543*20cf1dd8SToby Isaac ierr = PetscMalloc1(subNc,&subsp->u_alloc);CHKERRQ(ierr); 544*20cf1dd8SToby Isaac ierr = PetscMemcpy(subsp->u_alloc,u,subNc*sizeof(*subsp->u_alloc));CHKERRQ(ierr); 545*20cf1dd8SToby Isaac subsp->u = subsp->u_alloc; 546*20cf1dd8SToby Isaac } 547*20cf1dd8SToby Isaac if (Ju) { 548*20cf1dd8SToby Isaac ierr = PetscMalloc1(origNc * subNc,&subsp->Ju_alloc);CHKERRQ(ierr); 549*20cf1dd8SToby Isaac ierr = PetscMemcpy(subsp->Ju_alloc,Ju,origNc * subNc*sizeof(*subsp->Ju_alloc));CHKERRQ(ierr); 550*20cf1dd8SToby Isaac subsp->Ju = subsp->Ju_alloc; 551*20cf1dd8SToby Isaac } 552*20cf1dd8SToby Isaac break; 553*20cf1dd8SToby Isaac default: 554*20cf1dd8SToby Isaac SETERRQ(PetscObjectComm((PetscObject)origSpace),PETSC_ERR_ARG_OUTOFRANGE,"Unknown copy mode"); 555*20cf1dd8SToby Isaac } 556*20cf1dd8SToby Isaac ierr = PetscObjectReference((PetscObject)origSpace);CHKERRQ(ierr); 557*20cf1dd8SToby Isaac subsp->origSpace = origSpace; 558*20cf1dd8SToby Isaac ierr = PetscObjectReference((PetscObject)dualSubspace);CHKERRQ(ierr); 559*20cf1dd8SToby Isaac subsp->dualSubspace = dualSubspace; 560*20cf1dd8SToby Isaac ierr = PetscSpaceInitialize_Subspace(*subspace);CHKERRQ(ierr); 561*20cf1dd8SToby Isaac PetscFunctionReturn(0); 562*20cf1dd8SToby Isaac } 563*20cf1dd8SToby Isaac 564