Skip to content

Commit

Permalink
Cekees/mbd/chrono (#570)
Browse files Browse the repository at this point in the history
* re-enabled Nitsche adjoint for pinc

* updated optimized pressure calculation for projection scheme

* enabled optimized pressure solve for rans3p

* enabled optimized pressure solve for rans3p

* fixed mass flux for rotational form of pressure calc

* worked on isotropic interface size field

* fix 2D Darcy term and isotropic adaptive scheme

* update stack with fixed zoltan build
  • Loading branch information
cekees authored Jul 17, 2017
1 parent 0a00eaf commit 5bc5dcc
Show file tree
Hide file tree
Showing 8 changed files with 43 additions and 21 deletions.
2 changes: 1 addition & 1 deletion .hashstack_default
Original file line number Diff line number Diff line change
@@ -1 +1 @@
d3d27601f2de3c4538f11972f68973b142bdbc92
d88255e2a2741d2680b099135c13fa30dd8054fe
32 changes: 25 additions & 7 deletions proteus/MeshAdaptPUMI/SizeField.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -17,13 +17,9 @@ static void SmoothField(apf::Field* f);
thickness of refinement near the interface */
static double isotropicFormula(double phi, double hmin, double hmax)
{
static double const epsilon = 0.02;
phi = fabs(phi);
double size;
if (fabs(phi) < epsilon)
if (fabs(phi) < 5.0*hmin)
size = hmin;
else if (phi < 3 * epsilon)
size = (hmin + hmax) / 2;
else
size = hmax;
return size;
Expand All @@ -43,8 +39,30 @@ int MeshAdaptPUMIDrvr::calculateSizeField()
apf::setScalar(size_iso, v, 0, size);
}
m->end(it);
for(int i=0; i < 3; i++)
SmoothField(size_iso);
/*
If you just smooth then hmax will just diffuse into the hmin band
and you won't really get a band around phi=0 with uniform diameter
hmin. Instead, reset to hmin after each smooth within the band in
order to ensure the band uses hmin. Iterate on that process until
changes in the smoothed size are less than 50% of hmin.
*/
double err_h_max=hmax;
int its=0;
while (err_h_max > 0.5*hmin && its < 200)
{
its++;
SmoothField(size_iso);
err_h_max=0.0;
it = m->begin(0);
while ((v = m->iterate(it))) {
double phi = apf::getScalar(phif, v, 0);
double size_current = apf::getScalar(size_iso, v, 0);
double size = fmin(size_current,isotropicFormula(phi, hmin, hmax));
err_h_max = fmax(err_h_max,fabs(size_current-size));
apf::setScalar(size_iso, v, 0, size);
}
m->end(it);
}
return 0;
}

Expand Down
2 changes: 1 addition & 1 deletion proteus/MeshAdaptPUMI/cMeshAdaptPUMI.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -351,7 +351,7 @@ int MeshAdaptPUMIDrvr::adaptPUMIMesh()
size_iso = samSz::isoSize(m);
}
else if (size_field_config == "isotropic")
testIsotropicSizeField();
calculateSizeField();
else {
std::cerr << "unknown size field config " << size_field_config << '\n';
abort();
Expand Down
8 changes: 5 additions & 3 deletions proteus/NumericalSolution.py
Original file line number Diff line number Diff line change
Expand Up @@ -1239,7 +1239,7 @@ def calculateSolution(self,runName):
# print "Min / Max residual %s / %s" %(lr.min(),lr.max())

self.nSequenceSteps = 0
self.nSolveSteps = 0
self.nSolveSteps=self.nList[0].adaptMesh_nSteps-1
for (self.tn_last,self.tn) in zip(self.tnList[:-1],self.tnList[1:]):
logEvent("==============================================================",level=0)
logEvent("Solving over interval [%12.5e,%12.5e]" % (self.tn_last,self.tn),level=0)
Expand Down Expand Up @@ -1381,7 +1381,10 @@ def calculateSolution(self,runName):
self.tCount+=1
for index,model in enumerate(self.modelList):
self.archiveSolution(model,index,self.systemStepController.t_system_last)

#can only handle PUMIDomain's for now
self.nSolveSteps += 1
if(self.PUMI_estimateError()):
self.PUMI_adaptMesh()
#end system step iterations
if self.archiveFlag == ArchiveFlags.EVERY_USER_STEP:
self.tCount+=1
Expand All @@ -1398,7 +1401,6 @@ def calculateSolution(self,runName):
self.nSolveSteps += 1
if(self.PUMI_estimateError()):
self.PUMI_adaptMesh()

logEvent("Finished calculating solution",level=3)

for index,model in enumerate(self.modelList):
Expand Down
2 changes: 1 addition & 1 deletion proteus/mprans/MCorr.h
Original file line number Diff line number Diff line change
Expand Up @@ -1909,7 +1909,7 @@ namespace proteus
/* ck.calculateGScale(G,dir,h_phi); */
epsHeaviside=epsFactHeaviside*(useMetrics*h_phi+(1.0-useMetrics)*elementDiameter[eN]);
q_H[eN_k] = q_porosity[eN_k]*smoothedHeaviside(epsHeaviside,q_phi[eN_k]);
}//k
}//k
for (int i=0;i<nDOF_trial_element;i++)
{
int eN_i = eN*nDOF_trial_element + i;
Expand Down
5 changes: 3 additions & 2 deletions proteus/mprans/MCorr.py
Original file line number Diff line number Diff line change
Expand Up @@ -80,6 +80,7 @@ def attachModels(self,modelList):
self.ebq_H_vof = None
#correction
self.massCorrModel = modelList[self.me_model]
self.massCorrModel.setMassQuadrature()
self.vofModel.q[('m_last',0)][:] = self.vofModel.q[('m',0)]
if self.checkMass:
self.m_tmp = copy.deepcopy(self.massCorrModel.q[('r',0)])
Expand Down Expand Up @@ -419,6 +420,7 @@ def __init__(self,
self.ebqe={}
self.phi_ip={}
#mesh
#uncomment this to store q arrays, see calculateElementQuadrature below
#self.q['x'] = numpy.zeros((self.mesh.nElements_global,self.nQuadraturePoints_element,3),'d')
self.q[('u',0)] = numpy.zeros((self.mesh.nElements_global,self.nQuadraturePoints_element),'d')
self.q[('grad(u)',0)] = numpy.zeros((self.mesh.nElements_global,self.nQuadraturePoints_element,self.nSpace_global),'d')
Expand Down Expand Up @@ -863,6 +865,7 @@ def calculateElementQuadrature(self):
This function should be called only when the mesh changes.
"""
#uncomment this to store q arrays
#self.u[0].femSpace.elementMaps.getValues(self.elementQuadraturePoints,
# self.q['x'])
self.u[0].femSpace.elementMaps.getBasisValuesRef(self.elementQuadraturePoints)
Expand Down Expand Up @@ -1227,8 +1230,6 @@ def initializeTimeHistory(self):
for m,u,r in zip(self.model.levelModelList,
self.model.uList,
self.model.rList):
#pass
#m.setMassQuadrature()
u.flat[:]=0.0
m.getResidual(u,r)
m.coefficients.postStep(self.t_model)
Expand Down
11 changes: 5 additions & 6 deletions proteus/mprans/RANS2P2D.h
Original file line number Diff line number Diff line change
Expand Up @@ -862,11 +862,11 @@ namespace proteus
/* duc_du = u/(uc+1.0e-12); */
/* duc_dv = v/(uc+1.0e-12); */
//semi-implicit quadratic term
uc = sqrt(uStar*uStar+vStar*vStar*+wStar*wStar);
uc = sqrt(uStar*uStar+vStar*vStar);
duc_du = 0.0;
duc_dv = 0.0;
/* duc_dw = w/(uc+1.0e-12); */

/*
mom_u_source += H_s*viscosity*(alpha + beta*uc)*(u-u_s);
mom_v_source += H_s*viscosity*(alpha + beta*uc)*(v-v_s);
/* mom_w_source += H_s*viscosity*(alpha + beta*uc)*(w-w_s); */
Expand Down Expand Up @@ -2007,7 +2007,7 @@ namespace proteus
mom_v_ham,
dmom_v_ham_grad_p,
mom_w_ham,
dmom_w_ham_grad_p);
dmom_w_ham_grad_p);
//VRANS
mass_source = q_mass_source[eN_k];
//todo: decide if these should be lagged or not?
Expand All @@ -2031,7 +2031,7 @@ namespace proteus
w,
q_velocity_sge[eN_k_nSpace+0],
q_velocity_sge[eN_k_nSpace+1],
q_velocity_sge[eN_k_nSpace+2],
q_velocity_sge[eN_k_nSpace+1],//cek hack, should not be used
eps_solid[elementFlags[eN]],
phi_solid[eN_k],
q_velocity_solid[eN_k_nSpace+0],
Expand Down Expand Up @@ -2150,7 +2150,6 @@ namespace proteus
dmom_adv_sge[0] = dmom_u_acc_u*(q_velocity_sge[eN_k_nSpace+0] - MOVING_DOMAIN*xt);
dmom_adv_sge[1] = dmom_u_acc_u*(q_velocity_sge[eN_k_nSpace+1] - MOVING_DOMAIN*yt);
/* dmom_adv_sge[2] = dmom_u_acc_u*(q_velocity_sge[eN_k_nSpace+2] - MOVING_DOMAIN*zt); */

pdeResidual_u = ck.Mass_strong(mom_u_acc_t) +
ck.Advection_strong(dmom_adv_sge,grad_u) +
ck.Hamiltonian_strong(dmom_u_ham_grad_p,grad_p) +
Expand Down Expand Up @@ -3638,7 +3637,7 @@ namespace proteus
w,
q_velocity_sge[eN_k_nSpace+0],
q_velocity_sge[eN_k_nSpace+1],
q_velocity_sge[eN_k_nSpace+2],
q_velocity_sge[eN_k_nSpace+1],//cek hack, should not be used
eps_solid[elementFlags[eN]],
phi_solid[eN_k],
q_velocity_solid[eN_k_nSpace+0],
Expand Down
2 changes: 2 additions & 0 deletions proteus/mprans/VOF.py
Original file line number Diff line number Diff line change
Expand Up @@ -498,6 +498,7 @@ def __init__(self,
self.ebqe={}
self.phi_ip={}
#mesh
#uncomment this to store at q['x''] see also calculateElementQuadrature below
#self.q['x'] = numpy.zeros((self.mesh.nElements_global,self.nQuadraturePoints_element,3),'d')
self.ebqe['x'] = numpy.zeros((self.mesh.nExteriorElementBoundaries_global,self.nElementBoundaryQuadraturePoints_elementBoundary,3),'d')
self.q[('u',0)] = numpy.zeros((self.mesh.nElements_global,self.nQuadraturePoints_element),'d')
Expand Down Expand Up @@ -827,6 +828,7 @@ def calculateElementQuadrature(self):
This function should be called only when the mesh changes.
"""
#uncomment this to store at q['x']
#self.u[0].femSpace.elementMaps.getValues(self.elementQuadraturePoints,
# self.q['x'])
self.u[0].femSpace.elementMaps.getBasisValuesRef(self.elementQuadraturePoints)
Expand Down

0 comments on commit 5bc5dcc

Please sign in to comment.