Dear Martin, Thanks for your rapid response. 1. The MappingQ is set to be the same as the order of velocity, as is shown in the following code snippet:
> *template <int dim>* > *NS<dim>::NS (ParameterHandler &prm)* > *:* > * parameters (&prm),* > * degree (prm.get_integer("pressure degree")),* > * fe( FE_Q<dim>(QGaussLobatto<1>(degree+2)), dim,* > * FE_Q<dim>(QGaussLobatto<1>(degree+1)), 1),* > * fe_scalar (FE_Q<dim>(QGaussLobatto<1>(degree+2))),* > * dof_handler (triangulation),* > * dof_handler_scalar (triangulation),* > * mapping (degree+2),* > * computing_timer (std::cout,* > * TimerOutput::summary,** TimerOutput::wall_times)* I am not quite sure whether the computation of the lift/drag in my code is right, and my implementation is almost the same as the one in this post: https://groups.google.com/d/msg/dealii/rS6soTb69ig/C4QchAyEGwAJ The only change is the first line : *QGauss<dim-1> face_quadrature_formula(degree+2);* 2. I am using deal.II 8.4.0 now, and I think i am not using manifold description in my code 3. The mesh is shown as follows: <https://lh3.googleusercontent.com/-Vs0s6So28js/WYrCbyok9mI/AAAAAAAABSc/rfs-1D1TA8c71IQZpOEyOoU3C1idNr1DQCLcBGAs/s1600/Image.png> Best, Howe 在 2017年8月9日星期三 UTC+8下午3:25:20,Martin Kronbichler写道: > > Dear Howe, > > How did you run your simulation? From your picture, it appears that a > higher order method is worse at higher degrees than a lower order method, > which does not match with my experience. If that were the case, nobody > would use high orders. However, you need to bring many pieces in place to > really get to the benefit of the high order method for somewhat more > complicated examples such as the flow around a cylinder. Here is a list of > things to look at: > > - Do you use a high-order polynomial mapping MappingQ of the same or > higher degree as the interpolation space? Do you use this mapping in all > routines that evaluate quantities, such as the usual assembly, the > computation of the lift/drag, and so on? > - Do you use a manifold description that extends into the domain? (Look > into TransfiniteInterpolationManifold.) Without, you will not get more than > third order convergence. > - Do you have a good mesh around the area of interest? Flows around > cylinders tend to be really really sensitive to the mesh quality around the > cylinder. > > For the Navier-Stokes equations around the cylinder, if everything is done > right one gets significantly improved results in terms of accuracy over the > number of degrees of freedom up to degree (6,5) (velocity,pressure). Beyond > that picture is less clear. At least with the meshes that we tried in our > group it was not worth to go beyond. You can have a look a our results in > section 5.4 and Figs. 9 and 10 of this preprint: > https://arxiv.org/pdf/1706.09252.pdf > > Best, > Martin > > On 09.08.2017 09:01, Howe wrote: > > Dear Jaekwang > > Have you solved this problem? If yes, Could pls share your solution with > us? > I am simulating a steady state flow over a cylinder, and the drag/lift > coefficient shows an unexpected trend of change as i increase the > discretization order and refine the mesh. > > > <https://lh3.googleusercontent.com/-Gz5932Zt6e0/WYqvbO3P-4I/AAAAAAAABSM/6hBHn1FO5W0P7X3SPfQ4iyRaldDYzOt3QCLcBGAs/s1600/Image.png> > > As is shown in the figure, the Cd increased as the cells increased for all > the discretization orders, however, for a fixed cells, the Cd decreased as > the discretization order increased. > > In my opinion, to increase the order and refine the mesh should both make > the approximation more close to the exact solution, thus should have the > same trend of change. > > 在 2016年9月18日星期日 UTC+8下午11:57:16,Jaekwang Kim写道: >> >> >> Hello, I am a starter of dealii and am learning a lot these days with the >> help of video lectures and tutorial examples. >> >> I modified step-22 code (stokes flow code) into my own problem, the flow >> around sphere. >> >> and I intend to evaluate the drag force (which is analytically given by >> stokes equation) >> >> My code reached quite close to the value since the absolute error : >> abs(drag_calculated-drag_exact)/drag_exact is around 10^(-3) >> >> However, I expected that if I input higher 'degree' I will receive more >> accurate result, but it didn't >> >> Obviously Q2 is better than Q1. and Q3 is better than Q2. But Q4 or Q4 is >> not better than Q2 or Q3? >> >> Is there any reason on this? >> >> (To be specific, if i say degree 2 , that mean I use (2+1) for velocity, >> (2) for pressure, and (2+2) for Gauss integral.... >> >> >> Thank you >> >> Jaekwang Kim >> > -- > The deal.II project is located at http://www.dealii.org/ > For mailing list/forum options, see > https://groups.google.com/d/forum/dealii?hl=en > --- > You received this message because you are subscribed to the Google Groups > "deal.II User Group" group. > To unsubscribe from this group and stop receiving emails from it, send an > email to dealii+un...@googlegroups.com <javascript:>. > For more options, visit https://groups.google.com/d/optout. > > > -- The deal.II project is located at http://www.dealii.org/ For mailing list/forum options, see https://groups.google.com/d/forum/dealii?hl=en --- You received this message because you are subscribed to the Google Groups "deal.II User Group" group. To unsubscribe from this group and stop receiving emails from it, send an email to dealii+unsubscr...@googlegroups.com. For more options, visit https://groups.google.com/d/optout.