Deal all,
I want to use matrix-free method to implement a phase field simulation, the 
control equation consists of the allen-cahn equations and  mechanics 
equation.
Is there any tutorial dealing with vectored-problem using matrix-free? as I 
have looked at the step-37,48,59, they are all scalar problem.
As a implementation routine, I want to solve single mechanics equation by 
matrix-free at first and modify base on the step-37 code, and there are the 
problems:

1) as it is a vector-valued problem need to set FEEvaluation template 
n_component to dim
FEEvaluation<dim, fe_degree, fe_degree+1, dim, number> phi(data);
then
phi.submit_dof_value()  should use type: 
Tensor<1,dim,VectorizedArray<number>> 
as the dof index in phi.submit_dof_value() goes  through 
phi.dofs_per_component rather phi.dofs_per_cell.
is that right?

2) since need to extract the diagonal values of stiffness matrix, 
in step-37 use 
phi.submit_dof_value(make_vectorized_array<number>(1.), i);
to extract the diagonal values of stiffness matrix while set other dof 
values to zero by
phi.submit_dof_value(VectorizedArray<number> 
<https://dealii.org/developer/doxygen/deal.II/classVectorizedArray.html>(), 
j); 
So in vector-valued case I wanted to add a unit vector(Tensor) 
(1,0,0)/(0,1,0)/(0,0,1) to dof values at first
but as the dof index in phi.submit_dof_value() goes  through 
phi.dofs_per_component rather phi.dofs_per_cell
what I can to is to submit a Tensor<1,dim,VectoriedArray> (1,1,1) to the 
data field.
is that right?
then the local_compute_diagonal(...) function would be:
template <int dim, int fe_degree, typename number>
  void MechanicsOperator<dim, fe_degree, number>::local_compute_diagonal(
    const MatrixFree<dim, number> &             data,
    LinearAlgebra::distributed::Vector<number> &dst,
    const unsigned int &,
    const std::pair<unsigned int, unsigned int> &cell_range) const
  {
    FEEvaluation<dim, fe_degree, fe_degree + 1, dim, number> phi(data);
    AlignedVector<Tensor<1,dim,VectorizedArray<number>>>
      diagonal(phi.dofs_per_component);
    for (unsigned int cell = cell_range.first; cell < cell_range.second; ++
cell)
      {
        phi.reinit(cell);
        for (unsigned int i = 0; i < phi.dofs_per_component; ++i)
        {
          for (unsigned int j = 0; j < phi.dofs_per_component; ++j)
          {
            phi.submit_dof_value(Tensor<1,dim,VectorizedArray<number
>>(), j);
          }
          Tensor<1,dim,VectorizedArray<number>> diagonal_tensor;
          for (unsigned int d=0; d<dim; d++)
            diagonal_tensor[d] = 1.;
          phi.submit_dof_value(diagonal_tensor, i);
          phi.evaluate(false, true);
          for (unsigned int q = 0; q < phi.n_q_points; ++q)
            {
              
phi.submit_symmetric_gradient(Cmatx*phi.get_symmetric_gradient(q),q); 
// while Cmatx is the 4-rank elastic tensor
            }   
          phi.integrate(false,true);
          diagonal[i] = phi.get_dof_value(i);         
        }
        for (unsigned int i = 0; i < phi.dofs_per_component; ++i)
          phi.submit_dof_value(diagonal[i],i);               
        phi.distribute_local_to_global(dst);
      }
  }

3) is the Multigrid still work well in vector-valued problem? or in 
multi-physics couple problem? 
I want to use it to solve the couple equations(Allen-Cahn and mechanics 
equalibrium equations) 
because it might be the easy one to be build.
I am still searching others' answers in the mailing list but have not found 
the answer about this topic( 1) and 2) ) yet(maybe).
Any suggestion or any materials related to this implementation would be 
appreciated .

best,
m.

-- 
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 on the web visit 
https://groups.google.com/d/msgid/dealii/8e8569e3-c4fc-4c18-a574-e6e9c7bf1fdc%40googlegroups.com.

Reply via email to