Dear Jean-Paul,

I'm trying to solve a simpler problem of uniaxial response of a plate with 
hole using FE_Nothing as you recommended. I re-read steps 10, 27, and 46 to 
refresh my memory. I think I have the logic but I appreciate if you check 
some of the steps I have a little concern about. The two main step I took 
are:
1) I assign the cells inside the hole as FE_Nothing and the rest as regular 
FE.
2) During the stiffness matrix and residual force formation, I only go over 
the elements in regular FE part and exclude the FE_Nothing element from the 
calculation.

My only concern are the cells on the interface. Does performing assembly 
and residual formation just on the FE cells (and exclude FE_Nothing cells) 
guarantee the traction free surface of the cells on the surface of the 
hole?  Or I need to take extra action to guarantee the traction free 
surface of the hole?

Thanks
Reza
On Tuesday, April 20, 2021 at 5:19:13 PM UTC-4 Jean-Paul Pelteret wrote:

> Hi Reza,
>
> With no context as to how element deletion affects the geometry itself, 
> may I ask if you’ve considered element deactivation using the FE_Nothing 
> element? If you can afford to store the mesh but primarily want to ensure 
> that no DoFs are assigned to a region outside of the physical domain 
> boundaries, then this might be an alternate approach to adjusting the mesh 
> directly.
>
> Best,
> Jean-Paul
>
> On 20. Apr 2021, at 21:58, mohammadreza yaghoobi <mohammadre...@gmail.com> 
> wrote:
>
> Thanks Marc. 
>
> My biggest challenges are both time and memory. I'm trying to model a 
> crystal plasticity sample with ~20 million elements. In the case of no 
> element deletion, I can simply use global refinement and I won't face 
> memory issues while generating mesh and running the code (and of course no 
> time on element deletion).
>
> However, when I'm trying to add the element deletion feature, then I 
> cannot use global refinement,  and I'll face the memory issue, since all of 
> the triangulation will be on the coarsest level. On top of that, for every 
> processor, the code has to search the element deletion pattern for all of 
> the triangulation (~20 million elements). As I increase the sample size, 
> each processor still needs to handle all the elements which becomes a big 
> bottle-neck, and also uses the resources very inefficiently. Actually, 
> parallelization cannot help at all with this part of the code, which 
> actually serves as the serial with large memory requirements. Just to give 
> you an idea, it will take about 12 hrs in my example to bypass the element 
> deletion search part running on 1500 processors, each has 38gb memory. I 
> actually have to underutilize my node (instead of 24, I'm using 5) because 
> of memory issues.
>
> The ideal case would be each processor just handling their own part and 
> some reduction happened throughout the processors.
>
> I hope I was able to clarify the bottle neck.
>
> Thanks
> Reza
>
>
>
> On Tue, Apr 20, 2021 at 3:47 PM Marc Fehling <mafe...@gmail.com> wrote:
>
>> Hello Reza,
>>
>> `parallel::distributed::Triangulation` objects store the entire coarse 
>> mesh, see also the module documentation 
>> <https://www.dealii.org/current/doxygen/deal.II/group__distributed.html>. 
>> This allows for a quick execution of adaptive refinement. So from my 
>> understanding, you need to perform this operation on *all* active cells, 
>> including ghost and artificial ones. So the working code you have might be 
>> the only way to use 
>> `GridGenerator::create_triangulation_with_removed_cells()` in the parallel 
>> distributed context.
>>
>> Can you quantify how inefficient your code is? For example did you time 
>> how long your mesh creation takes and compare it to the overall runtime?
>>
>> Best,
>> Marc
>>
>> On Tuesday, April 20, 2021 at 12:14:10 PM UTC-6 mohammadre...@gmail.com 
>> wrote:
>>
>>> Hi,
>>>
>>> I'm trying to subtract some elements based on a user-defined criteria. 
>>> Ideally, I want to use each processor in a way that they are responsible 
>>> only for their partition (using is_locally_owned()) to remove the cells 
>>> from the triangulation 
>>> using GridGenerator::create_triangulation_with_removed_cells. However, when 
>>> I used  this, it will not proceed after 
>>> GridGenerator::create_triangulation_with_removed_cells. 
>>> Very inefficient solution is each processor will process the whole 
>>> triangulation space and remove the elements. This works. However, as the 
>>> simulation size increases, this becomes really a bottle neck.
>>>
>>> Is there any way I can tackle this smarter?
>>>
>>> Thanks
>>> Reza
>>>
>>> *Below is the version of the code which works but very inefficiently.*
>>> *I added the red lines so every processor be responsible just for their 
>>> own parts. It stops at 
>>> GridGenerator::create_triangulation_with_removed_cells.*
>>>
>>>
>>>
>>> GridGenerator::subdivided_hyper_rectangle (triangulation2, 
>>> userInputs.subdivisions, Point<dim>(), 
>>> Point<dim>(userInputs.span[0],userInputs.span[1],userInputs.span[2]), true);
>>>
>>> typename Triangulation< dim, dim >::active_cell_iterator cell = 
>>> triangulation2.begin_active(), endc = triangulation2.end();
>>>
>>> for (; cell!=endc; ++cell) {
>>> if (cell->is_locally_owned()){
>>> double pnt3[3];
>>> const Point<dim> pnt2=cell->center();
>>> for (unsigned int i=0; i<dim; ++i){
>>> pnt3[i]=pnt2[i];
>>> }
>>>
>>> gID=orientations_Mesh.getMaterialID(pnt3);
>>> cellOrientationMap_Mesh.push_back(gID);
>>> }
>>> }
>>>
>>> unsigned int materialID;
>>> unsigned int cell2=0;
>>> cell = triangulation2.begin_active(); endc = triangulation2.end();
>>> std::set<typename Triangulation< dim, dim >::active_cell_iterator> 
>>> cells_to_remove;
>>> unsigned int NumberOwned=0;
>>> for (; cell!=endc; ++cell) {
>>> if (cell->is_locally_owned()){
>>> materialID=cellOrientationMap_Mesh[cell2];
>>> NumberOwned=NumberOwned+1;
>>> if (materialID ==userInputs.deletionGrainID){
>>> cells_to_remove.insert (cell);
>>> }
>>> cell2=cell2+1;
>>> }
>>> }
>>>
>>> char buffer[200];
>>>
>>> if (cells_to_remove.size()>0){
>>>
>>> GridGenerator::create_triangulation_with_removed_cells(triangulation2,cells_to_remove,triangulation);
>>> }
>>> else{
>>> triangulation.copy_triangulation(triangulation2);
>>> }
>>> }
>>> else {
>>> triangulation.copy_triangulation(triangulation2);
>>> }
>>>
>>
>> -- 
>> 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/07d82b49-8eae-4eed-84b8-bbbb1f534baen%40googlegroups.com
>>  
>> <https://groups.google.com/d/msgid/dealii/07d82b49-8eae-4eed-84b8-bbbb1f534baen%40googlegroups.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+un...@googlegroups.com.
>
> To view this discussion on the web visit 
> https://groups.google.com/d/msgid/dealii/CAGxJ-Spk29Y2KNP_9W7%3DC0YkHrPpqLXiR56V7qiTUeCPLAu8DQ%40mail.gmail.com
>  
> <https://groups.google.com/d/msgid/dealii/CAGxJ-Spk29Y2KNP_9W7%3DC0YkHrPpqLXiR56V7qiTUeCPLAu8DQ%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/e496a2a0-ffa3-4038-975e-903401fec486n%40googlegroups.com.

Reply via email to