Subversion Repositories slepc-dev

Rev

Go to most recent revision | Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
2387 jroman 1
/*
2
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2575 eromero 4
   Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
2387 jroman 5
 
6
   This file is part of SLEPc.
7
 
8
   SLEPc is free software: you can redistribute it and/or modify it under  the
9
   terms of version 3 of the GNU Lesser General Public License as published by
10
   the Free Software Foundation.
11
 
12
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
13
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
14
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
15
   more details.
16
 
17
   You  should have received a copy of the GNU Lesser General  Public  License
18
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20
*/
21
 
22
static char help[] = "Tests multiple calls to EPSSolve with different matrix.\n\n";
23
 
24
#include <slepceps.h>
25
 
26
#undef __FUNCT__
27
#define __FUNCT__ "main"
28
int main(int argc,char **argv)
29
{
30
  Mat            A1,A2;       /* problem matrices */
31
  EPS            eps;         /* eigenproblem solver context */
2508 jroman 32
  PetscScalar    value[3];
2558 eromero 33
  PetscReal      tol=1000*PETSC_MACHINE_EPSILON,v;
2508 jroman 34
  Vec            d;
35
  PetscInt       n=30,i,Istart,Iend,col[3];
2387 jroman 36
  PetscBool      FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
37
  PetscRandom    myrand;
38
  PetscErrorCode ierr;
39
 
40
  SlepcInitialize(&argc,&argv,(char*)0,help);
41
 
42
  ierr = PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);CHKERRQ(ierr);
2508 jroman 43
  ierr = PetscPrintf(PETSC_COMM_WORLD,"\nTridiagonal with random diagonal, n=%D\n\n",n);CHKERRQ(ierr);
2387 jroman 44
 
45
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46
           Create matrix tridiag([-1 0 -1])
47
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
48
  ierr = MatCreate(PETSC_COMM_WORLD,&A1);CHKERRQ(ierr);
49
  ierr = MatSetSizes(A1,PETSC_DECIDE,PETSC_DECIDE,n,n);CHKERRQ(ierr);
50
  ierr = MatSetFromOptions(A1);CHKERRQ(ierr);
51
 
52
  ierr = MatGetOwnershipRange(A1,&Istart,&Iend);CHKERRQ(ierr);
53
  if (Istart==0) FirstBlock=PETSC_TRUE;
54
  if (Iend==n) LastBlock=PETSC_TRUE;
55
  value[0]=-1.0; value[1]=0.0; value[2]=-1.0;
56
  for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
57
    col[0]=i-1; col[1]=i; col[2]=i+1;
58
    ierr = MatSetValues(A1,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr);
59
  }
60
  if (LastBlock) {
61
    i=n-1; col[0]=n-2; col[1]=n-1;
62
    ierr = MatSetValues(A1,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
63
  }
64
  if (FirstBlock) {
65
    i=0; col[0]=0; col[1]=1; value[0]=0.0; value[1]=-1.0;
66
    ierr = MatSetValues(A1,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
67
  }
68
 
69
  ierr = MatAssemblyBegin(A1,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
70
  ierr = MatAssemblyEnd(A1,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
71
 
72
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73
       Create two matrices by filling the diagonal with rand values
74
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75
  ierr = MatDuplicate(A1,MAT_COPY_VALUES,&A2);CHKERRQ(ierr);
76
  ierr = MatGetVecs(A1,PETSC_NULL,&d);CHKERRQ(ierr);
77
  ierr = PetscRandomCreate(PETSC_COMM_WORLD,&myrand);CHKERRQ(ierr);
78
  ierr = PetscRandomSetFromOptions(myrand);CHKERRQ(ierr);
2558 eromero 79
  ierr = PetscRandomSetInterval(myrand,0.0,1.0);CHKERRQ(ierr);
80
  for (i=0; i<n; i++) {
81
    ierr = PetscRandomGetValueReal(myrand,&v);CHKERRQ(ierr);
82
    ierr = VecSetValue(d,i,v,INSERT_VALUES);CHKERRQ(ierr);
83
  }
84
  ierr = VecAssemblyBegin(d);CHKERRQ(ierr);
85
  ierr = VecAssemblyEnd(d);CHKERRQ(ierr);
2387 jroman 86
  ierr = MatDiagonalSet(A1,d,INSERT_VALUES);CHKERRQ(ierr);
2558 eromero 87
  for (i=0; i<n; i++) {
88
    ierr = PetscRandomGetValueReal(myrand,&v);CHKERRQ(ierr);
89
    ierr = VecSetValue(d,i,v,INSERT_VALUES);CHKERRQ(ierr);
90
  }
91
  ierr = VecAssemblyBegin(d);CHKERRQ(ierr);
92
  ierr = VecAssemblyEnd(d);CHKERRQ(ierr);
2387 jroman 93
  ierr = MatDiagonalSet(A2,d,INSERT_VALUES);CHKERRQ(ierr);
94
  ierr = VecDestroy(&d);CHKERRQ(ierr);
95
  ierr = PetscRandomDestroy(&myrand);CHKERRQ(ierr);
96
 
97
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98
                        Create the eigensolver
99
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100
  ierr = EPSCreate(PETSC_COMM_WORLD,&eps);CHKERRQ(ierr);
101
  ierr = EPSSetProblemType(eps,EPS_HEP);CHKERRQ(ierr);
2557 eromero 102
  ierr = EPSSetTolerances(eps,tol,PETSC_DECIDE);CHKERRQ(ierr);
2387 jroman 103
  ierr = EPSSetOperators(eps,A1,PETSC_NULL);CHKERRQ(ierr);
104
  ierr = EPSSetFromOptions(eps);CHKERRQ(ierr);
105
 
106
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107
                        Solve first eigenproblem
108
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109
  ierr = EPSSolve(eps);CHKERRQ(ierr);
2508 jroman 110
  ierr = PetscPrintf(PETSC_COMM_WORLD," - - - First matrix - - -\n");CHKERRQ(ierr);
111
  ierr = EPSPrintSolution(eps,PETSC_NULL);CHKERRQ(ierr);
2387 jroman 112
 
113
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114
                        Solve second eigenproblem
115
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116
  ierr = EPSSetOperators(eps,A2,PETSC_NULL);CHKERRQ(ierr);
117
  ierr = EPSSolve(eps);CHKERRQ(ierr);
2508 jroman 118
  ierr = PetscPrintf(PETSC_COMM_WORLD," - - - Second matrix - - -\n");CHKERRQ(ierr);
119
  ierr = EPSPrintSolution(eps,PETSC_NULL);CHKERRQ(ierr);
2387 jroman 120
 
121
  ierr = EPSDestroy(&eps);CHKERRQ(ierr);
122
  ierr = MatDestroy(&A1);CHKERRQ(ierr);
123
  ierr = MatDestroy(&A2);CHKERRQ(ierr);
124
  ierr = SlepcFinalize();CHKERRQ(ierr);
125
  return 0;
126
}