Dear Prof. Wolfganag,
                                       Thank you for the suggestion. I ran 
debugger with the ".prm" file and was successful to backtrace the error and 
resolved it. After this, using the step-19 particles tutorial refered by 
you, so far I am also successful in creating and exporting out the 
particles in ".vtu" file to visualize their ids and location (which is 
quadrature points location in my case).

As a next step, After creating particle on quadrature point, I am now 
trying to assign it the quadrature point value which is "*SymmetricTensor<2, 
dim>*" in my case. *I am using following simplified function similar to the 
one in particles tutorial file of step-19:*



































*template <int dim>    void  Problem<dim>::create_particles()    {      
next_unused_particle_id = 1;        FEValues<dim> fe_values (fe, 
quadrature_formula,                                     update_values | 
update_gradients |                                     
update_quadrature_points);        for (const auto &cell : 
dof_handler.active_cell_iterators())              {            
PointHistory<dim> *local_quadrature_points_history            
                   = reinterpret_cast<PointHistory<dim> 
*>(cell->user_pointer());                                      Assert 
(local_quadrature_points_history >=                                       
       &quadrature_point_history.front(),                                
                  ExcInternalError());                                   
Assert (local_quadrature_points_history <            
                               &quadrature_point_history.back(),            
                                      ExcInternalError());            
       fe_values.reinit(cell);                for (const unsigned int 
q_point :                        fe_values.quadrature_point_indices())    
              {                        const Point<dim> location =    
                            fe_values.quadrature_point(q_point);    
                    Particles::Particle<dim> new_particle;    
                    new_particle.set_location(location);    
                    new_particle.set_reference_location(    
                      mapping.transform_real_to_unit_cell(cell, 
location));                        
new_particle.set_id(next_unused_particle_id);                        
SymmetricTensor<2, dim> total_strain = 
local_quadrature_points_history[q_point].old_strain;*









*                        
new_particle.set_properties(make_array_view(total_strain));    
                    particle_handler.insert_particle(new_particle, 
cell);                        ++next_unused_particle_id;                  
}              }        particle_handler.update_cached_numbers();    }*

*But it gives the following error on running:*










*An error occurred in line <298> of file 
</home/muhammad/software/dealii-9.2.0/source/particles/particle.cc> in 
function    void dealii::Particles::Particle<dim, 
spacedim>::set_properties(const dealii::ArrayView<const double>&) [with int 
dim = 3; int spacedim = 3]The violated condition was:     property_pool != 
nullptrAdditional information:     This exception -- which is used in many 
places in the library -- usually indicates that some condition which the 
author of the code thought must be satisfied at a certain point in an 
algorithm, is not fulfilled. An example would be that the first part of an 
algorithm sorts elements of an array in ascending order, and a second part 
of the algorithm later encounters an element that is not larger than the 
previous one.There is usually not very much you can do if you encounter 
such an exception since it indicates an error in deal.II, not in your own 
program. Try to come up with the smallest possible program that still 
demonstrates the error and contact the deal.II mailing lists with it to 
obtain help.*I tried to fix it by initializing and setting the property 
pool to the new_particle in my current function :

*Particles::PropertyPool propertypool(6); // six entries present in the 
symmetric tensor of second order*
*new_particle.set_property_pool(propertypool);*

*But still it stops with the following error:*







*An error occurred in line <531> of file 
</home/muhammad/software/dealii-9.2.0/include/deal.II/base/array_view.h> in 
function    dealii::ArrayView<ElementType, MemorySpaceType>::value_type& 
dealii::ArrayView<ElementType, MemorySpaceType>::operator[](std::size_t) 
const [with ElementType = double; MemorySpaceType = 
dealii::MemorySpace::Host; dealii::ArrayView<ElementType, 
MemorySpaceType>::value_type = double; std::size_t = long unsigned int]The 
violated condition was:     static_cast<typename 
::dealii::internal::argument_type<void( typename 
std::common_type<decltype(i), decltype(n_elements)>::type)>:: type>(i) < 
static_cast<typename ::dealii::internal::argument_type<void( typename 
std::common_type<decltype(i), decltype(n_elements)>::type)>:: 
type>(n_elements)Additional information:     Index 3 is not in the 
half-open range [0,3).*
I explored the mailing list and found someone from the community already 
had similar problem of assigning the properties to the particles where it 
was suggested to use following way of assigning the properties:



*for( auto iter=particleHandler.begin(); iter!=particleHandler.end(); 
++iter, ++part ) *





* iter->set_properties 
(make_array_view(interpolatedParticleQuantities.begin()+ part*ncomponents, 
interpolatedParticleQuantities.begin()+ part*ncomponents + ncomponents)); *

If I have related it correctly to my case, I also then tried the similar 
way e.g.
   

*new_particle.set_properties(make_array_view(total_strain.begin_raw(),total_strain.end_raw()));*

instead of *new_particle.set_properties(make_array_view(total_strain));* 

It still stops with the same (2nd) error as mentioned above. Is there 
another way of doing it or am I not using the correct approach as per my 
limited knowledge so far? Sorry for the lengthy thread, I wanted to clear 
out the approaches which I tried so far as much as possible.

Thank you very much for your ongoing and productive support.

Regards,
Muhammad  
On Thursday, July 30, 2020 at 2:11:10 AM UTC+2 Wolfgang Bangerth wrote:

> On 7/28/20 11:05 AM, Muhammad Mashhood wrote:
> > 
> > /.~/working_dir/step-42$ /step-42 p1_chinese.prm
> > Segmentation fault (core dumped)/
> > 
> > I tried to comment out the line "  , particle_handler(triangulation, 
> mapping, 
> > /*n_properties=*/dim) " and by doing so the code starts running normally.
> > 
> > I have tried to run it in the debug mode but over there the program 
> exits with 
> > following error:
> > 
> > /(gdb) r
> > Starting program: ~/working_dir/step-42/step-42
> > [Thread debugging using libthread_db enabled]
> > Using host libthread_db library 
> "/lib/x86_64-linux-gnu/libthread_db.so.1".
> > *** Call this program as <./step-42 input.prm>
>
> You need to run the executable with
> r p1_chinese.prm
> to make sure it loads the appropriate input file, or you will not see the 
> segmentation fault in the debugger.
>
>
> > [Inferior 1 (process 5082) exited with code 01]/
> > 
> > Can there be a possible problem in using the " 
> > /particle_handler(triangulation, mapping, /*n_properties=*/dim) /" ? 
> Thank you!
>
> I don't know. Can you get a backtrace in the debugger?
>
> It's something that ought to work. If you can't figure out what is 
> happening, 
> send me the exact .cc file you're trying to run and I'll take a look!
>
> Best
> W.
>
> -- 
> ------------------------------------------------------------------------
> Wolfgang Bangerth email: bang...@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/80201761-2857-489e-9ac8-17f1471ac877n%40googlegroups.com.

Reply via email to