1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Time-dependent PDE in 2d for calculating joint PDF. \n"; 3*c4762a1bSJed Brown /* 4*c4762a1bSJed Brown p_t = -x_t*p_x -y_t*p_y + f(t)*p_yy 5*c4762a1bSJed Brown xmin < x < xmax, ymin < y < ymax; 6*c4762a1bSJed Brown 7*c4762a1bSJed Brown Boundary conditions Neumman using mirror values 8*c4762a1bSJed Brown 9*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. 10*c4762a1bSJed Brown x_t = (y - ws) y_t = (ws/2H)*(Pm - Pmax*sin(x)) 11*c4762a1bSJed Brown 12*c4762a1bSJed Brown */ 13*c4762a1bSJed Brown 14*c4762a1bSJed Brown #include <petscdm.h> 15*c4762a1bSJed Brown #include <petscdmda.h> 16*c4762a1bSJed Brown #include <petscts.h> 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown /* 19*c4762a1bSJed Brown User-defined data structures and routines 20*c4762a1bSJed Brown */ 21*c4762a1bSJed Brown typedef struct { 22*c4762a1bSJed Brown PetscScalar ws; /* Synchronous speed */ 23*c4762a1bSJed Brown PetscScalar H; /* Inertia constant */ 24*c4762a1bSJed Brown PetscScalar D; /* Damping constant */ 25*c4762a1bSJed Brown PetscScalar Pmax; /* Maximum power output of generator */ 26*c4762a1bSJed Brown PetscScalar PM_min; /* Mean mechanical power input */ 27*c4762a1bSJed Brown PetscScalar lambda; /* correlation time */ 28*c4762a1bSJed Brown PetscScalar q; /* noise strength */ 29*c4762a1bSJed Brown PetscScalar mux; /* Initial average angle */ 30*c4762a1bSJed Brown PetscScalar sigmax; /* Standard deviation of initial angle */ 31*c4762a1bSJed Brown PetscScalar muy; /* Average speed */ 32*c4762a1bSJed Brown PetscScalar sigmay; /* standard deviation of initial speed */ 33*c4762a1bSJed Brown PetscScalar rho; /* Cross-correlation coefficient */ 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 PetscScalar disper_coe; /* Dispersion coefficient */ 41*c4762a1bSJed Brown DM da; 42*c4762a1bSJed Brown PetscInt st_width; /* Stencil width */ 43*c4762a1bSJed Brown DMBoundaryType bx; /* x boundary type */ 44*c4762a1bSJed Brown DMBoundaryType by; /* y boundary type */ 45*c4762a1bSJed Brown PetscBool nonoiseinitial; 46*c4762a1bSJed Brown } AppCtx; 47*c4762a1bSJed Brown 48*c4762a1bSJed Brown PetscErrorCode Parameter_settings(AppCtx*); 49*c4762a1bSJed Brown PetscErrorCode ini_bou(Vec,AppCtx*); 50*c4762a1bSJed Brown PetscErrorCode IFunction(TS,PetscReal,Vec,Vec,Vec,void*); 51*c4762a1bSJed Brown PetscErrorCode IJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 52*c4762a1bSJed Brown PetscErrorCode PostStep(TS); 53*c4762a1bSJed Brown 54*c4762a1bSJed Brown int main(int argc, char **argv) 55*c4762a1bSJed Brown { 56*c4762a1bSJed Brown PetscErrorCode ierr; 57*c4762a1bSJed Brown Vec x; /* Solution vector */ 58*c4762a1bSJed Brown TS ts; /* Time-stepping context */ 59*c4762a1bSJed Brown AppCtx user; /* Application context */ 60*c4762a1bSJed Brown PetscViewer viewer; 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,"petscopt_ex7", help);if (ierr) return ierr; 63*c4762a1bSJed Brown 64*c4762a1bSJed Brown /* Get physics and time parameters */ 65*c4762a1bSJed Brown ierr = Parameter_settings(&user);CHKERRQ(ierr); 66*c4762a1bSJed Brown /* Create a 2D DA with dof = 1 */ 67*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); 68*c4762a1bSJed Brown ierr = DMSetFromOptions(user.da);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = DMSetUp(user.da);CHKERRQ(ierr); 70*c4762a1bSJed Brown /* Set x and y coordinates */ 71*c4762a1bSJed Brown ierr = DMDASetUniformCoordinates(user.da,user.xmin,user.xmax,user.ymin,user.ymax,0,0);CHKERRQ(ierr); 72*c4762a1bSJed Brown ierr = DMDASetCoordinateName(user.da,0,"X - the angle");CHKERRQ(ierr); 73*c4762a1bSJed Brown ierr = DMDASetCoordinateName(user.da,1,"Y - the speed");CHKERRQ(ierr); 74*c4762a1bSJed Brown 75*c4762a1bSJed Brown /* Get global vector x from DM */ 76*c4762a1bSJed Brown ierr = DMCreateGlobalVector(user.da,&x);CHKERRQ(ierr); 77*c4762a1bSJed Brown 78*c4762a1bSJed Brown ierr = ini_bou(x,&user);CHKERRQ(ierr); 79*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ini_x",FILE_MODE_WRITE,&viewer);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = VecView(x,viewer);CHKERRQ(ierr); 81*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 82*c4762a1bSJed Brown 83*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 84*c4762a1bSJed Brown ierr = TSSetDM(ts,user.da);CHKERRQ(ierr); 85*c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 86*c4762a1bSJed Brown ierr = TSSetType(ts,TSARKIMEX);CHKERRQ(ierr); 87*c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,IFunction,&user);CHKERRQ(ierr); 88*c4762a1bSJed Brown /* ierr = TSSetIJacobian(ts,NULL,NULL,IJacobian,&user);CHKERRQ(ierr); */ 89*c4762a1bSJed Brown ierr = TSSetApplicationContext(ts,&user);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 = DMDestroy(&user.da);CHKERRQ(ierr); 101*c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 102*c4762a1bSJed Brown ierr = PetscFinalize(); 103*c4762a1bSJed Brown return ierr; 104*c4762a1bSJed Brown } 105*c4762a1bSJed Brown 106*c4762a1bSJed Brown PetscErrorCode PostStep(TS ts) 107*c4762a1bSJed Brown { 108*c4762a1bSJed Brown PetscErrorCode ierr; 109*c4762a1bSJed Brown Vec X,gc; 110*c4762a1bSJed Brown AppCtx *user; 111*c4762a1bSJed Brown PetscScalar sum = 0,asum; 112*c4762a1bSJed Brown PetscReal t,**p; 113*c4762a1bSJed Brown DMDACoor2d **coors; 114*c4762a1bSJed Brown DM cda; 115*c4762a1bSJed Brown PetscInt i,j,xs,ys,xm,ym; 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown PetscFunctionBegin; 118*c4762a1bSJed Brown ierr = TSGetApplicationContext(ts,&user);CHKERRQ(ierr); 119*c4762a1bSJed Brown ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 120*c4762a1bSJed Brown ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 121*c4762a1bSJed Brown 122*c4762a1bSJed Brown ierr = DMGetCoordinateDM(user->da,&cda);CHKERRQ(ierr); 123*c4762a1bSJed Brown ierr = DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr); 124*c4762a1bSJed Brown ierr = DMGetCoordinates(user->da,&gc);CHKERRQ(ierr); 125*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(cda,gc,&coors);CHKERRQ(ierr); 126*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(user->da,X,&p);CHKERRQ(ierr); 127*c4762a1bSJed Brown for (i=xs; i < xs+xm; i++) { 128*c4762a1bSJed Brown for (j=ys; j < ys+ym; j++) { 129*c4762a1bSJed Brown if (coors[j][i].y < 5) sum += p[j][i]; 130*c4762a1bSJed Brown } 131*c4762a1bSJed Brown } 132*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(cda,gc,&coors);CHKERRQ(ierr); 133*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(user->da,X,&p);CHKERRQ(ierr); 134*c4762a1bSJed Brown ierr = MPI_Allreduce(&sum,&asum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)ts));CHKERRQ(ierr); 135*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"sum(p)*dw*dtheta at t = %f = %f\n",(double)t,(double)(asum));CHKERRQ(ierr); 136*c4762a1bSJed Brown if (sum < 1.0e-2) { 137*c4762a1bSJed Brown ierr = TSSetConvergedReason(ts,TS_CONVERGED_USER);CHKERRQ(ierr); 138*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"Exiting TS as the integral of PDF is almost zero\n");CHKERRQ(ierr); 139*c4762a1bSJed Brown } 140*c4762a1bSJed Brown PetscFunctionReturn(0); 141*c4762a1bSJed Brown } 142*c4762a1bSJed Brown 143*c4762a1bSJed Brown PetscErrorCode ini_bou(Vec X,AppCtx* user) 144*c4762a1bSJed Brown { 145*c4762a1bSJed Brown PetscErrorCode ierr; 146*c4762a1bSJed Brown DM cda; 147*c4762a1bSJed Brown DMDACoor2d **coors; 148*c4762a1bSJed Brown PetscScalar **p; 149*c4762a1bSJed Brown Vec gc; 150*c4762a1bSJed Brown PetscInt i,j; 151*c4762a1bSJed Brown PetscInt xs,ys,xm,ym,M,N; 152*c4762a1bSJed Brown PetscScalar xi,yi; 153*c4762a1bSJed Brown PetscScalar sigmax=user->sigmax,sigmay=user->sigmay; 154*c4762a1bSJed Brown PetscScalar rho =user->rho; 155*c4762a1bSJed Brown PetscScalar muy=user->muy,mux; 156*c4762a1bSJed Brown PetscMPIInt rank; 157*c4762a1bSJed Brown PetscScalar sum; 158*c4762a1bSJed Brown 159*c4762a1bSJed Brown PetscFunctionBeginUser; 160*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 161*c4762a1bSJed Brown ierr = DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 162*c4762a1bSJed Brown user->dx = (user->xmax - user->xmin)/(M-1); user->dy = (user->ymax - user->ymin)/(N-1); 163*c4762a1bSJed Brown ierr = DMGetCoordinateDM(user->da,&cda);CHKERRQ(ierr); 164*c4762a1bSJed Brown ierr = DMGetCoordinates(user->da,&gc);CHKERRQ(ierr); 165*c4762a1bSJed Brown ierr = DMDAVecGetArray(cda,gc,&coors);CHKERRQ(ierr); 166*c4762a1bSJed Brown ierr = DMDAVecGetArray(user->da,X,&p);CHKERRQ(ierr); 167*c4762a1bSJed Brown ierr = DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr); 168*c4762a1bSJed Brown 169*c4762a1bSJed Brown /* mux and muy need to be grid points in the x and y-direction otherwise the solution goes unstable 170*c4762a1bSJed Brown muy is set by choosing the y domain, no. of grid points along y-direction so that muy is a grid point 171*c4762a1bSJed Brown in the y-direction. We only modify mux here 172*c4762a1bSJed Brown */ 173*c4762a1bSJed Brown mux = user->mux = coors[0][M/2+10].x; /* For -pi < x < pi, this should be some angle between 0 and pi/2 */ 174*c4762a1bSJed Brown if (user->nonoiseinitial) { 175*c4762a1bSJed Brown for (i=xs; i < xs+xm; i++) { 176*c4762a1bSJed Brown for (j=ys; j < ys+ym; j++) { 177*c4762a1bSJed Brown xi = coors[j][i].x; yi = coors[j][i].y; 178*c4762a1bSJed Brown if ((xi == mux) && (yi == muy)) { 179*c4762a1bSJed Brown p[j][i] = 1.0; 180*c4762a1bSJed Brown } 181*c4762a1bSJed Brown } 182*c4762a1bSJed Brown } 183*c4762a1bSJed Brown } else { 184*c4762a1bSJed Brown /* Change PM_min accordingly */ 185*c4762a1bSJed Brown user->PM_min = user->Pmax*PetscSinScalar(mux); 186*c4762a1bSJed Brown for (i=xs; i < xs+xm; i++) { 187*c4762a1bSJed Brown for (j=ys; j < ys+ym; j++) { 188*c4762a1bSJed Brown xi = coors[j][i].x; yi = coors[j][i].y; 189*c4762a1bSJed Brown 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))); 190*c4762a1bSJed Brown } 191*c4762a1bSJed Brown } 192*c4762a1bSJed Brown } 193*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(cda,gc,&coors);CHKERRQ(ierr); 194*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(user->da,X,&p);CHKERRQ(ierr); 195*c4762a1bSJed Brown ierr = VecSum(X,&sum);CHKERRQ(ierr); 196*c4762a1bSJed Brown ierr = VecScale(X,1.0/sum);CHKERRQ(ierr); 197*c4762a1bSJed Brown PetscFunctionReturn(0); 198*c4762a1bSJed Brown } 199*c4762a1bSJed Brown 200*c4762a1bSJed Brown /* First advection term */ 201*c4762a1bSJed Brown PetscErrorCode adv1(PetscScalar **p,PetscScalar y,PetscInt i,PetscInt j,PetscInt M,PetscScalar *p1,AppCtx *user) 202*c4762a1bSJed Brown { 203*c4762a1bSJed Brown PetscScalar f,fpos,fneg; 204*c4762a1bSJed Brown PetscFunctionBegin; 205*c4762a1bSJed Brown f = (y - user->ws); 206*c4762a1bSJed Brown fpos = PetscMax(f,0); 207*c4762a1bSJed Brown fneg = PetscMin(f,0); 208*c4762a1bSJed Brown if (user->st_width == 1) { 209*c4762a1bSJed Brown *p1 = fpos*(p[j][i] - p[j][i-1])/user->dx + fneg*(p[j][i+1] - p[j][i])/user->dx; 210*c4762a1bSJed Brown } else if (user->st_width == 2) { 211*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); 212*c4762a1bSJed Brown } else if (user->st_width == 3) { 213*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); 214*c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No support for wider stencils"); 215*c4762a1bSJed Brown PetscFunctionReturn(0); 216*c4762a1bSJed Brown } 217*c4762a1bSJed Brown 218*c4762a1bSJed Brown /* Second advection term */ 219*c4762a1bSJed Brown PetscErrorCode adv2(PetscScalar **p,PetscScalar x,PetscInt i,PetscInt j,PetscInt N,PetscScalar *p2,AppCtx *user) 220*c4762a1bSJed Brown { 221*c4762a1bSJed Brown PetscScalar f,fpos,fneg; 222*c4762a1bSJed Brown PetscFunctionBegin; 223*c4762a1bSJed Brown f = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(x)); 224*c4762a1bSJed Brown fpos = PetscMax(f,0); 225*c4762a1bSJed Brown fneg = PetscMin(f,0); 226*c4762a1bSJed Brown if (user->st_width == 1) { 227*c4762a1bSJed Brown *p2 = fpos*(p[j][i] - p[j-1][i])/user->dy + fneg*(p[j+1][i] - p[j][i])/user->dy; 228*c4762a1bSJed Brown } else if (user->st_width ==2) { 229*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); 230*c4762a1bSJed Brown } else if (user->st_width == 3) { 231*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); 232*c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No support for wider stencils"); 233*c4762a1bSJed Brown PetscFunctionReturn(0); 234*c4762a1bSJed Brown } 235*c4762a1bSJed Brown 236*c4762a1bSJed Brown /* Diffusion term */ 237*c4762a1bSJed Brown PetscErrorCode diffuse(PetscScalar **p,PetscInt i,PetscInt j,PetscReal t,PetscScalar *p_diff,AppCtx * user) 238*c4762a1bSJed Brown { 239*c4762a1bSJed Brown PetscFunctionBeginUser; 240*c4762a1bSJed Brown if (user->st_width == 1) { 241*c4762a1bSJed Brown *p_diff = user->disper_coe*((p[j-1][i] - 2*p[j][i] + p[j+1][i])/(user->dy*user->dy)); 242*c4762a1bSJed Brown } else if (user->st_width == 2) { 243*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)); 244*c4762a1bSJed Brown } else if (user->st_width == 3) { 245*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)); 246*c4762a1bSJed Brown } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No support for wider stencils"); 247*c4762a1bSJed Brown PetscFunctionReturn(0); 248*c4762a1bSJed Brown } 249*c4762a1bSJed Brown 250*c4762a1bSJed Brown PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 251*c4762a1bSJed Brown { 252*c4762a1bSJed Brown PetscErrorCode ierr; 253*c4762a1bSJed Brown AppCtx *user=(AppCtx*)ctx; 254*c4762a1bSJed Brown DM cda; 255*c4762a1bSJed Brown DMDACoor2d **coors; 256*c4762a1bSJed Brown PetscScalar **p,**f,**pdot; 257*c4762a1bSJed Brown PetscInt i,j; 258*c4762a1bSJed Brown PetscInt xs,ys,xm,ym,M,N; 259*c4762a1bSJed Brown Vec localX,gc,localXdot; 260*c4762a1bSJed Brown PetscScalar p_adv1,p_adv2,p_diff; 261*c4762a1bSJed Brown 262*c4762a1bSJed Brown PetscFunctionBeginUser; 263*c4762a1bSJed Brown ierr = DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 264*c4762a1bSJed Brown ierr = DMGetCoordinateDM(user->da,&cda);CHKERRQ(ierr); 265*c4762a1bSJed Brown ierr = DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr); 266*c4762a1bSJed Brown 267*c4762a1bSJed Brown ierr = DMGetLocalVector(user->da,&localX);CHKERRQ(ierr); 268*c4762a1bSJed Brown ierr = DMGetLocalVector(user->da,&localXdot);CHKERRQ(ierr); 269*c4762a1bSJed Brown 270*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 271*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);CHKERRQ(ierr); 272*c4762a1bSJed Brown ierr = DMGlobalToLocalBegin(user->da,Xdot,INSERT_VALUES,localXdot);CHKERRQ(ierr); 273*c4762a1bSJed Brown ierr = DMGlobalToLocalEnd(user->da,Xdot,INSERT_VALUES,localXdot);CHKERRQ(ierr); 274*c4762a1bSJed Brown 275*c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(user->da,&gc);CHKERRQ(ierr); 276*c4762a1bSJed Brown 277*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(cda,gc,&coors);CHKERRQ(ierr); 278*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(user->da,localX,&p);CHKERRQ(ierr); 279*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(user->da,localXdot,&pdot);CHKERRQ(ierr); 280*c4762a1bSJed Brown ierr = DMDAVecGetArray(user->da,F,&f);CHKERRQ(ierr); 281*c4762a1bSJed Brown 282*c4762a1bSJed Brown user->disper_coe = PetscPowScalar((user->lambda*user->ws)/(2*user->H),2)*user->q*(1.0-PetscExpScalar(-t/user->lambda)); 283*c4762a1bSJed Brown for (i=xs; i < xs+xm; i++) { 284*c4762a1bSJed Brown for (j=ys; j < ys+ym; j++) { 285*c4762a1bSJed Brown ierr = adv1(p,coors[j][i].y,i,j,M,&p_adv1,user);CHKERRQ(ierr); 286*c4762a1bSJed Brown ierr = adv2(p,coors[j][i].x,i,j,N,&p_adv2,user);CHKERRQ(ierr); 287*c4762a1bSJed Brown ierr = diffuse(p,i,j,t,&p_diff,user);CHKERRQ(ierr); 288*c4762a1bSJed Brown f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i]; 289*c4762a1bSJed Brown } 290*c4762a1bSJed Brown } 291*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(user->da,localX,&p);CHKERRQ(ierr); 292*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(user->da,localX,&pdot);CHKERRQ(ierr); 293*c4762a1bSJed Brown ierr = DMRestoreLocalVector(user->da,&localX);CHKERRQ(ierr); 294*c4762a1bSJed Brown ierr = DMRestoreLocalVector(user->da,&localXdot);CHKERRQ(ierr); 295*c4762a1bSJed Brown ierr = DMDAVecRestoreArray(user->da,F,&f);CHKERRQ(ierr); 296*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(cda,gc,&coors);CHKERRQ(ierr); 297*c4762a1bSJed Brown 298*c4762a1bSJed Brown PetscFunctionReturn(0); 299*c4762a1bSJed Brown } 300*c4762a1bSJed Brown 301*c4762a1bSJed Brown PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ctx) 302*c4762a1bSJed Brown { 303*c4762a1bSJed Brown PetscErrorCode ierr; 304*c4762a1bSJed Brown AppCtx *user=(AppCtx*)ctx; 305*c4762a1bSJed Brown DM cda; 306*c4762a1bSJed Brown DMDACoor2d **coors; 307*c4762a1bSJed Brown PetscInt i,j; 308*c4762a1bSJed Brown PetscInt xs,ys,xm,ym,M,N; 309*c4762a1bSJed Brown Vec gc; 310*c4762a1bSJed Brown PetscScalar val[5],xi,yi; 311*c4762a1bSJed Brown MatStencil row,col[5]; 312*c4762a1bSJed Brown PetscScalar c1,c3,c5,c1pos,c1neg,c3pos,c3neg; 313*c4762a1bSJed Brown 314*c4762a1bSJed Brown PetscFunctionBeginUser; 315*c4762a1bSJed Brown ierr = DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr); 316*c4762a1bSJed Brown ierr = DMGetCoordinateDM(user->da,&cda);CHKERRQ(ierr); 317*c4762a1bSJed Brown ierr = DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);CHKERRQ(ierr); 318*c4762a1bSJed Brown 319*c4762a1bSJed Brown ierr = DMGetCoordinatesLocal(user->da,&gc);CHKERRQ(ierr); 320*c4762a1bSJed Brown ierr = DMDAVecGetArrayRead(cda,gc,&coors);CHKERRQ(ierr); 321*c4762a1bSJed Brown for (i=xs; i < xs+xm; i++) { 322*c4762a1bSJed Brown for (j=ys; j < ys+ym; j++) { 323*c4762a1bSJed Brown PetscInt nc = 0; 324*c4762a1bSJed Brown xi = coors[j][i].x; yi = coors[j][i].y; 325*c4762a1bSJed Brown row.i = i; row.j = j; 326*c4762a1bSJed Brown c1 = (yi-user->ws)/user->dx; 327*c4762a1bSJed Brown c1pos = PetscMax(c1,0); 328*c4762a1bSJed Brown c1neg = PetscMin(c1,0); 329*c4762a1bSJed Brown c3 = (user->ws/(2.0*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(xi))/user->dy; 330*c4762a1bSJed Brown c3pos = PetscMax(c3,0); 331*c4762a1bSJed Brown c3neg = PetscMin(c3,0); 332*c4762a1bSJed Brown c5 = (PetscPowScalar((user->lambda*user->ws)/(2*user->H),2)*user->q*(1.0-PetscExpScalar(-t/user->lambda)))/(user->dy*user->dy); 333*c4762a1bSJed Brown col[nc].i = i-1; col[nc].j = j; val[nc++] = c1pos; 334*c4762a1bSJed Brown col[nc].i = i+1; col[nc].j = j; val[nc++] = -c1neg; 335*c4762a1bSJed Brown col[nc].i = i; col[nc].j = j-1; val[nc++] = c3pos + c5; 336*c4762a1bSJed Brown col[nc].i = i; col[nc].j = j+1; val[nc++] = -c3neg + c5; 337*c4762a1bSJed Brown col[nc].i = i; col[nc].j = j; val[nc++] = -c1pos + c1neg -c3pos + c3neg -2*c5 -a; 338*c4762a1bSJed Brown ierr = MatSetValuesStencil(Jpre,1,&row,nc,col,val,INSERT_VALUES);CHKERRQ(ierr); 339*c4762a1bSJed Brown } 340*c4762a1bSJed Brown } 341*c4762a1bSJed Brown ierr = DMDAVecRestoreArrayRead(cda,gc,&coors);CHKERRQ(ierr); 342*c4762a1bSJed Brown 343*c4762a1bSJed Brown ierr = MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 344*c4762a1bSJed Brown ierr = MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 345*c4762a1bSJed Brown if (J != Jpre) { 346*c4762a1bSJed Brown ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 347*c4762a1bSJed Brown ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 348*c4762a1bSJed Brown } 349*c4762a1bSJed Brown PetscFunctionReturn(0); 350*c4762a1bSJed Brown } 351*c4762a1bSJed Brown 352*c4762a1bSJed Brown 353*c4762a1bSJed Brown 354*c4762a1bSJed Brown PetscErrorCode Parameter_settings(AppCtx *user) 355*c4762a1bSJed Brown { 356*c4762a1bSJed Brown PetscErrorCode ierr; 357*c4762a1bSJed Brown PetscBool flg; 358*c4762a1bSJed Brown 359*c4762a1bSJed Brown PetscFunctionBeginUser; 360*c4762a1bSJed Brown 361*c4762a1bSJed Brown /* Set default parameters */ 362*c4762a1bSJed Brown user->ws = 1.0; 363*c4762a1bSJed Brown user->H = 5.0; 364*c4762a1bSJed Brown user->Pmax = 2.1; 365*c4762a1bSJed Brown user->PM_min = 1.0; 366*c4762a1bSJed Brown user->lambda = 0.1; 367*c4762a1bSJed Brown user->q = 1.0; 368*c4762a1bSJed Brown user->mux = PetscAsinScalar(user->PM_min/user->Pmax); 369*c4762a1bSJed Brown user->sigmax = 0.1; 370*c4762a1bSJed Brown user->sigmay = 0.1; 371*c4762a1bSJed Brown user->rho = 0.0; 372*c4762a1bSJed Brown user->xmin = -PETSC_PI; 373*c4762a1bSJed Brown user->xmax = PETSC_PI; 374*c4762a1bSJed Brown user->bx = DM_BOUNDARY_PERIODIC; 375*c4762a1bSJed Brown user->by = DM_BOUNDARY_MIRROR; 376*c4762a1bSJed Brown user->nonoiseinitial = PETSC_FALSE; 377*c4762a1bSJed Brown 378*c4762a1bSJed Brown /* 379*c4762a1bSJed Brown ymin of -3 seems to let the unstable solution move up and leave a zero in its wake 380*c4762a1bSJed Brown with an ymin of -1 the wake is never exactly zero 381*c4762a1bSJed Brown */ 382*c4762a1bSJed Brown user->ymin = -3.0; 383*c4762a1bSJed Brown user->ymax = 10.0; 384*c4762a1bSJed Brown user->st_width = 1; 385*c4762a1bSJed Brown 386*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-ws",&user->ws,&flg);CHKERRQ(ierr); 387*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-Inertia",&user->H,&flg);CHKERRQ(ierr); 388*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-Pmax",&user->Pmax,&flg);CHKERRQ(ierr); 389*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-PM_min",&user->PM_min,&flg);CHKERRQ(ierr); 390*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-lambda",&user->lambda,&flg);CHKERRQ(ierr); 391*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-q",&user->q,&flg);CHKERRQ(ierr); 392*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-mux",&user->mux,&flg);CHKERRQ(ierr); 393*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-sigmax",&user->sigmax,&flg);CHKERRQ(ierr); 394*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-muy",&user->muy,&flg);CHKERRQ(ierr); 395*c4762a1bSJed Brown if (flg == 0) { 396*c4762a1bSJed Brown user->muy = user->ws; 397*c4762a1bSJed Brown } 398*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-sigmay",&user->sigmay,&flg);CHKERRQ(ierr); 399*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-rho",&user->rho,&flg);CHKERRQ(ierr); 400*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-xmin",&user->xmin,&flg);CHKERRQ(ierr); 401*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-xmax",&user->xmax,&flg);CHKERRQ(ierr); 402*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-ymin",&user->ymin,&flg);CHKERRQ(ierr); 403*c4762a1bSJed Brown ierr = PetscOptionsGetScalar(NULL,NULL,"-ymax",&user->ymax,&flg);CHKERRQ(ierr); 404*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-stencil_width",&user->st_width,&flg);CHKERRQ(ierr); 405*c4762a1bSJed Brown ierr = PetscOptionsGetEnum(NULL,NULL,"-bx",DMBoundaryTypes,(PetscEnum*)&user->bx,&flg);CHKERRQ(ierr); 406*c4762a1bSJed Brown ierr = PetscOptionsGetEnum(NULL,NULL,"-by",DMBoundaryTypes,(PetscEnum*)&user->by,&flg);CHKERRQ(ierr); 407*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-nonoiseinitial",&user->nonoiseinitial,&flg);CHKERRQ(ierr); 408*c4762a1bSJed Brown 409*c4762a1bSJed Brown PetscFunctionReturn(0); 410*c4762a1bSJed Brown } 411*c4762a1bSJed Brown 412*c4762a1bSJed Brown 413*c4762a1bSJed Brown /*TEST 414*c4762a1bSJed Brown 415*c4762a1bSJed Brown build: 416*c4762a1bSJed Brown requires: !complex !single 417*c4762a1bSJed Brown 418*c4762a1bSJed Brown test: 419*c4762a1bSJed Brown args: -ts_max_steps 2 420*c4762a1bSJed Brown localrunfiles: petscopt_ex7 421*c4762a1bSJed Brown 422*c4762a1bSJed Brown test: 423*c4762a1bSJed Brown suffix: 2 424*c4762a1bSJed Brown args: -ts_max_steps 2 -snes_mf_operator 425*c4762a1bSJed Brown output_file: output/ex7_1.out 426*c4762a1bSJed Brown localrunfiles: petscopt_ex7 427*c4762a1bSJed Brown timeoutfactor: 2 428*c4762a1bSJed Brown 429*c4762a1bSJed Brown test: 430*c4762a1bSJed Brown suffix: 3 431*c4762a1bSJed Brown args: -ts_max_steps 2 -snes_mf -pc_type none 432*c4762a1bSJed Brown output_file: output/ex7_1.out 433*c4762a1bSJed Brown localrunfiles: petscopt_ex7 434*c4762a1bSJed Brown timeoutfactor: 2 435*c4762a1bSJed Brown 436*c4762a1bSJed Brown TEST*/ 437