xref: /petsc/src/dm/impls/plex/tests/ex22.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown const char help[] = "Test DMPlexCoordinatesToReference().\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscdm.h>
5*c4762a1bSJed Brown #include <petscdmplex.h>
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown static PetscErrorCode testIdentity(DM dm, PetscBool dmIsSimplicial, PetscInt cell, PetscRandom randCtx, PetscInt numPoints, PetscReal tol)
8*c4762a1bSJed Brown {
9*c4762a1bSJed Brown   PetscInt       i, j, dimC, dimR;
10*c4762a1bSJed Brown   PetscReal      *preimage, *mapped, *inverted;
11*c4762a1bSJed Brown   PetscErrorCode ierr;
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown   PetscFunctionBegin;
14*c4762a1bSJed Brown   ierr = DMGetDimension(dm,&dimR);CHKERRQ(ierr);
15*c4762a1bSJed Brown   ierr = DMGetCoordinateDim(dm,&dimC);CHKERRQ(ierr);
16*c4762a1bSJed Brown 
17*c4762a1bSJed Brown   ierr = DMGetWorkArray(dm,dimR * numPoints,MPIU_REAL,&preimage);CHKERRQ(ierr);
18*c4762a1bSJed Brown   ierr = DMGetWorkArray(dm,dimC * numPoints,MPIU_REAL,&mapped);CHKERRQ(ierr);
19*c4762a1bSJed Brown   ierr = DMGetWorkArray(dm,dimR * numPoints,MPIU_REAL,&inverted);CHKERRQ(ierr);
20*c4762a1bSJed Brown 
21*c4762a1bSJed Brown   for (i = 0; i < dimR * numPoints; i++) {
22*c4762a1bSJed Brown     ierr = PetscRandomGetValueReal(randCtx, &preimage[i]);CHKERRQ(ierr);
23*c4762a1bSJed Brown   }
24*c4762a1bSJed Brown   if (dmIsSimplicial && dimR > 1) {
25*c4762a1bSJed Brown     for (i = 0; i < numPoints; i++) {
26*c4762a1bSJed Brown       for (j = 0; j < dimR; j++) {
27*c4762a1bSJed Brown         PetscReal x = preimage[i * dimR + j];
28*c4762a1bSJed Brown         PetscReal y = preimage[i * dimR + ((j + 1) % dimR)];
29*c4762a1bSJed Brown 
30*c4762a1bSJed Brown         preimage[i * dimR + ((j + 1) % dimR)] = -1. + (y + 1.) * 0.5 * (1. - x);
31*c4762a1bSJed Brown       }
32*c4762a1bSJed Brown     }
33*c4762a1bSJed Brown   }
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown   ierr = DMPlexReferenceToCoordinates(dm,cell,numPoints,preimage,mapped);CHKERRQ(ierr);
36*c4762a1bSJed Brown   ierr = DMPlexCoordinatesToReference(dm,cell,numPoints,mapped,inverted);CHKERRQ(ierr);
37*c4762a1bSJed Brown 
38*c4762a1bSJed Brown   for (i = 0; i < numPoints; i++) {
39*c4762a1bSJed Brown     PetscReal max = 0.;
40*c4762a1bSJed Brown 
41*c4762a1bSJed Brown     for (j = 0; j < dimR; j++) {
42*c4762a1bSJed Brown       max = PetscMax(max,PetscAbsReal(preimage[i * dimR + j] - inverted[i * dimR + j]));
43*c4762a1bSJed Brown     }
44*c4762a1bSJed Brown     if (max > tol) {
45*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"Bad inversion for cell %D with error %g (tol %g): (",cell,(double)max,(double)tol);CHKERRQ(ierr);
46*c4762a1bSJed Brown       for (j = 0; j < dimR; j++) {
47*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF,"%+f",(double) preimage[i * dimR + j]);CHKERRQ(ierr);
48*c4762a1bSJed Brown         if (j < dimR - 1) {
49*c4762a1bSJed Brown           ierr = PetscPrintf(PETSC_COMM_SELF,",");CHKERRQ(ierr);
50*c4762a1bSJed Brown         }
51*c4762a1bSJed Brown       }
52*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,") --> (");CHKERRQ(ierr);
53*c4762a1bSJed Brown       for (j = 0; j < dimC; j++) {
54*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF,"%+f",(double) mapped[i * dimC + j]);CHKERRQ(ierr);
55*c4762a1bSJed Brown         if (j < dimC - 1) {
56*c4762a1bSJed Brown           ierr = PetscPrintf(PETSC_COMM_SELF,",");CHKERRQ(ierr);
57*c4762a1bSJed Brown         }
58*c4762a1bSJed Brown       }
59*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,") --> (");CHKERRQ(ierr);
60*c4762a1bSJed Brown       for (j = 0; j < dimR; j++) {
61*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF,"%+f",(double) inverted[i * dimR + j]);CHKERRQ(ierr);
62*c4762a1bSJed Brown         if (j < dimR - 1) {
63*c4762a1bSJed Brown           ierr = PetscPrintf(PETSC_COMM_SELF,",");CHKERRQ(ierr);
64*c4762a1bSJed Brown         }
65*c4762a1bSJed Brown       }
66*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,")\n");CHKERRQ(ierr);
67*c4762a1bSJed Brown     } else {
68*c4762a1bSJed Brown       char   strBuf[BUFSIZ] = {'\0'};
69*c4762a1bSJed Brown       size_t offset = 0, count;
70*c4762a1bSJed Brown 
71*c4762a1bSJed Brown       ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,"Good inversion for cell %D: (",&count,cell);CHKERRQ(ierr);
72*c4762a1bSJed Brown       offset += count;
73*c4762a1bSJed Brown       for (j = 0; j < dimR; j++) {
74*c4762a1bSJed Brown         ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,"%+f",&count,(double) preimage[i * dimR + j]);CHKERRQ(ierr);
75*c4762a1bSJed Brown         offset += count;
76*c4762a1bSJed Brown         if (j < dimR - 1) {
77*c4762a1bSJed Brown           ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,",",&count);CHKERRQ(ierr);
78*c4762a1bSJed Brown           offset += count;
79*c4762a1bSJed Brown         }
80*c4762a1bSJed Brown       }
81*c4762a1bSJed Brown       ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,") --> (",&count);CHKERRQ(ierr);
82*c4762a1bSJed Brown       offset += count;
83*c4762a1bSJed Brown       for (j = 0; j < dimC; j++) {
84*c4762a1bSJed Brown         ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,"%+f",&count,(double) mapped[i * dimC + j]);CHKERRQ(ierr);
85*c4762a1bSJed Brown         offset += count;
86*c4762a1bSJed Brown         if (j < dimC - 1) {
87*c4762a1bSJed Brown           ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,",",&count);CHKERRQ(ierr);
88*c4762a1bSJed Brown           offset += count;
89*c4762a1bSJed Brown         }
90*c4762a1bSJed Brown       }
91*c4762a1bSJed Brown       ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,") --> (",&count);CHKERRQ(ierr);
92*c4762a1bSJed Brown       offset += count;
93*c4762a1bSJed Brown       for (j = 0; j < dimR; j++) {
94*c4762a1bSJed Brown         ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,"%+f",&count,(double) inverted[i * dimR + j]);CHKERRQ(ierr);
95*c4762a1bSJed Brown         offset += count;
96*c4762a1bSJed Brown         if (j < dimR - 1) {
97*c4762a1bSJed Brown           ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,",",&count);CHKERRQ(ierr);
98*c4762a1bSJed Brown           offset += count;
99*c4762a1bSJed Brown         }
100*c4762a1bSJed Brown       }
101*c4762a1bSJed Brown       ierr = PetscSNPrintfCount(strBuf + offset,BUFSIZ-offset,")\n",&count);CHKERRQ(ierr);
102*c4762a1bSJed Brown       ierr = PetscInfo1(dm,"%s",strBuf);CHKERRQ(ierr);
103*c4762a1bSJed Brown     }
104*c4762a1bSJed Brown   }
105*c4762a1bSJed Brown 
106*c4762a1bSJed Brown   ierr = DMRestoreWorkArray(dm,dimR * numPoints,MPIU_REAL,&inverted);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = DMRestoreWorkArray(dm,dimC * numPoints,MPIU_REAL,&mapped);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = DMRestoreWorkArray(dm,dimR * numPoints,MPIU_REAL,&preimage);CHKERRQ(ierr);
109*c4762a1bSJed Brown   PetscFunctionReturn(0);
110*c4762a1bSJed Brown }
111*c4762a1bSJed Brown 
112*c4762a1bSJed Brown static PetscErrorCode identityEmbedding(PetscInt dim, PetscReal time, const PetscReal *x, PetscInt Nf, PetscScalar *u, void *ctx)
113*c4762a1bSJed Brown {
114*c4762a1bSJed Brown   PetscInt i;
115*c4762a1bSJed Brown 
116*c4762a1bSJed Brown   for (i = 0; i < dim; i++) {u[i] = x[i];};
117*c4762a1bSJed Brown   return 0;
118*c4762a1bSJed Brown }
119*c4762a1bSJed Brown 
120*c4762a1bSJed Brown int main(int argc, char **argv)
121*c4762a1bSJed Brown {
122*c4762a1bSJed Brown   PetscRandom    randCtx;
123*c4762a1bSJed Brown   PetscInt       dim, dimC, isSimplex, isFE, numTests = 10;
124*c4762a1bSJed Brown   PetscReal      perturb = 0.1, tol = 10. * PETSC_SMALL;
125*c4762a1bSJed Brown   PetscErrorCode ierr;
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
128*c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_SELF,&randCtx);CHKERRQ(ierr);
129*c4762a1bSJed Brown   ierr = PetscRandomSetInterval(randCtx,-1.,1.);CHKERRQ(ierr);
130*c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(randCtx);CHKERRQ(ierr);
131*c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ex21",NULL);CHKERRQ(ierr);
132*c4762a1bSJed Brown   ierr = PetscOptionsReal("-vertex_perturbation","scale of random vertex distortion",NULL,perturb,&perturb,NULL);CHKERRQ(ierr);
133*c4762a1bSJed Brown   ierr = PetscOptionsBoundedInt("-num_test_points","number of points to test",NULL,numTests,&numTests,NULL,0);CHKERRQ(ierr);
134*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
135*c4762a1bSJed Brown   for (dim = 1; dim <= 3; dim++) {
136*c4762a1bSJed Brown     for (dimC = dim; dimC <= PetscMin(3,dim + 1); dimC++) {
137*c4762a1bSJed Brown       for (isSimplex = 0; isSimplex < 2; isSimplex++) {
138*c4762a1bSJed Brown         for (isFE = 0; isFE < 2; isFE++) {
139*c4762a1bSJed Brown           DM           dm;
140*c4762a1bSJed Brown           Vec          coords;
141*c4762a1bSJed Brown           PetscScalar *coordArray;
142*c4762a1bSJed Brown           PetscReal    noise;
143*c4762a1bSJed Brown           PetscInt     i, n, order = 1;
144*c4762a1bSJed Brown 
145*c4762a1bSJed Brown           ierr = DMPlexCreateReferenceCell(PETSC_COMM_SELF,dim,isSimplex ? PETSC_TRUE : PETSC_FALSE, &dm);CHKERRQ(ierr);
146*c4762a1bSJed Brown           if (isFE) {
147*c4762a1bSJed Brown             DM             dmCoord;
148*c4762a1bSJed Brown             PetscSpace     sp;
149*c4762a1bSJed Brown             PetscFE        fe;
150*c4762a1bSJed Brown             Vec            localCoords;
151*c4762a1bSJed Brown             PetscErrorCode (*funcs[1])(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *) = {identityEmbedding};
152*c4762a1bSJed Brown             void           *ctxs[1] = {NULL};
153*c4762a1bSJed Brown 
154*c4762a1bSJed Brown             ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm),dim,dim,isSimplex ? PETSC_TRUE : PETSC_FALSE,isSimplex ? NULL : "tensor_" ,PETSC_DEFAULT,&fe);CHKERRQ(ierr);
155*c4762a1bSJed Brown             ierr = PetscFEGetBasisSpace(fe,&sp);CHKERRQ(ierr);
156*c4762a1bSJed Brown             ierr = PetscSpaceGetDegree(sp,&order,NULL);CHKERRQ(ierr);
157*c4762a1bSJed Brown             ierr = DMSetField(dm,0,NULL,(PetscObject)fe);CHKERRQ(ierr);
158*c4762a1bSJed Brown             ierr = DMCreateDS(dm);CHKERRQ(ierr);
159*c4762a1bSJed Brown             ierr = DMCreateLocalVector(dm,&localCoords);CHKERRQ(ierr);
160*c4762a1bSJed Brown             ierr = DMProjectFunctionLocal(dm,0,funcs,ctxs,INSERT_VALUES,localCoords);CHKERRQ(ierr);
161*c4762a1bSJed Brown             ierr = VecSetDM(localCoords,NULL);CHKERRQ(ierr); /* This is necessary to prevent a reference loop */
162*c4762a1bSJed Brown             ierr = DMClone(dm,&dmCoord);CHKERRQ(ierr);
163*c4762a1bSJed Brown             ierr = DMSetField(dmCoord,0,NULL,(PetscObject)fe);CHKERRQ(ierr);
164*c4762a1bSJed Brown             ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
165*c4762a1bSJed Brown             ierr = DMCreateDS(dmCoord);CHKERRQ(ierr);
166*c4762a1bSJed Brown             ierr = DMSetCoordinateDM(dm,dmCoord);CHKERRQ(ierr);
167*c4762a1bSJed Brown             ierr = DMDestroy(&dmCoord);CHKERRQ(ierr);
168*c4762a1bSJed Brown             ierr = DMSetCoordinatesLocal(dm,localCoords);CHKERRQ(ierr);
169*c4762a1bSJed Brown             ierr = VecDestroy(&localCoords);CHKERRQ(ierr);
170*c4762a1bSJed Brown           }
171*c4762a1bSJed Brown           ierr = PetscInfo4(dm,"Testing %s%D %DD mesh embedded in %DD\n",isSimplex ? "P" : "Q" , order, dim, dimC);CHKERRQ(ierr);
172*c4762a1bSJed Brown           ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
173*c4762a1bSJed Brown           ierr = VecGetLocalSize(coords,&n);CHKERRQ(ierr);
174*c4762a1bSJed Brown           if (dimC > dim) { /* reembed in higher dimension */
175*c4762a1bSJed Brown             PetscSection sec, newSec;
176*c4762a1bSJed Brown             PetscInt     pStart, pEnd, p, i, newN;
177*c4762a1bSJed Brown             Vec          newVec;
178*c4762a1bSJed Brown             DM           coordDM;
179*c4762a1bSJed Brown             PetscScalar  *newCoordArray;
180*c4762a1bSJed Brown 
181*c4762a1bSJed Brown             ierr = DMGetCoordinateSection(dm,&sec);CHKERRQ(ierr);
182*c4762a1bSJed Brown             ierr = PetscSectionCreate(PetscObjectComm((PetscObject)sec),&newSec);CHKERRQ(ierr);
183*c4762a1bSJed Brown             ierr = PetscSectionSetNumFields(newSec,1);CHKERRQ(ierr);
184*c4762a1bSJed Brown             ierr = PetscSectionGetChart(sec,&pStart,&pEnd);CHKERRQ(ierr);
185*c4762a1bSJed Brown             ierr = PetscSectionSetChart(newSec,pStart,pEnd);CHKERRQ(ierr);
186*c4762a1bSJed Brown             for (p = pStart; p < pEnd; p++) {
187*c4762a1bSJed Brown               PetscInt nDof;
188*c4762a1bSJed Brown 
189*c4762a1bSJed Brown               ierr = PetscSectionGetDof(sec,p,&nDof);CHKERRQ(ierr);
190*c4762a1bSJed Brown               if (nDof % dim) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Coordinate section point %D has %D dofs != 0 mod %D",p,nDof,dim);
191*c4762a1bSJed Brown               ierr = PetscSectionSetDof(newSec,p,(nDof/dim)*dimC);CHKERRQ(ierr);
192*c4762a1bSJed Brown               ierr = PetscSectionSetFieldDof(newSec,p,0,(nDof/dim)*dimC);CHKERRQ(ierr);
193*c4762a1bSJed Brown             }
194*c4762a1bSJed Brown             ierr = PetscSectionSetUp(newSec);CHKERRQ(ierr);
195*c4762a1bSJed Brown             ierr = PetscSectionGetStorageSize(newSec,&newN);CHKERRQ(ierr);
196*c4762a1bSJed Brown             ierr = VecCreateSeq(PETSC_COMM_SELF,newN,&newVec);CHKERRQ(ierr);
197*c4762a1bSJed Brown             ierr = VecGetArray(newVec,&newCoordArray);CHKERRQ(ierr);
198*c4762a1bSJed Brown             ierr = VecGetArray(coords,&coordArray);CHKERRQ(ierr);
199*c4762a1bSJed Brown             for (i = 0; i < n / dim; i++) {
200*c4762a1bSJed Brown               PetscInt j;
201*c4762a1bSJed Brown 
202*c4762a1bSJed Brown               for (j = 0; j < dim; j++) {
203*c4762a1bSJed Brown                 newCoordArray[i * dimC + j] = coordArray[i * dim + j];
204*c4762a1bSJed Brown               }
205*c4762a1bSJed Brown               for (; j < dimC; j++) {
206*c4762a1bSJed Brown                 newCoordArray[i * dimC + j] = 0.;
207*c4762a1bSJed Brown               }
208*c4762a1bSJed Brown             }
209*c4762a1bSJed Brown             ierr = VecRestoreArray(coords,&coordArray);CHKERRQ(ierr);
210*c4762a1bSJed Brown             ierr = VecRestoreArray(newVec,&newCoordArray);CHKERRQ(ierr);
211*c4762a1bSJed Brown             ierr = DMSetCoordinateDim(dm,dimC);CHKERRQ(ierr);
212*c4762a1bSJed Brown             ierr = DMSetCoordinateSection(dm,dimC,newSec);CHKERRQ(ierr);
213*c4762a1bSJed Brown             ierr = DMSetCoordinatesLocal(dm,newVec);CHKERRQ(ierr);
214*c4762a1bSJed Brown             ierr = VecDestroy(&newVec);CHKERRQ(ierr);
215*c4762a1bSJed Brown             ierr = PetscSectionDestroy(&newSec);CHKERRQ(ierr);
216*c4762a1bSJed Brown             ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
217*c4762a1bSJed Brown             ierr = DMGetCoordinateDM(dm,&coordDM);CHKERRQ(ierr);
218*c4762a1bSJed Brown             if (isFE) {
219*c4762a1bSJed Brown               PetscFE fe;
220*c4762a1bSJed Brown 
221*c4762a1bSJed Brown               ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm),dim,dimC,isSimplex ? PETSC_TRUE : PETSC_FALSE,isSimplex ? NULL : "tensor_",PETSC_DEFAULT,&fe);CHKERRQ(ierr);
222*c4762a1bSJed Brown               ierr = DMSetField(coordDM,0,NULL,(PetscObject)fe);CHKERRQ(ierr);
223*c4762a1bSJed Brown               ierr = PetscFEDestroy(&fe);CHKERRQ(ierr);
224*c4762a1bSJed Brown               ierr = DMCreateDS(coordDM);CHKERRQ(ierr);
225*c4762a1bSJed Brown             }
226*c4762a1bSJed Brown             ierr = DMSetCoordinateDim(coordDM,dimC);CHKERRQ(ierr);
227*c4762a1bSJed Brown             ierr = VecGetLocalSize(coords,&n);CHKERRQ(ierr);
228*c4762a1bSJed Brown           }
229*c4762a1bSJed Brown           ierr = VecGetArray(coords,&coordArray);CHKERRQ(ierr);
230*c4762a1bSJed Brown           for (i = 0; i < n; i++) {
231*c4762a1bSJed Brown             ierr = PetscRandomGetValueReal(randCtx,&noise);CHKERRQ(ierr);
232*c4762a1bSJed Brown             coordArray[i] += noise * perturb;
233*c4762a1bSJed Brown           }
234*c4762a1bSJed Brown           ierr = VecRestoreArray(coords,&coordArray);CHKERRQ(ierr);
235*c4762a1bSJed Brown           ierr = DMSetCoordinatesLocal(dm, coords);CHKERRQ(ierr);
236*c4762a1bSJed Brown           ierr = testIdentity(dm, isSimplex ? PETSC_TRUE : PETSC_FALSE, 0, randCtx, numTests, tol);CHKERRQ(ierr);
237*c4762a1bSJed Brown           ierr = DMDestroy(&dm);CHKERRQ(ierr);
238*c4762a1bSJed Brown         }
239*c4762a1bSJed Brown       }
240*c4762a1bSJed Brown     }
241*c4762a1bSJed Brown   }
242*c4762a1bSJed Brown   ierr = PetscRandomDestroy(&randCtx);CHKERRQ(ierr);
243*c4762a1bSJed Brown   ierr = PetscFinalize();
244*c4762a1bSJed Brown   return ierr;
245*c4762a1bSJed Brown }
246*c4762a1bSJed Brown 
247*c4762a1bSJed Brown /*TEST
248*c4762a1bSJed Brown 
249*c4762a1bSJed Brown   test:
250*c4762a1bSJed Brown     suffix: 0
251*c4762a1bSJed Brown     args: -petscspace_degree 2 -tensor_petscspace_degree 2
252*c4762a1bSJed Brown 
253*c4762a1bSJed Brown TEST*/
254