1 static char help[] = "3D, tensor hexahedra (Q1-K), displacement finite element formulation\n\ 2 of linear elasticity. E=1.0, nu=1/3.\n\ 3 Unit cube domain with Dirichlet boundary\n\n"; 4 5 #include <petscdmplex.h> 6 #include <petscsnes.h> 7 #include <petscds.h> 8 #include <petscdmforest.h> 9 10 static PetscReal s_soft_alpha = 1.e-3; 11 static PetscReal s_mu = 0.4; 12 static PetscReal s_lambda = 0.4; 13 14 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[]) { 15 f0[0] = 1; /* x direction pull */ 16 f0[1] = -x[2]; /* add a twist around x-axis */ 17 f0[2] = x[1]; 18 } 19 20 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[]) { 21 const PetscInt Ncomp = dim; 22 PetscInt comp, d; 23 for (comp = 0; comp < Ncomp; ++comp) { 24 for (d = 0; d < dim; ++d) { f1[comp * dim + d] = 0.0; } 25 } 26 } 27 28 /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */ 29 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[]) { 30 PetscReal trace, mu = s_mu, lambda = s_lambda, rad; 31 PetscInt i, j; 32 for (i = 0, rad = 0.; i < dim; i++) { 33 PetscReal t = x[i]; 34 rad += t * t; 35 } 36 rad = PetscSqrtReal(rad); 37 if (rad > 0.25) { 38 mu *= s_soft_alpha; 39 lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */ 40 } 41 for (i = 0, trace = 0; i < dim; ++i) { trace += PetscRealPart(u_x[i * dim + i]); } 42 for (i = 0; i < dim; ++i) { 43 for (j = 0; j < dim; ++j) { f1[i * dim + j] = mu * (u_x[i * dim + j] + u_x[j * dim + i]); } 44 f1[i * dim + i] += lambda * trace; 45 } 46 } 47 48 /* gradU[comp*dim+d] = {u_x, u_y} or {u_x, u_y, u_z} */ 49 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[]) { 50 PetscReal trace, mu = s_mu, lambda = s_lambda; 51 PetscInt i, j; 52 for (i = 0, trace = 0; i < dim; ++i) { trace += PetscRealPart(u_x[i * dim + i]); } 53 for (i = 0; i < dim; ++i) { 54 for (j = 0; j < dim; ++j) { f1[i * dim + j] = mu * (u_x[i * dim + j] + u_x[j * dim + i]); } 55 f1[i * dim + i] += lambda * trace; 56 } 57 } 58 59 /* 3D elasticity */ 60 #define IDX(ii, jj, kk, ll) (27 * ii + 9 * jj + 3 * kk + ll) 61 62 void g3_uu_3d_private(PetscScalar g3[], const PetscReal mu, const PetscReal lambda) { 63 if (1) { 64 g3[0] += lambda; 65 g3[0] += mu; 66 g3[0] += mu; 67 g3[4] += lambda; 68 g3[8] += lambda; 69 g3[10] += mu; 70 g3[12] += mu; 71 g3[20] += mu; 72 g3[24] += mu; 73 g3[28] += mu; 74 g3[30] += mu; 75 g3[36] += lambda; 76 g3[40] += lambda; 77 g3[40] += mu; 78 g3[40] += mu; 79 g3[44] += lambda; 80 g3[50] += mu; 81 g3[52] += mu; 82 g3[56] += mu; 83 g3[60] += mu; 84 g3[68] += mu; 85 g3[70] += mu; 86 g3[72] += lambda; 87 g3[76] += lambda; 88 g3[80] += lambda; 89 g3[80] += mu; 90 g3[80] += mu; 91 } else { 92 int i, j, k, l; 93 static int cc = -1; 94 cc++; 95 for (i = 0; i < 3; ++i) { 96 for (j = 0; j < 3; ++j) { 97 for (k = 0; k < 3; ++k) { 98 for (l = 0; l < 3; ++l) { 99 if (k == l && i == j) g3[IDX(i, j, k, l)] += lambda; 100 if (i == k && j == l) g3[IDX(i, j, k, l)] += mu; 101 if (i == l && j == k) g3[IDX(i, j, k, l)] += mu; 102 if (k == l && i == j && !cc) (void)PetscPrintf(PETSC_COMM_WORLD, "g3[%d] += lambda;\n", IDX(i, j, k, l)); 103 if (i == k && j == l && !cc) (void)PetscPrintf(PETSC_COMM_WORLD, "g3[%d] += mu;\n", IDX(i, j, k, l)); 104 if (i == l && j == k && !cc) (void)PetscPrintf(PETSC_COMM_WORLD, "g3[%d] += mu;\n", IDX(i, j, k, l)); 105 } 106 } 107 } 108 } 109 } 110 } 111 112 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[]) { 113 PetscReal mu = s_mu, lambda = s_lambda, rad; 114 PetscInt i; 115 for (i = 0, rad = 0.; i < dim; i++) { 116 PetscReal t = x[i]; 117 rad += t * t; 118 } 119 rad = PetscSqrtReal(rad); 120 if (rad > 0.25) { 121 mu *= s_soft_alpha; 122 lambda *= s_soft_alpha; /* we could keep the bulk the same like rubberish */ 123 } 124 g3_uu_3d_private(g3, mu, lambda); 125 } 126 127 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[]) { 128 g3_uu_3d_private(g3, s_mu, s_lambda); 129 } 130 131 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[]) { 132 const PetscInt Ncomp = dim; 133 PetscInt comp; 134 135 for (comp = 0; comp < Ncomp; ++comp) f0[comp] = 0.0; 136 } 137 138 /* PI_i (x_i^4 - x_i^2) */ 139 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[]) { 140 const PetscInt Ncomp = dim; 141 PetscInt comp, i; 142 143 for (comp = 0; comp < Ncomp; ++comp) { 144 f0[comp] = 1e5; 145 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 */ } 146 } 147 } 148 149 PetscErrorCode zero(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nf, PetscScalar *u, void *ctx) { 150 const PetscInt Ncomp = dim; 151 PetscInt comp; 152 153 for (comp = 0; comp < Ncomp; ++comp) u[comp] = 0; 154 return 0; 155 } 156 157 int main(int argc, char **args) { 158 Mat Amat; 159 SNES snes; 160 KSP ksp; 161 MPI_Comm comm; 162 PetscMPIInt rank; 163 #if defined(PETSC_USE_LOG) 164 PetscLogStage stage[17]; 165 #endif 166 PetscBool test_nonzero_cols = PETSC_FALSE, use_nearnullspace = PETSC_TRUE, attach_nearnullspace = PETSC_FALSE; 167 Vec xx, bb; 168 PetscInt iter, i, N, dim = 3, cells[3] = {1, 1, 1}, max_conv_its, local_sizes[7], run_type = 1; 169 DM dm, distdm, basedm; 170 PetscBool flg; 171 char convType[256]; 172 PetscReal Lx, mdisp[10], err[10]; 173 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"}; 174 PetscFunctionBeginUser; 175 PetscFunctionBeginUser; 176 PetscCall(PetscInitialize(&argc, &args, (char *)0, help)); 177 comm = PETSC_COMM_WORLD; 178 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 179 /* options */ 180 PetscOptionsBegin(comm, NULL, "3D bilinear Q1 elasticity options", ""); 181 { 182 i = 3; 183 PetscCall(PetscOptionsIntArray("-cells", "Number of (flux tube) processor in each dimension", "ex56.c", cells, &i, NULL)); 184 185 Lx = 1.; /* or ne for rod */ 186 max_conv_its = 3; 187 PetscCall(PetscOptionsInt("-max_conv_its", "Number of iterations in convergence study", "", max_conv_its, &max_conv_its, NULL)); 188 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); 189 PetscCall(PetscOptionsReal("-lx", "Length of domain", "", Lx, &Lx, NULL)); 190 PetscCall(PetscOptionsReal("-alpha", "material coefficient inside circle", "", s_soft_alpha, &s_soft_alpha, NULL)); 191 PetscCall(PetscOptionsBool("-test_nonzero_cols", "nonzero test", "", test_nonzero_cols, &test_nonzero_cols, NULL)); 192 PetscCall(PetscOptionsBool("-use_mat_nearnullspace", "MatNearNullSpace API test", "", use_nearnullspace, &use_nearnullspace, NULL)); 193 PetscCall(PetscOptionsBool("-attach_mat_nearnullspace", "MatNearNullSpace API test (via MatSetNearNullSpace)", "", attach_nearnullspace, &attach_nearnullspace, NULL)); 194 PetscCall(PetscOptionsInt("-run_type", "0: twisting load on cantalever, 1: 3rd order accurate convergence test", "", run_type, &run_type, NULL)); 195 } 196 PetscOptionsEnd(); 197 PetscCall(PetscLogStageRegister("Mesh Setup", &stage[16])); 198 for (iter = 0; iter < max_conv_its; iter++) { 199 char str[] = "Solve 0"; 200 str[6] += iter; 201 PetscCall(PetscLogStageRegister(str, &stage[iter])); 202 } 203 /* create DM, Plex calls DMSetup */ 204 PetscCall(PetscLogStagePush(stage[16])); 205 PetscCall(DMPlexCreateBoxMesh(comm, dim, PETSC_FALSE, cells, NULL, NULL, NULL, PETSC_TRUE, &dm)); 206 { 207 DMLabel label; 208 IS is; 209 PetscCall(DMCreateLabel(dm, "boundary")); 210 PetscCall(DMGetLabel(dm, "boundary", &label)); 211 PetscCall(DMPlexMarkBoundaryFaces(dm, 1, label)); 212 if (run_type == 0) { 213 PetscCall(DMGetStratumIS(dm, "boundary", 1, &is)); 214 PetscCall(DMCreateLabel(dm, "Faces")); 215 if (is) { 216 PetscInt d, f, Nf; 217 const PetscInt *faces; 218 PetscInt csize; 219 PetscSection cs; 220 Vec coordinates; 221 DM cdm; 222 PetscCall(ISGetLocalSize(is, &Nf)); 223 PetscCall(ISGetIndices(is, &faces)); 224 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 225 PetscCall(DMGetCoordinateDM(dm, &cdm)); 226 PetscCall(DMGetLocalSection(cdm, &cs)); 227 /* Check for each boundary face if any component of its centroid is either 0.0 or 1.0 */ 228 for (f = 0; f < Nf; ++f) { 229 PetscReal faceCoord; 230 PetscInt b, v; 231 PetscScalar *coords = NULL; 232 PetscInt Nv; 233 PetscCall(DMPlexVecGetClosure(cdm, cs, coordinates, faces[f], &csize, &coords)); 234 Nv = csize / dim; /* Calculate mean coordinate vector */ 235 for (d = 0; d < dim; ++d) { 236 faceCoord = 0.0; 237 for (v = 0; v < Nv; ++v) faceCoord += PetscRealPart(coords[v * dim + d]); 238 faceCoord /= Nv; 239 for (b = 0; b < 2; ++b) { 240 if (PetscAbs(faceCoord - b) < PETSC_SMALL) { /* domain have not been set yet, still [0,1]^3 */ 241 PetscCall(DMSetLabelValue(dm, "Faces", faces[f], d * 2 + b + 1)); 242 } 243 } 244 } 245 PetscCall(DMPlexVecRestoreClosure(cdm, cs, coordinates, faces[f], &csize, &coords)); 246 } 247 PetscCall(ISRestoreIndices(is, &faces)); 248 } 249 PetscCall(ISDestroy(&is)); 250 PetscCall(DMGetLabel(dm, "Faces", &label)); 251 PetscCall(DMPlexLabelComplete(dm, label)); 252 } 253 } 254 { 255 PetscInt dimEmbed, i; 256 PetscInt nCoords; 257 PetscScalar *coords, bounds[] = { 258 0, 1, -.5, .5, -.5, .5, 259 }; /* x_min,x_max,y_min,y_max */ 260 Vec coordinates; 261 bounds[1] = Lx; 262 if (run_type == 1) { 263 for (i = 0; i < 2 * dim; i++) bounds[i] = (i % 2) ? 1 : 0; 264 } 265 PetscCall(DMGetCoordinatesLocal(dm, &coordinates)); 266 PetscCall(DMGetCoordinateDim(dm, &dimEmbed)); 267 PetscCheck(dimEmbed == dim, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "dimEmbed != dim %" PetscInt_FMT, dimEmbed); 268 PetscCall(VecGetLocalSize(coordinates, &nCoords)); 269 PetscCheck((nCoords % dimEmbed) == 0, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Coordinate vector the wrong size"); 270 PetscCall(VecGetArray(coordinates, &coords)); 271 for (i = 0; i < nCoords; i += dimEmbed) { 272 PetscInt j; 273 PetscScalar *coord = &coords[i]; 274 for (j = 0; j < dimEmbed; j++) { coord[j] = bounds[2 * j] + coord[j] * (bounds[2 * j + 1] - bounds[2 * j]); } 275 } 276 PetscCall(VecRestoreArray(coordinates, &coords)); 277 PetscCall(DMSetCoordinatesLocal(dm, coordinates)); 278 } 279 280 /* convert to p4est, and distribute */ 281 PetscOptionsBegin(comm, "", "Mesh conversion options", "DMPLEX"); 282 PetscCall(PetscOptionsFList("-dm_type", "Convert DMPlex to another format (should not be Plex!)", "ex56.c", DMList, DMPLEX, convType, 256, &flg)); 283 PetscOptionsEnd(); 284 if (flg) { 285 DM newdm; 286 PetscCall(DMConvert(dm, convType, &newdm)); 287 if (newdm) { 288 const char *prefix; 289 PetscBool isForest; 290 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 291 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)newdm, prefix)); 292 PetscCall(DMIsForest(newdm, &isForest)); 293 PetscCheck(isForest, PETSC_COMM_WORLD, PETSC_ERR_USER, "Converted to non Forest?"); 294 PetscCall(DMDestroy(&dm)); 295 dm = newdm; 296 } else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_USER, "Convert failed?"); 297 } else { 298 PetscPartitioner part; 299 /* Plex Distribute mesh over processes */ 300 PetscCall(DMPlexGetPartitioner(dm, &part)); 301 PetscCall(PetscPartitionerSetFromOptions(part)); 302 PetscCall(DMPlexDistribute(dm, 0, NULL, &distdm)); 303 if (distdm) { 304 const char *prefix; 305 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 306 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)distdm, prefix)); 307 PetscCall(DMDestroy(&dm)); 308 dm = distdm; 309 } 310 } 311 PetscCall(PetscLogStagePop()); 312 basedm = dm; 313 dm = NULL; 314 315 for (iter = 0; iter < max_conv_its; iter++) { 316 PetscCall(PetscLogStagePush(stage[16])); 317 /* make new DM */ 318 PetscCall(DMClone(basedm, &dm)); 319 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)dm, "ex56_")); 320 PetscCall(PetscObjectSetName((PetscObject)dm, "Mesh")); 321 if (max_conv_its > 1) { 322 /* If max_conv_its == 1, then we are not doing a convergence study. */ 323 PetscCall(PetscOptionsInsertString(NULL, options[iter])); 324 } 325 PetscCall(DMSetFromOptions(dm)); /* refinement done here in Plex, p4est */ 326 /* snes */ 327 PetscCall(SNESCreate(comm, &snes)); 328 PetscCall(SNESSetDM(snes, dm)); 329 /* fem */ 330 { 331 const PetscInt Ncomp = dim; 332 const PetscInt components[] = {0, 1, 2}; 333 const PetscInt Nfid = 1, Npid = 1; 334 const PetscInt fid[] = {1}; /* The fixed faces (x=0) */ 335 const PetscInt pid[] = {2}; /* The faces with loading (x=L_x) */ 336 PetscFE fe; 337 PetscDS prob; 338 DMLabel label; 339 DM cdm = dm; 340 341 PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, dim, PETSC_FALSE, NULL, PETSC_DECIDE, &fe)); /* elasticity */ 342 PetscCall(PetscObjectSetName((PetscObject)fe, "deformation")); 343 /* FEM prob */ 344 PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe)); 345 PetscCall(DMCreateDS(dm)); 346 PetscCall(DMGetDS(dm, &prob)); 347 /* setup problem */ 348 if (run_type == 1) { 349 PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d)); 350 PetscCall(PetscDSSetResidual(prob, 0, f0_u_x4, f1_u_3d)); 351 } else { 352 PetscWeakForm wf; 353 PetscInt bd, i; 354 355 PetscCall(PetscDSSetJacobian(prob, 0, 0, NULL, NULL, NULL, g3_uu_3d_alpha)); 356 PetscCall(PetscDSSetResidual(prob, 0, f0_u, f1_u_3d_alpha)); 357 358 PetscCall(DMGetLabel(dm, "Faces", &label)); 359 PetscCall(DMAddBoundary(dm, DM_BC_NATURAL, "traction", label, Npid, pid, 0, Ncomp, components, NULL, NULL, NULL, &bd)); 360 PetscCall(PetscDSGetBoundary(prob, bd, &wf, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL)); 361 for (i = 0; i < Npid; ++i) PetscCall(PetscWeakFormSetIndexBdResidual(wf, label, pid[i], 0, 0, 0, f0_bd_u_3d, 0, f1_bd_u)); 362 } 363 /* bcs */ 364 if (run_type == 1) { 365 PetscInt id = 1; 366 PetscCall(DMGetLabel(dm, "boundary", &label)); 367 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "wall", label, 1, &id, 0, 0, NULL, (void (*)(void))zero, NULL, NULL, NULL)); 368 } else { 369 PetscCall(DMGetLabel(dm, "Faces", &label)); 370 PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "fixed", label, Nfid, fid, 0, Ncomp, components, (void (*)(void))zero, NULL, NULL, NULL)); 371 } 372 while (cdm) { 373 PetscCall(DMCopyDisc(dm, cdm)); 374 PetscCall(DMGetCoarseDM(cdm, &cdm)); 375 } 376 PetscCall(PetscFEDestroy(&fe)); 377 } 378 /* vecs & mat */ 379 PetscCall(DMCreateGlobalVector(dm, &xx)); 380 PetscCall(VecDuplicate(xx, &bb)); 381 PetscCall(PetscObjectSetName((PetscObject)bb, "b")); 382 PetscCall(PetscObjectSetName((PetscObject)xx, "u")); 383 PetscCall(DMCreateMatrix(dm, &Amat)); 384 PetscCall(MatSetOption(Amat, MAT_SYMMETRIC, PETSC_TRUE)); /* Some matrix kernels can take advantage of symmetry if we set this. */ 385 PetscCall(MatSetOption(Amat, MAT_SYMMETRY_ETERNAL, PETSC_TRUE)); /* Inform PETSc that Amat is always symmetric, so info set above isn't lost. */ 386 PetscCall(MatSetBlockSize(Amat, 3)); 387 PetscCall(MatSetOption(Amat, MAT_SPD, PETSC_TRUE)); 388 PetscCall(MatSetOption(Amat, MAT_SPD_ETERNAL, PETSC_TRUE)); 389 PetscCall(VecGetSize(bb, &N)); 390 local_sizes[iter] = N; 391 PetscCall(PetscInfo(snes, "%" PetscInt_FMT " global equations, %" PetscInt_FMT " vertices\n", N, N / dim)); 392 if ((use_nearnullspace || attach_nearnullspace) && N / dim > 1) { 393 /* Set up the near null space (a.k.a. rigid body modes) that will be used by the multigrid preconditioner */ 394 DM subdm; 395 MatNullSpace nearNullSpace; 396 PetscInt fields = 0; 397 PetscObject deformation; 398 PetscCall(DMCreateSubDM(dm, 1, &fields, NULL, &subdm)); 399 PetscCall(DMPlexCreateRigidBody(subdm, 0, &nearNullSpace)); 400 PetscCall(DMGetField(dm, 0, NULL, &deformation)); 401 PetscCall(PetscObjectCompose(deformation, "nearnullspace", (PetscObject)nearNullSpace)); 402 PetscCall(DMDestroy(&subdm)); 403 if (attach_nearnullspace) PetscCall(MatSetNearNullSpace(Amat, nearNullSpace)); 404 PetscCall(MatNullSpaceDestroy(&nearNullSpace)); /* created by DM and destroyed by Mat */ 405 } 406 PetscCall(DMPlexSetSNESLocalFEM(dm, NULL, NULL, NULL)); 407 PetscCall(SNESSetJacobian(snes, Amat, Amat, NULL, NULL)); 408 PetscCall(SNESSetFromOptions(snes)); 409 PetscCall(DMSetUp(dm)); 410 PetscCall(PetscLogStagePop()); 411 PetscCall(PetscLogStagePush(stage[16])); 412 /* ksp */ 413 PetscCall(SNESGetKSP(snes, &ksp)); 414 PetscCall(KSPSetComputeSingularValues(ksp, PETSC_TRUE)); 415 /* test BCs */ 416 PetscCall(VecZeroEntries(xx)); 417 if (test_nonzero_cols) { 418 if (rank == 0) { PetscCall(VecSetValue(xx, 0, 1.0, INSERT_VALUES)); } 419 PetscCall(VecAssemblyBegin(xx)); 420 PetscCall(VecAssemblyEnd(xx)); 421 } 422 PetscCall(VecZeroEntries(bb)); 423 PetscCall(VecGetSize(bb, &i)); 424 local_sizes[iter] = i; 425 PetscCall(PetscInfo(snes, "%" PetscInt_FMT " equations in vector, %" PetscInt_FMT " vertices\n", i, i / dim)); 426 PetscCall(PetscLogStagePop()); 427 /* solve */ 428 PetscCall(PetscLogStagePush(stage[iter])); 429 PetscCall(SNESSolve(snes, bb, xx)); 430 PetscCall(PetscLogStagePop()); 431 PetscCall(VecNorm(xx, NORM_INFINITY, &mdisp[iter])); 432 PetscCall(DMViewFromOptions(dm, NULL, "-dm_view")); 433 { 434 PetscViewer viewer = NULL; 435 PetscViewerFormat fmt; 436 PetscCall(PetscOptionsGetViewer(comm, NULL, "ex56_", "-vec_view", &viewer, &fmt, &flg)); 437 if (flg) { 438 PetscCall(PetscViewerPushFormat(viewer, fmt)); 439 PetscCall(VecView(xx, viewer)); 440 PetscCall(VecView(bb, viewer)); 441 PetscCall(PetscViewerPopFormat(viewer)); 442 } 443 PetscCall(PetscViewerDestroy(&viewer)); 444 } 445 /* Free work space */ 446 PetscCall(DMDestroy(&dm)); 447 PetscCall(SNESDestroy(&snes)); 448 PetscCall(VecDestroy(&xx)); 449 PetscCall(VecDestroy(&bb)); 450 PetscCall(MatDestroy(&Amat)); 451 } 452 PetscCall(DMDestroy(&basedm)); 453 if (run_type == 1) err[0] = 59.975208 - mdisp[0]; /* error with what I think is the exact solution */ 454 else err[0] = 171.038 - mdisp[0]; 455 for (iter = 1; iter < max_conv_its; iter++) { 456 if (run_type == 1) err[iter] = 59.975208 - mdisp[iter]; 457 else err[iter] = 171.038 - mdisp[iter]; 458 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.)))); 459 } 460 461 PetscCall(PetscFinalize()); 462 return 0; 463 } 464 465 /*TEST 466 467 test: 468 suffix: 0 469 nsize: 4 470 requires: !single 471 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 472 timeoutfactor: 2 473 474 # HYPRE PtAP broken with complex numbers 475 test: 476 suffix: hypre 477 requires: hypre !single !complex !defined(PETSC_HAVE_HYPRE_DEVICE) 478 nsize: 4 479 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 480 481 test: 482 suffix: ml 483 requires: ml !single 484 nsize: 4 485 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 486 487 test: 488 suffix: hpddm 489 requires: hpddm slepc !single defined(PETSC_HAVE_DYNAMIC_LIBRARIES) defined(PETSC_USE_SHARED_LIBRARIES) 490 nsize: 4 491 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 492 493 test: 494 suffix: repart 495 nsize: 4 496 requires: parmetis !single 497 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 498 499 test: 500 suffix: bddc 501 nsize: 4 502 requires: !single 503 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 504 505 testset: 506 nsize: 4 507 requires: !single 508 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} 509 test: 510 suffix: bddc_approx_gamg 511 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 512 # HYPRE PtAP broken with complex numbers 513 test: 514 requires: hypre !complex !defined(PETSC_HAVE_HYPRE_DEVICE) 515 suffix: bddc_approx_hypre 516 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 517 test: 518 requires: ml 519 suffix: bddc_approx_ml 520 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 521 522 test: 523 suffix: fetidp 524 nsize: 4 525 requires: !single 526 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}} 527 528 test: 529 suffix: bddc_elast 530 nsize: 4 531 requires: !single 532 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 533 534 test: 535 suffix: fetidp_elast 536 nsize: 4 537 requires: !single 538 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 539 540 test: 541 suffix: gdsw 542 nsize: 4 543 requires: !single 544 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 \ 545 -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 546 547 testset: 548 nsize: 4 549 requires: !single 550 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 551 output_file: output/ex56_cuda.out 552 553 test: 554 suffix: cuda 555 requires: cuda 556 args: -ex56_dm_mat_type aijcusparse -ex56_dm_vec_type cuda 557 558 test: 559 suffix: viennacl 560 requires: viennacl 561 args: -ex56_dm_mat_type aijviennacl -ex56_dm_vec_type viennacl 562 563 test: 564 suffix: kokkos 565 requires: !sycl kokkos_kernels 566 args: -ex56_dm_mat_type aijkokkos -ex56_dm_vec_type kokkos 567 # Don't run AIJMKL caes with complex scalars because of convergence issues. 568 # 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. 569 test: 570 suffix: seqaijmkl 571 nsize: 1 572 requires: defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) !single !complex 573 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 574 timeoutfactor: 2 575 576 test: 577 suffix: mpiaijmkl 578 nsize: 2 579 requires: defined(PETSC_HAVE_MKL_SPARSE_OPTIMIZE) !single !complex 580 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 581 timeoutfactor: 2 582 583 TEST*/ 584