xref: /petsc/src/dm/dt/tests/ex4.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown static char help[] = "Tests dual space symmetry.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscfe.h>
4c4762a1bSJed Brown #include <petscdmplex.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown static PetscErrorCode CheckSymmetry(PetscInt dim, PetscInt order, PetscBool tensor)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   DM                dm;
9c4762a1bSJed Brown   PetscDualSpace    sp;
10c4762a1bSJed Brown   PetscInt          nFunc, *ids, *idsCopy, *idsCopy2, i, closureSize, *closure = NULL, offset, depth;
11c4762a1bSJed Brown   DMLabel           depthLabel;
12c4762a1bSJed Brown   PetscBool         printed = PETSC_FALSE;
13c4762a1bSJed Brown   PetscScalar       *vals, *valsCopy, *valsCopy2;
14c4762a1bSJed Brown   const PetscInt    *numDofs;
15c4762a1bSJed Brown   const PetscInt    ***perms = NULL;
16c4762a1bSJed Brown   const PetscScalar ***flips = NULL;
17c4762a1bSJed Brown 
18c4762a1bSJed Brown   PetscFunctionBegin;
19*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceCreate(PETSC_COMM_SELF,&sp));
20*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, tensor ? PETSC_FALSE : PETSC_TRUE), &dm));
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceSetType(sp,PETSCDUALSPACELAGRANGE));
22*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceSetDM(sp,dm));
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceSetOrder(sp,order));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceLagrangeSetContinuity(sp,PETSC_TRUE));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceLagrangeSetTensor(sp,tensor));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceSetFromOptions(sp));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceSetUp(sp));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceGetDimension(sp,&nFunc));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceGetSymmetries(sp,&perms,&flips));
30c4762a1bSJed Brown   if (!perms && !flips) {
31*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDualSpaceDestroy(&sp));
32*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDestroy(&dm));
33c4762a1bSJed Brown     PetscFunctionReturn(0);
34c4762a1bSJed Brown   }
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc6(nFunc,&ids,nFunc,&idsCopy,nFunc,&idsCopy2,nFunc*dim,&vals,nFunc*dim,&valsCopy,nFunc*dim,&valsCopy2));
36c4762a1bSJed Brown   for (i = 0; i < nFunc; i++) ids[i] = idsCopy2[i] = i;
37c4762a1bSJed Brown   for (i = 0; i < nFunc; i++) {
38c4762a1bSJed Brown     PetscQuadrature q;
39c4762a1bSJed Brown     PetscInt        numPoints, Nc, j;
40c4762a1bSJed Brown     const PetscReal *points;
41c4762a1bSJed Brown     const PetscReal *weights;
42c4762a1bSJed Brown 
43*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDualSpaceGetFunctional(sp,i,&q));
44*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscQuadratureGetData(q,NULL,&Nc,&numPoints,&points,&weights));
452c71b3e2SJacob Faibussowitsch     PetscCheckFalse(Nc != 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"Only support scalar quadrature, not %D components",Nc);
46c4762a1bSJed Brown     for (j = 0; j < dim; j++) vals[dim * i + j] = valsCopy2[dim * i + j] = (PetscScalar) points[j];
47c4762a1bSJed Brown   }
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceGetNumDof(sp,&numDofs));
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepth(dm,&depth));
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetTransitiveClosure(dm,0,PETSC_TRUE,&closureSize,&closure));
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexGetDepthLabel(dm,&depthLabel));
52c4762a1bSJed Brown   for (i = 0, offset = 0; i < closureSize; i++, offset += numDofs[depth]) {
53b5a892a1SMatthew G. Knepley     PetscInt          point = closure[2 * i], numFaces, j;
54c4762a1bSJed Brown     const PetscInt    **pointPerms = perms ? perms[i] : NULL;
55c4762a1bSJed Brown     const PetscScalar **pointFlips = flips ? flips[i] : NULL;
56c4762a1bSJed Brown     PetscBool         anyPrinted = PETSC_FALSE;
57c4762a1bSJed Brown 
58c4762a1bSJed Brown     if (!pointPerms && !pointFlips) continue;
59*5f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLabelGetValue(depthLabel,point,&depth));
60b5a892a1SMatthew G. Knepley     {
61b5a892a1SMatthew G. Knepley       DMPolytopeType ct;
62b5a892a1SMatthew G. Knepley       /* The number of arrangements is no longer based on the number of faces */
63*5f80ce2aSJacob Faibussowitsch       CHKERRQ(DMPlexGetCellType(dm, point, &ct));
64b5a892a1SMatthew G. Knepley       numFaces = DMPolytopeTypeGetNumArrangments(ct) / 2;
65b5a892a1SMatthew G. Knepley     }
66b5a892a1SMatthew G. Knepley     for (j = -numFaces; j < numFaces; j++) {
67c4762a1bSJed Brown       PetscInt          k, l;
68c4762a1bSJed Brown       const PetscInt    *perm = pointPerms ? pointPerms[j] : NULL;
69c4762a1bSJed Brown       const PetscScalar *flip = pointFlips ? pointFlips[j] : NULL;
70c4762a1bSJed Brown 
71c4762a1bSJed Brown       for (k = 0; k < numDofs[depth]; k++) {
72c4762a1bSJed Brown         PetscInt kLocal = perm ? perm[k] : k;
73c4762a1bSJed Brown 
74c4762a1bSJed Brown         idsCopy[kLocal] = ids[offset + k];
75c4762a1bSJed Brown         for (l = 0; l < dim; l++) {
76c4762a1bSJed Brown           valsCopy[kLocal * dim + l] = vals[(offset + k) * dim + l] * (flip ? flip[kLocal] : 1.);
77c4762a1bSJed Brown         }
78c4762a1bSJed Brown       }
79c4762a1bSJed Brown       if (!printed && numDofs[depth] > 1) {
80c4762a1bSJed Brown         IS   is;
81c4762a1bSJed Brown         Vec  vec;
82c4762a1bSJed Brown         char name[256];
83c4762a1bSJed Brown 
84c4762a1bSJed Brown         anyPrinted = PETSC_TRUE;
85*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscSNPrintf(name,256,"%DD, %s, Order %D, Point %D Symmetry %D",dim,tensor ? "Tensor" : "Simplex", order, point,j));
86*5f80ce2aSJacob Faibussowitsch         CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,numDofs[depth],idsCopy,PETSC_USE_POINTER,&is));
87*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscObjectSetName((PetscObject)is,name));
88*5f80ce2aSJacob Faibussowitsch         CHKERRQ(ISView(is,PETSC_VIEWER_STDOUT_SELF));
89*5f80ce2aSJacob Faibussowitsch         CHKERRQ(ISDestroy(&is));
90*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecCreateSeqWithArray(PETSC_COMM_SELF,dim,numDofs[depth]*dim,valsCopy,&vec));
91*5f80ce2aSJacob Faibussowitsch         CHKERRQ(PetscObjectSetName((PetscObject)vec,name));
92*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecView(vec,PETSC_VIEWER_STDOUT_SELF));
93*5f80ce2aSJacob Faibussowitsch         CHKERRQ(VecDestroy(&vec));
94c4762a1bSJed Brown       }
95c4762a1bSJed Brown       for (k = 0; k < numDofs[depth]; k++) {
96c4762a1bSJed Brown         PetscInt kLocal = perm ? perm[k] : k;
97c4762a1bSJed Brown 
98c4762a1bSJed Brown         idsCopy2[offset + k] = idsCopy[kLocal];
99c4762a1bSJed Brown         for (l = 0; l < dim; l++) {
100c4762a1bSJed Brown           valsCopy2[(offset + k) * dim + l] = valsCopy[kLocal * dim + l] * (flip ? PetscConj(flip[kLocal]) : 1.);
101c4762a1bSJed Brown         }
102c4762a1bSJed Brown       }
103c4762a1bSJed Brown       for (k = 0; k < nFunc; k++) {
1042c71b3e2SJacob Faibussowitsch         PetscCheckFalse(idsCopy2[k] != ids[k],PETSC_COMM_SELF,PETSC_ERR_PLIB,"Symmetry failure: %DD, %s, point %D, symmetry %D, order %D, functional %D: (%D != %D)",dim, tensor ? "Tensor" : "Simplex",point,j,order,k,ids[k],k);
105c4762a1bSJed Brown         for (l = 0; l < dim; l++) {
1062c71b3e2SJacob Faibussowitsch           PetscCheckFalse(valsCopy2[dim * k + l] != vals[dim * k + l],PETSC_COMM_SELF,PETSC_ERR_PLIB,"Symmetry failure: %DD, %s, point %D, symmetry %D, order %D, functional %D, component %D: (%D != %D)",dim, tensor ? "Tensor" : "Simplex",point,j,order,k,l);
107c4762a1bSJed Brown         }
108c4762a1bSJed Brown       }
109c4762a1bSJed Brown     }
110c4762a1bSJed Brown     if (anyPrinted && !printed) printed = PETSC_TRUE;
111c4762a1bSJed Brown   }
112*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexRestoreTransitiveClosure(dm,0,PETSC_TRUE,&closureSize,&closure));
113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree6(ids,idsCopy,idsCopy2,vals,valsCopy,valsCopy2));
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDualSpaceDestroy(&sp));
115*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dm));
116c4762a1bSJed Brown   PetscFunctionReturn(0);
117c4762a1bSJed Brown }
118c4762a1bSJed Brown 
119c4762a1bSJed Brown int main(int argc, char **argv)
120c4762a1bSJed Brown {
121c4762a1bSJed Brown   PetscInt       dim, order, tensor;
122c4762a1bSJed Brown   PetscErrorCode ierr;
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
125c4762a1bSJed Brown   for (tensor = 0; tensor < 2; tensor++) {
126c4762a1bSJed Brown     for (dim = 1; dim <= 3; dim++) {
127c4762a1bSJed Brown       if (dim == 1 && tensor) continue;
128c4762a1bSJed Brown       for (order = 0; order <= (tensor ? 5 : 6); order++) {
129*5f80ce2aSJacob Faibussowitsch         CHKERRQ(CheckSymmetry(dim,order,tensor ? PETSC_TRUE : PETSC_FALSE));
130c4762a1bSJed Brown       }
131c4762a1bSJed Brown     }
132c4762a1bSJed Brown   }
133c4762a1bSJed Brown   ierr = PetscFinalize();
134c4762a1bSJed Brown   return ierr;
135c4762a1bSJed Brown }
136c4762a1bSJed Brown 
137c4762a1bSJed Brown /*TEST
138c4762a1bSJed Brown   test:
139c4762a1bSJed Brown     suffix: 0
140c4762a1bSJed Brown TEST*/
141