Actual source code: test4.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  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.

 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.

 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: */

 22: static char help[] = "Test ST with four matrices.\n\n";

 24: #include <slepcst.h>

 28: int main(int argc,char **argv)
 29: {
 30:   Mat            A,B,C,D,mat[4];
 31:   ST             st;
 32:   Vec            v,w;
 33:   STType         type;
 34:   PetscScalar    value[3],sigma;
 35:   PetscInt       n=10,i,Istart,Iend,col[3];
 36:   PetscBool      FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;

 39:   SlepcInitialize(&argc,&argv,(char*)0,help);
 40:   PetscOptionsGetInt(NULL,"-n",&n,NULL);
 41:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Laplacian plus diagonal, n=%D\n\n",n);
 42:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 43:      Compute the operator matrix for the 1-D Laplacian
 44:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 46:   MatCreate(PETSC_COMM_WORLD,&A);
 47:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
 48:   MatSetFromOptions(A);
 49:   MatSetUp(A);

 51:   MatCreate(PETSC_COMM_WORLD,&B);
 52:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
 53:   MatSetFromOptions(B);
 54:   MatSetUp(B);

 56:   MatCreate(PETSC_COMM_WORLD,&C);
 57:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
 58:   MatSetFromOptions(C);
 59:   MatSetUp(C);

 61:   MatCreate(PETSC_COMM_WORLD,&D);
 62:   MatSetSizes(D,PETSC_DECIDE,PETSC_DECIDE,n,n);
 63:   MatSetFromOptions(D);
 64:   MatSetUp(D);

 66:   MatGetOwnershipRange(A,&Istart,&Iend);
 67:   if (Istart==0) FirstBlock=PETSC_TRUE;
 68:   if (Iend==n) LastBlock=PETSC_TRUE;
 69:   value[0]=-1.0; value[1]=2.0; value[2]=-1.0;
 70:   for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
 71:     col[0]=i-1; col[1]=i; col[2]=i+1;
 72:     MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);
 73:     MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES);
 74:   }
 75:   if (LastBlock) {
 76:     i=n-1; col[0]=n-2; col[1]=n-1;
 77:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 78:     MatSetValue(B,i,i,(PetscScalar)i,INSERT_VALUES);
 79:   }
 80:   if (FirstBlock) {
 81:     i=0; col[0]=0; col[1]=1; value[0]=2.0; value[1]=-1.0;
 82:     MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);
 83:     MatSetValue(B,i,i,-1.0,INSERT_VALUES);
 84:   }
 85:   for (i=Istart;i<Iend;i++) {
 86:     MatSetValue(C,i,n-i-1,1.0,INSERT_VALUES);
 87:     MatSetValue(D,i,i,i*.1,INSERT_VALUES);
 88:     if (i==0) {
 89:       MatSetValue(D,0,n-1,1.0,INSERT_VALUES);
 90:     }
 91:     if (i==n-1) {
 92:       MatSetValue(D,n-1,0,1.0,INSERT_VALUES);
 93:     }
 94:   }

 96:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 97:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 98:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 99:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
100:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
101:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
102:   MatAssemblyBegin(D,MAT_FINAL_ASSEMBLY);
103:   MatAssemblyEnd(D,MAT_FINAL_ASSEMBLY);
104:   MatGetVecs(A,&v,&w);
105:   VecSet(v,1.0);

107:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108:                 Create the spectral transformation object
109:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110:   STCreate(PETSC_COMM_WORLD,&st);
111:   mat[0] = A;
112:   mat[1] = B;
113:   mat[2] = C;
114:   mat[3] = D;
115:   STSetOperators(st,4,mat);
116:   STSetFromOptions(st);
117:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118:               Apply the transformed operator for several ST's
119:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
120:   /* shift, sigma=0.0 */
121:   STSetUp(st);
122:   STGetType(st,&type);
123:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
124:   for (i=0;i<4;i++) {
125:     STMatMult(st,i,v,w);
126:     PetscPrintf(PETSC_COMM_WORLD,"k= %D\n",i);
127:     VecView(w,NULL);
128:   }
129:   STMatSolve(st,v,w);
130:   PetscPrintf(PETSC_COMM_WORLD,"solve\n");
131:   VecView(w,NULL);

133:   /* shift, sigma=0.1 */
134:   sigma = 0.1;
135:   STSetShift(st,sigma);
136:   STGetShift(st,&sigma);
137:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
138:   for (i=0;i<4;i++) {
139:     STMatMult(st,i,v,w);
140:     PetscPrintf(PETSC_COMM_WORLD,"k= %D\n",i);
141:     VecView(w,NULL);
142:   }
143:   STMatSolve(st,v,w);
144:   PetscPrintf(PETSC_COMM_WORLD,"solve\n");
145:   VecView(w,NULL);

147:   /* sinvert, sigma=0.1 */
148:   STPostSolve(st);
149:   STSetType(st,STSINVERT);
150:   STGetType(st,&type);
151:   PetscPrintf(PETSC_COMM_WORLD,"ST type %s\n",type);
152:   STGetShift(st,&sigma);
153:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
154:   for (i=0;i<4;i++) {
155:     STMatMult(st,i,v,w);
156:     PetscPrintf(PETSC_COMM_WORLD,"k= %D\n",i);
157:     VecView(w,NULL);
158:   }
159:   STMatSolve(st,v,w);
160:   PetscPrintf(PETSC_COMM_WORLD,"solve\n");
161:   VecView(w,NULL);

163:   /* sinvert, sigma=-0.5 */
164:   sigma = -0.5;
165:   STSetShift(st,sigma);
166:   STGetShift(st,&sigma);
167:   PetscPrintf(PETSC_COMM_WORLD,"With shift=%g\n",(double)PetscRealPart(sigma));
168:   for (i=0;i<4;i++) {
169:     STMatMult(st,i,v,w);
170:     PetscPrintf(PETSC_COMM_WORLD,"k= %D\n",i);
171:     VecView(w,NULL);
172:   }
173:   STMatSolve(st,v,w);
174:   PetscPrintf(PETSC_COMM_WORLD,"solve\n");
175:   VecView(w,NULL);
176:   STDestroy(&st);
177:   MatDestroy(&A);
178:   MatDestroy(&B);
179:   MatDestroy(&C);
180:   MatDestroy(&D);
181:   VecDestroy(&v);
182:   VecDestroy(&w);
183:   SlepcFinalize();
184:   return 0;
185: }