xref: /petsc/src/snes/tutorials/ex47cu.cu (revision 499ce9560a63eddc2155c33efaa488a4e358d0ef)
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   PetscErrorCode ierr;
27c4762a1bSJed Brown   char           *tmp,typeName[256];
28c4762a1bSJed Brown   PetscBool      flg;
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
31c4762a1bSJed Brown   ierr = PetscOptionsGetString(NULL,NULL,"-dm_vec_type",typeName,256,&flg);CHKERRQ(ierr);
32c4762a1bSJed Brown   if (flg) {
33c4762a1bSJed Brown     ierr = PetscStrstr(typeName,"cuda",&tmp);CHKERRQ(ierr);
34c4762a1bSJed Brown     if (tmp) useCUDA = PETSC_TRUE;
35c4762a1bSJed Brown   }
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,8,1,1,NULL,&da);CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
39c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
40c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&x); VecDuplicate(x,&f);CHKERRQ(ierr);
41c4762a1bSJed Brown   ierr = DMCreateMatrix(da,&J);CHKERRQ(ierr);
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
44c4762a1bSJed Brown   ierr = SNESSetFunction(snes,f,ComputeFunction,da);CHKERRQ(ierr);
45c4762a1bSJed Brown   ierr = SNESSetJacobian(snes,J,J,ComputeJacobian,da);CHKERRQ(ierr);
46c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
47c4762a1bSJed Brown   ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr);
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
50c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
51c4762a1bSJed Brown   ierr = VecDestroy(&f);CHKERRQ(ierr);
52c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
53c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   ierr = PetscFinalize();
56c4762a1bSJed Brown   return ierr;
57c4762a1bSJed Brown }
58c4762a1bSJed Brown 
59c4762a1bSJed Brown struct ApplyStencil
60c4762a1bSJed Brown {
61c4762a1bSJed Brown   template <typename Tuple>
62c4762a1bSJed Brown   __host__ __device__
63c4762a1bSJed Brown   void operator()(Tuple t)
64c4762a1bSJed Brown   {
65c4762a1bSJed Brown     /* f = (2*x_i - x_(i+1) - x_(i-1))/h - h*exp(x_i) */
66c4762a1bSJed Brown     thrust::get<0>(t) = 1;
67c4762a1bSJed Brown     if ((thrust::get<4>(t) > 0) && (thrust::get<4>(t) < thrust::get<5>(t)-1)) {
68c4762a1bSJed Brown       thrust::get<0>(t) = (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));
69c4762a1bSJed Brown     } else if (thrust::get<4>(t) == 0) {
70c4762a1bSJed Brown       thrust::get<0>(t) = thrust::get<1>(t) / (thrust::get<6>(t));
71c4762a1bSJed Brown     } else if (thrust::get<4>(t) == thrust::get<5>(t)-1) {
72c4762a1bSJed Brown       thrust::get<0>(t) = thrust::get<1>(t) / (thrust::get<6>(t));
73c4762a1bSJed Brown     }
74c4762a1bSJed Brown   }
75c4762a1bSJed Brown };
76c4762a1bSJed Brown 
77c4762a1bSJed Brown PetscErrorCode ComputeFunction(SNES snes,Vec x,Vec f,void *ctx)
78c4762a1bSJed Brown {
79c4762a1bSJed Brown   PetscInt          i,Mx,xs,xm,xstartshift,xendshift,fstart,lsize;
80c4762a1bSJed Brown   PetscScalar       *xx,*ff,hx;
81c4762a1bSJed Brown   DM                da = (DM) ctx;
82c4762a1bSJed Brown   Vec               xlocal;
83c4762a1bSJed Brown   PetscErrorCode    ierr;
84c4762a1bSJed Brown   PetscMPIInt       rank,size;
85c4762a1bSJed Brown   MPI_Comm          comm;
86c4762a1bSJed Brown   PetscScalar const *xarray;
87c4762a1bSJed Brown   PetscScalar       *farray;
88c4762a1bSJed Brown 
89c4762a1bSJed Brown   ierr = 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);CHKERRQ(ierr);
90c4762a1bSJed Brown   hx   = 1.0/(PetscReal)(Mx-1);
91c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&xlocal);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal);CHKERRQ(ierr);
93c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal);CHKERRQ(ierr);
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   if (useCUDA) {
96c4762a1bSJed Brown     ierr = VecCUDAGetArrayRead(xlocal,&xarray);CHKERRQ(ierr);
97c4762a1bSJed Brown     ierr = VecCUDAGetArrayWrite(f,&farray);CHKERRQ(ierr);
98c4762a1bSJed Brown     ierr = PetscObjectGetComm((PetscObject)da,&comm);CHKERRQ(ierr);
99c4762a1bSJed Brown     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
100c4762a1bSJed Brown     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
101c4762a1bSJed Brown     if (rank) xstartshift = 1;
102c4762a1bSJed Brown     else xstartshift = 0;
103c4762a1bSJed Brown     if (rank != size-1) xendshift = 1;
104c4762a1bSJed Brown     else xendshift = 0;
105c4762a1bSJed Brown     ierr = VecGetOwnershipRange(f,&fstart,NULL);CHKERRQ(ierr);
106c4762a1bSJed Brown     ierr = VecGetLocalSize(x,&lsize);CHKERRQ(ierr);
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) {
130c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD, "Thrust is not working\n");CHKERRQ(ierr);
131c4762a1bSJed Brown     }
132c4762a1bSJed Brown     ierr = VecCUDARestoreArrayRead(xlocal,&xarray);CHKERRQ(ierr);
133c4762a1bSJed Brown     ierr = VecCUDARestoreArrayWrite(f,&farray);CHKERRQ(ierr);
134c4762a1bSJed Brown   } else {
135c4762a1bSJed Brown     ierr = DMDAVecGetArray(da,xlocal,&xx);CHKERRQ(ierr);
136c4762a1bSJed Brown     ierr = DMDAVecGetArray(da,f,&ff);CHKERRQ(ierr);
137c4762a1bSJed Brown     ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
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     }
143c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(da,xlocal,&xx);CHKERRQ(ierr);
144c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(da,f,&ff);CHKERRQ(ierr);
145c4762a1bSJed Brown   }
146c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&xlocal);CHKERRQ(ierr);
147c4762a1bSJed Brown   //  VecView(x,0);printf("f\n");
148c4762a1bSJed Brown   //  VecView(f,0);
149c4762a1bSJed Brown   return 0;
150c4762a1bSJed Brown 
151c4762a1bSJed Brown }
152c4762a1bSJed Brown PetscErrorCode ComputeJacobian(SNES snes,Vec x,Mat J,Mat B,void *ctx)
153c4762a1bSJed Brown {
154c4762a1bSJed Brown   DM             da = (DM) ctx;
155c4762a1bSJed Brown   PetscInt       i,Mx,xm,xs;
156c4762a1bSJed Brown   PetscScalar    hx,*xx;
157c4762a1bSJed Brown   Vec            xlocal;
158c4762a1bSJed Brown   PetscErrorCode ierr;
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   ierr = 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);CHKERRQ(ierr);
161c4762a1bSJed Brown   hx   = 1.0/(PetscReal)(Mx-1);
162c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&xlocal);DMGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal);CHKERRQ(ierr);
163c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal);CHKERRQ(ierr);
164c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,xlocal,&xx);CHKERRQ(ierr);
165c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
168c4762a1bSJed Brown     if (i == 0 || i == Mx-1) {
169c4762a1bSJed Brown       ierr = MatSetValue(J,i,i,1.0/hx,INSERT_VALUES);CHKERRQ(ierr);
170c4762a1bSJed Brown     } else {
171c4762a1bSJed Brown       ierr = MatSetValue(J,i,i-1,-1.0/hx,INSERT_VALUES);CHKERRQ(ierr);
172c4762a1bSJed Brown       ierr = MatSetValue(J,i,i,2.0/hx - hx*PetscExpScalar(xx[i]),INSERT_VALUES);CHKERRQ(ierr);
173c4762a1bSJed Brown       ierr = MatSetValue(J,i,i+1,-1.0/hx,INSERT_VALUES);CHKERRQ(ierr);
174c4762a1bSJed Brown     }
175c4762a1bSJed Brown   }
176c4762a1bSJed Brown   ierr  = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
177c4762a1bSJed Brown   ierr  = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
178c4762a1bSJed Brown   ierr  = DMDAVecRestoreArray(da,xlocal,&xx);CHKERRQ(ierr);
179c4762a1bSJed Brown   ierr  = DMRestoreLocalVector(da,&xlocal);CHKERRQ(ierr);
180c4762a1bSJed Brown   return 0;
181c4762a1bSJed Brown }
182c4762a1bSJed Brown 
183c4762a1bSJed Brown 
184c4762a1bSJed Brown 
185c4762a1bSJed Brown /*TEST
186c4762a1bSJed Brown 
187c4762a1bSJed Brown    build:
188c4762a1bSJed Brown       requires: cuda
189c4762a1bSJed Brown 
190*499ce956SJunchao Zhang    testset:
191*499ce956SJunchao Zhang       args: -snes_monitor_short -dm_mat_type aijcusparse -dm_vec_type cuda
192*499ce956SJunchao Zhang       output_file: output/ex47cu_1.out
193c4762a1bSJed Brown       test:
194*499ce956SJunchao Zhang         suffix: 1
195*499ce956SJunchao Zhang         nsize:  1
196*499ce956SJunchao Zhang       test:
197*499ce956SJunchao Zhang         suffix: 2
198*499ce956SJunchao Zhang         nsize:  2
199c4762a1bSJed Brown 
200c4762a1bSJed Brown TEST*/
201