1c6f61ee2SBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2c6f61ee2SBarry Smith #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/ 3c6f61ee2SBarry Smith #include <petscds.h> 4c6f61ee2SBarry Smith #include <petsc/private/petscimpl.h> 5c6f61ee2SBarry Smith #include <petsc/private/petscfeimpl.h> 6c6f61ee2SBarry Smith 7c6f61ee2SBarry Smith /*@C 8c6f61ee2SBarry Smith DMInterpolationCreate - Creates a `DMInterpolationInfo` context 9c6f61ee2SBarry Smith 10c6f61ee2SBarry Smith Collective 11c6f61ee2SBarry Smith 12c6f61ee2SBarry Smith Input Parameter: 13c6f61ee2SBarry Smith . comm - the communicator 14c6f61ee2SBarry Smith 15c6f61ee2SBarry Smith Output Parameter: 16c6f61ee2SBarry Smith . ctx - the context 17c6f61ee2SBarry Smith 18c6f61ee2SBarry Smith Level: beginner 19c6f61ee2SBarry Smith 20df514481SBarry Smith Developer Note: 21df514481SBarry Smith The naming is incorrect, either the object should be named `DMInterpolation` or all the routines should begin with `DMInterpolationInfo` 22df514481SBarry Smith 23df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationDestroy()` 24c6f61ee2SBarry Smith @*/ 25c6f61ee2SBarry Smith PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx) 26c6f61ee2SBarry Smith { 27c6f61ee2SBarry Smith PetscFunctionBegin; 28c6f61ee2SBarry Smith PetscAssertPointer(ctx, 2); 29c6f61ee2SBarry Smith PetscCall(PetscNew(ctx)); 30c6f61ee2SBarry Smith 31c6f61ee2SBarry Smith (*ctx)->comm = comm; 32c6f61ee2SBarry Smith (*ctx)->dim = -1; 33c6f61ee2SBarry Smith (*ctx)->nInput = 0; 34c6f61ee2SBarry Smith (*ctx)->points = NULL; 35c6f61ee2SBarry Smith (*ctx)->cells = NULL; 36c6f61ee2SBarry Smith (*ctx)->n = -1; 37c6f61ee2SBarry Smith (*ctx)->coords = NULL; 38c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 39c6f61ee2SBarry Smith } 40c6f61ee2SBarry Smith 41c6f61ee2SBarry Smith /*@C 42c6f61ee2SBarry Smith DMInterpolationSetDim - Sets the spatial dimension for the interpolation context 43c6f61ee2SBarry Smith 44c6f61ee2SBarry Smith Not Collective 45c6f61ee2SBarry Smith 46c6f61ee2SBarry Smith Input Parameters: 47c6f61ee2SBarry Smith + ctx - the context 48c6f61ee2SBarry Smith - dim - the spatial dimension 49c6f61ee2SBarry Smith 50c6f61ee2SBarry Smith Level: intermediate 51c6f61ee2SBarry Smith 52df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMInterpolationInfo`, `DMInterpolationGetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()` 53c6f61ee2SBarry Smith @*/ 54c6f61ee2SBarry Smith PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim) 55c6f61ee2SBarry Smith { 56c6f61ee2SBarry Smith PetscFunctionBegin; 57c6f61ee2SBarry Smith PetscCheck(!(dim < 1) && !(dim > 3), ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %" PetscInt_FMT, dim); 58c6f61ee2SBarry Smith ctx->dim = dim; 59c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 60c6f61ee2SBarry Smith } 61c6f61ee2SBarry Smith 62c6f61ee2SBarry Smith /*@C 63c6f61ee2SBarry Smith DMInterpolationGetDim - Gets the spatial dimension for the interpolation context 64c6f61ee2SBarry Smith 65c6f61ee2SBarry Smith Not Collective 66c6f61ee2SBarry Smith 67c6f61ee2SBarry Smith Input Parameter: 68c6f61ee2SBarry Smith . ctx - the context 69c6f61ee2SBarry Smith 70c6f61ee2SBarry Smith Output Parameter: 71c6f61ee2SBarry Smith . dim - the spatial dimension 72c6f61ee2SBarry Smith 73c6f61ee2SBarry Smith Level: intermediate 74c6f61ee2SBarry Smith 75df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMInterpolationInfo`, `DMInterpolationSetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()` 76c6f61ee2SBarry Smith @*/ 77c6f61ee2SBarry Smith PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim) 78c6f61ee2SBarry Smith { 79c6f61ee2SBarry Smith PetscFunctionBegin; 80c6f61ee2SBarry Smith PetscAssertPointer(dim, 2); 81c6f61ee2SBarry Smith *dim = ctx->dim; 82c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 83c6f61ee2SBarry Smith } 84c6f61ee2SBarry Smith 85c6f61ee2SBarry Smith /*@C 86c6f61ee2SBarry Smith DMInterpolationSetDof - Sets the number of fields interpolated at a point for the interpolation context 87c6f61ee2SBarry Smith 88c6f61ee2SBarry Smith Not Collective 89c6f61ee2SBarry Smith 90c6f61ee2SBarry Smith Input Parameters: 91c6f61ee2SBarry Smith + ctx - the context 92c6f61ee2SBarry Smith - dof - the number of fields 93c6f61ee2SBarry Smith 94c6f61ee2SBarry Smith Level: intermediate 95c6f61ee2SBarry Smith 96df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMInterpolationInfo`, `DMInterpolationGetDof()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()` 97c6f61ee2SBarry Smith @*/ 98c6f61ee2SBarry Smith PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof) 99c6f61ee2SBarry Smith { 100c6f61ee2SBarry Smith PetscFunctionBegin; 101c6f61ee2SBarry Smith PetscCheck(dof >= 1, ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %" PetscInt_FMT, dof); 102c6f61ee2SBarry Smith ctx->dof = dof; 103c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 104c6f61ee2SBarry Smith } 105c6f61ee2SBarry Smith 106c6f61ee2SBarry Smith /*@C 107c6f61ee2SBarry Smith DMInterpolationGetDof - Gets the number of fields interpolated at a point for the interpolation context 108c6f61ee2SBarry Smith 109c6f61ee2SBarry Smith Not Collective 110c6f61ee2SBarry Smith 111c6f61ee2SBarry Smith Input Parameter: 112c6f61ee2SBarry Smith . ctx - the context 113c6f61ee2SBarry Smith 114c6f61ee2SBarry Smith Output Parameter: 115c6f61ee2SBarry Smith . dof - the number of fields 116c6f61ee2SBarry Smith 117c6f61ee2SBarry Smith Level: intermediate 118c6f61ee2SBarry Smith 119df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMInterpolationInfo`, `DMInterpolationSetDof()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()` 120c6f61ee2SBarry Smith @*/ 121c6f61ee2SBarry Smith PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof) 122c6f61ee2SBarry Smith { 123c6f61ee2SBarry Smith PetscFunctionBegin; 124c6f61ee2SBarry Smith PetscAssertPointer(dof, 2); 125c6f61ee2SBarry Smith *dof = ctx->dof; 126c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 127c6f61ee2SBarry Smith } 128c6f61ee2SBarry Smith 129c6f61ee2SBarry Smith /*@C 130c6f61ee2SBarry Smith DMInterpolationAddPoints - Add points at which we will interpolate the fields 131c6f61ee2SBarry Smith 132c6f61ee2SBarry Smith Not Collective 133c6f61ee2SBarry Smith 134c6f61ee2SBarry Smith Input Parameters: 135c6f61ee2SBarry Smith + ctx - the context 136c6f61ee2SBarry Smith . n - the number of points 137df514481SBarry Smith - points - the coordinates for each point, an array of size `n` * dim 138c6f61ee2SBarry Smith 139c6f61ee2SBarry Smith Level: intermediate 140c6f61ee2SBarry Smith 141c6f61ee2SBarry Smith Note: 142df514481SBarry Smith The input coordinate information is copied into the object. 143c6f61ee2SBarry Smith 144df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMInterpolationInfo`, `DMInterpolationSetDim()`, `DMInterpolationEvaluate()`, `DMInterpolationCreate()` 145c6f61ee2SBarry Smith @*/ 146c6f61ee2SBarry Smith PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[]) 147c6f61ee2SBarry Smith { 148c6f61ee2SBarry Smith PetscFunctionBegin; 149c6f61ee2SBarry Smith PetscCheck(ctx->dim >= 0, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 150c6f61ee2SBarry Smith PetscCheck(!ctx->points, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet"); 151c6f61ee2SBarry Smith ctx->nInput = n; 152c6f61ee2SBarry Smith 153c6f61ee2SBarry Smith PetscCall(PetscMalloc1(n * ctx->dim, &ctx->points)); 154c6f61ee2SBarry Smith PetscCall(PetscArraycpy(ctx->points, points, n * ctx->dim)); 155c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 156c6f61ee2SBarry Smith } 157c6f61ee2SBarry Smith 158c6f61ee2SBarry Smith /*@C 159c6f61ee2SBarry Smith DMInterpolationSetUp - Compute spatial indices for point location during interpolation 160c6f61ee2SBarry Smith 161c6f61ee2SBarry Smith Collective 162c6f61ee2SBarry Smith 163c6f61ee2SBarry Smith Input Parameters: 164c6f61ee2SBarry Smith + ctx - the context 165c6f61ee2SBarry Smith . dm - the `DM` for the function space used for interpolation 166c6f61ee2SBarry Smith . redundantPoints - If `PETSC_TRUE`, all processes are passing in the same array of points. Otherwise, points need to be communicated among processes. 167c6f61ee2SBarry Smith - ignoreOutsideDomain - If `PETSC_TRUE`, ignore points outside the domain, otherwise return an error 168c6f61ee2SBarry Smith 169c6f61ee2SBarry Smith Level: intermediate 170c6f61ee2SBarry Smith 171df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 172c6f61ee2SBarry Smith @*/ 173c6f61ee2SBarry Smith PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints, PetscBool ignoreOutsideDomain) 174c6f61ee2SBarry Smith { 175c6f61ee2SBarry Smith MPI_Comm comm = ctx->comm; 176c6f61ee2SBarry Smith PetscScalar *a; 177c6f61ee2SBarry Smith PetscInt p, q, i; 178c6f61ee2SBarry Smith PetscMPIInt rank, size; 179c6f61ee2SBarry Smith Vec pointVec; 180c6f61ee2SBarry Smith PetscSF cellSF; 181c6f61ee2SBarry Smith PetscLayout layout; 182c6f61ee2SBarry Smith PetscReal *globalPoints; 183c6f61ee2SBarry Smith PetscScalar *globalPointsScalar; 184c6f61ee2SBarry Smith const PetscInt *ranges; 185*e91c04dfSPierre Jolivet PetscMPIInt *counts, *displs; 186c6f61ee2SBarry Smith const PetscSFNode *foundCells; 187c6f61ee2SBarry Smith const PetscInt *foundPoints; 1886497c311SBarry Smith PetscMPIInt *foundProcs, *globalProcs, in; 189c6f61ee2SBarry Smith PetscInt n, N, numFound; 190c6f61ee2SBarry Smith 191c6f61ee2SBarry Smith PetscFunctionBegin; 192c6f61ee2SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 193c6f61ee2SBarry Smith PetscCallMPI(MPI_Comm_size(comm, &size)); 194c6f61ee2SBarry Smith PetscCallMPI(MPI_Comm_rank(comm, &rank)); 195c6f61ee2SBarry Smith PetscCheck(ctx->dim >= 0, comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 196c6f61ee2SBarry Smith /* Locate points */ 197c6f61ee2SBarry Smith n = ctx->nInput; 198c6f61ee2SBarry Smith if (!redundantPoints) { 199c6f61ee2SBarry Smith PetscCall(PetscLayoutCreate(comm, &layout)); 200c6f61ee2SBarry Smith PetscCall(PetscLayoutSetBlockSize(layout, 1)); 201c6f61ee2SBarry Smith PetscCall(PetscLayoutSetLocalSize(layout, n)); 202c6f61ee2SBarry Smith PetscCall(PetscLayoutSetUp(layout)); 203c6f61ee2SBarry Smith PetscCall(PetscLayoutGetSize(layout, &N)); 204c6f61ee2SBarry Smith /* Communicate all points to all processes */ 205c6f61ee2SBarry Smith PetscCall(PetscMalloc3(N * ctx->dim, &globalPoints, size, &counts, size, &displs)); 206c6f61ee2SBarry Smith PetscCall(PetscLayoutGetRanges(layout, &ranges)); 207c6f61ee2SBarry Smith for (p = 0; p < size; ++p) { 2086497c311SBarry Smith PetscCall(PetscMPIIntCast((ranges[p + 1] - ranges[p]) * ctx->dim, &counts[p])); 2096497c311SBarry Smith PetscCall(PetscMPIIntCast(ranges[p] * ctx->dim, &displs[p])); 210c6f61ee2SBarry Smith } 2116497c311SBarry Smith PetscCall(PetscMPIIntCast(n * ctx->dim, &in)); 2126497c311SBarry Smith PetscCallMPI(MPI_Allgatherv(ctx->points, in, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm)); 213c6f61ee2SBarry Smith } else { 214c6f61ee2SBarry Smith N = n; 215c6f61ee2SBarry Smith globalPoints = ctx->points; 216c6f61ee2SBarry Smith counts = displs = NULL; 217c6f61ee2SBarry Smith layout = NULL; 218c6f61ee2SBarry Smith } 219c6f61ee2SBarry Smith #if 0 220c6f61ee2SBarry Smith PetscCall(PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs)); 221c6f61ee2SBarry Smith /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */ 222c6f61ee2SBarry Smith #else 223c6f61ee2SBarry Smith #if defined(PETSC_USE_COMPLEX) 224c6f61ee2SBarry Smith PetscCall(PetscMalloc1(N * ctx->dim, &globalPointsScalar)); 225c6f61ee2SBarry Smith for (i = 0; i < N * ctx->dim; i++) globalPointsScalar[i] = globalPoints[i]; 226c6f61ee2SBarry Smith #else 227c6f61ee2SBarry Smith globalPointsScalar = globalPoints; 228c6f61ee2SBarry Smith #endif 229c6f61ee2SBarry Smith PetscCall(VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N * ctx->dim, globalPointsScalar, &pointVec)); 230c6f61ee2SBarry Smith PetscCall(PetscMalloc2(N, &foundProcs, N, &globalProcs)); 231c6f61ee2SBarry Smith for (p = 0; p < N; ++p) foundProcs[p] = size; 232c6f61ee2SBarry Smith cellSF = NULL; 233c6f61ee2SBarry Smith PetscCall(DMLocatePoints(dm, pointVec, DM_POINTLOCATION_REMOVE, &cellSF)); 234c6f61ee2SBarry Smith PetscCall(PetscSFGetGraph(cellSF, NULL, &numFound, &foundPoints, &foundCells)); 235c6f61ee2SBarry Smith #endif 236c6f61ee2SBarry Smith for (p = 0; p < numFound; ++p) { 237c6f61ee2SBarry Smith if (foundCells[p].index >= 0) foundProcs[foundPoints ? foundPoints[p] : p] = rank; 238c6f61ee2SBarry Smith } 239c6f61ee2SBarry Smith /* Let the lowest rank process own each point */ 240*e91c04dfSPierre Jolivet PetscCallMPI(MPIU_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm)); 241c6f61ee2SBarry Smith ctx->n = 0; 242c6f61ee2SBarry Smith for (p = 0; p < N; ++p) { 243c6f61ee2SBarry Smith if (globalProcs[p] == size) { 244835f2295SStefano Zampini PetscCheck(ignoreOutsideDomain, comm, PETSC_ERR_PLIB, "Point %" PetscInt_FMT ": %g %g %g not located in mesh", p, (double)globalPoints[p * ctx->dim + 0], (double)(ctx->dim > 1 ? globalPoints[p * ctx->dim + 1] : 0), 245835f2295SStefano Zampini (double)(ctx->dim > 2 ? globalPoints[p * ctx->dim + 2] : 0)); 246c6f61ee2SBarry Smith if (rank == 0) ++ctx->n; 247c6f61ee2SBarry Smith } else if (globalProcs[p] == rank) ++ctx->n; 248c6f61ee2SBarry Smith } 249c6f61ee2SBarry Smith /* Create coordinates vector and array of owned cells */ 250c6f61ee2SBarry Smith PetscCall(PetscMalloc1(ctx->n, &ctx->cells)); 251c6f61ee2SBarry Smith PetscCall(VecCreate(comm, &ctx->coords)); 252c6f61ee2SBarry Smith PetscCall(VecSetSizes(ctx->coords, ctx->n * ctx->dim, PETSC_DECIDE)); 253c6f61ee2SBarry Smith PetscCall(VecSetBlockSize(ctx->coords, ctx->dim)); 254c6f61ee2SBarry Smith PetscCall(VecSetType(ctx->coords, VECSTANDARD)); 255c6f61ee2SBarry Smith PetscCall(VecGetArray(ctx->coords, &a)); 256c6f61ee2SBarry Smith for (p = 0, q = 0, i = 0; p < N; ++p) { 257c6f61ee2SBarry Smith if (globalProcs[p] == rank) { 258c6f61ee2SBarry Smith PetscInt d; 259c6f61ee2SBarry Smith 260c6f61ee2SBarry Smith for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p * ctx->dim + d]; 261c6f61ee2SBarry Smith ctx->cells[q] = foundCells[q].index; 262c6f61ee2SBarry Smith ++q; 263c6f61ee2SBarry Smith } 264c6f61ee2SBarry Smith if (globalProcs[p] == size && rank == 0) { 265c6f61ee2SBarry Smith PetscInt d; 266c6f61ee2SBarry Smith 267c6f61ee2SBarry Smith for (d = 0; d < ctx->dim; ++d, ++i) a[i] = 0.; 268c6f61ee2SBarry Smith ctx->cells[q] = -1; 269c6f61ee2SBarry Smith ++q; 270c6f61ee2SBarry Smith } 271c6f61ee2SBarry Smith } 272c6f61ee2SBarry Smith PetscCall(VecRestoreArray(ctx->coords, &a)); 273c6f61ee2SBarry Smith #if 0 274c6f61ee2SBarry Smith PetscCall(PetscFree3(foundCells,foundProcs,globalProcs)); 275c6f61ee2SBarry Smith #else 276c6f61ee2SBarry Smith PetscCall(PetscFree2(foundProcs, globalProcs)); 277c6f61ee2SBarry Smith PetscCall(PetscSFDestroy(&cellSF)); 278c6f61ee2SBarry Smith PetscCall(VecDestroy(&pointVec)); 279c6f61ee2SBarry Smith #endif 280c6f61ee2SBarry Smith if ((void *)globalPointsScalar != (void *)globalPoints) PetscCall(PetscFree(globalPointsScalar)); 281c6f61ee2SBarry Smith if (!redundantPoints) PetscCall(PetscFree3(globalPoints, counts, displs)); 282c6f61ee2SBarry Smith PetscCall(PetscLayoutDestroy(&layout)); 283c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 284c6f61ee2SBarry Smith } 285c6f61ee2SBarry Smith 286c6f61ee2SBarry Smith /*@C 287c6f61ee2SBarry Smith DMInterpolationGetCoordinates - Gets a `Vec` with the coordinates of each interpolation point 288c6f61ee2SBarry Smith 289c6f61ee2SBarry Smith Collective 290c6f61ee2SBarry Smith 291c6f61ee2SBarry Smith Input Parameter: 292c6f61ee2SBarry Smith . ctx - the context 293c6f61ee2SBarry Smith 294c6f61ee2SBarry Smith Output Parameter: 295c6f61ee2SBarry Smith . coordinates - the coordinates of interpolation points 296c6f61ee2SBarry Smith 297c6f61ee2SBarry Smith Level: intermediate 298c6f61ee2SBarry Smith 299c6f61ee2SBarry Smith Note: 300c6f61ee2SBarry Smith The local vector entries correspond to interpolation points lying on this process, according to the associated `DM`. 301c6f61ee2SBarry Smith This is a borrowed vector that the user should not destroy. 302c6f61ee2SBarry Smith 303df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 304c6f61ee2SBarry Smith @*/ 305c6f61ee2SBarry Smith PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates) 306c6f61ee2SBarry Smith { 307c6f61ee2SBarry Smith PetscFunctionBegin; 308c6f61ee2SBarry Smith PetscAssertPointer(coordinates, 2); 309c6f61ee2SBarry Smith PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 310c6f61ee2SBarry Smith *coordinates = ctx->coords; 311c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 312c6f61ee2SBarry Smith } 313c6f61ee2SBarry Smith 314c6f61ee2SBarry Smith /*@C 315c6f61ee2SBarry Smith DMInterpolationGetVector - Gets a `Vec` which can hold all the interpolated field values 316c6f61ee2SBarry Smith 317c6f61ee2SBarry Smith Collective 318c6f61ee2SBarry Smith 319c6f61ee2SBarry Smith Input Parameter: 320c6f61ee2SBarry Smith . ctx - the context 321c6f61ee2SBarry Smith 322c6f61ee2SBarry Smith Output Parameter: 323c6f61ee2SBarry Smith . v - a vector capable of holding the interpolated field values 324c6f61ee2SBarry Smith 325c6f61ee2SBarry Smith Level: intermediate 326c6f61ee2SBarry Smith 327c6f61ee2SBarry Smith Note: 328c6f61ee2SBarry Smith This vector should be returned using `DMInterpolationRestoreVector()`. 329c6f61ee2SBarry Smith 330df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMInterpolationInfo`, `DMInterpolationRestoreVector()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 331c6f61ee2SBarry Smith @*/ 332c6f61ee2SBarry Smith PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v) 333c6f61ee2SBarry Smith { 334c6f61ee2SBarry Smith PetscFunctionBegin; 335c6f61ee2SBarry Smith PetscAssertPointer(v, 2); 336c6f61ee2SBarry Smith PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 337c6f61ee2SBarry Smith PetscCall(VecCreate(ctx->comm, v)); 338c6f61ee2SBarry Smith PetscCall(VecSetSizes(*v, ctx->n * ctx->dof, PETSC_DECIDE)); 339c6f61ee2SBarry Smith PetscCall(VecSetBlockSize(*v, ctx->dof)); 340c6f61ee2SBarry Smith PetscCall(VecSetType(*v, VECSTANDARD)); 341c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 342c6f61ee2SBarry Smith } 343c6f61ee2SBarry Smith 344c6f61ee2SBarry Smith /*@C 345c6f61ee2SBarry Smith DMInterpolationRestoreVector - Returns a `Vec` which can hold all the interpolated field values 346c6f61ee2SBarry Smith 347c6f61ee2SBarry Smith Collective 348c6f61ee2SBarry Smith 349c6f61ee2SBarry Smith Input Parameters: 350c6f61ee2SBarry Smith + ctx - the context 351c6f61ee2SBarry Smith - v - a vector capable of holding the interpolated field values 352c6f61ee2SBarry Smith 353c6f61ee2SBarry Smith Level: intermediate 354c6f61ee2SBarry Smith 355df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMInterpolationInfo`, `DMInterpolationGetVector()`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 356c6f61ee2SBarry Smith @*/ 357c6f61ee2SBarry Smith PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v) 358c6f61ee2SBarry Smith { 359c6f61ee2SBarry Smith PetscFunctionBegin; 360c6f61ee2SBarry Smith PetscAssertPointer(v, 2); 361c6f61ee2SBarry Smith PetscCheck(ctx->coords, ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 362c6f61ee2SBarry Smith PetscCall(VecDestroy(v)); 363c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 364c6f61ee2SBarry Smith } 365c6f61ee2SBarry Smith 366c6f61ee2SBarry Smith static inline PetscErrorCode DMInterpolate_Segment_Private(DMInterpolationInfo ctx, DM dm, PetscInt p, Vec xLocal, Vec v) 367c6f61ee2SBarry Smith { 368c6f61ee2SBarry Smith const PetscInt c = ctx->cells[p]; 369c6f61ee2SBarry Smith const PetscInt dof = ctx->dof; 370c6f61ee2SBarry Smith PetscScalar *x = NULL; 371c6f61ee2SBarry Smith const PetscScalar *coords; 372c6f61ee2SBarry Smith PetscScalar *a; 373c6f61ee2SBarry Smith PetscReal v0, J, invJ, detJ, xir[1]; 374c6f61ee2SBarry Smith PetscInt xSize, comp; 375c6f61ee2SBarry Smith 376c6f61ee2SBarry Smith PetscFunctionBegin; 377c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(ctx->coords, &coords)); 378c6f61ee2SBarry Smith PetscCall(VecGetArray(v, &a)); 379c6f61ee2SBarry Smith PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, &v0, &J, &invJ, &detJ)); 380c6f61ee2SBarry Smith PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c); 381c6f61ee2SBarry Smith xir[0] = invJ * PetscRealPart(coords[p] - v0); 382c6f61ee2SBarry Smith PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x)); 383c6f61ee2SBarry Smith if (2 * dof == xSize) { 384c6f61ee2SBarry Smith for (comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp] * (1 - xir[0]) + x[1 * dof + comp] * xir[0]; 385c6f61ee2SBarry Smith } else if (dof == xSize) { 386c6f61ee2SBarry Smith for (comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp]; 387c6f61ee2SBarry Smith } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Input closure size %" PetscInt_FMT " must be either %" PetscInt_FMT " or %" PetscInt_FMT, xSize, 2 * dof, dof); 388c6f61ee2SBarry Smith PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x)); 389c6f61ee2SBarry Smith PetscCall(VecRestoreArray(v, &a)); 390c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 391c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 392c6f61ee2SBarry Smith } 393c6f61ee2SBarry Smith 394c6f61ee2SBarry Smith static inline PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, PetscInt p, Vec xLocal, Vec v) 395c6f61ee2SBarry Smith { 396c6f61ee2SBarry Smith const PetscInt c = ctx->cells[p]; 397c6f61ee2SBarry Smith PetscScalar *x = NULL; 398c6f61ee2SBarry Smith const PetscScalar *coords; 399c6f61ee2SBarry Smith PetscScalar *a; 400c6f61ee2SBarry Smith PetscReal *v0, *J, *invJ, detJ; 401c6f61ee2SBarry Smith PetscReal xi[4]; 402c6f61ee2SBarry Smith 403c6f61ee2SBarry Smith PetscFunctionBegin; 404c6f61ee2SBarry Smith PetscCall(PetscMalloc3(ctx->dim, &v0, ctx->dim * ctx->dim, &J, ctx->dim * ctx->dim, &invJ)); 405c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(ctx->coords, &coords)); 406c6f61ee2SBarry Smith PetscCall(VecGetArray(v, &a)); 407c6f61ee2SBarry Smith PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 408c6f61ee2SBarry Smith PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c); 409c6f61ee2SBarry Smith PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x)); 410c6f61ee2SBarry Smith for (PetscInt comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp]; 411c6f61ee2SBarry Smith for (PetscInt d = 0; d < ctx->dim; ++d) { 412c6f61ee2SBarry Smith xi[d] = 0.0; 413c6f61ee2SBarry Smith for (PetscInt f = 0; f < ctx->dim; ++f) xi[d] += invJ[d * ctx->dim + f] * 0.5 * PetscRealPart(coords[p * ctx->dim + f] - v0[f]); 414c6f61ee2SBarry Smith for (PetscInt comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] += PetscRealPart(x[(d + 1) * ctx->dof + comp] - x[0 * ctx->dof + comp]) * xi[d]; 415c6f61ee2SBarry Smith } 416c6f61ee2SBarry Smith PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x)); 417c6f61ee2SBarry Smith PetscCall(VecRestoreArray(v, &a)); 418c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 419c6f61ee2SBarry Smith PetscCall(PetscFree3(v0, J, invJ)); 420c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 421c6f61ee2SBarry Smith } 422c6f61ee2SBarry Smith 423c6f61ee2SBarry Smith static inline PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, PetscInt p, Vec xLocal, Vec v) 424c6f61ee2SBarry Smith { 425c6f61ee2SBarry Smith const PetscInt c = ctx->cells[p]; 426c6f61ee2SBarry Smith const PetscInt order[3] = {2, 1, 3}; 427c6f61ee2SBarry Smith PetscScalar *x = NULL; 428c6f61ee2SBarry Smith PetscReal *v0, *J, *invJ, detJ; 429c6f61ee2SBarry Smith const PetscScalar *coords; 430c6f61ee2SBarry Smith PetscScalar *a; 431c6f61ee2SBarry Smith PetscReal xi[4]; 432c6f61ee2SBarry Smith 433c6f61ee2SBarry Smith PetscFunctionBegin; 434c6f61ee2SBarry Smith PetscCall(PetscMalloc3(ctx->dim, &v0, ctx->dim * ctx->dim, &J, ctx->dim * ctx->dim, &invJ)); 435c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(ctx->coords, &coords)); 436c6f61ee2SBarry Smith PetscCall(VecGetArray(v, &a)); 437c6f61ee2SBarry Smith PetscCall(DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ)); 438c6f61ee2SBarry Smith PetscCheck(detJ > 0.0, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %" PetscInt_FMT, (double)detJ, c); 439c6f61ee2SBarry Smith PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x)); 440c6f61ee2SBarry Smith for (PetscInt comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp]; 441c6f61ee2SBarry Smith for (PetscInt d = 0; d < ctx->dim; ++d) { 442c6f61ee2SBarry Smith xi[d] = 0.0; 443c6f61ee2SBarry Smith for (PetscInt f = 0; f < ctx->dim; ++f) xi[d] += invJ[d * ctx->dim + f] * 0.5 * PetscRealPart(coords[p * ctx->dim + f] - v0[f]); 444c6f61ee2SBarry Smith for (PetscInt comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] += PetscRealPart(x[order[d] * ctx->dof + comp] - x[0 * ctx->dof + comp]) * xi[d]; 445c6f61ee2SBarry Smith } 446c6f61ee2SBarry Smith PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x)); 447c6f61ee2SBarry Smith PetscCall(VecRestoreArray(v, &a)); 448c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 449c6f61ee2SBarry Smith PetscCall(PetscFree3(v0, J, invJ)); 450c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 451c6f61ee2SBarry Smith } 452c6f61ee2SBarry Smith 453c6f61ee2SBarry Smith static inline PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) 454c6f61ee2SBarry Smith { 455c6f61ee2SBarry Smith const PetscScalar *vertices = (const PetscScalar *)ctx; 456c6f61ee2SBarry Smith const PetscScalar x0 = vertices[0]; 457c6f61ee2SBarry Smith const PetscScalar y0 = vertices[1]; 458c6f61ee2SBarry Smith const PetscScalar x1 = vertices[2]; 459c6f61ee2SBarry Smith const PetscScalar y1 = vertices[3]; 460c6f61ee2SBarry Smith const PetscScalar x2 = vertices[4]; 461c6f61ee2SBarry Smith const PetscScalar y2 = vertices[5]; 462c6f61ee2SBarry Smith const PetscScalar x3 = vertices[6]; 463c6f61ee2SBarry Smith const PetscScalar y3 = vertices[7]; 464c6f61ee2SBarry Smith const PetscScalar f_1 = x1 - x0; 465c6f61ee2SBarry Smith const PetscScalar g_1 = y1 - y0; 466c6f61ee2SBarry Smith const PetscScalar f_3 = x3 - x0; 467c6f61ee2SBarry Smith const PetscScalar g_3 = y3 - y0; 468c6f61ee2SBarry Smith const PetscScalar f_01 = x2 - x1 - x3 + x0; 469c6f61ee2SBarry Smith const PetscScalar g_01 = y2 - y1 - y3 + y0; 470c6f61ee2SBarry Smith const PetscScalar *ref; 471c6f61ee2SBarry Smith PetscScalar *real; 472c6f61ee2SBarry Smith 473c6f61ee2SBarry Smith PetscFunctionBegin; 474c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(Xref, &ref)); 475c6f61ee2SBarry Smith PetscCall(VecGetArray(Xreal, &real)); 476c6f61ee2SBarry Smith { 477c6f61ee2SBarry Smith const PetscScalar p0 = ref[0]; 478c6f61ee2SBarry Smith const PetscScalar p1 = ref[1]; 479c6f61ee2SBarry Smith 480c6f61ee2SBarry Smith real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1; 481c6f61ee2SBarry Smith real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1; 482c6f61ee2SBarry Smith } 483c6f61ee2SBarry Smith PetscCall(PetscLogFlops(28)); 484c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(Xref, &ref)); 485c6f61ee2SBarry Smith PetscCall(VecRestoreArray(Xreal, &real)); 486c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 487c6f61ee2SBarry Smith } 488c6f61ee2SBarry Smith 489c6f61ee2SBarry Smith #include <petsc/private/dmimpl.h> 490c6f61ee2SBarry Smith static inline PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx) 491c6f61ee2SBarry Smith { 492c6f61ee2SBarry Smith const PetscScalar *vertices = (const PetscScalar *)ctx; 493c6f61ee2SBarry Smith const PetscScalar x0 = vertices[0]; 494c6f61ee2SBarry Smith const PetscScalar y0 = vertices[1]; 495c6f61ee2SBarry Smith const PetscScalar x1 = vertices[2]; 496c6f61ee2SBarry Smith const PetscScalar y1 = vertices[3]; 497c6f61ee2SBarry Smith const PetscScalar x2 = vertices[4]; 498c6f61ee2SBarry Smith const PetscScalar y2 = vertices[5]; 499c6f61ee2SBarry Smith const PetscScalar x3 = vertices[6]; 500c6f61ee2SBarry Smith const PetscScalar y3 = vertices[7]; 501c6f61ee2SBarry Smith const PetscScalar f_01 = x2 - x1 - x3 + x0; 502c6f61ee2SBarry Smith const PetscScalar g_01 = y2 - y1 - y3 + y0; 503c6f61ee2SBarry Smith const PetscScalar *ref; 504c6f61ee2SBarry Smith 505c6f61ee2SBarry Smith PetscFunctionBegin; 506c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(Xref, &ref)); 507c6f61ee2SBarry Smith { 508c6f61ee2SBarry Smith const PetscScalar x = ref[0]; 509c6f61ee2SBarry Smith const PetscScalar y = ref[1]; 510c6f61ee2SBarry Smith const PetscInt rows[2] = {0, 1}; 511c6f61ee2SBarry Smith PetscScalar values[4]; 512c6f61ee2SBarry Smith 513c6f61ee2SBarry Smith values[0] = (x1 - x0 + f_01 * y) * 0.5; 514c6f61ee2SBarry Smith values[1] = (x3 - x0 + f_01 * x) * 0.5; 515c6f61ee2SBarry Smith values[2] = (y1 - y0 + g_01 * y) * 0.5; 516c6f61ee2SBarry Smith values[3] = (y3 - y0 + g_01 * x) * 0.5; 517c6f61ee2SBarry Smith PetscCall(MatSetValues(J, 2, rows, 2, rows, values, INSERT_VALUES)); 518c6f61ee2SBarry Smith } 519c6f61ee2SBarry Smith PetscCall(PetscLogFlops(30)); 520c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(Xref, &ref)); 521c6f61ee2SBarry Smith PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 522c6f61ee2SBarry Smith PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 523c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 524c6f61ee2SBarry Smith } 525c6f61ee2SBarry Smith 526c6f61ee2SBarry Smith static inline PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, PetscInt p, Vec xLocal, Vec v) 527c6f61ee2SBarry Smith { 528c6f61ee2SBarry Smith const PetscInt c = ctx->cells[p]; 529c6f61ee2SBarry Smith PetscFE fem = NULL; 530c6f61ee2SBarry Smith DM dmCoord; 531c6f61ee2SBarry Smith SNES snes; 532c6f61ee2SBarry Smith KSP ksp; 533c6f61ee2SBarry Smith PC pc; 534c6f61ee2SBarry Smith Vec coordsLocal, r, ref, real; 535c6f61ee2SBarry Smith Mat J; 536c6f61ee2SBarry Smith PetscTabulation T = NULL; 537c6f61ee2SBarry Smith const PetscScalar *coords; 538c6f61ee2SBarry Smith PetscScalar *a; 539c6f61ee2SBarry Smith PetscReal xir[2] = {0., 0.}; 540c6f61ee2SBarry Smith PetscInt Nf; 541c6f61ee2SBarry Smith const PetscInt dof = ctx->dof; 542c6f61ee2SBarry Smith PetscScalar *x = NULL, *vertices = NULL; 543c6f61ee2SBarry Smith PetscScalar *xi; 544c6f61ee2SBarry Smith PetscInt coordSize, xSize; 545c6f61ee2SBarry Smith 546c6f61ee2SBarry Smith PetscFunctionBegin; 547c6f61ee2SBarry Smith PetscCall(DMGetNumFields(dm, &Nf)); 548c6f61ee2SBarry Smith if (Nf) { 549c6f61ee2SBarry Smith PetscObject obj; 550c6f61ee2SBarry Smith PetscClassId id; 551c6f61ee2SBarry Smith 552c6f61ee2SBarry Smith PetscCall(DMGetField(dm, 0, NULL, &obj)); 553c6f61ee2SBarry Smith PetscCall(PetscObjectGetClassId(obj, &id)); 554c6f61ee2SBarry Smith if (id == PETSCFE_CLASSID) { 555c6f61ee2SBarry Smith fem = (PetscFE)obj; 556c6f61ee2SBarry Smith PetscCall(PetscFECreateTabulation(fem, 1, 1, xir, 0, &T)); 557c6f61ee2SBarry Smith } 558c6f61ee2SBarry Smith } 559c6f61ee2SBarry Smith PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal)); 560c6f61ee2SBarry Smith PetscCall(DMGetCoordinateDM(dm, &dmCoord)); 561c6f61ee2SBarry Smith PetscCall(SNESCreate(PETSC_COMM_SELF, &snes)); 562c6f61ee2SBarry Smith PetscCall(SNESSetOptionsPrefix(snes, "quad_interp_")); 563c6f61ee2SBarry Smith PetscCall(VecCreate(PETSC_COMM_SELF, &r)); 564c6f61ee2SBarry Smith PetscCall(VecSetSizes(r, 2, 2)); 565c6f61ee2SBarry Smith PetscCall(VecSetType(r, dm->vectype)); 566c6f61ee2SBarry Smith PetscCall(VecDuplicate(r, &ref)); 567c6f61ee2SBarry Smith PetscCall(VecDuplicate(r, &real)); 568c6f61ee2SBarry Smith PetscCall(MatCreate(PETSC_COMM_SELF, &J)); 569c6f61ee2SBarry Smith PetscCall(MatSetSizes(J, 2, 2, 2, 2)); 570c6f61ee2SBarry Smith PetscCall(MatSetType(J, MATSEQDENSE)); 571c6f61ee2SBarry Smith PetscCall(MatSetUp(J)); 572c6f61ee2SBarry Smith PetscCall(SNESSetFunction(snes, r, QuadMap_Private, NULL)); 573c6f61ee2SBarry Smith PetscCall(SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL)); 574c6f61ee2SBarry Smith PetscCall(SNESGetKSP(snes, &ksp)); 575c6f61ee2SBarry Smith PetscCall(KSPGetPC(ksp, &pc)); 576c6f61ee2SBarry Smith PetscCall(PCSetType(pc, PCLU)); 577c6f61ee2SBarry Smith PetscCall(SNESSetFromOptions(snes)); 578c6f61ee2SBarry Smith 579c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(ctx->coords, &coords)); 580c6f61ee2SBarry Smith PetscCall(VecGetArray(v, &a)); 581c6f61ee2SBarry Smith PetscCall(DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices)); 582c6f61ee2SBarry Smith PetscCheck(4 * 2 == coordSize, ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %" PetscInt_FMT " should be %d", coordSize, 4 * 2); 583c6f61ee2SBarry Smith PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x)); 584c6f61ee2SBarry Smith PetscCall(SNESSetFunction(snes, NULL, NULL, vertices)); 585c6f61ee2SBarry Smith PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, vertices)); 586c6f61ee2SBarry Smith PetscCall(VecGetArray(real, &xi)); 587c6f61ee2SBarry Smith xi[0] = coords[p * ctx->dim + 0]; 588c6f61ee2SBarry Smith xi[1] = coords[p * ctx->dim + 1]; 589c6f61ee2SBarry Smith PetscCall(VecRestoreArray(real, &xi)); 590c6f61ee2SBarry Smith PetscCall(SNESSolve(snes, real, ref)); 591c6f61ee2SBarry Smith PetscCall(VecGetArray(ref, &xi)); 592c6f61ee2SBarry Smith xir[0] = PetscRealPart(xi[0]); 593c6f61ee2SBarry Smith xir[1] = PetscRealPart(xi[1]); 594c6f61ee2SBarry Smith if (4 * dof == xSize) { 595c6f61ee2SBarry Smith for (PetscInt comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp] * (1 - xir[0]) * (1 - xir[1]) + x[1 * dof + comp] * xir[0] * (1 - xir[1]) + x[2 * dof + comp] * xir[0] * xir[1] + x[3 * dof + comp] * (1 - xir[0]) * xir[1]; 596c6f61ee2SBarry Smith } else if (dof == xSize) { 597c6f61ee2SBarry Smith for (PetscInt comp = 0; comp < dof; ++comp) a[p * dof + comp] = x[0 * dof + comp]; 598c6f61ee2SBarry Smith } else { 599c6f61ee2SBarry Smith PetscCheck(fem, ctx->comm, PETSC_ERR_ARG_WRONG, "Cannot have a higher order interpolant if the discretization is not PetscFE"); 600c6f61ee2SBarry Smith xir[0] = 2.0 * xir[0] - 1.0; 601c6f61ee2SBarry Smith xir[1] = 2.0 * xir[1] - 1.0; 602c6f61ee2SBarry Smith PetscCall(PetscFEComputeTabulation(fem, 1, xir, 0, T)); 603c6f61ee2SBarry Smith for (PetscInt comp = 0; comp < dof; ++comp) { 604c6f61ee2SBarry Smith a[p * dof + comp] = 0.0; 605c6f61ee2SBarry Smith for (PetscInt d = 0; d < xSize / dof; ++d) a[p * dof + comp] += x[d * dof + comp] * T->T[0][d * dof + comp]; 606c6f61ee2SBarry Smith } 607c6f61ee2SBarry Smith } 608c6f61ee2SBarry Smith PetscCall(VecRestoreArray(ref, &xi)); 609c6f61ee2SBarry Smith PetscCall(DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices)); 610c6f61ee2SBarry Smith PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x)); 611c6f61ee2SBarry Smith PetscCall(PetscTabulationDestroy(&T)); 612c6f61ee2SBarry Smith PetscCall(VecRestoreArray(v, &a)); 613c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 614c6f61ee2SBarry Smith 615c6f61ee2SBarry Smith PetscCall(SNESDestroy(&snes)); 616c6f61ee2SBarry Smith PetscCall(VecDestroy(&r)); 617c6f61ee2SBarry Smith PetscCall(VecDestroy(&ref)); 618c6f61ee2SBarry Smith PetscCall(VecDestroy(&real)); 619c6f61ee2SBarry Smith PetscCall(MatDestroy(&J)); 620c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 621c6f61ee2SBarry Smith } 622c6f61ee2SBarry Smith 623c6f61ee2SBarry Smith static inline PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) 624c6f61ee2SBarry Smith { 625c6f61ee2SBarry Smith const PetscScalar *vertices = (const PetscScalar *)ctx; 626c6f61ee2SBarry Smith const PetscScalar x0 = vertices[0]; 627c6f61ee2SBarry Smith const PetscScalar y0 = vertices[1]; 628c6f61ee2SBarry Smith const PetscScalar z0 = vertices[2]; 629c6f61ee2SBarry Smith const PetscScalar x1 = vertices[9]; 630c6f61ee2SBarry Smith const PetscScalar y1 = vertices[10]; 631c6f61ee2SBarry Smith const PetscScalar z1 = vertices[11]; 632c6f61ee2SBarry Smith const PetscScalar x2 = vertices[6]; 633c6f61ee2SBarry Smith const PetscScalar y2 = vertices[7]; 634c6f61ee2SBarry Smith const PetscScalar z2 = vertices[8]; 635c6f61ee2SBarry Smith const PetscScalar x3 = vertices[3]; 636c6f61ee2SBarry Smith const PetscScalar y3 = vertices[4]; 637c6f61ee2SBarry Smith const PetscScalar z3 = vertices[5]; 638c6f61ee2SBarry Smith const PetscScalar x4 = vertices[12]; 639c6f61ee2SBarry Smith const PetscScalar y4 = vertices[13]; 640c6f61ee2SBarry Smith const PetscScalar z4 = vertices[14]; 641c6f61ee2SBarry Smith const PetscScalar x5 = vertices[15]; 642c6f61ee2SBarry Smith const PetscScalar y5 = vertices[16]; 643c6f61ee2SBarry Smith const PetscScalar z5 = vertices[17]; 644c6f61ee2SBarry Smith const PetscScalar x6 = vertices[18]; 645c6f61ee2SBarry Smith const PetscScalar y6 = vertices[19]; 646c6f61ee2SBarry Smith const PetscScalar z6 = vertices[20]; 647c6f61ee2SBarry Smith const PetscScalar x7 = vertices[21]; 648c6f61ee2SBarry Smith const PetscScalar y7 = vertices[22]; 649c6f61ee2SBarry Smith const PetscScalar z7 = vertices[23]; 650c6f61ee2SBarry Smith const PetscScalar f_1 = x1 - x0; 651c6f61ee2SBarry Smith const PetscScalar g_1 = y1 - y0; 652c6f61ee2SBarry Smith const PetscScalar h_1 = z1 - z0; 653c6f61ee2SBarry Smith const PetscScalar f_3 = x3 - x0; 654c6f61ee2SBarry Smith const PetscScalar g_3 = y3 - y0; 655c6f61ee2SBarry Smith const PetscScalar h_3 = z3 - z0; 656c6f61ee2SBarry Smith const PetscScalar f_4 = x4 - x0; 657c6f61ee2SBarry Smith const PetscScalar g_4 = y4 - y0; 658c6f61ee2SBarry Smith const PetscScalar h_4 = z4 - z0; 659c6f61ee2SBarry Smith const PetscScalar f_01 = x2 - x1 - x3 + x0; 660c6f61ee2SBarry Smith const PetscScalar g_01 = y2 - y1 - y3 + y0; 661c6f61ee2SBarry Smith const PetscScalar h_01 = z2 - z1 - z3 + z0; 662c6f61ee2SBarry Smith const PetscScalar f_12 = x7 - x3 - x4 + x0; 663c6f61ee2SBarry Smith const PetscScalar g_12 = y7 - y3 - y4 + y0; 664c6f61ee2SBarry Smith const PetscScalar h_12 = z7 - z3 - z4 + z0; 665c6f61ee2SBarry Smith const PetscScalar f_02 = x5 - x1 - x4 + x0; 666c6f61ee2SBarry Smith const PetscScalar g_02 = y5 - y1 - y4 + y0; 667c6f61ee2SBarry Smith const PetscScalar h_02 = z5 - z1 - z4 + z0; 668c6f61ee2SBarry Smith const PetscScalar f_012 = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 669c6f61ee2SBarry Smith const PetscScalar g_012 = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 670c6f61ee2SBarry Smith const PetscScalar h_012 = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 671c6f61ee2SBarry Smith const PetscScalar *ref; 672c6f61ee2SBarry Smith PetscScalar *real; 673c6f61ee2SBarry Smith 674c6f61ee2SBarry Smith PetscFunctionBegin; 675c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(Xref, &ref)); 676c6f61ee2SBarry Smith PetscCall(VecGetArray(Xreal, &real)); 677c6f61ee2SBarry Smith { 678c6f61ee2SBarry Smith const PetscScalar p0 = ref[0]; 679c6f61ee2SBarry Smith const PetscScalar p1 = ref[1]; 680c6f61ee2SBarry Smith const PetscScalar p2 = ref[2]; 681c6f61ee2SBarry Smith 682c6f61ee2SBarry Smith real[0] = x0 + f_1 * p0 + f_3 * p1 + f_4 * p2 + f_01 * p0 * p1 + f_12 * p1 * p2 + f_02 * p0 * p2 + f_012 * p0 * p1 * p2; 683c6f61ee2SBarry Smith real[1] = y0 + g_1 * p0 + g_3 * p1 + g_4 * p2 + g_01 * p0 * p1 + g_01 * p0 * p1 + g_12 * p1 * p2 + g_02 * p0 * p2 + g_012 * p0 * p1 * p2; 684c6f61ee2SBarry Smith real[2] = z0 + h_1 * p0 + h_3 * p1 + h_4 * p2 + h_01 * p0 * p1 + h_01 * p0 * p1 + h_12 * p1 * p2 + h_02 * p0 * p2 + h_012 * p0 * p1 * p2; 685c6f61ee2SBarry Smith } 686c6f61ee2SBarry Smith PetscCall(PetscLogFlops(114)); 687c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(Xref, &ref)); 688c6f61ee2SBarry Smith PetscCall(VecRestoreArray(Xreal, &real)); 689c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 690c6f61ee2SBarry Smith } 691c6f61ee2SBarry Smith 692c6f61ee2SBarry Smith static inline PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat J, Mat M, void *ctx) 693c6f61ee2SBarry Smith { 694c6f61ee2SBarry Smith const PetscScalar *vertices = (const PetscScalar *)ctx; 695c6f61ee2SBarry Smith const PetscScalar x0 = vertices[0]; 696c6f61ee2SBarry Smith const PetscScalar y0 = vertices[1]; 697c6f61ee2SBarry Smith const PetscScalar z0 = vertices[2]; 698c6f61ee2SBarry Smith const PetscScalar x1 = vertices[9]; 699c6f61ee2SBarry Smith const PetscScalar y1 = vertices[10]; 700c6f61ee2SBarry Smith const PetscScalar z1 = vertices[11]; 701c6f61ee2SBarry Smith const PetscScalar x2 = vertices[6]; 702c6f61ee2SBarry Smith const PetscScalar y2 = vertices[7]; 703c6f61ee2SBarry Smith const PetscScalar z2 = vertices[8]; 704c6f61ee2SBarry Smith const PetscScalar x3 = vertices[3]; 705c6f61ee2SBarry Smith const PetscScalar y3 = vertices[4]; 706c6f61ee2SBarry Smith const PetscScalar z3 = vertices[5]; 707c6f61ee2SBarry Smith const PetscScalar x4 = vertices[12]; 708c6f61ee2SBarry Smith const PetscScalar y4 = vertices[13]; 709c6f61ee2SBarry Smith const PetscScalar z4 = vertices[14]; 710c6f61ee2SBarry Smith const PetscScalar x5 = vertices[15]; 711c6f61ee2SBarry Smith const PetscScalar y5 = vertices[16]; 712c6f61ee2SBarry Smith const PetscScalar z5 = vertices[17]; 713c6f61ee2SBarry Smith const PetscScalar x6 = vertices[18]; 714c6f61ee2SBarry Smith const PetscScalar y6 = vertices[19]; 715c6f61ee2SBarry Smith const PetscScalar z6 = vertices[20]; 716c6f61ee2SBarry Smith const PetscScalar x7 = vertices[21]; 717c6f61ee2SBarry Smith const PetscScalar y7 = vertices[22]; 718c6f61ee2SBarry Smith const PetscScalar z7 = vertices[23]; 719c6f61ee2SBarry Smith const PetscScalar f_xy = x2 - x1 - x3 + x0; 720c6f61ee2SBarry Smith const PetscScalar g_xy = y2 - y1 - y3 + y0; 721c6f61ee2SBarry Smith const PetscScalar h_xy = z2 - z1 - z3 + z0; 722c6f61ee2SBarry Smith const PetscScalar f_yz = x7 - x3 - x4 + x0; 723c6f61ee2SBarry Smith const PetscScalar g_yz = y7 - y3 - y4 + y0; 724c6f61ee2SBarry Smith const PetscScalar h_yz = z7 - z3 - z4 + z0; 725c6f61ee2SBarry Smith const PetscScalar f_xz = x5 - x1 - x4 + x0; 726c6f61ee2SBarry Smith const PetscScalar g_xz = y5 - y1 - y4 + y0; 727c6f61ee2SBarry Smith const PetscScalar h_xz = z5 - z1 - z4 + z0; 728c6f61ee2SBarry Smith const PetscScalar f_xyz = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 729c6f61ee2SBarry Smith const PetscScalar g_xyz = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 730c6f61ee2SBarry Smith const PetscScalar h_xyz = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 731c6f61ee2SBarry Smith const PetscScalar *ref; 732c6f61ee2SBarry Smith 733c6f61ee2SBarry Smith PetscFunctionBegin; 734c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(Xref, &ref)); 735c6f61ee2SBarry Smith { 736c6f61ee2SBarry Smith const PetscScalar x = ref[0]; 737c6f61ee2SBarry Smith const PetscScalar y = ref[1]; 738c6f61ee2SBarry Smith const PetscScalar z = ref[2]; 739c6f61ee2SBarry Smith const PetscInt rows[3] = {0, 1, 2}; 740c6f61ee2SBarry Smith PetscScalar values[9]; 741c6f61ee2SBarry Smith 742c6f61ee2SBarry Smith values[0] = (x1 - x0 + f_xy * y + f_xz * z + f_xyz * y * z) / 2.0; 743c6f61ee2SBarry Smith values[1] = (x3 - x0 + f_xy * x + f_yz * z + f_xyz * x * z) / 2.0; 744c6f61ee2SBarry Smith values[2] = (x4 - x0 + f_yz * y + f_xz * x + f_xyz * x * y) / 2.0; 745c6f61ee2SBarry Smith values[3] = (y1 - y0 + g_xy * y + g_xz * z + g_xyz * y * z) / 2.0; 746c6f61ee2SBarry Smith values[4] = (y3 - y0 + g_xy * x + g_yz * z + g_xyz * x * z) / 2.0; 747c6f61ee2SBarry Smith values[5] = (y4 - y0 + g_yz * y + g_xz * x + g_xyz * x * y) / 2.0; 748c6f61ee2SBarry Smith values[6] = (z1 - z0 + h_xy * y + h_xz * z + h_xyz * y * z) / 2.0; 749c6f61ee2SBarry Smith values[7] = (z3 - z0 + h_xy * x + h_yz * z + h_xyz * x * z) / 2.0; 750c6f61ee2SBarry Smith values[8] = (z4 - z0 + h_yz * y + h_xz * x + h_xyz * x * y) / 2.0; 751c6f61ee2SBarry Smith 752c6f61ee2SBarry Smith PetscCall(MatSetValues(J, 3, rows, 3, rows, values, INSERT_VALUES)); 753c6f61ee2SBarry Smith } 754c6f61ee2SBarry Smith PetscCall(PetscLogFlops(152)); 755c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(Xref, &ref)); 756c6f61ee2SBarry Smith PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 757c6f61ee2SBarry Smith PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 758c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 759c6f61ee2SBarry Smith } 760c6f61ee2SBarry Smith 761c6f61ee2SBarry Smith static inline PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, PetscInt p, Vec xLocal, Vec v) 762c6f61ee2SBarry Smith { 763c6f61ee2SBarry Smith const PetscInt c = ctx->cells[p]; 764c6f61ee2SBarry Smith DM dmCoord; 765c6f61ee2SBarry Smith SNES snes; 766c6f61ee2SBarry Smith KSP ksp; 767c6f61ee2SBarry Smith PC pc; 768c6f61ee2SBarry Smith Vec coordsLocal, r, ref, real; 769c6f61ee2SBarry Smith Mat J; 770c6f61ee2SBarry Smith const PetscScalar *coords; 771c6f61ee2SBarry Smith PetscScalar *a; 772c6f61ee2SBarry Smith PetscScalar *x = NULL, *vertices = NULL; 773c6f61ee2SBarry Smith PetscScalar *xi; 774c6f61ee2SBarry Smith PetscReal xir[3]; 775c6f61ee2SBarry Smith PetscInt coordSize, xSize; 776c6f61ee2SBarry Smith 777c6f61ee2SBarry Smith PetscFunctionBegin; 778c6f61ee2SBarry Smith PetscCall(DMGetCoordinatesLocal(dm, &coordsLocal)); 779c6f61ee2SBarry Smith PetscCall(DMGetCoordinateDM(dm, &dmCoord)); 780c6f61ee2SBarry Smith PetscCall(SNESCreate(PETSC_COMM_SELF, &snes)); 781c6f61ee2SBarry Smith PetscCall(SNESSetOptionsPrefix(snes, "hex_interp_")); 782c6f61ee2SBarry Smith PetscCall(VecCreate(PETSC_COMM_SELF, &r)); 783c6f61ee2SBarry Smith PetscCall(VecSetSizes(r, 3, 3)); 784c6f61ee2SBarry Smith PetscCall(VecSetType(r, dm->vectype)); 785c6f61ee2SBarry Smith PetscCall(VecDuplicate(r, &ref)); 786c6f61ee2SBarry Smith PetscCall(VecDuplicate(r, &real)); 787c6f61ee2SBarry Smith PetscCall(MatCreate(PETSC_COMM_SELF, &J)); 788c6f61ee2SBarry Smith PetscCall(MatSetSizes(J, 3, 3, 3, 3)); 789c6f61ee2SBarry Smith PetscCall(MatSetType(J, MATSEQDENSE)); 790c6f61ee2SBarry Smith PetscCall(MatSetUp(J)); 791c6f61ee2SBarry Smith PetscCall(SNESSetFunction(snes, r, HexMap_Private, NULL)); 792c6f61ee2SBarry Smith PetscCall(SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL)); 793c6f61ee2SBarry Smith PetscCall(SNESGetKSP(snes, &ksp)); 794c6f61ee2SBarry Smith PetscCall(KSPGetPC(ksp, &pc)); 795c6f61ee2SBarry Smith PetscCall(PCSetType(pc, PCLU)); 796c6f61ee2SBarry Smith PetscCall(SNESSetFromOptions(snes)); 797c6f61ee2SBarry Smith 798c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(ctx->coords, &coords)); 799c6f61ee2SBarry Smith PetscCall(VecGetArray(v, &a)); 800c6f61ee2SBarry Smith PetscCall(DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices)); 801c6f61ee2SBarry Smith PetscCheck(8 * 3 == coordSize, ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid coordinate closure size %" PetscInt_FMT " should be %d", coordSize, 8 * 3); 802c6f61ee2SBarry Smith PetscCall(DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x)); 803c6f61ee2SBarry Smith PetscCheck((8 * ctx->dof == xSize) || (ctx->dof == xSize), ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid input closure size %" PetscInt_FMT " should be %" PetscInt_FMT " or %" PetscInt_FMT, xSize, 8 * ctx->dof, ctx->dof); 804c6f61ee2SBarry Smith PetscCall(SNESSetFunction(snes, NULL, NULL, vertices)); 805c6f61ee2SBarry Smith PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, vertices)); 806c6f61ee2SBarry Smith PetscCall(VecGetArray(real, &xi)); 807c6f61ee2SBarry Smith xi[0] = coords[p * ctx->dim + 0]; 808c6f61ee2SBarry Smith xi[1] = coords[p * ctx->dim + 1]; 809c6f61ee2SBarry Smith xi[2] = coords[p * ctx->dim + 2]; 810c6f61ee2SBarry Smith PetscCall(VecRestoreArray(real, &xi)); 811c6f61ee2SBarry Smith PetscCall(SNESSolve(snes, real, ref)); 812c6f61ee2SBarry Smith PetscCall(VecGetArray(ref, &xi)); 813c6f61ee2SBarry Smith xir[0] = PetscRealPart(xi[0]); 814c6f61ee2SBarry Smith xir[1] = PetscRealPart(xi[1]); 815c6f61ee2SBarry Smith xir[2] = PetscRealPart(xi[2]); 816c6f61ee2SBarry Smith if (8 * ctx->dof == xSize) { 817c6f61ee2SBarry Smith for (PetscInt comp = 0; comp < ctx->dof; ++comp) { 818c6f61ee2SBarry Smith a[p * ctx->dof + comp] = x[0 * ctx->dof + comp] * (1 - xir[0]) * (1 - xir[1]) * (1 - xir[2]) + x[3 * ctx->dof + comp] * xir[0] * (1 - xir[1]) * (1 - xir[2]) + x[2 * ctx->dof + comp] * xir[0] * xir[1] * (1 - xir[2]) + x[1 * ctx->dof + comp] * (1 - xir[0]) * xir[1] * (1 - xir[2]) + 819c6f61ee2SBarry Smith x[4 * ctx->dof + comp] * (1 - xir[0]) * (1 - xir[1]) * xir[2] + x[5 * ctx->dof + comp] * xir[0] * (1 - xir[1]) * xir[2] + x[6 * ctx->dof + comp] * xir[0] * xir[1] * xir[2] + x[7 * ctx->dof + comp] * (1 - xir[0]) * xir[1] * xir[2]; 820c6f61ee2SBarry Smith } 821c6f61ee2SBarry Smith } else { 822c6f61ee2SBarry Smith for (PetscInt comp = 0; comp < ctx->dof; ++comp) a[p * ctx->dof + comp] = x[0 * ctx->dof + comp]; 823c6f61ee2SBarry Smith } 824c6f61ee2SBarry Smith PetscCall(VecRestoreArray(ref, &xi)); 825c6f61ee2SBarry Smith PetscCall(DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices)); 826c6f61ee2SBarry Smith PetscCall(DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x)); 827c6f61ee2SBarry Smith PetscCall(VecRestoreArray(v, &a)); 828c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 829c6f61ee2SBarry Smith 830c6f61ee2SBarry Smith PetscCall(SNESDestroy(&snes)); 831c6f61ee2SBarry Smith PetscCall(VecDestroy(&r)); 832c6f61ee2SBarry Smith PetscCall(VecDestroy(&ref)); 833c6f61ee2SBarry Smith PetscCall(VecDestroy(&real)); 834c6f61ee2SBarry Smith PetscCall(MatDestroy(&J)); 835c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 836c6f61ee2SBarry Smith } 837c6f61ee2SBarry Smith 838c6f61ee2SBarry Smith /*@C 839df514481SBarry Smith DMInterpolationEvaluate - Using the input from `dm` and `x`, calculates interpolated field values at the interpolation points. 840c6f61ee2SBarry Smith 841c6f61ee2SBarry Smith Input Parameters: 842df514481SBarry Smith + ctx - The `DMInterpolationInfo` context obtained with `DMInterpolationCreate()` 843c6f61ee2SBarry Smith . dm - The `DM` 844d447a1d4SMatthew Knepley - x - The local vector containing the field to be interpolated, can be created with `DMCreateGlobalVector()` 845c6f61ee2SBarry Smith 846c6f61ee2SBarry Smith Output Parameter: 847df514481SBarry Smith . v - The vector containing the interpolated values, obtained with `DMInterpolationGetVector()` 848c6f61ee2SBarry Smith 849c6f61ee2SBarry Smith Level: beginner 850c6f61ee2SBarry Smith 851df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMInterpolationInfo`, `DMInterpolationGetVector()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()`, `DMInterpolationGetCoordinates()` 852c6f61ee2SBarry Smith @*/ 853c6f61ee2SBarry Smith PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v) 854c6f61ee2SBarry Smith { 855c6f61ee2SBarry Smith PetscDS ds; 856c6f61ee2SBarry Smith PetscInt n, p, Nf, field; 857c6f61ee2SBarry Smith PetscBool useDS = PETSC_FALSE; 858c6f61ee2SBarry Smith 859c6f61ee2SBarry Smith PetscFunctionBegin; 860c6f61ee2SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 861c6f61ee2SBarry Smith PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 862c6f61ee2SBarry Smith PetscValidHeaderSpecific(v, VEC_CLASSID, 4); 863c6f61ee2SBarry Smith PetscCall(VecGetLocalSize(v, &n)); 864c6f61ee2SBarry Smith PetscCheck(n == ctx->n * ctx->dof, ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid input vector size %" PetscInt_FMT " should be %" PetscInt_FMT, n, ctx->n * ctx->dof); 865c6f61ee2SBarry Smith if (!n) PetscFunctionReturn(PETSC_SUCCESS); 866c6f61ee2SBarry Smith PetscCall(DMGetDS(dm, &ds)); 867c6f61ee2SBarry Smith if (ds) { 868c6f61ee2SBarry Smith useDS = PETSC_TRUE; 869c6f61ee2SBarry Smith PetscCall(PetscDSGetNumFields(ds, &Nf)); 870c6f61ee2SBarry Smith for (field = 0; field < Nf; ++field) { 871c6f61ee2SBarry Smith PetscObject obj; 872c6f61ee2SBarry Smith PetscClassId id; 873c6f61ee2SBarry Smith 874c6f61ee2SBarry Smith PetscCall(PetscDSGetDiscretization(ds, field, &obj)); 875c6f61ee2SBarry Smith PetscCall(PetscObjectGetClassId(obj, &id)); 876c6f61ee2SBarry Smith if (id != PETSCFE_CLASSID && id != PETSCFV_CLASSID) { 877c6f61ee2SBarry Smith useDS = PETSC_FALSE; 878c6f61ee2SBarry Smith break; 879c6f61ee2SBarry Smith } 880c6f61ee2SBarry Smith } 881c6f61ee2SBarry Smith } 882c6f61ee2SBarry Smith if (useDS) { 883c6f61ee2SBarry Smith const PetscScalar *coords; 884c6f61ee2SBarry Smith PetscScalar *interpolant; 885c6f61ee2SBarry Smith PetscInt cdim, d; 886c6f61ee2SBarry Smith 887c6f61ee2SBarry Smith PetscCall(DMGetCoordinateDim(dm, &cdim)); 888c6f61ee2SBarry Smith PetscCall(VecGetArrayRead(ctx->coords, &coords)); 889c6f61ee2SBarry Smith PetscCall(VecGetArrayWrite(v, &interpolant)); 890c6f61ee2SBarry Smith for (p = 0; p < ctx->n; ++p) { 891c6f61ee2SBarry Smith PetscReal pcoords[3], xi[3]; 892c6f61ee2SBarry Smith PetscScalar *xa = NULL; 893c6f61ee2SBarry Smith PetscInt coff = 0, foff = 0, clSize; 894c6f61ee2SBarry Smith 895c6f61ee2SBarry Smith if (ctx->cells[p] < 0) continue; 896c6f61ee2SBarry Smith for (d = 0; d < cdim; ++d) pcoords[d] = PetscRealPart(coords[p * cdim + d]); 897c6f61ee2SBarry Smith PetscCall(DMPlexCoordinatesToReference(dm, ctx->cells[p], 1, pcoords, xi)); 898c6f61ee2SBarry Smith PetscCall(DMPlexVecGetClosure(dm, NULL, x, ctx->cells[p], &clSize, &xa)); 899c6f61ee2SBarry Smith for (field = 0; field < Nf; ++field) { 900c6f61ee2SBarry Smith PetscTabulation T; 901c6f61ee2SBarry Smith PetscObject obj; 902c6f61ee2SBarry Smith PetscClassId id; 903c6f61ee2SBarry Smith 904c6f61ee2SBarry Smith PetscCall(PetscDSGetDiscretization(ds, field, &obj)); 905c6f61ee2SBarry Smith PetscCall(PetscObjectGetClassId(obj, &id)); 906c6f61ee2SBarry Smith if (id == PETSCFE_CLASSID) { 907c6f61ee2SBarry Smith PetscFE fe = (PetscFE)obj; 908c6f61ee2SBarry Smith 909c6f61ee2SBarry Smith PetscCall(PetscFECreateTabulation(fe, 1, 1, xi, 0, &T)); 910c6f61ee2SBarry Smith { 911c6f61ee2SBarry Smith const PetscReal *basis = T->T[0]; 912c6f61ee2SBarry Smith const PetscInt Nb = T->Nb; 913c6f61ee2SBarry Smith const PetscInt Nc = T->Nc; 914c6f61ee2SBarry Smith 915c6f61ee2SBarry Smith for (PetscInt fc = 0; fc < Nc; ++fc) { 916c6f61ee2SBarry Smith interpolant[p * ctx->dof + coff + fc] = 0.0; 917c6f61ee2SBarry Smith for (PetscInt f = 0; f < Nb; ++f) interpolant[p * ctx->dof + coff + fc] += xa[foff + f] * basis[(0 * Nb + f) * Nc + fc]; 918c6f61ee2SBarry Smith } 919c6f61ee2SBarry Smith coff += Nc; 920c6f61ee2SBarry Smith foff += Nb; 921c6f61ee2SBarry Smith } 922c6f61ee2SBarry Smith PetscCall(PetscTabulationDestroy(&T)); 923c6f61ee2SBarry Smith } else if (id == PETSCFV_CLASSID) { 924c6f61ee2SBarry Smith PetscFV fv = (PetscFV)obj; 925c6f61ee2SBarry Smith PetscInt Nc; 926c6f61ee2SBarry Smith 927c6f61ee2SBarry Smith // TODO Could use reconstruction if available 928c6f61ee2SBarry Smith PetscCall(PetscFVGetNumComponents(fv, &Nc)); 929c6f61ee2SBarry Smith for (PetscInt fc = 0; fc < Nc; ++fc) interpolant[p * ctx->dof + coff + fc] = xa[foff + fc]; 930c6f61ee2SBarry Smith coff += Nc; 931c6f61ee2SBarry Smith foff += Nc; 932c6f61ee2SBarry Smith } 933c6f61ee2SBarry Smith } 934c6f61ee2SBarry Smith PetscCall(DMPlexVecRestoreClosure(dm, NULL, x, ctx->cells[p], &clSize, &xa)); 935c6f61ee2SBarry Smith PetscCheck(coff == ctx->dof, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total components %" PetscInt_FMT " != %" PetscInt_FMT " dof specified for interpolation", coff, ctx->dof); 936c6f61ee2SBarry Smith PetscCheck(foff == clSize, PETSC_COMM_SELF, PETSC_ERR_PLIB, "Total FE/FV space size %" PetscInt_FMT " != %" PetscInt_FMT " closure size", foff, clSize); 937c6f61ee2SBarry Smith } 938c6f61ee2SBarry Smith PetscCall(VecRestoreArrayRead(ctx->coords, &coords)); 939c6f61ee2SBarry Smith PetscCall(VecRestoreArrayWrite(v, &interpolant)); 940c6f61ee2SBarry Smith } else { 941c6f61ee2SBarry Smith for (PetscInt p = 0; p < ctx->n; ++p) { 942c6f61ee2SBarry Smith const PetscInt cell = ctx->cells[p]; 943c6f61ee2SBarry Smith DMPolytopeType ct; 944c6f61ee2SBarry Smith 945c6f61ee2SBarry Smith PetscCall(DMPlexGetCellType(dm, cell, &ct)); 946c6f61ee2SBarry Smith switch (ct) { 947c6f61ee2SBarry Smith case DM_POLYTOPE_SEGMENT: 948c6f61ee2SBarry Smith PetscCall(DMInterpolate_Segment_Private(ctx, dm, p, x, v)); 949c6f61ee2SBarry Smith break; 950c6f61ee2SBarry Smith case DM_POLYTOPE_TRIANGLE: 951c6f61ee2SBarry Smith PetscCall(DMInterpolate_Triangle_Private(ctx, dm, p, x, v)); 952c6f61ee2SBarry Smith break; 953c6f61ee2SBarry Smith case DM_POLYTOPE_QUADRILATERAL: 954c6f61ee2SBarry Smith PetscCall(DMInterpolate_Quad_Private(ctx, dm, p, x, v)); 955c6f61ee2SBarry Smith break; 956c6f61ee2SBarry Smith case DM_POLYTOPE_TETRAHEDRON: 957c6f61ee2SBarry Smith PetscCall(DMInterpolate_Tetrahedron_Private(ctx, dm, p, x, v)); 958c6f61ee2SBarry Smith break; 959c6f61ee2SBarry Smith case DM_POLYTOPE_HEXAHEDRON: 960c6f61ee2SBarry Smith PetscCall(DMInterpolate_Hex_Private(ctx, dm, cell, x, v)); 961c6f61ee2SBarry Smith break; 962c6f61ee2SBarry Smith default: 963c6f61ee2SBarry Smith SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "No support for cell type %s", DMPolytopeTypes[PetscMax(0, PetscMin(ct, DM_NUM_POLYTOPES))]); 964c6f61ee2SBarry Smith } 965c6f61ee2SBarry Smith } 966c6f61ee2SBarry Smith } 967c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 968c6f61ee2SBarry Smith } 969c6f61ee2SBarry Smith 970c6f61ee2SBarry Smith /*@C 971c6f61ee2SBarry Smith DMInterpolationDestroy - Destroys a `DMInterpolationInfo` context 972c6f61ee2SBarry Smith 973c6f61ee2SBarry Smith Collective 974c6f61ee2SBarry Smith 975c6f61ee2SBarry Smith Input Parameter: 976c6f61ee2SBarry Smith . ctx - the context 977c6f61ee2SBarry Smith 978c6f61ee2SBarry Smith Level: beginner 979c6f61ee2SBarry Smith 980df514481SBarry Smith .seealso: [](ch_dmbase), `DM`, `DMInterpolationInfo`, `DMInterpolationEvaluate()`, `DMInterpolationAddPoints()`, `DMInterpolationCreate()` 981c6f61ee2SBarry Smith @*/ 982c6f61ee2SBarry Smith PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx) 983c6f61ee2SBarry Smith { 984c6f61ee2SBarry Smith PetscFunctionBegin; 985c6f61ee2SBarry Smith PetscAssertPointer(ctx, 1); 986c6f61ee2SBarry Smith PetscCall(VecDestroy(&(*ctx)->coords)); 987c6f61ee2SBarry Smith PetscCall(PetscFree((*ctx)->points)); 988c6f61ee2SBarry Smith PetscCall(PetscFree((*ctx)->cells)); 989c6f61ee2SBarry Smith PetscCall(PetscFree(*ctx)); 990c6f61ee2SBarry Smith *ctx = NULL; 991c6f61ee2SBarry Smith PetscFunctionReturn(PETSC_SUCCESS); 992c6f61ee2SBarry Smith } 993