Actual source code: ex47cu.cu
petsc-3.12.5 2020-03-29
1: static char help[] = "Solves -Laplacian u - exp(u) = 0, 0 < x < 1 using GPU\n\n";
2: /*
3: Same as ex47.c except it also uses the GPU to evaluate the function
4: */
6: #include <petscdm.h>
7: #include <petscdmda.h>
8: #include <petscsnes.h>
10: #include <thrust/device_ptr.h>
11: #include <thrust/for_each.h>
12: #include <thrust/tuple.h>
13: #include <thrust/iterator/constant_iterator.h>
14: #include <thrust/iterator/counting_iterator.h>
15: #include <thrust/iterator/zip_iterator.h>
17: extern PetscErrorCode ComputeFunction(SNES,Vec,Vec,void*), ComputeJacobian(SNES,Vec,Mat,Mat,void*);
18: PetscBool useCUDA = PETSC_FALSE;
20: int main(int argc,char **argv)
21: {
22: SNES snes;
23: Vec x,f;
24: Mat J;
25: DM da;
27: char *tmp,typeName[256];
28: PetscBool flg;
30: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
31: PetscOptionsGetString(NULL,NULL,"-dm_vec_type",typeName,256,&flg);
32: if (flg) {
33: PetscStrstr(typeName,"cuda",&tmp);
34: if (tmp) useCUDA = PETSC_TRUE;
35: }
37: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,8,1,1,NULL,&da);
38: DMSetFromOptions(da);
39: DMSetUp(da);
40: DMCreateGlobalVector(da,&x); VecDuplicate(x,&f);
41: DMSetMatType(da,MATAIJ);
42: DMCreateMatrix(da,&J);
44: SNESCreate(PETSC_COMM_WORLD,&snes);
45: SNESSetFunction(snes,f,ComputeFunction,da);
46: SNESSetJacobian(snes,J,J,ComputeJacobian,da);
47: SNESSetFromOptions(snes);
48: SNESSolve(snes,NULL,x);
50: MatDestroy(&J);
51: VecDestroy(&x);
52: VecDestroy(&f);
53: SNESDestroy(&snes);
54: DMDestroy(&da);
56: PetscFinalize();
57: return ierr;
58: }
60: struct ApplyStencil
61: {
62: template <typename Tuple>
63: __host__ __device__
64: void operator()(Tuple t)
65: {
66: /* f = (2*x_i - x_(i+1) - x_(i-1))/h - h*exp(x_i) */
67: thrust::get<0>(t) = 1;
68: if ((thrust::get<4>(t) > 0) && (thrust::get<4>(t) < thrust::get<5>(t)-1)) {
69: 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));
70: } else if (thrust::get<4>(t) == 0) {
71: thrust::get<0>(t) = thrust::get<1>(t) / (thrust::get<6>(t));
72: } else if (thrust::get<4>(t) == thrust::get<5>(t)-1) {
73: thrust::get<0>(t) = thrust::get<1>(t) / (thrust::get<6>(t));
74: }
75: }
76: };
78: PetscErrorCode ComputeFunction(SNES snes,Vec x,Vec f,void *ctx)
79: {
80: PetscInt i,Mx,xs,xm,xstartshift,xendshift,fstart,lsize;
81: PetscScalar *xx,*ff,hx;
82: DM da = (DM) ctx;
83: Vec xlocal;
84: PetscErrorCode ierr;
85: PetscMPIInt rank,size;
86: MPI_Comm comm;
87: PetscScalar const *xarray;
88: PetscScalar *farray;
90: 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);
91: hx = 1.0/(PetscReal)(Mx-1);
92: DMGetLocalVector(da,&xlocal);
93: DMGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal);
94: DMGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal);
96: if (useCUDA) {
97: VecCUDAGetArrayRead(xlocal,&xarray);
98: VecCUDAGetArrayWrite(f,&farray);
99: PetscObjectGetComm((PetscObject)da,&comm);
100: MPI_Comm_size(comm,&size);
101: MPI_Comm_rank(comm,&rank);
102: if (rank) xstartshift = 1;
103: else xstartshift = 0;
104: if (rank != size-1) xendshift = 1;
105: else xendshift = 0;
106: VecGetOwnershipRange(f,&fstart,NULL);
107: VecGetLocalSize(x,&lsize);
108: try {
109: thrust::for_each(
110: thrust::make_zip_iterator(
111: thrust::make_tuple(
112: thrust::device_ptr<PetscScalar>(farray),
113: thrust::device_ptr<const PetscScalar>(xarray + xstartshift),
114: thrust::device_ptr<const PetscScalar>(xarray + xstartshift + 1),
115: thrust::device_ptr<const PetscScalar>(xarray + xstartshift - 1),
116: thrust::counting_iterator<int>(fstart),
117: thrust::constant_iterator<int>(Mx),
118: thrust::constant_iterator<PetscScalar>(hx))),
119: thrust::make_zip_iterator(
120: thrust::make_tuple(
121: thrust::device_ptr<PetscScalar>(farray + lsize),
122: thrust::device_ptr<const PetscScalar>(xarray + lsize - xendshift),
123: thrust::device_ptr<const PetscScalar>(xarray + lsize - xendshift + 1),
124: thrust::device_ptr<const PetscScalar>(xarray + lsize - xendshift - 1),
125: thrust::counting_iterator<int>(fstart) + lsize,
126: thrust::constant_iterator<int>(Mx),
127: thrust::constant_iterator<PetscScalar>(hx))),
128: ApplyStencil());
129: }
130: catch (char *all) {
131: PetscPrintf(PETSC_COMM_WORLD, "Thrust is not working\n");
132: }
133: VecCUDARestoreArrayRead(xlocal,&xarray);
134: VecCUDARestoreArrayWrite(f,&farray);
135: } else {
136: DMDAVecGetArray(da,xlocal,&xx);
137: DMDAVecGetArray(da,f,&ff);
138: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
140: for (i=xs; i<xs+xm; i++) {
141: if (i == 0 || i == Mx-1) ff[i] = xx[i]/hx;
142: else ff[i] = (2.0*xx[i] - xx[i-1] - xx[i+1])/hx - hx*PetscExpScalar(xx[i]);
143: }
144: DMDAVecRestoreArray(da,xlocal,&xx);
145: DMDAVecRestoreArray(da,f,&ff);
146: }
147: DMRestoreLocalVector(da,&xlocal);
148: // VecView(x,0);printf("f\n");
149: // VecView(f,0);
150: return 0;
152: }
153: PetscErrorCode ComputeJacobian(SNES snes,Vec x,Mat J,Mat B,void *ctx)
154: {
155: DM da = (DM) ctx;
156: PetscInt i,Mx,xm,xs;
157: PetscScalar hx,*xx;
158: Vec xlocal;
161: 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);
162: hx = 1.0/(PetscReal)(Mx-1);
163: DMGetLocalVector(da,&xlocal);DMGlobalToLocalBegin(da,x,INSERT_VALUES,xlocal);
164: DMGlobalToLocalEnd(da,x,INSERT_VALUES,xlocal);
165: DMDAVecGetArray(da,xlocal,&xx);
166: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
168: for (i=xs; i<xs+xm; i++) {
169: if (i == 0 || i == Mx-1) {
170: MatSetValue(J,i,i,1.0/hx,INSERT_VALUES);
171: } else {
172: MatSetValue(J,i,i-1,-1.0/hx,INSERT_VALUES);
173: MatSetValue(J,i,i,2.0/hx - hx*PetscExpScalar(xx[i]),INSERT_VALUES);
174: MatSetValue(J,i,i+1,-1.0/hx,INSERT_VALUES);
175: }
176: }
177: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
178: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
179: DMDAVecRestoreArray(da,xlocal,&xx);
180: DMRestoreLocalVector(da,&xlocal);
181: return 0;
182: }
186: /*TEST
188: build:
189: requires: cuda
191: test:
192: args: -snes_monitor_short -dm_vec_type cuda
194: TEST*/