Actual source code: ex58.c
2: static char help[] = "Test VecCreate{Seq|MPI}CUDAWithArrays.\n\n";
4: #include "petsc.h"
6: int main(int argc, char **argv)
7: {
8: Vec x, y, z;
9: PetscMPIInt size;
10: PetscInt n = 5;
11: PetscScalar xHost[5] = {0., 1., 2., 3., 4.};
12: PetscScalar zHost[5];
13: PetscBool equal;
16: PetscInitialize(&argc, &argv, (char *)0, help);
17: MPI_Comm_size(PETSC_COMM_WORLD, &size);
19: PetscArraycpy(zHost, xHost, n);
20: VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, n, PETSC_DECIDE, zHost, &z); /* build z for comparison */
22: if (size == 1) VecCreateSeqCUDAWithArrays(PETSC_COMM_WORLD, 1, n, xHost, NULL, &x);
23: else VecCreateMPICUDAWithArrays(PETSC_COMM_WORLD, 1, n, PETSC_DECIDE, xHost, NULL, &x);
25: VecEqual(z, x, &equal);
28: VecSet(x, 42.0);
29: VecSet(z, 42.0);
30: VecEqual(z, x, &equal);
33: VecCreateMPIWithArray(PETSC_COMM_WORLD, 1, n, PETSC_DECIDE, xHost, &y);
34: VecSetFromOptions(y); /* changing y's type should not lose its value */
35: VecEqual(z, y, &equal);
38: VecDestroy(&y);
39: VecDestroy(&x);
40: VecDestroy(&z);
41: PetscFinalize();
42: return 0;
43: }
45: /*TEST
47: build:
48: requires: cuda
50: testset:
51: output_file: output/empty.out
52: nsize: {{1 2}}
54: test:
55: suffix: y_host
57: test:
58: TODO: we need something like VecConvert()
59: requires: kokkos_kernels
60: suffix: y_dev
61: args: -vec_type {{standard mpi cuda kokkos}}
62: TEST*/