1c4762a1bSJed Brown /* function subroutines used by power.c */ 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include "power.h" 4c4762a1bSJed Brown 59371c9d4SSatish Balay PetscErrorCode GetListofEdges_Power(PFDATA *pfdata, PetscInt *edgelist) { 6c4762a1bSJed Brown PetscInt i, fbus, tbus, nbranches = pfdata->nbranch; 7c4762a1bSJed Brown EDGE_Power branch = pfdata->branch; 8c4762a1bSJed Brown PetscBool netview = PETSC_FALSE; 9c4762a1bSJed Brown 10c4762a1bSJed Brown PetscFunctionBegin; 119566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-powernet_view", &netview)); 12c4762a1bSJed Brown for (i = 0; i < nbranches; i++) { 13c4762a1bSJed Brown fbus = branch[i].internal_i; 14c4762a1bSJed Brown tbus = branch[i].internal_j; 15c4762a1bSJed Brown edgelist[2 * i] = fbus; 16c4762a1bSJed Brown edgelist[2 * i + 1] = tbus; 17*48a46eb9SPierre Jolivet if (netview) PetscCall(PetscPrintf(PETSC_COMM_SELF, "branch %" PetscInt_FMT ", bus[%" PetscInt_FMT "] -> bus[%" PetscInt_FMT "]\n", i, fbus, tbus)); 18c4762a1bSJed Brown } 19c4762a1bSJed Brown if (netview) { 20c4762a1bSJed Brown for (i = 0; i < pfdata->nbus; i++) { 21c4762a1bSJed Brown if (pfdata->bus[i].ngen) { 2263a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " bus %" PetscInt_FMT ": gen\n", i)); 23c4762a1bSJed Brown } else if (pfdata->bus[i].nload) { 2463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " bus %" PetscInt_FMT ": load\n", i)); 25c4762a1bSJed Brown } 26c4762a1bSJed Brown } 27c4762a1bSJed Brown } 28c4762a1bSJed Brown PetscFunctionReturn(0); 29c4762a1bSJed Brown } 30c4762a1bSJed Brown 319371c9d4SSatish Balay PetscErrorCode FormJacobian_Power_private(DM networkdm, Vec localX, Mat J, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx) { 32c4762a1bSJed Brown const PetscScalar *xarr; 33c4762a1bSJed Brown PetscInt i, v, row[2], col[8], e, vfrom, vto; 34c4762a1bSJed Brown PetscInt offsetfrom, offsetto, goffsetfrom, goffsetto, numComps; 35c4762a1bSJed Brown PetscScalar values[8]; 36c4762a1bSJed Brown PetscInt j, key, offset, goffset; 37c4762a1bSJed Brown PetscScalar Vm; 38c4762a1bSJed Brown UserCtx_Power *user_power = (UserCtx_Power *)appctx; 39c4762a1bSJed Brown PetscScalar Sbase = user_power->Sbase; 40c4762a1bSJed Brown VERTEX_Power bus; 41c4762a1bSJed Brown PetscBool ghostvtex; 42c4762a1bSJed Brown void *component; 43c4762a1bSJed Brown 44c4762a1bSJed Brown PetscFunctionBegin; 459566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localX, &xarr)); 46c4762a1bSJed Brown 47c4762a1bSJed Brown for (v = 0; v < nv; v++) { 489566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[v], &ghostvtex)); 49c4762a1bSJed Brown 509566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm, vtx[v], &numComps)); 51c4762a1bSJed Brown for (j = 0; j < numComps; j++) { 529566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[v], ALL_COMPONENTS, &offset)); 539566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVecOffset(networkdm, vtx[v], ALL_COMPONENTS, &goffset)); 549566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vtx[v], j, &key, &component, NULL)); 55c4762a1bSJed Brown 56c4762a1bSJed Brown if (key == user_power->compkey_bus) { 57c4762a1bSJed Brown PetscInt nconnedges; 58c4762a1bSJed Brown const PetscInt *connedges; 59c4762a1bSJed Brown 60c4762a1bSJed Brown bus = (VERTEX_Power)(component); 61c4762a1bSJed Brown if (!ghostvtex) { 62c4762a1bSJed Brown /* Handle reference bus constrained dofs */ 63c4762a1bSJed Brown if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) { 649371c9d4SSatish Balay row[0] = goffset; 659371c9d4SSatish Balay row[1] = goffset + 1; 669371c9d4SSatish Balay col[0] = goffset; 679371c9d4SSatish Balay col[1] = goffset + 1; 689371c9d4SSatish Balay col[2] = goffset; 699371c9d4SSatish Balay col[3] = goffset + 1; 709371c9d4SSatish Balay values[0] = 1.0; 719371c9d4SSatish Balay values[1] = 0.0; 729371c9d4SSatish Balay values[2] = 0.0; 739371c9d4SSatish Balay values[3] = 1.0; 749566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 2, row, 2, col, values, ADD_VALUES)); 75c4762a1bSJed Brown break; 76c4762a1bSJed Brown } 77c4762a1bSJed Brown 78c4762a1bSJed Brown Vm = xarr[offset + 1]; 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* Shunt injections */ 819371c9d4SSatish Balay row[0] = goffset; 829371c9d4SSatish Balay row[1] = goffset + 1; 839371c9d4SSatish Balay col[0] = goffset; 849371c9d4SSatish Balay col[1] = goffset + 1; 85c4762a1bSJed Brown values[0] = values[1] = values[2] = values[3] = 0.0; 86c4762a1bSJed Brown if (bus->ide != PV_BUS) { 87c4762a1bSJed Brown values[1] = 2.0 * Vm * bus->gl / Sbase; 88c4762a1bSJed Brown values[3] = -2.0 * Vm * bus->bl / Sbase; 89c4762a1bSJed Brown } 909566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 2, row, 2, col, values, ADD_VALUES)); 91c4762a1bSJed Brown } 92c4762a1bSJed Brown 939566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSupportingEdges(networkdm, vtx[v], &nconnedges, &connedges)); 94c4762a1bSJed Brown 95c4762a1bSJed Brown for (i = 0; i < nconnedges; i++) { 96c4762a1bSJed Brown EDGE_Power branch; 97c4762a1bSJed Brown VERTEX_Power busf, bust; 98c4762a1bSJed Brown PetscInt keyf, keyt; 99c4762a1bSJed Brown PetscScalar Gff, Bff, Gft, Bft, Gtf, Btf, Gtt, Btt; 100c4762a1bSJed Brown const PetscInt *cone; 101c4762a1bSJed Brown PetscScalar Vmf, Vmt, thetaf, thetat, thetaft, thetatf; 102c4762a1bSJed Brown 103c4762a1bSJed Brown e = connedges[i]; 1049566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, e, 0, &key, (void **)&branch, NULL)); 105c4762a1bSJed Brown if (!branch->status) continue; 106c4762a1bSJed Brown 107c4762a1bSJed Brown Gff = branch->yff[0]; 108c4762a1bSJed Brown Bff = branch->yff[1]; 109c4762a1bSJed Brown Gft = branch->yft[0]; 110c4762a1bSJed Brown Bft = branch->yft[1]; 111c4762a1bSJed Brown Gtf = branch->ytf[0]; 112c4762a1bSJed Brown Btf = branch->ytf[1]; 113c4762a1bSJed Brown Gtt = branch->ytt[0]; 114c4762a1bSJed Brown Btt = branch->ytt[1]; 115c4762a1bSJed Brown 1169566063dSJacob Faibussowitsch PetscCall(DMNetworkGetConnectedVertices(networkdm, edges[e], &cone)); 117c4762a1bSJed Brown vfrom = cone[0]; 118c4762a1bSJed Brown vto = cone[1]; 119c4762a1bSJed Brown 1209566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vfrom, ALL_COMPONENTS, &offsetfrom)); 1219566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vto, ALL_COMPONENTS, &offsetto)); 1229566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVecOffset(networkdm, vfrom, ALL_COMPONENTS, &goffsetfrom)); 1239566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVecOffset(networkdm, vto, ALL_COMPONENTS, &goffsetto)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown if (goffsetto < 0) goffsetto = -goffsetto - 1; 126c4762a1bSJed Brown 127c4762a1bSJed Brown thetaf = xarr[offsetfrom]; 128c4762a1bSJed Brown Vmf = xarr[offsetfrom + 1]; 129c4762a1bSJed Brown thetat = xarr[offsetto]; 130c4762a1bSJed Brown Vmt = xarr[offsetto + 1]; 131c4762a1bSJed Brown thetaft = thetaf - thetat; 132c4762a1bSJed Brown thetatf = thetat - thetaf; 133c4762a1bSJed Brown 1349566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vfrom, 0, &keyf, (void **)&busf, NULL)); 1359566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vto, 0, &keyt, (void **)&bust, NULL)); 136c4762a1bSJed Brown 137c4762a1bSJed Brown if (vfrom == vtx[v]) { 138c4762a1bSJed Brown if (busf->ide != REF_BUS) { 139c4762a1bSJed Brown /* farr[offsetfrom] += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); */ 140c4762a1bSJed Brown row[0] = goffsetfrom; 1419371c9d4SSatish Balay col[0] = goffsetfrom; 1429371c9d4SSatish Balay col[1] = goffsetfrom + 1; 1439371c9d4SSatish Balay col[2] = goffsetto; 1449371c9d4SSatish Balay col[3] = goffsetto + 1; 145c4762a1bSJed Brown values[0] = Vmf * Vmt * (Gft * -PetscSinScalar(thetaft) + Bft * PetscCosScalar(thetaft)); /* df_dthetaf */ 146c4762a1bSJed Brown values[1] = 2.0 * Gff * Vmf + Vmt * (Gft * PetscCosScalar(thetaft) + Bft * PetscSinScalar(thetaft)); /* df_dVmf */ 147c4762a1bSJed Brown values[2] = Vmf * Vmt * (Gft * PetscSinScalar(thetaft) + Bft * -PetscCosScalar(thetaft)); /* df_dthetat */ 148c4762a1bSJed Brown values[3] = Vmf * (Gft * PetscCosScalar(thetaft) + Bft * PetscSinScalar(thetaft)); /* df_dVmt */ 149c4762a1bSJed Brown 1509566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 4, col, values, ADD_VALUES)); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown if (busf->ide != PV_BUS && busf->ide != REF_BUS) { 153c4762a1bSJed Brown row[0] = goffsetfrom + 1; 1549371c9d4SSatish Balay col[0] = goffsetfrom; 1559371c9d4SSatish Balay col[1] = goffsetfrom + 1; 1569371c9d4SSatish Balay col[2] = goffsetto; 1579371c9d4SSatish Balay col[3] = goffsetto + 1; 158c4762a1bSJed Brown /* farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); */ 159c4762a1bSJed Brown values[0] = Vmf * Vmt * (Bft * PetscSinScalar(thetaft) + Gft * PetscCosScalar(thetaft)); 160c4762a1bSJed Brown values[1] = -2.0 * Bff * Vmf + Vmt * (-Bft * PetscCosScalar(thetaft) + Gft * PetscSinScalar(thetaft)); 161c4762a1bSJed Brown values[2] = Vmf * Vmt * (-Bft * PetscSinScalar(thetaft) + Gft * -PetscCosScalar(thetaft)); 162c4762a1bSJed Brown values[3] = Vmf * (-Bft * PetscCosScalar(thetaft) + Gft * PetscSinScalar(thetaft)); 163c4762a1bSJed Brown 1649566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 4, col, values, ADD_VALUES)); 165c4762a1bSJed Brown } 166c4762a1bSJed Brown } else { 167c4762a1bSJed Brown if (bust->ide != REF_BUS) { 168c4762a1bSJed Brown row[0] = goffsetto; 1699371c9d4SSatish Balay col[0] = goffsetto; 1709371c9d4SSatish Balay col[1] = goffsetto + 1; 1719371c9d4SSatish Balay col[2] = goffsetfrom; 1729371c9d4SSatish Balay col[3] = goffsetfrom + 1; 173c4762a1bSJed Brown /* farr[offsetto] += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); */ 174c4762a1bSJed Brown values[0] = Vmt * Vmf * (Gtf * -PetscSinScalar(thetatf) + Btf * PetscCosScalar(thetaft)); /* df_dthetat */ 175c4762a1bSJed Brown values[1] = 2.0 * Gtt * Vmt + Vmf * (Gtf * PetscCosScalar(thetatf) + Btf * PetscSinScalar(thetatf)); /* df_dVmt */ 176c4762a1bSJed Brown values[2] = Vmt * Vmf * (Gtf * PetscSinScalar(thetatf) + Btf * -PetscCosScalar(thetatf)); /* df_dthetaf */ 177c4762a1bSJed Brown values[3] = Vmt * (Gtf * PetscCosScalar(thetatf) + Btf * PetscSinScalar(thetatf)); /* df_dVmf */ 178c4762a1bSJed Brown 1799566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 4, col, values, ADD_VALUES)); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown if (bust->ide != PV_BUS && bust->ide != REF_BUS) { 182c4762a1bSJed Brown row[0] = goffsetto + 1; 1839371c9d4SSatish Balay col[0] = goffsetto; 1849371c9d4SSatish Balay col[1] = goffsetto + 1; 1859371c9d4SSatish Balay col[2] = goffsetfrom; 1869371c9d4SSatish Balay col[3] = goffsetfrom + 1; 187c4762a1bSJed Brown /* farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); */ 188c4762a1bSJed Brown values[0] = Vmt * Vmf * (Btf * PetscSinScalar(thetatf) + Gtf * PetscCosScalar(thetatf)); 189c4762a1bSJed Brown values[1] = -2.0 * Btt * Vmt + Vmf * (-Btf * PetscCosScalar(thetatf) + Gtf * PetscSinScalar(thetatf)); 190c4762a1bSJed Brown values[2] = Vmt * Vmf * (-Btf * PetscSinScalar(thetatf) + Gtf * -PetscCosScalar(thetatf)); 191c4762a1bSJed Brown values[3] = Vmt * (-Btf * PetscCosScalar(thetatf) + Gtf * PetscSinScalar(thetatf)); 192c4762a1bSJed Brown 1939566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 4, col, values, ADD_VALUES)); 194c4762a1bSJed Brown } 195c4762a1bSJed Brown } 196c4762a1bSJed Brown } 197c4762a1bSJed Brown if (!ghostvtex && bus->ide == PV_BUS) { 1989371c9d4SSatish Balay row[0] = goffset + 1; 1999371c9d4SSatish Balay col[0] = goffset + 1; 2009371c9d4SSatish Balay values[0] = 1.0; 201c4762a1bSJed Brown if (user_power->jac_error) values[0] = 50.0; 2029566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 1, col, values, ADD_VALUES)); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown } 205c4762a1bSJed Brown } 206c4762a1bSJed Brown } 207c4762a1bSJed Brown 2089566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(localX, &xarr)); 209c4762a1bSJed Brown PetscFunctionReturn(0); 210c4762a1bSJed Brown } 211c4762a1bSJed Brown 2129371c9d4SSatish Balay PetscErrorCode FormJacobian_Power(SNES snes, Vec X, Mat J, Mat Jpre, void *appctx) { 213c4762a1bSJed Brown DM networkdm; 214c4762a1bSJed Brown Vec localX; 215c4762a1bSJed Brown PetscInt nv, ne; 216c4762a1bSJed Brown const PetscInt *vtx, *edges; 217c4762a1bSJed Brown 218c4762a1bSJed Brown PetscFunctionBegin; 2199566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(J)); 220c4762a1bSJed Brown 2219566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &networkdm)); 2229566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm, &localX)); 223c4762a1bSJed Brown 2249566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX)); 2259566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX)); 226c4762a1bSJed Brown 2279566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges)); 2289566063dSJacob Faibussowitsch PetscCall(FormJacobian_Power_private(networkdm, localX, J, nv, ne, vtx, edges, appctx)); 229c4762a1bSJed Brown 2309566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm, &localX)); 231c4762a1bSJed Brown 2329566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 2339566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 234c4762a1bSJed Brown PetscFunctionReturn(0); 235c4762a1bSJed Brown } 236c4762a1bSJed Brown 2379371c9d4SSatish Balay PetscErrorCode FormFunction_Power(DM networkdm, Vec localX, Vec localF, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx) { 238c4762a1bSJed Brown UserCtx_Power *User = (UserCtx_Power *)appctx; 239c4762a1bSJed Brown PetscInt e, v, vfrom, vto; 240c4762a1bSJed Brown const PetscScalar *xarr; 241c4762a1bSJed Brown PetscScalar *farr; 242c4762a1bSJed Brown PetscInt offsetfrom, offsetto, offset, i, j, key, numComps; 243c4762a1bSJed Brown PetscScalar Vm; 244c4762a1bSJed Brown PetscScalar Sbase = User->Sbase; 245c4762a1bSJed Brown VERTEX_Power bus = NULL; 246c4762a1bSJed Brown GEN gen; 247c4762a1bSJed Brown LOAD load; 248c4762a1bSJed Brown PetscBool ghostvtex; 249c4762a1bSJed Brown void *component; 250c4762a1bSJed Brown 251c4762a1bSJed Brown PetscFunctionBegin; 2529566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localX, &xarr)); 2539566063dSJacob Faibussowitsch PetscCall(VecGetArray(localF, &farr)); 254c4762a1bSJed Brown 255c4762a1bSJed Brown for (v = 0; v < nv; v++) { 2569566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[v], &ghostvtex)); 2579566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm, vtx[v], &numComps)); 2589566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[v], ALL_COMPONENTS, &offset)); 259c4762a1bSJed Brown 260c4762a1bSJed Brown for (j = 0; j < numComps; j++) { 2619566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vtx[v], j, &key, &component, NULL)); 262c4762a1bSJed Brown if (key == User->compkey_bus) { 263c4762a1bSJed Brown PetscInt nconnedges; 264c4762a1bSJed Brown const PetscInt *connedges; 265c4762a1bSJed Brown 266c4762a1bSJed Brown bus = (VERTEX_Power)(component); 267c4762a1bSJed Brown /* Handle reference bus constrained dofs */ 268c4762a1bSJed Brown if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) { 269c4762a1bSJed Brown farr[offset] = xarr[offset] - bus->va * PETSC_PI / 180.0; 270c4762a1bSJed Brown farr[offset + 1] = xarr[offset + 1] - bus->vm; 271c4762a1bSJed Brown break; 272c4762a1bSJed Brown } 273c4762a1bSJed Brown 274c4762a1bSJed Brown if (!ghostvtex) { 275c4762a1bSJed Brown Vm = xarr[offset + 1]; 276c4762a1bSJed Brown 277c4762a1bSJed Brown /* Shunt injections */ 278c4762a1bSJed Brown farr[offset] += Vm * Vm * bus->gl / Sbase; 279c4762a1bSJed Brown if (bus->ide != PV_BUS) farr[offset + 1] += -Vm * Vm * bus->bl / Sbase; 280c4762a1bSJed Brown } 281c4762a1bSJed Brown 2829566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSupportingEdges(networkdm, vtx[v], &nconnedges, &connedges)); 283c4762a1bSJed Brown for (i = 0; i < nconnedges; i++) { 284c4762a1bSJed Brown EDGE_Power branch; 285c4762a1bSJed Brown PetscInt keye; 286c4762a1bSJed Brown PetscScalar Gff, Bff, Gft, Bft, Gtf, Btf, Gtt, Btt; 287c4762a1bSJed Brown const PetscInt *cone; 288c4762a1bSJed Brown PetscScalar Vmf, Vmt, thetaf, thetat, thetaft, thetatf; 289c4762a1bSJed Brown 290c4762a1bSJed Brown e = connedges[i]; 2919566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, e, 0, &keye, (void **)&branch, NULL)); 292c4762a1bSJed Brown if (!branch->status) continue; 293c4762a1bSJed Brown Gff = branch->yff[0]; 294c4762a1bSJed Brown Bff = branch->yff[1]; 295c4762a1bSJed Brown Gft = branch->yft[0]; 296c4762a1bSJed Brown Bft = branch->yft[1]; 297c4762a1bSJed Brown Gtf = branch->ytf[0]; 298c4762a1bSJed Brown Btf = branch->ytf[1]; 299c4762a1bSJed Brown Gtt = branch->ytt[0]; 300c4762a1bSJed Brown Btt = branch->ytt[1]; 301c4762a1bSJed Brown 3029566063dSJacob Faibussowitsch PetscCall(DMNetworkGetConnectedVertices(networkdm, e, &cone)); 303c4762a1bSJed Brown vfrom = cone[0]; 304c4762a1bSJed Brown vto = cone[1]; 305c4762a1bSJed Brown 3069566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vfrom, ALL_COMPONENTS, &offsetfrom)); 3079566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vto, ALL_COMPONENTS, &offsetto)); 308c4762a1bSJed Brown 309c4762a1bSJed Brown thetaf = xarr[offsetfrom]; 310c4762a1bSJed Brown Vmf = xarr[offsetfrom + 1]; 311c4762a1bSJed Brown thetat = xarr[offsetto]; 312c4762a1bSJed Brown Vmt = xarr[offsetto + 1]; 313c4762a1bSJed Brown thetaft = thetaf - thetat; 314c4762a1bSJed Brown thetatf = thetat - thetaf; 315c4762a1bSJed Brown 316c4762a1bSJed Brown if (vfrom == vtx[v]) { 317c4762a1bSJed Brown farr[offsetfrom] += Gff * Vmf * Vmf + Vmf * Vmt * (Gft * PetscCosScalar(thetaft) + Bft * PetscSinScalar(thetaft)); 318c4762a1bSJed Brown farr[offsetfrom + 1] += -Bff * Vmf * Vmf + Vmf * Vmt * (-Bft * PetscCosScalar(thetaft) + Gft * PetscSinScalar(thetaft)); 319c4762a1bSJed Brown } else { 320c4762a1bSJed Brown farr[offsetto] += Gtt * Vmt * Vmt + Vmt * Vmf * (Gtf * PetscCosScalar(thetatf) + Btf * PetscSinScalar(thetatf)); 321c4762a1bSJed Brown farr[offsetto + 1] += -Btt * Vmt * Vmt + Vmt * Vmf * (-Btf * PetscCosScalar(thetatf) + Gtf * PetscSinScalar(thetatf)); 322c4762a1bSJed Brown } 323c4762a1bSJed Brown } 324c4762a1bSJed Brown } else if (key == User->compkey_gen) { 325c4762a1bSJed Brown if (!ghostvtex) { 326c4762a1bSJed Brown gen = (GEN)(component); 327c4762a1bSJed Brown if (!gen->status) continue; 328c4762a1bSJed Brown farr[offset] += -gen->pg / Sbase; 329c4762a1bSJed Brown farr[offset + 1] += -gen->qg / Sbase; 330c4762a1bSJed Brown } 331c4762a1bSJed Brown } else if (key == User->compkey_load) { 332c4762a1bSJed Brown if (!ghostvtex) { 333c4762a1bSJed Brown load = (LOAD)(component); 334c4762a1bSJed Brown farr[offset] += load->pl / Sbase; 335c4762a1bSJed Brown farr[offset + 1] += load->ql / Sbase; 336c4762a1bSJed Brown } 337c4762a1bSJed Brown } 338c4762a1bSJed Brown } 3399371c9d4SSatish Balay if (bus && bus->ide == PV_BUS) { farr[offset + 1] = xarr[offset + 1] - bus->vm; } 340c4762a1bSJed Brown } 3419566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(localX, &xarr)); 3429566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localF, &farr)); 343c4762a1bSJed Brown PetscFunctionReturn(0); 344c4762a1bSJed Brown } 345c4762a1bSJed Brown 3469371c9d4SSatish Balay PetscErrorCode SetInitialGuess_Power(DM networkdm, Vec localX, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx) { 347c4762a1bSJed Brown VERTEX_Power bus; 348c4762a1bSJed Brown PetscInt i; 349c4762a1bSJed Brown GEN gen; 3502bf73ac6SHong Zhang PetscBool ghostvtex, sharedv; 351c4762a1bSJed Brown PetscScalar *xarr; 352c4762a1bSJed Brown PetscInt key, numComps, j, offset; 353c4762a1bSJed Brown void *component; 354c4762a1bSJed Brown PetscMPIInt rank; 355c4762a1bSJed Brown MPI_Comm comm; 356c4762a1bSJed Brown UserCtx_Power *User = (UserCtx_Power *)appctx; 357c4762a1bSJed Brown 358c4762a1bSJed Brown PetscFunctionBegin; 3599566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)networkdm, &comm)); 3609566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 3619566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX, &xarr)); 362c4762a1bSJed Brown for (i = 0; i < nv; i++) { 3639566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[i], &ghostvtex)); 3649566063dSJacob Faibussowitsch PetscCall(DMNetworkIsSharedVertex(networkdm, vtx[i], &sharedv)); 3652bf73ac6SHong Zhang if (ghostvtex || sharedv) continue; 366c4762a1bSJed Brown 3679566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[i], ALL_COMPONENTS, &offset)); 3689566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm, vtx[i], &numComps)); 369c4762a1bSJed Brown for (j = 0; j < numComps; j++) { 3709566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vtx[i], j, &key, &component, NULL)); 371c4762a1bSJed Brown if (key == User->compkey_bus) { 372c4762a1bSJed Brown bus = (VERTEX_Power)(component); 373c4762a1bSJed Brown xarr[offset] = bus->va * PETSC_PI / 180.0; 374c4762a1bSJed Brown xarr[offset + 1] = bus->vm; 375c4762a1bSJed Brown } else if (key == User->compkey_gen) { 376c4762a1bSJed Brown gen = (GEN)(component); 377c4762a1bSJed Brown if (!gen->status) continue; 378c4762a1bSJed Brown xarr[offset + 1] = gen->vs; 379c4762a1bSJed Brown break; 380c4762a1bSJed Brown } 381c4762a1bSJed Brown } 382c4762a1bSJed Brown } 3839566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX, &xarr)); 384c4762a1bSJed Brown PetscFunctionReturn(0); 385c4762a1bSJed Brown } 386