xref: /petsc/src/snes/tutorials/ex47cu.cu (revision e09e07f2dd1b33d9446359ab9cf638f5f341c573)
1c4762a1bSJed Brown static char help[] = "Solves -Laplacian u - exp(u) = 0,  0 < x < 1 using GPU\n\n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown    Same as ex47.c except it also uses the GPU to evaluate the function
4c4762a1bSJed Brown */
5c4762a1bSJed Brown 
6c4762a1bSJed Brown #include <petscdm.h>
7c4762a1bSJed Brown #include <petscdmda.h>
8c4762a1bSJed Brown #include <petscsnes.h>
9c4762a1bSJed Brown 
10c4762a1bSJed Brown #include <thrust/device_ptr.h>
11c4762a1bSJed Brown #include <thrust/for_each.h>
12c4762a1bSJed Brown #include <thrust/tuple.h>
13c4762a1bSJed Brown #include <thrust/iterator/constant_iterator.h>
14c4762a1bSJed Brown #include <thrust/iterator/counting_iterator.h>
15c4762a1bSJed Brown #include <thrust/iterator/zip_iterator.h>
16c4762a1bSJed Brown 
17c4762a1bSJed Brown extern PetscErrorCode ComputeFunction(SNES,Vec,Vec,void*), ComputeJacobian(SNES,Vec,Mat,Mat,void*);
18c4762a1bSJed Brown PetscBool useCUDA = PETSC_FALSE;
19c4762a1bSJed Brown 
20c4762a1bSJed Brown int main(int argc,char **argv)
21c4762a1bSJed Brown {
22c4762a1bSJed Brown   SNES           snes;
23c4762a1bSJed Brown   Vec            x,f;
24c4762a1bSJed Brown   Mat            J;
25c4762a1bSJed Brown   DM             da;
26c4762a1bSJed Brown   char           *tmp,typeName[256];
27c4762a1bSJed Brown   PetscBool      flg;
28c4762a1bSJed Brown 
29327415f7SBarry Smith   PetscFunctionBeginUser;
309566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-dm_vec_type",typeName,sizeof(typeName),&flg));
32c4762a1bSJed Brown   if (flg) {
339566063dSJacob Faibussowitsch     PetscCall(PetscStrstr(typeName,"cuda",&tmp));
34c4762a1bSJed Brown     if (tmp) useCUDA = PETSC_TRUE;
35c4762a1bSJed Brown   }
36c4762a1bSJed Brown 
379566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,8,1,1,NULL,&da));
389566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
399566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
40*e09e07f2SBarry Smith   PetscCall(DMCreateGlobalVector(da,&x));
41*e09e07f2SBarry Smith   PetscCall(VecDuplicate(x,&f));
429566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(da,&J));
43c4762a1bSJed Brown 
449566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
459566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes,f,ComputeFunction,da));
469566063dSJacob Faibussowitsch   PetscCall(SNESSetJacobian(snes,J,J,ComputeJacobian,da));
479566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes));
489566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes,NULL,x));
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
519566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&f));
539566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes));
549566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
55c4762a1bSJed Brown 
569566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
57b122ec5aSJacob Faibussowitsch   return 0;
58c4762a1bSJed Brown }
59c4762a1bSJed Brown 
60c4762a1bSJed Brown struct ApplyStencil
61c4762a1bSJed Brown {
62c4762a1bSJed Brown   template <typename Tuple>
63c4762a1bSJed Brown   __host__ __device__
64c4762a1bSJed Brown   void operator()(Tuple t)
65c4762a1bSJed Brown   {
66c4762a1bSJed Brown     /* f = (2*x_i - x_(i+1) - x_(i-1))/h - h*exp(x_i) */
67c4762a1bSJed Brown     thrust::get<0>(t) = 1;
68c4762a1bSJed Brown     if ((thrust::get<4>(t) > 0) && (thrust::get<4>(t) < thrust::get<5>(t)-1)) {
6988e3df94SJunchao Zhang       thrust::get<0>(t) = (((PetscScalar)2.0)*thrust::get<1>(t) - thrust::get<2>(t) - thrust::get<3>(t)) / (thrust::get<6>(t)) - (thrust::get<6>(t))*exp(thrust::get<1>(t));
70c4762a1bSJed Brown     } else if (thrust::get<4>(t) == 0) {
71c4762a1bSJed Brown       thrust::get<0>(t) = thrust::get<1>(t) / (thrust::get<6>(t));
72c4762a1bSJed Brown     } else if (thrust::get<4>(t) == thrust::get<5>(t)-1) {
73c4762a1bSJed Brown       thrust::get<0>(t) = thrust::get<1>(t) / (thrust::get<6>(t));
74c4762a1bSJed Brown     }
75c4762a1bSJed Brown   }
76c4762a1bSJed Brown };
77c4762a1bSJed Brown 
78c4762a1bSJed Brown PetscErrorCode ComputeFunction(SNES snes,Vec x,Vec f,void *ctx)
79c4762a1bSJed Brown {
80c4762a1bSJed Brown   PetscInt          i,Mx,xs,xm,xstartshift,xendshift,fstart,lsize;
81c4762a1bSJed Brown   PetscScalar       *xx,*ff,hx;
82c4762a1bSJed Brown   DM                da = (DM) ctx;
83c4762a1bSJed Brown   Vec               xlocal;
84c4762a1bSJed Brown   PetscMPIInt       rank,size;
85c4762a1bSJed Brown   MPI_Comm          comm;
86c4762a1bSJed Brown   PetscScalar const *xarray;
87c4762a1bSJed Brown   PetscScalar       *farray;
88c4762a1bSJed Brown 
899566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
90c4762a1bSJed Brown   hx   = 1.0/(PetscReal)(Mx-1);
919566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&xlocal));
929566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal));
939566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal));
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   if (useCUDA) {
969566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetComm((PetscObject)da,&comm));
979566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_size(comm,&size));
989566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(comm,&rank));
9906cf5b03SBarry Smith     PetscCall(VecCUDAGetArrayRead(xlocal,&xarray));
10006cf5b03SBarry Smith     PetscCall(VecCUDAGetArrayWrite(f,&farray));
101c4762a1bSJed Brown     if (rank) xstartshift = 1;
102c4762a1bSJed Brown     else xstartshift = 0;
103c4762a1bSJed Brown     if (rank != size-1) xendshift = 1;
104c4762a1bSJed Brown     else xendshift = 0;
1059566063dSJacob Faibussowitsch     PetscCall(VecGetOwnershipRange(f,&fstart,NULL));
1069566063dSJacob Faibussowitsch     PetscCall(VecGetLocalSize(x,&lsize));
107c4762a1bSJed Brown     try {
108c4762a1bSJed Brown       thrust::for_each(
109c4762a1bSJed Brown         thrust::make_zip_iterator(
110c4762a1bSJed Brown           thrust::make_tuple(
111c4762a1bSJed Brown             thrust::device_ptr<PetscScalar>(farray),
112c4762a1bSJed Brown             thrust::device_ptr<const PetscScalar>(xarray + xstartshift),
113c4762a1bSJed Brown             thrust::device_ptr<const PetscScalar>(xarray + xstartshift + 1),
114c4762a1bSJed Brown             thrust::device_ptr<const PetscScalar>(xarray + xstartshift - 1),
115c4762a1bSJed Brown             thrust::counting_iterator<int>(fstart),
116c4762a1bSJed Brown             thrust::constant_iterator<int>(Mx),
117c4762a1bSJed Brown             thrust::constant_iterator<PetscScalar>(hx))),
118c4762a1bSJed Brown         thrust::make_zip_iterator(
119c4762a1bSJed Brown           thrust::make_tuple(
120c4762a1bSJed Brown             thrust::device_ptr<PetscScalar>(farray + lsize),
121c4762a1bSJed Brown             thrust::device_ptr<const PetscScalar>(xarray + lsize - xendshift),
122c4762a1bSJed Brown             thrust::device_ptr<const PetscScalar>(xarray + lsize - xendshift + 1),
123c4762a1bSJed Brown             thrust::device_ptr<const PetscScalar>(xarray + lsize - xendshift - 1),
124c4762a1bSJed Brown             thrust::counting_iterator<int>(fstart) + lsize,
125c4762a1bSJed Brown             thrust::constant_iterator<int>(Mx),
126c4762a1bSJed Brown             thrust::constant_iterator<PetscScalar>(hx))),
127c4762a1bSJed Brown         ApplyStencil());
128c4762a1bSJed Brown     }
129c4762a1bSJed Brown     catch (char *all) {
1309566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Thrust is not working\n"));
131c4762a1bSJed Brown     }
1329566063dSJacob Faibussowitsch     PetscCall(VecCUDARestoreArrayRead(xlocal,&xarray));
1339566063dSJacob Faibussowitsch     PetscCall(VecCUDARestoreArrayWrite(f,&farray));
134c4762a1bSJed Brown   } else {
1359566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(da,xlocal,&xx));
1369566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArray(da,f,&ff));
1379566063dSJacob Faibussowitsch     PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
138c4762a1bSJed Brown 
139c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
140c4762a1bSJed Brown       if (i == 0 || i == Mx-1) ff[i] = xx[i]/hx;
141c4762a1bSJed Brown       else ff[i] =  (2.0*xx[i] - xx[i-1] - xx[i+1])/hx - hx*PetscExpScalar(xx[i]);
142c4762a1bSJed Brown     }
1439566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(da,xlocal,&xx));
1449566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArray(da,f,&ff));
145c4762a1bSJed Brown   }
1469566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&xlocal));
147c4762a1bSJed Brown   return 0;
148c4762a1bSJed Brown 
149c4762a1bSJed Brown }
150c4762a1bSJed Brown PetscErrorCode ComputeJacobian(SNES snes,Vec x,Mat J,Mat B,void *ctx)
151c4762a1bSJed Brown {
152c4762a1bSJed Brown   DM             da = (DM) ctx;
153c4762a1bSJed Brown   PetscInt       i,Mx,xm,xs;
154c4762a1bSJed Brown   PetscScalar    hx,*xx;
155c4762a1bSJed Brown   Vec            xlocal;
156c4762a1bSJed Brown 
1579566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
158c4762a1bSJed Brown   hx   = 1.0/(PetscReal)(Mx-1);
1599566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&xlocal));
1609566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal));
1619566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal));
1629566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,xlocal,&xx));
1639566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
166c4762a1bSJed Brown     if (i == 0 || i == Mx-1) {
1679566063dSJacob Faibussowitsch       PetscCall(MatSetValue(J,i,i,1.0/hx,INSERT_VALUES));
168c4762a1bSJed Brown     } else {
1699566063dSJacob Faibussowitsch       PetscCall(MatSetValue(J,i,i-1,-1.0/hx,INSERT_VALUES));
1709566063dSJacob Faibussowitsch       PetscCall(MatSetValue(J,i,i,2.0/hx - hx*PetscExpScalar(xx[i]),INSERT_VALUES));
1719566063dSJacob Faibussowitsch       PetscCall(MatSetValue(J,i,i+1,-1.0/hx,INSERT_VALUES));
172c4762a1bSJed Brown     }
173c4762a1bSJed Brown   }
1749566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
1759566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
1769566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,xlocal,&xx));
1779566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&xlocal));
178c4762a1bSJed Brown   return 0;
179c4762a1bSJed Brown }
180c4762a1bSJed Brown 
181c4762a1bSJed Brown /*TEST
182c4762a1bSJed Brown 
183c4762a1bSJed Brown    build:
184c4762a1bSJed Brown       requires: cuda
185c4762a1bSJed Brown 
186499ce956SJunchao Zhang    testset:
187499ce956SJunchao Zhang       args: -snes_monitor_short -dm_mat_type aijcusparse -dm_vec_type cuda
188499ce956SJunchao Zhang       output_file: output/ex47cu_1.out
189c4762a1bSJed Brown       test:
190499ce956SJunchao Zhang         suffix: 1
191499ce956SJunchao Zhang         nsize:  1
192499ce956SJunchao Zhang       test:
193499ce956SJunchao Zhang         suffix: 2
194499ce956SJunchao Zhang         nsize:  2
195c4762a1bSJed Brown 
196c4762a1bSJed Brown TEST*/
197