Dear Franck,
> I am looking at step-43 using different method. I do a semi discretization
> in space and end up with an ode for the saturation and use different
> technique for time integration for the saturation.
>
> For this I will need to assemble the matrix H after solving for the
> velocity in the darcy system. Individual entries of H are given by:
>
> \mathbf{H}_{ij} & = - \Delta t^{(n)}_c \Big( F\left(S^{(n-1)}\right)
> \mathbf{v}_i,\nabla\phi_j\Big)_{\Omega}
>
>
> I am struggling with sparsity pattern that I have to define for this
> matrix. I try this:
>
> assemble_matrix_H.clear();
> DynamicSparsityPattern dsp (n_s, n_u);
> DoFTools::make_sparsity_pattern (saturation_dof_handler,
> darcy_dof_handler, dsp,
> saturation_constraints, false);
> ssemble_matrix_H.reinit(dsp);
>
>
>
> dealii is complaining and is saying
>
> /home/franckm/apps/candi/deal.II-toolchain/deal.II-v8.4.0/include/deal.II/dofs/dof_tools.h:563:3:
>
> note: template<class DoFHandlerType, class SparsityPatternType> void
> dealii::DoFTools::make_sparsity_pattern(const DoFHandlerType&, const
> DoFHandlerType&, SparsityPatternType&)
> make_sparsity_pattern (const DoFHandlerType &dof_row,
> ^
> /home/franckm/apps/candi/deal.II-toolchain/deal.II-v8.4.0/include/deal.II/dofs/dof_tools.h:563:3:
>
> note: template argument deduction/substitution failed:
> /home/franckm/candi-examples/step-43b/step-43b.cc:918:68: note:
> candidate expects 3 arguments, 5 provided
> saturation_constraints, false);
> ^
> make[3]: *** [CMakeFiles/step-43b.dir/step-43b.cc.o] Error 1
>
The error is not surprising as there isn't a variant of
DoFTools::make_sparsity_pattern that takes two DoFHandlers.
> How one define a sparsity pattern for assembling H in the above scenario?
>
What you can do, is to consider all the blocks of your system and use
\widetilde{H}=[0 H;0 0] instead. In this case your DoFHandler would have to
have two blocks and
you can use the variant [1] instead with an appropriate coupling.
A different approach is to build the SparsityPattern yourself by simulating
assembling the matrix. This means that instead of using
ConstraintMatrix.distribute_local_to_global
you would call SparsityPattern.add [2].
The question is if H really is the matrix you want to assemble. How does
the equation you are considering really look like?
Best,
Daniel
[1]
https://www.dealii.org/8.4.1/doxygen/deal.II/group__constraints.html#ga33ec1ee2ac7c86b16700022b704b69c8
[2]
https://www.dealii.org/8.4.0/doxygen/deal.II/classSparsityPattern.html#aa8dd1f7f393272cc5ae7789ac6be3413
--
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.