1 static const char help[] = "1D multiphysics prototype with analytic Jacobians to solve individual problems and a coupled problem.\n\n"; 2 3 /* Solve a PDE coupled to an algebraic system in 1D 4 * 5 * PDE (U): 6 * -(k u_x)_x = 1 on (0,1), subject to u(0) = 0, u(1) = 1 7 * Algebraic (K): 8 * exp(k-1) + k = 1/(1/(1+u) + 1/(1+u_x^2)) 9 * 10 * The discretization places k at staggered points, and a separate DMDA is used for each "physics". 11 * 12 * This example is a prototype for coupling in multi-physics problems, therefore residual evaluation and assembly for 13 * each problem (referred to as U and K) are written separately. This permits the same "physics" code to be used for 14 * solving each uncoupled problem as well as the coupled system. In particular, run with -problem_type 0 to solve only 15 * problem U (with K fixed), -problem_type 1 to solve only K (with U fixed), and -problem_type 2 to solve both at once. 16 * 17 * In all cases, a fully-assembled analytic Jacobian is available, so the systems can be solved with a direct solve or 18 * any other standard method. Additionally, by running with 19 * 20 * -pack_dm_mat_type nest 21 * 22 * The same code assembles a coupled matrix where each block is stored separately, which allows the use of PCFieldSplit 23 * without copying values to extract submatrices. 24 */ 25 26 #include <petscsnes.h> 27 #include <petscdm.h> 28 #include <petscdmda.h> 29 #include <petscdmcomposite.h> 30 31 typedef struct _UserCtx *User; 32 struct _UserCtx { 33 PetscInt ptype; 34 DM pack; 35 Vec Uloc, Kloc; 36 }; 37 38 static PetscErrorCode FormFunctionLocal_U(User user, DMDALocalInfo *info, const PetscScalar u[], const PetscScalar k[], PetscScalar f[]) { 39 PetscReal hx = 1. / info->mx; 40 PetscInt i; 41 42 PetscFunctionBeginUser; 43 for (i = info->xs; i < info->xs + info->xm; i++) { 44 if (i == 0) f[i] = 1. / hx * u[i]; 45 else if (i == info->mx - 1) f[i] = 1. / hx * (u[i] - 1.0); 46 else f[i] = hx * ((k[i - 1] * (u[i] - u[i - 1]) - k[i] * (u[i + 1] - u[i])) / (hx * hx) - 1.0); 47 } 48 PetscFunctionReturn(0); 49 } 50 51 static PetscErrorCode FormFunctionLocal_K(User user, DMDALocalInfo *info, const PetscScalar u[], const PetscScalar k[], PetscScalar f[]) { 52 PetscReal hx = 1. / info->mx; 53 PetscInt i; 54 55 PetscFunctionBeginUser; 56 for (i = info->xs; i < info->xs + info->xm; i++) { 57 const PetscScalar ubar = 0.5 * (u[i + 1] + u[i]), gradu = (u[i + 1] - u[i]) / hx, g = 1. + gradu * gradu, w = 1. / (1. + ubar) + 1. / g; 58 f[i] = hx * (PetscExpScalar(k[i] - 1.0) + k[i] - 1. / w); 59 } 60 PetscFunctionReturn(0); 61 } 62 63 static PetscErrorCode FormFunction_All(SNES snes, Vec X, Vec F, void *ctx) { 64 User user = (User)ctx; 65 DM dau, dak; 66 DMDALocalInfo infou, infok; 67 PetscScalar *u, *k; 68 PetscScalar *fu, *fk; 69 Vec Uloc, Kloc, Fu, Fk; 70 71 PetscFunctionBeginUser; 72 PetscCall(DMCompositeGetEntries(user->pack, &dau, &dak)); 73 PetscCall(DMDAGetLocalInfo(dau, &infou)); 74 PetscCall(DMDAGetLocalInfo(dak, &infok)); 75 PetscCall(DMCompositeGetLocalVectors(user->pack, &Uloc, &Kloc)); 76 switch (user->ptype) { 77 case 0: 78 PetscCall(DMGlobalToLocalBegin(dau, X, INSERT_VALUES, Uloc)); 79 PetscCall(DMGlobalToLocalEnd(dau, X, INSERT_VALUES, Uloc)); 80 PetscCall(DMDAVecGetArray(dau, Uloc, &u)); 81 PetscCall(DMDAVecGetArray(dak, user->Kloc, &k)); 82 PetscCall(DMDAVecGetArray(dau, F, &fu)); 83 PetscCall(FormFunctionLocal_U(user, &infou, u, k, fu)); 84 PetscCall(DMDAVecRestoreArray(dau, F, &fu)); 85 PetscCall(DMDAVecRestoreArray(dau, Uloc, &u)); 86 PetscCall(DMDAVecRestoreArray(dak, user->Kloc, &k)); 87 break; 88 case 1: 89 PetscCall(DMGlobalToLocalBegin(dak, X, INSERT_VALUES, Kloc)); 90 PetscCall(DMGlobalToLocalEnd(dak, X, INSERT_VALUES, Kloc)); 91 PetscCall(DMDAVecGetArray(dau, user->Uloc, &u)); 92 PetscCall(DMDAVecGetArray(dak, Kloc, &k)); 93 PetscCall(DMDAVecGetArray(dak, F, &fk)); 94 PetscCall(FormFunctionLocal_K(user, &infok, u, k, fk)); 95 PetscCall(DMDAVecRestoreArray(dak, F, &fk)); 96 PetscCall(DMDAVecRestoreArray(dau, user->Uloc, &u)); 97 PetscCall(DMDAVecRestoreArray(dak, Kloc, &k)); 98 break; 99 case 2: 100 PetscCall(DMCompositeScatter(user->pack, X, Uloc, Kloc)); 101 PetscCall(DMDAVecGetArray(dau, Uloc, &u)); 102 PetscCall(DMDAVecGetArray(dak, Kloc, &k)); 103 PetscCall(DMCompositeGetAccess(user->pack, F, &Fu, &Fk)); 104 PetscCall(DMDAVecGetArray(dau, Fu, &fu)); 105 PetscCall(DMDAVecGetArray(dak, Fk, &fk)); 106 PetscCall(FormFunctionLocal_U(user, &infou, u, k, fu)); 107 PetscCall(FormFunctionLocal_K(user, &infok, u, k, fk)); 108 PetscCall(DMDAVecRestoreArray(dau, Fu, &fu)); 109 PetscCall(DMDAVecRestoreArray(dak, Fk, &fk)); 110 PetscCall(DMCompositeRestoreAccess(user->pack, F, &Fu, &Fk)); 111 PetscCall(DMDAVecRestoreArray(dau, Uloc, &u)); 112 PetscCall(DMDAVecRestoreArray(dak, Kloc, &k)); 113 break; 114 } 115 PetscCall(DMCompositeRestoreLocalVectors(user->pack, &Uloc, &Kloc)); 116 PetscFunctionReturn(0); 117 } 118 119 static PetscErrorCode FormJacobianLocal_U(User user, DMDALocalInfo *info, const PetscScalar u[], const PetscScalar k[], Mat Buu) { 120 PetscReal hx = 1. / info->mx; 121 PetscInt i; 122 123 PetscFunctionBeginUser; 124 for (i = info->xs; i < info->xs + info->xm; i++) { 125 PetscInt row = i - info->gxs, cols[] = {row - 1, row, row + 1}; 126 PetscScalar val = 1. / hx; 127 if (i == 0) { 128 PetscCall(MatSetValuesLocal(Buu, 1, &row, 1, &row, &val, INSERT_VALUES)); 129 } else if (i == info->mx - 1) { 130 PetscCall(MatSetValuesLocal(Buu, 1, &row, 1, &row, &val, INSERT_VALUES)); 131 } else { 132 PetscScalar vals[] = {-k[i - 1] / hx, (k[i - 1] + k[i]) / hx, -k[i] / hx}; 133 PetscCall(MatSetValuesLocal(Buu, 1, &row, 3, cols, vals, INSERT_VALUES)); 134 } 135 } 136 PetscFunctionReturn(0); 137 } 138 139 static PetscErrorCode FormJacobianLocal_K(User user, DMDALocalInfo *info, const PetscScalar u[], const PetscScalar k[], Mat Bkk) { 140 PetscReal hx = 1. / info->mx; 141 PetscInt i; 142 143 PetscFunctionBeginUser; 144 for (i = info->xs; i < info->xs + info->xm; i++) { 145 PetscInt row = i - info->gxs; 146 PetscScalar vals[] = {hx * (PetscExpScalar(k[i] - 1.) + (PetscScalar)1.)}; 147 PetscCall(MatSetValuesLocal(Bkk, 1, &row, 1, &row, vals, INSERT_VALUES)); 148 } 149 PetscFunctionReturn(0); 150 } 151 152 static PetscErrorCode FormJacobianLocal_UK(User user, DMDALocalInfo *info, DMDALocalInfo *infok, const PetscScalar u[], const PetscScalar k[], Mat Buk) { 153 PetscReal hx = 1. / info->mx; 154 PetscInt i; 155 PetscInt row, cols[2]; 156 PetscScalar vals[2]; 157 158 PetscFunctionBeginUser; 159 if (!Buk) PetscFunctionReturn(0); /* Not assembling this block */ 160 for (i = info->xs; i < info->xs + info->xm; i++) { 161 if (i == 0 || i == info->mx - 1) continue; 162 row = i - info->gxs; 163 cols[0] = i - 1 - infok->gxs; 164 vals[0] = (u[i] - u[i - 1]) / hx; 165 cols[1] = i - infok->gxs; 166 vals[1] = (u[i] - u[i + 1]) / hx; 167 PetscCall(MatSetValuesLocal(Buk, 1, &row, 2, cols, vals, INSERT_VALUES)); 168 } 169 PetscFunctionReturn(0); 170 } 171 172 static PetscErrorCode FormJacobianLocal_KU(User user, DMDALocalInfo *info, DMDALocalInfo *infok, const PetscScalar u[], const PetscScalar k[], Mat Bku) { 173 PetscInt i; 174 PetscReal hx = 1. / (info->mx - 1); 175 176 PetscFunctionBeginUser; 177 if (!Bku) PetscFunctionReturn(0); /* Not assembling this block */ 178 for (i = infok->xs; i < infok->xs + infok->xm; i++) { 179 PetscInt row = i - infok->gxs, cols[2]; 180 PetscScalar vals[2]; 181 const PetscScalar ubar = 0.5 * (u[i] + u[i + 1]), ubar_L = 0.5, ubar_R = 0.5, gradu = (u[i + 1] - u[i]) / hx, gradu_L = -1. / hx, gradu_R = 1. / hx, g = 1. + PetscSqr(gradu), g_gradu = 2. * gradu, w = 1. / (1. + ubar) + 1. / g, w_ubar = -1. / PetscSqr(1. + ubar), w_gradu = -g_gradu / PetscSqr(g), iw = 1. / w, iw_ubar = -w_ubar * PetscSqr(iw), iw_gradu = -w_gradu * PetscSqr(iw); 182 cols[0] = i - info->gxs; 183 vals[0] = -hx * (iw_ubar * ubar_L + iw_gradu * gradu_L); 184 cols[1] = i + 1 - info->gxs; 185 vals[1] = -hx * (iw_ubar * ubar_R + iw_gradu * gradu_R); 186 PetscCall(MatSetValuesLocal(Bku, 1, &row, 2, cols, vals, INSERT_VALUES)); 187 } 188 PetscFunctionReturn(0); 189 } 190 191 static PetscErrorCode FormJacobian_All(SNES snes, Vec X, Mat J, Mat B, void *ctx) { 192 User user = (User)ctx; 193 DM dau, dak; 194 DMDALocalInfo infou, infok; 195 PetscScalar *u, *k; 196 Vec Uloc, Kloc; 197 198 PetscFunctionBeginUser; 199 PetscCall(DMCompositeGetEntries(user->pack, &dau, &dak)); 200 PetscCall(DMDAGetLocalInfo(dau, &infou)); 201 PetscCall(DMDAGetLocalInfo(dak, &infok)); 202 PetscCall(DMCompositeGetLocalVectors(user->pack, &Uloc, &Kloc)); 203 switch (user->ptype) { 204 case 0: 205 PetscCall(DMGlobalToLocalBegin(dau, X, INSERT_VALUES, Uloc)); 206 PetscCall(DMGlobalToLocalEnd(dau, X, INSERT_VALUES, Uloc)); 207 PetscCall(DMDAVecGetArray(dau, Uloc, &u)); 208 PetscCall(DMDAVecGetArray(dak, user->Kloc, &k)); 209 PetscCall(FormJacobianLocal_U(user, &infou, u, k, B)); 210 PetscCall(DMDAVecRestoreArray(dau, Uloc, &u)); 211 PetscCall(DMDAVecRestoreArray(dak, user->Kloc, &k)); 212 break; 213 case 1: 214 PetscCall(DMGlobalToLocalBegin(dak, X, INSERT_VALUES, Kloc)); 215 PetscCall(DMGlobalToLocalEnd(dak, X, INSERT_VALUES, Kloc)); 216 PetscCall(DMDAVecGetArray(dau, user->Uloc, &u)); 217 PetscCall(DMDAVecGetArray(dak, Kloc, &k)); 218 PetscCall(FormJacobianLocal_K(user, &infok, u, k, B)); 219 PetscCall(DMDAVecRestoreArray(dau, user->Uloc, &u)); 220 PetscCall(DMDAVecRestoreArray(dak, Kloc, &k)); 221 break; 222 case 2: { 223 Mat Buu, Buk, Bku, Bkk; 224 PetscBool nest; 225 IS *is; 226 PetscCall(DMCompositeScatter(user->pack, X, Uloc, Kloc)); 227 PetscCall(DMDAVecGetArray(dau, Uloc, &u)); 228 PetscCall(DMDAVecGetArray(dak, Kloc, &k)); 229 PetscCall(DMCompositeGetLocalISs(user->pack, &is)); 230 PetscCall(MatGetLocalSubMatrix(B, is[0], is[0], &Buu)); 231 PetscCall(MatGetLocalSubMatrix(B, is[0], is[1], &Buk)); 232 PetscCall(MatGetLocalSubMatrix(B, is[1], is[0], &Bku)); 233 PetscCall(MatGetLocalSubMatrix(B, is[1], is[1], &Bkk)); 234 PetscCall(FormJacobianLocal_U(user, &infou, u, k, Buu)); 235 PetscCall(PetscObjectTypeCompare((PetscObject)B, MATNEST, &nest)); 236 if (!nest) { 237 /* 238 DMCreateMatrix_Composite() for a nested matrix does not generate off-block matrices that one can call MatSetValuesLocal() on, it just creates dummy 239 matrices with no entries; there cannot insert values into them. Commit b6480e041dd2293a65f96222772d68cdb4ed6306 240 changed Mat_Nest() from returning NULL pointers for these submatrices to dummy matrices because PCFIELDSPLIT could not 241 handle the returned null matrices. 242 */ 243 PetscCall(FormJacobianLocal_UK(user, &infou, &infok, u, k, Buk)); 244 PetscCall(FormJacobianLocal_KU(user, &infou, &infok, u, k, Bku)); 245 } 246 PetscCall(FormJacobianLocal_K(user, &infok, u, k, Bkk)); 247 PetscCall(MatRestoreLocalSubMatrix(B, is[0], is[0], &Buu)); 248 PetscCall(MatRestoreLocalSubMatrix(B, is[0], is[1], &Buk)); 249 PetscCall(MatRestoreLocalSubMatrix(B, is[1], is[0], &Bku)); 250 PetscCall(MatRestoreLocalSubMatrix(B, is[1], is[1], &Bkk)); 251 PetscCall(DMDAVecRestoreArray(dau, Uloc, &u)); 252 PetscCall(DMDAVecRestoreArray(dak, Kloc, &k)); 253 254 PetscCall(ISDestroy(&is[0])); 255 PetscCall(ISDestroy(&is[1])); 256 PetscCall(PetscFree(is)); 257 } break; 258 } 259 PetscCall(DMCompositeRestoreLocalVectors(user->pack, &Uloc, &Kloc)); 260 PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 261 PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 262 if (J != B) { 263 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 264 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 265 } 266 PetscFunctionReturn(0); 267 } 268 269 static PetscErrorCode FormInitial_Coupled(User user, Vec X) { 270 DM dau, dak; 271 DMDALocalInfo infou, infok; 272 Vec Xu, Xk; 273 PetscScalar *u, *k, hx; 274 PetscInt i; 275 276 PetscFunctionBeginUser; 277 PetscCall(DMCompositeGetEntries(user->pack, &dau, &dak)); 278 PetscCall(DMCompositeGetAccess(user->pack, X, &Xu, &Xk)); 279 PetscCall(DMDAVecGetArray(dau, Xu, &u)); 280 PetscCall(DMDAVecGetArray(dak, Xk, &k)); 281 PetscCall(DMDAGetLocalInfo(dau, &infou)); 282 PetscCall(DMDAGetLocalInfo(dak, &infok)); 283 hx = 1. / (infok.mx); 284 for (i = infou.xs; i < infou.xs + infou.xm; i++) u[i] = (PetscScalar)i * hx * (1. - (PetscScalar)i * hx); 285 for (i = infok.xs; i < infok.xs + infok.xm; i++) k[i] = 1.0 + 0.5 * PetscSinScalar(2 * PETSC_PI * i * hx); 286 PetscCall(DMDAVecRestoreArray(dau, Xu, &u)); 287 PetscCall(DMDAVecRestoreArray(dak, Xk, &k)); 288 PetscCall(DMCompositeRestoreAccess(user->pack, X, &Xu, &Xk)); 289 PetscCall(DMCompositeScatter(user->pack, X, user->Uloc, user->Kloc)); 290 PetscFunctionReturn(0); 291 } 292 293 int main(int argc, char *argv[]) { 294 DM dau, dak, pack; 295 const PetscInt *lxu; 296 PetscInt *lxk, m, sizes; 297 User user; 298 SNES snes; 299 Vec X, F, Xu, Xk, Fu, Fk; 300 Mat B; 301 IS *isg; 302 PetscBool pass_dm; 303 304 PetscFunctionBeginUser; 305 PetscCall(PetscInitialize(&argc, &argv, 0, help)); 306 PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 10, 1, 1, NULL, &dau)); 307 PetscCall(DMSetOptionsPrefix(dau, "u_")); 308 PetscCall(DMSetFromOptions(dau)); 309 PetscCall(DMSetUp(dau)); 310 PetscCall(DMDAGetOwnershipRanges(dau, &lxu, 0, 0)); 311 PetscCall(DMDAGetInfo(dau, 0, &m, 0, 0, &sizes, 0, 0, 0, 0, 0, 0, 0, 0)); 312 PetscCall(PetscMalloc1(sizes, &lxk)); 313 PetscCall(PetscArraycpy(lxk, lxu, sizes)); 314 lxk[0]--; 315 PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, m - 1, 1, 1, lxk, &dak)); 316 PetscCall(DMSetOptionsPrefix(dak, "k_")); 317 PetscCall(DMSetFromOptions(dak)); 318 PetscCall(DMSetUp(dak)); 319 PetscCall(PetscFree(lxk)); 320 321 PetscCall(DMCompositeCreate(PETSC_COMM_WORLD, &pack)); 322 PetscCall(DMSetOptionsPrefix(pack, "pack_")); 323 PetscCall(DMCompositeAddDM(pack, dau)); 324 PetscCall(DMCompositeAddDM(pack, dak)); 325 PetscCall(DMDASetFieldName(dau, 0, "u")); 326 PetscCall(DMDASetFieldName(dak, 0, "k")); 327 PetscCall(DMSetFromOptions(pack)); 328 329 PetscCall(DMCreateGlobalVector(pack, &X)); 330 PetscCall(VecDuplicate(X, &F)); 331 332 PetscCall(PetscNew(&user)); 333 334 user->pack = pack; 335 336 PetscCall(DMCompositeGetGlobalISs(pack, &isg)); 337 PetscCall(DMCompositeGetLocalVectors(pack, &user->Uloc, &user->Kloc)); 338 PetscCall(DMCompositeScatter(pack, X, user->Uloc, user->Kloc)); 339 340 PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Coupled problem options", "SNES"); 341 { 342 user->ptype = 0; 343 pass_dm = PETSC_TRUE; 344 345 PetscCall(PetscOptionsInt("-problem_type", "0: solve for u only, 1: solve for k only, 2: solve for both", 0, user->ptype, &user->ptype, NULL)); 346 PetscCall(PetscOptionsBool("-pass_dm", "Pass the packed DM to SNES to use when determining splits and forward into splits", 0, pass_dm, &pass_dm, NULL)); 347 } 348 PetscOptionsEnd(); 349 350 PetscCall(FormInitial_Coupled(user, X)); 351 352 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 353 switch (user->ptype) { 354 case 0: 355 PetscCall(DMCompositeGetAccess(pack, X, &Xu, 0)); 356 PetscCall(DMCompositeGetAccess(pack, F, &Fu, 0)); 357 PetscCall(DMCreateMatrix(dau, &B)); 358 PetscCall(SNESSetFunction(snes, Fu, FormFunction_All, user)); 359 PetscCall(SNESSetJacobian(snes, B, B, FormJacobian_All, user)); 360 PetscCall(SNESSetFromOptions(snes)); 361 PetscCall(SNESSetDM(snes, dau)); 362 PetscCall(SNESSolve(snes, NULL, Xu)); 363 PetscCall(DMCompositeRestoreAccess(pack, X, &Xu, 0)); 364 PetscCall(DMCompositeRestoreAccess(pack, F, &Fu, 0)); 365 break; 366 case 1: 367 PetscCall(DMCompositeGetAccess(pack, X, 0, &Xk)); 368 PetscCall(DMCompositeGetAccess(pack, F, 0, &Fk)); 369 PetscCall(DMCreateMatrix(dak, &B)); 370 PetscCall(SNESSetFunction(snes, Fk, FormFunction_All, user)); 371 PetscCall(SNESSetJacobian(snes, B, B, FormJacobian_All, user)); 372 PetscCall(SNESSetFromOptions(snes)); 373 PetscCall(SNESSetDM(snes, dak)); 374 PetscCall(SNESSolve(snes, NULL, Xk)); 375 PetscCall(DMCompositeRestoreAccess(pack, X, 0, &Xk)); 376 PetscCall(DMCompositeRestoreAccess(pack, F, 0, &Fk)); 377 break; 378 case 2: 379 PetscCall(DMCreateMatrix(pack, &B)); 380 /* This example does not correctly allocate off-diagonal blocks. These options allows new nonzeros (slow). */ 381 PetscCall(MatSetOption(B, MAT_NEW_NONZERO_LOCATION_ERR, PETSC_FALSE)); 382 PetscCall(MatSetOption(B, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 383 PetscCall(SNESSetFunction(snes, F, FormFunction_All, user)); 384 PetscCall(SNESSetJacobian(snes, B, B, FormJacobian_All, user)); 385 PetscCall(SNESSetFromOptions(snes)); 386 if (!pass_dm) { /* Manually provide index sets and names for the splits */ 387 KSP ksp; 388 PC pc; 389 PetscCall(SNESGetKSP(snes, &ksp)); 390 PetscCall(KSPGetPC(ksp, &pc)); 391 PetscCall(PCFieldSplitSetIS(pc, "u", isg[0])); 392 PetscCall(PCFieldSplitSetIS(pc, "k", isg[1])); 393 } else { 394 /* The same names come from the options prefix for dau and dak. This option can support geometric multigrid inside 395 * of splits, but it requires using a DM (perhaps your own implementation). */ 396 PetscCall(SNESSetDM(snes, pack)); 397 } 398 PetscCall(SNESSolve(snes, NULL, X)); 399 break; 400 } 401 PetscCall(VecViewFromOptions(X, NULL, "-view_sol")); 402 403 if (0) { 404 PetscInt col = 0; 405 PetscBool mult_dup = PETSC_FALSE, view_dup = PETSC_FALSE; 406 Mat D; 407 Vec Y; 408 409 PetscCall(PetscOptionsGetInt(NULL, 0, "-col", &col, 0)); 410 PetscCall(PetscOptionsGetBool(NULL, 0, "-mult_dup", &mult_dup, 0)); 411 PetscCall(PetscOptionsGetBool(NULL, 0, "-view_dup", &view_dup, 0)); 412 413 PetscCall(VecDuplicate(X, &Y)); 414 /* PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); */ 415 /* PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); */ 416 PetscCall(MatConvert(B, MATAIJ, MAT_INITIAL_MATRIX, &D)); 417 PetscCall(VecZeroEntries(X)); 418 PetscCall(VecSetValue(X, col, 1.0, INSERT_VALUES)); 419 PetscCall(VecAssemblyBegin(X)); 420 PetscCall(VecAssemblyEnd(X)); 421 PetscCall(MatMult(mult_dup ? D : B, X, Y)); 422 PetscCall(MatView(view_dup ? D : B, PETSC_VIEWER_STDOUT_WORLD)); 423 /* PetscCall(VecView(X,PETSC_VIEWER_STDOUT_WORLD)); */ 424 PetscCall(VecView(Y, PETSC_VIEWER_STDOUT_WORLD)); 425 PetscCall(MatDestroy(&D)); 426 PetscCall(VecDestroy(&Y)); 427 } 428 429 PetscCall(DMCompositeRestoreLocalVectors(pack, &user->Uloc, &user->Kloc)); 430 PetscCall(PetscFree(user)); 431 432 PetscCall(ISDestroy(&isg[0])); 433 PetscCall(ISDestroy(&isg[1])); 434 PetscCall(PetscFree(isg)); 435 PetscCall(VecDestroy(&X)); 436 PetscCall(VecDestroy(&F)); 437 PetscCall(MatDestroy(&B)); 438 PetscCall(DMDestroy(&dau)); 439 PetscCall(DMDestroy(&dak)); 440 PetscCall(DMDestroy(&pack)); 441 PetscCall(SNESDestroy(&snes)); 442 PetscCall(PetscFinalize()); 443 return 0; 444 } 445 446 /*TEST 447 448 test: 449 suffix: 0 450 nsize: 3 451 args: -u_da_grid_x 20 -snes_converged_reason -snes_monitor_short -problem_type 0 452 453 test: 454 suffix: 1 455 nsize: 3 456 args: -u_da_grid_x 20 -snes_converged_reason -snes_monitor_short -problem_type 1 457 458 test: 459 suffix: 2 460 nsize: 3 461 args: -u_da_grid_x 20 -snes_converged_reason -snes_monitor_short -problem_type 2 462 463 test: 464 suffix: 3 465 nsize: 3 466 args: -u_da_grid_x 20 -snes_converged_reason -snes_monitor_short -ksp_monitor_short -problem_type 2 -snes_mf_operator -pack_dm_mat_type {{aij nest}} -pc_type fieldsplit -pc_fieldsplit_dm_splits -pc_fieldsplit_type additive -fieldsplit_u_ksp_type gmres -fieldsplit_k_pc_type jacobi 467 468 test: 469 suffix: 4 470 nsize: 6 471 args: -u_da_grid_x 257 -snes_converged_reason -snes_monitor_short -ksp_monitor_short -problem_type 2 -snes_mf_operator -pack_dm_mat_type aij -pc_type fieldsplit -pc_fieldsplit_type multiplicative -fieldsplit_u_ksp_type gmres -fieldsplit_u_ksp_pc_side right -fieldsplit_u_pc_type mg -fieldsplit_u_pc_mg_levels 4 -fieldsplit_u_mg_levels_ksp_type richardson -fieldsplit_u_mg_levels_ksp_max_it 1 -fieldsplit_u_mg_levels_pc_type sor -fieldsplit_u_pc_mg_galerkin pmat -fieldsplit_u_ksp_converged_reason -fieldsplit_k_pc_type jacobi 472 requires: !single 473 474 test: 475 suffix: glvis_composite_da_1d 476 args: -u_da_grid_x 20 -problem_type 0 -snes_monitor_solution glvis: 477 478 test: 479 suffix: cuda 480 nsize: 1 481 requires: cuda 482 args: -u_da_grid_x 20 -snes_converged_reason -snes_monitor_short -problem_type 2 -pack_dm_mat_type aijcusparse -pack_dm_vec_type cuda 483 484 test: 485 suffix: viennacl 486 nsize: 1 487 requires: viennacl 488 args: -u_da_grid_x 20 -snes_converged_reason -snes_monitor_short -problem_type 2 -pack_dm_mat_type aijviennacl -pack_dm_vec_type viennacl 489 490 TEST*/ 491