1cb1e1211SMatthew G Knepley #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2cb1e1211SMatthew G Knepley 3a0845e3aSMatthew G. Knepley #include <petscfe.h> 4a0845e3aSMatthew G. Knepley 5cb1e1211SMatthew G Knepley #undef __FUNCT__ 6cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexGetScale" 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 #undef __FUNCT__ 19cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexSetScale" 20cb1e1211SMatthew G Knepley PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale) 21cb1e1211SMatthew G Knepley { 22cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 23cb1e1211SMatthew G Knepley 24cb1e1211SMatthew G Knepley PetscFunctionBegin; 25cb1e1211SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 26cb1e1211SMatthew G Knepley mesh->scale[unit] = scale; 27cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 28cb1e1211SMatthew G Knepley } 29cb1e1211SMatthew G Knepley 30cb1e1211SMatthew G Knepley PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k) 31cb1e1211SMatthew G Knepley { 32cb1e1211SMatthew G Knepley switch (i) { 33cb1e1211SMatthew G Knepley case 0: 34cb1e1211SMatthew G Knepley switch (j) { 35cb1e1211SMatthew G Knepley case 0: return 0; 36cb1e1211SMatthew G Knepley case 1: 37cb1e1211SMatthew G Knepley switch (k) { 38cb1e1211SMatthew G Knepley case 0: return 0; 39cb1e1211SMatthew G Knepley case 1: return 0; 40cb1e1211SMatthew G Knepley case 2: return 1; 41cb1e1211SMatthew G Knepley } 42cb1e1211SMatthew G Knepley case 2: 43cb1e1211SMatthew G Knepley switch (k) { 44cb1e1211SMatthew G Knepley case 0: return 0; 45cb1e1211SMatthew G Knepley case 1: return -1; 46cb1e1211SMatthew G Knepley case 2: return 0; 47cb1e1211SMatthew G Knepley } 48cb1e1211SMatthew G Knepley } 49cb1e1211SMatthew G Knepley case 1: 50cb1e1211SMatthew G Knepley switch (j) { 51cb1e1211SMatthew G Knepley case 0: 52cb1e1211SMatthew G Knepley switch (k) { 53cb1e1211SMatthew G Knepley case 0: return 0; 54cb1e1211SMatthew G Knepley case 1: return 0; 55cb1e1211SMatthew G Knepley case 2: return -1; 56cb1e1211SMatthew G Knepley } 57cb1e1211SMatthew G Knepley case 1: return 0; 58cb1e1211SMatthew G Knepley case 2: 59cb1e1211SMatthew G Knepley switch (k) { 60cb1e1211SMatthew G Knepley case 0: return 1; 61cb1e1211SMatthew G Knepley case 1: return 0; 62cb1e1211SMatthew G Knepley case 2: return 0; 63cb1e1211SMatthew G Knepley } 64cb1e1211SMatthew G Knepley } 65cb1e1211SMatthew G Knepley case 2: 66cb1e1211SMatthew G Knepley switch (j) { 67cb1e1211SMatthew G Knepley case 0: 68cb1e1211SMatthew G Knepley switch (k) { 69cb1e1211SMatthew G Knepley case 0: return 0; 70cb1e1211SMatthew G Knepley case 1: return 1; 71cb1e1211SMatthew G Knepley case 2: return 0; 72cb1e1211SMatthew G Knepley } 73cb1e1211SMatthew G Knepley case 1: 74cb1e1211SMatthew G Knepley switch (k) { 75cb1e1211SMatthew G Knepley case 0: return -1; 76cb1e1211SMatthew G Knepley case 1: return 0; 77cb1e1211SMatthew G Knepley case 2: return 0; 78cb1e1211SMatthew G Knepley } 79cb1e1211SMatthew G Knepley case 2: return 0; 80cb1e1211SMatthew G Knepley } 81cb1e1211SMatthew G Knepley } 82cb1e1211SMatthew G Knepley return 0; 83cb1e1211SMatthew G Knepley } 84cb1e1211SMatthew G Knepley 85cb1e1211SMatthew G Knepley #undef __FUNCT__ 86cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexCreateRigidBody" 87cb1e1211SMatthew G Knepley /*@C 88cb1e1211SMatthew G Knepley DMPlexCreateRigidBody - create rigid body modes from coordinates 89cb1e1211SMatthew G Knepley 90cb1e1211SMatthew G Knepley Collective on DM 91cb1e1211SMatthew G Knepley 92cb1e1211SMatthew G Knepley Input Arguments: 93cb1e1211SMatthew G Knepley + dm - the DM 94cb1e1211SMatthew G Knepley . section - the local section associated with the rigid field, or NULL for the default section 95cb1e1211SMatthew G Knepley - globalSection - the global section associated with the rigid field, or NULL for the default section 96cb1e1211SMatthew G Knepley 97cb1e1211SMatthew G Knepley Output Argument: 98cb1e1211SMatthew G Knepley . sp - the null space 99cb1e1211SMatthew G Knepley 100cb1e1211SMatthew G Knepley Note: This is necessary to take account of Dirichlet conditions on the displacements 101cb1e1211SMatthew G Knepley 102cb1e1211SMatthew G Knepley Level: advanced 103cb1e1211SMatthew G Knepley 104cb1e1211SMatthew G Knepley .seealso: MatNullSpaceCreate() 105cb1e1211SMatthew G Knepley @*/ 106cb1e1211SMatthew G Knepley PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp) 107cb1e1211SMatthew G Knepley { 108cb1e1211SMatthew G Knepley MPI_Comm comm; 109cb1e1211SMatthew G Knepley Vec coordinates, localMode, mode[6]; 110cb1e1211SMatthew G Knepley PetscSection coordSection; 111cb1e1211SMatthew G Knepley PetscScalar *coords; 112cb1e1211SMatthew G Knepley PetscInt dim, vStart, vEnd, v, n, m, d, i, j; 113cb1e1211SMatthew G Knepley PetscErrorCode ierr; 114cb1e1211SMatthew G Knepley 115cb1e1211SMatthew G Knepley PetscFunctionBegin; 116cb1e1211SMatthew G Knepley ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 117cb1e1211SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 118cb1e1211SMatthew G Knepley if (dim == 1) { 119cb1e1211SMatthew G Knepley ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr); 120cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 121cb1e1211SMatthew G Knepley } 122cb1e1211SMatthew G Knepley if (!section) {ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr);} 123cb1e1211SMatthew G Knepley if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);} 124cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr); 125cb1e1211SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 12669d8a9ceSMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 127cb1e1211SMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 128cb1e1211SMatthew G Knepley m = (dim*(dim+1))/2; 129cb1e1211SMatthew G Knepley ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr); 130cb1e1211SMatthew G Knepley ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr); 131cb1e1211SMatthew G Knepley ierr = VecSetUp(mode[0]);CHKERRQ(ierr); 132cb1e1211SMatthew G Knepley for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);} 133cb1e1211SMatthew G Knepley /* Assume P1 */ 134cb1e1211SMatthew G Knepley ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr); 135cb1e1211SMatthew G Knepley for (d = 0; d < dim; ++d) { 136cb1e1211SMatthew G Knepley PetscScalar values[3] = {0.0, 0.0, 0.0}; 137cb1e1211SMatthew G Knepley 138cb1e1211SMatthew G Knepley values[d] = 1.0; 139cb1e1211SMatthew G Knepley ierr = VecSet(localMode, 0.0);CHKERRQ(ierr); 140cb1e1211SMatthew G Knepley for (v = vStart; v < vEnd; ++v) { 141cb1e1211SMatthew G Knepley ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr); 142cb1e1211SMatthew G Knepley } 143cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 144cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 145cb1e1211SMatthew G Knepley } 146cb1e1211SMatthew G Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 147cb1e1211SMatthew G Knepley for (d = dim; d < dim*(dim+1)/2; ++d) { 148cb1e1211SMatthew G Knepley PetscInt i, j, k = dim > 2 ? d - dim : d; 149cb1e1211SMatthew G Knepley 150cb1e1211SMatthew G Knepley ierr = VecSet(localMode, 0.0);CHKERRQ(ierr); 151cb1e1211SMatthew G Knepley for (v = vStart; v < vEnd; ++v) { 152cb1e1211SMatthew G Knepley PetscScalar values[3] = {0.0, 0.0, 0.0}; 153cb1e1211SMatthew G Knepley PetscInt off; 154cb1e1211SMatthew G Knepley 155cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 156cb1e1211SMatthew G Knepley for (i = 0; i < dim; ++i) { 157cb1e1211SMatthew G Knepley for (j = 0; j < dim; ++j) { 158cb1e1211SMatthew G Knepley values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]); 159cb1e1211SMatthew G Knepley } 160cb1e1211SMatthew G Knepley } 161cb1e1211SMatthew G Knepley ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr); 162cb1e1211SMatthew G Knepley } 163cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 164cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 165cb1e1211SMatthew G Knepley } 166cb1e1211SMatthew G Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 167cb1e1211SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr); 168cb1e1211SMatthew G Knepley for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);} 169cb1e1211SMatthew G Knepley /* Orthonormalize system */ 170cb1e1211SMatthew G Knepley for (i = dim; i < m; ++i) { 171cb1e1211SMatthew G Knepley PetscScalar dots[6]; 172cb1e1211SMatthew G Knepley 173cb1e1211SMatthew G Knepley ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr); 174cb1e1211SMatthew G Knepley for (j = 0; j < i; ++j) dots[j] *= -1.0; 175cb1e1211SMatthew G Knepley ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr); 176cb1e1211SMatthew G Knepley ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr); 177cb1e1211SMatthew G Knepley } 178cb1e1211SMatthew G Knepley ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr); 179cb1e1211SMatthew G Knepley for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);} 180cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 181cb1e1211SMatthew G Knepley } 182cb1e1211SMatthew G Knepley 183cb1e1211SMatthew G Knepley #undef __FUNCT__ 184cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunctionLocal" 185c110b1eeSGeoffrey Irving PetscErrorCode DMPlexProjectFunctionLocal(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec localX) 186cb1e1211SMatthew G Knepley { 18772f94c41SMatthew G. Knepley PetscDualSpace *sp; 18872f94c41SMatthew G. Knepley PetscSection section; 18972f94c41SMatthew G. Knepley PetscScalar *values; 19072f94c41SMatthew G. Knepley PetscReal *v0, *J, detJ; 19172f94c41SMatthew G. Knepley PetscInt numFields, numComp, dim, spDim, totDim = 0, numValues, cStart, cEnd, c, f, d, v; 192cb1e1211SMatthew G Knepley PetscErrorCode ierr; 193cb1e1211SMatthew G Knepley 194cb1e1211SMatthew G Knepley PetscFunctionBegin; 195cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 19672f94c41SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 197785e854fSJed Brown ierr = PetscMalloc1(numFields, &sp);CHKERRQ(ierr); 19872f94c41SMatthew G. Knepley for (f = 0; f < numFields; ++f) { 19972f94c41SMatthew G. Knepley ierr = PetscFEGetDualSpace(fe[f], &sp[f]);CHKERRQ(ierr); 20072f94c41SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr); 20172f94c41SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 20272f94c41SMatthew G. Knepley totDim += spDim*numComp; 203cb1e1211SMatthew G Knepley } 20472f94c41SMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 20572f94c41SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 20672f94c41SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, localX, cStart, &numValues, NULL);CHKERRQ(ierr); 20772f94c41SMatthew G. Knepley if (numValues != totDim) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "The section cell closure size %d != dual space dimension %d", numValues, totDim); 20872f94c41SMatthew G. Knepley ierr = DMGetWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 209dcca6d9dSJed Brown ierr = PetscMalloc2(dim,&v0,dim*dim,&J);CHKERRQ(ierr); 21072f94c41SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 21172f94c41SMatthew G. Knepley PetscCellGeometry geom; 212cb1e1211SMatthew G Knepley 213cb1e1211SMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr); 21472f94c41SMatthew G. Knepley geom.v0 = v0; 21572f94c41SMatthew G. Knepley geom.J = J; 21672f94c41SMatthew G. Knepley geom.detJ = &detJ; 21772f94c41SMatthew G. Knepley for (f = 0, v = 0; f < numFields; ++f) { 218c110b1eeSGeoffrey Irving void * const ctx = ctxs ? ctxs[f] : NULL; 21972f94c41SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr); 22072f94c41SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 22172f94c41SMatthew G. Knepley for (d = 0; d < spDim; ++d) { 222c110b1eeSGeoffrey Irving ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], ctx, &values[v]);CHKERRQ(ierr); 22372f94c41SMatthew G. Knepley v += numComp; 224cb1e1211SMatthew G Knepley } 225cb1e1211SMatthew G Knepley } 22672f94c41SMatthew G. Knepley ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr); 227cb1e1211SMatthew G Knepley } 22872f94c41SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 2291f2da991SMatthew G. Knepley ierr = PetscFree2(v0,J);CHKERRQ(ierr); 23072f94c41SMatthew G. Knepley ierr = PetscFree(sp);CHKERRQ(ierr); 231cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 232cb1e1211SMatthew G Knepley } 233cb1e1211SMatthew G Knepley 234cb1e1211SMatthew G Knepley #undef __FUNCT__ 235cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunction" 236cb1e1211SMatthew G Knepley /*@C 237cb1e1211SMatthew G Knepley DMPlexProjectFunction - This projects the given function into the function space provided. 238cb1e1211SMatthew G Knepley 239cb1e1211SMatthew G Knepley Input Parameters: 240cb1e1211SMatthew G Knepley + dm - The DM 24172f94c41SMatthew G. Knepley . fe - The PetscFE associated with the field 24272f94c41SMatthew G. Knepley . funcs - The coordinate functions to evaluate, one per field 243c110b1eeSGeoffrey Irving . ctxs - Optional array of contexts to pass to each coordinate function. ctxs itself may be null. 244cb1e1211SMatthew G Knepley - mode - The insertion mode for values 245cb1e1211SMatthew G Knepley 246cb1e1211SMatthew G Knepley Output Parameter: 247cb1e1211SMatthew G Knepley . X - vector 248cb1e1211SMatthew G Knepley 249cb1e1211SMatthew G Knepley Level: developer 250cb1e1211SMatthew G Knepley 251878cb397SSatish Balay .seealso: DMPlexComputeL2Diff() 252878cb397SSatish Balay @*/ 253c110b1eeSGeoffrey Irving PetscErrorCode DMPlexProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, InsertMode mode, Vec X) 254cb1e1211SMatthew G Knepley { 255cb1e1211SMatthew G Knepley Vec localX; 256cb1e1211SMatthew G Knepley PetscErrorCode ierr; 257cb1e1211SMatthew G Knepley 258cb1e1211SMatthew G Knepley PetscFunctionBegin; 2599a800dd8SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 260cb1e1211SMatthew G Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 261c110b1eeSGeoffrey Irving ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, mode, localX);CHKERRQ(ierr); 262cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr); 263cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr); 264cb1e1211SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 265cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 266cb1e1211SMatthew G Knepley } 267cb1e1211SMatthew G Knepley 268cb1e1211SMatthew G Knepley #undef __FUNCT__ 269cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff" 270cb1e1211SMatthew G Knepley /*@C 271cb1e1211SMatthew G Knepley DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 272cb1e1211SMatthew G Knepley 273cb1e1211SMatthew G Knepley Input Parameters: 274cb1e1211SMatthew G Knepley + dm - The DM 275c5bbbd5bSMatthew G. Knepley . fe - The PetscFE object for each field 276cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component 27751259fa3SMatthew G. Knepley . ctxs - Optional array of contexts to pass to each function, or NULL. 278cb1e1211SMatthew G Knepley - X - The coefficient vector u_h 279cb1e1211SMatthew G Knepley 280cb1e1211SMatthew G Knepley Output Parameter: 281cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2 282cb1e1211SMatthew G Knepley 283cb1e1211SMatthew G Knepley Level: developer 284cb1e1211SMatthew G Knepley 285cb1e1211SMatthew G Knepley .seealso: DMPlexProjectFunction() 286878cb397SSatish Balay @*/ 287c110b1eeSGeoffrey Irving PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, PetscReal *diff) 288cb1e1211SMatthew G Knepley { 289cb1e1211SMatthew G Knepley const PetscInt debug = 0; 290cb1e1211SMatthew G Knepley PetscSection section; 291c5bbbd5bSMatthew G. Knepley PetscQuadrature quad; 292cb1e1211SMatthew G Knepley Vec localX; 29372f94c41SMatthew G. Knepley PetscScalar *funcVal; 294cb1e1211SMatthew G Knepley PetscReal *coords, *v0, *J, *invJ, detJ; 295cb1e1211SMatthew G Knepley PetscReal localDiff = 0.0; 296cb1e1211SMatthew G Knepley PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 297cb1e1211SMatthew G Knepley PetscErrorCode ierr; 298cb1e1211SMatthew G Knepley 299cb1e1211SMatthew G Knepley PetscFunctionBegin; 300cb1e1211SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 301cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 302cb1e1211SMatthew G Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 303cb1e1211SMatthew G Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 304cb1e1211SMatthew G Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 305cb1e1211SMatthew G Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 306cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 307c5bbbd5bSMatthew G. Knepley PetscInt Nc; 308c5bbbd5bSMatthew G. Knepley 309c5bbbd5bSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 310c5bbbd5bSMatthew G. Knepley numComponents += Nc; 311cb1e1211SMatthew G Knepley } 312c110b1eeSGeoffrey Irving ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, ctxs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 313dcca6d9dSJed Brown ierr = PetscMalloc5(numComponents,&funcVal,dim,&coords,dim,&v0,dim*dim,&J,dim*dim,&invJ);CHKERRQ(ierr); 314cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 315c5bbbd5bSMatthew G. Knepley ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr); 316cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 317a1e44745SMatthew G. Knepley PetscScalar *x = NULL; 318cb1e1211SMatthew G Knepley PetscReal elemDiff = 0.0; 319cb1e1211SMatthew G Knepley 320cb1e1211SMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 321cb1e1211SMatthew G Knepley if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 322cb1e1211SMatthew G Knepley ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 323cb1e1211SMatthew G Knepley 324cb1e1211SMatthew G Knepley for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 325c110b1eeSGeoffrey Irving void * const ctx = ctxs ? ctxs[field] : NULL; 32621454ff5SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 327c5bbbd5bSMatthew G. Knepley PetscReal *basis; 32821454ff5SMatthew G. Knepley PetscInt numQuadPoints, numBasisFuncs, numBasisComps, q, d, e, fc, f; 329cb1e1211SMatthew G Knepley 33021454ff5SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 331c5bbbd5bSMatthew G. Knepley ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr); 332c5bbbd5bSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr); 333c5bbbd5bSMatthew G. Knepley ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr); 334cb1e1211SMatthew G Knepley if (debug) { 335cb1e1211SMatthew G Knepley char title[1024]; 336cb1e1211SMatthew G Knepley ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 337cb1e1211SMatthew G Knepley ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr); 338cb1e1211SMatthew G Knepley } 339cb1e1211SMatthew G Knepley for (q = 0; q < numQuadPoints; ++q) { 340cb1e1211SMatthew G Knepley for (d = 0; d < dim; d++) { 341cb1e1211SMatthew G Knepley coords[d] = v0[d]; 342cb1e1211SMatthew G Knepley for (e = 0; e < dim; e++) { 343cb1e1211SMatthew G Knepley coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 344cb1e1211SMatthew G Knepley } 345cb1e1211SMatthew G Knepley } 346c110b1eeSGeoffrey Irving (*funcs[field])(coords, funcVal, ctx); 347cb1e1211SMatthew G Knepley for (fc = 0; fc < numBasisComps; ++fc) { 348a1d24da5SMatthew G. Knepley PetscScalar interpolant = 0.0; 349a1d24da5SMatthew G. Knepley 350cb1e1211SMatthew G Knepley for (f = 0; f < numBasisFuncs; ++f) { 351cb1e1211SMatthew G Knepley const PetscInt fidx = f*numBasisComps+fc; 352a1d24da5SMatthew G. Knepley interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx]; 353cb1e1211SMatthew G Knepley } 35472f94c41SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);} 35572f94c41SMatthew G. Knepley elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 356cb1e1211SMatthew G Knepley } 357cb1e1211SMatthew G Knepley } 358cb1e1211SMatthew G Knepley comp += numBasisComps; 359cb1e1211SMatthew G Knepley fieldOffset += numBasisFuncs*numBasisComps; 360cb1e1211SMatthew G Knepley } 361cb1e1211SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 362cb1e1211SMatthew G Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 363cb1e1211SMatthew G Knepley localDiff += elemDiff; 364cb1e1211SMatthew G Knepley } 36572f94c41SMatthew G. Knepley ierr = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr); 366cb1e1211SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 36786a74ee0SMatthew G. Knepley ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 368cb1e1211SMatthew G Knepley *diff = PetscSqrtReal(*diff); 369cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 370cb1e1211SMatthew G Knepley } 371cb1e1211SMatthew G Knepley 372cb1e1211SMatthew G Knepley #undef __FUNCT__ 37340e14135SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeL2GradientDiff" 37440e14135SMatthew G. Knepley /*@C 37540e14135SMatthew G. Knepley DMPlexComputeL2GradientDiff - This function computes the L_2 difference between the gradient of a function u and an FEM interpolant solution grad u_h. 37640e14135SMatthew G. Knepley 37740e14135SMatthew G. Knepley Input Parameters: 37840e14135SMatthew G. Knepley + dm - The DM 37940e14135SMatthew G. Knepley . fe - The PetscFE object for each field 38040e14135SMatthew G. Knepley . funcs - The gradient functions to evaluate for each field component 38151259fa3SMatthew G. Knepley . ctxs - Optional array of contexts to pass to each function, or NULL. 38240e14135SMatthew G. Knepley . X - The coefficient vector u_h 38340e14135SMatthew G. Knepley - n - The vector to project along 38440e14135SMatthew G. Knepley 38540e14135SMatthew G. Knepley Output Parameter: 38640e14135SMatthew G. Knepley . diff - The diff ||(grad u - grad u_h) . n||_2 38740e14135SMatthew G. Knepley 38840e14135SMatthew G. Knepley Level: developer 38940e14135SMatthew G. Knepley 39040e14135SMatthew G. Knepley .seealso: DMPlexProjectFunction(), DMPlexComputeL2Diff() 39140e14135SMatthew G. Knepley @*/ 39251259fa3SMatthew G. Knepley PetscErrorCode DMPlexComputeL2GradientDiff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], const PetscReal [], PetscScalar *, void *), void **ctxs, Vec X, const PetscReal n[], PetscReal *diff) 393cb1e1211SMatthew G Knepley { 39440e14135SMatthew G. Knepley const PetscInt debug = 0; 395cb1e1211SMatthew G Knepley PetscSection section; 39640e14135SMatthew G. Knepley PetscQuadrature quad; 39740e14135SMatthew G. Knepley Vec localX; 39840e14135SMatthew G. Knepley PetscScalar *funcVal, *interpolantVec; 39940e14135SMatthew G. Knepley PetscReal *coords, *realSpaceDer, *v0, *J, *invJ, detJ; 40040e14135SMatthew G. Knepley PetscReal localDiff = 0.0; 40140e14135SMatthew G. Knepley PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 402cb1e1211SMatthew G Knepley PetscErrorCode ierr; 403cb1e1211SMatthew G Knepley 404cb1e1211SMatthew G Knepley PetscFunctionBegin; 40540e14135SMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 40640e14135SMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 40740e14135SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 40840e14135SMatthew G. Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 40940e14135SMatthew G. Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 41040e14135SMatthew G. Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 411652b88e8SMatthew G. Knepley for (field = 0; field < numFields; ++field) { 41240e14135SMatthew G. Knepley PetscInt Nc; 413652b88e8SMatthew G. Knepley 41440e14135SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 41540e14135SMatthew G. Knepley numComponents += Nc; 416652b88e8SMatthew G. Knepley } 41740e14135SMatthew G. Knepley /* ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); */ 41840e14135SMatthew G. Knepley ierr = PetscMalloc7(numComponents,&funcVal,dim,&coords,dim,&realSpaceDer,dim,&v0,dim*dim,&J,dim*dim,&invJ,dim,&interpolantVec);CHKERRQ(ierr); 41940e14135SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 42040e14135SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr); 42140e14135SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 42240e14135SMatthew G. Knepley PetscScalar *x = NULL; 42340e14135SMatthew G. Knepley PetscReal elemDiff = 0.0; 424652b88e8SMatthew G. Knepley 42540e14135SMatthew G. Knepley ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 42640e14135SMatthew G. Knepley if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 42740e14135SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 42840e14135SMatthew G. Knepley 42940e14135SMatthew G. Knepley for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 43051259fa3SMatthew G. Knepley void * const ctx = ctxs ? ctxs[field] : NULL; 43121454ff5SMatthew G. Knepley const PetscReal *quadPoints, *quadWeights; 43240e14135SMatthew G. Knepley PetscReal *basisDer; 43321454ff5SMatthew G. Knepley PetscInt numQuadPoints, Nb, Ncomp, q, d, e, fc, f, g; 43440e14135SMatthew G. Knepley 43521454ff5SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, &quadPoints, &quadWeights);CHKERRQ(ierr); 43640e14135SMatthew G. Knepley ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 43740e14135SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[field], &Ncomp);CHKERRQ(ierr); 43840e14135SMatthew G. Knepley ierr = PetscFEGetDefaultTabulation(fe[field], NULL, &basisDer, NULL);CHKERRQ(ierr); 43940e14135SMatthew G. Knepley if (debug) { 44040e14135SMatthew G. Knepley char title[1024]; 44140e14135SMatthew G. Knepley ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 44240e14135SMatthew G. Knepley ierr = DMPrintCellVector(c, title, Nb*Ncomp, &x[fieldOffset]);CHKERRQ(ierr); 443652b88e8SMatthew G. Knepley } 44440e14135SMatthew G. Knepley for (q = 0; q < numQuadPoints; ++q) { 44540e14135SMatthew G. Knepley for (d = 0; d < dim; d++) { 44640e14135SMatthew G. Knepley coords[d] = v0[d]; 44740e14135SMatthew G. Knepley for (e = 0; e < dim; e++) { 44840e14135SMatthew G. Knepley coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 449652b88e8SMatthew G. Knepley } 45040e14135SMatthew G. Knepley } 45151259fa3SMatthew G. Knepley (*funcs[field])(coords, n, funcVal, ctx); 45240e14135SMatthew G. Knepley for (fc = 0; fc < Ncomp; ++fc) { 45340e14135SMatthew G. Knepley PetscScalar interpolant = 0.0; 45440e14135SMatthew G. Knepley 45540e14135SMatthew G. Knepley for (d = 0; d < dim; ++d) interpolantVec[d] = 0.0; 45640e14135SMatthew G. Knepley for (f = 0; f < Nb; ++f) { 45740e14135SMatthew G. Knepley const PetscInt fidx = f*Ncomp+fc; 45840e14135SMatthew G. Knepley 45940e14135SMatthew G. Knepley for (d = 0; d < dim; ++d) { 46040e14135SMatthew G. Knepley realSpaceDer[d] = 0.0; 46140e14135SMatthew G. Knepley for (g = 0; g < dim; ++g) { 46240e14135SMatthew G. Knepley realSpaceDer[d] += invJ[g*dim+d]*basisDer[(q*Nb*Ncomp+fidx)*dim+g]; 46340e14135SMatthew G. Knepley } 46440e14135SMatthew G. Knepley interpolantVec[d] += x[fieldOffset+fidx]*realSpaceDer[d]; 46540e14135SMatthew G. Knepley } 46640e14135SMatthew G. Knepley } 46740e14135SMatthew G. Knepley for (d = 0; d < dim; ++d) interpolant += interpolantVec[d]*n[d]; 46840e14135SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d fieldDer %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ);CHKERRQ(ierr);} 46940e14135SMatthew G. Knepley elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 47040e14135SMatthew G. Knepley } 47140e14135SMatthew G. Knepley } 47240e14135SMatthew G. Knepley comp += Ncomp; 47340e14135SMatthew G. Knepley fieldOffset += Nb*Ncomp; 47440e14135SMatthew G. Knepley } 47540e14135SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 47640e14135SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 47740e14135SMatthew G. Knepley localDiff += elemDiff; 47840e14135SMatthew G. Knepley } 47940e14135SMatthew G. Knepley ierr = PetscFree7(funcVal,coords,realSpaceDer,v0,J,invJ,interpolantVec);CHKERRQ(ierr); 48040e14135SMatthew G. Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 48140e14135SMatthew G. Knepley ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 48240e14135SMatthew G. Knepley *diff = PetscSqrtReal(*diff); 483cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 484cb1e1211SMatthew G Knepley } 485cb1e1211SMatthew G Knepley 486a0845e3aSMatthew G. Knepley #undef __FUNCT__ 487a0845e3aSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeResidualFEM" 488a0845e3aSMatthew G. Knepley /*@ 489a0845e3aSMatthew G. Knepley DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user 490a0845e3aSMatthew G. Knepley 491a0845e3aSMatthew G. Knepley Input Parameters: 492a0845e3aSMatthew G. Knepley + dm - The mesh 493a0845e3aSMatthew G. Knepley . X - Local input vector 494a0845e3aSMatthew G. Knepley - user - The user context 495a0845e3aSMatthew G. Knepley 496a0845e3aSMatthew G. Knepley Output Parameter: 497a0845e3aSMatthew G. Knepley . F - Local output vector 498a0845e3aSMatthew G. Knepley 499a0845e3aSMatthew G. Knepley Note: 5008026896cSMatthew G. Knepley The first member of the user context must be an FEMContext. 501a0845e3aSMatthew G. Knepley 502a0845e3aSMatthew G. Knepley We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 503a0845e3aSMatthew G. Knepley like a GPU, or vectorize on a multicore machine. 504a0845e3aSMatthew G. Knepley 505a0845e3aSMatthew G. Knepley Level: developer 506a0845e3aSMatthew G. Knepley 507a0845e3aSMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM() 508a0845e3aSMatthew G. Knepley @*/ 509a0845e3aSMatthew G. Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 510a0845e3aSMatthew G. Knepley { 511a0845e3aSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 5129a559087SMatthew G. Knepley PetscFEM *fem = (PetscFEM *) user; 513a0845e3aSMatthew G. Knepley PetscFE *fe = fem->fe; 5149a559087SMatthew G. Knepley PetscFE *feAux = fem->feAux; 515f1ea0e2fSMatthew G. Knepley PetscFE *feBd = fem->feBd; 516a0845e3aSMatthew G. Knepley const char *name = "Residual"; 5179a559087SMatthew G. Knepley DM dmAux; 5189a559087SMatthew G. Knepley Vec A; 519a0845e3aSMatthew G. Knepley PetscQuadrature q; 520a0845e3aSMatthew G. Knepley PetscCellGeometry geom; 5219a559087SMatthew G. Knepley PetscSection section, sectionAux; 522a0845e3aSMatthew G. Knepley PetscReal *v0, *J, *invJ, *detJ; 52301599b20SMatthew G. Knepley PetscScalar *elemVec, *u, *a = NULL; 5249a559087SMatthew G. Knepley PetscInt dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c; 525a0845e3aSMatthew G. Knepley PetscInt cellDof = 0, numComponents = 0; 5269a559087SMatthew G. Knepley PetscInt cellDofAux = 0, numComponentsAux = 0; 527a0845e3aSMatthew G. Knepley PetscErrorCode ierr; 528a0845e3aSMatthew G. Knepley 529a0845e3aSMatthew G. Knepley PetscFunctionBegin; 530a0845e3aSMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 531a0845e3aSMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 532a0845e3aSMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 5339a559087SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 534a0845e3aSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 535a0845e3aSMatthew G. Knepley numCells = cEnd - cStart; 5369a559087SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 537a0845e3aSMatthew G. Knepley PetscInt Nb, Nc; 538a0845e3aSMatthew G. Knepley 539a0845e3aSMatthew G. Knepley ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 540a0845e3aSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 541a0845e3aSMatthew G. Knepley cellDof += Nb*Nc; 542a0845e3aSMatthew G. Knepley numComponents += Nc; 543a0845e3aSMatthew G. Knepley } 5449a559087SMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 5459a559087SMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 5469a559087SMatthew G. Knepley if (dmAux) { 5479a559087SMatthew G. Knepley ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 5489a559087SMatthew G. Knepley ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr); 5499a559087SMatthew G. Knepley } 5509a559087SMatthew G. Knepley for (f = 0; f < NfAux; ++f) { 5519a559087SMatthew G. Knepley PetscInt Nb, Nc; 5529a559087SMatthew G. Knepley 5539a559087SMatthew G. Knepley ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr); 5549a559087SMatthew G. Knepley ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr); 5559a559087SMatthew G. Knepley cellDofAux += Nb*Nc; 5569a559087SMatthew G. Knepley numComponentsAux += Nc; 5579a559087SMatthew G. Knepley } 558c110b1eeSGeoffrey Irving ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, fem->bcCtxs, INSERT_BC_VALUES, X);CHKERRQ(ierr); 559a0845e3aSMatthew G. Knepley ierr = VecSet(F, 0.0);CHKERRQ(ierr); 560dcca6d9dSJed Brown ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr); 561785e854fSJed Brown if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);} 562a0845e3aSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 563a0845e3aSMatthew G. Knepley PetscScalar *x = NULL; 564a0845e3aSMatthew G. Knepley PetscInt i; 565a0845e3aSMatthew G. Knepley 566a0845e3aSMatthew G. Knepley ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 567a0845e3aSMatthew G. Knepley if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 568a0845e3aSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 569a0845e3aSMatthew G. Knepley for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 570a0845e3aSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 5719a559087SMatthew G. Knepley if (dmAux) { 5729a559087SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 5739a559087SMatthew G. Knepley for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i]; 5749a559087SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 575a0845e3aSMatthew G. Knepley } 5769a559087SMatthew G. Knepley } 5779a559087SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 578c012ea0aSMatthew G. Knepley void (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0Funcs[f]; 579c012ea0aSMatthew G. Knepley void (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1Funcs[f]; 58021454ff5SMatthew G. Knepley PetscInt numQuadPoints, Nb; 581a0845e3aSMatthew G. Knepley /* Conforming batches */ 582f30c5766SMatthew G. Knepley PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 583a0845e3aSMatthew G. Knepley /* Remainder */ 584a0845e3aSMatthew G. Knepley PetscInt Nr, offset; 585a0845e3aSMatthew G. Knepley 586a0845e3aSMatthew G. Knepley ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr); 587a0845e3aSMatthew G. Knepley ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 588f30c5766SMatthew G. Knepley ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 58921454ff5SMatthew G. Knepley ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 59021454ff5SMatthew G. Knepley blockSize = Nb*numQuadPoints; 591a0845e3aSMatthew G. Knepley batchSize = numBlocks * blockSize; 592f30c5766SMatthew G. Knepley ierr = PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 593a0845e3aSMatthew G. Knepley numChunks = numCells / (numBatches*batchSize); 594a0845e3aSMatthew G. Knepley Ne = numChunks*numBatches*batchSize; 595a0845e3aSMatthew G. Knepley Nr = numCells % (numBatches*batchSize); 596a0845e3aSMatthew G. Knepley offset = numCells - Nr; 597a0845e3aSMatthew G. Knepley geom.v0 = v0; 598a0845e3aSMatthew G. Knepley geom.J = J; 599a0845e3aSMatthew G. Knepley geom.invJ = invJ; 600a0845e3aSMatthew G. Knepley geom.detJ = detJ; 6019a559087SMatthew G. Knepley ierr = PetscFEIntegrateResidual(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr); 602a0845e3aSMatthew G. Knepley geom.v0 = &v0[offset*dim]; 603a0845e3aSMatthew G. Knepley geom.J = &J[offset*dim*dim]; 604a0845e3aSMatthew G. Knepley geom.invJ = &invJ[offset*dim*dim]; 605a0845e3aSMatthew G. Knepley geom.detJ = &detJ[offset]; 6069a559087SMatthew G. Knepley ierr = PetscFEIntegrateResidual(fe[f], Nr, Nf, fe, f, geom, &u[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 607a0845e3aSMatthew G. Knepley } 608a0845e3aSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 609a0845e3aSMatthew G. Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 610a0845e3aSMatthew G. Knepley ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 611a0845e3aSMatthew G. Knepley } 612a0845e3aSMatthew G. Knepley ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 6139a559087SMatthew G. Knepley if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 614f1ea0e2fSMatthew G. Knepley if (feBd) { 615075da914SMatthew G. Knepley DMLabel label, depth; 616f1ea0e2fSMatthew G. Knepley IS pointIS; 617f1ea0e2fSMatthew G. Knepley const PetscInt *points; 618075da914SMatthew G. Knepley PetscInt dep, numPoints, p, numFaces; 619f1ea0e2fSMatthew G. Knepley PetscReal *n; 620f1ea0e2fSMatthew G. Knepley 621f1ea0e2fSMatthew G. Knepley ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr); 622075da914SMatthew G. Knepley ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr); 623f1ea0e2fSMatthew G. Knepley ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr); 624f1ea0e2fSMatthew G. Knepley ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr); 625f1ea0e2fSMatthew G. Knepley ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 626f1ea0e2fSMatthew G. Knepley for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) { 627f1ea0e2fSMatthew G. Knepley PetscInt Nb, Nc; 628f1ea0e2fSMatthew G. Knepley 629f1ea0e2fSMatthew G. Knepley ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr); 630f1ea0e2fSMatthew G. Knepley ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr); 631f1ea0e2fSMatthew G. Knepley cellDof += Nb*Nc; 632f1ea0e2fSMatthew G. Knepley numComponents += Nc; 633f1ea0e2fSMatthew G. Knepley } 634075da914SMatthew G. Knepley for (p = 0, numFaces = 0; p < numPoints; ++p) { 635075da914SMatthew G. Knepley ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 636075da914SMatthew G. Knepley if (dep == dim-1) ++numFaces; 637075da914SMatthew G. Knepley } 638dcca6d9dSJed Brown ierr = PetscMalloc7(numFaces*cellDof,&u,numFaces*dim,&v0,numFaces*dim,&n,numFaces*dim*dim,&J,numFaces*dim*dim,&invJ,numFaces,&detJ,numFaces*cellDof,&elemVec);CHKERRQ(ierr); 639075da914SMatthew G. Knepley for (p = 0, f = 0; p < numPoints; ++p) { 640f1ea0e2fSMatthew G. Knepley const PetscInt point = points[p]; 641f1ea0e2fSMatthew G. Knepley PetscScalar *x = NULL; 642f1ea0e2fSMatthew G. Knepley PetscInt i; 643f1ea0e2fSMatthew G. Knepley 644075da914SMatthew G. Knepley ierr = DMLabelGetValue(depth, points[p], &dep);CHKERRQ(ierr); 645075da914SMatthew G. Knepley if (dep != dim-1) continue; 646075da914SMatthew G. Knepley ierr = DMPlexComputeCellGeometry(dm, point, &v0[f*dim], &J[f*dim*dim], &invJ[f*dim*dim], &detJ[f]);CHKERRQ(ierr); 647a8007bbfSMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dm, point, NULL, NULL, &n[f*dim]); 648075da914SMatthew G. Knepley if (detJ[f] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[f], point); 649f1ea0e2fSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 650075da914SMatthew G. Knepley for (i = 0; i < cellDof; ++i) u[f*cellDof+i] = x[i]; 651f1ea0e2fSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 652075da914SMatthew G. Knepley ++f; 653f1ea0e2fSMatthew G. Knepley } 654f1ea0e2fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 655f1ea0e2fSMatthew G. Knepley void (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[f]; 656f1ea0e2fSMatthew G. Knepley void (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[f]; 65721454ff5SMatthew G. Knepley PetscInt numQuadPoints, Nb; 658f1ea0e2fSMatthew G. Knepley /* Conforming batches */ 659f1ea0e2fSMatthew G. Knepley PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 660f1ea0e2fSMatthew G. Knepley /* Remainder */ 661f1ea0e2fSMatthew G. Knepley PetscInt Nr, offset; 662f1ea0e2fSMatthew G. Knepley 663f1ea0e2fSMatthew G. Knepley ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr); 664f1ea0e2fSMatthew G. Knepley ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr); 665f1ea0e2fSMatthew G. Knepley ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 66621454ff5SMatthew G. Knepley ierr = PetscQuadratureGetData(q, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 66721454ff5SMatthew G. Knepley blockSize = Nb*numQuadPoints; 668f1ea0e2fSMatthew G. Knepley batchSize = numBlocks * blockSize; 669f1ea0e2fSMatthew G. Knepley ierr = PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 670075da914SMatthew G. Knepley numChunks = numFaces / (numBatches*batchSize); 671f1ea0e2fSMatthew G. Knepley Ne = numChunks*numBatches*batchSize; 672075da914SMatthew G. Knepley Nr = numFaces % (numBatches*batchSize); 673075da914SMatthew G. Knepley offset = numFaces - Nr; 674f1ea0e2fSMatthew G. Knepley geom.v0 = v0; 675f1ea0e2fSMatthew G. Knepley geom.n = n; 676f1ea0e2fSMatthew G. Knepley geom.J = J; 677f1ea0e2fSMatthew G. Knepley geom.invJ = invJ; 678f1ea0e2fSMatthew G. Knepley geom.detJ = detJ; 679f1ea0e2fSMatthew G. Knepley ierr = PetscFEIntegrateBdResidual(feBd[f], Ne, Nf, feBd, f, geom, u, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr); 680f1ea0e2fSMatthew G. Knepley geom.v0 = &v0[offset*dim]; 681f1ea0e2fSMatthew G. Knepley geom.n = &n[offset*dim]; 682f1ea0e2fSMatthew G. Knepley geom.J = &J[offset*dim*dim]; 683f1ea0e2fSMatthew G. Knepley geom.invJ = &invJ[offset*dim*dim]; 684f1ea0e2fSMatthew G. Knepley geom.detJ = &detJ[offset]; 685f1ea0e2fSMatthew G. Knepley ierr = PetscFEIntegrateBdResidual(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 686f1ea0e2fSMatthew G. Knepley } 687075da914SMatthew G. Knepley for (p = 0, f = 0; p < numPoints; ++p) { 688f1ea0e2fSMatthew G. Knepley const PetscInt point = points[p]; 689f1ea0e2fSMatthew G. Knepley 690075da914SMatthew G. Knepley ierr = DMLabelGetValue(depth, point, &dep);CHKERRQ(ierr); 691075da914SMatthew G. Knepley if (dep != dim-1) continue; 692075da914SMatthew G. Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[f*cellDof]);CHKERRQ(ierr);} 693075da914SMatthew G. Knepley ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[f*cellDof], ADD_VALUES);CHKERRQ(ierr); 694075da914SMatthew G. Knepley ++f; 695f1ea0e2fSMatthew G. Knepley } 696f1ea0e2fSMatthew G. Knepley ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 697f1ea0e2fSMatthew G. Knepley ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 698f1ea0e2fSMatthew G. Knepley ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr); 699f1ea0e2fSMatthew G. Knepley } 7006113b454SMatthew G. Knepley if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);} 701a0845e3aSMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 702a0845e3aSMatthew G. Knepley PetscFunctionReturn(0); 703a0845e3aSMatthew G. Knepley } 704a0845e3aSMatthew G. Knepley 705cb1e1211SMatthew G Knepley #undef __FUNCT__ 706cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianActionFEM" 707cb1e1211SMatthew G Knepley /*@C 708cb1e1211SMatthew G Knepley DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user 709cb1e1211SMatthew G Knepley 710cb1e1211SMatthew G Knepley Input Parameters: 711cb1e1211SMatthew G Knepley + dm - The mesh 712cb1e1211SMatthew G Knepley . J - The Jacobian shell matrix 713cb1e1211SMatthew G Knepley . X - Local input vector 714cb1e1211SMatthew G Knepley - user - The user context 715cb1e1211SMatthew G Knepley 716cb1e1211SMatthew G Knepley Output Parameter: 717cb1e1211SMatthew G Knepley . F - Local output vector 718cb1e1211SMatthew G Knepley 719cb1e1211SMatthew G Knepley Note: 7208026896cSMatthew G. Knepley The first member of the user context must be an FEMContext. 721cb1e1211SMatthew G Knepley 722cb1e1211SMatthew G Knepley We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 723cb1e1211SMatthew G Knepley like a GPU, or vectorize on a multicore machine. 724cb1e1211SMatthew G Knepley 7250059ad2aSSatish Balay Level: developer 7260059ad2aSSatish Balay 727cb1e1211SMatthew G Knepley .seealso: DMPlexComputeResidualFEM() 728878cb397SSatish Balay @*/ 729cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user) 730cb1e1211SMatthew G Knepley { 731cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 7329a559087SMatthew G. Knepley PetscFEM *fem = (PetscFEM *) user; 7330483ade4SMatthew G. Knepley PetscFE *fe = fem->fe; 7340483ade4SMatthew G. Knepley PetscQuadrature quad; 7350483ade4SMatthew G. Knepley PetscCellGeometry geom; 736cb1e1211SMatthew G Knepley PetscSection section; 737cb1e1211SMatthew G Knepley JacActionCtx *jctx; 738cb1e1211SMatthew G Knepley PetscReal *v0, *J, *invJ, *detJ; 739cb1e1211SMatthew G Knepley PetscScalar *elemVec, *u, *a; 7400483ade4SMatthew G. Knepley PetscInt dim, numFields, field, numCells, cStart, cEnd, c; 741cb1e1211SMatthew G Knepley PetscInt cellDof = 0; 742cb1e1211SMatthew G Knepley PetscErrorCode ierr; 743cb1e1211SMatthew G Knepley 744cb1e1211SMatthew G Knepley PetscFunctionBegin; 7450483ade4SMatthew G. Knepley /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */ 746cb1e1211SMatthew G Knepley ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 747cb1e1211SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 748cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 749cb1e1211SMatthew G Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 750cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 751cb1e1211SMatthew G Knepley numCells = cEnd - cStart; 752cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 7530483ade4SMatthew G. Knepley PetscInt Nb, Nc; 7540483ade4SMatthew G. Knepley 7550483ade4SMatthew G. Knepley ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 7560483ade4SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 7570483ade4SMatthew G. Knepley cellDof += Nb*Nc; 758cb1e1211SMatthew G Knepley } 759cb1e1211SMatthew G Knepley ierr = VecSet(F, 0.0);CHKERRQ(ierr); 760dcca6d9dSJed Brown ierr = PetscMalloc7(numCells*cellDof,&u,numCells*cellDof,&a,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof,&elemVec);CHKERRQ(ierr); 761cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 762a1e44745SMatthew G. Knepley PetscScalar *x = NULL; 763cb1e1211SMatthew G Knepley PetscInt i; 764cb1e1211SMatthew G Knepley 765cb1e1211SMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 766cb1e1211SMatthew G Knepley if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 767cb1e1211SMatthew G Knepley ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 768cb1e1211SMatthew G Knepley for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 769cb1e1211SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 770cb1e1211SMatthew G Knepley ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 771cb1e1211SMatthew G Knepley for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i]; 772cb1e1211SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 773cb1e1211SMatthew G Knepley } 774cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 77521454ff5SMatthew G. Knepley PetscInt numQuadPoints, Nb; 776cb1e1211SMatthew G Knepley /* Conforming batches */ 777cb1e1211SMatthew G Knepley PetscInt numBlocks = 1; 7780483ade4SMatthew G. Knepley PetscInt numBatches = 1; 7790483ade4SMatthew G. Knepley PetscInt numChunks, Ne, blockSize, batchSize; 780cb1e1211SMatthew G Knepley /* Remainder */ 7810483ade4SMatthew G. Knepley PetscInt Nr, offset; 782cb1e1211SMatthew G Knepley 7830483ade4SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr); 7840483ade4SMatthew G. Knepley ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 78521454ff5SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 78621454ff5SMatthew G. Knepley blockSize = Nb*numQuadPoints; 7870483ade4SMatthew G. Knepley batchSize = numBlocks * blockSize; 7880483ade4SMatthew G. Knepley numChunks = numCells / (numBatches*batchSize); 7890483ade4SMatthew G. Knepley Ne = numChunks*numBatches*batchSize; 7900483ade4SMatthew G. Knepley Nr = numCells % (numBatches*batchSize); 7910483ade4SMatthew G. Knepley offset = numCells - Nr; 7920483ade4SMatthew G. Knepley geom.v0 = v0; 7930483ade4SMatthew G. Knepley geom.J = J; 7940483ade4SMatthew G. Knepley geom.invJ = invJ; 7950483ade4SMatthew G. Knepley geom.detJ = detJ; 7960483ade4SMatthew G. Knepley ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr); 7970483ade4SMatthew G. Knepley geom.v0 = &v0[offset*dim]; 7980483ade4SMatthew G. Knepley geom.J = &J[offset*dim*dim]; 7990483ade4SMatthew G. Knepley geom.invJ = &invJ[offset*dim*dim]; 8000483ade4SMatthew G. Knepley geom.detJ = &detJ[offset]; 8010483ade4SMatthew G. Knepley ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof], 802cb1e1211SMatthew G Knepley fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr); 803cb1e1211SMatthew G Knepley } 804cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 805cb1e1211SMatthew G Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 806cb1e1211SMatthew G Knepley ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 807cb1e1211SMatthew G Knepley } 808cb1e1211SMatthew G Knepley ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 809cb1e1211SMatthew G Knepley if (mesh->printFEM) { 810cb1e1211SMatthew G Knepley PetscMPIInt rank, numProcs; 811cb1e1211SMatthew G Knepley PetscInt p; 812cb1e1211SMatthew G Knepley 813cb1e1211SMatthew G Knepley ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 814cb1e1211SMatthew G Knepley ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr); 81586a74ee0SMatthew G. Knepley ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr); 816cb1e1211SMatthew G Knepley for (p = 0; p < numProcs; ++p) { 817cb1e1211SMatthew G Knepley if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);} 818cb1e1211SMatthew G Knepley ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr); 819cb1e1211SMatthew G Knepley } 820cb1e1211SMatthew G Knepley } 8210483ade4SMatthew G. Knepley /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */ 822cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 823cb1e1211SMatthew G Knepley } 824cb1e1211SMatthew G Knepley 825cb1e1211SMatthew G Knepley #undef __FUNCT__ 826cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianFEM" 827cb1e1211SMatthew G Knepley /*@ 828cb1e1211SMatthew G Knepley DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user. 829cb1e1211SMatthew G Knepley 830cb1e1211SMatthew G Knepley Input Parameters: 831cb1e1211SMatthew G Knepley + dm - The mesh 832cb1e1211SMatthew G Knepley . X - Local input vector 833cb1e1211SMatthew G Knepley - user - The user context 834cb1e1211SMatthew G Knepley 835cb1e1211SMatthew G Knepley Output Parameter: 836cb1e1211SMatthew G Knepley . Jac - Jacobian matrix 837cb1e1211SMatthew G Knepley 838cb1e1211SMatthew G Knepley Note: 8398026896cSMatthew G. Knepley The first member of the user context must be an FEMContext. 840cb1e1211SMatthew G Knepley 841cb1e1211SMatthew G Knepley We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 842cb1e1211SMatthew G Knepley like a GPU, or vectorize on a multicore machine. 843cb1e1211SMatthew G Knepley 8440059ad2aSSatish Balay Level: developer 8450059ad2aSSatish Balay 846cb1e1211SMatthew G Knepley .seealso: FormFunctionLocal() 847878cb397SSatish Balay @*/ 8480405ed22SMatthew G. Knepley PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP,void *user) 849cb1e1211SMatthew G Knepley { 850cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 8519a559087SMatthew G. Knepley PetscFEM *fem = (PetscFEM *) user; 852a319912fSMatthew G. Knepley PetscFE *fe = fem->fe; 853754551f4SMatthew G. Knepley PetscFE *feAux = fem->feAux; 854a319912fSMatthew G. Knepley const char *name = "Jacobian"; 855754551f4SMatthew G. Knepley DM dmAux; 856754551f4SMatthew G. Knepley Vec A; 857a319912fSMatthew G. Knepley PetscQuadrature quad; 858a319912fSMatthew G. Knepley PetscCellGeometry geom; 859754551f4SMatthew G. Knepley PetscSection section, globalSection, sectionAux; 860cb1e1211SMatthew G Knepley PetscReal *v0, *J, *invJ, *detJ; 861754551f4SMatthew G. Knepley PetscScalar *elemMat, *u, *a; 862754551f4SMatthew G. Knepley PetscInt dim, Nf, NfAux = 0, f, fieldI, fieldJ, numCells, cStart, cEnd, c; 863cb1e1211SMatthew G Knepley PetscInt cellDof = 0, numComponents = 0; 864754551f4SMatthew G. Knepley PetscInt cellDofAux = 0, numComponentsAux = 0; 865cb1e1211SMatthew G Knepley PetscBool isShell; 866cb1e1211SMatthew G Knepley PetscErrorCode ierr; 867cb1e1211SMatthew G Knepley 868cb1e1211SMatthew G Knepley PetscFunctionBegin; 869a319912fSMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 870cb1e1211SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 871cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 872a319912fSMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr); 873754551f4SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 874cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 875cb1e1211SMatthew G Knepley numCells = cEnd - cStart; 876754551f4SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 877a319912fSMatthew G. Knepley PetscInt Nb, Nc; 878a319912fSMatthew G. Knepley 879a319912fSMatthew G. Knepley ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 880a319912fSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 881a319912fSMatthew G. Knepley cellDof += Nb*Nc; 882a319912fSMatthew G. Knepley numComponents += Nc; 883cb1e1211SMatthew G Knepley } 884754551f4SMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 885754551f4SMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 886754551f4SMatthew G. Knepley if (dmAux) { 887754551f4SMatthew G. Knepley ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 888754551f4SMatthew G. Knepley ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr); 889754551f4SMatthew G. Knepley } 890754551f4SMatthew G. Knepley for (f = 0; f < NfAux; ++f) { 891754551f4SMatthew G. Knepley PetscInt Nb, Nc; 892754551f4SMatthew G. Knepley 893754551f4SMatthew G. Knepley ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr); 894754551f4SMatthew G. Knepley ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr); 895754551f4SMatthew G. Knepley cellDofAux += Nb*Nc; 896754551f4SMatthew G. Knepley numComponentsAux += Nc; 897754551f4SMatthew G. Knepley } 898c110b1eeSGeoffrey Irving ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, fem->bcCtxs, INSERT_BC_VALUES, X);CHKERRQ(ierr); 899cb1e1211SMatthew G Knepley ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 900dcca6d9dSJed Brown ierr = PetscMalloc6(numCells*cellDof,&u,numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*cellDof*cellDof,&elemMat);CHKERRQ(ierr); 901785e854fSJed Brown if (dmAux) {ierr = PetscMalloc1(numCells*cellDofAux, &a);CHKERRQ(ierr);} 902cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 903a1e44745SMatthew G. Knepley PetscScalar *x = NULL; 904cb1e1211SMatthew G Knepley PetscInt i; 905cb1e1211SMatthew G Knepley 906cb1e1211SMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 907cb1e1211SMatthew G Knepley if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 908a319912fSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 909cb1e1211SMatthew G Knepley for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 910a319912fSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 911754551f4SMatthew G. Knepley if (dmAux) { 912754551f4SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 913754551f4SMatthew G. Knepley for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i]; 914754551f4SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 915754551f4SMatthew G. Knepley } 916cb1e1211SMatthew G Knepley } 917cb1e1211SMatthew G Knepley ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr); 918754551f4SMatthew G. Knepley for (fieldI = 0; fieldI < Nf; ++fieldI) { 91921454ff5SMatthew G. Knepley PetscInt numQuadPoints, Nb; 920cb1e1211SMatthew G Knepley /* Conforming batches */ 921754551f4SMatthew G. Knepley PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 922cb1e1211SMatthew G Knepley /* Remainder */ 923a319912fSMatthew G. Knepley PetscInt Nr, offset; 924cb1e1211SMatthew G Knepley 925754551f4SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr); 926754551f4SMatthew G. Knepley ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr); 927754551f4SMatthew G. Knepley ierr = PetscFEGetTileSizes(fe[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 92821454ff5SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 92921454ff5SMatthew G. Knepley blockSize = Nb*numQuadPoints; 930a319912fSMatthew G. Knepley batchSize = numBlocks * blockSize; 931754551f4SMatthew G. Knepley ierr = PetscFESetTileSizes(fe[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 932a319912fSMatthew G. Knepley numChunks = numCells / (numBatches*batchSize); 933a319912fSMatthew G. Knepley Ne = numChunks*numBatches*batchSize; 934a319912fSMatthew G. Knepley Nr = numCells % (numBatches*batchSize); 935a319912fSMatthew G. Knepley offset = numCells - Nr; 936754551f4SMatthew G. Knepley for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 937754551f4SMatthew G. Knepley void (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ]; 938754551f4SMatthew G. Knepley void (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*Nf+fieldJ]; 939754551f4SMatthew G. Knepley void (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*Nf+fieldJ]; 940754551f4SMatthew G. Knepley void (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*Nf+fieldJ]; 941754551f4SMatthew G. Knepley 942a319912fSMatthew G. Knepley geom.v0 = v0; 943a319912fSMatthew G. Knepley geom.J = J; 944a319912fSMatthew G. Knepley geom.invJ = invJ; 945a319912fSMatthew G. Knepley geom.detJ = detJ; 946754551f4SMatthew G. Knepley ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, Nf, fe, fieldI, fieldJ, geom, u, NfAux, feAux, a, g0, g1, g2, g3, elemMat);CHKERRQ(ierr); 947a319912fSMatthew G. Knepley geom.v0 = &v0[offset*dim]; 948a319912fSMatthew G. Knepley geom.J = &J[offset*dim*dim]; 949a319912fSMatthew G. Knepley geom.invJ = &invJ[offset*dim*dim]; 950a319912fSMatthew G. Knepley geom.detJ = &detJ[offset]; 951754551f4SMatthew G. Knepley ierr = PetscFEIntegrateJacobian(fe[fieldI], Nr, Nf, fe, fieldI, fieldJ, geom, &u[offset*cellDof], NfAux, feAux, &a[offset*cellDofAux], g0, g1, g2, g3, &elemMat[offset*cellDof*cellDof]);CHKERRQ(ierr); 952cb1e1211SMatthew G Knepley } 953cb1e1211SMatthew G Knepley } 954cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 955a319912fSMatthew G. Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);} 956a319912fSMatthew G. Knepley ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr); 957cb1e1211SMatthew G Knepley } 958cb1e1211SMatthew G Knepley ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr); 959754551f4SMatthew G. Knepley if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 960cb1e1211SMatthew G Knepley ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 961cb1e1211SMatthew G Knepley ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 962cb1e1211SMatthew G Knepley if (mesh->printFEM) { 963a319912fSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 964cb1e1211SMatthew G Knepley ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr); 965cb1e1211SMatthew G Knepley ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 966cb1e1211SMatthew G Knepley } 967a319912fSMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 968cb1e1211SMatthew G Knepley ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr); 969cb1e1211SMatthew G Knepley if (isShell) { 970cb1e1211SMatthew G Knepley JacActionCtx *jctx; 971cb1e1211SMatthew G Knepley 972cb1e1211SMatthew G Knepley ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 973cb1e1211SMatthew G Knepley ierr = VecCopy(X, jctx->u);CHKERRQ(ierr); 974cb1e1211SMatthew G Knepley } 975cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 976cb1e1211SMatthew G Knepley } 977bceba477SMatthew G. Knepley 978d69c5d34SMatthew G. Knepley #if 0 979d69c5d34SMatthew G. Knepley 98032356553SMatthew G. Knepley static void g0_identity_1d_static(const PetscScalar u[], const PetscScalar gradU[], const PetscScalar a[], const PetscScalar gradA[], const PetscReal x[], PetscScalar g3[]) 98132356553SMatthew G. Knepley { 98232356553SMatthew G. Knepley g3[0] = 1.0; 98332356553SMatthew G. Knepley } 98432356553SMatthew G. Knepley 98532356553SMatthew G. Knepley static void g0_identity_2d_static(const PetscScalar u[], const PetscScalar gradU[], const PetscScalar a[], const PetscScalar gradA[], const PetscReal x[], PetscScalar g3[]) 98632356553SMatthew G. Knepley { 98732356553SMatthew G. Knepley g3[0] = g3[3] = 1.0; 98832356553SMatthew G. Knepley } 98932356553SMatthew G. Knepley 99032356553SMatthew G. Knepley static void g0_identity_3d_static(const PetscScalar u[], const PetscScalar gradU[], const PetscScalar a[], const PetscScalar gradA[], const PetscReal x[], PetscScalar g3[]) 99132356553SMatthew G. Knepley { 99232356553SMatthew G. Knepley g3[0] = g3[4] = g3[8] = 1.0; 99332356553SMatthew G. Knepley } 99432356553SMatthew G. Knepley 995bceba477SMatthew G. Knepley #undef __FUNCT__ 996d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEMBroken" 997d69c5d34SMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEMBroken(DM dmc, DM dmf, Mat I, void *user) 998bceba477SMatthew G. Knepley { 999bceba477SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dmc->data; 1000bceba477SMatthew G. Knepley PetscFEM *fem = (PetscFEM *) user; 1001bceba477SMatthew G. Knepley PetscFE *fe = fem->fe; 1002bceba477SMatthew G. Knepley const char *name = "Interpolator"; 100332356553SMatthew G. Knepley PetscFE *feRef; 1004bceba477SMatthew G. Knepley PetscQuadrature quad, quadOld; 1005bceba477SMatthew G. Knepley PetscCellGeometry geom; 100632356553SMatthew G. Knepley PetscSection fsection, fglobalSection; 100732356553SMatthew G. Knepley PetscSection csection, cglobalSection; 1008bceba477SMatthew G. Knepley PetscReal *v0, *J, *invJ, *detJ; 1009bceba477SMatthew G. Knepley PetscScalar *elemMat; 1010bceba477SMatthew G. Knepley PetscInt dim, Nf, f, fieldI, fieldJ, numCells, cStart, cEnd, c; 101132356553SMatthew G. Knepley PetscInt rCellDof = 0, cCellDof = 0, numComponents = 0; 1012bceba477SMatthew G. Knepley PetscErrorCode ierr; 1013bceba477SMatthew G. Knepley 1014bceba477SMatthew G. Knepley PetscFunctionBegin; 1015bceba477SMatthew G. Knepley #if 0 1016bceba477SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1017bceba477SMatthew G. Knepley #endif 1018bceba477SMatthew G. Knepley ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr); 101932356553SMatthew G. Knepley ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 102032356553SMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 102132356553SMatthew G. Knepley ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 102232356553SMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 102332356553SMatthew G. Knepley ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 1024bceba477SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 1025bceba477SMatthew G. Knepley numCells = cEnd - cStart; 102632356553SMatthew G. Knepley ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr); 102732356553SMatthew G. Knepley for (fieldI = 0; fieldI < Nf; ++fieldI) { 102832356553SMatthew G. Knepley ierr = PetscFERefine(fe[fieldI], &feRef[fieldI]);CHKERRQ(ierr); 102932356553SMatthew G. Knepley } 1030bceba477SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 103132356553SMatthew G. Knepley PetscInt rNb, cNb, Nc; 1032bceba477SMatthew G. Knepley 103332356553SMatthew G. Knepley ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr); 103432356553SMatthew G. Knepley ierr = PetscFEGetDimension(fe[f], &cNb);CHKERRQ(ierr); 1035bceba477SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 1036bceba477SMatthew G. Knepley numComponents += Nc; 103732356553SMatthew G. Knepley rCellDof += rNb*Nc; 103832356553SMatthew G. Knepley cCellDof += cNb*Nc; 1039bceba477SMatthew G. Knepley } 1040bceba477SMatthew G. Knepley ierr = MatZeroEntries(I);CHKERRQ(ierr); 104132356553SMatthew G. Knepley ierr = PetscMalloc5(numCells*dim,&v0,numCells*dim*dim,&J,numCells*dim*dim,&invJ,numCells,&detJ,numCells*rCellDof*cCellDof,&elemMat);CHKERRQ(ierr); 1042bceba477SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 1043bceba477SMatthew G. Knepley ierr = DMPlexComputeCellGeometry(dmc, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 1044bceba477SMatthew G. Knepley if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 1045bceba477SMatthew G. Knepley } 104632356553SMatthew G. Knepley ierr = PetscMemzero(elemMat, numCells*rCellDof*cCellDof * sizeof(PetscScalar));CHKERRQ(ierr); 1047bceba477SMatthew G. Knepley for (fieldI = 0; fieldI < Nf; ++fieldI) { 1048bceba477SMatthew G. Knepley PetscInt numQuadPoints, Nb; 1049bceba477SMatthew G. Knepley /* Conforming batches */ 1050bceba477SMatthew G. Knepley PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 1051bceba477SMatthew G. Knepley /* Remainder */ 1052bceba477SMatthew G. Knepley PetscInt Nr, offset; 1053bceba477SMatthew G. Knepley 1054bceba477SMatthew G. Knepley /* Make new fine FE which refines the ref cell and the quadrature rule */ 105532356553SMatthew G. Knepley ierr = PetscFEGetQuadrature(feRef[fieldI], &quad);CHKERRQ(ierr); 105632356553SMatthew G. Knepley ierr = PetscFEGetDimension(feRef[fieldI], &Nb);CHKERRQ(ierr); 105732356553SMatthew G. Knepley ierr = PetscFEGetTileSizes(feRef[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 1058bceba477SMatthew G. Knepley ierr = PetscQuadratureGetData(quad, NULL, &numQuadPoints, NULL, NULL);CHKERRQ(ierr); 1059bceba477SMatthew G. Knepley blockSize = Nb*numQuadPoints; 1060bceba477SMatthew G. Knepley batchSize = numBlocks * blockSize; 106132356553SMatthew G. Knepley ierr = PetscFESetTileSizes(feRef[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 1062bceba477SMatthew G. Knepley numChunks = numCells / (numBatches*batchSize); 1063bceba477SMatthew G. Knepley Ne = numChunks*numBatches*batchSize; 1064bceba477SMatthew G. Knepley Nr = numCells % (numBatches*batchSize); 1065bceba477SMatthew G. Knepley offset = numCells - Nr; 1066d69c5d34SMatthew G. Knepley 1067bceba477SMatthew G. Knepley for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 106832356553SMatthew G. Knepley /* void (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ]; */ 106932356553SMatthew G. Knepley void (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = g0_identity_2d_static; 1070bceba477SMatthew G. Knepley 1071bceba477SMatthew G. Knepley /* Replace quadrature in coarse FE with refined quadrature */ 1072bceba477SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe[fieldJ], &quadOld);CHKERRQ(ierr); 1073bceba477SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) quadOld);CHKERRQ(ierr); 1074bceba477SMatthew G. Knepley ierr = PetscFESetQuadrature(fe[fieldJ], quad);CHKERRQ(ierr); 1075bceba477SMatthew G. Knepley geom.v0 = v0; 1076bceba477SMatthew G. Knepley geom.J = J; 1077bceba477SMatthew G. Knepley geom.invJ = invJ; 1078bceba477SMatthew G. Knepley geom.detJ = detJ; 107932356553SMatthew G. Knepley ierr = PetscFEIntegrateInterpolator_Basic(feRef[fieldI], Ne, Nf, feRef, fieldI, fe, fieldJ, geom, g0, elemMat);CHKERRQ(ierr); 1080bceba477SMatthew G. Knepley geom.v0 = &v0[offset*dim]; 1081bceba477SMatthew G. Knepley geom.J = &J[offset*dim*dim]; 1082bceba477SMatthew G. Knepley geom.invJ = &invJ[offset*dim*dim]; 1083bceba477SMatthew G. Knepley geom.detJ = &detJ[offset]; 108432356553SMatthew G. Knepley ierr = PetscFEIntegrateInterpolator_Basic(feRef[fieldI], Nr, Nf, feRef, fieldI, fe, fieldJ, geom, g0, &elemMat[offset*rCellDof*cCellDof]);CHKERRQ(ierr); 1085bceba477SMatthew G. Knepley ierr = PetscFESetQuadrature(fe[fieldJ], quadOld);CHKERRQ(ierr); 1086bceba477SMatthew G. Knepley } 1087bceba477SMatthew G. Knepley } 1088bceba477SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 108932356553SMatthew G. Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, rCellDof, cCellDof, &elemMat[c*rCellDof*cCellDof]);CHKERRQ(ierr);} 109032356553SMatthew G. Knepley ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, I, c, &elemMat[c*rCellDof*cCellDof], ADD_VALUES);CHKERRQ(ierr); 1091bceba477SMatthew G. Knepley } 1092bceba477SMatthew G. Knepley ierr = PetscFree5(v0,J,invJ,detJ,elemMat);CHKERRQ(ierr); 109332356553SMatthew G. Knepley ierr = PetscFree(feRef);CHKERRQ(ierr); 1094bceba477SMatthew G. Knepley ierr = MatAssemblyBegin(I, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1095bceba477SMatthew G. Knepley ierr = MatAssemblyEnd(I, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1096bceba477SMatthew G. Knepley if (mesh->printFEM) { 1097bceba477SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 1098bceba477SMatthew G. Knepley ierr = MatChop(I, 1.0e-10);CHKERRQ(ierr); 1099bceba477SMatthew G. Knepley ierr = MatView(I, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1100bceba477SMatthew G. Knepley } 1101bceba477SMatthew G. Knepley #if 0 1102bceba477SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1103bceba477SMatthew G. Knepley #endif 1104bceba477SMatthew G. Knepley PetscFunctionReturn(0); 1105bceba477SMatthew G. Knepley } 1106d69c5d34SMatthew G. Knepley #endif 1107d69c5d34SMatthew G. Knepley 1108d69c5d34SMatthew G. Knepley #undef __FUNCT__ 1109d69c5d34SMatthew G. Knepley #define __FUNCT__ "DMPlexComputeInterpolatorFEM" 1110d69c5d34SMatthew G. Knepley /*@ 1111d69c5d34SMatthew G. Knepley DMPlexComputeInterpolatorFEM - Form the local portion of the interpolation matrix I from the coarse DM to the uniformly refined DM. 1112d69c5d34SMatthew G. Knepley 1113d69c5d34SMatthew G. Knepley Input Parameters: 1114d69c5d34SMatthew G. Knepley + dmf - The fine mesh 1115d69c5d34SMatthew G. Knepley . dmc - The coarse mesh 1116d69c5d34SMatthew G. Knepley - user - The user context 1117d69c5d34SMatthew G. Knepley 1118d69c5d34SMatthew G. Knepley Output Parameter: 1119*934789fcSMatthew G. Knepley . In - The interpolation matrix 1120d69c5d34SMatthew G. Knepley 1121d69c5d34SMatthew G. Knepley Note: 1122d69c5d34SMatthew G. Knepley The first member of the user context must be an FEMContext. 1123d69c5d34SMatthew G. Knepley 1124d69c5d34SMatthew G. Knepley We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 1125d69c5d34SMatthew G. Knepley like a GPU, or vectorize on a multicore machine. 1126d69c5d34SMatthew G. Knepley 1127d69c5d34SMatthew G. Knepley Level: developer 1128d69c5d34SMatthew G. Knepley 1129d69c5d34SMatthew G. Knepley .seealso: DMPlexComputeJacobianFEM() 1130d69c5d34SMatthew G. Knepley @*/ 1131*934789fcSMatthew G. Knepley PetscErrorCode DMPlexComputeInterpolatorFEM(DM dmc, DM dmf, Mat In, void *user) 1132d69c5d34SMatthew G. Knepley { 1133d69c5d34SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dmc->data; 1134d69c5d34SMatthew G. Knepley PetscFEM *fem = (PetscFEM *) user; 1135d69c5d34SMatthew G. Knepley PetscFE *fe = fem->fe; 1136d69c5d34SMatthew G. Knepley const char *name = "Interpolator"; 1137d69c5d34SMatthew G. Knepley PetscFE *feRef; 1138d69c5d34SMatthew G. Knepley PetscSection fsection, fglobalSection; 1139d69c5d34SMatthew G. Knepley PetscSection csection, cglobalSection; 1140d69c5d34SMatthew G. Knepley PetscScalar *elemMat; 1141942a7a06SMatthew G. Knepley PetscInt dim, Nf, f, fieldI, fieldJ, offsetI, offsetJ, cStart, cEnd, c; 1142d69c5d34SMatthew G. Knepley PetscInt rCellDof = 0, cCellDof = 0; 1143d69c5d34SMatthew G. Knepley PetscErrorCode ierr; 1144d69c5d34SMatthew G. Knepley 1145d69c5d34SMatthew G. Knepley PetscFunctionBegin; 1146d69c5d34SMatthew G. Knepley #if 0 1147d69c5d34SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1148d69c5d34SMatthew G. Knepley #endif 1149d69c5d34SMatthew G. Knepley ierr = DMPlexGetDimension(dmf, &dim);CHKERRQ(ierr); 1150d69c5d34SMatthew G. Knepley ierr = DMGetDefaultSection(dmf, &fsection);CHKERRQ(ierr); 1151d69c5d34SMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dmf, &fglobalSection);CHKERRQ(ierr); 1152d69c5d34SMatthew G. Knepley ierr = DMGetDefaultSection(dmc, &csection);CHKERRQ(ierr); 1153d69c5d34SMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dmc, &cglobalSection);CHKERRQ(ierr); 1154d69c5d34SMatthew G. Knepley ierr = PetscSectionGetNumFields(fsection, &Nf);CHKERRQ(ierr); 1155d69c5d34SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dmc, 0, &cStart, &cEnd);CHKERRQ(ierr); 1156d69c5d34SMatthew G. Knepley ierr = PetscMalloc1(Nf,&feRef);CHKERRQ(ierr); 1157d69c5d34SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 1158d69c5d34SMatthew G. Knepley PetscInt rNb, cNb, Nc; 1159d69c5d34SMatthew G. Knepley 1160d69c5d34SMatthew G. Knepley ierr = PetscFERefine(fe[f], &feRef[f]);CHKERRQ(ierr); 1161d69c5d34SMatthew G. Knepley ierr = PetscFEGetDimension(feRef[f], &rNb);CHKERRQ(ierr); 1162d69c5d34SMatthew G. Knepley ierr = PetscFEGetDimension(fe[f], &cNb);CHKERRQ(ierr); 1163d69c5d34SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 1164d69c5d34SMatthew G. Knepley rCellDof += rNb*Nc; 1165d69c5d34SMatthew G. Knepley cCellDof += cNb*Nc; 1166d69c5d34SMatthew G. Knepley } 1167*934789fcSMatthew G. Knepley ierr = MatZeroEntries(In);CHKERRQ(ierr); 1168d69c5d34SMatthew G. Knepley ierr = PetscMalloc1(rCellDof*cCellDof,&elemMat);CHKERRQ(ierr); 1169d69c5d34SMatthew G. Knepley ierr = PetscMemzero(elemMat, rCellDof*cCellDof * sizeof(PetscScalar));CHKERRQ(ierr); 1170d69c5d34SMatthew G. Knepley for (fieldI = 0, offsetI = 0; fieldI < Nf; ++fieldI) { 1171d69c5d34SMatthew G. Knepley PetscDualSpace Qref; 1172d69c5d34SMatthew G. Knepley PetscQuadrature f; 1173d69c5d34SMatthew G. Knepley const PetscReal *qpoints, *qweights; 1174d69c5d34SMatthew G. Knepley PetscReal *points; 1175d69c5d34SMatthew G. Knepley PetscInt npoints = 0, Nc, Np, fpdim, i, k, p, d; 1176d69c5d34SMatthew G. Knepley 1177d69c5d34SMatthew G. Knepley /* Compose points from all dual basis functionals */ 1178d69c5d34SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[fieldI], &Nc);CHKERRQ(ierr); 1179d69c5d34SMatthew G. Knepley ierr = PetscFEGetDualSpace(feRef[fieldI], &Qref);CHKERRQ(ierr); 1180d69c5d34SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(Qref, &fpdim);CHKERRQ(ierr); 1181d69c5d34SMatthew G. Knepley for (i = 0; i < fpdim; ++i) { 1182d69c5d34SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1183d69c5d34SMatthew G. Knepley ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, NULL);CHKERRQ(ierr); 1184d69c5d34SMatthew G. Knepley npoints += Np; 1185d69c5d34SMatthew G. Knepley } 1186d69c5d34SMatthew G. Knepley ierr = PetscMalloc1(npoints*dim,&points);CHKERRQ(ierr); 1187d69c5d34SMatthew G. Knepley for (i = 0, k = 0; i < fpdim; ++i) { 1188d69c5d34SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1189d69c5d34SMatthew G. Knepley ierr = PetscQuadratureGetData(f, NULL, &Np, &qpoints, NULL);CHKERRQ(ierr); 1190d69c5d34SMatthew G. Knepley for (p = 0; p < Np; ++p, ++k) for (d = 0; d < dim; ++d) points[k*dim+d] = qpoints[p*dim+d]; 1191d69c5d34SMatthew G. Knepley } 1192d69c5d34SMatthew G. Knepley 1193d69c5d34SMatthew G. Knepley for (fieldJ = 0, offsetJ = 0; fieldJ < Nf; ++fieldJ) { 1194d69c5d34SMatthew G. Knepley PetscSpace P; 1195d69c5d34SMatthew G. Knepley PetscReal *B; 1196d69c5d34SMatthew G. Knepley PetscInt cpdim, j; 1197d69c5d34SMatthew G. Knepley 1198d69c5d34SMatthew G. Knepley /* Evaluate basis at points */ 1199d69c5d34SMatthew G. Knepley ierr = PetscFEGetBasisSpace(fe[fieldJ], &P);CHKERRQ(ierr); 1200d69c5d34SMatthew G. Knepley ierr = PetscFEGetDimension(fe[fieldJ], &cpdim);CHKERRQ(ierr); 1201d69c5d34SMatthew G. Knepley ierr = PetscFEGetTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr); 1202d69c5d34SMatthew G. Knepley for (i = 0, k = 0; i < fpdim; ++i) { 1203d69c5d34SMatthew G. Knepley ierr = PetscDualSpaceGetFunctional(Qref, i, &f);CHKERRQ(ierr); 1204d69c5d34SMatthew G. Knepley ierr = PetscQuadratureGetData(f, NULL, &Np, NULL, &qweights);CHKERRQ(ierr); 1205d69c5d34SMatthew G. Knepley for (p = 0; p < Np; ++p, ++k) { 1206d69c5d34SMatthew G. Knepley for (j = 0; j < cpdim; ++j) elemMat[(offsetI + i)*cCellDof + offsetJ+j] += B[k*cpdim+j]*qweights[p]; 1207d69c5d34SMatthew G. Knepley } 1208d69c5d34SMatthew G. Knepley } 1209d69c5d34SMatthew G. Knepley ierr = PetscFERestoreTabulation(fe[fieldJ], npoints, points, &B, NULL, NULL);CHKERRQ(ierr);CHKERRQ(ierr); 1210d69c5d34SMatthew G. Knepley offsetJ += cpdim; 1211d69c5d34SMatthew G. Knepley } 1212d69c5d34SMatthew G. Knepley offsetI += fpdim*Nc; 1213549a8adaSMatthew G. Knepley ierr = PetscFree(points);CHKERRQ(ierr); 1214d69c5d34SMatthew G. Knepley } 1215d69c5d34SMatthew G. Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(0, name, rCellDof, cCellDof, elemMat);CHKERRQ(ierr);} 1216d69c5d34SMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 1217*934789fcSMatthew G. Knepley ierr = DMPlexMatSetClosureRefined(dmf, fsection, fglobalSection, dmc, csection, cglobalSection, In, c, elemMat, INSERT_VALUES);CHKERRQ(ierr); 1218d69c5d34SMatthew G. Knepley } 1219549a8adaSMatthew G. Knepley for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&feRef[f]);CHKERRQ(ierr);} 1220d69c5d34SMatthew G. Knepley ierr = PetscFree(feRef);CHKERRQ(ierr); 1221549a8adaSMatthew G. Knepley ierr = PetscFree(elemMat);CHKERRQ(ierr); 1222*934789fcSMatthew G. Knepley ierr = MatAssemblyBegin(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1223*934789fcSMatthew G. Knepley ierr = MatAssemblyEnd(In, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1224d69c5d34SMatthew G. Knepley if (mesh->printFEM) { 1225d69c5d34SMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 1226*934789fcSMatthew G. Knepley ierr = MatChop(In, 1.0e-10);CHKERRQ(ierr); 1227*934789fcSMatthew G. Knepley ierr = MatView(In, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1228d69c5d34SMatthew G. Knepley } 1229d69c5d34SMatthew G. Knepley #if 0 1230d69c5d34SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_InterpolatorFEM,dmc,dmf,0,0);CHKERRQ(ierr); 1231d69c5d34SMatthew G. Knepley #endif 1232d69c5d34SMatthew G. Knepley PetscFunctionReturn(0); 1233d69c5d34SMatthew G. Knepley } 12346c73c22cSMatthew G. Knepley 12356c73c22cSMatthew G. Knepley #undef __FUNCT__ 12366c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexAddBoundary" 12376c73c22cSMatthew G. Knepley /* The ids can be overridden by the command line option -bc_<boundary name> */ 12386c73c22cSMatthew G. Knepley PetscErrorCode DMPlexAddBoundary(DM dm, PetscBool isEssential, const char name[], PetscInt field, void (*bcFunc)(), PetscInt numids, const PetscInt *ids, void *ctx) 12396c73c22cSMatthew G. Knepley { 12406c73c22cSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 12416c73c22cSMatthew G. Knepley DMBoundary b; 12426c73c22cSMatthew G. Knepley PetscErrorCode ierr; 12436c73c22cSMatthew G. Knepley 12446c73c22cSMatthew G. Knepley PetscFunctionBegin; 12456c73c22cSMatthew G. Knepley ierr = PetscNew(&b);CHKERRQ(ierr); 12466c73c22cSMatthew G. Knepley ierr = PetscStrallocpy(name, (char **) &b->name);CHKERRQ(ierr); 12476c73c22cSMatthew G. Knepley ierr = PetscMalloc1(numids, &b->ids);CHKERRQ(ierr); 12486c73c22cSMatthew G. Knepley ierr = PetscMemcpy(b->ids, ids, numids*sizeof(PetscInt));CHKERRQ(ierr); 12496c73c22cSMatthew G. Knepley b->essential = isEssential; 12506c73c22cSMatthew G. Knepley b->field = field; 12516c73c22cSMatthew G. Knepley b->func = bcFunc; 12526c73c22cSMatthew G. Knepley b->numids = numids; 12536c73c22cSMatthew G. Knepley b->ctx = ctx; 12546c73c22cSMatthew G. Knepley b->next = mesh->boundary; 12556c73c22cSMatthew G. Knepley mesh->boundary = b; 12566c73c22cSMatthew G. Knepley PetscFunctionReturn(0); 12576c73c22cSMatthew G. Knepley } 12586c73c22cSMatthew G. Knepley 12596c73c22cSMatthew G. Knepley #undef __FUNCT__ 12606c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexGetNumBoundary" 12616c73c22cSMatthew G. Knepley PetscErrorCode DMPlexGetNumBoundary(DM dm, PetscInt *numBd) 12626c73c22cSMatthew G. Knepley { 12636c73c22cSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 12646c73c22cSMatthew G. Knepley DMBoundary b = mesh->boundary; 12656c73c22cSMatthew G. Knepley 12666c73c22cSMatthew G. Knepley PetscFunctionBegin; 12676c73c22cSMatthew G. Knepley *numBd = 0; 12686c73c22cSMatthew G. Knepley while (b) {++(*numBd); b = b->next;} 12696c73c22cSMatthew G. Knepley PetscFunctionReturn(0); 12706c73c22cSMatthew G. Knepley } 12716c73c22cSMatthew G. Knepley 12726c73c22cSMatthew G. Knepley #undef __FUNCT__ 12736c73c22cSMatthew G. Knepley #define __FUNCT__ "DMPlexGetBoundary" 12746c73c22cSMatthew G. Knepley PetscErrorCode DMPlexGetBoundary(DM dm, PetscInt bd, PetscBool *isEssential, const char **name, PetscInt *field, void (**func)(), PetscInt *numids, const PetscInt **ids, void **ctx) 12756c73c22cSMatthew G. Knepley { 12766c73c22cSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 12776c73c22cSMatthew G. Knepley DMBoundary b = mesh->boundary; 12786c73c22cSMatthew G. Knepley PetscInt n = 0; 12796c73c22cSMatthew G. Knepley 12806c73c22cSMatthew G. Knepley PetscFunctionBegin; 12816c73c22cSMatthew G. Knepley while (b) { 12826c73c22cSMatthew G. Knepley if (n == bd) break; 12836c73c22cSMatthew G. Knepley b = b->next; 12846c73c22cSMatthew G. Knepley ++n; 12856c73c22cSMatthew G. Knepley } 12866c73c22cSMatthew G. Knepley if (n != bd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Boundary %d is not in [0, %d)", bd, n); 12876c73c22cSMatthew G. Knepley if (isEssential) { 12886c73c22cSMatthew G. Knepley PetscValidPointer(isEssential, 3); 12896c73c22cSMatthew G. Knepley *isEssential = b->essential; 12906c73c22cSMatthew G. Knepley } 12916c73c22cSMatthew G. Knepley if (name) { 12926c73c22cSMatthew G. Knepley PetscValidPointer(name, 4); 12936c73c22cSMatthew G. Knepley *name = b->name; 12946c73c22cSMatthew G. Knepley } 12956c73c22cSMatthew G. Knepley if (field) { 12966c73c22cSMatthew G. Knepley PetscValidPointer(field, 5); 12976c73c22cSMatthew G. Knepley *field = b->field; 12986c73c22cSMatthew G. Knepley } 12996c73c22cSMatthew G. Knepley if (func) { 13006c73c22cSMatthew G. Knepley PetscValidPointer(func, 6); 13016c73c22cSMatthew G. Knepley *func = b->func; 13026c73c22cSMatthew G. Knepley } 13036c73c22cSMatthew G. Knepley if (numids) { 13046c73c22cSMatthew G. Knepley PetscValidPointer(numids, 7); 13056c73c22cSMatthew G. Knepley *numids = b->numids; 13066c73c22cSMatthew G. Knepley } 13076c73c22cSMatthew G. Knepley if (ids) { 13086c73c22cSMatthew G. Knepley PetscValidPointer(ids, 8); 13096c73c22cSMatthew G. Knepley *ids = b->ids; 13106c73c22cSMatthew G. Knepley } 13116c73c22cSMatthew G. Knepley if (ctx) { 13126c73c22cSMatthew G. Knepley PetscValidPointer(ctx, 9); 13136c73c22cSMatthew G. Knepley *ctx = b->ctx; 13146c73c22cSMatthew G. Knepley } 13156c73c22cSMatthew G. Knepley PetscFunctionReturn(0); 13166c73c22cSMatthew G. Knepley } 1317