1c4762a1bSJed Brown static char help[] = "This example demonstrates the use of DMNetwork interface with subnetworks 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 The data file format used is from the MatPower package (http://www.pserc.cornell.edu//matpower/).\n\
4c4762a1bSJed Brown This example shows the use of subnetwork feature in DMNetwork. It creates duplicates of the same network which are treated as subnetworks.\n\
5c4762a1bSJed Brown Run this program: mpiexec -n <n> ./pf2\n\
6c4762a1bSJed Brown mpiexec -n <n> ./pf2 \n";
7c4762a1bSJed Brown
8c4762a1bSJed Brown #include "power.h"
9c4762a1bSJed Brown #include <petscdmnetwork.h>
10c4762a1bSJed Brown
FormFunction_Subnet(DM networkdm,Vec localX,Vec localF,PetscInt nv,PetscInt ne,const PetscInt * vtx,const PetscInt * edges,void * appctx)11d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction_Subnet(DM networkdm, Vec localX, Vec localF, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx)
12d71ae5a4SJacob Faibussowitsch {
13c4762a1bSJed Brown UserCtx_Power *User = (UserCtx_Power *)appctx;
14c4762a1bSJed Brown PetscInt e, v, vfrom, vto;
15c4762a1bSJed Brown const PetscScalar *xarr;
16c4762a1bSJed Brown PetscScalar *farr;
17c4762a1bSJed Brown PetscInt offsetfrom, offsetto, offset;
18c4762a1bSJed Brown
19c4762a1bSJed Brown PetscFunctionBegin;
209566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localX, &xarr));
219566063dSJacob Faibussowitsch PetscCall(VecGetArray(localF, &farr));
22c4762a1bSJed Brown
23c4762a1bSJed Brown for (v = 0; v < nv; v++) {
24c4762a1bSJed Brown PetscInt i, j, key;
25c4762a1bSJed Brown PetscScalar Vm;
26c4762a1bSJed Brown PetscScalar Sbase = User->Sbase;
27c4762a1bSJed Brown VERTEX_Power bus = NULL;
28c4762a1bSJed Brown GEN gen;
29c4762a1bSJed Brown LOAD load;
30c4762a1bSJed Brown PetscBool ghostvtex;
31c4762a1bSJed Brown PetscInt numComps;
32c4762a1bSJed Brown void *component;
33c4762a1bSJed Brown
349566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[v], &ghostvtex));
359566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm, vtx[v], &numComps));
369566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[v], ALL_COMPONENTS, &offset));
37c4762a1bSJed Brown for (j = 0; j < numComps; j++) {
389566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vtx[v], j, &key, &component, NULL));
39c4762a1bSJed Brown if (key == 1) {
40c4762a1bSJed Brown PetscInt nconnedges;
41c4762a1bSJed Brown const PetscInt *connedges;
42c4762a1bSJed Brown
43*57508eceSPierre Jolivet bus = (VERTEX_Power)component;
44c4762a1bSJed Brown /* Handle reference bus constrained dofs */
45c4762a1bSJed Brown if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) {
46c4762a1bSJed Brown farr[offset] = xarr[offset] - bus->va * PETSC_PI / 180.0;
47c4762a1bSJed Brown farr[offset + 1] = xarr[offset + 1] - bus->vm;
48c4762a1bSJed Brown break;
49c4762a1bSJed Brown }
50c4762a1bSJed Brown
51c4762a1bSJed Brown if (!ghostvtex) {
52c4762a1bSJed Brown Vm = xarr[offset + 1];
53c4762a1bSJed Brown
54c4762a1bSJed Brown /* Shunt injections */
55c4762a1bSJed Brown farr[offset] += Vm * Vm * bus->gl / Sbase;
56c4762a1bSJed Brown if (bus->ide != PV_BUS) farr[offset + 1] += -Vm * Vm * bus->bl / Sbase;
57c4762a1bSJed Brown }
58c4762a1bSJed Brown
599566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSupportingEdges(networkdm, vtx[v], &nconnedges, &connedges));
60c4762a1bSJed Brown for (i = 0; i < nconnedges; i++) {
61c4762a1bSJed Brown EDGE_Power branch;
62c4762a1bSJed Brown PetscInt keye;
63c4762a1bSJed Brown PetscScalar Gff, Bff, Gft, Bft, Gtf, Btf, Gtt, Btt;
64c4762a1bSJed Brown const PetscInt *cone;
65c4762a1bSJed Brown PetscScalar Vmf, Vmt, thetaf, thetat, thetaft, thetatf;
66c4762a1bSJed Brown
67c4762a1bSJed Brown e = connedges[i];
689566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, e, 0, &keye, (void **)&branch, NULL));
69c4762a1bSJed Brown if (!branch->status) continue;
70c4762a1bSJed Brown Gff = branch->yff[0];
71c4762a1bSJed Brown Bff = branch->yff[1];
72c4762a1bSJed Brown Gft = branch->yft[0];
73c4762a1bSJed Brown Bft = branch->yft[1];
74c4762a1bSJed Brown Gtf = branch->ytf[0];
75c4762a1bSJed Brown Btf = branch->ytf[1];
76c4762a1bSJed Brown Gtt = branch->ytt[0];
77c4762a1bSJed Brown Btt = branch->ytt[1];
78c4762a1bSJed Brown
799566063dSJacob Faibussowitsch PetscCall(DMNetworkGetConnectedVertices(networkdm, e, &cone));
80c4762a1bSJed Brown vfrom = cone[0];
81c4762a1bSJed Brown vto = cone[1];
82c4762a1bSJed Brown
839566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vfrom, ALL_COMPONENTS, &offsetfrom));
849566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vto, ALL_COMPONENTS, &offsetto));
85c4762a1bSJed Brown
86c4762a1bSJed Brown thetaf = xarr[offsetfrom];
87c4762a1bSJed Brown Vmf = xarr[offsetfrom + 1];
88c4762a1bSJed Brown thetat = xarr[offsetto];
89c4762a1bSJed Brown Vmt = xarr[offsetto + 1];
90c4762a1bSJed Brown thetaft = thetaf - thetat;
91c4762a1bSJed Brown thetatf = thetat - thetaf;
92c4762a1bSJed Brown
93c4762a1bSJed Brown if (vfrom == vtx[v]) {
94c4762a1bSJed Brown farr[offsetfrom] += Gff * Vmf * Vmf + Vmf * Vmt * (Gft * PetscCosScalar(thetaft) + Bft * PetscSinScalar(thetaft));
95c4762a1bSJed Brown farr[offsetfrom + 1] += -Bff * Vmf * Vmf + Vmf * Vmt * (-Bft * PetscCosScalar(thetaft) + Gft * PetscSinScalar(thetaft));
96c4762a1bSJed Brown } else {
97c4762a1bSJed Brown farr[offsetto] += Gtt * Vmt * Vmt + Vmt * Vmf * (Gtf * PetscCosScalar(thetatf) + Btf * PetscSinScalar(thetatf));
98c4762a1bSJed Brown farr[offsetto + 1] += -Btt * Vmt * Vmt + Vmt * Vmf * (-Btf * PetscCosScalar(thetatf) + Gtf * PetscSinScalar(thetatf));
99c4762a1bSJed Brown }
100c4762a1bSJed Brown }
101c4762a1bSJed Brown } else if (key == 2) {
102c4762a1bSJed Brown if (!ghostvtex) {
103*57508eceSPierre Jolivet gen = (GEN)component;
104c4762a1bSJed Brown if (!gen->status) continue;
105c4762a1bSJed Brown farr[offset] += -gen->pg / Sbase;
106c4762a1bSJed Brown farr[offset + 1] += -gen->qg / Sbase;
107c4762a1bSJed Brown }
108c4762a1bSJed Brown } else if (key == 3) {
109c4762a1bSJed Brown if (!ghostvtex) {
110*57508eceSPierre Jolivet load = (LOAD)component;
111c4762a1bSJed Brown farr[offset] += load->pl / Sbase;
112c4762a1bSJed Brown farr[offset + 1] += load->ql / Sbase;
113c4762a1bSJed Brown }
114c4762a1bSJed Brown }
115c4762a1bSJed Brown }
116ad540459SPierre Jolivet if (bus && bus->ide == PV_BUS) farr[offset + 1] = xarr[offset + 1] - bus->vm;
117c4762a1bSJed Brown }
1189566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(localX, &xarr));
1199566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localF, &farr));
1203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
121c4762a1bSJed Brown }
122c4762a1bSJed Brown
FormFunction(SNES snes,Vec X,Vec F,void * appctx)123d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *appctx)
124d71ae5a4SJacob Faibussowitsch {
125c4762a1bSJed Brown DM networkdm;
126c4762a1bSJed Brown Vec localX, localF;
127c4762a1bSJed Brown PetscInt nv, ne;
128c4762a1bSJed Brown const PetscInt *vtx, *edges;
129c4762a1bSJed Brown
130c4762a1bSJed Brown PetscFunctionBegin;
1319566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &networkdm));
1329566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm, &localX));
1339566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm, &localF));
1349566063dSJacob Faibussowitsch PetscCall(VecSet(F, 0.0));
135c4762a1bSJed Brown
1369566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
1379566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
138c4762a1bSJed Brown
1399566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm, F, INSERT_VALUES, localF));
1409566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm, F, INSERT_VALUES, localF));
141c4762a1bSJed Brown
142c4762a1bSJed Brown /* Form Function for first subnetwork */
1439566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges));
1449566063dSJacob Faibussowitsch PetscCall(FormFunction_Subnet(networkdm, localX, localF, nv, ne, vtx, edges, appctx));
145c4762a1bSJed Brown
146c4762a1bSJed Brown /* Form Function for second subnetwork */
1479566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 1, &nv, &ne, &vtx, &edges));
1489566063dSJacob Faibussowitsch PetscCall(FormFunction_Subnet(networkdm, localX, localF, nv, ne, vtx, edges, appctx));
149c4762a1bSJed Brown
1509566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm, &localX));
151c4762a1bSJed Brown
1529566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(networkdm, localF, ADD_VALUES, F));
1539566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(networkdm, localF, ADD_VALUES, F));
1549566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm, &localF));
1553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
156c4762a1bSJed Brown }
157c4762a1bSJed Brown
FormJacobian_Subnet(DM networkdm,Vec localX,Mat J,Mat Jpre,PetscInt nv,PetscInt ne,const PetscInt * vtx,const PetscInt * edges,void * appctx)158d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian_Subnet(DM networkdm, Vec localX, Mat J, Mat Jpre, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx)
159d71ae5a4SJacob Faibussowitsch {
160c4762a1bSJed Brown UserCtx_Power *User = (UserCtx_Power *)appctx;
161c4762a1bSJed Brown PetscInt e, v, vfrom, vto;
162c4762a1bSJed Brown const PetscScalar *xarr;
163c4762a1bSJed Brown PetscInt offsetfrom, offsetto, goffsetfrom, goffsetto;
164c4762a1bSJed Brown PetscInt row[2], col[8];
165c4762a1bSJed Brown PetscScalar values[8];
166c4762a1bSJed Brown
167c4762a1bSJed Brown PetscFunctionBegin;
1689566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localX, &xarr));
169c4762a1bSJed Brown
170c4762a1bSJed Brown for (v = 0; v < nv; v++) {
171c4762a1bSJed Brown PetscInt i, j, key;
172c4762a1bSJed Brown PetscInt offset, goffset;
173c4762a1bSJed Brown PetscScalar Vm;
174c4762a1bSJed Brown PetscScalar Sbase = User->Sbase;
175c4762a1bSJed Brown VERTEX_Power bus;
176c4762a1bSJed Brown PetscBool ghostvtex;
177c4762a1bSJed Brown PetscInt numComps;
178c4762a1bSJed Brown void *component;
179c4762a1bSJed Brown
1809566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[v], &ghostvtex));
1819566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm, vtx[v], &numComps));
182c4762a1bSJed Brown for (j = 0; j < numComps; j++) {
1839566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[v], ALL_COMPONENTS, &offset));
1849566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVecOffset(networkdm, vtx[v], ALL_COMPONENTS, &goffset));
1859566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vtx[v], j, &key, &component, NULL));
186c4762a1bSJed Brown if (key == 1) {
187c4762a1bSJed Brown PetscInt nconnedges;
188c4762a1bSJed Brown const PetscInt *connedges;
189c4762a1bSJed Brown
190*57508eceSPierre Jolivet bus = (VERTEX_Power)component;
191c4762a1bSJed Brown if (!ghostvtex) {
192c4762a1bSJed Brown /* Handle reference bus constrained dofs */
193c4762a1bSJed Brown if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) {
1949371c9d4SSatish Balay row[0] = goffset;
1959371c9d4SSatish Balay row[1] = goffset + 1;
1969371c9d4SSatish Balay col[0] = goffset;
1979371c9d4SSatish Balay col[1] = goffset + 1;
1989371c9d4SSatish Balay col[2] = goffset;
1999371c9d4SSatish Balay col[3] = goffset + 1;
2009371c9d4SSatish Balay values[0] = 1.0;
2019371c9d4SSatish Balay values[1] = 0.0;
2029371c9d4SSatish Balay values[2] = 0.0;
2039371c9d4SSatish Balay values[3] = 1.0;
2049566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 2, row, 2, col, values, ADD_VALUES));
205c4762a1bSJed Brown break;
206c4762a1bSJed Brown }
207c4762a1bSJed Brown
208c4762a1bSJed Brown Vm = xarr[offset + 1];
209c4762a1bSJed Brown
210c4762a1bSJed Brown /* Shunt injections */
2119371c9d4SSatish Balay row[0] = goffset;
2129371c9d4SSatish Balay row[1] = goffset + 1;
2139371c9d4SSatish Balay col[0] = goffset;
2149371c9d4SSatish Balay col[1] = goffset + 1;
215c4762a1bSJed Brown values[0] = values[1] = values[2] = values[3] = 0.0;
216c4762a1bSJed Brown if (bus->ide != PV_BUS) {
217c4762a1bSJed Brown values[1] = 2.0 * Vm * bus->gl / Sbase;
218c4762a1bSJed Brown values[3] = -2.0 * Vm * bus->bl / Sbase;
219c4762a1bSJed Brown }
2209566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 2, row, 2, col, values, ADD_VALUES));
221c4762a1bSJed Brown }
222c4762a1bSJed Brown
2239566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSupportingEdges(networkdm, vtx[v], &nconnedges, &connedges));
224c4762a1bSJed Brown for (i = 0; i < nconnedges; i++) {
225c4762a1bSJed Brown EDGE_Power branch;
226c4762a1bSJed Brown VERTEX_Power busf, bust;
227c4762a1bSJed Brown PetscInt keyf, keyt;
228c4762a1bSJed Brown PetscScalar Gff, Bff, Gft, Bft, Gtf, Btf, Gtt, Btt;
229c4762a1bSJed Brown const PetscInt *cone;
230c4762a1bSJed Brown PetscScalar Vmf, Vmt, thetaf, thetat, thetaft, thetatf;
231c4762a1bSJed Brown
232c4762a1bSJed Brown e = connedges[i];
2339566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, e, 0, &key, (void **)&branch, NULL));
234c4762a1bSJed Brown if (!branch->status) continue;
235c4762a1bSJed Brown
236c4762a1bSJed Brown Gff = branch->yff[0];
237c4762a1bSJed Brown Bff = branch->yff[1];
238c4762a1bSJed Brown Gft = branch->yft[0];
239c4762a1bSJed Brown Bft = branch->yft[1];
240c4762a1bSJed Brown Gtf = branch->ytf[0];
241c4762a1bSJed Brown Btf = branch->ytf[1];
242c4762a1bSJed Brown Gtt = branch->ytt[0];
243c4762a1bSJed Brown Btt = branch->ytt[1];
244c4762a1bSJed Brown
2459566063dSJacob Faibussowitsch PetscCall(DMNetworkGetConnectedVertices(networkdm, e, &cone));
246c4762a1bSJed Brown vfrom = cone[0];
247c4762a1bSJed Brown vto = cone[1];
248c4762a1bSJed Brown
2499566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vfrom, ALL_COMPONENTS, &offsetfrom));
2509566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vto, ALL_COMPONENTS, &offsetto));
2519566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVecOffset(networkdm, vfrom, ALL_COMPONENTS, &goffsetfrom));
2529566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVecOffset(networkdm, vto, ALL_COMPONENTS, &goffsetto));
253c4762a1bSJed Brown
254c4762a1bSJed Brown if (goffsetto < 0) goffsetto = -goffsetto - 1;
255c4762a1bSJed Brown
256c4762a1bSJed Brown thetaf = xarr[offsetfrom];
257c4762a1bSJed Brown Vmf = xarr[offsetfrom + 1];
258c4762a1bSJed Brown thetat = xarr[offsetto];
259c4762a1bSJed Brown Vmt = xarr[offsetto + 1];
260c4762a1bSJed Brown thetaft = thetaf - thetat;
261c4762a1bSJed Brown thetatf = thetat - thetaf;
262c4762a1bSJed Brown
2639566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vfrom, 0, &keyf, (void **)&busf, NULL));
2649566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vto, 0, &keyt, (void **)&bust, NULL));
265c4762a1bSJed Brown
266c4762a1bSJed Brown if (vfrom == vtx[v]) {
267c4762a1bSJed Brown if (busf->ide != REF_BUS) {
268c4762a1bSJed Brown /* farr[offsetfrom] += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); */
269c4762a1bSJed Brown row[0] = goffsetfrom;
2709371c9d4SSatish Balay col[0] = goffsetfrom;
2719371c9d4SSatish Balay col[1] = goffsetfrom + 1;
2729371c9d4SSatish Balay col[2] = goffsetto;
2739371c9d4SSatish Balay col[3] = goffsetto + 1;
274c4762a1bSJed Brown values[0] = Vmf * Vmt * (Gft * -PetscSinScalar(thetaft) + Bft * PetscCosScalar(thetaft)); /* df_dthetaf */
275c4762a1bSJed Brown values[1] = 2.0 * Gff * Vmf + Vmt * (Gft * PetscCosScalar(thetaft) + Bft * PetscSinScalar(thetaft)); /* df_dVmf */
276c4762a1bSJed Brown values[2] = Vmf * Vmt * (Gft * PetscSinScalar(thetaft) + Bft * -PetscCosScalar(thetaft)); /* df_dthetat */
277c4762a1bSJed Brown values[3] = Vmf * (Gft * PetscCosScalar(thetaft) + Bft * PetscSinScalar(thetaft)); /* df_dVmt */
278c4762a1bSJed Brown
2799566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 4, col, values, ADD_VALUES));
280c4762a1bSJed Brown }
281c4762a1bSJed Brown if (busf->ide != PV_BUS && busf->ide != REF_BUS) {
282c4762a1bSJed Brown row[0] = goffsetfrom + 1;
2839371c9d4SSatish Balay col[0] = goffsetfrom;
2849371c9d4SSatish Balay col[1] = goffsetfrom + 1;
2859371c9d4SSatish Balay col[2] = goffsetto;
2869371c9d4SSatish Balay col[3] = goffsetto + 1;
287c4762a1bSJed Brown /* farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); */
288c4762a1bSJed Brown values[0] = Vmf * Vmt * (Bft * PetscSinScalar(thetaft) + Gft * PetscCosScalar(thetaft));
289c4762a1bSJed Brown values[1] = -2.0 * Bff * Vmf + Vmt * (-Bft * PetscCosScalar(thetaft) + Gft * PetscSinScalar(thetaft));
290c4762a1bSJed Brown values[2] = Vmf * Vmt * (-Bft * PetscSinScalar(thetaft) + Gft * -PetscCosScalar(thetaft));
291c4762a1bSJed Brown values[3] = Vmf * (-Bft * PetscCosScalar(thetaft) + Gft * PetscSinScalar(thetaft));
292c4762a1bSJed Brown
2939566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 4, col, values, ADD_VALUES));
294c4762a1bSJed Brown }
295c4762a1bSJed Brown } else {
296c4762a1bSJed Brown if (bust->ide != REF_BUS) {
297c4762a1bSJed Brown row[0] = goffsetto;
2989371c9d4SSatish Balay col[0] = goffsetto;
2999371c9d4SSatish Balay col[1] = goffsetto + 1;
3009371c9d4SSatish Balay col[2] = goffsetfrom;
3019371c9d4SSatish Balay col[3] = goffsetfrom + 1;
302c4762a1bSJed Brown /* farr[offsetto] += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); */
303c4762a1bSJed Brown values[0] = Vmt * Vmf * (Gtf * -PetscSinScalar(thetatf) + Btf * PetscCosScalar(thetaft)); /* df_dthetat */
304c4762a1bSJed Brown values[1] = 2.0 * Gtt * Vmt + Vmf * (Gtf * PetscCosScalar(thetatf) + Btf * PetscSinScalar(thetatf)); /* df_dVmt */
305c4762a1bSJed Brown values[2] = Vmt * Vmf * (Gtf * PetscSinScalar(thetatf) + Btf * -PetscCosScalar(thetatf)); /* df_dthetaf */
306c4762a1bSJed Brown values[3] = Vmt * (Gtf * PetscCosScalar(thetatf) + Btf * PetscSinScalar(thetatf)); /* df_dVmf */
307c4762a1bSJed Brown
3089566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 4, col, values, ADD_VALUES));
309c4762a1bSJed Brown }
310c4762a1bSJed Brown if (bust->ide != PV_BUS && bust->ide != REF_BUS) {
311c4762a1bSJed Brown row[0] = goffsetto + 1;
3129371c9d4SSatish Balay col[0] = goffsetto;
3139371c9d4SSatish Balay col[1] = goffsetto + 1;
3149371c9d4SSatish Balay col[2] = goffsetfrom;
3159371c9d4SSatish Balay col[3] = goffsetfrom + 1;
316c4762a1bSJed Brown /* farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); */
317c4762a1bSJed Brown values[0] = Vmt * Vmf * (Btf * PetscSinScalar(thetatf) + Gtf * PetscCosScalar(thetatf));
318c4762a1bSJed Brown values[1] = -2.0 * Btt * Vmt + Vmf * (-Btf * PetscCosScalar(thetatf) + Gtf * PetscSinScalar(thetatf));
319c4762a1bSJed Brown values[2] = Vmt * Vmf * (-Btf * PetscSinScalar(thetatf) + Gtf * -PetscCosScalar(thetatf));
320c4762a1bSJed Brown values[3] = Vmt * (-Btf * PetscCosScalar(thetatf) + Gtf * PetscSinScalar(thetatf));
321c4762a1bSJed Brown
3229566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 4, col, values, ADD_VALUES));
323c4762a1bSJed Brown }
324c4762a1bSJed Brown }
325c4762a1bSJed Brown }
326c4762a1bSJed Brown if (!ghostvtex && bus->ide == PV_BUS) {
3279371c9d4SSatish Balay row[0] = goffset + 1;
3289371c9d4SSatish Balay col[0] = goffset + 1;
329c4762a1bSJed Brown values[0] = 1.0;
3309566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 1, col, values, ADD_VALUES));
331c4762a1bSJed Brown }
332c4762a1bSJed Brown }
333c4762a1bSJed Brown }
334c4762a1bSJed Brown }
3359566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(localX, &xarr));
3363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
337c4762a1bSJed Brown }
338c4762a1bSJed Brown
FormJacobian(SNES snes,Vec X,Mat J,Mat Jpre,void * appctx)339d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(SNES snes, Vec X, Mat J, Mat Jpre, void *appctx)
340d71ae5a4SJacob Faibussowitsch {
341c4762a1bSJed Brown DM networkdm;
342c4762a1bSJed Brown Vec localX;
343c4762a1bSJed Brown PetscInt ne, nv;
344c4762a1bSJed Brown const PetscInt *vtx, *edges;
345c4762a1bSJed Brown
346c4762a1bSJed Brown PetscFunctionBegin;
3479566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(J));
348c4762a1bSJed Brown
3499566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &networkdm));
3509566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm, &localX));
351c4762a1bSJed Brown
3529566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
3539566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
354c4762a1bSJed Brown
355c4762a1bSJed Brown /* Form Jacobian for first subnetwork */
3569566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges));
3579566063dSJacob Faibussowitsch PetscCall(FormJacobian_Subnet(networkdm, localX, J, Jpre, nv, ne, vtx, edges, appctx));
358c4762a1bSJed Brown
359c4762a1bSJed Brown /* Form Jacobian for second subnetwork */
3609566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 1, &nv, &ne, &vtx, &edges));
3619566063dSJacob Faibussowitsch PetscCall(FormJacobian_Subnet(networkdm, localX, J, Jpre, nv, ne, vtx, edges, appctx));
362c4762a1bSJed Brown
3639566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm, &localX));
364c4762a1bSJed Brown
3659566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
3669566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
3673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
368c4762a1bSJed Brown }
369c4762a1bSJed Brown
SetInitialValues_Subnet(DM networkdm,Vec localX,PetscInt nv,PetscInt ne,const PetscInt * vtx,const PetscInt * edges,void * appctx)370d71ae5a4SJacob Faibussowitsch PetscErrorCode SetInitialValues_Subnet(DM networkdm, Vec localX, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx)
371d71ae5a4SJacob Faibussowitsch {
372c4762a1bSJed Brown VERTEX_Power bus;
373c4762a1bSJed Brown PetscInt i;
374c4762a1bSJed Brown GEN gen;
375c4762a1bSJed Brown PetscBool ghostvtex;
376c4762a1bSJed Brown PetscScalar *xarr;
377c4762a1bSJed Brown PetscInt key, numComps, j, offset;
378c4762a1bSJed Brown void *component;
379c4762a1bSJed Brown
380c4762a1bSJed Brown PetscFunctionBegin;
3819566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX, &xarr));
382c4762a1bSJed Brown for (i = 0; i < nv; i++) {
3839566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[i], &ghostvtex));
384c4762a1bSJed Brown if (ghostvtex) continue;
385c4762a1bSJed Brown
3869566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[i], ALL_COMPONENTS, &offset));
3879566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm, vtx[i], &numComps));
388c4762a1bSJed Brown for (j = 0; j < numComps; j++) {
3899566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vtx[i], j, &key, &component, NULL));
390c4762a1bSJed Brown if (key == 1) {
391*57508eceSPierre Jolivet bus = (VERTEX_Power)component;
392c4762a1bSJed Brown xarr[offset] = bus->va * PETSC_PI / 180.0;
393c4762a1bSJed Brown xarr[offset + 1] = bus->vm;
394c4762a1bSJed Brown } else if (key == 2) {
395*57508eceSPierre Jolivet gen = (GEN)component;
396c4762a1bSJed Brown if (!gen->status) continue;
397c4762a1bSJed Brown xarr[offset + 1] = gen->vs;
398c4762a1bSJed Brown break;
399c4762a1bSJed Brown }
400c4762a1bSJed Brown }
401c4762a1bSJed Brown }
4029566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX, &xarr));
4033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
404c4762a1bSJed Brown }
405c4762a1bSJed Brown
SetInitialValues(DM networkdm,Vec X,void * appctx)406d71ae5a4SJacob Faibussowitsch PetscErrorCode SetInitialValues(DM networkdm, Vec X, void *appctx)
407d71ae5a4SJacob Faibussowitsch {
408c4762a1bSJed Brown PetscInt nv, ne;
409c4762a1bSJed Brown const PetscInt *vtx, *edges;
410c4762a1bSJed Brown Vec localX;
411c4762a1bSJed Brown
412c4762a1bSJed Brown PetscFunctionBegin;
4139566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm, &localX));
414c4762a1bSJed Brown
4159566063dSJacob Faibussowitsch PetscCall(VecSet(X, 0.0));
4169566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
4179566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
418c4762a1bSJed Brown
419c4762a1bSJed Brown /* Set initial guess for first subnetwork */
4209566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges));
4219566063dSJacob Faibussowitsch PetscCall(SetInitialValues_Subnet(networkdm, localX, nv, ne, vtx, edges, appctx));
422c4762a1bSJed Brown
423c4762a1bSJed Brown /* Set initial guess for second subnetwork */
4249566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 1, &nv, &ne, &vtx, &edges));
4259566063dSJacob Faibussowitsch PetscCall(SetInitialValues_Subnet(networkdm, localX, nv, ne, vtx, edges, appctx));
426c4762a1bSJed Brown
4279566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(networkdm, localX, ADD_VALUES, X));
4289566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(networkdm, localX, ADD_VALUES, X));
4299566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm, &localX));
4303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
431c4762a1bSJed Brown }
432c4762a1bSJed Brown
main(int argc,char ** argv)433d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
434d71ae5a4SJacob Faibussowitsch {
435c4762a1bSJed Brown char pfdata_file[PETSC_MAX_PATH_LEN] = "case9.m";
436c4762a1bSJed Brown PFDATA *pfdata1, *pfdata2;
437f11a936eSBarry Smith PetscInt numEdges1 = 0, numEdges2 = 0;
4382bf73ac6SHong Zhang PetscInt *edgelist1 = NULL, *edgelist2 = NULL, componentkey[4];
439c4762a1bSJed Brown DM networkdm;
440c4762a1bSJed Brown UserCtx_Power User;
441c4762a1bSJed Brown PetscLogStage stage1, stage2;
442c4762a1bSJed Brown PetscMPIInt rank;
4432bf73ac6SHong Zhang PetscInt nsubnet = 2, nv, ne, i, j, genj, loadj;
4442bf73ac6SHong Zhang const PetscInt *vtx, *edges;
445c4762a1bSJed Brown Vec X, F;
446c4762a1bSJed Brown Mat J;
447c4762a1bSJed Brown SNES snes;
448c4762a1bSJed Brown
449327415f7SBarry Smith PetscFunctionBeginUser;
4509566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, "poweroptions", help));
4519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
452c4762a1bSJed Brown {
453c4762a1bSJed 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 */
454c4762a1bSJed Brown /* this is an experiment to see how the analyzer reacts */
455c4762a1bSJed Brown const PetscMPIInt crank = rank;
456c4762a1bSJed Brown
457c4762a1bSJed Brown /* Create an empty network object */
4589566063dSJacob Faibussowitsch PetscCall(DMNetworkCreate(PETSC_COMM_WORLD, &networkdm));
459c4762a1bSJed Brown
460c4762a1bSJed Brown /* Register the components in the network */
4619566063dSJacob Faibussowitsch PetscCall(DMNetworkRegisterComponent(networkdm, "branchstruct", sizeof(struct _p_EDGE_Power), &componentkey[0]));
4629566063dSJacob Faibussowitsch PetscCall(DMNetworkRegisterComponent(networkdm, "busstruct", sizeof(struct _p_VERTEX_Power), &componentkey[1]));
4639566063dSJacob Faibussowitsch PetscCall(DMNetworkRegisterComponent(networkdm, "genstruct", sizeof(struct _p_GEN), &componentkey[2]));
4649566063dSJacob Faibussowitsch PetscCall(DMNetworkRegisterComponent(networkdm, "loadstruct", sizeof(struct _p_LOAD), &componentkey[3]));
465c4762a1bSJed Brown
4669566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Read Data", &stage1));
4673ba16761SJacob Faibussowitsch PetscCall(PetscLogStagePush(stage1));
468c4762a1bSJed Brown /* READ THE DATA */
469c4762a1bSJed Brown if (!crank) {
470c4762a1bSJed Brown /* Only rank 0 reads the data */
4719566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, NULL, "-pfdata", pfdata_file, sizeof(pfdata_file), NULL));
472c4762a1bSJed Brown /* HERE WE CREATE COPIES OF THE SAME NETWORK THAT WILL BE TREATED AS SUBNETWORKS */
473c4762a1bSJed Brown
474c4762a1bSJed Brown /* READ DATA FOR THE FIRST SUBNETWORK */
4759566063dSJacob Faibussowitsch PetscCall(PetscNew(&pfdata1));
4769566063dSJacob Faibussowitsch PetscCall(PFReadMatPowerData(pfdata1, pfdata_file));
477c4762a1bSJed Brown User.Sbase = pfdata1->sbase;
478c4762a1bSJed Brown
479c4762a1bSJed Brown numEdges1 = pfdata1->nbranch;
4809566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * numEdges1, &edgelist1));
4819566063dSJacob Faibussowitsch PetscCall(GetListofEdges_Power(pfdata1, edgelist1));
482c4762a1bSJed Brown
483c4762a1bSJed Brown /* READ DATA FOR THE SECOND SUBNETWORK */
4849566063dSJacob Faibussowitsch PetscCall(PetscNew(&pfdata2));
4859566063dSJacob Faibussowitsch PetscCall(PFReadMatPowerData(pfdata2, pfdata_file));
486c4762a1bSJed Brown User.Sbase = pfdata2->sbase;
487c4762a1bSJed Brown
488c4762a1bSJed Brown numEdges2 = pfdata2->nbranch;
4899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(2 * numEdges2, &edgelist2));
4909566063dSJacob Faibussowitsch PetscCall(GetListofEdges_Power(pfdata2, edgelist2));
491c4762a1bSJed Brown }
492c4762a1bSJed Brown
4933ba16761SJacob Faibussowitsch PetscCall(PetscLogStagePop());
4949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
4959566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("Create network", &stage2));
4963ba16761SJacob Faibussowitsch PetscCall(PetscLogStagePush(stage2));
497c4762a1bSJed Brown
4982bf73ac6SHong Zhang /* Set number of nodes/edges and edge connectivity */
4999566063dSJacob Faibussowitsch PetscCall(DMNetworkSetNumSubNetworks(networkdm, PETSC_DECIDE, nsubnet));
5009566063dSJacob Faibussowitsch PetscCall(DMNetworkAddSubnetwork(networkdm, "", numEdges1, edgelist1, NULL));
5019566063dSJacob Faibussowitsch PetscCall(DMNetworkAddSubnetwork(networkdm, "", numEdges2, edgelist2, NULL));
502c4762a1bSJed Brown
503c4762a1bSJed Brown /* Set up the network layout */
5049566063dSJacob Faibussowitsch PetscCall(DMNetworkLayoutSetUp(networkdm));
505c4762a1bSJed Brown
506c4762a1bSJed Brown /* Add network components only process 0 has any data to add*/
507c4762a1bSJed Brown if (!crank) {
5089371c9d4SSatish Balay genj = 0;
5099371c9d4SSatish Balay loadj = 0;
510c4762a1bSJed Brown
511c4762a1bSJed Brown /* ADD VARIABLES AND COMPONENTS FOR THE FIRST SUBNETWORK */
5129566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges));
513c4762a1bSJed Brown
51448a46eb9SPierre Jolivet for (i = 0; i < ne; i++) PetscCall(DMNetworkAddComponent(networkdm, edges[i], componentkey[0], &pfdata1->branch[i], 0));
515c4762a1bSJed Brown
516c4762a1bSJed Brown for (i = 0; i < nv; i++) {
5179566063dSJacob Faibussowitsch PetscCall(DMNetworkAddComponent(networkdm, vtx[i], componentkey[1], &pfdata1->bus[i], 2));
518c4762a1bSJed Brown if (pfdata1->bus[i].ngen) {
51948a46eb9SPierre Jolivet for (j = 0; j < pfdata1->bus[i].ngen; j++) PetscCall(DMNetworkAddComponent(networkdm, vtx[i], componentkey[2], &pfdata1->gen[genj++], 0));
520c4762a1bSJed Brown }
521c4762a1bSJed Brown if (pfdata1->bus[i].nload) {
52248a46eb9SPierre Jolivet for (j = 0; j < pfdata1->bus[i].nload; j++) PetscCall(DMNetworkAddComponent(networkdm, vtx[i], componentkey[3], &pfdata1->load[loadj++], 0));
523c4762a1bSJed Brown }
524c4762a1bSJed Brown }
525c4762a1bSJed Brown
5269371c9d4SSatish Balay genj = 0;
5279371c9d4SSatish Balay loadj = 0;
528c4762a1bSJed Brown
529c4762a1bSJed Brown /* ADD VARIABLES AND COMPONENTS FOR THE SECOND SUBNETWORK */
5309566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 1, &nv, &ne, &vtx, &edges));
531c4762a1bSJed Brown
53248a46eb9SPierre Jolivet for (i = 0; i < ne; i++) PetscCall(DMNetworkAddComponent(networkdm, edges[i], componentkey[0], &pfdata2->branch[i], 0));
533c4762a1bSJed Brown
534c4762a1bSJed Brown for (i = 0; i < nv; i++) {
5359566063dSJacob Faibussowitsch PetscCall(DMNetworkAddComponent(networkdm, vtx[i], componentkey[1], &pfdata2->bus[i], 2));
536c4762a1bSJed Brown if (pfdata2->bus[i].ngen) {
53748a46eb9SPierre Jolivet for (j = 0; j < pfdata2->bus[i].ngen; j++) PetscCall(DMNetworkAddComponent(networkdm, vtx[i], componentkey[2], &pfdata2->gen[genj++], 0));
538c4762a1bSJed Brown }
539c4762a1bSJed Brown if (pfdata2->bus[i].nload) {
54048a46eb9SPierre Jolivet for (j = 0; j < pfdata2->bus[i].nload; j++) PetscCall(DMNetworkAddComponent(networkdm, vtx[i], componentkey[3], &pfdata2->load[loadj++], 0));
541c4762a1bSJed Brown }
542c4762a1bSJed Brown }
543c4762a1bSJed Brown }
544c4762a1bSJed Brown
545c4762a1bSJed Brown /* Set up DM for use */
5469566063dSJacob Faibussowitsch PetscCall(DMSetUp(networkdm));
547c4762a1bSJed Brown
548c4762a1bSJed Brown if (!crank) {
5499566063dSJacob Faibussowitsch PetscCall(PetscFree(edgelist1));
5509566063dSJacob Faibussowitsch PetscCall(PetscFree(edgelist2));
551c4762a1bSJed Brown }
552c4762a1bSJed Brown
553c4762a1bSJed Brown if (!crank) {
5549566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata1->bus));
5559566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata1->gen));
5569566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata1->branch));
5579566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata1->load));
5589566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata1));
559c4762a1bSJed Brown
5609566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata2->bus));
5619566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata2->gen));
5629566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata2->branch));
5639566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata2->load));
5649566063dSJacob Faibussowitsch PetscCall(PetscFree(pfdata2));
565c4762a1bSJed Brown }
566c4762a1bSJed Brown
567c4762a1bSJed Brown /* Distribute networkdm to multiple processes */
5689566063dSJacob Faibussowitsch PetscCall(DMNetworkDistribute(&networkdm, 0));
569c4762a1bSJed Brown
5703ba16761SJacob Faibussowitsch PetscCall(PetscLogStagePop());
571c4762a1bSJed Brown
572c4762a1bSJed Brown /* Broadcast Sbase to all processors */
5739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Bcast(&User.Sbase, 1, MPIU_SCALAR, 0, PETSC_COMM_WORLD));
574c4762a1bSJed Brown
5759566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(networkdm, &X));
5769566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &F));
577c4762a1bSJed Brown
5789566063dSJacob Faibussowitsch PetscCall(DMCreateMatrix(networkdm, &J));
5799566063dSJacob Faibussowitsch PetscCall(MatSetOption(J, MAT_NEW_NONZERO_ALLOCATION_ERR, PETSC_FALSE));
580c4762a1bSJed Brown
5819566063dSJacob Faibussowitsch PetscCall(SetInitialValues(networkdm, X, &User));
582c4762a1bSJed Brown
583c4762a1bSJed Brown /* HOOK UP SOLVER */
5849566063dSJacob Faibussowitsch PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
5859566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, networkdm));
5869566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(snes, F, FormFunction, &User));
5879566063dSJacob Faibussowitsch PetscCall(SNESSetJacobian(snes, J, J, FormJacobian, &User));
5889566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes));
589c4762a1bSJed Brown
5909566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, X));
5919566063dSJacob Faibussowitsch /* PetscCall(VecView(X,PETSC_VIEWER_STDOUT_WORLD)); */
592c4762a1bSJed Brown
5939566063dSJacob Faibussowitsch PetscCall(VecDestroy(&X));
5949566063dSJacob Faibussowitsch PetscCall(VecDestroy(&F));
5959566063dSJacob Faibussowitsch PetscCall(MatDestroy(&J));
596c4762a1bSJed Brown
5979566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&snes));
5989566063dSJacob Faibussowitsch PetscCall(DMDestroy(&networkdm));
599c4762a1bSJed Brown }
6009566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
601b122ec5aSJacob Faibussowitsch return 0;
602c4762a1bSJed Brown }
603c4762a1bSJed Brown
604c4762a1bSJed Brown /*TEST
605c4762a1bSJed Brown
606c4762a1bSJed Brown build:
607c4762a1bSJed Brown depends: PFReadData.c pffunctions.c
608dfd57a17SPierre Jolivet requires: !complex double defined(PETSC_HAVE_ATTRIBUTEALIGNED)
609c4762a1bSJed Brown
610c4762a1bSJed Brown test:
611c4762a1bSJed Brown args: -snes_rtol 1.e-3
612c4762a1bSJed Brown localrunfiles: poweroptions case9.m
6132bf73ac6SHong Zhang output_file: output/power_1.out
614c4762a1bSJed Brown
615c4762a1bSJed Brown test:
616c4762a1bSJed Brown suffix: 2
617c4762a1bSJed Brown args: -snes_rtol 1.e-3 -petscpartitioner_type simple
618c4762a1bSJed Brown nsize: 4
619c4762a1bSJed Brown localrunfiles: poweroptions case9.m
6202bf73ac6SHong Zhang output_file: output/power_1.out
621c4762a1bSJed Brown
622c4762a1bSJed Brown TEST*/
623