xref: /petsc/src/ts/tutorials/network/pipeInterface.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
142dc13f1SHong Zhang #include "wash.h"
242dc13f1SHong Zhang 
342dc13f1SHong Zhang /* Subroutines for Pipe                                  */
442dc13f1SHong Zhang /* -------------------------------------------------------*/
542dc13f1SHong Zhang 
642dc13f1SHong Zhang /*
742dc13f1SHong Zhang    PipeCreate - Create Pipe object.
842dc13f1SHong Zhang 
942dc13f1SHong Zhang    Input Parameters:
1042dc13f1SHong Zhang    comm - MPI communicator
1142dc13f1SHong Zhang 
1242dc13f1SHong Zhang    Output Parameter:
1342dc13f1SHong Zhang .  pipe - location to put the PIPE context
1442dc13f1SHong Zhang */
15d71ae5a4SJacob Faibussowitsch PetscErrorCode PipeCreate(MPI_Comm comm, Pipe *pipe)
16d71ae5a4SJacob Faibussowitsch {
1742dc13f1SHong Zhang   PetscFunctionBegin;
189566063dSJacob Faibussowitsch   PetscCall(PetscNew(pipe));
19*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
2042dc13f1SHong Zhang }
2142dc13f1SHong Zhang 
2242dc13f1SHong Zhang /*
2342dc13f1SHong Zhang    PipeDestroy - Destroy Pipe object.
2442dc13f1SHong Zhang 
2542dc13f1SHong Zhang    Input Parameters:
2642dc13f1SHong Zhang    pipe - Reference to pipe intended to be destroyed.
2742dc13f1SHong Zhang */
28d71ae5a4SJacob Faibussowitsch PetscErrorCode PipeDestroy(Pipe *pipe)
29d71ae5a4SJacob Faibussowitsch {
3042dc13f1SHong Zhang   PetscFunctionBegin;
31*3ba16761SJacob Faibussowitsch   if (!*pipe) PetscFunctionReturn(PETSC_SUCCESS);
3242dc13f1SHong Zhang 
339566063dSJacob Faibussowitsch   PetscCall(PipeDestroyJacobian(*pipe));
349566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&(*pipe)->x));
359566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&(*pipe)->da));
36*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3742dc13f1SHong Zhang }
3842dc13f1SHong Zhang 
3942dc13f1SHong Zhang /*
4042dc13f1SHong Zhang    PipeSetParameters - Set parameters for Pipe context
4142dc13f1SHong Zhang 
4242dc13f1SHong Zhang    Input Parameter:
4342dc13f1SHong Zhang +  pipe - PIPE object
4442dc13f1SHong Zhang .  length -
4542dc13f1SHong Zhang .  nnodes -
4642dc13f1SHong Zhang .  D -
4742dc13f1SHong Zhang .  a -
4842dc13f1SHong Zhang -  fric -
4942dc13f1SHong Zhang */
50d71ae5a4SJacob Faibussowitsch PetscErrorCode PipeSetParameters(Pipe pipe, PetscReal length, PetscReal D, PetscReal a, PetscReal fric)
51d71ae5a4SJacob Faibussowitsch {
5242dc13f1SHong Zhang   PetscFunctionBegin;
5342dc13f1SHong Zhang   pipe->length = length;
5442dc13f1SHong Zhang   pipe->D      = D;
5542dc13f1SHong Zhang   pipe->a      = a;
5642dc13f1SHong Zhang   pipe->fric   = fric;
57*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5842dc13f1SHong Zhang }
5942dc13f1SHong Zhang 
6042dc13f1SHong Zhang /*
6142dc13f1SHong Zhang     PipeSetUp - Set up pipe based on set parameters.
6242dc13f1SHong Zhang */
63d71ae5a4SJacob Faibussowitsch PetscErrorCode PipeSetUp(Pipe pipe)
64d71ae5a4SJacob Faibussowitsch {
6542dc13f1SHong Zhang   DMDALocalInfo info;
6642dc13f1SHong Zhang 
6742dc13f1SHong Zhang   PetscFunctionBegin;
689566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_SELF, DM_BOUNDARY_GHOSTED, pipe->nnodes, 2, 1, NULL, &pipe->da));
699566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(pipe->da));
709566063dSJacob Faibussowitsch   PetscCall(DMSetUp(pipe->da));
719566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(pipe->da, 0, "Q"));
729566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(pipe->da, 1, "H"));
739566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(pipe->da, 0, pipe->length, 0, 0, 0, 0));
749566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(pipe->da, &(pipe->x)));
7542dc13f1SHong Zhang 
769566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(pipe->da, &info));
7742dc13f1SHong Zhang 
7842dc13f1SHong Zhang   pipe->rad = pipe->D / 2;
7942dc13f1SHong Zhang   pipe->A   = PETSC_PI * pipe->rad * pipe->rad;
8042dc13f1SHong Zhang   pipe->R   = pipe->fric / (2 * pipe->D * pipe->A);
81*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8242dc13f1SHong Zhang }
8342dc13f1SHong Zhang 
8442dc13f1SHong Zhang /*
8542dc13f1SHong Zhang     PipeCreateJacobian - Create Jacobian matrix structures for a Pipe.
8642dc13f1SHong Zhang 
87c3339decSBarry Smith     Collective
8842dc13f1SHong Zhang 
8942dc13f1SHong Zhang     Input Parameter:
9042dc13f1SHong Zhang +   pipe - the Pipe object
9142dc13f1SHong Zhang -   Jin - array of three constructed Jacobian matrices to be reused. Set NULL if it is not available
9242dc13f1SHong Zhang 
9342dc13f1SHong Zhang     Output Parameter:
9442dc13f1SHong Zhang .   J  - array of three empty Jacobian matrices
9542dc13f1SHong Zhang 
9642dc13f1SHong Zhang     Level: beginner
9742dc13f1SHong Zhang */
98d71ae5a4SJacob Faibussowitsch PetscErrorCode PipeCreateJacobian(Pipe pipe, Mat *Jin, Mat *J[])
99d71ae5a4SJacob Faibussowitsch {
10042dc13f1SHong Zhang   Mat         *Jpipe;
10142dc13f1SHong Zhang   PetscInt     M, rows[2], cols[2], *nz;
10242dc13f1SHong Zhang   PetscScalar *aa;
10342dc13f1SHong Zhang 
10442dc13f1SHong Zhang   PetscFunctionBegin;
10542dc13f1SHong Zhang   if (Jin) {
10642dc13f1SHong Zhang     *J             = Jin;
10742dc13f1SHong Zhang     pipe->jacobian = Jin;
1089566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject)(Jin[0])));
109*3ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
11042dc13f1SHong Zhang   }
1119566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(3, &Jpipe));
11242dc13f1SHong Zhang 
11342dc13f1SHong Zhang   /* Jacobian for this pipe */
1149566063dSJacob Faibussowitsch   PetscCall(DMSetMatrixStructureOnly(pipe->da, PETSC_TRUE));
1159566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(pipe->da, &Jpipe[0]));
1169566063dSJacob Faibussowitsch   PetscCall(DMSetMatrixStructureOnly(pipe->da, PETSC_FALSE));
11742dc13f1SHong Zhang 
11842dc13f1SHong Zhang   /* Jacobian for upstream vertex */
1199566063dSJacob Faibussowitsch   PetscCall(MatGetSize(Jpipe[0], &M, NULL));
1209566063dSJacob Faibussowitsch   PetscCall(PetscCalloc2(M, &nz, 4, &aa));
12142dc13f1SHong Zhang 
1229566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &Jpipe[1]));
1239566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Jpipe[1], PETSC_DECIDE, PETSC_DECIDE, M, 2));
1249566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Jpipe[1]));
1259566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Jpipe[1], MAT_STRUCTURE_ONLY, PETSC_TRUE));
1269371c9d4SSatish Balay   nz[0]   = 2;
1279371c9d4SSatish Balay   nz[1]   = 2;
1289371c9d4SSatish Balay   rows[0] = 0;
1299371c9d4SSatish Balay   rows[1] = 1;
1309371c9d4SSatish Balay   cols[0] = 0;
1319371c9d4SSatish Balay   cols[1] = 1;
1329566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Jpipe[1], 0, nz));
1339566063dSJacob Faibussowitsch   PetscCall(MatSetValues(Jpipe[1], 2, rows, 2, cols, aa, INSERT_VALUES));
1349566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpipe[1], MAT_FINAL_ASSEMBLY));
1359566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpipe[1], MAT_FINAL_ASSEMBLY));
13642dc13f1SHong Zhang 
13742dc13f1SHong Zhang   /* Jacobian for downstream vertex */
1389566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &Jpipe[2]));
1399566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Jpipe[2], PETSC_DECIDE, PETSC_DECIDE, M, 2));
1409566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Jpipe[2]));
1419566063dSJacob Faibussowitsch   PetscCall(MatSetOption(Jpipe[2], MAT_STRUCTURE_ONLY, PETSC_TRUE));
1429371c9d4SSatish Balay   nz[0]     = 0;
1439371c9d4SSatish Balay   nz[1]     = 0;
1449371c9d4SSatish Balay   nz[M - 2] = 2;
1459371c9d4SSatish Balay   nz[M - 1] = 2;
1469371c9d4SSatish Balay   rows[0]   = M - 2;
1479371c9d4SSatish Balay   rows[1]   = M - 1;
1489566063dSJacob Faibussowitsch   PetscCall(MatSeqAIJSetPreallocation(Jpipe[2], 0, nz));
1499566063dSJacob Faibussowitsch   PetscCall(MatSetValues(Jpipe[2], 2, rows, 2, cols, aa, INSERT_VALUES));
1509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpipe[2], MAT_FINAL_ASSEMBLY));
1519566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpipe[2], MAT_FINAL_ASSEMBLY));
15242dc13f1SHong Zhang 
1539566063dSJacob Faibussowitsch   PetscCall(PetscFree2(nz, aa));
15442dc13f1SHong Zhang 
15542dc13f1SHong Zhang   *J             = Jpipe;
15642dc13f1SHong Zhang   pipe->jacobian = Jpipe;
157*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15842dc13f1SHong Zhang }
15942dc13f1SHong Zhang 
160d71ae5a4SJacob Faibussowitsch PetscErrorCode PipeDestroyJacobian(Pipe pipe)
161d71ae5a4SJacob Faibussowitsch {
16242dc13f1SHong Zhang   Mat     *Jpipe = pipe->jacobian;
16342dc13f1SHong Zhang   PetscInt i;
16442dc13f1SHong Zhang 
16542dc13f1SHong Zhang   PetscFunctionBegin;
16642dc13f1SHong Zhang   if (Jpipe) {
16748a46eb9SPierre Jolivet     for (i = 0; i < 3; i++) PetscCall(MatDestroy(&Jpipe[i]));
16842dc13f1SHong Zhang   }
1699566063dSJacob Faibussowitsch   PetscCall(PetscFree(Jpipe));
170*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17142dc13f1SHong Zhang }
17242dc13f1SHong Zhang 
17342dc13f1SHong Zhang /*
17442dc13f1SHong Zhang     JunctionCreateJacobian - Create Jacobian matrices for a vertex.
17542dc13f1SHong Zhang 
176c3339decSBarry Smith     Collective
17742dc13f1SHong Zhang 
17842dc13f1SHong Zhang     Input Parameter:
17942dc13f1SHong Zhang +   dm - the DMNetwork object
18042dc13f1SHong Zhang .   v - vertex point
18142dc13f1SHong Zhang -   Jin - Jacobian patterns created by JunctionCreateJacobianSample() for reuse
18242dc13f1SHong Zhang 
18342dc13f1SHong Zhang     Output Parameter:
18442dc13f1SHong Zhang .   J  - array of Jacobian matrices (see dmnetworkimpl.h)
18542dc13f1SHong Zhang 
18642dc13f1SHong Zhang     Level: beginner
18742dc13f1SHong Zhang */
188d71ae5a4SJacob Faibussowitsch PetscErrorCode JunctionCreateJacobian(DM dm, PetscInt v, Mat *Jin, Mat *J[])
189d71ae5a4SJacob Faibussowitsch {
19042dc13f1SHong Zhang   Mat            *Jv;
19142dc13f1SHong Zhang   PetscInt        nedges, e, i, M, N, *rows, *cols;
19242dc13f1SHong Zhang   PetscBool       isSelf;
19342dc13f1SHong Zhang   const PetscInt *edges, *cone;
19442dc13f1SHong Zhang   PetscScalar    *zeros;
19542dc13f1SHong Zhang 
19642dc13f1SHong Zhang   PetscFunctionBegin;
197a5b23f4aSJose E. Roman   /* Get array size of Jv */
1989566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetSupportingEdges(dm, v, &nedges, &edges));
19963a3b9bcSJacob Faibussowitsch   PetscCheck(nedges > 0, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "%" PetscInt_FMT " vertex, nedges %" PetscInt_FMT, v, nedges);
20042dc13f1SHong Zhang 
20142dc13f1SHong Zhang   /* two Jacobians for each connected edge: J(v,e) and J(v,vc); adding J(v,v), total 2*nedges+1 Jacobians */
2029566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(2 * nedges + 1, &Jv));
20342dc13f1SHong Zhang 
20442dc13f1SHong Zhang   /* Create dense zero block for this vertex: J[0] = Jacobian(v,v) */
2059566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetComponent(dm, v, -1, NULL, NULL, &M));
20663a3b9bcSJacob Faibussowitsch   PetscCheck(M == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "%" PetscInt_FMT " != 2", M);
2079566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(M, &rows, M, &cols, M * M, &zeros));
2089566063dSJacob Faibussowitsch   PetscCall(PetscArrayzero(zeros, M * M));
20942dc13f1SHong Zhang   for (i = 0; i < M; i++) rows[i] = i;
21042dc13f1SHong Zhang 
21142dc13f1SHong Zhang   for (e = 0; e < nedges; e++) {
21242dc13f1SHong Zhang     /* create Jv[2*e+1] = Jacobian(v,e), e: supporting edge */
2139566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetConnectedVertices(dm, edges[e], &cone));
21442dc13f1SHong Zhang     isSelf = (v == cone[0]) ? PETSC_TRUE : PETSC_FALSE;
21542dc13f1SHong Zhang 
21642dc13f1SHong Zhang     if (Jin) {
21742dc13f1SHong Zhang       if (isSelf) {
21842dc13f1SHong Zhang         Jv[2 * e + 1] = Jin[0];
21942dc13f1SHong Zhang       } else {
22042dc13f1SHong Zhang         Jv[2 * e + 1] = Jin[1];
22142dc13f1SHong Zhang       }
22242dc13f1SHong Zhang       Jv[2 * e + 2] = Jin[2];
2239566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)(Jv[2 * e + 1])));
2249566063dSJacob Faibussowitsch       PetscCall(PetscObjectReference((PetscObject)(Jv[2 * e + 2])));
22542dc13f1SHong Zhang     } else {
22642dc13f1SHong Zhang       /* create J(v,e) */
2279566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_SELF, &Jv[2 * e + 1]));
2289566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetComponent(dm, edges[e], -1, NULL, NULL, &N));
2299566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(Jv[2 * e + 1], PETSC_DECIDE, PETSC_DECIDE, M, N));
2309566063dSJacob Faibussowitsch       PetscCall(MatSetFromOptions(Jv[2 * e + 1]));
2319566063dSJacob Faibussowitsch       PetscCall(MatSetOption(Jv[2 * e + 1], MAT_STRUCTURE_ONLY, PETSC_TRUE));
2329566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(Jv[2 * e + 1], 2, NULL));
23342dc13f1SHong Zhang       if (N) {
23442dc13f1SHong Zhang         if (isSelf) { /* coupling at upstream */
23542dc13f1SHong Zhang           for (i = 0; i < 2; i++) cols[i] = i;
23642dc13f1SHong Zhang         } else { /* coupling at downstream */
2379371c9d4SSatish Balay           cols[0] = N - 2;
2389371c9d4SSatish Balay           cols[1] = N - 1;
23942dc13f1SHong Zhang         }
2409566063dSJacob Faibussowitsch         PetscCall(MatSetValues(Jv[2 * e + 1], 2, rows, 2, cols, zeros, INSERT_VALUES));
24142dc13f1SHong Zhang       }
2429566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(Jv[2 * e + 1], MAT_FINAL_ASSEMBLY));
2439566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(Jv[2 * e + 1], MAT_FINAL_ASSEMBLY));
24442dc13f1SHong Zhang 
24542dc13f1SHong Zhang       /* create Jv[2*e+2] = Jacobian(v,vc), vc: connected vertex.
24642dc13f1SHong Zhang        In WashNetwork, v and vc are not connected, thus Jacobian(v,vc) is empty */
2479566063dSJacob Faibussowitsch       PetscCall(MatCreate(PETSC_COMM_SELF, &Jv[2 * e + 2]));
2489566063dSJacob Faibussowitsch       PetscCall(MatSetSizes(Jv[2 * e + 2], PETSC_DECIDE, PETSC_DECIDE, M, M)); /* empty matrix, sizes can be arbitrary */
2499566063dSJacob Faibussowitsch       PetscCall(MatSetFromOptions(Jv[2 * e + 2]));
2509566063dSJacob Faibussowitsch       PetscCall(MatSetOption(Jv[2 * e + 2], MAT_STRUCTURE_ONLY, PETSC_TRUE));
2519566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(Jv[2 * e + 2], 1, NULL));
2529566063dSJacob Faibussowitsch       PetscCall(MatAssemblyBegin(Jv[2 * e + 2], MAT_FINAL_ASSEMBLY));
2539566063dSJacob Faibussowitsch       PetscCall(MatAssemblyEnd(Jv[2 * e + 2], MAT_FINAL_ASSEMBLY));
25442dc13f1SHong Zhang     }
25542dc13f1SHong Zhang   }
2569566063dSJacob Faibussowitsch   PetscCall(PetscFree3(rows, cols, zeros));
25742dc13f1SHong Zhang 
25842dc13f1SHong Zhang   *J = Jv;
259*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
26042dc13f1SHong Zhang }
26142dc13f1SHong Zhang 
262d71ae5a4SJacob Faibussowitsch PetscErrorCode JunctionDestroyJacobian(DM dm, PetscInt v, Junction junc)
263d71ae5a4SJacob Faibussowitsch {
26442dc13f1SHong Zhang   Mat            *Jv = junc->jacobian;
26542dc13f1SHong Zhang   const PetscInt *edges;
26642dc13f1SHong Zhang   PetscInt        nedges, e;
26742dc13f1SHong Zhang 
26842dc13f1SHong Zhang   PetscFunctionBegin;
269*3ba16761SJacob Faibussowitsch   if (!Jv) PetscFunctionReturn(PETSC_SUCCESS);
27042dc13f1SHong Zhang 
2719566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetSupportingEdges(dm, v, &nedges, &edges));
27242dc13f1SHong Zhang   for (e = 0; e < nedges; e++) {
2739566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Jv[2 * e + 1]));
2749566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Jv[2 * e + 2]));
27542dc13f1SHong Zhang   }
2769566063dSJacob Faibussowitsch   PetscCall(PetscFree(Jv));
277*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
27842dc13f1SHong Zhang }
279