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 6*c4762a1bSJed Brown Boundary conditions: 7*c4762a1bSJed Brown Zero dirichlet in y using ghosted values 8*c4762a1bSJed Brown Periodic in x 9*c4762a1bSJed Brown 10*c4762a1bSJed 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. 11*c4762a1bSJed Brown x_t = (y - ws) 12*c4762a1bSJed Brown y_t = (ws/2H)*(Pm - Pmax*sin(x) - D*(w - ws)) 13*c4762a1bSJed Brown 14*c4762a1bSJed Brown In this example, we can see the effect of a fault, that zeroes the electrical power output 15*c4762a1bSJed Brown Pmax*sin(x), on the PDF. The fault on/off times can be controlled by options -tf and -tcl respectively. 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown */ 18*c4762a1bSJed Brown 19*c4762a1bSJed Brown #include <petscdm.h> 20*c4762a1bSJed Brown #include <petscdmda.h> 21*c4762a1bSJed Brown #include <petscts.h> 22*c4762a1bSJed Brown 23*c4762a1bSJed Brown /* 24*c4762a1bSJed Brown User-defined data structures and routines 25*c4762a1bSJed Brown */ 26*c4762a1bSJed Brown typedef struct { 27*c4762a1bSJed Brown PetscScalar ws; /* Synchronous speed */ 28*c4762a1bSJed Brown PetscScalar H; /* Inertia constant */ 29*c4762a1bSJed Brown PetscScalar D; /* Damping constant */ 30*c4762a1bSJed Brown PetscScalar Pmax,Pmax_s; /* Maximum power output of generator */ 31*c4762a1bSJed Brown PetscScalar PM_min; /* Mean mechanical power input */ 32*c4762a1bSJed Brown PetscScalar lambda; /* correlation time */ 33*c4762a1bSJed Brown PetscScalar q; /* noise strength */ 34*c4762a1bSJed Brown PetscScalar mux; /* Initial average angle */ 35*c4762a1bSJed Brown PetscScalar sigmax; /* Standard deviation of initial angle */ 36*c4762a1bSJed Brown PetscScalar muy; /* Average speed */ 37*c4762a1bSJed Brown PetscScalar sigmay; /* standard deviation of initial speed */ 38*c4762a1bSJed Brown PetscScalar rho; /* Cross-correlation coefficient */ 39*c4762a1bSJed Brown PetscScalar xmin; /* left boundary of angle */ 40*c4762a1bSJed Brown PetscScalar xmax; /* right boundary of angle */ 41*c4762a1bSJed Brown PetscScalar ymin; /* bottom boundary of speed */ 42*c4762a1bSJed Brown PetscScalar ymax; /* top boundary of speed */ 43*c4762a1bSJed Brown PetscScalar dx; /* x step size */ 44*c4762a1bSJed Brown PetscScalar dy; /* y step size */ 45*c4762a1bSJed Brown PetscScalar disper_coe; /* Dispersion coefficient */ 46*c4762a1bSJed Brown DM da; 47*c4762a1bSJed Brown PetscInt st_width; /* Stencil width */ 48*c4762a1bSJed Brown DMBoundaryType bx; /* x boundary type */ 49*c4762a1bSJed Brown DMBoundaryType by; /* y boundary type */ 50*c4762a1bSJed Brown PetscReal tf,tcl; /* Fault incidence and clearing times */ 51*c4762a1bSJed Brown } AppCtx; 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown PetscErrorCode Parameter_settings(AppCtx*); 54*c4762a1bSJed Brown PetscErrorCode ini_bou(Vec,AppCtx*); 55*c4762a1bSJed Brown PetscErrorCode IFunction(TS,PetscReal,Vec,Vec,Vec,void*); 56*c4762a1bSJed Brown PetscErrorCode IJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 57*c4762a1bSJed Brown PetscErrorCode PostStep(TS); 58*c4762a1bSJed Brown 59*c4762a1bSJed Brown int main(int argc, char **argv) 60*c4762a1bSJed Brown { 61*c4762a1bSJed Brown PetscErrorCode ierr; 62*c4762a1bSJed Brown Vec x; /* Solution vector */ 63*c4762a1bSJed Brown TS ts; /* Time-stepping context */ 64*c4762a1bSJed Brown AppCtx user; /* Application context */ 65*c4762a1bSJed Brown PetscViewer viewer; 66*c4762a1bSJed Brown 67*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,"petscopt_ex8", help);if (ierr) return ierr; 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown /* Get physics and time parameters */ 70*c4762a1bSJed Brown ierr = Parameter_settings(&user);CHKERRQ(ierr); 71*c4762a1bSJed Brown /* Create a 2D DA with dof = 1 */ 72*c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 73*c4762a1bSJed Brown ierr = DMSetFromOptions(user.da);CHKERRQ(ierr); 74*c4762a1bSJed Brown ierr = DMSetUp(user.da);CHKERRQ(ierr); 75*c4762a1bSJed Brown /* Set x and y coordinates */ 76*c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(user.da,user.xmin,user.xmax,user.ymin,user.ymax,0,0);CHKERRQ(ierr); 77*c4762a1bSJed Brown ierr = DMDASetCoordinateName(user.da,0,"X - the angle");CHKERRQ(ierr); 78*c4762a1bSJed Brown ierr = DMDASetCoordinateName(user.da,1,"Y - the speed");CHKERRQ(ierr); 79*c4762a1bSJed Brown 80*c4762a1bSJed Brown /* Get global vector x from DM */ 81*c4762a1bSJed Brown ierr = DMCreateGlobalVector(user.da,&x);CHKERRQ(ierr); 82*c4762a1bSJed Brown 83*c4762a1bSJed Brown ierr = ini_bou(x,&user);CHKERRQ(ierr); 84*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ini_x",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); 85*c4762a1bSJed Brown ierr = VecView(x,viewer);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 87*c4762a1bSJed Brown 88*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 89*c4762a1bSJed Brown ierr = TSSetDM(ts,user.da);CHKERRQ(ierr); 90*c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 91*c4762a1bSJed Brown ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr); 92*c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,IFunction,&user);CHKERRQ(ierr); 93*c4762a1bSJed Brown /* ierr = TSSetIJacobian(ts,NULL,NULL,IJacobian,&user);CHKERRQ(ierr); */ 94*c4762a1bSJed Brown ierr = TSSetApplicationContext(ts,&user);CHKERRQ(ierr); 95*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,.005);CHKERRQ(ierr); 96*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = TSSetPostStep(ts,PostStep);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = TSSolve(ts,x);CHKERRQ(ierr); 99*c4762a1bSJed Brown 100*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"fin_x",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); 101*c4762a1bSJed Brown ierr = VecView(x,viewer);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 103*c4762a1bSJed Brown 104*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 105*c4762a1bSJed Brown ierr = DMDestroy(&user.da);CHKERRQ(ierr); 106*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 107*c4762a1bSJed Brown ierr = PetscFinalize(); 108*c4762a1bSJed Brown return ierr; 109*c4762a1bSJed Brown } 110*c4762a1bSJed Brown 111*c4762a1bSJed Brown PetscErrorCode PostStep(TS ts) 112*c4762a1bSJed Brown { 113*c4762a1bSJed Brown PetscErrorCode ierr; 114*c4762a1bSJed Brown Vec X; 115*c4762a1bSJed Brown AppCtx *user; 116*c4762a1bSJed Brown PetscReal t; 117*c4762a1bSJed Brown PetscScalar asum; 118*c4762a1bSJed Brown 119*c4762a1bSJed Brown PetscFunctionBegin; 120*c4762a1bSJed Brown ierr = TSGetApplicationContext(ts,&user);CHKERRQ(ierr); 121*c4762a1bSJed Brown ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 122*c4762a1bSJed Brown ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 123*c4762a1bSJed Brown /* 124*c4762a1bSJed Brown if (t >= .2) { 125*c4762a1bSJed Brown ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 126*c4762a1bSJed Brown ierr = VecView(X,PETSC_VIEWER_BINARY_WORLD);CHKERRQ(ierr); 127*c4762a1bSJed Brown exit(0); 128*c4762a1bSJed Brown results in initial conditions after fault in binaryoutput 129*c4762a1bSJed Brown }*/ 130*c4762a1bSJed Brown 131*c4762a1bSJed 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 */ 132*c4762a1bSJed Brown else user->Pmax = user->Pmax_s; 133*c4762a1bSJed Brown 134*c4762a1bSJed Brown ierr = VecSum(X,&asum);CHKERRQ(ierr); 135*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"sum(p) at t = %f = %f\n",(double)t,(double)(asum));CHKERRQ(ierr); 136*c4762a1bSJed Brown PetscFunctionReturn(0); 137*c4762a1bSJed Brown } 138*c4762a1bSJed Brown 139*c4762a1bSJed Brown PetscErrorCode ini_bou(Vec X,AppCtx* user) 140*c4762a1bSJed Brown { 141*c4762a1bSJed Brown PetscErrorCode ierr; 142*c4762a1bSJed Brown DM cda; 143*c4762a1bSJed Brown DMDACoor2d **coors; 144*c4762a1bSJed Brown PetscScalar **p; 145*c4762a1bSJed Brown Vec gc; 146*c4762a1bSJed Brown PetscInt M,N,Ir,J; 147*c4762a1bSJed Brown PetscMPIInt rank; 148*c4762a1bSJed Brown 149*c4762a1bSJed Brown PetscFunctionBeginUser; 150*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 151*c4762a1bSJed Brown ierr = DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 152*c4762a1bSJed Brown user->dx = (user->xmax - user->xmin)/(M-1); user->dy = (user->ymax - user->ymin)/(N-1); 153*c4762a1bSJed Brown ierr = DMGetCoordinateDM(user->da,&cda);CHKERRQ(ierr); 154*c4762a1bSJed Brown ierr = DMGetCoordinates(user->da,&gc);CHKERRQ(ierr); 155*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(cda,gc,&coors);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = DMDAVecGetArray(user->da,X,&p);CHKERRQ(ierr); 157*c4762a1bSJed Brown 158*c4762a1bSJed Brown /* Point mass at (mux,muy) */ 159*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Original user->mux = %f, user->muy = %f\n",user->mux,user->muy);CHKERRQ(ierr); 160*c4762a1bSJed Brown ierr = DMDAGetLogicalCoordinate(user->da,user->mux,user->muy,0.0,&Ir,&J,NULL,&user->mux,&user->muy,NULL);CHKERRQ(ierr); 161*c4762a1bSJed Brown user->PM_min = user->Pmax*PetscSinScalar(user->mux); 162*c4762a1bSJed Brown ierr = 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);CHKERRQ(ierr); 163*c4762a1bSJed Brown if (Ir > -1 && J > -1) { 164*c4762a1bSJed Brown p[J][Ir] = 1.0; 165*c4762a1bSJed Brown } 166*c4762a1bSJed Brown 167*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(cda,gc,&coors);CHKERRQ(ierr); 168*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(user->da,X,&p);CHKERRQ(ierr); 169*c4762a1bSJed Brown PetscFunctionReturn(0); 170*c4762a1bSJed Brown } 171*c4762a1bSJed Brown 172*c4762a1bSJed Brown /* First advection term */ 173*c4762a1bSJed Brown PetscErrorCode adv1(PetscScalar **p,PetscScalar y,PetscInt i,PetscInt j,PetscInt M,PetscScalar *p1,AppCtx *user) 174*c4762a1bSJed Brown { 175*c4762a1bSJed Brown PetscScalar f,fpos,fneg; 176*c4762a1bSJed Brown PetscFunctionBegin; 177*c4762a1bSJed Brown f = (y - user->ws); 178*c4762a1bSJed Brown fpos = PetscMax(f,0); 179*c4762a1bSJed Brown fneg = PetscMin(f,0); 180*c4762a1bSJed Brown if (user->st_width == 1) { 181*c4762a1bSJed Brown *p1 = fpos*(p[j][i] - p[j][i-1])/user->dx + fneg*(p[j][i+1] - p[j][i])/user->dx; 182*c4762a1bSJed Brown } else if (user->st_width == 2) { 183*c4762a1bSJed 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); 184*c4762a1bSJed Brown } else if (user->st_width == 3) { 185*c4762a1bSJed 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); 186*c4762a1bSJed Brown } 187*c4762a1bSJed Brown /* *p1 = f*(p[j][i+1] - p[j][i-1])/user->dx;*/ 188*c4762a1bSJed Brown PetscFunctionReturn(0); 189*c4762a1bSJed Brown } 190*c4762a1bSJed Brown 191*c4762a1bSJed Brown /* Second advection term */ 192*c4762a1bSJed Brown PetscErrorCode adv2(PetscScalar **p,PetscScalar x,PetscScalar y,PetscInt i,PetscInt j,PetscInt N,PetscScalar *p2,AppCtx *user) 193*c4762a1bSJed Brown { 194*c4762a1bSJed Brown PetscScalar f,fpos,fneg; 195*c4762a1bSJed Brown PetscFunctionBegin; 196*c4762a1bSJed Brown f = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(x) - user->D*(y - user->ws)); 197*c4762a1bSJed Brown fpos = PetscMax(f,0); 198*c4762a1bSJed Brown fneg = PetscMin(f,0); 199*c4762a1bSJed Brown if (user->st_width == 1) { 200*c4762a1bSJed Brown *p2 = fpos*(p[j][i] - p[j-1][i])/user->dy + fneg*(p[j+1][i] - p[j][i])/user->dy; 201*c4762a1bSJed Brown } else if (user->st_width ==2) { 202*c4762a1bSJed 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); 203*c4762a1bSJed Brown } else if (user->st_width == 3) { 204*c4762a1bSJed 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); 205*c4762a1bSJed Brown } 206*c4762a1bSJed Brown 207*c4762a1bSJed Brown /* *p2 = f*(p[j+1][i] - p[j-1][i])/user->dy;*/ 208*c4762a1bSJed Brown PetscFunctionReturn(0); 209*c4762a1bSJed Brown } 210*c4762a1bSJed Brown 211*c4762a1bSJed Brown /* Diffusion term */ 212*c4762a1bSJed Brown PetscErrorCode diffuse(PetscScalar **p,PetscInt i,PetscInt j,PetscReal t,PetscScalar *p_diff,AppCtx * user) 213*c4762a1bSJed Brown { 214*c4762a1bSJed Brown PetscFunctionBeginUser; 215*c4762a1bSJed Brown if (user->st_width == 1) { 216*c4762a1bSJed Brown *p_diff = user->disper_coe*((p[j-1][i] - 2*p[j][i] + p[j+1][i])/(user->dy*user->dy)); 217*c4762a1bSJed Brown } else if (user->st_width == 2) { 218*c4762a1bSJed 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)); 219*c4762a1bSJed Brown } else if (user->st_width == 3) { 220*c4762a1bSJed 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)); 221*c4762a1bSJed Brown } 222*c4762a1bSJed Brown PetscFunctionReturn(0); 223*c4762a1bSJed Brown } 224*c4762a1bSJed Brown 225*c4762a1bSJed Brown PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 226*c4762a1bSJed Brown { 227*c4762a1bSJed Brown PetscErrorCode ierr; 228*c4762a1bSJed Brown AppCtx *user=(AppCtx*)ctx; 229*c4762a1bSJed Brown DM cda; 230*c4762a1bSJed Brown DMDACoor2d **coors; 231*c4762a1bSJed Brown PetscScalar **p,**f,**pdot; 232*c4762a1bSJed Brown PetscInt i,j; 233*c4762a1bSJed Brown PetscInt xs,ys,xm,ym,M,N; 234*c4762a1bSJed Brown Vec localX,gc,localXdot; 235*c4762a1bSJed Brown PetscScalar p_adv1 = 0.0,p_adv2 = 0.0,p_diff = 0; 236*c4762a1bSJed Brown PetscScalar diffuse1,gamma; 237*c4762a1bSJed Brown 238*c4762a1bSJed Brown PetscFunctionBeginUser; 239*c4762a1bSJed Brown ierr = DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 240*c4762a1bSJed Brown ierr = DMGetCoordinateDM(user->da,&cda);CHKERRQ(ierr); 241*c4762a1bSJed Brown ierr = DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr); 242*c4762a1bSJed Brown 243*c4762a1bSJed Brown ierr = DMGetLocalVector(user->da,&localX);CHKERRQ(ierr); 244*c4762a1bSJed Brown ierr = DMGetLocalVector(user->da,&localXdot);CHKERRQ(ierr); 245*c4762a1bSJed Brown 246*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 247*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 248*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(user->da,Xdot,INSERT_VALUES,localXdot);CHKERRQ(ierr); 249*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(user->da,Xdot,INSERT_VALUES,localXdot);CHKERRQ(ierr); 250*c4762a1bSJed Brown 251*c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(user->da,&gc);CHKERRQ(ierr); 252*c4762a1bSJed Brown 253*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(cda,gc,&coors);CHKERRQ(ierr); 254*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(user->da,localX,&p);CHKERRQ(ierr); 255*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(user->da,localXdot,&pdot);CHKERRQ(ierr); 256*c4762a1bSJed Brown ierr = DMDAVecGetArray(user->da,F,&f);CHKERRQ(ierr); 257*c4762a1bSJed Brown 258*c4762a1bSJed Brown gamma = user->D*user->ws/(2*user->H); 259*c4762a1bSJed Brown diffuse1 = user->lambda*user->lambda*user->q/(user->lambda*gamma+1)*(1.0 - PetscExpScalar(-t*(gamma+1.0)/user->lambda)); 260*c4762a1bSJed Brown user->disper_coe = user->ws*user->ws/(4*user->H*user->H)*diffuse1; 261*c4762a1bSJed Brown 262*c4762a1bSJed Brown for (i=xs; i < xs+xm; i++) { 263*c4762a1bSJed Brown for (j=ys; j < ys+ym; j++) { 264*c4762a1bSJed Brown ierr = adv1(p,coors[j][i].y,i,j,M,&p_adv1,user);CHKERRQ(ierr); 265*c4762a1bSJed Brown ierr = adv2(p,coors[j][i].x,coors[j][i].y,i,j,N,&p_adv2,user);CHKERRQ(ierr); 266*c4762a1bSJed Brown ierr = diffuse(p,i,j,t,&p_diff,user);CHKERRQ(ierr); 267*c4762a1bSJed Brown f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i]; 268*c4762a1bSJed Brown } 269*c4762a1bSJed Brown } 270*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(user->da,localX,&p);CHKERRQ(ierr); 271*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(user->da,localX,&pdot);CHKERRQ(ierr); 272*c4762a1bSJed Brown ierr = DMRestoreLocalVector(user->da,&localX);CHKERRQ(ierr); 273*c4762a1bSJed Brown ierr = DMRestoreLocalVector(user->da,&localXdot);CHKERRQ(ierr); 274*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(user->da,F,&f);CHKERRQ(ierr); 275*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(cda,gc,&coors);CHKERRQ(ierr); 276*c4762a1bSJed Brown 277*c4762a1bSJed Brown PetscFunctionReturn(0); 278*c4762a1bSJed Brown } 279*c4762a1bSJed Brown 280*c4762a1bSJed Brown PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ctx) 281*c4762a1bSJed Brown { 282*c4762a1bSJed Brown PetscErrorCode ierr; 283*c4762a1bSJed Brown AppCtx *user=(AppCtx*)ctx; 284*c4762a1bSJed Brown DM cda; 285*c4762a1bSJed Brown DMDACoor2d **coors; 286*c4762a1bSJed Brown PetscInt i,j; 287*c4762a1bSJed Brown PetscInt xs,ys,xm,ym,M,N; 288*c4762a1bSJed Brown Vec gc; 289*c4762a1bSJed Brown PetscScalar val[5],xi,yi; 290*c4762a1bSJed Brown MatStencil row,col[5]; 291*c4762a1bSJed Brown PetscScalar c1,c3,c5,c1pos,c1neg,c3pos,c3neg; 292*c4762a1bSJed Brown 293*c4762a1bSJed Brown PetscFunctionBeginUser; 294*c4762a1bSJed Brown ierr = DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 295*c4762a1bSJed Brown ierr = DMGetCoordinateDM(user->da,&cda);CHKERRQ(ierr); 296*c4762a1bSJed Brown ierr = DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr); 297*c4762a1bSJed Brown 298*c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(user->da,&gc);CHKERRQ(ierr); 299*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(cda,gc,&coors);CHKERRQ(ierr); 300*c4762a1bSJed Brown for (i=xs; i < xs+xm; i++) { 301*c4762a1bSJed Brown for (j=ys; j < ys+ym; j++) { 302*c4762a1bSJed Brown PetscInt nc = 0; 303*c4762a1bSJed Brown xi = coors[j][i].x; yi = coors[j][i].y; 304*c4762a1bSJed Brown row.i = i; row.j = j; 305*c4762a1bSJed Brown c1 = (yi-user->ws)/user->dx; 306*c4762a1bSJed Brown c1pos = PetscMax(c1,0); 307*c4762a1bSJed Brown c1neg = PetscMin(c1,0); 308*c4762a1bSJed Brown c3 = (user->ws/(2.0*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(xi) - user->D*(yi - user->ws))/user->dy; 309*c4762a1bSJed Brown c3pos = PetscMax(c3,0); 310*c4762a1bSJed Brown c3neg = PetscMin(c3,0); 311*c4762a1bSJed Brown c5 = (PetscPowScalar((user->lambda*user->ws)/(2*user->H),2)*user->q*(1.0-PetscExpScalar(-t/user->lambda)))/(user->dy*user->dy); 312*c4762a1bSJed Brown col[nc].i = i-1; col[nc].j = j; val[nc++] = c1pos; 313*c4762a1bSJed Brown col[nc].i = i+1; col[nc].j = j; val[nc++] = -c1neg; 314*c4762a1bSJed Brown col[nc].i = i; col[nc].j = j-1; val[nc++] = c3pos + c5; 315*c4762a1bSJed Brown col[nc].i = i; col[nc].j = j+1; val[nc++] = -c3neg + c5; 316*c4762a1bSJed Brown col[nc].i = i; col[nc].j = j; val[nc++] = -c1pos + c1neg -c3pos + c3neg -2*c5 -a; 317*c4762a1bSJed Brown ierr = MatSetValuesStencil(Jpre,1,&row,nc,col,val,INSERT_VALUES);CHKERRQ(ierr); 318*c4762a1bSJed Brown } 319*c4762a1bSJed Brown } 320*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(cda,gc,&coors);CHKERRQ(ierr); 321*c4762a1bSJed Brown 322*c4762a1bSJed Brown ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 323*c4762a1bSJed Brown ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 324*c4762a1bSJed Brown if (J != Jpre) { 325*c4762a1bSJed Brown ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 326*c4762a1bSJed Brown ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 327*c4762a1bSJed Brown } 328*c4762a1bSJed Brown PetscFunctionReturn(0); 329*c4762a1bSJed Brown } 330*c4762a1bSJed Brown 331*c4762a1bSJed Brown PetscErrorCode Parameter_settings(AppCtx *user) 332*c4762a1bSJed Brown { 333*c4762a1bSJed Brown PetscErrorCode ierr; 334*c4762a1bSJed Brown PetscBool flg; 335*c4762a1bSJed Brown 336*c4762a1bSJed Brown PetscFunctionBeginUser; 337*c4762a1bSJed Brown 338*c4762a1bSJed Brown /* Set default parameters */ 339*c4762a1bSJed Brown user->ws = 1.0; 340*c4762a1bSJed Brown user->H = 5.0; 341*c4762a1bSJed Brown user->D = 0.0; 342*c4762a1bSJed Brown user->Pmax = user->Pmax_s = 2.1; 343*c4762a1bSJed Brown user->PM_min = 1.0; 344*c4762a1bSJed Brown user->lambda = 0.1; 345*c4762a1bSJed Brown user->q = 1.0; 346*c4762a1bSJed Brown user->mux = PetscAsinScalar(user->PM_min/user->Pmax); 347*c4762a1bSJed Brown user->sigmax = 0.1; 348*c4762a1bSJed Brown user->sigmay = 0.1; 349*c4762a1bSJed Brown user->rho = 0.0; 350*c4762a1bSJed Brown user->xmin = -PETSC_PI; 351*c4762a1bSJed Brown user->xmax = PETSC_PI; 352*c4762a1bSJed Brown user->bx = DM_BOUNDARY_PERIODIC; 353*c4762a1bSJed Brown user->by = DM_BOUNDARY_GHOSTED; 354*c4762a1bSJed Brown user->tf = user->tcl = -1; 355*c4762a1bSJed Brown user->ymin = -2.0; 356*c4762a1bSJed Brown user->ymax = 2.0; 357*c4762a1bSJed Brown user->st_width = 1; 358*c4762a1bSJed Brown 359*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-ws",&user->ws,&flg);CHKERRQ(ierr); 360*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-Inertia",&user->H,&flg);CHKERRQ(ierr); 361*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-D",&user->D,&flg);CHKERRQ(ierr); 362*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-Pmax",&user->Pmax,&flg);CHKERRQ(ierr); 363*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-PM_min",&user->PM_min,&flg);CHKERRQ(ierr); 364*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-lambda",&user->lambda,&flg);CHKERRQ(ierr); 365*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-q",&user->q,&flg);CHKERRQ(ierr); 366*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-mux",&user->mux,&flg);CHKERRQ(ierr); 367*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-muy",&user->muy,&flg);CHKERRQ(ierr); 368*c4762a1bSJed Brown if (flg == 0) { 369*c4762a1bSJed Brown user->muy = user->ws; 370*c4762a1bSJed Brown } 371*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-xmin",&user->xmin,&flg);CHKERRQ(ierr); 372*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-xmax",&user->xmax,&flg);CHKERRQ(ierr); 373*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-ymin",&user->ymin,&flg);CHKERRQ(ierr); 374*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-ymax",&user->ymax,&flg);CHKERRQ(ierr); 375*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-stencil_width",&user->st_width,&flg);CHKERRQ(ierr); 376*c4762a1bSJed Brown ierr = PetscOptionsGetEnum(NULL,NULL,"-bx",DMBoundaryTypes,(PetscEnum*)&user->bx,&flg);CHKERRQ(ierr); 377*c4762a1bSJed Brown ierr = PetscOptionsGetEnum(NULL,NULL,"-by",DMBoundaryTypes,(PetscEnum*)&user->by,&flg);CHKERRQ(ierr); 378*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-tf",&user->tf,&flg);CHKERRQ(ierr); 379*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-tcl",&user->tcl,&flg);CHKERRQ(ierr); 380*c4762a1bSJed Brown PetscFunctionReturn(0); 381*c4762a1bSJed Brown } 382*c4762a1bSJed Brown 383*c4762a1bSJed Brown /*TEST 384*c4762a1bSJed Brown 385*c4762a1bSJed Brown build: 386*c4762a1bSJed Brown requires: !complex x 387*c4762a1bSJed Brown 388*c4762a1bSJed Brown test: 389*c4762a1bSJed Brown args: -ts_max_steps 1 390*c4762a1bSJed Brown localrunfiles: petscopt_ex8 391*c4762a1bSJed Brown timeoutfactor: 3 392*c4762a1bSJed Brown 393*c4762a1bSJed Brown TEST*/ 394