Actual source code: ex2.c
petsc-3.7.3 2016-08-01
1: /*
2: Formatted test for TS routines.
4: Solves U_t=F(t,u)
5: Where:
7: [2*u1+u2
8: F(t,u)= [u1+2*u2+u3
9: [ u2+2*u3
10: We can compare the solutions from euler, beuler and SUNDIALS to
11: see what is the difference.
13: */
15: static char help[] = "Solves a nonlinear ODE. \n\n";
17: #include <petscts.h>
18: #include <petscpc.h>
20: extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
21: extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
22: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
23: extern PetscErrorCode Initial(Vec,void*);
25: extern PetscReal solx(PetscReal);
26: extern PetscReal soly(PetscReal);
27: extern PetscReal solz(PetscReal);
31: int main(int argc,char **argv)
32: {
34: PetscInt time_steps = 100,steps;
35: PetscMPIInt size;
36: Vec global;
37: PetscReal dt,ftime;
38: TS ts;
39: Mat A = 0;
41: PetscInitialize(&argc,&argv,(char*)0,help);
42: MPI_Comm_size(PETSC_COMM_WORLD,&size);
44: PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL);
46: /* set initial conditions */
47: VecCreate(PETSC_COMM_WORLD,&global);
48: VecSetSizes(global,PETSC_DECIDE,3);
49: VecSetFromOptions(global);
50: Initial(global,NULL);
52: /* make timestep context */
53: TSCreate(PETSC_COMM_WORLD,&ts);
54: TSSetProblemType(ts,TS_NONLINEAR);
55: TSMonitorSet(ts,Monitor,NULL,NULL);
57: dt = 0.1;
59: /*
60: The user provides the RHS and Jacobian
61: */
62: TSSetRHSFunction(ts,NULL,RHSFunction,NULL);
63: MatCreate(PETSC_COMM_WORLD,&A);
64: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,3,3);
65: MatSetFromOptions(A);
66: MatSetUp(A);
67: RHSJacobian(ts,0.0,global,A,A,NULL);
68: TSSetRHSJacobian(ts,A,A,RHSJacobian,NULL);
70: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
71: TSSetFromOptions(ts);
73: TSSetInitialTimeStep(ts,0.0,dt);
74: TSSetDuration(ts,time_steps,1);
75: TSSetSolution(ts,global);
77: TSSolve(ts,global);
78: TSGetSolveTime(ts,&ftime);
79: TSGetTimeStepNumber(ts,&steps);
82: /* free the memories */
84: TSDestroy(&ts);
85: VecDestroy(&global);
86: MatDestroy(&A);
88: PetscFinalize();
89: return 0;
90: }
92: /* -------------------------------------------------------------------*/
95: /* this test problem has initial values (1,1,1). */
96: PetscErrorCode Initial(Vec global,void *ctx)
97: {
98: PetscScalar *localptr;
99: PetscInt i,mybase,myend,locsize;
102: /* determine starting point of each processor */
103: VecGetOwnershipRange(global,&mybase,&myend);
104: VecGetLocalSize(global,&locsize);
106: /* Initialize the array */
107: VecGetArray(global,&localptr);
108: for (i=0; i<locsize; i++) localptr[i] = 1.0;
110: if (mybase == 0) localptr[0]=1.0;
112: VecRestoreArray(global,&localptr);
113: return 0;
114: }
118: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec global,void *ctx)
119: {
120: VecScatter scatter;
121: IS from,to;
122: PetscInt i,n,*idx;
123: Vec tmp_vec;
124: PetscErrorCode ierr;
125: const PetscScalar *tmp;
127: /* Get the size of the vector */
128: VecGetSize(global,&n);
130: /* Set the index sets */
131: PetscMalloc1(n,&idx);
132: for (i=0; i<n; i++) idx[i]=i;
134: /* Create local sequential vectors */
135: VecCreateSeq(PETSC_COMM_SELF,n,&tmp_vec);
137: /* Create scatter context */
138: ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from);
139: ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to);
140: VecScatterCreate(global,from,tmp_vec,to,&scatter);
141: VecScatterBegin(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD);
142: VecScatterEnd(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD);
144: VecGetArrayRead(tmp_vec,&tmp);
145: PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e u = %14.6e %14.6e %14.6e \n",
146: time,PetscRealPart(tmp[0]),PetscRealPart(tmp[1]),PetscRealPart(tmp[2]));
147: PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e errors = %14.6e %14.6e %14.6e \n",
148: time,PetscRealPart(tmp[0]-solx(time)),PetscRealPart(tmp[1]-soly(time)),PetscRealPart(tmp[2]-solz(time)));
149: VecRestoreArrayRead(tmp_vec,&tmp);
150: VecScatterDestroy(&scatter);
151: ISDestroy(&from);
152: ISDestroy(&to);
153: PetscFree(idx);
154: VecDestroy(&tmp_vec);
155: return 0;
156: }
160: PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
161: {
162: PetscScalar *outptr;
163: const PetscScalar *inptr;
164: PetscInt i,n,*idx;
165: PetscErrorCode ierr;
166: IS from,to;
167: VecScatter scatter;
168: Vec tmp_in,tmp_out;
170: /* Get the length of parallel vector */
171: VecGetSize(globalin,&n);
173: /* Set the index sets */
174: PetscMalloc1(n,&idx);
175: for (i=0; i<n; i++) idx[i]=i;
177: /* Create local sequential vectors */
178: VecCreateSeq(PETSC_COMM_SELF,n,&tmp_in);
179: VecDuplicate(tmp_in,&tmp_out);
181: /* Create scatter context */
182: ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from);
183: ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to);
184: VecScatterCreate(globalin,from,tmp_in,to,&scatter);
185: VecScatterBegin(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD);
186: VecScatterEnd(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD);
187: VecScatterDestroy(&scatter);
189: /*Extract income array */
190: VecGetArrayRead(tmp_in,&inptr);
192: /* Extract outcome array*/
193: VecGetArray(tmp_out,&outptr);
195: outptr[0] = 2.0*inptr[0]+inptr[1];
196: outptr[1] = inptr[0]+2.0*inptr[1]+inptr[2];
197: outptr[2] = inptr[1]+2.0*inptr[2];
199: VecRestoreArrayRead(tmp_in,&inptr);
200: VecRestoreArray(tmp_out,&outptr);
202: VecScatterCreate(tmp_out,from,globalout,to,&scatter);
203: VecScatterBegin(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD);
204: VecScatterEnd(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD);
206: /* Destroy idx aand scatter */
207: ISDestroy(&from);
208: ISDestroy(&to);
209: VecScatterDestroy(&scatter);
210: VecDestroy(&tmp_in);
211: VecDestroy(&tmp_out);
212: PetscFree(idx);
213: return 0;
214: }
218: PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec x,Mat A,Mat BB,void *ctx)
219: {
220: PetscScalar v[3];
221: const PetscScalar *tmp;
222: PetscInt idx[3],i;
223: PetscErrorCode ierr;
225: idx[0]=0; idx[1]=1; idx[2]=2;
226: VecGetArrayRead(x,&tmp);
228: i = 0;
229: v[0] = 2.0; v[1] = 1.0; v[2] = 0.0;
230: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
232: i = 1;
233: v[0] = 1.0; v[1] = 2.0; v[2] = 1.0;
234: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
236: i = 2;
237: v[0] = 0.0; v[1] = 1.0; v[2] = 2.0;
238: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
240: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
241: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
243: VecRestoreArrayRead(x,&tmp);
245: return 0;
246: }
248: /*
249: The exact solutions
250: */
251: PetscReal solx(PetscReal t)
252: {
253: return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)) +
254: PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0));
255: }
257: PetscReal soly(PetscReal t)
258: {
259: return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/PetscSqrtReal(2.0) +
260: PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/PetscSqrtReal(2.0);
261: }
263: PetscReal solz(PetscReal t)
264: {
265: return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)) +
266: PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0));
267: }