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 PetscOptionsHeadBegin(PetscOptionsObject,"PetscSpace polynomial options"); 9 PetscCall(PetscOptionsInt("-petscspace_ptrimmed_form_degree", "form degree of trimmed space", "PetscSpacePTrimmedSetFormDegree", pt->formDegree, &(pt->formDegree), NULL)); 10 PetscOptionsHeadEnd(); 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 %" PetscInt_FMT "%s-forms of degree %" PetscInt_FMT " (P-%" PetscInt_FMT "/\\%" PetscInt_FMT ")\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 PetscCheck(pt->formDegree >= -sp->Nv && pt->formDegree <= sp->Nv,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Form degree %" PetscInt_FMT " not in valid range [%" PetscInt_FMT ",%" PetscInt_FMT "]", 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 PetscCheck(sp->Nc % Nf == 0,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_INCOMP, "Number of components %" PetscInt_FMT " is not a multiple of form dimension %" PetscInt_FMT, 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 PetscCall(PetscObjectSetName((PetscObject)subsp, "sum component")); 91 PetscCall(PetscSpaceSetType(subsp, PETSCSPACEPTRIMMED)); 92 PetscCall(PetscSpaceGetNumVariables(sp, &Nv)); 93 PetscCall(PetscSpaceSetNumVariables(subsp, Nv)); 94 PetscCall(PetscSpaceSetNumComponents(subsp, Nf)); 95 PetscCall(PetscSpaceGetDegree(sp, °, &maxDeg)); 96 PetscCall(PetscSpaceSetDegree(subsp, deg, maxDeg)); 97 PetscCall(PetscSpacePTrimmedSetFormDegree(subsp, formDegree)); 98 PetscCall(PetscSpaceSetUp(subsp)); 99 for (PetscInt i = 0; i < nCopies; i++) { 100 PetscCall(PetscSpaceSumSetSubspace(sp, i, subsp)); 101 } 102 PetscCall(PetscSpaceDestroy(&subsp)); 103 PetscCall(PetscSpaceSetUp(sp)); 104 PetscFunctionReturn(0); 105 } 106 if (sp->degree == PETSC_DEFAULT) { 107 sp->degree = 0; 108 } else if (sp->degree < 0) { 109 SETERRQ(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Invalid negative degree %" PetscInt_FMT, sp->degree); 110 } 111 sp->maxDegree = (pt->formDegree == 0 || PetscAbsInt(pt->formDegree) == sp->Nv) ? sp->degree : sp->degree + 1; 112 if (pt->formDegree == 0 || PetscAbsInt(pt->formDegree) == sp->Nv) { 113 // Convert to regular polynomial space 114 PetscCall(PetscSpaceSetType(sp, PETSCSPACEPOLYNOMIAL)); 115 PetscCall(PetscSpaceSetUp(sp)); 116 PetscFunctionReturn(0); 117 } 118 pt->setupCalled = PETSC_TRUE; 119 PetscFunctionReturn(0); 120 } 121 122 static PetscErrorCode PetscSpaceGetDimension_Ptrimmed(PetscSpace sp, PetscInt *dim) 123 { 124 PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 125 PetscInt f; 126 PetscInt Nf; 127 128 PetscFunctionBegin; 129 f = pt->formDegree; 130 // For PetscSpace, degree refers to the largest complete polynomial degree contained in the space which 131 // is equal to the index of a P trimmed space only for 0-forms: otherwise, the index is degree + 1 132 PetscCall(PetscDTPTrimmedSize(sp->Nv, f == 0 ? sp->degree : sp->degree + 1, pt->formDegree, dim)); 133 PetscCall(PetscDTBinomialInt(sp->Nv, PetscAbsInt(pt->formDegree), &Nf)); 134 *dim *= (sp->Nc / Nf); 135 PetscFunctionReturn(0); 136 } 137 138 /* 139 p in [0, npoints), i in [0, pdim), c in [0, Nc) 140 141 B[p][i][c] = B[p][i_scalar][c][c] 142 */ 143 static PetscErrorCode PetscSpaceEvaluate_Ptrimmed(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) 144 { 145 PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 146 DM dm = sp->dm; 147 PetscInt jet, degree, Nf, Ncopies, Njet; 148 PetscInt Nc = sp->Nc; 149 PetscInt f; 150 PetscInt dim = sp->Nv; 151 PetscReal *eval; 152 PetscInt Nb; 153 154 PetscFunctionBegin; 155 if (!pt->setupCalled) { 156 PetscCall(PetscSpaceSetUp(sp)); 157 PetscCall(PetscSpaceEvaluate(sp, npoints, points, B, D, H)); 158 PetscFunctionReturn(0); 159 } 160 if (H) { 161 jet = 2; 162 } else if (D) { 163 jet = 1; 164 } else { 165 jet = 0; 166 } 167 f = pt->formDegree; 168 degree = f == 0 ? sp->degree : sp->degree + 1; 169 PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(f), &Nf)); 170 Ncopies = Nc / Nf; 171 PetscCheck(Ncopies == 1,PetscObjectComm((PetscObject) sp), PETSC_ERR_PLIB, "Multicopy spaces should have been converted to PETSCSPACESUM"); 172 PetscCall(PetscDTBinomialInt(dim + jet, dim, &Njet)); 173 PetscCall(PetscDTPTrimmedSize(dim, degree, f, &Nb)); 174 PetscCall(DMGetWorkArray(dm, Nb * Nf * Njet * npoints, MPIU_REAL, &eval)); 175 PetscCall(PetscDTPTrimmedEvalJet(dim, npoints, points, degree, f, jet, eval)); 176 if (B) { 177 PetscInt p_strl = Nf*Nb; 178 PetscInt b_strl = Nf; 179 PetscInt v_strl = 1; 180 181 PetscInt b_strr = Nf*Njet*npoints; 182 PetscInt v_strr = Njet*npoints; 183 PetscInt p_strr = 1; 184 185 for (PetscInt v = 0; v < Nf; v++) { 186 for (PetscInt b = 0; b < Nb; b++) { 187 for (PetscInt p = 0; p < npoints; p++) { 188 B[p*p_strl + b*b_strl + v*v_strl] = eval[b*b_strr + v*v_strr + p*p_strr]; 189 } 190 } 191 } 192 } 193 if (D) { 194 PetscInt p_strl = dim*Nf*Nb; 195 PetscInt b_strl = dim*Nf; 196 PetscInt v_strl = dim; 197 PetscInt d_strl = 1; 198 199 PetscInt b_strr = Nf*Njet*npoints; 200 PetscInt v_strr = Njet*npoints; 201 PetscInt d_strr = npoints; 202 PetscInt p_strr = 1; 203 204 for (PetscInt v = 0; v < Nf; v++) { 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 + v*v_strl + d*d_strl] = eval[b*b_strr + v*v_strr + (1+d)*d_strr + p*p_strr]; 209 } 210 } 211 } 212 } 213 } 214 if (H) { 215 PetscInt p_strl = dim*dim*Nf*Nb; 216 PetscInt b_strl = dim*dim*Nf; 217 PetscInt v_strl = dim*dim; 218 PetscInt d1_strl = dim; 219 PetscInt d2_strl = 1; 220 221 PetscInt b_strr = Nf*Njet*npoints; 222 PetscInt v_strr = Njet*npoints; 223 PetscInt j_strr = npoints; 224 PetscInt p_strr = 1; 225 226 PetscInt *derivs; 227 PetscCall(PetscCalloc1(dim, &derivs)); 228 for (PetscInt d1 = 0; d1 < dim; d1++) { 229 for (PetscInt d2 = 0; d2 < dim; d2++) { 230 PetscInt j; 231 derivs[d1]++; 232 derivs[d2]++; 233 PetscCall(PetscDTGradedOrderToIndex(dim, derivs, &j)); 234 derivs[d1]--; 235 derivs[d2]--; 236 for (PetscInt v = 0; v < Nf; v++) { 237 for (PetscInt b = 0; b < Nb; b++) { 238 for (PetscInt p = 0; p < npoints; p++) { 239 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]; 240 } 241 } 242 } 243 } 244 } 245 PetscCall(PetscFree(derivs)); 246 } 247 PetscCall(DMRestoreWorkArray(dm, Nb * Nf * Njet * npoints, MPIU_REAL, &eval)); 248 PetscFunctionReturn(0); 249 } 250 251 /*@ 252 PetscSpacePTrimmedSetFormDegree - Set the form degree of the trimmed polynomials. 253 254 Input Parameters: 255 + sp - the function space object 256 - formDegree - the form degree 257 258 Options Database: 259 . -petscspace_ptrimmed_form_degree <int> - The trimmed polynomial form degree 260 261 Level: intermediate 262 263 .seealso: `PetscDTAltV`, `PetscDTPTrimmedEvalJet()`, `PetscSpacePTrimmedGetFormDegree()` 264 @*/ 265 PetscErrorCode PetscSpacePTrimmedSetFormDegree(PetscSpace sp, PetscInt formDegree) 266 { 267 PetscFunctionBegin; 268 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 269 PetscTryMethod(sp,"PetscSpacePTrimmedSetFormDegree_C",(PetscSpace,PetscInt),(sp,formDegree)); 270 PetscFunctionReturn(0); 271 } 272 273 /*@ 274 PetscSpacePTrimmedGetFormDegree - Get the form degree of the trimmed polynomials. 275 276 Input Parameters: 277 . sp - the function space object 278 279 Output Parameters: 280 . formDegee - the form degree 281 282 Level: intermediate 283 284 .seealso: `PetscDTAltV`, `PetscDTPTrimmedEvalJet()`, `PetscSpacePTrimmedSetFormDegree()` 285 @*/ 286 PetscErrorCode PetscSpacePTrimmedGetFormDegree(PetscSpace sp, PetscInt *formDegree) 287 { 288 PetscFunctionBegin; 289 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 290 PetscValidIntPointer(formDegree, 2); 291 PetscTryMethod(sp,"PetscSpacePTrimmedGetFormDegree_C",(PetscSpace,PetscInt*),(sp,formDegree)); 292 PetscFunctionReturn(0); 293 } 294 295 static PetscErrorCode PetscSpacePTrimmedSetFormDegree_Ptrimmed(PetscSpace sp, PetscInt formDegree) 296 { 297 PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 298 299 PetscFunctionBegin; 300 pt->formDegree = formDegree; 301 PetscFunctionReturn(0); 302 } 303 304 static PetscErrorCode PetscSpacePTrimmedGetFormDegree_Ptrimmed(PetscSpace sp, PetscInt *formDegree) 305 { 306 PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 307 308 PetscFunctionBegin; 309 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 310 PetscValidIntPointer(formDegree, 2); 311 *formDegree = pt->formDegree; 312 PetscFunctionReturn(0); 313 } 314 315 static PetscErrorCode PetscSpaceGetHeightSubspace_Ptrimmed(PetscSpace sp, PetscInt height, PetscSpace *subsp) 316 { 317 PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 318 PetscInt dim; 319 320 PetscFunctionBegin; 321 PetscCall(PetscSpaceGetNumVariables(sp, &dim)); 322 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); 323 if (!pt->subspaces) PetscCall(PetscCalloc1(dim, &(pt->subspaces))); 324 if ((dim - height) <= PetscAbsInt(pt->formDegree)) { 325 if (!pt->subspaces[height-1]) { 326 PetscInt Nc, degree, Nf, Ncopies, Nfsub; 327 PetscSpace sub; 328 const char *name; 329 330 PetscCall(PetscSpaceGetNumComponents(sp, &Nc)); 331 PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(pt->formDegree), &Nf)); 332 PetscCall(PetscDTBinomialInt((dim-height), PetscAbsInt(pt->formDegree), &Nfsub)); 333 Ncopies = Nf / Nc; 334 PetscCall(PetscSpaceGetDegree(sp, °ree, NULL)); 335 336 PetscCall(PetscSpaceCreate(PetscObjectComm((PetscObject) sp), &sub)); 337 PetscCall(PetscObjectGetName((PetscObject) sp, &name)); 338 PetscCall(PetscObjectSetName((PetscObject) sub, name)); 339 PetscCall(PetscSpaceSetType(sub, PETSCSPACEPTRIMMED)); 340 PetscCall(PetscSpaceSetNumComponents(sub, Nfsub * Ncopies)); 341 PetscCall(PetscSpaceSetDegree(sub, degree, PETSC_DETERMINE)); 342 PetscCall(PetscSpaceSetNumVariables(sub, dim-height)); 343 PetscCall(PetscSpacePTrimmedSetFormDegree(sub, pt->formDegree)); 344 PetscCall(PetscSpaceSetUp(sub)); 345 pt->subspaces[height-1] = sub; 346 } 347 *subsp = pt->subspaces[height-1]; 348 } else { 349 *subsp = NULL; 350 } 351 PetscFunctionReturn(0); 352 } 353 354 static PetscErrorCode PetscSpaceInitialize_Ptrimmed(PetscSpace sp) 355 { 356 PetscFunctionBegin; 357 PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePTrimmedGetFormDegree_C", PetscSpacePTrimmedGetFormDegree_Ptrimmed)); 358 PetscCall(PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePTrimmedSetFormDegree_C", PetscSpacePTrimmedSetFormDegree_Ptrimmed)); 359 sp->ops->setfromoptions = PetscSpaceSetFromOptions_Ptrimmed; 360 sp->ops->setup = PetscSpaceSetUp_Ptrimmed; 361 sp->ops->view = PetscSpaceView_Ptrimmed; 362 sp->ops->destroy = PetscSpaceDestroy_Ptrimmed; 363 sp->ops->getdimension = PetscSpaceGetDimension_Ptrimmed; 364 sp->ops->evaluate = PetscSpaceEvaluate_Ptrimmed; 365 sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Ptrimmed; 366 PetscFunctionReturn(0); 367 } 368 369 /*MC 370 PETSCSPACEPTRIMMED = "ptrimmed" - A PetscSpace object that encapsulates a trimmed polynomial space. 371 372 Level: intermediate 373 374 .seealso: `PetscSpaceType`, `PetscSpaceCreate()`, `PetscSpaceSetType()`, `PetscDTPTrimmedEvalJet()` 375 M*/ 376 377 PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Ptrimmed(PetscSpace sp) 378 { 379 PetscSpace_Ptrimmed *pt; 380 381 PetscFunctionBegin; 382 PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 383 PetscCall(PetscNewLog(sp,&pt)); 384 sp->data = pt; 385 386 pt->subspaces = NULL; 387 sp->Nc = PETSC_DETERMINE; 388 389 PetscCall(PetscSpaceInitialize_Ptrimmed(sp)); 390 PetscFunctionReturn(0); 391 } 392