Skip to content
GitLab
Explore
Sign in
Register
Primary navigation
Search or go to…
Project
dg
Manage
Activity
Members
Labels
Plan
Issues
Issue boards
Milestones
Wiki
Code
Merge requests
Repository
Branches
Commits
Tags
Repository graph
Compare revisions
Snippets
Build
Pipelines
Jobs
Pipeline schedules
Artifacts
Deploy
Releases
Package Registry
Model registry
Operate
Environments
Terraform modules
Monitor
Incidents
Analyze
Value stream analytics
Contributor analytics
CI/CD analytics
Repository analytics
Model experiments
Help
Help
Support
GitLab documentation
Compare GitLab plans
Community forum
Contribute to GitLab
Provide feedback
Terms and privacy
Keyboard shortcuts
?
Snippets
Groups
Projects
Show more breadcrumbs
gmsh
dg
Commits
29bed531
Commit
29bed531
authored
3 years ago
by
Matteo Cicuttin
Browse files
Options
Downloads
Patches
Plain Diff
USE_MPI cleanup, global DOF offsets on rank0.
parent
da525c94
No related branches found
No related tags found
No related merge requests found
Changes
4
Hide whitespace changes
Inline
Side-by-side
Showing
4 changed files
include/gmsh_io.h
+37
-11
37 additions, 11 deletions
include/gmsh_io.h
src/gmsh_io.cpp
+48
-14
48 additions, 14 deletions
src/gmsh_io.cpp
src/maxwell_solver.cpp
+12
-8
12 additions, 8 deletions
src/maxwell_solver.cpp
src/silo_io.cpp
+2
-2
2 additions, 2 deletions
src/silo_io.cpp
with
99 additions
and
35 deletions
include/gmsh_io.h
+
37
−
11
View file @
29bed531
...
...
@@ -45,9 +45,9 @@ struct boundary_descriptor
{
face_type
type
;
int
gmsh_entity
;
/* Actual GMSH entity */
#ifdef USE_MPI
int
parent_entity
;
/* Parent GMSH entity in partitioned model */
//auto operator <=> (const boundary_descriptor&) const = default;
#endif
/* USE_MPI */
bool
operator
<
(
const
boundary_descriptor
&
other
)
const
{
...
...
@@ -61,15 +61,19 @@ struct boundary_descriptor
}
boundary_descriptor
()
:
type
(
face_type
::
NONE
),
gmsh_entity
(
0
),
parent_entity
(
-
1
)
:
type
(
face_type
::
NONE
),
gmsh_entity
(
0
)
#ifdef USE_MPI
,
parent_entity
(
-
1
)
#endif
/* USE_MPI */
{}
int
material_tag
(
void
)
const
{
if
(
parent_entity
==
-
1
)
return
gmsh_entity
;
return
parent_entity
;
#ifdef USE_MPI
if
(
parent_entity
!=
-
1
)
return
parent_entity
;
#endif
return
gmsh_entity
;
}
};
...
...
@@ -78,17 +82,39 @@ struct boundary_descriptor
#define IMPLIES(a, b) ((not(a)) or (b))
#define IFF(a,b) (IMPLIES(a,b) and IMPLIES(b,a))
#ifdef USE_MPI
struct
dof_range
{
int
tag
;
int
dim
;
size_t
base
;
size_t
length
;
};
inline
std
::
ostream
&
operator
<<
(
std
::
ostream
&
os
,
const
dof_range
&
dr
)
{
os
<<
"Tag "
<<
dr
.
tag
<<
", dim "
<<
dr
.
dim
<<
": ("
;
os
<<
dr
.
base
<<
", "
<<
dr
.
length
<<
")"
;
return
os
;
}
#endif
class
model
{
int
geometric_order
;
int
approximation_order
;
#ifdef USE_MPI
int
num_partitions
;
#endif
/* USE_MPI */
/* Map from boundary tag to all the faces of the entity with that tag */
std
::
map
<
int
,
std
::
vector
<
element_key
>>
boundary_map
;
std
::
vector
<
entity
>
entities
;
std
::
vector
<
boundary_descriptor
>
bnd_descriptors
;
element_connectivity
<
element_key
>
conn
;
#ifdef USE_MPI
/* Map from boundary tag to all the partitions it belongs to */
std
::
map
<
int
,
std
::
vector
<
int
>>
comm_map
;
...
...
@@ -101,9 +127,9 @@ class model
/* Map from 3D entity tag to partition */
std
::
map
<
int
,
int
>
partition_inverse_map
;
std
::
vector
<
entity
>
entities
;
std
::
vector
<
boundary_descriptor
>
bnd_descriptor
s
;
element_connectivity
<
element_key
>
conn
;
/* Global dof range, only rank0 has this */
std
::
map
<
int
,
dof_range
>
global_dof_range
s
;
#endif
/* USE_MPI */
using
entofs_pair
=
std
::
pair
<
size_t
,
size_t
>
;
std
::
map
<
size_t
,
entofs_pair
>
etag_to_entity_offset
;
...
...
@@ -119,10 +145,10 @@ class model
void
populate_from_gmsh_rankN
(
void
);
void
import_gmsh_entities_rankN
(
void
);
int
my_partition
(
void
)
const
;
void
make_partition_to_entities_map
(
void
);
#endif
/* USE_MPI */
void
make_boundary_to_faces_map
(
void
);
void
make_partition_to_entities_map
(
void
);
public:
model
();
...
...
This diff is collapsed.
Click to expand it.
src/gmsh_io.cpp
+
48
−
14
View file @
29bed531
...
...
@@ -38,11 +38,17 @@ basis_grad_name(int order)
}
model
::
model
()
:
geometric_order
(
1
),
approximation_order
(
1
),
num_partitions
(
1
)
:
geometric_order
(
1
),
approximation_order
(
1
)
#ifdef USE_MPI
,
num_partitions
(
1
)
#endif
/* USE_MPI */
{}
model
::
model
(
int
pgo
,
int
pao
)
:
geometric_order
(
pgo
),
approximation_order
(
pao
),
num_partitions
(
1
)
:
geometric_order
(
pgo
),
approximation_order
(
pao
)
#ifdef USE_MPI
,
num_partitions
(
1
)
#endif
/* USE_MPI */
{
assert
(
geometric_order
>
0
);
assert
(
approximation_order
>
0
);
...
...
@@ -71,6 +77,7 @@ model::make_boundary_to_faces_map(void)
gm
::
getEntities
(
gmsh_entities
,
DIMENSION
(
2
));
for
(
auto
[
dim
,
tag
]
:
gmsh_entities
)
{
#ifdef USE_MPI
if
(
num_partitions
>
1
)
{
/* If the model is partitioned and entity has no parents,
...
...
@@ -90,6 +97,7 @@ model::make_boundary_to_faces_map(void)
if
(
pdim
==
dim
)
surface_to_parent_map
[
tag
]
=
ptag
;
}
#endif
/* USE_MPI */
std
::
vector
<
int
>
elemTypes
;
gmm
::
getElementTypes
(
elemTypes
,
dim
,
tag
);
...
...
@@ -103,19 +111,19 @@ model::make_boundary_to_faces_map(void)
}
}
for
(
auto
&
sp
:
surface_to_parent_map
)
std
::
cout
<<
sp
.
first
<<
" -> "
<<
sp
.
second
<<
std
::
endl
;
#ifdef USE_MPI
//
for (auto& sp : surface_to_parent_map)
//
std::cout << sp.first << " -> " << sp.second << std::endl;
assert
(
IFF
(
comm_map
.
size
()
==
0
,
num_partitions
==
1
)
);
#endif
/* USE_MPI */
}
#ifdef USE_MPI
void
model
::
make_partition_to_entities_map
()
{
#ifdef USE_MPI
ASSERT_MPI_RANK_0
;
#endif
/* USE_MPI */
gvp_t
gmsh_entities
;
...
...
@@ -144,6 +152,7 @@ model::make_partition_to_entities_map()
for
(
auto
&
[
pn
,
tags
]
:
partition_map
)
std
::
sort
(
tags
.
begin
(),
tags
.
end
());
/* For lookups */
}
#endif
/* USE_MPI */
void
model
::
update_connectivity
(
const
entity
&
e
,
size_t
entity_index
)
...
...
@@ -190,8 +199,13 @@ model::import_gmsh_entities_rank0(void)
size_t
dof_base
=
0
;
size_t
flux_base
=
0
;
size_t
index_base
=
0
;
#ifdef USE_MPI
size_t
global_dof_base
=
0
;
#endif
/* USE_MPI */
for
(
auto
[
dim
,
tag
]
:
gmsh_entities
)
{
#ifdef USE_MPI
if
(
num_partitions
>
1
)
{
/* If the model is partitioned and entity has no parents,
...
...
@@ -203,7 +217,7 @@ model::import_gmsh_entities_rank0(void)
if
(
ptag
==
-
1
)
continue
;
}
#endif
/* USE_MPI */
std
::
vector
<
int
>
elemTypes
;
gmm
::
getElementTypes
(
elemTypes
,
dim
,
tag
);
...
...
@@ -227,6 +241,12 @@ model::import_gmsh_entities_rank0(void)
if
(
entity_is_remote
)
{
dof_range
dr
;
dr
.
dim
=
dim
;
dr
.
tag
=
tag
;
dr
.
base
=
global_dof_base
;
dr
.
length
=
e
.
num_dofs
();
global_dof_ranges
[
tag
]
=
dr
;
global_dof_base
+=
e
.
num_dofs
();
assert
(
tag_partition
>
0
);
int
rank
=
tag_partition
-
1
;
e
.
mpi_send
(
rank
,
MPI_COMM_WORLD
);
...
...
@@ -252,13 +272,20 @@ model::import_gmsh_entities_rank0(void)
flux_base
+=
e
.
num_fluxes
();
index_base
+=
e
.
num_cells
();
#ifdef USE_MPI
dof_range
dr
;
dr
.
dim
=
dim
;
dr
.
tag
=
tag
;
dr
.
base
=
global_dof_base
;
dr
.
length
=
e
.
num_dofs
();
global_dof_ranges
[
tag
]
=
dr
;
global_dof_base
+=
e
.
num_dofs
();
#endif
/* USE_MPI */
update_connectivity
(
e
,
entity_index
);
entities
.
push_back
(
std
::
move
(
e
)
);
entity_index
++
;
}
}
std
::
cout
<<
dof_base
<<
" "
<<
flux_base
<<
std
::
endl
;
}
#ifdef USE_MPI
...
...
@@ -299,8 +326,6 @@ model::import_gmsh_entities_rankN(void)
entities
.
push_back
(
std
::
move
(
e
)
);
entity_index
++
;
}
std
::
cout
<<
dof_base
<<
" "
<<
flux_base
<<
std
::
endl
;
}
bool
...
...
@@ -370,19 +395,21 @@ model::map_boundaries(void)
}
/* and store the gmsh tag in the descriptor. */
bd
.
gmsh_entity
=
(
*
itor
).
second
;
#ifdef USE_MPI
if
(
num_partitions
>
1
)
{
auto
sitor
=
surface_to_parent_map
.
find
(
(
*
itor
).
second
);
if
(
sitor
!=
surface_to_parent_map
.
end
()
)
bd
.
parent_entity
=
surface_to_parent_map
.
at
(
(
*
itor
).
second
);
}
#endif
/* USE_MPI */
/* Finally, put the descriptor in the global array of faces. */
bnd_descriptors
.
at
(
fbase
+
iF
)
=
bd
;
}
fbase
+=
e
.
num_faces
();
}
#if 0
size_t normal = 0;
size_t interface = 0;
size_t boundary = 0;
...
...
@@ -402,6 +429,7 @@ model::map_boundaries(void)
std::cout << bfk.size() << " " << bnd_descriptors.size() << std::endl;
std::cout << normal << " " << interface << " " << boundary << " " << ipc_boundary << std::endl;
#endif
}
void
...
...
@@ -413,9 +441,10 @@ model::populate_from_gmsh_rank0()
gmm
::
setOrder
(
geometric_order
);
make_boundary_to_faces_map
();
make_partition_to_entities_map
();
#ifdef USE_MPI
make_partition_to_entities_map
();
/* At this point, the following members are populated:
* - geometric_order, approximation_order and num_partitions
* - boundary_map
...
...
@@ -432,6 +461,7 @@ model::populate_from_gmsh_rank0()
priv_MPI_Bcast
(
partition_map
,
0
,
MPI_COMM_WORLD
);
priv_MPI_Bcast
(
surface_to_parent_map
,
0
,
MPI_COMM_WORLD
);
/*
for (auto& [tag, vec] : partition_map)
{
std::cout << tag << ": " << vec.size() << " -> ";
...
...
@@ -439,6 +469,7 @@ model::populate_from_gmsh_rank0()
std::cout << v << " ";
std::cout << std::endl;
}
*/
#endif
/* USE_MPI */
...
...
@@ -460,6 +491,7 @@ model::populate_from_gmsh_rankN()
priv_MPI_Bcast
(
partition_map
,
0
,
MPI_COMM_WORLD
);
priv_MPI_Bcast
(
surface_to_parent_map
,
0
,
MPI_COMM_WORLD
);
/*
for (auto& [tag, vec] : partition_map)
{
std::cout << tag << ": " << vec.size() << " -> ";
...
...
@@ -467,7 +499,7 @@ model::populate_from_gmsh_rankN()
std::cout << v << " ";
std::cout << std::endl;
}
*/
entities
.
clear
();
import_gmsh_entities_rankN
();
map_boundaries
();
...
...
@@ -522,12 +554,14 @@ model::generate_mesh(double sf)
gmm
::
affineTransform
(
tr
);
}
#ifdef USE_MPI
void
model
::
partition_mesh
(
int
parts
)
{
gmm
::
partition
(
parts
);
num_partitions
=
parts
;
}
#endif
/* USE_MPI */
std
::
vector
<
entity
>::
const_iterator
model
::
begin
()
const
...
...
This diff is collapsed.
Click to expand it.
src/maxwell_solver.cpp
+
12
−
8
View file @
29bed531
...
...
@@ -220,10 +220,11 @@ void test_it(const model& mod, State& state, maxwell::parameter_loader& mpl)
timestep
(
state
,
mpl
,
ti
);
do_sources
(
mod
,
state
,
mpl
);
#ifdef USE_MPI
int
rank
;
MPI_Comm_rank
(
MPI_COMM_WORLD
,
&
rank
);
if
(
rank
==
0
)
#endif
/* USE_MPI */
if
(
(
silo_output_rate
!=
0
)
and
(
i
%
silo_output_rate
==
0
)
)
export_fields_to_silo
(
mod
,
state
,
mpl
);
...
...
@@ -309,14 +310,8 @@ int main(int argc, const char *argv[])
gmsh
::
open
(
mpl
.
sim_gmshmodel
()
);
model
mod
(
mpl
.
sim_geomorder
(),
mpl
.
sim_approxorder
()
);
/*
auto [scaled, sf] = mpl.mesh_scalefactor();
if (scaled)
mod.build(sf);
else
mod.build();
*/
#ifdef USE_MPI
if
(
rank
==
0
)
{
mod
.
generate_mesh
();
...
...
@@ -329,6 +324,15 @@ int main(int argc, const char *argv[])
{
mod
.
populate_from_gmsh
();
}
#else
/* USE_MPI */
auto
[
scaled
,
sf
]
=
mpl
.
mesh_scalefactor
();
if
(
scaled
)
mod
.
build
(
sf
);
else
mod
.
build
();
#endif
/* USE_MPI */
if
(
not
mpl
.
validate_model_params
(
mod
)
)
{
...
...
This diff is collapsed.
Click to expand it.
src/silo_io.cpp
+
2
−
2
View file @
29bed531
...
...
@@ -54,8 +54,8 @@ silo::gmsh_get_elements()
gmsh
::
model
::
getEntities
(
entities
,
3
);
for
(
auto
[
dim
,
tag
]
:
entities
)
{
if
(
tag
!=
4
and
tag
!=
5
)
continue
;
//
if (tag != 4 and tag != 5)
//
continue;
std
::
vector
<
int
>
elemTypes
;
gmsh
::
model
::
mesh
::
getElementTypes
(
elemTypes
,
dim
,
tag
);
...
...
This diff is collapsed.
Click to expand it.
Preview
0%
Loading
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Save comment
Cancel
Please
register
or
sign in
to comment