xref: /petsc/src/dm/impls/plex/tests/ex75.c (revision 62ed42821f84f33353546905c9271ed730d39ff0)
1605a06ccSksagiyam static char help[] = "Tests for submesh with both CG and DG coordinates\n\n";
2605a06ccSksagiyam 
3605a06ccSksagiyam #include <petscdmplex.h>
4605a06ccSksagiyam #include <petsc/private/dmimpl.h>
5605a06ccSksagiyam 
main(int argc,char ** argv)6605a06ccSksagiyam int main(int argc, char **argv)
7605a06ccSksagiyam {
8605a06ccSksagiyam   PetscInt       dim = 1, d, cStart, cEnd, c, q, degree = 2, coordSize, offset;
9605a06ccSksagiyam   PetscReal      R = 1.0;
10605a06ccSksagiyam   DM             dm, coordDM, subdm;
11605a06ccSksagiyam   PetscSection   coordSec;
12605a06ccSksagiyam   Vec            coordVec;
13605a06ccSksagiyam   PetscScalar   *coords;
14605a06ccSksagiyam   DMLabel        filter;
15605a06ccSksagiyam   const PetscInt filterValue = 1;
16605a06ccSksagiyam   MPI_Comm       comm;
17605a06ccSksagiyam   PetscMPIInt    size;
18605a06ccSksagiyam 
19605a06ccSksagiyam   PetscFunctionBeginUser;
20605a06ccSksagiyam   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
21605a06ccSksagiyam   comm = PETSC_COMM_WORLD;
22605a06ccSksagiyam   PetscCallMPI(MPI_Comm_size(comm, &size));
23605a06ccSksagiyam   if (size != 1) {
24605a06ccSksagiyam     PetscCall(PetscPrintf(comm, "This example is specifically designed for comm size == 1.\n"));
25605a06ccSksagiyam     PetscCall(PetscFinalize());
26605a06ccSksagiyam     return 0;
27605a06ccSksagiyam   }
28605a06ccSksagiyam   /* Make a unit circle with 16 elements. */
29605a06ccSksagiyam   PetscCall(DMPlexCreateSphereMesh(comm, dim, PETSC_TRUE, R, &dm));
30605a06ccSksagiyam   PetscUseTypeMethod(dm, createcellcoordinatedm, &coordDM);
31605a06ccSksagiyam   PetscCall(DMSetCellCoordinateDM(dm, coordDM));
32605a06ccSksagiyam   PetscCall(DMDestroy(&coordDM));
33605a06ccSksagiyam   PetscCall(DMGetCellCoordinateSection(dm, &coordSec));
34605a06ccSksagiyam   PetscCall(PetscSectionSetNumFields(coordSec, 1));
35605a06ccSksagiyam   PetscCall(PetscSectionSetFieldComponents(coordSec, 0, dim));
36605a06ccSksagiyam   PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd));
37605a06ccSksagiyam   PetscCall(PetscSectionSetChart(coordSec, cStart, cEnd));
38605a06ccSksagiyam   for (c = cStart; c < cEnd; ++c) {
39605a06ccSksagiyam     PetscCall(PetscSectionSetDof(coordSec, c, dim * (degree + 1)));
40605a06ccSksagiyam     PetscCall(PetscSectionSetFieldDof(coordSec, c, 0, dim * (degree + 1)));
41605a06ccSksagiyam   }
42605a06ccSksagiyam   PetscCall(PetscSectionSetUp(coordSec));
43605a06ccSksagiyam   PetscCall(PetscSectionGetStorageSize(coordSec, &coordSize));
44605a06ccSksagiyam   PetscCall(VecCreate(PETSC_COMM_SELF, &coordVec));
45605a06ccSksagiyam   PetscCall(PetscObjectSetName((PetscObject)coordVec, "cellcoordinates"));
46605a06ccSksagiyam   PetscCall(VecSetBlockSize(coordVec, dim));
47605a06ccSksagiyam   PetscCall(VecSetSizes(coordVec, coordSize, PETSC_DETERMINE));
48605a06ccSksagiyam   PetscCall(VecSetType(coordVec, VECSTANDARD));
49605a06ccSksagiyam   PetscCall(VecGetArray(coordVec, &coords));
50605a06ccSksagiyam   for (c = cStart; c < cEnd; ++c) {
51605a06ccSksagiyam     PetscCall(PetscSectionGetOffset(coordSec, c, &offset));
52605a06ccSksagiyam     for (q = 0; q < degree + 1; ++q) {
53605a06ccSksagiyam       /* Make some DG coordinates. Note that dim = 1.*/
54605a06ccSksagiyam       for (d = 0; d < dim; ++d) coords[offset + dim * q + d] = 100. + (PetscScalar)c + (1.0 / (PetscScalar)degree) * (PetscScalar)q;
55605a06ccSksagiyam     }
56605a06ccSksagiyam   }
57605a06ccSksagiyam   PetscCall(VecRestoreArray(coordVec, &coords));
58605a06ccSksagiyam   PetscCall(DMSetCellCoordinatesLocal(dm, coordVec));
59605a06ccSksagiyam   PetscCall(VecDestroy(&coordVec));
60605a06ccSksagiyam   /* Make submesh. */
61605a06ccSksagiyam   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "filter", &filter));
62605a06ccSksagiyam   PetscCall(DMLabelSetValue(filter, 15, filterValue)); /* last cell */
63*71f1c950SStefano Zampini   PetscCall(DMPlexFilter(dm, filter, filterValue, PETSC_FALSE, PETSC_FALSE, PetscObjectComm((PetscObject)dm), NULL, &subdm));
64605a06ccSksagiyam   PetscCall(DMLabelDestroy(&filter));
65605a06ccSksagiyam   PetscCall(PetscObjectSetName((PetscObject)subdm, "Example_SubDM"));
66605a06ccSksagiyam   PetscCall(DMViewFromOptions(subdm, NULL, "-dm_view"));
67605a06ccSksagiyam   PetscCall(DMDestroy(&subdm));
68605a06ccSksagiyam   PetscCall(DMDestroy(&dm));
69605a06ccSksagiyam   PetscCall(PetscFinalize());
70605a06ccSksagiyam   return 0;
71605a06ccSksagiyam }
72605a06ccSksagiyam 
73605a06ccSksagiyam /*TEST
74605a06ccSksagiyam 
75605a06ccSksagiyam   test:
76605a06ccSksagiyam     suffix: 0
77605a06ccSksagiyam     args: -dm_view ascii::ascii_info_detail
78605a06ccSksagiyam 
79605a06ccSksagiyam TEST*/
80