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