Hi,

It seems that some other problems (more serious) are present in
"KellyErrorEstimator" class (of course if I understood correctly), I just
mention some of them:


please refer to "error_estimator.cc", I focus on a pice of function
"integrate_over_regular_face":

first, it is obvious that each item of array:
"per_thread_data.psi[n][p][component]" is a vector (include gradient)

in line 1274: dot product of gradient at each qp by face-normal at each qp
is computed (i.e. du/dn), result for each qp is a scaler and this scaler is
cast to vector value gradient.

If we ignore performance leake for such excess operation. If implementation
of operator overloading in deallii vector algebra goes correct, We expect
that each component of gradient at qp include result of above mentioned dot
product (we should have duplicated storage).

But it's not happen in practice, if you print value of
"per_thread_data.psi[n][point][component]" prior and after dot product it is
not changed at all cases. e.g. effect of above dot product was not felt. I
do not know what is reason, maybe related to vector algebra, etc. (I use MT
version, i do not know can it cause or no).

to study this simply replace the following pice instead of line: 1271-1276
(of course add iostream, fstream headers to support output on stream) and
recompile library, i have used step-6 (though not important:i fixed its
coeff matrix to 1, zero dirichlet bc and rhs =4, by this exact solution is
known):

--------------------------

  std::cout<<"\n  prior to dot  = "<<per_thread_data.psi[0][0][0]<<"\n ";

  for (unsigned int n=0; n<n_solution_vectors; ++n)
    for (unsigned int component=0; component<n_components; ++component)
      for (unsigned int point=0; point<n_q_points; ++point)
 per_thread_data.phi[n][point][component]
   = (per_thread_data.psi[n][point][component] *
      per_thread_data.normal_vectors[point]);

std::cout<<"\n  after dot product  = "<<per_thread_data.psi[0][0][0]<<"\n ";


------------------------------

it seems that there are other problem besides to this, but lets to first
hear from experts.


Cheers

RT

On Thu, Nov 27, 2008 at 8:16 PM, Ruhollah Tavakoli <[EMAIL PROTECTED]> wrote:

> Hi,
>
>
> After being more than a year away from dealii, today start to use dealii
> egain, so i'm not very familiar with code and maybe my bug report be silly
> (apology in advance):
>
> my version: stable 6.1.0
>
> It seems that there is a logical bug in "KellyErrorEstimator" class, it
> takes a vector field (each field can be multi-component itself), so it
> should takes a vector boundary conditions (each one for each field), but it
> just take one neumann bc object and so use same bc for all fields.
>
> to see more about what happen refer to "error_estimator.cc":  line 1331 or
> line 1344, you can see that the same boundary condition is taken into
> account for all "n_solution_vectors"
>
> am i miss something?
>
>
> Cheers
>
> RT
>
>
> PS1: in this estimator for dirichlet faces, no action is considered, my
> question: does this mean that this estimator may not work as well for
> Dirichlet problems (in constrast to neumann ones)?
>
> PS2: notice that example folder is missed from "deal.nodoc-6.1.0.tar" and
> "deal.doc-6.1.0.tar", it's infeasible, so i suggest to fix it. Can someone
> just post zipped example folder of 6.1.0 version, i do not like to download
> a 60MB again!
>
>
_______________________________________________

Reply via email to