Dear Prof. Bangerth,
Thanks a lot for the explanation. I have used the setter function approach
to solve this problem but I think now I understand how the other approach
will work as well.
Regards,
Sabyasachi
On Friday, September 17, 2021 at 6:05:07 PM UTC-4 Wolfgang Bangerth wrote:
>
> Sabyasachi,
> it took me a while to understand that what you are asking is not actually
> a
> question about deal.II but about how to do what you are suggesting in C++.
>
> The relevant part of the setup_q_point_data() function in step-18 is this:
>
> quadrature_point_history.resize(
> triangulation.n_locally_owned_active_cells()
> * quadrature_formula.size());
>
> unsigned int history_index = 0;
> for (auto &cell : triangulation.active_cell_iterators())
> if (cell->is_locally_owned())
> {
> cell->set_user_pointer(&quadrature_point_history[history_index]);
> history_index += quadrature_formula.size();
> }
>
> What this is doing is simply allocating enough objects and then setting
> pointers from the cell to the individual objects. In other words, at the
> time
> where the objects are created (in the resize() call), the history objects
> are
> not yet associated with individual cells, and so you could not set a depth
> for
> each object.
>
> You could of course turn the order around and go over cells and on every
> cell
> you create an object that you then add to the collection via
> quadrature_point_history.push_back(...);
> Because you know which cell/q-point each object is associated when you
> create
> it, you can compute the depth and pass it to the constructor.
>
> Alternatively, you could leave the order in place and set the depth at a
> later
> time -- i.e., you'd add a "setter function" to your class, and in the body
> of
> the loop over all cells you'd loop over all quadrature points and set the
> depth for the corresponding history object.
>
> I hope this helps!
> Best
> W.
>
>
> On 9/8/21 2:04 PM, sabyasachi chatterjee wrote:
> > Hello all,
> >
> > I defined a class named SCDWrapper inside the point history class which
> needs
> > information about the depth of the gauss point corresponding to the
> point
> > history object's location from the surface when it is created.
> >
> > ```
> > template <int dim>
> > class PointHistory
> > {
> > public:
> > PointHistory(double&);
> > ~PointHistory();
> > SCDWrapper srscd;
> > };
> > ```
> >
> > So I added an argument storing the depth to the constructor of point
> history
> > class which then passes it to the member class.
> >
> > ```
> > template <int dim> PointHistory<dim>::PointHistory(double
> > &depth_in):srscd(depth_in)
> > {
> > }
> > ```
> >
> > My question is how do I initialize the point history objects with this
> > argument (which is variable across the gauss points) and pass it to the
> > constructor. For example, in Step 18 in the function
> > /setup_quadrature_point_history/ where the point history objects are
> > initialized, how do I pass the information about this depth variable.
> >
> > Thanks,
> > Sabyasachi
> >
> >
> >
> >
> > --
> > The deal.II project is located at http://www.dealii.org/
> > <
> https://nam10.safelinks.protection.outlook.com/?url=http%3A%2F%2Fwww.dealii.org%2F&data=04%7C01%7CWolfgang.Bangerth%40colostate.edu%7Cf76c21631c3a4524f70e08d97303efe6%7Cafb58802ff7a4bb1ab21367ff2ecfc8b%7C0%7C0%7C637667283758597366%7CUnknown%7CTWFpbGZsb3d8eyJWIjoiMC4wLjAwMDAiLCJQIjoiV2luMzIiLCJBTiI6Ik1haWwiLCJXVCI6Mn0%3D%7C3000&sdata=8pUEg2fir6ovEFyj%2BXUBO2XKLPi8pC9%2Fkn3%2BYSK8cVk%3D&reserved=0
> >
> > For mailing list/forum options, see
> > https://groups.google.com/d/forum/dealii?hl=en
> > <
> https://nam10.safelinks.protection.outlook.com/?url=https%3A%2F%2Fgroups.google.com%2Fd%2Fforum%2Fdealii%3Fhl%3Den&data=04%7C01%7CWolfgang.Bangerth%40colostate.edu%7Cf76c21631c3a4524f70e08d97303efe6%7Cafb58802ff7a4bb1ab21367ff2ecfc8b%7C0%7C0%7C637667283758607353%7CUnknown%7CTWFpbGZsb3d8eyJWIjoiMC4wLjAwMDAiLCJQIjoiV2luMzIiLCJBTiI6Ik1haWwiLCJXVCI6Mn0%3D%7C3000&sdata=J6zhwNNMmHk3%2F%2BMM%2FKIqP0PIEZSoda2zDHqTP2CAbC8%3D&reserved=0
> >
> > ---
> > 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 [email protected]
> > <mailto:[email protected]>.
> > To view this discussion on the web visit
> >
> https://groups.google.com/d/msgid/dealii/14ce9ef5-57d2-4bbd-bc07-52e5198a5490n%40googlegroups.com
>
> > <
> https://nam10.safelinks.protection.outlook.com/?url=https%3A%2F%2Fgroups.google.com%2Fd%2Fmsgid%2Fdealii%2F14ce9ef5-57d2-4bbd-bc07-52e5198a5490n%2540googlegroups.com%3Futm_medium%3Demail%26utm_source%3Dfooter&data=04%7C01%7CWolfgang.Bangerth%40colostate.edu%7Cf76c21631c3a4524f70e08d97303efe6%7Cafb58802ff7a4bb1ab21367ff2ecfc8b%7C0%7C0%7C637667283758617352%7CUnknown%7CTWFpbGZsb3d8eyJWIjoiMC4wLjAwMDAiLCJQIjoiV2luMzIiLCJBTiI6Ik1haWwiLCJXVCI6Mn0%3D%7C3000&sdata=SqhUO%2FsVvsLvv8R%2BCxUDNSfgGyPj4GLpUqW9qeJf3Lw%3D&reserved=0
> >.
>
>
> --
> ------------------------------------------------------------------------
> Wolfgang Bangerth email: [email protected]
> 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 [email protected].
To view this discussion on the web visit
https://groups.google.com/d/msgid/dealii/462b03fe-cc2a-45a5-af2f-82e058778e5an%40googlegroups.com.