Actual source code: neprefine.c
slepc-3.5.2 2014-10-10
1: /*
2: Newton refinement for NEP, 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/nepimpl.h>
25: #include <slepcblaslapack.h>
27: #define NREF_MAXIT 100
31: static PetscErrorCode NewtonSimpleRefSetUp(NEP nep,PetscInt nmat,Mat *A,PetscInt idx,Mat *M,Mat *T,PetscBool ini,Vec *t)
32: {
33: PetscErrorCode ierr;
34: PetscInt i,st,ml,m0,m1,mg;
35: PetscInt *dnz,*onz,ncols,*cols2,*nnz;
36: PetscScalar *array,zero=0.0,*coeffs;
37: PetscMPIInt rank,size;
38: MPI_Comm comm;
39: const PetscInt *cols;
40: const PetscScalar *vals;
41: Vec v,w=t[1],q=t[0];
44: comm = PetscObjectComm((PetscObject)A[0]);
45: PetscMalloc1(nmat,&coeffs);
46: MPI_Comm_rank(comm,&rank);
47: MPI_Comm_size(comm,&size);
48: if (ini) {
49: MatDuplicate(A[0],MAT_COPY_VALUES,T);
50: } else {
51: MatCopy(A[0],*T,SUBSET_NONZERO_PATTERN);
52: }
53: for (i=0;i<nmat;i++) {
54: FNEvaluateFunction(nep->f[i],nep->eigr[idx],coeffs+i);
55: }
56: if (coeffs[0]!=1.0) {
57: MatScale(*T,coeffs[0]);
58: }
59: for (i=1;i<nmat;i++) {
60: MatAXPY(*T,coeffs[i],A[i],(ini)?DIFFERENT_NONZERO_PATTERN:SUBSET_NONZERO_PATTERN);
61: }
62: MatGetSize(*T,&mg,NULL);
63: MatGetOwnershipRange(*T,&m0,&m1);
64: if (ini) {
65: MatCreate(comm,M);
66: MatGetLocalSize(*T,&ml,NULL);
67: if (rank==size-1) ml++;
68: MatSetSizes(*M,ml,ml,mg+1,mg+1);
69: MatSetFromOptions(*M);
70: MatSetUp(*M);
71: /* Preallocate M */
72: if (size>1) {
73: MatPreallocateInitialize(comm,ml,ml,dnz,onz);
74: for (i=m0;i<m1;i++) {
75: MatGetRow(*T,i,&ncols,&cols,NULL);
76: MatPreallocateSet(i,ncols,cols,dnz,onz);
77: MatPreallocateSet(i,1,&mg,dnz,onz);
78: MatRestoreRow(*T,i,&ncols,&cols,NULL);
79: }
80: if (rank==size-1) {
81: PetscCalloc1(mg+1,&cols2);
82: for (i=0;i<mg+1;i++) cols2[i]=i;
83: MatPreallocateSet(m1,mg+1,cols2,dnz,onz);
84: PetscFree(cols2);
85: }
86: MatMPIAIJSetPreallocation(*M,0,dnz,0,onz);
87: MatPreallocateFinalize(dnz,onz);
88: } else {
89: PetscCalloc1(mg+1,&nnz);
90: for (i=0;i<mg;i++) {
91: MatGetRow(*T,i,&ncols,NULL,NULL);
92: nnz[i] = ncols+1;
93: MatRestoreRow(*T,i,&ncols,NULL,NULL);
94: }
95: nnz[mg] = mg+1;
96: MatSeqAIJSetPreallocation(*M,0,nnz);
97: PetscFree(nnz);
98: }
99: }
100: for (i=0;i<nmat;i++) {
101: FNEvaluateDerivative(nep->f[i],nep->eigr[idx],coeffs+i);
102: }
103: st = 0;
104: for (i=0;i<nmat && PetscAbsScalar(coeffs[i])==0.0;i++) st++;
105: BVGetColumn(nep->V,idx,&v);
106: MatMult(A[st],v,w);
107: if (coeffs[st]!=1.0) {
108: VecScale(w,coeffs[st]);
109: }
110: for (i=st+1;i<nmat;i++) {
111: MatMult(A[i],v,q);
112: VecAXPY(w,coeffs[i],q);
113: }
114: BVRestoreColumn(nep->V,idx,&v);
115: /* Set values */
116: PetscMalloc1(m1-m0,&cols2);
117: for (i=0;i<m1-m0;i++) cols2[i]=m0+i;
118: VecGetArray(w,&array);
119: for (i=m0;i<m1;i++) {
120: MatGetRow(*T,i,&ncols,&cols,&vals);
121: MatSetValues(*M,1,&i,ncols,cols,vals,INSERT_VALUES);
122: MatRestoreRow(*T,i,&ncols,&cols,&vals);
123: MatSetValues(*M,1,&i,1,&mg,array+i-m0,INSERT_VALUES);
124: }
125: VecRestoreArray(w,&array);
126: BVGetColumn(nep->V,idx,&v);
127: VecConjugate(v);
128: VecGetArray(v,&array);
129: MatSetValues(*M,1,&mg,m1-m0,cols2,array,INSERT_VALUES);
130: MatSetValues(*M,1,&mg,1,&mg,&zero,INSERT_VALUES);
131: VecRestoreArray(v,&array);
132: VecConjugate(v);
133: BVRestoreColumn(nep->V,idx,&v);
134: MatAssemblyBegin(*M,MAT_FINAL_ASSEMBLY);
135: MatAssemblyEnd(*M,MAT_FINAL_ASSEMBLY);
136: PetscFree(cols2);
137: PetscFree(coeffs);
138: return(0);
139: }
143: PetscErrorCode NEPNewtonRefinementSimple(NEP nep,PetscInt *maxits,PetscReal *tol,PetscInt k)
144: {
146: PetscInt i,j,n,its;
147: PetscMPIInt rank,size;
148: KSP ksp;
149: Mat M,T;
150: MPI_Comm comm;
151: Vec r,v,dv,rr,dvv,t[2];
152: PetscScalar *array,*array2,dh;
153: PetscReal norm,error;
154: PetscBool ini=PETSC_TRUE;
157: PetscLogEventBegin(NEP_Refine,nep,0,0,0);
158: its = (maxits)?*maxits:NREF_MAXIT;
159: comm = PetscObjectComm((PetscObject)nep);
160: KSPCreate(comm,&ksp);
161: BVGetColumn(nep->V,0,&v);
162: VecDuplicate(v,&r);
163: VecDuplicate(v,&dv);
164: VecDuplicate(v,&t[0]);
165: VecDuplicate(v,&t[1]);
166: BVRestoreColumn(nep->V,0,&v);
167: MPI_Comm_size(comm,&size);
168: MPI_Comm_rank(comm,&rank);
169: VecGetLocalSize(r,&n);
170: /* Loop performing iterative refinements */
171: for (j=0;j<k;j++) {
172: #if !defined(PETSC_USE_COMPLEX)
173: if (nep->eigi[j]!=0.0) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Simple Refinement not implemented in real scalar for complex eigenvalues");
174: #endif
175: for (i=0;i<its;i++) {
176: if (tol) {
177: BVGetColumn(nep->V,j,&v);
178: NEPComputeRelativeError_Private(nep,nep->eigr[j],v,&error);
179: BVRestoreColumn(nep->V,j,&v);
180: if (error<=*tol) break;
181: }
182: NewtonSimpleRefSetUp(nep,nep->nt,nep->A,j,&M,&T,ini,t);
183: KSPSetOperators(ksp,M,M);
184: if (ini) {
185: KSPSetFromOptions(ksp);
186: MatGetVecs(M,&dvv,NULL);
187: VecDuplicate(dvv,&rr);
188: ini = PETSC_FALSE;
189: }
190: BVGetColumn(nep->V,j,&v);
191: MatMult(T,v,r);
192: BVRestoreColumn(nep->V,j,&v);
193: VecGetArray(r,&array);
194: if (rank==size-1) {
195: VecGetArray(rr,&array2);
196: PetscMemcpy(array2,array,n*sizeof(PetscScalar));
197: array2[n] = 0.0;
198: VecRestoreArray(rr,&array2);
199: } else {
200: VecPlaceArray(rr,array);
201: }
202: KSPSolve(ksp,rr,dvv);
203: if (rank != size-1) {
204: VecResetArray(rr);
205: }
206: VecRestoreArray(r,&array);
207: VecGetArray(dvv,&array);
208: VecPlaceArray(dv,array);
209: BVGetColumn(nep->V,j,&v);
210: VecAXPY(v,-1.0,dv);
211: VecNorm(v,NORM_2,&norm);
212: VecScale(v,1.0/norm);
213: BVRestoreColumn(nep->V,j,&v);
214: VecResetArray(dv);
215: if (rank==size-1) dh = array[n];
216: VecRestoreArray(dvv,&array);
217: MPI_Bcast(&dh,1,MPIU_SCALAR,size-1,comm);
218: nep->eigr[j] -= dh;
219: }
220: }
221: KSPDestroy(&ksp);
222: MatDestroy(&M);
223: MatDestroy(&T);
224: VecDestroy(&t[0]);
225: VecDestroy(&t[1]);
226: VecDestroy(&dv);
227: VecDestroy(&dvv);
228: VecDestroy(&r);
229: VecDestroy(&rr);
230: PetscLogEventEnd(NEP_Refine,nep,0,0,0);
231: return(0);
232: }