xref: /petsc/src/mat/tests/ex39.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests Elemental interface.\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscmat.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown int main(int argc,char **args)
7*c4762a1bSJed Brown {
8*c4762a1bSJed Brown   Mat            Cdense,B,C,Ct;
9*c4762a1bSJed Brown   Vec            d;
10*c4762a1bSJed Brown   PetscInt       i,j,m = 5,n,nrows,ncols;
11*c4762a1bSJed Brown   const PetscInt *rows,*cols;
12*c4762a1bSJed Brown   IS             isrows,iscols;
13*c4762a1bSJed Brown   PetscErrorCode ierr;
14*c4762a1bSJed Brown   PetscScalar    *v;
15*c4762a1bSJed Brown   PetscMPIInt    rank,size;
16*c4762a1bSJed Brown   PetscReal      Cnorm;
17*c4762a1bSJed Brown   PetscBool      flg,mats_view=PETSC_FALSE;
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
20*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
21*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
22*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
23*c4762a1bSJed Brown   n    = m;
24*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
25*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL,"-mats_view",&mats_view);CHKERRQ(ierr);
26*c4762a1bSJed Brown 
27*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
28*c4762a1bSJed Brown   ierr = MatSetSizes(C,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
29*c4762a1bSJed Brown   ierr = MatSetType(C,MATELEMENTAL);CHKERRQ(ierr);
30*c4762a1bSJed Brown   ierr = MatSetFromOptions(C);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = MatSetUp(C);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = MatGetOwnershipIS(C,&isrows,&iscols);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr);
34*c4762a1bSJed Brown   ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr);
35*c4762a1bSJed Brown   ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr);
36*c4762a1bSJed Brown   ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr);
37*c4762a1bSJed Brown   ierr = PetscMalloc1(nrows*ncols,&v);CHKERRQ(ierr);
38*c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
39*c4762a1bSJed Brown   PetscRandom rand;
40*c4762a1bSJed Brown   PetscScalar rval;
41*c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr);
43*c4762a1bSJed Brown   for (i=0; i<nrows; i++) {
44*c4762a1bSJed Brown     for (j=0; j<ncols; j++) {
45*c4762a1bSJed Brown       ierr         = PetscRandomGetValue(rand,&rval);CHKERRQ(ierr);
46*c4762a1bSJed Brown       v[i*ncols+j] = rval;
47*c4762a1bSJed Brown     }
48*c4762a1bSJed Brown   }
49*c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
50*c4762a1bSJed Brown #else
51*c4762a1bSJed Brown   for (i=0; i<nrows; i++) {
52*c4762a1bSJed Brown     for (j=0; j<ncols; j++) {
53*c4762a1bSJed Brown       v[i*ncols+j] = (PetscReal)(10000*rank+100*rows[i]+cols[j]);
54*c4762a1bSJed Brown     }
55*c4762a1bSJed Brown   }
56*c4762a1bSJed Brown #endif
57*c4762a1bSJed Brown   ierr = MatSetValues(C,nrows,rows,ncols,cols,v,INSERT_VALUES);CHKERRQ(ierr);
58*c4762a1bSJed Brown   ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr);
59*c4762a1bSJed Brown   ierr = ISRestoreIndices(iscols,&cols);CHKERRQ(ierr);
60*c4762a1bSJed Brown   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
61*c4762a1bSJed Brown   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
62*c4762a1bSJed Brown   ierr = ISDestroy(&isrows);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr = ISDestroy(&iscols);CHKERRQ(ierr);
64*c4762a1bSJed Brown 
65*c4762a1bSJed Brown   /* Test MatView(), MatDuplicate() and out-of-place MatConvert() */
66*c4762a1bSJed Brown   ierr = MatDuplicate(C,MAT_COPY_VALUES,&B);CHKERRQ(ierr);
67*c4762a1bSJed Brown   if (mats_view) {
68*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Duplicated C:\n");CHKERRQ(ierr);
69*c4762a1bSJed Brown     ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
70*c4762a1bSJed Brown   }
71*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
72*c4762a1bSJed Brown   ierr = MatConvert(C,MATMPIDENSE,MAT_INITIAL_MATRIX,&Cdense);CHKERRQ(ierr);
73*c4762a1bSJed Brown   ierr = MatMultEqual(C,Cdense,5,&flg);CHKERRQ(ierr);
74*c4762a1bSJed Brown   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NOTSAMETYPE,"Cdense != C. MatConvert() fails");
75*c4762a1bSJed Brown 
76*c4762a1bSJed Brown   /* Test MatNorm() */
77*c4762a1bSJed Brown   ierr = MatNorm(C,NORM_1,&Cnorm);CHKERRQ(ierr);
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown   /* Test MatTranspose(), MatZeroEntries() and MatGetDiagonal() */
80*c4762a1bSJed Brown   ierr = MatTranspose(C,MAT_INITIAL_MATRIX,&Ct);CHKERRQ(ierr);
81*c4762a1bSJed Brown   ierr = MatConjugate(Ct);CHKERRQ(ierr);
82*c4762a1bSJed Brown   if (mats_view) {
83*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"C's Transpose Conjugate:\n");CHKERRQ(ierr);
84*c4762a1bSJed Brown     ierr = MatView(Ct,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
85*c4762a1bSJed Brown   }
86*c4762a1bSJed Brown   ierr = MatZeroEntries(Ct);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&d);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = VecSetSizes(d,m>n ? n : m,PETSC_DECIDE);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = VecSetFromOptions(d);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = MatGetDiagonal(C,d);CHKERRQ(ierr);
91*c4762a1bSJed Brown   if (mats_view) {
92*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Diagonal of C:\n");CHKERRQ(ierr);
93*c4762a1bSJed Brown     ierr = VecView(d,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
94*c4762a1bSJed Brown   }
95*c4762a1bSJed Brown   if (m>n) {
96*c4762a1bSJed Brown     ierr = MatDiagonalScale(C,NULL,d);CHKERRQ(ierr);
97*c4762a1bSJed Brown   } else {
98*c4762a1bSJed Brown     ierr = MatDiagonalScale(C,d,NULL);CHKERRQ(ierr);
99*c4762a1bSJed Brown   }
100*c4762a1bSJed Brown   if (mats_view) {
101*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Diagonal Scaled C:\n");CHKERRQ(ierr);
102*c4762a1bSJed Brown     ierr = MatView(C,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
103*c4762a1bSJed Brown   }
104*c4762a1bSJed Brown 
105*c4762a1bSJed Brown   /* Test MatAXPY(), MatAYPX() and in-place MatConvert() */
106*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&B);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = MatSetSizes(B,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = MatSetType(B,MATELEMENTAL);CHKERRQ(ierr);
109*c4762a1bSJed Brown   ierr = MatSetFromOptions(B);CHKERRQ(ierr);
110*c4762a1bSJed Brown   ierr = MatSetUp(B);CHKERRQ(ierr);
111*c4762a1bSJed Brown   ierr = MatGetOwnershipIS(B,&isrows,&iscols);CHKERRQ(ierr);
112*c4762a1bSJed Brown   ierr = ISGetLocalSize(isrows,&nrows);CHKERRQ(ierr);
113*c4762a1bSJed Brown   ierr = ISGetIndices(isrows,&rows);CHKERRQ(ierr);
114*c4762a1bSJed Brown   ierr = ISGetLocalSize(iscols,&ncols);CHKERRQ(ierr);
115*c4762a1bSJed Brown   ierr = ISGetIndices(iscols,&cols);CHKERRQ(ierr);
116*c4762a1bSJed Brown   for (i=0; i<nrows; i++) {
117*c4762a1bSJed Brown     for (j=0; j<ncols; j++) {
118*c4762a1bSJed Brown       v[i*ncols+j] = (PetscReal)(1000*rows[i]+cols[j]);
119*c4762a1bSJed Brown     }
120*c4762a1bSJed Brown   }
121*c4762a1bSJed Brown   ierr = MatSetValues(B,nrows,rows,ncols,cols,v,INSERT_VALUES);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = PetscFree(v);CHKERRQ(ierr);
123*c4762a1bSJed Brown   ierr = ISRestoreIndices(isrows,&rows);CHKERRQ(ierr);
124*c4762a1bSJed Brown   ierr = ISRestoreIndices(iscols,&cols);CHKERRQ(ierr);
125*c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
126*c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
127*c4762a1bSJed Brown   ierr = MatAXPY(B,2.5,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
128*c4762a1bSJed Brown   ierr = MatAYPX(B,3.75,C,SAME_NONZERO_PATTERN);CHKERRQ(ierr);
129*c4762a1bSJed Brown   ierr = MatConvert(B,MATDENSE,MAT_INPLACE_MATRIX,&B);CHKERRQ(ierr);
130*c4762a1bSJed Brown   if (mats_view) {
131*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"B after MatAXPY and MatAYPX:\n");CHKERRQ(ierr);
132*c4762a1bSJed Brown     ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
133*c4762a1bSJed Brown   }
134*c4762a1bSJed Brown   ierr = ISDestroy(&isrows);CHKERRQ(ierr);
135*c4762a1bSJed Brown   ierr = ISDestroy(&iscols);CHKERRQ(ierr);
136*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
137*c4762a1bSJed Brown 
138*c4762a1bSJed Brown   /* Test MatMatTransposeMult(): B = C*C^T */
139*c4762a1bSJed Brown   ierr = MatMatTransposeMult(C,C,MAT_INITIAL_MATRIX,PETSC_DEFAULT,&B);CHKERRQ(ierr);
140*c4762a1bSJed Brown   if (mats_view) {
141*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"C MatMatTransposeMult C:\n");CHKERRQ(ierr);
142*c4762a1bSJed Brown     ierr = MatView(B,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
143*c4762a1bSJed Brown   }
144*c4762a1bSJed Brown 
145*c4762a1bSJed Brown   ierr = MatDestroy(&Cdense);CHKERRQ(ierr);
146*c4762a1bSJed Brown   ierr = PetscFree(v);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
148*c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
149*c4762a1bSJed Brown   ierr = MatDestroy(&Ct);CHKERRQ(ierr);
150*c4762a1bSJed Brown   ierr = VecDestroy(&d);CHKERRQ(ierr);
151*c4762a1bSJed Brown   ierr = PetscFinalize();
152*c4762a1bSJed Brown   return ierr;
153*c4762a1bSJed Brown }
154*c4762a1bSJed Brown 
155*c4762a1bSJed Brown 
156*c4762a1bSJed Brown 
157*c4762a1bSJed Brown 
158*c4762a1bSJed Brown /*TEST
159*c4762a1bSJed Brown 
160*c4762a1bSJed Brown    test:
161*c4762a1bSJed Brown       nsize: 2
162*c4762a1bSJed Brown       args: -m 3 -n 2
163*c4762a1bSJed Brown       requires: elemental
164*c4762a1bSJed Brown 
165*c4762a1bSJed Brown    test:
166*c4762a1bSJed Brown       suffix: 2
167*c4762a1bSJed Brown       nsize: 6
168*c4762a1bSJed Brown       args: -m 2 -n 3
169*c4762a1bSJed Brown       requires: elemental
170*c4762a1bSJed Brown 
171*c4762a1bSJed Brown TEST*/
172