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