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