Hello Isa,
perhabs you want to take look on the FE_Nothing element. With this you
can create a FECollection consisting of, for example, Taylor-Hood or
Raviart-Thomas elements (for Navier-Stokes), FE_Nothing (for the cells,
where nothing should happen) and FE_Q (for Laplace) and activate that
finite element on the cell, which you require for your current PDE.
Doing it this way, you do not have to set the degrees of freedom on the
unwanted part of the domain to zero explicitly, because there are just
no degrees of freedom anymore.
Best Regards,
Markus
Am 28.08.10 18:10, schrieb [email protected]:
Hello everybody!
I have a 2D geometry with different domains (each domain has its own
material value in order to be identified).
This geometry has associated its “Triangulation” and its “DoFHandler”.
The geometry and its material values is the following one:
-----------------------------------------------------------------
| | | | | | | |
| | | | | | | |
| 1 | 2 | 3 | 4 | 5 | 6 | 7 |
| | | | | | | |
-----------------------------------------------------------------
I want to calculate Navier-Stokes equations in domains 1,2,3,5,6 and 7.
So, I don't want to take into account the domain 4. That is to say, for
this system the line between domain 3 and 4 is a boundary and also the
line between domain 4 and 5.
To calculate the assemble of matrices I do:
...
for (; cell!=endc; ++cell)
{
if ((celltria->material_id()==1))
...
else if ((celltria->material_id()==2))
...
else if ((celltria->material_id()==3))
…
else if ((celltria->material_id()==5))
...
else if ((celltria->material_id()==6))
...
else if ((celltria->material_id()==7))
...
else if ((celltria->material_id()==4))
{
//Here I don't do anything because I don't want calculate nothing.
}
cell-> get_dof_indices(local_dof_indices);
for (unsigned int i=0;i<dofs_per_cell;++i)
{
for (unsigned int j=0;j<dofs_per_cell;j++)
general.add(local_dof_indices[i],
local_dof_indices[j],
cell_matrix(i,j));
vector_solution_rhs(local_dof_indices[i]) += cell_rhs(i);
}
celltria++;
}
//Here, I put all values inside domain 4 with a value of zero:
DoFTools::extract_subdomain_dofs(dof_handler,4, selected_dofsMEM);
for (unsigned int i=0;i<dof_handler.n_dofs();i++)
if ((selected_dofsMEM[i]==true))
boundary_values[i]=0.0;
//Apply the boundary values to my system:
MatrixTools::apply_boundary_values(boundary_values,general,
vector_solution,vector_solution_rhs,false);
I am not sure if this way is the correct way to proceed.
On the other hand, I want to calculate a Laplace equation in domains 3,4
and 5. So, I don't want to take into account the domains 1,2,6 and 7. So,
the line between 2 and 3 would be a boundary and the same to line between
5 and 6.
I was asking me if there are a way in other to calculate a particular
equation only in some domains belonging to the same geometry.
Thanks in advance.
Best.
Isa
_______________________________________________
dealii mailing list http://poisson.dealii.org/mailman/listinfo/dealii
_______________________________________________
dealii mailing list http://poisson.dealii.org/mailman/listinfo/dealii