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