Skip to content

Commit 148d51e

Browse files
Fix full parallel assembly for AMS solver and wave ports
1 parent f9f8872 commit 148d51e

File tree

2 files changed

+84
-70
lines changed

2 files changed

+84
-70
lines changed

palace/linalg/ams.cpp

+18-18
Original file line numberDiff line numberDiff line change
@@ -51,12 +51,13 @@ void HypreAmsSolver::ConstructAuxiliaryMatrices(mfem::ParFiniteElementSpace &nd_
5151
// HypreAMS:Init. Start with the discrete gradient matrix.
5252
{
5353
// XX TODO: Partial assembly option?
54-
auto grad = std::make_unique<mfem::DiscreteLinearOperator>(&h1_fespace, &nd_fespace);
55-
grad->AddDomainInterpolator(new mfem::GradientInterpolator);
56-
grad->SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
57-
grad->Assemble();
58-
grad->Finalize();
59-
ParOperator RAP_G(std::move(grad), h1_fespace, nd_fespace, true);
54+
mfem::DiscreteLinearOperator grad(&h1_fespace, &nd_fespace);
55+
grad.AddDomainInterpolator(new mfem::GradientInterpolator);
56+
grad.SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
57+
grad.Assemble();
58+
grad.Finalize();
59+
ParOperator RAP_G(std::unique_ptr<mfem::SparseMatrix>(grad.LoseMat()), h1_fespace,
60+
nd_fespace, true);
6061
G = RAP_G.StealParallelAssemble();
6162
}
6263

@@ -112,18 +113,17 @@ void HypreAmsSolver::ConstructAuxiliaryMatrices(mfem::ParFiniteElementSpace &nd_
112113
}
113114
else
114115
{
115-
{
116-
// XX TODO: Partial assembly option?
117-
mfem::ParFiniteElementSpace h1d_fespace(&mesh, h1_fespace.FEColl(), space_dim,
118-
mfem::Ordering::byVDIM);
119-
auto pi = std::make_unique<mfem::DiscreteLinearOperator>(&h1d_fespace, &nd_fespace);
120-
pi->AddDomainInterpolator(new mfem::IdentityInterpolator);
121-
pi->SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
122-
pi->Assemble();
123-
pi->Finalize();
124-
ParOperator RAP_Pi(std::move(pi), h1d_fespace, nd_fespace, true);
125-
Pi = RAP_Pi.StealParallelAssemble();
126-
}
116+
// XX TODO: Partial assembly option?
117+
mfem::ParFiniteElementSpace h1d_fespace(&mesh, h1_fespace.FEColl(), space_dim,
118+
mfem::Ordering::byVDIM);
119+
mfem::DiscreteLinearOperator pi(&h1d_fespace, &nd_fespace);
120+
pi.AddDomainInterpolator(new mfem::IdentityInterpolator);
121+
pi.SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
122+
pi.Assemble();
123+
pi.Finalize();
124+
ParOperator RAP_Pi(std::unique_ptr<mfem::SparseMatrix>(pi.LoseMat()), h1d_fespace,
125+
nd_fespace, true);
126+
Pi = RAP_Pi.StealParallelAssemble();
127127
if (cycle_type >= 10)
128128
{
129129
// Get blocks of Pi corresponding to each component, and free Pi.

palace/models/waveportoperator.cpp

+66-52
Original file line numberDiff line numberDiff line change
@@ -84,12 +84,13 @@ std::unique_ptr<ParOperator> GetBtt(const MaterialOperator &mat_op,
8484
constexpr MaterialPropertyType MatType = MaterialPropertyType::INV_PERMEABILITY;
8585
constexpr MeshElementType ElemType = MeshElementType::BDR_SUBMESH;
8686
MaterialPropertyCoefficient<MatType, ElemType> muinv_func(mat_op);
87-
auto btt = std::make_unique<mfem::SymmetricBilinearForm>(&nd_fespace);
88-
btt->AddDomainIntegrator(new mfem::VectorFEMassIntegrator(muinv_func));
89-
btt->SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
90-
btt->Assemble(skip_zeros);
91-
btt->Finalize(skip_zeros);
92-
return std::make_unique<ParOperator>(std::move(btt), nd_fespace);
87+
mfem::SymmetricBilinearForm btt(&nd_fespace);
88+
btt.AddDomainIntegrator(new mfem::VectorFEMassIntegrator(muinv_func));
89+
btt.SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
90+
btt.Assemble(skip_zeros);
91+
btt.Finalize(skip_zeros);
92+
return std::make_unique<ParOperator>(std::unique_ptr<mfem::SparseMatrix>(btt.LoseMat()),
93+
nd_fespace);
9394
}
9495

9596
std::unique_ptr<ParOperator> GetBtn(const MaterialOperator &mat_op,
@@ -100,12 +101,13 @@ std::unique_ptr<ParOperator> GetBtn(const MaterialOperator &mat_op,
100101
constexpr MaterialPropertyType MatType = MaterialPropertyType::INV_PERMEABILITY;
101102
constexpr MeshElementType ElemType = MeshElementType::BDR_SUBMESH;
102103
MaterialPropertyCoefficient<MatType, ElemType> muinv_func(mat_op);
103-
auto btn = std::make_unique<mfem::MixedBilinearForm>(&h1_fespace, &nd_fespace);
104-
btn->AddDomainIntegrator(new mfem::MixedVectorGradientIntegrator(muinv_func));
105-
btn->SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
106-
btn->Assemble(skip_zeros);
107-
btn->Finalize(skip_zeros);
108-
return std::make_unique<ParOperator>(std::move(btn), h1_fespace, nd_fespace, false);
104+
mfem::MixedBilinearForm btn(&h1_fespace, &nd_fespace);
105+
btn.AddDomainIntegrator(new mfem::MixedVectorGradientIntegrator(muinv_func));
106+
btn.SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
107+
btn.Assemble(skip_zeros);
108+
btn.Finalize(skip_zeros);
109+
return std::make_unique<ParOperator>(std::unique_ptr<mfem::SparseMatrix>(btn.LoseMat()),
110+
h1_fespace, nd_fespace, false);
109111
}
110112

111113
std::array<std::unique_ptr<ParOperator>, 3> GetBnn(const MaterialOperator &mat_op,
@@ -115,38 +117,44 @@ std::array<std::unique_ptr<ParOperator>, 3> GetBnn(const MaterialOperator &mat_o
115117
constexpr MaterialPropertyType MatTypeMuInv = MaterialPropertyType::INV_PERMEABILITY;
116118
constexpr MeshElementType ElemType = MeshElementType::BDR_SUBMESH;
117119
MaterialPropertyCoefficient<MatTypeMuInv, ElemType> muinv_func(mat_op);
118-
auto bnn1 = std::make_unique<mfem::SymmetricBilinearForm>(&h1_fespace);
119-
bnn1->AddDomainIntegrator(new mfem::DiffusionIntegrator(muinv_func));
120-
bnn1->SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
121-
bnn1->Assemble(skip_zeros);
122-
bnn1->Finalize(skip_zeros);
120+
mfem::SymmetricBilinearForm bnn1(&h1_fespace);
121+
bnn1.AddDomainIntegrator(new mfem::DiffusionIntegrator(muinv_func));
122+
bnn1.SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
123+
bnn1.Assemble(skip_zeros);
124+
bnn1.Finalize(skip_zeros);
123125

124126
constexpr MaterialPropertyType MatTypeEpsReal = MaterialPropertyType::PERMITTIVITY_REAL;
125127
NormalProjectedCoefficient epsilon_func(
126128
std::make_unique<MaterialPropertyCoefficient<MatTypeEpsReal, ElemType>>(mat_op));
127-
auto bnn2r = std::make_unique<mfem::SymmetricBilinearForm>(&h1_fespace);
128-
bnn2r->AddDomainIntegrator(new mfem::MixedScalarMassIntegrator(epsilon_func));
129-
bnn2r->SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
130-
bnn2r->Assemble(skip_zeros);
131-
bnn2r->Finalize(skip_zeros);
129+
mfem::SymmetricBilinearForm bnn2r(&h1_fespace);
130+
bnn2r.AddDomainIntegrator(new mfem::MixedScalarMassIntegrator(epsilon_func));
131+
bnn2r.SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
132+
bnn2r.Assemble(skip_zeros);
133+
bnn2r.Finalize(skip_zeros);
132134

133135
// Contribution for loss tangent: ε -> ε * (1 - i tan(δ)).
134136
if (!mat_op.HasLossTangent())
135137
{
136-
return {std::make_unique<ParOperator>(std::move(bnn1), h1_fespace),
137-
std::make_unique<ParOperator>(std::move(bnn2r), h1_fespace), nullptr};
138+
return {std::make_unique<ParOperator>(
139+
std::unique_ptr<mfem::SparseMatrix>(bnn1.LoseMat()), h1_fespace),
140+
std::make_unique<ParOperator>(
141+
std::unique_ptr<mfem::SparseMatrix>(bnn2r.LoseMat()), h1_fespace),
142+
nullptr};
138143
}
139144
constexpr MaterialPropertyType MatTypeEpsImag = MaterialPropertyType::PERMITTIVITY_IMAG;
140145
NormalProjectedCoefficient negepstandelta_func(
141146
std::make_unique<MaterialPropertyCoefficient<MatTypeEpsImag, ElemType>>(mat_op));
142-
auto bnn2i = std::make_unique<mfem::SymmetricBilinearForm>(&h1_fespace);
143-
bnn2i->AddDomainIntegrator(new mfem::MixedScalarMassIntegrator(negepstandelta_func));
144-
bnn2i->SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
145-
bnn2i->Assemble(skip_zeros);
146-
bnn2i->Finalize(skip_zeros);
147-
return {std::make_unique<ParOperator>(std::move(bnn1), h1_fespace),
148-
std::make_unique<ParOperator>(std::move(bnn2r), h1_fespace),
149-
std::make_unique<ParOperator>(std::move(bnn2i), h1_fespace)};
147+
mfem::SymmetricBilinearForm bnn2i(&h1_fespace);
148+
bnn2i.AddDomainIntegrator(new mfem::MixedScalarMassIntegrator(negepstandelta_func));
149+
bnn2i.SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
150+
bnn2i.Assemble(skip_zeros);
151+
bnn2i.Finalize(skip_zeros);
152+
return {std::make_unique<ParOperator>(std::unique_ptr<mfem::SparseMatrix>(bnn1.LoseMat()),
153+
h1_fespace),
154+
std::make_unique<ParOperator>(
155+
std::unique_ptr<mfem::SparseMatrix>(bnn2r.LoseMat()), h1_fespace),
156+
std::make_unique<ParOperator>(
157+
std::unique_ptr<mfem::SparseMatrix>(bnn2i.LoseMat()), h1_fespace)};
150158
}
151159

152160
std::array<std::unique_ptr<ParOperator>, 3> GetAtt(const MaterialOperator &mat_op,
@@ -157,36 +165,42 @@ std::array<std::unique_ptr<ParOperator>, 3> GetAtt(const MaterialOperator &mat_o
157165
constexpr MeshElementType ElemType = MeshElementType::BDR_SUBMESH;
158166
NormalProjectedCoefficient muinv_func(
159167
std::make_unique<MaterialPropertyCoefficient<MatTypeMuInv, ElemType>>(mat_op));
160-
auto att1 = std::make_unique<mfem::SymmetricBilinearForm>(&nd_fespace);
161-
att1->AddDomainIntegrator(new mfem::CurlCurlIntegrator(muinv_func));
162-
att1->SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
163-
att1->Assemble(skip_zeros);
164-
att1->Finalize(skip_zeros);
168+
mfem::SymmetricBilinearForm att1(&nd_fespace);
169+
att1.AddDomainIntegrator(new mfem::CurlCurlIntegrator(muinv_func));
170+
att1.SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
171+
att1.Assemble(skip_zeros);
172+
att1.Finalize(skip_zeros);
165173

166174
constexpr MaterialPropertyType MatTypeEpsReal = MaterialPropertyType::PERMITTIVITY_REAL;
167175
MaterialPropertyCoefficient<MatTypeEpsReal, ElemType> epsilon_func(mat_op);
168-
auto att2r = std::make_unique<mfem::SymmetricBilinearForm>(&nd_fespace);
169-
att2r->AddDomainIntegrator(new mfem::VectorFEMassIntegrator(epsilon_func));
170-
att2r->SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
171-
att2r->Assemble(skip_zeros);
172-
att2r->Finalize(skip_zeros);
176+
mfem::SymmetricBilinearForm att2r(&nd_fespace);
177+
att2r.AddDomainIntegrator(new mfem::VectorFEMassIntegrator(epsilon_func));
178+
att2r.SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
179+
att2r.Assemble(skip_zeros);
180+
att2r.Finalize(skip_zeros);
173181

174182
// Contribution for loss tangent: ε -> ε * (1 - i tan(δ)).
175183
if (!mat_op.HasLossTangent())
176184
{
177-
return {std::make_unique<ParOperator>(std::move(att1), nd_fespace),
178-
std::make_unique<ParOperator>(std::move(att2r), nd_fespace), nullptr};
185+
return {std::make_unique<ParOperator>(
186+
std::unique_ptr<mfem::SparseMatrix>(att1.LoseMat()), nd_fespace),
187+
std::make_unique<ParOperator>(
188+
std::unique_ptr<mfem::SparseMatrix>(att2r.LoseMat()), nd_fespace),
189+
nullptr};
179190
}
180191
constexpr MaterialPropertyType MatTypeEpsImag = MaterialPropertyType::PERMITTIVITY_IMAG;
181192
MaterialPropertyCoefficient<MatTypeEpsImag, ElemType> negepstandelta_func(mat_op);
182-
auto att2i = std::make_unique<mfem::SymmetricBilinearForm>(&nd_fespace);
183-
att2i->AddDomainIntegrator(new mfem::VectorFEMassIntegrator(negepstandelta_func));
184-
att2i->SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
185-
att2i->Assemble(skip_zeros);
186-
att2i->Finalize(skip_zeros);
187-
return {std::make_unique<ParOperator>(std::move(att1), nd_fespace),
188-
std::make_unique<ParOperator>(std::move(att2r), nd_fespace),
189-
std::make_unique<ParOperator>(std::move(att2i), nd_fespace)};
193+
mfem::SymmetricBilinearForm att2i(&nd_fespace);
194+
att2i.AddDomainIntegrator(new mfem::VectorFEMassIntegrator(negepstandelta_func));
195+
att2i.SetAssemblyLevel(mfem::AssemblyLevel::LEGACY);
196+
att2i.Assemble(skip_zeros);
197+
att2i.Finalize(skip_zeros);
198+
return {std::make_unique<ParOperator>(std::unique_ptr<mfem::SparseMatrix>(att1.LoseMat()),
199+
nd_fespace),
200+
std::make_unique<ParOperator>(
201+
std::unique_ptr<mfem::SparseMatrix>(att2r.LoseMat()), nd_fespace),
202+
std::make_unique<ParOperator>(
203+
std::unique_ptr<mfem::SparseMatrix>(att2i.LoseMat()), nd_fespace)};
190204
}
191205

192206
std::array<std::unique_ptr<mfem::HypreParMatrix>, 6>

0 commit comments

Comments
 (0)