Skip to content

Commit 1252657

Browse files
cmdougprj-
authored andcommitted
Fix possible segfault bug in SLEPc-code.hpp
This commit fixes a possible segfault in the `SLEPc-code.hpp` that can occur if the `larray` argument in `EPSSolve` is specified but `-eps_two_sided` is not since `pti` was not being properly allocated in that case.
1 parent 0363697 commit 1252657

File tree

1 file changed

+9
-10
lines changed

1 file changed

+9
-10
lines changed

plugin/mpi/SLEPc-code.hpp

Lines changed: 9 additions & 10 deletions
Original file line numberDiff line numberDiff line change
@@ -436,7 +436,7 @@ AnyType eigensolver<Type, K, SType>::E_eigensolver::operator()(Stack stack) cons
436436
otherarray->resize(nr, nconv);
437437
VecGetLocalSize(xr, &n);
438438
} else xr = xi = NULL;
439-
if(std::is_same<SType, EPS>::value && (othervectors || otherarray) && isTwoSided)
439+
if(std::is_same<SType, EPS>::value && (othervectors || otherarray))
440440
MatCreateVecs(ptA->_petsc, &yi, &yr);
441441
else yr = yi = NULL;
442442
for(PetscInt i = 0; i < nconv; ++i) {
@@ -468,23 +468,23 @@ AnyType eigensolver<Type, K, SType>::E_eigensolver::operator()(Stack stack) cons
468468
PetscScalar* tmp2r;
469469
PetscScalar* tmp2i;
470470
VecGetArray(xr, &tmpr);
471-
if(std::is_same<SType, EPS>::value && (othervectors || otherarray) && isTwoSided)
471+
if(std::is_same<SType, EPS>::value && (othervectors || otherarray))
472472
VecGetArray(yr, &tmp2r);
473473
K* pt, *pti;
474474
if(!std::is_same<SType, SVD>::value) {
475475
if(std::is_same<PetscScalar, double>::value && std::is_same<K, std::complex<double>>::value) {
476476
VecGetArray(xi, &tmpi);
477477
pt = new K[n];
478478
copy(pt, n, tmpr, tmpi);
479-
if(std::is_same<SType, EPS>::value && (othervectors || otherarray) && isTwoSided) {
479+
if(std::is_same<SType, EPS>::value && (othervectors || otherarray)) {
480480
VecGetArray(yi, &tmp2i);
481-
pti = new K[n];
482-
copy(pti, n, tmp2r, tmp2i);
481+
pti = new K[nr];
482+
copy(pti, nr, tmp2r, tmp2i);
483483
}
484484
}
485485
else {
486486
pt = reinterpret_cast<K*>(tmpr);
487-
if(std::is_same<SType, EPS>::value && (othervectors || otherarray) && isTwoSided)
487+
if(std::is_same<SType, EPS>::value && (othervectors || otherarray))
488488
pti = reinterpret_cast<K*>(tmp2r);
489489
}
490490
}
@@ -540,15 +540,14 @@ AnyType eigensolver<Type, K, SType>::E_eigensolver::operator()(Stack stack) cons
540540
}
541541
if(!std::is_same<SType, SVD>::value && std::is_same<PetscScalar, double>::value && std::is_same<K, std::complex<double>>::value) {
542542
delete [] pt;
543-
if(std::is_same<SType, EPS>::value && (othervectors || otherarray) && isTwoSided)
544-
delete [] pti;
543+
delete [] pti;
545544
}
546545
if(std::is_same<SType, SVD>::value || (std::is_same<PetscScalar, double>::value && std::is_same<K, std::complex<double>>::value)) {
547546
VecRestoreArray(xi, &tmpi);
548-
if(std::is_same<SType, EPS>::value && (othervectors || otherarray) && isTwoSided)
547+
if(std::is_same<SType, EPS>::value && (othervectors || otherarray))
549548
VecRestoreArray(yi, &tmp2i);
550549
}
551-
if(std::is_same<SType, EPS>::value && (othervectors || otherarray) && isTwoSided)
550+
if(std::is_same<SType, EPS>::value && (othervectors || otherarray))
552551
VecRestoreArray(yr, &tmp2r);
553552
VecRestoreArray(xr, &tmpr);
554553
}

0 commit comments

Comments
 (0)