xref: /petsc/src/vec/vec/tests/ex43.c (revision df2db7a0bab43fc882bb72877d20f2a85250a807)
1 static char help[] = "Tests VecMDot(),VecDot(),VecMTDot(), and VecTDot()\n";
2 
3 #include <petscvec.h>
4 
5 int main(int argc, char **argv)
6 {
7   Vec            *V,t;
8   PetscInt       i,j,reps,n=15,k=6;
9   PetscRandom    rctx;
10   PetscScalar    *val_dot,*val_mdot,*tval_dot,*tval_mdot;
11 
12   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
13   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
14   PetscCall(PetscOptionsGetInt(NULL,NULL,"-k",&k,NULL));
15   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test with %" PetscInt_FMT " random vectors of length %" PetscInt_FMT "\n",k,n));
16   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD,&rctx));
17   PetscCall(PetscRandomSetFromOptions(rctx));
18 #if defined(PETSC_USE_COMPLEX)
19   PetscCall(PetscRandomSetInterval(rctx,-1.+4.*PETSC_i,1.+5.*PETSC_i));
20 #else
21   PetscCall(PetscRandomSetInterval(rctx,-1.,1.));
22 #endif
23   PetscCall(VecCreate(PETSC_COMM_WORLD,&t));
24   PetscCall(VecSetSizes(t,n,PETSC_DECIDE));
25   PetscCall(VecSetFromOptions(t));
26   PetscCall(VecDuplicateVecs(t,k,&V));
27   PetscCall(VecSetRandom(t,rctx));
28   PetscCall(VecViewFromOptions(t,NULL,"-t_view"));
29   PetscCall(PetscMalloc1(k,&val_dot));
30   PetscCall(PetscMalloc1(k,&val_mdot));
31   PetscCall(PetscMalloc1(k,&tval_dot));
32   PetscCall(PetscMalloc1(k,&tval_mdot));
33   for (i=0; i<k; i++) PetscCall(VecSetRandom(V[i],rctx));
34   for (reps=0; reps<20; reps++) {
35     for (i=1; i<k; i++) {
36       PetscCall(VecMDot(t,i,V,val_mdot));
37       PetscCall(VecMTDot(t,i,V,tval_mdot));
38       for (j=0;j<i;j++) {
39         PetscCall(VecDot(t,V[j],&val_dot[j]));
40         PetscCall(VecTDot(t,V[j],&tval_dot[j]));
41       }
42       /* Check result */
43       for (j=0;j<i;j++) {
44         if (PetscAbsScalar(val_mdot[j] - val_dot[j])/PetscAbsScalar(val_dot[j]) > 1e-5) {
45           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[TEST FAILED] i=%" PetscInt_FMT ", j=%" PetscInt_FMT ", val_mdot[j]=%g, val_dot[j]=%g\n",i,j,(double)PetscAbsScalar(val_mdot[j]), (double)PetscAbsScalar(val_dot[j])));
46           break;
47         }
48         if (PetscAbsScalar(tval_mdot[j] - tval_dot[j])/PetscAbsScalar(tval_dot[j]) > 1e-5) {
49           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[TEST FAILED] i=%" PetscInt_FMT ", j=%" PetscInt_FMT ", tval_mdot[j]=%g, tval_dot[j]=%g\n",i,j,(double)PetscAbsScalar(tval_mdot[j]), (double)PetscAbsScalar(tval_dot[j])));
50           break;
51         }
52       }
53     }
54   }
55   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Test completed successfully!\n"));
56   PetscCall(PetscFree(val_dot));
57   PetscCall(PetscFree(val_mdot));
58   PetscCall(PetscFree(tval_dot));
59   PetscCall(PetscFree(tval_mdot));
60   PetscCall(VecDestroyVecs(k,&V));
61   PetscCall(VecDestroy(&t));
62   PetscCall(PetscRandomDestroy(&rctx));
63   PetscCall(PetscFinalize());
64   return 0;
65 }
66 
67 /*TEST
68 
69    test:
70 
71    testset:
72       output_file: output/ex43_1.out
73 
74       test:
75          suffix: cuda
76          args: -vec_type cuda -random_type curand
77          requires: cuda
78 
79       test:
80          suffix: kokkos
81          args: -vec_type kokkos
82          requires: kokkos_kernels
83 
84       test:
85          suffix: hip
86          args: -vec_type hip
87          requires: hip
88 TEST*/
89