MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
face_nbr_geom.cpp
Go to the documentation of this file.
1// Copyright (c) 2010-2024, Lawrence Livermore National Security, LLC. Produced
2// at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3// LICENSE and NOTICE for details. LLNL-CODE-806117.
4//
5// This file is part of the MFEM library. For more information and source code
6// availability visit https://mfem.org.
7//
8// MFEM is free software; you can redistribute it and/or modify it under the
9// terms of the BSD-3 license. We welcome feedback and contributions, see file
10// CONTRIBUTING.md for details.
11
12#include "face_nbr_geom.hpp"
13#include "../general/forall.hpp"
14
15namespace mfem
16{
17
19 const GeometricFactors &geom_) : num_neighbor_elems(0), geom(geom_)
20{
21#ifdef MFEM_USE_MPI
22 if (const ParMesh *par_mesh = dynamic_cast<const ParMesh*>(geom.mesh))
23 {
24 const int flags = geom.computed_factors;
25 const int dim = par_mesh->Dimension();
26 const int sdim = par_mesh->SpaceDimension();
27
28 const_cast<ParMesh*>(par_mesh)->ExchangeFaceNbrData();
29
31
33 {
35 }
37 {
39 }
41 {
43 }
44
45 // Free memory of work arrays.
49 }
50#endif
51}
52
54 const Vector &x_local, Vector &x_shared, const int vdim)
55{
56#ifdef MFEM_USE_MPI
57
58 const int nq = geom.IntRule->Size();
59 const int ndof_per_el = vdim * nq;
60
61 const ParMesh *mesh = static_cast<const ParMesh*>(geom.mesh);
62
63 const int n_face_nbr = mesh->GetNFaceNeighbors();
64 if (n_face_nbr == 0) { return; }
65
66 const int ne_send = mesh->send_face_nbr_elements.Size_of_connections();
67
68 x_shared.SetSize(ndof_per_el * num_neighbor_elems);
69 send_offsets.SetSize(n_face_nbr + 1);
70 send_data.SetSize(ndof_per_el * ne_send);
71 auto h_send_data = Reshape(send_data.HostWrite(), ndof_per_el, ne_send);
72 const auto h_x_local = Reshape(x_local.HostRead(), ndof_per_el, mesh->GetNE());
73 int idx = 0;
74 Array<int> row;
75 for (int i = 0; i < n_face_nbr; ++i)
76 {
77 send_offsets[i] = ndof_per_el*idx;
78 mesh->send_face_nbr_elements.GetRow(i, row);
79 for (const int el : row)
80 {
81 for (int j = 0; j < ndof_per_el; ++j)
82 {
83 h_send_data(j, idx) = h_x_local(j, el);
84 }
85 ++idx;
86 }
87 }
88 send_offsets[n_face_nbr] = ndof_per_el*idx;
89 MFEM_ASSERT(send_offsets[n_face_nbr] == send_data.Size(), "");
90
91 recv_offsets.SetSize(n_face_nbr + 1);
92 for (int i = 0; i < n_face_nbr + 1; ++i)
93 {
94 recv_offsets[i] = ndof_per_el * mesh->face_nbr_elements_offset[i];
95 }
96 MFEM_ASSERT(recv_offsets[n_face_nbr] == x_shared.Size(), "");
97
98 MPI_Comm comm = mesh->GetComm();
99 std::vector<MPI_Request> send_reqs(n_face_nbr);
100 std::vector<MPI_Request> recv_reqs(n_face_nbr);
101 std::vector<MPI_Status> statuses(n_face_nbr);
102
103 bool mpi_gpu_aware = Device::GetGPUAwareMPI();
104 const auto send_data_ptr = mpi_gpu_aware ? send_data.Read() :
106 auto x_shared_ptr = mpi_gpu_aware ? x_shared.Write() : x_shared.HostWrite();
107
108 for (int i = 0; i < n_face_nbr; ++i)
109 {
110 const int nbr_rank = mesh->GetFaceNbrRank(i);
111 const int tag = 0;
112
113 MPI_Isend(send_data_ptr + send_offsets[i],
114 send_offsets[i+1] - send_offsets[i],
115 MPITypeMap<real_t>::mpi_type, nbr_rank, tag, comm, &send_reqs[i]);
116
117 MPI_Irecv(x_shared_ptr + recv_offsets[i],
118 recv_offsets[i+1] - recv_offsets[i],
119 MPITypeMap<real_t>::mpi_type, nbr_rank, tag, comm, &recv_reqs[i]);
120 }
121
122 MPI_Waitall(n_face_nbr, send_reqs.data(), statuses.data());
123 MPI_Waitall(n_face_nbr, recv_reqs.data(), statuses.data());
124#endif
125}
126
127} // namespace mfem
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
Definition array.hpp:697
int Size() const
Return the logical size of the array.
Definition array.hpp:144
void DeleteAll()
Delete the whole array.
Definition array.hpp:864
static bool GetGPUAwareMPI()
Definition device.hpp:291
const GeometricFactors & geom
The GeometricFactors of the Mesh.
FaceNeighborGeometricFactors(const GeometricFactors &geom_)
Communicate (if needed) to gather the face neighbor geometric factors.
Vector X
Physical coordinates of the mesh.
Vector detJ
Jacobian determinants.
void ExchangeFaceNbrQVectors(const Vector &x_local, Vector &x_shared, const int vdim)
Given a Q-vector x_local with vdim components, fill the face-neighbor Q-vector x_shared by communicat...
int num_neighbor_elems
Number of face neighbor elements.
Structure for storing mesh geometric factors: coordinates, Jacobians, and determinants of the Jacobia...
Definition mesh.hpp:2790
Vector X
Mapped (physical) coordinates of all quadrature points.
Definition mesh.hpp:2820
const Mesh * mesh
Definition mesh.hpp:2796
const IntegrationRule * IntRule
Definition mesh.hpp:2797
Vector detJ
Determinants of the Jacobians at all quadrature points.
Definition mesh.hpp:2835
Vector J
Jacobians of the element transformations at all quadrature points.
Definition mesh.hpp:2829
int GetNE() const
Returns number of elements.
Definition mesh.hpp:1226
Class for parallel meshes.
Definition pmesh.hpp:34
Table send_face_nbr_elements
Definition pmesh.hpp:436
MPI_Comm GetComm() const
Definition pmesh.hpp:402
int GetNFaceNeighborElements() const
Definition pmesh.hpp:518
Array< int > face_nbr_elements_offset
Definition pmesh.hpp:431
int GetNFaceNeighbors() const
Definition pmesh.hpp:517
int GetFaceNbrRank(int fn) const
Definition pmesh.cpp:2796
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
Definition table.cpp:187
int Size_of_connections() const
Definition table.hpp:98
Vector data type.
Definition vector.hpp:80
virtual const real_t * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:478
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:474
void Destroy()
Destroy a vector.
Definition vector.hpp:615
int Size() const
Returns the size of the vector.
Definition vector.hpp:218
void SetSize(int s)
Resize the vector to size s.
Definition vector.hpp:538
virtual real_t * HostWrite()
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), false).
Definition vector.hpp:486
virtual real_t * Write(bool on_dev=true)
Shortcut for mfem::Write(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:482
int dim
Definition ex24.cpp:53
MFEM_HOST_DEVICE DeviceTensor< sizeof...(Dims), T > Reshape(T *ptr, Dims... dims)
Wrap a pointer as a DeviceTensor with automatically deduced template parameters.
Definition dtensor.hpp:131
Helper struct to convert a C++ type to an MPI type.