Skip to content

Commit dd1bb9a

Browse files
Pushback v to v_list, reset v to 0 for each face, make v_list a constant reference,
1 parent 9d4ba93 commit dd1bb9a

File tree

2 files changed

+8
-3
lines changed

2 files changed

+8
-3
lines changed

Code/Source/solver/trilinos_impl.cpp

Lines changed: 7 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -951,7 +951,7 @@ void constructJacobiScaling(const double *dirW, Epetra_Vector &diagonal)
951951
* \param v coupled boundary vector
952952
* \param isCoupledBC determines if coupled resistance BC is turned on
953953
*/
954-
void trilinos_bc_create_(const std::vector<Array<double>> v_list, bool &isCoupledBC)
954+
void trilinos_bc_create_(const std::vector<Array<double>> &v_list, bool &isCoupledBC)
955955
{
956956
//store as global to determine which matvec multiply to use in solver
957957
coupledBC = isCoupledBC;
@@ -1246,18 +1246,23 @@ void TrilinosLinearAlgebra::TrilinosImpl::init_dir_and_coup_neu(ComMod& com_mod,
12461246
bool isCoupledBC = false;
12471247

12481248
for (int faIn = 0; faIn < lhs.nFaces; faIn++) {
1249-
auto& face = lhs.face[faIn];
1249+
auto& face = lhs.face[faIn]; // Extract face
1250+
v = 0.0; // Initialize v to zero for each face
12501251
if (face.coupledFlag) {
12511252
isCoupledBC = true;
12521253
int faDof = std::min(face.dof,dof);
12531254

1255+
// Compute the coupled Neumann BC vector
12541256
for (int a = 0; a < face.nNo; a++) {
12551257
int Ac = face.glob(a);
12561258
for (int i = 0; i < faDof; i++) {
12571259
v(i,Ac) = v(i,Ac) + sqrt(fabs(res(faIn))) * face.val(i,a);
12581260
}
12591261
}
12601262
}
1263+
1264+
// Add the coupled Neumann BC vector to the list
1265+
v_list.push_back(v);
12611266
}
12621267

12631268
// Add the v vectors to global bdryVec_list

Code/Source/solver/trilinos_impl.h

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -201,7 +201,7 @@ class TrilinosMatVec: public virtual Epetra_Operator
201201
* \param v coeff in the scalar product
202202
* \param isCoupledBC determines if coupled resistance BC is turned on
203203
*/
204-
void trilinos_bc_create_(const std::vector<Array<double>> v_list, bool &isCoupledBC);
204+
void trilinos_bc_create_(const std::vector<Array<double>> &v_list, bool &isCoupledBC);
205205

206206
void trilinos_doassem_(int &numNodesPerElement, const int *eqN,
207207
const double *lK, double *lR);

0 commit comments

Comments
 (0)