1 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 2 3 static PetscErrorCode PetscSpaceSetFromOptions_Ptrimmed(PetscOptionItems *PetscOptionsObject,PetscSpace sp) 4 { 5 PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 6 7 PetscFunctionBegin; 8 PetscCall(PetscOptionsHead(PetscOptionsObject,"PetscSpace polynomial options")); 9 PetscCall(PetscOptionsInt("-petscspace_ptrimmed_form_degree", "form degree of trimmed space", "PetscSpacePTrimmedSetFormDegree", pt->formDegree, &(pt->formDegree), NULL)); 10 PetscCall(PetscOptionsTail()); 11 PetscFunctionReturn(0); 12 } 13 14 static PetscErrorCode PetscSpacePTrimmedView_Ascii(PetscSpace sp, PetscViewer v) 15 { 16 PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 17 PetscInt f, tdegree; 18 19 PetscFunctionBegin; 20 f = pt->formDegree; 21 tdegree = f == 0 ? sp->degree : sp->degree + 1; 22 PetscCall(PetscViewerASCIIPrintf(v, "Trimmed polynomials %D%s-forms of degree %D (P-%D/\\%D)\n", PetscAbsInt(f), f < 0 ? "*" : "", sp->degree, tdegree, PetscAbsInt(f))); 23 PetscFunctionReturn(0); 24 } 25 26 static PetscErrorCode PetscSpaceView_Ptrimmed(PetscSpace sp, PetscViewer viewer) 27 { 28 PetscBool iascii; 29 30 PetscFunctionBegin; 31 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 32 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 33 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 34 if (iascii) PetscCall(PetscSpacePTrimmedView_Ascii(sp, viewer)); 35 PetscFunctionReturn(0); 36 } 37 38 static PetscErrorCode PetscSpaceDestroy_Ptrimmed(PetscSpace sp) 39 { 40 PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 41 42 PetscFunctionBegin; 43 PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePTrimmedGetFormDegree_C", NULL)); 44 PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePTrimmedSetFormDegree_C", NULL)); 45 if (pt->subspaces) { 46 PetscInt d; 47 48 for (d = 0; d < sp->Nv; ++d) { 49 PetscCall(PetscSpaceDestroy(&pt->subspaces[d])); 50 } 51 } 52 PetscCall(PetscFree(pt->subspaces)); 53 PetscCall(PetscFree(pt)); 54 PetscFunctionReturn(0); 55 } 56 57 static PetscErrorCode PetscSpaceSetUp_Ptrimmed(PetscSpace sp) 58 { 59 PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 60 PetscInt Nf; 61 62 PetscFunctionBegin; 63 if (pt->setupCalled) PetscFunctionReturn(0); 64 PetscCheckFalse(pt->formDegree < -sp->Nv || pt->formDegree > sp->Nv,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Form degree %D not in valid range [%D,%D]", pt->formDegree, sp->Nv, sp->Nv); 65 PetscCall(PetscDTBinomialInt(sp->Nv, PetscAbsInt(pt->formDegree), &Nf)); 66 if (sp->Nc == PETSC_DETERMINE) { 67 sp->Nc = Nf; 68 } 69 PetscCheckFalse(sp->Nc % Nf,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_INCOMP, "Number of components %D is not a multiple of form dimension %D", sp->Nc, Nf); 70 if (sp->Nc != Nf) { 71 PetscSpace subsp; 72 PetscInt nCopies = sp->Nc / Nf; 73 PetscInt Nv, deg, maxDeg; 74 PetscInt formDegree = pt->formDegree; 75 const char *prefix; 76 const char *name; 77 char subname[PETSC_MAX_PATH_LEN]; 78 79 PetscCall(PetscSpaceSetType(sp, PETSCSPACESUM)); 80 PetscCall(PetscSpaceSumSetConcatenate(sp, PETSC_TRUE)); 81 PetscCall(PetscSpaceSumSetNumSubspaces(sp, nCopies)); 82 PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject)sp), &subsp)); 83 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)sp, &prefix)); 84 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)subsp, prefix)); 85 PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)subsp, "sumcomp_")); 86 if (((PetscObject)sp)->name) { 87 PetscCall(PetscObjectGetName((PetscObject)sp, &name)); 88 PetscCall(PetscSNPrintf(subname, PETSC_MAX_PATH_LEN-1, "%s sum component", name)); 89 PetscCall(PetscObjectSetName((PetscObject)subsp, subname)); 90 } else { 91 PetscCall(PetscObjectSetName((PetscObject)subsp, "sum component")); 92 } 93 PetscCall(PetscSpaceSetType(subsp, PETSCSPACEPTRIMMED)); 94 PetscCall(PetscSpaceGetNumVariables(sp, &Nv)); 95 PetscCall(PetscSpaceSetNumVariables(subsp, Nv)); 96 PetscCall(PetscSpaceSetNumComponents(subsp, Nf)); 97 PetscCall(PetscSpaceGetDegree(sp, °, &maxDeg)); 98 PetscCall(PetscSpaceSetDegree(subsp, deg, maxDeg)); 99 PetscCall(PetscSpacePTrimmedSetFormDegree(subsp, formDegree)); 100 PetscCall(PetscSpaceSetUp(subsp)); 101 for (PetscInt i = 0; i < nCopies; i++) { 102 PetscCall(PetscSpaceSumSetSubspace(sp, i, subsp)); 103 } 104 PetscCall(PetscSpaceDestroy(&subsp)); 105 PetscCall(PetscSpaceSetUp(sp)); 106 PetscFunctionReturn(0); 107 } 108 if (sp->degree == PETSC_DEFAULT) { 109 sp->degree = 0; 110 } else if (sp->degree < 0) { 111 SETERRQ(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Invalid negative degree %D", sp->degree); 112 } 113 sp->maxDegree = (pt->formDegree == 0 || PetscAbsInt(pt->formDegree) == sp->Nv) ? sp->degree : sp->degree + 1; 114 if (pt->formDegree == 0 || PetscAbsInt(pt->formDegree) == sp->Nv) { 115 // Convert to regular polynomial space 116 PetscCall(PetscSpaceSetType(sp, PETSCSPACEPOLYNOMIAL)); 117 PetscCall(PetscSpaceSetUp(sp)); 118 PetscFunctionReturn(0); 119 } 120 pt->setupCalled = PETSC_TRUE; 121 PetscFunctionReturn(0); 122 } 123 124 static PetscErrorCode PetscSpaceGetDimension_Ptrimmed(PetscSpace sp, PetscInt *dim) 125 { 126 PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 127 PetscInt f; 128 PetscInt Nf; 129 130 PetscFunctionBegin; 131 f = pt->formDegree; 132 // For PetscSpace, degree refers to the largest complete polynomial degree contained in the space which 133 // is equal to the index of a P trimmed space only for 0-forms: otherwise, the index is degree + 1 134 PetscCall(PetscDTPTrimmedSize(sp->Nv, f == 0 ? sp->degree : sp->degree + 1, pt->formDegree, dim)); 135 PetscCall(PetscDTBinomialInt(sp->Nv, PetscAbsInt(pt->formDegree), &Nf)); 136 *dim *= (sp->Nc / Nf); 137 PetscFunctionReturn(0); 138 } 139 140 /* 141 p in [0, npoints), i in [0, pdim), c in [0, Nc) 142 143 B[p][i][c] = B[p][i_scalar][c][c] 144 */ 145 static PetscErrorCode PetscSpaceEvaluate_Ptrimmed(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) 146 { 147 PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 148 DM dm = sp->dm; 149 PetscInt jet, degree, Nf, Ncopies, Njet; 150 PetscInt Nc = sp->Nc; 151 PetscInt f; 152 PetscInt dim = sp->Nv; 153 PetscReal *eval; 154 PetscInt Nb; 155 156 PetscFunctionBegin; 157 if (!pt->setupCalled) { 158 PetscCall(PetscSpaceSetUp(sp)); 159 PetscCall(PetscSpaceEvaluate(sp, npoints, points, B, D, H)); 160 PetscFunctionReturn(0); 161 } 162 if (H) { 163 jet = 2; 164 } else if (D) { 165 jet = 1; 166 } else { 167 jet = 0; 168 } 169 f = pt->formDegree; 170 degree = f == 0 ? sp->degree : sp->degree + 1; 171 PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(f), &Nf)); 172 Ncopies = Nc / Nf; 173 PetscCheckFalse(Ncopies != 1,PetscObjectComm((PetscObject) sp), PETSC_ERR_PLIB, "Multicopy spaces should have been converted to PETSCSPACESUM"); 174 PetscCall(PetscDTBinomialInt(dim + jet, dim, &Njet)); 175 PetscCall(PetscDTPTrimmedSize(dim, degree, f, &Nb)); 176 PetscCall(DMGetWorkArray(dm, Nb * Nf * Njet * npoints, MPIU_REAL, &eval)); 177 PetscCall(PetscDTPTrimmedEvalJet(dim, npoints, points, degree, f, jet, eval)); 178 if (B) { 179 PetscInt p_strl = Nf*Nb; 180 PetscInt b_strl = Nf; 181 PetscInt v_strl = 1; 182 183 PetscInt b_strr = Nf*Njet*npoints; 184 PetscInt v_strr = Njet*npoints; 185 PetscInt p_strr = 1; 186 187 for (PetscInt v = 0; v < Nf; v++) { 188 for (PetscInt b = 0; b < Nb; b++) { 189 for (PetscInt p = 0; p < npoints; p++) { 190 B[p*p_strl + b*b_strl + v*v_strl] = eval[b*b_strr + v*v_strr + p*p_strr]; 191 } 192 } 193 } 194 } 195 if (D) { 196 PetscInt p_strl = dim*Nf*Nb; 197 PetscInt b_strl = dim*Nf; 198 PetscInt v_strl = dim; 199 PetscInt d_strl = 1; 200 201 PetscInt b_strr = Nf*Njet*npoints; 202 PetscInt v_strr = Njet*npoints; 203 PetscInt d_strr = npoints; 204 PetscInt p_strr = 1; 205 206 for (PetscInt v = 0; v < Nf; v++) { 207 for (PetscInt d = 0; d < dim; d++) { 208 for (PetscInt b = 0; b < Nb; b++) { 209 for (PetscInt p = 0; p < npoints; p++) { 210 D[p*p_strl + b*b_strl + v*v_strl + d*d_strl] = eval[b*b_strr + v*v_strr + (1+d)*d_strr + p*p_strr]; 211 } 212 } 213 } 214 } 215 } 216 if (H) { 217 PetscInt p_strl = dim*dim*Nf*Nb; 218 PetscInt b_strl = dim*dim*Nf; 219 PetscInt v_strl = dim*dim; 220 PetscInt d1_strl = dim; 221 PetscInt d2_strl = 1; 222 223 PetscInt b_strr = Nf*Njet*npoints; 224 PetscInt v_strr = Njet*npoints; 225 PetscInt j_strr = npoints; 226 PetscInt p_strr = 1; 227 228 PetscInt *derivs; 229 PetscCall(PetscCalloc1(dim, &derivs)); 230 for (PetscInt d1 = 0; d1 < dim; d1++) { 231 for (PetscInt d2 = 0; d2 < dim; d2++) { 232 PetscInt j; 233 derivs[d1]++; 234 derivs[d2]++; 235 PetscCall(PetscDTGradedOrderToIndex(dim, derivs, &j)); 236 derivs[d1]--; 237 derivs[d2]--; 238 for (PetscInt v = 0; v < Nf; v++) { 239 for (PetscInt b = 0; b < Nb; b++) { 240 for (PetscInt p = 0; p < npoints; p++) { 241 H[p*p_strl + b*b_strl + v*v_strl + d1*d1_strl + d2*d2_strl] = eval[b*b_strr + v*v_strr + j*j_strr + p*p_strr]; 242 } 243 } 244 } 245 } 246 } 247 PetscCall(PetscFree(derivs)); 248 } 249 PetscCall(DMRestoreWorkArray(dm, Nb * Nf * Njet * npoints, MPIU_REAL, &eval)); 250 PetscFunctionReturn(0); 251 } 252 253 /*@ 254 PetscSpacePTrimmedSetFormDegree - Set the form degree of the trimmed polynomials. 255 256 Input Parameters: 257 + sp - the function space object 258 - formDegree - the form degree 259 260 Options Database: 261 . -petscspace_ptrimmed_form_degree <int> - The trimmed polynomial form degree 262 263 Level: intermediate 264 265 .seealso: PetscDTAltV, PetscDTPTrimmedEvalJet(), PetscSpacePTrimmedGetFormDegree() 266 @*/ 267 PetscErrorCode PetscSpacePTrimmedSetFormDegree(PetscSpace sp, PetscInt formDegree) 268 { 269 PetscFunctionBegin; 270 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 271 PetscTryMethod(sp,"PetscSpacePTrimmedSetFormDegree_C",(PetscSpace,PetscInt),(sp,formDegree)); 272 PetscFunctionReturn(0); 273 } 274 275 /*@ 276 PetscSpacePTrimmedGetFormDegree - Get the form degree of the trimmed polynomials. 277 278 Input Parameters: 279 . sp - the function space object 280 281 Output Parameters: 282 . formDegee - the form degree 283 284 Level: intermediate 285 286 .seealso: PetscDTAltV, PetscDTPTrimmedEvalJet(), PetscSpacePTrimmedSetFormDegree() 287 @*/ 288 PetscErrorCode PetscSpacePTrimmedGetFormDegree(PetscSpace sp, PetscInt *formDegree) 289 { 290 PetscFunctionBegin; 291 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 292 PetscValidIntPointer(formDegree, 2); 293 PetscTryMethod(sp,"PetscSpacePTrimmedGetFormDegree_C",(PetscSpace,PetscInt*),(sp,formDegree)); 294 PetscFunctionReturn(0); 295 } 296 297 static PetscErrorCode PetscSpacePTrimmedSetFormDegree_Ptrimmed(PetscSpace sp, PetscInt formDegree) 298 { 299 PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 300 301 PetscFunctionBegin; 302 pt->formDegree = formDegree; 303 PetscFunctionReturn(0); 304 } 305 306 static PetscErrorCode PetscSpacePTrimmedGetFormDegree_Ptrimmed(PetscSpace sp, PetscInt *formDegree) 307 { 308 PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 309 310 PetscFunctionBegin; 311 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 312 PetscValidIntPointer(formDegree, 2); 313 *formDegree = pt->formDegree; 314 PetscFunctionReturn(0); 315 } 316 317 static PetscErrorCode PetscSpaceGetHeightSubspace_Ptrimmed(PetscSpace sp, PetscInt height, PetscSpace *subsp) 318 { 319 PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 320 PetscInt dim; 321 322 PetscFunctionBegin; 323 PetscCall(PetscSpaceGetNumVariables(sp, &dim)); 324 PetscCheckFalse(height > dim || height < 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Asked for space at height %D for dimension %D space", height, dim); 325 if (!pt->subspaces) PetscCall(PetscCalloc1(dim, &(pt->subspaces))); 326 if ((dim - height) <= PetscAbsInt(pt->formDegree)) { 327 if (!pt->subspaces[height-1]) { 328 PetscInt Nc, degree, Nf, Ncopies, Nfsub; 329 PetscSpace sub; 330 const char *name; 331 332 PetscCall(PetscSpaceGetNumComponents(sp, &Nc)); 333 PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(pt->formDegree), &Nf)); 334 PetscCall(PetscDTBinomialInt((dim-height), PetscAbsInt(pt->formDegree), &Nfsub)); 335 Ncopies = Nf / Nc; 336 PetscCall(PetscSpaceGetDegree(sp, °ree, NULL)); 337 338 PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject) sp), &sub)); 339 PetscCall(PetscObjectGetName((PetscObject) sp, &name)); 340 PetscCall(PetscObjectSetName((PetscObject) sub, name)); 341 PetscCall(PetscSpaceSetType(sub, PETSCSPACEPTRIMMED)); 342 PetscCall(PetscSpaceSetNumComponents(sub, Nfsub * Ncopies)); 343 PetscCall(PetscSpaceSetDegree(sub, degree, PETSC_DETERMINE)); 344 PetscCall(PetscSpaceSetNumVariables(sub, dim-height)); 345 PetscCall(PetscSpacePTrimmedSetFormDegree(sub, pt->formDegree)); 346 PetscCall(PetscSpaceSetUp(sub)); 347 pt->subspaces[height-1] = sub; 348 } 349 *subsp = pt->subspaces[height-1]; 350 } else { 351 *subsp = NULL; 352 } 353 PetscFunctionReturn(0); 354 } 355 356 static PetscErrorCode PetscSpaceInitialize_Ptrimmed(PetscSpace sp) 357 { 358 PetscFunctionBegin; 359 PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePTrimmedGetFormDegree_C", PetscSpacePTrimmedGetFormDegree_Ptrimmed)); 360 PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePTrimmedSetFormDegree_C", PetscSpacePTrimmedSetFormDegree_Ptrimmed)); 361 sp->ops->setfromoptions = PetscSpaceSetFromOptions_Ptrimmed; 362 sp->ops->setup = PetscSpaceSetUp_Ptrimmed; 363 sp->ops->view = PetscSpaceView_Ptrimmed; 364 sp->ops->destroy = PetscSpaceDestroy_Ptrimmed; 365 sp->ops->getdimension = PetscSpaceGetDimension_Ptrimmed; 366 sp->ops->evaluate = PetscSpaceEvaluate_Ptrimmed; 367 sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Ptrimmed; 368 PetscFunctionReturn(0); 369 } 370 371 /*MC 372 PETSCSPACEPTRIMMED = "ptrimmed" - A PetscSpace object that encapsulates a trimmed polynomial space. 373 374 Level: intermediate 375 376 .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType(), PetscDTPTrimmedEvalJet() 377 M*/ 378 379 PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Ptrimmed(PetscSpace sp) 380 { 381 PetscSpace_Ptrimmed *pt; 382 383 PetscFunctionBegin; 384 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 385 PetscCall(PetscNewLog(sp,&pt)); 386 sp->data = pt; 387 388 pt->subspaces = NULL; 389 sp->Nc = PETSC_DETERMINE; 390 391 PetscCall(PetscSpaceInitialize_Ptrimmed(sp)); 392 PetscFunctionReturn(0); 393 } 394