10e2ec84fSDave May #include <petscdm.h> 20e2ec84fSDave May #include <petscdmda.h> 30e2ec84fSDave May #include <petscdmswarm.h> 41f52046fSDave May #include <petsc/private/dmswarmimpl.h> 5279f676cSBarry Smith #include "../src/dm/impls/swarm/data_bucket.h" 60e2ec84fSDave May 79371c9d4SSatish Balay PetscErrorCode private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(PetscInt dim, PetscInt np[], PetscInt *_npoints, PetscReal **_xi) { 80e2ec84fSDave May PetscReal *xi; 9f9368e01SSatish Balay PetscInt d, npoints = 0, cnt; 100e2ec84fSDave May PetscReal ds[] = {0.0, 0.0, 0.0}; 110e2ec84fSDave May PetscInt ii, jj, kk; 120e2ec84fSDave May 130e2ec84fSDave May PetscFunctionBegin; 140e2ec84fSDave May switch (dim) { 159371c9d4SSatish Balay case 1: npoints = np[0]; break; 169371c9d4SSatish Balay case 2: npoints = np[0] * np[1]; break; 179371c9d4SSatish Balay case 3: npoints = np[0] * np[1] * np[2]; break; 180e2ec84fSDave May } 199371c9d4SSatish Balay for (d = 0; d < dim; d++) { ds[d] = 2.0 / ((PetscReal)np[d]); } 200e2ec84fSDave May 219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * npoints, &xi)); 220e2ec84fSDave May switch (dim) { 230e2ec84fSDave May case 1: 240e2ec84fSDave May cnt = 0; 250e2ec84fSDave May for (ii = 0; ii < np[0]; ii++) { 260e2ec84fSDave May xi[dim * cnt + 0] = -1.0 + 0.5 * ds[d] + ii * ds[0]; 270e2ec84fSDave May cnt++; 280e2ec84fSDave May } 290e2ec84fSDave May break; 300e2ec84fSDave May 310e2ec84fSDave May case 2: 320e2ec84fSDave May cnt = 0; 330e2ec84fSDave May for (jj = 0; jj < np[1]; jj++) { 340e2ec84fSDave May for (ii = 0; ii < np[0]; ii++) { 350e2ec84fSDave May xi[dim * cnt + 0] = -1.0 + 0.5 * ds[0] + ii * ds[0]; 360e2ec84fSDave May xi[dim * cnt + 1] = -1.0 + 0.5 * ds[1] + jj * ds[1]; 370e2ec84fSDave May cnt++; 380e2ec84fSDave May } 390e2ec84fSDave May } 400e2ec84fSDave May break; 410e2ec84fSDave May 420e2ec84fSDave May case 3: 430e2ec84fSDave May cnt = 0; 440e2ec84fSDave May for (kk = 0; kk < np[2]; kk++) { 450e2ec84fSDave May for (jj = 0; jj < np[1]; jj++) { 460e2ec84fSDave May for (ii = 0; ii < np[0]; ii++) { 470e2ec84fSDave May xi[dim * cnt + 0] = -1.0 + 0.5 * ds[0] + ii * ds[0]; 480e2ec84fSDave May xi[dim * cnt + 1] = -1.0 + 0.5 * ds[1] + jj * ds[1]; 490e2ec84fSDave May xi[dim * cnt + 2] = -1.0 + 0.5 * ds[2] + kk * ds[2]; 500e2ec84fSDave May cnt++; 510e2ec84fSDave May } 520e2ec84fSDave May } 530e2ec84fSDave May } 540e2ec84fSDave May break; 550e2ec84fSDave May } 560e2ec84fSDave May *_npoints = npoints; 570e2ec84fSDave May *_xi = xi; 580e2ec84fSDave May PetscFunctionReturn(0); 590e2ec84fSDave May } 600e2ec84fSDave May 619371c9d4SSatish Balay PetscErrorCode private_DMSwarmCreateCellLocalCoords_DA_Q1_Gauss(PetscInt dim, PetscInt np_1d, PetscInt *_npoints, PetscReal **_xi) { 620e2ec84fSDave May PetscQuadrature quadrature; 630e2ec84fSDave May const PetscReal *quadrature_xi; 640e2ec84fSDave May PetscReal *xi; 650e2ec84fSDave May PetscInt d, q, npoints_q; 660e2ec84fSDave May 670e2ec84fSDave May PetscFunctionBegin; 689566063dSJacob Faibussowitsch PetscCall(PetscDTGaussTensorQuadrature(dim, 1, np_1d, -1.0, 1.0, &quadrature)); 699566063dSJacob Faibussowitsch PetscCall(PetscQuadratureGetData(quadrature, NULL, NULL, &npoints_q, &quadrature_xi, NULL)); 709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * npoints_q, &xi)); 710e2ec84fSDave May for (q = 0; q < npoints_q; q++) { 729371c9d4SSatish Balay for (d = 0; d < dim; d++) { xi[dim * q + d] = quadrature_xi[dim * q + d]; } 730e2ec84fSDave May } 749566063dSJacob Faibussowitsch PetscCall(PetscQuadratureDestroy(&quadrature)); 750e2ec84fSDave May *_npoints = npoints_q; 760e2ec84fSDave May *_xi = xi; 770e2ec84fSDave May PetscFunctionReturn(0); 780e2ec84fSDave May } 790e2ec84fSDave May 809371c9d4SSatish Balay PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA_Q1(DM dm, DM dmc, PetscInt npoints, DMSwarmPICLayoutType layout) { 810e2ec84fSDave May PetscInt dim, npoints_q; 820e2ec84fSDave May PetscInt nel, npe, e, q, k, d; 830e2ec84fSDave May const PetscInt *element_list; 840e2ec84fSDave May PetscReal **basis; 850e2ec84fSDave May PetscReal *xi; 860e2ec84fSDave May Vec coor; 870e2ec84fSDave May const PetscScalar *_coor; 880e2ec84fSDave May PetscReal *elcoor; 890e2ec84fSDave May PetscReal *swarm_coor; 900e2ec84fSDave May PetscInt *swarm_cellid; 910e2ec84fSDave May PetscInt pcnt; 920e2ec84fSDave May 930e2ec84fSDave May PetscFunctionBegin; 949566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 950e2ec84fSDave May switch (layout) { 969371c9d4SSatish Balay case DMSWARMPIC_LAYOUT_REGULAR: { 970e2ec84fSDave May PetscInt np_dir[3]; 980e2ec84fSDave May np_dir[0] = np_dir[1] = np_dir[2] = npoints; 999566063dSJacob Faibussowitsch PetscCall(private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(dim, np_dir, &npoints_q, &xi)); 1009371c9d4SSatish Balay } break; 1019371c9d4SSatish Balay case DMSWARMPIC_LAYOUT_GAUSS: PetscCall(private_DMSwarmCreateCellLocalCoords_DA_Q1_Gauss(dim, npoints, &npoints_q, &xi)); break; 102d3393ee1SDave May 1039371c9d4SSatish Balay case DMSWARMPIC_LAYOUT_SUBDIVISION: { 104d3393ee1SDave May PetscInt s, nsub; 105d3393ee1SDave May PetscInt np_dir[3]; 106d3393ee1SDave May nsub = npoints; 107d3393ee1SDave May np_dir[0] = 1; 1089371c9d4SSatish Balay for (s = 0; s < nsub; s++) { np_dir[0] *= 2; } 109d3393ee1SDave May np_dir[1] = np_dir[0]; 110d3393ee1SDave May np_dir[2] = np_dir[0]; 1119566063dSJacob Faibussowitsch PetscCall(private_DMSwarmCreateCellLocalCoords_DA_Q1_Regular(dim, np_dir, &npoints_q, &xi)); 1129371c9d4SSatish Balay } break; 1139371c9d4SSatish Balay default: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "A valid DMSwarmPIC layout must be provided"); 1140e2ec84fSDave May } 1150e2ec84fSDave May 1169566063dSJacob Faibussowitsch PetscCall(DMDAGetElements(dmc, &nel, &npe, &element_list)); 1179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(dim * npe, &elcoor)); 1189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npoints_q, &basis)); 1190e2ec84fSDave May for (q = 0; q < npoints_q; q++) { 1209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(npe, &basis[q])); 1210e2ec84fSDave May 1220e2ec84fSDave May switch (dim) { 1230e2ec84fSDave May case 1: 1240e2ec84fSDave May basis[q][0] = 0.5 * (1.0 - xi[dim * q + 0]); 1250e2ec84fSDave May basis[q][1] = 0.5 * (1.0 + xi[dim * q + 0]); 1260e2ec84fSDave May break; 1270e2ec84fSDave May case 2: 1280e2ec84fSDave May basis[q][0] = 0.25 * (1.0 - xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]); 1290e2ec84fSDave May basis[q][1] = 0.25 * (1.0 + xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]); 1300e2ec84fSDave May basis[q][2] = 0.25 * (1.0 + xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]); 1310e2ec84fSDave May basis[q][3] = 0.25 * (1.0 - xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]); 1320e2ec84fSDave May break; 1330e2ec84fSDave May 1340e2ec84fSDave May case 3: 1350e2ec84fSDave May basis[q][0] = 0.125 * (1.0 - xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]) * (1.0 - xi[dim * q + 2]); 1360e2ec84fSDave May basis[q][1] = 0.125 * (1.0 + xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]) * (1.0 - xi[dim * q + 2]); 1370e2ec84fSDave May basis[q][2] = 0.125 * (1.0 + xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]) * (1.0 - xi[dim * q + 2]); 1380e2ec84fSDave May basis[q][3] = 0.125 * (1.0 - xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]) * (1.0 - xi[dim * q + 2]); 1390e2ec84fSDave May basis[q][4] = 0.125 * (1.0 - xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]) * (1.0 + xi[dim * q + 2]); 1400e2ec84fSDave May basis[q][5] = 0.125 * (1.0 + xi[dim * q + 0]) * (1.0 - xi[dim * q + 1]) * (1.0 + xi[dim * q + 2]); 1410e2ec84fSDave May basis[q][6] = 0.125 * (1.0 + xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]) * (1.0 + xi[dim * q + 2]); 1420e2ec84fSDave May basis[q][7] = 0.125 * (1.0 - xi[dim * q + 0]) * (1.0 + xi[dim * q + 1]) * (1.0 + xi[dim * q + 2]); 1430e2ec84fSDave May break; 1440e2ec84fSDave May } 1450e2ec84fSDave May } 1460e2ec84fSDave May 1479566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dm, npoints_q * nel, -1)); 1489566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 1499566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 1500e2ec84fSDave May 1519566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dmc, &coor)); 1529566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coor, &_coor)); 1530e2ec84fSDave May pcnt = 0; 1540e2ec84fSDave May for (e = 0; e < nel; e++) { 1550e2ec84fSDave May const PetscInt *element = &element_list[npe * e]; 1560e2ec84fSDave May 1570e2ec84fSDave May for (k = 0; k < npe; k++) { 1589371c9d4SSatish Balay for (d = 0; d < dim; d++) { elcoor[dim * k + d] = PetscRealPart(_coor[dim * element[k] + d]); } 1590e2ec84fSDave May } 1600e2ec84fSDave May 1610e2ec84fSDave May for (q = 0; q < npoints_q; q++) { 1629371c9d4SSatish Balay for (d = 0; d < dim; d++) { swarm_coor[dim * pcnt + d] = 0.0; } 1630e2ec84fSDave May for (k = 0; k < npe; k++) { 1649371c9d4SSatish Balay for (d = 0; d < dim; d++) { swarm_coor[dim * pcnt + d] += basis[q][k] * elcoor[dim * k + d]; } 1650e2ec84fSDave May } 1660e2ec84fSDave May swarm_cellid[pcnt] = e; 1670e2ec84fSDave May pcnt++; 1680e2ec84fSDave May } 1690e2ec84fSDave May } 1709566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coor, &_coor)); 1719566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_cellid, NULL, NULL, (void **)&swarm_cellid)); 1729566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dm, DMSwarmPICField_coor, NULL, NULL, (void **)&swarm_coor)); 1739566063dSJacob Faibussowitsch PetscCall(DMDARestoreElements(dmc, &nel, &npe, &element_list)); 1740e2ec84fSDave May 1759566063dSJacob Faibussowitsch PetscCall(PetscFree(xi)); 1769566063dSJacob Faibussowitsch PetscCall(PetscFree(elcoor)); 177*48a46eb9SPierre Jolivet for (q = 0; q < npoints_q; q++) PetscCall(PetscFree(basis[q])); 1789566063dSJacob Faibussowitsch PetscCall(PetscFree(basis)); 1790e2ec84fSDave May PetscFunctionReturn(0); 1800e2ec84fSDave May } 1810e2ec84fSDave May 1829371c9d4SSatish Balay PetscErrorCode private_DMSwarmInsertPointsUsingCellDM_DA(DM dm, DM celldm, DMSwarmPICLayoutType layout, PetscInt layout_param) { 1830e2ec84fSDave May DMDAElementType etype; 1840e2ec84fSDave May PetscInt dim; 1850e2ec84fSDave May 1860e2ec84fSDave May PetscFunctionBegin; 1879566063dSJacob Faibussowitsch PetscCall(DMDAGetElementType(celldm, &etype)); 1889566063dSJacob Faibussowitsch PetscCall(DMGetDimension(celldm, &dim)); 1890e2ec84fSDave May switch (etype) { 1909371c9d4SSatish Balay case DMDA_ELEMENT_P1: SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DA support is not currently available for DMDA_ELEMENT_P1"); 1910e2ec84fSDave May case DMDA_ELEMENT_Q1: 19208401ef6SPierre Jolivet PetscCheck(dim != 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Support only available for dim = 2, 3"); 1939566063dSJacob Faibussowitsch PetscCall(private_DMSwarmInsertPointsUsingCellDM_DA_Q1(dm, celldm, layout_param, layout)); 1940e2ec84fSDave May break; 1950e2ec84fSDave May } 1960e2ec84fSDave May PetscFunctionReturn(0); 1970e2ec84fSDave May } 1981f52046fSDave May 1999371c9d4SSatish Balay PetscErrorCode DMSwarmProjectField_ApproxQ1_DA_2D(DM swarm, PetscReal *swarm_field, DM dm, Vec v_field) { 2001f52046fSDave May Vec v_field_l, denom_l, coor_l, denom; 201a5f152d1SDave May PetscScalar *_field_l, *_denom_l; 2021f52046fSDave May PetscInt k, p, e, npoints, nel, npe; 2031f52046fSDave May PetscInt *mpfield_cell; 2041f52046fSDave May PetscReal *mpfield_coor; 2051f52046fSDave May const PetscInt *element_list; 2061f52046fSDave May const PetscInt *element; 207a5f152d1SDave May PetscScalar xi_p[2], Ni[4]; 208a5f152d1SDave May const PetscScalar *_coor; 2091f52046fSDave May 2101f52046fSDave May PetscFunctionBegin; 2119566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(v_field)); 2121f52046fSDave May 2139566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &v_field_l)); 2149566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &denom)); 2159566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &denom_l)); 2169566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(v_field_l)); 2179566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(denom)); 2189566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(denom_l)); 2191f52046fSDave May 2209566063dSJacob Faibussowitsch PetscCall(VecGetArray(v_field_l, &_field_l)); 2219566063dSJacob Faibussowitsch PetscCall(VecGetArray(denom_l, &_denom_l)); 2221f52046fSDave May 2239566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coor_l)); 2249566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(coor_l, &_coor)); 2251f52046fSDave May 2269566063dSJacob Faibussowitsch PetscCall(DMDAGetElements(dm, &nel, &npe, &element_list)); 2279566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(swarm, &npoints)); 2289566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&mpfield_coor)); 2299566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(swarm, DMSwarmPICField_cellid, NULL, NULL, (void **)&mpfield_cell)); 2301f52046fSDave May 2311f52046fSDave May for (p = 0; p < npoints; p++) { 2321f52046fSDave May PetscReal *coor_p; 233a5f152d1SDave May const PetscScalar *x0; 234a5f152d1SDave May const PetscScalar *x2; 235a5f152d1SDave May PetscScalar dx[2]; 2361f52046fSDave May 2371f52046fSDave May e = mpfield_cell[p]; 2381f52046fSDave May coor_p = &mpfield_coor[2 * p]; 2391f52046fSDave May element = &element_list[npe * e]; 2401f52046fSDave May 2411f52046fSDave May /* compute local coordinates: (xp-x0)/dx = (xip+1)/2 */ 2421f52046fSDave May x0 = &_coor[2 * element[0]]; 2431f52046fSDave May x2 = &_coor[2 * element[2]]; 2441f52046fSDave May 2451f52046fSDave May dx[0] = x2[0] - x0[0]; 2461f52046fSDave May dx[1] = x2[1] - x0[1]; 2471f52046fSDave May 2481f52046fSDave May xi_p[0] = 2.0 * (coor_p[0] - x0[0]) / dx[0] - 1.0; 2491f52046fSDave May xi_p[1] = 2.0 * (coor_p[1] - x0[1]) / dx[1] - 1.0; 2501f52046fSDave May 2511f52046fSDave May /* evaluate basis functions */ 2521f52046fSDave May Ni[0] = 0.25 * (1.0 - xi_p[0]) * (1.0 - xi_p[1]); 2531f52046fSDave May Ni[1] = 0.25 * (1.0 + xi_p[0]) * (1.0 - xi_p[1]); 2541f52046fSDave May Ni[2] = 0.25 * (1.0 + xi_p[0]) * (1.0 + xi_p[1]); 2551f52046fSDave May Ni[3] = 0.25 * (1.0 - xi_p[0]) * (1.0 + xi_p[1]); 2561f52046fSDave May 2571f52046fSDave May for (k = 0; k < npe; k++) { 2581f52046fSDave May _field_l[element[k]] += Ni[k] * swarm_field[p]; 2591f52046fSDave May _denom_l[element[k]] += Ni[k]; 2601f52046fSDave May } 2611f52046fSDave May } 2621f52046fSDave May 2639566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(swarm, DMSwarmPICField_cellid, NULL, NULL, (void **)&mpfield_cell)); 2649566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(swarm, DMSwarmPICField_coor, NULL, NULL, (void **)&mpfield_coor)); 2659566063dSJacob Faibussowitsch PetscCall(DMDARestoreElements(dm, &nel, &npe, &element_list)); 2669566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(coor_l, &_coor)); 2679566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v_field_l, &_field_l)); 2689566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(denom_l, &_denom_l)); 2691f52046fSDave May 2709566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, v_field_l, ADD_VALUES, v_field)); 2719566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, v_field_l, ADD_VALUES, v_field)); 2729566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, denom_l, ADD_VALUES, denom)); 2739566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, denom_l, ADD_VALUES, denom)); 2741f52046fSDave May 2759566063dSJacob Faibussowitsch PetscCall(VecPointwiseDivide(v_field, v_field, denom)); 2761f52046fSDave May 2779566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &v_field_l)); 2789566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &denom_l)); 2799566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &denom)); 2801f52046fSDave May PetscFunctionReturn(0); 2811f52046fSDave May } 2821f52046fSDave May 2839371c9d4SSatish Balay PetscErrorCode private_DMSwarmProjectFields_DA(DM swarm, DM celldm, PetscInt project_type, PetscInt nfields, DMSwarmDataField dfield[], Vec vecs[]) { 2841f52046fSDave May PetscInt f, dim; 2851f52046fSDave May DMDAElementType etype; 2861f52046fSDave May 2871f52046fSDave May PetscFunctionBegin; 2889566063dSJacob Faibussowitsch PetscCall(DMDAGetElementType(celldm, &etype)); 28908401ef6SPierre Jolivet PetscCheck(etype != DMDA_ELEMENT_P1, PetscObjectComm((PetscObject)swarm), PETSC_ERR_SUP, "Only Q1 DMDA supported"); 2901f52046fSDave May 2919566063dSJacob Faibussowitsch PetscCall(DMGetDimension(swarm, &dim)); 2921f52046fSDave May switch (dim) { 2931f52046fSDave May case 2: 2941f52046fSDave May for (f = 0; f < nfields; f++) { 2951f52046fSDave May PetscReal *swarm_field; 2961f52046fSDave May 2979566063dSJacob Faibussowitsch PetscCall(DMSwarmDataFieldGetEntries(dfield[f], (void **)&swarm_field)); 2989566063dSJacob Faibussowitsch PetscCall(DMSwarmProjectField_ApproxQ1_DA_2D(swarm, swarm_field, celldm, vecs[f])); 2991f52046fSDave May } 3001f52046fSDave May break; 3019371c9d4SSatish Balay case 3: SETERRQ(PetscObjectComm((PetscObject)swarm), PETSC_ERR_SUP, "No support for 3D"); 3029371c9d4SSatish Balay default: break; 3031f52046fSDave May } 3041f52046fSDave May PetscFunctionReturn(0); 3051f52046fSDave May } 306