Skip to content

Commit 12eaa25

Browse files
committed
Removed debug comments
1 parent d8260f4 commit 12eaa25

3 files changed

Lines changed: 61 additions & 36 deletions

File tree

src/pmpo_MPMesh.cpp

Lines changed: 46 additions & 31 deletions
Original file line numberDiff line numberDiff line change
@@ -52,10 +52,10 @@ void MPMesh::calculateStrain(){
5252
v22 = v22 + MPsBasisGrads(mp, i*2 + 1) * velField(iVertex, 1);
5353
uTanOverR = uTanOverR + MPsBasis(mp, i) * tanLatVertexRotatedOverRadius(iVertex, 0) * velField(iVertex, 0);
5454
vTanOverR = vTanOverR + MPsBasis(mp, i) * tanLatVertexRotatedOverRadius(iVertex, 0) * velField(iVertex, 1);
55-
//Debugging
56-
if(MPsAppID(mp)==0){
57-
//printf("Strain Calc: iVertex %d vel field %.15e %.15e \n", iVertex, velField(iVertex, 0), velField(iVertex, 1));
58-
}
55+
/*
56+
if(MPsAppID(mp)==0)
57+
printf("Strain Calc: iVertex %d vel field %.15e %.15e \n", iVertex, velField(iVertex, 0), velField(iVertex, 1));
58+
*/
5959
}
6060
MPsStrainRate(mp, 0) = v11 - vTanOverR;
6161
MPsStrainRate(mp, 1) = v22;
@@ -84,21 +84,28 @@ void MPMesh::calculateStress(){
8484

8585
Vec3d strain_rate (MPsStrainRate(mp, 0), MPsStrainRate(mp, 1), MPsStrainRate(mp, 2));
8686
Vec3d stress(MPsStress(mp, 0), MPsStress(mp, 1), MPsStress(mp, 2));
87-
constitutive_evp(strain_rate, stress, MPsIcePressure(mp, 0), MPsRepPressure(mp, 0), MPsArea(mp, 0), elasticTimeStep, dampingTimescale);
87+
//constitutive_evp(strain_rate, stress, MPsIcePressure(mp, 0), MPsRepPressure(mp, 0), MPsArea(mp, 0), elasticTimeStep, dampingTimescale);
88+
constitutive_linear(strain_rate, stress);
8889
for (int m=0 ; m<3; m++)
8990
MPsStress(mp, m) = stress[m];
90-
//Debugging
91+
/*
9192
if(MPsAppID(mp)==0){
92-
//printf("Strain in GPU: %.15e %.15e %.15e\n", MPsStrainRate(mp, 0), MPsStrainRate(mp, 1), MPsStrainRate(mp, 2));
93-
//printf("Stress in GPU: %.15e %.15e %.15e\n", MPsStress(mp, 0), MPsStress(mp, 1), MPsStress(mp, 2));
93+
printf("Strain in GPU: %.15e %.15e %.15e\n", MPsStrainRate(mp, 0), MPsStrainRate(mp, 1), MPsStrainRate(mp, 2));
94+
printf("Stress in GPU: %.15e %.15e %.15e\n", MPsStress(mp, 0), MPsStress(mp, 1), MPsStress(mp, 2));
9495
}
96+
*/
9597
}
9698
};
9799
p_MPs->parallel_for(setMPStress, "setMPStress");
98100
}
99101

100102
void MPMesh::calculateStressDivergence(){
101-
103+
104+
int self, numProcsTot;
105+
MPI_Comm comm = p_MPs->getMPIComm();
106+
MPI_Comm_rank(comm, &self);
107+
MPI_Comm_size(comm, &numProcsTot);
108+
102109
//Mesh Information
103110
auto elm2VtxConn = p_mesh->getElm2VtxConn();
104111
int numVtx = p_mesh->getNumVertices();
@@ -125,11 +132,11 @@ void MPMesh::calculateStressDivergence(){
125132
radius=p_mesh->getSphereRadius();
126133

127134
//Reconstructed the stress
128-
Kokkos::View<double*> stress_divU("stress_divu", p_mesh->getNumVertices());
129-
Kokkos::View<double*> stress_divV("stress_divv", p_mesh->getNumVertices());
135+
Kokkos::View<doubleSclr_t*> stress_divU("stress_divu", p_mesh->getNumVertices());
136+
Kokkos::View<doubleSclr_t*> stress_divV("stress_divv", p_mesh->getNumVertices());
130137

131-
Kokkos::View<double*> divU_edge("divUedge", p_mesh->getNumVertices());
132-
Kokkos::View<double*> divV_edge("divVedge", p_mesh->getNumVertices());
138+
Kokkos::View<doubleSclr_t*> divU_edge("divUedge", p_mesh->getNumVertices());
139+
Kokkos::View<doubleSclr_t*> divV_edge("divVedge", p_mesh->getNumVertices());
133140

134141
//Assemble fields for Stress Divergence
135142
auto stress_div = PS_LAMBDA(const int& elm, const int& mp, const int& mask) {
@@ -153,45 +160,53 @@ void MPMesh::calculateStressDivergence(){
153160
VtxCoeffs_new(vID,2, 2)*CoordDiffs[2] +
154161
VtxCoeffs_new(vID,2, 3)*CoordDiffs[3]);
155162

156-
Kokkos::atomic_add(&stress_divU(vID), factor1 * MPsStress(mp, 0) + factor2 * MPsStress(mp, 2) -
163+
Kokkos::atomic_add(&stress_divU(vID, 0), factor1 * MPsStress(mp, 0) + factor2 * MPsStress(mp, 2) -
157164
2 * tanLatVertexRotatedOverRadius(vID, 0) * factor * MPsStress(mp, 2));
158-
Kokkos::atomic_add(&stress_divV(vID), factor2 * MPsStress(mp, 1) + factor1*MPsStress(mp, 2) +
165+
Kokkos::atomic_add(&stress_divV(vID, 0), factor2 * MPsStress(mp, 1) + factor1*MPsStress(mp, 2) +
159166
factor * tanLatVertexRotatedOverRadius(vID, 0) * (MPsStress(mp, 0)-MPsStress(mp, 1)));
160167

161168

162-
Kokkos::atomic_add(&divU_edge(vID), - weight_grads(mp, i*2 + 0) * MPsStress(mp, 0) - weight_grads(mp, i*2 + 1) * MPsStress(mp, 2) -
169+
Kokkos::atomic_add(&divU_edge(vID, 0), - weight_grads(mp, i*2 + 0) * MPsStress(mp, 0) - weight_grads(mp, i*2 + 1) * MPsStress(mp, 2) -
163170
2 * tanLatVertexRotatedOverRadius(vID, 0) * w_vtx * MPsStress(mp, 2));
164-
Kokkos::atomic_add(&divV_edge(vID), - weight_grads(mp, i*2 + 1) * MPsStress(mp, 1) - weight_grads(mp, i*2 + 0) * MPsStress(mp, 2) +
171+
Kokkos::atomic_add(&divV_edge(vID, 0), - weight_grads(mp, i*2 + 1) * MPsStress(mp, 1) - weight_grads(mp, i*2 + 0) * MPsStress(mp, 2) +
165172
w_vtx * tanLatVertexRotatedOverRadius(vID, 0) * (MPsStress(mp, 0)- MPsStress(mp, 1)));
166173

167174
}
168175
}
169176
};
170-
p_MPs->parallel_for(stress_div, "assembly");
171-
172-
//TODO COMMUNICATE THE VERTEX FIELDS
173-
177+
p_MPs->parallel_for(stress_div, " stress_div_assembly");
178+
179+
//TO DO make to 2 fields instrad of 4 to reduce latency
180+
if(numProcsTot>1){
181+
//Takes contribution of halo vertices and adds it in halo cells
182+
communicate_and_take_halo_contributions(stress_divU, numVertices, 1, 0, 0);
183+
communicate_and_take_halo_contributions(stress_divV, numVertices, 1, 0, 0);
184+
communicate_and_take_halo_contributions(divU_edge, numVertices, 1, 0, 0);
185+
communicate_and_take_halo_contributions(divV_edge, numVertices, 1, 0, 0);
186+
//Does it need to transfer the correct values at owned vertices to halo vertices
187+
communicate_and_take_halo_contributions(stress_divU, numVertices, 1, 1, 1);
188+
communicate_and_take_halo_contributions(stress_divV, numVertices, 1, 1, 1);
189+
communicate_and_take_halo_contributions(divU_edge, numVertices, 1, 1, 1);
190+
communicate_and_take_halo_contributions(divV_edge, numVertices, 1, 1, 1);
191+
}
192+
174193
auto stressDivergence = p_mesh->getMeshField<MeshF_StressDivergence>();
175-
176194
Kokkos::parallel_for("calculate_divergence", numVtx, KOKKOS_LAMBDA(const int vtx){
177195

178196
double ramp = nearAnEdge_l(vtx);
179197
double invM = 1.0/vtxMatrixMass_l(vtx);
180198
invM = vtxMatrixMass_l(vtx) >1e-4 ? invM : 0;
181199

182-
stressDivergence(vtx, 0) = ramp * stress_divU(vtx) + (1 - ramp) * divU_edge(vtx) * invM ;
183-
stressDivergence(vtx, 1) = ramp * stress_divV(vtx) + (1 - ramp) * divV_edge(vtx) * invM ;
184-
//Debugging
185-
if (vtx == 10 || vtx == 11 || vtx == 1476 || vtx == 1481) {
186-
//printf("Vtx %d Divergence %.15e %.15e \n", vtx, stressDivergence(vtx, 0), stressDivergence(vtx, 1));
200+
stressDivergence(vtx, 0) = ramp * stress_divU(vtx, 0) + (1 - ramp) * divU_edge(vtx, 0) * invM ;
201+
stressDivergence(vtx, 1) = ramp * stress_divV(vtx, 0) + (1 - ramp) * divV_edge(vtx, 0) * invM ;
202+
/*
203+
if (ent2global(vtx) == 1714) {
204+
printf("Vtx %d Divergence %.15e %.15e %.15e %.15e \n", vtx, stressDivergence(vtx, 0), stressDivergence(vtx, 1));
187205
}
206+
*/
188207
});
189-
190-
//TODO COMMUNICATE THE VERTEX FIELDS
191-
192208
}
193209

194-
195210
void MPMesh::calcBasis() {
196211
assert(p_mesh->getGeomType() == geom_spherical_surf);
197212

src/pmpo_MPMesh_assembly.hpp

Lines changed: 7 additions & 5 deletions
Original file line numberDiff line numberDiff line change
@@ -96,13 +96,15 @@ void MPMesh::assemblyElm0() {
9696
}
9797

9898
void MPMesh::reconstruct_coeff_full(){
99-
std::cout<<__FUNCTION__<<std::endl;
10099
Kokkos::Timer timer;
101100
int self, numProcsTot;
102101
MPI_Comm comm = p_MPs->getMPIComm();
103102
MPI_Comm_rank(comm, &self);
104103
MPI_Comm_size(comm, &numProcsTot);
105-
104+
105+
static int coeff_count=0;
106+
if(!self) std::cout<<"===="<<__FUNCTION__<<" "<<coeff_count<<"===="<<std::endl;
107+
coeff_count++;
106108
//Mesh Information
107109
auto elm2VtxConn = p_mesh->getElm2VtxConn();
108110
int numVtx = p_mesh->getNumVertices();
@@ -177,7 +179,6 @@ void MPMesh::reconstruct_coeff_full(){
177179
}
178180

179181
void MPMesh::invertMatrix(const Kokkos::View<double**>& vtxMatrices, const double& radius){
180-
std::cout<<__FUNCTION__<<std::endl;
181182

182183
int nVertices = p_mesh->getNumVertices();
183184
auto vtxCoords = p_mesh->getMeshField<polyMPO::MeshF_VtxCoords>();
@@ -320,7 +321,6 @@ void MPMesh::invertMatrix(const Kokkos::View<double**>& vtxMatrices, const doubl
320321

321322
template <MeshFieldIndex meshFieldIndex>
322323
void MPMesh::assemblyVtx1(){
323-
std::cout<<__FUNCTION__<<std::endl;
324324
Kokkos::Timer timer;
325325

326326
int self, numProcsTot;
@@ -379,8 +379,10 @@ void MPMesh::assemblyVtx1(){
379379
pumipic::RecordTime("Assemble Field per process" + std::to_string(self), timer.seconds());
380380

381381
timer.reset();
382-
if(numProcsTot>1)
382+
if(numProcsTot>1){
383383
communicate_and_take_halo_contributions(meshField, numVertices, numEntries, 0, 0);
384+
communicate_and_take_halo_contributions(meshField, numVertices, numEntries, 1, 1);
385+
}
384386
pumipic::RecordTime("Communicate Field Values" + std::to_string(self), timer.seconds());
385387
}
386388

src/pmpo_const_relation.hpp

Lines changed: 8 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -37,6 +37,14 @@ void constitutive_evp(const Vec3d& strain, Vec3d& stress, const double& icePress
3737
stress[1] = 0.5 * (stress1 - stress2);
3838
}
3939

40+
KOKKOS_INLINE_FUNCTION
41+
void constitutive_linear(const Vec3d& strain, Vec3d& stress){
42+
double lambda = 1.0;
43+
stress[0] = lambda * strain[0];
44+
stress[1] = lambda * strain[1];
45+
stress[2] = lambda * strain[2];
46+
}
47+
4048
}
4149
#endif
4250

0 commit comments

Comments
 (0)