1*cb1e1211SMatthew G Knepley #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2*cb1e1211SMatthew G Knepley 3*cb1e1211SMatthew G Knepley #undef __FUNCT__ 4*cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexGetScale" 5*cb1e1211SMatthew G Knepley PetscErrorCode DMPlexGetScale(DM dm, PetscUnit unit, PetscReal *scale) 6*cb1e1211SMatthew G Knepley { 7*cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 8*cb1e1211SMatthew G Knepley 9*cb1e1211SMatthew G Knepley PetscFunctionBegin; 10*cb1e1211SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11*cb1e1211SMatthew G Knepley PetscValidPointer(scale, 3); 12*cb1e1211SMatthew G Knepley *scale = mesh->scale[unit]; 13*cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 14*cb1e1211SMatthew G Knepley } 15*cb1e1211SMatthew G Knepley 16*cb1e1211SMatthew G Knepley #undef __FUNCT__ 17*cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexSetScale" 18*cb1e1211SMatthew G Knepley PetscErrorCode DMPlexSetScale(DM dm, PetscUnit unit, PetscReal scale) 19*cb1e1211SMatthew G Knepley { 20*cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 21*cb1e1211SMatthew G Knepley 22*cb1e1211SMatthew G Knepley PetscFunctionBegin; 23*cb1e1211SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 24*cb1e1211SMatthew G Knepley mesh->scale[unit] = scale; 25*cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 26*cb1e1211SMatthew G Knepley } 27*cb1e1211SMatthew G Knepley 28*cb1e1211SMatthew G Knepley PETSC_STATIC_INLINE PetscInt epsilon(PetscInt i, PetscInt j, PetscInt k) 29*cb1e1211SMatthew G Knepley { 30*cb1e1211SMatthew G Knepley switch (i) { 31*cb1e1211SMatthew G Knepley case 0: 32*cb1e1211SMatthew G Knepley switch (j) { 33*cb1e1211SMatthew G Knepley case 0: return 0; 34*cb1e1211SMatthew G Knepley case 1: 35*cb1e1211SMatthew G Knepley switch (k) { 36*cb1e1211SMatthew G Knepley case 0: return 0; 37*cb1e1211SMatthew G Knepley case 1: return 0; 38*cb1e1211SMatthew G Knepley case 2: return 1; 39*cb1e1211SMatthew G Knepley } 40*cb1e1211SMatthew G Knepley case 2: 41*cb1e1211SMatthew G Knepley switch (k) { 42*cb1e1211SMatthew G Knepley case 0: return 0; 43*cb1e1211SMatthew G Knepley case 1: return -1; 44*cb1e1211SMatthew G Knepley case 2: return 0; 45*cb1e1211SMatthew G Knepley } 46*cb1e1211SMatthew G Knepley } 47*cb1e1211SMatthew G Knepley case 1: 48*cb1e1211SMatthew G Knepley switch (j) { 49*cb1e1211SMatthew G Knepley case 0: 50*cb1e1211SMatthew G Knepley switch (k) { 51*cb1e1211SMatthew G Knepley case 0: return 0; 52*cb1e1211SMatthew G Knepley case 1: return 0; 53*cb1e1211SMatthew G Knepley case 2: return -1; 54*cb1e1211SMatthew G Knepley } 55*cb1e1211SMatthew G Knepley case 1: return 0; 56*cb1e1211SMatthew G Knepley case 2: 57*cb1e1211SMatthew G Knepley switch (k) { 58*cb1e1211SMatthew G Knepley case 0: return 1; 59*cb1e1211SMatthew G Knepley case 1: return 0; 60*cb1e1211SMatthew G Knepley case 2: return 0; 61*cb1e1211SMatthew G Knepley } 62*cb1e1211SMatthew G Knepley } 63*cb1e1211SMatthew G Knepley case 2: 64*cb1e1211SMatthew G Knepley switch (j) { 65*cb1e1211SMatthew G Knepley case 0: 66*cb1e1211SMatthew G Knepley switch (k) { 67*cb1e1211SMatthew G Knepley case 0: return 0; 68*cb1e1211SMatthew G Knepley case 1: return 1; 69*cb1e1211SMatthew G Knepley case 2: return 0; 70*cb1e1211SMatthew G Knepley } 71*cb1e1211SMatthew G Knepley case 1: 72*cb1e1211SMatthew G Knepley switch (k) { 73*cb1e1211SMatthew G Knepley case 0: return -1; 74*cb1e1211SMatthew G Knepley case 1: return 0; 75*cb1e1211SMatthew G Knepley case 2: return 0; 76*cb1e1211SMatthew G Knepley } 77*cb1e1211SMatthew G Knepley case 2: return 0; 78*cb1e1211SMatthew G Knepley } 79*cb1e1211SMatthew G Knepley } 80*cb1e1211SMatthew G Knepley return 0; 81*cb1e1211SMatthew G Knepley } 82*cb1e1211SMatthew G Knepley 83*cb1e1211SMatthew G Knepley #undef __FUNCT__ 84*cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexCreateRigidBody" 85*cb1e1211SMatthew G Knepley /*@C 86*cb1e1211SMatthew G Knepley DMPlexCreateRigidBody - create rigid body modes from coordinates 87*cb1e1211SMatthew G Knepley 88*cb1e1211SMatthew G Knepley Collective on DM 89*cb1e1211SMatthew G Knepley 90*cb1e1211SMatthew G Knepley Input Arguments: 91*cb1e1211SMatthew G Knepley + dm - the DM 92*cb1e1211SMatthew G Knepley . section - the local section associated with the rigid field, or NULL for the default section 93*cb1e1211SMatthew G Knepley - globalSection - the global section associated with the rigid field, or NULL for the default section 94*cb1e1211SMatthew G Knepley 95*cb1e1211SMatthew G Knepley Output Argument: 96*cb1e1211SMatthew G Knepley . sp - the null space 97*cb1e1211SMatthew G Knepley 98*cb1e1211SMatthew G Knepley Note: This is necessary to take account of Dirichlet conditions on the displacements 99*cb1e1211SMatthew G Knepley 100*cb1e1211SMatthew G Knepley Level: advanced 101*cb1e1211SMatthew G Knepley 102*cb1e1211SMatthew G Knepley .seealso: MatNullSpaceCreate() 103*cb1e1211SMatthew G Knepley @*/ 104*cb1e1211SMatthew G Knepley PetscErrorCode DMPlexCreateRigidBody(DM dm, PetscSection section, PetscSection globalSection, MatNullSpace *sp) 105*cb1e1211SMatthew G Knepley { 106*cb1e1211SMatthew G Knepley MPI_Comm comm; 107*cb1e1211SMatthew G Knepley Vec coordinates, localMode, mode[6]; 108*cb1e1211SMatthew G Knepley PetscSection coordSection; 109*cb1e1211SMatthew G Knepley PetscScalar *coords; 110*cb1e1211SMatthew G Knepley PetscInt dim, vStart, vEnd, v, n, m, d, i, j; 111*cb1e1211SMatthew G Knepley PetscErrorCode ierr; 112*cb1e1211SMatthew G Knepley 113*cb1e1211SMatthew G Knepley PetscFunctionBegin; 114*cb1e1211SMatthew G Knepley ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 115*cb1e1211SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 116*cb1e1211SMatthew G Knepley if (dim == 1) { 117*cb1e1211SMatthew G Knepley ierr = MatNullSpaceCreate(comm, PETSC_TRUE, 0, NULL, sp);CHKERRQ(ierr); 118*cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 119*cb1e1211SMatthew G Knepley } 120*cb1e1211SMatthew G Knepley if (!section) {ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr);} 121*cb1e1211SMatthew G Knepley if (!globalSection) {ierr = DMGetDefaultGlobalSection(dm, &globalSection);CHKERRQ(ierr);} 122*cb1e1211SMatthew G Knepley ierr = PetscSectionGetConstrainedStorageSize(globalSection, &n);CHKERRQ(ierr); 123*cb1e1211SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 124*cb1e1211SMatthew G Knepley ierr = DMPlexGetCoordinateSection(dm, &coordSection);CHKERRQ(ierr); 125*cb1e1211SMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 126*cb1e1211SMatthew G Knepley m = (dim*(dim+1))/2; 127*cb1e1211SMatthew G Knepley ierr = VecCreate(comm, &mode[0]);CHKERRQ(ierr); 128*cb1e1211SMatthew G Knepley ierr = VecSetSizes(mode[0], n, PETSC_DETERMINE);CHKERRQ(ierr); 129*cb1e1211SMatthew G Knepley ierr = VecSetUp(mode[0]);CHKERRQ(ierr); 130*cb1e1211SMatthew G Knepley for (i = 1; i < m; ++i) {ierr = VecDuplicate(mode[0], &mode[i]);CHKERRQ(ierr);} 131*cb1e1211SMatthew G Knepley /* Assume P1 */ 132*cb1e1211SMatthew G Knepley ierr = DMGetLocalVector(dm, &localMode);CHKERRQ(ierr); 133*cb1e1211SMatthew G Knepley for (d = 0; d < dim; ++d) { 134*cb1e1211SMatthew G Knepley PetscScalar values[3] = {0.0, 0.0, 0.0}; 135*cb1e1211SMatthew G Knepley 136*cb1e1211SMatthew G Knepley values[d] = 1.0; 137*cb1e1211SMatthew G Knepley ierr = VecSet(localMode, 0.0);CHKERRQ(ierr); 138*cb1e1211SMatthew G Knepley for (v = vStart; v < vEnd; ++v) { 139*cb1e1211SMatthew G Knepley ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr); 140*cb1e1211SMatthew G Knepley } 141*cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 142*cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 143*cb1e1211SMatthew G Knepley } 144*cb1e1211SMatthew G Knepley ierr = VecGetArray(coordinates, &coords);CHKERRQ(ierr); 145*cb1e1211SMatthew G Knepley for (d = dim; d < dim*(dim+1)/2; ++d) { 146*cb1e1211SMatthew G Knepley PetscInt i, j, k = dim > 2 ? d - dim : d; 147*cb1e1211SMatthew G Knepley 148*cb1e1211SMatthew G Knepley ierr = VecSet(localMode, 0.0);CHKERRQ(ierr); 149*cb1e1211SMatthew G Knepley for (v = vStart; v < vEnd; ++v) { 150*cb1e1211SMatthew G Knepley PetscScalar values[3] = {0.0, 0.0, 0.0}; 151*cb1e1211SMatthew G Knepley PetscInt off; 152*cb1e1211SMatthew G Knepley 153*cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 154*cb1e1211SMatthew G Knepley for (i = 0; i < dim; ++i) { 155*cb1e1211SMatthew G Knepley for (j = 0; j < dim; ++j) { 156*cb1e1211SMatthew G Knepley values[j] += epsilon(i, j, k)*PetscRealPart(coords[off+i]); 157*cb1e1211SMatthew G Knepley } 158*cb1e1211SMatthew G Knepley } 159*cb1e1211SMatthew G Knepley ierr = DMPlexVecSetClosure(dm, section, localMode, v, values, INSERT_VALUES);CHKERRQ(ierr); 160*cb1e1211SMatthew G Knepley } 161*cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalBegin(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 162*cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalEnd(dm, localMode, INSERT_VALUES, mode[d]);CHKERRQ(ierr); 163*cb1e1211SMatthew G Knepley } 164*cb1e1211SMatthew G Knepley ierr = VecRestoreArray(coordinates, &coords);CHKERRQ(ierr); 165*cb1e1211SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localMode);CHKERRQ(ierr); 166*cb1e1211SMatthew G Knepley for (i = 0; i < dim; ++i) {ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr);} 167*cb1e1211SMatthew G Knepley /* Orthonormalize system */ 168*cb1e1211SMatthew G Knepley for (i = dim; i < m; ++i) { 169*cb1e1211SMatthew G Knepley PetscScalar dots[6]; 170*cb1e1211SMatthew G Knepley 171*cb1e1211SMatthew G Knepley ierr = VecMDot(mode[i], i, mode, dots);CHKERRQ(ierr); 172*cb1e1211SMatthew G Knepley for (j = 0; j < i; ++j) dots[j] *= -1.0; 173*cb1e1211SMatthew G Knepley ierr = VecMAXPY(mode[i], i, dots, mode);CHKERRQ(ierr); 174*cb1e1211SMatthew G Knepley ierr = VecNormalize(mode[i], NULL);CHKERRQ(ierr); 175*cb1e1211SMatthew G Knepley } 176*cb1e1211SMatthew G Knepley ierr = MatNullSpaceCreate(comm, PETSC_FALSE, m, mode, sp);CHKERRQ(ierr); 177*cb1e1211SMatthew G Knepley for (i = 0; i< m; ++i) {ierr = VecDestroy(&mode[i]);CHKERRQ(ierr);} 178*cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 179*cb1e1211SMatthew G Knepley } 180*cb1e1211SMatthew G Knepley /******************************************************************************* 181*cb1e1211SMatthew G Knepley This should be in a separate Discretization object, but I am not sure how to lay 182*cb1e1211SMatthew G Knepley it out yet, so I am stuffing things here while I experiment. 183*cb1e1211SMatthew G Knepley *******************************************************************************/ 184*cb1e1211SMatthew G Knepley #undef __FUNCT__ 185*cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexSetFEMIntegration" 186*cb1e1211SMatthew G Knepley PetscErrorCode DMPlexSetFEMIntegration(DM dm, 187*cb1e1211SMatthew G Knepley PetscErrorCode (*integrateResidualFEM)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], 188*cb1e1211SMatthew G Knepley const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], 189*cb1e1211SMatthew G Knepley void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 190*cb1e1211SMatthew G Knepley void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]), 191*cb1e1211SMatthew G Knepley PetscErrorCode (*integrateJacobianActionFEM)(PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], const PetscScalar[], 192*cb1e1211SMatthew G Knepley const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], 193*cb1e1211SMatthew G Knepley void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 194*cb1e1211SMatthew G Knepley void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 195*cb1e1211SMatthew G Knepley void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 196*cb1e1211SMatthew G Knepley void (**)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[]), 197*cb1e1211SMatthew G Knepley PetscErrorCode (*integrateJacobianFEM)(PetscInt, PetscInt, PetscInt, PetscInt, PetscQuadrature[], const PetscScalar[], 198*cb1e1211SMatthew G Knepley const PetscReal[], const PetscReal[], const PetscReal[], const PetscReal[], 199*cb1e1211SMatthew G Knepley void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 200*cb1e1211SMatthew G Knepley void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 201*cb1e1211SMatthew G Knepley void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 202*cb1e1211SMatthew G Knepley void (*)(const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), PetscScalar[])) 203*cb1e1211SMatthew G Knepley { 204*cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 205*cb1e1211SMatthew G Knepley 206*cb1e1211SMatthew G Knepley PetscFunctionBegin; 207*cb1e1211SMatthew G Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 208*cb1e1211SMatthew G Knepley mesh->integrateResidualFEM = integrateResidualFEM; 209*cb1e1211SMatthew G Knepley mesh->integrateJacobianActionFEM = integrateJacobianActionFEM; 210*cb1e1211SMatthew G Knepley mesh->integrateJacobianFEM = integrateJacobianFEM; 211*cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 212*cb1e1211SMatthew G Knepley } 213*cb1e1211SMatthew G Knepley 214*cb1e1211SMatthew G Knepley #undef __FUNCT__ 215*cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunctionLocal" 216*cb1e1211SMatthew G Knepley PetscErrorCode DMPlexProjectFunctionLocal(DM dm, PetscInt numComp, PetscScalar (**funcs)(const PetscReal []), InsertMode mode, Vec localX) 217*cb1e1211SMatthew G Knepley { 218*cb1e1211SMatthew G Knepley Vec coordinates; 219*cb1e1211SMatthew G Knepley PetscSection section, cSection; 220*cb1e1211SMatthew G Knepley PetscInt dim, vStart, vEnd, v, c, d; 221*cb1e1211SMatthew G Knepley PetscScalar *values, *cArray; 222*cb1e1211SMatthew G Knepley PetscReal *coords; 223*cb1e1211SMatthew G Knepley PetscErrorCode ierr; 224*cb1e1211SMatthew G Knepley 225*cb1e1211SMatthew G Knepley PetscFunctionBegin; 226*cb1e1211SMatthew G Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 227*cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 228*cb1e1211SMatthew G Knepley ierr = DMPlexGetCoordinateSection(dm, &cSection);CHKERRQ(ierr); 229*cb1e1211SMatthew G Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 230*cb1e1211SMatthew G Knepley ierr = PetscMalloc(numComp * sizeof(PetscScalar), &values);CHKERRQ(ierr); 231*cb1e1211SMatthew G Knepley ierr = VecGetArray(coordinates, &cArray);CHKERRQ(ierr); 232*cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(cSection, vStart, &dim);CHKERRQ(ierr); 233*cb1e1211SMatthew G Knepley ierr = PetscMalloc(dim * sizeof(PetscReal),&coords);CHKERRQ(ierr); 234*cb1e1211SMatthew G Knepley for (v = vStart; v < vEnd; ++v) { 235*cb1e1211SMatthew G Knepley PetscInt dof, off; 236*cb1e1211SMatthew G Knepley 237*cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(cSection, v, &dof);CHKERRQ(ierr); 238*cb1e1211SMatthew G Knepley ierr = PetscSectionGetOffset(cSection, v, &off);CHKERRQ(ierr); 239*cb1e1211SMatthew G Knepley if (dof > dim) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Cannot have more coordinates %d then dimensions %d", dof, dim); 240*cb1e1211SMatthew G Knepley for (d = 0; d < dof; ++d) coords[d] = PetscRealPart(cArray[off+d]); 241*cb1e1211SMatthew G Knepley for (c = 0; c < numComp; ++c) values[c] = (*funcs[c])(coords); 242*cb1e1211SMatthew G Knepley ierr = VecSetValuesSection(localX, section, v, values, mode);CHKERRQ(ierr); 243*cb1e1211SMatthew G Knepley } 244*cb1e1211SMatthew G Knepley ierr = VecRestoreArray(coordinates, &cArray);CHKERRQ(ierr); 245*cb1e1211SMatthew G Knepley /* Temporary, must be replaced by a projection on the finite element basis */ 246*cb1e1211SMatthew G Knepley { 247*cb1e1211SMatthew G Knepley PetscInt eStart = 0, eEnd = 0, e, depth; 248*cb1e1211SMatthew G Knepley 249*cb1e1211SMatthew G Knepley ierr = DMPlexGetLabelSize(dm, "depth", &depth);CHKERRQ(ierr); 250*cb1e1211SMatthew G Knepley --depth; 251*cb1e1211SMatthew G Knepley if (depth > 1) {ierr = DMPlexGetDepthStratum(dm, 1, &eStart, &eEnd);CHKERRQ(ierr);} 252*cb1e1211SMatthew G Knepley for (e = eStart; e < eEnd; ++e) { 253*cb1e1211SMatthew G Knepley const PetscInt *cone = NULL; 254*cb1e1211SMatthew G Knepley PetscInt coneSize, d; 255*cb1e1211SMatthew G Knepley PetscScalar *coordsA, *coordsB; 256*cb1e1211SMatthew G Knepley 257*cb1e1211SMatthew G Knepley ierr = DMPlexGetConeSize(dm, e, &coneSize);CHKERRQ(ierr); 258*cb1e1211SMatthew G Knepley ierr = DMPlexGetCone(dm, e, &cone);CHKERRQ(ierr); 259*cb1e1211SMatthew G Knepley if (coneSize != 2) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_SIZ, "Cone size %d for point %d should be 2", coneSize, e); 260*cb1e1211SMatthew G Knepley ierr = VecGetValuesSection(coordinates, cSection, cone[0], &coordsA);CHKERRQ(ierr); 261*cb1e1211SMatthew G Knepley ierr = VecGetValuesSection(coordinates, cSection, cone[1], &coordsB);CHKERRQ(ierr); 262*cb1e1211SMatthew G Knepley for (d = 0; d < dim; ++d) { 263*cb1e1211SMatthew G Knepley coords[d] = 0.5*(PetscRealPart(coordsA[d]) + PetscRealPart(coordsB[d])); 264*cb1e1211SMatthew G Knepley } 265*cb1e1211SMatthew G Knepley for (c = 0; c < numComp; ++c) values[c] = (*funcs[c])(coords); 266*cb1e1211SMatthew G Knepley ierr = VecSetValuesSection(localX, section, e, values, mode);CHKERRQ(ierr); 267*cb1e1211SMatthew G Knepley } 268*cb1e1211SMatthew G Knepley } 269*cb1e1211SMatthew G Knepley 270*cb1e1211SMatthew G Knepley ierr = PetscFree(coords);CHKERRQ(ierr); 271*cb1e1211SMatthew G Knepley ierr = PetscFree(values);CHKERRQ(ierr); 272*cb1e1211SMatthew G Knepley #if 0 273*cb1e1211SMatthew G Knepley const PetscInt localDof = this->_mesh->sizeWithBC(s, *cells->begin()); 274*cb1e1211SMatthew G Knepley PetscReal detJ; 275*cb1e1211SMatthew G Knepley 276*cb1e1211SMatthew G Knepley ierr = PetscMalloc(localDof * sizeof(PetscScalar), &values);CHKERRQ(ierr); 277*cb1e1211SMatthew G Knepley ierr = PetscMalloc2(dim,PetscReal,&v0,dim*dim,PetscReal,&J);CHKERRQ(ierr); 278*cb1e1211SMatthew G Knepley ALE::ISieveVisitor::PointRetriever<PETSC_MESH_TYPE::sieve_type> pV(PetscPowInt(this->_mesh->getSieve()->getMaxConeSize(),dim+1), true); 279*cb1e1211SMatthew G Knepley 280*cb1e1211SMatthew G Knepley for (PetscInt c = cStart; c < cEnd; ++c) { 281*cb1e1211SMatthew G Knepley ALE::ISieveTraversal<PETSC_MESH_TYPE::sieve_type>::orientedClosure(*this->_mesh->getSieve(), c, pV); 282*cb1e1211SMatthew G Knepley const PETSC_MESH_TYPE::point_type *oPoints = pV.getPoints(); 283*cb1e1211SMatthew G Knepley const int oSize = pV.getSize(); 284*cb1e1211SMatthew G Knepley int v = 0; 285*cb1e1211SMatthew G Knepley 286*cb1e1211SMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr); 287*cb1e1211SMatthew G Knepley for (PetscInt cl = 0; cl < oSize; ++cl) { 288*cb1e1211SMatthew G Knepley const PetscInt fDim; 289*cb1e1211SMatthew G Knepley 290*cb1e1211SMatthew G Knepley ierr = PetscSectionGetDof(oPoints[cl], &fDim);CHKERRQ(ierr); 291*cb1e1211SMatthew G Knepley if (pointDim) { 292*cb1e1211SMatthew G Knepley for (PetscInt d = 0; d < fDim; ++d, ++v) { 293*cb1e1211SMatthew G Knepley values[v] = (*this->_options.integrate)(v0, J, v, initFunc); 294*cb1e1211SMatthew G Knepley } 295*cb1e1211SMatthew G Knepley } 296*cb1e1211SMatthew G Knepley } 297*cb1e1211SMatthew G Knepley ierr = DMPlexVecSetClosure(dm, NULL, localX, c, values);CHKERRQ(ierr); 298*cb1e1211SMatthew G Knepley pV.clear(); 299*cb1e1211SMatthew G Knepley } 300*cb1e1211SMatthew G Knepley ierr = PetscFree2(v0,J);CHKERRQ(ierr); 301*cb1e1211SMatthew G Knepley ierr = PetscFree(values);CHKERRQ(ierr); 302*cb1e1211SMatthew G Knepley #endif 303*cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 304*cb1e1211SMatthew G Knepley } 305*cb1e1211SMatthew G Knepley 306*cb1e1211SMatthew G Knepley #undef __FUNCT__ 307*cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexProjectFunction" 308*cb1e1211SMatthew G Knepley /*@C 309*cb1e1211SMatthew G Knepley DMPlexProjectFunction - This projects the given function into the function space provided. 310*cb1e1211SMatthew G Knepley 311*cb1e1211SMatthew G Knepley Input Parameters: 312*cb1e1211SMatthew G Knepley + dm - The DM 313*cb1e1211SMatthew G Knepley . numComp - The number of components (functions) 314*cb1e1211SMatthew G Knepley . funcs - The coordinate functions to evaluate 315*cb1e1211SMatthew G Knepley - mode - The insertion mode for values 316*cb1e1211SMatthew G Knepley 317*cb1e1211SMatthew G Knepley Output Parameter: 318*cb1e1211SMatthew G Knepley . X - vector 319*cb1e1211SMatthew G Knepley 320*cb1e1211SMatthew G Knepley Level: developer 321*cb1e1211SMatthew G Knepley 322*cb1e1211SMatthew G Knepley Note: 323*cb1e1211SMatthew G Knepley This currently just calls the function with the coordinates of each vertex and edge midpoint, and stores the result in a vector. 324*cb1e1211SMatthew G Knepley We will eventually fix it. 325*cb1e1211SMatthew G Knepley 326*cb1e1211SMatthew G Knepley ,seealso: DMPlexComputeL2Diff() 327*cb1e1211SMatthew G Knepley */ 328*cb1e1211SMatthew G Knepley PetscErrorCode DMPlexProjectFunction(DM dm, PetscInt numComp, PetscScalar (**funcs)(const PetscReal []), InsertMode mode, Vec X) 329*cb1e1211SMatthew G Knepley { 330*cb1e1211SMatthew G Knepley Vec localX; 331*cb1e1211SMatthew G Knepley PetscErrorCode ierr; 332*cb1e1211SMatthew G Knepley 333*cb1e1211SMatthew G Knepley PetscFunctionBegin; 334*cb1e1211SMatthew G Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 335*cb1e1211SMatthew G Knepley ierr = DMPlexProjectFunctionLocal(dm, numComp, funcs, mode, localX);CHKERRQ(ierr); 336*cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalBegin(dm, localX, mode, X);CHKERRQ(ierr); 337*cb1e1211SMatthew G Knepley ierr = DMLocalToGlobalEnd(dm, localX, mode, X);CHKERRQ(ierr); 338*cb1e1211SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 339*cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 340*cb1e1211SMatthew G Knepley } 341*cb1e1211SMatthew G Knepley 342*cb1e1211SMatthew G Knepley #undef __FUNCT__ 343*cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeL2Diff" 344*cb1e1211SMatthew G Knepley /*@C 345*cb1e1211SMatthew G Knepley DMPlexComputeL2Diff - This function computes the L_2 difference between a function u and an FEM interpolant solution u_h. 346*cb1e1211SMatthew G Knepley 347*cb1e1211SMatthew G Knepley Input Parameters: 348*cb1e1211SMatthew G Knepley + dm - The DM 349*cb1e1211SMatthew G Knepley . quad - The PetscQuadrature object for each field 350*cb1e1211SMatthew G Knepley . funcs - The functions to evaluate for each field component 351*cb1e1211SMatthew G Knepley - X - The coefficient vector u_h 352*cb1e1211SMatthew G Knepley 353*cb1e1211SMatthew G Knepley Output Parameter: 354*cb1e1211SMatthew G Knepley . diff - The diff ||u - u_h||_2 355*cb1e1211SMatthew G Knepley 356*cb1e1211SMatthew G Knepley Level: developer 357*cb1e1211SMatthew G Knepley 358*cb1e1211SMatthew G Knepley .seealso: DMPlexProjectFunction() 359*cb1e1211SMatthew G Knepley */ 360*cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeL2Diff(DM dm, PetscQuadrature quad[], PetscScalar (**funcs)(const PetscReal []), Vec X, PetscReal *diff) 361*cb1e1211SMatthew G Knepley { 362*cb1e1211SMatthew G Knepley const PetscInt debug = 0; 363*cb1e1211SMatthew G Knepley PetscSection section; 364*cb1e1211SMatthew G Knepley Vec localX; 365*cb1e1211SMatthew G Knepley PetscReal *coords, *v0, *J, *invJ, detJ; 366*cb1e1211SMatthew G Knepley PetscReal localDiff = 0.0; 367*cb1e1211SMatthew G Knepley PetscInt dim, numFields, numComponents = 0, cStart, cEnd, c, field, fieldOffset, comp; 368*cb1e1211SMatthew G Knepley PetscErrorCode ierr; 369*cb1e1211SMatthew G Knepley 370*cb1e1211SMatthew G Knepley PetscFunctionBegin; 371*cb1e1211SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 372*cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 373*cb1e1211SMatthew G Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 374*cb1e1211SMatthew G Knepley ierr = DMGetLocalVector(dm, &localX);CHKERRQ(ierr); 375*cb1e1211SMatthew G Knepley ierr = DMGlobalToLocalBegin(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 376*cb1e1211SMatthew G Knepley ierr = DMGlobalToLocalEnd(dm, X, INSERT_VALUES, localX);CHKERRQ(ierr); 377*cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 378*cb1e1211SMatthew G Knepley numComponents += quad[field].numComponents; 379*cb1e1211SMatthew G Knepley } 380*cb1e1211SMatthew G Knepley ierr = DMPlexProjectFunctionLocal(dm, numComponents, funcs, INSERT_BC_VALUES, localX);CHKERRQ(ierr); 381*cb1e1211SMatthew G Knepley ierr = PetscMalloc4(dim,PetscReal,&coords,dim,PetscReal,&v0,dim*dim,PetscReal,&J,dim*dim,PetscReal,&invJ);CHKERRQ(ierr); 382*cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 383*cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 384*cb1e1211SMatthew G Knepley PetscScalar *x; 385*cb1e1211SMatthew G Knepley PetscReal elemDiff = 0.0; 386*cb1e1211SMatthew G Knepley 387*cb1e1211SMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, v0, J, invJ, &detJ);CHKERRQ(ierr); 388*cb1e1211SMatthew G Knepley if (detJ <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ, c); 389*cb1e1211SMatthew G Knepley ierr = DMPlexVecGetClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 390*cb1e1211SMatthew G Knepley 391*cb1e1211SMatthew G Knepley for (field = 0, comp = 0, fieldOffset = 0; field < numFields; ++field) { 392*cb1e1211SMatthew G Knepley const PetscInt numQuadPoints = quad[field].numQuadPoints; 393*cb1e1211SMatthew G Knepley const PetscReal *quadPoints = quad[field].quadPoints; 394*cb1e1211SMatthew G Knepley const PetscReal *quadWeights = quad[field].quadWeights; 395*cb1e1211SMatthew G Knepley const PetscInt numBasisFuncs = quad[field].numBasisFuncs; 396*cb1e1211SMatthew G Knepley const PetscInt numBasisComps = quad[field].numComponents; 397*cb1e1211SMatthew G Knepley const PetscReal *basis = quad[field].basis; 398*cb1e1211SMatthew G Knepley PetscInt q, d, e, fc, f; 399*cb1e1211SMatthew G Knepley 400*cb1e1211SMatthew G Knepley if (debug) { 401*cb1e1211SMatthew G Knepley char title[1024]; 402*cb1e1211SMatthew G Knepley ierr = PetscSNPrintf(title, 1023, "Solution for Field %d", field);CHKERRQ(ierr); 403*cb1e1211SMatthew G Knepley ierr = DMPrintCellVector(c, title, numBasisFuncs*numBasisComps, &x[fieldOffset]);CHKERRQ(ierr); 404*cb1e1211SMatthew G Knepley } 405*cb1e1211SMatthew G Knepley for (q = 0; q < numQuadPoints; ++q) { 406*cb1e1211SMatthew G Knepley for (d = 0; d < dim; d++) { 407*cb1e1211SMatthew G Knepley coords[d] = v0[d]; 408*cb1e1211SMatthew G Knepley for (e = 0; e < dim; e++) { 409*cb1e1211SMatthew G Knepley coords[d] += J[d*dim+e]*(quadPoints[q*dim+e] + 1.0); 410*cb1e1211SMatthew G Knepley } 411*cb1e1211SMatthew G Knepley } 412*cb1e1211SMatthew G Knepley for (fc = 0; fc < numBasisComps; ++fc) { 413*cb1e1211SMatthew G Knepley const PetscReal funcVal = PetscRealPart((*funcs[comp+fc])(coords)); 414*cb1e1211SMatthew G Knepley PetscReal interpolant = 0.0; 415*cb1e1211SMatthew G Knepley for (f = 0; f < numBasisFuncs; ++f) { 416*cb1e1211SMatthew G Knepley const PetscInt fidx = f*numBasisComps+fc; 417*cb1e1211SMatthew G Knepley interpolant += PetscRealPart(x[fieldOffset+fidx])*basis[q*numBasisFuncs*numBasisComps+fidx]; 418*cb1e1211SMatthew G Knepley } 419*cb1e1211SMatthew G Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d field %d diff %g\n", c, field, PetscSqr(interpolant - funcVal)*quadWeights[q]*detJ);CHKERRQ(ierr);} 420*cb1e1211SMatthew G Knepley elemDiff += PetscSqr(interpolant - funcVal)*quadWeights[q]*detJ; 421*cb1e1211SMatthew G Knepley } 422*cb1e1211SMatthew G Knepley } 423*cb1e1211SMatthew G Knepley comp += numBasisComps; 424*cb1e1211SMatthew G Knepley fieldOffset += numBasisFuncs*numBasisComps; 425*cb1e1211SMatthew G Knepley } 426*cb1e1211SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, localX, c, NULL, &x);CHKERRQ(ierr); 427*cb1e1211SMatthew G Knepley if (debug) {ierr = PetscPrintf(PETSC_COMM_SELF, " elem %d diff %g\n", c, elemDiff);CHKERRQ(ierr);} 428*cb1e1211SMatthew G Knepley localDiff += elemDiff; 429*cb1e1211SMatthew G Knepley } 430*cb1e1211SMatthew G Knepley ierr = PetscFree4(coords,v0,J,invJ);CHKERRQ(ierr); 431*cb1e1211SMatthew G Knepley ierr = DMRestoreLocalVector(dm, &localX);CHKERRQ(ierr); 432*cb1e1211SMatthew G Knepley ierr = MPI_Allreduce(&localDiff, diff, 1, MPIU_REAL, MPI_SUM, PETSC_COMM_WORLD);CHKERRQ(ierr); 433*cb1e1211SMatthew G Knepley *diff = PetscSqrtReal(*diff); 434*cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 435*cb1e1211SMatthew G Knepley } 436*cb1e1211SMatthew G Knepley 437*cb1e1211SMatthew G Knepley #undef __FUNCT__ 438*cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeResidualFEM" 439*cb1e1211SMatthew G Knepley /*@ 440*cb1e1211SMatthew G Knepley DMPlexComputeResidualFEM - Form the local residual F from the local input X using pointwise functions specified by the user 441*cb1e1211SMatthew G Knepley 442*cb1e1211SMatthew G Knepley Input Parameters: 443*cb1e1211SMatthew G Knepley + dm - The mesh 444*cb1e1211SMatthew G Knepley . X - Local input vector 445*cb1e1211SMatthew G Knepley - user - The user context 446*cb1e1211SMatthew G Knepley 447*cb1e1211SMatthew G Knepley Output Parameter: 448*cb1e1211SMatthew G Knepley . F - Local output vector 449*cb1e1211SMatthew G Knepley 450*cb1e1211SMatthew G Knepley Note: 451*cb1e1211SMatthew G Knepley The second member of the user context must be an FEMContext. 452*cb1e1211SMatthew G Knepley 453*cb1e1211SMatthew G Knepley We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 454*cb1e1211SMatthew G Knepley like a GPU, or vectorize on a multicore machine. 455*cb1e1211SMatthew G Knepley 456*cb1e1211SMatthew G Knepley .seealso: DMPlexComputeJacobianActionFEM() 457*cb1e1211SMatthew G Knepley */ 458*cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeResidualFEM(DM dm, Vec X, Vec F, void *user) 459*cb1e1211SMatthew G Knepley { 460*cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 461*cb1e1211SMatthew G Knepley PetscFEM *fem = (PetscFEM*) &((DM*) user)[1]; 462*cb1e1211SMatthew G Knepley PetscQuadrature *quad = fem->quad; 463*cb1e1211SMatthew G Knepley PetscSection section; 464*cb1e1211SMatthew G Knepley PetscReal *v0, *J, *invJ, *detJ; 465*cb1e1211SMatthew G Knepley PetscScalar *elemVec, *u; 466*cb1e1211SMatthew G Knepley PetscInt dim, numFields, field, numBatchesTmp = 1, numCells, cStart, cEnd, c; 467*cb1e1211SMatthew G Knepley PetscInt cellDof = 0, numComponents = 0; 468*cb1e1211SMatthew G Knepley PetscErrorCode ierr; 469*cb1e1211SMatthew G Knepley 470*cb1e1211SMatthew G Knepley PetscFunctionBegin; 471*cb1e1211SMatthew G Knepley /* ierr = PetscLogEventBegin(ResidualFEMEvent,0,0,0,0);CHKERRQ(ierr); */ 472*cb1e1211SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 473*cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 474*cb1e1211SMatthew G Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 475*cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 476*cb1e1211SMatthew G Knepley numCells = cEnd - cStart; 477*cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 478*cb1e1211SMatthew G Knepley cellDof += quad[field].numBasisFuncs*quad[field].numComponents; 479*cb1e1211SMatthew G Knepley numComponents += quad[field].numComponents; 480*cb1e1211SMatthew G Knepley } 481*cb1e1211SMatthew G Knepley ierr = DMPlexProjectFunctionLocal(dm, numComponents, fem->bcFuncs, INSERT_BC_VALUES, X);CHKERRQ(ierr); 482*cb1e1211SMatthew G Knepley ierr = VecSet(F, 0.0);CHKERRQ(ierr); 483*cb1e1211SMatthew 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); 484*cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 485*cb1e1211SMatthew G Knepley PetscScalar *x; 486*cb1e1211SMatthew G Knepley PetscInt i; 487*cb1e1211SMatthew G Knepley 488*cb1e1211SMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 489*cb1e1211SMatthew G Knepley if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 490*cb1e1211SMatthew G Knepley ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 491*cb1e1211SMatthew G Knepley 492*cb1e1211SMatthew G Knepley for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 493*cb1e1211SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 494*cb1e1211SMatthew G Knepley } 495*cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 496*cb1e1211SMatthew G Knepley const PetscInt numQuadPoints = quad[field].numQuadPoints; 497*cb1e1211SMatthew G Knepley const PetscInt numBasisFuncs = quad[field].numBasisFuncs; 498*cb1e1211SMatthew G Knepley void (*f0)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f0[]) = fem->f0Funcs[field]; 499*cb1e1211SMatthew G Knepley void (*f1)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar f1[]) = fem->f1Funcs[field]; 500*cb1e1211SMatthew G Knepley /* Conforming batches */ 501*cb1e1211SMatthew G Knepley PetscInt blockSize = numBasisFuncs*numQuadPoints; 502*cb1e1211SMatthew G Knepley PetscInt numBlocks = 1; 503*cb1e1211SMatthew G Knepley PetscInt batchSize = numBlocks * blockSize; 504*cb1e1211SMatthew G Knepley PetscInt numBatches = numBatchesTmp; 505*cb1e1211SMatthew G Knepley PetscInt numChunks = numCells / (numBatches*batchSize); 506*cb1e1211SMatthew G Knepley /* Remainder */ 507*cb1e1211SMatthew G Knepley PetscInt numRemainder = numCells % (numBatches * batchSize); 508*cb1e1211SMatthew G Knepley PetscInt offset = numCells - numRemainder; 509*cb1e1211SMatthew G Knepley 510*cb1e1211SMatthew G Knepley ierr = (*mesh->integrateResidualFEM)(numChunks*numBatches*batchSize, numFields, field, quad, u, v0, J, invJ, detJ, f0, f1, elemVec);CHKERRQ(ierr); 511*cb1e1211SMatthew G Knepley ierr = (*mesh->integrateResidualFEM)(numRemainder, numFields, field, quad, &u[offset*cellDof], &v0[offset*dim], &J[offset*dim*dim], &invJ[offset*dim*dim], &detJ[offset], 512*cb1e1211SMatthew G Knepley f0, f1, &elemVec[offset*cellDof]);CHKERRQ(ierr); 513*cb1e1211SMatthew G Knepley } 514*cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 515*cb1e1211SMatthew G Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Residual", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 516*cb1e1211SMatthew G Knepley ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 517*cb1e1211SMatthew G Knepley } 518*cb1e1211SMatthew G Knepley ierr = PetscFree6(u,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 519*cb1e1211SMatthew G Knepley if (mesh->printFEM) { 520*cb1e1211SMatthew G Knepley PetscMPIInt rank, numProcs; 521*cb1e1211SMatthew G Knepley PetscInt p; 522*cb1e1211SMatthew G Knepley 523*cb1e1211SMatthew G Knepley ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 524*cb1e1211SMatthew G Knepley ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr); 525*cb1e1211SMatthew G Knepley ierr = PetscPrintf(PETSC_COMM_WORLD, "Residual:\n");CHKERRQ(ierr); 526*cb1e1211SMatthew G Knepley for (p = 0; p < numProcs; ++p) { 527*cb1e1211SMatthew G Knepley if (p == rank) { 528*cb1e1211SMatthew G Knepley Vec f; 529*cb1e1211SMatthew G Knepley 530*cb1e1211SMatthew G Knepley ierr = VecDuplicate(F, &f);CHKERRQ(ierr); 531*cb1e1211SMatthew G Knepley ierr = VecCopy(F, f);CHKERRQ(ierr); 532*cb1e1211SMatthew G Knepley ierr = VecChop(f, 1.0e-10);CHKERRQ(ierr); 533*cb1e1211SMatthew G Knepley ierr = VecView(f, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 534*cb1e1211SMatthew G Knepley ierr = VecDestroy(&f);CHKERRQ(ierr); 535*cb1e1211SMatthew G Knepley ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 536*cb1e1211SMatthew G Knepley } 537*cb1e1211SMatthew G Knepley ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr); 538*cb1e1211SMatthew G Knepley } 539*cb1e1211SMatthew G Knepley } 540*cb1e1211SMatthew G Knepley /* ierr = PetscLogEventEnd(ResidualFEMEvent,0,0,0,0);CHKERRQ(ierr); */ 541*cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 542*cb1e1211SMatthew G Knepley } 543*cb1e1211SMatthew G Knepley 544*cb1e1211SMatthew G Knepley #undef __FUNCT__ 545*cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianActionFEM" 546*cb1e1211SMatthew G Knepley /*@C 547*cb1e1211SMatthew G Knepley DMPlexComputeJacobianActionFEM - Form the local action of Jacobian J(u) on the local input X using pointwise functions specified by the user 548*cb1e1211SMatthew G Knepley 549*cb1e1211SMatthew G Knepley Input Parameters: 550*cb1e1211SMatthew G Knepley + dm - The mesh 551*cb1e1211SMatthew G Knepley . J - The Jacobian shell matrix 552*cb1e1211SMatthew G Knepley . X - Local input vector 553*cb1e1211SMatthew G Knepley - user - The user context 554*cb1e1211SMatthew G Knepley 555*cb1e1211SMatthew G Knepley Output Parameter: 556*cb1e1211SMatthew G Knepley . F - Local output vector 557*cb1e1211SMatthew G Knepley 558*cb1e1211SMatthew G Knepley Note: 559*cb1e1211SMatthew G Knepley The second member of the user context must be an FEMContext. 560*cb1e1211SMatthew G Knepley 561*cb1e1211SMatthew G Knepley We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 562*cb1e1211SMatthew G Knepley like a GPU, or vectorize on a multicore machine. 563*cb1e1211SMatthew G Knepley 564*cb1e1211SMatthew G Knepley .seealso: DMPlexComputeResidualFEM() 565*cb1e1211SMatthew G Knepley */ 566*cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianActionFEM(DM dm, Mat Jac, Vec X, Vec F, void *user) 567*cb1e1211SMatthew G Knepley { 568*cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 569*cb1e1211SMatthew G Knepley PetscFEM *fem = (PetscFEM*) &((DM*) user)[1]; 570*cb1e1211SMatthew G Knepley PetscQuadrature *quad = fem->quad; 571*cb1e1211SMatthew G Knepley PetscSection section; 572*cb1e1211SMatthew G Knepley JacActionCtx *jctx; 573*cb1e1211SMatthew G Knepley PetscReal *v0, *J, *invJ, *detJ; 574*cb1e1211SMatthew G Knepley PetscScalar *elemVec, *u, *a; 575*cb1e1211SMatthew G Knepley PetscInt dim, numFields, field, numBatchesTmp = 1, numCells, cStart, cEnd, c; 576*cb1e1211SMatthew G Knepley PetscInt cellDof = 0; 577*cb1e1211SMatthew G Knepley PetscErrorCode ierr; 578*cb1e1211SMatthew G Knepley 579*cb1e1211SMatthew G Knepley PetscFunctionBegin; 580*cb1e1211SMatthew G Knepley /* ierr = PetscLogEventBegin(JacobianActionFEMEvent,0,0,0,0);CHKERRQ(ierr); */ 581*cb1e1211SMatthew G Knepley ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 582*cb1e1211SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 583*cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 584*cb1e1211SMatthew G Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 585*cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 586*cb1e1211SMatthew G Knepley numCells = cEnd - cStart; 587*cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 588*cb1e1211SMatthew G Knepley cellDof += quad[field].numBasisFuncs*quad[field].numComponents; 589*cb1e1211SMatthew G Knepley } 590*cb1e1211SMatthew G Knepley ierr = VecSet(F, 0.0);CHKERRQ(ierr); 591*cb1e1211SMatthew 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); 592*cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 593*cb1e1211SMatthew G Knepley PetscScalar *x; 594*cb1e1211SMatthew G Knepley PetscInt i; 595*cb1e1211SMatthew G Knepley 596*cb1e1211SMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 597*cb1e1211SMatthew G Knepley if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 598*cb1e1211SMatthew G Knepley ierr = DMPlexVecGetClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 599*cb1e1211SMatthew G Knepley for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 600*cb1e1211SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, jctx->u, c, NULL, &x);CHKERRQ(ierr); 601*cb1e1211SMatthew G Knepley ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 602*cb1e1211SMatthew G Knepley for (i = 0; i < cellDof; ++i) a[c*cellDof+i] = x[i]; 603*cb1e1211SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 604*cb1e1211SMatthew G Knepley } 605*cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 606*cb1e1211SMatthew G Knepley const PetscInt numQuadPoints = quad[field].numQuadPoints; 607*cb1e1211SMatthew G Knepley const PetscInt numBasisFuncs = quad[field].numBasisFuncs; 608*cb1e1211SMatthew G Knepley /* Conforming batches */ 609*cb1e1211SMatthew G Knepley PetscInt blockSize = numBasisFuncs*numQuadPoints; 610*cb1e1211SMatthew G Knepley PetscInt numBlocks = 1; 611*cb1e1211SMatthew G Knepley PetscInt batchSize = numBlocks * blockSize; 612*cb1e1211SMatthew G Knepley PetscInt numBatches = numBatchesTmp; 613*cb1e1211SMatthew G Knepley PetscInt numChunks = numCells / (numBatches*batchSize); 614*cb1e1211SMatthew G Knepley /* Remainder */ 615*cb1e1211SMatthew G Knepley PetscInt numRemainder = numCells % (numBatches * batchSize); 616*cb1e1211SMatthew G Knepley PetscInt offset = numCells - numRemainder; 617*cb1e1211SMatthew G Knepley 618*cb1e1211SMatthew G Knepley ierr = (*mesh->integrateJacobianActionFEM)(numChunks*numBatches*batchSize, numFields, field, quad, u, a, v0, J, invJ, detJ, fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, elemVec);CHKERRQ(ierr); 619*cb1e1211SMatthew G Knepley ierr = (*mesh->integrateJacobianActionFEM)(numRemainder, numFields, field, quad, &u[offset*cellDof], &a[offset*cellDof], &v0[offset*dim], &J[offset*dim*dim], &invJ[offset*dim*dim], &detJ[offset], 620*cb1e1211SMatthew G Knepley fem->g0Funcs, fem->g1Funcs, fem->g2Funcs, fem->g3Funcs, &elemVec[offset*cellDof]);CHKERRQ(ierr); 621*cb1e1211SMatthew G Knepley } 622*cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 623*cb1e1211SMatthew G Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellVector(c, "Jacobian Action", cellDof, &elemVec[c*cellDof]);CHKERRQ(ierr);} 624*cb1e1211SMatthew G Knepley ierr = DMPlexVecSetClosure(dm, NULL, F, c, &elemVec[c*cellDof], ADD_VALUES);CHKERRQ(ierr); 625*cb1e1211SMatthew G Knepley } 626*cb1e1211SMatthew G Knepley ierr = PetscFree7(u,a,v0,J,invJ,detJ,elemVec);CHKERRQ(ierr); 627*cb1e1211SMatthew G Knepley if (mesh->printFEM) { 628*cb1e1211SMatthew G Knepley PetscMPIInt rank, numProcs; 629*cb1e1211SMatthew G Knepley PetscInt p; 630*cb1e1211SMatthew G Knepley 631*cb1e1211SMatthew G Knepley ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm), &rank);CHKERRQ(ierr); 632*cb1e1211SMatthew G Knepley ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &numProcs);CHKERRQ(ierr); 633*cb1e1211SMatthew G Knepley ierr = PetscPrintf(PETSC_COMM_WORLD, "Jacobian Action:\n");CHKERRQ(ierr); 634*cb1e1211SMatthew G Knepley for (p = 0; p < numProcs; ++p) { 635*cb1e1211SMatthew G Knepley if (p == rank) {ierr = VecView(F, PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);} 636*cb1e1211SMatthew G Knepley ierr = PetscBarrier((PetscObject) dm);CHKERRQ(ierr); 637*cb1e1211SMatthew G Knepley } 638*cb1e1211SMatthew G Knepley } 639*cb1e1211SMatthew G Knepley /* ierr = PetscLogEventEnd(JacobianActionFEMEvent,0,0,0,0);CHKERRQ(ierr); */ 640*cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 641*cb1e1211SMatthew G Knepley } 642*cb1e1211SMatthew G Knepley 643*cb1e1211SMatthew G Knepley #undef __FUNCT__ 644*cb1e1211SMatthew G Knepley #define __FUNCT__ "DMPlexComputeJacobianFEM" 645*cb1e1211SMatthew G Knepley /*@ 646*cb1e1211SMatthew G Knepley DMPlexComputeJacobianFEM - Form the local portion of the Jacobian matrix J at the local solution X using pointwise functions specified by the user. 647*cb1e1211SMatthew G Knepley 648*cb1e1211SMatthew G Knepley Input Parameters: 649*cb1e1211SMatthew G Knepley + dm - The mesh 650*cb1e1211SMatthew G Knepley . X - Local input vector 651*cb1e1211SMatthew G Knepley - user - The user context 652*cb1e1211SMatthew G Knepley 653*cb1e1211SMatthew G Knepley Output Parameter: 654*cb1e1211SMatthew G Knepley . Jac - Jacobian matrix 655*cb1e1211SMatthew G Knepley 656*cb1e1211SMatthew G Knepley Note: 657*cb1e1211SMatthew G Knepley The second member of the user context must be an FEMContext. 658*cb1e1211SMatthew G Knepley 659*cb1e1211SMatthew G Knepley We form the residual one batch of elements at a time. This allows us to offload work onto an accelerator, 660*cb1e1211SMatthew G Knepley like a GPU, or vectorize on a multicore machine. 661*cb1e1211SMatthew G Knepley 662*cb1e1211SMatthew G Knepley .seealso: FormFunctionLocal() 663*cb1e1211SMatthew G Knepley */ 664*cb1e1211SMatthew G Knepley PetscErrorCode DMPlexComputeJacobianFEM(DM dm, Vec X, Mat Jac, Mat JacP, MatStructure *str,void *user) 665*cb1e1211SMatthew G Knepley { 666*cb1e1211SMatthew G Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 667*cb1e1211SMatthew G Knepley PetscFEM *fem = (PetscFEM*) &((DM*) user)[1]; 668*cb1e1211SMatthew G Knepley PetscQuadrature *quad = fem->quad; 669*cb1e1211SMatthew G Knepley PetscSection section; 670*cb1e1211SMatthew G Knepley PetscReal *v0, *J, *invJ, *detJ; 671*cb1e1211SMatthew G Knepley PetscScalar *elemMat, *u; 672*cb1e1211SMatthew G Knepley PetscInt dim, numFields, field, fieldI, numBatchesTmp = 1, numCells, cStart, cEnd, c; 673*cb1e1211SMatthew G Knepley PetscInt cellDof = 0, numComponents = 0; 674*cb1e1211SMatthew G Knepley PetscBool isShell; 675*cb1e1211SMatthew G Knepley PetscErrorCode ierr; 676*cb1e1211SMatthew G Knepley 677*cb1e1211SMatthew G Knepley PetscFunctionBegin; 678*cb1e1211SMatthew G Knepley /* ierr = PetscLogEventBegin(JacobianFEMEvent,0,0,0,0);CHKERRQ(ierr); */ 679*cb1e1211SMatthew G Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 680*cb1e1211SMatthew G Knepley ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 681*cb1e1211SMatthew G Knepley ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 682*cb1e1211SMatthew G Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 683*cb1e1211SMatthew G Knepley numCells = cEnd - cStart; 684*cb1e1211SMatthew G Knepley for (field = 0; field < numFields; ++field) { 685*cb1e1211SMatthew G Knepley cellDof += quad[field].numBasisFuncs*quad[field].numComponents; 686*cb1e1211SMatthew G Knepley numComponents += quad[field].numComponents; 687*cb1e1211SMatthew G Knepley } 688*cb1e1211SMatthew G Knepley ierr = DMPlexProjectFunctionLocal(dm, numComponents, fem->bcFuncs, INSERT_BC_VALUES, X);CHKERRQ(ierr); 689*cb1e1211SMatthew G Knepley ierr = MatZeroEntries(JacP);CHKERRQ(ierr); 690*cb1e1211SMatthew 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); 691*cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 692*cb1e1211SMatthew G Knepley PetscScalar *x; 693*cb1e1211SMatthew G Knepley PetscInt i; 694*cb1e1211SMatthew G Knepley 695*cb1e1211SMatthew G Knepley ierr = DMPlexComputeCellGeometry(dm, c, &v0[c*dim], &J[c*dim*dim], &invJ[c*dim*dim], &detJ[c]);CHKERRQ(ierr); 696*cb1e1211SMatthew G Knepley if (detJ[c] <= 0.0) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid determinant %g for element %d", detJ[c], c); 697*cb1e1211SMatthew G Knepley ierr = DMPlexVecGetClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 698*cb1e1211SMatthew G Knepley 699*cb1e1211SMatthew G Knepley for (i = 0; i < cellDof; ++i) u[c*cellDof+i] = x[i]; 700*cb1e1211SMatthew G Knepley ierr = DMPlexVecRestoreClosure(dm, NULL, X, c, NULL, &x);CHKERRQ(ierr); 701*cb1e1211SMatthew G Knepley } 702*cb1e1211SMatthew G Knepley ierr = PetscMemzero(elemMat, numCells*cellDof*cellDof * sizeof(PetscScalar));CHKERRQ(ierr); 703*cb1e1211SMatthew G Knepley for (fieldI = 0; fieldI < numFields; ++fieldI) { 704*cb1e1211SMatthew G Knepley const PetscInt numQuadPoints = quad[fieldI].numQuadPoints; 705*cb1e1211SMatthew G Knepley const PetscInt numBasisFuncs = quad[fieldI].numBasisFuncs; 706*cb1e1211SMatthew G Knepley PetscInt fieldJ; 707*cb1e1211SMatthew G Knepley 708*cb1e1211SMatthew G Knepley for (fieldJ = 0; fieldJ < numFields; ++fieldJ) { 709*cb1e1211SMatthew G Knepley void (*g0)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar g0[]) = fem->g0Funcs[fieldI*numFields+fieldJ]; 710*cb1e1211SMatthew G Knepley void (*g1)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar g1[]) = fem->g1Funcs[fieldI*numFields+fieldJ]; 711*cb1e1211SMatthew G Knepley void (*g2)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar g2[]) = fem->g2Funcs[fieldI*numFields+fieldJ]; 712*cb1e1211SMatthew G Knepley void (*g3)(const PetscScalar u[], const PetscScalar gradU[], const PetscReal x[], PetscScalar g3[]) = fem->g3Funcs[fieldI*numFields+fieldJ]; 713*cb1e1211SMatthew G Knepley /* Conforming batches */ 714*cb1e1211SMatthew G Knepley PetscInt blockSize = numBasisFuncs*numQuadPoints; 715*cb1e1211SMatthew G Knepley PetscInt numBlocks = 1; 716*cb1e1211SMatthew G Knepley PetscInt batchSize = numBlocks * blockSize; 717*cb1e1211SMatthew G Knepley PetscInt numBatches = numBatchesTmp; 718*cb1e1211SMatthew G Knepley PetscInt numChunks = numCells / (numBatches*batchSize); 719*cb1e1211SMatthew G Knepley /* Remainder */ 720*cb1e1211SMatthew G Knepley PetscInt numRemainder = numCells % (numBatches * batchSize); 721*cb1e1211SMatthew G Knepley PetscInt offset = numCells - numRemainder; 722*cb1e1211SMatthew G Knepley 723*cb1e1211SMatthew G Knepley ierr = (*mesh->integrateJacobianFEM)(numChunks*numBatches*batchSize, numFields, fieldI, fieldJ, quad, u, v0, J, invJ, detJ, g0, g1, g2, g3, elemMat);CHKERRQ(ierr); 724*cb1e1211SMatthew G Knepley ierr = (*mesh->integrateJacobianFEM)(numRemainder, numFields, fieldI, fieldJ, quad, &u[offset*cellDof], &v0[offset*dim], &J[offset*dim*dim], &invJ[offset*dim*dim], &detJ[offset], 725*cb1e1211SMatthew G Knepley g0, g1, g2, g3, &elemMat[offset*cellDof*cellDof]);CHKERRQ(ierr); 726*cb1e1211SMatthew G Knepley } 727*cb1e1211SMatthew G Knepley } 728*cb1e1211SMatthew G Knepley for (c = cStart; c < cEnd; ++c) { 729*cb1e1211SMatthew G Knepley if (mesh->printFEM > 1) {ierr = DMPrintCellMatrix(c, "Jacobian", cellDof, cellDof, &elemMat[c*cellDof*cellDof]);CHKERRQ(ierr);} 730*cb1e1211SMatthew G Knepley ierr = DMPlexMatSetClosure(dm, NULL, NULL, JacP, c, &elemMat[c*cellDof*cellDof], ADD_VALUES);CHKERRQ(ierr); 731*cb1e1211SMatthew G Knepley } 732*cb1e1211SMatthew G Knepley ierr = PetscFree6(u,v0,J,invJ,detJ,elemMat);CHKERRQ(ierr); 733*cb1e1211SMatthew G Knepley 734*cb1e1211SMatthew G Knepley /* Assemble matrix, using the 2-step process: 735*cb1e1211SMatthew G Knepley MatAssemblyBegin(), MatAssemblyEnd(). */ 736*cb1e1211SMatthew G Knepley ierr = MatAssemblyBegin(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 737*cb1e1211SMatthew G Knepley ierr = MatAssemblyEnd(JacP, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 738*cb1e1211SMatthew G Knepley 739*cb1e1211SMatthew G Knepley if (mesh->printFEM) { 740*cb1e1211SMatthew G Knepley ierr = PetscPrintf(PETSC_COMM_WORLD, "Jacobian:\n");CHKERRQ(ierr); 741*cb1e1211SMatthew G Knepley ierr = MatChop(JacP, 1.0e-10);CHKERRQ(ierr); 742*cb1e1211SMatthew G Knepley ierr = MatView(JacP, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 743*cb1e1211SMatthew G Knepley } 744*cb1e1211SMatthew G Knepley /* ierr = PetscLogEventEnd(JacobianFEMEvent,0,0,0,0);CHKERRQ(ierr); */ 745*cb1e1211SMatthew G Knepley ierr = PetscObjectTypeCompare((PetscObject)Jac, MATSHELL, &isShell);CHKERRQ(ierr); 746*cb1e1211SMatthew G Knepley if (isShell) { 747*cb1e1211SMatthew G Knepley JacActionCtx *jctx; 748*cb1e1211SMatthew G Knepley 749*cb1e1211SMatthew G Knepley ierr = MatShellGetContext(Jac, &jctx);CHKERRQ(ierr); 750*cb1e1211SMatthew G Knepley ierr = VecCopy(X, jctx->u);CHKERRQ(ierr); 751*cb1e1211SMatthew G Knepley } 752*cb1e1211SMatthew G Knepley *str = SAME_NONZERO_PATTERN; 753*cb1e1211SMatthew G Knepley PetscFunctionReturn(0); 754*cb1e1211SMatthew G Knepley } 755