[Dune] Handling checkinside for curvilinear elements
Andreas Dedner
a.s.dedner at warwick.ac.uk
Mon Jul 28 18:10:01 CEST 2014
PS: for local I would have though you can use a Newton method.to find a
root of
F(y) = global(y)-x
for some given global point x. That should be uniquely solvable
otherwise your reference mapping is not
one-on-one and that would be "!unfortunate".
Andreas
On 28/07/14 17:03, Andreas Dedner wrote:
> Hi.
> I think all your problems boil down to the last one: implementing the
> method local.
> If you have local then a point is a global point x is within an
> element if local(x) is within the
> reference element. That is what we use for the hoerarchicSearch class
> (a utility class in
> dune grid which should then be directly useable for your second point).
> I don't think I understand your first problem since you are still
> using the same reference element, right?
> So a condition for a point being inside the reference element has
> nothing to do with your global mapping,
> be it straight faced or not.
> Best
> Andreas
>
> On 28/07/14 16:55, Aleksejs Fomins wrote:
>> Dear Dune,
>>
>> I am currently implementing curvilinear tetrahedral grid for dune.
>>
>> The following functionality will be required:
>>
>> 1) dune/geometry/genericreferenceelements.hh has function isinside()
>> which checks whether the point is inside of the tetrahedron. This
>> functionality will not be valid for curvilinear tetrahedrons.
>>
>> Where should I implement the new functionality? Should I directly edit
>> the genericreferenceelements.hh or should I define a new method inside
>> my new version of multilineargeometry.hh, or sth else?
>>
>> Where is this functionality used?
>>
>> Also, if you have suggestions on how to mathematically check this,
>> please advise me. I have a conjecture that the point is inside if it is
>> on the correct side of each face, and it is on the correct side of the
>> face if det(d, v1, v2) has the right sign, where d is the distance
>> vector between this point and the closest point on the surface, and v1
>> and v2 are the surface tangential unit vectors given by the parametric
>> coordinates.
>>
>> 2) We personally would really need the functionality of finding the
>> tetrahedron to which a given global coordinate belongs. Extending this
>> to curvilinear is trivial - first find the element to which it would
>> have belonged to if the tetrahedrons were straight-sided, then check
>> whether it is actually inside this tetrahedron or one of its neighbors
>> using the isinside() method.
>>
>> Is this (straight-sided) functionality already implemented in dune
>> somewhere? Is it efficient? I mean, whether it loops over all elements
>> and checks isinside(), or whether it uses something more sophisticated
>> like octrees.
>>
>> 3) Another mathematics-related question. To map from local to global
>> coordinates within an element I simply evaluate the linear combination
>> of lagrange polynomials.
>>
>> How do I map from global to local?
>>
>> The best thing that comes to my mind is to write an iterative
>> optimization algorithm minimizing the distance between given local
>> coordinate and the expected one. It will have to be random like MCMC
>> because the problem will have local minima. Please comment if you have
>> good ideas :)
>>
>> Thank you,
>> Aleksejs
>>
>>
>>
>> _______________________________________________
>> Dune mailing list
>> Dune at dune-project.org
>> http://lists.dune-project.org/mailman/listinfo/dune
>
>
>
> _______________________________________________
> Dune mailing list
> Dune at dune-project.org
> http://lists.dune-project.org/mailman/listinfo/dune
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <https://lists.dune-project.org/pipermail/dune/attachments/20140728/93934824/attachment.htm>
More information about the Dune
mailing list