[Trilinos-Users] problem computing eigenvalues with Sundance
Andrea Penza
pepe.aero at alice.it
Fri Feb 3 09:38:22 MST 2012
Dear Janos,
thanks to your suggestion I've been able to fix that problem. Anyway the simulation does not end fine.
The problem I'm trying to solve is a generalized non-hermitian eigenvalue problem, and I'm using these parameters:
<ParameterList>
<ParameterList name="Eigensolver">
<Parameter name="Method" type="string" value="Block Krylov Schur"/>
<Parameter name="Number of Eigenvalues" type="int" value="10"/>
<Parameter name="Block Size" type="int" value="4"/>
<Parameter name="Num Blocks" type="int" value="4"/>
<Parameter name="Verbosity" type="int" value="0"/>
<Parameter name="Which" type="string" value="SM"/>
<Parameter name="Use Preconditioner" type="bool" value="true"/>
<Parameter name="Is Hermitian" type="bool" value="false"/>
<Parameter name="Maximum Restarts" type="int" value="100"/>
<Parameter name="Use Locking" type="bool" value="false"/>
<Parameter name="Max Locked" type="int" value="1"/>
<Parameter name="Convergence Tolerance" type="double" value="1.0e-12"/>
<ParameterList name="Preconditioner">
<Parameter name="Type" type="string" value="ML"/>
<Parameter name="Problem Type" type="string" value="SA"/>
<ParameterList name="ML Settings">
<Parameter name="output" type="int" value="0"/>
</ParameterList>
</ParameterList>
</ParameterList>
</ParameterList>
I don't know if there is something wrong in this....
It seems to fail reaching convergence and the system reports this error:
p=0 | return code = 1
Sundance detected exception:
/usr/local/include/PlayaAnasaziEigensolverImpl.hpp:181:
Throw number = 1
Throw test that evaluated to true: returnCode != Anasazi::Converged
Anasazi did not converge!
test FAILED
Hope this sounds familiar to you...thanks to all for your help!
Best regards,
Andrea
Il giorno 02/feb/2012, alle ore 23.28, benk ha scritto:
> Dear Andrea,
>
> Try
> Expr massExpr = vx*ux + vy*uy + 0.0*q*p;
> instead of
> Expr massExpr = Integral( interior, vx*ux + vy*uy + 0.0*q*p, quad2);
> This should cure the problem.
> Sundance internally creates the Integral object.
>
> Best Regards,
> Janos
>
>
> On Thu, 2012-02-02 at 21:39 +0100, Andrea Penza wrote:
>> Dear all,
>>
>> I have been trying for days to solve the eigenvalue problem related to
>> Navier-Stokes equations in primitive variables using Sundance. Here is
>> the formulation that brings me troubles.
>>
>> // --------------------------------------------------------
>> bool lumpedMass = false;
>>
>> BasisFamily L1 = new Lagrange(1);
>> BasisFamily L2 = new Lagrange(2);
>>
>> Expr ux = new UnknownFunction(L2, "ux");
>> Expr vx = new TestFunction(L2, "vx");
>>
>> Expr uy = new UnknownFunction(L2, "uy");
>> Expr vy = new TestFunction(L2, "vy");
>>
>> Expr p = new UnknownFunction(L1, "p");
>> Expr q = new TestFunction(L1, "q");
>>
>> Expr u = List(ux, uy);
>> Expr v = List(vx, vy);
>>
>> // U is the solution of main flow, previously computed:
>> Expr U = List(u0[0], u0[1]);
>>
>> // Define the new weak form for the operator A:
>> Expr eqn3 = Integral(interior, (1/reynolds)*((grad*vx)*(grad*ux) +
>> (grad*vy)*(grad*uy)) + vx*(U*grad)*ux + vy*(U*grad)*uy +
>> vx*(u*grad)*u0[0] + vy*(u*grad)*u0[1] - p*(dx*vx+dy*vy), quad2 );
>>
>> Expr eqn4 = Integral( interior, (dx*ux + dy*uy) * q, quad2 );
>>
>> Expr eqnEig = eqn3 + eqn4;
>>
>> // Define the weak form to build the operator M:
>> Expr massExpr = Integral( interior, vx*ux + vy*uy + 0.0*q*p, quad2);
>>
>> LinearEigenproblem probEig(mesh, eqnEig, massExpr, List(vx, vy, q),
>> List(ux, uy, p), vecType, lumpedMass);
>> // -------------------------------------------------------------
>>
>> The compilation of the whole code ends fine.
>> Trying to run the executable the system reports this error:
>>
>>
>> Sundance detected exception:
>> /home/andrea/Projects/trilinos-10.8.3-Source/packages/Sundance/src-core/InternalExprs/SundanceSymbPreprocessor.cpp:221:
>>
>> Throw number = 1
>>
>> Throw test that evaluated to true: e==0
>>
>> Non-evaluatable expr Sum of Integrals[
>> Integral[
>> rqc=Integration Region
>> [6] cell filter=MaximalCellFilter
>> [6] quadrature rule=GaussianQuadrature[order=2]
>> [6] watchpoint=[]
>>
>> expr=vx*ux+vy*uy
>> ]
>> ]
>> given to SymbPreprocessor::setupExpr()
>> test FAILED
>>
>>
>> Has anyone the idea of what could be the problem?
>> Thanks in advance for your support.
>> Kind regards,
>>
>> Andrea
>>
>>
>>
>>
>>
>> _______________________________________________
>> Trilinos-Users mailing list
>> Trilinos-Users at software.sandia.gov
>> http://software.sandia.gov/mailman/listinfo/trilinos-users
>>
>
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: https://software.sandia.gov/pipermail/trilinos-users/attachments/20120203/3e0c5d47/attachment.html
More information about the Trilinos-Users
mailing list