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