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