1*20cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 2*20cf1dd8SToby Isaac 3*20cf1dd8SToby Isaac PetscErrorCode PetscSpaceSetFromOptions_Polynomial(PetscOptionItems *PetscOptionsObject,PetscSpace sp) 4*20cf1dd8SToby Isaac { 5*20cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 6*20cf1dd8SToby Isaac PetscErrorCode ierr; 7*20cf1dd8SToby Isaac 8*20cf1dd8SToby Isaac PetscFunctionBegin; 9*20cf1dd8SToby Isaac ierr = PetscOptionsHead(PetscOptionsObject,"PetscSpace polynomial options");CHKERRQ(ierr); 10*20cf1dd8SToby Isaac ierr = PetscOptionsBool("-petscspace_poly_sym", "Use only symmetric polynomials", "PetscSpacePolynomialSetSymmetric", poly->symmetric, &poly->symmetric, NULL);CHKERRQ(ierr); 11*20cf1dd8SToby Isaac ierr = PetscOptionsBool("-petscspace_poly_tensor", "Use the tensor product polynomials", "PetscSpacePolynomialSetTensor", poly->tensor, &poly->tensor, NULL);CHKERRQ(ierr); 12*20cf1dd8SToby Isaac ierr = PetscOptionsTail();CHKERRQ(ierr); 13*20cf1dd8SToby Isaac PetscFunctionReturn(0); 14*20cf1dd8SToby Isaac } 15*20cf1dd8SToby Isaac 16*20cf1dd8SToby Isaac static PetscErrorCode PetscSpacePolynomialView_Ascii(PetscSpace sp, PetscViewer viewer) 17*20cf1dd8SToby Isaac { 18*20cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 19*20cf1dd8SToby Isaac PetscErrorCode ierr; 20*20cf1dd8SToby Isaac 21*20cf1dd8SToby Isaac PetscFunctionBegin; 22*20cf1dd8SToby Isaac if (poly->tensor) {ierr = PetscViewerASCIIPrintf(viewer, "Tensor polynomial space of degree %D\n", sp->degree);CHKERRQ(ierr);} 23*20cf1dd8SToby Isaac else {ierr = PetscViewerASCIIPrintf(viewer, "Polynomial space of degree %D\n", sp->degree);CHKERRQ(ierr);} 24*20cf1dd8SToby Isaac PetscFunctionReturn(0); 25*20cf1dd8SToby Isaac } 26*20cf1dd8SToby Isaac 27*20cf1dd8SToby Isaac PetscErrorCode PetscSpaceView_Polynomial(PetscSpace sp, PetscViewer viewer) 28*20cf1dd8SToby Isaac { 29*20cf1dd8SToby Isaac PetscBool iascii; 30*20cf1dd8SToby Isaac PetscErrorCode ierr; 31*20cf1dd8SToby Isaac 32*20cf1dd8SToby Isaac PetscFunctionBegin; 33*20cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 34*20cf1dd8SToby Isaac PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 35*20cf1dd8SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 36*20cf1dd8SToby Isaac if (iascii) {ierr = PetscSpacePolynomialView_Ascii(sp, viewer);CHKERRQ(ierr);} 37*20cf1dd8SToby Isaac PetscFunctionReturn(0); 38*20cf1dd8SToby Isaac } 39*20cf1dd8SToby Isaac 40*20cf1dd8SToby Isaac PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp) 41*20cf1dd8SToby Isaac { 42*20cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 43*20cf1dd8SToby Isaac PetscInt ndegree = sp->degree+1; 44*20cf1dd8SToby Isaac PetscInt deg; 45*20cf1dd8SToby Isaac PetscErrorCode ierr; 46*20cf1dd8SToby Isaac 47*20cf1dd8SToby Isaac PetscFunctionBegin; 48*20cf1dd8SToby Isaac ierr = PetscMalloc1(ndegree, &poly->degrees);CHKERRQ(ierr); 49*20cf1dd8SToby Isaac for (deg = 0; deg < ndegree; ++deg) poly->degrees[deg] = deg; 50*20cf1dd8SToby Isaac if (poly->tensor) { 51*20cf1dd8SToby Isaac sp->maxDegree = sp->degree + PetscMax(sp->Nv - 1,0); 52*20cf1dd8SToby Isaac } else { 53*20cf1dd8SToby Isaac sp->maxDegree = sp->degree; 54*20cf1dd8SToby Isaac } 55*20cf1dd8SToby Isaac PetscFunctionReturn(0); 56*20cf1dd8SToby Isaac } 57*20cf1dd8SToby Isaac 58*20cf1dd8SToby Isaac PetscErrorCode PetscSpaceDestroy_Polynomial(PetscSpace sp) 59*20cf1dd8SToby Isaac { 60*20cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 61*20cf1dd8SToby Isaac PetscErrorCode ierr; 62*20cf1dd8SToby Isaac 63*20cf1dd8SToby Isaac PetscFunctionBegin; 64*20cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", NULL);CHKERRQ(ierr); 65*20cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", NULL);CHKERRQ(ierr); 66*20cf1dd8SToby Isaac ierr = PetscFree(poly->degrees);CHKERRQ(ierr); 67*20cf1dd8SToby Isaac if (poly->subspaces) { 68*20cf1dd8SToby Isaac PetscInt d; 69*20cf1dd8SToby Isaac 70*20cf1dd8SToby Isaac for (d = 0; d < sp->Nv; ++d) { 71*20cf1dd8SToby Isaac ierr = PetscSpaceDestroy(&poly->subspaces[d]);CHKERRQ(ierr); 72*20cf1dd8SToby Isaac } 73*20cf1dd8SToby Isaac } 74*20cf1dd8SToby Isaac ierr = PetscFree(poly->subspaces);CHKERRQ(ierr); 75*20cf1dd8SToby Isaac ierr = PetscFree(poly);CHKERRQ(ierr); 76*20cf1dd8SToby Isaac PetscFunctionReturn(0); 77*20cf1dd8SToby Isaac } 78*20cf1dd8SToby Isaac 79*20cf1dd8SToby Isaac /* We treat the space as a tensor product of scalar polynomial spaces, so the dimension is multiplied by Nc */ 80*20cf1dd8SToby Isaac PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim) 81*20cf1dd8SToby Isaac { 82*20cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 83*20cf1dd8SToby Isaac PetscInt deg = sp->degree; 84*20cf1dd8SToby Isaac PetscInt n = sp->Nv, i; 85*20cf1dd8SToby Isaac PetscReal D = 1.0; 86*20cf1dd8SToby Isaac 87*20cf1dd8SToby Isaac PetscFunctionBegin; 88*20cf1dd8SToby Isaac if (poly->tensor) { 89*20cf1dd8SToby Isaac *dim = 1; 90*20cf1dd8SToby Isaac for (i = 0; i < n; ++i) *dim *= (deg+1); 91*20cf1dd8SToby Isaac } else { 92*20cf1dd8SToby Isaac for (i = 1; i <= n; ++i) { 93*20cf1dd8SToby Isaac D *= ((PetscReal) (deg+i))/i; 94*20cf1dd8SToby Isaac } 95*20cf1dd8SToby Isaac *dim = (PetscInt) (D + 0.5); 96*20cf1dd8SToby Isaac } 97*20cf1dd8SToby Isaac *dim *= sp->Nc; 98*20cf1dd8SToby Isaac PetscFunctionReturn(0); 99*20cf1dd8SToby Isaac } 100*20cf1dd8SToby Isaac 101*20cf1dd8SToby Isaac /* 102*20cf1dd8SToby Isaac LatticePoint_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to 'sum'. 103*20cf1dd8SToby Isaac 104*20cf1dd8SToby Isaac Input Parameters: 105*20cf1dd8SToby Isaac + len - The length of the tuple 106*20cf1dd8SToby Isaac . sum - The sum of all entries in the tuple 107*20cf1dd8SToby Isaac - ind - The current multi-index of the tuple, initialized to the 0 tuple 108*20cf1dd8SToby Isaac 109*20cf1dd8SToby Isaac Output Parameter: 110*20cf1dd8SToby Isaac + ind - The multi-index of the tuple, -1 indicates the iteration has terminated 111*20cf1dd8SToby Isaac . tup - A tuple of len integers addig to sum 112*20cf1dd8SToby Isaac 113*20cf1dd8SToby Isaac Level: developer 114*20cf1dd8SToby Isaac 115*20cf1dd8SToby Isaac .seealso: 116*20cf1dd8SToby Isaac */ 117*20cf1dd8SToby Isaac static PetscErrorCode LatticePoint_Internal(PetscInt len, PetscInt sum, PetscInt ind[], PetscInt tup[]) 118*20cf1dd8SToby Isaac { 119*20cf1dd8SToby Isaac PetscInt i; 120*20cf1dd8SToby Isaac PetscErrorCode ierr; 121*20cf1dd8SToby Isaac 122*20cf1dd8SToby Isaac PetscFunctionBegin; 123*20cf1dd8SToby Isaac if (len == 1) { 124*20cf1dd8SToby Isaac ind[0] = -1; 125*20cf1dd8SToby Isaac tup[0] = sum; 126*20cf1dd8SToby Isaac } else if (sum == 0) { 127*20cf1dd8SToby Isaac for (i = 0; i < len; ++i) {ind[0] = -1; tup[i] = 0;} 128*20cf1dd8SToby Isaac } else { 129*20cf1dd8SToby Isaac tup[0] = sum - ind[0]; 130*20cf1dd8SToby Isaac ierr = LatticePoint_Internal(len-1, ind[0], &ind[1], &tup[1]);CHKERRQ(ierr); 131*20cf1dd8SToby Isaac if (ind[1] < 0) { 132*20cf1dd8SToby Isaac if (ind[0] == sum) {ind[0] = -1;} 133*20cf1dd8SToby Isaac else {ind[1] = 0; ++ind[0];} 134*20cf1dd8SToby Isaac } 135*20cf1dd8SToby Isaac } 136*20cf1dd8SToby Isaac PetscFunctionReturn(0); 137*20cf1dd8SToby Isaac } 138*20cf1dd8SToby Isaac 139*20cf1dd8SToby Isaac /* 140*20cf1dd8SToby Isaac TensorPoint_Internal - Returns all tuples of size 'len' with nonnegative integers that are less than 'max'. 141*20cf1dd8SToby Isaac 142*20cf1dd8SToby Isaac Input Parameters: 143*20cf1dd8SToby Isaac + len - The length of the tuple 144*20cf1dd8SToby Isaac . max - The max for all entries in the tuple 145*20cf1dd8SToby Isaac - ind - The current multi-index of the tuple, initialized to the 0 tuple 146*20cf1dd8SToby Isaac 147*20cf1dd8SToby Isaac Output Parameter: 148*20cf1dd8SToby Isaac + ind - The multi-index of the tuple, -1 indicates the iteration has terminated 149*20cf1dd8SToby Isaac . tup - A tuple of len integers less than max 150*20cf1dd8SToby Isaac 151*20cf1dd8SToby Isaac Level: developer 152*20cf1dd8SToby Isaac 153*20cf1dd8SToby Isaac .seealso: 154*20cf1dd8SToby Isaac */ 155*20cf1dd8SToby Isaac static PetscErrorCode TensorPoint_Internal(PetscInt len, PetscInt max, PetscInt ind[], PetscInt tup[]) 156*20cf1dd8SToby Isaac { 157*20cf1dd8SToby Isaac PetscInt i; 158*20cf1dd8SToby Isaac PetscErrorCode ierr; 159*20cf1dd8SToby Isaac 160*20cf1dd8SToby Isaac PetscFunctionBegin; 161*20cf1dd8SToby Isaac if (len == 1) { 162*20cf1dd8SToby Isaac tup[0] = ind[0]++; 163*20cf1dd8SToby Isaac ind[0] = ind[0] >= max ? -1 : ind[0]; 164*20cf1dd8SToby Isaac } else if (max == 0) { 165*20cf1dd8SToby Isaac for (i = 0; i < len; ++i) {ind[0] = -1; tup[i] = 0;} 166*20cf1dd8SToby Isaac } else { 167*20cf1dd8SToby Isaac tup[0] = ind[0]; 168*20cf1dd8SToby Isaac ierr = TensorPoint_Internal(len-1, max, &ind[1], &tup[1]);CHKERRQ(ierr); 169*20cf1dd8SToby Isaac if (ind[1] < 0) { 170*20cf1dd8SToby Isaac ind[1] = 0; 171*20cf1dd8SToby Isaac if (ind[0] == max-1) {ind[0] = -1;} 172*20cf1dd8SToby Isaac else {++ind[0];} 173*20cf1dd8SToby Isaac } 174*20cf1dd8SToby Isaac } 175*20cf1dd8SToby Isaac PetscFunctionReturn(0); 176*20cf1dd8SToby Isaac } 177*20cf1dd8SToby Isaac 178*20cf1dd8SToby Isaac /* 179*20cf1dd8SToby Isaac p in [0, npoints), i in [0, pdim), c in [0, Nc) 180*20cf1dd8SToby Isaac 181*20cf1dd8SToby Isaac B[p][i][c] = B[p][i_scalar][c][c] 182*20cf1dd8SToby Isaac */ 183*20cf1dd8SToby Isaac PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) 184*20cf1dd8SToby Isaac { 185*20cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 186*20cf1dd8SToby Isaac DM dm = sp->dm; 187*20cf1dd8SToby Isaac PetscInt Nc = sp->Nc; 188*20cf1dd8SToby Isaac PetscInt ndegree = sp->degree+1; 189*20cf1dd8SToby Isaac PetscInt *degrees = poly->degrees; 190*20cf1dd8SToby Isaac PetscInt dim = sp->Nv; 191*20cf1dd8SToby Isaac PetscReal *lpoints, *tmp, *LB, *LD, *LH; 192*20cf1dd8SToby Isaac PetscInt *ind, *tup; 193*20cf1dd8SToby Isaac PetscInt c, pdim, d, e, der, der2, i, p, deg, o; 194*20cf1dd8SToby Isaac PetscErrorCode ierr; 195*20cf1dd8SToby Isaac 196*20cf1dd8SToby Isaac PetscFunctionBegin; 197*20cf1dd8SToby Isaac ierr = PetscSpaceGetDimension(sp, &pdim);CHKERRQ(ierr); 198*20cf1dd8SToby Isaac pdim /= Nc; 199*20cf1dd8SToby Isaac ierr = DMGetWorkArray(dm, npoints, MPIU_REAL, &lpoints);CHKERRQ(ierr); 200*20cf1dd8SToby Isaac ierr = DMGetWorkArray(dm, npoints*ndegree*3, MPIU_REAL, &tmp);CHKERRQ(ierr); 201*20cf1dd8SToby Isaac if (B) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LB);CHKERRQ(ierr);} 202*20cf1dd8SToby Isaac if (D) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LD);CHKERRQ(ierr);} 203*20cf1dd8SToby Isaac if (H) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LH);CHKERRQ(ierr);} 204*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 205*20cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 206*20cf1dd8SToby Isaac lpoints[p] = points[p*dim+d]; 207*20cf1dd8SToby Isaac } 208*20cf1dd8SToby Isaac ierr = PetscDTLegendreEval(npoints, lpoints, ndegree, degrees, tmp, &tmp[1*npoints*ndegree], &tmp[2*npoints*ndegree]);CHKERRQ(ierr); 209*20cf1dd8SToby Isaac /* LB, LD, LH (ndegree * dim x npoints) */ 210*20cf1dd8SToby Isaac for (deg = 0; deg < ndegree; ++deg) { 211*20cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 212*20cf1dd8SToby Isaac if (B) LB[(deg*dim + d)*npoints + p] = tmp[(0*npoints + p)*ndegree+deg]; 213*20cf1dd8SToby Isaac if (D) LD[(deg*dim + d)*npoints + p] = tmp[(1*npoints + p)*ndegree+deg]; 214*20cf1dd8SToby Isaac if (H) LH[(deg*dim + d)*npoints + p] = tmp[(2*npoints + p)*ndegree+deg]; 215*20cf1dd8SToby Isaac } 216*20cf1dd8SToby Isaac } 217*20cf1dd8SToby Isaac } 218*20cf1dd8SToby Isaac /* Multiply by A (pdim x ndegree * dim) */ 219*20cf1dd8SToby Isaac ierr = PetscMalloc2(dim,&ind,dim,&tup);CHKERRQ(ierr); 220*20cf1dd8SToby Isaac if (B) { 221*20cf1dd8SToby Isaac /* B (npoints x pdim x Nc) */ 222*20cf1dd8SToby Isaac ierr = PetscMemzero(B, npoints*pdim*Nc*Nc * sizeof(PetscReal));CHKERRQ(ierr); 223*20cf1dd8SToby Isaac if (poly->tensor) { 224*20cf1dd8SToby Isaac i = 0; 225*20cf1dd8SToby Isaac ierr = PetscMemzero(ind, dim * sizeof(PetscInt));CHKERRQ(ierr); 226*20cf1dd8SToby Isaac while (ind[0] >= 0) { 227*20cf1dd8SToby Isaac ierr = TensorPoint_Internal(dim, sp->degree+1, ind, tup);CHKERRQ(ierr); 228*20cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 229*20cf1dd8SToby Isaac B[(p*pdim + i)*Nc*Nc] = 1.0; 230*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 231*20cf1dd8SToby Isaac B[(p*pdim + i)*Nc*Nc] *= LB[(tup[d]*dim + d)*npoints + p]; 232*20cf1dd8SToby Isaac } 233*20cf1dd8SToby Isaac } 234*20cf1dd8SToby Isaac ++i; 235*20cf1dd8SToby Isaac } 236*20cf1dd8SToby Isaac } else { 237*20cf1dd8SToby Isaac i = 0; 238*20cf1dd8SToby Isaac for (o = 0; o <= sp->degree; ++o) { 239*20cf1dd8SToby Isaac ierr = PetscMemzero(ind, dim * sizeof(PetscInt));CHKERRQ(ierr); 240*20cf1dd8SToby Isaac while (ind[0] >= 0) { 241*20cf1dd8SToby Isaac ierr = LatticePoint_Internal(dim, o, ind, tup);CHKERRQ(ierr); 242*20cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 243*20cf1dd8SToby Isaac B[(p*pdim + i)*Nc*Nc] = 1.0; 244*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 245*20cf1dd8SToby Isaac B[(p*pdim + i)*Nc*Nc] *= LB[(tup[d]*dim + d)*npoints + p]; 246*20cf1dd8SToby Isaac } 247*20cf1dd8SToby Isaac } 248*20cf1dd8SToby Isaac ++i; 249*20cf1dd8SToby Isaac } 250*20cf1dd8SToby Isaac } 251*20cf1dd8SToby Isaac } 252*20cf1dd8SToby Isaac /* Make direct sum basis for multicomponent space */ 253*20cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 254*20cf1dd8SToby Isaac for (i = 0; i < pdim; ++i) { 255*20cf1dd8SToby Isaac for (c = 1; c < Nc; ++c) { 256*20cf1dd8SToby Isaac B[(p*pdim*Nc + i*Nc + c)*Nc + c] = B[(p*pdim + i)*Nc*Nc]; 257*20cf1dd8SToby Isaac } 258*20cf1dd8SToby Isaac } 259*20cf1dd8SToby Isaac } 260*20cf1dd8SToby Isaac } 261*20cf1dd8SToby Isaac if (D) { 262*20cf1dd8SToby Isaac /* D (npoints x pdim x Nc x dim) */ 263*20cf1dd8SToby Isaac ierr = PetscMemzero(D, npoints*pdim*Nc*Nc*dim * sizeof(PetscReal));CHKERRQ(ierr); 264*20cf1dd8SToby Isaac if (poly->tensor) { 265*20cf1dd8SToby Isaac i = 0; 266*20cf1dd8SToby Isaac ierr = PetscMemzero(ind, dim * sizeof(PetscInt));CHKERRQ(ierr); 267*20cf1dd8SToby Isaac while (ind[0] >= 0) { 268*20cf1dd8SToby Isaac ierr = TensorPoint_Internal(dim, sp->degree+1, ind, tup);CHKERRQ(ierr); 269*20cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 270*20cf1dd8SToby Isaac for (der = 0; der < dim; ++der) { 271*20cf1dd8SToby Isaac D[(p*pdim + i)*Nc*Nc*dim + der] = 1.0; 272*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 273*20cf1dd8SToby Isaac if (d == der) { 274*20cf1dd8SToby Isaac D[(p*pdim + i)*Nc*Nc*dim + der] *= LD[(tup[d]*dim + d)*npoints + p]; 275*20cf1dd8SToby Isaac } else { 276*20cf1dd8SToby Isaac D[(p*pdim + i)*Nc*Nc*dim + der] *= LB[(tup[d]*dim + d)*npoints + p]; 277*20cf1dd8SToby Isaac } 278*20cf1dd8SToby Isaac } 279*20cf1dd8SToby Isaac } 280*20cf1dd8SToby Isaac } 281*20cf1dd8SToby Isaac ++i; 282*20cf1dd8SToby Isaac } 283*20cf1dd8SToby Isaac } else { 284*20cf1dd8SToby Isaac i = 0; 285*20cf1dd8SToby Isaac for (o = 0; o <= sp->degree; ++o) { 286*20cf1dd8SToby Isaac ierr = PetscMemzero(ind, dim * sizeof(PetscInt));CHKERRQ(ierr); 287*20cf1dd8SToby Isaac while (ind[0] >= 0) { 288*20cf1dd8SToby Isaac ierr = LatticePoint_Internal(dim, o, ind, tup);CHKERRQ(ierr); 289*20cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 290*20cf1dd8SToby Isaac for (der = 0; der < dim; ++der) { 291*20cf1dd8SToby Isaac D[(p*pdim + i)*Nc*Nc*dim + der] = 1.0; 292*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 293*20cf1dd8SToby Isaac if (d == der) { 294*20cf1dd8SToby Isaac D[(p*pdim + i)*Nc*Nc*dim + der] *= LD[(tup[d]*dim + d)*npoints + p]; 295*20cf1dd8SToby Isaac } else { 296*20cf1dd8SToby Isaac D[(p*pdim + i)*Nc*Nc*dim + der] *= LB[(tup[d]*dim + d)*npoints + p]; 297*20cf1dd8SToby Isaac } 298*20cf1dd8SToby Isaac } 299*20cf1dd8SToby Isaac } 300*20cf1dd8SToby Isaac } 301*20cf1dd8SToby Isaac ++i; 302*20cf1dd8SToby Isaac } 303*20cf1dd8SToby Isaac } 304*20cf1dd8SToby Isaac } 305*20cf1dd8SToby Isaac /* Make direct sum basis for multicomponent space */ 306*20cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 307*20cf1dd8SToby Isaac for (i = 0; i < pdim; ++i) { 308*20cf1dd8SToby Isaac for (c = 1; c < Nc; ++c) { 309*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 310*20cf1dd8SToby Isaac D[((p*pdim*Nc + i*Nc + c)*Nc + c)*dim + d] = D[(p*pdim + i)*Nc*Nc*dim + d]; 311*20cf1dd8SToby Isaac } 312*20cf1dd8SToby Isaac } 313*20cf1dd8SToby Isaac } 314*20cf1dd8SToby Isaac } 315*20cf1dd8SToby Isaac } 316*20cf1dd8SToby Isaac if (H) { 317*20cf1dd8SToby Isaac /* H (npoints x pdim x Nc x Nc x dim x dim) */ 318*20cf1dd8SToby Isaac ierr = PetscMemzero(H, npoints*pdim*Nc*Nc*dim*dim * sizeof(PetscReal));CHKERRQ(ierr); 319*20cf1dd8SToby Isaac if (poly->tensor) { 320*20cf1dd8SToby Isaac i = 0; 321*20cf1dd8SToby Isaac ierr = PetscMemzero(ind, dim * sizeof(PetscInt));CHKERRQ(ierr); 322*20cf1dd8SToby Isaac while (ind[0] >= 0) { 323*20cf1dd8SToby Isaac ierr = TensorPoint_Internal(dim, sp->degree+1, ind, tup);CHKERRQ(ierr); 324*20cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 325*20cf1dd8SToby Isaac for (der = 0; der < dim; ++der) { 326*20cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der] = 1.0; 327*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 328*20cf1dd8SToby Isaac if (d == der) { 329*20cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der] *= LH[(tup[d]*dim + d)*npoints + p]; 330*20cf1dd8SToby Isaac } else { 331*20cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der] *= LB[(tup[d]*dim + d)*npoints + p]; 332*20cf1dd8SToby Isaac } 333*20cf1dd8SToby Isaac } 334*20cf1dd8SToby Isaac for (der2 = der + 1; der2 < dim; ++der2) { 335*20cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] = 1.0; 336*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 337*20cf1dd8SToby Isaac if (d == der || d == der2) { 338*20cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] *= LD[(tup[d]*dim + d)*npoints + p]; 339*20cf1dd8SToby Isaac } else { 340*20cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] *= LB[(tup[d]*dim + d)*npoints + p]; 341*20cf1dd8SToby Isaac } 342*20cf1dd8SToby Isaac } 343*20cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der2) * dim + der] = H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2]; 344*20cf1dd8SToby Isaac } 345*20cf1dd8SToby Isaac } 346*20cf1dd8SToby Isaac } 347*20cf1dd8SToby Isaac ++i; 348*20cf1dd8SToby Isaac } 349*20cf1dd8SToby Isaac } else { 350*20cf1dd8SToby Isaac i = 0; 351*20cf1dd8SToby Isaac for (o = 0; o <= sp->degree; ++o) { 352*20cf1dd8SToby Isaac ierr = PetscMemzero(ind, dim * sizeof(PetscInt));CHKERRQ(ierr); 353*20cf1dd8SToby Isaac while (ind[0] >= 0) { 354*20cf1dd8SToby Isaac ierr = LatticePoint_Internal(dim, o, ind, tup);CHKERRQ(ierr); 355*20cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 356*20cf1dd8SToby Isaac for (der = 0; der < dim; ++der) { 357*20cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der)*dim + der] = 1.0; 358*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 359*20cf1dd8SToby Isaac if (d == der) { 360*20cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der)*dim + der] *= LH[(tup[d]*dim + d)*npoints + p]; 361*20cf1dd8SToby Isaac } else { 362*20cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der)*dim + der] *= LB[(tup[d]*dim + d)*npoints + p]; 363*20cf1dd8SToby Isaac } 364*20cf1dd8SToby Isaac } 365*20cf1dd8SToby Isaac for (der2 = der + 1; der2 < dim; ++der2) { 366*20cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] = 1.0; 367*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 368*20cf1dd8SToby Isaac if (d == der || d == der2) { 369*20cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] *= LD[(tup[d]*dim + d)*npoints + p]; 370*20cf1dd8SToby Isaac } else { 371*20cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] *= LB[(tup[d]*dim + d)*npoints + p]; 372*20cf1dd8SToby Isaac } 373*20cf1dd8SToby Isaac } 374*20cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der2) * dim + der] = H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2]; 375*20cf1dd8SToby Isaac } 376*20cf1dd8SToby Isaac } 377*20cf1dd8SToby Isaac } 378*20cf1dd8SToby Isaac ++i; 379*20cf1dd8SToby Isaac } 380*20cf1dd8SToby Isaac } 381*20cf1dd8SToby Isaac } 382*20cf1dd8SToby Isaac /* Make direct sum basis for multicomponent space */ 383*20cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 384*20cf1dd8SToby Isaac for (i = 0; i < pdim; ++i) { 385*20cf1dd8SToby Isaac for (c = 1; c < Nc; ++c) { 386*20cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 387*20cf1dd8SToby Isaac for (e = 0; e < dim; ++e) { 388*20cf1dd8SToby Isaac H[(((p*pdim*Nc + i*Nc + c)*Nc + c)*dim + d)*dim + e] = H[((p*pdim + i)*Nc*Nc*dim + d)*dim + e]; 389*20cf1dd8SToby Isaac } 390*20cf1dd8SToby Isaac } 391*20cf1dd8SToby Isaac } 392*20cf1dd8SToby Isaac } 393*20cf1dd8SToby Isaac } 394*20cf1dd8SToby Isaac } 395*20cf1dd8SToby Isaac ierr = PetscFree2(ind,tup);CHKERRQ(ierr); 396*20cf1dd8SToby Isaac if (B) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LB);CHKERRQ(ierr);} 397*20cf1dd8SToby Isaac if (D) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LD);CHKERRQ(ierr);} 398*20cf1dd8SToby Isaac if (H) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LH);CHKERRQ(ierr);} 399*20cf1dd8SToby Isaac ierr = DMRestoreWorkArray(dm, npoints*ndegree*3, MPIU_REAL, &tmp);CHKERRQ(ierr); 400*20cf1dd8SToby Isaac ierr = DMRestoreWorkArray(dm, npoints, MPIU_REAL, &lpoints);CHKERRQ(ierr); 401*20cf1dd8SToby Isaac PetscFunctionReturn(0); 402*20cf1dd8SToby Isaac } 403*20cf1dd8SToby Isaac 404*20cf1dd8SToby Isaac /*@ 405*20cf1dd8SToby Isaac PetscSpacePolynomialSetTensor - Set whether a function space is a space of tensor polynomials (the space is spanned 406*20cf1dd8SToby Isaac by polynomials whose degree in each variabl is bounded by the given order), as opposed to polynomials (the space is 407*20cf1dd8SToby Isaac spanned by polynomials whose total degree---summing over all variables---is bounded by the given order). 408*20cf1dd8SToby Isaac 409*20cf1dd8SToby Isaac Input Parameters: 410*20cf1dd8SToby Isaac + sp - the function space object 411*20cf1dd8SToby Isaac - tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space 412*20cf1dd8SToby Isaac 413*20cf1dd8SToby Isaac Level: beginner 414*20cf1dd8SToby Isaac 415*20cf1dd8SToby Isaac .seealso: PetscSpacePolynomialGetTensor(), PetscSpaceSetDegree(), PetscSpaceSetNumVariables() 416*20cf1dd8SToby Isaac @*/ 417*20cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace sp, PetscBool tensor) 418*20cf1dd8SToby Isaac { 419*20cf1dd8SToby Isaac PetscErrorCode ierr; 420*20cf1dd8SToby Isaac 421*20cf1dd8SToby Isaac PetscFunctionBegin; 422*20cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 423*20cf1dd8SToby Isaac ierr = PetscTryMethod(sp,"PetscSpacePolynomialSetTensor_C",(PetscSpace,PetscBool),(sp,tensor));CHKERRQ(ierr); 424*20cf1dd8SToby Isaac PetscFunctionReturn(0); 425*20cf1dd8SToby Isaac } 426*20cf1dd8SToby Isaac 427*20cf1dd8SToby Isaac /*@ 428*20cf1dd8SToby Isaac PetscSpacePolynomialGetTensor - Get whether a function space is a space of tensor polynomials (the space is spanned 429*20cf1dd8SToby Isaac by polynomials whose degree in each variabl is bounded by the given order), as opposed to polynomials (the space is 430*20cf1dd8SToby Isaac spanned by polynomials whose total degree---summing over all variables---is bounded by the given order). 431*20cf1dd8SToby Isaac 432*20cf1dd8SToby Isaac Input Parameters: 433*20cf1dd8SToby Isaac . sp - the function space object 434*20cf1dd8SToby Isaac 435*20cf1dd8SToby Isaac Output Parameters: 436*20cf1dd8SToby Isaac . tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space 437*20cf1dd8SToby Isaac 438*20cf1dd8SToby Isaac Level: beginner 439*20cf1dd8SToby Isaac 440*20cf1dd8SToby Isaac .seealso: PetscSpacePolynomialSetTensor(), PetscSpaceSetDegree(), PetscSpaceSetNumVariables() 441*20cf1dd8SToby Isaac @*/ 442*20cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace sp, PetscBool *tensor) 443*20cf1dd8SToby Isaac { 444*20cf1dd8SToby Isaac PetscErrorCode ierr; 445*20cf1dd8SToby Isaac 446*20cf1dd8SToby Isaac PetscFunctionBegin; 447*20cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 448*20cf1dd8SToby Isaac PetscValidPointer(tensor, 2); 449*20cf1dd8SToby Isaac ierr = PetscTryMethod(sp,"PetscSpacePolynomialGetTensor_C",(PetscSpace,PetscBool*),(sp,tensor));CHKERRQ(ierr); 450*20cf1dd8SToby Isaac PetscFunctionReturn(0); 451*20cf1dd8SToby Isaac } 452*20cf1dd8SToby Isaac 453*20cf1dd8SToby Isaac static PetscErrorCode PetscSpacePolynomialSetTensor_Polynomial(PetscSpace sp, PetscBool tensor) 454*20cf1dd8SToby Isaac { 455*20cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 456*20cf1dd8SToby Isaac 457*20cf1dd8SToby Isaac PetscFunctionBegin; 458*20cf1dd8SToby Isaac poly->tensor = tensor; 459*20cf1dd8SToby Isaac PetscFunctionReturn(0); 460*20cf1dd8SToby Isaac } 461*20cf1dd8SToby Isaac 462*20cf1dd8SToby Isaac static PetscErrorCode PetscSpacePolynomialGetTensor_Polynomial(PetscSpace sp, PetscBool *tensor) 463*20cf1dd8SToby Isaac { 464*20cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 465*20cf1dd8SToby Isaac 466*20cf1dd8SToby Isaac PetscFunctionBegin; 467*20cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 468*20cf1dd8SToby Isaac PetscValidPointer(tensor, 2); 469*20cf1dd8SToby Isaac *tensor = poly->tensor; 470*20cf1dd8SToby Isaac PetscFunctionReturn(0); 471*20cf1dd8SToby Isaac } 472*20cf1dd8SToby Isaac 473*20cf1dd8SToby Isaac static PetscErrorCode PetscSpaceGetHeightSubspace_Polynomial(PetscSpace sp, PetscInt height, PetscSpace *subsp) 474*20cf1dd8SToby Isaac { 475*20cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 476*20cf1dd8SToby Isaac PetscInt Nc, dim, order; 477*20cf1dd8SToby Isaac PetscBool tensor; 478*20cf1dd8SToby Isaac PetscErrorCode ierr; 479*20cf1dd8SToby Isaac 480*20cf1dd8SToby Isaac PetscFunctionBegin; 481*20cf1dd8SToby Isaac ierr = PetscSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr); 482*20cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(sp, &dim);CHKERRQ(ierr); 483*20cf1dd8SToby Isaac ierr = PetscSpaceGetDegree(sp, &order, NULL);CHKERRQ(ierr); 484*20cf1dd8SToby Isaac ierr = PetscSpacePolynomialGetTensor(sp, &tensor);CHKERRQ(ierr); 485*20cf1dd8SToby Isaac if (height > dim || height < 0) {SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Asked for space at height %D for dimension %D space", height, dim);} 486*20cf1dd8SToby Isaac if (!poly->subspaces) {ierr = PetscCalloc1(dim, &poly->subspaces);CHKERRQ(ierr);} 487*20cf1dd8SToby Isaac if (height <= dim) { 488*20cf1dd8SToby Isaac if (!poly->subspaces[height-1]) { 489*20cf1dd8SToby Isaac PetscSpace sub; 490*20cf1dd8SToby Isaac 491*20cf1dd8SToby Isaac ierr = PetscSpaceCreate(PetscObjectComm((PetscObject) sp), &sub);CHKERRQ(ierr); 492*20cf1dd8SToby Isaac ierr = PetscSpaceSetNumComponents(sub, Nc);CHKERRQ(ierr); 493*20cf1dd8SToby Isaac ierr = PetscSpaceSetDegree(sub, order);CHKERRQ(ierr); 494*20cf1dd8SToby Isaac ierr = PetscSpaceSetType(sub, PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr); 495*20cf1dd8SToby Isaac ierr = PetscSpaceSetNumVariables(sub, dim-height);CHKERRQ(ierr); 496*20cf1dd8SToby Isaac ierr = PetscSpacePolynomialSetTensor(sub, tensor);CHKERRQ(ierr); 497*20cf1dd8SToby Isaac ierr = PetscSpaceSetUp(sub);CHKERRQ(ierr); 498*20cf1dd8SToby Isaac poly->subspaces[height-1] = sub; 499*20cf1dd8SToby Isaac } 500*20cf1dd8SToby Isaac *subsp = poly->subspaces[height-1]; 501*20cf1dd8SToby Isaac } else { 502*20cf1dd8SToby Isaac *subsp = NULL; 503*20cf1dd8SToby Isaac } 504*20cf1dd8SToby Isaac PetscFunctionReturn(0); 505*20cf1dd8SToby Isaac } 506*20cf1dd8SToby Isaac 507*20cf1dd8SToby Isaac PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp) 508*20cf1dd8SToby Isaac { 509*20cf1dd8SToby Isaac PetscErrorCode ierr; 510*20cf1dd8SToby Isaac 511*20cf1dd8SToby Isaac PetscFunctionBegin; 512*20cf1dd8SToby Isaac sp->ops->setfromoptions = PetscSpaceSetFromOptions_Polynomial; 513*20cf1dd8SToby Isaac sp->ops->setup = PetscSpaceSetUp_Polynomial; 514*20cf1dd8SToby Isaac sp->ops->view = PetscSpaceView_Polynomial; 515*20cf1dd8SToby Isaac sp->ops->destroy = PetscSpaceDestroy_Polynomial; 516*20cf1dd8SToby Isaac sp->ops->getdimension = PetscSpaceGetDimension_Polynomial; 517*20cf1dd8SToby Isaac sp->ops->evaluate = PetscSpaceEvaluate_Polynomial; 518*20cf1dd8SToby Isaac sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Polynomial; 519*20cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Polynomial);CHKERRQ(ierr); 520*20cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", PetscSpacePolynomialSetTensor_Polynomial);CHKERRQ(ierr); 521*20cf1dd8SToby Isaac PetscFunctionReturn(0); 522*20cf1dd8SToby Isaac } 523*20cf1dd8SToby Isaac 524*20cf1dd8SToby Isaac /*MC 525*20cf1dd8SToby Isaac PETSCSPACEPOLYNOMIAL = "poly" - A PetscSpace object that encapsulates a polynomial space, e.g. P1 is the space of 526*20cf1dd8SToby Isaac linear polynomials. The space is replicated for each component. 527*20cf1dd8SToby Isaac 528*20cf1dd8SToby Isaac Level: intermediate 529*20cf1dd8SToby Isaac 530*20cf1dd8SToby Isaac .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType() 531*20cf1dd8SToby Isaac M*/ 532*20cf1dd8SToby Isaac 533*20cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp) 534*20cf1dd8SToby Isaac { 535*20cf1dd8SToby Isaac PetscSpace_Poly *poly; 536*20cf1dd8SToby Isaac PetscErrorCode ierr; 537*20cf1dd8SToby Isaac 538*20cf1dd8SToby Isaac PetscFunctionBegin; 539*20cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 540*20cf1dd8SToby Isaac ierr = PetscNewLog(sp,&poly);CHKERRQ(ierr); 541*20cf1dd8SToby Isaac sp->data = poly; 542*20cf1dd8SToby Isaac 543*20cf1dd8SToby Isaac poly->symmetric = PETSC_FALSE; 544*20cf1dd8SToby Isaac poly->tensor = PETSC_FALSE; 545*20cf1dd8SToby Isaac poly->degrees = NULL; 546*20cf1dd8SToby Isaac poly->subspaces = NULL; 547*20cf1dd8SToby Isaac 548*20cf1dd8SToby Isaac ierr = PetscSpaceInitialize_Polynomial(sp);CHKERRQ(ierr); 549*20cf1dd8SToby Isaac PetscFunctionReturn(0); 550*20cf1dd8SToby Isaac } 551*20cf1dd8SToby Isaac 552*20cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace sp, PetscBool sym) 553*20cf1dd8SToby Isaac { 554*20cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 555*20cf1dd8SToby Isaac 556*20cf1dd8SToby Isaac PetscFunctionBegin; 557*20cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 558*20cf1dd8SToby Isaac poly->symmetric = sym; 559*20cf1dd8SToby Isaac PetscFunctionReturn(0); 560*20cf1dd8SToby Isaac } 561*20cf1dd8SToby Isaac 562*20cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialGetSymmetric(PetscSpace sp, PetscBool *sym) 563*20cf1dd8SToby Isaac { 564*20cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 565*20cf1dd8SToby Isaac 566*20cf1dd8SToby Isaac PetscFunctionBegin; 567*20cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 568*20cf1dd8SToby Isaac PetscValidPointer(sym, 2); 569*20cf1dd8SToby Isaac *sym = poly->symmetric; 570*20cf1dd8SToby Isaac PetscFunctionReturn(0); 571*20cf1dd8SToby Isaac } 572*20cf1dd8SToby Isaac 573