Ted,
I assume your dirichlet_value=2.0 and neumann_value=0.5 on your
assembly routines for the boundary ? If yes, the boundary assembly is
not correct. A generic code for bc specification would look like :
for (unsigned int side=0; side<elem->n_sides(); side++)
if (elem->neighbor(side) == NULL)
{
const Real penalty = 1.e10;
const std::vector<std::vector<Real> >& phi_face =
fe_face->get_phi();
const std::vector<Real>& JxW_face = fe_face->get_JxW();
const std::vector<Point >& qface_point = fe_face->get_xyz();
fe_face->reinit(elem, side);
for (unsigned int qp=0; qp<qface.n_points(); qp++)
{
const Real xf = qface_point[qp](0);
const Real yf = qface_point[qp](1);
const Real zf = qface_point[qp](2);
if (bc_id == 0) {
const Real value = exact_solution(xf, yf, zf);
for (unsigned int i=0; i<phi_face.size(); i++)
for (unsigned int j=0; j<phi_face.size(); j++)
Ke(i,j) +=
JxW_face[qp]*penalty*phi_face[i][qp]*phi_face[j][qp];
for (unsigned int i=0; i<phi_face.size(); i++)
Fe(i) +=
JxW_face[qp]*penalty*value*phi_face[i][qp];
}
else if (bc_id == 1)
{
const Real value = neumann_value(xf, yf, zf);
for (unsigned int i=0; i<phi_face.size(); i++)
Fe(i) += JxW_face[qp]*value*phi_face[i][qp];
}
}
Most part of it is taken from ex4.C. Substitute in the neumann_value
for the side and you should get the right solution. Just look at the
variational form for your problem and this should be easy to see.
Vijay
On Mon, Aug 31, 2009 at 6:57 PM, Ted Kord<[email protected]> wrote:
>>
>> Well, to impose a Neumann BC on part of the boundary you should only add
>> Neumann BC terms to the entries of the rhs vector corresponding to dofs on
>> the relevant part of the boundary. So you'd need to detect which part of the
>> boundary you're on and add the Neumann BC terms accordingly. Have a look at
>> ex13 for an example usage of boundary IDs for this kind of purpose. For more
>> complicated domains, I typically generate a mesh using gmsh, which enables
>> you to specify boundary IDs (using a gmsh "physical").
>>
>> In principle 1D is the same except that in that simpler case the Neumann
>> term is not an integral so you don't need a boundary quadrature rule etc...
>> Also there are only 2 boundary elements in a 1D mesh so you could easily
>> determine which boundary element you're on by testing the value of
>> elem->centroid() or something.
>>
>> If you can't get it to work, then just send your BC code fragment to the
>> list...
>>
>> - Dave
>>
>
>
> This is what I've tried. The Dirichlet BC works fine but the plot of the
> solution implies that something's terriby wrong with the Neumann B.C :
>
> for(unsigned int s=0; s<elem->n_sides(); s++)
> {
> {
> Ke(s,s) += penalty;
>
> short int bc_id = mesh.boundary_info->boundary_id (elem,s);
>
> if (bc_id == 0) //Left end : Dirichlet BC
> Fe(s) += 2*penalty;
> else if (bc_id == 1) //Right end : Neumann BC
> Fe(s) += - 0.5;
>
> }
> }
>
> What am I doing incorrectly?
>
> Thanks
>
> Ted
> ------------------------------------------------------------------------------
> Let Crystal Reports handle the reporting - Free Crystal Reports 2008 30-Day
> trial. Simplify your report design, integration and deployment - and focus on
> what you do best, core application coding. Discover what's new with
> Crystal Reports now. http://p.sf.net/sfu/bobj-july
> _______________________________________________
> Libmesh-users mailing list
> [email protected]
> https://lists.sourceforge.net/lists/listinfo/libmesh-users
>
------------------------------------------------------------------------------
Let Crystal Reports handle the reporting - Free Crystal Reports 2008 30-Day
trial. Simplify your report design, integration and deployment - and focus on
what you do best, core application coding. Discover what's new with
Crystal Reports now. http://p.sf.net/sfu/bobj-july
_______________________________________________
Libmesh-users mailing list
[email protected]
https://lists.sourceforge.net/lists/listinfo/libmesh-users