[Dune] finding local coordinate in parent entity

Aleksejs Fomins aleksejs.fomins at lspr.ch
Thu Apr 16 19:51:47 CEST 2015


-----BEGIN PGP SIGNED MESSAGE-----
Hash: SHA1

Damn, now I feel really ashamed, this is so natural and obvious

Thanks a lot Martin :)

On 16/04/15 19:11, Martin Nolte wrote:
> Hi Aleksejs,
> 
> the answer is:
> 
> local3D = intersection.geometryInInside().global( local2D );
> 
> Best,
> 
> Martin
> 
> On 04/16/2015 06:48 PM, Aleksejs Fomins wrote:
> Dear Dune,
> 
> I am now working on our electromagnetic code to make it compatible with the curvilinear grid. When it finally runs we will have good evidence to believe that the grid manager is operational.
> 
> 
> When reading through the code, I realise that we quite frequently execute the following sequence
> 
> -------------------------------------------------------------------------------------------
> 
> 1) Iterate over elements
> 2) Iterate over intersections of an element
> 3) Integrate a certain function over the intersection using Quadrature
> 3.1) The integral is done over the intersection, so quadrature gives 2D local coordinate
> 3.2) Function is defined over the element, so it requires 3D local coordinate. Currently we do it this way:
> 
> FieldVector<ctype, dim-1> local2D = quadratureIterator->position();
> FieldVector<ctype, dim>   global3D = intersection.geometry().global(local2D);
> FieldVector<ctype, dim>   local3D = elementInside.geometry().local(global3D);
> 
> -------------------------------------------------------------------------------------------
> 
> I know that for linear geometries this is probably not so relevant, but for curvilinear geometries calling the local method is very expensive (Newton method), especially if the point is on a face of an element.
> 
> If I am missing it, could anyone suggest how the above sequence can be run cheaper using the current grid functionality.
> 
> 
> If this is indeed the fastest way this can be currently done, then I propose to extend the functionality of the intersection class with the following 2 functions:
> 
> FieldVector<ctype, dim>  localInside(FieldVector<ctype, dim-1> local)
> FieldVector<ctype, dim>  localOutside(FieldVector<ctype, dim-1> local)
> 
> which would find the local coordinate inside of one of the neighbouring elements based on the local coordinate of the intersection. For the simplex geometries the implementation of these methods would boil down to inserting the face coordinate into the element coordinate and permuting based on the face subentity number. In case of non-conformal refinement one would also then have to shift and scale the coordinate to the correct neighbour, but nothing impossible.
> 
> Also, linear grid implementers could simply do it using the above local->global->local transformation, because it is not too expensive in the linear case.
> 
> 
> Please tell me what you think,
> Regards,
> Aleksejs
> 
> 
> 
> 
>>
>> _______________________________________________
>> Dune mailing list
>> Dune at dune-project.org
>> http://lists.dune-project.org/mailman/listinfo/dune
>>
> 
-----BEGIN PGP SIGNATURE-----
Version: GnuPG v1.4.12 (GNU/Linux)

iQIcBAEBAgAGBQJVL/azAAoJEDkNM7UOEwMZ/lEQAInQ1NJYaTMz2I394VJ22ZTU
MhwSjKDuVCnvXzeme2xR73JkVVeEkPRr/LFDtJLBjexhMX+W0mknfKV1MLR7s+8X
PBiHIwb74/scXbHU2wseYJl96FPuEPeAGMUsbL9e63lyvLifvcQCot2IaYPqrZk3
ZfBJKTGSpksx0XZ7RCpDIc1+OWWz+td43O7AaueeGduE2CDdp279Le7PkXUrFrGu
gBPWD6Xgaks67fYKNaqJQSKMTuxnbd4RAiCobItvdjA65ERqAejfInujtEkGw2qU
lzIFyt5HmfHO1v14c1RaXFsXBZSOdInAzMJGsUzmdcX6jvr/pvh+TOgp05gxtzNZ
MPmHQjOgCGx5i3PGOda2z816uKtvx4z/RTJFCYk4kW74vBodeswogh+1MAtYJz3Y
epjXRW0HDvH2u7h69YIWchfJTNiH/1PnNQYD+5DrS7f12WzWDGxFdXXdYjLMGhg8
KsppHzH11VPL8Fx9wEDzi5Eo/sk1/ATMZLPRmaPt066/FBUtrbXHRd/1KdihLnc3
bSb94rZ2Qx4pwid4xBXxBkuqlwRnNe/ARZpw5Oz+m6GLYOFrY5KkcpHmnT3IdGZ4
wVVzVnavLMoyPsddsltIyodAvTEssF2fh4vTy7ZzX89yr9itHG1InZve/BBjmJfD
C2jJgaj1z+XFyA/uexQS
=Kd2m
-----END PGP SIGNATURE-----




More information about the Dune mailing list