Albedo boundary conditions

154 views
Skip to first unread message

Antoni Vidal

unread,
Apr 22, 2013, 7:09:07 AM4/22/13
to dea...@googlegroups.com
Dear Deal.ii's users,
I'm solving an eigenvalue problem with a version of the step -36 using SLEPC to solve the eigenvalue problem.

And I would like to implement the following boundary condition :

d phi/dn = factor * phi 


I would be really grateful if you could help me with some idea in order to carry it out.


Thanks in advance.

Antoni Vidal

Guido Kanschat

unread,
Apr 22, 2013, 9:34:05 AM4/22/13
to deal.II user group

Dear Antoni,

this is what is usually called Robin boundary condition when you consider elliptic partial differential equations (which I assume you do).

The implementation is as follows: in your weak formulation, treat the boundary where you want to impose this condition like a Neumann boundary, but add the boundary integral

\int_{\Gamma_N} factor * u * v * ds

Best,
Guido




Antoni Vidal

unread,
May 8, 2013, 6:37:51 AM5/8/13
to dea...@googlegroups.com
Thank you Guido it worked as you said simply adding to my quadrature integral.

   for (unsigned int face=0; face<GeometryInfo<dim>::faces_per_cell; ++face)
  if (cell->face(face)->at_boundary())
  {
   fe_face_values.reinit (cell, face);
   for (unsigned int i=0; i<dofs_per_cell; ++i)
   for (unsigned int j=0; j<dofs_per_cell; ++j)
  for (unsigned int q_point=0; q_point<n_face_q_points; ++q_point)
  {
  cell_A_matrix(i, j) += (factor*
  fe_face_values.shape_value(i, q_point) *
  fe_face_values.shape_value(j, q_point)
  )*fe_face_values.JxW(q_point);
  }
  }

However, now , the error of the eigenvalues is 10 times greater (compared with the analytical ones) using the same grid and refination. Am I forgetting something? Should I made something with the mass_matrix (B matrix) ?

Thanks,
Antoni Vidal

Guido Kanschat

unread,
May 8, 2013, 12:12:02 PM5/8/13
to deal. II user group

DearAntoni,

what happens if you refine? Did the eigenvalues converge to the right value before? It is difficult to deduce from results on a single mesh. Did you try the boundary term with opposite sign? Is your analysis consistent?

Best, Guido

Antoni Vidal

unread,
May 10, 2013, 6:59:01 AM5/10/13
to dea...@googlegroups.com
Dear Guido,

what happens if you refine? Did the eigenvalues converge to the right value before? It is difficult to deduce from results on a single mesh. Did you try the boundary term with opposite sign? Is your analysis consistent?

The eigenvalues do not converge to the right values but they only differ in less than 0.01 whatever the factor is, so the difference is not really important. I've tried changing the FE degree, the refination and the quadrature formula. The eigenfunction also differ around 1%. 

Any idea?

Antoni Vidal
Reply all
Reply to author
Forward
0 new messages