Skip to content

Commit a720e33

Browse files
committed
Eliminate some dead _GPU_ code.
Specifically, parts of GradNonLinearForm and GradFaceIntegrator had _GPU_ specific code that was superseded by device-specific code in the Gradients class and thus was never called. This code has been eliminated to improve hygiene in preparation for changes to address #198.
1 parent f02d085 commit a720e33

File tree

3 files changed

+12
-103
lines changed

3 files changed

+12
-103
lines changed

src/faceGradientIntegration.cpp

+7-24
Original file line numberDiff line numberDiff line change
@@ -60,18 +60,17 @@ void GradFaceIntegrator::AssembleFaceVector(const FiniteElement &el1, const Fini
6060
elvect.SetSize((dof1 + dof2) * num_equation * dim);
6161
elvect = 0.0;
6262

63-
// #ifdef _GPU_
64-
// DenseMatrix elfun1_mat(elfun.GetData(), dof1, num_equation);
65-
// DenseMatrix elfun2_mat(elfun.GetData()+dof1*num_equation,
66-
// dof2, num_equation);
67-
#ifndef _GPU_
68-
DenseMatrix elfun1_mat(elfun.GetData(), dof1, num_equation * dim);
69-
DenseMatrix elfun2_mat(elfun.GetData() + dof1 * num_equation * dim, dof2, num_equation * dim);
63+
// NB: Incoming Vector &elfun is in gradient space. The first
64+
// component is the state, and the rest is zero. See
65+
// GradNonLinearForm::Mult in gradNonLinearForm.cpp for details.
66+
DenseMatrix elfun1_mat(elfun.GetData(), dof1, num_equation);
67+
DenseMatrix elfun2_mat(elfun.GetData() + dof1 * num_equation * dim, dof2, num_equation);
68+
69+
7070
DenseMatrix elvect1_mat(elvect.GetData(), dof1, num_equation * dim);
7171
DenseMatrix elvect2_mat(elvect.GetData() + dof1 * num_equation * dim, dof2, num_equation * dim);
7272
elvect1_mat = 0.;
7373
elvect2_mat = 0.;
74-
#endif
7574

7675
// Integration order calculation from DGTraceIntegrator
7776
int intorder;
@@ -104,20 +103,12 @@ void GradFaceIntegrator::AssembleFaceVector(const FiniteElement &el1, const Fini
104103
for (int eq = 0; eq < num_equation; eq++) {
105104
double sum = 0.;
106105
for (int k = 0; k < dof1; k++) {
107-
#ifdef _GPU_
108-
sum += shape1(k) * elfun(k + eq * dof1);
109-
#else
110106
sum += shape1[k] * elfun1_mat(k, eq);
111-
#endif
112107
}
113108
iUp1(eq) = sum;
114109
sum = 0.;
115110
for (int k = 0; k < dof2; k++) {
116-
#ifdef _GPU_
117-
sum += shape2(k) * elfun(k + eq * dof2 + dof1 * num_equation);
118-
#else
119111
sum += shape2[k] * elfun2_mat(k, eq);
120-
#endif
121112
}
122113
iUp2(eq) = sum;
123114
mean(eq) = (iUp1(eq) + iUp2(eq)) / 2.;
@@ -133,18 +124,10 @@ void GradFaceIntegrator::AssembleFaceVector(const FiniteElement &el1, const Fini
133124
const double du1n = (mean(eq) - iUp1(eq)) * nor(d);
134125
const double du2n = (mean(eq) - iUp2(eq)) * nor(d);
135126
for (int k = 0; k < dof1; k++) {
136-
#ifdef _GPU_
137-
elvect(k + eq + d * num_equation) += du1n * shape1(k);
138-
#else
139127
elvect1_mat(k, eq + d * num_equation) += du1n * shape1(k);
140-
#endif
141128
}
142129
for (int k = 0; k < dof2; k++) {
143-
#ifdef _GPU_
144-
elvect(k + eq + d * num_equation + dof1 * num_equation * dim) -= du2n * shape2(k);
145-
#else
146130
elvect2_mat(k, eq + d * num_equation) -= du2n * shape2(k);
147-
#endif
148131
}
149132
}
150133
}

src/gradNonLinearForm.cpp

+5-73
Original file line numberDiff line numberDiff line change
@@ -39,82 +39,14 @@ GradNonLinearForm::GradNonLinearForm(ParFiniteElementSpace *_vfes, IntegrationRu
3939

4040
void GradNonLinearForm::Mult(const ParGridFunction *Up, Vector &y) {
4141
Vector x;
42-
43-
#ifdef _GPU_
44-
// TestNonLinearForm::setToZero_gpu(y,y.Size()); already done in gradient.cpd
45-
x.UseDevice(true);
46-
if (fnfi.Size()) {
47-
x.NewMemoryAndSize(Up->GetMemory(), vfes->GetNDofs() * num_equation, false);
48-
Mult_gpu(x, y);
49-
}
50-
#else
5142
const double *data = Up->GetData();
43+
44+
// NB: Setting x here accounts for the fact that the space Up lives
45+
// in is different from the space of the gradient.
5246
x.SetSize(vfes->GetNDofs() * num_equation * dim);
5347
x = 0.;
5448
for (int i = 0; i < vfes->GetNDofs() * num_equation; i++) x(i) = data[i];
55-
ParNonlinearForm::Mult(x, y);
56-
#endif
57-
}
58-
59-
#ifdef _GPU_
60-
void GradNonLinearForm::Mult_gpu(const Vector &x, Vector &y) {
61-
// INTEGRATION SHARED FACES
62-
const Vector &px = Prolongate(x);
63-
64-
const int dofs = vfes->GetNDofs();
65-
66-
// Terms over shared interior faces in parallel.
67-
ParFiniteElementSpace *pfes = ParFESpace();
68-
ParMesh *pmesh = pfes->GetParMesh();
69-
FaceElementTransformations *tr;
70-
const FiniteElement *fe1, *fe2;
71-
Array<int> vdofs1, vdofs2;
72-
Vector el_x, el_y;
73-
el_x.UseDevice(true);
74-
el_y.UseDevice(true);
75-
76-
// if (bfnfi.Size()==0) y.HostReadWrite();
77-
X.MakeRef(aux1, 0); // aux1 contains P.x
78-
X.ExchangeFaceNbrData();
79-
const int n_shared_faces = pmesh->GetNSharedFaces();
80-
for (int i = 0; i < n_shared_faces; i++) {
81-
tr = pmesh->GetSharedFaceTransformations(i, true);
82-
int Elem2NbrNo = tr->Elem2No - pmesh->GetNE();
8349

84-
fe1 = pfes->GetFE(tr->Elem1No);
85-
fe2 = pfes->GetFaceNbrFE(Elem2NbrNo);
86-
87-
pfes->GetElementVDofs(tr->Elem1No, vdofs1);
88-
pfes->GetFaceNbrElementVDofs(Elem2NbrNo, vdofs2);
89-
90-
el_x.SetSize(vdofs1.Size() + vdofs2.Size());
91-
X.GetSubVector(vdofs1, el_x.GetData());
92-
X.FaceNbrData().GetSubVector(vdofs2, el_x.GetData() + vdofs1.Size());
93-
94-
const int elDof1 = fe1->GetDof();
95-
const int elDof2 = fe2->GetDof();
96-
for (int k = 0; k < fnfi.Size(); k++) {
97-
fnfi[k]->AssembleFaceVector(*fe1, *fe2, *tr, el_x, el_y);
98-
// y.AddElementVector(vdofs1, el_y.GetData());
99-
GradNonLinearForm::IndexedAddToGlobalMemory(vdofs1, y, el_y, dofs, elDof1, num_equation, dim);
100-
}
101-
}
102-
}
103-
104-
void GradNonLinearForm::IndexedAddToGlobalMemory(const Array<int> &vdofs, Vector &y, const Vector &el_y, const int &dof,
105-
const int &elDof, const int &num_equation, const int &dim) {
106-
double *dy = y.ReadWrite();
107-
const double *d_ely = el_y.Read();
108-
auto d_vdofs = vdofs.Read();
109-
110-
MFEM_FORALL(i, elDof, {
111-
const int index = d_vdofs[i];
112-
for (int eq = 0; eq < num_equation; eq++) {
113-
for (int d = 0; d < dim; d++) {
114-
dy[index + eq * dof + d * num_equation * dof] += d_ely[i + eq * elDof + d * num_equation * elDof];
115-
}
116-
}
117-
});
50+
// After setting x as above, base class Mult does the right thing
51+
ParNonlinearForm::Mult(x, y);
11852
}
119-
120-
#endif // _GPU_

src/gradNonLinearForm.hpp

-6
Original file line numberDiff line numberDiff line change
@@ -51,12 +51,6 @@ class GradNonLinearForm : public ParNonlinearForm {
5151
GradNonLinearForm(ParFiniteElementSpace *f, IntegrationRules *intRules, const int dim, const int num_equation);
5252

5353
void Mult(const ParGridFunction *Up, Vector &y);
54-
55-
#ifdef _GPU_
56-
void Mult_gpu(const Vector &x, Vector &y);
57-
static void IndexedAddToGlobalMemory(const Array<int> &vdofs, Vector &y, const Vector &el_y, const int &dof,
58-
const int &elDof, const int &num_equation, const int &dim);
59-
#endif
6054
};
6155

6256
#endif // GRADNONLINEARFORM_HPP_

0 commit comments

Comments
 (0)