1c4762a1bSJed Brown /* 2c4762a1bSJed Brown Formatted test for TS routines. 3c4762a1bSJed Brown 4c4762a1bSJed Brown Solves U_t=F(t,u) 5c4762a1bSJed Brown Where: 6c4762a1bSJed Brown 7c4762a1bSJed Brown [2*u1+u2 8c4762a1bSJed Brown F(t,u)= [u1+2*u2+u3 9c4762a1bSJed Brown [ u2+2*u3 10c4762a1bSJed Brown We can compare the solutions from euler, beuler and SUNDIALS to 11c4762a1bSJed Brown see what is the difference. 12c4762a1bSJed Brown 13c4762a1bSJed Brown */ 14c4762a1bSJed Brown 15c4762a1bSJed Brown static char help[] = "Solves a linear ODE. \n\n"; 16c4762a1bSJed Brown 17c4762a1bSJed Brown #include <petscts.h> 18c4762a1bSJed Brown #include <petscpc.h> 19c4762a1bSJed Brown 20c4762a1bSJed Brown extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*); 21c4762a1bSJed Brown extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*); 22c4762a1bSJed Brown extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*); 23c4762a1bSJed Brown extern PetscErrorCode Initial(Vec,void*); 24c4762a1bSJed Brown extern PetscErrorCode MyMatMult(Mat,Vec,Vec); 25c4762a1bSJed Brown 26c4762a1bSJed Brown extern PetscReal solx(PetscReal); 27c4762a1bSJed Brown extern PetscReal soly(PetscReal); 28c4762a1bSJed Brown extern PetscReal solz(PetscReal); 29c4762a1bSJed Brown 30c4762a1bSJed Brown int main(int argc,char **argv) 31c4762a1bSJed Brown { 32c4762a1bSJed Brown PetscInt time_steps = 100,steps; 33c4762a1bSJed Brown Vec global; 34c4762a1bSJed Brown PetscReal dt,ftime; 35c4762a1bSJed Brown TS ts; 36c4762a1bSJed Brown Mat A = 0,S; 37c4762a1bSJed Brown 389566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 399566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL)); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* set initial conditions */ 429566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_WORLD,&global)); 439566063dSJacob Faibussowitsch PetscCall(VecSetSizes(global,PETSC_DECIDE,3)); 449566063dSJacob Faibussowitsch PetscCall(VecSetFromOptions(global)); 459566063dSJacob Faibussowitsch PetscCall(Initial(global,NULL)); 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* make timestep context */ 489566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 499566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts,TS_NONLINEAR)); 509566063dSJacob Faibussowitsch PetscCall(TSMonitorSet(ts,Monitor,NULL,NULL)); 51c4762a1bSJed Brown dt = 0.001; 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* 54c4762a1bSJed Brown The user provides the RHS and Jacobian 55c4762a1bSJed Brown */ 569566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts,NULL,RHSFunction,NULL)); 579566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 589566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,3,3)); 599566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 609566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 619566063dSJacob Faibussowitsch PetscCall(RHSJacobian(ts,0.0,global,A,A,NULL)); 629566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts,A,A,RHSJacobian,NULL)); 63c4762a1bSJed Brown 649566063dSJacob Faibussowitsch PetscCall(MatCreateShell(PETSC_COMM_WORLD,3,3,3,3,NULL,&S)); 659566063dSJacob Faibussowitsch PetscCall(MatShellSetOperation(S,MATOP_MULT,(void (*)(void))MyMatMult)); 669566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts,S,A,RHSJacobian,NULL)); 67c4762a1bSJed Brown 689566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 699566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 70c4762a1bSJed Brown 719566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,dt)); 729566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts,time_steps)); 739566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,1)); 749566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,global)); 75c4762a1bSJed Brown 769566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,global)); 779566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&ftime)); 789566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&steps)); 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* free the memories */ 81c4762a1bSJed Brown 829566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&global)); 849566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&S)); 86c4762a1bSJed Brown 879566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 88b122ec5aSJacob Faibussowitsch return 0; 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat S,Vec x,Vec y) 92c4762a1bSJed Brown { 93c4762a1bSJed Brown const PetscScalar *inptr; 94c4762a1bSJed Brown PetscScalar *outptr; 95c4762a1bSJed Brown 96c4762a1bSJed Brown PetscFunctionBegin; 979566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&inptr)); 989566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(y,&outptr)); 99c4762a1bSJed Brown 100c4762a1bSJed Brown outptr[0] = 2.0*inptr[0]+inptr[1]; 101c4762a1bSJed Brown outptr[1] = inptr[0]+2.0*inptr[1]+inptr[2]; 102c4762a1bSJed Brown outptr[2] = inptr[1]+2.0*inptr[2]; 103c4762a1bSJed Brown 1049566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&inptr)); 1059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(y,&outptr)); 106c4762a1bSJed Brown PetscFunctionReturn(0); 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 109c4762a1bSJed Brown /* this test problem has initial values (1,1,1). */ 110c4762a1bSJed Brown PetscErrorCode Initial(Vec global,void *ctx) 111c4762a1bSJed Brown { 112c4762a1bSJed Brown PetscScalar *localptr; 113c4762a1bSJed Brown PetscInt i,mybase,myend,locsize; 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* determine starting point of each processor */ 1169566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(global,&mybase,&myend)); 1179566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(global,&locsize)); 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* Initialize the array */ 1209566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(global,&localptr)); 121c4762a1bSJed Brown for (i=0; i<locsize; i++) localptr[i] = 1.0; 122c4762a1bSJed Brown 123c4762a1bSJed Brown if (mybase == 0) localptr[0]=1.0; 124c4762a1bSJed Brown 1259566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(global,&localptr)); 126c4762a1bSJed Brown return 0; 127c4762a1bSJed Brown } 128c4762a1bSJed Brown 129c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec global,void *ctx) 130c4762a1bSJed Brown { 131c4762a1bSJed Brown VecScatter scatter; 132c4762a1bSJed Brown IS from,to; 133c4762a1bSJed Brown PetscInt i,n,*idx; 134c4762a1bSJed Brown Vec tmp_vec; 135c4762a1bSJed Brown const PetscScalar *tmp; 136c4762a1bSJed Brown 137c4762a1bSJed Brown /* Get the size of the vector */ 1389566063dSJacob Faibussowitsch PetscCall(VecGetSize(global,&n)); 139c4762a1bSJed Brown 140c4762a1bSJed Brown /* Set the index sets */ 1419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&idx)); 142c4762a1bSJed Brown for (i=0; i<n; i++) idx[i]=i; 143c4762a1bSJed Brown 144c4762a1bSJed Brown /* Create local sequential vectors */ 1459566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&tmp_vec)); 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* Create scatter context */ 1489566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from)); 1499566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to)); 1509566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(global,from,tmp_vec,to,&scatter)); 1519566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD)); 1529566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD)); 153c4762a1bSJed Brown 1549566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tmp_vec,&tmp)); 155*d0609cedSBarry Smith PetscCall(PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e u = %14.6e %14.6e %14.6e \n", 156*d0609cedSBarry Smith (double)time,(double)PetscRealPart(tmp[0]),(double)PetscRealPart(tmp[1]),(double)PetscRealPart(tmp[2]))); 157*d0609cedSBarry Smith PetscCall(PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e errors = %14.6e %14.6e %14.6e \n", 158*d0609cedSBarry Smith (double)time,(double)PetscRealPart(tmp[0]-solx(time)),(double)PetscRealPart(tmp[1]-soly(time)),(double)PetscRealPart(tmp[2]-solz(time)))); 1599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tmp_vec,&tmp)); 1609566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&scatter)); 1619566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 1629566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 1639566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 1649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tmp_vec)); 165c4762a1bSJed Brown return 0; 166c4762a1bSJed Brown } 167c4762a1bSJed Brown 168c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx) 169c4762a1bSJed Brown { 170c4762a1bSJed Brown PetscScalar *outptr; 171c4762a1bSJed Brown const PetscScalar *inptr; 172c4762a1bSJed Brown PetscInt i,n,*idx; 173c4762a1bSJed Brown IS from,to; 174c4762a1bSJed Brown VecScatter scatter; 175c4762a1bSJed Brown Vec tmp_in,tmp_out; 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* Get the length of parallel vector */ 1789566063dSJacob Faibussowitsch PetscCall(VecGetSize(globalin,&n)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* Set the index sets */ 1819566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&idx)); 182c4762a1bSJed Brown for (i=0; i<n; i++) idx[i]=i; 183c4762a1bSJed Brown 184c4762a1bSJed Brown /* Create local sequential vectors */ 1859566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,n,&tmp_in)); 1869566063dSJacob Faibussowitsch PetscCall(VecDuplicate(tmp_in,&tmp_out)); 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* Create scatter context */ 1899566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from)); 1909566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to)); 1919566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(globalin,from,tmp_in,to,&scatter)); 1929566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD)); 1939566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD)); 1949566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&scatter)); 195c4762a1bSJed Brown 196c4762a1bSJed Brown /*Extract income array */ 1979566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(tmp_in,&inptr)); 198c4762a1bSJed Brown 199c4762a1bSJed Brown /* Extract outcome array*/ 2009566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(tmp_out,&outptr)); 201c4762a1bSJed Brown 202c4762a1bSJed Brown outptr[0] = 2.0*inptr[0]+inptr[1]; 203c4762a1bSJed Brown outptr[1] = inptr[0]+2.0*inptr[1]+inptr[2]; 204c4762a1bSJed Brown outptr[2] = inptr[1]+2.0*inptr[2]; 205c4762a1bSJed Brown 2069566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(tmp_in,&inptr)); 2079566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(tmp_out,&outptr)); 208c4762a1bSJed Brown 2099566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(tmp_out,from,globalout,to,&scatter)); 2109566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD)); 2119566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD)); 212c4762a1bSJed Brown 213c4762a1bSJed Brown /* Destroy idx aand scatter */ 2149566063dSJacob Faibussowitsch PetscCall(ISDestroy(&from)); 2159566063dSJacob Faibussowitsch PetscCall(ISDestroy(&to)); 2169566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&scatter)); 2179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tmp_in)); 2189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&tmp_out)); 2199566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 220c4762a1bSJed Brown return 0; 221c4762a1bSJed Brown } 222c4762a1bSJed Brown 223c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec x,Mat A,Mat BB,void *ctx) 224c4762a1bSJed Brown { 225c4762a1bSJed Brown PetscScalar v[3]; 226c4762a1bSJed Brown const PetscScalar *tmp; 227c4762a1bSJed Brown PetscInt idx[3],i; 228c4762a1bSJed Brown 229c4762a1bSJed Brown idx[0]=0; idx[1]=1; idx[2]=2; 2309566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(x,&tmp)); 231c4762a1bSJed Brown 232c4762a1bSJed Brown i = 0; 233c4762a1bSJed Brown v[0] = 2.0; v[1] = 1.0; v[2] = 0.0; 2349566063dSJacob Faibussowitsch PetscCall(MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES)); 235c4762a1bSJed Brown 236c4762a1bSJed Brown i = 1; 237c4762a1bSJed Brown v[0] = 1.0; v[1] = 2.0; v[2] = 1.0; 2389566063dSJacob Faibussowitsch PetscCall(MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES)); 239c4762a1bSJed Brown 240c4762a1bSJed Brown i = 2; 241c4762a1bSJed Brown v[0] = 0.0; v[1] = 1.0; v[2] = 2.0; 2429566063dSJacob Faibussowitsch PetscCall(MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES)); 243c4762a1bSJed Brown 2449566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(BB,MAT_FINAL_ASSEMBLY)); 2459566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(BB,MAT_FINAL_ASSEMBLY)); 246c4762a1bSJed Brown 247c4762a1bSJed Brown if (A != BB) { 2489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 2499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 250c4762a1bSJed Brown } 2519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(x,&tmp)); 252c4762a1bSJed Brown 253c4762a1bSJed Brown return 0; 254c4762a1bSJed Brown } 255c4762a1bSJed Brown 256c4762a1bSJed Brown /* 257c4762a1bSJed Brown The exact solutions 258c4762a1bSJed Brown */ 259c4762a1bSJed Brown PetscReal solx(PetscReal t) 260c4762a1bSJed Brown { 261c4762a1bSJed Brown return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)) + 262c4762a1bSJed Brown PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)); 263c4762a1bSJed Brown } 264c4762a1bSJed Brown 265c4762a1bSJed Brown PetscReal soly(PetscReal t) 266c4762a1bSJed Brown { 267c4762a1bSJed Brown return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/PetscSqrtReal(2.0) + 268c4762a1bSJed Brown PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/PetscSqrtReal(2.0); 269c4762a1bSJed Brown } 270c4762a1bSJed Brown 271c4762a1bSJed Brown PetscReal solz(PetscReal t) 272c4762a1bSJed Brown { 273c4762a1bSJed Brown return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)) + 274c4762a1bSJed Brown PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)); 275c4762a1bSJed Brown } 276c4762a1bSJed Brown 277c4762a1bSJed Brown /*TEST 278c4762a1bSJed Brown 279c4762a1bSJed Brown test: 280c4762a1bSJed Brown suffix: euler 281c4762a1bSJed Brown args: -ts_type euler 282c4762a1bSJed Brown requires: !single 283c4762a1bSJed Brown 284c4762a1bSJed Brown test: 285c4762a1bSJed Brown suffix: beuler 286c4762a1bSJed Brown args: -ts_type beuler 287c4762a1bSJed Brown requires: !single 288c4762a1bSJed Brown 289c4762a1bSJed Brown TEST*/ 290