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