Multi-physics implementation

87 views
Skip to first unread message

Ray Mclaren

unread,
Sep 21, 2020, 8:31:01 PM9/21/20
to deal.II User Group
Hi,
I have a question regrading multi-physics problem implementation.

I am developing a program that will solve the steady-state thermal conduction (laplace) first and then the allen-Cahn equation.

both of the problems are solved on the same mesh, however due to the nature of my specific implementation they solved separately. To solve the allen-Cahn equation, I need to evaluate temperature dependent functions at each node in the domain.

Hence, they only share the same triangulation, but everything else is separate.

My question is, how to use the solution vector from the thermal conduction (the temperature) in the assembly of the allen-cahn equation ?

 how to access the temperature at these nodes and make sure that they are synchronized ? 

thank you

simon...@gmail.com

unread,
Sep 22, 2020, 4:01:42 AM9/22/20
to deal.II User Group
Hi,
I would do it like this:

  // Setup two DoFHandlers, one for the heat and one for Allen-Cahn
  DoFHandler<dim> heat_dof_handler(triangulation);
  DoFHandler<dim> allen_cahn_dof_handler(triangulation);

   // Solve for this first
  Vector<double> heat_solution;

  // Lots of things in between here...

  // When assembling for allen_cahn, create two FEValues objects.
  FEValues<dim> fe_values(allen_cahn_element, quadrature, update_flags);
  FEValues<dim> heat_fe_values(heat_element, quadrature, update_values);

  for (auto cell : allen_cahn_dof_handler.active_cell_iterators())
    {
      // Cast to an iterator on the heat_dof_handler
      typename DoFHandler<dim>::active_cell_iterator heat_cell(
        &triangulation, cell->level(), cell->index(), &heat_dof_handler);

      // Call reinit with the corresponding cells.
      fe_values.reinit(cell);
      heat_fe_values.reinit(heat_cell);

      // Get the heat values at the quadrature points using heat_fe_values
      std::vector<double> heat_values(heat_element.dofs_per_cell);
      heat_fe_values.get_function_values(heat_solution, heat_values);

      // Assemble the system using fe_values here ...
    }

Best,
Simon

Ray Mclaren

unread,
Sep 22, 2020, 10:03:34 AM9/22/20
to deal.II User Group
Mr. Simon

thank you for your response, It was exactly what I was looking for. 

Jean-Paul Pelteret

unread,
Sep 22, 2020, 3:55:04 PM9/22/20
to dea...@googlegroups.com
Hi Ray,

I’m happy to see that Simon has already given a great explanation as to how to solve your problem. It’s probably more on point than the code that mentioned that I would direct you towards. So it really was worthwhile asking the question again here :-)

Best,
Jean-Paul

-- 
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+un...@googlegroups.com.
To view this discussion on the web visit https://groups.google.com/d/msgid/dealii/391eb569-24d4-4ea4-9dc6-b9eada95f9edn%40googlegroups.com.

Ray Mclaren

unread,
Sep 23, 2020, 6:21:37 AM9/23/20
to deal.II User Group
Mr. Pelteret 

Sir, thank you for your guidance. 

I am new to using dealii and I am very thankful for your support

Thanks 

Reply all
Reply to author
Forward
0 new messages