Hello,
I've been working on a phase field code using the matrix free approach from
Step 48. When I use adaptive meshing, I notice some artifacts at the
boundary between levels of refinement and I'm trying to understand their
origin.
For example, if I make the following modification to the "local_apply"
function in Step 48, the solution should remain unchanged from time step to
time step:
template <int dim, int fe_degree>
void SineGordonOperation<dim, fe_degree>::
local_apply (const MatrixFree<dim> &data,
parallel::distributed::Vector<double> &dst,
const std::vector<parallel::distributed::Vector<double>*> &
src,
const std::pair<unsigned int,unsigned int> &cell_range) const
{
AssertDimension (src.size(), 2);
FEEvaluation<dim,fe_degree> current (data), old (data);
for (unsigned int cell=cell_range.first; cell<cell_range.second; ++cell)
{
current.reinit (cell);
old.reinit (cell);
current.read_dof_values (*src[0]);
old.read_dof_values (*src[1]);
current.evaluate (true, true, false);
old.evaluate (true, false, false);
for (unsigned int q=0; q<current.n_q_points; ++q)
{
const VectorizedArray<double> current_value = current.get_value(
q);
const VectorizedArray<double> old_value = old.get_value(q);
// ORIGINAL
// current.submit_value (2.*current_value - old_value -
// delta_t_sqr *
std::sin(current_value),q);
// current.submit_gradient (- delta_t_sqr *
// current.get_gradient(q), q);
// CHANGES
current.submit_value (current_value,q);
current.submit_gradient (- delta_t_sqr *
current.get_gradient(q)*0.0, q);
}
current.integrate (true,true);
current.distribute_local_to_global (dst);
}
}
However, there is a non-trivial change in the solution (same plot with and
without the mesh):
<https://lh3.googleusercontent.com/-3PmzXmsz0WM/WbKjJ3tOnZI/AAAAAAAAA50/jPaQpAd8uW82BTfgtN6MiO3ZuLkycer2QCLcBGAs/s1600/step-48_artifact.tiff>
<https://lh3.googleusercontent.com/-p_c2uy9fmno/WbKjYO-dJyI/AAAAAAAAA54/1wmh__EsxOc2WYJ3smlia6qUEdd9gfvZQCLcBGAs/s1600/step-48_artifact_n0_mesh.tiff>
Has anyone noticed this before? In this example, the magnitude of the
artifacts is low, but in some of my other tests, the magnitude has been on
the order of the "real" change in the solution. When doing implicit time
stepping (still in the matrix free framework), the artifacts don't show up,
which suggests that they aren't just a consequence of having an adapted
mesh. In the Step 48 description, it's noted that the inverse of the mass
matrix isn't computed exactly. Could the artifacts be related to that?
For Step 48, the artifacts only appear from the initial time step to the
first time step (in my application, I'm getting the artifacts every time
step and I'm still in the midst of finding what the source of that
difference is).
Thank you for any insight you can give.
Best,
Steve
Stephen DeWitt
Asst. Research Scientist
PRISMS Center
Dept. of Materials Science and Engineering
University of Michigan
--
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].
For more options, visit https://groups.google.com/d/optout.