1 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 #include <petsc-private/petscfeimpl.h> 4 #include <petsc-private/petscfvimpl.h> 5 6 #undef __FUNCT__ 7 #define __FUNCT__ "DMPlexGetScale" 8 PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale) 9 { 10 DM_Plex *mesh = (DM_Plex*) dm->data; 11 12 PetscFunctionBegin; 13 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 14 PetscValidPointer(scale, 3); 15 *scale = mesh->scale[unit]; 16 PetscFunctionReturn(0); 17 } 18 19 #undef __FUNCT__ 20 #define __FUNCT__ "DMPlexSetScale" 21 PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale) 22 { 23 DM_Plex *mesh = (DM_Plex*) dm->data; 24 25 PetscFunctionBegin; 26 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 27 mesh->scale[unit] = scale; 28 PetscFunctionReturn(0); 29 } 30 31 PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k) 32 { 33 switch (i) { 34 case 0: 35 switch (j) { 36 case 0: return 0; 37 case 1: 38 switch (k) { 39 case 0: return 0; 40 case 1: return 0; 41 case 2: return 1; 42 } 43 case 2: 44 switch (k) { 45 case 0: return 0; 46 case 1: return -1; 47 case 2: return 0; 48 } 49 } 50 case 1: 51 switch (j) { 52 case 0: 53 switch (k) { 54 case 0: return 0; 55 case 1: return 0; 56 case 2: return -1; 57 } 58 case 1: return 0; 59 case 2: 60 switch (k) { 61 case 0: return 1; 62 case 1: return 0; 63 case 2: return 0; 64 } 65 } 66 case 2: 67 switch (j) { 68 case 0: 69 switch (k) { 70 case 0: return 0; 71 case 1: return 1; 72 case 2: return 0; 73 } 74 case 1: 75 switch (k) { 76 case 0: return -1; 77 case 1: return 0; 78 case 2: return 0; 79 } 80 case 2: return 0; 81 } 82 } 83 return 0; 84 } 85 86 #undef __FUNCT__ 87 #define __FUNCT__ "DMPlexCreateRigidBody" 88 /*@C 89 DMPlexCreateRigidBody - create rigid body modes from coordinates 90 91 Collective on DM 92 93 Input Arguments: 94 + dm - the DM 95 . section - the local section associated with the rigid field, or NULL for the default section 96 - globalSection - the global section associated with the rigid field, or NULL for the default section 97 98 Output Argument: 99 . sp - the null space 100 101 Note: This is necessary to take account of Dirichlet conditions on the displacements 102 103 Level: advanced 104 105 .seealso: MatNullSpaceCreate() 106 @*/ 107 PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp) 108 { 109 MPI_Comm comm; 110 Vec coordinates, localMode, mode[6]; 111 PetscSection coordSection; 112 PetscScalar *coords; 113 PetscInt dim, vStart, vEnd, v, n, m, d, i, j; 114 PetscErrorCode ierr; 115 116 PetscFunctionBegin; 117 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 118 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 119 if (dim == 1) { 120 ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr); 121 PetscFunctionReturn(0); 122 } 123 if (!section) {ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr);} 124 if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);} 125 ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr); 126 ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 127 ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 128 ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 129 m = (dim*(dim+1))/2; 130 ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr); 131 ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr); 132 ierr = VecSetUp(mode[0]);CHKERRQ(ierr); 133 for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);} 134 /* Assume P1 */ 135 ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr); 136 for (d = 0; d < dim; ++d) { 137 PetscScalar values[3] = {0.0, 0.0, 0.0}; 138 139 values[d] = 1.0; 140 ierr = VecSet(localMode, 0.0);CHKERRQ(ierr); 141 for (v = vStart; v < vEnd; ++v) { 142 ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr); 143 } 144 ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 145 ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 146 } 147 ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 148 for (d = dim; d < dim*(dim+1)/2; ++d) { 149 PetscInt i, j, k = dim > 2 ? d - dim : d; 150 151 ierr = VecSet(localMode, 0.0);CHKERRQ(ierr); 152 for (v = vStart; v < vEnd; ++v) { 153 PetscScalar values[3] = {0.0, 0.0, 0.0}; 154 PetscInt off; 155 156 ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 157 for (i = 0; i < dim; ++i) { 158 for (j = 0; j < dim; ++j) { 159 values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]); 160 } 161 } 162 ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr); 163 } 164 ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 165 ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 166 } 167 ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 168 ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr); 169 for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);} 170 /* Orthonormalize system */ 171 for (i = dim; i < m; ++i) { 172 PetscScalar dots[6]; 173 174 ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr); 175 for (j = 0; j < i; ++j) dots[j] *= -1.0; 176 ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr); 177 ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr); 178 } 179 ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr); 180 for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);} 181 PetscFunctionReturn(0); 182 } 183 184 #undef __FUNCT__ 185 #define __FUNCT__ "DMPlexProjectFunctionLabelLocal" 186 PetscErrorCode DMPlexProjectFunctionLabelLocal(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 187 { 188 PetscFE fe; 189 PetscDualSpace *sp; 190 PetscSection section; 191 PetscScalar *values; 192 PetscBool *fieldActive; 193 PetscInt numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, f, d, v, i, comp; 194 PetscErrorCode ierr; 195 196 PetscFunctionBegin; 197 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 198 if (cEnd <= cStart) PetscFunctionReturn(0); 199 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 200 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 201 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 202 ierr = PetscMalloc1(numFields,&sp);CHKERRQ(ierr); 203 for (f = 0; f < numFields; ++f) { 204 ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr); 205 ierr = PetscFEGetDualSpace(fe, &sp[f]);CHKERRQ(ierr); 206 ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr); 207 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 208 totDim += spDim*numComp; 209 } 210 ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 211 if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 212 ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 213 ierr = DMGetWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr); 214 for (f = 0; f < numFields; ++f) fieldActive[f] = funcs[f] ? PETSC_TRUE : PETSC_FALSE; 215 for (i = 0; i < numIds; ++i) { 216 IS pointIS; 217 const PetscInt *points; 218 PetscInt n, p; 219 220 ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 221 ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr); 222 ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 223 for (p = 0; p < n; ++p) { 224 const PetscInt point = points[p]; 225 PetscFECellGeom geom; 226 227 if ((point < cStart) || (point >= cEnd)) continue; 228 ierr = DMPlexComputeCellGeometryFEM(dm, point, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr); 229 for (f = 0, v = 0; f < numFields; ++f) { 230 void * const ctx = ctxs ? ctxs[f] : NULL; 231 232 ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr); 233 ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr); 234 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 235 for (d = 0; d < spDim; ++d) { 236 if (funcs[f]) { 237 ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr); 238 } else { 239 for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0; 240 } 241 v += numComp; 242 } 243 } 244 ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr); 245 } 246 ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 247 ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 248 } 249 ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 250 ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr); 251 ierr = PetscFree(sp);CHKERRQ(ierr); 252 PetscFunctionReturn(0); 253 } 254 255 #undef __FUNCT__ 256 #define __FUNCT__ "DMPlexProjectFunctionLocal" 257 PetscErrorCode DMPlexProjectFunctionLocal(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 258 { 259 PetscDualSpace *sp; 260 PetscInt *numComp; 261 PetscSection section; 262 PetscScalar *values; 263 PetscInt numFields, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v, comp; 264 PetscErrorCode ierr; 265 266 PetscFunctionBegin; 267 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 268 if (cEnd <= cStart) PetscFunctionReturn(0); 269 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 270 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 271 ierr = PetscMalloc2(numFields, &sp, numFields, &numComp);CHKERRQ(ierr); 272 for (f = 0; f < numFields; ++f) { 273 PetscObject obj; 274 PetscClassId id; 275 276 ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr); 277 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 278 if (id == PETSCFE_CLASSID) { 279 PetscFE fe = (PetscFE) obj; 280 281 ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr); 282 ierr = PetscFEGetDualSpace(fe, &sp[f]);CHKERRQ(ierr); 283 ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr); 284 } else if (id == PETSCFV_CLASSID) { 285 PetscFV fv = (PetscFV) obj; 286 PetscQuadrature q; 287 288 ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr); 289 ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr); 290 ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fv), &sp[f]);CHKERRQ(ierr); 291 ierr = PetscDualSpaceSetDM(sp[f], dm);CHKERRQ(ierr); 292 ierr = PetscDualSpaceSetType(sp[f], PETSCDUALSPACESIMPLE);CHKERRQ(ierr); 293 ierr = PetscDualSpaceSimpleSetDimension(sp[f], 1);CHKERRQ(ierr); 294 ierr = PetscDualSpaceSimpleSetFunctional(sp[f], 0, q);CHKERRQ(ierr); 295 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 296 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 297 totDim += spDim*numComp[f]; 298 } 299 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 300 ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 301 if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 302 ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 303 for (c = cStart; c < cEnd; ++c) { 304 PetscFECellGeom geom; 305 306 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr); 307 for (f = 0, v = 0; f < numFields; ++f) { 308 void *const ctx = ctxs ? ctxs[f] : NULL; 309 310 ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 311 for (d = 0; d < spDim; ++d) { 312 if (funcs[f]) { 313 ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr); 314 } else { 315 for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0; 316 } 317 v += numComp[f]; 318 } 319 } 320 ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr); 321 } 322 ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 323 for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);} 324 ierr = PetscFree2(sp, numComp);CHKERRQ(ierr); 325 PetscFunctionReturn(0); 326 } 327 328 #undef __FUNCT__ 329 #define __FUNCT__ "DMPlexProjectFieldLocal" 330 PetscErrorCode DMPlexProjectFieldLocal(DM dm, Vec localU, void (**funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode mode, Vec localX) 331 { 332 DM dmAux; 333 PetscDS prob, probAux; 334 Vec A; 335 PetscSection section, sectionAux; 336 PetscScalar *values, *u, *u_x, *a, *a_x; 337 PetscReal *x, *v0, *J, *invJ, detJ, **basisField, **basisFieldDer, **basisFieldAux, **basisFieldDerAux; 338 PetscInt Nf, dim, spDim, totDim, numValues, cStart, cEnd, c, f, d, v, comp; 339 PetscErrorCode ierr; 340 341 PetscFunctionBegin; 342 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 343 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 344 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 345 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 346 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 347 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 348 ierr = PetscDSGetTabulation(prob, &basisField, &basisFieldDer);CHKERRQ(ierr); 349 ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr); 350 ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr); 351 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 352 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 353 if (dmAux) { 354 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 355 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 356 ierr = PetscDSGetTabulation(prob, &basisFieldAux, &basisFieldDerAux);CHKERRQ(ierr); 357 ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr); 358 } 359 ierr = DMPlexInsertBoundaryValues(dm, localU, 0.0, NULL, NULL, NULL);CHKERRQ(ierr); 360 ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 361 if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 362 ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 363 ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 364 for (c = cStart; c < cEnd; ++c) { 365 PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 366 367 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 368 ierr = DMPlexVecGetClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr); 369 if (dmAux) {ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);} 370 for (f = 0, v = 0; f < Nf; ++f) { 371 PetscFE fe; 372 PetscDualSpace sp; 373 PetscInt Ncf; 374 375 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 376 ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr); 377 ierr = PetscFEGetNumComponents(fe, &Ncf);CHKERRQ(ierr); 378 ierr = PetscDualSpaceGetDimension(sp, &spDim);CHKERRQ(ierr); 379 for (d = 0; d < spDim; ++d) { 380 PetscQuadrature quad; 381 const PetscReal *points, *weights; 382 PetscInt numPoints, q; 383 384 if (funcs[f]) { 385 ierr = PetscDualSpaceGetFunctional(sp, d, &quad);CHKERRQ(ierr); 386 ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr); 387 ierr = PetscFEGetTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr); 388 for (q = 0; q < numPoints; ++q) { 389 CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x); 390 ierr = EvaluateFieldJets(prob, PETSC_FALSE, q, invJ, coefficients, NULL, u, u_x, NULL);CHKERRQ(ierr); 391 ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, coefficientsAux, NULL, a, a_x, NULL);CHKERRQ(ierr); 392 (*funcs[f])(u, NULL, u_x, a, NULL, a_x, x, &values[v]); 393 } 394 ierr = PetscFERestoreTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr); 395 } else { 396 for (comp = 0; comp < Ncf; ++comp) values[v+comp] = 0.0; 397 } 398 v += Ncf; 399 } 400 } 401 ierr = DMPlexVecRestoreClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr); 402 if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);} 403 ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr); 404 } 405 ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 406 ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr); 407 PetscFunctionReturn(0); 408 } 409 410 #undef __FUNCT__ 411 #define __FUNCT__ "DMPlexInsertBoundaryValues_FEM_Internal" 412 static PetscErrorCode DMPlexInsertBoundaryValues_FEM_Internal(DM dm, PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], void (*func)(const PetscReal[], PetscScalar *, void *), void *ctx, Vec locX) 413 { 414 void (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx); 415 void **ctxs; 416 PetscInt numFields; 417 PetscErrorCode ierr; 418 419 PetscFunctionBegin; 420 ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr); 421 ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr); 422 funcs[field] = func; 423 ctxs[field] = ctx; 424 ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr); 425 ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr); 426 PetscFunctionReturn(0); 427 } 428 429 #undef __FUNCT__ 430 #define __FUNCT__ "DMPlexInsertBoundaryValues_FVM_Internal" 431 static PetscErrorCode DMPlexInsertBoundaryValues_FVM_Internal(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad, 432 PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *ctx, Vec locX) 433 { 434 PetscFV fv; 435 DM dmFace, dmCell, dmGrad; 436 const PetscScalar *facegeom, *cellgeom, *grad; 437 PetscScalar *x, *fx; 438 PetscInt dim, fStart, fEnd, pdim, i; 439 PetscErrorCode ierr; 440 441 PetscFunctionBegin; 442 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 443 ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 444 ierr = DMGetField(dm, field, (PetscObject *) &fv);CHKERRQ(ierr); 445 ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 446 ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 447 ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 448 ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 449 if (Grad) { 450 ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr); 451 ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr); 452 ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr); 453 ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr); 454 } 455 ierr = VecGetArray(locX, &x);CHKERRQ(ierr); 456 for (i = 0; i < numids; ++i) { 457 IS faceIS; 458 const PetscInt *faces; 459 PetscInt numFaces, f; 460 461 ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr); 462 if (!faceIS) continue; /* No points with that id on this process */ 463 ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr); 464 ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr); 465 for (f = 0; f < numFaces; ++f) { 466 const PetscInt face = faces[f], *cells; 467 const PetscFVFaceGeom *fg; 468 469 if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */ 470 ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 471 ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 472 if (Grad) { 473 const PetscFVCellGeom *cg; 474 const PetscScalar *cx, *cgrad; 475 PetscScalar *xG; 476 PetscReal dx[3]; 477 PetscInt d; 478 479 ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr); 480 ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr); 481 ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr); 482 ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr); 483 DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx); 484 for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx); 485 ierr = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx);CHKERRQ(ierr); 486 } else { 487 const PetscScalar *xI; 488 PetscScalar *xG; 489 490 ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr); 491 ierr = DMPlexPointLocalRef(dm, cells[1], x, &xG);CHKERRQ(ierr); 492 ierr = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx);CHKERRQ(ierr); 493 } 494 } 495 ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 496 ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 497 } 498 ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr); 499 if (Grad) { 500 ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr); 501 ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr); 502 } 503 ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 504 ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 505 PetscFunctionReturn(0); 506 } 507 508 #undef __FUNCT__ 509 #define __FUNCT__ "DMPlexInsertBoundaryValues" 510 PetscErrorCode DMPlexInsertBoundaryValues(DM dm, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM) 511 { 512 PetscInt numBd, b; 513 PetscErrorCode ierr; 514 515 PetscFunctionBegin; 516 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 517 PetscValidHeaderSpecific(locX, VEC_CLASSID, 2); 518 if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);} 519 if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);} 520 if (gradFVM) {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);} 521 ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 522 for (b = 0; b < numBd; ++b) { 523 PetscBool isEssential; 524 const char *labelname; 525 DMLabel label; 526 PetscInt field; 527 PetscObject obj; 528 PetscClassId id; 529 void (*func)(); 530 PetscInt numids; 531 const PetscInt *ids; 532 void *ctx; 533 534 ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr); 535 if (!isEssential) continue; 536 ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr); 537 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 538 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 539 if (id == PETSCFE_CLASSID) { 540 ierr = DMPlexInsertBoundaryValues_FEM_Internal(dm, field, label, numids, ids, (void (*)(const PetscReal[], PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr); 541 } else if (id == PETSCFV_CLASSID) { 542 ierr = DMPlexInsertBoundaryValues_FVM_Internal(dm, time, faceGeomFVM, cellGeomFVM, gradFVM, 543 field, label, numids, ids, (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr); 544 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 545 } 546 PetscFunctionReturn(0); 547 } 548 549 #undef __FUNCT__ 550 #define __FUNCT__ "DMPlexComputeL2Diff" 551 /*@C 552 DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 553 554 Input Parameters: 555 + dm - The DM 556 . funcs - The functions to evaluate for each field component 557 . ctxs - Optional array of contexts to pass to each function, or NULL. 558 - X - The coefficient vector u_h 559 560 Output Parameter: 561 . diff - The diff ||u - u_h||_2 562 563 Level: developer 564 565 .seealso: DMPlexProjectFunction(), DMPlexComputeL2FieldDiff(), DMPlexComputeL2GradientDiff() 566 @*/ 567 PetscErrorCode DMPlexComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff) 568 { 569 const PetscInt debug = 0; 570 PetscSection section; 571 PetscQuadrature quad; 572 Vec localX; 573 PetscScalar *funcVal, *interpolant; 574 PetscReal *coords, *v0, *J, *invJ, detJ; 575 PetscReal localDiff = 0.0; 576 const PetscReal *quadPoints, *quadWeights; 577 PetscInt dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, c, field, fieldOffset; 578 PetscErrorCode ierr; 579 580 PetscFunctionBegin; 581 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 582 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 583 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 584 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 585 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 586 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 587 for (field = 0; field < numFields; ++field) { 588 PetscObject obj; 589 PetscClassId id; 590 PetscInt Nc; 591 592 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 593 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 594 if (id == PETSCFE_CLASSID) { 595 PetscFE fe = (PetscFE) obj; 596 597 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 598 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 599 } else if (id == PETSCFV_CLASSID) { 600 PetscFV fv = (PetscFV) obj; 601 602 ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 603 ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 604 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 605 numComponents += Nc; 606 } 607 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 608 ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 609 ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 610 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 611 for (c = cStart; c < cEnd; ++c) { 612 PetscScalar *x = NULL; 613 PetscReal elemDiff = 0.0; 614 615 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 616 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 617 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 618 619 for (field = 0, fieldOffset = 0; field < numFields; ++field) { 620 PetscObject obj; 621 PetscClassId id; 622 void * const ctx = ctxs ? ctxs[field] : NULL; 623 PetscInt Nb, Nc, q, fc; 624 625 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 626 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 627 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 628 else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 629 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 630 if (debug) { 631 char title[1024]; 632 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 633 ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr); 634 } 635 for (q = 0; q < numQuadPoints; ++q) { 636 CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords); 637 (*funcs[field])(coords, funcVal, ctx); 638 if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 639 else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 640 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 641 for (fc = 0; fc < Nc; ++fc) { 642 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);} 643 elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ; 644 } 645 } 646 fieldOffset += Nb*Nc; 647 } 648 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 649 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 650 localDiff += elemDiff; 651 } 652 ierr = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr); 653 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 654 ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 655 *diff = PetscSqrtReal(*diff); 656 PetscFunctionReturn(0); 657 } 658 659 #undef __FUNCT__ 660 #define __FUNCT__ "DMPlexComputeL2GradientDiff" 661 /*@C 662 DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h. 663 664 Input Parameters: 665 + dm - The DM 666 . funcs - The gradient functions to evaluate for each field component 667 . ctxs - Optional array of contexts to pass to each function, or NULL. 668 . X - The coefficient vector u_h 669 - n - The vector to project along 670 671 Output Parameter: 672 . diff - The diff ||(grad u - grad u_h) . n||_2 673 674 Level: developer 675 676 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff() 677 @*/ 678 PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff) 679 { 680 const PetscInt debug = 0; 681 PetscSection section; 682 PetscQuadrature quad; 683 Vec localX; 684 PetscScalar *funcVal, *interpolantVec; 685 PetscReal *coords, *realSpaceDer, *v0, *J, *invJ, detJ; 686 PetscReal localDiff = 0.0; 687 PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 688 PetscErrorCode ierr; 689 690 PetscFunctionBegin; 691 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 692 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 693 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 694 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 695 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 696 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 697 for (field = 0; field < numFields; ++field) { 698 PetscFE fe; 699 PetscInt Nc; 700 701 ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 702 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 703 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 704 numComponents += Nc; 705 } 706 /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */ 707 ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr); 708 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 709 for (c = cStart; c < cEnd; ++c) { 710 PetscScalar *x = NULL; 711 PetscReal elemDiff = 0.0; 712 713 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 714 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 715 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 716 717 for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 718 PetscFE fe; 719 void * const ctx = ctxs ? ctxs[field] : NULL; 720 const PetscReal *quadPoints, *quadWeights; 721 PetscReal *basisDer; 722 PetscInt numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g; 723 724 ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 725 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 726 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 727 ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr); 728 ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr); 729 if (debug) { 730 char title[1024]; 731 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 732 ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr); 733 } 734 for (q = 0; q < numQuadPoints; ++q) { 735 for (d = 0; d < dim; d++) { 736 coords[d] = v0[d]; 737 for (e = 0; e < dim; e++) { 738 coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 739 } 740 } 741 (*funcs[field])(coords, n, funcVal, ctx); 742 for (fc = 0; fc < Ncomp; ++fc) { 743 PetscScalar interpolant = 0.0; 744 745 for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0; 746 for (f = 0; f < Nb; ++f) { 747 const PetscInt fidx = f*Ncomp+fc; 748 749 for (d = 0; d < dim; ++d) { 750 realSpaceDer[d] = 0.0; 751 for (g = 0; g < dim; ++g) { 752 realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g]; 753 } 754 interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d]; 755 } 756 } 757 for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d]; 758 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d fieldDer %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);} 759 elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 760 } 761 } 762 comp += Ncomp; 763 fieldOffset += Nb*Ncomp; 764 } 765 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 766 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 767 localDiff += elemDiff; 768 } 769 ierr = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr); 770 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 771 ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 772 *diff = PetscSqrtReal(*diff); 773 PetscFunctionReturn(0); 774 } 775 776 #undef __FUNCT__ 777 #define __FUNCT__ "DMPlexComputeL2FieldDiff" 778 /*@C 779 DMPlexComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components. 780 781 Input Parameters: 782 + dm - The DM 783 . funcs - The functions to evaluate for each field component 784 . ctxs - Optional array of contexts to pass to each function, or NULL. 785 - X - The coefficient vector u_h 786 787 Output Parameter: 788 . diff - The array of differences, ||u^f - u^f_h||_2 789 790 Level: developer 791 792 .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff(), DMPlexComputeL2GradientDiff() 793 @*/ 794 PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[]) 795 { 796 const PetscInt debug = 0; 797 PetscSection section; 798 PetscQuadrature quad; 799 Vec localX; 800 PetscScalar *funcVal, *interpolant; 801 PetscReal *coords, *v0, *J, *invJ, detJ; 802 PetscReal *localDiff; 803 const PetscReal *quadPoints, *quadWeights; 804 PetscInt dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, c, field, fieldOffset; 805 PetscErrorCode ierr; 806 807 PetscFunctionBegin; 808 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 809 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 810 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 811 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 812 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 813 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 814 for (field = 0; field < numFields; ++field) { 815 PetscObject obj; 816 PetscClassId id; 817 PetscInt Nc; 818 819 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 820 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 821 if (id == PETSCFE_CLASSID) { 822 PetscFE fe = (PetscFE) obj; 823 824 ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 825 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 826 } else if (id == PETSCFV_CLASSID) { 827 PetscFV fv = (PetscFV) obj; 828 829 ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 830 ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 831 } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 832 numComponents += Nc; 833 } 834 ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 835 ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 836 ierr = PetscCalloc7(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 837 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 838 for (c = cStart; c < cEnd; ++c) { 839 PetscScalar *x = NULL; 840 841 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 842 if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 843 ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 844 845 for (field = 0, fieldOffset = 0; field < numFields; ++field) { 846 PetscObject obj; 847 PetscClassId id; 848 void * const ctx = ctxs ? ctxs[field] : NULL; 849 PetscInt Nb, Nc, q, fc; 850 851 PetscReal elemDiff = 0.0; 852 853 ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 854 ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 855 if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 856 else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 857 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 858 if (debug) { 859 char title[1024]; 860 ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 861 ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr); 862 } 863 for (q = 0; q < numQuadPoints; ++q) { 864 CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords); 865 (*funcs[field])(coords, funcVal, ctx); 866 if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 867 else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 868 else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 869 for (fc = 0; fc < Nc; ++fc) { 870 if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);} 871 elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ; 872 } 873 } 874 fieldOffset += Nb*Nc; 875 localDiff[field] += elemDiff; 876 } 877 ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 878 } 879 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 880 ierr = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 881 for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]); 882 ierr = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr); 883 PetscFunctionReturn(0); 884 } 885 886 #undef __FUNCT__ 887 #define __FUNCT__ "DMPlexComputeIntegralFEM" 888 /*@ 889 DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user 890 891 Input Parameters: 892 + dm - The mesh 893 . X - Local input vector 894 - user - The user context 895 896 Output Parameter: 897 . integral - Local integral for each field 898 899 Level: developer 900 901 .seealso: DMPlexComputeResidualFEM() 902 @*/ 903 PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user) 904 { 905 DM_Plex *mesh = (DM_Plex *) dm->data; 906 DM dmAux; 907 Vec localX, A; 908 PetscDS prob, probAux = NULL; 909 PetscQuadrature q; 910 PetscSection section, sectionAux; 911 PetscFECellGeom *cgeom; 912 PetscScalar *u, *a = NULL; 913 PetscInt dim, Nf, f, numCells, cStart, cEnd, c; 914 PetscInt totDim, totDimAux; 915 PetscErrorCode ierr; 916 917 PetscFunctionBegin; 918 /*ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/ 919 ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 920 ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 921 ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 922 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 923 ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 924 ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 925 ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 926 ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 927 ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 928 numCells = cEnd - cStart; 929 for (f = 0; f < Nf; ++f) {integral[f] = 0.0;} 930 ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 931 ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 932 if (dmAux) { 933 ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 934 ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 935 ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 936 } 937 ierr = DMPlexInsertBoundaryValues(dm, localX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr); 938 ierr = PetscMalloc2(numCells*totDim,&u,numCells,&cgeom);CHKERRQ(ierr); 939 if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);} 940 for (c = cStart; c < cEnd; ++c) { 941 PetscScalar *x = NULL; 942 PetscInt i; 943 944 ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeom[c].v0, cgeom[c].J, cgeom[c].invJ, &cgeom[c].detJ);CHKERRQ(ierr); 945 if (cgeom[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeom[c].detJ, c); 946 ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr); 947 for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i]; 948 ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr); 949 if (dmAux) { 950 ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 951 for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i]; 952 ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 953 } 954 } 955 for (f = 0; f < Nf; ++f) { 956 PetscFE fe; 957 PetscInt numQuadPoints, Nb; 958 /* Conforming batches */ 959 PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 960 /* Remainder */ 961 PetscInt Nr, offset; 962 963 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 964 ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 965 ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 966 ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 967 ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 968 blockSize = Nb*numQuadPoints; 969 batchSize = numBlocks * blockSize; 970 ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 971 numChunks = numCells / (numBatches*batchSize); 972 Ne = numChunks*numBatches*batchSize; 973 Nr = numCells % (numBatches*batchSize); 974 offset = numCells - Nr; 975 ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeom, u, probAux, a, integral);CHKERRQ(ierr); 976 ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeom[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], integral);CHKERRQ(ierr); 977 } 978 ierr = PetscFree2(u,cgeom);CHKERRQ(ierr); 979 if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 980 if (mesh->printFEM) { 981 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr); 982 for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", integral[f]);CHKERRQ(ierr);} 983 ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr); 984 } 985 ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 986 /* TODO: Allreduce for integral */ 987 /*ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/ 988 PetscFunctionReturn(0); 989 } 990 991 #undef __FUNCT__ 992 #define __FUNCT__ "DMPlexComputeInterpolatorFEM" 993 /*@ 994 DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM. 995 996 Input Parameters: 997 + dmf - The fine mesh 998 . dmc - The coarse mesh 999 - user - The user context 1000 1001 Output Parameter: 1002 . In - The interpolation matrix 1003 1004 Note: 1005 The first member of the user context must be an FEMContext. 1006 1007 We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 1008 like a GPU, or vectorize on a multicore machine. 1009 1010 Level: developer 1011 1012 .seealso: DMPlexComputeJacobianFEM() 1013 @*/ 1014 PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user) 1015 { 1016 DM_Plex *mesh = (DM_Plex *) dmc->data; 1017 const char *name = "Interpolator"; 1018 PetscDS prob; 1019 PetscFE *feRef; 1020 PetscSection fsection, fglobalSection; 1021 PetscSection csection, cglobalSection; 1022 PetscScalar *elemMat; 1023 PetscInt dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c; 1024 PetscInt cTotDim, rTotDim = 0; 1025 PetscErrorCode ierr; 1026 1027 PetscFunctionBegin; 1028 ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1029 ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr); 1030 ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 1031 ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 1032 ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 1033 ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 1034 ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 1035 ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 1036 ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr); 1037 ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr); 1038 for (f = 0; f < Nf; ++f) { 1039 PetscFE fe; 1040 PetscInt rNb, Nc; 1041 1042 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1043 ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr); 1044 ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr); 1045 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1046 rTotDim += rNb*Nc; 1047 } 1048 ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr); 1049 ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr); 1050 ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr); 1051 for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) { 1052 PetscDualSpace Qref; 1053 PetscQuadrature f; 1054 const PetscReal *qpoints, *qweights; 1055 PetscReal *points; 1056 PetscInt npoints = 0, Nc, Np, fpdim, i, k, p, d; 1057 1058 /* Compose points from all dual basis functionals */ 1059 ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr); 1060 ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr); 1061 ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr); 1062 for (i = 0; i < fpdim; ++i) { 1063 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1064 ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr); 1065 npoints += Np; 1066 } 1067 ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr); 1068 for (i = 0, k = 0; i < fpdim; ++i) { 1069 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1070 ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr); 1071 for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d]; 1072 } 1073 1074 for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) { 1075 PetscFE fe; 1076 PetscReal *B; 1077 PetscInt NcJ, cpdim, j; 1078 1079 /* Evaluate basis at points */ 1080 ierr = PetscDSGetDiscretization(prob, fieldJ, (PetscObject *) &fe);CHKERRQ(ierr); 1081 ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr); 1082 ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr); 1083 /* For now, fields only interpolate themselves */ 1084 if (fieldI == fieldJ) { 1085 if (Nc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %d does not match coarse field %d", Nc, NcJ); 1086 ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr); 1087 for (i = 0, k = 0; i < fpdim; ++i) { 1088 ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1089 ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr); 1090 for (p = 0; p < Np; ++p, ++k) { 1091 for (j = 0; j < cpdim; ++j) { 1092 for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cTotDim + offsetJ + j*NcJ + c] += B[k*cpdim*NcJ+j*Nc+c]*qweights[p]; 1093 } 1094 } 1095 } 1096 ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 1097 } 1098 offsetJ += cpdim*NcJ; 1099 } 1100 offsetI += fpdim*Nc; 1101 ierr = PetscFree(points);CHKERRQ(ierr); 1102 } 1103 if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);} 1104 /* Preallocate matrix */ 1105 { 1106 PetscHashJK ht; 1107 PetscLayout rLayout; 1108 PetscInt *dnz, *onz, *cellCIndices, *cellFIndices; 1109 PetscInt locRows, rStart, rEnd, cell, r; 1110 1111 ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr); 1112 ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr); 1113 ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr); 1114 ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr); 1115 ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr); 1116 ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr); 1117 ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr); 1118 ierr = PetscCalloc4(locRows,&dnz,locRows,&onz,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr); 1119 ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr); 1120 for (cell = cStart; cell < cEnd; ++cell) { 1121 ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr); 1122 for (r = 0; r < rTotDim; ++r) { 1123 PetscHashJKKey key; 1124 PetscHashJKIter missing, iter; 1125 1126 key.j = cellFIndices[r]; 1127 if (key.j < 0) continue; 1128 for (c = 0; c < cTotDim; ++c) { 1129 key.k = cellCIndices[c]; 1130 if (key.k < 0) continue; 1131 ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr); 1132 if (missing) { 1133 ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr); 1134 if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart]; 1135 else ++onz[key.j-rStart]; 1136 } 1137 } 1138 } 1139 } 1140 ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr); 1141 ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr); 1142 ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 1143 ierr = PetscFree4(dnz,onz,cellCIndices,cellFIndices);CHKERRQ(ierr); 1144 } 1145 /* Fill matrix */ 1146 ierr = MatZeroEntries(In);CHKERRQ(ierr); 1147 for (c = cStart; c < cEnd; ++c) { 1148 ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr); 1149 } 1150 for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);} 1151 ierr = PetscFree(feRef);CHKERRQ(ierr); 1152 ierr = PetscFree(elemMat);CHKERRQ(ierr); 1153 ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1154 ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1155 if (mesh->printFEM) { 1156 ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 1157 ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr); 1158 ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1159 } 1160 ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1161 PetscFunctionReturn(0); 1162 } 1163 1164 #undef __FUNCT__ 1165 #define __FUNCT__ "DMPlexComputeInjectorFEM" 1166 PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user) 1167 { 1168 PetscDS prob; 1169 PetscFE *feRef; 1170 Vec fv, cv; 1171 IS fis, cis; 1172 PetscSection fsection, fglobalSection, csection, cglobalSection; 1173 PetscInt *cmap, *cellCIndices, *cellFIndices, *cindices, *findices; 1174 PetscInt cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, c, dim, d, startC, offsetC, offsetF, m; 1175 PetscErrorCode ierr; 1176 1177 PetscFunctionBegin; 1178 ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1179 ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr); 1180 ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 1181 ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 1182 ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 1183 ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 1184 ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 1185 ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 1186 ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr); 1187 ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr); 1188 for (f = 0; f < Nf; ++f) { 1189 PetscFE fe; 1190 PetscInt fNb, Nc; 1191 1192 ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 1193 ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr); 1194 ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr); 1195 ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 1196 fTotDim += fNb*Nc; 1197 } 1198 ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr); 1199 ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr); 1200 for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) { 1201 PetscFE feC; 1202 PetscDualSpace QF, QC; 1203 PetscInt NcF, NcC, fpdim, cpdim; 1204 1205 ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr); 1206 ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr); 1207 ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr); 1208 if (NcF != NcC) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in fine space field %d does not match coarse field %d", NcF, NcC); 1209 ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr); 1210 ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr); 1211 ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr); 1212 ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr); 1213 for (c = 0; c < cpdim; ++c) { 1214 PetscQuadrature cfunc; 1215 const PetscReal *cqpoints; 1216 PetscInt NpC; 1217 1218 ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr); 1219 ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr); 1220 if (NpC != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments"); 1221 for (f = 0; f < fpdim; ++f) { 1222 PetscQuadrature ffunc; 1223 const PetscReal *fqpoints; 1224 PetscReal sum = 0.0; 1225 PetscInt NpF, comp; 1226 1227 ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr); 1228 ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr); 1229 if (NpC != NpF) continue; 1230 for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]); 1231 if (sum > 1.0e-9) continue; 1232 for (comp = 0; comp < NcC; ++comp) { 1233 cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp; 1234 } 1235 break; 1236 } 1237 } 1238 offsetC += cpdim*NcC; 1239 offsetF += fpdim*NcF; 1240 } 1241 for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);} 1242 ierr = PetscFree(feRef);CHKERRQ(ierr); 1243 1244 ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr); 1245 ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr); 1246 ierr = VecGetOwnershipRange(cv, &startC, NULL);CHKERRQ(ierr); 1247 ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr); 1248 ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr); 1249 ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr); 1250 ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr); 1251 for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1; 1252 for (c = cStart; c < cEnd; ++c) { 1253 ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr); 1254 for (d = 0; d < cTotDim; ++d) { 1255 if (cellCIndices[d] < 0) continue; 1256 if ((findices[cellCIndices[d]-startC] >= 0) && (findices[cellCIndices[d]-startC] != cellFIndices[cmap[d]])) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Coarse dof %d maps to both %d and %d", cindices[cellCIndices[d]-startC], findices[cellCIndices[d]-startC], cellFIndices[cmap[d]]); 1257 cindices[cellCIndices[d]-startC] = cellCIndices[d]; 1258 findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]]; 1259 } 1260 } 1261 ierr = PetscFree(cmap);CHKERRQ(ierr); 1262 ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr); 1263 1264 ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr); 1265 ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr); 1266 ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr); 1267 ierr = ISDestroy(&cis);CHKERRQ(ierr); 1268 ierr = ISDestroy(&fis);CHKERRQ(ierr); 1269 ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr); 1270 ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr); 1271 ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1272 PetscFunctionReturn(0); 1273 } 1274