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