Actual source code: power.c

petsc-3.10.3 2018-12-18
Report Typos and Errors
  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:                       The data file format used is from the MatPower package (http://www.pserc.cornell.edu//matpower/).\n\
  4:                       Run this program: mpiexec -n <n> ./pf\n\
  5:                       mpiexec -n <n> ./pfc \n";

  7: /* T
  8:    Concepts: DMNetwork
  9:    Concepts: PETSc SNES solver
 10: */

 12: #include "power.h"
 13:  #include <petscdmnetwork.h>

 15: PetscErrorCode FormFunction(SNES snes,Vec X, Vec F,void *appctx)
 16: {
 18:   DM             networkdm;
 19:   UserCtx_Power  *User=(UserCtx_Power*)appctx;
 20:   Vec            localX,localF;
 21:   PetscInt       nv,ne;
 22:   const PetscInt *vtx,*edges;

 25:   SNESGetDM(snes,&networkdm);
 26:   DMGetLocalVector(networkdm,&localX);
 27:   DMGetLocalVector(networkdm,&localF);
 28:   VecSet(F,0.0);
 29:   VecSet(localF,0.0);

 31:   DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);
 32:   DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);

 34:   DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);
 35:   FormFunction_Power(networkdm,localX,localF,nv,ne,vtx,edges,User);

 37:   DMRestoreLocalVector(networkdm,&localX);

 39:   DMLocalToGlobalBegin(networkdm,localF,ADD_VALUES,F);
 40:   DMLocalToGlobalEnd(networkdm,localF,ADD_VALUES,F);
 41:   DMRestoreLocalVector(networkdm,&localF);
 42:   return(0);
 43: }

 45: PetscErrorCode SetInitialValues(DM networkdm,Vec X,void* appctx)
 46: {
 48:   PetscInt       vStart,vEnd,nv,ne;
 49:   const PetscInt *vtx,*edges;
 50:   Vec            localX;
 51:   UserCtx_Power  *user_power=(UserCtx_Power*)appctx;

 54:   DMNetworkGetVertexRange(networkdm,&vStart, &vEnd);

 56:   DMGetLocalVector(networkdm,&localX);

 58:   VecSet(X,0.0);
 59:   DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);
 60:   DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);

 62:   DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);
 63:   SetInitialGuess_Power(networkdm,localX,nv,ne,vtx,edges,user_power);

 65:   DMLocalToGlobalBegin(networkdm,localX,ADD_VALUES,X);
 66:   DMLocalToGlobalEnd(networkdm,localX,ADD_VALUES,X);
 67:   DMRestoreLocalVector(networkdm,&localX);
 68:   return(0);
 69: }

 71: int main(int argc,char ** argv)
 72: {
 73:   PetscErrorCode   ierr;
 74:   char             pfdata_file[PETSC_MAX_PATH_LEN]="case9.m";
 75:   PFDATA           *pfdata;
 76:   PetscInt         numEdges=0,numVertices=0,NumEdges=PETSC_DETERMINE,NumVertices=PETSC_DETERMINE;
 77:   PetscInt         *edges = NULL;
 78:   PetscInt         i;
 79:   DM               networkdm;
 80:   UserCtx_Power    User;
 81:   PetscLogStage    stage1,stage2;
 82:   PetscMPIInt      rank;
 83:   PetscInt         eStart, eEnd, vStart, vEnd,j;
 84:   PetscInt         genj,loadj;
 85:   Vec              X,F;
 86:   Mat              J;
 87:   SNES             snes;

 89:   PetscInitialize(&argc,&argv,"poweroptions",help);
 90:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 91:   {
 92:     /* 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 */
 93:     /* this is an experiment to see how the analyzer reacts */
 94:     const PetscMPIInt crank = rank;

 96:     /* Create an empty network object */
 97:     DMNetworkCreate(PETSC_COMM_WORLD,&networkdm);
 98:     /* Register the components in the network */
 99:     DMNetworkRegisterComponent(networkdm,"branchstruct",sizeof(struct _p_EDGE_Power),&User.compkey_branch);
100:     DMNetworkRegisterComponent(networkdm,"busstruct",sizeof(struct _p_VERTEX_Power),&User.compkey_bus);
101:     DMNetworkRegisterComponent(networkdm,"genstruct",sizeof(struct _p_GEN),&User.compkey_gen);
102:     DMNetworkRegisterComponent(networkdm,"loadstruct",sizeof(struct _p_LOAD),&User.compkey_load);

104:     PetscLogStageRegister("Read Data",&stage1);
105:     PetscLogStagePush(stage1);
106:     /* READ THE DATA */
107:     if (!crank) {
108:       /*    READ DATA */
109:       /* Only rank 0 reads the data */
110:       PetscOptionsGetString(NULL,NULL,"-pfdata",pfdata_file,PETSC_MAX_PATH_LEN-1,NULL);
111:       PetscNew(&pfdata);
112:       PFReadMatPowerData(pfdata,pfdata_file);
113:       User.Sbase = pfdata->sbase;

115:       numEdges = pfdata->nbranch;
116:       numVertices = pfdata->nbus;

118:       PetscMalloc1(2*numEdges,&edges);
119:       GetListofEdges_Power(pfdata,edges);
120:     }

122:     /* If external option activated. Introduce error in jacobian */
123:     PetscOptionsHasName(NULL,NULL, "-jac_error", &User.jac_error);

125:     PetscLogStagePop();
126:     MPI_Barrier(PETSC_COMM_WORLD);
127:     PetscLogStageRegister("Create network",&stage2);
128:     PetscLogStagePush(stage2);
129:     /* Set number of nodes/edges */
130:     DMNetworkSetSizes(networkdm,1,0,&numVertices,&numEdges,&NumVertices,&NumEdges);
131:     /* Add edge connectivity */
132:     DMNetworkSetEdgeList(networkdm,&edges,NULL);
133:     /* Set up the network layout */
134:     DMNetworkLayoutSetUp(networkdm);

136:     if (!crank) {
137:       PetscFree(edges);
138:     }

140:     /* Add network components only process 0 has any data to add*/
141:     if (!crank) {
142:       genj=0; loadj=0;
143:       DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);
144:       for (i = eStart; i < eEnd; i++) {
145:         DMNetworkAddComponent(networkdm,i,User.compkey_branch,&pfdata->branch[i-eStart]);
146:       }
147:       DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);
148:       for (i = vStart; i < vEnd; i++) {
149:         DMNetworkAddComponent(networkdm,i,User.compkey_bus,&pfdata->bus[i-vStart]);
150:         if (pfdata->bus[i-vStart].ngen) {
151:           for (j = 0; j < pfdata->bus[i-vStart].ngen; j++) {
152:             DMNetworkAddComponent(networkdm,i,User.compkey_gen,&pfdata->gen[genj++]);
153:           }
154:         }
155:         if (pfdata->bus[i-vStart].nload) {
156:           for (j=0; j < pfdata->bus[i-vStart].nload; j++) {
157:             DMNetworkAddComponent(networkdm,i,User.compkey_load,&pfdata->load[loadj++]);
158:           }
159:         }
160:         /* Add number of variables */
161:         DMNetworkAddNumVariables(networkdm,i,2);
162:       }
163:     }

165:     /* Set up DM for use */
166:     DMSetUp(networkdm);

168:     if (!crank) {
169:       PetscFree(pfdata->bus);
170:       PetscFree(pfdata->gen);
171:       PetscFree(pfdata->branch);
172:       PetscFree(pfdata->load);
173:       PetscFree(pfdata);
174:     }

176:     /* Distribute networkdm to multiple processes */
177:     DMNetworkDistribute(&networkdm,0);

179:     PetscLogStagePop();
180:     DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);
181:     DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);

183: #if 0
184:     EDGE_Power     edge;
185:     PetscInt       key,kk,numComponents;
186:     VERTEX_Power   bus;
187:     GEN            gen;
188:     LOAD           load;

190:     for (i = eStart; i < eEnd; i++) {
191:       DMNetworkGetComponent(networkdm,i,0,&key,(void**)&edge);
192:       DMNetworkGetNumComponents(networkdm,i,&numComponents);
193:       PetscPrintf(PETSC_COMM_SELF,"Rank %d ncomps = %d Line %d ---- %d\n",crank,numComponents,edge->internal_i,edge->internal_j);
194:     }

196:     for (i = vStart; i < vEnd; i++) {
197:       DMNetworkGetNumComponents(networkdm,i,&numComponents);
198:       for (kk=0; kk < numComponents; kk++) {
199:         DMNetworkGetComponent(networkdm,i,kk,&key,&component);
200:         if (key == 1) {
201:           bus = (VERTEX_Power)(component);
202:           PetscPrintf(PETSC_COMM_SELF,"Rank %d ncomps = %d Bus %d\n",crank,numComponents,bus->internal_i);
203:         } else if (key == 2) {
204:           gen = (GEN)(component);
205:           PetscPrintf(PETSC_COMM_SELF,"Rank %d Gen pg = %f qg = %f\n",crank,gen->pg,gen->qg);
206:         } else if (key == 3) {
207:           load = (LOAD)(component);
208:           PetscPrintf(PETSC_COMM_SELF,"Rank %d Load pl = %f ql = %f\n",crank,load->pl,load->ql);
209:         }
210:       }
211:     }
212: #endif
213:     /* Broadcast Sbase to all processors */
214:     MPI_Bcast(&User.Sbase,1,MPIU_SCALAR,0,PETSC_COMM_WORLD);

216:     DMCreateGlobalVector(networkdm,&X);
217:     VecDuplicate(X,&F);

219:     DMCreateMatrix(networkdm,&J);
220:     MatSetOption(J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);

222:     SetInitialValues(networkdm,X,&User);

224:     /* HOOK UP SOLVER */
225:     SNESCreate(PETSC_COMM_WORLD,&snes);
226:     SNESSetDM(snes,networkdm);
227:     SNESSetFunction(snes,F,FormFunction,&User);
228:     SNESSetJacobian(snes,J,J,FormJacobian_Power,&User);
229:     SNESSetFromOptions(snes);

231:     SNESSolve(snes,NULL,X);
232:     /* VecView(X,PETSC_VIEWER_STDOUT_WORLD); */

234:     VecDestroy(&X);
235:     VecDestroy(&F);
236:     MatDestroy(&J);

238:     SNESDestroy(&snes);
239:     DMDestroy(&networkdm);
240:   }
241:   PetscFinalize();
242:   return ierr;
243: }

245: /*TEST

247:    build:
248:      depends: PFReadData.c pffunctions.c
249:      requires: !complex double define(PETSC_HAVE_ATTRIBUTEALIGNED)


252:    test:
253:      args: -snes_rtol 1.e-3
254:      localrunfiles: poweroptions case9.m
255:      output_file: output/power_1.out
256:      requires: double !complex define(PETSC_HAVE_ATTRIBUTEALIGNED)

258:    test:
259:      suffix: 2
260:      args: -snes_rtol 1.e-3 -petscpartitioner_type simple
261:      nsize: 4
262:      localrunfiles: poweroptions case9.m
263:      output_file: output/power_1.out
264:      requires: double !complex define(PETSC_HAVE_ATTRIBUTEALIGNED)

266: TEST*/