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 UserCtx_Power *User = (UserCtx_Power*)appctx; 19c4762a1bSJed Brown PetscInt e,v,vfrom,vto; 20c4762a1bSJed Brown const PetscScalar *xarr; 21c4762a1bSJed Brown PetscScalar *farr; 22c4762a1bSJed Brown PetscInt offsetfrom,offsetto,offset; 23c4762a1bSJed Brown 24c4762a1bSJed Brown PetscFunctionBegin; 255f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(localX,&xarr)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(localF,&farr)); 27c4762a1bSJed Brown 28c4762a1bSJed Brown for (v=0; v<nv; v++) { 29c4762a1bSJed Brown PetscInt i,j,key; 30c4762a1bSJed Brown PetscScalar Vm; 31c4762a1bSJed Brown PetscScalar Sbase = User->Sbase; 32c4762a1bSJed Brown VERTEX_Power bus = NULL; 33c4762a1bSJed Brown GEN gen; 34c4762a1bSJed Brown LOAD load; 35c4762a1bSJed Brown PetscBool ghostvtex; 36c4762a1bSJed Brown PetscInt numComps; 37c4762a1bSJed Brown void* component; 38c4762a1bSJed Brown 395f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkIsGhostVertex(networkdm,vtx[v],&ghostvtex)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetNumComponents(networkdm,vtx[v],&numComps)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetLocalVecOffset(networkdm,vtx[v],ALL_COMPONENTS,&offset)); 42c4762a1bSJed Brown for (j = 0; j < numComps; j++) { 435f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetComponent(networkdm,vtx[v],j,&key,&component,NULL)); 44c4762a1bSJed Brown if (key == 1) { 45c4762a1bSJed Brown PetscInt nconnedges; 46c4762a1bSJed Brown const PetscInt *connedges; 47c4762a1bSJed Brown 48c4762a1bSJed Brown bus = (VERTEX_Power)(component); 49c4762a1bSJed Brown /* Handle reference bus constrained dofs */ 50c4762a1bSJed Brown if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) { 51c4762a1bSJed Brown farr[offset] = xarr[offset] - bus->va*PETSC_PI/180.0; 52c4762a1bSJed Brown farr[offset+1] = xarr[offset+1] - bus->vm; 53c4762a1bSJed Brown break; 54c4762a1bSJed Brown } 55c4762a1bSJed Brown 56c4762a1bSJed Brown if (!ghostvtex) { 57c4762a1bSJed Brown Vm = xarr[offset+1]; 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* Shunt injections */ 60c4762a1bSJed Brown farr[offset] += Vm*Vm*bus->gl/Sbase; 61c4762a1bSJed Brown if (bus->ide != PV_BUS) farr[offset+1] += -Vm*Vm*bus->bl/Sbase; 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 645f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetSupportingEdges(networkdm,vtx[v],&nconnedges,&connedges)); 65c4762a1bSJed Brown for (i=0; i < nconnedges; i++) { 66c4762a1bSJed Brown EDGE_Power branch; 67c4762a1bSJed Brown PetscInt keye; 68c4762a1bSJed Brown PetscScalar Gff,Bff,Gft,Bft,Gtf,Btf,Gtt,Btt; 69c4762a1bSJed Brown const PetscInt *cone; 70c4762a1bSJed Brown PetscScalar Vmf,Vmt,thetaf,thetat,thetaft,thetatf; 71c4762a1bSJed Brown 72c4762a1bSJed Brown e = connedges[i]; 735f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetComponent(networkdm,e,0,&keye,(void**)&branch,NULL)); 74c4762a1bSJed Brown if (!branch->status) continue; 75c4762a1bSJed Brown Gff = branch->yff[0]; 76c4762a1bSJed Brown Bff = branch->yff[1]; 77c4762a1bSJed Brown Gft = branch->yft[0]; 78c4762a1bSJed Brown Bft = branch->yft[1]; 79c4762a1bSJed Brown Gtf = branch->ytf[0]; 80c4762a1bSJed Brown Btf = branch->ytf[1]; 81c4762a1bSJed Brown Gtt = branch->ytt[0]; 82c4762a1bSJed Brown Btt = branch->ytt[1]; 83c4762a1bSJed Brown 845f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetConnectedVertices(networkdm,e,&cone)); 85c4762a1bSJed Brown vfrom = cone[0]; 86c4762a1bSJed Brown vto = cone[1]; 87c4762a1bSJed Brown 885f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetLocalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&offsetfrom)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetLocalVecOffset(networkdm,vto,ALL_COMPONENTS,&offsetto)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown thetaf = xarr[offsetfrom]; 92c4762a1bSJed Brown Vmf = xarr[offsetfrom+1]; 93c4762a1bSJed Brown thetat = xarr[offsetto]; 94c4762a1bSJed Brown Vmt = xarr[offsetto+1]; 95c4762a1bSJed Brown thetaft = thetaf - thetat; 96c4762a1bSJed Brown thetatf = thetat - thetaf; 97c4762a1bSJed Brown 98c4762a1bSJed Brown if (vfrom == vtx[v]) { 99c4762a1bSJed Brown farr[offsetfrom] += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); 100c4762a1bSJed Brown farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); 101c4762a1bSJed Brown } else { 102c4762a1bSJed Brown farr[offsetto] += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); 103c4762a1bSJed Brown farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); 104c4762a1bSJed Brown } 105c4762a1bSJed Brown } 106c4762a1bSJed Brown } else if (key == 2) { 107c4762a1bSJed Brown if (!ghostvtex) { 108c4762a1bSJed Brown gen = (GEN)(component); 109c4762a1bSJed Brown if (!gen->status) continue; 110c4762a1bSJed Brown farr[offset] += -gen->pg/Sbase; 111c4762a1bSJed Brown farr[offset+1] += -gen->qg/Sbase; 112c4762a1bSJed Brown } 113c4762a1bSJed Brown } else if (key == 3) { 114c4762a1bSJed Brown if (!ghostvtex) { 115c4762a1bSJed Brown load = (LOAD)(component); 116c4762a1bSJed Brown farr[offset] += load->pl/Sbase; 117c4762a1bSJed Brown farr[offset+1] += load->ql/Sbase; 118c4762a1bSJed Brown } 119c4762a1bSJed Brown } 120c4762a1bSJed Brown } 121c4762a1bSJed Brown if (bus && bus->ide == PV_BUS) { 122c4762a1bSJed Brown farr[offset+1] = xarr[offset+1] - bus->vm; 123c4762a1bSJed Brown } 124c4762a1bSJed Brown } 1255f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(localX,&xarr)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(localF,&farr)); 127c4762a1bSJed Brown PetscFunctionReturn(0); 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X, Vec F,void *appctx) 131c4762a1bSJed Brown { 132c4762a1bSJed Brown DM networkdm; 133c4762a1bSJed Brown Vec localX,localF; 134c4762a1bSJed Brown PetscInt nv,ne; 135c4762a1bSJed Brown const PetscInt *vtx,*edges; 136c4762a1bSJed Brown 137c4762a1bSJed Brown PetscFunctionBegin; 1385f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetDM(snes,&networkdm)); 1395f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(networkdm,&localX)); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(networkdm,&localF)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(F,0.0)); 142c4762a1bSJed Brown 1435f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX)); 145c4762a1bSJed Brown 1465f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(networkdm,F,INSERT_VALUES,localF)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(networkdm,F,INSERT_VALUES,localF)); 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* Form Function for first subnetwork */ 1505f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(FormFunction_Subnet(networkdm,localX,localF,nv,ne,vtx,edges,appctx)); 152c4762a1bSJed Brown 153c4762a1bSJed Brown /* Form Function for second subnetwork */ 1545f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetSubnetwork(networkdm,1,&nv,&ne,&vtx,&edges)); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(FormFunction_Subnet(networkdm,localX,localF,nv,ne,vtx,edges,appctx)); 156c4762a1bSJed Brown 1575f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(networkdm,&localX)); 158c4762a1bSJed Brown 1595f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalBegin(networkdm,localF,ADD_VALUES,F)); 1605f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalEnd(networkdm,localF,ADD_VALUES,F)); 1615f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(networkdm,&localF)); 162c4762a1bSJed Brown PetscFunctionReturn(0); 163c4762a1bSJed Brown } 164c4762a1bSJed Brown 165c4762a1bSJed Brown PetscErrorCode FormJacobian_Subnet(DM networkdm,Vec localX, Mat J, Mat Jpre, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx) 166c4762a1bSJed Brown { 167c4762a1bSJed Brown UserCtx_Power *User=(UserCtx_Power*)appctx; 168c4762a1bSJed Brown PetscInt e,v,vfrom,vto; 169c4762a1bSJed Brown const PetscScalar *xarr; 170c4762a1bSJed Brown PetscInt offsetfrom,offsetto,goffsetfrom,goffsetto; 171c4762a1bSJed Brown PetscInt row[2],col[8]; 172c4762a1bSJed Brown PetscScalar values[8]; 173c4762a1bSJed Brown 174c4762a1bSJed Brown PetscFunctionBegin; 1755f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(localX,&xarr)); 176c4762a1bSJed Brown 177c4762a1bSJed Brown for (v=0; v<nv; v++) { 178c4762a1bSJed Brown PetscInt i,j,key; 179c4762a1bSJed Brown PetscInt offset,goffset; 180c4762a1bSJed Brown PetscScalar Vm; 181c4762a1bSJed Brown PetscScalar Sbase=User->Sbase; 182c4762a1bSJed Brown VERTEX_Power bus; 183c4762a1bSJed Brown PetscBool ghostvtex; 184c4762a1bSJed Brown PetscInt numComps; 185c4762a1bSJed Brown void* component; 186c4762a1bSJed Brown 1875f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkIsGhostVertex(networkdm,vtx[v],&ghostvtex)); 1885f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetNumComponents(networkdm,vtx[v],&numComps)); 189c4762a1bSJed Brown for (j = 0; j < numComps; j++) { 1905f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetLocalVecOffset(networkdm,vtx[v],ALL_COMPONENTS,&offset)); 1915f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetGlobalVecOffset(networkdm,vtx[v],ALL_COMPONENTS,&goffset)); 1925f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetComponent(networkdm,vtx[v],j,&key,&component,NULL)); 193c4762a1bSJed Brown if (key == 1) { 194c4762a1bSJed Brown PetscInt nconnedges; 195c4762a1bSJed Brown const PetscInt *connedges; 196c4762a1bSJed Brown 197c4762a1bSJed Brown bus = (VERTEX_Power)(component); 198c4762a1bSJed Brown if (!ghostvtex) { 199c4762a1bSJed Brown /* Handle reference bus constrained dofs */ 200c4762a1bSJed Brown if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) { 201c4762a1bSJed Brown row[0] = goffset; row[1] = goffset+1; 202c4762a1bSJed Brown col[0] = goffset; col[1] = goffset+1; col[2] = goffset; col[3] = goffset+1; 203c4762a1bSJed Brown values[0] = 1.0; values[1] = 0.0; values[2] = 0.0; values[3] = 1.0; 2045f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(J,2,row,2,col,values,ADD_VALUES)); 205c4762a1bSJed Brown break; 206c4762a1bSJed Brown } 207c4762a1bSJed Brown 208c4762a1bSJed Brown Vm = xarr[offset+1]; 209c4762a1bSJed Brown 210c4762a1bSJed Brown /* Shunt injections */ 211c4762a1bSJed Brown row[0] = goffset; row[1] = goffset+1; 212c4762a1bSJed Brown col[0] = goffset; col[1] = goffset+1; 213c4762a1bSJed Brown values[0] = values[1] = values[2] = values[3] = 0.0; 214c4762a1bSJed Brown if (bus->ide != PV_BUS) { 215c4762a1bSJed Brown values[1] = 2.0*Vm*bus->gl/Sbase; 216c4762a1bSJed Brown values[3] = -2.0*Vm*bus->bl/Sbase; 217c4762a1bSJed Brown } 2185f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(J,2,row,2,col,values,ADD_VALUES)); 219c4762a1bSJed Brown } 220c4762a1bSJed Brown 2215f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetSupportingEdges(networkdm,vtx[v],&nconnedges,&connedges)); 222c4762a1bSJed Brown for (i=0; i < nconnedges; i++) { 223c4762a1bSJed Brown EDGE_Power branch; 224c4762a1bSJed Brown VERTEX_Power busf,bust; 225c4762a1bSJed Brown PetscInt keyf,keyt; 226c4762a1bSJed Brown PetscScalar Gff,Bff,Gft,Bft,Gtf,Btf,Gtt,Btt; 227c4762a1bSJed Brown const PetscInt *cone; 228c4762a1bSJed Brown PetscScalar Vmf,Vmt,thetaf,thetat,thetaft,thetatf; 229c4762a1bSJed Brown 230c4762a1bSJed Brown e = connedges[i]; 2315f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetComponent(networkdm,e,0,&key,(void**)&branch,NULL)); 232c4762a1bSJed Brown if (!branch->status) continue; 233c4762a1bSJed Brown 234c4762a1bSJed Brown Gff = branch->yff[0]; 235c4762a1bSJed Brown Bff = branch->yff[1]; 236c4762a1bSJed Brown Gft = branch->yft[0]; 237c4762a1bSJed Brown Bft = branch->yft[1]; 238c4762a1bSJed Brown Gtf = branch->ytf[0]; 239c4762a1bSJed Brown Btf = branch->ytf[1]; 240c4762a1bSJed Brown Gtt = branch->ytt[0]; 241c4762a1bSJed Brown Btt = branch->ytt[1]; 242c4762a1bSJed Brown 2435f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetConnectedVertices(networkdm,e,&cone)); 244c4762a1bSJed Brown vfrom = cone[0]; 245c4762a1bSJed Brown vto = cone[1]; 246c4762a1bSJed Brown 2475f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetLocalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&offsetfrom)); 2485f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetLocalVecOffset(networkdm,vto,ALL_COMPONENTS,&offsetto)); 2495f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetGlobalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&goffsetfrom)); 2505f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetGlobalVecOffset(networkdm,vto,ALL_COMPONENTS,&goffsetto)); 251c4762a1bSJed Brown 252c4762a1bSJed Brown if (goffsetto < 0) goffsetto = -goffsetto - 1; 253c4762a1bSJed Brown 254c4762a1bSJed Brown thetaf = xarr[offsetfrom]; 255c4762a1bSJed Brown Vmf = xarr[offsetfrom+1]; 256c4762a1bSJed Brown thetat = xarr[offsetto]; 257c4762a1bSJed Brown Vmt = xarr[offsetto+1]; 258c4762a1bSJed Brown thetaft = thetaf - thetat; 259c4762a1bSJed Brown thetatf = thetat - thetaf; 260c4762a1bSJed Brown 2615f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetComponent(networkdm,vfrom,0,&keyf,(void**)&busf,NULL)); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetComponent(networkdm,vto,0,&keyt,(void**)&bust,NULL)); 263c4762a1bSJed Brown 264c4762a1bSJed Brown if (vfrom == vtx[v]) { 265c4762a1bSJed Brown if (busf->ide != REF_BUS) { 266c4762a1bSJed Brown /* farr[offsetfrom] += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); */ 267c4762a1bSJed Brown row[0] = goffsetfrom; 268c4762a1bSJed Brown col[0] = goffsetfrom; col[1] = goffsetfrom+1; col[2] = goffsetto; col[3] = goffsetto+1; 269c4762a1bSJed Brown values[0] = Vmf*Vmt*(Gft*-PetscSinScalar(thetaft) + Bft*PetscCosScalar(thetaft)); /* df_dthetaf */ 270c4762a1bSJed Brown values[1] = 2.0*Gff*Vmf + Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); /* df_dVmf */ 271c4762a1bSJed Brown values[2] = Vmf*Vmt*(Gft*PetscSinScalar(thetaft) + Bft*-PetscCosScalar(thetaft)); /* df_dthetat */ 272c4762a1bSJed Brown values[3] = Vmf*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); /* df_dVmt */ 273c4762a1bSJed Brown 2745f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(J,1,row,4,col,values,ADD_VALUES)); 275c4762a1bSJed Brown } 276c4762a1bSJed Brown if (busf->ide != PV_BUS && busf->ide != REF_BUS) { 277c4762a1bSJed Brown row[0] = goffsetfrom+1; 278c4762a1bSJed Brown col[0] = goffsetfrom; col[1] = goffsetfrom+1; col[2] = goffsetto; col[3] = goffsetto+1; 279c4762a1bSJed Brown /* farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); */ 280c4762a1bSJed Brown values[0] = Vmf*Vmt*(Bft*PetscSinScalar(thetaft) + Gft*PetscCosScalar(thetaft)); 281c4762a1bSJed Brown values[1] = -2.0*Bff*Vmf + Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); 282c4762a1bSJed Brown values[2] = Vmf*Vmt*(-Bft*PetscSinScalar(thetaft) + Gft*-PetscCosScalar(thetaft)); 283c4762a1bSJed Brown values[3] = Vmf*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); 284c4762a1bSJed Brown 2855f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(J,1,row,4,col,values,ADD_VALUES)); 286c4762a1bSJed Brown } 287c4762a1bSJed Brown } else { 288c4762a1bSJed Brown if (bust->ide != REF_BUS) { 289c4762a1bSJed Brown row[0] = goffsetto; 290c4762a1bSJed Brown col[0] = goffsetto; col[1] = goffsetto+1; col[2] = goffsetfrom; col[3] = goffsetfrom+1; 291c4762a1bSJed Brown /* farr[offsetto] += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); */ 292c4762a1bSJed Brown values[0] = Vmt*Vmf*(Gtf*-PetscSinScalar(thetatf) + Btf*PetscCosScalar(thetaft)); /* df_dthetat */ 293c4762a1bSJed Brown values[1] = 2.0*Gtt*Vmt + Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); /* df_dVmt */ 294c4762a1bSJed Brown values[2] = Vmt*Vmf*(Gtf*PetscSinScalar(thetatf) + Btf*-PetscCosScalar(thetatf)); /* df_dthetaf */ 295c4762a1bSJed Brown values[3] = Vmt*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); /* df_dVmf */ 296c4762a1bSJed Brown 2975f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(J,1,row,4,col,values,ADD_VALUES)); 298c4762a1bSJed Brown } 299c4762a1bSJed Brown if (bust->ide != PV_BUS && bust->ide != REF_BUS) { 300c4762a1bSJed Brown row[0] = goffsetto+1; 301c4762a1bSJed Brown col[0] = goffsetto; col[1] = goffsetto+1; col[2] = goffsetfrom; col[3] = goffsetfrom+1; 302c4762a1bSJed Brown /* farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); */ 303c4762a1bSJed Brown values[0] = Vmt*Vmf*(Btf*PetscSinScalar(thetatf) + Gtf*PetscCosScalar(thetatf)); 304c4762a1bSJed Brown values[1] = -2.0*Btt*Vmt + Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); 305c4762a1bSJed Brown values[2] = Vmt*Vmf*(-Btf*PetscSinScalar(thetatf) + Gtf*-PetscCosScalar(thetatf)); 306c4762a1bSJed Brown values[3] = Vmt*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); 307c4762a1bSJed Brown 3085f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(J,1,row,4,col,values,ADD_VALUES)); 309c4762a1bSJed Brown } 310c4762a1bSJed Brown } 311c4762a1bSJed Brown } 312c4762a1bSJed Brown if (!ghostvtex && bus->ide == PV_BUS) { 313c4762a1bSJed Brown row[0] = goffset+1; col[0] = goffset+1; 314c4762a1bSJed Brown values[0] = 1.0; 3155f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(J,1,row,1,col,values,ADD_VALUES)); 316c4762a1bSJed Brown } 317c4762a1bSJed Brown } 318c4762a1bSJed Brown } 319c4762a1bSJed Brown } 3205f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(localX,&xarr)); 321c4762a1bSJed Brown PetscFunctionReturn(0); 322c4762a1bSJed Brown } 323c4762a1bSJed Brown 324c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec X, Mat J,Mat Jpre,void *appctx) 325c4762a1bSJed Brown { 326c4762a1bSJed Brown DM networkdm; 327c4762a1bSJed Brown Vec localX; 328c4762a1bSJed Brown PetscInt ne,nv; 329c4762a1bSJed Brown const PetscInt *vtx,*edges; 330c4762a1bSJed Brown 331c4762a1bSJed Brown PetscFunctionBegin; 3325f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroEntries(J)); 333c4762a1bSJed Brown 3345f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetDM(snes,&networkdm)); 3355f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(networkdm,&localX)); 336c4762a1bSJed Brown 3375f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX)); 3385f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX)); 339c4762a1bSJed Brown 340c4762a1bSJed Brown /* Form Jacobian for first subnetwork */ 3415f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges)); 3425f80ce2aSJacob Faibussowitsch CHKERRQ(FormJacobian_Subnet(networkdm,localX,J,Jpre,nv,ne,vtx,edges,appctx)); 343c4762a1bSJed Brown 344c4762a1bSJed Brown /* Form Jacobian for second subnetwork */ 3455f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetSubnetwork(networkdm,1,&nv,&ne,&vtx,&edges)); 3465f80ce2aSJacob Faibussowitsch CHKERRQ(FormJacobian_Subnet(networkdm,localX,J,Jpre,nv,ne,vtx,edges,appctx)); 347c4762a1bSJed Brown 3485f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(networkdm,&localX)); 349c4762a1bSJed Brown 3505f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 3515f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 352c4762a1bSJed Brown PetscFunctionReturn(0); 353c4762a1bSJed Brown } 354c4762a1bSJed Brown 355c4762a1bSJed Brown PetscErrorCode SetInitialValues_Subnet(DM networkdm,Vec localX,PetscInt nv,PetscInt ne, const PetscInt *vtx, const PetscInt *edges,void* appctx) 356c4762a1bSJed Brown { 357c4762a1bSJed Brown VERTEX_Power bus; 358c4762a1bSJed Brown PetscInt i; 359c4762a1bSJed Brown GEN gen; 360c4762a1bSJed Brown PetscBool ghostvtex; 361c4762a1bSJed Brown PetscScalar *xarr; 362c4762a1bSJed Brown PetscInt key,numComps,j,offset; 363c4762a1bSJed Brown void* component; 364c4762a1bSJed Brown 365c4762a1bSJed Brown PetscFunctionBegin; 3665f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(localX,&xarr)); 367c4762a1bSJed Brown for (i = 0; i < nv; i++) { 3685f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkIsGhostVertex(networkdm,vtx[i],&ghostvtex)); 369c4762a1bSJed Brown if (ghostvtex) continue; 370c4762a1bSJed Brown 3715f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetLocalVecOffset(networkdm,vtx[i],ALL_COMPONENTS,&offset)); 3725f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetNumComponents(networkdm,vtx[i],&numComps)); 373c4762a1bSJed Brown for (j=0; j < numComps; j++) { 3745f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetComponent(networkdm,vtx[i],j,&key,&component,NULL)); 375c4762a1bSJed Brown if (key == 1) { 376c4762a1bSJed Brown bus = (VERTEX_Power)(component); 377c4762a1bSJed Brown xarr[offset] = bus->va*PETSC_PI/180.0; 378c4762a1bSJed Brown xarr[offset+1] = bus->vm; 379c4762a1bSJed Brown } else if (key == 2) { 380c4762a1bSJed Brown gen = (GEN)(component); 381c4762a1bSJed Brown if (!gen->status) continue; 382c4762a1bSJed Brown xarr[offset+1] = gen->vs; 383c4762a1bSJed Brown break; 384c4762a1bSJed Brown } 385c4762a1bSJed Brown } 386c4762a1bSJed Brown } 3875f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(localX,&xarr)); 388c4762a1bSJed Brown PetscFunctionReturn(0); 389c4762a1bSJed Brown } 390c4762a1bSJed Brown 391c4762a1bSJed Brown PetscErrorCode SetInitialValues(DM networkdm, Vec X,void* appctx) 392c4762a1bSJed Brown { 393c4762a1bSJed Brown PetscInt nv,ne; 394c4762a1bSJed Brown const PetscInt *vtx,*edges; 395c4762a1bSJed Brown Vec localX; 396c4762a1bSJed Brown 397c4762a1bSJed Brown PetscFunctionBegin; 3985f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(networkdm,&localX)); 399c4762a1bSJed Brown 4005f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(X,0.0)); 4015f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX)); 4025f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX)); 403c4762a1bSJed Brown 404c4762a1bSJed Brown /* Set initial guess for first subnetwork */ 4055f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges)); 4065f80ce2aSJacob Faibussowitsch CHKERRQ(SetInitialValues_Subnet(networkdm,localX,nv,ne,vtx,edges,appctx)); 407c4762a1bSJed Brown 408c4762a1bSJed Brown /* Set initial guess for second subnetwork */ 4095f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetSubnetwork(networkdm,1,&nv,&ne,&vtx,&edges)); 4105f80ce2aSJacob Faibussowitsch CHKERRQ(SetInitialValues_Subnet(networkdm,localX,nv,ne,vtx,edges,appctx)); 411c4762a1bSJed Brown 4125f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalBegin(networkdm,localX,ADD_VALUES,X)); 4135f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalEnd(networkdm,localX,ADD_VALUES,X)); 4145f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(networkdm,&localX)); 415c4762a1bSJed Brown PetscFunctionReturn(0); 416c4762a1bSJed Brown } 417c4762a1bSJed Brown 418c4762a1bSJed Brown int main(int argc,char ** argv) 419c4762a1bSJed Brown { 420c4762a1bSJed Brown char pfdata_file[PETSC_MAX_PATH_LEN]="case9.m"; 421c4762a1bSJed Brown PFDATA *pfdata1,*pfdata2; 422f11a936eSBarry Smith PetscInt numEdges1=0,numEdges2=0; 4232bf73ac6SHong Zhang PetscInt *edgelist1 = NULL,*edgelist2 = NULL,componentkey[4]; 424c4762a1bSJed Brown DM networkdm; 425c4762a1bSJed Brown UserCtx_Power User; 426956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 427c4762a1bSJed Brown PetscLogStage stage1,stage2; 428956f8c0dSBarry Smith #endif 429c4762a1bSJed Brown PetscMPIInt rank; 4302bf73ac6SHong Zhang PetscInt nsubnet = 2,nv,ne,i,j,genj,loadj; 4312bf73ac6SHong Zhang const PetscInt *vtx,*edges; 432c4762a1bSJed Brown Vec X,F; 433c4762a1bSJed Brown Mat J; 434c4762a1bSJed Brown SNES snes; 435c4762a1bSJed Brown 436*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,"poweroptions",help)); 4375f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 438c4762a1bSJed Brown { 439c4762a1bSJed 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 */ 440c4762a1bSJed Brown /* this is an experiment to see how the analyzer reacts */ 441c4762a1bSJed Brown const PetscMPIInt crank = rank; 442c4762a1bSJed Brown 443c4762a1bSJed Brown /* Create an empty network object */ 4445f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkCreate(PETSC_COMM_WORLD,&networkdm)); 445c4762a1bSJed Brown 446c4762a1bSJed Brown /* Register the components in the network */ 4475f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkRegisterComponent(networkdm,"branchstruct",sizeof(struct _p_EDGE_Power),&componentkey[0])); 4485f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkRegisterComponent(networkdm,"busstruct",sizeof(struct _p_VERTEX_Power),&componentkey[1])); 4495f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkRegisterComponent(networkdm,"genstruct",sizeof(struct _p_GEN),&componentkey[2])); 4505f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkRegisterComponent(networkdm,"loadstruct",sizeof(struct _p_LOAD),&componentkey[3])); 451c4762a1bSJed Brown 4525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister("Read Data",&stage1)); 453c4762a1bSJed Brown PetscLogStagePush(stage1); 454c4762a1bSJed Brown /* READ THE DATA */ 455c4762a1bSJed Brown if (!crank) { 456c4762a1bSJed Brown /* Only rank 0 reads the data */ 4575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-pfdata",pfdata_file,sizeof(pfdata_file),NULL)); 458c4762a1bSJed Brown /* HERE WE CREATE COPIES OF THE SAME NETWORK THAT WILL BE TREATED AS SUBNETWORKS */ 459c4762a1bSJed Brown 460c4762a1bSJed Brown /* READ DATA FOR THE FIRST SUBNETWORK */ 4615f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&pfdata1)); 4625f80ce2aSJacob Faibussowitsch CHKERRQ(PFReadMatPowerData(pfdata1,pfdata_file)); 463c4762a1bSJed Brown User.Sbase = pfdata1->sbase; 464c4762a1bSJed Brown 465c4762a1bSJed Brown numEdges1 = pfdata1->nbranch; 4665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(2*numEdges1,&edgelist1)); 4675f80ce2aSJacob Faibussowitsch CHKERRQ(GetListofEdges_Power(pfdata1,edgelist1)); 468c4762a1bSJed Brown 469c4762a1bSJed Brown /* READ DATA FOR THE SECOND SUBNETWORK */ 4705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&pfdata2)); 4715f80ce2aSJacob Faibussowitsch CHKERRQ(PFReadMatPowerData(pfdata2,pfdata_file)); 472c4762a1bSJed Brown User.Sbase = pfdata2->sbase; 473c4762a1bSJed Brown 474c4762a1bSJed Brown numEdges2 = pfdata2->nbranch; 4755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(2*numEdges2,&edgelist2)); 4765f80ce2aSJacob Faibussowitsch CHKERRQ(GetListofEdges_Power(pfdata2,edgelist2)); 477c4762a1bSJed Brown } 478c4762a1bSJed Brown 479c4762a1bSJed Brown PetscLogStagePop(); 4805f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD)); 4815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister("Create network",&stage2)); 482c4762a1bSJed Brown PetscLogStagePush(stage2); 483c4762a1bSJed Brown 4842bf73ac6SHong Zhang /* Set number of nodes/edges and edge connectivity */ 4855f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkSetNumSubNetworks(networkdm,PETSC_DECIDE,nsubnet)); 4865f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddSubnetwork(networkdm,"",numEdges1,edgelist1,NULL)); 4875f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddSubnetwork(networkdm,"",numEdges2,edgelist2,NULL)); 488c4762a1bSJed Brown 489c4762a1bSJed Brown /* Set up the network layout */ 4905f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkLayoutSetUp(networkdm)); 491c4762a1bSJed Brown 492c4762a1bSJed Brown /* Add network components only process 0 has any data to add*/ 493c4762a1bSJed Brown if (!crank) { 494c4762a1bSJed Brown genj=0; loadj=0; 495c4762a1bSJed Brown 496c4762a1bSJed Brown /* ADD VARIABLES AND COMPONENTS FOR THE FIRST SUBNETWORK */ 4975f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges)); 498c4762a1bSJed Brown 499c4762a1bSJed Brown for (i = 0; i < ne; i++) { 5005f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddComponent(networkdm,edges[i],componentkey[0],&pfdata1->branch[i],0)); 501c4762a1bSJed Brown } 502c4762a1bSJed Brown 503c4762a1bSJed Brown for (i = 0; i < nv; i++) { 5045f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddComponent(networkdm,vtx[i],componentkey[1],&pfdata1->bus[i],2)); 505c4762a1bSJed Brown if (pfdata1->bus[i].ngen) { 506c4762a1bSJed Brown for (j = 0; j < pfdata1->bus[i].ngen; j++) { 5075f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddComponent(networkdm,vtx[i],componentkey[2],&pfdata1->gen[genj++],0)); 508c4762a1bSJed Brown } 509c4762a1bSJed Brown } 510c4762a1bSJed Brown if (pfdata1->bus[i].nload) { 511c4762a1bSJed Brown for (j=0; j < pfdata1->bus[i].nload; j++) { 5125f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddComponent(networkdm,vtx[i],componentkey[3],&pfdata1->load[loadj++],0)); 513c4762a1bSJed Brown } 514c4762a1bSJed Brown } 515c4762a1bSJed Brown } 516c4762a1bSJed Brown 517c4762a1bSJed Brown genj=0; loadj=0; 518c4762a1bSJed Brown 519c4762a1bSJed Brown /* ADD VARIABLES AND COMPONENTS FOR THE SECOND SUBNETWORK */ 5205f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetSubnetwork(networkdm,1,&nv,&ne,&vtx,&edges)); 521c4762a1bSJed Brown 522c4762a1bSJed Brown for (i = 0; i < ne; i++) { 5235f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddComponent(networkdm,edges[i],componentkey[0],&pfdata2->branch[i],0)); 524c4762a1bSJed Brown } 525c4762a1bSJed Brown 526c4762a1bSJed Brown for (i = 0; i < nv; i++) { 5275f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddComponent(networkdm,vtx[i],componentkey[1],&pfdata2->bus[i],2)); 528c4762a1bSJed Brown if (pfdata2->bus[i].ngen) { 529c4762a1bSJed Brown for (j = 0; j < pfdata2->bus[i].ngen; j++) { 5305f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddComponent(networkdm,vtx[i],componentkey[2],&pfdata2->gen[genj++],0)); 531c4762a1bSJed Brown } 532c4762a1bSJed Brown } 533c4762a1bSJed Brown if (pfdata2->bus[i].nload) { 534c4762a1bSJed Brown for (j=0; j < pfdata2->bus[i].nload; j++) { 5355f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddComponent(networkdm,vtx[i],componentkey[3],&pfdata2->load[loadj++],0)); 536c4762a1bSJed Brown } 537c4762a1bSJed Brown } 538c4762a1bSJed Brown } 539c4762a1bSJed Brown } 540c4762a1bSJed Brown 541c4762a1bSJed Brown /* Set up DM for use */ 5425f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(networkdm)); 543c4762a1bSJed Brown 544c4762a1bSJed Brown if (!crank) { 5455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(edgelist1)); 5465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(edgelist2)); 547c4762a1bSJed Brown } 548c4762a1bSJed Brown 549c4762a1bSJed Brown if (!crank) { 5505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pfdata1->bus)); 5515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pfdata1->gen)); 5525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pfdata1->branch)); 5535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pfdata1->load)); 5545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pfdata1)); 555c4762a1bSJed Brown 5565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pfdata2->bus)); 5575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pfdata2->gen)); 5585f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pfdata2->branch)); 5595f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pfdata2->load)); 5605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pfdata2)); 561c4762a1bSJed Brown } 562c4762a1bSJed Brown 563c4762a1bSJed Brown /* Distribute networkdm to multiple processes */ 5645f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkDistribute(&networkdm,0)); 565c4762a1bSJed Brown 566c4762a1bSJed Brown PetscLogStagePop(); 567c4762a1bSJed Brown 568c4762a1bSJed Brown /* Broadcast Sbase to all processors */ 5695f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Bcast(&User.Sbase,1,MPIU_SCALAR,0,PETSC_COMM_WORLD)); 570c4762a1bSJed Brown 5715f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(networkdm,&X)); 5725f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(X,&F)); 573c4762a1bSJed Brown 5745f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(networkdm,&J)); 5755f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 576c4762a1bSJed Brown 5775f80ce2aSJacob Faibussowitsch CHKERRQ(SetInitialValues(networkdm,X,&User)); 578c4762a1bSJed Brown 579c4762a1bSJed Brown /* HOOK UP SOLVER */ 5805f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes)); 5815f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetDM(snes,networkdm)); 5825f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(snes,F,FormFunction,&User)); 5835f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian,&User)); 5845f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 585c4762a1bSJed Brown 5865f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes,NULL,X)); 5875f80ce2aSJacob Faibussowitsch /* CHKERRQ(VecView(X,PETSC_VIEWER_STDOUT_WORLD)); */ 588c4762a1bSJed Brown 5895f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X)); 5905f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&F)); 5915f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); 592c4762a1bSJed Brown 5935f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&snes)); 5945f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&networkdm)); 595c4762a1bSJed Brown } 596*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 597*b122ec5aSJacob Faibussowitsch return 0; 598c4762a1bSJed Brown } 599c4762a1bSJed Brown 600c4762a1bSJed Brown /*TEST 601c4762a1bSJed Brown 602c4762a1bSJed Brown build: 603c4762a1bSJed Brown depends: PFReadData.c pffunctions.c 604dfd57a17SPierre Jolivet requires: !complex double defined(PETSC_HAVE_ATTRIBUTEALIGNED) 605c4762a1bSJed Brown 606c4762a1bSJed Brown test: 607c4762a1bSJed Brown args: -snes_rtol 1.e-3 608c4762a1bSJed Brown localrunfiles: poweroptions case9.m 6092bf73ac6SHong Zhang output_file: output/power_1.out 610c4762a1bSJed Brown 611c4762a1bSJed Brown test: 612c4762a1bSJed Brown suffix: 2 613c4762a1bSJed Brown args: -snes_rtol 1.e-3 -petscpartitioner_type simple 614c4762a1bSJed Brown nsize: 4 615c4762a1bSJed Brown localrunfiles: poweroptions case9.m 6162bf73ac6SHong Zhang output_file: output/power_1.out 617c4762a1bSJed Brown 618c4762a1bSJed Brown TEST*/ 619