[Trilinos-Users] problem computing eigenvalues with Sundance

benk benk at in.tum.de
Thu Feb 2 15:28:36 MST 2012


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
> 





More information about the Trilinos-Users mailing list