155a74a43SLisandro Dalcin #include <petsc.h>
255a74a43SLisandro Dalcin
355a74a43SLisandro Dalcin EXTERN_C_BEGIN
455a74a43SLisandro Dalcin extern void formInitial(int*,int*,int*,double*,
555a74a43SLisandro Dalcin double*,double*);
655a74a43SLisandro Dalcin extern void formFunction(const int*,const int*,const int*,const double*,
755a74a43SLisandro Dalcin const double*,const double[],const double[],double[]);
855a74a43SLisandro Dalcin EXTERN_C_END
955a74a43SLisandro Dalcin
1055a74a43SLisandro Dalcin typedef struct AppCtx {
1155a74a43SLisandro Dalcin PetscInt nx,ny,nz;
1255a74a43SLisandro Dalcin PetscScalar h[3];
1355a74a43SLisandro Dalcin } AppCtx;
1455a74a43SLisandro Dalcin
FormInitial(PetscReal t,Vec X,PetscCtx ctx)15*2a8381b2SBarry Smith PetscErrorCode FormInitial(PetscReal t, Vec X, PetscCtx ctx)
1655a74a43SLisandro Dalcin {
1755a74a43SLisandro Dalcin PetscScalar *x;
1855a74a43SLisandro Dalcin AppCtx *app = (AppCtx*) ctx;
194d86920dSPierre Jolivet
2055a74a43SLisandro Dalcin PetscFunctionBegin;
2155a74a43SLisandro Dalcin PetscCall(VecGetArray(X,&x));
2255a74a43SLisandro Dalcin /**/
2355a74a43SLisandro Dalcin formInitial(&app->nx,&app->ny,&app->nz,app->h,&t,x);
2455a74a43SLisandro Dalcin /**/
2555a74a43SLisandro Dalcin PetscCall(VecRestoreArray(X,&x));
2655a74a43SLisandro Dalcin PetscFunctionReturn(PETSC_SUCCESS);
2755a74a43SLisandro Dalcin }
2855a74a43SLisandro Dalcin
FormFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,PetscCtx ctx)29*2a8381b2SBarry Smith PetscErrorCode FormFunction(TS ts, PetscReal t, Vec X, Vec Xdot,Vec F, PetscCtx ctx)
3055a74a43SLisandro Dalcin {
3155a74a43SLisandro Dalcin const PetscScalar *x;
3255a74a43SLisandro Dalcin const PetscScalar *xdot;
3355a74a43SLisandro Dalcin PetscScalar *f;
3455a74a43SLisandro Dalcin AppCtx *app = (AppCtx*) ctx;
354d86920dSPierre Jolivet
3655a74a43SLisandro Dalcin PetscFunctionBegin;
3755a74a43SLisandro Dalcin PetscCall(VecGetArrayRead(X,&x));
3855a74a43SLisandro Dalcin PetscCall(VecGetArrayRead(Xdot,&xdot));
3955a74a43SLisandro Dalcin PetscCall(VecGetArray(F,&f));
4055a74a43SLisandro Dalcin /**/
4155a74a43SLisandro Dalcin formFunction(&app->nx,&app->ny,&app->nz,app->h,&t,x,xdot,f);
4255a74a43SLisandro Dalcin /**/
4355a74a43SLisandro Dalcin PetscCall(VecRestoreArrayRead(X,&x));
4455a74a43SLisandro Dalcin PetscCall(VecRestoreArrayRead(Xdot,&xdot));
4555a74a43SLisandro Dalcin PetscCall(VecRestoreArray(F,&f));
4655a74a43SLisandro Dalcin PetscFunctionReturn(PETSC_SUCCESS);
4755a74a43SLisandro Dalcin }
4855a74a43SLisandro Dalcin
RunTest(int nx,int ny,int nz,int loops,double * wt)4955a74a43SLisandro Dalcin PetscErrorCode RunTest(int nx, int ny, int nz, int loops, double *wt)
5055a74a43SLisandro Dalcin {
5155a74a43SLisandro Dalcin Vec x,f;
5255a74a43SLisandro Dalcin TS ts;
5355a74a43SLisandro Dalcin AppCtx _app,*app=&_app;
5455a74a43SLisandro Dalcin double t1,t2;
5555a74a43SLisandro Dalcin
564d86920dSPierre Jolivet PetscFunctionBegin;
5755a74a43SLisandro Dalcin app->nx = nx; app->h[0] = 1./(nx-1);
5855a74a43SLisandro Dalcin app->ny = ny; app->h[1] = 1./(ny-1);
5955a74a43SLisandro Dalcin app->nz = nz; app->h[2] = 1./(nz-1);
6055a74a43SLisandro Dalcin
6155a74a43SLisandro Dalcin PetscCall(VecCreate(PETSC_COMM_SELF,&x));
6255a74a43SLisandro Dalcin PetscCall(VecSetSizes(x,nx*ny*nz,nx*ny*nz));
6355a74a43SLisandro Dalcin PetscCall(VecSetUp(x));
6455a74a43SLisandro Dalcin PetscCall(VecDuplicate(x,&f));
6555a74a43SLisandro Dalcin
6655a74a43SLisandro Dalcin PetscCall(TSCreate(PETSC_COMM_SELF,&ts));
6755a74a43SLisandro Dalcin PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
6855a74a43SLisandro Dalcin PetscCall(TSSetType(ts,TSTHETA));
6955a74a43SLisandro Dalcin PetscCall(TSThetaSetTheta(ts,1.0));
7055a74a43SLisandro Dalcin PetscCall(TSSetTimeStep(ts,0.01));
7155a74a43SLisandro Dalcin PetscCall(TSSetTime(ts,0.0));
7255a74a43SLisandro Dalcin PetscCall(TSSetMaxTime(ts,1.0));
7355a74a43SLisandro Dalcin PetscCall(TSSetMaxSteps(ts,10));
7455a74a43SLisandro Dalcin PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
7555a74a43SLisandro Dalcin
7655a74a43SLisandro Dalcin PetscCall(TSSetSolution(ts,x));
7755a74a43SLisandro Dalcin PetscCall(TSSetIFunction(ts,f,FormFunction,app));
7855a74a43SLisandro Dalcin PetscCall(PetscOptionsSetValue(NULL,"-snes_mf","1"));
7955a74a43SLisandro Dalcin {
8055a74a43SLisandro Dalcin SNES snes;
8155a74a43SLisandro Dalcin KSP ksp;
8255a74a43SLisandro Dalcin PetscCall(TSGetSNES(ts,&snes));
8355a74a43SLisandro Dalcin PetscCall(SNESGetKSP(snes,&ksp));
8455a74a43SLisandro Dalcin PetscCall(KSPSetType(ksp,KSPCG));
8555a74a43SLisandro Dalcin }
8655a74a43SLisandro Dalcin PetscCall(TSSetFromOptions(ts));
8755a74a43SLisandro Dalcin PetscCall(TSSetUp(ts));
8855a74a43SLisandro Dalcin
8955a74a43SLisandro Dalcin *wt = 1e300;
9055a74a43SLisandro Dalcin while (loops-- > 0) {
9155a74a43SLisandro Dalcin PetscCall(FormInitial(0.0,x,app));
9255a74a43SLisandro Dalcin PetscCall(PetscTime(&t1));
9355a74a43SLisandro Dalcin PetscCall(TSSolve(ts,x));
9455a74a43SLisandro Dalcin PetscCall(PetscTime(&t2));
9555a74a43SLisandro Dalcin *wt = PetscMin(*wt,t2-t1);
9655a74a43SLisandro Dalcin }
9755a74a43SLisandro Dalcin
9855a74a43SLisandro Dalcin PetscCall(VecDestroy(&x));
9955a74a43SLisandro Dalcin PetscCall(VecDestroy(&f));
10055a74a43SLisandro Dalcin PetscCall(TSDestroy(&ts));
10155a74a43SLisandro Dalcin PetscFunctionReturn(PETSC_SUCCESS);
10255a74a43SLisandro Dalcin }
10355a74a43SLisandro Dalcin
GetInt(const char * name,PetscInt * v,PetscInt defv)10455a74a43SLisandro Dalcin PetscErrorCode GetInt(const char* name, PetscInt *v, PetscInt defv)
10555a74a43SLisandro Dalcin {
10655a74a43SLisandro Dalcin PetscFunctionBegin;
10755a74a43SLisandro Dalcin *v = defv;
10855a74a43SLisandro Dalcin PetscCall(PetscOptionsGetInt(NULL,NULL,name,v,NULL));
10955a74a43SLisandro Dalcin PetscFunctionReturn(PETSC_SUCCESS);
11055a74a43SLisandro Dalcin }
11155a74a43SLisandro Dalcin
main(int argc,char * argv[])11255a74a43SLisandro Dalcin int main(int argc, char *argv[])
11355a74a43SLisandro Dalcin {
11455a74a43SLisandro Dalcin double wt;
11555a74a43SLisandro Dalcin PetscInt n,start,step,stop,samples;
11655a74a43SLisandro Dalcin
11755a74a43SLisandro Dalcin PetscCall(PetscInitialize(&argc,&argv,NULL,NULL));
11855a74a43SLisandro Dalcin
11955a74a43SLisandro Dalcin PetscCall(GetInt("-start", &start, 12));
12055a74a43SLisandro Dalcin PetscCall(GetInt("-step", &step, 4));
12155a74a43SLisandro Dalcin PetscCall(GetInt("-stop", &stop, start));
12255a74a43SLisandro Dalcin PetscCall(GetInt("-samples", &samples, 1));
12355a74a43SLisandro Dalcin
12455a74a43SLisandro Dalcin for (n=start; n<=stop; n+=step) {
12555a74a43SLisandro Dalcin int nx=n+1, ny=n+1, nz=n+1;
12655a74a43SLisandro Dalcin PetscCall(RunTest(nx,ny,nz,samples,&wt));
12755a74a43SLisandro Dalcin PetscCall(PetscPrintf(PETSC_COMM_SELF,"Grid %3d x %3d x %3d -> %f seconds (%2d samples)\n",nx,ny,nz,wt,samples));
12855a74a43SLisandro Dalcin }
12955a74a43SLisandro Dalcin PetscCall(PetscFinalize());
13055a74a43SLisandro Dalcin return 0;
13155a74a43SLisandro Dalcin }
132