Actual source code: ex1_nest.c

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1: static char help[] = "This example is based on ex1 using MATNEST. \n";


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

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

 12: /* The topology looks like:

 14:             (1)
 15:             /|\
 16:            / | \
 17:           /  |  \
 18:          R   R   V
 19:         /    |b4  \
 20:     b1 /    (4)    \ b2
 21:       /    /   \    R
 22:      /   R       R   \
 23:     /  /           \  \
 24:    / / b5        b6  \ \
 25:   //                   \\
 26: (2)--------- R -------- (3)
 27:              b3

 29:   Nodes:          (1), ... (4)
 30:   Branches:       b1, ... b6
 31:   Resistances:    R
 32:   Voltage source: V

 34:   Additionally, there is a current source from (2) to (1).
 35: */

 37: /*
 38:   Structures containing physical data of circuit.
 39:   Note that no topology is defined
 40: */

 42: typedef struct {
 43:   PetscInt    id;  /* node id */
 44:   PetscScalar inj; /* current injection (A) */
 45:   PetscBool   gr;  /* grounded ? */
 46: } Node;

 48: typedef struct {
 49:   PetscInt    id;  /* branch id */
 50:   PetscScalar r;   /* resistance (ohms) */
 51:   PetscScalar bat; /* battery (V) */
 52: } Branch;

 54: /*
 55:   read_data: this routine fills data structures with problem data.
 56:   This can be substituted by an external parser.
 57: */

 59: PetscErrorCode read_data(PetscInt *pnnode,PetscInt *pnbranch,Node **pnode,Branch **pbranch,int **pedgelist)
 60: {
 61:   PetscErrorCode    ierr;
 62:   PetscInt          nnode, nbranch, i;
 63:   Branch            *branch;
 64:   Node              *node;
 65:   int               *edgelist;

 67:   nnode   = 4;
 68:   nbranch = 6;

 70:   PetscCalloc1(nnode*sizeof(Node),&node);
 71:   PetscCalloc1(nbranch*sizeof(Branch),&branch);

 73:   for (i = 0; i < nnode; i++) {
 74:     node[i].id  = i;
 75:     node[i].inj = 0;
 76:     node[i].gr = PETSC_FALSE;
 77:   }

 79:   for (i = 0; i < nbranch; i++) {
 80:     branch[i].id  = i;
 81:     branch[i].r   = 1.0;
 82:     branch[i].bat = 0;
 83:   }

 85:   /*
 86:     Branch 1 contains a voltage source of 12.0 V
 87:     From node 0 to 1 there exists a current source of 4.0 A
 88:     Node 3 is grounded, hence the voltage is 0.
 89:   */
 90:   branch[1].bat = 12.0;
 91:   node[0].inj   = -4.0;
 92:   node[1].inj   =  4.0;
 93:   node[3].gr    =  PETSC_TRUE;

 95:   /*
 96:     edgelist is an array of len = 2*nbranch. that defines the
 97:     topology of the network. For branch i we would have that:
 98:       edgelist[2*i]     = from node
 99:       edgelist[2*i + 1] = to node
100:   */

102:   PetscCalloc1(2*nbranch*sizeof(int), &edgelist);

104:   for (i = 0; i < nbranch; i++) {
105:     switch (i) {
106:       case 0:
107:         edgelist[2*i]     = 0;
108:         edgelist[2*i + 1] = 1;
109:         break;
110:       case 1:
111:         edgelist[2*i]     = 0;
112:         edgelist[2*i + 1] = 2;
113:         break;
114:       case 2:
115:         edgelist[2*i]     = 1;
116:         edgelist[2*i + 1] = 2;
117:         break;
118:       case 3:
119:         edgelist[2*i]     = 0;
120:         edgelist[2*i + 1] = 3;
121:         break;
122:       case 4:
123:         edgelist[2*i]     = 1;
124:         edgelist[2*i + 1] = 3;
125:         break;
126:       case 5:
127:         edgelist[2*i]     = 2;
128:         edgelist[2*i + 1] = 3;
129:         break;
130:       default:
131:         break;
132:     }
133:   }

135:   /* assign pointers */
136:   *pnnode    = nnode;
137:   *pnbranch  = nbranch;
138:   *pedgelist = edgelist;
139:   *pbranch   = branch;
140:   *pnode     = node;
141:   return(0);
142: }

144: PetscErrorCode FormOperator(DM networkdm,Mat A,Vec b)
145: {
146:   PetscErrorCode    ierr;
147:   Vec               localb;
148:   Branch            *branch;
149:   Node              *node;
150:   PetscInt          e;
151:   PetscInt          v,vStart,vEnd;
152:   PetscInt          eStart, eEnd;
153:   PetscBool         ghost;
154:   const PetscInt    *cone;
155:   PetscScalar       *barr;
156:   PetscInt          lofst, lofst_to, lofst_fr;
157:   PetscInt          compoffset, key;
158:   PetscInt          row[2],col[6];
159:   PetscScalar       val[6];
160:   Mat               e11, c12, c21, v22;
161:   Mat               **mats;
162:   DMNetworkComponentGenericDataType *arr;

165:   DMGetLocalVector(networkdm,&localb);
166:   VecSet(b,0.0);
167:   VecSet(localb,0.0);

169:   VecGetArray(localb,&barr);

171:   /* Get nested submatrices */
172:   MatNestGetSubMats(A,NULL,NULL,&mats);
173:   e11 = mats[0][0];  /* edges */
174:   c12 = mats[0][1];  /* couplings */
175:   c21 = mats[1][0];  /* couplings */
176:   v22 = mats[1][1];  /* vertices */

178:   DMNetworkGetComponentDataArray(networkdm,&arr);

180:   /* Get vertices and edge range */
181:   DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);
182:   DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);

184:   for (e = 0; e < eEnd; e++) {
185:     DMNetworkGetComponentKeyOffset(networkdm,e,0,&key,&compoffset);
186:     DMNetworkGetEdgeOffset(networkdm,e,&lofst);

188:     DMNetworkGetConnectedVertices(networkdm,e,&cone);
189:     DMNetworkGetVertexOffset(networkdm,cone[0],&lofst_fr);
190:     DMNetworkGetVertexOffset(networkdm,cone[1],&lofst_to);
191:     branch = (Branch*)(arr + compoffset);

193:     /* These are edge-edge and go to e11 */
194:     row[0] = lofst;
195:     col[0] = lofst;     val[0] =  1;
196:     MatSetValuesLocal(e11,1,row,1,col,val,INSERT_VALUES);

198:     /* These are edge-vertex and go to c12 */
199:     col[0] = lofst_to;  val[0] =  1;
200:     col[1] = lofst_fr;  val[1] = -1;
201:     MatSetValuesLocal(c12,1,row,2,col,val,INSERT_VALUES);

203:     /* These are edge-vertex and go to c12 */
204:     /* from node */
205:     DMNetworkGetComponentKeyOffset(networkdm,cone[0],0,&key,&compoffset);
206:     node = (Node*)(arr + compoffset);

208:     if (!node->gr) {
209:       row[0] = lofst_fr;
210:       col[0] = lofst;   val[0] =  1;
211:       MatSetValuesLocal(c21,1,row,1,col,val,INSERT_VALUES);
212:     }

214:     /* to node */
215:     DMNetworkGetComponentKeyOffset(networkdm,cone[1],0,&key,&compoffset);
216:     node = (Node*)(arr + compoffset);

218:     if (!node->gr) {
219:       row[0] = lofst_to;
220:       col[0] = lofst;   val[0] =  -1;
221:       MatSetValuesLocal(c21,1,row,1,col,val,INSERT_VALUES);
222:     }

224:     /* TODO: this is not a nested vector. Need to implement nested vector */
225:     DMNetworkGetVariableOffset(networkdm,e,&lofst);
226:     barr[lofst] = branch->bat;
227:   }

229:   for (v = vStart; v < vEnd; v++) {
230:     DMNetworkIsGhostVertex(networkdm,v,&ghost);
231:     if (!ghost) {
232:       DMNetworkGetComponentKeyOffset(networkdm,v,0,&key,&compoffset);
233:       DMNetworkGetVertexOffset(networkdm,v,&lofst);
234:       node = (Node*)(arr + compoffset);

236:       if (node->gr) {
237:         row[0] = lofst;
238:         col[0] = lofst;   val[0] =  1;
239:         MatSetValuesLocal(v22,1,row,1,col,val,INSERT_VALUES);
240:       } else {
241:         /* TODO: this is not a nested vector. Need to implement nested vector */
242:         DMNetworkGetVariableOffset(networkdm,v,&lofst);
243:         barr[lofst] -= node->inj;
244:       }
245:     }
246:   }

248:   VecRestoreArray(localb,&barr);

250:   DMLocalToGlobalBegin(networkdm,localb,ADD_VALUES,b);
251:   DMLocalToGlobalEnd(networkdm,localb,ADD_VALUES,b);
252:   DMRestoreLocalVector(networkdm,&localb);

254:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
255:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
256:   return(0);
257: }

259: int main(int argc,char ** argv)
260: {
261:   PetscErrorCode    ierr;
262:   PetscInt          i, nnode = 0, nbranch = 0;
263:   PetscInt          eStart, eEnd, vStart, vEnd;
264:   PetscMPIInt       size, rank;
265:   DM                networkdm;
266:   Vec               x, b;
267:   Mat               A;
268:   KSP               ksp;
269:   int               *edgelist = NULL;
270:   PetscInt          componentkey[2];
271:   Node              *node;
272:   Branch            *branch;

274:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
275:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
276:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

278:   /* "read" data only for processor 0 */
279:   if (!rank) {
280:     read_data(&nnode, &nbranch, &node, &branch, &edgelist);
281:   }

283:   DMNetworkCreate(PETSC_COMM_WORLD,&networkdm);

285:   DMNetworkRegisterComponent(networkdm,"nstr",sizeof(Node),&componentkey[0]);
286:   DMNetworkRegisterComponent(networkdm,"bsrt",sizeof(Branch),&componentkey[1]);


289:   /* Set number of nodes/edges */
290:   DMNetworkSetSizes(networkdm,nnode,nbranch,PETSC_DETERMINE,PETSC_DETERMINE);
291:   /* Add edge connectivity */
292:   DMNetworkSetEdgeList(networkdm,edgelist);
293:   /* Set up the network layout */
294:   DMNetworkLayoutSetUp(networkdm);

296:   /* Add network components: physical parameters of nodes and branches*/
297:   if (!rank) {
298:     DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);
299:     for (i = eStart; i < eEnd; i++) {
300:       DMNetworkAddComponent(networkdm,i,componentkey[1],&branch[i-eStart]);
301:       DMNetworkAddNumVariables(networkdm,i,1);
302:     }

304:     DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);
305:     for (i = vStart; i < vEnd; i++) {
306:       DMNetworkAddComponent(networkdm,i,componentkey[0],&node[i-vStart]);
307:       /* Add number of variables */
308:       DMNetworkAddNumVariables(networkdm,i,1);
309:     }
310:   }

312:   /* Network partitioning and distribution of data */
313:   DMSetUp(networkdm);
314:   DMNetworkDistribute(&networkdm,0);

316:   DMNetworkAssembleGraphStructures(networkdm);

318:   /* Print some info */
319: #if 0
320:   PetscInt offset, goffset;
321:   DMNetworkGetEdgeRange(networkdm,&eStart,&eEnd);
322:   DMNetworkGetVertexRange(networkdm,&vStart,&vEnd);

324:   for (i = eStart; i < eEnd; i++) {
325:     DMNetworkGetVariableOffset(networkdm,i,&offset);
326:     DMNetworkGetVariableGlobalOffset(networkdm,i,&goffset);
327:     PetscPrintf(PETSC_COMM_SELF,"rank[%d] edge %d - loff: %d, goff: %d .\n",rank,i,offset,goffset);
328:   }
329:   for (i = vStart; i < vEnd; i++) {
330:     DMNetworkGetVariableOffset(networkdm,i,&offset);
331:     DMNetworkGetVariableGlobalOffset(networkdm,i,&goffset);
332:     PetscPrintf("rank[%d] vertex %d - loff: %d, goff: %d .\n",rank,i,offset,goffset);
333:   }
334: #endif

336:   /* We don't use these data structures anymore since they have been copied to networkdm */
337:   if (!rank) {
338:     PetscFree(edgelist);
339:     PetscFree(node);
340:     PetscFree(branch);
341:   }

343:   DMCreateGlobalVector(networkdm,&x);
344:   VecDuplicate(x,&b);

346:   DMSetMatType(networkdm, MATNEST);
347:   DMCreateMatrix(networkdm,&A);

349:   /* Assembly system of equations */
350:   FormOperator(networkdm,A,b);

352:   KSPCreate(PETSC_COMM_WORLD, &ksp);
353:   KSPSetOperators(ksp, A, A);
354:   KSPSetFromOptions(ksp);
355:   KSPSolve(ksp, b, x);
356:   VecView(x, 0);

358:   VecDestroy(&x);
359:   VecDestroy(&b);
360:   MatDestroy(&A);
361:   KSPDestroy(&ksp);
362:   DMDestroy(&networkdm);
363:   PetscFinalize();
364:   return ierr;
365: }