Actual source code: test1.c
slepc-3.21.1 2024-04-26
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: static char help[] = "Test MatCreateTile.\n\n";
13: #include <slepcsys.h>
15: int main(int argc,char **argv)
16: {
17: Mat T,E,A;
18: PetscInt i,Istart,Iend,n=10;
20: PetscFunctionBeginUser;
21: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
22: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
23: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"MatCreateTile test, n=%" PetscInt_FMT "\n",n));
25: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
26: Create T=tridiag([-1 2 -1],n,n) and E=eye(n)
27: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
29: PetscCall(MatCreate(PETSC_COMM_WORLD,&T));
30: PetscCall(MatSetSizes(T,PETSC_DECIDE,PETSC_DECIDE,n,n));
31: PetscCall(MatSetFromOptions(T));
33: PetscCall(MatGetOwnershipRange(T,&Istart,&Iend));
34: for (i=Istart;i<Iend;i++) {
35: if (i>0) PetscCall(MatSetValue(T,i,i-1,-1.0,INSERT_VALUES));
36: if (i<n-1) PetscCall(MatSetValue(T,i,i+1,-1.0,INSERT_VALUES));
37: PetscCall(MatSetValue(T,i,i,2.0,INSERT_VALUES));
38: }
39: PetscCall(MatAssemblyBegin(T,MAT_FINAL_ASSEMBLY));
40: PetscCall(MatAssemblyEnd(T,MAT_FINAL_ASSEMBLY));
42: PetscCall(MatCreate(PETSC_COMM_WORLD,&E));
43: PetscCall(MatSetSizes(E,PETSC_DECIDE,PETSC_DECIDE,n,n));
44: PetscCall(MatSetFromOptions(E));
46: PetscCall(MatGetOwnershipRange(E,&Istart,&Iend));
47: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(E,i,i,1.0,INSERT_VALUES));
48: PetscCall(MatAssemblyBegin(E,MAT_FINAL_ASSEMBLY));
49: PetscCall(MatAssemblyEnd(E,MAT_FINAL_ASSEMBLY));
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Create tiled matrix A = [ 2*T -E; 0 3*T ]
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54: PetscCall(MatCreateTile(2.0,T,-1.0,E,0.0,E,3.0,T,&A));
55: PetscCall(MatView(A,NULL));
57: PetscCall(MatDestroy(&T));
58: PetscCall(MatDestroy(&E));
59: PetscCall(MatDestroy(&A));
60: PetscCall(SlepcFinalize());
61: return 0;
62: }
64: /*TEST
66: test:
67: suffix: 1
68: nsize: 1
70: test:
71: suffix: 2
72: nsize: 2
74: TEST*/