xref: /petsc/src/dm/dt/space/impls/poly/spacepoly.c (revision 4ab777544fba258ca449a465d0202acf6f77424b)
120cf1dd8SToby Isaac #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/
220cf1dd8SToby Isaac 
320cf1dd8SToby Isaac PetscErrorCode PetscSpaceSetFromOptions_Polynomial(PetscOptionItems *PetscOptionsObject,PetscSpace sp)
420cf1dd8SToby Isaac {
520cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
620cf1dd8SToby Isaac   PetscErrorCode   ierr;
720cf1dd8SToby Isaac 
820cf1dd8SToby Isaac   PetscFunctionBegin;
920cf1dd8SToby Isaac   ierr = PetscOptionsHead(PetscOptionsObject,"PetscSpace polynomial options");CHKERRQ(ierr);
1020cf1dd8SToby Isaac   ierr = PetscOptionsBool("-petscspace_poly_sym", "Use only symmetric polynomials", "PetscSpacePolynomialSetSymmetric", poly->symmetric, &poly->symmetric, NULL);CHKERRQ(ierr);
1120cf1dd8SToby Isaac   ierr = PetscOptionsBool("-petscspace_poly_tensor", "Use the tensor product polynomials", "PetscSpacePolynomialSetTensor", poly->tensor, &poly->tensor, NULL);CHKERRQ(ierr);
1220cf1dd8SToby Isaac   ierr = PetscOptionsTail();CHKERRQ(ierr);
1320cf1dd8SToby Isaac   PetscFunctionReturn(0);
1420cf1dd8SToby Isaac }
1520cf1dd8SToby Isaac 
1620cf1dd8SToby Isaac static PetscErrorCode PetscSpacePolynomialView_Ascii(PetscSpace sp, PetscViewer viewer)
1720cf1dd8SToby Isaac {
1820cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
1920cf1dd8SToby Isaac   PetscErrorCode   ierr;
2020cf1dd8SToby Isaac 
2120cf1dd8SToby Isaac   PetscFunctionBegin;
2220cf1dd8SToby Isaac   if (poly->tensor) {ierr = PetscViewerASCIIPrintf(viewer, "Tensor polynomial space of degree %D\n", sp->degree);CHKERRQ(ierr);}
2320cf1dd8SToby Isaac   else              {ierr = PetscViewerASCIIPrintf(viewer, "Polynomial space of degree %D\n", sp->degree);CHKERRQ(ierr);}
2420cf1dd8SToby Isaac   PetscFunctionReturn(0);
2520cf1dd8SToby Isaac }
2620cf1dd8SToby Isaac 
2720cf1dd8SToby Isaac PetscErrorCode PetscSpaceView_Polynomial(PetscSpace sp, PetscViewer viewer)
2820cf1dd8SToby Isaac {
2920cf1dd8SToby Isaac   PetscBool      iascii;
3020cf1dd8SToby Isaac   PetscErrorCode ierr;
3120cf1dd8SToby Isaac 
3220cf1dd8SToby Isaac   PetscFunctionBegin;
3320cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
3420cf1dd8SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
3520cf1dd8SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
3620cf1dd8SToby Isaac   if (iascii) {ierr = PetscSpacePolynomialView_Ascii(sp, viewer);CHKERRQ(ierr);}
3720cf1dd8SToby Isaac   PetscFunctionReturn(0);
3820cf1dd8SToby Isaac }
3920cf1dd8SToby Isaac 
4020cf1dd8SToby Isaac PetscErrorCode PetscSpaceSetUp_Polynomial(PetscSpace sp)
4120cf1dd8SToby Isaac {
4220cf1dd8SToby Isaac   PetscSpace_Poly *poly    = (PetscSpace_Poly *) sp->data;
4320cf1dd8SToby Isaac   PetscInt         ndegree = sp->degree+1;
4420cf1dd8SToby Isaac   PetscInt         deg;
4520cf1dd8SToby Isaac   PetscErrorCode   ierr;
4620cf1dd8SToby Isaac 
4720cf1dd8SToby Isaac   PetscFunctionBegin;
48982162a4SToby Isaac   if (poly->setupCalled) PetscFunctionReturn(0);
4920cf1dd8SToby Isaac   ierr = PetscMalloc1(ndegree, &poly->degrees);CHKERRQ(ierr);
5020cf1dd8SToby Isaac   for (deg = 0; deg < ndegree; ++deg) poly->degrees[deg] = deg;
5120cf1dd8SToby Isaac   if (poly->tensor) {
5220cf1dd8SToby Isaac     sp->maxDegree = sp->degree + PetscMax(sp->Nv - 1,0);
5320cf1dd8SToby Isaac   } else {
5420cf1dd8SToby Isaac     sp->maxDegree = sp->degree;
5520cf1dd8SToby Isaac   }
56982162a4SToby Isaac   poly->setupCalled = PETSC_TRUE;
5720cf1dd8SToby Isaac   PetscFunctionReturn(0);
5820cf1dd8SToby Isaac }
5920cf1dd8SToby Isaac 
6020cf1dd8SToby Isaac PetscErrorCode PetscSpaceDestroy_Polynomial(PetscSpace sp)
6120cf1dd8SToby Isaac {
6220cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
6320cf1dd8SToby Isaac   PetscErrorCode   ierr;
6420cf1dd8SToby Isaac 
6520cf1dd8SToby Isaac   PetscFunctionBegin;
6620cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", NULL);CHKERRQ(ierr);
6720cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", NULL);CHKERRQ(ierr);
6820cf1dd8SToby Isaac   ierr = PetscFree(poly->degrees);CHKERRQ(ierr);
6920cf1dd8SToby Isaac   if (poly->subspaces) {
7020cf1dd8SToby Isaac     PetscInt d;
7120cf1dd8SToby Isaac 
7220cf1dd8SToby Isaac     for (d = 0; d < sp->Nv; ++d) {
7320cf1dd8SToby Isaac       ierr = PetscSpaceDestroy(&poly->subspaces[d]);CHKERRQ(ierr);
7420cf1dd8SToby Isaac     }
7520cf1dd8SToby Isaac   }
7620cf1dd8SToby Isaac   ierr = PetscFree(poly->subspaces);CHKERRQ(ierr);
7720cf1dd8SToby Isaac   ierr = PetscFree(poly);CHKERRQ(ierr);
7820cf1dd8SToby Isaac   PetscFunctionReturn(0);
7920cf1dd8SToby Isaac }
8020cf1dd8SToby Isaac 
8120cf1dd8SToby Isaac /* We treat the space as a tensor product of scalar polynomial spaces, so the dimension is multiplied by Nc */
8220cf1dd8SToby Isaac PetscErrorCode PetscSpaceGetDimension_Polynomial(PetscSpace sp, PetscInt *dim)
8320cf1dd8SToby Isaac {
8420cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
8520cf1dd8SToby Isaac   PetscInt         deg  = sp->degree;
8620cf1dd8SToby Isaac   PetscInt         n    = sp->Nv, i;
8720cf1dd8SToby Isaac   PetscReal        D    = 1.0;
8820cf1dd8SToby Isaac 
8920cf1dd8SToby Isaac   PetscFunctionBegin;
9020cf1dd8SToby Isaac   if (poly->tensor) {
9120cf1dd8SToby Isaac     *dim = 1;
9220cf1dd8SToby Isaac     for (i = 0; i < n; ++i) *dim *= (deg+1);
9320cf1dd8SToby Isaac   } else {
9420cf1dd8SToby Isaac     for (i = 1; i <= n; ++i) {
9520cf1dd8SToby Isaac       D *= ((PetscReal) (deg+i))/i;
9620cf1dd8SToby Isaac     }
9720cf1dd8SToby Isaac     *dim = (PetscInt) (D + 0.5);
9820cf1dd8SToby Isaac   }
9920cf1dd8SToby Isaac   *dim *= sp->Nc;
10020cf1dd8SToby Isaac   PetscFunctionReturn(0);
10120cf1dd8SToby Isaac }
10220cf1dd8SToby Isaac 
10320cf1dd8SToby Isaac /*
10420cf1dd8SToby Isaac   LatticePoint_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to 'sum'.
10520cf1dd8SToby Isaac 
10620cf1dd8SToby Isaac   Input Parameters:
10720cf1dd8SToby Isaac + len - The length of the tuple
10820cf1dd8SToby Isaac . sum - The sum of all entries in the tuple
10920cf1dd8SToby Isaac - ind - The current multi-index of the tuple, initialized to the 0 tuple
11020cf1dd8SToby Isaac 
11120cf1dd8SToby Isaac   Output Parameter:
11220cf1dd8SToby Isaac + ind - The multi-index of the tuple, -1 indicates the iteration has terminated
11320cf1dd8SToby Isaac . tup - A tuple of len integers addig to sum
11420cf1dd8SToby Isaac 
11520cf1dd8SToby Isaac   Level: developer
11620cf1dd8SToby Isaac 
11720cf1dd8SToby Isaac .seealso:
11820cf1dd8SToby Isaac */
11920cf1dd8SToby Isaac static PetscErrorCode LatticePoint_Internal(PetscInt len, PetscInt sum, PetscInt ind[], PetscInt tup[])
12020cf1dd8SToby Isaac {
12120cf1dd8SToby Isaac   PetscInt       i;
12220cf1dd8SToby Isaac   PetscErrorCode ierr;
12320cf1dd8SToby Isaac 
12420cf1dd8SToby Isaac   PetscFunctionBegin;
12520cf1dd8SToby Isaac   if (len == 1) {
12620cf1dd8SToby Isaac     ind[0] = -1;
12720cf1dd8SToby Isaac     tup[0] = sum;
12820cf1dd8SToby Isaac   } else if (sum == 0) {
12920cf1dd8SToby Isaac     for (i = 0; i < len; ++i) {ind[0] = -1; tup[i] = 0;}
13020cf1dd8SToby Isaac   } else {
13120cf1dd8SToby Isaac     tup[0] = sum - ind[0];
13220cf1dd8SToby Isaac     ierr = LatticePoint_Internal(len-1, ind[0], &ind[1], &tup[1]);CHKERRQ(ierr);
13320cf1dd8SToby Isaac     if (ind[1] < 0) {
13420cf1dd8SToby Isaac       if (ind[0] == sum) {ind[0] = -1;}
13520cf1dd8SToby Isaac       else               {ind[1] = 0; ++ind[0];}
13620cf1dd8SToby Isaac     }
13720cf1dd8SToby Isaac   }
13820cf1dd8SToby Isaac   PetscFunctionReturn(0);
13920cf1dd8SToby Isaac }
14020cf1dd8SToby Isaac 
14120cf1dd8SToby Isaac /*
14220cf1dd8SToby Isaac   TensorPoint_Internal - Returns all tuples of size 'len' with nonnegative integers that are less than 'max'.
14320cf1dd8SToby Isaac 
14420cf1dd8SToby Isaac   Input Parameters:
14520cf1dd8SToby Isaac + len - The length of the tuple
14620cf1dd8SToby Isaac . max - The max for all entries in the tuple
14720cf1dd8SToby Isaac - ind - The current multi-index of the tuple, initialized to the 0 tuple
14820cf1dd8SToby Isaac 
14920cf1dd8SToby Isaac   Output Parameter:
15020cf1dd8SToby Isaac + ind - The multi-index of the tuple, -1 indicates the iteration has terminated
15120cf1dd8SToby Isaac . tup - A tuple of len integers less than max
15220cf1dd8SToby Isaac 
15320cf1dd8SToby Isaac   Level: developer
15420cf1dd8SToby Isaac 
15520cf1dd8SToby Isaac .seealso:
15620cf1dd8SToby Isaac */
15720cf1dd8SToby Isaac static PetscErrorCode TensorPoint_Internal(PetscInt len, PetscInt max, PetscInt ind[], PetscInt tup[])
15820cf1dd8SToby Isaac {
15920cf1dd8SToby Isaac   PetscInt       i;
16020cf1dd8SToby Isaac   PetscErrorCode ierr;
16120cf1dd8SToby Isaac 
16220cf1dd8SToby Isaac   PetscFunctionBegin;
16320cf1dd8SToby Isaac   if (len == 1) {
16420cf1dd8SToby Isaac     tup[0] = ind[0]++;
16520cf1dd8SToby Isaac     ind[0] = ind[0] >= max ? -1 : ind[0];
16620cf1dd8SToby Isaac   } else if (max == 0) {
16720cf1dd8SToby Isaac     for (i = 0; i < len; ++i) {ind[0] = -1; tup[i] = 0;}
16820cf1dd8SToby Isaac   } else {
16920cf1dd8SToby Isaac     tup[0] = ind[0];
17020cf1dd8SToby Isaac     ierr = TensorPoint_Internal(len-1, max, &ind[1], &tup[1]);CHKERRQ(ierr);
17120cf1dd8SToby Isaac     if (ind[1] < 0) {
17220cf1dd8SToby Isaac       ind[1] = 0;
17320cf1dd8SToby Isaac       if (ind[0] == max-1) {ind[0] = -1;}
17420cf1dd8SToby Isaac       else                 {++ind[0];}
17520cf1dd8SToby Isaac     }
17620cf1dd8SToby Isaac   }
17720cf1dd8SToby Isaac   PetscFunctionReturn(0);
17820cf1dd8SToby Isaac }
17920cf1dd8SToby Isaac 
18020cf1dd8SToby Isaac /*
18120cf1dd8SToby Isaac   p in [0, npoints), i in [0, pdim), c in [0, Nc)
18220cf1dd8SToby Isaac 
18320cf1dd8SToby Isaac   B[p][i][c] = B[p][i_scalar][c][c]
18420cf1dd8SToby Isaac */
18520cf1dd8SToby Isaac PetscErrorCode PetscSpaceEvaluate_Polynomial(PetscSpace sp, PetscInt npoints, const PetscReal points[], PetscReal B[], PetscReal D[], PetscReal H[])
18620cf1dd8SToby Isaac {
18720cf1dd8SToby Isaac   PetscSpace_Poly *poly    = (PetscSpace_Poly *) sp->data;
18820cf1dd8SToby Isaac   DM               dm      = sp->dm;
18920cf1dd8SToby Isaac   PetscInt         Nc      = sp->Nc;
19020cf1dd8SToby Isaac   PetscInt         ndegree = sp->degree+1;
19120cf1dd8SToby Isaac   PetscInt        *degrees = poly->degrees;
19220cf1dd8SToby Isaac   PetscInt         dim     = sp->Nv;
19320cf1dd8SToby Isaac   PetscReal       *lpoints, *tmp, *LB, *LD, *LH;
19420cf1dd8SToby Isaac   PetscInt        *ind, *tup;
19520cf1dd8SToby Isaac   PetscInt         c, pdim, d, e, der, der2, i, p, deg, o;
19620cf1dd8SToby Isaac   PetscErrorCode   ierr;
19720cf1dd8SToby Isaac 
19820cf1dd8SToby Isaac   PetscFunctionBegin;
19920cf1dd8SToby Isaac   ierr = PetscSpaceGetDimension(sp, &pdim);CHKERRQ(ierr);
20020cf1dd8SToby Isaac   pdim /= Nc;
20120cf1dd8SToby Isaac   ierr = DMGetWorkArray(dm, npoints, MPIU_REAL, &lpoints);CHKERRQ(ierr);
20220cf1dd8SToby Isaac   ierr = DMGetWorkArray(dm, npoints*ndegree*3, MPIU_REAL, &tmp);CHKERRQ(ierr);
2038c7916f5SToby Isaac   if (B || D || H) {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LB);CHKERRQ(ierr);}
2048c7916f5SToby Isaac   if (D || H)      {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LD);CHKERRQ(ierr);}
20520cf1dd8SToby Isaac   if (H)           {ierr = DMGetWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LH);CHKERRQ(ierr);}
20620cf1dd8SToby Isaac   for (d = 0; d < dim; ++d) {
20720cf1dd8SToby Isaac     for (p = 0; p < npoints; ++p) {
20820cf1dd8SToby Isaac       lpoints[p] = points[p*dim+d];
20920cf1dd8SToby Isaac     }
21020cf1dd8SToby Isaac     ierr = PetscDTLegendreEval(npoints, lpoints, ndegree, degrees, tmp, &tmp[1*npoints*ndegree], &tmp[2*npoints*ndegree]);CHKERRQ(ierr);
21120cf1dd8SToby Isaac     /* LB, LD, LH (ndegree * dim x npoints) */
21220cf1dd8SToby Isaac     for (deg = 0; deg < ndegree; ++deg) {
21320cf1dd8SToby Isaac       for (p = 0; p < npoints; ++p) {
2148c7916f5SToby Isaac         if (B || D || H) LB[(deg*dim + d)*npoints + p] = tmp[(0*npoints + p)*ndegree+deg];
2158c7916f5SToby Isaac         if (D || H)      LD[(deg*dim + d)*npoints + p] = tmp[(1*npoints + p)*ndegree+deg];
21620cf1dd8SToby Isaac         if (H)           LH[(deg*dim + d)*npoints + p] = tmp[(2*npoints + p)*ndegree+deg];
21720cf1dd8SToby Isaac       }
21820cf1dd8SToby Isaac     }
21920cf1dd8SToby Isaac   }
22020cf1dd8SToby Isaac   /* Multiply by A (pdim x ndegree * dim) */
22120cf1dd8SToby Isaac   ierr = PetscMalloc2(dim,&ind,dim,&tup);CHKERRQ(ierr);
22220cf1dd8SToby Isaac   if (B) {
22320cf1dd8SToby Isaac     /* B (npoints x pdim x Nc) */
22420cf1dd8SToby Isaac     ierr = PetscMemzero(B, npoints*pdim*Nc*Nc * sizeof(PetscReal));CHKERRQ(ierr);
22520cf1dd8SToby Isaac     if (poly->tensor) {
22620cf1dd8SToby Isaac       i = 0;
22720cf1dd8SToby Isaac       ierr = PetscMemzero(ind, dim * sizeof(PetscInt));CHKERRQ(ierr);
22820cf1dd8SToby Isaac       while (ind[0] >= 0) {
22920cf1dd8SToby Isaac         ierr = TensorPoint_Internal(dim, sp->degree+1, ind, tup);CHKERRQ(ierr);
23020cf1dd8SToby Isaac         for (p = 0; p < npoints; ++p) {
23120cf1dd8SToby Isaac           B[(p*pdim + i)*Nc*Nc] = 1.0;
23220cf1dd8SToby Isaac           for (d = 0; d < dim; ++d) {
23320cf1dd8SToby Isaac             B[(p*pdim + i)*Nc*Nc] *= LB[(tup[d]*dim + d)*npoints + p];
23420cf1dd8SToby Isaac           }
23520cf1dd8SToby Isaac         }
23620cf1dd8SToby Isaac         ++i;
23720cf1dd8SToby Isaac       }
23820cf1dd8SToby Isaac     } else {
23920cf1dd8SToby Isaac       i = 0;
24020cf1dd8SToby Isaac       for (o = 0; o <= sp->degree; ++o) {
24120cf1dd8SToby Isaac         ierr = PetscMemzero(ind, dim * sizeof(PetscInt));CHKERRQ(ierr);
24220cf1dd8SToby Isaac         while (ind[0] >= 0) {
24320cf1dd8SToby Isaac           ierr = LatticePoint_Internal(dim, o, 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       }
25320cf1dd8SToby Isaac     }
25420cf1dd8SToby Isaac     /* Make direct sum basis for multicomponent space */
25520cf1dd8SToby Isaac     for (p = 0; p < npoints; ++p) {
25620cf1dd8SToby Isaac       for (i = 0; i < pdim; ++i) {
25720cf1dd8SToby Isaac         for (c = 1; c < Nc; ++c) {
25820cf1dd8SToby Isaac           B[(p*pdim*Nc + i*Nc + c)*Nc + c] = B[(p*pdim + i)*Nc*Nc];
25920cf1dd8SToby Isaac         }
26020cf1dd8SToby Isaac       }
26120cf1dd8SToby Isaac     }
26220cf1dd8SToby Isaac   }
26320cf1dd8SToby Isaac   if (D) {
26420cf1dd8SToby Isaac     /* D (npoints x pdim x Nc x dim) */
26520cf1dd8SToby Isaac     ierr = PetscMemzero(D, npoints*pdim*Nc*Nc*dim * sizeof(PetscReal));CHKERRQ(ierr);
26620cf1dd8SToby Isaac     if (poly->tensor) {
26720cf1dd8SToby Isaac       i = 0;
26820cf1dd8SToby Isaac       ierr = PetscMemzero(ind, dim * sizeof(PetscInt));CHKERRQ(ierr);
26920cf1dd8SToby Isaac       while (ind[0] >= 0) {
27020cf1dd8SToby Isaac         ierr = TensorPoint_Internal(dim, sp->degree+1, ind, tup);CHKERRQ(ierr);
27120cf1dd8SToby Isaac         for (p = 0; p < npoints; ++p) {
27220cf1dd8SToby Isaac           for (der = 0; der < dim; ++der) {
27320cf1dd8SToby Isaac             D[(p*pdim + i)*Nc*Nc*dim + der] = 1.0;
27420cf1dd8SToby Isaac             for (d = 0; d < dim; ++d) {
27520cf1dd8SToby Isaac               if (d == der) {
27620cf1dd8SToby Isaac                 D[(p*pdim + i)*Nc*Nc*dim + der] *= LD[(tup[d]*dim + d)*npoints + p];
27720cf1dd8SToby Isaac               } else {
27820cf1dd8SToby Isaac                 D[(p*pdim + i)*Nc*Nc*dim + der] *= LB[(tup[d]*dim + d)*npoints + p];
27920cf1dd8SToby Isaac               }
28020cf1dd8SToby Isaac             }
28120cf1dd8SToby Isaac           }
28220cf1dd8SToby Isaac         }
28320cf1dd8SToby Isaac         ++i;
28420cf1dd8SToby Isaac       }
28520cf1dd8SToby Isaac     } else {
28620cf1dd8SToby Isaac       i = 0;
28720cf1dd8SToby Isaac       for (o = 0; o <= sp->degree; ++o) {
28820cf1dd8SToby Isaac         ierr = PetscMemzero(ind, dim * sizeof(PetscInt));CHKERRQ(ierr);
28920cf1dd8SToby Isaac         while (ind[0] >= 0) {
29020cf1dd8SToby Isaac           ierr = LatticePoint_Internal(dim, o, ind, tup);CHKERRQ(ierr);
29120cf1dd8SToby Isaac           for (p = 0; p < npoints; ++p) {
29220cf1dd8SToby Isaac             for (der = 0; der < dim; ++der) {
29320cf1dd8SToby Isaac               D[(p*pdim + i)*Nc*Nc*dim + der] = 1.0;
29420cf1dd8SToby Isaac               for (d = 0; d < dim; ++d) {
29520cf1dd8SToby Isaac                 if (d == der) {
29620cf1dd8SToby Isaac                   D[(p*pdim + i)*Nc*Nc*dim + der] *= LD[(tup[d]*dim + d)*npoints + p];
29720cf1dd8SToby Isaac                 } else {
29820cf1dd8SToby Isaac                   D[(p*pdim + i)*Nc*Nc*dim + der] *= LB[(tup[d]*dim + d)*npoints + p];
29920cf1dd8SToby Isaac                 }
30020cf1dd8SToby Isaac               }
30120cf1dd8SToby Isaac             }
30220cf1dd8SToby Isaac           }
30320cf1dd8SToby Isaac           ++i;
30420cf1dd8SToby Isaac         }
30520cf1dd8SToby Isaac       }
30620cf1dd8SToby Isaac     }
30720cf1dd8SToby Isaac     /* Make direct sum basis for multicomponent space */
30820cf1dd8SToby Isaac     for (p = 0; p < npoints; ++p) {
30920cf1dd8SToby Isaac       for (i = 0; i < pdim; ++i) {
31020cf1dd8SToby Isaac         for (c = 1; c < Nc; ++c) {
31120cf1dd8SToby Isaac           for (d = 0; d < dim; ++d) {
31220cf1dd8SToby Isaac             D[((p*pdim*Nc + i*Nc + c)*Nc + c)*dim + d] = D[(p*pdim + i)*Nc*Nc*dim + d];
31320cf1dd8SToby Isaac           }
31420cf1dd8SToby Isaac         }
31520cf1dd8SToby Isaac       }
31620cf1dd8SToby Isaac     }
31720cf1dd8SToby Isaac   }
31820cf1dd8SToby Isaac   if (H) {
31920cf1dd8SToby Isaac     /* H (npoints x pdim x Nc x Nc x dim x dim) */
32020cf1dd8SToby Isaac     ierr = PetscMemzero(H, npoints*pdim*Nc*Nc*dim*dim * sizeof(PetscReal));CHKERRQ(ierr);
32120cf1dd8SToby Isaac     if (poly->tensor) {
32220cf1dd8SToby Isaac       i = 0;
32320cf1dd8SToby Isaac       ierr = PetscMemzero(ind, dim * sizeof(PetscInt));CHKERRQ(ierr);
32420cf1dd8SToby Isaac       while (ind[0] >= 0) {
32520cf1dd8SToby Isaac         ierr = TensorPoint_Internal(dim, sp->degree+1, ind, tup);CHKERRQ(ierr);
32620cf1dd8SToby Isaac         for (p = 0; p < npoints; ++p) {
32720cf1dd8SToby Isaac           for (der = 0; der < dim; ++der) {
32820cf1dd8SToby Isaac             H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der] = 1.0;
32920cf1dd8SToby Isaac             for (d = 0; d < dim; ++d) {
33020cf1dd8SToby Isaac               if (d == der) {
33120cf1dd8SToby Isaac                 H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der] *= LH[(tup[d]*dim + d)*npoints + p];
33220cf1dd8SToby Isaac               } else {
33320cf1dd8SToby Isaac                 H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der] *= LB[(tup[d]*dim + d)*npoints + p];
33420cf1dd8SToby Isaac               }
33520cf1dd8SToby Isaac             }
33620cf1dd8SToby Isaac             for (der2 = der + 1; der2 < dim; ++der2) {
33720cf1dd8SToby Isaac               H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] = 1.0;
33820cf1dd8SToby Isaac               for (d = 0; d < dim; ++d) {
33920cf1dd8SToby Isaac                 if (d == der || d == der2) {
34020cf1dd8SToby Isaac                   H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] *= LD[(tup[d]*dim + d)*npoints + p];
34120cf1dd8SToby Isaac                 } else {
34220cf1dd8SToby Isaac                   H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] *= LB[(tup[d]*dim + d)*npoints + p];
34320cf1dd8SToby Isaac                 }
34420cf1dd8SToby Isaac               }
34520cf1dd8SToby Isaac               H[((p*pdim + i)*Nc*Nc*dim + der2) * dim + der] = H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2];
34620cf1dd8SToby Isaac             }
34720cf1dd8SToby Isaac           }
34820cf1dd8SToby Isaac         }
34920cf1dd8SToby Isaac         ++i;
35020cf1dd8SToby Isaac       }
35120cf1dd8SToby Isaac     } else {
35220cf1dd8SToby Isaac       i = 0;
35320cf1dd8SToby Isaac       for (o = 0; o <= sp->degree; ++o) {
35420cf1dd8SToby Isaac         ierr = PetscMemzero(ind, dim * sizeof(PetscInt));CHKERRQ(ierr);
35520cf1dd8SToby Isaac         while (ind[0] >= 0) {
35620cf1dd8SToby Isaac           ierr = LatticePoint_Internal(dim, o, ind, tup);CHKERRQ(ierr);
35720cf1dd8SToby Isaac           for (p = 0; p < npoints; ++p) {
35820cf1dd8SToby Isaac             for (der = 0; der < dim; ++der) {
35920cf1dd8SToby Isaac               H[((p*pdim + i)*Nc*Nc*dim + der)*dim + der] = 1.0;
36020cf1dd8SToby Isaac               for (d = 0; d < dim; ++d) {
36120cf1dd8SToby Isaac                 if (d == der) {
36220cf1dd8SToby Isaac                   H[((p*pdim + i)*Nc*Nc*dim + der)*dim + der] *= LH[(tup[d]*dim + d)*npoints + p];
36320cf1dd8SToby Isaac                 } else {
36420cf1dd8SToby Isaac                   H[((p*pdim + i)*Nc*Nc*dim + der)*dim + der] *= LB[(tup[d]*dim + d)*npoints + p];
36520cf1dd8SToby Isaac                 }
36620cf1dd8SToby Isaac               }
36720cf1dd8SToby Isaac               for (der2 = der + 1; der2 < dim; ++der2) {
36820cf1dd8SToby Isaac                 H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] = 1.0;
36920cf1dd8SToby Isaac                 for (d = 0; d < dim; ++d) {
37020cf1dd8SToby Isaac                   if (d == der || d == der2) {
37120cf1dd8SToby Isaac                     H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] *= LD[(tup[d]*dim + d)*npoints + p];
37220cf1dd8SToby Isaac                   } else {
37320cf1dd8SToby Isaac                     H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2] *= LB[(tup[d]*dim + d)*npoints + p];
37420cf1dd8SToby Isaac                   }
37520cf1dd8SToby Isaac                 }
37620cf1dd8SToby Isaac                 H[((p*pdim + i)*Nc*Nc*dim + der2) * dim + der] = H[((p*pdim + i)*Nc*Nc*dim + der) * dim + der2];
37720cf1dd8SToby Isaac               }
37820cf1dd8SToby Isaac             }
37920cf1dd8SToby Isaac           }
38020cf1dd8SToby Isaac           ++i;
38120cf1dd8SToby Isaac         }
38220cf1dd8SToby Isaac       }
38320cf1dd8SToby Isaac     }
38420cf1dd8SToby Isaac     /* Make direct sum basis for multicomponent space */
38520cf1dd8SToby Isaac     for (p = 0; p < npoints; ++p) {
38620cf1dd8SToby Isaac       for (i = 0; i < pdim; ++i) {
38720cf1dd8SToby Isaac         for (c = 1; c < Nc; ++c) {
38820cf1dd8SToby Isaac           for (d = 0; d < dim; ++d) {
38920cf1dd8SToby Isaac             for (e = 0; e < dim; ++e) {
39020cf1dd8SToby Isaac               H[(((p*pdim*Nc + i*Nc + c)*Nc + c)*dim + d)*dim + e] = H[((p*pdim + i)*Nc*Nc*dim + d)*dim + e];
39120cf1dd8SToby Isaac             }
39220cf1dd8SToby Isaac           }
39320cf1dd8SToby Isaac         }
39420cf1dd8SToby Isaac       }
39520cf1dd8SToby Isaac     }
39620cf1dd8SToby Isaac   }
39720cf1dd8SToby Isaac   ierr = PetscFree2(ind,tup);CHKERRQ(ierr);
39820cf1dd8SToby Isaac   if (H)           {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LH);CHKERRQ(ierr);}
3998c7916f5SToby Isaac   if (D || H)      {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LD);CHKERRQ(ierr);}
4008c7916f5SToby Isaac   if (B || D || H) {ierr = DMRestoreWorkArray(dm, npoints*dim*ndegree, MPIU_REAL, &LB);CHKERRQ(ierr);}
40120cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(dm, npoints*ndegree*3, MPIU_REAL, &tmp);CHKERRQ(ierr);
40220cf1dd8SToby Isaac   ierr = DMRestoreWorkArray(dm, npoints, MPIU_REAL, &lpoints);CHKERRQ(ierr);
40320cf1dd8SToby Isaac   PetscFunctionReturn(0);
40420cf1dd8SToby Isaac }
40520cf1dd8SToby Isaac 
40620cf1dd8SToby Isaac /*@
40720cf1dd8SToby Isaac   PetscSpacePolynomialSetTensor - Set whether a function space is a space of tensor polynomials (the space is spanned
40820cf1dd8SToby Isaac   by polynomials whose degree in each variabl is bounded by the given order), as opposed to polynomials (the space is
40920cf1dd8SToby Isaac   spanned by polynomials whose total degree---summing over all variables---is bounded by the given order).
41020cf1dd8SToby Isaac 
41120cf1dd8SToby Isaac   Input Parameters:
41220cf1dd8SToby Isaac + sp     - the function space object
41320cf1dd8SToby Isaac - tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space
41420cf1dd8SToby Isaac 
415*4ab77754SMatthew G. Knepley   Options Database:
416*4ab77754SMatthew G. Knepley . -petscspace_poly_tensor <bool> - Whether to use tensor product polynomials in higher dimension
417*4ab77754SMatthew G. Knepley 
41820cf1dd8SToby Isaac   Level: beginner
41920cf1dd8SToby Isaac 
42020cf1dd8SToby Isaac .seealso: PetscSpacePolynomialGetTensor(), PetscSpaceSetDegree(), PetscSpaceSetNumVariables()
42120cf1dd8SToby Isaac @*/
42220cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace sp, PetscBool tensor)
42320cf1dd8SToby Isaac {
42420cf1dd8SToby Isaac   PetscErrorCode ierr;
42520cf1dd8SToby Isaac 
42620cf1dd8SToby Isaac   PetscFunctionBegin;
42720cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
42820cf1dd8SToby Isaac   ierr = PetscTryMethod(sp,"PetscSpacePolynomialSetTensor_C",(PetscSpace,PetscBool),(sp,tensor));CHKERRQ(ierr);
42920cf1dd8SToby Isaac   PetscFunctionReturn(0);
43020cf1dd8SToby Isaac }
43120cf1dd8SToby Isaac 
43220cf1dd8SToby Isaac /*@
43320cf1dd8SToby Isaac   PetscSpacePolynomialGetTensor - Get whether a function space is a space of tensor polynomials (the space is spanned
43420cf1dd8SToby Isaac   by polynomials whose degree in each variabl is bounded by the given order), as opposed to polynomials (the space is
43520cf1dd8SToby Isaac   spanned by polynomials whose total degree---summing over all variables---is bounded by the given order).
43620cf1dd8SToby Isaac 
43720cf1dd8SToby Isaac   Input Parameters:
43820cf1dd8SToby Isaac . sp     - the function space object
43920cf1dd8SToby Isaac 
44020cf1dd8SToby Isaac   Output Parameters:
44120cf1dd8SToby Isaac . tensor - PETSC_TRUE for a tensor polynomial space, PETSC_FALSE for a polynomial space
44220cf1dd8SToby Isaac 
44320cf1dd8SToby Isaac   Level: beginner
44420cf1dd8SToby Isaac 
44520cf1dd8SToby Isaac .seealso: PetscSpacePolynomialSetTensor(), PetscSpaceSetDegree(), PetscSpaceSetNumVariables()
44620cf1dd8SToby Isaac @*/
44720cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace sp, PetscBool *tensor)
44820cf1dd8SToby Isaac {
44920cf1dd8SToby Isaac   PetscErrorCode ierr;
45020cf1dd8SToby Isaac 
45120cf1dd8SToby Isaac   PetscFunctionBegin;
45220cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
45320cf1dd8SToby Isaac   PetscValidPointer(tensor, 2);
45420cf1dd8SToby Isaac   ierr = PetscTryMethod(sp,"PetscSpacePolynomialGetTensor_C",(PetscSpace,PetscBool*),(sp,tensor));CHKERRQ(ierr);
45520cf1dd8SToby Isaac   PetscFunctionReturn(0);
45620cf1dd8SToby Isaac }
45720cf1dd8SToby Isaac 
45820cf1dd8SToby Isaac static PetscErrorCode PetscSpacePolynomialSetTensor_Polynomial(PetscSpace sp, PetscBool tensor)
45920cf1dd8SToby Isaac {
46020cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
46120cf1dd8SToby Isaac 
46220cf1dd8SToby Isaac   PetscFunctionBegin;
46320cf1dd8SToby Isaac   poly->tensor = tensor;
46420cf1dd8SToby Isaac   PetscFunctionReturn(0);
46520cf1dd8SToby Isaac }
46620cf1dd8SToby Isaac 
46720cf1dd8SToby Isaac static PetscErrorCode PetscSpacePolynomialGetTensor_Polynomial(PetscSpace sp, PetscBool *tensor)
46820cf1dd8SToby Isaac {
46920cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
47020cf1dd8SToby Isaac 
47120cf1dd8SToby Isaac   PetscFunctionBegin;
47220cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
47320cf1dd8SToby Isaac   PetscValidPointer(tensor, 2);
47420cf1dd8SToby Isaac   *tensor = poly->tensor;
47520cf1dd8SToby Isaac   PetscFunctionReturn(0);
47620cf1dd8SToby Isaac }
47720cf1dd8SToby Isaac 
47820cf1dd8SToby Isaac static PetscErrorCode PetscSpaceGetHeightSubspace_Polynomial(PetscSpace sp, PetscInt height, PetscSpace *subsp)
47920cf1dd8SToby Isaac {
48020cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
48120cf1dd8SToby Isaac   PetscInt         Nc, dim, order;
48220cf1dd8SToby Isaac   PetscBool        tensor;
48320cf1dd8SToby Isaac   PetscErrorCode   ierr;
48420cf1dd8SToby Isaac 
48520cf1dd8SToby Isaac   PetscFunctionBegin;
48620cf1dd8SToby Isaac   ierr = PetscSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr);
48720cf1dd8SToby Isaac   ierr = PetscSpaceGetNumVariables(sp, &dim);CHKERRQ(ierr);
48820cf1dd8SToby Isaac   ierr = PetscSpaceGetDegree(sp, &order, NULL);CHKERRQ(ierr);
48920cf1dd8SToby Isaac   ierr = PetscSpacePolynomialGetTensor(sp, &tensor);CHKERRQ(ierr);
49020cf1dd8SToby 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);}
49120cf1dd8SToby Isaac   if (!poly->subspaces) {ierr = PetscCalloc1(dim, &poly->subspaces);CHKERRQ(ierr);}
49220cf1dd8SToby Isaac   if (height <= dim) {
49320cf1dd8SToby Isaac     if (!poly->subspaces[height-1]) {
49420cf1dd8SToby Isaac       PetscSpace sub;
49520cf1dd8SToby Isaac 
49620cf1dd8SToby Isaac       ierr = PetscSpaceCreate(PetscObjectComm((PetscObject) sp), &sub);CHKERRQ(ierr);
49720cf1dd8SToby Isaac       ierr = PetscSpaceSetType(sub, PETSCSPACEPOLYNOMIAL);CHKERRQ(ierr);
498d39dd5f5SToby Isaac       ierr = PetscSpaceSetNumComponents(sub, Nc);CHKERRQ(ierr);
499d39dd5f5SToby Isaac       ierr = PetscSpaceSetDegree(sub, order, PETSC_DETERMINE);CHKERRQ(ierr);
50020cf1dd8SToby Isaac       ierr = PetscSpaceSetNumVariables(sub, dim-height);CHKERRQ(ierr);
50120cf1dd8SToby Isaac       ierr = PetscSpacePolynomialSetTensor(sub, tensor);CHKERRQ(ierr);
50220cf1dd8SToby Isaac       ierr = PetscSpaceSetUp(sub);CHKERRQ(ierr);
50320cf1dd8SToby Isaac       poly->subspaces[height-1] = sub;
50420cf1dd8SToby Isaac     }
50520cf1dd8SToby Isaac     *subsp = poly->subspaces[height-1];
50620cf1dd8SToby Isaac   } else {
50720cf1dd8SToby Isaac     *subsp = NULL;
50820cf1dd8SToby Isaac   }
50920cf1dd8SToby Isaac   PetscFunctionReturn(0);
51020cf1dd8SToby Isaac }
51120cf1dd8SToby Isaac 
51220cf1dd8SToby Isaac PetscErrorCode PetscSpaceInitialize_Polynomial(PetscSpace sp)
51320cf1dd8SToby Isaac {
51420cf1dd8SToby Isaac   PetscErrorCode ierr;
51520cf1dd8SToby Isaac 
51620cf1dd8SToby Isaac   PetscFunctionBegin;
51720cf1dd8SToby Isaac   sp->ops->setfromoptions    = PetscSpaceSetFromOptions_Polynomial;
51820cf1dd8SToby Isaac   sp->ops->setup             = PetscSpaceSetUp_Polynomial;
51920cf1dd8SToby Isaac   sp->ops->view              = PetscSpaceView_Polynomial;
52020cf1dd8SToby Isaac   sp->ops->destroy           = PetscSpaceDestroy_Polynomial;
52120cf1dd8SToby Isaac   sp->ops->getdimension      = PetscSpaceGetDimension_Polynomial;
52220cf1dd8SToby Isaac   sp->ops->evaluate          = PetscSpaceEvaluate_Polynomial;
52320cf1dd8SToby Isaac   sp->ops->getheightsubspace = PetscSpaceGetHeightSubspace_Polynomial;
52420cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialGetTensor_C", PetscSpacePolynomialGetTensor_Polynomial);CHKERRQ(ierr);
52520cf1dd8SToby Isaac   ierr = PetscObjectComposeFunction((PetscObject) sp, "PetscSpacePolynomialSetTensor_C", PetscSpacePolynomialSetTensor_Polynomial);CHKERRQ(ierr);
52620cf1dd8SToby Isaac   PetscFunctionReturn(0);
52720cf1dd8SToby Isaac }
52820cf1dd8SToby Isaac 
52920cf1dd8SToby Isaac /*MC
53020cf1dd8SToby Isaac   PETSCSPACEPOLYNOMIAL = "poly" - A PetscSpace object that encapsulates a polynomial space, e.g. P1 is the space of
53120cf1dd8SToby Isaac   linear polynomials. The space is replicated for each component.
53220cf1dd8SToby Isaac 
53320cf1dd8SToby Isaac   Level: intermediate
53420cf1dd8SToby Isaac 
53520cf1dd8SToby Isaac .seealso: PetscSpaceType, PetscSpaceCreate(), PetscSpaceSetType()
53620cf1dd8SToby Isaac M*/
53720cf1dd8SToby Isaac 
53820cf1dd8SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceCreate_Polynomial(PetscSpace sp)
53920cf1dd8SToby Isaac {
54020cf1dd8SToby Isaac   PetscSpace_Poly *poly;
54120cf1dd8SToby Isaac   PetscErrorCode   ierr;
54220cf1dd8SToby Isaac 
54320cf1dd8SToby Isaac   PetscFunctionBegin;
54420cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
54520cf1dd8SToby Isaac   ierr     = PetscNewLog(sp,&poly);CHKERRQ(ierr);
54620cf1dd8SToby Isaac   sp->data = poly;
54720cf1dd8SToby Isaac 
54820cf1dd8SToby Isaac   poly->symmetric    = PETSC_FALSE;
54920cf1dd8SToby Isaac   poly->tensor       = PETSC_FALSE;
55020cf1dd8SToby Isaac   poly->degrees      = NULL;
55120cf1dd8SToby Isaac   poly->subspaces    = NULL;
55220cf1dd8SToby Isaac 
55320cf1dd8SToby Isaac   ierr = PetscSpaceInitialize_Polynomial(sp);CHKERRQ(ierr);
55420cf1dd8SToby Isaac   PetscFunctionReturn(0);
55520cf1dd8SToby Isaac }
55620cf1dd8SToby Isaac 
55720cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace sp, PetscBool sym)
55820cf1dd8SToby Isaac {
55920cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
56020cf1dd8SToby Isaac 
56120cf1dd8SToby Isaac   PetscFunctionBegin;
56220cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
56320cf1dd8SToby Isaac   poly->symmetric = sym;
56420cf1dd8SToby Isaac   PetscFunctionReturn(0);
56520cf1dd8SToby Isaac }
56620cf1dd8SToby Isaac 
56720cf1dd8SToby Isaac PetscErrorCode PetscSpacePolynomialGetSymmetric(PetscSpace sp, PetscBool *sym)
56820cf1dd8SToby Isaac {
56920cf1dd8SToby Isaac   PetscSpace_Poly *poly = (PetscSpace_Poly *) sp->data;
57020cf1dd8SToby Isaac 
57120cf1dd8SToby Isaac   PetscFunctionBegin;
57220cf1dd8SToby Isaac   PetscValidHeaderSpecific(sp, PETSCSPACE_CLASSID, 1);
57320cf1dd8SToby Isaac   PetscValidPointer(sym, 2);
57420cf1dd8SToby Isaac   *sym = poly->symmetric;
57520cf1dd8SToby Isaac   PetscFunctionReturn(0);
57620cf1dd8SToby Isaac }
57720cf1dd8SToby Isaac 
578