1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2da97024aSMatthew G. Knepley #include <petscsf.h> 3cb1e1211SMatthew G Knepley 4af0996ceSBarry Smith #include <petsc/private/petscfeimpl.h> 5af0996ceSBarry Smith #include <petsc/private/petscfvimpl.h> 6a0845e3aSMatthew G. Knepley 7cb1e1211SMatthew G Knepley PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale) 8cb1e1211SMatthew G Knepley { 9cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 10cb1e1211SMatthew G Knepley 11cb1e1211SMatthew G Knepley PetscFunctionBegin; 12cb1e1211SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 13cb1e1211SMatthew G Knepley PetscValidPointer(scale, 3); 14cb1e1211SMatthew G Knepley *scale = mesh->scale[unit]; 15cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 16cb1e1211SMatthew G Knepley } 17cb1e1211SMatthew G Knepley 18cb1e1211SMatthew G Knepley PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale) 19cb1e1211SMatthew G Knepley { 20cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 21cb1e1211SMatthew G Knepley 22cb1e1211SMatthew G Knepley PetscFunctionBegin; 23cb1e1211SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 24cb1e1211SMatthew G Knepley mesh->scale[unit] = scale; 25cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 26cb1e1211SMatthew G Knepley } 27cb1e1211SMatthew G Knepley 28cb1e1211SMatthew G Knepley PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k) 29cb1e1211SMatthew G Knepley { 30cb1e1211SMatthew G Knepley switch (i) { 31cb1e1211SMatthew G Knepley case 0: 32cb1e1211SMatthew G Knepley switch (j) { 33cb1e1211SMatthew G Knepley case 0: return 0; 34cb1e1211SMatthew G Knepley case 1: 35cb1e1211SMatthew G Knepley switch (k) { 36cb1e1211SMatthew G Knepley case 0: return 0; 37cb1e1211SMatthew G Knepley case 1: return 0; 38cb1e1211SMatthew G Knepley case 2: return 1; 39cb1e1211SMatthew G Knepley } 40cb1e1211SMatthew G Knepley case 2: 41cb1e1211SMatthew G Knepley switch (k) { 42cb1e1211SMatthew G Knepley case 0: return 0; 43cb1e1211SMatthew G Knepley case 1: return -1; 44cb1e1211SMatthew G Knepley case 2: return 0; 45cb1e1211SMatthew G Knepley } 46cb1e1211SMatthew G Knepley } 47cb1e1211SMatthew G Knepley case 1: 48cb1e1211SMatthew G Knepley switch (j) { 49cb1e1211SMatthew G Knepley case 0: 50cb1e1211SMatthew G Knepley switch (k) { 51cb1e1211SMatthew G Knepley case 0: return 0; 52cb1e1211SMatthew G Knepley case 1: return 0; 53cb1e1211SMatthew G Knepley case 2: return -1; 54cb1e1211SMatthew G Knepley } 55cb1e1211SMatthew G Knepley case 1: return 0; 56cb1e1211SMatthew G Knepley case 2: 57cb1e1211SMatthew G Knepley switch (k) { 58cb1e1211SMatthew G Knepley case 0: return 1; 59cb1e1211SMatthew G Knepley case 1: return 0; 60cb1e1211SMatthew G Knepley case 2: return 0; 61cb1e1211SMatthew G Knepley } 62cb1e1211SMatthew G Knepley } 63cb1e1211SMatthew G Knepley case 2: 64cb1e1211SMatthew G Knepley switch (j) { 65cb1e1211SMatthew G Knepley case 0: 66cb1e1211SMatthew G Knepley switch (k) { 67cb1e1211SMatthew G Knepley case 0: return 0; 68cb1e1211SMatthew G Knepley case 1: return 1; 69cb1e1211SMatthew G Knepley case 2: return 0; 70cb1e1211SMatthew G Knepley } 71cb1e1211SMatthew G Knepley case 1: 72cb1e1211SMatthew G Knepley switch (k) { 73cb1e1211SMatthew G Knepley case 0: return -1; 74cb1e1211SMatthew G Knepley case 1: return 0; 75cb1e1211SMatthew G Knepley case 2: return 0; 76cb1e1211SMatthew G Knepley } 77cb1e1211SMatthew G Knepley case 2: return 0; 78cb1e1211SMatthew G Knepley } 79cb1e1211SMatthew G Knepley } 80cb1e1211SMatthew G Knepley return 0; 81cb1e1211SMatthew G Knepley } 82cb1e1211SMatthew G Knepley 83d1828a1cSMatthew G. Knepley static PetscErrorCode DMPlexProjectRigidBody_Private(PetscInt dim, PetscReal t, const PetscReal X[], PetscInt Nf, PetscScalar *mode, void *ctx) 84026175e5SToby Isaac { 85026175e5SToby Isaac PetscInt *ctxInt = (PetscInt *) ctx; 86ad917190SMatthew G. Knepley PetscInt dim2 = ctxInt[0]; 87026175e5SToby Isaac PetscInt d = ctxInt[1]; 88026175e5SToby Isaac PetscInt i, j, k = dim > 2 ? d - dim : d; 89026175e5SToby Isaac 90ad917190SMatthew G. Knepley PetscFunctionBegin; 91ad917190SMatthew G. Knepley if (dim != dim2) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Input dimension %d does not match context dimension %d", dim, dim2); 92026175e5SToby Isaac for (i = 0; i < dim; i++) mode[i] = 0.; 93026175e5SToby Isaac if (d < dim) { 94026175e5SToby Isaac mode[d] = 1.; 95ad917190SMatthew G. Knepley } else { 96026175e5SToby Isaac for (i = 0; i < dim; i++) { 97026175e5SToby Isaac for (j = 0; j < dim; j++) { 98026175e5SToby Isaac mode[j] += epsilon(i, j, k)*X[i]; 99026175e5SToby Isaac } 100026175e5SToby Isaac } 101026175e5SToby Isaac } 102ad917190SMatthew G. Knepley PetscFunctionReturn(0); 103026175e5SToby Isaac } 104026175e5SToby Isaac 105cb1e1211SMatthew G Knepley /*@C 106026175e5SToby Isaac DMPlexCreateRigidBody - for the default global section, create rigid body modes from coordinates 107cb1e1211SMatthew G Knepley 108cb1e1211SMatthew G Knepley Collective on DM 109cb1e1211SMatthew G Knepley 110cb1e1211SMatthew G Knepley Input Arguments: 111026175e5SToby Isaac . dm - the DM 112cb1e1211SMatthew G Knepley 113cb1e1211SMatthew G Knepley Output Argument: 114cb1e1211SMatthew G Knepley . sp - the null space 115cb1e1211SMatthew G Knepley 116cb1e1211SMatthew G Knepley Note: This is necessary to take account of Dirichlet conditions on the displacements 117cb1e1211SMatthew G Knepley 118cb1e1211SMatthew G Knepley Level: advanced 119cb1e1211SMatthew G Knepley 120cb1e1211SMatthew G Knepley .seealso: MatNullSpaceCreate() 121cb1e1211SMatthew G Knepley @*/ 122026175e5SToby Isaac PetscErrorCode DMPlexCreateRigidBody(DM dm, MatNullSpace *sp) 123cb1e1211SMatthew G Knepley { 124cb1e1211SMatthew G Knepley MPI_Comm comm; 125026175e5SToby Isaac Vec mode[6]; 126026175e5SToby Isaac PetscSection section, globalSection; 1279d8fbdeaSMatthew G. Knepley PetscInt dim, dimEmbed, n, m, d, i, j; 128cb1e1211SMatthew G Knepley PetscErrorCode ierr; 129cb1e1211SMatthew G Knepley 130cb1e1211SMatthew G Knepley PetscFunctionBegin; 131cb1e1211SMatthew G Knepley ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 132c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1339d8fbdeaSMatthew G. Knepley ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 134cb1e1211SMatthew G Knepley if (dim == 1) { 135cb1e1211SMatthew G Knepley ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr); 136cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 137cb1e1211SMatthew G Knepley } 138026175e5SToby Isaac ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 139026175e5SToby Isaac ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr); 140cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr); 141cb1e1211SMatthew G Knepley m = (dim*(dim+1))/2; 142cb1e1211SMatthew G Knepley ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr); 143cb1e1211SMatthew G Knepley ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr); 144cb1e1211SMatthew G Knepley ierr = VecSetUp(mode[0]);CHKERRQ(ierr); 145cb1e1211SMatthew G Knepley for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);} 146026175e5SToby Isaac for (d = 0; d < m; d++) { 147330485fdSToby Isaac PetscInt ctx[2]; 148d1828a1cSMatthew G. Knepley PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal *, PetscInt, PetscScalar *, void *) = DMPlexProjectRigidBody_Private; 149026175e5SToby Isaac void *voidctx = (void *) (&ctx[0]); 150cb1e1211SMatthew G Knepley 1519d8fbdeaSMatthew G. Knepley ctx[0] = dimEmbed; 152330485fdSToby Isaac ctx[1] = d; 1530709b2feSToby Isaac ierr = DMProjectFunction(dm, 0.0, &func, &voidctx, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 154cb1e1211SMatthew G Knepley } 155cb1e1211SMatthew G Knepley for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);} 156cb1e1211SMatthew G Knepley /* Orthonormalize system */ 157cb1e1211SMatthew G Knepley for (i = dim; i < m; ++i) { 158cb1e1211SMatthew G Knepley PetscScalar dots[6]; 159cb1e1211SMatthew G Knepley 160cb1e1211SMatthew G Knepley ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr); 161cb1e1211SMatthew G Knepley for (j = 0; j < i; ++j) dots[j] *= -1.0; 162cb1e1211SMatthew G Knepley ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr); 163cb1e1211SMatthew G Knepley ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr); 164cb1e1211SMatthew G Knepley } 165cb1e1211SMatthew G Knepley ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr); 166cb1e1211SMatthew G Knepley for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);} 167cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 168cb1e1211SMatthew G Knepley } 169cb1e1211SMatthew G Knepley 170b29cfa1cSToby Isaac /*@ 171b29cfa1cSToby Isaac DMPlexSetMaxProjectionHeight - In DMPlexProjectXXXLocal() functions, the projected values of a basis function's dofs 172b29cfa1cSToby Isaac are computed by associating the basis function with one of the mesh points in its transitively-closed support, and 173b29cfa1cSToby Isaac evaluating the dual space basis of that point. A basis function is associated with the point in its 174b29cfa1cSToby Isaac transitively-closed support whose mesh height is highest (w.r.t. DAG height), but not greater than the maximum 175b29cfa1cSToby Isaac projection height, which is set with this function. By default, the maximum projection height is zero, which means 176b29cfa1cSToby Isaac that only mesh cells are used to project basis functions. A height of one, for example, evaluates a cell-interior 177b29cfa1cSToby Isaac basis functions using its cells dual space basis, but all other basis functions with the dual space basis of a face. 178b29cfa1cSToby Isaac 179b29cfa1cSToby Isaac Input Parameters: 180b29cfa1cSToby Isaac + dm - the DMPlex object 181b29cfa1cSToby Isaac - height - the maximum projection height >= 0 182b29cfa1cSToby Isaac 183b29cfa1cSToby Isaac Level: advanced 184b29cfa1cSToby Isaac 1854d6f44ffSToby Isaac .seealso: DMPlexGetMaxProjectionHeight(), DMProjectFunctionLocal(), DMProjectFunctionLabelLocal() 186b29cfa1cSToby Isaac @*/ 187b29cfa1cSToby Isaac PetscErrorCode DMPlexSetMaxProjectionHeight(DM dm, PetscInt height) 188b29cfa1cSToby Isaac { 189b29cfa1cSToby Isaac DM_Plex *plex = (DM_Plex *) dm->data; 190b29cfa1cSToby Isaac 191b29cfa1cSToby Isaac PetscFunctionBegin; 192b29cfa1cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 193b29cfa1cSToby Isaac plex->maxProjectionHeight = height; 194b29cfa1cSToby Isaac PetscFunctionReturn(0); 195b29cfa1cSToby Isaac } 196b29cfa1cSToby Isaac 197b29cfa1cSToby Isaac /*@ 198b29cfa1cSToby Isaac DMPlexGetMaxProjectionHeight - Get the maximum height (w.r.t. DAG) of mesh points used to evaluate dual bases in 199b29cfa1cSToby Isaac DMPlexProjectXXXLocal() functions. 200b29cfa1cSToby Isaac 201b29cfa1cSToby Isaac Input Parameters: 202b29cfa1cSToby Isaac . dm - the DMPlex object 203b29cfa1cSToby Isaac 204b29cfa1cSToby Isaac Output Parameters: 205b29cfa1cSToby Isaac . height - the maximum projection height 206b29cfa1cSToby Isaac 207b29cfa1cSToby Isaac Level: intermediate 208b29cfa1cSToby Isaac 2094d6f44ffSToby Isaac .seealso: DMPlexSetMaxProjectionHeight(), DMProjectFunctionLocal(), DMProjectFunctionLabelLocal() 210b29cfa1cSToby Isaac @*/ 211b29cfa1cSToby Isaac PetscErrorCode DMPlexGetMaxProjectionHeight(DM dm, PetscInt *height) 212b29cfa1cSToby Isaac { 213b29cfa1cSToby Isaac DM_Plex *plex = (DM_Plex *) dm->data; 214b29cfa1cSToby Isaac 215b29cfa1cSToby Isaac PetscFunctionBegin; 216b29cfa1cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 217b29cfa1cSToby Isaac *height = plex->maxProjectionHeight; 218b29cfa1cSToby Isaac PetscFunctionReturn(0); 219b29cfa1cSToby Isaac } 220b29cfa1cSToby Isaac 2210163fd50SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FEM_Internal(DM dm, PetscReal time, PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *ctx, Vec locX) 22255f2e967SMatthew G. Knepley { 2230163fd50SMatthew G. Knepley PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal x[], PetscInt, PetscScalar *u, void *ctx); 22455f2e967SMatthew G. Knepley void **ctxs; 225d7ddef95SMatthew G. Knepley PetscInt numFields; 226d7ddef95SMatthew G. Knepley PetscErrorCode ierr; 227d7ddef95SMatthew G. Knepley 228d7ddef95SMatthew G. Knepley PetscFunctionBegin; 229d7ddef95SMatthew G. Knepley ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr); 230d7ddef95SMatthew G. Knepley ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr); 231d7ddef95SMatthew G. Knepley funcs[field] = func; 232d7ddef95SMatthew G. Knepley ctxs[field] = ctx; 2330709b2feSToby Isaac ierr = DMProjectFunctionLabelLocal(dm, time, label, numids, ids, funcs, ctxs, INSERT_BC_VALUES, locX);CHKERRQ(ierr); 234d7ddef95SMatthew G. Knepley ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr); 235d7ddef95SMatthew G. Knepley PetscFunctionReturn(0); 236d7ddef95SMatthew G. Knepley } 237d7ddef95SMatthew G. Knepley 238c60e475cSMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FEM_AuxField_Internal(DM dm, PetscReal time, Vec locU, PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], 239c60e475cSMatthew G. Knepley void (*func)(PetscInt, PetscInt, PetscInt, 240c60e475cSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 241c60e475cSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 242c60e475cSMatthew G. Knepley PetscReal, const PetscReal[], PetscScalar[]), void *ctx, Vec locX) 243c60e475cSMatthew G. Knepley { 244c60e475cSMatthew G. Knepley void (**funcs)(PetscInt, PetscInt, PetscInt, 245c60e475cSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 246c60e475cSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 247c60e475cSMatthew G. Knepley PetscReal, const PetscReal[], PetscScalar[]); 248c60e475cSMatthew G. Knepley void **ctxs; 249c60e475cSMatthew G. Knepley PetscInt numFields; 250c60e475cSMatthew G. Knepley PetscErrorCode ierr; 251c60e475cSMatthew G. Knepley 252c60e475cSMatthew G. Knepley PetscFunctionBegin; 253c60e475cSMatthew G. Knepley ierr = DMGetNumFields(dm, &numFields);CHKERRQ(ierr); 254c60e475cSMatthew G. Knepley ierr = PetscCalloc2(numFields,&funcs,numFields,&ctxs);CHKERRQ(ierr); 255c60e475cSMatthew G. Knepley funcs[field] = func; 256c60e475cSMatthew G. Knepley ctxs[field] = ctx; 257c60e475cSMatthew G. Knepley ierr = DMProjectFieldLabelLocal(dm, time, label, numids, ids, locU, funcs, INSERT_BC_VALUES, locX);CHKERRQ(ierr); 258c60e475cSMatthew G. Knepley ierr = PetscFree2(funcs,ctxs);CHKERRQ(ierr); 259c60e475cSMatthew G. Knepley PetscFunctionReturn(0); 260c60e475cSMatthew G. Knepley } 261c60e475cSMatthew G. Knepley 2620e0f25f2SMatthew G. Knepley /* This ignores numcomps/comps */ 263d7ddef95SMatthew G. Knepley static PetscErrorCode DMPlexInsertBoundaryValues_FVM_Internal(DM dm, PetscReal time, Vec faceGeometry, Vec cellGeometry, Vec Grad, 264d7ddef95SMatthew G. Knepley PetscInt field, DMLabel label, PetscInt numids, const PetscInt ids[], PetscErrorCode (*func)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*), void *ctx, Vec locX) 265d7ddef95SMatthew G. Knepley { 26661f58d28SMatthew G. Knepley PetscDS prob; 267da97024aSMatthew G. Knepley PetscSF sf; 268d7ddef95SMatthew G. Knepley DM dmFace, dmCell, dmGrad; 26920369375SToby Isaac const PetscScalar *facegeom, *cellgeom = NULL, *grad; 270da97024aSMatthew G. Knepley const PetscInt *leaves; 271d7ddef95SMatthew G. Knepley PetscScalar *x, *fx; 272520b3818SMatthew G. Knepley PetscInt dim, nleaves, loc, fStart, fEnd, pdim, i; 273e735a8a9SMatthew G. Knepley PetscErrorCode ierr, ierru = 0; 274d7ddef95SMatthew G. Knepley 275d7ddef95SMatthew G. Knepley PetscFunctionBegin; 276da97024aSMatthew G. Knepley ierr = DMGetPointSF(dm, &sf);CHKERRQ(ierr); 277da97024aSMatthew G. Knepley ierr = PetscSFGetGraph(sf, NULL, &nleaves, &leaves, NULL);CHKERRQ(ierr); 278da97024aSMatthew G. Knepley nleaves = PetscMax(0, nleaves); 279d7ddef95SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 280d7ddef95SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 28161f58d28SMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 282d7ddef95SMatthew G. Knepley ierr = VecGetDM(faceGeometry, &dmFace);CHKERRQ(ierr); 283d7ddef95SMatthew G. Knepley ierr = VecGetArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 28420369375SToby Isaac if (cellGeometry) { 285d7ddef95SMatthew G. Knepley ierr = VecGetDM(cellGeometry, &dmCell);CHKERRQ(ierr); 286d7ddef95SMatthew G. Knepley ierr = VecGetArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr); 28720369375SToby Isaac } 288d7ddef95SMatthew G. Knepley if (Grad) { 289c0a6632aSMatthew G. Knepley PetscFV fv; 290c0a6632aSMatthew G. Knepley 291c0a6632aSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fv);CHKERRQ(ierr); 292d7ddef95SMatthew G. Knepley ierr = VecGetDM(Grad, &dmGrad);CHKERRQ(ierr); 293d7ddef95SMatthew G. Knepley ierr = VecGetArrayRead(Grad, &grad);CHKERRQ(ierr); 294d7ddef95SMatthew G. Knepley ierr = PetscFVGetNumComponents(fv, &pdim);CHKERRQ(ierr); 295d7ddef95SMatthew G. Knepley ierr = DMGetWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr); 296d7ddef95SMatthew G. Knepley } 297d7ddef95SMatthew G. Knepley ierr = VecGetArray(locX, &x);CHKERRQ(ierr); 298d7ddef95SMatthew G. Knepley for (i = 0; i < numids; ++i) { 299d7ddef95SMatthew G. Knepley IS faceIS; 300d7ddef95SMatthew G. Knepley const PetscInt *faces; 301d7ddef95SMatthew G. Knepley PetscInt numFaces, f; 302d7ddef95SMatthew G. Knepley 303d7ddef95SMatthew G. Knepley ierr = DMLabelGetStratumIS(label, ids[i], &faceIS);CHKERRQ(ierr); 304d7ddef95SMatthew G. Knepley if (!faceIS) continue; /* No points with that id on this process */ 305d7ddef95SMatthew G. Knepley ierr = ISGetLocalSize(faceIS, &numFaces);CHKERRQ(ierr); 306d7ddef95SMatthew G. Knepley ierr = ISGetIndices(faceIS, &faces);CHKERRQ(ierr); 307d7ddef95SMatthew G. Knepley for (f = 0; f < numFaces; ++f) { 308d7ddef95SMatthew G. Knepley const PetscInt face = faces[f], *cells; 309640bce14SSatish Balay PetscFVFaceGeom *fg; 310d7ddef95SMatthew G. Knepley 311d7ddef95SMatthew G. Knepley if ((face < fStart) || (face >= fEnd)) continue; /* Refinement adds non-faces to labels */ 312da97024aSMatthew G. Knepley ierr = PetscFindInt(face, nleaves, (PetscInt *) leaves, &loc);CHKERRQ(ierr); 313da97024aSMatthew G. Knepley if (loc >= 0) continue; 314d7ddef95SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmFace, face, facegeom, &fg);CHKERRQ(ierr); 315d7ddef95SMatthew G. Knepley ierr = DMPlexGetSupport(dm, face, &cells);CHKERRQ(ierr); 316d7ddef95SMatthew G. Knepley if (Grad) { 317640bce14SSatish Balay PetscFVCellGeom *cg; 318640bce14SSatish Balay PetscScalar *cx, *cgrad; 319d7ddef95SMatthew G. Knepley PetscScalar *xG; 320d7ddef95SMatthew G. Knepley PetscReal dx[3]; 321d7ddef95SMatthew G. Knepley PetscInt d; 322d7ddef95SMatthew G. Knepley 323d7ddef95SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmCell, cells[0], cellgeom, &cg);CHKERRQ(ierr); 324d7ddef95SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[0], x, &cx);CHKERRQ(ierr); 325d7ddef95SMatthew G. Knepley ierr = DMPlexPointLocalRead(dmGrad, cells[0], grad, &cgrad);CHKERRQ(ierr); 326520b3818SMatthew G. Knepley ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr); 327d7ddef95SMatthew G. Knepley DMPlex_WaxpyD_Internal(dim, -1, cg->centroid, fg->centroid, dx); 328d7ddef95SMatthew G. Knepley for (d = 0; d < pdim; ++d) fx[d] = cx[d] + DMPlex_DotD_Internal(dim, &cgrad[d*dim], dx); 329e735a8a9SMatthew G. Knepley ierru = (*func)(time, fg->centroid, fg->normal, fx, xG, ctx); 330e735a8a9SMatthew G. Knepley if (ierru) { 331e735a8a9SMatthew G. Knepley ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 332e735a8a9SMatthew G. Knepley ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 333e735a8a9SMatthew G. Knepley goto cleanup; 334e735a8a9SMatthew G. Knepley } 335d7ddef95SMatthew G. Knepley } else { 336640bce14SSatish Balay PetscScalar *xI; 337d7ddef95SMatthew G. Knepley PetscScalar *xG; 338d7ddef95SMatthew G. Knepley 339d7ddef95SMatthew G. Knepley ierr = DMPlexPointLocalRead(dm, cells[0], x, &xI);CHKERRQ(ierr); 340520b3818SMatthew G. Knepley ierr = DMPlexPointLocalFieldRef(dm, cells[1], field, x, &xG);CHKERRQ(ierr); 341e735a8a9SMatthew G. Knepley ierru = (*func)(time, fg->centroid, fg->normal, xI, xG, ctx); 342e735a8a9SMatthew G. Knepley if (ierru) { 343e735a8a9SMatthew G. Knepley ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 344e735a8a9SMatthew G. Knepley ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 345e735a8a9SMatthew G. Knepley goto cleanup; 346e735a8a9SMatthew G. Knepley } 347d7ddef95SMatthew G. Knepley } 348d7ddef95SMatthew G. Knepley } 349d7ddef95SMatthew G. Knepley ierr = ISRestoreIndices(faceIS, &faces);CHKERRQ(ierr); 350d7ddef95SMatthew G. Knepley ierr = ISDestroy(&faceIS);CHKERRQ(ierr); 351d7ddef95SMatthew G. Knepley } 352e735a8a9SMatthew G. Knepley cleanup: 353d7ddef95SMatthew G. Knepley ierr = VecRestoreArray(locX, &x);CHKERRQ(ierr); 354d7ddef95SMatthew G. Knepley if (Grad) { 355d7ddef95SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, pdim, PETSC_SCALAR, &fx);CHKERRQ(ierr); 356d7ddef95SMatthew G. Knepley ierr = VecRestoreArrayRead(Grad, &grad);CHKERRQ(ierr); 357d7ddef95SMatthew G. Knepley } 358e735a8a9SMatthew G. Knepley if (cellGeometry) {ierr = VecRestoreArrayRead(cellGeometry, &cellgeom);CHKERRQ(ierr);} 359d7ddef95SMatthew G. Knepley ierr = VecRestoreArrayRead(faceGeometry, &facegeom);CHKERRQ(ierr); 360e735a8a9SMatthew G. Knepley CHKERRQ(ierru); 361d7ddef95SMatthew G. Knepley PetscFunctionReturn(0); 362d7ddef95SMatthew G. Knepley } 363d7ddef95SMatthew G. Knepley 3649d24f7bbSMatthew G. Knepley /*@ 3659d24f7bbSMatthew G. Knepley DMPlexInsertBoundaryValues - Puts coefficients which represent boundary values into the local solution vector 3669d24f7bbSMatthew G. Knepley 3679d24f7bbSMatthew G. Knepley Input Parameters: 3689d24f7bbSMatthew G. Knepley + dm - The DM 3699d24f7bbSMatthew G. Knepley . insertEssential - Should I insert essential (e.g. Dirichlet) or inessential (e.g. Neumann) boundary conditions 3709d24f7bbSMatthew G. Knepley . time - The time 3719d24f7bbSMatthew G. Knepley . faceGeomFVM - Face geometry data for FV discretizations 3729d24f7bbSMatthew G. Knepley . cellGeomFVM - Cell geometry data for FV discretizations 3739d24f7bbSMatthew G. Knepley - gradFVM - Gradient reconstruction data for FV discretizations 3749d24f7bbSMatthew G. Knepley 3759d24f7bbSMatthew G. Knepley Output Parameters: 3769d24f7bbSMatthew G. Knepley . locX - Solution updated with boundary values 3779d24f7bbSMatthew G. Knepley 3789d24f7bbSMatthew G. Knepley Level: developer 3799d24f7bbSMatthew G. Knepley 3809d24f7bbSMatthew G. Knepley .seealso: DMProjectFunctionLabelLocal() 3819d24f7bbSMatthew G. Knepley @*/ 382bdd6f66aSToby Isaac PetscErrorCode DMPlexInsertBoundaryValues(DM dm, PetscBool insertEssential, Vec locX, PetscReal time, Vec faceGeomFVM, Vec cellGeomFVM, Vec gradFVM) 383d7ddef95SMatthew G. Knepley { 384d7ddef95SMatthew G. Knepley PetscInt numBd, b; 38555f2e967SMatthew G. Knepley PetscErrorCode ierr; 38655f2e967SMatthew G. Knepley 38755f2e967SMatthew G. Knepley PetscFunctionBegin; 38855f2e967SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 389d7ddef95SMatthew G. Knepley PetscValidHeaderSpecific(locX, VEC_CLASSID, 2); 390d7ddef95SMatthew G. Knepley if (faceGeomFVM) {PetscValidHeaderSpecific(faceGeomFVM, VEC_CLASSID, 4);} 391d7ddef95SMatthew G. Knepley if (cellGeomFVM) {PetscValidHeaderSpecific(cellGeomFVM, VEC_CLASSID, 5);} 392d7ddef95SMatthew G. Knepley if (gradFVM) {PetscValidHeaderSpecific(gradFVM, VEC_CLASSID, 6);} 393dff059c6SToby Isaac ierr = PetscDSGetNumBoundary(dm->prob, &numBd);CHKERRQ(ierr); 39455f2e967SMatthew G. Knepley for (b = 0; b < numBd; ++b) { 395f971fd6bSMatthew G. Knepley DMBoundaryConditionType type; 396d7ddef95SMatthew G. Knepley const char *labelname; 397d7ddef95SMatthew G. Knepley DMLabel label; 398d7ddef95SMatthew G. Knepley PetscInt field; 399d7ddef95SMatthew G. Knepley PetscObject obj; 400d7ddef95SMatthew G. Knepley PetscClassId id; 40155f2e967SMatthew G. Knepley void (*func)(); 402d7ddef95SMatthew G. Knepley PetscInt numids; 403d7ddef95SMatthew G. Knepley const PetscInt *ids; 40455f2e967SMatthew G. Knepley void *ctx; 40555f2e967SMatthew G. Knepley 406f971fd6bSMatthew G. Knepley ierr = DMGetBoundary(dm, b, &type, NULL, &labelname, &field, NULL, NULL, &func, &numids, &ids, &ctx);CHKERRQ(ierr); 407f971fd6bSMatthew G. Knepley if (insertEssential != (type & DM_BC_ESSENTIAL)) continue; 408c58f1c22SToby Isaac ierr = DMGetLabel(dm, labelname, &label);CHKERRQ(ierr); 409d7ddef95SMatthew G. Knepley ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 410d7ddef95SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 411d7ddef95SMatthew G. Knepley if (id == PETSCFE_CLASSID) { 412c60e475cSMatthew G. Knepley switch (type) { 413c60e475cSMatthew G. Knepley /* for FEM, there is no insertion to be done for non-essential boundary conditions */ 414c60e475cSMatthew G. Knepley case DM_BC_ESSENTIAL: 415092e5057SToby Isaac ierr = DMPlexLabelAddCells(dm,label);CHKERRQ(ierr); 4160163fd50SMatthew G. Knepley ierr = DMPlexInsertBoundaryValues_FEM_Internal(dm, time, field, label, numids, ids, (PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *)) func, ctx, locX);CHKERRQ(ierr); 417092e5057SToby Isaac ierr = DMPlexLabelClearCells(dm,label);CHKERRQ(ierr); 418c60e475cSMatthew G. Knepley break; 419c60e475cSMatthew G. Knepley case DM_BC_ESSENTIAL_FIELD: 420c60e475cSMatthew G. Knepley ierr = DMPlexLabelAddCells(dm,label);CHKERRQ(ierr); 421c60e475cSMatthew G. Knepley ierr = DMPlexInsertBoundaryValues_FEM_AuxField_Internal(dm, time, locX, field, label, numids, ids, (void (*)(PetscInt, PetscInt, PetscInt, 422c60e475cSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 423c60e475cSMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 424c60e475cSMatthew G. Knepley PetscReal, const PetscReal[], PetscScalar[])) func, ctx, locX);CHKERRQ(ierr); 425c60e475cSMatthew G. Knepley ierr = DMPlexLabelClearCells(dm,label);CHKERRQ(ierr); 426c60e475cSMatthew G. Knepley break; 427c60e475cSMatthew G. Knepley default: break; 428c60e475cSMatthew G. Knepley } 429d7ddef95SMatthew G. Knepley } else if (id == PETSCFV_CLASSID) { 43043ea7facSMatthew G. Knepley if (!faceGeomFVM) continue; 431d7ddef95SMatthew G. Knepley ierr = DMPlexInsertBoundaryValues_FVM_Internal(dm, time, faceGeomFVM, cellGeomFVM, gradFVM, 432d7ddef95SMatthew G. Knepley field, label, numids, ids, (PetscErrorCode (*)(PetscReal,const PetscReal*,const PetscReal*,const PetscScalar*,PetscScalar*,void*)) func, ctx, locX);CHKERRQ(ierr); 433d7ddef95SMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 43455f2e967SMatthew G. Knepley } 43555f2e967SMatthew G. Knepley PetscFunctionReturn(0); 43655f2e967SMatthew G. Knepley } 43755f2e967SMatthew G. Knepley 4380709b2feSToby Isaac PetscErrorCode DMComputeL2Diff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff) 439cb1e1211SMatthew G Knepley { 440cb1e1211SMatthew G Knepley const PetscInt debug = 0; 441cb1e1211SMatthew G Knepley PetscSection section; 442c5bbbd5bSMatthew G. Knepley PetscQuadrature quad; 443cb1e1211SMatthew G Knepley Vec localX; 44415496722SMatthew G. Knepley PetscScalar *funcVal, *interpolant; 4457318780aSToby Isaac PetscReal *coords, *detJ, *J; 446cb1e1211SMatthew G Knepley PetscReal localDiff = 0.0; 4477318780aSToby Isaac const PetscReal *quadWeights; 4489c3cf19fSMatthew G. Knepley PetscInt dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, cEndInterior, c, field, fieldOffset; 449cb1e1211SMatthew G Knepley PetscErrorCode ierr; 450cb1e1211SMatthew G Knepley 451cb1e1211SMatthew G Knepley PetscFunctionBegin; 452c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 4537318780aSToby Isaac ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr); 454cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 455cb1e1211SMatthew G Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 456cb1e1211SMatthew G Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 457bdd6f66aSToby Isaac ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 458cb1e1211SMatthew G Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 459cb1e1211SMatthew G Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 460cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 46115496722SMatthew G. Knepley PetscObject obj; 46215496722SMatthew G. Knepley PetscClassId id; 463c5bbbd5bSMatthew G. Knepley PetscInt Nc; 464c5bbbd5bSMatthew G. Knepley 46515496722SMatthew G. Knepley ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 46615496722SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 46715496722SMatthew G. Knepley if (id == PETSCFE_CLASSID) { 46815496722SMatthew G. Knepley PetscFE fe = (PetscFE) obj; 46915496722SMatthew G. Knepley 4700f2d7e86SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 4710f2d7e86SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 47215496722SMatthew G. Knepley } else if (id == PETSCFV_CLASSID) { 47315496722SMatthew G. Knepley PetscFV fv = (PetscFV) obj; 47415496722SMatthew G. Knepley 47515496722SMatthew G. Knepley ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 47615496722SMatthew G. Knepley ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 47715496722SMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 478c5bbbd5bSMatthew G. Knepley numComponents += Nc; 479cb1e1211SMatthew G Knepley } 4809c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, NULL, &quadWeights);CHKERRQ(ierr); 481*beaa55a6SMatthew G. Knepley if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents); 4827318780aSToby Isaac ierr = PetscMalloc5(numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J);CHKERRQ(ierr); 483cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 4849ac3fadcSMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 4859ac3fadcSMatthew G. Knepley cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 486cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 487a1e44745SMatthew G. Knepley PetscScalar *x = NULL; 488cb1e1211SMatthew G Knepley PetscReal elemDiff = 0.0; 4899c3cf19fSMatthew G. Knepley PetscInt qc = 0; 490cb1e1211SMatthew G Knepley 4917318780aSToby Isaac ierr = DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, NULL, detJ);CHKERRQ(ierr); 492cb1e1211SMatthew G Knepley ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 493cb1e1211SMatthew G Knepley 49415496722SMatthew G. Knepley for (field = 0, fieldOffset = 0; field < numFields; ++field) { 49515496722SMatthew G. Knepley PetscObject obj; 49615496722SMatthew G. Knepley PetscClassId id; 497c110b1eeSGeoffrey Irving void * const ctx = ctxs ? ctxs[field] : NULL; 49815496722SMatthew G. Knepley PetscInt Nb, Nc, q, fc; 499cb1e1211SMatthew G Knepley 50015496722SMatthew G. Knepley ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 50115496722SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 50215496722SMatthew G. Knepley if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 50315496722SMatthew G. Knepley else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 50415496722SMatthew G. Knepley else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 505cb1e1211SMatthew G Knepley if (debug) { 506cb1e1211SMatthew G Knepley char title[1024]; 507cb1e1211SMatthew G Knepley ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 50815496722SMatthew G. Knepley ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr); 509cb1e1211SMatthew G Knepley } 5107318780aSToby Isaac for (q = 0; q < Nq; ++q) { 5117318780aSToby Isaac if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, point %D", detJ[q], c, q); 5127318780aSToby Isaac ierr = (*funcs[field])(coordDim, time, &coords[coordDim * q], Nc, funcVal, ctx); 513e735a8a9SMatthew G. Knepley if (ierr) { 514e735a8a9SMatthew G. Knepley PetscErrorCode ierr2; 515e735a8a9SMatthew G. Knepley ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2); 516e735a8a9SMatthew G. Knepley ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2); 5177318780aSToby Isaac ierr2 = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr2); 518e735a8a9SMatthew G. Knepley CHKERRQ(ierr); 519e735a8a9SMatthew G. Knepley } 52015496722SMatthew G. Knepley if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 52115496722SMatthew G. Knepley else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 52215496722SMatthew G. Knepley else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 52315496722SMatthew G. Knepley for (fc = 0; fc < Nc; ++fc) { 524*beaa55a6SMatthew G. Knepley const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)]; 525*beaa55a6SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]);CHKERRQ(ierr);} 526*beaa55a6SMatthew G. Knepley elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]; 527cb1e1211SMatthew G Knepley } 528cb1e1211SMatthew G Knepley } 5299c3cf19fSMatthew G. Knepley fieldOffset += Nb; 530*beaa55a6SMatthew G. Knepley qc += Nc; 531cb1e1211SMatthew G Knepley } 532cb1e1211SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 533cb1e1211SMatthew G Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 534cb1e1211SMatthew G Knepley localDiff += elemDiff; 535cb1e1211SMatthew G Knepley } 5367318780aSToby Isaac ierr = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr); 537cb1e1211SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 538b2566f29SBarry Smith ierr = MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 539cb1e1211SMatthew G Knepley *diff = PetscSqrtReal(*diff); 540cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 541cb1e1211SMatthew G Knepley } 542cb1e1211SMatthew G Knepley 543b698f381SToby Isaac PetscErrorCode DMComputeL2GradientDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff) 544cb1e1211SMatthew G Knepley { 54540e14135SMatthew G. Knepley const PetscInt debug = 0; 546cb1e1211SMatthew G Knepley PetscSection section; 54740e14135SMatthew G. Knepley PetscQuadrature quad; 54840e14135SMatthew G. Knepley Vec localX; 5499c3cf19fSMatthew G. Knepley PetscScalar *funcVal, *interpolant; 5509c3cf19fSMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 5517318780aSToby Isaac PetscReal *coords, *realSpaceDer, *J, *invJ, *detJ; 55240e14135SMatthew G. Knepley PetscReal localDiff = 0.0; 5539c3cf19fSMatthew G. Knepley PetscInt dim, coordDim, qNc = 0, Nq = 0, numFields, numComponents = 0, cStart, cEnd, cEndInterior, c, field, fieldOffset; 554cb1e1211SMatthew G Knepley PetscErrorCode ierr; 555cb1e1211SMatthew G Knepley 556cb1e1211SMatthew G Knepley PetscFunctionBegin; 557c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 5587318780aSToby Isaac ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr); 55940e14135SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 56040e14135SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 56140e14135SMatthew G. Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 56240e14135SMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 56340e14135SMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 564652b88e8SMatthew G. Knepley for (field = 0; field < numFields; ++field) { 5650f2d7e86SMatthew G. Knepley PetscFE fe; 56640e14135SMatthew G. Knepley PetscInt Nc; 567652b88e8SMatthew G. Knepley 5680f2d7e86SMatthew G. Knepley ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 5690f2d7e86SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 5700f2d7e86SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 57140e14135SMatthew G. Knepley numComponents += Nc; 572652b88e8SMatthew G. Knepley } 5739c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 574*beaa55a6SMatthew G. Knepley if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents); 5754d6f44ffSToby Isaac /* ierr = DMProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */ 5769c3cf19fSMatthew G. Knepley ierr = PetscMalloc7(numComponents,&funcVal,coordDim*Nq,&coords,coordDim*Nq,&realSpaceDer,coordDim*coordDim*Nq,&J,coordDim*coordDim*Nq,&invJ,numComponents,&interpolant,Nq,&detJ);CHKERRQ(ierr); 57740e14135SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 5789ac3fadcSMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 5799ac3fadcSMatthew G. Knepley cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 58040e14135SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 58140e14135SMatthew G. Knepley PetscScalar *x = NULL; 58240e14135SMatthew G. Knepley PetscReal elemDiff = 0.0; 5839c3cf19fSMatthew G. Knepley PetscInt qc = 0; 584652b88e8SMatthew G. Knepley 5857318780aSToby Isaac ierr = DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, invJ, detJ);CHKERRQ(ierr); 58640e14135SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 58740e14135SMatthew G. Knepley 5889c3cf19fSMatthew G. Knepley for (field = 0, fieldOffset = 0; field < numFields; ++field) { 5890f2d7e86SMatthew G. Knepley PetscFE fe; 59051259fa3SMatthew G. Knepley void * const ctx = ctxs ? ctxs[field] : NULL; 59140e14135SMatthew G. Knepley PetscReal *basisDer; 5929c3cf19fSMatthew G. Knepley PetscInt Nb, Nc, q, fc; 59340e14135SMatthew G. Knepley 5940f2d7e86SMatthew G. Knepley ierr = DMGetField(dm, field, (PetscObject *) &fe);CHKERRQ(ierr); 5950f2d7e86SMatthew G. Knepley ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 5969c3cf19fSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 5970f2d7e86SMatthew G. Knepley ierr = PetscFEGetDefaultTabulation(fe, NULL, &basisDer, NULL);CHKERRQ(ierr); 59840e14135SMatthew G. Knepley if (debug) { 59940e14135SMatthew G. Knepley char title[1024]; 60040e14135SMatthew G. Knepley ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 6019c3cf19fSMatthew G. Knepley ierr = DMPrintCellVector(c, title, Nb, &x[fieldOffset]);CHKERRQ(ierr); 602652b88e8SMatthew G. Knepley } 6039c3cf19fSMatthew G. Knepley for (q = 0; q < Nq; ++q) { 6047318780aSToby Isaac if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", detJ[q], c, q); 6057318780aSToby Isaac ierr = (*funcs[field])(coordDim, time, &coords[q*coordDim], n, numFields, funcVal, ctx); 606e735a8a9SMatthew G. Knepley if (ierr) { 607e735a8a9SMatthew G. Knepley PetscErrorCode ierr2; 608e735a8a9SMatthew G. Knepley ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2); 609e735a8a9SMatthew G. Knepley ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2); 6109c3cf19fSMatthew G. Knepley ierr2 = PetscFree7(funcVal,coords,realSpaceDer,J,invJ,interpolant,detJ);CHKERRQ(ierr2); 611e735a8a9SMatthew G. Knepley CHKERRQ(ierr); 612e735a8a9SMatthew G. Knepley } 6139c3cf19fSMatthew G. Knepley ierr = PetscFEInterpolateGradient_Static(fe, &x[fieldOffset], coordDim, invJ, n, q, interpolant);CHKERRQ(ierr); 6149c3cf19fSMatthew G. Knepley for (fc = 0; fc < Nc; ++fc) { 615*beaa55a6SMatthew G. Knepley const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)]; 616*beaa55a6SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d fieldDer %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]);CHKERRQ(ierr);} 617*beaa55a6SMatthew G. Knepley elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]; 61840e14135SMatthew G. Knepley } 61940e14135SMatthew G. Knepley } 6209c3cf19fSMatthew G. Knepley fieldOffset += Nb; 6219c3cf19fSMatthew G. Knepley qc += Nc; 62240e14135SMatthew G. Knepley } 62340e14135SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 62440e14135SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 62540e14135SMatthew G. Knepley localDiff += elemDiff; 62640e14135SMatthew G. Knepley } 6279c3cf19fSMatthew G. Knepley ierr = PetscFree7(funcVal,coords,realSpaceDer,J,invJ,interpolant,detJ);CHKERRQ(ierr); 62840e14135SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 629b2566f29SBarry Smith ierr = MPIU_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 63040e14135SMatthew G. Knepley *diff = PetscSqrtReal(*diff); 631cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 632cb1e1211SMatthew G Knepley } 633cb1e1211SMatthew G Knepley 634c6eecec3SToby Isaac PetscErrorCode DMComputeL2FieldDiff_Plex(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff) 63573d901b8SMatthew G. Knepley { 63673d901b8SMatthew G. Knepley const PetscInt debug = 0; 63773d901b8SMatthew G. Knepley PetscSection section; 63873d901b8SMatthew G. Knepley PetscQuadrature quad; 63973d901b8SMatthew G. Knepley Vec localX; 64015496722SMatthew G. Knepley PetscScalar *funcVal, *interpolant; 6417318780aSToby Isaac PetscReal *coords, *detJ, *J; 64273d901b8SMatthew G. Knepley PetscReal *localDiff; 64315496722SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 6449c3cf19fSMatthew G. Knepley PetscInt dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, cEndInterior, c, field, fieldOffset; 64573d901b8SMatthew G. Knepley PetscErrorCode ierr; 64673d901b8SMatthew G. Knepley 64773d901b8SMatthew G. Knepley PetscFunctionBegin; 648c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 6497318780aSToby Isaac ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr); 65073d901b8SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 65173d901b8SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 65273d901b8SMatthew G. Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 653bdd6f66aSToby Isaac ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 65473d901b8SMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 65573d901b8SMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 65673d901b8SMatthew G. Knepley for (field = 0; field < numFields; ++field) { 65715496722SMatthew G. Knepley PetscObject obj; 65815496722SMatthew G. Knepley PetscClassId id; 65973d901b8SMatthew G. Knepley PetscInt Nc; 66073d901b8SMatthew G. Knepley 66115496722SMatthew G. Knepley ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 66215496722SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 66315496722SMatthew G. Knepley if (id == PETSCFE_CLASSID) { 66415496722SMatthew G. Knepley PetscFE fe = (PetscFE) obj; 66515496722SMatthew G. Knepley 6660f2d7e86SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 6670f2d7e86SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 66815496722SMatthew G. Knepley } else if (id == PETSCFV_CLASSID) { 66915496722SMatthew G. Knepley PetscFV fv = (PetscFV) obj; 67015496722SMatthew G. Knepley 67115496722SMatthew G. Knepley ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 67215496722SMatthew G. Knepley ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 67315496722SMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 67473d901b8SMatthew G. Knepley numComponents += Nc; 67573d901b8SMatthew G. Knepley } 6769c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 677*beaa55a6SMatthew G. Knepley if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents); 6787318780aSToby Isaac ierr = PetscCalloc6(numFields,&localDiff,numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J);CHKERRQ(ierr); 67973d901b8SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 6809ac3fadcSMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 6819ac3fadcSMatthew G. Knepley cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 68273d901b8SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 68373d901b8SMatthew G. Knepley PetscScalar *x = NULL; 6849c3cf19fSMatthew G. Knepley PetscInt qc = 0; 68573d901b8SMatthew G. Knepley 6867318780aSToby Isaac ierr = DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, NULL, detJ);CHKERRQ(ierr); 68773d901b8SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 68873d901b8SMatthew G. Knepley 68915496722SMatthew G. Knepley for (field = 0, fieldOffset = 0; field < numFields; ++field) { 69015496722SMatthew G. Knepley PetscObject obj; 69115496722SMatthew G. Knepley PetscClassId id; 69273d901b8SMatthew G. Knepley void * const ctx = ctxs ? ctxs[field] : NULL; 69315496722SMatthew G. Knepley PetscInt Nb, Nc, q, fc; 69473d901b8SMatthew G. Knepley 69515496722SMatthew G. Knepley PetscReal elemDiff = 0.0; 69615496722SMatthew G. Knepley 69715496722SMatthew G. Knepley ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 69815496722SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 69915496722SMatthew G. Knepley if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 70015496722SMatthew G. Knepley else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 70115496722SMatthew G. Knepley else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 70273d901b8SMatthew G. Knepley if (debug) { 70373d901b8SMatthew G. Knepley char title[1024]; 70473d901b8SMatthew G. Knepley ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 70515496722SMatthew G. Knepley ierr = DMPrintCellVector(c, title, Nb*Nc, &x[fieldOffset]);CHKERRQ(ierr); 70673d901b8SMatthew G. Knepley } 7077318780aSToby Isaac for (q = 0; q < Nq; ++q) { 7087318780aSToby Isaac if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature point %D", detJ, c, q); 7097318780aSToby Isaac ierr = (*funcs[field])(coordDim, time, &coords[coordDim*q], numFields, funcVal, ctx); 710e735a8a9SMatthew G. Knepley if (ierr) { 711e735a8a9SMatthew G. Knepley PetscErrorCode ierr2; 712e735a8a9SMatthew G. Knepley ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2); 713e735a8a9SMatthew G. Knepley ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2); 7147318780aSToby Isaac ierr2 = PetscFree6(localDiff,funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr2); 715e735a8a9SMatthew G. Knepley CHKERRQ(ierr); 716e735a8a9SMatthew G. Knepley } 71715496722SMatthew G. Knepley if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 71815496722SMatthew G. Knepley else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 71915496722SMatthew G. Knepley else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 72015496722SMatthew G. Knepley for (fc = 0; fc < Nc; ++fc) { 721*beaa55a6SMatthew G. Knepley const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)]; 722*beaa55a6SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d field %d point %g %g %g diff %g\n", c, field, coordDim > 0 ? coords[0] : 0., coordDim > 1 ? coords[1] : 0., coordDim > 2 ? coords[2] : 0., PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]);CHKERRQ(ierr);} 723*beaa55a6SMatthew G. Knepley elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]; 72473d901b8SMatthew G. Knepley } 72573d901b8SMatthew G. Knepley } 726*beaa55a6SMatthew G. Knepley fieldOffset += Nb; 7279c3cf19fSMatthew G. Knepley qc += Nc; 72873d901b8SMatthew G. Knepley localDiff[field] += elemDiff; 72973d901b8SMatthew G. Knepley } 73073d901b8SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 73173d901b8SMatthew G. Knepley } 73273d901b8SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 733b2566f29SBarry Smith ierr = MPIU_Allreduce(localDiff, diff, numFields, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 73473d901b8SMatthew G. Knepley for (field = 0; field < numFields; ++field) diff[field] = PetscSqrtReal(diff[field]); 7357318780aSToby Isaac ierr = PetscFree6(localDiff,funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr); 73673d901b8SMatthew G. Knepley PetscFunctionReturn(0); 73773d901b8SMatthew G. Knepley } 73873d901b8SMatthew G. Knepley 739e729f68cSMatthew G. Knepley /*@C 740e729f68cSMatthew G. Knepley DMPlexComputeL2DiffVec - This function computes the cellwise L_2 difference between a function u and an FEM interpolant solution u_h, and stores it in a Vec. 741e729f68cSMatthew G. Knepley 742e729f68cSMatthew G. Knepley Input Parameters: 743e729f68cSMatthew G. Knepley + dm - The DM 7440163fd50SMatthew G. Knepley . time - The time 745ca3eba1bSToby Isaac . funcs - The functions to evaluate for each field component: NULL means that component does not contribute to error calculation 746e729f68cSMatthew G. Knepley . ctxs - Optional array of contexts to pass to each function, or NULL. 747e729f68cSMatthew G. Knepley - X - The coefficient vector u_h 748e729f68cSMatthew G. Knepley 749e729f68cSMatthew G. Knepley Output Parameter: 750e729f68cSMatthew G. Knepley . D - A Vec which holds the difference ||u - u_h||_2 for each cell 751e729f68cSMatthew G. Knepley 752e729f68cSMatthew G. Knepley Level: developer 753e729f68cSMatthew G. Knepley 754b698f381SToby Isaac .seealso: DMProjectFunction(), DMComputeL2Diff(), DMPlexComputeL2FieldDiff(), DMComputeL2GradientDiff() 755e729f68cSMatthew G. Knepley @*/ 7560163fd50SMatthew G. Knepley PetscErrorCode DMPlexComputeL2DiffVec(DM dm, PetscReal time, PetscErrorCode (**funcs)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void **ctxs, Vec X, Vec D) 757e729f68cSMatthew G. Knepley { 758e729f68cSMatthew G. Knepley PetscSection section; 759e729f68cSMatthew G. Knepley PetscQuadrature quad; 760e729f68cSMatthew G. Knepley Vec localX; 761e729f68cSMatthew G. Knepley PetscScalar *funcVal, *interpolant; 7627318780aSToby Isaac PetscReal *coords, *detJ, *J; 763e729f68cSMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 7649c3cf19fSMatthew G. Knepley PetscInt dim, coordDim, numFields, numComponents = 0, qNc, Nq, cStart, cEnd, cEndInterior, c, field, fieldOffset; 765e729f68cSMatthew G. Knepley PetscErrorCode ierr; 766e729f68cSMatthew G. Knepley 767e729f68cSMatthew G. Knepley PetscFunctionBegin; 768e729f68cSMatthew G. Knepley ierr = VecSet(D, 0.0);CHKERRQ(ierr); 769e729f68cSMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 7707318780aSToby Isaac ierr = DMGetCoordinateDim(dm, &coordDim);CHKERRQ(ierr); 771e729f68cSMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 772e729f68cSMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 773e729f68cSMatthew G. Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 774bdd6f66aSToby Isaac ierr = DMProjectFunctionLocal(dm, time, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 775e729f68cSMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 776e729f68cSMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 777e729f68cSMatthew G. Knepley for (field = 0; field < numFields; ++field) { 778e729f68cSMatthew G. Knepley PetscObject obj; 779e729f68cSMatthew G. Knepley PetscClassId id; 780e729f68cSMatthew G. Knepley PetscInt Nc; 781e729f68cSMatthew G. Knepley 782e729f68cSMatthew G. Knepley ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 783e729f68cSMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 784e729f68cSMatthew G. Knepley if (id == PETSCFE_CLASSID) { 785e729f68cSMatthew G. Knepley PetscFE fe = (PetscFE) obj; 786e729f68cSMatthew G. Knepley 787e729f68cSMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &quad);CHKERRQ(ierr); 788e729f68cSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 789e729f68cSMatthew G. Knepley } else if (id == PETSCFV_CLASSID) { 790e729f68cSMatthew G. Knepley PetscFV fv = (PetscFV) obj; 791e729f68cSMatthew G. Knepley 792e729f68cSMatthew G. Knepley ierr = PetscFVGetQuadrature(fv, &quad);CHKERRQ(ierr); 793e729f68cSMatthew G. Knepley ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 794e729f68cSMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 795e729f68cSMatthew G. Knepley numComponents += Nc; 796e729f68cSMatthew G. Knepley } 7979c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &qNc, &Nq, &quadPoints, &quadWeights);CHKERRQ(ierr); 798*beaa55a6SMatthew G. Knepley if ((qNc != 1) && (qNc != numComponents)) SETERRQ2(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_SIZ, "Quadrature components %D != %D field components", qNc, numComponents); 7997318780aSToby Isaac ierr = PetscMalloc5(numComponents,&funcVal,numComponents,&interpolant,coordDim*Nq,&coords,Nq,&detJ,coordDim*coordDim*Nq,&J);CHKERRQ(ierr); 800e729f68cSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 801e729f68cSMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 802e729f68cSMatthew G. Knepley cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 803e729f68cSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 804e729f68cSMatthew G. Knepley PetscScalar *x = NULL; 8056f288a59SMatthew G. Knepley PetscScalar elemDiff = 0.0; 8069c3cf19fSMatthew G. Knepley PetscInt qc = 0; 807e729f68cSMatthew G. Knepley 8087318780aSToby Isaac ierr = DMPlexComputeCellGeometryFEM(dm, c, quad, coords, J, NULL, detJ);CHKERRQ(ierr); 809e729f68cSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 810e729f68cSMatthew G. Knepley 811e729f68cSMatthew G. Knepley for (field = 0, fieldOffset = 0; field < numFields; ++field) { 812e729f68cSMatthew G. Knepley PetscObject obj; 813e729f68cSMatthew G. Knepley PetscClassId id; 814e729f68cSMatthew G. Knepley void * const ctx = ctxs ? ctxs[field] : NULL; 815e729f68cSMatthew G. Knepley PetscInt Nb, Nc, q, fc; 816e729f68cSMatthew G. Knepley 817e729f68cSMatthew G. Knepley ierr = DMGetField(dm, field, &obj);CHKERRQ(ierr); 818e729f68cSMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 819e729f68cSMatthew G. Knepley if (id == PETSCFE_CLASSID) {ierr = PetscFEGetNumComponents((PetscFE) obj, &Nc);CHKERRQ(ierr);ierr = PetscFEGetDimension((PetscFE) obj, &Nb);CHKERRQ(ierr);} 820e729f68cSMatthew G. Knepley else if (id == PETSCFV_CLASSID) {ierr = PetscFVGetNumComponents((PetscFV) obj, &Nc);CHKERRQ(ierr);Nb = 1;} 821e729f68cSMatthew G. Knepley else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 82223f34ed2SToby Isaac if (funcs[field]) { 8237318780aSToby Isaac for (q = 0; q < Nq; ++q) { 8247318780aSToby Isaac if (detJ[q] <= 0.0) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %D, quadrature points %D", detJ[q], c, q); 8257318780aSToby Isaac ierr = (*funcs[field])(dim, time, &coords[q*coordDim], Nc, funcVal, ctx); 826e735a8a9SMatthew G. Knepley if (ierr) { 827e735a8a9SMatthew G. Knepley PetscErrorCode ierr2; 828e735a8a9SMatthew G. Knepley ierr2 = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr2); 8297318780aSToby Isaac ierr2 = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr); 830e735a8a9SMatthew G. Knepley ierr2 = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr2); 831e735a8a9SMatthew G. Knepley CHKERRQ(ierr); 832e735a8a9SMatthew G. Knepley } 833e729f68cSMatthew G. Knepley if (id == PETSCFE_CLASSID) {ierr = PetscFEInterpolate_Static((PetscFE) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 834e729f68cSMatthew G. Knepley else if (id == PETSCFV_CLASSID) {ierr = PetscFVInterpolate_Static((PetscFV) obj, &x[fieldOffset], q, interpolant);CHKERRQ(ierr);} 835e729f68cSMatthew G. Knepley else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", field); 836e729f68cSMatthew G. Knepley for (fc = 0; fc < Nc; ++fc) { 837*beaa55a6SMatthew G. Knepley const PetscReal wt = quadWeights[q*qNc+(qNc == 1 ? 0 : qc+fc)]; 838*beaa55a6SMatthew G. Knepley elemDiff += PetscSqr(PetscRealPart(interpolant[fc] - funcVal[fc]))*wt*detJ[q]; 839e729f68cSMatthew G. Knepley } 840e729f68cSMatthew G. Knepley } 84123f34ed2SToby Isaac } 842*beaa55a6SMatthew G. Knepley fieldOffset += Nb; 8439c3cf19fSMatthew G. Knepley qc += Nc; 844e729f68cSMatthew G. Knepley } 845e729f68cSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 84623f34ed2SToby Isaac ierr = VecSetValue(D, c - cStart, elemDiff, INSERT_VALUES);CHKERRQ(ierr); 847e729f68cSMatthew G. Knepley } 8487318780aSToby Isaac ierr = PetscFree5(funcVal,interpolant,coords,detJ,J);CHKERRQ(ierr); 849e729f68cSMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 850e729f68cSMatthew G. Knepley ierr = VecSqrtAbs(D);CHKERRQ(ierr); 851e729f68cSMatthew G. Knepley PetscFunctionReturn(0); 852e729f68cSMatthew G. Knepley } 853e729f68cSMatthew G. Knepley 85473d901b8SMatthew G. Knepley /*@ 85573d901b8SMatthew G. Knepley DMPlexComputeIntegralFEM - Form the local integral F from the local input X using pointwise functions specified by the user 85673d901b8SMatthew G. Knepley 85773d901b8SMatthew G. Knepley Input Parameters: 85873d901b8SMatthew G. Knepley + dm - The mesh 85973d901b8SMatthew G. Knepley . X - Local input vector 86073d901b8SMatthew G. Knepley - user - The user context 86173d901b8SMatthew G. Knepley 86273d901b8SMatthew G. Knepley Output Parameter: 86373d901b8SMatthew G. Knepley . integral - Local integral for each field 86473d901b8SMatthew G. Knepley 86573d901b8SMatthew G. Knepley Level: developer 86673d901b8SMatthew G. Knepley 86773d901b8SMatthew G. Knepley .seealso: DMPlexComputeResidualFEM() 86873d901b8SMatthew G. Knepley @*/ 8690f2d7e86SMatthew G. Knepley PetscErrorCode DMPlexComputeIntegralFEM(DM dm, Vec X, PetscReal *integral, void *user) 87073d901b8SMatthew G. Knepley { 87173d901b8SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 872b5a3613cSMatthew G. Knepley DM dmAux, dmGrad; 873425f1808SMatthew G. Knepley Vec localX, A, cellGeometryFVM = NULL, faceGeometryFVM = NULL, locGrad = NULL; 8742764a2aaSMatthew G. Knepley PetscDS prob, probAux = NULL; 87573d901b8SMatthew G. Knepley PetscSection section, sectionAux; 876b5a3613cSMatthew G. Knepley PetscFV fvm = NULL; 877b5a3613cSMatthew G. Knepley PetscFECellGeom *cgeomFEM; 878b5a3613cSMatthew G. Knepley PetscFVFaceGeom *fgeomFVM; 879b5a3613cSMatthew G. Knepley PetscFVCellGeom *cgeomFVM; 88073d901b8SMatthew G. Knepley PetscScalar *u, *a = NULL; 881b5a3613cSMatthew G. Knepley const PetscScalar *lgrad; 882b5a3613cSMatthew G. Knepley PetscReal *lintegral; 883b5a3613cSMatthew G. Knepley PetscInt *uOff, *uOff_x, *aOff = NULL; 884b5a3613cSMatthew G. Knepley PetscInt dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, cEndInterior, c; 8850f2d7e86SMatthew G. Knepley PetscInt totDim, totDimAux; 886425f1808SMatthew G. Knepley PetscBool useFVM = PETSC_FALSE; 88773d901b8SMatthew G. Knepley PetscErrorCode ierr; 88873d901b8SMatthew G. Knepley 88973d901b8SMatthew G. Knepley PetscFunctionBegin; 890c1f031eeSMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr); 89173d901b8SMatthew G. Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 892bdd6f66aSToby Isaac ierr = DMPlexInsertBoundaryValues(dm, PETSC_TRUE, localX, 0.0, NULL, NULL, NULL);CHKERRQ(ierr); 89373d901b8SMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 89473d901b8SMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 895c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 89673d901b8SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 8972764a2aaSMatthew G. Knepley ierr = DMGetDS(dm, &prob);CHKERRQ(ierr); 8982764a2aaSMatthew G. Knepley ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 899b5a3613cSMatthew G. Knepley ierr = PetscDSGetComponentOffsets(prob, &uOff);CHKERRQ(ierr); 900b5a3613cSMatthew G. Knepley ierr = PetscDSGetComponentDerivativeOffsets(prob, &uOff_x);CHKERRQ(ierr); 90173d901b8SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 90273d901b8SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 9039ac3fadcSMatthew G. Knepley ierr = DMPlexGetHybridBounds(dm, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 9049ac3fadcSMatthew G. Knepley cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 90573d901b8SMatthew G. Knepley numCells = cEnd - cStart; 90673d901b8SMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 90773d901b8SMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 90873d901b8SMatthew G. Knepley if (dmAux) { 9092764a2aaSMatthew G. Knepley ierr = DMGetDS(dmAux, &probAux);CHKERRQ(ierr); 910b5a3613cSMatthew G. Knepley ierr = PetscDSGetNumFields(probAux, &NfAux);CHKERRQ(ierr); 911b5a3613cSMatthew G. Knepley ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 9122764a2aaSMatthew G. Knepley ierr = PetscDSGetTotalDimension(probAux, &totDimAux);CHKERRQ(ierr); 913b5a3613cSMatthew G. Knepley ierr = PetscDSGetComponentOffsets(probAux, &aOff);CHKERRQ(ierr); 914b5a3613cSMatthew G. Knepley ierr = PetscDSGetEvaluationArrays(probAux, &a, NULL, NULL);CHKERRQ(ierr); 91573d901b8SMatthew G. Knepley } 916b5a3613cSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 917b5a3613cSMatthew G. Knepley PetscObject obj; 918b5a3613cSMatthew G. Knepley PetscClassId id; 919b5a3613cSMatthew G. Knepley 920b5a3613cSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 921b5a3613cSMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 922b5a3613cSMatthew G. Knepley if (id == PETSCFV_CLASSID) {useFVM = PETSC_TRUE; fvm = (PetscFV) obj;} 923b5a3613cSMatthew G. Knepley } 924b5a3613cSMatthew G. Knepley if (useFVM) { 925b5a3613cSMatthew G. Knepley Vec grad; 926b5a3613cSMatthew G. Knepley PetscInt fStart, fEnd; 927b5a3613cSMatthew G. Knepley PetscBool compGrad; 928b5a3613cSMatthew G. Knepley 929b5a3613cSMatthew G. Knepley ierr = PetscFVGetComputeGradients(fvm, &compGrad);CHKERRQ(ierr); 930b5a3613cSMatthew G. Knepley ierr = PetscFVSetComputeGradients(fvm, PETSC_TRUE);CHKERRQ(ierr); 931b5a3613cSMatthew G. Knepley ierr = DMPlexComputeGeometryFVM(dm, &cellGeometryFVM, &faceGeometryFVM);CHKERRQ(ierr); 932b5a3613cSMatthew G. Knepley ierr = DMPlexComputeGradientFVM(dm, fvm, faceGeometryFVM, cellGeometryFVM, &dmGrad);CHKERRQ(ierr); 933b5a3613cSMatthew G. Knepley ierr = PetscFVSetComputeGradients(fvm, compGrad);CHKERRQ(ierr); 934b5a3613cSMatthew G. Knepley ierr = VecGetArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);CHKERRQ(ierr); 935b5a3613cSMatthew G. Knepley ierr = VecGetArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr); 936b5a3613cSMatthew G. Knepley /* Reconstruct and limit cell gradients */ 937b5a3613cSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 1, &fStart, &fEnd);CHKERRQ(ierr); 938b5a3613cSMatthew G. Knepley ierr = DMGetGlobalVector(dmGrad, &grad);CHKERRQ(ierr); 939b5a3613cSMatthew G. Knepley ierr = DMPlexReconstructGradients_Internal(dm, fvm, fStart, fEnd, faceGeometryFVM, cellGeometryFVM, localX, grad);CHKERRQ(ierr); 940b5a3613cSMatthew G. Knepley /* Communicate gradient values */ 941b5a3613cSMatthew G. Knepley ierr = DMGetLocalVector(dmGrad, &locGrad);CHKERRQ(ierr); 942b5a3613cSMatthew G. Knepley ierr = DMGlobalToLocalBegin(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 943b5a3613cSMatthew G. Knepley ierr = DMGlobalToLocalEnd(dmGrad, grad, INSERT_VALUES, locGrad);CHKERRQ(ierr); 944b5a3613cSMatthew G. Knepley ierr = DMRestoreGlobalVector(dmGrad, &grad);CHKERRQ(ierr); 945b5a3613cSMatthew G. Knepley /* Handle non-essential (e.g. outflow) boundary values */ 946b5a3613cSMatthew G. Knepley ierr = DMPlexInsertBoundaryValues(dm, PETSC_FALSE, localX, 0.0, faceGeometryFVM, cellGeometryFVM, locGrad);CHKERRQ(ierr); 947b5a3613cSMatthew G. Knepley ierr = VecGetArrayRead(locGrad, &lgrad);CHKERRQ(ierr); 948b5a3613cSMatthew G. Knepley } 949b5a3613cSMatthew G. Knepley ierr = PetscMalloc3(Nf,&lintegral,numCells*totDim,&u,numCells,&cgeomFEM);CHKERRQ(ierr); 9500f2d7e86SMatthew G. Knepley if (dmAux) {ierr = PetscMalloc1(numCells*totDimAux, &a);CHKERRQ(ierr);} 951c1f031eeSMatthew G. Knepley for (f = 0; f < Nf; ++f) {lintegral[f] = 0.0;} 95273d901b8SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 95373d901b8SMatthew G. Knepley PetscScalar *x = NULL; 95473d901b8SMatthew G. Knepley PetscInt i; 95573d901b8SMatthew G. Knepley 956b5a3613cSMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, cgeomFEM[c].v0, cgeomFEM[c].J, cgeomFEM[c].invJ, &cgeomFEM[c].detJ);CHKERRQ(ierr); 957b5a3613cSMatthew G. Knepley if (cgeomFEM[c].detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", cgeomFEM[c].detJ, c); 95873d901b8SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr); 9590f2d7e86SMatthew G. Knepley for (i = 0; i < totDim; ++i) u[c*totDim+i] = x[i]; 96073d901b8SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, section, localX, c, NULL, &x);CHKERRQ(ierr); 96173d901b8SMatthew G. Knepley if (dmAux) { 96273d901b8SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 9630f2d7e86SMatthew G. Knepley for (i = 0; i < totDimAux; ++i) a[c*totDimAux+i] = x[i]; 96473d901b8SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 96573d901b8SMatthew G. Knepley } 96673d901b8SMatthew G. Knepley } 96773d901b8SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 968c1f031eeSMatthew G. Knepley PetscObject obj; 969c1f031eeSMatthew G. Knepley PetscClassId id; 970c1f031eeSMatthew G. Knepley PetscInt numChunks, numBatches, batchSize, numBlocks, blockSize, Ne, Nr, offset; 97173d901b8SMatthew G. Knepley 972c1f031eeSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 973c1f031eeSMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 974c1f031eeSMatthew G. Knepley if (id == PETSCFE_CLASSID) { 975c1f031eeSMatthew G. Knepley PetscFE fe = (PetscFE) obj; 976c1f031eeSMatthew G. Knepley PetscQuadrature q; 977c1f031eeSMatthew G. Knepley PetscInt Nq, Nb; 978c1f031eeSMatthew G. Knepley 9790f2d7e86SMatthew G. Knepley ierr = PetscFEGetTileSizes(fe, NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 980c1f031eeSMatthew G. Knepley ierr = PetscFEGetQuadrature(fe, &q);CHKERRQ(ierr); 9819c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(q, NULL, NULL, &Nq, NULL, NULL);CHKERRQ(ierr); 982c1f031eeSMatthew G. Knepley ierr = PetscFEGetDimension(fe, &Nb);CHKERRQ(ierr); 983c1f031eeSMatthew G. Knepley blockSize = Nb*Nq; 98473d901b8SMatthew G. Knepley batchSize = numBlocks * blockSize; 9850f2d7e86SMatthew G. Knepley ierr = PetscFESetTileSizes(fe, blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 98673d901b8SMatthew G. Knepley numChunks = numCells / (numBatches*batchSize); 98773d901b8SMatthew G. Knepley Ne = numChunks*numBatches*batchSize; 98873d901b8SMatthew G. Knepley Nr = numCells % (numBatches*batchSize); 98973d901b8SMatthew G. Knepley offset = numCells - Nr; 990b5a3613cSMatthew G. Knepley ierr = PetscFEIntegrate(fe, prob, f, Ne, cgeomFEM, u, probAux, a, lintegral);CHKERRQ(ierr); 991b5a3613cSMatthew G. Knepley ierr = PetscFEIntegrate(fe, prob, f, Nr, &cgeomFEM[offset], &u[offset*totDim], probAux, &a[offset*totDimAux], lintegral);CHKERRQ(ierr); 992c1f031eeSMatthew G. Knepley } else if (id == PETSCFV_CLASSID) { 993b69edc29SMatthew G. Knepley /* PetscFV fv = (PetscFV) obj; */ 994c1f031eeSMatthew G. Knepley PetscInt foff; 995420e96edSMatthew G. Knepley PetscPointFunc obj_func; 996b69edc29SMatthew G. Knepley PetscScalar lint; 997c1f031eeSMatthew G. Knepley 998c1f031eeSMatthew G. Knepley ierr = PetscDSGetObjective(prob, f, &obj_func);CHKERRQ(ierr); 999c1f031eeSMatthew G. Knepley ierr = PetscDSGetFieldOffset(prob, f, &foff);CHKERRQ(ierr); 1000c1f031eeSMatthew G. Knepley if (obj_func) { 1001c1f031eeSMatthew G. Knepley for (c = 0; c < numCells; ++c) { 1002b5a3613cSMatthew G. Knepley PetscScalar *u_x; 1003b5a3613cSMatthew G. Knepley 1004b5a3613cSMatthew G. Knepley ierr = DMPlexPointLocalRead(dmGrad, c, lgrad, &u_x);CHKERRQ(ierr); 1005b5a3613cSMatthew G. Knepley obj_func(dim, Nf, NfAux, uOff, uOff_x, &u[totDim*c+foff], NULL, u_x, aOff, NULL, a, NULL, NULL, 0.0, cgeomFVM[c].centroid, &lint); 1006b5a3613cSMatthew G. Knepley lintegral[f] += PetscRealPart(lint)*cgeomFVM[c].volume; 100773d901b8SMatthew G. Knepley } 1008c1f031eeSMatthew G. Knepley } 1009c1f031eeSMatthew G. Knepley } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Unknown discretization type for field %d", f); 1010c1f031eeSMatthew G. Knepley } 1011b5a3613cSMatthew G. Knepley if (useFVM) { 1012b5a3613cSMatthew G. Knepley ierr = VecRestoreArrayRead(locGrad, &lgrad);CHKERRQ(ierr); 1013b5a3613cSMatthew G. Knepley ierr = VecRestoreArrayRead(faceGeometryFVM, (const PetscScalar **) &fgeomFVM);CHKERRQ(ierr); 1014b5a3613cSMatthew G. Knepley ierr = VecRestoreArrayRead(cellGeometryFVM, (const PetscScalar **) &cgeomFVM);CHKERRQ(ierr); 1015b5a3613cSMatthew G. Knepley ierr = DMRestoreLocalVector(dmGrad, &locGrad);CHKERRQ(ierr); 1016b5a3613cSMatthew G. Knepley ierr = VecDestroy(&faceGeometryFVM);CHKERRQ(ierr); 1017b5a3613cSMatthew G. Knepley ierr = VecDestroy(&cellGeometryFVM);CHKERRQ(ierr); 1018b5a3613cSMatthew G. Knepley ierr = DMDestroy(&dmGrad);CHKERRQ(ierr); 1019b5a3613cSMatthew G. Knepley } 102073d901b8SMatthew G. Knepley if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 102173d901b8SMatthew G. Knepley if (mesh->printFEM) { 102273d901b8SMatthew G. Knepley ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "Local integral:");CHKERRQ(ierr); 1023c1f031eeSMatthew G. Knepley for (f = 0; f < Nf; ++f) {ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), " %g", lintegral[f]);CHKERRQ(ierr);} 102473d901b8SMatthew G. Knepley ierr = PetscPrintf(PetscObjectComm((PetscObject) dm), "\n");CHKERRQ(ierr); 102573d901b8SMatthew G. Knepley } 102673d901b8SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 1027b2566f29SBarry Smith ierr = MPIU_Allreduce(lintegral, integral, Nf, MPIU_REAL, MPIU_SUM, PetscObjectComm((PetscObject) dm));CHKERRQ(ierr); 1028b5a3613cSMatthew G. Knepley ierr = PetscFree3(lintegral,u,cgeomFEM);CHKERRQ(ierr); 1029c1f031eeSMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_IntegralFEM,dm,0,0,0);CHKERRQ(ierr); 103073d901b8SMatthew G. Knepley PetscFunctionReturn(0); 103173d901b8SMatthew G. Knepley } 103273d901b8SMatthew G. Knepley 1033d69c5d34SMatthew G. Knepley /*@ 103468132eb9SMatthew G. Knepley DMPlexComputeInterpolatorNested - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM. 1035d69c5d34SMatthew G. Knepley 1036d69c5d34SMatthew G. Knepley Input Parameters: 1037d69c5d34SMatthew G. Knepley + dmf - The fine mesh 1038d69c5d34SMatthew G. Knepley . dmc - The coarse mesh 1039d69c5d34SMatthew G. Knepley - user - The user context 1040d69c5d34SMatthew G. Knepley 1041d69c5d34SMatthew G. Knepley Output Parameter: 1042934789fcSMatthew G. Knepley . In - The interpolation matrix 1043d69c5d34SMatthew G. Knepley 1044d69c5d34SMatthew G. Knepley Level: developer 1045d69c5d34SMatthew G. Knepley 104668132eb9SMatthew G. Knepley .seealso: DMPlexComputeInterpolatorGeneral(), DMPlexComputeJacobianFEM() 1047d69c5d34SMatthew G. Knepley @*/ 104868132eb9SMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorNested(DM dmc, DM dmf, Mat In, void *user) 1049d69c5d34SMatthew G. Knepley { 1050d69c5d34SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dmc->data; 1051d69c5d34SMatthew G. Knepley const char *name = "Interpolator"; 10522764a2aaSMatthew G. Knepley PetscDS prob; 1053d69c5d34SMatthew G. Knepley PetscFE *feRef; 105497c42addSMatthew G. Knepley PetscFV *fvRef; 1055d69c5d34SMatthew G. Knepley PetscSection fsection, fglobalSection; 1056d69c5d34SMatthew G. Knepley PetscSection csection, cglobalSection; 1057d69c5d34SMatthew G. Knepley PetscScalar *elemMat; 10589ac3fadcSMatthew G. Knepley PetscInt dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, cEndInterior, c; 10590f2d7e86SMatthew G. Knepley PetscInt cTotDim, rTotDim = 0; 1060d69c5d34SMatthew G. Knepley PetscErrorCode ierr; 1061d69c5d34SMatthew G. Knepley 1062d69c5d34SMatthew G. Knepley PetscFunctionBegin; 1063d69c5d34SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1064c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr); 1065d69c5d34SMatthew G. Knepley ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 1066d69c5d34SMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 1067d69c5d34SMatthew G. Knepley ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 1068d69c5d34SMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 1069d69c5d34SMatthew G. Knepley ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 1070d69c5d34SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 10719ac3fadcSMatthew G. Knepley ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 10729ac3fadcSMatthew G. Knepley cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 10732764a2aaSMatthew G. Knepley ierr = DMGetDS(dmf, &prob);CHKERRQ(ierr); 107497c42addSMatthew G. Knepley ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr); 1075d69c5d34SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 107697c42addSMatthew G. Knepley PetscObject obj; 107797c42addSMatthew G. Knepley PetscClassId id; 1078aa7890ccSMatthew G. Knepley PetscInt rNb = 0, Nc = 0; 1079d69c5d34SMatthew G. Knepley 108097c42addSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 108197c42addSMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 108297c42addSMatthew G. Knepley if (id == PETSCFE_CLASSID) { 108397c42addSMatthew G. Knepley PetscFE fe = (PetscFE) obj; 108497c42addSMatthew G. Knepley 10850f2d7e86SMatthew G. Knepley ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr); 1086d69c5d34SMatthew G. Knepley ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr); 10870f2d7e86SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 108897c42addSMatthew G. Knepley } else if (id == PETSCFV_CLASSID) { 108997c42addSMatthew G. Knepley PetscFV fv = (PetscFV) obj; 109097c42addSMatthew G. Knepley PetscDualSpace Q; 109197c42addSMatthew G. Knepley 109297c42addSMatthew G. Knepley ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr); 109397c42addSMatthew G. Knepley ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr); 109497c42addSMatthew G. Knepley ierr = PetscDualSpaceGetDimension(Q, &rNb);CHKERRQ(ierr); 109597c42addSMatthew G. Knepley ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 109697c42addSMatthew G. Knepley } 10979c3cf19fSMatthew G. Knepley rTotDim += rNb; 1098d69c5d34SMatthew G. Knepley } 10992764a2aaSMatthew G. Knepley ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr); 11000f2d7e86SMatthew G. Knepley ierr = PetscMalloc1(rTotDim*cTotDim,&elemMat);CHKERRQ(ierr); 11010f2d7e86SMatthew G. Knepley ierr = PetscMemzero(elemMat, rTotDim*cTotDim * sizeof(PetscScalar));CHKERRQ(ierr); 1102d69c5d34SMatthew G. Knepley for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) { 1103d69c5d34SMatthew G. Knepley PetscDualSpace Qref; 1104d69c5d34SMatthew G. Knepley PetscQuadrature f; 1105d69c5d34SMatthew G. Knepley const PetscReal *qpoints, *qweights; 1106d69c5d34SMatthew G. Knepley PetscReal *points; 1107d69c5d34SMatthew G. Knepley PetscInt npoints = 0, Nc, Np, fpdim, i, k, p, d; 1108d69c5d34SMatthew G. Knepley 1109d69c5d34SMatthew G. Knepley /* Compose points from all dual basis functionals */ 111097c42addSMatthew G. Knepley if (feRef[fieldI]) { 1111d69c5d34SMatthew G. Knepley ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr); 11120f2d7e86SMatthew G. Knepley ierr = PetscFEGetNumComponents(feRef[fieldI], &Nc);CHKERRQ(ierr); 111397c42addSMatthew G. Knepley } else { 111497c42addSMatthew G. Knepley ierr = PetscFVGetDualSpace(fvRef[fieldI], &Qref);CHKERRQ(ierr); 111597c42addSMatthew G. Knepley ierr = PetscFVGetNumComponents(fvRef[fieldI], &Nc);CHKERRQ(ierr); 111697c42addSMatthew G. Knepley } 1117d69c5d34SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr); 1118d69c5d34SMatthew G. Knepley for (i = 0; i < fpdim; ++i) { 1119d69c5d34SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 11209c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(f, NULL, NULL, &Np, NULL, NULL);CHKERRQ(ierr); 1121d69c5d34SMatthew G. Knepley npoints += Np; 1122d69c5d34SMatthew G. Knepley } 1123d69c5d34SMatthew G. Knepley ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr); 1124d69c5d34SMatthew G. Knepley for (i = 0, k = 0; i < fpdim; ++i) { 1125d69c5d34SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 11269c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(f, NULL, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr); 1127d69c5d34SMatthew G. Knepley for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d]; 1128d69c5d34SMatthew G. Knepley } 1129d69c5d34SMatthew G. Knepley 1130d69c5d34SMatthew G. Knepley for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) { 113197c42addSMatthew G. Knepley PetscObject obj; 113297c42addSMatthew G. Knepley PetscClassId id; 1133d69c5d34SMatthew G. Knepley PetscReal *B; 11349c3cf19fSMatthew G. Knepley PetscInt NcJ = 0, cpdim = 0, j, qNc; 1135d69c5d34SMatthew G. Knepley 113697c42addSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, fieldJ, &obj);CHKERRQ(ierr); 113797c42addSMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 113897c42addSMatthew G. Knepley if (id == PETSCFE_CLASSID) { 113997c42addSMatthew G. Knepley PetscFE fe = (PetscFE) obj; 1140d69c5d34SMatthew G. Knepley 1141d69c5d34SMatthew G. Knepley /* Evaluate basis at points */ 11420f2d7e86SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &NcJ);CHKERRQ(ierr); 11430f2d7e86SMatthew G. Knepley ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr); 1144ffe73a53SMatthew G. Knepley /* For now, fields only interpolate themselves */ 1145ffe73a53SMatthew G. Knepley if (fieldI == fieldJ) { 11469c3cf19fSMatthew 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); 11470f2d7e86SMatthew G. Knepley ierr = PetscFEGetTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr); 1148d69c5d34SMatthew G. Knepley for (i = 0, k = 0; i < fpdim; ++i) { 1149d69c5d34SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 11509c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(f, NULL, &qNc, &Np, NULL, &qweights);CHKERRQ(ierr); 11519c3cf19fSMatthew G. Knepley if (qNc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %D does not match coarse field %D", qNc, NcJ); 1152d69c5d34SMatthew G. Knepley for (p = 0; p < Np; ++p, ++k) { 115336a6d9c0SMatthew G. Knepley for (j = 0; j < cpdim; ++j) { 11549c3cf19fSMatthew G. Knepley /* NOTE: This is not quite right, unless fpdim == number of fine grid functional quad points */ 11559c3cf19fSMatthew G. Knepley for (c = 0; c < Nc; ++c) elemMat[(offsetI + i)*cTotDim + offsetJ + j] += B[k*cpdim*NcJ+j*Nc+c]*qweights[p*qNc+c]; 115636a6d9c0SMatthew G. Knepley } 1157d69c5d34SMatthew G. Knepley } 1158d69c5d34SMatthew G. Knepley } 11590f2d7e86SMatthew G. Knepley ierr = PetscFERestoreTabulation(fe, npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 1160ffe73a53SMatthew G. Knepley } 116197c42addSMatthew G. Knepley } else if (id == PETSCFV_CLASSID) { 116297c42addSMatthew G. Knepley PetscFV fv = (PetscFV) obj; 116397c42addSMatthew G. Knepley 116497c42addSMatthew G. Knepley /* Evaluate constant function at points */ 116597c42addSMatthew G. Knepley ierr = PetscFVGetNumComponents(fv, &NcJ);CHKERRQ(ierr); 116697c42addSMatthew G. Knepley cpdim = 1; 116797c42addSMatthew G. Knepley /* For now, fields only interpolate themselves */ 116897c42addSMatthew G. Knepley if (fieldI == fieldJ) { 116997c42addSMatthew 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); 117097c42addSMatthew G. Knepley for (i = 0, k = 0; i < fpdim; ++i) { 117197c42addSMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 11729c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(f, NULL, &qNc, &Np, NULL, &qweights);CHKERRQ(ierr); 11739c3cf19fSMatthew G. Knepley if (qNc != NcJ) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %D does not match coarse field %D", qNc, NcJ); 117497c42addSMatthew G. Knepley for (p = 0; p < Np; ++p, ++k) { 117597c42addSMatthew G. Knepley for (j = 0; j < cpdim; ++j) { 11769c3cf19fSMatthew G. Knepley for (c = 0; c < Nc; ++c) elemMat[(offsetI + i*Nc + c)*cTotDim + offsetJ + j*NcJ + c] += 1.0*qweights[p*qNc+c]; 117797c42addSMatthew G. Knepley } 117897c42addSMatthew G. Knepley } 117997c42addSMatthew G. Knepley } 118097c42addSMatthew G. Knepley } 118197c42addSMatthew G. Knepley } 118236a6d9c0SMatthew G. Knepley offsetJ += cpdim*NcJ; 1183d69c5d34SMatthew G. Knepley } 1184d69c5d34SMatthew G. Knepley offsetI += fpdim*Nc; 1185549a8adaSMatthew G. Knepley ierr = PetscFree(points);CHKERRQ(ierr); 1186d69c5d34SMatthew G. Knepley } 11870f2d7e86SMatthew G. Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rTotDim, cTotDim, elemMat);CHKERRQ(ierr);} 11887f5b169aSMatthew G. Knepley /* Preallocate matrix */ 11897f5b169aSMatthew G. Knepley { 1190c094ef40SMatthew G. Knepley Mat preallocator; 1191c094ef40SMatthew G. Knepley PetscScalar *vals; 1192c094ef40SMatthew G. Knepley PetscInt *cellCIndices, *cellFIndices; 1193c094ef40SMatthew G. Knepley PetscInt locRows, locCols, cell; 11947f5b169aSMatthew G. Knepley 1195c094ef40SMatthew G. Knepley ierr = MatGetLocalSize(In, &locRows, &locCols);CHKERRQ(ierr); 1196c094ef40SMatthew G. Knepley ierr = MatCreate(PetscObjectComm((PetscObject) In), &preallocator);CHKERRQ(ierr); 1197c094ef40SMatthew G. Knepley ierr = MatSetType(preallocator, MATPREALLOCATOR);CHKERRQ(ierr); 1198c094ef40SMatthew G. Knepley ierr = MatSetSizes(preallocator, locRows, locCols, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 1199c094ef40SMatthew G. Knepley ierr = MatSetUp(preallocator);CHKERRQ(ierr); 1200c094ef40SMatthew G. Knepley ierr = PetscCalloc3(rTotDim*cTotDim, &vals,cTotDim,&cellCIndices,rTotDim,&cellFIndices);CHKERRQ(ierr); 12017f5b169aSMatthew G. Knepley for (cell = cStart; cell < cEnd; ++cell) { 12027f5b169aSMatthew G. Knepley ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, cell, cellCIndices, cellFIndices);CHKERRQ(ierr); 1203c094ef40SMatthew G. Knepley ierr = MatSetValues(preallocator, rTotDim, cellFIndices, cTotDim, cellCIndices, vals, INSERT_VALUES);CHKERRQ(ierr); 12047f5b169aSMatthew G. Knepley } 1205c094ef40SMatthew G. Knepley ierr = PetscFree3(vals,cellCIndices,cellFIndices);CHKERRQ(ierr); 1206c094ef40SMatthew G. Knepley ierr = MatAssemblyBegin(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1207c094ef40SMatthew G. Knepley ierr = MatAssemblyEnd(preallocator, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1208c094ef40SMatthew G. Knepley ierr = MatPreallocatorPreallocate(preallocator, PETSC_TRUE, In);CHKERRQ(ierr); 1209c094ef40SMatthew G. Knepley ierr = MatDestroy(&preallocator);CHKERRQ(ierr); 12107f5b169aSMatthew G. Knepley } 12117f5b169aSMatthew G. Knepley /* Fill matrix */ 12127f5b169aSMatthew G. Knepley ierr = MatZeroEntries(In);CHKERRQ(ierr); 1213d69c5d34SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 1214934789fcSMatthew G. Knepley ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr); 1215d69c5d34SMatthew G. Knepley } 1216549a8adaSMatthew G. Knepley for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);} 121797c42addSMatthew G. Knepley ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr); 1218549a8adaSMatthew G. Knepley ierr = PetscFree(elemMat);CHKERRQ(ierr); 1219934789fcSMatthew G. Knepley ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1220934789fcSMatthew G. Knepley ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1221d69c5d34SMatthew G. Knepley if (mesh->printFEM) { 1222d69c5d34SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 1223934789fcSMatthew G. Knepley ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr); 1224934789fcSMatthew G. Knepley ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1225d69c5d34SMatthew G. Knepley } 1226d69c5d34SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1227d69c5d34SMatthew G. Knepley PetscFunctionReturn(0); 1228d69c5d34SMatthew G. Knepley } 12296c73c22cSMatthew G. Knepley 123068132eb9SMatthew G. Knepley /*@ 123168132eb9SMatthew G. Knepley DMPlexComputeInterpolatorGeneral - Form the local portion of the interpolation matrix I from the coarse DM to a non-nested fine DM. 123268132eb9SMatthew G. Knepley 123368132eb9SMatthew G. Knepley Input Parameters: 123468132eb9SMatthew G. Knepley + dmf - The fine mesh 123568132eb9SMatthew G. Knepley . dmc - The coarse mesh 123668132eb9SMatthew G. Knepley - user - The user context 123768132eb9SMatthew G. Knepley 123868132eb9SMatthew G. Knepley Output Parameter: 123968132eb9SMatthew G. Knepley . In - The interpolation matrix 124068132eb9SMatthew G. Knepley 124168132eb9SMatthew G. Knepley Level: developer 124268132eb9SMatthew G. Knepley 124368132eb9SMatthew G. Knepley .seealso: DMPlexComputeInterpolatorNested(), DMPlexComputeJacobianFEM() 124468132eb9SMatthew G. Knepley @*/ 124568132eb9SMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorGeneral(DM dmc, DM dmf, Mat In, void *user) 12464ef9d792SMatthew G. Knepley { 124764e98e1dSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dmf->data; 124864e98e1dSMatthew G. Knepley const char *name = "Interpolator"; 12494ef9d792SMatthew G. Knepley PetscDS prob; 12504ef9d792SMatthew G. Knepley PetscSection fsection, csection, globalFSection, globalCSection; 12514ef9d792SMatthew G. Knepley PetscHashJK ht; 12524ef9d792SMatthew G. Knepley PetscLayout rLayout; 12534ef9d792SMatthew G. Knepley PetscInt *dnz, *onz; 12544ef9d792SMatthew G. Knepley PetscInt locRows, rStart, rEnd; 12554ef9d792SMatthew G. Knepley PetscReal *x, *v0, *J, *invJ, detJ; 12564ef9d792SMatthew G. Knepley PetscReal *v0c, *Jc, *invJc, detJc; 12574ef9d792SMatthew G. Knepley PetscScalar *elemMat; 12584ef9d792SMatthew G. Knepley PetscInt dim, Nf, field, totDim, cStart, cEnd, cell, ccell; 12594ef9d792SMatthew G. Knepley PetscErrorCode ierr; 12604ef9d792SMatthew G. Knepley 12614ef9d792SMatthew G. Knepley PetscFunctionBegin; 126277711781SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 12634ef9d792SMatthew G. Knepley ierr = DMGetCoordinateDim(dmc, &dim);CHKERRQ(ierr); 12644ef9d792SMatthew G. Knepley ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr); 12654ef9d792SMatthew G. Knepley ierr = PetscDSGetRefCoordArrays(prob, &x, NULL);CHKERRQ(ierr); 12664ef9d792SMatthew G. Knepley ierr = PetscDSGetNumFields(prob, &Nf);CHKERRQ(ierr); 12674ef9d792SMatthew G. Knepley ierr = PetscMalloc3(dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 12684ef9d792SMatthew G. Knepley ierr = PetscMalloc3(dim,&v0c,dim*dim,&Jc,dim*dim,&invJc);CHKERRQ(ierr); 12694ef9d792SMatthew G. Knepley ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 12704ef9d792SMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dmf, &globalFSection);CHKERRQ(ierr); 12714ef9d792SMatthew G. Knepley ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 12724ef9d792SMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dmc, &globalCSection);CHKERRQ(ierr); 12734ef9d792SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dmf, 0, &cStart, &cEnd);CHKERRQ(ierr); 12744ef9d792SMatthew G. Knepley ierr = PetscDSGetTotalDimension(prob, &totDim);CHKERRQ(ierr); 12759c3cf19fSMatthew G. Knepley ierr = PetscMalloc1(totDim, &elemMat);CHKERRQ(ierr); 12764ef9d792SMatthew G. Knepley 12774ef9d792SMatthew G. Knepley ierr = MatGetLocalSize(In, &locRows, NULL);CHKERRQ(ierr); 12784ef9d792SMatthew G. Knepley ierr = PetscLayoutCreate(PetscObjectComm((PetscObject) In), &rLayout);CHKERRQ(ierr); 12794ef9d792SMatthew G. Knepley ierr = PetscLayoutSetLocalSize(rLayout, locRows);CHKERRQ(ierr); 12804ef9d792SMatthew G. Knepley ierr = PetscLayoutSetBlockSize(rLayout, 1);CHKERRQ(ierr); 12814ef9d792SMatthew G. Knepley ierr = PetscLayoutSetUp(rLayout);CHKERRQ(ierr); 12824ef9d792SMatthew G. Knepley ierr = PetscLayoutGetRange(rLayout, &rStart, &rEnd);CHKERRQ(ierr); 12834ef9d792SMatthew G. Knepley ierr = PetscLayoutDestroy(&rLayout);CHKERRQ(ierr); 12844ef9d792SMatthew G. Knepley ierr = PetscCalloc2(locRows,&dnz,locRows,&onz);CHKERRQ(ierr); 12854ef9d792SMatthew G. Knepley ierr = PetscHashJKCreate(&ht);CHKERRQ(ierr); 12864ef9d792SMatthew G. Knepley for (field = 0; field < Nf; ++field) { 12874ef9d792SMatthew G. Knepley PetscObject obj; 12884ef9d792SMatthew G. Knepley PetscClassId id; 1289c0d7054bSMatthew G. Knepley PetscDualSpace Q = NULL; 12904ef9d792SMatthew G. Knepley PetscQuadrature f; 129117f047d8SMatthew G. Knepley const PetscReal *qpoints; 129217f047d8SMatthew G. Knepley PetscInt Nc, Np, fpdim, i, d; 12934ef9d792SMatthew G. Knepley 12944ef9d792SMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr); 12954ef9d792SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 12964ef9d792SMatthew G. Knepley if (id == PETSCFE_CLASSID) { 12974ef9d792SMatthew G. Knepley PetscFE fe = (PetscFE) obj; 12984ef9d792SMatthew G. Knepley 12994ef9d792SMatthew G. Knepley ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr); 13004ef9d792SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 13014ef9d792SMatthew G. Knepley } else if (id == PETSCFV_CLASSID) { 13024ef9d792SMatthew G. Knepley PetscFV fv = (PetscFV) obj; 13034ef9d792SMatthew G. Knepley 13044ef9d792SMatthew G. Knepley ierr = PetscFVGetDualSpace(fv, &Q);CHKERRQ(ierr); 13054ef9d792SMatthew G. Knepley Nc = 1; 13064ef9d792SMatthew G. Knepley } 13074ef9d792SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(Q, &fpdim);CHKERRQ(ierr); 13084ef9d792SMatthew G. Knepley /* For each fine grid cell */ 13094ef9d792SMatthew G. Knepley for (cell = cStart; cell < cEnd; ++cell) { 13104ef9d792SMatthew G. Knepley PetscInt *findices, *cindices; 13114ef9d792SMatthew G. Knepley PetscInt numFIndices, numCIndices; 13124ef9d792SMatthew G. Knepley 13136ecaa68aSToby Isaac ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr); 13144ef9d792SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 13159c3cf19fSMatthew G. Knepley if (numFIndices != fpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %d != %d dual basis vecs", numFIndices, fpdim); 13164ef9d792SMatthew G. Knepley for (i = 0; i < fpdim; ++i) { 13174ef9d792SMatthew G. Knepley Vec pointVec; 13184ef9d792SMatthew G. Knepley PetscScalar *pV; 13193a93e3b7SToby Isaac PetscSF coarseCellSF = NULL; 13203a93e3b7SToby Isaac const PetscSFNode *coarseCells; 13219c3cf19fSMatthew G. Knepley PetscInt numCoarseCells, q, c; 13224ef9d792SMatthew G. Knepley 13234ef9d792SMatthew G. Knepley /* Get points from the dual basis functional quadrature */ 13244ef9d792SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(Q, i, &f);CHKERRQ(ierr); 13259c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(f, NULL, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr); 13264ef9d792SMatthew G. Knepley ierr = VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);CHKERRQ(ierr); 13274ef9d792SMatthew G. Knepley ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr); 13284ef9d792SMatthew G. Knepley ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr); 13294ef9d792SMatthew G. Knepley for (q = 0; q < Np; ++q) { 13304ef9d792SMatthew G. Knepley /* Transform point to real space */ 13314ef9d792SMatthew G. Knepley CoordinatesRefToReal(dim, dim, v0, J, &qpoints[q*dim], x); 13324ef9d792SMatthew G. Knepley for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d]; 13334ef9d792SMatthew G. Knepley } 13344ef9d792SMatthew G. Knepley ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr); 13354ef9d792SMatthew G. Knepley /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */ 133662a38674SMatthew G. Knepley ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr); 13373a93e3b7SToby Isaac ierr = PetscSFViewFromOptions(coarseCellSF, NULL, "-interp_sf_view");CHKERRQ(ierr); 13384ef9d792SMatthew G. Knepley /* Update preallocation info */ 13393a93e3b7SToby Isaac ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr); 13403a93e3b7SToby Isaac if (numCoarseCells != Np) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located"); 13419c3cf19fSMatthew G. Knepley { 13424ef9d792SMatthew G. Knepley PetscHashJKKey key; 13434ef9d792SMatthew G. Knepley PetscHashJKIter missing, iter; 13444ef9d792SMatthew G. Knepley 13459c3cf19fSMatthew G. Knepley key.j = findices[i]; 13468c543595SMatthew G. Knepley if (key.j >= 0) { 13474ef9d792SMatthew G. Knepley /* Get indices for coarse elements */ 13484ef9d792SMatthew G. Knepley for (ccell = 0; ccell < numCoarseCells; ++ccell) { 13493a93e3b7SToby Isaac ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr); 13504ef9d792SMatthew G. Knepley for (c = 0; c < numCIndices; ++c) { 13514ef9d792SMatthew G. Knepley key.k = cindices[c]; 13524ef9d792SMatthew G. Knepley if (key.k < 0) continue; 13534ef9d792SMatthew G. Knepley ierr = PetscHashJKPut(ht, key, &missing, &iter);CHKERRQ(ierr); 13544ef9d792SMatthew G. Knepley if (missing) { 13554ef9d792SMatthew G. Knepley ierr = PetscHashJKSet(ht, iter, 1);CHKERRQ(ierr); 13564ef9d792SMatthew G. Knepley if ((key.k >= rStart) && (key.k < rEnd)) ++dnz[key.j-rStart]; 13574ef9d792SMatthew G. Knepley else ++onz[key.j-rStart]; 13584ef9d792SMatthew G. Knepley } 13594ef9d792SMatthew G. Knepley } 13603a93e3b7SToby Isaac ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr); 13614ef9d792SMatthew G. Knepley } 13624ef9d792SMatthew G. Knepley } 13638c543595SMatthew G. Knepley } 13643a93e3b7SToby Isaac ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr); 13654ef9d792SMatthew G. Knepley ierr = VecDestroy(&pointVec);CHKERRQ(ierr); 13664ef9d792SMatthew G. Knepley } 136746bdb399SToby Isaac ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr); 13684ef9d792SMatthew G. Knepley } 13694ef9d792SMatthew G. Knepley } 13704ef9d792SMatthew G. Knepley ierr = PetscHashJKDestroy(&ht);CHKERRQ(ierr); 13714ef9d792SMatthew G. Knepley ierr = MatXAIJSetPreallocation(In, 1, dnz, onz, NULL, NULL);CHKERRQ(ierr); 13724ef9d792SMatthew G. Knepley ierr = MatSetOption(In, MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 13734ef9d792SMatthew G. Knepley ierr = PetscFree2(dnz,onz);CHKERRQ(ierr); 13744ef9d792SMatthew G. Knepley for (field = 0; field < Nf; ++field) { 13754ef9d792SMatthew G. Knepley PetscObject obj; 13764ef9d792SMatthew G. Knepley PetscClassId id; 1377c0d7054bSMatthew G. Knepley PetscDualSpace Q = NULL; 13784ef9d792SMatthew G. Knepley PetscQuadrature f; 13794ef9d792SMatthew G. Knepley const PetscReal *qpoints, *qweights; 13809c3cf19fSMatthew G. Knepley PetscInt Nc, qNc, Np, fpdim, i, d; 13814ef9d792SMatthew G. Knepley 13824ef9d792SMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, field, &obj);CHKERRQ(ierr); 13834ef9d792SMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 13844ef9d792SMatthew G. Knepley if (id == PETSCFE_CLASSID) { 13854ef9d792SMatthew G. Knepley PetscFE fe = (PetscFE) obj; 13864ef9d792SMatthew G. Knepley 13874ef9d792SMatthew G. Knepley ierr = PetscFEGetDualSpace(fe, &Q);CHKERRQ(ierr); 13884ef9d792SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 13894ef9d792SMatthew G. Knepley } else if (id == PETSCFV_CLASSID) { 13904ef9d792SMatthew G. Knepley PetscFV fv = (PetscFV) obj; 13914ef9d792SMatthew G. Knepley 13924ef9d792SMatthew G. Knepley ierr = PetscFVGetDualSpace(fv, &Q);CHKERRQ(ierr); 13934ef9d792SMatthew G. Knepley Nc = 1; 13944ef9d792SMatthew G. Knepley } 13954ef9d792SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(Q, &fpdim);CHKERRQ(ierr); 13964ef9d792SMatthew G. Knepley /* For each fine grid cell */ 13974ef9d792SMatthew G. Knepley for (cell = cStart; cell < cEnd; ++cell) { 13984ef9d792SMatthew G. Knepley PetscInt *findices, *cindices; 13994ef9d792SMatthew G. Knepley PetscInt numFIndices, numCIndices; 14004ef9d792SMatthew G. Knepley 14016ecaa68aSToby Isaac ierr = DMPlexGetClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr); 14024ef9d792SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dmf, cell, NULL, v0, J, invJ, &detJ);CHKERRQ(ierr); 14039c3cf19fSMatthew G. Knepley if (numFIndices != fpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of fine indices %d != %d dual basis vecs", numFIndices, fpdim); 14044ef9d792SMatthew G. Knepley for (i = 0; i < fpdim; ++i) { 14054ef9d792SMatthew G. Knepley Vec pointVec; 14064ef9d792SMatthew G. Knepley PetscScalar *pV; 140712111d7cSToby Isaac PetscSF coarseCellSF = NULL; 14083a93e3b7SToby Isaac const PetscSFNode *coarseCells; 140917f047d8SMatthew G. Knepley PetscInt numCoarseCells, cpdim, q, c, j; 14104ef9d792SMatthew G. Knepley 14114ef9d792SMatthew G. Knepley /* Get points from the dual basis functional quadrature */ 14124ef9d792SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(Q, i, &f);CHKERRQ(ierr); 14139c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(f, NULL, &qNc, &Np, &qpoints, &qweights);CHKERRQ(ierr); 14149c3cf19fSMatthew G. Knepley if (qNc != Nc) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Number of components in quadrature %D does not match coarse field %D", qNc, Nc); 14154ef9d792SMatthew G. Knepley ierr = VecCreateSeq(PETSC_COMM_SELF, Np*dim, &pointVec);CHKERRQ(ierr); 14164ef9d792SMatthew G. Knepley ierr = VecSetBlockSize(pointVec, dim);CHKERRQ(ierr); 14174ef9d792SMatthew G. Knepley ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr); 14184ef9d792SMatthew G. Knepley for (q = 0; q < Np; ++q) { 14194ef9d792SMatthew G. Knepley /* Transform point to real space */ 14204ef9d792SMatthew G. Knepley CoordinatesRefToReal(dim, dim, v0, J, &qpoints[q*dim], x); 14214ef9d792SMatthew G. Knepley for (d = 0; d < dim; ++d) pV[q*dim+d] = x[d]; 14224ef9d792SMatthew G. Knepley } 14234ef9d792SMatthew G. Knepley ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr); 14244ef9d792SMatthew G. Knepley /* Get set of coarse cells that overlap points (would like to group points by coarse cell) */ 142562a38674SMatthew G. Knepley ierr = DMLocatePoints(dmc, pointVec, DM_POINTLOCATION_NEAREST, &coarseCellSF);CHKERRQ(ierr); 14264ef9d792SMatthew G. Knepley /* Update preallocation info */ 14273a93e3b7SToby Isaac ierr = PetscSFGetGraph(coarseCellSF, NULL, &numCoarseCells, NULL, &coarseCells);CHKERRQ(ierr); 14283a93e3b7SToby Isaac if (numCoarseCells != Np) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Not all closure points located"); 14294ef9d792SMatthew G. Knepley ierr = VecGetArray(pointVec, &pV);CHKERRQ(ierr); 14304ef9d792SMatthew G. Knepley for (ccell = 0; ccell < numCoarseCells; ++ccell) { 1431826eb36dSMatthew G. Knepley PetscReal pVReal[3]; 1432826eb36dSMatthew G. Knepley 14333a93e3b7SToby Isaac ierr = DMPlexGetClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr); 14344ef9d792SMatthew G. Knepley /* Transform points from real space to coarse reference space */ 14353a93e3b7SToby Isaac ierr = DMPlexComputeCellGeometryFEM(dmc, coarseCells[ccell].index, NULL, v0c, Jc, invJc, &detJc);CHKERRQ(ierr); 1436e2d86523SMatthew G. Knepley for (d = 0; d < dim; ++d) pVReal[d] = PetscRealPart(pV[ccell*dim+d]); 1437826eb36dSMatthew G. Knepley CoordinatesRealToRef(dim, dim, v0c, invJc, pVReal, x); 14384ef9d792SMatthew G. Knepley 14394ef9d792SMatthew G. Knepley if (id == PETSCFE_CLASSID) { 14404ef9d792SMatthew G. Knepley PetscFE fe = (PetscFE) obj; 14414ef9d792SMatthew G. Knepley PetscReal *B; 14424ef9d792SMatthew G. Knepley 14434ef9d792SMatthew G. Knepley /* Evaluate coarse basis on contained point */ 14444ef9d792SMatthew G. Knepley ierr = PetscFEGetDimension(fe, &cpdim);CHKERRQ(ierr); 14454ef9d792SMatthew G. Knepley ierr = PetscFEGetTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr); 14469c3cf19fSMatthew G. Knepley ierr = PetscMemzero(elemMat, cpdim * sizeof(PetscScalar));CHKERRQ(ierr); 14474ef9d792SMatthew G. Knepley /* Get elemMat entries by multiplying by weight */ 14484ef9d792SMatthew G. Knepley for (j = 0; j < cpdim; ++j) { 14499c3cf19fSMatthew G. Knepley for (c = 0; c < Nc; ++c) elemMat[j] += B[j*Nc + c]*qweights[ccell*qNc + c]; 14504ef9d792SMatthew G. Knepley } 14514ef9d792SMatthew G. Knepley ierr = PetscFERestoreTabulation(fe, 1, x, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 14524ef9d792SMatthew G. Knepley } else { 14534ef9d792SMatthew G. Knepley cpdim = 1; 14544ef9d792SMatthew G. Knepley for (j = 0; j < cpdim; ++j) { 14559c3cf19fSMatthew G. Knepley for (c = 0; c < Nc; ++c) elemMat[j] += 1.0*qweights[ccell*qNc + c]; 14564ef9d792SMatthew G. Knepley } 14574ef9d792SMatthew G. Knepley } 14584ef9d792SMatthew G. Knepley /* Update interpolator */ 14599c3cf19fSMatthew G. Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(cell, name, 1, numCIndices, elemMat);CHKERRQ(ierr);} 14609c3cf19fSMatthew G. Knepley if (numCIndices != cpdim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Number of element matrix columns %D != %D", numCIndices, cpdim); 14619c3cf19fSMatthew G. Knepley ierr = MatSetValues(In, 1, &findices[i], numCIndices, cindices, elemMat, INSERT_VALUES);CHKERRQ(ierr); 14623a93e3b7SToby Isaac ierr = DMPlexRestoreClosureIndices(dmc, csection, globalCSection, coarseCells[ccell].index, &numCIndices, &cindices, NULL);CHKERRQ(ierr); 14634ef9d792SMatthew G. Knepley } 14644ef9d792SMatthew G. Knepley ierr = VecRestoreArray(pointVec, &pV);CHKERRQ(ierr); 14653a93e3b7SToby Isaac ierr = PetscSFDestroy(&coarseCellSF);CHKERRQ(ierr); 14664ef9d792SMatthew G. Knepley ierr = VecDestroy(&pointVec);CHKERRQ(ierr); 14674ef9d792SMatthew G. Knepley } 146846bdb399SToby Isaac ierr = DMPlexRestoreClosureIndices(dmf, fsection, globalFSection, cell, &numFIndices, &findices, NULL);CHKERRQ(ierr); 14694ef9d792SMatthew G. Knepley } 14704ef9d792SMatthew G. Knepley } 14714ef9d792SMatthew G. Knepley ierr = PetscFree3(v0,J,invJ);CHKERRQ(ierr); 14724ef9d792SMatthew G. Knepley ierr = PetscFree3(v0c,Jc,invJc);CHKERRQ(ierr); 14734ef9d792SMatthew G. Knepley ierr = PetscFree(elemMat);CHKERRQ(ierr); 14744ef9d792SMatthew G. Knepley ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14754ef9d792SMatthew G. Knepley ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 147677711781SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 14774ef9d792SMatthew G. Knepley PetscFunctionReturn(0); 14784ef9d792SMatthew G. Knepley } 14794ef9d792SMatthew G. Knepley 14807c927364SMatthew G. Knepley PetscErrorCode DMPlexComputeInjectorFEM(DM dmc, DM dmf, VecScatter *sc, void *user) 14817c927364SMatthew G. Knepley { 1482e9d4ef1bSMatthew G. Knepley PetscDS prob; 14837c927364SMatthew G. Knepley PetscFE *feRef; 148497c42addSMatthew G. Knepley PetscFV *fvRef; 14857c927364SMatthew G. Knepley Vec fv, cv; 14867c927364SMatthew G. Knepley IS fis, cis; 14877c927364SMatthew G. Knepley PetscSection fsection, fglobalSection, csection, cglobalSection; 14887c927364SMatthew G. Knepley PetscInt *cmap, *cellCIndices, *cellFIndices, *cindices, *findices; 14890bd915a7SMatthew G. Knepley PetscInt cTotDim, fTotDim = 0, Nf, f, field, cStart, cEnd, cEndInterior, c, dim, d, startC, endC, offsetC, offsetF, m; 14907c927364SMatthew G. Knepley PetscErrorCode ierr; 14917c927364SMatthew G. Knepley 14927c927364SMatthew G. Knepley PetscFunctionBegin; 149375a69067SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1494c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dmf, &dim);CHKERRQ(ierr); 14957c927364SMatthew G. Knepley ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 14967c927364SMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 14977c927364SMatthew G. Knepley ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 14987c927364SMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 14997c927364SMatthew G. Knepley ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 15007c927364SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 15019ac3fadcSMatthew G. Knepley ierr = DMPlexGetHybridBounds(dmc, &cEndInterior, NULL, NULL, NULL);CHKERRQ(ierr); 15029ac3fadcSMatthew G. Knepley cEnd = cEndInterior < 0 ? cEnd : cEndInterior; 1503e9d4ef1bSMatthew G. Knepley ierr = DMGetDS(dmc, &prob);CHKERRQ(ierr); 150497c42addSMatthew G. Knepley ierr = PetscCalloc2(Nf,&feRef,Nf,&fvRef);CHKERRQ(ierr); 15057c927364SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 150697c42addSMatthew G. Knepley PetscObject obj; 150797c42addSMatthew G. Knepley PetscClassId id; 1508aa7890ccSMatthew G. Knepley PetscInt fNb = 0, Nc = 0; 15097c927364SMatthew G. Knepley 151097c42addSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, f, &obj);CHKERRQ(ierr); 151197c42addSMatthew G. Knepley ierr = PetscObjectGetClassId(obj, &id);CHKERRQ(ierr); 151297c42addSMatthew G. Knepley if (id == PETSCFE_CLASSID) { 151397c42addSMatthew G. Knepley PetscFE fe = (PetscFE) obj; 151497c42addSMatthew G. Knepley 15157c927364SMatthew G. Knepley ierr = PetscFERefine(fe, &feRef[f]);CHKERRQ(ierr); 15167c927364SMatthew G. Knepley ierr = PetscFEGetDimension(feRef[f], &fNb);CHKERRQ(ierr); 15177c927364SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe, &Nc);CHKERRQ(ierr); 151897c42addSMatthew G. Knepley } else if (id == PETSCFV_CLASSID) { 151997c42addSMatthew G. Knepley PetscFV fv = (PetscFV) obj; 152097c42addSMatthew G. Knepley PetscDualSpace Q; 152197c42addSMatthew G. Knepley 152297c42addSMatthew G. Knepley ierr = PetscFVRefine(fv, &fvRef[f]);CHKERRQ(ierr); 152397c42addSMatthew G. Knepley ierr = PetscFVGetDualSpace(fvRef[f], &Q);CHKERRQ(ierr); 152497c42addSMatthew G. Knepley ierr = PetscDualSpaceGetDimension(Q, &fNb);CHKERRQ(ierr); 152597c42addSMatthew G. Knepley ierr = PetscFVGetNumComponents(fv, &Nc);CHKERRQ(ierr); 152697c42addSMatthew G. Knepley } 15277c927364SMatthew G. Knepley fTotDim += fNb*Nc; 15287c927364SMatthew G. Knepley } 1529e9d4ef1bSMatthew G. Knepley ierr = PetscDSGetTotalDimension(prob, &cTotDim);CHKERRQ(ierr); 15307c927364SMatthew G. Knepley ierr = PetscMalloc1(cTotDim,&cmap);CHKERRQ(ierr); 15317c927364SMatthew G. Knepley for (field = 0, offsetC = 0, offsetF = 0; field < Nf; ++field) { 15327c927364SMatthew G. Knepley PetscFE feC; 153397c42addSMatthew G. Knepley PetscFV fvC; 15347c927364SMatthew G. Knepley PetscDualSpace QF, QC; 15357c927364SMatthew G. Knepley PetscInt NcF, NcC, fpdim, cpdim; 15367c927364SMatthew G. Knepley 153797c42addSMatthew G. Knepley if (feRef[field]) { 1538e9d4ef1bSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &feC);CHKERRQ(ierr); 15397c927364SMatthew G. Knepley ierr = PetscFEGetNumComponents(feC, &NcC);CHKERRQ(ierr); 15407c927364SMatthew G. Knepley ierr = PetscFEGetNumComponents(feRef[field], &NcF);CHKERRQ(ierr); 15417c927364SMatthew G. Knepley ierr = PetscFEGetDualSpace(feRef[field], &QF);CHKERRQ(ierr); 15427c927364SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr); 15437c927364SMatthew G. Knepley ierr = PetscFEGetDualSpace(feC, &QC);CHKERRQ(ierr); 15447c927364SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr); 154597c42addSMatthew G. Knepley } else { 154697c42addSMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, field, (PetscObject *) &fvC);CHKERRQ(ierr); 154797c42addSMatthew G. Knepley ierr = PetscFVGetNumComponents(fvC, &NcC);CHKERRQ(ierr); 154897c42addSMatthew G. Knepley ierr = PetscFVGetNumComponents(fvRef[field], &NcF);CHKERRQ(ierr); 154997c42addSMatthew G. Knepley ierr = PetscFVGetDualSpace(fvRef[field], &QF);CHKERRQ(ierr); 155097c42addSMatthew G. Knepley ierr = PetscDualSpaceGetDimension(QF, &fpdim);CHKERRQ(ierr); 155197c42addSMatthew G. Knepley ierr = PetscFVGetDualSpace(fvC, &QC);CHKERRQ(ierr); 155297c42addSMatthew G. Knepley ierr = PetscDualSpaceGetDimension(QC, &cpdim);CHKERRQ(ierr); 155397c42addSMatthew G. Knepley } 155497c42addSMatthew 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); 15557c927364SMatthew G. Knepley for (c = 0; c < cpdim; ++c) { 15567c927364SMatthew G. Knepley PetscQuadrature cfunc; 15577c927364SMatthew G. Knepley const PetscReal *cqpoints; 15587c927364SMatthew G. Knepley PetscInt NpC; 155997c42addSMatthew G. Knepley PetscBool found = PETSC_FALSE; 15607c927364SMatthew G. Knepley 15617c927364SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(QC, c, &cfunc);CHKERRQ(ierr); 15629c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(cfunc, NULL, NULL, &NpC, &cqpoints, NULL);CHKERRQ(ierr); 156397c42addSMatthew G. Knepley if (NpC != 1 && feRef[field]) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Do not know how to do injection for moments"); 15647c927364SMatthew G. Knepley for (f = 0; f < fpdim; ++f) { 15657c927364SMatthew G. Knepley PetscQuadrature ffunc; 15667c927364SMatthew G. Knepley const PetscReal *fqpoints; 15677c927364SMatthew G. Knepley PetscReal sum = 0.0; 15687c927364SMatthew G. Knepley PetscInt NpF, comp; 15697c927364SMatthew G. Knepley 15707c927364SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(QF, f, &ffunc);CHKERRQ(ierr); 15719c3cf19fSMatthew G. Knepley ierr = PetscQuadratureGetData(ffunc, NULL, NULL, &NpF, &fqpoints, NULL);CHKERRQ(ierr); 15727c927364SMatthew G. Knepley if (NpC != NpF) continue; 15737c927364SMatthew G. Knepley for (d = 0; d < dim; ++d) sum += PetscAbsReal(cqpoints[d] - fqpoints[d]); 15747c927364SMatthew G. Knepley if (sum > 1.0e-9) continue; 15757c927364SMatthew G. Knepley for (comp = 0; comp < NcC; ++comp) { 15767c927364SMatthew G. Knepley cmap[(offsetC+c)*NcC+comp] = (offsetF+f)*NcF+comp; 15777c927364SMatthew G. Knepley } 157897c42addSMatthew G. Knepley found = PETSC_TRUE; 15797c927364SMatthew G. Knepley break; 15807c927364SMatthew G. Knepley } 158197c42addSMatthew G. Knepley if (!found) { 158297c42addSMatthew G. Knepley /* TODO We really want the average here, but some asshole put VecScatter in the interface */ 158397c42addSMatthew G. Knepley if (fvRef[field]) { 158497c42addSMatthew G. Knepley PetscInt comp; 158597c42addSMatthew G. Knepley for (comp = 0; comp < NcC; ++comp) { 158697c42addSMatthew G. Knepley cmap[(offsetC+c)*NcC+comp] = (offsetF+0)*NcF+comp; 158797c42addSMatthew G. Knepley } 158897c42addSMatthew G. Knepley } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Could not locate matching functional for injection"); 158997c42addSMatthew G. Knepley } 15907c927364SMatthew G. Knepley } 15917c927364SMatthew G. Knepley offsetC += cpdim*NcC; 15927c927364SMatthew G. Knepley offsetF += fpdim*NcF; 15937c927364SMatthew G. Knepley } 159497c42addSMatthew G. Knepley for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);ierr = PetscFVDestroy(&fvRef[f]);CHKERRQ(ierr);} 159597c42addSMatthew G. Knepley ierr = PetscFree2(feRef,fvRef);CHKERRQ(ierr); 15967c927364SMatthew G. Knepley 15977c927364SMatthew G. Knepley ierr = DMGetGlobalVector(dmf, &fv);CHKERRQ(ierr); 15987c927364SMatthew G. Knepley ierr = DMGetGlobalVector(dmc, &cv);CHKERRQ(ierr); 15990bd915a7SMatthew G. Knepley ierr = VecGetOwnershipRange(cv, &startC, &endC);CHKERRQ(ierr); 16007c927364SMatthew G. Knepley ierr = PetscSectionGetConstrainedStorageSize(cglobalSection, &m);CHKERRQ(ierr); 16017c927364SMatthew G. Knepley ierr = PetscMalloc2(cTotDim,&cellCIndices,fTotDim,&cellFIndices);CHKERRQ(ierr); 1602aa7da3c4SMatthew G. Knepley ierr = PetscMalloc1(m,&cindices);CHKERRQ(ierr); 1603aa7da3c4SMatthew G. Knepley ierr = PetscMalloc1(m,&findices);CHKERRQ(ierr); 16047c927364SMatthew G. Knepley for (d = 0; d < m; ++d) cindices[d] = findices[d] = -1; 16057c927364SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 16067c927364SMatthew G. Knepley ierr = DMPlexMatGetClosureIndicesRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, c, cellCIndices, cellFIndices);CHKERRQ(ierr); 16077c927364SMatthew G. Knepley for (d = 0; d < cTotDim; ++d) { 16080bd915a7SMatthew G. Knepley if ((cellCIndices[d] < startC) || (cellCIndices[d] >= endC)) continue; 16097c927364SMatthew 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]]); 16107c927364SMatthew G. Knepley cindices[cellCIndices[d]-startC] = cellCIndices[d]; 16117c927364SMatthew G. Knepley findices[cellCIndices[d]-startC] = cellFIndices[cmap[d]]; 16127c927364SMatthew G. Knepley } 16137c927364SMatthew G. Knepley } 16147c927364SMatthew G. Knepley ierr = PetscFree(cmap);CHKERRQ(ierr); 16157c927364SMatthew G. Knepley ierr = PetscFree2(cellCIndices,cellFIndices);CHKERRQ(ierr); 16167c927364SMatthew G. Knepley 16177c927364SMatthew G. Knepley ierr = ISCreateGeneral(PETSC_COMM_SELF, m, cindices, PETSC_OWN_POINTER, &cis);CHKERRQ(ierr); 16187c927364SMatthew G. Knepley ierr = ISCreateGeneral(PETSC_COMM_SELF, m, findices, PETSC_OWN_POINTER, &fis);CHKERRQ(ierr); 16197c927364SMatthew G. Knepley ierr = VecScatterCreate(cv, cis, fv, fis, sc);CHKERRQ(ierr); 16207c927364SMatthew G. Knepley ierr = ISDestroy(&cis);CHKERRQ(ierr); 16217c927364SMatthew G. Knepley ierr = ISDestroy(&fis);CHKERRQ(ierr); 16227c927364SMatthew G. Knepley ierr = DMRestoreGlobalVector(dmf, &fv);CHKERRQ(ierr); 16237c927364SMatthew G. Knepley ierr = DMRestoreGlobalVector(dmc, &cv);CHKERRQ(ierr); 162475a69067SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_InjectorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1625cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 1626cb1e1211SMatthew G Knepley } 1627