1552f7358SJed Brown #include <petscdmplex.h> /*I "petscdmplex.h" I*/ 2552f7358SJed Brown #include <petscsnes.h> /*I "petscsnes.h" I*/ 3552f7358SJed Brown 4552f7358SJed Brown #undef __FUNCT__ 5552f7358SJed Brown #define __FUNCT__ "DMInterpolationCreate" 60adebc6cSBarry Smith PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx) 70adebc6cSBarry Smith { 8552f7358SJed Brown PetscErrorCode ierr; 9552f7358SJed Brown 10552f7358SJed Brown PetscFunctionBegin; 11552f7358SJed Brown PetscValidPointer(ctx, 2); 12552f7358SJed Brown ierr = PetscMalloc(sizeof(struct _DMInterpolationInfo), ctx);CHKERRQ(ierr); 131aa26658SKarl Rupp 14552f7358SJed Brown (*ctx)->comm = comm; 15552f7358SJed Brown (*ctx)->dim = -1; 16552f7358SJed Brown (*ctx)->nInput = 0; 17552f7358SJed Brown (*ctx)->points = PETSC_NULL; 18552f7358SJed Brown (*ctx)->cells = PETSC_NULL; 19552f7358SJed Brown (*ctx)->n = -1; 20552f7358SJed Brown (*ctx)->coords = PETSC_NULL; 21552f7358SJed Brown PetscFunctionReturn(0); 22552f7358SJed Brown } 23552f7358SJed Brown 24552f7358SJed Brown #undef __FUNCT__ 25552f7358SJed Brown #define __FUNCT__ "DMInterpolationSetDim" 260adebc6cSBarry Smith PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim) 270adebc6cSBarry Smith { 28552f7358SJed Brown PetscFunctionBegin; 290adebc6cSBarry Smith if ((dim < 1) || (dim > 3)) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %d", dim); 30552f7358SJed Brown ctx->dim = dim; 31552f7358SJed Brown PetscFunctionReturn(0); 32552f7358SJed Brown } 33552f7358SJed Brown 34552f7358SJed Brown #undef __FUNCT__ 35552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetDim" 360adebc6cSBarry Smith PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim) 370adebc6cSBarry Smith { 38552f7358SJed Brown PetscFunctionBegin; 39552f7358SJed Brown PetscValidIntPointer(dim, 2); 40552f7358SJed Brown *dim = ctx->dim; 41552f7358SJed Brown PetscFunctionReturn(0); 42552f7358SJed Brown } 43552f7358SJed Brown 44552f7358SJed Brown #undef __FUNCT__ 45552f7358SJed Brown #define __FUNCT__ "DMInterpolationSetDof" 460adebc6cSBarry Smith PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof) 470adebc6cSBarry Smith { 48552f7358SJed Brown PetscFunctionBegin; 490adebc6cSBarry Smith if (dof < 1) SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %d", dof); 50552f7358SJed Brown ctx->dof = dof; 51552f7358SJed Brown PetscFunctionReturn(0); 52552f7358SJed Brown } 53552f7358SJed Brown 54552f7358SJed Brown #undef __FUNCT__ 55552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetDof" 560adebc6cSBarry Smith PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof) 570adebc6cSBarry Smith { 58552f7358SJed Brown PetscFunctionBegin; 59552f7358SJed Brown PetscValidIntPointer(dof, 2); 60552f7358SJed Brown *dof = ctx->dof; 61552f7358SJed Brown PetscFunctionReturn(0); 62552f7358SJed Brown } 63552f7358SJed Brown 64552f7358SJed Brown #undef __FUNCT__ 65552f7358SJed Brown #define __FUNCT__ "DMInterpolationAddPoints" 660adebc6cSBarry Smith PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[]) 670adebc6cSBarry Smith { 68552f7358SJed Brown PetscErrorCode ierr; 69552f7358SJed Brown 70552f7358SJed Brown PetscFunctionBegin; 710adebc6cSBarry Smith if (ctx->dim < 0) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 720adebc6cSBarry Smith if (ctx->points) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet"); 73552f7358SJed Brown ctx->nInput = n; 741aa26658SKarl Rupp 75552f7358SJed Brown ierr = PetscMalloc(n*ctx->dim * sizeof(PetscReal), &ctx->points);CHKERRQ(ierr); 76552f7358SJed Brown ierr = PetscMemcpy(ctx->points, points, n*ctx->dim * sizeof(PetscReal));CHKERRQ(ierr); 77552f7358SJed Brown PetscFunctionReturn(0); 78552f7358SJed Brown } 79552f7358SJed Brown 80552f7358SJed Brown #undef __FUNCT__ 81552f7358SJed Brown #define __FUNCT__ "DMInterpolationSetUp" 820adebc6cSBarry Smith PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints) 830adebc6cSBarry Smith { 84552f7358SJed Brown MPI_Comm comm = ctx->comm; 85552f7358SJed Brown PetscScalar *a; 86552f7358SJed Brown PetscInt p, q, i; 87552f7358SJed Brown PetscMPIInt rank, size; 88552f7358SJed Brown PetscErrorCode ierr; 89552f7358SJed Brown Vec pointVec; 90552f7358SJed Brown IS cellIS; 91552f7358SJed Brown PetscLayout layout; 92552f7358SJed Brown PetscReal *globalPoints; 93cb313848SJed Brown PetscScalar *globalPointsScalar; 94552f7358SJed Brown const PetscInt *ranges; 95552f7358SJed Brown PetscMPIInt *counts, *displs; 96552f7358SJed Brown const PetscInt *foundCells; 97552f7358SJed Brown PetscMPIInt *foundProcs, *globalProcs; 9819436ca2SJed Brown PetscInt n, N; 99552f7358SJed Brown 10019436ca2SJed Brown PetscFunctionBegin; 10119436ca2SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 10219436ca2SJed Brown ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 10319436ca2SJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1040adebc6cSBarry Smith if (ctx->dim < 0) SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 10519436ca2SJed Brown /* Locate points */ 10619436ca2SJed Brown n = ctx->nInput; 107552f7358SJed Brown if (!redundantPoints) { 108552f7358SJed Brown ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr); 109552f7358SJed Brown ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr); 110552f7358SJed Brown ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr); 111552f7358SJed Brown ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr); 112552f7358SJed Brown ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr); 113552f7358SJed Brown /* Communicate all points to all processes */ 114552f7358SJed Brown ierr = PetscMalloc3(N*ctx->dim,PetscReal,&globalPoints,size,PetscMPIInt,&counts,size,PetscMPIInt,&displs);CHKERRQ(ierr); 115552f7358SJed Brown ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr); 116552f7358SJed Brown for (p = 0; p < size; ++p) { 117552f7358SJed Brown counts[p] = (ranges[p+1] - ranges[p])*ctx->dim; 118552f7358SJed Brown displs[p] = ranges[p]*ctx->dim; 119552f7358SJed Brown } 120552f7358SJed Brown ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRQ(ierr); 121552f7358SJed Brown } else { 122552f7358SJed Brown N = n; 123*f5af7f23SKarl Rupp 124552f7358SJed Brown globalPoints = ctx->points; 125552f7358SJed Brown } 126552f7358SJed Brown #if 0 127552f7358SJed Brown ierr = PetscMalloc3(N,PetscInt,&foundCells,N,PetscMPIInt,&foundProcs,N,PetscMPIInt,&globalProcs);CHKERRQ(ierr); 12819436ca2SJed Brown /* foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); */ 129552f7358SJed Brown #else 130cb313848SJed Brown #if defined(PETSC_USE_COMPLEX) 131cb313848SJed Brown ierr = PetscMalloc(N*sizeof(PetscScalar),&globalPointsScalar);CHKERRQ(ierr); 132cb313848SJed Brown for (i=0; i<N; i++) globalPointsScalar[i] = globalPoints[i]; 133cb313848SJed Brown #else 134cb313848SJed Brown globalPointsScalar = globalPoints; 135cb313848SJed Brown #endif 13604706141SMatthew G Knepley ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N*ctx->dim, globalPointsScalar, &pointVec);CHKERRQ(ierr); 137552f7358SJed Brown ierr = PetscMalloc2(N,PetscMPIInt,&foundProcs,N,PetscMPIInt,&globalProcs);CHKERRQ(ierr); 138552f7358SJed Brown ierr = DMLocatePoints(dm, pointVec, &cellIS);CHKERRQ(ierr); 139552f7358SJed Brown ierr = ISGetIndices(cellIS, &foundCells);CHKERRQ(ierr); 140552f7358SJed Brown #endif 141552f7358SJed Brown for (p = 0; p < N; ++p) { 1421aa26658SKarl Rupp if (foundCells[p] >= 0) foundProcs[p] = rank; 1431aa26658SKarl Rupp else foundProcs[p] = size; 144552f7358SJed Brown } 145552f7358SJed Brown /* Let the lowest rank process own each point */ 146efab3cc2SJed Brown ierr = MPI_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);CHKERRQ(ierr); 147552f7358SJed Brown ctx->n = 0; 148552f7358SJed Brown for (p = 0; p < N; ++p) { 1490adebc6cSBarry 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); 1501aa26658SKarl Rupp else if (globalProcs[p] == rank) ctx->n++; 151552f7358SJed Brown } 152552f7358SJed Brown /* Create coordinates vector and array of owned cells */ 153552f7358SJed Brown ierr = PetscMalloc(ctx->n * sizeof(PetscInt), &ctx->cells);CHKERRQ(ierr); 154552f7358SJed Brown ierr = VecCreate(comm, &ctx->coords);CHKERRQ(ierr); 155552f7358SJed Brown ierr = VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);CHKERRQ(ierr); 156552f7358SJed Brown ierr = VecSetBlockSize(ctx->coords, ctx->dim);CHKERRQ(ierr); 157552f7358SJed Brown ierr = VecSetFromOptions(ctx->coords);CHKERRQ(ierr); 158552f7358SJed Brown ierr = VecGetArray(ctx->coords, &a);CHKERRQ(ierr); 159552f7358SJed Brown for (p = 0, q = 0, i = 0; p < N; ++p) { 160552f7358SJed Brown if (globalProcs[p] == rank) { 161552f7358SJed Brown PetscInt d; 162552f7358SJed Brown 1631aa26658SKarl Rupp for (d = 0; d < ctx->dim; ++d, ++i) a[i] = globalPoints[p*ctx->dim+d]; 164552f7358SJed Brown ctx->cells[q++] = foundCells[p]; 165552f7358SJed Brown } 166552f7358SJed Brown } 167552f7358SJed Brown ierr = VecRestoreArray(ctx->coords, &a);CHKERRQ(ierr); 168552f7358SJed Brown #if 0 169552f7358SJed Brown ierr = PetscFree3(foundCells,foundProcs,globalProcs);CHKERRQ(ierr); 170552f7358SJed Brown #else 171552f7358SJed Brown ierr = PetscFree2(foundProcs,globalProcs);CHKERRQ(ierr); 172552f7358SJed Brown ierr = ISRestoreIndices(cellIS, &foundCells);CHKERRQ(ierr); 173552f7358SJed Brown ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 174552f7358SJed Brown ierr = VecDestroy(&pointVec);CHKERRQ(ierr); 175552f7358SJed Brown #endif 176cb313848SJed Brown if ((void*)globalPointsScalar != (void*)globalPoints) {ierr = PetscFree(globalPointsScalar);CHKERRQ(ierr);} 177552f7358SJed Brown if (!redundantPoints) { 178552f7358SJed Brown ierr = PetscFree3(globalPoints,counts,displs);CHKERRQ(ierr); 179552f7358SJed Brown ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr); 180552f7358SJed Brown } 181552f7358SJed Brown PetscFunctionReturn(0); 182552f7358SJed Brown } 183552f7358SJed Brown 184552f7358SJed Brown #undef __FUNCT__ 185552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetCoordinates" 1860adebc6cSBarry Smith PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates) 1870adebc6cSBarry Smith { 188552f7358SJed Brown PetscFunctionBegin; 189552f7358SJed Brown PetscValidPointer(coordinates, 2); 1900adebc6cSBarry Smith if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 191552f7358SJed Brown *coordinates = ctx->coords; 192552f7358SJed Brown PetscFunctionReturn(0); 193552f7358SJed Brown } 194552f7358SJed Brown 195552f7358SJed Brown #undef __FUNCT__ 196552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetVector" 1970adebc6cSBarry Smith PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v) 1980adebc6cSBarry Smith { 199552f7358SJed Brown PetscErrorCode ierr; 200552f7358SJed Brown 201552f7358SJed Brown PetscFunctionBegin; 202552f7358SJed Brown PetscValidPointer(v, 2); 2030adebc6cSBarry Smith if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 204552f7358SJed Brown ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr); 205552f7358SJed Brown ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr); 206552f7358SJed Brown ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr); 207552f7358SJed Brown ierr = VecSetFromOptions(*v);CHKERRQ(ierr); 208552f7358SJed Brown PetscFunctionReturn(0); 209552f7358SJed Brown } 210552f7358SJed Brown 211552f7358SJed Brown #undef __FUNCT__ 212552f7358SJed Brown #define __FUNCT__ "DMInterpolationRestoreVector" 2130adebc6cSBarry Smith PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v) 2140adebc6cSBarry Smith { 215552f7358SJed Brown PetscErrorCode ierr; 216552f7358SJed Brown 217552f7358SJed Brown PetscFunctionBegin; 218552f7358SJed Brown PetscValidPointer(v, 2); 2190adebc6cSBarry Smith if (!ctx->coords) SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup."); 220552f7358SJed Brown ierr = VecDestroy(v);CHKERRQ(ierr); 221552f7358SJed Brown PetscFunctionReturn(0); 222552f7358SJed Brown } 223552f7358SJed Brown 224552f7358SJed Brown #undef __FUNCT__ 225552f7358SJed Brown #define __FUNCT__ "DMInterpolate_Simplex_Private" 226a6dfd86eSKarl Rupp PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Simplex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 227a6dfd86eSKarl Rupp { 228552f7358SJed Brown PetscReal *v0, *J, *invJ, detJ; 229552f7358SJed Brown PetscScalar *a, *coords; 230552f7358SJed Brown PetscInt p; 231552f7358SJed Brown PetscErrorCode ierr; 232552f7358SJed Brown 233552f7358SJed Brown PetscFunctionBegin; 234552f7358SJed Brown ierr = PetscMalloc3(ctx->dim,PetscReal,&v0,ctx->dim*ctx->dim,PetscReal,&J,ctx->dim*ctx->dim,PetscReal,&invJ);CHKERRQ(ierr); 235552f7358SJed Brown ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr); 236552f7358SJed Brown ierr = VecGetArray(v, &a);CHKERRQ(ierr); 237552f7358SJed Brown for (p = 0; p < ctx->n; ++p) { 238552f7358SJed Brown PetscInt c = ctx->cells[p]; 239552f7358SJed Brown const PetscScalar *x; 240552f7358SJed Brown PetscReal xi[4]; 241552f7358SJed Brown PetscInt d, f, comp; 242552f7358SJed Brown 243552f7358SJed Brown ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 244552f7358SJed Brown if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 245552f7358SJed Brown ierr = DMPlexVecGetClosure(dm, PETSC_NULL, xLocal, c, PETSC_NULL, &x);CHKERRQ(ierr); 2461aa26658SKarl Rupp for (comp = 0; comp < ctx->dof; ++comp) a[p*ctx->dof+comp] = x[0*ctx->dof+comp]; 2471aa26658SKarl Rupp 248552f7358SJed Brown for (d = 0; d < ctx->dim; ++d) { 249552f7358SJed Brown xi[d] = 0.0; 2501aa26658SKarl 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]); 2511aa26658SKarl 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]; 252552f7358SJed Brown } 253552f7358SJed Brown ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, xLocal, c, PETSC_NULL, &x);CHKERRQ(ierr); 254552f7358SJed Brown } 255552f7358SJed Brown ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 256552f7358SJed Brown ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr); 257552f7358SJed Brown ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr); 258552f7358SJed Brown PetscFunctionReturn(0); 259552f7358SJed Brown } 260552f7358SJed Brown 261552f7358SJed Brown #undef __FUNCT__ 262552f7358SJed Brown #define __FUNCT__ "QuadMap_Private" 2635820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) 264552f7358SJed Brown { 265552f7358SJed Brown const PetscScalar *vertices = (const PetscScalar*) ctx; 266552f7358SJed Brown const PetscScalar x0 = vertices[0]; 267552f7358SJed Brown const PetscScalar y0 = vertices[1]; 268552f7358SJed Brown const PetscScalar x1 = vertices[2]; 269552f7358SJed Brown const PetscScalar y1 = vertices[3]; 270552f7358SJed Brown const PetscScalar x2 = vertices[4]; 271552f7358SJed Brown const PetscScalar y2 = vertices[5]; 272552f7358SJed Brown const PetscScalar x3 = vertices[6]; 273552f7358SJed Brown const PetscScalar y3 = vertices[7]; 274552f7358SJed Brown const PetscScalar f_1 = x1 - x0; 275552f7358SJed Brown const PetscScalar g_1 = y1 - y0; 276552f7358SJed Brown const PetscScalar f_3 = x3 - x0; 277552f7358SJed Brown const PetscScalar g_3 = y3 - y0; 278552f7358SJed Brown const PetscScalar f_01 = x2 - x1 - x3 + x0; 279552f7358SJed Brown const PetscScalar g_01 = y2 - y1 - y3 + y0; 280552f7358SJed Brown PetscScalar *ref, *real; 281552f7358SJed Brown PetscErrorCode ierr; 282552f7358SJed Brown 283552f7358SJed Brown PetscFunctionBegin; 284552f7358SJed Brown ierr = VecGetArray(Xref, &ref);CHKERRQ(ierr); 285552f7358SJed Brown ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr); 286552f7358SJed Brown { 287552f7358SJed Brown const PetscScalar p0 = ref[0]; 288552f7358SJed Brown const PetscScalar p1 = ref[1]; 289552f7358SJed Brown 290552f7358SJed Brown real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1; 291552f7358SJed Brown real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1; 292552f7358SJed Brown } 293552f7358SJed Brown ierr = PetscLogFlops(28);CHKERRQ(ierr); 294552f7358SJed Brown ierr = VecRestoreArray(Xref, &ref);CHKERRQ(ierr); 295552f7358SJed Brown ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr); 296552f7358SJed Brown PetscFunctionReturn(0); 297552f7358SJed Brown } 298552f7358SJed Brown 299552f7358SJed Brown #undef __FUNCT__ 300552f7358SJed Brown #define __FUNCT__ "QuadJacobian_Private" 3015820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx) 302552f7358SJed Brown { 303552f7358SJed Brown const PetscScalar *vertices = (const PetscScalar*) ctx; 304552f7358SJed Brown const PetscScalar x0 = vertices[0]; 305552f7358SJed Brown const PetscScalar y0 = vertices[1]; 306552f7358SJed Brown const PetscScalar x1 = vertices[2]; 307552f7358SJed Brown const PetscScalar y1 = vertices[3]; 308552f7358SJed Brown const PetscScalar x2 = vertices[4]; 309552f7358SJed Brown const PetscScalar y2 = vertices[5]; 310552f7358SJed Brown const PetscScalar x3 = vertices[6]; 311552f7358SJed Brown const PetscScalar y3 = vertices[7]; 312552f7358SJed Brown const PetscScalar f_01 = x2 - x1 - x3 + x0; 313552f7358SJed Brown const PetscScalar g_01 = y2 - y1 - y3 + y0; 314552f7358SJed Brown PetscScalar *ref; 315552f7358SJed Brown PetscErrorCode ierr; 316552f7358SJed Brown 317552f7358SJed Brown PetscFunctionBegin; 318552f7358SJed Brown ierr = VecGetArray(Xref, &ref);CHKERRQ(ierr); 319552f7358SJed Brown { 320552f7358SJed Brown const PetscScalar x = ref[0]; 321552f7358SJed Brown const PetscScalar y = ref[1]; 322552f7358SJed Brown const PetscInt rows[2] = {0, 1}; 323da80777bSKarl Rupp PetscScalar values[4]; 324da80777bSKarl Rupp 325da80777bSKarl Rupp values[0] = (x1 - x0 + f_01*y) * 0.5; values[1] = (x3 - x0 + f_01*x) * 0.5; 326da80777bSKarl Rupp values[2] = (y1 - y0 + g_01*y) * 0.5; values[3] = (y3 - y0 + g_01*x) * 0.5; 327552f7358SJed Brown ierr = MatSetValues(*J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr); 328552f7358SJed Brown } 329552f7358SJed Brown ierr = PetscLogFlops(30);CHKERRQ(ierr); 330552f7358SJed Brown ierr = VecRestoreArray(Xref, &ref);CHKERRQ(ierr); 331552f7358SJed Brown ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 332552f7358SJed Brown ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 333552f7358SJed Brown PetscFunctionReturn(0); 334552f7358SJed Brown } 335552f7358SJed Brown 336552f7358SJed Brown #undef __FUNCT__ 337552f7358SJed Brown #define __FUNCT__ "DMInterpolate_Quad_Private" 338a6dfd86eSKarl Rupp PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 339a6dfd86eSKarl Rupp { 340fafc0619SMatthew G Knepley DM dmCoord; 341552f7358SJed Brown SNES snes; 342552f7358SJed Brown KSP ksp; 343552f7358SJed Brown PC pc; 344552f7358SJed Brown Vec coordsLocal, r, ref, real; 345552f7358SJed Brown Mat J; 346552f7358SJed Brown PetscScalar *a, *coords; 347552f7358SJed Brown PetscInt p; 348552f7358SJed Brown PetscErrorCode ierr; 349552f7358SJed Brown 350552f7358SJed Brown PetscFunctionBegin; 351552f7358SJed Brown ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 352fafc0619SMatthew G Knepley ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr); 353552f7358SJed Brown ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr); 354552f7358SJed Brown ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr); 355552f7358SJed Brown ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr); 356552f7358SJed Brown ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr); 357552f7358SJed Brown ierr = VecSetFromOptions(r);CHKERRQ(ierr); 358552f7358SJed Brown ierr = VecDuplicate(r, &ref);CHKERRQ(ierr); 359552f7358SJed Brown ierr = VecDuplicate(r, &real);CHKERRQ(ierr); 360552f7358SJed Brown ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr); 361552f7358SJed Brown ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr); 362552f7358SJed Brown ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr); 363552f7358SJed Brown ierr = MatSetUp(J);CHKERRQ(ierr); 364552f7358SJed Brown ierr = SNESSetFunction(snes, r, QuadMap_Private, PETSC_NULL);CHKERRQ(ierr); 365552f7358SJed Brown ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, PETSC_NULL);CHKERRQ(ierr); 366552f7358SJed Brown ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 367552f7358SJed Brown ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 368552f7358SJed Brown ierr = PCSetType(pc, PCLU);CHKERRQ(ierr); 369552f7358SJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 370552f7358SJed Brown 371552f7358SJed Brown ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr); 372552f7358SJed Brown ierr = VecGetArray(v, &a);CHKERRQ(ierr); 373552f7358SJed Brown for (p = 0; p < ctx->n; ++p) { 374552f7358SJed Brown const PetscScalar *x, *vertices; 375552f7358SJed Brown PetscScalar *xi; 376cb313848SJed Brown PetscReal xir[2]; 377552f7358SJed Brown PetscInt c = ctx->cells[p], comp, coordSize, xSize; 378552f7358SJed Brown 379552f7358SJed Brown /* Can make this do all points at once */ 380fafc0619SMatthew G Knepley ierr = DMPlexVecGetClosure(dmCoord, PETSC_NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 3810adebc6cSBarry Smith if (4*2 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 4*2); 382552f7358SJed Brown ierr = DMPlexVecGetClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 3830adebc6cSBarry Smith if (4*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 4*ctx->dof); 38488c71bb6SMatthew G Knepley ierr = SNESSetFunction(snes, PETSC_NULL, PETSC_NULL, (void*) vertices);CHKERRQ(ierr); 38588c71bb6SMatthew G Knepley ierr = SNESSetJacobian(snes, PETSC_NULL, PETSC_NULL, PETSC_NULL, (void*) vertices);CHKERRQ(ierr); 386552f7358SJed Brown ierr = VecGetArray(real, &xi);CHKERRQ(ierr); 387552f7358SJed Brown xi[0] = coords[p*ctx->dim+0]; 388552f7358SJed Brown xi[1] = coords[p*ctx->dim+1]; 389552f7358SJed Brown ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr); 390552f7358SJed Brown ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr); 391552f7358SJed Brown ierr = VecGetArray(ref, &xi);CHKERRQ(ierr); 392cb313848SJed Brown xir[0] = PetscRealPart(xi[0]); 393cb313848SJed Brown xir[1] = PetscRealPart(xi[1]); 3941aa26658SKarl 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]; 3951aa26658SKarl Rupp 396552f7358SJed Brown ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr); 397fafc0619SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dmCoord, PETSC_NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 398552f7358SJed Brown ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 399552f7358SJed Brown } 400552f7358SJed Brown ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 401552f7358SJed Brown ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr); 402552f7358SJed Brown 403552f7358SJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 404552f7358SJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 405552f7358SJed Brown ierr = VecDestroy(&ref);CHKERRQ(ierr); 406552f7358SJed Brown ierr = VecDestroy(&real);CHKERRQ(ierr); 407552f7358SJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 408552f7358SJed Brown PetscFunctionReturn(0); 409552f7358SJed Brown } 410552f7358SJed Brown 411552f7358SJed Brown #undef __FUNCT__ 412552f7358SJed Brown #define __FUNCT__ "HexMap_Private" 4135820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) 414552f7358SJed Brown { 415552f7358SJed Brown const PetscScalar *vertices = (const PetscScalar*) ctx; 416552f7358SJed Brown const PetscScalar x0 = vertices[0]; 417552f7358SJed Brown const PetscScalar y0 = vertices[1]; 418552f7358SJed Brown const PetscScalar z0 = vertices[2]; 419552f7358SJed Brown const PetscScalar x1 = vertices[3]; 420552f7358SJed Brown const PetscScalar y1 = vertices[4]; 421552f7358SJed Brown const PetscScalar z1 = vertices[5]; 422552f7358SJed Brown const PetscScalar x2 = vertices[6]; 423552f7358SJed Brown const PetscScalar y2 = vertices[7]; 424552f7358SJed Brown const PetscScalar z2 = vertices[8]; 425552f7358SJed Brown const PetscScalar x3 = vertices[9]; 426552f7358SJed Brown const PetscScalar y3 = vertices[10]; 427552f7358SJed Brown const PetscScalar z3 = vertices[11]; 428552f7358SJed Brown const PetscScalar x4 = vertices[12]; 429552f7358SJed Brown const PetscScalar y4 = vertices[13]; 430552f7358SJed Brown const PetscScalar z4 = vertices[14]; 431552f7358SJed Brown const PetscScalar x5 = vertices[15]; 432552f7358SJed Brown const PetscScalar y5 = vertices[16]; 433552f7358SJed Brown const PetscScalar z5 = vertices[17]; 434552f7358SJed Brown const PetscScalar x6 = vertices[18]; 435552f7358SJed Brown const PetscScalar y6 = vertices[19]; 436552f7358SJed Brown const PetscScalar z6 = vertices[20]; 437552f7358SJed Brown const PetscScalar x7 = vertices[21]; 438552f7358SJed Brown const PetscScalar y7 = vertices[22]; 439552f7358SJed Brown const PetscScalar z7 = vertices[23]; 440552f7358SJed Brown const PetscScalar f_1 = x1 - x0; 441552f7358SJed Brown const PetscScalar g_1 = y1 - y0; 442552f7358SJed Brown const PetscScalar h_1 = z1 - z0; 443552f7358SJed Brown const PetscScalar f_3 = x3 - x0; 444552f7358SJed Brown const PetscScalar g_3 = y3 - y0; 445552f7358SJed Brown const PetscScalar h_3 = z3 - z0; 446552f7358SJed Brown const PetscScalar f_4 = x4 - x0; 447552f7358SJed Brown const PetscScalar g_4 = y4 - y0; 448552f7358SJed Brown const PetscScalar h_4 = z4 - z0; 449552f7358SJed Brown const PetscScalar f_01 = x2 - x1 - x3 + x0; 450552f7358SJed Brown const PetscScalar g_01 = y2 - y1 - y3 + y0; 451552f7358SJed Brown const PetscScalar h_01 = z2 - z1 - z3 + z0; 452552f7358SJed Brown const PetscScalar f_12 = x7 - x3 - x4 + x0; 453552f7358SJed Brown const PetscScalar g_12 = y7 - y3 - y4 + y0; 454552f7358SJed Brown const PetscScalar h_12 = z7 - z3 - z4 + z0; 455552f7358SJed Brown const PetscScalar f_02 = x5 - x1 - x4 + x0; 456552f7358SJed Brown const PetscScalar g_02 = y5 - y1 - y4 + y0; 457552f7358SJed Brown const PetscScalar h_02 = z5 - z1 - z4 + z0; 458552f7358SJed Brown const PetscScalar f_012 = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 459552f7358SJed Brown const PetscScalar g_012 = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 460552f7358SJed Brown const PetscScalar h_012 = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 461552f7358SJed Brown PetscScalar *ref, *real; 462552f7358SJed Brown PetscErrorCode ierr; 463552f7358SJed Brown 464552f7358SJed Brown PetscFunctionBegin; 465552f7358SJed Brown ierr = VecGetArray(Xref, &ref);CHKERRQ(ierr); 466552f7358SJed Brown ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr); 467552f7358SJed Brown { 468552f7358SJed Brown const PetscScalar p0 = ref[0]; 469552f7358SJed Brown const PetscScalar p1 = ref[1]; 470552f7358SJed Brown const PetscScalar p2 = ref[2]; 471552f7358SJed Brown 472552f7358SJed 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; 473552f7358SJed 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; 474552f7358SJed 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; 475552f7358SJed Brown } 476552f7358SJed Brown ierr = PetscLogFlops(114);CHKERRQ(ierr); 477552f7358SJed Brown ierr = VecRestoreArray(Xref, &ref);CHKERRQ(ierr); 478552f7358SJed Brown ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr); 479552f7358SJed Brown PetscFunctionReturn(0); 480552f7358SJed Brown } 481552f7358SJed Brown 482552f7358SJed Brown #undef __FUNCT__ 483552f7358SJed Brown #define __FUNCT__ "HexJacobian_Private" 4845820edbdSMatthew G Knepley PETSC_STATIC_INLINE PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx) 485552f7358SJed Brown { 486552f7358SJed Brown const PetscScalar *vertices = (const PetscScalar*) ctx; 487552f7358SJed Brown const PetscScalar x0 = vertices[0]; 488552f7358SJed Brown const PetscScalar y0 = vertices[1]; 489552f7358SJed Brown const PetscScalar z0 = vertices[2]; 490552f7358SJed Brown const PetscScalar x1 = vertices[3]; 491552f7358SJed Brown const PetscScalar y1 = vertices[4]; 492552f7358SJed Brown const PetscScalar z1 = vertices[5]; 493552f7358SJed Brown const PetscScalar x2 = vertices[6]; 494552f7358SJed Brown const PetscScalar y2 = vertices[7]; 495552f7358SJed Brown const PetscScalar z2 = vertices[8]; 496552f7358SJed Brown const PetscScalar x3 = vertices[9]; 497552f7358SJed Brown const PetscScalar y3 = vertices[10]; 498552f7358SJed Brown const PetscScalar z3 = vertices[11]; 499552f7358SJed Brown const PetscScalar x4 = vertices[12]; 500552f7358SJed Brown const PetscScalar y4 = vertices[13]; 501552f7358SJed Brown const PetscScalar z4 = vertices[14]; 502552f7358SJed Brown const PetscScalar x5 = vertices[15]; 503552f7358SJed Brown const PetscScalar y5 = vertices[16]; 504552f7358SJed Brown const PetscScalar z5 = vertices[17]; 505552f7358SJed Brown const PetscScalar x6 = vertices[18]; 506552f7358SJed Brown const PetscScalar y6 = vertices[19]; 507552f7358SJed Brown const PetscScalar z6 = vertices[20]; 508552f7358SJed Brown const PetscScalar x7 = vertices[21]; 509552f7358SJed Brown const PetscScalar y7 = vertices[22]; 510552f7358SJed Brown const PetscScalar z7 = vertices[23]; 511552f7358SJed Brown const PetscScalar f_xy = x2 - x1 - x3 + x0; 512552f7358SJed Brown const PetscScalar g_xy = y2 - y1 - y3 + y0; 513552f7358SJed Brown const PetscScalar h_xy = z2 - z1 - z3 + z0; 514552f7358SJed Brown const PetscScalar f_yz = x7 - x3 - x4 + x0; 515552f7358SJed Brown const PetscScalar g_yz = y7 - y3 - y4 + y0; 516552f7358SJed Brown const PetscScalar h_yz = z7 - z3 - z4 + z0; 517552f7358SJed Brown const PetscScalar f_xz = x5 - x1 - x4 + x0; 518552f7358SJed Brown const PetscScalar g_xz = y5 - y1 - y4 + y0; 519552f7358SJed Brown const PetscScalar h_xz = z5 - z1 - z4 + z0; 520552f7358SJed Brown const PetscScalar f_xyz = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 521552f7358SJed Brown const PetscScalar g_xyz = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 522552f7358SJed Brown const PetscScalar h_xyz = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 523552f7358SJed Brown PetscScalar *ref; 524552f7358SJed Brown PetscErrorCode ierr; 525552f7358SJed Brown 526552f7358SJed Brown PetscFunctionBegin; 527552f7358SJed Brown ierr = VecGetArray(Xref, &ref);CHKERRQ(ierr); 528552f7358SJed Brown { 529552f7358SJed Brown const PetscScalar x = ref[0]; 530552f7358SJed Brown const PetscScalar y = ref[1]; 531552f7358SJed Brown const PetscScalar z = ref[2]; 532552f7358SJed Brown const PetscInt rows[3] = {0, 1, 2}; 533da80777bSKarl Rupp PetscScalar values[9]; 534da80777bSKarl Rupp 535da80777bSKarl Rupp values[0] = (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0; 536da80777bSKarl Rupp values[1] = (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0; 537da80777bSKarl Rupp values[2] = (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0; 538da80777bSKarl Rupp values[3] = (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0; 539da80777bSKarl Rupp values[4] = (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0; 540da80777bSKarl Rupp values[5] = (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0; 541da80777bSKarl Rupp values[6] = (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0; 542da80777bSKarl Rupp values[7] = (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0; 543da80777bSKarl Rupp values[8] = (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0; 5441aa26658SKarl Rupp 545552f7358SJed Brown ierr = MatSetValues(*J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr); 546552f7358SJed Brown } 547552f7358SJed Brown ierr = PetscLogFlops(152);CHKERRQ(ierr); 548552f7358SJed Brown ierr = VecRestoreArray(Xref, &ref);CHKERRQ(ierr); 549552f7358SJed Brown ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 550552f7358SJed Brown ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 551552f7358SJed Brown PetscFunctionReturn(0); 552552f7358SJed Brown } 553552f7358SJed Brown 554552f7358SJed Brown #undef __FUNCT__ 555552f7358SJed Brown #define __FUNCT__ "DMInterpolate_Hex_Private" 556a6dfd86eSKarl Rupp PETSC_STATIC_INLINE PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) 557a6dfd86eSKarl Rupp { 558fafc0619SMatthew G Knepley DM dmCoord; 559552f7358SJed Brown SNES snes; 560552f7358SJed Brown KSP ksp; 561552f7358SJed Brown PC pc; 562552f7358SJed Brown Vec coordsLocal, r, ref, real; 563552f7358SJed Brown Mat J; 564552f7358SJed Brown PetscScalar *a, *coords; 565552f7358SJed Brown PetscInt p; 566552f7358SJed Brown PetscErrorCode ierr; 567552f7358SJed Brown 568552f7358SJed Brown PetscFunctionBegin; 569552f7358SJed Brown ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 570fafc0619SMatthew G Knepley ierr = DMGetCoordinateDM(dm, &dmCoord);CHKERRQ(ierr); 571552f7358SJed Brown ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr); 572552f7358SJed Brown ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr); 573552f7358SJed Brown ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr); 574552f7358SJed Brown ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr); 575552f7358SJed Brown ierr = VecSetFromOptions(r);CHKERRQ(ierr); 576552f7358SJed Brown ierr = VecDuplicate(r, &ref);CHKERRQ(ierr); 577552f7358SJed Brown ierr = VecDuplicate(r, &real);CHKERRQ(ierr); 578552f7358SJed Brown ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr); 579552f7358SJed Brown ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr); 580552f7358SJed Brown ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr); 581552f7358SJed Brown ierr = MatSetUp(J);CHKERRQ(ierr); 582552f7358SJed Brown ierr = SNESSetFunction(snes, r, HexMap_Private, PETSC_NULL);CHKERRQ(ierr); 583552f7358SJed Brown ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, PETSC_NULL);CHKERRQ(ierr); 584552f7358SJed Brown ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 585552f7358SJed Brown ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 586552f7358SJed Brown ierr = PCSetType(pc, PCLU);CHKERRQ(ierr); 587552f7358SJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 588552f7358SJed Brown 589552f7358SJed Brown ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr); 590552f7358SJed Brown ierr = VecGetArray(v, &a);CHKERRQ(ierr); 591552f7358SJed Brown for (p = 0; p < ctx->n; ++p) { 592552f7358SJed Brown const PetscScalar *x, *vertices; 593552f7358SJed Brown PetscScalar *xi; 594cb313848SJed Brown PetscReal xir[3]; 595552f7358SJed Brown PetscInt c = ctx->cells[p], comp, coordSize, xSize; 596552f7358SJed Brown 597552f7358SJed Brown /* Can make this do all points at once */ 598fafc0619SMatthew G Knepley ierr = DMPlexVecGetClosure(dmCoord, PETSC_NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 5990adebc6cSBarry Smith if (8*3 != coordSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 8*3); 600552f7358SJed Brown ierr = DMPlexVecGetClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 6010adebc6cSBarry Smith if (8*ctx->dof != xSize) SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 8*ctx->dof); 60288c71bb6SMatthew G Knepley ierr = SNESSetFunction(snes, PETSC_NULL, PETSC_NULL, (void*) vertices);CHKERRQ(ierr); 60388c71bb6SMatthew G Knepley ierr = SNESSetJacobian(snes, PETSC_NULL, PETSC_NULL, PETSC_NULL, (void*) vertices);CHKERRQ(ierr); 604552f7358SJed Brown ierr = VecGetArray(real, &xi);CHKERRQ(ierr); 605552f7358SJed Brown xi[0] = coords[p*ctx->dim+0]; 606552f7358SJed Brown xi[1] = coords[p*ctx->dim+1]; 607552f7358SJed Brown xi[2] = coords[p*ctx->dim+2]; 608552f7358SJed Brown ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr); 609552f7358SJed Brown ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr); 610552f7358SJed Brown ierr = VecGetArray(ref, &xi);CHKERRQ(ierr); 611cb313848SJed Brown xir[0] = PetscRealPart(xi[0]); 612cb313848SJed Brown xir[1] = PetscRealPart(xi[1]); 613cb313848SJed Brown xir[2] = PetscRealPart(xi[2]); 614552f7358SJed Brown for (comp = 0; comp < ctx->dof; ++comp) { 615552f7358SJed Brown a[p*ctx->dof+comp] = 616cb313848SJed Brown x[0*ctx->dof+comp]*(1-xir[0])*(1-xir[1])*(1-xir[2]) + 617cb313848SJed Brown x[1*ctx->dof+comp]* xir[0]*(1-xir[1])*(1-xir[2]) + 618cb313848SJed Brown x[2*ctx->dof+comp]* xir[0]* xir[1]*(1-xir[2]) + 619cb313848SJed Brown x[3*ctx->dof+comp]*(1-xir[0])* xir[1]*(1-xir[2]) + 620cb313848SJed Brown x[4*ctx->dof+comp]*(1-xir[0])*(1-xir[1])* xir[2] + 621cb313848SJed Brown x[5*ctx->dof+comp]* xir[0]*(1-xir[1])* xir[2] + 622cb313848SJed Brown x[6*ctx->dof+comp]* xir[0]* xir[1]* xir[2] + 623cb313848SJed Brown x[7*ctx->dof+comp]*(1-xir[0])* xir[1]* xir[2]; 624552f7358SJed Brown } 625552f7358SJed Brown ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr); 626fafc0619SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dmCoord, PETSC_NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 627552f7358SJed Brown ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 628552f7358SJed Brown } 629552f7358SJed Brown ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 630552f7358SJed Brown ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr); 631552f7358SJed Brown 632552f7358SJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 633552f7358SJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 634552f7358SJed Brown ierr = VecDestroy(&ref);CHKERRQ(ierr); 635552f7358SJed Brown ierr = VecDestroy(&real);CHKERRQ(ierr); 636552f7358SJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 637552f7358SJed Brown PetscFunctionReturn(0); 638552f7358SJed Brown } 639552f7358SJed Brown 640552f7358SJed Brown #undef __FUNCT__ 641552f7358SJed Brown #define __FUNCT__ "DMInterpolationEvaluate" 642552f7358SJed Brown /* 643552f7358SJed Brown Input Parameters: 644552f7358SJed Brown + ctx - The DMInterpolationInfo context 645552f7358SJed Brown . dm - The DM 646552f7358SJed Brown - x - The local vector containing the field to be interpolated 647552f7358SJed Brown 648552f7358SJed Brown Output Parameters: 649552f7358SJed Brown . v - The vector containing the interpolated values 650552f7358SJed Brown */ 6510adebc6cSBarry Smith PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v) 6520adebc6cSBarry Smith { 653552f7358SJed Brown PetscInt dim, coneSize, n; 654552f7358SJed Brown PetscErrorCode ierr; 655552f7358SJed Brown 656552f7358SJed Brown PetscFunctionBegin; 657552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 658552f7358SJed Brown PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 659552f7358SJed Brown PetscValidHeaderSpecific(v, VEC_CLASSID, 4); 660552f7358SJed Brown ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr); 6610adebc6cSBarry 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); 662552f7358SJed Brown if (n) { 663552f7358SJed Brown ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 664552f7358SJed Brown ierr = DMPlexGetConeSize(dm, ctx->cells[0], &coneSize);CHKERRQ(ierr); 665552f7358SJed Brown if (dim == 2) { 666552f7358SJed Brown if (coneSize == 3) { 667552f7358SJed Brown ierr = DMInterpolate_Simplex_Private(ctx, dm, x, v);CHKERRQ(ierr); 668552f7358SJed Brown } else if (coneSize == 4) { 669552f7358SJed Brown ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr); 6700adebc6cSBarry Smith } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim); 671552f7358SJed Brown } else if (dim == 3) { 672552f7358SJed Brown if (coneSize == 4) { 673552f7358SJed Brown ierr = DMInterpolate_Simplex_Private(ctx, dm, x, v);CHKERRQ(ierr); 674552f7358SJed Brown } else { 675552f7358SJed Brown ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr); 676552f7358SJed Brown } 6770adebc6cSBarry Smith } else SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim); 678552f7358SJed Brown } 679552f7358SJed Brown PetscFunctionReturn(0); 680552f7358SJed Brown } 681552f7358SJed Brown 682552f7358SJed Brown #undef __FUNCT__ 683552f7358SJed Brown #define __FUNCT__ "DMInterpolationDestroy" 6840adebc6cSBarry Smith PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx) 6850adebc6cSBarry Smith { 686552f7358SJed Brown PetscErrorCode ierr; 687552f7358SJed Brown 688552f7358SJed Brown PetscFunctionBegin; 689552f7358SJed Brown PetscValidPointer(ctx, 2); 690552f7358SJed Brown ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr); 691552f7358SJed Brown ierr = PetscFree((*ctx)->points);CHKERRQ(ierr); 692552f7358SJed Brown ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr); 693552f7358SJed Brown ierr = PetscFree(*ctx);CHKERRQ(ierr); 694552f7358SJed Brown *ctx = PETSC_NULL; 695552f7358SJed Brown PetscFunctionReturn(0); 696552f7358SJed Brown } 697