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