xref: /petsc/src/dm/impls/plex/plexadapt.c (revision 2ee120b006de7f0acd33a3aec6d625988d480ddd)
10bbe5d31Sbarral #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
20bbe5d31Sbarral #ifdef PETSC_HAVE_PRAGMATIC
30bbe5d31Sbarral #include <pragmatic/cpragmatic.h>
40bbe5d31Sbarral #endif
50bbe5d31Sbarral 
6f7bcbcbbSMatthew G. Knepley #undef __FUNCT__
7f7bcbcbbSMatthew G. Knepley #define __FUNCT__ "DMPlexRemesh_Internal"
8f7bcbcbbSMatthew G. Knepley /*
9f7bcbcbbSMatthew G. Knepley   DMPlexRemesh_Internal - Generates a new mesh conforming to a metric field.
100bbe5d31Sbarral 
11f7bcbcbbSMatthew G. Knepley   Input Parameters:
12f7bcbcbbSMatthew G. Knepley + dm - The DM object
13f7bcbcbbSMatthew G. Knepley . vertexMetric - The metric to which the mesh is adapted, defined vertex-wise.
14f7bcbcbbSMatthew G. Knepley - bdLabelName - Label name for boundary tags which are preserved in dmNew, or NULL. Should not be "_boundary_".
15f7bcbcbbSMatthew G. Knepley 
16f7bcbcbbSMatthew G. Knepley   Output Parameter:
17f7bcbcbbSMatthew G. Knepley . dmNew  - the new DM
18f7bcbcbbSMatthew G. Knepley 
19f7bcbcbbSMatthew G. Knepley   Level: advanced
20f7bcbcbbSMatthew G. Knepley 
21f7bcbcbbSMatthew G. Knepley .seealso: DMCoarsen(), DMRefine()
22f7bcbcbbSMatthew G. Knepley */
23f7bcbcbbSMatthew G. Knepley PetscErrorCode DMPlexRemesh_Internal(DM dm, Vec vertexMetric, const char bdLabelName[], DM *dmNew)
24f7bcbcbbSMatthew G. Knepley {
25f7bcbcbbSMatthew G. Knepley   const char        *bdName = "_boundary_";
26f7bcbcbbSMatthew G. Knepley   DM                 udm, cdm;
27*2ee120b0SMatthew G. Knepley   DMLabel            bdLabel = NULL, bdLabelFull;
28f7bcbcbbSMatthew G. Knepley   IS                 bdIS;
29f7bcbcbbSMatthew G. Knepley   PetscSection       coordSection;
30f7bcbcbbSMatthew G. Knepley   Vec                coordinates;
31f7bcbcbbSMatthew G. Knepley   const PetscScalar *coords, *met;
32f7bcbcbbSMatthew G. Knepley   const PetscInt    *bdFacesFull;
33*2ee120b0SMatthew G. Knepley   PetscInt          *bdFaces, *bdFaceIds;
34f7bcbcbbSMatthew G. Knepley   PetscReal         *x, *y, *z, *metric;
35*2ee120b0SMatthew G. Knepley   PetscInt          *cells;
36*2ee120b0SMatthew G. Knepley   PetscInt           dim, cStart, cEnd, numCells, c, coff, vStart, vEnd, numVertices, v;
37*2ee120b0SMatthew G. Knepley   PetscInt           off, maxConeSize, numBdFaces, f, bdSize;
38f7bcbcbbSMatthew G. Knepley   PetscBool          flg;
39*2ee120b0SMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC
40*2ee120b0SMatthew G. Knepley   DMLabel            bdLabelNew;
41*2ee120b0SMatthew G. Knepley   PetscScalar       *coordsNew;
42*2ee120b0SMatthew G. Knepley   PetscInt          *bdTags;
43*2ee120b0SMatthew G. Knepley   PetscReal         *xNew[3] = {NULL, NULL, NULL};
44*2ee120b0SMatthew G. Knepley   PetscInt          *cellsNew;
45*2ee120b0SMatthew G. Knepley   PetscInt           d, numCellsNew, numVerticesNew;
46*2ee120b0SMatthew G. Knepley   PetscInt           numCornersNew, fStart, fEnd;
47*2ee120b0SMatthew G. Knepley #endif
48f7bcbcbbSMatthew G. Knepley   PetscErrorCode     ierr;
49f7bcbcbbSMatthew G. Knepley 
50f7bcbcbbSMatthew G. Knepley   PetscFunctionBegin;
51f7bcbcbbSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
52f7bcbcbbSMatthew G. Knepley   PetscValidHeaderSpecific(vertexMetric, VEC_CLASSID, 2);
53f7bcbcbbSMatthew G. Knepley   PetscValidCharPointer(bdLabelName, 3);
54f7bcbcbbSMatthew G. Knepley   PetscValidPointer(dmNew, 4);
55f7bcbcbbSMatthew G. Knepley   if (bdLabelName) {
56f7bcbcbbSMatthew G. Knepley     size_t len;
57f7bcbcbbSMatthew G. Knepley 
58f7bcbcbbSMatthew G. Knepley     ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr);
59f7bcbcbbSMatthew G. Knepley     if (flg) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName);
60f7bcbcbbSMatthew G. Knepley     ierr = PetscStrlen(bdLabelName, &len);CHKERRQ(ierr);
61f7bcbcbbSMatthew G. Knepley     if (len) {
62f7bcbcbbSMatthew G. Knepley       ierr = DMGetLabel(dm, bdLabelName, &bdLabel);CHKERRQ(ierr);
63f7bcbcbbSMatthew G. Knepley       if (!bdLabel) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Label \"%s\" does not exist in DM", bdLabelName);
64f7bcbcbbSMatthew G. Knepley     }
65f7bcbcbbSMatthew G. Knepley   }
66f7bcbcbbSMatthew G. Knepley   /* Get mesh information */
67f7bcbcbbSMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
68f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
69f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
70f7bcbcbbSMatthew G. Knepley   ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
71f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr);
72f7bcbcbbSMatthew G. Knepley   numCells    = cEnd - cStart;
73f7bcbcbbSMatthew G. Knepley   numVertices = vEnd - vStart;
74f7bcbcbbSMatthew G. Knepley   ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr);
75f7bcbcbbSMatthew G. Knepley   for (c = 0, coff = 0; c < numCells; ++c) {
76f7bcbcbbSMatthew G. Knepley     const PetscInt *cone;
77f7bcbcbbSMatthew G. Knepley     PetscInt        coneSize, cl;
78f7bcbcbbSMatthew G. Knepley 
79f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
80f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
81f7bcbcbbSMatthew G. Knepley     for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart;
82f7bcbcbbSMatthew G. Knepley   }
83f7bcbcbbSMatthew G. Knepley   ierr = DMDestroy(&udm);CHKERRQ(ierr);
84f7bcbcbbSMatthew G. Knepley   ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr);
85f7bcbcbbSMatthew G. Knepley   ierr = DMGetDefaultSection(cdm, &coordSection);CHKERRQ(ierr);
86f7bcbcbbSMatthew G. Knepley   ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
87f7bcbcbbSMatthew G. Knepley   ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
88f7bcbcbbSMatthew G. Knepley   for (v = vStart; v < vEnd; ++v) {
89f7bcbcbbSMatthew G. Knepley     ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
905f2e139eSMatthew G. Knepley     x[v-vStart] = PetscRealPart(coords[off+0]);
915f2e139eSMatthew G. Knepley     if (dim > 1) y[v-vStart] = PetscRealPart(coords[off+1]);
925f2e139eSMatthew G. Knepley     if (dim > 2) z[v-vStart] = PetscRealPart(coords[off+2]);
93f7bcbcbbSMatthew G. Knepley   }
94f7bcbcbbSMatthew G. Knepley   ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
95f7bcbcbbSMatthew G. Knepley   /* Get boundary mesh */
96f7bcbcbbSMatthew G. Knepley   ierr = DMLabelCreate(bdName, &bdLabelFull);CHKERRQ(ierr);
97f7bcbcbbSMatthew G. Knepley   ierr = DMPlexMarkBoundaryFaces(dm, bdLabelFull);CHKERRQ(ierr);
98f7bcbcbbSMatthew G. Knepley   ierr = DMLabelGetStratumIS(bdLabelFull, 1, &bdIS);CHKERRQ(ierr);
99f7bcbcbbSMatthew G. Knepley   ierr = DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces);CHKERRQ(ierr);
100f7bcbcbbSMatthew G. Knepley   ierr = ISGetIndices(bdIS, &bdFacesFull);CHKERRQ(ierr);
101f7bcbcbbSMatthew G. Knepley   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
102f7bcbcbbSMatthew G. Knepley     PetscInt *closure = NULL;
103f7bcbcbbSMatthew G. Knepley     PetscInt  closureSize, cl;
104f7bcbcbbSMatthew G. Knepley 
105f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
106f7bcbcbbSMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
107f7bcbcbbSMatthew G. Knepley       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
108f7bcbcbbSMatthew G. Knepley     }
109f7bcbcbbSMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
110f7bcbcbbSMatthew G. Knepley   }
111f7bcbcbbSMatthew G. Knepley   ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr);
112f7bcbcbbSMatthew G. Knepley   for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
113f7bcbcbbSMatthew G. Knepley     PetscInt *closure = NULL;
114f7bcbcbbSMatthew G. Knepley     PetscInt  closureSize, cl;
115f7bcbcbbSMatthew G. Knepley 
116f7bcbcbbSMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
117f7bcbcbbSMatthew G. Knepley     for (cl = 0; cl < closureSize*2; cl += 2) {
118f7bcbcbbSMatthew G. Knepley       if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart;
119f7bcbcbbSMatthew G. Knepley     }
120f7bcbcbbSMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
121f7bcbcbbSMatthew G. Knepley     if (bdLabel) {ierr = DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]);CHKERRQ(ierr);}
122f7bcbcbbSMatthew G. Knepley     else         {bdFaceIds[f] = 1;}
123f7bcbcbbSMatthew G. Knepley   }
124f7bcbcbbSMatthew G. Knepley   ierr = ISDestroy(&bdIS);CHKERRQ(ierr);
125f7bcbcbbSMatthew G. Knepley   ierr = DMLabelDestroy(&bdLabelFull);CHKERRQ(ierr);
126f7bcbcbbSMatthew G. Knepley   /* Get metric */
127f7bcbcbbSMatthew G. Knepley   ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr);
1285f2e139eSMatthew G. Knepley   for (v = 0; v < (vEnd-vStart)*PetscSqr(dim); ++v) metric[v] = PetscRealPart(met[v]);
129f7bcbcbbSMatthew G. Knepley   ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr);
130f7bcbcbbSMatthew G. Knepley   /* Create new mesh */
131f7bcbcbbSMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC
132f7bcbcbbSMatthew G. Knepley   switch (dim) {
133f7bcbcbbSMatthew G. Knepley   case 2:
134f7bcbcbbSMatthew G. Knepley     pragmatic_2d_init(&numVertices, &numCells, cells, x, y);break;
135f7bcbcbbSMatthew G. Knepley   case 3:
136f7bcbcbbSMatthew G. Knepley     pragmatic_3d_init(&numVertices, &numCells, cells, x, y, z);break;
137f7bcbcbbSMatthew G. Knepley   default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
138f7bcbcbbSMatthew G. Knepley   }
139f7bcbcbbSMatthew G. Knepley   pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds);
140f7bcbcbbSMatthew G. Knepley   pragmatic_set_metric(metric);
141f7bcbcbbSMatthew G. Knepley   pragmatic_adapt();
142f7bcbcbbSMatthew G. Knepley   /* Read out mesh */
143f7bcbcbbSMatthew G. Knepley   pragmatic_get_info(&numVerticesNew, &numCellsNew);
144f7bcbcbbSMatthew G. Knepley   ierr = PetscMalloc1(numVerticesNew*dim, &coordsNew);CHKERRQ(ierr);
145f7bcbcbbSMatthew G. Knepley   switch (dim) {
146f7bcbcbbSMatthew G. Knepley   case 2:
147f7bcbcbbSMatthew G. Knepley     numCornersNew = 3;
148f7bcbcbbSMatthew G. Knepley     ierr = PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]);CHKERRQ(ierr);
149f7bcbcbbSMatthew G. Knepley     pragmatic_get_coords_2d(xNew[0], xNew[1]);
150f7bcbcbbSMatthew G. Knepley     break;
151f7bcbcbbSMatthew G. Knepley   case 3:
152f7bcbcbbSMatthew G. Knepley     numCornersNew = 4;
153f7bcbcbbSMatthew G. Knepley     ierr = PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]);CHKERRQ(ierr);
154f7bcbcbbSMatthew G. Knepley     pragmatic_get_coords_3d(xNew[0], xNew[1], xNew[2]);
155f7bcbcbbSMatthew G. Knepley     break;
156f7bcbcbbSMatthew G. Knepley   default:
157f7bcbcbbSMatthew G. Knepley     SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
158f7bcbcbbSMatthew G. Knepley   }
159f7bcbcbbSMatthew G. Knepley   for (v = 0; v < numVerticesNew; ++v) {for (d = 0; d < dim; ++d) coordsNew[v*dim+d] = xNew[d][v];}
160f7bcbcbbSMatthew G. Knepley   ierr = PetscMalloc1(numCellsNew*(dim+1), &cellsNew);CHKERRQ(ierr);
161f7bcbcbbSMatthew G. Knepley   pragmatic_get_elements(cellsNew);
162f7bcbcbbSMatthew G. Knepley   ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject) dm), dim, numCellsNew, numVerticesNew, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, dmNew);CHKERRQ(ierr);
163f7bcbcbbSMatthew G. Knepley   /* Read out boundary label */
164f7bcbcbbSMatthew G. Knepley   pragmatic_get_boundaryTags(&bdTags);
165f7bcbcbbSMatthew G. Knepley   ierr = DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName);CHKERRQ(ierr);
166f7bcbcbbSMatthew G. Knepley   ierr = DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew);CHKERRQ(ierr);
167f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr);
168f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr);
169f7bcbcbbSMatthew G. Knepley   ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr);
170f7bcbcbbSMatthew G. Knepley   for (c = cStart; c < cEnd; ++c) {
171f7bcbcbbSMatthew G. Knepley     /* Only for simplicial meshes */
172f7bcbcbbSMatthew G. Knepley     coff = (c-cStart)*(dim+1);
173f7bcbcbbSMatthew G. Knepley     for (d = 0; d < dim+1; ++d) {
174f7bcbcbbSMatthew G. Knepley       if (bdTags[coff+d]) {
175f7bcbcbbSMatthew G. Knepley         const PetscInt *faces;
176f7bcbcbbSMatthew G. Knepley         PetscInt        numFaces, vertices[3], p;
177f7bcbcbbSMatthew G. Knepley 
178f7bcbcbbSMatthew G. Knepley         /* Mark face opposite to this vertex */
179f7bcbcbbSMatthew G. Knepley         for (p = 0, v = 0; p < dim+1; ++p) if (p != d) {vertices[v] = cellsNew[coff+p] + vStart; ++v;}
180f7bcbcbbSMatthew G. Knepley         ierr = DMPlexGetFullJoin(*dmNew, dim, vertices, &numFaces, &faces);CHKERRQ(ierr);
181f7bcbcbbSMatthew G. Knepley         if (numFaces != 1) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Cannot find boundary face in new cell %D for vertex %D(%D) with tag %D", c, cellsNew[coff+d], d, bdTags[coff+d]);
182f7bcbcbbSMatthew G. Knepley         if (faces[0] < fStart || faces[0] >= fEnd) SETERRQ5(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Boundary point %D in new cell %D for vertex %D(%D) with tag %D is not a face", faces[0], c, cellsNew[coff+d], d, bdTags[coff+d]);
183f7bcbcbbSMatthew G. Knepley         ierr = DMLabelSetValue(bdLabelNew, faces[0], bdTags[coff+d]);CHKERRQ(ierr);
184f7bcbcbbSMatthew G. Knepley         ierr = DMPlexRestoreJoin(*dmNew, dim, vertices, &numFaces, &faces);CHKERRQ(ierr);
185f7bcbcbbSMatthew G. Knepley       }
186f7bcbcbbSMatthew G. Knepley     }
187f7bcbcbbSMatthew G. Knepley   }
188f7bcbcbbSMatthew G. Knepley   /* Cleanup */
189f7bcbcbbSMatthew G. Knepley   switch (dim) {
190f7bcbcbbSMatthew G. Knepley   case 2: ierr = PetscFree2(xNew[0], xNew[1]);CHKERRQ(ierr);break;
191f7bcbcbbSMatthew G. Knepley   case 3: ierr = PetscFree3(xNew[0], xNew[1], xNew[2]);CHKERRQ(ierr);break;
192f7bcbcbbSMatthew G. Knepley   }
193f7bcbcbbSMatthew G. Knepley   ierr = PetscFree(cellsNew);CHKERRQ(ierr);
194f7bcbcbbSMatthew G. Knepley   ierr = PetscFree5(x, y, z, metric, cells);CHKERRQ(ierr);
195f7bcbcbbSMatthew G. Knepley   ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr);
196f7bcbcbbSMatthew G. Knepley   ierr = PetscFree(coordsNew);CHKERRQ(ierr);
197f7bcbcbbSMatthew G. Knepley   pragmatic_finalize();
198f7bcbcbbSMatthew G. Knepley #else
199f7bcbcbbSMatthew G. Knepley   SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Remeshing needs external package support.\nPlease reconfigure with --download-pragmatic.");
200f7bcbcbbSMatthew G. Knepley #endif
201f7bcbcbbSMatthew G. Knepley   PetscFunctionReturn(0);
202f7bcbcbbSMatthew G. Knepley }
2030bbe5d31Sbarral 
2040bbe5d31Sbarral #undef __FUNCT__
20552b40ca0Sbarral #define __FUNCT__ "DMPlexAdapt"
2060f7e110dSbarral /*@
2070f7e110dSbarral   DMPlexAdapt - Generates a mesh adapted to the specified metric field using the pragmatic library.
2080f7e110dSbarral 
2090f7e110dSbarral   Input Parameters:
2100f7e110dSbarral + dm - The DM object
2110f7e110dSbarral . metric - The metric to which the mesh is adapted, defined vertex-wise.
212379fadc5Sbarral - bdyLabelName - Label name for boundary tags. These will be preserved in the output mesh. bdyLabelName should be "" (empty string) if there is no such label, and should be different from "boundary".
2130f7e110dSbarral 
2140f7e110dSbarral   Output Parameter:
215f7bcbcbbSMatthew G. Knepley . dmAdapt  - Pointer to the DM object containing the adapted mesh
2160f7e110dSbarral 
2170f7e110dSbarral   Level: advanced
2180f7e110dSbarral 
2190f7e110dSbarral .seealso: DMCoarsen(), DMRefine()
2200f7e110dSbarral @*/
221f7bcbcbbSMatthew G. Knepley PetscErrorCode DMPlexAdapt(DM dm, Vec metric, const char bdLabelName[], DM *dmAdapt)
2220bbe5d31Sbarral {
2230bbe5d31Sbarral   PetscErrorCode ierr;
2240bbe5d31Sbarral 
2250bbe5d31Sbarral   PetscFunctionBegin;
226f7bcbcbbSMatthew G. Knepley   ierr = DMPlexRemesh_Internal(dm, metric, bdLabelName, dmAdapt);CHKERRQ(ierr);
2270bbe5d31Sbarral   PetscFunctionReturn(0);
2280bbe5d31Sbarral }
229