1130d5748SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 2130d5748SToby Isaac 3130d5748SToby Isaac static PetscErrorCode PetscSpaceSetFromOptions_Ptrimmed(PetscOptionItems *PetscOptionsObject,PetscSpace sp) 4130d5748SToby Isaac { 5130d5748SToby Isaac PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 6130d5748SToby Isaac PetscErrorCode ierr; 7130d5748SToby Isaac 8130d5748SToby Isaac PetscFunctionBegin; 9130d5748SToby Isaac ierr = PetscOptionsHead(PetscOptionsObject,"PetscSpace polynomial options");CHKERRQ(ierr); 10130d5748SToby Isaac ierr = PetscOptionsInt("-petscspace_ptrimmed_form_degree", "form degree of trimmed space", "PetscSpacePTrimmedSetFormDegree", pt->formDegree, &(pt->formDegree), NULL);CHKERRQ(ierr); 11130d5748SToby Isaac ierr = PetscOptionsTail();CHKERRQ(ierr); 12130d5748SToby Isaac PetscFunctionReturn(0); 13130d5748SToby Isaac } 14130d5748SToby Isaac 15130d5748SToby Isaac static PetscErrorCode PetscSpacePTrimmedView_Ascii(PetscSpace sp, PetscViewer v) 16130d5748SToby Isaac { 17130d5748SToby Isaac PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 18*ccb4e88aSToby Isaac PetscInt f, tdegree; 19130d5748SToby Isaac PetscErrorCode ierr; 20130d5748SToby Isaac 21130d5748SToby Isaac PetscFunctionBegin; 22130d5748SToby Isaac f = pt->formDegree; 23*ccb4e88aSToby Isaac tdegree = f == 0 ? sp->degree : sp->degree + 1; 24*ccb4e88aSToby Isaac ierr = PetscViewerASCIIPrintf(v, "Trimmed polynomials %D%s-forms of degree %D (P-%D/\\%D)\n", PetscAbsInt(f), f < 0 ? "*" : "", sp->degree, tdegree, PetscAbsInt(f));CHKERRQ(ierr); 25130d5748SToby Isaac PetscFunctionReturn(0); 26130d5748SToby Isaac } 27130d5748SToby Isaac 28130d5748SToby Isaac static PetscErrorCode PetscSpaceView_Ptrimmed(PetscSpace sp, PetscViewer viewer) 29130d5748SToby Isaac { 30130d5748SToby Isaac PetscBool iascii; 31130d5748SToby Isaac PetscErrorCode ierr; 32130d5748SToby Isaac 33130d5748SToby Isaac PetscFunctionBegin; 34130d5748SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 35130d5748SToby Isaac PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 36130d5748SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 37130d5748SToby Isaac if (iascii) {ierr = PetscSpacePTrimmedView_Ascii(sp, viewer);CHKERRQ(ierr);} 38130d5748SToby Isaac PetscFunctionReturn(0); 39130d5748SToby Isaac } 40130d5748SToby Isaac 41130d5748SToby Isaac static PetscErrorCode PetscSpaceDestroy_Ptrimmed(PetscSpace sp) 42130d5748SToby Isaac { 43130d5748SToby Isaac PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 44130d5748SToby Isaac PetscErrorCode ierr; 45130d5748SToby Isaac 46130d5748SToby Isaac PetscFunctionBegin; 47130d5748SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePTrimmedGetFormDegree_C", NULL);CHKERRQ(ierr); 48130d5748SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePTrimmedSetFormDegree_C", NULL);CHKERRQ(ierr); 49130d5748SToby Isaac if (pt->subspaces) { 50130d5748SToby Isaac PetscInt d; 51130d5748SToby Isaac 52130d5748SToby Isaac for (d = 0; d < sp->Nv; ++d) { 53130d5748SToby Isaac ierr = PetscSpaceDestroy(&pt->subspaces[d]);CHKERRQ(ierr); 54130d5748SToby Isaac } 55130d5748SToby Isaac } 56130d5748SToby Isaac ierr = PetscFree(pt->subspaces);CHKERRQ(ierr); 57130d5748SToby Isaac ierr = PetscFree(pt);CHKERRQ(ierr); 58130d5748SToby Isaac PetscFunctionReturn(0); 59130d5748SToby Isaac } 60130d5748SToby Isaac 61130d5748SToby Isaac static PetscErrorCode PetscSpaceSetUp_Ptrimmed(PetscSpace sp) 62130d5748SToby Isaac { 63130d5748SToby Isaac PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 64130d5748SToby Isaac PetscInt Nf; 65130d5748SToby Isaac PetscErrorCode ierr; 66130d5748SToby Isaac 67130d5748SToby Isaac PetscFunctionBegin; 68130d5748SToby Isaac if (pt->setupCalled) PetscFunctionReturn(0); 69130d5748SToby Isaac if (pt->formDegree < -sp->Nv || pt->formDegree > sp->Nv) SETERRQ3(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Form degree %D not in valid range [%D,%D]", pt->formDegree, sp->Nv, sp->Nv); 70130d5748SToby Isaac ierr = PetscDTBinomialInt(sp->Nv, PetscAbsInt(pt->formDegree), &Nf);CHKERRQ(ierr); 71130d5748SToby Isaac if (sp->Nc == PETSC_DETERMINE) { 72130d5748SToby Isaac sp->Nc = Nf; 73130d5748SToby Isaac } 74130d5748SToby Isaac if (sp->Nc % Nf) SETERRQ2(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_INCOMP, "Number of components %D is not a multiple of form dimension %D\n", sp->Nc, Nf); 75130d5748SToby Isaac if (sp->Nc != Nf) { 76130d5748SToby Isaac PetscSpace subsp; 77130d5748SToby Isaac PetscInt nCopies = sp->Nc / Nf; 78130d5748SToby Isaac PetscInt Nv, deg, maxDeg; 79130d5748SToby Isaac PetscInt formDegree = pt->formDegree; 80130d5748SToby Isaac const char *prefix; 81130d5748SToby Isaac const char *name; 82130d5748SToby Isaac char subname[PETSC_MAX_PATH_LEN]; 83130d5748SToby Isaac 84130d5748SToby Isaac ierr = PetscSpaceSetType(sp, PETSCSPACESUM);CHKERRQ(ierr); 85130d5748SToby Isaac ierr = PetscSpaceSumSetConcatenate(sp, PETSC_TRUE);CHKERRQ(ierr); 86130d5748SToby Isaac ierr = PetscSpaceSumSetNumSubspaces(sp, nCopies);CHKERRQ(ierr); 87130d5748SToby Isaac ierr = PetscSpaceCreate(PetscObjectComm((PetscObject)sp), &subsp);CHKERRQ(ierr); 88130d5748SToby Isaac ierr = PetscObjectGetOptionsPrefix((PetscObject)sp, &prefix);CHKERRQ(ierr); 89130d5748SToby Isaac ierr = PetscObjectSetOptionsPrefix((PetscObject)subsp, prefix);CHKERRQ(ierr); 90130d5748SToby Isaac ierr = PetscObjectAppendOptionsPrefix((PetscObject)subsp, "sumcomp_");CHKERRQ(ierr); 91130d5748SToby Isaac if (((PetscObject)sp)->name) { 92130d5748SToby Isaac ierr = PetscObjectGetName((PetscObject)sp, &name);CHKERRQ(ierr); 93130d5748SToby Isaac ierr = PetscSNPrintf(subname, PETSC_MAX_PATH_LEN-1, "%s sum component", name);CHKERRQ(ierr); 94130d5748SToby Isaac ierr = PetscObjectSetName((PetscObject)subsp, subname);CHKERRQ(ierr); 95130d5748SToby Isaac } else { 96130d5748SToby Isaac ierr = PetscObjectSetName((PetscObject)subsp, "sum component");CHKERRQ(ierr); 97130d5748SToby Isaac } 98130d5748SToby Isaac ierr = PetscSpaceSetType(subsp, PETSCSPACEPTRIMMED);CHKERRQ(ierr); 99130d5748SToby Isaac ierr = PetscSpaceGetNumVariables(sp, &Nv);CHKERRQ(ierr); 100130d5748SToby Isaac ierr = PetscSpaceSetNumVariables(subsp, Nv);CHKERRQ(ierr); 101130d5748SToby Isaac ierr = PetscSpaceSetNumComponents(subsp, Nf);CHKERRQ(ierr); 102130d5748SToby Isaac ierr = PetscSpaceGetDegree(sp, °, &maxDeg);CHKERRQ(ierr); 103130d5748SToby Isaac ierr = PetscSpaceSetDegree(subsp, deg, maxDeg);CHKERRQ(ierr); 104130d5748SToby Isaac ierr = PetscSpacePTrimmedSetFormDegree(subsp, formDegree);CHKERRQ(ierr); 105130d5748SToby Isaac ierr = PetscSpaceSetUp(subsp);CHKERRQ(ierr); 106130d5748SToby Isaac for (PetscInt i = 0; i < nCopies; i++) { 107130d5748SToby Isaac ierr = PetscSpaceSumSetSubspace(sp, i, subsp);CHKERRQ(ierr); 108130d5748SToby Isaac } 109130d5748SToby Isaac ierr = PetscSpaceDestroy(&subsp);CHKERRQ(ierr); 110130d5748SToby Isaac ierr = PetscSpaceSetUp(sp);CHKERRQ(ierr); 111130d5748SToby Isaac PetscFunctionReturn(0); 112130d5748SToby Isaac } 113130d5748SToby Isaac if (sp->degree == PETSC_DEFAULT) { 114130d5748SToby Isaac sp->degree = 0; 115130d5748SToby Isaac } else if (sp->degree < 0) { 116130d5748SToby Isaac SETERRQ1(PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_OUTOFRANGE, "Invalid negative degree %D\n", sp->degree); 117130d5748SToby Isaac } 118*ccb4e88aSToby Isaac sp->maxDegree = (pt->formDegree == 0 || PetscAbsInt(pt->formDegree) == sp->Nv) ? sp->degree : sp->degree + 1; 119130d5748SToby Isaac if (pt->formDegree == 0 || PetscAbsInt(pt->formDegree) == sp->Nv) { 120130d5748SToby Isaac // Convert to regular polynomial space 121130d5748SToby Isaac ierr = PetscSpaceSetType(sp, PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr); 122130d5748SToby Isaac ierr = PetscSpaceSetUp(sp);CHKERRQ(ierr); 123130d5748SToby Isaac PetscFunctionReturn(0); 124130d5748SToby Isaac } 125130d5748SToby Isaac pt->setupCalled = PETSC_TRUE; 126130d5748SToby Isaac PetscFunctionReturn(0); 127130d5748SToby Isaac } 128130d5748SToby Isaac 129130d5748SToby Isaac static PetscErrorCode PetscSpaceGetDimension_Ptrimmed(PetscSpace sp, PetscInt *dim) 130130d5748SToby Isaac { 131130d5748SToby Isaac PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 132130d5748SToby Isaac PetscInt f; 133130d5748SToby Isaac PetscInt Nf; 134130d5748SToby Isaac PetscErrorCode ierr; 135130d5748SToby Isaac 136130d5748SToby Isaac PetscFunctionBegin; 137130d5748SToby Isaac f = pt->formDegree; 138130d5748SToby Isaac // For PetscSpace, degree refers to the largest complete polynomial degree contained in the space which 139130d5748SToby Isaac // is equal to the index of a P trimmed space only for 0-forms: otherwise, the index is degree + 1 140130d5748SToby Isaac ierr = PetscDTPTrimmedSize(sp->Nv, f == 0 ? sp->degree : sp->degree + 1, pt->formDegree, dim);CHKERRQ(ierr); 141130d5748SToby Isaac ierr = PetscDTBinomialInt(sp->Nv, PetscAbsInt(pt->formDegree), &Nf);CHKERRQ(ierr); 142130d5748SToby Isaac *dim *= (sp->Nc / Nf); 143130d5748SToby Isaac PetscFunctionReturn(0); 144130d5748SToby Isaac } 145130d5748SToby Isaac 146130d5748SToby Isaac /* 147130d5748SToby Isaac p in [0, npoints), i in [0, pdim), c in [0, Nc) 148130d5748SToby Isaac 149130d5748SToby Isaac B[p][i][c] = B[p][i_scalar][c][c] 150130d5748SToby Isaac */ 151130d5748SToby Isaac static PetscErrorCode PetscSpaceEvaluate_Ptrimmed(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) 152130d5748SToby Isaac { 153130d5748SToby Isaac PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 154130d5748SToby Isaac DM dm = sp->dm; 155130d5748SToby Isaac PetscInt jet, degree, Nf, Ncopies, Njet; 156130d5748SToby Isaac PetscInt Nc = sp->Nc; 157130d5748SToby Isaac PetscInt f; 158130d5748SToby Isaac PetscInt dim = sp->Nv; 159130d5748SToby Isaac PetscReal *eval; 160130d5748SToby Isaac PetscInt Nb; 161130d5748SToby Isaac PetscErrorCode ierr; 162130d5748SToby Isaac 163130d5748SToby Isaac PetscFunctionBegin; 164130d5748SToby Isaac if (!pt->setupCalled) { 165130d5748SToby Isaac ierr = PetscSpaceSetUp(sp);CHKERRQ(ierr); 166130d5748SToby Isaac ierr = PetscSpaceEvaluate(sp, npoints, points, B, D, H);CHKERRQ(ierr); 167130d5748SToby Isaac PetscFunctionReturn(0); 168130d5748SToby Isaac } 169130d5748SToby Isaac if (H) { 170130d5748SToby Isaac jet = 2; 171130d5748SToby Isaac } else if (D) { 172130d5748SToby Isaac jet = 1; 173130d5748SToby Isaac } else { 174130d5748SToby Isaac jet = 0; 175130d5748SToby Isaac } 176130d5748SToby Isaac f = pt->formDegree; 177130d5748SToby Isaac degree = f == 0 ? sp->degree : sp->degree + 1; 178130d5748SToby Isaac ierr = PetscDTBinomialInt(dim, PetscAbsInt(f), &Nf);CHKERRQ(ierr); 179130d5748SToby Isaac Ncopies = Nc / Nf; 180130d5748SToby Isaac if (Ncopies != 1) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_PLIB, "Multicopy spaces should have been converted to PETSCSPACESUM"); 181130d5748SToby Isaac ierr = PetscDTBinomialInt(dim + jet, dim, &Njet);CHKERRQ(ierr); 182130d5748SToby Isaac ierr = PetscDTPTrimmedSize(dim, degree, f, &Nb);CHKERRQ(ierr); 183130d5748SToby Isaac ierr = DMGetWorkArray(dm, Nb * Nf * Njet * npoints, MPIU_REAL, &eval);CHKERRQ(ierr); 184130d5748SToby Isaac ierr = PetscDTPTrimmedEvalJet(dim, npoints, points, degree, f, jet, eval);CHKERRQ(ierr); 185130d5748SToby Isaac if (B) { 186130d5748SToby Isaac PetscInt p_strl = Nf*Nb; 187130d5748SToby Isaac PetscInt b_strl = Nf; 188130d5748SToby Isaac PetscInt v_strl = 1; 189130d5748SToby Isaac 190130d5748SToby Isaac PetscInt b_strr = Nf*Njet*npoints; 191130d5748SToby Isaac PetscInt v_strr = Njet*npoints; 192130d5748SToby Isaac PetscInt p_strr = 1; 193130d5748SToby Isaac 194130d5748SToby Isaac for (PetscInt v = 0; v < Nf; v++) { 195130d5748SToby Isaac for (PetscInt b = 0; b < Nb; b++) { 196130d5748SToby Isaac for (PetscInt p = 0; p < npoints; p++) { 197130d5748SToby Isaac B[p*p_strl + b*b_strl + v*v_strl] = eval[b*b_strr + v*v_strr + p*p_strr]; 198130d5748SToby Isaac } 199130d5748SToby Isaac } 200130d5748SToby Isaac } 201130d5748SToby Isaac } 202130d5748SToby Isaac if (D) { 203*ccb4e88aSToby Isaac PetscInt p_strl = dim*Nf*Nb; 204*ccb4e88aSToby Isaac PetscInt b_strl = dim*Nf; 205130d5748SToby Isaac PetscInt v_strl = dim; 206130d5748SToby Isaac PetscInt d_strl = 1; 207130d5748SToby Isaac 208130d5748SToby Isaac PetscInt b_strr = Nf*Njet*npoints; 209130d5748SToby Isaac PetscInt v_strr = Njet*npoints; 210130d5748SToby Isaac PetscInt d_strr = npoints; 211130d5748SToby Isaac PetscInt p_strr = 1; 212130d5748SToby Isaac 213130d5748SToby Isaac for (PetscInt v = 0; v < Nf; v++) { 214130d5748SToby Isaac for (PetscInt d = 0; d < dim; d++) { 215130d5748SToby Isaac for (PetscInt b = 0; b < Nb; b++) { 216130d5748SToby Isaac for (PetscInt p = 0; p < npoints; p++) { 217*ccb4e88aSToby Isaac 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]; 218130d5748SToby Isaac } 219130d5748SToby Isaac } 220130d5748SToby Isaac } 221130d5748SToby Isaac } 222130d5748SToby Isaac } 223130d5748SToby Isaac if (H) { 224*ccb4e88aSToby Isaac PetscInt p_strl = dim*dim*Nf*Nb; 225*ccb4e88aSToby Isaac PetscInt b_strl = dim*dim*Nf; 226*ccb4e88aSToby Isaac PetscInt v_strl = dim*dim; 227130d5748SToby Isaac PetscInt d1_strl = dim; 228130d5748SToby Isaac PetscInt d2_strl = 1; 229130d5748SToby Isaac 230130d5748SToby Isaac PetscInt b_strr = Nf*Njet*npoints; 231130d5748SToby Isaac PetscInt v_strr = Njet*npoints; 232130d5748SToby Isaac PetscInt j_strr = npoints; 233130d5748SToby Isaac PetscInt p_strr = 1; 234130d5748SToby Isaac 235130d5748SToby Isaac PetscInt *derivs; 236130d5748SToby Isaac ierr = PetscCalloc1(dim, &derivs);CHKERRQ(ierr); 237130d5748SToby Isaac for (PetscInt d1 = 0; d1 < dim; d1++) { 238130d5748SToby Isaac for (PetscInt d2 = 0; d2 < dim; d2++) { 239130d5748SToby Isaac PetscInt j; 240130d5748SToby Isaac derivs[d1]++; 241130d5748SToby Isaac derivs[d2]++; 242130d5748SToby Isaac ierr = PetscDTGradedOrderToIndex(dim, derivs, &j);CHKERRQ(ierr); 243130d5748SToby Isaac derivs[d1]--; 244130d5748SToby Isaac derivs[d2]--; 245130d5748SToby Isaac for (PetscInt v = 0; v < Nf; v++) { 246130d5748SToby Isaac for (PetscInt b = 0; b < Nb; b++) { 247130d5748SToby Isaac for (PetscInt p = 0; p < npoints; p++) { 248*ccb4e88aSToby Isaac 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]; 249130d5748SToby Isaac } 250130d5748SToby Isaac } 251130d5748SToby Isaac } 252130d5748SToby Isaac } 253130d5748SToby Isaac } 254130d5748SToby Isaac ierr = PetscFree(derivs);CHKERRQ(ierr); 255130d5748SToby Isaac } 256130d5748SToby Isaac ierr = DMRestoreWorkArray(dm, Nb * Nf * Njet * npoints, MPIU_REAL, &eval);CHKERRQ(ierr); 257130d5748SToby Isaac PetscFunctionReturn(0); 258130d5748SToby Isaac } 259130d5748SToby Isaac 260130d5748SToby Isaac /*@ 261130d5748SToby Isaac PetscSpacePTrimmedSetFormDegree - Set the form degree of the trimmed polynomials. 262130d5748SToby Isaac 263130d5748SToby Isaac Input Parameters: 264130d5748SToby Isaac + sp - the function space object 265130d5748SToby Isaac - formDegree - the form degree 266130d5748SToby Isaac 267130d5748SToby Isaac Options Database: 268130d5748SToby Isaac . -petscspace_ptrimmed_form_degree <int> - The trimmed polynomial form degree 269130d5748SToby Isaac 270130d5748SToby Isaac Level: intermediate 271130d5748SToby Isaac 272130d5748SToby Isaac .seealso: PetscDTAltV, PetscDTPTrimmedEvalJet(), PetscSpacePTrimmedGetFormDegree() 273130d5748SToby Isaac @*/ 274130d5748SToby Isaac PetscErrorCode PetscSpacePTrimmedSetFormDegree(PetscSpace sp, PetscInt formDegree) 275130d5748SToby Isaac { 276130d5748SToby Isaac PetscErrorCode ierr; 277130d5748SToby Isaac 278130d5748SToby Isaac PetscFunctionBegin; 279130d5748SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 280130d5748SToby Isaac ierr = PetscTryMethod(sp,"PetscSpacePTrimmedSetFormDegree_C",(PetscSpace,PetscInt),(sp,formDegree));CHKERRQ(ierr); 281130d5748SToby Isaac PetscFunctionReturn(0); 282130d5748SToby Isaac } 283130d5748SToby Isaac 284130d5748SToby Isaac /*@ 285130d5748SToby Isaac PetscSpacePTrimmedGetFormDegree - Get the form degree of the trimmed polynomials. 286130d5748SToby Isaac 287130d5748SToby Isaac Input Parameters: 288130d5748SToby Isaac . sp - the function space object 289130d5748SToby Isaac 290130d5748SToby Isaac Output Parameters: 291130d5748SToby Isaac . formDegee - the form degree 292130d5748SToby Isaac 293130d5748SToby Isaac Level: intermediate 294130d5748SToby Isaac 295130d5748SToby Isaac .seealso: PetscDTAltV, PetscDTPTrimmedEvalJet(), PetscSpacePTrimmedSetFormDegree() 296130d5748SToby Isaac @*/ 297130d5748SToby Isaac PetscErrorCode PetscSpacePTrimmedGetFormDegree(PetscSpace sp, PetscInt *formDegree) 298130d5748SToby Isaac { 299130d5748SToby Isaac PetscErrorCode ierr; 300130d5748SToby Isaac 301130d5748SToby Isaac PetscFunctionBegin; 302130d5748SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 303130d5748SToby Isaac PetscValidPointer(formDegree, 2); 304130d5748SToby Isaac ierr = PetscTryMethod(sp,"PetscSpacePTrimmedGetFormDegree_C",(PetscSpace,PetscInt*),(sp,formDegree));CHKERRQ(ierr); 305130d5748SToby Isaac PetscFunctionReturn(0); 306130d5748SToby Isaac } 307130d5748SToby Isaac 308130d5748SToby Isaac static PetscErrorCode PetscSpacePTrimmedSetFormDegree_Ptrimmed(PetscSpace sp, PetscInt formDegree) 309130d5748SToby Isaac { 310130d5748SToby Isaac PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 311130d5748SToby Isaac 312130d5748SToby Isaac PetscFunctionBegin; 313130d5748SToby Isaac pt->formDegree = formDegree; 314130d5748SToby Isaac PetscFunctionReturn(0); 315130d5748SToby Isaac } 316130d5748SToby Isaac 317130d5748SToby Isaac static PetscErrorCode PetscSpacePTrimmedGetFormDegree_Ptrimmed(PetscSpace sp, PetscInt *formDegree) 318130d5748SToby Isaac { 319130d5748SToby Isaac PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 320130d5748SToby Isaac 321130d5748SToby Isaac PetscFunctionBegin; 322130d5748SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 323130d5748SToby Isaac PetscValidPointer(formDegree, 2); 324130d5748SToby Isaac *formDegree = pt->formDegree; 325130d5748SToby Isaac PetscFunctionReturn(0); 326130d5748SToby Isaac } 327130d5748SToby Isaac 328130d5748SToby Isaac static PetscErrorCode PetscSpaceGetHeightSubspace_Ptrimmed(PetscSpace sp, PetscInt height, PetscSpace *subsp) 329130d5748SToby Isaac { 330130d5748SToby Isaac PetscSpace_Ptrimmed *pt = (PetscSpace_Ptrimmed *) sp->data; 331130d5748SToby Isaac PetscInt dim; 332130d5748SToby Isaac PetscErrorCode ierr; 333130d5748SToby Isaac 334130d5748SToby Isaac PetscFunctionBegin; 335130d5748SToby Isaac ierr = PetscSpaceGetNumVariables(sp, &dim);CHKERRQ(ierr); 336130d5748SToby 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); 337130d5748SToby Isaac if (!pt->subspaces) {ierr = PetscCalloc1(dim, &(pt->subspaces));CHKERRQ(ierr);} 338130d5748SToby Isaac if ((dim - height) <= PetscAbsInt(pt->formDegree)) { 339130d5748SToby Isaac if (!pt->subspaces[height-1]) { 340130d5748SToby Isaac PetscInt Nc, degree, Nf, Ncopies, Nfsub; 341130d5748SToby Isaac PetscSpace sub; 342130d5748SToby Isaac const char *name; 343130d5748SToby Isaac 344130d5748SToby Isaac ierr = PetscSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr); 345130d5748SToby Isaac ierr = PetscDTBinomialInt(dim, PetscAbsInt(pt->formDegree), &Nf);CHKERRQ(ierr); 346130d5748SToby Isaac ierr = PetscDTBinomialInt((dim-height), PetscAbsInt(pt->formDegree), &Nfsub);CHKERRQ(ierr); 347130d5748SToby Isaac Ncopies = Nf / Nc; 348130d5748SToby Isaac ierr = PetscSpaceGetDegree(sp, °ree, NULL);CHKERRQ(ierr); 349130d5748SToby Isaac 350130d5748SToby Isaac ierr = PetscSpaceCreate(PetscObjectComm((PetscObject) sp), &sub);CHKERRQ(ierr); 351130d5748SToby Isaac ierr = PetscObjectGetName((PetscObject) sp, &name);CHKERRQ(ierr); 352130d5748SToby Isaac ierr = PetscObjectSetName((PetscObject) sub, name);CHKERRQ(ierr); 353130d5748SToby Isaac ierr = PetscSpaceSetType(sub, PETSCSPACEPTRIMMED);CHKERRQ(ierr); 354130d5748SToby Isaac ierr = PetscSpaceSetNumComponents(sub, Nfsub * Ncopies);CHKERRQ(ierr); 355130d5748SToby Isaac ierr = PetscSpaceSetDegree(sub, degree, PETSC_DETERMINE);CHKERRQ(ierr); 356130d5748SToby Isaac ierr = PetscSpaceSetNumVariables(sub, dim-height);CHKERRQ(ierr); 357130d5748SToby Isaac ierr = PetscSpacePTrimmedSetFormDegree(sub, pt->formDegree);CHKERRQ(ierr); 358130d5748SToby Isaac ierr = PetscSpaceSetUp(sub);CHKERRQ(ierr); 359130d5748SToby Isaac pt->subspaces[height-1] = sub; 360130d5748SToby Isaac } 361130d5748SToby Isaac *subsp = pt->subspaces[height-1]; 362130d5748SToby Isaac } else { 363130d5748SToby Isaac *subsp = NULL; 364130d5748SToby Isaac } 365130d5748SToby Isaac PetscFunctionReturn(0); 366130d5748SToby Isaac } 367130d5748SToby Isaac 368130d5748SToby Isaac static PetscErrorCode PetscSpaceInitialize_Ptrimmed(PetscSpace sp) 369130d5748SToby Isaac { 370130d5748SToby Isaac PetscErrorCode ierr; 371130d5748SToby Isaac 372130d5748SToby Isaac PetscFunctionBegin; 373130d5748SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePTrimmedGetFormDegree_C", PetscSpacePTrimmedGetFormDegree_Ptrimmed);CHKERRQ(ierr); 374130d5748SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePTrimmedSetFormDegree_C", PetscSpacePTrimmedSetFormDegree_Ptrimmed);CHKERRQ(ierr); 375130d5748SToby Isaac sp->ops->setfromoptions = PetscSpaceSetFromOptions_Ptrimmed; 376130d5748SToby Isaac sp->ops->setup = PetscSpaceSetUp_Ptrimmed; 377130d5748SToby Isaac sp->ops->view = PetscSpaceView_Ptrimmed; 378130d5748SToby Isaac sp->ops->destroy = PetscSpaceDestroy_Ptrimmed; 379130d5748SToby Isaac sp->ops->getdimension = PetscSpaceGetDimension_Ptrimmed; 380130d5748SToby Isaac sp->ops->evaluate = PetscSpaceEvaluate_Ptrimmed; 381130d5748SToby Isaac sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Ptrimmed; 382130d5748SToby Isaac PetscFunctionReturn(0); 383130d5748SToby Isaac } 384130d5748SToby Isaac 385130d5748SToby Isaac /*MC 386130d5748SToby Isaac PETSCSPACEPTRIMMED = "ptrimmed" - A PetscSpace object that encapsulates a trimmed polynomial space. 387130d5748SToby Isaac 388130d5748SToby Isaac Level: intermediate 389130d5748SToby Isaac 390130d5748SToby Isaac .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType(), PetscDTPTrimmedEvalJet() 391130d5748SToby Isaac M*/ 392130d5748SToby Isaac 393130d5748SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Ptrimmed(PetscSpace sp) 394130d5748SToby Isaac { 395130d5748SToby Isaac PetscSpace_Ptrimmed *pt; 396130d5748SToby Isaac PetscErrorCode ierr; 397130d5748SToby Isaac 398130d5748SToby Isaac PetscFunctionBegin; 399130d5748SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 400130d5748SToby Isaac ierr = PetscNewLog(sp,&pt);CHKERRQ(ierr); 401130d5748SToby Isaac sp->data = pt; 402130d5748SToby Isaac 403130d5748SToby Isaac pt->subspaces = NULL; 404130d5748SToby Isaac sp->Nc = PETSC_DETERMINE; 405130d5748SToby Isaac 406130d5748SToby Isaac ierr = PetscSpaceInitialize_Ptrimmed(sp);CHKERRQ(ierr); 407130d5748SToby Isaac PetscFunctionReturn(0); 408130d5748SToby Isaac } 409130d5748SToby Isaac 410