1*552f7358SJed Brown #include <petscdmplex.h> /*I "petscdmplex.h" I*/ 2*552f7358SJed Brown #include <petscsnes.h> /*I "petscsnes.h" I*/ 3*552f7358SJed Brown 4*552f7358SJed Brown #undef __FUNCT__ 5*552f7358SJed Brown #define __FUNCT__ "DMInterpolationCreate" 6*552f7358SJed Brown PetscErrorCode DMInterpolationCreate(MPI_Comm comm, DMInterpolationInfo *ctx) { 7*552f7358SJed Brown PetscErrorCode ierr; 8*552f7358SJed Brown 9*552f7358SJed Brown PetscFunctionBegin; 10*552f7358SJed Brown PetscValidPointer(ctx, 2); 11*552f7358SJed Brown ierr = PetscMalloc(sizeof(struct _DMInterpolationInfo), ctx);CHKERRQ(ierr); 12*552f7358SJed Brown (*ctx)->comm = comm; 13*552f7358SJed Brown (*ctx)->dim = -1; 14*552f7358SJed Brown (*ctx)->nInput = 0; 15*552f7358SJed Brown (*ctx)->points = PETSC_NULL; 16*552f7358SJed Brown (*ctx)->cells = PETSC_NULL; 17*552f7358SJed Brown (*ctx)->n = -1; 18*552f7358SJed Brown (*ctx)->coords = PETSC_NULL; 19*552f7358SJed Brown PetscFunctionReturn(0); 20*552f7358SJed Brown } 21*552f7358SJed Brown 22*552f7358SJed Brown #undef __FUNCT__ 23*552f7358SJed Brown #define __FUNCT__ "DMInterpolationSetDim" 24*552f7358SJed Brown PetscErrorCode DMInterpolationSetDim(DMInterpolationInfo ctx, PetscInt dim) { 25*552f7358SJed Brown PetscFunctionBegin; 26*552f7358SJed Brown if ((dim < 1) || (dim > 3)) {SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid dimension for points: %d", dim);} 27*552f7358SJed Brown ctx->dim = dim; 28*552f7358SJed Brown PetscFunctionReturn(0); 29*552f7358SJed Brown } 30*552f7358SJed Brown 31*552f7358SJed Brown #undef __FUNCT__ 32*552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetDim" 33*552f7358SJed Brown PetscErrorCode DMInterpolationGetDim(DMInterpolationInfo ctx, PetscInt *dim) { 34*552f7358SJed Brown PetscFunctionBegin; 35*552f7358SJed Brown PetscValidIntPointer(dim, 2); 36*552f7358SJed Brown *dim = ctx->dim; 37*552f7358SJed Brown PetscFunctionReturn(0); 38*552f7358SJed Brown } 39*552f7358SJed Brown 40*552f7358SJed Brown #undef __FUNCT__ 41*552f7358SJed Brown #define __FUNCT__ "DMInterpolationSetDof" 42*552f7358SJed Brown PetscErrorCode DMInterpolationSetDof(DMInterpolationInfo ctx, PetscInt dof) { 43*552f7358SJed Brown PetscFunctionBegin; 44*552f7358SJed Brown if (dof < 1) {SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Invalid number of components: %d", dof);} 45*552f7358SJed Brown ctx->dof = dof; 46*552f7358SJed Brown PetscFunctionReturn(0); 47*552f7358SJed Brown } 48*552f7358SJed Brown 49*552f7358SJed Brown #undef __FUNCT__ 50*552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetDof" 51*552f7358SJed Brown PetscErrorCode DMInterpolationGetDof(DMInterpolationInfo ctx, PetscInt *dof) { 52*552f7358SJed Brown PetscFunctionBegin; 53*552f7358SJed Brown PetscValidIntPointer(dof, 2); 54*552f7358SJed Brown *dof = ctx->dof; 55*552f7358SJed Brown PetscFunctionReturn(0); 56*552f7358SJed Brown } 57*552f7358SJed Brown 58*552f7358SJed Brown #undef __FUNCT__ 59*552f7358SJed Brown #define __FUNCT__ "DMInterpolationAddPoints" 60*552f7358SJed Brown PetscErrorCode DMInterpolationAddPoints(DMInterpolationInfo ctx, PetscInt n, PetscReal points[]) { 61*552f7358SJed Brown PetscErrorCode ierr; 62*552f7358SJed Brown 63*552f7358SJed Brown PetscFunctionBegin; 64*552f7358SJed Brown if (ctx->dim < 0) { 65*552f7358SJed Brown SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 66*552f7358SJed Brown } 67*552f7358SJed Brown if (ctx->points) { 68*552f7358SJed Brown SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "Cannot add points multiple times yet"); 69*552f7358SJed Brown } 70*552f7358SJed Brown ctx->nInput = n; 71*552f7358SJed Brown ierr = PetscMalloc(n*ctx->dim * sizeof(PetscReal), &ctx->points);CHKERRQ(ierr); 72*552f7358SJed Brown ierr = PetscMemcpy(ctx->points, points, n*ctx->dim * sizeof(PetscReal));CHKERRQ(ierr); 73*552f7358SJed Brown PetscFunctionReturn(0); 74*552f7358SJed Brown } 75*552f7358SJed Brown 76*552f7358SJed Brown #undef __FUNCT__ 77*552f7358SJed Brown #define __FUNCT__ "DMInterpolationSetUp" 78*552f7358SJed Brown PetscErrorCode DMInterpolationSetUp(DMInterpolationInfo ctx, DM dm, PetscBool redundantPoints) { 79*552f7358SJed Brown MPI_Comm comm = ctx->comm; 80*552f7358SJed Brown PetscScalar *a; 81*552f7358SJed Brown PetscInt p, q, i; 82*552f7358SJed Brown PetscMPIInt rank, size; 83*552f7358SJed Brown PetscErrorCode ierr; 84*552f7358SJed Brown 85*552f7358SJed Brown PetscFunctionBegin; 86*552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 87*552f7358SJed Brown ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 88*552f7358SJed Brown ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 89*552f7358SJed Brown if (ctx->dim < 0) { 90*552f7358SJed Brown SETERRQ(comm, PETSC_ERR_ARG_WRONGSTATE, "The spatial dimension has not been set"); 91*552f7358SJed Brown } 92*552f7358SJed Brown // Locate points 93*552f7358SJed Brown Vec pointVec; 94*552f7358SJed Brown IS cellIS; 95*552f7358SJed Brown PetscLayout layout; 96*552f7358SJed Brown PetscReal *globalPoints; 97*552f7358SJed Brown const PetscInt *ranges; 98*552f7358SJed Brown PetscMPIInt *counts, *displs; 99*552f7358SJed Brown const PetscInt *foundCells; 100*552f7358SJed Brown PetscMPIInt *foundProcs, *globalProcs; 101*552f7358SJed Brown PetscInt n = ctx->nInput, N; 102*552f7358SJed Brown 103*552f7358SJed Brown if (!redundantPoints) { 104*552f7358SJed Brown ierr = PetscLayoutCreate(comm, &layout);CHKERRQ(ierr); 105*552f7358SJed Brown ierr = PetscLayoutSetBlockSize(layout, 1);CHKERRQ(ierr); 106*552f7358SJed Brown ierr = PetscLayoutSetLocalSize(layout, n);CHKERRQ(ierr); 107*552f7358SJed Brown ierr = PetscLayoutSetUp(layout);CHKERRQ(ierr); 108*552f7358SJed Brown ierr = PetscLayoutGetSize(layout, &N);CHKERRQ(ierr); 109*552f7358SJed Brown /* Communicate all points to all processes */ 110*552f7358SJed Brown ierr = PetscMalloc3(N*ctx->dim,PetscReal,&globalPoints,size,PetscMPIInt,&counts,size,PetscMPIInt,&displs);CHKERRQ(ierr); 111*552f7358SJed Brown ierr = PetscLayoutGetRanges(layout, &ranges);CHKERRQ(ierr); 112*552f7358SJed Brown for (p = 0; p < size; ++p) { 113*552f7358SJed Brown counts[p] = (ranges[p+1] - ranges[p])*ctx->dim; 114*552f7358SJed Brown displs[p] = ranges[p]*ctx->dim; 115*552f7358SJed Brown } 116*552f7358SJed Brown ierr = MPI_Allgatherv(ctx->points, n*ctx->dim, MPIU_REAL, globalPoints, counts, displs, MPIU_REAL, comm);CHKERRQ(ierr); 117*552f7358SJed Brown } else { 118*552f7358SJed Brown N = n; 119*552f7358SJed Brown globalPoints = ctx->points; 120*552f7358SJed Brown } 121*552f7358SJed Brown #if 0 122*552f7358SJed Brown ierr = PetscMalloc3(N,PetscInt,&foundCells,N,PetscMPIInt,&foundProcs,N,PetscMPIInt,&globalProcs);CHKERRQ(ierr); 123*552f7358SJed Brown //foundCells[p] = m->locatePoint(&globalPoints[p*ctx->dim]); 124*552f7358SJed Brown #else 125*552f7358SJed Brown ierr = VecCreateSeqWithArray(PETSC_COMM_SELF, ctx->dim, N, globalPoints, &pointVec);CHKERRQ(ierr); 126*552f7358SJed Brown ierr = PetscMalloc2(N,PetscMPIInt,&foundProcs,N,PetscMPIInt,&globalProcs);CHKERRQ(ierr); 127*552f7358SJed Brown ierr = DMLocatePoints(dm, pointVec, &cellIS);CHKERRQ(ierr); 128*552f7358SJed Brown ierr = ISGetIndices(cellIS, &foundCells);CHKERRQ(ierr); 129*552f7358SJed Brown #endif 130*552f7358SJed Brown for (p = 0; p < N; ++p) { 131*552f7358SJed Brown if (foundCells[p] >= 0) { 132*552f7358SJed Brown foundProcs[p] = rank; 133*552f7358SJed Brown } else { 134*552f7358SJed Brown foundProcs[p] = size; 135*552f7358SJed Brown } 136*552f7358SJed Brown } 137*552f7358SJed Brown /* Let the lowest rank process own each point */ 138*552f7358SJed Brown ierr = MPI_Allreduce(foundProcs, globalProcs, N, MPI_INT, MPI_MIN, comm);CHKERRQ(ierr); 139*552f7358SJed Brown ctx->n = 0; 140*552f7358SJed Brown for (p = 0; p < N; ++p) { 141*552f7358SJed Brown if (globalProcs[p] == size) { 142*552f7358SJed Brown 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); 143*552f7358SJed Brown } else if (globalProcs[p] == rank) { 144*552f7358SJed Brown ctx->n++; 145*552f7358SJed Brown } 146*552f7358SJed Brown } 147*552f7358SJed Brown /* Create coordinates vector and array of owned cells */ 148*552f7358SJed Brown ierr = PetscMalloc(ctx->n * sizeof(PetscInt), &ctx->cells);CHKERRQ(ierr); 149*552f7358SJed Brown ierr = VecCreate(comm, &ctx->coords);CHKERRQ(ierr); 150*552f7358SJed Brown ierr = VecSetSizes(ctx->coords, ctx->n*ctx->dim, PETSC_DECIDE);CHKERRQ(ierr); 151*552f7358SJed Brown ierr = VecSetBlockSize(ctx->coords, ctx->dim);CHKERRQ(ierr); 152*552f7358SJed Brown ierr = VecSetFromOptions(ctx->coords);CHKERRQ(ierr); 153*552f7358SJed Brown ierr = VecGetArray(ctx->coords, &a);CHKERRQ(ierr); 154*552f7358SJed Brown for (p = 0, q = 0, i = 0; p < N; ++p) { 155*552f7358SJed Brown if (globalProcs[p] == rank) { 156*552f7358SJed Brown PetscInt d; 157*552f7358SJed Brown 158*552f7358SJed Brown for (d = 0; d < ctx->dim; ++d, ++i) { 159*552f7358SJed Brown a[i] = globalPoints[p*ctx->dim+d]; 160*552f7358SJed Brown } 161*552f7358SJed Brown ctx->cells[q++] = foundCells[p]; 162*552f7358SJed Brown } 163*552f7358SJed Brown } 164*552f7358SJed Brown ierr = VecRestoreArray(ctx->coords, &a);CHKERRQ(ierr); 165*552f7358SJed Brown #if 0 166*552f7358SJed Brown ierr = PetscFree3(foundCells,foundProcs,globalProcs);CHKERRQ(ierr); 167*552f7358SJed Brown #else 168*552f7358SJed Brown ierr = PetscFree2(foundProcs,globalProcs);CHKERRQ(ierr); 169*552f7358SJed Brown ierr = ISRestoreIndices(cellIS, &foundCells);CHKERRQ(ierr); 170*552f7358SJed Brown ierr = ISDestroy(&cellIS);CHKERRQ(ierr); 171*552f7358SJed Brown ierr = VecDestroy(&pointVec);CHKERRQ(ierr); 172*552f7358SJed Brown #endif 173*552f7358SJed Brown if (!redundantPoints) { 174*552f7358SJed Brown ierr = PetscFree3(globalPoints,counts,displs);CHKERRQ(ierr); 175*552f7358SJed Brown ierr = PetscLayoutDestroy(&layout);CHKERRQ(ierr); 176*552f7358SJed Brown } 177*552f7358SJed Brown PetscFunctionReturn(0); 178*552f7358SJed Brown } 179*552f7358SJed Brown 180*552f7358SJed Brown #undef __FUNCT__ 181*552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetCoordinates" 182*552f7358SJed Brown PetscErrorCode DMInterpolationGetCoordinates(DMInterpolationInfo ctx, Vec *coordinates) { 183*552f7358SJed Brown PetscFunctionBegin; 184*552f7358SJed Brown PetscValidPointer(coordinates, 2); 185*552f7358SJed Brown if (!ctx->coords) {SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");} 186*552f7358SJed Brown *coordinates = ctx->coords; 187*552f7358SJed Brown PetscFunctionReturn(0); 188*552f7358SJed Brown } 189*552f7358SJed Brown 190*552f7358SJed Brown #undef __FUNCT__ 191*552f7358SJed Brown #define __FUNCT__ "DMInterpolationGetVector" 192*552f7358SJed Brown PetscErrorCode DMInterpolationGetVector(DMInterpolationInfo ctx, Vec *v) { 193*552f7358SJed Brown PetscErrorCode ierr; 194*552f7358SJed Brown 195*552f7358SJed Brown PetscFunctionBegin; 196*552f7358SJed Brown PetscValidPointer(v, 2); 197*552f7358SJed Brown if (!ctx->coords) {SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");} 198*552f7358SJed Brown ierr = VecCreate(ctx->comm, v);CHKERRQ(ierr); 199*552f7358SJed Brown ierr = VecSetSizes(*v, ctx->n*ctx->dof, PETSC_DECIDE);CHKERRQ(ierr); 200*552f7358SJed Brown ierr = VecSetBlockSize(*v, ctx->dof);CHKERRQ(ierr); 201*552f7358SJed Brown ierr = VecSetFromOptions(*v);CHKERRQ(ierr); 202*552f7358SJed Brown PetscFunctionReturn(0); 203*552f7358SJed Brown } 204*552f7358SJed Brown 205*552f7358SJed Brown #undef __FUNCT__ 206*552f7358SJed Brown #define __FUNCT__ "DMInterpolationRestoreVector" 207*552f7358SJed Brown PetscErrorCode DMInterpolationRestoreVector(DMInterpolationInfo ctx, Vec *v) { 208*552f7358SJed Brown PetscErrorCode ierr; 209*552f7358SJed Brown 210*552f7358SJed Brown PetscFunctionBegin; 211*552f7358SJed Brown PetscValidPointer(v, 2); 212*552f7358SJed Brown if (!ctx->coords) {SETERRQ(ctx->comm, PETSC_ERR_ARG_WRONGSTATE, "The interpolation context has not been setup.");} 213*552f7358SJed Brown ierr = VecDestroy(v);CHKERRQ(ierr); 214*552f7358SJed Brown PetscFunctionReturn(0); 215*552f7358SJed Brown } 216*552f7358SJed Brown 217*552f7358SJed Brown #undef __FUNCT__ 218*552f7358SJed Brown #define __FUNCT__ "DMInterpolate_Simplex_Private" 219*552f7358SJed Brown PetscErrorCode DMInterpolate_Simplex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) { 220*552f7358SJed Brown PetscReal *v0, *J, *invJ, detJ; 221*552f7358SJed Brown PetscScalar *a, *coords; 222*552f7358SJed Brown PetscInt p; 223*552f7358SJed Brown PetscErrorCode ierr; 224*552f7358SJed Brown 225*552f7358SJed Brown PetscFunctionBegin; 226*552f7358SJed Brown ierr = PetscMalloc3(ctx->dim,PetscReal,&v0,ctx->dim*ctx->dim,PetscReal,&J,ctx->dim*ctx->dim,PetscReal,&invJ);CHKERRQ(ierr); 227*552f7358SJed Brown ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr); 228*552f7358SJed Brown ierr = VecGetArray(v, &a);CHKERRQ(ierr); 229*552f7358SJed Brown for(p = 0; p < ctx->n; ++p) { 230*552f7358SJed Brown PetscInt c = ctx->cells[p]; 231*552f7358SJed Brown const PetscScalar *x; 232*552f7358SJed Brown PetscReal xi[4]; 233*552f7358SJed Brown PetscInt d, f, comp; 234*552f7358SJed Brown 235*552f7358SJed Brown ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 236*552f7358SJed Brown if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 237*552f7358SJed Brown ierr = DMPlexVecGetClosure(dm, PETSC_NULL, xLocal, c, PETSC_NULL, &x);CHKERRQ(ierr); 238*552f7358SJed Brown for(comp = 0; comp < ctx->dof; ++comp) { 239*552f7358SJed Brown a[p*ctx->dof+comp] = x[0*ctx->dof+comp]; 240*552f7358SJed Brown } 241*552f7358SJed Brown for(d = 0; d < ctx->dim; ++d) { 242*552f7358SJed Brown xi[d] = 0.0; 243*552f7358SJed Brown for(f = 0; f < ctx->dim; ++f) { 244*552f7358SJed Brown xi[d] += invJ[d*ctx->dim+f]*0.5*(coords[p*ctx->dim+f] - v0[f]); 245*552f7358SJed Brown } 246*552f7358SJed Brown for(comp = 0; comp < ctx->dof; ++comp) { 247*552f7358SJed Brown a[p*ctx->dof+comp] += (x[(d+1)*ctx->dof+comp] - x[0*ctx->dof+comp])*xi[d]; 248*552f7358SJed Brown } 249*552f7358SJed Brown } 250*552f7358SJed Brown ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, xLocal, c, PETSC_NULL, &x);CHKERRQ(ierr); 251*552f7358SJed Brown } 252*552f7358SJed Brown ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 253*552f7358SJed Brown ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr); 254*552f7358SJed Brown ierr = PetscFree3(v0, J, invJ);CHKERRQ(ierr); 255*552f7358SJed Brown PetscFunctionReturn(0); 256*552f7358SJed Brown } 257*552f7358SJed Brown 258*552f7358SJed Brown #undef __FUNCT__ 259*552f7358SJed Brown #define __FUNCT__ "QuadMap_Private" 260*552f7358SJed Brown PetscErrorCode QuadMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) 261*552f7358SJed Brown { 262*552f7358SJed Brown const PetscScalar*vertices = (const PetscScalar *) ctx; 263*552f7358SJed Brown const PetscScalar x0 = vertices[0]; 264*552f7358SJed Brown const PetscScalar y0 = vertices[1]; 265*552f7358SJed Brown const PetscScalar x1 = vertices[2]; 266*552f7358SJed Brown const PetscScalar y1 = vertices[3]; 267*552f7358SJed Brown const PetscScalar x2 = vertices[4]; 268*552f7358SJed Brown const PetscScalar y2 = vertices[5]; 269*552f7358SJed Brown const PetscScalar x3 = vertices[6]; 270*552f7358SJed Brown const PetscScalar y3 = vertices[7]; 271*552f7358SJed Brown const PetscScalar f_1 = x1 - x0; 272*552f7358SJed Brown const PetscScalar g_1 = y1 - y0; 273*552f7358SJed Brown const PetscScalar f_3 = x3 - x0; 274*552f7358SJed Brown const PetscScalar g_3 = y3 - y0; 275*552f7358SJed Brown const PetscScalar f_01 = x2 - x1 - x3 + x0; 276*552f7358SJed Brown const PetscScalar g_01 = y2 - y1 - y3 + y0; 277*552f7358SJed Brown PetscScalar *ref, *real; 278*552f7358SJed Brown PetscErrorCode ierr; 279*552f7358SJed Brown 280*552f7358SJed Brown PetscFunctionBegin; 281*552f7358SJed Brown ierr = VecGetArray(Xref, &ref);CHKERRQ(ierr); 282*552f7358SJed Brown ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr); 283*552f7358SJed Brown { 284*552f7358SJed Brown const PetscScalar p0 = ref[0]; 285*552f7358SJed Brown const PetscScalar p1 = ref[1]; 286*552f7358SJed Brown 287*552f7358SJed Brown real[0] = x0 + f_1 * p0 + f_3 * p1 + f_01 * p0 * p1; 288*552f7358SJed Brown real[1] = y0 + g_1 * p0 + g_3 * p1 + g_01 * p0 * p1; 289*552f7358SJed Brown } 290*552f7358SJed Brown ierr = PetscLogFlops(28);CHKERRQ(ierr); 291*552f7358SJed Brown ierr = VecRestoreArray(Xref, &ref);CHKERRQ(ierr); 292*552f7358SJed Brown ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr); 293*552f7358SJed Brown PetscFunctionReturn(0); 294*552f7358SJed Brown } 295*552f7358SJed Brown 296*552f7358SJed Brown #undef __FUNCT__ 297*552f7358SJed Brown #define __FUNCT__ "QuadJacobian_Private" 298*552f7358SJed Brown PetscErrorCode QuadJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx) 299*552f7358SJed Brown { 300*552f7358SJed Brown const PetscScalar*vertices = (const PetscScalar *) ctx; 301*552f7358SJed Brown const PetscScalar x0 = vertices[0]; 302*552f7358SJed Brown const PetscScalar y0 = vertices[1]; 303*552f7358SJed Brown const PetscScalar x1 = vertices[2]; 304*552f7358SJed Brown const PetscScalar y1 = vertices[3]; 305*552f7358SJed Brown const PetscScalar x2 = vertices[4]; 306*552f7358SJed Brown const PetscScalar y2 = vertices[5]; 307*552f7358SJed Brown const PetscScalar x3 = vertices[6]; 308*552f7358SJed Brown const PetscScalar y3 = vertices[7]; 309*552f7358SJed Brown const PetscScalar f_01 = x2 - x1 - x3 + x0; 310*552f7358SJed Brown const PetscScalar g_01 = y2 - y1 - y3 + y0; 311*552f7358SJed Brown PetscScalar *ref; 312*552f7358SJed Brown PetscErrorCode ierr; 313*552f7358SJed Brown 314*552f7358SJed Brown PetscFunctionBegin; 315*552f7358SJed Brown ierr = VecGetArray(Xref, &ref);CHKERRQ(ierr); 316*552f7358SJed Brown { 317*552f7358SJed Brown const PetscScalar x = ref[0]; 318*552f7358SJed Brown const PetscScalar y = ref[1]; 319*552f7358SJed Brown const PetscInt rows[2] = {0, 1}; 320*552f7358SJed Brown const PetscScalar values[4] = {(x1 - x0 + f_01*y) * 0.5, (x3 - x0 + f_01*x) * 0.5, 321*552f7358SJed Brown (y1 - y0 + g_01*y) * 0.5, (y3 - y0 + g_01*x) * 0.5}; 322*552f7358SJed Brown ierr = MatSetValues(*J, 2, rows, 2, rows, values, INSERT_VALUES);CHKERRQ(ierr); 323*552f7358SJed Brown } 324*552f7358SJed Brown ierr = PetscLogFlops(30);CHKERRQ(ierr); 325*552f7358SJed Brown ierr = VecRestoreArray(Xref, &ref);CHKERRQ(ierr); 326*552f7358SJed Brown ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 327*552f7358SJed Brown ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 328*552f7358SJed Brown PetscFunctionReturn(0); 329*552f7358SJed Brown } 330*552f7358SJed Brown 331*552f7358SJed Brown #undef __FUNCT__ 332*552f7358SJed Brown #define __FUNCT__ "DMInterpolate_Quad_Private" 333*552f7358SJed Brown PetscErrorCode DMInterpolate_Quad_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) { 334*552f7358SJed Brown SNES snes; 335*552f7358SJed Brown KSP ksp; 336*552f7358SJed Brown PC pc; 337*552f7358SJed Brown Vec coordsLocal, r, ref, real; 338*552f7358SJed Brown Mat J; 339*552f7358SJed Brown PetscScalar *a, *coords; 340*552f7358SJed Brown PetscInt p; 341*552f7358SJed Brown PetscErrorCode ierr; 342*552f7358SJed Brown 343*552f7358SJed Brown PetscFunctionBegin; 344*552f7358SJed Brown ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 345*552f7358SJed Brown ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr); 346*552f7358SJed Brown ierr = SNESSetOptionsPrefix(snes, "quad_interp_");CHKERRQ(ierr); 347*552f7358SJed Brown ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr); 348*552f7358SJed Brown ierr = VecSetSizes(r, 2, 2);CHKERRQ(ierr); 349*552f7358SJed Brown ierr = VecSetFromOptions(r);CHKERRQ(ierr); 350*552f7358SJed Brown ierr = VecDuplicate(r, &ref);CHKERRQ(ierr); 351*552f7358SJed Brown ierr = VecDuplicate(r, &real);CHKERRQ(ierr); 352*552f7358SJed Brown ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr); 353*552f7358SJed Brown ierr = MatSetSizes(J, 2, 2, 2, 2);CHKERRQ(ierr); 354*552f7358SJed Brown ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr); 355*552f7358SJed Brown ierr = MatSetUp(J);CHKERRQ(ierr); 356*552f7358SJed Brown ierr = SNESSetFunction(snes, r, QuadMap_Private, PETSC_NULL);CHKERRQ(ierr); 357*552f7358SJed Brown ierr = SNESSetJacobian(snes, J, J, QuadJacobian_Private, PETSC_NULL);CHKERRQ(ierr); 358*552f7358SJed Brown ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 359*552f7358SJed Brown ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 360*552f7358SJed Brown ierr = PCSetType(pc, PCLU);CHKERRQ(ierr); 361*552f7358SJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 362*552f7358SJed Brown 363*552f7358SJed Brown ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr); 364*552f7358SJed Brown ierr = VecGetArray(v, &a);CHKERRQ(ierr); 365*552f7358SJed Brown for (p = 0; p < ctx->n; ++p) { 366*552f7358SJed Brown const PetscScalar *x, *vertices; 367*552f7358SJed Brown PetscScalar *xi; 368*552f7358SJed Brown PetscInt c = ctx->cells[p], comp, coordSize, xSize; 369*552f7358SJed Brown 370*552f7358SJed Brown /* Can make this do all points at once */ 371*552f7358SJed Brown ierr = DMPlexVecGetClosure(dm, PETSC_NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 372*552f7358SJed Brown if (4*2 != coordSize) {SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 4*2);} 373*552f7358SJed Brown ierr = DMPlexVecGetClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 374*552f7358SJed Brown if (4*ctx->dof != xSize) {SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 4*ctx->dof);} 375*552f7358SJed Brown ierr = SNESSetFunction(snes, PETSC_NULL, PETSC_NULL, vertices);CHKERRQ(ierr); 376*552f7358SJed Brown ierr = SNESSetJacobian(snes, PETSC_NULL, PETSC_NULL, PETSC_NULL, vertices);CHKERRQ(ierr); 377*552f7358SJed Brown ierr = VecGetArray(real, &xi);CHKERRQ(ierr); 378*552f7358SJed Brown xi[0] = coords[p*ctx->dim+0]; 379*552f7358SJed Brown xi[1] = coords[p*ctx->dim+1]; 380*552f7358SJed Brown ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr); 381*552f7358SJed Brown ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr); 382*552f7358SJed Brown ierr = VecGetArray(ref, &xi);CHKERRQ(ierr); 383*552f7358SJed Brown for (comp = 0; comp < ctx->dof; ++comp) { 384*552f7358SJed Brown a[p*ctx->dof+comp] = x[0*ctx->dof+comp]*(1 - xi[0])*(1 - xi[1]) + x[1*ctx->dof+comp]*xi[0]*(1 - xi[1]) + x[2*ctx->dof+comp]*xi[0]*xi[1] + x[3*ctx->dof+comp]*(1 - xi[0])*xi[1]; 385*552f7358SJed Brown } 386*552f7358SJed Brown ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr); 387*552f7358SJed Brown ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 388*552f7358SJed Brown } 389*552f7358SJed Brown ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 390*552f7358SJed Brown ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr); 391*552f7358SJed Brown 392*552f7358SJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 393*552f7358SJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 394*552f7358SJed Brown ierr = VecDestroy(&ref);CHKERRQ(ierr); 395*552f7358SJed Brown ierr = VecDestroy(&real);CHKERRQ(ierr); 396*552f7358SJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 397*552f7358SJed Brown PetscFunctionReturn(0); 398*552f7358SJed Brown } 399*552f7358SJed Brown 400*552f7358SJed Brown #undef __FUNCT__ 401*552f7358SJed Brown #define __FUNCT__ "HexMap_Private" 402*552f7358SJed Brown PetscErrorCode HexMap_Private(SNES snes, Vec Xref, Vec Xreal, void *ctx) 403*552f7358SJed Brown { 404*552f7358SJed Brown const PetscScalar*vertices = (const PetscScalar *) ctx; 405*552f7358SJed Brown const PetscScalar x0 = vertices[0]; 406*552f7358SJed Brown const PetscScalar y0 = vertices[1]; 407*552f7358SJed Brown const PetscScalar z0 = vertices[2]; 408*552f7358SJed Brown const PetscScalar x1 = vertices[3]; 409*552f7358SJed Brown const PetscScalar y1 = vertices[4]; 410*552f7358SJed Brown const PetscScalar z1 = vertices[5]; 411*552f7358SJed Brown const PetscScalar x2 = vertices[6]; 412*552f7358SJed Brown const PetscScalar y2 = vertices[7]; 413*552f7358SJed Brown const PetscScalar z2 = vertices[8]; 414*552f7358SJed Brown const PetscScalar x3 = vertices[9]; 415*552f7358SJed Brown const PetscScalar y3 = vertices[10]; 416*552f7358SJed Brown const PetscScalar z3 = vertices[11]; 417*552f7358SJed Brown const PetscScalar x4 = vertices[12]; 418*552f7358SJed Brown const PetscScalar y4 = vertices[13]; 419*552f7358SJed Brown const PetscScalar z4 = vertices[14]; 420*552f7358SJed Brown const PetscScalar x5 = vertices[15]; 421*552f7358SJed Brown const PetscScalar y5 = vertices[16]; 422*552f7358SJed Brown const PetscScalar z5 = vertices[17]; 423*552f7358SJed Brown const PetscScalar x6 = vertices[18]; 424*552f7358SJed Brown const PetscScalar y6 = vertices[19]; 425*552f7358SJed Brown const PetscScalar z6 = vertices[20]; 426*552f7358SJed Brown const PetscScalar x7 = vertices[21]; 427*552f7358SJed Brown const PetscScalar y7 = vertices[22]; 428*552f7358SJed Brown const PetscScalar z7 = vertices[23]; 429*552f7358SJed Brown const PetscScalar f_1 = x1 - x0; 430*552f7358SJed Brown const PetscScalar g_1 = y1 - y0; 431*552f7358SJed Brown const PetscScalar h_1 = z1 - z0; 432*552f7358SJed Brown const PetscScalar f_3 = x3 - x0; 433*552f7358SJed Brown const PetscScalar g_3 = y3 - y0; 434*552f7358SJed Brown const PetscScalar h_3 = z3 - z0; 435*552f7358SJed Brown const PetscScalar f_4 = x4 - x0; 436*552f7358SJed Brown const PetscScalar g_4 = y4 - y0; 437*552f7358SJed Brown const PetscScalar h_4 = z4 - z0; 438*552f7358SJed Brown const PetscScalar f_01 = x2 - x1 - x3 + x0; 439*552f7358SJed Brown const PetscScalar g_01 = y2 - y1 - y3 + y0; 440*552f7358SJed Brown const PetscScalar h_01 = z2 - z1 - z3 + z0; 441*552f7358SJed Brown const PetscScalar f_12 = x7 - x3 - x4 + x0; 442*552f7358SJed Brown const PetscScalar g_12 = y7 - y3 - y4 + y0; 443*552f7358SJed Brown const PetscScalar h_12 = z7 - z3 - z4 + z0; 444*552f7358SJed Brown const PetscScalar f_02 = x5 - x1 - x4 + x0; 445*552f7358SJed Brown const PetscScalar g_02 = y5 - y1 - y4 + y0; 446*552f7358SJed Brown const PetscScalar h_02 = z5 - z1 - z4 + z0; 447*552f7358SJed Brown const PetscScalar f_012 = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 448*552f7358SJed Brown const PetscScalar g_012 = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 449*552f7358SJed Brown const PetscScalar h_012 = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 450*552f7358SJed Brown PetscScalar *ref, *real; 451*552f7358SJed Brown PetscErrorCode ierr; 452*552f7358SJed Brown 453*552f7358SJed Brown PetscFunctionBegin; 454*552f7358SJed Brown ierr = VecGetArray(Xref, &ref);CHKERRQ(ierr); 455*552f7358SJed Brown ierr = VecGetArray(Xreal, &real);CHKERRQ(ierr); 456*552f7358SJed Brown { 457*552f7358SJed Brown const PetscScalar p0 = ref[0]; 458*552f7358SJed Brown const PetscScalar p1 = ref[1]; 459*552f7358SJed Brown const PetscScalar p2 = ref[2]; 460*552f7358SJed Brown 461*552f7358SJed 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; 462*552f7358SJed 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; 463*552f7358SJed 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; 464*552f7358SJed Brown } 465*552f7358SJed Brown ierr = PetscLogFlops(114);CHKERRQ(ierr); 466*552f7358SJed Brown ierr = VecRestoreArray(Xref, &ref);CHKERRQ(ierr); 467*552f7358SJed Brown ierr = VecRestoreArray(Xreal, &real);CHKERRQ(ierr); 468*552f7358SJed Brown PetscFunctionReturn(0); 469*552f7358SJed Brown } 470*552f7358SJed Brown 471*552f7358SJed Brown #undef __FUNCT__ 472*552f7358SJed Brown #define __FUNCT__ "HexJacobian_Private" 473*552f7358SJed Brown PetscErrorCode HexJacobian_Private(SNES snes, Vec Xref, Mat *J, Mat *M, MatStructure *flag, void *ctx) 474*552f7358SJed Brown { 475*552f7358SJed Brown const PetscScalar*vertices = (const PetscScalar *) ctx; 476*552f7358SJed Brown const PetscScalar x0 = vertices[0]; 477*552f7358SJed Brown const PetscScalar y0 = vertices[1]; 478*552f7358SJed Brown const PetscScalar z0 = vertices[2]; 479*552f7358SJed Brown const PetscScalar x1 = vertices[3]; 480*552f7358SJed Brown const PetscScalar y1 = vertices[4]; 481*552f7358SJed Brown const PetscScalar z1 = vertices[5]; 482*552f7358SJed Brown const PetscScalar x2 = vertices[6]; 483*552f7358SJed Brown const PetscScalar y2 = vertices[7]; 484*552f7358SJed Brown const PetscScalar z2 = vertices[8]; 485*552f7358SJed Brown const PetscScalar x3 = vertices[9]; 486*552f7358SJed Brown const PetscScalar y3 = vertices[10]; 487*552f7358SJed Brown const PetscScalar z3 = vertices[11]; 488*552f7358SJed Brown const PetscScalar x4 = vertices[12]; 489*552f7358SJed Brown const PetscScalar y4 = vertices[13]; 490*552f7358SJed Brown const PetscScalar z4 = vertices[14]; 491*552f7358SJed Brown const PetscScalar x5 = vertices[15]; 492*552f7358SJed Brown const PetscScalar y5 = vertices[16]; 493*552f7358SJed Brown const PetscScalar z5 = vertices[17]; 494*552f7358SJed Brown const PetscScalar x6 = vertices[18]; 495*552f7358SJed Brown const PetscScalar y6 = vertices[19]; 496*552f7358SJed Brown const PetscScalar z6 = vertices[20]; 497*552f7358SJed Brown const PetscScalar x7 = vertices[21]; 498*552f7358SJed Brown const PetscScalar y7 = vertices[22]; 499*552f7358SJed Brown const PetscScalar z7 = vertices[23]; 500*552f7358SJed Brown const PetscScalar f_xy = x2 - x1 - x3 + x0; 501*552f7358SJed Brown const PetscScalar g_xy = y2 - y1 - y3 + y0; 502*552f7358SJed Brown const PetscScalar h_xy = z2 - z1 - z3 + z0; 503*552f7358SJed Brown const PetscScalar f_yz = x7 - x3 - x4 + x0; 504*552f7358SJed Brown const PetscScalar g_yz = y7 - y3 - y4 + y0; 505*552f7358SJed Brown const PetscScalar h_yz = z7 - z3 - z4 + z0; 506*552f7358SJed Brown const PetscScalar f_xz = x5 - x1 - x4 + x0; 507*552f7358SJed Brown const PetscScalar g_xz = y5 - y1 - y4 + y0; 508*552f7358SJed Brown const PetscScalar h_xz = z5 - z1 - z4 + z0; 509*552f7358SJed Brown const PetscScalar f_xyz = x6 - x0 + x1 - x2 + x3 + x4 - x5 - x7; 510*552f7358SJed Brown const PetscScalar g_xyz = y6 - y0 + y1 - y2 + y3 + y4 - y5 - y7; 511*552f7358SJed Brown const PetscScalar h_xyz = z6 - z0 + z1 - z2 + z3 + z4 - z5 - z7; 512*552f7358SJed Brown PetscScalar *ref; 513*552f7358SJed Brown PetscErrorCode ierr; 514*552f7358SJed Brown 515*552f7358SJed Brown PetscFunctionBegin; 516*552f7358SJed Brown ierr = VecGetArray(Xref, &ref);CHKERRQ(ierr); 517*552f7358SJed Brown { 518*552f7358SJed Brown const PetscScalar x = ref[0]; 519*552f7358SJed Brown const PetscScalar y = ref[1]; 520*552f7358SJed Brown const PetscScalar z = ref[2]; 521*552f7358SJed Brown const PetscInt rows[3] = {0, 1, 2}; 522*552f7358SJed Brown const PetscScalar values[9] = { 523*552f7358SJed Brown (x1 - x0 + f_xy*y + f_xz*z + f_xyz*y*z) / 2.0, 524*552f7358SJed Brown (x3 - x0 + f_xy*x + f_yz*z + f_xyz*x*z) / 2.0, 525*552f7358SJed Brown (x4 - x0 + f_yz*y + f_xz*x + f_xyz*x*y) / 2.0, 526*552f7358SJed Brown (y1 - y0 + g_xy*y + g_xz*z + g_xyz*y*z) / 2.0, 527*552f7358SJed Brown (y3 - y0 + g_xy*x + g_yz*z + g_xyz*x*z) / 2.0, 528*552f7358SJed Brown (y4 - y0 + g_yz*y + g_xz*x + g_xyz*x*y) / 2.0, 529*552f7358SJed Brown (z1 - z0 + h_xy*y + h_xz*z + h_xyz*y*z) / 2.0, 530*552f7358SJed Brown (z3 - z0 + h_xy*x + h_yz*z + h_xyz*x*z) / 2.0, 531*552f7358SJed Brown (z4 - z0 + h_yz*y + h_xz*x + h_xyz*x*y) / 2.0}; 532*552f7358SJed Brown ierr = MatSetValues(*J, 3, rows, 3, rows, values, INSERT_VALUES);CHKERRQ(ierr); 533*552f7358SJed Brown } 534*552f7358SJed Brown ierr = PetscLogFlops(152);CHKERRQ(ierr); 535*552f7358SJed Brown ierr = VecRestoreArray(Xref, &ref);CHKERRQ(ierr); 536*552f7358SJed Brown ierr = MatAssemblyBegin(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 537*552f7358SJed Brown ierr = MatAssemblyEnd(*J, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 538*552f7358SJed Brown PetscFunctionReturn(0); 539*552f7358SJed Brown } 540*552f7358SJed Brown 541*552f7358SJed Brown #undef __FUNCT__ 542*552f7358SJed Brown #define __FUNCT__ "DMInterpolate_Hex_Private" 543*552f7358SJed Brown PetscErrorCode DMInterpolate_Hex_Private(DMInterpolationInfo ctx, DM dm, Vec xLocal, Vec v) { 544*552f7358SJed Brown SNES snes; 545*552f7358SJed Brown KSP ksp; 546*552f7358SJed Brown PC pc; 547*552f7358SJed Brown Vec coordsLocal, r, ref, real; 548*552f7358SJed Brown Mat J; 549*552f7358SJed Brown PetscScalar *a, *coords; 550*552f7358SJed Brown PetscInt p; 551*552f7358SJed Brown PetscErrorCode ierr; 552*552f7358SJed Brown 553*552f7358SJed Brown PetscFunctionBegin; 554*552f7358SJed Brown ierr = DMGetCoordinatesLocal(dm, &coordsLocal);CHKERRQ(ierr); 555*552f7358SJed Brown ierr = SNESCreate(PETSC_COMM_SELF, &snes);CHKERRQ(ierr); 556*552f7358SJed Brown ierr = SNESSetOptionsPrefix(snes, "hex_interp_");CHKERRQ(ierr); 557*552f7358SJed Brown ierr = VecCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr); 558*552f7358SJed Brown ierr = VecSetSizes(r, 3, 3);CHKERRQ(ierr); 559*552f7358SJed Brown ierr = VecSetFromOptions(r);CHKERRQ(ierr); 560*552f7358SJed Brown ierr = VecDuplicate(r, &ref);CHKERRQ(ierr); 561*552f7358SJed Brown ierr = VecDuplicate(r, &real);CHKERRQ(ierr); 562*552f7358SJed Brown ierr = MatCreate(PETSC_COMM_SELF, &J);CHKERRQ(ierr); 563*552f7358SJed Brown ierr = MatSetSizes(J, 3, 3, 3, 3);CHKERRQ(ierr); 564*552f7358SJed Brown ierr = MatSetType(J, MATSEQDENSE);CHKERRQ(ierr); 565*552f7358SJed Brown ierr = MatSetUp(J);CHKERRQ(ierr); 566*552f7358SJed Brown ierr = SNESSetFunction(snes, r, HexMap_Private, PETSC_NULL);CHKERRQ(ierr); 567*552f7358SJed Brown ierr = SNESSetJacobian(snes, J, J, HexJacobian_Private, PETSC_NULL);CHKERRQ(ierr); 568*552f7358SJed Brown ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 569*552f7358SJed Brown ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 570*552f7358SJed Brown ierr = PCSetType(pc, PCLU);CHKERRQ(ierr); 571*552f7358SJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 572*552f7358SJed Brown 573*552f7358SJed Brown ierr = VecGetArray(ctx->coords, &coords);CHKERRQ(ierr); 574*552f7358SJed Brown ierr = VecGetArray(v, &a);CHKERRQ(ierr); 575*552f7358SJed Brown for(p = 0; p < ctx->n; ++p) { 576*552f7358SJed Brown const PetscScalar *x, *vertices; 577*552f7358SJed Brown PetscScalar *xi; 578*552f7358SJed Brown PetscInt c = ctx->cells[p], comp, coordSize, xSize; 579*552f7358SJed Brown 580*552f7358SJed Brown /* Can make this do all points at once */ 581*552f7358SJed Brown ierr = DMPlexVecGetClosure(dm, PETSC_NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 582*552f7358SJed Brown if (8*3 != coordSize) {SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", coordSize, 8*3);} 583*552f7358SJed Brown ierr = DMPlexVecGetClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 584*552f7358SJed Brown if (8*ctx->dof != xSize) {SETERRQ2(ctx->comm, PETSC_ERR_ARG_SIZ, "Invalid closure size %d should be %d", xSize, 8*ctx->dof);} 585*552f7358SJed Brown ierr = SNESSetFunction(snes, PETSC_NULL, PETSC_NULL, vertices);CHKERRQ(ierr); 586*552f7358SJed Brown ierr = SNESSetJacobian(snes, PETSC_NULL, PETSC_NULL, PETSC_NULL, vertices);CHKERRQ(ierr); 587*552f7358SJed Brown ierr = VecGetArray(real, &xi);CHKERRQ(ierr); 588*552f7358SJed Brown xi[0] = coords[p*ctx->dim+0]; 589*552f7358SJed Brown xi[1] = coords[p*ctx->dim+1]; 590*552f7358SJed Brown xi[2] = coords[p*ctx->dim+2]; 591*552f7358SJed Brown ierr = VecRestoreArray(real, &xi);CHKERRQ(ierr); 592*552f7358SJed Brown ierr = SNESSolve(snes, real, ref);CHKERRQ(ierr); 593*552f7358SJed Brown ierr = VecGetArray(ref, &xi);CHKERRQ(ierr); 594*552f7358SJed Brown for (comp = 0; comp < ctx->dof; ++comp) { 595*552f7358SJed Brown a[p*ctx->dof+comp] = 596*552f7358SJed Brown x[0*ctx->dof+comp]*(1-xi[0])*(1-xi[1])*(1-xi[2]) + 597*552f7358SJed Brown x[1*ctx->dof+comp]* xi[0]*(1-xi[1])*(1-xi[2]) + 598*552f7358SJed Brown x[2*ctx->dof+comp]* xi[0]* xi[1]*(1-xi[2]) + 599*552f7358SJed Brown x[3*ctx->dof+comp]*(1-xi[0])* xi[1]*(1-xi[2]) + 600*552f7358SJed Brown x[4*ctx->dof+comp]*(1-xi[0])*(1-xi[1])* xi[2] + 601*552f7358SJed Brown x[5*ctx->dof+comp]* xi[0]*(1-xi[1])* xi[2] + 602*552f7358SJed Brown x[6*ctx->dof+comp]* xi[0]* xi[1]* xi[2] + 603*552f7358SJed Brown x[7*ctx->dof+comp]*(1-xi[0])* xi[1]* xi[2]; 604*552f7358SJed Brown } 605*552f7358SJed Brown ierr = VecRestoreArray(ref, &xi);CHKERRQ(ierr); 606*552f7358SJed Brown ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, coordsLocal, c, &coordSize, &vertices);CHKERRQ(ierr); 607*552f7358SJed Brown ierr = DMPlexVecRestoreClosure(dm, PETSC_NULL, xLocal, c, &xSize, &x);CHKERRQ(ierr); 608*552f7358SJed Brown } 609*552f7358SJed Brown ierr = VecRestoreArray(v, &a);CHKERRQ(ierr); 610*552f7358SJed Brown ierr = VecRestoreArray(ctx->coords, &coords);CHKERRQ(ierr); 611*552f7358SJed Brown 612*552f7358SJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 613*552f7358SJed Brown ierr = VecDestroy(&r);CHKERRQ(ierr); 614*552f7358SJed Brown ierr = VecDestroy(&ref);CHKERRQ(ierr); 615*552f7358SJed Brown ierr = VecDestroy(&real);CHKERRQ(ierr); 616*552f7358SJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 617*552f7358SJed Brown PetscFunctionReturn(0); 618*552f7358SJed Brown } 619*552f7358SJed Brown 620*552f7358SJed Brown #undef __FUNCT__ 621*552f7358SJed Brown #define __FUNCT__ "DMInterpolationEvaluate" 622*552f7358SJed Brown /* 623*552f7358SJed Brown Input Parameters: 624*552f7358SJed Brown + ctx - The DMInterpolationInfo context 625*552f7358SJed Brown . dm - The DM 626*552f7358SJed Brown - x - The local vector containing the field to be interpolated 627*552f7358SJed Brown 628*552f7358SJed Brown Output Parameters: 629*552f7358SJed Brown . v - The vector containing the interpolated values 630*552f7358SJed Brown */ 631*552f7358SJed Brown PetscErrorCode DMInterpolationEvaluate(DMInterpolationInfo ctx, DM dm, Vec x, Vec v) { 632*552f7358SJed Brown PetscInt dim, coneSize, n; 633*552f7358SJed Brown PetscErrorCode ierr; 634*552f7358SJed Brown 635*552f7358SJed Brown PetscFunctionBegin; 636*552f7358SJed Brown PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 637*552f7358SJed Brown PetscValidHeaderSpecific(x, VEC_CLASSID, 3); 638*552f7358SJed Brown PetscValidHeaderSpecific(v, VEC_CLASSID, 4); 639*552f7358SJed Brown ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr); 640*552f7358SJed Brown 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);} 641*552f7358SJed Brown if (n) { 642*552f7358SJed Brown ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 643*552f7358SJed Brown ierr = DMPlexGetConeSize(dm, ctx->cells[0], &coneSize);CHKERRQ(ierr); 644*552f7358SJed Brown if (dim == 2) { 645*552f7358SJed Brown if (coneSize == 3) { 646*552f7358SJed Brown ierr = DMInterpolate_Simplex_Private(ctx, dm, x, v);CHKERRQ(ierr); 647*552f7358SJed Brown } else if (coneSize == 4) { 648*552f7358SJed Brown ierr = DMInterpolate_Quad_Private(ctx, dm, x, v);CHKERRQ(ierr); 649*552f7358SJed Brown } else { 650*552f7358SJed Brown SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim); 651*552f7358SJed Brown } 652*552f7358SJed Brown } else if (dim == 3) { 653*552f7358SJed Brown if (coneSize == 4) { 654*552f7358SJed Brown ierr = DMInterpolate_Simplex_Private(ctx, dm, x, v);CHKERRQ(ierr); 655*552f7358SJed Brown } else { 656*552f7358SJed Brown ierr = DMInterpolate_Hex_Private(ctx, dm, x, v);CHKERRQ(ierr); 657*552f7358SJed Brown } 658*552f7358SJed Brown } else { 659*552f7358SJed Brown SETERRQ1(ctx->comm, PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dimension %d for point interpolation", dim); 660*552f7358SJed Brown } 661*552f7358SJed Brown } 662*552f7358SJed Brown PetscFunctionReturn(0); 663*552f7358SJed Brown } 664*552f7358SJed Brown 665*552f7358SJed Brown #undef __FUNCT__ 666*552f7358SJed Brown #define __FUNCT__ "DMInterpolationDestroy" 667*552f7358SJed Brown PetscErrorCode DMInterpolationDestroy(DMInterpolationInfo *ctx) { 668*552f7358SJed Brown PetscErrorCode ierr; 669*552f7358SJed Brown 670*552f7358SJed Brown PetscFunctionBegin; 671*552f7358SJed Brown PetscValidPointer(ctx, 2); 672*552f7358SJed Brown ierr = VecDestroy(&(*ctx)->coords);CHKERRQ(ierr); 673*552f7358SJed Brown ierr = PetscFree((*ctx)->points);CHKERRQ(ierr); 674*552f7358SJed Brown ierr = PetscFree((*ctx)->cells);CHKERRQ(ierr); 675*552f7358SJed Brown ierr = PetscFree(*ctx);CHKERRQ(ierr); 676*552f7358SJed Brown *ctx = PETSC_NULL; 677*552f7358SJed Brown PetscFunctionReturn(0); 678*552f7358SJed Brown } 679