Re: [deal.II] Project cell quadrature point data to boundary

2024-02-23 Thread Abbas Ballout
Almonaci, 

I think you'll need the compute_projection_from_face_quadrature_points_matrix 
tfunction. 
Here's how you might use it:

std::vector distributed_field;
distributed_field.reinit(dof_handler.locally_owned_dofs(), 
mpi_communicator);

FullMatrix qpoint_to_dof_matrix(fe.dofs_per_face,
face_quadrature.size());

std::vector> local_values_at_qpoints, local_fe_values;
local_values_at_qpoints.reinit(face_quadrature.size());
local_fe_values.reinit(fe.dofs_per_cell);

for (const auto cell :dof_handler.active_cell_iterators())
for (unsigned int face = 0; face < GeometryInfo::faces_per_cell; ++
face)
if(face_at_boundary)
{
FETools::compute_projection_from_face_quadrature_points_matrix(fe,
face_quadrature,
face_quadrature,
cell,
face,
qpoint_to_dof_matrix);

for (unsigned int q = 0; q < face_quadrature->size(); ++q)
{
local_values_at_qpoints[q] = kirchof_stress_or_something;
}

qpoint_to_dof_matrix.vmult_add(local_fe_values,
local_values_at_qpoints);


cell->set_dof_values(local_fe_values,
distributed_field);

local_fe_values = 0;
}


The  distributed_field should have what you want in the end. Hope this 
helps. 

Abbas 

-- 
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 on the web visit 
https://groups.google.com/d/msgid/dealii/e83497f4-4336-487d-aebc-41875a7d36dbn%40googlegroups.com.


Re: [deal.II] Project cell quadrature point data to boundary

2024-02-23 Thread Wolfgang Bangerth

On 2/21/24 15:10, Javier Almonacid wrote:


Is there any method already implemented in deal.II to project cell-based 
quadrature point data onto the boundary? The end goal is to compute some 
forces on the surface in step-44 and for that, I need to project the Kirchhoff 
stresses (tau) onto the boundary.


An idea that occurred to me was to first project the quadrature point data 
into a more global element with DOFs (a finite element function) using one of 
the recommended extensions in step-18. However, I don't know if this is a good 
start because then the problem is how to project this object (coming from say, 
FE_DGQ) onto the boundary (these are all cell-wise operations by the way, as 
in step-18).


Javier,
do I understand correctly that you have some quantity (in your case the 
Kirchhoff stress) stored at quadrature points of each cell only, and you 
wonder how you could get an approximation of that quantity at points on the 
faces of that cell?


I can think of a number of ways of doing this, including interpolating or 
projecting onto a finite element field and evaluating it, but also just taking 
the cell-q_point that is closest to the face-q_point. In the end, can I ask 
you for a concise mathematical definition of what you want to do first, and 
then one can think of how to implement it?


Best
 W.

--

Wolfgang Bangerth  email: bange...@colostate.edu
   www: http://www.math.colostate.edu/~bangerth/


--
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 on the web visit 
https://groups.google.com/d/msgid/dealii/f17bd770-d481-4029-9007-0fdb9db66e9a%40colostate.edu.