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