xref: /petsc/src/snes/tutorials/network/power/power.c (revision d0e6bf2ad94dcc89b258ce16c7987200a4714786)
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 
FormFunction(SNES snes,Vec X,Vec F,void * appctx)12d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *appctx)
13d71ae5a4SJacob Faibussowitsch {
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));
383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39c4762a1bSJed Brown }
40c4762a1bSJed Brown 
SetInitialValues(DM networkdm,Vec X,void * appctx)41d71ae5a4SJacob Faibussowitsch PetscErrorCode SetInitialValues(DM networkdm, Vec X, void *appctx)
42d71ae5a4SJacob Faibussowitsch {
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));
633ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
main(int argc,char ** argv)66d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
67d71ae5a4SJacob Faibussowitsch {
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;
75c4762a1bSJed Brown   PetscLogStage stage1, stage2;
76c4762a1bSJed Brown   PetscMPIInt   rank;
77c4762a1bSJed Brown   PetscInt      eStart, eEnd, vStart, vEnd, j;
78c4762a1bSJed Brown   PetscInt      genj, loadj;
79c4762a1bSJed Brown   Vec           X, F;
80c4762a1bSJed Brown   Mat           J;
81c4762a1bSJed Brown   SNES          snes;
82c4762a1bSJed Brown 
83327415f7SBarry Smith   PetscFunctionBeginUser;
849566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, "poweroptions", help));
859566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
86c4762a1bSJed Brown   {
87c4762a1bSJed 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 */
88c4762a1bSJed Brown     /* this is an experiment to see how the analyzer reacts */
89c4762a1bSJed Brown     const PetscMPIInt crank = rank;
90c4762a1bSJed Brown 
91c4762a1bSJed Brown     /* Create an empty network object */
929566063dSJacob Faibussowitsch     PetscCall(DMNetworkCreate(PETSC_COMM_WORLD, &networkdm));
93c4762a1bSJed Brown     /* Register the components in the network */
949566063dSJacob Faibussowitsch     PetscCall(DMNetworkRegisterComponent(networkdm, "branchstruct", sizeof(struct _p_EDGE_Power), &User.compkey_branch));
959566063dSJacob Faibussowitsch     PetscCall(DMNetworkRegisterComponent(networkdm, "busstruct", sizeof(struct _p_VERTEX_Power), &User.compkey_bus));
969566063dSJacob Faibussowitsch     PetscCall(DMNetworkRegisterComponent(networkdm, "genstruct", sizeof(struct _p_GEN), &User.compkey_gen));
979566063dSJacob Faibussowitsch     PetscCall(DMNetworkRegisterComponent(networkdm, "loadstruct", sizeof(struct _p_LOAD), &User.compkey_load));
98c4762a1bSJed Brown 
999566063dSJacob Faibussowitsch     PetscCall(PetscLogStageRegister("Read Data", &stage1));
1003ba16761SJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage1));
101c4762a1bSJed Brown     /* READ THE DATA */
102c4762a1bSJed Brown     if (!crank) {
103c4762a1bSJed Brown       /*    READ DATA */
104c4762a1bSJed Brown       /* Only rank 0 reads the data */
1059566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL, NULL, "-pfdata", pfdata_file, sizeof(pfdata_file), NULL));
1069566063dSJacob Faibussowitsch       PetscCall(PetscNew(&pfdata));
1079566063dSJacob Faibussowitsch       PetscCall(PFReadMatPowerData(pfdata, pfdata_file));
108c4762a1bSJed Brown       User.Sbase = pfdata->sbase;
109c4762a1bSJed Brown 
110c4762a1bSJed Brown       numEdges = pfdata->nbranch;
1119566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(2 * numEdges, &edges));
1129566063dSJacob Faibussowitsch       PetscCall(GetListofEdges_Power(pfdata, edges));
113c4762a1bSJed Brown     }
114c4762a1bSJed Brown 
115c4762a1bSJed Brown     /* If external option activated. Introduce error in jacobian */
1169566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(NULL, NULL, "-jac_error", &User.jac_error));
117c4762a1bSJed Brown 
1183ba16761SJacob Faibussowitsch     PetscCall(PetscLogStagePop());
1199566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
1209566063dSJacob Faibussowitsch     PetscCall(PetscLogStageRegister("Create network", &stage2));
1213ba16761SJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage2));
122c4762a1bSJed Brown     /* Set number of nodes/edges */
1239566063dSJacob Faibussowitsch     PetscCall(DMNetworkSetNumSubNetworks(networkdm, PETSC_DECIDE, 1));
1249566063dSJacob Faibussowitsch     PetscCall(DMNetworkAddSubnetwork(networkdm, "", numEdges, edges, NULL));
1252bf73ac6SHong Zhang 
126c4762a1bSJed Brown     /* Set up the network layout */
1279566063dSJacob Faibussowitsch     PetscCall(DMNetworkLayoutSetUp(networkdm));
128c4762a1bSJed Brown 
12948a46eb9SPierre Jolivet     if (!crank) PetscCall(PetscFree(edges));
130c4762a1bSJed Brown 
131c4762a1bSJed Brown     /* Add network components only process 0 has any data to add */
132c4762a1bSJed Brown     if (!crank) {
1339371c9d4SSatish Balay       genj  = 0;
1349371c9d4SSatish Balay       loadj = 0;
1359566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetEdgeRange(networkdm, &eStart, &eEnd));
13648a46eb9SPierre Jolivet       for (i = eStart; i < eEnd; i++) PetscCall(DMNetworkAddComponent(networkdm, i, User.compkey_branch, &pfdata->branch[i - eStart], 0));
1379566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd));
138c4762a1bSJed Brown       for (i = vStart; i < vEnd; i++) {
1399566063dSJacob Faibussowitsch         PetscCall(DMNetworkAddComponent(networkdm, i, User.compkey_bus, &pfdata->bus[i - vStart], 2));
140c4762a1bSJed Brown         if (pfdata->bus[i - vStart].ngen) {
14148a46eb9SPierre Jolivet           for (j = 0; j < pfdata->bus[i - vStart].ngen; j++) PetscCall(DMNetworkAddComponent(networkdm, i, User.compkey_gen, &pfdata->gen[genj++], 0));
142c4762a1bSJed Brown         }
143c4762a1bSJed Brown         if (pfdata->bus[i - vStart].nload) {
14448a46eb9SPierre Jolivet           for (j = 0; j < pfdata->bus[i - vStart].nload; j++) PetscCall(DMNetworkAddComponent(networkdm, i, User.compkey_load, &pfdata->load[loadj++], 0));
145c4762a1bSJed Brown         }
146c4762a1bSJed Brown       }
147c4762a1bSJed Brown     }
148c4762a1bSJed Brown 
149c4762a1bSJed Brown     /* Set up DM for use */
1509566063dSJacob Faibussowitsch     PetscCall(DMSetUp(networkdm));
151c4762a1bSJed Brown 
152c4762a1bSJed Brown     if (!crank) {
1539566063dSJacob Faibussowitsch       PetscCall(PetscFree(pfdata->bus));
1549566063dSJacob Faibussowitsch       PetscCall(PetscFree(pfdata->gen));
1559566063dSJacob Faibussowitsch       PetscCall(PetscFree(pfdata->branch));
1569566063dSJacob Faibussowitsch       PetscCall(PetscFree(pfdata->load));
1579566063dSJacob Faibussowitsch       PetscCall(PetscFree(pfdata));
158c4762a1bSJed Brown     }
159c4762a1bSJed Brown 
160c4762a1bSJed Brown     /* Distribute networkdm to multiple processes */
1619566063dSJacob Faibussowitsch     PetscCall(DMNetworkDistribute(&networkdm, 0));
162c4762a1bSJed Brown 
1633ba16761SJacob Faibussowitsch     PetscCall(PetscLogStagePop());
1649566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetEdgeRange(networkdm, &eStart, &eEnd));
1659566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown #if 0
168c4762a1bSJed Brown     EDGE_Power     edge;
169c4762a1bSJed Brown     PetscInt       key,kk,numComponents;
170c4762a1bSJed Brown     VERTEX_Power   bus;
171c4762a1bSJed Brown     GEN            gen;
172c4762a1bSJed Brown     LOAD           load;
173c4762a1bSJed Brown 
174c4762a1bSJed Brown     for (i = eStart; i < eEnd; i++) {
1759566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetComponent(networkdm,i,0,&key,(void**)&edge));
1769566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetNumComponents(networkdm,i,&numComponents));
1779566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_SELF,"Rank %d ncomps = %d Line %d ---- %d\n",crank,numComponents,edge->internal_i,edge->internal_j));
178c4762a1bSJed Brown     }
179c4762a1bSJed Brown 
180c4762a1bSJed Brown     for (i = vStart; i < vEnd; i++) {
1819566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetNumComponents(networkdm,i,&numComponents));
182c4762a1bSJed Brown       for (kk=0; kk < numComponents; kk++) {
1839566063dSJacob Faibussowitsch         PetscCall(DMNetworkGetComponent(networkdm,i,kk,&key,&component));
184c4762a1bSJed Brown         if (key == 1) {
185*57508eceSPierre Jolivet           bus = (VERTEX_Power)component;
1869566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF,"Rank %d ncomps = %d Bus %d\n",crank,numComponents,bus->internal_i));
187c4762a1bSJed Brown         } else if (key == 2) {
188*57508eceSPierre Jolivet           gen = (GEN)component;
18963a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF,"Rank %d Gen pg = %f qg = %f\n",crank,(double)gen->pg,(double)gen->qg));
190c4762a1bSJed Brown         } else if (key == 3) {
191*57508eceSPierre Jolivet           load = (LOAD)component;
19263a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF,"Rank %d Load pl = %f ql = %f\n",crank,(double)load->pl,(double)load->ql));
193c4762a1bSJed Brown         }
194c4762a1bSJed Brown       }
195c4762a1bSJed Brown     }
196c4762a1bSJed Brown #endif
197c4762a1bSJed Brown     /* Broadcast Sbase to all processors */
1989566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&User.Sbase, 1, MPIU_SCALAR, 0, PETSC_COMM_WORLD));
199c4762a1bSJed Brown 
2009566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(networkdm, &X));
2019566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X, &F));
202c4762a1bSJed Brown 
2039566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(networkdm, &J));
2049566063dSJacob Faibussowitsch     PetscCall(MatSetOption(J, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
205c4762a1bSJed Brown 
2069566063dSJacob Faibussowitsch     PetscCall(SetInitialValues(networkdm, X, &User));
207c4762a1bSJed Brown 
208c4762a1bSJed Brown     /* HOOK UP SOLVER */
2099566063dSJacob Faibussowitsch     PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
2109566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(snes, networkdm));
2119566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes, F, FormFunction, &User));
2129566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes, J, J, FormJacobian_Power, &User));
2139566063dSJacob Faibussowitsch     PetscCall(SNESSetFromOptions(snes));
214c4762a1bSJed Brown 
2159566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes, NULL, X));
2169566063dSJacob Faibussowitsch     /* PetscCall(VecView(X,PETSC_VIEWER_STDOUT_WORLD)); */
217c4762a1bSJed Brown 
2189566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&X));
2199566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&F));
2209566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&J));
221c4762a1bSJed Brown 
2229566063dSJacob Faibussowitsch     PetscCall(SNESDestroy(&snes));
2239566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&networkdm));
224c4762a1bSJed Brown   }
2259566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
226b122ec5aSJacob Faibussowitsch   return 0;
227c4762a1bSJed Brown }
228c4762a1bSJed Brown 
229c4762a1bSJed Brown /*TEST
230c4762a1bSJed Brown 
231c4762a1bSJed Brown    build:
232c4762a1bSJed Brown      depends: PFReadData.c pffunctions.c
233dfd57a17SPierre Jolivet      requires: !complex double defined(PETSC_HAVE_ATTRIBUTEALIGNED)
234c4762a1bSJed Brown 
235c4762a1bSJed Brown    test:
236c4762a1bSJed Brown      args: -snes_rtol 1.e-3
237c4762a1bSJed Brown      localrunfiles: poweroptions case9.m
238c4762a1bSJed Brown      output_file: output/power_1.out
239c4762a1bSJed Brown 
240c4762a1bSJed Brown    test:
241c4762a1bSJed Brown      suffix: 2
242c4762a1bSJed Brown      args: -snes_rtol 1.e-3 -petscpartitioner_type simple
243c4762a1bSJed Brown      nsize: 4
244c4762a1bSJed Brown      localrunfiles: poweroptions case9.m
245c4762a1bSJed Brown      output_file: output/power_1.out
246c4762a1bSJed Brown 
247c4762a1bSJed Brown TEST*/
248