1 static char help[] = "Basic problem for multi-rate method.\n"; 2 3 /*F 4 5 \begin{eqnarray} 6 ys' = -2.0*\frac{-1.0+ys^2.0-\cos(t)}{2.0*ys}+0.05*\frac{-2.0+yf^2-\cos(5.0*t)}{2.0*yf}-\frac{\sin(t)}{2.0*ys}\\ 7 yf' = 0.05*\frac{-1.0+ys^2-\cos(t)}{2.0*ys}-\frac{-2.0+yf^2-\cos(5.0*t)}{2.0*yf}-5.0*\frac{\sin(5.0*t)}{2.0*yf}\\ 8 \end{eqnarray} 9 10 F*/ 11 12 #include <petscts.h> 13 14 typedef struct { 15 PetscReal Tf,dt; 16 } AppCtx; 17 18 static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx) 19 { 20 const PetscScalar *u; 21 PetscScalar *f; 22 23 PetscFunctionBegin; 24 CHKERRQ(VecGetArrayRead(U,&u)); 25 CHKERRQ(VecGetArray(F,&f)); 26 f[0] = -2.0*(-1.0+u[0]*u[0]-PetscCosScalar(t))/(2.0*u[0])+0.05*(-2.0+u[1]*u[1]-PetscCosScalar(5.0*t))/(2.0*u[1])-PetscSinScalar(t)/(2.0*u[0]); 27 f[1] = 0.05*(-1.0+u[0]*u[0]-PetscCosScalar(t))/(2.0*u[0])-(-2.0+u[1]*u[1]-PetscCosScalar(5.0*t))/(2.0*u[1])-5.0*PetscSinScalar(5.0*t)/(2.0*u[1]); 28 CHKERRQ(VecRestoreArrayRead(U,&u)); 29 CHKERRQ(VecRestoreArray(F,&f)); 30 PetscFunctionReturn(0); 31 } 32 33 static PetscErrorCode RHSFunctionslow(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx) 34 { 35 const PetscScalar *u; 36 PetscScalar *f; 37 38 PetscFunctionBegin; 39 CHKERRQ(VecGetArrayRead(U,&u)); 40 CHKERRQ(VecGetArray(F,&f)); 41 f[0] = -2.0*(-1.0+u[0]*u[0]-PetscCosScalar(t))/(2.0*u[0])+0.05*(-2.0+u[1]*u[1]-PetscCosScalar(5.0*t))/(2.0*u[1])-PetscSinScalar(t)/(2.0*u[0]); 42 CHKERRQ(VecRestoreArrayRead(U,&u)); 43 CHKERRQ(VecRestoreArray(F,&f)); 44 PetscFunctionReturn(0); 45 } 46 47 static PetscErrorCode RHSFunctionfast(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx) 48 { 49 const PetscScalar *u; 50 PetscScalar *f; 51 52 PetscFunctionBegin; 53 CHKERRQ(VecGetArrayRead(U,&u)); 54 CHKERRQ(VecGetArray(F,&f)); 55 f[0] = 0.05*(-1.0+u[0]*u[0]-PetscCosScalar(t))/(2.0*u[0])-(-2.0+u[1]*u[1]-PetscCosScalar(5.0*t))/(2.0*u[1])-5.0*PetscSinScalar(5.0*t)/(2.0*u[1]); 56 CHKERRQ(VecRestoreArrayRead(U,&u)); 57 CHKERRQ(VecRestoreArray(F,&f)); 58 PetscFunctionReturn(0); 59 } 60 61 static PetscErrorCode sol_true(PetscReal t,Vec U) 62 { 63 PetscScalar *u; 64 65 PetscFunctionBegin; 66 CHKERRQ(VecGetArray(U,&u)); 67 u[0] = PetscSqrtScalar(1.0+PetscCosScalar(t)); 68 u[1] = PetscSqrtScalar(2.0+PetscCosScalar(5.0*t)); 69 CHKERRQ(VecRestoreArray(U,&u)); 70 PetscFunctionReturn(0); 71 } 72 73 int main(int argc,char **argv) 74 { 75 TS ts; /* ODE integrator */ 76 Vec U; /* solution will be stored here */ 77 Vec Utrue; 78 PetscErrorCode ierr; 79 PetscMPIInt size; 80 AppCtx ctx; 81 PetscScalar *u; 82 IS iss; 83 IS isf; 84 PetscInt *indicess; 85 PetscInt *indicesf; 86 PetscInt n=2; 87 PetscReal error,tt; 88 89 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 90 Initialize program 91 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 92 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 93 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 94 PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs"); 95 96 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 97 Create index for slow part and fast part 98 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 99 CHKERRQ(PetscMalloc1(1,&indicess)); 100 indicess[0]=0; 101 CHKERRQ(PetscMalloc1(1,&indicesf)); 102 indicesf[0]=1; 103 CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,1,indicess,PETSC_COPY_VALUES,&iss)); 104 CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,1,indicesf,PETSC_COPY_VALUES,&isf)); 105 106 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 107 Create necesary vector 108 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 109 CHKERRQ(VecCreate(PETSC_COMM_WORLD,&U)); 110 CHKERRQ(VecSetSizes(U,n,PETSC_DETERMINE)); 111 CHKERRQ(VecSetFromOptions(U)); 112 CHKERRQ(VecDuplicate(U,&Utrue)); 113 CHKERRQ(VecCopy(U,Utrue)); 114 115 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 116 Set initial condition 117 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 118 CHKERRQ(VecGetArray(U,&u)); 119 u[0] = PetscSqrtScalar(2.0); 120 u[1] = PetscSqrtScalar(3.0); 121 CHKERRQ(VecRestoreArray(U,&u)); 122 123 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 124 Create timestepping solver context 125 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 126 CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 127 CHKERRQ(TSSetType(ts,TSMPRK)); 128 129 CHKERRQ(TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx)); 130 CHKERRQ(TSRHSSplitSetIS(ts,"slow",iss)); 131 CHKERRQ(TSRHSSplitSetIS(ts,"fast",isf)); 132 CHKERRQ(TSRHSSplitSetRHSFunction(ts,"slow",NULL,(TSRHSFunction)RHSFunctionslow,&ctx)); 133 CHKERRQ(TSRHSSplitSetRHSFunction(ts,"fast",NULL,(TSRHSFunction)RHSFunctionfast,&ctx)); 134 135 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 136 Set initial conditions 137 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 138 CHKERRQ(TSSetSolution(ts,U)); 139 140 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 141 Set solver options 142 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 143 ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"ODE options","");CHKERRQ(ierr); 144 { 145 ctx.Tf = 0.3; 146 ctx.dt = 0.01; 147 CHKERRQ(PetscOptionsScalar("-Tf","","",ctx.Tf,&ctx.Tf,NULL)); 148 CHKERRQ(PetscOptionsScalar("-dt","","",ctx.dt,&ctx.dt,NULL)); 149 } 150 ierr = PetscOptionsEnd();CHKERRQ(ierr); 151 CHKERRQ(TSSetMaxTime(ts,ctx.Tf)); 152 CHKERRQ(TSSetTimeStep(ts,ctx.dt)); 153 CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 154 CHKERRQ(TSSetFromOptions(ts)); 155 156 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 157 Solve linear system 158 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 159 CHKERRQ(TSSolve(ts,U)); 160 CHKERRQ(VecView(U,PETSC_VIEWER_STDOUT_WORLD)); 161 162 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 163 Check the error of the Petsc solution 164 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 165 CHKERRQ(TSGetTime(ts,&tt)); 166 CHKERRQ(sol_true(tt,Utrue)); 167 CHKERRQ(VecAXPY(Utrue,-1.0,U)); 168 CHKERRQ(VecNorm(Utrue,NORM_2,&error)); 169 170 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 171 Print norm2 error 172 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 173 CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"l2 error norm: %g\n", error)); 174 175 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 176 Free work space. All PETSc objects should be destroyed when they are no longer needed. 177 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 178 CHKERRQ(VecDestroy(&U)); 179 CHKERRQ(TSDestroy(&ts)); 180 CHKERRQ(VecDestroy(&Utrue)); 181 CHKERRQ(ISDestroy(&iss)); 182 CHKERRQ(ISDestroy(&isf)); 183 CHKERRQ(PetscFree(indicess)); 184 CHKERRQ(PetscFree(indicesf)); 185 ierr = PetscFinalize(); 186 return ierr; 187 } 188 189 /*TEST 190 build: 191 requires: !complex 192 193 test: 194 195 TEST*/ 196