1 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 2 3 static PetscErrorCode PetscSpaceSetFromOptions_Polynomial(PetscOptionItems *PetscOptionsObject,PetscSpace sp) 4 { 5 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 6 7 PetscFunctionBegin; 8 PetscCall(PetscOptionsHead(PetscOptionsObject,"PetscSpace polynomial options")); 9 PetscCall(PetscOptionsBool("-petscspace_poly_tensor", "Use the tensor product polynomials", "PetscSpacePolynomialSetTensor", poly->tensor, &poly->tensor, NULL)); 10 PetscCall(PetscOptionsTail()); 11 PetscFunctionReturn(0); 12 } 13 14 static PetscErrorCode PetscSpacePolynomialView_Ascii(PetscSpace sp, PetscViewer v) 15 { 16 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 17 18 PetscFunctionBegin; 19 PetscCall(PetscViewerASCIIPrintf(v, "%s space of degree %D\n", poly->tensor ? "Tensor polynomial" : "Polynomial", sp->degree)); 20 PetscFunctionReturn(0); 21 } 22 23 static PetscErrorCode PetscSpaceView_Polynomial(PetscSpace sp, PetscViewer viewer) 24 { 25 PetscBool iascii; 26 27 PetscFunctionBegin; 28 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 29 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 30 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 31 if (iascii) PetscCall(PetscSpacePolynomialView_Ascii(sp, viewer)); 32 PetscFunctionReturn(0); 33 } 34 35 static PetscErrorCode PetscSpaceDestroy_Polynomial(PetscSpace sp) 36 { 37 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 38 39 PetscFunctionBegin; 40 PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", NULL)); 41 PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", NULL)); 42 if (poly->subspaces) { 43 PetscInt d; 44 45 for (d = 0; d < sp->Nv; ++d) { 46 PetscCall(PetscSpaceDestroy(&poly->subspaces[d])); 47 } 48 } 49 PetscCall(PetscFree(poly->subspaces)); 50 PetscCall(PetscFree(poly)); 51 PetscFunctionReturn(0); 52 } 53 54 static PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp) 55 { 56 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 57 58 PetscFunctionBegin; 59 if (poly->setupCalled) PetscFunctionReturn(0); 60 if (sp->Nv <= 1) { 61 poly->tensor = PETSC_FALSE; 62 } 63 if (sp->Nc != 1) { 64 PetscInt Nc = sp->Nc; 65 PetscBool tensor = poly->tensor; 66 PetscInt Nv = sp->Nv; 67 PetscInt degree = sp->degree; 68 const char *prefix; 69 const char *name; 70 char subname[PETSC_MAX_PATH_LEN]; 71 PetscSpace subsp; 72 73 PetscCall(PetscSpaceSetType(sp, PETSCSPACESUM)); 74 PetscCall(PetscSpaceSumSetNumSubspaces(sp, Nc)); 75 PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)sp), &subsp)); 76 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sp, &prefix)); 77 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)subsp, prefix)); 78 PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)subsp, "sumcomp_")); 79 if (((PetscObject)sp)->name) { 80 PetscCall(PetscObjectGetName((PetscObject)sp, &name)); 81 PetscCall(PetscSNPrintf(subname, PETSC_MAX_PATH_LEN-1, "%s sum component", name)); 82 PetscCall(PetscObjectSetName((PetscObject)subsp, subname)); 83 } else { 84 PetscCall(PetscObjectSetName((PetscObject)subsp, "sum component")); 85 } 86 PetscCall(PetscSpaceSetType(subsp, PETSCSPACEPOLYNOMIAL)); 87 PetscCall(PetscSpaceSetDegree(subsp, degree, PETSC_DETERMINE)); 88 PetscCall(PetscSpaceSetNumComponents(subsp, 1)); 89 PetscCall(PetscSpaceSetNumVariables(subsp, Nv)); 90 PetscCall(PetscSpacePolynomialSetTensor(subsp, tensor)); 91 PetscCall(PetscSpaceSetUp(subsp)); 92 for (PetscInt i = 0; i < Nc; i++) { 93 PetscCall(PetscSpaceSumSetSubspace(sp, i, subsp)); 94 } 95 PetscCall(PetscSpaceDestroy(&subsp)); 96 PetscCall(PetscSpaceSetUp(sp)); 97 PetscFunctionReturn(0); 98 } 99 if (poly->tensor) { 100 sp->maxDegree = PETSC_DETERMINE; 101 PetscCall(PetscSpaceSetType(sp, PETSCSPACETENSOR)); 102 PetscCall(PetscSpaceSetUp(sp)); 103 PetscFunctionReturn(0); 104 } 105 PetscCheckFalse(sp->degree < 0,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Negative degree %D invalid", sp->degree); 106 sp->maxDegree = sp->degree; 107 poly->setupCalled = PETSC_TRUE; 108 PetscFunctionReturn(0); 109 } 110 111 static PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim) 112 { 113 PetscInt deg = sp->degree; 114 PetscInt n = sp->Nv; 115 116 PetscFunctionBegin; 117 PetscCall(PetscDTBinomialInt(n + deg, n, dim)); 118 *dim *= sp->Nc; 119 PetscFunctionReturn(0); 120 } 121 122 static PetscErrorCode CoordinateBasis(PetscInt dim, PetscInt npoints, const PetscReal points[], PetscInt jet, PetscInt Njet, PetscReal pScalar[]) 123 { 124 PetscFunctionBegin; 125 PetscCall(PetscArrayzero(pScalar, (1 + dim) * Njet * npoints)); 126 for (PetscInt b = 0; b < 1 + dim; b++) { 127 for (PetscInt j = 0; j < PetscMin(1 + dim, Njet); j++) { 128 if (j == 0) { 129 if (b == 0) { 130 for (PetscInt pt = 0; pt < npoints; pt++) { 131 pScalar[b * Njet * npoints + j * npoints + pt] = 1.; 132 } 133 } else { 134 for (PetscInt pt = 0; pt < npoints; pt++) { 135 pScalar[b * Njet * npoints + j * npoints + pt] = points[pt * dim + (b-1)]; 136 } 137 } 138 } else if (j == b) { 139 for (PetscInt pt = 0; pt < npoints; pt++) { 140 pScalar[b * Njet * npoints + j * npoints + pt] = 1.; 141 } 142 } 143 } 144 } 145 PetscFunctionReturn(0); 146 } 147 148 static PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) 149 { 150 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 151 DM dm = sp->dm; 152 PetscInt dim = sp->Nv; 153 PetscInt Nb, jet, Njet; 154 PetscReal *pScalar; 155 156 PetscFunctionBegin; 157 if (!poly->setupCalled) { 158 PetscCall(PetscSpaceSetUp(sp)); 159 PetscCall(PetscSpaceEvaluate(sp, npoints, points, B, D, H)); 160 PetscFunctionReturn(0); 161 } 162 PetscCheckFalse(poly->tensor || sp->Nc != 1,PETSC_COMM_SELF, PETSC_ERR_PLIB, "tensor and multicomponent spaces should have been converted"); 163 PetscCall(PetscDTBinomialInt(dim + sp->degree, dim, &Nb)); 164 if (H) { 165 jet = 2; 166 } else if (D) { 167 jet = 1; 168 } else { 169 jet = 0; 170 } 171 PetscCall(PetscDTBinomialInt(dim + jet, dim, &Njet)); 172 PetscCall(DMGetWorkArray(dm, Nb * Njet * npoints, MPIU_REAL, &pScalar)); 173 // Why are we handling the case degree == 1 specially? Because we don't want numerical noise when we evaluate hat 174 // functions at the vertices of a simplex, which happens when we invert the Vandermonde matrix of the PKD basis. 175 // We don't make any promise about which basis is used. 176 if (sp->degree == 1) { 177 PetscCall(CoordinateBasis(dim, npoints, points, jet, Njet, pScalar)); 178 } else { 179 PetscCall(PetscDTPKDEvalJet(dim, npoints, points, sp->degree, jet, pScalar)); 180 } 181 if (B) { 182 PetscInt p_strl = Nb; 183 PetscInt b_strl = 1; 184 185 PetscInt b_strr = Njet*npoints; 186 PetscInt p_strr = 1; 187 188 PetscCall(PetscArrayzero(B, npoints * Nb)); 189 for (PetscInt b = 0; b < Nb; b++) { 190 for (PetscInt p = 0; p < npoints; p++) { 191 B[p*p_strl + b*b_strl] = pScalar[b*b_strr + p*p_strr]; 192 } 193 } 194 } 195 if (D) { 196 PetscInt p_strl = dim*Nb; 197 PetscInt b_strl = dim; 198 PetscInt d_strl = 1; 199 200 PetscInt b_strr = Njet*npoints; 201 PetscInt d_strr = npoints; 202 PetscInt p_strr = 1; 203 204 PetscCall(PetscArrayzero(D, npoints * Nb * dim)); 205 for (PetscInt d = 0; d < dim; d++) { 206 for (PetscInt b = 0; b < Nb; b++) { 207 for (PetscInt p = 0; p < npoints; p++) { 208 D[p*p_strl + b*b_strl + d*d_strl] = pScalar[b*b_strr + (1+d)*d_strr + p*p_strr]; 209 } 210 } 211 } 212 } 213 if (H) { 214 PetscInt p_strl = dim*dim*Nb; 215 PetscInt b_strl = dim*dim; 216 PetscInt d1_strl = dim; 217 PetscInt d2_strl = 1; 218 219 PetscInt b_strr = Njet*npoints; 220 PetscInt j_strr = npoints; 221 PetscInt p_strr = 1; 222 223 PetscInt *derivs; 224 PetscCall(PetscCalloc1(dim, &derivs)); 225 PetscCall(PetscArrayzero(H, npoints * Nb * dim * dim)); 226 for (PetscInt d1 = 0; d1 < dim; d1++) { 227 for (PetscInt d2 = 0; d2 < dim; d2++) { 228 PetscInt j; 229 derivs[d1]++; 230 derivs[d2]++; 231 PetscCall(PetscDTGradedOrderToIndex(dim, derivs, &j)); 232 derivs[d1]--; 233 derivs[d2]--; 234 for (PetscInt b = 0; b < Nb; b++) { 235 for (PetscInt p = 0; p < npoints; p++) { 236 H[p*p_strl + b*b_strl + d1*d1_strl + d2*d2_strl] = pScalar[b*b_strr + j*j_strr + p*p_strr]; 237 } 238 } 239 } 240 } 241 PetscCall(PetscFree(derivs)); 242 } 243 PetscCall(DMRestoreWorkArray(dm, Nb * Njet * npoints, MPIU_REAL, &pScalar)); 244 PetscFunctionReturn(0); 245 } 246 247 /*@ 248 PetscSpacePolynomialSetTensor - Set whether a function space is a space of tensor polynomials (the space is spanned 249 by polynomials whose degree in each variable is bounded by the given order), as opposed to polynomials (the space is 250 spanned by polynomials whose total degree---summing over all variables---is bounded by the given order). 251 252 Input Parameters: 253 + sp - the function space object 254 - tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space 255 256 Options Database: 257 . -petscspace_poly_tensor <bool> - Whether to use tensor product polynomials in higher dimension 258 259 Level: intermediate 260 261 .seealso: PetscSpacePolynomialGetTensor(), PetscSpaceSetDegree(), PetscSpaceSetNumVariables() 262 @*/ 263 PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace sp, PetscBool tensor) 264 { 265 PetscFunctionBegin; 266 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 267 PetscTryMethod(sp,"PetscSpacePolynomialSetTensor_C",(PetscSpace,PetscBool),(sp,tensor)); 268 PetscFunctionReturn(0); 269 } 270 271 /*@ 272 PetscSpacePolynomialGetTensor - Get whether a function space is a space of tensor polynomials (the space is spanned 273 by polynomials whose degree in each variabl is bounded by the given order), as opposed to polynomials (the space is 274 spanned by polynomials whose total degree---summing over all variables---is bounded by the given order). 275 276 Input Parameters: 277 . sp - the function space object 278 279 Output Parameters: 280 . tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space 281 282 Level: intermediate 283 284 .seealso: PetscSpacePolynomialSetTensor(), PetscSpaceSetDegree(), PetscSpaceSetNumVariables() 285 @*/ 286 PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace sp, PetscBool *tensor) 287 { 288 PetscFunctionBegin; 289 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 290 PetscValidBoolPointer(tensor, 2); 291 PetscTryMethod(sp,"PetscSpacePolynomialGetTensor_C",(PetscSpace,PetscBool*),(sp,tensor)); 292 PetscFunctionReturn(0); 293 } 294 295 static PetscErrorCode PetscSpacePolynomialSetTensor_Polynomial(PetscSpace sp, PetscBool tensor) 296 { 297 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 298 299 PetscFunctionBegin; 300 poly->tensor = tensor; 301 PetscFunctionReturn(0); 302 } 303 304 static PetscErrorCode PetscSpacePolynomialGetTensor_Polynomial(PetscSpace sp, PetscBool *tensor) 305 { 306 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 307 308 PetscFunctionBegin; 309 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 310 PetscValidBoolPointer(tensor, 2); 311 *tensor = poly->tensor; 312 PetscFunctionReturn(0); 313 } 314 315 static PetscErrorCode PetscSpaceGetHeightSubspace_Polynomial(PetscSpace sp, PetscInt height, PetscSpace *subsp) 316 { 317 PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 318 PetscInt Nc, dim, order; 319 PetscBool tensor; 320 321 PetscFunctionBegin; 322 PetscCall(PetscSpaceGetNumComponents(sp, &Nc)); 323 PetscCall(PetscSpaceGetNumVariables(sp, &dim)); 324 PetscCall(PetscSpaceGetDegree(sp, &order, NULL)); 325 PetscCall(PetscSpacePolynomialGetTensor(sp, &tensor)); 326 PetscCheckFalse(height > dim || height < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Asked for space at height %D for dimension %D space", height, dim); 327 if (!poly->subspaces) PetscCall(PetscCalloc1(dim, &poly->subspaces)); 328 if (height <= dim) { 329 if (!poly->subspaces[height-1]) { 330 PetscSpace sub; 331 const char *name; 332 333 PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject) sp), &sub)); 334 PetscCall(PetscObjectGetName((PetscObject) sp, &name)); 335 PetscCall(PetscObjectSetName((PetscObject) sub, name)); 336 PetscCall(PetscSpaceSetType(sub, PETSCSPACEPOLYNOMIAL)); 337 PetscCall(PetscSpaceSetNumComponents(sub, Nc)); 338 PetscCall(PetscSpaceSetDegree(sub, order, PETSC_DETERMINE)); 339 PetscCall(PetscSpaceSetNumVariables(sub, dim-height)); 340 PetscCall(PetscSpacePolynomialSetTensor(sub, tensor)); 341 PetscCall(PetscSpaceSetUp(sub)); 342 poly->subspaces[height-1] = sub; 343 } 344 *subsp = poly->subspaces[height-1]; 345 } else { 346 *subsp = NULL; 347 } 348 PetscFunctionReturn(0); 349 } 350 351 static PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp) 352 { 353 PetscFunctionBegin; 354 sp->ops->setfromoptions = PetscSpaceSetFromOptions_Polynomial; 355 sp->ops->setup = PetscSpaceSetUp_Polynomial; 356 sp->ops->view = PetscSpaceView_Polynomial; 357 sp->ops->destroy = PetscSpaceDestroy_Polynomial; 358 sp->ops->getdimension = PetscSpaceGetDimension_Polynomial; 359 sp->ops->evaluate = PetscSpaceEvaluate_Polynomial; 360 sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Polynomial; 361 PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Polynomial)); 362 PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", PetscSpacePolynomialSetTensor_Polynomial)); 363 PetscFunctionReturn(0); 364 } 365 366 /*MC 367 PETSCSPACEPOLYNOMIAL = "poly" - A PetscSpace object that encapsulates a polynomial space, e.g. P1 is the space of 368 linear polynomials. The space is replicated for each component. 369 370 Level: intermediate 371 372 .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType() 373 M*/ 374 375 PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp) 376 { 377 PetscSpace_Poly *poly; 378 379 PetscFunctionBegin; 380 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 381 PetscCall(PetscNewLog(sp,&poly)); 382 sp->data = poly; 383 384 poly->tensor = PETSC_FALSE; 385 poly->subspaces = NULL; 386 387 PetscCall(PetscSpaceInitialize_Polynomial(sp)); 388 PetscFunctionReturn(0); 389 } 390