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 #include "power.h" 10c4762a1bSJed Brown #include <petscdmnetwork.h> 11c4762a1bSJed Brown 12c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X, Vec F,void *appctx) 13c4762a1bSJed Brown { 14c4762a1bSJed Brown DM networkdm; 15c4762a1bSJed Brown UserCtx_Power *User=(UserCtx_Power*)appctx; 16c4762a1bSJed Brown Vec localX,localF; 17c4762a1bSJed Brown PetscInt nv,ne; 18c4762a1bSJed Brown const PetscInt *vtx,*edges; 19c4762a1bSJed Brown 20c4762a1bSJed Brown PetscFunctionBegin; 219566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&networkdm)); 229566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm,&localX)); 239566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm,&localF)); 249566063dSJacob Faibussowitsch PetscCall(VecSet(F,0.0)); 259566063dSJacob Faibussowitsch PetscCall(VecSet(localF,0.0)); 26c4762a1bSJed Brown 279566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX)); 289566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX)); 29c4762a1bSJed Brown 309566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges)); 319566063dSJacob Faibussowitsch PetscCall(FormFunction_Power(networkdm,localX,localF,nv,ne,vtx,edges,User)); 32c4762a1bSJed Brown 339566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm,&localX)); 34c4762a1bSJed Brown 359566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(networkdm,localF,ADD_VALUES,F)); 369566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(networkdm,localF,ADD_VALUES,F)); 379566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm,&localF)); 38c4762a1bSJed Brown PetscFunctionReturn(0); 39c4762a1bSJed Brown } 40c4762a1bSJed Brown 41c4762a1bSJed Brown PetscErrorCode SetInitialValues(DM networkdm,Vec X,void* appctx) 42c4762a1bSJed Brown { 43c4762a1bSJed Brown PetscInt vStart,vEnd,nv,ne; 44c4762a1bSJed Brown const PetscInt *vtx,*edges; 45c4762a1bSJed Brown Vec localX; 46c4762a1bSJed Brown UserCtx_Power *user_power=(UserCtx_Power*)appctx; 47c4762a1bSJed Brown 48c4762a1bSJed Brown PetscFunctionBegin; 499566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(networkdm,&vStart, &vEnd)); 50c4762a1bSJed Brown 519566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm,&localX)); 52c4762a1bSJed Brown 539566063dSJacob Faibussowitsch PetscCall(VecSet(X,0.0)); 549566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX)); 559566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX)); 56c4762a1bSJed Brown 579566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges)); 589566063dSJacob Faibussowitsch PetscCall(SetInitialGuess_Power(networkdm,localX,nv,ne,vtx,edges,user_power)); 59c4762a1bSJed Brown 609566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(networkdm,localX,ADD_VALUES,X)); 619566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(networkdm,localX,ADD_VALUES,X)); 629566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm,&localX)); 63c4762a1bSJed Brown PetscFunctionReturn(0); 64c4762a1bSJed Brown } 65c4762a1bSJed Brown 66c4762a1bSJed Brown int main(int argc,char ** argv) 67c4762a1bSJed Brown { 68c4762a1bSJed Brown char pfdata_file[PETSC_MAX_PATH_LEN]="case9.m"; 69c4762a1bSJed Brown PFDATA *pfdata; 70f11a936eSBarry Smith PetscInt numEdges=0; 71c4762a1bSJed Brown PetscInt *edges = NULL; 72c4762a1bSJed Brown PetscInt i; 73c4762a1bSJed Brown DM networkdm; 74c4762a1bSJed Brown UserCtx_Power User; 75956f8c0dSBarry Smith #if defined(PETSC_USE_LOG) 76c4762a1bSJed Brown PetscLogStage stage1,stage2; 77956f8c0dSBarry Smith #endif 78c4762a1bSJed Brown PetscMPIInt rank; 79c4762a1bSJed Brown PetscInt eStart, eEnd, vStart, vEnd,j; 80c4762a1bSJed Brown PetscInt genj,loadj; 81c4762a1bSJed Brown Vec X,F; 82c4762a1bSJed Brown Mat J; 83c4762a1bSJed Brown SNES snes; 84c4762a1bSJed Brown 85*327415f7SBarry Smith PetscFunctionBeginUser; 869566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,"poweroptions",help)); 879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 88c4762a1bSJed Brown { 89c4762a1bSJed 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 */ 90c4762a1bSJed Brown /* this is an experiment to see how the analyzer reacts */ 91c4762a1bSJed Brown const PetscMPIInt crank = rank; 92c4762a1bSJed Brown 93c4762a1bSJed Brown /* Create an empty network object */ 949566063dSJacob Faibussowitsch PetscCall(DMNetworkCreate(PETSC_COMM_WORLD,&networkdm)); 95c4762a1bSJed Brown /* Register the components in the network */ 969566063dSJacob Faibussowitsch PetscCall(DMNetworkRegisterComponent(networkdm,"branchstruct",sizeof(struct _p_EDGE_Power),&User.compkey_branch)); 979566063dSJacob Faibussowitsch PetscCall(DMNetworkRegisterComponent(networkdm,"busstruct",sizeof(struct _p_VERTEX_Power),&User.compkey_bus)); 989566063dSJacob Faibussowitsch PetscCall(DMNetworkRegisterComponent(networkdm,"genstruct",sizeof(struct _p_GEN),&User.compkey_gen)); 999566063dSJacob Faibussowitsch PetscCall(DMNetworkRegisterComponent(networkdm,"loadstruct",sizeof(struct _p_LOAD),&User.compkey_load)); 100c4762a1bSJed Brown 1019566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Read Data",&stage1)); 102c4762a1bSJed Brown PetscLogStagePush(stage1); 103c4762a1bSJed Brown /* READ THE DATA */ 104c4762a1bSJed Brown if (!crank) { 105c4762a1bSJed Brown /* READ DATA */ 106c4762a1bSJed Brown /* Only rank 0 reads the data */ 1079566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-pfdata",pfdata_file,sizeof(pfdata_file),NULL)); 1089566063dSJacob Faibussowitsch PetscCall(PetscNew(&pfdata)); 1099566063dSJacob Faibussowitsch PetscCall(PFReadMatPowerData(pfdata,pfdata_file)); 110c4762a1bSJed Brown User.Sbase = pfdata->sbase; 111c4762a1bSJed Brown 112c4762a1bSJed Brown numEdges = pfdata->nbranch; 1139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2*numEdges,&edges)); 1149566063dSJacob Faibussowitsch PetscCall(GetListofEdges_Power(pfdata,edges)); 115c4762a1bSJed Brown } 116c4762a1bSJed Brown 117c4762a1bSJed Brown /* If external option activated. Introduce error in jacobian */ 1189566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL, "-jac_error", &User.jac_error)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown PetscLogStagePop(); 1219566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD)); 1229566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Create network",&stage2)); 123c4762a1bSJed Brown PetscLogStagePush(stage2); 124c4762a1bSJed Brown /* Set number of nodes/edges */ 1259566063dSJacob Faibussowitsch PetscCall(DMNetworkSetNumSubNetworks(networkdm,PETSC_DECIDE,1)); 1269566063dSJacob Faibussowitsch PetscCall(DMNetworkAddSubnetwork(networkdm,"",numEdges,edges,NULL)); 1272bf73ac6SHong Zhang 128c4762a1bSJed Brown /* Set up the network layout */ 1299566063dSJacob Faibussowitsch PetscCall(DMNetworkLayoutSetUp(networkdm)); 130c4762a1bSJed Brown 131c4762a1bSJed Brown if (!crank) { 1329566063dSJacob Faibussowitsch PetscCall(PetscFree(edges)); 133c4762a1bSJed Brown } 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* Add network components only process 0 has any data to add */ 136c4762a1bSJed Brown if (!crank) { 137c4762a1bSJed Brown genj=0; loadj=0; 1389566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd)); 139c4762a1bSJed Brown for (i = eStart; i < eEnd; i++) { 1409566063dSJacob Faibussowitsch PetscCall(DMNetworkAddComponent(networkdm,i,User.compkey_branch,&pfdata->branch[i-eStart],0)); 141c4762a1bSJed Brown } 1429566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(networkdm,&vStart,&vEnd)); 143c4762a1bSJed Brown for (i = vStart; i < vEnd; i++) { 1449566063dSJacob Faibussowitsch PetscCall(DMNetworkAddComponent(networkdm,i,User.compkey_bus,&pfdata->bus[i-vStart],2)); 145c4762a1bSJed Brown if (pfdata->bus[i-vStart].ngen) { 146c4762a1bSJed Brown for (j = 0; j < pfdata->bus[i-vStart].ngen; j++) { 1479566063dSJacob Faibussowitsch PetscCall(DMNetworkAddComponent(networkdm,i,User.compkey_gen,&pfdata->gen[genj++],0)); 148c4762a1bSJed Brown } 149c4762a1bSJed Brown } 150c4762a1bSJed Brown if (pfdata->bus[i-vStart].nload) { 151c4762a1bSJed Brown for (j=0; j < pfdata->bus[i-vStart].nload; j++) { 1529566063dSJacob Faibussowitsch PetscCall(DMNetworkAddComponent(networkdm,i,User.compkey_load,&pfdata->load[loadj++],0)); 153c4762a1bSJed Brown } 154c4762a1bSJed Brown } 155c4762a1bSJed Brown } 156c4762a1bSJed Brown } 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* Set up DM for use */ 1599566063dSJacob Faibussowitsch PetscCall(DMSetUp(networkdm)); 160c4762a1bSJed Brown 161c4762a1bSJed Brown if (!crank) { 1629566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata->bus)); 1639566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata->gen)); 1649566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata->branch)); 1659566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata->load)); 1669566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata)); 167c4762a1bSJed Brown } 168c4762a1bSJed Brown 169c4762a1bSJed Brown /* Distribute networkdm to multiple processes */ 1709566063dSJacob Faibussowitsch PetscCall(DMNetworkDistribute(&networkdm,0)); 171c4762a1bSJed Brown 172c4762a1bSJed Brown PetscLogStagePop(); 1739566063dSJacob Faibussowitsch PetscCall(DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd)); 1749566063dSJacob Faibussowitsch PetscCall(DMNetworkGetVertexRange(networkdm,&vStart,&vEnd)); 175c4762a1bSJed Brown 176c4762a1bSJed Brown #if 0 177c4762a1bSJed Brown EDGE_Power edge; 178c4762a1bSJed Brown PetscInt key,kk,numComponents; 179c4762a1bSJed Brown VERTEX_Power bus; 180c4762a1bSJed Brown GEN gen; 181c4762a1bSJed Brown LOAD load; 182c4762a1bSJed Brown 183c4762a1bSJed Brown for (i = eStart; i < eEnd; i++) { 1849566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm,i,0,&key,(void**)&edge)); 1859566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm,i,&numComponents)); 1869566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Rank %d ncomps = %d Line %d ---- %d\n",crank,numComponents,edge->internal_i,edge->internal_j)); 187c4762a1bSJed Brown } 188c4762a1bSJed Brown 189c4762a1bSJed Brown for (i = vStart; i < vEnd; i++) { 1909566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm,i,&numComponents)); 191c4762a1bSJed Brown for (kk=0; kk < numComponents; kk++) { 1929566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm,i,kk,&key,&component)); 193c4762a1bSJed Brown if (key == 1) { 194c4762a1bSJed Brown bus = (VERTEX_Power)(component); 1959566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Rank %d ncomps = %d Bus %d\n",crank,numComponents,bus->internal_i)); 196c4762a1bSJed Brown } else if (key == 2) { 197c4762a1bSJed Brown gen = (GEN)(component); 19863a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Rank %d Gen pg = %f qg = %f\n",crank,(double)gen->pg,(double)gen->qg)); 199c4762a1bSJed Brown } else if (key == 3) { 200c4762a1bSJed Brown load = (LOAD)(component); 20163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF,"Rank %d Load pl = %f ql = %f\n",crank,(double)load->pl,(double)load->ql)); 202c4762a1bSJed Brown } 203c4762a1bSJed Brown } 204c4762a1bSJed Brown } 205c4762a1bSJed Brown #endif 206c4762a1bSJed Brown /* Broadcast Sbase to all processors */ 2079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&User.Sbase,1,MPIU_SCALAR,0,PETSC_COMM_WORLD)); 208c4762a1bSJed Brown 2099566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(networkdm,&X)); 2109566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X,&F)); 211c4762a1bSJed Brown 2129566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(networkdm,&J)); 2139566063dSJacob Faibussowitsch PetscCall(MatSetOption(J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE)); 214c4762a1bSJed Brown 2159566063dSJacob Faibussowitsch PetscCall(SetInitialValues(networkdm,X,&User)); 216c4762a1bSJed Brown 217c4762a1bSJed Brown /* HOOK UP SOLVER */ 2189566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes)); 2199566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes,networkdm)); 2209566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes,F,FormFunction,&User)); 2219566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes,J,J,FormJacobian_Power,&User)); 2229566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 223c4762a1bSJed Brown 2249566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes,NULL,X)); 2259566063dSJacob Faibussowitsch /* PetscCall(VecView(X,PETSC_VIEWER_STDOUT_WORLD)); */ 226c4762a1bSJed Brown 2279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X)); 2289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&F)); 2299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J)); 230c4762a1bSJed Brown 2319566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes)); 2329566063dSJacob Faibussowitsch PetscCall(DMDestroy(&networkdm)); 233c4762a1bSJed Brown } 2349566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 235b122ec5aSJacob Faibussowitsch return 0; 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 238c4762a1bSJed Brown /*TEST 239c4762a1bSJed Brown 240c4762a1bSJed Brown build: 241c4762a1bSJed Brown depends: PFReadData.c pffunctions.c 242dfd57a17SPierre Jolivet requires: !complex double defined(PETSC_HAVE_ATTRIBUTEALIGNED) 243c4762a1bSJed Brown 244c4762a1bSJed Brown test: 245c4762a1bSJed Brown args: -snes_rtol 1.e-3 246c4762a1bSJed Brown localrunfiles: poweroptions case9.m 247c4762a1bSJed Brown output_file: output/power_1.out 248c4762a1bSJed Brown 249c4762a1bSJed Brown test: 250c4762a1bSJed Brown suffix: 2 251c4762a1bSJed Brown args: -snes_rtol 1.e-3 -petscpartitioner_type simple 252c4762a1bSJed Brown nsize: 4 253c4762a1bSJed Brown localrunfiles: poweroptions case9.m 254c4762a1bSJed Brown output_file: output/power_1.out 255c4762a1bSJed Brown 256c4762a1bSJed Brown TEST*/ 257