xref: /petsc/src/ts/tutorials/power_grid/ex8.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "Time-dependent PDE in 2d for calculating joint PDF. \n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown    p_t = -x_t*p_x -y_t*p_y + f(t)*p_yy
4c4762a1bSJed Brown    xmin < x < xmax, ymin < y < ymax;
5c4762a1bSJed Brown 
6c4762a1bSJed Brown    Boundary conditions:
7c4762a1bSJed Brown    Zero dirichlet in y using ghosted values
8c4762a1bSJed Brown    Periodic in x
9c4762a1bSJed Brown 
10c4762a1bSJed Brown    Note that x_t and y_t in the above are given functions of x and y; they are not derivatives of x and y.
11c4762a1bSJed Brown    x_t = (y - ws)
12c4762a1bSJed Brown    y_t = (ws/2H)*(Pm - Pmax*sin(x) - D*(w - ws))
13c4762a1bSJed Brown 
14c4762a1bSJed Brown    In this example, we can see the effect of a fault, that zeroes the electrical power output
15c4762a1bSJed Brown    Pmax*sin(x), on the PDF. The fault on/off times can be controlled by options -tf and -tcl respectively.
16c4762a1bSJed Brown 
17c4762a1bSJed Brown */
18c4762a1bSJed Brown 
19c4762a1bSJed Brown #include <petscdm.h>
20c4762a1bSJed Brown #include <petscdmda.h>
21c4762a1bSJed Brown #include <petscts.h>
22c4762a1bSJed Brown 
23c4762a1bSJed Brown /*
24c4762a1bSJed Brown    User-defined data structures and routines
25c4762a1bSJed Brown */
26c4762a1bSJed Brown typedef struct {
27c4762a1bSJed Brown   PetscScalar ws;   /* Synchronous speed */
28c4762a1bSJed Brown   PetscScalar H;    /* Inertia constant */
29c4762a1bSJed Brown   PetscScalar D;    /* Damping constant */
30c4762a1bSJed Brown   PetscScalar Pmax,Pmax_s; /* Maximum power output of generator */
31c4762a1bSJed Brown   PetscScalar PM_min; /* Mean mechanical power input */
32c4762a1bSJed Brown   PetscScalar lambda; /* correlation time */
33c4762a1bSJed Brown   PetscScalar q;      /* noise strength */
34c4762a1bSJed Brown   PetscScalar mux;    /* Initial average angle */
35c4762a1bSJed Brown   PetscScalar sigmax; /* Standard deviation of initial angle */
36c4762a1bSJed Brown   PetscScalar muy;    /* Average speed */
37c4762a1bSJed Brown   PetscScalar sigmay; /* standard deviation of initial speed */
38c4762a1bSJed Brown   PetscScalar rho;    /* Cross-correlation coefficient */
39c4762a1bSJed Brown   PetscScalar xmin;   /* left boundary of angle */
40c4762a1bSJed Brown   PetscScalar xmax;   /* right boundary of angle */
41c4762a1bSJed Brown   PetscScalar ymin;   /* bottom boundary of speed */
42c4762a1bSJed Brown   PetscScalar ymax;   /* top boundary of speed */
43c4762a1bSJed Brown   PetscScalar dx;     /* x step size */
44c4762a1bSJed Brown   PetscScalar dy;     /* y step size */
45c4762a1bSJed Brown   PetscScalar disper_coe; /* Dispersion coefficient */
46c4762a1bSJed Brown   DM          da;
47c4762a1bSJed Brown   PetscInt    st_width; /* Stencil width */
48c4762a1bSJed Brown   DMBoundaryType bx; /* x boundary type */
49c4762a1bSJed Brown   DMBoundaryType by; /* y boundary type */
50c4762a1bSJed Brown   PetscReal        tf,tcl; /* Fault incidence and clearing times */
51c4762a1bSJed Brown } AppCtx;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown PetscErrorCode Parameter_settings(AppCtx*);
54c4762a1bSJed Brown PetscErrorCode ini_bou(Vec,AppCtx*);
55c4762a1bSJed Brown PetscErrorCode IFunction(TS,PetscReal,Vec,Vec,Vec,void*);
56c4762a1bSJed Brown PetscErrorCode IJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
57c4762a1bSJed Brown PetscErrorCode PostStep(TS);
58c4762a1bSJed Brown 
59c4762a1bSJed Brown int main(int argc, char **argv)
60c4762a1bSJed Brown {
61c4762a1bSJed Brown   Vec         x;                /* Solution vector */
62c4762a1bSJed Brown   TS          ts;               /* Time-stepping context */
63c4762a1bSJed Brown   AppCtx      user;             /* Application context */
64c4762a1bSJed Brown   PetscViewer viewer;
65c4762a1bSJed Brown 
66*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,"petscopt_ex8", help));
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* Get physics and time parameters */
695f80ce2aSJacob Faibussowitsch   CHKERRQ(Parameter_settings(&user));
70c4762a1bSJed Brown   /* Create a 2D DA with dof = 1 */
715f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,user.bx,user.by,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,user.st_width,NULL,NULL,&user.da));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(user.da));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(user.da));
74c4762a1bSJed Brown   /* Set x and y coordinates */
755f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(user.da,user.xmin,user.xmax,user.ymin,user.ymax,0,0));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetCoordinateName(user.da,0,"X - the angle"));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetCoordinateName(user.da,1,"Y - the speed"));
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* Get global vector x from DM  */
805f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(user.da,&x));
81c4762a1bSJed Brown 
825f80ce2aSJacob Faibussowitsch   CHKERRQ(ini_bou(x,&user));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ini_x",FILE_MODE_WRITE,&viewer));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,viewer));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
86c4762a1bSJed Brown 
875f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts,user.da));
895f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSARKIMEX));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIFunction(ts,NULL,IFunction,&user));
925f80ce2aSJacob Faibussowitsch   /*  CHKERRQ(TSSetIJacobian(ts,NULL,NULL,IJacobian,&user)); */
935f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetApplicationContext(ts,&user));
945f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,.005));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetPostStep(ts,PostStep));
975f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,x));
98c4762a1bSJed Brown 
995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"fin_x",FILE_MODE_WRITE,&viewer));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,viewer));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&viewer));
102c4762a1bSJed Brown 
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&user.da));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
106*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
107*b122ec5aSJacob Faibussowitsch   return 0;
108c4762a1bSJed Brown }
109c4762a1bSJed Brown 
110c4762a1bSJed Brown PetscErrorCode PostStep(TS ts)
111c4762a1bSJed Brown {
112c4762a1bSJed Brown   Vec          X;
113c4762a1bSJed Brown   AppCtx      *user;
114c4762a1bSJed Brown   PetscReal    t;
115c4762a1bSJed Brown   PetscScalar  asum;
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   PetscFunctionBegin;
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetApplicationContext(ts,&user));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTime(ts,&t));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolution(ts,&X));
121c4762a1bSJed Brown   /*
122c4762a1bSJed Brown   if (t >= .2) {
1235f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetSolution(ts,&X));
1245f80ce2aSJacob Faibussowitsch     CHKERRQ(VecView(X,PETSC_VIEWER_BINARY_WORLD));
125c4762a1bSJed Brown     exit(0);
126c4762a1bSJed Brown      results in initial conditions after fault in binaryoutput
127c4762a1bSJed Brown   }*/
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   if ((t > user->tf) && (t < user->tcl)) user->Pmax = 0.0; /* A short-circuit that drives the electrical power output (Pmax*sin(delta)) to zero */
130c4762a1bSJed Brown   else user->Pmax = user->Pmax_s;
131c4762a1bSJed Brown 
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSum(X,&asum));
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"sum(p) at t = %f = %f\n",(double)t,(double)(asum)));
134c4762a1bSJed Brown   PetscFunctionReturn(0);
135c4762a1bSJed Brown }
136c4762a1bSJed Brown 
137c4762a1bSJed Brown PetscErrorCode ini_bou(Vec X,AppCtx* user)
138c4762a1bSJed Brown {
139c4762a1bSJed Brown   DM            cda;
140c4762a1bSJed Brown   DMDACoor2d  **coors;
141c4762a1bSJed Brown   PetscScalar **p;
142c4762a1bSJed Brown   Vec           gc;
143c4762a1bSJed Brown   PetscInt      M,N,Ir,J;
144c4762a1bSJed Brown   PetscMPIInt   rank;
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   PetscFunctionBeginUser;
1475f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
149c4762a1bSJed Brown   user->dx = (user->xmax - user->xmin)/(M-1); user->dy = (user->ymax - user->ymin)/(N-1);
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(user->da,&cda));
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinates(user->da,&gc));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(cda,gc,&coors));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(user->da,X,&p));
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   /* Point mass at (mux,muy) */
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Original user->mux = %f, user->muy = %f\n",user->mux,user->muy));
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLogicalCoordinate(user->da,user->mux,user->muy,0.0,&Ir,&J,NULL,&user->mux,&user->muy,NULL));
158c4762a1bSJed Brown   user->PM_min = user->Pmax*PetscSinScalar(user->mux);
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Corrected user->mux = %f, user->muy = %f user->PM_min = %f,user->dx = %f\n",user->mux,user->muy,user->PM_min,user->dx));
160c4762a1bSJed Brown   if (Ir > -1 && J > -1) {
161c4762a1bSJed Brown     p[J][Ir] = 1.0;
162c4762a1bSJed Brown   }
163c4762a1bSJed Brown 
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(cda,gc,&coors));
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(user->da,X,&p));
166c4762a1bSJed Brown   PetscFunctionReturn(0);
167c4762a1bSJed Brown }
168c4762a1bSJed Brown 
169c4762a1bSJed Brown /* First advection term */
170c4762a1bSJed Brown PetscErrorCode adv1(PetscScalar **p,PetscScalar y,PetscInt i,PetscInt j,PetscInt M,PetscScalar *p1,AppCtx *user)
171c4762a1bSJed Brown {
172c4762a1bSJed Brown   PetscScalar f,fpos,fneg;
173c4762a1bSJed Brown   PetscFunctionBegin;
174c4762a1bSJed Brown   f   =  (y - user->ws);
175c4762a1bSJed Brown   fpos = PetscMax(f,0);
176c4762a1bSJed Brown   fneg = PetscMin(f,0);
177c4762a1bSJed Brown   if (user->st_width == 1) {
178c4762a1bSJed Brown     *p1 = fpos*(p[j][i] - p[j][i-1])/user->dx + fneg*(p[j][i+1] - p[j][i])/user->dx;
179c4762a1bSJed Brown   } else if (user->st_width == 2) {
180c4762a1bSJed Brown     *p1 = fpos*(3*p[j][i] - 4*p[j][i-1] + p[j][i-2])/(2*user->dx) + fneg*(-p[j][i+2] + 4*p[j][i+1] - 3*p[j][i])/(2*user->dx);
181c4762a1bSJed Brown   } else if (user->st_width == 3) {
182c4762a1bSJed Brown     *p1 = fpos*(2*p[j][i+1] + 3*p[j][i] - 6*p[j][i-1] + p[j][i-2])/(6*user->dx) + fneg*(-p[j][i+2] + 6*p[j][i+1] - 3*p[j][i] - 2*p[j][i-1])/(6*user->dx);
183c4762a1bSJed Brown   }
184c4762a1bSJed Brown   /* *p1 = f*(p[j][i+1] - p[j][i-1])/user->dx;*/
185c4762a1bSJed Brown   PetscFunctionReturn(0);
186c4762a1bSJed Brown }
187c4762a1bSJed Brown 
188c4762a1bSJed Brown /* Second advection term */
189c4762a1bSJed Brown PetscErrorCode adv2(PetscScalar **p,PetscScalar x,PetscScalar y,PetscInt i,PetscInt j,PetscInt N,PetscScalar *p2,AppCtx *user)
190c4762a1bSJed Brown {
191c4762a1bSJed Brown   PetscScalar f,fpos,fneg;
192c4762a1bSJed Brown   PetscFunctionBegin;
193c4762a1bSJed Brown   f   = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(x) - user->D*(y - user->ws));
194c4762a1bSJed Brown   fpos = PetscMax(f,0);
195c4762a1bSJed Brown   fneg = PetscMin(f,0);
196c4762a1bSJed Brown   if (user->st_width == 1) {
197c4762a1bSJed Brown     *p2 = fpos*(p[j][i] - p[j-1][i])/user->dy + fneg*(p[j+1][i] - p[j][i])/user->dy;
198c4762a1bSJed Brown   } else if (user->st_width ==2) {
199c4762a1bSJed Brown     *p2 = fpos*(3*p[j][i] - 4*p[j-1][i] + p[j-2][i])/(2*user->dy) + fneg*(-p[j+2][i] + 4*p[j+1][i] - 3*p[j][i])/(2*user->dy);
200c4762a1bSJed Brown   } else if (user->st_width == 3) {
201c4762a1bSJed Brown     *p2 = fpos*(2*p[j+1][i] + 3*p[j][i] - 6*p[j-1][i] + p[j-2][i])/(6*user->dy) + fneg*(-p[j+2][i] + 6*p[j+1][i] - 3*p[j][i] - 2*p[j-1][i])/(6*user->dy);
202c4762a1bSJed Brown   }
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   /* *p2 = f*(p[j+1][i] - p[j-1][i])/user->dy;*/
205c4762a1bSJed Brown   PetscFunctionReturn(0);
206c4762a1bSJed Brown }
207c4762a1bSJed Brown 
208c4762a1bSJed Brown /* Diffusion term */
209c4762a1bSJed Brown PetscErrorCode diffuse(PetscScalar **p,PetscInt i,PetscInt j,PetscReal t,PetscScalar *p_diff,AppCtx * user)
210c4762a1bSJed Brown {
211c4762a1bSJed Brown   PetscFunctionBeginUser;
212c4762a1bSJed Brown   if (user->st_width == 1) {
213c4762a1bSJed Brown     *p_diff = user->disper_coe*((p[j-1][i] - 2*p[j][i] + p[j+1][i])/(user->dy*user->dy));
214c4762a1bSJed Brown   } else if (user->st_width == 2) {
215c4762a1bSJed Brown     *p_diff = user->disper_coe*((-p[j-2][i] + 16*p[j-1][i] - 30*p[j][i] + 16*p[j+1][i] - p[j+2][i])/(12.0*user->dy*user->dy));
216c4762a1bSJed Brown   } else if (user->st_width == 3) {
217c4762a1bSJed Brown     *p_diff = user->disper_coe*((2*p[j-3][i] - 27*p[j-2][i] + 270*p[j-1][i] - 490*p[j][i] + 270*p[j+1][i] - 27*p[j+2][i] + 2*p[j+3][i])/(180.0*user->dy*user->dy));
218c4762a1bSJed Brown   }
219c4762a1bSJed Brown   PetscFunctionReturn(0);
220c4762a1bSJed Brown }
221c4762a1bSJed Brown 
222c4762a1bSJed Brown PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
223c4762a1bSJed Brown {
224c4762a1bSJed Brown   AppCtx       *user   = (AppCtx*)ctx;
225c4762a1bSJed Brown   DM            cda;
226c4762a1bSJed Brown   DMDACoor2d  **coors;
227c4762a1bSJed Brown   PetscScalar **p,**f,**pdot;
228c4762a1bSJed Brown   PetscInt      i,j;
229c4762a1bSJed Brown   PetscInt      xs,ys,xm,ym,M,N;
230c4762a1bSJed Brown   Vec           localX,gc,localXdot;
231c4762a1bSJed Brown   PetscScalar   p_adv1 = 0.0,p_adv2 = 0.0,p_diff = 0;
232c4762a1bSJed Brown   PetscScalar   diffuse1,gamma;
233c4762a1bSJed Brown 
234c4762a1bSJed Brown   PetscFunctionBeginUser;
2355f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
2365f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(user->da,&cda));
2375f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0));
238c4762a1bSJed Brown 
2395f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(user->da,&localX));
2405f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(user->da,&localXdot));
241c4762a1bSJed Brown 
2425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX));
2435f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX));
2445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(user->da,Xdot,INSERT_VALUES,localXdot));
2455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(user->da,Xdot,INSERT_VALUES,localXdot));
246c4762a1bSJed Brown 
2475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(user->da,&gc));
248c4762a1bSJed Brown 
2495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(cda,gc,&coors));
2505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(user->da,localX,&p));
2515f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(user->da,localXdot,&pdot));
2525f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(user->da,F,&f));
253c4762a1bSJed Brown 
254c4762a1bSJed Brown   gamma = user->D*user->ws/(2*user->H);
255c4762a1bSJed Brown   diffuse1 = user->lambda*user->lambda*user->q/(user->lambda*gamma+1)*(1.0 - PetscExpScalar(-t*(gamma+1.0)/user->lambda));
256c4762a1bSJed Brown   user->disper_coe = user->ws*user->ws/(4*user->H*user->H)*diffuse1;
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   for (i=xs; i < xs+xm; i++) {
259c4762a1bSJed Brown     for (j=ys; j < ys+ym; j++) {
2605f80ce2aSJacob Faibussowitsch       CHKERRQ(adv1(p,coors[j][i].y,i,j,M,&p_adv1,user));
2615f80ce2aSJacob Faibussowitsch       CHKERRQ(adv2(p,coors[j][i].x,coors[j][i].y,i,j,N,&p_adv2,user));
2625f80ce2aSJacob Faibussowitsch       CHKERRQ(diffuse(p,i,j,t,&p_diff,user));
263c4762a1bSJed Brown       f[j][i] = -p_adv1 - p_adv2  + p_diff - pdot[j][i];
264c4762a1bSJed Brown     }
265c4762a1bSJed Brown   }
2665f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(user->da,localX,&p));
2675f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(user->da,localX,&pdot));
2685f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(user->da,&localX));
2695f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(user->da,&localXdot));
2705f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(user->da,F,&f));
2715f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(cda,gc,&coors));
272c4762a1bSJed Brown 
273c4762a1bSJed Brown   PetscFunctionReturn(0);
274c4762a1bSJed Brown }
275c4762a1bSJed Brown 
276c4762a1bSJed Brown PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ctx)
277c4762a1bSJed Brown {
278c4762a1bSJed Brown   AppCtx         *user=(AppCtx*)ctx;
279c4762a1bSJed Brown   DM             cda;
280c4762a1bSJed Brown   DMDACoor2d     **coors;
281c4762a1bSJed Brown   PetscInt       i,j;
282c4762a1bSJed Brown   PetscInt       xs,ys,xm,ym,M,N;
283c4762a1bSJed Brown   Vec            gc;
284c4762a1bSJed Brown   PetscScalar    val[5],xi,yi;
285c4762a1bSJed Brown   MatStencil     row,col[5];
286c4762a1bSJed Brown   PetscScalar    c1,c3,c5,c1pos,c1neg,c3pos,c3neg;
287c4762a1bSJed Brown 
288c4762a1bSJed Brown   PetscFunctionBeginUser;
2895f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
2905f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(user->da,&cda));
2915f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0));
292c4762a1bSJed Brown 
2935f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(user->da,&gc));
2945f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(cda,gc,&coors));
295c4762a1bSJed Brown   for (i=xs; i < xs+xm; i++) {
296c4762a1bSJed Brown     for (j=ys; j < ys+ym; j++) {
297c4762a1bSJed Brown       PetscInt nc = 0;
298c4762a1bSJed Brown       xi = coors[j][i].x; yi = coors[j][i].y;
299c4762a1bSJed Brown       row.i = i; row.j = j;
300c4762a1bSJed Brown       c1        = (yi-user->ws)/user->dx;
301c4762a1bSJed Brown       c1pos    = PetscMax(c1,0);
302c4762a1bSJed Brown       c1neg    = PetscMin(c1,0);
303c4762a1bSJed Brown       c3        = (user->ws/(2.0*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(xi) - user->D*(yi - user->ws))/user->dy;
304c4762a1bSJed Brown       c3pos    = PetscMax(c3,0);
305c4762a1bSJed Brown       c3neg    = PetscMin(c3,0);
306c4762a1bSJed Brown       c5        = (PetscPowScalar((user->lambda*user->ws)/(2*user->H),2)*user->q*(1.0-PetscExpScalar(-t/user->lambda)))/(user->dy*user->dy);
307c4762a1bSJed Brown       col[nc].i = i-1; col[nc].j = j;   val[nc++] = c1pos;
308c4762a1bSJed Brown       col[nc].i = i+1; col[nc].j = j;   val[nc++] = -c1neg;
309c4762a1bSJed Brown       col[nc].i = i;   col[nc].j = j-1; val[nc++] = c3pos + c5;
310c4762a1bSJed Brown       col[nc].i = i;   col[nc].j = j+1; val[nc++] = -c3neg + c5;
311c4762a1bSJed Brown       col[nc].i = i;   col[nc].j = j;   val[nc++] = -c1pos + c1neg -c3pos + c3neg -2*c5 -a;
3125f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValuesStencil(Jpre,1,&row,nc,col,val,INSERT_VALUES));
313c4762a1bSJed Brown     }
314c4762a1bSJed Brown   }
3155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(cda,gc,&coors));
316c4762a1bSJed Brown 
3175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY));
3185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY));
319c4762a1bSJed Brown   if (J != Jpre) {
3205f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
3215f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
322c4762a1bSJed Brown   }
323c4762a1bSJed Brown   PetscFunctionReturn(0);
324c4762a1bSJed Brown }
325c4762a1bSJed Brown 
326c4762a1bSJed Brown PetscErrorCode Parameter_settings(AppCtx *user)
327c4762a1bSJed Brown {
328c4762a1bSJed Brown   PetscBool flg;
329c4762a1bSJed Brown 
330c4762a1bSJed Brown   PetscFunctionBeginUser;
331c4762a1bSJed Brown   /* Set default parameters */
332c4762a1bSJed Brown   user->ws     = 1.0;
333c4762a1bSJed Brown   user->H      = 5.0;
334c4762a1bSJed Brown   user->D      = 0.0;
335c4762a1bSJed Brown   user->Pmax = user->Pmax_s  = 2.1;
336c4762a1bSJed Brown   user->PM_min = 1.0;
337c4762a1bSJed Brown   user->lambda = 0.1;
338c4762a1bSJed Brown   user->q      = 1.0;
339c4762a1bSJed Brown   user->mux    = PetscAsinScalar(user->PM_min/user->Pmax);
340c4762a1bSJed Brown   user->sigmax = 0.1;
341c4762a1bSJed Brown   user->sigmay = 0.1;
342c4762a1bSJed Brown   user->rho    = 0.0;
343c4762a1bSJed Brown   user->xmin   = -PETSC_PI;
344c4762a1bSJed Brown   user->xmax   = PETSC_PI;
345c4762a1bSJed Brown   user->bx     = DM_BOUNDARY_PERIODIC;
346c4762a1bSJed Brown   user->by     = DM_BOUNDARY_GHOSTED;
347c4762a1bSJed Brown   user->tf = user->tcl = -1;
348c4762a1bSJed Brown   user->ymin   = -2.0;
349c4762a1bSJed Brown   user->ymax   = 2.0;
350c4762a1bSJed Brown   user->st_width = 1;
351c4762a1bSJed Brown 
3525f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-ws",&user->ws,&flg));
3535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-Inertia",&user->H,&flg));
3545f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-D",&user->D,&flg));
3555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-Pmax",&user->Pmax,&flg));
3565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-PM_min",&user->PM_min,&flg));
3575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-lambda",&user->lambda,&flg));
3585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-q",&user->q,&flg));
3595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-mux",&user->mux,&flg));
3605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-muy",&user->muy,&flg));
3615f80ce2aSJacob Faibussowitsch   if (flg == 0) user->muy = user->ws;
3625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-xmin",&user->xmin,&flg));
3635f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-xmax",&user->xmax,&flg));
3645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-ymin",&user->ymin,&flg));
3655f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-ymax",&user->ymax,&flg));
3665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-stencil_width",&user->st_width,&flg));
3675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetEnum(NULL,NULL,"-bx",DMBoundaryTypes,(PetscEnum*)&user->bx,&flg));
3685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetEnum(NULL,NULL,"-by",DMBoundaryTypes,(PetscEnum*)&user->by,&flg));
3695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-tf",&user->tf,&flg));
3705f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-tcl",&user->tcl,&flg));
371c4762a1bSJed Brown   PetscFunctionReturn(0);
372c4762a1bSJed Brown }
373c4762a1bSJed Brown 
374c4762a1bSJed Brown /*TEST
375c4762a1bSJed Brown 
376c4762a1bSJed Brown    build:
377c4762a1bSJed Brown       requires: !complex x
378c4762a1bSJed Brown 
379c4762a1bSJed Brown    test:
380c4762a1bSJed Brown       args: -ts_max_steps 1
381c4762a1bSJed Brown       localrunfiles: petscopt_ex8
382c4762a1bSJed Brown       timeoutfactor: 3
383c4762a1bSJed Brown 
384c4762a1bSJed Brown TEST*/
385