xref: /petsc/src/ts/tutorials/network/pipeImpls.c (revision 31d78bcd2b98084dc1368b20eb1129c8b9fb39fe)
142dc13f1SHong Zhang #include "pipe.h"
242dc13f1SHong Zhang 
342dc13f1SHong Zhang /* Initial Function for PIPE       */
442dc13f1SHong Zhang /*-------------------------------- */
542dc13f1SHong Zhang /*
642dc13f1SHong Zhang      Q(x) = Q0 (constant)
742dc13f1SHong Zhang      H(x) = H0 - (R/gA) Q0*|Q0|* x
842dc13f1SHong Zhang  */
942dc13f1SHong Zhang /* ----------------------------------- */
PipeComputeSteadyState(Pipe pipe,PetscScalar Q0,PetscScalar H0)10d71ae5a4SJacob Faibussowitsch PetscErrorCode PipeComputeSteadyState(Pipe pipe, PetscScalar Q0, PetscScalar H0)
11d71ae5a4SJacob Faibussowitsch {
1242dc13f1SHong Zhang   DM           cda;
1342dc13f1SHong Zhang   PipeField   *x;
1442dc13f1SHong Zhang   PetscInt     i, start, n;
1542dc13f1SHong Zhang   Vec          local;
1642dc13f1SHong Zhang   PetscScalar *coords, c = pipe->R / (GRAV * pipe->A);
1742dc13f1SHong Zhang 
1842dc13f1SHong Zhang   PetscFunctionBegin;
199566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(pipe->da, &cda));
209566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(pipe->da, &local));
219566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(pipe->da, pipe->x, &x));
229566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cda, local, &coords));
239566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(pipe->da, &start, 0, 0, &n, 0, 0));
2442dc13f1SHong Zhang 
2542dc13f1SHong Zhang   for (i = start; i < start + n; i++) {
2642dc13f1SHong Zhang     x[i].q = Q0;
2742dc13f1SHong Zhang     x[i].h = H0 - c * Q0 * PetscAbsScalar(Q0) * coords[i];
2842dc13f1SHong Zhang   }
2942dc13f1SHong Zhang 
309566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(pipe->da, pipe->x, &x));
319566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cda, local, &coords));
32*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
3342dc13f1SHong Zhang }
3442dc13f1SHong Zhang 
35d5b43468SJose E. Roman /* Function evaluations for PIPE    */
3642dc13f1SHong Zhang /*-------------------------------- */
3742dc13f1SHong Zhang /* consider using a one-sided higher order fd derivative at boundary. */
dqdx(PipeField * x,PetscInt i,PetscInt ilast,PetscReal dx)38d71ae5a4SJacob Faibussowitsch static inline PetscScalar dqdx(PipeField *x, PetscInt i, PetscInt ilast, PetscReal dx)
39d71ae5a4SJacob Faibussowitsch {
4042dc13f1SHong Zhang   if (i == 0) {
4142dc13f1SHong Zhang     return (x[i + 1].q - x[i].q) / dx;
4242dc13f1SHong Zhang   } else if (i == ilast) {
4342dc13f1SHong Zhang     return (x[i].q - x[i - 1].q) / dx;
4442dc13f1SHong Zhang   } else {
4542dc13f1SHong Zhang     return (x[i + 1].q - x[i - 1].q) / (2 * dx);
4642dc13f1SHong Zhang   }
4742dc13f1SHong Zhang }
4842dc13f1SHong Zhang 
dhdx(PipeField * x,PetscInt i,PetscInt ilast,PetscReal dx)49d71ae5a4SJacob Faibussowitsch static inline PetscScalar dhdx(PipeField *x, PetscInt i, PetscInt ilast, PetscReal dx)
50d71ae5a4SJacob Faibussowitsch {
5142dc13f1SHong Zhang   if (i == 0) {
5242dc13f1SHong Zhang     return (x[i + 1].h - x[i].h) / dx;
5342dc13f1SHong Zhang   } else if (i == ilast) {
5442dc13f1SHong Zhang     return (x[i].h - x[i - 1].h) / dx;
5542dc13f1SHong Zhang   } else {
5642dc13f1SHong Zhang     return (x[i + 1].h - x[i - 1].h) / (2 * dx);
5742dc13f1SHong Zhang   }
5842dc13f1SHong Zhang }
5942dc13f1SHong Zhang 
PipeIFunctionLocal_Lax(DMDALocalInfo * info,PetscReal ptime,PipeField * x,PipeField * xdot,PetscScalar * f,Pipe pipe)60d71ae5a4SJacob Faibussowitsch PetscErrorCode PipeIFunctionLocal_Lax(DMDALocalInfo *info, PetscReal ptime, PipeField *x, PipeField *xdot, PetscScalar *f, Pipe pipe)
61d71ae5a4SJacob Faibussowitsch {
6242dc13f1SHong Zhang   PetscInt    i, start, n, ilast;
6342dc13f1SHong Zhang   PetscReal   a = pipe->a, A = pipe->A, R = pipe->R, c = a * a / (GRAV * A);
6442dc13f1SHong Zhang   PetscReal   dx = pipe->length / (info->mx - 1), dt = pipe->dt;
6542dc13f1SHong Zhang   PetscScalar qavg, xold_i, ha, hb, qa, qb;
6642dc13f1SHong Zhang   PipeField  *xold = pipe->xold;
6742dc13f1SHong Zhang 
6842dc13f1SHong Zhang   PetscFunctionBegin;
699566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(pipe->da, &start, 0, 0, &n, 0, 0));
7042dc13f1SHong Zhang 
7142dc13f1SHong Zhang   /* interior and boundary */
7242dc13f1SHong Zhang   ilast = start + n - 1;
7342dc13f1SHong Zhang   for (i = start + 1; i < start + n - 1; i++) {
7442dc13f1SHong Zhang     qavg = (xold[i + 1].q + xold[i - 1].q) / 2.0;
759371c9d4SSatish Balay     qa   = xold[i - 1].q;
769371c9d4SSatish Balay     qb   = xold[i + 1].q;
779371c9d4SSatish Balay     ha   = xold[i - 1].h;
789371c9d4SSatish Balay     hb   = xold[i + 1].h;
7942dc13f1SHong Zhang 
8042dc13f1SHong Zhang     /* xdot[i].q = (x[i].q - old_i)/dt */
8142dc13f1SHong Zhang     xold_i             = 0.5 * (qa + qb);
8242dc13f1SHong Zhang     f[2 * (i - 1) + 2] = (x[i].q - xold_i) + dt * (GRAV * pipe->A * dhdx(xold, i, ilast, dx) + pipe->R * qavg * PetscAbsScalar(qavg));
8342dc13f1SHong Zhang 
8442dc13f1SHong Zhang     /* xdot[i].h = (x[i].h - xold_i)/dt */
8542dc13f1SHong Zhang     xold_i             = 0.5 * (ha + hb);
8642dc13f1SHong Zhang     f[2 * (i - 1) + 3] = (x[i].h - xold_i) + dt * c * dqdx(xold, i, ilast, dx);
8742dc13f1SHong Zhang   }
8842dc13f1SHong Zhang 
8942dc13f1SHong Zhang   /* Characteristic equations */
9042dc13f1SHong Zhang   f[start + 1] = x[start].q - xold[start + 1].q - ((GRAV * A) / a) * (x[start].h - xold[start + 1].h) + dt * R * xold[start + 1].q * PetscAbsScalar(xold[start + 1].q);
9142dc13f1SHong Zhang   f[2 * ilast] = x[ilast].q - xold[ilast - 1].q + ((GRAV * A) / a) * (x[ilast].h - xold[ilast - 1].h) + dt * R * xold[ilast - 1].q * PetscAbsScalar(xold[ilast - 1].q);
92*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9342dc13f1SHong Zhang }
94