Queries on calculation on average stress.

조회수 49회
읽지 않은 첫 메시지로 건너뛰기

Deepika Kushwah

읽지 않음,
2023. 3. 1. 오후 12:38:2123. 3. 1.
받는사람 dea...@googlegroups.com
Hello Everyone,

I have the following queries :
1. How can I calculate the average stress and strain values for the following geometry (Where element size is not uniform) using the Data PostProcedure?

image.png
2. Is it possible to make a uniform mesh so the size of elements in both (inclusion and matrix) will be the same?

Thanks & Regards
Deepika

**************************************************************************
This e-mail is for the sole use of the intended recipient(s) and may
contain confidential and privileged information. If you are not the
intended recipient, please contact the sender by reply e-mail and destroy
all copies and the original message. Any unauthorized review, use,
disclosure, dissemination, forwarding, printing or copying of this email
is strictly prohibited and appropriate legal action will be taken.
************************************************************************************************

Daniel Arndt

읽지 않음,
2023. 3. 2. 오전 6:42:4623. 3. 2.
받는사람 dea...@googlegroups.com
Deepika,

> How can I calculate the average stress and strain values for the following geometry (Where element size is not uniform) using the Data PostProcedure?

Can you describe what you are trying to do in a formula? DataPostProcessor typically gives you a result per grid point in the output step.

> 2. Is it possible to make a uniform mesh so the size of elements in both (inclusion and matrix) will be the same?

If you can draw a mesh that you would like to use, you should also be able to create that in deal.II to use. Of course, it will be very difficult to achieve the same size for all elements since you will need more cells to describe the inclusions anyway.

Best,
Daniel

--
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/CAKTdrpp2Ckswjqi8RsU9DDwm6AkT-d05Cin-D33bWmfLv20Ytw%40mail.gmail.com.

Abbas

읽지 않음,
2023. 3. 2. 오전 9:09:4423. 3. 2.
받는사람 deal.II User Group
Daniel, 

Side question: 
If I have to output gradients of a variable, the straightforward approach would be to evaluate the gradients over quadrature points and then average. This gives me an average value of my gradient per cell.  
Do you know how does data DataPostProcessor evaluates the value of the gradient at the mesh nodes? How does it it work under the hood? Do you mind explaining?  

Deepika, 

This is a strain post processor. Create a function that does this: 
template <int dim>
class StrainPostprocessor : public DataPostprocessorTensor<dim>
{
public:
StrainPostprocessor()
: DataPostprocessorTensor<dim>("strain",
update_gradients)
{
}

virtual void evaluate_vector_field(const DataPostprocessorInputs::Vector<dim> &input_data, std::vector<Vector<double>> &computed_quantities) const override
{
AssertDimension(input_data.solution_gradients.size(),
computed_quantities.size());

for (unsigned int p = 0; p < input_data.solution_gradients.size(); ++p)
{
AssertDimension(computed_quantities[p].size(),
(Tensor<2, dim>::n_independent_components));

for (unsigned int d = 0; d < dim; ++d)
for (unsigned int e = 0; e < dim; ++e)
computed_quantities[p][Tensor<2, dim>::component_to_unrolled_index(TableIndices<2>(d, e))] = (input_data.solution_gradients[p][d][e] +
input_data.solution_gradients[p][e][d]) /
2;
}
}
};

Then in data output you can do this: 
void CG_Elasticity_Linear<dim>::output_results(unsigned int refinement_level, unsigned int cycle) const
{

std::vector<DataComponentInterpretation::DataComponentInterpretation>
interpretation(dim,
DataComponentInterpretation::component_is_part_of_vector);
std::vector<std::string> solution_names(dim, "u");

DataOut<dim> data_out;

data_out.add_data_vector(dof_handler, solution, solution_names, interpretation);
const StrainPostprocessor<dim> strain;
data_out.add_data_vector(dof_handler, solution, strain);
data_out.build_patches();

std::ofstream output("solution" + std::to_string(refinement_level) + "_" + std::to_string(cycle) + ".vtu");

data_out.write_vtu(output);

}
For averages you'll have to loop in the DataPostProcessor derived class and average your quantity.  


Best,
Abbas

Deepika Kushwah

읽지 않음,
2023. 3. 3. 오전 12:58:2123. 3. 3.
받는사람 dea...@googlegroups.com
Hello Everyone,

Thanks for your response.

The problem is how to calculate the average strain and stress values using Area Average Method (for the geometry attached in previous mail)?

I am using this formula Avg Strain = (A1E1+A2E2.....AnEn/A) for the calculation. I have calculated the values of E1,E2.. En using DataPostProcedure Class but now the problem is how to calculate A1,A2....An in the same class as active.cell_iterator showing the error?

What I have done so far:
1. I have calculated strain and stress using DataPostProcedure Class

 template <int dim>
class StrainPostprocessor : public DataPostprocessorTensor<dim>
{
public:
  StrainPostprocessor ()
    :
    DataPostprocessorTensor<dim> ("strain",
                                  update_gradients)
  {}
 
  virtual
  void
  evaluate_vector_field
  (const DataPostprocessorInputs::Vector<dim> &input_data,
   std::vector<Vector<double> > &computed_quantities) const override
 
  {
    AssertDimension (input_data.solution_gradients.size(),
                     computed_quantities.size());
         
   
         
    for (unsigned int p=0; p<input_data.solution_gradients.size(); ++p)
      {
        AssertDimension (computed_quantities[p].size(),
                         (Tensor<2,dim>::n_independent_components));
     
        for (unsigned int d=0; d<dim; ++d)
          for (unsigned int e=0; e<dim; ++e)

            computed_quantities[p][Tensor<2,dim>::component_to_unrolled_index(TableIndices<2>(d,e))]
              = (input_data.solution_gradients[p][d][e]
                 +
                 input_data.solution_gradients[p][e][d]) / 2;
     }
   }          
};

2. I have calculated area of each cell using active.cell_iterator separately 

  Triangulation<2> triangulation;
 
  std::vector<double>  x(4);
  std::vector<double>  y(4);
  std::vector<double>  A(4);
  GridGenerator::hyper_cube(triangulation,-1,1);
 //double A;
// double B;
 // double C;
 
triangulation.refine_global(1);

int j =0;
  for (const auto &cell : triangulation.active_cell_iterators())
    {
       int i =0;
       
          for (const auto v : cell->vertex_indices())
           
             {
             x[i] = cell->vertex(v)(0);
             y[i] = cell->vertex(v)(1);
             
              std::cout<< "cord_x" << x[i]<<std::endl;
              std::cout<< "cord_y" << y[i]<<std::endl;
               i = i+1 ;
              }
             
             A[j] = (1./2.)*((x[0]*y[1] + x[1]*y[3] + x[3]*y[2] + x[2]*y[0]) - (x[1]*y[0] + x[3]*y[1] + x[2]*y[3] + x[0]*y[2]));
               
             std::cout<<"area"<<""<< A[j] <<std::endl;
             j=j+1;
         }

what I want to do:
Now I want to calculate both strain/stress and area in a single go so that we can get the product of these two.

Kindly suggest if there is any other method to calculate avg strain/stress using composite material where the size of elements in matrix and inclusion are different.



Thanks,
Deepika

Wolfgang Bangerth

읽지 않음,
2023. 3. 3. 오후 1:33:0823. 3. 3.
받는사람 dea...@googlegroups.com
On 3/2/23 22:57, Deepika Kushwah wrote:
>
> Kindly suggest if there is any other method to calculate avg strain/stress
> using composite material where the size of elements in matrix and inclusion
> are different.

DataPostprocessor is the wrong approach to the problem. It is used to generate
graphical output, but that is not what you want.

The average strain is

s_avg = \int_\Omega s(x) dx / \int_\Omega 1 dx

where the second integral is simply the volume. You can approximate both
integrals by a sum over all cells and then applying quadrature. For example,
the first integral is computed via a loop that will look roughly like this:

double s_integral = 0;
for (cell=....)
{
fe_values.reinit(cell);
fe_values.get_function_gradients (solution, solution_grads);
for (q=...)
s_integral += (compute strain/stress from solution_grads[q]) *
fe_values.JxW(q);
}

The volume of the domain is computed similarly, but you can also use
GridTools::volume() for that.

Best
W.

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


전체답장
작성자에게 답글
전달
새 메시지 0개