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); 126cb1e1211SMatthew G Knepley ierr = DMPlexGetCoordinateSection(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" 18572f94c41SMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLocal(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), 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); 19772f94c41SMatthew G. Knepley ierr = PetscMalloc(numFields * sizeof(PetscDualSpace), &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); 209cb1e1211SMatthew G Knepley ierr = PetscMalloc2(dim,PetscReal,&v0,dim*dim,PetscReal,&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) { 21872f94c41SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[f], &numComp);CHKERRQ(ierr); 21972f94c41SMatthew G. Knepley ierr = PetscDualSpaceGetDimension(sp[f], &spDim);CHKERRQ(ierr); 22072f94c41SMatthew G. Knepley for (d = 0; d < spDim; ++d) { 22172f94c41SMatthew G. Knepley ierr = PetscDualSpaceApply(sp[f], d, geom, numComp, funcs[f], &values[v]);CHKERRQ(ierr); 22272f94c41SMatthew G. Knepley v += numComp; 223cb1e1211SMatthew G Knepley } 224cb1e1211SMatthew G Knepley } 22572f94c41SMatthew G. Knepley ierr = DMPlexVecSetClosure(dm, section, localX, c, values, mode);CHKERRQ(ierr); 226cb1e1211SMatthew G Knepley } 22772f94c41SMatthew G. Knepley ierr = DMRestoreWorkArray(dm, numValues, PETSC_SCALAR, &values);CHKERRQ(ierr); 228*1f2da991SMatthew G. Knepley ierr = PetscFree2(v0,J);CHKERRQ(ierr); 22972f94c41SMatthew G. Knepley ierr = PetscFree(sp);CHKERRQ(ierr); 230cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 231cb1e1211SMatthew G Knepley } 232cb1e1211SMatthew G Knepley 233cb1e1211SMatthew G Knepley #undef __FUNCT__ 234cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunction" 235cb1e1211SMatthew G Knepley /*@C 236cb1e1211SMatthew G Knepley DMPlexProjectFunction - This projects the given function into the function space provided. 237cb1e1211SMatthew G Knepley 238cb1e1211SMatthew G Knepley Input Parameters: 239cb1e1211SMatthew G Knepley + dm - The DM 24072f94c41SMatthew G. Knepley . fe - The PetscFE associated with the field 24172f94c41SMatthew G. Knepley . funcs - The coordinate functions to evaluate, one per field 242cb1e1211SMatthew G Knepley - mode - The insertion mode for values 243cb1e1211SMatthew G Knepley 244cb1e1211SMatthew G Knepley Output Parameter: 245cb1e1211SMatthew G Knepley . X - vector 246cb1e1211SMatthew G Knepley 247cb1e1211SMatthew G Knepley Level: developer 248cb1e1211SMatthew G Knepley 249878cb397SSatish Balay .seealso: DMPlexComputeL2Diff() 250878cb397SSatish Balay @*/ 25172f94c41SMatthew G. Knepley PetscErrorCode DMPlexProjectFunction(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec X) 252cb1e1211SMatthew G Knepley { 253cb1e1211SMatthew G Knepley Vec localX; 254cb1e1211SMatthew G Knepley PetscErrorCode ierr; 255cb1e1211SMatthew G Knepley 256cb1e1211SMatthew G Knepley PetscFunctionBegin; 257cb1e1211SMatthew G Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 25872f94c41SMatthew G. Knepley ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, mode, localX);CHKERRQ(ierr); 259cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr); 260cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr); 261cb1e1211SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 262cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 263cb1e1211SMatthew G Knepley } 264cb1e1211SMatthew G Knepley 265cb1e1211SMatthew G Knepley #undef __FUNCT__ 266cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff" 267cb1e1211SMatthew G Knepley /*@C 268cb1e1211SMatthew G Knepley DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 269cb1e1211SMatthew G Knepley 270cb1e1211SMatthew G Knepley Input Parameters: 271cb1e1211SMatthew G Knepley + dm - The DM 272c5bbbd5bSMatthew G. Knepley . fe - The PetscFE object for each field 273cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component 274cb1e1211SMatthew G Knepley - X - The coefficient vector u_h 275cb1e1211SMatthew G Knepley 276cb1e1211SMatthew G Knepley Output Parameter: 277cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2 278cb1e1211SMatthew G Knepley 279cb1e1211SMatthew G Knepley Level: developer 280cb1e1211SMatthew G Knepley 281cb1e1211SMatthew G Knepley .seealso: DMPlexProjectFunction() 282878cb397SSatish Balay @*/ 283c5bbbd5bSMatthew G. Knepley PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), Vec X, PetscReal *diff) 284cb1e1211SMatthew G Knepley { 285cb1e1211SMatthew G Knepley const PetscInt debug = 0; 286cb1e1211SMatthew G Knepley PetscSection section; 287c5bbbd5bSMatthew G. Knepley PetscQuadrature quad; 288cb1e1211SMatthew G Knepley Vec localX; 28972f94c41SMatthew G. Knepley PetscScalar *funcVal; 290cb1e1211SMatthew G Knepley PetscReal *coords, *v0, *J, *invJ, detJ; 291cb1e1211SMatthew G Knepley PetscReal localDiff = 0.0; 292cb1e1211SMatthew G Knepley PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 293cb1e1211SMatthew G Knepley PetscErrorCode ierr; 294cb1e1211SMatthew G Knepley 295cb1e1211SMatthew G Knepley PetscFunctionBegin; 296cb1e1211SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 297cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 298cb1e1211SMatthew G Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 299cb1e1211SMatthew G Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 300cb1e1211SMatthew G Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 301cb1e1211SMatthew G Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 302cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 303c5bbbd5bSMatthew G. Knepley PetscInt Nc; 304c5bbbd5bSMatthew G. Knepley 305c5bbbd5bSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 306c5bbbd5bSMatthew G. Knepley numComponents += Nc; 307cb1e1211SMatthew G Knepley } 30872f94c41SMatthew G. Knepley ierr = DMPlexProjectFunctionLocal(dm, fe, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 30972f94c41SMatthew G. Knepley ierr = PetscMalloc5(numComponents,PetscScalar,&funcVal,dim,PetscReal,&coords,dim,PetscReal,&v0,dim*dim,PetscReal,&J,dim*dim,PetscReal,&invJ);CHKERRQ(ierr); 310cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 311c5bbbd5bSMatthew G. Knepley ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr); 312cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 313a1e44745SMatthew G. Knepley PetscScalar *x = NULL; 314cb1e1211SMatthew G Knepley PetscReal elemDiff = 0.0; 315cb1e1211SMatthew G Knepley 316cb1e1211SMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 317cb1e1211SMatthew G Knepley if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 318cb1e1211SMatthew G Knepley ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 319cb1e1211SMatthew G Knepley 320cb1e1211SMatthew G Knepley for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 321f9fd7fdbSMatthew G. Knepley const PetscInt numQuadPoints = quad.numPoints; 322f9fd7fdbSMatthew G. Knepley const PetscReal *quadPoints = quad.points; 323f9fd7fdbSMatthew G. Knepley const PetscReal *quadWeights = quad.weights; 324c5bbbd5bSMatthew G. Knepley PetscReal *basis; 325c5bbbd5bSMatthew G. Knepley PetscInt numBasisFuncs, numBasisComps, q, d, e, fc, f; 326cb1e1211SMatthew G Knepley 327c5bbbd5bSMatthew G. Knepley ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr); 328c5bbbd5bSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr); 329c5bbbd5bSMatthew G. Knepley ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr); 330cb1e1211SMatthew G Knepley if (debug) { 331cb1e1211SMatthew G Knepley char title[1024]; 332cb1e1211SMatthew G Knepley ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 333cb1e1211SMatthew G Knepley ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr); 334cb1e1211SMatthew G Knepley } 335cb1e1211SMatthew G Knepley for (q = 0; q < numQuadPoints; ++q) { 336cb1e1211SMatthew G Knepley for (d = 0; d < dim; d++) { 337cb1e1211SMatthew G Knepley coords[d] = v0[d]; 338cb1e1211SMatthew G Knepley for (e = 0; e < dim; e++) { 339cb1e1211SMatthew G Knepley coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 340cb1e1211SMatthew G Knepley } 341cb1e1211SMatthew G Knepley } 34272f94c41SMatthew G. Knepley (*funcs[field])(coords, funcVal); 343cb1e1211SMatthew G Knepley for (fc = 0; fc < numBasisComps; ++fc) { 344a1d24da5SMatthew G. Knepley PetscScalar interpolant = 0.0; 345a1d24da5SMatthew G. Knepley 346cb1e1211SMatthew G Knepley for (f = 0; f < numBasisFuncs; ++f) { 347cb1e1211SMatthew G Knepley const PetscInt fidx = f*numBasisComps+fc; 348a1d24da5SMatthew G. Knepley interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx]; 349cb1e1211SMatthew G Knepley } 35072f94c41SMatthew 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);} 35172f94c41SMatthew G. Knepley elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal[fc]))*quadWeights[q]*detJ; 352cb1e1211SMatthew G Knepley } 353cb1e1211SMatthew G Knepley } 354cb1e1211SMatthew G Knepley comp += numBasisComps; 355cb1e1211SMatthew G Knepley fieldOffset += numBasisFuncs*numBasisComps; 356cb1e1211SMatthew G Knepley } 357cb1e1211SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 358cb1e1211SMatthew G Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 359cb1e1211SMatthew G Knepley localDiff += elemDiff; 360cb1e1211SMatthew G Knepley } 36172f94c41SMatthew G. Knepley ierr = PetscFree5(funcVal,coords,v0,J,invJ);CHKERRQ(ierr); 362cb1e1211SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 36386a74ee0SMatthew G. Knepley ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 364cb1e1211SMatthew G Knepley *diff = PetscSqrtReal(*diff); 365cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 366cb1e1211SMatthew G Knepley } 367cb1e1211SMatthew G Knepley 368a0845e3aSMatthew G. Knepley #if 0 369a0845e3aSMatthew G. Knepley 370cb1e1211SMatthew G Knepley #undef __FUNCT__ 371cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeResidualFEM" 372cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 373cb1e1211SMatthew G Knepley { 374cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 3759a559087SMatthew G. Knepley PetscFEM *fem = (PetscFEM *) user; 376cb1e1211SMatthew G Knepley PetscQuadrature *quad = fem->quad; 377652b88e8SMatthew G. Knepley PetscQuadrature *quadBd = fem->quadBd; 378cb1e1211SMatthew G Knepley PetscSection section; 37902a80efeSMatthew G. Knepley PetscReal *v0, *n, *J, *invJ, *detJ; 380cb1e1211SMatthew G Knepley PetscScalar *elemVec, *u; 381cb1e1211SMatthew G Knepley PetscInt dim, numFields, field, numBatchesTmp = 1, numCells, cStart, cEnd, c; 382652b88e8SMatthew G. Knepley PetscInt cellDof, numComponents; 383652b88e8SMatthew G. Knepley PetscBool has; 384cb1e1211SMatthew G Knepley PetscErrorCode ierr; 385cb1e1211SMatthew G Knepley 386cb1e1211SMatthew G Knepley PetscFunctionBegin; 387652b88e8SMatthew G. Knepley if (has && quadBd) { 388652b88e8SMatthew G. Knepley DMLabel label; 389652b88e8SMatthew G. Knepley IS pointIS; 390652b88e8SMatthew G. Knepley const PetscInt *points; 391652b88e8SMatthew G. Knepley PetscInt numPoints, p; 392652b88e8SMatthew G. Knepley 393652b88e8SMatthew G. Knepley ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr); 394652b88e8SMatthew G. Knepley ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr); 395652b88e8SMatthew G. Knepley ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr); 396652b88e8SMatthew G. Knepley ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 397652b88e8SMatthew G. Knepley for (field = 0, cellDof = 0, numComponents = 0; field < numFields; ++field) { 398652b88e8SMatthew G. Knepley cellDof += quadBd[field].numBasisFuncs*quadBd[field].numComponents; 399652b88e8SMatthew G. Knepley numComponents += quadBd[field].numComponents; 400652b88e8SMatthew G. Knepley } 40102a80efeSMatthew G. Knepley ierr = PetscMalloc7(numPoints*cellDof,PetscScalar,&u,numPoints*dim,PetscReal,&v0,numPoints*dim,PetscReal,&n,numPoints*dim*dim,PetscReal,&J,numPoints*dim*dim,PetscReal,&invJ,numPoints,PetscReal,&detJ,numPoints*cellDof,PetscScalar,&elemVec);CHKERRQ(ierr); 402652b88e8SMatthew G. Knepley for (p = 0; p < numPoints; ++p) { 403652b88e8SMatthew G. Knepley const PetscInt point = points[p]; 404652b88e8SMatthew G. Knepley PetscScalar *x; 405652b88e8SMatthew G. Knepley PetscInt i; 406652b88e8SMatthew G. Knepley 40702a80efeSMatthew G. Knepley /* TODO: Add normal determination here */ 408652b88e8SMatthew G. Knepley ierr = DMPlexComputeCellGeometry(dm, point, &v0[p*dim], &J[p*dim*dim], &invJ[p*dim*dim], &detJ[p]);CHKERRQ(ierr); 4091d930511SMatthew G. Knepley if (detJ[p] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[p], point); 410652b88e8SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr); 411652b88e8SMatthew G. Knepley 412652b88e8SMatthew G. Knepley for (i = 0; i < cellDof; ++i) u[p*cellDof+i] = x[i]; 413652b88e8SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr); 414652b88e8SMatthew G. Knepley } 415652b88e8SMatthew G. Knepley for (field = 0; field < numFields; ++field) { 416652b88e8SMatthew G. Knepley const PetscInt numQuadPoints = quadBd[field].numQuadPoints; 417652b88e8SMatthew G. Knepley const PetscInt numBasisFuncs = quadBd[field].numBasisFuncs; 418a9dc2124SMatthew G. Knepley void (*f0)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[field]; 419a9dc2124SMatthew G. Knepley void (*f1)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[field]; 420652b88e8SMatthew G. Knepley /* Conforming batches */ 421652b88e8SMatthew G. Knepley PetscInt blockSize = numBasisFuncs*numQuadPoints; 422652b88e8SMatthew G. Knepley PetscInt numBlocks = 1; 423652b88e8SMatthew G. Knepley PetscInt batchSize = numBlocks * blockSize; 424652b88e8SMatthew G. Knepley PetscInt numBatches = numBatchesTmp; 425652b88e8SMatthew G. Knepley PetscInt numChunks = numPoints / (numBatches*batchSize); 426652b88e8SMatthew G. Knepley /* Remainder */ 427652b88e8SMatthew G. Knepley PetscInt numRemainder = numPoints % (numBatches * batchSize); 428652b88e8SMatthew G. Knepley PetscInt offset = numPoints - numRemainder; 429652b88e8SMatthew G. Knepley 43002a80efeSMatthew G. Knepley ierr = (*mesh->integrateBdResidualFEM)(numChunks*numBatches*batchSize, numFields, field, quadBd, u, v0, n, J, invJ, detJ, f0, f1, elemVec);CHKERRQ(ierr); 43102a80efeSMatthew G. Knepley ierr = (*mesh->integrateBdResidualFEM)(numRemainder, numFields, field, quadBd, &u[offset*cellDof], &v0[offset*dim], &n[offset*dim], &J[offset*dim*dim], &invJ[offset*dim*dim], &detJ[offset], 432652b88e8SMatthew G. Knepley f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 433652b88e8SMatthew G. Knepley } 434652b88e8SMatthew G. Knepley for (p = 0; p < numPoints; ++p) { 435652b88e8SMatthew G. Knepley const PetscInt point = points[p]; 436652b88e8SMatthew G. Knepley 437652b88e8SMatthew G. Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "Residual", cellDof, &elemVec[p*cellDof]);CHKERRQ(ierr);} 438652b88e8SMatthew G. Knepley ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[p*cellDof], ADD_VALUES);CHKERRQ(ierr); 439652b88e8SMatthew G. Knepley } 440652b88e8SMatthew G. Knepley ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 441652b88e8SMatthew G. Knepley ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 44202a80efeSMatthew G. Knepley ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr); 443652b88e8SMatthew G. Knepley } 444cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 445cb1e1211SMatthew G Knepley } 446cb1e1211SMatthew G Knepley 447a0845e3aSMatthew G. Knepley #else 448a0845e3aSMatthew G. Knepley 449a0845e3aSMatthew G. Knepley #undef __FUNCT__ 450a0845e3aSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeResidualFEM" 451a0845e3aSMatthew G. Knepley /*@ 452a0845e3aSMatthew G. Knepley DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user 453a0845e3aSMatthew G. Knepley 454a0845e3aSMatthew G. Knepley Input Parameters: 455a0845e3aSMatthew G. Knepley + dm - The mesh 456a0845e3aSMatthew G. Knepley . X - Local input vector 457a0845e3aSMatthew G. Knepley - user - The user context 458a0845e3aSMatthew G. Knepley 459a0845e3aSMatthew G. Knepley Output Parameter: 460a0845e3aSMatthew G. Knepley . F - Local output vector 461a0845e3aSMatthew G. Knepley 462a0845e3aSMatthew G. Knepley Note: 463a0845e3aSMatthew G. Knepley The second member of the user context must be an FEMContext. 464a0845e3aSMatthew G. Knepley 465a0845e3aSMatthew G. Knepley We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 466a0845e3aSMatthew G. Knepley like a GPU, or vectorize on a multicore machine. 467a0845e3aSMatthew G. Knepley 468a0845e3aSMatthew G. Knepley Level: developer 469a0845e3aSMatthew G. Knepley 470a0845e3aSMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM() 471a0845e3aSMatthew G. Knepley @*/ 472a0845e3aSMatthew G. Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 473a0845e3aSMatthew G. Knepley { 474a0845e3aSMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 4759a559087SMatthew G. Knepley PetscFEM *fem = (PetscFEM *) user; 476a0845e3aSMatthew G. Knepley PetscFE *fe = fem->fe; 4779a559087SMatthew G. Knepley PetscFE *feAux = fem->feAux; 478f1ea0e2fSMatthew G. Knepley PetscFE *feBd = fem->feBd; 479a0845e3aSMatthew G. Knepley const char *name = "Residual"; 4809a559087SMatthew G. Knepley DM dmAux; 4819a559087SMatthew G. Knepley Vec A; 482a0845e3aSMatthew G. Knepley PetscQuadrature q; 483a0845e3aSMatthew G. Knepley PetscCellGeometry geom; 4849a559087SMatthew G. Knepley PetscSection section, sectionAux; 485a0845e3aSMatthew G. Knepley PetscReal *v0, *J, *invJ, *detJ; 4869a559087SMatthew G. Knepley PetscScalar *elemVec, *u, *a; 4879a559087SMatthew G. Knepley PetscInt dim, Nf, NfAux = 0, f, numCells, cStart, cEnd, c; 488a0845e3aSMatthew G. Knepley PetscInt cellDof = 0, numComponents = 0; 4899a559087SMatthew G. Knepley PetscInt cellDofAux = 0, numComponentsAux = 0; 490a0845e3aSMatthew G. Knepley PetscErrorCode ierr; 491a0845e3aSMatthew G. Knepley 492a0845e3aSMatthew G. Knepley PetscFunctionBegin; 493a0845e3aSMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 494a0845e3aSMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 495a0845e3aSMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 4969a559087SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 497a0845e3aSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 498a0845e3aSMatthew G. Knepley numCells = cEnd - cStart; 4999a559087SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 500a0845e3aSMatthew G. Knepley PetscInt Nb, Nc; 501a0845e3aSMatthew G. Knepley 502a0845e3aSMatthew G. Knepley ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 503a0845e3aSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 504a0845e3aSMatthew G. Knepley cellDof += Nb*Nc; 505a0845e3aSMatthew G. Knepley numComponents += Nc; 506a0845e3aSMatthew G. Knepley } 5079a559087SMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 5089a559087SMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 5099a559087SMatthew G. Knepley if (dmAux) { 5109a559087SMatthew G. Knepley ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 5119a559087SMatthew G. Knepley ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr); 5129a559087SMatthew G. Knepley } 5139a559087SMatthew G. Knepley for (f = 0; f < NfAux; ++f) { 5149a559087SMatthew G. Knepley PetscInt Nb, Nc; 5159a559087SMatthew G. Knepley 5169a559087SMatthew G. Knepley ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr); 5179a559087SMatthew G. Knepley ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr); 5189a559087SMatthew G. Knepley cellDofAux += Nb*Nc; 5199a559087SMatthew G. Knepley numComponentsAux += Nc; 5209a559087SMatthew G. Knepley } 52172f94c41SMatthew G. Knepley ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, INSERT_BC_VALUES, X);CHKERRQ(ierr); 522a0845e3aSMatthew G. Knepley ierr = VecSet(F, 0.0);CHKERRQ(ierr); 523a0845e3aSMatthew G. Knepley ierr = PetscMalloc6(numCells*cellDof,PetscScalar,&u,numCells*dim,PetscReal,&v0,numCells*dim*dim,PetscReal,&J,numCells*dim*dim,PetscReal,&invJ,numCells,PetscReal,&detJ,numCells*cellDof,PetscScalar,&elemVec);CHKERRQ(ierr); 5249a559087SMatthew G. Knepley if (dmAux) {ierr = PetscMalloc(numCells*cellDofAux * sizeof(PetscScalar), &a);CHKERRQ(ierr);} 525a0845e3aSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 526a0845e3aSMatthew G. Knepley PetscScalar *x = NULL; 527a0845e3aSMatthew G. Knepley PetscInt i; 528a0845e3aSMatthew G. Knepley 529a0845e3aSMatthew G. Knepley ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 530a0845e3aSMatthew G. Knepley if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 531a0845e3aSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 532a0845e3aSMatthew G. Knepley for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 533a0845e3aSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 5349a559087SMatthew G. Knepley if (dmAux) { 5359a559087SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 5369a559087SMatthew G. Knepley for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i]; 5379a559087SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 538a0845e3aSMatthew G. Knepley } 5399a559087SMatthew G. Knepley } 5409a559087SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 541c012ea0aSMatthew G. Knepley void (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f0Funcs[f]; 542c012ea0aSMatthew G. Knepley void (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->f1Funcs[f]; 543a0845e3aSMatthew G. Knepley PetscInt Nb; 544a0845e3aSMatthew G. Knepley /* Conforming batches */ 545f30c5766SMatthew G. Knepley PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 546a0845e3aSMatthew G. Knepley /* Remainder */ 547a0845e3aSMatthew G. Knepley PetscInt Nr, offset; 548a0845e3aSMatthew G. Knepley 549a0845e3aSMatthew G. Knepley ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr); 550a0845e3aSMatthew G. Knepley ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 551f30c5766SMatthew G. Knepley ierr = PetscFEGetTileSizes(fe[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 552f9fd7fdbSMatthew G. Knepley blockSize = Nb*q.numPoints; 553a0845e3aSMatthew G. Knepley batchSize = numBlocks * blockSize; 554f30c5766SMatthew G. Knepley ierr = PetscFESetTileSizes(fe[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 555a0845e3aSMatthew G. Knepley numChunks = numCells / (numBatches*batchSize); 556a0845e3aSMatthew G. Knepley Ne = numChunks*numBatches*batchSize; 557a0845e3aSMatthew G. Knepley Nr = numCells % (numBatches*batchSize); 558a0845e3aSMatthew G. Knepley offset = numCells - Nr; 559a0845e3aSMatthew G. Knepley geom.v0 = v0; 560a0845e3aSMatthew G. Knepley geom.J = J; 561a0845e3aSMatthew G. Knepley geom.invJ = invJ; 562a0845e3aSMatthew G. Knepley geom.detJ = detJ; 5639a559087SMatthew G. Knepley ierr = PetscFEIntegrateResidual(fe[f], Ne, Nf, fe, f, geom, u, NfAux, feAux, a, f0, f1, elemVec);CHKERRQ(ierr); 564a0845e3aSMatthew G. Knepley geom.v0 = &v0[offset*dim]; 565a0845e3aSMatthew G. Knepley geom.J = &J[offset*dim*dim]; 566a0845e3aSMatthew G. Knepley geom.invJ = &invJ[offset*dim*dim]; 567a0845e3aSMatthew G. Knepley geom.detJ = &detJ[offset]; 5689a559087SMatthew 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); 569a0845e3aSMatthew G. Knepley } 570a0845e3aSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 571a0845e3aSMatthew G. Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 572a0845e3aSMatthew G. Knepley ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 573a0845e3aSMatthew G. Knepley } 574a0845e3aSMatthew G. Knepley ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 5759a559087SMatthew G. Knepley if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 576f1ea0e2fSMatthew G. Knepley if (feBd) { 577f1ea0e2fSMatthew G. Knepley DMLabel label; 578f1ea0e2fSMatthew G. Knepley IS pointIS; 579f1ea0e2fSMatthew G. Knepley const PetscInt *points; 580f1ea0e2fSMatthew G. Knepley PetscInt numPoints, p; 581f1ea0e2fSMatthew G. Knepley PetscReal *n; 582f1ea0e2fSMatthew G. Knepley 583f1ea0e2fSMatthew G. Knepley ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr); 584f1ea0e2fSMatthew G. Knepley ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr); 585f1ea0e2fSMatthew G. Knepley ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr); 586f1ea0e2fSMatthew G. Knepley ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 587f1ea0e2fSMatthew G. Knepley for (f = 0, cellDof = 0, numComponents = 0; f < Nf; ++f) { 588f1ea0e2fSMatthew G. Knepley PetscInt Nb, Nc; 589f1ea0e2fSMatthew G. Knepley 590f1ea0e2fSMatthew G. Knepley ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr); 591f1ea0e2fSMatthew G. Knepley ierr = PetscFEGetNumComponents(feBd[f], &Nc);CHKERRQ(ierr); 592f1ea0e2fSMatthew G. Knepley cellDof += Nb*Nc; 593f1ea0e2fSMatthew G. Knepley numComponents += Nc; 594f1ea0e2fSMatthew G. Knepley } 595f1ea0e2fSMatthew G. Knepley ierr = PetscMalloc7(numPoints*cellDof,PetscScalar,&u,numPoints*dim,PetscReal,&v0,numPoints*dim,PetscReal,&n,numPoints*dim*dim,PetscReal,&J,numPoints*dim*dim,PetscReal,&invJ,numPoints,PetscReal,&detJ,numPoints*cellDof,PetscScalar,&elemVec);CHKERRQ(ierr); 596f1ea0e2fSMatthew G. Knepley for (p = 0; p < numPoints; ++p) { 597f1ea0e2fSMatthew G. Knepley const PetscInt point = points[p]; 598f1ea0e2fSMatthew G. Knepley PetscScalar *x = NULL; 599f1ea0e2fSMatthew G. Knepley PetscInt i; 600f1ea0e2fSMatthew G. Knepley 601f1ea0e2fSMatthew G. Knepley /* TODO: Add normal determination here */ 602f1ea0e2fSMatthew G. Knepley ierr = DMPlexComputeCellGeometry(dm, point, &v0[p*dim], &J[p*dim*dim], &invJ[p*dim*dim], &detJ[p]);CHKERRQ(ierr); 603f1ea0e2fSMatthew G. Knepley if (detJ[p] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[p], point); 604f1ea0e2fSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 605f1ea0e2fSMatthew G. Knepley for (i = 0; i < cellDof; ++i) u[p*cellDof+i] = x[i]; 606f1ea0e2fSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, section, X, point, NULL, &x);CHKERRQ(ierr); 607f1ea0e2fSMatthew G. Knepley } 608f1ea0e2fSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 609f1ea0e2fSMatthew G. Knepley void (*f0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[f]; 610f1ea0e2fSMatthew G. Knepley void (*f1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[f]; 611f1ea0e2fSMatthew G. Knepley PetscInt Nb; 612f1ea0e2fSMatthew G. Knepley /* Conforming batches */ 613f1ea0e2fSMatthew G. Knepley PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 614f1ea0e2fSMatthew G. Knepley /* Remainder */ 615f1ea0e2fSMatthew G. Knepley PetscInt Nr, offset; 616f1ea0e2fSMatthew G. Knepley 617f1ea0e2fSMatthew G. Knepley ierr = PetscFEGetQuadrature(feBd[f], &q);CHKERRQ(ierr); 618f1ea0e2fSMatthew G. Knepley ierr = PetscFEGetDimension(feBd[f], &Nb);CHKERRQ(ierr); 619f1ea0e2fSMatthew G. Knepley ierr = PetscFEGetTileSizes(feBd[f], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 620f9fd7fdbSMatthew G. Knepley blockSize = Nb*q.numPoints; 621f1ea0e2fSMatthew G. Knepley batchSize = numBlocks * blockSize; 622f1ea0e2fSMatthew G. Knepley ierr = PetscFESetTileSizes(feBd[f], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 623f1ea0e2fSMatthew G. Knepley numChunks = numPoints / (numBatches*batchSize); 624f1ea0e2fSMatthew G. Knepley Ne = numChunks*numBatches*batchSize; 625f1ea0e2fSMatthew G. Knepley Nr = numPoints % (numBatches*batchSize); 626f1ea0e2fSMatthew G. Knepley offset = numPoints - Nr; 627f1ea0e2fSMatthew G. Knepley geom.v0 = v0; 628f1ea0e2fSMatthew G. Knepley geom.n = n; 629f1ea0e2fSMatthew G. Knepley geom.J = J; 630f1ea0e2fSMatthew G. Knepley geom.invJ = invJ; 631f1ea0e2fSMatthew G. Knepley geom.detJ = detJ; 632f1ea0e2fSMatthew G. Knepley ierr = PetscFEIntegrateBdResidual(feBd[f], Ne, Nf, feBd, f, geom, u, 0, NULL, NULL, f0, f1, elemVec);CHKERRQ(ierr); 633f1ea0e2fSMatthew G. Knepley geom.v0 = &v0[offset*dim]; 634f1ea0e2fSMatthew G. Knepley geom.n = &n[offset*dim]; 635f1ea0e2fSMatthew G. Knepley geom.J = &J[offset*dim*dim]; 636f1ea0e2fSMatthew G. Knepley geom.invJ = &invJ[offset*dim*dim]; 637f1ea0e2fSMatthew G. Knepley geom.detJ = &detJ[offset]; 638f1ea0e2fSMatthew G. Knepley ierr = PetscFEIntegrateBdResidual(feBd[f], Nr, Nf, feBd, f, geom, &u[offset*cellDof], 0, NULL, NULL, f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 639f1ea0e2fSMatthew G. Knepley } 640f1ea0e2fSMatthew G. Knepley for (p = 0; p < numPoints; ++p) { 641f1ea0e2fSMatthew G. Knepley const PetscInt point = points[p]; 642f1ea0e2fSMatthew G. Knepley 643f1ea0e2fSMatthew G. Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "BdResidual", cellDof, &elemVec[p*cellDof]);CHKERRQ(ierr);} 644f1ea0e2fSMatthew G. Knepley ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[p*cellDof], ADD_VALUES);CHKERRQ(ierr); 645f1ea0e2fSMatthew G. Knepley } 646f1ea0e2fSMatthew G. Knepley ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 647f1ea0e2fSMatthew G. Knepley ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 648f1ea0e2fSMatthew G. Knepley ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr); 649f1ea0e2fSMatthew G. Knepley } 6506113b454SMatthew G. Knepley if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, mesh->printTol, F);CHKERRQ(ierr);} 651a0845e3aSMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 652a0845e3aSMatthew G. Knepley PetscFunctionReturn(0); 653a0845e3aSMatthew G. Knepley } 654a0845e3aSMatthew G. Knepley 655a0845e3aSMatthew G. Knepley #endif 656a0845e3aSMatthew G. Knepley 657cb1e1211SMatthew G Knepley #undef __FUNCT__ 658cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianActionFEM" 659cb1e1211SMatthew G Knepley /*@C 660cb1e1211SMatthew G Knepley DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user 661cb1e1211SMatthew G Knepley 662cb1e1211SMatthew G Knepley Input Parameters: 663cb1e1211SMatthew G Knepley + dm - The mesh 664cb1e1211SMatthew G Knepley . J - The Jacobian shell matrix 665cb1e1211SMatthew G Knepley . X - Local input vector 666cb1e1211SMatthew G Knepley - user - The user context 667cb1e1211SMatthew G Knepley 668cb1e1211SMatthew G Knepley Output Parameter: 669cb1e1211SMatthew G Knepley . F - Local output vector 670cb1e1211SMatthew G Knepley 671cb1e1211SMatthew G Knepley Note: 672cb1e1211SMatthew G Knepley The second member of the user context must be an FEMContext. 673cb1e1211SMatthew G Knepley 674cb1e1211SMatthew G Knepley We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 675cb1e1211SMatthew G Knepley like a GPU, or vectorize on a multicore machine. 676cb1e1211SMatthew G Knepley 6770059ad2aSSatish Balay Level: developer 6780059ad2aSSatish Balay 679cb1e1211SMatthew G Knepley .seealso: DMPlexComputeResidualFEM() 680878cb397SSatish Balay @*/ 681cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user) 682cb1e1211SMatthew G Knepley { 683cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 6849a559087SMatthew G. Knepley PetscFEM *fem = (PetscFEM *) user; 6850483ade4SMatthew G. Knepley PetscFE *fe = fem->fe; 6860483ade4SMatthew G. Knepley PetscQuadrature quad; 6870483ade4SMatthew G. Knepley PetscCellGeometry geom; 688cb1e1211SMatthew G Knepley PetscSection section; 689cb1e1211SMatthew G Knepley JacActionCtx *jctx; 690cb1e1211SMatthew G Knepley PetscReal *v0, *J, *invJ, *detJ; 691cb1e1211SMatthew G Knepley PetscScalar *elemVec, *u, *a; 6920483ade4SMatthew G. Knepley PetscInt dim, numFields, field, numCells, cStart, cEnd, c; 693cb1e1211SMatthew G Knepley PetscInt cellDof = 0; 694cb1e1211SMatthew G Knepley PetscErrorCode ierr; 695cb1e1211SMatthew G Knepley 696cb1e1211SMatthew G Knepley PetscFunctionBegin; 6970483ade4SMatthew G. Knepley /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */ 698cb1e1211SMatthew G Knepley ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 699cb1e1211SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 700cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 701cb1e1211SMatthew G Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 702cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 703cb1e1211SMatthew G Knepley numCells = cEnd - cStart; 704cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 7050483ade4SMatthew G. Knepley PetscInt Nb, Nc; 7060483ade4SMatthew G. Knepley 7070483ade4SMatthew G. Knepley ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 7080483ade4SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 7090483ade4SMatthew G. Knepley cellDof += Nb*Nc; 710cb1e1211SMatthew G Knepley } 711cb1e1211SMatthew G Knepley ierr = VecSet(F, 0.0);CHKERRQ(ierr); 712cb1e1211SMatthew G Knepley ierr = PetscMalloc7(numCells*cellDof,PetscScalar,&u,numCells*cellDof,PetscScalar,&a,numCells*dim,PetscReal,&v0,numCells*dim*dim,PetscReal,&J,numCells*dim*dim,PetscReal,&invJ,numCells,PetscReal,&detJ,numCells*cellDof,PetscScalar,&elemVec);CHKERRQ(ierr); 713cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 714a1e44745SMatthew G. Knepley PetscScalar *x = NULL; 715cb1e1211SMatthew G Knepley PetscInt i; 716cb1e1211SMatthew G Knepley 717cb1e1211SMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 718cb1e1211SMatthew G Knepley if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 719cb1e1211SMatthew G Knepley ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 720cb1e1211SMatthew G Knepley for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 721cb1e1211SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 722cb1e1211SMatthew G Knepley ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 723cb1e1211SMatthew G Knepley for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i]; 724cb1e1211SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 725cb1e1211SMatthew G Knepley } 726cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 7270483ade4SMatthew G. Knepley PetscInt Nb; 728cb1e1211SMatthew G Knepley /* Conforming batches */ 729cb1e1211SMatthew G Knepley PetscInt numBlocks = 1; 7300483ade4SMatthew G. Knepley PetscInt numBatches = 1; 7310483ade4SMatthew G. Knepley PetscInt numChunks, Ne, blockSize, batchSize; 732cb1e1211SMatthew G Knepley /* Remainder */ 7330483ade4SMatthew G. Knepley PetscInt Nr, offset; 734cb1e1211SMatthew G Knepley 7350483ade4SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr); 7360483ade4SMatthew G. Knepley ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 737f9fd7fdbSMatthew G. Knepley blockSize = Nb*quad.numPoints; 7380483ade4SMatthew G. Knepley batchSize = numBlocks * blockSize; 7390483ade4SMatthew G. Knepley numChunks = numCells / (numBatches*batchSize); 7400483ade4SMatthew G. Knepley Ne = numChunks*numBatches*batchSize; 7410483ade4SMatthew G. Knepley Nr = numCells % (numBatches*batchSize); 7420483ade4SMatthew G. Knepley offset = numCells - Nr; 7430483ade4SMatthew G. Knepley geom.v0 = v0; 7440483ade4SMatthew G. Knepley geom.J = J; 7450483ade4SMatthew G. Knepley geom.invJ = invJ; 7460483ade4SMatthew G. Knepley geom.detJ = detJ; 7470483ade4SMatthew G. Knepley ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr); 7480483ade4SMatthew G. Knepley geom.v0 = &v0[offset*dim]; 7490483ade4SMatthew G. Knepley geom.J = &J[offset*dim*dim]; 7500483ade4SMatthew G. Knepley geom.invJ = &invJ[offset*dim*dim]; 7510483ade4SMatthew G. Knepley geom.detJ = &detJ[offset]; 7520483ade4SMatthew G. Knepley ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof], 753cb1e1211SMatthew G Knepley fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr); 754cb1e1211SMatthew G Knepley } 755cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 756cb1e1211SMatthew G Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 757cb1e1211SMatthew G Knepley ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 758cb1e1211SMatthew G Knepley } 759cb1e1211SMatthew G Knepley ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 760cb1e1211SMatthew G Knepley if (mesh->printFEM) { 761cb1e1211SMatthew G Knepley PetscMPIInt rank, numProcs; 762cb1e1211SMatthew G Knepley PetscInt p; 763cb1e1211SMatthew G Knepley 764cb1e1211SMatthew G Knepley ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 765cb1e1211SMatthew G Knepley ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr); 76686a74ee0SMatthew G. Knepley ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr); 767cb1e1211SMatthew G Knepley for (p = 0; p < numProcs; ++p) { 768cb1e1211SMatthew G Knepley if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);} 769cb1e1211SMatthew G Knepley ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr); 770cb1e1211SMatthew G Knepley } 771cb1e1211SMatthew G Knepley } 7720483ade4SMatthew G. Knepley /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */ 773cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 774cb1e1211SMatthew G Knepley } 775cb1e1211SMatthew G Knepley 776cb1e1211SMatthew G Knepley #undef __FUNCT__ 777cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianFEM" 778cb1e1211SMatthew G Knepley /*@ 779cb1e1211SMatthew G Knepley DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user. 780cb1e1211SMatthew G Knepley 781cb1e1211SMatthew G Knepley Input Parameters: 782cb1e1211SMatthew G Knepley + dm - The mesh 783cb1e1211SMatthew G Knepley . X - Local input vector 784cb1e1211SMatthew G Knepley - user - The user context 785cb1e1211SMatthew G Knepley 786cb1e1211SMatthew G Knepley Output Parameter: 787cb1e1211SMatthew G Knepley . Jac - Jacobian matrix 788cb1e1211SMatthew G Knepley 789cb1e1211SMatthew G Knepley Note: 790cb1e1211SMatthew G Knepley The second member of the user context must be an FEMContext. 791cb1e1211SMatthew G Knepley 792cb1e1211SMatthew G Knepley We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 793cb1e1211SMatthew G Knepley like a GPU, or vectorize on a multicore machine. 794cb1e1211SMatthew G Knepley 7950059ad2aSSatish Balay Level: developer 7960059ad2aSSatish Balay 797cb1e1211SMatthew G Knepley .seealso: FormFunctionLocal() 798878cb397SSatish Balay @*/ 799cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, MatStructure *str,void *user) 800cb1e1211SMatthew G Knepley { 801cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 8029a559087SMatthew G. Knepley PetscFEM *fem = (PetscFEM *) user; 803a319912fSMatthew G. Knepley PetscFE *fe = fem->fe; 804754551f4SMatthew G. Knepley PetscFE *feAux = fem->feAux; 805a319912fSMatthew G. Knepley const char *name = "Jacobian"; 806754551f4SMatthew G. Knepley DM dmAux; 807754551f4SMatthew G. Knepley Vec A; 808a319912fSMatthew G. Knepley PetscQuadrature quad; 809a319912fSMatthew G. Knepley PetscCellGeometry geom; 810754551f4SMatthew G. Knepley PetscSection section, globalSection, sectionAux; 811cb1e1211SMatthew G Knepley PetscReal *v0, *J, *invJ, *detJ; 812754551f4SMatthew G. Knepley PetscScalar *elemMat, *u, *a; 813754551f4SMatthew G. Knepley PetscInt dim, Nf, NfAux = 0, f, fieldI, fieldJ, numCells, cStart, cEnd, c; 814cb1e1211SMatthew G Knepley PetscInt cellDof = 0, numComponents = 0; 815754551f4SMatthew G. Knepley PetscInt cellDofAux = 0, numComponentsAux = 0; 816cb1e1211SMatthew G Knepley PetscBool isShell; 817cb1e1211SMatthew G Knepley PetscErrorCode ierr; 818cb1e1211SMatthew G Knepley 819cb1e1211SMatthew G Knepley PetscFunctionBegin; 820a319912fSMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 821cb1e1211SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 822cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 823a319912fSMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr); 824754551f4SMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &Nf);CHKERRQ(ierr); 825cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 826cb1e1211SMatthew G Knepley numCells = cEnd - cStart; 827754551f4SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 828a319912fSMatthew G. Knepley PetscInt Nb, Nc; 829a319912fSMatthew G. Knepley 830a319912fSMatthew G. Knepley ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 831a319912fSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 832a319912fSMatthew G. Knepley cellDof += Nb*Nc; 833a319912fSMatthew G. Knepley numComponents += Nc; 834cb1e1211SMatthew G Knepley } 835754551f4SMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "dmAux", (PetscObject *) &dmAux);CHKERRQ(ierr); 836754551f4SMatthew G. Knepley ierr = PetscObjectQuery((PetscObject) dm, "A", (PetscObject *) &A);CHKERRQ(ierr); 837754551f4SMatthew G. Knepley if (dmAux) { 838754551f4SMatthew G. Knepley ierr = DMGetDefaultSection(dmAux, §ionAux);CHKERRQ(ierr); 839754551f4SMatthew G. Knepley ierr = PetscSectionGetNumFields(sectionAux, &NfAux);CHKERRQ(ierr); 840754551f4SMatthew G. Knepley } 841754551f4SMatthew G. Knepley for (f = 0; f < NfAux; ++f) { 842754551f4SMatthew G. Knepley PetscInt Nb, Nc; 843754551f4SMatthew G. Knepley 844754551f4SMatthew G. Knepley ierr = PetscFEGetDimension(feAux[f], &Nb);CHKERRQ(ierr); 845754551f4SMatthew G. Knepley ierr = PetscFEGetNumComponents(feAux[f], &Nc);CHKERRQ(ierr); 846754551f4SMatthew G. Knepley cellDofAux += Nb*Nc; 847754551f4SMatthew G. Knepley numComponentsAux += Nc; 848754551f4SMatthew G. Knepley } 84972f94c41SMatthew G. Knepley ierr = DMPlexProjectFunctionLocal(dm, fe, fem->bcFuncs, INSERT_BC_VALUES, X);CHKERRQ(ierr); 850cb1e1211SMatthew G Knepley ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 851cb1e1211SMatthew G Knepley ierr = PetscMalloc6(numCells*cellDof,PetscScalar,&u,numCells*dim,PetscReal,&v0,numCells*dim*dim,PetscReal,&J,numCells*dim*dim,PetscReal,&invJ,numCells,PetscReal,&detJ,numCells*cellDof*cellDof,PetscScalar,&elemMat);CHKERRQ(ierr); 852754551f4SMatthew G. Knepley if (dmAux) {ierr = PetscMalloc(numCells*cellDofAux * sizeof(PetscScalar), &a);CHKERRQ(ierr);} 853cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 854a1e44745SMatthew G. Knepley PetscScalar *x = NULL; 855cb1e1211SMatthew G Knepley PetscInt i; 856cb1e1211SMatthew G Knepley 857cb1e1211SMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 858cb1e1211SMatthew G Knepley if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 859a319912fSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 860cb1e1211SMatthew G Knepley for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 861a319912fSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 862754551f4SMatthew G. Knepley if (dmAux) { 863754551f4SMatthew G. Knepley ierr = DMPlexVecGetClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 864754551f4SMatthew G. Knepley for (i = 0; i < cellDofAux; ++i) a[c*cellDofAux+i] = x[i]; 865754551f4SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dmAux, sectionAux, A, c, NULL, &x);CHKERRQ(ierr); 866754551f4SMatthew G. Knepley } 867cb1e1211SMatthew G Knepley } 868cb1e1211SMatthew G Knepley ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr); 869754551f4SMatthew G. Knepley for (fieldI = 0; fieldI < Nf; ++fieldI) { 870a319912fSMatthew G. Knepley PetscInt Nb; 871cb1e1211SMatthew G Knepley /* Conforming batches */ 872754551f4SMatthew G. Knepley PetscInt numChunks, numBatches, numBlocks, Ne, blockSize, batchSize; 873cb1e1211SMatthew G Knepley /* Remainder */ 874a319912fSMatthew G. Knepley PetscInt Nr, offset; 875cb1e1211SMatthew G Knepley 876754551f4SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr); 877754551f4SMatthew G. Knepley ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr); 878754551f4SMatthew G. Knepley ierr = PetscFEGetTileSizes(fe[fieldI], NULL, &numBlocks, NULL, &numBatches);CHKERRQ(ierr); 879f9fd7fdbSMatthew G. Knepley blockSize = Nb*quad.numPoints; 880a319912fSMatthew G. Knepley batchSize = numBlocks * blockSize; 881754551f4SMatthew G. Knepley ierr = PetscFESetTileSizes(fe[fieldI], blockSize, numBlocks, batchSize, numBatches);CHKERRQ(ierr); 882a319912fSMatthew G. Knepley numChunks = numCells / (numBatches*batchSize); 883a319912fSMatthew G. Knepley Ne = numChunks*numBatches*batchSize; 884a319912fSMatthew G. Knepley Nr = numCells % (numBatches*batchSize); 885a319912fSMatthew G. Knepley offset = numCells - Nr; 886754551f4SMatthew G. Knepley for (fieldJ = 0; fieldJ < Nf; ++fieldJ) { 887754551f4SMatthew G. Knepley void (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*Nf+fieldJ]; 888754551f4SMatthew G. Knepley void (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*Nf+fieldJ]; 889754551f4SMatthew G. Knepley void (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*Nf+fieldJ]; 890754551f4SMatthew G. Knepley void (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*Nf+fieldJ]; 891754551f4SMatthew G. Knepley 892a319912fSMatthew G. Knepley geom.v0 = v0; 893a319912fSMatthew G. Knepley geom.J = J; 894a319912fSMatthew G. Knepley geom.invJ = invJ; 895a319912fSMatthew G. Knepley geom.detJ = detJ; 896754551f4SMatthew G. Knepley ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, Nf, fe, fieldI, fieldJ, geom, u, NfAux, feAux, a, g0, g1, g2, g3, elemMat);CHKERRQ(ierr); 897a319912fSMatthew G. Knepley geom.v0 = &v0[offset*dim]; 898a319912fSMatthew G. Knepley geom.J = &J[offset*dim*dim]; 899a319912fSMatthew G. Knepley geom.invJ = &invJ[offset*dim*dim]; 900a319912fSMatthew G. Knepley geom.detJ = &detJ[offset]; 901754551f4SMatthew 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); 902cb1e1211SMatthew G Knepley } 903cb1e1211SMatthew G Knepley } 904cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 905a319912fSMatthew G. Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);} 906a319912fSMatthew G. Knepley ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr); 907cb1e1211SMatthew G Knepley } 908cb1e1211SMatthew G Knepley ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr); 909754551f4SMatthew G. Knepley if (dmAux) {ierr = PetscFree(a);CHKERRQ(ierr);} 910cb1e1211SMatthew G Knepley ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 911cb1e1211SMatthew G Knepley ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 912cb1e1211SMatthew G Knepley if (mesh->printFEM) { 913a319912fSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 914cb1e1211SMatthew G Knepley ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr); 915cb1e1211SMatthew G Knepley ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 916cb1e1211SMatthew G Knepley } 917a319912fSMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 918cb1e1211SMatthew G Knepley ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr); 919cb1e1211SMatthew G Knepley if (isShell) { 920cb1e1211SMatthew G Knepley JacActionCtx *jctx; 921cb1e1211SMatthew G Knepley 922cb1e1211SMatthew G Knepley ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 923cb1e1211SMatthew G Knepley ierr = VecCopy(X, jctx->u);CHKERRQ(ierr); 924cb1e1211SMatthew G Knepley } 925cb1e1211SMatthew G Knepley *str = SAME_NONZERO_PATTERN; 926cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 927cb1e1211SMatthew G Knepley } 928