Actual source code: ex2.c
1: /*
2: Formatted test for TS routines.
4: Solves U_t=F(t,u)
5: Where:
6:
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>
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: MatStructure A_structure;
40: Mat A = 0;
42: PetscInitialize(&argc,&argv,(char*)0,help);
43: MPI_Comm_size(PETSC_COMM_WORLD,&size);
45: PetscOptionsGetInt(PETSC_NULL,"-time",&time_steps,PETSC_NULL);
47: /* set initial conditions */
48: VecCreate(PETSC_COMM_WORLD,&global);
49: VecSetSizes(global,PETSC_DECIDE,3);
50: VecSetFromOptions(global);
51: Initial(global,PETSC_NULL);
53: /* make timestep context */
54: TSCreate(PETSC_COMM_WORLD,&ts);
55: TSSetProblemType(ts,TS_NONLINEAR);
56: TSMonitorSet(ts,Monitor,PETSC_NULL,PETSC_NULL);
58: dt = 0.1;
60: /*
61: The user provides the RHS and Jacobian
62: */
63: TSSetRHSFunction(ts,PETSC_NULL,RHSFunction,NULL);
64: MatCreate(PETSC_COMM_WORLD,&A);
65: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,3,3);
66: MatSetFromOptions(A);
67: RHSJacobian(ts,0.0,global,&A,&A,&A_structure,NULL);
68: TSSetRHSJacobian(ts,A,A,RHSJacobian,NULL);
70: TSSetFromOptions(ts);
72: TSSetInitialTimeStep(ts,0.0,dt);
73: TSSetDuration(ts,time_steps,1);
74: TSSetSolution(ts,global);
76: TSSolve(ts,global,&ftime);
77: TSGetTimeStepNumber(ts,&steps);
80: /* free the memories */
82: TSDestroy(&ts);
83: VecDestroy(&global);
84: ierr= MatDestroy(&A);
86: PetscFinalize();
87: return 0;
88: }
90: /* -------------------------------------------------------------------*/
93: /* this test problem has initial values (1,1,1). */
94: PetscErrorCode Initial(Vec global,void *ctx)
95: {
96: PetscScalar *localptr;
97: PetscInt i,mybase,myend,locsize;
100: /* determine starting point of each processor */
101: VecGetOwnershipRange(global,&mybase,&myend);
102: VecGetLocalSize(global,&locsize);
104: /* Initialize the array */
105: VecGetArray(global,&localptr);
106: for (i=0; i<locsize; i++) {
107: localptr[i] = 1.0;
108: }
109:
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;
125: PetscScalar *tmp;
127: /* Get the size of the vector */
128: VecGetSize(global,&n);
130: /* Set the index sets */
131: PetscMalloc(n*sizeof(PetscInt),&idx);
132: for(i=0; i<n; i++) idx[i]=i;
133:
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: VecGetArray(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: VecRestoreArray(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 *inptr,*outptr;
163: PetscInt i,n,*idx;
165: IS from,to;
166: VecScatter scatter;
167: Vec tmp_in,tmp_out;
169: /* Get the length of parallel vector */
170: VecGetSize(globalin,&n);
172: /* Set the index sets */
173: PetscMalloc(n*sizeof(PetscInt),&idx);
174: for(i=0; i<n; i++) idx[i]=i;
175:
176: /* Create local sequential vectors */
177: VecCreateSeq(PETSC_COMM_SELF,n,&tmp_in);
178: VecDuplicate(tmp_in,&tmp_out);
180: /* Create scatter context */
181: ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from);
182: ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to);
183: VecScatterCreate(globalin,from,tmp_in,to,&scatter);
184: VecScatterBegin(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD);
185: VecScatterEnd(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD);
186: VecScatterDestroy(&scatter);
188: /*Extract income array */
189: VecGetArray(tmp_in,&inptr);
191: /* Extract outcome array*/
192: VecGetArray(tmp_out,&outptr);
194: outptr[0] = 2.0*inptr[0]+inptr[1];
195: outptr[1] = inptr[0]+2.0*inptr[1]+inptr[2];
196: outptr[2] = inptr[1]+2.0*inptr[2];
198: VecRestoreArray(tmp_in,&inptr);
199: VecRestoreArray(tmp_out,&outptr);
201: VecScatterCreate(tmp_out,from,globalout,to,&scatter);
202: VecScatterBegin(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD);
203: VecScatterEnd(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD);
205: /* Destroy idx aand scatter */
206: ISDestroy(&from);
207: ISDestroy(&to);
208: VecScatterDestroy(&scatter);
209: VecDestroy(&tmp_in);
210: VecDestroy(&tmp_out);
211: PetscFree(idx);
212: return 0;
213: }
217: PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
218: {
219: Mat A = *AA;
220: PetscScalar v[3],*tmp;
221: PetscInt idx[3],i;
223:
224: *str = SAME_NONZERO_PATTERN;
226: idx[0]=0; idx[1]=1; idx[2]=2;
227: VecGetArray(x,&tmp);
229: i = 0;
230: v[0] = 2.0; v[1] = 1.0; v[2] = 0.0;
231: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
233: i = 1;
234: v[0] = 1.0; v[1] = 2.0; v[2] = 1.0;
235: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
236:
237: i = 2;
238: v[0]= 0.0; v[1] = 1.0; v[2] = 2.0;
239: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
241: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
242: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
244: VecRestoreArray(x,&tmp);
246: return 0;
247: }
249: /*
250: The exact solutions
251: */
252: PetscReal solx(PetscReal t)
253: {
254: return exp((2.0 - sqrt(2.0))*t)/2.0 - exp((2.0 - sqrt(2.0))*t)/(2.0*sqrt(2.0)) +
255: exp((2.0 + sqrt(2.0))*t)/2.0 + exp((2.0 + sqrt(2.0))*t)/(2.0*sqrt(2.0));
256: }
258: PetscReal soly(PetscReal t)
259: {
260: return exp((2.0 - sqrt(2.0))*t)/2.0 - exp((2.0 - sqrt(2.0))*t)/sqrt(2.0) +
261: exp((2.0 + sqrt(2.0))*t)/2.0 + exp((2.0 + sqrt(2.0))*t)/sqrt(2.0);
262: }
263:
264: PetscReal solz(PetscReal t)
265: {
266: return exp((2.0 - sqrt(2.0))*t)/2.0 - exp((2.0 - sqrt(2.0))*t)/(2.0*sqrt(2.0)) +
267: exp((2.0 + sqrt(2.0))*t)/2.0 + exp((2.0 + sqrt(2.0))*t)/(2.0*sqrt(2.0));
268: }