1c4762a1bSJed Brown static char help[] = "Time-dependent PDE in 2d for calculating joint PDF. \n"; 2c4762a1bSJed Brown /* 3c4762a1bSJed Brown p_t = -x_t*p_x -y_t*p_y + f(t)*p_yy 4c4762a1bSJed Brown xmin < x < xmax, ymin < y < ymax; 5c4762a1bSJed Brown x_t = (y - ws) y_t = (ws/2H)*(Pm - Pmax*sin(x)) 6c4762a1bSJed Brown 7c4762a1bSJed Brown Boundary conditions: -bc_type 0 => Zero dirichlet boundary 8c4762a1bSJed Brown -bc_type 1 => Steady state boundary condition 9c4762a1bSJed Brown Steady state boundary condition found by setting p_t = 0 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown 12c4762a1bSJed Brown #include <petscdm.h> 13c4762a1bSJed Brown #include <petscdmda.h> 14c4762a1bSJed Brown #include <petscts.h> 15c4762a1bSJed Brown 16c4762a1bSJed Brown /* 17c4762a1bSJed Brown User-defined data structures and routines 18c4762a1bSJed Brown */ 19c4762a1bSJed Brown typedef struct { 20c4762a1bSJed Brown PetscScalar ws; /* Synchronous speed */ 21c4762a1bSJed Brown PetscScalar H; /* Inertia constant */ 22c4762a1bSJed Brown PetscScalar D; /* Damping constant */ 23c4762a1bSJed Brown PetscScalar Pmax; /* Maximum power output of generator */ 24c4762a1bSJed Brown PetscScalar PM_min; /* Mean mechanical power input */ 25c4762a1bSJed Brown PetscScalar lambda; /* correlation time */ 26c4762a1bSJed Brown PetscScalar q; /* noise strength */ 27c4762a1bSJed Brown PetscScalar mux; /* Initial average angle */ 28c4762a1bSJed Brown PetscScalar sigmax; /* Standard deviation of initial angle */ 29c4762a1bSJed Brown PetscScalar muy; /* Average speed */ 30c4762a1bSJed Brown PetscScalar sigmay; /* standard deviation of initial speed */ 31c4762a1bSJed Brown PetscScalar rho; /* Cross-correlation coefficient */ 32c4762a1bSJed Brown PetscScalar t0; /* Initial time */ 33c4762a1bSJed Brown PetscScalar tmax; /* Final time */ 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 PetscInt bc; /* Boundary conditions */ 41c4762a1bSJed Brown PetscScalar disper_coe; /* Dispersion coefficient */ 42c4762a1bSJed Brown DM da; 43c4762a1bSJed Brown } AppCtx; 44c4762a1bSJed Brown 45c4762a1bSJed Brown PetscErrorCode Parameter_settings(AppCtx*); 46c4762a1bSJed Brown PetscErrorCode ini_bou(Vec,AppCtx*); 47c4762a1bSJed Brown PetscErrorCode IFunction(TS,PetscReal,Vec,Vec,Vec,void*); 48c4762a1bSJed Brown PetscErrorCode IJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*); 49c4762a1bSJed Brown PetscErrorCode PostStep(TS); 50c4762a1bSJed Brown 51c4762a1bSJed Brown int main(int argc, char **argv) 52c4762a1bSJed Brown { 53c4762a1bSJed Brown Vec x; /* Solution vector */ 54c4762a1bSJed Brown TS ts; /* Time-stepping context */ 55c4762a1bSJed Brown AppCtx user; /* Application context */ 56c4762a1bSJed Brown Mat J; 57c4762a1bSJed Brown PetscViewer viewer; 58c4762a1bSJed Brown 59*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,"petscopt_ex6", help)); 60c4762a1bSJed Brown /* Get physics and time parameters */ 615f80ce2aSJacob Faibussowitsch CHKERRQ(Parameter_settings(&user)); 62c4762a1bSJed Brown /* Create a 2D DA with dof = 1 */ 635f80ce2aSJacob Faibussowitsch CHKERRQ(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)); 645f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(user.da)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(user.da)); 66c4762a1bSJed Brown /* Set x and y coordinates */ 675f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(user.da,user.xmin,user.xmax,user.ymin,user.ymax,0.0,1.0)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Get global vector x from DM */ 705f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(user.da,&x)); 71c4762a1bSJed Brown 725f80ce2aSJacob Faibussowitsch CHKERRQ(ini_bou(x,&user)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ini_x",FILE_MODE_WRITE,&viewer)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(x,viewer)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 76c4762a1bSJed Brown 77c4762a1bSJed Brown /* Get Jacobian matrix structure from the da */ 785f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetMatType(user.da,MATAIJ)); 795f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateMatrix(user.da,&J)); 80c4762a1bSJed Brown 815f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 825f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(ts,NULL,IFunction,&user)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIJacobian(ts,J,J,IJacobian,&user)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetApplicationContext(ts,&user)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,user.tmax)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 885f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTime(ts,user.t0)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,.005)); 905f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 915f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetPostStep(ts,PostStep)); 925f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,x)); 93c4762a1bSJed Brown 945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"fin_x",FILE_MODE_WRITE,&viewer)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(x,viewer)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 97c4762a1bSJed Brown 985f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&J)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&user.da)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 102*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 103*b122ec5aSJacob Faibussowitsch return 0; 104c4762a1bSJed Brown } 105c4762a1bSJed Brown 106c4762a1bSJed Brown PetscErrorCode PostStep(TS ts) 107c4762a1bSJed Brown { 108c4762a1bSJed Brown Vec X; 109c4762a1bSJed Brown AppCtx *user; 110c4762a1bSJed Brown PetscScalar sum; 111c4762a1bSJed Brown PetscReal t; 112c4762a1bSJed Brown 113c4762a1bSJed Brown PetscFunctionBegin; 1145f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetApplicationContext(ts,&user)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTime(ts,&t)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolution(ts,&X)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(VecSum(X,&sum)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"sum(p)*dw*dtheta at t = %3.2f = %3.6f\n",(double)t,(double)sum*user->dx*user->dy)); 119c4762a1bSJed Brown PetscFunctionReturn(0); 120c4762a1bSJed Brown } 121c4762a1bSJed Brown 122c4762a1bSJed Brown PetscErrorCode ini_bou(Vec X,AppCtx* user) 123c4762a1bSJed Brown { 124c4762a1bSJed Brown DM cda; 125c4762a1bSJed Brown DMDACoor2d **coors; 126c4762a1bSJed Brown PetscScalar **p; 127c4762a1bSJed Brown Vec gc; 128c4762a1bSJed Brown PetscInt i,j; 129c4762a1bSJed Brown PetscInt xs,ys,xm,ym,M,N; 130c4762a1bSJed Brown PetscScalar xi,yi; 131c4762a1bSJed Brown PetscScalar sigmax=user->sigmax,sigmay=user->sigmay; 132c4762a1bSJed Brown PetscScalar rho =user->rho; 133c4762a1bSJed Brown PetscScalar mux =user->mux,muy=user->muy; 134c4762a1bSJed Brown PetscMPIInt rank; 135c4762a1bSJed Brown 136c4762a1bSJed Brown PetscFunctionBeginUser; 1375f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 1385f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 139c4762a1bSJed Brown user->dx = (user->xmax - user->xmin)/(M-1); user->dy = (user->ymax - user->ymin)/(N-1); 1405f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(user->da,&cda)); 1415f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinates(user->da,&gc)); 1425f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(cda,gc,&coors)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArray(user->da,X,&p)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0)); 145c4762a1bSJed Brown for (i=xs; i < xs+xm; i++) { 146c4762a1bSJed Brown for (j=ys; j < ys+ym; j++) { 147c4762a1bSJed Brown xi = coors[j][i].x; yi = coors[j][i].y; 148c4762a1bSJed Brown if (i == 0 || j == 0 || i == M-1 || j == N-1) p[j][i] = 0.0; 149c4762a1bSJed 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))); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown } 152c4762a1bSJed Brown /* p[N/2+N%2][M/2+M%2] = 1/(user->dx*user->dy); */ 153c4762a1bSJed Brown 1545f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(cda,gc,&coors)); 1555f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(user->da,X,&p)); 156c4762a1bSJed Brown PetscFunctionReturn(0); 157c4762a1bSJed Brown } 158c4762a1bSJed Brown 159c4762a1bSJed Brown /* First advection term */ 160c4762a1bSJed Brown PetscErrorCode adv1(PetscScalar **p,PetscScalar y,PetscInt i,PetscInt j,PetscInt M,PetscScalar *p1,AppCtx *user) 161c4762a1bSJed Brown { 162c4762a1bSJed Brown PetscScalar f; 163c4762a1bSJed Brown /* PetscScalar v1,v2,v3,v4,v5,s1,s2,s3; */ 164c4762a1bSJed Brown PetscFunctionBegin; 165c4762a1bSJed Brown /* if (i > 2 && i < M-2) { 166c4762a1bSJed Brown v1 = (y-user->ws)*(p[j][i-2] - p[j][i-3])/user->dx; 167c4762a1bSJed Brown v2 = (y-user->ws)*(p[j][i-1] - p[j][i-2])/user->dx; 168c4762a1bSJed Brown v3 = (y-user->ws)*(p[j][i] - p[j][i-1])/user->dx; 169c4762a1bSJed Brown v4 = (y-user->ws)*(p[j][i+1] - p[j][i])/user->dx; 170c4762a1bSJed Brown v5 = (y-user->ws)*(p[j][i+1] - p[j][i+2])/user->dx; 171c4762a1bSJed Brown 172c4762a1bSJed Brown s1 = v1/3.0 - (7.0/6.0)*v2 + (11.0/6.0)*v3; 173c4762a1bSJed Brown s2 =-v2/6.0 + (5.0/6.0)*v3 + (1.0/3.0)*v4; 174c4762a1bSJed Brown s3 = v3/3.0 + (5.0/6.0)*v4 - (1.0/6.0)*v5; 175c4762a1bSJed Brown 176c4762a1bSJed Brown *p1 = 0.1*s1 + 0.6*s2 + 0.3*s3; 177c4762a1bSJed Brown } else *p1 = 0.0; */ 178c4762a1bSJed Brown f = (y - user->ws); 179c4762a1bSJed Brown *p1 = f*(p[j][i+1] - p[j][i-1])/(2*user->dx); 180c4762a1bSJed Brown PetscFunctionReturn(0); 181c4762a1bSJed Brown } 182c4762a1bSJed Brown 183c4762a1bSJed Brown /* Second advection term */ 184c4762a1bSJed Brown PetscErrorCode adv2(PetscScalar **p,PetscScalar x,PetscInt i,PetscInt j,PetscInt N,PetscScalar *p2,AppCtx *user) 185c4762a1bSJed Brown { 186c4762a1bSJed Brown PetscScalar f; 187c4762a1bSJed Brown /* PetscScalar v1,v2,v3,v4,v5,s1,s2,s3; */ 188c4762a1bSJed Brown PetscFunctionBegin; 189c4762a1bSJed Brown /* if (j > 2 && j < N-2) { 190c4762a1bSJed Brown v1 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j-2][i] - p[j-3][i])/user->dy; 191c4762a1bSJed Brown v2 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j-1][i] - p[j-2][i])/user->dy; 192c4762a1bSJed Brown v3 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j][i] - p[j-1][i])/user->dy; 193c4762a1bSJed Brown v4 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j+1][i] - p[j][i])/user->dy; 194c4762a1bSJed Brown v5 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j+2][i] - p[j+1][i])/user->dy; 195c4762a1bSJed Brown 196c4762a1bSJed Brown s1 = v1/3.0 - (7.0/6.0)*v2 + (11.0/6.0)*v3; 197c4762a1bSJed Brown s2 =-v2/6.0 + (5.0/6.0)*v3 + (1.0/3.0)*v4; 198c4762a1bSJed Brown s3 = v3/3.0 + (5.0/6.0)*v4 - (1.0/6.0)*v5; 199c4762a1bSJed Brown 200c4762a1bSJed Brown *p2 = 0.1*s1 + 0.6*s2 + 0.3*s3; 201c4762a1bSJed Brown } else *p2 = 0.0; */ 202c4762a1bSJed Brown f = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(x)); 203c4762a1bSJed Brown *p2 = f*(p[j+1][i] - p[j-1][i])/(2*user->dy); 204c4762a1bSJed Brown PetscFunctionReturn(0); 205c4762a1bSJed Brown } 206c4762a1bSJed Brown 207c4762a1bSJed Brown /* Diffusion term */ 208c4762a1bSJed Brown PetscErrorCode diffuse(PetscScalar **p,PetscInt i,PetscInt j,PetscReal t,PetscScalar *p_diff,AppCtx * user) 209c4762a1bSJed Brown { 210c4762a1bSJed Brown PetscFunctionBeginUser; 211c4762a1bSJed Brown 212c4762a1bSJed Brown *p_diff = user->disper_coe*((p[j-1][i] - 2*p[j][i] + p[j+1][i])/(user->dy*user->dy)); 213c4762a1bSJed Brown PetscFunctionReturn(0); 214c4762a1bSJed Brown } 215c4762a1bSJed Brown 216c4762a1bSJed Brown PetscErrorCode BoundaryConditions(PetscScalar **p,DMDACoor2d **coors,PetscInt i,PetscInt j,PetscInt M, PetscInt N,PetscScalar **f,AppCtx *user) 217c4762a1bSJed Brown { 218c4762a1bSJed Brown PetscScalar fwc,fthetac; 219c4762a1bSJed Brown PetscScalar w=coors[j][i].y,theta=coors[j][i].x; 220c4762a1bSJed Brown 221c4762a1bSJed Brown PetscFunctionBeginUser; 222c4762a1bSJed Brown if (user->bc == 0) { /* Natural boundary condition */ 223c4762a1bSJed Brown f[j][i] = p[j][i]; 224c4762a1bSJed Brown } else { /* Steady state boundary condition */ 225c4762a1bSJed Brown fthetac = user->ws/(2*user->H)*(user->PM_min - user->Pmax*PetscSinScalar(theta)); 226c4762a1bSJed Brown fwc = (w*w/2.0 - user->ws*w); 227c4762a1bSJed Brown if (i == 0 && j == 0) { /* left bottom corner */ 228c4762a1bSJed 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; 229c4762a1bSJed Brown } else if (i == 0 && j == N-1) { /* right bottom corner */ 230c4762a1bSJed 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; 231c4762a1bSJed Brown } else if (i == M-1 && j == 0) { /* left top corner */ 232c4762a1bSJed 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; 233c4762a1bSJed Brown } else if (i == M-1 && j == N-1) { /* right top corner */ 234c4762a1bSJed 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; 235c4762a1bSJed Brown } else if (i == 0) { /* Bottom edge */ 236c4762a1bSJed 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); 237c4762a1bSJed Brown } else if (i == M-1) { /* Top edge */ 238c4762a1bSJed 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); 239c4762a1bSJed Brown } else if (j == 0) { /* Left edge */ 240c4762a1bSJed 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); 241c4762a1bSJed Brown } else if (j == N-1) { /* Right edge */ 242c4762a1bSJed 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); 243c4762a1bSJed Brown } 244c4762a1bSJed Brown } 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; 2605f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 2615f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(user->da,&cda)); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0)); 263c4762a1bSJed Brown 2645f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(user->da,&localX)); 2655f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetLocalVector(user->da,&localXdot)); 266c4762a1bSJed Brown 2675f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX)); 2685f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX)); 2695f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalBegin(user->da,Xdot,INSERT_VALUES,localXdot)); 2705f80ce2aSJacob Faibussowitsch CHKERRQ(DMGlobalToLocalEnd(user->da,Xdot,INSERT_VALUES,localXdot)); 271c4762a1bSJed Brown 2725f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(user->da,&gc)); 273c4762a1bSJed Brown 2745f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(cda,gc,&coors)); 2755f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(user->da,localX,&p)); 2765f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(user->da,localXdot,&pdot)); 2775f80ce2aSJacob 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++) { 282c4762a1bSJed Brown if (i == 0 || j == 0 || i == M-1 || j == N-1) { 2835f80ce2aSJacob Faibussowitsch CHKERRQ(BoundaryConditions(p,coors,i,j,M,N,f,user)); 284c4762a1bSJed Brown } else { 2855f80ce2aSJacob Faibussowitsch CHKERRQ(adv1(p,coors[j][i].y,i,j,M,&p_adv1,user)); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(adv2(p,coors[j][i].x,i,j,N,&p_adv2,user)); 2875f80ce2aSJacob Faibussowitsch CHKERRQ(diffuse(p,i,j,t,&p_diff,user)); 288c4762a1bSJed Brown f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i]; 289c4762a1bSJed Brown } 290c4762a1bSJed Brown } 291c4762a1bSJed Brown } 2925f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(user->da,localX,&p)); 2935f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(user->da,localX,&pdot)); 2945f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(user->da,&localX)); 2955f80ce2aSJacob Faibussowitsch CHKERRQ(DMRestoreLocalVector(user->da,&localXdot)); 2965f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(user->da,F,&f)); 2975f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(cda,gc,&coors)); 298c4762a1bSJed Brown 299c4762a1bSJed Brown PetscFunctionReturn(0); 300c4762a1bSJed Brown } 301c4762a1bSJed Brown 302c4762a1bSJed Brown PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ctx) 303c4762a1bSJed Brown { 304c4762a1bSJed Brown AppCtx *user=(AppCtx*)ctx; 305c4762a1bSJed Brown DM cda; 306c4762a1bSJed Brown DMDACoor2d **coors; 307c4762a1bSJed Brown PetscInt i,j; 308c4762a1bSJed Brown PetscInt xs,ys,xm,ym,M,N; 309c4762a1bSJed Brown Vec gc; 310c4762a1bSJed Brown PetscScalar val[5],xi,yi; 311c4762a1bSJed Brown MatStencil row,col[5]; 312c4762a1bSJed Brown PetscScalar c1,c3,c5; 313c4762a1bSJed Brown 314c4762a1bSJed Brown PetscFunctionBeginUser; 3155f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL)); 3165f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinateDM(user->da,&cda)); 3175f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0)); 318c4762a1bSJed Brown 3195f80ce2aSJacob Faibussowitsch CHKERRQ(DMGetCoordinatesLocal(user->da,&gc)); 3205f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecGetArrayRead(cda,gc,&coors)); 321c4762a1bSJed Brown for (i=xs; i < xs+xm; i++) { 322c4762a1bSJed Brown for (j=ys; j < ys+ym; j++) { 323c4762a1bSJed Brown PetscInt nc = 0; 324c4762a1bSJed Brown xi = coors[j][i].x; yi = coors[j][i].y; 325c4762a1bSJed Brown row.i = i; row.j = j; 326c4762a1bSJed Brown if (i == 0 || j == 0 || i == M-1 || j == N-1) { 327c4762a1bSJed Brown if (user->bc == 0) { 328c4762a1bSJed Brown col[nc].i = i; col[nc].j = j; val[nc++] = 1.0; 329c4762a1bSJed Brown } else { 330c4762a1bSJed Brown PetscScalar fthetac,fwc; 331c4762a1bSJed Brown fthetac = user->ws/(2*user->H)*(user->PM_min - user->Pmax*PetscSinScalar(xi)); 332c4762a1bSJed Brown fwc = (yi*yi/2.0 - user->ws*yi); 333c4762a1bSJed Brown if (i==0 && j==0) { 334c4762a1bSJed Brown col[nc].i = i+1; col[nc].j = j; val[nc++] = fwc/user->dx; 335c4762a1bSJed Brown col[nc].i = i; col[nc].j = j+1; val[nc++] = -user->disper_coe/user->dy; 336c4762a1bSJed Brown col[nc].i = i; col[nc].j = j; val[nc++] = -fwc/user->dx + fthetac + user->disper_coe/user->dy; 337c4762a1bSJed Brown } else if (i==0 && j == N-1) { 338c4762a1bSJed Brown col[nc].i = i+1; col[nc].j = j; val[nc++] = fwc/user->dx; 339c4762a1bSJed Brown col[nc].i = i; col[nc].j = j-1; val[nc++] = user->disper_coe/user->dy; 340c4762a1bSJed Brown col[nc].i = i; col[nc].j = j; val[nc++] = -fwc/user->dx + fthetac - user->disper_coe/user->dy; 341c4762a1bSJed Brown } else if (i== M-1 && j == 0) { 342c4762a1bSJed Brown col[nc].i = i-1; col[nc].j = j; val[nc++] = -fwc/user->dx; 343c4762a1bSJed Brown col[nc].i = i; col[nc].j = j+1; val[nc++] = -user->disper_coe/user->dy; 344c4762a1bSJed Brown col[nc].i = i; col[nc].j = j; val[nc++] = fwc/user->dx + fthetac + user->disper_coe/user->dy; 345c4762a1bSJed Brown } else if (i == M-1 && j == N-1) { 346c4762a1bSJed Brown col[nc].i = i-1; col[nc].j = j; val[nc++] = -fwc/user->dx; 347c4762a1bSJed Brown col[nc].i = i; col[nc].j = j-1; val[nc++] = user->disper_coe/user->dy; 348c4762a1bSJed Brown col[nc].i = i; col[nc].j = j; val[nc++] = fwc/user->dx + fthetac - user->disper_coe/user->dy; 349c4762a1bSJed Brown } else if (i==0) { 350c4762a1bSJed Brown col[nc].i = i+1; col[nc].j = j; val[nc++] = fwc/user->dx; 351c4762a1bSJed Brown col[nc].i = i; col[nc].j = j+1; val[nc++] = -user->disper_coe/(2*user->dy); 352c4762a1bSJed Brown col[nc].i = i; col[nc].j = j-1; val[nc++] = user->disper_coe/(2*user->dy); 353c4762a1bSJed Brown col[nc].i = i; col[nc].j = j; val[nc++] = -fwc/user->dx + fthetac; 354c4762a1bSJed Brown } else if (i == M-1) { 355c4762a1bSJed Brown col[nc].i = i-1; col[nc].j = j; val[nc++] = -fwc/user->dx; 356c4762a1bSJed Brown col[nc].i = i; col[nc].j = j+1; val[nc++] = -user->disper_coe/(2*user->dy); 357c4762a1bSJed Brown col[nc].i = i; col[nc].j = j-1; val[nc++] = user->disper_coe/(2*user->dy); 358c4762a1bSJed Brown col[nc].i = i; col[nc].j = j; val[nc++] = fwc/user->dx + fthetac; 359c4762a1bSJed Brown } else if (j==0) { 360c4762a1bSJed Brown col[nc].i = i+1; col[nc].j = j; val[nc++] = fwc/(2*user->dx); 361c4762a1bSJed Brown col[nc].i = i-1; col[nc].j = j; val[nc++] = -fwc/(2*user->dx); 362c4762a1bSJed Brown col[nc].i = i; col[nc].j = j+1; val[nc++] = -user->disper_coe/user->dy; 363c4762a1bSJed Brown col[nc].i = i; col[nc].j = j; val[nc++] = user->disper_coe/user->dy + fthetac; 364c4762a1bSJed Brown } else if (j == N-1) { 365c4762a1bSJed Brown col[nc].i = i+1; col[nc].j = j; val[nc++] = fwc/(2*user->dx); 366c4762a1bSJed Brown col[nc].i = i-1; col[nc].j = j; val[nc++] = -fwc/(2*user->dx); 367c4762a1bSJed Brown col[nc].i = i; col[nc].j = j-1; val[nc++] = user->disper_coe/user->dy; 368c4762a1bSJed Brown col[nc].i = i; col[nc].j = j; val[nc++] = -user->disper_coe/user->dy + fthetac; 369c4762a1bSJed Brown } 370c4762a1bSJed Brown } 371c4762a1bSJed Brown } else { 372c4762a1bSJed Brown c1 = (yi-user->ws)/(2*user->dx); 373c4762a1bSJed Brown c3 = (user->ws/(2.0*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(xi))/(2*user->dy); 374c4762a1bSJed Brown c5 = (PetscPowScalar((user->lambda*user->ws)/(2*user->H),2)*user->q*(1.0-PetscExpScalar(-t/user->lambda)))/(user->dy*user->dy); 375c4762a1bSJed Brown col[nc].i = i-1; col[nc].j = j; val[nc++] = c1; 376c4762a1bSJed Brown col[nc].i = i+1; col[nc].j = j; val[nc++] = -c1; 377c4762a1bSJed Brown col[nc].i = i; col[nc].j = j-1; val[nc++] = c3 + c5; 378c4762a1bSJed Brown col[nc].i = i; col[nc].j = j+1; val[nc++] = -c3 + c5; 379c4762a1bSJed Brown col[nc].i = i; col[nc].j = j; val[nc++] = -2*c5 -a; 380c4762a1bSJed Brown } 3815f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValuesStencil(Jpre,1,&row,nc,col,val,INSERT_VALUES)); 382c4762a1bSJed Brown } 383c4762a1bSJed Brown } 3845f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArrayRead(cda,gc,&coors)); 385c4762a1bSJed Brown 3865f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY)); 3875f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY)); 388c4762a1bSJed Brown if (J != Jpre) { 3895f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 3905f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 391c4762a1bSJed Brown } 392c4762a1bSJed Brown PetscFunctionReturn(0); 393c4762a1bSJed Brown } 394c4762a1bSJed Brown 395c4762a1bSJed Brown PetscErrorCode Parameter_settings(AppCtx *user) 396c4762a1bSJed Brown { 397c4762a1bSJed Brown PetscBool flg; 398c4762a1bSJed Brown 399c4762a1bSJed Brown PetscFunctionBeginUser; 400c4762a1bSJed Brown 401c4762a1bSJed Brown /* Set default parameters */ 402c4762a1bSJed Brown user->ws = 1.0; 403c4762a1bSJed Brown user->H = 5.0; user->Pmax = 2.1; 404c4762a1bSJed Brown user->PM_min = 1.0; user->lambda = 0.1; 405c4762a1bSJed Brown user->q = 1.0; user->mux = PetscAsinScalar(user->PM_min/user->Pmax); 406c4762a1bSJed Brown user->sigmax = 0.1; 407c4762a1bSJed Brown user->sigmay = 0.1; user->rho = 0.0; 408c4762a1bSJed Brown user->t0 = 0.0; user->tmax = 2.0; 409c4762a1bSJed Brown user->xmin = -1.0; user->xmax = 10.0; 410c4762a1bSJed Brown user->ymin = -1.0; user->ymax = 10.0; 411c4762a1bSJed Brown user->bc = 0; 412c4762a1bSJed Brown 4135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-ws",&user->ws,&flg)); 4145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-Inertia",&user->H,&flg)); 4155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-Pmax",&user->Pmax,&flg)); 4165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-PM_min",&user->PM_min,&flg)); 4175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-lambda",&user->lambda,&flg)); 4185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-q",&user->q,&flg)); 4195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-mux",&user->mux,&flg)); 4205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-sigmax",&user->sigmax,&flg)); 4215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-muy",&user->muy,&flg)); 4225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-sigmay",&user->sigmay,&flg)); 4235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-rho",&user->rho,&flg)); 4245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-t0",&user->t0,&flg)); 4255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-tmax",&user->tmax,&flg)); 4265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-xmin",&user->xmin,&flg)); 4275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-xmax",&user->xmax,&flg)); 4285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-ymin",&user->ymin,&flg)); 4295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetScalar(NULL,NULL,"-ymax",&user->ymax,&flg)); 4305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-bc_type",&user->bc,&flg)); 431c4762a1bSJed Brown user->muy = user->ws; 432c4762a1bSJed Brown PetscFunctionReturn(0); 433c4762a1bSJed Brown } 434c4762a1bSJed Brown 435c4762a1bSJed Brown /*TEST 436c4762a1bSJed Brown 437c4762a1bSJed Brown build: 438c4762a1bSJed Brown requires: !complex 439c4762a1bSJed Brown 440c4762a1bSJed Brown test: 441c4762a1bSJed Brown args: -nox -ts_max_steps 2 442c4762a1bSJed Brown localrunfiles: petscopt_ex6 443c4762a1bSJed Brown timeoutfactor: 4 444c4762a1bSJed Brown 445c4762a1bSJed Brown TEST*/ 446