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