Hi Wolfgang,


Thank you for the quick response.


I am fine with writing the displacement and the crack phase-field in 
separate files, as long as the objective of visualizing the phase-field on 
the deformed mesh is achieved.


Listed below are two possible ways of doing this (along with the 
issues/queries regarding them).


1. Write the phase-field in a separate file using a MappingQEulerian type 
object based on the displacement with nSubdivisions equaling the polynomial 
degree for FE ansatz of phase-field. The problem with this approach is that 
the command


dataOutDisp.add_data_vector(solution.block(0),

dispNames,

dealii::DataOut<dim>::type_dof_data,

dataCompIntepretationDisp);

dataOutDisp.build_patches(qMapping, polyDegrees[0]);


results in the following error:


An error occurred in line <1069> of file 
</calculate/iwtm99_lokal/Paras_HomeOffice/spack/var/spack/stage/dealii-9.0.1-jem5nejch2egpiislts7pg4o767wpuhu/dealii-9.0.1/include/deal.II/numerics/data_out_dof_data.templates.h>
 
in function

void dealii::DataOut_DoFData<DoFHandlerType, patch_dim, 
patch_space_dim>::add_data_vector_internal(const DoFHandlerType*, const 
VectorType&, const std::vector<std::__cxx11::basic_string<char> >&, 
dealii::DataOut_DoFData<DoFHandlerType, patch_dim, 
patch_space_dim>::DataVectorType, const 
std::vector<dealii::DataComponentInterpretation::DataComponentInterpretation>&, 
bool) [with VectorType = dealii::Vector<double>; DoFHandlerType = 
dealii::DoFHandler<2, 2>; int patch_dim = 2; int patch_space_dim = 2]

The violated condition was:

data_vector.size() == dof_handler->n_dofs()

Additional information:

The vector has size 132372 but the DoFHandler object says that there are 
198558 degrees of freedom and there are 65731 active cells. The size of 
your vector needs to be either equal to the number of degrees of freedom, 
or equal to the number of active cells.


A probable solution could be to create a new DOfHandler during each output 
call for both the displacement as well as the phase-field and attach the 
respective solution blocks to their DofHandlers. Not sure this would work 
since the mapping which would be used for the phase-field would be based on 
a different DoFHandler. Also, this would probably be computationally 
expensive (not sure?)


2. Another idea could be to write both of them in the same file as done in 
Step-44, but use nSubdividisions equaling either min or max of the 
polynomial degrees of the FE ansatz for the two variables. Since, I 
probably not yet fully understand the implication of the parameter 
nSubdivisions in the definition of build_patches() 
[https://www.dealii.org/9.0.0/doxygen/deal.II/classDataOut.html#a5eb51872b8736849bb7e8d2007fae086],
 
could you please explain what happens in the below example.


Let pDisp=2 and pPF=3 (it is known that pPF >= pDisp always) be the 
polynomial degrees used for the FE ansatz of displacement and phase-field 
respectively. 


a) What information would we loose for phase-field, if nSubDivisions=pDisp, 
since, there are 16 values available per cell, and we only visualize 9 in 
PARAVIEW? Are the non-vertex values used for interpolation of the 
phase-field on the 4 sub-cells or are they completely neglected?


b) If we use nSubDivisions=pPF, then we would have the best “possible” 
visualization for the phase-field, but how would the values of the 
displacement field be interpolated for the 27 sub-cells since in reality 
only values for 9 points per cell are known. Could this lead to incorrect 
visualization of displacement?


Best regards,

Paras 


On Friday, August 7, 2020 at 8:15:06 PM UTC+2, Wolfgang Bangerth wrote:
>
> On 8/7/20 11:20 AM, Paras Kumar wrote: 
> > 
> > I have a query regarding writing the solution to a VTK file for 
> visualization 
> > with PARAVIEW.  If one wanted to do this via a mapping based on the 
> > displacement field, then one could simply follow the output_results() 
> function 
> > of step-44. 
> > My question is that, Is it possible to write to the same vtk file, 
> fields 
> > based on different MappingQEulerian 
> > <
> https://nam01.safelinks.protection.outlook.com/?url=https%3A%2F%2Fwww.dealii.org%2F9.0.0%2Fdoxygen%2Fdeal.II%2FclassMappingQEulerian.html&data=02%7C01%7CWolfgang.Bangerth%40colostate.edu%7Ccbd1a9a17b8540f0018708d83af62a7b%7Cafb58802ff7a4bb1ab21367ff2ecfc8b%7C0%7C0%7C637324176247361942&sdata=0fDQN1511kXggwdxcB8jf%2FO05Wn8ghyd3QzReWpRrqc%3D&reserved=0>
>  
>
> > objects and with different values of n_subdivsions for the 
> > data_out.build_patches 
> > <
> https://nam01.safelinks.protection.outlook.com/?url=https%3A%2F%2Fwww.dealii.org%2F9.0.0%2Fdoxygen%2Fdeal.II%2FclassDataOut.html%23a5eb51872b8736849bb7e8d2007fae086&data=02%7C01%7CWolfgang.Bangerth%40colostate.edu%7Ccbd1a9a17b8540f0018708d83af62a7b%7Cafb58802ff7a4bb1ab21367ff2ecfc8b%7C0%7C0%7C637324176247371938&sdata=d1P5xlbOdyblD2C38dcCkBaFpm2sCqSI8HWRFvKNYtw%3D&reserved=0>()
>  
>
> > function? 
>
> Paras: 
> No, you can't. But there is nothing that prevents you from creating two 
> .vtu 
> files that contain different things, and visualizing the at the same time. 
> Both paraview and visit allow you to open several .vtu files at the same 
> time 
> and showing their contents in the same window. 
>
> Best 
>   W. 
>
>
> -- 
> ------------------------------------------------------------------------ 
> Wolfgang Bangerth          email:                 bang...@colostate.edu 
> <javascript:> 
>                             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/cc27883c-3c4e-4d9d-afca-a5f47d9ba15eo%40googlegroups.com.

Reply via email to