1c4762a1bSJed Brown static char help[] = "This example demonstrates the use of DMNetwork interface with subnetworks for solving a nonlinear electric power grid problem.\n\ 2c4762a1bSJed Brown The available solver options are in the poweroptions file and the data files are in the datafiles directory.\n\ 3c4762a1bSJed Brown The data file format used is from the MatPower package (http://www.pserc.cornell.edu//matpower/).\n\ 4c4762a1bSJed Brown This example shows the use of subnetwork feature in DMNetwork. It creates duplicates of the same network which are treated as subnetworks.\n\ 5c4762a1bSJed Brown Run this program: mpiexec -n <n> ./pf2\n\ 6c4762a1bSJed Brown mpiexec -n <n> ./pf2 \n"; 7c4762a1bSJed Brown 8c4762a1bSJed Brown #include "power.h" 9c4762a1bSJed Brown #include <petscdmnetwork.h> 10c4762a1bSJed Brown 119371c9d4SSatish Balay PetscErrorCode FormFunction_Subnet(DM networkdm, Vec localX, Vec localF, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx) { 12c4762a1bSJed Brown UserCtx_Power *User = (UserCtx_Power *)appctx; 13c4762a1bSJed Brown PetscInt e, v, vfrom, vto; 14c4762a1bSJed Brown const PetscScalar *xarr; 15c4762a1bSJed Brown PetscScalar *farr; 16c4762a1bSJed Brown PetscInt offsetfrom, offsetto, offset; 17c4762a1bSJed Brown 18c4762a1bSJed Brown PetscFunctionBegin; 199566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localX, &xarr)); 209566063dSJacob Faibussowitsch PetscCall(VecGetArray(localF, &farr)); 21c4762a1bSJed Brown 22c4762a1bSJed Brown for (v = 0; v < nv; v++) { 23c4762a1bSJed Brown PetscInt i, j, key; 24c4762a1bSJed Brown PetscScalar Vm; 25c4762a1bSJed Brown PetscScalar Sbase = User->Sbase; 26c4762a1bSJed Brown VERTEX_Power bus = NULL; 27c4762a1bSJed Brown GEN gen; 28c4762a1bSJed Brown LOAD load; 29c4762a1bSJed Brown PetscBool ghostvtex; 30c4762a1bSJed Brown PetscInt numComps; 31c4762a1bSJed Brown void *component; 32c4762a1bSJed Brown 339566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[v], &ghostvtex)); 349566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm, vtx[v], &numComps)); 359566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[v], ALL_COMPONENTS, &offset)); 36c4762a1bSJed Brown for (j = 0; j < numComps; j++) { 379566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vtx[v], j, &key, &component, NULL)); 38c4762a1bSJed Brown if (key == 1) { 39c4762a1bSJed Brown PetscInt nconnedges; 40c4762a1bSJed Brown const PetscInt *connedges; 41c4762a1bSJed Brown 42c4762a1bSJed Brown bus = (VERTEX_Power)(component); 43c4762a1bSJed Brown /* Handle reference bus constrained dofs */ 44c4762a1bSJed Brown if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) { 45c4762a1bSJed Brown farr[offset] = xarr[offset] - bus->va * PETSC_PI / 180.0; 46c4762a1bSJed Brown farr[offset + 1] = xarr[offset + 1] - bus->vm; 47c4762a1bSJed Brown break; 48c4762a1bSJed Brown } 49c4762a1bSJed Brown 50c4762a1bSJed Brown if (!ghostvtex) { 51c4762a1bSJed Brown Vm = xarr[offset + 1]; 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* Shunt injections */ 54c4762a1bSJed Brown farr[offset] += Vm * Vm * bus->gl / Sbase; 55c4762a1bSJed Brown if (bus->ide != PV_BUS) farr[offset + 1] += -Vm * Vm * bus->bl / Sbase; 56c4762a1bSJed Brown } 57c4762a1bSJed Brown 589566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSupportingEdges(networkdm, vtx[v], &nconnedges, &connedges)); 59c4762a1bSJed Brown for (i = 0; i < nconnedges; i++) { 60c4762a1bSJed Brown EDGE_Power branch; 61c4762a1bSJed Brown PetscInt keye; 62c4762a1bSJed Brown PetscScalar Gff, Bff, Gft, Bft, Gtf, Btf, Gtt, Btt; 63c4762a1bSJed Brown const PetscInt *cone; 64c4762a1bSJed Brown PetscScalar Vmf, Vmt, thetaf, thetat, thetaft, thetatf; 65c4762a1bSJed Brown 66c4762a1bSJed Brown e = connedges[i]; 679566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, e, 0, &keye, (void **)&branch, NULL)); 68c4762a1bSJed Brown if (!branch->status) continue; 69c4762a1bSJed Brown Gff = branch->yff[0]; 70c4762a1bSJed Brown Bff = branch->yff[1]; 71c4762a1bSJed Brown Gft = branch->yft[0]; 72c4762a1bSJed Brown Bft = branch->yft[1]; 73c4762a1bSJed Brown Gtf = branch->ytf[0]; 74c4762a1bSJed Brown Btf = branch->ytf[1]; 75c4762a1bSJed Brown Gtt = branch->ytt[0]; 76c4762a1bSJed Brown Btt = branch->ytt[1]; 77c4762a1bSJed Brown 789566063dSJacob Faibussowitsch PetscCall(DMNetworkGetConnectedVertices(networkdm, e, &cone)); 79c4762a1bSJed Brown vfrom = cone[0]; 80c4762a1bSJed Brown vto = cone[1]; 81c4762a1bSJed Brown 829566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vfrom, ALL_COMPONENTS, &offsetfrom)); 839566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vto, ALL_COMPONENTS, &offsetto)); 84c4762a1bSJed Brown 85c4762a1bSJed Brown thetaf = xarr[offsetfrom]; 86c4762a1bSJed Brown Vmf = xarr[offsetfrom + 1]; 87c4762a1bSJed Brown thetat = xarr[offsetto]; 88c4762a1bSJed Brown Vmt = xarr[offsetto + 1]; 89c4762a1bSJed Brown thetaft = thetaf - thetat; 90c4762a1bSJed Brown thetatf = thetat - thetaf; 91c4762a1bSJed Brown 92c4762a1bSJed Brown if (vfrom == vtx[v]) { 93c4762a1bSJed Brown farr[offsetfrom] += Gff * Vmf * Vmf + Vmf * Vmt * (Gft * PetscCosScalar(thetaft) + Bft * PetscSinScalar(thetaft)); 94c4762a1bSJed Brown farr[offsetfrom + 1] += -Bff * Vmf * Vmf + Vmf * Vmt * (-Bft * PetscCosScalar(thetaft) + Gft * PetscSinScalar(thetaft)); 95c4762a1bSJed Brown } else { 96c4762a1bSJed Brown farr[offsetto] += Gtt * Vmt * Vmt + Vmt * Vmf * (Gtf * PetscCosScalar(thetatf) + Btf * PetscSinScalar(thetatf)); 97c4762a1bSJed Brown farr[offsetto + 1] += -Btt * Vmt * Vmt + Vmt * Vmf * (-Btf * PetscCosScalar(thetatf) + Gtf * PetscSinScalar(thetatf)); 98c4762a1bSJed Brown } 99c4762a1bSJed Brown } 100c4762a1bSJed Brown } else if (key == 2) { 101c4762a1bSJed Brown if (!ghostvtex) { 102c4762a1bSJed Brown gen = (GEN)(component); 103c4762a1bSJed Brown if (!gen->status) continue; 104c4762a1bSJed Brown farr[offset] += -gen->pg / Sbase; 105c4762a1bSJed Brown farr[offset + 1] += -gen->qg / Sbase; 106c4762a1bSJed Brown } 107c4762a1bSJed Brown } else if (key == 3) { 108c4762a1bSJed Brown if (!ghostvtex) { 109c4762a1bSJed Brown load = (LOAD)(component); 110c4762a1bSJed Brown farr[offset] += load->pl / Sbase; 111c4762a1bSJed Brown farr[offset + 1] += load->ql / Sbase; 112c4762a1bSJed Brown } 113c4762a1bSJed Brown } 114c4762a1bSJed Brown } 1159371c9d4SSatish Balay if (bus && bus->ide == PV_BUS) { farr[offset + 1] = xarr[offset + 1] - bus->vm; } 116c4762a1bSJed Brown } 1179566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(localX, &xarr)); 1189566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localF, &farr)); 119c4762a1bSJed Brown PetscFunctionReturn(0); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 1229371c9d4SSatish Balay PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *appctx) { 123c4762a1bSJed Brown DM networkdm; 124c4762a1bSJed Brown Vec localX, localF; 125c4762a1bSJed Brown PetscInt nv, ne; 126c4762a1bSJed Brown const PetscInt *vtx, *edges; 127c4762a1bSJed Brown 128c4762a1bSJed Brown PetscFunctionBegin; 1299566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &networkdm)); 1309566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm, &localX)); 1319566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm, &localF)); 1329566063dSJacob Faibussowitsch PetscCall(VecSet(F, 0.0)); 133c4762a1bSJed Brown 1349566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX)); 1359566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX)); 136c4762a1bSJed Brown 1379566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm, F, INSERT_VALUES, localF)); 1389566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm, F, INSERT_VALUES, localF)); 139c4762a1bSJed Brown 140c4762a1bSJed Brown /* Form Function for first subnetwork */ 1419566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges)); 1429566063dSJacob Faibussowitsch PetscCall(FormFunction_Subnet(networkdm, localX, localF, nv, ne, vtx, edges, appctx)); 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* Form Function for second subnetwork */ 1459566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 1, &nv, &ne, &vtx, &edges)); 1469566063dSJacob Faibussowitsch PetscCall(FormFunction_Subnet(networkdm, localX, localF, nv, ne, vtx, edges, appctx)); 147c4762a1bSJed Brown 1489566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm, &localX)); 149c4762a1bSJed Brown 1509566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(networkdm, localF, ADD_VALUES, F)); 1519566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(networkdm, localF, ADD_VALUES, F)); 1529566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm, &localF)); 153c4762a1bSJed Brown PetscFunctionReturn(0); 154c4762a1bSJed Brown } 155c4762a1bSJed Brown 1569371c9d4SSatish Balay PetscErrorCode FormJacobian_Subnet(DM networkdm, Vec localX, Mat J, Mat Jpre, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx) { 157c4762a1bSJed Brown UserCtx_Power *User = (UserCtx_Power *)appctx; 158c4762a1bSJed Brown PetscInt e, v, vfrom, vto; 159c4762a1bSJed Brown const PetscScalar *xarr; 160c4762a1bSJed Brown PetscInt offsetfrom, offsetto, goffsetfrom, goffsetto; 161c4762a1bSJed Brown PetscInt row[2], col[8]; 162c4762a1bSJed Brown PetscScalar values[8]; 163c4762a1bSJed Brown 164c4762a1bSJed Brown PetscFunctionBegin; 1659566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localX, &xarr)); 166c4762a1bSJed Brown 167c4762a1bSJed Brown for (v = 0; v < nv; v++) { 168c4762a1bSJed Brown PetscInt i, j, key; 169c4762a1bSJed Brown PetscInt offset, goffset; 170c4762a1bSJed Brown PetscScalar Vm; 171c4762a1bSJed Brown PetscScalar Sbase = User->Sbase; 172c4762a1bSJed Brown VERTEX_Power bus; 173c4762a1bSJed Brown PetscBool ghostvtex; 174c4762a1bSJed Brown PetscInt numComps; 175c4762a1bSJed Brown void *component; 176c4762a1bSJed Brown 1779566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[v], &ghostvtex)); 1789566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm, vtx[v], &numComps)); 179c4762a1bSJed Brown for (j = 0; j < numComps; j++) { 1809566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[v], ALL_COMPONENTS, &offset)); 1819566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVecOffset(networkdm, vtx[v], ALL_COMPONENTS, &goffset)); 1829566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vtx[v], j, &key, &component, NULL)); 183c4762a1bSJed Brown if (key == 1) { 184c4762a1bSJed Brown PetscInt nconnedges; 185c4762a1bSJed Brown const PetscInt *connedges; 186c4762a1bSJed Brown 187c4762a1bSJed Brown bus = (VERTEX_Power)(component); 188c4762a1bSJed Brown if (!ghostvtex) { 189c4762a1bSJed Brown /* Handle reference bus constrained dofs */ 190c4762a1bSJed Brown if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) { 1919371c9d4SSatish Balay row[0] = goffset; 1929371c9d4SSatish Balay row[1] = goffset + 1; 1939371c9d4SSatish Balay col[0] = goffset; 1949371c9d4SSatish Balay col[1] = goffset + 1; 1959371c9d4SSatish Balay col[2] = goffset; 1969371c9d4SSatish Balay col[3] = goffset + 1; 1979371c9d4SSatish Balay values[0] = 1.0; 1989371c9d4SSatish Balay values[1] = 0.0; 1999371c9d4SSatish Balay values[2] = 0.0; 2009371c9d4SSatish Balay values[3] = 1.0; 2019566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 2, row, 2, col, values, ADD_VALUES)); 202c4762a1bSJed Brown break; 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205c4762a1bSJed Brown Vm = xarr[offset + 1]; 206c4762a1bSJed Brown 207c4762a1bSJed Brown /* Shunt injections */ 2089371c9d4SSatish Balay row[0] = goffset; 2099371c9d4SSatish Balay row[1] = goffset + 1; 2109371c9d4SSatish Balay col[0] = goffset; 2119371c9d4SSatish Balay col[1] = goffset + 1; 212c4762a1bSJed Brown values[0] = values[1] = values[2] = values[3] = 0.0; 213c4762a1bSJed Brown if (bus->ide != PV_BUS) { 214c4762a1bSJed Brown values[1] = 2.0 * Vm * bus->gl / Sbase; 215c4762a1bSJed Brown values[3] = -2.0 * Vm * bus->bl / Sbase; 216c4762a1bSJed Brown } 2179566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 2, row, 2, col, values, ADD_VALUES)); 218c4762a1bSJed Brown } 219c4762a1bSJed Brown 2209566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSupportingEdges(networkdm, vtx[v], &nconnedges, &connedges)); 221c4762a1bSJed Brown for (i = 0; i < nconnedges; i++) { 222c4762a1bSJed Brown EDGE_Power branch; 223c4762a1bSJed Brown VERTEX_Power busf, bust; 224c4762a1bSJed Brown PetscInt keyf, keyt; 225c4762a1bSJed Brown PetscScalar Gff, Bff, Gft, Bft, Gtf, Btf, Gtt, Btt; 226c4762a1bSJed Brown const PetscInt *cone; 227c4762a1bSJed Brown PetscScalar Vmf, Vmt, thetaf, thetat, thetaft, thetatf; 228c4762a1bSJed Brown 229c4762a1bSJed Brown e = connedges[i]; 2309566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, e, 0, &key, (void **)&branch, NULL)); 231c4762a1bSJed Brown if (!branch->status) continue; 232c4762a1bSJed Brown 233c4762a1bSJed Brown Gff = branch->yff[0]; 234c4762a1bSJed Brown Bff = branch->yff[1]; 235c4762a1bSJed Brown Gft = branch->yft[0]; 236c4762a1bSJed Brown Bft = branch->yft[1]; 237c4762a1bSJed Brown Gtf = branch->ytf[0]; 238c4762a1bSJed Brown Btf = branch->ytf[1]; 239c4762a1bSJed Brown Gtt = branch->ytt[0]; 240c4762a1bSJed Brown Btt = branch->ytt[1]; 241c4762a1bSJed Brown 2429566063dSJacob Faibussowitsch PetscCall(DMNetworkGetConnectedVertices(networkdm, e, &cone)); 243c4762a1bSJed Brown vfrom = cone[0]; 244c4762a1bSJed Brown vto = cone[1]; 245c4762a1bSJed Brown 2469566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vfrom, ALL_COMPONENTS, &offsetfrom)); 2479566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vto, ALL_COMPONENTS, &offsetto)); 2489566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVecOffset(networkdm, vfrom, ALL_COMPONENTS, &goffsetfrom)); 2499566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVecOffset(networkdm, vto, ALL_COMPONENTS, &goffsetto)); 250c4762a1bSJed Brown 251c4762a1bSJed Brown if (goffsetto < 0) goffsetto = -goffsetto - 1; 252c4762a1bSJed Brown 253c4762a1bSJed Brown thetaf = xarr[offsetfrom]; 254c4762a1bSJed Brown Vmf = xarr[offsetfrom + 1]; 255c4762a1bSJed Brown thetat = xarr[offsetto]; 256c4762a1bSJed Brown Vmt = xarr[offsetto + 1]; 257c4762a1bSJed Brown thetaft = thetaf - thetat; 258c4762a1bSJed Brown thetatf = thetat - thetaf; 259c4762a1bSJed Brown 2609566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vfrom, 0, &keyf, (void **)&busf, NULL)); 2619566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vto, 0, &keyt, (void **)&bust, NULL)); 262c4762a1bSJed Brown 263c4762a1bSJed Brown if (vfrom == vtx[v]) { 264c4762a1bSJed Brown if (busf->ide != REF_BUS) { 265c4762a1bSJed Brown /* farr[offsetfrom] += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); */ 266c4762a1bSJed Brown row[0] = goffsetfrom; 2679371c9d4SSatish Balay col[0] = goffsetfrom; 2689371c9d4SSatish Balay col[1] = goffsetfrom + 1; 2699371c9d4SSatish Balay col[2] = goffsetto; 2709371c9d4SSatish Balay col[3] = goffsetto + 1; 271c4762a1bSJed Brown values[0] = Vmf * Vmt * (Gft * -PetscSinScalar(thetaft) + Bft * PetscCosScalar(thetaft)); /* df_dthetaf */ 272c4762a1bSJed Brown values[1] = 2.0 * Gff * Vmf + Vmt * (Gft * PetscCosScalar(thetaft) + Bft * PetscSinScalar(thetaft)); /* df_dVmf */ 273c4762a1bSJed Brown values[2] = Vmf * Vmt * (Gft * PetscSinScalar(thetaft) + Bft * -PetscCosScalar(thetaft)); /* df_dthetat */ 274c4762a1bSJed Brown values[3] = Vmf * (Gft * PetscCosScalar(thetaft) + Bft * PetscSinScalar(thetaft)); /* df_dVmt */ 275c4762a1bSJed Brown 2769566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 4, col, values, ADD_VALUES)); 277c4762a1bSJed Brown } 278c4762a1bSJed Brown if (busf->ide != PV_BUS && busf->ide != REF_BUS) { 279c4762a1bSJed Brown row[0] = goffsetfrom + 1; 2809371c9d4SSatish Balay col[0] = goffsetfrom; 2819371c9d4SSatish Balay col[1] = goffsetfrom + 1; 2829371c9d4SSatish Balay col[2] = goffsetto; 2839371c9d4SSatish Balay col[3] = goffsetto + 1; 284c4762a1bSJed Brown /* farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); */ 285c4762a1bSJed Brown values[0] = Vmf * Vmt * (Bft * PetscSinScalar(thetaft) + Gft * PetscCosScalar(thetaft)); 286c4762a1bSJed Brown values[1] = -2.0 * Bff * Vmf + Vmt * (-Bft * PetscCosScalar(thetaft) + Gft * PetscSinScalar(thetaft)); 287c4762a1bSJed Brown values[2] = Vmf * Vmt * (-Bft * PetscSinScalar(thetaft) + Gft * -PetscCosScalar(thetaft)); 288c4762a1bSJed Brown values[3] = Vmf * (-Bft * PetscCosScalar(thetaft) + Gft * PetscSinScalar(thetaft)); 289c4762a1bSJed Brown 2909566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 4, col, values, ADD_VALUES)); 291c4762a1bSJed Brown } 292c4762a1bSJed Brown } else { 293c4762a1bSJed Brown if (bust->ide != REF_BUS) { 294c4762a1bSJed Brown row[0] = goffsetto; 2959371c9d4SSatish Balay col[0] = goffsetto; 2969371c9d4SSatish Balay col[1] = goffsetto + 1; 2979371c9d4SSatish Balay col[2] = goffsetfrom; 2989371c9d4SSatish Balay col[3] = goffsetfrom + 1; 299c4762a1bSJed Brown /* farr[offsetto] += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); */ 300c4762a1bSJed Brown values[0] = Vmt * Vmf * (Gtf * -PetscSinScalar(thetatf) + Btf * PetscCosScalar(thetaft)); /* df_dthetat */ 301c4762a1bSJed Brown values[1] = 2.0 * Gtt * Vmt + Vmf * (Gtf * PetscCosScalar(thetatf) + Btf * PetscSinScalar(thetatf)); /* df_dVmt */ 302c4762a1bSJed Brown values[2] = Vmt * Vmf * (Gtf * PetscSinScalar(thetatf) + Btf * -PetscCosScalar(thetatf)); /* df_dthetaf */ 303c4762a1bSJed Brown values[3] = Vmt * (Gtf * PetscCosScalar(thetatf) + Btf * PetscSinScalar(thetatf)); /* df_dVmf */ 304c4762a1bSJed Brown 3059566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 4, col, values, ADD_VALUES)); 306c4762a1bSJed Brown } 307c4762a1bSJed Brown if (bust->ide != PV_BUS && bust->ide != REF_BUS) { 308c4762a1bSJed Brown row[0] = goffsetto + 1; 3099371c9d4SSatish Balay col[0] = goffsetto; 3109371c9d4SSatish Balay col[1] = goffsetto + 1; 3119371c9d4SSatish Balay col[2] = goffsetfrom; 3129371c9d4SSatish Balay col[3] = goffsetfrom + 1; 313c4762a1bSJed Brown /* farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); */ 314c4762a1bSJed Brown values[0] = Vmt * Vmf * (Btf * PetscSinScalar(thetatf) + Gtf * PetscCosScalar(thetatf)); 315c4762a1bSJed Brown values[1] = -2.0 * Btt * Vmt + Vmf * (-Btf * PetscCosScalar(thetatf) + Gtf * PetscSinScalar(thetatf)); 316c4762a1bSJed Brown values[2] = Vmt * Vmf * (-Btf * PetscSinScalar(thetatf) + Gtf * -PetscCosScalar(thetatf)); 317c4762a1bSJed Brown values[3] = Vmt * (-Btf * PetscCosScalar(thetatf) + Gtf * PetscSinScalar(thetatf)); 318c4762a1bSJed Brown 3199566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 4, col, values, ADD_VALUES)); 320c4762a1bSJed Brown } 321c4762a1bSJed Brown } 322c4762a1bSJed Brown } 323c4762a1bSJed Brown if (!ghostvtex && bus->ide == PV_BUS) { 3249371c9d4SSatish Balay row[0] = goffset + 1; 3259371c9d4SSatish Balay col[0] = goffset + 1; 326c4762a1bSJed Brown values[0] = 1.0; 3279566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 1, col, values, ADD_VALUES)); 328c4762a1bSJed Brown } 329c4762a1bSJed Brown } 330c4762a1bSJed Brown } 331c4762a1bSJed Brown } 3329566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(localX, &xarr)); 333c4762a1bSJed Brown PetscFunctionReturn(0); 334c4762a1bSJed Brown } 335c4762a1bSJed Brown 3369371c9d4SSatish Balay PetscErrorCode FormJacobian(SNES snes, Vec X, Mat J, Mat Jpre, void *appctx) { 337c4762a1bSJed Brown DM networkdm; 338c4762a1bSJed Brown Vec localX; 339c4762a1bSJed Brown PetscInt ne, nv; 340c4762a1bSJed Brown const PetscInt *vtx, *edges; 341c4762a1bSJed Brown 342c4762a1bSJed Brown PetscFunctionBegin; 3439566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(J)); 344c4762a1bSJed Brown 3459566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &networkdm)); 3469566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm, &localX)); 347c4762a1bSJed Brown 3489566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX)); 3499566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX)); 350c4762a1bSJed Brown 351c4762a1bSJed Brown /* Form Jacobian for first subnetwork */ 3529566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges)); 3539566063dSJacob Faibussowitsch PetscCall(FormJacobian_Subnet(networkdm, localX, J, Jpre, nv, ne, vtx, edges, appctx)); 354c4762a1bSJed Brown 355c4762a1bSJed Brown /* Form Jacobian for second subnetwork */ 3569566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 1, &nv, &ne, &vtx, &edges)); 3579566063dSJacob Faibussowitsch PetscCall(FormJacobian_Subnet(networkdm, localX, J, Jpre, nv, ne, vtx, edges, appctx)); 358c4762a1bSJed Brown 3599566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm, &localX)); 360c4762a1bSJed Brown 3619566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 3629566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 363c4762a1bSJed Brown PetscFunctionReturn(0); 364c4762a1bSJed Brown } 365c4762a1bSJed Brown 3669371c9d4SSatish Balay PetscErrorCode SetInitialValues_Subnet(DM networkdm, Vec localX, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx) { 367c4762a1bSJed Brown VERTEX_Power bus; 368c4762a1bSJed Brown PetscInt i; 369c4762a1bSJed Brown GEN gen; 370c4762a1bSJed Brown PetscBool ghostvtex; 371c4762a1bSJed Brown PetscScalar *xarr; 372c4762a1bSJed Brown PetscInt key, numComps, j, offset; 373c4762a1bSJed Brown void *component; 374c4762a1bSJed Brown 375c4762a1bSJed Brown PetscFunctionBegin; 3769566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX, &xarr)); 377c4762a1bSJed Brown for (i = 0; i < nv; i++) { 3789566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[i], &ghostvtex)); 379c4762a1bSJed Brown if (ghostvtex) continue; 380c4762a1bSJed Brown 3819566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[i], ALL_COMPONENTS, &offset)); 3829566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm, vtx[i], &numComps)); 383c4762a1bSJed Brown for (j = 0; j < numComps; j++) { 3849566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vtx[i], j, &key, &component, NULL)); 385c4762a1bSJed Brown if (key == 1) { 386c4762a1bSJed Brown bus = (VERTEX_Power)(component); 387c4762a1bSJed Brown xarr[offset] = bus->va * PETSC_PI / 180.0; 388c4762a1bSJed Brown xarr[offset + 1] = bus->vm; 389c4762a1bSJed Brown } else if (key == 2) { 390c4762a1bSJed Brown gen = (GEN)(component); 391c4762a1bSJed Brown if (!gen->status) continue; 392c4762a1bSJed Brown xarr[offset + 1] = gen->vs; 393c4762a1bSJed Brown break; 394c4762a1bSJed Brown } 395c4762a1bSJed Brown } 396c4762a1bSJed Brown } 3979566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX, &xarr)); 398c4762a1bSJed Brown PetscFunctionReturn(0); 399c4762a1bSJed Brown } 400c4762a1bSJed Brown 4019371c9d4SSatish Balay PetscErrorCode SetInitialValues(DM networkdm, Vec X, void *appctx) { 402c4762a1bSJed Brown PetscInt nv, ne; 403c4762a1bSJed Brown const PetscInt *vtx, *edges; 404c4762a1bSJed Brown Vec localX; 405c4762a1bSJed Brown 406c4762a1bSJed Brown PetscFunctionBegin; 4079566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm, &localX)); 408c4762a1bSJed Brown 4099566063dSJacob Faibussowitsch PetscCall(VecSet(X, 0.0)); 4109566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX)); 4119566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX)); 412c4762a1bSJed Brown 413c4762a1bSJed Brown /* Set initial guess for first subnetwork */ 4149566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges)); 4159566063dSJacob Faibussowitsch PetscCall(SetInitialValues_Subnet(networkdm, localX, nv, ne, vtx, edges, appctx)); 416c4762a1bSJed Brown 417c4762a1bSJed Brown /* Set initial guess for second subnetwork */ 4189566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 1, &nv, &ne, &vtx, &edges)); 4199566063dSJacob Faibussowitsch PetscCall(SetInitialValues_Subnet(networkdm, localX, nv, ne, vtx, edges, appctx)); 420c4762a1bSJed Brown 4219566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(networkdm, localX, ADD_VALUES, X)); 4229566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(networkdm, localX, ADD_VALUES, X)); 4239566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm, &localX)); 424c4762a1bSJed Brown PetscFunctionReturn(0); 425c4762a1bSJed Brown } 426c4762a1bSJed Brown 4279371c9d4SSatish Balay int main(int argc, char **argv) { 428c4762a1bSJed Brown char pfdata_file[PETSC_MAX_PATH_LEN] = "case9.m"; 429c4762a1bSJed Brown PFDATA *pfdata1, *pfdata2; 430f11a936eSBarry Smith PetscInt numEdges1 = 0, numEdges2 = 0; 4312bf73ac6SHong Zhang PetscInt *edgelist1 = NULL, *edgelist2 = NULL, componentkey[4]; 432c4762a1bSJed Brown DM networkdm; 433c4762a1bSJed Brown UserCtx_Power User; 434956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 435c4762a1bSJed Brown PetscLogStage stage1, stage2; 436956f8c0dSBarry Smith #endif 437c4762a1bSJed Brown PetscMPIInt rank; 4382bf73ac6SHong Zhang PetscInt nsubnet = 2, nv, ne, i, j, genj, loadj; 4392bf73ac6SHong Zhang const PetscInt *vtx, *edges; 440c4762a1bSJed Brown Vec X, F; 441c4762a1bSJed Brown Mat J; 442c4762a1bSJed Brown SNES snes; 443c4762a1bSJed Brown 444327415f7SBarry Smith PetscFunctionBeginUser; 4459566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, "poweroptions", help)); 4469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank)); 447c4762a1bSJed Brown { 448c4762a1bSJed Brown /* introduce the const crank so the clang static analyzer realizes that if it enters any of the if (crank) then it must have entered the first */ 449c4762a1bSJed Brown /* this is an experiment to see how the analyzer reacts */ 450c4762a1bSJed Brown const PetscMPIInt crank = rank; 451c4762a1bSJed Brown 452c4762a1bSJed Brown /* Create an empty network object */ 4539566063dSJacob Faibussowitsch PetscCall(DMNetworkCreate(PETSC_COMM_WORLD, &networkdm)); 454c4762a1bSJed Brown 455c4762a1bSJed Brown /* Register the components in the network */ 4569566063dSJacob Faibussowitsch PetscCall(DMNetworkRegisterComponent(networkdm, "branchstruct", sizeof(struct _p_EDGE_Power), &componentkey[0])); 4579566063dSJacob Faibussowitsch PetscCall(DMNetworkRegisterComponent(networkdm, "busstruct", sizeof(struct _p_VERTEX_Power), &componentkey[1])); 4589566063dSJacob Faibussowitsch PetscCall(DMNetworkRegisterComponent(networkdm, "genstruct", sizeof(struct _p_GEN), &componentkey[2])); 4599566063dSJacob Faibussowitsch PetscCall(DMNetworkRegisterComponent(networkdm, "loadstruct", sizeof(struct _p_LOAD), &componentkey[3])); 460c4762a1bSJed Brown 4619566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Read Data", &stage1)); 462c4762a1bSJed Brown PetscLogStagePush(stage1); 463c4762a1bSJed Brown /* READ THE DATA */ 464c4762a1bSJed Brown if (!crank) { 465c4762a1bSJed Brown /* Only rank 0 reads the data */ 4669566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-pfdata", pfdata_file, sizeof(pfdata_file), NULL)); 467c4762a1bSJed Brown /* HERE WE CREATE COPIES OF THE SAME NETWORK THAT WILL BE TREATED AS SUBNETWORKS */ 468c4762a1bSJed Brown 469c4762a1bSJed Brown /* READ DATA FOR THE FIRST SUBNETWORK */ 4709566063dSJacob Faibussowitsch PetscCall(PetscNew(&pfdata1)); 4719566063dSJacob Faibussowitsch PetscCall(PFReadMatPowerData(pfdata1, pfdata_file)); 472c4762a1bSJed Brown User.Sbase = pfdata1->sbase; 473c4762a1bSJed Brown 474c4762a1bSJed Brown numEdges1 = pfdata1->nbranch; 4759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * numEdges1, &edgelist1)); 4769566063dSJacob Faibussowitsch PetscCall(GetListofEdges_Power(pfdata1, edgelist1)); 477c4762a1bSJed Brown 478c4762a1bSJed Brown /* READ DATA FOR THE SECOND SUBNETWORK */ 4799566063dSJacob Faibussowitsch PetscCall(PetscNew(&pfdata2)); 4809566063dSJacob Faibussowitsch PetscCall(PFReadMatPowerData(pfdata2, pfdata_file)); 481c4762a1bSJed Brown User.Sbase = pfdata2->sbase; 482c4762a1bSJed Brown 483c4762a1bSJed Brown numEdges2 = pfdata2->nbranch; 4849566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * numEdges2, &edgelist2)); 4859566063dSJacob Faibussowitsch PetscCall(GetListofEdges_Power(pfdata2, edgelist2)); 486c4762a1bSJed Brown } 487c4762a1bSJed Brown 488c4762a1bSJed Brown PetscLogStagePop(); 4899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 4909566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Create network", &stage2)); 491c4762a1bSJed Brown PetscLogStagePush(stage2); 492c4762a1bSJed Brown 4932bf73ac6SHong Zhang /* Set number of nodes/edges and edge connectivity */ 4949566063dSJacob Faibussowitsch PetscCall(DMNetworkSetNumSubNetworks(networkdm, PETSC_DECIDE, nsubnet)); 4959566063dSJacob Faibussowitsch PetscCall(DMNetworkAddSubnetwork(networkdm, "", numEdges1, edgelist1, NULL)); 4969566063dSJacob Faibussowitsch PetscCall(DMNetworkAddSubnetwork(networkdm, "", numEdges2, edgelist2, NULL)); 497c4762a1bSJed Brown 498c4762a1bSJed Brown /* Set up the network layout */ 4999566063dSJacob Faibussowitsch PetscCall(DMNetworkLayoutSetUp(networkdm)); 500c4762a1bSJed Brown 501c4762a1bSJed Brown /* Add network components only process 0 has any data to add*/ 502c4762a1bSJed Brown if (!crank) { 5039371c9d4SSatish Balay genj = 0; 5049371c9d4SSatish Balay loadj = 0; 505c4762a1bSJed Brown 506c4762a1bSJed Brown /* ADD VARIABLES AND COMPONENTS FOR THE FIRST SUBNETWORK */ 5079566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges)); 508c4762a1bSJed Brown 509*48a46eb9SPierre Jolivet for (i = 0; i < ne; i++) PetscCall(DMNetworkAddComponent(networkdm, edges[i], componentkey[0], &pfdata1->branch[i], 0)); 510c4762a1bSJed Brown 511c4762a1bSJed Brown for (i = 0; i < nv; i++) { 5129566063dSJacob Faibussowitsch PetscCall(DMNetworkAddComponent(networkdm, vtx[i], componentkey[1], &pfdata1->bus[i], 2)); 513c4762a1bSJed Brown if (pfdata1->bus[i].ngen) { 514*48a46eb9SPierre Jolivet for (j = 0; j < pfdata1->bus[i].ngen; j++) PetscCall(DMNetworkAddComponent(networkdm, vtx[i], componentkey[2], &pfdata1->gen[genj++], 0)); 515c4762a1bSJed Brown } 516c4762a1bSJed Brown if (pfdata1->bus[i].nload) { 517*48a46eb9SPierre Jolivet for (j = 0; j < pfdata1->bus[i].nload; j++) PetscCall(DMNetworkAddComponent(networkdm, vtx[i], componentkey[3], &pfdata1->load[loadj++], 0)); 518c4762a1bSJed Brown } 519c4762a1bSJed Brown } 520c4762a1bSJed Brown 5219371c9d4SSatish Balay genj = 0; 5229371c9d4SSatish Balay loadj = 0; 523c4762a1bSJed Brown 524c4762a1bSJed Brown /* ADD VARIABLES AND COMPONENTS FOR THE SECOND SUBNETWORK */ 5259566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 1, &nv, &ne, &vtx, &edges)); 526c4762a1bSJed Brown 527*48a46eb9SPierre Jolivet for (i = 0; i < ne; i++) PetscCall(DMNetworkAddComponent(networkdm, edges[i], componentkey[0], &pfdata2->branch[i], 0)); 528c4762a1bSJed Brown 529c4762a1bSJed Brown for (i = 0; i < nv; i++) { 5309566063dSJacob Faibussowitsch PetscCall(DMNetworkAddComponent(networkdm, vtx[i], componentkey[1], &pfdata2->bus[i], 2)); 531c4762a1bSJed Brown if (pfdata2->bus[i].ngen) { 532*48a46eb9SPierre Jolivet for (j = 0; j < pfdata2->bus[i].ngen; j++) PetscCall(DMNetworkAddComponent(networkdm, vtx[i], componentkey[2], &pfdata2->gen[genj++], 0)); 533c4762a1bSJed Brown } 534c4762a1bSJed Brown if (pfdata2->bus[i].nload) { 535*48a46eb9SPierre Jolivet for (j = 0; j < pfdata2->bus[i].nload; j++) PetscCall(DMNetworkAddComponent(networkdm, vtx[i], componentkey[3], &pfdata2->load[loadj++], 0)); 536c4762a1bSJed Brown } 537c4762a1bSJed Brown } 538c4762a1bSJed Brown } 539c4762a1bSJed Brown 540c4762a1bSJed Brown /* Set up DM for use */ 5419566063dSJacob Faibussowitsch PetscCall(DMSetUp(networkdm)); 542c4762a1bSJed Brown 543c4762a1bSJed Brown if (!crank) { 5449566063dSJacob Faibussowitsch PetscCall(PetscFree(edgelist1)); 5459566063dSJacob Faibussowitsch PetscCall(PetscFree(edgelist2)); 546c4762a1bSJed Brown } 547c4762a1bSJed Brown 548c4762a1bSJed Brown if (!crank) { 5499566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata1->bus)); 5509566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata1->gen)); 5519566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata1->branch)); 5529566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata1->load)); 5539566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata1)); 554c4762a1bSJed Brown 5559566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata2->bus)); 5569566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata2->gen)); 5579566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata2->branch)); 5589566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata2->load)); 5599566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata2)); 560c4762a1bSJed Brown } 561c4762a1bSJed Brown 562c4762a1bSJed Brown /* Distribute networkdm to multiple processes */ 5639566063dSJacob Faibussowitsch PetscCall(DMNetworkDistribute(&networkdm, 0)); 564c4762a1bSJed Brown 565c4762a1bSJed Brown PetscLogStagePop(); 566c4762a1bSJed Brown 567c4762a1bSJed Brown /* Broadcast Sbase to all processors */ 5689566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&User.Sbase, 1, MPIU_SCALAR, 0, PETSC_COMM_WORLD)); 569c4762a1bSJed Brown 5709566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(networkdm, &X)); 5719566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &F)); 572c4762a1bSJed Brown 5739566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(networkdm, &J)); 5749566063dSJacob Faibussowitsch PetscCall(MatSetOption(J, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE)); 575c4762a1bSJed Brown 5769566063dSJacob Faibussowitsch PetscCall(SetInitialValues(networkdm, X, &User)); 577c4762a1bSJed Brown 578c4762a1bSJed Brown /* HOOK UP SOLVER */ 5799566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes)); 5809566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, networkdm)); 5819566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, F, FormFunction, &User)); 5829566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, &User)); 5839566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 584c4762a1bSJed Brown 5859566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, X)); 5869566063dSJacob Faibussowitsch /* PetscCall(VecView(X,PETSC_VIEWER_STDOUT_WORLD)); */ 587c4762a1bSJed Brown 5889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 5899566063dSJacob Faibussowitsch PetscCall(VecDestroy(&F)); 5909566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 591c4762a1bSJed Brown 5929566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 5939566063dSJacob Faibussowitsch PetscCall(DMDestroy(&networkdm)); 594c4762a1bSJed Brown } 5959566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 596b122ec5aSJacob Faibussowitsch return 0; 597c4762a1bSJed Brown } 598c4762a1bSJed Brown 599c4762a1bSJed Brown /*TEST 600c4762a1bSJed Brown 601c4762a1bSJed Brown build: 602c4762a1bSJed Brown depends: PFReadData.c pffunctions.c 603dfd57a17SPierre Jolivet requires: !complex double defined(PETSC_HAVE_ATTRIBUTEALIGNED) 604c4762a1bSJed Brown 605c4762a1bSJed Brown test: 606c4762a1bSJed Brown args: -snes_rtol 1.e-3 607c4762a1bSJed Brown localrunfiles: poweroptions case9.m 6082bf73ac6SHong Zhang output_file: output/power_1.out 609c4762a1bSJed Brown 610c4762a1bSJed Brown test: 611c4762a1bSJed Brown suffix: 2 612c4762a1bSJed Brown args: -snes_rtol 1.e-3 -petscpartitioner_type simple 613c4762a1bSJed Brown nsize: 4 614c4762a1bSJed Brown localrunfiles: poweroptions case9.m 6152bf73ac6SHong Zhang output_file: output/power_1.out 616c4762a1bSJed Brown 617c4762a1bSJed Brown TEST*/ 618