Actual source code: ex2.c

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1: static char help[] = "This example is based on ex1.c, but generates a random network of chosen sizes and parameters. \n\
  2:   Usage: -n determines number of nodes. The nonnegative seed can be specified with the flag -seed, otherwise the program generates a random seed.\n\n";

  4: /* T
  5:   Concepts: DMNetwork
  6:   Concepts: KSP
  7: */

  9:  #include <petscdmnetwork.h>
 10:  #include <petscksp.h>
 11: #include <time.h>

 13: typedef struct {
 14:   PetscInt    id; /* node id */
 15:   PetscScalar inj; /* current injection (A) */
 16:   PetscBool   gr; /* grounded ? */
 17: } Node;

 19: typedef struct {
 20:   PetscInt    id;  /* branch id */
 21:   PetscScalar r;   /* resistance (ohms) */
 22:   PetscScalar bat; /* battery (V) */
 23: } Branch;

 25: typedef struct Edge {
 26:   PetscInt      n;
 27:   PetscInt      i;
 28:   PetscInt      j;
 29:   struct Edge   *next;
 30: } Edge;

 32: PetscReal distance(PetscReal x1, PetscReal x2, PetscReal y1, PetscReal y2)
 33: {
 34:   return PetscSqrtReal(PetscPowReal(x2-x1,2.0) + PetscPowReal(y2-y1,2.0));
 35: }

 37: /*
 38:   The algorithm for network formation is based on the paper:
 39:   Routing of Multipoint Connections, Bernard M. Waxman. 1988
 40: */

 42: PetscErrorCode random_network(PetscInt nvertex,PetscInt *pnbranch,Node **pnode,Branch **pbranch,int **pedgelist,PetscInt seed)
 43: {
 45:   PetscInt       i, j, nedges = 0;
 46:   int            *edgelist;
 47:   PetscInt       nbat, ncurr, fr, to;
 48:   PetscReal      *x, *y, value, xmax = 10.0; /* generate points in square */
 49:   PetscReal    maxdist = 0.0, dist, alpha, beta, prob;
 50:   PetscRandom    rnd;
 51:   Branch         *branch;
 52:   Node           *node;
 53:   Edge           *head = NULL, *nnew= NULL, *aux= NULL;

 56:   PetscRandomCreate(PETSC_COMM_SELF,&rnd);
 57:   PetscRandomSetFromOptions(rnd);

 59:   PetscRandomSetSeed(rnd, seed);
 60:   PetscRandomSeed(rnd);

 62:   /* These parameters might be modified for experimentation */
 63:   nbat  = (PetscInt)(0.1*nvertex);
 64:   ncurr = (PetscInt)(0.1*nvertex);
 65:   alpha = 0.6;
 66:   beta  = 0.2;

 68:   PetscMalloc2(nvertex,&x,nvertex,&y);

 70:   PetscRandomSetInterval(rnd,0.0,xmax);
 71:   for (i=0; i<nvertex; i++) {
 72:     PetscRandomGetValueReal(rnd,&x[i]);
 73:     PetscRandomGetValueReal(rnd,&y[i]);
 74:   }

 76:   /* find maximum distance */
 77:   for (i=0; i<nvertex; i++) {
 78:     for (j=0; j<nvertex; j++) {
 79:       dist = distance(x[i],x[j],y[i],y[j]);
 80:       if (dist >= maxdist) maxdist = dist;
 81:     }
 82:   }

 84:   PetscRandomSetInterval(rnd,0.0,1.0);
 85:   for (i=0; i<nvertex; i++) {
 86:     for (j=0; j<nvertex; j++) {
 87:       if (j != i) {
 88:         dist = distance(x[i],x[j],y[i],y[j]);
 89:         prob = beta*PetscExpScalar(-dist/(maxdist*alpha));
 90:         PetscRandomGetValueReal(rnd,&value);
 91:         if (value <= prob) {
 92:           PetscMalloc1(1,&nnew);
 93:           if (head == NULL) {
 94:             head       = nnew;
 95:             head->next = NULL;
 96:             head->n    = nedges;
 97:             head->i    = i;
 98:             head->j    = j;
 99:           } else {
100:             aux = head;
101:             head = nnew;
102:             head->n    = nedges;
103:             head->next = aux;
104:             head->i    = i;
105:             head->j    = j;
106:           }
107:           nedges += 1;
108:         }
109:       }
110:     }
111:   }

113:   PetscMalloc1(2*nedges,&edgelist);

115:   for (aux = head; aux; aux = aux->next) {
116:     edgelist[(aux->n)*2]     = aux->i;
117:     edgelist[(aux->n)*2 + 1] = aux->j;
118:   }

120:   aux = head;
121:   while (aux != NULL) {
122:     nnew = aux;
123:     aux = aux->next;
124:     PetscFree(nnew);
125:   }

127:   PetscCalloc2(nvertex,&node,nedges,&branch);
128: 
129:   for (i = 0; i < nvertex; i++) {
130:     node[i].id  = i;
131:     node[i].inj = 0;
132:     node[i].gr = PETSC_FALSE;
133:   }

135:   for (i = 0; i < nedges; i++) {
136:     branch[i].id  = i;
137:     branch[i].r   = 1.0;
138:     branch[i].bat = 0;
139:   }
140: 
141:   /* Chose random node as ground voltage */
142:   PetscRandomSetInterval(rnd,0.0,nvertex);
143:   PetscRandomGetValueReal(rnd,&value);
144:   node[(int)value].gr = PETSC_TRUE;
145: 
146:   /* Create random current and battery injectionsa */
147:   for (i=0; i<ncurr; i++) {
148:     PetscRandomSetInterval(rnd,0.0,nvertex);
149:     PetscRandomGetValueReal(rnd,&value);
150:     fr   = edgelist[(int)value*2];
151:     to   = edgelist[(int)value*2 + 1];
152:     node[fr].inj += 1.0;
153:     node[to].inj -= 1.0;
154:   }

156:   for (i=0; i<nbat; i++) {
157:     PetscRandomSetInterval(rnd,0.0,nedges);
158:     PetscRandomGetValueReal(rnd,&value);
159:     branch[(int)value].bat += 1.0;
160:   }

162:   PetscFree2(x,y);
163:   PetscRandomDestroy(&rnd);

165:   /* assign pointers */
166:   *pnbranch  = nedges;
167:   *pedgelist = edgelist;
168:   *pbranch   = branch;
169:   *pnode     = node;
170:   PetscFunctionReturn(ierr);
171: }

173: PetscErrorCode FormOperator(DM networkdm,Mat A,Vec b)
174: {
175:   PetscErrorCode    ierr;
176:   Vec               localb;
177:   Branch            *branch;
178:   Node              *node;
179:   PetscInt          e,v,vStart,vEnd,eStart, eEnd;
180:   PetscInt          lofst,lofst_to,lofst_fr,compoffset,row[2],col[6];
181:   PetscBool         ghost;
182:   const PetscInt    *cone;
183:   PetscScalar       *barr,val[6];
184:   DMNetworkComponentGenericDataType *arr;

187:   DMGetLocalVector(networkdm,&localb);
188:   VecSet(b,0.0);
189:   VecSet(localb,0.0);
190:   MatZeroEntries(A);

192:   VecGetArray(localb,&barr);

194:   /*
195:     The component data array stores the information which we had in the
196:     node and branch data structures. We access the correct element  with
197:     a variable offset that the DMNetwork provides.
198:   */
199:   DMNetworkGetComponentDataArray(networkdm,&arr);

201:   /*
202:     We can define the current as a "edge characteristic" and the voltage
203:     and the voltage as a "vertex characteristic". With that, we can iterate
204:     the list of edges and vertices, query the associated voltages and currents
205:     and use them to write the Kirchoff equations.
206:   */

208:   /* Branch equations: i/r + uj - ui = battery */
209:   DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);
210:   for (e = 0; e < eEnd; e++) {
211:     DMNetworkGetComponentKeyOffset(networkdm,e,0,NULL,&compoffset);
212:     DMNetworkGetVariableOffset(networkdm,e,&lofst);

214:     DMNetworkGetConnectedVertices(networkdm,e,&cone);
215:     DMNetworkGetVariableOffset(networkdm,cone[0],&lofst_fr);
216:     DMNetworkGetVariableOffset(networkdm,cone[1],&lofst_to);

218:     branch = (Branch*)(arr + compoffset);

220:     barr[lofst] = branch->bat;

222:     row[0] = lofst;
223:     col[0] = lofst;     val[0] =  1;
224:     col[1] = lofst_to;  val[1] =  1;
225:     col[2] = lofst_fr;  val[2] = -1;
226:     MatSetValuesLocal(A,1,row,3,col,val,ADD_VALUES);

228:     /* from node */
229:     DMNetworkGetComponentKeyOffset(networkdm,cone[0],0,NULL,&compoffset);
230:     node = (Node*)(arr + compoffset);

232:     if (!node->gr) {
233:       row[0] = lofst_fr;
234:       col[0] = lofst;   val[0] =  1;
235:       MatSetValuesLocal(A,1,row,1,col,val,ADD_VALUES);
236:     }

238:     /* to node */
239:     DMNetworkGetComponentKeyOffset(networkdm,cone[1],0,NULL,&compoffset);
240:     node = (Node*)(arr + compoffset);

242:     if (!node->gr) {
243:       row[0] = lofst_to;
244:       col[0] = lofst;   val[0] =  -1;
245:       MatSetValuesLocal(A,1,row,1,col,val,ADD_VALUES);
246:     }
247:   }

249:   DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);
250:   for (v = vStart; v < vEnd; v++) {
251:     DMNetworkIsGhostVertex(networkdm,v,&ghost);
252:     if (!ghost) {
253:       DMNetworkGetComponentKeyOffset(networkdm,v,0,NULL,&compoffset);
254:       DMNetworkGetVariableOffset(networkdm,v,&lofst);
255:       node = (Node*)(arr + compoffset);

257:       if (node->gr) {
258:         row[0] = lofst;
259:         col[0] = lofst;   val[0] =  1;
260:         MatSetValuesLocal(A,1,row,1,col,val,ADD_VALUES);
261:       } else {
262:         barr[lofst] -= node->inj;
263:       }
264:     }
265:   }

267:   VecRestoreArray(localb,&barr);

269:   DMLocalToGlobalBegin(networkdm,localb,ADD_VALUES,b);
270:   DMLocalToGlobalEnd(networkdm,localb,ADD_VALUES,b);
271:   DMRestoreLocalVector(networkdm,&localb);

273:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
274:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
275:   return(0);
276: }

278: int main(int argc,char ** argv)
279: {
280:   PetscErrorCode    ierr;
281:   PetscInt          i, nbranch = 0, eStart, eEnd, vStart, vEnd;
282:   PetscInt          seed = 0, nnode = 0;
283:   PetscMPIInt       size, rank;
284:   DM                networkdm;
285:   Vec               x, b;
286:   Mat               A;
287:   KSP               ksp;
288:   int               *edgelist = NULL;
289:   PetscInt          componentkey[2];
290:   Node              *node;
291:   Branch            *branch;
292: #if defined(PETSC_USE_LOG)
293:   PetscLogStage stage[3];
294: #endif

296:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
297:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
298:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

300:   PetscOptionsGetInt(NULL,NULL,"-seed",&seed,NULL);

302:   PetscLogStageRegister("Network Creation", &stage[0]);
303:   PetscLogStageRegister("DMNetwork data structures", &stage[1]);
304:   PetscLogStageRegister("KSP", &stage[2]);

306:   PetscLogStagePush(stage[0]);
307:   /* "read" data only for processor 0 */
308:   if (!rank) {
309:     nnode = 100;
310:     PetscOptionsGetInt(NULL,NULL,"-n",&nnode,NULL);
311:     random_network(nnode, &nbranch, &node, &branch, &edgelist, seed);
312:   }
313:   PetscLogStagePop();

315:   PetscLogStagePush(stage[1]);
316:   DMNetworkCreate(PETSC_COMM_WORLD,&networkdm);
317:   DMNetworkRegisterComponent(networkdm,"nstr",sizeof(Node),&componentkey[0]);
318:   DMNetworkRegisterComponent(networkdm,"bsrt",sizeof(Branch),&componentkey[1]);

320:   /* Set number of nodes/edges */
321:   DMNetworkSetSizes(networkdm,nnode,nbranch,PETSC_DETERMINE,PETSC_DETERMINE);
322:   /* Add edge connectivity */
323:   DMNetworkSetEdgeList(networkdm,edgelist);
324:   /* Set up the network layout */
325:   DMNetworkLayoutSetUp(networkdm);

327:   /* Add network components: physical parameters of nodes and branches*/
328:   if (!rank) {
329:     DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);
330:     for (i = eStart; i < eEnd; i++) {
331:       DMNetworkAddComponent(networkdm,i,componentkey[1],&branch[i-eStart]);
332:       DMNetworkAddNumVariables(networkdm,i,1);
333:     }

335:     DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);
336:     for (i = vStart; i < vEnd; i++) {
337:       DMNetworkAddComponent(networkdm,i,componentkey[0],&node[i-vStart]);
338:       /* Add number of variables */
339:       DMNetworkAddNumVariables(networkdm,i,1);
340:     }
341:   }

343:   /* Network partitioning and distribution of data */
344:   DMSetUp(networkdm);
345:   DMNetworkDistribute(&networkdm,0);
346:   DMNetworkAssembleGraphStructures(networkdm);

348:   /* We don't use these data structures anymore since they have been copied to networkdm */
349:   if (!rank) {
350:     PetscFree(edgelist);
351:     PetscFree2(node,branch);
352:   }

354:   /* Create vectors and matrix */
355:   DMCreateGlobalVector(networkdm,&x);
356:   VecDuplicate(x,&b);
357:   DMCreateMatrix(networkdm,&A);

359:   PetscLogStagePop();

361:   PetscLogStagePush(stage[2]);
362:   /* Assembly system of equations */
363:   FormOperator(networkdm,A,b);

365:   /* Solve linear system: A x = b */
366:   KSPCreate(PETSC_COMM_WORLD, &ksp);
367:   KSPSetOperators(ksp, A, A);
368:   KSPSetFromOptions(ksp);
369:   KSPSolve(ksp, b, x);

371:   PetscLogStagePop();
372: 
373:   /* Free work space */
374:   VecDestroy(&x);
375:   VecDestroy(&b);
376:   MatDestroy(&A);
377:   KSPDestroy(&ksp);
378:   DMDestroy(&networkdm);
379:   PetscFinalize();
380:   return ierr;
381: }