1 static char help[] = "This example demonstrates the use of DMNetwork with subnetworks for solving a coupled nonlinear \n\ 2 electric power grid and water pipe problem.\n\ 3 The available solver options are in the ex1options file \n\ 4 and the data files are in the datafiles of subdirectories.\n\ 5 Run this program: mpiexec -n <n> ./ex1 \n\\n"; 6 /* 7 Example: 8 mpiexec -n 3 ./ex1 -petscpartitioner_type parmetis -dmnetwork_view draw -dmnetwork_view_distributed draw -dmnetwork_view_rank_range 0,1,2 9 mpiexec -n <n> ./ex1 -monitorIteration -monitorColor -power_snes_max_it 0 -water_snes_max_it 0 -coupled_snes_max_it 10 -draw_pause 5.0 10 */ 11 12 #include "power/power.h" 13 #include "water/water.h" 14 15 typedef struct { 16 UserCtx_Power appctx_power; 17 AppCtx_Water appctx_water; 18 PetscInt subsnes_id; /* snes solver id */ 19 PetscInt it; /* iteration number */ 20 Vec localXold; /* store previous solution, used by FormFunction_Dummy() */ 21 PetscBool monitorColor; 22 } UserCtx; 23 24 /* 25 UserMonitor -- called at the end of every SNES iteration via option `-monitorIteration' or `-monitorColor' 26 */ 27 PetscErrorCode UserMonitor(SNES snes, PetscInt its, PetscReal fnorm, void *appctx) 28 { 29 UserCtx *user = (UserCtx *)appctx; 30 PetscMPIInt rank; 31 MPI_Comm comm; 32 PetscInt it; 33 34 PetscFunctionBegin; 35 PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 36 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 37 38 PetscCall(SNESGetIterationNumber(snes, &it)); 39 if (rank == 0) { 40 PetscCall(SNESGetIterationNumber(snes, &it)); 41 if (user->subsnes_id == 0 || user->subsnes_id == 1) { 42 PetscCall(PetscPrintf(PETSC_COMM_SELF, " subsnes_id %" PetscInt_FMT ", it %" PetscInt_FMT ", fnorm %g\n", user->subsnes_id, it, (double)fnorm)); 43 } else { 44 PetscCall(PetscPrintf(PETSC_COMM_SELF, " coupled_snes_it %" PetscInt_FMT ", total_snes_it %" PetscInt_FMT ", fnorm %g\n", it, user->it, (double)fnorm)); 45 } 46 } 47 48 if (user->monitorColor) { 49 DM networkdm, dmcoords; 50 Vec F; 51 PetscInt v, vStart, vEnd, offset, gidx, rstart; 52 PetscReal *color; 53 PetscScalar *farr; 54 PetscBool ghost; 55 56 PetscCall(SNESGetDM(snes, &networkdm)); 57 PetscCall(DMGetCoordinateDM(networkdm, &dmcoords)); 58 59 PetscCall(SNESGetFunction(snes, &F, NULL, NULL)); 60 PetscCall(VecGetOwnershipRange(F, &rstart, NULL)); 61 62 PetscCall(VecGetArray(F, &farr)); 63 PetscCall(DMNetworkGetVertexRange(dmcoords, &vStart, &vEnd)); 64 65 PetscCall(PetscPrintf(MPI_COMM_WORLD, "\nColorPrint:\n")); 66 for (v = vStart; v < vEnd; v++) { 67 PetscCall(DMNetworkIsGhostVertex(networkdm, v, &ghost)); 68 PetscCall(DMNetworkGetComponent(dmcoords, v, 0, NULL, (void **)&color, NULL)); 69 PetscCall(DMNetworkGetGlobalVertexIndex(networkdm, v, &gidx)); 70 if (!ghost) { 71 PetscCall(DMNetworkGetGlobalVecOffset(networkdm, v, 0, &offset)); 72 *color = (PetscRealPart(farr[offset - rstart])); 73 } 74 PetscCall(PetscSynchronizedPrintf(MPI_COMM_WORLD, "[%d] v %" PetscInt_FMT ": color[%" PetscInt_FMT "] = %g\n", rank, gidx, offset - rstart, *color)); 75 } 76 PetscCall(PetscSynchronizedFlush(MPI_COMM_WORLD, NULL)); 77 PetscCall(VecRestoreArray(F, &farr)); 78 79 PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_CSV)); 80 PetscCall(DMView(networkdm, PETSC_VIEWER_DRAW_WORLD)); 81 PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 82 } 83 PetscFunctionReturn(PETSC_SUCCESS); 84 } 85 86 PetscErrorCode FormJacobian_subPower(SNES snes, Vec X, Mat J, Mat Jpre, void *appctx) 87 { 88 DM networkdm; 89 Vec localX; 90 PetscInt nv, ne, i, j, offset, nvar, row; 91 const PetscInt *vtx, *edges; 92 PetscBool ghostvtex; 93 PetscScalar one = 1.0; 94 PetscMPIInt rank; 95 MPI_Comm comm; 96 97 PetscFunctionBegin; 98 PetscCall(SNESGetDM(snes, &networkdm)); 99 PetscCall(DMGetLocalVector(networkdm, &localX)); 100 101 PetscCall(PetscObjectGetComm((PetscObject)networkdm, &comm)); 102 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 103 104 PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX)); 105 PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX)); 106 107 PetscCall(MatZeroEntries(J)); 108 109 /* Power subnetwork: copied from power/FormJacobian_Power() */ 110 PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges)); 111 PetscCall(FormJacobian_Power_private(networkdm, localX, J, nv, ne, vtx, edges, appctx)); 112 113 /* Water subnetwork: Identity */ 114 PetscCall(DMNetworkGetSubnetwork(networkdm, 1, &nv, &ne, &vtx, &edges)); 115 for (i = 0; i < nv; i++) { 116 PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[i], &ghostvtex)); 117 if (ghostvtex) continue; 118 119 PetscCall(DMNetworkGetGlobalVecOffset(networkdm, vtx[i], ALL_COMPONENTS, &offset)); 120 PetscCall(DMNetworkGetComponent(networkdm, vtx[i], ALL_COMPONENTS, NULL, NULL, &nvar)); 121 for (j = 0; j < nvar; j++) { 122 row = offset + j; 123 PetscCall(MatSetValues(J, 1, &row, 1, &row, &one, ADD_VALUES)); 124 } 125 } 126 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 127 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 128 129 PetscCall(DMRestoreLocalVector(networkdm, &localX)); 130 PetscFunctionReturn(PETSC_SUCCESS); 131 } 132 133 /* Dummy equation localF(X) = localX - localXold */ 134 PetscErrorCode FormFunction_Dummy(DM networkdm, Vec localX, Vec localF, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx) 135 { 136 const PetscScalar *xarr, *xoldarr; 137 PetscScalar *farr; 138 PetscInt i, j, offset, nvar; 139 PetscBool ghostvtex; 140 UserCtx *user = (UserCtx *)appctx; 141 Vec localXold = user->localXold; 142 143 PetscFunctionBegin; 144 PetscCall(VecGetArrayRead(localX, &xarr)); 145 PetscCall(VecGetArrayRead(localXold, &xoldarr)); 146 PetscCall(VecGetArray(localF, &farr)); 147 148 for (i = 0; i < nv; i++) { 149 PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[i], &ghostvtex)); 150 if (ghostvtex) continue; 151 152 PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[i], ALL_COMPONENTS, &offset)); 153 PetscCall(DMNetworkGetComponent(networkdm, vtx[i], ALL_COMPONENTS, NULL, NULL, &nvar)); 154 for (j = 0; j < nvar; j++) farr[offset + j] = xarr[offset + j] - xoldarr[offset + j]; 155 } 156 157 PetscCall(VecRestoreArrayRead(localX, &xarr)); 158 PetscCall(VecRestoreArrayRead(localXold, &xoldarr)); 159 PetscCall(VecRestoreArray(localF, &farr)); 160 PetscFunctionReturn(PETSC_SUCCESS); 161 } 162 163 PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *appctx) 164 { 165 DM networkdm; 166 Vec localX, localF; 167 PetscInt nv, ne, v; 168 const PetscInt *vtx, *edges; 169 PetscMPIInt rank; 170 MPI_Comm comm; 171 UserCtx *user = (UserCtx *)appctx; 172 UserCtx_Power appctx_power = (*user).appctx_power; 173 AppCtx_Water appctx_water = (*user).appctx_water; 174 175 PetscFunctionBegin; 176 PetscCall(SNESGetDM(snes, &networkdm)); 177 PetscCall(PetscObjectGetComm((PetscObject)networkdm, &comm)); 178 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 179 180 PetscCall(DMGetLocalVector(networkdm, &localX)); 181 PetscCall(DMGetLocalVector(networkdm, &localF)); 182 PetscCall(VecSet(F, 0.0)); 183 PetscCall(VecSet(localF, 0.0)); 184 185 PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX)); 186 PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX)); 187 188 /* Form Function for power subnetwork */ 189 PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges)); 190 if (user->subsnes_id == 1) { /* snes_water only */ 191 PetscCall(FormFunction_Dummy(networkdm, localX, localF, nv, ne, vtx, edges, user)); 192 } else { 193 PetscCall(FormFunction_Power(networkdm, localX, localF, nv, ne, vtx, edges, &appctx_power)); 194 } 195 196 /* Form Function for water subnetwork */ 197 PetscCall(DMNetworkGetSubnetwork(networkdm, 1, &nv, &ne, &vtx, &edges)); 198 if (user->subsnes_id == 0) { /* snes_power only */ 199 PetscCall(FormFunction_Dummy(networkdm, localX, localF, nv, ne, vtx, edges, user)); 200 } else { 201 PetscCall(FormFunction_Water(networkdm, localX, localF, nv, ne, vtx, edges, NULL)); 202 } 203 204 /* Illustrate how to access the coupling vertex of the subnetworks without doing anything to F yet */ 205 PetscCall(DMNetworkGetSharedVertices(networkdm, &nv, &vtx)); 206 for (v = 0; v < nv; v++) { 207 PetscInt key, ncomp, nvar, nconnedges, k, e, keye, goffset[3]; 208 void *component; 209 const PetscInt *connedges; 210 211 PetscCall(DMNetworkGetComponent(networkdm, vtx[v], ALL_COMPONENTS, NULL, NULL, &nvar)); 212 PetscCall(DMNetworkGetNumComponents(networkdm, vtx[v], &ncomp)); 213 /* printf(" [%d] coupling vertex[%" PetscInt_FMT "]: v %" PetscInt_FMT ", ncomp %" PetscInt_FMT "; nvar %" PetscInt_FMT "\n",rank,v,vtx[v], ncomp,nvar); */ 214 215 for (k = 0; k < ncomp; k++) { 216 PetscCall(DMNetworkGetComponent(networkdm, vtx[v], k, &key, &component, &nvar)); 217 PetscCall(DMNetworkGetGlobalVecOffset(networkdm, vtx[v], k, &goffset[k])); 218 219 /* Verify the coupling vertex is a powernet load vertex or a water vertex */ 220 switch (k) { 221 case 0: 222 PetscCheck(key == appctx_power.compkey_bus && nvar == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "key %" PetscInt_FMT " not a power bus vertex or nvar %" PetscInt_FMT " != 2", key, nvar); 223 break; 224 case 1: 225 PetscCheck(key == appctx_power.compkey_load && nvar == 0 && goffset[1] == goffset[0] + 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Not a power load vertex"); 226 break; 227 case 2: 228 PetscCheck(key == appctx_water.compkey_vtx && nvar == 1 && goffset[2] == goffset[1], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Not a water vertex"); 229 break; 230 default: 231 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "k %" PetscInt_FMT " is wrong", k); 232 } 233 /* printf(" [%d] coupling vertex[%" PetscInt_FMT "]: key %" PetscInt_FMT "; nvar %" PetscInt_FMT ", goffset %" PetscInt_FMT "\n",rank,v,key,nvar,goffset[k]); */ 234 } 235 236 /* Get its supporting edges */ 237 PetscCall(DMNetworkGetSupportingEdges(networkdm, vtx[v], &nconnedges, &connedges)); 238 /* printf("\n[%d] coupling vertex: nconnedges %" PetscInt_FMT "\n",rank,nconnedges); */ 239 for (k = 0; k < nconnedges; k++) { 240 e = connedges[k]; 241 PetscCall(DMNetworkGetNumComponents(networkdm, e, &ncomp)); 242 /* printf("\n [%d] connected edge[%" PetscInt_FMT "]=%" PetscInt_FMT " has ncomp %" PetscInt_FMT "\n",rank,k,e,ncomp); */ 243 PetscCall(DMNetworkGetComponent(networkdm, e, 0, &keye, &component, NULL)); 244 if (keye != appctx_water.compkey_edge) PetscCheck(keye == appctx_power.compkey_branch, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Not a power branch"); 245 } 246 } 247 248 PetscCall(DMRestoreLocalVector(networkdm, &localX)); 249 250 PetscCall(DMLocalToGlobalBegin(networkdm, localF, ADD_VALUES, F)); 251 PetscCall(DMLocalToGlobalEnd(networkdm, localF, ADD_VALUES, F)); 252 PetscCall(DMRestoreLocalVector(networkdm, &localF)); 253 #if 0 254 if (rank == 0) printf("F:\n"); 255 PetscCall(VecView(F,PETSC_VIEWER_STDOUT_WORLD)); 256 #endif 257 PetscFunctionReturn(PETSC_SUCCESS); 258 } 259 260 PetscErrorCode SetInitialGuess(DM networkdm, Vec X, void *appctx) 261 { 262 PetscInt nv, ne, i, j, ncomp, offset, key; 263 const PetscInt *vtx, *edges; 264 UserCtx *user = (UserCtx *)appctx; 265 Vec localX = user->localXold; 266 UserCtx_Power appctx_power = (*user).appctx_power; 267 AppCtx_Water appctx_water = (*user).appctx_water; 268 PetscBool ghost; 269 PetscScalar *xarr; 270 VERTEX_Power bus; 271 VERTEX_Water vertex; 272 void *component; 273 GEN gen; 274 275 PetscFunctionBegin; 276 PetscCall(VecSet(X, 0.0)); 277 PetscCall(VecSet(localX, 0.0)); 278 279 /* Set initial guess for power subnetwork */ 280 PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges)); 281 PetscCall(SetInitialGuess_Power(networkdm, localX, nv, ne, vtx, edges, &appctx_power)); 282 283 /* Set initial guess for water subnetwork */ 284 PetscCall(DMNetworkGetSubnetwork(networkdm, 1, &nv, &ne, &vtx, &edges)); 285 PetscCall(SetInitialGuess_Water(networkdm, localX, nv, ne, vtx, edges, NULL)); 286 287 /* Set initial guess at the coupling vertex */ 288 PetscCall(VecGetArray(localX, &xarr)); 289 PetscCall(DMNetworkGetSharedVertices(networkdm, &nv, &vtx)); 290 for (i = 0; i < nv; i++) { 291 PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[i], &ghost)); 292 if (ghost) continue; 293 294 PetscCall(DMNetworkGetNumComponents(networkdm, vtx[i], &ncomp)); 295 for (j = 0; j < ncomp; j++) { 296 PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[i], j, &offset)); 297 PetscCall(DMNetworkGetComponent(networkdm, vtx[i], j, &key, (void **)&component, NULL)); 298 if (key == appctx_power.compkey_bus) { 299 bus = (VERTEX_Power)(component); 300 xarr[offset] = bus->va * PETSC_PI / 180.0; 301 xarr[offset + 1] = bus->vm; 302 } else if (key == appctx_power.compkey_gen) { 303 gen = (GEN)(component); 304 if (!gen->status) continue; 305 xarr[offset + 1] = gen->vs; 306 } else if (key == appctx_water.compkey_vtx) { 307 vertex = (VERTEX_Water)(component); 308 if (vertex->type == VERTEX_TYPE_JUNCTION) { 309 xarr[offset] = 100; 310 } else if (vertex->type == VERTEX_TYPE_RESERVOIR) { 311 xarr[offset] = vertex->res.head; 312 } else { 313 xarr[offset] = vertex->tank.initlvl + vertex->tank.elev; 314 } 315 } 316 } 317 } 318 PetscCall(VecRestoreArray(localX, &xarr)); 319 320 PetscCall(DMLocalToGlobalBegin(networkdm, localX, ADD_VALUES, X)); 321 PetscCall(DMLocalToGlobalEnd(networkdm, localX, ADD_VALUES, X)); 322 PetscFunctionReturn(PETSC_SUCCESS); 323 } 324 325 /* Set coordinates */ 326 static PetscErrorCode CoordinateVecSetUp(DM dmcoords, Vec coords) 327 { 328 PetscInt i, gidx, offset, v, nv, Nsubnet; 329 const PetscInt *vtx; 330 PetscScalar *carray; 331 PetscReal *color; 332 333 PetscFunctionBeginUser; 334 PetscCall(VecGetArrayWrite(coords, &carray)); 335 PetscCall(DMNetworkGetNumSubNetworks(dmcoords, NULL, &Nsubnet)); 336 for (i = 0; i < Nsubnet; i++) { 337 PetscCall(DMNetworkGetSubnetwork(dmcoords, i, &nv, NULL, &vtx, NULL)); 338 for (v = 0; v < nv; v++) { 339 PetscCall(DMNetworkGetGlobalVertexIndex(dmcoords, vtx[v], &gidx)); 340 PetscCall(DMNetworkGetLocalVecOffset(dmcoords, vtx[v], 0, &offset)); 341 PetscCall(DMNetworkGetComponent(dmcoords, vtx[v], 0, NULL, (void **)&color, NULL)); 342 *color = 0.0; 343 switch (gidx) { 344 case 0: 345 carray[offset] = -1.0; 346 carray[offset + 1] = -1.0; 347 break; 348 case 1: 349 carray[offset] = -2.0; 350 carray[offset + 1] = 2.0; 351 break; 352 case 2: 353 carray[offset] = 0.0; 354 carray[offset + 1] = 2.0; 355 break; 356 case 3: 357 carray[offset] = -1.0; 358 carray[offset + 1] = 0.0; 359 break; 360 case 4: 361 carray[offset] = 0.0; 362 carray[offset + 1] = 0.0; 363 break; 364 case 5: 365 carray[offset] = 0.0; 366 carray[offset + 1] = 1.0; 367 break; 368 case 6: 369 carray[offset] = -1.0; 370 carray[offset + 1] = 1.0; 371 break; 372 case 7: 373 carray[offset] = -2.0; 374 carray[offset + 1] = 1.0; 375 break; 376 case 8: 377 carray[offset] = -2.0; 378 carray[offset + 1] = 0.0; 379 break; 380 case 9: 381 carray[offset] = 1.0; 382 carray[offset + 1] = 0.0; 383 break; 384 case 10: 385 carray[offset] = 1.0; 386 carray[offset + 1] = -1.0; 387 break; 388 case 11: 389 carray[offset] = 2.0; 390 carray[offset + 1] = -1.0; 391 break; 392 case 12: 393 carray[offset] = 2.0; 394 carray[offset + 1] = 0.0; 395 break; 396 case 13: 397 carray[offset] = 0.0; 398 carray[offset + 1] = -1.0; 399 break; 400 case 14: 401 carray[offset] = 2.0; 402 carray[offset + 1] = 1.0; 403 break; 404 default: 405 PetscCheck(gidx < 15 && gidx > -1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "gidx %" PetscInt_FMT "must between 0 and 14", gidx); 406 } 407 } 408 } 409 PetscCall(VecRestoreArrayWrite(coords, &carray)); 410 PetscFunctionReturn(PETSC_SUCCESS); 411 } 412 413 static PetscErrorCode CoordinatePrint(DM dm) 414 { 415 DM dmcoords; 416 PetscInt cdim, v, off, vglobal, vStart, vEnd; 417 const PetscScalar *carray; 418 Vec coords; 419 MPI_Comm comm; 420 PetscMPIInt rank; 421 422 PetscFunctionBegin; 423 /* get info from dm */ 424 PetscCall(DMGetCoordinateDim(dm, &cdim)); 425 PetscCall(DMGetCoordinatesLocal(dm, &coords)); 426 427 PetscCall(DMGetCoordinateDM(dm, &dmcoords)); 428 PetscCall(PetscObjectGetComm((PetscObject)dmcoords, &comm)); 429 PetscCallMPI(MPI_Comm_rank(comm, &rank)); 430 431 /* print coordinates from dmcoords */ 432 PetscCall(PetscPrintf(MPI_COMM_WORLD, "\nCoordinatePrint, cdim %" PetscInt_FMT ":\n", cdim)); 433 PetscCall(PetscSynchronizedPrintf(MPI_COMM_WORLD, "[%d]\n", rank)); 434 435 PetscCall(DMNetworkGetVertexRange(dmcoords, &vStart, &vEnd)); 436 PetscCall(VecGetArrayRead(coords, &carray)); 437 for (v = vStart; v < vEnd; v++) { 438 PetscCall(DMNetworkGetLocalVecOffset(dmcoords, v, 0, &off)); 439 PetscCall(DMNetworkGetGlobalVertexIndex(dmcoords, v, &vglobal)); 440 switch (cdim) { 441 case 2: 442 PetscCall(PetscSynchronizedPrintf(MPI_COMM_WORLD, "Vertex: %" PetscInt_FMT ", x = %f y = %f \n", vglobal, (double)PetscRealPart(carray[off]), (double)PetscRealPart(carray[off + 1]))); 443 break; 444 default: 445 PetscCheck(cdim == 2, MPI_COMM_WORLD, PETSC_ERR_SUP, "Only supports Network embedding dimension of 2, not supplied %" PetscInt_FMT, cdim); 446 break; 447 } 448 } 449 PetscCall(PetscSynchronizedFlush(MPI_COMM_WORLD, NULL)); 450 PetscCall(VecRestoreArrayRead(coords, &carray)); 451 PetscFunctionReturn(PETSC_SUCCESS); 452 } 453 454 int main(int argc, char **argv) 455 { 456 DM networkdm, dmcoords; 457 PetscMPIInt rank, size; 458 PetscInt Nsubnet = 2, numVertices[2], numEdges[2], i, j, nv, ne, it_max = 10; 459 PetscInt vStart, vEnd, compkey; 460 const PetscInt *vtx, *edges; 461 PetscReal *color; 462 Vec X, F, coords; 463 SNES snes, snes_power, snes_water; 464 Mat Jac; 465 PetscBool ghost, viewJ = PETSC_FALSE, viewX = PETSC_FALSE, test = PETSC_FALSE, distribute = PETSC_TRUE, flg, printCoord = PETSC_FALSE, viewCSV = PETSC_FALSE, monitorIt = PETSC_FALSE; 466 UserCtx user; 467 SNESConvergedReason reason; 468 #if defined(PETSC_USE_LOG) 469 PetscLogStage stage[4]; 470 #endif 471 472 /* Power subnetwork */ 473 UserCtx_Power *appctx_power = &user.appctx_power; 474 char pfdata_file[PETSC_MAX_PATH_LEN] = "power/case9.m"; 475 PFDATA *pfdata = NULL; 476 PetscInt genj, loadj, *edgelist_power = NULL, power_netnum; 477 PetscScalar Sbase = 0.0; 478 479 /* Water subnetwork */ 480 AppCtx_Water *appctx_water = &user.appctx_water; 481 WATERDATA *waterdata = NULL; 482 char waterdata_file[PETSC_MAX_PATH_LEN] = "water/sample1.inp"; 483 PetscInt *edgelist_water = NULL, water_netnum; 484 485 /* Shared vertices between subnetworks */ 486 PetscInt power_svtx, water_svtx; 487 488 PetscFunctionBeginUser; 489 PetscCall(PetscInitialize(&argc, &argv, "ex1options", help)); 490 PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 491 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 492 493 /* (1) Read Data - Only rank 0 reads the data */ 494 PetscCall(PetscLogStageRegister("Read Data", &stage[0])); 495 PetscCall(PetscLogStagePush(stage[0])); 496 497 for (i = 0; i < Nsubnet; i++) { 498 numVertices[i] = 0; 499 numEdges[i] = 0; 500 } 501 502 /* All processes READ THE DATA FOR THE FIRST SUBNETWORK: Electric Power Grid */ 503 /* Used for shared vertex, because currently the coupling info must be available in all processes!!! */ 504 PetscCall(PetscOptionsGetString(NULL, NULL, "-pfdata", pfdata_file, PETSC_MAX_PATH_LEN - 1, NULL)); 505 PetscCall(PetscNew(&pfdata)); 506 PetscCall(PFReadMatPowerData(pfdata, pfdata_file)); 507 if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Power network: nb = %" PetscInt_FMT ", ngen = %" PetscInt_FMT ", nload = %" PetscInt_FMT ", nbranch = %" PetscInt_FMT "\n", pfdata->nbus, pfdata->ngen, pfdata->nload, pfdata->nbranch)); 508 Sbase = pfdata->sbase; 509 if (rank == 0) { /* proc[0] will create Electric Power Grid */ 510 numEdges[0] = pfdata->nbranch; 511 numVertices[0] = pfdata->nbus; 512 513 PetscCall(PetscMalloc1(2 * numEdges[0], &edgelist_power)); 514 PetscCall(GetListofEdges_Power(pfdata, edgelist_power)); 515 } 516 /* Broadcast power Sbase to all processors */ 517 PetscCallMPI(MPI_Bcast(&Sbase, 1, MPIU_SCALAR, 0, PETSC_COMM_WORLD)); 518 appctx_power->Sbase = Sbase; 519 appctx_power->jac_error = PETSC_FALSE; 520 /* If external option activated. Introduce error in jacobian */ 521 PetscCall(PetscOptionsHasName(NULL, NULL, "-jac_error", &appctx_power->jac_error)); 522 523 /* All processes READ THE DATA FOR THE SECOND SUBNETWORK: Water */ 524 /* Used for shared vertex, because currently the coupling info must be available in all processes!!! */ 525 PetscCall(PetscNew(&waterdata)); 526 PetscCall(PetscOptionsGetString(NULL, NULL, "-waterdata", waterdata_file, PETSC_MAX_PATH_LEN - 1, NULL)); 527 PetscCall(WaterReadData(waterdata, waterdata_file)); 528 if (size == 1 || (size > 1 && rank == 1)) { 529 PetscCall(PetscCalloc1(2 * waterdata->nedge, &edgelist_water)); 530 PetscCall(GetListofEdges_Water(waterdata, edgelist_water)); 531 numEdges[1] = waterdata->nedge; 532 numVertices[1] = waterdata->nvertex; 533 } 534 PetscCall(PetscLogStagePop()); 535 536 /* (2) Create a network consist of two subnetworks */ 537 PetscCall(PetscLogStageRegister("Net Setup", &stage[1])); 538 PetscCall(PetscLogStagePush(stage[1])); 539 540 PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewCSV", &viewCSV, NULL)); 541 PetscCall(PetscOptionsGetBool(NULL, NULL, "-printCoord", &printCoord, NULL)); 542 PetscCall(PetscOptionsGetBool(NULL, NULL, "-test", &test, NULL)); 543 PetscCall(PetscOptionsGetBool(NULL, NULL, "-distribute", &distribute, NULL)); 544 545 /* Create an empty network object */ 546 PetscCall(DMNetworkCreate(PETSC_COMM_WORLD, &networkdm)); 547 548 /* Register the components in the network */ 549 PetscCall(DMNetworkRegisterComponent(networkdm, "branchstruct", sizeof(struct _p_EDGE_Power), &appctx_power->compkey_branch)); 550 PetscCall(DMNetworkRegisterComponent(networkdm, "busstruct", sizeof(struct _p_VERTEX_Power), &appctx_power->compkey_bus)); 551 PetscCall(DMNetworkRegisterComponent(networkdm, "genstruct", sizeof(struct _p_GEN), &appctx_power->compkey_gen)); 552 PetscCall(DMNetworkRegisterComponent(networkdm, "loadstruct", sizeof(struct _p_LOAD), &appctx_power->compkey_load)); 553 554 PetscCall(DMNetworkRegisterComponent(networkdm, "edge_water", sizeof(struct _p_EDGE_Water), &appctx_water->compkey_edge)); 555 PetscCall(DMNetworkRegisterComponent(networkdm, "vertex_water", sizeof(struct _p_VERTEX_Water), &appctx_water->compkey_vtx)); 556 557 PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] Total local nvertices %" PetscInt_FMT " + %" PetscInt_FMT " = %" PetscInt_FMT ", nedges %" PetscInt_FMT " + %" PetscInt_FMT " = %" PetscInt_FMT "\n", rank, numVertices[0], numVertices[1], numVertices[0] + numVertices[1], numEdges[0], numEdges[1], numEdges[0] + numEdges[1])); 558 PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT)); 559 560 PetscCall(DMNetworkSetNumSubNetworks(networkdm, PETSC_DECIDE, Nsubnet)); 561 PetscCall(DMNetworkAddSubnetwork(networkdm, "power", numEdges[0], edgelist_power, &power_netnum)); 562 PetscCall(DMNetworkAddSubnetwork(networkdm, "water", numEdges[1], edgelist_water, &water_netnum)); 563 564 /* vertex subnet[0].4 shares with vertex subnet[1].0 */ 565 power_svtx = 4; 566 water_svtx = 0; 567 PetscCall(DMNetworkAddSharedVertices(networkdm, power_netnum, water_netnum, 1, &power_svtx, &water_svtx)); 568 569 /* Set up the network layout */ 570 PetscCall(DMNetworkLayoutSetUp(networkdm)); 571 572 /* ADD VARIABLES AND COMPONENTS FOR THE POWER SUBNETWORK */ 573 genj = 0; 574 loadj = 0; 575 PetscCall(DMNetworkGetSubnetwork(networkdm, power_netnum, &nv, &ne, &vtx, &edges)); 576 577 for (i = 0; i < ne; i++) PetscCall(DMNetworkAddComponent(networkdm, edges[i], appctx_power->compkey_branch, &pfdata->branch[i], 0)); 578 579 for (i = 0; i < nv; i++) { 580 PetscCall(DMNetworkIsSharedVertex(networkdm, vtx[i], &flg)); 581 if (flg) continue; 582 583 PetscCall(DMNetworkAddComponent(networkdm, vtx[i], appctx_power->compkey_bus, &pfdata->bus[i], 2)); 584 if (pfdata->bus[i].ngen) { 585 for (j = 0; j < pfdata->bus[i].ngen; j++) PetscCall(DMNetworkAddComponent(networkdm, vtx[i], appctx_power->compkey_gen, &pfdata->gen[genj++], 0)); 586 } 587 if (pfdata->bus[i].nload) { 588 for (j = 0; j < pfdata->bus[i].nload; j++) PetscCall(DMNetworkAddComponent(networkdm, vtx[i], appctx_power->compkey_load, &pfdata->load[loadj++], 0)); 589 } 590 } 591 592 /* ADD VARIABLES AND COMPONENTS FOR THE WATER SUBNETWORK */ 593 PetscCall(DMNetworkGetSubnetwork(networkdm, water_netnum, &nv, &ne, &vtx, &edges)); 594 for (i = 0; i < ne; i++) PetscCall(DMNetworkAddComponent(networkdm, edges[i], appctx_water->compkey_edge, &waterdata->edge[i], 0)); 595 596 for (i = 0; i < nv; i++) { 597 PetscCall(DMNetworkIsSharedVertex(networkdm, vtx[i], &flg)); 598 if (flg) continue; 599 600 PetscCall(DMNetworkAddComponent(networkdm, vtx[i], appctx_water->compkey_vtx, &waterdata->vertex[i], 1)); 601 } 602 603 /* ADD VARIABLES AND COMPONENTS AT THE SHARED VERTEX: net[0].4 coupls with net[1].0 -- owning and all ghost ranks of the vertex do this */ 604 PetscCall(DMNetworkGetSharedVertices(networkdm, &nv, &vtx)); 605 for (i = 0; i < nv; i++) { 606 /* power */ 607 PetscCall(DMNetworkAddComponent(networkdm, vtx[i], appctx_power->compkey_bus, &pfdata->bus[4], 2)); 608 /* bus[4] is a load, add its component */ 609 PetscCall(DMNetworkAddComponent(networkdm, vtx[i], appctx_power->compkey_load, &pfdata->load[0], 0)); 610 611 /* water */ 612 PetscCall(DMNetworkAddComponent(networkdm, vtx[i], appctx_water->compkey_vtx, &waterdata->vertex[0], 1)); 613 } 614 615 /* Set coordinates for visualization */ 616 PetscCall(DMSetCoordinateDim(networkdm, 2)); 617 PetscCall(DMGetCoordinateDM(networkdm, &dmcoords)); 618 PetscCall(DMNetworkGetVertexRange(dmcoords, &vStart, &vEnd)); 619 620 PetscCall(PetscCalloc1(vEnd - vStart, &color)); 621 PetscCall(DMNetworkRegisterComponent(dmcoords, "coordinate&color", sizeof(PetscReal), &compkey)); 622 for (i = vStart; i < vEnd; i++) PetscCall(DMNetworkAddComponent(dmcoords, i, compkey, &color[i - vStart], 2)); 623 PetscCall(DMNetworkFinalizeComponents(dmcoords)); 624 625 PetscCall(DMCreateLocalVector(dmcoords, &coords)); 626 PetscCall(DMSetCoordinatesLocal(networkdm, coords)); /* set/get coords to/from networkdm */ 627 PetscCall(CoordinateVecSetUp(dmcoords, coords)); 628 if (printCoord) PetscCall(CoordinatePrint(networkdm)); 629 630 /* Set up DM for use */ 631 PetscCall(DMSetUp(networkdm)); 632 633 /* Free user objects */ 634 PetscCall(PetscFree(edgelist_power)); 635 PetscCall(PetscFree(pfdata->bus)); 636 PetscCall(PetscFree(pfdata->gen)); 637 PetscCall(PetscFree(pfdata->branch)); 638 PetscCall(PetscFree(pfdata->load)); 639 PetscCall(PetscFree(pfdata)); 640 641 PetscCall(PetscFree(edgelist_water)); 642 PetscCall(PetscFree(waterdata->vertex)); 643 PetscCall(PetscFree(waterdata->edge)); 644 PetscCall(PetscFree(waterdata)); 645 646 /* Re-distribute networkdm to multiple processes for better job balance */ 647 if (distribute) { 648 PetscCall(DMNetworkDistribute(&networkdm, 0)); 649 650 if (printCoord) PetscCall(CoordinatePrint(networkdm)); 651 if (viewCSV) { /* CSV View of network with coordinates */ 652 PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_CSV)); 653 PetscCall(DMView(networkdm, PETSC_VIEWER_STDOUT_WORLD)); 654 PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 655 } 656 } 657 PetscCall(VecDestroy(&coords)); 658 659 /* Test DMNetworkGetSubnetwork() and DMNetworkGetSubnetworkSharedVertices() */ 660 if (test) { 661 PetscInt v, gidx; 662 PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 663 for (i = 0; i < Nsubnet; i++) { 664 PetscCall(DMNetworkGetSubnetwork(networkdm, i, &nv, &ne, &vtx, &edges)); 665 PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] After distribute, subnet[%" PetscInt_FMT "] ne %" PetscInt_FMT ", nv %" PetscInt_FMT "\n", rank, i, ne, nv)); 666 PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 667 668 for (v = 0; v < nv; v++) { 669 PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[v], &ghost)); 670 PetscCall(DMNetworkGetGlobalVertexIndex(networkdm, vtx[v], &gidx)); 671 PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] subnet[%" PetscInt_FMT "] v %" PetscInt_FMT " %" PetscInt_FMT "; ghost %d\n", rank, i, vtx[v], gidx, ghost)); 672 } 673 PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 674 } 675 PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 676 677 PetscCall(DMNetworkGetSharedVertices(networkdm, &nv, &vtx)); 678 PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] After distribute, num of shared vertices nsv = %" PetscInt_FMT "\n", rank, nv)); 679 for (v = 0; v < nv; v++) { 680 PetscCall(DMNetworkGetGlobalVertexIndex(networkdm, vtx[v], &gidx)); 681 PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] sv %" PetscInt_FMT ", gidx=%" PetscInt_FMT "\n", rank, vtx[v], gidx)); 682 } 683 PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 684 } 685 686 /* Create solution vector X */ 687 PetscCall(DMCreateGlobalVector(networkdm, &X)); 688 PetscCall(VecDuplicate(X, &F)); 689 PetscCall(DMGetLocalVector(networkdm, &user.localXold)); 690 PetscCall(PetscLogStagePop()); 691 692 /* (3) Setup Solvers */ 693 PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewJ", &viewJ, NULL)); 694 PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewX", &viewX, NULL)); 695 696 PetscCall(PetscLogStageRegister("SNES Setup", &stage[2])); 697 PetscCall(PetscLogStagePush(stage[2])); 698 699 PetscCall(SetInitialGuess(networkdm, X, &user)); 700 701 /* Create coupled snes */ 702 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "SNES_coupled setup ......\n")); 703 user.subsnes_id = Nsubnet; 704 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 705 PetscCall(SNESSetDM(snes, networkdm)); 706 PetscCall(SNESSetOptionsPrefix(snes, "coupled_")); 707 PetscCall(SNESSetFunction(snes, F, FormFunction, &user)); 708 /* set maxit=1 which can be changed via option '-coupled_snes_max_it <>', see ex1options */ 709 PetscCall(SNESSetTolerances(snes, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 1, PETSC_DEFAULT)); 710 PetscCall(SNESSetFromOptions(snes)); 711 712 if (viewJ) { 713 /* View Jac structure */ 714 PetscCall(SNESGetJacobian(snes, &Jac, NULL, NULL, NULL)); 715 PetscCall(MatView(Jac, PETSC_VIEWER_DRAW_WORLD)); 716 } 717 718 if (viewX) { 719 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Solution:\n")); 720 PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD)); 721 } 722 723 if (viewJ) { 724 /* View assembled Jac */ 725 PetscCall(MatView(Jac, PETSC_VIEWER_DRAW_WORLD)); 726 } 727 728 /* Create snes_power */ 729 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "SNES_power setup ......\n")); 730 user.subsnes_id = 0; 731 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes_power)); 732 PetscCall(SNESSetDM(snes_power, networkdm)); 733 PetscCall(SNESSetOptionsPrefix(snes_power, "power_")); 734 PetscCall(SNESSetFunction(snes_power, F, FormFunction, &user)); 735 /* set maxit=1 which can be changed via option '-power_snes_max_it <>', see ex1options */ 736 PetscCall(SNESSetTolerances(snes_power, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 1, PETSC_DEFAULT)); 737 738 /* Use user-provide Jacobian */ 739 PetscCall(DMCreateMatrix(networkdm, &Jac)); 740 PetscCall(SNESSetJacobian(snes_power, Jac, Jac, FormJacobian_subPower, &user)); 741 PetscCall(SNESSetFromOptions(snes_power)); 742 743 if (viewX) { 744 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Power Solution:\n")); 745 PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD)); 746 } 747 if (viewJ) { 748 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Power Jac:\n")); 749 PetscCall(SNESGetJacobian(snes_power, &Jac, NULL, NULL, NULL)); 750 PetscCall(MatView(Jac, PETSC_VIEWER_DRAW_WORLD)); 751 } 752 753 /* Create snes_water */ 754 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "SNES_water setup......\n")); 755 756 user.subsnes_id = 1; 757 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes_water)); 758 PetscCall(SNESSetDM(snes_water, networkdm)); 759 PetscCall(SNESSetOptionsPrefix(snes_water, "water_")); 760 PetscCall(SNESSetFunction(snes_water, F, FormFunction, &user)); 761 /* set maxit=1 which can be changed via option '-water_snes_max_it <>', see ex1options */ 762 PetscCall(SNESSetTolerances(snes_water, PETSC_DEFAULT, PETSC_DEFAULT, PETSC_DEFAULT, 1, PETSC_DEFAULT)); 763 PetscCall(SNESSetFromOptions(snes_water)); 764 765 if (viewX) { 766 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Water Solution:\n")); 767 PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD)); 768 } 769 770 /* Monitor snes, snes_power and snes_water iterations */ 771 PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitorIteration", &monitorIt, NULL)); 772 user.monitorColor = PETSC_FALSE; 773 PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitorColor", &user.monitorColor, NULL)); 774 if (user.monitorColor) monitorIt = PETSC_TRUE; /* require installation of pandas and matplotlib */ 775 if (monitorIt) { 776 PetscCall(SNESMonitorSet(snes_power, UserMonitor, &user, NULL)); 777 PetscCall(SNESMonitorSet(snes_water, UserMonitor, &user, NULL)); 778 PetscCall(SNESMonitorSet(snes, UserMonitor, &user, NULL)); 779 } 780 PetscCall(PetscLogStagePop()); 781 782 /* (4) Solve: we must update user.localXold after each call of SNESSolve(). 783 See "PETSc DMNetwork: A Library for Scalable Network PDE-Based Multiphysics Simulations", 784 https://dl.acm.org/doi/10.1145/3344587 785 */ 786 PetscCall(PetscLogStageRegister("SNES Solve", &stage[3])); 787 PetscCall(PetscLogStagePush(stage[3])); 788 user.it = 0; /* total_snes_it */ 789 reason = SNES_DIVERGED_DTOL; 790 while (user.it < it_max && (PetscInt)reason < 0) { 791 user.subsnes_id = 0; 792 PetscCall(SNESSolve(snes_power, NULL, X)); 793 PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, user.localXold)); 794 PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, user.localXold)); 795 796 user.subsnes_id = 1; 797 PetscCall(SNESSolve(snes_water, NULL, X)); 798 PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, user.localXold)); 799 PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, user.localXold)); 800 801 user.subsnes_id = Nsubnet; 802 PetscCall(SNESSolve(snes, NULL, X)); 803 PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, user.localXold)); 804 PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, user.localXold)); 805 806 PetscCall(SNESGetConvergedReason(snes, &reason)); 807 user.it++; 808 } 809 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Coupled_SNES converged in %" PetscInt_FMT " iterations\n", user.it)); 810 if (viewX) { 811 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final Solution:\n")); 812 PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD)); 813 } 814 PetscCall(PetscLogStagePop()); 815 816 /* Free objects */ 817 /* -------------*/ 818 PetscCall(PetscFree(color)); 819 PetscCall(VecDestroy(&X)); 820 PetscCall(VecDestroy(&F)); 821 PetscCall(DMRestoreLocalVector(networkdm, &user.localXold)); 822 823 PetscCall(SNESDestroy(&snes)); 824 PetscCall(MatDestroy(&Jac)); 825 PetscCall(SNESDestroy(&snes_power)); 826 PetscCall(SNESDestroy(&snes_water)); 827 828 PetscCall(DMDestroy(&networkdm)); 829 PetscCall(PetscFinalize()); 830 return 0; 831 } 832 833 /*TEST 834 835 build: 836 requires: !complex double defined(PETSC_HAVE_ATTRIBUTEALIGNED) 837 depends: power/PFReadData.c power/pffunctions.c water/waterreaddata.c water/waterfunctions.c 838 839 test: 840 args: -options_left no -dmnetwork_view -fp_trap 0 841 localrunfiles: ex1options power/case9.m water/sample1.inp 842 output_file: output/ex1.out 843 844 test: 845 suffix: 2 846 nsize: 3 847 args: -options_left no -petscpartitioner_type parmetis -fp_trap 0 848 localrunfiles: ex1options power/case9.m water/sample1.inp 849 output_file: output/ex1_2.out 850 requires: parmetis 851 852 test: 853 suffix: 3 854 nsize: 3 855 args: -options_left no -distribute false -fp_trap 0 856 localrunfiles: ex1options power/case9.m water/sample1.inp 857 output_file: output/ex1_2.out 858 859 test: 860 suffix: 4 861 nsize: 4 862 args: -options_left no -petscpartitioner_type simple -dmnetwork_view -dmnetwork_view_distributed -fp_trap 0 863 localrunfiles: ex1options power/case9.m water/sample1.inp 864 output_file: output/ex1_4.out 865 866 test: 867 suffix: 5 868 args: -options_left no -viewCSV -fp_trap 0 869 localrunfiles: ex1options power/case9.m water/sample1.inp 870 output_file: output/ex1_5.out 871 872 test: 873 suffix: 6 874 nsize: 3 875 args: -options_left no -petscpartitioner_type parmetis -dmnetwork_view_distributed draw:null -fp_trap 0 876 localrunfiles: ex1options power/case9.m water/sample1.inp 877 output_file: output/ex1_2.out 878 requires: parmetis 879 880 TEST*/ 881