Actual source code: ex1.c
petsc-3.10.5 2019-03-28
1: static char help[] = "This example demonstrates the use of DMNetwork interface with subnetworks for solving a coupled nonlinear \n\
2: electric power grid and water pipe problem.\n\
3: The available solver options are in the ex1options file \n\
4: and the data files are in the datafiles of subdirectories.\n\
5: This example shows the use of subnetwork feature in DMNetwork. \n\
6: Run this program: mpiexec -n <n> ./ex1 \n\\n";
8: /* T
9: Concepts: DMNetwork
10: Concepts: PETSc SNES solver
11: */
13: #include "power/power.h"
14: #include "water/water.h"
16: typedef struct{
17: UserCtx_Power appctx_power;
18: AppCtx_Water appctx_water;
19: PetscInt subsnes_id; /* snes solver id */
20: PetscInt it; /* iteration number */
21: Vec localXold; /* store previous solution, used by FormFunction_Dummy() */
22: } UserCtx;
24: PetscErrorCode UserMonitor(SNES snes,PetscInt its,PetscReal fnorm ,void *appctx)
25: {
27: UserCtx *user = (UserCtx*)appctx;
28: Vec X,localXold=user->localXold;
29: DM networkdm;
30: PetscMPIInt rank;
31: MPI_Comm comm;
34: PetscObjectGetComm((PetscObject)snes,&comm);
35: MPI_Comm_rank(comm,&rank);
36: #if 0
37: if (!rank) {
38: PetscInt subsnes_id=user->subsnes_id;
39: if (subsnes_id == 2) {
40: PetscPrintf(PETSC_COMM_SELF," it %d, subsnes_id %d, fnorm %g\n",user->it,user->subsnes_id,fnorm);
41: } else {
42: PetscPrintf(PETSC_COMM_SELF," subsnes_id %d, fnorm %g\n",user->subsnes_id,fnorm);
43: }
44: }
45: #endif
46: SNESGetSolution(snes,&X);
47: SNESGetDM(snes,&networkdm);
48: DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localXold);
49: DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localXold);
50: return(0);
51: }
53: PetscErrorCode FormJacobian_subPower(SNES snes,Vec X, Mat J,Mat Jpre,void *appctx)
54: {
56: DM networkdm;
57: Vec localX;
58: PetscInt nv,ne,i,j,offset,nvar,row;
59: const PetscInt *vtx,*edges;
60: PetscBool ghostvtex;
61: PetscScalar one = 1.0;
62: PetscMPIInt rank;
63: MPI_Comm comm;
66: SNESGetDM(snes,&networkdm);
67: DMGetLocalVector(networkdm,&localX);
69: PetscObjectGetComm((PetscObject)networkdm,&comm);
70: MPI_Comm_rank(comm,&rank);
72: DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);
73: DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);
75: MatZeroEntries(J);
77: /* Power subnetwork: copied from power/FormJacobian_Power() */
78: DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);
79: FormJacobian_Power_private(networkdm,localX,J,nv,ne,vtx,edges,appctx);
81: /* Water subnetwork: Identity */
82: DMNetworkGetSubnetworkInfo(networkdm,1,&nv,&ne,&vtx,&edges);
83: for (i=0; i<nv; i++) {
84: DMNetworkIsGhostVertex(networkdm,vtx[i],&ghostvtex);
85: if (ghostvtex) continue;
87: DMNetworkGetVariableGlobalOffset(networkdm,vtx[i],&offset);
88: DMNetworkGetNumVariables(networkdm,vtx[i],&nvar);
89: for (j=0; j<nvar; j++) {
90: row = offset + j;
91: MatSetValues(J,1,&row,1,&row,&one,ADD_VALUES);
92: }
93: }
94: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
95: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
97: DMRestoreLocalVector(networkdm,&localX);
98: return(0);
99: }
101: /* Dummy equation localF(X) = localX - localXold */
102: PetscErrorCode FormFunction_Dummy(DM networkdm,Vec localX, Vec localF,PetscInt nv,PetscInt ne,const PetscInt* vtx,const PetscInt* edges,void* appctx)
103: {
104: PetscErrorCode ierr;
105: const PetscScalar *xarr,*xoldarr;
106: PetscScalar *farr;
107: PetscInt i,j,offset,nvar;
108: PetscBool ghostvtex;
109: UserCtx *user = (UserCtx*)appctx;
110: Vec localXold = user->localXold;
113: VecGetArrayRead(localX,&xarr);
114: VecGetArrayRead(localXold,&xoldarr);
115: VecGetArray(localF,&farr);
117: for (i=0; i<nv; i++) {
118: DMNetworkIsGhostVertex(networkdm,vtx[i],&ghostvtex);
119: if(ghostvtex) continue;
121: DMNetworkGetVariableOffset(networkdm,vtx[i],&offset);
122: DMNetworkGetNumVariables(networkdm,vtx[i],&nvar);
123: for (j=0; j<nvar; j++) {
124: farr[offset+j] = xarr[offset+j] - xoldarr[offset+j];
125: }
126: }
128: VecRestoreArrayRead(localX,&xarr);
129: VecRestoreArrayRead(localXold,&xoldarr);
130: VecRestoreArray(localF,&farr);
131: return(0);
132: }
134: PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *appctx)
135: {
137: DM networkdm;
138: Vec localX,localF;
139: PetscInt nv,ne;
140: const PetscInt *vtx,*edges;
141: PetscMPIInt rank;
142: MPI_Comm comm;
143: UserCtx *user = (UserCtx*)appctx;
144: UserCtx_Power appctx_power = (*user).appctx_power;
145: AppCtx_Water appctx_water = (*user).appctx_water;
148: SNESGetDM(snes,&networkdm);
149: PetscObjectGetComm((PetscObject)networkdm,&comm);
150: MPI_Comm_rank(comm,&rank);
152: DMGetLocalVector(networkdm,&localX);
153: DMGetLocalVector(networkdm,&localF);
154: VecSet(F,0.0);
155: VecSet(localF,0.0);
157: DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);
158: DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);
160: /* Form Function for power subnetwork */
161: DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);
162: if (user->subsnes_id == 1) { /* snes_water only */
163: FormFunction_Dummy(networkdm,localX,localF,nv,ne,vtx,edges,user);
164: } else {
165: FormFunction_Power(networkdm,localX,localF,nv,ne,vtx,edges,&appctx_power);
166: }
168: /* Form Function for water subnetwork */
169: DMNetworkGetSubnetworkInfo(networkdm,1,&nv,&ne,&vtx,&edges);
170: if (user->subsnes_id == 0) { /* snes_power only */
171: FormFunction_Dummy(networkdm,localX,localF,nv,ne,vtx,edges,user);
172: } else {
173: FormFunction_Water(networkdm,localX,localF,nv,ne,vtx,edges,NULL);
174: }
176: /* Form Function for the coupling subnetwork -- experimental, not done yet */
177: DMNetworkGetSubnetworkInfo(networkdm,2,&nv,&ne,&vtx,&edges);
178: if (user->subsnes_id != 0 && ne) {
179: const PetscInt *cone,*connedges,*econe;
180: PetscInt key,vid[2],nconnedges,k,e,keye;
181: void* component;
183: DMNetworkGetConnectedVertices(networkdm,edges[0],&cone);
185: DMNetworkGetGlobalVertexIndex(networkdm,cone[0],&vid[0]);
186: DMNetworkGetGlobalVertexIndex(networkdm,cone[1],&vid[1]);
187: /* PetscPrintf(PETSC_COMM_SELF,"[%d] Formfunction, coupling subnetwork: nv %d, ne %d; connected vertices %d %d\n",rank,nv,ne,vid[0],vid[1]); */
189: /* Get coupling powernet load vertex */
190: DMNetworkGetComponent(networkdm,cone[0],1,&key,&component);
191: if (key != appctx_power.compkey_load) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not a power load vertex");
193: /* Get coupling water vertex and pump edge */
194: DMNetworkGetComponent(networkdm,cone[1],0,&key,&component);
195: if (key != appctx_water.compkey_vtx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not a water vertex");
197: /* get its supporting edges */
198: DMNetworkGetSupportingEdges(networkdm,cone[1],&nconnedges,&connedges);
199: /* printf("nconnedges %d\n",nconnedges); */
201: for (k=0; k<nconnedges; k++) {
202: e = connedges[k];
203: DMNetworkGetComponent(networkdm,e,0,&keye,&component);
205: if (keye == appctx_water.compkey_edge) { /* water_compkey_edge */
206: EDGE_Water edge=(EDGE_Water)component;
207: if (edge->type == EDGE_TYPE_PUMP) {
208: /* compute flow_pump */
209: PetscInt offsetnode1,offsetnode2,key_0,key_1;
210: const PetscScalar *xarr;
211: PetscScalar *farr;
212: VERTEX_Water vertexnode1,vertexnode2;
214: DMNetworkGetConnectedVertices(networkdm,e,&econe);
215: DMNetworkGetGlobalVertexIndex(networkdm,econe[0],&vid[0]);
216: DMNetworkGetGlobalVertexIndex(networkdm,econe[1],&vid[1]);
218: VecGetArray(localF,&farr);
219: DMNetworkGetVariableOffset(networkdm,econe[0],&offsetnode1);
220: DMNetworkGetVariableOffset(networkdm,econe[1],&offsetnode2);
222: VecGetArrayRead(localX,&xarr);
223: #if 0
224: PetscScalar hf,ht;
225: Pump *pump;
226: pump = &edge->pump;
227: hf = xarr[offsetnode1];
228: ht = xarr[offsetnode2];
230: PetscScalar flow = Flow_Pump(pump,hf,ht);
231: PetscScalar Hp = 0.1; /* load->pl */
232: PetscScalar flow_couple = 8.81*Hp*1.e6/(ht-hf); /* pump->h0; */
233: /* printf("pump %d: connected vtx %d %d; flow_pump %g flow_couple %g; offset %d %d\n",e,vid[0],vid[1],flow,flow_couple,offsetnode1,offsetnode2); */
234: #endif
235: /* Get the components at the two vertices */
236: DMNetworkGetComponent(networkdm,econe[0],0,&key_0,(void**)&vertexnode1);
237: if (key_0 != appctx_water.compkey_vtx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not a water vertex");
238: DMNetworkGetComponent(networkdm,econe[1],0,&key_1,(void**)&vertexnode2);
239: if (key_1 != appctx_water.compkey_vtx) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Not a water vertex");
240: #if 0
241: /* subtract flow_pump computed in FormFunction_Water() and add flow_couple to connected nodes */
242: if (vertexnode1->type == VERTEX_TYPE_JUNCTION) {
243: farr[offsetnode1] += flow;
244: farr[offsetnode1] -= flow_couple;
245: }
246: if (vertexnode2->type == VERTEX_TYPE_JUNCTION) {
247: farr[offsetnode2] -= flow;
248: farr[offsetnode2] += flow_couple;
249: }
250: #endif
251: VecRestoreArrayRead(localX,&xarr);
252: VecRestoreArray(localF,&farr);
253: }
254: }
255: }
256: }
258: DMRestoreLocalVector(networkdm,&localX);
260: DMLocalToGlobalBegin(networkdm,localF,ADD_VALUES,F);
261: DMLocalToGlobalEnd(networkdm,localF,ADD_VALUES,F);
262: DMRestoreLocalVector(networkdm,&localF);
263: /* VecView(F,PETSC_VIEWER_STDOUT_WORLD); */
264: return(0);
265: }
267: PetscErrorCode SetInitialGuess(DM networkdm,Vec X,void* appctx)
268: {
270: PetscInt nv,ne;
271: const PetscInt *vtx,*edges;
272: UserCtx *user = (UserCtx*)appctx;
273: Vec localX = user->localXold;
274: UserCtx_Power appctx_power = (*user).appctx_power;
277: VecSet(X,0.0);
278: VecSet(localX,0.0);
280: /* Set initial guess for power subnetwork */
281: DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);
282: SetInitialGuess_Power(networkdm,localX,nv,ne,vtx,edges,&appctx_power);
284: /* Set initial guess for water subnetwork */
285: DMNetworkGetSubnetworkInfo(networkdm,1,&nv,&ne,&vtx,&edges);
286: SetInitialGuess_Water(networkdm,localX,nv,ne,vtx,edges,NULL);
288: #if 0
289: /* Set initial guess for the coupling subnet */
290: DMNetworkGetSubnetworkInfo(networkdm,2,&nv,&ne,&vtx,&edges);
291: if (ne) {
292: const PetscInt *cone;
293: DMNetworkGetConnectedVertices(networkdm,edges[0],&cone);
294: }
295: #endif
297: DMLocalToGlobalBegin(networkdm,localX,ADD_VALUES,X);
298: DMLocalToGlobalEnd(networkdm,localX,ADD_VALUES,X);
299: return(0);
300: }
302: int main(int argc,char **argv)
303: {
304: PetscErrorCode ierr;
305: DM networkdm;
306: PetscLogStage stage[4];
307: PetscMPIInt rank;
308: PetscInt nsubnet = 3,numVertices[3],NumVertices[3],numEdges[3],NumEdges[3];
309: PetscInt i,j,nv,ne;
310: PetscInt *edgelist[3];
311: const PetscInt *vtx,*edges;
312: Vec X,F;
313: SNES snes,snes_power,snes_water;
314: Mat Jac;
315: PetscBool viewJ=PETSC_FALSE,viewX=PETSC_FALSE;
316: UserCtx user;
317: PetscInt it_max=10;
318: SNESConvergedReason reason;
320: /* Power subnetwork */
321: UserCtx_Power *appctx_power = &user.appctx_power;
322: char pfdata_file[PETSC_MAX_PATH_LEN]="power/case9.m";
323: PFDATA *pfdata=NULL;
324: PetscInt genj,loadj;
325: PetscInt *edgelist_power=NULL;
326: PetscScalar Sbase;
328: /* Water subnetwork */
329: AppCtx_Water *appctx_water = &user.appctx_water;
330: WATERDATA *waterdata=NULL;
331: char waterdata_file[PETSC_MAX_PATH_LEN]="water/sample1.inp";
332: PetscInt *edgelist_water=NULL;
334: /* Coupling subnetwork */
335: PetscInt *edgelist_couple=NULL;
337: PetscInitialize(&argc,&argv,"ex1options",help);
338: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
340: /* (1) Read Data - Only rank 0 reads the data */
341: /*--------------------------------------------*/
342: PetscLogStageRegister("Read Data",&stage[0]);
343: PetscLogStagePush(stage[0]);
345: for (i=0; i<nsubnet; i++) {
346: numVertices[i] = 0; NumVertices[i] = PETSC_DETERMINE;
347: numEdges[i] = 0; NumEdges[i] = PETSC_DETERMINE;
348: }
350: /* READ THE DATA FOR THE FIRST SUBNETWORK: Electric Power Grid */
351: if (!rank) {
352: PetscOptionsGetString(NULL,NULL,"-pfdata",pfdata_file,PETSC_MAX_PATH_LEN-1,NULL);
353: PetscNew(&pfdata);
354: PFReadMatPowerData(pfdata,pfdata_file);
355: Sbase = pfdata->sbase;
357: numEdges[0] = pfdata->nbranch;
358: numVertices[0] = pfdata->nbus;
360: PetscMalloc1(2*numEdges[0],&edgelist_power);
361: GetListofEdges_Power(pfdata,edgelist_power);
362: #if 0
363: printf("edgelist_power:\n");
364: for (i=0; i<numEdges[0]; i++) {
365: PetscPrintf(PETSC_COMM_SELF,"[%D %D]",edgelist_power[2*i],edgelist_power[2*i+1]);
366: }
367: printf("\n");
368: #endif
369: }
370: /* Broadcast power Sbase to all processors */
371: MPI_Bcast(&Sbase,1,MPIU_SCALAR,0,PETSC_COMM_WORLD);
372: appctx_power->Sbase = Sbase;
373: appctx_power->jac_error = PETSC_FALSE;
374: /* If external option activated. Introduce error in jacobian */
375: PetscOptionsHasName(NULL,NULL, "-jac_error", &appctx_power->jac_error);
377: /* GET DATA FOR THE SECOND SUBNETWORK: Water */
378: PetscNew(&waterdata);
379: if (!rank) {
380: PetscOptionsGetString(NULL,NULL,"-waterdata",waterdata_file,PETSC_MAX_PATH_LEN-1,NULL);
381: WaterReadData(waterdata,waterdata_file);
383: PetscCalloc1(2*waterdata->nedge,&edgelist_water);
384: GetListofEdges_Water(waterdata,edgelist_water);
385: numEdges[1] = waterdata->nedge;
386: numVertices[1] = waterdata->nvertex;
387: #if 0
388: printf("edgelist_water:\n");
389: for (i=0; i<numEdges[1]; i++) {
390: PetscPrintf(PETSC_COMM_SELF,"[%D %D]",edgelist_water[2*i],edgelist_water[2*i+1]);
391: }
392: printf("\n");
393: #endif
394: }
396: /* Get data for the coupling subnetwork */
397: if (!rank) {
398: numEdges[2] = 1; numVertices[2] = 0;
399: PetscMalloc1(4*numEdges[2],&edgelist_couple);
400: edgelist_couple[0] = 0; edgelist_couple[1] = 4; /* from node: net[0] vertex[4] */
401: edgelist_couple[2] = 1; edgelist_couple[3] = 0; /* to node: net[1] vertex[0] */
402: }
403: PetscLogStagePop();
405: /* (2) Create network */
406: /*--------------------*/
407: MPI_Barrier(PETSC_COMM_WORLD);
408: PetscLogStageRegister("Net Setup",&stage[1]);
409: PetscLogStagePush(stage[1]);
411: /* Create an empty network object */
412: DMNetworkCreate(PETSC_COMM_WORLD,&networkdm);
414: /* Register the components in the network */
415: DMNetworkRegisterComponent(networkdm,"branchstruct",sizeof(struct _p_EDGE_Power),&appctx_power->compkey_branch);
416: DMNetworkRegisterComponent(networkdm,"busstruct",sizeof(struct _p_VERTEX_Power),&appctx_power->compkey_bus);
417: DMNetworkRegisterComponent(networkdm,"genstruct",sizeof(struct _p_GEN),&appctx_power->compkey_gen);
418: DMNetworkRegisterComponent(networkdm,"loadstruct",sizeof(struct _p_LOAD),&appctx_power->compkey_load);
420: DMNetworkRegisterComponent(networkdm,"edge_water",sizeof(struct _p_EDGE_Water),&appctx_water->compkey_edge);
421: DMNetworkRegisterComponent(networkdm,"vertex_water",sizeof(struct _p_VERTEX_Water),&appctx_water->compkey_vtx);
423: if (!rank) {
424: PetscPrintf(PETSC_COMM_SELF,"[%d] Total local nvertices %D + %D = %D, nedges %D + %D = %D\n",rank,numVertices[0],numVertices[1],numVertices[0]+numVertices[1],numEdges[0],numEdges[1],numEdges[0]+numEdges[1]);
425: }
426: DMNetworkSetSizes(networkdm,nsubnet-1,1,numVertices,numEdges,NumVertices,NumEdges);
428: /* Add edge connectivity */
429: edgelist[0] = edgelist_power;
430: edgelist[1] = edgelist_water;
431: edgelist[2] = edgelist_couple;
432: DMNetworkSetEdgeList(networkdm,edgelist,edgelist+2);
434: /* Set up the network layout */
435: DMNetworkLayoutSetUp(networkdm);
437: /* Add network components - only process[0] has any data to add */
438: /* ADD VARIABLES AND COMPONENTS FOR THE POWER SUBNETWORK */
439: DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);
440: PetscPrintf(PETSC_COMM_WORLD,"Power network: nv %D, ne %D\n",nv,ne);
441: genj = 0; loadj = 0;
442: if (!rank) {
443: for (i = 0; i < ne; i++) {
444: DMNetworkAddComponent(networkdm,edges[i],appctx_power->compkey_branch,&pfdata->branch[i]);
445: }
447: for (i = 0; i < nv; i++) {
448: DMNetworkAddComponent(networkdm,vtx[i],appctx_power->compkey_bus,&pfdata->bus[i]);
449: if (pfdata->bus[i].ngen) {
450: for (j = 0; j < pfdata->bus[i].ngen; j++) {
451: DMNetworkAddComponent(networkdm,vtx[i],appctx_power->compkey_gen,&pfdata->gen[genj++]);
452: }
453: }
454: if (pfdata->bus[i].nload) {
455: for (j=0; j < pfdata->bus[i].nload; j++) {
456: DMNetworkAddComponent(networkdm,vtx[i],appctx_power->compkey_load,&pfdata->load[loadj++]);
457: }
458: }
459: /* Add number of variables */
460: DMNetworkAddNumVariables(networkdm,vtx[i],2);
461: }
462: }
464: /* ADD VARIABLES AND COMPONENTS FOR THE WATER SUBNETWORK */
465: DMNetworkGetSubnetworkInfo(networkdm,1,&nv,&ne,&vtx,&edges);
466: PetscPrintf(PETSC_COMM_WORLD,"Water network: nv %D, ne %D\n",nv,ne);
467: if (!rank) {
468: for (i = 0; i < ne; i++) {
469: DMNetworkAddComponent(networkdm,edges[i],appctx_water->compkey_edge,&waterdata->edge[i]);
470: }
472: for (i = 0; i < nv; i++) {
473: DMNetworkAddComponent(networkdm,vtx[i],appctx_water->compkey_vtx,&waterdata->vertex[i]);
474: /* Add number of variables */
475: DMNetworkAddNumVariables(networkdm,vtx[i],1);
476: }
477: }
479: /* Set up DM for use */
480: DMSetUp(networkdm);
482: /* Distribute networkdm to multiple processes */
483: DMNetworkDistribute(&networkdm,0);
485: DMCreateGlobalVector(networkdm,&X);
486: VecDuplicate(X,&F);
487: DMGetLocalVector(networkdm,&user.localXold);
489: PetscLogStagePop();
491: /* (3) Setup Solvers */
492: /*-------------------*/
493: PetscOptionsGetBool(NULL,NULL,"-viewJ",&viewJ,NULL);
494: PetscOptionsGetBool(NULL,NULL,"-viewX",&viewX,NULL);
496: PetscLogStageRegister("SNES Setup",&stage[2]);
497: PetscLogStagePush(stage[2]);
499: SetInitialGuess(networkdm,X,&user);
500: /* VecView(X,PETSC_VIEWER_STDOUT_WORLD); */
502: /* Create coupled snes */
503: /*-------------------- */
504: PetscPrintf(PETSC_COMM_WORLD,"SNES_coupled setup ......\n");
505: user.subsnes_id = nsubnet-1;
506: SNESCreate(PETSC_COMM_WORLD,&snes);
507: SNESSetDM(snes,networkdm);
508: SNESSetOptionsPrefix(snes,"coupled_");
509: SNESSetFunction(snes,F,FormFunction,&user);
510: SNESMonitorSet(snes,UserMonitor,&user,NULL);
511: SNESSetFromOptions(snes);
513: if (viewJ) {
514: /* View Jac structure */
515: SNESGetJacobian(snes,&Jac,NULL,NULL,NULL);
516: MatView(Jac,PETSC_VIEWER_DRAW_WORLD);
517: }
519: if (viewX) {
520: PetscPrintf(PETSC_COMM_WORLD,"Solution:\n");
521: VecView(X,PETSC_VIEWER_STDOUT_WORLD);
522: }
524: if (viewJ) {
525: /* View assembled Jac */
526: MatView(Jac,PETSC_VIEWER_DRAW_WORLD);
527: }
529: /* Create snes_power */
530: /*-------------------*/
531: PetscPrintf(PETSC_COMM_WORLD,"SNES_power setup ......\n");
532: SetInitialGuess(networkdm,X,&user);
534: user.subsnes_id = 0;
535: SNESCreate(PETSC_COMM_WORLD,&snes_power);
536: SNESSetDM(snes_power,networkdm);
537: SNESSetOptionsPrefix(snes_power,"power_");
538: SNESSetFunction(snes_power,F,FormFunction,&user);
539: SNESMonitorSet(snes_power,UserMonitor,&user,NULL);
541: /* Use user-provide Jacobian */
542: DMCreateMatrix(networkdm,&Jac);
543: SNESSetJacobian(snes_power,Jac,Jac,FormJacobian_subPower,&user);
544: SNESSetFromOptions(snes_power);
546: if (viewX) {
547: PetscPrintf(PETSC_COMM_WORLD,"Power Solution:\n");
548: VecView(X,PETSC_VIEWER_STDOUT_WORLD);
549: }
550: if (viewJ) {
551: PetscPrintf(PETSC_COMM_WORLD,"Power Jac:\n");
552: SNESGetJacobian(snes_power,&Jac,NULL,NULL,NULL);
553: MatView(Jac,PETSC_VIEWER_DRAW_WORLD);
554: /* MatView(Jac,PETSC_VIEWER_STDOUT_WORLD); */
555: }
557: /* Create snes_water */
558: /*-------------------*/
559: PetscPrintf(PETSC_COMM_WORLD,"SNES_water setup......\n");
560: SetInitialGuess(networkdm,X,&user);
562: user.subsnes_id = 1;
563: SNESCreate(PETSC_COMM_WORLD,&snes_water);
564: SNESSetDM(snes_water,networkdm);
565: SNESSetOptionsPrefix(snes_water,"water_");
566: SNESSetFunction(snes_water,F,FormFunction,&user);
567: SNESMonitorSet(snes_water,UserMonitor,&user,NULL);
568: SNESSetFromOptions(snes_water);
570: if (viewX) {
571: PetscPrintf(PETSC_COMM_WORLD,"Water Solution:\n");
572: VecView(X,PETSC_VIEWER_STDOUT_WORLD);
573: }
574: PetscLogStagePop();
576: /* (4) Solve */
577: /*-----------*/
578: PetscLogStageRegister("SNES Solve",&stage[3]);
579: PetscLogStagePush(stage[3]);
580: user.it = 0;
581: reason = SNES_DIVERGED_DTOL;
582: while (user.it < it_max && (PetscInt)reason<0) {
583: #if 0
584: user.subsnes_id = 0;
585: SNESSolve(snes_power,NULL,X);
587: user.subsnes_id = 1;
588: SNESSolve(snes_water,NULL,X);
589: #endif
590: user.subsnes_id = nsubnet-1;
591: SNESSolve(snes,NULL,X);
593: SNESGetConvergedReason(snes,&reason);
594: user.it++;
595: }
596: PetscPrintf(PETSC_COMM_WORLD,"Coupled_SNES converged in %D iterations\n",user.it);
597: if (viewX) {
598: PetscPrintf(PETSC_COMM_WORLD,"Final Solution:\n");
599: VecView(X,PETSC_VIEWER_STDOUT_WORLD);
600: }
601: PetscLogStagePop();
603: /* Free objects */
604: /* -------------*/
605: VecDestroy(&X);
606: VecDestroy(&F);
607: DMRestoreLocalVector(networkdm,&user.localXold);
609: SNESDestroy(&snes);
610: MatDestroy(&Jac);
611: SNESDestroy(&snes_power);
612: SNESDestroy(&snes_water);
614: if (!rank) {
615: PetscFree(edgelist_power);
616: PetscFree(pfdata->bus);
617: PetscFree(pfdata->gen);
618: PetscFree(pfdata->branch);
619: PetscFree(pfdata->load);
620: PetscFree(pfdata);
622: PetscFree(edgelist_water);
623: PetscFree(waterdata->vertex);
624: PetscFree(waterdata->edge);
626: PetscFree(edgelist_couple);
627: }
628: PetscFree(waterdata);
629: DMDestroy(&networkdm);
631: PetscFinalize();
632: return ierr;
633: }
635: /*TEST
637: build:
638: requires: !complex double define(PETSC_HAVE_ATTRIBUTEALIGNED)
639: depends: power/PFReadData.c power/pffunctions.c water/waterreaddata.c water/waterfunctions.c
641: test:
642: args: -coupled_snes_converged_reason -options_left no
643: localrunfiles: ex1options power/case9.m water/sample1.inp
644: output_file: output/ex1.out
645: requires: double !complex define(PETSC_HAVE_ATTRIBUTEALIGNED)
647: test:
648: suffix: 2
649: nsize: 3
650: args: -coupled_snes_converged_reason -options_left no
651: localrunfiles: ex1options power/case9.m water/sample1.inp
652: output_file: output/ex1.out
653: requires: double !complex define(PETSC_HAVE_ATTRIBUTEALIGNED)
655: TEST*/