From: Steffen P. <ste...@tu...> - 2006-01-17 14:15:57
|
> Hi, I also use the generalized eigenvalue problem and it=20 > works fine. Many > thanks. I have a question: I have a real symmetric problem,=20 > so I get only real > eigenvalues, but I need to get the lowest ones. So I use: >=20 > eigen_system.eigen_solver->set_position_of_spectrum(SMALLEST_M AGNITUDE); >=20 > The problem is, that all the solvers from slepc package=20 > (ARNOLDI, SUBSPACE, > POWER) don't support the "SMALLEST_MAGNITUDE" option, see=20 > their sources, for > example slepc-2.3.0/src/eps/impls/arnoldi/arnoldi.c, line 38: >=20 > if (eps->which!=3DEPS_LARGEST_MAGNITUDE) > SETERRQ(1,"Wrong value of eps->which"); >=20 > Except the "LAPACK" solver. Lapack works fine, but it=20 > retrieves all the > eigenvalues and is really slow, in fact for a bigger mesh=20 > it's completely > unusable (but works). Is there a way how to get the lowest=20 > eigenvalues with > slepc (maybe using some trick)? The actual slepc version (2.3.0) also offers a Lanczos solver (for Hermitian problems) that can be used to compute the smallest values in the spectrum (since I have been using an older slepc version so far I have not tested this solver yet). I have just commited some changes to also include the Lanczos algorithm that comes with the current version. For computing the smallest eigenvalues I have used slepc (2.2.1) together with the arpack package, but it seems that I need to change some cunfiguration stuff to get everything working with the 2.3.0 version. Steffen =20 > In the meantime, I export the element matrices to a file and=20 > use external > solver, pysparse >=20 > http://people.web.psi.ch/geus/pyfemax/pysparse.html >=20 > so I have a small (88 lines) python script, which assembles=20 > the global matrices > (they are called A and B in libmesh), and calls Jacobi=20 > Davidson method from > pysparse module. The solver is written in C (exported to=20 > python), it's fast and > it works.=20 >=20 > Ondrej >=20 > On Mon, Jan 16, 2006 at 12:25:50AM +0100, Steffen Petersen wrote: > > Hello Michael, > >=20 > > a few weeks ago you have asked for > > solving generalized eigenvalue problems with libMesh. > > Meanwhile, I have extended the slepc interface > > to also support generalized problems > > (with a few modifaications example 16 > > may be used for solving a generalized eigenvalue > > problem). > > Feel free to ask if you are interested in some > > more details. > >=20 > > Steffen > >=20 > >=20 > >=20 > > Michael Povolotskyi schrieb: > > >Dear Libmesh developers, > > >I have to solve an eigenvalue problem for a quantum=20 > physics calculations. > > >My questions are: > > >1) Is it possible to solve a general eigenvalue problem=20 > like Hx =3D Sx,=20 > > >where H and S are known matrixes. > > >2) Is it possible to have complex (self-conjugated) matrixes? > > > > > >Thank you, > > >Michael. > > >=20 > > >-------------------------------------------------------=20 > This SF.net=20 > > >email is sponsored by: Splunk Inc. Do you grep through log=20 > files for=20 > > >problems? Stop! Download the new AJAX search engine that=20 > makes searching=20 > > >your log files as easy as surfing the web. DOWNLOAD SPLUNK!=20 > > >http://ads.osdn.com/?ad_idv37&alloc_id=16865&op=3Dclick=20 > > >_______________________________________________=20 > Libmesh-users mailing=20 > > >list Lib...@li...=20 > > >https://lists.sourceforge.net/lists/listinfo/libmesh-users > >=20 > >=20 > >=20 > > ------------------------------------------------------- > > This SF.net email is sponsored by: Splunk Inc. Do you grep=20 > through log files > > for problems? Stop! Download the new AJAX search engine that makes > > searching your log files as easy as surfing the web. =20 > DOWNLOAD SPLUNK! > > http://ads.osdn.com/?ad_id=3D7637&alloc_id=3D16865&op=3Dclick > > _______________________________________________ > > Libmesh-users mailing list > > Lib...@li... > > https://lists.sourceforge.net/lists/listinfo/libmesh-users >=20 |