Thank you a lot! I`ve tried to play around and my conclusion is trivial and 
exactly what you have said - we do really need to solve a full 3d equation 
for 2D geometry too, it comes from physics, not from coding or math. Just 
imagine a typical 1D (or 2D) problem - an infinite subwavelength plate made 
from uniaxial crystal, with a wave passing it at normal incidence. Such a 
1D setup (geometrically) we need all field components to describe a 
phenomenon, an extraordinary ray for sure will move away form initial 1D 
line and in general will move  out of 2D plane of consideration.

I am still a bit confused. E.g. to simulate  Fresnel equation in 1D should 
I need to consider a 1D stack of hexahedral? I do actually need only 2 
components for field out of six (one normal for electric field, and one 
normal for magnetic for 1D case). However, to keep compatibility with 3D 
case I need to evaluate all six components, so it will be easy to move to a 
2D as a plane of hexahedra or to a 3d case with a volume.  Probably I 
should simply put a constraint for components that I am sure are not needed 
to be equal zero. This way the number of DoF per element should be 
reduced...




понедельник, 26 декабря 2016 г., 23:44:49 UTC+3 пользователь Timo Heister 
написал:
>
> The problem is that .curl returns a scalar in 2d and a tensor in 3d. A 
> typical definition of maxwell curl notation is to treat the 2d case as 
> a 3d solution that is constant in the (implicit) 3rd dimension, right? 
> Maybe you can do something like 
>
> Tensor<1,3> make_3d_curl(double curl_2d) 
> { 
> Tensor<1,3> r; r[2] = curl_2d; return r; 
> } 
>
> Tensor<1,3> make_3d_curl(Tensor<1,3> curl_3d) 
> { 
> return curl_3d; 
> } 
>
> and then do all the math in 3d: 
>
> Tensor<1,dim> curl_W_j = make_3d_curl(fe_values[E_field].curl  (j, q)); 
>
>   local_rhs_KT (i, j) += curl_W_j * ... 
>
> etc. 
>
>
>
>
>
> On Mon, Dec 26, 2016 at 1:00 PM, Konstantin Ladutenko 
> <kosty...@gmail.com <javascript:>> wrote: 
> > Dear all, 
> > 
> > Any ideas, what is the best way to go when implementing Maxwell 
> equations to 
> > make them dimension independent? Or at least to have the same code for 
> 2d 
> > and 3d cases? 
> > 
> > Simply using dim template parameter doesn`t work by definition - when 
> > switching from 3D to 2D  you need to specify, should it be TE or TM mode 
> > (which field, electric of magnetic, should be normal to the 2D plane). 
> > Introducing some boolean variable  `isTEmode` to be True or False can 
> help, 
> > but this will lead to many if\then codes... As soon as it is still the 
> same 
> > Maxwell equation I hope there is a better way to do this. Any ideas? 
> > 
> > From the code point of view it does not work either at the moment. Using 
> > 
> >      fe (FE_RaviartThomas<dim> (degree), 1, FE_Nedelec<dim> (degree), 
> 1), 
> > 
> > where Raviart-Thomas face element is needed for magnetic field, Nedelec 
> edge 
> > element is needed for electric field and 
> > 
> >     const FEValuesExtractors::Vector B_field (0); 
> >     const FEValuesExtractors::Vector E_field (dim); 
> >     const auto &curl_W_j = fe_values[E_field].curl  (j, q); 
> >     const auto &F_i      = fe_values[B_field].value (i, q); 
> > 
> > works fine for 3D (some wave is traveling in free space) and breaks down 
> in 
> > 2D during compilation of a weak-form integral local component. 
> > 
> > 
> /home/tig/symlink/dealii/examples/step-maxwell-nonliner/step-maxwell-linear.cc:
>  
>
> > In instantiation of ‘void Maxwell::MaxwellTD<dim>::assemble_system() 
> [with 
> > int dim = 2]’: 
> > 
> /home/tig/symlink/dealii/examples/step-maxwell-nonliner/step-maxwell-linear.cc:549:20:
>  
>
> > required from ‘void Maxwell::MaxwellTD<dim>::run() [with int dim = 2]’ 
> > 
> /home/tig/symlink/dealii/examples/step-maxwell-nonliner/step-maxwell-linear.cc:648:32:
>  
>
> > required from here 
> > 
> /home/tig/symlink/dealii/examples/step-maxwell-nonliner/step-maxwell-linear.cc:394:51:
>  
>
> > error: no match for ‘operator*’ (operand types are ‘const 
> dealii::Tensor<1, 
> > 1, double>’ and ‘const dealii::Tensor<1, 2, double>’) 
> >                    local_rhs_KT (i, j) += curl_W_j * F_i * fe_values.JxW 
> (q) 
> > 
> > Full example is available online at GitHub 
> > 
> > 
> > 
> https://urldefense.proofpoint.com/v2/url?u=https-3A__github.com_kostyfisik_dealii_blob_0188064ee4dd0cdd53d515f33932732b0a2d5d66_examples_step-2Dmaxwell-2Dnonliner_step-2Dmaxwell-2Dlinear.cc&d=CwIFaQ&c=Ngd-ta5yRYsqeUsEDgxhcqsYYY1Xs5ogLxWPA_2Wlc4&r=4k7iKXbjGC8LfYxVJJXiaYVu6FRWmEjX38S7JmlS9Vw&m=DjvKZYOEy46NTfpFc92OVXJsfNl6mHgZuzQYoA80LvQ&s=3pcGss3gRopgWpsIOGUaJCr3t5KJooO1atsssf2KC4w&e=
>  
> > 
> > Best regards, 
> > Konstantin Ladutenko 
> > 
> > 
> > 
> > -- 
> > The deal.II project is located at 
> https://urldefense.proofpoint.com/v2/url?u=http-3A__www.dealii.org_&d=CwIFaQ&c=Ngd-ta5yRYsqeUsEDgxhcqsYYY1Xs5ogLxWPA_2Wlc4&r=4k7iKXbjGC8LfYxVJJXiaYVu6FRWmEjX38S7JmlS9Vw&m=DjvKZYOEy46NTfpFc92OVXJsfNl6mHgZuzQYoA80LvQ&s=7bS84-mMZwrlyNuve21BTEi-eEz4iDqOpRq00ZCjQMg&e=
>  
> > For mailing list/forum options, see 
> > 
> https://urldefense.proofpoint.com/v2/url?u=https-3A__groups.google.com_d_forum_dealii-3Fhl-3Den&d=CwIFaQ&c=Ngd-ta5yRYsqeUsEDgxhcqsYYY1Xs5ogLxWPA_2Wlc4&r=4k7iKXbjGC8LfYxVJJXiaYVu6FRWmEjX38S7JmlS9Vw&m=DjvKZYOEy46NTfpFc92OVXJsfNl6mHgZuzQYoA80LvQ&s=H8P3DRu9PSgiI5yEqB9YN_F4LexIDdvjPYDm06dbW20&e=
>  
> > --- 
> > 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://urldefense.proofpoint.com/v2/url?u=https-3A__groups.google.com_d_optout&d=CwIFaQ&c=Ngd-ta5yRYsqeUsEDgxhcqsYYY1Xs5ogLxWPA_2Wlc4&r=4k7iKXbjGC8LfYxVJJXiaYVu6FRWmEjX38S7JmlS9Vw&m=DjvKZYOEy46NTfpFc92OVXJsfNl6mHgZuzQYoA80LvQ&s=UUGQFwTwiWuIW__4tomLbDzJbxwcubyUjW53IWswUeQ&e=
>  
> . 
>
>
>
> -- 
> Timo Heister 
> http://www.math.clemson.edu/~heister/ 
>

-- 
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.

Reply via email to