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