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 PetscErrorCode ierr; 33c4762a1bSJed Brown PetscInt time_steps = 100,steps; 34c4762a1bSJed Brown Vec global; 35c4762a1bSJed Brown PetscReal dt,ftime; 36c4762a1bSJed Brown TS ts; 37c4762a1bSJed Brown Mat A = 0,S; 38c4762a1bSJed Brown 39c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL)); 41c4762a1bSJed Brown 42c4762a1bSJed Brown /* set initial conditions */ 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreate(PETSC_COMM_WORLD,&global)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetSizes(global,PETSC_DECIDE,3)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetFromOptions(global)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(Initial(global,NULL)); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* make timestep context */ 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSMonitorSet(ts,Monitor,NULL,NULL)); 52c4762a1bSJed Brown dt = 0.001; 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* 55c4762a1bSJed Brown The user provides the RHS and Jacobian 56c4762a1bSJed Brown */ 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunction,NULL)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 59*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,3,3)); 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(RHSJacobian(ts,0.0,global,A,A,NULL)); 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,A,A,RHSJacobian,NULL)); 64c4762a1bSJed Brown 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateShell(PETSC_COMM_WORLD,3,3,3,3,NULL,&S)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(S,MATOP_MULT,(void (*)(void))MyMatMult)); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,S,A,RHSJacobian,NULL)); 68c4762a1bSJed Brown 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 71c4762a1bSJed Brown 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,dt)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxSteps(ts,time_steps)); 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,1)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,global)); 76c4762a1bSJed Brown 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,global)); 78*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetSolveTime(ts,&ftime)); 79*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetStepNumber(ts,&steps)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* free the memories */ 82c4762a1bSJed Brown 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&global)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&S)); 87c4762a1bSJed Brown 88c4762a1bSJed Brown ierr = PetscFinalize(); 89c4762a1bSJed Brown return ierr; 90c4762a1bSJed Brown } 91c4762a1bSJed Brown 92c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat S,Vec x,Vec y) 93c4762a1bSJed Brown { 94c4762a1bSJed Brown const PetscScalar *inptr; 95c4762a1bSJed Brown PetscScalar *outptr; 96c4762a1bSJed Brown 97c4762a1bSJed Brown PetscFunctionBegin; 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&inptr)); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayWrite(y,&outptr)); 100c4762a1bSJed Brown 101c4762a1bSJed Brown outptr[0] = 2.0*inptr[0]+inptr[1]; 102c4762a1bSJed Brown outptr[1] = inptr[0]+2.0*inptr[1]+inptr[2]; 103c4762a1bSJed Brown outptr[2] = inptr[1]+2.0*inptr[2]; 104c4762a1bSJed Brown 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&inptr)); 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayWrite(y,&outptr)); 107c4762a1bSJed Brown PetscFunctionReturn(0); 108c4762a1bSJed Brown } 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* this test problem has initial values (1,1,1). */ 111c4762a1bSJed Brown PetscErrorCode Initial(Vec global,void *ctx) 112c4762a1bSJed Brown { 113c4762a1bSJed Brown PetscScalar *localptr; 114c4762a1bSJed Brown PetscInt i,mybase,myend,locsize; 115c4762a1bSJed Brown 116c4762a1bSJed Brown /* determine starting point of each processor */ 117*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(global,&mybase,&myend)); 118*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(global,&locsize)); 119c4762a1bSJed Brown 120c4762a1bSJed Brown /* Initialize the array */ 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayWrite(global,&localptr)); 122c4762a1bSJed Brown for (i=0; i<locsize; i++) localptr[i] = 1.0; 123c4762a1bSJed Brown 124c4762a1bSJed Brown if (mybase == 0) localptr[0]=1.0; 125c4762a1bSJed Brown 126*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayWrite(global,&localptr)); 127c4762a1bSJed Brown return 0; 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130c4762a1bSJed Brown PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec global,void *ctx) 131c4762a1bSJed Brown { 132c4762a1bSJed Brown VecScatter scatter; 133c4762a1bSJed Brown IS from,to; 134c4762a1bSJed Brown PetscInt i,n,*idx; 135c4762a1bSJed Brown Vec tmp_vec; 136c4762a1bSJed Brown PetscErrorCode ierr; 137c4762a1bSJed Brown const PetscScalar *tmp; 138c4762a1bSJed Brown 139c4762a1bSJed Brown /* Get the size of the vector */ 140*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(global,&n)); 141c4762a1bSJed Brown 142c4762a1bSJed Brown /* Set the index sets */ 143*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n,&idx)); 144c4762a1bSJed Brown for (i=0; i<n; i++) idx[i]=i; 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* Create local sequential vectors */ 147*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&tmp_vec)); 148c4762a1bSJed Brown 149c4762a1bSJed Brown /* Create scatter context */ 150*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from)); 151*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to)); 152*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(global,from,tmp_vec,to,&scatter)); 153*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD)); 154*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD)); 155c4762a1bSJed Brown 156*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(tmp_vec,&tmp)); 157c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e u = %14.6e %14.6e %14.6e \n", 158c4762a1bSJed Brown (double)time,(double)PetscRealPart(tmp[0]),(double)PetscRealPart(tmp[1]),(double)PetscRealPart(tmp[2]));CHKERRQ(ierr); 159c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e errors = %14.6e %14.6e %14.6e \n", 160c4762a1bSJed Brown (double)time,(double)PetscRealPart(tmp[0]-solx(time)),(double)PetscRealPart(tmp[1]-soly(time)),(double)PetscRealPart(tmp[2]-solz(time)));CHKERRQ(ierr); 161*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(tmp_vec,&tmp)); 162*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&scatter)); 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&from)); 164*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&to)); 165*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(idx)); 166*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&tmp_vec)); 167c4762a1bSJed Brown return 0; 168c4762a1bSJed Brown } 169c4762a1bSJed Brown 170c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx) 171c4762a1bSJed Brown { 172c4762a1bSJed Brown PetscScalar *outptr; 173c4762a1bSJed Brown const PetscScalar *inptr; 174c4762a1bSJed Brown PetscInt i,n,*idx; 175c4762a1bSJed Brown IS from,to; 176c4762a1bSJed Brown VecScatter scatter; 177c4762a1bSJed Brown Vec tmp_in,tmp_out; 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* Get the length of parallel vector */ 180*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(globalin,&n)); 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* Set the index sets */ 183*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n,&idx)); 184c4762a1bSJed Brown for (i=0; i<n; i++) idx[i]=i; 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* Create local sequential vectors */ 187*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,n,&tmp_in)); 188*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(tmp_in,&tmp_out)); 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* Create scatter context */ 191*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from)); 192*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to)); 193*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(globalin,from,tmp_in,to,&scatter)); 194*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD)); 195*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD)); 196*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&scatter)); 197c4762a1bSJed Brown 198c4762a1bSJed Brown /*Extract income array */ 199*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(tmp_in,&inptr)); 200c4762a1bSJed Brown 201c4762a1bSJed Brown /* Extract outcome array*/ 202*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayWrite(tmp_out,&outptr)); 203c4762a1bSJed Brown 204c4762a1bSJed Brown outptr[0] = 2.0*inptr[0]+inptr[1]; 205c4762a1bSJed Brown outptr[1] = inptr[0]+2.0*inptr[1]+inptr[2]; 206c4762a1bSJed Brown outptr[2] = inptr[1]+2.0*inptr[2]; 207c4762a1bSJed Brown 208*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(tmp_in,&inptr)); 209*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayWrite(tmp_out,&outptr)); 210c4762a1bSJed Brown 211*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(tmp_out,from,globalout,to,&scatter)); 212*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD)); 213*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD)); 214c4762a1bSJed Brown 215c4762a1bSJed Brown /* Destroy idx aand scatter */ 216*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&from)); 217*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&to)); 218*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&scatter)); 219*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&tmp_in)); 220*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&tmp_out)); 221*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(idx)); 222c4762a1bSJed Brown return 0; 223c4762a1bSJed Brown } 224c4762a1bSJed Brown 225c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec x,Mat A,Mat BB,void *ctx) 226c4762a1bSJed Brown { 227c4762a1bSJed Brown PetscScalar v[3]; 228c4762a1bSJed Brown const PetscScalar *tmp; 229c4762a1bSJed Brown PetscInt idx[3],i; 230c4762a1bSJed Brown 231c4762a1bSJed Brown idx[0]=0; idx[1]=1; idx[2]=2; 232*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(x,&tmp)); 233c4762a1bSJed Brown 234c4762a1bSJed Brown i = 0; 235c4762a1bSJed Brown v[0] = 2.0; v[1] = 1.0; v[2] = 0.0; 236*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES)); 237c4762a1bSJed Brown 238c4762a1bSJed Brown i = 1; 239c4762a1bSJed Brown v[0] = 1.0; v[1] = 2.0; v[2] = 1.0; 240*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES)); 241c4762a1bSJed Brown 242c4762a1bSJed Brown i = 2; 243c4762a1bSJed Brown v[0] = 0.0; v[1] = 1.0; v[2] = 2.0; 244*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(BB,1,&i,3,idx,v,INSERT_VALUES)); 245c4762a1bSJed Brown 246*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(BB,MAT_FINAL_ASSEMBLY)); 247*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(BB,MAT_FINAL_ASSEMBLY)); 248c4762a1bSJed Brown 249c4762a1bSJed Brown if (A != BB) { 250*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 251*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 252c4762a1bSJed Brown } 253*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(x,&tmp)); 254c4762a1bSJed Brown 255c4762a1bSJed Brown return 0; 256c4762a1bSJed Brown } 257c4762a1bSJed Brown 258c4762a1bSJed Brown /* 259c4762a1bSJed Brown The exact solutions 260c4762a1bSJed Brown */ 261c4762a1bSJed Brown PetscReal solx(PetscReal t) 262c4762a1bSJed Brown { 263c4762a1bSJed Brown return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)) + 264c4762a1bSJed Brown PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)); 265c4762a1bSJed Brown } 266c4762a1bSJed Brown 267c4762a1bSJed Brown PetscReal soly(PetscReal t) 268c4762a1bSJed Brown { 269c4762a1bSJed Brown return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/PetscSqrtReal(2.0) + 270c4762a1bSJed Brown PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/PetscSqrtReal(2.0); 271c4762a1bSJed Brown } 272c4762a1bSJed Brown 273c4762a1bSJed Brown PetscReal solz(PetscReal t) 274c4762a1bSJed Brown { 275c4762a1bSJed Brown return PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/2.0 - PetscExpReal((2.0 - PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)) + 276c4762a1bSJed Brown PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/2.0 + PetscExpReal((2.0 + PetscSqrtReal(2.0))*t)/(2.0*PetscSqrtReal(2.0)); 277c4762a1bSJed Brown } 278c4762a1bSJed Brown 279c4762a1bSJed Brown /*TEST 280c4762a1bSJed Brown 281c4762a1bSJed Brown test: 282c4762a1bSJed Brown suffix: euler 283c4762a1bSJed Brown args: -ts_type euler 284c4762a1bSJed Brown requires: !single 285c4762a1bSJed Brown 286c4762a1bSJed Brown test: 287c4762a1bSJed Brown suffix: beuler 288c4762a1bSJed Brown args: -ts_type beuler 289c4762a1bSJed Brown requires: !single 290c4762a1bSJed Brown 291c4762a1bSJed Brown TEST*/ 292