120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 220cf1dd8SToby Isaac 33596293dSMatthew G. Knepley const char *const PetscSpacePolynomialTypes[] = {"P", "PMINUS_HDIV", "PMINUS_HCURL", "PetscSpacePolynomialType", "PETSCSPACE_POLYNOMIALTYPE_",0}; 43596293dSMatthew G. Knepley 5*29b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceSetFromOptions_Polynomial(PetscOptionItems *PetscOptionsObject,PetscSpace sp) 620cf1dd8SToby Isaac { 720cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 820cf1dd8SToby Isaac PetscErrorCode ierr; 920cf1dd8SToby Isaac 1020cf1dd8SToby Isaac PetscFunctionBegin; 1120cf1dd8SToby Isaac ierr = PetscOptionsHead(PetscOptionsObject,"PetscSpace polynomial options");CHKERRQ(ierr); 1220cf1dd8SToby Isaac ierr = PetscOptionsBool("-petscspace_poly_sym", "Use only symmetric polynomials", "PetscSpacePolynomialSetSymmetric", poly->symmetric, &poly->symmetric, NULL);CHKERRQ(ierr); 1320cf1dd8SToby Isaac ierr = PetscOptionsBool("-petscspace_poly_tensor", "Use the tensor product polynomials", "PetscSpacePolynomialSetTensor", poly->tensor, &poly->tensor, NULL);CHKERRQ(ierr); 143596293dSMatthew G. Knepley ierr = PetscOptionsEnum("-petscspace_poly_type", "Type of polynomial space", "PetscSpacePolynomialSetType", PetscSpacePolynomialTypes, (PetscEnum)poly->ptype, (PetscEnum*)&poly->ptype, NULL);CHKERRQ(ierr); 1520cf1dd8SToby Isaac ierr = PetscOptionsTail();CHKERRQ(ierr); 1620cf1dd8SToby Isaac PetscFunctionReturn(0); 1720cf1dd8SToby Isaac } 1820cf1dd8SToby Isaac 19d9bac1caSLisandro Dalcin static PetscErrorCode PetscSpacePolynomialView_Ascii(PetscSpace sp, PetscViewer v) 2020cf1dd8SToby Isaac { 2120cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 2220cf1dd8SToby Isaac PetscErrorCode ierr; 2320cf1dd8SToby Isaac 2420cf1dd8SToby Isaac PetscFunctionBegin; 253596293dSMatthew G. Knepley ierr = PetscViewerASCIIPrintf(v, "%s%s%s space of degree %D\n", poly->ptype ? PetscSpacePolynomialTypes[poly->ptype] : "", poly->ptype ? " " : "", poly->tensor ? "Tensor polynomial" : "Polynomial", sp->degree);CHKERRQ(ierr); 2620cf1dd8SToby Isaac PetscFunctionReturn(0); 2720cf1dd8SToby Isaac } 2820cf1dd8SToby Isaac 29*29b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceView_Polynomial(PetscSpace sp, PetscViewer viewer) 3020cf1dd8SToby Isaac { 3120cf1dd8SToby Isaac PetscBool iascii; 3220cf1dd8SToby Isaac PetscErrorCode ierr; 3320cf1dd8SToby Isaac 3420cf1dd8SToby Isaac PetscFunctionBegin; 3520cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 3620cf1dd8SToby Isaac PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 3720cf1dd8SToby Isaac ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 3820cf1dd8SToby Isaac if (iascii) {ierr = PetscSpacePolynomialView_Ascii(sp, viewer);CHKERRQ(ierr);} 3920cf1dd8SToby Isaac PetscFunctionReturn(0); 4020cf1dd8SToby Isaac } 4120cf1dd8SToby Isaac 42*29b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp) 4320cf1dd8SToby Isaac { 4420cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 4520cf1dd8SToby Isaac PetscInt ndegree = sp->degree+1; 4620cf1dd8SToby Isaac PetscInt deg; 4720cf1dd8SToby Isaac PetscErrorCode ierr; 4820cf1dd8SToby Isaac 4920cf1dd8SToby Isaac PetscFunctionBegin; 50982162a4SToby Isaac if (poly->setupCalled) PetscFunctionReturn(0); 5120cf1dd8SToby Isaac ierr = PetscMalloc1(ndegree, &poly->degrees);CHKERRQ(ierr); 5220cf1dd8SToby Isaac for (deg = 0; deg < ndegree; ++deg) poly->degrees[deg] = deg; 5320cf1dd8SToby Isaac if (poly->tensor) { 5420cf1dd8SToby Isaac sp->maxDegree = sp->degree + PetscMax(sp->Nv - 1,0); 5520cf1dd8SToby Isaac } else { 5620cf1dd8SToby Isaac sp->maxDegree = sp->degree; 5720cf1dd8SToby Isaac } 58982162a4SToby Isaac poly->setupCalled = PETSC_TRUE; 5920cf1dd8SToby Isaac PetscFunctionReturn(0); 6020cf1dd8SToby Isaac } 6120cf1dd8SToby Isaac 62*29b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceDestroy_Polynomial(PetscSpace sp) 6320cf1dd8SToby Isaac { 6420cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 6520cf1dd8SToby Isaac PetscErrorCode ierr; 6620cf1dd8SToby Isaac 6720cf1dd8SToby Isaac PetscFunctionBegin; 6820cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", NULL);CHKERRQ(ierr); 6920cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", NULL);CHKERRQ(ierr); 7020cf1dd8SToby Isaac ierr = PetscFree(poly->degrees);CHKERRQ(ierr); 7120cf1dd8SToby Isaac if (poly->subspaces) { 7220cf1dd8SToby Isaac PetscInt d; 7320cf1dd8SToby Isaac 7420cf1dd8SToby Isaac for (d = 0; d < sp->Nv; ++d) { 7520cf1dd8SToby Isaac ierr = PetscSpaceDestroy(&poly->subspaces[d]);CHKERRQ(ierr); 7620cf1dd8SToby Isaac } 7720cf1dd8SToby Isaac } 7820cf1dd8SToby Isaac ierr = PetscFree(poly->subspaces);CHKERRQ(ierr); 7920cf1dd8SToby Isaac ierr = PetscFree(poly);CHKERRQ(ierr); 8020cf1dd8SToby Isaac PetscFunctionReturn(0); 8120cf1dd8SToby Isaac } 8220cf1dd8SToby Isaac 8320cf1dd8SToby Isaac /* We treat the space as a tensor product of scalar polynomial spaces, so the dimension is multiplied by Nc */ 84*29b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim) 8520cf1dd8SToby Isaac { 8620cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 8720cf1dd8SToby Isaac PetscInt deg = sp->degree; 883596293dSMatthew G. Knepley PetscInt n = sp->Nv, N, i; 8920cf1dd8SToby Isaac PetscReal D = 1.0; 9020cf1dd8SToby Isaac 9120cf1dd8SToby Isaac PetscFunctionBegin; 923596293dSMatthew G. Knepley if ((poly->ptype == PETSCSPACE_POLYNOMIALTYPE_PMINUS_HDIV) || (poly->ptype == PETSCSPACE_POLYNOMIALTYPE_PMINUS_HCURL)) --deg; 9320cf1dd8SToby Isaac if (poly->tensor) { 943596293dSMatthew G. Knepley N = 1; 953596293dSMatthew G. Knepley for (i = 0; i < n; ++i) N *= (deg+1); 9620cf1dd8SToby Isaac } else { 9720cf1dd8SToby Isaac for (i = 1; i <= n; ++i) { 9820cf1dd8SToby Isaac D *= ((PetscReal) (deg+i))/i; 9920cf1dd8SToby Isaac } 1003596293dSMatthew G. Knepley N = (PetscInt) (D + 0.5); 10120cf1dd8SToby Isaac } 1023596293dSMatthew G. Knepley if ((poly->ptype == PETSCSPACE_POLYNOMIALTYPE_PMINUS_HDIV) || (poly->ptype == PETSCSPACE_POLYNOMIALTYPE_PMINUS_HCURL)) { 1033596293dSMatthew G. Knepley N *= sp->Nc + 1; 1043596293dSMatthew G. Knepley } else { 1053596293dSMatthew G. Knepley N *= sp->Nc; 1063596293dSMatthew G. Knepley } 1073596293dSMatthew G. Knepley *dim = N; 10820cf1dd8SToby Isaac PetscFunctionReturn(0); 10920cf1dd8SToby Isaac } 11020cf1dd8SToby Isaac 11120cf1dd8SToby Isaac /* 11220cf1dd8SToby Isaac LatticePoint_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to 'sum'. 11320cf1dd8SToby Isaac 11420cf1dd8SToby Isaac Input Parameters: 11520cf1dd8SToby Isaac + len - The length of the tuple 11620cf1dd8SToby Isaac . sum - The sum of all entries in the tuple 11720cf1dd8SToby Isaac - ind - The current multi-index of the tuple, initialized to the 0 tuple 11820cf1dd8SToby Isaac 11920cf1dd8SToby Isaac Output Parameter: 12020cf1dd8SToby Isaac + ind - The multi-index of the tuple, -1 indicates the iteration has terminated 12120cf1dd8SToby Isaac . tup - A tuple of len integers addig to sum 12220cf1dd8SToby Isaac 12320cf1dd8SToby Isaac Level: developer 12420cf1dd8SToby Isaac 12520cf1dd8SToby Isaac .seealso: 12620cf1dd8SToby Isaac */ 12720cf1dd8SToby Isaac static PetscErrorCode LatticePoint_Internal(PetscInt len, PetscInt sum, PetscInt ind[], PetscInt tup[]) 12820cf1dd8SToby Isaac { 12920cf1dd8SToby Isaac PetscInt i; 13020cf1dd8SToby Isaac PetscErrorCode ierr; 13120cf1dd8SToby Isaac 13220cf1dd8SToby Isaac PetscFunctionBegin; 13320cf1dd8SToby Isaac if (len == 1) { 13420cf1dd8SToby Isaac ind[0] = -1; 13520cf1dd8SToby Isaac tup[0] = sum; 13620cf1dd8SToby Isaac } else if (sum == 0) { 13720cf1dd8SToby Isaac for (i = 0; i < len; ++i) {ind[0] = -1; tup[i] = 0;} 13820cf1dd8SToby Isaac } else { 13920cf1dd8SToby Isaac tup[0] = sum - ind[0]; 14020cf1dd8SToby Isaac ierr = LatticePoint_Internal(len-1, ind[0], &ind[1], &tup[1]);CHKERRQ(ierr); 14120cf1dd8SToby Isaac if (ind[1] < 0) { 14220cf1dd8SToby Isaac if (ind[0] == sum) {ind[0] = -1;} 14320cf1dd8SToby Isaac else {ind[1] = 0; ++ind[0];} 14420cf1dd8SToby Isaac } 14520cf1dd8SToby Isaac } 14620cf1dd8SToby Isaac PetscFunctionReturn(0); 14720cf1dd8SToby Isaac } 14820cf1dd8SToby Isaac 14920cf1dd8SToby Isaac /* 15020cf1dd8SToby Isaac TensorPoint_Internal - Returns all tuples of size 'len' with nonnegative integers that are less than 'max'. 15120cf1dd8SToby Isaac 15220cf1dd8SToby Isaac Input Parameters: 15320cf1dd8SToby Isaac + len - The length of the tuple 15420cf1dd8SToby Isaac . max - The max for all entries in the tuple 15520cf1dd8SToby Isaac - ind - The current multi-index of the tuple, initialized to the 0 tuple 15620cf1dd8SToby Isaac 15720cf1dd8SToby Isaac Output Parameter: 15820cf1dd8SToby Isaac + ind - The multi-index of the tuple, -1 indicates the iteration has terminated 15920cf1dd8SToby Isaac . tup - A tuple of len integers less than max 16020cf1dd8SToby Isaac 16120cf1dd8SToby Isaac Level: developer 16220cf1dd8SToby Isaac 16320cf1dd8SToby Isaac .seealso: 16420cf1dd8SToby Isaac */ 16520cf1dd8SToby Isaac static PetscErrorCode TensorPoint_Internal(PetscInt len, PetscInt max, PetscInt ind[], PetscInt tup[]) 16620cf1dd8SToby Isaac { 16720cf1dd8SToby Isaac PetscInt i; 16820cf1dd8SToby Isaac PetscErrorCode ierr; 16920cf1dd8SToby Isaac 17020cf1dd8SToby Isaac PetscFunctionBegin; 17120cf1dd8SToby Isaac if (len == 1) { 17220cf1dd8SToby Isaac tup[0] = ind[0]++; 17320cf1dd8SToby Isaac ind[0] = ind[0] >= max ? -1 : ind[0]; 17420cf1dd8SToby Isaac } else if (max == 0) { 17520cf1dd8SToby Isaac for (i = 0; i < len; ++i) {ind[0] = -1; tup[i] = 0;} 17620cf1dd8SToby Isaac } else { 17720cf1dd8SToby Isaac tup[0] = ind[0]; 17820cf1dd8SToby Isaac ierr = TensorPoint_Internal(len-1, max, &ind[1], &tup[1]);CHKERRQ(ierr); 17920cf1dd8SToby Isaac if (ind[1] < 0) { 18020cf1dd8SToby Isaac ind[1] = 0; 18120cf1dd8SToby Isaac if (ind[0] == max-1) {ind[0] = -1;} 18220cf1dd8SToby Isaac else {++ind[0];} 18320cf1dd8SToby Isaac } 18420cf1dd8SToby Isaac } 18520cf1dd8SToby Isaac PetscFunctionReturn(0); 18620cf1dd8SToby Isaac } 18720cf1dd8SToby Isaac 18820cf1dd8SToby Isaac /* 18920cf1dd8SToby Isaac p in [0, npoints), i in [0, pdim), c in [0, Nc) 19020cf1dd8SToby Isaac 19120cf1dd8SToby Isaac B[p][i][c] = B[p][i_scalar][c][c] 19220cf1dd8SToby Isaac */ 193*29b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[]) 19420cf1dd8SToby Isaac { 1953596293dSMatthew G. Knepley const PetscInt eps[3][3][3] = {{{0, 0, 0}, {0, 0, 1}, {0, -1, 0}}, {{0, 0, -1}, {0, 0, 0}, {1, 0, 0}}, {{0, 1, 0}, {-1, 0, 0}, {0, 0, 0}}}; 19620cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 19720cf1dd8SToby Isaac DM dm = sp->dm; 19820cf1dd8SToby Isaac PetscInt Nc = sp->Nc; 19920cf1dd8SToby Isaac PetscInt ndegree = sp->degree+1; 20020cf1dd8SToby Isaac PetscInt *degrees = poly->degrees; 20120cf1dd8SToby Isaac PetscInt dim = sp->Nv; 20220cf1dd8SToby Isaac PetscReal *lpoints, *tmp, *LB, *LD, *LH; 20320cf1dd8SToby Isaac PetscInt *ind, *tup; 2043596293dSMatthew G. Knepley PetscInt c, pdim, pdimRed, d, e, der, der2, i, p, deg, o; 20520cf1dd8SToby Isaac PetscErrorCode ierr; 20620cf1dd8SToby Isaac 20720cf1dd8SToby Isaac PetscFunctionBegin; 20820cf1dd8SToby Isaac ierr = PetscSpaceGetDimension(sp, &pdim);CHKERRQ(ierr); 20920cf1dd8SToby Isaac pdim /= Nc; 21020cf1dd8SToby Isaac ierr = DMGetWorkArray(dm, npoints, MPIU_REAL, &lpoints);CHKERRQ(ierr); 21120cf1dd8SToby Isaac ierr = DMGetWorkArray(dm, npoints*ndegree*3, MPIU_REAL, &tmp);CHKERRQ(ierr); 2128c7916f5SToby Isaac if (B || D || H) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LB);CHKERRQ(ierr);} 2138c7916f5SToby Isaac if (D || H) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LD);CHKERRQ(ierr);} 21420cf1dd8SToby Isaac if (H) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LH);CHKERRQ(ierr);} 21520cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 21620cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 21720cf1dd8SToby Isaac lpoints[p] = points[p*dim+d]; 21820cf1dd8SToby Isaac } 21920cf1dd8SToby Isaac ierr = PetscDTLegendreEval(npoints, lpoints, ndegree, degrees, tmp, &tmp[1*npoints*ndegree], &tmp[2*npoints*ndegree]);CHKERRQ(ierr); 22020cf1dd8SToby Isaac /* LB, LD, LH (ndegree * dim x npoints) */ 22120cf1dd8SToby Isaac for (deg = 0; deg < ndegree; ++deg) { 22220cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 2238c7916f5SToby Isaac if (B || D || H) LB[(deg*dim + d)*npoints + p] = tmp[(0*npoints + p)*ndegree+deg]; 2248c7916f5SToby Isaac if (D || H) LD[(deg*dim + d)*npoints + p] = tmp[(1*npoints + p)*ndegree+deg]; 22520cf1dd8SToby Isaac if (H) LH[(deg*dim + d)*npoints + p] = tmp[(2*npoints + p)*ndegree+deg]; 22620cf1dd8SToby Isaac } 22720cf1dd8SToby Isaac } 22820cf1dd8SToby Isaac } 22920cf1dd8SToby Isaac /* Multiply by A (pdim x ndegree * dim) */ 23020cf1dd8SToby Isaac ierr = PetscMalloc2(dim,&ind,dim,&tup);CHKERRQ(ierr); 23120cf1dd8SToby Isaac if (B) { 2323596293dSMatthew G. Knepley PetscInt topDegree = sp->degree; 2333596293dSMatthew G. Knepley 23420cf1dd8SToby Isaac /* B (npoints x pdim x Nc) */ 235580bdb30SBarry Smith ierr = PetscArrayzero(B, npoints*pdim*Nc*Nc);CHKERRQ(ierr); 2363596293dSMatthew G. Knepley if ((poly->ptype == PETSCSPACE_POLYNOMIALTYPE_PMINUS_HDIV) || (poly->ptype == PETSCSPACE_POLYNOMIALTYPE_PMINUS_HCURL)) topDegree--; 2373596293dSMatthew G. Knepley /* Make complete space portion */ 23820cf1dd8SToby Isaac if (poly->tensor) { 2393596293dSMatthew G. Knepley if (poly->ptype != PETSCSPACE_POLYNOMIALTYPE_P) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_SUP, "Tensor spaces not supported for P^- spaces (%s)", PetscSpacePolynomialTypes[poly->ptype]); 24020cf1dd8SToby Isaac i = 0; 241580bdb30SBarry Smith ierr = PetscArrayzero(ind, dim);CHKERRQ(ierr); 24220cf1dd8SToby Isaac while (ind[0] >= 0) { 24320cf1dd8SToby Isaac ierr = TensorPoint_Internal(dim, sp->degree+1, ind, tup);CHKERRQ(ierr); 24420cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 24520cf1dd8SToby Isaac B[(p*pdim + i)*Nc*Nc] = 1.0; 24620cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 24720cf1dd8SToby Isaac B[(p*pdim + i)*Nc*Nc] *= LB[(tup[d]*dim + d)*npoints + p]; 24820cf1dd8SToby Isaac } 24920cf1dd8SToby Isaac } 25020cf1dd8SToby Isaac ++i; 25120cf1dd8SToby Isaac } 25220cf1dd8SToby Isaac } else { 25320cf1dd8SToby Isaac i = 0; 2543596293dSMatthew G. Knepley for (o = 0; o <= topDegree; ++o) { 255580bdb30SBarry Smith ierr = PetscArrayzero(ind, dim);CHKERRQ(ierr); 25620cf1dd8SToby Isaac while (ind[0] >= 0) { 25720cf1dd8SToby Isaac ierr = LatticePoint_Internal(dim, o, ind, tup);CHKERRQ(ierr); 25820cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 25920cf1dd8SToby Isaac B[(p*pdim + i)*Nc*Nc] = 1.0; 26020cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 26120cf1dd8SToby Isaac B[(p*pdim + i)*Nc*Nc] *= LB[(tup[d]*dim + d)*npoints + p]; 26220cf1dd8SToby Isaac } 26320cf1dd8SToby Isaac } 26420cf1dd8SToby Isaac ++i; 26520cf1dd8SToby Isaac } 26620cf1dd8SToby Isaac } 26720cf1dd8SToby Isaac } 2683596293dSMatthew G. Knepley pdimRed = i; 26920cf1dd8SToby Isaac /* Make direct sum basis for multicomponent space */ 27020cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 2713596293dSMatthew G. Knepley for (i = 0; i < pdimRed; ++i) { 27220cf1dd8SToby Isaac for (c = 1; c < Nc; ++c) { 27320cf1dd8SToby Isaac B[(p*pdim*Nc + i*Nc + c)*Nc + c] = B[(p*pdim + i)*Nc*Nc]; 27420cf1dd8SToby Isaac } 27520cf1dd8SToby Isaac } 27620cf1dd8SToby Isaac } 2773596293dSMatthew G. Knepley /* Make homogeneous part */ 2783596293dSMatthew G. Knepley if (topDegree < sp->degree) { 2793596293dSMatthew G. Knepley if (poly->tensor) { 2803596293dSMatthew G. Knepley } else { 2813596293dSMatthew G. Knepley i = pdimRed; 282580bdb30SBarry Smith ierr = PetscArrayzero(ind, dim);CHKERRQ(ierr); 2833596293dSMatthew G. Knepley while (ind[0] >= 0) { 2843596293dSMatthew G. Knepley ierr = LatticePoint_Internal(dim, topDegree, ind, tup);CHKERRQ(ierr); 2853596293dSMatthew G. Knepley for (p = 0; p < npoints; ++p) { 2863596293dSMatthew G. Knepley for (c = 0; c < Nc; ++c) { 2873596293dSMatthew G. Knepley B[(p*pdim*Nc + i*Nc + c)*Nc + c] = 1.0; 2883596293dSMatthew G. Knepley for (d = 0; d < dim; ++d) { 2893596293dSMatthew G. Knepley B[(p*pdim*Nc + i*Nc + c)*Nc + c] *= LB[(tup[d]*dim + d)*npoints + p]; 2903596293dSMatthew G. Knepley } 2913596293dSMatthew G. Knepley switch (poly->ptype) { 2923596293dSMatthew G. Knepley case PETSCSPACE_POLYNOMIALTYPE_PMINUS_HDIV: 2933596293dSMatthew G. Knepley B[(p*pdim*Nc + i*Nc + c)*Nc + c] *= LB[(c*dim + d)*npoints + p];break; 2943596293dSMatthew G. Knepley case PETSCSPACE_POLYNOMIALTYPE_PMINUS_HCURL: 2953596293dSMatthew G. Knepley { 2963596293dSMatthew G. Knepley PetscReal sum = 0.0; 2973596293dSMatthew G. Knepley for (d = 0; d < dim; ++d) for (e = 0; e < dim; ++e) sum += eps[c][d][e]*LB[(d*dim + d)*npoints + p]; 2983596293dSMatthew G. Knepley B[(p*pdim*Nc + i*Nc + c)*Nc + c] *= sum; 2993596293dSMatthew G. Knepley break; 3003596293dSMatthew G. Knepley } 3013596293dSMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_SUP, "Invalid polynomial type %s", PetscSpacePolynomialTypes[poly->ptype]); 3023596293dSMatthew G. Knepley } 3033596293dSMatthew G. Knepley } 3043596293dSMatthew G. Knepley } 3053596293dSMatthew G. Knepley ++i; 3063596293dSMatthew G. Knepley } 3073596293dSMatthew G. Knepley } 3083596293dSMatthew G. Knepley } 30920cf1dd8SToby Isaac } 31020cf1dd8SToby Isaac if (D) { 3113596293dSMatthew G. Knepley if (poly->ptype != PETSCSPACE_POLYNOMIALTYPE_P) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_SUP, "Derivatives not supported for P^- spaces (%s)", PetscSpacePolynomialTypes[poly->ptype]); 31220cf1dd8SToby Isaac /* D (npoints x pdim x Nc x dim) */ 313580bdb30SBarry Smith ierr = PetscArrayzero(D, npoints*pdim*Nc*Nc*dim);CHKERRQ(ierr); 31420cf1dd8SToby Isaac if (poly->tensor) { 31520cf1dd8SToby Isaac i = 0; 316580bdb30SBarry Smith ierr = PetscArrayzero(ind, dim);CHKERRQ(ierr); 31720cf1dd8SToby Isaac while (ind[0] >= 0) { 31820cf1dd8SToby Isaac ierr = TensorPoint_Internal(dim, sp->degree+1, ind, tup);CHKERRQ(ierr); 31920cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 32020cf1dd8SToby Isaac for (der = 0; der < dim; ++der) { 32120cf1dd8SToby Isaac D[(p*pdim + i)*Nc*Nc*dim + der] = 1.0; 32220cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 32320cf1dd8SToby Isaac if (d == der) { 32420cf1dd8SToby Isaac D[(p*pdim + i)*Nc*Nc*dim + der] *= LD[(tup[d]*dim + d)*npoints + p]; 32520cf1dd8SToby Isaac } else { 32620cf1dd8SToby Isaac D[(p*pdim + i)*Nc*Nc*dim + der] *= LB[(tup[d]*dim + d)*npoints + p]; 32720cf1dd8SToby Isaac } 32820cf1dd8SToby Isaac } 32920cf1dd8SToby Isaac } 33020cf1dd8SToby Isaac } 33120cf1dd8SToby Isaac ++i; 33220cf1dd8SToby Isaac } 33320cf1dd8SToby Isaac } else { 33420cf1dd8SToby Isaac i = 0; 33520cf1dd8SToby Isaac for (o = 0; o <= sp->degree; ++o) { 336580bdb30SBarry Smith ierr = PetscArrayzero(ind, dim);CHKERRQ(ierr); 33720cf1dd8SToby Isaac while (ind[0] >= 0) { 33820cf1dd8SToby Isaac ierr = LatticePoint_Internal(dim, o, ind, tup);CHKERRQ(ierr); 33920cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 34020cf1dd8SToby Isaac for (der = 0; der < dim; ++der) { 34120cf1dd8SToby Isaac D[(p*pdim + i)*Nc*Nc*dim + der] = 1.0; 34220cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 34320cf1dd8SToby Isaac if (d == der) { 34420cf1dd8SToby Isaac D[(p*pdim + i)*Nc*Nc*dim + der] *= LD[(tup[d]*dim + d)*npoints + p]; 34520cf1dd8SToby Isaac } else { 34620cf1dd8SToby Isaac D[(p*pdim + i)*Nc*Nc*dim + der] *= LB[(tup[d]*dim + d)*npoints + p]; 34720cf1dd8SToby Isaac } 34820cf1dd8SToby Isaac } 34920cf1dd8SToby Isaac } 35020cf1dd8SToby Isaac } 35120cf1dd8SToby Isaac ++i; 35220cf1dd8SToby Isaac } 35320cf1dd8SToby Isaac } 35420cf1dd8SToby Isaac } 35520cf1dd8SToby Isaac /* Make direct sum basis for multicomponent space */ 35620cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 35720cf1dd8SToby Isaac for (i = 0; i < pdim; ++i) { 35820cf1dd8SToby Isaac for (c = 1; c < Nc; ++c) { 35920cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 36020cf1dd8SToby Isaac D[((p*pdim*Nc + i*Nc + c)*Nc + c)*dim + d] = D[(p*pdim + i)*Nc*Nc*dim + d]; 36120cf1dd8SToby Isaac } 36220cf1dd8SToby Isaac } 36320cf1dd8SToby Isaac } 36420cf1dd8SToby Isaac } 36520cf1dd8SToby Isaac } 36620cf1dd8SToby Isaac if (H) { 3673596293dSMatthew G. Knepley if (poly->ptype != PETSCSPACE_POLYNOMIALTYPE_P) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_SUP, "Hessians not supported for P^- spaces (%s)", PetscSpacePolynomialTypes[poly->ptype]); 36820cf1dd8SToby Isaac /* H (npoints x pdim x Nc x Nc x dim x dim) */ 369580bdb30SBarry Smith ierr = PetscArrayzero(H, npoints*pdim*Nc*Nc*dim*dim);CHKERRQ(ierr); 37020cf1dd8SToby Isaac if (poly->tensor) { 37120cf1dd8SToby Isaac i = 0; 372580bdb30SBarry Smith ierr = PetscArrayzero(ind, dim);CHKERRQ(ierr); 37320cf1dd8SToby Isaac while (ind[0] >= 0) { 37420cf1dd8SToby Isaac ierr = TensorPoint_Internal(dim, sp->degree+1, ind, tup);CHKERRQ(ierr); 37520cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 37620cf1dd8SToby Isaac for (der = 0; der < dim; ++der) { 37720cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der] = 1.0; 37820cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 37920cf1dd8SToby Isaac if (d == der) { 38020cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der] *= LH[(tup[d]*dim + d)*npoints + p]; 38120cf1dd8SToby Isaac } else { 38220cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der] *= LB[(tup[d]*dim + d)*npoints + p]; 38320cf1dd8SToby Isaac } 38420cf1dd8SToby Isaac } 38520cf1dd8SToby Isaac for (der2 = der + 1; der2 < dim; ++der2) { 38620cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] = 1.0; 38720cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 38820cf1dd8SToby Isaac if (d == der || d == der2) { 38920cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] *= LD[(tup[d]*dim + d)*npoints + p]; 39020cf1dd8SToby Isaac } else { 39120cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] *= LB[(tup[d]*dim + d)*npoints + p]; 39220cf1dd8SToby Isaac } 39320cf1dd8SToby Isaac } 39420cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der2) * dim + der] = H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2]; 39520cf1dd8SToby Isaac } 39620cf1dd8SToby Isaac } 39720cf1dd8SToby Isaac } 39820cf1dd8SToby Isaac ++i; 39920cf1dd8SToby Isaac } 40020cf1dd8SToby Isaac } else { 40120cf1dd8SToby Isaac i = 0; 40220cf1dd8SToby Isaac for (o = 0; o <= sp->degree; ++o) { 403580bdb30SBarry Smith ierr = PetscArrayzero(ind, dim);CHKERRQ(ierr); 40420cf1dd8SToby Isaac while (ind[0] >= 0) { 40520cf1dd8SToby Isaac ierr = LatticePoint_Internal(dim, o, ind, tup);CHKERRQ(ierr); 40620cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 40720cf1dd8SToby Isaac for (der = 0; der < dim; ++der) { 40820cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der)*dim + der] = 1.0; 40920cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 41020cf1dd8SToby Isaac if (d == der) { 41120cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der)*dim + der] *= LH[(tup[d]*dim + d)*npoints + p]; 41220cf1dd8SToby Isaac } else { 41320cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der)*dim + der] *= LB[(tup[d]*dim + d)*npoints + p]; 41420cf1dd8SToby Isaac } 41520cf1dd8SToby Isaac } 41620cf1dd8SToby Isaac for (der2 = der + 1; der2 < dim; ++der2) { 41720cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] = 1.0; 41820cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 41920cf1dd8SToby Isaac if (d == der || d == der2) { 42020cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] *= LD[(tup[d]*dim + d)*npoints + p]; 42120cf1dd8SToby Isaac } else { 42220cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] *= LB[(tup[d]*dim + d)*npoints + p]; 42320cf1dd8SToby Isaac } 42420cf1dd8SToby Isaac } 42520cf1dd8SToby Isaac H[((p*pdim + i)*Nc*Nc*dim + der2) * dim + der] = H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2]; 42620cf1dd8SToby Isaac } 42720cf1dd8SToby Isaac } 42820cf1dd8SToby Isaac } 42920cf1dd8SToby Isaac ++i; 43020cf1dd8SToby Isaac } 43120cf1dd8SToby Isaac } 43220cf1dd8SToby Isaac } 43320cf1dd8SToby Isaac /* Make direct sum basis for multicomponent space */ 43420cf1dd8SToby Isaac for (p = 0; p < npoints; ++p) { 43520cf1dd8SToby Isaac for (i = 0; i < pdim; ++i) { 43620cf1dd8SToby Isaac for (c = 1; c < Nc; ++c) { 43720cf1dd8SToby Isaac for (d = 0; d < dim; ++d) { 43820cf1dd8SToby Isaac for (e = 0; e < dim; ++e) { 43920cf1dd8SToby Isaac H[(((p*pdim*Nc + i*Nc + c)*Nc + c)*dim + d)*dim + e] = H[((p*pdim + i)*Nc*Nc*dim + d)*dim + e]; 44020cf1dd8SToby Isaac } 44120cf1dd8SToby Isaac } 44220cf1dd8SToby Isaac } 44320cf1dd8SToby Isaac } 44420cf1dd8SToby Isaac } 44520cf1dd8SToby Isaac } 44620cf1dd8SToby Isaac ierr = PetscFree2(ind,tup);CHKERRQ(ierr); 44720cf1dd8SToby Isaac if (H) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LH);CHKERRQ(ierr);} 4488c7916f5SToby Isaac if (D || H) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LD);CHKERRQ(ierr);} 4498c7916f5SToby Isaac if (B || D || H) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LB);CHKERRQ(ierr);} 45020cf1dd8SToby Isaac ierr = DMRestoreWorkArray(dm, npoints*ndegree*3, MPIU_REAL, &tmp);CHKERRQ(ierr); 45120cf1dd8SToby Isaac ierr = DMRestoreWorkArray(dm, npoints, MPIU_REAL, &lpoints);CHKERRQ(ierr); 45220cf1dd8SToby Isaac PetscFunctionReturn(0); 45320cf1dd8SToby Isaac } 45420cf1dd8SToby Isaac 45520cf1dd8SToby Isaac /*@ 45620cf1dd8SToby Isaac PetscSpacePolynomialSetTensor - Set whether a function space is a space of tensor polynomials (the space is spanned 45720cf1dd8SToby Isaac by polynomials whose degree in each variabl is bounded by the given order), as opposed to polynomials (the space is 45820cf1dd8SToby Isaac spanned by polynomials whose total degree---summing over all variables---is bounded by the given order). 45920cf1dd8SToby Isaac 46020cf1dd8SToby Isaac Input Parameters: 46120cf1dd8SToby Isaac + sp - the function space object 46220cf1dd8SToby Isaac - tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space 46320cf1dd8SToby Isaac 4644ab77754SMatthew G. Knepley Options Database: 4654ab77754SMatthew G. Knepley . -petscspace_poly_tensor <bool> - Whether to use tensor product polynomials in higher dimension 4664ab77754SMatthew G. Knepley 467*29b5c115SMatthew G. Knepley Level: intermediate 46820cf1dd8SToby Isaac 46920cf1dd8SToby Isaac .seealso: PetscSpacePolynomialGetTensor(), PetscSpaceSetDegree(), PetscSpaceSetNumVariables() 47020cf1dd8SToby Isaac @*/ 47120cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace sp, PetscBool tensor) 47220cf1dd8SToby Isaac { 47320cf1dd8SToby Isaac PetscErrorCode ierr; 47420cf1dd8SToby Isaac 47520cf1dd8SToby Isaac PetscFunctionBegin; 47620cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 47720cf1dd8SToby Isaac ierr = PetscTryMethod(sp,"PetscSpacePolynomialSetTensor_C",(PetscSpace,PetscBool),(sp,tensor));CHKERRQ(ierr); 47820cf1dd8SToby Isaac PetscFunctionReturn(0); 47920cf1dd8SToby Isaac } 48020cf1dd8SToby Isaac 48120cf1dd8SToby Isaac /*@ 48220cf1dd8SToby Isaac PetscSpacePolynomialGetTensor - Get whether a function space is a space of tensor polynomials (the space is spanned 48320cf1dd8SToby Isaac by polynomials whose degree in each variabl is bounded by the given order), as opposed to polynomials (the space is 48420cf1dd8SToby Isaac spanned by polynomials whose total degree---summing over all variables---is bounded by the given order). 48520cf1dd8SToby Isaac 48620cf1dd8SToby Isaac Input Parameters: 48720cf1dd8SToby Isaac . sp - the function space object 48820cf1dd8SToby Isaac 48920cf1dd8SToby Isaac Output Parameters: 49020cf1dd8SToby Isaac . tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space 49120cf1dd8SToby Isaac 492*29b5c115SMatthew G. Knepley Level: intermediate 49320cf1dd8SToby Isaac 49420cf1dd8SToby Isaac .seealso: PetscSpacePolynomialSetTensor(), PetscSpaceSetDegree(), PetscSpaceSetNumVariables() 49520cf1dd8SToby Isaac @*/ 49620cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace sp, PetscBool *tensor) 49720cf1dd8SToby Isaac { 49820cf1dd8SToby Isaac PetscErrorCode ierr; 49920cf1dd8SToby Isaac 50020cf1dd8SToby Isaac PetscFunctionBegin; 50120cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 50220cf1dd8SToby Isaac PetscValidPointer(tensor, 2); 50320cf1dd8SToby Isaac ierr = PetscTryMethod(sp,"PetscSpacePolynomialGetTensor_C",(PetscSpace,PetscBool*),(sp,tensor));CHKERRQ(ierr); 50420cf1dd8SToby Isaac PetscFunctionReturn(0); 50520cf1dd8SToby Isaac } 50620cf1dd8SToby Isaac 50720cf1dd8SToby Isaac static PetscErrorCode PetscSpacePolynomialSetTensor_Polynomial(PetscSpace sp, PetscBool tensor) 50820cf1dd8SToby Isaac { 50920cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 51020cf1dd8SToby Isaac 51120cf1dd8SToby Isaac PetscFunctionBegin; 51220cf1dd8SToby Isaac poly->tensor = tensor; 51320cf1dd8SToby Isaac PetscFunctionReturn(0); 51420cf1dd8SToby Isaac } 51520cf1dd8SToby Isaac 51620cf1dd8SToby Isaac static PetscErrorCode PetscSpacePolynomialGetTensor_Polynomial(PetscSpace sp, PetscBool *tensor) 51720cf1dd8SToby Isaac { 51820cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 51920cf1dd8SToby Isaac 52020cf1dd8SToby Isaac PetscFunctionBegin; 52120cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 52220cf1dd8SToby Isaac PetscValidPointer(tensor, 2); 52320cf1dd8SToby Isaac *tensor = poly->tensor; 52420cf1dd8SToby Isaac PetscFunctionReturn(0); 52520cf1dd8SToby Isaac } 52620cf1dd8SToby Isaac 52720cf1dd8SToby Isaac static PetscErrorCode PetscSpaceGetHeightSubspace_Polynomial(PetscSpace sp, PetscInt height, PetscSpace *subsp) 52820cf1dd8SToby Isaac { 52920cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 53020cf1dd8SToby Isaac PetscInt Nc, dim, order; 53120cf1dd8SToby Isaac PetscBool tensor; 53220cf1dd8SToby Isaac PetscErrorCode ierr; 53320cf1dd8SToby Isaac 53420cf1dd8SToby Isaac PetscFunctionBegin; 53520cf1dd8SToby Isaac ierr = PetscSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr); 53620cf1dd8SToby Isaac ierr = PetscSpaceGetNumVariables(sp, &dim);CHKERRQ(ierr); 53720cf1dd8SToby Isaac ierr = PetscSpaceGetDegree(sp, &order, NULL);CHKERRQ(ierr); 53820cf1dd8SToby Isaac ierr = PetscSpacePolynomialGetTensor(sp, &tensor);CHKERRQ(ierr); 53920cf1dd8SToby 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);} 54020cf1dd8SToby Isaac if (!poly->subspaces) {ierr = PetscCalloc1(dim, &poly->subspaces);CHKERRQ(ierr);} 54120cf1dd8SToby Isaac if (height <= dim) { 54220cf1dd8SToby Isaac if (!poly->subspaces[height-1]) { 54320cf1dd8SToby Isaac PetscSpace sub; 5443f6b16c7SMatthew G. Knepley const char *name; 54520cf1dd8SToby Isaac 54620cf1dd8SToby Isaac ierr = PetscSpaceCreate(PetscObjectComm((PetscObject) sp), &sub);CHKERRQ(ierr); 5473f6b16c7SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) sp, &name);CHKERRQ(ierr); 5483f6b16c7SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) sub, name);CHKERRQ(ierr); 54920cf1dd8SToby Isaac ierr = PetscSpaceSetType(sub, PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr); 550d39dd5f5SToby Isaac ierr = PetscSpaceSetNumComponents(sub, Nc);CHKERRQ(ierr); 551d39dd5f5SToby Isaac ierr = PetscSpaceSetDegree(sub, order, PETSC_DETERMINE);CHKERRQ(ierr); 55220cf1dd8SToby Isaac ierr = PetscSpaceSetNumVariables(sub, dim-height);CHKERRQ(ierr); 55320cf1dd8SToby Isaac ierr = PetscSpacePolynomialSetTensor(sub, tensor);CHKERRQ(ierr); 55420cf1dd8SToby Isaac ierr = PetscSpaceSetUp(sub);CHKERRQ(ierr); 55520cf1dd8SToby Isaac poly->subspaces[height-1] = sub; 55620cf1dd8SToby Isaac } 55720cf1dd8SToby Isaac *subsp = poly->subspaces[height-1]; 55820cf1dd8SToby Isaac } else { 55920cf1dd8SToby Isaac *subsp = NULL; 56020cf1dd8SToby Isaac } 56120cf1dd8SToby Isaac PetscFunctionReturn(0); 56220cf1dd8SToby Isaac } 56320cf1dd8SToby Isaac 564*29b5c115SMatthew G. Knepley static PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp) 56520cf1dd8SToby Isaac { 56620cf1dd8SToby Isaac PetscErrorCode ierr; 56720cf1dd8SToby Isaac 56820cf1dd8SToby Isaac PetscFunctionBegin; 56920cf1dd8SToby Isaac sp->ops->setfromoptions = PetscSpaceSetFromOptions_Polynomial; 57020cf1dd8SToby Isaac sp->ops->setup = PetscSpaceSetUp_Polynomial; 57120cf1dd8SToby Isaac sp->ops->view = PetscSpaceView_Polynomial; 57220cf1dd8SToby Isaac sp->ops->destroy = PetscSpaceDestroy_Polynomial; 57320cf1dd8SToby Isaac sp->ops->getdimension = PetscSpaceGetDimension_Polynomial; 57420cf1dd8SToby Isaac sp->ops->evaluate = PetscSpaceEvaluate_Polynomial; 57520cf1dd8SToby Isaac sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Polynomial; 57620cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Polynomial);CHKERRQ(ierr); 57720cf1dd8SToby Isaac ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", PetscSpacePolynomialSetTensor_Polynomial);CHKERRQ(ierr); 57820cf1dd8SToby Isaac PetscFunctionReturn(0); 57920cf1dd8SToby Isaac } 58020cf1dd8SToby Isaac 58120cf1dd8SToby Isaac /*MC 58220cf1dd8SToby Isaac PETSCSPACEPOLYNOMIAL = "poly" - A PetscSpace object that encapsulates a polynomial space, e.g. P1 is the space of 58320cf1dd8SToby Isaac linear polynomials. The space is replicated for each component. 58420cf1dd8SToby Isaac 58520cf1dd8SToby Isaac Level: intermediate 58620cf1dd8SToby Isaac 58720cf1dd8SToby Isaac .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType() 58820cf1dd8SToby Isaac M*/ 58920cf1dd8SToby Isaac 59020cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp) 59120cf1dd8SToby Isaac { 59220cf1dd8SToby Isaac PetscSpace_Poly *poly; 59320cf1dd8SToby Isaac PetscErrorCode ierr; 59420cf1dd8SToby Isaac 59520cf1dd8SToby Isaac PetscFunctionBegin; 59620cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 59720cf1dd8SToby Isaac ierr = PetscNewLog(sp,&poly);CHKERRQ(ierr); 59820cf1dd8SToby Isaac sp->data = poly; 59920cf1dd8SToby Isaac 60020cf1dd8SToby Isaac poly->symmetric = PETSC_FALSE; 60120cf1dd8SToby Isaac poly->tensor = PETSC_FALSE; 60220cf1dd8SToby Isaac poly->degrees = NULL; 60320cf1dd8SToby Isaac poly->subspaces = NULL; 60420cf1dd8SToby Isaac 60520cf1dd8SToby Isaac ierr = PetscSpaceInitialize_Polynomial(sp);CHKERRQ(ierr); 60620cf1dd8SToby Isaac PetscFunctionReturn(0); 60720cf1dd8SToby Isaac } 60820cf1dd8SToby Isaac 609*29b5c115SMatthew G. Knepley /*@ 610*29b5c115SMatthew G. Knepley PetscSpacePolynomialSetSymmetric - Set whether a function space is a space of symmetric polynomials 611*29b5c115SMatthew G. Knepley 612*29b5c115SMatthew G. Knepley Input Parameters: 613*29b5c115SMatthew G. Knepley + sp - the function space object 614*29b5c115SMatthew G. Knepley - sym - flag for symmetric polynomials 615*29b5c115SMatthew G. Knepley 616*29b5c115SMatthew G. Knepley Options Database: 617*29b5c115SMatthew G. Knepley . -petscspace_poly_sym <bool> - Whether to use symmetric polynomials 618*29b5c115SMatthew G. Knepley 619*29b5c115SMatthew G. Knepley Level: intermediate 620*29b5c115SMatthew G. Knepley 621*29b5c115SMatthew G. Knepley .seealso: PetscSpacePolynomialGetSymmetric(), PetscSpacePolynomialGetTensor(), PetscSpaceSetDegree(), PetscSpaceSetNumVariables() 622*29b5c115SMatthew G. Knepley @*/ 62320cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace sp, PetscBool sym) 62420cf1dd8SToby Isaac { 62520cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 62620cf1dd8SToby Isaac 62720cf1dd8SToby Isaac PetscFunctionBegin; 62820cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 62920cf1dd8SToby Isaac poly->symmetric = sym; 63020cf1dd8SToby Isaac PetscFunctionReturn(0); 63120cf1dd8SToby Isaac } 63220cf1dd8SToby Isaac 633*29b5c115SMatthew G. Knepley /*@ 634*29b5c115SMatthew G. Knepley PetscSpacePolynomialGetSymmetric - Get whether a function space is a space of symmetric polynomials 635*29b5c115SMatthew G. Knepley 636*29b5c115SMatthew G. Knepley Input Parameter: 637*29b5c115SMatthew G. Knepley . sp - the function space object 638*29b5c115SMatthew G. Knepley 639*29b5c115SMatthew G. Knepley Output Parameter: 640*29b5c115SMatthew G. Knepley . sym - flag for symmetric polynomials 641*29b5c115SMatthew G. Knepley 642*29b5c115SMatthew G. Knepley Level: intermediate 643*29b5c115SMatthew G. Knepley 644*29b5c115SMatthew G. Knepley .seealso: PetscSpacePolynomialSetSymmetric(), PetscSpacePolynomialGetTensor(), PetscSpaceSetDegree(), PetscSpaceSetNumVariables() 645*29b5c115SMatthew G. Knepley @*/ 64620cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialGetSymmetric(PetscSpace sp, PetscBool *sym) 64720cf1dd8SToby Isaac { 64820cf1dd8SToby Isaac PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data; 64920cf1dd8SToby Isaac 65020cf1dd8SToby Isaac PetscFunctionBegin; 65120cf1dd8SToby Isaac PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1); 65220cf1dd8SToby Isaac PetscValidPointer(sym, 2); 65320cf1dd8SToby Isaac *sym = poly->symmetric; 65420cf1dd8SToby Isaac PetscFunctionReturn(0); 65520cf1dd8SToby Isaac } 656