147d993e7Ssuyashtn /* Portions of this code are under: 247d993e7Ssuyashtn Copyright (C) 2022 Advanced Micro Devices, Inc. All rights reserved. 347d993e7Ssuyashtn */ 4*3667b0a5SMark Adams static char help[] = "3D tensor hexahedra & 3D Laplacian 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 13*3667b0a5SMark Adams static PetscReal s_soft_alpha = 0.01; 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 66*3667b0a5SMark Adams static void f1_u_lap(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[]) 67*3667b0a5SMark Adams { 68*3667b0a5SMark Adams PetscInt d; 69*3667b0a5SMark Adams for (d = 0; d < dim; ++d) f1[d] = u_x[d]; 70*3667b0a5SMark Adams } 71*3667b0a5SMark Adams 72c4762a1bSJed Brown /* 3D elasticity */ 73c4762a1bSJed Brown #define IDX(ii, jj, kk, ll) (27 * ii + 9 * jj + 3 * kk + ll) 74c4762a1bSJed Brown 75d71ae5a4SJacob Faibussowitsch void g3_uu_3d_private(PetscScalar g3[], const PetscReal mu, const PetscReal lambda) 76d71ae5a4SJacob Faibussowitsch { 77c4762a1bSJed Brown if (1) { 78c4762a1bSJed Brown g3[0] += lambda; 79c4762a1bSJed Brown g3[0] += mu; 80c4762a1bSJed Brown g3[0] += mu; 81c4762a1bSJed Brown g3[4] += lambda; 82c4762a1bSJed Brown g3[8] += lambda; 83c4762a1bSJed Brown g3[10] += mu; 84c4762a1bSJed Brown g3[12] += mu; 85c4762a1bSJed Brown g3[20] += mu; 86c4762a1bSJed Brown g3[24] += mu; 87c4762a1bSJed Brown g3[28] += mu; 88c4762a1bSJed Brown g3[30] += mu; 89c4762a1bSJed Brown g3[36] += lambda; 90c4762a1bSJed Brown g3[40] += lambda; 91c4762a1bSJed Brown g3[40] += mu; 92c4762a1bSJed Brown g3[40] += mu; 93c4762a1bSJed Brown g3[44] += lambda; 94c4762a1bSJed Brown g3[50] += mu; 95c4762a1bSJed Brown g3[52] += mu; 96c4762a1bSJed Brown g3[56] += mu; 97c4762a1bSJed Brown g3[60] += mu; 98c4762a1bSJed Brown g3[68] += mu; 99c4762a1bSJed Brown g3[70] += mu; 100c4762a1bSJed Brown g3[72] += lambda; 101c4762a1bSJed Brown g3[76] += lambda; 102c4762a1bSJed Brown g3[80] += lambda; 103c4762a1bSJed Brown g3[80] += mu; 104c4762a1bSJed Brown g3[80] += mu; 105c4762a1bSJed Brown } else { 106c4762a1bSJed Brown int i, j, k, l; 107c4762a1bSJed Brown static int cc = -1; 108c4762a1bSJed Brown cc++; 109c4762a1bSJed Brown for (i = 0; i < 3; ++i) { 110c4762a1bSJed Brown for (j = 0; j < 3; ++j) { 111c4762a1bSJed Brown for (k = 0; k < 3; ++k) { 112c4762a1bSJed Brown for (l = 0; l < 3; ++l) { 113c4762a1bSJed Brown if (k == l && i == j) g3[IDX(i, j, k, l)] += lambda; 114c4762a1bSJed Brown if (i == k && j == l) g3[IDX(i, j, k, l)] += mu; 115c4762a1bSJed Brown if (i == l && j == k) g3[IDX(i, j, k, l)] += mu; 116c4762a1bSJed Brown if (k == l && i == j && !cc) (void)PetscPrintf(PETSC_COMM_WORLD, "g3[%d] += lambda;\n", IDX(i, j, k, l)); 117c4762a1bSJed Brown if (i == k && j == l && !cc) (void)PetscPrintf(PETSC_COMM_WORLD, "g3[%d] += mu;\n", IDX(i, j, k, l)); 118c4762a1bSJed Brown if (i == l && j == k && !cc) (void)PetscPrintf(PETSC_COMM_WORLD, "g3[%d] += mu;\n", IDX(i, j, k, l)); 119c4762a1bSJed Brown } 120c4762a1bSJed Brown } 121c4762a1bSJed Brown } 122c4762a1bSJed Brown } 123c4762a1bSJed Brown } 124c4762a1bSJed Brown } 125c4762a1bSJed Brown 126d71ae5a4SJacob 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[]) 127d71ae5a4SJacob Faibussowitsch { 128c4762a1bSJed Brown PetscReal mu = s_mu, lambda = s_lambda, rad; 129c4762a1bSJed Brown PetscInt i; 130c4762a1bSJed Brown for (i = 0, rad = 0.; i < dim; i++) { 131c4762a1bSJed Brown PetscReal t = x[i]; 132c4762a1bSJed Brown rad += t * t; 133c4762a1bSJed Brown } 134c4762a1bSJed Brown rad = PetscSqrtReal(rad); 135c4762a1bSJed Brown if (rad > 0.25) { 136c4762a1bSJed Brown mu *= s_soft_alpha; 137c4762a1bSJed Brown lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */ 138c4762a1bSJed Brown } 139c4762a1bSJed Brown g3_uu_3d_private(g3, mu, lambda); 140c4762a1bSJed Brown } 141c4762a1bSJed Brown 142d71ae5a4SJacob 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[]) 143d71ae5a4SJacob Faibussowitsch { 144c4762a1bSJed Brown g3_uu_3d_private(g3, s_mu, s_lambda); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147*3667b0a5SMark Adams static void g3_lap(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[]) 148*3667b0a5SMark Adams { 149*3667b0a5SMark Adams PetscInt d; 150*3667b0a5SMark Adams for (d = 0; d < dim; ++d) g3[d * dim + d] = 1.0; 151*3667b0a5SMark Adams } 152*3667b0a5SMark Adams 153*3667b0a5SMark Adams static void g3_lap_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[]) 154*3667b0a5SMark Adams { 155*3667b0a5SMark Adams PetscReal lambda = 1, rad; 156*3667b0a5SMark Adams PetscInt i; 157*3667b0a5SMark Adams for (i = 0, rad = 0.; i < dim; i++) { 158*3667b0a5SMark Adams PetscReal t = x[i]; 159*3667b0a5SMark Adams rad += t * t; 160*3667b0a5SMark Adams } 161*3667b0a5SMark Adams rad = PetscSqrtReal(rad); 162*3667b0a5SMark Adams if (rad > 0.25) { lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */ } 163*3667b0a5SMark Adams for (int d = 0; d < dim; ++d) g3[d * dim + d] = lambda; 164*3667b0a5SMark Adams } 165*3667b0a5SMark Adams 166d71ae5a4SJacob 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[]) 167d71ae5a4SJacob Faibussowitsch { 168c4762a1bSJed Brown const PetscInt Ncomp = dim; 169c4762a1bSJed Brown PetscInt comp; 170c4762a1bSJed Brown 171c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) f0[comp] = 0.0; 172c4762a1bSJed Brown } 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* PI_i (x_i^4 - x_i^2) */ 175d71ae5a4SJacob 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[]) 176d71ae5a4SJacob Faibussowitsch { 177*3667b0a5SMark Adams for (int comp = 0; comp < Nf; ++comp) { 178c4762a1bSJed Brown f0[comp] = 1e5; 179*3667b0a5SMark Adams for (int i = 0; i < dim; ++i) { f0[comp] *= /* (comp+1)* */ (x[i] * x[i] * x[i] * x[i] - x[i] * x[i]); /* assumes (0,1]^D domain */ } 180c4762a1bSJed Brown } 181c4762a1bSJed Brown } 182c4762a1bSJed Brown 183d71ae5a4SJacob Faibussowitsch PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) 184d71ae5a4SJacob Faibussowitsch { 185c4762a1bSJed Brown const PetscInt Ncomp = dim; 186c4762a1bSJed Brown PetscInt comp; 187c4762a1bSJed Brown 188c4762a1bSJed Brown for (comp = 0; comp < Ncomp; ++comp) u[comp] = 0; 189c4762a1bSJed Brown return 0; 190c4762a1bSJed Brown } 191c4762a1bSJed Brown 192d71ae5a4SJacob Faibussowitsch int main(int argc, char **args) 193d71ae5a4SJacob Faibussowitsch { 194c4762a1bSJed Brown Mat Amat; 195c4762a1bSJed Brown SNES snes; 196c4762a1bSJed Brown KSP ksp; 197c4762a1bSJed Brown MPI_Comm comm; 198c4762a1bSJed Brown PetscMPIInt rank; 199956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 200c4762a1bSJed Brown PetscLogStage stage[17]; 201956f8c0dSBarry Smith #endif 202c4762a1bSJed Brown PetscBool test_nonzero_cols = PETSC_FALSE, use_nearnullspace = PETSC_TRUE, attach_nearnullspace = PETSC_FALSE; 203c4762a1bSJed Brown Vec xx, bb; 204*3667b0a5SMark Adams PetscInt iter, i, N, dim = 3, max_conv_its, sizes[7], run_type = 1, Ncomp = dim; 205*3667b0a5SMark Adams DM dm; 206c4762a1bSJed Brown PetscBool flg; 207c4762a1bSJed Brown PetscReal Lx, mdisp[10], err[10]; 208c4762a1bSJed Brown PetscFunctionBeginUser; 209327415f7SBarry Smith PetscFunctionBeginUser; 2109566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 211c4762a1bSJed Brown comm = PETSC_COMM_WORLD; 2129566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 213c4762a1bSJed Brown /* options */ 214d0609cedSBarry Smith PetscOptionsBegin(comm, NULL, "3D bilinear Q1 elasticity options", ""); 215c4762a1bSJed Brown { 216c4762a1bSJed Brown Lx = 1.; /* or ne for rod */ 217c4762a1bSJed Brown max_conv_its = 3; 2189566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-max_conv_its", "Number of iterations in convergence study", "", max_conv_its, &max_conv_its, NULL)); 219*3667b0a5SMark Adams PetscCheck(max_conv_its > 0 && max_conv_its < 8, PETSC_COMM_WORLD, PETSC_ERR_USER, "Bad number of iterations for convergence test (%" PetscInt_FMT ")", max_conv_its); 2209566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-lx", "Length of domain", "", Lx, &Lx, NULL)); 2219566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-alpha", "material coefficient inside circle", "", s_soft_alpha, &s_soft_alpha, NULL)); 2229566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-test_nonzero_cols", "nonzero test", "", test_nonzero_cols, &test_nonzero_cols, NULL)); 2239566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-use_mat_nearnullspace", "MatNearNullSpace API test", "", use_nearnullspace, &use_nearnullspace, NULL)); 2249566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-attach_mat_nearnullspace", "MatNearNullSpace API test (via MatSetNearNullSpace)", "", attach_nearnullspace, &attach_nearnullspace, NULL)); 225*3667b0a5SMark Adams PetscCall(PetscOptionsInt("-run_type", "0: twisting load on cantalever, 1: Elasticty convergence test on cube, 2: Laplacian, 3: soft core Laplacian", "", run_type, &run_type, NULL)); 226c4762a1bSJed Brown } 227d0609cedSBarry Smith PetscOptionsEnd(); 2289566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Mesh Setup", &stage[16])); 229c4762a1bSJed Brown for (iter = 0; iter < max_conv_its; iter++) { 230c4762a1bSJed Brown char str[] = "Solve 0"; 231c4762a1bSJed Brown str[6] += iter; 2329566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister(str, &stage[iter])); 233c4762a1bSJed Brown } 234c4762a1bSJed Brown /* create DM, Plex calls DMSetup */ 2359566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[16])); 236*3667b0a5SMark Adams PetscCall(DMCreate(comm, &dm)); 237*3667b0a5SMark Adams PetscCall(DMSetType(dm, DMPLEX)); 238*3667b0a5SMark Adams PetscCall(PetscObjectSetName((PetscObject)dm, "Mesh")); 239*3667b0a5SMark Adams PetscCall(DMSetFromOptions(dm)); 240*3667b0a5SMark Adams PetscCall(DMPlexDistributeSetDefault(dm, PETSC_FALSE)); 241*3667b0a5SMark Adams PetscCall(DMGetDimension(dm, &dim)); 242c4762a1bSJed Brown { 243c4762a1bSJed Brown DMLabel label; 244c4762a1bSJed Brown IS is; 2459566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "boundary")); 2469566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "boundary", &label)); 2479566063dSJacob Faibussowitsch PetscCall(DMPlexMarkBoundaryFaces(dm, 1, label)); 2482a5e0493SMark Adams if (run_type == 0) { 2499566063dSJacob Faibussowitsch PetscCall(DMGetStratumIS(dm, "boundary", 1, &is)); 2509566063dSJacob Faibussowitsch PetscCall(DMCreateLabel(dm, "Faces")); 251c4762a1bSJed Brown if (is) { 252c4762a1bSJed Brown PetscInt d, f, Nf; 253c4762a1bSJed Brown const PetscInt *faces; 254c4762a1bSJed Brown PetscInt csize; 255c4762a1bSJed Brown PetscSection cs; 256c4762a1bSJed Brown Vec coordinates; 257c4762a1bSJed Brown DM cdm; 2589566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(is, &Nf)); 2599566063dSJacob Faibussowitsch PetscCall(ISGetIndices(is, &faces)); 2609566063dSJacob Faibussowitsch PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 2619566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm, &cdm)); 2629566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(cdm, &cs)); 263c4762a1bSJed Brown /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */ 264c4762a1bSJed Brown for (f = 0; f < Nf; ++f) { 265c4762a1bSJed Brown PetscReal faceCoord; 266c4762a1bSJed Brown PetscInt b, v; 267c4762a1bSJed Brown PetscScalar *coords = NULL; 268c4762a1bSJed Brown PetscInt Nv; 2699566063dSJacob Faibussowitsch PetscCall(DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords)); 270c4762a1bSJed Brown Nv = csize / dim; /* Calculate mean coordinate vector */ 271c4762a1bSJed Brown for (d = 0; d < dim; ++d) { 272c4762a1bSJed Brown faceCoord = 0.0; 273c4762a1bSJed Brown for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v * dim + d]); 274c4762a1bSJed Brown faceCoord /= Nv; 275c4762a1bSJed Brown for (b = 0; b < 2; ++b) { 276c4762a1bSJed Brown if (PetscAbs(faceCoord - b) < PETSC_SMALL) { /* domain have not been set yet, still [0,1]^3 */ 2779566063dSJacob Faibussowitsch PetscCall(DMSetLabelValue(dm, "Faces", faces[f], d * 2 + b + 1)); 278c4762a1bSJed Brown } 279c4762a1bSJed Brown } 280c4762a1bSJed Brown } 2819566063dSJacob Faibussowitsch PetscCall(DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords)); 282c4762a1bSJed Brown } 2839566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(is, &faces)); 284c4762a1bSJed Brown } 2859566063dSJacob Faibussowitsch PetscCall(ISDestroy(&is)); 2869566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Faces", &label)); 2879566063dSJacob Faibussowitsch PetscCall(DMPlexLabelComplete(dm, label)); 288c4762a1bSJed Brown } 289c4762a1bSJed Brown } 2909566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 291c4762a1bSJed Brown for (iter = 0; iter < max_conv_its; iter++) { 2929566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[16])); 293c4762a1bSJed Brown /* snes */ 2949566063dSJacob Faibussowitsch PetscCall(SNESCreate(comm, &snes)); 2959566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm)); 296*3667b0a5SMark Adams PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 297c4762a1bSJed Brown /* fem */ 298c4762a1bSJed Brown { 299c4762a1bSJed Brown const PetscInt components[] = {0, 1, 2}; 300c4762a1bSJed Brown const PetscInt Nfid = 1, Npid = 1; 301*3667b0a5SMark Adams PetscInt fid[] = {1}; /* The fixed faces (x=0) */ 302c4762a1bSJed Brown const PetscInt pid[] = {2}; /* The faces with loading (x=L_x) */ 303c4762a1bSJed Brown PetscFE fe; 304c4762a1bSJed Brown PetscDS prob; 30545480ffeSMatthew G. Knepley DMLabel label; 306c4762a1bSJed Brown 307*3667b0a5SMark Adams if (run_type == 2 || run_type == 3) Ncomp = 1; 308*3667b0a5SMark Adams else Ncomp = dim; 309*3667b0a5SMark Adams PetscCall(PetscFECreateDefault(PETSC_COMM_SELF, dim, Ncomp, PETSC_FALSE, NULL, PETSC_DECIDE, &fe)); 3109566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)fe, "deformation")); 311c4762a1bSJed Brown /* FEM prob */ 3129566063dSJacob Faibussowitsch PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 3139566063dSJacob Faibussowitsch PetscCall(DMCreateDS(dm)); 3149566063dSJacob Faibussowitsch PetscCall(DMGetDS(dm, &prob)); 315c4762a1bSJed Brown /* setup problem */ 316*3667b0a5SMark Adams if (run_type == 1) { // elast 3179566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d)); 3189566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_u_x4, f1_u_3d)); 319*3667b0a5SMark Adams } else if (run_type == 0) { //twisted not maintained 32045480ffeSMatthew G. Knepley PetscWeakForm wf; 32145480ffeSMatthew G. Knepley PetscInt bd, i; 3229566063dSJacob Faibussowitsch PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d_alpha)); 3239566063dSJacob Faibussowitsch PetscCall(PetscDSSetResidual(prob, 0, f0_u, f1_u_3d_alpha)); 3249566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Faces", &label)); 3259566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "traction", label, Npid, pid, 0, Ncomp, components, NULL, NULL, NULL, &bd)); 3269566063dSJacob Faibussowitsch PetscCall(PetscDSGetBoundary(prob, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 3279566063dSJacob Faibussowitsch for (i = 0; i < Npid; ++i) PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, pid[i], 0, 0, 0, f0_bd_u_3d, 0, f1_bd_u)); 328*3667b0a5SMark Adams } else if (run_type == 2) { // Laplacian 329*3667b0a5SMark Adams PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_lap)); 330*3667b0a5SMark Adams PetscCall(PetscDSSetResidual(prob, 0, f0_u_x4, f1_u_lap)); 331*3667b0a5SMark Adams } else if (run_type == 3) { // soft core Laplacian 332*3667b0a5SMark Adams PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_lap_alpha)); 333*3667b0a5SMark Adams PetscCall(PetscDSSetResidual(prob, 0, f0_u_x4, f1_u_lap)); 334c4762a1bSJed Brown } 335c4762a1bSJed Brown /* bcs */ 336*3667b0a5SMark Adams if (run_type != 0) { 337c4762a1bSJed Brown PetscInt id = 1; 3389566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "boundary", &label)); 3399566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))zero, NULL, NULL, NULL)); 340c4762a1bSJed Brown } else { 3419566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "Faces", &label)); 3429566063dSJacob Faibussowitsch PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed", label, Nfid, fid, 0, Ncomp, components, (void (*)(void))zero, NULL, NULL, NULL)); 343c4762a1bSJed Brown } 3449566063dSJacob Faibussowitsch PetscCall(PetscFEDestroy(&fe)); 345c4762a1bSJed Brown } 346c4762a1bSJed Brown /* vecs & mat */ 3479566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm, &xx)); 3489566063dSJacob Faibussowitsch PetscCall(VecDuplicate(xx, &bb)); 3499566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)bb, "b")); 3509566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)xx, "u")); 3519566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(dm, &Amat)); 3529566063dSJacob Faibussowitsch PetscCall(MatSetOption(Amat, MAT_SYMMETRIC, PETSC_TRUE)); /* Some matrix kernels can take advantage of symmetry if we set this. */ 3539566063dSJacob Faibussowitsch PetscCall(MatSetOption(Amat, MAT_SYMMETRY_ETERNAL, PETSC_TRUE)); /* Inform PETSc that Amat is always symmetric, so info set above isn't lost. */ 354*3667b0a5SMark Adams PetscCall(MatSetBlockSize(Amat, Ncomp)); 3559566063dSJacob Faibussowitsch PetscCall(MatSetOption(Amat, MAT_SPD, PETSC_TRUE)); 356b94d7dedSBarry Smith PetscCall(MatSetOption(Amat, MAT_SPD_ETERNAL, PETSC_TRUE)); 3579566063dSJacob Faibussowitsch PetscCall(VecGetSize(bb, &N)); 358*3667b0a5SMark Adams sizes[iter] = N; 35963a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "%" PetscInt_FMT " global equations, %" PetscInt_FMT " vertices\n", N, N / dim)); 360*3667b0a5SMark Adams if ((use_nearnullspace || attach_nearnullspace) && N / dim > 1 && Ncomp > 1) { 361c4762a1bSJed Brown /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */ 362c4762a1bSJed Brown DM subdm; 363c4762a1bSJed Brown MatNullSpace nearNullSpace; 364c4762a1bSJed Brown PetscInt fields = 0; 365c4762a1bSJed Brown PetscObject deformation; 3669566063dSJacob Faibussowitsch PetscCall(DMCreateSubDM(dm, 1, &fields, NULL, &subdm)); 3679566063dSJacob Faibussowitsch PetscCall(DMPlexCreateRigidBody(subdm, 0, &nearNullSpace)); 3689566063dSJacob Faibussowitsch PetscCall(DMGetField(dm, 0, NULL, &deformation)); 3699566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose(deformation, "nearnullspace", (PetscObject)nearNullSpace)); 3709566063dSJacob Faibussowitsch PetscCall(DMDestroy(&subdm)); 3711baa6e33SBarry Smith if (attach_nearnullspace) PetscCall(MatSetNearNullSpace(Amat, nearNullSpace)); 3729566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nearNullSpace)); /* created by DM and destroyed by Mat */ 373c4762a1bSJed Brown } 3749566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm, NULL, NULL, NULL)); 3759566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, Amat, Amat, NULL, NULL)); 3769566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 3779566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 3789566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 3799566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[16])); 380c4762a1bSJed Brown /* ksp */ 3819566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 3829566063dSJacob Faibussowitsch PetscCall(KSPSetComputeSingularValues(ksp, PETSC_TRUE)); 383c4762a1bSJed Brown /* test BCs */ 3849566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(xx)); 385c4762a1bSJed Brown if (test_nonzero_cols) { 38648a46eb9SPierre Jolivet if (rank == 0) PetscCall(VecSetValue(xx, 0, 1.0, INSERT_VALUES)); 3879566063dSJacob Faibussowitsch PetscCall(VecAssemblyBegin(xx)); 3889566063dSJacob Faibussowitsch PetscCall(VecAssemblyEnd(xx)); 389c4762a1bSJed Brown } 3909566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(bb)); 3919566063dSJacob Faibussowitsch PetscCall(VecGetSize(bb, &i)); 392*3667b0a5SMark Adams sizes[iter] = i; 39363a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "%" PetscInt_FMT " equations in vector, %" PetscInt_FMT " vertices\n", i, i / dim)); 3949566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 395c4762a1bSJed Brown /* solve */ 396*3667b0a5SMark Adams PetscCall(SNESComputeJacobian(snes, xx, Amat, Amat)); 397*3667b0a5SMark Adams PetscCall(MatViewFromOptions(Amat, NULL, "-my_mat_view")); 3989566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage[iter])); 3999566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, bb, xx)); 4009566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 4019566063dSJacob Faibussowitsch PetscCall(VecNorm(xx, NORM_INFINITY, &mdisp[iter])); 402c4762a1bSJed Brown { 403c4762a1bSJed Brown PetscViewer viewer = NULL; 404c4762a1bSJed Brown PetscViewerFormat fmt; 405*3667b0a5SMark Adams PetscCall(PetscOptionsGetViewer(comm, NULL, "", "-vec_view", &viewer, &fmt, &flg)); 406c4762a1bSJed Brown if (flg) { 4079566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, fmt)); 4089566063dSJacob Faibussowitsch PetscCall(VecView(xx, viewer)); 4099566063dSJacob Faibussowitsch PetscCall(VecView(bb, viewer)); 4109566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 411c4762a1bSJed Brown } 4129566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 413c4762a1bSJed Brown } 414c4762a1bSJed Brown /* Free work space */ 4159566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 4169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&xx)); 4179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&bb)); 4189566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Amat)); 419*3667b0a5SMark Adams if (iter + 1 < max_conv_its) { 420*3667b0a5SMark Adams DM newdm; 421*3667b0a5SMark Adams PetscCall(DMViewFromOptions(dm, NULL, "-my_dm_view")); 422*3667b0a5SMark Adams PetscCall(DMRefine(dm, comm, &newdm)); 423*3667b0a5SMark Adams if (rank == -1) { 424*3667b0a5SMark Adams PetscDS prob; 425*3667b0a5SMark Adams PetscCall(DMGetDS(dm, &prob)); 426*3667b0a5SMark Adams PetscCall(PetscDSViewFromOptions(prob, NULL, "-ds_view")); 427*3667b0a5SMark Adams PetscCall(DMGetDS(newdm, &prob)); 428*3667b0a5SMark Adams PetscCall(PetscDSViewFromOptions(prob, NULL, "-ds_view")); 429c4762a1bSJed Brown } 430*3667b0a5SMark Adams PetscCall(DMDestroy(&dm)); 431*3667b0a5SMark Adams dm = newdm; 432*3667b0a5SMark Adams PetscCall(PetscObjectSetName((PetscObject)dm, "Mesh")); 433*3667b0a5SMark Adams PetscCall(DMViewFromOptions(dm, NULL, "-my_dm_view")); 434*3667b0a5SMark Adams PetscCall(DMSetFromOptions(dm)); 435*3667b0a5SMark Adams } 436*3667b0a5SMark Adams } 437*3667b0a5SMark Adams PetscCall(DMDestroy(&dm)); 438*3667b0a5SMark Adams if (run_type == 1) err[0] = 5.9731627e+01 - mdisp[0]; /* error with what I think is the exact solution */ 439*3667b0a5SMark Adams else if (run_type == 0) err[0] = 0; 440*3667b0a5SMark Adams else if (run_type == 2) err[0] = 3.527795e+01 - mdisp[0]; 441*3667b0a5SMark Adams else err[0] = 0; 442*3667b0a5SMark Adams PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%d] %d) N=%12" PetscInt_FMT ", max displ=%9.7e, error=%4.3e\n", rank, 0, sizes[0], (double)mdisp[0], (double)err[0])); 443c4762a1bSJed Brown for (iter = 1; iter < max_conv_its; iter++) { 444*3667b0a5SMark Adams if (run_type == 1) err[iter] = 5.9731627e+01 - mdisp[iter]; 445*3667b0a5SMark Adams else if (run_type == 0) err[iter] = 0; 446*3667b0a5SMark Adams else if (run_type == 2) err[iter] = 3.527795e+01 - mdisp[iter]; 447*3667b0a5SMark Adams else err[iter] = 0; 448*3667b0a5SMark Adams 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, sizes[iter], (double)mdisp[iter], (double)(mdisp[iter] - mdisp[iter - 1]), (double)err[iter], (double)(PetscLogReal(PetscAbs(err[iter - 1] / err[iter])) / PetscLogReal(2.)))); 449c4762a1bSJed Brown } 450c4762a1bSJed Brown 4519566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 452b122ec5aSJacob Faibussowitsch return 0; 453c4762a1bSJed Brown } 454c4762a1bSJed Brown 455c4762a1bSJed Brown /*TEST 456c4762a1bSJed Brown 457*3667b0a5SMark Adams testset: 458c4762a1bSJed Brown nsize: 4 459c4762a1bSJed Brown requires: !single 460*3667b0a5SMark Adams args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -max_conv_its 2 -petscspace_degree 2 -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 1 -pc_gamg_threshold 0.001 -ksp_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 -my_dm_view -snes_lag_jacobian -2 -snes_type ksponly 461c4762a1bSJed Brown timeoutfactor: 2 462*3667b0a5SMark Adams test: 463*3667b0a5SMark Adams suffix: 0 464*3667b0a5SMark Adams args: -run_type 1 465*3667b0a5SMark Adams test: 466*3667b0a5SMark Adams suffix: 1 467*3667b0a5SMark Adams args: -run_type 2 468c4762a1bSJed Brown 469c4762a1bSJed Brown # HYPRE PtAP broken with complex numbers 470c4762a1bSJed Brown test: 471c4762a1bSJed Brown suffix: hypre 472263f2b91SStefano Zampini requires: hypre !single !complex !defined(PETSC_HAVE_HYPRE_DEVICE) 473c4762a1bSJed Brown nsize: 4 474*3667b0a5SMark Adams args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 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 475c4762a1bSJed Brown 476c4762a1bSJed Brown test: 477c4762a1bSJed Brown suffix: ml 478c4762a1bSJed Brown requires: ml !single 479c4762a1bSJed Brown nsize: 4 480*3667b0a5SMark Adams args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 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 481c4762a1bSJed Brown 482c4762a1bSJed Brown test: 483c4762a1bSJed Brown suffix: hpddm 484dfd57a17SPierre Jolivet requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 485c4762a1bSJed Brown nsize: 4 486*3667b0a5SMark Adams args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 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 487c4762a1bSJed Brown 488c4762a1bSJed Brown test: 48963b77682SMark Adams suffix: repart 490c4762a1bSJed Brown nsize: 4 491c4762a1bSJed Brown requires: parmetis !single 492*3667b0a5SMark Adams args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 2,2,1 -max_conv_its 2 -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 -pc_gamg_process_eq_limit 20 -pc_gamg_coarse_eq_limit 10 -ksp_converged_reason -pc_gamg_reuse_interpolation true -petscpartitioner_type simple 493c4762a1bSJed Brown 494c4762a1bSJed Brown test: 495c4762a1bSJed Brown suffix: bddc 496c4762a1bSJed Brown nsize: 4 497c4762a1bSJed Brown requires: !single 498*3667b0a5SMark Adams args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 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 -dm_mat_type is -matis_localmat_type {{sbaij baij aij}} -pc_type bddc 499c4762a1bSJed Brown 500c4762a1bSJed Brown testset: 501c4762a1bSJed Brown nsize: 4 502c4762a1bSJed Brown requires: !single 503*3667b0a5SMark Adams args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 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 -dm_mat_type is -matis_localmat_type aij -pc_type bddc -attach_mat_nearnullspace {{0 1}separate output} 504c4762a1bSJed Brown test: 505c4762a1bSJed Brown suffix: bddc_approx_gamg 506bae903cbSmarkadams4 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 507c4762a1bSJed Brown # HYPRE PtAP broken with complex numbers 508c4762a1bSJed Brown test: 509263f2b91SStefano Zampini requires: hypre !complex !defined(PETSC_HAVE_HYPRE_DEVICE) 510c4762a1bSJed Brown suffix: bddc_approx_hypre 511c4762a1bSJed 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 512c4762a1bSJed Brown test: 513c4762a1bSJed Brown requires: ml 514c4762a1bSJed Brown suffix: bddc_approx_ml 5150e75e42fSMark 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 516c4762a1bSJed Brown 517c4762a1bSJed Brown test: 518c4762a1bSJed Brown suffix: fetidp 519c4762a1bSJed Brown nsize: 4 520c4762a1bSJed Brown requires: !single 521*3667b0a5SMark Adams args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 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 -dm_mat_type is -matis_localmat_type {{sbaij baij aij}} 522c4762a1bSJed Brown 523c4762a1bSJed Brown test: 524c4762a1bSJed Brown suffix: bddc_elast 525c4762a1bSJed Brown nsize: 4 526c4762a1bSJed Brown requires: !single 527*3667b0a5SMark Adams args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 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 -dm_mat_type is -matis_localmat_type sbaij -pc_type bddc -pc_bddc_monolithic -attach_mat_nearnullspace 528c4762a1bSJed Brown 529c4762a1bSJed Brown test: 530c4762a1bSJed Brown suffix: fetidp_elast 531c4762a1bSJed Brown nsize: 4 532c4762a1bSJed Brown requires: !single 533*3667b0a5SMark Adams args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 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 -dm_mat_type is -matis_localmat_type sbaij -fetidp_bddc_pc_bddc_monolithic -attach_mat_nearnullspace 534c4762a1bSJed Brown 535908b9b43SStefano Zampini test: 536908b9b43SStefano Zampini suffix: gdsw 537908b9b43SStefano Zampini nsize: 4 538908b9b43SStefano Zampini requires: !single 539*3667b0a5SMark Adams args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 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 -dm_mat_type is -attach_mat_nearnullspace \ 540908b9b43SStefano 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 541908b9b43SStefano Zampini 54235990778SJunchao Zhang testset: 54335990778SJunchao Zhang nsize: 4 54435990778SJunchao Zhang requires: !single 545*3667b0a5SMark Adams args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 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_monitor_short -dm_view -petscpartitioner_type simple -pc_gamg_process_eq_limit 20 54635990778SJunchao Zhang output_file: output/ex56_cuda.out 54735990778SJunchao Zhang 548c4762a1bSJed Brown test: 549c4762a1bSJed Brown suffix: cuda 55035990778SJunchao Zhang requires: cuda 551*3667b0a5SMark Adams args: -dm_mat_type aijcusparse -dm_vec_type cuda 5526cb74228SMark Adams 5536cb74228SMark Adams test: 55447d993e7Ssuyashtn suffix: hip 55547d993e7Ssuyashtn requires: hip 556*3667b0a5SMark Adams args: -dm_mat_type aijhipsparse -dm_vec_type hip 55747d993e7Ssuyashtn 55847d993e7Ssuyashtn test: 5596cb74228SMark Adams suffix: viennacl 56035990778SJunchao Zhang requires: viennacl 561*3667b0a5SMark Adams args: -dm_mat_type aijviennacl -dm_vec_type viennacl 5626cb74228SMark Adams 56335990778SJunchao Zhang test: 56435990778SJunchao Zhang suffix: kokkos 565dcfd994dSJunchao Zhang requires: kokkos_kernels 566*3667b0a5SMark Adams args: -dm_mat_type aijkokkos -dm_vec_type kokkos 567dea3b165SRichard Tran Mills # Don't run AIJMKL caes with complex scalars because of convergence issues. 568dea3b165SRichard 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. 569a8736bf8SRichard Tran Mills test: 570a8736bf8SRichard Tran Mills suffix: seqaijmkl 571a8736bf8SRichard Tran Mills nsize: 1 572dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) !single !complex 573*3667b0a5SMark Adams args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 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 -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 -dm_view -mat_seqaij_type seqaijmkl 574a8736bf8SRichard Tran Mills timeoutfactor: 2 575a8736bf8SRichard Tran Mills 576dea3b165SRichard Tran Mills test: 577dea3b165SRichard Tran Mills suffix: mpiaijmkl 578*3667b0a5SMark Adams nsize: 4 579dfd57a17SPierre Jolivet requires: defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) !single !complex 580*3667b0a5SMark Adams args: -dm_plex_dim 3 -dm_plex_simplex 0 -dm_plex_box_lower 0,0,0 -dm_plex_box_upper 1,1,1 -run_type 1 -dm_plex_box_faces 2,2,1 -petscpartitioner_simple_process_grid 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 -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 -dm_view -mat_seqaij_type seqaijmkl 581dea3b165SRichard Tran Mills timeoutfactor: 2 582dea3b165SRichard Tran Mills 583c4762a1bSJed Brown TEST*/ 584