xref: /petsc/src/ts/tutorials/power_grid/ex6.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Time-dependent PDE in 2d for calculating joint PDF. \n";
2*c4762a1bSJed Brown /*
3*c4762a1bSJed Brown    p_t = -x_t*p_x -y_t*p_y + f(t)*p_yy
4*c4762a1bSJed Brown    xmin < x < xmax, ymin < y < ymax;
5*c4762a1bSJed Brown    x_t = (y - ws)  y_t = (ws/2H)*(Pm - Pmax*sin(x))
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown    Boundary conditions: -bc_type 0 => Zero dirichlet boundary
8*c4762a1bSJed Brown                         -bc_type 1 => Steady state boundary condition
9*c4762a1bSJed Brown    Steady state boundary condition found by setting p_t = 0
10*c4762a1bSJed Brown */
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown #include <petscdm.h>
13*c4762a1bSJed Brown #include <petscdmda.h>
14*c4762a1bSJed Brown #include <petscts.h>
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown /*
17*c4762a1bSJed Brown    User-defined data structures and routines
18*c4762a1bSJed Brown */
19*c4762a1bSJed Brown typedef struct {
20*c4762a1bSJed Brown   PetscScalar ws;   /* Synchronous speed */
21*c4762a1bSJed Brown   PetscScalar H;    /* Inertia constant */
22*c4762a1bSJed Brown   PetscScalar D;    /* Damping constant */
23*c4762a1bSJed Brown   PetscScalar Pmax; /* Maximum power output of generator */
24*c4762a1bSJed Brown   PetscScalar PM_min; /* Mean mechanical power input */
25*c4762a1bSJed Brown   PetscScalar lambda; /* correlation time */
26*c4762a1bSJed Brown   PetscScalar q;      /* noise strength */
27*c4762a1bSJed Brown   PetscScalar mux;    /* Initial average angle */
28*c4762a1bSJed Brown   PetscScalar sigmax; /* Standard deviation of initial angle */
29*c4762a1bSJed Brown   PetscScalar muy;    /* Average speed */
30*c4762a1bSJed Brown   PetscScalar sigmay; /* standard deviation of initial speed */
31*c4762a1bSJed Brown   PetscScalar rho;    /* Cross-correlation coefficient */
32*c4762a1bSJed Brown   PetscScalar t0;     /* Initial time */
33*c4762a1bSJed Brown   PetscScalar tmax;   /* Final time */
34*c4762a1bSJed Brown   PetscScalar xmin;   /* left boundary of angle */
35*c4762a1bSJed Brown   PetscScalar xmax;   /* right boundary of angle */
36*c4762a1bSJed Brown   PetscScalar ymin;   /* bottom boundary of speed */
37*c4762a1bSJed Brown   PetscScalar ymax;   /* top boundary of speed */
38*c4762a1bSJed Brown   PetscScalar dx;     /* x step size */
39*c4762a1bSJed Brown   PetscScalar dy;     /* y step size */
40*c4762a1bSJed Brown   PetscInt    bc; /* Boundary conditions */
41*c4762a1bSJed Brown   PetscScalar disper_coe; /* Dispersion coefficient */
42*c4762a1bSJed Brown   DM          da;
43*c4762a1bSJed Brown } AppCtx;
44*c4762a1bSJed Brown 
45*c4762a1bSJed Brown PetscErrorCode Parameter_settings(AppCtx*);
46*c4762a1bSJed Brown PetscErrorCode ini_bou(Vec,AppCtx*);
47*c4762a1bSJed Brown PetscErrorCode IFunction(TS,PetscReal,Vec,Vec,Vec,void*);
48*c4762a1bSJed Brown PetscErrorCode IJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
49*c4762a1bSJed Brown PetscErrorCode PostStep(TS);
50*c4762a1bSJed Brown 
51*c4762a1bSJed Brown int main(int argc, char **argv)
52*c4762a1bSJed Brown {
53*c4762a1bSJed Brown   PetscErrorCode ierr;
54*c4762a1bSJed Brown   Vec            x;  /* Solution vector */
55*c4762a1bSJed Brown   TS             ts;   /* Time-stepping context */
56*c4762a1bSJed Brown   AppCtx         user; /* Application context */
57*c4762a1bSJed Brown   Mat            J;
58*c4762a1bSJed Brown   PetscViewer    viewer;
59*c4762a1bSJed Brown 
60*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,"petscopt_ex6", help);if (ierr) return ierr;
61*c4762a1bSJed Brown   /* Get physics and time parameters */
62*c4762a1bSJed Brown   ierr = Parameter_settings(&user);CHKERRQ(ierr);
63*c4762a1bSJed Brown   /* Create a 2D DA with dof = 1 */
64*c4762a1bSJed Brown   ierr = DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&user.da);CHKERRQ(ierr);
65*c4762a1bSJed Brown   ierr = DMSetFromOptions(user.da);CHKERRQ(ierr);
66*c4762a1bSJed Brown   ierr = DMSetUp(user.da);CHKERRQ(ierr);
67*c4762a1bSJed Brown   /* Set x and y coordinates */
68*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(user.da,user.xmin,user.xmax,user.ymin,user.ymax,0.0,1.0);CHKERRQ(ierr);
69*c4762a1bSJed Brown 
70*c4762a1bSJed Brown   /* Get global vector x from DM  */
71*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(user.da,&x);CHKERRQ(ierr);
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown   ierr = ini_bou(x,&user);CHKERRQ(ierr);
74*c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ini_x",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
75*c4762a1bSJed Brown   ierr = VecView(x,viewer);CHKERRQ(ierr);
76*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
77*c4762a1bSJed Brown 
78*c4762a1bSJed Brown   /* Get Jacobian matrix structure from the da */
79*c4762a1bSJed Brown   ierr = DMSetMatType(user.da,MATAIJ);CHKERRQ(ierr);
80*c4762a1bSJed Brown   ierr = DMCreateMatrix(user.da,&J);CHKERRQ(ierr);
81*c4762a1bSJed Brown 
82*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
84*c4762a1bSJed Brown   ierr = TSSetIFunction(ts,NULL,IFunction,&user);CHKERRQ(ierr);
85*c4762a1bSJed Brown   ierr = TSSetIJacobian(ts,J,J,IJacobian,&user);CHKERRQ(ierr);
86*c4762a1bSJed Brown   ierr = TSSetApplicationContext(ts,&user);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,user.tmax);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = TSSetTime(ts,user.t0);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,.005);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
92*c4762a1bSJed Brown   ierr = TSSetPostStep(ts,PostStep);CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr = TSSolve(ts,x);CHKERRQ(ierr);
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"fin_x",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr);
96*c4762a1bSJed Brown   ierr = VecView(x,viewer);CHKERRQ(ierr);
97*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
98*c4762a1bSJed Brown 
99*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
100*c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
101*c4762a1bSJed Brown   ierr = DMDestroy(&user.da);CHKERRQ(ierr);
102*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
103*c4762a1bSJed Brown   ierr = PetscFinalize();
104*c4762a1bSJed Brown   return ierr;
105*c4762a1bSJed Brown }
106*c4762a1bSJed Brown 
107*c4762a1bSJed Brown PetscErrorCode PostStep(TS ts)
108*c4762a1bSJed Brown {
109*c4762a1bSJed Brown   PetscErrorCode ierr;
110*c4762a1bSJed Brown   Vec            X;
111*c4762a1bSJed Brown   AppCtx         *user;
112*c4762a1bSJed Brown   PetscScalar    sum;
113*c4762a1bSJed Brown   PetscReal      t;
114*c4762a1bSJed Brown 
115*c4762a1bSJed Brown   PetscFunctionBegin;
116*c4762a1bSJed Brown   ierr = TSGetApplicationContext(ts,&user);CHKERRQ(ierr);
117*c4762a1bSJed Brown   ierr = TSGetTime(ts,&t);CHKERRQ(ierr);
118*c4762a1bSJed Brown   ierr = TSGetSolution(ts,&X);CHKERRQ(ierr);
119*c4762a1bSJed Brown   ierr = VecSum(X,&sum);CHKERRQ(ierr);
120*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"sum(p)*dw*dtheta at t = %3.2f = %3.6f\n",(double)t,(double)sum*user->dx*user->dy);CHKERRQ(ierr);
121*c4762a1bSJed Brown   PetscFunctionReturn(0);
122*c4762a1bSJed Brown }
123*c4762a1bSJed Brown 
124*c4762a1bSJed Brown PetscErrorCode ini_bou(Vec X,AppCtx* user)
125*c4762a1bSJed Brown {
126*c4762a1bSJed Brown   PetscErrorCode ierr;
127*c4762a1bSJed Brown   DM             cda;
128*c4762a1bSJed Brown   DMDACoor2d     **coors;
129*c4762a1bSJed Brown   PetscScalar    **p;
130*c4762a1bSJed Brown   Vec            gc;
131*c4762a1bSJed Brown   PetscInt       i,j;
132*c4762a1bSJed Brown   PetscInt       xs,ys,xm,ym,M,N;
133*c4762a1bSJed Brown   PetscScalar    xi,yi;
134*c4762a1bSJed Brown   PetscScalar    sigmax=user->sigmax,sigmay=user->sigmay;
135*c4762a1bSJed Brown   PetscScalar    rho   =user->rho;
136*c4762a1bSJed Brown   PetscScalar    mux   =user->mux,muy=user->muy;
137*c4762a1bSJed Brown   PetscMPIInt    rank;
138*c4762a1bSJed Brown 
139*c4762a1bSJed Brown   PetscFunctionBeginUser;
140*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
141*c4762a1bSJed Brown   ierr = DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
142*c4762a1bSJed Brown   user->dx = (user->xmax - user->xmin)/(M-1); user->dy = (user->ymax - user->ymin)/(N-1);
143*c4762a1bSJed Brown   ierr = DMGetCoordinateDM(user->da,&cda);CHKERRQ(ierr);
144*c4762a1bSJed Brown   ierr = DMGetCoordinates(user->da,&gc);CHKERRQ(ierr);
145*c4762a1bSJed Brown   ierr = DMDAVecGetArray(cda,gc,&coors);CHKERRQ(ierr);
146*c4762a1bSJed Brown   ierr = DMDAVecGetArray(user->da,X,&p);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr);
148*c4762a1bSJed Brown   for (i=xs; i < xs+xm; i++) {
149*c4762a1bSJed Brown     for (j=ys; j < ys+ym; j++) {
150*c4762a1bSJed Brown       xi = coors[j][i].x; yi = coors[j][i].y;
151*c4762a1bSJed Brown       if (i == 0 || j == 0 || i == M-1 || j == N-1) p[j][i] = 0.0;
152*c4762a1bSJed Brown       else p[j][i] = (0.5/(PETSC_PI*sigmax*sigmay*PetscSqrtScalar(1.0-rho*rho)))*PetscExpScalar(-0.5/(1-rho*rho)*(PetscPowScalar((xi-mux)/sigmax,2) + PetscPowScalar((yi-muy)/sigmay,2) - 2*rho*(xi-mux)*(yi-muy)/(sigmax*sigmay)));
153*c4762a1bSJed Brown     }
154*c4762a1bSJed Brown   }
155*c4762a1bSJed Brown   /*  p[N/2+N%2][M/2+M%2] = 1/(user->dx*user->dy); */
156*c4762a1bSJed Brown 
157*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(cda,gc,&coors);CHKERRQ(ierr);
158*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(user->da,X,&p);CHKERRQ(ierr);
159*c4762a1bSJed Brown   PetscFunctionReturn(0);
160*c4762a1bSJed Brown }
161*c4762a1bSJed Brown 
162*c4762a1bSJed Brown /* First advection term */
163*c4762a1bSJed Brown PetscErrorCode adv1(PetscScalar **p,PetscScalar y,PetscInt i,PetscInt j,PetscInt M,PetscScalar *p1,AppCtx *user)
164*c4762a1bSJed Brown {
165*c4762a1bSJed Brown   PetscScalar f;
166*c4762a1bSJed Brown   /*  PetscScalar v1,v2,v3,v4,v5,s1,s2,s3; */
167*c4762a1bSJed Brown   PetscFunctionBegin;
168*c4762a1bSJed Brown   /*  if (i > 2 && i < M-2) {
169*c4762a1bSJed Brown     v1 = (y-user->ws)*(p[j][i-2] - p[j][i-3])/user->dx;
170*c4762a1bSJed Brown     v2 = (y-user->ws)*(p[j][i-1] - p[j][i-2])/user->dx;
171*c4762a1bSJed Brown     v3 = (y-user->ws)*(p[j][i] - p[j][i-1])/user->dx;
172*c4762a1bSJed Brown     v4 = (y-user->ws)*(p[j][i+1] - p[j][i])/user->dx;
173*c4762a1bSJed Brown     v5 = (y-user->ws)*(p[j][i+1] - p[j][i+2])/user->dx;
174*c4762a1bSJed Brown 
175*c4762a1bSJed Brown     s1 = v1/3.0 - (7.0/6.0)*v2 + (11.0/6.0)*v3;
176*c4762a1bSJed Brown     s2 =-v2/6.0 + (5.0/6.0)*v3 + (1.0/3.0)*v4;
177*c4762a1bSJed Brown     s3 = v3/3.0 + (5.0/6.0)*v4 - (1.0/6.0)*v5;
178*c4762a1bSJed Brown 
179*c4762a1bSJed Brown     *p1 = 0.1*s1 + 0.6*s2 + 0.3*s3;
180*c4762a1bSJed Brown     } else *p1 = 0.0; */
181*c4762a1bSJed Brown   f   =  (y - user->ws);
182*c4762a1bSJed Brown   *p1 = f*(p[j][i+1] - p[j][i-1])/(2*user->dx);
183*c4762a1bSJed Brown   PetscFunctionReturn(0);
184*c4762a1bSJed Brown }
185*c4762a1bSJed Brown 
186*c4762a1bSJed Brown /* Second advection term */
187*c4762a1bSJed Brown PetscErrorCode adv2(PetscScalar **p,PetscScalar x,PetscInt i,PetscInt j,PetscInt N,PetscScalar *p2,AppCtx *user)
188*c4762a1bSJed Brown {
189*c4762a1bSJed Brown   PetscScalar f;
190*c4762a1bSJed Brown   /*  PetscScalar v1,v2,v3,v4,v5,s1,s2,s3; */
191*c4762a1bSJed Brown   PetscFunctionBegin;
192*c4762a1bSJed Brown   /*  if (j > 2 && j < N-2) {
193*c4762a1bSJed Brown     v1 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j-2][i] - p[j-3][i])/user->dy;
194*c4762a1bSJed Brown     v2 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j-1][i] - p[j-2][i])/user->dy;
195*c4762a1bSJed Brown     v3 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j][i] - p[j-1][i])/user->dy;
196*c4762a1bSJed Brown     v4 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j+1][i] - p[j][i])/user->dy;
197*c4762a1bSJed Brown     v5 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j+2][i] - p[j+1][i])/user->dy;
198*c4762a1bSJed Brown 
199*c4762a1bSJed Brown     s1 = v1/3.0 - (7.0/6.0)*v2 + (11.0/6.0)*v3;
200*c4762a1bSJed Brown     s2 =-v2/6.0 + (5.0/6.0)*v3 + (1.0/3.0)*v4;
201*c4762a1bSJed Brown     s3 = v3/3.0 + (5.0/6.0)*v4 - (1.0/6.0)*v5;
202*c4762a1bSJed Brown 
203*c4762a1bSJed Brown     *p2 = 0.1*s1 + 0.6*s2 + 0.3*s3;
204*c4762a1bSJed Brown     } else *p2 = 0.0; */
205*c4762a1bSJed Brown   f   = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(x));
206*c4762a1bSJed Brown   *p2 = f*(p[j+1][i] - p[j-1][i])/(2*user->dy);
207*c4762a1bSJed Brown   PetscFunctionReturn(0);
208*c4762a1bSJed Brown }
209*c4762a1bSJed Brown 
210*c4762a1bSJed Brown /* Diffusion term */
211*c4762a1bSJed Brown PetscErrorCode diffuse(PetscScalar **p,PetscInt i,PetscInt j,PetscReal t,PetscScalar *p_diff,AppCtx * user)
212*c4762a1bSJed Brown {
213*c4762a1bSJed Brown   PetscFunctionBeginUser;
214*c4762a1bSJed Brown 
215*c4762a1bSJed Brown   *p_diff = user->disper_coe*((p[j-1][i] - 2*p[j][i] + p[j+1][i])/(user->dy*user->dy));
216*c4762a1bSJed Brown   PetscFunctionReturn(0);
217*c4762a1bSJed Brown }
218*c4762a1bSJed Brown 
219*c4762a1bSJed Brown PetscErrorCode BoundaryConditions(PetscScalar **p,DMDACoor2d **coors,PetscInt i,PetscInt j,PetscInt M, PetscInt N,PetscScalar **f,AppCtx *user)
220*c4762a1bSJed Brown {
221*c4762a1bSJed Brown   PetscScalar fwc,fthetac;
222*c4762a1bSJed Brown   PetscScalar w=coors[j][i].y,theta=coors[j][i].x;
223*c4762a1bSJed Brown 
224*c4762a1bSJed Brown   PetscFunctionBeginUser;
225*c4762a1bSJed Brown   if (user->bc == 0) { /* Natural boundary condition */
226*c4762a1bSJed Brown     f[j][i] = p[j][i];
227*c4762a1bSJed Brown   } else { /* Steady state boundary condition */
228*c4762a1bSJed Brown     fthetac = user->ws/(2*user->H)*(user->PM_min - user->Pmax*PetscSinScalar(theta));
229*c4762a1bSJed Brown     fwc = (w*w/2.0 - user->ws*w);
230*c4762a1bSJed Brown     if (i == 0 && j == 0) { /* left bottom corner */
231*c4762a1bSJed Brown       f[j][i] = fwc*(p[j][i+1] - p[j][i])/user->dx + fthetac*p[j][i] - user->disper_coe*(p[j+1][i] - p[j][i])/user->dy;
232*c4762a1bSJed Brown     } else if (i == 0 && j == N-1) { /* right bottom corner */
233*c4762a1bSJed Brown       f[j][i] = fwc*(p[j][i+1] - p[j][i])/user->dx + fthetac*p[j][i] - user->disper_coe*(p[j][i] - p[j-1][i])/user->dy;
234*c4762a1bSJed Brown     } else if (i == M-1 && j == 0) { /* left top corner */
235*c4762a1bSJed Brown       f[j][i] = fwc*(p[j][i] - p[j][i-1])/user->dx + fthetac*p[j][i] - user->disper_coe*(p[j+1][i] - p[j][i])/user->dy;
236*c4762a1bSJed Brown     } else if (i == M-1 && j == N-1) { /* right top corner */
237*c4762a1bSJed Brown       f[j][i] = fwc*(p[j][i] - p[j][i-1])/user->dx + fthetac*p[j][i] - user->disper_coe*(p[j][i] - p[j-1][i])/user->dy;
238*c4762a1bSJed Brown     } else if (i == 0) { /* Bottom edge */
239*c4762a1bSJed Brown       f[j][i] = fwc*(p[j][i+1] - p[j][i])/(user->dx) + fthetac*p[j][i] - user->disper_coe*(p[j+1][i] - p[j-1][i])/(2*user->dy);
240*c4762a1bSJed Brown     } else if (i == M-1) { /* Top edge */
241*c4762a1bSJed Brown       f[j][i] = fwc*(p[j][i] - p[j][i-1])/(user->dx) + fthetac*p[j][i] - user->disper_coe*(p[j+1][i] - p[j-1][i])/(2*user->dy);
242*c4762a1bSJed Brown     } else if (j == 0) { /* Left edge */
243*c4762a1bSJed Brown       f[j][i] = fwc*(p[j][i+1] - p[j][i-1])/(2*user->dx) + fthetac*p[j][i] - user->disper_coe*(p[j+1][i] - p[j][i])/(user->dy);
244*c4762a1bSJed Brown     } else if (j == N-1) { /* Right edge */
245*c4762a1bSJed Brown       f[j][i] = fwc*(p[j][i+1] - p[j][i-1])/(2*user->dx) + fthetac*p[j][i] - user->disper_coe*(p[j][i] - p[j-1][i])/(user->dy);
246*c4762a1bSJed Brown     }
247*c4762a1bSJed Brown   }
248*c4762a1bSJed Brown   PetscFunctionReturn(0);
249*c4762a1bSJed Brown }
250*c4762a1bSJed Brown 
251*c4762a1bSJed Brown PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
252*c4762a1bSJed Brown {
253*c4762a1bSJed Brown   PetscErrorCode ierr;
254*c4762a1bSJed Brown   AppCtx         *user=(AppCtx*)ctx;
255*c4762a1bSJed Brown   DM             cda;
256*c4762a1bSJed Brown   DMDACoor2d     **coors;
257*c4762a1bSJed Brown   PetscScalar    **p,**f,**pdot;
258*c4762a1bSJed Brown   PetscInt       i,j;
259*c4762a1bSJed Brown   PetscInt       xs,ys,xm,ym,M,N;
260*c4762a1bSJed Brown   Vec            localX,gc,localXdot;
261*c4762a1bSJed Brown   PetscScalar    p_adv1,p_adv2,p_diff;
262*c4762a1bSJed Brown 
263*c4762a1bSJed Brown   PetscFunctionBeginUser;
264*c4762a1bSJed Brown   ierr = DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
265*c4762a1bSJed Brown   ierr = DMGetCoordinateDM(user->da,&cda);CHKERRQ(ierr);
266*c4762a1bSJed Brown   ierr = DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr);
267*c4762a1bSJed Brown 
268*c4762a1bSJed Brown   ierr = DMGetLocalVector(user->da,&localX);CHKERRQ(ierr);
269*c4762a1bSJed Brown   ierr = DMGetLocalVector(user->da,&localXdot);CHKERRQ(ierr);
270*c4762a1bSJed Brown 
271*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
272*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
273*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(user->da,Xdot,INSERT_VALUES,localXdot);CHKERRQ(ierr);
274*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(user->da,Xdot,INSERT_VALUES,localXdot);CHKERRQ(ierr);
275*c4762a1bSJed Brown 
276*c4762a1bSJed Brown   ierr = DMGetCoordinatesLocal(user->da,&gc);CHKERRQ(ierr);
277*c4762a1bSJed Brown 
278*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(cda,gc,&coors);CHKERRQ(ierr);
279*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(user->da,localX,&p);CHKERRQ(ierr);
280*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(user->da,localXdot,&pdot);CHKERRQ(ierr);
281*c4762a1bSJed Brown   ierr = DMDAVecGetArray(user->da,F,&f);CHKERRQ(ierr);
282*c4762a1bSJed Brown 
283*c4762a1bSJed Brown   user->disper_coe = PetscPowScalar((user->lambda*user->ws)/(2*user->H),2)*user->q*(1.0-PetscExpScalar(-t/user->lambda));
284*c4762a1bSJed Brown   for (i=xs; i < xs+xm; i++) {
285*c4762a1bSJed Brown     for (j=ys; j < ys+ym; j++) {
286*c4762a1bSJed Brown       if (i == 0 || j == 0 || i == M-1 || j == N-1) {
287*c4762a1bSJed Brown         ierr = BoundaryConditions(p,coors,i,j,M,N,f,user);CHKERRQ(ierr);
288*c4762a1bSJed Brown       } else {
289*c4762a1bSJed Brown         ierr = adv1(p,coors[j][i].y,i,j,M,&p_adv1,user);CHKERRQ(ierr);
290*c4762a1bSJed Brown         ierr = adv2(p,coors[j][i].x,i,j,N,&p_adv2,user);CHKERRQ(ierr);
291*c4762a1bSJed Brown         ierr = diffuse(p,i,j,t,&p_diff,user);CHKERRQ(ierr);
292*c4762a1bSJed Brown         f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i];
293*c4762a1bSJed Brown       }
294*c4762a1bSJed Brown     }
295*c4762a1bSJed Brown   }
296*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(user->da,localX,&p);CHKERRQ(ierr);
297*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(user->da,localX,&pdot);CHKERRQ(ierr);
298*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(user->da,&localX);CHKERRQ(ierr);
299*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(user->da,&localXdot);CHKERRQ(ierr);
300*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(user->da,F,&f);CHKERRQ(ierr);
301*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(cda,gc,&coors);CHKERRQ(ierr);
302*c4762a1bSJed Brown 
303*c4762a1bSJed Brown   PetscFunctionReturn(0);
304*c4762a1bSJed Brown }
305*c4762a1bSJed Brown 
306*c4762a1bSJed Brown PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ctx)
307*c4762a1bSJed Brown {
308*c4762a1bSJed Brown   PetscErrorCode ierr;
309*c4762a1bSJed Brown   AppCtx         *user=(AppCtx*)ctx;
310*c4762a1bSJed Brown   DM             cda;
311*c4762a1bSJed Brown   DMDACoor2d     **coors;
312*c4762a1bSJed Brown   PetscInt       i,j;
313*c4762a1bSJed Brown   PetscInt       xs,ys,xm,ym,M,N;
314*c4762a1bSJed Brown   Vec            gc;
315*c4762a1bSJed Brown   PetscScalar    val[5],xi,yi;
316*c4762a1bSJed Brown   MatStencil     row,col[5];
317*c4762a1bSJed Brown   PetscScalar    c1,c3,c5;
318*c4762a1bSJed Brown 
319*c4762a1bSJed Brown   PetscFunctionBeginUser;
320*c4762a1bSJed Brown   ierr = DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
321*c4762a1bSJed Brown   ierr = DMGetCoordinateDM(user->da,&cda);CHKERRQ(ierr);
322*c4762a1bSJed Brown   ierr = DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr);
323*c4762a1bSJed Brown 
324*c4762a1bSJed Brown   ierr = DMGetCoordinatesLocal(user->da,&gc);CHKERRQ(ierr);
325*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(cda,gc,&coors);CHKERRQ(ierr);
326*c4762a1bSJed Brown   for (i=xs; i < xs+xm; i++) {
327*c4762a1bSJed Brown     for (j=ys; j < ys+ym; j++) {
328*c4762a1bSJed Brown       PetscInt nc = 0;
329*c4762a1bSJed Brown       xi = coors[j][i].x; yi = coors[j][i].y;
330*c4762a1bSJed Brown       row.i = i; row.j = j;
331*c4762a1bSJed Brown       if (i == 0 || j == 0 || i == M-1 || j == N-1) {
332*c4762a1bSJed Brown         if (user->bc == 0) {
333*c4762a1bSJed Brown           col[nc].i = i; col[nc].j = j; val[nc++] = 1.0;
334*c4762a1bSJed Brown         } else {
335*c4762a1bSJed Brown           PetscScalar fthetac,fwc;
336*c4762a1bSJed Brown           fthetac = user->ws/(2*user->H)*(user->PM_min - user->Pmax*PetscSinScalar(xi));
337*c4762a1bSJed Brown           fwc     = (yi*yi/2.0 - user->ws*yi);
338*c4762a1bSJed Brown           if (i==0 && j==0) {
339*c4762a1bSJed Brown             col[nc].i = i+1; col[nc].j = j;   val[nc++] = fwc/user->dx;
340*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j+1; val[nc++] = -user->disper_coe/user->dy;
341*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j;   val[nc++] = -fwc/user->dx + fthetac + user->disper_coe/user->dy;
342*c4762a1bSJed Brown           } else if (i==0 && j == N-1) {
343*c4762a1bSJed Brown             col[nc].i = i+1; col[nc].j = j;   val[nc++] = fwc/user->dx;
344*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j-1; val[nc++] = user->disper_coe/user->dy;
345*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j;   val[nc++] = -fwc/user->dx + fthetac - user->disper_coe/user->dy;
346*c4762a1bSJed Brown           } else if (i== M-1 && j == 0) {
347*c4762a1bSJed Brown             col[nc].i = i-1; col[nc].j = j;   val[nc++] = -fwc/user->dx;
348*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j+1; val[nc++] = -user->disper_coe/user->dy;
349*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j;   val[nc++] =  fwc/user->dx + fthetac + user->disper_coe/user->dy;
350*c4762a1bSJed Brown           } else if (i == M-1 && j == N-1) {
351*c4762a1bSJed Brown             col[nc].i = i-1; col[nc].j = j;   val[nc++] = -fwc/user->dx;
352*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j-1; val[nc++] =  user->disper_coe/user->dy;
353*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j;   val[nc++] =  fwc/user->dx + fthetac - user->disper_coe/user->dy;
354*c4762a1bSJed Brown           } else if (i==0) {
355*c4762a1bSJed Brown             col[nc].i = i+1; col[nc].j = j;   val[nc++] = fwc/user->dx;
356*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j+1; val[nc++] = -user->disper_coe/(2*user->dy);
357*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j-1; val[nc++] =  user->disper_coe/(2*user->dy);
358*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j;   val[nc++] = -fwc/user->dx + fthetac;
359*c4762a1bSJed Brown           } else if (i == M-1) {
360*c4762a1bSJed Brown             col[nc].i = i-1; col[nc].j = j;   val[nc++] = -fwc/user->dx;
361*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j+1; val[nc++] = -user->disper_coe/(2*user->dy);
362*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j-1; val[nc++] =  user->disper_coe/(2*user->dy);
363*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j;   val[nc++] = fwc/user->dx + fthetac;
364*c4762a1bSJed Brown           } else if (j==0) {
365*c4762a1bSJed Brown             col[nc].i = i+1; col[nc].j = j;   val[nc++] = fwc/(2*user->dx);
366*c4762a1bSJed Brown             col[nc].i = i-1; col[nc].j = j;   val[nc++] = -fwc/(2*user->dx);
367*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j+1; val[nc++] = -user->disper_coe/user->dy;
368*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j;   val[nc++] = user->disper_coe/user->dy + fthetac;
369*c4762a1bSJed Brown           } else if (j == N-1) {
370*c4762a1bSJed Brown             col[nc].i = i+1; col[nc].j = j;   val[nc++] = fwc/(2*user->dx);
371*c4762a1bSJed Brown             col[nc].i = i-1; col[nc].j = j;   val[nc++] = -fwc/(2*user->dx);
372*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j-1; val[nc++] = user->disper_coe/user->dy;
373*c4762a1bSJed Brown             col[nc].i = i;   col[nc].j = j;   val[nc++] = -user->disper_coe/user->dy + fthetac;
374*c4762a1bSJed Brown           }
375*c4762a1bSJed Brown         }
376*c4762a1bSJed Brown       } else {
377*c4762a1bSJed Brown         c1        = (yi-user->ws)/(2*user->dx);
378*c4762a1bSJed Brown         c3        = (user->ws/(2.0*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(xi))/(2*user->dy);
379*c4762a1bSJed Brown         c5        = (PetscPowScalar((user->lambda*user->ws)/(2*user->H),2)*user->q*(1.0-PetscExpScalar(-t/user->lambda)))/(user->dy*user->dy);
380*c4762a1bSJed Brown         col[nc].i = i-1; col[nc].j = j;   val[nc++] = c1;
381*c4762a1bSJed Brown         col[nc].i = i+1; col[nc].j = j;   val[nc++] = -c1;
382*c4762a1bSJed Brown         col[nc].i = i;   col[nc].j = j-1; val[nc++] = c3 + c5;
383*c4762a1bSJed Brown         col[nc].i = i;   col[nc].j = j+1; val[nc++] = -c3 + c5;
384*c4762a1bSJed Brown         col[nc].i = i;   col[nc].j = j;   val[nc++] = -2*c5 -a;
385*c4762a1bSJed Brown       }
386*c4762a1bSJed Brown       ierr = MatSetValuesStencil(Jpre,1,&row,nc,col,val,INSERT_VALUES);CHKERRQ(ierr);
387*c4762a1bSJed Brown     }
388*c4762a1bSJed Brown   }
389*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(cda,gc,&coors);CHKERRQ(ierr);
390*c4762a1bSJed Brown 
391*c4762a1bSJed Brown   ierr =  MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
392*c4762a1bSJed Brown   ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
393*c4762a1bSJed Brown   if (J != Jpre) {
394*c4762a1bSJed Brown     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
395*c4762a1bSJed Brown     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
396*c4762a1bSJed Brown   }
397*c4762a1bSJed Brown   PetscFunctionReturn(0);
398*c4762a1bSJed Brown }
399*c4762a1bSJed Brown 
400*c4762a1bSJed Brown PetscErrorCode Parameter_settings(AppCtx *user)
401*c4762a1bSJed Brown {
402*c4762a1bSJed Brown   PetscErrorCode ierr;
403*c4762a1bSJed Brown   PetscBool      flg;
404*c4762a1bSJed Brown 
405*c4762a1bSJed Brown   PetscFunctionBeginUser;
406*c4762a1bSJed Brown 
407*c4762a1bSJed Brown   /* Set default parameters */
408*c4762a1bSJed Brown   user->ws     = 1.0;
409*c4762a1bSJed Brown   user->H      = 5.0;  user->Pmax   = 2.1;
410*c4762a1bSJed Brown   user->PM_min = 1.0;  user->lambda = 0.1;
411*c4762a1bSJed Brown   user->q      = 1.0;  user->mux    = PetscAsinScalar(user->PM_min/user->Pmax);
412*c4762a1bSJed Brown   user->sigmax = 0.1;
413*c4762a1bSJed Brown   user->sigmay = 0.1;  user->rho  = 0.0;
414*c4762a1bSJed Brown   user->t0     = 0.0;  user->tmax = 2.0;
415*c4762a1bSJed Brown   user->xmin   = -1.0; user->xmax = 10.0;
416*c4762a1bSJed Brown   user->ymin   = -1.0; user->ymax = 10.0;
417*c4762a1bSJed Brown   user->bc     = 0;
418*c4762a1bSJed Brown 
419*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-ws",&user->ws,&flg);CHKERRQ(ierr);
420*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-Inertia",&user->H,&flg);CHKERRQ(ierr);
421*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-Pmax",&user->Pmax,&flg);CHKERRQ(ierr);
422*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-PM_min",&user->PM_min,&flg);CHKERRQ(ierr);
423*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-lambda",&user->lambda,&flg);CHKERRQ(ierr);
424*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-q",&user->q,&flg);CHKERRQ(ierr);
425*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-mux",&user->mux,&flg);CHKERRQ(ierr);
426*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-sigmax",&user->sigmax,&flg);CHKERRQ(ierr);
427*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-muy",&user->muy,&flg);CHKERRQ(ierr);
428*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-sigmay",&user->sigmay,&flg);CHKERRQ(ierr);
429*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-rho",&user->rho,&flg);CHKERRQ(ierr);
430*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-t0",&user->t0,&flg);CHKERRQ(ierr);
431*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-tmax",&user->tmax,&flg);CHKERRQ(ierr);
432*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-xmin",&user->xmin,&flg);CHKERRQ(ierr);
433*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-xmax",&user->xmax,&flg);CHKERRQ(ierr);
434*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-ymin",&user->ymin,&flg);CHKERRQ(ierr);
435*c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-ymax",&user->ymax,&flg);CHKERRQ(ierr);
436*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-bc_type",&user->bc,&flg);CHKERRQ(ierr);
437*c4762a1bSJed Brown   user->muy = user->ws;
438*c4762a1bSJed Brown   PetscFunctionReturn(0);
439*c4762a1bSJed Brown }
440*c4762a1bSJed Brown 
441*c4762a1bSJed Brown 
442*c4762a1bSJed Brown /*TEST
443*c4762a1bSJed Brown 
444*c4762a1bSJed Brown    build:
445*c4762a1bSJed Brown       requires: !complex
446*c4762a1bSJed Brown 
447*c4762a1bSJed Brown    test:
448*c4762a1bSJed Brown       args: -nox -ts_max_steps 2
449*c4762a1bSJed Brown       localrunfiles: petscopt_ex6
450*c4762a1bSJed Brown       timeoutfactor: 4
451*c4762a1bSJed Brown 
452*c4762a1bSJed Brown TEST*/
453