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