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 PetscErrorCode ierr; 8c4762a1bSJed Brown PetscInt i,fbus,tbus,nbranches=pfdata->nbranch; 9c4762a1bSJed Brown EDGE_Power branch=pfdata->branch; 10c4762a1bSJed Brown PetscBool netview=PETSC_FALSE; 11c4762a1bSJed Brown 12c4762a1bSJed Brown PetscFunctionBegin; 13c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-powernet_view",&netview);CHKERRQ(ierr); 14c4762a1bSJed Brown for (i=0; i<nbranches; i++) { 15c4762a1bSJed Brown fbus = branch[i].internal_i; 16c4762a1bSJed Brown tbus = branch[i].internal_j; 17c4762a1bSJed Brown edgelist[2*i] = fbus; 18c4762a1bSJed Brown edgelist[2*i+1] = tbus; 19c4762a1bSJed Brown if (netview) { 20c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"branch %d, bus[%d] -> bus[%d]\n",i,fbus,tbus);CHKERRQ(ierr); 21c4762a1bSJed Brown } 22c4762a1bSJed Brown } 23c4762a1bSJed Brown if (netview) { 24c4762a1bSJed Brown for (i=0; i<pfdata->nbus; i++) { 25c4762a1bSJed Brown if (pfdata->bus[i].ngen) { 26c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," bus %D: gen\n",i);CHKERRQ(ierr); 27c4762a1bSJed Brown } else if (pfdata->bus[i].nload) { 28c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF," bus %D: load\n",i);CHKERRQ(ierr); 29c4762a1bSJed Brown } 30c4762a1bSJed Brown } 31c4762a1bSJed Brown } 32c4762a1bSJed Brown PetscFunctionReturn(0); 33c4762a1bSJed Brown } 34c4762a1bSJed Brown 35c4762a1bSJed Brown PetscErrorCode FormJacobian_Power_private(DM networkdm,Vec localX,Mat J,PetscInt nv,PetscInt ne,const PetscInt* vtx,const PetscInt* edges,void* appctx) 36c4762a1bSJed Brown { 37c4762a1bSJed Brown PetscErrorCode ierr; 38c4762a1bSJed Brown const PetscScalar *xarr; 39c4762a1bSJed Brown PetscInt i,v,row[2],col[8],e,vfrom,vto; 40c4762a1bSJed Brown PetscInt offsetfrom,offsetto,goffsetfrom,goffsetto,numComps; 41c4762a1bSJed Brown PetscScalar values[8]; 42c4762a1bSJed Brown PetscInt j,key,offset,goffset; 43c4762a1bSJed Brown PetscScalar Vm; 44c4762a1bSJed Brown UserCtx_Power *user_power=(UserCtx_Power*)appctx; 45c4762a1bSJed Brown PetscScalar Sbase=user_power->Sbase; 46c4762a1bSJed Brown VERTEX_Power bus; 47c4762a1bSJed Brown PetscBool ghostvtex; 48c4762a1bSJed Brown void* component; 49c4762a1bSJed Brown 50c4762a1bSJed Brown PetscFunctionBegin; 51c4762a1bSJed Brown ierr = VecGetArrayRead(localX,&xarr);CHKERRQ(ierr); 52c4762a1bSJed Brown 53c4762a1bSJed Brown for (v=0; v<nv; v++) { 54c4762a1bSJed Brown ierr = DMNetworkIsGhostVertex(networkdm,vtx[v],&ghostvtex);CHKERRQ(ierr); 55c4762a1bSJed Brown 56c4762a1bSJed Brown ierr = DMNetworkGetNumComponents(networkdm,vtx[v],&numComps);CHKERRQ(ierr); 57c4762a1bSJed Brown for (j = 0; j < numComps; j++) { 58*2bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vtx[v],ALL_COMPONENTS,&offset);CHKERRQ(ierr); 59*2bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVecOffset(networkdm,vtx[v],ALL_COMPONENTS,&goffset);CHKERRQ(ierr); 60*2bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,vtx[v],j,&key,&component,NULL);CHKERRQ(ierr); 61c4762a1bSJed Brown 62c4762a1bSJed Brown if (key == user_power->compkey_bus) { 63c4762a1bSJed Brown PetscInt nconnedges; 64c4762a1bSJed Brown const PetscInt *connedges; 65c4762a1bSJed Brown 66c4762a1bSJed Brown bus = (VERTEX_Power)(component); 67c4762a1bSJed Brown if (!ghostvtex) { 68c4762a1bSJed Brown /* Handle reference bus constrained dofs */ 69c4762a1bSJed Brown if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) { 70c4762a1bSJed Brown row[0] = goffset; row[1] = goffset+1; 71c4762a1bSJed Brown col[0] = goffset; col[1] = goffset+1; col[2] = goffset; col[3] = goffset+1; 72c4762a1bSJed Brown values[0] = 1.0; values[1] = 0.0; values[2] = 0.0; values[3] = 1.0; 73c4762a1bSJed Brown ierr = MatSetValues(J,2,row,2,col,values,ADD_VALUES);CHKERRQ(ierr); 74c4762a1bSJed Brown break; 75c4762a1bSJed Brown } 76c4762a1bSJed Brown 77c4762a1bSJed Brown Vm = xarr[offset+1]; 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Shunt injections */ 80c4762a1bSJed Brown row[0] = goffset; row[1] = goffset+1; 81c4762a1bSJed Brown col[0] = goffset; col[1] = goffset+1; 82c4762a1bSJed Brown values[0] = values[1] = values[2] = values[3] = 0.0; 83c4762a1bSJed Brown if (bus->ide != PV_BUS) { 84c4762a1bSJed Brown values[1] = 2.0*Vm*bus->gl/Sbase; 85c4762a1bSJed Brown values[3] = -2.0*Vm*bus->bl/Sbase; 86c4762a1bSJed Brown } 87c4762a1bSJed Brown ierr = MatSetValues(J,2,row,2,col,values,ADD_VALUES);CHKERRQ(ierr); 88c4762a1bSJed Brown } 89c4762a1bSJed Brown 90c4762a1bSJed Brown ierr = DMNetworkGetSupportingEdges(networkdm,vtx[v],&nconnedges,&connedges);CHKERRQ(ierr); 91c4762a1bSJed Brown 92c4762a1bSJed Brown for (i=0; i < nconnedges; i++) { 93c4762a1bSJed Brown EDGE_Power branch; 94c4762a1bSJed Brown VERTEX_Power busf,bust; 95c4762a1bSJed Brown PetscInt keyf,keyt; 96c4762a1bSJed Brown PetscScalar Gff,Bff,Gft,Bft,Gtf,Btf,Gtt,Btt; 97c4762a1bSJed Brown const PetscInt *cone; 98c4762a1bSJed Brown PetscScalar Vmf,Vmt,thetaf,thetat,thetaft,thetatf; 99c4762a1bSJed Brown 100c4762a1bSJed Brown e = connedges[i]; 101*2bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,e,0,&key,(void**)&branch,NULL);CHKERRQ(ierr); 102c4762a1bSJed Brown if (!branch->status) continue; 103c4762a1bSJed Brown 104c4762a1bSJed Brown Gff = branch->yff[0]; 105c4762a1bSJed Brown Bff = branch->yff[1]; 106c4762a1bSJed Brown Gft = branch->yft[0]; 107c4762a1bSJed Brown Bft = branch->yft[1]; 108c4762a1bSJed Brown Gtf = branch->ytf[0]; 109c4762a1bSJed Brown Btf = branch->ytf[1]; 110c4762a1bSJed Brown Gtt = branch->ytt[0]; 111c4762a1bSJed Brown Btt = branch->ytt[1]; 112c4762a1bSJed Brown 113c4762a1bSJed Brown ierr = DMNetworkGetConnectedVertices(networkdm,edges[e],&cone);CHKERRQ(ierr); 114c4762a1bSJed Brown vfrom = cone[0]; 115c4762a1bSJed Brown vto = cone[1]; 116c4762a1bSJed Brown 117*2bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&offsetfrom);CHKERRQ(ierr); 118*2bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vto,ALL_COMPONENTS,&offsetto);CHKERRQ(ierr); 119*2bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&goffsetfrom);CHKERRQ(ierr); 120*2bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVecOffset(networkdm,vto,ALL_COMPONENTS,&goffsetto);CHKERRQ(ierr); 121c4762a1bSJed Brown 122c4762a1bSJed Brown if (goffsetto < 0) goffsetto = -goffsetto - 1; 123c4762a1bSJed Brown 124c4762a1bSJed Brown thetaf = xarr[offsetfrom]; 125c4762a1bSJed Brown Vmf = xarr[offsetfrom+1]; 126c4762a1bSJed Brown thetat = xarr[offsetto]; 127c4762a1bSJed Brown Vmt = xarr[offsetto+1]; 128c4762a1bSJed Brown thetaft = thetaf - thetat; 129c4762a1bSJed Brown thetatf = thetat - thetaf; 130c4762a1bSJed Brown 131*2bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,vfrom,0,&keyf,(void**)&busf,NULL);CHKERRQ(ierr); 132*2bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,vto,0,&keyt,(void**)&bust,NULL);CHKERRQ(ierr); 133c4762a1bSJed Brown 134c4762a1bSJed Brown if (vfrom == vtx[v]) { 135c4762a1bSJed Brown if (busf->ide != REF_BUS) { 136c4762a1bSJed Brown /* farr[offsetfrom] += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); */ 137c4762a1bSJed Brown row[0] = goffsetfrom; 138c4762a1bSJed Brown col[0] = goffsetfrom; col[1] = goffsetfrom+1; col[2] = goffsetto; col[3] = goffsetto+1; 139c4762a1bSJed Brown values[0] = Vmf*Vmt*(Gft*-PetscSinScalar(thetaft) + Bft*PetscCosScalar(thetaft)); /* df_dthetaf */ 140c4762a1bSJed Brown values[1] = 2.0*Gff*Vmf + Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); /* df_dVmf */ 141c4762a1bSJed Brown values[2] = Vmf*Vmt*(Gft*PetscSinScalar(thetaft) + Bft*-PetscCosScalar(thetaft)); /* df_dthetat */ 142c4762a1bSJed Brown values[3] = Vmf*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); /* df_dVmt */ 143c4762a1bSJed Brown 144c4762a1bSJed Brown ierr = MatSetValues(J,1,row,4,col,values,ADD_VALUES);CHKERRQ(ierr); 145c4762a1bSJed Brown } 146c4762a1bSJed Brown if (busf->ide != PV_BUS && busf->ide != REF_BUS) { 147c4762a1bSJed Brown row[0] = goffsetfrom+1; 148c4762a1bSJed Brown col[0] = goffsetfrom; col[1] = goffsetfrom+1; col[2] = goffsetto; col[3] = goffsetto+1; 149c4762a1bSJed Brown /* farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); */ 150c4762a1bSJed Brown values[0] = Vmf*Vmt*(Bft*PetscSinScalar(thetaft) + Gft*PetscCosScalar(thetaft)); 151c4762a1bSJed Brown values[1] = -2.0*Bff*Vmf + Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); 152c4762a1bSJed Brown values[2] = Vmf*Vmt*(-Bft*PetscSinScalar(thetaft) + Gft*-PetscCosScalar(thetaft)); 153c4762a1bSJed Brown values[3] = Vmf*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); 154c4762a1bSJed Brown 155c4762a1bSJed Brown ierr = MatSetValues(J,1,row,4,col,values,ADD_VALUES);CHKERRQ(ierr); 156c4762a1bSJed Brown } 157c4762a1bSJed Brown } else { 158c4762a1bSJed Brown if (bust->ide != REF_BUS) { 159c4762a1bSJed Brown row[0] = goffsetto; 160c4762a1bSJed Brown col[0] = goffsetto; col[1] = goffsetto+1; col[2] = goffsetfrom; col[3] = goffsetfrom+1; 161c4762a1bSJed Brown /* farr[offsetto] += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); */ 162c4762a1bSJed Brown values[0] = Vmt*Vmf*(Gtf*-PetscSinScalar(thetatf) + Btf*PetscCosScalar(thetaft)); /* df_dthetat */ 163c4762a1bSJed Brown values[1] = 2.0*Gtt*Vmt + Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); /* df_dVmt */ 164c4762a1bSJed Brown values[2] = Vmt*Vmf*(Gtf*PetscSinScalar(thetatf) + Btf*-PetscCosScalar(thetatf)); /* df_dthetaf */ 165c4762a1bSJed Brown values[3] = Vmt*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); /* df_dVmf */ 166c4762a1bSJed Brown 167c4762a1bSJed Brown ierr = MatSetValues(J,1,row,4,col,values,ADD_VALUES);CHKERRQ(ierr); 168c4762a1bSJed Brown } 169c4762a1bSJed Brown if (bust->ide != PV_BUS && bust->ide != REF_BUS) { 170c4762a1bSJed Brown row[0] = goffsetto+1; 171c4762a1bSJed Brown col[0] = goffsetto; col[1] = goffsetto+1; col[2] = goffsetfrom; col[3] = goffsetfrom+1; 172c4762a1bSJed Brown /* farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); */ 173c4762a1bSJed Brown values[0] = Vmt*Vmf*(Btf*PetscSinScalar(thetatf) + Gtf*PetscCosScalar(thetatf)); 174c4762a1bSJed Brown values[1] = -2.0*Btt*Vmt + Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); 175c4762a1bSJed Brown values[2] = Vmt*Vmf*(-Btf*PetscSinScalar(thetatf) + Gtf*-PetscCosScalar(thetatf)); 176c4762a1bSJed Brown values[3] = Vmt*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); 177c4762a1bSJed Brown 178c4762a1bSJed Brown ierr = MatSetValues(J,1,row,4,col,values,ADD_VALUES);CHKERRQ(ierr); 179c4762a1bSJed Brown } 180c4762a1bSJed Brown } 181c4762a1bSJed Brown } 182c4762a1bSJed Brown if (!ghostvtex && bus->ide == PV_BUS) { 183c4762a1bSJed Brown row[0] = goffset+1; col[0] = goffset+1; values[0] = 1.0; 184c4762a1bSJed Brown if (user_power->jac_error) values[0] = 50.0; 185c4762a1bSJed Brown ierr = MatSetValues(J,1,row,1,col,values,ADD_VALUES);CHKERRQ(ierr); 186c4762a1bSJed Brown } 187c4762a1bSJed Brown } 188c4762a1bSJed Brown } 189c4762a1bSJed Brown } 190c4762a1bSJed Brown 191c4762a1bSJed Brown ierr = VecRestoreArrayRead(localX,&xarr);CHKERRQ(ierr); 192c4762a1bSJed Brown PetscFunctionReturn(0); 193c4762a1bSJed Brown } 194c4762a1bSJed Brown 195c4762a1bSJed Brown PetscErrorCode FormJacobian_Power(SNES snes,Vec X, Mat J,Mat Jpre,void *appctx) 196c4762a1bSJed Brown { 197c4762a1bSJed Brown PetscErrorCode ierr; 198c4762a1bSJed Brown DM networkdm; 199c4762a1bSJed Brown Vec localX; 200c4762a1bSJed Brown PetscInt nv,ne; 201c4762a1bSJed Brown const PetscInt *vtx,*edges; 202c4762a1bSJed Brown 203c4762a1bSJed Brown PetscFunctionBegin; 204c4762a1bSJed Brown ierr = MatZeroEntries(J);CHKERRQ(ierr); 205c4762a1bSJed Brown 206c4762a1bSJed Brown ierr = SNESGetDM(snes,&networkdm);CHKERRQ(ierr); 207c4762a1bSJed Brown ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr); 208c4762a1bSJed Brown 209c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 210c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 211c4762a1bSJed Brown 212*2bf73ac6SHong Zhang ierr = DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 213c4762a1bSJed Brown ierr = FormJacobian_Power_private(networkdm,localX,J,nv,ne,vtx,edges,appctx);CHKERRQ(ierr); 214c4762a1bSJed Brown 215c4762a1bSJed Brown ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr); 216c4762a1bSJed Brown 217c4762a1bSJed Brown ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 218c4762a1bSJed Brown ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 219c4762a1bSJed Brown PetscFunctionReturn(0); 220c4762a1bSJed Brown } 221c4762a1bSJed Brown 222c4762a1bSJed Brown PetscErrorCode FormFunction_Power(DM networkdm,Vec localX, Vec localF,PetscInt nv,PetscInt ne,const PetscInt* vtx,const PetscInt* edges,void* appctx) 223c4762a1bSJed Brown { 224c4762a1bSJed Brown PetscErrorCode ierr; 225c4762a1bSJed Brown UserCtx_Power *User=(UserCtx_Power*)appctx; 226c4762a1bSJed Brown PetscInt e,v,vfrom,vto; 227c4762a1bSJed Brown const PetscScalar *xarr; 228c4762a1bSJed Brown PetscScalar *farr; 229c4762a1bSJed Brown PetscInt offsetfrom,offsetto,offset,i,j,key,numComps; 230c4762a1bSJed Brown PetscScalar Vm; 231c4762a1bSJed Brown PetscScalar Sbase=User->Sbase; 232c4762a1bSJed Brown VERTEX_Power bus=NULL; 233c4762a1bSJed Brown GEN gen; 234c4762a1bSJed Brown LOAD load; 235c4762a1bSJed Brown PetscBool ghostvtex; 236c4762a1bSJed Brown void* component; 237c4762a1bSJed Brown 238c4762a1bSJed Brown PetscFunctionBegin; 239c4762a1bSJed Brown ierr = VecGetArrayRead(localX,&xarr);CHKERRQ(ierr); 240c4762a1bSJed Brown ierr = VecGetArray(localF,&farr);CHKERRQ(ierr); 241c4762a1bSJed Brown 242c4762a1bSJed Brown for (v=0; v<nv; v++) { 243c4762a1bSJed Brown ierr = DMNetworkIsGhostVertex(networkdm,vtx[v],&ghostvtex);CHKERRQ(ierr); 244c4762a1bSJed Brown ierr = DMNetworkGetNumComponents(networkdm,vtx[v],&numComps);CHKERRQ(ierr); 245*2bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vtx[v],ALL_COMPONENTS,&offset);CHKERRQ(ierr); 246c4762a1bSJed Brown 247c4762a1bSJed Brown for (j = 0; j < numComps; j++) { 248*2bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,vtx[v],j,&key,&component,NULL);CHKERRQ(ierr); 249c4762a1bSJed Brown if (key == User->compkey_bus) { 250c4762a1bSJed Brown PetscInt nconnedges; 251c4762a1bSJed Brown const PetscInt *connedges; 252c4762a1bSJed Brown 253c4762a1bSJed Brown bus = (VERTEX_Power)(component); 254c4762a1bSJed Brown /* Handle reference bus constrained dofs */ 255c4762a1bSJed Brown if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) { 256c4762a1bSJed Brown farr[offset] = xarr[offset] - bus->va*PETSC_PI/180.0; 257c4762a1bSJed Brown farr[offset+1] = xarr[offset+1] - bus->vm; 258c4762a1bSJed Brown break; 259c4762a1bSJed Brown } 260c4762a1bSJed Brown 261c4762a1bSJed Brown if (!ghostvtex) { 262c4762a1bSJed Brown Vm = xarr[offset+1]; 263c4762a1bSJed Brown 264c4762a1bSJed Brown /* Shunt injections */ 265c4762a1bSJed Brown farr[offset] += Vm*Vm*bus->gl/Sbase; 266c4762a1bSJed Brown if (bus->ide != PV_BUS) farr[offset+1] += -Vm*Vm*bus->bl/Sbase; 267c4762a1bSJed Brown } 268c4762a1bSJed Brown 269c4762a1bSJed Brown ierr = DMNetworkGetSupportingEdges(networkdm,vtx[v],&nconnedges,&connedges);CHKERRQ(ierr); 270c4762a1bSJed Brown for (i=0; i < nconnedges; i++) { 271c4762a1bSJed Brown EDGE_Power branch; 272c4762a1bSJed Brown PetscInt keye; 273c4762a1bSJed Brown PetscScalar Gff,Bff,Gft,Bft,Gtf,Btf,Gtt,Btt; 274c4762a1bSJed Brown const PetscInt *cone; 275c4762a1bSJed Brown PetscScalar Vmf,Vmt,thetaf,thetat,thetaft,thetatf; 276c4762a1bSJed Brown 277c4762a1bSJed Brown e = connedges[i]; 278*2bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,e,0,&keye,(void**)&branch,NULL);CHKERRQ(ierr); 279c4762a1bSJed Brown if (!branch->status) continue; 280c4762a1bSJed Brown Gff = branch->yff[0]; 281c4762a1bSJed Brown Bff = branch->yff[1]; 282c4762a1bSJed Brown Gft = branch->yft[0]; 283c4762a1bSJed Brown Bft = branch->yft[1]; 284c4762a1bSJed Brown Gtf = branch->ytf[0]; 285c4762a1bSJed Brown Btf = branch->ytf[1]; 286c4762a1bSJed Brown Gtt = branch->ytt[0]; 287c4762a1bSJed Brown Btt = branch->ytt[1]; 288c4762a1bSJed Brown 289c4762a1bSJed Brown ierr = DMNetworkGetConnectedVertices(networkdm,e,&cone);CHKERRQ(ierr); 290c4762a1bSJed Brown vfrom = cone[0]; 291c4762a1bSJed Brown vto = cone[1]; 292c4762a1bSJed Brown 293*2bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&offsetfrom);CHKERRQ(ierr); 294*2bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vto,ALL_COMPONENTS,&offsetto);CHKERRQ(ierr); 295c4762a1bSJed Brown 296c4762a1bSJed Brown thetaf = xarr[offsetfrom]; 297c4762a1bSJed Brown Vmf = xarr[offsetfrom+1]; 298c4762a1bSJed Brown thetat = xarr[offsetto]; 299c4762a1bSJed Brown Vmt = xarr[offsetto+1]; 300c4762a1bSJed Brown thetaft = thetaf - thetat; 301c4762a1bSJed Brown thetatf = thetat - thetaf; 302c4762a1bSJed Brown 303c4762a1bSJed Brown if (vfrom == vtx[v]) { 304c4762a1bSJed Brown farr[offsetfrom] += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); 305c4762a1bSJed Brown farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); 306c4762a1bSJed Brown } else { 307c4762a1bSJed Brown farr[offsetto] += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); 308c4762a1bSJed Brown farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); 309c4762a1bSJed Brown } 310c4762a1bSJed Brown } 311c4762a1bSJed Brown } else if (key == User->compkey_gen) { 312c4762a1bSJed Brown if (!ghostvtex) { 313c4762a1bSJed Brown gen = (GEN)(component); 314c4762a1bSJed Brown if (!gen->status) continue; 315c4762a1bSJed Brown farr[offset] += -gen->pg/Sbase; 316c4762a1bSJed Brown farr[offset+1] += -gen->qg/Sbase; 317c4762a1bSJed Brown } 318c4762a1bSJed Brown } else if (key == User->compkey_load) { 319c4762a1bSJed Brown if (!ghostvtex) { 320c4762a1bSJed Brown load = (LOAD)(component); 321c4762a1bSJed Brown farr[offset] += load->pl/Sbase; 322c4762a1bSJed Brown farr[offset+1] += load->ql/Sbase; 323c4762a1bSJed Brown } 324c4762a1bSJed Brown } 325c4762a1bSJed Brown } 326c4762a1bSJed Brown if (bus && bus->ide == PV_BUS) { 327c4762a1bSJed Brown farr[offset+1] = xarr[offset+1] - bus->vm; 328c4762a1bSJed Brown } 329c4762a1bSJed Brown } 330c4762a1bSJed Brown ierr = VecRestoreArrayRead(localX,&xarr);CHKERRQ(ierr); 331c4762a1bSJed Brown ierr = VecRestoreArray(localF,&farr);CHKERRQ(ierr); 332c4762a1bSJed Brown PetscFunctionReturn(0); 333c4762a1bSJed Brown } 334c4762a1bSJed Brown 335c4762a1bSJed Brown PetscErrorCode SetInitialGuess_Power(DM networkdm,Vec localX,PetscInt nv,PetscInt ne,const PetscInt *vtx,const PetscInt *edges,void* appctx) 336c4762a1bSJed Brown { 337c4762a1bSJed Brown PetscErrorCode ierr; 338c4762a1bSJed Brown VERTEX_Power bus; 339c4762a1bSJed Brown PetscInt i; 340c4762a1bSJed Brown GEN gen; 341*2bf73ac6SHong Zhang PetscBool ghostvtex,sharedv; 342c4762a1bSJed Brown PetscScalar *xarr; 343c4762a1bSJed Brown PetscInt key,numComps,j,offset; 344c4762a1bSJed Brown void* component; 345c4762a1bSJed Brown PetscMPIInt rank; 346c4762a1bSJed Brown MPI_Comm comm; 347c4762a1bSJed Brown UserCtx_Power *User=(UserCtx_Power*)appctx; 348c4762a1bSJed Brown 349c4762a1bSJed Brown PetscFunctionBegin; 350c4762a1bSJed Brown ierr = PetscObjectGetComm((PetscObject)networkdm,&comm);CHKERRQ(ierr); 351ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 352c4762a1bSJed Brown ierr = VecGetArray(localX,&xarr);CHKERRQ(ierr); 353c4762a1bSJed Brown for (i = 0; i < nv; i++) { 354c4762a1bSJed Brown ierr = DMNetworkIsGhostVertex(networkdm,vtx[i],&ghostvtex);CHKERRQ(ierr); 355*2bf73ac6SHong Zhang ierr = DMNetworkIsSharedVertex(networkdm,vtx[i],&sharedv);CHKERRQ(ierr); 356*2bf73ac6SHong Zhang if (ghostvtex ||sharedv) continue; 357c4762a1bSJed Brown 358*2bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vtx[i],ALL_COMPONENTS,&offset);CHKERRQ(ierr); 359c4762a1bSJed Brown ierr = DMNetworkGetNumComponents(networkdm,vtx[i],&numComps);CHKERRQ(ierr); 360c4762a1bSJed Brown for (j=0; j < numComps; j++) { 361*2bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,vtx[i],j,&key,&component,NULL);CHKERRQ(ierr); 362c4762a1bSJed Brown if (key == User->compkey_bus) { 363c4762a1bSJed Brown bus = (VERTEX_Power)(component); 364c4762a1bSJed Brown xarr[offset] = bus->va*PETSC_PI/180.0; 365c4762a1bSJed Brown xarr[offset+1] = bus->vm; 366c4762a1bSJed Brown } else if (key == User->compkey_gen) { 367c4762a1bSJed Brown gen = (GEN)(component); 368c4762a1bSJed Brown if (!gen->status) continue; 369c4762a1bSJed Brown xarr[offset+1] = gen->vs; 370c4762a1bSJed Brown break; 371c4762a1bSJed Brown } 372c4762a1bSJed Brown } 373c4762a1bSJed Brown } 374c4762a1bSJed Brown ierr = VecRestoreArray(localX,&xarr);CHKERRQ(ierr); 375c4762a1bSJed Brown PetscFunctionReturn(0); 376c4762a1bSJed Brown } 377