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 38*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL)); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* set initial conditions */ 425f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&global)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(global,PETSC_DECIDE,3)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(global)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(Initial(global,NULL)); 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* make timestep context */ 485f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 565f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunction,NULL)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,3,3)); 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(RHSJacobian(ts,0.0,global,A,A,NULL)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,A,A,RHSJacobian,NULL)); 63c4762a1bSJed Brown 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,3,3,3,3,NULL,&S)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S,MATOP_MULT,(void (*)(void))MyMatMult)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,S,A,RHSJacobian,NULL)); 67c4762a1bSJed Brown 685f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 70c4762a1bSJed Brown 715f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,dt)); 725f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts,time_steps)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,1)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,global)); 75c4762a1bSJed Brown 765f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,global)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolveTime(ts,&ftime)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetStepNumber(ts,&steps)); 79c4762a1bSJed Brown 80c4762a1bSJed Brown /* free the memories */ 81c4762a1bSJed Brown 825f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&global)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&S)); 86c4762a1bSJed Brown 87*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 88*b122ec5aSJacob 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; 975f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&inptr)); 985f80ce2aSJacob Faibussowitsch CHKERRQ(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 1045f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&inptr)); 1055f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 1165f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(global,&mybase,&myend)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(global,&locsize)); 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* Initialize the array */ 1205f80ce2aSJacob Faibussowitsch CHKERRQ(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 1255f80ce2aSJacob Faibussowitsch CHKERRQ(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 PetscErrorCode ierr; 136c4762a1bSJed Brown const PetscScalar *tmp; 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* Get the size of the vector */ 1395f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(global,&n)); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* Set the index sets */ 1425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n,&idx)); 143c4762a1bSJed Brown for (i=0; i<n; i++) idx[i]=i; 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* Create local sequential vectors */ 1465f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&tmp_vec)); 147c4762a1bSJed Brown 148c4762a1bSJed Brown /* Create scatter context */ 1495f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from)); 1505f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(global,from,tmp_vec,to,&scatter)); 1525f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD)); 1535f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD)); 154c4762a1bSJed Brown 1555f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(tmp_vec,&tmp)); 156c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e u = %14.6e %14.6e %14.6e \n", 157c4762a1bSJed Brown (double)time,(double)PetscRealPart(tmp[0]),(double)PetscRealPart(tmp[1]),(double)PetscRealPart(tmp[2]));CHKERRQ(ierr); 158c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e errors = %14.6e %14.6e %14.6e \n", 159c4762a1bSJed Brown (double)time,(double)PetscRealPart(tmp[0]-solx(time)),(double)PetscRealPart(tmp[1]-soly(time)),(double)PetscRealPart(tmp[2]-solz(time)));CHKERRQ(ierr); 1605f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(tmp_vec,&tmp)); 1615f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&scatter)); 1625f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&from)); 1635f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&to)); 1645f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(idx)); 1655f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&tmp_vec)); 166c4762a1bSJed Brown return 0; 167c4762a1bSJed Brown } 168c4762a1bSJed Brown 169c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx) 170c4762a1bSJed Brown { 171c4762a1bSJed Brown PetscScalar *outptr; 172c4762a1bSJed Brown const PetscScalar *inptr; 173c4762a1bSJed Brown PetscInt i,n,*idx; 174c4762a1bSJed Brown IS from,to; 175c4762a1bSJed Brown VecScatter scatter; 176c4762a1bSJed Brown Vec tmp_in,tmp_out; 177c4762a1bSJed Brown 178c4762a1bSJed Brown /* Get the length of parallel vector */ 1795f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(globalin,&n)); 180c4762a1bSJed Brown 181c4762a1bSJed Brown /* Set the index sets */ 1825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n,&idx)); 183c4762a1bSJed Brown for (i=0; i<n; i++) idx[i]=i; 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* Create local sequential vectors */ 1865f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&tmp_in)); 1875f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(tmp_in,&tmp_out)); 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* Create scatter context */ 1905f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from)); 1915f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to)); 1925f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(globalin,from,tmp_in,to,&scatter)); 1935f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD)); 1945f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD)); 1955f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&scatter)); 196c4762a1bSJed Brown 197c4762a1bSJed Brown /*Extract income array */ 1985f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(tmp_in,&inptr)); 199c4762a1bSJed Brown 200c4762a1bSJed Brown /* Extract outcome array*/ 2015f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayWrite(tmp_out,&outptr)); 202c4762a1bSJed Brown 203c4762a1bSJed Brown outptr[0] = 2.0*inptr[0]+inptr[1]; 204c4762a1bSJed Brown outptr[1] = inptr[0]+2.0*inptr[1]+inptr[2]; 205c4762a1bSJed Brown outptr[2] = inptr[1]+2.0*inptr[2]; 206c4762a1bSJed Brown 2075f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(tmp_in,&inptr)); 2085f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayWrite(tmp_out,&outptr)); 209c4762a1bSJed Brown 2105f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(tmp_out,from,globalout,to,&scatter)); 2115f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD)); 213c4762a1bSJed Brown 214c4762a1bSJed Brown /* Destroy idx aand scatter */ 2155f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&from)); 2165f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&to)); 2175f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&scatter)); 2185f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&tmp_in)); 2195f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&tmp_out)); 2205f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(idx)); 221c4762a1bSJed Brown return 0; 222c4762a1bSJed Brown } 223c4762a1bSJed Brown 224c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec x,Mat A,Mat BB,void *ctx) 225c4762a1bSJed Brown { 226c4762a1bSJed Brown PetscScalar v[3]; 227c4762a1bSJed Brown const PetscScalar *tmp; 228c4762a1bSJed Brown PetscInt idx[3],i; 229c4762a1bSJed Brown 230c4762a1bSJed Brown idx[0]=0; idx[1]=1; idx[2]=2; 2315f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&tmp)); 232c4762a1bSJed Brown 233c4762a1bSJed Brown i = 0; 234c4762a1bSJed Brown v[0] = 2.0; v[1] = 1.0; v[2] = 0.0; 2355f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES)); 236c4762a1bSJed Brown 237c4762a1bSJed Brown i = 1; 238c4762a1bSJed Brown v[0] = 1.0; v[1] = 2.0; v[2] = 1.0; 2395f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES)); 240c4762a1bSJed Brown 241c4762a1bSJed Brown i = 2; 242c4762a1bSJed Brown v[0] = 0.0; v[1] = 1.0; v[2] = 2.0; 2435f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES)); 244c4762a1bSJed Brown 2455f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(BB,MAT_FINAL_ASSEMBLY)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(BB,MAT_FINAL_ASSEMBLY)); 247c4762a1bSJed Brown 248c4762a1bSJed Brown if (A != BB) { 2495f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 2505f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 251c4762a1bSJed Brown } 2525f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&tmp)); 253c4762a1bSJed Brown 254c4762a1bSJed Brown return 0; 255c4762a1bSJed Brown } 256c4762a1bSJed Brown 257c4762a1bSJed Brown /* 258c4762a1bSJed Brown The exact solutions 259c4762a1bSJed Brown */ 260c4762a1bSJed Brown PetscReal solx(PetscReal t) 261c4762a1bSJed Brown { 262c4762a1bSJed Brown return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)) + 263c4762a1bSJed Brown PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)); 264c4762a1bSJed Brown } 265c4762a1bSJed Brown 266c4762a1bSJed Brown PetscReal soly(PetscReal t) 267c4762a1bSJed Brown { 268c4762a1bSJed Brown return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/PetscSqrtReal(2.0) + 269c4762a1bSJed Brown PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/PetscSqrtReal(2.0); 270c4762a1bSJed Brown } 271c4762a1bSJed Brown 272c4762a1bSJed Brown PetscReal solz(PetscReal t) 273c4762a1bSJed Brown { 274c4762a1bSJed Brown return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)) + 275c4762a1bSJed Brown PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)); 276c4762a1bSJed Brown } 277c4762a1bSJed Brown 278c4762a1bSJed Brown /*TEST 279c4762a1bSJed Brown 280c4762a1bSJed Brown test: 281c4762a1bSJed Brown suffix: euler 282c4762a1bSJed Brown args: -ts_type euler 283c4762a1bSJed Brown requires: !single 284c4762a1bSJed Brown 285c4762a1bSJed Brown test: 286c4762a1bSJed Brown suffix: beuler 287c4762a1bSJed Brown args: -ts_type beuler 288c4762a1bSJed Brown requires: !single 289c4762a1bSJed Brown 290c4762a1bSJed Brown TEST*/ 291