Added methods projectSolution and injectNodeVec

git-svn-id: http://svn.sintef.no/trondheim/IFEM/trunk@868 e10b68d5-8a6e-419e-a041-bce267b0401d
This commit is contained in:
kmo 2011-03-22 08:21:04 +00:00 committed by Knut Morten Okstad
parent 36148701c2
commit cbc130117a
9 changed files with 234 additions and 51 deletions

View File

@ -479,6 +479,20 @@ void ASMbase::extractNodeVec (const Vector& globRes, Vector& nodeVec,
}
void ASMbase::injectNodeVec (const Vector& nodeVec, Vector& globRes,
unsigned char nndof) const
{
if (nndof == 0) nndof = nf;
for (size_t i = 0; i < MLGN.size(); i++)
{
int n = MLGN[i]-1;
for (unsigned char j = 0; j < nndof; j++)
globRes[nndof*n+j] = nodeVec[nndof*i+j];
}
}
bool ASMbase::tesselate (ElementBlock&, const int*) const
{
std::cerr <<" *** ASMBase::tesselate: Must be implemented in sub-class."

View File

@ -245,8 +245,11 @@ public:
//!
//! \details The secondary solution is derived from the primary solution,
//! which is assumed to be stored in the \a integrand for current patch.
//! If \a npe is NULL, the solution is evaluated at the Greville points and
//! then projected onto the spline basis to obtain the control point values,
//! which then are returned through \a sField.
virtual bool evalSolution(Matrix& sField, const Integrand& integrand,
const int* npe) const;
const int* npe = 0) const;
// Methods for result extraction
@ -264,6 +267,13 @@ public:
virtual void extractNodeVec(const Vector& globVec, Vector& nodeVec,
unsigned char nndof = 0) const;
//! \brief Inject nodal results for this patch into the global vector.
//! \param[in] nodeVec Nodal result vector for this patch
//! \param globVec Global solution vector in DOF-order
//! \param[in] nndof Number of DOFs per node (the default is \a nf)
virtual void injectNodeVec(const Vector& nodeVec, Vector& globVec,
unsigned char nndof = 0) const;
protected:
// Internal methods for preprocessing of boundary conditions

View File

@ -729,12 +729,47 @@ bool ASMs1D::evalSolution (Matrix& sField, const Vector& locSol,
bool ASMs1D::evalSolution (Matrix& sField, const Integrand& integrand,
const int* npe) const
{
sField.resize(0,0);
if (npe)
{
// Compute parameter values of the result sampling points
DoubleVec gpar;
if (this->getGridParameters(gpar,npe[0]-1))
// Evaluate the secondary solution at all sampling points
return this->evalSolution(sField,integrand,gpar);
}
else
{
Go::SplineCurve* c = this->projectSolution(integrand);
if (c)
{
sField.resize(c->dimension(),c->numCoefs());
sField.fill(&(*c->coefs_begin()));
delete c;
return true;
}
}
// Compute parameter values of the result sampling points
DoubleVec gpar;
if (!this->getGridParameters(gpar,npe[0]-1))
return false;
return false;
}
Go::GeomObject* ASMs1D::evalSolution (const Integrand& integrand) const
{
return this->projectSolution(integrand);
}
Go::SplineCurve* ASMs1D::projectSolution (const Integrand& integrand) const
{
std::cerr <<"ASMs1D::projectSolution: Not implemented yet!"<< std::endl;
return 0;
}
bool ASMs1D::evalSolution (Matrix& sField, const Integrand& integrand,
const RealArray& gpar) const
{
sField.resize(0,0);
const int p1 = curv->order();

View File

@ -139,14 +139,31 @@ public:
//! \param[out] sField Solution field
//! \param[in] integrand Object with problem-specific data and methods
//! \param[in] npe Number of visualization nodes over each knot span
//!
//! \details If \a npe is NULL, the solution is evaluated at the Greville
//! points and then projected onto the spline basis to obtain the control
//! point values, which then are returned through \a sField.
virtual bool evalSolution(Matrix& sField, const Integrand& integrand,
const int* npe) const;
const int* npe = 0) const;
//! \brief Projects the secondary solution field onto the primary basis.
//! \param[in] integrand Object with problem-specific data and methods
Go::SplineCurve* projectSolution(const Integrand& integrand) const;
//! \brief Projects the secondary solution field onto the primary basis.
virtual Go::GeomObject* evalSolution(const Integrand& integrand) const;
protected:
// Internal utility methods
// ========================
//! \brief Evaluates the secondary solution field at the given points.
//! \param[out] sField Solution field
//! \param[in] integrand Object with problem-specific data and methods
//! \param[in] gpar Parameter values of the result sampling points
bool evalSolution(Matrix& sField, const Integrand& integrand,
const RealArray& gpar) const;
//! \brief Calculates parameter values for the visualization nodal points.
//! \param[out] prm Parameter values for all points
//! \param[in] nSegSpan Number of visualization segments over each knot-span

View File

@ -754,6 +754,19 @@ void ASMs2D::extractBasis (const Go::BasisDerivsSf2& spline,
}
#if SP_DEBUG > 4
std::ostream& operator<<(std::ostream& os, const Go::BasisDerivsSf& bder)
{
os <<" : (u,v) = "<< bder.param[0] <<" "<< bder.param[1]
<<" left_idx = "<< bder.left_idx[0] <<" "<< bder.left_idx[1] << std::endl;
for (size_t i = 0; i < bder.basisValues.size(); i++)
os << 1+i <<'\t'<< bder.basisValues[i] <<'\t'
<< bder.basisDerivs_u[i] <<'\t'<< bder.basisDerivs_v[i] << std::endl;
return os;
}
#endif
bool ASMs2D::integrate (Integrand& integrand,
GlobalIntegral& glInt,
const TimeDomain& time,
@ -795,6 +808,11 @@ bool ASMs2D::integrate (Integrand& integrand,
else
surf->computeBasisGrid(gpar[0],gpar[1],spline);
#if SP_DEBUG > 4
for (size_t i = 0; i < spline.size(); i++)
std::cout <<"\nBasis functions at integration point "<< 1+i << spline[i];
#endif
const int p1 = surf->order_u();
const int p2 = surf->order_v();
const int n1 = surf->numCoefs_u();
@ -899,6 +917,11 @@ bool ASMs2D::integrate (Integrand& integrand,
if (!utl::Hessian(Hess,d2NdX2,Jac,Xnod,d2Ndu2,dNdu))
return false;
#if SP_DEBUG > 4
std::cout <<"\niel, ip = "<< iel <<" "<< ip
<<"\nN ="<< N <<"dNdX ="<< dNdX << std::endl;
#endif
// Cartesian coordinates of current integration point
X = Xnod * N;
X.t = time.t;
@ -1208,20 +1231,54 @@ bool ASMs2D::evalSolution (Matrix& sField, const Vector& locSol,
bool ASMs2D::evalSolution (Matrix& sField, const Integrand& integrand,
const int* npe) const
{
sField.resize(0,0);
// Compute parameter values of the result sampling points
DoubleVec gpar[2];
bool retVal = true;
for (int dir = 0; dir < 2; dir++)
if (npe)
if (npe)
{
// Compute parameter values of the result sampling points
DoubleVec gpar[2];
for (int dir = 0; dir < 2 && retVal; dir++)
retVal = this->getGridParameters(gpar[dir],dir,npe[dir]-1);
// Evaluate the secondary solution at all sampling points
if (retVal)
retVal = this->evalSolution(sField,integrand,gpar);
}
else
{
Go::SplineSurface* s = this->projectSolution(integrand);
if (s)
{
sField.resize(s->dimension(),s->numCoefs_u()*s->numCoefs_v());
sField.fill(&(*s->coefs_begin()));
delete s;
}
else
retVal = this->getGrevilleParameters(gpar[dir],dir);
retVal = false;
}
return retVal;
}
Go::GeomObject* ASMs2D::evalSolution (const Integrand& integrand) const
{
return this->projectSolution(integrand);
}
Go::SplineSurface* ASMs2D::projectSolution (const Integrand& integrand) const
{
// Compute parameter values of the result sampling points (Greville points)
DoubleVec gpar[2];
for (int dir = 0; dir < 2; dir++)
if (!this->getGrevilleParameters(gpar[dir],dir))
return 0;
// Evaluate the secondary solution at all sampling points
if (retVal) retVal = this->evalSolution(sField,integrand,gpar);
if (!retVal || npe) return retVal;
Matrix sValues;
if (!this->evalSolution(sValues,integrand,gpar))
return 0;
// Project the results onto the spline basis to find control point
// values based on the result values evaluated at the Greville points.
@ -1230,22 +1287,18 @@ bool ASMs2D::evalSolution (Matrix& sField, const Integrand& integrand,
// the result array. Think that is always the case, but beware if trying
// other projection schemes later.
Vector sValues = sField;
DoubleVec weights;
if (surf->rational())
surf->getWeights(weights);
Go::SplineSurface* s =
Go::SurfaceInterpolator::regularInterpolation(surf->basis(0),
surf->basis(1),
gpar[0], gpar[1],
sValues, sField.rows(),
surf->rational(), weights);
sField.fill(&(*s->coefs_begin()));
delete s;
return retVal;
const Vector& vec = sValues;
return Go::SurfaceInterpolator::regularInterpolation(surf->basis(0),
surf->basis(1),
gpar[0], gpar[1],
const_cast<Vector&>(vec),
sValues.rows(),
surf->rational(),
weights);
}

View File

@ -217,7 +217,13 @@ public:
//! points and then projected onto the spline basis to obtain the control
//! point values, which then are returned through \a sField.
virtual bool evalSolution(Matrix& sField, const Integrand& integrand,
const int* npe) const;
const int* npe = 0) const;
//! \brief Projects the secondary solution field onto the primary basis.
//! \param[in] integrand Object with problem-specific data and methods
Go::SplineSurface* projectSolution(const Integrand& integrand) const;
//! \brief Projects the secondary solution field onto the primary basis.
virtual Go::GeomObject* evalSolution(const Integrand& integrand) const;
protected:

View File

@ -1236,6 +1236,11 @@ bool ASMs3D::integrate (Integrand& integrand,
if (!utl::Hessian(Hess,d2NdX2,Jac,Xnod,d2Ndu2,dNdu))
return false;
#if SP_DEBUG > 4
std::cout <<"\niel, ip = "<< iel <<" "<< ip
<<"\nN ="<< N <<"dNdX ="<< dNdX << std::endl;
#endif
// Cartesian coordinates of current integration point
X = Xnod * N;
X.t = time.t;
@ -1734,19 +1739,56 @@ bool ASMs3D::evalSolution (Matrix& sField, const Vector& locSol,
bool ASMs3D::evalSolution (Matrix& sField, const Integrand& integrand,
const int* npe) const
{
bool retVal = true;
if (npe)
{
// Compute parameter values of the result sampling points
DoubleVec gpar[3];
for (int dir = 0; dir < 3 && retVal; dir++)
retVal = this->getGridParameters(gpar[dir],dir,npe[dir]-1);
// Evaluate the secondary solution at all sampling points
if (retVal)
retVal = this->evalSolution(sField,integrand,gpar);
}
else
{
Go::SplineVolume* v = this->projectSolution(integrand);
if (v)
{
sField.resize(v->dimension(),
v->numCoefs(0)*v->numCoefs(1)*v->numCoefs(2));
sField.fill(&(*v->coefs_begin()));
delete v;
}
else
retVal = false;
}
return retVal;
}
Go::GeomObject* ASMs3D::evalSolution (const Integrand& integrand) const
{
return this->projectSolution(integrand);
}
Go::SplineVolume* ASMs3D::projectSolution (const Integrand& integrand) const
{
// Compute parameter values of the result sampling points
DoubleVec gpar[3];
bool retVal = true;
for (int dir = 0; dir < 3 && retVal; dir++)
if (npe)
retVal = this->getGridParameters(gpar[dir],dir,npe[dir]-1);
else
retVal = this->getGrevilleParameters(gpar[dir],dir);
for (int dir = 0; dir < 3; dir++)
if (!this->getGrevilleParameters(gpar[dir],dir))
return 0;
// Evaluate the secondary solution at all sampling points
if (retVal) retVal = this->evalSolution(sField,integrand,gpar);
if (!retVal || npe) return retVal;
Matrix sValues;
if (!this->evalSolution(sValues,integrand,gpar))
return false;
// Project the results onto the spline basis to find control point
// values based on the result values evaluated at the Greville points.
@ -1755,23 +1797,19 @@ bool ASMs3D::evalSolution (Matrix& sField, const Integrand& integrand,
// the result array. Think that is always the case, but beware if trying
// other projection schemes later.
Vector sValues = sField;
DoubleVec weights;
if (svol->rational())
svol->getWeights(weights);
Go::SplineVolume* s =
Go::VolumeInterpolator::regularInterpolation(svol->basis(0),
svol->basis(1),
svol->basis(2),
gpar[0], gpar[1], gpar[2],
sValues, sField.rows(),
svol->rational(), weights);
sField.fill(&(*s->coefs_begin()));
delete s;
return retVal;
const Vector& vec = sValues;
return Go::VolumeInterpolator::regularInterpolation(svol->basis(0),
svol->basis(1),
svol->basis(2),
gpar[0], gpar[1], gpar[2],
const_cast<Vector&>(vec),
sValues.rows(),
svol->rational(),
weights);
}

View File

@ -282,7 +282,13 @@ public:
//! points and then projected onto the spline basis to obtain the control
//! point values, which then are returned through \a sField.
virtual bool evalSolution(Matrix& sField, const Integrand& integrand,
const int* npe) const;
const int* npe = 0) const;
//! \brief Projects the secondary solution field onto the primary basis.
//! \param[in] integrand Object with problem-specific data and methods
Go::SplineVolume* projectSolution(const Integrand& integrand) const;
//! \brief Projects the secondary solution field onto the primary basis.
virtual Go::GeomObject* evalSolution(const Integrand& integrand) const;
protected:

View File

@ -48,6 +48,10 @@ public:
//! \brief Resets the global element and node counters.
static void resetNumbering() { gEl = gNod = 0; }
//! \brief Projects the secondary solution field onto the primary basis.
//! \param[in] integrand Object with problem-specific data and methods
virtual Go::GeomObject* evalSolution(const Integrand& integrand) const = 0;
protected:
Go::GeomObject* geo; //!< Pointer to the actual spline geometry object