1552f7358SJed Brown #include <petscdmplex.h> /*I "petscdmplex.h" I*/ 2552f7358SJed Brown #include <petscsnes.h> /*I "petscsnes.h" I*/ 3afcb2eb5SJed Brown #include <petsc-private/petscimpl.h> 4552f7358SJed Brown 5552f7358SJed Brown #undef __FUNCT__ 6552f7358SJed Brown #define __FUNCT__ "DMInterpolationCreate" 70adebc6cSBarry Smith PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx) 80adebc6cSBarry Smith { 9552f7358SJed Brown PetscErrorCode ierr; 10552f7358SJed Brown 11552f7358SJed Brown PetscFunctionBegin; 12552f7358SJed Brown PetscValidPointer(ctx, 2); 13552f7358SJed Brown ierr = PetscMalloc(sizeof(struct _DMInterpolationInfo), ctx);CHKERRQ(ierr); 141aa26658SKarl Rupp 15552f7358SJed Brown (*ctx)->comm = comm; 16552f7358SJed Brown (*ctx)->dim = -1; 17552f7358SJed Brown (*ctx)->nInput = 0; 180298fd71SBarry Smith (*ctx)->points = NULL; 190298fd71SBarry Smith (*ctx)->cells = NULL; 20552f7358SJed Brown (*ctx)->n = -1; 210298fd71SBarry Smith (*ctx)->coords = NULL; 22552f7358SJed Brown PetscFunctionReturn(0); 23552f7358SJed Brown } 24552f7358SJed Brown 25552f7358SJed Brown #undef __FUNCT__ 26552f7358SJed Brown #define __FUNCT__ "DMInterpolationSetDim" 270adebc6cSBarry Smith PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim) 280adebc6cSBarry Smith { 29552f7358SJed Brown PetscFunctionBegin; 300adebc6cSBarry Smith if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %d", dim); 31552f7358SJed Brown ctx->dim = dim; 32552f7358SJed Brown PetscFunctionReturn(0); 33552f7358SJed Brown } 34552f7358SJed Brown 35552f7358SJed Brown #undef __FUNCT__ 36552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetDim" 370adebc6cSBarry Smith PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim) 380adebc6cSBarry Smith { 39552f7358SJed Brown PetscFunctionBegin; 40552f7358SJed Brown PetscValidIntPointer(dim, 2); 41552f7358SJed Brown *dim = ctx->dim; 42552f7358SJed Brown PetscFunctionReturn(0); 43552f7358SJed Brown } 44552f7358SJed Brown 45552f7358SJed Brown #undef __FUNCT__ 46552f7358SJed Brown #define __FUNCT__ "DMInterpolationSetDof" 470adebc6cSBarry Smith PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof) 480adebc6cSBarry Smith { 49552f7358SJed Brown PetscFunctionBegin; 500adebc6cSBarry Smith if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %d", dof); 51552f7358SJed Brown ctx->dof = dof; 52552f7358SJed Brown PetscFunctionReturn(0); 53552f7358SJed Brown } 54552f7358SJed Brown 55552f7358SJed Brown #undef __FUNCT__ 56552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetDof" 570adebc6cSBarry Smith PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof) 580adebc6cSBarry Smith { 59552f7358SJed Brown PetscFunctionBegin; 60552f7358SJed Brown PetscValidIntPointer(dof, 2); 61552f7358SJed Brown *dof = ctx->dof; 62552f7358SJed Brown PetscFunctionReturn(0); 63552f7358SJed Brown } 64552f7358SJed Brown 65552f7358SJed Brown #undef __FUNCT__ 66552f7358SJed Brown #define __FUNCT__ "DMInterpolationAddPoints" 670adebc6cSBarry Smith PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[]) 680adebc6cSBarry Smith { 69552f7358SJed Brown PetscErrorCode ierr; 70552f7358SJed Brown 71552f7358SJed Brown PetscFunctionBegin; 720adebc6cSBarry Smith if (ctx->dim < 0) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 730adebc6cSBarry Smith if (ctx->points) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet"); 74552f7358SJed Brown ctx->nInput = n; 751aa26658SKarl Rupp 76552f7358SJed Brown ierr = PetscMalloc(n*ctx->dim * sizeof(PetscReal), &ctx->points);CHKERRQ(ierr); 77552f7358SJed Brown ierr = PetscMemcpy(ctx->points, points, n*ctx->dim * sizeof(PetscReal));CHKERRQ(ierr); 78552f7358SJed Brown PetscFunctionReturn(0); 79552f7358SJed Brown } 80552f7358SJed Brown 81552f7358SJed Brown #undef __FUNCT__ 82552f7358SJed Brown #define __FUNCT__ "DMInterpolationSetUp" 830adebc6cSBarry Smith PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints) 840adebc6cSBarry Smith { 85552f7358SJed Brown MPI_Comm comm = ctx->comm; 86552f7358SJed Brown PetscScalar *a; 87552f7358SJed Brown PetscInt p, q, i; 88552f7358SJed Brown PetscMPIInt rank, size; 89552f7358SJed Brown PetscErrorCode ierr; 90552f7358SJed Brown Vec pointVec; 91552f7358SJed Brown IS cellIS; 92552f7358SJed Brown PetscLayout layout; 93552f7358SJed Brown PetscReal *globalPoints; 94cb313848SJed Brown PetscScalar *globalPointsScalar; 95552f7358SJed Brown const PetscInt *ranges; 96552f7358SJed Brown PetscMPIInt *counts, *displs; 97552f7358SJed Brown const PetscInt *foundCells; 98552f7358SJed Brown PetscMPIInt *foundProcs, *globalProcs; 9919436ca2SJed Brown PetscInt n, N; 100552f7358SJed Brown 10119436ca2SJed Brown PetscFunctionBegin; 10219436ca2SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 10319436ca2SJed Brown ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 10419436ca2SJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1050adebc6cSBarry Smith if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 10619436ca2SJed Brown /* Locate points */ 10719436ca2SJed Brown n = ctx->nInput; 108552f7358SJed Brown if (!redundantPoints) { 109552f7358SJed Brown ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr); 110552f7358SJed Brown ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr); 111552f7358SJed Brown ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr); 112552f7358SJed Brown ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr); 113552f7358SJed Brown ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr); 114552f7358SJed Brown /* Communicate all points to all processes */ 115*dcca6d9dSJed Brown ierr = PetscMalloc3(N*ctx->dim,&globalPoints,size,&counts,size,&displs);CHKERRQ(ierr); 116552f7358SJed Brown ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr); 117552f7358SJed Brown for (p = 0; p < size; ++p) { 118552f7358SJed Brown counts[p] = (ranges[p+1] - ranges[p])*ctx->dim; 119552f7358SJed Brown displs[p] = ranges[p]*ctx->dim; 120552f7358SJed Brown } 121552f7358SJed Brown ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRQ(ierr); 122552f7358SJed Brown } else { 123552f7358SJed Brown N = n; 124f5af7f23SKarl Rupp 125552f7358SJed Brown globalPoints = ctx->points; 126552f7358SJed Brown } 127552f7358SJed Brown #if 0 128*dcca6d9dSJed Brown ierr = PetscMalloc3(N,&foundCells,N,&foundProcs,N,&globalProcs);CHKERRQ(ierr); 12919436ca2SJed Brown /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */ 130552f7358SJed Brown #else 131cb313848SJed Brown #if defined(PETSC_USE_COMPLEX) 132cb313848SJed Brown ierr = PetscMalloc(N*sizeof(PetscScalar),&globalPointsScalar);CHKERRQ(ierr); 133cb313848SJed Brown for (i=0; i<N; i++) globalPointsScalar[i] = globalPoints[i]; 134cb313848SJed Brown #else 135cb313848SJed Brown globalPointsScalar = globalPoints; 136cb313848SJed Brown #endif 13704706141SMatthew G Knepley ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);CHKERRQ(ierr); 138*dcca6d9dSJed Brown ierr = PetscMalloc2(N,&foundProcs,N,&globalProcs);CHKERRQ(ierr); 139552f7358SJed Brown ierr = DMLocatePoints(dm, pointVec, &cellIS);CHKERRQ(ierr); 140552f7358SJed Brown ierr = ISGetIndices(cellIS, &foundCells);CHKERRQ(ierr); 141552f7358SJed Brown #endif 142552f7358SJed Brown for (p = 0; p < N; ++p) { 1431aa26658SKarl Rupp if (foundCells[p] >= 0) foundProcs[p] = rank; 1441aa26658SKarl Rupp else foundProcs[p] = size; 145552f7358SJed Brown } 146552f7358SJed Brown /* Let the lowest rank process own each point */ 147efab3cc2SJed Brown ierr = MPI_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);CHKERRQ(ierr); 148552f7358SJed Brown ctx->n = 0; 149552f7358SJed Brown for (p = 0; p < N; ++p) { 1500adebc6cSBarry Smith if (globalProcs[p] == size) SETERRQ4(comm, PETSC_ERR_PLIB, "Point %d: %g %g %g not located in mesh", p, globalPoints[p*ctx->dim+0], ctx->dim > 1 ? globalPoints[p*ctx->dim+1] : 0.0, ctx->dim > 2 ? globalPoints[p*ctx->dim+2] : 0.0); 1511aa26658SKarl Rupp else if (globalProcs[p] == rank) ctx->n++; 152552f7358SJed Brown } 153552f7358SJed Brown /* Create coordinates vector and array of owned cells */ 154552f7358SJed Brown ierr = PetscMalloc(ctx->n * sizeof(PetscInt), &ctx->cells);CHKERRQ(ierr); 155552f7358SJed Brown ierr = VecCreate(comm, &ctx->coords);CHKERRQ(ierr); 156552f7358SJed Brown ierr = VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);CHKERRQ(ierr); 157552f7358SJed Brown ierr = VecSetBlockSize(ctx->coords, ctx->dim);CHKERRQ(ierr); 158c0dedaeaSBarry Smith ierr = VecSetType(ctx->coords,VECSTANDARD);CHKERRQ(ierr); 159552f7358SJed Brown ierr = VecGetArray(ctx->coords, &a);CHKERRQ(ierr); 160552f7358SJed Brown for (p = 0, q = 0, i = 0; p < N; ++p) { 161552f7358SJed Brown if (globalProcs[p] == rank) { 162552f7358SJed Brown PetscInt d; 163552f7358SJed Brown 1641aa26658SKarl Rupp for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d]; 165552f7358SJed Brown ctx->cells[q++] = foundCells[p]; 166552f7358SJed Brown } 167552f7358SJed Brown } 168552f7358SJed Brown ierr = VecRestoreArray(ctx->coords, &a);CHKERRQ(ierr); 169552f7358SJed Brown #if 0 170552f7358SJed Brown ierr = PetscFree3(foundCells,foundProcs,globalProcs);CHKERRQ(ierr); 171552f7358SJed Brown #else 172552f7358SJed Brown ierr = PetscFree2(foundProcs,globalProcs);CHKERRQ(ierr); 173552f7358SJed Brown ierr = ISRestoreIndices(cellIS, &foundCells);CHKERRQ(ierr); 174552f7358SJed Brown ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 175552f7358SJed Brown ierr = VecDestroy(&pointVec);CHKERRQ(ierr); 176552f7358SJed Brown #endif 177cb313848SJed Brown if ((void*)globalPointsScalar != (void*)globalPoints) {ierr = PetscFree(globalPointsScalar);CHKERRQ(ierr);} 178552f7358SJed Brown if (!redundantPoints) { 179552f7358SJed Brown ierr = PetscFree3(globalPoints,counts,displs);CHKERRQ(ierr); 180552f7358SJed Brown ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr); 181552f7358SJed Brown } 182552f7358SJed Brown PetscFunctionReturn(0); 183552f7358SJed Brown } 184552f7358SJed Brown 185552f7358SJed Brown #undef __FUNCT__ 186552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetCoordinates" 1870adebc6cSBarry Smith PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates) 1880adebc6cSBarry Smith { 189552f7358SJed Brown PetscFunctionBegin; 190552f7358SJed Brown PetscValidPointer(coordinates, 2); 1910adebc6cSBarry Smith if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 192552f7358SJed Brown *coordinates = ctx->coords; 193552f7358SJed Brown PetscFunctionReturn(0); 194552f7358SJed Brown } 195552f7358SJed Brown 196552f7358SJed Brown #undef __FUNCT__ 197552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetVector" 1980adebc6cSBarry Smith PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v) 1990adebc6cSBarry Smith { 200552f7358SJed Brown PetscErrorCode ierr; 201552f7358SJed Brown 202552f7358SJed Brown PetscFunctionBegin; 203552f7358SJed Brown PetscValidPointer(v, 2); 2040adebc6cSBarry Smith if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 205552f7358SJed Brown ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr); 206552f7358SJed Brown ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr); 207552f7358SJed Brown ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr); 208c0dedaeaSBarry Smith ierr = VecSetType(*v,VECSTANDARD);CHKERRQ(ierr); 209552f7358SJed Brown PetscFunctionReturn(0); 210552f7358SJed Brown } 211552f7358SJed Brown 212552f7358SJed Brown #undef __FUNCT__ 213552f7358SJed Brown #define __FUNCT__ "DMInterpolationRestoreVector" 2140adebc6cSBarry Smith PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v) 2150adebc6cSBarry Smith { 216552f7358SJed Brown PetscErrorCode ierr; 217552f7358SJed Brown 218552f7358SJed Brown PetscFunctionBegin; 219552f7358SJed Brown PetscValidPointer(v, 2); 2200adebc6cSBarry Smith if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 221552f7358SJed Brown ierr = VecDestroy(v);CHKERRQ(ierr); 222552f7358SJed Brown PetscFunctionReturn(0); 223552f7358SJed Brown } 224552f7358SJed Brown 225552f7358SJed Brown #undef __FUNCT__ 2267a1931ceSMatthew G. Knepley #define __FUNCT__ "DMInterpolate_Triangle_Private" 2277a1931ceSMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Triangle_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 228a6dfd86eSKarl Rupp { 229552f7358SJed Brown PetscReal *v0, *J, *invJ, detJ; 230552f7358SJed Brown PetscScalar *a, *coords; 231552f7358SJed Brown PetscInt p; 232552f7358SJed Brown PetscErrorCode ierr; 233552f7358SJed Brown 234552f7358SJed Brown PetscFunctionBegin; 235*dcca6d9dSJed Brown ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr); 236552f7358SJed Brown ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr); 237552f7358SJed Brown ierr = VecGetArray(v, &a);CHKERRQ(ierr); 238552f7358SJed Brown for (p = 0; p < ctx->n; ++p) { 239552f7358SJed Brown PetscInt c = ctx->cells[p]; 240a1e44745SMatthew G. Knepley PetscScalar *x = NULL; 241552f7358SJed Brown PetscReal xi[4]; 242552f7358SJed Brown PetscInt d, f, comp; 243552f7358SJed Brown 244552f7358SJed Brown ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 245552f7358SJed Brown if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 2460298fd71SBarry Smith ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr); 2471aa26658SKarl Rupp for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp]; 2481aa26658SKarl Rupp 249552f7358SJed Brown for (d = 0; d < ctx->dim; ++d) { 250552f7358SJed Brown xi[d] = 0.0; 2511aa26658SKarl Rupp for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]); 2521aa26658SKarl Rupp for (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]; 253552f7358SJed Brown } 2540298fd71SBarry Smith ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr); 255552f7358SJed Brown } 256552f7358SJed Brown ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 257552f7358SJed Brown ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr); 258552f7358SJed Brown ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr); 259552f7358SJed Brown PetscFunctionReturn(0); 260552f7358SJed Brown } 261552f7358SJed Brown 262552f7358SJed Brown #undef __FUNCT__ 2637a1931ceSMatthew G. Knepley #define __FUNCT__ "DMInterpolate_Tetrahedron_Private" 2647a1931ceSMatthew G. Knepley PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Tetrahedron_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 2657a1931ceSMatthew G. Knepley { 2667a1931ceSMatthew G. Knepley PetscReal *v0, *J, *invJ, detJ; 2677a1931ceSMatthew G. Knepley PetscScalar *a, *coords; 2687a1931ceSMatthew G. Knepley PetscInt p; 2697a1931ceSMatthew G. Knepley PetscErrorCode ierr; 2707a1931ceSMatthew G. Knepley 2717a1931ceSMatthew G. Knepley PetscFunctionBegin; 272*dcca6d9dSJed Brown ierr = PetscMalloc3(ctx->dim,&v0,ctx->dim*ctx->dim,&J,ctx->dim*ctx->dim,&invJ);CHKERRQ(ierr); 2737a1931ceSMatthew G. Knepley ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr); 2747a1931ceSMatthew G. Knepley ierr = VecGetArray(v, &a);CHKERRQ(ierr); 2757a1931ceSMatthew G. Knepley for (p = 0; p < ctx->n; ++p) { 2767a1931ceSMatthew G. Knepley PetscInt c = ctx->cells[p]; 2777a1931ceSMatthew G. Knepley const PetscInt order[3] = {2, 1, 3}; 2782584bbe8SMatthew G. Knepley PetscScalar *x = NULL; 2797a1931ceSMatthew G. Knepley PetscReal xi[4]; 2807a1931ceSMatthew G. Knepley PetscInt d, f, comp; 2817a1931ceSMatthew G. Knepley 2827a1931ceSMatthew G. Knepley ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 2837a1931ceSMatthew G. Knepley if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 2847a1931ceSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr); 2857a1931ceSMatthew G. Knepley for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp]; 2867a1931ceSMatthew G. Knepley 2877a1931ceSMatthew G. Knepley for (d = 0; d < ctx->dim; ++d) { 2887a1931ceSMatthew G. Knepley xi[d] = 0.0; 2897a1931ceSMatthew G. Knepley for (f = 0; f < ctx->dim; ++f) xi[d] += invJ[d*ctx->dim+f]*0.5*PetscRealPart(coords[p*ctx->dim+f] - v0[f]); 2907a1931ceSMatthew G. Knepley for (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]; 2917a1931ceSMatthew G. Knepley } 2927a1931ceSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, NULL, &x);CHKERRQ(ierr); 2937a1931ceSMatthew G. Knepley } 2947a1931ceSMatthew G. Knepley ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 2957a1931ceSMatthew G. Knepley ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr); 2967a1931ceSMatthew G. Knepley ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr); 2977a1931ceSMatthew G. Knepley PetscFunctionReturn(0); 2987a1931ceSMatthew G. Knepley } 2997a1931ceSMatthew G. Knepley 3007a1931ceSMatthew G. Knepley #undef __FUNCT__ 301552f7358SJed Brown #define __FUNCT__ "QuadMap_Private" 3025820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) 303552f7358SJed Brown { 304552f7358SJed Brown const PetscScalar *vertices = (const PetscScalar*) ctx; 305552f7358SJed Brown const PetscScalar x0 = vertices[0]; 306552f7358SJed Brown const PetscScalar y0 = vertices[1]; 307552f7358SJed Brown const PetscScalar x1 = vertices[2]; 308552f7358SJed Brown const PetscScalar y1 = vertices[3]; 309552f7358SJed Brown const PetscScalar x2 = vertices[4]; 310552f7358SJed Brown const PetscScalar y2 = vertices[5]; 311552f7358SJed Brown const PetscScalar x3 = vertices[6]; 312552f7358SJed Brown const PetscScalar y3 = vertices[7]; 313552f7358SJed Brown const PetscScalar f_1 = x1 - x0; 314552f7358SJed Brown const PetscScalar g_1 = y1 - y0; 315552f7358SJed Brown const PetscScalar f_3 = x3 - x0; 316552f7358SJed Brown const PetscScalar g_3 = y3 - y0; 317552f7358SJed Brown const PetscScalar f_01 = x2 - x1 - x3 + x0; 318552f7358SJed Brown const PetscScalar g_01 = y2 - y1 - y3 + y0; 319552f7358SJed Brown PetscScalar *ref, *real; 320552f7358SJed Brown PetscErrorCode ierr; 321552f7358SJed Brown 322552f7358SJed Brown PetscFunctionBegin; 323552f7358SJed Brown ierr = VecGetArray(Xref, &ref);CHKERRQ(ierr); 324552f7358SJed Brown ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr); 325552f7358SJed Brown { 326552f7358SJed Brown const PetscScalar p0 = ref[0]; 327552f7358SJed Brown const PetscScalar p1 = ref[1]; 328552f7358SJed Brown 329552f7358SJed Brown real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1; 330552f7358SJed Brown real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1; 331552f7358SJed Brown } 332552f7358SJed Brown ierr = PetscLogFlops(28);CHKERRQ(ierr); 333552f7358SJed Brown ierr = VecRestoreArray(Xref, &ref);CHKERRQ(ierr); 334552f7358SJed Brown ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr); 335552f7358SJed Brown PetscFunctionReturn(0); 336552f7358SJed Brown } 337552f7358SJed Brown 338c0dedaeaSBarry Smith #include <petsc-private/dmimpl.h> 339552f7358SJed Brown #undef __FUNCT__ 340552f7358SJed Brown #define __FUNCT__ "QuadJacobian_Private" 3415820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx) 342552f7358SJed Brown { 343552f7358SJed Brown const PetscScalar *vertices = (const PetscScalar*) ctx; 344552f7358SJed Brown const PetscScalar x0 = vertices[0]; 345552f7358SJed Brown const PetscScalar y0 = vertices[1]; 346552f7358SJed Brown const PetscScalar x1 = vertices[2]; 347552f7358SJed Brown const PetscScalar y1 = vertices[3]; 348552f7358SJed Brown const PetscScalar x2 = vertices[4]; 349552f7358SJed Brown const PetscScalar y2 = vertices[5]; 350552f7358SJed Brown const PetscScalar x3 = vertices[6]; 351552f7358SJed Brown const PetscScalar y3 = vertices[7]; 352552f7358SJed Brown const PetscScalar f_01 = x2 - x1 - x3 + x0; 353552f7358SJed Brown const PetscScalar g_01 = y2 - y1 - y3 + y0; 354552f7358SJed Brown PetscScalar *ref; 355552f7358SJed Brown PetscErrorCode ierr; 356552f7358SJed Brown 357552f7358SJed Brown PetscFunctionBegin; 358552f7358SJed Brown ierr = VecGetArray(Xref, &ref);CHKERRQ(ierr); 359552f7358SJed Brown { 360552f7358SJed Brown const PetscScalar x = ref[0]; 361552f7358SJed Brown const PetscScalar y = ref[1]; 362552f7358SJed Brown const PetscInt rows[2] = {0, 1}; 363da80777bSKarl Rupp PetscScalar values[4]; 364da80777bSKarl Rupp 365da80777bSKarl Rupp values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5; 366da80777bSKarl Rupp values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5; 367552f7358SJed Brown ierr = MatSetValues(*J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr); 368552f7358SJed Brown } 369552f7358SJed Brown ierr = PetscLogFlops(30);CHKERRQ(ierr); 370552f7358SJed Brown ierr = VecRestoreArray(Xref, &ref);CHKERRQ(ierr); 371552f7358SJed Brown ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 372552f7358SJed Brown ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 373552f7358SJed Brown PetscFunctionReturn(0); 374552f7358SJed Brown } 375552f7358SJed Brown 376552f7358SJed Brown #undef __FUNCT__ 377552f7358SJed Brown #define __FUNCT__ "DMInterpolate_Quad_Private" 378a6dfd86eSKarl Rupp PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 379a6dfd86eSKarl Rupp { 380fafc0619SMatthew G Knepley DM dmCoord; 381552f7358SJed Brown SNES snes; 382552f7358SJed Brown KSP ksp; 383552f7358SJed Brown PC pc; 384552f7358SJed Brown Vec coordsLocal, r, ref, real; 385552f7358SJed Brown Mat J; 386552f7358SJed Brown PetscScalar *a, *coords; 387552f7358SJed Brown PetscInt p; 388552f7358SJed Brown PetscErrorCode ierr; 389552f7358SJed Brown 390552f7358SJed Brown PetscFunctionBegin; 391552f7358SJed Brown ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 392fafc0619SMatthew G Knepley ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr); 393552f7358SJed Brown ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr); 394552f7358SJed Brown ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr); 395552f7358SJed Brown ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr); 396552f7358SJed Brown ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr); 397c0dedaeaSBarry Smith ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr); 398552f7358SJed Brown ierr = VecDuplicate(r, &ref);CHKERRQ(ierr); 399552f7358SJed Brown ierr = VecDuplicate(r, &real);CHKERRQ(ierr); 400552f7358SJed Brown ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr); 401552f7358SJed Brown ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr); 402552f7358SJed Brown ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr); 403552f7358SJed Brown ierr = MatSetUp(J);CHKERRQ(ierr); 4040298fd71SBarry Smith ierr = SNESSetFunction(snes, r, QuadMap_Private, NULL);CHKERRQ(ierr); 4050298fd71SBarry Smith ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, NULL);CHKERRQ(ierr); 406552f7358SJed Brown ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 407552f7358SJed Brown ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 408552f7358SJed Brown ierr = PCSetType(pc, PCLU);CHKERRQ(ierr); 409552f7358SJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 410552f7358SJed Brown 411552f7358SJed Brown ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr); 412552f7358SJed Brown ierr = VecGetArray(v, &a);CHKERRQ(ierr); 413552f7358SJed Brown for (p = 0; p < ctx->n; ++p) { 414a1e44745SMatthew G. Knepley PetscScalar *x = NULL, *vertices = NULL; 415552f7358SJed Brown PetscScalar *xi; 416cb313848SJed Brown PetscReal xir[2]; 417552f7358SJed Brown PetscInt c = ctx->cells[p], comp, coordSize, xSize; 418552f7358SJed Brown 419552f7358SJed Brown /* Can make this do all points at once */ 4200298fd71SBarry Smith ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 4210adebc6cSBarry Smith if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 4*2); 4220298fd71SBarry Smith ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 4230adebc6cSBarry Smith if (4*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 4*ctx->dof); 4240298fd71SBarry Smith ierr = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr); 4250298fd71SBarry Smith ierr = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr); 426552f7358SJed Brown ierr = VecGetArray(real, &xi);CHKERRQ(ierr); 427552f7358SJed Brown xi[0] = coords[p*ctx->dim+0]; 428552f7358SJed Brown xi[1] = coords[p*ctx->dim+1]; 429552f7358SJed Brown ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr); 430552f7358SJed Brown ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr); 431552f7358SJed Brown ierr = VecGetArray(ref, &xi);CHKERRQ(ierr); 432cb313848SJed Brown xir[0] = PetscRealPart(xi[0]); 433cb313848SJed Brown xir[1] = PetscRealPart(xi[1]); 4341aa26658SKarl Rupp for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp]*(1 - xir[0])*(1 - xir[1]) + x[1*ctx->dof+comp]*xir[0]*(1 - xir[1]) + x[2*ctx->dof+comp]*xir[0]*xir[1] + x[3*ctx->dof+comp]*(1 - xir[0])*xir[1]; 4351aa26658SKarl Rupp 436552f7358SJed Brown ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr); 4370298fd71SBarry Smith ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 4380298fd71SBarry Smith ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 439552f7358SJed Brown } 440552f7358SJed Brown ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 441552f7358SJed Brown ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr); 442552f7358SJed Brown 443552f7358SJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 444552f7358SJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 445552f7358SJed Brown ierr = VecDestroy(&ref);CHKERRQ(ierr); 446552f7358SJed Brown ierr = VecDestroy(&real);CHKERRQ(ierr); 447552f7358SJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 448552f7358SJed Brown PetscFunctionReturn(0); 449552f7358SJed Brown } 450552f7358SJed Brown 451552f7358SJed Brown #undef __FUNCT__ 452552f7358SJed Brown #define __FUNCT__ "HexMap_Private" 4535820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) 454552f7358SJed Brown { 455552f7358SJed Brown const PetscScalar *vertices = (const PetscScalar*) ctx; 456552f7358SJed Brown const PetscScalar x0 = vertices[0]; 457552f7358SJed Brown const PetscScalar y0 = vertices[1]; 458552f7358SJed Brown const PetscScalar z0 = vertices[2]; 4597a1931ceSMatthew G. Knepley const PetscScalar x1 = vertices[9]; 4607a1931ceSMatthew G. Knepley const PetscScalar y1 = vertices[10]; 4617a1931ceSMatthew G. Knepley const PetscScalar z1 = vertices[11]; 462552f7358SJed Brown const PetscScalar x2 = vertices[6]; 463552f7358SJed Brown const PetscScalar y2 = vertices[7]; 464552f7358SJed Brown const PetscScalar z2 = vertices[8]; 4657a1931ceSMatthew G. Knepley const PetscScalar x3 = vertices[3]; 4667a1931ceSMatthew G. Knepley const PetscScalar y3 = vertices[4]; 4677a1931ceSMatthew G. Knepley const PetscScalar z3 = vertices[5]; 468552f7358SJed Brown const PetscScalar x4 = vertices[12]; 469552f7358SJed Brown const PetscScalar y4 = vertices[13]; 470552f7358SJed Brown const PetscScalar z4 = vertices[14]; 471552f7358SJed Brown const PetscScalar x5 = vertices[15]; 472552f7358SJed Brown const PetscScalar y5 = vertices[16]; 473552f7358SJed Brown const PetscScalar z5 = vertices[17]; 474552f7358SJed Brown const PetscScalar x6 = vertices[18]; 475552f7358SJed Brown const PetscScalar y6 = vertices[19]; 476552f7358SJed Brown const PetscScalar z6 = vertices[20]; 477552f7358SJed Brown const PetscScalar x7 = vertices[21]; 478552f7358SJed Brown const PetscScalar y7 = vertices[22]; 479552f7358SJed Brown const PetscScalar z7 = vertices[23]; 480552f7358SJed Brown const PetscScalar f_1 = x1 - x0; 481552f7358SJed Brown const PetscScalar g_1 = y1 - y0; 482552f7358SJed Brown const PetscScalar h_1 = z1 - z0; 483552f7358SJed Brown const PetscScalar f_3 = x3 - x0; 484552f7358SJed Brown const PetscScalar g_3 = y3 - y0; 485552f7358SJed Brown const PetscScalar h_3 = z3 - z0; 486552f7358SJed Brown const PetscScalar f_4 = x4 - x0; 487552f7358SJed Brown const PetscScalar g_4 = y4 - y0; 488552f7358SJed Brown const PetscScalar h_4 = z4 - z0; 489552f7358SJed Brown const PetscScalar f_01 = x2 - x1 - x3 + x0; 490552f7358SJed Brown const PetscScalar g_01 = y2 - y1 - y3 + y0; 491552f7358SJed Brown const PetscScalar h_01 = z2 - z1 - z3 + z0; 492552f7358SJed Brown const PetscScalar f_12 = x7 - x3 - x4 + x0; 493552f7358SJed Brown const PetscScalar g_12 = y7 - y3 - y4 + y0; 494552f7358SJed Brown const PetscScalar h_12 = z7 - z3 - z4 + z0; 495552f7358SJed Brown const PetscScalar f_02 = x5 - x1 - x4 + x0; 496552f7358SJed Brown const PetscScalar g_02 = y5 - y1 - y4 + y0; 497552f7358SJed Brown const PetscScalar h_02 = z5 - z1 - z4 + z0; 498552f7358SJed Brown const PetscScalar f_012 = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 499552f7358SJed Brown const PetscScalar g_012 = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 500552f7358SJed Brown const PetscScalar h_012 = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 501552f7358SJed Brown PetscScalar *ref, *real; 502552f7358SJed Brown PetscErrorCode ierr; 503552f7358SJed Brown 504552f7358SJed Brown PetscFunctionBegin; 505552f7358SJed Brown ierr = VecGetArray(Xref, &ref);CHKERRQ(ierr); 506552f7358SJed Brown ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr); 507552f7358SJed Brown { 508552f7358SJed Brown const PetscScalar p0 = ref[0]; 509552f7358SJed Brown const PetscScalar p1 = ref[1]; 510552f7358SJed Brown const PetscScalar p2 = ref[2]; 511552f7358SJed Brown 512552f7358SJed Brown 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; 513552f7358SJed Brown 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; 514552f7358SJed Brown 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; 515552f7358SJed Brown } 516552f7358SJed Brown ierr = PetscLogFlops(114);CHKERRQ(ierr); 517552f7358SJed Brown ierr = VecRestoreArray(Xref, &ref);CHKERRQ(ierr); 518552f7358SJed Brown ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr); 519552f7358SJed Brown PetscFunctionReturn(0); 520552f7358SJed Brown } 521552f7358SJed Brown 522552f7358SJed Brown #undef __FUNCT__ 523552f7358SJed Brown #define __FUNCT__ "HexJacobian_Private" 5245820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx) 525552f7358SJed Brown { 526552f7358SJed Brown const PetscScalar *vertices = (const PetscScalar*) ctx; 527552f7358SJed Brown const PetscScalar x0 = vertices[0]; 528552f7358SJed Brown const PetscScalar y0 = vertices[1]; 529552f7358SJed Brown const PetscScalar z0 = vertices[2]; 5307a1931ceSMatthew G. Knepley const PetscScalar x1 = vertices[9]; 5317a1931ceSMatthew G. Knepley const PetscScalar y1 = vertices[10]; 5327a1931ceSMatthew G. Knepley const PetscScalar z1 = vertices[11]; 533552f7358SJed Brown const PetscScalar x2 = vertices[6]; 534552f7358SJed Brown const PetscScalar y2 = vertices[7]; 535552f7358SJed Brown const PetscScalar z2 = vertices[8]; 5367a1931ceSMatthew G. Knepley const PetscScalar x3 = vertices[3]; 5377a1931ceSMatthew G. Knepley const PetscScalar y3 = vertices[4]; 5387a1931ceSMatthew G. Knepley const PetscScalar z3 = vertices[5]; 539552f7358SJed Brown const PetscScalar x4 = vertices[12]; 540552f7358SJed Brown const PetscScalar y4 = vertices[13]; 541552f7358SJed Brown const PetscScalar z4 = vertices[14]; 542552f7358SJed Brown const PetscScalar x5 = vertices[15]; 543552f7358SJed Brown const PetscScalar y5 = vertices[16]; 544552f7358SJed Brown const PetscScalar z5 = vertices[17]; 545552f7358SJed Brown const PetscScalar x6 = vertices[18]; 546552f7358SJed Brown const PetscScalar y6 = vertices[19]; 547552f7358SJed Brown const PetscScalar z6 = vertices[20]; 548552f7358SJed Brown const PetscScalar x7 = vertices[21]; 549552f7358SJed Brown const PetscScalar y7 = vertices[22]; 550552f7358SJed Brown const PetscScalar z7 = vertices[23]; 551552f7358SJed Brown const PetscScalar f_xy = x2 - x1 - x3 + x0; 552552f7358SJed Brown const PetscScalar g_xy = y2 - y1 - y3 + y0; 553552f7358SJed Brown const PetscScalar h_xy = z2 - z1 - z3 + z0; 554552f7358SJed Brown const PetscScalar f_yz = x7 - x3 - x4 + x0; 555552f7358SJed Brown const PetscScalar g_yz = y7 - y3 - y4 + y0; 556552f7358SJed Brown const PetscScalar h_yz = z7 - z3 - z4 + z0; 557552f7358SJed Brown const PetscScalar f_xz = x5 - x1 - x4 + x0; 558552f7358SJed Brown const PetscScalar g_xz = y5 - y1 - y4 + y0; 559552f7358SJed Brown const PetscScalar h_xz = z5 - z1 - z4 + z0; 560552f7358SJed Brown const PetscScalar f_xyz = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 561552f7358SJed Brown const PetscScalar g_xyz = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 562552f7358SJed Brown const PetscScalar h_xyz = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 563552f7358SJed Brown PetscScalar *ref; 564552f7358SJed Brown PetscErrorCode ierr; 565552f7358SJed Brown 566552f7358SJed Brown PetscFunctionBegin; 567552f7358SJed Brown ierr = VecGetArray(Xref, &ref);CHKERRQ(ierr); 568552f7358SJed Brown { 569552f7358SJed Brown const PetscScalar x = ref[0]; 570552f7358SJed Brown const PetscScalar y = ref[1]; 571552f7358SJed Brown const PetscScalar z = ref[2]; 572552f7358SJed Brown const PetscInt rows[3] = {0, 1, 2}; 573da80777bSKarl Rupp PetscScalar values[9]; 574da80777bSKarl Rupp 575da80777bSKarl Rupp values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0; 576da80777bSKarl Rupp values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0; 577da80777bSKarl Rupp values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0; 578da80777bSKarl Rupp values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0; 579da80777bSKarl Rupp values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0; 580da80777bSKarl Rupp values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0; 581da80777bSKarl Rupp values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0; 582da80777bSKarl Rupp values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0; 583da80777bSKarl Rupp values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0; 5841aa26658SKarl Rupp 585552f7358SJed Brown ierr = MatSetValues(*J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr); 586552f7358SJed Brown } 587552f7358SJed Brown ierr = PetscLogFlops(152);CHKERRQ(ierr); 588552f7358SJed Brown ierr = VecRestoreArray(Xref, &ref);CHKERRQ(ierr); 589552f7358SJed Brown ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 590552f7358SJed Brown ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 591552f7358SJed Brown PetscFunctionReturn(0); 592552f7358SJed Brown } 593552f7358SJed Brown 594552f7358SJed Brown #undef __FUNCT__ 595552f7358SJed Brown #define __FUNCT__ "DMInterpolate_Hex_Private" 596a6dfd86eSKarl Rupp PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 597a6dfd86eSKarl Rupp { 598fafc0619SMatthew G Knepley DM dmCoord; 599552f7358SJed Brown SNES snes; 600552f7358SJed Brown KSP ksp; 601552f7358SJed Brown PC pc; 602552f7358SJed Brown Vec coordsLocal, r, ref, real; 603552f7358SJed Brown Mat J; 604552f7358SJed Brown PetscScalar *a, *coords; 605552f7358SJed Brown PetscInt p; 606552f7358SJed Brown PetscErrorCode ierr; 607552f7358SJed Brown 608552f7358SJed Brown PetscFunctionBegin; 609552f7358SJed Brown ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 610fafc0619SMatthew G Knepley ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr); 611552f7358SJed Brown ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr); 612552f7358SJed Brown ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr); 613552f7358SJed Brown ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr); 614552f7358SJed Brown ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr); 615c0dedaeaSBarry Smith ierr = VecSetType(r,dm->vectype);CHKERRQ(ierr); 616552f7358SJed Brown ierr = VecDuplicate(r, &ref);CHKERRQ(ierr); 617552f7358SJed Brown ierr = VecDuplicate(r, &real);CHKERRQ(ierr); 618552f7358SJed Brown ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr); 619552f7358SJed Brown ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr); 620552f7358SJed Brown ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr); 621552f7358SJed Brown ierr = MatSetUp(J);CHKERRQ(ierr); 6220298fd71SBarry Smith ierr = SNESSetFunction(snes, r, HexMap_Private, NULL);CHKERRQ(ierr); 6230298fd71SBarry Smith ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, NULL);CHKERRQ(ierr); 624552f7358SJed Brown ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 625552f7358SJed Brown ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 626552f7358SJed Brown ierr = PCSetType(pc, PCLU);CHKERRQ(ierr); 627552f7358SJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 628552f7358SJed Brown 629552f7358SJed Brown ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr); 630552f7358SJed Brown ierr = VecGetArray(v, &a);CHKERRQ(ierr); 631552f7358SJed Brown for (p = 0; p < ctx->n; ++p) { 632a1e44745SMatthew G. Knepley PetscScalar *x = NULL, *vertices = NULL; 633552f7358SJed Brown PetscScalar *xi; 634cb313848SJed Brown PetscReal xir[3]; 635552f7358SJed Brown PetscInt c = ctx->cells[p], comp, coordSize, xSize; 636552f7358SJed Brown 637552f7358SJed Brown /* Can make this do all points at once */ 6380298fd71SBarry Smith ierr = DMPlexVecGetClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 6390adebc6cSBarry Smith if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 8*3); 6400298fd71SBarry Smith ierr = DMPlexVecGetClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 6410adebc6cSBarry Smith if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 8*ctx->dof); 6420298fd71SBarry Smith ierr = SNESSetFunction(snes, NULL, NULL, (void*) vertices);CHKERRQ(ierr); 6430298fd71SBarry Smith ierr = SNESSetJacobian(snes, NULL, NULL, NULL, (void*) vertices);CHKERRQ(ierr); 644552f7358SJed Brown ierr = VecGetArray(real, &xi);CHKERRQ(ierr); 645552f7358SJed Brown xi[0] = coords[p*ctx->dim+0]; 646552f7358SJed Brown xi[1] = coords[p*ctx->dim+1]; 647552f7358SJed Brown xi[2] = coords[p*ctx->dim+2]; 648552f7358SJed Brown ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr); 649552f7358SJed Brown ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr); 650552f7358SJed Brown ierr = VecGetArray(ref, &xi);CHKERRQ(ierr); 651cb313848SJed Brown xir[0] = PetscRealPart(xi[0]); 652cb313848SJed Brown xir[1] = PetscRealPart(xi[1]); 653cb313848SJed Brown xir[2] = PetscRealPart(xi[2]); 654552f7358SJed Brown for (comp = 0; comp < ctx->dof; ++comp) { 655552f7358SJed Brown a[p*ctx->dof+comp] = 656cb313848SJed Brown x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) + 6577a1931ceSMatthew G. Knepley x[3*ctx->dof+comp]* xir[0]*(1-xir[1])*(1-xir[2]) + 658cb313848SJed Brown x[2*ctx->dof+comp]* xir[0]* xir[1]*(1-xir[2]) + 6597a1931ceSMatthew G. Knepley x[1*ctx->dof+comp]*(1-xir[0])* xir[1]*(1-xir[2]) + 660cb313848SJed Brown x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])* xir[2] + 661cb313848SJed Brown x[5*ctx->dof+comp]* xir[0]*(1-xir[1])* xir[2] + 662cb313848SJed Brown x[6*ctx->dof+comp]* xir[0]* xir[1]* xir[2] + 663cb313848SJed Brown x[7*ctx->dof+comp]*(1-xir[0])* xir[1]* xir[2]; 664552f7358SJed Brown } 665552f7358SJed Brown ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr); 6660298fd71SBarry Smith ierr = DMPlexVecRestoreClosure(dmCoord, NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 6670298fd71SBarry Smith ierr = DMPlexVecRestoreClosure(dm, NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 668552f7358SJed Brown } 669552f7358SJed Brown ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 670552f7358SJed Brown ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr); 671552f7358SJed Brown 672552f7358SJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 673552f7358SJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 674552f7358SJed Brown ierr = VecDestroy(&ref);CHKERRQ(ierr); 675552f7358SJed Brown ierr = VecDestroy(&real);CHKERRQ(ierr); 676552f7358SJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 677552f7358SJed Brown PetscFunctionReturn(0); 678552f7358SJed Brown } 679552f7358SJed Brown 680552f7358SJed Brown #undef __FUNCT__ 681552f7358SJed Brown #define __FUNCT__ "DMInterpolationEvaluate" 682552f7358SJed Brown /* 683552f7358SJed Brown Input Parameters: 684552f7358SJed Brown + ctx - The DMInterpolationInfo context 685552f7358SJed Brown . dm - The DM 686552f7358SJed Brown - x - The local vector containing the field to be interpolated 687552f7358SJed Brown 688552f7358SJed Brown Output Parameters: 689552f7358SJed Brown . v - The vector containing the interpolated values 690552f7358SJed Brown */ 6910adebc6cSBarry Smith PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v) 6920adebc6cSBarry Smith { 693552f7358SJed Brown PetscInt dim, coneSize, n; 694552f7358SJed Brown PetscErrorCode ierr; 695552f7358SJed Brown 696552f7358SJed Brown PetscFunctionBegin; 697552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 698552f7358SJed Brown PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 699552f7358SJed Brown PetscValidHeaderSpecific(v, VEC_CLASSID, 4); 700552f7358SJed Brown ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr); 7010adebc6cSBarry Smith if (n != ctx->n*ctx->dof) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid input vector size %d should be %d", n, ctx->n*ctx->dof); 702552f7358SJed Brown if (n) { 703552f7358SJed Brown ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 704552f7358SJed Brown ierr = DMPlexGetConeSize(dm, ctx->cells[0], &coneSize);CHKERRQ(ierr); 705552f7358SJed Brown if (dim == 2) { 706552f7358SJed Brown if (coneSize == 3) { 7077a1931ceSMatthew G. Knepley ierr = DMInterpolate_Triangle_Private(ctx, dm, x, v);CHKERRQ(ierr); 708552f7358SJed Brown } else if (coneSize == 4) { 709552f7358SJed Brown ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr); 7100adebc6cSBarry Smith } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim); 711552f7358SJed Brown } else if (dim == 3) { 712552f7358SJed Brown if (coneSize == 4) { 7137a1931ceSMatthew G. Knepley ierr = DMInterpolate_Tetrahedron_Private(ctx, dm, x, v);CHKERRQ(ierr); 714552f7358SJed Brown } else { 715552f7358SJed Brown ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr); 716552f7358SJed Brown } 7170adebc6cSBarry Smith } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim); 718552f7358SJed Brown } 719552f7358SJed Brown PetscFunctionReturn(0); 720552f7358SJed Brown } 721552f7358SJed Brown 722552f7358SJed Brown #undef __FUNCT__ 723552f7358SJed Brown #define __FUNCT__ "DMInterpolationDestroy" 7240adebc6cSBarry Smith PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx) 7250adebc6cSBarry Smith { 726552f7358SJed Brown PetscErrorCode ierr; 727552f7358SJed Brown 728552f7358SJed Brown PetscFunctionBegin; 729552f7358SJed Brown PetscValidPointer(ctx, 2); 730552f7358SJed Brown ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr); 731552f7358SJed Brown ierr = PetscFree((*ctx)->points);CHKERRQ(ierr); 732552f7358SJed Brown ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr); 733552f7358SJed Brown ierr = PetscFree(*ctx);CHKERRQ(ierr); 7340298fd71SBarry Smith *ctx = NULL; 735552f7358SJed Brown PetscFunctionReturn(0); 736552f7358SJed Brown } 737