Dear Subramanya,

The actual cells associated with the integer indices 'cell' (that are batches in case of SIMD vectorization) are not necessarily the same if you create two separate MatrixFree objects and initialize each with the respective DoFHandler for the two fields. If you intend to use coupled problems like the thermal-stress problem you describe, you should use the MatrixFree::reinit function that takes both DoFHandler objects as the arguments to an std::vector<DoFHandler*> and respectively for the affine constraints. Then, you can initialize FEEvaluation objects using the 'dof_no' argument in the constructor: https://www.dealii.org/developer/doxygen/deal.II/classFEEvaluation.html#ade40e452949cbefd4f7e6d5580af81c8 - it is the second argument next to MatrixFree.

You can also try to run with separate MatrixFree objects, but then you need to manually link the two objects and one FEEvaluation needs to be initialized with an array of cells (that you constructed manually). I have done that in some applications, and there is at least one code inside deal.II that does this: https://github.com/dealii/dealii/blob/8f97fd0dca1eb878bc7ae9922ac89430cab800db/include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h#L4041-L4068 - but it is much less obvious, so I really recommend the first variant.

I hope this helps!

Best regards,
Martin

Am 03.01.25 um 17:14 schrieb subramanya gautam:
I am trying to implement a staggered thermal-elastic solver using the matrix-free multigrid solver approach from step-37. I got an elasticity solver and a thermal solver working and now I am trying to associate the two. I have a shared triangulation, with two separate dof_handlers for the  thermal and stress problems. The question I have is when I do the assembly loop, does the cell refer to the same cell for both the thermal and  elastic problems? So from step-37,  can I just call phi_temp.reinit(cell), and phi_stress.reinit(cell), and have cell refer to the same physical cell in the triangulation?


for (unsigned int cell = 0;
  cell < system_matrix.get_matrix_free()->n_cell_batches();
  ++cell)
  {
  phi.reinit(cell);
  for (const unsigned int q : phi.quadrature_point_indices())
  phi.submit_value(make_vectorized_array <https://www.dealii.org/current/doxygen/deal.II/vectorization_8h.html#ad7d7e08942faeecf438c75a254e06cbe><double>(1.0), q);   phi.integrate(EvaluationFlags::values <https://www.dealii.org/current/doxygen/deal.II/namespaceEvaluationFlags.html#a9b7c6d689cb76386839d0d13640f59aeaf9825c682f693a6a200094641a0d6a58>);
  phi.distribute_local_to_global(system_rhs);
  }
Thanks!
Subramanya.
--
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 dealii+unsubscr...@googlegroups.com. To view this discussion visit https://groups.google.com/d/msgid/dealii/d8f59b24-2a30-4ade-83e4-75bc72112afen%40googlegroups.com <https://groups.google.com/d/msgid/dealii/d8f59b24-2a30-4ade-83e4-75bc72112afen%40googlegroups.com?utm_medium=email&utm_source=footer>.

--
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 dealii+unsubscr...@googlegroups.com.
To view this discussion visit 
https://groups.google.com/d/msgid/dealii/71ea7f4c-19c8-43f6-9037-775aa504b470%40rub.de.

Reply via email to