Najwa:
if you have a discontinuous viscosity in the Stokes equations, I think you
want to read this paper:
https://www.math.colostate.edu/~bangerth/publications/2017-boussinesq.pdf
Specifically, section 3.3 talks about this situation, and it shows in figures
2 and 4 that you do get spikes in the pressure.
Best
W.
On 4/24/25 05:17, Najwa Alshehri wrote:
*** Caution: EXTERNAL Sender ***
Dear team,
I am seeking your help in understanding how few things work in dealii.
I am solving Stokes interface problem where viscosity differs in two parts of
the domain ( outer part and immersed circle). I am solving using unfitted
boundary approach ( boundary between domains do not align with the mesh). I am
using Q2-Q1 Finite element.
I use the following data_out:
std::vector<std::string> solution_names(dim, "velocity");
solution_names.emplace_back("pressure");
std::vector<DataComponentInterpretation::DataComponentInterpretation>
data_component_interpretation(
dim, DataComponentInterpretation::component_is_part_of_vector);
data_component_interpretation.push_back(
DataComponentInterpretation::component_is_scalar);
const MappingQ<dim> mapping(1);
DataOut<dim> data_out;
data_out.attach_dof_handler(dof_handler);
data_out.add_data_vector(error_u, "error_u");
data_out.add_data_vector(error_p, "error_p");
data_out.add_data_vector(solution,
solution_names,
DataOut<dim>::type_dof_data,
data_component_interpretation);
data_out.add_data_vector(error_per_cell_omega, "error_indicator_omega");
data_out.build_patches(mapping, fe.degree + 1,
DataOut<dim>::curved_inner_cells);
I removed the mean value of the approximated pressure in different ways. One
of them is just like in step-56
const double mean_pressure = VectorTools::compute_mean_value(
dof_handler, QGauss<dim>(2 * pressure_degree + 10), solution, dim);
solution.block(1).add(-mean_pressure);
I noticed that the solution of the pressure have some (unexpected)
oscillations at the interface even though the exact solution of the pressure
is a smooth function (x^2-y^2). In fact, if I just look at the solution in 2d
using heat map, those oscillations are not clear.
2d_p.png
However, in ParaView, when I use Warp by a scalar filter, they are so clear.
3d_p.png
My questions are,
Is this "compute mean value" works with any order of finite element?
Am I using "build batches" wrong?
Is this a common behaviour if the coefficients jumps?
Any other suggestions are appreciated.
Thank you so much,
Najwa
--
The deal.II project is located at http://www.dealii.org/
<http://www.dealii.org/>
For mailing list/forum options, see
https://groups.google.com/d/forum/dealii?hl=en
<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 [email protected]
<mailto:[email protected]>.
To view this discussion visit
https://groups.google.com/d/msgid/dealii/d9a8c115-f363-42b5-bef2-9b26b62ae616n%40googlegroups.com <https://groups.google.com/d/msgid/dealii/d9a8c115-f363-42b5-bef2-9b26b62ae616n%40googlegroups.com?utm_medium=email&utm_source=footer>.
--
------------------------------------------------------------------------
Wolfgang Bangerth email: [email protected]
www: http://www.math.colostate.edu/~bangerth/
--
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 [email protected].
To view this discussion visit
https://groups.google.com/d/msgid/dealii/a422fc45-ff8a-423d-9c04-613343d0dbb2%40colostate.edu.