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" 185a1d24da5SMatthew G. Knepley PetscErrorCode DMPlexProjectFunctionLocal(DM dm, PetscInt numComp, void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec localX) 186cb1e1211SMatthew G Knepley { 187cb1e1211SMatthew G Knepley Vec coordinates; 188cb1e1211SMatthew G Knepley PetscSection section, cSection; 189cb1e1211SMatthew G Knepley PetscInt dim, vStart, vEnd, v, c, d; 190cb1e1211SMatthew G Knepley PetscScalar *values, *cArray; 191cb1e1211SMatthew G Knepley PetscReal *coords; 192cb1e1211SMatthew G Knepley PetscErrorCode ierr; 193cb1e1211SMatthew G Knepley 194cb1e1211SMatthew G Knepley PetscFunctionBegin; 195cb1e1211SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 196cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 197cb1e1211SMatthew G Knepley ierr = DMPlexGetCoordinateSection(dm, &cSection);CHKERRQ(ierr); 198cb1e1211SMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 199cb1e1211SMatthew G Knepley ierr = PetscMalloc(numComp * sizeof(PetscScalar), &values);CHKERRQ(ierr); 200cb1e1211SMatthew G Knepley ierr = VecGetArray(coordinates, &cArray);CHKERRQ(ierr); 201cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(cSection, vStart, &dim);CHKERRQ(ierr); 202cb1e1211SMatthew G Knepley ierr = PetscMalloc(dim * sizeof(PetscReal),&coords);CHKERRQ(ierr); 203cb1e1211SMatthew G Knepley for (v = vStart; v < vEnd; ++v) { 204cb1e1211SMatthew G Knepley PetscInt dof, off; 205cb1e1211SMatthew G Knepley 206cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(cSection, v, &dof);CHKERRQ(ierr); 207cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(cSection, v, &off);CHKERRQ(ierr); 208cb1e1211SMatthew G Knepley if (dof > dim) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Cannot have more coordinates %d then dimensions %d", dof, dim); 209cb1e1211SMatthew G Knepley for (d = 0; d < dof; ++d) coords[d] = PetscRealPart(cArray[off+d]); 210a1d24da5SMatthew G. Knepley for (c = 0; c < numComp; ++c) (*funcs[c])(coords, &values[c]); 211cb1e1211SMatthew G Knepley ierr = VecSetValuesSection(localX, section, v, values, mode);CHKERRQ(ierr); 212cb1e1211SMatthew G Knepley } 213cb1e1211SMatthew G Knepley ierr = VecRestoreArray(coordinates, &cArray);CHKERRQ(ierr); 214cb1e1211SMatthew G Knepley /* Temporary, must be replaced by a projection on the finite element basis */ 215cb1e1211SMatthew G Knepley { 216cb1e1211SMatthew G Knepley PetscInt eStart = 0, eEnd = 0, e, depth; 217cb1e1211SMatthew G Knepley 218cb1e1211SMatthew G Knepley ierr = DMPlexGetLabelSize(dm, "depth", &depth);CHKERRQ(ierr); 219cb1e1211SMatthew G Knepley --depth; 220cb1e1211SMatthew G Knepley if (depth > 1) {ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr);} 221cb1e1211SMatthew G Knepley for (e = eStart; e < eEnd; ++e) { 222cb1e1211SMatthew G Knepley const PetscInt *cone = NULL; 223cb1e1211SMatthew G Knepley PetscInt coneSize, d; 224cb1e1211SMatthew G Knepley PetscScalar *coordsA, *coordsB; 225cb1e1211SMatthew G Knepley 226cb1e1211SMatthew G Knepley ierr = DMPlexGetConeSize(dm, e, &coneSize);CHKERRQ(ierr); 227cb1e1211SMatthew G Knepley ierr = DMPlexGetCone(dm, e, &cone);CHKERRQ(ierr); 228cb1e1211SMatthew G Knepley if (coneSize != 2) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Cone size %d for point %d should be 2", coneSize, e); 229cb1e1211SMatthew G Knepley ierr = VecGetValuesSection(coordinates, cSection, cone[0], &coordsA);CHKERRQ(ierr); 230cb1e1211SMatthew G Knepley ierr = VecGetValuesSection(coordinates, cSection, cone[1], &coordsB);CHKERRQ(ierr); 231cb1e1211SMatthew G Knepley for (d = 0; d < dim; ++d) { 232cb1e1211SMatthew G Knepley coords[d] = 0.5*(PetscRealPart(coordsA[d]) + PetscRealPart(coordsB[d])); 233cb1e1211SMatthew G Knepley } 234a1d24da5SMatthew G. Knepley for (c = 0; c < numComp; ++c) (*funcs[c])(coords, &values[c]); 235cb1e1211SMatthew G Knepley ierr = VecSetValuesSection(localX, section, e, values, mode);CHKERRQ(ierr); 236cb1e1211SMatthew G Knepley } 237cb1e1211SMatthew G Knepley } 238cb1e1211SMatthew G Knepley 239cb1e1211SMatthew G Knepley ierr = PetscFree(coords);CHKERRQ(ierr); 240cb1e1211SMatthew G Knepley ierr = PetscFree(values);CHKERRQ(ierr); 241cb1e1211SMatthew G Knepley #if 0 242cb1e1211SMatthew G Knepley const PetscInt localDof = this->_mesh->sizeWithBC(s, *cells->begin()); 243cb1e1211SMatthew G Knepley PetscReal detJ; 244cb1e1211SMatthew G Knepley 245cb1e1211SMatthew G Knepley ierr = PetscMalloc(localDof * sizeof(PetscScalar), &values);CHKERRQ(ierr); 246cb1e1211SMatthew G Knepley ierr = PetscMalloc2(dim,PetscReal,&v0,dim*dim,PetscReal,&J);CHKERRQ(ierr); 247cb1e1211SMatthew G Knepley ALE::ISieveVisitor::PointRetriever<PETSC_MESH_TYPE::sieve_type> pV(PetscPowInt(this->_mesh->getSieve()->getMaxConeSize(),dim+1), true); 248cb1e1211SMatthew G Knepley 249cb1e1211SMatthew G Knepley for (PetscInt c = cStart; c < cEnd; ++c) { 250cb1e1211SMatthew G Knepley ALE::ISieveTraversal<PETSC_MESH_TYPE::sieve_type>::orientedClosure(*this->_mesh->getSieve(), c, pV); 251cb1e1211SMatthew G Knepley const PETSC_MESH_TYPE::point_type *oPoints = pV.getPoints(); 252cb1e1211SMatthew G Knepley const int oSize = pV.getSize(); 253cb1e1211SMatthew G Knepley int v = 0; 254cb1e1211SMatthew G Knepley 255cb1e1211SMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr); 256cb1e1211SMatthew G Knepley for (PetscInt cl = 0; cl < oSize; ++cl) { 257cb1e1211SMatthew G Knepley const PetscInt fDim; 258cb1e1211SMatthew G Knepley 259cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(oPoints[cl], &fDim);CHKERRQ(ierr); 260cb1e1211SMatthew G Knepley if (pointDim) { 261cb1e1211SMatthew G Knepley for (PetscInt d = 0; d < fDim; ++d, ++v) { 262cb1e1211SMatthew G Knepley values[v] = (*this->_options.integrate)(v0, J, v, initFunc); 263cb1e1211SMatthew G Knepley } 264cb1e1211SMatthew G Knepley } 265cb1e1211SMatthew G Knepley } 266cb1e1211SMatthew G Knepley ierr = DMPlexVecSetClosure(dm, NULL, localX, c, values);CHKERRQ(ierr); 267cb1e1211SMatthew G Knepley pV.clear(); 268cb1e1211SMatthew G Knepley } 269cb1e1211SMatthew G Knepley ierr = PetscFree2(v0,J);CHKERRQ(ierr); 270cb1e1211SMatthew G Knepley ierr = PetscFree(values);CHKERRQ(ierr); 271cb1e1211SMatthew G Knepley #endif 272cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 273cb1e1211SMatthew G Knepley } 274cb1e1211SMatthew G Knepley 275cb1e1211SMatthew G Knepley #undef __FUNCT__ 276cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunction" 277cb1e1211SMatthew G Knepley /*@C 278cb1e1211SMatthew G Knepley DMPlexProjectFunction - This projects the given function into the function space provided. 279cb1e1211SMatthew G Knepley 280cb1e1211SMatthew G Knepley Input Parameters: 281cb1e1211SMatthew G Knepley + dm - The DM 282cb1e1211SMatthew G Knepley . numComp - The number of components (functions) 283cb1e1211SMatthew G Knepley . funcs - The coordinate functions to evaluate 284cb1e1211SMatthew G Knepley - mode - The insertion mode for values 285cb1e1211SMatthew G Knepley 286cb1e1211SMatthew G Knepley Output Parameter: 287cb1e1211SMatthew G Knepley . X - vector 288cb1e1211SMatthew G Knepley 289cb1e1211SMatthew G Knepley Level: developer 290cb1e1211SMatthew G Knepley 291cb1e1211SMatthew G Knepley Note: 292cb1e1211SMatthew G Knepley This currently just calls the function with the coordinates of each vertex and edge midpoint, and stores the result in a vector. 293cb1e1211SMatthew G Knepley We will eventually fix it. 294cb1e1211SMatthew G Knepley 295878cb397SSatish Balay .seealso: DMPlexComputeL2Diff() 296878cb397SSatish Balay @*/ 297a1d24da5SMatthew G. Knepley PetscErrorCode DMPlexProjectFunction(DM dm, PetscInt numComp, void (**funcs)(const PetscReal [], PetscScalar *), InsertMode mode, Vec X) 298cb1e1211SMatthew G Knepley { 299cb1e1211SMatthew G Knepley Vec localX; 300cb1e1211SMatthew G Knepley PetscErrorCode ierr; 301cb1e1211SMatthew G Knepley 302cb1e1211SMatthew G Knepley PetscFunctionBegin; 303cb1e1211SMatthew G Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 304cb1e1211SMatthew G Knepley ierr = DMPlexProjectFunctionLocal(dm, numComp, funcs, mode, localX);CHKERRQ(ierr); 305cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr); 306cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr); 307cb1e1211SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 308cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 309cb1e1211SMatthew G Knepley } 310cb1e1211SMatthew G Knepley 311cb1e1211SMatthew G Knepley #undef __FUNCT__ 312cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff" 313cb1e1211SMatthew G Knepley /*@C 314cb1e1211SMatthew G Knepley DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 315cb1e1211SMatthew G Knepley 316cb1e1211SMatthew G Knepley Input Parameters: 317cb1e1211SMatthew G Knepley + dm - The DM 318*c5bbbd5bSMatthew G. Knepley . fe - The PetscFE object for each field 319cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component 320cb1e1211SMatthew G Knepley - X - The coefficient vector u_h 321cb1e1211SMatthew G Knepley 322cb1e1211SMatthew G Knepley Output Parameter: 323cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2 324cb1e1211SMatthew G Knepley 325cb1e1211SMatthew G Knepley Level: developer 326cb1e1211SMatthew G Knepley 327cb1e1211SMatthew G Knepley .seealso: DMPlexProjectFunction() 328878cb397SSatish Balay @*/ 329*c5bbbd5bSMatthew G. Knepley PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscFE fe[], void (**funcs)(const PetscReal [], PetscScalar *), Vec X, PetscReal *diff) 330cb1e1211SMatthew G Knepley { 331cb1e1211SMatthew G Knepley const PetscInt debug = 0; 332cb1e1211SMatthew G Knepley PetscSection section; 333*c5bbbd5bSMatthew G. Knepley PetscQuadrature quad; 334cb1e1211SMatthew G Knepley Vec localX; 335cb1e1211SMatthew G Knepley PetscReal *coords, *v0, *J, *invJ, detJ; 336cb1e1211SMatthew G Knepley PetscReal localDiff = 0.0; 337cb1e1211SMatthew G Knepley PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 338cb1e1211SMatthew G Knepley PetscErrorCode ierr; 339cb1e1211SMatthew G Knepley 340cb1e1211SMatthew G Knepley PetscFunctionBegin; 341cb1e1211SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 342cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 343cb1e1211SMatthew G Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 344cb1e1211SMatthew G Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 345cb1e1211SMatthew G Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 346cb1e1211SMatthew G Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 347cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 348*c5bbbd5bSMatthew G. Knepley PetscInt Nc; 349*c5bbbd5bSMatthew G. Knepley 350*c5bbbd5bSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 351*c5bbbd5bSMatthew G. Knepley numComponents += Nc; 352cb1e1211SMatthew G Knepley } 353cb1e1211SMatthew G Knepley ierr = DMPlexProjectFunctionLocal(dm, numComponents, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 354cb1e1211SMatthew G Knepley ierr = PetscMalloc4(dim,PetscReal,&coords,dim,PetscReal,&v0,dim*dim,PetscReal,&J,dim*dim,PetscReal,&invJ);CHKERRQ(ierr); 355cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 356*c5bbbd5bSMatthew G. Knepley ierr = PetscFEGetQuadrature(fe[0], &quad);CHKERRQ(ierr); 357cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 358cb1e1211SMatthew G Knepley PetscScalar *x; 359cb1e1211SMatthew G Knepley PetscReal elemDiff = 0.0; 360cb1e1211SMatthew G Knepley 361cb1e1211SMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 362cb1e1211SMatthew G Knepley if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 363cb1e1211SMatthew G Knepley ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 364cb1e1211SMatthew G Knepley 365cb1e1211SMatthew G Knepley for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 366*c5bbbd5bSMatthew G. Knepley const PetscInt numQuadPoints = quad.numQuadPoints; 367*c5bbbd5bSMatthew G. Knepley const PetscReal *quadPoints = quad.quadPoints; 368*c5bbbd5bSMatthew G. Knepley const PetscReal *quadWeights = quad.quadWeights; 369*c5bbbd5bSMatthew G. Knepley PetscReal *basis; 370*c5bbbd5bSMatthew G. Knepley PetscInt numBasisFuncs, numBasisComps, q, d, e, fc, f; 371cb1e1211SMatthew G Knepley 372*c5bbbd5bSMatthew G. Knepley ierr = PetscFEGetDimension(fe[field], &numBasisFuncs);CHKERRQ(ierr); 373*c5bbbd5bSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[field], &numBasisComps);CHKERRQ(ierr); 374*c5bbbd5bSMatthew G. Knepley ierr = PetscFEGetDefaultTabulation(fe[field], &basis, NULL, NULL);CHKERRQ(ierr); 375cb1e1211SMatthew G Knepley if (debug) { 376cb1e1211SMatthew G Knepley char title[1024]; 377cb1e1211SMatthew G Knepley ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 378cb1e1211SMatthew G Knepley ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr); 379cb1e1211SMatthew G Knepley } 380cb1e1211SMatthew G Knepley for (q = 0; q < numQuadPoints; ++q) { 381cb1e1211SMatthew G Knepley for (d = 0; d < dim; d++) { 382cb1e1211SMatthew G Knepley coords[d] = v0[d]; 383cb1e1211SMatthew G Knepley for (e = 0; e < dim; e++) { 384cb1e1211SMatthew G Knepley coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 385cb1e1211SMatthew G Knepley } 386cb1e1211SMatthew G Knepley } 387cb1e1211SMatthew G Knepley for (fc = 0; fc < numBasisComps; ++fc) { 388a1d24da5SMatthew G. Knepley PetscScalar funcVal; 389a1d24da5SMatthew G. Knepley PetscScalar interpolant = 0.0; 390a1d24da5SMatthew G. Knepley 391a1d24da5SMatthew G. Knepley (*funcs[comp+fc])(coords, &funcVal); 392cb1e1211SMatthew G Knepley for (f = 0; f < numBasisFuncs; ++f) { 393cb1e1211SMatthew G Knepley const PetscInt fidx = f*numBasisComps+fc; 394a1d24da5SMatthew G. Knepley interpolant += x[fieldOffset+fidx]*basis[q*numBasisFuncs*numBasisComps+fidx]; 395cb1e1211SMatthew G Knepley } 396a1d24da5SMatthew G. Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d field %d diff %g\n", c, field, PetscSqr(PetscRealPart(interpolant - funcVal))*quadWeights[q]*detJ);CHKERRQ(ierr);} 397a1d24da5SMatthew G. Knepley elemDiff += PetscSqr(PetscRealPart(interpolant - funcVal))*quadWeights[q]*detJ; 398cb1e1211SMatthew G Knepley } 399cb1e1211SMatthew G Knepley } 400cb1e1211SMatthew G Knepley comp += numBasisComps; 401cb1e1211SMatthew G Knepley fieldOffset += numBasisFuncs*numBasisComps; 402cb1e1211SMatthew G Knepley } 403cb1e1211SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 404cb1e1211SMatthew G Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 405cb1e1211SMatthew G Knepley localDiff += elemDiff; 406cb1e1211SMatthew G Knepley } 407cb1e1211SMatthew G Knepley ierr = PetscFree4(coords,v0,J,invJ);CHKERRQ(ierr); 408cb1e1211SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 40986a74ee0SMatthew G. Knepley ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 410cb1e1211SMatthew G Knepley *diff = PetscSqrtReal(*diff); 411cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 412cb1e1211SMatthew G Knepley } 413cb1e1211SMatthew G Knepley 414a0845e3aSMatthew G. Knepley #if 0 415a0845e3aSMatthew G. Knepley 416cb1e1211SMatthew G Knepley #undef __FUNCT__ 417cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeResidualFEM" 418cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 419cb1e1211SMatthew G Knepley { 420cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 421cb1e1211SMatthew G Knepley PetscFEM *fem = (PetscFEM*) &((DM*) user)[1]; 422cb1e1211SMatthew G Knepley PetscQuadrature *quad = fem->quad; 423652b88e8SMatthew G. Knepley PetscQuadrature *quadBd = fem->quadBd; 424cb1e1211SMatthew G Knepley PetscSection section; 42502a80efeSMatthew G. Knepley PetscReal *v0, *n, *J, *invJ, *detJ; 426cb1e1211SMatthew G Knepley PetscScalar *elemVec, *u; 427cb1e1211SMatthew G Knepley PetscInt dim, numFields, field, numBatchesTmp = 1, numCells, cStart, cEnd, c; 428652b88e8SMatthew G. Knepley PetscInt cellDof, numComponents; 429652b88e8SMatthew G. Knepley PetscBool has; 430cb1e1211SMatthew G Knepley PetscErrorCode ierr; 431cb1e1211SMatthew G Knepley 432cb1e1211SMatthew G Knepley PetscFunctionBegin; 433652b88e8SMatthew G. Knepley if (has && quadBd) { 434652b88e8SMatthew G. Knepley DMLabel label; 435652b88e8SMatthew G. Knepley IS pointIS; 436652b88e8SMatthew G. Knepley const PetscInt *points; 437652b88e8SMatthew G. Knepley PetscInt numPoints, p; 438652b88e8SMatthew G. Knepley 439652b88e8SMatthew G. Knepley ierr = DMPlexGetLabel(dm, "boundary", &label);CHKERRQ(ierr); 440652b88e8SMatthew G. Knepley ierr = DMLabelGetStratumSize(label, 1, &numPoints);CHKERRQ(ierr); 441652b88e8SMatthew G. Knepley ierr = DMLabelGetStratumIS(label, 1, &pointIS);CHKERRQ(ierr); 442652b88e8SMatthew G. Knepley ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 443652b88e8SMatthew G. Knepley for (field = 0, cellDof = 0, numComponents = 0; field < numFields; ++field) { 444652b88e8SMatthew G. Knepley cellDof += quadBd[field].numBasisFuncs*quadBd[field].numComponents; 445652b88e8SMatthew G. Knepley numComponents += quadBd[field].numComponents; 446652b88e8SMatthew G. Knepley } 44702a80efeSMatthew 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); 448652b88e8SMatthew G. Knepley for (p = 0; p < numPoints; ++p) { 449652b88e8SMatthew G. Knepley const PetscInt point = points[p]; 450652b88e8SMatthew G. Knepley PetscScalar *x; 451652b88e8SMatthew G. Knepley PetscInt i; 452652b88e8SMatthew G. Knepley 45302a80efeSMatthew G. Knepley /* TODO: Add normal determination here */ 454652b88e8SMatthew G. Knepley ierr = DMPlexComputeCellGeometry(dm, point, &v0[p*dim], &J[p*dim*dim], &invJ[p*dim*dim], &detJ[p]);CHKERRQ(ierr); 4551d930511SMatthew G. Knepley if (detJ[p] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for face %d", detJ[p], point); 456652b88e8SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr); 457652b88e8SMatthew G. Knepley 458652b88e8SMatthew G. Knepley for (i = 0; i < cellDof; ++i) u[p*cellDof+i] = x[i]; 459652b88e8SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, X, point, NULL, &x);CHKERRQ(ierr); 460652b88e8SMatthew G. Knepley } 461652b88e8SMatthew G. Knepley for (field = 0; field < numFields; ++field) { 462652b88e8SMatthew G. Knepley const PetscInt numQuadPoints = quadBd[field].numQuadPoints; 463652b88e8SMatthew G. Knepley const PetscInt numBasisFuncs = quadBd[field].numBasisFuncs; 464a9dc2124SMatthew G. Knepley void (*f0)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f0BdFuncs[field]; 465a9dc2124SMatthew G. Knepley void (*f1)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]) = fem->f1BdFuncs[field]; 466652b88e8SMatthew G. Knepley /* Conforming batches */ 467652b88e8SMatthew G. Knepley PetscInt blockSize = numBasisFuncs*numQuadPoints; 468652b88e8SMatthew G. Knepley PetscInt numBlocks = 1; 469652b88e8SMatthew G. Knepley PetscInt batchSize = numBlocks * blockSize; 470652b88e8SMatthew G. Knepley PetscInt numBatches = numBatchesTmp; 471652b88e8SMatthew G. Knepley PetscInt numChunks = numPoints / (numBatches*batchSize); 472652b88e8SMatthew G. Knepley /* Remainder */ 473652b88e8SMatthew G. Knepley PetscInt numRemainder = numPoints % (numBatches * batchSize); 474652b88e8SMatthew G. Knepley PetscInt offset = numPoints - numRemainder; 475652b88e8SMatthew G. Knepley 47602a80efeSMatthew G. Knepley ierr = (*mesh->integrateBdResidualFEM)(numChunks*numBatches*batchSize, numFields, field, quadBd, u, v0, n, J, invJ, detJ, f0, f1, elemVec);CHKERRQ(ierr); 47702a80efeSMatthew 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], 478652b88e8SMatthew G. Knepley f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 479652b88e8SMatthew G. Knepley } 480652b88e8SMatthew G. Knepley for (p = 0; p < numPoints; ++p) { 481652b88e8SMatthew G. Knepley const PetscInt point = points[p]; 482652b88e8SMatthew G. Knepley 483652b88e8SMatthew G. Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellVector(point, "Residual", cellDof, &elemVec[p*cellDof]);CHKERRQ(ierr);} 484652b88e8SMatthew G. Knepley ierr = DMPlexVecSetClosure(dm, NULL, F, point, &elemVec[p*cellDof], ADD_VALUES);CHKERRQ(ierr); 485652b88e8SMatthew G. Knepley } 486652b88e8SMatthew G. Knepley ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 487652b88e8SMatthew G. Knepley ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 48802a80efeSMatthew G. Knepley ierr = PetscFree7(u,v0,n,J,invJ,detJ,elemVec);CHKERRQ(ierr); 489652b88e8SMatthew G. Knepley } 490cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 491cb1e1211SMatthew G Knepley } 492cb1e1211SMatthew G Knepley 493a0845e3aSMatthew G. Knepley #else 494a0845e3aSMatthew G. Knepley 495a0845e3aSMatthew G. Knepley #undef __FUNCT__ 496a0845e3aSMatthew G. Knepley #define __FUNCT__ "DMPlexComputeResidualFEM" 497a0845e3aSMatthew G. Knepley /*@ 498a0845e3aSMatthew G. Knepley DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user 499a0845e3aSMatthew G. Knepley 500a0845e3aSMatthew G. Knepley Input Parameters: 501a0845e3aSMatthew G. Knepley + dm - The mesh 502a0845e3aSMatthew G. Knepley . X - Local input vector 503a0845e3aSMatthew G. Knepley - user - The user context 504a0845e3aSMatthew G. Knepley 505a0845e3aSMatthew G. Knepley Output Parameter: 506a0845e3aSMatthew G. Knepley . F - Local output vector 507a0845e3aSMatthew G. Knepley 508a0845e3aSMatthew G. Knepley Note: 509a0845e3aSMatthew G. Knepley The second member of the user context must be an FEMContext. 510a0845e3aSMatthew G. Knepley 511a0845e3aSMatthew G. Knepley We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 512a0845e3aSMatthew G. Knepley like a GPU, or vectorize on a multicore machine. 513a0845e3aSMatthew G. Knepley 514a0845e3aSMatthew G. Knepley Level: developer 515a0845e3aSMatthew G. Knepley 516a0845e3aSMatthew G. Knepley .seealso: DMPlexComputeJacobianActionFEM() 517a0845e3aSMatthew G. Knepley @*/ 518a0845e3aSMatthew G. Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 519a0845e3aSMatthew G. Knepley { 520a0845e3aSMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 521a0845e3aSMatthew G. Knepley PetscFEM *fem = (PetscFEM*) &((DM*) user)[1]; 522a0845e3aSMatthew G. Knepley PetscFE *fe = fem->fe; 523a0845e3aSMatthew G. Knepley const char *name = "Residual"; 524a0845e3aSMatthew G. Knepley PetscQuadrature q; 525a0845e3aSMatthew G. Knepley PetscCellGeometry geom; 526a0845e3aSMatthew G. Knepley PetscSection section; 527a0845e3aSMatthew G. Knepley PetscReal *v0, *J, *invJ, *detJ; 528a0845e3aSMatthew G. Knepley PetscScalar *elemVec, *u; 529a0845e3aSMatthew G. Knepley PetscInt dim, numFields, f, numCells, cStart, cEnd, c; 530a0845e3aSMatthew G. Knepley PetscInt cellDof = 0, numComponents = 0; 531a0845e3aSMatthew G. Knepley PetscErrorCode ierr; 532a0845e3aSMatthew G. Knepley 533a0845e3aSMatthew G. Knepley PetscFunctionBegin; 534a0845e3aSMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 535a0845e3aSMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 536a0845e3aSMatthew G. Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 537a0845e3aSMatthew G. Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 538a0845e3aSMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 539a0845e3aSMatthew G. Knepley numCells = cEnd - cStart; 540a0845e3aSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 541a0845e3aSMatthew G. Knepley PetscInt Nb, Nc; 542a0845e3aSMatthew G. Knepley 543a0845e3aSMatthew G. Knepley ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 544a0845e3aSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 545a0845e3aSMatthew G. Knepley cellDof += Nb*Nc; 546a0845e3aSMatthew G. Knepley numComponents += Nc; 547a0845e3aSMatthew G. Knepley } 548a0845e3aSMatthew G. Knepley ierr = DMPlexProjectFunctionLocal(dm, numComponents, fem->bcFuncs, INSERT_BC_VALUES, X);CHKERRQ(ierr); 549a0845e3aSMatthew G. Knepley ierr = VecSet(F, 0.0);CHKERRQ(ierr); 550a0845e3aSMatthew 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); 551a0845e3aSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 552a0845e3aSMatthew G. Knepley PetscScalar *x = NULL; 553a0845e3aSMatthew G. Knepley PetscInt i; 554a0845e3aSMatthew G. Knepley 555a0845e3aSMatthew G. Knepley ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 556a0845e3aSMatthew G. Knepley if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 557a0845e3aSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 558a0845e3aSMatthew G. Knepley for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 559a0845e3aSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 560a0845e3aSMatthew G. Knepley } 561a0845e3aSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 562a0845e3aSMatthew G. Knepley void (*f0)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[]) = fem->f0Funcs[f]; 563a0845e3aSMatthew G. Knepley void (*f1)(const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscScalar[]) = fem->f1Funcs[f]; 564a0845e3aSMatthew G. Knepley PetscInt Nb; 565a0845e3aSMatthew G. Knepley /* Conforming batches */ 566a0845e3aSMatthew G. Knepley PetscInt numBlocks = 1; 567a0845e3aSMatthew G. Knepley PetscInt numBatches = 1; 568a0845e3aSMatthew G. Knepley PetscInt numChunks, Ne, blockSize, batchSize; 569a0845e3aSMatthew G. Knepley /* Remainder */ 570a0845e3aSMatthew G. Knepley PetscInt Nr, offset; 571a0845e3aSMatthew G. Knepley 572a0845e3aSMatthew G. Knepley ierr = PetscFEGetQuadrature(fe[f], &q);CHKERRQ(ierr); 573a0845e3aSMatthew G. Knepley ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 574a0845e3aSMatthew G. Knepley blockSize = Nb*q.numQuadPoints; 575a0845e3aSMatthew G. Knepley batchSize = numBlocks * blockSize; 576a0845e3aSMatthew G. Knepley numChunks = numCells / (numBatches*batchSize); 577a0845e3aSMatthew G. Knepley Ne = numChunks*numBatches*batchSize; 578a0845e3aSMatthew G. Knepley Nr = numCells % (numBatches*batchSize); 579a0845e3aSMatthew G. Knepley offset = numCells - Nr; 580a0845e3aSMatthew G. Knepley geom.v0 = v0; 581a0845e3aSMatthew G. Knepley geom.J = J; 582a0845e3aSMatthew G. Knepley geom.invJ = invJ; 583a0845e3aSMatthew G. Knepley geom.detJ = detJ; 5840483ade4SMatthew G. Knepley ierr = PetscFEIntegrateResidual(fe[f], Ne, numFields, fe, f, geom, u, f0, f1, elemVec);CHKERRQ(ierr); 585a0845e3aSMatthew G. Knepley geom.v0 = &v0[offset*dim]; 586a0845e3aSMatthew G. Knepley geom.J = &J[offset*dim*dim]; 587a0845e3aSMatthew G. Knepley geom.invJ = &invJ[offset*dim*dim]; 588a0845e3aSMatthew G. Knepley geom.detJ = &detJ[offset]; 5890483ade4SMatthew G. Knepley ierr = PetscFEIntegrateResidual(fe[f], Nr, numFields, fe, f, geom, &u[offset*cellDof], f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 590a0845e3aSMatthew G. Knepley } 591a0845e3aSMatthew G. Knepley for (c = cStart; c < cEnd; ++c) { 592a0845e3aSMatthew G. Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, name, cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 593a0845e3aSMatthew G. Knepley ierr = DMPlexVecSetClosure(dm, section, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 594a0845e3aSMatthew G. Knepley } 595a0845e3aSMatthew G. Knepley ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 596a0845e3aSMatthew G. Knepley if (mesh->printFEM) {ierr = DMPrintLocalVec(dm, name, F);CHKERRQ(ierr);} 597a0845e3aSMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_ResidualFEM,dm,0,0,0);CHKERRQ(ierr); 598a0845e3aSMatthew G. Knepley PetscFunctionReturn(0); 599a0845e3aSMatthew G. Knepley } 600a0845e3aSMatthew G. Knepley 601a0845e3aSMatthew G. Knepley #endif 602a0845e3aSMatthew G. Knepley 603cb1e1211SMatthew G Knepley #undef __FUNCT__ 604cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianActionFEM" 605cb1e1211SMatthew G Knepley /*@C 606cb1e1211SMatthew G Knepley DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user 607cb1e1211SMatthew G Knepley 608cb1e1211SMatthew G Knepley Input Parameters: 609cb1e1211SMatthew G Knepley + dm - The mesh 610cb1e1211SMatthew G Knepley . J - The Jacobian shell matrix 611cb1e1211SMatthew G Knepley . X - Local input vector 612cb1e1211SMatthew G Knepley - user - The user context 613cb1e1211SMatthew G Knepley 614cb1e1211SMatthew G Knepley Output Parameter: 615cb1e1211SMatthew G Knepley . F - Local output vector 616cb1e1211SMatthew G Knepley 617cb1e1211SMatthew G Knepley Note: 618cb1e1211SMatthew G Knepley The second member of the user context must be an FEMContext. 619cb1e1211SMatthew G Knepley 620cb1e1211SMatthew G Knepley We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 621cb1e1211SMatthew G Knepley like a GPU, or vectorize on a multicore machine. 622cb1e1211SMatthew G Knepley 6230059ad2aSSatish Balay Level: developer 6240059ad2aSSatish Balay 625cb1e1211SMatthew G Knepley .seealso: DMPlexComputeResidualFEM() 626878cb397SSatish Balay @*/ 627cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user) 628cb1e1211SMatthew G Knepley { 629cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 630cb1e1211SMatthew G Knepley PetscFEM *fem = (PetscFEM*) &((DM*) user)[1]; 6310483ade4SMatthew G. Knepley PetscFE *fe = fem->fe; 6320483ade4SMatthew G. Knepley PetscQuadrature quad; 6330483ade4SMatthew G. Knepley PetscCellGeometry geom; 634cb1e1211SMatthew G Knepley PetscSection section; 635cb1e1211SMatthew G Knepley JacActionCtx *jctx; 636cb1e1211SMatthew G Knepley PetscReal *v0, *J, *invJ, *detJ; 637cb1e1211SMatthew G Knepley PetscScalar *elemVec, *u, *a; 6380483ade4SMatthew G. Knepley PetscInt dim, numFields, field, numCells, cStart, cEnd, c; 639cb1e1211SMatthew G Knepley PetscInt cellDof = 0; 640cb1e1211SMatthew G Knepley PetscErrorCode ierr; 641cb1e1211SMatthew G Knepley 642cb1e1211SMatthew G Knepley PetscFunctionBegin; 6430483ade4SMatthew G. Knepley /* ierr = PetscLogEventBegin(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */ 644cb1e1211SMatthew G Knepley ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 645cb1e1211SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 646cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 647cb1e1211SMatthew G Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 648cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 649cb1e1211SMatthew G Knepley numCells = cEnd - cStart; 650cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 6510483ade4SMatthew G. Knepley PetscInt Nb, Nc; 6520483ade4SMatthew G. Knepley 6530483ade4SMatthew G. Knepley ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 6540483ade4SMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[field], &Nc);CHKERRQ(ierr); 6550483ade4SMatthew G. Knepley cellDof += Nb*Nc; 656cb1e1211SMatthew G Knepley } 657cb1e1211SMatthew G Knepley ierr = VecSet(F, 0.0);CHKERRQ(ierr); 658cb1e1211SMatthew 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); 659cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 660cb1e1211SMatthew G Knepley PetscScalar *x; 661cb1e1211SMatthew G Knepley PetscInt i; 662cb1e1211SMatthew G Knepley 663cb1e1211SMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 664cb1e1211SMatthew G Knepley if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 665cb1e1211SMatthew G Knepley ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 666cb1e1211SMatthew G Knepley for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 667cb1e1211SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 668cb1e1211SMatthew G Knepley ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 669cb1e1211SMatthew G Knepley for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i]; 670cb1e1211SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 671cb1e1211SMatthew G Knepley } 672cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 6730483ade4SMatthew G. Knepley PetscInt Nb; 674cb1e1211SMatthew G Knepley /* Conforming batches */ 675cb1e1211SMatthew G Knepley PetscInt numBlocks = 1; 6760483ade4SMatthew G. Knepley PetscInt numBatches = 1; 6770483ade4SMatthew G. Knepley PetscInt numChunks, Ne, blockSize, batchSize; 678cb1e1211SMatthew G Knepley /* Remainder */ 6790483ade4SMatthew G. Knepley PetscInt Nr, offset; 680cb1e1211SMatthew G Knepley 6810483ade4SMatthew G. Knepley ierr = PetscFEGetQuadrature(fe[field], &quad);CHKERRQ(ierr); 6820483ade4SMatthew G. Knepley ierr = PetscFEGetDimension(fe[field], &Nb);CHKERRQ(ierr); 6830483ade4SMatthew G. Knepley blockSize = Nb*quad.numQuadPoints; 6840483ade4SMatthew G. Knepley batchSize = numBlocks * blockSize; 6850483ade4SMatthew G. Knepley numChunks = numCells / (numBatches*batchSize); 6860483ade4SMatthew G. Knepley Ne = numChunks*numBatches*batchSize; 6870483ade4SMatthew G. Knepley Nr = numCells % (numBatches*batchSize); 6880483ade4SMatthew G. Knepley offset = numCells - Nr; 6890483ade4SMatthew G. Knepley geom.v0 = v0; 6900483ade4SMatthew G. Knepley geom.J = J; 6910483ade4SMatthew G. Knepley geom.invJ = invJ; 6920483ade4SMatthew G. Knepley geom.detJ = detJ; 6930483ade4SMatthew G. Knepley ierr = PetscFEIntegrateJacobianAction(fe[field], Ne, numFields, fe, field, geom, u, a, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr); 6940483ade4SMatthew G. Knepley geom.v0 = &v0[offset*dim]; 6950483ade4SMatthew G. Knepley geom.J = &J[offset*dim*dim]; 6960483ade4SMatthew G. Knepley geom.invJ = &invJ[offset*dim*dim]; 6970483ade4SMatthew G. Knepley geom.detJ = &detJ[offset]; 6980483ade4SMatthew G. Knepley ierr = PetscFEIntegrateJacobianAction(fe[field], Nr, numFields, fe, field, geom, &u[offset*cellDof], &a[offset*cellDof], 699cb1e1211SMatthew G Knepley fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr); 700cb1e1211SMatthew G Knepley } 701cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 702cb1e1211SMatthew G Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 703cb1e1211SMatthew G Knepley ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 704cb1e1211SMatthew G Knepley } 705cb1e1211SMatthew G Knepley ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 706cb1e1211SMatthew G Knepley if (mesh->printFEM) { 707cb1e1211SMatthew G Knepley PetscMPIInt rank, numProcs; 708cb1e1211SMatthew G Knepley PetscInt p; 709cb1e1211SMatthew G Knepley 710cb1e1211SMatthew G Knepley ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 711cb1e1211SMatthew G Knepley ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr); 71286a74ee0SMatthew G. Knepley ierr = PetscPrintf(PetscObjectComm((PetscObject)dm), "Jacobian Action:\n");CHKERRQ(ierr); 713cb1e1211SMatthew G Knepley for (p = 0; p < numProcs; ++p) { 714cb1e1211SMatthew G Knepley if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);} 715cb1e1211SMatthew G Knepley ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr); 716cb1e1211SMatthew G Knepley } 717cb1e1211SMatthew G Knepley } 7180483ade4SMatthew G. Knepley /* ierr = PetscLogEventEnd(DMPLEX_JacobianActionFEM,dm,0,0,0);CHKERRQ(ierr); */ 719cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 720cb1e1211SMatthew G Knepley } 721cb1e1211SMatthew G Knepley 722cb1e1211SMatthew G Knepley #undef __FUNCT__ 723cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianFEM" 724cb1e1211SMatthew G Knepley /*@ 725cb1e1211SMatthew G Knepley DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user. 726cb1e1211SMatthew G Knepley 727cb1e1211SMatthew G Knepley Input Parameters: 728cb1e1211SMatthew G Knepley + dm - The mesh 729cb1e1211SMatthew G Knepley . X - Local input vector 730cb1e1211SMatthew G Knepley - user - The user context 731cb1e1211SMatthew G Knepley 732cb1e1211SMatthew G Knepley Output Parameter: 733cb1e1211SMatthew G Knepley . Jac - Jacobian matrix 734cb1e1211SMatthew G Knepley 735cb1e1211SMatthew G Knepley Note: 736cb1e1211SMatthew G Knepley The second member of the user context must be an FEMContext. 737cb1e1211SMatthew G Knepley 738cb1e1211SMatthew G Knepley We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 739cb1e1211SMatthew G Knepley like a GPU, or vectorize on a multicore machine. 740cb1e1211SMatthew G Knepley 7410059ad2aSSatish Balay Level: developer 7420059ad2aSSatish Balay 743cb1e1211SMatthew G Knepley .seealso: FormFunctionLocal() 744878cb397SSatish Balay @*/ 745cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, MatStructure *str,void *user) 746cb1e1211SMatthew G Knepley { 747cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 748cb1e1211SMatthew G Knepley PetscFEM *fem = (PetscFEM*) &((DM*) user)[1]; 749a319912fSMatthew G. Knepley PetscFE *fe = fem->fe; 750a319912fSMatthew G. Knepley const char *name = "Jacobian"; 751a319912fSMatthew G. Knepley PetscQuadrature quad; 752a319912fSMatthew G. Knepley PetscCellGeometry geom; 753a319912fSMatthew G. Knepley PetscSection section, globalSection; 754cb1e1211SMatthew G Knepley PetscReal *v0, *J, *invJ, *detJ; 755cb1e1211SMatthew G Knepley PetscScalar *elemMat, *u; 756a319912fSMatthew G. Knepley PetscInt dim, numFields, f, fieldI, fieldJ, numCells, cStart, cEnd, c; 757cb1e1211SMatthew G Knepley PetscInt cellDof = 0, numComponents = 0; 758cb1e1211SMatthew G Knepley PetscBool isShell; 759cb1e1211SMatthew G Knepley PetscErrorCode ierr; 760cb1e1211SMatthew G Knepley 761cb1e1211SMatthew G Knepley PetscFunctionBegin; 762a319912fSMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 763cb1e1211SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 764cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 765a319912fSMatthew G. Knepley ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr); 766cb1e1211SMatthew G Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 767cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 768cb1e1211SMatthew G Knepley numCells = cEnd - cStart; 769a319912fSMatthew G. Knepley for (f = 0; f < numFields; ++f) { 770a319912fSMatthew G. Knepley PetscInt Nb, Nc; 771a319912fSMatthew G. Knepley 772a319912fSMatthew G. Knepley ierr = PetscFEGetDimension(fe[f], &Nb);CHKERRQ(ierr); 773a319912fSMatthew G. Knepley ierr = PetscFEGetNumComponents(fe[f], &Nc);CHKERRQ(ierr); 774a319912fSMatthew G. Knepley cellDof += Nb*Nc; 775a319912fSMatthew G. Knepley numComponents += Nc; 776cb1e1211SMatthew G Knepley } 777cb1e1211SMatthew G Knepley ierr = DMPlexProjectFunctionLocal(dm, numComponents, fem->bcFuncs, INSERT_BC_VALUES, X);CHKERRQ(ierr); 778cb1e1211SMatthew G Knepley ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 779cb1e1211SMatthew 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); 780cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 781cb1e1211SMatthew G Knepley PetscScalar *x; 782cb1e1211SMatthew G Knepley PetscInt i; 783cb1e1211SMatthew G Knepley 784cb1e1211SMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 785cb1e1211SMatthew G Knepley if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 786a319912fSMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 787cb1e1211SMatthew G Knepley for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 788a319912fSMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, section, X, c, NULL, &x);CHKERRQ(ierr); 789cb1e1211SMatthew G Knepley } 790cb1e1211SMatthew G Knepley ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr); 791cb1e1211SMatthew G Knepley for (fieldI = 0; fieldI < numFields; ++fieldI) { 792a319912fSMatthew G. Knepley PetscInt Nb; 793a319912fSMatthew G. Knepley ierr = PetscFEGetQuadrature(fe[fieldI], &quad);CHKERRQ(ierr); 794a319912fSMatthew G. Knepley ierr = PetscFEGetDimension(fe[fieldI], &Nb);CHKERRQ(ierr); 795cb1e1211SMatthew G Knepley for (fieldJ = 0; fieldJ < numFields; ++fieldJ) { 796a0845e3aSMatthew G. Knepley void (*g0)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g0Funcs[fieldI*numFields+fieldJ]; 797a0845e3aSMatthew G. Knepley void (*g1)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g1Funcs[fieldI*numFields+fieldJ]; 798a0845e3aSMatthew G. Knepley void (*g2)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g2Funcs[fieldI*numFields+fieldJ]; 799a0845e3aSMatthew G. Knepley void (*g3)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]) = fem->g3Funcs[fieldI*numFields+fieldJ]; 800cb1e1211SMatthew G Knepley /* Conforming batches */ 801cb1e1211SMatthew G Knepley PetscInt numBlocks = 1; 802a319912fSMatthew G. Knepley PetscInt numBatches = 1; 803a319912fSMatthew G. Knepley PetscInt numChunks, Ne, blockSize, batchSize; 804cb1e1211SMatthew G Knepley /* Remainder */ 805a319912fSMatthew G. Knepley PetscInt Nr, offset; 806cb1e1211SMatthew G Knepley 807a319912fSMatthew G. Knepley blockSize = Nb*quad.numQuadPoints; 808a319912fSMatthew G. Knepley batchSize = numBlocks * blockSize; 809a319912fSMatthew G. Knepley numChunks = numCells / (numBatches*batchSize); 810a319912fSMatthew G. Knepley Ne = numChunks*numBatches*batchSize; 811a319912fSMatthew G. Knepley Nr = numCells % (numBatches*batchSize); 812a319912fSMatthew G. Knepley offset = numCells - Nr; 813a319912fSMatthew G. Knepley geom.v0 = v0; 814a319912fSMatthew G. Knepley geom.J = J; 815a319912fSMatthew G. Knepley geom.invJ = invJ; 816a319912fSMatthew G. Knepley geom.detJ = detJ; 8170483ade4SMatthew G. Knepley ierr = PetscFEIntegrateJacobian(fe[fieldI], Ne, numFields, fe, fieldI, fieldJ, geom, u, g0, g1, g2, g3, elemMat);CHKERRQ(ierr); 818a319912fSMatthew G. Knepley geom.v0 = &v0[offset*dim]; 819a319912fSMatthew G. Knepley geom.J = &J[offset*dim*dim]; 820a319912fSMatthew G. Knepley geom.invJ = &invJ[offset*dim*dim]; 821a319912fSMatthew G. Knepley geom.detJ = &detJ[offset]; 8220483ade4SMatthew G. Knepley ierr = PetscFEIntegrateJacobian(fe[fieldI], Nr, numFields, fe, fieldI, fieldJ, geom, &u[offset*cellDof], g0, g1, g2, g3, &elemMat[offset*cellDof*cellDof]);CHKERRQ(ierr); 823cb1e1211SMatthew G Knepley } 824cb1e1211SMatthew G Knepley } 825cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 826a319912fSMatthew G. Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, name, cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);} 827a319912fSMatthew G. Knepley ierr = DMPlexMatSetClosure(dm, section, globalSection, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr); 828cb1e1211SMatthew G Knepley } 829cb1e1211SMatthew G Knepley ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr); 830cb1e1211SMatthew G Knepley ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 831cb1e1211SMatthew G Knepley ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 832cb1e1211SMatthew G Knepley if (mesh->printFEM) { 833a319912fSMatthew G. Knepley ierr = PetscPrintf(PETSC_COMM_WORLD, "%s:\n", name);CHKERRQ(ierr); 834cb1e1211SMatthew G Knepley ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr); 835cb1e1211SMatthew G Knepley ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 836cb1e1211SMatthew G Knepley } 837a319912fSMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_JacobianFEM,dm,0,0,0);CHKERRQ(ierr); 838cb1e1211SMatthew G Knepley ierr = PetscObjectTypeCompare((PetscObject) Jac, MATSHELL, &isShell);CHKERRQ(ierr); 839cb1e1211SMatthew G Knepley if (isShell) { 840cb1e1211SMatthew G Knepley JacActionCtx *jctx; 841cb1e1211SMatthew G Knepley 842cb1e1211SMatthew G Knepley ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 843cb1e1211SMatthew G Knepley ierr = VecCopy(X, jctx->u);CHKERRQ(ierr); 844cb1e1211SMatthew G Knepley } 845cb1e1211SMatthew G Knepley *str = SAME_NONZERO_PATTERN; 846cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 847cb1e1211SMatthew G Knepley } 848