Actual source code: peprefine.c
slepc-3.5.2 2014-10-10
1: /*
2: Newton refinement for PEP, simple version.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc-private/pepimpl.h>
25: #include <slepcblaslapack.h>
27: #define NREF_MAXIT 100
31: static PetscErrorCode PEPEvaluateFunctionDerivatives(PEP pep,PetscScalar alpha,PetscScalar *vals)
32: {
33: PetscInt i,nmat=pep->nmat;
34: PetscScalar a0,a1,a2;
35: PetscReal *a=pep->pbc,*b=a+nmat,*g=b+nmat;
38: a0 = 0.0;
39: a1 = 1.0;
40: vals[0] = 0.0;
41: if (nmat>1) vals[1] = 1/a[0];
42: for (i=2;i<nmat;i++) {
43: a2 = ((alpha-b[i-2])*a1-g[i-2]*a0)/a[i-2];
44: vals[i] = (a2+(alpha-b[i-1])*vals[i-1]-g[i-1]*vals[i-2])/a[i-1];
45: a0 = a1; a1 = a2;
46: }
47: return(0);
48: }
52: PetscErrorCode PEPNSimpleRefSetUp(PEP pep,Mat *A,PetscInt idx,Mat *M,Mat *T,PetscBool ini,Vec *t)
53: {
54: PetscErrorCode ierr;
55: PetscInt i,nmat=pep->nmat,ml,m0,m1,mg;
56: PetscInt *dnz,*onz,ncols,*cols2,*nnz;
57: PetscScalar *array,zero=0.0,*coeffs;
58: PetscMPIInt rank,size;
59: MPI_Comm comm;
60: const PetscInt *cols;
61: const PetscScalar *vals;
62: MatStructure str;
63: Vec v,w=t[1],q=t[0];
66: comm = PetscObjectComm((PetscObject)A[0]);
67: STGetMatStructure(pep->st,&str);
68: PetscMalloc1(nmat,&coeffs);
69: MPI_Comm_rank(comm,&rank);
70: MPI_Comm_size(comm,&size);
71: if (ini) {
72: MatDuplicate(A[0],MAT_COPY_VALUES,T);
73: } else {
74: MatCopy(A[0],*T,SUBSET_NONZERO_PATTERN);
75: }
76: PEPEvaluateBasis(pep,pep->eigr[idx],0,coeffs,NULL);
77: if (coeffs[0]!=1.0) {
78: MatScale(*T,coeffs[0]);
79: }
80: for (i=1;i<nmat;i++) {
81: MatAXPY(*T,coeffs[i],A[i],(ini)?str:SUBSET_NONZERO_PATTERN);
82: }
83: MatGetSize(*T,&mg,NULL);
84: MatGetOwnershipRange(*T,&m0,&m1);
85: if (ini) {
86: MatCreate(comm,M);
87: MatGetLocalSize(*T,&ml,NULL);
88: if (rank==size-1) ml++;
89: MatSetSizes(*M,ml,ml,mg+1,mg+1);
90: MatSetFromOptions(*M);
91: MatSetUp(*M);
92: /* Preallocate M */
93: if (size>1) {
94: MatPreallocateInitialize(comm,ml,ml,dnz,onz);
95: for (i=m0;i<m1;i++) {
96: MatGetRow(*T,i,&ncols,&cols,NULL);
97: MatPreallocateSet(i,ncols,cols,dnz,onz);
98: MatPreallocateSet(i,1,&mg,dnz,onz);
99: MatRestoreRow(*T,i,&ncols,&cols,NULL);
100: }
101: if (rank==size-1) {
102: PetscCalloc1(mg+1,&cols2);
103: for (i=0;i<mg+1;i++) cols2[i]=i;
104: MatPreallocateSet(m1,mg+1,cols2,dnz,onz);
105: PetscFree(cols2);
106: }
107: MatMPIAIJSetPreallocation(*M,0,dnz,0,onz);
108: MatPreallocateFinalize(dnz,onz);
109: } else {
110: PetscCalloc1(mg+1,&nnz);
111: for (i=0;i<mg;i++) {
112: MatGetRow(*T,i,&ncols,NULL,NULL);
113: nnz[i] = ncols+1;
114: MatRestoreRow(*T,i,&ncols,NULL,NULL);
115: }
116: nnz[mg] = mg+1;
117: MatSeqAIJSetPreallocation(*M,0,nnz);
118: PetscFree(nnz);
119: }
120: }
121: PEPEvaluateFunctionDerivatives(pep,pep->eigr[idx],coeffs);
122: for (i=0;i<nmat && PetscAbsScalar(coeffs[i])==0.0;i++);
123: BVGetColumn(pep->V,idx,&v);
124: MatMult(A[i],v,w);
125: if (coeffs[i]!=1.0) {
126: VecScale(w,coeffs[i]);
127: }
128: for (i++;i<nmat;i++) {
129: MatMult(A[i],v,q);
130: VecAXPY(w,coeffs[i],q);
131: }
132: BVRestoreColumn(pep->V,idx,&v);
133:
134: /* Set values */
135: PetscMalloc1(m1-m0,&cols2);
136: for (i=0;i<m1-m0;i++) cols2[i]=m0+i;
137: VecGetArray(w,&array);
138: for (i=m0;i<m1;i++) {
139: MatGetRow(*T,i,&ncols,&cols,&vals);
140: MatSetValues(*M,1,&i,ncols,cols,vals,INSERT_VALUES);
141: MatRestoreRow(*T,i,&ncols,&cols,&vals);
142: MatSetValues(*M,1,&i,1,&mg,array+i-m0,INSERT_VALUES);
143: }
144: VecRestoreArray(w,&array);
145: BVGetColumn(pep->V,idx,&v);
146: VecConjugate(v);
147: VecGetArray(v,&array);
148: MatSetValues(*M,1,&mg,m1-m0,cols2,array,INSERT_VALUES);
149: MatSetValues(*M,1,&mg,1,&mg,&zero,INSERT_VALUES);
150: VecRestoreArray(v,&array);
151: VecConjugate(v);
152: BVRestoreColumn(pep->V,idx,&v);
153: MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);
154: MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);
155: PetscFree(cols2);
156: PetscFree(coeffs);
157: return(0);
158: }
162: PetscErrorCode PEPNewtonRefinementSimple(PEP pep,PetscInt *maxits,PetscReal *tol,PetscInt k)
163: {
165: PetscInt i,j,n,its;
166: PetscMPIInt rank,size;
167: KSP ksp;
168: Mat M=NULL,T=NULL;
169: MPI_Comm comm;
170: Vec r,v,dv,rr=NULL,dvv=NULL,t[2];
171: PetscScalar *array,*array2,dh;
172: PetscReal norm,error;
173: PetscBool ini=PETSC_TRUE;
176: PetscLogEventBegin(PEP_Refine,pep,0,0,0);
177: its = (maxits)?*maxits:NREF_MAXIT;
178: comm = PetscObjectComm((PetscObject)pep);
179: KSPCreate(comm,&ksp);
180: BVGetColumn(pep->V,0,&v);
181: VecDuplicate(v,&r);
182: VecDuplicate(v,&dv);
183: VecDuplicate(v,&t[0]);
184: VecDuplicate(v,&t[1]);
185: BVRestoreColumn(pep->V,0,&v);
186: MPI_Comm_size(comm,&size);
187: MPI_Comm_rank(comm,&rank);
188: VecGetLocalSize(r,&n);
189: /* Loop performing iterative refinements */
190: for (j=0;j<k;j++) {
191: #if !defined(PETSC_USE_COMPLEX)
192: if (pep->eigi[j]!=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),1,"Simple Refinement not implemented in real scalar for complex eigenvalues");
193: #endif
194: for (i=0;i<its;i++) {
195: if (tol) {
196: BVGetColumn(pep->V,j,&v);
197: PEPComputeRelativeError_Private(pep,pep->eigr[j],0.0,v,NULL,&error);
198: BVRestoreColumn(pep->V,j,&v);
199: if (error<=*tol) break;
200: }
201: PEPNSimpleRefSetUp(pep,pep->A,j,&M,&T,ini,t);
202: KSPSetOperators(ksp,M,M);
203: if (ini) {
204: KSPSetFromOptions(ksp);
205: MatGetVecs(M,&dvv,NULL);
206: VecDuplicate(dvv,&rr);
207: ini = PETSC_FALSE;
208: }
209: BVGetColumn(pep->V,j,&v);
210: MatMult(T,v,r);
211: BVRestoreColumn(pep->V,j,&v);
212: VecGetArray(r,&array);
213: if (rank==size-1) {
214: VecGetArray(rr,&array2);
215: PetscMemcpy(array2,array,n*sizeof(PetscScalar));
216: array2[n] = 0.0;
217: VecRestoreArray(rr,&array2);
218: } else {
219: VecPlaceArray(rr,array);
220: }
221: KSPSolve(ksp,rr,dvv);
222: if (rank != size-1) {
223: VecResetArray(rr);
224: }
225: VecRestoreArray(r,&array);
226: VecGetArray(dvv,&array);
227: VecPlaceArray(dv,array);
228: BVGetColumn(pep->V,j,&v);
229: VecAXPY(v,-1.0,dv);
230: VecNorm(v,NORM_2,&norm);
231: VecScale(v,1.0/norm);
232: BVRestoreColumn(pep->V,j,&v);
233: VecResetArray(dv);
234: if (rank==size-1) dh = array[n];
235: VecRestoreArray(dvv,&array);
236: MPI_Bcast(&dh,1,MPIU_SCALAR,size-1,comm);
237: pep->eigr[j] -= dh;
238: }
239: }
240: KSPDestroy(&ksp);
241: MatDestroy(&M);
242: MatDestroy(&T);
243: VecDestroy(&t[0]);
244: VecDestroy(&t[1]);
245: VecDestroy(&dv);
246: VecDestroy(&dvv);
247: VecDestroy(&r);
248: VecDestroy(&rr);
249: PetscLogEventEnd(PEP_Refine,pep,0,0,0);
250: return(0);
251: }