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

On Thursday, March 2, 2023 at 12:42:46 PM UTC+1 d.arnd...@gmail.com wrote:

> 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 
>
> On Wed, Mar 1, 2023 at 6:38 PM Deepika Kushwah <deepika...@iitgoa.ac.in> 
> wrote:
>
>> 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: 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. 
>> ************************************************************
>> ************************************ 
>>
>> -- 
>> 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
>>  
>> <https://groups.google.com/d/msgid/dealii/CAKTdrpp2Ckswjqi8RsU9DDwm6AkT-d05Cin-D33bWmfLv20Ytw%40mail.gmail.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 on the web visit 
https://groups.google.com/d/msgid/dealii/c5336533-06ac-460a-a2e4-9f3f575b545dn%40googlegroups.com.

Reply via email to