[Nek5000-users] multiple runs and poisson for velocity

nek5000-users at lists.mcs.anl.gov nek5000-users at lists.mcs.anl.gov
Mon Apr 7 01:54:28 CDT 2014


Thank you Paul, I managed to obtain the Poisson solver I needed.
Best,
Giuseppe



Il giorno 05/apr/2014, alle ore 20:34, <nek5000-users at lists.mcs.anl.gov>
 <nek5000-users at lists.mcs.anl.gov> ha scritto:

> 
> Hi Giuseppe,
> 
> Your understanding is essentially correct.  One minor point is that the rhs should be
> 
> 
>    rhs = B*f
> 
> i.e., as given by:
> 
>       call col3(rhs,bm1,f,n)
> 
> or (more typically)
> 
>      call col2(rhs,bm1,n)
> 
> where rhs() would be filled with the desired rhs.   The idea here is that the rhs should
> be multiplied by the mass matrix --- that is, when nek solves the Poisson equation it
> is expressed as:
> 
> 
>        A u = B f
> 
> for  -del^2 u = f(x,y), say.
> 
> I hope this helps.
> 
> Best, Paul
> 
> ________________________________________
> From: nek5000-users-bounces at lists.mcs.anl.gov [nek5000-users-bounces at lists.mcs.anl.gov] on behalf of nek5000-users at lists.mcs.anl.gov [nek5000-users at lists.mcs.anl.gov]
> Sent: Saturday, April 05, 2014 8:36 AM
> To: <nek5000-users at lists.mcs.anl.gov>
> Subject: Re: [Nek5000-users] multiple runs and poisson for velocity
> 
> Hi Paul,
> thank you for referring me to the psi-omega example. This was indeed useful, if I have understood correctly the way it works.
> If I am right, to solve a Poisson problem for a variable psi for instance, the steps are the following:
> 
> 1) write the desired values to the rhs vector
> 2) call rone(h1,n)  ! this sets to one the coefficient multiplying the stiffness matrix
> 3) call rzero(h2,n) ! this sets to zero the coefficient multiplying the mass matrix
> 4) call hmholtz('psi ',psi,rhs,h1,h2,tmask,tmult,imsh,tol,maxiter,isd) ! call the poisson solver for the psi vector, imsh specifies which mesh is used, tol is the tolerance for the CG, maxiter the maximum number of iterations, isd should be set to two for axisymmetric problems. If i am right tmult is the vector with the quadrature weights for the temperature mesh, should I replace tmult with vmult is psi is defined on the velocity mesh? tmask if I am right is related to the Dirichlet b.c. on the temperature mesh, should I replace this with v1mask if I work on velocity components?
> 5) then I should be done, the example continues computing the gradient of psi to get the velocity components from the velocity potential, but I don't need this.
> 
> Thanks also to Aleks for his help with the first part of my question, I will write again if I have other questions.
> Best regards,
> Giuseppe
> 
> 
> 
> Il giorno 04/apr/2014, alle ore 20:01, <nek5000-users at lists.mcs.anl.gov>
> <nek5000-users at lists.mcs.anl.gov> ha scritto:
> 
>> 
>> Hi Giuseppe,
>> 
>> If you look at the .usr file in examples/eddy_psi_omega/
>> you'll see a call for a Poisson solve for a little psi-omega formulation
>> of Navier-Stokes.
>> 
>> Please let me know if this provides the insight you seek for setting
>> up your Poisson solves.
>> 
>> Best regards,
>> 
>> Paul
>> 
>> ________________________________________
>> From: nek5000-users-bounces at lists.mcs.anl.gov [nek5000-users-bounces at lists.mcs.anl.gov] on behalf of nek5000-users at lists.mcs.anl.gov [nek5000-users at lists.mcs.anl.gov]
>> Sent: Friday, April 04, 2014 7:09 AM
>> To: nek5000-users at lists.mcs.anl.gov
>> Subject: [Nek5000-users] multiple runs and poisson for velocity
>> 
>> Dear Neks,
>> I have two issues that I would like to discuss with you.
>> Firstly, I would like to perform many steady state simulations with different values of the viscosity, each time using the previous solution as starting condition.
>> For this reason, I have built the following subroutine:
>> 
>> c-----------------------------------------------------------------------
>>     subroutine userchk
>>     include 'SIZE'
>>     include 'TOTAL'
>> 
>> 
>>     do 120 ire=1,10
>>     param(2) = -1*ire
>>     call setprop
>>     igeom=1
>>     call fluid(igeom)
>> 120 continue
>> 
>>     ifxyo = .true.
>>     call outpost(vx,vy,vz,pr,t,'ste')
>>     call exitt
>> 
>>     return
>>     end
>> c-----------------------------------------------------------------------
>> 
>> but this does not work as I expected, could somebody suggest me what I should do instead?
>> 
>> Secondly, I would like to solve (in postprocessing mode) a Poisson problem for the velocity components, in the form -lapl(vx) = fx and -lapl(vy)=fy, with fx and fy known.
>> Is it possible to do this within NEK? I think this could be done by adapting the Uzawa pressure solver in navier1.f, but I am not sure if this could work.
>> Best regards,
>> Giuseppe
>> _______________________________________________
>> Nek5000-users mailing list
>> Nek5000-users at lists.mcs.anl.gov
>> https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users
>> _______________________________________________
>> Nek5000-users mailing list
>> Nek5000-users at lists.mcs.anl.gov
>> https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users
> 
> _______________________________________________
> Nek5000-users mailing list
> Nek5000-users at lists.mcs.anl.gov
> https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users
> _______________________________________________
> Nek5000-users mailing list
> Nek5000-users at lists.mcs.anl.gov
> https://lists.mcs.anl.gov/mailman/listinfo/nek5000-users



More information about the Nek5000-users mailing list