1c4762a1bSJed Brown /* function subroutines used by power.c */ 2c4762a1bSJed Brown 3c4762a1bSJed Brown #include "power.h" 4c4762a1bSJed Brown 5c4762a1bSJed Brown PetscErrorCode GetListofEdges_Power(PFDATA *pfdata,PetscInt *edgelist) 6c4762a1bSJed Brown { 7c4762a1bSJed Brown PetscInt i,fbus,tbus,nbranches=pfdata->nbranch; 8c4762a1bSJed Brown EDGE_Power branch=pfdata->branch; 9c4762a1bSJed Brown PetscBool netview=PETSC_FALSE; 10c4762a1bSJed Brown 11c4762a1bSJed Brown PetscFunctionBegin; 129566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL, "-powernet_view",&netview)); 13c4762a1bSJed Brown for (i=0; i<nbranches; i++) { 14c4762a1bSJed Brown fbus = branch[i].internal_i; 15c4762a1bSJed Brown tbus = branch[i].internal_j; 16c4762a1bSJed Brown edgelist[2*i] = fbus; 17c4762a1bSJed Brown edgelist[2*i+1] = tbus; 18c4762a1bSJed Brown if (netview) { 19*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"branch %" PetscInt_FMT ", bus[%" PetscInt_FMT "] -> bus[%" PetscInt_FMT "]\n",i,fbus,tbus)); 20c4762a1bSJed Brown } 21c4762a1bSJed Brown } 22c4762a1bSJed Brown if (netview) { 23c4762a1bSJed Brown for (i=0; i<pfdata->nbus; i++) { 24c4762a1bSJed Brown if (pfdata->bus[i].ngen) { 25*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," bus %" PetscInt_FMT ": gen\n",i)); 26c4762a1bSJed Brown } else if (pfdata->bus[i].nload) { 27*63a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF," bus %" PetscInt_FMT ": load\n",i)); 28c4762a1bSJed Brown } 29c4762a1bSJed Brown } 30c4762a1bSJed Brown } 31c4762a1bSJed Brown PetscFunctionReturn(0); 32c4762a1bSJed Brown } 33c4762a1bSJed Brown 34c4762a1bSJed Brown PetscErrorCode FormJacobian_Power_private(DM networkdm,Vec localX,Mat J,PetscInt nv,PetscInt ne,const PetscInt* vtx,const PetscInt* edges,void* appctx) 35c4762a1bSJed Brown { 36c4762a1bSJed Brown const PetscScalar *xarr; 37c4762a1bSJed Brown PetscInt i,v,row[2],col[8],e,vfrom,vto; 38c4762a1bSJed Brown PetscInt offsetfrom,offsetto,goffsetfrom,goffsetto,numComps; 39c4762a1bSJed Brown PetscScalar values[8]; 40c4762a1bSJed Brown PetscInt j,key,offset,goffset; 41c4762a1bSJed Brown PetscScalar Vm; 42c4762a1bSJed Brown UserCtx_Power *user_power=(UserCtx_Power*)appctx; 43c4762a1bSJed Brown PetscScalar Sbase=user_power->Sbase; 44c4762a1bSJed Brown VERTEX_Power bus; 45c4762a1bSJed Brown PetscBool ghostvtex; 46c4762a1bSJed Brown void* component; 47c4762a1bSJed Brown 48c4762a1bSJed Brown PetscFunctionBegin; 499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localX,&xarr)); 50c4762a1bSJed Brown 51c4762a1bSJed Brown for (v=0; v<nv; v++) { 529566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm,vtx[v],&ghostvtex)); 53c4762a1bSJed Brown 549566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm,vtx[v],&numComps)); 55c4762a1bSJed Brown for (j = 0; j < numComps; j++) { 569566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm,vtx[v],ALL_COMPONENTS,&offset)); 579566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVecOffset(networkdm,vtx[v],ALL_COMPONENTS,&goffset)); 589566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm,vtx[v],j,&key,&component,NULL)); 59c4762a1bSJed Brown 60c4762a1bSJed Brown if (key == user_power->compkey_bus) { 61c4762a1bSJed Brown PetscInt nconnedges; 62c4762a1bSJed Brown const PetscInt *connedges; 63c4762a1bSJed Brown 64c4762a1bSJed Brown bus = (VERTEX_Power)(component); 65c4762a1bSJed Brown if (!ghostvtex) { 66c4762a1bSJed Brown /* Handle reference bus constrained dofs */ 67c4762a1bSJed Brown if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) { 68c4762a1bSJed Brown row[0] = goffset; row[1] = goffset+1; 69c4762a1bSJed Brown col[0] = goffset; col[1] = goffset+1; col[2] = goffset; col[3] = goffset+1; 70c4762a1bSJed Brown values[0] = 1.0; values[1] = 0.0; values[2] = 0.0; values[3] = 1.0; 719566063dSJacob Faibussowitsch PetscCall(MatSetValues(J,2,row,2,col,values,ADD_VALUES)); 72c4762a1bSJed Brown break; 73c4762a1bSJed Brown } 74c4762a1bSJed Brown 75c4762a1bSJed Brown Vm = xarr[offset+1]; 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* Shunt injections */ 78c4762a1bSJed Brown row[0] = goffset; row[1] = goffset+1; 79c4762a1bSJed Brown col[0] = goffset; col[1] = goffset+1; 80c4762a1bSJed Brown values[0] = values[1] = values[2] = values[3] = 0.0; 81c4762a1bSJed Brown if (bus->ide != PV_BUS) { 82c4762a1bSJed Brown values[1] = 2.0*Vm*bus->gl/Sbase; 83c4762a1bSJed Brown values[3] = -2.0*Vm*bus->bl/Sbase; 84c4762a1bSJed Brown } 859566063dSJacob Faibussowitsch PetscCall(MatSetValues(J,2,row,2,col,values,ADD_VALUES)); 86c4762a1bSJed Brown } 87c4762a1bSJed Brown 889566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSupportingEdges(networkdm,vtx[v],&nconnedges,&connedges)); 89c4762a1bSJed Brown 90c4762a1bSJed Brown for (i=0; i < nconnedges; i++) { 91c4762a1bSJed Brown EDGE_Power branch; 92c4762a1bSJed Brown VERTEX_Power busf,bust; 93c4762a1bSJed Brown PetscInt keyf,keyt; 94c4762a1bSJed Brown PetscScalar Gff,Bff,Gft,Bft,Gtf,Btf,Gtt,Btt; 95c4762a1bSJed Brown const PetscInt *cone; 96c4762a1bSJed Brown PetscScalar Vmf,Vmt,thetaf,thetat,thetaft,thetatf; 97c4762a1bSJed Brown 98c4762a1bSJed Brown e = connedges[i]; 999566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm,e,0,&key,(void**)&branch,NULL)); 100c4762a1bSJed Brown if (!branch->status) continue; 101c4762a1bSJed Brown 102c4762a1bSJed Brown Gff = branch->yff[0]; 103c4762a1bSJed Brown Bff = branch->yff[1]; 104c4762a1bSJed Brown Gft = branch->yft[0]; 105c4762a1bSJed Brown Bft = branch->yft[1]; 106c4762a1bSJed Brown Gtf = branch->ytf[0]; 107c4762a1bSJed Brown Btf = branch->ytf[1]; 108c4762a1bSJed Brown Gtt = branch->ytt[0]; 109c4762a1bSJed Brown Btt = branch->ytt[1]; 110c4762a1bSJed Brown 1119566063dSJacob Faibussowitsch PetscCall(DMNetworkGetConnectedVertices(networkdm,edges[e],&cone)); 112c4762a1bSJed Brown vfrom = cone[0]; 113c4762a1bSJed Brown vto = cone[1]; 114c4762a1bSJed Brown 1159566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&offsetfrom)); 1169566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm,vto,ALL_COMPONENTS,&offsetto)); 1179566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&goffsetfrom)); 1189566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVecOffset(networkdm,vto,ALL_COMPONENTS,&goffsetto)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown if (goffsetto < 0) goffsetto = -goffsetto - 1; 121c4762a1bSJed Brown 122c4762a1bSJed Brown thetaf = xarr[offsetfrom]; 123c4762a1bSJed Brown Vmf = xarr[offsetfrom+1]; 124c4762a1bSJed Brown thetat = xarr[offsetto]; 125c4762a1bSJed Brown Vmt = xarr[offsetto+1]; 126c4762a1bSJed Brown thetaft = thetaf - thetat; 127c4762a1bSJed Brown thetatf = thetat - thetaf; 128c4762a1bSJed Brown 1299566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm,vfrom,0,&keyf,(void**)&busf,NULL)); 1309566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm,vto,0,&keyt,(void**)&bust,NULL)); 131c4762a1bSJed Brown 132c4762a1bSJed Brown if (vfrom == vtx[v]) { 133c4762a1bSJed Brown if (busf->ide != REF_BUS) { 134c4762a1bSJed Brown /* farr[offsetfrom] += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); */ 135c4762a1bSJed Brown row[0] = goffsetfrom; 136c4762a1bSJed Brown col[0] = goffsetfrom; col[1] = goffsetfrom+1; col[2] = goffsetto; col[3] = goffsetto+1; 137c4762a1bSJed Brown values[0] = Vmf*Vmt*(Gft*-PetscSinScalar(thetaft) + Bft*PetscCosScalar(thetaft)); /* df_dthetaf */ 138c4762a1bSJed Brown values[1] = 2.0*Gff*Vmf + Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); /* df_dVmf */ 139c4762a1bSJed Brown values[2] = Vmf*Vmt*(Gft*PetscSinScalar(thetaft) + Bft*-PetscCosScalar(thetaft)); /* df_dthetat */ 140c4762a1bSJed Brown values[3] = Vmf*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); /* df_dVmt */ 141c4762a1bSJed Brown 1429566063dSJacob Faibussowitsch PetscCall(MatSetValues(J,1,row,4,col,values,ADD_VALUES)); 143c4762a1bSJed Brown } 144c4762a1bSJed Brown if (busf->ide != PV_BUS && busf->ide != REF_BUS) { 145c4762a1bSJed Brown row[0] = goffsetfrom+1; 146c4762a1bSJed Brown col[0] = goffsetfrom; col[1] = goffsetfrom+1; col[2] = goffsetto; col[3] = goffsetto+1; 147c4762a1bSJed Brown /* farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); */ 148c4762a1bSJed Brown values[0] = Vmf*Vmt*(Bft*PetscSinScalar(thetaft) + Gft*PetscCosScalar(thetaft)); 149c4762a1bSJed Brown values[1] = -2.0*Bff*Vmf + Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); 150c4762a1bSJed Brown values[2] = Vmf*Vmt*(-Bft*PetscSinScalar(thetaft) + Gft*-PetscCosScalar(thetaft)); 151c4762a1bSJed Brown values[3] = Vmf*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); 152c4762a1bSJed Brown 1539566063dSJacob Faibussowitsch PetscCall(MatSetValues(J,1,row,4,col,values,ADD_VALUES)); 154c4762a1bSJed Brown } 155c4762a1bSJed Brown } else { 156c4762a1bSJed Brown if (bust->ide != REF_BUS) { 157c4762a1bSJed Brown row[0] = goffsetto; 158c4762a1bSJed Brown col[0] = goffsetto; col[1] = goffsetto+1; col[2] = goffsetfrom; col[3] = goffsetfrom+1; 159c4762a1bSJed Brown /* farr[offsetto] += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); */ 160c4762a1bSJed Brown values[0] = Vmt*Vmf*(Gtf*-PetscSinScalar(thetatf) + Btf*PetscCosScalar(thetaft)); /* df_dthetat */ 161c4762a1bSJed Brown values[1] = 2.0*Gtt*Vmt + Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); /* df_dVmt */ 162c4762a1bSJed Brown values[2] = Vmt*Vmf*(Gtf*PetscSinScalar(thetatf) + Btf*-PetscCosScalar(thetatf)); /* df_dthetaf */ 163c4762a1bSJed Brown values[3] = Vmt*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); /* df_dVmf */ 164c4762a1bSJed Brown 1659566063dSJacob Faibussowitsch PetscCall(MatSetValues(J,1,row,4,col,values,ADD_VALUES)); 166c4762a1bSJed Brown } 167c4762a1bSJed Brown if (bust->ide != PV_BUS && bust->ide != REF_BUS) { 168c4762a1bSJed Brown row[0] = goffsetto+1; 169c4762a1bSJed Brown col[0] = goffsetto; col[1] = goffsetto+1; col[2] = goffsetfrom; col[3] = goffsetfrom+1; 170c4762a1bSJed Brown /* farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); */ 171c4762a1bSJed Brown values[0] = Vmt*Vmf*(Btf*PetscSinScalar(thetatf) + Gtf*PetscCosScalar(thetatf)); 172c4762a1bSJed Brown values[1] = -2.0*Btt*Vmt + Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); 173c4762a1bSJed Brown values[2] = Vmt*Vmf*(-Btf*PetscSinScalar(thetatf) + Gtf*-PetscCosScalar(thetatf)); 174c4762a1bSJed Brown values[3] = Vmt*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); 175c4762a1bSJed Brown 1769566063dSJacob Faibussowitsch PetscCall(MatSetValues(J,1,row,4,col,values,ADD_VALUES)); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown } 179c4762a1bSJed Brown } 180c4762a1bSJed Brown if (!ghostvtex && bus->ide == PV_BUS) { 181c4762a1bSJed Brown row[0] = goffset+1; col[0] = goffset+1; values[0] = 1.0; 182c4762a1bSJed Brown if (user_power->jac_error) values[0] = 50.0; 1839566063dSJacob Faibussowitsch PetscCall(MatSetValues(J,1,row,1,col,values,ADD_VALUES)); 184c4762a1bSJed Brown } 185c4762a1bSJed Brown } 186c4762a1bSJed Brown } 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 1899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(localX,&xarr)); 190c4762a1bSJed Brown PetscFunctionReturn(0); 191c4762a1bSJed Brown } 192c4762a1bSJed Brown 193c4762a1bSJed Brown PetscErrorCode FormJacobian_Power(SNES snes,Vec X, Mat J,Mat Jpre,void *appctx) 194c4762a1bSJed Brown { 195c4762a1bSJed Brown DM networkdm; 196c4762a1bSJed Brown Vec localX; 197c4762a1bSJed Brown PetscInt nv,ne; 198c4762a1bSJed Brown const PetscInt *vtx,*edges; 199c4762a1bSJed Brown 200c4762a1bSJed Brown PetscFunctionBegin; 2019566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(J)); 202c4762a1bSJed Brown 2039566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&networkdm)); 2049566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm,&localX)); 205c4762a1bSJed Brown 2069566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX)); 2079566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX)); 208c4762a1bSJed Brown 2099566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges)); 2109566063dSJacob Faibussowitsch PetscCall(FormJacobian_Power_private(networkdm,localX,J,nv,ne,vtx,edges,appctx)); 211c4762a1bSJed Brown 2129566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm,&localX)); 213c4762a1bSJed Brown 2149566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 2159566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 216c4762a1bSJed Brown PetscFunctionReturn(0); 217c4762a1bSJed Brown } 218c4762a1bSJed Brown 219c4762a1bSJed Brown PetscErrorCode FormFunction_Power(DM networkdm,Vec localX, Vec localF,PetscInt nv,PetscInt ne,const PetscInt* vtx,const PetscInt* edges,void* appctx) 220c4762a1bSJed Brown { 221c4762a1bSJed Brown UserCtx_Power *User=(UserCtx_Power*)appctx; 222c4762a1bSJed Brown PetscInt e,v,vfrom,vto; 223c4762a1bSJed Brown const PetscScalar *xarr; 224c4762a1bSJed Brown PetscScalar *farr; 225c4762a1bSJed Brown PetscInt offsetfrom,offsetto,offset,i,j,key,numComps; 226c4762a1bSJed Brown PetscScalar Vm; 227c4762a1bSJed Brown PetscScalar Sbase=User->Sbase; 228c4762a1bSJed Brown VERTEX_Power bus=NULL; 229c4762a1bSJed Brown GEN gen; 230c4762a1bSJed Brown LOAD load; 231c4762a1bSJed Brown PetscBool ghostvtex; 232c4762a1bSJed Brown void* component; 233c4762a1bSJed Brown 234c4762a1bSJed Brown PetscFunctionBegin; 2359566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localX,&xarr)); 2369566063dSJacob Faibussowitsch PetscCall(VecGetArray(localF,&farr)); 237c4762a1bSJed Brown 238c4762a1bSJed Brown for (v=0; v<nv; v++) { 2399566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm,vtx[v],&ghostvtex)); 2409566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm,vtx[v],&numComps)); 2419566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm,vtx[v],ALL_COMPONENTS,&offset)); 242c4762a1bSJed Brown 243c4762a1bSJed Brown for (j = 0; j < numComps; j++) { 2449566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm,vtx[v],j,&key,&component,NULL)); 245c4762a1bSJed Brown if (key == User->compkey_bus) { 246c4762a1bSJed Brown PetscInt nconnedges; 247c4762a1bSJed Brown const PetscInt *connedges; 248c4762a1bSJed Brown 249c4762a1bSJed Brown bus = (VERTEX_Power)(component); 250c4762a1bSJed Brown /* Handle reference bus constrained dofs */ 251c4762a1bSJed Brown if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) { 252c4762a1bSJed Brown farr[offset] = xarr[offset] - bus->va*PETSC_PI/180.0; 253c4762a1bSJed Brown farr[offset+1] = xarr[offset+1] - bus->vm; 254c4762a1bSJed Brown break; 255c4762a1bSJed Brown } 256c4762a1bSJed Brown 257c4762a1bSJed Brown if (!ghostvtex) { 258c4762a1bSJed Brown Vm = xarr[offset+1]; 259c4762a1bSJed Brown 260c4762a1bSJed Brown /* Shunt injections */ 261c4762a1bSJed Brown farr[offset] += Vm*Vm*bus->gl/Sbase; 262c4762a1bSJed Brown if (bus->ide != PV_BUS) farr[offset+1] += -Vm*Vm*bus->bl/Sbase; 263c4762a1bSJed Brown } 264c4762a1bSJed Brown 2659566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSupportingEdges(networkdm,vtx[v],&nconnedges,&connedges)); 266c4762a1bSJed Brown for (i=0; i < nconnedges; i++) { 267c4762a1bSJed Brown EDGE_Power branch; 268c4762a1bSJed Brown PetscInt keye; 269c4762a1bSJed Brown PetscScalar Gff,Bff,Gft,Bft,Gtf,Btf,Gtt,Btt; 270c4762a1bSJed Brown const PetscInt *cone; 271c4762a1bSJed Brown PetscScalar Vmf,Vmt,thetaf,thetat,thetaft,thetatf; 272c4762a1bSJed Brown 273c4762a1bSJed Brown e = connedges[i]; 2749566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm,e,0,&keye,(void**)&branch,NULL)); 275c4762a1bSJed Brown if (!branch->status) continue; 276c4762a1bSJed Brown Gff = branch->yff[0]; 277c4762a1bSJed Brown Bff = branch->yff[1]; 278c4762a1bSJed Brown Gft = branch->yft[0]; 279c4762a1bSJed Brown Bft = branch->yft[1]; 280c4762a1bSJed Brown Gtf = branch->ytf[0]; 281c4762a1bSJed Brown Btf = branch->ytf[1]; 282c4762a1bSJed Brown Gtt = branch->ytt[0]; 283c4762a1bSJed Brown Btt = branch->ytt[1]; 284c4762a1bSJed Brown 2859566063dSJacob Faibussowitsch PetscCall(DMNetworkGetConnectedVertices(networkdm,e,&cone)); 286c4762a1bSJed Brown vfrom = cone[0]; 287c4762a1bSJed Brown vto = cone[1]; 288c4762a1bSJed Brown 2899566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&offsetfrom)); 2909566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm,vto,ALL_COMPONENTS,&offsetto)); 291c4762a1bSJed Brown 292c4762a1bSJed Brown thetaf = xarr[offsetfrom]; 293c4762a1bSJed Brown Vmf = xarr[offsetfrom+1]; 294c4762a1bSJed Brown thetat = xarr[offsetto]; 295c4762a1bSJed Brown Vmt = xarr[offsetto+1]; 296c4762a1bSJed Brown thetaft = thetaf - thetat; 297c4762a1bSJed Brown thetatf = thetat - thetaf; 298c4762a1bSJed Brown 299c4762a1bSJed Brown if (vfrom == vtx[v]) { 300c4762a1bSJed Brown farr[offsetfrom] += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); 301c4762a1bSJed Brown farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); 302c4762a1bSJed Brown } else { 303c4762a1bSJed Brown farr[offsetto] += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); 304c4762a1bSJed Brown farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); 305c4762a1bSJed Brown } 306c4762a1bSJed Brown } 307c4762a1bSJed Brown } else if (key == User->compkey_gen) { 308c4762a1bSJed Brown if (!ghostvtex) { 309c4762a1bSJed Brown gen = (GEN)(component); 310c4762a1bSJed Brown if (!gen->status) continue; 311c4762a1bSJed Brown farr[offset] += -gen->pg/Sbase; 312c4762a1bSJed Brown farr[offset+1] += -gen->qg/Sbase; 313c4762a1bSJed Brown } 314c4762a1bSJed Brown } else if (key == User->compkey_load) { 315c4762a1bSJed Brown if (!ghostvtex) { 316c4762a1bSJed Brown load = (LOAD)(component); 317c4762a1bSJed Brown farr[offset] += load->pl/Sbase; 318c4762a1bSJed Brown farr[offset+1] += load->ql/Sbase; 319c4762a1bSJed Brown } 320c4762a1bSJed Brown } 321c4762a1bSJed Brown } 322c4762a1bSJed Brown if (bus && bus->ide == PV_BUS) { 323c4762a1bSJed Brown farr[offset+1] = xarr[offset+1] - bus->vm; 324c4762a1bSJed Brown } 325c4762a1bSJed Brown } 3269566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(localX,&xarr)); 3279566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localF,&farr)); 328c4762a1bSJed Brown PetscFunctionReturn(0); 329c4762a1bSJed Brown } 330c4762a1bSJed Brown 331c4762a1bSJed Brown PetscErrorCode SetInitialGuess_Power(DM networkdm,Vec localX,PetscInt nv,PetscInt ne,const PetscInt *vtx,const PetscInt *edges,void* appctx) 332c4762a1bSJed Brown { 333c4762a1bSJed Brown VERTEX_Power bus; 334c4762a1bSJed Brown PetscInt i; 335c4762a1bSJed Brown GEN gen; 3362bf73ac6SHong Zhang PetscBool ghostvtex,sharedv; 337c4762a1bSJed Brown PetscScalar *xarr; 338c4762a1bSJed Brown PetscInt key,numComps,j,offset; 339c4762a1bSJed Brown void* component; 340c4762a1bSJed Brown PetscMPIInt rank; 341c4762a1bSJed Brown MPI_Comm comm; 342c4762a1bSJed Brown UserCtx_Power *User=(UserCtx_Power*)appctx; 343c4762a1bSJed Brown 344c4762a1bSJed Brown PetscFunctionBegin; 3459566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)networkdm,&comm)); 3469566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm,&rank)); 3479566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX,&xarr)); 348c4762a1bSJed Brown for (i = 0; i < nv; i++) { 3499566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm,vtx[i],&ghostvtex)); 3509566063dSJacob Faibussowitsch PetscCall(DMNetworkIsSharedVertex(networkdm,vtx[i],&sharedv)); 3512bf73ac6SHong Zhang if (ghostvtex ||sharedv) continue; 352c4762a1bSJed Brown 3539566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm,vtx[i],ALL_COMPONENTS,&offset)); 3549566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm,vtx[i],&numComps)); 355c4762a1bSJed Brown for (j=0; j < numComps; j++) { 3569566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm,vtx[i],j,&key,&component,NULL)); 357c4762a1bSJed Brown if (key == User->compkey_bus) { 358c4762a1bSJed Brown bus = (VERTEX_Power)(component); 359c4762a1bSJed Brown xarr[offset] = bus->va*PETSC_PI/180.0; 360c4762a1bSJed Brown xarr[offset+1] = bus->vm; 361c4762a1bSJed Brown } else if (key == User->compkey_gen) { 362c4762a1bSJed Brown gen = (GEN)(component); 363c4762a1bSJed Brown if (!gen->status) continue; 364c4762a1bSJed Brown xarr[offset+1] = gen->vs; 365c4762a1bSJed Brown break; 366c4762a1bSJed Brown } 367c4762a1bSJed Brown } 368c4762a1bSJed Brown } 3699566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX,&xarr)); 370c4762a1bSJed Brown PetscFunctionReturn(0); 371c4762a1bSJed Brown } 372