1c4762a1bSJed Brown static char help[] = "This example demonstrates the use of DMNetwork interface 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 See 'Evaluation of overlapping restricted additive schwarz preconditioning for parallel solution \n\ 4c4762a1bSJed Brown of very large power flow problems' https://dl.acm.org/citation.cfm?id=2536784).\n\ 5c4762a1bSJed Brown The data file format used is from the MatPower package (http://www.pserc.cornell.edu//matpower/).\n\ 6c4762a1bSJed Brown Run this program: mpiexec -n <n> ./pf\n\ 7c4762a1bSJed Brown mpiexec -n <n> ./pfc \n"; 8c4762a1bSJed Brown 9c4762a1bSJed Brown /* T 10c4762a1bSJed Brown Concepts: DMNetwork 11c4762a1bSJed Brown Concepts: PETSc SNES solver 12c4762a1bSJed Brown */ 13c4762a1bSJed Brown 14c4762a1bSJed Brown #include "power.h" 15c4762a1bSJed Brown #include <petscdmnetwork.h> 16c4762a1bSJed Brown 17c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X, Vec F,void *appctx) 18c4762a1bSJed Brown { 19c4762a1bSJed Brown DM networkdm; 20c4762a1bSJed Brown UserCtx_Power *User=(UserCtx_Power*)appctx; 21c4762a1bSJed Brown Vec localX,localF; 22c4762a1bSJed Brown PetscInt nv,ne; 23c4762a1bSJed Brown const PetscInt *vtx,*edges; 24c4762a1bSJed Brown 25c4762a1bSJed Brown PetscFunctionBegin; 265f80ce2aSJacob Faibussowitsch CHKERRQ(SNESGetDM(snes,&networkdm)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(networkdm,&localX)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(networkdm,&localF)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(F,0.0)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(localF,0.0)); 31c4762a1bSJed Brown 325f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX)); 335f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX)); 34c4762a1bSJed Brown 355f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(FormFunction_Power(networkdm,localX,localF,nv,ne,vtx,edges,User)); 37c4762a1bSJed Brown 385f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(networkdm,&localX)); 39c4762a1bSJed Brown 405f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalBegin(networkdm,localF,ADD_VALUES,F)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalEnd(networkdm,localF,ADD_VALUES,F)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(networkdm,&localF)); 43c4762a1bSJed Brown PetscFunctionReturn(0); 44c4762a1bSJed Brown } 45c4762a1bSJed Brown 46c4762a1bSJed Brown PetscErrorCode SetInitialValues(DM networkdm,Vec X,void* appctx) 47c4762a1bSJed Brown { 48c4762a1bSJed Brown PetscInt vStart,vEnd,nv,ne; 49c4762a1bSJed Brown const PetscInt *vtx,*edges; 50c4762a1bSJed Brown Vec localX; 51c4762a1bSJed Brown UserCtx_Power *user_power=(UserCtx_Power*)appctx; 52c4762a1bSJed Brown 53c4762a1bSJed Brown PetscFunctionBegin; 545f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetVertexRange(networkdm,&vStart, &vEnd)); 55c4762a1bSJed Brown 565f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(networkdm,&localX)); 57c4762a1bSJed Brown 585f80ce2aSJacob Faibussowitsch CHKERRQ(VecSet(X,0.0)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX)); 61c4762a1bSJed Brown 625f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges)); 635f80ce2aSJacob Faibussowitsch CHKERRQ(SetInitialGuess_Power(networkdm,localX,nv,ne,vtx,edges,user_power)); 64c4762a1bSJed Brown 655f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalBegin(networkdm,localX,ADD_VALUES,X)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(DMLocalToGlobalEnd(networkdm,localX,ADD_VALUES,X)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(networkdm,&localX)); 68c4762a1bSJed Brown PetscFunctionReturn(0); 69c4762a1bSJed Brown } 70c4762a1bSJed Brown 71c4762a1bSJed Brown int main(int argc,char ** argv) 72c4762a1bSJed Brown { 73c4762a1bSJed Brown char pfdata_file[PETSC_MAX_PATH_LEN]="case9.m"; 74c4762a1bSJed Brown PFDATA *pfdata; 75f11a936eSBarry Smith PetscInt numEdges=0; 76c4762a1bSJed Brown PetscInt *edges = NULL; 77c4762a1bSJed Brown PetscInt i; 78c4762a1bSJed Brown DM networkdm; 79c4762a1bSJed Brown UserCtx_Power User; 80956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 81c4762a1bSJed Brown PetscLogStage stage1,stage2; 82956f8c0dSBarry Smith #endif 83c4762a1bSJed Brown PetscMPIInt rank; 84c4762a1bSJed Brown PetscInt eStart, eEnd, vStart, vEnd,j; 85c4762a1bSJed Brown PetscInt genj,loadj; 86c4762a1bSJed Brown Vec X,F; 87c4762a1bSJed Brown Mat J; 88c4762a1bSJed Brown SNES snes; 89c4762a1bSJed Brown 90*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,"poweroptions",help)); 915f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 92c4762a1bSJed Brown { 93c4762a1bSJed 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 */ 94c4762a1bSJed Brown /* this is an experiment to see how the analyzer reacts */ 95c4762a1bSJed Brown const PetscMPIInt crank = rank; 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* Create an empty network object */ 985f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkCreate(PETSC_COMM_WORLD,&networkdm)); 99c4762a1bSJed Brown /* Register the components in the network */ 1005f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkRegisterComponent(networkdm,"branchstruct",sizeof(struct _p_EDGE_Power),&User.compkey_branch)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkRegisterComponent(networkdm,"busstruct",sizeof(struct _p_VERTEX_Power),&User.compkey_bus)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkRegisterComponent(networkdm,"genstruct",sizeof(struct _p_GEN),&User.compkey_gen)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkRegisterComponent(networkdm,"loadstruct",sizeof(struct _p_LOAD),&User.compkey_load)); 104c4762a1bSJed Brown 1055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister("Read Data",&stage1)); 106c4762a1bSJed Brown PetscLogStagePush(stage1); 107c4762a1bSJed Brown /* READ THE DATA */ 108c4762a1bSJed Brown if (!crank) { 109c4762a1bSJed Brown /* READ DATA */ 110c4762a1bSJed Brown /* Only rank 0 reads the data */ 1115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-pfdata",pfdata_file,sizeof(pfdata_file),NULL)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNew(&pfdata)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(PFReadMatPowerData(pfdata,pfdata_file)); 114c4762a1bSJed Brown User.Sbase = pfdata->sbase; 115c4762a1bSJed Brown 116c4762a1bSJed Brown numEdges = pfdata->nbranch; 1175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(2*numEdges,&edges)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(GetListofEdges_Power(pfdata,edges)); 119c4762a1bSJed Brown } 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* If external option activated. Introduce error in jacobian */ 1225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL, "-jac_error", &User.jac_error)); 123c4762a1bSJed Brown 124c4762a1bSJed Brown PetscLogStagePop(); 1255f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Barrier(PETSC_COMM_WORLD)); 1265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogStageRegister("Create network",&stage2)); 127c4762a1bSJed Brown PetscLogStagePush(stage2); 128c4762a1bSJed Brown /* Set number of nodes/edges */ 1295f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkSetNumSubNetworks(networkdm,PETSC_DECIDE,1)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddSubnetwork(networkdm,"",numEdges,edges,NULL)); 1312bf73ac6SHong Zhang 132c4762a1bSJed Brown /* Set up the network layout */ 1335f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkLayoutSetUp(networkdm)); 134c4762a1bSJed Brown 135c4762a1bSJed Brown if (!crank) { 1365f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(edges)); 137c4762a1bSJed Brown } 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* Add network components only process 0 has any data to add */ 140c4762a1bSJed Brown if (!crank) { 141c4762a1bSJed Brown genj=0; loadj=0; 1425f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd)); 143c4762a1bSJed Brown for (i = eStart; i < eEnd; i++) { 1445f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddComponent(networkdm,i,User.compkey_branch,&pfdata->branch[i-eStart],0)); 145c4762a1bSJed Brown } 1465f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetVertexRange(networkdm,&vStart,&vEnd)); 147c4762a1bSJed Brown for (i = vStart; i < vEnd; i++) { 1485f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddComponent(networkdm,i,User.compkey_bus,&pfdata->bus[i-vStart],2)); 149c4762a1bSJed Brown if (pfdata->bus[i-vStart].ngen) { 150c4762a1bSJed Brown for (j = 0; j < pfdata->bus[i-vStart].ngen; j++) { 1515f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddComponent(networkdm,i,User.compkey_gen,&pfdata->gen[genj++],0)); 152c4762a1bSJed Brown } 153c4762a1bSJed Brown } 154c4762a1bSJed Brown if (pfdata->bus[i-vStart].nload) { 155c4762a1bSJed Brown for (j=0; j < pfdata->bus[i-vStart].nload; j++) { 1565f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkAddComponent(networkdm,i,User.compkey_load,&pfdata->load[loadj++],0)); 157c4762a1bSJed Brown } 158c4762a1bSJed Brown } 159c4762a1bSJed Brown } 160c4762a1bSJed Brown } 161c4762a1bSJed Brown 162c4762a1bSJed Brown /* Set up DM for use */ 1635f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(networkdm)); 164c4762a1bSJed Brown 165c4762a1bSJed Brown if (!crank) { 1665f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pfdata->bus)); 1675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pfdata->gen)); 1685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pfdata->branch)); 1695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pfdata->load)); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pfdata)); 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173c4762a1bSJed Brown /* Distribute networkdm to multiple processes */ 1745f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkDistribute(&networkdm,0)); 175c4762a1bSJed Brown 176c4762a1bSJed Brown PetscLogStagePop(); 1775f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd)); 1785f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetVertexRange(networkdm,&vStart,&vEnd)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown #if 0 181c4762a1bSJed Brown EDGE_Power edge; 182c4762a1bSJed Brown PetscInt key,kk,numComponents; 183c4762a1bSJed Brown VERTEX_Power bus; 184c4762a1bSJed Brown GEN gen; 185c4762a1bSJed Brown LOAD load; 186c4762a1bSJed Brown 187c4762a1bSJed Brown for (i = eStart; i < eEnd; i++) { 1885f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetComponent(networkdm,i,0,&key,(void**)&edge)); 1895f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetNumComponents(networkdm,i,&numComponents)); 1905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Rank %d ncomps = %d Line %d ---- %d\n",crank,numComponents,edge->internal_i,edge->internal_j)); 191c4762a1bSJed Brown } 192c4762a1bSJed Brown 193c4762a1bSJed Brown for (i = vStart; i < vEnd; i++) { 1945f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetNumComponents(networkdm,i,&numComponents)); 195c4762a1bSJed Brown for (kk=0; kk < numComponents; kk++) { 1965f80ce2aSJacob Faibussowitsch CHKERRQ(DMNetworkGetComponent(networkdm,i,kk,&key,&component)); 197c4762a1bSJed Brown if (key == 1) { 198c4762a1bSJed Brown bus = (VERTEX_Power)(component); 1995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Rank %d ncomps = %d Bus %d\n",crank,numComponents,bus->internal_i)); 200c4762a1bSJed Brown } else if (key == 2) { 201c4762a1bSJed Brown gen = (GEN)(component); 2025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Rank %d Gen pg = %f qg = %f\n",crank,gen->pg,gen->qg)); 203c4762a1bSJed Brown } else if (key == 3) { 204c4762a1bSJed Brown load = (LOAD)(component); 2055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"Rank %d Load pl = %f ql = %f\n",crank,load->pl,load->ql)); 206c4762a1bSJed Brown } 207c4762a1bSJed Brown } 208c4762a1bSJed Brown } 209c4762a1bSJed Brown #endif 210c4762a1bSJed Brown /* Broadcast Sbase to all processors */ 2115f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Bcast(&User.Sbase,1,MPIU_SCALAR,0,PETSC_COMM_WORLD)); 212c4762a1bSJed Brown 2135f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(networkdm,&X)); 2145f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(X,&F)); 215c4762a1bSJed Brown 2165f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(networkdm,&J)); 2175f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 218c4762a1bSJed Brown 2195f80ce2aSJacob Faibussowitsch CHKERRQ(SetInitialValues(networkdm,X,&User)); 220c4762a1bSJed Brown 221c4762a1bSJed Brown /* HOOK UP SOLVER */ 2225f80ce2aSJacob Faibussowitsch CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes)); 2235f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetDM(snes,networkdm)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFunction(snes,F,FormFunction,&User)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetJacobian(snes,J,J,FormJacobian_Power,&User)); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSetFromOptions(snes)); 227c4762a1bSJed Brown 2285f80ce2aSJacob Faibussowitsch CHKERRQ(SNESSolve(snes,NULL,X)); 2295f80ce2aSJacob Faibussowitsch /* CHKERRQ(VecView(X,PETSC_VIEWER_STDOUT_WORLD)); */ 230c4762a1bSJed Brown 2315f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&X)); 2325f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&F)); 2335f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); 234c4762a1bSJed Brown 2355f80ce2aSJacob Faibussowitsch CHKERRQ(SNESDestroy(&snes)); 2365f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&networkdm)); 237c4762a1bSJed Brown } 238*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 239*b122ec5aSJacob Faibussowitsch return 0; 240c4762a1bSJed Brown } 241c4762a1bSJed Brown 242c4762a1bSJed Brown /*TEST 243c4762a1bSJed Brown 244c4762a1bSJed Brown build: 245c4762a1bSJed Brown depends: PFReadData.c pffunctions.c 246dfd57a17SPierre Jolivet requires: !complex double defined(PETSC_HAVE_ATTRIBUTEALIGNED) 247c4762a1bSJed Brown 248c4762a1bSJed Brown test: 249c4762a1bSJed Brown args: -snes_rtol 1.e-3 250c4762a1bSJed Brown localrunfiles: poweroptions case9.m 251c4762a1bSJed Brown output_file: output/power_1.out 252c4762a1bSJed Brown 253c4762a1bSJed Brown test: 254c4762a1bSJed Brown suffix: 2 255c4762a1bSJed Brown args: -snes_rtol 1.e-3 -petscpartitioner_type simple 256c4762a1bSJed Brown nsize: 4 257c4762a1bSJed Brown localrunfiles: poweroptions case9.m 258c4762a1bSJed Brown output_file: output/power_1.out 259c4762a1bSJed Brown 260c4762a1bSJed Brown TEST*/ 261