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