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 /* T 9c4762a1bSJed Brown Concepts: DMNetwork 10c4762a1bSJed Brown Concepts: PETSc SNES solver 11c4762a1bSJed Brown */ 12c4762a1bSJed Brown 13c4762a1bSJed Brown #include "power.h" 14c4762a1bSJed Brown #include <petscdmnetwork.h> 15c4762a1bSJed Brown 16c4762a1bSJed Brown PetscErrorCode FormFunction_Subnet(DM networkdm,Vec localX, Vec localF,PetscInt nv,PetscInt ne,const PetscInt* vtx,const PetscInt* edges,void* appctx) 17c4762a1bSJed Brown { 18c4762a1bSJed Brown PetscErrorCode ierr; 19c4762a1bSJed Brown UserCtx_Power *User = (UserCtx_Power*)appctx; 20c4762a1bSJed Brown PetscInt e,v,vfrom,vto; 21c4762a1bSJed Brown const PetscScalar *xarr; 22c4762a1bSJed Brown PetscScalar *farr; 23c4762a1bSJed Brown PetscInt offsetfrom,offsetto,offset; 24c4762a1bSJed Brown 25c4762a1bSJed Brown PetscFunctionBegin; 26c4762a1bSJed Brown ierr = VecGetArrayRead(localX,&xarr);CHKERRQ(ierr); 27c4762a1bSJed Brown ierr = VecGetArray(localF,&farr);CHKERRQ(ierr); 28c4762a1bSJed Brown 29c4762a1bSJed Brown for (v=0; v<nv; v++) { 30c4762a1bSJed Brown PetscInt i,j,key; 31c4762a1bSJed Brown PetscScalar Vm; 32c4762a1bSJed Brown PetscScalar Sbase = User->Sbase; 33c4762a1bSJed Brown VERTEX_Power bus = NULL; 34c4762a1bSJed Brown GEN gen; 35c4762a1bSJed Brown LOAD load; 36c4762a1bSJed Brown PetscBool ghostvtex; 37c4762a1bSJed Brown PetscInt numComps; 38c4762a1bSJed Brown void* component; 39c4762a1bSJed Brown 40c4762a1bSJed Brown ierr = DMNetworkIsGhostVertex(networkdm,vtx[v],&ghostvtex);CHKERRQ(ierr); 41c4762a1bSJed Brown ierr = DMNetworkGetNumComponents(networkdm,vtx[v],&numComps);CHKERRQ(ierr); 422bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vtx[v],ALL_COMPONENTS,&offset);CHKERRQ(ierr); 43c4762a1bSJed Brown for (j = 0; j < numComps; j++) { 442bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,vtx[v],j,&key,&component,NULL);CHKERRQ(ierr); 45c4762a1bSJed Brown if (key == 1) { 46c4762a1bSJed Brown PetscInt nconnedges; 47c4762a1bSJed Brown const PetscInt *connedges; 48c4762a1bSJed Brown 49c4762a1bSJed Brown bus = (VERTEX_Power)(component); 50c4762a1bSJed Brown /* Handle reference bus constrained dofs */ 51c4762a1bSJed Brown if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) { 52c4762a1bSJed Brown farr[offset] = xarr[offset] - bus->va*PETSC_PI/180.0; 53c4762a1bSJed Brown farr[offset+1] = xarr[offset+1] - bus->vm; 54c4762a1bSJed Brown break; 55c4762a1bSJed Brown } 56c4762a1bSJed Brown 57c4762a1bSJed Brown if (!ghostvtex) { 58c4762a1bSJed Brown Vm = xarr[offset+1]; 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Shunt injections */ 61c4762a1bSJed Brown farr[offset] += Vm*Vm*bus->gl/Sbase; 62c4762a1bSJed Brown if (bus->ide != PV_BUS) farr[offset+1] += -Vm*Vm*bus->bl/Sbase; 63c4762a1bSJed Brown } 64c4762a1bSJed Brown 65c4762a1bSJed Brown ierr = DMNetworkGetSupportingEdges(networkdm,vtx[v],&nconnedges,&connedges);CHKERRQ(ierr); 66c4762a1bSJed Brown for (i=0; i < nconnedges; i++) { 67c4762a1bSJed Brown EDGE_Power branch; 68c4762a1bSJed Brown PetscInt keye; 69c4762a1bSJed Brown PetscScalar Gff,Bff,Gft,Bft,Gtf,Btf,Gtt,Btt; 70c4762a1bSJed Brown const PetscInt *cone; 71c4762a1bSJed Brown PetscScalar Vmf,Vmt,thetaf,thetat,thetaft,thetatf; 72c4762a1bSJed Brown 73c4762a1bSJed Brown e = connedges[i]; 742bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,e,0,&keye,(void**)&branch,NULL);CHKERRQ(ierr); 75c4762a1bSJed Brown if (!branch->status) continue; 76c4762a1bSJed Brown Gff = branch->yff[0]; 77c4762a1bSJed Brown Bff = branch->yff[1]; 78c4762a1bSJed Brown Gft = branch->yft[0]; 79c4762a1bSJed Brown Bft = branch->yft[1]; 80c4762a1bSJed Brown Gtf = branch->ytf[0]; 81c4762a1bSJed Brown Btf = branch->ytf[1]; 82c4762a1bSJed Brown Gtt = branch->ytt[0]; 83c4762a1bSJed Brown Btt = branch->ytt[1]; 84c4762a1bSJed Brown 85c4762a1bSJed Brown ierr = DMNetworkGetConnectedVertices(networkdm,e,&cone);CHKERRQ(ierr); 86c4762a1bSJed Brown vfrom = cone[0]; 87c4762a1bSJed Brown vto = cone[1]; 88c4762a1bSJed Brown 892bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&offsetfrom);CHKERRQ(ierr); 902bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vto,ALL_COMPONENTS,&offsetto);CHKERRQ(ierr); 91c4762a1bSJed Brown 92c4762a1bSJed Brown thetaf = xarr[offsetfrom]; 93c4762a1bSJed Brown Vmf = xarr[offsetfrom+1]; 94c4762a1bSJed Brown thetat = xarr[offsetto]; 95c4762a1bSJed Brown Vmt = xarr[offsetto+1]; 96c4762a1bSJed Brown thetaft = thetaf - thetat; 97c4762a1bSJed Brown thetatf = thetat - thetaf; 98c4762a1bSJed Brown 99c4762a1bSJed Brown if (vfrom == vtx[v]) { 100c4762a1bSJed Brown farr[offsetfrom] += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); 101c4762a1bSJed Brown farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); 102c4762a1bSJed Brown } else { 103c4762a1bSJed Brown farr[offsetto] += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); 104c4762a1bSJed Brown farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); 105c4762a1bSJed Brown } 106c4762a1bSJed Brown } 107c4762a1bSJed Brown } else if (key == 2) { 108c4762a1bSJed Brown if (!ghostvtex) { 109c4762a1bSJed Brown gen = (GEN)(component); 110c4762a1bSJed Brown if (!gen->status) continue; 111c4762a1bSJed Brown farr[offset] += -gen->pg/Sbase; 112c4762a1bSJed Brown farr[offset+1] += -gen->qg/Sbase; 113c4762a1bSJed Brown } 114c4762a1bSJed Brown } else if (key == 3) { 115c4762a1bSJed Brown if (!ghostvtex) { 116c4762a1bSJed Brown load = (LOAD)(component); 117c4762a1bSJed Brown farr[offset] += load->pl/Sbase; 118c4762a1bSJed Brown farr[offset+1] += load->ql/Sbase; 119c4762a1bSJed Brown } 120c4762a1bSJed Brown } 121c4762a1bSJed Brown } 122c4762a1bSJed Brown if (bus && bus->ide == PV_BUS) { 123c4762a1bSJed Brown farr[offset+1] = xarr[offset+1] - bus->vm; 124c4762a1bSJed Brown } 125c4762a1bSJed Brown } 126c4762a1bSJed Brown ierr = VecRestoreArrayRead(localX,&xarr);CHKERRQ(ierr); 127c4762a1bSJed Brown ierr = VecRestoreArray(localF,&farr);CHKERRQ(ierr); 128c4762a1bSJed Brown PetscFunctionReturn(0); 129c4762a1bSJed Brown } 130c4762a1bSJed Brown 131c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X, Vec F,void *appctx) 132c4762a1bSJed Brown { 133c4762a1bSJed Brown PetscErrorCode ierr; 134c4762a1bSJed Brown DM networkdm; 135c4762a1bSJed Brown Vec localX,localF; 136c4762a1bSJed Brown PetscInt nv,ne; 137c4762a1bSJed Brown const PetscInt *vtx,*edges; 138c4762a1bSJed Brown 139c4762a1bSJed Brown PetscFunctionBegin; 140c4762a1bSJed Brown ierr = SNESGetDM(snes,&networkdm);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = DMGetLocalVector(networkdm,&localF);CHKERRQ(ierr); 143c4762a1bSJed Brown ierr = VecSet(F,0.0);CHKERRQ(ierr); 144c4762a1bSJed Brown 145c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 146c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 147c4762a1bSJed Brown 148c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(networkdm,F,INSERT_VALUES,localF);CHKERRQ(ierr); 149c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(networkdm,F,INSERT_VALUES,localF);CHKERRQ(ierr); 150c4762a1bSJed Brown 151c4762a1bSJed Brown /* Form Function for first subnetwork */ 1522bf73ac6SHong Zhang ierr = DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 153c4762a1bSJed Brown ierr = FormFunction_Subnet(networkdm,localX,localF,nv,ne,vtx,edges,appctx);CHKERRQ(ierr); 154c4762a1bSJed Brown 155c4762a1bSJed Brown /* Form Function for second subnetwork */ 1562bf73ac6SHong Zhang ierr = DMNetworkGetSubnetwork(networkdm,1,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 157c4762a1bSJed Brown ierr = FormFunction_Subnet(networkdm,localX,localF,nv,ne,vtx,edges,appctx);CHKERRQ(ierr); 158c4762a1bSJed Brown 159c4762a1bSJed Brown ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr); 160c4762a1bSJed Brown 161c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(networkdm,localF,ADD_VALUES,F);CHKERRQ(ierr); 162c4762a1bSJed Brown ierr = DMLocalToGlobalEnd(networkdm,localF,ADD_VALUES,F);CHKERRQ(ierr); 163c4762a1bSJed Brown ierr = DMRestoreLocalVector(networkdm,&localF);CHKERRQ(ierr); 164c4762a1bSJed Brown PetscFunctionReturn(0); 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167c4762a1bSJed Brown PetscErrorCode FormJacobian_Subnet(DM networkdm,Vec localX, Mat J, Mat Jpre, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx) 168c4762a1bSJed Brown { 169c4762a1bSJed Brown PetscErrorCode ierr; 170c4762a1bSJed Brown UserCtx_Power *User=(UserCtx_Power*)appctx; 171c4762a1bSJed Brown PetscInt e,v,vfrom,vto; 172c4762a1bSJed Brown const PetscScalar *xarr; 173c4762a1bSJed Brown PetscInt offsetfrom,offsetto,goffsetfrom,goffsetto; 174c4762a1bSJed Brown PetscInt row[2],col[8]; 175c4762a1bSJed Brown PetscScalar values[8]; 176c4762a1bSJed Brown 177c4762a1bSJed Brown PetscFunctionBegin; 178c4762a1bSJed Brown ierr = VecGetArrayRead(localX,&xarr);CHKERRQ(ierr); 179c4762a1bSJed Brown 180c4762a1bSJed Brown for (v=0; v<nv; v++) { 181c4762a1bSJed Brown PetscInt i,j,key; 182c4762a1bSJed Brown PetscInt offset,goffset; 183c4762a1bSJed Brown PetscScalar Vm; 184c4762a1bSJed Brown PetscScalar Sbase=User->Sbase; 185c4762a1bSJed Brown VERTEX_Power bus; 186c4762a1bSJed Brown PetscBool ghostvtex; 187c4762a1bSJed Brown PetscInt numComps; 188c4762a1bSJed Brown void* component; 189c4762a1bSJed Brown 190c4762a1bSJed Brown ierr = DMNetworkIsGhostVertex(networkdm,vtx[v],&ghostvtex);CHKERRQ(ierr); 191c4762a1bSJed Brown ierr = DMNetworkGetNumComponents(networkdm,vtx[v],&numComps);CHKERRQ(ierr); 192c4762a1bSJed Brown for (j = 0; j < numComps; j++) { 1932bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vtx[v],ALL_COMPONENTS,&offset);CHKERRQ(ierr); 1942bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVecOffset(networkdm,vtx[v],ALL_COMPONENTS,&goffset);CHKERRQ(ierr); 1952bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,vtx[v],j,&key,&component,NULL);CHKERRQ(ierr); 196c4762a1bSJed Brown if (key == 1) { 197c4762a1bSJed Brown PetscInt nconnedges; 198c4762a1bSJed Brown const PetscInt *connedges; 199c4762a1bSJed Brown 200c4762a1bSJed Brown bus = (VERTEX_Power)(component); 201c4762a1bSJed Brown if (!ghostvtex) { 202c4762a1bSJed Brown /* Handle reference bus constrained dofs */ 203c4762a1bSJed Brown if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) { 204c4762a1bSJed Brown row[0] = goffset; row[1] = goffset+1; 205c4762a1bSJed Brown col[0] = goffset; col[1] = goffset+1; col[2] = goffset; col[3] = goffset+1; 206c4762a1bSJed Brown values[0] = 1.0; values[1] = 0.0; values[2] = 0.0; values[3] = 1.0; 207c4762a1bSJed Brown ierr = MatSetValues(J,2,row,2,col,values,ADD_VALUES);CHKERRQ(ierr); 208c4762a1bSJed Brown break; 209c4762a1bSJed Brown } 210c4762a1bSJed Brown 211c4762a1bSJed Brown Vm = xarr[offset+1]; 212c4762a1bSJed Brown 213c4762a1bSJed Brown /* Shunt injections */ 214c4762a1bSJed Brown row[0] = goffset; row[1] = goffset+1; 215c4762a1bSJed Brown col[0] = goffset; col[1] = goffset+1; 216c4762a1bSJed Brown values[0] = values[1] = values[2] = values[3] = 0.0; 217c4762a1bSJed Brown if (bus->ide != PV_BUS) { 218c4762a1bSJed Brown values[1] = 2.0*Vm*bus->gl/Sbase; 219c4762a1bSJed Brown values[3] = -2.0*Vm*bus->bl/Sbase; 220c4762a1bSJed Brown } 221c4762a1bSJed Brown ierr = MatSetValues(J,2,row,2,col,values,ADD_VALUES);CHKERRQ(ierr); 222c4762a1bSJed Brown } 223c4762a1bSJed Brown 224c4762a1bSJed Brown ierr = DMNetworkGetSupportingEdges(networkdm,vtx[v],&nconnedges,&connedges);CHKERRQ(ierr); 225c4762a1bSJed Brown for (i=0; i < nconnedges; i++) { 226c4762a1bSJed Brown EDGE_Power branch; 227c4762a1bSJed Brown VERTEX_Power busf,bust; 228c4762a1bSJed Brown PetscInt keyf,keyt; 229c4762a1bSJed Brown PetscScalar Gff,Bff,Gft,Bft,Gtf,Btf,Gtt,Btt; 230c4762a1bSJed Brown const PetscInt *cone; 231c4762a1bSJed Brown PetscScalar Vmf,Vmt,thetaf,thetat,thetaft,thetatf; 232c4762a1bSJed Brown 233c4762a1bSJed Brown e = connedges[i]; 2342bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,e,0,&key,(void**)&branch,NULL);CHKERRQ(ierr); 235c4762a1bSJed Brown if (!branch->status) continue; 236c4762a1bSJed Brown 237c4762a1bSJed Brown Gff = branch->yff[0]; 238c4762a1bSJed Brown Bff = branch->yff[1]; 239c4762a1bSJed Brown Gft = branch->yft[0]; 240c4762a1bSJed Brown Bft = branch->yft[1]; 241c4762a1bSJed Brown Gtf = branch->ytf[0]; 242c4762a1bSJed Brown Btf = branch->ytf[1]; 243c4762a1bSJed Brown Gtt = branch->ytt[0]; 244c4762a1bSJed Brown Btt = branch->ytt[1]; 245c4762a1bSJed Brown 246c4762a1bSJed Brown ierr = DMNetworkGetConnectedVertices(networkdm,e,&cone);CHKERRQ(ierr); 247c4762a1bSJed Brown vfrom = cone[0]; 248c4762a1bSJed Brown vto = cone[1]; 249c4762a1bSJed Brown 2502bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&offsetfrom);CHKERRQ(ierr); 2512bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vto,ALL_COMPONENTS,&offsetto);CHKERRQ(ierr); 2522bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&goffsetfrom);CHKERRQ(ierr); 2532bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVecOffset(networkdm,vto,ALL_COMPONENTS,&goffsetto);CHKERRQ(ierr); 254c4762a1bSJed Brown 255c4762a1bSJed Brown if (goffsetto < 0) goffsetto = -goffsetto - 1; 256c4762a1bSJed Brown 257c4762a1bSJed Brown thetaf = xarr[offsetfrom]; 258c4762a1bSJed Brown Vmf = xarr[offsetfrom+1]; 259c4762a1bSJed Brown thetat = xarr[offsetto]; 260c4762a1bSJed Brown Vmt = xarr[offsetto+1]; 261c4762a1bSJed Brown thetaft = thetaf - thetat; 262c4762a1bSJed Brown thetatf = thetat - thetaf; 263c4762a1bSJed Brown 2642bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,vfrom,0,&keyf,(void**)&busf,NULL);CHKERRQ(ierr); 2652bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,vto,0,&keyt,(void**)&bust,NULL);CHKERRQ(ierr); 266c4762a1bSJed Brown 267c4762a1bSJed Brown if (vfrom == vtx[v]) { 268c4762a1bSJed Brown if (busf->ide != REF_BUS) { 269c4762a1bSJed Brown /* farr[offsetfrom] += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); */ 270c4762a1bSJed Brown row[0] = goffsetfrom; 271c4762a1bSJed Brown col[0] = goffsetfrom; col[1] = goffsetfrom+1; col[2] = goffsetto; col[3] = goffsetto+1; 272c4762a1bSJed Brown values[0] = Vmf*Vmt*(Gft*-PetscSinScalar(thetaft) + Bft*PetscCosScalar(thetaft)); /* df_dthetaf */ 273c4762a1bSJed Brown values[1] = 2.0*Gff*Vmf + Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); /* df_dVmf */ 274c4762a1bSJed Brown values[2] = Vmf*Vmt*(Gft*PetscSinScalar(thetaft) + Bft*-PetscCosScalar(thetaft)); /* df_dthetat */ 275c4762a1bSJed Brown values[3] = Vmf*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); /* df_dVmt */ 276c4762a1bSJed Brown 277c4762a1bSJed Brown ierr = MatSetValues(J,1,row,4,col,values,ADD_VALUES);CHKERRQ(ierr); 278c4762a1bSJed Brown } 279c4762a1bSJed Brown if (busf->ide != PV_BUS && busf->ide != REF_BUS) { 280c4762a1bSJed Brown row[0] = goffsetfrom+1; 281c4762a1bSJed Brown col[0] = goffsetfrom; col[1] = goffsetfrom+1; col[2] = goffsetto; col[3] = goffsetto+1; 282c4762a1bSJed Brown /* farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); */ 283c4762a1bSJed Brown values[0] = Vmf*Vmt*(Bft*PetscSinScalar(thetaft) + Gft*PetscCosScalar(thetaft)); 284c4762a1bSJed Brown values[1] = -2.0*Bff*Vmf + Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); 285c4762a1bSJed Brown values[2] = Vmf*Vmt*(-Bft*PetscSinScalar(thetaft) + Gft*-PetscCosScalar(thetaft)); 286c4762a1bSJed Brown values[3] = Vmf*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); 287c4762a1bSJed Brown 288c4762a1bSJed Brown ierr = MatSetValues(J,1,row,4,col,values,ADD_VALUES);CHKERRQ(ierr); 289c4762a1bSJed Brown } 290c4762a1bSJed Brown } else { 291c4762a1bSJed Brown if (bust->ide != REF_BUS) { 292c4762a1bSJed Brown row[0] = goffsetto; 293c4762a1bSJed Brown col[0] = goffsetto; col[1] = goffsetto+1; col[2] = goffsetfrom; col[3] = goffsetfrom+1; 294c4762a1bSJed Brown /* farr[offsetto] += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); */ 295c4762a1bSJed Brown values[0] = Vmt*Vmf*(Gtf*-PetscSinScalar(thetatf) + Btf*PetscCosScalar(thetaft)); /* df_dthetat */ 296c4762a1bSJed Brown values[1] = 2.0*Gtt*Vmt + Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); /* df_dVmt */ 297c4762a1bSJed Brown values[2] = Vmt*Vmf*(Gtf*PetscSinScalar(thetatf) + Btf*-PetscCosScalar(thetatf)); /* df_dthetaf */ 298c4762a1bSJed Brown values[3] = Vmt*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); /* df_dVmf */ 299c4762a1bSJed Brown 300c4762a1bSJed Brown ierr = MatSetValues(J,1,row,4,col,values,ADD_VALUES);CHKERRQ(ierr); 301c4762a1bSJed Brown } 302c4762a1bSJed Brown if (bust->ide != PV_BUS && bust->ide != REF_BUS) { 303c4762a1bSJed Brown row[0] = goffsetto+1; 304c4762a1bSJed Brown col[0] = goffsetto; col[1] = goffsetto+1; col[2] = goffsetfrom; col[3] = goffsetfrom+1; 305c4762a1bSJed Brown /* farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); */ 306c4762a1bSJed Brown values[0] = Vmt*Vmf*(Btf*PetscSinScalar(thetatf) + Gtf*PetscCosScalar(thetatf)); 307c4762a1bSJed Brown values[1] = -2.0*Btt*Vmt + Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); 308c4762a1bSJed Brown values[2] = Vmt*Vmf*(-Btf*PetscSinScalar(thetatf) + Gtf*-PetscCosScalar(thetatf)); 309c4762a1bSJed Brown values[3] = Vmt*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); 310c4762a1bSJed Brown 311c4762a1bSJed Brown ierr = MatSetValues(J,1,row,4,col,values,ADD_VALUES);CHKERRQ(ierr); 312c4762a1bSJed Brown } 313c4762a1bSJed Brown } 314c4762a1bSJed Brown } 315c4762a1bSJed Brown if (!ghostvtex && bus->ide == PV_BUS) { 316c4762a1bSJed Brown row[0] = goffset+1; col[0] = goffset+1; 317c4762a1bSJed Brown values[0] = 1.0; 318c4762a1bSJed Brown ierr = MatSetValues(J,1,row,1,col,values,ADD_VALUES);CHKERRQ(ierr); 319c4762a1bSJed Brown } 320c4762a1bSJed Brown } 321c4762a1bSJed Brown } 322c4762a1bSJed Brown } 323c4762a1bSJed Brown ierr = VecRestoreArrayRead(localX,&xarr);CHKERRQ(ierr); 324c4762a1bSJed Brown PetscFunctionReturn(0); 325c4762a1bSJed Brown } 326c4762a1bSJed Brown 327c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec X, Mat J,Mat Jpre,void *appctx) 328c4762a1bSJed Brown { 329c4762a1bSJed Brown PetscErrorCode ierr; 330c4762a1bSJed Brown DM networkdm; 331c4762a1bSJed Brown Vec localX; 332c4762a1bSJed Brown PetscInt ne,nv; 333c4762a1bSJed Brown const PetscInt *vtx,*edges; 334c4762a1bSJed Brown 335c4762a1bSJed Brown PetscFunctionBegin; 336c4762a1bSJed Brown ierr = MatZeroEntries(J);CHKERRQ(ierr); 337c4762a1bSJed Brown 338c4762a1bSJed Brown ierr = SNESGetDM(snes,&networkdm);CHKERRQ(ierr); 339c4762a1bSJed Brown ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr); 340c4762a1bSJed Brown 341c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 342c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 343c4762a1bSJed Brown 344c4762a1bSJed Brown /* Form Jacobian for first subnetwork */ 3452bf73ac6SHong Zhang ierr = DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 346c4762a1bSJed Brown ierr = FormJacobian_Subnet(networkdm,localX,J,Jpre,nv,ne,vtx,edges,appctx);CHKERRQ(ierr); 347c4762a1bSJed Brown 348c4762a1bSJed Brown /* Form Jacobian for second subnetwork */ 3492bf73ac6SHong Zhang ierr = DMNetworkGetSubnetwork(networkdm,1,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 350c4762a1bSJed Brown ierr = FormJacobian_Subnet(networkdm,localX,J,Jpre,nv,ne,vtx,edges,appctx);CHKERRQ(ierr); 351c4762a1bSJed Brown 352c4762a1bSJed Brown ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr); 353c4762a1bSJed Brown 354c4762a1bSJed Brown ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 355c4762a1bSJed Brown ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 356c4762a1bSJed Brown PetscFunctionReturn(0); 357c4762a1bSJed Brown } 358c4762a1bSJed Brown 359c4762a1bSJed Brown PetscErrorCode SetInitialValues_Subnet(DM networkdm,Vec localX,PetscInt nv,PetscInt ne, const PetscInt *vtx, const PetscInt *edges,void* appctx) 360c4762a1bSJed Brown { 361c4762a1bSJed Brown PetscErrorCode ierr; 362c4762a1bSJed Brown VERTEX_Power bus; 363c4762a1bSJed Brown PetscInt i; 364c4762a1bSJed Brown GEN gen; 365c4762a1bSJed Brown PetscBool ghostvtex; 366c4762a1bSJed Brown PetscScalar *xarr; 367c4762a1bSJed Brown PetscInt key,numComps,j,offset; 368c4762a1bSJed Brown void* component; 369c4762a1bSJed Brown 370c4762a1bSJed Brown PetscFunctionBegin; 371c4762a1bSJed Brown ierr = VecGetArray(localX,&xarr);CHKERRQ(ierr); 372c4762a1bSJed Brown for (i = 0; i < nv; i++) { 373c4762a1bSJed Brown ierr = DMNetworkIsGhostVertex(networkdm,vtx[i],&ghostvtex);CHKERRQ(ierr); 374c4762a1bSJed Brown if (ghostvtex) continue; 375c4762a1bSJed Brown 3762bf73ac6SHong Zhang ierr = DMNetworkGetLocalVecOffset(networkdm,vtx[i],ALL_COMPONENTS,&offset);CHKERRQ(ierr); 377c4762a1bSJed Brown ierr = DMNetworkGetNumComponents(networkdm,vtx[i],&numComps);CHKERRQ(ierr); 378c4762a1bSJed Brown for (j=0; j < numComps; j++) { 3792bf73ac6SHong Zhang ierr = DMNetworkGetComponent(networkdm,vtx[i],j,&key,&component,NULL);CHKERRQ(ierr); 380c4762a1bSJed Brown if (key == 1) { 381c4762a1bSJed Brown bus = (VERTEX_Power)(component); 382c4762a1bSJed Brown xarr[offset] = bus->va*PETSC_PI/180.0; 383c4762a1bSJed Brown xarr[offset+1] = bus->vm; 384c4762a1bSJed Brown } else if (key == 2) { 385c4762a1bSJed Brown gen = (GEN)(component); 386c4762a1bSJed Brown if (!gen->status) continue; 387c4762a1bSJed Brown xarr[offset+1] = gen->vs; 388c4762a1bSJed Brown break; 389c4762a1bSJed Brown } 390c4762a1bSJed Brown } 391c4762a1bSJed Brown } 392c4762a1bSJed Brown ierr = VecRestoreArray(localX,&xarr);CHKERRQ(ierr); 393c4762a1bSJed Brown PetscFunctionReturn(0); 394c4762a1bSJed Brown } 395c4762a1bSJed Brown 396c4762a1bSJed Brown PetscErrorCode SetInitialValues(DM networkdm, Vec X,void* appctx) 397c4762a1bSJed Brown { 398c4762a1bSJed Brown PetscErrorCode ierr; 399c4762a1bSJed Brown PetscInt nv,ne; 400c4762a1bSJed Brown const PetscInt *vtx,*edges; 401c4762a1bSJed Brown Vec localX; 402c4762a1bSJed Brown 403c4762a1bSJed Brown PetscFunctionBegin; 404c4762a1bSJed Brown ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr); 405c4762a1bSJed Brown 406c4762a1bSJed Brown ierr = VecSet(X,0.0);CHKERRQ(ierr); 407c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 408c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 409c4762a1bSJed Brown 410c4762a1bSJed Brown /* Set initial guess for first subnetwork */ 4112bf73ac6SHong Zhang ierr = DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 412c4762a1bSJed Brown ierr = SetInitialValues_Subnet(networkdm,localX,nv,ne,vtx,edges,appctx);CHKERRQ(ierr); 413c4762a1bSJed Brown 414c4762a1bSJed Brown /* Set initial guess for second subnetwork */ 4152bf73ac6SHong Zhang ierr = DMNetworkGetSubnetwork(networkdm,1,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 416c4762a1bSJed Brown ierr = SetInitialValues_Subnet(networkdm,localX,nv,ne,vtx,edges,appctx);CHKERRQ(ierr); 417c4762a1bSJed Brown 418c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(networkdm,localX,ADD_VALUES,X);CHKERRQ(ierr); 419c4762a1bSJed Brown ierr = DMLocalToGlobalEnd(networkdm,localX,ADD_VALUES,X);CHKERRQ(ierr); 420c4762a1bSJed Brown ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr); 421c4762a1bSJed Brown PetscFunctionReturn(0); 422c4762a1bSJed Brown } 423c4762a1bSJed Brown 424c4762a1bSJed Brown int main(int argc,char ** argv) 425c4762a1bSJed Brown { 426c4762a1bSJed Brown PetscErrorCode ierr; 427c4762a1bSJed Brown char pfdata_file[PETSC_MAX_PATH_LEN]="case9.m"; 428c4762a1bSJed Brown PFDATA *pfdata1,*pfdata2; 429c4762a1bSJed Brown PetscInt numEdges1=0,numVertices1=0,numEdges2=0,numVertices2=0; 4302bf73ac6SHong Zhang PetscInt *edgelist1 = NULL,*edgelist2 = NULL,componentkey[4]; 431c4762a1bSJed Brown DM networkdm; 432c4762a1bSJed Brown UserCtx_Power User; 433956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 434c4762a1bSJed Brown PetscLogStage stage1,stage2; 435956f8c0dSBarry Smith #endif 436c4762a1bSJed Brown PetscMPIInt rank; 4372bf73ac6SHong Zhang PetscInt nsubnet = 2,nv,ne,i,j,genj,loadj; 4382bf73ac6SHong Zhang const PetscInt *vtx,*edges; 439c4762a1bSJed Brown Vec X,F; 440c4762a1bSJed Brown Mat J; 441c4762a1bSJed Brown SNES snes; 442c4762a1bSJed Brown 443c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,"poweroptions",help);if (ierr) return ierr; 444ffc4695bSBarry Smith ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr); 445c4762a1bSJed Brown { 446c4762a1bSJed 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 */ 447c4762a1bSJed Brown /* this is an experiment to see how the analyzer reacts */ 448c4762a1bSJed Brown const PetscMPIInt crank = rank; 449c4762a1bSJed Brown 450c4762a1bSJed Brown /* Create an empty network object */ 451c4762a1bSJed Brown ierr = DMNetworkCreate(PETSC_COMM_WORLD,&networkdm);CHKERRQ(ierr); 452c4762a1bSJed Brown 453c4762a1bSJed Brown /* Register the components in the network */ 454c4762a1bSJed Brown ierr = DMNetworkRegisterComponent(networkdm,"branchstruct",sizeof(struct _p_EDGE_Power),&componentkey[0]);CHKERRQ(ierr); 455c4762a1bSJed Brown ierr = DMNetworkRegisterComponent(networkdm,"busstruct",sizeof(struct _p_VERTEX_Power),&componentkey[1]);CHKERRQ(ierr); 456c4762a1bSJed Brown ierr = DMNetworkRegisterComponent(networkdm,"genstruct",sizeof(struct _p_GEN),&componentkey[2]);CHKERRQ(ierr); 457c4762a1bSJed Brown ierr = DMNetworkRegisterComponent(networkdm,"loadstruct",sizeof(struct _p_LOAD),&componentkey[3]);CHKERRQ(ierr); 458c4762a1bSJed Brown 459c4762a1bSJed Brown ierr = PetscLogStageRegister("Read Data",&stage1);CHKERRQ(ierr); 460c4762a1bSJed Brown PetscLogStagePush(stage1); 461c4762a1bSJed Brown /* READ THE DATA */ 462c4762a1bSJed Brown if (!crank) { 463c4762a1bSJed Brown /* Only rank 0 reads the data */ 464589a23caSBarry Smith ierr = PetscOptionsGetString(NULL,NULL,"-pfdata",pfdata_file,sizeof(pfdata_file),NULL);CHKERRQ(ierr); 465c4762a1bSJed Brown /* HERE WE CREATE COPIES OF THE SAME NETWORK THAT WILL BE TREATED AS SUBNETWORKS */ 466c4762a1bSJed Brown 467c4762a1bSJed Brown /* READ DATA FOR THE FIRST SUBNETWORK */ 468c4762a1bSJed Brown ierr = PetscNew(&pfdata1);CHKERRQ(ierr); 469c4762a1bSJed Brown ierr = PFReadMatPowerData(pfdata1,pfdata_file);CHKERRQ(ierr); 470c4762a1bSJed Brown User.Sbase = pfdata1->sbase; 471c4762a1bSJed Brown 472c4762a1bSJed Brown numEdges1 = pfdata1->nbranch; 473c4762a1bSJed Brown numVertices1 = pfdata1->nbus; 474c4762a1bSJed Brown 475c4762a1bSJed Brown ierr = PetscMalloc1(2*numEdges1,&edgelist1);CHKERRQ(ierr); 476c4762a1bSJed Brown ierr = GetListofEdges_Power(pfdata1,edgelist1);CHKERRQ(ierr); 477c4762a1bSJed Brown 478c4762a1bSJed Brown /* READ DATA FOR THE SECOND SUBNETWORK */ 479c4762a1bSJed Brown ierr = PetscNew(&pfdata2);CHKERRQ(ierr); 480c4762a1bSJed Brown ierr = PFReadMatPowerData(pfdata2,pfdata_file);CHKERRQ(ierr); 481c4762a1bSJed Brown User.Sbase = pfdata2->sbase; 482c4762a1bSJed Brown 483c4762a1bSJed Brown numEdges2 = pfdata2->nbranch; 484c4762a1bSJed Brown numVertices2 = pfdata2->nbus; 485c4762a1bSJed Brown 486c4762a1bSJed Brown ierr = PetscMalloc1(2*numEdges2,&edgelist2);CHKERRQ(ierr); 487c4762a1bSJed Brown ierr = GetListofEdges_Power(pfdata2,edgelist2);CHKERRQ(ierr); 488c4762a1bSJed Brown } 489c4762a1bSJed Brown 490c4762a1bSJed Brown PetscLogStagePop(); 491ffc4695bSBarry Smith ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRMPI(ierr); 492c4762a1bSJed Brown ierr = PetscLogStageRegister("Create network",&stage2);CHKERRQ(ierr); 493c4762a1bSJed Brown PetscLogStagePush(stage2); 494c4762a1bSJed Brown 4952bf73ac6SHong Zhang /* Set number of nodes/edges and edge connectivity */ 4962bf73ac6SHong Zhang ierr = DMNetworkSetNumSubNetworks(networkdm,PETSC_DECIDE,nsubnet);CHKERRQ(ierr); 4972bf73ac6SHong Zhang ierr = DMNetworkAddSubnetwork(networkdm,"",numVertices1,numEdges1,edgelist1,NULL);CHKERRQ(ierr); 4982bf73ac6SHong Zhang ierr = DMNetworkAddSubnetwork(networkdm,"",numVertices2,numEdges2,edgelist2,NULL);CHKERRQ(ierr); 499c4762a1bSJed Brown 500c4762a1bSJed Brown /* Set up the network layout */ 501c4762a1bSJed Brown ierr = DMNetworkLayoutSetUp(networkdm);CHKERRQ(ierr); 502c4762a1bSJed Brown 503c4762a1bSJed Brown /* Add network components only process 0 has any data to add*/ 504c4762a1bSJed Brown if (!crank) { 505c4762a1bSJed Brown genj=0; loadj=0; 506c4762a1bSJed Brown 507c4762a1bSJed Brown /* ADD VARIABLES AND COMPONENTS FOR THE FIRST SUBNETWORK */ 5082bf73ac6SHong Zhang ierr = DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 509c4762a1bSJed Brown 510c4762a1bSJed Brown for (i = 0; i < ne; i++) { 5112bf73ac6SHong Zhang ierr = DMNetworkAddComponent(networkdm,edges[i],componentkey[0],&pfdata1->branch[i],0);CHKERRQ(ierr); 512c4762a1bSJed Brown } 513c4762a1bSJed Brown 514c4762a1bSJed Brown for (i = 0; i < nv; i++) { 5152bf73ac6SHong Zhang ierr = DMNetworkAddComponent(networkdm,vtx[i],componentkey[1],&pfdata1->bus[i],2);CHKERRQ(ierr); 516c4762a1bSJed Brown if (pfdata1->bus[i].ngen) { 517c4762a1bSJed Brown for (j = 0; j < pfdata1->bus[i].ngen; j++) { 5182bf73ac6SHong Zhang ierr = DMNetworkAddComponent(networkdm,vtx[i],componentkey[2],&pfdata1->gen[genj++],0);CHKERRQ(ierr); 519c4762a1bSJed Brown } 520c4762a1bSJed Brown } 521c4762a1bSJed Brown if (pfdata1->bus[i].nload) { 522c4762a1bSJed Brown for (j=0; j < pfdata1->bus[i].nload; j++) { 5232bf73ac6SHong Zhang ierr = DMNetworkAddComponent(networkdm,vtx[i],componentkey[3],&pfdata1->load[loadj++],0);CHKERRQ(ierr); 524c4762a1bSJed Brown } 525c4762a1bSJed Brown } 526c4762a1bSJed Brown } 527c4762a1bSJed Brown 528c4762a1bSJed Brown genj=0; loadj=0; 529c4762a1bSJed Brown 530c4762a1bSJed Brown /* ADD VARIABLES AND COMPONENTS FOR THE SECOND SUBNETWORK */ 5312bf73ac6SHong Zhang ierr = DMNetworkGetSubnetwork(networkdm,1,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 532c4762a1bSJed Brown 533c4762a1bSJed Brown for (i = 0; i < ne; i++) { 5342bf73ac6SHong Zhang ierr = DMNetworkAddComponent(networkdm,edges[i],componentkey[0],&pfdata2->branch[i],0);CHKERRQ(ierr); 535c4762a1bSJed Brown } 536c4762a1bSJed Brown 537c4762a1bSJed Brown for (i = 0; i < nv; i++) { 5382bf73ac6SHong Zhang ierr = DMNetworkAddComponent(networkdm,vtx[i],componentkey[1],&pfdata2->bus[i],2);CHKERRQ(ierr); 539c4762a1bSJed Brown if (pfdata2->bus[i].ngen) { 540c4762a1bSJed Brown for (j = 0; j < pfdata2->bus[i].ngen; j++) { 5412bf73ac6SHong Zhang ierr = DMNetworkAddComponent(networkdm,vtx[i],componentkey[2],&pfdata2->gen[genj++],0);CHKERRQ(ierr); 542c4762a1bSJed Brown } 543c4762a1bSJed Brown } 544c4762a1bSJed Brown if (pfdata2->bus[i].nload) { 545c4762a1bSJed Brown for (j=0; j < pfdata2->bus[i].nload; j++) { 5462bf73ac6SHong Zhang ierr = DMNetworkAddComponent(networkdm,vtx[i],componentkey[3],&pfdata2->load[loadj++],0);CHKERRQ(ierr); 547c4762a1bSJed Brown } 548c4762a1bSJed Brown } 549c4762a1bSJed Brown } 550c4762a1bSJed Brown } 551c4762a1bSJed Brown 552c4762a1bSJed Brown /* Set up DM for use */ 553c4762a1bSJed Brown ierr = DMSetUp(networkdm);CHKERRQ(ierr); 554c4762a1bSJed Brown 555c4762a1bSJed Brown if (!crank) { 556c4762a1bSJed Brown ierr = PetscFree(edgelist1);CHKERRQ(ierr); 557c4762a1bSJed Brown ierr = PetscFree(edgelist2);CHKERRQ(ierr); 558c4762a1bSJed Brown } 559c4762a1bSJed Brown 560c4762a1bSJed Brown if (!crank) { 561c4762a1bSJed Brown ierr = PetscFree(pfdata1->bus);CHKERRQ(ierr); 562c4762a1bSJed Brown ierr = PetscFree(pfdata1->gen);CHKERRQ(ierr); 563c4762a1bSJed Brown ierr = PetscFree(pfdata1->branch);CHKERRQ(ierr); 564c4762a1bSJed Brown ierr = PetscFree(pfdata1->load);CHKERRQ(ierr); 565c4762a1bSJed Brown ierr = PetscFree(pfdata1);CHKERRQ(ierr); 566c4762a1bSJed Brown 567c4762a1bSJed Brown ierr = PetscFree(pfdata2->bus);CHKERRQ(ierr); 568c4762a1bSJed Brown ierr = PetscFree(pfdata2->gen);CHKERRQ(ierr); 569c4762a1bSJed Brown ierr = PetscFree(pfdata2->branch);CHKERRQ(ierr); 570c4762a1bSJed Brown ierr = PetscFree(pfdata2->load);CHKERRQ(ierr); 571c4762a1bSJed Brown ierr = PetscFree(pfdata2);CHKERRQ(ierr); 572c4762a1bSJed Brown } 573c4762a1bSJed Brown 574c4762a1bSJed Brown /* Distribute networkdm to multiple processes */ 575c4762a1bSJed Brown ierr = DMNetworkDistribute(&networkdm,0);CHKERRQ(ierr); 576c4762a1bSJed Brown 577c4762a1bSJed Brown PetscLogStagePop(); 578c4762a1bSJed Brown 579c4762a1bSJed Brown /* Broadcast Sbase to all processors */ 580ffc4695bSBarry Smith ierr = MPI_Bcast(&User.Sbase,1,MPIU_SCALAR,0,PETSC_COMM_WORLD);CHKERRMPI(ierr); 581c4762a1bSJed Brown 582c4762a1bSJed Brown ierr = DMCreateGlobalVector(networkdm,&X);CHKERRQ(ierr); 583c4762a1bSJed Brown ierr = VecDuplicate(X,&F);CHKERRQ(ierr); 584c4762a1bSJed Brown 585c4762a1bSJed Brown ierr = DMCreateMatrix(networkdm,&J);CHKERRQ(ierr); 586c4762a1bSJed Brown ierr = MatSetOption(J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 587c4762a1bSJed Brown 588c4762a1bSJed Brown ierr = SetInitialValues(networkdm,X,&User);CHKERRQ(ierr); 589c4762a1bSJed Brown 590c4762a1bSJed Brown /* HOOK UP SOLVER */ 591c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 592c4762a1bSJed Brown ierr = SNESSetDM(snes,networkdm);CHKERRQ(ierr); 593c4762a1bSJed Brown ierr = SNESSetFunction(snes,F,FormFunction,&User);CHKERRQ(ierr); 594c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,FormJacobian,&User);CHKERRQ(ierr); 595c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 596c4762a1bSJed Brown 597c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,X);CHKERRQ(ierr); 598c4762a1bSJed Brown /* ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 599c4762a1bSJed Brown 600c4762a1bSJed Brown ierr = VecDestroy(&X);CHKERRQ(ierr); 601c4762a1bSJed Brown ierr = VecDestroy(&F);CHKERRQ(ierr); 602c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 603c4762a1bSJed Brown 604c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 605c4762a1bSJed Brown ierr = DMDestroy(&networkdm);CHKERRQ(ierr); 606c4762a1bSJed Brown } 607c4762a1bSJed Brown ierr = PetscFinalize(); 608c4762a1bSJed Brown return ierr; 609c4762a1bSJed Brown } 610c4762a1bSJed Brown 611c4762a1bSJed Brown /*TEST 612c4762a1bSJed Brown 613c4762a1bSJed Brown build: 614c4762a1bSJed Brown depends: PFReadData.c pffunctions.c 615*dfd57a17SPierre Jolivet requires: !complex double defined(PETSC_HAVE_ATTRIBUTEALIGNED) 616c4762a1bSJed Brown 617c4762a1bSJed Brown test: 618c4762a1bSJed Brown args: -snes_rtol 1.e-3 619c4762a1bSJed Brown localrunfiles: poweroptions case9.m 6202bf73ac6SHong Zhang output_file: output/power_1.out 621c4762a1bSJed Brown 622c4762a1bSJed Brown test: 623c4762a1bSJed Brown suffix: 2 624c4762a1bSJed Brown args: -snes_rtol 1.e-3 -petscpartitioner_type simple 625c4762a1bSJed Brown nsize: 4 626c4762a1bSJed Brown localrunfiles: poweroptions case9.m 6272bf73ac6SHong Zhang output_file: output/power_1.out 628c4762a1bSJed Brown 629c4762a1bSJed Brown TEST*/ 630