[Dune] finding local coordinate in parent entity

Martin Nolte nolte at mathematik.uni-freiburg.de
Thu Apr 16 19:11:45 CEST 2015


Hi Aleksejs,

the answer is:

local3D = intersection.geometryInInside().global( local2D );

Best,

Martin

On 04/16/2015 06:48 PM, Aleksejs Fomins wrote:
> -----BEGIN PGP SIGNED MESSAGE-----
> Hash: SHA1
>
> 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
>
>
>
>
> -----BEGIN PGP SIGNATURE-----
> Version: GnuPG v1.4.12 (GNU/Linux)
>
> iQIcBAEBAgAGBQJVL+fDAAoJEDkNM7UOEwMZo/wP/1E03LXwcqoCCBEeT3RaKnCO
> ZSb59LUTNcGD90kFV1ggz7NUZDq+OxCqhEmvQf5GYw7pA+xPgjtFvoQ6N+stExZc
> UbMPNpHPJIClMF+r2rXDyFEVKsAtoEpNADU98sMrKpjzotZD2bYjCB/UFBZsVufV
> IAqtoB4Sr2T0O/+pnLoYszOR20yGeWlqUa6hZqNi61y03bOvlPaQCwwRwT4Rwavd
> 6TOuuRrqY0S1S8rjTdTmNFj1acFK9IOR3M8UlS/mtDnw9CFU6GVmHf5Xi6V0ZsVm
> 5P7m6JwAnOQSbscijjyEQJKnxN1yJWTL82P56ZXxrZkUB9OESBhFqnEwHB8/AEod
> nL5RS2vOVncS7DM+c/eCX2c6ftEYTWDXYKLpRpo1v2Ienk783pvbu2Lf2ILQvXcn
> hzFgzopOiFt5SNMlYVgp6pDJSI9x7HLONGO6jakDQmPeqfEXkKmPc/P0U2j3kDBB
> NU9QtwhOdTzifzzhCL+2kGJUjyKGPq/XnjGh2qfgCGRxwUa+mwaxHOkwcxWAqO1R
> RS2j1blAL1fj/5TqZVzW/Gg1Dm57LClGOA096bxv2b+41lZ0206nQZrSAowg7ijJ
> TR5d1crvjVg+O9S1MYV7yprpFxdilGnwgEl3CIYJnEC89oKk93wGWlhhhEBhtpsS
> X02GMwUMJMzXZOp52aK8
> =o4iI
> -----END PGP SIGNATURE-----
>
> _______________________________________________
> Dune mailing list
> Dune at dune-project.org
> http://lists.dune-project.org/mailman/listinfo/dune
>

-- 
Dr. Martin Nolte <nolte at mathematik.uni-freiburg.de>

Universität Freiburg                                   phone: +49-761-203-5630
Abteilung für angewandte Mathematik                    fax:   +49-761-203-5632
Hermann-Herder-Straße 10
79104 Freiburg, Germany




More information about the Dune mailing list