1cb1e1211SMatthew G Knepley #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2cb1e1211SMatthew G Knepley 30f2d7e86SMatthew G. Knepley #include <petsc-private/petscfeimpl.h> 415496722SMatthew G. Knepley #include <petsc-private/petscfvimpl.h> 5a0845e3aSMatthew G. Knepley 6cb1e1211SMatthew G Knepley #undef __FUNCT__ 7cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexGetScale" 8cb1e1211SMatthew G Knepley PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale) 9cb1e1211SMatthew G Knepley { 10cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 11cb1e1211SMatthew G Knepley 12cb1e1211SMatthew G Knepley PetscFunctionBegin; 13cb1e1211SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 14cb1e1211SMatthew G Knepley PetscValidPointer(scale, 3); 15cb1e1211SMatthew G Knepley *scale = mesh->scale[unit]; 16cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 17cb1e1211SMatthew G Knepley } 18cb1e1211SMatthew G Knepley 19cb1e1211SMatthew G Knepley #undef __FUNCT__ 20cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexSetScale" 21cb1e1211SMatthew G Knepley PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale) 22cb1e1211SMatthew G Knepley { 23cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 24cb1e1211SMatthew G Knepley 25cb1e1211SMatthew G Knepley PetscFunctionBegin; 26cb1e1211SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 27cb1e1211SMatthew G Knepley mesh->scale[unit] = scale; 28cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 29cb1e1211SMatthew G Knepley } 30cb1e1211SMatthew G Knepley 31cb1e1211SMatthew G Knepley PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k) 32cb1e1211SMatthew G Knepley { 33cb1e1211SMatthew G Knepley switch (i) { 34cb1e1211SMatthew G Knepley case 0: 35cb1e1211SMatthew G Knepley switch (j) { 36cb1e1211SMatthew G Knepley case 0: return 0; 37cb1e1211SMatthew G Knepley case 1: 38cb1e1211SMatthew G Knepley switch (k) { 39cb1e1211SMatthew G Knepley case 0: return 0; 40cb1e1211SMatthew G Knepley case 1: return 0; 41cb1e1211SMatthew G Knepley case 2: return 1; 42cb1e1211SMatthew G Knepley } 43cb1e1211SMatthew G Knepley case 2: 44cb1e1211SMatthew G Knepley switch (k) { 45cb1e1211SMatthew G Knepley case 0: return 0; 46cb1e1211SMatthew G Knepley case 1: return -1; 47cb1e1211SMatthew G Knepley case 2: return 0; 48cb1e1211SMatthew G Knepley } 49cb1e1211SMatthew G Knepley } 50cb1e1211SMatthew G Knepley case 1: 51cb1e1211SMatthew G Knepley switch (j) { 52cb1e1211SMatthew G Knepley case 0: 53cb1e1211SMatthew G Knepley switch (k) { 54cb1e1211SMatthew G Knepley case 0: return 0; 55cb1e1211SMatthew G Knepley case 1: return 0; 56cb1e1211SMatthew G Knepley case 2: return -1; 57cb1e1211SMatthew G Knepley } 58cb1e1211SMatthew G Knepley case 1: return 0; 59cb1e1211SMatthew G Knepley case 2: 60cb1e1211SMatthew G Knepley switch (k) { 61cb1e1211SMatthew G Knepley case 0: return 1; 62cb1e1211SMatthew G Knepley case 1: return 0; 63cb1e1211SMatthew G Knepley case 2: return 0; 64cb1e1211SMatthew G Knepley } 65cb1e1211SMatthew G Knepley } 66cb1e1211SMatthew G Knepley case 2: 67cb1e1211SMatthew G Knepley switch (j) { 68cb1e1211SMatthew G Knepley case 0: 69cb1e1211SMatthew G Knepley switch (k) { 70cb1e1211SMatthew G Knepley case 0: return 0; 71cb1e1211SMatthew G Knepley case 1: return 1; 72cb1e1211SMatthew G Knepley case 2: return 0; 73cb1e1211SMatthew G Knepley } 74cb1e1211SMatthew G Knepley case 1: 75cb1e1211SMatthew G Knepley switch (k) { 76cb1e1211SMatthew G Knepley case 0: return -1; 77cb1e1211SMatthew G Knepley case 1: return 0; 78cb1e1211SMatthew G Knepley case 2: return 0; 79cb1e1211SMatthew G Knepley } 80cb1e1211SMatthew G Knepley case 2: return 0; 81cb1e1211SMatthew G Knepley } 82cb1e1211SMatthew G Knepley } 83cb1e1211SMatthew G Knepley return 0; 84cb1e1211SMatthew G Knepley } 85cb1e1211SMatthew G Knepley 86cb1e1211SMatthew G Knepley #undef __FUNCT__ 87cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexCreateRigidBody" 88cb1e1211SMatthew G Knepley /*@C 89cb1e1211SMatthew G Knepley DMPlexCreateRigidBody - create rigid body modes from coordinates 90cb1e1211SMatthew G Knepley 91cb1e1211SMatthew G Knepley Collective on DM 92cb1e1211SMatthew G Knepley 93cb1e1211SMatthew G Knepley Input Arguments: 94cb1e1211SMatthew G Knepley + dm - the DM 95cb1e1211SMatthew G Knepley . section - the local section associated with the rigid field, or NULL for the default section 96cb1e1211SMatthew G Knepley - globalSection - the global section associated with the rigid field, or NULL for the default section 97cb1e1211SMatthew G Knepley 98cb1e1211SMatthew G Knepley Output Argument: 99cb1e1211SMatthew G Knepley . sp - the null space 100cb1e1211SMatthew G Knepley 101cb1e1211SMatthew G Knepley Note: This is necessary to take account of Dirichlet conditions on the displacements 102cb1e1211SMatthew G Knepley 103cb1e1211SMatthew G Knepley Level: advanced 104cb1e1211SMatthew G Knepley 105cb1e1211SMatthew G Knepley .seealso: MatNullSpaceCreate() 106cb1e1211SMatthew G Knepley @*/ 107cb1e1211SMatthew G Knepley PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp) 108cb1e1211SMatthew G Knepley { 109cb1e1211SMatthew G Knepley MPI_Comm comm; 110cb1e1211SMatthew G Knepley Vec coordinates, localMode, mode[6]; 111cb1e1211SMatthew G Knepley PetscSection coordSection; 112cb1e1211SMatthew G Knepley PetscScalar *coords; 113cb1e1211SMatthew G Knepley PetscInt dim, vStart, vEnd, v, n, m, d, i, j; 114cb1e1211SMatthew G Knepley PetscErrorCode ierr; 115cb1e1211SMatthew G Knepley 116cb1e1211SMatthew G Knepley PetscFunctionBegin; 117cb1e1211SMatthew G Knepley ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 118c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 119cb1e1211SMatthew G Knepley if (dim == 1) { 120cb1e1211SMatthew G Knepley ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr); 121cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 122cb1e1211SMatthew G Knepley } 123cb1e1211SMatthew G Knepley if (!section) {ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr);} 124cb1e1211SMatthew G Knepley if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);} 125cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr); 126cb1e1211SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 12769d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 128cb1e1211SMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 129cb1e1211SMatthew G Knepley m = (dim*(dim+1))/2; 130cb1e1211SMatthew G Knepley ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr); 131cb1e1211SMatthew G Knepley ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr); 132cb1e1211SMatthew G Knepley ierr = VecSetUp(mode[0]);CHKERRQ(ierr); 133cb1e1211SMatthew G Knepley for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);} 134cb1e1211SMatthew G Knepley /* Assume P1 */ 135cb1e1211SMatthew G Knepley ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr); 136cb1e1211SMatthew G Knepley for (d = 0; d < dim; ++d) { 137cb1e1211SMatthew G Knepley PetscScalar values[3] = {0.0, 0.0, 0.0}; 138cb1e1211SMatthew G Knepley 139cb1e1211SMatthew G Knepley values[d] = 1.0; 140cb1e1211SMatthew G Knepley ierr = VecSet(localMode, 0.0);CHKERRQ(ierr); 141cb1e1211SMatthew G Knepley for (v = vStart; v < vEnd; ++v) { 142cb1e1211SMatthew G Knepley ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr); 143cb1e1211SMatthew G Knepley } 144cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 145cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 146cb1e1211SMatthew G Knepley } 147cb1e1211SMatthew G Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 148cb1e1211SMatthew G Knepley for (d = dim; d < dim*(dim+1)/2; ++d) { 149cb1e1211SMatthew G Knepley PetscInt i, j, k = dim > 2 ? d - dim : d; 150cb1e1211SMatthew G Knepley 151cb1e1211SMatthew G Knepley ierr = VecSet(localMode, 0.0);CHKERRQ(ierr); 152cb1e1211SMatthew G Knepley for (v = vStart; v < vEnd; ++v) { 153cb1e1211SMatthew G Knepley PetscScalar values[3] = {0.0, 0.0, 0.0}; 154cb1e1211SMatthew G Knepley PetscInt off; 155cb1e1211SMatthew G Knepley 156cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 157cb1e1211SMatthew G Knepley for (i = 0; i < dim; ++i) { 158cb1e1211SMatthew G Knepley for (j = 0; j < dim; ++j) { 159cb1e1211SMatthew G Knepley values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]); 160cb1e1211SMatthew G Knepley } 161cb1e1211SMatthew G Knepley } 162cb1e1211SMatthew G Knepley ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr); 163cb1e1211SMatthew G Knepley } 164cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 165cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 166cb1e1211SMatthew G Knepley } 167cb1e1211SMatthew G Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 168cb1e1211SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr); 169cb1e1211SMatthew G Knepley for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);} 170cb1e1211SMatthew G Knepley /* Orthonormalize system */ 171cb1e1211SMatthew G Knepley for (i = dim; i < m; ++i) { 172cb1e1211SMatthew G Knepley PetscScalar dots[6]; 173cb1e1211SMatthew G Knepley 174cb1e1211SMatthew G Knepley ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr); 175cb1e1211SMatthew G Knepley for (j = 0; j < i; ++j) dots[j] *= -1.0; 176cb1e1211SMatthew G Knepley ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr); 177cb1e1211SMatthew G Knepley ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr); 178cb1e1211SMatthew G Knepley } 179cb1e1211SMatthew G Knepley ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr); 180cb1e1211SMatthew G Knepley for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);} 181cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 182cb1e1211SMatthew G Knepley } 183cb1e1211SMatthew G Knepley 184cb1e1211SMatthew G Knepley #undef __FUNCT__ 185a18a7fb9SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFunctionLabelLocal" 186*bf3434eeSMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLabelLocal(DM dm, DMLabel label, PetscInt numIds, const PetscInt ids[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 187a18a7fb9SMatthew G. Knepley { 188*bf3434eeSMatthew G. Knepley PetscFE fe; 189a18a7fb9SMatthew G. Knepley PetscDualSpace *sp; 190a18a7fb9SMatthew G. Knepley PetscSection section; 191a18a7fb9SMatthew G. Knepley PetscScalar *values; 192ad96f515SMatthew G. Knepley PetscBool *fieldActive; 193a18a7fb9SMatthew G. Knepley PetscInt numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, f, d, v, i, comp; 194a18a7fb9SMatthew G. Knepley PetscErrorCode ierr; 195a18a7fb9SMatthew G. Knepley 196a18a7fb9SMatthew G. Knepley PetscFunctionBegin; 197c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 198a18a7fb9SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 199a18a7fb9SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 200e1d0b1adSMatthew G. Knepley ierr = PetscMalloc1(numFields,&sp);CHKERRQ(ierr); 201a18a7fb9SMatthew G. Knepley for (f = 0; f < numFields; ++f) { 202*bf3434eeSMatthew G. Knepley ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr); 203*bf3434eeSMatthew G. Knepley ierr = PetscFEGetDualSpace(fe, &sp[f]);CHKERRQ(ierr); 204*bf3434eeSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr); 205a18a7fb9SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 206a18a7fb9SMatthew G. Knepley totDim += spDim*numComp; 207a18a7fb9SMatthew G. Knepley } 208a18a7fb9SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 209a18a7fb9SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 210a18a7fb9SMatthew G. Knepley if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 211a18a7fb9SMatthew G. Knepley ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 212ad96f515SMatthew G. Knepley ierr = DMGetWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr); 213ad96f515SMatthew G. Knepley for (f = 0; f < numFields; ++f) fieldActive[f] = funcs[f] ? PETSC_TRUE : PETSC_FALSE; 214a18a7fb9SMatthew G. Knepley for (i = 0; i < numIds; ++i) { 215a18a7fb9SMatthew G. Knepley IS pointIS; 216a18a7fb9SMatthew G. Knepley const PetscInt *points; 217a18a7fb9SMatthew G. Knepley PetscInt n, p; 218a18a7fb9SMatthew G. Knepley 219a18a7fb9SMatthew G. Knepley ierr = DMLabelGetStratumIS(label, ids[i], &pointIS);CHKERRQ(ierr); 220a18a7fb9SMatthew G. Knepley ierr = ISGetLocalSize(pointIS, &n);CHKERRQ(ierr); 221a18a7fb9SMatthew G. Knepley ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 222a18a7fb9SMatthew G. Knepley for (p = 0; p < n; ++p) { 223a18a7fb9SMatthew G. Knepley const PetscInt point = points[p]; 224e1d0b1adSMatthew G. Knepley PetscFECellGeom geom; 225a18a7fb9SMatthew G. Knepley 226a18a7fb9SMatthew G. Knepley if ((point < cStart) || (point >= cEnd)) continue; 227e1d0b1adSMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, point, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr); 228a18a7fb9SMatthew G. Knepley for (f = 0, v = 0; f < numFields; ++f) { 229a18a7fb9SMatthew G. Knepley void * const ctx = ctxs ? ctxs[f] : NULL; 230*bf3434eeSMatthew G. Knepley 231*bf3434eeSMatthew G. Knepley ierr = DMGetField(dm, f, (PetscObject *) &fe);CHKERRQ(ierr); 232*bf3434eeSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &numComp);CHKERRQ(ierr); 233a18a7fb9SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 234a18a7fb9SMatthew G. Knepley for (d = 0; d < spDim; ++d) { 235a18a7fb9SMatthew G. Knepley if (funcs[f]) { 236e1d0b1adSMatthew G. Knepley ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr); 237a18a7fb9SMatthew G. Knepley } else { 238a18a7fb9SMatthew G. Knepley for (comp = 0; comp < numComp; ++comp) values[v+comp] = 0.0; 239a18a7fb9SMatthew G. Knepley } 240a18a7fb9SMatthew G. Knepley v += numComp; 241a18a7fb9SMatthew G. Knepley } 242a18a7fb9SMatthew G. Knepley } 243ad96f515SMatthew G. Knepley ierr = DMPlexVecSetFieldClosure_Internal(dm, section, localX, fieldActive, point, values, mode);CHKERRQ(ierr); 244a18a7fb9SMatthew G. Knepley } 245a18a7fb9SMatthew G. Knepley ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 246a18a7fb9SMatthew G. Knepley ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 247a18a7fb9SMatthew G. Knepley } 248a18a7fb9SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 249ad96f515SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, numFields, PETSC_BOOL, &fieldActive);CHKERRQ(ierr); 250e1d0b1adSMatthew G. Knepley ierr = PetscFree(sp);CHKERRQ(ierr); 251a18a7fb9SMatthew G. Knepley PetscFunctionReturn(0); 252a18a7fb9SMatthew G. Knepley } 253a18a7fb9SMatthew G. Knepley 254a18a7fb9SMatthew G. Knepley #undef __FUNCT__ 255cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunctionLocal" 2560f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLocal(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 257cb1e1211SMatthew G Knepley { 25872f94c41SMatthew G. Knepley PetscDualSpace *sp; 25915496722SMatthew G. Knepley PetscInt *numComp; 26072f94c41SMatthew G. Knepley PetscSection section; 26172f94c41SMatthew G. Knepley PetscScalar *values; 26215496722SMatthew G. Knepley PetscInt numFields, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v, comp; 263cb1e1211SMatthew G Knepley PetscErrorCode ierr; 264cb1e1211SMatthew G Knepley 265cb1e1211SMatthew G Knepley PetscFunctionBegin; 266cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 26772f94c41SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 26815496722SMatthew G. Knepley ierr = PetscMalloc2(numFields, &sp, numFields, &numComp);CHKERRQ(ierr); 26972f94c41SMatthew G. Knepley for (f = 0; f < numFields; ++f) { 27015496722SMatthew G. Knepley PetscObject obj; 27115496722SMatthew G. Knepley PetscClassId id; 2720f2d7e86SMatthew G. Knepley 27315496722SMatthew G. Knepley ierr = DMGetField(dm, f, &obj);CHKERRQ(ierr); 27415496722SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 27515496722SMatthew G. Knepley if (id == PETSCFE_CLASSID) { 27615496722SMatthew G. Knepley PetscFE fe = (PetscFE) obj; 27715496722SMatthew G. Knepley 27815496722SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &numComp[f]);CHKERRQ(ierr); 2790f2d7e86SMatthew G. Knepley ierr = PetscFEGetDualSpace(fe, &sp[f]);CHKERRQ(ierr); 28015496722SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) sp[f]);CHKERRQ(ierr); 28115496722SMatthew G. Knepley } else if (id == PETSCFV_CLASSID) { 28215496722SMatthew G. Knepley PetscFV fv = (PetscFV) obj; 28315496722SMatthew G. Knepley PetscQuadrature q; 28415496722SMatthew G. Knepley 28515496722SMatthew G. Knepley ierr = PetscFVGetNumComponents(fv, &numComp[f]);CHKERRQ(ierr); 28615496722SMatthew G. Knepley ierr = PetscFVGetQuadrature(fv, &q);CHKERRQ(ierr); 28715496722SMatthew G. Knepley ierr = PetscDualSpaceCreate(PetscObjectComm((PetscObject) fv), &sp[f]);CHKERRQ(ierr); 28815496722SMatthew G. Knepley ierr = PetscDualSpaceSetDM(sp[f], dm);CHKERRQ(ierr); 28915496722SMatthew G. Knepley ierr = PetscDualSpaceSetType(sp[f], PETSCDUALSPACESIMPLE);CHKERRQ(ierr); 29015496722SMatthew G. Knepley ierr = PetscDualSpaceSimpleSetDimension(sp[f], 1);CHKERRQ(ierr); 29115496722SMatthew G. Knepley ierr = PetscDualSpaceSimpleSetFunctional(sp[f], 0, q);CHKERRQ(ierr); 29215496722SMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 29372f94c41SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 29415496722SMatthew G. Knepley totDim += spDim*numComp[f]; 295cb1e1211SMatthew G Knepley } 296c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 29772f94c41SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 29872f94c41SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 29972f94c41SMatthew G. Knepley if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 30072f94c41SMatthew G. Knepley ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 30172f94c41SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 302e1d0b1adSMatthew G. Knepley PetscFECellGeom geom; 303cb1e1211SMatthew G Knepley 304e1d0b1adSMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, geom.v0, geom.J, NULL, &geom.detJ);CHKERRQ(ierr); 30572f94c41SMatthew G. Knepley for (f = 0, v = 0; f < numFields; ++f) { 306c110b1eeSGeoffrey Irving void *const ctx = ctxs ? ctxs[f] : NULL; 3070f2d7e86SMatthew G. Knepley 30872f94c41SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 30972f94c41SMatthew G. Knepley for (d = 0; d < spDim; ++d) { 310120386c5SMatthew G. Knepley if (funcs[f]) { 311e1d0b1adSMatthew G. Knepley ierr = PetscDualSpaceApply(sp[f], d, &geom, numComp[f], funcs[f], ctx, &values[v]);CHKERRQ(ierr); 312120386c5SMatthew G. Knepley } else { 31315496722SMatthew G. Knepley for (comp = 0; comp < numComp[f]; ++comp) values[v+comp] = 0.0; 314120386c5SMatthew G. Knepley } 31515496722SMatthew G. Knepley v += numComp[f]; 316cb1e1211SMatthew G Knepley } 317cb1e1211SMatthew G Knepley } 31872f94c41SMatthew G. Knepley ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr); 319cb1e1211SMatthew G Knepley } 32072f94c41SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 32115496722SMatthew G. Knepley for (f = 0; f < numFields; ++f) {ierr = PetscDualSpaceDestroy(&sp[f]);CHKERRQ(ierr);} 32215496722SMatthew G. Knepley ierr = PetscFree2(sp, numComp);CHKERRQ(ierr); 323cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 324cb1e1211SMatthew G Knepley } 325cb1e1211SMatthew G Knepley 326cb1e1211SMatthew G Knepley #undef __FUNCT__ 327cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunction" 328cb1e1211SMatthew G Knepley /*@C 329cb1e1211SMatthew G Knepley DMPlexProjectFunction - This projects the given function into the function space provided. 330cb1e1211SMatthew G Knepley 331cb1e1211SMatthew G Knepley Input Parameters: 332cb1e1211SMatthew G Knepley + dm - The DM 33372f94c41SMatthew G. Knepley . funcs - The coordinate functions to evaluate, one per field 334c110b1eeSGeoffrey Irving . ctxs - Optional array of contexts to pass to each coordinate function. ctxs itself may be null. 335cb1e1211SMatthew G Knepley - mode - The insertion mode for values 336cb1e1211SMatthew G Knepley 337cb1e1211SMatthew G Knepley Output Parameter: 338cb1e1211SMatthew G Knepley . X - vector 339cb1e1211SMatthew G Knepley 340cb1e1211SMatthew G Knepley Level: developer 341cb1e1211SMatthew G Knepley 342878cb397SSatish Balay .seealso: DMPlexComputeL2Diff() 343878cb397SSatish Balay @*/ 3440f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexProjectFunction(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X) 345cb1e1211SMatthew G Knepley { 346cb1e1211SMatthew G Knepley Vec localX; 347cb1e1211SMatthew G Knepley PetscErrorCode ierr; 348cb1e1211SMatthew G Knepley 349cb1e1211SMatthew G Knepley PetscFunctionBegin; 3509a800dd8SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 351cb1e1211SMatthew G Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 3520f2d7e86SMatthew G. Knepley ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, mode, localX);CHKERRQ(ierr); 353cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr); 354cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr); 355cb1e1211SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 356cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 357cb1e1211SMatthew G Knepley } 358cb1e1211SMatthew G Knepley 35955f2e967SMatthew G. Knepley #undef __FUNCT__ 3600f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectFieldLocal" 3613bc3b0a0SMatthew G. Knepley 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) 3620f2d7e86SMatthew G. Knepley { 3630f2d7e86SMatthew G. Knepley DM dmAux; 3642764a2aaSMatthew G. Knepley PetscDS prob, probAux; 3650f2d7e86SMatthew G. Knepley Vec A; 366326413afSMatthew G. Knepley PetscSection section, sectionAux; 367326413afSMatthew G. Knepley PetscScalar *values, *u, *u_x, *a, *a_x; 3680f2d7e86SMatthew G. Knepley PetscReal *x, *v0, *J, *invJ, detJ, **basisField, **basisFieldDer, **basisFieldAux, **basisFieldDerAux; 369326413afSMatthew G. Knepley PetscInt Nf, dim, spDim, totDim, numValues, cStart, cEnd, c, f, d, v, comp; 3700f2d7e86SMatthew G. Knepley PetscErrorCode ierr; 3710f2d7e86SMatthew G. Knepley 3720f2d7e86SMatthew G. Knepley PetscFunctionBegin; 3732764a2aaSMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 374c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 3750f2d7e86SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 3760f2d7e86SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 3770f2d7e86SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 3782764a2aaSMatthew G. Knepley ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 3792764a2aaSMatthew G. Knepley ierr = PetscDSGetTabulation(prob, &basisField, &basisFieldDer);CHKERRQ(ierr); 3802764a2aaSMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(prob, &u, NULL, &u_x);CHKERRQ(ierr); 3812764a2aaSMatthew G. Knepley ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr); 3820f2d7e86SMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 3830f2d7e86SMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 3840f2d7e86SMatthew G. Knepley if (dmAux) { 3852764a2aaSMatthew G. Knepley ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 386326413afSMatthew G. Knepley ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 3872764a2aaSMatthew G. Knepley ierr = PetscDSGetTabulation(prob, &basisFieldAux, &basisFieldDerAux);CHKERRQ(ierr); 3882764a2aaSMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, &a_x);CHKERRQ(ierr); 3890f2d7e86SMatthew G. Knepley } 3900f2d7e86SMatthew G. Knepley ierr = DMPlexInsertBoundaryValuesFEM(dm, localU);CHKERRQ(ierr); 3910f2d7e86SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 3920f2d7e86SMatthew G. Knepley if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 3930f2d7e86SMatthew G. Knepley ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 3940f2d7e86SMatthew G. Knepley ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 3950f2d7e86SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 396326413afSMatthew G. Knepley PetscScalar *coefficients = NULL, *coefficientsAux = NULL; 397326413afSMatthew G. Knepley 3988e0841e0SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 399326413afSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr); 400326413afSMatthew G. Knepley if (dmAux) {ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);} 4010f2d7e86SMatthew G. Knepley for (f = 0, v = 0; f < Nf; ++f) { 4023113607cSMatthew G. Knepley PetscFE fe; 4033113607cSMatthew G. Knepley PetscDualSpace sp; 4043113607cSMatthew G. Knepley PetscInt Ncf; 4053113607cSMatthew G. Knepley 4062764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 4073113607cSMatthew G. Knepley ierr = PetscFEGetDualSpace(fe, &sp);CHKERRQ(ierr); 4083113607cSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Ncf);CHKERRQ(ierr); 4093113607cSMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp, &spDim);CHKERRQ(ierr); 4100f2d7e86SMatthew G. Knepley for (d = 0; d < spDim; ++d) { 4110f2d7e86SMatthew G. Knepley PetscQuadrature quad; 4120f2d7e86SMatthew G. Knepley const PetscReal *points, *weights; 4130f2d7e86SMatthew G. Knepley PetscInt numPoints, q; 4140f2d7e86SMatthew G. Knepley 4150f2d7e86SMatthew G. Knepley if (funcs[f]) { 4163113607cSMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(sp, d, &quad);CHKERRQ(ierr); 4170f2d7e86SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &numPoints, &points, &weights);CHKERRQ(ierr); 4183113607cSMatthew G. Knepley ierr = PetscFEGetTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr); 4190f2d7e86SMatthew G. Knepley for (q = 0; q < numPoints; ++q) { 4200f2d7e86SMatthew G. Knepley CoordinatesRefToReal(dim, dim, v0, J, &points[q*dim], x); 421326413afSMatthew G. Knepley ierr = EvaluateFieldJets(prob, PETSC_FALSE, q, invJ, coefficients, NULL, u, u_x, NULL);CHKERRQ(ierr); 422326413afSMatthew G. Knepley ierr = EvaluateFieldJets(probAux, PETSC_FALSE, q, invJ, coefficientsAux, NULL, a, a_x, NULL);CHKERRQ(ierr); 4233bc3b0a0SMatthew G. Knepley (*funcs[f])(u, NULL, u_x, a, NULL, a_x, x, &values[v]); 4240f2d7e86SMatthew G. Knepley } 4253113607cSMatthew G. Knepley ierr = PetscFERestoreTabulation(fe, numPoints, points, &basisField[f], &basisFieldDer[f], NULL);CHKERRQ(ierr); 4260f2d7e86SMatthew G. Knepley } else { 4270f2d7e86SMatthew G. Knepley for (comp = 0; comp < Ncf; ++comp) values[v+comp] = 0.0; 4280f2d7e86SMatthew G. Knepley } 4290f2d7e86SMatthew G. Knepley v += Ncf; 4300f2d7e86SMatthew G. Knepley } 4310f2d7e86SMatthew G. Knepley } 432326413afSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, section, localU, c, NULL, &coefficients);CHKERRQ(ierr); 433326413afSMatthew G. Knepley if (dmAux) {ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &coefficientsAux);CHKERRQ(ierr);} 4340f2d7e86SMatthew G. Knepley ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr); 4350f2d7e86SMatthew G. Knepley } 4360f2d7e86SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 4370f2d7e86SMatthew G. Knepley ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr); 4380f2d7e86SMatthew G. Knepley PetscFunctionReturn(0); 4390f2d7e86SMatthew G. Knepley } 4400f2d7e86SMatthew G. Knepley 4410f2d7e86SMatthew G. Knepley #undef __FUNCT__ 4420f2d7e86SMatthew G. Knepley #define __FUNCT__ "DMPlexProjectField" 4430f2d7e86SMatthew G. Knepley /*@C 4440f2d7e86SMatthew G. Knepley DMPlexProjectField - This projects the given function of the fields into the function space provided. 4450f2d7e86SMatthew G. Knepley 4460f2d7e86SMatthew G. Knepley Input Parameters: 4470f2d7e86SMatthew G. Knepley + dm - The DM 4480f2d7e86SMatthew G. Knepley . U - The input field vector 4490f2d7e86SMatthew G. Knepley . funcs - The functions to evaluate, one per field 4500f2d7e86SMatthew G. Knepley - mode - The insertion mode for values 4510f2d7e86SMatthew G. Knepley 4520f2d7e86SMatthew G. Knepley Output Parameter: 4530f2d7e86SMatthew G. Knepley . X - The output vector 4540f2d7e86SMatthew G. Knepley 4550f2d7e86SMatthew G. Knepley Level: developer 4560f2d7e86SMatthew G. Knepley 4570f2d7e86SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff() 4580f2d7e86SMatthew G. Knepley @*/ 4593bc3b0a0SMatthew G. Knepley 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) 4600f2d7e86SMatthew G. Knepley { 4610f2d7e86SMatthew G. Knepley Vec localX, localU; 4620f2d7e86SMatthew G. Knepley PetscErrorCode ierr; 4630f2d7e86SMatthew G. Knepley 4640f2d7e86SMatthew G. Knepley PetscFunctionBegin; 4650f2d7e86SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4660f2d7e86SMatthew G. Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 4670f2d7e86SMatthew G. Knepley ierr = DMGetLocalVector(dm, &localU);CHKERRQ(ierr); 4680f2d7e86SMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, U, INSERT_VALUES, localU);CHKERRQ(ierr); 4690f2d7e86SMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, U, INSERT_VALUES, localU);CHKERRQ(ierr); 4703113607cSMatthew G. Knepley ierr = DMPlexProjectFieldLocal(dm, localU, funcs, mode, localX);CHKERRQ(ierr); 4710f2d7e86SMatthew G. Knepley ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr); 4720f2d7e86SMatthew G. Knepley ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr); 4730f2d7e86SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 4740f2d7e86SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localU);CHKERRQ(ierr); 4750f2d7e86SMatthew G. Knepley PetscFunctionReturn(0); 4760f2d7e86SMatthew G. Knepley } 4770f2d7e86SMatthew G. Knepley 4780f2d7e86SMatthew G. Knepley #undef __FUNCT__ 4793351dd3dSMatthew G. Knepley #define __FUNCT__ "DMPlexInsertBoundaryValuesFEM" 4803351dd3dSMatthew G. Knepley PetscErrorCode DMPlexInsertBoundaryValuesFEM(DM dm, Vec localX) 48155f2e967SMatthew G. Knepley { 48255f2e967SMatthew G. Knepley void (**funcs)(const PetscReal x[], PetscScalar *u, void *ctx); 48355f2e967SMatthew G. Knepley void **ctxs; 48455f2e967SMatthew G. Knepley PetscFE *fe; 48555f2e967SMatthew G. Knepley PetscInt numFields, f, numBd, b; 48655f2e967SMatthew G. Knepley PetscErrorCode ierr; 48755f2e967SMatthew G. Knepley 48855f2e967SMatthew G. Knepley PetscFunctionBegin; 48955f2e967SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 49055f2e967SMatthew G. Knepley PetscValidHeaderSpecific(localX, VEC_CLASSID, 2); 49155f2e967SMatthew G. Knepley ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr); 49255f2e967SMatthew G. Knepley ierr = PetscMalloc3(numFields,&fe,numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr); 49355f2e967SMatthew G. Knepley for (f = 0; f < numFields; ++f) {ierr = DMGetField(dm, f, (PetscObject *) &fe[f]);CHKERRQ(ierr);} 49455f2e967SMatthew G. Knepley /* OPT: Could attempt to do multiple BCs at once */ 49555f2e967SMatthew G. Knepley ierr = DMPlexGetNumBoundary(dm, &numBd);CHKERRQ(ierr); 49655f2e967SMatthew G. Knepley for (b = 0; b < numBd; ++b) { 497a18a7fb9SMatthew G. Knepley DMLabel label; 49855f2e967SMatthew G. Knepley const PetscInt *ids; 49963d5297fSMatthew G. Knepley const char *labelname; 50055f2e967SMatthew G. Knepley PetscInt numids, field; 50155f2e967SMatthew G. Knepley PetscBool isEssential; 50255f2e967SMatthew G. Knepley void (*func)(); 50355f2e967SMatthew G. Knepley void *ctx; 50455f2e967SMatthew G. Knepley 50563d5297fSMatthew G. Knepley ierr = DMPlexGetBoundary(dm, b, &isEssential, NULL, &labelname, &field, &func, &numids, &ids, &ctx);CHKERRQ(ierr); 50663d5297fSMatthew G. Knepley ierr = DMPlexGetLabel(dm, labelname, &label);CHKERRQ(ierr); 50755f2e967SMatthew G. Knepley for (f = 0; f < numFields; ++f) { 50855f2e967SMatthew G. Knepley funcs[f] = field == f ? (void (*)(const PetscReal[], PetscScalar *, void *)) func : NULL; 50955f2e967SMatthew G. Knepley ctxs[f] = field == f ? ctx : NULL; 51055f2e967SMatthew G. Knepley } 511a18a7fb9SMatthew G. Knepley ierr = DMPlexProjectFunctionLabelLocal(dm, label, numids, ids, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 51255f2e967SMatthew G. Knepley } 51355f2e967SMatthew G. Knepley ierr = PetscFree3(fe,funcs,ctxs);CHKERRQ(ierr); 51455f2e967SMatthew G. Knepley PetscFunctionReturn(0); 51555f2e967SMatthew G. Knepley } 51655f2e967SMatthew G. Knepley 517cb1e1211SMatthew G Knepley #undef __FUNCT__ 518cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff" 519cb1e1211SMatthew G Knepley /*@C 520cb1e1211SMatthew G Knepley DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 521cb1e1211SMatthew G Knepley 522cb1e1211SMatthew G Knepley Input Parameters: 523cb1e1211SMatthew G Knepley + dm - The DM 524cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component 52551259fa3SMatthew G. Knepley . ctxs - Optional array of contexts to pass to each function, or NULL. 526cb1e1211SMatthew G Knepley - X - The coefficient vector u_h 527cb1e1211SMatthew G Knepley 528cb1e1211SMatthew G Knepley Output Parameter: 529cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2 530cb1e1211SMatthew G Knepley 531cb1e1211SMatthew G Knepley Level: developer 532cb1e1211SMatthew G Knepley 53315496722SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2FieldDiff(), DMPlexComputeL2GradientDiff() 534878cb397SSatish Balay @*/ 5350f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2Diff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff) 536cb1e1211SMatthew G Knepley { 537cb1e1211SMatthew G Knepley const PetscInt debug = 0; 538cb1e1211SMatthew G Knepley PetscSection section; 539c5bbbd5bSMatthew G. Knepley PetscQuadrature quad; 540cb1e1211SMatthew G Knepley Vec localX; 54115496722SMatthew G. Knepley PetscScalar *funcVal, *interpolant; 542cb1e1211SMatthew G Knepley PetscReal *coords, *v0, *J, *invJ, detJ; 543cb1e1211SMatthew G Knepley PetscReal localDiff = 0.0; 54415496722SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 54515496722SMatthew G. Knepley PetscInt dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, c, field, fieldOffset; 546cb1e1211SMatthew G Knepley PetscErrorCode ierr; 547cb1e1211SMatthew G Knepley 548cb1e1211SMatthew G Knepley PetscFunctionBegin; 549c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 550cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 551cb1e1211SMatthew G Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 552cb1e1211SMatthew G Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 553cb1e1211SMatthew G Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 554cb1e1211SMatthew G Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 555cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 55615496722SMatthew G. Knepley PetscObject obj; 55715496722SMatthew G. Knepley PetscClassId id; 558c5bbbd5bSMatthew G. Knepley PetscInt Nc; 559c5bbbd5bSMatthew G. Knepley 56015496722SMatthew G. Knepley ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 56115496722SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 56215496722SMatthew G. Knepley if (id == PETSCFE_CLASSID) { 56315496722SMatthew G. Knepley PetscFE fe = (PetscFE) obj; 56415496722SMatthew G. Knepley 5650f2d7e86SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 5660f2d7e86SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 56715496722SMatthew G. Knepley } else if (id == PETSCFV_CLASSID) { 56815496722SMatthew G. Knepley PetscFV fv = (PetscFV) obj; 56915496722SMatthew G. Knepley 57015496722SMatthew G. Knepley ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 57115496722SMatthew G. Knepley ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 57215496722SMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 573c5bbbd5bSMatthew G. Knepley numComponents += Nc; 574cb1e1211SMatthew G Knepley } 57515496722SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 5760f2d7e86SMatthew G. Knepley ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 57715496722SMatthew G. Knepley ierr = PetscMalloc6(numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 578cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 579cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 580a1e44745SMatthew G. Knepley PetscScalar *x = NULL; 581cb1e1211SMatthew G Knepley PetscReal elemDiff = 0.0; 582cb1e1211SMatthew G Knepley 5838e0841e0SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 584cb1e1211SMatthew G Knepley if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 585cb1e1211SMatthew G Knepley ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 586cb1e1211SMatthew G Knepley 58715496722SMatthew G. Knepley for (field = 0, fieldOffset = 0; field < numFields; ++field) { 58815496722SMatthew G. Knepley PetscObject obj; 58915496722SMatthew G. Knepley PetscClassId id; 590c110b1eeSGeoffrey Irving void * const ctx = ctxs ? ctxs[field] : NULL; 59115496722SMatthew G. Knepley PetscInt Nb, Nc, q, fc; 592cb1e1211SMatthew G Knepley 59315496722SMatthew G. Knepley ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 59415496722SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 59515496722SMatthew G. Knepley if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 59615496722SMatthew G. Knepley else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 59715496722SMatthew G. Knepley else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 598cb1e1211SMatthew G Knepley if (debug) { 599cb1e1211SMatthew G Knepley char title[1024]; 600cb1e1211SMatthew G Knepley ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 60115496722SMatthew G. Knepley ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr); 602cb1e1211SMatthew G Knepley } 603cb1e1211SMatthew G Knepley for (q = 0; q < numQuadPoints; ++q) { 60415496722SMatthew G. Knepley CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords); 605c110b1eeSGeoffrey Irving (*funcs[field])(coords, funcVal, ctx); 60615496722SMatthew G. Knepley if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 60715496722SMatthew G. Knepley else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 60815496722SMatthew G. Knepley else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 60915496722SMatthew G. Knepley for (fc = 0; fc < Nc; ++fc) { 61015496722SMatthew G. Knepley 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);} 61115496722SMatthew G. Knepley elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ; 612cb1e1211SMatthew G Knepley } 613cb1e1211SMatthew G Knepley } 61415496722SMatthew G. Knepley fieldOffset += Nb*Nc; 615cb1e1211SMatthew G Knepley } 616cb1e1211SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 617cb1e1211SMatthew G Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 618cb1e1211SMatthew G Knepley localDiff += elemDiff; 619cb1e1211SMatthew G Knepley } 62015496722SMatthew G. Knepley ierr = PetscFree6(funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr); 621cb1e1211SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 62286a74ee0SMatthew G. Knepley ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 623cb1e1211SMatthew G Knepley *diff = PetscSqrtReal(*diff); 624cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 625cb1e1211SMatthew G Knepley } 626cb1e1211SMatthew G Knepley 627cb1e1211SMatthew G Knepley #undef __FUNCT__ 62840e14135SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2GradientDiff" 62940e14135SMatthew G. Knepley /*@C 63040e14135SMatthew G. Knepley DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h. 63140e14135SMatthew G. Knepley 63240e14135SMatthew G. Knepley Input Parameters: 63340e14135SMatthew G. Knepley + dm - The DM 63440e14135SMatthew G. Knepley . funcs - The gradient functions to evaluate for each field component 63551259fa3SMatthew G. Knepley . ctxs - Optional array of contexts to pass to each function, or NULL. 63640e14135SMatthew G. Knepley . X - The coefficient vector u_h 63740e14135SMatthew G. Knepley - n - The vector to project along 63840e14135SMatthew G. Knepley 63940e14135SMatthew G. Knepley Output Parameter: 64040e14135SMatthew G. Knepley . diff - The diff ||(grad u - grad u_h) . n||_2 64140e14135SMatthew G. Knepley 64240e14135SMatthew G. Knepley Level: developer 64340e14135SMatthew G. Knepley 64440e14135SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff() 64540e14135SMatthew G. Knepley @*/ 6460f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff) 647cb1e1211SMatthew G Knepley { 64840e14135SMatthew G. Knepley const PetscInt debug = 0; 649cb1e1211SMatthew G Knepley PetscSection section; 65040e14135SMatthew G. Knepley PetscQuadrature quad; 65140e14135SMatthew G. Knepley Vec localX; 65240e14135SMatthew G. Knepley PetscScalar *funcVal, *interpolantVec; 65340e14135SMatthew G. Knepley PetscReal *coords, *realSpaceDer, *v0, *J, *invJ, detJ; 65440e14135SMatthew G. Knepley PetscReal localDiff = 0.0; 65540e14135SMatthew G. Knepley PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 656cb1e1211SMatthew G Knepley PetscErrorCode ierr; 657cb1e1211SMatthew G Knepley 658cb1e1211SMatthew G Knepley PetscFunctionBegin; 659c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 66040e14135SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 66140e14135SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 66240e14135SMatthew G. Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 66340e14135SMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 66440e14135SMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 665652b88e8SMatthew G. Knepley for (field = 0; field < numFields; ++field) { 6660f2d7e86SMatthew G. Knepley PetscFE fe; 66740e14135SMatthew G. Knepley PetscInt Nc; 668652b88e8SMatthew G. Knepley 6690f2d7e86SMatthew G. Knepley ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 6700f2d7e86SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 6710f2d7e86SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 67240e14135SMatthew G. Knepley numComponents += Nc; 673652b88e8SMatthew G. Knepley } 67440e14135SMatthew G. Knepley /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */ 67540e14135SMatthew G. Knepley ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr); 67640e14135SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 67740e14135SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 67840e14135SMatthew G. Knepley PetscScalar *x = NULL; 67940e14135SMatthew G. Knepley PetscReal elemDiff = 0.0; 680652b88e8SMatthew G. Knepley 6818e0841e0SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 68240e14135SMatthew G. Knepley if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 68340e14135SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 68440e14135SMatthew G. Knepley 68540e14135SMatthew G. Knepley for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 6860f2d7e86SMatthew G. Knepley PetscFE fe; 68751259fa3SMatthew G. Knepley void * const ctx = ctxs ? ctxs[field] : NULL; 68821454ff5SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 68940e14135SMatthew G. Knepley PetscReal *basisDer; 69021454ff5SMatthew G. Knepley PetscInt numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g; 69140e14135SMatthew G. Knepley 6920f2d7e86SMatthew G. Knepley ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 69321454ff5SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 6940f2d7e86SMatthew G. Knepley ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 6950f2d7e86SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Ncomp);CHKERRQ(ierr); 6960f2d7e86SMatthew G. Knepley ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr); 69740e14135SMatthew G. Knepley if (debug) { 69840e14135SMatthew G. Knepley char title[1024]; 69940e14135SMatthew G. Knepley ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 70040e14135SMatthew G. Knepley ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr); 701652b88e8SMatthew G. Knepley } 70240e14135SMatthew G. Knepley for (q = 0; q < numQuadPoints; ++q) { 70340e14135SMatthew G. Knepley for (d = 0; d < dim; d++) { 70440e14135SMatthew G. Knepley coords[d] = v0[d]; 70540e14135SMatthew G. Knepley for (e = 0; e < dim; e++) { 70640e14135SMatthew G. Knepley coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 707652b88e8SMatthew G. Knepley } 70840e14135SMatthew G. Knepley } 70951259fa3SMatthew G. Knepley (*funcs[field])(coords, n, funcVal, ctx); 71040e14135SMatthew G. Knepley for (fc = 0; fc < Ncomp; ++fc) { 71140e14135SMatthew G. Knepley PetscScalar interpolant = 0.0; 71240e14135SMatthew G. Knepley 71340e14135SMatthew G. Knepley for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0; 71440e14135SMatthew G. Knepley for (f = 0; f < Nb; ++f) { 71540e14135SMatthew G. Knepley const PetscInt fidx = f*Ncomp+fc; 71640e14135SMatthew G. Knepley 71740e14135SMatthew G. Knepley for (d = 0; d < dim; ++d) { 71840e14135SMatthew G. Knepley realSpaceDer[d] = 0.0; 71940e14135SMatthew G. Knepley for (g = 0; g < dim; ++g) { 72040e14135SMatthew G. Knepley realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g]; 72140e14135SMatthew G. Knepley } 72240e14135SMatthew G. Knepley interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d]; 72340e14135SMatthew G. Knepley } 72440e14135SMatthew G. Knepley } 72540e14135SMatthew G. Knepley for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d]; 72640e14135SMatthew G. Knepley 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);} 72740e14135SMatthew G. Knepley elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 72840e14135SMatthew G. Knepley } 72940e14135SMatthew G. Knepley } 73040e14135SMatthew G. Knepley comp += Ncomp; 73140e14135SMatthew G. Knepley fieldOffset += Nb*Ncomp; 73240e14135SMatthew G. Knepley } 73340e14135SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 73440e14135SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 73540e14135SMatthew G. Knepley localDiff += elemDiff; 73640e14135SMatthew G. Knepley } 73740e14135SMatthew G. Knepley ierr = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr); 73840e14135SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 73940e14135SMatthew G. Knepley ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 74040e14135SMatthew G. Knepley *diff = PetscSqrtReal(*diff); 741cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 742cb1e1211SMatthew G Knepley } 743cb1e1211SMatthew G Knepley 744a0845e3aSMatthew G. Knepley #undef __FUNCT__ 74573d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2FieldDiff" 74615496722SMatthew G. Knepley /*@C 74715496722SMatthew G. Knepley DMPlexComputeL2FieldDiff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h, separated into field components. 74815496722SMatthew G. Knepley 74915496722SMatthew G. Knepley Input Parameters: 75015496722SMatthew G. Knepley + dm - The DM 75115496722SMatthew G. Knepley . funcs - The functions to evaluate for each field component 75215496722SMatthew G. Knepley . ctxs - Optional array of contexts to pass to each function, or NULL. 75315496722SMatthew G. Knepley - X - The coefficient vector u_h 75415496722SMatthew G. Knepley 75515496722SMatthew G. Knepley Output Parameter: 75615496722SMatthew G. Knepley . diff - The array of differences, ||u^f - u^f_h||_2 75715496722SMatthew G. Knepley 75815496722SMatthew G. Knepley Level: developer 75915496722SMatthew G. Knepley 76015496722SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff(), DMPlexComputeL2GradientDiff() 76115496722SMatthew G. Knepley @*/ 7620f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeL2FieldDiff(DM dm, void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal diff[]) 76373d901b8SMatthew G. Knepley { 76473d901b8SMatthew G. Knepley const PetscInt debug = 0; 76573d901b8SMatthew G. Knepley PetscSection section; 76673d901b8SMatthew G. Knepley PetscQuadrature quad; 76773d901b8SMatthew G. Knepley Vec localX; 76815496722SMatthew G. Knepley PetscScalar *funcVal, *interpolant; 76973d901b8SMatthew G. Knepley PetscReal *coords, *v0, *J, *invJ, detJ; 77073d901b8SMatthew G. Knepley PetscReal *localDiff; 77115496722SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 77215496722SMatthew G. Knepley PetscInt dim, numFields, numComponents = 0, numQuadPoints, cStart, cEnd, c, field, fieldOffset; 77373d901b8SMatthew G. Knepley PetscErrorCode ierr; 77473d901b8SMatthew G. Knepley 77573d901b8SMatthew G. Knepley PetscFunctionBegin; 776c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 77773d901b8SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 77873d901b8SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 77973d901b8SMatthew G. Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 78073d901b8SMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 78173d901b8SMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 78273d901b8SMatthew G. Knepley for (field = 0; field < numFields; ++field) { 78315496722SMatthew G. Knepley PetscObject obj; 78415496722SMatthew G. Knepley PetscClassId id; 78573d901b8SMatthew G. Knepley PetscInt Nc; 78673d901b8SMatthew G. Knepley 78715496722SMatthew G. Knepley ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 78815496722SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 78915496722SMatthew G. Knepley if (id == PETSCFE_CLASSID) { 79015496722SMatthew G. Knepley PetscFE fe = (PetscFE) obj; 79115496722SMatthew G. Knepley 7920f2d7e86SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 7930f2d7e86SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 79415496722SMatthew G. Knepley } else if (id == PETSCFV_CLASSID) { 79515496722SMatthew G. Knepley PetscFV fv = (PetscFV) obj; 79615496722SMatthew G. Knepley 79715496722SMatthew G. Knepley ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 79815496722SMatthew G. Knepley ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 79915496722SMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 80073d901b8SMatthew G. Knepley numComponents += Nc; 80173d901b8SMatthew G. Knepley } 80215496722SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 8030f2d7e86SMatthew G. Knepley ierr = DMPlexProjectFunctionLocal(dm, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 80415496722SMatthew G. Knepley ierr = PetscCalloc7(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 80573d901b8SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 80673d901b8SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 80773d901b8SMatthew G. Knepley PetscScalar *x = NULL; 80873d901b8SMatthew G. Knepley 8098e0841e0SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 81073d901b8SMatthew G. Knepley if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 81173d901b8SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 81273d901b8SMatthew G. Knepley 81315496722SMatthew G. Knepley for (field = 0, fieldOffset = 0; field < numFields; ++field) { 81415496722SMatthew G. Knepley PetscObject obj; 81515496722SMatthew G. Knepley PetscClassId id; 81673d901b8SMatthew G. Knepley void * const ctx = ctxs ? ctxs[field] : NULL; 81715496722SMatthew G. Knepley PetscInt Nb, Nc, q, fc; 81873d901b8SMatthew G. Knepley 81915496722SMatthew G. Knepley PetscReal elemDiff = 0.0; 82015496722SMatthew G. Knepley 82115496722SMatthew G. Knepley ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 82215496722SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 82315496722SMatthew G. Knepley if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 82415496722SMatthew G. Knepley else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 82515496722SMatthew G. Knepley else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 82673d901b8SMatthew G. Knepley if (debug) { 82773d901b8SMatthew G. Knepley char title[1024]; 82873d901b8SMatthew G. Knepley ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 82915496722SMatthew G. Knepley ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr); 83073d901b8SMatthew G. Knepley } 83173d901b8SMatthew G. Knepley for (q = 0; q < numQuadPoints; ++q) { 83215496722SMatthew G. Knepley CoordinatesRefToReal(dim, dim, v0, J, &quadPoints[q*dim], coords); 83373d901b8SMatthew G. Knepley (*funcs[field])(coords, funcVal, ctx); 83415496722SMatthew G. Knepley if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 83515496722SMatthew G. Knepley else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 83615496722SMatthew G. Knepley else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 83715496722SMatthew G. Knepley for (fc = 0; fc < Nc; ++fc) { 83815496722SMatthew G. Knepley 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);} 83915496722SMatthew G. Knepley elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*quadWeights[q]*detJ; 84073d901b8SMatthew G. Knepley } 84173d901b8SMatthew G. Knepley } 84215496722SMatthew G. Knepley fieldOffset += Nb*Nc; 84373d901b8SMatthew G. Knepley localDiff[field] += elemDiff; 84473d901b8SMatthew G. Knepley } 84573d901b8SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 84673d901b8SMatthew G. Knepley } 84773d901b8SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 84873d901b8SMatthew G. Knepley ierr = MPI_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 84973d901b8SMatthew G. Knepley for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]); 85015496722SMatthew G. Knepley ierr = PetscFree7(localDiff,funcVal,interpolant,coords,v0,J,invJ);CHKERRQ(ierr); 85173d901b8SMatthew G. Knepley PetscFunctionReturn(0); 85273d901b8SMatthew G. Knepley } 85373d901b8SMatthew G. Knepley 85473d901b8SMatthew G. Knepley #undef __FUNCT__ 85573d901b8SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeIntegralFEM" 85673d901b8SMatthew G. Knepley /*@ 85773d901b8SMatthew G. Knepley DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user 85873d901b8SMatthew G. Knepley 85973d901b8SMatthew G. Knepley Input Parameters: 86073d901b8SMatthew G. Knepley + dm - The mesh 86173d901b8SMatthew G. Knepley . X - Local input vector 86273d901b8SMatthew G. Knepley - user - The user context 86373d901b8SMatthew G. Knepley 86473d901b8SMatthew G. Knepley Output Parameter: 86573d901b8SMatthew G. Knepley . integral - Local integral for each field 86673d901b8SMatthew G. Knepley 86773d901b8SMatthew G. Knepley Level: developer 86873d901b8SMatthew G. Knepley 86973d901b8SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM() 87073d901b8SMatthew G. Knepley @*/ 8710f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user) 87273d901b8SMatthew G. Knepley { 87373d901b8SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 87473d901b8SMatthew G. Knepley DM dmAux; 87573d901b8SMatthew G. Knepley Vec localX, A; 8762764a2aaSMatthew G. Knepley PetscDS prob, probAux = NULL; 87773d901b8SMatthew G. Knepley PetscQuadrature q; 87873d901b8SMatthew G. Knepley PetscSection section, sectionAux; 879bbce034cSMatthew G. Knepley PetscFECellGeom *cgeom; 88073d901b8SMatthew G. Knepley PetscScalar *u, *a = NULL; 8810f2d7e86SMatthew G. Knepley PetscInt dim, Nf, f, numCells, cStart, cEnd, c; 8820f2d7e86SMatthew G. Knepley PetscInt totDim, totDimAux; 88373d901b8SMatthew G. Knepley PetscErrorCode ierr; 88473d901b8SMatthew G. Knepley 88573d901b8SMatthew G. Knepley PetscFunctionBegin; 88673d901b8SMatthew G. Knepley /*ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/ 88773d901b8SMatthew G. Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 88873d901b8SMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 88973d901b8SMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 890c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 89173d901b8SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 8922764a2aaSMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 8932764a2aaSMatthew G. Knepley ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 89473d901b8SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 89573d901b8SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 89673d901b8SMatthew G. Knepley numCells = cEnd - cStart; 8970f2d7e86SMatthew G. Knepley for (f = 0; f < Nf; ++f) {integral[f] = 0.0;} 89873d901b8SMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 89973d901b8SMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 90073d901b8SMatthew G. Knepley if (dmAux) { 90173d901b8SMatthew G. Knepley ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 9022764a2aaSMatthew G. Knepley ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 9032764a2aaSMatthew G. Knepley ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 90473d901b8SMatthew G. Knepley } 90573d901b8SMatthew G. Knepley ierr = DMPlexInsertBoundaryValuesFEM(dm, localX);CHKERRQ(ierr); 906bbce034cSMatthew G. Knepley ierr = PetscMalloc2(numCells*totDim,&u,numCells,&cgeom);CHKERRQ(ierr); 9070f2d7e86SMatthew G. Knepley if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);} 90873d901b8SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 90973d901b8SMatthew G. Knepley PetscScalar *x = NULL; 91073d901b8SMatthew G. Knepley PetscInt i; 91173d901b8SMatthew G. Knepley 912bbce034cSMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeom[c].v0, cgeom[c].J, cgeom[c].invJ, &cgeom[c].detJ);CHKERRQ(ierr); 913bbce034cSMatthew G. Knepley if (cgeom[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeom[c].detJ, c); 91473d901b8SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr); 9150f2d7e86SMatthew G. Knepley for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i]; 91673d901b8SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr); 91773d901b8SMatthew G. Knepley if (dmAux) { 91873d901b8SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 9190f2d7e86SMatthew G. Knepley for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i]; 92073d901b8SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 92173d901b8SMatthew G. Knepley } 92273d901b8SMatthew G. Knepley } 92373d901b8SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 9240f2d7e86SMatthew G. Knepley PetscFE fe; 92573d901b8SMatthew G. Knepley PetscInt numQuadPoints, Nb; 92673d901b8SMatthew G. Knepley /* Conforming batches */ 92773d901b8SMatthew G. Knepley PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 92873d901b8SMatthew G. Knepley /* Remainder */ 92973d901b8SMatthew G. Knepley PetscInt Nr, offset; 93073d901b8SMatthew G. Knepley 9312764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 9320f2d7e86SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 9330f2d7e86SMatthew G. Knepley ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 9340f2d7e86SMatthew G. Knepley ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 93573d901b8SMatthew G. Knepley ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 93673d901b8SMatthew G. Knepley blockSize = Nb*numQuadPoints; 93773d901b8SMatthew G. Knepley batchSize = numBlocks * blockSize; 9380f2d7e86SMatthew G. Knepley ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 93973d901b8SMatthew G. Knepley numChunks = numCells / (numBatches*batchSize); 94073d901b8SMatthew G. Knepley Ne = numChunks*numBatches*batchSize; 94173d901b8SMatthew G. Knepley Nr = numCells % (numBatches*batchSize); 94273d901b8SMatthew G. Knepley offset = numCells - Nr; 943bbce034cSMatthew G. Knepley ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeom, u, probAux, a, integral);CHKERRQ(ierr); 944bbce034cSMatthew G. Knepley ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeom[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], integral);CHKERRQ(ierr); 94573d901b8SMatthew G. Knepley } 946bbce034cSMatthew G. Knepley ierr = PetscFree2(u,cgeom);CHKERRQ(ierr); 94773d901b8SMatthew G. Knepley if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 94873d901b8SMatthew G. Knepley if (mesh->printFEM) { 94973d901b8SMatthew G. Knepley ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr); 95073d901b8SMatthew G. Knepley for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", integral[f]);CHKERRQ(ierr);} 95173d901b8SMatthew G. Knepley ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr); 95273d901b8SMatthew G. Knepley } 95373d901b8SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 95473d901b8SMatthew G. Knepley /* TODO: Allreduce for integral */ 95573d901b8SMatthew G. Knepley /*ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr);*/ 95673d901b8SMatthew G. Knepley PetscFunctionReturn(0); 95773d901b8SMatthew G. Knepley } 95873d901b8SMatthew G. Knepley 95973d901b8SMatthew G. Knepley #undef __FUNCT__ 960d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEM" 961d69c5d34SMatthew G. Knepley /*@ 962d69c5d34SMatthew G. Knepley DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM. 963d69c5d34SMatthew G. Knepley 964d69c5d34SMatthew G. Knepley Input Parameters: 965d69c5d34SMatthew G. Knepley + dmf - The fine mesh 966d69c5d34SMatthew G. Knepley . dmc - The coarse mesh 967d69c5d34SMatthew G. Knepley - user - The user context 968d69c5d34SMatthew G. Knepley 969d69c5d34SMatthew G. Knepley Output Parameter: 970934789fcSMatthew G. Knepley . In - The interpolation matrix 971d69c5d34SMatthew G. Knepley 972d69c5d34SMatthew G. Knepley Note: 973d69c5d34SMatthew G. Knepley The first member of the user context must be an FEMContext. 974d69c5d34SMatthew G. Knepley 975d69c5d34SMatthew G. Knepley We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 976d69c5d34SMatthew G. Knepley like a GPU, or vectorize on a multicore machine. 977d69c5d34SMatthew G. Knepley 978d69c5d34SMatthew G. Knepley Level: developer 979d69c5d34SMatthew G. Knepley 980d69c5d34SMatthew G. Knepley .seealso: DMPlexComputeJacobianFEM() 981d69c5d34SMatthew G. Knepley @*/ 982934789fcSMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user) 983d69c5d34SMatthew G. Knepley { 984d69c5d34SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dmc->data; 985d69c5d34SMatthew G. Knepley const char *name = "Interpolator"; 9862764a2aaSMatthew G. Knepley PetscDS prob; 987d69c5d34SMatthew G. Knepley PetscFE *feRef; 988d69c5d34SMatthew G. Knepley PetscSection fsection, fglobalSection; 989d69c5d34SMatthew G. Knepley PetscSection csection, cglobalSection; 990d69c5d34SMatthew G. Knepley PetscScalar *elemMat; 991942a7a06SMatthew G. Knepley PetscInt dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c; 9920f2d7e86SMatthew G. Knepley PetscInt cTotDim, rTotDim = 0; 993d69c5d34SMatthew G. Knepley PetscErrorCode ierr; 994d69c5d34SMatthew G. Knepley 995d69c5d34SMatthew G. Knepley PetscFunctionBegin; 996d69c5d34SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 997c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr); 998d69c5d34SMatthew G. Knepley ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 999d69c5d34SMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 1000d69c5d34SMatthew G. Knepley ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 1001d69c5d34SMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 1002d69c5d34SMatthew G. Knepley ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 1003d69c5d34SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 10042764a2aaSMatthew G. Knepley ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr); 1005d69c5d34SMatthew G. Knepley ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr); 1006d69c5d34SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 10070f2d7e86SMatthew G. Knepley PetscFE fe; 10080f2d7e86SMatthew G. Knepley PetscInt rNb, Nc; 1009d69c5d34SMatthew G. Knepley 10102764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 10110f2d7e86SMatthew G. Knepley ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr); 1012d69c5d34SMatthew G. Knepley ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr); 10130f2d7e86SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 10140f2d7e86SMatthew G. Knepley rTotDim += rNb*Nc; 1015d69c5d34SMatthew G. Knepley } 10162764a2aaSMatthew G. Knepley ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr); 10170f2d7e86SMatthew G. Knepley ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr); 10180f2d7e86SMatthew G. Knepley ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr); 1019d69c5d34SMatthew G. Knepley for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) { 1020d69c5d34SMatthew G. Knepley PetscDualSpace Qref; 1021d69c5d34SMatthew G. Knepley PetscQuadrature f; 1022d69c5d34SMatthew G. Knepley const PetscReal *qpoints, *qweights; 1023d69c5d34SMatthew G. Knepley PetscReal *points; 1024d69c5d34SMatthew G. Knepley PetscInt npoints = 0, Nc, Np, fpdim, i, k, p, d; 1025d69c5d34SMatthew G. Knepley 1026d69c5d34SMatthew G. Knepley /* Compose points from all dual basis functionals */ 1027d69c5d34SMatthew G. Knepley ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr); 10280f2d7e86SMatthew G. Knepley ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr); 1029d69c5d34SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr); 1030d69c5d34SMatthew G. Knepley for (i = 0; i < fpdim; ++i) { 1031d69c5d34SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1032d69c5d34SMatthew G. Knepley ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr); 1033d69c5d34SMatthew G. Knepley npoints += Np; 1034d69c5d34SMatthew G. Knepley } 1035d69c5d34SMatthew G. Knepley ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr); 1036d69c5d34SMatthew G. Knepley for (i = 0, k = 0; i < fpdim; ++i) { 1037d69c5d34SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1038d69c5d34SMatthew G. Knepley ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr); 1039d69c5d34SMatthew G. Knepley for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d]; 1040d69c5d34SMatthew G. Knepley } 1041d69c5d34SMatthew G. Knepley 1042d69c5d34SMatthew G. Knepley for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) { 10430f2d7e86SMatthew G. Knepley PetscFE fe; 1044d69c5d34SMatthew G. Knepley PetscReal *B; 104536a6d9c0SMatthew G. Knepley PetscInt NcJ, cpdim, j; 1046d69c5d34SMatthew G. Knepley 1047d69c5d34SMatthew G. Knepley /* Evaluate basis at points */ 10482764a2aaSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, fieldJ, (PetscObject *) &fe);CHKERRQ(ierr); 10490f2d7e86SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr); 10500f2d7e86SMatthew G. Knepley ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr); 1051ffe73a53SMatthew G. Knepley /* For now, fields only interpolate themselves */ 1052ffe73a53SMatthew G. Knepley if (fieldI == fieldJ) { 10537c927364SMatthew G. Knepley 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); 10540f2d7e86SMatthew G. Knepley ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr); 1055d69c5d34SMatthew G. Knepley for (i = 0, k = 0; i < fpdim; ++i) { 1056d69c5d34SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1057d69c5d34SMatthew G. Knepley ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr); 1058d69c5d34SMatthew G. Knepley for (p = 0; p < Np; ++p, ++k) { 105936a6d9c0SMatthew G. Knepley for (j = 0; j < cpdim; ++j) { 10600f2d7e86SMatthew G. Knepley 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]; 106136a6d9c0SMatthew G. Knepley } 1062d69c5d34SMatthew G. Knepley } 1063d69c5d34SMatthew G. Knepley } 10640f2d7e86SMatthew G. Knepley ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 1065ffe73a53SMatthew G. Knepley } 106636a6d9c0SMatthew G. Knepley offsetJ += cpdim*NcJ; 1067d69c5d34SMatthew G. Knepley } 1068d69c5d34SMatthew G. Knepley offsetI += fpdim*Nc; 1069549a8adaSMatthew G. Knepley ierr = PetscFree(points);CHKERRQ(ierr); 1070d69c5d34SMatthew G. Knepley } 10710f2d7e86SMatthew G. Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);} 10727f5b169aSMatthew G. Knepley /* Preallocate matrix */ 10737f5b169aSMatthew G. Knepley { 10747f5b169aSMatthew G. Knepley PetscHashJK ht; 10757f5b169aSMatthew G. Knepley PetscLayout rLayout; 10767f5b169aSMatthew G. Knepley PetscInt *dnz, *onz, *cellCIndices, *cellFIndices; 10777f5b169aSMatthew G. Knepley PetscInt locRows, rStart, rEnd, cell, r; 10787f5b169aSMatthew G. Knepley 10797f5b169aSMatthew G. Knepley ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr); 10807f5b169aSMatthew G. Knepley ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr); 10817f5b169aSMatthew G. Knepley ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr); 10827f5b169aSMatthew G. Knepley ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr); 10837f5b169aSMatthew G. Knepley ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr); 10847f5b169aSMatthew G. Knepley ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr); 10857f5b169aSMatthew G. Knepley ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr); 10867f5b169aSMatthew G. Knepley ierr = PetscCalloc4(locRows,&dnz,locRows,&onz,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr); 10877f5b169aSMatthew G. Knepley ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr); 10887f5b169aSMatthew G. Knepley for (cell = cStart; cell < cEnd; ++cell) { 10897f5b169aSMatthew G. Knepley ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr); 10907f5b169aSMatthew G. Knepley for (r = 0; r < rTotDim; ++r) { 10917f5b169aSMatthew G. Knepley PetscHashJKKey key; 10927f5b169aSMatthew G. Knepley PetscHashJKIter missing, iter; 10937f5b169aSMatthew G. Knepley 10947f5b169aSMatthew G. Knepley key.j = cellFIndices[r]; 10957f5b169aSMatthew G. Knepley if (key.j < 0) continue; 10967f5b169aSMatthew G. Knepley for (c = 0; c < cTotDim; ++c) { 10977f5b169aSMatthew G. Knepley key.k = cellCIndices[c]; 10987f5b169aSMatthew G. Knepley if (key.k < 0) continue; 10997f5b169aSMatthew G. Knepley ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr); 11007f5b169aSMatthew G. Knepley if (missing) { 11017f5b169aSMatthew G. Knepley ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr); 11027f5b169aSMatthew G. Knepley if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart]; 11037f5b169aSMatthew G. Knepley else ++onz[key.j-rStart]; 11047f5b169aSMatthew G. Knepley } 11057f5b169aSMatthew G. Knepley } 11067f5b169aSMatthew G. Knepley } 11077f5b169aSMatthew G. Knepley } 11087f5b169aSMatthew G. Knepley ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr); 11097f5b169aSMatthew G. Knepley ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr); 11107f5b169aSMatthew G. Knepley ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 11117f5b169aSMatthew G. Knepley ierr = PetscFree4(dnz,onz,cellCIndices,cellFIndices);CHKERRQ(ierr); 11127f5b169aSMatthew G. Knepley } 11137f5b169aSMatthew G. Knepley /* Fill matrix */ 11147f5b169aSMatthew G. Knepley ierr = MatZeroEntries(In);CHKERRQ(ierr); 1115d69c5d34SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 1116934789fcSMatthew G. Knepley ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr); 1117d69c5d34SMatthew G. Knepley } 1118549a8adaSMatthew G. Knepley for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);} 1119d69c5d34SMatthew G. Knepley ierr = PetscFree(feRef);CHKERRQ(ierr); 1120549a8adaSMatthew G. Knepley ierr = PetscFree(elemMat);CHKERRQ(ierr); 1121934789fcSMatthew G. Knepley ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1122934789fcSMatthew G. Knepley ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1123d69c5d34SMatthew G. Knepley if (mesh->printFEM) { 1124d69c5d34SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 1125934789fcSMatthew G. Knepley ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr); 1126934789fcSMatthew G. Knepley ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1127d69c5d34SMatthew G. Knepley } 1128d69c5d34SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1129d69c5d34SMatthew G. Knepley PetscFunctionReturn(0); 1130d69c5d34SMatthew G. Knepley } 11316c73c22cSMatthew G. Knepley 11326c73c22cSMatthew G. Knepley #undef __FUNCT__ 11337c927364SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInjectorFEM" 11347c927364SMatthew G. Knepley PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user) 11357c927364SMatthew G. Knepley { 1136e9d4ef1bSMatthew G. Knepley PetscDS prob; 11377c927364SMatthew G. Knepley PetscFE *feRef; 11387c927364SMatthew G. Knepley Vec fv, cv; 11397c927364SMatthew G. Knepley IS fis, cis; 11407c927364SMatthew G. Knepley PetscSection fsection, fglobalSection, csection, cglobalSection; 11417c927364SMatthew G. Knepley PetscInt *cmap, *cellCIndices, *cellFIndices, *cindices, *findices; 11427c927364SMatthew G. Knepley PetscInt cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, c, dim, d, startC, offsetC, offsetF, m; 11437c927364SMatthew G. Knepley PetscErrorCode ierr; 11447c927364SMatthew G. Knepley 11457c927364SMatthew G. Knepley PetscFunctionBegin; 114675a69067SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1147c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr); 11487c927364SMatthew G. Knepley ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 11497c927364SMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 11507c927364SMatthew G. Knepley ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 11517c927364SMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 11527c927364SMatthew G. Knepley ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 11537c927364SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 1154e9d4ef1bSMatthew G. Knepley ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr); 11557c927364SMatthew G. Knepley ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr); 11567c927364SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 11577c927364SMatthew G. Knepley PetscFE fe; 11587c927364SMatthew G. Knepley PetscInt fNb, Nc; 11597c927364SMatthew G. Knepley 1160e9d4ef1bSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, f, (PetscObject *) &fe);CHKERRQ(ierr); 11617c927364SMatthew G. Knepley ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr); 11627c927364SMatthew G. Knepley ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr); 11637c927364SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 11647c927364SMatthew G. Knepley fTotDim += fNb*Nc; 11657c927364SMatthew G. Knepley } 1166e9d4ef1bSMatthew G. Knepley ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr); 11677c927364SMatthew G. Knepley ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr); 11687c927364SMatthew G. Knepley for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) { 11697c927364SMatthew G. Knepley PetscFE feC; 11707c927364SMatthew G. Knepley PetscDualSpace QF, QC; 11717c927364SMatthew G. Knepley PetscInt NcF, NcC, fpdim, cpdim; 11727c927364SMatthew G. Knepley 1173e9d4ef1bSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr); 11747c927364SMatthew G. Knepley ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr); 11757c927364SMatthew G. Knepley ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr); 11767c927364SMatthew G. Knepley 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); 11777c927364SMatthew G. Knepley ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr); 11787c927364SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr); 11797c927364SMatthew G. Knepley ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr); 11807c927364SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr); 11817c927364SMatthew G. Knepley for (c = 0; c < cpdim; ++c) { 11827c927364SMatthew G. Knepley PetscQuadrature cfunc; 11837c927364SMatthew G. Knepley const PetscReal *cqpoints; 11847c927364SMatthew G. Knepley PetscInt NpC; 11857c927364SMatthew G. Knepley 11867c927364SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr); 11877c927364SMatthew G. Knepley ierr = PetscQuadratureGetData(cfunc, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr); 11887c927364SMatthew G. Knepley if (NpC != 1) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments"); 11897c927364SMatthew G. Knepley for (f = 0; f < fpdim; ++f) { 11907c927364SMatthew G. Knepley PetscQuadrature ffunc; 11917c927364SMatthew G. Knepley const PetscReal *fqpoints; 11927c927364SMatthew G. Knepley PetscReal sum = 0.0; 11937c927364SMatthew G. Knepley PetscInt NpF, comp; 11947c927364SMatthew G. Knepley 11957c927364SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr); 11967c927364SMatthew G. Knepley ierr = PetscQuadratureGetData(ffunc, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr); 11977c927364SMatthew G. Knepley if (NpC != NpF) continue; 11987c927364SMatthew G. Knepley for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]); 11997c927364SMatthew G. Knepley if (sum > 1.0e-9) continue; 12007c927364SMatthew G. Knepley for (comp = 0; comp < NcC; ++comp) { 12017c927364SMatthew G. Knepley cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp; 12027c927364SMatthew G. Knepley } 12037c927364SMatthew G. Knepley break; 12047c927364SMatthew G. Knepley } 12057c927364SMatthew G. Knepley } 12067c927364SMatthew G. Knepley offsetC += cpdim*NcC; 12077c927364SMatthew G. Knepley offsetF += fpdim*NcF; 12087c927364SMatthew G. Knepley } 12097c927364SMatthew G. Knepley for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);} 12107c927364SMatthew G. Knepley ierr = PetscFree(feRef);CHKERRQ(ierr); 12117c927364SMatthew G. Knepley 12127c927364SMatthew G. Knepley ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr); 12137c927364SMatthew G. Knepley ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr); 12147c927364SMatthew G. Knepley ierr = VecGetOwnershipRange(cv, &startC, NULL);CHKERRQ(ierr); 12157c927364SMatthew G. Knepley ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr); 12167c927364SMatthew G. Knepley ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr); 1217aa7da3c4SMatthew G. Knepley ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr); 1218aa7da3c4SMatthew G. Knepley ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr); 12197c927364SMatthew G. Knepley for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1; 12207c927364SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 12217c927364SMatthew G. Knepley ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr); 12227c927364SMatthew G. Knepley for (d = 0; d < cTotDim; ++d) { 12237c927364SMatthew G. Knepley if (cellCIndices[d] < 0) continue; 12247c927364SMatthew G. Knepley 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]]); 12257c927364SMatthew G. Knepley cindices[cellCIndices[d]-startC] = cellCIndices[d]; 12267c927364SMatthew G. Knepley findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]]; 12277c927364SMatthew G. Knepley } 12287c927364SMatthew G. Knepley } 12297c927364SMatthew G. Knepley ierr = PetscFree(cmap);CHKERRQ(ierr); 12307c927364SMatthew G. Knepley ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr); 12317c927364SMatthew G. Knepley 12327c927364SMatthew G. Knepley ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr); 12337c927364SMatthew G. Knepley ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr); 12347c927364SMatthew G. Knepley ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr); 12357c927364SMatthew G. Knepley ierr = ISDestroy(&cis);CHKERRQ(ierr); 12367c927364SMatthew G. Knepley ierr = ISDestroy(&fis);CHKERRQ(ierr); 12377c927364SMatthew G. Knepley ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr); 12387c927364SMatthew G. Knepley ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr); 123975a69067SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1240cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 1241cb1e1211SMatthew G Knepley } 1242