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