W dniu sobota, 28 października 2017 18:52:57 UTC+2 użytkownik Martin
Kronbichler napisał:
>
> Dear Michal,
>
> You mean because once you apply an inhomogeneous Dirichlet condition on
> the velocity you also get a contribution because the divergence is not
> zero?
Exactly.
> You could work around that by applying a function in the whole
> domain that satisfies the Dirichlet constraints and is divergence free,
> in which case you should not get a contribution to right right hand side
> of the divergence equation but only the momentum equation. On the other
> hand, I'm not sure how easy it would be to construct such a function...
>
Constructing such a function is not an option (for now), it is too
expensive.
One question about that way of imposing Dirichlet boundary conditions. I've
modified step-37 following the same idea, but in different way:
template <int dim>
void LaplaceProblem<dim>::solve ()
{
SystemMatrixType system_matrix_no_dirichlet;
{
typename MatrixFree<dim,double>::AdditionalData additional_data;
additional_data.tasks_parallel_scheme =
MatrixFree<dim,double>::AdditionalData::none;
additional_data.mapping_update_flags = (update_gradients |
update_JxW_values |
update_quadrature_points);
std::shared_ptr<MatrixFree<dim,double> >
system_mf_storage(new MatrixFree<dim,double>());
system_mf_storage->reinit (dof_handler,
constraints_without_dirichlet, QGauss<1>(fe.degree+1),
additional_data);
system_matrix_no_dirichlet.initialize (system_mf_storage);
}
system_matrix_no_dirichlet.evaluate_coefficient(Coefficient<dim>());
//
constraints.distribute(solution);
LinearAlgebra::distributed::Vector<double> residual=system_rhs;
system_matrix_no_dirichlet.vmult(residual, solution);
system_rhs -= residual;
LinearAlgebra::distributed::Vector<double> solution_update=system_rhs;
solution_update=0;
......
cg.solve (system_matrix, solution_update, system_rhs,
preconditioner);
solution+= solution_update;
constraints.distribute(solution);
In this way I don't have to modify my operator, bum I'm not sure if it will
work for all cases. For now it works and the results looks OK.
>
> Regarding the workaround: It could work by the modification you
> suggested, but we cannot really fix that in the library because the
> operators are meant for linear operators and not affine ones. The
> problematic issue is that there are some cases where you want to apply
> an operator with a non-zero constraint and others where you don't. For
> example, all Krylov solvers expect you to provide a matrix-vector for
> the linear operator, not the affine one as you would get when you have
> inhomogeneous constraints inside the mat-vec. In particular, I do not
> really understand the part in your second post where you talk about the
> vmult_interface_{up,down} functions because those are only used for
> multigrid where the library always assumes to work on homogeneous
> problems.
>
> You're right, forget about what I written previously. I've totally
misunderstood a few things.
> I have not gone down this route and as far as I know, nobody else has
> previously - including the matrix-based solvers we have available in
> deal.II whose way of operation resembles the setup I described in my
> previous post, but do it inside the
> ConstraintMatrix::distribute_local_to_global call - so I cannot say how
> much work that would be. We would be happy to provide a solution that is
> generic if you find one, but I'm not sure I will be able to help much
> along this direction.
>
> Best,
> Martin
>
--
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.