Actual source code: ex2.c
1: /*$Id: ex2.c,v 1.36 2001/08/07 21:31:27 bsmith Exp $*/
2: /*
3: Formatted test for TS routines.
5: Solves U_t=F(t,u)
6: Where:
7:
8: [2*u1+u2
9: F(t,u)= [u1+2*u2+u3
10: [ u2+2*u3
11: We can compare the solutions from euler, beuler and PVODE to
12: see what is the difference.
14: */
16: static char help[] = "Solves a nonlinear ODE. nn";
18: #include petscsys.h
19: #include petscts.h
20: #include petscpc.h
22: extern int RHSFunction(TS,PetscReal,Vec,Vec,void*);
23: extern int RHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure *,void*);
24: extern int Monitor(TS,int,PetscReal,Vec,void *);
25: extern int Initial(Vec,void *);
27: extern PetscReal solx(PetscReal);
28: extern PetscReal soly(PetscReal);
29: extern PetscReal solz(PetscReal);
31: int main(int argc,char **argv)
32: {
33: int ierr,time_steps = 100,steps,size;
34: Vec global;
35: PetscReal dt,ftime;
36: TS ts;
37: MatStructure A_structure;
38: Mat A = 0;
39:
40: PetscInitialize(&argc,&argv,(char*)0,help);
41: MPI_Comm_size(PETSC_COMM_WORLD,&size);
42:
43: PetscOptionsGetInt(PETSC_NULL,"-time",&time_steps,PETSC_NULL);
44:
45: /* set initial conditions */
46: VecCreate(PETSC_COMM_WORLD,&global);
47: VecSetSizes(global,PETSC_DECIDE,3);
48: VecSetFromOptions(global);
49: Initial(global,PETSC_NULL);
50:
51: /* make timestep context */
52: TSCreate(PETSC_COMM_WORLD,&ts);
53: TSSetProblemType(ts,TS_NONLINEAR);
54: TSSetMonitor(ts,Monitor,PETSC_NULL,PETSC_NULL);
56: dt = 0.1;
58: /*
59: The user provides the RHS and Jacobian
60: */
61: TSSetRHSFunction(ts,RHSFunction,NULL);
62: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,3,3,&A);
63: MatSetFromOptions(A);
64: RHSJacobian(ts,0.0,global,&A,&A,&A_structure,NULL);
65: TSSetRHSJacobian(ts,A,A,RHSJacobian,NULL);
66:
67: TSSetFromOptions(ts);
69: TSSetInitialTimeStep(ts,0.0,dt);
70: TSSetDuration(ts,time_steps,1);
71: TSSetSolution(ts,global);
74: TSSetUp(ts);
75: TSStep(ts,&steps,&ftime);
78: /* free the memories */
79:
80: TSDestroy(ts);
81: VecDestroy(global);
82: ierr= MatDestroy(A);
84: PetscFinalize();
85: return 0;
86: }
88: /* -------------------------------------------------------------------*/
89: /* this test problem has initial values (1,1,1). */
90: int Initial(Vec global,void *ctx)
91: {
92: PetscScalar *localptr;
93: int i,mybase,myend,ierr,locsize;
95: /* determine starting point of each processor */
96: VecGetOwnershipRange(global,&mybase,&myend);
97: VecGetLocalSize(global,&locsize);
99: /* Initialize the array */
100: VecGetArray(global,&localptr);
101: for (i=0; i<locsize; i++) {
102: localptr[i] = 1.0;
103: }
104:
105: if (mybase == 0) localptr[0]=1.0;
107: VecRestoreArray(global,&localptr);
108: return 0;
109: }
111: int Monitor(TS ts,int step,PetscReal time,Vec global,void *ctx)
112: {
113: VecScatter scatter;
114: IS from,to;
115: int i,n,*idx;
116: Vec tmp_vec;
117: int ierr;
118: PetscScalar *tmp;
120: /* Get the size of the vector */
121: VecGetSize(global,&n);
123: /* Set the index sets */
124: PetscMalloc(n*sizeof(int),&idx);
125: for(i=0; i<n; i++) idx[i]=i;
126:
127: /* Create local sequential vectors */
128: VecCreateSeq(PETSC_COMM_SELF,n,&tmp_vec);
130: /* Create scatter context */
131: ISCreateGeneral(PETSC_COMM_SELF,n,idx,&from);
132: ISCreateGeneral(PETSC_COMM_SELF,n,idx,&to);
133: VecScatterCreate(global,from,tmp_vec,to,&scatter);
134: VecScatterBegin(global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD,scatter);
135: VecScatterEnd(global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD,scatter);
137: VecGetArray(tmp_vec,&tmp);
138: PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e u = %14.6e %14.6e %14.6e n",
139: time,PetscRealPart(tmp[0]),PetscRealPart(tmp[1]),PetscRealPart(tmp[2]));
140: PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e errors = %14.6e %14.6e %14.6e n",
141: time,PetscRealPart(tmp[0]-solx(time)),PetscRealPart(tmp[1]-soly(time)),PetscRealPart(tmp[2]-solz(time)));
142: VecRestoreArray(tmp_vec,&tmp);
143: VecScatterDestroy(scatter);
144: ISDestroy(from);
145: ISDestroy(to);
146: PetscFree(idx);
147: VecDestroy(tmp_vec);
148: return 0;
149: }
151: int RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
152: {
153: PetscScalar *inptr,*outptr;
154: int i,n,ierr,*idx;
155: IS from,to;
156: VecScatter scatter;
157: Vec tmp_in,tmp_out;
159: /* Get the length of parallel vector */
160: VecGetSize(globalin,&n);
162: /* Set the index sets */
163: PetscMalloc(n*sizeof(int),&idx);
164: for(i=0; i<n; i++) idx[i]=i;
165:
166: /* Create local sequential vectors */
167: VecCreateSeq(PETSC_COMM_SELF,n,&tmp_in);
168: VecDuplicate(tmp_in,&tmp_out);
170: /* Create scatter context */
171: ISCreateGeneral(PETSC_COMM_SELF,n,idx,&from);
172: ISCreateGeneral(PETSC_COMM_SELF,n,idx,&to);
173: VecScatterCreate(globalin,from,tmp_in,to,&scatter);
174: VecScatterBegin(globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD,scatter);
175: VecScatterEnd(globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD,scatter);
176: VecScatterDestroy(scatter);
178: /*Extract income array */
179: VecGetArray(tmp_in,&inptr);
181: /* Extract outcome array*/
182: VecGetArray(tmp_out,&outptr);
184: outptr[0] = 2.0*inptr[0]+inptr[1];
185: outptr[1] = inptr[0]+2.0*inptr[1]+inptr[2];
186: outptr[2] = inptr[1]+2.0*inptr[2];
188: VecRestoreArray(tmp_in,&inptr);
189: VecRestoreArray(tmp_out,&outptr);
191: VecScatterCreate(tmp_out,from,globalout,to,&scatter);
192: VecScatterBegin(tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD,scatter);
193: VecScatterEnd(tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD,scatter);
195: /* Destroy idx aand scatter */
196: ISDestroy(from);
197: ISDestroy(to);
198: VecScatterDestroy(scatter);
199: VecDestroy(tmp_in);
200: VecDestroy(tmp_out);
201: PetscFree(idx);
202: return 0;
203: }
205: int RHSJacobian(TS ts,PetscReal t,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
206: {
207: Mat A = *AA;
208: PetscScalar v[3],*tmp;
209: int idx[3],i,ierr;
210:
211: *str = SAME_NONZERO_PATTERN;
213: idx[0]=0; idx[1]=1; idx[2]=2;
214: VecGetArray(x,&tmp);
216: i = 0;
217: v[0] = 2.0; v[1] = 1.0; v[2] = 0.0;
218: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
220: i = 1;
221: v[0] = 1.0; v[1] = 2.0; v[2] = 1.0;
222: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
223:
224: i = 2;
225: v[0]= 0.0; v[1] = 1.0; v[2] = 2.0;
226: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
228: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
229: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
231: VecRestoreArray(x,&tmp);
233: return 0;
234: }
236: /*
237: The exact solutions
238: */
239: PetscReal solx(PetscReal t)
240: {
241: return exp((2.0 - sqrt(2.0))*t)/2.0 - exp((2.0 - sqrt(2.0))*t)/(2.0*sqrt(2.0)) +
242: exp((2.0 + sqrt(2.0))*t)/2.0 + exp((2.0 + sqrt(2.0))*t)/(2.0*sqrt(2.0));
243: }
245: PetscReal soly(PetscReal t)
246: {
247: return exp((2.0 - sqrt(2.0))*t)/2.0 - exp((2.0 - sqrt(2.0))*t)/sqrt(2.0) +
248: exp((2.0 + sqrt(2.0))*t)/2.0 + exp((2.0 + sqrt(2.0))*t)/sqrt(2.0);
249: }
250:
251: PetscReal solz(PetscReal t)
252: {
253: return exp((2.0 - sqrt(2.0))*t)/2.0 - exp((2.0 - sqrt(2.0))*t)/(2.0*sqrt(2.0)) +
254: exp((2.0 + sqrt(2.0))*t)/2.0 + exp((2.0 + sqrt(2.0))*t)/(2.0*sqrt(2.0));
255: }