Actual source code: ex2.c

petsc-3.8.3 2017-12-09
Report Typos and Errors
  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);

 29: int main(int argc,char **argv)
 30: {
 32:   PetscInt       time_steps = 100,steps;
 33:   PetscMPIInt    size;
 34:   Vec            global;
 35:   PetscReal      dt,ftime;
 36:   TS             ts;
 37:   Mat            A = 0;

 39:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 40:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 42:   PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL);

 44:   /* set initial conditions */
 45:   VecCreate(PETSC_COMM_WORLD,&global);
 46:   VecSetSizes(global,PETSC_DECIDE,3);
 47:   VecSetFromOptions(global);
 48:   Initial(global,NULL);

 50:   /* make timestep context */
 51:   TSCreate(PETSC_COMM_WORLD,&ts);
 52:   TSSetProblemType(ts,TS_NONLINEAR);
 53:   TSMonitorSet(ts,Monitor,NULL,NULL);

 55:   dt = 0.1;

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

 68:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
 69:   TSSetFromOptions(ts);

 71:   TSSetTimeStep(ts,dt);
 72:   TSSetMaxSteps(ts,time_steps);
 73:   TSSetMaxTime(ts,1);
 74:   TSSetSolution(ts,global);

 76:   TSSolve(ts,global);
 77:   TSGetSolveTime(ts,&ftime);
 78:   TSGetStepNumber(ts,&steps);


 81:   /* free the memories */

 83:   TSDestroy(&ts);
 84:   VecDestroy(&global);
 85:   MatDestroy(&A);

 87:   PetscFinalize();
 88:   return ierr;
 89: }

 91: /* -------------------------------------------------------------------*/
 92: /* this test problem has initial values (1,1,1).                      */
 93: PetscErrorCode Initial(Vec global,void *ctx)
 94: {
 95:   PetscScalar    *localptr;
 96:   PetscInt       i,mybase,myend,locsize;

 99:   /* determine starting point of each processor */
100:   VecGetOwnershipRange(global,&mybase,&myend);
101:   VecGetLocalSize(global,&locsize);

103:   /* Initialize the array */
104:   VecGetArray(global,&localptr);
105:   for (i=0; i<locsize; i++) localptr[i] = 1.0;

107:   if (mybase == 0) localptr[0]=1.0;

109:   VecRestoreArray(global,&localptr);
110:   return 0;
111: }

113: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec global,void *ctx)
114: {
115:   VecScatter        scatter;
116:   IS                from,to;
117:   PetscInt          i,n,*idx;
118:   Vec               tmp_vec;
119:   PetscErrorCode    ierr;
120:   const PetscScalar *tmp;

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

125:   /* Set the index sets */
126:   PetscMalloc1(n,&idx);
127:   for (i=0; i<n; i++) idx[i]=i;

129:   /* Create local sequential vectors */
130:   VecCreateSeq(PETSC_COMM_SELF,n,&tmp_vec);

132:   /* Create scatter context */
133:   ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from);
134:   ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to);
135:   VecScatterCreate(global,from,tmp_vec,to,&scatter);
136:   VecScatterBegin(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD);
137:   VecScatterEnd(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD);

139:   VecGetArrayRead(tmp_vec,&tmp);
140:   PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e u = %14.6e  %14.6e  %14.6e \n",
141:                      time,PetscRealPart(tmp[0]),PetscRealPart(tmp[1]),PetscRealPart(tmp[2]));
142:   PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e errors = %14.6e  %14.6e  %14.6e \n",
143:                      time,PetscRealPart(tmp[0]-solx(time)),PetscRealPart(tmp[1]-soly(time)),PetscRealPart(tmp[2]-solz(time)));
144:   VecRestoreArrayRead(tmp_vec,&tmp);
145:   VecScatterDestroy(&scatter);
146:   ISDestroy(&from);
147:   ISDestroy(&to);
148:   PetscFree(idx);
149:   VecDestroy(&tmp_vec);
150:   return 0;
151: }

153: PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
154: {
155:   PetscScalar       *outptr;
156:   const PetscScalar *inptr;
157:   PetscInt          i,n,*idx;
158:   PetscErrorCode    ierr;
159:   IS                from,to;
160:   VecScatter        scatter;
161:   Vec               tmp_in,tmp_out;

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

166:   /* Set the index sets */
167:   PetscMalloc1(n,&idx);
168:   for (i=0; i<n; i++) idx[i]=i;

170:   /* Create local sequential vectors */
171:   VecCreateSeq(PETSC_COMM_SELF,n,&tmp_in);
172:   VecDuplicate(tmp_in,&tmp_out);

174:   /* Create scatter context */
175:   ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from);
176:   ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to);
177:   VecScatterCreate(globalin,from,tmp_in,to,&scatter);
178:   VecScatterBegin(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD);
179:   VecScatterEnd(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD);
180:   VecScatterDestroy(&scatter);

182:   /*Extract income array */
183:   VecGetArrayRead(tmp_in,&inptr);

185:   /* Extract outcome array*/
186:   VecGetArray(tmp_out,&outptr);

188:   outptr[0] = 2.0*inptr[0]+inptr[1];
189:   outptr[1] = inptr[0]+2.0*inptr[1]+inptr[2];
190:   outptr[2] = inptr[1]+2.0*inptr[2];

192:   VecRestoreArrayRead(tmp_in,&inptr);
193:   VecRestoreArray(tmp_out,&outptr);

195:   VecScatterCreate(tmp_out,from,globalout,to,&scatter);
196:   VecScatterBegin(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD);
197:   VecScatterEnd(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD);

199:   /* Destroy idx aand scatter */
200:   ISDestroy(&from);
201:   ISDestroy(&to);
202:   VecScatterDestroy(&scatter);
203:   VecDestroy(&tmp_in);
204:   VecDestroy(&tmp_out);
205:   PetscFree(idx);
206:   return 0;
207: }

209: PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec x,Mat A,Mat BB,void *ctx)
210: {
211:   PetscScalar       v[3];
212:   const PetscScalar *tmp;
213:   PetscInt          idx[3],i;
214:   PetscErrorCode    ierr;

216:   idx[0]=0; idx[1]=1; idx[2]=2;
217:   VecGetArrayRead(x,&tmp);

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

223:   i    = 1;
224:   v[0] = 1.0; v[1] = 2.0; v[2] = 1.0;
225:   MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);

227:   i    = 2;
228:   v[0] = 0.0; v[1] = 1.0; v[2] = 2.0;
229:   MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);

231:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
232:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

234:   VecRestoreArrayRead(x,&tmp);

236:   return 0;
237: }

239: /*
240:       The exact solutions
241: */
242: PetscReal solx(PetscReal t)
243: {
244:   return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)) +
245:          PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0));
246: }

248: PetscReal soly(PetscReal t)
249: {
250:   return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/PetscSqrtReal(2.0) +
251:          PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/PetscSqrtReal(2.0);
252: }

254: PetscReal solz(PetscReal t)
255: {
256:   return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)) +
257:          PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0));
258: }