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. 14*2e62ab5aSMatthew G. Knepley . remeshBd - Flag to allow boundary changes 15f7bcbcbbSMatthew G. Knepley - bdLabelName - Label name for boundary tags which are preserved in dmNew, or NULL. Should not be "_boundary_". 16f7bcbcbbSMatthew G. Knepley 17f7bcbcbbSMatthew G. Knepley Output Parameter: 18f7bcbcbbSMatthew G. Knepley . dmNew - the new DM 19f7bcbcbbSMatthew G. Knepley 20f7bcbcbbSMatthew G. Knepley Level: advanced 21f7bcbcbbSMatthew G. Knepley 22f7bcbcbbSMatthew G. Knepley .seealso: DMCoarsen(), DMRefine() 23f7bcbcbbSMatthew G. Knepley */ 24*2e62ab5aSMatthew G. Knepley PetscErrorCode DMPlexRemesh_Internal(DM dm, Vec vertexMetric, const char bdLabelName[], PetscBool remeshBd, DM *dmNew) 25f7bcbcbbSMatthew G. Knepley { 26f7bcbcbbSMatthew G. Knepley const char *bdName = "_boundary_"; 27f7bcbcbbSMatthew G. Knepley DM udm, cdm; 282ee120b0SMatthew G. Knepley DMLabel bdLabel = NULL, bdLabelFull; 29f7bcbcbbSMatthew G. Knepley IS bdIS; 30f7bcbcbbSMatthew G. Knepley PetscSection coordSection; 31f7bcbcbbSMatthew G. Knepley Vec coordinates; 32f7bcbcbbSMatthew G. Knepley const PetscScalar *coords, *met; 33f7bcbcbbSMatthew G. Knepley const PetscInt *bdFacesFull; 342ee120b0SMatthew G. Knepley PetscInt *bdFaces, *bdFaceIds; 35f7bcbcbbSMatthew G. Knepley PetscReal *x, *y, *z, *metric; 362ee120b0SMatthew G. Knepley PetscInt *cells; 372ee120b0SMatthew G. Knepley PetscInt dim, cStart, cEnd, numCells, c, coff, vStart, vEnd, numVertices, v; 382ee120b0SMatthew G. Knepley PetscInt off, maxConeSize, numBdFaces, f, bdSize; 39f7bcbcbbSMatthew G. Knepley PetscBool flg; 402ee120b0SMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC 412ee120b0SMatthew G. Knepley DMLabel bdLabelNew; 422ee120b0SMatthew G. Knepley PetscScalar *coordsNew; 432ee120b0SMatthew G. Knepley PetscInt *bdTags; 442ee120b0SMatthew G. Knepley PetscReal *xNew[3] = {NULL, NULL, NULL}; 452ee120b0SMatthew G. Knepley PetscInt *cellsNew; 462ee120b0SMatthew G. Knepley PetscInt d, numCellsNew, numVerticesNew; 472ee120b0SMatthew G. Knepley PetscInt numCornersNew, fStart, fEnd; 482ee120b0SMatthew G. Knepley #endif 49f7bcbcbbSMatthew G. Knepley PetscErrorCode ierr; 50f7bcbcbbSMatthew G. Knepley 51f7bcbcbbSMatthew G. Knepley PetscFunctionBegin; 52f7bcbcbbSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 53f7bcbcbbSMatthew G. Knepley PetscValidHeaderSpecific(vertexMetric, VEC_CLASSID, 2); 54f7bcbcbbSMatthew G. Knepley PetscValidCharPointer(bdLabelName, 3); 55*2e62ab5aSMatthew G. Knepley PetscValidPointer(dmNew, 5); 56f7bcbcbbSMatthew G. Knepley if (bdLabelName) { 57f7bcbcbbSMatthew G. Knepley size_t len; 58f7bcbcbbSMatthew G. Knepley 59f7bcbcbbSMatthew G. Knepley ierr = PetscStrcmp(bdLabelName, bdName, &flg);CHKERRQ(ierr); 60f7bcbcbbSMatthew G. Knepley if (flg) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "\"%s\" cannot be used as label for boundary facets", bdLabelName); 61f7bcbcbbSMatthew G. Knepley ierr = PetscStrlen(bdLabelName, &len);CHKERRQ(ierr); 62f7bcbcbbSMatthew G. Knepley if (len) { 63f7bcbcbbSMatthew G. Knepley ierr = DMGetLabel(dm, bdLabelName, &bdLabel);CHKERRQ(ierr); 64f7bcbcbbSMatthew G. Knepley if (!bdLabel) SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Label \"%s\" does not exist in DM", bdLabelName); 65f7bcbcbbSMatthew G. Knepley } 66f7bcbcbbSMatthew G. Knepley } 67f7bcbcbbSMatthew G. Knepley /* Get mesh information */ 68f7bcbcbbSMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 69f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 70f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 71f7bcbcbbSMatthew G. Knepley ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr); 72f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr); 73f7bcbcbbSMatthew G. Knepley numCells = cEnd - cStart; 74f7bcbcbbSMatthew G. Knepley numVertices = vEnd - vStart; 75f7bcbcbbSMatthew G. Knepley ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr); 76f7bcbcbbSMatthew G. Knepley for (c = 0, coff = 0; c < numCells; ++c) { 77f7bcbcbbSMatthew G. Knepley const PetscInt *cone; 78f7bcbcbbSMatthew G. Knepley PetscInt coneSize, cl; 79f7bcbcbbSMatthew G. Knepley 80f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr); 81f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr); 82f7bcbcbbSMatthew G. Knepley for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart; 83f7bcbcbbSMatthew G. Knepley } 84f7bcbcbbSMatthew G. Knepley ierr = DMDestroy(&udm);CHKERRQ(ierr); 85f7bcbcbbSMatthew G. Knepley ierr = DMGetCoordinateDM(dm, &cdm);CHKERRQ(ierr); 86f7bcbcbbSMatthew G. Knepley ierr = DMGetDefaultSection(cdm, &coordSection);CHKERRQ(ierr); 87f7bcbcbbSMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 88f7bcbcbbSMatthew G. Knepley ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 89f7bcbcbbSMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 90f7bcbcbbSMatthew G. Knepley ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 915f2e139eSMatthew G. Knepley x[v-vStart] = PetscRealPart(coords[off+0]); 925f2e139eSMatthew G. Knepley if (dim > 1) y[v-vStart] = PetscRealPart(coords[off+1]); 935f2e139eSMatthew G. Knepley if (dim > 2) z[v-vStart] = PetscRealPart(coords[off+2]); 94f7bcbcbbSMatthew G. Knepley } 95f7bcbcbbSMatthew G. Knepley ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 96f7bcbcbbSMatthew G. Knepley /* Get boundary mesh */ 97f7bcbcbbSMatthew G. Knepley ierr = DMLabelCreate(bdName, &bdLabelFull);CHKERRQ(ierr); 98f7bcbcbbSMatthew G. Knepley ierr = DMPlexMarkBoundaryFaces(dm, bdLabelFull);CHKERRQ(ierr); 99f7bcbcbbSMatthew G. Knepley ierr = DMLabelGetStratumIS(bdLabelFull, 1, &bdIS);CHKERRQ(ierr); 100f7bcbcbbSMatthew G. Knepley ierr = DMLabelGetStratumSize(bdLabelFull, 1, &numBdFaces);CHKERRQ(ierr); 101f7bcbcbbSMatthew G. Knepley ierr = ISGetIndices(bdIS, &bdFacesFull);CHKERRQ(ierr); 102f7bcbcbbSMatthew G. Knepley for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 103f7bcbcbbSMatthew G. Knepley PetscInt *closure = NULL; 104f7bcbcbbSMatthew G. Knepley PetscInt closureSize, cl; 105f7bcbcbbSMatthew G. Knepley 106f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 107f7bcbcbbSMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 108f7bcbcbbSMatthew G. Knepley if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize; 109f7bcbcbbSMatthew G. Knepley } 110f7bcbcbbSMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 111f7bcbcbbSMatthew G. Knepley } 112f7bcbcbbSMatthew G. Knepley ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr); 113f7bcbcbbSMatthew G. Knepley for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 114f7bcbcbbSMatthew G. Knepley PetscInt *closure = NULL; 115f7bcbcbbSMatthew G. Knepley PetscInt closureSize, cl; 116f7bcbcbbSMatthew G. Knepley 117f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 118f7bcbcbbSMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 119f7bcbcbbSMatthew G. Knepley if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart; 120f7bcbcbbSMatthew G. Knepley } 121f7bcbcbbSMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, bdFacesFull[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 122f7bcbcbbSMatthew G. Knepley if (bdLabel) {ierr = DMLabelGetValue(bdLabel, bdFacesFull[f], &bdFaceIds[f]);CHKERRQ(ierr);} 123f7bcbcbbSMatthew G. Knepley else {bdFaceIds[f] = 1;} 124f7bcbcbbSMatthew G. Knepley } 125f7bcbcbbSMatthew G. Knepley ierr = ISDestroy(&bdIS);CHKERRQ(ierr); 126f7bcbcbbSMatthew G. Knepley ierr = DMLabelDestroy(&bdLabelFull);CHKERRQ(ierr); 127f7bcbcbbSMatthew G. Knepley /* Get metric */ 128f7bcbcbbSMatthew G. Knepley ierr = VecGetArrayRead(vertexMetric, &met);CHKERRQ(ierr); 1295f2e139eSMatthew G. Knepley for (v = 0; v < (vEnd-vStart)*PetscSqr(dim); ++v) metric[v] = PetscRealPart(met[v]); 130f7bcbcbbSMatthew G. Knepley ierr = VecRestoreArrayRead(vertexMetric, &met);CHKERRQ(ierr); 131f7bcbcbbSMatthew G. Knepley /* Create new mesh */ 132f7bcbcbbSMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC 133f7bcbcbbSMatthew G. Knepley switch (dim) { 134f7bcbcbbSMatthew G. Knepley case 2: 135f7bcbcbbSMatthew G. Knepley pragmatic_2d_init(&numVertices, &numCells, cells, x, y);break; 136f7bcbcbbSMatthew G. Knepley case 3: 137f7bcbcbbSMatthew G. Knepley pragmatic_3d_init(&numVertices, &numCells, cells, x, y, z);break; 138f7bcbcbbSMatthew G. Knepley default: SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim); 139f7bcbcbbSMatthew G. Knepley } 140f7bcbcbbSMatthew G. Knepley pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds); 141f7bcbcbbSMatthew G. Knepley pragmatic_set_metric(metric); 142*2e62ab5aSMatthew G. Knepley pragmatic_adapt(remeshBd ? 1 : 0); 143f7bcbcbbSMatthew G. Knepley /* Read out mesh */ 144f7bcbcbbSMatthew G. Knepley pragmatic_get_info(&numVerticesNew, &numCellsNew); 145f7bcbcbbSMatthew G. Knepley ierr = PetscMalloc1(numVerticesNew*dim, &coordsNew);CHKERRQ(ierr); 146f7bcbcbbSMatthew G. Knepley switch (dim) { 147f7bcbcbbSMatthew G. Knepley case 2: 148f7bcbcbbSMatthew G. Knepley numCornersNew = 3; 149f7bcbcbbSMatthew G. Knepley ierr = PetscMalloc2(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1]);CHKERRQ(ierr); 150f7bcbcbbSMatthew G. Knepley pragmatic_get_coords_2d(xNew[0], xNew[1]); 151f7bcbcbbSMatthew G. Knepley break; 152f7bcbcbbSMatthew G. Knepley case 3: 153f7bcbcbbSMatthew G. Knepley numCornersNew = 4; 154f7bcbcbbSMatthew G. Knepley ierr = PetscMalloc3(numVerticesNew, &xNew[0], numVerticesNew, &xNew[1], numVerticesNew, &xNew[2]);CHKERRQ(ierr); 155f7bcbcbbSMatthew G. Knepley pragmatic_get_coords_3d(xNew[0], xNew[1], xNew[2]); 156f7bcbcbbSMatthew G. Knepley break; 157f7bcbcbbSMatthew G. Knepley default: 158f7bcbcbbSMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim); 159f7bcbcbbSMatthew G. Knepley } 160f7bcbcbbSMatthew G. Knepley for (v = 0; v < numVerticesNew; ++v) {for (d = 0; d < dim; ++d) coordsNew[v*dim+d] = xNew[d][v];} 161f7bcbcbbSMatthew G. Knepley ierr = PetscMalloc1(numCellsNew*(dim+1), &cellsNew);CHKERRQ(ierr); 162f7bcbcbbSMatthew G. Knepley pragmatic_get_elements(cellsNew); 163f7bcbcbbSMatthew G. Knepley ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject) dm), dim, numCellsNew, numVerticesNew, numCornersNew, PETSC_TRUE, cellsNew, dim, coordsNew, dmNew);CHKERRQ(ierr); 164f7bcbcbbSMatthew G. Knepley /* Read out boundary label */ 165f7bcbcbbSMatthew G. Knepley pragmatic_get_boundaryTags(&bdTags); 166f7bcbcbbSMatthew G. Knepley ierr = DMCreateLabel(*dmNew, bdLabel ? bdLabelName : bdName);CHKERRQ(ierr); 167f7bcbcbbSMatthew G. Knepley ierr = DMGetLabel(*dmNew, bdLabel ? bdLabelName : bdName, &bdLabelNew);CHKERRQ(ierr); 168f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetHeightStratum(*dmNew, 0, &cStart, &cEnd);CHKERRQ(ierr); 169f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetHeightStratum(*dmNew, 1, &fStart, &fEnd);CHKERRQ(ierr); 170f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetDepthStratum(*dmNew, 0, &vStart, &vEnd);CHKERRQ(ierr); 171f7bcbcbbSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 172f7bcbcbbSMatthew G. Knepley /* Only for simplicial meshes */ 173f7bcbcbbSMatthew G. Knepley coff = (c-cStart)*(dim+1); 174f7bcbcbbSMatthew G. Knepley for (d = 0; d < dim+1; ++d) { 175f7bcbcbbSMatthew G. Knepley if (bdTags[coff+d]) { 176f7bcbcbbSMatthew G. Knepley const PetscInt *faces; 177f7bcbcbbSMatthew G. Knepley PetscInt numFaces, vertices[3], p; 178f7bcbcbbSMatthew G. Knepley 179f7bcbcbbSMatthew G. Knepley /* Mark face opposite to this vertex */ 180f7bcbcbbSMatthew G. Knepley for (p = 0, v = 0; p < dim+1; ++p) if (p != d) {vertices[v] = cellsNew[coff+p] + vStart; ++v;} 181f7bcbcbbSMatthew G. Knepley ierr = DMPlexGetFullJoin(*dmNew, dim, vertices, &numFaces, &faces);CHKERRQ(ierr); 182f7bcbcbbSMatthew 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]); 183f7bcbcbbSMatthew 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]); 184f7bcbcbbSMatthew G. Knepley ierr = DMLabelSetValue(bdLabelNew, faces[0], bdTags[coff+d]);CHKERRQ(ierr); 185f7bcbcbbSMatthew G. Knepley ierr = DMPlexRestoreJoin(*dmNew, dim, vertices, &numFaces, &faces);CHKERRQ(ierr); 186f7bcbcbbSMatthew G. Knepley } 187f7bcbcbbSMatthew G. Knepley } 188f7bcbcbbSMatthew G. Knepley } 189f7bcbcbbSMatthew G. Knepley /* Cleanup */ 190f7bcbcbbSMatthew G. Knepley switch (dim) { 191f7bcbcbbSMatthew G. Knepley case 2: ierr = PetscFree2(xNew[0], xNew[1]);CHKERRQ(ierr);break; 192f7bcbcbbSMatthew G. Knepley case 3: ierr = PetscFree3(xNew[0], xNew[1], xNew[2]);CHKERRQ(ierr);break; 193f7bcbcbbSMatthew G. Knepley } 194f7bcbcbbSMatthew G. Knepley ierr = PetscFree(cellsNew);CHKERRQ(ierr); 195f7bcbcbbSMatthew G. Knepley ierr = PetscFree5(x, y, z, metric, cells);CHKERRQ(ierr); 196f7bcbcbbSMatthew G. Knepley ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr); 197f7bcbcbbSMatthew G. Knepley ierr = PetscFree(coordsNew);CHKERRQ(ierr); 198f7bcbcbbSMatthew G. Knepley pragmatic_finalize(); 199f7bcbcbbSMatthew G. Knepley #else 200f7bcbcbbSMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Remeshing needs external package support.\nPlease reconfigure with --download-pragmatic."); 201f7bcbcbbSMatthew G. Knepley #endif 202f7bcbcbbSMatthew G. Knepley PetscFunctionReturn(0); 203f7bcbcbbSMatthew G. Knepley } 2040bbe5d31Sbarral 2050bbe5d31Sbarral #undef __FUNCT__ 20652b40ca0Sbarral #define __FUNCT__ "DMPlexAdapt" 2070f7e110dSbarral /*@ 2080f7e110dSbarral DMPlexAdapt - Generates a mesh adapted to the specified metric field using the pragmatic library. 2090f7e110dSbarral 2100f7e110dSbarral Input Parameters: 2110f7e110dSbarral + dm - The DM object 2120f7e110dSbarral . metric - The metric to which the mesh is adapted, defined vertex-wise. 213379fadc5Sbarral - 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". 2140f7e110dSbarral 2150f7e110dSbarral Output Parameter: 216f7bcbcbbSMatthew G. Knepley . dmAdapt - Pointer to the DM object containing the adapted mesh 2170f7e110dSbarral 2180f7e110dSbarral Level: advanced 2190f7e110dSbarral 2200f7e110dSbarral .seealso: DMCoarsen(), DMRefine() 2210f7e110dSbarral @*/ 222f7bcbcbbSMatthew G. Knepley PetscErrorCode DMPlexAdapt(DM dm, Vec metric, const char bdLabelName[], DM *dmAdapt) 2230bbe5d31Sbarral { 224*2e62ab5aSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 2250bbe5d31Sbarral PetscErrorCode ierr; 2260bbe5d31Sbarral 2270bbe5d31Sbarral PetscFunctionBegin; 228*2e62ab5aSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 229*2e62ab5aSMatthew G. Knepley PetscValidHeaderSpecific(metric, VEC_CLASSID, 2); 230*2e62ab5aSMatthew G. Knepley PetscValidCharPointer(bdLabelName, 3); 231*2e62ab5aSMatthew G. Knepley PetscValidPointer(dmAdapt, 4); 232*2e62ab5aSMatthew G. Knepley ierr = DMPlexRemesh_Internal(dm, metric, bdLabelName, mesh->remeshBd, dmAdapt);CHKERRQ(ierr); 2330bbe5d31Sbarral PetscFunctionReturn(0); 2340bbe5d31Sbarral } 235