Actual source code: ex21.c
2: static char help[] = "Tests VecMax() with index.\n\
3: -n <length> : vector length\n\n";
5: #include <petscvec.h>
7: int main(int argc,char **argv)
8: {
10: PetscInt n = 5,idx;
11: PetscReal value,value2;
12: Vec x;
13: PetscScalar one = 1.0;
15: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
16: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
18: /* create vector */
19: VecCreate(PETSC_COMM_WORLD,&x);
20: VecSetSizes(x,PETSC_DECIDE,n);
21: VecSetFromOptions(x);
23: VecSet(x,one);
24: VecSetValue(x,0,0.0,INSERT_VALUES);
25: VecSetValue(x,n-1,2.0,INSERT_VALUES);
26: VecAssemblyBegin(x);
27: VecAssemblyEnd(x);
28: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
29: VecMax(x,&idx,&value);
30: VecMax(x,NULL,&value2);
31: PetscPrintf(PETSC_COMM_WORLD,"Maximum value %g index %D (no index %g)\n",(double)value,idx,(double)value2);
32: VecMin(x,&idx,&value);
33: VecMin(x,NULL,&value2);
34: PetscPrintf(PETSC_COMM_WORLD,"Minimum value %g index %D (no index %g)\n",(double)value,idx,(double)value2);
36: VecDestroy(&x);
38: PetscFinalize();
39: return ierr;
40: }
42: /*TEST
44: testset:
45: diff_args: -j
46: filter: grep -v type | grep -v "MPI processes" | grep -v Process
47: output_file: output/ex21_1.out
49: test:
50: suffix: 1
51: args: -vec_type {{seq mpi}}
53: test:
54: requires: cuda
55: suffix: 1_cuda
56: args: -vec_type {{cuda mpicuda}}
58: test:
59: requires: kokkos_kernels
60: suffix: 1_kokkos
61: args: -vec_type {{kokkos mpikokkos}}
63: test:
64: requires: hip
65: suffix: 1_hip
66: args: -vec_type {{hip mpihip}}
68: testset:
69: diff_args: -j
70: filter: grep -v type
71: output_file: output/ex21_2.out
72: nsize: 2
74: test:
75: suffix: 2
77: test:
78: requires: cuda
79: suffix: 2_cuda
80: args: -vec_type cuda
82: test:
83: requires: kokkos_kernels
84: suffix: 2_kokkos
85: args: -vec_type kokkos
87: test:
88: requires: hip
89: suffix: 2_hip
90: args: -vec_type hip
92: TEST*/