Actual source code: ex2.c

  1: /*$Id: ex2.c,v 1.34 2001/04/10 19:37:11 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,double,Vec,Vec,void*);
 23: extern int RHSJacobian(TS,double,Vec,Mat*,Mat*,MatStructure *,void*);
 24: extern int Monitor(TS,int,double,Vec,void *);
 25: extern int Initial(Vec,void *);

 27: extern double solx(double);
 28: extern double soly(double);
 29: extern double solz(double);

 31: int main(int argc,char **argv)
 32: {
 33:   int           ierr,time_steps = 100,steps,size;
 34:   Vec           global;
 35:   double        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,PETSC_DECIDE,3,&global);
 47:   VecSetFromOptions(global);
 48:   Initial(global,PETSC_NULL);
 49: 
 50:   /* make timestep context */
 51:   TSCreate(PETSC_COMM_WORLD,TS_NONLINEAR,&ts);
 52:   TSSetMonitor(ts,Monitor,PETSC_NULL,PETSC_NULL);

 54:   dt = 0.1;

 56:   /*
 57:     The user provides the RHS and Jacobian
 58:   */
 59:   TSSetRHSFunction(ts,RHSFunction,NULL);
 60:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,3,3,&A);
 61:   MatSetFromOptions(A);
 62:   RHSJacobian(ts,0.0,global,&A,&A,&A_structure,NULL);
 63:   TSSetRHSJacobian(ts,A,A,RHSJacobian,NULL);
 64: 
 65:   TSSetFromOptions(ts);

 67:   TSSetInitialTimeStep(ts,0.0,dt);
 68:   TSSetDuration(ts,time_steps,1);
 69:   TSSetSolution(ts,global);


 72:   TSSetUp(ts);
 73:   TSStep(ts,&steps,&ftime);


 76:   /* free the memories */
 77: 
 78:   TSDestroy(ts);
 79:   VecDestroy(global);
 80:   ierr= MatDestroy(A);

 82:   PetscFinalize();
 83:   return 0;
 84: }

 86: /* -------------------------------------------------------------------*/
 87: /* this test problem has initial values (1,1,1).                      */
 88: int Initial(Vec global,void *ctx)
 89: {
 90:   Scalar *localptr;
 91:   int    i,mybase,myend,ierr,locsize;

 93:   /* determine starting point of each processor */
 94:   VecGetOwnershipRange(global,&mybase,&myend);
 95:   VecGetLocalSize(global,&locsize);

 97:   /* Initialize the array */
 98:   VecGetArray(global,&localptr);
 99:   for (i=0; i<locsize; i++) {
100:     localptr[i] = 1.0;
101:   }
102: 
103:   if (mybase == 0) localptr[0]=1.0;

105:   VecRestoreArray(global,&localptr);
106:   return 0;
107: }

109: int Monitor(TS ts,int step,double time,Vec global,void *ctx)
110: {
111:   VecScatter scatter;
112:   IS         from,to;
113:   int        i,n,*idx;
114:   Vec        tmp_vec;
115:   int        ierr;
116:   Scalar     *tmp;

118:   /* Get the size of the vector */
119:   VecGetSize(global,&n);

121:   /* Set the index sets */
122:   PetscMalloc(n*sizeof(int),&idx);
123:   for(i=0; i<n; i++) idx[i]=i;
124: 
125:   /* Create local sequential vectors */
126:   VecCreateSeq(PETSC_COMM_SELF,n,&tmp_vec);

128:   /* Create scatter context */
129:   ISCreateGeneral(PETSC_COMM_SELF,n,idx,&from);
130:   ISCreateGeneral(PETSC_COMM_SELF,n,idx,&to);
131:   VecScatterCreate(global,from,tmp_vec,to,&scatter);
132:   VecScatterBegin(global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD,scatter);
133:   VecScatterEnd(global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD,scatter);

135:   VecGetArray(tmp_vec,&tmp);
136:   PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e u = %14.6e  %14.6e  %14.6e n",
137:                      time,PetscRealPart(tmp[0]),PetscRealPart(tmp[1]),PetscRealPart(tmp[2]));
138:   PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e errors = %14.6e  %14.6e  %14.6e n",
139:                      time,PetscRealPart(tmp[0]-solx(time)),PetscRealPart(tmp[1]-soly(time)),PetscRealPart(tmp[2]-solz(time)));
140:   VecRestoreArray(tmp_vec,&tmp);
141:   VecScatterDestroy(scatter);
142:   ISDestroy(from);
143:   ISDestroy(to);
144:   PetscFree(idx);
145:   VecDestroy(tmp_vec);
146:   return 0;
147: }

149: int RHSFunction(TS ts,double t,Vec globalin,Vec globalout,void *ctx)
150: {
151:   Scalar     *inptr,*outptr;
152:   int        i,n,ierr,*idx;
153:   IS         from,to;
154:   VecScatter scatter;
155:   Vec        tmp_in,tmp_out;

157:   /* Get the length of parallel vector */
158:   VecGetSize(globalin,&n);

160:   /* Set the index sets */
161:   PetscMalloc(n*sizeof(int),&idx);
162:   for(i=0; i<n; i++) idx[i]=i;
163: 
164:   /* Create local sequential vectors */
165:   VecCreateSeq(PETSC_COMM_SELF,n,&tmp_in);
166:   VecDuplicate(tmp_in,&tmp_out);

168:   /* Create scatter context */
169:   ISCreateGeneral(PETSC_COMM_SELF,n,idx,&from);
170:   ISCreateGeneral(PETSC_COMM_SELF,n,idx,&to);
171:   VecScatterCreate(globalin,from,tmp_in,to,&scatter);
172:   VecScatterBegin(globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD,scatter);
173:   VecScatterEnd(globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD,scatter);
174:   VecScatterDestroy(scatter);

176:   /*Extract income array */
177:   VecGetArray(tmp_in,&inptr);

179:   /* Extract outcome array*/
180:   VecGetArray(tmp_out,&outptr);

182:   outptr[0] = 2.0*inptr[0]+inptr[1];
183:   outptr[1] = inptr[0]+2.0*inptr[1]+inptr[2];
184:   outptr[2] = inptr[1]+2.0*inptr[2];

186:   VecRestoreArray(tmp_in,&inptr);
187:   VecRestoreArray(tmp_out,&outptr);

189:   VecScatterCreate(tmp_out,from,globalout,to,&scatter);
190:   VecScatterBegin(tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD,scatter);
191:   VecScatterEnd(tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD,scatter);

193:   /* Destroy idx aand scatter */
194:   ISDestroy(from);
195:   ISDestroy(to);
196:   VecScatterDestroy(scatter);
197:   VecDestroy(tmp_in);
198:   VecDestroy(tmp_out);
199:   PetscFree(idx);
200:   return 0;
201: }

203: int RHSJacobian(TS ts,double t,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
204: {
205:   Mat    A = *AA;
206:   Scalar v[3],*tmp;
207:   int    idx[3],i,ierr;
208: 
209:   *str = SAME_NONZERO_PATTERN;

211:   idx[0]=0; idx[1]=1; idx[2]=2;
212:   VecGetArray(x,&tmp);

214:   i = 0;
215:   v[0] = 2.0; v[1] = 1.0; v[2] = 0.0;
216:   MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);

218:   i = 1;
219:   v[0] = 1.0; v[1] = 2.0; v[2] = 1.0;
220:   MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
221: 
222:   i = 2;
223:   v[0]= 0.0; v[1] = 1.0; v[2] = 2.0;
224:   MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);

226:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
227:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

229:   VecRestoreArray(x,&tmp);

231:   return 0;
232: }

234: /*
235:       The exact solutions 
236: */
237: double solx(double t)
238: {
239:   return exp((2.0 - sqrt(2.0))*t)/2.0 - exp((2.0 - sqrt(2.0))*t)/(2.0*sqrt(2.0)) +
240:          exp((2.0 + sqrt(2.0))*t)/2.0 + exp((2.0 + sqrt(2.0))*t)/(2.0*sqrt(2.0));
241: }

243: double soly(double t)
244: {
245:   return exp((2.0 - sqrt(2.0))*t)/2.0 - exp((2.0 - sqrt(2.0))*t)/sqrt(2.0) +
246:          exp((2.0 + sqrt(2.0))*t)/2.0 + exp((2.0 + sqrt(2.0))*t)/sqrt(2.0);
247: }
248: 
249: double solz(double t)
250: {
251:   return exp((2.0 - sqrt(2.0))*t)/2.0 - exp((2.0 - sqrt(2.0))*t)/(2.0*sqrt(2.0)) +
252:          exp((2.0 + sqrt(2.0))*t)/2.0 + exp((2.0 + sqrt(2.0))*t)/(2.0*sqrt(2.0));
253: }