1*c4762a1bSJed Brown static char help[] = "This example demonstrates the use of DMNetwork interface for solving a nonlinear electric power grid problem.\n\ 2*c4762a1bSJed Brown The available solver options are in the poweroptions file and the data files are in the datafiles directory.\n\ 3*c4762a1bSJed Brown See 'Evaluation of overlapping restricted additive schwarz preconditioning for parallel solution \n\ 4*c4762a1bSJed Brown of very large power flow problems' https://dl.acm.org/citation.cfm?id=2536784).\n\ 5*c4762a1bSJed Brown The data file format used is from the MatPower package (http://www.pserc.cornell.edu//matpower/).\n\ 6*c4762a1bSJed Brown Run this program: mpiexec -n <n> ./pf\n\ 7*c4762a1bSJed Brown mpiexec -n <n> ./pfc \n"; 8*c4762a1bSJed Brown 9*c4762a1bSJed Brown /* T 10*c4762a1bSJed Brown Concepts: DMNetwork 11*c4762a1bSJed Brown Concepts: PETSc SNES solver 12*c4762a1bSJed Brown */ 13*c4762a1bSJed Brown 14*c4762a1bSJed Brown #include "power.h" 15*c4762a1bSJed Brown #include <petscdmnetwork.h> 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X, Vec F,void *appctx) 18*c4762a1bSJed Brown { 19*c4762a1bSJed Brown PetscErrorCode ierr; 20*c4762a1bSJed Brown DM networkdm; 21*c4762a1bSJed Brown UserCtx_Power *User=(UserCtx_Power*)appctx; 22*c4762a1bSJed Brown Vec localX,localF; 23*c4762a1bSJed Brown PetscInt nv,ne; 24*c4762a1bSJed Brown const PetscInt *vtx,*edges; 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown PetscFunctionBegin; 27*c4762a1bSJed Brown ierr = SNESGetDM(snes,&networkdm);CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = DMGetLocalVector(networkdm,&localF);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = VecSet(F,0.0);CHKERRQ(ierr); 31*c4762a1bSJed Brown ierr = VecSet(localF,0.0);CHKERRQ(ierr); 32*c4762a1bSJed Brown 33*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 35*c4762a1bSJed Brown 36*c4762a1bSJed Brown ierr = DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 37*c4762a1bSJed Brown ierr = FormFunction_Power(networkdm,localX,localF,nv,ne,vtx,edges,User);CHKERRQ(ierr); 38*c4762a1bSJed Brown 39*c4762a1bSJed Brown ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr); 40*c4762a1bSJed Brown 41*c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(networkdm,localF,ADD_VALUES,F);CHKERRQ(ierr); 42*c4762a1bSJed Brown ierr = DMLocalToGlobalEnd(networkdm,localF,ADD_VALUES,F);CHKERRQ(ierr); 43*c4762a1bSJed Brown ierr = DMRestoreLocalVector(networkdm,&localF);CHKERRQ(ierr); 44*c4762a1bSJed Brown PetscFunctionReturn(0); 45*c4762a1bSJed Brown } 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown PetscErrorCode SetInitialValues(DM networkdm,Vec X,void* appctx) 48*c4762a1bSJed Brown { 49*c4762a1bSJed Brown PetscErrorCode ierr; 50*c4762a1bSJed Brown PetscInt vStart,vEnd,nv,ne; 51*c4762a1bSJed Brown const PetscInt *vtx,*edges; 52*c4762a1bSJed Brown Vec localX; 53*c4762a1bSJed Brown UserCtx_Power *user_power=(UserCtx_Power*)appctx; 54*c4762a1bSJed Brown 55*c4762a1bSJed Brown PetscFunctionBegin; 56*c4762a1bSJed Brown ierr = DMNetworkGetVertexRange(networkdm,&vStart, &vEnd);CHKERRQ(ierr); 57*c4762a1bSJed Brown 58*c4762a1bSJed Brown ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr); 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown ierr = VecSet(X,0.0);CHKERRQ(ierr); 61*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 62*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr); 63*c4762a1bSJed Brown 64*c4762a1bSJed Brown ierr = DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 65*c4762a1bSJed Brown ierr = SetInitialGuess_Power(networkdm,localX,nv,ne,vtx,edges,user_power);CHKERRQ(ierr); 66*c4762a1bSJed Brown 67*c4762a1bSJed Brown ierr = DMLocalToGlobalBegin(networkdm,localX,ADD_VALUES,X);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = DMLocalToGlobalEnd(networkdm,localX,ADD_VALUES,X);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr); 70*c4762a1bSJed Brown PetscFunctionReturn(0); 71*c4762a1bSJed Brown } 72*c4762a1bSJed Brown 73*c4762a1bSJed Brown int main(int argc,char ** argv) 74*c4762a1bSJed Brown { 75*c4762a1bSJed Brown PetscErrorCode ierr; 76*c4762a1bSJed Brown char pfdata_file[PETSC_MAX_PATH_LEN]="case9.m"; 77*c4762a1bSJed Brown PFDATA *pfdata; 78*c4762a1bSJed Brown PetscInt numEdges=0,numVertices=0; 79*c4762a1bSJed Brown PetscInt *edges = NULL; 80*c4762a1bSJed Brown PetscInt i; 81*c4762a1bSJed Brown DM networkdm; 82*c4762a1bSJed Brown UserCtx_Power User; 83*c4762a1bSJed Brown PetscLogStage stage1,stage2; 84*c4762a1bSJed Brown PetscMPIInt rank; 85*c4762a1bSJed Brown PetscInt eStart, eEnd, vStart, vEnd,j; 86*c4762a1bSJed Brown PetscInt genj,loadj; 87*c4762a1bSJed Brown Vec X,F; 88*c4762a1bSJed Brown Mat J; 89*c4762a1bSJed Brown SNES snes; 90*c4762a1bSJed Brown 91*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,"poweroptions",help);if (ierr) return ierr; 92*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 93*c4762a1bSJed Brown { 94*c4762a1bSJed 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 */ 95*c4762a1bSJed Brown /* this is an experiment to see how the analyzer reacts */ 96*c4762a1bSJed Brown const PetscMPIInt crank = rank; 97*c4762a1bSJed Brown 98*c4762a1bSJed Brown /* Create an empty network object */ 99*c4762a1bSJed Brown ierr = DMNetworkCreate(PETSC_COMM_WORLD,&networkdm);CHKERRQ(ierr); 100*c4762a1bSJed Brown /* Register the components in the network */ 101*c4762a1bSJed Brown ierr = DMNetworkRegisterComponent(networkdm,"branchstruct",sizeof(struct _p_EDGE_Power),&User.compkey_branch);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = DMNetworkRegisterComponent(networkdm,"busstruct",sizeof(struct _p_VERTEX_Power),&User.compkey_bus);CHKERRQ(ierr); 103*c4762a1bSJed Brown ierr = DMNetworkRegisterComponent(networkdm,"genstruct",sizeof(struct _p_GEN),&User.compkey_gen);CHKERRQ(ierr); 104*c4762a1bSJed Brown ierr = DMNetworkRegisterComponent(networkdm,"loadstruct",sizeof(struct _p_LOAD),&User.compkey_load);CHKERRQ(ierr); 105*c4762a1bSJed Brown 106*c4762a1bSJed Brown ierr = PetscLogStageRegister("Read Data",&stage1);CHKERRQ(ierr); 107*c4762a1bSJed Brown PetscLogStagePush(stage1); 108*c4762a1bSJed Brown /* READ THE DATA */ 109*c4762a1bSJed Brown if (!crank) { 110*c4762a1bSJed Brown /* READ DATA */ 111*c4762a1bSJed Brown /* Only rank 0 reads the data */ 112*c4762a1bSJed Brown ierr = PetscOptionsGetString(NULL,NULL,"-pfdata",pfdata_file,PETSC_MAX_PATH_LEN-1,NULL);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = PetscNew(&pfdata);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = PFReadMatPowerData(pfdata,pfdata_file);CHKERRQ(ierr); 115*c4762a1bSJed Brown User.Sbase = pfdata->sbase; 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown numEdges = pfdata->nbranch; 118*c4762a1bSJed Brown numVertices = pfdata->nbus; 119*c4762a1bSJed Brown 120*c4762a1bSJed Brown ierr = PetscMalloc1(2*numEdges,&edges);CHKERRQ(ierr); 121*c4762a1bSJed Brown ierr = GetListofEdges_Power(pfdata,edges);CHKERRQ(ierr); 122*c4762a1bSJed Brown } 123*c4762a1bSJed Brown 124*c4762a1bSJed Brown /* If external option activated. Introduce error in jacobian */ 125*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL, "-jac_error", &User.jac_error);CHKERRQ(ierr); 126*c4762a1bSJed Brown 127*c4762a1bSJed Brown PetscLogStagePop(); 128*c4762a1bSJed Brown ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRQ(ierr); 129*c4762a1bSJed Brown ierr = PetscLogStageRegister("Create network",&stage2);CHKERRQ(ierr); 130*c4762a1bSJed Brown PetscLogStagePush(stage2); 131*c4762a1bSJed Brown /* Set number of nodes/edges */ 132*c4762a1bSJed Brown ierr = DMNetworkSetSizes(networkdm,1,&numVertices,&numEdges,0,NULL);CHKERRQ(ierr); 133*c4762a1bSJed Brown /* Add edge connectivity */ 134*c4762a1bSJed Brown ierr = DMNetworkSetEdgeList(networkdm,&edges,NULL);CHKERRQ(ierr); 135*c4762a1bSJed Brown /* Set up the network layout */ 136*c4762a1bSJed Brown ierr = DMNetworkLayoutSetUp(networkdm);CHKERRQ(ierr); 137*c4762a1bSJed Brown 138*c4762a1bSJed Brown if (!crank) { 139*c4762a1bSJed Brown ierr = PetscFree(edges);CHKERRQ(ierr); 140*c4762a1bSJed Brown } 141*c4762a1bSJed Brown 142*c4762a1bSJed Brown /* Add network components only process 0 has any data to add*/ 143*c4762a1bSJed Brown if (!crank) { 144*c4762a1bSJed Brown genj=0; loadj=0; 145*c4762a1bSJed Brown ierr = DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);CHKERRQ(ierr); 146*c4762a1bSJed Brown for (i = eStart; i < eEnd; i++) { 147*c4762a1bSJed Brown ierr = DMNetworkAddComponent(networkdm,i,User.compkey_branch,&pfdata->branch[i-eStart]);CHKERRQ(ierr); 148*c4762a1bSJed Brown } 149*c4762a1bSJed Brown ierr = DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);CHKERRQ(ierr); 150*c4762a1bSJed Brown for (i = vStart; i < vEnd; i++) { 151*c4762a1bSJed Brown ierr = DMNetworkAddComponent(networkdm,i,User.compkey_bus,&pfdata->bus[i-vStart]);CHKERRQ(ierr); 152*c4762a1bSJed Brown if (pfdata->bus[i-vStart].ngen) { 153*c4762a1bSJed Brown for (j = 0; j < pfdata->bus[i-vStart].ngen; j++) { 154*c4762a1bSJed Brown ierr = DMNetworkAddComponent(networkdm,i,User.compkey_gen,&pfdata->gen[genj++]);CHKERRQ(ierr); 155*c4762a1bSJed Brown } 156*c4762a1bSJed Brown } 157*c4762a1bSJed Brown if (pfdata->bus[i-vStart].nload) { 158*c4762a1bSJed Brown for (j=0; j < pfdata->bus[i-vStart].nload; j++) { 159*c4762a1bSJed Brown ierr = DMNetworkAddComponent(networkdm,i,User.compkey_load,&pfdata->load[loadj++]);CHKERRQ(ierr); 160*c4762a1bSJed Brown } 161*c4762a1bSJed Brown } 162*c4762a1bSJed Brown /* Add number of variables */ 163*c4762a1bSJed Brown ierr = DMNetworkAddNumVariables(networkdm,i,2);CHKERRQ(ierr); 164*c4762a1bSJed Brown } 165*c4762a1bSJed Brown } 166*c4762a1bSJed Brown 167*c4762a1bSJed Brown /* Set up DM for use */ 168*c4762a1bSJed Brown ierr = DMSetUp(networkdm);CHKERRQ(ierr); 169*c4762a1bSJed Brown 170*c4762a1bSJed Brown if (!crank) { 171*c4762a1bSJed Brown ierr = PetscFree(pfdata->bus);CHKERRQ(ierr); 172*c4762a1bSJed Brown ierr = PetscFree(pfdata->gen);CHKERRQ(ierr); 173*c4762a1bSJed Brown ierr = PetscFree(pfdata->branch);CHKERRQ(ierr); 174*c4762a1bSJed Brown ierr = PetscFree(pfdata->load);CHKERRQ(ierr); 175*c4762a1bSJed Brown ierr = PetscFree(pfdata);CHKERRQ(ierr); 176*c4762a1bSJed Brown } 177*c4762a1bSJed Brown 178*c4762a1bSJed Brown /* Distribute networkdm to multiple processes */ 179*c4762a1bSJed Brown ierr = DMNetworkDistribute(&networkdm,0);CHKERRQ(ierr); 180*c4762a1bSJed Brown 181*c4762a1bSJed Brown PetscLogStagePop(); 182*c4762a1bSJed Brown ierr = DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);CHKERRQ(ierr); 183*c4762a1bSJed Brown ierr = DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);CHKERRQ(ierr); 184*c4762a1bSJed Brown 185*c4762a1bSJed Brown #if 0 186*c4762a1bSJed Brown EDGE_Power edge; 187*c4762a1bSJed Brown PetscInt key,kk,numComponents; 188*c4762a1bSJed Brown VERTEX_Power bus; 189*c4762a1bSJed Brown GEN gen; 190*c4762a1bSJed Brown LOAD load; 191*c4762a1bSJed Brown 192*c4762a1bSJed Brown for (i = eStart; i < eEnd; i++) { 193*c4762a1bSJed Brown ierr = DMNetworkGetComponent(networkdm,i,0,&key,(void**)&edge);CHKERRQ(ierr); 194*c4762a1bSJed Brown ierr = DMNetworkGetNumComponents(networkdm,i,&numComponents);CHKERRQ(ierr); 195*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Rank %d ncomps = %d Line %d ---- %d\n",crank,numComponents,edge->internal_i,edge->internal_j);CHKERRQ(ierr); 196*c4762a1bSJed Brown } 197*c4762a1bSJed Brown 198*c4762a1bSJed Brown for (i = vStart; i < vEnd; i++) { 199*c4762a1bSJed Brown ierr = DMNetworkGetNumComponents(networkdm,i,&numComponents);CHKERRQ(ierr); 200*c4762a1bSJed Brown for (kk=0; kk < numComponents; kk++) { 201*c4762a1bSJed Brown ierr = DMNetworkGetComponent(networkdm,i,kk,&key,&component);CHKERRQ(ierr); 202*c4762a1bSJed Brown if (key == 1) { 203*c4762a1bSJed Brown bus = (VERTEX_Power)(component); 204*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Rank %d ncomps = %d Bus %d\n",crank,numComponents,bus->internal_i);CHKERRQ(ierr); 205*c4762a1bSJed Brown } else if (key == 2) { 206*c4762a1bSJed Brown gen = (GEN)(component); 207*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Rank %d Gen pg = %f qg = %f\n",crank,gen->pg,gen->qg);CHKERRQ(ierr); 208*c4762a1bSJed Brown } else if (key == 3) { 209*c4762a1bSJed Brown load = (LOAD)(component); 210*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_SELF,"Rank %d Load pl = %f ql = %f\n",crank,load->pl,load->ql);CHKERRQ(ierr); 211*c4762a1bSJed Brown } 212*c4762a1bSJed Brown } 213*c4762a1bSJed Brown } 214*c4762a1bSJed Brown #endif 215*c4762a1bSJed Brown /* Broadcast Sbase to all processors */ 216*c4762a1bSJed Brown ierr = MPI_Bcast(&User.Sbase,1,MPIU_SCALAR,0,PETSC_COMM_WORLD);CHKERRQ(ierr); 217*c4762a1bSJed Brown 218*c4762a1bSJed Brown ierr = DMCreateGlobalVector(networkdm,&X);CHKERRQ(ierr); 219*c4762a1bSJed Brown ierr = VecDuplicate(X,&F);CHKERRQ(ierr); 220*c4762a1bSJed Brown 221*c4762a1bSJed Brown ierr = DMCreateMatrix(networkdm,&J);CHKERRQ(ierr); 222*c4762a1bSJed Brown ierr = MatSetOption(J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr); 223*c4762a1bSJed Brown 224*c4762a1bSJed Brown ierr = SetInitialValues(networkdm,X,&User);CHKERRQ(ierr); 225*c4762a1bSJed Brown 226*c4762a1bSJed Brown /* HOOK UP SOLVER */ 227*c4762a1bSJed Brown ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr); 228*c4762a1bSJed Brown ierr = SNESSetDM(snes,networkdm);CHKERRQ(ierr); 229*c4762a1bSJed Brown ierr = SNESSetFunction(snes,F,FormFunction,&User);CHKERRQ(ierr); 230*c4762a1bSJed Brown ierr = SNESSetJacobian(snes,J,J,FormJacobian_Power,&User);CHKERRQ(ierr); 231*c4762a1bSJed Brown ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 232*c4762a1bSJed Brown 233*c4762a1bSJed Brown ierr = SNESSolve(snes,NULL,X);CHKERRQ(ierr); 234*c4762a1bSJed Brown /* ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 235*c4762a1bSJed Brown 236*c4762a1bSJed Brown ierr = VecDestroy(&X);CHKERRQ(ierr); 237*c4762a1bSJed Brown ierr = VecDestroy(&F);CHKERRQ(ierr); 238*c4762a1bSJed Brown ierr = MatDestroy(&J);CHKERRQ(ierr); 239*c4762a1bSJed Brown 240*c4762a1bSJed Brown ierr = SNESDestroy(&snes);CHKERRQ(ierr); 241*c4762a1bSJed Brown ierr = DMDestroy(&networkdm);CHKERRQ(ierr); 242*c4762a1bSJed Brown } 243*c4762a1bSJed Brown ierr = PetscFinalize(); 244*c4762a1bSJed Brown return ierr; 245*c4762a1bSJed Brown } 246*c4762a1bSJed Brown 247*c4762a1bSJed Brown /*TEST 248*c4762a1bSJed Brown 249*c4762a1bSJed Brown build: 250*c4762a1bSJed Brown depends: PFReadData.c pffunctions.c 251*c4762a1bSJed Brown requires: !complex double define(PETSC_HAVE_ATTRIBUTEALIGNED) 252*c4762a1bSJed Brown 253*c4762a1bSJed Brown 254*c4762a1bSJed Brown test: 255*c4762a1bSJed Brown args: -snes_rtol 1.e-3 256*c4762a1bSJed Brown localrunfiles: poweroptions case9.m 257*c4762a1bSJed Brown output_file: output/power_1.out 258*c4762a1bSJed Brown requires: double !complex define(PETSC_HAVE_ATTRIBUTEALIGNED) 259*c4762a1bSJed Brown 260*c4762a1bSJed Brown test: 261*c4762a1bSJed Brown suffix: 2 262*c4762a1bSJed Brown args: -snes_rtol 1.e-3 -petscpartitioner_type simple 263*c4762a1bSJed Brown nsize: 4 264*c4762a1bSJed Brown localrunfiles: poweroptions case9.m 265*c4762a1bSJed Brown output_file: output/power_1.out 266*c4762a1bSJed Brown requires: double !complex define(PETSC_HAVE_ATTRIBUTEALIGNED) 267*c4762a1bSJed Brown 268*c4762a1bSJed Brown TEST*/ 269