1*47d993e7Ssuyashtn /* Portions of this code are under: 2*47d993e7Ssuyashtn Copyright (C) 2022 Advanced Micro Devices, Inc. All rights reserved. 3*47d993e7Ssuyashtn */ 42a5e0493SMark Adams static char help[] = "3D, tensor hexahedra (Q1-K), displacement finite element formulation\n\ 5c4762a1bSJed Brown of linear elasticity. E=1.0, nu=1/3.\n\ 6c4762a1bSJed Brown Unit cube domain with Dirichlet boundary\n\n"; 7c4762a1bSJed Brown 8c4762a1bSJed Brown #include <petscdmplex.h> 9c4762a1bSJed Brown #include <petscsnes.h> 10c4762a1bSJed Brown #include <petscds.h> 11c4762a1bSJed Brown #include <petscdmforest.h> 12c4762a1bSJed Brown 13c4762a1bSJed Brown static PetscReal s_soft_alpha = 1.e-3; 14c4762a1bSJed Brown static PetscReal s_mu = 0.4; 15c4762a1bSJed Brown static PetscReal s_lambda = 0.4; 16c4762a1bSJed Brown 17d71ae5a4SJacob Faibussowitsch static void f0_bd_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 18d71ae5a4SJacob Faibussowitsch { 19c4762a1bSJed Brown f0[0] = 1; /* x direction pull */ 20c4762a1bSJed Brown f0[1] = -x[2]; /* add a twist around x-axis */ 21c4762a1bSJed Brown f0[2] = x[1]; 22c4762a1bSJed Brown } 23c4762a1bSJed Brown 24d71ae5a4SJacob Faibussowitsch static void f1_bd_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], const PetscReal n[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 25d71ae5a4SJacob Faibussowitsch { 26c4762a1bSJed Brown const PetscInt Ncomp = dim; 27c4762a1bSJed Brown PetscInt comp, d; 28c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) { 29ad540459SPierre Jolivet for (d = 0; d < dim; ++d) f1[comp * dim + d] = 0.0; 30c4762a1bSJed Brown } 31c4762a1bSJed Brown } 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */ 34d71ae5a4SJacob Faibussowitsch static void f1_u_3d_alpha(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 35d71ae5a4SJacob Faibussowitsch { 36c4762a1bSJed Brown PetscReal trace, mu = s_mu, lambda = s_lambda, rad; 37c4762a1bSJed Brown PetscInt i, j; 38c4762a1bSJed Brown for (i = 0, rad = 0.; i < dim; i++) { 39c4762a1bSJed Brown PetscReal t = x[i]; 40c4762a1bSJed Brown rad += t * t; 41c4762a1bSJed Brown } 42c4762a1bSJed Brown rad = PetscSqrtReal(rad); 43c4762a1bSJed Brown if (rad > 0.25) { 44c4762a1bSJed Brown mu *= s_soft_alpha; 45c4762a1bSJed Brown lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */ 46c4762a1bSJed Brown } 47ad540459SPierre Jolivet for (i = 0, trace = 0; i < dim; ++i) trace += PetscRealPart(u_x[i * dim + i]); 48c4762a1bSJed Brown for (i = 0; i < dim; ++i) { 49ad540459SPierre Jolivet for (j = 0; j < dim; ++j) f1[i * dim + j] = mu * (u_x[i * dim + j] + u_x[j * dim + i]); 50c4762a1bSJed Brown f1[i * dim + i] += lambda * trace; 51c4762a1bSJed Brown } 52c4762a1bSJed Brown } 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */ 55d71ae5a4SJacob Faibussowitsch static void f1_u_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[]) 56d71ae5a4SJacob Faibussowitsch { 57c4762a1bSJed Brown PetscReal trace, mu = s_mu, lambda = s_lambda; 58c4762a1bSJed Brown PetscInt i, j; 59ad540459SPierre Jolivet for (i = 0, trace = 0; i < dim; ++i) trace += PetscRealPart(u_x[i * dim + i]); 60c4762a1bSJed Brown for (i = 0; i < dim; ++i) { 61ad540459SPierre Jolivet for (j = 0; j < dim; ++j) f1[i * dim + j] = mu * (u_x[i * dim + j] + u_x[j * dim + i]); 62c4762a1bSJed Brown f1[i * dim + i] += lambda * trace; 63c4762a1bSJed Brown } 64c4762a1bSJed Brown } 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* 3D elasticity */ 67c4762a1bSJed Brown #define IDX(ii, jj, kk, ll) (27 * ii + 9 * jj + 3 * kk + ll) 68c4762a1bSJed Brown 69d71ae5a4SJacob Faibussowitsch void g3_uu_3d_private(PetscScalar g3[], const PetscReal mu, const PetscReal lambda) 70d71ae5a4SJacob Faibussowitsch { 71c4762a1bSJed Brown if (1) { 72c4762a1bSJed Brown g3[0] += lambda; 73c4762a1bSJed Brown g3[0] += mu; 74c4762a1bSJed Brown g3[0] += mu; 75c4762a1bSJed Brown g3[4] += lambda; 76c4762a1bSJed Brown g3[8] += lambda; 77c4762a1bSJed Brown g3[10] += mu; 78c4762a1bSJed Brown g3[12] += mu; 79c4762a1bSJed Brown g3[20] += mu; 80c4762a1bSJed Brown g3[24] += mu; 81c4762a1bSJed Brown g3[28] += mu; 82c4762a1bSJed Brown g3[30] += mu; 83c4762a1bSJed Brown g3[36] += lambda; 84c4762a1bSJed Brown g3[40] += lambda; 85c4762a1bSJed Brown g3[40] += mu; 86c4762a1bSJed Brown g3[40] += mu; 87c4762a1bSJed Brown g3[44] += lambda; 88c4762a1bSJed Brown g3[50] += mu; 89c4762a1bSJed Brown g3[52] += mu; 90c4762a1bSJed Brown g3[56] += mu; 91c4762a1bSJed Brown g3[60] += mu; 92c4762a1bSJed Brown g3[68] += mu; 93c4762a1bSJed Brown g3[70] += mu; 94c4762a1bSJed Brown g3[72] += lambda; 95c4762a1bSJed Brown g3[76] += lambda; 96c4762a1bSJed Brown g3[80] += lambda; 97c4762a1bSJed Brown g3[80] += mu; 98c4762a1bSJed Brown g3[80] += mu; 99c4762a1bSJed Brown } else { 100c4762a1bSJed Brown int i, j, k, l; 101c4762a1bSJed Brown static int cc = -1; 102c4762a1bSJed Brown cc++; 103c4762a1bSJed Brown for (i = 0; i < 3; ++i) { 104c4762a1bSJed Brown for (j = 0; j < 3; ++j) { 105c4762a1bSJed Brown for (k = 0; k < 3; ++k) { 106c4762a1bSJed Brown for (l = 0; l < 3; ++l) { 107c4762a1bSJed Brown if (k == l && i == j) g3[IDX(i, j, k, l)] += lambda; 108c4762a1bSJed Brown if (i == k && j == l) g3[IDX(i, j, k, l)] += mu; 109c4762a1bSJed Brown if (i == l && j == k) g3[IDX(i, j, k, l)] += mu; 110c4762a1bSJed Brown if (k == l && i == j && !cc) (void)PetscPrintf(PETSC_COMM_WORLD, "g3[%d] += lambda;\n", IDX(i, j, k, l)); 111c4762a1bSJed Brown if (i == k && j == l && !cc) (void)PetscPrintf(PETSC_COMM_WORLD, "g3[%d] += mu;\n", IDX(i, j, k, l)); 112c4762a1bSJed Brown if (i == l && j == k && !cc) (void)PetscPrintf(PETSC_COMM_WORLD, "g3[%d] += mu;\n", IDX(i, j, k, l)); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown } 115c4762a1bSJed Brown } 116c4762a1bSJed Brown } 117c4762a1bSJed Brown } 118c4762a1bSJed Brown } 119c4762a1bSJed Brown 120d71ae5a4SJacob Faibussowitsch static void g3_uu_3d_alpha(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 121d71ae5a4SJacob Faibussowitsch { 122c4762a1bSJed Brown PetscReal mu = s_mu, lambda = s_lambda, rad; 123c4762a1bSJed Brown PetscInt i; 124c4762a1bSJed Brown for (i = 0, rad = 0.; i < dim; i++) { 125c4762a1bSJed Brown PetscReal t = x[i]; 126c4762a1bSJed Brown rad += t * t; 127c4762a1bSJed Brown } 128c4762a1bSJed Brown rad = PetscSqrtReal(rad); 129c4762a1bSJed Brown if (rad > 0.25) { 130c4762a1bSJed Brown mu *= s_soft_alpha; 131c4762a1bSJed Brown lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */ 132c4762a1bSJed Brown } 133c4762a1bSJed Brown g3_uu_3d_private(g3, mu, lambda); 134c4762a1bSJed Brown } 135c4762a1bSJed Brown 136d71ae5a4SJacob Faibussowitsch static void g3_uu_3d(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[]) 137d71ae5a4SJacob Faibussowitsch { 138c4762a1bSJed Brown g3_uu_3d_private(g3, s_mu, s_lambda); 139c4762a1bSJed Brown } 140c4762a1bSJed Brown 141d71ae5a4SJacob Faibussowitsch static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 142d71ae5a4SJacob Faibussowitsch { 143c4762a1bSJed Brown const PetscInt Ncomp = dim; 144c4762a1bSJed Brown PetscInt comp; 145c4762a1bSJed Brown 146c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) f0[comp] = 0.0; 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* PI_i (x_i^4 - x_i^2) */ 150d71ae5a4SJacob Faibussowitsch static void f0_u_x4(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[]) 151d71ae5a4SJacob Faibussowitsch { 152c4762a1bSJed Brown const PetscInt Ncomp = dim; 153c4762a1bSJed Brown PetscInt comp, i; 154c4762a1bSJed Brown 155c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) { 156c4762a1bSJed Brown f0[comp] = 1e5; 1579371c9d4SSatish Balay for (i = 0; i < Ncomp; ++i) { f0[comp] *= /* (comp+1)* */ (x[i] * x[i] * x[i] * x[i] - x[i] * x[i]); /* assumes (0,1]^D domain */ } 158c4762a1bSJed Brown } 159c4762a1bSJed Brown } 160c4762a1bSJed Brown 161d71ae5a4SJacob Faibussowitsch PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 162d71ae5a4SJacob Faibussowitsch { 163c4762a1bSJed Brown const PetscInt Ncomp = dim; 164c4762a1bSJed Brown PetscInt comp; 165c4762a1bSJed Brown 166c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) u[comp] = 0; 167c4762a1bSJed Brown return 0; 168c4762a1bSJed Brown } 169c4762a1bSJed Brown 170d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 171d71ae5a4SJacob Faibussowitsch { 172c4762a1bSJed Brown Mat Amat; 173c4762a1bSJed Brown SNES snes; 174c4762a1bSJed Brown KSP ksp; 175c4762a1bSJed Brown MPI_Comm comm; 176c4762a1bSJed Brown PetscMPIInt rank; 177956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 178c4762a1bSJed Brown PetscLogStage stage[17]; 179956f8c0dSBarry Smith #endif 180c4762a1bSJed Brown PetscBool test_nonzero_cols = PETSC_FALSE, use_nearnullspace = PETSC_TRUE, attach_nearnullspace = PETSC_FALSE; 181c4762a1bSJed Brown Vec xx, bb; 182c4762a1bSJed Brown PetscInt iter, i, N, dim = 3, cells[3] = {1, 1, 1}, max_conv_its, local_sizes[7], run_type = 1; 183c4762a1bSJed Brown DM dm, distdm, basedm; 184c4762a1bSJed Brown PetscBool flg; 185c4762a1bSJed Brown char convType[256]; 186c4762a1bSJed Brown PetscReal Lx, mdisp[10], err[10]; 1879371c9d4SSatish Balay const char *const options[10] = {"-ex56_dm_refine 0", "-ex56_dm_refine 1", "-ex56_dm_refine 2", "-ex56_dm_refine 3", "-ex56_dm_refine 4", "-ex56_dm_refine 5", "-ex56_dm_refine 6", "-ex56_dm_refine 7", "-ex56_dm_refine 8", "-ex56_dm_refine 9"}; 188c4762a1bSJed Brown PetscFunctionBeginUser; 189327415f7SBarry Smith PetscFunctionBeginUser; 1909566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 191c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 1929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 193c4762a1bSJed Brown /* options */ 194d0609cedSBarry Smith PetscOptionsBegin(comm, NULL, "3D bilinear Q1 elasticity options", ""); 195c4762a1bSJed Brown { 196c4762a1bSJed Brown i = 3; 1979566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-cells", "Number of (flux tube) processor in each dimension", "ex56.c", cells, &i, NULL)); 198c4762a1bSJed Brown 199c4762a1bSJed Brown Lx = 1.; /* or ne for rod */ 200c4762a1bSJed Brown max_conv_its = 3; 2019566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-max_conv_its", "Number of iterations in convergence study", "", max_conv_its, &max_conv_its, NULL)); 20263a3b9bcSJacob Faibussowitsch PetscCheck(max_conv_its > 0 && max_conv_its < 7, PETSC_COMM_WORLD, PETSC_ERR_USER, "Bad number of iterations for convergence test (%" PetscInt_FMT ")", max_conv_its); 2039566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-lx", "Length of domain", "", Lx, &Lx, NULL)); 2049566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-alpha", "material coefficient inside circle", "", s_soft_alpha, &s_soft_alpha, NULL)); 2059566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_nonzero_cols", "nonzero test", "", test_nonzero_cols, &test_nonzero_cols, NULL)); 2069566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_mat_nearnullspace", "MatNearNullSpace API test", "", use_nearnullspace, &use_nearnullspace, NULL)); 2079566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-attach_mat_nearnullspace", "MatNearNullSpace API test (via MatSetNearNullSpace)", "", attach_nearnullspace, &attach_nearnullspace, NULL)); 2089566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-run_type", "0: twisting load on cantalever, 1: 3rd order accurate convergence test", "", run_type, &run_type, NULL)); 209c4762a1bSJed Brown } 210d0609cedSBarry Smith PetscOptionsEnd(); 2119566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Mesh Setup", &stage[16])); 212c4762a1bSJed Brown for (iter = 0; iter < max_conv_its; iter++) { 213c4762a1bSJed Brown char str[] = "Solve 0"; 214c4762a1bSJed Brown str[6] += iter; 2159566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister(str, &stage[iter])); 216c4762a1bSJed Brown } 217c4762a1bSJed Brown /* create DM, Plex calls DMSetup */ 2189566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[16])); 2199566063dSJacob Faibussowitsch PetscCall(DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, cells, NULL, NULL, NULL, PETSC_TRUE, &dm)); 220c4762a1bSJed Brown { 221c4762a1bSJed Brown DMLabel label; 222c4762a1bSJed Brown IS is; 2239566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "boundary")); 2249566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "boundary", &label)); 2259566063dSJacob Faibussowitsch PetscCall(DMPlexMarkBoundaryFaces(dm, 1, label)); 2262a5e0493SMark Adams if (run_type == 0) { 2279566063dSJacob Faibussowitsch PetscCall(DMGetStratumIS(dm, "boundary", 1, &is)); 2289566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "Faces")); 229c4762a1bSJed Brown if (is) { 230c4762a1bSJed Brown PetscInt d, f, Nf; 231c4762a1bSJed Brown const PetscInt *faces; 232c4762a1bSJed Brown PetscInt csize; 233c4762a1bSJed Brown PetscSection cs; 234c4762a1bSJed Brown Vec coordinates; 235c4762a1bSJed Brown DM cdm; 2369566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &Nf)); 2379566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &faces)); 2389566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 2399566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 2409566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(cdm, &cs)); 241c4762a1bSJed Brown /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */ 242c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 243c4762a1bSJed Brown PetscReal faceCoord; 244c4762a1bSJed Brown PetscInt b, v; 245c4762a1bSJed Brown PetscScalar *coords = NULL; 246c4762a1bSJed Brown PetscInt Nv; 2479566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords)); 248c4762a1bSJed Brown Nv = csize / dim; /* Calculate mean coordinate vector */ 249c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 250c4762a1bSJed Brown faceCoord = 0.0; 251c4762a1bSJed Brown for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v * dim + d]); 252c4762a1bSJed Brown faceCoord /= Nv; 253c4762a1bSJed Brown for (b = 0; b < 2; ++b) { 254c4762a1bSJed Brown if (PetscAbs(faceCoord - b) < PETSC_SMALL) { /* domain have not been set yet, still [0,1]^3 */ 2559566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(dm, "Faces", faces[f], d * 2 + b + 1)); 256c4762a1bSJed Brown } 257c4762a1bSJed Brown } 258c4762a1bSJed Brown } 2599566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords)); 260c4762a1bSJed Brown } 2619566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &faces)); 262c4762a1bSJed Brown } 2639566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 2649566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Faces", &label)); 2659566063dSJacob Faibussowitsch PetscCall(DMPlexLabelComplete(dm, label)); 266c4762a1bSJed Brown } 267c4762a1bSJed Brown } 268c4762a1bSJed Brown { 269c4762a1bSJed Brown PetscInt dimEmbed, i; 270c4762a1bSJed Brown PetscInt nCoords; 2719371c9d4SSatish Balay PetscScalar *coords, bounds[] = { 2729371c9d4SSatish Balay 0, 1, -.5, .5, -.5, .5, 2739371c9d4SSatish Balay }; /* x_min,x_max,y_min,y_max */ 274c4762a1bSJed Brown Vec coordinates; 275c4762a1bSJed Brown bounds[1] = Lx; 276c4762a1bSJed Brown if (run_type == 1) { 277c4762a1bSJed Brown for (i = 0; i < 2 * dim; i++) bounds[i] = (i % 2) ? 1 : 0; 278c4762a1bSJed Brown } 2799566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 2809566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(dm, &dimEmbed)); 28163a3b9bcSJacob Faibussowitsch PetscCheck(dimEmbed == dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "dimEmbed != dim %" PetscInt_FMT, dimEmbed); 2829566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(coordinates, &nCoords)); 283e00437b9SBarry Smith PetscCheck((nCoords % dimEmbed) == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Coordinate vector the wrong size"); 2849566063dSJacob Faibussowitsch PetscCall(VecGetArray(coordinates, &coords)); 285c4762a1bSJed Brown for (i = 0; i < nCoords; i += dimEmbed) { 286c4762a1bSJed Brown PetscInt j; 287c4762a1bSJed Brown PetscScalar *coord = &coords[i]; 288ad540459SPierre Jolivet for (j = 0; j < dimEmbed; j++) coord[j] = bounds[2 * j] + coord[j] * (bounds[2 * j + 1] - bounds[2 * j]); 289c4762a1bSJed Brown } 2909566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(coordinates, &coords)); 2919566063dSJacob Faibussowitsch PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 292c4762a1bSJed Brown } 293c4762a1bSJed Brown 294c4762a1bSJed Brown /* convert to p4est, and distribute */ 295d0609cedSBarry Smith PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX"); 2969566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-dm_type", "Convert DMPlex to another format (should not be Plex!)", "ex56.c", DMList, DMPLEX, convType, 256, &flg)); 297d0609cedSBarry Smith PetscOptionsEnd(); 298c4762a1bSJed Brown if (flg) { 299c4762a1bSJed Brown DM newdm; 3009566063dSJacob Faibussowitsch PetscCall(DMConvert(dm, convType, &newdm)); 301c4762a1bSJed Brown if (newdm) { 302c4762a1bSJed Brown const char *prefix; 303c4762a1bSJed Brown PetscBool isForest; 3049566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 3059566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)newdm, prefix)); 3069566063dSJacob Faibussowitsch PetscCall(DMIsForest(newdm, &isForest)); 30728b400f6SJacob Faibussowitsch PetscCheck(isForest, PETSC_COMM_WORLD, PETSC_ERR_USER, "Converted to non Forest?"); 3089566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 309c4762a1bSJed Brown dm = newdm; 310c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Convert failed?"); 311c4762a1bSJed Brown } else { 312c4762a1bSJed Brown PetscPartitioner part; 313c4762a1bSJed Brown /* Plex Distribute mesh over processes */ 3149566063dSJacob Faibussowitsch PetscCall(DMPlexGetPartitioner(dm, &part)); 3159566063dSJacob Faibussowitsch PetscCall(PetscPartitionerSetFromOptions(part)); 3169566063dSJacob Faibussowitsch PetscCall(DMPlexDistribute(dm, 0, NULL, &distdm)); 317c4762a1bSJed Brown if (distdm) { 318c4762a1bSJed Brown const char *prefix; 3199566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 3209566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)distdm, prefix)); 3219566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 322c4762a1bSJed Brown dm = distdm; 323c4762a1bSJed Brown } 324c4762a1bSJed Brown } 3259566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 3269371c9d4SSatish Balay basedm = dm; 3279371c9d4SSatish Balay dm = NULL; 328c4762a1bSJed Brown 329c4762a1bSJed Brown for (iter = 0; iter < max_conv_its; iter++) { 3309566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[16])); 331c4762a1bSJed Brown /* make new DM */ 3329566063dSJacob Faibussowitsch PetscCall(DMClone(basedm, &dm)); 3339566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, "ex56_")); 3349566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dm, "Mesh")); 335c4762a1bSJed Brown if (max_conv_its > 1) { 3360e75e42fSMark /* If max_conv_its == 1, then we are not doing a convergence study. */ 3379566063dSJacob Faibussowitsch PetscCall(PetscOptionsInsertString(NULL, options[iter])); 338c4762a1bSJed Brown } 3399566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(dm)); /* refinement done here in Plex, p4est */ 340c4762a1bSJed Brown /* snes */ 3419566063dSJacob Faibussowitsch PetscCall(SNESCreate(comm, &snes)); 3429566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 343c4762a1bSJed Brown /* fem */ 344c4762a1bSJed Brown { 345c4762a1bSJed Brown const PetscInt Ncomp = dim; 346c4762a1bSJed Brown const PetscInt components[] = {0, 1, 2}; 347c4762a1bSJed Brown const PetscInt Nfid = 1, Npid = 1; 348c4762a1bSJed Brown const PetscInt fid[] = {1}; /* The fixed faces (x=0) */ 349c4762a1bSJed Brown const PetscInt pid[] = {2}; /* The faces with loading (x=L_x) */ 350c4762a1bSJed Brown PetscFE fe; 351c4762a1bSJed Brown PetscDS prob; 35245480ffeSMatthew G. Knepley DMLabel label; 353c4762a1bSJed Brown DM cdm = dm; 354c4762a1bSJed Brown 3559566063dSJacob Faibussowitsch PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, dim, PETSC_FALSE, NULL, PETSC_DECIDE, &fe)); /* elasticity */ 3569566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, "deformation")); 357c4762a1bSJed Brown /* FEM prob */ 3589566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 3599566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 3609566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 361c4762a1bSJed Brown /* setup problem */ 362c4762a1bSJed Brown if (run_type == 1) { 3639566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d)); 3649566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_u_x4, f1_u_3d)); 365c4762a1bSJed Brown } else { 36645480ffeSMatthew G. Knepley PetscWeakForm wf; 36745480ffeSMatthew G. Knepley PetscInt bd, i; 36845480ffeSMatthew G. Knepley 3699566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d_alpha)); 3709566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_u, f1_u_3d_alpha)); 37145480ffeSMatthew G. Knepley 3729566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Faces", &label)); 3739566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "traction", label, Npid, pid, 0, Ncomp, components, NULL, NULL, NULL, &bd)); 3749566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(prob, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 3759566063dSJacob Faibussowitsch for (i = 0; i < Npid; ++i) PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, pid[i], 0, 0, 0, f0_bd_u_3d, 0, f1_bd_u)); 376c4762a1bSJed Brown } 377c4762a1bSJed Brown /* bcs */ 378c4762a1bSJed Brown if (run_type == 1) { 379c4762a1bSJed Brown PetscInt id = 1; 3809566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "boundary", &label)); 3819566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))zero, NULL, NULL, NULL)); 382c4762a1bSJed Brown } else { 3839566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Faces", &label)); 3849566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed", label, Nfid, fid, 0, Ncomp, components, (void (*)(void))zero, NULL, NULL, NULL)); 385c4762a1bSJed Brown } 386c4762a1bSJed Brown while (cdm) { 3879566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, cdm)); 3889566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(cdm, &cdm)); 389c4762a1bSJed Brown } 3909566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 391c4762a1bSJed Brown } 392c4762a1bSJed Brown /* vecs & mat */ 3939566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &xx)); 3949566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &bb)); 3959566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)bb, "b")); 3969566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)xx, "u")); 3979566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &Amat)); 3989566063dSJacob Faibussowitsch PetscCall(MatSetOption(Amat, MAT_SYMMETRIC, PETSC_TRUE)); /* Some matrix kernels can take advantage of symmetry if we set this. */ 3999566063dSJacob Faibussowitsch PetscCall(MatSetOption(Amat, MAT_SYMMETRY_ETERNAL, PETSC_TRUE)); /* Inform PETSc that Amat is always symmetric, so info set above isn't lost. */ 4009566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(Amat, 3)); 4019566063dSJacob Faibussowitsch PetscCall(MatSetOption(Amat, MAT_SPD, PETSC_TRUE)); 402b94d7dedSBarry Smith PetscCall(MatSetOption(Amat, MAT_SPD_ETERNAL, PETSC_TRUE)); 4039566063dSJacob Faibussowitsch PetscCall(VecGetSize(bb, &N)); 404c4762a1bSJed Brown local_sizes[iter] = N; 40563a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "%" PetscInt_FMT " global equations, %" PetscInt_FMT " vertices\n", N, N / dim)); 406c4762a1bSJed Brown if ((use_nearnullspace || attach_nearnullspace) && N / dim > 1) { 407c4762a1bSJed Brown /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */ 408c4762a1bSJed Brown DM subdm; 409c4762a1bSJed Brown MatNullSpace nearNullSpace; 410c4762a1bSJed Brown PetscInt fields = 0; 411c4762a1bSJed Brown PetscObject deformation; 4129566063dSJacob Faibussowitsch PetscCall(DMCreateSubDM(dm, 1, &fields, NULL, &subdm)); 4139566063dSJacob Faibussowitsch PetscCall(DMPlexCreateRigidBody(subdm, 0, &nearNullSpace)); 4149566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, 0, NULL, &deformation)); 4159566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose(deformation, "nearnullspace", (PetscObject)nearNullSpace)); 4169566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdm)); 4171baa6e33SBarry Smith if (attach_nearnullspace) PetscCall(MatSetNearNullSpace(Amat, nearNullSpace)); 4189566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nearNullSpace)); /* created by DM and destroyed by Mat */ 419c4762a1bSJed Brown } 4209566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, NULL, NULL, NULL)); 4219566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, Amat, Amat, NULL, NULL)); 4229566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 4239566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 4249566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 4259566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[16])); 426c4762a1bSJed Brown /* ksp */ 4279566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 4289566063dSJacob Faibussowitsch PetscCall(KSPSetComputeSingularValues(ksp, PETSC_TRUE)); 429c4762a1bSJed Brown /* test BCs */ 4309566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(xx)); 431c4762a1bSJed Brown if (test_nonzero_cols) { 43248a46eb9SPierre Jolivet if (rank == 0) PetscCall(VecSetValue(xx, 0, 1.0, INSERT_VALUES)); 4339566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(xx)); 4349566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(xx)); 435c4762a1bSJed Brown } 4369566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(bb)); 4379566063dSJacob Faibussowitsch PetscCall(VecGetSize(bb, &i)); 438c4762a1bSJed Brown local_sizes[iter] = i; 43963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "%" PetscInt_FMT " equations in vector, %" PetscInt_FMT " vertices\n", i, i / dim)); 4409566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 441c4762a1bSJed Brown /* solve */ 4429566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[iter])); 4439566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, bb, xx)); 4449566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 4459566063dSJacob Faibussowitsch PetscCall(VecNorm(xx, NORM_INFINITY, &mdisp[iter])); 4469566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 447c4762a1bSJed Brown { 448c4762a1bSJed Brown PetscViewer viewer = NULL; 449c4762a1bSJed Brown PetscViewerFormat fmt; 4509566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(comm, NULL, "ex56_", "-vec_view", &viewer, &fmt, &flg)); 451c4762a1bSJed Brown if (flg) { 4529566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, fmt)); 4539566063dSJacob Faibussowitsch PetscCall(VecView(xx, viewer)); 4549566063dSJacob Faibussowitsch PetscCall(VecView(bb, viewer)); 4559566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 456c4762a1bSJed Brown } 4579566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 458c4762a1bSJed Brown } 459c4762a1bSJed Brown /* Free work space */ 4609566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dm)); 4619566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 4629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 4639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bb)); 4649566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Amat)); 465c4762a1bSJed Brown } 4669566063dSJacob Faibussowitsch PetscCall(DMDestroy(&basedm)); 467c4762a1bSJed Brown if (run_type == 1) err[0] = 59.975208 - mdisp[0]; /* error with what I think is the exact solution */ 468c4762a1bSJed Brown else err[0] = 171.038 - mdisp[0]; 469c4762a1bSJed Brown for (iter = 1; iter < max_conv_its; iter++) { 470c4762a1bSJed Brown if (run_type == 1) err[iter] = 59.975208 - mdisp[iter]; 471c4762a1bSJed Brown else err[iter] = 171.038 - mdisp[iter]; 4729371c9d4SSatish Balay PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%d] %" PetscInt_FMT ") N=%12" PetscInt_FMT ", max displ=%9.7e, disp diff=%9.2e, error=%4.3e, rate=%3.2g\n", rank, iter, local_sizes[iter], (double)mdisp[iter], (double)(mdisp[iter] - mdisp[iter - 1]), (double)err[iter], (double)(PetscLogReal(err[iter - 1] / err[iter]) / PetscLogReal(2.)))); 473c4762a1bSJed Brown } 474c4762a1bSJed Brown 4759566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 476b122ec5aSJacob Faibussowitsch return 0; 477c4762a1bSJed Brown } 478c4762a1bSJed Brown 479c4762a1bSJed Brown /*TEST 480c4762a1bSJed Brown 481c4762a1bSJed Brown test: 482c4762a1bSJed Brown suffix: 0 483c4762a1bSJed Brown nsize: 4 484c4762a1bSJed Brown requires: !single 485bae903cbSmarkadams4 args: -cells 2,2,1 -max_conv_its 2 -petscspace_degree 3 -snes_max_it 1 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-10 -ksp_norm_type unpreconditioned -pc_type gamg -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true -pc_gamg_aggressive_coarsening 0 -pc_gamg_threshold 0.001 -ksp_converged_reason -snes_converged_reason -use_mat_nearnullspace true -mg_levels_ksp_max_it 2 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.2,0,1.1 -mg_levels_pc_type jacobi -petscpartitioner_type simple -ex56_dm_view -snes_lag_jacobian -2 -snes_type ksponly -use_gpu_aware_mpi true 486c4762a1bSJed Brown timeoutfactor: 2 487c4762a1bSJed Brown 488c4762a1bSJed Brown # HYPRE PtAP broken with complex numbers 489c4762a1bSJed Brown test: 490c4762a1bSJed Brown suffix: hypre 491263f2b91SStefano Zampini requires: hypre !single !complex !defined(PETSC_HAVE_HYPRE_DEVICE) 492c4762a1bSJed Brown nsize: 4 4930e75e42fSMark args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -pc_type hypre -pc_hypre_type boomeramg -pc_hypre_boomeramg_no_CF true -pc_hypre_boomeramg_agg_nl 1 -pc_hypre_boomeramg_coarsen_type HMIS -pc_hypre_boomeramg_interp_type ext+i -ksp_converged_reason -use_mat_nearnullspace true -petscpartitioner_type simple 494c4762a1bSJed Brown 495c4762a1bSJed Brown test: 496c4762a1bSJed Brown suffix: ml 497c4762a1bSJed Brown requires: ml !single 498c4762a1bSJed Brown nsize: 4 4990e75e42fSMark args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_converged_reason -ksp_rtol 1.e-8 -pc_type ml -mg_levels_ksp_type chebyshev -mg_levels_ksp_max_it 3 -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type sor -petscpartitioner_type simple -use_mat_nearnullspace 500c4762a1bSJed Brown 501c4762a1bSJed Brown test: 502c4762a1bSJed Brown suffix: hpddm 503dfd57a17SPierre Jolivet requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 504c4762a1bSJed Brown nsize: 4 5050e75e42fSMark args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type fgmres -ksp_monitor_short -ksp_converged_reason -ksp_rtol 1.e-8 -pc_type hpddm -petscpartitioner_type simple -pc_hpddm_levels_1_sub_pc_type lu -pc_hpddm_levels_1_eps_nev 6 -pc_hpddm_coarse_p 1 -pc_hpddm_coarse_pc_type svd 506c4762a1bSJed Brown 507c4762a1bSJed Brown test: 50863b77682SMark Adams suffix: repart 509c4762a1bSJed Brown nsize: 4 510c4762a1bSJed Brown requires: parmetis !single 511bae903cbSmarkadams4 args: -cells 8,2,2 -max_conv_its 1 -petscspace_degree 2 -snes_max_it 4 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-2 -ksp_norm_type unpreconditioned -snes_rtol 1.e-3 -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_aggressive_coarsening 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -use_mat_nearnullspace true -mg_levels_ksp_max_it 2 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type jacobi -pc_gamg_mat_partitioning_type parmetis -pc_gamg_repartition true -snes_converged_reason -pc_gamg_process_eq_limit 20 -pc_gamg_coarse_eq_limit 10 -ksp_converged_reason -snes_converged_reason -pc_gamg_reuse_interpolation true 512c4762a1bSJed Brown 513c4762a1bSJed Brown test: 514c4762a1bSJed Brown suffix: bddc 515c4762a1bSJed Brown nsize: 4 516c4762a1bSJed Brown requires: !single 5170e75e42fSMark args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -petscpartitioner_type simple -ex56_dm_mat_type is -matis_localmat_type {{sbaij baij aij}} -pc_type bddc 518c4762a1bSJed Brown 519c4762a1bSJed Brown testset: 520c4762a1bSJed Brown nsize: 4 521c4762a1bSJed Brown requires: !single 5220e75e42fSMark args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-10 -ksp_converged_reason -petscpartitioner_type simple -ex56_dm_mat_type is -matis_localmat_type aij -pc_type bddc -attach_mat_nearnullspace {{0 1}separate output} 523c4762a1bSJed Brown test: 524c4762a1bSJed Brown suffix: bddc_approx_gamg 525bae903cbSmarkadams4 args: -pc_bddc_switch_static -prefix_push pc_bddc_dirichlet_ -approximate -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_reuse_interpolation true -pc_gamg_aggressive_coarsening 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -prefix_pop -prefix_push pc_bddc_neumann_ -approximate -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true -pc_gamg_aggressive_coarsening 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -prefix_pop 526c4762a1bSJed Brown # HYPRE PtAP broken with complex numbers 527c4762a1bSJed Brown test: 528263f2b91SStefano Zampini requires: hypre !complex !defined(PETSC_HAVE_HYPRE_DEVICE) 529c4762a1bSJed Brown suffix: bddc_approx_hypre 530c4762a1bSJed Brown args: -pc_bddc_switch_static -prefix_push pc_bddc_dirichlet_ -pc_type hypre -pc_hypre_boomeramg_no_CF true -pc_hypre_boomeramg_strong_threshold 0.75 -pc_hypre_boomeramg_agg_nl 1 -pc_hypre_boomeramg_coarsen_type HMIS -pc_hypre_boomeramg_interp_type ext+i -prefix_pop -prefix_push pc_bddc_neumann_ -pc_type hypre -pc_hypre_boomeramg_no_CF true -pc_hypre_boomeramg_strong_threshold 0.75 -pc_hypre_boomeramg_agg_nl 1 -pc_hypre_boomeramg_coarsen_type HMIS -pc_hypre_boomeramg_interp_type ext+i -prefix_pop 531c4762a1bSJed Brown test: 532c4762a1bSJed Brown requires: ml 533c4762a1bSJed Brown suffix: bddc_approx_ml 5340e75e42fSMark args: -pc_bddc_switch_static -prefix_push pc_bddc_dirichlet_ -approximate -pc_type ml -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -prefix_pop -prefix_push pc_bddc_neumann_ -approximate -pc_type ml -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -prefix_pop 535c4762a1bSJed Brown 536c4762a1bSJed Brown test: 537c4762a1bSJed Brown suffix: fetidp 538c4762a1bSJed Brown nsize: 4 539c4762a1bSJed Brown requires: !single 5400e75e42fSMark args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type fetidp -fetidp_ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -petscpartitioner_type simple -ex56_dm_mat_type is -matis_localmat_type {{sbaij baij aij}} 541c4762a1bSJed Brown 542c4762a1bSJed Brown test: 543c4762a1bSJed Brown suffix: bddc_elast 544c4762a1bSJed Brown nsize: 4 545c4762a1bSJed Brown requires: !single 5460e75e42fSMark args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -petscpartitioner_type simple -ex56_dm_mat_type is -matis_localmat_type sbaij -pc_type bddc -pc_bddc_monolithic -attach_mat_nearnullspace 547c4762a1bSJed Brown 548c4762a1bSJed Brown test: 549c4762a1bSJed Brown suffix: fetidp_elast 550c4762a1bSJed Brown nsize: 4 551c4762a1bSJed Brown requires: !single 5520e75e42fSMark args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type fetidp -fetidp_ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -petscpartitioner_type simple -ex56_dm_mat_type is -matis_localmat_type sbaij -fetidp_bddc_pc_bddc_monolithic -attach_mat_nearnullspace 553c4762a1bSJed Brown 554908b9b43SStefano Zampini test: 555908b9b43SStefano Zampini suffix: gdsw 556908b9b43SStefano Zampini nsize: 4 557908b9b43SStefano Zampini requires: !single 558908b9b43SStefano Zampini args: -cells 2,2,1 -max_conv_its 2 -lx 1. -alpha .01 -petscspace_degree 2 -ksp_type cg -ksp_monitor_short -ksp_rtol 1.e-8 -ksp_converged_reason -petscpartitioner_type simple -ex56_dm_mat_type is -attach_mat_nearnullspace \ 559908b9b43SStefano Zampini -pc_type mg -pc_mg_galerkin -pc_mg_adapt_interp_coarse_space gdsw -pc_mg_levels 2 -mg_levels_pc_type bjacobi -mg_levels_sub_pc_type icc 560908b9b43SStefano Zampini 56135990778SJunchao Zhang testset: 56235990778SJunchao Zhang nsize: 4 56335990778SJunchao Zhang requires: !single 564bae903cbSmarkadams4 args: -cells 2,2,1 -max_conv_its 2 -petscspace_degree 2 -snes_max_it 2 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-10 -ksp_norm_type unpreconditioned -snes_rtol 1.e-10 -pc_type gamg -pc_gamg_esteig_ksp_max_it 10 -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 10 -pc_gamg_reuse_interpolation true -pc_gamg_aggressive_coarsening 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -use_mat_nearnullspace true -mg_levels_ksp_max_it 2 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.05 -mg_levels_pc_type jacobi -ksp_monitor_short -ksp_converged_reason -snes_converged_reason -snes_monitor_short -ex56_dm_view -petscpartitioner_type simple -pc_gamg_process_eq_limit 20 56535990778SJunchao Zhang output_file: output/ex56_cuda.out 56635990778SJunchao Zhang 567c4762a1bSJed Brown test: 568c4762a1bSJed Brown suffix: cuda 56935990778SJunchao Zhang requires: cuda 57035990778SJunchao Zhang args: -ex56_dm_mat_type aijcusparse -ex56_dm_vec_type cuda 5716cb74228SMark Adams 5726cb74228SMark Adams test: 573*47d993e7Ssuyashtn suffix: hip 574*47d993e7Ssuyashtn requires: hip 575*47d993e7Ssuyashtn args: -ex56_dm_mat_type aijhipsparse -ex56_dm_vec_type hip 576*47d993e7Ssuyashtn 577*47d993e7Ssuyashtn test: 5786cb74228SMark Adams suffix: viennacl 57935990778SJunchao Zhang requires: viennacl 58035990778SJunchao Zhang args: -ex56_dm_mat_type aijviennacl -ex56_dm_vec_type viennacl 5816cb74228SMark Adams 58235990778SJunchao Zhang test: 58335990778SJunchao Zhang suffix: kokkos 584dcfd994dSJunchao Zhang requires: kokkos_kernels 58535990778SJunchao Zhang args: -ex56_dm_mat_type aijkokkos -ex56_dm_vec_type kokkos 586dea3b165SRichard Tran Mills # Don't run AIJMKL caes with complex scalars because of convergence issues. 587dea3b165SRichard Tran Mills # Note that we need to test both single and multiple MPI rank cases, because these use different sparse MKL routines to implement the PtAP operation. 588a8736bf8SRichard Tran Mills test: 589a8736bf8SRichard Tran Mills suffix: seqaijmkl 590a8736bf8SRichard Tran Mills nsize: 1 591dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) !single !complex 592bae903cbSmarkadams4 args: -cells 2,2,1 -max_conv_its 2 -petscspace_degree 2 -snes_max_it 2 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned -snes_rtol 1.e-10 -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 1000 -pc_gamg_reuse_interpolation true -pc_gamg_aggressive_coarsening 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -ksp_converged_reason -snes_monitor_short -ksp_monitor_short -snes_converged_reason -use_mat_nearnullspace true -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.1 -mg_levels_pc_type jacobi -petscpartitioner_type simple -mat_block_size 3 -ex56_dm_view -run_type 1 -mat_seqaij_type seqaijmkl 593a8736bf8SRichard Tran Mills timeoutfactor: 2 594a8736bf8SRichard Tran Mills 595dea3b165SRichard Tran Mills test: 596dea3b165SRichard Tran Mills suffix: mpiaijmkl 597dea3b165SRichard Tran Mills nsize: 2 598dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) !single !complex 599bae903cbSmarkadams4 args: -cells 2,2,1 -max_conv_its 2 -petscspace_degree 2 -snes_max_it 2 -ksp_max_it 100 -ksp_type cg -ksp_rtol 1.e-11 -ksp_norm_type unpreconditioned -snes_rtol 1.e-10 -pc_type gamg -pc_gamg_type agg -pc_gamg_agg_nsmooths 1 -pc_gamg_coarse_eq_limit 1000 -pc_gamg_reuse_interpolation true -pc_gamg_aggressive_coarsening 1 -pc_gamg_threshold 0.05 -pc_gamg_threshold_scale .0 -ksp_converged_reason -snes_monitor_short -ksp_monitor_short -snes_converged_reason -use_mat_nearnullspace true -mg_levels_ksp_max_it 1 -mg_levels_ksp_type chebyshev -pc_gamg_esteig_ksp_type cg -pc_gamg_esteig_ksp_max_it 10 -mg_levels_ksp_chebyshev_esteig 0,0.05,0,1.1 -mg_levels_pc_type jacobi -petscpartitioner_type simple -mat_block_size 3 -ex56_dm_view -run_type 1 -mat_seqaij_type seqaijmkl 600dea3b165SRichard Tran Mills timeoutfactor: 2 601dea3b165SRichard Tran Mills 602c4762a1bSJed Brown TEST*/ 603